Belle II Software  release-06-00-14
DQMHistAnalysisCDCDedx.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 #include <dqm/analysis/modules/DQMHistAnalysisCDCDedx.h>
10 
11 using namespace std;
12 using namespace Belle2;
13 using boost::format;
14 
15 REG_MODULE(DQMHistAnalysisCDCDedx)
16 
17 //-----------------------------------------------------------------
20 {
21  //Parameter definition here
22  B2DEBUG(20, "DQMHistAnalysisCDCDedx: Constructor done.");
23  addParam("mmode", mmode, "default monitoring mode is basic", std::string("basic"));
24 }
25 
26 //--------------------------------------------------------------
27 DQMHistAnalysisCDCDedxModule::~DQMHistAnalysisCDCDedxModule()
28 {
29 
30 }
31 
32 //---------------------------------------------
33 void DQMHistAnalysisCDCDedxModule::initialize()
34 {
35 
36  gStyle->SetOptStat(0);
37  gROOT->cd();
38 
39  B2DEBUG(20, "DQMHistAnalysisCDCDedx: initialized.");
40 
41  //per run canvas
42  c_pr_dedx = new TCanvas("CDCDedx/c_CDCdedxMean", "", 700, 600);
43 
44  //2D plots
45  c_pr_bands = new TCanvas("CDCDedx/c_pr_bands", "", 700, 600);
46  c_pr_dedxcos = new TCanvas("CDCDedx/c_pr_dedxcos", "", 700, 500);
47  c_pr_dedxphi = new TCanvas("CDCDedx/c_pr_dedxphi", "", 700, 500);
48 
49  //intra-run cavnas 1D
50  c_ir_mean = new TCanvas("CDCDedx/c_ir_mean", "", 700, 500);
51  c_ir_reso = new TCanvas("CDCDedx/c_ir_reso", "", 700, 500);
52 
53  if (mmode != "basic") {
54  c_pr_wires = new TCanvas("CDCDedx/c_pr_wires", "", 700, 700);
55  c_ir_dedx = new TCanvas("CDCDedx/c_ir_dedx", "", 700, 500);
56  }
57 
58  f_gaus = new TF1("f_gaus", "gaus", 0.0, 2.5);
59  f_gaus->SetLineColor(kRed);
60 
61  l_line = new TLine();
62  l_line->SetLineColor(kRed);
63  l_line->SetLineStyle(1);
64  l_line->SetLineWidth(2);
65 
66  //Monitoring object for run dependency at mirabelle
67  m_monObj = getMonitoringObject("cdcdedx");
68  m_monObj->addCanvas(c_pr_dedx);
69  m_monObj->addCanvas(c_pr_bands);
70  m_monObj->addCanvas(c_pr_dedxcos);
71  m_monObj->addCanvas(c_ir_mean);
72  m_monObj->addCanvas(c_ir_reso);
73  m_monObj->addCanvas(c_pr_dedxphi);
74  if (mmode != "basic") {
75  m_monObj->addCanvas(c_pr_wires);
76  m_monObj->addCanvas(c_ir_dedx);
77  }
78 
79 }
80 
81 
82 //-------------------------------------------
83 void DQMHistAnalysisCDCDedxModule::beginRun()
84 {
85  B2DEBUG(20, "DQMHistAnalysisCDCDedx: beginRun called.");
86 }
87 
88 
89 //----------------------------------------
90 void DQMHistAnalysisCDCDedxModule::event()
91 {
92 
93  B2DEBUG(20, "DQMHistAnalysisCDCDedx: event called.");
94 
95  //Plot 0 getmeta those are useful for cosmetics
96  getMetadata();
97 
98  //Plot 1 dE/dx per run gain/reso
99  drawDedx();
100 
101  //Plot 2 dE/dx intra-run gain/reso
102  drawDedxIR();
103 
104  //Plot 3 wire status
105  if (mmode != "basic")drawWireStatus();
106 
107  //Plot 4 dE/dx bands vs p
108  drawBandPlot();
109 
110  //Plot 5 dE/dx vs phi and costh
111  drawDedxCosPhi();
112 
113 }
114 
115 //-----------------------------------------
116 void DQMHistAnalysisCDCDedxModule::endRun()
117 {
118  B2DEBUG(20, "DQMHistAnalysisCDCDedx : endRun called");
119 }
120 
121 //--------------------------------------------
122 void DQMHistAnalysisCDCDedxModule::terminate()
123 {
124 
125  B2DEBUG(20, "DQMHistAnalysisCDCDedx : terminate called");
126 
127  delete c_pr_dedx;
128  delete c_pr_dedxcos;
129  delete c_pr_bands;
130  delete c_ir_mean;
131  delete c_ir_reso;
132  delete c_pr_dedxphi;
133  if (mmode != "basic") {
134  delete c_pr_wires;
135  delete c_ir_dedx;
136  }
137 
138 }
139 
140 //-------------------------------------------
141 void DQMHistAnalysisCDCDedxModule::getMetadata()
142 {
143 
144  TH1D* h_Meta = (TH1D*)findHist("CDCDedx/hMeta");
145 
146  if (h_Meta != nullptr) {
147  m_nallevt = int(h_Meta->GetBinContent(1));
148  m_nbhabhaevt = int(h_Meta->GetBinContent(2));
149  m_nhadevt = int(h_Meta->GetBinContent(3));
150 
151  std::string m_title = h_Meta->GetTitle();
152 
153  first = m_title.find("Exp:");
154  last = m_title.find(", Run:");
155  std::string expN = m_title.substr(first + 4, last - first - 4);
156 
157  first = m_title.find(", Run:");
158  last = m_title.find(", RG:");
159  std::string runN = m_title.substr(first + 6, last - first - 6);
160 
161  first = m_title.find(", RG:");
162  last = m_title.find(")");
163  std::string runGain = m_title.substr(first + 5, last - first - 5);
164 
165  m_exp = std::stoi(expN.c_str());
166  m_run = std::stoi(runN.c_str());
167  m_dbrg = std::stof(runGain.c_str());
168  }
169 
170 }
171 
172 //-------------------------------------------
173 void DQMHistAnalysisCDCDedxModule::drawDedx()
174 {
175 
176  m_mean = 0.0;
177  m_sigma = 0.0;
178 
179  TH1* h_dEdx = findHist("CDCDedx/hdEdx");
180  if (h_dEdx != nullptr) {
181 
182  double m_meanerr = 0.0;
183  double m_sigmaerr = 0.0;
184 
185  c_pr_dedx->Clear();
186  l_line->Clear();
187 
188  m_status = "LowStats";
189 
190  TH1D* h_dEdxClone = (TH1D*)h_dEdx->Clone("hdEdx_clone");
191  if (h_dEdxClone != nullptr && h_dEdxClone->Integral() > 250) {
192  fitHistogram(h_dEdxClone, m_status);
193  if (m_status == "OK") {
194  m_mean = h_dEdxClone->GetFunction("f_gaus")->GetParameter(1);
195  m_sigma = h_dEdxClone->GetFunction("f_gaus")->GetParameter(2);
196  m_meanerr = h_dEdxClone->GetFunction("f_gaus")->GetParError(1);
197  m_sigmaerr = h_dEdxClone->GetFunction("f_gaus")->GetParError(2);
198  }
199  }
200 
201  m_monObj->setVariable("CDCDedxMean", m_mean);
202  m_monObj->setVariable("CDCDedxReso", m_sigma);
203  m_monObj->setVariable("CDCDedxMeanErr", m_meanerr);
204  m_monObj->setVariable("CDCDedxResoErr", m_sigmaerr);
205 
206  c_pr_dedx->cd();
207  set_Pad_Style(0.120, 0.046, 0.071, 0.0);
208  set_Hist_Style(h_dEdxClone);
209  h_dEdxClone->SetTitle("CDC-dEdx");
210  h_dEdxClone->DrawCopy("");
211 
212  //Draw line for dE/dx mean
213  l_line->DrawLine(m_mean, 0, m_mean, h_dEdxClone->GetMaximum());
214 
215  TF1* f_gausC = (TF1*)f_gaus->Clone("f_gausC");
216  f_gausC->SetLineColor(kRed);
217  f_gausC->DrawClone("same");
218 
219 
220  TPaveText* pinfo0 = new TPaveText(0.15, 0.69, 0.37, 0.89, "NBNDC");
221  set_Text_Style(pinfo0);
222  pinfo0->AddText(Form("CDC dE/dx (e^{-}e^{+})"));
223  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
224  pinfo0->AddText(Form("-------------"));
225  pinfo0->AddText(Form("Fit Status: %s", m_status.data()));
226  pinfo0->AddText(Form("Fit #mu^{dE/dx}: %0.3f ", m_mean));
227  pinfo0->AddText(Form("Fit #sigma^{dE/dx}: %0.3f ", m_sigma));
228  pinfo0->SetTextColor(kRed);
229  pinfo0->Draw("same");
230 
231 
232  TPaveText* pinfo1 = new TPaveText(0.68, 0.77, 0.90, 0.89, "NBNDC");
233  set_Text_Style(pinfo1);
234  pinfo1->AddText(Form("-- Expert info"));
235  pinfo1->AddText(Form("Prev Gain: %0.03f", m_dbrg));
236  if (m_nbhabhaevt > 1e5)
237  pinfo1->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
238  if (m_nbhabhaevt > 1e3)
239  pinfo1->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
240  else
241  pinfo1->AddText(Form("Events: %d", m_nbhabhaevt));
242  pinfo1->Draw("same");
243 
244  m_status.clear();
245  c_pr_dedx->Modified();
246  c_pr_dedx->Update();
247  }
248 
249 }
250 
251 //--------------------------------------------
252 void DQMHistAnalysisCDCDedxModule::drawDedxIR()
253 {
254 
255  //1. Draw Scattered plot
256  if (mmode != "basic") {
257 
258  TH2D* hdEdxIRScat = (TH2D*)findHist("CDCDedx/hdEdxvsEvt");
259  if (hdEdxIRScat != nullptr) {
260 
261  c_ir_dedx->Clear();
262  c_ir_dedx->cd();
263  set_Pad_Style(0.143, 0.045, 0.077, 0.0);
264 
265  if (hdEdxIRScat->GetEntries() > 0) {
266  hdEdxIRScat->GetXaxis()->SetRange(hdEdxIRScat->FindFirstBinAbove(0, 1), hdEdxIRScat->FindLastBinAbove(0, 1));
267  }
268 
269  set_Hist_Style(hdEdxIRScat);
270  hdEdxIRScat->Draw("");
271  TPaveText* pinfo = new TPaveText(0.689, 0.690, 0.942, 0.911, "NBNDC");
272  set_Text_Style(pinfo);
273  pinfo->AddText("CDC-dE/dx Intra-run");
274  pinfo->AddText("Electrons (e^{+}e^{-})");
275  pinfo->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
276  if (m_nbhabhaevt > 1e5)
277  pinfo->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
278  if (m_nbhabhaevt > 1e3)
279  pinfo->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
280  else
281  pinfo->AddText(Form("Events: %d", m_nbhabhaevt));
282  pinfo->Draw("same");
283 
284  c_ir_dedx->Modified();
285  c_ir_dedx->Update();
286  }
287  }
288 
289 
290  //Intra rungain/reso variation
291  TH2D* hdEdxIRScatC = (TH2D*)findHist("CDCDedx/hdEdxvsEvt");
292  if (hdEdxIRScatC != nullptr) {
293 
294  c_ir_mean->Clear();
295  c_ir_mean->cd();
296  gPad->SetGridy(1);
297 
298  int fbin = hdEdxIRScatC->FindFirstBinAbove(0, 1);
299  int lbin = hdEdxIRScatC->FindLastBinAbove(0, 1);
300  int nbin = lbin - fbin + 1;
301 
302  TH1D* hdEdxIRInd[nbin];
303  for (int ibin = 0; ibin < nbin; ibin++) {
304  int localbin = ibin + fbin;
305  hdEdxIRInd[ibin] = (TH1D*)hdEdxIRScatC->ProjectionY(Form("htemp_%d", localbin), localbin, localbin);
306  }
307 
308  if (nbin <= 0)nbin = 1;
309  TH1F* hdEdxIRMean = new TH1F("hdEdxIRMean", "", nbin, 0.5, nbin + 0.5);
310  TH1F* hdEdxIRSigma = new TH1F("hdEdxIRSigma", "", nbin, 0.5, nbin + 0.5);
311 
312  for (int ibin = 0; ibin < nbin; ibin++) {
313 
314  double m_imean = 0.0, m_imeanerr = 0.0;
315  double m_isigma = 0.0, m_isigmaerr = 0.0;
316  fitHistogram(hdEdxIRInd[ibin], m_status);
317 
318  if (m_status == "OK") {
319  m_imean = hdEdxIRInd[ibin]->GetFunction("f_gaus")->GetParameter(1);
320  m_imeanerr = hdEdxIRInd[ibin]->GetFunction("f_gaus")->GetParError(1);
321  m_isigma = hdEdxIRInd[ibin]->GetFunction("f_gaus")->GetParameter(2);
322  m_isigmaerr = hdEdxIRInd[ibin]->GetFunction("f_gaus")->GetParError(2);
323  }
324  hdEdxIRMean->SetBinContent(ibin + 1, m_imean);
325  hdEdxIRMean->SetBinError(ibin + 1, m_imeanerr);
326  hdEdxIRSigma->SetBinContent(ibin + 1, m_isigma);
327  hdEdxIRSigma->SetBinError(ibin + 1, m_isigmaerr);
328  }
329 
330  //2 intra-gain trend
331  set_Pad_Style(0.143, 0.045, 0.077, 0.0);
332  set_Hist_Style(hdEdxIRMean);
333  hdEdxIRMean->SetMarkerColor(kRed);
334  hdEdxIRMean->SetMarkerStyle(20);
335  hdEdxIRMean->SetMarkerSize(1.10);
336  hdEdxIRMean->GetYaxis()->SetRangeUser(0.80, 1.20);
337  hdEdxIRMean->GetYaxis()->SetTitle("dE/dx (#mu_{fit})");
338  hdEdxIRMean->GetXaxis()->SetTitle("Events(M)");
339  hdEdxIRMean->SetTitle("CDC-dE/dx gain(#mu): intra-run variation");
340  hdEdxIRMean->Draw("");
341 
342  l_line->DrawLine(0.5, m_mean, hdEdxIRMean->GetXaxis()->GetBinUpEdge(nbin), m_mean);
343 
344  TPaveText* pinfo0 = new TPaveText(0.689, 0.680, 0.942, 0.911, "NBNDC");
345  set_Text_Style(pinfo0);
346  pinfo0->AddText("Intra-run variation");
347  pinfo0->AddText("Electrons (e^{+}e^{-})");
348  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
349  pinfo0->AddText(Form("Avg #mu_{fit}: %0.3f", m_mean));
350  if (m_nbhabhaevt > 1e5)
351  pinfo0->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
352  if (m_nbhabhaevt > 1e3)
353  pinfo0->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
354  else
355  pinfo0->AddText(Form("Events: %d", m_nbhabhaevt));
356  pinfo0->Draw("same");
357  c_ir_mean->Modified();
358  c_ir_mean->Update();
359 
360  //3 intra-resolution trend
361  c_ir_reso->Clear();
362  c_ir_reso->cd();
363  gPad->SetGridy(1);
364  set_Pad_Style(0.143, 0.045, 0.077, 0.0);
365  set_Hist_Style(hdEdxIRSigma);
366  hdEdxIRSigma->SetMarkerColor(kRed);
367  hdEdxIRSigma->SetMarkerStyle(20);
368  hdEdxIRSigma->SetMarkerSize(1.10);
369  hdEdxIRSigma->GetYaxis()->SetRangeUser(0.03, 0.20);
370  hdEdxIRSigma->GetYaxis()->SetTitle("dE/dx (#sigma_{fit})");
371  hdEdxIRSigma->GetXaxis()->SetTitle("Events(M)");
372  hdEdxIRSigma->SetTitle("CDC-dE/dx reso.(#sigma): intra-run variation");
373  hdEdxIRSigma->Draw("");
374 
375  l_line->DrawLine(0.5, m_sigma, hdEdxIRSigma->GetXaxis()->GetBinUpEdge(nbin), m_sigma);
376 
377  TPaveText* pinfo1 = new TPaveText(0.689, 0.680, 0.942, 0.911, "NBNDC");
378  set_Text_Style(pinfo1);
379  pinfo1->AddText("Intra-run variation");
380  pinfo1->AddText("Electrons (e^{+}e^{-})");
381  pinfo1->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
382  pinfo1->AddText(Form("Avg #sigma_{fit}: %0.3f", m_sigma));
383  if (m_nbhabhaevt > 1e5)
384  pinfo1->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
385  if (m_nbhabhaevt > 1e3)
386  pinfo1->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
387  else
388  pinfo1->AddText(Form("Events: %d", m_nbhabhaevt));
389  pinfo1->Draw("same");
390  c_ir_reso->Modified();
391  c_ir_reso->Update();
392 
393  //Let's reset histogram here
394  for (int ibin = 0; ibin < nbin; ibin++) hdEdxIRInd[ibin]->Reset();
395  }
396 
397 }
398 
399 //------------------------------------------------
400 void DQMHistAnalysisCDCDedxModule::drawWireStatus()
401 {
402 
403  //Draw Scattered plot
404  TH2D* hWires = (TH2D*)findHist("CDCDedx/hWires");
405  TH2D* hWireStatus = (TH2D*)findHist("CDCDedx/hWireStatus");
406  if (hWires != nullptr && hWireStatus != nullptr) {
407  c_pr_wires->Clear();
408  c_pr_wires->cd();
409  set_Hist_Style(hWires);
410  hWires->SetMarkerColor(kGray + 1);
411  hWires->Draw("");
412 
413  std::string s_ndead = hWireStatus->GetTitle();
414  int m_ndead = atof(s_ndead.c_str());
415  m_monObj->setVariable("CDCDedxDeadWires", m_ndead);
416 
417  set_Hist_Style(hWireStatus);
418  hWireStatus->SetMarkerColor(kRed);
419  hWireStatus->SetMarkerStyle(7);
420  hWireStatus->Draw("same");
421 
422  TPaveText* pinfo0 = new TPaveText(0.117, 0.832, 0.148, 0.976, "NBNDC");
423  set_Text_Style(pinfo0);
424  pinfo0->AddText(Form("CDC Wire Status"));
425  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
426  pinfo0->AddText(Form("Dead: %d (%0.02f%%)", m_ndead, (100.0 * m_ndead / 14336.0)));
427  if (m_nallevt > 1e5)
428  pinfo0->AddText(Form("Events: %0.02fM", double(m_nallevt / 1e6)));
429  if (m_nallevt > 1e3)
430  pinfo0->AddText(Form("Events: %0.02fK", double(m_nallevt / 1e3)));
431  else
432  pinfo0->AddText(Form("Events: %d", m_nallevt));
433  pinfo0->Draw("same");
434 
435  c_pr_wires->Modified();
436  c_pr_wires->Update();
437  }
438 }
439 
440 //-----------------------------------------------
441 void DQMHistAnalysisCDCDedxModule::drawBandPlot()
442 {
443  //Draw Scattered plot
444  TH2D* hdEdxVsP = (TH2D*)findHist("CDCDedx/hdEdxVsP");
445  if (hdEdxVsP != nullptr) {
446  c_pr_bands->Clear();
447  c_pr_bands->cd();
448  c_pr_bands->SetLogx();
449  c_pr_bands->SetLogy();
450 
451  set_Plot_Style();
452  set_Hist_Style(hdEdxVsP);
453  hdEdxVsP->SetTitle("CDC-dEdx band plots");
454  hdEdxVsP->Draw("col");
455 
456  TPaveText* pinfo0 = new TPaveText(0.59, 0.75, 0.80, 0.88, "NBNDC");
457  set_Text_Style(pinfo0);
458  pinfo0->AddText(Form("IP tracks (hadron)"));
459  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
460  if (m_nhadevt > 1e5)
461  pinfo0->AddText(Form("Events: %0.02fM", double(m_nhadevt / 1e6)));
462  if (m_nhadevt > 1e3)
463  pinfo0->AddText(Form("Events: %0.02fK", double(m_nhadevt / 1e3)));
464  else
465  pinfo0->AddText(Form("Events: %d", m_nhadevt));
466  pinfo0->DrawClone("same");
467 
468  c_pr_bands->Modified();
469  c_pr_bands->Update();
470  }
471 
472 }
473 
474 //---------------------------------------------
475 void DQMHistAnalysisCDCDedxModule::drawDedxCosPhi()
476 {
477 
478  TH2D* hdEdxvsPhi = (TH2D*)findHist("CDCDedx/hdEdxvsPhi");
479  if (hdEdxvsPhi != nullptr) {
480  c_pr_dedxphi->Clear();
481  c_pr_dedxphi->cd();
482  set_Pad_Style(0.143, 0.045, 0.077, 0.0);
483 
484  set_Hist_Style(hdEdxvsPhi);
485  hdEdxvsPhi->SetTitle("CDC-dEdx vs Phi");
486  hdEdxvsPhi->Draw("col");
487 
488  l_line->DrawLine(-3.20, 1.0, 3.20, 1.0);
489 
490  TPaveText* pinfo0 = new TPaveText(0.689, 0.751, 0.942, 0.89, "NBNDC");
491  set_Text_Style(pinfo0);
492  pinfo0->AddText(Form("Electrons (e^{+}e^{-})"));
493  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
494  if (m_nbhabhaevt > 1e5)
495  pinfo0->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
496  if (m_nbhabhaevt > 1e3)
497  pinfo0->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
498  else
499  pinfo0->AddText(Form("Events: %d", m_nbhabhaevt));
500  pinfo0->DrawClone("same");
501 
502  c_pr_dedxphi->Modified();
503  c_pr_dedxphi->Update();
504  }
505 
506  //plot # 2
507  TH2D* hdEdxvsCosth = (TH2D*)findHist("CDCDedx/hdEdxvsCosth");
508  if (hdEdxvsCosth != nullptr) {
509  c_pr_dedxcos->Clear();
510  c_pr_dedxcos->cd();
511  set_Pad_Style(0.143, 0.045, 0.077, 0.0);
512 
513  set_Hist_Style(hdEdxvsCosth);
514  hdEdxvsCosth->SetTitle("CDC-dEdx vs Costh");
515  hdEdxvsCosth->Draw("col");
516 
517  l_line->DrawLine(-1.0, 1.0, 1.0, 1.0);
518 
519  TPaveText* pinfo1 = new TPaveText(0.689, 0.751, 0.942, 0.89, "NBNDC");
520  set_Text_Style(pinfo1);
521  pinfo1->AddText(Form("Electrons (e^{+}e^{-})"));
522  pinfo1->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
523  if (m_nbhabhaevt > 1e5)
524  pinfo1->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
525  if (m_nbhabhaevt > 1e3)
526  pinfo1->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
527  else
528  pinfo1->AddText(Form("Events: %d", m_nbhabhaevt));
529  pinfo1->DrawClone("same");
530 
531  c_pr_dedxcos->Modified();
532  c_pr_dedxcos->Update();
533  }
534 
535 }
536 
537 //----------------------------------------------------------------------------------------
538 void DQMHistAnalysisCDCDedxModule::fitHistogram(TH1D*& temphist, std::string& status)
539 {
540 
541  if (temphist != nullptr) {
542  temphist->GetXaxis()->SetRange(temphist->FindFirstBinAbove(0, 1), temphist->FindLastBinAbove(0, 1));
543  int fs = temphist->Fit(f_gaus, "QR");
544  if (fs != 0) status = "Failed";
545  else {
546  double mean = temphist->GetFunction("f_gaus")->GetParameter(1);
547  double width = temphist->GetFunction("f_gaus")->GetParameter(2);
548  temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
549  fs = temphist->Fit(f_gaus, "QR", "", mean - 2.5 * width, mean + 2.5 * width);
550  if (fs != 0)status = "Failed";
551  else {
552  temphist->GetXaxis()->SetRangeUser(mean - 6.0 * width, mean + 6.0 * width);
553  status = "OK";
554  }
555  }
556  }
557 
558 }
559 
560 //------------------------------------------------
561 void DQMHistAnalysisCDCDedxModule::set_Plot_Style()
562 {
563 
564  const Int_t NRGBs = 6;
565  const Int_t NCont = 255;
566  Double_t stops[NRGBs] = { 0.00, 0.00, 0.34, 0.61, 0.84, 1.00 };
567  Double_t red[NRGBs] = { 0.00, 0.00, 0.00, 0.87, 1.00, 0.51 };
568  Double_t green[NRGBs] = { 0.00, 0.00, 0.81, 1.00, 0.20, 0.00 };
569  Double_t blue[NRGBs] = { 0.00, 0.51, 1.00, 0.12, 0.00, 0.00 };
570  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
571  gStyle->SetNumberContours(NCont);
572 
573 }
574 
575 
576 //------------------------------------------------------------------
577 void DQMHistAnalysisCDCDedxModule::set_Text_Style(TPaveText*& obj)
578 {
579 
580  obj->SetFillColor(0);
581  obj->SetFillStyle(0);
582 
583  obj->SetLineColor(TColor::GetColor("#000000"));
584  obj->SetLineWidth(0);
585  obj->SetTextAlign(12);
586 
587  obj->SetTextColor(kGray + 3);
588  obj->SetTextFont(82);
589  obj->SetTextSize(0.03157895);
590  obj->SetTextAlign(12);
591 
592 }
593 
594 
595 //------------------------------------------------------------------
596 void DQMHistAnalysisCDCDedxModule::set_Hist_Style(TH1* obj)
597 {
598 
599  obj->SetStats(0);
600  obj->SetTitle("");
601  obj->SetFillColor(kYellow);
602 
603  obj->SetTitleOffset(1.15, "x");
604  obj->SetTitleSize(.040, "x");
605  obj->SetTitleOffset(1.15, "y");
606  obj->SetTitleSize(.040, "y");
607 
608  obj->SetLabelOffset(0.015, "x");
609  obj->SetLabelSize(.040, "x");
610  obj->SetLabelOffset(0.015, "y");
611  obj->SetLabelSize(.040, "y");
612 
613  obj->SetTickLength(0.03, "x");
614  obj->SetTickLength(0.02, "y");
615 
616 }
617 
618 
619 //-------------------------------------------------------------------------------------
620 void DQMHistAnalysisCDCDedxModule::set_Pad_Style(double l, double r, double t, double b)
621 {
622 
623  if (l != 0)gPad->SetLeftMargin(l);
624  if (r != 0)gPad->SetRightMargin(r);
625  if (t != 0)gPad->SetTopMargin(t);
626  if (b != 0)gPad->SetBottomMargin(b);
627  gPad->SetTickx(1);
628  gPad->SetTicky(1);
629 
630 }
DQM analysis module grab canvases from DQM module and perform higher level operation like histogram f...
The base class for the histogram analysis module.
#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.