CMS 3D CMS Logo

Phase2OTHarvestTrackingParticles.cc
Go to the documentation of this file.
9 
11 public:
14  void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override;
15  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
16 
17 private:
18  // ----------member data ---------------------------
21 };
22 
24  : topFolderName_(iConfig.getParameter<std::string>("TopFolderName")) {}
25 
27 
28 // ------------ method called once each job just after ending the event loop
29 // ------------
31  using namespace edm;
32 
33  float eta_bins[] = {0.0, 0.7, 1.0, 1.2, 1.6, 2.0, 2.4};
34  int eta_binnum = 6;
35 
36  dbe = nullptr;
38 
39  if (dbe) {
40  // Find all monitor elements for histograms
41  MonitorElement *meN_eta = dbe->get(topFolderName_ + "/Efficiency/match_tp_eta");
42  MonitorElement *meD_eta = dbe->get(topFolderName_ + "/Efficiency/tp_eta");
43  MonitorElement *meN_pt = dbe->get(topFolderName_ + "/Efficiency/match_tp_pt");
44  MonitorElement *meD_pt = dbe->get(topFolderName_ + "/Efficiency/tp_pt");
45  MonitorElement *meN_pt_zoom = dbe->get(topFolderName_ + "/Efficiency/match_tp_pt_zoom");
46  MonitorElement *meD_pt_zoom = dbe->get(topFolderName_ + "/Efficiency/tp_pt_zoom");
47  MonitorElement *meN_d0 = dbe->get(topFolderName_ + "/Efficiency/match_tp_d0");
48  MonitorElement *meD_d0 = dbe->get(topFolderName_ + "/Efficiency/tp_d0");
49  MonitorElement *meN_VtxR = dbe->get(topFolderName_ + "/Efficiency/match_tp_VtxR");
50  MonitorElement *meD_VtxR = dbe->get(topFolderName_ + "/Efficiency/tp_VtxR");
51  MonitorElement *meN_VtxZ = dbe->get(topFolderName_ + "/Efficiency/match_tp_VtxZ");
52  MonitorElement *meD_VtxZ = dbe->get(topFolderName_ + "/Efficiency/tp_VtxZ");
53 
54  MonitorElement *merespt_eta0to0p7_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta0to0p7_pt2to3");
55  MonitorElement *merespt_eta0p7to1_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta0p7to1_pt2to3");
56  MonitorElement *merespt_eta1to1p2_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta1to1p2_pt2to3");
57  MonitorElement *merespt_eta1p2to1p6_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p2to1p6_pt2to3");
58  MonitorElement *merespt_eta1p6to2_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p6to2_pt2to3");
59  MonitorElement *merespt_eta2to2p4_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta2to2p4_pt2to3");
60  MonitorElement *merespt_eta0to0p7_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta0to0p7_pt3to8");
61  MonitorElement *merespt_eta0p7to1_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta0p7to1_pt3to8");
62  MonitorElement *merespt_eta1to1p2_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta1to1p2_pt3to8");
63  MonitorElement *merespt_eta1p2to1p6_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p2to1p6_pt3to8");
64  MonitorElement *merespt_eta1p6to2_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p6to2_pt3to8");
65  MonitorElement *merespt_eta2to2p4_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta2to2p4_pt3to8");
66  MonitorElement *merespt_eta0to0p7_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta0to0p7_pt8toInf");
67  MonitorElement *merespt_eta0p7to1_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta0p7to1_pt8toInf");
68  MonitorElement *merespt_eta1to1p2_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta1to1p2_pt8toInf");
69  MonitorElement *merespt_eta1p2to1p6_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta1p2to1p6_pt8toInf");
70  MonitorElement *merespt_eta1p6to2_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta1p6to2_pt8toInf");
71  MonitorElement *merespt_eta2to2p4_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta2to2p4_pt8toInf");
72 
73  MonitorElement *mereseta_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/reseta_eta0to0p7");
74  MonitorElement *mereseta_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/reseta_eta0p7to1");
75  MonitorElement *mereseta_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/reseta_eta1to1p2");
76  MonitorElement *mereseta_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/reseta_eta1p2to1p6");
77  MonitorElement *mereseta_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/reseta_eta1p6to2");
78  MonitorElement *mereseta_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/reseta_eta2to2p4");
79 
80  MonitorElement *meresphi_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/resphi_eta0to0p7");
81  MonitorElement *meresphi_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/resphi_eta0p7to1");
82  MonitorElement *meresphi_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/resphi_eta1to1p2");
83  MonitorElement *meresphi_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/resphi_eta1p2to1p6");
84  MonitorElement *meresphi_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/resphi_eta1p6to2");
85  MonitorElement *meresphi_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/resphi_eta2to2p4");
86 
87  MonitorElement *meresVtxZ_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta0to0p7");
88  MonitorElement *meresVtxZ_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta0p7to1");
89  MonitorElement *meresVtxZ_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta1to1p2");
90  MonitorElement *meresVtxZ_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta1p2to1p6");
91  MonitorElement *meresVtxZ_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta1p6to2");
92  MonitorElement *meresVtxZ_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta2to2p4");
93 
94  MonitorElement *meresd0_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/resd0_eta0to0p7");
95  MonitorElement *meresd0_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/resd0_eta0p7to1");
96  MonitorElement *meresd0_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/resd0_eta1to1p2");
97  MonitorElement *meresd0_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/resd0_eta1p2to1p6");
98  MonitorElement *meresd0_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/resd0_eta1p6to2");
99  MonitorElement *meresd0_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/resd0_eta2to2p4");
100 
101  if (meN_eta && meD_eta) {
102  // Get the numerator and denominator histograms
103  TH1F *numerator = meN_eta->getTH1F();
104  TH1F *denominator = meD_eta->getTH1F();
105  numerator->Sumw2();
106  denominator->Sumw2();
107 
108  // Set the current directory
109  dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
110 
111  // Book the new histogram to contain the results
112  MonitorElement *me_effic_eta = ibooker.book1D("EtaEfficiency",
113  "#eta efficiency",
114  numerator->GetNbinsX(),
115  numerator->GetXaxis()->GetXmin(),
116  numerator->GetXaxis()->GetXmax());
117 
118  // Calculate the efficiency
119  me_effic_eta->getTH1F()->Divide(numerator, denominator, 1., 1., "B");
120  me_effic_eta->setAxisTitle("tracking particle #eta");
121  me_effic_eta->getTH1F()->GetYaxis()->SetTitle("Efficiency");
122  me_effic_eta->getTH1F()->SetMaximum(1.0);
123  me_effic_eta->getTH1F()->SetMinimum(0.0);
124  me_effic_eta->getTH1F()->SetStats(false);
125  } // if ME found
126  else {
127  edm::LogWarning("DataNotFound") << "Monitor elements for eta efficiency cannot be found!\n";
128  }
129 
130  if (meN_pt && meD_pt) {
131  // Get the numerator and denominator histograms
132  TH1F *numerator2 = meN_pt->getTH1F();
133  numerator2->Sumw2();
134  TH1F *denominator2 = meD_pt->getTH1F();
135  denominator2->Sumw2();
136 
137  // Set the current directory
138  dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
139 
140  // Book the new histogram to contain the results
141  MonitorElement *me_effic_pt = ibooker.book1D("PtEfficiency",
142  "p_{T} efficiency",
143  numerator2->GetNbinsX(),
144  numerator2->GetXaxis()->GetXmin(),
145  numerator2->GetXaxis()->GetXmax());
146 
147  // Calculate the efficiency
148  me_effic_pt->getTH1F()->Divide(numerator2, denominator2, 1., 1., "B");
149  me_effic_pt->setAxisTitle("Tracking particle p_{T} [GeV]");
150  me_effic_pt->getTH1F()->GetYaxis()->SetTitle("Efficiency");
151  me_effic_pt->getTH1F()->SetMaximum(1.0);
152  me_effic_pt->getTH1F()->SetMinimum(0.0);
153  me_effic_pt->getTH1F()->SetStats(false);
154  } // if ME found
155  else {
156  edm::LogWarning("DataNotFound") << "Monitor elements for pT efficiency cannot be found!\n";
157  }
158 
159  if (meN_pt_zoom && meD_pt_zoom) {
160  // Get the numerator and denominator histograms
161  TH1F *numerator2_zoom = meN_pt_zoom->getTH1F();
162  numerator2_zoom->Sumw2();
163  TH1F *denominator2_zoom = meD_pt_zoom->getTH1F();
164  denominator2_zoom->Sumw2();
165 
166  // Set the current directory
167  dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
168 
169  // Book the new histogram to contain the results
170  MonitorElement *me_effic_pt_zoom = ibooker.book1D("PtEfficiency_zoom",
171  "p_{T} efficiency",
172  numerator2_zoom->GetNbinsX(),
173  numerator2_zoom->GetXaxis()->GetXmin(),
174  numerator2_zoom->GetXaxis()->GetXmax());
175 
176  // Calculate the efficiency
177  me_effic_pt_zoom->getTH1F()->Divide(numerator2_zoom, denominator2_zoom, 1., 1., "B");
178  me_effic_pt_zoom->setAxisTitle("Tracking particle p_{T} [GeV]");
179  me_effic_pt_zoom->getTH1F()->GetYaxis()->SetTitle("Efficiency");
180  me_effic_pt_zoom->getTH1F()->SetMaximum(1.0);
181  me_effic_pt_zoom->getTH1F()->SetMinimum(0.0);
182  me_effic_pt_zoom->getTH1F()->SetStats(false);
183  } // if ME found
184  else {
185  edm::LogWarning("DataNotFound") << "Monitor elements for zoom pT efficiency cannot be found!\n";
186  }
187 
188  if (meN_d0 && meD_d0) {
189  // Get the numerator and denominator histograms
190  TH1F *numerator5 = meN_d0->getTH1F();
191  numerator5->Sumw2();
192  TH1F *denominator5 = meD_d0->getTH1F();
193  denominator5->Sumw2();
194 
195  // Set the current directory
196  dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
197 
198  // Book the new histogram to contain the results
199  MonitorElement *me_effic_d0 = ibooker.book1D("d0Efficiency",
200  "d_{0} efficiency",
201  numerator5->GetNbinsX(),
202  numerator5->GetXaxis()->GetXmin(),
203  numerator5->GetXaxis()->GetXmax());
204 
205  // Calculate the efficiency
206  me_effic_d0->getTH1F()->Divide(numerator5, denominator5, 1., 1., "B");
207  me_effic_d0->setAxisTitle("Tracking particle d_{0} [cm]");
208  me_effic_d0->getTH1F()->GetYaxis()->SetTitle("Efficiency");
209  me_effic_d0->getTH1F()->SetMaximum(1.0);
210  me_effic_d0->getTH1F()->SetMinimum(0.0);
211  me_effic_d0->getTH1F()->SetStats(false);
212  } // if ME found
213  else {
214  edm::LogWarning("DataNotFound") << "Monitor elements for d0 efficiency cannot be found!\n";
215  }
216 
217  if (meN_VtxR && meD_VtxR) {
218  // Get the numerator and denominator histograms
219  TH1F *numerator6 = meN_VtxR->getTH1F();
220  numerator6->Sumw2();
221  TH1F *denominator6 = meD_VtxR->getTH1F();
222  denominator6->Sumw2();
223 
224  // Set the current directory
225  dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
226 
227  // Book the new histogram to contain the results
228  MonitorElement *me_effic_VtxR = ibooker.book1D("VtxREfficiency",
229  "Vtx R efficiency",
230  numerator6->GetNbinsX(),
231  numerator6->GetXaxis()->GetXmin(),
232  numerator6->GetXaxis()->GetXmax());
233 
234  // Calculate the efficiency
235  me_effic_VtxR->getTH1F()->Divide(numerator6, denominator6, 1., 1., "B");
236  me_effic_VtxR->setAxisTitle("Tracking particle VtxR [cm]");
237  me_effic_VtxR->getTH1F()->GetYaxis()->SetTitle("Efficiency");
238  me_effic_VtxR->getTH1F()->SetMaximum(1.0);
239  me_effic_VtxR->getTH1F()->SetMinimum(0.0);
240  me_effic_VtxR->getTH1F()->SetStats(false);
241  } // if ME found
242  else {
243  edm::LogWarning("DataNotFound") << "Monitor elements for VtxR efficiency cannot be found!\n";
244  }
245 
246  if (meN_VtxZ && meD_VtxZ) {
247  // Get the numerator and denominator histograms
248  TH1F *numerator7 = meN_VtxZ->getTH1F();
249  numerator7->Sumw2();
250  TH1F *denominator7 = meD_VtxZ->getTH1F();
251  denominator7->Sumw2();
252 
253  // Set the current directory
254  dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
255 
256  // Book the new histogram to contain the results
257  MonitorElement *me_effic_VtxZ = ibooker.book1D("VtxZEfficiency",
258  "Vtx Z efficiency",
259  numerator7->GetNbinsX(),
260  numerator7->GetXaxis()->GetXmin(),
261  numerator7->GetXaxis()->GetXmax());
262 
263  // Calculate the efficiency
264  me_effic_VtxZ->getTH1F()->Divide(numerator7, denominator7, 1., 1., "B");
265  me_effic_VtxZ->setAxisTitle("Tracking particle VtxZ [cm]");
266  me_effic_VtxZ->getTH1F()->GetYaxis()->SetTitle("Efficiency");
267  me_effic_VtxZ->getTH1F()->SetMaximum(1.0);
268  me_effic_VtxZ->getTH1F()->SetMinimum(0.0);
269  me_effic_VtxZ->getTH1F()->SetStats(false);
270  } // if ME found
271  else {
272  edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ efficiency cannot be found!\n";
273  }
274 
275  if (merespt_eta0to0p7_pt2to3 && merespt_eta0p7to1_pt2to3 && merespt_eta1to1p2_pt2to3 &&
276  merespt_eta1p2to1p6_pt2to3 && merespt_eta1p6to2_pt2to3 && merespt_eta2to2p4_pt2to3) {
277  // Set the current directoy
278  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
279 
280  // Grab the histograms
281  TH1F *resPt1a = merespt_eta0to0p7_pt2to3->getTH1F();
282  TH1F *resPt2a = merespt_eta0p7to1_pt2to3->getTH1F();
283  TH1F *resPt3a = merespt_eta1to1p2_pt2to3->getTH1F();
284  TH1F *resPt4a = merespt_eta1p2to1p6_pt2to3->getTH1F();
285  TH1F *resPt5a = merespt_eta1p6to2_pt2to3->getTH1F();
286  TH1F *resPt6a = merespt_eta2to2p4_pt2to3->getTH1F();
287 
288  // Book the new histogram to contain the results
289  MonitorElement *me_res_pt1 =
290  ibooker.book1D("pTResVsEta_2-3", "p_{T} resolution vs |#eta|, for p_{T}: 2-3 GeV", eta_binnum, eta_bins);
291  TH1F *resPt1 = me_res_pt1->getTH1F();
292  resPt1->GetXaxis()->SetTitle("tracking particle |#eta|");
293  resPt1->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
294  resPt1->SetMinimum(0.0);
295  resPt1->SetStats(false);
296 
297  std::vector<TH1F *> vResPt1 = {resPt1a, resPt2a, resPt3a, resPt4a, resPt5a, resPt6a};
298  for (int i = 0; i < 6; i++) {
299  resPt1->SetBinContent(i + 1, vResPt1[i]->GetStdDev());
300  resPt1->SetBinError(i + 1, vResPt1[i]->GetStdDevError());
301  }
302  } // if ME found
303  else {
304  edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (2-3) cannot be found!\n";
305  }
306 
307  if (merespt_eta0to0p7_pt3to8 && merespt_eta0p7to1_pt3to8 && merespt_eta1to1p2_pt3to8 &&
308  merespt_eta1p2to1p6_pt3to8 && merespt_eta1p6to2_pt3to8 && merespt_eta2to2p4_pt3to8) {
309  // Set the current directoy
310  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
311 
312  // Grab the histograms
313  TH1F *resPt1b = merespt_eta0to0p7_pt3to8->getTH1F();
314  TH1F *resPt2b = merespt_eta0p7to1_pt3to8->getTH1F();
315  TH1F *resPt3b = merespt_eta1to1p2_pt3to8->getTH1F();
316  TH1F *resPt4b = merespt_eta1p2to1p6_pt3to8->getTH1F();
317  TH1F *resPt5b = merespt_eta1p6to2_pt3to8->getTH1F();
318  TH1F *resPt6b = merespt_eta2to2p4_pt3to8->getTH1F();
319 
320  // Book the new histogram to contain the results
321  MonitorElement *me_res_pt2 =
322  ibooker.book1D("pTResVsEta_3-8", "p_{T} resolution vs |#eta|, for p_{T}: 3-8 GeV", eta_binnum, eta_bins);
323  TH1F *resPt2 = me_res_pt2->getTH1F();
324  resPt2->GetXaxis()->SetTitle("tracking particle |#eta|");
325  resPt2->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
326  resPt2->SetMinimum(0.0);
327  resPt2->SetStats(false);
328 
329  std::vector<TH1F *> vResPt2 = {resPt1b, resPt2b, resPt3b, resPt4b, resPt5b, resPt6b};
330  for (int i = 0; i < 6; i++) {
331  resPt2->SetBinContent(i + 1, vResPt2[i]->GetStdDev());
332  resPt2->SetBinError(i + 1, vResPt2[i]->GetStdDevError());
333  }
334  } // if ME found
335  else {
336  edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (3-8) cannot be found!\n";
337  }
338 
339  if (merespt_eta0to0p7_pt8toInf && merespt_eta0p7to1_pt8toInf && merespt_eta1to1p2_pt8toInf &&
340  merespt_eta1p2to1p6_pt8toInf && merespt_eta1p6to2_pt8toInf && merespt_eta2to2p4_pt8toInf) {
341  // Set the current directoy
342  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
343 
344  // Grab the histograms
345  TH1F *resPt1c = merespt_eta0to0p7_pt8toInf->getTH1F();
346  TH1F *resPt2c = merespt_eta0p7to1_pt8toInf->getTH1F();
347  TH1F *resPt3c = merespt_eta1to1p2_pt8toInf->getTH1F();
348  TH1F *resPt4c = merespt_eta1p2to1p6_pt8toInf->getTH1F();
349  TH1F *resPt5c = merespt_eta1p6to2_pt8toInf->getTH1F();
350  TH1F *resPt6c = merespt_eta2to2p4_pt8toInf->getTH1F();
351 
352  // Book the new histogram to contain the results
353  MonitorElement *me_res_pt3 =
354  ibooker.book1D("pTResVsEta_8-inf", "p_{T} resolution vs |#eta|, for p_{T}: >8 GeV", eta_binnum, eta_bins);
355  TH1F *resPt3 = me_res_pt3->getTH1F();
356  resPt3->GetXaxis()->SetTitle("tracking particle |#eta|");
357  resPt3->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
358  resPt3->SetMinimum(0.0);
359  resPt3->SetStats(false);
360 
361  std::vector<TH1F *> vResPt3 = {resPt1c, resPt2c, resPt3c, resPt4c, resPt5c, resPt6c};
362  for (int i = 0; i < 6; i++) {
363  resPt3->SetBinContent(i + 1, vResPt3[i]->GetStdDev());
364  resPt3->SetBinError(i + 1, vResPt3[i]->GetStdDevError());
365  }
366  } // if ME found
367  else {
368  edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (8-inf) cannot be found!\n";
369  }
370 
371  if (mereseta_eta0to0p7 && mereseta_eta0p7to1 && mereseta_eta1to1p2 && mereseta_eta1p2to1p6 && mereseta_eta1p6to2 &&
372  mereseta_eta2to2p4) {
373  // Set the current directoy
374  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
375 
376  // Grab the histograms
377  TH1F *resEta1 = mereseta_eta0to0p7->getTH1F();
378  TH1F *resEta2 = mereseta_eta0p7to1->getTH1F();
379  TH1F *resEta3 = mereseta_eta1to1p2->getTH1F();
380  TH1F *resEta4 = mereseta_eta1p2to1p6->getTH1F();
381  TH1F *resEta5 = mereseta_eta1p6to2->getTH1F();
382  TH1F *resEta6 = mereseta_eta2to2p4->getTH1F();
383 
384  // Book the new histogram to contain the results
385  MonitorElement *me_res_eta = ibooker.book1D("EtaResolution", "#eta resolution vs |#eta|", eta_binnum, eta_bins);
386  TH1F *resEta = me_res_eta->getTH1F();
387  resEta->GetXaxis()->SetTitle("tracking particle |#eta|");
388  resEta->GetYaxis()->SetTitle("#sigma(#Delta#eta)");
389  resEta->SetMinimum(0.0);
390  resEta->SetStats(false);
391 
392  std::vector<TH1F *> vResEta = {resEta1, resEta2, resEta3, resEta4, resEta5, resEta6};
393  for (int i = 0; i < 6; i++) {
394  resEta->SetBinContent(i + 1, vResEta[i]->GetStdDev());
395  resEta->SetBinError(i + 1, vResEta[i]->GetStdDevError());
396  }
397  } // if ME found
398  else {
399  edm::LogWarning("DataNotFound") << "Monitor elements for eta resolution cannot be found!\n";
400  }
401 
402  if (meresphi_eta0to0p7 && meresphi_eta0p7to1 && meresphi_eta1to1p2 && meresphi_eta1p2to1p6 && meresphi_eta1p6to2 &&
403  meresphi_eta2to2p4) {
404  // Set the current directoy
405  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
406 
407  // Grab the histograms
408  TH1F *resPhi1 = meresphi_eta0to0p7->getTH1F();
409  TH1F *resPhi2 = meresphi_eta0p7to1->getTH1F();
410  TH1F *resPhi3 = meresphi_eta1to1p2->getTH1F();
411  TH1F *resPhi4 = meresphi_eta1p2to1p6->getTH1F();
412  TH1F *resPhi5 = meresphi_eta1p6to2->getTH1F();
413  TH1F *resPhi6 = meresphi_eta2to2p4->getTH1F();
414 
415  // Book the new histogram to contain the results
416  MonitorElement *me_res_phi = ibooker.book1D("PhiResolution", "#phi resolution vs |#eta|", eta_binnum, eta_bins);
417  TH1F *resPhi = me_res_phi->getTH1F();
418  resPhi->GetXaxis()->SetTitle("tracking particle |#eta|");
419  resPhi->GetYaxis()->SetTitle("#sigma(#Delta#phi)");
420  resPhi->SetMinimum(0.0);
421  resPhi->SetStats(false);
422 
423  std::vector<TH1F *> vResPhi = {resPhi1, resPhi2, resPhi3, resPhi4, resPhi5, resPhi6};
424  for (int i = 0; i < 6; i++) {
425  resPhi->SetBinContent(i + 1, vResPhi[i]->GetStdDev());
426  resPhi->SetBinError(i + 1, vResPhi[i]->GetStdDevError());
427  }
428  } // if ME found
429  else {
430  edm::LogWarning("DataNotFound") << "Monitor elements for phi resolution cannot be found!\n";
431  }
432 
433  if (meresVtxZ_eta0to0p7 && meresVtxZ_eta0p7to1 && meresVtxZ_eta1to1p2 && meresVtxZ_eta1p2to1p6 &&
434  meresVtxZ_eta1p6to2 && meresVtxZ_eta2to2p4) {
435  // Set the current directoy
436  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
437 
438  // Grab the histograms
439  TH1F *resVtxZ_1 = meresVtxZ_eta0to0p7->getTH1F();
440  TH1F *resVtxZ_2 = meresVtxZ_eta0p7to1->getTH1F();
441  TH1F *resVtxZ_3 = meresVtxZ_eta1to1p2->getTH1F();
442  TH1F *resVtxZ_4 = meresVtxZ_eta1p2to1p6->getTH1F();
443  TH1F *resVtxZ_5 = meresVtxZ_eta1p6to2->getTH1F();
444  TH1F *resVtxZ_6 = meresVtxZ_eta2to2p4->getTH1F();
445 
446  // Book the new histogram to contain the results
447  MonitorElement *me_res_VtxZ = ibooker.book1D("VtxZResolution", "VtxZ resolution vs |#eta|", eta_binnum, eta_bins);
448  TH1F *resVtxZ = me_res_VtxZ->getTH1F();
449  resVtxZ->GetXaxis()->SetTitle("tracking particle |#eta|");
450  resVtxZ->GetYaxis()->SetTitle("#sigma(#DeltaVtxZ) [cm]");
451  resVtxZ->SetMinimum(0.0);
452  resVtxZ->SetStats(false);
453 
454  std::vector<TH1F *> vResVtxZ = {resVtxZ_1, resVtxZ_2, resVtxZ_3, resVtxZ_4, resVtxZ_5, resVtxZ_6};
455  for (int i = 0; i < 6; i++) {
456  resVtxZ->SetBinContent(i + 1, vResVtxZ[i]->GetStdDev());
457  resVtxZ->SetBinError(i + 1, vResVtxZ[i]->GetStdDevError());
458  }
459  } // if ME found
460  else {
461  edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ resolution cannot be found!\n";
462  }
463 
464  if (meresd0_eta0to0p7 && meresd0_eta0p7to1 && meresd0_eta1to1p2 && meresd0_eta1p2to1p6 && meresd0_eta1p6to2 &&
465  meresd0_eta2to2p4) {
466  // Set the current directoy
467  dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
468 
469  // Grab the histograms
470  TH1F *resd0_1 = meresd0_eta0to0p7->getTH1F();
471  TH1F *resd0_2 = meresd0_eta0p7to1->getTH1F();
472  TH1F *resd0_3 = meresd0_eta1to1p2->getTH1F();
473  TH1F *resd0_4 = meresd0_eta1p2to1p6->getTH1F();
474  TH1F *resd0_5 = meresd0_eta1p6to2->getTH1F();
475  TH1F *resd0_6 = meresd0_eta2to2p4->getTH1F();
476 
477  // Book the new histogram to contain the results
478  MonitorElement *me_res_d0 = ibooker.book1D("d0Resolution", "d_{0} resolution vs |#eta|", eta_binnum, eta_bins);
479  TH1F *resd0 = me_res_d0->getTH1F();
480  resd0->GetXaxis()->SetTitle("tracking particle |#eta|");
481  resd0->GetYaxis()->SetTitle("#sigma(#Deltad_{0}) [cm]");
482  resd0->SetMinimum(0.0);
483  resd0->SetStats(false);
484 
485  std::vector<TH1F *> vResD0 = {resd0_1, resd0_2, resd0_3, resd0_4, resd0_5, resd0_6};
486  for (int i = 0; i < 6; i++) {
487  resd0->SetBinContent(i + 1, vResD0[i]->GetStdDev());
488  resd0->SetBinError(i + 1, vResD0[i]->GetStdDevError());
489  }
490  } // if ME found
491  else {
492  edm::LogWarning("DataNotFound") << "Monitor elements for d0 resolution cannot be found!\n";
493  }
494 
495  } // if dbe found
496  else {
497  edm::LogWarning("DataNotFound") << "Cannot find valid DQM back end \n";
498  }
499 } // end dqmEndJob
500 
505  desc.add<std::string>("TopFolderName", "TrackerPhase2OTL1TrackV");
506  descriptions.add("Phase2OTHarvestTrackingParticles", desc);
507 }
void setCurrentFolder(std::string const &fullpath) override
Definition: DQMStore.h:647
void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Phase2OTHarvestTrackingParticles(const edm::ParameterSet &)
virtual TH1F * getTH1F() const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:712
HLT enums.
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)