Belle II Software  release-06-02-00
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 
23 using 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
69  feaBS = (feaUE - feaLE) / fnEntaBinG;
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 = 0;
164  ibinEA = (ieaHit - feaLE) / feaBS ; //from 0
165 
166  ibinDOCA = 0;
167  ibinDOCA = (idocaHit - fdocaLE) / fdocaBS;
168 
169  if (ibinEA >= fnEntaBinG || ibinDOCA >= fnDocaBinG) continue; //bin stats from 0
170 
171  if (IsLocalBin) ibinEA = fEntaBinNums.at(ibinEA);
172 
173  if (layer->at(j) < 8) {
174  hILDocaEntaG->Fill(idocaHit, ieaHit);
175  if (IsLocalBin)hILDocaEntaL->Fill(idocaHit, ieaHit);
176  hILdEdx_all->Fill(dedxhit->at(j));
177  hILdEdxhitInEntaDocaBin[ibinEA][ibinDOCA]->Fill(dedxhit->at(j));
178  } else {
179  hOLDocaEntaG->Fill(idocaHit, ieaHit);
180  if (IsLocalBin)hOLDocaEntaL->Fill(idocaHit, ieaHit);
181  hOLdEdx_all->Fill(dedxhit->at(j));
182  hOLdEdxhitInEntaDocaBin[ibinEA][ibinDOCA]->Fill(dedxhit->at(j));
183  }
184  }
185  }
186 
187 
188  if (IsMakePlots) {
189  TCanvas* ctmpde = new TCanvas("DocavsEnta", "Doca vs Enta distributions", 800, 400);
190  if (IsLocalBin) {
191  ctmpde->SetCanvasSize(800, 800);
192  ctmpde->Divide(2, 2);
193  ctmpde->cd(1); hILDocaEntaG->Draw("colz");
194  ctmpde->cd(2); hOLDocaEntaG->Draw("colz");
195  ctmpde->cd(3); hILDocaEntaL->Draw("colz");
196  ctmpde->cd(4); hOLDocaEntaL->Draw("colz");
197  } else {
198  ctmpde->Divide(2, 1);
199  ctmpde->cd(1); hILDocaEntaG->Draw("colz");
200  ctmpde->cd(2); hOLDocaEntaG->Draw("colz");
201  }
202  ctmpde->SaveAs(Form("DocavsEnta_%s.pdf", fSetPrefix.data()));
203  ctmpde->SaveAs(Form("DocavsEnta_%s.root", fSetPrefix.data()));
204 
205  TCanvas* ctem = new TCanvas("Layerhisto", "Inner and Outer Histo", 600, 600);
206  hOLdEdx_all->Draw("histo");
207  hILdEdx_all->SetMarkerColor(kRed);
208  hILdEdx_all->Draw("same histo");
209  ctem->SaveAs(Form("Layerhistodedxhit_TwoDCorr_%s.pdf", fSetPrefix.data()));
210  }
211 
212 
213  //Calculationg 5-75% global truncation mean
214  double InsumPer5 = 0.0, InsumPer75 = 0.0;
215  double OutsumPer5 = 0.0, OutsumPer75 = 0.0;
216  double InLayInt = hILdEdx_all->Integral();
217  double OutLayInt = hOLdEdx_all->Integral();
218 
219  Int_t lBinInLayer = 1, hBinInLayer = 1;
220  Int_t lBinOutLayer = 1, hBinOutLayer = 1;
221 
222  for (int ibin = 1; ibin <= hILdEdx_all->GetNbinsX(); ibin++) {
223 
224  if (InsumPer5 <= 0.05 * InLayInt) {
225  InsumPer5 += hILdEdx_all->GetBinContent(ibin);
226  lBinInLayer = ibin;
227  }
228 
229  if (InsumPer75 <= 0.75 * InLayInt) {
230  InsumPer75 += hILdEdx_all->GetBinContent(ibin);
231  hBinInLayer = ibin;
232  }
233 
234  if (OutsumPer5 <= 0.05 * OutLayInt) {
235  OutsumPer5 += hOLdEdx_all->GetBinContent(ibin);
236  lBinOutLayer = ibin;
237  }
238 
239  if (OutsumPer75 <= 0.75 * OutLayInt) {
240  OutsumPer75 += hOLdEdx_all->GetBinContent(ibin);
241  hBinOutLayer = ibin;
242  }
243  }
244 
245  short version = 1;
246  TCanvas* ctmp = new TCanvas("tmp", "tmp", 1200, 1200);
247  ctmp->Divide(4, 4);
248  std::stringstream psname; psname << Form("dedx_2dcell_%s.pdf[", fSetPrefix.data());
249  TLine* tl = new TLine();
250  tl->SetLineColor(kRed);
251 
252  TBox* tb = new TBox();
253 
254  if (IsMakePlots) {
255  ctmp->Print(psname.str().c_str());
256  psname.str(""); psname << Form("dedx_2dcell_%s.pdf", fSetPrefix.data());
257  }
258 
259  gStyle->SetOptStat("nemriou");
260  TH1D* htemp = NULL;
261  std::vector<TH2F> twodcors;
262  TH2F tempTwoD = TH2F("tempTwoD", "dE/dx in bins of DOCA/Enta;DOCA;Entrance Angle", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
263  RmapEntaValue);
264  TH2F twodcor = TH2F("twodcorrection", "dE/dx in bins of DOCA/Enta;DOCA;Entrance Angle", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG,
265  feaLE, feaUE);
266 
267  std::ofstream file2DILout;
268 
269  for (int iIOLayer = 0; iIOLayer <= 1; iIOLayer++) {
270 
271  file2DILout.open(Form("f2D_Constants_Local_Layer%d_%s.txt", iIOLayer, fSetPrefix.data()));
272 
273  int startfrom = 1, endat = 1;
274  if (iIOLayer == 0) {
275  startfrom = lBinInLayer; endat = hBinInLayer;
276  } else {
277  startfrom = lBinOutLayer; endat = hBinOutLayer;
278  }
279 
280  //std::cout << "Layer I/O # = " << iIOLayer << std::endl;
281  for (int iea = 1; iea <= fnEntaBinL; iea++) {
282 
283  Int_t ieaprime = iea; //rotation symmtery for 1<->3 and 4<->2
284  if (IsRS)ieaprime = GetRotationSymmericBin(fnEntaBinL, iea);
285 
286  for (int idoca = 1; idoca <= fnDocaBinL; idoca++) {
287 
288  if (iIOLayer == 0)
289  htemp = (TH1D*)hILdEdxhitInEntaDocaBin[ieaprime - 1][idoca - 1]->Clone(Form("hL%d_Ea%d_Doca%d", iIOLayer, iea, idoca));
290  else if (iIOLayer == 1)
291  htemp = (TH1D*)hOLdEdxhitInEntaDocaBin[ieaprime - 1][idoca - 1]->Clone(Form("hL%d_Ea%d_Doca%d", iIOLayer, iea, idoca));
292  else continue;
293 
294  double truncMean = 1.0;
295  if (htemp->GetEntries() < 400) truncMean = 1.0; //low stats
296  else {
297  double binweights = 0.0;
298  int sumofbc = 0;
299  for (int ibin = startfrom; ibin <= endat; ibin++) {
300  //std::cout << " dedxhit bin = " << ibin << ", Entries =" << htemp->GetBinContent(ibin) << std::endl;
301  if (htemp->GetBinContent(ibin) > 0) {
302  binweights += (htemp->GetBinContent(ibin) * htemp->GetBinCenter(ibin));
303  sumofbc += htemp->GetBinContent(ibin);
304  }
305  }
306  if (sumofbc > 0)truncMean = (double)(binweights / sumofbc);
307  else truncMean = 1.0;
308  }
309 
310  if (truncMean <= 0)truncMean = 1.0; //protection only
311  tempTwoD.SetBinContent(idoca, iea, truncMean); //binning starts at 1
312 
313  if (IsMakePlots) {
314  ctmp->cd(((idoca - 1) % 16) + 1);
315  htemp->SetTitle(Form("%s, mean = %0.4f", htemp->GetTitle(), truncMean));
316  htemp->SetFillColorAlpha(kGreen, 0.30);
317  if (truncMean >= 1.02 || truncMean <= 0.98)htemp->SetFillColor(kYellow);
318  if (truncMean >= 1.05 || truncMean <= 0.95)htemp->SetFillColor(kOrange);
319  if (truncMean >= 1.10 || truncMean <= 0.90)htemp->SetFillColor(kRed);
320  htemp->DrawClone(); //clone is nessesory for pointer survival
321  tl->SetLineColor(kRed);
322  tl->SetX1(truncMean); tl->SetX2(truncMean);
323  tl->SetY1(0); tl->SetY2(htemp->GetMaximum());
324  tl->DrawClone("same");
325 
326  tb->SetLineColor(kPink);
327  tb->SetLineStyle(6);
328  tb->SetFillColorAlpha(kPink, 0.15);
329  tb->SetX1(htemp->GetBinLowEdge(startfrom));
330  tb->SetY1(0);
331  tb->SetX2(htemp->GetBinLowEdge(endat));
332  tb->SetY2(htemp->GetMaximum());
333  tb->DrawClone("same");
334  if (idoca % 16 == 0)ctmp->Print(psname.str().c_str());
335  }
336 
337  file2DILout << iea << ", \t" << idoca << ", \t" << truncMean << "\n";
338  htemp->Reset();
339  }
340  }
341 
342  //conversion from local to global
343  ibinEA = 0, ibinDOCA = 0;
344  for (int iea = 0; iea < fnEntaBinG; iea++) {
345  ibinEA = iea;
346  if (IsLocalBin)ibinEA = fEntaBinNums.at(iea);
347  for (int idoca = 0; idoca < fnDocaBinG; idoca++) {
348  ibinDOCA = idoca;
349  twodcor.SetBinContent(idoca + 1, iea + 1, tempTwoD.GetBinContent(ibinDOCA + 1, ibinEA + 1));
350  }
351  }
352 
353  if (iIOLayer == 0)twodcor.SetTitle("InnerLayer: dE/dx in bins of DOCA/Enta");
354  else if (iIOLayer == 1)twodcor.SetTitle("OuterLayer: dE/dx in bins of DOCA/Enta");
355  twodcors.push_back(twodcor);
356  twodcor.Reset();
357  tempTwoD.Reset();
358  file2DILout.close();
359  }
360 
361  if (IsMakePlots) {
362  psname.str(""); psname << Form("dedx_2dcell_%s.pdf]", fSetPrefix.data());
363  ctmp->Print(psname.str().c_str());
364 
365  // //Drawing final constants
366  TCanvas* cconst = new TCanvas("FinalConstantHistoMap", "Inner and Outer Histo", 800, 400);
367  cconst->Divide(2, 1);
368  cconst->cd(1); twodcors.at(0).Draw("colz");
369  cconst->cd(2); twodcors.at(1).Draw("colz");
370  cconst->SaveAs(Form("Final2DConstantMap_wGlobalBin_%s.pdf", fSetPrefix.data()));
371  cconst->SaveAs(Form("Final2DConstantMap_wGlobalBin_%s.root", fSetPrefix.data()));
372  }
373 
374 
375  B2INFO("dE/dx calibration done for 2D correction");
376  CDCDedx2DCell* gain = new CDCDedx2DCell(version, twodcors);
377  saveCalibration(gain, "CDCDedx2DCell");
378 
379  delete htemp;
380  delete ctmp;
381  delete tl;
382  delete tb;
383  return c_OK;
384 
385 }
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)
funtion 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)
funtion 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.