9#include <reconstruction/calibration/CDCdEdx/CDCDedxInjectTimeAlgorithm.h>
30 m_prefix(
"cdcdedx_injcal"),
34 setDescription(
"A calibration algorithm for CDC dE/dx injection time gain and reso");
47 B2FATAL(
"There is no valid payload for Injection time");
50 auto ttree = getObjectPtr<TTree>(
"tree");
53 double dedx = 0.0, injtime = 0.0, injring = 1.0, costh, mom;
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);
73 std::array<std::vector<TH1D*>, numdedx::nrings> hdedx, htime, hchi, hdedx_corr;
79 const double tzedges[4] = {0, 2.0e3, 0.1e6, 20e6};
80 std::array<std::array<TH1D*, 3>, numdedx::nrings> hztime;
84 for (
int i = 0; i < ttree->GetEntries(); ++i) {
87 if (dedx <= 0 || injtime < 0 || injring < 0)
continue;
94 if (injring > 0.5) wr = 1;
97 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
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);
106 hdedx[wr][tb]->Fill(dedx);
107 htime[wr][tb]->Fill(injtime);
124 std::map<int, std::vector<double>> vmeans, vresos, vtimes;
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);
140 std::map<int, std::vector<double>> vmeanscorr;
144 std::array<double, numdedx::nrings> scale;
145 std::map<int, std::vector<double>> vmeanscal;
154 for (
int i = 0; i < ttree->GetEntries(); ++i) {
157 if (dedx <= 0 || injtime < 0 || injring < 0)
continue;
159 double corrcetion =
getCorrection(injring, injtime, vmeanscal);
160 double old_cor =
m_DBInjectTime->getCorrection(
"mean", injring, injtime);
162 dedx = (dedx * old_cor) / corrcetion;
168 if (injring > 0.5) wr = 1;
171 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
176 double predsigma = s.nhitPrediction(nhits) * s.ionzPrediction(dedx) * s.cosPrediction(costh);
178 double chi = (dedx - predmean) / predsigma;
179 hdedx_corr[wr][tb]->Fill(dedx);
180 hchi[wr][tb]->Fill(chi);
184 std::map<int, std::vector<double>> vmeans_chi, vresos_chi;
188 std::map<int, std::vector<double>> vresoscorr;
192 std::map<int, std::vector<double>> vresoscal;
193 std::array<double, numdedx::nrings> scale_reso;
197 std::map<int, std::vector<double>> vmeans_corr, vresos_corr;
205 for (
int ir = 0; ir < 2; ir++) {
247 B2INFO(
"dE/dx Injection time calibration done");
257 B2INFO(
"Saving calibration for: " <<
m_suffix <<
"");
263 std::map<
int, std::vector<double>>& vmeans, std::map<
int, std::vector<double>>& vresos)
265 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
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) {
273 if (status != fitOK) {
274 mean = hvar[ir][it]->GetMean();
275 hvar[ir][it]->SetTitle(Form(
"%s, (%d)", hvar[ir][it]->GetTitle(), status));
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()));
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);
297 double histmean = temphist->GetMean();
298 double histrms = temphist->GetRMS();
299 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
301 int fs = temphist->Fit(
"gaus",
"Q0");
303 B2INFO(Form(
"\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
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);
312 B2INFO(Form(
"\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
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));
328 if (cruns == 0)B2INFO(
"CDCDedxInjectTimeAlgorithm: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
333 int estart = erStart.first;
334 int rstart = erStart.second;
337 int rend = erEnd.second;
342 m_suffix = Form(
"e%dr%dto%d_nruns%d", estart, rstart, rend, cruns);
343 B2INFO(
"\t+ run = " << rend <<
", m_suffix = " <<
m_suffix <<
"");
345 m_suffix = Form(
"e%dr%d", estart, rstart);
346 B2INFO(
"tool run = " << estart <<
", exp = " << estart <<
", m_suffix = " <<
m_suffix <<
"");
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;
368 for (
unsigned int ib = 0; ib <
m_vtedges.size(); ib++)
376 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
378 for (
unsigned int it = 0; it <
m_tbins; 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);
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()));
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()));
411 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
412 for (
unsigned int it = 3; it <
m_tbins; it++) {
424 std::map<
int, std::vector<double>>& var,
425 std::map<
int, std::vector<double>>& time, TH1D*& htimes)
430 for (
int ir = 0; ir < 2; ir++) {
432 for (
int ix = varcorr[ir * 2].size(); ix -- > 0;) {
434 double var_thisbin = 1.0;
435 var_thisbin = var[ir * 2].at(ix);
437 double atime_thisbin = time[ir * 2].at(ix);
438 double ctime_thisbin = htimes->GetBinCenter(ix + 1);
440 if (atime_thisbin > 0 && atime_thisbin < 4e4 * 0.99) {
442 double var_nextbin = 1.0;
443 var_nextbin = var[ir * 2].at(ix + 1);
444 double var_diff = var_nextbin - var_thisbin;
446 double atime_nextbin = time[ir * 2].at(ix + 1);
447 double atime_diff = atime_nextbin - atime_thisbin;
449 double slope = (atime_diff > 0) ? var_diff / atime_diff : -1.0;
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));
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));
465 std::map<
int, std::vector<double>>& var,
466 std::map<
int, std::vector<double>>& varscal, std::string svar)
470 B2INFO(
"CDCDedxInjectTimeAlgorithm: normalising constants with plateau");
471 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
473 unsigned int msize = varscal[ir * 2].size();
476 for (
unsigned int im = 0; im < msize; im++) {
478 double mean = varscal[ir * 2].at(im);
479 if (time > 4e4 && mean > 0) {
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];
493 bool incomp_bin =
false;
494 std::vector<std::vector<double>> oldvectors;
496 int vsize = oldvectors.size();
497 if (vsize != 6) incomp_bin =
true;
499 for (
int iv = 0; iv < 2; iv++) {
500 if (oldvectors[iv * 3 + 1].size() != varscal[iv * 2].size()) incomp_bin =
true;
504 B2INFO(
"CDCDedxInjectTimeAlgorithm: started merging relative constants");
505 for (
int ir = 0; ir < 2; ir++) {
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) ;
515 }
else B2ERROR(
"CDCDedxInjectTimeAlgorithm: found incompatible bins for merging");
516 }
else B2INFO(
"CDCDedxInjectTimeAlgorithm: saving final (abs) calibration");
522 TCanvas cfit(
"cfit",
"cfit", 1000, 500);
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());
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++) {
532 hvar[ir][it]->SetFillColorAlpha(ir + 5, 0.25);
533 hvar[ir][it]->Draw();
535 cfit.Print(psname_fit.str().c_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());
546 TCanvas cestat(
"cestat",
"cestat", 1000, 500);
547 cestat.SetBatch(kTRUE);
551 auto hestats = getObjectPtr<TH1I>(
"hestats");
553 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
554 hestats->SetStats(0);
555 hestats->Draw(
"hist text");
558 auto htstats = getObjectPtr<TH1I>(
"htstats");
560 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
561 htstats->SetStats(0);
562 htstats->Draw(
"hist text");
564 cestat.Print(Form(
"%s_eventstat_%s.pdf",
m_prefix.data(),
m_suffix.data()));
570 TCanvas ctzoom(
"ctzoom",
"ctzoom", 1500, 450);
571 ctzoom.SetBatch(kTRUE);
573 for (
int wt = 0; wt < 3; wt++) {
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);
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");
587 ctzoom.Print(Form(
"%s_timezoom_%s.pdf]",
m_prefix.data(),
m_suffix.data()));
592 std::map<
int, std::vector<double>>& vresos, std::map<
int, std::vector<double>>& corr, std::string svar)
594 std::string sname[3] = {
"mean",
"reso"};
595 const int lcolors[
c_rings] = {2, 4};
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);
606 TLegend* mleg =
new TLegend(0.50, 0.54, 0.80, 0.75, NULL,
"brNDC");
607 mleg->SetBorderSize(0);
608 mleg->SetFillStyle(0);
610 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
612 std::string mtitle = Form(
"#mu(dedx), relative const. compare, (%s)",
m_suffix.data());
615 std::string rtitle = Form(
"#sigma(#chi), relative const. compare, (%s)",
m_suffix.data());
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()));
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()));
631 for (
unsigned int it = 0; it <
m_tbins; it++) {
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());
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());
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());
648 mleg->AddEntry(hmean[ir], Form(
"%s",
m_sring[ir].data()),
"lep");
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);
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");
661 if (ir == 1)mleg->Draw(
"same");
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");
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");
673 if (ir == 1)mleg->Draw(
"same");
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()));
681 for (
int ic = 0; ic < 2; ic++) {
692 const int lcolors[
c_rings] = {2, 4};
696 TCanvas* cconst =
new TCanvas(
"ctimestatconst",
"", 900, 500);
699 TLegend* rleg =
new TLegend(0.40, 0.60, 0.80, 0.72, NULL,
"brNDC");
700 rleg->SetBorderSize(0);
701 rleg->SetFillStyle(0);
703 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
705 std::string title = Form(
"injection time, her-ler comparison, (%s)",
m_suffix.data());
707 htimestat[ir]->SetTitle(Form(
"%s;injection time(#mu-second);norm. entries", title.data()));
709 for (
unsigned int it = 0; it <
m_tbins; 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());
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");
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()));
734 std::map<
int, std::vector<double>>& vresos,
735 std::array<double, numdedx::nrings>& scale, std::array<double, numdedx::nrings>& scale_reso)
738 std::vector<std::vector<double>> oldvectors;
741 const int c_type = 2;
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};
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);
755 TLegend* mleg =
new TLegend(0.50, 0.54, 0.80, 0.75, NULL,
"brNDC");
756 mleg->SetBorderSize(0);
757 mleg->SetFillStyle(0);
759 TLegend* rleg =
new TLegend(0.50, 0.54, 0.80, 0.75, NULL,
"brNDC");
760 rleg->SetBorderSize(0);
761 rleg->SetFillStyle(0);
763 for (
unsigned int ip = 0; ip < c_type; ip++) {
765 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
767 std::string hname = Form(
"hfmean_%s_%s_%s",
m_sring[ir].data(), stype[ip].data(),
m_suffix.data());
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()));
772 hname = Form(
"hfreso_%s_%s_%s",
m_sring[ir].data(), stype[ip].data(),
m_suffix.data());
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()));
777 for (
unsigned int it = 0; it <
m_tbins; it++) {
780 double mean = 0.0, reso = 0.0;
787 mean = oldvectors[ir * 3 + 1].at(it);
788 reso = oldvectors[ir * 3 + 2].at(it);
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());
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());
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");
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");
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()));
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];
836 std::map<
int, std::vector<double>>& vmeans)
839 unsigned int iv = ring * 2;
843 std::vector<unsigned int> tedges(sizet);
846 if (time >= 5e6) time = 5e6 - 10;
854 int thisbin = it, nextbin = it;
855 if (center != time && it > 0) {
860 if (it < sizet - 2)nextbin = it + 1;
861 else thisbin = it - 1;
865 double diff = vmeans[iv].at(2) - vmeans[iv].at(1) ;
868 if (it == 1) nextbin = it;
869 else nextbin = it + 1;
879 double thisdedx = vmeans[iv].at(thisbin);
880 double nextdedx = vmeans[iv].at(nextbin);
885 double newdedx = vmeans[iv].at(it);
886 if (thisbin != nextbin)
887 newdedx = thisdedx + ((nextdedx - thisdedx) / (nexttime - thistime)) * (time - thistime);
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
double m_chiMax
max range of chi
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)
double m_chiMin
min range of chi
bool m_ismakePlots
produce plots for monitoring
std::string m_suffix
string suffix for object names
int m_chiBins
bins for chi histogram
static const int c_rings
injection ring constants
int m_dedxBins
bins for dedx histogram
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
double m_dedxMax
max range of dedx
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.
double m_dedxMin
min range of dedx
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
Abstract base class for different kinds of events.