Belle II Software development
TimeWalkCalibrationAlgorithm.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#include <cdc/calibration/TimeWalkCalibrationAlgorithm.h>
9#include <cdc/dbobjects/CDCTimeWalks.h>
10#include <cdc/geometry/CDCGeometryPar.h>
11#include <cdc/dataobjects/WireID.h>
12
13#include <TF1.h>
14#include <TFile.h>
15#include <TDirectory.h>
16#include <TROOT.h>
17#include <TTree.h>
18#include <TStopwatch.h>
19
20#include <framework/database/DBObjPtr.h>
21#include <framework/database/IntervalOfValidity.h>
22#include <framework/logging/Logger.h>
23
24using namespace std;
25using namespace Belle2;
26using namespace CDC;
27
29{
31 " -------------------------- Time Walk Calibration Algorithm -------------------------\n"
32 );
33}
34
36{
37 B2INFO("Creating and filling histograms");
38
40
41 double halfCSize[56];
42 for (int i = 0; i < 56; ++i) {
43 double R = cdcgeo.senseWireR(i);
44 double nW = cdcgeo.nWiresInLayer(i);
45 halfCSize[i] = M_PI * R / nW;
46 }
47
48 // Histogram for each board
49 for (int i = 0; i < 300; ++i) {
50 m_h2[i] = new TH2D(Form("board_%d", i), Form("board_%d", i), 50, 0., 500, 60, -30, 30);
51 }
52
53 // Read data
54
55 auto tree = getObjectPtr<TTree>("tree");
56
57 Float_t x;
58 Float_t t_mea;
59 Float_t t_fit;
60 Float_t ndf;
61 Float_t Pval;
62 UShort_t adc;
63 UShort_t IWire;
64 UChar_t lay;
65
66 tree->SetBranchAddress("lay", &lay);
67 tree->SetBranchAddress("IWire", &IWire);
68 tree->SetBranchAddress("x_u", &x);
69 tree->SetBranchAddress("t", &t_mea);
70 tree->SetBranchAddress("t_fit", &t_fit);
71 tree->SetBranchAddress("ndf", &ndf);
72 tree->SetBranchAddress("Pval", &Pval);
73 tree->SetBranchAddress("adc", &adc);
74
75 /* Disable unused branch */
76 std::vector<TString> list_vars = {"lay", "IWire", "x_u", "t", "t_fit", "Pval", "ndf", "adc"};
77 tree->SetBranchStatus("*", 0);
78
79 for (TString brname : list_vars) {
80 tree->SetBranchStatus(brname, 1);
81 }
82
83 const Long64_t nEntries = tree->GetEntries();
84 B2INFO("Number of entries: " << nEntries);
85 TStopwatch time;
86 time.Start();
87 for (Long64_t i = 0; i < nEntries; ++i) {
88 tree->GetEntry(i);
89 const double xmax = halfCSize[lay] - 0.12;
90 if ((fabs(x) < m_xmin) || (fabs(x) > xmax)
91 || (ndf < m_minNdf)
92 || (Pval < m_minPval)) continue; /*select good region*/
93
94 m_h2[cdcgeo.getBoardID(WireID(lay, IWire))]->Fill(adc, fabs(t_mea) - fabs(t_fit));
95 }
96 time.Stop();
97 B2INFO("Time to fill histograms: " << time.RealTime() << "s");
98
99}
100
102{
103 B2INFO("Start calibration");
104 gROOT->SetBatch(1);
105
106 const auto exprun = getRunList()[0];
107 B2INFO("ExpRun used for DB Geometry : " << exprun.first << " " << exprun.second);
108 updateDBObjPtrs(1, exprun.second, exprun.first);
109 B2INFO("Creating CDCGeometryPar object");
111
112 prepare();
113 createHisto();
114
115 TF1* fold = nullptr;
116 if (m_twParamMode == 0)
117 fold = new TF1("fold", "[0]/sqrt(x)", 0, 500);
118 else if (m_twParamMode == 1)
119 fold = new TF1("fold", "[0]*exp(-1*[1]*x)", 0, 500);
120
121 if (fold == nullptr) {
122 B2FATAL("Old fitting function is not defined.");
123 }
124
125
126 B2INFO("time walk formular: ");
127 [](TF1 * f) { auto title = f->GetTitle(); B2INFO("Expression : " << title);}(fold);
128 // B2INFO("New time walk mode : " << m_twParamMode_new << " with " << m_nTwParams_new << " parameters");
129
130 for (int ib = 1; ib < 300; ++ib) {
131 m_flag[ib] = 1;
132 B2DEBUG(21, "Board ID:" << ib);
133 m_h2[ib]->SetDirectory(0);
134
135 // Ignore if histogram has low stat. (<500 entries)
136 if (m_h2[ib]->GetEntries() < 500) {
137 m_flag[ib] = 0;
138 continue;
139 }
140 m_h2[ib]->FitSlicesY(0, 1, -1, 10);
141 TString name = m_h2[ib]->GetName();
142 TString hm_name = name + "_1";
143 m_h1[ib] = (TH1D*)gDirectory->Get(hm_name);
144 if (!m_h1[ib]) {m_flag[ib] = 0; continue;}
145 m_h1[ib]->SetDirectory(0);
146 if (m_h1[ib]->GetEntries() < 5) {
147 m_flag[ib] = 0;
148 B2WARNING("Low statistic, number of points after slice fit: " << m_h1[ib]->GetEntries());
149 continue;
150 }
151
152 // Add previous correction to this
153 for (int p = 0; p < m_nTwParams; ++p) {
154 fold->SetParameter(p, m_tw_old[ib][p]);
155 }
156
157
158 if (m_twParamMode == 0) {
159 TF1* func = new TF1("func", "[0]+[1]/sqrt(x)", 0, 500);
160 func->SetParameters(-4, 28);
161 m_h1[ib]->Fit("func", "MQ", "", 20, 150);
162 } else if (m_twParamMode == 1) {
163 m_h1[ib]->Add(fold);
165 } else {
166 B2FATAL("Mode " << m_twParamMode << " is not available, please check again");
167 }
168
169 // Read fitted parameters
170 TF1* f1 = m_h1[ib]->GetFunction("func");
171 if (!f1) { m_flag[ib] = 0; continue;}
172 m_constTerm[ib] = f1->GetParameter(0);
173 for (int i = 1; i <= m_nTwParams; ++i) {
174 m_tw_new[ib][i - 1] = f1->GetParameter(i);
175 }
176
177
178 B2DEBUG(21, "Prob of fitting:" << f1->GetProb());
179 B2DEBUG(21, "Fitting Param 0-1:" << f1->GetParameter(0) << " - " << f1->GetParameter(1));
180
181 }
182
183 //Write histogram to file
184 if (m_storeHisto) {
185 storeHist();
186 }
187
188 write();
189 updateT0();
190
191 return checkConvergence();
192}
193
195{
196 B2INFO("Storing histogram");
197 B2DEBUG(21, "Store 1D histogram");
198 TFile* fhist = new TFile(m_histName.c_str(), "RECREATE");
199 auto hNDF = getObjectPtr<TH1F>("hNDF");
200 auto hPval = getObjectPtr<TH1F>("hPval");
201 auto hEvtT0 = getObjectPtr<TH1F>("hEventT0");
202 //store NDF, P-val. EventT0 histogram for monitoring during calibration
203 if (hNDF && hPval && hEvtT0) {
204 hEvtT0->Write();
205 hPval->Write();
206 hNDF->Write();
207 }
208
209 TDirectory* old = gDirectory;
210 TDirectory* h1D = old->mkdir("h1D");
211 TDirectory* h2D = old->mkdir("h2D");
212 h1D->cd();
213 for (int ib = 1; ib < 300; ++ib) {
214 if (!m_h1[ib] || m_flag[ib] != 1) continue;
215 if (m_h1[ib]->GetEntries() < 5) continue;
216 m_h1[ib]->SetMinimum(-5);
217 m_h1[ib]->SetMaximum(15);
218 m_h1[ib]->Write();
219 }
220
221 B2DEBUG(21, "Store 2D histogram");
222 h2D->cd();
223 for (int ib = 1; ib < 300; ++ib) {
224 if (m_h2[ib] == nullptr) continue;
225 if (m_h2[ib]->GetEntries() < 5) continue;
226 m_h2[ib]->Write();
227
228 }
229
230 fhist->Close();
231 B2INFO("Hitograms were stored");
232}
233
235{
236 TH1F* hDtw = new TH1F("hDtw", "", 100, -1, 1);
237 for (int ib = 0; ib < 300; ++ib) {
238 float dtw = (m_tw_new[ib][0] - m_tw_old[ib][0]) / m_tw_old[ib][0];
239 if (std::isnan(dtw) == 0) {
240 hDtw->Fill(dtw);
241 }
242 }
243
244 B2INFO(hDtw->GetRMS());
245
246 if (hDtw->GetRMS() < 0.02) {
247 return c_OK;
248 } else {
249 return c_Iterate;
250 }
251 delete hDtw;
252}
253
255{
256 B2INFO("Save to the local DB");
257 CDCTimeWalks* dbTw = new CDCTimeWalks();
258 int nfailure = 0;
259
261 for (int ib = 0; ib < 300; ++ib) {
262 if (m_flag[ib] != 1) {
263 nfailure += 1;
264 }
265 if (m_twParamMode == 0) {
266 const int num = static_cast<int>(m_tw_old[ib].size());
267 for (int i = 0; i < num; ++i) {
268 m_tw_new[ib][i] += m_tw_old[ib][i];
269 }
270 }
271
272 dbTw->setTimeWalkParams(ib, m_tw_new[ib]);
273 }
274
275 if (m_textOutput == true) {
277 B2RESULT("Time-walk coeff. table has been written to " << m_outputTWFileName.c_str());
278 }
279
280 saveCalibration(dbTw, "CDCTimeWalks");
281 B2RESULT("Failure to calibrate time-walk for " << nfailure << " board");
282
283}
284
286{
287 B2INFO("Prepare calibration");
288
290 m_twParamMode = dbTw->getTwParamMode();
291
292 if (!(m_twParamMode == 0 || m_twParamMode == 1)) {
293 B2FATAL("Mode " << m_twParamMode << " haven't implemented yet.");
294 }
295
296 B2INFO("tw param mode " << m_twParamMode);
297 const int nEntries = dbTw->getEntries();
298 for (int ib = 0; ib < nEntries; ++ib) {
299 m_tw_old[ib] = dbTw->getTimeWalkParams(ib);
300 m_tw_new[ib].resize(m_nTwParams, 0.0);
301 }
302}
303
305{
306 B2INFO("Add constant term into T0 database");
308 CDCTimeZeros* tz = new CDCTimeZeros();
309 double T0;
310 for (int ilay = 0; ilay < 56; ++ilay) {
311 const unsigned int nW = cdcgeo.nWiresInLayer(ilay);
312 for (unsigned int iwire = 0; iwire < nW; ++iwire) {
313 WireID wireid(ilay, iwire);
314 int bID = cdcgeo.getBoardID(wireid);
315 T0 = cdcgeo.getT0(wireid);
316 if (m_flag[bID] == 1) {
317 tz->setT0(wireid, T0 - m_constTerm[bID]);
318 } else {
319 tz->setT0(wireid, T0);
320 }
321 }
322 }
323 if (m_textOutput == true) {
325 }
326 saveCalibration(tz, "CDCTimeZeros");
327}
328
330{
331 h1->SetDirectory(0);
332 int max = h1->GetMaximumBin();
333 double maxX = h1->GetBinCenter(max);
334 double maxY = h1->GetBinContent(max);
335 B2DEBUG(21, "Max: id - x - y : " << max << " " << maxX << " " << maxY);
336
337 //search for p0
338 double p0 = -1;
339 h1->Fit("pol0", "MQE", "", maxX + 125, 500);
340 if (h1->GetFunction("pol0")) {
341 p0 = h1->GetFunction("pol0")->GetParameter(0);
342 }
343 //create histo = old-p0;
344 // fit with expo function to find intial parameters
345 TH1D* hshift = new TH1D("hshift", "shift", h1->GetNbinsX(), 0, 500);
346 hshift->SetDirectory(0);
347 for (int i = 0; i <= h1->GetNbinsX(); ++i) {
348 hshift->SetBinContent(i, h1->GetBinContent(i) - p0);
349 }
350 hshift->Fit("expo", "MQE", "", 0, 400);
351 double p1 = maxY + p0;
352 double p2 = -0.04;
353 if (hshift->GetFunction("expo")) {
354 p1 = exp(hshift->GetFunction("expo")->GetParameter(0));
355 p2 = hshift->GetFunction("expo")->GetParameter(1);
356 }
357
358 // fit with function
359 TF1* func = new TF1("func", "[0]+ [1]*exp(-1*[2]*x)", 0, 500);
360 func->SetParameters(p0, p1, -1 * p2);
361 func->SetParLimits(0, -5, 5);
362 func->SetParLimits(1, -5, 500);
363 func->SetParLimits(2, -0.01, 0.1);
364 h1->Fit("func", "MQE", "", 10, 400);
365}
double R
typedef autogenerated by FFTW
Database object for time-walk.
Definition: CDCTimeWalks.h:25
void outputToFile(std::string fileName) const
Output the contents in text file format.
Definition: CDCTimeWalks.h:134
void setTwParamMode(unsigned short mode)
Set tw parameterization mode mode=0: tw (in ns) = p0/sqrt(FADCsum); mode=1: tw( in ns) = p0*exp(-p1*F...
Definition: CDCTimeWalks.h:38
void setTimeWalkParams(unsigned short boardID, const std::vector< float > &params)
Set the time-walk coefficients in the list.
Definition: CDCTimeWalks.h:48
Database object for timing offset (t0).
Definition: CDCTimeZeros.h:26
void outputToFile(std::string fileName) const
Output the contents in text file format.
Definition: CDCTimeZeros.h:123
void setT0(unsigned short iCLayer, unsigned short iWire, float t0)
Set t0 in the list.
Definition: CDCTimeZeros.h:40
The Class for CDC Geometry Parameters.
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
float getT0(const WireID &wireID) const
Returns t0 parameter of the specified sense wire.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
double m_xmin
minimum value cut of drift length.
std::string m_outputT0FileName
t0 file name after calibration.
TH2D * m_h2[300]
2D histogram of residual vs ADC for each board.
unsigned short m_twParamMode
=0 for P0/Sqrt(ADC); =1 for P0*Exp(-P1*ADC).
double m_constTerm[300]
const term in fitting, it will be added to T0 instead tw
double m_minNdf
minimum number of degree of freedom required for track.
unsigned short m_flag[300]
flag for fit status
double m_minPval
minimum number of Prob(chi2) of fitted track.
void updateT0()
update constant term to t0 database.
std::vector< float > m_tw_old[300]
tw list before calibration.
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
std::string m_outputTWFileName
tw file name after calibration.
std::vector< float > m_tw_new[300]
tw list after calibration.
TH1D * m_h1[300]
Mean of residual as function of ADC of each board.
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 successfully =0 in Python.
@ c_Iterate
Needs iteration =1 in Python.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Class to identify a wire inside the CDC.
Definition: WireID.h:34
Abstract base class for different kinds of events.
STL namespace.