CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
OuterTrackerMonitorTTTrack.cc
Go to the documentation of this file.
1 // Package: SiOuterTracker
2 // Class: SiOuterTracker
3 //
4 // Original Author: Isis Marina Van Parijs
5 // Modified by: Emily MacDonald (emily.kaelyn.macdonald@cern.ch)
6 
7 // system include files
8 #include <memory>
9 #include <numeric>
10 #include <string>
11 #include <vector>
12 #include <bitset>
13 
14 // user include files
23 
36 
39 
41 public:
43  ~OuterTrackerMonitorTTTrack() override;
44  void analyze(const edm::Event &, const edm::EventSetup &) override;
45  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
47  MonitorElement *Track_All_N = nullptr; // Number of tracks per event
48  MonitorElement *Track_All_NStubs = nullptr; // Number of stubs per track
49  MonitorElement *Track_All_NLayersMissed = nullptr; // Number of layers missed per track
50  MonitorElement *Track_All_Eta_NStubs = nullptr; // Number of stubs per track vs eta
51  MonitorElement *Track_All_Pt = nullptr; // pT distrubtion for tracks
52  MonitorElement *Track_All_Eta = nullptr; // eta distrubtion for tracks
53  MonitorElement *Track_All_Phi = nullptr; // phi distrubtion for tracks
54  MonitorElement *Track_All_D0 = nullptr; // d0 distrubtion for tracks
55  MonitorElement *Track_All_VtxZ = nullptr; // z0 distrubtion for tracks
56  MonitorElement *Track_All_BendChi2 = nullptr; // Bendchi2 distrubtion for tracks
57  MonitorElement *Track_All_Chi2 = nullptr; // chi2 distrubtion for tracks
58  MonitorElement *Track_All_Chi2Red = nullptr; // chi2/dof distrubtion for tracks
59  MonitorElement *Track_All_Chi2RZ = nullptr; // chi2 r-phi distrubtion for tracks
60  MonitorElement *Track_All_Chi2RPhi = nullptr; // chi2 r-z distrubtion for tracks
61  MonitorElement *Track_All_Chi2Red_NStubs = nullptr; // chi2/dof vs number of stubs
62  MonitorElement *Track_All_Chi2Red_Eta = nullptr; // chi2/dof vs eta of track
63  MonitorElement *Track_All_Eta_BarrelStubs = nullptr; // eta vs number of stubs in barrel
64  MonitorElement *Track_All_Eta_ECStubs = nullptr; // eta vs number of stubs in end caps
65  MonitorElement *Track_All_Chi2_Probability = nullptr; // chi2 probability
66 
68  // Quality cuts: chi2/dof<10, bendchi2<2.2 (Prompt), default in config
69  // Quality cuts: chi2/dof<40, bendchi2<2.4 (Extended/Displaced tracks)
70  MonitorElement *Track_HQ_N = nullptr; // Number of tracks per event
71  MonitorElement *Track_HQ_NStubs = nullptr; // Number of stubs per track
72  MonitorElement *Track_HQ_NLayersMissed = nullptr; // Number of layers missed per track
73  MonitorElement *Track_HQ_Eta_NStubs = nullptr; // Number of stubs per track vs eta
74  MonitorElement *Track_HQ_Pt = nullptr; // pT distrubtion for tracks
75  MonitorElement *Track_HQ_Eta = nullptr; // eta distrubtion for tracks
76  MonitorElement *Track_HQ_Phi = nullptr; // phi distrubtion for tracks
77  MonitorElement *Track_HQ_D0 = nullptr; // d0 distrubtion for tracks
78  MonitorElement *Track_HQ_VtxZ = nullptr; // z0 distrubtion for tracks
79  MonitorElement *Track_HQ_BendChi2 = nullptr; // Bendchi2 distrubtion for tracks
80  MonitorElement *Track_HQ_Chi2 = nullptr; // chi2 distrubtion for tracks
81  MonitorElement *Track_HQ_Chi2Red = nullptr; // chi2/dof distrubtion for tracks
82  MonitorElement *Track_HQ_Chi2RZ = nullptr; // chi2 r-z distrubtion for tracks
83  MonitorElement *Track_HQ_Chi2RPhi = nullptr; // chi2 r-phi distrubtion for tracks
84  MonitorElement *Track_HQ_Chi2Red_NStubs = nullptr; // chi2/dof vs number of stubs
85  MonitorElement *Track_HQ_Chi2Red_Eta = nullptr; // chi2/dof vs eta of track
86  MonitorElement *Track_HQ_Eta_BarrelStubs = nullptr; // eta vs number of stubs in barrel
87  MonitorElement *Track_HQ_Eta_ECStubs = nullptr; // eta vs number of stubs in end caps
88  MonitorElement *Track_HQ_Chi2_Probability = nullptr; // chi2 probability
89 
90 private:
93 
94  unsigned int HQNStubs_;
95  double HQChi2dof_;
96  double HQBendChi2_;
98 };
99 
100 // constructors and destructor
102  topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
103  ttTrackToken_ =
104  consumes<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>>(conf_.getParameter<edm::InputTag>("TTTracksTag"));
105  HQNStubs_ = conf_.getParameter<int>("HQNStubs");
106  HQChi2dof_ = conf_.getParameter<double>("HQChi2dof");
107  HQBendChi2_ = conf_.getParameter<double>("HQBendChi2");
108 }
109 
111  // do anything here that needs to be done at desctruction time
112  // (e.g. close files, deallocate resources etc.)
113 }
114 
115 // ------------ method called for each event ------------
117  // L1 Primaries
119  iEvent.getByToken(ttTrackToken_, TTTrackHandle);
120 
122  unsigned int numAllTracks = 0;
123  unsigned int numHQTracks = 0;
124 
125  // Adding protection
126  if (!TTTrackHandle.isValid())
127  return;
128 
130  unsigned int tkCnt = 0;
131  for (const auto &iterTTTrack : *TTTrackHandle) {
132  edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_>> tempTrackPtr(TTTrackHandle, tkCnt++);
133 
134  unsigned int nStubs = tempTrackPtr->getStubRefs().size();
135  int nBarrelStubs = 0;
136  int nECStubs = 0;
137 
138  float track_eta = tempTrackPtr->momentum().eta();
139  float track_d0 = tempTrackPtr->d0();
140  float track_bendchi2 = tempTrackPtr->stubPtConsistency();
141  float track_chi2 = tempTrackPtr->chi2();
142  float track_chi2dof = tempTrackPtr->chi2Red();
143  float track_chi2rz = tempTrackPtr->chi2Z();
144  float track_chi2rphi = tempTrackPtr->chi2XY();
145  int nLayersMissed = 0;
146  unsigned int hitPattern_ = (unsigned int)tempTrackPtr->hitPattern();
147 
148  int nbits = floor(log2(hitPattern_)) + 1;
149  int lay_i = 0;
150  bool seq = false;
151  for (int i = 0; i < nbits; i++) {
152  lay_i = ((1 << i) & hitPattern_) >> i; //0 or 1 in ith bit (right to left)
153  if (lay_i && !seq)
154  seq = true; //sequence starts when first 1 found
155  if (!lay_i && seq)
156  nLayersMissed++;
157  }
158 
159  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
160  theStubs = iterTTTrack.getStubRefs();
161  for (const auto &istub : theStubs) {
162  bool inBarrel = false;
163  bool inEC = false;
164 
165  if (istub->getDetId().subdetId() == StripSubdetector::TOB)
166  inBarrel = true;
167  else if (istub->getDetId().subdetId() == StripSubdetector::TID)
168  inEC = true;
169  if (inBarrel)
170  nBarrelStubs++;
171  else if (inEC)
172  nECStubs++;
173  } // end loop over stubs
174 
175  // HQ tracks: bendchi2<2.2 and chi2/dof<10
176  if (nStubs >= HQNStubs_ && track_chi2dof <= HQChi2dof_ && track_bendchi2 <= HQBendChi2_) {
177  numHQTracks++;
178 
179  Track_HQ_NStubs->Fill(nStubs);
180  Track_HQ_NLayersMissed->Fill(nLayersMissed);
181  Track_HQ_Eta_NStubs->Fill(track_eta, nStubs);
182  Track_HQ_Pt->Fill(tempTrackPtr->momentum().perp());
183  Track_HQ_Eta->Fill(track_eta);
184  Track_HQ_Phi->Fill(tempTrackPtr->momentum().phi());
185  Track_HQ_D0->Fill(track_d0);
186  Track_HQ_VtxZ->Fill(tempTrackPtr->z0());
187  Track_HQ_BendChi2->Fill(track_bendchi2);
188  Track_HQ_Chi2->Fill(track_chi2);
189  Track_HQ_Chi2RZ->Fill(track_chi2rz);
190  Track_HQ_Chi2RPhi->Fill(track_chi2rphi);
191  Track_HQ_Chi2Red->Fill(track_chi2dof);
192  Track_HQ_Chi2Red_NStubs->Fill(nStubs, track_chi2dof);
193  Track_HQ_Chi2Red_Eta->Fill(track_eta, track_chi2dof);
194  Track_HQ_Eta_BarrelStubs->Fill(track_eta, nBarrelStubs);
195  Track_HQ_Eta_ECStubs->Fill(track_eta, nECStubs);
197  }
198 
199  // All tracks (including HQ tracks)
200  numAllTracks++;
201  Track_All_NStubs->Fill(nStubs);
202  Track_All_NLayersMissed->Fill(nLayersMissed);
203  Track_All_Eta_NStubs->Fill(track_eta, nStubs);
204  Track_All_Pt->Fill(tempTrackPtr->momentum().perp());
205  Track_All_Eta->Fill(track_eta);
206  Track_All_Phi->Fill(tempTrackPtr->momentum().phi());
207  Track_All_D0->Fill(track_d0);
208  Track_All_VtxZ->Fill(tempTrackPtr->z0());
209  Track_All_BendChi2->Fill(track_bendchi2);
210  Track_All_Chi2->Fill(track_chi2);
211  Track_All_Chi2RZ->Fill(track_chi2rz);
212  Track_All_Chi2RPhi->Fill(track_chi2rphi);
213  Track_All_Chi2Red->Fill(track_chi2dof);
214  Track_All_Chi2Red_NStubs->Fill(nStubs, track_chi2dof);
215  Track_All_Chi2Red_Eta->Fill(track_eta, track_chi2dof);
216  Track_All_Eta_BarrelStubs->Fill(track_eta, nBarrelStubs);
217  Track_All_Eta_ECStubs->Fill(track_eta, nECStubs);
219 
220  } // End of loop over TTTracks
221 
222  Track_HQ_N->Fill(numHQTracks);
223  Track_All_N->Fill(numAllTracks);
224 } // end of method
225 
226 // ------------ method called once each job just before starting event loop
227 // ------------
228 // Creating all histograms for DQM file output
230  edm::Run const &run,
231  edm::EventSetup const &es) {
233 
235  iBooker.setCurrentFolder(topFolderName_ + "/Tracks/All");
236  // Nb of L1Tracks
237  HistoName = "Track_All_N";
238  edm::ParameterSet psTrack_N = conf_.getParameter<edm::ParameterSet>("TH1_NTracks");
239  Track_All_N = iBooker.book1D(HistoName,
240  HistoName,
241  psTrack_N.getParameter<int32_t>("Nbinsx"),
242  psTrack_N.getParameter<double>("xmin"),
243  psTrack_N.getParameter<double>("xmax"));
244  Track_All_N->setAxisTitle("# L1 Tracks", 1);
245  Track_All_N->setAxisTitle("# Events", 2);
246 
247  // Number of stubs
248  edm::ParameterSet psTrack_NStubs = conf_.getParameter<edm::ParameterSet>("TH1_NStubs");
249  HistoName = "Track_All_NStubs";
250  Track_All_NStubs = iBooker.book1D(HistoName,
251  HistoName,
252  psTrack_NStubs.getParameter<int32_t>("Nbinsx"),
253  psTrack_NStubs.getParameter<double>("xmin"),
254  psTrack_NStubs.getParameter<double>("xmax"));
255  Track_All_NStubs->setAxisTitle("# L1 Stubs per L1 Track", 1);
256  Track_All_NStubs->setAxisTitle("# L1 Tracks", 2);
257 
258  // Number of layers missed
259  HistoName = "Track_All_NLayersMissed";
260  Track_All_NLayersMissed = iBooker.book1D(HistoName,
261  HistoName,
262  psTrack_NStubs.getParameter<int32_t>("Nbinsx"),
263  psTrack_NStubs.getParameter<double>("xmin"),
264  psTrack_NStubs.getParameter<double>("xmax"));
265  Track_All_NLayersMissed->setAxisTitle("# Layers missed", 1);
266  Track_All_NLayersMissed->setAxisTitle("# L1 Tracks", 2);
267 
268  edm::ParameterSet psTrack_Eta_NStubs = conf_.getParameter<edm::ParameterSet>("TH2_Track_Eta_NStubs");
269  HistoName = "Track_All_Eta_NStubs";
270  Track_All_Eta_NStubs = iBooker.book2D(HistoName,
271  HistoName,
272  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsx"),
273  psTrack_Eta_NStubs.getParameter<double>("xmin"),
274  psTrack_Eta_NStubs.getParameter<double>("xmax"),
275  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsy"),
276  psTrack_Eta_NStubs.getParameter<double>("ymin"),
277  psTrack_Eta_NStubs.getParameter<double>("ymax"));
279  Track_All_Eta_NStubs->setAxisTitle("# L1 Stubs", 2);
280 
281  // Pt of the tracks
282  edm::ParameterSet psTrack_Pt = conf_.getParameter<edm::ParameterSet>("TH1_Track_Pt");
283  HistoName = "Track_All_Pt";
284  Track_All_Pt = iBooker.book1D(HistoName,
285  HistoName,
286  psTrack_Pt.getParameter<int32_t>("Nbinsx"),
287  psTrack_Pt.getParameter<double>("xmin"),
288  psTrack_Pt.getParameter<double>("xmax"));
289  Track_All_Pt->setAxisTitle("p_{T} [GeV]", 1);
290  Track_All_Pt->setAxisTitle("# L1 Tracks", 2);
291 
292  // Phi
293  edm::ParameterSet psTrack_Phi = conf_.getParameter<edm::ParameterSet>("TH1_Track_Phi");
294  HistoName = "Track_All_Phi";
295  Track_All_Phi = iBooker.book1D(HistoName,
296  HistoName,
297  psTrack_Phi.getParameter<int32_t>("Nbinsx"),
298  psTrack_Phi.getParameter<double>("xmin"),
299  psTrack_Phi.getParameter<double>("xmax"));
300  Track_All_Phi->setAxisTitle("#phi", 1);
301  Track_All_Phi->setAxisTitle("# L1 Tracks", 2);
302 
303  // D0
304  edm::ParameterSet psTrack_D0 = conf_.getParameter<edm::ParameterSet>("TH1_Track_D0");
305  HistoName = "Track_All_D0";
306  Track_All_D0 = iBooker.book1D(HistoName,
307  HistoName,
308  psTrack_D0.getParameter<int32_t>("Nbinsx"),
309  psTrack_D0.getParameter<double>("xmin"),
310  psTrack_D0.getParameter<double>("xmax"));
311  Track_All_D0->setAxisTitle("Track D0", 1);
312  Track_All_D0->setAxisTitle("# L1 Tracks", 2);
313 
314  // Eta
315  edm::ParameterSet psTrack_Eta = conf_.getParameter<edm::ParameterSet>("TH1_Track_Eta");
316  HistoName = "Track_All_Eta";
317  Track_All_Eta = iBooker.book1D(HistoName,
318  HistoName,
319  psTrack_Eta.getParameter<int32_t>("Nbinsx"),
320  psTrack_Eta.getParameter<double>("xmin"),
321  psTrack_Eta.getParameter<double>("xmax"));
322  Track_All_Eta->setAxisTitle("#eta", 1);
323  Track_All_Eta->setAxisTitle("# L1 Tracks", 2);
324 
325  // VtxZ
326  edm::ParameterSet psTrack_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1_Track_VtxZ");
327  HistoName = "Track_All_VtxZ";
328  Track_All_VtxZ = iBooker.book1D(HistoName,
329  HistoName,
330  psTrack_VtxZ.getParameter<int32_t>("Nbinsx"),
331  psTrack_VtxZ.getParameter<double>("xmin"),
332  psTrack_VtxZ.getParameter<double>("xmax"));
333  Track_All_VtxZ->setAxisTitle("L1 Track vertex position z [cm]", 1);
334  Track_All_VtxZ->setAxisTitle("# L1 Tracks", 2);
335 
336  // chi2
337  edm::ParameterSet psTrack_Chi2 = conf_.getParameter<edm::ParameterSet>("TH1_Track_Chi2");
338  HistoName = "Track_All_Chi2";
339  Track_All_Chi2 = iBooker.book1D(HistoName,
340  HistoName,
341  psTrack_Chi2.getParameter<int32_t>("Nbinsx"),
342  psTrack_Chi2.getParameter<double>("xmin"),
343  psTrack_Chi2.getParameter<double>("xmax"));
344  Track_All_Chi2->setAxisTitle("L1 Track #chi^{2}", 1);
345  Track_All_Chi2->setAxisTitle("# L1 Tracks", 2);
346 
347  // chi2 r-z
348  HistoName = "Track_All_Chi2RZ";
349  Track_All_Chi2RZ = iBooker.book1D(HistoName,
350  HistoName,
351  psTrack_Chi2.getParameter<int32_t>("Nbinsx"),
352  psTrack_Chi2.getParameter<double>("xmin"),
353  psTrack_Chi2.getParameter<double>("xmax"));
354  Track_All_Chi2RZ->setAxisTitle("L1 Track #chi^{2} r-z", 1);
355  Track_All_Chi2RZ->setAxisTitle("# L1 Tracks", 2);
356 
357  // chi2 r-phi
358  HistoName = "Track_All_Chi2RPhi";
359  Track_All_Chi2RPhi = iBooker.book1D(HistoName,
360  HistoName,
361  psTrack_Chi2.getParameter<int32_t>("Nbinsx"),
362  psTrack_Chi2.getParameter<double>("xmin"),
363  psTrack_Chi2.getParameter<double>("xmax"));
364  Track_All_Chi2RPhi->setAxisTitle("L1 Track #chi^{2}", 1);
365  Track_All_Chi2RPhi->setAxisTitle("# L1 Tracks", 2);
366 
367  // Bendchi2
368  edm::ParameterSet psTrack_Chi2R = conf_.getParameter<edm::ParameterSet>("TH1_Track_Chi2R");
369  HistoName = "Track_All_BendChi2";
370  Track_All_BendChi2 = iBooker.book1D(HistoName,
371  HistoName,
372  psTrack_Chi2R.getParameter<int32_t>("Nbinsx"),
373  psTrack_Chi2R.getParameter<double>("xmin"),
374  psTrack_Chi2R.getParameter<double>("xmax"));
375  Track_All_BendChi2->setAxisTitle("L1 Track Bend #chi^{2}", 1);
376  Track_All_BendChi2->setAxisTitle("# L1 Tracks", 2);
377 
378  // chi2Red
379  edm::ParameterSet psTrack_Chi2Red = conf_.getParameter<edm::ParameterSet>("TH1_Track_Chi2R");
380  HistoName = "Track_All_Chi2Red";
381  Track_All_Chi2Red = iBooker.book1D(HistoName,
382  HistoName,
383  psTrack_Chi2R.getParameter<int32_t>("Nbinsx"),
384  psTrack_Chi2R.getParameter<double>("xmin"),
385  psTrack_Chi2R.getParameter<double>("xmax"));
386  Track_All_Chi2Red->setAxisTitle("L1 Track #chi^{2}/ndf", 1);
387  Track_All_Chi2Red->setAxisTitle("# L1 Tracks", 2);
388 
389  // Chi2 prob
390  edm::ParameterSet psTrack_Chi2_Probability = conf_.getParameter<edm::ParameterSet>("TH1_Track_Chi2_Probability");
391  HistoName = "Track_All_Chi2_Probability";
392  Track_All_Chi2_Probability = iBooker.book1D(HistoName,
393  HistoName,
394  psTrack_Chi2_Probability.getParameter<int32_t>("Nbinsx"),
395  psTrack_Chi2_Probability.getParameter<double>("xmin"),
396  psTrack_Chi2_Probability.getParameter<double>("xmax"));
397  Track_All_Chi2_Probability->setAxisTitle("#chi^{2} probability", 1);
398  Track_All_Chi2_Probability->setAxisTitle("# L1 Tracks", 2);
399 
400  // Reduced chi2 vs #stubs
401  edm::ParameterSet psTrack_Chi2R_NStubs = conf_.getParameter<edm::ParameterSet>("TH2_Track_Chi2R_NStubs");
402  HistoName = "Track_All_Chi2Red_NStubs";
403  Track_All_Chi2Red_NStubs = iBooker.book2D(HistoName,
404  HistoName,
405  psTrack_Chi2R_NStubs.getParameter<int32_t>("Nbinsx"),
406  psTrack_Chi2R_NStubs.getParameter<double>("xmin"),
407  psTrack_Chi2R_NStubs.getParameter<double>("xmax"),
408  psTrack_Chi2R_NStubs.getParameter<int32_t>("Nbinsy"),
409  psTrack_Chi2R_NStubs.getParameter<double>("ymin"),
410  psTrack_Chi2R_NStubs.getParameter<double>("ymax"));
411  Track_All_Chi2Red_NStubs->setAxisTitle("# L1 Stubs", 1);
412  Track_All_Chi2Red_NStubs->setAxisTitle("L1 Track #chi^{2}/ndf", 2);
413 
414  // chi2/dof vs eta
415  edm::ParameterSet psTrack_Chi2R_Eta = conf_.getParameter<edm::ParameterSet>("TH2_Track_Chi2R_Eta");
416  HistoName = "Track_All_Chi2Red_Eta";
417  Track_All_Chi2Red_Eta = iBooker.book2D(HistoName,
418  HistoName,
419  psTrack_Chi2R_Eta.getParameter<int32_t>("Nbinsx"),
420  psTrack_Chi2R_Eta.getParameter<double>("xmin"),
421  psTrack_Chi2R_Eta.getParameter<double>("xmax"),
422  psTrack_Chi2R_Eta.getParameter<int32_t>("Nbinsy"),
423  psTrack_Chi2R_Eta.getParameter<double>("ymin"),
424  psTrack_Chi2R_Eta.getParameter<double>("ymax"));
426  Track_All_Chi2Red_Eta->setAxisTitle("L1 Track #chi^{2}/ndf", 2);
427 
428  // Eta vs #stubs in barrel
429  HistoName = "Track_All_Eta_BarrelStubs";
430  Track_All_Eta_BarrelStubs = iBooker.book2D(HistoName,
431  HistoName,
432  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsx"),
433  psTrack_Eta_NStubs.getParameter<double>("xmin"),
434  psTrack_Eta_NStubs.getParameter<double>("xmax"),
435  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsy"),
436  psTrack_Eta_NStubs.getParameter<double>("ymin"),
437  psTrack_Eta_NStubs.getParameter<double>("ymax"));
439  Track_All_Eta_BarrelStubs->setAxisTitle("# L1 Barrel Stubs", 2);
440 
441  // Eta vs #stubs in EC
442  HistoName = "Track_LQ_Eta_ECStubs";
443  Track_All_Eta_ECStubs = iBooker.book2D(HistoName,
444  HistoName,
445  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsx"),
446  psTrack_Eta_NStubs.getParameter<double>("xmin"),
447  psTrack_Eta_NStubs.getParameter<double>("xmax"),
448  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsy"),
449  psTrack_Eta_NStubs.getParameter<double>("ymin"),
450  psTrack_Eta_NStubs.getParameter<double>("ymax"));
452  Track_All_Eta_ECStubs->setAxisTitle("# L1 EC Stubs", 2);
453 
455  iBooker.setCurrentFolder(topFolderName_ + "/Tracks/HQ");
456  // Nb of L1Tracks
457  HistoName = "Track_HQ_N";
458  Track_HQ_N = iBooker.book1D(HistoName,
459  HistoName,
460  psTrack_N.getParameter<int32_t>("Nbinsx"),
461  psTrack_N.getParameter<double>("xmin"),
462  psTrack_N.getParameter<double>("xmax"));
463  Track_HQ_N->setAxisTitle("# L1 Tracks", 1);
464  Track_HQ_N->setAxisTitle("# Events", 2);
465 
466  // Number of stubs
467  HistoName = "Track_HQ_NStubs";
468  Track_HQ_NStubs = iBooker.book1D(HistoName,
469  HistoName,
470  psTrack_NStubs.getParameter<int32_t>("Nbinsx"),
471  psTrack_NStubs.getParameter<double>("xmin"),
472  psTrack_NStubs.getParameter<double>("xmax"));
473  Track_HQ_NStubs->setAxisTitle("# L1 Stubs per L1 Track", 1);
474  Track_HQ_NStubs->setAxisTitle("# L1 Tracks", 2);
475 
476  // Number of layers missed
477  HistoName = "Track_HQ_NLayersMissed";
478  Track_HQ_NLayersMissed = iBooker.book1D(HistoName,
479  HistoName,
480  psTrack_NStubs.getParameter<int32_t>("Nbinsx"),
481  psTrack_NStubs.getParameter<double>("xmin"),
482  psTrack_NStubs.getParameter<double>("xmax"));
483  Track_HQ_NLayersMissed->setAxisTitle("# Layers missed", 1);
484  Track_HQ_NLayersMissed->setAxisTitle("# L1 Tracks", 2);
485 
486  // Track eta vs #stubs
487  HistoName = "Track_HQ_Eta_NStubs";
488  Track_HQ_Eta_NStubs = iBooker.book2D(HistoName,
489  HistoName,
490  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsx"),
491  psTrack_Eta_NStubs.getParameter<double>("xmin"),
492  psTrack_Eta_NStubs.getParameter<double>("xmax"),
493  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsy"),
494  psTrack_Eta_NStubs.getParameter<double>("ymin"),
495  psTrack_Eta_NStubs.getParameter<double>("ymax"));
496  Track_HQ_Eta_NStubs->setAxisTitle("#eta", 1);
497  Track_HQ_Eta_NStubs->setAxisTitle("# L1 Stubs", 2);
498 
499  // Pt of the tracks
500  HistoName = "Track_HQ_Pt";
501  Track_HQ_Pt = iBooker.book1D(HistoName,
502  HistoName,
503  psTrack_Pt.getParameter<int32_t>("Nbinsx"),
504  psTrack_Pt.getParameter<double>("xmin"),
505  psTrack_Pt.getParameter<double>("xmax"));
506  Track_HQ_Pt->setAxisTitle("p_{T} [GeV]", 1);
507  Track_HQ_Pt->setAxisTitle("# L1 Tracks", 2);
508 
509  // Phi
510  HistoName = "Track_HQ_Phi";
511  Track_HQ_Phi = iBooker.book1D(HistoName,
512  HistoName,
513  psTrack_Phi.getParameter<int32_t>("Nbinsx"),
514  psTrack_Phi.getParameter<double>("xmin"),
515  psTrack_Phi.getParameter<double>("xmax"));
516  Track_HQ_Phi->setAxisTitle("#phi", 1);
517  Track_HQ_Phi->setAxisTitle("# L1 Tracks", 2);
518 
519  // D0
520  HistoName = "Track_HQ_D0";
521  Track_HQ_D0 = iBooker.book1D(HistoName,
522  HistoName,
523  psTrack_D0.getParameter<int32_t>("Nbinsx"),
524  psTrack_D0.getParameter<double>("xmin"),
525  psTrack_D0.getParameter<double>("xmax"));
526  Track_HQ_D0->setAxisTitle("Track D0", 1);
527  Track_HQ_D0->setAxisTitle("# L1 Tracks", 2);
528 
529  // Eta
530  HistoName = "Track_HQ_Eta";
531  Track_HQ_Eta = iBooker.book1D(HistoName,
532  HistoName,
533  psTrack_Eta.getParameter<int32_t>("Nbinsx"),
534  psTrack_Eta.getParameter<double>("xmin"),
535  psTrack_Eta.getParameter<double>("xmax"));
536  Track_HQ_Eta->setAxisTitle("#eta", 1);
537  Track_HQ_Eta->setAxisTitle("# L1 Tracks", 2);
538 
539  // VtxZ
540  HistoName = "Track_HQ_VtxZ";
541  Track_HQ_VtxZ = iBooker.book1D(HistoName,
542  HistoName,
543  psTrack_VtxZ.getParameter<int32_t>("Nbinsx"),
544  psTrack_VtxZ.getParameter<double>("xmin"),
545  psTrack_VtxZ.getParameter<double>("xmax"));
546  Track_HQ_VtxZ->setAxisTitle("L1 Track vertex position z [cm]", 1);
547  Track_HQ_VtxZ->setAxisTitle("# L1 Tracks", 2);
548 
549  // chi2
550  HistoName = "Track_HQ_Chi2";
551  Track_HQ_Chi2 = iBooker.book1D(HistoName,
552  HistoName,
553  psTrack_Chi2.getParameter<int32_t>("Nbinsx"),
554  psTrack_Chi2.getParameter<double>("xmin"),
555  psTrack_Chi2.getParameter<double>("xmax"));
556  Track_HQ_Chi2->setAxisTitle("L1 Track #chi^{2}", 1);
557  Track_HQ_Chi2->setAxisTitle("# L1 Tracks", 2);
558 
559  // Bendchi2
560  HistoName = "Track_HQ_BendChi2";
561  Track_HQ_BendChi2 = iBooker.book1D(HistoName,
562  HistoName,
563  psTrack_Chi2R.getParameter<int32_t>("Nbinsx"),
564  psTrack_Chi2R.getParameter<double>("xmin"),
565  psTrack_Chi2R.getParameter<double>("xmax"));
566  Track_HQ_BendChi2->setAxisTitle("L1 Track Bend #chi^{2}", 1);
567  Track_HQ_BendChi2->setAxisTitle("# L1 Tracks", 2);
568 
569  // chi2 r-z
570  HistoName = "Track_HQ_Chi2RZ";
571  Track_HQ_Chi2RZ = iBooker.book1D(HistoName,
572  HistoName,
573  psTrack_Chi2.getParameter<int32_t>("Nbinsx"),
574  psTrack_Chi2.getParameter<double>("xmin"),
575  psTrack_Chi2.getParameter<double>("xmax"));
576  Track_HQ_Chi2RZ->setAxisTitle("L1 Track #chi^{2} r-z", 1);
577  Track_HQ_Chi2RZ->setAxisTitle("# L1 Tracks", 2);
578 
579  HistoName = "Track_HQ_Chi2RPhi";
580  Track_HQ_Chi2RPhi = iBooker.book1D(HistoName,
581  HistoName,
582  psTrack_Chi2.getParameter<int32_t>("Nbinsx"),
583  psTrack_Chi2.getParameter<double>("xmin"),
584  psTrack_Chi2.getParameter<double>("xmax"));
585  Track_HQ_Chi2RPhi->setAxisTitle("L1 Track #chi^{2} r-phi", 1);
586  Track_HQ_Chi2RPhi->setAxisTitle("# L1 Tracks", 2);
587 
588  // chi2Red
589  HistoName = "Track_HQ_Chi2Red";
590  Track_HQ_Chi2Red = iBooker.book1D(HistoName,
591  HistoName,
592  psTrack_Chi2R.getParameter<int32_t>("Nbinsx"),
593  psTrack_Chi2R.getParameter<double>("xmin"),
594  psTrack_Chi2R.getParameter<double>("xmax"));
595  Track_HQ_Chi2Red->setAxisTitle("L1 Track #chi^{2}/ndf", 1);
596  Track_HQ_Chi2Red->setAxisTitle("# L1 Tracks", 2);
597 
598  // Chi2 prob
599  HistoName = "Track_HQ_Chi2_Probability";
600  Track_HQ_Chi2_Probability = iBooker.book1D(HistoName,
601  HistoName,
602  psTrack_Chi2_Probability.getParameter<int32_t>("Nbinsx"),
603  psTrack_Chi2_Probability.getParameter<double>("xmin"),
604  psTrack_Chi2_Probability.getParameter<double>("xmax"));
605  Track_HQ_Chi2_Probability->setAxisTitle("#chi^{2} probability", 1);
606  Track_HQ_Chi2_Probability->setAxisTitle("# L1 Tracks", 2);
607 
608  // Reduced chi2 vs #stubs
609  HistoName = "Track_HQ_Chi2Red_NStubs";
610  Track_HQ_Chi2Red_NStubs = iBooker.book2D(HistoName,
611  HistoName,
612  psTrack_Chi2R_NStubs.getParameter<int32_t>("Nbinsx"),
613  psTrack_Chi2R_NStubs.getParameter<double>("xmin"),
614  psTrack_Chi2R_NStubs.getParameter<double>("xmax"),
615  psTrack_Chi2R_NStubs.getParameter<int32_t>("Nbinsy"),
616  psTrack_Chi2R_NStubs.getParameter<double>("ymin"),
617  psTrack_Chi2R_NStubs.getParameter<double>("ymax"));
618  Track_HQ_Chi2Red_NStubs->setAxisTitle("# L1 Stubs", 1);
619  Track_HQ_Chi2Red_NStubs->setAxisTitle("L1 Track #chi^{2}/ndf", 2);
620 
621  // chi2/dof vs eta
622  HistoName = "Track_HQ_Chi2Red_Eta";
623  Track_HQ_Chi2Red_Eta = iBooker.book2D(HistoName,
624  HistoName,
625  psTrack_Chi2R_Eta.getParameter<int32_t>("Nbinsx"),
626  psTrack_Chi2R_Eta.getParameter<double>("xmin"),
627  psTrack_Chi2R_Eta.getParameter<double>("xmax"),
628  psTrack_Chi2R_Eta.getParameter<int32_t>("Nbinsy"),
629  psTrack_Chi2R_Eta.getParameter<double>("ymin"),
630  psTrack_Chi2R_Eta.getParameter<double>("ymax"));
632  Track_HQ_Chi2Red_Eta->setAxisTitle("L1 Track #chi^{2}/ndf", 2);
633 
634  // eta vs #stubs in barrel
635  HistoName = "Track_HQ_Eta_BarrelStubs";
636  Track_HQ_Eta_BarrelStubs = iBooker.book2D(HistoName,
637  HistoName,
638  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsx"),
639  psTrack_Eta_NStubs.getParameter<double>("xmin"),
640  psTrack_Eta_NStubs.getParameter<double>("xmax"),
641  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsy"),
642  psTrack_Eta_NStubs.getParameter<double>("ymin"),
643  psTrack_Eta_NStubs.getParameter<double>("ymax"));
645  Track_HQ_Eta_BarrelStubs->setAxisTitle("# L1 Barrel Stubs", 2);
646 
647  // eta vs #stubs in EC
648  HistoName = "Track_HQ_Eta_ECStubs";
649  Track_HQ_Eta_ECStubs = iBooker.book2D(HistoName,
650  HistoName,
651  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsx"),
652  psTrack_Eta_NStubs.getParameter<double>("xmin"),
653  psTrack_Eta_NStubs.getParameter<double>("xmax"),
654  psTrack_Eta_NStubs.getParameter<int32_t>("Nbinsy"),
655  psTrack_Eta_NStubs.getParameter<double>("ymin"),
656  psTrack_Eta_NStubs.getParameter<double>("ymax"));
658  Track_HQ_Eta_ECStubs->setAxisTitle("# L1 EC Stubs", 2);
659 
660 } // end of method
661 
MonitorElement * Track_All_N
Low-quality TTTracks (All tracks)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void Fill(long long x)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > ttTrackToken_
MonitorElement * Track_HQ_N
High-quality TTTracks; different depending on prompt vs displaced tracks.
float ChiSquaredProbability(double chiSquared, double nrDOF)
static constexpr auto TOB
bool isValid() const
Definition: HandleBase.h:70
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
OuterTrackerMonitorTTTrack(const edm::ParameterSet &)
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::string HistoName
void analyze(const edm::Event &, const edm::EventSetup &) override
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 constexpr auto TID
Definition: Run.h:45
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)