CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GEMRecHitValidation.cc
Go to the documentation of this file.
5 
7  : GEMBaseValidation(pset, "GEMRecHitsValidation") {
8  const auto& rechit_pset = pset.getParameterSet("gemRecHit");
9  const auto& rechit_tag = rechit_pset.getParameter<edm::InputTag>("inputTag");
10  rechit_token_ = consumes<GEMRecHitCollection>(rechit_tag);
11 
12  const auto& simhit_pset = pset.getParameterSet("gemSimHit");
13  const auto& simhit_tag = simhit_pset.getParameter<edm::InputTag>("inputTag");
14 
15  const auto& digisimlink_tag = pset.getParameter<edm::InputTag>("gemDigiSimLink");
16  digisimlink_token_ = consumes<edm::DetSetVector<GEMDigiSimLink>>(digisimlink_tag);
17 
18  simhit_token_ = consumes<edm::PSimHitContainer>(simhit_tag);
19  geomToken_ = esConsumes<GEMGeometry, MuonGeometryRecord>();
20  geomTokenBeginRun_ = esConsumes<GEMGeometry, MuonGeometryRecord, edm::Transition::BeginRun>();
21 }
22 
24 
26  const GEMGeometry* gem = &setup.getData(geomTokenBeginRun_);
27 
28  // NOTE Cluster Size
29  booker.setCurrentFolder("GEM/RecHits");
30 
31  TString cls_title = "Cluster Size Distribution";
32  TString cls_x_title = "Cluster size";
33 
34  if (detail_plot_) {
35  me_detail_cls_total_ = booker.book1D("cls", cls_title + ";" + cls_x_title + ";" + "Entries", 10, 0.5, 10.5);
36 
37  for (const auto& station : gem->regions()[0]->stations()) {
38  Int_t station_id = station->station();
39  for (const auto& roll : station->superChambers()[0]->chambers()[0]->etaPartitions()) {
40  Int_t ieta = roll->id().ieta();
41  ME2IdsKey key{station_id, ieta};
42  me_detail_cls_roll_[key] = booker.book1D(Form("cls_GE%d1-E%d", station_id, ieta),
43  Form("Cluster Size Distribution : GE%d1-E%d", station_id, ieta),
44  10,
45  0.5,
46  10.5);
47  }
48  }
49 
50  for (const auto& region : gem->regions()) {
51  Int_t region_id = region->region();
52 
53  for (const auto& station : region->stations()) {
54  Int_t station_id = station->station();
55 
56  const auto& superChamberVec = station->superChambers();
57  if (superChamberVec.empty() || superChamberVec[0] == nullptr) {
58  edm::LogError(kLogCategory_) << "Super chambers missing or null for region = " << region_id
59  << " and station = " << station_id;
60  } else {
61  for (const auto& chamber : superChamberVec[0]->chambers()) {
62  Int_t layer_id = chamber->id().layer();
63 
64  for (const auto& roll : chamber->etaPartitions()) {
65  Int_t ieta = roll->id().ieta();
66  ME4IdsKey key4{region_id, station_id, layer_id, ieta};
67 
68  me_detail_cls_[key4] = bookHist1D(booker, key4, "cls", "Cluster Size Distribution", 11, -0.5, 10.5);
69  } // roll loop
70  } // chamber loop
71  } // end else
72  } // station loop
73  } // region loop
74  } // detail plot
75 
76  // NOTE Residual
77  for (const auto& station : gem->regions()[0]->stations()) {
78  Int_t station_id = station->station();
79  for (const auto& roll : station->superChambers()[0]->chambers()[0]->etaPartitions()) {
80  Int_t ieta = roll->id().ieta();
81  ME2IdsKey key{station_id, ieta};
82 
83  me_residual_y_[key] = booker.book1D(Form("residual_y_GE%d1-E%d", station_id, ieta),
84  Form("Residual in Y : GE%d1-E%d; Residual in Y [cm]", station_id, ieta),
85  60,
86  -15,
87  15);
88 
90  booker.book1D(Form("residual_rphi_GE%d1-E%d", station_id, ieta),
91  Form("Residual in R #times #phi : GE%d1-E%d; Residual in r #times #phi [cm]", station_id, ieta),
92  60,
93  -15,
94  15);
95  }
96  }
97 
98  if (detail_plot_) {
99  for (const auto& region : gem->regions()) {
100  Int_t region_id = region->region();
101 
102  for (const auto& station : region->stations()) {
103  Int_t station_id = station->station();
104 
105  const auto& superChamberVec = station->superChambers();
106  if (!superChamberVec.empty() && superChamberVec[0] != nullptr) {
107  for (const auto& chamber : superChamberVec[0]->chambers()) {
108  Int_t layer_id = chamber->id().layer();
109 
110  for (const auto& roll : chamber->etaPartitions()) {
111  Int_t ieta = roll->id().ieta();
112  ME4IdsKey key4{region_id, station_id, layer_id, ieta};
113 
114  me_detail_residual_y_[key4] =
115  bookHist1D(booker, key4, "residual_y", "Residual in y", 60, -15, 15, "Residual in y [cm]");
116 
117  me_detail_residual_rphi_[key4] = bookHist1D(booker,
118  key4,
119  "residual_rphi",
120  "Residual in r #times #phi",
121  60,
122  -15,
123  15,
124  "Residual in r #times #phi [cm]");
125  } // roll loop
126  } // chamber loop
127  } // end if
128  } // station loop
129  } // region loop
130  } // detail plot
131 
132  // NOTE Pull
133  if (detail_plot_) {
134  for (const auto& station : gem->regions()[0]->stations()) {
135  Int_t station_id = station->station();
136  for (const auto& roll : station->superChambers()[0]->chambers()[0]->etaPartitions()) {
137  Int_t ieta = roll->id().ieta();
138  ME2IdsKey key{station_id, ieta};
139 
140  me_detail_pull_x_[key] = booker.book1D(
141  Form("pull_x_GE%d1-E%d", station_id, ieta), Form("Pull in X : GE%d1-E%d", station_id, ieta), 60, -3, 3);
142 
143  me_detail_pull_y_[key] = booker.book1D(
144  Form("pull_y_GE%d1-E%d", station_id, ieta), Form("Pull in Y : GE%d1-E%d", station_id, ieta), 60, -3, 3);
145  }
146  }
147 
148  for (const auto& region : gem->regions()) {
149  Int_t region_id = region->region();
150 
151  for (const auto& station : region->stations()) {
152  Int_t station_id = station->station();
153 
154  const auto& superChamberVec = station->superChambers();
155  if (!superChamberVec.empty() && superChamberVec[0] != nullptr) {
156  for (const auto& chamber : superChamberVec[0]->chambers()) {
157  Int_t layer_id = chamber->id().layer();
158 
159  for (const auto& roll : chamber->etaPartitions()) {
160  Int_t ieta = roll->id().ieta();
161  ME4IdsKey key4{region_id, station_id, layer_id, ieta};
162 
163  me_detail_pull_x_la_[key4] = bookHist1D(booker, key4, "pull_x", "Pull in x", 60, -3, 3);
164 
165  me_detail_pull_y_la_[key4] = bookHist1D(booker, key4, "pull_y", "Pull in y", 60, -3, 3);
166  } // roll loop
167  } // chamber loop
168  } // end if
169  } // station loop
170  } // region loop
171  } // detail plot
172 
173  // NOTE Occupancy
174  for (const auto& region : gem->regions()) {
175  Int_t region_id = region->region();
176 
177  if (detail_plot_)
178  me_detail_occ_zr_[region_id] = bookZROccupancy(booker, region_id, "rechit", "RecHit");
179 
180  for (const auto& station : region->stations()) {
181  Int_t station_id = station->station();
182  ME2IdsKey key2{region_id, station_id};
183 
184  if (detail_plot_)
185  me_detail_rechit_occ_det_[key2] = bookDetectorOccupancy(booker, key2, station, "sim_matched", "Matched RecHit");
186 
187  const auto& superChamberVec = station->superChambers();
188  if (!superChamberVec.empty() && superChamberVec[0] != nullptr) {
189  for (const auto& chamber : superChamberVec[0]->chambers()) {
190  Int_t layer_id = chamber->id().layer();
191  ME3IdsKey key3{region_id, station_id, layer_id};
192 
193  Int_t num_eta_partitions = chamber->nEtaPartitions();
194 
195  me_rechit_occ_eta_[key3] = bookHist1D(booker,
196  key3,
197  "sim_matched_occ_eta",
198  "Matched RecHit Eta Occupancy",
199  16,
200  eta_range_[station_id * 2 + 0],
201  eta_range_[station_id * 2 + 1],
202  "|#eta|");
203 
204  me_rechit_occ_phi_[key3] =
205  bookHist1D(booker, key3, "sim_matched_occ_phi", "Matched RecHit Phi Occupancy", 36, -5, 355, "#phi");
206 
207  if (detail_plot_) {
209  bookHist1D(booker, key3, "total_rechit", "Number of rec hits per event", 25, -0.5, 24.5);
210 
211  me_detail_occ_pid_[key3] = bookPIDHist(booker, key3, "sim_occ_pid", "Number of entreis for each particle");
212 
213  me_detail_occ_ieta_[key3] = bookHist1D(booker,
214  key3,
215  "occ_ieta",
216  "Rechit Occupancy per eta partition",
217  num_eta_partitions,
218  0.5,
219  num_eta_partitions + 0.5);
220 
221  me_detail_occ_phi_[key3] = bookHist1D(booker, key3, "occ_phi", "Rechit Phi Occupancy", 108, -5, 355);
222 
223  me_detail_occ_xy_[key3] = bookXYOccupancy(booker, key3, "rechit", "RecHit");
224 
225  me_detail_occ_polar_[key3] = bookPolarOccupancy(booker, key3, "rechit", "RecHit");
226  }
227  } // chamber loop
228  } // end if
229  } // station loop
230  } // region_loop
231 }
232 
234  Int_t cls = rechit->clusterSize();
235  Int_t rechit_first_strip = rechit->firstClusterStrip();
236 
237  if (cls == 1) {
238  return simhit_strip == rechit_first_strip;
239  } else {
240  Int_t rechit_last_strip = rechit_first_strip + cls - 1;
241  return (simhit_strip >= rechit_first_strip) and (simhit_strip <= rechit_last_strip);
242  }
243 }
244 
246  const GEMGeometry* gem = &setup.getData(geomToken_);
247 
249  event.getByToken(digisimlink_token_, digiSimLink);
250  if (not digiSimLink.isValid()) {
251  edm::LogError(kLogCategory_) << "Failed to get GEMDigiSimLink." << std::endl;
252  return;
253  }
254 
255  edm::Handle<edm::PSimHitContainer> simhit_container;
256  event.getByToken(simhit_token_, simhit_container);
257  if (not simhit_container.isValid()) {
258  edm::LogError(kLogCategory_) << "Failed to get PSimHitContainer." << std::endl;
259  return;
260  }
261 
262  edm::Handle<GEMRecHitCollection> rechit_collection;
263  event.getByToken(rechit_token_, rechit_collection);
264  if (not rechit_collection.isValid()) {
265  edm::LogError(kLogCategory_) << "Failed to get GEMRecHitCollection" << std::endl;
266  return;
267  }
268 
269  std::map<ME3IdsKey, Int_t> total_rechit;
270  for (const auto& rechit : *rechit_collection) {
271  GEMDetId gem_id{rechit.gemId()};
272  Int_t region_id = gem_id.region();
273  Int_t station_id = gem_id.station();
274  Int_t layer_id = gem_id.layer();
275  Int_t ieta = gem_id.ieta();
276  ME4IdsKey key4{region_id, station_id, layer_id, ieta};
277 
278  ME2IdsKey key{station_id, ieta};
279  ME2IdsKey key2{region_id, station_id};
280  ME3IdsKey key3{region_id, station_id, layer_id};
281 
282  const BoundPlane& surface = gem->idToDet(gem_id)->surface();
283  GlobalPoint&& rechit_global_pos = surface.toGlobal(rechit.localPosition());
284 
285  Float_t rechit_g_x = rechit_global_pos.x();
286  Float_t rechit_g_y = rechit_global_pos.y();
287  Float_t rechit_g_abs_z = std::fabs(rechit_global_pos.z());
288  Float_t rechit_g_r = rechit_global_pos.perp();
289  Float_t rechit_g_phi = toDegree(rechit_global_pos.phi());
290 
291  Int_t first_strip = rechit.firstClusterStrip();
292  Int_t cls = rechit.clusterSize();
293  cls = cls > 10 ? 10 : cls;
294 
295  total_rechit[key3]++;
296 
297  if (detail_plot_) {
299  me_detail_cls_[key4]->Fill(cls);
300  me_detail_cls_roll_[key]->Fill(cls);
301  me_detail_occ_ieta_[key3]->Fill(ieta);
302  me_detail_occ_phi_[key3]->Fill(rechit_g_phi);
303  me_detail_occ_zr_[region_id]->Fill(rechit_g_abs_z, rechit_g_r);
304  me_detail_occ_xy_[key3]->Fill(rechit_g_x, rechit_g_y);
305  me_detail_occ_polar_[key3]->Fill(rechit_g_phi, rechit_g_r);
306  } // detail plot
307 
308  auto links = digiSimLink->find(gem_id);
309 
310  if (links == digiSimLink->end())
311  continue;
312  std::map<Int_t, Int_t> pid_count;
313 
314  for (Int_t strip = first_strip; strip < first_strip + cls; strip++) {
315  for (const auto& link : *links) {
316  Int_t link_strip = link.getStrip();
317  if (link_strip == strip) {
318  Int_t pid = link.getParticleType();
319  pid_count[pid]++;
320  break;
321  }
322  }
323  }
324  if (detail_plot_) {
325  Int_t max_pid = 0;
326  Int_t max_count = 0;
327  for (auto& [pid, count] : pid_count) {
328  if (max_count < count) {
329  max_pid = pid;
330  max_count = count;
331  }
332  }
333  Int_t pid_idx = getPidIdx(max_pid);
334  me_detail_occ_pid_[key3]->Fill(pid_idx);
335  }
336  }
337 
338  if (detail_plot_) {
339  for (auto [key, num_total_rechit] : total_rechit) {
340  me_detail_total_rechit_[key]->Fill(num_total_rechit);
341  }
342  }
343 
344  // NOTE
345  for (const auto& simhit : *simhit_container.product()) {
346  if (gem->idToDet(simhit.detUnitId()) == nullptr) {
347  edm::LogError(kLogCategory_) << "MuonGEMHit didn't matched with GEMGeometry." << std::endl;
348  continue;
349  }
350 
351  GEMDetId simhit_gemid{simhit.detUnitId()};
352  const BoundPlane& surface = gem->idToDet(simhit_gemid)->surface();
353 
354  Int_t region_id = simhit_gemid.region();
355  Int_t station_id = simhit_gemid.station();
356  Int_t layer_id = simhit_gemid.layer();
357  Int_t chamber_id = simhit_gemid.chamber();
358  Int_t ieta = simhit_gemid.ieta();
359  Int_t num_layers = simhit_gemid.nlayers();
360 
361  ME2IdsKey key{station_id, ieta};
362  ME2IdsKey key2{region_id, station_id};
363  ME3IdsKey key3{region_id, station_id, layer_id};
364  ME4IdsKey key4{region_id, station_id, layer_id, ieta};
365 
366  const LocalPoint& simhit_local_pos = simhit.localPosition();
367  const GlobalPoint& simhit_global_pos = surface.toGlobal(simhit_local_pos);
368 
369  Float_t simhit_g_abs_eta = std::fabs(simhit_global_pos.eta());
370  Float_t simhit_g_phi = toDegree(simhit_global_pos.phi());
371 
372  Int_t det_occ_bin_x = getDetOccBinX(num_layers, chamber_id, layer_id);
373 
374  auto simhit_trackId = simhit.trackId();
375 
376  auto links = digiSimLink->find(simhit_gemid);
377  if (links == digiSimLink->end())
378  continue;
379 
380  Int_t simhit_strip = -1;
381  for (const auto& link : *links) {
382  if (simhit_trackId == link.getTrackId()) {
383  simhit_strip = link.getStrip();
384  break;
385  }
386  }
387 
388  GEMRecHitCollection::range range = rechit_collection->get(simhit_gemid);
389  for (auto rechit = range.first; rechit != range.second; ++rechit) {
390  if (gem->idToDet(rechit->gemId()) == nullptr) {
391  edm::LogError(kLogCategory_) << "GEMRecHit didn't matched with GEMGeometry." << std::endl;
392  continue;
393  }
394 
395  if (not isMuonSimHit(simhit))
396  continue;
397 
398  if (matchRecHitAgainstSimHit(rechit, simhit_strip)) {
399  const LocalPoint& rechit_local_pos = rechit->localPosition();
400 
401  Float_t resolution_x = std::sqrt(rechit->localPositionError().xx());
402  Float_t resolution_y = std::sqrt(rechit->localPositionError().yy());
403 
404  Float_t residual_x = rechit_local_pos.x() - simhit_local_pos.x();
405  Float_t residual_y = rechit_local_pos.y() - simhit_local_pos.y();
406  Float_t residual_r = sqrt(pow(residual_x, 2) + pow(residual_y, 2));
407  Float_t residual_phi = rechit_local_pos.phi() - simhit_local_pos.phi();
408  Float_t residual_rphi = residual_r * residual_phi;
409 
410  Float_t pull_x = residual_x / resolution_x;
411  Float_t pull_y = residual_y / resolution_y;
412 
413  me_residual_y_[key]->Fill(residual_y);
414  me_residual_rphi_[key]->Fill(residual_rphi);
415 
416  me_rechit_occ_eta_[key3]->Fill(simhit_g_abs_eta);
417  me_rechit_occ_phi_[key3]->Fill(simhit_g_phi);
418 
419  if (detail_plot_) {
420  me_detail_rechit_occ_det_[key2]->Fill(det_occ_bin_x, ieta);
421 
422  me_detail_residual_y_[key4]->Fill(residual_y);
423  me_detail_residual_rphi_[key4]->Fill(residual_rphi);
424 
425  me_detail_pull_x_[key]->Fill(pull_x);
426  me_detail_pull_y_[key]->Fill(pull_y);
427  me_detail_pull_x_la_[key4]->Fill(pull_x);
428  me_detail_pull_y_la_[key4]->Fill(pull_y);
429  } // detail_plot
430 
431  // If we find GEMRecHit that matches PSimHit, then exit
432  // GEMRecHitCollection loop.
433  break;
434 
435  } // if rechit matches against simhit
436  } // rechit loop
437  } // simhit loop
438 }
dqm::impl::MonitorElement * bookHist1D(DQMStore::IBooker &booker, const T &key, const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, const char *x_title="", const char *y_title="Entries")
Bool_t matchRecHitAgainstSimHit(GEMRecHitCollection::const_iterator, Int_t)
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
MEMap2Ids me_detail_rechit_occ_det_
T perp() const
Definition: PV3DBase.h:69
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
Bool_t isMuonSimHit(const PSimHit &)
std::tuple< Int_t, Int_t, Int_t, Int_t > ME4IdsKey
const GeomDet * idToDet(DetId) const override
Definition: GEMGeometry.cc:25
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T y() const
Definition: PV3DBase.h:60
MonitorElement * me_detail_cls_total_
Log< level::Error, false > LogError
dqm::impl::MonitorElement * bookZROccupancy(DQMStore::IBooker &booker, Int_t region_id, const char *name_prfix, const char *title_prefix)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const uint16_t range(const Frame &aFrame)
void Fill(long long x)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
edm::ESGetToken< GEMGeometry, MuonGeometryRecord > geomTokenBeginRun_
dqm::impl::MonitorElement * bookXYOccupancy(DQMStore::IBooker &booker, const T &key, const char *name_prefix, const char *title_prefix)
T sqrt(T t)
Definition: SSEVec.h:19
constexpr int region() const
Definition: GEMDetId.h:171
const std::vector< const GEMRegion * > & regions() const
Return a vector of all GEM regions.
Definition: GEMGeometry.cc:30
T z() const
Definition: PV3DBase.h:61
Int_t getPidIdx(Int_t pid)
tuple key
prepare the HTCondor submission files and eventually submit them
std::vector< Double_t > eta_range_
edm::EDGetTokenT< edm::PSimHitContainer > simhit_token_
void analyze(const edm::Event &, const edm::EventSetup &) override
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< edm::DetSetVector< GEMDigiSimLink > > digisimlink_token_
T const * product() const
Definition: Handle.h:70
ParameterSet const & getParameterSet(std::string const &) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< GEMRecHitCollection > rechit_token_
MEMap4Ids me_detail_residual_rphi_
T eta() const
Definition: PV3DBase.h:73
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
GEMRecHitValidation(const edm::ParameterSet &)
std::tuple< Int_t, Int_t > ME2IdsKey
Float_t toDegree(Float_t radian)
std::tuple< Int_t, Int_t, Int_t > ME3IdsKey
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
static char chambers[264][20]
Definition: ReadPGInfo.cc:243
Int_t getDetOccBinX(Int_t num_layers, Int_t chamber_id, Int_t layer_id)
T x() const
Definition: PV3DBase.h:59
dqm::impl::MonitorElement * bookPIDHist(DQMStore::IBooker &booker, const T &key, const char *name, const char *title)
const std::string kLogCategory_
dqm::impl::MonitorElement * bookPolarOccupancy(DQMStore::IBooker &booker, const T &key, const char *name_prefix, const char *title_prefix)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Definition: Run.h:45
edm::ESGetToken< GEMGeometry, MuonGeometryRecord > geomToken_
dqm::impl::MonitorElement * bookDetectorOccupancy(DQMStore::IBooker &booker, const T &key, const GEMStation *station, const char *name_prfix, const char *title_prefix)