1 #ifndef ROOT_TEfficiency_cxx 2 #define ROOT_TEfficiency_cxx 683 Info(
"TEfficiency",
"given histograms are filled with weights");
688 Error(
"TEfficiency(const TH1&,const TH1&)",
"histograms are not consistent -> results are useless");
689 Warning(
"TEfficiency(const TH1&,const TH1&)",
"using two empty TH1D('h1','h1',10,0,10)");
909 fTotalHistogram =
new TH3D(
"total",
"total",nbinsx,xlow,xup,nbinsy,ylow,yup,nbinsz,zlow,zup);
910 fPassedHistogram =
new TH3D(
"passed",
"passed",nbinsx,xlow,xup,nbinsy,ylow,yup,nbinsz,zlow,zup);
1016 rEff.TAttLine::Copy(*
this);
1017 rEff.TAttFill::Copy(*
this);
1018 rEff.TAttMarker::Copy(*
this);
1079 Double_t mode = (passed + 0.5 * kappa * kappa) / (total + kappa * kappa);
1080 Double_t delta = kappa *
std::sqrt(mode * (1 - mode) / (total + kappa * kappa));
1083 return ((mode + delta) > 1) ? 1.0 : (mode + delta);
1085 return ((mode - delta) < 0) ? 0.0 : (mode - delta);
1102 ::Error(
"FeldmanCousins",
"Error running FC method - return 0 or 1");
1104 return (bUpper) ? upper : lower;
1133 double alpha = 1.-level;
1157 const double alpha = 1. - level;
1158 const bool equal_tailed =
true;
1159 const double alpha_min = equal_tailed ? alpha/2 : alpha;
1160 const double tol = 1
e-9;
1170 if ( passed > 0 && passed < 1) {
1173 p = (p1 - p0) * passed + p0;
1177 while (std::abs(pmax - pmin) > tol) {
1178 p = (pmin + pmax)/2;
1186 double vmin = (bUpper) ? alpha_min : 1.- alpha_min;
1258 return (bUpper) ? upper : lower;
1276 if((a > 0) && (b > 0))
1279 gROOT->Error(
"TEfficiency::BayesianCentral",
"Invalid input parameters - return 1");
1284 if((a > 0) && (b > 0))
1287 gROOT->Error(
"TEfficiency::BayesianCentral",
"Invalid input parameters - return 0");
1293 struct Beta_interval_length {
1295 fCL(level), fAlpha(alpha), fBeta(beta)
1331 if (a <= 0 || b <= 0) {
1332 lower = 0; upper = 1;
1333 gROOT->Error(
"TEfficiency::BayesianShortest",
"Invalid input parameters - return [0,1]");
1352 if ( a==b && a<=1.0) {
1360 Beta_interval_length intervalLength(level,a,b);
1364 minim.
SetFunction(func, 0, intervalLength.LowerMax() );
1366 bool ret = minim.
Minimize(100, 1.
E-10,1.
E-10);
1368 gROOT->Error(
"TEfficiency::BayesianShortes",
"Error finding the shortest interval");
1385 if (a <= 0 || b <= 0 ) {
1386 gROOT->Error(
"TEfficiency::BayesianMean",
"Invalid input parameters - return 0");
1408 if (a <= 0 || b <= 0 ) {
1409 gROOT->Error(
"TEfficiency::BayesianMode",
"Invalid input parameters - return 0");
1412 if ( a <= 1 || b <= 1) {
1413 if ( a < b)
return 0;
1414 if ( a > b)
return 1;
1415 if (a == b)
return 0.5;
1419 Double_t mode = (a - 1.0) / (a + b -2.0);
1455 const TAxis* ax1 = 0;
1456 const TAxis* ax2 = 0;
1476 gROOT->Info(
"TEfficiency::CheckBinning",
"Histograms are not consistent: they have different number of bins");
1482 gROOT->Info(
"TEfficiency::CheckBinning",
"Histograms are not consistent: they have different bin edges");
1504 gROOT->Error(
"TEfficiency::CheckConsistency",
"passed TEfficiency objects have different dimensions");
1509 gROOT->Error(
"TEfficiency::CheckConsistency",
"passed TEfficiency objects have different binning");
1514 gROOT->Error(
"TEfficiency::CheckConsistency",
"passed TEfficiency objects do not have consistent bin contents");
1537 Int_t nbinsx, nbinsy, nbinsz, nbins;
1544 case 1: nbins = nbinsx + 2;
break;
1545 case 2: nbins = (nbinsx + 2) * (nbinsy + 2);
break;
1546 case 3: nbins = (nbinsx + 2) * (nbinsy + 2) * (nbinsz + 2);
break;
1550 for(
Int_t i = 0; i < nbins; ++i) {
1552 gROOT->Info(
"TEfficiency::CheckEntries",
"Histograms are not consistent: passed bin content > total bin content");
1574 double tolerance = (total.IsA() ==
TH1F::Class() ) ? 1.
E-5 : 1.
E-12;
1595 Error(
"CreatePaintingGraph",
"Call this function only for dimension == 1");
1618 Bool_t plot0Bins =
false;
1619 if (option.
Contains(
"e0") ) plot0Bins =
true;
1628 double * px = graph->
GetX();
1629 double * py = graph->
GetY();
1635 for (
Int_t i = 0; i < npoints; ++i) {
1644 if (j >= graph->
GetN() ) {
1664 if (oldTitle != newTitle ) {
1692 Error(
"CreatePaintingistogram",
"Call this function only for dimension == 2");
1741 for(
Int_t i = 0; i < nbinsx + 2; ++i) {
1742 for(
Int_t j = 0; j < nbinsy + 2; ++j) {
1803 Double_t alpha = (1.0 - level) / 2;
1914 for (
int i = 0; i <
n ; ++i) {
1915 if(pass[i] > total[i]) {
1916 ::Error(
"TEfficiency::Combine",
"total events = %i < passed events %i",total[i],pass[i]);
1917 ::Info(
"TEfficiency::Combine",
"stop combining");
1921 ntot += w[i] * total[i];
1922 ktot += w[i] * pass[i];
1927 double norm = sumw/sumw2;
1931 ::Error(
"TEfficiency::Combine",
"total = %f < passed %f",ntot,ktot);
1932 ::Info(
"TEfficiency::Combine",
"stop combining");
1936 double a = ktot + alpha;
1937 double b = ntot - ktot +
beta;
1939 double mean = a/(a+
b);
1945 if (shortestInterval)
1952 if (option.
Contains(
"mode"))
return mode;
2003 std::vector<TH1*> vTotal; vTotal.reserve(n);
2004 std::vector<TH1*> vPassed; vPassed.reserve(n);
2005 std::vector<Double_t> vWeights; vWeights.reserve(n);
2021 level = atof( opt(pos,opt.
Length() ).Data() );
2022 if((level <= 0) || (level >= 1))
2030 for(
Int_t k = 0; k <
n; ++k) {
2032 vWeights.push_back(w[k]);
2034 gROOT->Error(
"TEfficiency::Combine",
"invalid custom weight found w = %.2lf",w[k]);
2035 gROOT->Info(
"TEfficiency::Combine",
"stop combining");
2044 while((obj = next())) {
2070 vWeights.push_back(pEff->
fWeight);
2085 if(vTotal.empty()) {
2086 gROOT->Error(
"TEfficiency::Combine",
"no TEfficiency objects in given list");
2087 gROOT->Info(
"TEfficiency::Combine",
"stop combining");
2092 if(bWeights && (n != (
Int_t)vTotal.size())) {
2093 gROOT->Error(
"TEfficiency::Combine",
"number of weights n=%i differs from number of TEfficiency objects k=%i which should be combined",n,(
Int_t)vTotal.size());
2094 gROOT->Info(
"TEfficiency::Combine",
"stop combining");
2098 Int_t nbins_max = vTotal.at(0)->GetNbinsX();
2100 for(
UInt_t i=0; i<vTotal.size(); ++i) {
2102 gROOT->Warning(
"TEfficiency::Combine",
"histograms have not the same binning -> results may be useless");
2103 if(vTotal.at(i)->GetNbinsX() < nbins_max) nbins_max = vTotal.at(i)->GetNbinsX();
2108 gROOT->Info(
"TEfficiency::Combine",
"combining %i TEfficiency objects",(
Int_t)vTotal.size());
2110 gROOT->Info(
"TEfficiency::Combine",
"using custom weights");
2112 gROOT->Info(
"TEfficiency::Combine",
"using the following prior probability for the efficiency: P(e) ~ Beta(e,%.3lf,%.3lf)",alpha,beta);
2115 gROOT->Info(
"TEfficiency::Combine",
"using individual priors of each TEfficiency object");
2116 gROOT->Info(
"TEfficiency::Combine",
"confidence level = %.2lf",level);
2120 std::vector<Double_t>
x(nbins_max);
2121 std::vector<Double_t> xlow(nbins_max);
2122 std::vector<Double_t> xhigh(nbins_max);
2123 std::vector<Double_t> eff(nbins_max);
2124 std::vector<Double_t> efflow(nbins_max);
2125 std::vector<Double_t> effhigh(nbins_max);
2129 Int_t num = vTotal.size();
2130 std::vector<Int_t> pass(num);
2131 std::vector<Int_t>
total(num);
2136 for(
Int_t i=1; i <= nbins_max; ++i) {
2138 x[i-1] = vTotal.at(0)->GetBinCenter(i);
2139 xlow[i-1] = x[i-1] - vTotal.at(0)->GetBinLowEdge(i);
2140 xhigh[i-1] = vTotal.at(0)->GetBinWidth(i) - xlow[i-1];
2142 for(
Int_t j = 0; j < num; ++j) {
2143 pass[j] = (
Int_t)(vPassed.at(j)->GetBinContent(i) + 0.5);
2144 total[j] = (
Int_t)(vTotal.at(j)->GetBinContent(i) + 0.5);
2148 eff[i-1] =
Combine(up,low,num,&pass[0],&total[0],alpha,beta,level,&vWeights[0],opt.
Data());
2150 if(eff[i-1] == -1) {
2151 gROOT->Error(
"TEfficiency::Combine",
"error occurred during combining");
2152 gROOT->Info(
"TEfficiency::Combine",
"stop combining");
2155 efflow[i-1]= eff[i-1] - low;
2156 effhigh[i-1]= up - eff[i-1];
2199 if (option.
IsNull() ) option =
"ap";
2205 if (!option.
Contains(
"a") ) option +=
"a";
2209 if (!option.
Contains(
"p") ) option +=
"p";
2342 Bool_t bDeleteOld =
true;
2353 TF1* pFunc =
new TF1(*f1);
2358 while((obj = next())) {
2480 if (tw2 <= 0 )
return pw/tw;
2483 double norm = tw/tw2;
2484 aa = pw * norm + alpha;
2485 bb = (tw - pw) * norm + beta;
2489 aa = passed + alpha;
2490 bb = total - passed +
beta;
2533 if (tw2 <= 0)
return 0;
2538 Double_t bb = (tw - pw) * norm + beta;
2554 Warning(
"GetEfficiencyErrorLow",
"frequentist confidence intervals for weights are only supported by the normal approximation");
2555 Info(
"GetEfficiencyErrorLow",
"setting statistic option to kFNormal");
2559 Double_t variance = ( pw2 * (1. - 2 * eff) + tw2 * eff *eff ) / ( tw * tw) ;
2566 return (eff - delta < 0) ? eff : delta;
2613 if (tw2 <= 0)
return 0;
2618 Double_t bb = (tw - pw) * norm + beta;
2634 Warning(
"GetEfficiencyErrorUp",
"frequentist confidence intervals for weights are only supported by the normal approximation");
2635 Info(
"GetEfficiencyErrorUp",
"setting statistic option to kFNormal");
2639 Double_t variance = ( pw2 * (1. - 2 * eff) + tw2 * eff *eff ) / ( tw * tw) ;
2645 return (eff + delta > 1) ? 1.-eff : delta;
2702 while((obj = next())) {
2736 if (total == 0)
return (bUpper) ? 1 : 0;
2742 return ((average + delta) > 1) ? 1.0 : (average + delta);
2744 return ((average - delta) < 0) ? 0.0 : (average - delta);
2769 Fatal(
"operator+=",
"Adding to a non consistent TEfficiency object which has not a total or a passed histogram ");
2774 Warning(
"operator+=",
"no operation: adding an empty object");
2778 Fatal(
"operator+=",
"Adding a non consistent TEfficiency object which has not a total or a passed histogram ");
2834 rhs.TAttLine::Copy(*
this);
2835 rhs.TAttFill::Copy(*
this);
2836 rhs.TAttMarker::Copy(*
this);
2888 while((obj = next())) {
2891 ((
TF1*)obj)->Paint(
"sameC");
2911 Warning(
"Paint",
"Painting 3D efficiency is not implemented");
2925 static Int_t naxis = 0;
2926 TString sxaxis=
"xAxis",syaxis=
"yAxis",szaxis=
"zAxis";
2949 out << indent <<
"Double_t " << sxaxis <<
"[" 2952 if (i != 0) out <<
", ";
2955 out <<
"}; " << std::endl;
2958 out << indent <<
"Double_t " << syaxis <<
"[" 2961 if (i != 0) out <<
", ";
2964 out <<
"}; " << std::endl;
2968 out << indent <<
"Double_t " << szaxis <<
"[" 2971 if (i != 0) out <<
", ";
2974 out <<
"}; " << std::endl;
2979 static Int_t eff_count = 0;
2982 eff_name += eff_count;
2984 const char*
name = eff_name.
Data();
2987 const char quote =
'"';
2988 out << indent << std::endl;
2990 <<
"(" << quote <<
GetName() << quote <<
"," << quote
3018 out <<
");" << std::endl;
3019 out << indent << std::endl;
3022 out << indent << name <<
"->SetConfidenceLevel(" <<
fConfLevel <<
");" 3024 out << indent << name <<
"->SetBetaAlpha(" <<
fBeta_alpha <<
");" 3026 out << indent << name <<
"->SetBetaBeta(" <<
fBeta_beta <<
");" << std::endl;
3027 out << indent << name <<
"->SetWeight(" <<
fWeight <<
");" << std::endl;
3028 out << indent << name <<
"->SetStatisticOption(" <<
fStatisticOption <<
");" 3033 out << indent << name <<
"->SetUseWeightedEvents();" << std::endl;
3038 out << indent << name <<
"->SetBetaBinParameters(" << i <<
"," <<
fBeta_bin_params.at(i).first
3050 for(
Int_t i = 0; i < nbins; ++i) {
3051 out << indent << name <<
"->SetTotalEvents(" << i <<
"," <<
3053 out << indent << name <<
"->SetPassedEvents(" << i <<
"," <<
3060 while((obj = next())) {
3063 out << indent << name <<
"->GetListOfFunctions()->Add(" 3064 << obj->
GetName() <<
");" << std::endl;
3077 out<< indent << name<<
"->Draw(" << quote << opt << quote <<
");" 3096 Warning(
"SetBetaAlpha(Double_t)",
"invalid shape parameter %.2lf",alpha);
3114 Warning(
"SetBetaBeta(Double_t)",
"invalid shape parameter %.2lf",beta);
3154 Error(
"SetBins",
"Using wrong SetBins function for a %d-d histogram",
GetDimension());
3158 Warning(
"SetBins",
"Histogram entries will be lost after SetBins");
3174 Error(
"SetBins",
"Using wrong SetBins function for a %d-d histogram",
GetDimension());
3178 Warning(
"SetBins",
"Histogram entries will be lost after SetBins");
3194 Error(
"SetBins",
"Using wrong SetBins function for a %d-d histogram",
GetDimension());
3198 Warning(
"SetBins",
"Histogram entries will be lost after SetBins");
3214 Error(
"SetBins",
"Using wrong SetBins function for a %d-d histogram",
GetDimension());
3218 Warning(
"SetBins",
"Histogram entries will be lost after SetBins");
3235 Error(
"SetBins",
"Using wrong SetBins function for a %d-d histogram",
GetDimension());
3239 Warning(
"SetBins",
"Histogram entries will be lost after SetBins");
3256 Error(
"SetBins",
"Using wrong SetBins function for a %d-d histogram",
GetDimension());
3260 Warning(
"SetBins",
"Histogram entries will be lost after SetBins");
3275 if((level > 0) && (level < 1))
3278 Warning(
"SetConfidenceLevel(Double_t)",
"invalid confidence level %.2lf",level);
3315 TString name_passed = name + TString(
"_passed");
3330 if(events <= fTotalHistogram->GetBinContent(bin)) {
3335 Error(
"SetPassedEvents(Int_t,Int_t)",
"total number of events (%.1lf) in bin %i is less than given number of passed events %i",
fTotalHistogram->
GetBinContent(bin),bin,events);
3495 title_passed.
Insert(pos,
" (passed)");
3496 title_total.
Insert(pos,
" (total)");
3499 title_passed.
Append(
" (passed)");
3500 title_total.
Append(
" (total)");
3528 Error(
"SetTotalEvents(Int_t,Int_t)",
"passed number of events (%.1lf) in bin %i is bigger than given number of total events %i",
fPassedHistogram->
GetBinContent(bin),bin,events);
3588 gROOT->Info(
"TEfficiency::SetUseWeightedEvents",
"Handle weighted events for computing efficiency");
3608 Warning(
"SetWeight",
"invalid weight %.2lf",weight);
3636 if (total == 0)
return (bUpper) ? 1 : 0;
3640 Double_t mode = (passed + 0.5 * kappa * kappa) / (total + kappa * kappa);
3642 * (1 - average) + kappa * kappa / 4);
3644 return ((mode + delta) > 1) ? 1.0 : (mode + delta);
3646 return ((mode - delta) < 0) ? 0.0 : (mode - delta);
TGraphAsymmErrors * CreateGraph(Option_t *opt="") const
Create the graph used be painted (for dim=1 TEfficiency) The return object is managed by the caller...
Double_t * GetEXlow() const
double normal_quantile(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the lower tail of the normal (Gaussian) distri...
void SetConfidenceLevel(Double_t level)
Sets the confidence level (0 < level < 1) The default value is 1-sigma :~ 0.683.
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 ...
virtual const char * GetName() const
Returns name of object.
virtual Int_t GetNcells() const
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
virtual void Paint(Option_t *option="")
Control routine to paint any kind of histograms.
static Bool_t CheckEntries(const TH1 &pass, const TH1 &total, Option_t *opt="")
Checks whether bin contents are compatible with binomial statistics.
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
void Build(const char *name, const char *title)
Building standard data structure of a TEfficiency object.
Class to handle efficiency histograms.
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
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.
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
double normal_quantile_c(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the upper tail of the normal (Gaussian) distri...
Binomial fitter for the division of two histograms.
TString & ReplaceAll(const TString &s1, const TString &s2)
const Double_t * GetArray() const
virtual void SetBins(Int_t nx, Double_t xmin, Double_t xmax)
Redefine x axis parameters.
static Bool_t CheckBinning(const TH1 &pass, const TH1 &total)
Checks binning for each axis.
virtual void SetName(const char *name)
Set the name of the TNamed.
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
TEfficiency & operator+=(const TEfficiency &rhs)
Adds the histograms of another TEfficiency object to current histograms.
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
double beta_pdf(double x, double a, double b)
Probability density function of the beta distribution.
TAxis * GetYaxis() const
Get y axis of the graph.
Double_t GetBetaBeta(Int_t bin=-1) const
virtual Int_t GetNbinsZ() const
const TEfficiency operator+(const TEfficiency &lhs, const TEfficiency &rhs)
Addition operator.
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
const Double_t kDefBetaBeta
static Bool_t AddDirectoryStatus()
Static function: cannot be inlined on Windows/NT.
void FillGraph(TGraphAsymmErrors *graph, Option_t *opt) const
Fill the graph to be painted with information from TEfficiency Internal method called by TEfficiency:...
void ToLower()
Change string to lower-case.
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
virtual void SetTitle(const char *title="")
Set graph title.
void SetBetaBinParameters(Int_t bin, Double_t alpha, Double_t beta)
Sets different shape parameter α and β for the prior distribution for each bin.
TFitResultPtr Fit(TF1 *f1, Option_t *opt="")
Fits the efficiency using the TBinomialEfficiencyFitter class.
~TEfficiency()
default destructor
virtual Double_t GetBinLowEdge(Int_t bin) const
Return bin lower edge for 1D histogram.
EStatOption fStatisticOption
Bool_t SetPassedEvents(Int_t bin, Int_t events)
Sets the number of passed events in the given global bin.
TString & Insert(Ssiz_t pos, const char *s)
Double_t(* fBoundary)(Double_t, Double_t, Double_t, Bool_t)
TRObject operator()(const T1 &t1) const
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
virtual void GetStats(Double_t *stats) const
fill the array stats from the contents of this histogram The array stats must be correctly dimensione...
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
double beta(double x, double y)
Calculates the beta function.
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Double_t GetWeight() const
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Template class to wrap any C++ callable object which takes one argument i.e.
TH1 * GetCopyPassedHisto() const
Returns a cloned version of fPassedHistogram.
virtual Int_t GetDimension() const
static struct mg_connection * fc(struct mg_context *ctx)
TGraph with asymmetric error bars.
EStatOption GetStatisticOption() const
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Double_t GetEfficiency(Int_t bin) const
Returns the efficiency in the given global bin.
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="")
Calculates the combined efficiency and its uncertainties.
Fill Area Attributes class.
virtual void Paint(Option_t *chopt="")
Draw this graph with its current attributes.
virtual void SetName(const char *name="")
Set graph name.
TGraphAsymmErrors * fPaintGraph
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.
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
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.
The TNamed class is the base class for all named ROOT classes.
static Bool_t CheckWeights(const TH1 &pass, const TH1 &total)
Check if both histogram are weighted.
TH1 * fPassedHistogram
temporary histogram for painting
TH2 * fPaintHisto
temporary graph for painting
TList * GetListOfFunctions()
TString & Append(const char *cs)
TList * fFunctions
pointer to directory holding this TEfficiency object
Ssiz_t First(char c) const
Find first occurrence of a character c.
Bool_t UsesBayesianStat() const
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.
virtual TArrayD * GetSumw2()
Bool_t SetPassedHistogram(const TH1 &rPassed, Option_t *opt)
Sets the histogram containing the passed events.
void SetTitle(const char *title)
Sets the title.
User class for performing function minimization.
static Double_t BetaMean(Double_t alpha, Double_t beta)
Compute the mean (average) of the beta distribution.
double beta_cdf(double x, double a, double b)
Cumulative distribution function of the beta distribution Upper tail of the integral of the beta_pdf...
The 3-D histogram classes derived from the 1-D histogram classes.
object has not been deleted
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...
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...
Double_t fConfLevel
pointer to a method calculating the boundaries of confidence intervals
void SetStatisticOption(EStatOption option)
Sets the statistic option which affects the calculation of the confidence interval.
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
TH1 * GetCopyTotalHisto() const
Returns a cloned version of fTotalHistogram.
Bool_t SetTotalEvents(Int_t bin, Int_t events)
Sets the number of total events in the given global bin.
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a line.
void SetWeight(Double_t weight)
Sets the global weight for this TEfficiency object.
THist< 3, double, THistStatContent, THistStatUncertainty > TH3D
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
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.
const char * GetTitle() const
Returns title of object.
Service class for 2-Dim histogram classes.
Class to manage histogram axis.
void SetBetaAlpha(Double_t alpha)
Sets the shape parameter α.
if object ctor succeeded but object should not be used
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
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 distr...
void SetFunction(const ROOT::Math::IGenFunction &f, double xlow, double xup)
Sets function to be minimized.
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Collection abstract base class.
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.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
virtual void Append(TObject *obj, Bool_t replace=kFALSE)
Append object to this directory.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a graph.
const Double_t kDefWeight
void FillHistogram(TH2 *h2) const
Fill the 2d histogram to be painted with information from TEfficiency 2D Internal method called by TE...
static double p1(double t, double a, double b)
void SavePrimitive(std::ostream &out, Option_t *opt="")
Have histograms fixed bins along each axis?
const Double_t kDefBetaAlpha
virtual Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Int_t FindFixBin(Double_t x, Double_t y=0, Double_t z=0) const
Returns the global bin number containing the given values.
constexpr Double_t E()
Base of natural log: .
TAxis * GetXaxis() const
Get x axis of the graph.
double beta_quantile_c(double x, double a, double b)
Inverse ( ) of the cumulative distribution function of the lower tail of the beta distribution (beta_...
virtual void SetName(const char *name)
Change the name of this histogram.
Double_t At(Int_t i) const
static unsigned int total
virtual Int_t GetSumw2N() const
Bin contents are average (used by Add)
virtual Bool_t IsEmpty() const
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width for 1D histogram.
virtual TObject * Remove(TObject *)
Remove an object from the in-memory list.
const TEfficiency::EStatOption kDefStatOpt
void SetDirectory(TDirectory *dir)
Sets the directory holding this TEfficiency object.
TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases...
void SetUseWeightedEvents(Bool_t on=kTRUE)
Bool_t SetTotalHistogram(const TH1 &rTotal, Option_t *opt)
Sets the histogram containing all events.
virtual void PaintStats(TF1 *fit)
Draw the stats.
Describe directory structure in memory.
Double_t GetEfficiencyErrorUp(Int_t bin) const
Returns the upper error on the efficiency in the given global bin.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a graph.
Double_t GetConfidenceLevel() const
Double_t GetBetaAlpha(Int_t bin=-1) const
virtual double FValMinimum() const
Return function value at current estimate of the minimum.
TEfficiency & operator=(const TEfficiency &rhs)
Assignment operator.
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
void SetBetaBeta(Double_t beta)
Sets the shape parameter β.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
virtual Double_t GetEntries() const
Return the current number of entries.
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Long64_t Merge(TCollection *list)
Merges the TEfficiency objects in the given list to the given TEfficiency object using the operator+=...
Bool_t IsVariableBinSize() const
virtual Bool_t Add(TF1 *h1, Double_t c1=1, Option_t *option="")
Performs the operation: this = this + c1*f1 if errors are defined (see TH1::Sumw2), errors are also recalculated.
TH2 * CreateHistogram(Option_t *opt="") const
Create the histogram used to be painted (for dim=2 TEfficiency) The return object is managed by the c...
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.
double beta_cdf_c(double x, double a, double b)
Complement of the cumulative distribution function of the beta distribution.
virtual double XMinimum() const
Return current estimate of the position of the minimum.
Mother of all ROOT objects.
virtual Int_t FindFixBin(Double_t x) const
Find bin number corresponding to abscissa x.
you should not use this method at all Int_t Int_t z
Double_t GetEfficiencyErrorLow(Int_t bin) const
Returns the lower error on the efficiency in the given global bin.
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: ...
void SetName(const char *name)
Sets the name.
TFitResultPtr Fit(TF1 *f1, Option_t *option="")
Carry out the fit of the given function to the given histograms.
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
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.
void Draw(Option_t *opt="")
Draws the current TEfficiency object.
void Paint(Option_t *opt)
Paints this TEfficiency object.
virtual void Add(TObject *obj)
Double_t * GetEYlow() const
Double_t * GetEYhigh() const
void SetNpx(int npx)
Set the number of point used to bracket root using a grid.
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.
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.
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
const Double_t kDefConfLevel
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
std::vector< std::pair< Double_t, Double_t > > fBeta_bin_params
Int_t GetDimension() const
returns the dimension of the current TEfficiency object
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
virtual bool Minimize(int maxIter, double absTol=1.E-8, double relTol=1.E-10)
Find minimum position iterating until convergence specified by the absolute and relative tolerance or...
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
virtual Int_t GetNbinsX() const
Double_t * GetEXhigh() const
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
virtual const char * GetName() const
Returns name of object.
virtual void Set(Int_t n)
Set number of points in the graph Existing coordinates are preserved New coordinates above fNpoints a...
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
THist< 2, float, THistStatContent, THistStatUncertainty > TH2F
virtual void SetNormFactor(Double_t factor=1)
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
TEfficiency()
default constructor
const TArrayD * GetXbins() const
static Double_t BetaMode(Double_t alpha, Double_t beta)
Compute the mode of the beta distribution.
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.
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
static Double_t Wilson(Double_t total, Double_t passed, Double_t level, Bool_t bUpper)
Calculates the boundaries for the frequentist Wilson interval.
static Bool_t CheckConsistency(const TH1 &pass, const TH1 &total, Option_t *opt="")
Checks the consistence of the given histograms.
virtual Int_t GetNbinsY() const
virtual const char * GetTitle() const
Returns title of object.
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
double beta_quantile(double x, double a, double b)
Inverse ( ) of the cumulative distribution function of the upper tail of the beta distribution (beta_...
const char * Data() const
void Calculate(const double X, const double n)
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
virtual void SetPointError(Double_t exl, Double_t exh, Double_t eyl, Double_t eyh)
Set ex and ey values for point pointed by the mouse.