CMS 3D CMS Logo

TotemRPDQMSource.cc
Go to the documentation of this file.
1 /****************************************************************************
2 *
3 * This is a part of TotemDQM and TOTEM offline software.
4 * Authors:
5 * Jan Kašpar (jan.kaspar@gmail.com)
6 * Rafał Leszko (rafal.leszko@gmail.com)
7 *
8 ****************************************************************************/
9 
16 
20 
29 
33 
34 #include <string>
35 
36 //----------------------------------------------------------------------------------------------------
37 
39 {
40  public:
42  ~TotemRPDQMSource() override;
43 
44  protected:
45  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
46  void analyze(edm::Event const& e, edm::EventSetup const& eSetup) override;
47 
48  private:
49  unsigned int verbosity;
50 
57 
59  struct PotPlots
60  {
62 
63  MonitorElement *activity=nullptr, *activity_u=nullptr, *activity_v=nullptr;
66  MonitorElement *patterns_u=nullptr, *patterns_v=nullptr;
72 
73  PotPlots() {}
74  PotPlots(DQMStore::IBooker &ibooker, unsigned int id);
75  };
76 
77  std::map<unsigned int, PotPlots> potPlots;
78 
80  struct PlanePlots
81  {
82  MonitorElement *digi_profile_cumulative = nullptr;
83  MonitorElement *cluster_profile_cumulative = nullptr;
84  MonitorElement *hit_multiplicity = nullptr;
85  MonitorElement *cluster_size = nullptr;
86  MonitorElement *efficiency_num = nullptr, *efficiency_den = nullptr;
87 
89  PlanePlots(DQMStore::IBooker &ibooker, unsigned int id);
90  };
91 
92  std::map<unsigned int, PlanePlots> planePlots;
93 };
94 
95 //----------------------------------------------------------------------------------------------------
96 //----------------------------------------------------------------------------------------------------
97 
98 using namespace std;
99 using namespace edm;
100 
101 //----------------------------------------------------------------------------------------------------
102 
104 {
105  string path;
107  ibooker.setCurrentFolder(path);
108 
109  string title;
111 
112  vfat_problem = ibooker.book2D("vfats with any problem", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
113  vfat_missing = ibooker.book2D("vfats missing", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
114  vfat_ec_bc_error = ibooker.book2D("vfats with EC or BC error", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
115  vfat_corruption = ibooker.book2D("vfats with data corruption", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
116 
117  activity = ibooker.book1D("active planes", title+";number of active planes", 11, -0.5, 10.5);
118  activity_u = ibooker.book1D("active planes U", title+";number of active U planes", 11, -0.5, 10.5);
119  activity_v = ibooker.book1D("active planes V", title+";number of active V planes", 11, -0.5, 10.5);
120 
121  activity_per_bx = ibooker.book1D("activity per BX", title+";Event.BX", 4002, -1.5, 4000. + 0.5);
122  activity_per_bx_short = ibooker.book1D("activity per BX (short)", title+";Event.BX", 102, -1.5, 100. + 0.5);
123 
124  hit_plane_hist = ibooker.book2D("activity in planes (2D)", title+";plane number;strip number", 10, -0.5, 9.5, 32, -0.5, 511.5);
125 
126  patterns_u = ibooker.book1D("recognized patterns U", title+";number of recognized U patterns", 11, -0.5, 10.5);
127  patterns_v = ibooker.book1D("recognized patterns V", title+";number of recognized V patterns", 11, -0.5, 10.5);
128 
129  h_planes_fit_u = ibooker.book1D("planes contributing to fit U", title+";number of planes contributing to U fit", 6, -0.5, 5.5);
130  h_planes_fit_v = ibooker.book1D("planes contributing to fit V", title+";number of planes contributing to V fit", 6, -0.5, 5.5);
131 
132  event_category = ibooker.book1D("event category", title+";event category", 5, -0.5, 4.5);
133  TH1F *event_category_h = event_category->getTH1F();
134  event_category_h->GetXaxis()->SetBinLabel(1, "empty");
135  event_category_h->GetXaxis()->SetBinLabel(2, "insufficient");
136  event_category_h->GetXaxis()->SetBinLabel(3, "single-track");
137  event_category_h->GetXaxis()->SetBinLabel(4, "multi-track");
138  event_category_h->GetXaxis()->SetBinLabel(5, "shower");
139 
140  trackHitsCumulativeHist = ibooker.book2D("track XY profile", title+";x (mm);y (mm)", 100, -18., +18., 100, -18., +18.);
141 
142  track_u_profile = ibooker.book1D("track profile U", title+"; U (mm)", 512, -256*66E-3, +256*66E-3);
143  track_v_profile = ibooker.book1D("track profile V", title+"; V (mm)", 512, -256*66E-3, +256*66E-3);
144 
145  triggerSectorUVCorrelation_all = ibooker.book2D("trigger sector UV correlation (no cond)", title+";V sector;U sector", 16, -0.5, 15.5, 16, -0.5, 15.5);
146  triggerSectorUVCorrelation_mult2 = ibooker.book2D("trigger sector UV correlation (max mult 2)", title+";V sector;U sector", 16, -0.5, 15.5, 16, -0.5, 15.5);
147  triggerSectorUVCorrelation_track = ibooker.book2D("trigger sector UV correlation (track)", title+";V sector;U sector", 16, -0.5, 15.5, 16, -0.5, 15.5);
148 }
149 
150 //----------------------------------------------------------------------------------------------------
151 
153 {
154  string path;
156  ibooker.setCurrentFolder(path);
157 
158  string title;
160 
161  digi_profile_cumulative = ibooker.book1D("digi profile", title+";strip number", 512, -0.5, 511.5);
162  cluster_profile_cumulative = ibooker.book1D("cluster profile", title+";cluster center", 1024, -0.25, 511.75);
163  hit_multiplicity = ibooker.book1D("hit multiplicity", title+";hits/detector/event", 6, -0.5, 5.5);
164  cluster_size = ibooker.book1D("cluster size", title+";hits per cluster", 5, 0.5, 5.5);
165 
166  efficiency_num = ibooker.book1D("efficiency num", title+";track position (mm)", 30, -15., 0.);
167  efficiency_den = ibooker.book1D("efficiency den", title+";track position (mm)", 30, -15., 0.);
168 }
169 
170 //----------------------------------------------------------------------------------------------------
171 //----------------------------------------------------------------------------------------------------
172 
174  verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0))
175 {
176  tokenStatus = consumes<DetSetVector<TotemVFATStatus>>(ps.getParameter<edm::InputTag>("tagStatus"));
177 
178  tokenDigi = consumes< DetSetVector<TotemRPDigi> >(ps.getParameter<edm::InputTag>("tagDigi"));
179  tokenCluster = consumes< edm::DetSetVector<TotemRPCluster> >(ps.getParameter<edm::InputTag>("tagCluster"));
180  tokenRecHit = consumes< edm::DetSetVector<TotemRPRecHit> >(ps.getParameter<edm::InputTag>("tagRecHit"));
181  tokenUVPattern = consumes< DetSetVector<TotemRPUVPattern> >(ps.getParameter<edm::InputTag>("tagUVPattern"));
182  tokenLocalTrack = consumes< DetSetVector<TotemRPLocalTrack> >(ps.getParameter<edm::InputTag>("tagLocalTrack"));
183 }
184 
185 //----------------------------------------------------------------------------------------------------
186 
188 {
189 }
190 
191 //----------------------------------------------------------------------------------------------------
192 
194 {
195  ibooker.cd();
196  ibooker.setCurrentFolder("CTPPS");
197 
198  // loop over arms
199  for (unsigned int arm : {0, 1})
200  {
201  // loop over RPs
202  for (unsigned int st_rp : {2, 3, 4, 5, 24, 25})
203  {
204  const unsigned int st = st_rp / 10;
205  const unsigned int rp = st_rp % 10;
206 
207  TotemRPDetId rpId(arm, st, rp);
208  potPlots[rpId] = PotPlots(ibooker, rpId);
209 
210  // loop over planes
211  for (unsigned int pl = 0; pl < 10; ++pl)
212  {
213  TotemRPDetId plId(arm, st, rp, pl);
214  planePlots[plId] = PlanePlots(ibooker, plId);
215  }
216  }
217  }
218 }
219 
220 //----------------------------------------------------------------------------------------------------
221 
223 {
224  // get event setup data
226  eventSetup.get<VeryForwardRealGeometryRecord>().get(geometry);
227 
228  // get event data
230  event.getByToken(tokenStatus, status);
231 
233  event.getByToken(tokenDigi, digi);
234 
236  event.getByToken(tokenCluster, digCluster);
237 
239  event.getByToken(tokenRecHit, hits);
240 
242  event.getByToken(tokenUVPattern, patterns);
243 
245  event.getByToken(tokenLocalTrack, tracks);
246 
247  // check validity
248  bool valid = true;
249  valid &= status.isValid();
250  valid &= digi.isValid();
251  valid &= digCluster.isValid();
252  valid &= hits.isValid();
253  valid &= patterns.isValid();
254  valid &= tracks.isValid();
255 
256  if (!valid)
257  {
258  if (verbosity)
259  {
260  LogProblem("TotemRPDQMSource") <<
261  "ERROR in TotemDQMModuleRP::analyze > some of the required inputs are not valid. Skipping this event.\n"
262  << " status.isValid = " << status.isValid() << "\n"
263  << " digi.isValid = " << digi.isValid() << "\n"
264  << " digCluster.isValid = " << digCluster.isValid() << "\n"
265  << " hits.isValid = " << hits.isValid() << "\n"
266  << " patterns.isValid = " << patterns.isValid() << "\n"
267  << " tracks.isValid = " << tracks.isValid();
268  }
269 
270  return;
271  }
272 
273  //------------------------------
274  // Status Plots
275 
276  for (auto &ds : *status)
277  {
278  TotemRPDetId detId(ds.detId());
279  unsigned int plNum = detId.plane();
280  CTPPSDetId rpId = detId.getRPId();
281 
282  auto it = potPlots.find(rpId);
283  if (it == potPlots.end())
284  continue;
285  auto &plots = it->second;
286 
287  for (auto &s : ds)
288  {
289  if (s.isMissing())
290  {
291  plots.vfat_problem->Fill(plNum, s.getChipPosition());
292  plots.vfat_missing->Fill(plNum, s.getChipPosition());
293  }
294 
295  if (s.isECProgressError() || s.isBCProgressError())
296  {
297  plots.vfat_problem->Fill(plNum, s.getChipPosition());
298  plots.vfat_ec_bc_error->Fill(plNum, s.getChipPosition());
299  }
300 
301  if (s.isIDMismatch() || s.isFootprintError() || s.isCRCError())
302  {
303  plots.vfat_problem->Fill(plNum, s.getChipPosition());
304  plots.vfat_corruption->Fill(plNum, s.getChipPosition());
305  }
306  }
307  }
308 
309  //------------------------------
310  // Plane Plots
311 
312  // digi profile cumulative
313  for (DetSetVector<TotemRPDigi>::const_iterator it = digi->begin(); it != digi->end(); ++it)
314  {
315  TotemRPDetId detId(it->detId());
316 
317  auto plIt = planePlots.find(detId);
318  if (plIt == planePlots.end())
319  continue;
320  auto &plots = plIt->second;
321 
322  for (DetSet<TotemRPDigi>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
323  plots.digi_profile_cumulative->Fill(dit->getStripNumber());
324  }
325 
326  // cluster plots
327  for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++)
328  {
329  TotemRPDetId detId(it->detId());
330 
331  auto plIt = planePlots.find(detId);
332  if (plIt == planePlots.end())
333  continue;
334  auto &plots = plIt->second;
335 
336  // hit multiplicity
337  plots.hit_multiplicity->Fill(it->size());
338 
339  for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
340  {
341  // profile cumulative
342  plots.cluster_profile_cumulative->Fill(dit->getCenterStripPosition());
343 
344  // cluster size
345  plots.cluster_size->Fill(dit->getNumberOfStrips());
346  }
347  }
348 
349  // plane efficiency plots
350  for (auto &ds : *tracks)
351  {
352  CTPPSDetId rpId(ds.detId());
353 
354  for (auto &ft : ds)
355  {
356  if (!ft.isValid())
357  continue;
358 
359  double rp_z = geometry->getRPTranslation(rpId).z();
360 
361  for (unsigned int plNum = 0; plNum < 10; ++plNum)
362  {
363  TotemRPDetId plId = rpId;
364  plId.setPlane(plNum);
365 
366  auto plIt = planePlots.find(plId);
367  if (plIt == planePlots.end())
368  continue;
369  auto &plots = plIt->second;
370 
371  double ft_z = ft.getZ0();
372  double ft_x = ft.getX0() + ft.getTx() * (ft_z - rp_z);
373  double ft_y = ft.getY0() + ft.getTy() * (ft_z - rp_z);
374 
375  double ft_v = geometry->globalToLocal(plId, CLHEP::Hep3Vector(ft_x, ft_y, ft_z)).y();
376 
377  bool hasMatchingHit = false;
378  const auto &hit_ds_it = hits->find(plId);
379  if (hit_ds_it != hits->end())
380  {
381  for (const auto &h : *hit_ds_it)
382  {
383  bool match = (fabs(ft_v - h.getPosition()) < 2.*0.066);
384  if (match)
385  {
386  hasMatchingHit = true;
387  break;
388  }
389  }
390  }
391 
392  plots.efficiency_den->Fill(ft_v);
393  if (hasMatchingHit)
394  plots.efficiency_num->Fill(ft_v);
395  }
396  }
397  }
398 
399 
400  //------------------------------
401  // Roman Pots Plots
402 
403  // determine active planes (from RecHits and VFATStatus)
404  map<unsigned int, set<unsigned int> > planes;
405  map<unsigned int, set<unsigned int> > planes_u;
406  map<unsigned int, set<unsigned int> > planes_v;
407  for (const auto &ds : *hits)
408  {
409  if (ds.empty())
410  continue;
411 
412  TotemRPDetId detId(ds.detId());
413  unsigned int planeNum = detId.plane();
414  CTPPSDetId rpId = detId.getRPId();
415 
416  planes[rpId].insert(planeNum);
417  if (detId.isStripsCoordinateUDirection())
418  planes_u[rpId].insert(planeNum);
419  else
420  planes_v[rpId].insert(planeNum);
421  }
422 
423  for (const auto &ds : *status)
424  {
425  bool activity = false;
426  for (const auto &s : ds)
427  {
428  if (s.isNumberOfClustersSpecified() && s.getNumberOfClusters() > 0)
429  {
430  activity = true;
431  break;
432  }
433  }
434 
435  if (!activity)
436  continue;
437 
438  TotemRPDetId detId(ds.detId());
439  unsigned int planeNum = detId.plane();
440  CTPPSDetId rpId = detId.getRPId();
441 
442  planes[rpId].insert(planeNum);
443  if (detId.isStripsCoordinateUDirection())
444  planes_u[rpId].insert(planeNum);
445  else
446  planes_v[rpId].insert(planeNum);
447  }
448 
449  // plane activity histogram
450  for (std::map<unsigned int, PotPlots>::iterator it = potPlots.begin(); it != potPlots.end(); it++)
451  {
452  it->second.activity->Fill(planes[it->first].size());
453  it->second.activity_u->Fill(planes_u[it->first].size());
454  it->second.activity_v->Fill(planes_v[it->first].size());
455 
456  if (planes[it->first].size() >= 6)
457  {
458  it->second.activity_per_bx->Fill(event.bunchCrossing());
459  it->second.activity_per_bx_short->Fill(event.bunchCrossing());
460  }
461  }
462 
463  for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++)
464  {
465  TotemRPDetId detId(it->detId());
466  unsigned int planeNum = detId.plane();
467  CTPPSDetId rpId = detId.getRPId();
468 
469  auto plIt = potPlots.find(rpId);
470  if (plIt == potPlots.end())
471  continue;
472  auto &plots = plIt->second;
473 
474  for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
475  plots.hit_plane_hist->Fill(planeNum, dit->getCenterStripPosition());
476  }
477 
478  // recognized pattern histograms
479  for (auto &ds : *patterns)
480  {
481  CTPPSDetId rpId(ds.detId());
482 
483  auto plIt = potPlots.find(rpId);
484  if (plIt == potPlots.end())
485  continue;
486  auto &plots = plIt->second;
487 
488  // count U and V patterns
489  unsigned int u = 0, v = 0;
490  for (auto &p : ds)
491  {
492  if (! p.getFittable())
493  continue;
494 
495  if (p.getProjection() == TotemRPUVPattern::projU)
496  u++;
497 
498  if (p.getProjection() == TotemRPUVPattern::projV)
499  v++;
500  }
501 
502  plots.patterns_u->Fill(u);
503  plots.patterns_v->Fill(v);
504  }
505 
506  // event-category histogram
507  for (auto &it : potPlots)
508  {
509  TotemRPDetId rpId(it.first);
510  auto &pp = it.second;
511 
512  // process hit data for this plot
513  unsigned int pl_u = planes_u[rpId].size();
514  unsigned int pl_v = planes_v[rpId].size();
515 
516  // process pattern data for this pot
517  const auto &rp_pat_it = patterns->find(rpId);
518 
519  unsigned int pat_u = 0, pat_v = 0;
520  if (rp_pat_it != patterns->end())
521  {
522  for (auto &p : *rp_pat_it)
523  {
524  if (! p.getFittable())
525  continue;
526 
527  if (p.getProjection() == TotemRPUVPattern::projU)
528  pat_u++;
529 
530  if (p.getProjection() == TotemRPUVPattern::projV)
531  pat_v++;
532  }
533  }
534 
535  // determine category
536  signed int category = -1;
537 
538  if (pl_u == 0 && pl_v == 0) category = 0; // empty
539 
540  if (category == -1 && pat_u + pat_v <= 1)
541  {
542  if (pl_u + pl_v < 6)
543  category = 1; // insuff
544  else
545  category = 4; // shower
546  }
547 
548  if (pat_u == 1 && pat_v == 1) category = 2; // 1-track
549 
550  if (category == -1) category = 3; // multi-track
551 
552  pp.event_category->Fill(category);
553  }
554 
555  // RP track-fit plots
556  set<unsigned int> rps_with_tracks;
557 
558  for (auto &ds : *tracks)
559  {
560  CTPPSDetId rpId(ds.detId());
561 
562  rps_with_tracks.insert(rpId);
563 
564  auto plIt = potPlots.find(rpId);
565  if (plIt == potPlots.end())
566  continue;
567  auto &plots = plIt->second;
568 
569  for (auto &ft : ds)
570  {
571  if (!ft.isValid())
572  continue;
573 
574  // number of planes contributing to (valid) fits
575  unsigned int n_pl_in_fit_u = 0, n_pl_in_fit_v = 0;
576  for (auto &hds : ft.getHits())
577  {
578  TotemRPDetId plId(hds.detId());
579  bool uProj = plId.isStripsCoordinateUDirection();
580 
581  for (auto &h : hds)
582  {
583  h.getPosition(); // just to keep compiler silent
584  if (uProj)
585  n_pl_in_fit_u++;
586  else
587  n_pl_in_fit_v++;
588  }
589  }
590 
591  plots.h_planes_fit_u->Fill(n_pl_in_fit_u);
592  plots.h_planes_fit_v->Fill(n_pl_in_fit_v);
593 
594  // mean position of U and V planes
595  TotemRPDetId plId_V(rpId); plId_V.setPlane(0);
596  TotemRPDetId plId_U(rpId); plId_U.setPlane(1);
597 
598  double rp_x = ( geometry->getSensor(plId_V)->translation().x() +
599  geometry->getSensor(plId_U)->translation().x() ) / 2.;
600  double rp_y = ( geometry->getSensor(plId_V)->translation().y() +
601  geometry->getSensor(plId_U)->translation().y() ) / 2.;
602 
603  // mean read-out direction of U and V planes
604  CLHEP::Hep3Vector rod_U = geometry->localToGlobalDirection(plId_U, CLHEP::Hep3Vector(0., 1., 0.));
605  CLHEP::Hep3Vector rod_V = geometry->localToGlobalDirection(plId_V, CLHEP::Hep3Vector(0., 1., 0.));
606 
607  double x = ft.getX0() - rp_x;
608  double y = ft.getY0() - rp_y;
609 
610  plots.trackHitsCumulativeHist->Fill(x, y);
611 
612  double U = x * rod_U.x() + y * rod_U.y();
613  double V = x * rod_V.x() + y * rod_V.y();
614 
615  plots.track_u_profile->Fill(U);
616  plots.track_v_profile->Fill(V);
617  }
618  }
619 
620  // restore trigger-sector map from digis
621  map<unsigned int, map<unsigned int, map<unsigned int, unsigned int>>> triggerSectorMap; // [rpId, U/V flag, sector] --> number of planes
622  for (const auto &dp : *digi)
623  {
624  TotemRPDetId plId(dp.detId());
625  CTPPSDetId rpId = plId.getRPId();
626  unsigned int uvFlag = (plId.isStripsCoordinateUDirection()) ? 0 : 1;
627 
628  set<unsigned int> sectors;
629  for (const auto &d : dp)
630  {
631  unsigned int sector = d.getStripNumber() / 32;
632  sectors.insert(sector);
633  }
634 
635  for (const auto &sector : sectors)
636  triggerSectorMap[rpId][uvFlag][sector]++;
637  }
638 
639  for (auto &rpp : triggerSectorMap)
640  {
641  const unsigned int rpId = rpp.first;
642 
643  // trigger sector is counted as active if at least 3 planes report activity
644 
645  set<unsigned int> triggerSectorsU;
646  for (const auto sp : rpp.second[0])
647  {
648  if (sp.second >= 3)
649  triggerSectorsU.insert(sp.first);
650  }
651 
652  set<unsigned int> triggerSectorsV;
653  for (const auto sp : rpp.second[1])
654  {
655  if (sp.second >= 3)
656  triggerSectorsV.insert(sp.first);
657  }
658 
659  auto plIt = potPlots.find(rpId);
660  if (plIt == potPlots.end())
661  continue;
662  auto &plots = plIt->second;
663 
664  const bool high_mult = (triggerSectorsU.size() > 2 && triggerSectorsV.size() > 2);
665 
666  const bool has_track = (rps_with_tracks.find(rpId) != rps_with_tracks.end());
667 
668  for (const auto &secU : triggerSectorsU)
669  {
670  for (const auto &secV : triggerSectorsV)
671  {
672  plots.triggerSectorUVCorrelation_all->Fill(secV, secU);
673 
674  if (!high_mult)
675  plots.triggerSectorUVCorrelation_mult2->Fill(secV, secU);
676 
677  if (has_track)
678  plots.triggerSectorUVCorrelation_track->Fill(secV, secU);
679  }
680  }
681  }
682 }
683 
684 //----------------------------------------------------------------------------------------------------
685 
Detector ID class for TOTEM Si strip detectors.
Definition: TotemRPDetId.h:30
MonitorElement * event_category
T getParameter(std::string const &) const
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * triggerSectorUVCorrelation_all
Translation translation() const
Definition: DetGeomDesc.h:66
MonitorElement * vfat_corruption
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
TH1F * getTH1F() const
MonitorElement * h_planes_fit_u
MonitorElement * activity_per_bx
MonitorElement * triggerSectorUVCorrelation_track
int bunchCrossing() const
Definition: EventBase.h:64
uint32_t plane() const
Definition: TotemRPDetId.h:46
MonitorElement * track_u_profile
bool isStripsCoordinateUDirection() const
Definition: TotemRPDetId.h:66
edm::EDGetTokenT< edm::DetSetVector< TotemRPDigi > > tokenDigi
Event setup record containing the real (actual) geometry information.
plots related to one RP plane
std::map< unsigned int, PlanePlots > planePlots
edm::EDGetTokenT< edm::DetSetVector< TotemRPRecHit > > tokenRecHit
const DetGeomDesc * getSensor(unsigned int id) const
returns geometry of a detector performs necessary checks, returns NULL if fails
unsigned int verbosity
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
~TotemRPDQMSource() override
edm::EDGetTokenT< edm::DetSetVector< TotemRPUVPattern > > tokenUVPattern
plots related to one RP
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
void rpName(std::string &name, NameFlag flag=nFull) const
Definition: CTPPSDetId.h:128
bool insert(Storage &iStorage, ItemType *iItem, const IdTag &iIdTag)
Definition: HCMethods.h:50
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * track_v_profile
CLHEP::Hep3Vector getRPTranslation(unsigned int id) const
edm::EDGetTokenT< edm::DetSetVector< TotemRPCluster > > tokenCluster
MonitorElement * hit_plane_hist
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
CLHEP::Hep3Vector localToGlobalDirection(unsigned int id, const CLHEP::Hep3Vector &) const
void setPlane(uint32_t det)
Definition: TotemRPDetId.h:48
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
TotemRPDQMSource(const edm::ParameterSet &ps)
MonitorElement * vfat_ec_bc_error
ESHandle< TrackerGeometry > geometry
MonitorElement * triggerSectorUVCorrelation_mult2
HLT enums.
T get() const
Definition: EventSetup.h:71
std::map< unsigned int, PotPlots > potPlots
void planeName(std::string &name, NameFlag flag=nFull) const
Definition: TotemRPDetId.h:79
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
MonitorElement * h_planes_fit_v
CLHEP::Hep3Vector globalToLocal(const DetGeomDesc *, const CLHEP::Hep3Vector &) const
edm::EDGetTokenT< edm::DetSetVector< TotemRPLocalTrack > > tokenLocalTrack
MonitorElement * activity_per_bx_short
edm::EDGetTokenT< edm::DetSetVector< TotemVFATStatus > > tokenStatus
Definition: event.py:1
Definition: Run.h:45
MonitorElement * trackHitsCumulativeHist