Belle II Software  release-06-00-14
PerformanceEvaluationBaseClass.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 <tracking/modules/trackingPerformanceEvaluation/PerformanceEvaluationBaseClass.h>
10 
11 #include <framework/logging/Logger.h>
12 
13 #include <root/TAxis.h>
14 
15 using namespace Belle2;
16 
17 PerformanceEvaluationBaseClass:: PerformanceEvaluationBaseClass()
18 {
19 }
20 
21 PerformanceEvaluationBaseClass::~ PerformanceEvaluationBaseClass()
22 {
23 
24 }
25 
26 TH1F* PerformanceEvaluationBaseClass::createHistogram1D(const char* name, const char* title,
27  Int_t nbins, Double_t min, Double_t max,
28  const char* xtitle, TList* histoList)
29 {
30 
31  TH1F* h = new TH1F(name, title, nbins, min, max);
32 
33  h->GetXaxis()->SetTitle(xtitle);
34 
35  if (histoList)
36  histoList->Add(h);
37 
38  return h;
39 }
40 
41 TH1F* PerformanceEvaluationBaseClass::createHistogram1D(const char* name, const char* title,
42  Int_t nbins, Double_t* bins,
43  const char* xtitle, TList* histoList)
44 {
45 
46  TH1F* h = new TH1F(name, title, nbins, bins);
47 
48  h->GetXaxis()->SetTitle(xtitle);
49 
50  if (histoList)
51  histoList->Add(h);
52 
53  return h;
54 }
55 
56 TH2F* PerformanceEvaluationBaseClass::createHistogram2D(const char* name, const char* title,
57  Int_t nbinsX, Double_t minX, Double_t maxX,
58  const char* titleX,
59  Int_t nbinsY, Double_t minY, Double_t maxY,
60  const char* titleY, TList* histoList)
61 {
62 
63  TH2F* h = new TH2F(name, title, nbinsX, minX, maxX, nbinsY, minY, maxY);
64 
65  h->GetXaxis()->SetTitle(titleX);
66  h->GetYaxis()->SetTitle(titleY);
67 
68  if (histoList)
69  histoList->Add(h);
70 
71  return h;
72 }
73 
74 TH2F* PerformanceEvaluationBaseClass::createHistogram2D(const char* name, const char* title,
75  Int_t nbinsX, Double_t* binsX,
76  const char* titleX,
77  Int_t nbinsY, Double_t* binsY,
78  const char* titleY,
79  TList* histoList)
80 {
81 
82  TH2F* h = new TH2F(name, title, nbinsX, binsX, nbinsY, binsY);
83 
84  h->GetXaxis()->SetTitle(titleX);
85  h->GetYaxis()->SetTitle(titleY);
86 
87  if (histoList)
88  histoList->Add(h);
89 
90  return h;
91 }
92 
93 TH3F* PerformanceEvaluationBaseClass::createHistogram3D(const char* name, const char* title,
94  Int_t nbinsX, Double_t minX, Double_t maxX,
95  const char* titleX,
96  Int_t nbinsY, Double_t minY, Double_t maxY,
97  const char* titleY,
98  Int_t nbinsZ, Double_t minZ, Double_t maxZ,
99  const char* titleZ,
100  TList* histoList)
101 {
102 
103  TH3F* h = new TH3F(name, title, nbinsX, minX, maxX, nbinsY, minY, maxY, nbinsZ, minZ, maxZ);
104 
105  h->GetXaxis()->SetTitle(titleX);
106  h->GetYaxis()->SetTitle(titleY);
107  h->GetZaxis()->SetTitle(titleZ);
108 
109  if (histoList)
110  histoList->Add(h);
111 
112  return h;
113 }
114 
115 TH3F* PerformanceEvaluationBaseClass::createHistogram3D(const char* name, const char* title,
116  Int_t nbinsX, Double_t* binsX,
117  const char* titleX,
118  Int_t nbinsY, Double_t* binsY,
119  const char* titleY,
120  Int_t nbinsZ, Double_t* binsZ,
121  const char* titleZ,
122  TList* histoList)
123 {
124 
125  TH3F* h = new TH3F(name, title, nbinsX, binsX, nbinsY, binsY, nbinsZ, binsZ);
126 
127  h->GetXaxis()->SetTitle(titleX);
128  h->GetYaxis()->SetTitle(titleY);
129  h->GetZaxis()->SetTitle(titleZ);
130 
131  if (histoList)
132  histoList->Add(h);
133 
134  return h;
135 }
136 
137 TH1* PerformanceEvaluationBaseClass::duplicateHistogram(const char* newname, const char* newtitle,
138  TH1* h, TList* histoList)
139 {
140 
141  TH1F* h1 = dynamic_cast<TH1F*>(h);
142  TH2F* h2 = dynamic_cast<TH2F*>(h);
143  TH3F* h3 = dynamic_cast<TH3F*>(h);
144 
145  TH1* newh = 0;
146 
147  if (h1)
148  newh = new TH1F(*h1);
149  if (h2)
150  newh = new TH2F(*h2);
151  if (h3)
152  newh = new TH3F(*h3);
153 
154  newh->SetName(newname);
155  newh->SetTitle(newtitle);
156 
157  if (histoList)
158  histoList->Add(newh);
159 
160 
161  return newh;
162 }
163 
164 TH1F* PerformanceEvaluationBaseClass::createHistogramsRatio(const char* name, const char* title,
165  TH1* hNum, TH1* hDen, bool isEffPlot,
166  int axisRef)
167 {
168 
169  TH1F* h1den = dynamic_cast<TH1F*>(hDen);
170  TH1F* h1num = dynamic_cast<TH1F*>(hNum);
171  TH2F* h2den = dynamic_cast<TH2F*>(hDen);
172  TH2F* h2num = dynamic_cast<TH2F*>(hNum);
173  TH3F* h3den = dynamic_cast<TH3F*>(hDen);
174  TH3F* h3num = dynamic_cast<TH3F*>(hNum);
175 
176  TH1* hden = 0;
177  TH1* hnum = 0;
178 
179  if (h1den) {
180  hden = new TH1F(*h1den);
181  hnum = new TH1F(*h1num);
182  }
183  if (h2den) {
184  hden = new TH2F(*h2den);
185  hnum = new TH2F(*h2num);
186  }
187  if (h3den) {
188  hden = new TH3F(*h3den);
189  hnum = new TH3F(*h3num);
190  }
191 
192  TAxis* the_axis;
193  TAxis* the_other1;
194  TAxis* the_other2;
195 
196  if (axisRef == 0) {
197  the_axis = hden->GetXaxis();
198  the_other1 = hden->GetYaxis();
199  the_other2 = hden->GetZaxis();
200  } else if (axisRef == 1) {
201  the_axis = hden->GetYaxis();
202  the_other1 = hden->GetXaxis();
203  the_other2 = hden->GetZaxis();
204  } else if (axisRef == 2) {
205  the_axis = hden->GetZaxis();
206  the_other1 = hden->GetXaxis();
207  the_other2 = hden->GetYaxis();
208  } else
209  return nullptr;
210 
211 
212  TH1F* h;
213  if (the_axis->GetXbins()->GetSize())
214  h = new TH1F(name, title, the_axis->GetNbins(), (the_axis->GetXbins())->GetArray());
215  else
216  h = new TH1F(name, title, the_axis->GetNbins(), the_axis->GetXmin(), the_axis->GetXmax());
217  h->GetXaxis()->SetTitle(the_axis->GetTitle());
218 
219  h->GetYaxis()->SetRangeUser(0.00001, 1);
220 
221  Int_t bin = 0;
222  Int_t nBins = 0;
223 
224  for (int the_bin = 1; the_bin < the_axis->GetNbins() + 1; the_bin++) {
225 
226  double num = 0;
227  double den = 0;
228 
229  for (int other1_bin = 1; other1_bin < the_other1->GetNbins() + 1; other1_bin++)
230  for (int other2_bin = 1; other2_bin < the_other2->GetNbins() + 1; other2_bin++) {
231 
232  if (axisRef == 0) bin = hden->GetBin(the_bin, other1_bin, other2_bin);
233  else if (axisRef == 1) bin = hden->GetBin(other1_bin, the_bin, other2_bin);
234  else if (axisRef == 2) bin = hden->GetBin(other1_bin, other2_bin, the_bin);
235 
236  if (hden->IsBinUnderflow(bin))
237  B2INFO(" bin = " << bin << "(" << the_bin << "," << other1_bin << "," << other2_bin << "), UNDERFLOW");
238  if (hden->IsBinOverflow(bin))
239  B2INFO(" bin = " << bin << "(" << the_bin << "," << other1_bin << "," << other2_bin << "), OVERFLOW");
240 
241  num += hnum->GetBinContent(bin);
242  den += hden->GetBinContent(bin);
243 
244  nBins++;
245  }
246 
247  double eff = 0;
248  double err = 0;
249 
250  if (den > 0) {
251  eff = (double)num / den;
252  err = sqrt(eff * (1 - eff)) / sqrt(den);
253  }
254 
255  if (isEffPlot) {
256  h->SetBinContent(the_bin, eff);
257  h->SetBinError(the_bin, err);
258  } else {
259  h->SetBinContent(the_bin, 1 - eff);
260  h->SetBinError(the_bin, err);
261  }
262  }
263 
264  return h;
265 
266 }
267 
268 
269 void PerformanceEvaluationBaseClass::addInefficiencyPlots(TList* histoList, TH3F* h3_xPerMCParticle, TH3F* h3_MCParticle)
270 {
271 
272  if ((h3_xPerMCParticle == nullptr) || (h3_MCParticle == nullptr))
273  return;
274 
275  //normalized to MCParticles
276  TH1F* h_ineff_pt = createHistogramsRatio("hineffpt", "inefficiency VS pt, normalized to MCParticles", h3_xPerMCParticle,
277  h3_MCParticle, false, 0);
278  histoList->Add(h_ineff_pt);
279 
280  TH1F* h_ineff_theta = createHistogramsRatio("hinefftheta", "inefficiency VS #theta, normalized to MCParticles",
281  h3_xPerMCParticle, h3_MCParticle, false, 1);
282  histoList->Add(h_ineff_theta);
283 
284  TH1F* h_ineff_phi = createHistogramsRatio("hineffphi", "inefficiency VS #phi, normalized to MCParticles", h3_xPerMCParticle,
285  h3_MCParticle, false, 2);
286  histoList->Add(h_ineff_phi);
287 
288 }
289 
290 void PerformanceEvaluationBaseClass::addEfficiencyPlots(TList* histoList, TH3F* h3_xPerMCParticle, TH3F* h3_MCParticle)
291 {
292  if ((h3_xPerMCParticle == nullptr) || (h3_MCParticle == nullptr))
293  return;
294 
295  //normalized to MCParticles
296  TH1F* h_eff_pt = createHistogramsRatio("heffpt", "efficiency VS pt, normalized to MCParticles", h3_xPerMCParticle,
297  h3_MCParticle, true, 0);
298  histoList->Add(h_eff_pt);
299 
300  TH1F* h_eff_theta = createHistogramsRatio("hefftheta", "efficiency VS #theta, normalized to MCParticles", h3_xPerMCParticle,
301  h3_MCParticle, true, 1);
302  histoList->Add(h_eff_theta);
303 
304  TH1F* h_eff_phi = createHistogramsRatio("heffphi", "efficiency VS #phi, normalized to MCParticles", h3_xPerMCParticle,
305  h3_MCParticle, true, 2);
306  histoList->Add(h_eff_phi);
307 
308 }
309 
310 
311 
312 void PerformanceEvaluationBaseClass::addPurityPlots(TList* histoList, TH3F* h3_MCParticlesPerX, TH3F* h3_X)
313 {
314  if ((h3_X == nullptr) || (h3_MCParticlesPerX == nullptr))
315  return;
316 
317  //purity histograms
318  TH1F* h_pur_pt = createHistogramsRatio("hpurpt", "purity VS pt", h3_MCParticlesPerX, h3_X, true, 0);
319  histoList->Add(h_pur_pt);
320 
321  TH1F* h_pur_theta = createHistogramsRatio("hpurtheta", "purity VS #theta", h3_MCParticlesPerX, h3_X, true, 1);
322  histoList->Add(h_pur_theta);
323 
324  TH1F* h_pur_phi = createHistogramsRatio("hpurphi", "purity VS #phi", h3_MCParticlesPerX, h3_X, true, 2);
325  histoList->Add(h_pur_phi);
326 
327 }
328 
329 TH1F* PerformanceEvaluationBaseClass::effPlot1D(TH1F* h1_den, TH1F* h1_num, const char* name, const char* title,
330  bool geo_accettance, TList* histoList)
331 {
332 
333  std::string total;
334  std::string trueTitle;
335 
336  if (geo_accettance == false) {
337  std::string name1 = "_noGeoAcc";
338  total = std::string(name) + name1;
339  trueTitle = std::string(title) + name1;
340  } else {
341  std::string name2 = "_withGeoAcc";
342  total = std::string(name) + name2;
343  trueTitle = std::string(title) + name2;
344  }
345 
346  TH1F* h = (TH1F*)duplicateHistogram(total.c_str(), trueTitle.c_str(), h1_den, histoList);
347  h->GetYaxis()->SetRangeUser(0., 1);
348 
349  for (int bin = 0; bin < h->GetXaxis()->GetNbins(); bin++) {
350  float num = h1_num->GetBinContent(bin + 1);
351  float den = h1_den->GetBinContent(bin + 1);
352  double eff = 0.;
353  double err = 0.;
354 
355  if (den > 0) {
356  eff = (double)num / den;
357  err = sqrt(eff * (1 - eff)) / sqrt(den);
358  }
359  h->SetBinContent(bin + 1, eff);
360  h->SetBinError(bin + 1, err);
361  }
362 
363 
364  if (histoList)
365  histoList->Add(h);
366 
367  return h;
368 }
369 
370 TH1F* PerformanceEvaluationBaseClass::effPlot1D(TH1F* h1_MC, TH1F* h1_RecoTrack, TH1F* h1_Track, const char* name,
371  const char* title, TList* histoList)
372 {
373  if (h1_Track == nullptr) B2INFO("h_Track missing");
374 
375  std::string name1 = "_noGeoAcc";
376  std::string name2 = "_withGeoAcc";
377 
378  std::string total1 = std::string(name) + name1;
379  std::string total2 = std::string(name) + name2;
380 
381  std::string title1 = std::string(title) + name1;
382  std::string title2 = std::string(title) + name2;
383 
384  TH1F* h[2];
385 
386  h[0] = (TH1F*)duplicateHistogram(total2.c_str(), title2.c_str(), h1_RecoTrack, histoList);
387  h[0]->GetYaxis()->SetRangeUser(0., 1);
388 
389  for (int bin = 0; bin < h[0]->GetXaxis()->GetNbins(); bin++) {
390  float num = h1_Track->GetBinContent(bin + 1);
391  float den = h1_RecoTrack->GetBinContent(bin + 1);
392  double eff = 0.;
393  double err = 0.;
394 
395  if (den > 0) {
396  eff = (double)num / den;
397  err = sqrt(eff * (1 - eff)) / sqrt(den);
398  }
399  h[0]->SetBinContent(bin + 1, eff);
400  h[0]->SetBinError(bin + 1, err);
401  }
402 
403  h[1] = (TH1F*)duplicateHistogram(total1.c_str(), title1.c_str(), h1_MC, histoList);
404  h[1]->GetYaxis()->SetRangeUser(0., 1);
405 
406  for (int bin = 0; bin < h[1]->GetXaxis()->GetNbins(); bin++) {
407  float num = h1_Track->GetBinContent(bin + 1);
408  float den = h1_MC->GetBinContent(bin + 1);
409  double eff = 0.;
410  double err = 0.;
411 
412  if (den > 0) {
413  eff = (double)num / den;
414  err = sqrt(eff * (1 - eff)) / sqrt(den);
415  }
416  h[1]->SetBinContent(bin + 1, eff);
417  h[1]->SetBinError(bin + 1, err);
418  }
419 
420  if (histoList) {
421  histoList->Add(h[0]);
422  histoList->Add(h[1]);
423  }
424 
425  return *h;
426 }
427 
428 TH2F* PerformanceEvaluationBaseClass::effPlot2D(TH2F* h2_den, TH2F* h2_num,
429  const char* name, const char* title, bool geo_accettance, TList* histoList)
430 {
431  std::string err_char = "_error_";
432  std::string addTitle = "Errors, ";
433 
434  std::string total;
435  std::string error;
436  std::string trueTitle;
437 
438  std::string titleErr = addTitle + std::string(title);
439 
440  if (geo_accettance == false) {
441  std::string name1 = "_noGeoAcc";
442  total = std::string(name) + name1;
443  trueTitle = std::string(title) + name1;
444  error = std::string(name) + err_char + std::string(name1);
445  } else {
446  std::string name2 = "_withGeoAcc";
447  total = std::string(name) + name2;
448  trueTitle = std::string(title) + name2;
449  error = std::string(name) + err_char + std::string(name2);
450  }
451 
452  TH2F* h2[2];
453  h2[0] = (TH2F*)duplicateHistogram(total.c_str(), trueTitle.c_str(), h2_den, histoList);
454  h2[1] = (TH2F*)duplicateHistogram(error.c_str(), titleErr.c_str(), h2_den, histoList);
455 
456  for (int binX = 0; binX < h2[0]->GetXaxis()->GetNbins(); binX++) {
457  for (int binY = 0; binY < h2[0]->GetYaxis()->GetNbins(); binY++) {
458  float num = h2_num->GetBinContent(binX + 1, binY + 1);
459  float den = h2_den->GetBinContent(binX + 1, binY + 1);
460  double eff = 0.;
461  double err = 0.;
462 
463  if (den > 0) {
464  eff = (double)num / den;
465  err = sqrt(eff * (1 - eff)) / sqrt(den);
466  }
467 
468  h2[0]->SetBinContent(binX + 1, binY + 1, eff);
469  h2[0]->SetBinError(binX + 1, binY + 1, err);
470  h2[1]->SetBinContent(binX + 1, binY + 1, err);
471  }
472  }
473 
474 
475  if (histoList) {
476  histoList->Add(h2[0]);
477  histoList->Add(h2[1]);
478  }
479 
480  return *h2;
481 
482 }
483 
484 TH2F* PerformanceEvaluationBaseClass::effPlot2D(TH2F* h2_MC, TH2F* h2_RecoTrack, TH2F* h2_Track, const char* name,
485  const char* title, TList* histoList)
486 {
487  if (h2_Track == nullptr) B2INFO("h_Track missing");
488 
489  std::string name1 = "_noGeoAcc";
490  std::string name2 = "_withGeoAcc";
491  std::string err_char = "_error_";
492  std::string addTitle = "Errors, ";
493 
494  std::string total1 = std::string(name) + name1;
495  std::string total2 = std::string(name) + name2;
496 
497  std::string title1 = std::string(title) + name1;
498  std::string title2 = std::string(title) + name2;
499 
500  std::string error1 = std::string(name) + err_char + name1;
501  std::string error2 = std::string(name) + err_char + name2;
502  std::string titleErr = addTitle + std::string(title);
503 
504  TH2F* h2[4];
505 
506  h2[0] = (TH2F*)duplicateHistogram(total2.c_str(), title2.c_str(), h2_RecoTrack, histoList);
507  h2[1] = (TH2F*)duplicateHistogram(error2.c_str(), titleErr.c_str(), h2_RecoTrack, histoList);
508 
509  for (int binX = 0; binX < h2[0]->GetXaxis()->GetNbins(); binX++) {
510  for (int binY = 0; binY < h2[0]->GetYaxis()->GetNbins(); binY++) {
511  float num = h2_Track->GetBinContent(binX + 1, binY + 1);
512  float den = h2_RecoTrack->GetBinContent(binX + 1, binY + 1);
513  double eff = 0.;
514  double err = 0.;
515 
516  if (den > 0) {
517  eff = num / den;
518  err = sqrt(eff * (1 - eff)) / sqrt(den);
519  }
520 
521  h2[0]->SetBinContent(binX + 1, binY + 1, eff);
522  h2[0]->SetBinError(binX + 1, binY + 1, err);
523  h2[1]->SetBinContent(binX + 1, binY + 1, err);
524  }
525  }
526 
527  h2[2] = (TH2F*)duplicateHistogram(total1.c_str(), title1.c_str(), h2_MC, histoList);
528  h2[3] = (TH2F*)duplicateHistogram(error1.c_str(), titleErr.c_str(), h2_MC, histoList);
529 
530  for (int binX = 0; binX < h2[2]->GetXaxis()->GetNbins(); binX++) {
531  for (int binY = 0; binY < h2[2]->GetYaxis()->GetNbins(); binY++) {
532  float num = h2_Track->GetBinContent(binX + 1, binY + 1);
533  float den = h2_MC->GetBinContent(binX + 1, binY + 1);
534  double eff = 0.;
535  double err = 0.;
536 
537  if (den > 0) {
538  eff = num / den;
539  err = sqrt(eff * (1 - eff)) / sqrt(den);
540  }
541 
542  h2[2]->SetBinContent(binX + 1, binY + 1, eff);
543  h2[2]->SetBinError(binX + 1, binY + 1, err);
544  h2[3]->SetBinContent(binX + 1, binY + 1, err);
545  }
546  }
547 
548  if (histoList) {
549  histoList->Add(h2[0]);
550  histoList->Add(h2[1]);
551  histoList->Add(h2[2]);
552  histoList->Add(h2[3]);
553  }
554 
555  return *h2;
556 }
557 
558 TH1F* PerformanceEvaluationBaseClass::geoAcc1D(TH1F* h1_den, TH1F* h1_num, const char* name, const char* title, TList* histoList)
559 {
560  TH1F* h = (TH1F*)duplicateHistogram(name, title, h1_den, histoList);
561  h->GetYaxis()->SetRangeUser(0., 1);
562 
563  for (int bin = 0; bin < h->GetXaxis()->GetNbins(); bin++) {
564  float num = h1_num->GetBinContent(bin + 1);
565  float den = h1_den->GetBinContent(bin + 1);
566  double eff = 0.;
567  double err = 0.;
568 
569  if (den > 0) {
570  eff = (double)num / den;
571  err = sqrt(eff * (1 - eff)) / sqrt(den);
572  }
573  h->SetBinContent(bin + 1, eff);
574  h->SetBinError(bin + 1, err);
575  }
576 
577  if (histoList)
578  histoList->Add(h);
579 
580  return h;
581 }
582 
583 TH2F* PerformanceEvaluationBaseClass::geoAcc2D(TH2F* h2_den, TH2F* h2_num,
584  const char* name, const char* title, TList* histoList)
585 {
586  std::string err_char = "_err";
587  std::string addTitle = "Errors, ";
588 
589  std::string error = std::string(name) + err_char;
590  std::string titleErr = addTitle + std::string(title);
591 
592  TH2F* h2[2];
593  h2[0] = (TH2F*)duplicateHistogram(name, title, h2_den, histoList);
594  h2[1] = (TH2F*)duplicateHistogram(error.c_str(), titleErr.c_str(), h2_den, histoList);
595 
596  for (int binX = 0; binX < h2[0]->GetXaxis()->GetNbins(); binX++) {
597  for (int binY = 0; binY < h2[0]->GetYaxis()->GetNbins(); binY++) {
598  float num = h2_num->GetBinContent(binX + 1, binY + 1);
599  float den = h2_den->GetBinContent(binX + 1, binY + 1);
600  double eff = 0.;
601  double err = 0.;
602 
603  if (den > 0) {
604  eff = (double)num / den;
605  err = sqrt(eff * (1 - eff)) / sqrt(den);
606  }
607  h2[0]->SetBinContent(binX + 1, binY + 1, eff);
608  h2[0]->SetBinError(binX + 1, binY + 1, err);
609  h2[1]->SetBinContent(binX + 1, binY + 1, err);
610  }
611  }
612 
613  if (histoList) {
614  histoList->Add(h2[0]);
615  histoList->Add(h2[1]);
616  }
617 
618  return *h2;
619 
620 }
621 
622 TH1F* PerformanceEvaluationBaseClass::V0FinderEff(TH1F* h1_dau0, TH1F* h1_dau1, TH1F* h1_Mother,
623  const char* name, const char* title, TList* histoList)
624 {
625 
626  TH1F* h = (TH1F*)duplicateHistogram(name, title, h1_Mother, histoList);
627  h->GetYaxis()->SetRangeUser(0., 1);
628 
629  for (int bin = 0; bin < h->GetXaxis()->GetNbins(); bin++) {
630  double dau0 = h1_dau0->GetBinContent(bin);
631  double dau1 = h1_dau1->GetBinContent(bin);
632  double Mother = h1_Mother->GetBinContent(bin);
633  double dau0Err = h1_dau0->GetBinError(bin);
634  double dau1Err = h1_dau1->GetBinError(bin);
635  double MotherErr = h1_Mother->GetBinError(bin);
636 
637  double binCont = 1. * Mother / dau0 / dau1;
638  double binErr = binCont * sqrt((dau0Err / dau0) * (dau0Err / dau0) + (dau1Err / dau1) * (dau1Err / dau1) * (MotherErr / Mother) *
639  (MotherErr / Mother));
640 
641  h->SetBinContent(bin, binCont);
642  h->SetBinError(bin, binErr);
643  }
644 
645  if (histoList)
646  histoList->Add(h);
647 
648  return h;
649 }
650 
TH1F * createHistogram1D(const char *name, const char *title, Int_t nbins, Double_t min, Double_t max, const char *xtitle, TList *histoList=nullptr)
Create a 1D histogram and add it to the TList of 1D-histograms.
TH2F * effPlot2D(TH2F *h2_den, TH2F *h2_num, const char *name, const char *title, bool geo_accettance, TList *histoList=nullptr)
Create a 2D efficiency histogram and add it to the TList of 2D-histograms.
TH1F * effPlot1D(TH1F *h1_den, TH1F *h1_num, const char *name, const char *title, bool geo_accettance, TList *histoList=nullptr)
Create a 1D efficiency histogram and add it to the TList of 1D-histograms.
TH1 * duplicateHistogram(const char *newname, const char *newtitle, TH1 *h, TList *histoList=nullptr)
Make a copy of a 1D histogram and add it to the TList of 1D-histograms.
TH1F * V0FinderEff(TH1F *h1_dau0, TH1F *h1_dau1, TH1F *h1_Mother, const char *name, const char *title, TList *histoList=nullptr)
Create a 1D efficiency histogram for V0 finding and add it to the TList of 1D-histograms.
void addInefficiencyPlots(TList *graphList=nullptr, TH3F *h3_xPerMCParticle=nullptr, TH3F *h3_MCParticle=nullptr)
Create pt-, theta- and phi-inefficiency 1D histograms and add them to the TList of 1D-histograms.
void addPurityPlots(TList *graphList=nullptr, TH3F *h3_xPerMCParticle=nullptr, TH3F *h3_MCParticle=nullptr)
Create pt-, theta- and phi-purity 1D histograms and add them to the TList of 1D-histograms.
void addEfficiencyPlots(TList *graphList=nullptr, TH3F *h3_xPerMCParticle=nullptr, TH3F *h3_MCParticle=nullptr)
Create pt-, theta- and phi-efficiency 1D histograms and add them to the TList of 1D-histograms.
TH1F * geoAcc1D(TH1F *h1_den, TH1F *h1_num, const char *name, const char *title, TList *histoList=nullptr)
Create a 1D efficiency histogram for geometric acceptance and add it to the TList of 1D-histograms.
TH3F * createHistogram3D(const char *name, const char *title, Int_t nbinsX, Double_t minX, Double_t maxX, const char *titleX, Int_t nbinsY, Double_t minY, Double_t maxY, const char *titleY, Int_t nbinsZ, Double_t minZ, Double_t maxZ, const char *titleZ, TList *histoList=nullptr)
Create a 3D histogram and add it to the TList of 3D-histograms.
TH2F * createHistogram2D(const char *name, const char *title, Int_t nbinsX, Double_t minX, Double_t maxX, const char *titleX, Int_t nbinsY, Double_t minY, Double_t maxY, const char *titleY, TList *histoList=nullptr)
Create a 2D histogram and add it to the TList of 2D-histograms.
TH2F * geoAcc2D(TH2F *h2_den, TH2F *h2_num, const char *name, const char *title, TList *histoList=nullptr)
Create a 2D efficiency histogram for geometric acceptance and add it to the TList of 2D-histograms.
TH1F * createHistogramsRatio(const char *name, const char *title, TH1 *hNum, TH1 *hDen, bool isEffPlot, int axisRef)
Make a new 1D histogram from the ratio of two others and add it to the TList of 1D-histograms.
Abstract base class for different kinds of events.