Belle II Software  release-08-02-06
CDCDedx1DCellAlgorithm.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/CDCDedx1DCellAlgorithm.h>
10 using namespace Belle2;
11 
12 //-----------------------------------------------------------------
13 // Implementation
14 //-----------------------------------------------------------------
16  CalibrationAlgorithm("CDCDedxElectronCollector"),
17  m_eaMin(-TMath::Pi() / 2),
18  m_eaMax(+TMath::Pi() / 2),
19  m_eaBin(316),
20  m_dedxMin(0.0),
21  m_dedxMax(5.0),
22  m_dedxBin(250),
23  m_ptMax(8.0),
24  m_cosMax(1.0),
25  m_truncMin(0.05),
26  m_truncMax(0.75),
27  m_binSplit(3),
28  m_chargeType(0),
29  m_adjustFac(1.00),
30  isFixTrunc(false),
31  isVarBins(true),
32  isRotSymm(false),
33  isMakePlots(true),
34  isPrintLog(false),
35  isMerge(true),
36  m_suffix("")
37 {
38  // Set module properties
39  setDescription("A calibration algorithm for the CDC dE/dx entrance angle cleanup correction");
40 }
41 
42 //-----------------------------------------------------------------
43 // Run the calibration
44 //-----------------------------------------------------------------
46 {
47 
48  getExpRunInfo();
49 
50  if (!m_DBOneDCell.isValid())
51  B2FATAL("There is no valid previous payload for CDCDedx1DCell");
52 
53  //reading radiative electron collector TREE
54  auto ttree = getObjectPtr<TTree>("tree");
55  if (!ttree) return c_NotEnoughData;
56 
57  std::vector<double>* dedxhit = 0, *enta = 0;
58  std::vector<int>* layer = 0;
59  double pt = 0, costh = 0;
60 
61  ttree->SetBranchAddress("dedxhit", &dedxhit);
62  ttree->SetBranchAddress("entaRS", &enta);
63  ttree->SetBranchAddress("layer", &layer);
64  ttree->SetBranchAddress("pt", &pt);
65  ttree->SetBranchAddress("costh", &costh);
66 
67  //repair nbins if they are not divisible accordingly
69  m_eaBW = (m_eaMax - m_eaMin) / m_eaBin;
70 
71  //Settings of variables bins
73 
74  if (isPrintLog) {
75  B2INFO("inner layers bins: " << m_eaBinLocal[0]);
76  B2INFO("outer layers bins: " << m_eaBinLocal[1]);
77  }
78 
79  // dedxhit vector to store dE/dx values for each enta bin
80  std::vector<TH1D*> hdedxhit[2];
81  TH1D* hdedxlay[2];
82  TH1D* hentalay[2];
83 
84  TH2D* hptcosth = new TH2D("hptcosth", "pt vs costh dist;pt;costh", 1000, -8.0, 8.0, 1000, -1.0, 1.0);
85 
86  defineHisto(hdedxhit, hdedxlay, hentalay);
87 
88  //Star filling histogram defined above
89  for (int i = 0; i < ttree->GetEntries(); ++i) {
90 
91  ttree->GetEvent(i);
92 
93  if (abs(costh) > m_cosMax) continue;
94 
95  // remove wide angle bhabha tracks
96  // double mom = pt/sqrt(1-costh*costh);
97  // if(abs(pt)<2.4 && abs(mom)>3.6)continue;
98 
99  if (abs(pt) > m_ptMax) continue;
100 
101  //change to random 10%
102  int rand = gRandom->Integer(100);
103  if (rand < 10) hptcosth->Fill(pt, costh);
104 
105  for (unsigned int j = 0; j < dedxhit->size(); ++j) {
106 
107  if (dedxhit->at(j) == 0) continue;
108 
109  double entaval = enta->at(j);
110  //Mapped bin corresponds to entaval
111  int ibin = std::floor((entaval - m_eaMin) / m_eaBW);
112  if (ibin < 0 || ibin > m_eaBin) continue;
113 
114  int mL = -1;
115  if (layer->at(j) < 8)mL = 0;
116  else mL = 1;
117 
118  hdedxlay[mL]->Fill(dedxhit->at(j));
119  if (rand < 10) hentalay[mL]->Fill(entaval);
120 
121  int jbinea = ibin;
122  if (isVarBins) jbinea = m_binIndex[mL].at(ibin);
123  hdedxhit[mL][jbinea]->Fill(dedxhit->at(j));
124  }
125  }
126 
127  for (int il = 0; il < 2; il++) {
128 
129  int minlay = 0, maxlay = 0;
130 
131  if (isFixTrunc) {
132  getTruncatedBins(hdedxlay[il], minlay, maxlay);
133  hdedxlay[il]->SetTitle(Form("%s;%d;%d", hdedxlay[il]->GetTitle(), minlay, maxlay));
134  }
135 
136  std::vector<double>tempconst;
137  tempconst.reserve(m_eaBinLocal[il]);
138 
139  for (int iea = 0; iea < m_eaBinLocal[il]; iea++) {
140 
141  int jea = iea;
142 
143  // rotation symmtery for 1<->3 and 4<->2 but only symmetric bin
144  if (!isVarBins && isRotSymm) jea = rotationalBin(m_eaBinLocal[il], jea);
145 
146  TH1D* htemp = (TH1D*)hdedxhit[il][jea]->Clone(Form("h_%s_b%d_c", m_label[il].data(), jea));
147 
148  int minbin = 1, maxbin = 1;
149  if (isFixTrunc) {
150  minbin = minlay;
151  maxbin = maxlay;
152  } else {
153  //extract truncation window per bin
154  getTruncatedBins(htemp, minbin, maxbin);
155  }
156 
157  double dedxmean;
158  dedxmean = getTruncationMean(htemp, minbin, maxbin);
159  tempconst.push_back(dedxmean);
160 
161  hdedxhit[il][iea]->SetTitle(Form("%s, #mu_{truc} = %0.5f;%d;%d", hdedxhit[il][iea]->GetTitle(), dedxmean, minbin, maxbin));
162  }
163 
164  //Expending constants
165  std::vector<double>layerconst;
166  layerconst.reserve(m_eaBin);
167 
168  for (int iea = 0; iea < m_eaBin; iea++) {
169  int jea = iea;
170  if (isVarBins) jea = m_binIndex[il].at(iea);
171  layerconst.push_back(tempconst.at(jea));
172  }
173 
174  // plot the rel constants var/sym bins
175  if (isMakePlots) plotRelConst(tempconst, layerconst, il);
176  m_onedcors.push_back(layerconst);
177 
178  layerconst.clear();
179  tempconst.clear();
180  }
181 
182  //Saving final constants
183  createPayload();
184 
185  if (isMakePlots) {
186 
187  //1. dE/dx dist. for entrance angle bins
188  plotdedxHist(hdedxhit);
189 
190  //3. Inner and Outer layer dE/dx distributions
191  plotLayerDist(hdedxlay);
192 
193  //4. entrance angle distribution sym/var bins
194  plotQaPars(hentalay, hptcosth);
195 
196  //6. draw the final constants
197  plotConstants();
198 
199  //7. plot statistics related plots here
200  plotEventStats();
201  }
202 
203  for (int il = 0; il < 2; il++) {
204  delete hentalay[il];
205  delete hdedxlay[il];
206  for (int iea = 0; iea < m_eaBinLocal[il]; iea++)
207  delete hdedxhit[il][iea];
208  }
209 
210  delete hptcosth;
211  return c_OK;
212 }
213 
214 //--------------------------------------------------
216 {
217 
218  int cruns = 0;
219  for (auto expRun : getRunList()) {
220  if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
221  cruns++;
222  }
223 
224  const auto erStart = getRunList()[0];
225  int estart = erStart.first;
226  int rstart = erStart.second;
227 
228  const auto erEnd = getRunList()[cruns - 1];
229  int eend = erEnd.first;
230  int rend = erEnd.second;
231 
232  updateDBObjPtrs(1, rstart, estart);
233 
234  m_runExp = Form("Range (%d:%d,%d:%d)", estart, rstart, eend, rend);
235  if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
236  else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
237 }
238 
239 //--------------------------------------------------
241 {
242 
243  std::map<int, std::vector<double>> bounds;
244  std::map<int, std::vector<int>> steps;
245 
246  const std::array<int, 2> nDev{8, 4};
247  bounds[0] = {0, 108, 123, 133, 158, 183, 193, 208, 316}; //il boundries
248  steps[0] = {9, 3, 2, 1, 1, 2, 3, 9}; //il steps
249  bounds[1] = {0, 38, 158, 278, 316}; //OL boundries
250  steps[1] = {2, 1, 1, 2}; //OL steps
251 
252  for (int il = 0; il < 2; il++) {
253 
254  for (int ibin = 0; ibin <= nDev[il]; ibin++) bounds[il][ibin] = bounds[il][ibin] * m_binSplit;
255 
256  int ieaprime = -1, temp = -99, ibin = 0;
257 
258  double pastbin = m_eaMin;
259  m_binValue[il].push_back(pastbin);
260 
261  for (int iea = 0; iea < m_eaBin; iea++) {
262 
263  if (isVarBins) {
264  if (iea % int(bounds[il][ibin + 1]) == 0 && iea > 0) ibin++;
265  int diff = iea - int(bounds[il][ibin]);
266  if (diff % steps[il][ibin] == 0) ieaprime++;
267  } else ieaprime = iea;
268 
269  m_binIndex[il].push_back(ieaprime);
270 
271  if (ieaprime != temp) {
272  double binwidth = m_eaBW;
273  if (isVarBins) binwidth = m_eaBW * steps[il][ibin];
274  double binvalue = pastbin + binwidth;
275  pastbin = binvalue;
276  if (abs(binvalue) < 1e-5)binvalue = 0;
277  m_binValue[il].push_back(binvalue);
278  }
279  temp = ieaprime;
280  }
281  m_eaBinLocal.push_back(int(m_binValue[il].size()) - 1) ;
282  }
283  if (isMakePlots) plotMergeFactor(bounds, nDev, steps);
284 }
285 
286 //--------------------------------------------------
287 void CDCDedx1DCellAlgorithm::defineHisto(std::vector<TH1D*> hdedxhit[2], TH1D* hdedxlay[2], TH1D* hentalay[2])
288 {
289  for (int il = 0; il < 2; il++) {
290 
291  std::string title = Form("dedxhit dist (%s): %s ; dedxhit;entries", m_label[il].data(), m_runExp.data());
292  hdedxlay[il] = new TH1D(Form("hdedxlay%s", m_label[il].data()), "", m_dedxBin, m_dedxMin, m_dedxMax);
293  hdedxlay[il]->SetTitle(Form("%s", title.data()));
294 
295  Double_t* nvarBins;
296  nvarBins = &m_binValue[il][0];
297 
298  if (isVarBins) title = Form("entaRS dist (variable bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
299  else title = Form("entaRS dist (sym. bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
300 
301  hentalay[il] = new TH1D(Form("hentalay%s", m_label[il].data()), "", m_eaBinLocal[il], nvarBins);
302  hentalay[il]->SetTitle(Form("%s", title.data()));
303 
304  for (int iea = 0; iea < m_eaBinLocal[il]; iea++) {
305 
306  double min = m_binValue[il].at(iea);
307  double max = m_binValue[il].at(iea + 1);
308  double width = max - min;
309 
310  if (isPrintLog) B2INFO("bin: " << iea << " ], min:" << min << " , max: " << max << " , width: " << width);
311 
312  title = Form("%s: entaRS = (%0.03f to %0.03f)", m_label[il].data(), min, max);
313  hdedxhit[il].push_back(new TH1D(Form("hdedxhit_%s_bin%d", m_label[il].data(), iea), "", m_dedxBin, m_dedxMin, m_dedxMax));
314  hdedxhit[il][iea]->SetTitle(Form("%s", title.data()));
315  }
316  }
317 }
318 
319 //--------------------------------------------------
320 void CDCDedx1DCellAlgorithm::getTruncatedBins(TH1D* hist, int& binlow, int& binhigh)
321 {
322 
323  //calculating truncation average
324  double sum = hist->Integral();
325  if (sum <= 0 || hist->GetNbinsX() <= 0) {
326  binlow = 1; binhigh = 1;
327  return ;
328  }
329 
330  binlow = 1.0; binhigh = 1.0;
331  double sumPer5 = 0.0, sumPer75 = 0.0;
332  for (int ibin = 1; ibin <= hist->GetNbinsX(); ibin++) {
333  double bcdedx = hist->GetBinContent(ibin);
334  if (sumPer5 <= m_truncMin * sum) {
335  sumPer5 += bcdedx;
336  binlow = ibin;
337  }
338  if (sumPer75 <= m_truncMax * sum) {
339  sumPer75 += bcdedx;
340  binhigh = ibin;
341  }
342  }
343  return;
344 }
345 
346 //--------------------------
347 double CDCDedx1DCellAlgorithm::getTruncationMean(TH1D* hist, int binlow, int binhigh)
348 {
349 
350  //calculating truncation average
351  if (hist->Integral() < 100) return 1.0;
352 
353  if (binlow <= 0 || binhigh > hist->GetNbinsX())return 1.0;
354 
355  double binweights = 0., sumofbc = 0.;
356  for (int ibin = binlow; ibin <= binhigh; ibin++) {
357  double bcdedx = hist->GetBinContent(ibin);
358  if (bcdedx > 0) {
359  binweights += (bcdedx * hist->GetBinCenter(ibin));
360  sumofbc += bcdedx;
361  }
362  }
363  if (sumofbc > 0) return binweights / sumofbc;
364  else return 1.0;
365 }
366 
367 //--------------------------------------------------
369 {
370 
371  B2INFO("dE/dx one cell calibration: Generating payloads");
372 
373  for (unsigned int il = 0; il < 2; il++) {
374  if (isMerge) {
375  unsigned int nbins = m_DBOneDCell->getNBins(il);
376 
377  if (int(nbins) != m_eaBin)
378  B2ERROR("merging failed because of unmatch bins (old " << m_eaBin << " new " << nbins << ")");
379 
380  for (unsigned int iea = 0; iea < nbins; iea++) {
381  double prev = m_DBOneDCell->getMean(8 * il + 1, iea);
382  m_onedcors[il][iea] *= prev;
383  // m_onedcors[il][iea] /= 0.98;
384  }
385  }
386 
387  if (m_chargeType > 0)
388  for (int ie = 0; ie < m_eaBin / 2; ie++) m_onedcors[il][ie] *= m_adjustFac;
389  if (m_chargeType < 0)
390  for (int ie = m_eaBin / 2; ie < m_eaBin; ie++) m_onedcors[il][ie] *= m_adjustFac;
391 
392  }
393  //Saving constants
394  B2INFO("dE/dx Calibration done for CDCDedx1DCell");
395  CDCDedx1DCell* gain = new CDCDedx1DCell(0, m_onedcors);
396  saveCalibration(gain, "CDCDedx1DCell");
397 }
398 
399 //--------------------------------------------------
400 void CDCDedx1DCellAlgorithm::plotMergeFactor(std::map<int, std::vector<double>> bounds, const std::array<int, 2> nDev,
401  std::map<int, std::vector<int>> steps)
402 {
403 
404  TCanvas cmfactor("cmfactor", "Merging factors", 800, 400);
405  cmfactor.Divide(2, 1);
406 
407  for (int il = 0; il < 2; il++) {
408  Double_t* nvarBins;
409  nvarBins = &bounds[il][0];
410 
411  TH1I* hist = new TH1I(Form("hist_%s", m_label[il].data()), "", nDev[il], nvarBins);
412  hist->SetTitle(Form("Merging factor for %s bins;binindex;merge-factors", m_label[il].data()));
413 
414  for (int ibin = 0; ibin < nDev[il]; ibin++) hist->SetBinContent(ibin + 1, steps[il][ibin]);
415 
416  cmfactor.cd(il + 1);
417  hist->SetFillColor(kYellow);
418  hist->Draw("hist");
419  delete hist;
420  }
421 
422  cmfactor.SaveAs(Form("cdcdedx_1dcell_mergefactor%s.pdf", m_suffix.data()));
423  cmfactor.SaveAs(Form("cdcdedx_1dcell_mergefactor%s.root", m_suffix.data()));
424 }
425 
426 //--------------------------------------------------
427 void CDCDedx1DCellAlgorithm::plotdedxHist(std::vector<TH1D*> hdedxhit[2])
428 {
429 
430  TCanvas ctmp("tmp", "tmp", 1200, 1200);
431  ctmp.Divide(4, 4);
432  std::stringstream psname;
433 
434  psname << Form("cdcdedx_1dcell_dedxhit%s.pdf[", m_suffix.data());
435  ctmp.Print(psname.str().c_str());
436  psname.str("");
437  psname << Form("cdcdedx_1dcell_dedxhit%s.pdf", m_suffix.data());
438 
439  for (int il = 0; il < 2; il++) {
440 
441  for (int jea = 0; jea < m_eaBinLocal[il]; jea++) {
442 
443  int minbin = std::stoi(hdedxhit[il][jea]->GetXaxis()->GetTitle());
444  int maxbin = std::stoi(hdedxhit[il][jea]->GetYaxis()->GetTitle());
445 
446  ctmp.cd(jea % 16 + 1);
447  hdedxhit[il][jea]->SetFillColor(4 + il);
448 
449  hdedxhit[il][jea]->SetTitle(Form("%s;dedxhit;entries", hdedxhit[il][jea]->GetTitle()));
450  hdedxhit[il][jea]->DrawClone("hist");
451  TH1D* htempC = (TH1D*)hdedxhit[il][jea]->Clone(Form("%sc2", hdedxhit[il][jea]->GetName()));
452  htempC->GetXaxis()->SetRange(minbin, maxbin);
453  htempC->SetFillColor(kGray);
454  htempC->DrawClone("same hist");
455 
456  if (jea % 16 == 15 || (jea == m_eaBinLocal[il] - 1)) {
457  ctmp.Print(psname.str().c_str());
458  gPad->Clear("D");
459  ctmp.Clear("D");
460  }
461  delete htempC;
462  }
463  }
464  psname.str("");
465  psname << Form("cdcdedx_1dcell_dedxhit%s.pdf]", m_suffix.data());
466  ctmp.Print(psname.str().c_str());
467 }
468 
469 //--------------------------------------------------
471 {
472 
473  TCanvas cdedxlayer("layerdedxhit", "Inner and Outer Layer dedxhit dist", 900, 400);
474  cdedxlayer.Divide(2, 1);
475 
476  for (int il = 0; il < 2; il++) {
477  int minlay = 0, maxlay = 0;
478  if (isFixTrunc) {
479  minlay = std::stoi(hdedxlay[il]->GetXaxis()->GetTitle());
480  maxlay = std::stoi(hdedxlay[il]->GetYaxis()->GetTitle());
481  double lowedge = hdedxlay[il]->GetXaxis()->GetBinLowEdge(minlay);
482  double upedge = hdedxlay[il]->GetXaxis()->GetBinUpEdge(maxlay);
483  hdedxlay[il]->SetTitle(Form("%s, trunc #rightarrow: %0.02f - %0.02f;dedxhit;entries", hdedxlay[il]->GetTitle(), lowedge, upedge));
484  }
485 
486  cdedxlayer.cd(il + 1);
487  hdedxlay[il]->SetFillColor(kYellow);
488  hdedxlay[il]->Draw("histo");
489 
490  if (isFixTrunc) {
491  TH1D* hdedxlayC = (TH1D*)hdedxlay[il]->Clone(Form("hdedxlayC%d", il));
492  hdedxlayC->GetXaxis()->SetRange(minlay, maxlay);
493  hdedxlayC->SetFillColor(kAzure + 1);
494  hdedxlayC->Draw("same histo");
495  }
496  }
497 
498  cdedxlayer.SaveAs(Form("cdcdedx_1dcell_dedxlay%s.pdf", m_suffix.data()));
499  cdedxlayer.SaveAs(Form("cdcdedx_1dcell_dedxlay%s.root", m_suffix.data()));
500 }
501 
502 //--------------------------------------------------
503 void CDCDedx1DCellAlgorithm::plotQaPars(TH1D* hentalay[2], TH2D* hptcosth)
504 {
505 
506  TCanvas ceadist("ceadist", "Enta distributions", 800, 400);
507  ceadist.Divide(2, 1);
508 
509  for (int il = 0; il < 2; il++) {
510 
511  ceadist.cd(il + 1);
512  gPad->SetLogy();
513  hentalay[il]->SetFillColor(kYellow);
514  hentalay[il]->Draw("hist");
515  }
516 
517  TCanvas cptcos("cptcos", "pt vs costh dist.", 400, 400);
518  cptcos.cd();
519  hptcosth->Draw("colz");
520 
521  cptcos.SaveAs(Form("cdcdedx_ptcosth_%s.pdf", m_suffix.data()));
522  ceadist.SaveAs(Form("cdcdedx_1dcell_enta%s.pdf", m_suffix.data()));
523  ceadist.SaveAs(Form("cdcdedx_1dcell_enta%s.root", m_suffix.data()));
524 }
525 
526 //--------------------------------------------------
527 void CDCDedx1DCellAlgorithm::plotRelConst(std::vector<double>tempconst, std::vector<double>layerconst, int il)
528 {
529 
530  TH1D* hconst, *hconstvar;
531 
532  Double_t* nvarBins;
533  nvarBins = &m_binValue[il][0];
534 
535  hconst = new TH1D(Form("hconst%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
536  std::string title = Form("calibration const dist: %s: (%s); entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
537  hconst->SetTitle(Form("%s", title.data()));
538 
539  hconstvar = new TH1D(Form("hconstvar%s", m_label[il].data()), "", m_eaBinLocal[il], nvarBins);
540  title = Form("calibration const dist (var bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
541  hconstvar->SetTitle(Form("%s", title.data()));
542 
543  if (isVarBins) {
544  for (int iea = 0; iea < m_eaBinLocal[il]; iea++)
545  hconstvar->SetBinContent(iea + 1, tempconst.at(iea));
546  }
547 
548  for (int jea = 0; jea < m_eaBin; jea++) hconst->SetBinContent(jea + 1, layerconst.at(jea));
549 
550  gStyle->SetOptStat("ne");
551  TCanvas cconst("cconst", "Calirbation Constants", 800, 400);
552  if (isVarBins) {
553  cconst.Divide(2, 1);
554  cconst.SetWindowSize(1000, 800);
555  }
556 
557  cconst.cd(1);
558  hconst->SetFillColor(kYellow);
559  hconst->Draw("histo");
560  if (isVarBins) {
561  cconst.cd(2);
562  hconstvar->SetFillColor(kBlue);
563  hconstvar->Draw("hist");
564  }
565  cconst.SaveAs(Form("cdcdedx_1dcell_relconst%s_%s.pdf", m_label[il].data(), m_suffix.data()));
566  cconst.SaveAs(Form("cdcdedx_1dcell_relconst%s_%s.root", m_label[il].data(), m_suffix.data()));
567 
568  delete hconst;
569  delete hconstvar;
570 }
571 
572 //--------------------------------------------------
574 {
575 
576  //Draw New/Old final constants
577  TH1D* hnewconst[2], *holdconst[2];
578  double min[2], max[2];
579 
580  for (unsigned int il = 0; il < 2; il++) {
581  unsigned int nbins = m_DBOneDCell->getNBins(il);
582 
583  std::string title = Form("final calibration const dist (%s): %s; entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
584  hnewconst[il] = new TH1D(Form("hnewconst_%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
585  hnewconst[il]->SetTitle(Form("%s", title.data()));
586 
587  title = Form("old calibration const dist (%s): %s; entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
588  holdconst[il] = new TH1D(Form("holdconst_%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
589  holdconst[il]->SetTitle(Form("%s", title.data()));
590 
591  for (unsigned int iea = 0; iea < nbins; iea++) {
592  double prev = m_DBOneDCell->getMean(8 * il + 1, iea);
593  holdconst[il]->SetBinContent(iea + 1, prev);
594  hnewconst[il]->SetBinContent(iea + 1, m_onedcors[il][iea]);
595  }
596  min[il] = hnewconst[il]->GetMinimum();
597  max[il] = hnewconst[il]->GetMaximum();
598  }
599 
600  //Ploting final constants
601  if (max[1] < max[0])max[1] = max[0];
602  if (min[1] > min[0])min[1] = min[0];
603 
604  gStyle->SetOptStat("ne");
605  TCanvas cfconst("cfconst", "Final calirbation constants", 800, 400);
606  cfconst.Divide(2, 1);
607 
608  for (int il = 0; il < 2; il++) {
609  cfconst.cd(il + 1);
610  hnewconst[il]->GetYaxis()->SetRangeUser(min[1] * 0.95, max[1] * 1.05);
611  hnewconst[il]->SetLineColor(kBlack);
612  hnewconst[il]->Draw("histo");
613  holdconst[il]->SetLineColor(kRed);
614  holdconst[il]->Draw("histo same");
615 
616  auto legend = new TLegend(0.4, 0.75, 0.56, 0.85);
617  legend->AddEntry(holdconst[il], "Old", "lep");
618  legend->AddEntry(hnewconst[il], "New", "lep");
619  legend->Draw();
620  }
621 
622  cfconst.SaveAs(Form("cdcdedx_1dcell_fconsts%s.pdf", m_suffix.data()));
623  cfconst.SaveAs(Form("cdcdedx_1dcell_fconsts%s.root", m_suffix.data()));
624 
625  for (int il = 0; il < 2; il++) {
626  delete hnewconst[il];
627  delete holdconst[il];
628  }
629 }
630 
631 //------------------------------------
633 {
634 
635  TCanvas cstats("cstats", "cstats", 1000, 500);
636  cstats.SetBatch(kTRUE);
637  cstats.Divide(2, 1);
638 
639  cstats.cd(1);
640  auto hestats = getObjectPtr<TH1I>("hestats");
641  if (hestats) {
642  hestats->SetName(Form("hestats_%s", m_suffix.data()));
643  hestats->SetStats(0);
644  hestats->DrawCopy("");
645  }
646 
647  cstats.cd(2);
648  auto htstats = getObjectPtr<TH1I>("htstats");
649  if (htstats) {
650  hestats->SetName(Form("htstats_%s", m_suffix.data()));
651  htstats->DrawCopy("");
652  hestats->SetStats(0);
653  }
654  cstats.Print(Form("cdcdedx_1dcell_stats_%s.pdf", m_suffix.data()));
655 }
CDCDedx1DCellAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
std::string m_label[2]
add inner/outer layer label
double m_eaMax
upper edge of enta angle
double m_truncMax
uppr thershold on truncation
int m_binSplit
multiply nbins by this factor in full range
std::array< std::vector< int >, 2 > m_binIndex
symm/Var bin numebrs
double m_truncMin
lower thershold on truncation
void plotMergeFactor(std::map< int, std::vector< double >> bounds, const std::array< int, 2 > nDev, std::map< int, std::vector< int >> steps)
funtion to plot merging factor
double m_adjustFac
faactor with that one what to adjust baseline
void getTruncatedBins(TH1D *hist, int &binlow, int &binhigh)
function to get bins of trunction from histogram
void CreateBinMapping()
class function to create vectors for bin mappping (Var->symm)
double m_chargeType
charge type for baseline adj
void getExpRunInfo()
function to get extract calibration run/exp
DBObjPtr< CDCDedx1DCell > m_DBOneDCell
One cell correction DB object.
double m_cosMax
a limit on cos theta
bool isPrintLog
print more debug information
std::array< std::vector< double >, 2 > m_binValue
etna Var bin values
std::string m_suffix
add suffix to all plot name
double getTruncationMean(TH1D *hist, int binlow, int binhigh)
function to get truncated mean
double m_ptMax
a limit on transverse momentum
void plotConstants()
function to draw the old/new final constants
bool isFixTrunc
true = fix window for all out/inner layers
bool isVarBins
true: if variable bin size is requested
void plotdedxHist(std::vector< TH1D * > hdedxhit[2])
function to draw the dE/dx histrogram in enta bins
void defineHisto(std::vector< TH1D * > hdedxhit[2], TH1D *hdedxlay[2], TH1D *hentalay[2])
function to define histograms
double m_eaBW
binwdith of enta angle bin
bool isRotSymm
if rotation symmtery requested
std::string m_runExp
add suffix to all plot name
void plotEventStats()
function to draw the stats plots
int rotationalBin(int nbin, int ibin)
class funtion to set rotation symmetry
virtual EResult calibrate() override
1D cell algorithm
void plotQaPars(TH1D *hentalay[2], TH2D *hptcosth)
funtion to draw pt vs costh and entrance angle distribution for Inner/Outer layer
double m_dedxMax
upper edge of dedxhit
void createPayload()
funtion to generate final constants
bool isMakePlots
produce plots for status
void plotRelConst(std::vector< double >tempconst, std::vector< double >layerconst, int il)
function to draw symm/Var layer constant
std::vector< std::vector< double > > m_onedcors
final vectors of calibration
bool isMerge
print more debug information
double m_dedxMin
lower edge of dedxhit
void plotLayerDist(TH1D *hdedxL[2])
funtion to draw dedx dist.
double m_eaMin
lower edge of enta angle
dE/dx wire gain calibration constants
Definition: CDCDedx1DCell.h:27
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
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)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Abstract base class for different kinds of events.