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