Belle II Software  release-05-01-25
DQMHistAnalysisCDCDedx.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2018 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Jitendra Kumar *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <dqm/analysis/modules/DQMHistAnalysisCDCDedx.h>
12 #include <boost/format.hpp>
13 #include <TROOT.h>
14 
15 #include <TF1.h>
16 #include <TCanvas.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 using boost::format;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(DQMHistAnalysisCDCDedx)
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
33  dedxmean(0.0),
34  dedxsigma(0.0)
35 {
36  //Parameter definition
37  B2DEBUG(20, "DQMHistAnalysisCDCDedx: Constructor done.");
38 }
39 
40 
41 DQMHistAnalysisCDCDedxModule::~DQMHistAnalysisCDCDedxModule() { }
42 
43 void DQMHistAnalysisCDCDedxModule::initialize()
44 {
45  gROOT->cd(); // this seems to be important, or strange things happen
46 
47  B2DEBUG(1, "DQMHistAnalysisCDCDedx: initialized.");
48 
49  runstatus = "";
50  dedxmean = 0.0;
51  dedxsigma = 0.0;
52 
53  tLine = new TLine();
54  tLine->SetLineColor(4);
55  tLine->SetLineStyle(9);
56 
57  c_CDCdedxSigma = new TCanvas("CDCDedx/c_CDCdedxSigma");
58 
59  c_CDCdedxMean = new TCanvas("CDCDedx/c_CDCdedxMean");
60  h_CDCdedxMean = new TH1F("CDCDedx/h_CDCdedxMean", "dEdx distribution", 200, 0.0, 2.0);
61  h_CDCdedxMean->SetDirectory(0);
62  h_CDCdedxMean->SetStats(false);
63 
64  f_fGaus = new TF1("f_Gaus", "gaus", 0.0, 2.0);
65  f_fGaus->SetParameter(1, 1.00);
66  f_fGaus->SetParameter(2, 0.06);
67  f_fGaus->SetLineColor(kRed);
68 }
69 
70 
71 void DQMHistAnalysisCDCDedxModule::beginRun()
72 {
73  B2DEBUG(1, "DQMHistAnalysisCDCDedx: beginRun called.");
74 }
75 
76 void DQMHistAnalysisCDCDedxModule::event()
77 {
78  computedEdxMeanSigma();
79 }
80 
81 
82 void DQMHistAnalysisCDCDedxModule::computedEdxMeanSigma()
83 {
84 
85  runstatus.clear();
86  tLine->Clear();
87 
88  TH1* hh1 = findHist("CDCDedx/hdEdx_PerRun");
89  if (hh1 != NULL) {
90  c_CDCdedxMean->Clear();
91  c_CDCdedxMean->cd();
92  //Copy bins and X ranges as well as bin-content (Clone and Copy histogram not working here)
93  h_CDCdedxMean->SetBins(int(hh1->GetXaxis()->GetNbins()), double(hh1->GetXaxis()->GetXmin()), double(hh1->GetXaxis()->GetXmax()));
94  h_CDCdedxMean->SetTitle(hh1->GetTitle());
95  for (Int_t i = 1; i <= hh1->GetXaxis()->GetNbins(); i++) {
96  h_CDCdedxMean->SetBinContent(i, hh1->GetBinContent(i));
97  h_CDCdedxMean->SetBinError(i, hh1->GetBinError(i));
98  }
99 
100  h_CDCdedxMean->Fit(f_fGaus, "Q");
101 
102  runnumber = h_CDCdedxMean->GetTitle();
103  runnumber = runnumber.substr(12, 5);
104 
105  if (h_CDCdedxMean->GetEntries() < 100) {
106  dedxmean = 0.0; dedxsigma = 0.0;
107  runstatus = "Low Stats";
108  } else {
109  if (h_CDCdedxMean->GetFunction("f_Gaus") == NULL || !h_CDCdedxMean->GetFunction("f_Gaus")->IsValid()) {
110  dedxmean = 0.0; dedxsigma = 0.0;
111  runstatus = "Fit Failed";
112  } else {
113  dedxmean = h_CDCdedxMean->GetFunction("f_Gaus")->GetParameter(1);
114  dedxsigma = h_CDCdedxMean->GetFunction("f_Gaus")->GetParameter(2);
115  runstatus = "Good";
116  }
117  }
118 
119  h_CDCdedxMean->GetXaxis()->SetRangeUser(0.0, 2.0);
120  h_CDCdedxMean->SetTitle(Form("Run ## %s, Status: %s, mean: %0.04f, sigma: %0.04f", runnumber.data(), runstatus.data(), dedxmean,
121  dedxsigma));
122  h_CDCdedxMean->SetFillColor(kYellow);
123  h_CDCdedxMean->Draw("hist");
124 
125  f_fGaus->Draw("same");
126 
127  tLine->SetX1(dedxmean);
128  tLine->SetX2(dedxmean);
129  tLine->SetY1(0);
130  tLine->SetY2(h_CDCdedxMean->GetMaximum());
131  tLine->DrawClone("same");
132 
133  c_CDCdedxMean->Modified();
134  }
135 }
136 
137 
138 void DQMHistAnalysisCDCDedxModule::endRun()
139 {
140  B2DEBUG(20, "DQMHistAnalysisCDCDedx : endRun called");
141 }
142 
143 
144 void DQMHistAnalysisCDCDedxModule::terminate()
145 {
146  B2DEBUG(20, "terminate called");
147 }
148 
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::DQMHistAnalysisCDCDedxModule
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisCDCDedx.h:32
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27