Belle II Software development
CDCDedxInjectTimeAlgorithm.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 <cdc/calibration/CDCdEdx/CDCDedxInjectTimeAlgorithm.h>
10
11using namespace Belle2;
12
13//-----------------------------------------------------------------
14// Implementation
15//-----------------------------------------------------------------
17 CalibrationAlgorithm("CDCDedxElectronCollector"),
18 m_sigmaR(2.0),
19 m_dedxBins(250),
20 m_dedxMin(0.0),
21 m_dedxMax(2.5),
22 m_chiBins(250),
23 m_chiMin(-10.0),
24 m_chiMax(10.0),
25 m_countR(0),
26 m_thersE(1000),
27 m_isminStat(false),
28 m_ismakePlots(true),
29 m_isMerge(true),
30 m_prefix("cdcdedx_injcal"),
31 m_suffix("")
32{
33 // Set module properties
34 setDescription("A calibration algorithm for CDC dE/dx injection time gain and reso");
35}
36
37//-----------------------------------------------------------------
38// Run the calibration
39//-----------------------------------------------------------------
41{
42
44
45 //existing inject time payload for merging
46 if (!m_DBInjectTime.isValid())
47 B2FATAL("There is no valid payload for Injection time");
48
49 // Get data objects
50 auto ttree = getObjectPtr<TTree>("tree");
51 if (!ttree) return c_NotEnoughData;
52
53 double dedx = 0.0, injtime = 0.0, injring = 1.0, costh, mom;
54 int nhits;
55 ttree->SetBranchAddress("dedx", &dedx);
56 ttree->SetBranchAddress("injtime", &injtime);
57 ttree->SetBranchAddress("injring", &injring);
58 ttree->SetBranchAddress("costh", &costh);
59 ttree->SetBranchAddress("p", &mom);
60 ttree->SetBranchAddress("nhits", &nhits);
61
62 //way to define time bins/edges only once
63 if (m_countR == 0) {
64 defineTimeBins(); //returns m_vtlocaledges
65 m_tbins = m_vtlocaledges.size() - 1;
67 m_countR++;
68 }
69
70 TH1D* htimes = new TH1D(Form("htimes_%s", m_suffix.data()), "", m_tbins, m_tedges);
71
72 //time bins are changeable from out so vector is used
73 std::array<std::vector<TH1D*>, numdedx::nrings> hdedx, htime, hchi, hdedx_corr;
74 defineHisto(hdedx, "dedx");
75 defineHisto(hdedx_corr, "dedx_corr");
76 defineHisto(htime, "timeinj");
77 defineHisto(hchi, "chi");
78
79 const double tzedges[4] = {0, 2.0e3, 0.1e6, 20e6};
80 std::array<std::array<TH1D*, 3>, numdedx::nrings> hztime;
81 defineTimeHisto(hztime);
82
83 //fill histograms
84 for (int i = 0; i < ttree->GetEntries(); ++i) {
85
86 ttree->GetEvent(i);
87 if (dedx <= 0 || injtime < 0 || injring < 0) continue;
88
89 //add larger times to the last bin
90 if (injtime > m_tedges[m_tbins]) injtime = m_tedges[m_tbins] - 10.0;
91
92 //injection ring
93 int wr = 0;
94 if (injring > 0.5) wr = 1;
95
96 //injection time bin
97 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
98 if (tb > m_tbins)tb = m_tbins; //overflow
99 tb = tb - 1;
100
101 htimes->Fill(injtime);
102 if (injtime < tzedges[1]) hztime[wr][0]->Fill(injtime);
103 else if (injtime < tzedges[2]) hztime[wr][1]->Fill(injtime);
104 else hztime[wr][2]->Fill(injtime);
105
106 hdedx[wr][tb]->Fill(dedx);
107 htime[wr][tb]->Fill(injtime);
108 }
109
110 //keep merging runs to achieve enough stats
111 m_isminStat = false;
112 checkStatistics(hdedx);
113 if (m_isminStat) {
114 deleteHisto(htime);
115 deleteHisto(hdedx);
116 deleteHisto(hdedx_corr);
117 deleteHisto(hchi);
118 deleteTimeHisto(hztime);
119 delete htimes;
120 return c_NotEnoughData;
121 }
122
123 //clear vector of existing constants
124 std::map<int, std::vector<double>> vmeans, vresos, vtimes;
125
126 // get time vector
127 for (unsigned int ir = 0; ir < c_rings; ir++) {
128 for (unsigned int it = 0; it < m_tbins; it++) {
129 double avgtime = htime[ir][it]->GetMean();
130 double avgtimeerr = htime[ir][it]->GetMeanError();
131 vtimes[ir * 2].push_back(avgtime);
132 vtimes[ir * 2 + 1].push_back(avgtimeerr);
133 }
134 }
135
136 //Fit dedx to get mean and resolution
137 getMeanReso(hdedx, vmeans, vresos);
138
139 //Bin-bias correction to mean
140 std::map<int, std::vector<double>> vmeanscorr;
141 correctBinBias(vmeanscorr, vmeans, vtimes, htimes);
142
143 //scale the mean and merge with old constants
144 std::array<double, numdedx::nrings> scale;
145 std::map<int, std::vector<double>> vmeanscal;
146 createPayload(scale, vmeanscorr, vmeanscal, "mean");
147
148 //................................................
149 // Do the calibration for resolution
150 //................................................
151 CDCDedxMeanPred mbg;
153 mbg.setParameters();
154 sbg.setParameters();
155
156 //fill histograms for the resolution
157 for (int i = 0; i < ttree->GetEntries(); ++i) {
158
159 ttree->GetEvent(i);
160 if (dedx <= 0 || injtime < 0 || injring < 0) continue;
161
162 double corrcetion = getCorrection(injring, injtime, vmeanscal);
163 double old_cor = m_DBInjectTime->getCorrection("mean", injring, injtime);
164
165 dedx = (dedx * old_cor) / corrcetion;
166 //add larger times to the last bin
167 if (injtime > m_tedges[m_tbins]) injtime = m_tedges[m_tbins] - 10.0;
168
169 //injection ring
170 int wr = 0;
171 if (injring > 0.5) wr = 1;
172
173 //injection time bin
174 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
175 if (tb > m_tbins)tb = m_tbins; //overflow
176 tb = tb - 1;
177
178 double predmean = mbg.getMean(mom / Const::electron.getMass());
179 double predsigma = sbg.nhitPrediction(nhits) * sbg.ionzPrediction(dedx) * sbg.cosPrediction(costh);
180
181 double chi = (dedx - predmean) / predsigma;
182 hdedx_corr[wr][tb]->Fill(dedx);
183 hchi[wr][tb]->Fill(chi);
184 }
185
186 // fit chi to get mean and resolution
187 std::map<int, std::vector<double>> vmeans_chi, vresos_chi;
188 getMeanReso(hchi, vmeans_chi, vresos_chi);
189
190 //bin-bias correction to the resolution
191 std::map<int, std::vector<double>> vresoscorr;
192 correctBinBias(vresoscorr, vresos_chi, vtimes, htimes);
193
194 // scale the resolution
195 std::map<int, std::vector<double>> vresoscal;
196 std::array<double, numdedx::nrings> scale_reso;
197 createPayload(scale_reso, vresoscorr, vresoscal, "reso");
198
199 //Fit the corrected mean to check for consistency
200 std::map<int, std::vector<double>> vmeans_corr, vresos_corr;
201 getMeanReso(hdedx_corr, vmeans_corr, vresos_corr);
202
203 //................................................
204 //preparing final payload
205 //................................................
206 m_vinjPayload.clear();
207 m_vinjPayload.reserve(6);
208 for (int ir = 0; ir < 2; ir++) {
209 m_vinjPayload.push_back(m_vtlocaledges);
210 m_vinjPayload.push_back(vmeanscal[ir * 2]);
211 m_vinjPayload.push_back(vresoscal[ir * 2]);
212 }
213
214 if (m_ismakePlots) {
215
216 //0 plot event track statistics
218
219 //1 plot injection time
220 plotInjectionTime(hztime);
221
222 //2. Draw dedxfits
223 plotBinLevelDist(hdedx, "dedxfits");
224
225 //3. Draw chifits
226 plotBinLevelDist(hchi, "chifits");
227
228 //4. Draw timedist
229 plotBinLevelDist(htime, "timedist");
230
231 //5. plot relative const., bias-bias corrected for dedx
232 plotRelConstants(vmeans, vresos, vmeanscorr, "dedx");
233
234 //6. plot relative const., bias-bias corrected for chi
235 plotRelConstants(vmeans_chi, vresos_chi, vresoscorr, "chi");
236
237 //7. plot mean and resolution of corrected dedx to check for consistency
238 plotRelConstants(vmeans_corr, vresos_corr, vresoscorr, "dedx_corr");
239
240 //8. plot time statistics dist
241 plotTimeStat(htime);
242
243 //9. plot final merged const. and comparison to old
244 plotFinalConstants(vmeanscal, vresoscal, scale, scale_reso);
245 }
246
247 //saving payloads;
249 saveCalibration(gains, "CDCDedxInjectionTime");
250 B2INFO("dE/dx Injection time calibration done");
251
252 //delete all histograms
253 deleteHisto(htime);
254 deleteHisto(hdedx);
255 deleteHisto(hdedx_corr);
256 deleteHisto(hchi);
257 deleteTimeHisto(hztime);
258 delete htimes;
259
260 B2INFO("Saving calibration for: " << m_suffix << "");
261 return c_OK;
262}
263
264//------------------------------------
265void CDCDedxInjectTimeAlgorithm::getMeanReso(std::array<std::vector<TH1D*>, numdedx::nrings>& hvar,
266 std::map<int, std::vector<double>>& vmeans, std::map<int, std::vector<double>>& vresos)
267{
268 for (unsigned int ir = 0; ir < c_rings; ir++) {
269
270 for (unsigned int it = 0; it < m_tbins; it++) {
271 double mean = 1.00, meanerr = 0.0;
272 double reso = 1.00, resoerr = 0.0;
273 if (hvar[ir][it]->Integral() > 250) {
274 fstatus status;
275 fitGaussianWRange(hvar[ir][it], status);
276 if (status != fitOK) {
277 mean = hvar[ir][it]->GetMean();
278 hvar[ir][it]->SetTitle(Form("%s, (%d)", hvar[ir][it]->GetTitle(), status));
279 } else {
280 mean = hvar[ir][it]->GetFunction("gaus")->GetParameter(1);
281 meanerr = hvar[ir][it]->GetFunction("gaus")->GetParError(1);
282 reso = hvar[ir][it]->GetFunction("gaus")->GetParameter(2);
283 resoerr = hvar[ir][it]->GetFunction("gaus")->GetParError(2);
284 std::string title = Form("#mu_{fit}: %0.03f, #sigma_{fit}: %0.03f", mean, reso);
285 hvar[ir][it]->SetTitle(Form("%s, %s", hvar[ir][it]->GetTitle(), title.data()));
286 }
287 }
288
289 vmeans[ir * 2].push_back(mean);
290 vresos[ir * 2].push_back(reso);
291 vmeans[ir * 2 + 1].push_back(meanerr);
292 vresos[ir * 2 + 1].push_back(resoerr);
293 }
294 }
295}
296
297//----------------------------------------
298void CDCDedxInjectTimeAlgorithm::fitGaussianWRange(TH1D*& temphist, fstatus& status)
299{
300 double histmean = temphist->GetMean();
301 double histrms = temphist->GetRMS();
302 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
303
304 int fs = temphist->Fit("gaus", "Q0");
305 if (fs != 0) {
306 B2INFO(Form("\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
307 status = fitFailed;
308 return;
309 } else {
310 double mean = temphist->GetFunction("gaus")->GetParameter(1);
311 double width = temphist->GetFunction("gaus")->GetParameter(2);
312 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
313 fs = temphist->Fit("gaus", "QR", "", mean - m_sigmaR * width, mean + m_sigmaR * width);
314 if (fs != 0) {
315 B2INFO(Form("\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
316 status = fitFailed;
317 return;
318 } else {
319 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
320 B2INFO(Form("\tFit for hist (%s) successful (status = %d)", temphist->GetName(), fs));
321 status = fitOK;
322 }
323 }
324}
325
326//------------------------------------
328{
329 int cruns = 0;
330 for (auto expRun : getRunList()) {
331 if (cruns == 0)B2INFO("CDCDedxInjectTimeAlgorithm: start exp " << expRun.first << " and run " << expRun.second << "");
332 cruns++;
333 }
334
335 const auto erStart = getRunList()[0];
336 int estart = erStart.first;
337 int rstart = erStart.second;
338
339 const auto erEnd = getRunList()[cruns - 1];
340 int rend = erEnd.second;
341
342 updateDBObjPtrs(1, erStart.second, erStart.first);
343
344 if (m_isminStat) {
345 m_suffix = Form("e%dr%dto%d_nruns%d", estart, rstart, rend, cruns);
346 B2INFO("\t+ run = " << rend << ", m_suffix = " << m_suffix << "");
347 } else {
348 m_suffix = Form("e%dr%d", estart, rstart);
349 B2INFO("tool run = " << estart << ", exp = " << estart << ", m_suffix = " << m_suffix << "");
350 }
351}
352
353//------------------------------------
355{
356 //empty local vector or find a way to execulate this function
357 //only once
358 if (!m_vtlocaledges.empty()) m_vtlocaledges.clear();
359 if (m_vtedges.empty()) {
360 double fixedges[69];
361 for (int ib = 0; ib < 69; ib++) {
362 fixedges[ib] = ib * 0.5 * 1e3;
363 if (ib > 40 && ib <= 60) fixedges[ib] = fixedges[ib - 1] + 1.0 * 1e3;
364 else if (ib > 60 && ib <= 64) fixedges[ib] = fixedges[ib - 1] + 10.0 * 1e3;
365 else if (ib > 64 && ib <= 65) fixedges[ib] = fixedges[ib - 1] + 420.0 * 1e3;
366 else if (ib > 65 && ib <= 66) fixedges[ib] = fixedges[ib - 1] + 500.0 * 1e3;
367 else if (ib > 66) fixedges[ib] = fixedges[ib - 1] + 2e6;
368 m_vtlocaledges.push_back(fixedges[ib]);
369 }
370 } else {
371 for (unsigned int ib = 0; ib < m_vtedges.size(); ib++)
372 m_vtlocaledges.push_back(m_vtedges.at(ib));
373 }
374}
375
376//------------------------------------
377void CDCDedxInjectTimeAlgorithm::defineHisto(std::array<std::vector<TH1D*>, numdedx::nrings>& htemp, std::string var)
378{
379 for (unsigned int ir = 0; ir < c_rings; ir++) {
380 htemp[ir].resize(m_tbins);
381 for (unsigned int it = 0; it < m_tbins; it++) {
382 std::string label = getTimeBinLabel(m_tedges[it], it);
383 std::string title = Form("%s(%s), time(%s)", m_suffix.data(), m_sring[ir].data(), label.data());
384 std::string hname = Form("h%s_%s_%s_t%d", var.data(), m_sring[ir].data(), m_suffix.data(), it);
385 if (var == "dedx" or var == "dedx_corr") htemp[ir][it] = new TH1D(hname.data(), "", m_dedxBins, m_dedxMin, m_dedxMax);
386 else if (var == "chi") htemp[ir][it] = new TH1D(hname.data(), "", m_chiBins, m_chiMin, m_chiMax);
387 else htemp[ir][it] = new TH1D(hname.data(), "", 50, m_tedges[it], m_tedges[it + 1]);
388 htemp[ir][it]->SetTitle(Form("%s;%s;entries", title.data(), var.data()));
389 }
390 }
391}
392
393//------------------------------------
394void CDCDedxInjectTimeAlgorithm::defineTimeHisto(std::array<std::array<TH1D*, 3>, numdedx::nrings>& htemp)
395{
396 const int tzoom = 3;
397 const int nt[tzoom] = {50, 500, 1000};
398 double tzedges[tzoom + 1] = {0, 2.0e3, 0.1e6, 20e6};
399 std::string stname[tzoom] = {"early", "mid", "later"};
400 std::string stlabel[tzoom] = {"zoom <2ms", "early time <= 100ms", "later time >100ms"};
401 for (unsigned int ir = 0; ir < c_rings; ir++) {
402 for (int wt = 0; wt < tzoom; wt++) {
403 std::string title = Form("inject time (%s), %s (%s)", stlabel[wt].data(), m_sring[ir].data(), m_suffix.data());
404 std::string hname = Form("htimezoom_%s_%s_%s", m_sring[ir].data(), stname[wt].data(), m_suffix.data());
405 htemp[ir][wt] = new TH1D(Form("%s", hname.data()), "", nt[wt], tzedges[wt], tzedges[wt + 1]);
406 htemp[ir][wt]->SetTitle(Form("%s;injection time(#mu-sec);entries", title.data()));
407 }
408 }
409}
410
411//------------------------------------
412void CDCDedxInjectTimeAlgorithm::checkStatistics(std::array<std::vector<TH1D*>, numdedx::nrings>& hvar)
413{
414 for (unsigned int ir = 0; ir < c_rings; ir++) {
415 for (unsigned int it = 3; it < m_tbins; it++) {
416 //check statiscs from 1-40ms
417 if (m_tedges[it] < 4e4 && hvar[ir][it]->Integral() < m_thersE) {
418 m_isminStat = true;
419 break;
420 } else continue;
421 }
422 }
423}
424
425//------------------------------------
426void CDCDedxInjectTimeAlgorithm::correctBinBias(std::map<int, std::vector<double>>& varcorr,
427 std::map<int, std::vector<double>>& var,
428 std::map<int, std::vector<double>>& time, TH1D*& htimes)
429{
430 //Deep copy OK
431 varcorr = var;
432
433 for (int ir = 0; ir < 2; ir++) {
434
435 for (int ix = varcorr[ir * 2].size(); ix -- > 0;) {
436
437 double var_thisbin = 1.0;
438 var_thisbin = var[ir * 2].at(ix);
439
440 double atime_thisbin = time[ir * 2].at(ix);
441 double ctime_thisbin = htimes->GetBinCenter(ix + 1);
442
443 if (atime_thisbin > 0 && atime_thisbin < 4e4 * 0.99) {
444
445 double var_nextbin = 1.0;
446 var_nextbin = var[ir * 2].at(ix + 1);
447 double var_diff = var_nextbin - var_thisbin;
448
449 double atime_nextbin = time[ir * 2].at(ix + 1);
450 double atime_diff = atime_nextbin - atime_thisbin;
451
452 double slope = (atime_diff > 0) ? var_diff / atime_diff : -1.0;
453
454 //extrapolation after veto only
455 if (var_diff > 0 && slope > 0)varcorr[ir * 2].at(ix) = var_thisbin + (ctime_thisbin - atime_thisbin) * (slope);
456 printf("\t %s ix = %d, center = %0.2f(%0.3f), var = %0.5f(%0.5f) \n", m_sring[ir].data(), ix, ctime_thisbin, atime_thisbin,
457 var_thisbin, varcorr[ir * 2].at(ix));
458 } else {
459 printf("\t %s --> ix = %d, center = %0.2f(%0.3f), var = %0.5f(%0.5f) \n", m_sring[ir].data(), ix, ctime_thisbin, atime_thisbin,
460 var_thisbin, varcorr[ir * 2].at(ix));
461 }
462 }
463 }
464}
465
466//-------------------------------------
467void CDCDedxInjectTimeAlgorithm::createPayload(std::array<double, numdedx::nrings>& scale,
468 std::map<int, std::vector<double>>& var,
469 std::map<int, std::vector<double>>& varscal, std::string svar)
470{
471 varscal = var;
472
473 B2INFO("CDCDedxInjectTimeAlgorithm: normalising constants with plateau");
474 for (unsigned int ir = 0; ir < c_rings; ir++) {
475 //scaling means with time >40ms
476 unsigned int msize = varscal[ir * 2].size();
477 int countsum = 0;
478 scale[ir] = 0;
479 for (unsigned int im = 0; im < msize; im++) {
480 double time = m_vtlocaledges.at(im);
481 double mean = varscal[ir * 2].at(im);
482 if (time > 4e4 && mean > 0) {
483 scale[ir] += mean;
484 countsum++;
485 }
486 }
487 if (countsum > 0 && scale[ir] > 0) {
488 scale[ir] /= countsum;
489 for (unsigned int im = 0; im < msize; im++) {
490 varscal[ir * 2].at(im) /= scale[ir];
491 }
492 }
493 }
494 if (m_isMerge && svar == "mean") {
495 //merge only no change in payload structure
496 bool incomp_bin = false;
497 std::vector<std::vector<double>> oldvectors;
498 if (m_DBInjectTime) oldvectors = m_DBInjectTime->getConstVector();
499 int vsize = oldvectors.size();
500 if (vsize != 6) incomp_bin = true;
501 else {
502 for (int iv = 0; iv < 2; iv++) {
503 if (oldvectors[iv * 3 + 1].size() != varscal[iv * 2].size()) incomp_bin = true;
504 }
505 }
506 if (!incomp_bin) {
507 B2INFO("CDCDedxInjectTimeAlgorithm: started merging relative constants");
508 for (int ir = 0; ir < 2; ir++) {//merging only means
509 unsigned int msize = varscal[ir * 2].size();
510 for (unsigned int im = 0; im < msize; im++) {
511 double relvalue = varscal[ir * 2].at(im);
512 double oldvalue = oldvectors[ir * 3 + 1].at(im);
513 double merged = oldvalue * relvalue;
514 printf("%s: rel %0.03f, old %0.03f, merged %0.03f\n", m_suffix.data(), relvalue, oldvalue, merged);
515 varscal[ir * 2].at(im) *= oldvectors[ir * 3 + 1].at(im) ;
516 }
517 }
518 } else B2ERROR("CDCDedxInjectTimeAlgorithm: found incompatible bins for merging");
519 } else B2INFO("CDCDedxInjectTimeAlgorithm: saving final (abs) calibration");
520}
521
522//------------------------------------
523void CDCDedxInjectTimeAlgorithm::plotBinLevelDist(std::array<std::vector<TH1D*>, numdedx::nrings>& hvar, std::string var)
524{
525 TCanvas cfit("cfit", "cfit", 1000, 500);
526 cfit.Divide(2, 1);
527 std::stringstream psname_fit;
528 psname_fit << Form("%s_%s_%s.pdf[", m_prefix.data(), var.data(), m_suffix.data());
529 cfit.Print(psname_fit.str().c_str());
530 psname_fit.str("");
531 psname_fit << Form("%s_%s_%s.pdf", m_prefix.data(), var.data(), m_suffix.data());
532 for (unsigned int it = 0; it < m_tbins; it++) {
533 for (unsigned int ir = 0; ir < c_rings; ir++) {
534 cfit.cd(ir + 1);
535 hvar[ir][it]->SetFillColorAlpha(ir + 5, 0.25);
536 hvar[ir][it]->Draw();
537 }
538 cfit.Print(psname_fit.str().c_str());
539 }
540 psname_fit.str("");
541 psname_fit << Form("%s_%s_%s.pdf]", m_prefix.data(), var.data(), m_suffix.data());
542 cfit.Print(psname_fit.str().c_str());
543}
544
545//----------------------------------------
547{
548 // draw event and track statistics
549 TCanvas cestat("cestat", "cestat", 1000, 500);
550 cestat.SetBatch(kTRUE);
551 cestat.Divide(2, 1);
552
553 cestat.cd(1);
554 auto hestats = getObjectPtr<TH1I>("hestats");
555 if (hestats) {
556 hestats->SetName(Form("hestats_%s", m_suffix.data()));
557 hestats->SetStats(0);
558 hestats->Draw("hist text");
559 }
560 cestat.cd(2);
561 auto htstats = getObjectPtr<TH1I>("htstats");
562 if (htstats) {
563 htstats->SetName(Form("htstats_%s", m_suffix.data()));
564 htstats->SetStats(0);
565 htstats->Draw("hist text");
566 }
567 cestat.Print(Form("%s_eventstat_%s.pdf", m_prefix.data(), m_suffix.data()));
568}
569
570//------------------------------------
571void CDCDedxInjectTimeAlgorithm::plotInjectionTime(std::array<std::array<TH1D*, 3>, numdedx::nrings>& hvar)
572{
573 TCanvas ctzoom("ctzoom", "ctzoom", 1500, 450);
574 ctzoom.SetBatch(kTRUE);
575 ctzoom.Divide(3, 1);
576 for (int wt = 0; wt < 3; wt++) {
577 ctzoom.cd(wt + 1);
578 if (wt == 2) gPad->SetLogy();
579 for (unsigned int ir = 0; ir < c_rings; ir++) {
580 hvar[ir][wt]->SetStats(0);
581 hvar[ir][wt]->SetFillColorAlpha(5 + ir, 0.20);
582 if (ir == 0) {
583 double max1 = hvar[ir][wt]->GetMaximum();
584 double max2 = hvar[c_rings - 1][wt]->GetMaximum();
585 if (max2 > max1) hvar[ir][wt]->SetMaximum(max2 * 1.05);
586 hvar[ir][wt]->Draw("");
587 } else hvar[ir][wt]->Draw("same");
588 }
589 }
590 ctzoom.Print(Form("%s_timezoom_%s.pdf]", m_prefix.data(), m_suffix.data()));
591}
592
593//------------------------------------
594void CDCDedxInjectTimeAlgorithm::plotRelConstants(std::map<int, std::vector<double>>& vmeans,
595 std::map<int, std::vector<double>>& vresos, std::map<int, std::vector<double>>& corr, std::string svar)
596{
597 std::string sname[3] = {"mean", "reso"};
598 const int lcolors[c_rings] = {2, 4};
599
600 TH1D* hmean[c_rings], *hcorr[c_rings];
601 TH1D* hreso[c_rings];
602
603 TCanvas* cconst[2];
604 for (int ic = 0; ic < 2; ic++) {
605 cconst[ic] = new TCanvas(Form("c%sconst", sname[ic].data()), "", 900, 500);
606 cconst[ic]->SetGridy(1);
607 }
608
609 TLegend* mleg = new TLegend(0.50, 0.54, 0.80, 0.75, NULL, "brNDC");
610 mleg->SetBorderSize(0);
611 mleg->SetFillStyle(0);
612
613 for (unsigned int ir = 0; ir < c_rings; ir++) {
614
615 std::string mtitle = Form("#mu(dedx), relative const. compare, (%s)", m_suffix.data());
616 hmean[ir] = new TH1D(Form("hmean_%s_%s", m_suffix.data(), m_sring[ir].data()), "", m_tbins, 0, m_tbins);
617
618 std::string rtitle = Form("#sigma(#chi), relative const. compare, (%s)", m_suffix.data());
619 hreso[ir] = new TH1D(Form("hreso_%s_%s", m_suffix.data(), m_sring[ir].data()), "", m_tbins, 0, m_tbins);
620
621 hcorr[ir] = new TH1D(Form("hcorr_%s_%s", m_suffix.data(), m_sring[ir].data()), "", m_tbins, 0, m_tbins);
622 if (svar == "chi") {
623 hmean[ir]->SetTitle(Form("%s;injection time(#mu-second);#mu (#chi-fit)", rtitle.data()));
624 hreso[ir]->SetTitle(Form("%s;injection time(#mu-second);#sigma (#chi-fit)", rtitle.data()));
625 hcorr[ir]->SetTitle(Form("%s;injection time(#mu-second);#sigma (#chi-fit bin-bais-corr)", rtitle.data()));
626
627 } else {
628 hmean[ir]->SetTitle(Form("%s;injection time(#mu-second);#mu (dedx-fit)", mtitle.data()));
629 hreso[ir]->SetTitle(Form("%s;injection time(#mu-second);#sigma (dedx-fit)", mtitle.data()));
630 hcorr[ir]->SetTitle(Form("%s;injection time(#mu-second);#mu (dedx-fit, bin-bais-corr)", mtitle.data()));
631 }
632
633
634 for (unsigned int it = 0; it < m_tbins; it++) {
635
636 std::string label = getTimeBinLabel(m_tedges[it], it);
637
638 hmean[ir]->SetBinContent(it + 1, vmeans[ir * 2].at(it));
639 hmean[ir]->SetBinError(it + 1, vmeans[ir * 2 + 1].at(it));
640 hmean[ir]->GetXaxis()->SetBinLabel(it + 1, label.data());
641
642 hreso[ir]->SetBinContent(it + 1, vresos[ir * 2].at(it));
643 hreso[ir]->SetBinError(it + 1, vresos[ir * 2 + 1].at(it));
644 hreso[ir]->GetXaxis()->SetBinLabel(it + 1, label.data());
645
646 hcorr[ir]->SetBinContent(it + 1, corr[ir * 2].at(it));
647 hcorr[ir]->SetBinError(it + 1, corr[ir * 2 + 1].at(it));
648 hcorr[ir]->GetXaxis()->SetBinLabel(it + 1, label.data());
649 }
650
651 mleg->AddEntry(hmean[ir], Form("%s", m_sring[ir].data()), "lep");
652 cconst[0]->cd();
653 if (svar == "chi") setHistStyle(hmean[ir], lcolors[ir], ir + 24, -0.60, 0.60);
654 else if (svar == "dedx") setHistStyle(hmean[ir], lcolors[ir], ir + 24, 0.60, 1.10);
655 else setHistStyle(hmean[ir], lcolors[ir], ir + 24, 0.9, 1.10);
656
657 if (ir == 0)hmean[ir]->Draw("");
658 else hmean[ir]->Draw("same");
659 if (svar == "dedx") {
660 mleg->AddEntry(hcorr[ir], Form("%s (bin-bias-corr)", m_sring[ir].data()), "lep");
661 setHistStyle(hcorr[ir], lcolors[ir], ir + 20, 0.60, 1.10);
662 hcorr[ir]->Draw("same");
663 }
664 if (ir == 1)mleg->Draw("same");
665
666 cconst[1]->cd();
667 if (svar == "chi") setHistStyle(hreso[ir], lcolors[ir], ir + 24, 0.5, 1.50);
668 else setHistStyle(hreso[ir], lcolors[ir], ir + 24, 0.0, 0.15);
669 if (ir == 0)hreso[ir]->Draw("");
670 else hreso[ir]->Draw("same");
671 if (svar == "chi") {
672 mleg->AddEntry(hcorr[ir], Form("%s (bin-bias-corr)", m_sring[ir].data()), "lep");
673 setHistStyle(hcorr[ir], lcolors[ir], ir + 20, 0.5, 1.50);
674 hcorr[ir]->Draw("same");
675 }
676 if (ir == 1)mleg->Draw("same");
677 }
678
679 for (int ic = 0; ic < 2; ic++) {
680 cconst[ic]->SaveAs(Form("%s_relconst_%s_%s_%s.pdf", m_prefix.data(), sname[ic].data(), svar.data(), m_suffix.data()));
681 cconst[ic]->SaveAs(Form("%s_relconst_%s_%s_%s.root", m_prefix.data(), sname[ic].data(), svar.data(), m_suffix.data()));
682 delete cconst[ic];
683 }
684 for (int ic = 0; ic < 2; ic++) {
685
686 delete hmean[ic];
687 delete hreso[ic];
688 delete hcorr[ic];
689 }
690}
691
692//------------------------------------
693void CDCDedxInjectTimeAlgorithm::plotTimeStat(std::array<std::vector<TH1D*>, numdedx::nrings>& htime)
694{
695 const int lcolors[c_rings] = {2, 4};
696
697 TH1D* htimestat[c_rings];
698
699 TCanvas* cconst = new TCanvas("ctimestatconst", "", 900, 500);
700 cconst->SetGridy(1);
701
702 TLegend* rleg = new TLegend(0.40, 0.60, 0.80, 0.72, NULL, "brNDC");
703 rleg->SetBorderSize(0);
704 rleg->SetFillStyle(0);
705
706 for (unsigned int ir = 0; ir < c_rings; ir++) {
707
708 std::string title = Form("injection time, her-ler comparison, (%s)", m_suffix.data());
709 htimestat[ir] = new TH1D(Form("htimestat_%s_%s", m_suffix.data(), m_sring[ir].data()), "", m_tbins, 0, m_tbins);
710 htimestat[ir]->SetTitle(Form("%s;injection time(#mu-second);norm. entries", title.data()));
711
712 for (unsigned int it = 0; it < m_tbins; it++) {
713 std::string label = getTimeBinLabel(m_tedges[it], it);
714 htimestat[ir]->SetBinContent(it + 1, htime[ir][it]->Integral());
715 htimestat[ir]->SetBinError(it + 1, 0);
716 htimestat[ir]->GetXaxis()->SetBinLabel(it + 1, label.data());
717 }
718
719 cconst->cd();
720 double norm = htimestat[ir]->GetMaximum();
721 rleg->AddEntry(htimestat[ir], Form("%s (scaled with %0.02f)", m_sring[ir].data(), norm), "lep");
722 htimestat[ir]->Scale(1.0 / norm);
723 setHistStyle(htimestat[ir], lcolors[ir], ir + 20, 0.0, 1.10);
724 htimestat[ir]->SetFillColorAlpha(lcolors[ir], 0.30);
725 if (ir == 0) htimestat[ir]->Draw("hist");
726 else htimestat[ir]->Draw("hist same");
727 if (ir == 1)rleg->Draw("same");
728 }
729
730 cconst->SaveAs(Form("%s_relconst_timestat_%s.pdf", m_prefix.data(), m_suffix.data()));
731 cconst->SaveAs(Form("%s_relconst_timestat_%s.root", m_prefix.data(), m_suffix.data()));
732 delete cconst;
733}
734
735//------------------------------------
736void CDCDedxInjectTimeAlgorithm::plotFinalConstants(std::map<int, std::vector<double>>& vmeans,
737 std::map<int, std::vector<double>>& vresos,
738 std::array<double, numdedx::nrings>& scale, std::array<double, numdedx::nrings>& scale_reso)
739{
740
741 std::vector<std::vector<double>> oldvectors;
742 if (m_DBInjectTime)oldvectors = m_DBInjectTime->getConstVector();
743
744 const int c_type = 2; //old and new
745 std::string sname[c_rings] = {"mean", "reso"};
746 std::string stype[c_type] = {"new", "old"};
747 const int lcolors[c_rings] = {2, 4};
748 const int lmarker[c_type] = {20, 24}; //+2 for different rings
749
750 TH1D* hmean[c_rings][c_type], *hreso[c_rings][c_type];
751
752 TCanvas* cconst[c_rings];
753 for (int ic = 0; ic < 2; ic++) {
754 cconst[ic] = new TCanvas(Form("c%sconst", sname[ic].data()), "", 900, 500);
755 cconst[ic]->SetGridy(1);
756 }
757
758 TLegend* mleg = new TLegend(0.50, 0.54, 0.80, 0.75, NULL, "brNDC");
759 mleg->SetBorderSize(0);
760 mleg->SetFillStyle(0);
761
762 TLegend* rleg = new TLegend(0.50, 0.54, 0.80, 0.75, NULL, "brNDC");
763 rleg->SetBorderSize(0);
764 rleg->SetFillStyle(0);
765
766 for (unsigned int ip = 0; ip < c_type; ip++) {
767
768 for (unsigned int ir = 0; ir < c_rings; ir++) {
769
770 std::string hname = Form("hfmean_%s_%s_%s", m_sring[ir].data(), stype[ip].data(), m_suffix.data());
771 hmean[ir][ip] = new TH1D(hname.data(), "", m_tbins, 0, m_tbins);
772 std::string title = Form("#mu(dedx), final-mean-compare (%s)", m_suffix.data());
773 hmean[ir][ip]->SetTitle(Form("%s;injection time(#mu-second);#mu (dedx-fit)", title.data()));
774
775 hname = Form("hfreso_%s_%s_%s", m_sring[ir].data(), stype[ip].data(), m_suffix.data());
776 hreso[ir][ip] = new TH1D(hname.data(), "", m_tbins, 0, m_tbins);
777 title = Form("#sigma(#chi), final-reso-compare (%s)", m_suffix.data());
778 hreso[ir][ip]->SetTitle(Form("%s;injection time(#mu-second);#sigma (#chi-fit)", title.data()));
779
780 for (unsigned int it = 0; it < m_tbins; it++) {
781
782 std::string label = getTimeBinLabel(m_tedges[it], it);
783 double mean = 0.0, reso = 0.0;
784 if (ip == 0) {
785 mean = m_vinjPayload[ir * 3 + 1].at(it);
786 //reso is reso/mu (reso is relative so mean needs to be relative)
787 reso = m_vinjPayload[ir * 3 + 2].at(it);
788
789 } else {
790 mean = oldvectors[ir * 3 + 1].at(it);
791 reso = oldvectors[ir * 3 + 2].at(it);
792
793 }
794
795 //old payloads
796 hmean[ir][ip]->SetBinContent(it + 1, mean);
797 hmean[ir][ip]->SetBinError(it + 1, vmeans[ir * 2 + 1].at(it));
798 hmean[ir][ip]->GetXaxis()->SetBinLabel(it + 1, label.data());
799
800 hreso[ir][ip]->SetBinContent(it + 1, reso);
801 hreso[ir][ip]->SetBinError(it + 1, vresos[ir * 2 + 1].at(it));
802 hreso[ir][ip]->GetXaxis()->SetBinLabel(it + 1, label.data());
803 }
804
805 cconst[0]->cd();
806 if (ip == 1)mleg->AddEntry(hmean[ir][ip], Form("%s, %s", m_sring[ir].data(), stype[ip].data()), "lep");
807 else mleg->AddEntry(hmean[ir][ip], Form("%s, %s (scaled by %0.03f)", m_sring[ir].data(), stype[ip].data(), scale[ir]), "lep");
808 setHistStyle(hmean[ir][ip], lcolors[ir], lmarker[ip] + ir * 2, 0.60, 1.05);
809 if (ir == 0 && ip == 0)hmean[ir][ip]->Draw("");
810 else hmean[ir][ip]->Draw("same");
811 if (ir == 1 && ip == 1)mleg->Draw("same");
812
813 cconst[1]->cd();
814 if (ip == 1)rleg->AddEntry(hreso[ir][ip], Form("%s, %s", m_sring[ir].data(), stype[ip].data()), "lep");
815 else rleg->AddEntry(hreso[ir][ip], Form("%s, %s (scaled by %0.03f)", m_sring[ir].data(), stype[ip].data(), scale_reso[ir]), "lep");
816 setHistStyle(hreso[ir][ip], lcolors[ir], lmarker[ip] + ir * 2, 0.6, 1.9);
817 if (ir == 0 && ip == 0)hreso[ir][ip]->Draw("");
818 else hreso[ir][ip]->Draw("same");
819 if (ir == 1 && ip == 1)rleg->Draw("same");
820 }
821 }
822
823 for (int ic = 0; ic < 2; ic++) {
824 cconst[ic]->SaveAs(Form("%s_finalconst_%s_%s.pdf", m_prefix.data(), sname[ic].data(), m_suffix.data()));
825 cconst[ic]->SaveAs(Form("%s_finalconst_%s_%s.root", m_prefix.data(), sname[ic].data(), m_suffix.data()));
826 delete cconst[ic];
827 }
828
829 for (unsigned int ir = 0; ir < c_rings; ir++) {
830 for (unsigned int ip = 0; ip < c_type; ip++) {
831 delete hmean[ir][ip];
832 delete hreso[ir][ip];
833 }
834 }
835}
836
837//------------------------------------
838double CDCDedxInjectTimeAlgorithm::getCorrection(unsigned int ring, unsigned int time,
839 std::map<int, std::vector<double>>& vmeans)
840{
841
842 unsigned int iv = ring * 2;
843
844 unsigned int sizet = m_vtlocaledges.size(); //time
845
846 std::vector<unsigned int> tedges(sizet); //time edges array
847 std::copy(m_vtlocaledges.begin(), m_vtlocaledges.end(), tedges.begin());
848
849 if (time >= 5e6) time = 5e6 - 10;
850 unsigned int it = m_DBInjectTime->getTimeBin(tedges, time);
851
852 double center = 0.5 * (m_vtlocaledges.at(it) + m_vtlocaledges.at(it + 1));
853
854 //no corr before veto bin (usually one or two starting bin)
855 //intrapolation for entire range except
856 //--extrapolation (for first half and last half of intended bin)
857 int thisbin = it, nextbin = it;
858 if (center != time && it > 0) {
859
860 if (time < center) {
861 thisbin = it - 1;
862 } else {
863 if (it < sizet - 2)nextbin = it + 1;
864 else thisbin = it - 1;
865 }
866
867 if (it <= 2) {
868 double diff = vmeans[iv].at(2) - vmeans[iv].at(1) ;
869 if (diff < -0.015) { //difference above 1.0%
870 thisbin = it;
871 if (it == 1) nextbin = it;
872 else nextbin = it + 1;
873 } else {
874 if (it == 1) {
875 thisbin = it;
876 nextbin = it + 1;
877 }
878 }
879 }
880 }
881
882 double thisdedx = vmeans[iv].at(thisbin);
883 double nextdedx = vmeans[iv].at(nextbin);
884
885 double thistime = 0.5 * (m_vtlocaledges.at(thisbin) + m_vtlocaledges.at(thisbin + 1));
886 double nexttime = 0.5 * (m_vtlocaledges.at(nextbin) + m_vtlocaledges.at(nextbin + 1));
887
888 double newdedx = vmeans[iv].at(it);
889 if (thisbin != nextbin)
890 newdedx = thisdedx + ((nextdedx - thisdedx) / (nexttime - thistime)) * (time - thistime);
891
892 return newdedx;
893}
bool m_isMerge
merge payload when rel constant
void plotBinLevelDist(std::array< std::vector< TH1D * >, numdedx::nrings > &hvar, std::string var)
function to draw dedx, chi and time dist.
std::vector< double > m_vtlocaledges
internal time vector
double m_sigmaR
fit dedx dist in sigma range
void correctBinBias(std::map< int, std::vector< double > > &varcorr, std::map< int, std::vector< double > > &var, std::map< int, std::vector< double > > &time, TH1D *&htimes)
function to correct dedx mean/reso and return corrected vector map
void defineHisto(std::array< std::vector< TH1D * >, numdedx::nrings > &htemp, std::string var)
function to define histograms for dedx and time dist.
double getCorrection(unsigned int ring, unsigned int time, std::map< int, std::vector< double > > &vmeans)
function to get the correction factor of mean
std::string m_prefix
string prefix for plot names
void getExpRunInfo()
function to get exp/run information (payload object, plotting)
void setHistStyle(TH1D *&htemp, const int ic, const int is, const double min, const double max)
function to set histogram cosmetics
void plotFinalConstants(std::map< int, std::vector< double > > &vmeans, std::map< int, std::vector< double > > &vresos, std::array< double, numdedx::nrings > &scale, std::array< double, numdedx::nrings > &scale_reso)
function to final constant from merging or abs fits
int m_countR
a hack for running functions once
double * m_tedges
internal time array (copy of vtlocaledges)
bool m_ismakePlots
produce plots for monitoring
std::string m_suffix
string suffix for object names
static const int c_rings
injection ring constants
void createPayload(std::array< double, numdedx::nrings > &scale, std::map< int, std::vector< double > > &vmeans, std::map< int, std::vector< double > > &varscal, std::string svar)
function to store payloads after full calibration
void plotInjectionTime(std::array< std::array< TH1D *, 3 >, numdedx::nrings > &hvar)
function to injection time distributions (HER/LER in three bins)
void plotTimeStat(std::array< std::vector< TH1D * >, numdedx::nrings > &htime)
function to draw time stats
void fitGaussianWRange(TH1D *&temphist, fstatus &status)
function to perform gauss fit for input histogram
void defineTimeBins()
function to set/reset time bins
void plotRelConstants(std::map< int, std::vector< double > > &vmeans, std::map< int, std::vector< double > > &vresos, std::map< int, std::vector< double > > &corr, std::string svar)
function to relative constant from dedx fit mean and chi fit reso
void defineTimeHisto(std::array< std::array< TH1D *, 3 >, numdedx::nrings > &htemp)
function to define injection time bins histograms (monitoring only)
void deleteHisto(std::array< std::vector< TH1D * >, numdedx::nrings > &htemp)
function to delete histograms for dedx and time dist.
void deleteTimeHisto(std::array< std::array< TH1D *, 3 >, numdedx::nrings > &htemp)
function to define injection time bins histograms (monitoring only)
void plotEventStats()
function to draw event/track statistics plots
virtual EResult calibrate() override
CDC dE/dx Injection time algorithm.
void checkStatistics(std::array< std::vector< TH1D * >, numdedx::nrings > &hvar)
check statistics for obtaining calibration const.
std::vector< double > m_vtedges
external time vector
std::vector< std::vector< double > > m_vinjPayload
vector to store payload values
bool m_isminStat
flag to merge runs for statistics thershold
CDCDedxInjectTimeAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
std::array< std::string, numdedx::nrings > m_sring
injection ring name
DBObjPtr< CDCDedxInjectionTime > m_DBInjectTime
Injection time DB object.
unsigned int m_tbins
internal time bins
int m_thersE
min tracks to start calibration
void getMeanReso(std::array< std::vector< TH1D * >, numdedx::nrings > &hvar, std::map< int, std::vector< double > > &vmeans, std::map< int, std::vector< double > > &vresos)
function to get mean and reso of histogram
std::string getTimeBinLabel(const double &tedges, const int &it)
function to return time label for histograms labeling
dE/dx injection time calibration constants
Class to hold the prediction of mean as a function of beta-gamma (bg)
double getMean(double bg)
Return the predicted mean value as a function of beta-gamma (bg)
void setParameters(std::string infile)
set the parameters from file
Class to hold the prediction of resolution depending dE/dx, nhit, and cos(theta)
double ionzPrediction(double dedx)
Return sigma from the ionization parameterization.
double cosPrediction(double cos)
Return sigma from the cos parameterization.
double nhitPrediction(double nhit)
Return sigma from the nhit parameterization.
void setParameters(std::string infile)
set the parameters from file
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_NotEnoughData
Needs more data =2 in Python.
static const ChargedStable electron
electron particle
Definition: Const.h:659
Abstract base class for different kinds of events.