9#include <cdc/calibration/CDCdEdx/CDCDedxCosineAlgorithm.h>
38 setDescription(
"A calibration algorithm for CDC dE/dx electron cos(theta) dependence");
51 B2FATAL(
"There is no valid previous payload for CDCDedxCosineCor");
53 B2INFO(
"Preparing dE/dx calibration for CDC dE/dx electron saturation");
58 B2ERROR(
"Input tree 'tree' not found");
63 double dedx, costh;
int charge;
64 ttree->SetBranchAddress(
"dedx", &dedx);
65 ttree->SetBranchAddress(
"costh", &costh);
66 ttree->SetBranchAddress(
"charge", &charge);
71 std::vector<TH1D*> hDedxCos_neg, hDedxCos_pos, hDedxCos_all;
84 for (
int i = 0; i < ttree->GetEntries(); ++i) {
89 if (dedx <= 0 || charge == 0)
continue;
92 if (costh < TMath::Cos(150 * TMath::DegToRad()) || costh > TMath::Cos(17 * TMath::DegToRad()))
continue;
94 int bin = int((costh -
m_cosMin) / binW);
95 if (bin < 0 || bin >=
static_cast<int>(
m_cosBin))
continue;
99 hCosth_neg->Fill(costh);
100 hDedxCos_neg[bin]->Fill(dedx);
101 }
else if (charge > 0) {
102 hCosth_pos->Fill(costh);
103 hDedxCos_pos[bin]->Fill(dedx);
106 hCosth_all->Fill(costh);
107 hDedxCos_all[bin]->Fill(dedx);
112 std::vector<double> cosine;
114 std::vector<std::vector<double>> dedxAll(4);
115 std::vector<std::vector<double>> dedxNeg(4);
116 std::vector<std::vector<double>> dedxPos(4);
118 for (
unsigned int i = 0; i <
m_cosBin; ++i) {
120 double meanDedx = 1.0;
121 double meanDedxErr = 0.0;
126 meanDedx = fitAll.
mean;
128 dedxAll[0].push_back(fitAll.
mean);
129 dedxAll[1].push_back(fitAll.
meanErr);
130 dedxAll[2].push_back(fitAll.
sigma);
131 dedxAll[3].push_back(fitAll.
sigmaErr);
141 if (fitPos.
status !=
"FitOK" && fitNeg.
status ==
"FitOK") {
143 hDedxCos_pos[i]->SetTitle(Form(
"%s, mean (manual) = elec left", hDedxCos_pos[i]->GetTitle()));
144 }
else if (fitNeg.
status !=
"FitOK" && fitPos.
status ==
"FitOK") {
146 hDedxCos_neg[i]->SetTitle(Form(
"%s, mean (manual) = posi right", hDedxCos_neg[i]->GetTitle()));
147 }
else if (fitNeg.
status !=
"FitOK" && fitPos.
status !=
"FitOK") {
152 dedxNeg[0].push_back(fitNeg.
mean);
153 dedxNeg[1].push_back(fitNeg.
meanErr);
154 dedxNeg[2].push_back(fitNeg.
sigma);
155 dedxNeg[3].push_back(fitNeg.
sigmaErr);
157 dedxPos[0].push_back(fitPos.
mean);
158 dedxPos[1].push_back(fitPos.
meanErr);
159 dedxPos[2].push_back(fitPos.
sigma);
160 dedxPos[3].push_back(fitPos.
sigmaErr);
162 meanDedx = 0.5 * (fitNeg.
mean + fitPos.
mean);
163 if (meanDedx <= 0.0) meanDedx = 1.0;
168 dedxAll[0].push_back(meanDedx);
169 dedxAll[1].push_back(meanDedxErr);
173 cosine.push_back(meanDedx);
208 if (cruns == 0) B2INFO(
"CDCDedxCosineCor: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
213 int estart = erStart.first;
214 int rstart = erStart.second;
217 int rend = erEnd.second;
222 else m_suffix = Form(
"e%d_r%dr%d", estart, rstart, rend);
227 const std::string& chargeLabel)
234 for (
unsigned int i = 0; i <
m_cosBin; ++i) {
235 double coslow = i * binW +
m_cosMin;
236 double coshigh = coslow + binW;
240 hdedx[i]->SetTitle(Form(
"dE/dx dist (%s) in costh (%0.02f, %0.02f);dE/dx (no had sat, for %s);Entries", chargeLabel.c_str(), coslow,
241 coshigh, chargeLabel.c_str()));
251 hist->SetTitle(
"cos(#theta) dist (e- and e+); cos(#theta); Entries");
263 hist->SetFillColorAlpha(kAzure + 1, 0.30);
265 if (fitValues.
status ==
"FitOK") {
266 TF1* fitFunc = hist->GetFunction(
"gaus");
268 fitValues.
mean = fitFunc->GetParameter(1);
269 fitValues.
meanErr = fitFunc->GetParError(1);
270 fitValues.
sigma = fitFunc->GetParameter(2);
271 fitValues.
sigmaErr = fitFunc->GetParError(2);
273 std::string fitSummary = Form(
"#mu_{fit}: %0.03f #pm %0.03f, #sigma_{fit}: %0.03f",
276 hist->SetTitle(Form(
"%s, %s", hist->GetTitle(), fitSummary.data()));
286 if (temphist->Integral() < 2000) {
287 B2INFO(Form(
"\tThis hist (%s) have insufficient entries to perform fit (%0.03f)", temphist->GetName(), temphist->Integral()));
291 temphist->GetXaxis()->SetRange(temphist->FindFirstBinAbove(0, 1), temphist->FindLastBinAbove(0, 1));
292 int fs = temphist->Fit(
"gaus",
"QR");
294 B2INFO(Form(
"\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
295 status =
"FitFailed";
298 double meanDedx = temphist->GetFunction(
"gaus")->GetParameter(1);
299 double width = temphist->GetFunction(
"gaus")->GetParameter(2);
300 temphist->GetXaxis()->SetRangeUser(meanDedx - 5.0 * width, meanDedx + 5.0 * width);
301 fs = temphist->Fit(
"gaus",
"QR",
"", meanDedx -
m_sigLim * width, meanDedx +
m_sigLim * width);
303 B2INFO(Form(
"\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
304 status =
"FitFailed";
307 temphist->GetXaxis()->SetRangeUser(meanDedx - 5.0 * width, meanDedx + 5.0 * width);
308 B2INFO(Form(
"\tFit for hist (%s) sucessfull (status = %d)", temphist->GetName(), fs));
320 for (
unsigned int il = 0; il <
m_kNGroups; il++) {
325 B2ERROR(
"merging failed because of unmatch bins (old "
326 << nbins <<
" new " <<
m_cosBin <<
")");
330 for (
unsigned int ibin = 0; ibin < nbins; ibin++) {
332 double value = cosine[ibin];
339 B2INFO(
"Cosine Corr for " <<
m_label[il]
341 <<
", Previous = " << prev
342 <<
", Relative = " << cosine[ibin]
343 <<
", Merged = " << value);
351 B2INFO(
"dE/dx calibration done for CDC dE/dx electron saturation");
353 std::vector<unsigned int> layerToGroup(56);
355 for (
unsigned int layer = 0; layer < 56; layer++) {
356 if (layer < 8) layerToGroup[layer] = 0;
357 else if (layer < 14) layerToGroup[layer] = 1;
358 else layerToGroup[layer] = 2;
367 std::vector<TH1D*>& hDedxCos_neg,
368 std::vector<TH1D*>& hDedxCos_pos)
371 TCanvas ctmp(
"tmp",
"tmp", 1200, 1200);
374 unsigned int nPads = nx * ny;
377 std::stringstream psname;
379 psname << Form(
"cdcdedx_coscorr_dedx_%s.pdf[",
m_suffix.data());
380 ctmp.Print(psname.str().c_str());
382 psname << Form(
"cdcdedx_coscorr_dedx_%s.pdf",
m_suffix.data());
384 for (
unsigned int ic = 0; ic <
m_cosBin; ic++) {
386 ctmp.cd(ic % nPads + 1);
387 hDedxCos_all[ic]->SetStats(0);
388 hDedxCos_all[ic]->SetFillColorAlpha(kYellow, 0.25);
389 hDedxCos_all[ic]->DrawCopy();
391 if (ic % nPads == nPads - 1 || ic ==
m_cosBin - 1) {
392 ctmp.Print(psname.str().c_str());
400 hDedxCos_neg[ic]->SetFillColorAlpha(kRed, 0.25);
401 hDedxCos_neg[ic]->DrawCopy();
406 hDedxCos_pos[ic]->SetFillColorAlpha(kBlue, 0.25);
407 hDedxCos_pos[ic]->DrawCopy();
409 ctmp.Print(psname.str().c_str());
417 psname << Form(
"cdcdedx_coscorr_dedx_%s.pdf]",
m_suffix.data());
418 ctmp.Print(psname.str().c_str());
425 TCanvas ceadist(
"ceadist",
"Cosine distributions", 800, 600);
431 TLegend* leg =
new TLegend(0.6, 0.7, 0.8, 0.9);
434 hCosth_neg->SetLineColor(kRed);
435 hCosth_neg->SetFillColorAlpha(kYellow, 0.55);
436 hCosth_neg->SetStats(0);
437 hCosth_neg->Draw(
"hist");
438 leg->AddEntry(hCosth_neg,
"neg",
"f");
441 hCosth_pos->SetLineColor(kBlue);
442 hCosth_pos->SetFillColorAlpha(kGray, 0.35);
443 hCosth_pos->SetStats(0);
444 hCosth_pos->Draw(
"hist same");
445 leg->AddEntry(hCosth_pos,
"pos",
"f");
452 hCosth_all->SetFillColorAlpha(kGray, 0.25);
453 hCosth_all->SetLineColor(kGray);
454 hCosth_all->SetStats(0);
455 hCosth_all->Draw(
"hist");
459 ceadist.SaveAs(Form(
"cdcdedx_coscorr_cosine_%s.pdf",
m_suffix.data()));
460 ceadist.SaveAs(Form(
"cdcdedx_coscorr_cosine_%s.root",
m_suffix.data()));
466 const std::vector<std::vector<double>>& dedxNeg,
467 const std::vector<std::vector<double>>& dedxPos)
485 for (
unsigned int i = 0; i <
m_cosBin; i++) {
486 hMean_all->SetBinContent(i + 1, dedxAll[0][i]);
487 hMean_all->SetBinError(i + 1, dedxAll[1][i]);
490 hMean_el->SetBinContent(i + 1, dedxNeg[0][i]);
491 hMean_el->SetBinError(i + 1, dedxNeg[1][i]);
492 hSig_el->SetBinContent(i + 1, dedxNeg[2][i]);
493 hSig_el->SetBinError(i + 1, dedxNeg[3][i]);
495 hMean_po->SetBinContent(i + 1, dedxPos[0][i]);
496 hMean_po->SetBinError(i + 1, dedxPos[1][i]);
497 hSig_po->SetBinContent(i + 1, dedxPos[2][i]);
498 hSig_po->SetBinError(i + 1, dedxPos[3][i]);
500 hSig_all->SetBinContent(i + 1, dedxAll[2][i]);
501 hSig_all->SetBinError(i + 1, dedxAll[3][i]);
506 TCanvas* ctmp =
new TCanvas(
"c_fit",
"Mean & Sigma", 1000, 500);
511 setHist(hMean_all, kBlack,
"dedx rel(#mu_{fit}) for e- and e+ combined", 0.97, 1.04);
515 setHist(hMean_el, kRed,
"comparison of dedx #mu_{fit}^{rel}", 0.96, 1.04);
516 setHist(hMean_po, kBlue,
"", 0.96, 1.04);
518 hMean_el->Draw(
"E1");
519 hMean_po->Draw(
"E1 same");
520 hMean_all->Draw(
"E1 same");
523 hMean_all->Draw(
"E1");
529 setHist(hSig_all, kBlack,
"dedx rel(#sigma_{fit}) for e- and e+ combined", 0.05, 0.23);
533 setHist(hSig_el, kRed,
"comparison of dedx #sigma_{fit}^{rel}", 0.05, 0.23, 24);
534 setHist(hSig_po, kBlue,
"", 0.05, 0.23, 25);
537 hSig_po->Draw(
"E1 same");
540 hSig_all->Draw(
"E1");
542 B2INFO(
"Plotting finished ");
545 ctmp->SaveAs(Form(
"cdcdedx_coscorr_fit_%s.pdf",
m_suffix.data()));
564 TH1D* hnew =
new TH1D(Form(
"hnew_%s",
m_label[il].data()), Form(
"Final const: %s;cos(#theta);dedx #mu_{fit}",
m_label[il].data()),
568 TH1D* hold =
new TH1D(Form(
"hold_%s",
m_label[il].data()), Form(
"Final const: %s;cos(#theta);dedx #mu_{fit}",
m_label[il].data()),
572 for (
unsigned int iea = 0; iea < nbins; iea++) {
576 hold->SetBinContent(iea + 1, oldv);
577 hnew->SetBinContent(iea + 1, newv);
581 TH1D* hratio = (TH1D*)hnew->Clone(Form(
"hratio_%s",
m_label[il].data()));
582 hratio->Divide(hold);
584 TCanvas c(Form(
"c_%s",
m_label[il].data()), Form(
"Final constants %s",
m_label[il].data()), 1000, 500);
590 hnew->SetLineColor(kBlack);
592 hold->SetLineColor(kRed);
595 double min = std::min(hnew->GetMinimum(), hold->GetMinimum());
596 double max = std::max(hnew->GetMaximum(), hold->GetMaximum());
597 hnew->GetYaxis()->SetRangeUser(min * 0.95, max * 1.05);
600 hold->Draw(
"hist same");
602 auto leg =
new TLegend(0.6, 0.75, 0.85, 0.88);
603 leg->SetBorderSize(0);
604 leg->SetFillStyle(0);
605 leg->AddEntry(hnew,
"New",
"l");
606 leg->AddEntry(hold,
"Old",
"l");
611 hratio->SetLineColor(kBlue);
613 hratio->SetTitle(Form(
"Ratio: new/old, %s;cos(#theta); New / Old",
m_label[il].data()));
614 hratio->GetYaxis()->SetRangeUser(0.8, 1.2);
615 hratio->Draw(
"hist");
618 line->SetLineStyle(2);
621 c.SaveAs(Form(
"cdcdedx_coscorr_fconsts_%s_%s.pdf",
m_label[il].data(),
m_suffix.data()));
622 c.SaveAs(Form(
"cdcdedx_coscorr_fconsts_%s_%s.root",
m_label[il].data(),
m_suffix.data()));
636 TCanvas cstats(
"cstats",
"cstats", 1000, 500);
637 cstats.SetBatch(kTRUE);
643 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
644 hestats->SetStats(0);
645 hestats->DrawCopy(
"");
651 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
652 htstats->SetStats(0);
653 htstats->DrawCopy(
"");
655 cstats.Print(Form(
"cdcdedx_coscorr_stats_%s.pdf",
m_suffix.data()));
FitValues fitHistogram(TH1D *&hist)
Fit histogram with Gaussian and return mean, error, and width.
void getExpRunInfo()
function to extract calibration run/exp
TH1D * defineCosthHist(const std::string &tag)
function to define cosine histograms
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 plotCosThetaDist(TH1D *hCosth_all, TH1D *hCosth_pos, TH1D *hCosth_neg)
Plot cos(theta) distributions for all, positive, and negative tracks.
bool isMergePayload
merge payload at the time of calibration
unsigned int getRepresentativeLayer(unsigned int igroup) const
Representative CDC layer for each SL group (used to access group-wise constants): SL0 => 1,...
std::string m_suffix
add suffix to all plot name
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
void plotConstants()
function to draw the old/new final constants
CDCDedxCosineAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
void fitGaussianWithRange(TH1D *&temphist, TString &status)
function to fit histogram in each cosine bin
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.
double m_sigLim
gaussian fit sigma limit
void plotFitResults(const std::vector< std::vector< double > > &dedxAll, const std::vector< std::vector< double > > &dedxNeg, const std::vector< std::vector< double > > &dedxPos)
Plot dE/dx fit results for all, negative, and positive tracks.
void plotdedxHist(std::vector< TH1D * > &hDedxCos_all, std::vector< TH1D * > &hDedxCos_neg, std::vector< TH1D * > &hDedxCos_pos)
function to draw the dE/dx histogram in costh bins
bool isMethodSep
if e+ e- need to be consider sep
void createPayload(std::vector< double > cosine)
function to store new payload after full calibration
double m_dedxMax
max dedx range for gain cal
bool isMakePlots
produce plots for status
double m_dedxMin
min dedx range for gain cal
void setHist(TH1D *h, int color, const char *title, double ymin, double ymax, int marker=20)
Set basic style (color, marker, title, y-range) for a TH1D histogram.
unsigned int m_cosBin
number of bins across cosine range
void defineHisto(std::vector< TH1D * > &hdedx, const std::string &tag, const std::string &chargeLabel)
function to define dE/dx histograms
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.
@ c_Failure
Failed =3 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.
Container for Gaussian fit results of a histogram.
double sigma
sigma : fitted width
double meanErr
meanErr : uncertainty on the mean
double sigmaErr
sigmaErr : uncertainty on the width
double mean
mean : fitted mean value
TString status
status : fit status flag (e.g.