9#include <cdc/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;
157 for (
int i = 0; i < ttree->GetEntries(); ++i) {
160 if (dedx <= 0 || injtime < 0 || injring < 0)
continue;
162 double corrcetion =
getCorrection(injring, injtime, vmeanscal);
163 double old_cor =
m_DBInjectTime->getCorrection(
"mean", injring, injtime);
165 dedx = (dedx * old_cor) / corrcetion;
171 if (injring > 0.5) wr = 1;
174 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
181 double chi = (dedx - predmean) / predsigma;
182 hdedx_corr[wr][tb]->Fill(dedx);
183 hchi[wr][tb]->Fill(chi);
187 std::map<int, std::vector<double>> vmeans_chi, vresos_chi;
191 std::map<int, std::vector<double>> vresoscorr;
195 std::map<int, std::vector<double>> vresoscal;
196 std::array<double, numdedx::nrings> scale_reso;
200 std::map<int, std::vector<double>> vmeans_corr, vresos_corr;
208 for (
int ir = 0; ir < 2; ir++) {
250 B2INFO(
"dE/dx Injection time calibration done");
260 B2INFO(
"Saving calibration for: " <<
m_suffix <<
"");
266 std::map<
int, std::vector<double>>& vmeans, std::map<
int, std::vector<double>>& vresos)
268 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
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) {
276 if (status != fitOK) {
277 mean = hvar[ir][it]->GetMean();
278 hvar[ir][it]->SetTitle(Form(
"%s, (%d)", hvar[ir][it]->GetTitle(), status));
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()));
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);
300 double histmean = temphist->GetMean();
301 double histrms = temphist->GetRMS();
302 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
304 int fs = temphist->Fit(
"gaus",
"Q0");
306 B2INFO(Form(
"\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
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);
315 B2INFO(Form(
"\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
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));
331 if (cruns == 0)B2INFO(
"CDCDedxInjectTimeAlgorithm: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
336 int estart = erStart.first;
337 int rstart = erStart.second;
340 int rend = erEnd.second;
345 m_suffix = Form(
"e%dr%dto%d_nruns%d", estart, rstart, rend, cruns);
346 B2INFO(
"\t+ run = " << rend <<
", m_suffix = " <<
m_suffix <<
"");
348 m_suffix = Form(
"e%dr%d", estart, rstart);
349 B2INFO(
"tool run = " << estart <<
", exp = " << estart <<
", m_suffix = " <<
m_suffix <<
"");
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;
371 for (
unsigned int ib = 0; ib <
m_vtedges.size(); ib++)
379 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
381 for (
unsigned int it = 0; it <
m_tbins; 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);
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()));
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()));
414 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
415 for (
unsigned int it = 3; it <
m_tbins; it++) {
427 std::map<
int, std::vector<double>>& var,
428 std::map<
int, std::vector<double>>& time, TH1D*& htimes)
433 for (
int ir = 0; ir < 2; ir++) {
435 for (
int ix = varcorr[ir * 2].size(); ix -- > 0;) {
437 double var_thisbin = 1.0;
438 var_thisbin = var[ir * 2].at(ix);
440 double atime_thisbin = time[ir * 2].at(ix);
441 double ctime_thisbin = htimes->GetBinCenter(ix + 1);
443 if (atime_thisbin > 0 && atime_thisbin < 4e4 * 0.99) {
445 double var_nextbin = 1.0;
446 var_nextbin = var[ir * 2].at(ix + 1);
447 double var_diff = var_nextbin - var_thisbin;
449 double atime_nextbin = time[ir * 2].at(ix + 1);
450 double atime_diff = atime_nextbin - atime_thisbin;
452 double slope = (atime_diff > 0) ? var_diff / atime_diff : -1.0;
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));
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));
468 std::map<
int, std::vector<double>>& var,
469 std::map<
int, std::vector<double>>& varscal, std::string svar)
473 B2INFO(
"CDCDedxInjectTimeAlgorithm: normalising constants with plateau");
474 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
476 unsigned int msize = varscal[ir * 2].size();
479 for (
unsigned int im = 0; im < msize; im++) {
481 double mean = varscal[ir * 2].at(im);
482 if (time > 4e4 && mean > 0) {
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];
496 bool incomp_bin =
false;
497 std::vector<std::vector<double>> oldvectors;
499 int vsize = oldvectors.size();
500 if (vsize != 6) incomp_bin =
true;
502 for (
int iv = 0; iv < 2; iv++) {
503 if (oldvectors[iv * 3 + 1].size() != varscal[iv * 2].size()) incomp_bin =
true;
507 B2INFO(
"CDCDedxInjectTimeAlgorithm: started merging relative constants");
508 for (
int ir = 0; ir < 2; ir++) {
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) ;
518 }
else B2ERROR(
"CDCDedxInjectTimeAlgorithm: found incompatible bins for merging");
519 }
else B2INFO(
"CDCDedxInjectTimeAlgorithm: saving final (abs) calibration");
525 TCanvas cfit(
"cfit",
"cfit", 1000, 500);
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());
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++) {
535 hvar[ir][it]->SetFillColorAlpha(ir + 5, 0.25);
536 hvar[ir][it]->Draw();
538 cfit.Print(psname_fit.str().c_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());
549 TCanvas cestat(
"cestat",
"cestat", 1000, 500);
550 cestat.SetBatch(kTRUE);
554 auto hestats = getObjectPtr<TH1I>(
"hestats");
556 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
557 hestats->SetStats(0);
558 hestats->Draw(
"hist text");
561 auto htstats = getObjectPtr<TH1I>(
"htstats");
563 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
564 htstats->SetStats(0);
565 htstats->Draw(
"hist text");
567 cestat.Print(Form(
"%s_eventstat_%s.pdf",
m_prefix.data(),
m_suffix.data()));
573 TCanvas ctzoom(
"ctzoom",
"ctzoom", 1500, 450);
574 ctzoom.SetBatch(kTRUE);
576 for (
int wt = 0; wt < 3; wt++) {
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);
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");
590 ctzoom.Print(Form(
"%s_timezoom_%s.pdf]",
m_prefix.data(),
m_suffix.data()));
595 std::map<
int, std::vector<double>>& vresos, std::map<
int, std::vector<double>>& corr, std::string svar)
597 std::string sname[3] = {
"mean",
"reso"};
598 const int lcolors[
c_rings] = {2, 4};
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);
609 TLegend* mleg =
new TLegend(0.50, 0.54, 0.80, 0.75, NULL,
"brNDC");
610 mleg->SetBorderSize(0);
611 mleg->SetFillStyle(0);
613 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
615 std::string mtitle = Form(
"#mu(dedx), relative const. compare, (%s)",
m_suffix.data());
618 std::string rtitle = Form(
"#sigma(#chi), relative const. compare, (%s)",
m_suffix.data());
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()));
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()));
634 for (
unsigned int it = 0; it <
m_tbins; it++) {
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());
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());
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());
651 mleg->AddEntry(hmean[ir], Form(
"%s",
m_sring[ir].data()),
"lep");
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);
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");
664 if (ir == 1)mleg->Draw(
"same");
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");
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");
676 if (ir == 1)mleg->Draw(
"same");
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()));
684 for (
int ic = 0; ic < 2; ic++) {
695 const int lcolors[
c_rings] = {2, 4};
699 TCanvas* cconst =
new TCanvas(
"ctimestatconst",
"", 900, 500);
702 TLegend* rleg =
new TLegend(0.40, 0.60, 0.80, 0.72, NULL,
"brNDC");
703 rleg->SetBorderSize(0);
704 rleg->SetFillStyle(0);
706 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
708 std::string title = Form(
"injection time, her-ler comparison, (%s)",
m_suffix.data());
710 htimestat[ir]->SetTitle(Form(
"%s;injection time(#mu-second);norm. entries", title.data()));
712 for (
unsigned int it = 0; it <
m_tbins; 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());
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");
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()));
737 std::map<
int, std::vector<double>>& vresos,
738 std::array<double, numdedx::nrings>& scale, std::array<double, numdedx::nrings>& scale_reso)
741 std::vector<std::vector<double>> oldvectors;
744 const int c_type = 2;
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};
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);
758 TLegend* mleg =
new TLegend(0.50, 0.54, 0.80, 0.75, NULL,
"brNDC");
759 mleg->SetBorderSize(0);
760 mleg->SetFillStyle(0);
762 TLegend* rleg =
new TLegend(0.50, 0.54, 0.80, 0.75, NULL,
"brNDC");
763 rleg->SetBorderSize(0);
764 rleg->SetFillStyle(0);
766 for (
unsigned int ip = 0; ip < c_type; ip++) {
768 for (
unsigned int ir = 0; ir <
c_rings; ir++) {
770 std::string hname = Form(
"hfmean_%s_%s_%s",
m_sring[ir].data(), stype[ip].data(),
m_suffix.data());
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()));
775 hname = Form(
"hfreso_%s_%s_%s",
m_sring[ir].data(), stype[ip].data(),
m_suffix.data());
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()));
780 for (
unsigned int it = 0; it <
m_tbins; it++) {
783 double mean = 0.0, reso = 0.0;
790 mean = oldvectors[ir * 3 + 1].at(it);
791 reso = oldvectors[ir * 3 + 2].at(it);
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());
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());
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");
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");
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()));
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];
839 std::map<
int, std::vector<double>>& vmeans)
842 unsigned int iv = ring * 2;
846 std::vector<unsigned int> tedges(sizet);
849 if (time >= 5e6) time = 5e6 - 10;
857 int thisbin = it, nextbin = it;
858 if (center != time && it > 0) {
863 if (it < sizet - 2)nextbin = it + 1;
864 else thisbin = it - 1;
868 double diff = vmeans[iv].at(2) - vmeans[iv].at(1) ;
871 if (it == 1) nextbin = it;
872 else nextbin = it + 1;
882 double thisdedx = vmeans[iv].at(thisbin);
883 double nextdedx = vmeans[iv].at(nextbin);
888 double newdedx = vmeans[iv].at(it);
889 if (thisbin != nextbin)
890 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)
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
Abstract base class for different kinds of events.