Belle II Software  release-08-01-10
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 #include <vector>
11 
12 using namespace std;
13 using namespace Belle2;
14 using boost::format;
15 
16 REG_MODULE(DQMHistAnalysisCDCDedx);
17 
18 //-----------------------------------------------------------------
19 DQMHistAnalysisCDCDedxModule::DQMHistAnalysisCDCDedxModule()
21 {
22  B2DEBUG(20, "DQMHistAnalysisCDCDedx: Constructor done.");
23  setDescription("Module to draw and compute values related to dEdx for CDC. ");
24 
25  //Parameter definition here
26  addParam("mmode", mmode, "default monitoring mode is basic", std::string("basic"));
27 }
28 
29 
30 //---------------------------------------------
32 {
33 
34  gStyle->SetOptStat(0);
35  gROOT->cd();
36 
37  B2DEBUG(20, "DQMHistAnalysisCDCDedx: initialized.");
38 
39  //intra-run cavnas 1D
40  c_ir_dedx = new TCanvas("CDCDedx/c_ir_dedx", "", 1400, 500);
41 
42  //per run canvas
43  c_pr_dedx = new TCanvas("CDCDedx/c_CDCdedxMean", "", 800, 800);
44 
45  f_gaus = new TF1("f_gaus", "gaus", 0.0, 2.5);
46  f_gaus->SetLineColor(kRed);
47 
48  l_line = new TLine();
49  l_line->SetLineColor(kRed);
50  l_line->SetLineStyle(1);
51  l_line->SetLineWidth(2);
52 
53  //Monitoring object for run dependency at mirabelle
54  m_monObj = getMonitoringObject("cdcdedx");
57 
58 }
59 
60 
61 //-------------------------------------------
63 {
64  B2DEBUG(20, "DQMHistAnalysisCDCDedx: beginRun called.");
65 }
66 
67 
68 //----------------------------------------
70 {
71 
72  B2DEBUG(20, "DQMHistAnalysisCDCDedx: event called.");
73 
74  //Plot 0 getmeta those are useful for cosmetics
75  getMetadata();
76 
77  c_pr_dedx->Clear();
78  if (mmode != "basic") {
79  c_pr_dedx->SetWindowSize(1400, 800);
80  c_pr_dedx->Divide(3, 3);
81  } else {
82  c_pr_dedx->SetWindowSize(1000, 800);
83  c_pr_dedx->Divide(3, 2);
84  }
85 
86  //Plot 1 dE/dx per run gain/reso
87  drawDedxPR();
88 
89  // Plot 2 dE/dx bands vs p
90  drawBandPlot();
91 
92  // Plot 3/4 dE/dx vs phi and costh
94 
95  // Plot 5/6 dE/dx mean/reso vs inject time
97 
98  c_pr_dedx->Modified();
99  c_pr_dedx->Update();
100 
101  c_ir_dedx->Clear();
102  if (mmode != "basic") {
103  c_ir_dedx->SetWindowSize(1400, 400);
104  c_ir_dedx->Divide(3, 1);
105  } else {
106  c_ir_dedx->SetWindowSize(900, 400);
107  c_ir_dedx->Divide(2, 1);
108  }
109 
110  //Plot 1 dE/dx intra-run gain/reso
111  drawDedxIR();
112 
113  c_ir_dedx->Modified();
114  c_ir_dedx->Update();
115 
116  //Plot 7/8 wire status/ injection time
117  if (mmode != "basic") {
118  drawDedxInjTime();
119  drawWireStatus();
120  }
121 }
122 
123 //-----------------------------------------
125 {
126  B2DEBUG(20, "DQMHistAnalysisCDCDedx : endRun called");
127 }
128 
129 //--------------------------------------------
131 {
132 
133  B2DEBUG(20, "DQMHistAnalysisCDCDedx : terminate called");
134 
135  delete c_pr_dedx;
136  delete c_ir_dedx;
137 
138 }
139 
140 //-------------------------------------------
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 //-------------------------------------------
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->cd(1);
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) {
192  if (h_dEdxClone->Integral() > 250) {
193  fitHistogram(h_dEdxClone, m_status);
194  if (m_status == "OK") {
195  m_mean = h_dEdxClone->GetFunction("f_gaus")->GetParameter(1);
196  m_sigma = h_dEdxClone->GetFunction("f_gaus")->GetParameter(2);
197  m_meanerr = h_dEdxClone->GetFunction("f_gaus")->GetParError(1);
198  m_sigmaerr = h_dEdxClone->GetFunction("f_gaus")->GetParError(2);
199  }
200  }
201  setHistStyle(h_dEdxClone);
202  h_dEdxClone->SetTitle("CDC-dEdx");
203  h_dEdxClone->DrawCopy("");
204  //Draw line for dE/dx mean
205  l_line->DrawLine(m_mean, 0, m_mean, h_dEdxClone->GetMaximum());
206  }
207 
208  m_monObj->setVariable("CDCDedxMean", m_mean);
209  m_monObj->setVariable("CDCDedxReso", m_sigma);
210  m_monObj->setVariable("CDCDedxMeanErr", m_meanerr);
211  m_monObj->setVariable("CDCDedxResoErr", m_sigmaerr);
212 
213  TF1* f_gausC = (TF1*)f_gaus->Clone("f_gausC");
214  f_gausC->SetLineColor(kRed);
215  f_gausC->DrawClone("same");
216 
217  TPaveText* pinfo0 = new TPaveText(0.12, 0.72, 0.37, 0.89, "NBNDC");
218  setTextStyle(pinfo0);
219  pinfo0->AddText(Form("Fit Status: %s", m_status.data()));
220  pinfo0->AddText(Form("Fit #mu^{dE/dx}: %0.3f ", m_mean));
221  pinfo0->AddText(Form("Fit #sigma^{dE/dx}: %0.3f ", m_sigma));
222  pinfo0->SetTextColor(kRed);
223  pinfo0->Draw("same");
224 
225  TPaveText* pinfo1 = new TPaveText(0.60, 0.69, 0.85, 0.89, "NBNDC");
226  setTextStyle(pinfo1);
227  pinfo1->AddText(Form("-- Expert info"));
228  pinfo1->AddText(Form("-------------"));
229  setBEvtInfo(pinfo1);
230  pinfo1->AddText(Form("Prev Gain: %0.03f", m_dbrg));
231  pinfo1->Draw("same");
232 
233  m_status.clear();
234  }
235 
236 }
237 
238 //--------------------------------------------
240 {
241 
242  //1. Draw Scattered plot
243  if (mmode != "basic") {
244 
245  c_ir_dedx->cd(3);
246 
247  TH2D* hdEdxIRScat = (TH2D*)findHist("CDCDedx/hdEdxvsEvt");
248  if (hdEdxIRScat != nullptr) {
249 
250  setPadStyle(0.143, 0.045, 0.077, 0.0);
251 
252  if (hdEdxIRScat->GetEntries() > 0) {
253  hdEdxIRScat->GetXaxis()->SetRange(hdEdxIRScat->FindFirstBinAbove(0, 1), hdEdxIRScat->FindLastBinAbove(0, 1));
254  }
255 
256  setHistStyle(hdEdxIRScat);
257  hdEdxIRScat->Draw("");
258  TPaveText* pinfo = new TPaveText(0.609, 0.710, 0.942, 0.911, "NBNDC");
259  setTextStyle(pinfo);
260  pinfo->AddText("CDC-dE/dx Intra-run");
261  setBEvtInfo(pinfo);
262  pinfo->Draw("same");
263  }
264  }
265 
266  //Intra rungain/reso variation
267  TH2D* hdEdxIRScatC = (TH2D*)findHist("CDCDedx/hdEdxvsEvt");
268 
269  if (hdEdxIRScatC != nullptr) {
270 
271  c_ir_dedx->cd(1);
272  gPad->SetGridy(1);
273 
274  int fbin = hdEdxIRScatC->FindFirstBinAbove(0, 1);
275  int lbin = hdEdxIRScatC->FindLastBinAbove(0, 1);
276  int nbin = lbin - fbin + 1;
277 
278  if (nbin <= 0) nbin = 1;
279 
280  TH1F* hdEdxIRMean = new TH1F("hdEdxIRMean", "", nbin, 0.5, nbin + 0.5);
281  hdEdxIRMean->SetTitle("CDC-dE/dx gain(#mu): intra-run variation;Events(M);dE/dx (#mu_{fit})");
282 
283  TH1F* hdEdxIRSigma = new TH1F("hdEdxIRSigma", "", nbin, 0.5, nbin + 0.5);
284  hdEdxIRSigma->SetTitle("CDC-dE/dx reso.(#sigma): intra-run variation;Events(M);dE/dx (#sigma_{fit})");
285 
286  setHistPars(hdEdxIRScatC, hdEdxIRMean, hdEdxIRSigma, nbin);
287 
288  //2 intra-gain trend
289 
290  drawHistPars(hdEdxIRMean, nbin, m_mean, 0.10, "#mu_{fit}");
291 
292  //3 intra-resolution trend
293  c_ir_dedx->cd(2);
294  gPad->SetGridy(1);
295 
296  drawHistPars(hdEdxIRSigma, nbin, m_sigma, 0.04, "#sigma_{fit}");
297  }
298 }
299 
300 //-----------------------------------------------
302 {
303  //Draw Scattered plot
304  TH2D* hdEdxVsP = (TH2D*)findHist("CDCDedx/hdEdxVsP");
305  if (hdEdxVsP != nullptr) {
306 
307  c_pr_dedx->cd(2);
308  gPad->SetLogx();
309  gPad->SetLogy();
310 
311  setPlotStyle();
312  setHistStyle(hdEdxVsP);
313  hdEdxVsP->SetTitle("CDC-dEdx band plot");
314  hdEdxVsP->SetMinimum(0.10);
315  hdEdxVsP->Draw("col");
316 
317  TPaveText* pinfo0 = new TPaveText(0.60, 0.77, 0.85, 0.89, "NBNDC");
318  setTextStyle(pinfo0);
319  pinfo0->AddText(Form("IP tracks (hadron)"));
320  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
321  if (m_nhadevt > 1e5)
322  pinfo0->AddText(Form("Events: %0.02fM", double(m_nhadevt / 1e6)));
323  if (m_nhadevt > 1e3)
324  pinfo0->AddText(Form("Events: %0.02fK", double(m_nhadevt / 1e3)));
325  else
326  pinfo0->AddText(Form("Events: %d", m_nhadevt));
327  pinfo0->DrawClone("same");
328  }
329 
330 }
331 
332 
333 //---------------------------------------------
335 {
336 
337  TH2D* hdEdxvsPhi = (TH2D*)findHist("CDCDedx/hdEdxvsPhi");
338  if (hdEdxvsPhi != nullptr) {
339 
340  c_pr_dedx->cd(3);
341 
342  setHistStyle(hdEdxvsPhi);
343  hdEdxvsPhi->SetTitle("CDC-dEdx vs Phi");
344  hdEdxvsPhi->Draw("col");
345 
346  l_line->DrawLine(-3.20, m_mean, 3.20, m_mean);
347 
348  TPaveText* pinfo0 = new TPaveText(0.60, 0.77, 0.85, 0.89, "NBNDC");
349  setTextStyle(pinfo0);
350  setBEvtInfo(pinfo0);
351  pinfo0->DrawClone("same");
352  }
353 
354  //plot # 2
355  TH2D* hdEdxvsCosth = (TH2D*)findHist("CDCDedx/hdEdxvsCosth");
356  if (hdEdxvsCosth != nullptr) {
357 
358  c_pr_dedx->cd(4);
359 
360  setHistStyle(hdEdxvsCosth);
361  hdEdxvsCosth->SetTitle("CDC-dEdx vs Costh");
362  hdEdxvsCosth->Draw("col");
363 
364  l_line->DrawLine(-1.0, m_mean, 1.0, m_mean);
365 
366  TPaveText* pinfo1 = new TPaveText(0.60, 0.77, 0.85, 0.89, "NBNDC");
367  setTextStyle(pinfo1);
368  setBEvtInfo(pinfo1);
369  pinfo1->DrawClone("same");
370  }
371 
372 }
373 
374 //-----------------------------------------------
376 {
377 
378  TH2D* hinjtimeHer = (TH2D*)findHist("CDCDedx/hinjtimeHer");
379  TH2D* hinjtimeLer = (TH2D*)findHist("CDCDedx/hinjtimeLer");
380 
381  if (hinjtimeHer != nullptr && hinjtimeLer != nullptr) {
382 
383  c_pr_dedx->cd(7);
384 
385  setPlotStyle();
386  setHistStyle(hinjtimeHer);
387  hinjtimeHer->SetTitle("Time since last injection (HER)");
388  hinjtimeHer->Draw("box");
389 
390  setHistStyle(hinjtimeLer);
391  hinjtimeLer->SetFillColor(kRed);
392  hinjtimeLer->Draw("same box");
393 
394  l_line->DrawLine(0, m_mean, 80e3, m_mean);
395 
396  TLegend* lego = new TLegend(0.50, 0.77, 0.60, 0.89);
397  lego->AddEntry(hinjtimeHer, "HER", "f");
398  lego->AddEntry(hinjtimeLer, "LER", "f");
399  lego->Draw("same");
400 
401  TPaveText* pinfo0 = new TPaveText(0.60, 0.77, 0.85, 0.89, "NBNDC");
402  setTextStyle(pinfo0);
403  setBEvtInfo(pinfo0);
404  pinfo0->DrawClone("same");
405  }
406 }
407 
408 //---------------------------------------------
410 {
411 
412  //Injection time variation
413  TH2D* hdEdxITHer = (TH2D*)findHist("CDCDedx/hinjtimeHer");
414  TH2D* hdEdxITLer = (TH2D*)findHist("CDCDedx/hinjtimeLer");
415 
416  if (hdEdxITHer != nullptr && hdEdxITLer != nullptr) {
417 
418  c_pr_dedx->cd(5);
419 
420  int fbin = hdEdxITHer->FindFirstBinAbove(0, 1);
421  int lbin = hdEdxITHer->FindLastBinAbove(0, 1);
422  int nbin = (lbin - fbin + 1) / 2;
423  if (nbin <= 0)nbin = 1;
424 
425  TH1F* hinjectmean[2];
426  TH1F* hinjectsigma[2];
427 
428  for (int i = 0; i < 2; i++) {
429  hinjectmean[i] = new TH1F(Form("hinjectmean%d", i), "", nbin, 0.5, nbin + 0.5);
430  hinjectmean[i]->SetTitle("CDC-dE/dx gain(#mu);Injection time (ms);dE/dx (#mu_{fit})");
431  hinjectsigma[i] = new TH1F(Form("hinjectsigma%d", i), "", nbin, 0.5, nbin + 0.5);
432  hinjectsigma[i]->SetTitle("CDC-dE/dx reso.(#sigma);Injection time (ms);dE/dx (#sigma_{fit})");
433  }
434 
435  setHistPars(hdEdxITHer, hinjectmean[0], hinjectsigma[0], nbin);
436  setHistPars(hdEdxITLer, hinjectmean[1], hinjectsigma[1], nbin);
437 
438  //2 intra-gain trend
439 
440  gPad->SetGridy(1);
441 
442  drawHistPars(hinjectmean[0], nbin, m_mean, 0.40, "#mu_{fit}");
443 
444  hinjectmean[1]->SetMarkerColor(kBlue);
445  hinjectmean[1]->SetMarkerStyle(20);
446  hinjectmean[1]->SetMarkerSize(1.10);
447  hinjectmean[1]->Draw("same");
448 
449  TLegend* lego = new TLegend(0.45, 0.77, 0.60, 0.89);
450  lego->AddEntry(hinjectmean[0], "HER", "p");
451  lego->AddEntry(hinjectmean[1], "LER", "p");
452  lego->Draw("same");
453 
454  //3 intra-resolution trend
455  c_pr_dedx->cd(6);
456  gPad->SetGridy(1);
457  drawHistPars(hinjectsigma[0], nbin, m_sigma, 0.15, "#sigma_{fit}");
458 
459  hinjectsigma[1]->SetMarkerColor(kBlue);
460  hinjectsigma[1]->SetMarkerStyle(20);
461  hinjectsigma[1]->SetMarkerSize(1.10);
462  hinjectsigma[1]->Draw("same");
463 
464  TLegend* lego1 = new TLegend(0.45, 0.77, 0.60, 0.89);
465  lego1->AddEntry(hinjectsigma[0], "HER", "p");
466  lego1->AddEntry(hinjectsigma[1], "LER", "p");
467  lego1->Draw("same");
468  }
469 }
470 
471 //------------------------------------------------
473 {
474 
475  //Draw Scattered plot
476  TH2D* hWires = (TH2D*)findHist("CDCDedx/hWires");
477  TH2D* hWireStatus = (TH2D*)findHist("CDCDedx/hWireStatus");
478  if (hWires != nullptr && hWireStatus != nullptr) {
479 
480  c_pr_dedx->cd(8);
481  setHistStyle(hWires);
482  hWires->SetMarkerColor(kGray);
483  hWires->Draw("");
484 
485  std::string s_ndead = hWireStatus->GetTitle();
486  int m_ndead = atof(s_ndead.c_str());
487  m_monObj->setVariable("CDCDedxDeadWires", m_ndead);
488 
489  setHistStyle(hWireStatus);
490  hWireStatus->SetMarkerColor(kRed);
491  hWireStatus->SetMarkerStyle(7);
492  hWireStatus->Draw("same");
493 
494  TPaveText* pinfo0 = new TPaveText(0.117, 0.832, 0.148, 0.976, "NBNDC");
495  setTextStyle(pinfo0);
496  pinfo0->AddText(Form("CDC Wire Status"));
497  pinfo0->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
498  pinfo0->AddText(Form("Dead: %d (%0.02f%%)", m_ndead, (100.0 * m_ndead / c_nSenseWires)));
499  if (m_nallevt > 1e5)
500  pinfo0->AddText(Form("Events: %0.02fM", double(m_nallevt / 1e6)));
501  if (m_nallevt > 1e3)
502  pinfo0->AddText(Form("Events: %0.02fK", double(m_nallevt / 1e3)));
503  else
504  pinfo0->AddText(Form("Events: %d", m_nallevt));
505  pinfo0->Draw("same");
506  }
507 }
508 
509 //-----------------------------------------------
510 void DQMHistAnalysisCDCDedxModule::setHistPars(TH2D* hdEdx, TH1F* hmean, TH1F* hsigma, int nbin)
511 {
512 
513  int fbin = hdEdx->FindFirstBinAbove(0, 1);
514 
515  for (int ibin = 0; ibin < nbin; ibin++) {
516  int localbin = ibin + fbin;
517  TH1D* hdEdxIRInd = (TH1D*)hdEdx->ProjectionY(Form("htemp_%d", localbin), localbin, localbin);
518 
519  double mean = 0.0, meanerr = 0.0;
520  double sigma = 0.0, sigmaerr = 0.0;
521 
522  fitHistogram(hdEdxIRInd, m_status);
523 
524  if (m_status == "OK") {
525  mean = hdEdxIRInd->GetFunction("f_gaus")->GetParameter(1);
526  meanerr = hdEdxIRInd->GetFunction("f_gaus")->GetParError(1);
527  sigma = hdEdxIRInd->GetFunction("f_gaus")->GetParameter(2);
528  sigmaerr = hdEdxIRInd->GetFunction("f_gaus")->GetParError(2);
529  }
530 
531  hmean->SetBinContent(ibin + 1, mean);
532  hmean->SetBinError(ibin + 1, meanerr);
533  hsigma->SetBinContent(ibin + 1, sigma);
534  hsigma->SetBinError(ibin + 1, sigmaerr);
535  }
536 
537 }
538 
539 //-----------------------------------------------
540 void DQMHistAnalysisCDCDedxModule::drawHistPars(TH1F* hist, int nbin, double pars, double fac, std::string var)
541 {
542  std::string hname = hist->GetName();
543  setPadStyle(0.143, 0.045, 0.077, 0.0);
544  hist->SetMarkerColor(kRed);
545  hist->SetMarkerStyle(20);
546  hist->SetMarkerSize(1.10);
547  hist->GetYaxis()->SetRangeUser(pars - fac, pars + fac); //m_mean - 0.10, m_mean + 0.10);
548  hist->Draw("");
549 
550  l_line->DrawLine(0.5, pars, hist->GetXaxis()->GetBinUpEdge(nbin), pars);
551 
552  TPaveText* pinfo0 = new TPaveText(0.609, 0.680, 0.942, 0.911, "NBNDC");
553  setTextStyle(pinfo0);
554  if (hname == "hdEdxIRMean" || hname == "hdEdxIRSigma") pinfo0->AddText("Intra-run variation");
555  setBEvtInfo(pinfo0);
556  pinfo0->AddText(Form("Avg %s: %0.3f", var.data(), pars));
557  pinfo0->Draw("same");
558 }
559 
560 //-----------------------------------------------
562 {
563 
564  pt->AddText("CDC dE/dx (e^{+}e^{-})");
565  pt->AddText(Form("Exp/Run: %d/%d", m_exp, m_run));
566  if (m_nbhabhaevt > 1e5)
567  pt->AddText(Form("Events: %0.02fM", double(m_nbhabhaevt / 1e6)));
568  if (m_nbhabhaevt > 1e3)
569  pt->AddText(Form("Events: %0.02fK", double(m_nbhabhaevt / 1e3)));
570  else
571  pt->AddText(Form("Events: %d", m_nbhabhaevt));
572 }
573 
574 //----------------------------------------------------------------------------------------
575 void DQMHistAnalysisCDCDedxModule::fitHistogram(TH1D*& temphist, std::string& status)
576 {
577 
578  if (temphist != nullptr) {
579  temphist->GetXaxis()->SetRange(temphist->FindFirstBinAbove(0, 1), temphist->FindLastBinAbove(0, 1));
580  int fs = temphist->Fit(f_gaus, "QR");
581  if (fs != 0) {
582  status = "Failed";
583  } else {
584  double mean = temphist->GetFunction("f_gaus")->GetParameter(1);
585  double width = temphist->GetFunction("f_gaus")->GetParameter(2);
586  temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
587  fs = temphist->Fit(f_gaus, "QR", "", mean - 3.0 * width, mean + 3.0 * width);
588  if (fs != 0)status = "Failed";
589  else {
590  temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
591  status = "OK";
592  }
593  }
594  }
595 
596 }
597 
598 //------------------------------------------------
600 {
601 
602  const Int_t NRGBs = 6;
603  const Int_t NCont = 255;
604  Double_t stops[NRGBs] = { 0.00, 0.00, 0.34, 0.61, 0.84, 1.00 };
605  Double_t red[NRGBs] = { 0.00, 0.00, 0.00, 0.87, 1.00, 0.51 };
606  Double_t green[NRGBs] = { 0.00, 0.00, 0.81, 1.00, 0.20, 0.00 };
607  Double_t blue[NRGBs] = { 0.00, 0.51, 1.00, 0.12, 0.00, 0.00 };
608  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
609  gStyle->SetNumberContours(NCont);
610 
611 }
612 
613 
614 //------------------------------------------------------------------
616 {
617 
618  obj->SetFillColor(0);
619  obj->SetFillStyle(0);
620 
621  obj->SetLineColor(TColor::GetColor("#000000"));
622  obj->SetLineWidth(0);
623  obj->SetTextAlign(12);
624 
625  obj->SetTextColor(kGray + 3);
626  obj->SetTextFont(82);
627  obj->SetTextSize(0.03157895);
628  obj->SetTextAlign(12);
629 
630 }
631 
632 //------------------------------------------------------------------
634 {
635 
636  obj->SetStats(0);
637  obj->SetTitle("");
638  obj->SetFillColor(kYellow);
639 
640  obj->SetTitleOffset(1.15, "x");
641  obj->SetTitleSize(.040, "x");
642  obj->SetTitleOffset(1.15, "y");
643  obj->SetTitleSize(.040, "y");
644 
645  obj->SetLabelOffset(0.015, "x");
646  obj->SetLabelSize(.040, "x");
647  obj->SetLabelOffset(0.015, "y");
648  obj->SetLabelSize(.040, "y");
649 
650  obj->SetTickLength(0.03, "x");
651  obj->SetTickLength(0.02, "y");
652 
653 }
654 
655 //-------------------------------------------------------------------------------------
656 void DQMHistAnalysisCDCDedxModule::setPadStyle(double l, double r, double t, double b)
657 {
658 
659  if (l != 0)gPad->SetLeftMargin(l);
660  if (r != 0)gPad->SetRightMargin(r);
661  if (t != 0)gPad->SetTopMargin(t);
662  if (b != 0)gPad->SetBottomMargin(b);
663  gPad->SetTickx(1);
664  gPad->SetTicky(1);
665 
666 }
void drawDedxInjTime()
function to draw the dEdx vs injection time
void drawDedxCosPhi()
funtion to draw dEdx vs costh and phi
void initialize() override final
init function for default values
void drawBandPlot()
funtion to dedx bands P
void drawWireStatus()
funtion to draw ADC-based dead wire status of CDC
void drawHistPars(TH1F *hist, int nbin, double pars, double fac, std::string var)
funtion to draw the histograms
void setHistStyle(TH1 *obj)
funtion to set the style of histogram
MonitoringObject * m_monObj
MonitoringObject for mirabelle.
void terminate() override final
terminating at the end of last run
void drawDedxInjTimeBin()
funtion to draw the mean/reso of dedx vs injection time
void event() override final
event by event function
TCanvas * c_ir_dedx
intra-run dedx status
void setHistPars(TH2D *hist, TH1F *hmean, TH1F *hsigma, int nbin)
funtion to set the mean and sigma histograms
void endRun() override final
end of each run
void setPadStyle(double l, double r, double t, double b)
funtion to reset pad margins
void getMetadata()
funtion to get metadata from histogram
void setBEvtInfo(TPaveText *pt)
funtion to set the bhabha event info
void beginRun() override final
begin each run
void drawDedxPR()
funtion to draw dEdx+Fit
void fitHistogram(TH1D *&temphist, std::string &status)
funtion to fit gaussian dist.
void drawDedxIR()
funtion to draw dEdx+Fit for run variation
void setTextStyle(TPaveText *&obj)
funtion to add text style
std::string mmode
monitoring mode all/basic
void setPlotStyle()
funtion to add plot style
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.