CMS 3D CMS Logo

OuterTrackerMCHarvester.cc
Go to the documentation of this file.
2 
4 {
5 }
6 
8 {
9 }
10 
11 // ------------ method called once each job just after ending the event loop ------------
13  using namespace edm;
14 
15  // Global variables
16  TF1* fit = new TF1("fit", "gaus", -0.01,0.01);
17  TF1* fit2 = new TF1("fit2", "gaus", -0.1,0.1);
18  TF1* fit3 = new TF1("fit3", "gaus", -1,1);
19 
20  std::vector<double> sigma_pt1;
21  std::vector<double> error_pt1;
22  std::vector<double> sigma_pt2;
23  std::vector<double> error_pt2;
24  std::vector<double> sigma_pt3;
25  std::vector<double> error_pt3;
26  std::vector<double> sigma_eta;
27  std::vector<double> error_eta;
28  std::vector<double> sigma_phi;
29  std::vector<double> error_phi;
30  std::vector<double> sigma_VtxZ;
31  std::vector<double> error_VtxZ;
32  std::vector<double> sigma_d0;
33  std::vector<double> error_d0;
34 
35  float eta_bins[]={0.0,0.7,1.0,1.2,1.6,2.0,2.4};
36  int eta_binnum = 6;
37 
38  dbe = nullptr;
40 
41  if (dbe) {
42  // Find all monitor elements for histograms
43  MonitorElement *meN_eta = dbe->get("SiOuterTrackerV/Tracks/Efficiency/match_tp_eta");
44  MonitorElement *meD_eta = dbe->get("SiOuterTrackerV/Tracks/Efficiency/tp_eta");
45  MonitorElement *meN_pt = dbe->get("SiOuterTrackerV/Tracks/Efficiency/match_tp_pt");
46  MonitorElement *meD_pt = dbe->get("SiOuterTrackerV/Tracks/Efficiency/tp_pt");
47  MonitorElement *meN_pt_zoom = dbe->get("SiOuterTrackerV/Tracks/Efficiency/match_tp_pt_zoom");
48  MonitorElement *meD_pt_zoom = dbe->get("SiOuterTrackerV/Tracks/Efficiency/tp_pt_zoom");
49  MonitorElement *meN_d0 = dbe->get("SiOuterTrackerV/Tracks/Efficiency/match_tp_d0");
50  MonitorElement *meD_d0 = dbe->get("SiOuterTrackerV/Tracks/Efficiency/tp_d0");
51  MonitorElement *meN_VtxR = dbe->get("SiOuterTrackerV/Tracks/Efficiency/match_tp_VtxR");
52  MonitorElement *meD_VtxR = dbe->get("SiOuterTrackerV/Tracks/Efficiency/tp_VtxR");
53  MonitorElement *meN_VtxZ = dbe->get("SiOuterTrackerV/Tracks/Efficiency/match_tp_VtxZ");
54  MonitorElement *meD_VtxZ = dbe->get("SiOuterTrackerV/Tracks/Efficiency/tp_VtxZ");
55 
56  MonitorElement *merespt_eta0to0p7_pt2to3 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta0to0p7_pt2to3");
57  MonitorElement *merespt_eta0p7to1_pt2to3 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta0p7to1_pt2to3");
58  MonitorElement *merespt_eta1to1p2_pt2to3 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1to1p2_pt2to3");
59  MonitorElement *merespt_eta1p2to1p6_pt2to3 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1p2to1p6_pt2to3");
60  MonitorElement *merespt_eta1p6to2_pt2to3 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1p6to2_pt2to3");
61  MonitorElement *merespt_eta2to2p4_pt2to3 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta2to2p4_pt2to3");
62  MonitorElement *merespt_eta0to0p7_pt3to8 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta0to0p7_pt3to8");
63  MonitorElement *merespt_eta0p7to1_pt3to8 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta0p7to1_pt3to8");
64  MonitorElement *merespt_eta1to1p2_pt3to8 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1to1p2_pt3to8");
65  MonitorElement *merespt_eta1p2to1p6_pt3to8 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1p2to1p6_pt3to8");
66  MonitorElement *merespt_eta1p6to2_pt3to8 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1p6to2_pt3to8");
67  MonitorElement *merespt_eta2to2p4_pt3to8 = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta2to2p4_pt3to8");
68  MonitorElement *merespt_eta0to0p7_pt8toInf = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta0to0p7_pt8toInf");
69  MonitorElement *merespt_eta0p7to1_pt8toInf = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta0p7to1_pt8toInf");
70  MonitorElement *merespt_eta1to1p2_pt8toInf = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1to1p2_pt8toInf");
71  MonitorElement *merespt_eta1p2to1p6_pt8toInf = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1p2to1p6_pt8toInf");
72  MonitorElement *merespt_eta1p6to2_pt8toInf = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta1p6to2_pt8toInf");
73  MonitorElement *merespt_eta2to2p4_pt8toInf = dbe->get("SiOuterTrackerV/Tracks/Resolution/respt_eta2to2p4_pt8toInf");
74 
75  MonitorElement *mereseta_eta0to0p7 = dbe->get("SiOuterTrackerV/Tracks/Resolution/reseta_eta0to0p7");
76  MonitorElement *mereseta_eta0p7to1 = dbe->get("SiOuterTrackerV/Tracks/Resolution/reseta_eta0p7to1");
77  MonitorElement *mereseta_eta1to1p2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/reseta_eta1to1p2");
78  MonitorElement *mereseta_eta1p2to1p6 = dbe->get("SiOuterTrackerV/Tracks/Resolution/reseta_eta1p2to1p6");
79  MonitorElement *mereseta_eta1p6to2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/reseta_eta1p6to2");
80  MonitorElement *mereseta_eta2to2p4 = dbe->get("SiOuterTrackerV/Tracks/Resolution/reseta_eta2to2p4");
81 
82  MonitorElement *meresphi_eta0to0p7 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resphi_eta0to0p7");
83  MonitorElement *meresphi_eta0p7to1 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resphi_eta0p7to1");
84  MonitorElement *meresphi_eta1to1p2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resphi_eta1to1p2");
85  MonitorElement *meresphi_eta1p2to1p6 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resphi_eta1p2to1p6");
86  MonitorElement *meresphi_eta1p6to2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resphi_eta1p6to2");
87  MonitorElement *meresphi_eta2to2p4 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resphi_eta2to2p4");
88 
89  MonitorElement *meresVtxZ_eta0to0p7 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resVtxZ_eta0to0p7");
90  MonitorElement *meresVtxZ_eta0p7to1 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resVtxZ_eta0p7to1");
91  MonitorElement *meresVtxZ_eta1to1p2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resVtxZ_eta1to1p2");
92  MonitorElement *meresVtxZ_eta1p2to1p6 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resVtxZ_eta1p2to1p6");
93  MonitorElement *meresVtxZ_eta1p6to2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resVtxZ_eta1p6to2");
94  MonitorElement *meresVtxZ_eta2to2p4 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resVtxZ_eta2to2p4");
95 
96  MonitorElement *meresd0_eta0to0p7 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resd0_eta0to0p7");
97  MonitorElement *meresd0_eta0p7to1 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resd0_eta0p7to1");
98  MonitorElement *meresd0_eta1to1p2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resd0_eta1to1p2");
99  MonitorElement *meresd0_eta1p2to1p6 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resd0_eta1p2to1p6");
100  MonitorElement *meresd0_eta1p6to2 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resd0_eta1p6to2");
101  MonitorElement *meresd0_eta2to2p4 = dbe->get("SiOuterTrackerV/Tracks/Resolution/resd0_eta2to2p4");
102 
103  if (meN_eta && meD_eta) {
104  // Get the numerator and denominator histograms
105  TH1F *numerator = meN_eta->getTH1F();
106  TH1F *denominator = meD_eta->getTH1F();
107  numerator->Sumw2();
108  denominator->Sumw2();
109 
110  // Set the current directory
111  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalEfficiency");
112 
113  // Book the new histogram to contain the results
114  MonitorElement *me_effic_eta =
115  ibooker.book1D("EtaEfficiency","#eta efficiency",
116  numerator->GetNbinsX(),
117  numerator->GetXaxis()->GetXmin(),
118  numerator->GetXaxis()->GetXmax());
119 
120  // Calculate the efficiency
121  me_effic_eta->getTH1F()->Divide(numerator, denominator, 1., 1., "B");
122  me_effic_eta->getTH1F()->GetXaxis()->SetTitle("tracking particle #eta");
123  me_effic_eta->getTH1F()->GetYaxis()->SetTitle("Efficiency");
124  me_effic_eta->getTH1F()->SetMaximum(1.0);
125  me_effic_eta->getTH1F()->SetMinimum(0.0);
126  me_effic_eta->getTH1F()->SetStats(false);
127  } //if ME found
128  else {edm::LogWarning("DataNotFound") << "Monitor elements for eta efficiency cannot be found!\n";}
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("SiOuterTrackerV/Tracks/FinalEfficiency");
139 
140  // Book the new histogram to contain the results
141  MonitorElement *me_effic_pt =
142  ibooker.book1D("PtEfficiency","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->getTH1F()->GetXaxis()->SetTitle("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 {edm::LogWarning("DataNotFound") << "Monitor elements for pT efficiency cannot be found!\n";}
156 
157  if (meN_pt_zoom && meD_pt_zoom) {
158  // Get the numerator and denominator histograms
159  TH1F *numerator2_zoom = meN_pt_zoom->getTH1F();
160  numerator2_zoom->Sumw2();
161  TH1F *denominator2_zoom = meD_pt_zoom->getTH1F();
162  denominator2_zoom->Sumw2();
163 
164  // Set the current directory
165  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalEfficiency");
166 
167  // Book the new histogram to contain the results
168  MonitorElement *me_effic_pt_zoom =
169  ibooker.book1D("PtEfficiency_zoom","p_{T} efficiency",
170  numerator2_zoom->GetNbinsX(),
171  numerator2_zoom->GetXaxis()->GetXmin(),
172  numerator2_zoom->GetXaxis()->GetXmax());
173 
174  // Calculate the efficiency
175  me_effic_pt_zoom->getTH1F()->Divide(numerator2_zoom, denominator2_zoom, 1., 1., "B");
176  me_effic_pt_zoom->getTH1F()->GetXaxis()->SetTitle("Tracking particle p_{T} [GeV]");
177  me_effic_pt_zoom->getTH1F()->GetYaxis()->SetTitle("Efficiency");
178  me_effic_pt_zoom->getTH1F()->SetMaximum(1.0);
179  me_effic_pt_zoom->getTH1F()->SetMinimum(0.0);
180  me_effic_pt_zoom->getTH1F()->SetStats(false);
181  } //if ME found
182  else {edm::LogWarning("DataNotFound") << "Monitor elements for zoom pT efficiency cannot be found!\n";}
183 
184  if (meN_d0 && meD_d0) {
185  // Get the numerator and denominator histograms
186  TH1F *numerator5 = meN_d0->getTH1F();
187  numerator5->Sumw2();
188  TH1F *denominator5 = meD_d0->getTH1F();
189  denominator5->Sumw2();
190 
191  // Set the current directory
192  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalEfficiency");
193 
194  // Book the new histogram to contain the results
195  MonitorElement *me_effic_d0 =
196  ibooker.book1D("d0Efficiency","d_{0} efficiency",
197  numerator5->GetNbinsX(),
198  numerator5->GetXaxis()->GetXmin(),
199  numerator5->GetXaxis()->GetXmax());
200 
201  // Calculate the efficiency
202  me_effic_d0->getTH1F()->Divide(numerator5, denominator5, 1., 1., "B");
203  me_effic_d0->getTH1F()->GetXaxis()->SetTitle("Tracking particle d_{0} [cm]");
204  me_effic_d0->getTH1F()->GetYaxis()->SetTitle("Efficiency");
205  me_effic_d0->getTH1F()->SetMaximum(1.0);
206  me_effic_d0->getTH1F()->SetMinimum(0.0);
207  me_effic_d0->getTH1F()->SetStats(false);
208  } //if ME found
209  else {edm::LogWarning("DataNotFound") << "Monitor elements for d0 efficiency cannot be found!\n";}
210 
211  if (meN_VtxR && meD_VtxR) {
212  // Get the numerator and denominator histograms
213  TH1F *numerator6 = meN_VtxR->getTH1F();
214  numerator6->Sumw2();
215  TH1F *denominator6 = meD_VtxR->getTH1F();
216  denominator6->Sumw2();
217 
218  // Set the current directory
219  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalEfficiency");
220 
221  // Book the new histogram to contain the results
222  MonitorElement *me_effic_VtxR =
223  ibooker.book1D("VtxREfficiency","Vtx R efficiency",
224  numerator6->GetNbinsX(),
225  numerator6->GetXaxis()->GetXmin(),
226  numerator6->GetXaxis()->GetXmax());
227 
228  // Calculate the efficiency
229  me_effic_VtxR->getTH1F()->Divide(numerator6, denominator6, 1., 1., "B");
230  me_effic_VtxR->getTH1F()->GetXaxis()->SetTitle("Tracking particle VtxR [cm]");
231  me_effic_VtxR->getTH1F()->GetYaxis()->SetTitle("Efficiency");
232  me_effic_VtxR->getTH1F()->SetMaximum(1.0);
233  me_effic_VtxR->getTH1F()->SetMinimum(0.0);
234  me_effic_VtxR->getTH1F()->SetStats(false);
235  } //if ME found
236  else {edm::LogWarning("DataNotFound") << "Monitor elements for VtxR efficiency cannot be found!\n";}
237 
238  if (meN_VtxZ && meD_VtxZ) {
239  // Get the numerator and denominator histograms
240  TH1F *numerator7 = meN_VtxZ->getTH1F();
241  numerator7->Sumw2();
242  TH1F *denominator7 = meD_VtxZ->getTH1F();
243  denominator7->Sumw2();
244 
245  // Set the current directory
246  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalEfficiency");
247 
248  // Book the new histogram to contain the results
249  MonitorElement *me_effic_VtxZ =
250  ibooker.book1D("VtxZEfficiency","Vtx Z efficiency",
251  numerator7->GetNbinsX(),
252  numerator7->GetXaxis()->GetXmin(),
253  numerator7->GetXaxis()->GetXmax());
254 
255  // Calculate the efficiency
256  me_effic_VtxZ->getTH1F()->Divide(numerator7, denominator7, 1., 1., "B");
257  me_effic_VtxZ->getTH1F()->GetXaxis()->SetTitle("Tracking particle VtxZ [cm]");
258  me_effic_VtxZ->getTH1F()->GetYaxis()->SetTitle("Efficiency");
259  me_effic_VtxZ->getTH1F()->SetMaximum(1.0);
260  me_effic_VtxZ->getTH1F()->SetMinimum(0.0);
261  me_effic_VtxZ->getTH1F()->SetStats(false);
262  } //if ME found
263  else {edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ efficiency cannot be found!\n";}
264 
265  if (merespt_eta0to0p7_pt2to3 && merespt_eta0p7to1_pt2to3 && merespt_eta1to1p2_pt2to3
266  && merespt_eta1p2to1p6_pt2to3 && merespt_eta1p6to2_pt2to3 && merespt_eta2to2p4_pt2to3) {
267  // Set the current directoy
268  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
269 
270  // Grab the histograms
271  TH1F *resPt1a = merespt_eta0to0p7_pt2to3->getTH1F();
272  TH1F *resPt2a = merespt_eta0p7to1_pt2to3->getTH1F();
273  TH1F *resPt3a = merespt_eta1to1p2_pt2to3->getTH1F();
274  TH1F *resPt4a = merespt_eta1p2to1p6_pt2to3->getTH1F();
275  TH1F *resPt5a = merespt_eta1p6to2_pt2to3->getTH1F();
276  TH1F *resPt6a = merespt_eta2to2p4_pt2to3->getTH1F();
277 
278  // Book the new histogram to contain the results
279  MonitorElement *me_res_pt1 = ibooker.book1D("pTResVsEta_2-3","p_{T} resolution vs |#eta|, for p_{T}: 2-3 GeV",eta_binnum,eta_bins);
280  TH1F *resPt1 = me_res_pt1->getTH1F();
281  resPt1->GetXaxis()->SetTitle("tracking particle |#eta|");
282  resPt1->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
283  resPt1->SetMinimum(0.0);
284  resPt1->SetStats(false);
285 
286  int testNumEntries1 = resPt1a->GetEntries();
287  if (testNumEntries1>0) {
288  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
289  resPt1a->Fit(fit2,"R"); resPt2a->Fit(fit2,"R"); resPt3a->Fit(fit2,"R");
290  resPt4a->Fit(fit2,"R"); resPt5a->Fit(fit2,"R"); resPt6a->Fit(fit2,"R");
291  sigma_pt1.push_back(resPt1a->GetFunction("fit2")->GetParameter(2));
292  sigma_pt1.push_back(resPt2a->GetFunction("fit2")->GetParameter(2));
293  sigma_pt1.push_back(resPt3a->GetFunction("fit2")->GetParameter(2));
294  sigma_pt1.push_back(resPt4a->GetFunction("fit2")->GetParameter(2));
295  sigma_pt1.push_back(resPt5a->GetFunction("fit2")->GetParameter(2));
296  sigma_pt1.push_back(resPt6a->GetFunction("fit2")->GetParameter(2));
297  error_pt1.push_back(resPt1a->GetFunction("fit2")->GetParError(2));
298  error_pt1.push_back(resPt2a->GetFunction("fit2")->GetParError(2));
299  error_pt1.push_back(resPt3a->GetFunction("fit2")->GetParError(2));
300  error_pt1.push_back(resPt4a->GetFunction("fit2")->GetParError(2));
301  error_pt1.push_back(resPt5a->GetFunction("fit2")->GetParError(2));
302  error_pt1.push_back(resPt6a->GetFunction("fit2")->GetParError(2));
303 
304  // Fill the new histogram to create resolution plot
305  for(int i=0;i<6;i++) {
306  resPt1->SetBinContent(i+1,sigma_pt1[i]);
307  resPt1->SetBinError(i+1,error_pt1[i]);
308  }
309  }
310  else {
311  edm::LogWarning("DataNotFound")<<"L1 tracks not found for pT resolution (2-3)!\n";
312  for(int i=0;i<6;i++) {
313  resPt1->SetBinContent(i+1,-1);
314  resPt1->SetBinError(i+1,-1);
315  }
316  }
317  } //if ME found
318  else {edm::LogWarning("DataNotFound")<<"Monitor elements for pT resolution (2-3) cannot be found!\n";}
319 
320  if (merespt_eta0to0p7_pt3to8 && merespt_eta0p7to1_pt3to8 && merespt_eta1to1p2_pt3to8
321  && merespt_eta1p2to1p6_pt3to8 && merespt_eta1p6to2_pt3to8 && merespt_eta2to2p4_pt3to8) {
322  // Set the current directoy
323  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
324 
325  // Grab the histograms
326  TH1F *resPt1b = merespt_eta0to0p7_pt3to8->getTH1F();
327  TH1F *resPt2b = merespt_eta0p7to1_pt3to8->getTH1F();
328  TH1F *resPt3b = merespt_eta1to1p2_pt3to8->getTH1F();
329  TH1F *resPt4b = merespt_eta1p2to1p6_pt3to8->getTH1F();
330  TH1F *resPt5b = merespt_eta1p6to2_pt3to8->getTH1F();
331  TH1F *resPt6b = merespt_eta2to2p4_pt3to8->getTH1F();
332 
333  // Book the new histogram to contain the results
334  MonitorElement *me_res_pt2 = ibooker.book1D("pTResVsEta_3-8","p_{T} resolution vs |#eta|, for p_{T}: 3-8 GeV",eta_binnum,eta_bins);
335  TH1F *resPt2 = me_res_pt2->getTH1F();
336  resPt2->GetXaxis()->SetTitle("tracking particle |#eta|");
337  resPt2->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
338  resPt2->SetMinimum(0.0);
339  resPt2->SetStats(false);
340 
341  int testNumEntries2 = resPt1b->GetEntries();
342  if (testNumEntries2>0) {
343  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
344  resPt1b->Fit(fit2,"R"); resPt2b->Fit(fit2,"R"); resPt3b->Fit(fit2,"R");
345  resPt4b->Fit(fit2,"R"); resPt5b->Fit(fit2,"R"); resPt6b->Fit(fit2,"R");
346  sigma_pt2.push_back(resPt1b->GetFunction("fit2")->GetParameter(2));
347  sigma_pt2.push_back(resPt2b->GetFunction("fit2")->GetParameter(2));
348  sigma_pt2.push_back(resPt3b->GetFunction("fit2")->GetParameter(2));
349  sigma_pt2.push_back(resPt4b->GetFunction("fit2")->GetParameter(2));
350  sigma_pt2.push_back(resPt5b->GetFunction("fit2")->GetParameter(2));
351  sigma_pt2.push_back(resPt6b->GetFunction("fit2")->GetParameter(2));
352  error_pt2.push_back(resPt1b->GetFunction("fit2")->GetParError(2));
353  error_pt2.push_back(resPt2b->GetFunction("fit2")->GetParError(2));
354  error_pt2.push_back(resPt3b->GetFunction("fit2")->GetParError(2));
355  error_pt2.push_back(resPt4b->GetFunction("fit2")->GetParError(2));
356  error_pt2.push_back(resPt5b->GetFunction("fit2")->GetParError(2));
357  error_pt2.push_back(resPt6b->GetFunction("fit2")->GetParError(2));
358 
359  // Fill the new histogram to create resolution plot
360  for(int i=0;i<6;i++) {
361  resPt2->SetBinContent(i+1,sigma_pt2[i]);
362  resPt2->SetBinError(i+1,error_pt2[i]);
363  }
364  }
365  else {
366  edm::LogWarning("DataNotFound")<<"L1 tracks not found for pT resolution (3-8)!\n";
367  for(int i=0;i<6;i++) {
368  resPt2->SetBinContent(i+1,-1);
369  resPt2->SetBinError(i+1,-1);
370  }
371  }
372  } //if ME found
373  else {edm::LogWarning("DataNotFound")<<"Monitor elements for pT resolution (3-8) cannot be found!\n";}
374 
375  if (merespt_eta0to0p7_pt8toInf && merespt_eta0p7to1_pt8toInf && merespt_eta1to1p2_pt8toInf
376  && merespt_eta1p2to1p6_pt8toInf && merespt_eta1p6to2_pt8toInf && merespt_eta2to2p4_pt8toInf) {
377  // Set the current directoy
378  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
379 
380  // Grab the histograms
381  TH1F *resPt1c = merespt_eta0to0p7_pt8toInf->getTH1F();
382  TH1F *resPt2c = merespt_eta0p7to1_pt8toInf->getTH1F();
383  TH1F *resPt3c = merespt_eta1to1p2_pt8toInf->getTH1F();
384  TH1F *resPt4c = merespt_eta1p2to1p6_pt8toInf->getTH1F();
385  TH1F *resPt5c = merespt_eta1p6to2_pt8toInf->getTH1F();
386  TH1F *resPt6c = merespt_eta2to2p4_pt8toInf->getTH1F();
387 
388  // Book the new histogram to contain the results
389  MonitorElement *me_res_pt3 = ibooker.book1D("pTResVsEta_8-inf","p_{T} resolution vs |#eta|, for p_{T}: >8 GeV",eta_binnum,eta_bins);
390  TH1F *resPt3 = me_res_pt3->getTH1F();
391  resPt3->GetXaxis()->SetTitle("tracking particle |#eta|");
392  resPt3->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
393  resPt3->SetMinimum(0.0);
394  resPt3->SetStats(false);
395 
396  int testNumEntries3 = resPt1c->GetEntries();
397  if (testNumEntries3>0) {
398  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
399  resPt1c->Fit(fit2,"R"); resPt2c->Fit(fit2,"R"); resPt3c->Fit(fit2,"R");
400  resPt4c->Fit(fit2,"R"); resPt5c->Fit(fit2,"R"); resPt6c->Fit(fit2,"R");
401  sigma_pt3.push_back(resPt1c->GetFunction("fit2")->GetParameter(2));
402  sigma_pt3.push_back(resPt2c->GetFunction("fit2")->GetParameter(2));
403  sigma_pt3.push_back(resPt3c->GetFunction("fit2")->GetParameter(2));
404  sigma_pt3.push_back(resPt4c->GetFunction("fit2")->GetParameter(2));
405  sigma_pt3.push_back(resPt5c->GetFunction("fit2")->GetParameter(2));
406  sigma_pt3.push_back(resPt6c->GetFunction("fit2")->GetParameter(2));
407  error_pt3.push_back(resPt1c->GetFunction("fit2")->GetParError(2));
408  error_pt3.push_back(resPt2c->GetFunction("fit2")->GetParError(2));
409  error_pt3.push_back(resPt3c->GetFunction("fit2")->GetParError(2));
410  error_pt3.push_back(resPt4c->GetFunction("fit2")->GetParError(2));
411  error_pt3.push_back(resPt5c->GetFunction("fit2")->GetParError(2));
412  error_pt3.push_back(resPt6c->GetFunction("fit2")->GetParError(2));
413 
414  // Fill the new histogram to create resolution plot
415  for(int i=0;i<6;i++) {
416  resPt3->SetBinContent(i+1,sigma_pt3[i]);
417  resPt3->SetBinError(i+1,error_pt3[i]);
418  }
419  }
420  else {
421  edm::LogWarning("DataNotFound")<<"L1 tracks not found for pT resolution (8-inf)!\n";
422  for(int i=0;i<6;i++) {
423  resPt3->SetBinContent(i+1,-1);
424  resPt3->SetBinError(i+1,-1);
425  }
426  }
427  } //if ME found
428  else {edm::LogWarning("DataNotFound")<<"Monitor elements for pT resolution (8-inf) cannot be found!\n";}
429 
430 
431  if (mereseta_eta0to0p7 && mereseta_eta0p7to1 && mereseta_eta1to1p2 && mereseta_eta1p2to1p6 && mereseta_eta1p6to2 && mereseta_eta2to2p4) {
432  // Set the current directoy
433  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
434 
435  // Grab the histograms
436  TH1F *resEta1 = mereseta_eta0to0p7->getTH1F();
437  TH1F *resEta2 = mereseta_eta0p7to1->getTH1F();
438  TH1F *resEta3 = mereseta_eta1to1p2->getTH1F();
439  TH1F *resEta4 = mereseta_eta1p2to1p6->getTH1F();
440  TH1F *resEta5 = mereseta_eta1p6to2->getTH1F();
441  TH1F *resEta6 = mereseta_eta2to2p4->getTH1F();
442 
443  // Book the new histogram to contain the results
444  MonitorElement *me_res_eta = ibooker.book1D("EtaResolution","#eta resolution vs |#eta|",eta_binnum,eta_bins);
445  TH1F *resEta = me_res_eta->getTH1F();
446  resEta->GetXaxis()->SetTitle("tracking particle |#eta|");
447  resEta->GetYaxis()->SetTitle("#sigma(#Delta#eta)");
448  resEta->SetMinimum(0.0);
449  resEta->SetStats(false);
450 
451  int testNumEntries4 = resEta1->GetEntries();
452  if (testNumEntries4>0) {
453 
454  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
455  resEta1->Fit(fit,"R"); resEta2->Fit(fit,"R"); resEta3->Fit(fit,"R");
456  resEta4->Fit(fit,"R"); resEta5->Fit(fit,"R"); resEta6->Fit(fit,"R");
457  sigma_eta.push_back(resEta1->GetFunction("fit")->GetParameter(2));
458  sigma_eta.push_back(resEta2->GetFunction("fit")->GetParameter(2));
459  sigma_eta.push_back(resEta3->GetFunction("fit")->GetParameter(2));
460  sigma_eta.push_back(resEta4->GetFunction("fit")->GetParameter(2));
461  sigma_eta.push_back(resEta5->GetFunction("fit")->GetParameter(2));
462  sigma_eta.push_back(resEta6->GetFunction("fit")->GetParameter(2));
463  error_eta.push_back(resEta1->GetFunction("fit")->GetParError(2));
464  error_eta.push_back(resEta2->GetFunction("fit")->GetParError(2));
465  error_eta.push_back(resEta3->GetFunction("fit")->GetParError(2));
466  error_eta.push_back(resEta4->GetFunction("fit")->GetParError(2));
467  error_eta.push_back(resEta5->GetFunction("fit")->GetParError(2));
468  error_eta.push_back(resEta6->GetFunction("fit")->GetParError(2));
469 
470  // Fill the new histogram to create resolution plot
471  for(int i=0;i<6;i++) {
472  resEta->SetBinContent(i+1,sigma_eta[i]);
473  resEta->SetBinError(i+1,error_eta[i]);
474  }
475  }
476  else {
477  edm::LogWarning("DataNotFound")<<"L1 tracks not found for eta resolution!\n";
478  for(int i=0;i<6;i++) {
479  resEta->SetBinContent(i+1,-1);
480  resEta->SetBinError(i+1,-1);
481  }
482  }
483  } //if ME found
484  else {edm::LogWarning("DataNotFound")<<"Monitor elements for eta resolution cannot be found!\n";}
485 
486  if (meresphi_eta0to0p7 && meresphi_eta0p7to1 && meresphi_eta1to1p2 && meresphi_eta1p2to1p6 && meresphi_eta1p6to2 && meresphi_eta2to2p4) {
487  // Set the current directoy
488  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
489 
490  // Grab the histograms
491  TH1F *resPhi1 = meresphi_eta0to0p7->getTH1F();
492  TH1F *resPhi2 = meresphi_eta0p7to1->getTH1F();
493  TH1F *resPhi3 = meresphi_eta1to1p2->getTH1F();
494  TH1F *resPhi4 = meresphi_eta1p2to1p6->getTH1F();
495  TH1F *resPhi5 = meresphi_eta1p6to2->getTH1F();
496  TH1F *resPhi6 = meresphi_eta2to2p4->getTH1F();
497 
498  // Book the new histogram to contain the results
499  MonitorElement *me_res_phi = ibooker.book1D("PhiResolution","#phi resolution vs |#eta|",eta_binnum,eta_bins);
500  TH1F *resPhi = me_res_phi->getTH1F();
501  resPhi->GetXaxis()->SetTitle("tracking particle |#eta|");
502  resPhi->GetYaxis()->SetTitle("#sigma(#Delta#phi)");
503  resPhi->SetMinimum(0.0);
504  resPhi->SetStats(false);
505 
506  int testNumEntries5 = resPhi1->GetEntries();
507  if (testNumEntries5>0) {
508  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
509  resPhi1->Fit(fit,"R"); resPhi2->Fit(fit,"R"); resPhi3->Fit(fit,"R");
510  resPhi4->Fit(fit,"R"); resPhi5->Fit(fit,"R"); resPhi6->Fit(fit,"R");
511  sigma_phi.push_back(resPhi1->GetFunction("fit")->GetParameter(2));
512  sigma_phi.push_back(resPhi2->GetFunction("fit")->GetParameter(2));
513  sigma_phi.push_back(resPhi3->GetFunction("fit")->GetParameter(2));
514  sigma_phi.push_back(resPhi4->GetFunction("fit")->GetParameter(2));
515  sigma_phi.push_back(resPhi5->GetFunction("fit")->GetParameter(2));
516  sigma_phi.push_back(resPhi6->GetFunction("fit")->GetParameter(2));
517  error_phi.push_back(resPhi1->GetFunction("fit")->GetParError(2));
518  error_phi.push_back(resPhi2->GetFunction("fit")->GetParError(2));
519  error_phi.push_back(resPhi3->GetFunction("fit")->GetParError(2));
520  error_phi.push_back(resPhi4->GetFunction("fit")->GetParError(2));
521  error_phi.push_back(resPhi5->GetFunction("fit")->GetParError(2));
522  error_phi.push_back(resPhi6->GetFunction("fit")->GetParError(2));
523 
524  // Fill the new histogram to create resolution plot
525  for(int i=0;i<6;i++) {
526  resPhi->SetBinContent(i+1,sigma_phi[i]);
527  resPhi->SetBinError(i+1,error_phi[i]);
528  }
529  }
530  else {
531  edm::LogWarning("DataNotFound")<<"L1 tracks not found for phi resolution!\n";
532  for(int i=0;i<6;i++) {
533  resPhi->SetBinContent(i+1,-1);
534  resPhi->SetBinError(i+1,-1);
535  }
536  }
537  } //if ME found
538  else {edm::LogWarning("DataNotFound")<<"Monitor elements for phi resolution cannot be found!\n";}
539 
540  if (meresVtxZ_eta0to0p7 && meresVtxZ_eta0p7to1 && meresVtxZ_eta1to1p2 && meresVtxZ_eta1p2to1p6 && meresVtxZ_eta1p6to2 && meresVtxZ_eta2to2p4) {
541  // Set the current directoy
542  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
543 
544  // Grab the histograms
545  TH1F *resVtxZ_1 = meresVtxZ_eta0to0p7->getTH1F();
546  TH1F *resVtxZ_2 = meresVtxZ_eta0p7to1->getTH1F();
547  TH1F *resVtxZ_3 = meresVtxZ_eta1to1p2->getTH1F();
548  TH1F *resVtxZ_4 = meresVtxZ_eta1p2to1p6->getTH1F();
549  TH1F *resVtxZ_5 = meresVtxZ_eta1p6to2->getTH1F();
550  TH1F *resVtxZ_6 = meresVtxZ_eta2to2p4->getTH1F();
551 
552  // Book the new histogram to contain the results
553  MonitorElement *me_res_VtxZ = ibooker.book1D("VtxZResolution","VtxZ resolution vs |#eta|",eta_binnum,eta_bins);
554  TH1F *resVtxZ = me_res_VtxZ->getTH1F();
555  resVtxZ->GetXaxis()->SetTitle("tracking particle |#eta|");
556  resVtxZ->GetYaxis()->SetTitle("#sigma(#DeltaVtxZ) [cm]");
557  resVtxZ->SetMinimum(0.0);
558  resVtxZ->SetStats(false);
559 
560  int testNumEntries6 = resVtxZ_1->GetEntries();
561  if (testNumEntries6>0) {
562  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
563  resVtxZ_1->Fit(fit3,"R"); resVtxZ_2->Fit(fit3,"R"); resVtxZ_3->Fit(fit3,"R");
564  resVtxZ_4->Fit(fit3,"R"); resVtxZ_5->Fit(fit3,"R"); resVtxZ_6->Fit(fit3,"R");
565  sigma_VtxZ.push_back(resVtxZ_1->GetFunction("fit3")->GetParameter(2));
566  sigma_VtxZ.push_back(resVtxZ_2->GetFunction("fit3")->GetParameter(2));
567  sigma_VtxZ.push_back(resVtxZ_3->GetFunction("fit3")->GetParameter(2));
568  sigma_VtxZ.push_back(resVtxZ_4->GetFunction("fit3")->GetParameter(2));
569  sigma_VtxZ.push_back(resVtxZ_5->GetFunction("fit3")->GetParameter(2));
570  sigma_VtxZ.push_back(resVtxZ_6->GetFunction("fit3")->GetParameter(2));
571  error_VtxZ.push_back(resVtxZ_1->GetFunction("fit3")->GetParError(2));
572  error_VtxZ.push_back(resVtxZ_2->GetFunction("fit3")->GetParError(2));
573  error_VtxZ.push_back(resVtxZ_3->GetFunction("fit3")->GetParError(2));
574  error_VtxZ.push_back(resVtxZ_4->GetFunction("fit3")->GetParError(2));
575  error_VtxZ.push_back(resVtxZ_5->GetFunction("fit3")->GetParError(2));
576  error_VtxZ.push_back(resVtxZ_6->GetFunction("fit3")->GetParError(2));
577 
578  // Fill the new histogram to create resolution plot
579  for(int i=0;i<6;i++) {
580  resVtxZ->SetBinContent(i+1,sigma_VtxZ[i]);
581  resVtxZ->SetBinError(i+1,error_VtxZ[i]);
582  }
583  }
584  else {
585  edm::LogWarning("DataNotFound")<<"L1 tracks not found for VtxZ resolution!\n";
586  for(int i=0;i<6;i++) {
587  resVtxZ->SetBinContent(i+1,-1);
588  resVtxZ->SetBinError(i+1,-1);
589  }
590  }
591  } //if ME found
592  else {edm::LogWarning("DataNotFound")<<"Monitor elements for VtxZ resolution cannot be found!\n";}
593 
594  if (meresd0_eta0to0p7 && meresd0_eta0p7to1 && meresd0_eta1to1p2 && meresd0_eta1p2to1p6 && meresd0_eta1p6to2 && meresd0_eta2to2p4) {
595  // Set the current directoy
596  dbe->setCurrentFolder("SiOuterTrackerV/Tracks/FinalResolution");
597 
598  // Grab the histograms
599  TH1F *resd0_1 = meresd0_eta0to0p7->getTH1F();
600  TH1F *resd0_2 = meresd0_eta0p7to1->getTH1F();
601  TH1F *resd0_3 = meresd0_eta1to1p2->getTH1F();
602  TH1F *resd0_4 = meresd0_eta1p2to1p6->getTH1F();
603  TH1F *resd0_5 = meresd0_eta1p6to2->getTH1F();
604  TH1F *resd0_6 = meresd0_eta2to2p4->getTH1F();
605 
606  // Book the new histogram to contain the results
607  MonitorElement *me_res_d0 = ibooker.book1D("d0Resolution","d_{0} resolution vs |#eta|",eta_binnum,eta_bins);
608  TH1F *resd0 = me_res_d0->getTH1F();
609  resd0->GetXaxis()->SetTitle("tracking particle |#eta|");
610  resd0->GetYaxis()->SetTitle("#sigma(#Deltad_{0}) [cm]");
611  resd0->SetMinimum(0.0);
612  resd0->SetStats(false);
613 
614  int testNumEntries7 = resd0_1->GetEntries();
615  if (testNumEntries7>0) {
616  // Fit the histograms with a gaussian curve - take sigma and the error from the fit
617  resd0_1->Fit(fit,"R"); resd0_2->Fit(fit,"R"); resd0_3->Fit(fit,"R");
618  resd0_4->Fit(fit,"R"); resd0_5->Fit(fit,"R"); resd0_6->Fit(fit,"R");
619  sigma_d0.push_back(resd0_1->GetFunction("fit")->GetParameter(2));
620  sigma_d0.push_back(resd0_2->GetFunction("fit")->GetParameter(2));
621  sigma_d0.push_back(resd0_3->GetFunction("fit")->GetParameter(2));
622  sigma_d0.push_back(resd0_4->GetFunction("fit")->GetParameter(2));
623  sigma_d0.push_back(resd0_5->GetFunction("fit")->GetParameter(2));
624  sigma_d0.push_back(resd0_6->GetFunction("fit")->GetParameter(2));
625  error_d0.push_back(resd0_1->GetFunction("fit")->GetParError(2));
626  error_d0.push_back(resd0_2->GetFunction("fit")->GetParError(2));
627  error_d0.push_back(resd0_3->GetFunction("fit")->GetParError(2));
628  error_d0.push_back(resd0_4->GetFunction("fit")->GetParError(2));
629  error_d0.push_back(resd0_5->GetFunction("fit")->GetParError(2));
630  error_d0.push_back(resd0_6->GetFunction("fit")->GetParError(2));
631 
632  // Fill the new histogram to create resolution plot
633  for(int i=0;i<6;i++) {
634  resd0->SetBinContent(i+1,sigma_d0[i]);
635  resd0->SetBinError(i+1,error_d0[i]);
636  }
637  }
638  else {
639  edm::LogWarning("DataNotFound")<<"L1 tracks not found for d0 resolution!\n";
640  for(int i=0;i<6;i++) {
641  resd0->SetBinContent(i+1,-1);
642  resd0->SetBinError(i+1,-1);
643  }
644  }
645  } //if ME found
646  else {edm::LogWarning("DataNotFound")<<"Monitor elements for d0 resolution cannot be found!\n";}
647 
648  } //if dbe found
649  else {edm::LogWarning("DataNotFound")<<"Cannot find valid DQM back end \n";}
650  delete fit;
651  delete fit2;
652  delete fit3;
653 } //end dqmEndJob
654 
TH1F * getTH1F() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
OuterTrackerMCHarvester(const edm::ParameterSet &)
HLT enums.
void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override