Belle II Software development
PerformanceEvaluationBaseClass Class Reference

This module takes the MCParticles, the genfit Tracks, the genfit TrackCand, and the MCTrackCands input and produce a root file containing various histograms showing the performance of the tracking package: fitter, pattern recognition algorithms. More...

#include <PerformanceEvaluationBaseClass.h>

Inheritance diagram for PerformanceEvaluationBaseClass:
EffPlotsModule FillTrackFitNtupleModule TrackingPerformanceEvaluationModule V0findingPerformanceEvaluationModule

Public Member Functions

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.
 
TH1F * createHistogram1D (const char *name, const char *title, Int_t nbins, Double_t *bins, const char *xtitle, TList *histoList=nullptr)
 Create a 1D histogram and add it to the TList of 1D-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 * createHistogram2D (const char *name, const char *title, Int_t nbinsX, Double_t *binsX, const char *titleX, Int_t nbinsY, Double_t *binsY, const char *titleY, TList *histoList=nullptr)
 Create a 2D histogram and add it to the TList of 2D-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.
 
TH3F * createHistogram3D (const char *name, const char *title, Int_t nbinsX, Double_t *binsX, const char *titleX, Int_t nbinsY, Double_t *binsY, const char *titleY, Int_t nbinsZ, Double_t *binsZ, const char *titleZ, TList *histoList=nullptr)
 Create a 3D histogram and add it to the TList of 3D-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 * 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.
 
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.
 
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.
 
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.
 
TH1F * effPlot1D (TH1F *h1_MC, TH1F *h1_RecoTrack, TH1F *h1_Track, const char *name, const char *title, TList *histoList=nullptr)
 Create a 1D efficiency 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.
 
TH2F * effPlot2D (TH2F *h2_MC, TH2F *h2_RecoTrack, TH2F *h2_Track, const char *name, const char *title, TList *histoList)
 Create a 2D efficiency histogram and add it to the TList of 2D-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.
 
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 * 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.
 

Public Attributes

TList * m_histoList = nullptr
 List of performance-evaluation histograms.
 
TList * m_histoList_multiplicity = nullptr
 List of multiplicity histograms.
 
TList * m_histoList_evtCharacterization = nullptr
 List of event-characterization histograms.
 
TList * m_histoList_trkQuality = nullptr
 List of track-quality histograms.
 
TList * m_histoList_firstHit = nullptr
 List of first-hit-position histograms.
 
TList * m_histoList_pr = nullptr
 List of pattern-recognition histograms.
 
TList * m_histoList_fit = nullptr
 List of track-fit histograms.
 
TList * m_histoList_efficiency = nullptr
 List of efficiency histograms.
 
TList * m_histoList_purity = nullptr
 List of purity histograms.
 
TList * m_histoList_others = nullptr
 List of other performance-evaluation histograms.
 
std::string m_rootFileName
 root file name
 
TFile * m_rootFilePtr = nullptr
 pointer at root file used for storing histograms
 

Detailed Description

This module takes the MCParticles, the genfit Tracks, the genfit TrackCand, and the MCTrackCands input and produce a root file containing various histograms showing the performance of the tracking package: fitter, pattern recognition algorithms.

Definition at line 30 of file PerformanceEvaluationBaseClass.h.

Constructor & Destructor Documentation

◆ PerformanceEvaluationBaseClass()

Definition at line 17 of file PerformanceEvaluationBaseClass.cc.

18{
19}

◆ ~PerformanceEvaluationBaseClass()

Definition at line 21 of file PerformanceEvaluationBaseClass.cc.

22{
23
24}

Member Function Documentation

◆ addEfficiencyPlots()

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.

Definition at line 293 of file PerformanceEvaluationBaseClass.cc.

294{
295 if ((h3_xPerMCParticle == nullptr) || (h3_MCParticle == nullptr))
296 return;
297
298 //normalized to MCParticles
299 TH1F* h_eff_pt = createHistogramsRatio("heffpt", "efficiency VS pt, normalized to MCParticles", h3_xPerMCParticle,
300 h3_MCParticle, true, 0);
301 histoList->Add(h_eff_pt);
302
303 TH1F* h_eff_theta = createHistogramsRatio("hefftheta", "efficiency VS #theta, normalized to MCParticles", h3_xPerMCParticle,
304 h3_MCParticle, true, 1);
305 histoList->Add(h_eff_theta);
306
307 TH1F* h_eff_phi = createHistogramsRatio("heffphi", "efficiency VS #phi, normalized to MCParticles", h3_xPerMCParticle,
308 h3_MCParticle, true, 2);
309 histoList->Add(h_eff_phi);
310
311}
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.

◆ addInefficiencyPlots()

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.

Definition at line 272 of file PerformanceEvaluationBaseClass.cc.

273{
274
275 if ((h3_xPerMCParticle == nullptr) || (h3_MCParticle == nullptr))
276 return;
277
278 //normalized to MCParticles
279 TH1F* h_ineff_pt = createHistogramsRatio("hineffpt", "inefficiency VS pt, normalized to MCParticles", h3_xPerMCParticle,
280 h3_MCParticle, false, 0);
281 histoList->Add(h_ineff_pt);
282
283 TH1F* h_ineff_theta = createHistogramsRatio("hinefftheta", "inefficiency VS #theta, normalized to MCParticles",
284 h3_xPerMCParticle, h3_MCParticle, false, 1);
285 histoList->Add(h_ineff_theta);
286
287 TH1F* h_ineff_phi = createHistogramsRatio("hineffphi", "inefficiency VS #phi, normalized to MCParticles", h3_xPerMCParticle,
288 h3_MCParticle, false, 2);
289 histoList->Add(h_ineff_phi);
290
291}

◆ addPurityPlots()

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.

Definition at line 315 of file PerformanceEvaluationBaseClass.cc.

316{
317 if ((h3_X == nullptr) || (h3_MCParticlesPerX == nullptr))
318 return;
319
320 //purity histograms
321 TH1F* h_pur_pt = createHistogramsRatio("hpurpt", "purity VS pt", h3_MCParticlesPerX, h3_X, true, 0);
322 histoList->Add(h_pur_pt);
323
324 TH1F* h_pur_theta = createHistogramsRatio("hpurtheta", "purity VS #theta", h3_MCParticlesPerX, h3_X, true, 1);
325 histoList->Add(h_pur_theta);
326
327 TH1F* h_pur_phi = createHistogramsRatio("hpurphi", "purity VS #phi", h3_MCParticlesPerX, h3_X, true, 2);
328 histoList->Add(h_pur_phi);
329
330}

◆ createHistogram1D() [1/2]

TH1F * createHistogram1D ( const char *  name,
const char *  title,
Int_t  nbins,
Double_t *  bins,
const char *  xtitle,
TList *  histoList = nullptr 
)

Create a 1D histogram and add it to the TList of 1D-histograms.

Definition at line 41 of file PerformanceEvaluationBaseClass.cc.

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}

◆ createHistogram1D() [2/2]

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.

Definition at line 26 of file PerformanceEvaluationBaseClass.cc.

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}

◆ createHistogram2D() [1/2]

TH2F * createHistogram2D ( const char *  name,
const char *  title,
Int_t  nbinsX,
Double_t *  binsX,
const char *  titleX,
Int_t  nbinsY,
Double_t *  binsY,
const char *  titleY,
TList *  histoList = nullptr 
)

Create a 2D histogram and add it to the TList of 2D-histograms.

Definition at line 74 of file PerformanceEvaluationBaseClass.cc.

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}

◆ createHistogram2D() [2/2]

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.

Create 2D histogram

Definition at line 56 of file PerformanceEvaluationBaseClass.cc.

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}

◆ createHistogram3D() [1/2]

TH3F * createHistogram3D ( const char *  name,
const char *  title,
Int_t  nbinsX,
Double_t *  binsX,
const char *  titleX,
Int_t  nbinsY,
Double_t *  binsY,
const char *  titleY,
Int_t  nbinsZ,
Double_t *  binsZ,
const char *  titleZ,
TList *  histoList = nullptr 
)

Create a 3D histogram and add it to the TList of 3D-histograms.

Definition at line 115 of file PerformanceEvaluationBaseClass.cc.

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}

◆ createHistogram3D() [2/2]

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.

Definition at line 93 of file PerformanceEvaluationBaseClass.cc.

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}

◆ createHistogramsRatio()

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.

Definition at line 170 of file PerformanceEvaluationBaseClass.cc.

173{
174
175 TH1F* h1den = dynamic_cast<TH1F*>(hDen);
176 TH1F* h1num = dynamic_cast<TH1F*>(hNum);
177 TH2F* h2den = dynamic_cast<TH2F*>(hDen);
178 TH2F* h2num = dynamic_cast<TH2F*>(hNum);
179 TH3F* h3den = dynamic_cast<TH3F*>(hDen);
180 TH3F* h3num = dynamic_cast<TH3F*>(hNum);
181
182 TH1* hden = 0;
183 TH1* hnum = 0;
184
185 if (h1den) {
186 hden = new TH1F(*h1den);
187 hnum = new TH1F(*h1num);
188 }
189 if (h2den) {
190 hden = new TH2F(*h2den);
191 hnum = new TH2F(*h2num);
192 }
193 if (h3den) {
194 hden = new TH3F(*h3den);
195 hnum = new TH3F(*h3num);
196 }
197
198 TAxis* the_axis;
199 TAxis* the_other1;
200 TAxis* the_other2;
201
202 if (axisRef == 0) {
203 the_axis = hden->GetXaxis();
204 the_other1 = hden->GetYaxis();
205 the_other2 = hden->GetZaxis();
206 } else if (axisRef == 1) {
207 the_axis = hden->GetYaxis();
208 the_other1 = hden->GetXaxis();
209 the_other2 = hden->GetZaxis();
210 } else if (axisRef == 2) {
211 the_axis = hden->GetZaxis();
212 the_other1 = hden->GetXaxis();
213 the_other2 = hden->GetYaxis();
214 } else
215 return nullptr;
216
217
218 TH1F* h;
219 if (the_axis->GetXbins()->GetSize())
220 h = new TH1F(name, title, the_axis->GetNbins(), (the_axis->GetXbins())->GetArray());
221 else
222 h = new TH1F(name, title, the_axis->GetNbins(), the_axis->GetXmin(), the_axis->GetXmax());
223 h->GetXaxis()->SetTitle(the_axis->GetTitle());
224
225 h->GetYaxis()->SetRangeUser(0.00001, 1);
226
227 Int_t bin = 0;
228
229 for (int the_bin = 1; the_bin < the_axis->GetNbins() + 1; the_bin++) {
230
231 double num = 0;
232 double den = 0;
233
234 for (int other1_bin = 1; other1_bin < the_other1->GetNbins() + 1; other1_bin++)
235 for (int other2_bin = 1; other2_bin < the_other2->GetNbins() + 1; other2_bin++) {
236
237 if (axisRef == 0) bin = hden->GetBin(the_bin, other1_bin, other2_bin);
238 else if (axisRef == 1) bin = hden->GetBin(other1_bin, the_bin, other2_bin);
239 else if (axisRef == 2) bin = hden->GetBin(other1_bin, other2_bin, the_bin);
240
241 if (hden->IsBinUnderflow(bin))
242 B2INFO(" bin = " << bin << "(" << the_bin << "," << other1_bin << "," << other2_bin << "), UNDERFLOW");
243 if (hden->IsBinOverflow(bin))
244 B2INFO(" bin = " << bin << "(" << the_bin << "," << other1_bin << "," << other2_bin << "), OVERFLOW");
245
246 num += hnum->GetBinContent(bin);
247 den += hden->GetBinContent(bin);
248 }
249
250 double eff = 0;
251 double err = 0;
252
253 if (den > 0) {
254 eff = (double)num / den;
255 err = sqrt(eff * (1 - eff)) / sqrt(den);
256 }
257
258 if (isEffPlot) {
259 h->SetBinContent(the_bin, eff);
260 h->SetBinError(the_bin, err);
261 } else {
262 h->SetBinContent(the_bin, 1 - eff);
263 h->SetBinError(the_bin, err);
264 }
265 }
266
267 return h;
268
269}
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ duplicateHistogram()

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.

Definition at line 137 of file PerformanceEvaluationBaseClass.cc.

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 = nullptr;
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 if (newh == nullptr) {
155 B2ERROR("In function duplicateHistogram: newh is a nullptr. This shouldn't happen."\
156 "Don't continue creation of duplicate histogram in this case and return nullptr.");
157 return nullptr;
158 }
159
160 newh->SetName(newname);
161 newh->SetTitle(newtitle);
162
163 if (histoList)
164 histoList->Add(newh);
165
166
167 return newh;
168}

◆ effPlot1D() [1/2]

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.

Definition at line 332 of file PerformanceEvaluationBaseClass.cc.

334{
335 if (h1_den == nullptr or h1_num == nullptr) {
336 B2ERROR("One of the input histograms for function effPlot1D is a nullptr, "\
337 "can't create new histograms from this. Returning a nullptr.");
338 return nullptr;
339 }
340
341 std::string total;
342 std::string trueTitle;
343
344 if (geo_accettance == false) {
345 std::string name1 = "_noGeoAcc";
346 total = std::string(name) + name1;
347 trueTitle = std::string(title) + name1;
348 } else {
349 std::string name2 = "_withGeoAcc";
350 total = std::string(name) + name2;
351 trueTitle = std::string(title) + name2;
352 }
353
354 TH1F* h = (TH1F*)duplicateHistogram(total.c_str(), trueTitle.c_str(), h1_den, histoList);
355 h->GetYaxis()->SetRangeUser(0., 1);
356
357 for (int bin = 0; bin < h->GetXaxis()->GetNbins(); bin++) {
358 float num = h1_num->GetBinContent(bin + 1);
359 float den = h1_den->GetBinContent(bin + 1);
360 double eff = 0.;
361 double err = 0.;
362
363 if (den > 0) {
364 eff = (double)num / den;
365 err = sqrt(eff * (1 - eff)) / sqrt(den);
366 }
367 h->SetBinContent(bin + 1, eff);
368 h->SetBinError(bin + 1, err);
369 }
370
371
372 if (histoList)
373 histoList->Add(h);
374
375 return h;
376}
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.

◆ effPlot1D() [2/2]

TH1F * effPlot1D ( TH1F *  h1_MC,
TH1F *  h1_RecoTrack,
TH1F *  h1_Track,
const char *  name,
const char *  title,
TList *  histoList = nullptr 
)

Create a 1D efficiency histogram and add it to the TList of 1D-histograms.

Definition at line 378 of file PerformanceEvaluationBaseClass.cc.

380{
381 if (h1_MC == nullptr or h1_RecoTrack == nullptr or h1_Track == nullptr) {
382 B2ERROR("One of the input histograms for function effPlot1D is a nullptr, "\
383 "can't create new histograms from this. Returning a nullptr.");
384 return nullptr;
385 }
386
387 std::string name1 = "_noGeoAcc";
388 std::string name2 = "_withGeoAcc";
389
390 std::string total1 = std::string(name) + name1;
391 std::string total2 = std::string(name) + name2;
392
393 std::string title1 = std::string(title) + name1;
394 std::string title2 = std::string(title) + name2;
395
396 TH1F* h[2];
397
398 h[0] = (TH1F*)duplicateHistogram(total2.c_str(), title2.c_str(), h1_RecoTrack, histoList);
399 h[0]->GetYaxis()->SetRangeUser(0., 1);
400
401 for (int bin = 0; bin < h[0]->GetXaxis()->GetNbins(); bin++) {
402 float num = h1_Track->GetBinContent(bin + 1);
403 float den = h1_RecoTrack->GetBinContent(bin + 1);
404 double eff = 0.;
405 double err = 0.;
406
407 if (den > 0) {
408 eff = (double)num / den;
409 err = sqrt(eff * (1 - eff)) / sqrt(den);
410 }
411 h[0]->SetBinContent(bin + 1, eff);
412 h[0]->SetBinError(bin + 1, err);
413 }
414
415 h[1] = (TH1F*)duplicateHistogram(total1.c_str(), title1.c_str(), h1_MC, histoList);
416 h[1]->GetYaxis()->SetRangeUser(0., 1);
417
418 for (int bin = 0; bin < h[1]->GetXaxis()->GetNbins(); bin++) {
419 float num = h1_Track->GetBinContent(bin + 1);
420 float den = h1_MC->GetBinContent(bin + 1);
421 double eff = 0.;
422 double err = 0.;
423
424 if (den > 0) {
425 eff = (double)num / den;
426 err = sqrt(eff * (1 - eff)) / sqrt(den);
427 }
428 h[1]->SetBinContent(bin + 1, eff);
429 h[1]->SetBinError(bin + 1, err);
430 }
431
432 if (histoList) {
433 histoList->Add(h[0]);
434 histoList->Add(h[1]);
435 }
436
437 return *h;
438}

◆ effPlot2D() [1/2]

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.

Definition at line 440 of file PerformanceEvaluationBaseClass.cc.

442{
443 if (h2_den == nullptr or h2_num == nullptr) {
444 B2ERROR("One of the input histograms for function effPlot1D is a nullptr, "\
445 "can't create new histograms from this. Returning a nullptr.");
446 return nullptr;
447 }
448
449 std::string err_char = "_error_";
450 std::string addTitle = "Errors, ";
451
452 std::string total;
453 std::string error;
454 std::string trueTitle;
455
456 std::string titleErr = addTitle + std::string(title);
457
458 if (geo_accettance == false) {
459 std::string name1 = "_noGeoAcc";
460 total = std::string(name) + name1;
461 trueTitle = std::string(title) + name1;
462 error = std::string(name) + err_char + std::string(name1);
463 } else {
464 std::string name2 = "_withGeoAcc";
465 total = std::string(name) + name2;
466 trueTitle = std::string(title) + name2;
467 error = std::string(name) + err_char + std::string(name2);
468 }
469
470 TH2F* h2[2];
471 h2[0] = (TH2F*)duplicateHistogram(total.c_str(), trueTitle.c_str(), h2_den, histoList);
472 h2[1] = (TH2F*)duplicateHistogram(error.c_str(), titleErr.c_str(), h2_den, histoList);
473
474 for (int binX = 0; binX < h2[0]->GetXaxis()->GetNbins(); binX++) {
475 for (int binY = 0; binY < h2[0]->GetYaxis()->GetNbins(); binY++) {
476 float num = h2_num->GetBinContent(binX + 1, binY + 1);
477 float den = h2_den->GetBinContent(binX + 1, binY + 1);
478 double eff = 0.;
479 double err = 0.;
480
481 if (den > 0) {
482 eff = (double)num / den;
483 err = sqrt(eff * (1 - eff)) / sqrt(den);
484 }
485
486 h2[0]->SetBinContent(binX + 1, binY + 1, eff);
487 h2[0]->SetBinError(binX + 1, binY + 1, err);
488 h2[1]->SetBinContent(binX + 1, binY + 1, err);
489 }
490 }
491
492
493 if (histoList) {
494 histoList->Add(h2[0]);
495 histoList->Add(h2[1]);
496 }
497
498 return *h2;
499
500}

◆ effPlot2D() [2/2]

TH2F * effPlot2D ( TH2F *  h2_MC,
TH2F *  h2_RecoTrack,
TH2F *  h2_Track,
const char *  name,
const char *  title,
TList *  histoList 
)

Create a 2D efficiency histogram and add it to the TList of 2D-histograms.

Definition at line 502 of file PerformanceEvaluationBaseClass.cc.

504{
505 if (h2_MC == nullptr or h2_RecoTrack == nullptr or h2_Track == nullptr) {
506 B2ERROR("One of the input histograms for function effPlot1D is a nullptr, "\
507 "can't create new histograms from this. Returning a nullptr.");
508 return nullptr;
509 }
510
511 std::string name1 = "_noGeoAcc";
512 std::string name2 = "_withGeoAcc";
513 std::string err_char = "_error_";
514 std::string addTitle = "Errors, ";
515
516 std::string total1 = std::string(name) + name1;
517 std::string total2 = std::string(name) + name2;
518
519 std::string title1 = std::string(title) + name1;
520 std::string title2 = std::string(title) + name2;
521
522 std::string error1 = std::string(name) + err_char + name1;
523 std::string error2 = std::string(name) + err_char + name2;
524 std::string titleErr = addTitle + std::string(title);
525
526 TH2F* h2[4];
527
528 h2[0] = (TH2F*)duplicateHistogram(total2.c_str(), title2.c_str(), h2_RecoTrack, histoList);
529 h2[1] = (TH2F*)duplicateHistogram(error2.c_str(), titleErr.c_str(), h2_RecoTrack, histoList);
530
531 for (int binX = 0; binX < h2[0]->GetXaxis()->GetNbins(); binX++) {
532 for (int binY = 0; binY < h2[0]->GetYaxis()->GetNbins(); binY++) {
533 float num = h2_Track->GetBinContent(binX + 1, binY + 1);
534 float den = h2_RecoTrack->GetBinContent(binX + 1, binY + 1);
535 double eff = 0.;
536 double err = 0.;
537
538 if (den > 0) {
539 eff = num / den;
540 err = sqrt(eff * (1 - eff)) / sqrt(den);
541 }
542
543 h2[0]->SetBinContent(binX + 1, binY + 1, eff);
544 h2[0]->SetBinError(binX + 1, binY + 1, err);
545 h2[1]->SetBinContent(binX + 1, binY + 1, err);
546 }
547 }
548
549 h2[2] = (TH2F*)duplicateHistogram(total1.c_str(), title1.c_str(), h2_MC, histoList);
550 h2[3] = (TH2F*)duplicateHistogram(error1.c_str(), titleErr.c_str(), h2_MC, histoList);
551
552 for (int binX = 0; binX < h2[2]->GetXaxis()->GetNbins(); binX++) {
553 for (int binY = 0; binY < h2[2]->GetYaxis()->GetNbins(); binY++) {
554 float num = h2_Track->GetBinContent(binX + 1, binY + 1);
555 float den = h2_MC->GetBinContent(binX + 1, binY + 1);
556 double eff = 0.;
557 double err = 0.;
558
559 if (den > 0) {
560 eff = num / den;
561 err = sqrt(eff * (1 - eff)) / sqrt(den);
562 }
563
564 h2[2]->SetBinContent(binX + 1, binY + 1, eff);
565 h2[2]->SetBinError(binX + 1, binY + 1, err);
566 h2[3]->SetBinContent(binX + 1, binY + 1, err);
567 }
568 }
569
570 if (histoList) {
571 histoList->Add(h2[0]);
572 histoList->Add(h2[1]);
573 histoList->Add(h2[2]);
574 histoList->Add(h2[3]);
575 }
576
577 return *h2;
578}

◆ geoAcc1D()

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.

Definition at line 580 of file PerformanceEvaluationBaseClass.cc.

581{
582 TH1F* h = (TH1F*)duplicateHistogram(name, title, h1_den, histoList);
583 h->GetYaxis()->SetRangeUser(0., 1);
584
585 for (int bin = 0; bin < h->GetXaxis()->GetNbins(); bin++) {
586 float num = h1_num->GetBinContent(bin + 1);
587 float den = h1_den->GetBinContent(bin + 1);
588 double eff = 0.;
589 double err = 0.;
590
591 if (den > 0) {
592 eff = (double)num / den;
593 err = sqrt(eff * (1 - eff)) / sqrt(den);
594 }
595 h->SetBinContent(bin + 1, eff);
596 h->SetBinError(bin + 1, err);
597 }
598
599 if (histoList)
600 histoList->Add(h);
601
602 return h;
603}

◆ geoAcc2D()

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.

Definition at line 605 of file PerformanceEvaluationBaseClass.cc.

607{
608 std::string err_char = "_err";
609 std::string addTitle = "Errors, ";
610
611 std::string error = std::string(name) + err_char;
612 std::string titleErr = addTitle + std::string(title);
613
614 TH2F* h2[2];
615 h2[0] = (TH2F*)duplicateHistogram(name, title, h2_den, histoList);
616 h2[1] = (TH2F*)duplicateHistogram(error.c_str(), titleErr.c_str(), h2_den, histoList);
617
618 for (int binX = 0; binX < h2[0]->GetXaxis()->GetNbins(); binX++) {
619 for (int binY = 0; binY < h2[0]->GetYaxis()->GetNbins(); binY++) {
620 float num = h2_num->GetBinContent(binX + 1, binY + 1);
621 float den = h2_den->GetBinContent(binX + 1, binY + 1);
622 double eff = 0.;
623 double err = 0.;
624
625 if (den > 0) {
626 eff = (double)num / den;
627 err = sqrt(eff * (1 - eff)) / sqrt(den);
628 }
629 h2[0]->SetBinContent(binX + 1, binY + 1, eff);
630 h2[0]->SetBinError(binX + 1, binY + 1, err);
631 h2[1]->SetBinContent(binX + 1, binY + 1, err);
632 }
633 }
634
635 if (histoList) {
636 histoList->Add(h2[0]);
637 histoList->Add(h2[1]);
638 }
639
640 return *h2;
641
642}

◆ V0FinderEff()

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.

Definition at line 644 of file PerformanceEvaluationBaseClass.cc.

646{
647
648 TH1F* h = (TH1F*)duplicateHistogram(name, title, h1_Mother, histoList);
649 h->GetYaxis()->SetRangeUser(0., 1);
650
651 for (int bin = 0; bin < h->GetXaxis()->GetNbins(); bin++) {
652 double dau0 = h1_dau0->GetBinContent(bin);
653 double dau1 = h1_dau1->GetBinContent(bin);
654 double Mother = h1_Mother->GetBinContent(bin);
655 double dau0Err = h1_dau0->GetBinError(bin);
656 double dau1Err = h1_dau1->GetBinError(bin);
657 double MotherErr = h1_Mother->GetBinError(bin);
658
659 double binCont = 1. * Mother / dau0 / dau1;
660 double binErr = binCont * sqrt((dau0Err / dau0) * (dau0Err / dau0) + (dau1Err / dau1) * (dau1Err / dau1) * (MotherErr / Mother) *
661 (MotherErr / Mother));
662
663 h->SetBinContent(bin, binCont);
664 h->SetBinError(bin, binErr);
665 }
666
667 if (histoList)
668 histoList->Add(h);
669
670 return h;
671}

Member Data Documentation

◆ m_histoList

TList* m_histoList = nullptr

List of performance-evaluation histograms.

Definition at line 38 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_efficiency

TList* m_histoList_efficiency = nullptr

List of efficiency histograms.

Definition at line 46 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_evtCharacterization

TList* m_histoList_evtCharacterization = nullptr

List of event-characterization histograms.

Definition at line 41 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_firstHit

TList* m_histoList_firstHit = nullptr

List of first-hit-position histograms.

Definition at line 43 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_fit

TList* m_histoList_fit = nullptr

List of track-fit histograms.

Definition at line 45 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_multiplicity

TList* m_histoList_multiplicity = nullptr

List of multiplicity histograms.

Definition at line 40 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_others

TList* m_histoList_others = nullptr

List of other performance-evaluation histograms.

Definition at line 48 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_pr

TList* m_histoList_pr = nullptr

List of pattern-recognition histograms.

Definition at line 44 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_purity

TList* m_histoList_purity = nullptr

List of purity histograms.

Definition at line 47 of file PerformanceEvaluationBaseClass.h.

◆ m_histoList_trkQuality

TList* m_histoList_trkQuality = nullptr

List of track-quality histograms.

Definition at line 42 of file PerformanceEvaluationBaseClass.h.

◆ m_rootFileName

std::string m_rootFileName

root file name

Definition at line 137 of file PerformanceEvaluationBaseClass.h.

◆ m_rootFilePtr

TFile* m_rootFilePtr = nullptr

pointer at root file used for storing histograms

Definition at line 140 of file PerformanceEvaluationBaseClass.h.


The documentation for this class was generated from the following files: