9#include <reconstruction/calibration/CDCdEdx/CDCDedx2DCellAlgorithm.h>
11#include <reconstruction/dbobjects/CDCDedx2DCell.h>
36 feaLE(-TMath::Pi() / 2),
37 feaUE(+TMath::Pi() / 2),
46 setDescription(
"A calibration algorithm for the CDC dE/dx two dimensional correction");
57 auto ttree = getObjectPtr<TTree>(
"tree");
60 std::vector<double>* dedxhit = 0, *doca = 0, *enta = 0;
61 std::vector<int>* layer = 0;
63 ttree->SetBranchAddress(
"dedxhit", &dedxhit);
64 ttree->SetBranchAddress(
"layer", &layer);
65 ttree->SetBranchAddress(
"docaRS", &doca);
66 ttree->SetBranchAddress(
"entaRS", &enta);
72 std::vector<int> globalbinsEnta, globalbinsDoca;
73 for (
int ibin = 0; ibin <
fnEntaBinG; ibin++)globalbinsEnta.push_back(ibin);
74 for (
int ibin = 0; ibin <
fnDocaBinG; ibin++)globalbinsDoca.push_back(ibin);
87 std::vector<std::vector<TH1F*>> hILdEdxhitInEntaDocaBin(
fnEntaBinL, std::vector<TH1F*>(
fnDocaBinL, 0));
88 std::vector<std::vector<TH1F*>> hOLdEdxhitInEntaDocaBin(
fnEntaBinL, std::vector<TH1F*>(
fnDocaBinL, 0));
89 Double_t ifeaLE = 0, ifeaUE = 0, ifdocaLE = 0, ifdocaUE = 0;
97 ifeaUE = ifeaLE +
feaBS;
100 for (
int idoca = 0; idoca <
fnDocaBinL; idoca++) {
105 hILdEdxhitInEntaDocaBin[iea][idoca] =
new TH1F(Form(
"hILdEdxhitInEntaBin%dDocaBin%d", iea, idoca),
"bla-bla", 500, 0., 5.);
106 hILdEdxhitInEntaDocaBin[iea][idoca]->SetTitle(Form(
"IL: EntA = (%0.03f to %0.03f) and Doca = (%0.03f to %0.03f)", ifeaLE,
107 ifeaUE, ifdocaLE, ifdocaUE));
108 hILdEdxhitInEntaDocaBin[iea][idoca]->GetXaxis()->SetTitle(
"dedxhits in Inner Layer");
109 hILdEdxhitInEntaDocaBin[iea][idoca]->GetYaxis()->SetTitle(
"Entries");
111 hOLdEdxhitInEntaDocaBin[iea][idoca] =
new TH1F(Form(
"hOLdEdxhitInEntaBin%dDocaBin%d", iea, idoca),
"bla-bla", 500, 0., 5.);
112 hOLdEdxhitInEntaDocaBin[iea][idoca]->SetTitle(Form(
"OL: EntA = (%0.03f to %0.03f) and Doca = (%0.03f to %0.03f)", ifeaLE,
113 ifeaUE, ifdocaLE, ifdocaUE));
114 hOLdEdxhitInEntaDocaBin[iea][idoca]->GetXaxis()->SetTitle(
"dedxhits in Outer Layer");
115 hOLdEdxhitInEntaDocaBin[iea][idoca]->GetYaxis()->SetTitle(
"Entries");
122 hILDocaEntaG->GetXaxis()->SetTitle(
"Doca");
123 hILDocaEntaG->GetYaxis()->SetTitle(
"Entrance angle (#theta)");
126 hOLDocaEntaG->GetXaxis()->SetTitle(
"Doca");
127 hOLDocaEntaG->GetYaxis()->SetTitle(
"Entrance angle (#theta)");
133 hILDocaEntaL->GetXaxis()->SetTitle(
"Doca");
134 hILDocaEntaL->GetYaxis()->SetTitle(
"Entrance angle (#theta)");
138 hOLDocaEntaL->GetXaxis()->SetTitle(
"Doca");
139 hOLDocaEntaL->GetYaxis()->SetTitle(
"Entrance angle (#theta)");
141 TH1D* hILdEdx_all =
new TH1D(
"hILdEdx_all",
"", 500, 0., 5.);
142 TH1D* hOLdEdx_all =
new TH1D(
"hOLdEdx_all",
"", 500, 0., 5.);
145 Int_t ibinEA = 0, ibinDOCA = 0;
146 for (
int i = 0; i < ttree->GetEntries(); ++i) {
150 for (
unsigned int j = 0; j < dedxhit->size(); ++j) {
152 if (dedxhit->at(j) == 0)
continue;
154 Double_t ieaHit = enta->at(j);
155 if (ieaHit < -TMath::Pi() / 2.0) ieaHit += TMath::Pi() / 2.0;
156 else if (ieaHit > TMath::Pi() / 2.0) ieaHit -= TMath::Pi() / 2.0;
157 if (abs(ieaHit) > TMath::Pi() / 2.0)
continue;
159 Double_t idocaHit = doca->at(j);
160 if (abs(idocaHit) > 1.50)
continue;
171 if (layer->at(j) < 8) {
172 hILDocaEntaG->Fill(idocaHit, ieaHit);
173 if (
IsLocalBin)hILDocaEntaL->Fill(idocaHit, ieaHit);
174 hILdEdx_all->Fill(dedxhit->at(j));
175 hILdEdxhitInEntaDocaBin[ibinEA][ibinDOCA]->Fill(dedxhit->at(j));
177 hOLDocaEntaG->Fill(idocaHit, ieaHit);
178 if (
IsLocalBin)hOLDocaEntaL->Fill(idocaHit, ieaHit);
179 hOLdEdx_all->Fill(dedxhit->at(j));
180 hOLdEdxhitInEntaDocaBin[ibinEA][ibinDOCA]->Fill(dedxhit->at(j));
187 TCanvas* ctmpde =
new TCanvas(
"DocavsEnta",
"Doca vs Enta distributions", 800, 400);
189 ctmpde->SetCanvasSize(800, 800);
190 ctmpde->Divide(2, 2);
191 ctmpde->cd(1); hILDocaEntaG->Draw(
"colz");
192 ctmpde->cd(2); hOLDocaEntaG->Draw(
"colz");
193 ctmpde->cd(3); hILDocaEntaL->Draw(
"colz");
194 ctmpde->cd(4); hOLDocaEntaL->Draw(
"colz");
196 ctmpde->Divide(2, 1);
197 ctmpde->cd(1); hILDocaEntaG->Draw(
"colz");
198 ctmpde->cd(2); hOLDocaEntaG->Draw(
"colz");
200 ctmpde->SaveAs(Form(
"DocavsEnta_%s.pdf",
fSetPrefix.data()));
201 ctmpde->SaveAs(Form(
"DocavsEnta_%s.root",
fSetPrefix.data()));
203 TCanvas* ctem =
new TCanvas(
"Layerhisto",
"Inner and Outer Histo", 600, 600);
204 hOLdEdx_all->Draw(
"histo");
205 hILdEdx_all->SetMarkerColor(kRed);
206 hILdEdx_all->Draw(
"same histo");
207 ctem->SaveAs(Form(
"Layerhistodedxhit_TwoDCorr_%s.pdf",
fSetPrefix.data()));
212 double InsumPer5 = 0.0, InsumPer75 = 0.0;
213 double OutsumPer5 = 0.0, OutsumPer75 = 0.0;
214 double InLayInt = hILdEdx_all->Integral();
215 double OutLayInt = hOLdEdx_all->Integral();
217 Int_t lBinInLayer = 1, hBinInLayer = 1;
218 Int_t lBinOutLayer = 1, hBinOutLayer = 1;
220 for (
int ibin = 1; ibin <= hILdEdx_all->GetNbinsX(); ibin++) {
222 if (InsumPer5 <= 0.05 * InLayInt) {
223 InsumPer5 += hILdEdx_all->GetBinContent(ibin);
227 if (InsumPer75 <= 0.75 * InLayInt) {
228 InsumPer75 += hILdEdx_all->GetBinContent(ibin);
232 if (OutsumPer5 <= 0.05 * OutLayInt) {
233 OutsumPer5 += hOLdEdx_all->GetBinContent(ibin);
237 if (OutsumPer75 <= 0.75 * OutLayInt) {
238 OutsumPer75 += hOLdEdx_all->GetBinContent(ibin);
244 TCanvas* ctmp =
new TCanvas(
"tmp",
"tmp", 1200, 1200);
246 std::stringstream psname; psname << Form(
"dedx_2dcell_%s.pdf[",
fSetPrefix.data());
247 TLine* tl =
new TLine();
248 tl->SetLineColor(kRed);
250 TBox* tb =
new TBox();
253 ctmp->Print(psname.str().c_str());
254 psname.str(
""); psname << Form(
"dedx_2dcell_%s.pdf",
fSetPrefix.data());
257 gStyle->SetOptStat(
"nemriou");
259 std::vector<TH2F> twodcors;
265 std::ofstream file2DILout;
267 for (
int iIOLayer = 0; iIOLayer <= 1; iIOLayer++) {
269 file2DILout.open(Form(
"f2D_Constants_Local_Layer%d_%s.txt", iIOLayer,
fSetPrefix.data()));
271 int startfrom = 1, endat = 1;
273 startfrom = lBinInLayer; endat = hBinInLayer;
275 startfrom = lBinOutLayer; endat = hBinOutLayer;
281 Int_t ieaprime = iea;
284 for (
int idoca = 1; idoca <=
fnDocaBinL; idoca++) {
287 htemp = (TH1D*)hILdEdxhitInEntaDocaBin[ieaprime - 1][idoca - 1]->Clone(Form(
"hL%d_Ea%d_Doca%d", iIOLayer, iea, idoca));
288 else if (iIOLayer == 1)
289 htemp = (TH1D*)hOLdEdxhitInEntaDocaBin[ieaprime - 1][idoca - 1]->Clone(Form(
"hL%d_Ea%d_Doca%d", iIOLayer, iea, idoca));
292 double truncMean = 1.0;
293 if (htemp->GetEntries() < 400) truncMean = 1.0;
295 double binweights = 0.0;
297 for (
int ibin = startfrom; ibin <= endat; ibin++) {
299 if (htemp->GetBinContent(ibin) > 0) {
300 binweights += (htemp->GetBinContent(ibin) * htemp->GetBinCenter(ibin));
301 sumofbc += htemp->GetBinContent(ibin);
304 if (sumofbc > 0)truncMean = (double)(binweights / sumofbc);
305 else truncMean = 1.0;
308 if (truncMean <= 0)truncMean = 1.0;
309 tempTwoD.SetBinContent(idoca, iea, truncMean);
312 ctmp->cd(((idoca - 1) % 16) + 1);
313 htemp->SetTitle(Form(
"%s, mean = %0.4f", htemp->GetTitle(), truncMean));
314 htemp->SetFillColorAlpha(kGreen, 0.30);
315 if (truncMean >= 1.02 || truncMean <= 0.98)htemp->SetFillColor(kYellow);
316 if (truncMean >= 1.05 || truncMean <= 0.95)htemp->SetFillColor(kOrange);
317 if (truncMean >= 1.10 || truncMean <= 0.90)htemp->SetFillColor(kRed);
319 tl->SetLineColor(kRed);
320 tl->SetX1(truncMean); tl->SetX2(truncMean);
321 tl->SetY1(0); tl->SetY2(htemp->GetMaximum());
322 tl->DrawClone(
"same");
324 tb->SetLineColor(kPink);
326 tb->SetFillColorAlpha(kPink, 0.15);
327 tb->SetX1(htemp->GetBinLowEdge(startfrom));
329 tb->SetX2(htemp->GetBinLowEdge(endat));
330 tb->SetY2(htemp->GetMaximum());
331 tb->DrawClone(
"same");
332 if (idoca % 16 == 0)ctmp->Print(psname.str().c_str());
335 file2DILout << iea <<
", \t" << idoca <<
", \t" << truncMean <<
"\n";
344 for (
int idoca = 0; idoca <
fnDocaBinG; idoca++) {
346 twodcor.SetBinContent(idoca + 1, iea + 1, tempTwoD.GetBinContent(ibinDOCA + 1, ibinEA + 1));
350 if (iIOLayer == 0)twodcor.SetTitle(
"InnerLayer: dE/dx in bins of DOCA/Enta");
351 else if (iIOLayer == 1)twodcor.SetTitle(
"OuterLayer: dE/dx in bins of DOCA/Enta");
352 twodcors.push_back(twodcor);
359 psname.str(
""); psname << Form(
"dedx_2dcell_%s.pdf]",
fSetPrefix.data());
360 ctmp->Print(psname.str().c_str());
363 TCanvas* cconst =
new TCanvas(
"FinalConstantHistoMap",
"Inner and Outer Histo", 800, 400);
364 cconst->Divide(2, 1);
365 cconst->cd(1); twodcors.at(0).Draw(
"colz");
366 cconst->cd(2); twodcors.at(1).Draw(
"colz");
367 cconst->SaveAs(Form(
"Final2DConstantMap_wGlobalBin_%s.pdf",
fSetPrefix.data()));
368 cconst->SaveAs(Form(
"Final2DConstantMap_wGlobalBin_%s.root",
fSetPrefix.data()));
372 B2INFO(
"dE/dx calibration done for 2D correction");
bool IsRS
if rotation symmtery requested
bool IsLocalBin
if local variable bin requested
std::vector< double > fEntaBinValues
Vector for doca asym bin values.
Double_t fdocaUE
Upper edge of doca.
Double_t fdocaLE
Lower edge of doca.
int fnDocaBinG
doca angle bins, Global
Double_t feaBS
Binwidth edge of enta angle.
std::vector< int > fEntaBinNums
Vector for enta asym bin values.
int fnDocaBinL
doca angle bins, Local
CDCDedx2DCellAlgorithm()
Constructor: Sets the description the properties and the parameters of the algorithm.
int fnEntaBinL
etna angle bins, Local
std::string fSetPrefix
prefix to filename
Double_t fdocaBS
Binwidth edge of doca.
int GetRotationSymmericBin(int nbin, int ibin)
function to set rotation symmetry
bool IsMakePlots
produce plots for status
virtual EResult calibrate() override
2D Cell algorithm algorithm
int fnEntaBinG
Save arithmetic and truncated mean for the 'dedx' values.
std::vector< int > fDocaBinNums
Vector for enta asym bin #.
Double_t feaLE
Lower edge of enta angle.
void GetVariableBin(int nbin, std::vector< int > &nBinEnta0to100Per)
function to set variable bins
Double_t feaUE
Upper edge of enta angle.
dE/dx wire gain calibration constants
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
Abstract base class for different kinds of events.