Belle II Software  release-08-01-10
DQMHistAnalysisEventT0Efficiency.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 : DQMHistAnalysisEventT0Efficiency.h
10 // Description : module for DQM histogram analysis of EventT0 algorithm efficiencies
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisEventT0Efficiency.h>
15 
16 #include <TROOT.h>
17 #include <TGraphAsymmErrors.h>
18 #include <TStyle.h>
19 
20 using namespace Belle2;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(DQMHistAnalysisEventT0Efficiency);
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
33 {
34  setDescription("Calculate EventT0 algorithm efficiencies for the different trigger types, EventT0 algorithm sources, and event types (hadron, bhabha, ยตยต).");
35 
36  //Parameter definition
37  addParam("min_nEntries", m_nEntriesMin, "Minimum number of entries to process the histogram.", m_nEntriesMin);
38  addParam("prefixCanvas", m_prefixCanvas, "Prefix to be added to canvas filename when saved as pdf.", std::string("c"));
39  addParam("printCanvas", m_printCanvas, "If true, prints pdf of the analysis canvas.", bool(false));
40 }
41 
42 
44 
46 {
47  gROOT->cd();
48 
49  // EventT0 source fractions
50  m_cT0FractionsHadronECLTRG = new TCanvas("EventT0/c_HadronECLTRG", "Fractions ECLTRG for hadron events");
51  m_cT0FractionsHadronCDCTRG = new TCanvas("EventT0/c_HadronCDCTRG", "Fractions CDCTRG for hadron events");
52  m_cT0FractionsHadronTOPTRG = new TCanvas("EventT0/c_HadronTOPTRG", "Fractions TOPTRG for hadron events");
53 
54  m_cT0FractionsBhaBhaECLTRG = new TCanvas("EventT0/c_BhaBhaECLTRG", "Fractions ECLTRG for BhaBha events");
55  m_cT0FractionsBhaBhaCDCTRG = new TCanvas("EventT0/c_BhaBhaCDCTRG", "Fractions CDCTRG for BhaBha events");
56  m_cT0FractionsBhaBhaTOPTRG = new TCanvas("EventT0/c_BhaBhaTOPTRG", "Fractions TOPTRG for BhaBha events");
57 
58  m_cT0FractionsMuMuECLTRG = new TCanvas("EventT0/c_MuMuECLTRG", "Fractions ECLTRG for #mu#mu events");
59  m_cT0FractionsMuMuCDCTRG = new TCanvas("EventT0/c_MuMuCDCTRG", "Fractions CDCTRG for #mu#mu events");
60  m_cT0FractionsMuMuTOPTRG = new TCanvas("EventT0/c_MuMuTOPTRG", "Fractions TOPTRG for #mu#mu events");
61 
63  new TEfficiency("effAlgorithmSourceFractionsHadronL1ECLTRG",
64  "EventT0 source fractions, hadronic events, L1TRG timing from ECL;Algorithm;Fraction #epsilon",
65  6, 0, 6);
67  new TEfficiency("effAlgorithmSourceFractionsHadronL1CDCTRG",
68  "EventT0 source fractions, hadronic events, L1TRG timing from CDC;Algorithm;Fraction #epsilon",
69  6, 0, 6);
71  new TEfficiency("effAlgorithmSourceFractionsHadronL1TOPTRG",
72  "EventT0 source fractions, hadronic events, L1TRG timing from TOP;Algorithm;Fraction #epsilon",
73  6, 0, 6);
75  new TEfficiency("effAlgorithmSourceFractionsBhaBhaL1ECLTRG",
76  "EventT0 source fractions, Bhabha events, L1TRG timing from ECL;Algorithm;Fraction #epsilon",
77  6, 0, 6);
79  new TEfficiency("effAlgorithmSourceFractionsBhaBhaL1CDCTRG",
80  "EventT0 source fractions, Bhabha events, L1TRG timing from CDC;Algorithm;Fraction #epsilon",
81  6, 0, 6);
83  new TEfficiency("effAlgorithmSourceFractionsBhaBhaL1TOPTRG",
84  "EventT0 source fractions, Bhabha events, L1TRG timing from TOP;Algorithm;Fraction #epsilon",
85  6, 0, 6);
87  new TEfficiency("effAlgorithmSourceFractionsMuMuL1ECLTRG",
88  "EventT0 source fractions, #mu#mu events, L1TRG timing from ECL;Algorithm;Fraction #epsilon",
89  6, 0, 6);
91  new TEfficiency("effAlgorithmSourceFractionsMuMuL1CDCTRG",
92  "EventT0 source fractions, #mu#mu events, L1TRG timing from CDC;Algorithm;Fraction #epsilon",
93  6, 0, 6);
95  new TEfficiency("effAlgorithmSourceFractionsMuMuL1TOPTRG",
96  "EventT0 source fractions, #mu#mu events, L1TRG timing from TOP;Algorithm;Fraction #epsilon",
97  6, 0, 6);
98 
99  m_monObj = getMonitoringObject("eventT0");
100 }
101 
102 
104 {
111  m_cT0FractionsMuMuECLTRG->Clear();
112  m_cT0FractionsMuMuCDCTRG->Clear();
113  m_cT0FractionsMuMuTOPTRG->Clear();
114 }
115 
117 {
118  std::string histname = "AlgorithmSourceFractionsHadronL1ECLTRG";
121  m_cT0FractionsHadronECLTRG->SetFillColor(0);
122  m_cT0FractionsHadronECLTRG->Modified();
123  m_cT0FractionsHadronECLTRG->Update();
124  } else {
125  B2DEBUG(29, "Histogram EventT0 source fractions for hadrons from ECLTRG events (" << histname <<
126  ") from EventT0 DQM not processed!");
127  m_cT0FractionsHadronECLTRG->SetFillColor(kGray);
128  }
129 
130  histname = "AlgorithmSourceFractionsHadronL1CDCTRG";
133  m_cT0FractionsHadronCDCTRG->SetFillColor(0);
134  m_cT0FractionsHadronCDCTRG->Modified();
135  m_cT0FractionsHadronCDCTRG->Update();
136  } else {
137  B2DEBUG(29, "Histogram EventT0 source fractions for hadrons from CDCTRG events (" << histname <<
138  ") from EventT0 DQM not processed!");
139  m_cT0FractionsHadronCDCTRG->SetFillColor(kGray);
140  }
141 
142  histname = "AlgorithmSourceFractionsHadronL1TOPTRG";
145  m_cT0FractionsHadronTOPTRG->SetFillColor(0);
146  m_cT0FractionsHadronTOPTRG->Modified();
147  m_cT0FractionsHadronTOPTRG->Update();
148  } else {
149  B2DEBUG(29, "Histogram EventT0 source fractions for hadrons from TOPTRG events (" << histname <<
150  ") from EventT0 DQM not processed!");
151  m_cT0FractionsHadronTOPTRG->SetFillColor(kGray);
152  }
153 
154 
155  histname = "AlgorithmSourceFractionsBhaBhaL1ECLTRG";
158  m_cT0FractionsBhaBhaECLTRG->SetFillColor(0);
159  m_cT0FractionsBhaBhaECLTRG->Modified();
160  m_cT0FractionsBhaBhaECLTRG->Update();
161  } else {
162  B2DEBUG(29, "Histogram EventT0 source fractions for BhaBha from ECLTRG events (" << histname <<
163  ") from EventT0 DQM not processed!");
164  m_cT0FractionsBhaBhaECLTRG->SetFillColor(kGray);
165  }
166 
167  histname = "AlgorithmSourceFractionsBhaBhaL1CDCTRG";
170  m_cT0FractionsBhaBhaCDCTRG->SetFillColor(0);
171  m_cT0FractionsBhaBhaCDCTRG->Modified();
172  m_cT0FractionsBhaBhaCDCTRG->Update();
173  } else {
174  B2DEBUG(29, "Histogram EventT0 source fractions for BhaBha from CDCTRG events (" << histname <<
175  ") from EventT0 DQM not processed!");
176  m_cT0FractionsBhaBhaCDCTRG->SetFillColor(kGray);
177  }
178 
179  histname = "AlgorithmSourceFractionsBhaBhaL1TOPTRG";
182  m_cT0FractionsBhaBhaTOPTRG->SetFillColor(0);
183  m_cT0FractionsBhaBhaTOPTRG->Modified();
184  m_cT0FractionsBhaBhaTOPTRG->Update();
185  } else {
186  B2DEBUG(29, "Histogram EventT0 source fractions for BhaBha from TOPTRG events (" << histname <<
187  ") from EventT0 DQM not processed!");
188  m_cT0FractionsBhaBhaTOPTRG->SetFillColor(kGray);
189  }
190 
191 
192  histname = "AlgorithmSourceFractionsMuMuL1ECLTRG";
195  m_cT0FractionsMuMuECLTRG->SetFillColor(0);
196  m_cT0FractionsMuMuECLTRG->Modified();
197  m_cT0FractionsMuMuECLTRG->Update();
198  } else {
199  B2DEBUG(29, "Histogram EventT0 source fractions for MuMu from ECLTRG events (" << histname << ") from EventT0 DQM not processed!");
200  m_cT0FractionsMuMuECLTRG->SetFillColor(kGray);
201  }
202 
203  histname = "AlgorithmSourceFractionsMuMuL1CDCTRG";
206  m_cT0FractionsMuMuCDCTRG->SetFillColor(0);
207  m_cT0FractionsMuMuCDCTRG->Modified();
208  m_cT0FractionsMuMuCDCTRG->Update();
209  } else {
210  B2DEBUG(29, "Histogram EventT0 source fractions for MuMu from CDCTRG events (" << histname << ") from EventT0 DQM not processed!");
211  m_cT0FractionsMuMuCDCTRG->SetFillColor(kGray);
212  }
213 
214  histname = "AlgorithmSourceFractionsMuMuL1TOPTRG";
217  m_cT0FractionsMuMuTOPTRG->SetFillColor(0);
218  m_cT0FractionsMuMuTOPTRG->Modified();
219  m_cT0FractionsMuMuTOPTRG->Update();
220  } else {
221  B2DEBUG(29, "Histogram EventT0 source fractions for MuMu from TOPTRG events (" << histname << ") from EventT0 DQM not processed!");
222  m_cT0FractionsMuMuTOPTRG->SetFillColor(kGray);
223  }
224 
225  if (m_printCanvas) {
226  m_cT0FractionsHadronECLTRG->Print("EventT0_Algorithm_Efficiency.pdf(");
227  m_cT0FractionsHadronCDCTRG->Print("EventT0_Algorithm_Efficiency.pdf");
228  m_cT0FractionsHadronTOPTRG->Print("EventT0_Algorithm_Efficiency.pdf");
229 
230  m_cT0FractionsBhaBhaECLTRG->Print("EventT0_Algorithm_Efficiency.pdf");
231  m_cT0FractionsBhaBhaCDCTRG->Print("EventT0_Algorithm_Efficiency.pdf");
232  m_cT0FractionsBhaBhaTOPTRG->Print("EventT0_Algorithm_Efficiency.pdf");
233 
234  m_cT0FractionsMuMuECLTRG->Print("EventT0_Algorithm_Efficiency.pdf");
235  m_cT0FractionsMuMuCDCTRG->Print("EventT0_Algorithm_Efficiency.pdf");
236  m_cT0FractionsMuMuTOPTRG->Print("EventT0_Algorithm_Efficiency.pdf)");
237  }
238 
239 }
240 
241 
243 {
253 }
254 
255 
256 bool DQMHistAnalysisEventT0EfficiencyModule::FillEfficiencyHistogram(const std::string& histname, TEfficiency* eff)
257 {
258  B2DEBUG(20, "Begin processing histogram " << histname << " ...");
259  TH1* h = findHist("EventT0/" + histname);
260  if (not h) {
261  return false;
262  }
263 
264  // Admittedly quite a hacky way to obtain the normalisation values: Create a new histogram and fill each of the bins with
265  // 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.
266  const auto totalEntries = h->GetBinContent(-1);
267  const auto nBins = h->GetNbinsX();
268  TH1D* totalHist = new TH1D("total", "total;Algorithm;Fraction #epsilon", nBins, 0, nBins);
269  for (int i = 0; i < nBins; i++) {
270  totalHist->SetBinContent(i + 1, totalEntries);
271  }
272  eff->SetPassedHistogram(*h, "f");
273  eff->SetTotalHistogram(*totalHist, "f");
274 
275  eff->Paint("AP");
276 
277  TGraphAsymmErrors* graph = eff->GetPaintedGraph();
278  if (not graph) {
279  return false;
280  }
281 
282  auto ax = graph->GetXaxis();
283  if (not ax) {
284  return false;
285  }
286  // Print x-axis bin labels horizontally
287  ax->SetTitleOffset(1.0);
288  ax->CenterTitle(kTRUE);
289  ax->Set(nBins, 0, nBins);
290  for (int i = 0; i < nBins; i++) {
291  ax->SetBinLabel(i + 1, c_eventT0Algorithms[i]);
292  }
293 
294  auto ay = graph->GetYaxis();
295  if (not ay) {
296  return false;
297  }
298  ay->SetTitleOffset(1.0);
299  ay->SetRangeUser(0, 1.05);
300 
301  graph->Draw("AP");
302 
303  B2DEBUG(20, "Finished processing histogram " << histname << "!");
304 
305  return true;
306 }
307 
TEfficiency * m_eAlgorithmSourceFractionsMuMuL1ECLTRG
Fraction of events with EventT0 from a given algorithm, HLT mumu events, L1 time by ECL trigger.
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.
const char * c_eventT0Algorithms[6]
EventT0 algorithms for which to calculate fractions of abundance.
TCanvas * m_cT0FractionsMuMuECLTRG
Canvas for time fractions for ECLTRG mumu.
bool FillEfficiencyHistogram(const std::string &histname, TEfficiency *eff)
Fill the TEfficiency plots.
std::string m_prefixCanvas
prefix to be added to canvas name when saved as pdf
MonitoringObject * m_monObj
MonitoringObject to be produced by this module.
TCanvas * m_cT0FractionsBhaBhaCDCTRG
Canvas for time fractions for CDCTRG bhabhas.
TCanvas * m_cT0FractionsHadronCDCTRG
Canvas for time fractions for CDCTRG hadrons.
TEfficiency * m_eAlgorithmSourceFractionsBhaBhaL1CDCTRG
Fraction of events with EventT0 from a given algorithm, HLT bhabha events, L1 time by CDC trigger.
bool m_printCanvas
if true print the pdf of the canvases
TCanvas * m_cT0FractionsBhaBhaECLTRG
Canvas for time fractions for ECLTRG bhabhas.
uint m_nEntriesMin
minimum number of entries to process the histogram
TCanvas * m_cT0FractionsHadronECLTRG
Canvas for time fractions for ECLTRG hadrons.
TEfficiency * m_eAlgorithmSourceFractionsHadronL1ECLTRG
Fraction of events with EventT0 from a given algorithm, HLT hadronic events, L1 time by ECL trigger.
TCanvas * m_cT0FractionsMuMuTOPTRG
Canvas for time fractions for TOPTRG mumu.
TCanvas * m_cT0FractionsMuMuCDCTRG
Canvas for time fractions for CDCTRG mumu.
TEfficiency * m_eAlgorithmSourceFractionsHadronL1CDCTRG
Fraction of events with EventT0 from a given algorithm, HLT hadronic events, L1 time by CDC trigger.
TEfficiency * m_eAlgorithmSourceFractionsMuMuL1TOPTRG
Fraction of events with EventT0 from a given algorithm, HLT mumu events, L1 time by TOP trigger.
TCanvas * m_cT0FractionsHadronTOPTRG
Canvas for time fractions for TOPTRG 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.
TEfficiency * m_eAlgorithmSourceFractionsHadronL1TOPTRG
Fraction of events with EventT0 from a given algorithm, HLT hadronic events, L1 time by TOP trigger.
TCanvas * m_cT0FractionsBhaBhaTOPTRG
Canvas for time fractions for TOPTRG 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
REG_MODULE(arichBtest)
Register the Module.
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
Abstract base class for different kinds of events.