9#include <cdc/calibration/CDCdEdx/CDCDedxCosLayerAlgorithm.h>
39 setDescription(
"A calibration algorithm for CDC dE/dx electron cos(theta) dependence");
52 B2FATAL(
"There is no valid previous payload for CDCDedxCosineCor");
54 B2INFO(
"Preparing dE/dx calibration for CDC dE/dx electron saturation");
60 std::vector<double>* lDedx =
nullptr;
61 std::vector<int>* lLayer =
nullptr;
65 ttree->SetBranchAddress(
"ldedx", &lDedx);
66 ttree->SetBranchAddress(
"lLayer", &lLayer);
67 ttree->SetBranchAddress(
"costh", &costh);
68 ttree->SetBranchAddress(
"charge", &charge);
75 std::array<std::vector<double>,
m_kNGroups> corrFactor;
78 corrFactor[ig].assign(
m_cosBin, 1.0);
81 for (
int iter = 0; iter < 3; iter++) {
84 std::array<std::vector<TH1D*>,
m_kNGroups> hDedxCos_neg;
85 std::array<std::vector<TH1D*>,
m_kNGroups> hDedxCos_pos;
86 std::array<std::vector<TH1D*>,
m_kNGroups> hDedxCos_all;
88 defineHisto(hDedxCos_neg, Form(
"neg_iter%d", iter),
"e-");
89 defineHisto(hDedxCos_pos, Form(
"pos_iter%d", iter),
"e+");
90 defineHisto(hDedxCos_all, Form(
"all_iter%d", iter),
"e-,e+");
92 std::array<TH1D*, m_kNGroups> hDedxGroup{};
94 std::string title = Form(
"dedxhit dist (%s); dedxhit;entries",
m_label[il].data());
97 hDedxGroup[il]->SetTitle(title.c_str());
102 for (
int i = 0; i < ttree->GetEntries(); ++i) {
106 if (!lDedx || !lLayer)
continue;
107 if (lDedx->size() != lLayer->size())
continue;
108 if (charge == 0)
continue;
110 if (costh < TMath::Cos(150 * TMath::DegToRad()) ||
111 costh > TMath::Cos(17 * TMath::DegToRad()))
continue;
113 int bin = int((costh -
m_cosMin) / binW);
114 if (bin < 0 || bin >=
int(
m_cosBin))
continue;
116 for (
size_t j = 0; j < lDedx->size(); ++j) {
118 double val = lDedx->at(j);
119 int lay = lLayer->at(j);
121 if (val <= 0)
continue;
123 int ig = (lay < 8) ? 0 : ((lay < 14) ? 1 : 2);
125 val /= corrFactor[ig][bin];
129 hDedxCos_neg[ig][bin]->Fill(val);
131 hDedxCos_pos[ig][bin]->Fill(val);
134 hDedxCos_all[ig][bin]->Fill(val);
136 hDedxGroup[ig]->Fill(val);
142 if (charge < 0) hCosth_neg->Fill(costh);
143 else if (charge > 0) hCosth_pos->Fill(costh);
145 hCosth_all->Fill(costh);
149 std::array<std::vector<double>,
m_kNGroups> cosine;
153 int minGroup = 0, maxGroup = 0;
154 std::vector<double> vmean_neg, vmean_pos;
158 hDedxGroup[il]->SetTitle(
159 Form(
"%s;%d;%d", hDedxGroup[il]->GetTitle(), minGroup, maxGroup));
164 for (
unsigned int ibin = 0; ibin <
m_cosBin; ++ibin) {
170 double mean_neg =
extractCosMean(hDedxCos_neg[il][ibin], minGroup, maxGroup);
171 double mean_pos =
extractCosMean(hDedxCos_pos[il][ibin], minGroup, maxGroup);
173 bool has_neg = (hDedxCos_neg[il][ibin]->Integral() > 0);
174 bool has_pos = (hDedxCos_pos[il][ibin]->Integral() > 0);
176 if (has_neg && has_pos) mean = 0.5 * (mean_neg + mean_pos);
177 else if (has_neg) mean = mean_neg;
178 else if (has_pos) mean = mean_pos;
179 vmean_neg.push_back(mean_neg);
180 vmean_pos.push_back(mean_pos);
183 mean =
extractCosMean(hDedxCos_all[il][ibin], minGroup, maxGroup);
186 cosine[il].push_back(mean);
187 if (mean > 0) corrFactor[il][ibin] *= (mean / 1.25);
190 std::array<std::vector<double>, 3> cosMeanSets = {vmean_neg, vmean_pos, cosine[il]};
222 plotQaPars(hCosth_all, hCosth_pos, hCosth_neg);
243 if (cruns == 0) B2INFO(
"CDCDedxBadWires: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
248 int estart = erStart.first;
249 int rstart = erStart.second;
252 int eend = erEnd.first;
253 int rend = erEnd.second;
257 m_runExp = Form(
"Range (%d:%d,%d:%d)", estart, rstart, eend, rend);
259 else m_suffix = Form(
"e%d_r%dr%d", estart, rstart, rend);
267 hist->SetTitle(Form(
"cos(#theta) dist (%s); cos(#theta); Entries", chargeLabel.c_str()));
274 const std::string& chargeLabel)
282 for (
unsigned int i = 0; i <
m_cosBin; ++i) {
283 double coslow = i * binW +
m_cosMin;
284 double coshigh = coslow + binW;
286 hdedx[il].push_back(
new TH1D(Form(
"hDedxCos_%s_g%d_bin%d_%s", tag.c_str(), il, i,
m_suffix.data()),
289 hdedx[il][i]->SetTitle(Form(
"%s dE/dx dist (%s) in costh (%0.02f, %0.02f)",
290 m_label[il].c_str(), chargeLabel.c_str(), coslow, coshigh));
302 for (
unsigned int il = 0; il <
m_kNGroups; il++) {
306 B2ERROR(
"merging failed because of unmatch bins (old " <<
m_cosBin <<
" new " << nbins <<
")");
308 for (
unsigned int ibin = 0; ibin < nbins; ibin++) {
310 B2INFO(
"Cosine Corr for " <<
m_label[il] <<
" Bin # " << ibin <<
", Previous = " << prev <<
", Relative = " <<
m_coscors[il][ibin]
311 <<
", Merged = " << prev *
m_coscors[il][ibin]);
319 B2INFO(
"dE/dx calibration done for CDC dE/dx _eltron saturation");
321 std::vector<unsigned int> layerToGroup(56);
323 for (
unsigned int layer = 0; layer < 56; layer++) {
324 if (layer < 8) layerToGroup[layer] = 0;
325 else if (layer < 14) layerToGroup[layer] = 1;
326 else layerToGroup[layer] = 2;
337 TCanvas ctmp(
"tmp",
"tmp", 1200, 1200);
340 unsigned int nPads = nx * ny;
343 std::stringstream psname;
345 psname << Form(
"cdcdedx_coscorr_ldedx_%s_%s.pdf[", tag.c_str(),
m_suffix.data());
346 ctmp.Print(psname.str().c_str());
348 psname << Form(
"cdcdedx_coscorr_ldedx_%s_%s.pdf", tag.c_str(),
m_suffix.data());
352 for (
unsigned int ic = 0; ic <
m_cosBin; ic++) {
354 ctmp.cd(ic % nPads + 1);
355 hdedx[il][ic]->SetFillColor(4 + il);
357 hdedx[il][ic]->SetTitle(Form(
"%s;ldedx;entries", hdedx[il][ic]->GetTitle()));
358 hdedx[il][ic]->DrawClone(
"hist");
360 if (ic % nPads == nPads - 1 || ic ==
m_cosBin - 1) {
361 ctmp.Print(psname.str().c_str());
368 psname << Form(
"cdcdedx_coscorr_ldedx_%s_%s.pdf]", tag.c_str(),
m_suffix.data());
369 ctmp.Print(psname.str().c_str());
376 TCanvas cdedxlayer(Form(
"layerdedxhit_iter%d", iter),
"Inner and Outer Layer dedxhit dist", 2400, 800);
377 cdedxlayer.Divide(3, 1);
380 int minlay = 0, maxlay = 0;
382 minlay = std::stoi(hDedxGroup[il]->GetXaxis()->GetTitle());
383 maxlay = std::stoi(hDedxGroup[il]->GetYaxis()->GetTitle());
384 double lowedge = hDedxGroup[il]->GetXaxis()->GetBinLowEdge(minlay);
385 double upedge = hDedxGroup[il]->GetXaxis()->GetBinUpEdge(maxlay);
386 hDedxGroup[il]->SetTitle(Form(
"%s, trunc #rightarrow: %0.02f - %0.02f;dedxhit;entries", hDedxGroup[il]->GetTitle(), lowedge,
390 cdedxlayer.cd(il + 1);
391 hDedxGroup[il]->SetFillColor(kYellow);
392 hDedxGroup[il]->Draw(
"histo");
395 TH1D* hDedxGroupC = (TH1D*)hDedxGroup[il]->Clone(Form(
"hDedxGroupC%d", il));
396 hDedxGroupC->GetXaxis()->SetRange(minlay, maxlay);
397 hDedxGroupC->SetFillColor(kAzure + 1);
398 hDedxGroupC->Draw(
"same histo");
402 cdedxlayer.SaveAs(Form(
"cdcdedx_coscorr_dedxlay%s_iter%d.pdf",
m_suffix.data(), iter));
403 cdedxlayer.SaveAs(Form(
"cdcdedx_coscorr_dedxlay%s_iter%d.root",
m_suffix.data(), iter));
410 TCanvas ceadist(
"ceadist",
"Cosine distributions", 800, 600);
413 TLegend* leg =
new TLegend(0.6, 0.7, 0.8, 0.9);
417 hCosth_all->SetFillColor(kYellow);
418 hCosth_all->SetLineColor(kBlack);
419 hCosth_all->SetStats(0);
420 hCosth_all->Draw(
"hist");
421 leg->AddEntry(hCosth_all,
"all",
"f");
428 hCosth_pos->SetLineColor(kRed);
429 hCosth_pos->SetFillStyle(0);
430 hCosth_pos->SetStats(0);
431 hCosth_pos->Draw(
"hist same");
432 leg->AddEntry(hCosth_pos,
"pos",
"l");
436 hCosth_neg->SetLineColor(kBlue);
437 hCosth_neg->SetFillStyle(0);
438 hCosth_neg->SetStats(0);
439 hCosth_neg->Draw(
"hist same");
440 leg->AddEntry(hCosth_neg,
"neg",
"l");
446 ceadist.SaveAs(Form(
"cdcdedx_coscorr_cosine_%s.pdf",
m_suffix.data()));
447 ceadist.SaveAs(Form(
"cdcdedx_coscorr_cosine_%s.root",
m_suffix.data()));
453 TCanvas cconst(
"cconst",
"calibration Constants", 800, 600);
456 TLegend* leg =
new TLegend(0.6, 0.8, 0.9, 0.9);
457 leg->SetBorderSize(0);
458 leg->SetFillStyle(0);
460 std::vector<TH1D*> hists;
461 std::vector<int> colors = {kRed, kBlue, kGreen + 2};
466 TH1D* h =
new TH1D(Form(
"hconst_%d_%s", il,
m_suffix.data()),
"Relative constants; cos(#theta); constant",
m_cosBin,
m_cosMin,
470 for (
unsigned int jea = 0; jea <
m_cosBin; jea++) {
471 if (jea < cosine[il].size())
472 h->SetBinContent(jea + 1, cosine[il].at(jea));
475 double hmax = h->GetMaximum();
476 if (hmax > ymax) ymax = hmax;
484 hists[il]->SetLineColor(colors[il]);
485 hists[il]->SetStats(0);
486 hists[il]->SetLineWidth(2);
490 hists[il]->SetMaximum(ymax + 0.01);
491 hists[il]->Draw(
"hist");
493 hists[il]->Draw(
"hist same");
496 leg->AddEntry(hists[il],
m_label[il].data(),
"l");
501 cconst.SaveAs(Form(
"cdcdedx_coscorr_relconst_%s_iter%d.pdf",
m_suffix.data(), iter));
502 cconst.SaveAs(Form(
"cdcdedx_coscorr_relconst_%s_iter%d.root",
m_suffix.data(), iter));
505 for (
auto h : hists)
delete h;
523 for (
unsigned int iea = 0; iea < nbins; iea++) {
527 hold->SetBinContent(iea + 1, oldv);
528 hnew->SetBinContent(iea + 1, newv);
532 TH1D* hratio = (TH1D*)hnew->Clone(Form(
"hratio_%s",
m_label[il].data()));
533 hratio->Divide(hold);
535 TCanvas c(Form(
"c_%s",
m_label[il].data()), Form(
"Final constants %s",
m_label[il].data()), 1000, 500);
539 hnew->SetLineColor(kBlack);
541 hold->SetLineColor(kRed);
544 double min = std::min(hnew->GetMinimum(), hold->GetMinimum());
545 double max = std::max(hnew->GetMaximum(), hold->GetMaximum());
546 hnew->GetYaxis()->SetRangeUser(min * 0.95, max * 1.05);
549 hold->Draw(
"hist same");
551 auto leg =
new TLegend(0.6, 0.75, 0.85, 0.88);
552 leg->SetBorderSize(0);
553 leg->SetFillStyle(0);
554 leg->AddEntry(hnew,
"New",
"l");
555 leg->AddEntry(hold,
"Old",
"l");
560 hratio->SetLineColor(kBlue);
562 hratio->SetTitle(Form(
"Ratio: new/old, %s;cos(#theta); New / Old",
m_label[il].data()));
563 hratio->GetYaxis()->SetRangeUser(0.2, 1.2);
564 hratio->Draw(
"hist");
567 line->SetLineStyle(2);
570 c.SaveAs(Form(
"cdcdedx_coscorr_fconsts_%s_%s.pdf",
m_label[il].data(),
m_suffix.data()));
571 c.SaveAs(Form(
"cdcdedx_coscorr_fconsts_%s_%s.root",
m_label[il].data(),
m_suffix.data()));
585 TCanvas cconst(Form(
"cconst_%s_iter%d", sltag.c_str(), iter),
"calibration Constants", 800, 600);
588 TLegend* leg =
new TLegend(0.6, 0.8, 0.9, 0.9);
589 leg->SetBorderSize(0);
590 leg->SetFillStyle(0);
592 std::vector<TH1D*> hists;
594 std::array<std::string, 3> labels = {
"e^{+}",
"e^{-}",
"Average"};
595 std::vector<int> colors = {kRed, kBlue, kBlack + 2};
597 for (
int il = 0; il < 3; il++) {
599 TH1D* h =
new TH1D(Form(
"hconst_%d_%s_%s_iter%d", il, sltag.c_str(),
m_suffix.data(), iter),
600 Form(
"Relative mean %s, iter %d; cos(#theta); constant", sltag.c_str(), iter),
603 for (
unsigned int jea = 0; jea <
m_cosBin; jea++) {
604 if (jea < mean[il].size())
605 h->SetBinContent(jea + 1, mean[il].at(jea));
608 h->SetLineColor(colors[il]);
612 if (il == 0) h->Draw(
"hist");
613 else h->Draw(
"hist same");
615 leg->AddEntry(h, labels[il].c_str(),
"l");
621 cconst.SaveAs(Form(
"cdcdedx_coscorr_relmean_iter%d_%s_%s.pdf", iter, sltag.c_str(),
m_suffix.data()));
622 cconst.SaveAs(Form(
"cdcdedx_coscorr_relmean_iter%d_%s_%s.root", iter, sltag.c_str(),
m_suffix.data()));
624 for (
auto h : hists)
delete h;
632 TCanvas cstats(
"cstats",
"cstats", 1000, 500);
633 cstats.SetBatch(kTRUE);
639 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
640 hestats->SetStats(0);
641 hestats->DrawCopy(
"");
647 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
648 htstats->SetStats(0);
649 htstats->DrawCopy(
"");
651 cstats.Print(Form(
"cdcdedx_coscorr_stats_%s.pdf",
m_suffix.data()));
659 double sum = hist->Integral();
660 if (sum <= 0 || hist->GetNbinsX() <= 0) {
661 binlow = 1; binhigh = 1;
665 binlow = 1.0; binhigh = 1.0;
666 double sumPer5 = 0.0, sumPer75 = 0.0;
667 for (
int ibin = 1; ibin <= hist->GetNbinsX(); ibin++) {
668 double bcdedx = hist->GetBinContent(ibin);
686 if (hist->Integral() < 100)
return 1.0;
688 if (binlow <= 0 || binhigh > hist->GetNbinsX())
return 1.0;
690 double binweights = 0., sumofbc = 0.;
691 for (
int ibin = binlow; ibin <= binhigh; ibin++) {
692 double bcdedx = hist->GetBinContent(ibin);
694 binweights += (bcdedx * hist->GetBinCenter(ibin));
698 if (sumofbc > 0)
return binweights / sumofbc;
704 if (!hist || hist->Integral() <= 0)
return 1.0;
709 hist->SetTitle(Form(
"%s, mean = %0.5f", hist->GetTitle(), hist->GetMean()));
710 return hist->GetMean();
713 int minbin = 1, maxbin = 1;
723 hist->SetTitle(Form(
"%s, mean = %0.5f;%d;%d", hist->GetTitle(), mean, minbin, maxbin));
void plotdedxHist(std::array< std::vector< TH1D * >, m_kNGroups > &hdedx, const std::string &tag)
function to draw the dE/dx histogram in costh bins
void plotRelConst(const std::array< std::vector< double >, m_kNGroups > &cosine, int iter)
Plot relative calibration constants vs costh for all SL groups (overlayed)
void plotQaPars(TH1D *hCosth_all, TH1D *hCosth_pos, TH1D *hCosth_neg)
function to costh distribution for Inner/Outer layer
double m_truncMax
upper threshold on truncation
void plotLayerDist(std::array< TH1D *, m_kNGroups > &hdedxlay, int iter)
function to draw dedx dist.
double m_truncMin
lower threshold on truncation
void getTruncatedBins(TH1D *hist, int &binlow, int &binhigh)
function to get bins of truncation from histogram
void defineHisto(std::array< std::vector< TH1D * >, m_kNGroups > &hdedx, const std::string &tag, const std::string &chargeLabel)
function to define dE/dx histograms
void getExpRunInfo()
function to get extract calibration run/exp
static constexpr int m_kNGroups
SL grouping: inner (SL0), middle (SL1), outer (SL2–8)
int m_dedxBin
number of bins for dedx histogram
double m_cosMax
max cosine angle for cal
void plotmeanChargeOverlay(const std::array< std::vector< double >, 3 > &cosine_pos, const std::string &sltag, int iter)
Plot overlay of positive, negative, and average cosine means for one SL group.
unsigned int getRepresentativeLayer(unsigned int igroup) const
Representative CDC layer for each SL group (used to access group-wise constants): SL0 => 1,...
CDCDedxCosLayerAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
TH1D * defineCosthHist(const std::string &tag, const std::string &chargeLabel)
function to define cosine histograms
std::string m_suffix
add suffix to all plot name
double extractCosMean(TH1D *&hist, int fixedLow=1, int fixedHigh=1)
Extract mean dE/dx vs costh for a given group from the histogram.
double getTruncationMean(TH1D *hist, int binlow, int binhigh)
function to get truncated mean
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
void plotConstants()
function to draw the old/new final constants
bool isFixTrunc
true = fix window for all out/inner layers
std::string m_runExp
add run and exp to title of plot
const std::array< std::string, m_kNGroups > m_label
add inner/outer superlayer label
double m_cosMin
min cosine angle for cal
void plotEventStats()
function to draw the stats plots
virtual EResult calibrate() override
Cosine algorithm.
bool isMethodSep
if e+e- need to be consider sep
double m_dedxMax
max dedx range for gain cal
void createPayload()
function to store new payload after full calibration
bool isMakePlots
produce plots for status
bool isMerge
merge payload at the of calibration
bool isUseTrunc
true if truncated mean for SL0,1
double m_dedxMin
min dedx range for gain cal
unsigned int m_cosBin
number of bins across cosine range
std::vector< std::vector< double > > m_coscors
final vectors of calibration
dE/dx cosine gain calibration constants
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
static 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.
CalibrationAlgorithm(const std::string &collectorModuleName)
Constructor - sets the prefix for collected objects (won't be accesses until execute(....
std::shared_ptr< T > getObjectPtr(const std::string &name, const std::vector< Calibration::ExpRun > &requestedRuns)
Get calibration data object by name and list of runs, the Merge function will be called to generate t...
Abstract base class for different kinds of events.