Belle II Software development
CDCDedxRunGainAlgorithm.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 <reconstruction/calibration/CDCDedxRunGainAlgorithm.h>
10
11#include <TF1.h>
12#include <TCanvas.h>
13#include <TLine.h>
14#include <iostream>
15#include <fstream>
16#include <TH1I.h>
17
18using namespace Belle2;
19
20//-----------------------------------------------------------------
21// Implementation
22//-----------------------------------------------------------------
24 CalibrationAlgorithm("CDCDedxElectronCollector"),
25 isMakePlots(true),
26 isMergePayload(true),
27 fsrun(""),
28 fSigLim(2.5),
29 fdEdxBins(1000),
30 fdEdxMin(0.0),
31 fdEdxMax(5.0),
32 fAdjust(1.00)
33{
34 // Set module properties
35 setDescription("A calibration algorithm for CDC dE/dx run gains");
36}
37
38//-----------------------------------------------------------------
39// Run the calibration
40//-----------------------------------------------------------------
42{
43 B2INFO("Preparing dE/dx calibration for CDC dE/dx run gain");
44
45 const auto expRun = getRunList()[0];
46 updateDBObjPtrs(1, expRun.second, expRun.first);
47 double ExistingRG = m_DBRunGain->getRunGain();
48
49 // Get data objects
50 auto ttree = getObjectPtr<TTree>("tree");
51 if (!ttree)return c_NotEnoughData;
52
53 double dedx;
54 int run = 0, runtemp = 0;
55 TString status = "";
56 ttree->SetBranchAddress("dedx", &dedx);
57 ttree->SetBranchAddress("run", &run);
58
59 TH1D* hDedx = new TH1D("hDedx", "hDedx;dEdx rel (no hadsat);entries", fdEdxBins, fdEdxMin, fdEdxMax);
60 TH1D* hDedxAbs = new TH1D("hDedxAbs", "hDedxAbs;dEdx abs(no hadsat);entries", fdEdxBins, fdEdxMin, fdEdxMax);
61
62 for (int i = 0; i < ttree->GetEntries(); ++i) {
63 ttree->GetEvent(i);
64 if (dedx <= 0) continue;
65 hDedx->Fill(dedx);
66 hDedxAbs->Fill(dedx * ExistingRG);
67 runtemp = run;
68 }
69
70 //min 500 good events otherwise merge to next run
71 if (hDedx->Integral() < 1000) {
72 B2INFO("Not enough data for run (" << runtemp << ") so far only " << hDedx->Integral() << " therefore adding next run");
73 fsrun += Form("%d_", runtemp);
74 delete hDedx;
75 delete hDedxAbs;
76 return c_NotEnoughData;
77 }
78
79 B2INFO("Total track for this run: " << runtemp << " = " << hDedx->Integral());
80 fsrun += Form("%d", runtemp);
81 hDedx->SetName(Form("%s_%s", hDedx->GetName(), fsrun.Data()));
82 hDedxAbs->SetName(Form("%s_%s", hDedxAbs->GetName(), fsrun.Data()));
83
84 //Fit absolute run for the reference only
85 double RGabs_Mean = 1.0;
86 FitGaussianWRange(hDedxAbs, status);
87 if (status != "FitOK") {
88 hDedxAbs->SetTitle(Form("dE/dx abs, Run %s (%s)", fsrun.Data(), status.Data()));
89 } else {
90 RGabs_Mean = hDedxAbs->GetFunction("gaus")->GetParameter(1);
91 double RGabs_MeanErr = hDedxAbs->GetFunction("gaus")->GetParError(1);
92 double RGabs_Sigma = hDedxAbs->GetFunction("gaus")->GetParameter(2);
93 hDedxAbs->SetTitle(Form("dE/dx abs, Run %s (%s), #mu_{fit}: %0.04f#pm%0.04f,, #sigma_{fit}: %0.04f", fsrun.Data(), status.Data(),
94 RGabs_Mean, RGabs_MeanErr, RGabs_Sigma));
95 }
96
97 //Fit relative run for actual calibration
98 double RGrel_Mean = 1.0, RGrel_MeanErr = -99.0, RGrel_Sigma = -99.0;
99 FitGaussianWRange(hDedx, status);
100 if (status != "FitOK") {
101 RGrel_Mean = hDedx->GetMean(); //mean for hist for fit failure
102 hDedx->SetTitle(Form("dE/dx abs, Run %s (%s)", fsrun.Data(), status.Data()));
103 } else {
104 RGrel_Mean = hDedx->GetFunction("gaus")->GetParameter(1);
105 RGrel_MeanErr = hDedx->GetFunction("gaus")->GetParError(1);
106 RGrel_Sigma = hDedx->GetFunction("gaus")->GetParameter(2);
107 hDedx->SetTitle(Form("dE/dx rel, Run %s (%s), #mu_{fit}: %0.04f#pm%0.04f,, #sigma_{fit}: %0.04f", fsrun.Data(), status.Data(),
108 RGrel_Mean, RGrel_MeanErr, RGrel_Sigma));
109 }
110
111 if (RGrel_Mean <= 0.0)RGrel_Mean = 1.0; //safe button
112
113 //for validation purpose
114 if (isMakePlots) {
115
116 std::ofstream fileRGout;
117 fileRGout.open(Form("fRG_Constants_%s.txt", fsrun.Data()));
118 fileRGout << expRun.second << " " << ExistingRG << " " << RGrel_Mean << " " << RGrel_MeanErr << " " << RGrel_Sigma << "\n";
119 fileRGout.close();
120
121 TCanvas* ctmp = new TCanvas("tmp", "tmp", 1000, 500);
122 ctmp->SetBatch(kTRUE);
123 ctmp->Divide(2, 1);
124
125 ctmp->cd(1);
126 hDedx->SetFillColorAlpha(kYellow, 0.30);
127 hDedx->DrawCopy("");
128
129 TLine* tl = new TLine();
130 tl->SetLineColor(kRed);
131 tl->SetX1(RGrel_Mean); tl->SetX2(RGrel_Mean);
132 tl->SetY1(0); tl->SetY2(hDedx->GetMaximum());
133 tl->DrawClone("same");
134
135 ctmp->cd(2);
136 hDedxAbs->SetFillColorAlpha(kBlue, 0.10);
137 hDedxAbs->DrawCopy("");
138
139 TLine* tlAbs = new TLine();
140 tlAbs->SetLineColor(kRed);
141 tlAbs->SetX1(RGabs_Mean); tlAbs->SetX2(RGabs_Mean);
142 tlAbs->SetY1(0); tlAbs->SetY2(hDedxAbs->GetMaximum());
143 tlAbs->DrawClone("same");
144
145 ctmp->Print(Form("cdcdedx_rungain_fits_%s_%s.pdf", fsrun.Data(), status.Data()));
146
147 TCanvas* cstats = new TCanvas("cstats", "cstats", 1000, 500);
148 cstats->SetBatch(kTRUE);
149 cstats->Divide(2, 1);
150 cstats->cd(1);
151 auto hestats = getObjectPtr<TH1I>("hestats");
152 if (hestats)hestats->DrawCopy("");
153 cstats->cd(2);
154 auto htstats = getObjectPtr<TH1I>("htstats");
155 if (htstats)htstats->DrawCopy("");
156 cstats->Print(Form("cdcdedx_rungain_stats_%s.pdf", fsrun.Data()));
157
158 delete ctmp;
159 delete tl;
160 delete tlAbs;
161 delete cstats;
162 }
163
164 B2INFO("Saving new rung for (Exp, Run) : (" << expRun.first << "," << expRun.second << ")");
165 generateNewPayloads(RGrel_Mean, ExistingRG);
166
167 fsrun.Resize(0);
168 delete hDedx;
169 delete hDedxAbs;
170 return c_OK;
171}
172
173
174void CDCDedxRunGainAlgorithm::generateNewPayloads(double RGrel_Mean, double ExistingRG)
175{
176 if (isMergePayload) {
177 B2INFO("--> RunGain: Previous = " << ExistingRG << ", Relative = " << RGrel_Mean << ", Adjustment = " << fAdjust <<
178 ", Merged (saved) = " << RGrel_Mean*
179 ExistingRG * fAdjust);
180 RGrel_Mean *= ExistingRG;
181 RGrel_Mean *= fAdjust; //default is 1.0
182 } else {
183 B2INFO("--> RunGain: Previous = " << ExistingRG << ", Relative (saved) = " << RGrel_Mean);
184 }
185
186 CDCDedxRunGain* gain = new CDCDedxRunGain(RGrel_Mean);
187 saveCalibration(gain, "CDCDedxRunGain");
188}
189
190
191void CDCDedxRunGainAlgorithm::FitGaussianWRange(TH1D*& temphist, TString& status)
192{
193
194 double histmean = temphist->GetMean();
195 double histrms = temphist->GetRMS();
196 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
197
198 int fs = temphist->Fit("gaus", "Q0");
199 if (fs != 0) {
200 B2INFO(Form("\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
201 status = "FitFailed";
202 return;
203 } else {
204 double mean = temphist->GetFunction("gaus")->GetParameter(1);
205 double width = temphist->GetFunction("gaus")->GetParameter(2);
206 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
207 fs = temphist->Fit("gaus", "QR", "", mean - fSigLim * width, mean + fSigLim * width);
208 if (fs != 0) {
209 B2INFO(Form("\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
210 status = "FitFailed";
211 return;
212 } else {
213 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
214 B2INFO(Form("\tFit for hist (%s) successful (status = %d)", temphist->GetName(), fs));
215 status = "FitOK";
216 }
217 }
218}
CDCDedxRunGainAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
DBObjPtr< CDCDedxRunGain > m_DBRunGain
Run gain DB object.
void generateNewPayloads(double RunGainConst, double ExistingRG)
function to store new payload after full calibration
double fdEdxMax
max dedx range for gain cal
bool isMergePayload
merge payload at the of calibration
double fAdjust
factor to adjust dedx gain
virtual EResult calibrate() override
Run algorithm.
double fdEdxMin
min dedx range for gain cal
void FitGaussianWRange(TH1D *&temphist, TString &status)
function to perform fit run by run
TString fsrun
flag to identify low stats runs
int fdEdxBins
number of bins for dedx histogram
bool isMakePlots
produce plots for status
double fSigLim
fit range limit based on sigma
dE/dx run gain calibration constants
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
Abstract base class for different kinds of events.