Belle II Software  release-05-02-19
CDCDedx2DCellAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: jvbennett *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <reconstruction/calibration/CDCDedx2DCellAlgorithm.h>
12 
13 #include <reconstruction/dbobjects/CDCDedx2DCell.h>
14 
15 #include <TBox.h>
16 #include <TCanvas.h>
17 #include <TH2F.h>
18 #include <TLine.h>
19 #include <TMath.h>
20 #include <TStyle.h>
21 
22 #include <algorithm>
23 #include <fstream>
24 
25 using namespace Belle2;
26 
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
33  CalibrationAlgorithm("CDCDedxElectronCollector"),
34  fnEntaBinG(128),
35  fnDocaBinG(64),
36  fnEntaBinL(64),
37  fnDocaBinL(28),
38  feaLE(-TMath::Pi() / 2),
39  feaUE(+TMath::Pi() / 2),
40  fdocaLE(-1.50),
41  fdocaUE(1.50),
42  fSetPrefix("_it0"),
43  IsLocalBin(true),
44  IsMakePlots(false),
45  IsRS(true)
46 {
47  // Set module properties
48  setDescription("A calibration algorithm for the CDC dE/dx two dimensional correction");
49 }
50 
51 //-----------------------------------------------------------------
52 // Run the calibration
53 //-----------------------------------------------------------------
54 
56 {
57 
58  //reading electron collector TREE
59  auto ttree = getObjectPtr<TTree>("tree");
60  if (ttree->GetEntries() < 100)return c_NotEnoughData;
61 
62  std::vector<double>* dedxhit = 0, *doca = 0, *enta = 0;
63  std::vector<int>* layer = 0;
64 
65  ttree->SetBranchAddress("dedxhit", &dedxhit);
66  ttree->SetBranchAddress("layer", &layer);
67  ttree->SetBranchAddress("docaRS", &doca);
68  ttree->SetBranchAddress("entaRS", &enta);
69 
70  // Setting up bins for doca and entra angle
71  feaBS = (feaUE - feaLE) / fnEntaBinG;
73 
74  std::vector<int> globalbinsEnta, globalbinsDoca;
75  for (int ibin = 0; ibin < fnEntaBinG; ibin++)globalbinsEnta.push_back(ibin);
76  for (int ibin = 0; ibin < fnDocaBinG; ibin++)globalbinsDoca.push_back(ibin);
77 
78  if (!IsLocalBin) {
79  fEntaBinNums = globalbinsEnta;
81  } else {
83  fnEntaBinL = fEntaBinNums.at(fEntaBinNums.size() - 1) + 1;
84  }
85 
86  fDocaBinNums = globalbinsDoca;
88 
89  std::vector<std::vector<TH1F*>> hILdEdxhitInEntaDocaBin(fnEntaBinL, std::vector<TH1F*>(fnDocaBinL, 0));
90  std::vector<std::vector<TH1F*>> hOLdEdxhitInEntaDocaBin(fnEntaBinL, std::vector<TH1F*>(fnDocaBinL, 0));
91  Double_t ifeaLE = 0, ifeaUE = 0, ifdocaLE = 0, ifdocaUE = 0;
92 
93  for (int iea = 0; iea < fnEntaBinL; iea++) {
94  if (IsLocalBin) {
95  ifeaLE = fEntaBinValues.at(iea);
96  ifeaUE = fEntaBinValues.at(iea + 1);
97  } else {
98  ifeaLE = iea * feaBS - feaUE; //- because of -ive range shifting
99  ifeaUE = ifeaLE + feaBS;
100  }
101 
102  for (int idoca = 0; idoca < fnDocaBinL; idoca++) {
103 
104  ifdocaLE = idoca * fdocaBS - fdocaUE;
105  ifdocaUE = ifdocaLE + fdocaBS;
106 
107  hILdEdxhitInEntaDocaBin[iea][idoca] = new TH1F(Form("hILdEdxhitInEntaBin%dDocaBin%d", iea, idoca), "bla-bla", 500, 0., 5.);
108  hILdEdxhitInEntaDocaBin[iea][idoca]->SetTitle(Form("IL: EntA = (%0.03f to %0.03f) and Doca = (%0.03f to %0.03f)", ifeaLE,
109  ifeaUE, ifdocaLE, ifdocaUE));
110  hILdEdxhitInEntaDocaBin[iea][idoca]->GetXaxis()->SetTitle("dedxhits in Inner Layer");
111  hILdEdxhitInEntaDocaBin[iea][idoca]->GetYaxis()->SetTitle("Entries");
112 
113  hOLdEdxhitInEntaDocaBin[iea][idoca] = new TH1F(Form("hOLdEdxhitInEntaBin%dDocaBin%d", iea, idoca), "bla-bla", 500, 0., 5.);
114  hOLdEdxhitInEntaDocaBin[iea][idoca]->SetTitle(Form("OL: EntA = (%0.03f to %0.03f) and Doca = (%0.03f to %0.03f)", ifeaLE,
115  ifeaUE, ifdocaLE, ifdocaUE));
116  hOLdEdxhitInEntaDocaBin[iea][idoca]->GetXaxis()->SetTitle("dedxhits in Outer Layer");
117  hOLdEdxhitInEntaDocaBin[iea][idoca]->GetYaxis()->SetTitle("Entries");
118 
119  }
120  }
121 
122  //Doca vs Enta stats
123  TH2D* hILDocaEntaG = new TH2D("hILDocaEntaG", "Doca vs EntA: Inner Layer", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG, feaLE, feaUE);
124  hILDocaEntaG->GetXaxis()->SetTitle("Doca");
125  hILDocaEntaG->GetYaxis()->SetTitle("Entrance angle (#theta)");
126 
127  TH2D* hOLDocaEntaG = new TH2D("hOLDocaEntaG", "Doca vs EntA: Outer Layer", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG, feaLE, feaUE);
128  hOLDocaEntaG->GetXaxis()->SetTitle("Doca");
129  hOLDocaEntaG->GetYaxis()->SetTitle("Entrance angle (#theta)");
130 
131  //when local enta angle is demanded
132  Double_t* RmapEntaValue = &fEntaBinValues[0];
133  TH2D* hILDocaEntaL = new TH2D("hILDocaEntaL", "Doca vs EntA: Inner Layer (rebin)", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
134  RmapEntaValue);
135  hILDocaEntaL->GetXaxis()->SetTitle("Doca");
136  hILDocaEntaL->GetYaxis()->SetTitle("Entrance angle (#theta)");
137 
138  TH2D* hOLDocaEntaL = new TH2D("hOLDocaEntaL", "Doca vs EntA: Outer Layer (rebin)", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
139  RmapEntaValue);
140  hOLDocaEntaL->GetXaxis()->SetTitle("Doca");
141  hOLDocaEntaL->GetYaxis()->SetTitle("Entrance angle (#theta)");
142 
143  TH1D* hILdEdx_all = new TH1D("hILdEdx_all", "", 500, 0., 5.);
144  TH1D* hOLdEdx_all = new TH1D("hOLdEdx_all", "", 500, 0., 5.);
145 
146 
147  Int_t ibinEA = 0, ibinDOCA = 0;
148  for (int i = 0; i < ttree->GetEntries(); ++i) {
149 
150  ttree->GetEvent(i);
151 
152  for (unsigned int j = 0; j < dedxhit->size(); ++j) {
153 
154  if (dedxhit->at(j) == 0) continue;
155 
156  Double_t ieaHit = enta->at(j);
157  if (ieaHit < -TMath::Pi() / 2.0) ieaHit += TMath::Pi() / 2.0;
158  else if (ieaHit > TMath::Pi() / 2.0) ieaHit -= TMath::Pi() / 2.0;
159  if (abs(ieaHit) > TMath::Pi() / 2.0) continue;
160 
161  Double_t idocaHit = doca->at(j);
162  if (abs(idocaHit) > 1.50) continue;
163 
164  //Bin corresponds to enta and doca value
165  ibinEA = 0;
166  ibinEA = (ieaHit - feaLE) / feaBS ; //from 0
167 
168  ibinDOCA = 0;
169  ibinDOCA = (idocaHit - fdocaLE) / fdocaBS;
170 
171  if (ibinEA >= fnEntaBinG || ibinDOCA >= fnDocaBinG) continue; //bin stats from 0
172 
173  if (IsLocalBin) ibinEA = fEntaBinNums.at(ibinEA);
174 
175  if (layer->at(j) < 8) {
176  hILDocaEntaG->Fill(idocaHit, ieaHit);
177  if (IsLocalBin)hILDocaEntaL->Fill(idocaHit, ieaHit);
178  hILdEdx_all->Fill(dedxhit->at(j));
179  hILdEdxhitInEntaDocaBin[ibinEA][ibinDOCA]->Fill(dedxhit->at(j));
180  } else {
181  hOLDocaEntaG->Fill(idocaHit, ieaHit);
182  if (IsLocalBin)hOLDocaEntaL->Fill(idocaHit, ieaHit);
183  hOLdEdx_all->Fill(dedxhit->at(j));
184  hOLdEdxhitInEntaDocaBin[ibinEA][ibinDOCA]->Fill(dedxhit->at(j));
185  }
186  }
187  }
188 
189 
190  if (IsMakePlots) {
191  TCanvas* ctmpde = new TCanvas("DocavsEnta", "Doca vs Enta distributions", 800, 400);
192  if (IsLocalBin) {
193  ctmpde->SetCanvasSize(800, 800);
194  ctmpde->Divide(2, 2);
195  ctmpde->cd(1); hILDocaEntaG->Draw("colz");
196  ctmpde->cd(2); hOLDocaEntaG->Draw("colz");
197  ctmpde->cd(3); hILDocaEntaL->Draw("colz");
198  ctmpde->cd(4); hOLDocaEntaL->Draw("colz");
199  } else {
200  ctmpde->Divide(2, 1);
201  ctmpde->cd(1); hILDocaEntaG->Draw("colz");
202  ctmpde->cd(2); hOLDocaEntaG->Draw("colz");
203  }
204  ctmpde->SaveAs(Form("DocavsEnta_%s.pdf", fSetPrefix.data()));
205  ctmpde->SaveAs(Form("DocavsEnta_%s.root", fSetPrefix.data()));
206 
207  TCanvas* ctem = new TCanvas("Layerhisto", "Inner and Outer Histo", 600, 600);
208  hOLdEdx_all->Draw("histo");
209  hILdEdx_all->SetMarkerColor(kRed);
210  hILdEdx_all->Draw("same histo");
211  ctem->SaveAs(Form("Layerhistodedxhit_TwoDCorr_%s.pdf", fSetPrefix.data()));
212  }
213 
214 
215  //Calculationg 5-75% global truncation mean
216  double InsumPer5 = 0.0, InsumPer75 = 0.0;
217  double OutsumPer5 = 0.0, OutsumPer75 = 0.0;
218  double InLayInt = hILdEdx_all->Integral();
219  double OutLayInt = hOLdEdx_all->Integral();
220 
221  Int_t lBinInLayer = 1, hBinInLayer = 1;
222  Int_t lBinOutLayer = 1, hBinOutLayer = 1;
223 
224  for (int ibin = 1; ibin <= hILdEdx_all->GetNbinsX(); ibin++) {
225 
226  if (InsumPer5 <= 0.05 * InLayInt) {
227  InsumPer5 += hILdEdx_all->GetBinContent(ibin);
228  lBinInLayer = ibin;
229  }
230 
231  if (InsumPer75 <= 0.75 * InLayInt) {
232  InsumPer75 += hILdEdx_all->GetBinContent(ibin);
233  hBinInLayer = ibin;
234  }
235 
236  if (OutsumPer5 <= 0.05 * OutLayInt) {
237  OutsumPer5 += hOLdEdx_all->GetBinContent(ibin);
238  lBinOutLayer = ibin;
239  }
240 
241  if (OutsumPer75 <= 0.75 * OutLayInt) {
242  OutsumPer75 += hOLdEdx_all->GetBinContent(ibin);
243  hBinOutLayer = ibin;
244  }
245  }
246 
247  short version = 1;
248  TCanvas* ctmp = new TCanvas("tmp", "tmp", 1200, 1200);
249  ctmp->Divide(4, 4);
250  std::stringstream psname; psname << Form("dedx_2dcell_%s.pdf[", fSetPrefix.data());
251  TLine* tl = new TLine();
252  tl->SetLineColor(kRed);
253 
254  TBox* tb = new TBox();
255 
256  if (IsMakePlots) {
257  ctmp->Print(psname.str().c_str());
258  psname.str(""); psname << Form("dedx_2dcell_%s.pdf", fSetPrefix.data());
259  }
260 
261  gStyle->SetOptStat("nemriou");
262  TH1D* htemp = NULL;
263  std::vector<TH2F> twodcors;
264  TH2F tempTwoD = TH2F("tempTwoD", "dE/dx in bins of DOCA/Enta;DOCA;Entrance Angle", fnDocaBinL, fdocaLE, fdocaUE, fnEntaBinL,
265  RmapEntaValue);
266  TH2F twodcor = TH2F("twodcorrection", "dE/dx in bins of DOCA/Enta;DOCA;Entrance Angle", fnDocaBinG, fdocaLE, fdocaUE, fnEntaBinG,
267  feaLE, feaUE);
268 
269  std::ofstream file2DILout;
270 
271  for (int iIOLayer = 0; iIOLayer <= 1; iIOLayer++) {
272 
273  file2DILout.open(Form("f2D_Constants_Local_Layer%d_%s.txt", iIOLayer, fSetPrefix.data()));
274 
275  int startfrom = 1, endat = 1;
276  if (iIOLayer == 0) {
277  startfrom = lBinInLayer; endat = hBinInLayer;
278  } else {
279  startfrom = lBinOutLayer; endat = hBinOutLayer;
280  }
281 
282  //std::cout << "Layer I/O # = " << iIOLayer << std::endl;
283  for (int iea = 1; iea <= fnEntaBinL; iea++) {
284 
285  Int_t ieaprime = iea; //rotation symmtery for 1<->3 and 4<->2
286  if (IsRS)ieaprime = GetRotationSymmericBin(fnEntaBinL, iea);
287 
288  for (int idoca = 1; idoca <= fnDocaBinL; idoca++) {
289 
290  if (iIOLayer == 0)
291  htemp = (TH1D*)hILdEdxhitInEntaDocaBin[ieaprime - 1][idoca - 1]->Clone(Form("hL%d_Ea%d_Doca%d", iIOLayer, iea, idoca));
292  else if (iIOLayer == 1)
293  htemp = (TH1D*)hOLdEdxhitInEntaDocaBin[ieaprime - 1][idoca - 1]->Clone(Form("hL%d_Ea%d_Doca%d", iIOLayer, iea, idoca));
294  else continue;
295 
296  double truncMean = 1.0;
297  if (htemp->GetEntries() < 400) truncMean = 1.0; //low stats
298  else {
299  double binweights = 0.0;
300  int sumofbc = 0;
301  for (int ibin = startfrom; ibin <= endat; ibin++) {
302  //std::cout << " dedxhit bin = " << ibin << ", Entries =" << htemp->GetBinContent(ibin) << std::endl;
303  if (htemp->GetBinContent(ibin) > 0) {
304  binweights += (htemp->GetBinContent(ibin) * htemp->GetBinCenter(ibin));
305  sumofbc += htemp->GetBinContent(ibin);
306  }
307  }
308  if (sumofbc > 0)truncMean = (double)(binweights / sumofbc);
309  else truncMean = 1.0;
310  }
311 
312  if (truncMean <= 0)truncMean = 1.0; //protection only
313  tempTwoD.SetBinContent(idoca, iea, truncMean); //binning starts at 1
314 
315  if (IsMakePlots) {
316  ctmp->cd(((idoca - 1) % 16) + 1);
317  htemp->SetTitle(Form("%s, mean = %0.4f", htemp->GetTitle(), truncMean));
318  htemp->SetFillColorAlpha(kGreen, 0.30);
319  if (truncMean >= 1.02 || truncMean <= 0.98)htemp->SetFillColor(kYellow);
320  if (truncMean >= 1.05 || truncMean <= 0.95)htemp->SetFillColor(kOrange);
321  if (truncMean >= 1.10 || truncMean <= 0.90)htemp->SetFillColor(kRed);
322  htemp->DrawClone(); //clone is nessesory for pointer survival
323  tl->SetLineColor(kRed);
324  tl->SetX1(truncMean); tl->SetX2(truncMean);
325  tl->SetY1(0); tl->SetY2(htemp->GetMaximum());
326  tl->DrawClone("same");
327 
328  tb->SetLineColor(kPink);
329  tb->SetLineStyle(6);
330  tb->SetFillColorAlpha(kPink, 0.15);
331  tb->SetX1(htemp->GetBinLowEdge(startfrom));
332  tb->SetY1(0);
333  tb->SetX2(htemp->GetBinLowEdge(endat));
334  tb->SetY2(htemp->GetMaximum());
335  tb->DrawClone("same");
336  if (idoca % 16 == 0)ctmp->Print(psname.str().c_str());
337  }
338 
339  file2DILout << iea << ", \t" << idoca << ", \t" << truncMean << "\n";
340  htemp->Reset();
341  }
342  }
343 
344  //conversion from local to global
345  ibinEA = 0, ibinDOCA = 0;
346  for (int iea = 0; iea < fnEntaBinG; iea++) {
347  ibinEA = iea;
348  if (IsLocalBin)ibinEA = fEntaBinNums.at(iea);
349  for (int idoca = 0; idoca < fnDocaBinG; idoca++) {
350  ibinDOCA = idoca;
351  twodcor.SetBinContent(idoca + 1, iea + 1, tempTwoD.GetBinContent(ibinDOCA + 1, ibinEA + 1));
352  }
353  }
354 
355  if (iIOLayer == 0)twodcor.SetTitle("InnerLayer: dE/dx in bins of DOCA/Enta");
356  else if (iIOLayer == 1)twodcor.SetTitle("OuterLayer: dE/dx in bins of DOCA/Enta");
357  twodcors.push_back(twodcor);
358  twodcor.Reset();
359  tempTwoD.Reset();
360  file2DILout.close();
361  }
362 
363  if (IsMakePlots) {
364  psname.str(""); psname << Form("dedx_2dcell_%s.pdf]", fSetPrefix.data());
365  ctmp->Print(psname.str().c_str());
366 
367  // //Drawing final constants
368  TCanvas* cconst = new TCanvas("FinalConstantHistoMap", "Inner and Outer Histo", 800, 400);
369  cconst->Divide(2, 1);
370  cconst->cd(1); twodcors.at(0).Draw("colz");
371  cconst->cd(2); twodcors.at(1).Draw("colz");
372  cconst->SaveAs(Form("Final2DConstantMap_wGlobalBin_%s.pdf", fSetPrefix.data()));
373  cconst->SaveAs(Form("Final2DConstantMap_wGlobalBin_%s.root", fSetPrefix.data()));
374  }
375 
376 
377  B2INFO("dE/dx calibration done for 2D correction");
378  CDCDedx2DCell* gain = new CDCDedx2DCell(version, twodcors);
379  saveCalibration(gain, "CDCDedx2DCell");
380 
381  delete htemp;
382  delete ctmp;
383  delete tl;
384  delete tb;
385  return c_OK;
386 
387 }
Belle2::CDCDedx2DCellAlgorithm::fdocaUE
Double_t fdocaUE
Upper edge of doca.
Definition: CDCDedx2DCellAlgorithm.h:168
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::CDCDedx2DCellAlgorithm::IsRS
bool IsRS
if rotation symmtery requested
Definition: CDCDedx2DCellAlgorithm.h:180
Belle2::CDCDedx2DCellAlgorithm::IsMakePlots
bool IsMakePlots
produce plots for status
Definition: CDCDedx2DCellAlgorithm.h:179
Belle2::CDCDedx2DCellAlgorithm::GetRotationSymmericBin
int GetRotationSymmericBin(int nbin, int ibin)
funtion to set rotation symmetry
Definition: CDCDedx2DCellAlgorithm.h:82
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::CDCDedx2DCellAlgorithm::fnDocaBinG
int fnDocaBinG
doca angle bins, Global
Definition: CDCDedx2DCellAlgorithm.h:159
Belle2::CDCDedx2DCellAlgorithm::fDocaBinNums
std::vector< int > fDocaBinNums
Vector for enta asym bin #.
Definition: CDCDedx2DCellAlgorithm.h:172
Belle2::CDCDedx2DCellAlgorithm::fnDocaBinL
int fnDocaBinL
doca angle bins, Local
Definition: CDCDedx2DCellAlgorithm.h:161
Belle2::CDCDedx2DCellAlgorithm::calibrate
virtual EResult calibrate() override
2D Cell algorithm algorithm
Definition: CDCDedx2DCellAlgorithm.cc:55
Belle2::CDCDedx2DCellAlgorithm::fSetPrefix
std::string fSetPrefix
prefix to filename
Definition: CDCDedx2DCellAlgorithm.h:171
Belle2::CDCDedx2DCellAlgorithm::fdocaLE
Double_t fdocaLE
Lower edge of doca.
Definition: CDCDedx2DCellAlgorithm.h:167
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CDCDedx2DCellAlgorithm::fEntaBinValues
std::vector< double > fEntaBinValues
Vector for doca asym bin values.
Definition: CDCDedx2DCellAlgorithm.h:176
Belle2::CDCDedx2DCell
dE/dx wire gain calibration constants
Definition: CDCDedx2DCell.h:36
Belle2::CDCDedx2DCellAlgorithm::fEntaBinNums
std::vector< int > fEntaBinNums
Vector for enta asym bin values.
Definition: CDCDedx2DCellAlgorithm.h:175
Belle2::CDCDedx2DCellAlgorithm::feaUE
Double_t feaUE
Upper edge of enta angle.
Definition: CDCDedx2DCellAlgorithm.h:164
Belle2::CDCDedx2DCellAlgorithm::feaLE
Double_t feaLE
Lower edge of enta angle.
Definition: CDCDedx2DCellAlgorithm.h:163
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::CDCDedx2DCellAlgorithm::fdocaBS
Double_t fdocaBS
Binwidth edge of doca.
Definition: CDCDedx2DCellAlgorithm.h:169
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::CDCDedx2DCellAlgorithm::fnEntaBinL
int fnEntaBinL
etna angle bins, Local
Definition: CDCDedx2DCellAlgorithm.h:160
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::CDCDedx2DCellAlgorithm::CDCDedx2DCellAlgorithm
CDCDedx2DCellAlgorithm()
Constructor: Sets the description the properties and the parameters of the algorithm.
Definition: CDCDedx2DCellAlgorithm.cc:32
Belle2::CDCDedx2DCellAlgorithm::IsLocalBin
bool IsLocalBin
if local variable bin requested
Definition: CDCDedx2DCellAlgorithm.h:178
Belle2::CDCDedx2DCellAlgorithm::feaBS
Double_t feaBS
Binwidth edge of enta angle.
Definition: CDCDedx2DCellAlgorithm.h:165
Belle2::CDCDedx2DCellAlgorithm::GetVariableBin
void GetVariableBin(int nbin, std::vector< int > &nBinEnta0to100Per)
funtion to set variable bins
Definition: CDCDedx2DCellAlgorithm.h:97
Belle2::CDCDedx2DCellAlgorithm::fnEntaBinG
int fnEntaBinG
Save arithmetic and truncated mean for the 'dedx' values.
Definition: CDCDedx2DCellAlgorithm.h:158