CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 //#include "RecoTotemRP/RPRecoDataFormats/interface/RPMulFittedTrackCollection.h"
30 
34 
35 #include <string>
36 
37 //----------------------------------------------------------------------------------------------------
38 
40 {
41  public:
43  virtual ~TotemRPDQMSource();
44 
45  protected:
46  void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
47  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
48  void analyze(edm::Event const& e, edm::EventSetup const& eSetup);
50  void endLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const& eSetup);
51  void endRun(edm::Run const& run, edm::EventSetup const& eSetup);
52 
53  private:
54  unsigned int verbosity;
55 
62  //edm::EDGetTokenT< RPMulFittedTrackCollection > tokenMultiTrackColl;
63 
65  struct GlobalPlots
66  {
69 
70  void Init(DQMStore::IBooker &ibooker);
71  };
72 
74 
77  {
78  int id;
79 
82 
84 
85  DiagonalPlots(DQMStore::IBooker &ibooker, int _id);
86  };
87 
88  std::map<unsigned int, DiagonalPlots> diagonalPlots;
89 
91  struct ArmPlots
92  {
93  int id;
94 
97 
99 
100  ArmPlots(DQMStore::IBooker &ibooker, int _id);
101  };
102 
103  std::map<unsigned int, ArmPlots> armPlots;
104 
107  {
109  StationPlots(DQMStore::IBooker &ibooker, int _id);
110  };
111 
112  std::map<unsigned int, StationPlots> stationPlots;
113 
115  struct PotPlots
116  {
118 
127 
128  PotPlots() {}
129  PotPlots(DQMStore::IBooker &ibooker, unsigned int id);
130  };
131 
132  std::map<unsigned int, PotPlots> potPlots;
133 
135  struct PlanePlots
136  {
142 
144  PlanePlots(DQMStore::IBooker &ibooker, unsigned int id);
145  };
146 
147  std::map<unsigned int, PlanePlots> planePlots;
148 };
149 
150 //----------------------------------------------------------------------------------------------------
151 //----------------------------------------------------------------------------------------------------
152 
153 using namespace std;
154 using namespace edm;
155 
156 //----------------------------------------------------------------------------------------------------
157 
159 {
160  ibooker.setCurrentFolder("CTPPS/TrackingStrip");
161 
162  events_per_bx = ibooker.book1D("events per BX", "rp;Event.BX", 4002, -1.5, 4000. + 0.5);
163  events_per_bx_short = ibooker.book1D("events per BX (short)", "rp;Event.BX", 102, -1.5, 100. + 0.5);
164 
165  h_trackCorr_hor = ibooker.book2D("track correlation RP-210-hor", "rp, 210, hor", 4, -0.5, 3.5, 4, -0.5, 3.5);
166  TH2F *hist = h_trackCorr_hor->getTH2F();
167  TAxis *xa = hist->GetXaxis(), *ya = hist->GetYaxis();
168  xa->SetBinLabel(1, "45, 210, near"); ya->SetBinLabel(1, "45, 210, near");
169  xa->SetBinLabel(2, "45, 210, far"); ya->SetBinLabel(2, "45, 210, far");
170  xa->SetBinLabel(3, "56, 210, near"); ya->SetBinLabel(3, "56, 210, near");
171  xa->SetBinLabel(4, "56, 210, far"); ya->SetBinLabel(4, "56, 210, far");
172 }
173 
174 //----------------------------------------------------------------------------------------------------
175 
177 {
178  bool top45 = id & 2;
179  bool top56 = id & 1;
180  bool diag = (top45 != top56);
181 
182  char name[50];
183  sprintf(name, "%s 45%s - 56%s",
184  (diag) ? "diagonal" : "antidiagonal",
185  (top45) ? "top" : "bot",
186  (top56) ? "top" : "bot"
187  );
188 
189  ibooker.setCurrentFolder(string("CTPPS/TrackingStrip/") + name);
190 
191  // TODO: define ranges! If defined automatically, can lead to problems when histograms are merged from several instances of the module.
192  h_lrc_x_d = ibooker.book2D("dx left vs right", string(name) + " : dx left vs. right, histogram;#Delta x_{45};#Delta x_{56}", 50, 0., 0., 50, 0., 0.);
193  h_lrc_x_n = ibooker.book2D("xn left vs right", string(name) + " : xn left vs. right, histogram;x^{N}_{45};x^{N}_{56}", 50, 0., 0., 50, 0., 0.);
194  h_lrc_x_f = ibooker.book2D("xf left vs right", string(name) + " : xf left vs. right, histogram;x^{F}_{45};x^{F}_{56}", 50, 0., 0., 50, 0., 0.);
195 
196  h_lrc_y_d = ibooker.book2D("dy left vs right", string(name) + " : dy left vs. right, histogram;#Delta y_{45};#Delta y_{56}", 50, 0., 0., 50, 0., 0.);
197  h_lrc_y_n = ibooker.book2D("yn left vs right", string(name) + " : yn left vs. right, histogram;y^{N}_{45};y^{N}_{56}", 50, 0., 0., 50, 0., 0.);
198  h_lrc_y_f = ibooker.book2D("yf left vs right", string(name) + " : yf left vs. right, histogram;y^{F}_{45};y^{F}_{56}", 50, 0., 0., 50, 0., 0.);
199 }
200 
201 //----------------------------------------------------------------------------------------------------
202 
204 {
205  string path;
207  ibooker.setCurrentFolder(path);
208 
209  string title;
211 
212  h_numRPWithTrack_top = ibooker.book1D("number of top RPs with tracks", title+";number of top RPs with tracks", 5, -0.5, 4.5);
213  h_numRPWithTrack_hor = ibooker.book1D("number of hor RPs with tracks", title+";number of hor RPs with tracks", 5, -0.5, 4.5);
214  h_numRPWithTrack_bot = ibooker.book1D("number of bot RPs with tracks", title+";number of bot RPs with tracks", 5, -0.5, 4.5);
215 
216  h_trackCorr = ibooker.book2D("track RP correlation", title, 13, -0.5, 12.5, 13, -0.5, 12.5);
217  TH2F *h_trackCorr_h = h_trackCorr->getTH2F();
218  TAxis *xa = h_trackCorr_h->GetXaxis(), *ya = h_trackCorr_h->GetYaxis();
219  xa->SetBinLabel(1, "210, near, top"); ya->SetBinLabel(1, "210, near, top");
220  xa->SetBinLabel(2, "bot"); ya->SetBinLabel(2, "bot");
221  xa->SetBinLabel(3, "hor"); ya->SetBinLabel(3, "hor");
222  xa->SetBinLabel(4, "far, hor"); ya->SetBinLabel(4, "far, hor");
223  xa->SetBinLabel(5, "top"); ya->SetBinLabel(5, "top");
224  xa->SetBinLabel(6, "bot"); ya->SetBinLabel(6, "bot");
225  xa->SetBinLabel(8, "220, near, top"); ya->SetBinLabel(8, "220, near, top");
226  xa->SetBinLabel(9, "bot"); ya->SetBinLabel(9, "bot");
227  xa->SetBinLabel(10, "hor"); ya->SetBinLabel(10, "hor");
228  xa->SetBinLabel(11, "far, hor"); ya->SetBinLabel(11, "far, hor");
229  xa->SetBinLabel(12, "top"); ya->SetBinLabel(12, "top");
230  xa->SetBinLabel(13, "bot"); ya->SetBinLabel(13, "bot");
231 
232  h_trackCorr_overlap = ibooker.book2D("track RP correlation hor-vert overlaps", title, 13, -0.5, 12.5, 13, -0.5, 12.5);
233  h_trackCorr_h = h_trackCorr_overlap->getTH2F();
234  xa = h_trackCorr_h->GetXaxis(); ya = h_trackCorr_h->GetYaxis();
235  xa->SetBinLabel(1, "210, near, top"); ya->SetBinLabel(1, "210, near, top");
236  xa->SetBinLabel(2, "bot"); ya->SetBinLabel(2, "bot");
237  xa->SetBinLabel(3, "hor"); ya->SetBinLabel(3, "hor");
238  xa->SetBinLabel(4, "far, hor"); ya->SetBinLabel(4, "far, hor");
239  xa->SetBinLabel(5, "top"); ya->SetBinLabel(5, "top");
240  xa->SetBinLabel(6, "bot"); ya->SetBinLabel(6, "bot");
241  xa->SetBinLabel(8, "220, near, top"); ya->SetBinLabel(8, "220, near, top");
242  xa->SetBinLabel(9, "bot"); ya->SetBinLabel(9, "bot");
243  xa->SetBinLabel(10, "hor"); ya->SetBinLabel(10, "hor");
244  xa->SetBinLabel(11, "far, hor"); ya->SetBinLabel(11, "far, hor");
245  xa->SetBinLabel(12, "top"); ya->SetBinLabel(12, "top");
246  xa->SetBinLabel(13, "bot"); ya->SetBinLabel(13, "bot");
247 }
248 
249 //----------------------------------------------------------------------------------------------------
250 
252 {
253  string path;
255  ibooker.setCurrentFolder(path);
256 }
257 
258 //----------------------------------------------------------------------------------------------------
259 //----------------------------------------------------------------------------------------------------
260 
262 {
263  string path;
265  ibooker.setCurrentFolder(path);
266 
267  string title;
269 
270  vfat_problem = ibooker.book2D("vfats with any problem", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
271  vfat_missing = ibooker.book2D("vfats missing", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
272  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);
273  vfat_corruption = ibooker.book2D("vfats with data corruption", title+";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
274 
275  activity = ibooker.book1D("active planes", title+";number of active planes", 11, -0.5, 10.5);
276  activity_u = ibooker.book1D("active planes U", title+";number of active U planes", 11, -0.5, 10.5);
277  activity_v = ibooker.book1D("active planes V", title+";number of active V planes", 11, -0.5, 10.5);
278 
279  activity_per_bx = ibooker.book1D("activity per BX", title+";Event.BX", 4002, -1.5, 4000. + 0.5);
280  activity_per_bx_short = ibooker.book1D("activity per BX (short)", title+";Event.BX", 102, -1.5, 100. + 0.5);
281 
282  hit_plane_hist = ibooker.book2D("activity in planes (2D)", title+";plane number;strip number", 10, -0.5, 9.5, 32, -0.5, 511.5);
283 
284  patterns_u = ibooker.book1D("recognized patterns U", title+";number of recognized U patterns", 11, -0.5, 10.5);
285  patterns_v = ibooker.book1D("recognized patterns V", title+";number of recognized V patterns", 11, -0.5, 10.5);
286 
287  h_planes_fit_u = ibooker.book1D("planes contributing to fit U", title+";number of planes contributing to U fit", 6, -0.5, 5.5);
288  h_planes_fit_v = ibooker.book1D("planes contributing to fit V", title+";number of planes contributing to V fit", 6, -0.5, 5.5);
289 
290  event_category = ibooker.book1D("event category", title+";event category", 5, -0.5, 4.5);
291  TH1F *event_category_h = event_category->getTH1F();
292  event_category_h->GetXaxis()->SetBinLabel(1, "empty");
293  event_category_h->GetXaxis()->SetBinLabel(2, "insufficient");
294  event_category_h->GetXaxis()->SetBinLabel(3, "single-track");
295  event_category_h->GetXaxis()->SetBinLabel(4, "multi-track");
296  event_category_h->GetXaxis()->SetBinLabel(5, "shower");
297 
298  trackHitsCumulativeHist = ibooker.book2D("track XY profile", title+";x (mm);y (mm)", 100, -18., +18., 100, -18., +18.);
299 
300  track_u_profile = ibooker.book1D("track profile U", title+"; U (mm)", 512, -256*66E-3, +256*66E-3);
301  track_v_profile = ibooker.book1D("track profile V", title+"; V (mm)", 512, -256*66E-3, +256*66E-3);
302 }
303 
304 //----------------------------------------------------------------------------------------------------
305 //----------------------------------------------------------------------------------------------------
306 
308 {
309  string path;
311  ibooker.setCurrentFolder(path);
312 
313  string title;
315 
316  digi_profile_cumulative = ibooker.book1D("digi profile", title+";strip number", 512, -0.5, 511.5);
317  cluster_profile_cumulative = ibooker.book1D("cluster profile", title+";cluster center", 1024, -0.25, 511.75);
318  hit_multiplicity = ibooker.book1D("hit multiplicity", title+";hits/detector/event", 6, -0.5, 5.5);
319  cluster_size = ibooker.book1D("cluster size", title+";hits per cluster", 5, 0.5, 5.5);
320 
321  efficiency_num = ibooker.book1D("efficiency num", title+";track position (mm)", 30, -15., 0.);
322  efficiency_den = ibooker.book1D("efficiency den", title+";track position (mm)", 30, -15., 0.);
323 }
324 
325 //----------------------------------------------------------------------------------------------------
326 //----------------------------------------------------------------------------------------------------
327 
329  verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0))
330 {
331  tokenStatus = consumes<DetSetVector<TotemVFATStatus>>(ps.getParameter<edm::InputTag>("tagStatus"));
332 
333  tokenDigi = consumes< DetSetVector<TotemRPDigi> >(ps.getParameter<edm::InputTag>("tagDigi"));
334  tokenCluster = consumes< edm::DetSetVector<TotemRPCluster> >(ps.getParameter<edm::InputTag>("tagCluster"));
335  tokenRecHit = consumes< edm::DetSetVector<TotemRPRecHit> >(ps.getParameter<edm::InputTag>("tagRecHit"));
336  tokenUVPattern = consumes< DetSetVector<TotemRPUVPattern> >(ps.getParameter<edm::InputTag>("tagUVPattern"));
337  tokenLocalTrack = consumes< DetSetVector<TotemRPLocalTrack> >(ps.getParameter<edm::InputTag>("tagLocalTrack"));
338  //tokenMultiTrackColl = consumes< RPMulFittedTrackCollection >(ps.getParameter<edm::InputTag>("tagMultiTrackColl"));
339 }
340 
341 //----------------------------------------------------------------------------------------------------
342 
344 {
345 }
346 
347 //----------------------------------------------------------------------------------------------------
348 
350 {
351 }
352 
353 //----------------------------------------------------------------------------------------------------
354 
356 {
357  ibooker.cd();
358  ibooker.setCurrentFolder("CTPPS");
359 
360  // global plots
361  globalPlots.Init(ibooker);
362 
363  // temporarily disabled
364  /*
365  // initialize diagonals
366  diagonalPlots[1] = DiagonalPlots(ibooker, 1); // 45 bot - 56 top
367  diagonalPlots[2] = DiagonalPlots(ibooker, 2); // 45 top - 45 bot
368 
369  // initialize anti-diagonals
370  diagonalPlots[0] = DiagonalPlots(ibooker, 0); // 45 bot - 56 bot
371  diagonalPlots[3] = DiagonalPlots(ibooker, 3); // 45 top - 56 top
372  */
373 
374  // loop over arms
375  for (unsigned int arm = 0; arm < 2; arm++)
376  {
377  TotemRPDetId armId(arm, 0);
378  armPlots[armId] = ArmPlots(ibooker, armId);
379 
380  // loop over stations
381  for (unsigned int st = 0; st < 3; st += 2)
382  {
383  TotemRPDetId stId(arm, st);
384  stationPlots[stId] = StationPlots(ibooker, stId);
385 
386  // loop over RPs
387  for (unsigned int rp = 0; rp < 6; ++rp)
388  {
389  TotemRPDetId rpId(arm, st, rp);
390  potPlots[rpId] = PotPlots(ibooker, rpId);
391 
392  // loop over planes
393  for (unsigned int pl = 0; pl < 10; ++pl)
394  {
395  TotemRPDetId plId(arm, st, rp, pl);
396  planePlots[plId] = PlanePlots(ibooker, plId);
397  }
398  }
399  }
400  }
401 }
402 
403 //----------------------------------------------------------------------------------------------------
404 
406  edm::EventSetup const& context)
407 {
408 }
409 
410 //----------------------------------------------------------------------------------------------------
411 
413 {
414  // get event setup data
416  eventSetup.get<VeryForwardRealGeometryRecord>().get(geometry);
417 
418  // get event data
420  event.getByToken(tokenStatus, status);
421 
423  event.getByToken(tokenDigi, digi);
424 
426  event.getByToken(tokenCluster, digCluster);
427 
429  event.getByToken(tokenRecHit, hits);
430 
432  event.getByToken(tokenUVPattern, patterns);
433 
435  event.getByToken(tokenLocalTrack, tracks);
436 
437  //Handle< RPMulFittedTrackCollection > multiTracks;
438  //event.getByToken(tokenMultiTrackColl, multiTracks);
439 
440  // check validity
441  bool valid = true;
442  valid &= status.isValid();
443  valid &= digi.isValid();
444  valid &= digCluster.isValid();
445  valid &= hits.isValid();
446  valid &= patterns.isValid();
447  valid &= tracks.isValid();
448  //valid &= multiTracks.isValid();
449 
450  if (!valid)
451  {
452  if (verbosity)
453  {
454  LogProblem("TotemRPDQMSource") <<
455  "ERROR in TotemDQMModuleRP::analyze > some of the required inputs are not valid. Skipping this event.\n"
456  << " status.isValid = " << status.isValid() << "\n"
457  << " digi.isValid = " << digi.isValid() << "\n"
458  << " digCluster.isValid = " << digCluster.isValid() << "\n"
459  << " hits.isValid = " << hits.isValid() << "\n"
460  << " patterns.isValid = " << patterns.isValid() << "\n"
461  << " tracks.isValid = " << tracks.isValid();
462  //<< " multiTracks.isValid = %i\n", multiTracks.isValid()
463  }
464 
465  return;
466  }
467 
468  //------------------------------
469  // Global Plots
470 
473 
474  for (auto &ds1 : *tracks)
475  {
476  for (auto &tr1 : ds1)
477  {
478  if (! tr1.isValid())
479  continue;
480 
481  CTPPSDetId rpId1(ds1.detId());
482  unsigned int arm1 = rpId1.arm();
483  unsigned int stNum1 = rpId1.station();
484  unsigned int rpNum1 = rpId1.rp();
485  if (stNum1 != 0 || (rpNum1 != 2 && rpNum1 != 3))
486  continue;
487  unsigned int idx1 = arm1*2 + rpNum1-2;
488 
489  for (auto &ds2 : *tracks)
490  {
491  for (auto &tr2 : ds2)
492  {
493  if (! tr2.isValid())
494  continue;
495 
496  CTPPSDetId rpId2(ds2.detId());
497  unsigned int arm2 = rpId2.arm();
498  unsigned int stNum2 = rpId2.station();
499  unsigned int rpNum2 = rpId2.rp();
500  if (stNum2 != 0 || (rpNum2 != 2 && rpNum2 != 3))
501  continue;
502  unsigned int idx2 = arm2*2 + rpNum2-2;
503 
504  globalPlots.h_trackCorr_hor->Fill(idx1, idx2);
505  }
506  }
507  }
508  }
509 
510  //------------------------------
511  // Status Plots
512 
513  for (auto &ds : *status)
514  {
515  TotemRPDetId detId(ds.detId());
516  unsigned int plNum = detId.plane();
517  CTPPSDetId rpId = detId.getRPId();
518 
519  auto &plots = potPlots[rpId];
520 
521  for (auto &s : ds)
522  {
523  if (s.isMissing())
524  {
525  plots.vfat_problem->Fill(plNum, s.getChipPosition());
526  plots.vfat_missing->Fill(plNum, s.getChipPosition());
527  }
528 
529  if (s.isECProgressError() || s.isBCProgressError())
530  {
531  plots.vfat_problem->Fill(plNum, s.getChipPosition());
532  plots.vfat_ec_bc_error->Fill(plNum, s.getChipPosition());
533  }
534 
535  if (s.isIDMismatch() || s.isFootprintError() || s.isCRCError())
536  {
537  plots.vfat_problem->Fill(plNum, s.getChipPosition());
538  plots.vfat_corruption->Fill(plNum, s.getChipPosition());
539  }
540  }
541  }
542 
543  //------------------------------
544  // Plane Plots
545 
546  // digi profile cumulative
547  for (DetSetVector<TotemRPDigi>::const_iterator it = digi->begin(); it != digi->end(); ++it)
548  {
549  TotemRPDetId detId(it->detId());
550  for (DetSet<TotemRPDigi>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
551  planePlots[detId].digi_profile_cumulative->Fill(dit->getStripNumber());
552  }
553 
554  // cluster profile cumulative
555  for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++)
556  {
557  TotemRPDetId detId(it->detId());
558  for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
559  planePlots[detId].cluster_profile_cumulative->Fill(dit->getCenterStripPosition());
560  }
561 
562  // hit multiplicity
563  for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++)
564  {
565  TotemRPDetId detId(it->detId());
566  planePlots[detId].hit_multiplicity->Fill(it->size());
567  }
568 
569  // cluster size
570  for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++)
571  {
572  TotemRPDetId detId(it->detId());
573  for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
574  planePlots[detId].cluster_size->Fill(dit->getNumberOfStrips());
575  }
576 
577  // plane efficiency plots
578  for (auto &ds : *tracks)
579  {
580  CTPPSDetId rpId(ds.detId());
581 
582  for (auto &ft : ds)
583  {
584  if (!ft.isValid())
585  continue;
586 
587  double rp_z = geometry->GetRPGlobalTranslation(rpId).z();
588 
589  for (unsigned int plNum = 0; plNum < 10; ++plNum)
590  {
591  TotemRPDetId plId = rpId;
592  plId.setPlane(plNum);
593 
594  double ft_z = ft.getZ0();
595  double ft_x = ft.getX0() + ft.getTx() * (ft_z - rp_z);
596  double ft_y = ft.getY0() + ft.getTy() * (ft_z - rp_z);
597 
598  double ft_v = geometry->GlobalToLocal(plId, CLHEP::Hep3Vector(ft_x, ft_y, ft_z)).y();
599 
600  bool hasMatchingHit = false;
601  const auto &hit_ds_it = hits->find(plId);
602  if (hit_ds_it != hits->end())
603  {
604  for (const auto &h : *hit_ds_it)
605  {
606  bool match = (fabs(ft_v - h.getPosition()) < 2.*0.066);
607  if (match)
608  {
609  hasMatchingHit = true;
610  break;
611  }
612  }
613  }
614 
615  auto &pp = planePlots[plId];
616 
617  pp.efficiency_den->Fill(ft_v);
618  if (hasMatchingHit)
619  pp.efficiency_num->Fill(ft_v);
620  }
621  }
622  }
623 
624 
625  //------------------------------
626  // Roman Pots Plots
627 
628  // determine active planes (from RecHits and VFATStatus)
629  map<unsigned int, set<unsigned int> > planes;
630  map<unsigned int, set<unsigned int> > planes_u;
631  map<unsigned int, set<unsigned int> > planes_v;
632  for (const auto &ds : *hits)
633  {
634  if (ds.empty())
635  continue;
636 
637  TotemRPDetId detId(ds.detId());
638  unsigned int planeNum = detId.plane();
639  CTPPSDetId rpId = detId.getRPId();
640 
641  planes[rpId].insert(planeNum);
642  if (detId.isStripsCoordinateUDirection())
643  planes_u[rpId].insert(planeNum);
644  else
645  planes_v[rpId].insert(planeNum);
646  }
647 
648  for (const auto &ds : *status)
649  {
650  bool activity = false;
651  for (const auto &s : ds)
652  {
653  if (s.isNumberOfClustersSpecified() && s.getNumberOfClusters() > 0)
654  {
655  activity = true;
656  break;
657  }
658  }
659 
660  if (!activity)
661  continue;
662 
663  TotemRPDetId detId(ds.detId());
664  unsigned int planeNum = detId.plane();
665  CTPPSDetId rpId = detId.getRPId();
666 
667  planes[rpId].insert(planeNum);
668  if (detId.isStripsCoordinateUDirection())
669  planes_u[rpId].insert(planeNum);
670  else
671  planes_v[rpId].insert(planeNum);
672  }
673 
674  // plane activity histogram
675  for (std::map<unsigned int, PotPlots>::iterator it = potPlots.begin(); it != potPlots.end(); it++)
676  {
677  it->second.activity->Fill(planes[it->first].size());
678  it->second.activity_u->Fill(planes_u[it->first].size());
679  it->second.activity_v->Fill(planes_v[it->first].size());
680 
681  if (planes[it->first].size() >= 6)
682  {
683  it->second.activity_per_bx->Fill(event.bunchCrossing());
684  it->second.activity_per_bx_short->Fill(event.bunchCrossing());
685  }
686  }
687 
688  for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++)
689  {
690  TotemRPDetId detId(it->detId());
691  unsigned int planeNum = detId.plane();
692  CTPPSDetId rpId = detId.getRPId();
693 
694  PotPlots &pp = potPlots[rpId];
695  for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
696  pp.hit_plane_hist->Fill(planeNum, dit->getCenterStripPosition());
697  }
698 
699  // recognized pattern histograms
700  for (auto &ds : *patterns)
701  {
702  CTPPSDetId rpId(ds.detId());
703 
704  PotPlots &pp = potPlots[rpId];
705 
706  // count U and V patterns
707  unsigned int u = 0, v = 0;
708  for (auto &p : ds)
709  {
710  if (! p.getFittable())
711  continue;
712 
713  if (p.getProjection() == TotemRPUVPattern::projU)
714  u++;
715 
716  if (p.getProjection() == TotemRPUVPattern::projV)
717  v++;
718  }
719 
720  pp.patterns_u->Fill(u);
721  pp.patterns_v->Fill(v);
722  }
723 
724  // event-category histogram
725  for (auto &it : potPlots)
726  {
727  TotemRPDetId rpId(it.first);
728  auto &pp = it.second;
729 
730  // process hit data for this plot
731  unsigned int pl_u = planes_u[rpId].size();
732  unsigned int pl_v = planes_v[rpId].size();
733 
734  // process pattern data for this pot
735  const auto &rp_pat_it = patterns->find(rpId);
736 
737  unsigned int pat_u = 0, pat_v = 0;
738  if (rp_pat_it != patterns->end())
739  {
740  for (auto &p : *rp_pat_it)
741  {
742  if (! p.getFittable())
743  continue;
744 
745  if (p.getProjection() == TotemRPUVPattern::projU)
746  pat_u++;
747 
748  if (p.getProjection() == TotemRPUVPattern::projV)
749  pat_v++;
750  }
751  }
752 
753  // determine category
754  signed int category = -1;
755 
756  if (pl_u == 0 && pl_v == 0) category = 0; // empty
757 
758  if (category == -1 && pat_u + pat_v <= 1)
759  {
760  if (pl_u + pl_v < 6)
761  category = 1; // insuff
762  else
763  category = 4; // shower
764  }
765 
766  if (pat_u == 1 && pat_v == 1) category = 2; // 1-track
767 
768  if (category == -1) category = 3; // multi-track
769 
770  pp.event_category->Fill(category);
771  }
772 
773  // RP track-fit plots
774  for (auto &ds : *tracks)
775  {
776  CTPPSDetId rpId(ds.detId());
777 
778  PotPlots &pp = potPlots[rpId];
779 
780  for (auto &ft : ds)
781  {
782  if (!ft.isValid())
783  continue;
784 
785  // number of planes contributing to (valid) fits
786  unsigned int n_pl_in_fit_u = 0, n_pl_in_fit_v = 0;
787  for (auto &hds : ft.getHits())
788  {
789  TotemRPDetId plId(hds.detId());
790  bool uProj = plId.isStripsCoordinateUDirection();
791 
792  for (auto &h : hds)
793  {
794  h.getPosition(); // just to keep compiler silent
795  if (uProj)
796  n_pl_in_fit_u++;
797  else
798  n_pl_in_fit_v++;
799  }
800  }
801 
802  pp.h_planes_fit_u->Fill(n_pl_in_fit_u);
803  pp.h_planes_fit_v->Fill(n_pl_in_fit_v);
804 
805  // mean position of U and V planes
806  TotemRPDetId plId_V(rpId); plId_V.setPlane(0);
807  TotemRPDetId plId_U(rpId); plId_U.setPlane(1);
808 
809  double rp_x = ( geometry->GetDetector(plId_V)->translation().x() +
810  geometry->GetDetector(plId_U)->translation().x() ) / 2.;
811  double rp_y = ( geometry->GetDetector(plId_V)->translation().y() +
812  geometry->GetDetector(plId_U)->translation().y() ) / 2.;
813 
814  // mean read-out direction of U and V planes
815  CLHEP::Hep3Vector rod_U = geometry->LocalToGlobalDirection(plId_U, CLHEP::Hep3Vector(0., 1., 0.));
816  CLHEP::Hep3Vector rod_V = geometry->LocalToGlobalDirection(plId_V, CLHEP::Hep3Vector(0., 1., 0.));
817 
818  double x = ft.getX0() - rp_x;
819  double y = ft.getY0() - rp_y;
820 
821  pp.trackHitsCumulativeHist->Fill(x, y);
822 
823  double U = x * rod_U.x() + y * rod_U.y();
824  double V = x * rod_V.x() + y * rod_V.y();
825 
826  pp.track_u_profile->Fill(U);
827  pp.track_v_profile->Fill(V);
828  }
829  }
830 
831  //------------------------------
832  // Station Plots
833 
834 
835  //------------------------------
836  // Arm Plots
837  {
838  map<unsigned int, unsigned int> mTop, mHor, mBot;
839 
840  for (auto p : armPlots)
841  {
842  mTop[p.first] = 0;
843  mHor[p.first] = 0;
844  mBot[p.first] = 0;
845  }
846 
847  for (auto &ds : *tracks)
848  {
849  CTPPSDetId rpId(ds.detId());
850  unsigned int rpNum = rpId.rp();
851  CTPPSDetId armId = rpId.getArmId();
852 
853  for (auto &tr : ds)
854  {
855  if (! tr.isValid())
856  continue;
857 
858  if (rpNum == 0 || rpNum == 4)
859  mTop[armId]++;
860  if (rpNum == 2 || rpNum == 3)
861  mHor[armId]++;
862  if (rpNum == 1 || rpNum == 5)
863  mBot[armId]++;
864  }
865  }
866 
867  for (auto &p : armPlots)
868  {
869  p.second.h_numRPWithTrack_top->Fill(mTop[p.first]);
870  p.second.h_numRPWithTrack_hor->Fill(mHor[p.first]);
871  p.second.h_numRPWithTrack_bot->Fill(mBot[p.first]);
872  }
873 
874  // track RP correlation
875  for (auto &ds1 : *tracks)
876  {
877  for (auto &tr1 : ds1)
878  {
879  if (! tr1.isValid())
880  continue;
881 
882  CTPPSDetId rpId1(ds1.detId());
883  unsigned int arm1 = rpId1.arm();
884  unsigned int stNum1 = rpId1.station();
885  unsigned int rpNum1 = rpId1.rp();
886  unsigned int idx1 = stNum1/2 * 7 + rpNum1;
887  bool hor1 = (rpNum1 == 2 || rpNum1 == 3);
888 
889  CTPPSDetId armId = rpId1.getArmId();
890  ArmPlots &ap = armPlots[armId];
891 
892  for (auto &ds2 : *tracks)
893  {
894  for (auto &tr2 : ds2)
895  {
896  if (! tr2.isValid())
897  continue;
898 
899  CTPPSDetId rpId2(ds2.detId());
900  unsigned int arm2 = rpId2.arm();
901  unsigned int stNum2 = rpId2.station();
902  unsigned int rpNum2 = rpId2.rp();
903  unsigned int idx2 = stNum2/2 * 7 + rpNum2;
904  bool hor2 = (rpNum2 == 2 || rpNum2 == 3);
905 
906  if (arm1 != arm2)
907  continue;
908 
909  ap.h_trackCorr->Fill(idx1, idx2);
910 
911  if (hor1 != hor2)
912  ap.h_trackCorr_overlap->Fill(idx1, idx2);
913  }
914  }
915  }
916  }
917  }
918 
919  //------------------------------
920  // RP-system plots
921  // TODO: this code needs
922  // * generalization for more than two RPs per arm
923  // * updating for tracks as DetSetVector
924  /*
925  for (auto &dp : diagonalPlots)
926  {
927  unsigned int id = dp.first;
928  bool top45 = id & 2;
929  bool top56 = id & 1;
930 
931  unsigned int id_45_n = (top45) ? 20 : 21;
932  unsigned int id_45_f = (top45) ? 24 : 25;
933  unsigned int id_56_n = (top56) ? 120 : 121;
934  unsigned int id_56_f = (top56) ? 124 : 125;
935 
936  bool h_45_n = (tracks->find(id_45_n) != tracks->end() && tracks->find(id_45_n)->second.IsValid());
937  bool h_45_f = (tracks->find(id_45_f) != tracks->end() && tracks->find(id_45_f)->second.IsValid());
938  bool h_56_n = (tracks->find(id_56_n) != tracks->end() && tracks->find(id_56_n)->second.IsValid());
939  bool h_56_f = (tracks->find(id_56_f) != tracks->end() && tracks->find(id_56_f)->second.IsValid());
940 
941  if (! (h_45_n && h_45_f && h_56_n && h_56_f) )
942  continue;
943 
944  double x_45_n = tracks->find(id_45_n)->second.X0(), y_45_n = tracks->find(id_45_n)->second.Y0();
945  double x_45_f = tracks->find(id_45_f)->second.X0(), y_45_f = tracks->find(id_45_f)->second.Y0();
946  double x_56_n = tracks->find(id_56_n)->second.X0(), y_56_n = tracks->find(id_56_n)->second.Y0();
947  double x_56_f = tracks->find(id_56_f)->second.X0(), y_56_f = tracks->find(id_56_f)->second.Y0();
948 
949  double dx_45 = x_45_f - x_45_n;
950  double dy_45 = y_45_f - y_45_n;
951  double dx_56 = x_56_f - x_56_n;
952  double dy_56 = y_56_f - y_56_n;
953 
954  DiagonalPlots &pl = dp.second;
955 
956  pl.h_lrc_x_d->Fill(dx_45, dx_56);
957  pl.h_lrc_y_d->Fill(dy_45, dy_56);
958 
959  pl.h_lrc_x_n->Fill(x_45_n, x_56_n);
960  pl.h_lrc_y_n->Fill(y_45_n, y_56_n);
961 
962  pl.h_lrc_x_f->Fill(x_45_f, x_56_f);
963  pl.h_lrc_y_f->Fill(y_45_f, y_56_f);
964  }
965  */
966 }
967 
968 //----------------------------------------------------------------------------------------------------
969 
971 {
972 }
973 
974 //----------------------------------------------------------------------------------------------------
975 
977 {
978 }
979 
980 //----------------------------------------------------------------------------------------------------
981 
Detector ID class for TOTEM Si strip detectors.
Definition: TotemRPDetId.h:30
T getParameter(std::string const &) const
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
GlobalPlots globalPlots
MonitorElement * cluster_profile_cumulative
PhiMemoryImage patterns[9]
void Init(DQMStore::IBooker &ibooker)
MonitorElement * h_trackCorr_overlap
MonitorElement * vfat_corruption
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
tuple pp
Definition: createTree.py:15
void cd(void)
Definition: DQMStore.cc:268
void beginLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &eSetup)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
tuple lumi
Definition: fjr2json.py:35
MonitorElement * activity_per_bx
std::map< unsigned int, ArmPlots > armPlots
void analyze(edm::Event const &e, edm::EventSetup const &eSetup)
int bunchCrossing() const
Definition: EventBase.h:65
#define NULL
Definition: scimark2.h:8
plots related to one arm
uint32_t plane() const
Definition: TotemRPDetId.h:49
MonitorElement * track_u_profile
bool isStripsCoordinateUDirection() const
Definition: TotemRPDetId.h:80
edm::EDGetTokenT< edm::DetSetVector< TotemRPDigi > > tokenDigi
Event setup record containing the real (actual) geometry information.
MonitorElement * h_numRPWithTrack_top
void Fill(long long x)
plots related to one RP plane
std::map< unsigned int, PlanePlots > planePlots
MonitorElement * h_numRPWithTrack_hor
edm::EDGetTokenT< edm::DetSetVector< TotemRPRecHit > > tokenRecHit
unsigned int verbosity
plots related to one (anti)diagonal
virtual ~TotemRPDQMSource()
edm::EDGetTokenT< edm::DetSetVector< TotemRPUVPattern > > tokenUVPattern
plots related to one RP
MonitorElement * h_numRPWithTrack_bot
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
void rpName(std::string &name, NameFlag flag=nFull) const
Definition: CTPPSDetId.h:139
CTPPSDetId getRPId() const
Definition: CTPPSDetId.h:97
bool insert(Storage &iStorage, ItemType *iItem, const IdTag &iIdTag)
Definition: HCMethods.h:49
uint32_t arm() const
Definition: CTPPSDetId.h:52
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:75
MonitorElement * events_per_bx_short
void armName(std::string &name, NameFlag flag=nFull) const
Definition: CTPPSDetId.h:115
MonitorElement * track_v_profile
void endLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &eSetup)
plots related to one station
iterator begin()
Definition: DetSet.h:59
edm::EDGetTokenT< edm::DetSetVector< TotemRPCluster > > tokenCluster
std::map< unsigned int, StationPlots > stationPlots
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
CTPPSDetId getArmId() const
Definition: CTPPSDetId.h:87
void endRun(edm::Run const &run, edm::EventSetup const &eSetup)
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:56
void setPlane(uint32_t det)
Definition: TotemRPDetId.h:54
TH1F * getTH1F(void) const
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
TotemRPDQMSource(const edm::ParameterSet &ps)
MonitorElement * vfat_ec_bc_error
ESHandle< TrackerGeometry > geometry
void stationName(std::string &name, NameFlag flag=nFull) const
Definition: CTPPSDetId.h:127
std::map< unsigned int, PotPlots > potPlots
void planeName(std::string &name, NameFlag flag=nFull) const
Definition: TotemRPDetId.h:105
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
plots related to the whole system
std::map< unsigned int, DiagonalPlots > diagonalPlots
TH2F * getTH2F(void) const
MonitorElement * digi_profile_cumulative
edm::EDGetTokenT< edm::DetSetVector< TotemRPLocalTrack > > tokenLocalTrack
MonitorElement * activity_per_bx_short
edm::EDGetTokenT< edm::DetSetVector< TotemVFATStatus > > tokenStatus
Definition: Run.h:43
uint32_t rp() const
Definition: CTPPSDetId.h:74
tuple status
Definition: mps_update.py:57
MonitorElement * trackHitsCumulativeHist