Belle II Software  release-08-01-08
DQMHistAnalysisEventT0.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 //+
9 // File : DQMHistAnalysisEventT0.cc
10 // Description : module for trigger jitter/EventT0 DQM histogram analysis
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisEventT0.h>
15 
16 #include <TROOT.h>
17 #include <TGraphAsymmErrors.h>
18 #include <TStyle.h>
19 #include <TF1.h>
20 #include <TMath.h>
21 
22 using namespace std;
23 using namespace Belle2;
24 
25 //-----------------------------------------------------------------
26 // Register the Module
27 //-----------------------------------------------------------------
28 REG_MODULE(DQMHistAnalysisEventT0);
29 
30 //-----------------------------------------------------------------
31 // Implementation
32 //-----------------------------------------------------------------
33 
34 DQMHistAnalysisEventT0Module::DQMHistAnalysisEventT0Module()
36 {
37  setDescription("Determining and processing EventT0s from different subdetector (TOP and SVD for now) for different L1 trigger sources (ECL and CDC for now) to estimate trigger jitter information.");
38 
39  //Parameter definition
40  addParam("min_nEntries", m_nEntriesMin, "Minimum number of entries to process the histogram.", m_nEntriesMin);
41  addParam("prefixCanvas", m_prefixCanvas, "Prefix to be added to canvas filename when saved as pdf.", std::string("c"));
42  addParam("printCanvas", m_printCanvas, "If true, prints pdf of the analysis canvas.", bool(false));
43 }
44 
45 
47 
49 {
50  gROOT->cd();
51 
52  //ECLTRG canvas
53  m_cTOPTimeForECLTRG = new TCanvas("EventT0/c_TOPTimeECLTRGjitter", "ECLTRG jitter based on TOP time", 1200, 400);
54  m_topPad1ECLTRG = new TPad("topPad1ECLTRG", "TOP pad1 ECLTRG", 0.03, 0.02, 0.33, 0.98);
55  m_topPad2ECLTRG = new TPad("topPad2ECLTRG", "TOP pad2 ECLTRG", 0.35, 0.02, 0.65, 0.98);
56  m_topPad3ECLTRG = new TPad("topPad3ECLTRG", "TOP pad3 ECLTRG", 0.67, 0.02, 0.97, 0.98);
57 
58  //CDCTRG canvases
59  m_cTOPTimeForCDCTRG = new TCanvas("EventT0/c_TOPTimeCDCTRGjitter", "CDCTRG jitter based on TOP time", 1200, 400);
60  m_topPad1CDCTRG = new TPad("topPad1CDCTRG", "TOP pad1 CDCTRG", 0.03, 0.02, 0.33, 0.98);
61  m_topPad2CDCTRG = new TPad("topPad2CDCTRG", "TOP pad2 CDCTRG", 0.35, 0.02, 0.65, 0.98);
62  m_topPad3CDCTRG = new TPad("topPad3CDCTRG", "TOP pad3 CDCTRG", 0.67, 0.02, 0.97, 0.98);
63 
64  //SVD canvases
65  //ECLTRG canvas
66  m_cSVDTimeForECLTRG = new TCanvas("EventT0/c_SVDTimeECLTRGjitter", "ECLTRG jitter based on SVD time", 1200, 400);
67  m_svdPad1ECLTRG = new TPad("svdPad1ECLTRG", "SVD pad1 ECLTRG", 0.03, 0.02, 0.33, 0.98);
68  m_svdPad2ECLTRG = new TPad("svdPad2ECLTRG", "SVD pad2 ECLTRG", 0.35, 0.02, 0.65, 0.98);
69  m_svdPad3ECLTRG = new TPad("svdPad3ECLTRG", "SVD pad3 ECLTRG", 0.67, 0.02, 0.97, 0.98);
70 
71  // CDC TRG canvas
72  m_cSVDTimeForCDCTRG = new TCanvas("EventT0/c_SVDTimeCDCTRGjitter", "CDCTRG jitter based on SVD time", 1200, 400);
73  m_svdPad1CDCTRG = new TPad("svdPad1CDCTRG", "SVD pad1 CDCTRG", 0.03, 0.02, 0.33, 0.98);
74  m_svdPad2CDCTRG = new TPad("svdPad2CDCTRG", "SVD pad2 CDCTRG", 0.35, 0.02, 0.65, 0.98);
75  m_svdPad3CDCTRG = new TPad("svdPad3CDCTRG", "SVD pad3 CDCTRG", 0.67, 0.02, 0.97, 0.98);
76 
77  // EventT0 source fractions
78  m_cT0FractionsForHadrons = new TCanvas("EventT0/c_T0FractionsForHadrons", "EventT0 source fractions for hadron events", 1200, 400);
79  m_pHadronECLTRG = new TPad("pHadronECLTRG", "Fractions ECLTRG", 0.03, 0.02, 0.33, 0.98);
80  m_pHadronCDCTRG = new TPad("pHadronCDCTRG", "Fractions CDCTRG", 0.35, 0.02, 0.65, 0.98);
81  m_pHadronTOPTRG = new TPad("pHadronTOPTRG", "Fractions TOPTRG", 0.67, 0.02, 0.97, 0.98);
82 
83  m_cT0FractionsForBhaBhas = new TCanvas("EventT0/c_cT0FractionsForBhaBhas", "EventT0 source fractions for bhabha events", 1200, 400);
84  m_pBhaBhaECLTRG = new TPad("pBhaBhaECLTRG", "Fractions ECLTRG", 0.03, 0.02, 0.33, 0.98);
85  m_pBhaBhaCDCTRG = new TPad("pBhaBhaCDCTRG", "Fractions CDCTRG", 0.35, 0.02, 0.65, 0.98);
86  m_pBhaBhaTOPTRG = new TPad("pBhaBhaTOPTRG", "Fractions TOPTRG", 0.67, 0.02, 0.97, 0.98);
87 
88  m_cT0FractionsForMuMus = new TCanvas("EventT0/c_cT0FractionsForMuMus", "EventT0 source fractions for #mu#mu events", 1200, 400);
89  m_pMuMuECLTRG = new TPad("pMuMuECLTRG", "Fractions ECLTRG", 0.03, 0.02, 0.33, 0.98);
90  m_pMuMuCDCTRG = new TPad("pMuMuCDCTRG", "Fractions CDCTRG", 0.35, 0.02, 0.65, 0.98);
91  m_pMuMuTOPTRG = new TPad("pMuMuTOPTRG", "Fractions TOPTRG", 0.67, 0.02, 0.97, 0.98);
92 
94  new TEfficiency("effAlgorithmSourceFractionsHadronL1ECLTRG",
95  "EventT0 source fractions, hadronic events, L1ECLTRG;Algorithm;Fraction #epsilon",
96  6, 0, 6);
98  new TEfficiency("effAlgorithmSourceFractionsHadronL1CDCTRG",
99  "EventT0 source fractions, hadronic events, L1CDCTRG;Algorithm;Fraction #epsilon",
100  6, 0, 6);
102  new TEfficiency("effAlgorithmSourceFractionsHadronL1TOPTRG",
103  "EventT0 source fractions, hadronic events, L1TOPTRG;Algorithm;Fraction #epsilon",
104  6, 0, 6);
106  new TEfficiency("effAlgorithmSourceFractionsBhaBhaL1ECLTRG",
107  "EventT0 source fractions, Bhabha events, L1ECLTRG;Algorithm;Fraction #epsilon",
108  6, 0, 6);
110  new TEfficiency("effAlgorithmSourceFractionsBhaBhaL1CDCTRG",
111  "EventT0 source fractions, Bhabha events, L1CDCTRG;Algorithm;Fraction #epsilon",
112  6, 0, 6);
114  new TEfficiency("effAlgorithmSourceFractionsBhaBhaL1TOPTRG",
115  "EventT0 source fractions, Bhabha events, L1TOPTRG;Algorithm;Fraction #epsilon",
116  6, 0, 6);
118  new TEfficiency("effAlgorithmSourceFractionsMuMuL1ECLTRG",
119  "EventT0 source fractions, #mu#mu events, L1ECLTRG;Algorithm;Fraction #epsilon",
120  6, 0, 6);
122  new TEfficiency("effAlgorithmSourceFractionsMuMuL1CDCTRG",
123  "EventT0 source fractions, #mu#mu events, L1CDCTRG;Algorithm;Fraction #epsilon",
124  6, 0, 6);
126  new TEfficiency("effAlgorithmSourceFractionsMuMuL1TOPTRG",
127  "EventT0 source fractions, #mu#mu events, L1TOPTRG;Algorithm;Fraction #epsilon",
128  6, 0, 6);
129 
130  m_monObj = getMonitoringObject("eventT0");
135 }
136 
137 
139 {
140  m_cTOPTimeForECLTRG->Clear();
141  m_cTOPTimeForCDCTRG->Clear();
142  m_cSVDTimeForECLTRG->Clear();
143  m_cSVDTimeForCDCTRG->Clear();
144 
145  m_cT0FractionsForHadrons->Clear();
146  m_cT0FractionsForBhaBhas->Clear();
147  m_cT0FractionsForMuMus->Clear();
148 }
149 
151 {
152  std::string histname = "AlgorithmSourceFractionsHadronL1ECLTRG";
153  m_pHadronECLTRG->cd();
155  m_pHadronECLTRG->SetFillColor(0);
156  m_pHadronECLTRG->Modified();
157  m_pHadronECLTRG->Update();
158  } else {
159  B2WARNING("Histogram EventT0 source fractions for hadrons from ECLTRG events (" << histname << ") from EventT0 DQM not processed!");
160  m_pHadronECLTRG->SetFillColor(kGray);
161  }
162 
163  histname = "AlgorithmSourceFractionsHadronL1CDCTRG";
164  m_pHadronCDCTRG->cd();
166  m_pHadronCDCTRG->SetFillColor(0);
167  m_pHadronCDCTRG->Modified();
168  m_pHadronCDCTRG->Update();
169  } else {
170  B2WARNING("Histogram EventT0 source fractions for hadrons from CDCTRG events (" << histname << ") from EventT0 DQM not processed!");
171  m_pHadronCDCTRG->SetFillColor(kGray);
172  }
173 
174  histname = "AlgorithmSourceFractionsHadronL1TOPTRG";
175  m_pHadronTOPTRG->cd();
177  m_pHadronTOPTRG->SetFillColor(0);
178  m_pHadronTOPTRG->Modified();
179  m_pHadronTOPTRG->Update();
180  } else {
181  B2WARNING("Histogram EventT0 source fractions for hadrons from TOPTRG events (" << histname << ") from EventT0 DQM not processed!");
182  m_pHadronTOPTRG->SetFillColor(kGray);
183  }
184 
185  // anonymous helper to draw the pads consistently without repeating oneself
186  const auto drawPad = [](TPad * pad) {
187  pad->SetMargin(0.07, 0.0, 0.07, 0.07);
188  pad->Draw();
189  };
191  drawPad(m_pHadronECLTRG);
192  drawPad(m_pHadronCDCTRG);
193  drawPad(m_pHadronTOPTRG);
194 
195  if (m_printCanvas)
196  m_cT0FractionsForHadrons->Print(Form("%s_EventT0FractionsHadron.pdf", m_prefixCanvas.c_str()));
197 
198 
199  histname = "AlgorithmSourceFractionsBhaBhaL1ECLTRG";
200  m_pBhaBhaECLTRG->cd();
202  m_pBhaBhaECLTRG->SetFillColor(0);
203  m_pBhaBhaECLTRG->Modified();
204  m_pBhaBhaECLTRG->Update();
205  } else {
206  B2WARNING("Histogram EventT0 source fractions for BhaBha from ECLTRG events (" << histname << ") from EventT0 DQM not processed!");
207  m_pBhaBhaECLTRG->SetFillColor(kGray);
208  }
209 
210  histname = "AlgorithmSourceFractionsBhaBhaL1CDCTRG";
211  m_pBhaBhaCDCTRG->cd();
213  m_pBhaBhaCDCTRG->SetFillColor(0);
214  m_pBhaBhaCDCTRG->Modified();
215  m_pBhaBhaCDCTRG->Update();
216  } else {
217  B2WARNING("Histogram EventT0 source fractions for BhaBha from CDCTRG events (" << histname << ") from EventT0 DQM not processed!");
218  m_pBhaBhaCDCTRG->SetFillColor(kGray);
219  }
220 
221  histname = "AlgorithmSourceFractionsBhaBhaL1TOPTRG";
222  m_pBhaBhaTOPTRG->cd();
224  m_pBhaBhaTOPTRG->SetFillColor(0);
225  m_pBhaBhaTOPTRG->Modified();
226  m_pBhaBhaTOPTRG->Update();
227  } else {
228  B2WARNING("Histogram EventT0 source fractions for BhaBha from TOPTRG events (" << histname << ") from EventT0 DQM not processed!");
229  m_pBhaBhaTOPTRG->SetFillColor(kGray);
230  }
231 
233  drawPad(m_pBhaBhaECLTRG);
234  drawPad(m_pBhaBhaCDCTRG);
235  drawPad(m_pBhaBhaTOPTRG);
236 
237  if (m_printCanvas)
238  m_cT0FractionsForBhaBhas->Print(Form("%s_EventT0FractionsBhaBha.pdf", m_prefixCanvas.c_str()));
239 
240 
241  histname = "AlgorithmSourceFractionsMuMuL1ECLTRG";
242  m_pMuMuECLTRG->cd();
244  m_pMuMuECLTRG->SetFillColor(0);
245  m_pMuMuECLTRG->Modified();
246  m_pMuMuECLTRG->Update();
247  } else {
248  B2WARNING("Histogram EventT0 source fractions for MuMu from ECLTRG events (" << histname << ") from EventT0 DQM not processed!");
249  m_pMuMuECLTRG->SetFillColor(kGray);
250  }
251 
252  histname = "AlgorithmSourceFractionsMuMuL1CDCTRG";
253  m_pMuMuCDCTRG->cd();
255  m_pMuMuCDCTRG->SetFillColor(0);
256  m_pMuMuCDCTRG->Modified();
257  m_pMuMuCDCTRG->Update();
258  } else {
259  B2WARNING("Histogram EventT0 source fractions for MuMu from CDCTRG events (" << histname << ") from EventT0 DQM not processed!");
260  m_pMuMuCDCTRG->SetFillColor(kGray);
261  }
262 
263  histname = "AlgorithmSourceFractionsMuMuL1TOPTRG";
264  m_pMuMuTOPTRG->cd();
266  m_pMuMuTOPTRG->SetFillColor(0);
267  m_pMuMuTOPTRG->Modified();
268  m_pMuMuTOPTRG->Update();
269  } else {
270  B2WARNING("Histogram EventT0 source fractions for MuMu from TOPTRG events (" << histname << ") from EventT0 DQM not processed!");
271  m_pMuMuTOPTRG->SetFillColor(kGray);
272  }
273 
275  drawPad(m_pMuMuECLTRG);
276  drawPad(m_pMuMuCDCTRG);
277  drawPad(m_pMuMuTOPTRG);
278 
279  if (m_printCanvas)
280  m_cT0FractionsForMuMus->Print(Form("%s_EventT0FractionsMuMu.pdf", m_prefixCanvas.c_str()));
281 
282 }
283 
285 {
286  // --- TOP EventT0 plots for ECLTRG ---
287 
288  // find TOP EventT0 Hadrons ECLTRG histogram and process it
289  TH1* h = findHist("EventT0/m_histEventT0_TOP_hadron_L1_ECLTRG");
290  TString tag = "hadronECLTRG";
291  m_topPad1ECLTRG->cd();
292  if (processHistogram(h, tag)) {
293  m_topPad1ECLTRG->SetFillColor(0);
294  m_topPad1ECLTRG->Modified();
295  m_topPad1ECLTRG->Update();
296  } else {
297  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
298  h->Draw();
299  m_topPad1ECLTRG->SetFillColor(kGray);
300  }
301 
302  // find TOP EventT0 Bhabhas ECLTRG histogram and process it
303  h = findHist("EventT0/m_histEventT0_TOP_bhabha_L1_ECLTRG");
304  tag = "bhabhaECLTRG";
305  m_topPad2ECLTRG->cd();
306  if (processHistogram(h, tag)) {
307  m_topPad2ECLTRG->SetFillColor(0);
308  m_topPad2ECLTRG->Modified();
309  m_topPad2ECLTRG->Update();
310  } else {
311  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
312  h->Draw();
313  m_topPad2ECLTRG->SetFillColor(kGray);
314  }
315 
316  // find TOP EventT0 Mumus ECLTRG histogram and process it
317  h = findHist("EventT0/m_histEventT0_TOP_mumu_L1_ECLTRG");
318  tag = "mumuECLTRG";
319  m_topPad3ECLTRG->cd();
320  if (processHistogram(h, tag)) {
321  m_topPad3ECLTRG->SetFillColor(0);
322  m_topPad3ECLTRG->Modified();
323  m_topPad3ECLTRG->Update();
324  } else {
325  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
326  h->Draw();
327  m_topPad3ECLTRG->SetFillColor(kGray);
328  }
329 
330  m_cTOPTimeForECLTRG->cd();
331  m_topPad1ECLTRG->Draw();
332  m_topPad2ECLTRG->Draw();
333  m_topPad3ECLTRG->Draw();
334 
335  if (m_printCanvas)
336  m_cTOPTimeForECLTRG->Print(Form("%s_ECLTRG.pdf", m_prefixCanvas.c_str()));
337 
338 
339  // --- TOP EventT0 plots for CDCTRG ---
340 
341  // find TOP EventT0 Hadrons CDCTRG histogram and process it
342  h = findHist("EventT0/m_histEventT0_TOP_hadron_L1_CDCTRG");
343  tag = "hadronCDCTRG";
344  m_topPad1CDCTRG->cd();
345  if (processHistogram(h, tag)) {
346  m_topPad1CDCTRG->SetFillColor(0);
347  m_topPad1CDCTRG->Modified();
348  m_topPad1CDCTRG->Update();
349  m_cTOPTimeForCDCTRG->cd();
350  m_topPad1CDCTRG->Draw();
351  } else {
352  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
353  h->Draw();
354  m_topPad1CDCTRG->SetFillColor(kGray);
355  m_cTOPTimeForCDCTRG->cd();
356  m_topPad1CDCTRG->Draw();
357  }
358 
359  // find TOP EventT0 Bhabhas CDCTRG histogram and process it
360  h = findHist("EventT0/m_histEventT0_TOP_bhabha_L1_CDCTRG");
361  tag = "bhabhaCDCTRG";
362  m_topPad2CDCTRG->cd();
363  if (processHistogram(h, tag)) {
364  m_topPad2CDCTRG->SetFillColor(0);
365  m_topPad2CDCTRG->Modified();
366  m_topPad2CDCTRG->Update();
367  m_cTOPTimeForCDCTRG->cd();
368  m_topPad2CDCTRG->Draw();
369  } else {
370  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
371  h->Draw();
372  m_topPad2CDCTRG->SetFillColor(kGray);
373  m_cTOPTimeForCDCTRG->cd();
374  m_topPad2CDCTRG->Draw();
375  }
376 
377  // find TOP EventT0 Mumus CDCTRG histogram and process it
378  h = findHist("EventT0/m_histEventT0_TOP_mumu_L1_CDCTRG");
379  tag = "mumuCDCTRG";
380  m_topPad3CDCTRG->cd();
381  if (processHistogram(h, tag)) {
382  m_topPad3CDCTRG->SetFillColor(0);
383  m_topPad3CDCTRG->Modified();
384  m_topPad3CDCTRG->Update();
385  } else {
386  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
387  h->Draw();
388  m_topPad3CDCTRG->SetFillColor(kGray);
389  }
390 
391  m_cTOPTimeForCDCTRG->cd();
392  m_topPad1CDCTRG->Draw();
393  m_topPad2CDCTRG->Draw();
394  m_topPad3CDCTRG->Draw();
395 
396  if (m_printCanvas)
397  m_cTOPTimeForCDCTRG->Print(Form("%s_CDCTRG.pdf", m_prefixCanvas.c_str()));
398 
399 
400 
401  // --- SVD EventT0 plots for ECLTRG ---
402 
403  // find SVD EventT0 Hadrons ECLTRG histogram and process it
404  h = findHist("EventT0/m_histEventT0_SVD_hadron_L1_ECLTRG");
405  tag = "hadronECLTRG";
406  m_svdPad1ECLTRG->cd();
407  if (processHistogram(h, tag)) {
408  m_svdPad1ECLTRG->SetFillColor(0);
409  m_svdPad1ECLTRG->Modified();
410  m_svdPad1ECLTRG->Update();
411  } else {
412  B2WARNING(Form("Histogram SVD EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
413  h->Draw();
414  m_svdPad1ECLTRG->SetFillColor(kGray);
415  }
416 
417  // find SVD EventT0 Bhabhas ECLTRG histogram and process it
418  h = findHist("EventT0/m_histEventT0_SVD_bhabha_L1_ECLTRG");
419  tag = "bhabhaECLTRG";
420  m_svdPad2ECLTRG->cd();
421  if (processHistogram(h, tag)) {
422  m_svdPad2ECLTRG->SetFillColor(0);
423  m_svdPad2ECLTRG->Modified();
424  m_svdPad2ECLTRG->Update();
425  } else {
426  B2WARNING(Form("Histogram SVD EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
427  h->Draw();
428  m_svdPad2ECLTRG->SetFillColor(kGray);
429  }
430 
431  // find SVD EventT0 Mumus ECLTRG histogram and process it
432  h = findHist("EventT0/m_histEventT0_SVD_mumu_L1_ECLTRG");
433  tag = "mumuECLTRG";
434  m_svdPad3ECLTRG->cd();
435  if (processHistogram(h, tag)) {
436  m_svdPad3ECLTRG->SetFillColor(0);
437  m_svdPad3ECLTRG->Modified();
438  m_svdPad3ECLTRG->Update();
439  } else {
440  B2WARNING(Form("Histogram SVD EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
441  h->Draw();
442  m_svdPad3ECLTRG->SetFillColor(kGray);
443  }
444 
445  m_cSVDTimeForECLTRG->cd();
446  m_svdPad1ECLTRG->Draw();
447  m_svdPad2ECLTRG->Draw();
448  m_svdPad3ECLTRG->Draw();
449 
450  if (m_printCanvas)
451  m_cSVDTimeForECLTRG->Print(Form("%s_SVDECLTRG.pdf", m_prefixCanvas.c_str()));
452 
453 
454  // --- SVD EventT0 plots for CDCTRG ---
455 
456  // find SVD EventT0 Hadrons CDCTRG histogram and process it
457  h = findHist("EventT0/m_histEventT0_SVD_hadron_L1_CDCTRG");
458  tag = "hadronCDCTRG";
459  m_svdPad1CDCTRG->cd();
460  if (processHistogram(h, tag)) {
461  m_svdPad1CDCTRG->SetFillColor(0);
462  m_svdPad1CDCTRG->Modified();
463  m_svdPad1CDCTRG->Update();
464  m_cSVDTimeForCDCTRG->cd();
465  m_svdPad1CDCTRG->Draw();
466  } else {
467  B2WARNING(Form("Histogram SVD EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
468  h->Draw();
469  m_svdPad1CDCTRG->SetFillColor(kGray);
470  m_cSVDTimeForCDCTRG->cd();
471  m_svdPad1CDCTRG->Draw();
472  }
473 
474  // find SVD EventT0 Bhabhas CDCTRG histogram and process it
475  h = findHist("EventT0/m_histEventT0_SVD_bhabha_L1_CDCTRG");
476  tag = "bhabhaCDCTRG";
477  m_svdPad2CDCTRG->cd();
478  if (processHistogram(h, tag)) {
479  m_svdPad2CDCTRG->SetFillColor(0);
480  m_svdPad2CDCTRG->Modified();
481  m_svdPad2CDCTRG->Update();
482  m_cSVDTimeForCDCTRG->cd();
483  m_svdPad2CDCTRG->Draw();
484  } else {
485  B2WARNING(Form("Histogram SVD EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
486  h->Draw();
487  m_svdPad2CDCTRG->SetFillColor(kGray);
488  m_cSVDTimeForCDCTRG->cd();
489  m_svdPad2CDCTRG->Draw();
490  }
491 
492 
493  // find SVD EventT0 Mumus CDCTRG histogram and process it
494  h = findHist("EventT0/m_histEventT0_SVD_mumu_L1_CDCTRG");
495  tag = "mumuCDCTRG";
496  m_svdPad3CDCTRG->cd();
497  if (processHistogram(h, tag)) {
498  m_svdPad3CDCTRG->SetFillColor(0);
499  m_svdPad3CDCTRG->Modified();
500  m_svdPad3CDCTRG->Update();
501  } else {
502  B2WARNING(Form("Histogram SVD EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
503  h->Draw();
504  m_svdPad3CDCTRG->SetFillColor(kGray);
505  }
506 
507  m_cSVDTimeForCDCTRG->cd();
508  m_svdPad1CDCTRG->Draw();
509  m_svdPad2CDCTRG->Draw();
510  m_svdPad3CDCTRG->Draw();
511 
512  if (m_printCanvas)
513  m_cSVDTimeForCDCTRG->Print(Form("%s_SVDCDCTRG.pdf", m_prefixCanvas.c_str()));
514 
515 }
516 
518 {
519 
520  delete m_cTOPTimeForECLTRG;
521  delete m_cTOPTimeForCDCTRG;
522  delete m_cSVDTimeForECLTRG;
523  delete m_cSVDTimeForCDCTRG;
524 
527  delete m_cT0FractionsForMuMus;
528 
529 }
530 
531 double DQMHistAnalysisEventT0Module::fDoubleGaus(double* x, double* par)
532 {
533  double N = par[0];
534  double frac = par[1];
535  double mean = par[2];
536  double sigma = par[3];
537  double mean2 = par[4];
538  double sigma2 = par[5];
539 
540  return N * frac * TMath::Gaus(x[0], mean, sigma) + N * (1 - frac) * TMath::Gaus(x[0], mean2, sigma2);
541 }
542 
544 {
545 
546  if (h == nullptr) {
547  B2DEBUG(20, "h == nullptr");
548  m_monObj->setVariable(Form("fit_%s", tag.Data()), 0);
549  return false;
550  }
551 
552  // The default value for the EventT0 value is -1000, but bins start at -100, so we might mostly fill the underflow bin if
553  // EventT0 for a detector is not present. And also the nominal EventT0 might be too big or too small. Only use the content
554  // of the actually useful bins to decide whether or not to fit the histogram.
555  auto nValidEntries = h->GetEntries() - h->GetBinContent(0) - h->GetBinContent(h->GetNbinsX() + 1);
556  if (static_cast<uint>(nValidEntries) < m_nEntriesMin) {
557  B2DEBUG(20, "not enough entries");
558  m_monObj->setVariable(Form("fit_%s", tag.Data()), 0);
559  return false;
560  }
561 
562 
563  //scale the histogram only with content of valid bins, ignore over and underflow bins
564  h->Scale(1. / nValidEntries);
565  h->GetXaxis()->SetRangeUser(-50, 50);
566 
567  //define the fitting function
568  TF1 fitf("fit", DQMHistAnalysisEventT0Module::fDoubleGaus, -50, 50, 6);
569  fitf.SetParNames("N", "f_{1}", "#mu_{1}", "#sigma_{1}", "#mu_{2}", "#sigma_{2}");
570  fitf.SetParameters(0.1, 0.8, 0, 5, 0, 15);
571  fitf.SetParLimits(1, 0, 1); //fraction
572  fitf.SetParLimits(3, 0, 100); //sigma1
573  fitf.SetParLimits(5, 0, 100); //sigma2
574 
575  if (h->Fit(&fitf, "SR+") != 0) {
576  B2DEBUG(20, "failed fit");
577  m_monObj->setVariable(Form("fit_%s", tag.Data()), 0);
578  return false;
579  }
580 
581  Double_t par[6];
582  fitf.GetParameters(&par[0]);
583  Double_t parErr[6];
584  for (int i = 0; i < 6; i++)
585  parErr[i] = fitf.GetParError(i) ;
586 
587 
588  //define gaussian components
589  TF1 gauss1("gauss1", "gaus", -100, 100);
590  TF1 gauss2("gauss2", "gaus", -100, 100);
591 
592  gauss1.SetLineColor(kBlue);
593  gauss1.SetLineStyle(kDashed);
594  gauss1.SetParameters(par[0]*par[1], par[2], par[3]);
595 
596  gauss2.SetLineColor(kRed);
597  gauss2.SetLineStyle(kDashed);
598  gauss2.SetParameters(par[0] * (1 - par[1]), par[4], par[5]);
599 
600  m_monObj->setVariable(Form("fit_%s", tag.Data()), 1);
601  m_monObj->setVariable(Form("N_%s", tag.Data()), nValidEntries, TMath::Sqrt(nValidEntries));
602  m_monObj->setVariable(Form("f_%s", tag.Data()), par[1], parErr[1]);
603  m_monObj->setVariable(Form("mean1_%s", tag.Data()), par[2], parErr[2]);
604  m_monObj->setVariable(Form("sigma1_%s", tag.Data()), par[3], parErr[3]);
605  m_monObj->setVariable(Form("mean2_%s", tag.Data()), par[4], parErr[4]);
606  m_monObj->setVariable(Form("sigma2_%s", tag.Data()), par[5], parErr[5]);
607 
608  //SETUP gSTYLE - all plots
609  gStyle->SetOptFit(1111);
610 
611  h->DrawClone();
612  fitf.DrawClone("same");
613  gauss1.DrawClone("same");
614  gauss2.DrawClone("same");
615 
616  return true;
617 
618 }
619 
620 bool DQMHistAnalysisEventT0Module::FillEfficiencyHistogram(const std::string& histname, TEfficiency* eff)
621 {
622  B2DEBUG(20, "Begin processing histogram " << histname << " ...");
623  TH1* h = findHist("EventT0/" + histname);
624  if (not h) {
625  return false;
626  }
627 
628  // Admittedly quite a hacky way to obtain the normalisation values: Create a new histogram and fill each of the bins with
629  // the bin content of the -1 bin of h which is used for bin counting, and at the same time set the corresponding bin label.
630  const auto totalEntries = h->GetBinContent(-1);
631  const auto nBins = h->GetNbinsX();
632  TH1D* totalHist = new TH1D("total", "total;Algorithm;Fraction #epsilon", nBins, 0, nBins);
633  for (int i = 0; i < nBins; i++) {
634  totalHist->SetBinContent(i + 1, totalEntries);
635  }
636  eff->SetPassedHistogram(*h, "f");
637  eff->SetTotalHistogram(*totalHist, "f");
638 
639  eff->Paint("AP");
640 
641  TGraphAsymmErrors* graph = eff->GetPaintedGraph();
642  if (not graph) {
643  return false;
644  }
645 
646  auto ax = graph->GetXaxis();
647  if (not ax) {
648  return false;
649  }
650  // Print x-axis bin labels horizontally
651  ax->SetTitleOffset(1.0);
652  ax->CenterTitle(kTRUE);
653  ax->Set(nBins, 0, nBins);
654  for (int i = 0; i < nBins; i++) {
655  ax->SetBinLabel(i + 1, c_eventT0Algorithms[i]);
656  }
657 
658  auto ay = graph->GetYaxis();
659  if (not ay) {
660  return false;
661  }
662  ay->SetTitleOffset(1.0);
663  ay->SetRangeUser(0, 1.05);
664 
665  graph->Draw("AP");
666 
667  B2DEBUG(20, "Finished processing histogram " << histname << "!");
668 
669  return true;
670 }
671 
TCanvas * m_cTOPTimeForECLTRG
TOP EventT0 for ECLTRG plots canvas.
TEfficiency * m_eAlgorithmSourceFractionsMuMuL1ECLTRG
Fraction of events with EventT0 from a given algorithm, HLT mumu events, L1 time by ECL trigger.
TCanvas * m_cSVDTimeForECLTRG
SVD EventT0 for ECLTRG plots canvas.
void initialize() override final
create TCanvas and MonitoringObject
TEfficiency * m_eAlgorithmSourceFractionsBhaBhaL1ECLTRG
Fraction of events with EventT0 from a given algorithm, HLT bhabha events, L1 time by ECL trigger.
TCanvas * m_cT0FractionsForMuMus
EventT0 fractions plots canvas for MuMu events.
static double fDoubleGaus(double *x, double *par)
double gaussian fitting function for the jitter distribution
const char * c_eventT0Algorithms[6]
EventT0 algorithms for which to calculate fractions of abundance.
TPad * m_pBhaBhaTOPTRG
pad for time fractions for TOPTRG bhabhas
TPad * m_pMuMuTOPTRG
pad for time fractions for TOPTRG mumu
bool FillEfficiencyHistogram(const std::string &histname, TEfficiency *eff)
Fill the TEfficiency plots.
TPad * m_svdPad1CDCTRG
pad for SVD time CDCTRG hadrons
TPad * m_svdPad3CDCTRG
pad for SVD time CDCTRG mumu
TPad * m_pMuMuCDCTRG
pad for time fractions for CDCTRG mumu
TPad * m_pHadronTOPTRG
pad for time fractions for TOPTRG hadrons
TCanvas * m_cT0FractionsForBhaBhas
EventT0 fractions plots canvas for BhaBha events.
std::string m_prefixCanvas
prefix to be added to canvas name when saved as pdf
TPad * m_pMuMuECLTRG
pad for time fractions for ECLTRG mumu
TPad * m_svdPad2ECLTRG
pad for SVD time ECLTRG bhabhas
MonitoringObject * m_monObj
MonitoringObject to be produced by this module.
void terminate() override final
delete pointers
TEfficiency * m_eAlgorithmSourceFractionsBhaBhaL1CDCTRG
Fraction of events with EventT0 from a given algorithm, HLT bhabha events, L1 time by CDC trigger.
TCanvas * m_cT0FractionsForHadrons
EventT0 fractions plots canvas for hadron events.
bool m_printCanvas
if true print the pdf of the canvases
TPad * m_svdPad2CDCTRG
pad for SVD time CDCTRG bhabhas
TPad * m_pHadronECLTRG
pad for time fractions for ECLTRG hadrons
TPad * m_topPad1CDCTRG
pad for TOP time CDCTRG hadrons
uint m_nEntriesMin
minimum number of entries to process the histogram
TEfficiency * m_eAlgorithmSourceFractionsHadronL1ECLTRG
Fraction of events with EventT0 from a given algorithm, HLT hadronic events, L1 time by ECL trigger.
void endRun() override final
fit the histograms
TPad * m_topPad3CDCTRG
pad for TOP time CDCTRG mumu
TPad * m_svdPad1ECLTRG
pad for SVD time ECLTRG hadrons
TEfficiency * m_eAlgorithmSourceFractionsHadronL1CDCTRG
Fraction of events with EventT0 from a given algorithm, HLT hadronic events, L1 time by CDC trigger.
bool processHistogram(TH1 *h, TString tag)
process the EventT0 distribution fitting with two gaussians filling the MonitoringObject
void beginRun() override final
clear TCanvas
TPad * m_topPad2CDCTRG
pad for TOP time CDCTRG bhabhas
TPad * m_pHadronCDCTRG
pad for time fractions for CDCTRG hadrons
TPad * m_topPad3ECLTRG
pad for TOP time ECLTRG mumu
TPad * m_topPad2ECLTRG
pad for TOP time ECLTRG bhabhas
TEfficiency * m_eAlgorithmSourceFractionsMuMuL1TOPTRG
Fraction of events with EventT0 from a given algorithm, HLT mumu events, L1 time by TOP trigger.
TCanvas * m_cTOPTimeForCDCTRG
TOP EventT0 for CDCTRG plots canvas.
TCanvas * m_cSVDTimeForCDCTRG
SVD EventT0 for CDCTRG plots canvas.
TPad * m_topPad1ECLTRG
pad for TOP time ECLTRG hadrons
TEfficiency * m_eAlgorithmSourceFractionsMuMuL1CDCTRG
Fraction of events with EventT0 from a given algorithm, HLT mumu events, L1 time by CDC trigger.
TEfficiency * m_eAlgorithmSourceFractionsBhaBhaL1TOPTRG
Fraction of events with EventT0 from a given algorithm, HLT bhabha events, L1 time by TOP trigger.
TPad * m_svdPad3ECLTRG
pad for SVD time ECLTRG mumu
TPad * m_pBhaBhaCDCTRG
pad for time fractions for CDCTRG bhabhas
TEfficiency * m_eAlgorithmSourceFractionsHadronL1TOPTRG
Fraction of events with EventT0 from a given algorithm, HLT hadronic events, L1 time by TOP trigger.
TPad * m_pBhaBhaECLTRG
pad for time fractions for ECLTRG bhabhas
The base class for the histogram analysis module.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
void addCanvas(TCanvas *canv)
Add Canvas to monitoring object.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.