Belle II Software  release-05-01-25
DQMHistAnalysisEventT0.cc
1 //+
2 // File : DQMHistAnalysisEventT0.cc
3 // Description : module for DQM histogram analysis of SVD sensors occupancies
4 //
5 // Author : Giulia Casarosa (PI)
6 // Date : 20210210
7 //-
8 
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisEventT0.h>
11 
12 #include <TROOT.h>
13 #include <TStyle.h>
14 #include <TString.h>
15 #include <TFitResult.h>
16 
17 #include <TMath.h>
18 #include <iostream>
19 
20 using namespace std;
21 using namespace Belle2;
22 
23 //-----------------------------------------------------------------
24 // Register the Module
25 //-----------------------------------------------------------------
26 REG_MODULE(DQMHistAnalysisEventT0)
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
34 {
35  //Parameter definition
36  addParam("min_nEntries", m_nEntriesMin, "minimum numeber of entries to process the histogram", m_nEntriesMin);
37  addParam("prefixCanvas", m_prefixCanvas, "prefix to be added to canvas filename when saved as pdf", std::string("c"));
38  addParam("printCanvas", m_printCanvas, "if True prints pdf of the analysis canvas", bool(false));
39 }
40 
41 
42 DQMHistAnalysisEventT0Module::~DQMHistAnalysisEventT0Module() { }
43 
44 void DQMHistAnalysisEventT0Module::initialize()
45 {
46  gROOT->cd();
47 
48  //ECLTRG canvas
49  m_cECLTRG = new TCanvas("EventT0/c_ECLTRGjitter", "ECLTRG jitter", 1200, 400);
50  m_pad1ECLTRG = new TPad("pad1ECLTRG", "pad1 ECLTRG", 0.03, 0.02, 0.33, 0.98);
51  m_pad2ECLTRG = new TPad("pad2ECLTRG", "pad2 ECLTRG", 0.35, 0.02, 0.65, 0.98);
52  m_pad3ECLTRG = new TPad("pad3ECLTRG", "pad3 ECLTRG", 0.67, 0.02, 0.97, 0.98);
53 
54  //CDCTRG canvases
55  m_cCDCTRG = new TCanvas("EventT0/c_CDCTRGjitter", "CDCTRG jitter", 1200, 400);
56  m_pad1CDCTRG = new TPad("pad1CDCTRG", "pad1 CDCTRG", 0.03, 0.02, 0.33, 0.98);
57  m_pad2CDCTRG = new TPad("pad2CDCTRG", "pad2 CDCTRG", 0.35, 0.02, 0.65, 0.98);
58  m_pad3CDCTRG = new TPad("pad3CDCTRG", "pad3 CDCTRG", 0.67, 0.02, 0.97, 0.98);
59 
60  m_monObj = getMonitoringObject("eventT0");
61  m_monObj->addCanvas(m_cECLTRG);
62  m_monObj->addCanvas(m_cCDCTRG);
63 }
64 
65 
66 void DQMHistAnalysisEventT0Module::beginRun()
67 {
68  m_cECLTRG->Clear();
69  m_cCDCTRG->Clear();
70 }
71 
72 void DQMHistAnalysisEventT0Module::endRun()
73 {
74 
75 
76  // --- ECLTRG ---
77 
78  //find TOP EventT0 Hadrons ECLTRG histogram and process it
79  TH1* h = findHist("EventT0DQMdir/m_histEventT0_TOP_hadron_L1_ECLTRG");
80  TString tag = "hadronECLTRG";
81  m_pad1ECLTRG->cd();
82  if (processHistogram(h, tag)) {
83  m_pad1ECLTRG->SetFillColor(0);
84  m_pad1ECLTRG->Modified();
85  m_pad1ECLTRG->Update();
86  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
87  h->Draw();
88  m_pad1ECLTRG->SetFillColor(kGray);
89  }
90 
91  //find TOP EventT0 Bhabhas ECLTRG histogram and process it
92  h = findHist("EventT0DQMdir/m_histEventT0_TOP_bhabha_L1_ECLTRG");
93  tag = "bhabhaECLTRG";
94  m_pad2ECLTRG->cd();
95  if (processHistogram(h, tag)) {
96  m_pad2ECLTRG->SetFillColor(0);
97  m_pad2ECLTRG->Modified();
98  m_pad2ECLTRG->Update();
99  } else {
100  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
101  h->Draw();
102  m_pad2ECLTRG->SetFillColor(kGray);
103  }
104 
105 
106  //find TOP EventT0 Mumus ECLTRG histogram and process it
107  h = findHist("EventT0DQMdir/m_histEventT0_TOP_mumu_L1_ECLTRG");
108  tag = "mumuECLTRG";
109  m_pad3ECLTRG->cd();
110  if (processHistogram(h, tag)) {
111  m_pad3ECLTRG->SetFillColor(0);
112  m_pad3ECLTRG->Modified();
113  m_pad3ECLTRG->Update();
114  } else {
115  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
116  h->Draw();
117  m_pad3ECLTRG->SetFillColor(kGray);
118  }
119 
120  m_cECLTRG->cd();
121  m_pad1ECLTRG->Draw();
122  m_pad2ECLTRG->Draw();
123  m_pad3ECLTRG->Draw();
124 
125  if (m_printCanvas)
126  m_cECLTRG->Print(Form("%s_ECLTRG.pdf", m_prefixCanvas.c_str()));
127 
128 
129  // --- CDCTRG ---
130 
131  //find TOP EventT0 Hadrons CDCTRG histogram and process it
132  h = findHist("EventT0DQMdir/m_histEventT0_TOP_hadron_L1_CDCTRG");
133  tag = "hadronCDCTRG";
134  m_pad1CDCTRG->cd();
135  if (processHistogram(h, tag)) {
136  m_pad1CDCTRG->SetFillColor(0);
137  m_pad1CDCTRG->Modified();
138  m_pad1CDCTRG->Update();
139  m_cCDCTRG->cd();
140  m_pad1CDCTRG->Draw();
141  } else {
142  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
143  h->Draw();
144  m_pad1CDCTRG->SetFillColor(kGray);
145  m_cCDCTRG->cd();
146  m_pad1CDCTRG->Draw();
147  }
148 
149  //find TOP EventT0 Bhabhas CDCTRG histogram and process it
150  h = findHist("EventT0DQMdir/m_histEventT0_TOP_bhabha_L1_CDCTRG");
151  tag = "bhabhaCDCTRG";
152  m_pad2CDCTRG->cd();
153  if (processHistogram(h, tag)) {
154  m_pad2CDCTRG->SetFillColor(0);
155  m_pad2CDCTRG->Modified();
156  m_pad2CDCTRG->Update();
157  m_cCDCTRG->cd();
158  m_pad2CDCTRG->Draw();
159  } else {
160  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
161  h->Draw();
162  m_pad2CDCTRG->SetFillColor(kGray);
163  m_cCDCTRG->cd();
164  m_pad2CDCTRG->Draw();
165  }
166 
167 
168  //find TOP EventT0 Mumus CDCTRG histogram and process it
169  h = findHist("EventT0DQMdir/m_histEventT0_TOP_mumu_L1_CDCTRG");
170  tag = "mumuCDCTRG";
171  m_pad3CDCTRG->cd();
172  if (processHistogram(h, tag)) {
173  m_pad3CDCTRG->SetFillColor(0);
174  m_pad3CDCTRG->Modified();
175  m_pad3CDCTRG->Update();
176  } else {
177  B2WARNING(Form("Histogram TOP EventT0 for %s from EventT0 DQM not processed!", tag.Data()));
178  h->Draw();
179  m_pad3CDCTRG->SetFillColor(kGray);
180  }
181 
182  m_cCDCTRG->cd();
183  m_pad1CDCTRG->Draw();
184  m_pad2CDCTRG->Draw();
185  m_pad3CDCTRG->Draw();
186 
187  if (m_printCanvas)
188  m_cCDCTRG->Print(Form("%s_CDCTRG.pdf", m_prefixCanvas.c_str()));
189 
190 
191 }
192 
193 void DQMHistAnalysisEventT0Module::terminate()
194 {
195 
196  delete m_cECLTRG;
197  delete m_cCDCTRG;
198 
199 }
200 
201 
202 double DQMHistAnalysisEventT0Module::fDoubleGaus(double* x, double* par)
203 {
204  double N = par[0];
205  double frac = par[1];
206  double mean = par[2];
207  double sigma = par[3];
208  double mean2 = par[4];
209  double sigma2 = par[5];
210 
211  return N * frac * TMath::Gaus(x[0], mean, sigma) + N * (1 - frac) * TMath::Gaus(x[0], mean2, sigma2);
212 }
213 
214 bool DQMHistAnalysisEventT0Module::processHistogram(TH1* h, TString tag)
215 {
216 
217  if (h == nullptr) {
218  B2DEBUG(20, "h == nullptr");
219  m_monObj->setVariable(Form("fit_%s", tag.Data()), 0);
220  return false;
221  }
222 
223  int nToFit = h->GetEntries();// - h->GetBinContent(0) - h->GetBinContent(h->GetNbinsX()+1);
224  if (nToFit < m_nEntriesMin) {
225  B2DEBUG(20, "not enough entries");
226  m_monObj->setVariable(Form("fit_%s", tag.Data()), 0);
227  return false;
228  }
229 
230 
231  //scale the histogram
232  h->Scale(1. / h->GetEntries());
233  h->GetXaxis()->SetRangeUser(-50, 50);
234 
235  //define the fitting function
236  TF1 fitf("fit", DQMHistAnalysisEventT0Module::fDoubleGaus, -50, 50, 6);
237  fitf.SetParNames("N", "f_{1}", "#mu_{1}", "#sigma_{1}", "#mu_{2}", "#sigma_{2}");
238  fitf.SetParameters(0.1, 0.8, 0, 5, 0, 15);
239  fitf.SetParLimits(1, 0, 1); //fraction
240  fitf.SetParLimits(3, 0, 100); //sigma1
241  fitf.SetParLimits(5, 0, 100); //sigma2
242 
243  if (h->Fit(&fitf, "SR+") != 0) {
244  B2DEBUG(20, "failed fit");
245  m_monObj->setVariable(Form("fit_%s", tag.Data()), 0);
246  return false;
247  }
248 
249  Double_t par[6];
250  fitf.GetParameters(&par[0]);
251  Double_t parErr[6];
252  for (int i = 0; i < 6; i++)
253  parErr[i] = fitf.GetParError(i) ;
254 
255 
256  //define gaussian components
257  TF1 gauss1("gauss1", "gaus", -100, 100);
258  TF1 gauss2("gauss2", "gaus", -100, 100);
259 
260  gauss1.SetLineColor(kBlue);
261  gauss1.SetLineStyle(kDashed);
262  gauss1.SetParameters(par[0]*par[1], par[2], par[3]);
263 
264  gauss2.SetLineColor(kBlue);
265  gauss2.SetLineStyle(kDashed);
266  gauss2.SetParameters(par[0] * (1 - par[1]), par[4], par[5]);
267 
268  m_monObj->setVariable(Form("fit_%s", tag.Data()), 1);
269  m_monObj->setVariable(Form("N_%s", tag.Data()), h->GetEntries(), TMath::Sqrt(h->GetEntries()));
270  m_monObj->setVariable(Form("f_%s", tag.Data()), par[1], parErr[1]);
271  m_monObj->setVariable(Form("mean1_%s", tag.Data()), par[2], parErr[2]);
272  m_monObj->setVariable(Form("sigma1_%s", tag.Data()), par[3], parErr[3]);
273  m_monObj->setVariable(Form("mean2_%s", tag.Data()), par[4], parErr[4]);
274  m_monObj->setVariable(Form("sigma2_%s", tag.Data()), par[5], parErr[5]);
275 
276  //SETUP gSTYLE - all plots
277  gStyle->SetOptFit(1111);
278 
279  h->DrawClone();
280  fitf.DrawClone("same");
281  gauss1.DrawClone("same");
282  gauss2.DrawClone("same");
283 
284  return true;
285 
286 }
287 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27
Belle2::DQMHistAnalysisEventT0Module
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisEventT0.h:26