Belle II Software development
CDCDedx2DCellAlgorithm.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <reconstruction/calibration/CDCDedx2DCellAlgorithm.h>
10
11#include <reconstruction/dbobjects/CDCDedx2DCell.h>
12
13#include <TBox.h>
14#include <TCanvas.h>
15#include <TH2F.h>
16#include <TLine.h>
17#include <TMath.h>
18#include <TStyle.h>
19
20#include <algorithm>
21#include <fstream>
22
23using namespace Belle2;
24
25
26//-----------------------------------------------------------------
27// Implementation
28//-----------------------------------------------------------------
29
31 CalibrationAlgorithm("CDCDedxElectronCollector"),
32 fnEntaBinG(128),
33 fnDocaBinG(64),
34 fnEntaBinL(64),
35 fnDocaBinL(28),
36 feaLE(-TMath::Pi() / 2),
37 feaUE(+TMath::Pi() / 2),
38 fdocaLE(-1.50),
39 fdocaUE(1.50),
40 fSetPrefix("_it0"),
41 IsLocalBin(true),
42 IsMakePlots(false),
43 IsRS(true)
44{
45 // Set module properties
46 setDescription("A calibration algorithm for the CDC dE/dx two dimensional correction");
47}
48
49//-----------------------------------------------------------------
50// Run the calibration
51//-----------------------------------------------------------------
52
54{
55
56 //reading electron collector TREE
57 auto ttree = getObjectPtr<TTree>("tree");
58 if (ttree->GetEntries() < 100)return c_NotEnoughData;
59
60 std::vector<double>* dedxhit = 0, *doca = 0, *enta = 0;
61 std::vector<int>* layer = 0;
62
63 ttree->SetBranchAddress("dedxhit", &dedxhit);
64 ttree->SetBranchAddress("layer", &layer);
65 ttree->SetBranchAddress("docaRS", &doca);
66 ttree->SetBranchAddress("entaRS", &enta);
67
68 // Setting up bins for doca and entra angle
71
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);
75
76 if (!IsLocalBin) {
77 fEntaBinNums = globalbinsEnta;
79 } else {
81 fnEntaBinL = fEntaBinNums.at(fEntaBinNums.size() - 1) + 1;
82 }
83
84 fDocaBinNums = globalbinsDoca;
86
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;
90
91 for (int iea = 0; iea < fnEntaBinL; iea++) {
92 if (IsLocalBin) {
93 ifeaLE = fEntaBinValues.at(iea);
94 ifeaUE = fEntaBinValues.at(iea + 1);
95 } else {
96 ifeaLE = iea * feaBS - feaUE; //- because of -ive range shifting
97 ifeaUE = ifeaLE + feaBS;
98 }
99
100 for (int idoca = 0; idoca < fnDocaBinL; idoca++) {
101
102 ifdocaLE = idoca * fdocaBS - fdocaUE;
103 ifdocaUE = ifdocaLE + fdocaBS;
104
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");
110
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");
116
117 }
118 }
119
120 //Doca vs Enta stats
121 TH2D* hILDocaEntaG = new TH2D("hILDocaEntaG", "Doca vs EntA: Inner Layer", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG, feaLE, feaUE);
122 hILDocaEntaG->GetXaxis()->SetTitle("Doca");
123 hILDocaEntaG->GetYaxis()->SetTitle("Entrance angle (#theta)");
124
125 TH2D* hOLDocaEntaG = new TH2D("hOLDocaEntaG", "Doca vs EntA: Outer Layer", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG, feaLE, feaUE);
126 hOLDocaEntaG->GetXaxis()->SetTitle("Doca");
127 hOLDocaEntaG->GetYaxis()->SetTitle("Entrance angle (#theta)");
128
129 //when local enta angle is demanded
130 Double_t* RmapEntaValue = &fEntaBinValues[0];
131 TH2D* hILDocaEntaL = new TH2D("hILDocaEntaL", "Doca vs EntA: Inner Layer (rebin)", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
132 RmapEntaValue);
133 hILDocaEntaL->GetXaxis()->SetTitle("Doca");
134 hILDocaEntaL->GetYaxis()->SetTitle("Entrance angle (#theta)");
135
136 TH2D* hOLDocaEntaL = new TH2D("hOLDocaEntaL", "Doca vs EntA: Outer Layer (rebin)", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
137 RmapEntaValue);
138 hOLDocaEntaL->GetXaxis()->SetTitle("Doca");
139 hOLDocaEntaL->GetYaxis()->SetTitle("Entrance angle (#theta)");
140
141 TH1D* hILdEdx_all = new TH1D("hILdEdx_all", "", 500, 0., 5.);
142 TH1D* hOLdEdx_all = new TH1D("hOLdEdx_all", "", 500, 0., 5.);
143
144
145 Int_t ibinEA = 0, ibinDOCA = 0;
146 for (int i = 0; i < ttree->GetEntries(); ++i) {
147
148 ttree->GetEvent(i);
149
150 for (unsigned int j = 0; j < dedxhit->size(); ++j) {
151
152 if (dedxhit->at(j) == 0) continue;
153
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;
158
159 Double_t idocaHit = doca->at(j);
160 if (abs(idocaHit) > 1.50) continue;
161
162 //Bin corresponds to enta and doca value
163 ibinEA = (ieaHit - feaLE) / feaBS ; //from 0
164
165 ibinDOCA = (idocaHit - fdocaLE) / fdocaBS;
166
167 if (ibinEA >= fnEntaBinG || ibinDOCA >= fnDocaBinG) continue; //bin stats from 0
168
169 if (IsLocalBin) ibinEA = fEntaBinNums.at(ibinEA);
170
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));
176 } else {
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));
181 }
182 }
183 }
184
185
186 if (IsMakePlots) {
187 TCanvas* ctmpde = new TCanvas("DocavsEnta", "Doca vs Enta distributions", 800, 400);
188 if (IsLocalBin) {
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");
195 } else {
196 ctmpde->Divide(2, 1);
197 ctmpde->cd(1); hILDocaEntaG->Draw("colz");
198 ctmpde->cd(2); hOLDocaEntaG->Draw("colz");
199 }
200 ctmpde->SaveAs(Form("DocavsEnta_%s.pdf", fSetPrefix.data()));
201 ctmpde->SaveAs(Form("DocavsEnta_%s.root", fSetPrefix.data()));
202
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()));
208 }
209
210
211 //Calculating 5-75% global truncation mean
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();
216
217 Int_t lBinInLayer = 1, hBinInLayer = 1;
218 Int_t lBinOutLayer = 1, hBinOutLayer = 1;
219
220 for (int ibin = 1; ibin <= hILdEdx_all->GetNbinsX(); ibin++) {
221
222 if (InsumPer5 <= 0.05 * InLayInt) {
223 InsumPer5 += hILdEdx_all->GetBinContent(ibin);
224 lBinInLayer = ibin;
225 }
226
227 if (InsumPer75 <= 0.75 * InLayInt) {
228 InsumPer75 += hILdEdx_all->GetBinContent(ibin);
229 hBinInLayer = ibin;
230 }
231
232 if (OutsumPer5 <= 0.05 * OutLayInt) {
233 OutsumPer5 += hOLdEdx_all->GetBinContent(ibin);
234 lBinOutLayer = ibin;
235 }
236
237 if (OutsumPer75 <= 0.75 * OutLayInt) {
238 OutsumPer75 += hOLdEdx_all->GetBinContent(ibin);
239 hBinOutLayer = ibin;
240 }
241 }
242
243 short version = 1;
244 TCanvas* ctmp = new TCanvas("tmp", "tmp", 1200, 1200);
245 ctmp->Divide(4, 4);
246 std::stringstream psname; psname << Form("dedx_2dcell_%s.pdf[", fSetPrefix.data());
247 TLine* tl = new TLine();
248 tl->SetLineColor(kRed);
249
250 TBox* tb = new TBox();
251
252 if (IsMakePlots) {
253 ctmp->Print(psname.str().c_str());
254 psname.str(""); psname << Form("dedx_2dcell_%s.pdf", fSetPrefix.data());
255 }
256
257 gStyle->SetOptStat("nemriou");
258 TH1D* htemp = NULL;
259 std::vector<TH2F> twodcors;
260 TH2F tempTwoD = TH2F("tempTwoD", "dE/dx in bins of DOCA/Enta;DOCA;Entrance Angle", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
261 RmapEntaValue);
262 TH2F twodcor = TH2F("twodcorrection", "dE/dx in bins of DOCA/Enta;DOCA;Entrance Angle", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG,
263 feaLE, feaUE);
264
265 std::ofstream file2DILout;
266
267 for (int iIOLayer = 0; iIOLayer <= 1; iIOLayer++) {
268
269 file2DILout.open(Form("f2D_Constants_Local_Layer%d_%s.txt", iIOLayer, fSetPrefix.data()));
270
271 int startfrom = 1, endat = 1;
272 if (iIOLayer == 0) {
273 startfrom = lBinInLayer; endat = hBinInLayer;
274 } else {
275 startfrom = lBinOutLayer; endat = hBinOutLayer;
276 }
277
278 //std::cout << "Layer I/O # = " << iIOLayer << std::endl;
279 for (int iea = 1; iea <= fnEntaBinL; iea++) {
280
281 Int_t ieaprime = iea; //rotation symmtery for 1<->3 and 4<->2
282 if (IsRS)ieaprime = GetRotationSymmericBin(fnEntaBinL, iea);
283
284 for (int idoca = 1; idoca <= fnDocaBinL; idoca++) {
285
286 if (iIOLayer == 0)
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));
290 else continue;
291
292 double truncMean = 1.0;
293 if (htemp->GetEntries() < 400) truncMean = 1.0; //low stats
294 else {
295 double binweights = 0.0;
296 int sumofbc = 0;
297 for (int ibin = startfrom; ibin <= endat; ibin++) {
298 //std::cout << " dedxhit bin = " << ibin << ", Entries =" << htemp->GetBinContent(ibin) << std::endl;
299 if (htemp->GetBinContent(ibin) > 0) {
300 binweights += (htemp->GetBinContent(ibin) * htemp->GetBinCenter(ibin));
301 sumofbc += htemp->GetBinContent(ibin);
302 }
303 }
304 if (sumofbc > 0)truncMean = (double)(binweights / sumofbc);
305 else truncMean = 1.0;
306 }
307
308 if (truncMean <= 0)truncMean = 1.0; //protection only
309 tempTwoD.SetBinContent(idoca, iea, truncMean); //binning starts at 1
310
311 if (IsMakePlots) {
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);
318 htemp->DrawClone(); //clone is nessesory for pointer survival
319 tl->SetLineColor(kRed);
320 tl->SetX1(truncMean); tl->SetX2(truncMean);
321 tl->SetY1(0); tl->SetY2(htemp->GetMaximum());
322 tl->DrawClone("same");
323
324 tb->SetLineColor(kPink);
325 tb->SetLineStyle(6);
326 tb->SetFillColorAlpha(kPink, 0.15);
327 tb->SetX1(htemp->GetBinLowEdge(startfrom));
328 tb->SetY1(0);
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());
333 }
334
335 file2DILout << iea << ", \t" << idoca << ", \t" << truncMean << "\n";
336 htemp->Reset();
337 }
338 }
339
340 //conversion from local to global
341 for (int iea = 0; iea < fnEntaBinG; iea++) {
342 ibinEA = iea;
343 if (IsLocalBin)ibinEA = fEntaBinNums.at(iea);
344 for (int idoca = 0; idoca < fnDocaBinG; idoca++) {
345 ibinDOCA = idoca;
346 twodcor.SetBinContent(idoca + 1, iea + 1, tempTwoD.GetBinContent(ibinDOCA + 1, ibinEA + 1));
347 }
348 }
349
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);
353 twodcor.Reset();
354 tempTwoD.Reset();
355 file2DILout.close();
356 }
357
358 if (IsMakePlots) {
359 psname.str(""); psname << Form("dedx_2dcell_%s.pdf]", fSetPrefix.data());
360 ctmp->Print(psname.str().c_str());
361
362 // //Drawing final constants
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()));
369 }
370
371
372 B2INFO("dE/dx calibration done for 2D correction");
373 CDCDedx2DCell* gain = new CDCDedx2DCell(version, twodcors);
374 saveCalibration(gain, "CDCDedx2DCell");
375
376 delete htemp;
377 delete ctmp;
378 delete tl;
379 delete tb;
380 return c_OK;
381
382}
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
Definition: CDCDedx2DCell.h:26
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 successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
Abstract base class for different kinds of events.