Logo ROOT  
Reference Guide
TEfficiency Class Reference

Class to handle efficiency histograms.

I. Overview

This class handles the calculation of efficiencies and their uncertainties. It provides several statistical methods for calculating frequentist and Bayesian confidence intervals as well as a function for combining several efficiencies.

Efficiencies have a lot of applications and meanings but in principle, they can be described by the fraction of good/passed events k out of sample containing N events. One is usually interested in the dependency of the efficiency on other (binned) variables. The number of passed and total events is therefore stored internally in two histograms (TEfficiency::fTotalHistogram and TEfficiency::fPassedHistogram). Then the efficiency, as well as its upper and lower error, can be calculated for each bin individually.

As the efficiency can be regarded as a parameter of a binomial distribution, the number of passed and total events must always be integer numbers. Therefore a filling with weights is not possible. However, you can assign a global weight to each TEfficiency object (TEfficiency::SetWeight). It is necessary to create one TEfficiency object for each weight if you investigate a process involving different weights. This procedure needs more effort but enables you to re-use the filled object in cases where you want to change one or more weights. This would not be possible if all events with different weights were filled in the same histogram.

II. Creating a TEfficiency object

If you start a new analysis, it is highly recommended to use the TEfficiency class from the beginning. You can then use one of the constructors for fixed or variable bin size and your desired dimension. These constructors append the created TEfficiency object to the current directory. So it will be written automatically to a file during the next TFile::Write command.

Example: create a two-dimensional TEfficiency object with

  • name = "eff"
  • title = "my efficiency"
  • axis titles: x, y and LaTeX-formatted epsilon as a label for Z axis
  • 10 bins with constant bin width (= 1) along X axis starting at 0 (lower edge from the first bin) up to 10 (upper edge of last bin)
  • 20 bins with constant bin width (= 0.5) along Y axis starting at -5 (lower edge from the first bin) up to 5 (upper edge of last bin)
      TEfficiency* pEff = new TEfficiency("eff","my efficiency;x;y;#epsilon",10,0,10,20,-5,5);
    
    If you already have two histograms filled with the number of passed and total events, you will use the constructor TEfficiency(const TH1& passed,const TH1& total) to construct the TEfficiency object. The histograms "passed" and "total" have to fulfill the conditions mentioned in TEfficiency::CheckConsistency, otherwise the construction will fail. As the histograms already exist, the new TEfficiency is by default not attached to the current directory to avoid duplication of data. If you want to store the new object anyway, you can either write it directly by calling TObject::Write or attach it to a directory using TEfficiency::SetDirectory. This also applies to TEfficiency objects created by the copy constructor TEfficiency::TEfficiency(const TEfficiency& rEff).

Example 1

TEfficiency* pEff = 0;
TFile* pFile = new TFile("myfile.root","recreate");
//h_pass and h_total are valid and consistent histograms
if(TEfficiency::CheckConsistency(h_pass,h_total))
{
pEff = new TEfficiency(h_pass,h_total);
// this will write the TEfficiency object to "myfile.root"
// AND pEff will be attached to the current directory
pEff->Write();
}
Class to handle efficiency histograms.
Definition: TEfficiency.h:28
TEfficiency()
default constructor
static Bool_t CheckConsistency(const TH1 &pass, const TH1 &total, Option_t *opt="")
Checks the consistence of the given histograms.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:53
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:796

Example 2

TEfficiency* pEff = 0;
TFile* pFile = new TFile("myfile.root","recreate");
//h_pass and h_total are valid and consistent histograms
if(TEfficiency::CheckConsistency(h_pass,h_total))
{
pEff = new TEfficiency(h_pass,h_total);
//this will attach the TEfficiency object to the current directory
//now all objects in gDirectory will be written to "myfile.root"
pFile->Write();
}
#define gDirectory
Definition: TDirectory.h:229
void SetDirectory(TDirectory *dir)
Sets the directory holding this TEfficiency object.
Int_t Write(const char *name=nullptr, Int_t opt=0, Int_t bufsiz=0) override
Write memory objects to this file.
Definition: TFile.cxx:2297

In case you already have two filled histograms and you only want to plot them as a graph, you should rather use TGraphAsymmErrors::TGraphAsymmErrors(const TH1* pass,const TH1* total,Option_t* opt) to create a graph object.

III. Filling with events

You can fill the TEfficiency object by calling the TEfficiency::Fill(Bool_t bPassed,Double_t x,Double_t y,Double_t z) method. The "bPassed" boolean flag indicates whether the current event is good (both histograms are filled) or not (only TEfficiency::fTotalHistogram is filled). The x, y and z variables determine the bin which is filled. For lower dimensions, the z- or even the y-value may be omitted.

{
//canvas only needed for this documentation
TCanvas* c1 = new TCanvas("example","",600,400);
c1->SetFillStyle(1001);
c1->SetFillColor(kWhite);
//create one-dimensional TEfficiency object with fixed bin size
TEfficiency* pEff = new TEfficiency("eff","my efficiency;x;#epsilon",20,0,10);
TRandom3 rand3;
bool bPassed;
double x;
for(int i=0; i<10000; ++i)
{
//simulate events with variable under investigation
x = rand3.Uniform(10);
//check selection: bPassed = DoesEventPassSelection(x)
bPassed = rand3.Rndm() < TMath::Gaus(x,5,4);
pEff->Fill(bPassed,x);
}
pEff->Draw("AP");
//only for this documentation
return c1;
}
@ kWhite
Definition: Rtypes.h:63
The Canvas class.
Definition: TCanvas.h:27
void Draw(Option_t *opt="")
Draws the current TEfficiency object.
void Fill(Bool_t bPassed, Double_t x, Double_t y=0, Double_t z=0)
This function is used for filling the two histograms.
Random number generator class based on M.
Definition: TRandom3.h:27
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom3.cxx:100
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:635
return c1
Definition: legend1.C:41
Double_t x[n]
Definition: legend1.C:17
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculate a gaussian function with mean and sigma.
Definition: TMath.cxx:448

You can also set the number of passed or total events for a bin directly by using the TEfficiency::SetPassedEvents or TEfficiency::SetTotalEvents method.

IV. Statistic options

The calculation of the estimated efficiency depends on the chosen statistic option. Let k denotes the number of passed events and N the number of total events.

Frequentist methods

The expectation value of the number of passed events is given by the true efficiency times the total number of events. One can estimate the efficiency by replacing the expected number of passed events by the observed number of passed events.

\[ k = \epsilon \times N \Rightarrow \hat{\varepsilon} = \frac{k}{N} \]

Bayesian methods

In Bayesian statistics a likelihood-function (how probable is it to get the observed data assuming a true efficiency) and a prior probability (what is the probability that a certain true efficiency is actually realised) are used to determine a posterior probability by using Bayes theorem. At the moment, only beta distributions (have 2 free parameters) are supported as prior probabilities.

\begin{eqnarray*} P(\epsilon | k ; N) &=& \frac{1}{norm} \times P(k | \epsilon ; N) \times Prior(\epsilon) \\ P(k | \epsilon ; N) &=& Binomial(N,k) \times \epsilon^{k} \times (1 - \epsilon)^{N - k} ...\ binomial\ distribution \\ Prior(\epsilon) &=& \frac{1}{B(\alpha,\beta)} \times \epsilon ^{\alpha - 1} \times (1 - \epsilon)^{\beta - 1} \equiv Beta(\epsilon; \alpha,\beta) \\ \Rightarrow P(\epsilon | k ; N) &=& \frac{1}{norm'} \times \epsilon^{k + \alpha - 1} \times (1 - \epsilon)^{N - k + \beta - 1} \equiv Beta(\epsilon; k + \alpha, N - k + \beta) \end{eqnarray*}

By default the expectation value of this posterior distribution is used as an estimator for the efficiency:

\[ \hat{\varepsilon} = \frac{k + \alpha}{N + \alpha + \beta} \]

Optionally the mode can also be used as a value for the estimated efficiency. This can be done by calling SetBit(kPosteriorMode) or TEfficiency::SetPosteriorMode. In this case, the estimated efficiency is:

\[ \hat{\varepsilon} = \frac{k + \alpha -1}{N + \alpha + \beta - 2} \]

In the case of a uniform prior distribution, B(x,1,1), the posterior mode is k/n, equivalent to the frequentist estimate (the maximum likelihood value).

The statistic options also specify which confidence interval is used for calculating the uncertainties of the efficiency. The following properties define the error calculation:

In the following table, the implemented confidence intervals are listed with their corresponding statistic option. For more details on the calculation, please have a look at the mentioned functions.

name statistic option function kIsBayesian parameters
Clopper-Pearson kFCP TEfficiency::ClopperPearson false total events, passed events, confidence level
normal approximation kFNormal TEfficiency::Normal false total events, passed events, confidence level
Wilson kFWilson TEfficiency::Wilson false total events, passed events, confidence level
Agresti-Coull kFAC TEfficiency::AgrestiCoull false total events, passed events. confidence level
Feldman-Cousins kFFC TEfficiency::FeldmanCousins false total events, passed events, confidence level
Mid-P Lancaster kMidP TEfficiency::MidPInterval false total events, passed events, confidence level
Jeffrey kBJeffrey TEfficiency::Bayesian true total events, passed events, confidence level, fBeta_alpha = 0.5, fBeta_beta = 0.5
Uniform prior kBUniform TEfficiency::Bayesian true total events, passed events, confidence level, fBeta_alpha = 1, fBeta_beta = 1
custom prior kBBayesian TEfficiency::Bayesian true total events, passed events, confidence level, fBeta_alpha, fBeta_beta

The following example demonstrates the effect of different statistic options and confidence levels.

{
//canvas only needed for the documentation
TCanvas* c1 = new TCanvas("c1","",600,400);
c1->Divide(2);
c1->SetFillStyle(1001);
c1->SetFillColor(kWhite);
//create one-dimensional TEfficiency object with fixed bin size
TEfficiency* pEff = new TEfficiency("eff","different confidence levels;x;#epsilon",20,0,10);
TRandom3 rand3;
bool bPassed;
double x;
for(int i=0; i<1000; ++i)
{
//simulate events with variable under investigation
x = rand3.Uniform(10);
//check selection: bPassed = DoesEventPassSelection(x)
bPassed = rand3.Rndm() < TMath::Gaus(x,5,4);
pEff->Fill(bPassed,x);
}
//set style attributes
pEff->SetFillStyle(3004);
//copy current TEfficiency object and set new confidence level
TEfficiency* pCopy = new TEfficiency(*pEff);
pCopy->SetConfidenceLevel(0.90);
//set style attributes
pCopy->SetFillStyle(3005);
c1->cd(1);
//add legend
TLegend* leg1 = new TLegend(0.3,0.1,0.7,0.5);
leg1->AddEntry(pEff,"95%","F");
leg1->AddEntry(pCopy,"68.3%","F");
pEff->Draw("A4");
pCopy->Draw("same4");
leg1->Draw("same");
//use same confidence level but different statistic methods
TEfficiency* pEff2 = new TEfficiency(*pEff);
TEfficiency* pCopy2 = new TEfficiency(*pEff);
pEff2->SetTitle("different statistic options;x;#epsilon");
//set style attributes
pCopy2->SetFillStyle(3005);
pCopy2->SetFillColor(kBlue);
c1->cd(2);
//add legend
TLegend* leg2 = new TLegend(0.3,0.1,0.7,0.5);
leg2->AddEntry(pEff2,"kFNormal","F");
leg2->AddEntry(pCopy2,"kFAC","F");
pEff2->Draw("a4");
pCopy2->Draw("same4");
leg2->Draw("same");
//only for this documentation
c1->cd(0);
return c1;
}
@ kRed
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
void SetStatisticOption(EStatOption option)
Sets the statistic option which affects the calculation of the confidence interval.
void SetConfidenceLevel(Double_t level)
Sets the confidence level (0 < level < 1) The default value is 1-sigma :~ 0.683.
void SetTitle(const char *title)
Sets the title.
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:423

The prior probability of the efficiency in Bayesian statistics can be given in terms of a beta distribution. The beta distribution has two positive shape parameters. The resulting priors for different combinations of these shape parameters are shown in the plot below.

{
//canvas only needed for the documentation
TCanvas* c1 = new TCanvas("c1","",600,400);
c1->SetFillStyle(1001);
c1->SetFillColor(kWhite);
//create different beta distributions
TF1* f1 = new TF1("f1","TMath::BetaDist(x,1,1)",0,1);
TF1* f2 = new TF1("f2","TMath::BetaDist(x,0.5,0.5)",0,1);
TF1* f3 = new TF1("f3","TMath::BetaDist(x,1,5)",0,1);
f3->SetTitle("Beta distributions as priors;#epsilon;P(#epsilon)");
TF1* f4 = new TF1("f4","TMath::BetaDist(x,4,3)",0,1);
//add legend
TLegend* leg = new TLegend(0.25,0.5,0.85,0.89);
leg->SetFillColor(kWhite);
leg->SetFillStyle(1001);
leg->AddEntry(f1,"a=1, b=1","L");
leg->AddEntry(f2,"a=0.5, b=0.5","L");
leg->AddEntry(f3,"a=1, b=5","L");
leg->AddEntry(f4,"a=4, b=3","L");
f3->Draw();
f1->Draw("same");
f2->Draw("Same");
f4->Draw("same");
leg->Draw("same");
//only for this documentation
return c1;
}
@ kGreen
Definition: Rtypes.h:64
@ kViolet
Definition: Rtypes.h:65
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
1-Dim function class
Definition: TF1.h:210
virtual void SetTitle(const char *title="")
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition: TF1.cxx:3551
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1320
TF1 * f1
Definition: legend1.C:11
leg
Definition: legend1.C:34

IV.1 Coverage probabilities for different methods

The following pictures illustrate the actual coverage probability for the different values of the true efficiency and the total number of events when a confidence level of 95% is desired.

Normal Approximation
Wilson
Agresti Coull
Clopper Pearson
Bayesian with Uniform Prior
Bayesian with Jeffrey Prior

The average (over all possible true efficiencies) coverage probability for different number of total events is shown in the next picture.

Average Coverage

V. Merging and combining TEfficiency objects

In many applications, the efficiency should be calculated for an inhomogeneous sample in the sense that it contains events with different weights. In order to be able to determine the correct overall efficiency, it is necessary to use for each subsample (= all events with the same weight) a different TEfficiency object. After finishing your analysis you can then construct the overall efficiency with its uncertainty.

This procedure has the advantage that you can change the weight of one subsample easily without rerunning the whole analysis. On the other hand, more effort is needed to handle several TEfficiency objects instead of one histogram. In the case of many different or even continuously distributed weights, this approach becomes cumbersome. One possibility to overcome this problem is the usage of binned weights.

Example

In particle physics weights arises from the fact that you want to normalise your results to a certain reference value. A very common formula for calculating weights is

\begin{eqnarray*} w &=& \frac{\sigma L}{N_{gen} \epsilon_{trig}} \\ &-& \sigma ...\ cross\ section \\ &-& L ...\ luminosity \\ &-& N_{gen}\ ... number\ of\ generated\ events \\ &-& \epsilon_{trig}\ ...\ (known)\ trigger\ efficiency \\ \end{eqnarray*}

The reason for different weights can therefore be:

  • different processes
  • other integrated luminosity
  • varying trigger efficiency
  • different sample sizes
  • ...
  • or even combination of them

Depending on the actual meaning of different weights in your case, you should either merge or combine them to get the overall efficiency.

V.1 When should I use merging?

If the weights are artificial and do not represent real alternative hypotheses, you should merge the different TEfficiency objects. That means especially for the Bayesian case that the prior probability should be the same for all merged TEfficiency objects. The merging can be done by invoking one of the following operations:

  • eff1.Add(eff2)
  • eff1 += eff2
  • eff1 = eff1 + eff2

The result of the merging is stored in the TEfficiency object which is marked bold above. The contents of the internal histograms of both TEfficiency objects are added and a new weight is assigned. The statistic options are not changed.

\[ \frac{1}{w_{new}} = \frac{1}{w_{1}} + \frac{1}{w_{2}} \]

Example:

If you use two samples with different numbers of generated events for the same process and you want to normalise both to the same integrated luminosity and trigger efficiency, the different weights then arise just from the fact that you have different numbers of events. The TEfficiency objects should be merged because the samples do not represent true alternatives. You expect the same result as if you would have a big sample with all events in it.

\[ w_{1} = \frac{\sigma L}{\epsilon N_{1}}, w_{2} = \frac{\sigma L}{\epsilon N_{2}} \Rightarrow w_{new} = \frac{\sigma L}{\epsilon (N_{1} + N_{2})} = \frac{1}{\frac{1}{w_{1}} + \frac{1}{w_{2}}} \]

V.2 When should I use combining?

You should combine TEfficiency objects whenever the weights represent alternatives processes for the efficiency. As the combination of two TEfficiency objects is not always consistent with the representation by two internal histograms, the result is not stored in a TEfficiency object but a TGraphAsymmErrors is returned which shows the estimated combined efficiency and its uncertainty for each bin. At the moment the combination method TEfficiency::Combine only supports a combination of 1-dimensional efficiencies in a Bayesian approach.

For calculating the combined efficiency and its uncertainty for each bin only Bayesian statistics is used. No frequentists methods are presently supported for computing the combined efficiency and its confidence interval. In the case of the Bayesian statistics, a combined posterior is constructed taking into account the weight of each TEfficiency object. The same prior is used for all the TEfficiency objects.

\begin{eqnarray*} P_{comb}(\epsilon | {w_{i}}, {k_{i}} , {N_{i}}) = \frac{1}{norm} \prod_{i}{L(k_{i} | N_{i}, \epsilon)}^{w_{i}} \Pi( \epsilon )\\ L(k_{i} | N_{i}, \epsilon)\ is\ the\ likelihood\ function\ for\ the\ sample\ i\ (a\ Binomial\ distribution)\\ \Pi( \epsilon)\ is\ the\ prior,\ a\ beta\ distribution\ B(\epsilon, \alpha, \beta).\\ The\ resulting\ combined\ posterior\ is \\ P_{comb}(\epsilon |{w_{i}}; {k_{i}}; {N_{i}}) = B(\epsilon, \sum_{i}{ w_{i} k_{i}} + \alpha, \sum_{i}{ w_{i}(n_{i}-k_{i})}+\beta) \\ \hat{\varepsilon} = \int_{0}^{1} \epsilon \times P_{comb}(\epsilon | {k_{i}} , {N_{i}}) d\epsilon \\ confidence\ level = 1 - \alpha \\ \frac{\alpha}{2} = \int_{0}^{\epsilon_{low}} P_{comb}(\epsilon | {k_{i}} , {N_{i}}) d\epsilon ...\ defines\ lower\ boundary \\ 1- \frac{\alpha}{2} = \int_{0}^{\epsilon_{up}} P_{comb}(\epsilon | {k_{i}} , {N_{i}}) d\epsilon ...\ defines\ upper\ boundary \end{eqnarray*}

Example:

If you use cuts to select electrons which can originate from two different processes, you can determine the selection efficiency for each process. The overall selection efficiency is then the combined efficiency. The weights to be used in the combination should be the probability that an electron comes from the corresponding process.

\[ p_{1} = \frac{\sigma_{1}}{\sigma_{1} + \sigma_{2}} = \frac{N_{1}w_{1}}{N_{1}w_{1} + N_{2}w_{2}}\\ p_{2} = \frac{\sigma_{2}}{\sigma_{1} + \sigma_{2}} = \frac{N_{2}w_{2}}{N_{1}w_{1} + N_{2}w_{2}} \]

VI. Further operations

VI.Information about the internal histograms

The methods TEfficiency::GetPassedHistogram and TEfficiency::GetTotalHistogram return a constant pointer to the internal histograms. They can be used to obtain information about the internal histograms (e.g., the binning, number of passed / total events in a bin, mean values...). One can obtain a clone of the internal histograms by calling TEfficiency::GetCopyPassedHisto or TEfficiency::GetCopyTotalHisto. The returned histograms are completely independent from the current TEfficiency object. By default, they are not attached to a directory to avoid the duplication of data and the user is responsible for deleting them.

//open a root file which contains a TEfficiency object
TFile* pFile = new TFile("myfile.root","update");
//get TEfficiency object with name "my_eff"
TEfficiency* pEff = (TEfficiency*)pFile->Get("my_eff");
//get clone of total histogram
TH1* clone = pEff->GetCopyTotalHisto();
//change clone...
//save changes of clone directly
clone->Write();
//or append it to the current directory and write the file
//clone->SetDirectory(gDirectory);
//pFile->Write();
//delete histogram object
delete clone;
clone = 0;
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
TH1 * GetCopyTotalHisto() const
Returns a cloned version of fTotalHistogram.
The TH1 histogram class.
Definition: TH1.h:56

It is also possible to set the internal total or passed histogram by using the methods TEfficiency::SetPassedHistogram or TEfficiency::SetTotalHistogram.

In order to ensure the validity of the TEfficiency object, the consistency of the new histogram and the stored histogram is checked. It might be impossible sometimes to change the histograms in a consistent way. Therefore one can force the replacement by passing the "f" option. Then the user has to ensure that the other internal histogram is replaced as well and that the TEfficiency object is in a valid state.

VI.2 Fitting

The efficiency can be fitted using the TEfficiency::Fit function which internally uses the TBinomialEfficiencyFitter::Fit method. As this method is using a maximum-likelihood-fit, it is necessary to initialise the given fit function with reasonable start values. The resulting fit function is attached to the list of associated functions and will be drawn automatically during the next TEfficiency::Draw command. The list of associated function can be modified by using the pointer returned by TEfficiency::GetListOfFunctions.

{
//canvas only needed for this documentation
TCanvas* c1 = new TCanvas("example","",600,400);
c1->SetFillStyle(1001);
c1->SetFillColor(kWhite);
//create one-dimensional TEfficiency object with fixed bin size
TEfficiency* pEff = new TEfficiency("eff","my efficiency;x;#epsilon",20,0,10);
TRandom3 rand3;
bool bPassed;
double x;
for(int i=0; i<10000; ++i)
{
//simulate events with variable under investigation
x = rand3.Uniform(10);
//check selection: bPassed = DoesEventPassSelection(x)
bPassed = rand3.Rndm() < TMath::Gaus(x,5,4);
pEff->Fill(bPassed,x);
}
//create a function for fitting and do the fit
TF1* f1 = new TF1("f1","gaus",0,10);
f1->SetParameters(1,5,2);
pEff->Fit(f1);
//create a threshold function
TF1* f2 = new TF1("thres","0.8",0,10);
//add it to the list of functions
//use add first because the parameters of the last function will be displayed
pEff->Draw("AP");
//only for this documentation
return c1;
}
TList * GetListOfFunctions()
TFitResultPtr Fit(TF1 *f1, Option_t *opt="")
Fits the efficiency using the TBinomialEfficiencyFitter class.
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
virtual void AddFirst(TObject *obj)
Add object at the beginning of the list.
Definition: TList.cxx:99

VI.3 Draw a TEfficiency object

A TEfficiency object can be drawn by calling the usual TEfficiency::Draw method. At the moment drawing is only supported for 1- and 2-dimensional TEfficiency objects. In the 1-dimensional case, you can use the same options as for the TGraphAsymmErrors::Draw method. For 2-dimensional TEfficiency objects, you can pass the same options as for a TH2::Draw object.

Definition at line 27 of file TEfficiency.h.

Public Types

enum  EStatOption {
  kFCP = 0 , kFNormal , kFWilson , kFAC ,
  kFFC , kBJeffrey , kBUniform , kBBayesian ,
  kMidP
}
 
- Public Types inherited from TObject
enum  {
  kIsOnHeap = 0x01000000 , kNotDeleted = 0x02000000 , kZombie = 0x04000000 , kInconsistent = 0x08000000 ,
  kBitMask = 0x00ffffff
}
 
enum  { kSingleKey = BIT(0) , kOverwrite = BIT(1) , kWriteDelete = BIT(2) }
 
enum  EDeprecatedStatusBits { kObjInCanvas = BIT(3) }
 
enum  EStatusBits {
  kCanDelete = BIT(0) , kMustCleanup = BIT(3) , kIsReferenced = BIT(4) , kHasUUID = BIT(5) ,
  kCannotPick = BIT(6) , kNoContextMenu = BIT(8) , kInvalidObject = BIT(13)
}
 

Public Member Functions

 TEfficiency ()
 default constructor More...
 
 TEfficiency (const char *name, const char *title, Int_t nbins, const Double_t *xbins)
 Create 1-dimensional TEfficiency object with variable bin size. More...
 
 TEfficiency (const char *name, const char *title, Int_t nbins, Double_t xlow, Double_t xup)
 Create 1-dimensional TEfficiency object with fixed bins size. More...
 
 TEfficiency (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, const Double_t *ybins)
 Create 2-dimensional TEfficiency object with variable bin size. More...
 
 TEfficiency (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, const Double_t *ybins, Int_t nbinsz, const Double_t *zbins)
 Create 3-dimensional TEfficiency object with variable bin size. More...
 
 TEfficiency (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup)
 Create 2-dimensional TEfficiency object with fixed bin size. More...
 
 TEfficiency (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup, Int_t nbinsz, Double_t zlow, Double_t zup)
 Create 3-dimensional TEfficiency object with fixed bin size. More...
 
 TEfficiency (const TEfficiency &heff)
 Copy constructor. More...
 
 TEfficiency (const TH1 &passed, const TH1 &total)
 constructor using two existing histograms as input More...
 
 ~TEfficiency ()
 default destructor More...
 
void Add (const TEfficiency &rEff)
 
void Browse (TBrowser *)
 Browse object. May be overridden for another default action. More...
 
TGraphAsymmErrorsCreateGraph (Option_t *opt="") const
 Create the graph used be painted (for dim=1 TEfficiency) The return object is managed by the caller. More...
 
TH2CreateHistogram (Option_t *opt="") const
 Create the histogram used to be painted (for dim=2 TEfficiency) The return object is managed by the caller. More...
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 Compute distance from point px,py to a graph. More...
 
void Draw (Option_t *opt="")
 Draws the current TEfficiency object. More...
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 Execute action corresponding to one event. More...
 
void Fill (Bool_t bPassed, Double_t x, Double_t y=0, Double_t z=0)
 This function is used for filling the two histograms. More...
 
void FillWeighted (Bool_t bPassed, Double_t weight, Double_t x, Double_t y=0, Double_t z=0)
 This function is used for filling the two histograms with a weight. More...
 
Int_t FindFixBin (Double_t x, Double_t y=0, Double_t z=0) const
 Returns the global bin number containing the given values. More...
 
TFitResultPtr Fit (TF1 *f1, Option_t *opt="")
 Fits the efficiency using the TBinomialEfficiencyFitter class. More...
 
Double_t GetBetaAlpha (Int_t bin=-1) const
 
Double_t GetBetaBeta (Int_t bin=-1) const
 
Double_t GetConfidenceLevel () const
 
TH1GetCopyPassedHisto () const
 Returns a cloned version of fPassedHistogram. More...
 
TH1GetCopyTotalHisto () const
 Returns a cloned version of fTotalHistogram. More...
 
Int_t GetDimension () const
 returns the dimension of the current TEfficiency object More...
 
TDirectoryGetDirectory () const
 
Double_t GetEfficiency (Int_t bin) const
 Returns the efficiency in the given global bin. More...
 
Double_t GetEfficiencyErrorLow (Int_t bin) const
 Returns the lower error on the efficiency in the given global bin. More...
 
Double_t GetEfficiencyErrorUp (Int_t bin) const
 Returns the upper error on the efficiency in the given global bin. More...
 
Int_t GetGlobalBin (Int_t binx, Int_t biny=0, Int_t binz=0) const
 Returns the global bin number which can be used as argument for the following functions: More...
 
TListGetListOfFunctions ()
 
TGraphAsymmErrorsGetPaintedGraph () const
 
TH2GetPaintedHistogram () const
 
const TH1GetPassedHistogram () const
 
EStatOption GetStatisticOption () const
 
const TH1GetTotalHistogram () const
 
Double_t GetWeight () const
 
Long64_t Merge (TCollection *list)
 Merges the TEfficiency objects in the given list to the given TEfficiency object using the operator+=(TEfficiency&) More...
 
TEfficiencyoperator+= (const TEfficiency &rhs)
 Adds the histograms of another TEfficiency object to current histograms. More...
 
TEfficiencyoperator= (const TEfficiency &rhs)
 Assignment operator. More...
 
void Paint (Option_t *opt)
 Paints this TEfficiency object. More...
 
void SavePrimitive (std::ostream &out, Option_t *opt="")
 Have histograms fixed bins along each axis? More...
 
void SetBetaAlpha (Double_t alpha)
 Sets the shape parameter α. More...
 
void SetBetaBeta (Double_t beta)
 Sets the shape parameter β. More...
 
void SetBetaBinParameters (Int_t bin, Double_t alpha, Double_t beta)
 Sets different shape parameter α and β for the prior distribution for each bin. More...
 
Bool_t SetBins (Int_t nx, const Double_t *xBins)
 Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost. More...
 
Bool_t SetBins (Int_t nx, const Double_t *xBins, Int_t ny, const Double_t *yBins)
 Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost. More...
 
Bool_t SetBins (Int_t nx, const Double_t *xBins, Int_t ny, const Double_t *yBins, Int_t nz, const Double_t *zBins)
 Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost. More...
 
Bool_t SetBins (Int_t nx, Double_t xmin, Double_t xmax)
 Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost. More...
 
Bool_t SetBins (Int_t nx, Double_t xmin, Double_t xmax, Int_t ny, Double_t ymin, Double_t ymax)
 Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost. More...
 
Bool_t SetBins (Int_t nx, Double_t xmin, Double_t xmax, Int_t ny, Double_t ymin, Double_t ymax, Int_t nz, Double_t zmin, Double_t zmax)
 Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost. More...
 
void SetCentralInterval (Bool_t on=true)
 
void SetConfidenceLevel (Double_t level)
 Sets the confidence level (0 < level < 1) The default value is 1-sigma :~ 0.683. More...
 
void SetDirectory (TDirectory *dir)
 Sets the directory holding this TEfficiency object. More...
 
void SetName (const char *name)
 Sets the name. More...
 
Bool_t SetPassedEvents (Int_t bin, Int_t events)
 Sets the number of passed events in the given global bin. More...
 
Bool_t SetPassedHistogram (const TH1 &rPassed, Option_t *opt)
 Sets the histogram containing the passed events. More...
 
void SetPosteriorAverage (Bool_t on=true)
 
void SetPosteriorMode (Bool_t on=true)
 
void SetShortestInterval (Bool_t on=true)
 
void SetStatisticOption (EStatOption option)
 Sets the statistic option which affects the calculation of the confidence interval. More...
 
void SetTitle (const char *title)
 Sets the title. More...
 
Bool_t SetTotalEvents (Int_t bin, Int_t events)
 Sets the number of total events in the given global bin. More...
 
Bool_t SetTotalHistogram (const TH1 &rTotal, Option_t *opt)
 Sets the histogram containing all events. More...
 
void SetUseWeightedEvents (Bool_t on=kTRUE)
 
void SetWeight (Double_t weight)
 Sets the global weight for this TEfficiency object. More...
 
Bool_t UsesBayesianStat () const
 
Bool_t UsesCentralInterval () const
 
Bool_t UsesPosteriorAverage () const
 
Bool_t UsesPosteriorMode () const
 
Bool_t UsesShortestInterval () const
 
Bool_t UsesWeights () const
 
- Public Member Functions inherited from TNamed
 TNamed ()
 
 TNamed (const char *name, const char *title)
 
 TNamed (const TNamed &named)
 TNamed copy ctor. More...
 
 TNamed (const TString &name, const TString &title)
 
virtual ~TNamed ()
 TNamed destructor. More...
 
virtual void Clear (Option_t *option="")
 Set name and title to empty strings (""). More...
 
virtual TObjectClone (const char *newname="") const
 Make a clone of an object using the Streamer facility. More...
 
virtual Int_t Compare (const TObject *obj) const
 Compare two TNamed objects. More...
 
virtual void Copy (TObject &named) const
 Copy this to obj. More...
 
virtual void FillBuffer (char *&buffer)
 Encode TNamed into output buffer. More...
 
virtual const char * GetName () const
 Returns name of object. More...
 
virtual const char * GetTitle () const
 Returns title of object. More...
 
virtual ULong_t Hash () const
 Return hash value for this object. More...
 
virtual Bool_t IsSortable () const
 
virtual void ls (Option_t *option="") const
 List TNamed name and title. More...
 
TNamedoperator= (const TNamed &rhs)
 TNamed assignment operator. More...
 
virtual void Print (Option_t *option="") const
 Print TNamed name and title. More...
 
virtual void SetName (const char *name)
 Set the name of the TNamed. More...
 
virtual void SetNameTitle (const char *name, const char *title)
 Set all the TNamed parameters (name and title). More...
 
virtual void SetTitle (const char *title="")
 Set the title of the TNamed. More...
 
virtual Int_t Sizeof () const
 Return size of the TNamed part of the TObject. More...
 
- Public Member Functions inherited from TObject
 TObject ()
 TObject constructor. More...
 
 TObject (const TObject &object)
 TObject copy ctor. More...
 
virtual ~TObject ()
 TObject destructor. More...
 
void AbstractMethod (const char *method) const
 Use this method to implement an "abstract" method that you don't want to leave purely abstract. More...
 
virtual void AppendPad (Option_t *option="")
 Append graphics object to current pad. More...
 
virtual void Browse (TBrowser *b)
 Browse object. May be overridden for another default action. More...
 
ULong_t CheckedHash ()
 Check and record whether this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object. More...
 
virtual const char * ClassName () const
 Returns name of class to which the object belongs. More...
 
virtual void Clear (Option_t *="")
 
virtual TObjectClone (const char *newname="") const
 Make a clone of an object using the Streamer facility. More...
 
virtual Int_t Compare (const TObject *obj) const
 Compare abstract method. More...
 
virtual void Copy (TObject &object) const
 Copy this to obj. More...
 
virtual void Delete (Option_t *option="")
 Delete this object. More...
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 Computes distance from point (px,py) to the object. More...
 
virtual void Draw (Option_t *option="")
 Default Draw method for all objects. More...
 
virtual void DrawClass () const
 Draw class inheritance tree of the class to which this object belongs. More...
 
virtual TObjectDrawClone (Option_t *option="") const
 Draw a clone of this object in the current selected pad for instance with: gROOT->SetSelectedPad(gPad). More...
 
virtual void Dump () const
 Dump contents of object on stdout. More...
 
virtual void Error (const char *method, const char *msgfmt,...) const
 Issue error message. More...
 
virtual void Execute (const char *method, const char *params, Int_t *error=0)
 Execute method on this object with the given parameter string, e.g. More...
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=0)
 Execute method on this object with parameters stored in the TObjArray. More...
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 Execute action corresponding to an event at (px,py). More...
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 Issue fatal error message. More...
 
virtual TObjectFindObject (const char *name) const
 Must be redefined in derived classes. More...
 
virtual TObjectFindObject (const TObject *obj) const
 Must be redefined in derived classes. More...
 
virtual Option_tGetDrawOption () const
 Get option used by the graphics system to draw this object. More...
 
virtual const char * GetIconName () const
 Returns mime type name of object. More...
 
virtual const char * GetName () const
 Returns name of object. More...
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const
 Returns string containing info about the object at position (px,py). More...
 
virtual Option_tGetOption () const
 
virtual const char * GetTitle () const
 Returns title of object. More...
 
virtual UInt_t GetUniqueID () const
 Return the unique object id. More...
 
virtual Bool_t HandleTimer (TTimer *timer)
 Execute action in response of a timer timing out. More...
 
virtual ULong_t Hash () const
 Return hash value for this object. More...
 
Bool_t HasInconsistentHash () const
 Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e. More...
 
virtual void Info (const char *method, const char *msgfmt,...) const
 Issue info message. More...
 
virtual Bool_t InheritsFrom (const char *classname) const
 Returns kTRUE if object inherits from class "classname". More...
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 Returns kTRUE if object inherits from TClass cl. More...
 
virtual void Inspect () const
 Dump contents of this object in a graphics canvas. More...
 
void InvertBit (UInt_t f)
 
virtual Bool_t IsEqual (const TObject *obj) const
 Default equal comparison (objects are equal if they have the same address in memory). More...
 
virtual Bool_t IsFolder () const
 Returns kTRUE in case object contains browsable objects (like containers or lists of other objects). More...
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
virtual Bool_t IsSortable () const
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
virtual void ls (Option_t *option="") const
 The ls function lists the contents of a class on stdout. More...
 
void MayNotUse (const char *method) const
 Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary). More...
 
virtual Bool_t Notify ()
 This method must be overridden to handle object notification. More...
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 Use this method to declare a method obsolete. More...
 
void operator delete (void *ptr)
 Operator delete. More...
 
void operator delete[] (void *ptr)
 Operator delete []. More...
 
voidoperator new (size_t sz)
 
voidoperator new (size_t sz, void *vp)
 
voidoperator new[] (size_t sz)
 
voidoperator new[] (size_t sz, void *vp)
 
TObjectoperator= (const TObject &rhs)
 TObject assignment operator. More...
 
virtual void Paint (Option_t *option="")
 This method must be overridden if a class wants to paint itself. More...
 
virtual void Pop ()
 Pop on object drawn in a pad to the top of the display list. More...
 
virtual void Print (Option_t *option="") const
 This method must be overridden when a class wants to print itself. More...
 
virtual Int_t Read (const char *name)
 Read contents of object with specified name from the current directory. More...
 
virtual void RecursiveRemove (TObject *obj)
 Recursively remove this object from a list. More...
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 Save this object in the file specified by filename. More...
 
virtual void SavePrimitive (std::ostream &out, Option_t *option="")
 Save a primitive as a C++ statement(s) on output stream "out". More...
 
void SetBit (UInt_t f)
 
void SetBit (UInt_t f, Bool_t set)
 Set or unset the user status bits as specified in f. More...
 
virtual void SetDrawOption (Option_t *option="")
 Set drawing option for object. More...
 
virtual void SetUniqueID (UInt_t uid)
 Set the unique object id. More...
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 Issue system error message. More...
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
virtual void UseCurrentStyle ()
 Set current style settings in this object This function is called when either TCanvas::UseCurrentStyle or TROOT::ForceStyle have been invoked. More...
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 Issue warning message. More...
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0)
 Write this object to the current directory. More...
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0) const
 Write this object to the current directory. More...
 
- Public Member Functions inherited from TAttLine
 TAttLine ()
 AttLine default constructor. More...
 
 TAttLine (Color_t lcolor, Style_t lstyle, Width_t lwidth)
 AttLine normal constructor. More...
 
virtual ~TAttLine ()
 AttLine destructor. More...
 
void Copy (TAttLine &attline) const
 Copy this line attributes to a new TAttLine. More...
 
Int_t DistancetoLine (Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
 Compute distance from point px,py to a line. More...
 
virtual Color_t GetLineColor () const
 Return the line color. More...
 
virtual Style_t GetLineStyle () const
 Return the line style. More...
 
virtual Width_t GetLineWidth () const
 Return the line width. More...
 
virtual void Modify ()
 Change current line attributes if necessary. More...
 
virtual void ResetAttLine (Option_t *option="")
 Reset this line attributes to default values. More...
 
virtual void SaveLineAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
 Save line attributes as C++ statement(s) on output stream out. More...
 
virtual void SetLineAttributes ()
 Invoke the DialogCanvas Line attributes. More...
 
virtual void SetLineColor (Color_t lcolor)
 Set the line color. More...
 
virtual void SetLineColorAlpha (Color_t lcolor, Float_t lalpha)
 Set a transparent line color. More...
 
virtual void SetLineStyle (Style_t lstyle)
 Set the line style. More...
 
virtual void SetLineWidth (Width_t lwidth)
 Set the line width. More...
 
- Public Member Functions inherited from TAttFill
 TAttFill ()
 AttFill default constructor. More...
 
 TAttFill (Color_t fcolor, Style_t fstyle)
 AttFill normal constructor. More...
 
virtual ~TAttFill ()
 AttFill destructor. More...
 
void Copy (TAttFill &attfill) const
 Copy this fill attributes to a new TAttFill. More...
 
virtual Color_t GetFillColor () const
 Return the fill area color. More...
 
virtual Style_t GetFillStyle () const
 Return the fill area style. More...
 
virtual Bool_t IsTransparent () const
 
virtual void Modify ()
 Change current fill area attributes if necessary. More...
 
virtual void ResetAttFill (Option_t *option="")
 Reset this fill attributes to default values. More...
 
virtual void SaveFillAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
 Save fill attributes as C++ statement(s) on output stream out. More...
 
virtual void SetFillAttributes ()
 Invoke the DialogCanvas Fill attributes. More...
 
virtual void SetFillColor (Color_t fcolor)
 Set the fill area color. More...
 
virtual void SetFillColorAlpha (Color_t fcolor, Float_t falpha)
 Set a transparent fill color. More...
 
virtual void SetFillStyle (Style_t fstyle)
 Set the fill area style. More...
 
- Public Member Functions inherited from TAttMarker
 TAttMarker ()
 TAttMarker default constructor. More...
 
 TAttMarker (Color_t color, Style_t style, Size_t msize)
 TAttMarker normal constructor. More...
 
virtual ~TAttMarker ()
 TAttMarker destructor. More...
 
void Copy (TAttMarker &attmarker) const
 Copy this marker attributes to a new TAttMarker. More...
 
virtual Color_t GetMarkerColor () const
 Return the marker color. More...
 
virtual Size_t GetMarkerSize () const
 Return the marker size. More...
 
virtual Style_t GetMarkerStyle () const
 Return the marker style. More...
 
virtual void Modify ()
 Change current marker attributes if necessary. More...
 
virtual void ResetAttMarker (Option_t *toption="")
 Reset this marker attributes to the default values. More...
 
virtual void SaveMarkerAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
 Save line attributes as C++ statement(s) on output stream out. More...
 
virtual void SetMarkerAttributes ()
 Invoke the DialogCanvas Marker attributes. More...
 
virtual void SetMarkerColor (Color_t mcolor=1)
 Set the marker color. More...
 
virtual void SetMarkerColorAlpha (Color_t mcolor, Float_t malpha)
 Set a transparent marker color. More...
 
virtual void SetMarkerSize (Size_t msize=1)
 Set the marker size. More...
 
virtual void SetMarkerStyle (Style_t mstyle=1)
 Set the marker style. More...
 

Static Public Member Functions

static Double_t AgrestiCoull (Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
 Calculates the boundaries for the frequentist Agresti-Coull interval. More...
 
static Double_t Bayesian (Double_t total, Double_t passed, Double_t level, Double_t alpha, Double_t beta, Bool_t bUpper, Bool_t bShortest=false)
 Calculates the boundaries for a Bayesian confidence interval (shortest or central interval depending on the option) More...
 
static Double_t BetaCentralInterval (Double_t level, Double_t alpha, Double_t beta, Bool_t bUpper)
 Calculates the boundaries for a central confidence interval for a Beta distribution. More...
 
static Double_t BetaMean (Double_t alpha, Double_t beta)
 Compute the mean (average) of the beta distribution. More...
 
static Double_t BetaMode (Double_t alpha, Double_t beta)
 Compute the mode of the beta distribution. More...
 
static Bool_t BetaShortestInterval (Double_t level, Double_t alpha, Double_t beta, Double_t &lower, Double_t &upper)
 Calculates the boundaries for a shortest confidence interval for a Beta distribution. More...
 
static Bool_t CheckBinning (const TH1 &pass, const TH1 &total)
 Checks binning for each axis. More...
 
static Bool_t CheckConsistency (const TH1 &pass, const TH1 &total, Option_t *opt="")
 Checks the consistence of the given histograms. More...
 
static Bool_t CheckEntries (const TH1 &pass, const TH1 &total, Option_t *opt="")
 Checks whether bin contents are compatible with binomial statistics. More...
 
static Bool_t CheckWeights (const TH1 &pass, const TH1 &total)
 Check if both histogram are weighted. More...
 
static Double_t ClopperPearson (Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
 Calculates the boundaries for the frequentist Clopper-Pearson interval. More...
 
static Double_t Combine (Double_t &up, Double_t &low, Int_t n, const Int_t *pass, const Int_t *total, Double_t alpha, Double_t beta, Double_t level=0.683, const Double_t *w=0, Option_t *opt="")
 
static TGraphAsymmErrorsCombine (TCollection *pList, Option_t *opt="", Int_t n=0, const Double_t *w=0)
 Combines a list of 1-dimensional TEfficiency objects. More...
 
static Double_t FeldmanCousins (Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
 Calculates the boundaries for the frequentist Feldman-Cousins interval. More...
 
static Bool_t FeldmanCousinsInterval (Double_t total, Double_t passed, Double_t level, Double_t &lower, Double_t &upper)
 Calculates the interval boundaries using the frequentist methods of Feldman-Cousins. More...
 
static Double_t MidPInterval (Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
 Calculates the boundaries using the mid-P binomial interval (Lancaster method) from B. More...
 
static Double_t Normal (Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
 Returns the confidence limits for the efficiency supposing that the efficiency follows a normal distribution with the rms below. More...
 
static Double_t Wilson (Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
 Calculates the boundaries for the frequentist Wilson interval. More...
 
- Static Public Member Functions inherited from TObject
static Long_t GetDtorOnly ()
 Return destructor only flag. More...
 
static Bool_t GetObjectStat ()
 Get status of object stat flag. More...
 
static void SetDtorOnly (void *obj)
 Set destructor only flag. More...
 
static void SetObjectStat (Bool_t stat)
 Turn on/off tracking of objects in the TObjectTable. More...
 
- Static Public Member Functions inherited from TAttMarker
static Width_t GetMarkerLineWidth (Style_t style)
 Internal helper function that returns the line width of the given marker style (0 = filled marker) More...
 
static Style_t GetMarkerStyleBase (Style_t style)
 Internal helper function that returns the corresponding marker style with line width 1 for the given style. More...
 

Protected Types

enum  EStatusBits {
  kIsBayesian = BIT(14) , kPosteriorMode = BIT(15) , kShortestInterval = BIT(16) , kUseBinPrior = BIT(17) ,
  kUseWeights = BIT(18)
}
 
- Protected Types inherited from TObject
enum  { kOnlyPrepStep = BIT(3) }
 

Protected Member Functions

void Build (const char *name, const char *title)
 Building standard data structure of a TEfficiency object. More...
 
void FillGraph (TGraphAsymmErrors *graph, Option_t *opt) const
 Fill the graph to be painted with information from TEfficiency Internal method called by TEfficiency::Paint or TEfficiency::CreateGraph. More...
 
void FillHistogram (TH2 *h2) const
 Fill the 2d histogram to be painted with information from TEfficiency 2D Internal method called by TEfficiency::Paint or TEfficiency::CreatePaintingGraph. More...
 
- Protected Member Functions inherited from TObject
virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 Interface to ErrorHandler (protected). More...
 
void MakeZombie ()
 

Protected Attributes

Double_t fBeta_alpha
 
Double_t fBeta_beta
 
std::vector< std::pair< Double_t, Double_t > > fBeta_bin_params
 
Double_t(* fBoundary )(Double_t, Double_t, Double_t, Bool_t)
 
Double_t fConfLevel
 pointer to a method calculating the boundaries of confidence intervals More...
 
TDirectoryfDirectory
 
TListfFunctions
 pointer to directory holding this TEfficiency object More...
 
TGraphAsymmErrorsfPaintGraph
 
TH2fPaintHisto
 temporary graph for painting More...
 
TH1fPassedHistogram
 temporary histogram for painting More...
 
EStatOption fStatisticOption
 
TH1fTotalHistogram
 
Double_t fWeight
 
- Protected Attributes inherited from TNamed
TString fName
 
TString fTitle
 
- Protected Attributes inherited from TAttLine
Color_t fLineColor
 Line color. More...
 
Style_t fLineStyle
 Line style. More...
 
Width_t fLineWidth
 Line width. More...
 
- Protected Attributes inherited from TAttFill
Color_t fFillColor
 Fill area color. More...
 
Style_t fFillStyle
 Fill area style. More...
 
- Protected Attributes inherited from TAttMarker
Color_t fMarkerColor
 Marker color. More...
 
Size_t fMarkerSize
 Marker size. More...
 
Style_t fMarkerStyle
 Marker style. More...
 

#include <TEfficiency.h>

Inheritance diagram for TEfficiency:
[legend]

Member Enumeration Documentation

◆ EStatOption

Enumerator
kFCP 
kFNormal 
kFWilson 
kFAC 
kFFC 
kBJeffrey 
kBUniform 
kBBayesian 
kMidP 

Definition at line 32 of file TEfficiency.h.

◆ EStatusBits

enum TEfficiency::EStatusBits
protected
Enumerator
kIsBayesian 
kPosteriorMode 
kShortestInterval 
kUseBinPrior 
kUseWeights 

Definition at line 61 of file TEfficiency.h.

Constructor & Destructor Documentation

◆ TEfficiency() [1/9]

TEfficiency::TEfficiency ( )

default constructor

should not be used explicitly

Definition at line 619 of file TEfficiency.cxx.

◆ TEfficiency() [2/9]

TEfficiency::TEfficiency ( const TH1 passed,
const TH1 total 
)

constructor using two existing histograms as input

Input: passed - contains the events fulfilling some criteria total - contains all investigated events

Notes: - both histograms have to fulfill the conditions of CheckConsistency

  • dimension of the resulting efficiency object depends on the dimension of the given histograms
  • Clones of both histograms are stored internally
  • The function SetName(total.GetName() + "_clone") is called to set the names of the new object and the internal histograms..
  • The created TEfficiency object is NOT appended to a directory. It will not be written to disk during the next TFile::Write() command in order to prevent duplication of data. If you want to save this TEfficiency object anyway, you can either append it to a directory by calling SetDirectory(TDirectory*) or write it explicitly to disk by calling Write().

Definition at line 658 of file TEfficiency.cxx.

◆ TEfficiency() [3/9]

TEfficiency::TEfficiency ( const char *  name,
const char *  title,
Int_t  nbins,
const Double_t xbins 
)

Create 1-dimensional TEfficiency object with variable bin size.

Constructor creates two new and empty histograms with a given binning

Input:

  • name: the common part of the name for both histograms (no blanks) fTotalHistogram has name: name + "_total" fPassedHistogram has name: name + "_passed"
  • title: the common part of the title for both histogram fTotalHistogram has title: title + " (total)" fPassedHistogram has title: title + " (passed)" It is possible to label the axis by passing a title with the following format: "title;xlabel;ylabel".
  • nbins: number of bins on the x-axis
  • xbins: array of length (nbins + 1) with low-edges for each bin xbins[nbinsx] ... lower edge for overflow bin

Definition at line 724 of file TEfficiency.cxx.

◆ TEfficiency() [4/9]

TEfficiency::TEfficiency ( const char *  name,
const char *  title,
Int_t  nbinsx,
Double_t  xlow,
Double_t  xup 
)

Create 1-dimensional TEfficiency object with fixed bins size.

Constructor creates two new and empty histograms with a fixed binning.

Input:

  • name: the common part of the name for both histograms(no blanks) fTotalHistogram has name: name + "_total" fPassedHistogram has name: name + "_passed"
  • title: the common part of the title for both histogram fTotalHistogram has title: title + " (total)" fPassedHistogram has title: title + " (passed)" It is possible to label the axis by passing a title with the following format: "title;xlabel;ylabel".
  • nbinsx: number of bins on the x-axis
  • xlow: lower edge of first bin
  • xup: upper edge of last bin

Definition at line 763 of file TEfficiency.cxx.

◆ TEfficiency() [5/9]

TEfficiency::TEfficiency ( const char *  name,
const char *  title,
Int_t  nbinsx,
Double_t  xlow,
Double_t  xup,
Int_t  nbinsy,
Double_t  ylow,
Double_t  yup 
)

Create 2-dimensional TEfficiency object with fixed bin size.

Constructor creates two new and empty histograms with a fixed binning.

Input:

  • name: the common part of the name for both histograms(no blanks) fTotalHistogram has name: name + "_total" fPassedHistogram has name: name + "_passed"
  • title: the common part of the title for both histogram fTotalHistogram has title: title + " (total)" fPassedHistogram has title: title + " (passed)" It is possible to label the axis by passing a title with the following format: "title;xlabel;ylabel;zlabel".
  • nbinsx: number of bins on the x-axis
  • xlow: lower edge of first x-bin
  • xup: upper edge of last x-bin
  • nbinsy: number of bins on the y-axis
  • ylow: lower edge of first y-bin
  • yup: upper edge of last y-bin

Definition at line 805 of file TEfficiency.cxx.

◆ TEfficiency() [6/9]

TEfficiency::TEfficiency ( const char *  name,
const char *  title,
Int_t  nbinsx,
const Double_t xbins,
Int_t  nbinsy,
const Double_t ybins 
)

Create 2-dimensional TEfficiency object with variable bin size.

Constructor creates two new and empty histograms with a given binning.

Input:

  • name: the common part of the name for both histograms(no blanks) fTotalHistogram has name: name + "_total" fPassedHistogram has name: name + "_passed"
  • title: the common part of the title for both histogram fTotalHistogram has title: title + " (total)" fPassedHistogram has title: title + " (passed)" It is possible to label the axis by passing a title with the following format: "title;xlabel;ylabel;zlabel".
  • nbinsx: number of bins on the x-axis
  • xbins: array of length (nbins + 1) with low-edges for each bin xbins[nbinsx] ... lower edge for overflow x-bin
  • nbinsy: number of bins on the y-axis
  • ybins: array of length (nbins + 1) with low-edges for each bin ybins[nbinsy] ... lower edge for overflow y-bin

Definition at line 848 of file TEfficiency.cxx.

◆ TEfficiency() [7/9]

TEfficiency::TEfficiency ( const char *  name,
const char *  title,
Int_t  nbinsx,
Double_t  xlow,
Double_t  xup,
Int_t  nbinsy,
Double_t  ylow,
Double_t  yup,
Int_t  nbinsz,
Double_t  zlow,
Double_t  zup 
)

Create 3-dimensional TEfficiency object with fixed bin size.

Constructor creates two new and empty histograms with a fixed binning.

Input:

  • name: the common part of the name for both histograms(no blanks) fTotalHistogram has name: name + "_total" fPassedHistogram has name: name + "_passed"
  • title: the common part of the title for both histogram fTotalHistogram has title: title + " (total)" fPassedHistogram has title: title + " (passed)" It is possible to label the axis by passing a title with the following format: "title;xlabel;ylabel;zlabel".
  • nbinsx: number of bins on the x-axis
  • xlow: lower edge of first x-bin
  • xup: upper edge of last x-bin
  • nbinsy: number of bins on the y-axis
  • ylow: lower edge of first y-bin
  • yup: upper edge of last y-bin
  • nbinsz: number of bins on the z-axis
  • zlow: lower edge of first z-bin
  • zup: upper edge of last z-bin

Definition at line 894 of file TEfficiency.cxx.

◆ TEfficiency() [8/9]

TEfficiency::TEfficiency ( const char *  name,
const char *  title,
Int_t  nbinsx,
const Double_t xbins,
Int_t  nbinsy,
const Double_t ybins,
Int_t  nbinsz,
const Double_t zbins 
)

Create 3-dimensional TEfficiency object with variable bin size.

Constructor creates two new and empty histograms with a given binning.

Input:

  • name: the common part of the name for both histograms(no blanks) fTotalHistogram has name: name + "_total" fPassedHistogram has name: name + "_passed"
  • title: the common part of the title for both histogram fTotalHistogram has title: title + " (total)" fPassedHistogram has title: title + " (passed)" It is possible to label the axis by passing a title with the following format: "title;xlabel;ylabel;zlabel".
  • nbinsx: number of bins on the x-axis
  • xbins: array of length (nbins + 1) with low-edges for each bin xbins[nbinsx] ... lower edge for overflow x-bin
  • nbinsy: number of bins on the y-axis
  • ybins: array of length (nbins + 1) with low-edges for each bin xbins[nbinsx] ... lower edge for overflow y-bin
  • nbinsz: number of bins on the z-axis
  • zbins: array of length (nbins + 1) with low-edges for each bin xbins[nbinsx] ... lower edge for overflow z-bin

Definition at line 941 of file TEfficiency.cxx.

◆ TEfficiency() [9/9]

TEfficiency::TEfficiency ( const TEfficiency rEff)

Copy constructor.

The list of associated objects (e.g. fitted functions) is not copied.

Note:

  • SetName(rEff.GetName() + "_copy") is called to set the names of the object and the histograms.
  • The titles are set by calling SetTitle("[copy] " + rEff.GetTitle()).
  • The copied TEfficiency object is NOT appended to a directory. It will not be written to disk during the next TFile::Write() command in order to prevent duplication of data. If you want to save this TEfficiency object anyway, you can either append it to a directory by calling SetDirectory(TDirectory*) or write it explicitly to disk by calling Write().

Definition at line 980 of file TEfficiency.cxx.

◆ ~TEfficiency()

TEfficiency::~TEfficiency ( )

default destructor

Definition at line 1024 of file TEfficiency.cxx.

Member Function Documentation

◆ Add()

void TEfficiency::Add ( const TEfficiency rEff)
inline

Definition at line 94 of file TEfficiency.h.

◆ AgrestiCoull()

Double_t TEfficiency::AgrestiCoull ( Double_t  total,
Double_t  passed,
Double_t  level,
Bool_t  bUpper 
)
static

Calculates the boundaries for the frequentist Agresti-Coull interval.

Parameters
totalnumber of total events
passed0 <= number of passed events <= total
levelconfidence level
bUppertrue - upper boundary is returned false - lower boundary is returned

\begin{eqnarray*} \alpha &=& 1 - \frac{level}{2} \\ \kappa &=& \Phi^{-1}(1 - \alpha,1)\ ... normal\ quantile\ function\\ mode &=& \frac{passed + \frac{\kappa^{2}}{2}}{total + \kappa^{2}}\\ \Delta &=& \kappa * \sqrt{\frac{mode * (1 - mode)}{total + \kappa^{2}}}\\ return &=& max(0,mode - \Delta)\ or\ min(1,mode + \Delta) \end{eqnarray*}

Definition at line 1074 of file TEfficiency.cxx.

◆ Bayesian()

Double_t TEfficiency::Bayesian ( Double_t  total,
Double_t  passed,
Double_t  level,
Double_t  alpha,
Double_t  beta,
Bool_t  bUpper,
Bool_t  bShortest = false 
)
static

Calculates the boundaries for a Bayesian confidence interval (shortest or central interval depending on the option)

Parameters
[in]totalnumber of total events
[in]passed0 <= number of passed events <= total
[in]levelconfidence level
[in]alphashape parameter > 0 for the prior distribution (fBeta_alpha)
[in]betashape parameter > 0 for the prior distribution (fBeta_beta)
[in]bUpper
  • true - upper boundary is returned
  • false - lower boundary is returned
[in]bShortest??

Note: In the case central confidence interval is calculated. when passed = 0 (or passed = total) the lower (or upper) interval values will be larger than 0 (or smaller than 1).

Calculation:

The posterior probability in bayesian statistics is given by:

\[ P(\varepsilon |k,N) \propto L(\varepsilon|k,N) \times Prior(\varepsilon) \]

As an efficiency can be interpreted as probability of a positive outcome of a Bernoullli trial the likelihood function is given by the binomial distribution:

\[ L(\varepsilon|k,N) = Binomial(N,k) \varepsilon ^{k} (1 - \varepsilon)^{N-k} \]

At the moment only beta distributions are supported as prior probabilities of the efficiency ( \( B(\alpha,\beta)\) is the beta function):

\[ Prior(\varepsilon) = \frac{1}{B(\alpha,\beta)} \varepsilon ^{\alpha - 1} (1 - \varepsilon)^{\beta - 1} \]

The posterior probability is therefore again given by a beta distribution:

\[ P(\varepsilon |k,N) \propto \varepsilon ^{k + \alpha - 1} (1 - \varepsilon)^{N - k + \beta - 1} \]

In case of central intervals the lower boundary for the equal-tailed confidence interval is given by the inverse cumulative (= quantile) function for the quantile \( \frac{1 - level}{2} \). The upper boundary for the equal-tailed confidence interval is given by the inverse cumulative (= quantile) function for the quantile \( \frac{1 + level}{2} \). Hence it is the solution \( \varepsilon \) of the following equation:

\[ I_{\varepsilon}(k + \alpha,N - k + \beta) = \frac{1}{norm} \int_{0}^{\varepsilon} dt t^{k + \alpha - 1} (1 - t)^{N - k + \beta - 1} = \frac{1 \pm level}{2} \]

In the case of shortest interval the minimum interval around the mode is found by minimizing the length of all intervals width the given probability content. See TEfficiency::BetaShortestInterval

Definition at line 1249 of file TEfficiency.cxx.

◆ BetaCentralInterval()

Double_t TEfficiency::BetaCentralInterval ( Double_t  level,
Double_t  a,
Double_t  b,
Bool_t  bUpper 
)
static

Calculates the boundaries for a central confidence interval for a Beta distribution.

Parameters
[in]levelconfidence level
[in]aparameter > 0 for the beta distribution (for a posterior is passed + prior_alpha
[in]bparameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta
[in]bUppertrue - upper boundary is returned false - lower boundary is returned

Definition at line 1273 of file TEfficiency.cxx.

◆ BetaMean()

Double_t TEfficiency::BetaMean ( Double_t  a,
Double_t  b 
)
static

Compute the mean (average) of the beta distribution.

Parameters
[in]aparameter > 0 for the beta distribution (for a posterior is passed + prior_alpha
[in]bparameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta

Definition at line 1383 of file TEfficiency.cxx.

◆ BetaMode()

Double_t TEfficiency::BetaMode ( Double_t  a,
Double_t  b 
)
static

Compute the mode of the beta distribution.

Parameters
[in]aparameter > 0 for the beta distribution (for a posterior is passed + prior_alpha
[in]bparameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta

note the mode is defined for a Beta(a,b) only if (a,b)>1 (a = passed+alpha; b = total-passed+beta) return then the following in case (a,b) < 1:

  • if (a==b) return 0.5 (it is really undefined)
  • if (a < b) return 0;
  • if (a > b) return 1;

Definition at line 1406 of file TEfficiency.cxx.

◆ BetaShortestInterval()

Bool_t TEfficiency::BetaShortestInterval ( Double_t  level,
Double_t  a,
Double_t  b,
Double_t lower,
Double_t upper 
)
static

Calculates the boundaries for a shortest confidence interval for a Beta distribution.

Parameters
[in]levelconfidence level
[in]aparameter > 0 for the beta distribution (for a posterior is passed + prior_alpha
[in]bparameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta
[out]upperupper boundary is returned
[out]lowerlower boundary is returned

The lower/upper boundary are then obtained by finding the shortest interval of the beta distribution contained the desired probability level. The length of all possible intervals is minimized in order to find the shortest one

Definition at line 1329 of file TEfficiency.cxx.

◆ Browse()

void TEfficiency::Browse ( TBrowser b)
inlinevirtual

Browse object. May be overridden for another default action.

Reimplemented from TObject.

Definition at line 95 of file TEfficiency.h.

◆ Build()

void TEfficiency::Build ( const char *  name,
const char *  title 
)
protected

Building standard data structure of a TEfficiency object.

Notes:

  • calls: SetName(name), SetTitle(title)
  • set the statistic option to the default (kFCP)
  • appends this object to the current directory SetDirectory(gDirectory)

Definition at line 1430 of file TEfficiency.cxx.

◆ CheckBinning()

Bool_t TEfficiency::CheckBinning ( const TH1 pass,
const TH1 total 
)
static

Checks binning for each axis.

It is assumed that the passed histograms have the same dimension.

Definition at line 1452 of file TEfficiency.cxx.

◆ CheckConsistency()

Bool_t TEfficiency::CheckConsistency ( const TH1 pass,
const TH1 total,
Option_t opt = "" 
)
static

Checks the consistence of the given histograms.

The histograms are considered as consistent if:

  • both have the same dimension
  • both have the same binning
  • pass.GetBinContent(i) <= total.GetBinContent(i) for each bin i

Definition at line 1501 of file TEfficiency.cxx.

◆ CheckEntries()

Bool_t TEfficiency::CheckEntries ( const TH1 pass,
const TH1 total,
Option_t opt = "" 
)
static

Checks whether bin contents are compatible with binomial statistics.

The following inequality has to be valid for each bin i: total.GetBinContent(i) >= pass.GetBinContent(i)

Note:

  • It is assumed that both histograms have the same dimension and binning.

Definition at line 1533 of file TEfficiency.cxx.

◆ CheckWeights()

Bool_t TEfficiency::CheckWeights ( const TH1 pass,
const TH1 total 
)
static

Check if both histogram are weighted.

If they are weighted a true is returned

Definition at line 1563 of file TEfficiency.cxx.

◆ ClopperPearson()

Double_t TEfficiency::ClopperPearson ( Double_t  total,
Double_t  passed,
Double_t  level,
Bool_t  bUpper 
)
static

Calculates the boundaries for the frequentist Clopper-Pearson interval.

This interval is recommended by the PDG.

Parameters
[in]totalnumber of total events
[in]passed0 <= number of passed events <= total
[in]levelconfidence level
[in]bUppertrue - upper boundary is returned ;false - lower boundary is returned

Calculation:

The lower boundary of the Clopper-Pearson interval is the "exact" inversion of the test:

\begin{eqnarray*} P(x \geq passed; total) &=& \frac{1 - level}{2}\\ P(x \geq passed; total) &=& 1 - P(x \leq passed - 1; total)\\ &=& 1 - \frac{1}{norm} * \int_{0}^{1 - \varepsilon} t^{total - passed} (1 - t)^{passed - 1} dt\\ &=& 1 - \frac{1}{norm} * \int_{\varepsilon}^{1} t^{passed - 1} (1 - t)^{total - passed} dt\\ &=& \frac{1}{norm} * \int_{0}^{\varepsilon} t^{passed - 1} (1 - t)^{total - passed} dt\\ &=& I_{\varepsilon}(passed,total - passed + 1) \end{eqnarray*}

The lower boundary is therefore given by the \( \frac{1 - level}{2}\) quantile of the beta distribution.

The upper boundary of the Clopper-Pearson interval is the "exact" inversion of the test:

\begin{eqnarray*} P(x \leq passed; total) &=& \frac{1 - level}{2}\\ P(x \leq passed; total) &=& \frac{1}{norm} * \int_{0}^{1 - \varepsilon} t^{total - passed - 1} (1 - t)^{passed} dt\\ &=& \frac{1}{norm} * \int_{\varepsilon}^{1} t^{passed} (1 - t)^{total - passed - 1} dt\\ &=& 1 - \frac{1}{norm} * \int_{0}^{\varepsilon} t^{passed} (1 - t)^{total - passed - 1} dt\\ \Rightarrow 1 - \frac{1 - level}{2} &=& \frac{1}{norm} * \int_{0}^{\varepsilon} t^{passed} (1 - t)^{total - passed -1} dt\\ \frac{1 + level}{2} &=& I_{\varepsilon}(passed + 1,total - passed) \end{eqnarray*}

The upper boundary is therefore given by the \(\frac{1 + level}{2}\) quantile of the beta distribution.

Note: The connection between the binomial distribution and the regularized incomplete beta function \( I_{\varepsilon}(\alpha,\beta)\) has been used.

Definition at line 1821 of file TEfficiency.cxx.

◆ Combine() [1/2]

Double_t TEfficiency::Combine ( Double_t up,
Double_t low,
Int_t  n,
const Int_t pass,
const Int_t total,
Double_t  alpha,
Double_t  beta,
Double_t  level = 0.683,
const Double_t w = 0,
Option_t opt = "" 
)
static
Calculates the combined efficiency and its uncertainties

This method does a bayesian combination of the given samples.

\param[in] up  contains the upper limit of the confidence interval afterwards
\param[in] low  contains the lower limit of the confidence interval afterwards
\param[in] n    number of samples which are combined
\param[in] pass array of length n containing the number of passed events
\param[in] total array of length n containing the corresponding numbers of total events
\param[in] alpha  shape parameters for the beta distribution as prior
\param[in] beta   shape parameters for the beta distribution as prior
\param[in] level  desired confidence level
\param[in] w weights for each sample; if not given, all samples get the weight 1
          The weights do not need to be normalized, since they are internally renormalized
          to the number of effective entries.
\param[in] opt
  -  mode : The mode is returned instead of the mean of the posterior as best value
            When using the mode the shortest interval is also computed instead of the central one
  -  shortest: compute shortest interval (done by default if mode option is set)
  -  central: compute central interval (done by default if mode option is NOT set)

Calculation:

The combined posterior distributions is calculated from the Bayes theorem assuming a common prior Beta distribution.
    It is easy to proof that the combined posterior is then:

\begin{eqnarray*} P_{comb}(\epsilon |{w_{i}}; {k_{i}}; {N_{i}}) &=& B(\epsilon, \sum_{i}{ w_{i} k_{i}} + \alpha, \sum_{i}{ w_{i}(n_{i}-k_{i})}+\beta)\\ w_{i} &=& weight\ for\ each\ sample\ renormalized\ to\ the\ effective\ entries\\ w^{'}_{i} &=& w_{i} \frac{ \sum_{i} {w_{i} } } { \sum_{i} {w_{i}^{2} } } \end{eqnarray*}

The estimated efficiency is the mode (or the mean) of the obtained posterior distribution

The boundaries of the confidence interval for a confidence level (1 - a) are given by the a/2 and 1-a/2 quantiles of the resulting cumulative distribution.

Example (uniform prior distribution):

{
TCanvas* c1 = new TCanvas("c1","",600,800);
c1->Divide(1,2);
c1->SetFillStyle(1001);
c1->SetFillColor(kWhite);
TF1* p1 = new TF1("p1","TMath::BetaDist(x,19,9)",0,1);
TF1* p2 = new TF1("p2","TMath::BetaDist(x,4,8)",0,1);
TF1* comb = new TF1("comb2","TMath::BetaDist(x,[0],[1])",0,1);
double nrm = 1./(0.6*0.6+0.4*0.4); // weight normalization
double a = 0.6*18.0 + 0.4*3.0 + 1.0; // new alpha parameter of combined beta dist.
double b = 0.6*10+0.4*7+1.0; // new beta parameter of combined beta dist.
comb->SetParameters(nrm*a ,nrm *b );
TF1* const1 = new TF1("const1","0.05",0,1);
TF1* const2 = new TF1("const2","0.95",0,1);
p1->SetTitle("combined posteriors;#epsilon;P(#epsilon|k,N)");
comb->SetLineColor(kGreen+2);
TLegend* leg1 = new TLegend(0.12,0.65,0.5,0.85);
leg1->AddEntry(p1,"k1 = 18, N1 = 26","l");
leg1->AddEntry(p2,"k2 = 3, N2 = 10","l");
leg1->AddEntry(comb,"combined: p1 = 0.6, p2=0.4","l");
c1->cd(1);
comb->Draw();
p1->Draw("same");
p2->Draw("same");
leg1->Draw("same");
c1->cd(2);
const1->SetLineWidth(1);
const2->SetLineWidth(1);
TGraph* gr = (TGraph*)comb->DrawIntegral();
gr->SetTitle("cumulative function of combined posterior with boundaries for cl = 95%;#epsilon;CDF");
const1->Draw("same");
const2->Draw("same");
c1->cd(0);
return c1;
}
#define b(i)
Definition: RSha256.hxx:100
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
virtual TObject * DrawIntegral(Option_t *option="al")
Draw integral of this function.
Definition: TF1.cxx:1397
A TGraph is an object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
virtual void SetTitle(const char *title="")
Change (i.e.
Definition: TGraph.cxx:2324
TGraphErrors * gr
Definition: legend1.C:25
auto * a
Definition: textangle.C:12

Definition at line 1917 of file TEfficiency.cxx.

◆ Combine() [2/2]

TGraphAsymmErrors * TEfficiency::Combine ( TCollection pList,
Option_t option = "",
Int_t  n = 0,
const Double_t w = 0 
)
static

Combines a list of 1-dimensional TEfficiency objects.

A TGraphAsymmErrors object is returned which contains the estimated efficiency and its uncertainty for each bin. If the combination fails, a zero pointer is returned.

At the moment the combining is only implemented for bayesian statistics.

Parameters
[in]pListlist containing TEfficiency objects which should be combined only one-dimensional efficiencies are taken into account
[in]option
  • s : strict combining; only TEfficiency objects with the same beta prior and the flag kIsBayesian == true are combined If not specified the prior parameter of the first TEfficiency object is used
  • v : verbose mode; print information about combining
  • cl=x : set confidence level (0 < cl < 1). If not specified, the confidence level of the first TEfficiency object is used.
  • mode Use mode of combined posterior as estimated value for the efficiency
  • shortest: compute shortest interval (done by default if mode option is set)
  • central: compute central interval (done by default if mode option is NOT set)
[in]nnumber of weights (has to be the number of one-dimensional TEfficiency objects in pList) If no weights are passed, the internal weights GetWeight() of the given TEfficiency objects are used.
[in]warray of length n with weights for each TEfficiency object in pList (w[0] correspond to pList->First ... w[n-1] -> pList->Last) The weights do not have to be normalised.

For each bin the calculation is done by the Combine(double&, double& ...) method.

Definition at line 2007 of file TEfficiency.cxx.

◆ CreateGraph()

TGraphAsymmErrors * TEfficiency::CreateGraph ( Option_t opt = "") const

Create the graph used be painted (for dim=1 TEfficiency) The return object is managed by the caller.

Definition at line 1592 of file TEfficiency.cxx.

◆ CreateHistogram()

TH2 * TEfficiency::CreateHistogram ( Option_t opt = "") const

Create the histogram used to be painted (for dim=2 TEfficiency) The return object is managed by the caller.

Definition at line 1697 of file TEfficiency.cxx.

◆ DistancetoPrimitive()

Int_t TEfficiency::DistancetoPrimitive ( Int_t  px,
Int_t  py 
)
virtual

Compute distance from point px,py to a graph.

Compute the closest distance of approach from point px,py to this line. The distance is computed in pixels units.

Forward the call to the painted graph

Reimplemented from TObject.

Definition at line 2192 of file TEfficiency.cxx.

◆ Draw()

void TEfficiency::Draw ( Option_t opt = "")
virtual

Draws the current TEfficiency object.

Parameters
[in]opt

Specific TEfficiency drawing options:

  • E0 - plot bins where the total number of passed events is zero (the error interval will be [0,1] )

Reimplemented from TObject.

Definition at line 2213 of file TEfficiency.cxx.

◆ ExecuteEvent()

void TEfficiency::ExecuteEvent ( Int_t  event,
Int_t  px,
Int_t  py 
)
virtual

Execute action corresponding to one event.

This member function is called when the drawn class is clicked with the locator If Left button clicked on one of the line end points, this point follows the cursor until button is released.

if Middle button clicked, the line is moved parallel to itself until the button is released. Forward the call to the underlying graph

Reimplemented from TObject.

Definition at line 2246 of file TEfficiency.cxx.

◆ FeldmanCousins()

Double_t TEfficiency::FeldmanCousins ( Double_t  total,
Double_t  passed,
Double_t  level,
Bool_t  bUpper 
)
static

Calculates the boundaries for the frequentist Feldman-Cousins interval.

Parameters
totalnumber of total events
passed0 <= number of passed events <= total
levelconfidence level
bUppertrue - upper boundary is returned false - lower boundary is returned

Definition at line 1097 of file TEfficiency.cxx.

◆ FeldmanCousinsInterval()

Bool_t TEfficiency::FeldmanCousinsInterval ( Double_t  total,
Double_t  passed,
Double_t  level,
Double_t lower,
Double_t upper 
)
static

Calculates the interval boundaries using the frequentist methods of Feldman-Cousins.

Parameters
[in]totalnumber of total events
[in]passed0 <= number of passed events <= total
[in]levelconfidence level
[out]lowerlower boundary returned on exit
[out]upperlower boundary returned on exit
Returns
a flag with the status of the calculation

Calculation:

The Feldman-Cousins is a frequentist method where the interval is estimated using a Neyman construction where the ordering is based on the likelihood ratio:

\[ LR = \frac{Binomial(k | N, \epsilon)}{Binomial(k | N, \hat{\epsilon} ) } \]

See G. J. Feldman and R. D. Cousins, Phys. Rev. D57 (1998) 3873 and R. D. Cousins, K. E. Hymes, J. Tucker, Nuclear Instruments and Methods in Physics Research A 612 (2010) 388

Implemented using classes developed by Jordan Tucker and Luca Lista See File hist/hist/src/TEfficiencyHelper.h

Definition at line 1130 of file TEfficiency.cxx.

◆ Fill()

void TEfficiency::Fill ( Bool_t  bPassed,
Double_t  x,
Double_t  y = 0,
Double_t  z = 0 
)

This function is used for filling the two histograms.

Parameters
[in]bPassedflag whether the current event passed the selection
  • true: both histograms are filled
  • false: only the total histogram is filled
[in]xx-value
[in]yy-value (use default=0 for 1-D efficiencies)
[in]zz-value (use default=0 for 2-D or 1-D efficiencies)

Definition at line 2262 of file TEfficiency.cxx.

◆ FillGraph()

void TEfficiency::FillGraph ( TGraphAsymmErrors graph,
Option_t opt 
) const
protected

Fill the graph to be painted with information from TEfficiency Internal method called by TEfficiency::Paint or TEfficiency::CreateGraph.

Definition at line 1613 of file TEfficiency.cxx.

◆ FillHistogram()

void TEfficiency::FillHistogram ( TH2 h2) const
protected

Fill the 2d histogram to be painted with information from TEfficiency 2D Internal method called by TEfficiency::Paint or TEfficiency::CreatePaintingGraph.

Definition at line 1735 of file TEfficiency.cxx.

◆ FillWeighted()

void TEfficiency::FillWeighted ( Bool_t  bPassed,
Double_t  weight,
Double_t  x,
Double_t  y = 0,
Double_t  z = 0 
)

This function is used for filling the two histograms with a weight.

Parameters
[in]bPassedflag whether the current event passed the selection
  • true: both histograms are filled
  • false: only the total histogram is filled
[in]weightweight for the event
[in]xx-value
[in]yy-value (use default=0 for 1-D efficiencies)
[in]zz-value (use default=0 for 2-D or 1-D efficiencies)

Note: - this function will call SetUseWeightedEvents if it was not called by the user before

Definition at line 2296 of file TEfficiency.cxx.

◆ FindFixBin()

Int_t TEfficiency::FindFixBin ( Double_t  x,
Double_t  y = 0,
Double_t  z = 0 
) const

Returns the global bin number containing the given values.

Note:

  • values which belong to dimensions higher than the current dimension of the TEfficiency object are ignored (i.e. for 1-dimensional efficiencies only the x-value is considered)

Definition at line 2332 of file TEfficiency.cxx.

◆ Fit()

TFitResultPtr TEfficiency::Fit ( TF1 f1,
Option_t opt = "" 
)

Fits the efficiency using the TBinomialEfficiencyFitter class.

The resulting fit function is added to the list of associated functions.

Options:

  • "+": previous fitted functions in the list are kept, by default all functions in the list are deleted
  • for more fitting options see TBinomialEfficiencyFitter::Fit

Definition at line 2356 of file TEfficiency.cxx.

◆ GetBetaAlpha()

Double_t TEfficiency::GetBetaAlpha ( Int_t  bin = -1) const
inline

Definition at line 106 of file TEfficiency.h.

◆ GetBetaBeta()

Double_t TEfficiency::GetBetaBeta ( Int_t  bin = -1) const
inline

Definition at line 107 of file TEfficiency.h.

◆ GetConfidenceLevel()

Double_t TEfficiency::GetConfidenceLevel ( ) const
inline

Definition at line 108 of file TEfficiency.h.

◆ GetCopyPassedHisto()

TH1 * TEfficiency::GetCopyPassedHisto ( ) const

Returns a cloned version of fPassedHistogram.

Notes:

  • The histogram is filled with unit weights. You might want to scale it with the global weight GetWeight().
  • The returned object is owned by the user who has to care about the deletion of the new TH1 object.
  • This histogram is by default NOT attached to the current directory to avoid duplication of data. If you want to store it automatically during the next TFile::Write() command, you have to attach it to the corresponding directory.
TFile* pFile = new TFile("passed.root","update");
TEfficiency* pEff = (TEfficiency*)gDirectory->Get("my_eff");
TH1* copy = pEff->GetCopyPassedHisto();
pFile->Write();
TH1 * GetCopyPassedHisto() const
Returns a cloned version of fPassedHistogram.
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8393

Definition at line 2416 of file TEfficiency.cxx.

◆ GetCopyTotalHisto()

TH1 * TEfficiency::GetCopyTotalHisto ( ) const

Returns a cloned version of fTotalHistogram.

Notes:

  • The histogram is filled with unit weights. You might want to scale it with the global weight GetWeight().
  • The returned object is owned by the user who has to care about the deletion of the new TH1 object.
  • This histogram is by default NOT attached to the current directory to avoid duplication of data. If you want to store it automatically during the next TFile::Write() command, you have to attach it to the corresponding directory.
TFile* pFile = new TFile("total.root","update");
TEfficiency* pEff = (TEfficiency*)gDirectory->Get("my_eff");
TH1* copy = pEff->GetCopyTotalHisto();
pFile->Write();

Definition at line 2447 of file TEfficiency.cxx.

◆ GetDimension()

Int_t TEfficiency::GetDimension ( ) const

returns the dimension of the current TEfficiency object

Definition at line 2460 of file TEfficiency.cxx.

◆ GetDirectory()

TDirectory * TEfficiency::GetDirectory ( ) const
inline

Definition at line 112 of file TEfficiency.h.

◆ GetEfficiency()

Double_t TEfficiency::GetEfficiency ( Int_t  bin) const

Returns the efficiency in the given global bin.

Note:

  • The estimated efficiency depends on the chosen statistic option: for frequentist ones: \( \hat{\varepsilon} = \frac{passed}{total} \) for bayesian ones the expectation value of the resulting posterior distribution is returned: \( \hat{\varepsilon} = \frac{passed + \alpha}{total + \alpha + \beta} \) If the bit kPosteriorMode is set (or the method TEfficiency::UsePosteriorMode() has been called ) the mode (most probable value) of the posterior is returned: \( \hat{\varepsilon} = \frac{passed + \alpha -1}{total + \alpha + \beta -2} \)
    • If the denominator is equal to 0, an efficiency of 0 is returned.
    • When \( passed + \alpha < 1 \) or \( total - passed + \beta < 1 \) the above formula for the mode is not valid. In these cases values the estimated efficiency is 0 or 1.

Definition at line 2482 of file TEfficiency.cxx.

◆ GetEfficiencyErrorLow()

Double_t TEfficiency::GetEfficiencyErrorLow ( Int_t  bin) const

Returns the lower error on the efficiency in the given global bin.

The result depends on the current confidence level fConfLevel and the chosen statistic option fStatisticOption. See SetStatisticOption(Int_t) for more details.

Note: If the histograms are filled with weights, only bayesian methods and the normal approximation are supported.

Definition at line 2533 of file TEfficiency.cxx.

◆ GetEfficiencyErrorUp()

Double_t TEfficiency::GetEfficiencyErrorUp ( Int_t  bin) const

Returns the upper error on the efficiency in the given global bin.

The result depends on the current confidence level fConfLevel and the chosen statistic option fStatisticOption. See SetStatisticOption(Int_t) for more details.

Note: If the histograms are filled with weights, only bayesian methods and the normal approximation are supported.

Definition at line 2613 of file TEfficiency.cxx.

◆ GetGlobalBin()

Int_t TEfficiency::GetGlobalBin ( Int_t  binx,
Int_t  biny = 0,
Int_t  binz = 0 
) const

Returns the global bin number which can be used as argument for the following functions:

  • GetEfficiency(bin), GetEfficiencyErrorLow(bin), GetEfficiencyErrorUp(bin)
  • SetPassedEvents(bin), SetTotalEvents(bin)

see TH1::GetBin() for conventions on numbering bins

Definition at line 2691 of file TEfficiency.cxx.

◆ GetListOfFunctions()

TList * TEfficiency::GetListOfFunctions ( )

Definition at line 2698 of file TEfficiency.cxx.

◆ GetPaintedGraph()

TGraphAsymmErrors * TEfficiency::GetPaintedGraph ( ) const
inline

Definition at line 117 of file TEfficiency.h.

◆ GetPaintedHistogram()

TH2 * TEfficiency::GetPaintedHistogram ( ) const
inline

Definition at line 118 of file TEfficiency.h.

◆ GetPassedHistogram()

const TH1 * TEfficiency::GetPassedHistogram ( ) const
inline

Definition at line 120 of file TEfficiency.h.

◆ GetStatisticOption()

EStatOption TEfficiency::GetStatisticOption ( ) const
inline

Definition at line 121 of file TEfficiency.h.

◆ GetTotalHistogram()

const TH1 * TEfficiency::GetTotalHistogram ( ) const
inline

Definition at line 122 of file TEfficiency.h.

◆ GetWeight()

Double_t TEfficiency::GetWeight ( ) const
inline

Definition at line 123 of file TEfficiency.h.

◆ Merge()

Long64_t TEfficiency::Merge ( TCollection pList)

Merges the TEfficiency objects in the given list to the given TEfficiency object using the operator+=(TEfficiency&)

The merged result is stored in the current object. The statistic options and the confidence level are taken from the current object.

This function should be used when all TEfficiency objects correspond to the same process.

The new weight is set according to: \( \frac{1}{w_{new}} = \sum_{i} \frac{1}{w_{i}} \)

Definition at line 2716 of file TEfficiency.cxx.

◆ MidPInterval()

Double_t TEfficiency::MidPInterval ( Double_t  total,
Double_t  passed,
Double_t  level,
Bool_t  bUpper 
)
static

Calculates the boundaries using the mid-P binomial interval (Lancaster method) from B.

Cousing and J. Tucker. See http://arxiv.org/abs/0905.3831 for a description and references for the method

Modify equal_tailed to get the kind of interval you want. Can also be converted to interval on ratio of poisson means X/Y by the substitutions

X = passed
total = X + Y
lower_poisson = lower/(1 - lower)
upper_poisson = upper/(1 - upper)
static unsigned int total

Definition at line 1155 of file TEfficiency.cxx.

◆ Normal()

Double_t TEfficiency::Normal ( Double_t  total,
Double_t  passed,
Double_t  level,
Bool_t  bUpper 
)
static

Returns the confidence limits for the efficiency supposing that the efficiency follows a normal distribution with the rms below.

Parameters
[in]totalnumber of total events
[in]passed0 <= number of passed events <= total
[in]levelconfidence level
[in]bUpper
  • true - upper boundary is returned
  • false - lower boundary is returned

Calculation:

\begin{eqnarray*} \hat{\varepsilon} &=& \frac{passed}{total}\\ \sigma_{\varepsilon} &=& \sqrt{\frac{\hat{\varepsilon} (1 - \hat{\varepsilon})}{total}}\\ \varepsilon_{low} &=& \hat{\varepsilon} \pm \Phi^{-1}(\frac{level}{2},\sigma_{\varepsilon}) \end{eqnarray*}

Definition at line 2753 of file TEfficiency.cxx.

◆ operator+=()

TEfficiency & TEfficiency::operator+= ( const TEfficiency rhs)

Adds the histograms of another TEfficiency object to current histograms.

The statistic options and the confidence level remain unchanged.

fTotalHistogram += rhs.fTotalHistogram; fPassedHistogram += rhs.fPassedHistogram;

calculates a new weight: current weight of this TEfficiency object = \( w_{1} \) weight of rhs = \( w_{2} \) \( w_{new} = \frac{w_{1} \times w_{2}}{w_{1} + w_{2}} \)

Definition at line 2780 of file TEfficiency.cxx.

◆ operator=()

TEfficiency & TEfficiency::operator= ( const TEfficiency rhs)

Assignment operator.

The histograms, statistic option, confidence level, weight and paint styles of rhs are copied to the this TEfficiency object.

Note: - The list of associated functions is not copied. After this operation the list of associated functions is empty.

Definition at line 2822 of file TEfficiency.cxx.

◆ Paint()

void TEfficiency::Paint ( Option_t opt)
virtual

Paints this TEfficiency object.

For details on the possible option see Draw(Option_t*)

Note for 1D classes In 1D the TEfficiency uses a TGraphAsymmErrors for drawing The TGraph is created only the first time Paint is used. The user can manipulate the TGraph via the method TEfficiency::GetPaintedGraph() The TGraph creates behing an histogram for the axis. The histogram is created also only the first time. If the axis needs to be updated because in the meantime the class changed use this trick which will trigger a re-calculation of the axis of the graph TEfficiency::GetPaintedGraph()->Set(0)

Note that in order to access the painted graph via GetPaintedGraph() you need either to call Paint or better gPad->Update();

Reimplemented from TObject.

Definition at line 2880 of file TEfficiency.cxx.

◆ SavePrimitive()

void TEfficiency::SavePrimitive ( std::ostream &  out,
Option_t opt = "" 
)
virtual

Have histograms fixed bins along each axis?

Reimplemented from TObject.

Definition at line 2937 of file TEfficiency.cxx.

◆ SetBetaAlpha()

void TEfficiency::SetBetaAlpha ( Double_t  alpha)

Sets the shape parameter α.

The prior probability of the efficiency is given by the beta distribution:

\[ f(\varepsilon;\alpha;\beta) = \frac{1}{B(\alpha,\beta)} \varepsilon^{\alpha-1} (1 - \varepsilon)^{\beta-1} \]

Note: - both shape parameters have to be positive (i.e. > 0)

Definition at line 3111 of file TEfficiency.cxx.

◆ SetBetaBeta()

void TEfficiency::SetBetaBeta ( Double_t  beta)

Sets the shape parameter β.

The prior probability of the efficiency is given by the beta distribution:

\[ f(\varepsilon;\alpha,\beta) = \frac{1}{B(\alpha,\beta)} \varepsilon^{\alpha-1} (1 - \varepsilon)^{\beta-1} \]

Note: - both shape parameters have to be positive (i.e. > 0)

Definition at line 3129 of file TEfficiency.cxx.

◆ SetBetaBinParameters()

void TEfficiency::SetBetaBinParameters ( Int_t  bin,
Double_t  alpha,
Double_t  beta 
)

Sets different shape parameter α and β for the prior distribution for each bin.

By default the global parameter are used if they are not set for the specific bin The prior probability of the efficiency is given by the beta distribution:

\[ f(\varepsilon;\alpha;\beta) = \frac{1}{B(\alpha,\beta)} \varepsilon^{\alpha-1} (1 - \varepsilon)^{\beta-1} \]

Note:

  • both shape parameters have to be positive (i.e. > 0)
  • bin gives the global bin number (cf. GetGlobalBin)

Definition at line 3150 of file TEfficiency.cxx.

◆ SetBins() [1/6]

Bool_t TEfficiency::SetBins ( Int_t  nx,
const Double_t xBins 
)

Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost.

Definition at line 3191 of file TEfficiency.cxx.

◆ SetBins() [2/6]

Bool_t TEfficiency::SetBins ( Int_t  nx,
const Double_t xBins,
Int_t  ny,
const Double_t yBins 
)

Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost.

Definition at line 3231 of file TEfficiency.cxx.

◆ SetBins() [3/6]

Bool_t TEfficiency::SetBins ( Int_t  nx,
const Double_t xBins,
Int_t  ny,
const Double_t yBins,
Int_t  nz,
const Double_t zBins 
)

Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost.

Definition at line 3272 of file TEfficiency.cxx.

◆ SetBins() [4/6]

Bool_t TEfficiency::SetBins ( Int_t  nx,
Double_t  xmin,
Double_t  xmax 
)

Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost.

Definition at line 3171 of file TEfficiency.cxx.

◆ SetBins() [5/6]

Bool_t TEfficiency::SetBins ( Int_t  nx,
Double_t  xmin,
Double_t  xmax,
Int_t  ny,
Double_t  ymin,
Double_t  ymax 
)

Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost.

Definition at line 3211 of file TEfficiency.cxx.

◆ SetBins() [6/6]

Bool_t TEfficiency::SetBins ( Int_t  nx,
Double_t  xmin,
Double_t  xmax,
Int_t  ny,
Double_t  ymin,
Double_t  ymax,
Int_t  nz,
Double_t  zmin,
Double_t  zmax 
)

Set the bins for the underlined passed and total histograms If the class have been already filled the previous contents will be lost.

Definition at line 3251 of file TEfficiency.cxx.

◆ SetCentralInterval()

void TEfficiency::SetCentralInterval ( Bool_t  on = true)
inline

Definition at line 140 of file TEfficiency.h.

◆ SetConfidenceLevel()

void TEfficiency::SetConfidenceLevel ( Double_t  level)

Sets the confidence level (0 < level < 1) The default value is 1-sigma :~ 0.683.

Definition at line 3293 of file TEfficiency.cxx.

◆ SetDirectory()

void TEfficiency::SetDirectory ( TDirectory dir)

Sets the directory holding this TEfficiency object.

A reference to this TEfficiency object is removed from the current directory (if it exists) and a new reference to this TEfficiency object is added to the given directory.

Notes: - If the given directory is 0, the TEfficiency object does not belong to any directory and will not be written to file during the next TFile::Write() command.

Definition at line 3312 of file TEfficiency.cxx.

◆ SetName()

void TEfficiency::SetName ( const char *  name)
virtual

Sets the name.

Note: The names of the internal histograms are set to "name + _total" and "name + _passed" respectively.

Reimplemented from TNamed.

Definition at line 3329 of file TEfficiency.cxx.

◆ SetPassedEvents()

Bool_t TEfficiency::SetPassedEvents ( Int_t  bin,
Int_t  events 
)

Sets the number of passed events in the given global bin.

returns "true" if the number of passed events has been updated otherwise "false" ist returned

Note: - requires: 0 <= events <= fTotalHistogram->GetBinContent(bin)

Definition at line 3348 of file TEfficiency.cxx.

◆ SetPassedHistogram()

Bool_t TEfficiency::SetPassedHistogram ( const TH1 rPassed,
Option_t opt 
)

Sets the histogram containing the passed events.

The given histogram is cloned and stored internally as histogram containing the passed events. The given histogram has to be consistent with the current fTotalHistogram (see CheckConsistency(const TH1&,const TH1&)). The method returns whether the fPassedHistogram has been replaced (true) or not (false).

Note: The list of associated functions fFunctions is cleared.

Option:

  • "f": force the replacement without checking the consistency This can lead to inconsistent histograms and useless results or unexpected behaviour. But sometimes it might be the only way to change the histograms. If you use this option, you should ensure that the fTotalHistogram is replaced by a consistent one (with respect to rPassed) as well.

Definition at line 3379 of file TEfficiency.cxx.

◆ SetPosteriorAverage()

void TEfficiency::SetPosteriorAverage ( Bool_t  on = true)
inline

Definition at line 138 of file TEfficiency.h.

◆ SetPosteriorMode()

void TEfficiency::SetPosteriorMode ( Bool_t  on = true)
inline

Definition at line 137 of file TEfficiency.h.

◆ SetShortestInterval()

void TEfficiency::SetShortestInterval ( Bool_t  on = true)
inline

Definition at line 139 of file TEfficiency.h.

◆ SetStatisticOption()

void TEfficiency::SetStatisticOption ( EStatOption  option)

Sets the statistic option which affects the calculation of the confidence interval.

Options:

  • kFCP (=0)(default): using the Clopper-Pearson interval (recommended by PDG) sets kIsBayesian = false see also ClopperPearson
  • kFNormal (=1) : using the normal approximation sets kIsBayesian = false see also Normal
  • kFWilson (=2) : using the Wilson interval sets kIsBayesian = false see also Wilson
  • kFAC (=3) : using the Agresti-Coull interval sets kIsBayesian = false see also AgrestiCoull
  • kFFC (=4) : using the Feldman-Cousins frequentist method sets kIsBayesian = false see also FeldmanCousins
  • kBJeffrey (=5) : using the Jeffrey interval sets kIsBayesian = true, fBeta_alpha = 0.5 and fBeta_beta = 0.5 see also Bayesian
  • kBUniform (=6) : using a uniform prior sets kIsBayesian = true, fBeta_alpha = 1 and fBeta_beta = 1 see also Bayesian
  • kBBayesian (=7) : using a custom prior defined by fBeta_alpha and fBeta_beta sets kIsBayesian = true see also Bayesian
  • kMidP (=8) : using the Lancaster Mid-P method sets kIsBayesian = false

Definition at line 3443 of file TEfficiency.cxx.

◆ SetTitle()

void TEfficiency::SetTitle ( const char *  title)
virtual

Sets the title.

Notes:

  • The titles of the internal histograms are set to "title + (total)" or "title + (passed)" respectively.
  • It is possible to label the axis of the histograms as usual (see TH1::SetTitle).

Example: Setting the title to "My Efficiency" and label the axis pEff->SetTitle("My Efficiency;x label;eff");

Reimplemented from TNamed.

Definition at line 3507 of file TEfficiency.cxx.

◆ SetTotalEvents()

Bool_t TEfficiency::SetTotalEvents ( Int_t  bin,
Int_t  events 
)

Sets the number of total events in the given global bin.

returns "true" if the number of total events has been updated otherwise "false" ist returned

Note: - requires: fPassedHistogram->GetBinContent(bin) <= events

Definition at line 3541 of file TEfficiency.cxx.

◆ SetTotalHistogram()

Bool_t TEfficiency::SetTotalHistogram ( const TH1 rTotal,
Option_t opt 
)

Sets the histogram containing all events.

The given histogram is cloned and stored internally as histogram containing all events. The given histogram has to be consistent with the current fPassedHistogram (see CheckConsistency(const TH1&,const TH1&)). The method returns whether the fTotalHistogram has been replaced (true) or not (false).

Note: The list of associated functions fFunctions is cleared.

Option:

  • "f": force the replacement without checking the consistency This can lead to inconsistent histograms and useless results or unexpected behaviour. But sometimes it might be the only way to change the histograms. If you use this option, you should ensure that the fPassedHistogram is replaced by a consistent one (with respect to rTotal) as well.

Definition at line 3572 of file TEfficiency.cxx.

◆ SetUseWeightedEvents()

void TEfficiency::SetUseWeightedEvents ( Bool_t  on = kTRUE)

Definition at line 3605 of file TEfficiency.cxx.

◆ SetWeight()

void TEfficiency::SetWeight ( Double_t  weight)

Sets the global weight for this TEfficiency object.

Note: - weight has to be positive ( > 0)

Definition at line 3623 of file TEfficiency.cxx.

◆ UsesBayesianStat()

Bool_t TEfficiency::UsesBayesianStat ( ) const
inline

Definition at line 156 of file TEfficiency.h.

◆ UsesCentralInterval()

Bool_t TEfficiency::UsesCentralInterval ( ) const
inline

Definition at line 160 of file TEfficiency.h.

◆ UsesPosteriorAverage()

Bool_t TEfficiency::UsesPosteriorAverage ( ) const
inline

Definition at line 159 of file TEfficiency.h.

◆ UsesPosteriorMode()

Bool_t TEfficiency::UsesPosteriorMode ( ) const
inline

Definition at line 157 of file TEfficiency.h.

◆ UsesShortestInterval()

Bool_t TEfficiency::UsesShortestInterval ( ) const
inline

Definition at line 158 of file TEfficiency.h.

◆ UsesWeights()

Bool_t TEfficiency::UsesWeights ( ) const
inline

Definition at line 161 of file TEfficiency.h.

◆ Wilson()

Double_t TEfficiency::Wilson ( Double_t  total,
Double_t  passed,
Double_t  level,
Bool_t  bUpper 
)
static

Calculates the boundaries for the frequentist Wilson interval.

Parameters
[in]totalnumber of total events
[in]passed0 <= number of passed events <= total
[in]levelconfidence level
[in]bUpper
  • true - upper boundary is returned
  • false - lower boundary is returned

Calculation:

\begin{eqnarray*} \alpha &=& 1 - \frac{level}{2}\\ \kappa &=& \Phi^{-1}(1 - \alpha,1) ...\ normal\ quantile\ function\\ mode &=& \frac{passed + \frac{\kappa^{2}}{2}}{total + \kappa^{2}}\\ \Delta &=& \frac{\kappa}{total + \kappa^{2}} * \sqrt{passed (1 - \frac{passed}{total}) + \frac{\kappa^{2}}{4}}\\ return &=& max(0,mode - \Delta)\ or\ min(1,mode + \Delta) \end{eqnarray*}

Definition at line 3653 of file TEfficiency.cxx.

Member Data Documentation

◆ fBeta_alpha

Double_t TEfficiency::fBeta_alpha
protected

Definition at line 46 of file TEfficiency.h.

◆ fBeta_beta

Double_t TEfficiency::fBeta_beta
protected

Definition at line 47 of file TEfficiency.h.

◆ fBeta_bin_params

std::vector<std::pair<Double_t, Double_t> > TEfficiency::fBeta_bin_params
protected

Definition at line 48 of file TEfficiency.h.

◆ fBoundary

Double_t(* TEfficiency::fBoundary) (Double_t, Double_t, Double_t, Bool_t)
protected

Definition at line 50 of file TEfficiency.h.

◆ fConfLevel

Double_t TEfficiency::fConfLevel
protected

pointer to a method calculating the boundaries of confidence intervals

Definition at line 51 of file TEfficiency.h.

◆ fDirectory

TDirectory* TEfficiency::fDirectory
protected

Definition at line 52 of file TEfficiency.h.

◆ fFunctions

TList* TEfficiency::fFunctions
protected

pointer to directory holding this TEfficiency object

Definition at line 53 of file TEfficiency.h.

◆ fPaintGraph

TGraphAsymmErrors* TEfficiency::fPaintGraph
protected

Definition at line 54 of file TEfficiency.h.

◆ fPaintHisto

TH2* TEfficiency::fPaintHisto
protected

temporary graph for painting

Definition at line 55 of file TEfficiency.h.

◆ fPassedHistogram

TH1* TEfficiency::fPassedHistogram
protected

temporary histogram for painting

Definition at line 56 of file TEfficiency.h.

◆ fStatisticOption

EStatOption TEfficiency::fStatisticOption
protected

Definition at line 57 of file TEfficiency.h.

◆ fTotalHistogram

TH1* TEfficiency::fTotalHistogram
protected

Definition at line 58 of file TEfficiency.h.

◆ fWeight

Double_t TEfficiency::fWeight
protected

Definition at line 59 of file TEfficiency.h.

Libraries for TEfficiency:
[legend]

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