00001 #ifndef ROOT_TEfficiency_cxx 00002 #define ROOT_TEfficiency_cxx 00003 00004 //standard header 00005 #include <vector> 00006 #include <string> 00007 #include <cmath> 00008 #include <stdlib.h> 00009 00010 //ROOT headers 00011 #include "Math/ProbFunc.h" 00012 #include "Math/QuantFunc.h" 00013 #include "TBinomialEfficiencyFitter.h" 00014 #include "TDirectory.h" 00015 #include "TF1.h" 00016 #include "TGraphAsymmErrors.h" 00017 #include "TH1.h" 00018 #include "TH2.h" 00019 #include "TH3.h" 00020 #include "TList.h" 00021 #include "TMath.h" 00022 #include "TROOT.h" 00023 #include "TStyle.h" 00024 #include "TVirtualPad.h" 00025 #include "TError.h" 00026 #include "Math/BrentMinimizer1D.h" 00027 #include "Math/WrappedFunction.h" 00028 00029 //custom headers 00030 #include "TEfficiency.h" 00031 00032 // file with extra class for FC method 00033 #include "TEfficiencyHelper.h" 00034 00035 //default values 00036 const Double_t kDefBetaAlpha = 1; 00037 const Double_t kDefBetaBeta = 1; 00038 const Double_t kDefConfLevel = 0.682689492137; // 1 sigma 00039 const TEfficiency::EStatOption kDefStatOpt = TEfficiency::kFCP; 00040 const Double_t kDefWeight = 1; 00041 00042 ClassImp(TEfficiency) 00043 00044 //______________________________________________________________________________ 00045 // Begin_Html <center><h2>TEfficiency - a class to handle efficiency 00046 // histograms</h2></center> 00047 // <ol style="list-style-type: upper-roman;"> 00048 // <li><a href="#over">Overview</a></li> 00049 // <li><a href="#create">Creating a TEfficiency object</a></li> 00050 // <li><a href="#fill">Fill in events</a></li> 00051 // <li><a href="#stat">Statistic options</a></li> 00052 // <ol><li><a href="#compare">Coverage probabilities for different methods</a></li></ol> 00053 // <li><a href="#cm">Merging and combining TEfficiency objects</a></li> 00054 // <ol> 00055 // <li><a href="#merge">When should I use merging?</a></li> 00056 // <li><a href="#comb">When should I use combining?</a></li> 00057 // </ol> 00058 // <li><a href="#other">Further operations</a> 00059 // <ol> 00060 // <li><a href="#histo">Information about the internal histograms</a></li> 00061 // <li><a href="#fit">Fitting</a></li> 00062 // <li><a href="#draw">Draw a TEfficiency object</a></li> 00063 // </ol> 00064 // <li><a href="#class">TEfficiency class</a></li> 00065 // </ol> 00066 // 00067 // <h3><a name="over">I. Overview</a></h3> 00068 // This class handles the calculation of efficiencies and their uncertainties. It 00069 // provides several statistical methods for calculating frequentist and bayesian 00070 // confidence intervals as well as a function for combining several efficiencies. 00071 // <br /> 00072 // Efficiencies have a lot of applications and meanings but in principle they can 00073 // be described by the fraction of good/passed events k out of sample containing 00074 // N events. One is usually interested in the dependency of the efficiency on other 00075 // (binned) variables. The number of passed and total events is therefore stored 00076 // internally in two histograms (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:fTotalHistogram">fTotalHistogram</a> and 00077 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:fPassedHistogram">fPassedHistogram</a>). Then the 00078 // efficiency as well as its upper and lower error an be calculated for each bin 00079 // individually.<br /> 00080 // As the efficiency can be regarded as a parameter of a binomial distribution, the 00081 // number of pass ed and total events must always be integer numbers. Therefore a 00082 // filling with weights is not possible however you can assign a global weight to each 00083 // TEfficiency object (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetWeight">SetWeight</a>). It is necessary to create one TEfficiency object 00084 // for each weight if you investigate a process involving different weights. This 00085 // procedure needs more effort but enables you to re-use the filled object in cases 00086 // where you want to change one or more weights. This would not be possible if all 00087 // events with different weights were filled in the same histogram. 00088 // 00089 // <h3><a name="create">II. Creating a TEfficiency object</a></h3> 00090 // If you start a new analysis, it is highly recommended to use the TEfficiency class 00091 // from the beginning. You can then use one of the constructors for fixed or 00092 // variable bin size and your desired dimension. These constructors append the 00093 // created TEfficiency object to the current directory. So it will be written 00094 // automatically to a file during the next <a href="http://root.cern.ch/root/html/TFile.html#TFile:Write">TFile::Write</a> command. 00095 // <div class="code"><pre> 00096 // <b>Example:</b> 00097 //create a twodimensional TEfficiency object with 00098 //- name = "eff" 00099 //- title = "my efficiency" 00100 //- axistitles: x, y and LaTeX formated epsilon as label for Z axis 00101 //- 10 bins with constant bin width (= 1) along X axis starting at 0 (lower edge 00102 // from first bin) upto 10 (upper edge of last bin) 00103 //- 20 bins with constant bin width (= 0.5) along Y axis starting at -5 (lower 00104 // edge from first bin) upto 5 (upper edge of last bin) 00105 // TEfficiency* pEff = new TEfficiency("eff","my efficiency;x;y;#epsilon",10,0,10,20,-5,5); 00106 // </pre></div><div class="clear" /> 00107 // If you already have two histograms filled with the number of passed and total 00108 // events, you will use the constructor <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:TEfficiency%1">TEfficiency(const TH1& passed,const TH1& total)</a> 00109 // to construct the TEfficiency object. The histograms "passed" and "total" have 00110 // to fullfill the conditions mentioned in <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:CheckConsistency">CheckConsistency</a>, otherwise the construction will fail. 00111 // As the histograms already exist, the new TEfficiency is by default <b>not</b> attached 00112 // to the current directory to avoid duplication of data. If you want to store the 00113 // new object anyway, you can either write it directly by calling <a href="http://root.cern.ch/root/html/TObject.html#TObject:Write">Write</a> or attach it to a directory using <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetDirectory">SetDirectory</a>. 00114 // This also applies for TEfficiency objects created by the copy constructor <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:TEfficiency%8">TEfficiency(const TEfficiency& rEff)</a>. 00115 // <div class="code"> 00116 // <pre> 00117 // <b>Example 1:</b> 00118 // TEfficiency* pEff = 0; 00119 // TFile* pFile = new TFile("myfile.root","recreate"); 00120 // 00121 // //h_pass and h_total are valid and consistent histograms 00122 // if(TEfficiency::CheckConsistency(h_pass,h_total)) 00123 // { 00124 // pEff = new TEfficiency(h_pass,h_total); 00125 // // this will write the TEfficiency object to "myfile.root" 00126 // // AND pEff will be attached to the current directory 00127 // pEff->Write(); 00128 // } 00129 // 00130 // <b>Example 2:</b> 00131 // TEfficiency* pEff = 0; 00132 // TFile* pFile = new TFile("myfile.root","recreate"); 00133 // 00134 // //h_pass and h_total are valid and consistent histograms 00135 // if(TEfficiency::CheckConsistency(h_pass,h_total)) 00136 // { 00137 // pEff = new TEfficiency(h_pass,h_total); 00138 // //this will attach the TEfficiency object to the current directory 00139 // pEff->SetDirectory(gDirectory); 00140 // //now all objects in gDirectory will be written to "myfile.root" 00141 // pFile->Write(); 00142 // } 00143 // </pre> 00144 // </div><div class="clear" /> 00145 // In the case that you already have two filled histograms and you only want to 00146 // plot them as a graph, you should rather use <a href="http://root.cern.ch/root/html/TGraphAsymmErrors.html#TGraphAsymmErrors:TGraphAsymmErrors%8">TGraphAsymmErrors::TGraphAsymmErrors(const TH1* pass,const TH1* total,Option_t* opt)</a> 00147 // to create a graph object. 00148 // 00149 // <h3><a name="fill">III. Filling with events</a></h3> 00150 // You can fill the TEfficiency object by calling the <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Fill">Fill(Bool_t bPassed,Double_t x,Double_t y,Double_t z)</a> method. 00151 // The boolean flag "bPassed" indicates whether the current event is a good 00152 // (both histograms are filled) or not (only <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:fTotalHistogram">fTotalHistogram</a> is filled). 00153 // The variables x,y and z determine the bin which is filled. For lower 00154 // dimensions the z- or even the y-value may be omitted. 00155 // End_Html 00156 // Begin_Macro(source) 00157 // { 00158 // //canvas only needed for this documentation 00159 // TCanvas* c1 = new TCanvas("example","",600,400); 00160 // c1->SetFillStyle(1001); 00161 // c1->SetFillColor(kWhite); 00162 // 00163 // //create one-dimensional TEfficiency object with fixed bin size 00164 // TEfficiency* pEff = new TEfficiency("eff","my efficiency;x;#epsilon",20,0,10); 00165 // TRandom3 rand; 00166 // 00167 // bool bPassed; 00168 // double x; 00169 // for(int i=0; i<10000; ++i) 00170 // { 00171 // //simulate events with variable under investigation 00172 // x = rand.Uniform(10); 00173 // //check selection: bPassed = DoesEventPassSelection(x) 00174 // bPassed = rand.Rndm() < TMath::Gaus(x,5,4); 00175 // pEff->Fill(bPassed,x); 00176 // } 00177 // 00178 // pEff->Draw("AP"); 00179 // 00180 // //only for this documentation 00181 // return c1; 00182 // } 00183 // End_Macro 00184 // Begin_Html 00185 // You can also set the number of passed or total events for a bin directly by 00186 // using the <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetPassedEvents">SetPassedEvents</a> or <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetTotalEvents">SetTotalEvents</a> method. 00187 // 00188 // <h3><a name="stat">IV. Statistic options</a></h3> 00189 // The calculation of the estimated efficiency depends on the chosen statistic 00190 // option. Let k denotes the number of passed events and N the number of total 00191 // events.<br /> 00192 // <b>Frequentist methods</b><br /> 00193 // The expectation value of the number of passed events is given by the true 00194 // efficiency times the total number of events. One can estimate the efficiency 00195 // by replacing the expected number of passed events by the observed number of 00196 // passed events. End_Html 00197 // Begin_Latex 00198 // #LT k #GT = #epsilon #times N #Rightarrow #hat{#varepsilon} = #frac{k}{N} 00199 // End_Latex 00200 // Begin_Html 00201 // <b>Bayesian methods</b><br /> 00202 // In bayesian statistics a likelihood-function (how probable is it to get the 00203 // observed data assuming a true efficiency) and a prior probability (what is the 00204 // probability that a certain true efficiency is actually realised) are used to 00205 // determine a posterior probability by using Bayes theorem. At the moment, only 00206 // beta distributions (have 2 free parameters) are supported as prior 00207 // probabilities. 00208 // End_Html 00209 // Begin_Latex(separator='=',align='rl') 00210 // P(#epsilon | k ; N) = #frac{1}{norm} #times P(k | #epsilon ; N) #times Prior(#epsilon) 00211 // P(k | #epsilon ; N) = Binomial(N,k) #times #epsilon^{k} #times (1 - #epsilon)^{N - k} ... binomial distribution 00212 // Prior(#epsilon) = #frac{1}{B(#alpha,#beta)} #times #epsilon ^{#alpha - 1} #times (1 - #epsilon)^{#beta - 1} #equiv Beta(#epsilon; #alpha,#beta) 00213 // #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) 00214 // End_Latex 00215 // Begin_Html 00216 // By default the expectation value of this posterior distribution is used as estimator for the efficiency: 00217 // End_Html 00218 // Begin_Latex 00219 // #hat{#varepsilon} = #frac{k + #alpha}{N + #alpha + #beta} 00220 // End_Latex 00221 // Begin_Html 00222 // Optionally the mode can also be used as value for the estimated efficiency. This can be done by calling SetBit(kPosteriorMode) or 00223 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetPosteriorMode">SetPosteriorMode</a>. In this case the estimated efficiency is: 00224 // End_Html 00225 // Begin_Latex 00226 // #hat{#varepsilon} = #frac{k + #alpha -1}{N + #alpha + #beta - 2} 00227 // End_Latex 00228 // Begin_Html 00229 // 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). 00230 // 00231 // The statistic options also specifiy which confidence interval is used for calculating 00232 // the uncertainties of the efficiency. The following properties define the error 00233 // calculation: 00234 // <ul> 00235 // <li><b>fConfLevel:</b> desired confidence level: 0 < fConfLevel < 1 (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetConfidenceLevel">GetConfidenceLevel</a> / <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetConfidenceLevel">SetConfidenceLevel</a>)</li> 00236 // <li><b>fStatisticOption:</b> defines which method is used to calculate the boundaries of the confidence interval (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetStatisticOption">SetStatisticOption</a>)</li> 00237 // <li><b>fBeta_alpha, fBeta_beta:</b> parameters for the prior distribution which is only used in the bayesian case (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetBetaAlpha">GetBetaAlpha</a> / <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetBetaBeta">GetBetaBeta</a> / <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetBetaAlpha">SetBetaAlpha</a> / <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetBetaBeta">SetBetaBeta</a>)</li> 00238 // <li><b>kIsBayesian:</b> flag whether bayesian statistics are used or not (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:UsesBayesianStat">UsesBayesianStat</a>)</li> 00239 // <li><b>kShortestInterval:</b> flag whether shortest interval (instead of central one) are used in case of Bayesian statistics (<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:UsesShortestInterval">UsesShortestInterval</a>). Normally shortest interval should be used in combination with the mode (see <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:UsesPosteriorMode">UsesPosteriorMode</a>)</li> 00240 // <li><b>fWeight:</b> global weight for this TEfficiency object which is used during combining or merging with other TEfficiency objects(<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetWeight">GetWeight</a> / <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetWeight">SetWeight</a>)</li> 00241 // </ul> 00242 // In the following table the implemented confidence intervals are listed 00243 // with their corresponding statistic option. For more details on the calculation, 00244 // please have a look at the the mentioned functions.<br /><br /> 00245 // <table align="center" border="1" cellpadding="5" rules="rows" vspace="10"> 00246 // <caption align="bottom">implemented confidence intervals and their options</caption> 00247 // <tr> 00248 // <th>name</th><th>statistic option</th><th>function</th><th>kIsBayesian</th><th>parameters</th> 00249 // </tr> 00250 // <tr> 00251 // <td>Clopper-Pearson</td><td>kFCP</td> 00252 // <td> 00253 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:ClopperPearson">ClopperPearson</a> 00254 // </td> 00255 // <td>false</td> 00256 // <td> 00257 // <ul><li>total events</li><li>passed events</li><li>confidence level</li></ul> 00258 // </td> 00259 // </tr> 00260 // <tr> 00261 // <td>normal approximation</td><td>kFNormal</td> 00262 // <td> 00263 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Normal">Normal</a> 00264 // </td> 00265 // <td>false</td> 00266 // <td> 00267 // <ul><li>total events</li><li>passed events</li><li>confidence level</li></ul> 00268 // </td> 00269 // </tr> 00270 // <tr> 00271 // <td>Wilson</td><td>kFWilson</td> 00272 // <td> 00273 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Wilson">Wilson</a> 00274 // </td> 00275 // <td>false</td> 00276 // <td> 00277 // <ul><li>total events</li><li>passed events</li><li>confidence level</li></ul> 00278 // </td> 00279 // </tr> 00280 // <tr> 00281 // <td>Agresti-Coull</td><td>kFAC</td> 00282 // <td> 00283 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:AgrestiCoull">AgrestiCoull</a> 00284 // </td> 00285 // <td>false</td> 00286 // <td> 00287 // <ul><li>total events</li><li>passed events</li><li>confidence level</li></ul> 00288 // </td> 00289 // </tr> 00290 // <tr> 00291 // <td>Feldman-Cousins</td><td>kFAC</td> 00292 // <td> 00293 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:FeldmanCousins">FeldmanCousins</a> 00294 // </td> 00295 // <td>false</td> 00296 // <td> 00297 // <ul><li>total events</li><li>passed events</li><li>confidence level</li></ul> 00298 // </td> 00299 // </tr> 00300 // <tr> 00301 // <td>Jeffrey</td><td>kBJeffrey</td> 00302 // <td> 00303 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Bayesian">Bayesian</a> 00304 // </td> 00305 // <td>true</td> 00306 // <td> 00307 // <ul><li>total events</li><li>passed events</li><li>confidence level</li><li>fBeta_alpha = 0.5</li><li>fBeta_beta = 0.5</li></ul> 00308 // </td> 00309 // </tr> 00310 // <td>Uniform prior</td><td>kBUniform</td> 00311 // <td> 00312 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Bayesian">Bayesian</a> 00313 // </td> 00314 // <td>true</td> 00315 // <td> 00316 // <ul><li>total events</li><li>passed events</li><li>confidence level</li><li>fBeta_alpha = 1</li><li>fBeta_beta = 1</li></ul> 00317 // </td> 00318 // </tr> 00319 // <td>custom prior</td><td>kBBayesian</td> 00320 // <td> 00321 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Bayesian">Bayesian</a> 00322 // </td> 00323 // <td>true</td> 00324 // <td> 00325 // <ul><li>total events</li><li>passed events</li><li>confidence level</li><li>fBeta_alpha</li><li>fBeta_beta</li></ul> 00326 // </td> 00327 // </tr> 00328 // </table> 00329 // <br /> 00330 // The following example demonstrates the effect of different statistic options and 00331 // confidence levels. 00332 // End_Html 00333 // Begin_Macro(source) 00334 // { 00335 // //canvas only needed for the documentation 00336 // TCanvas* c1 = new TCanvas("c1","",600,400); 00337 // c1->Divide(2); 00338 // c1->SetFillStyle(1001); 00339 // c1->SetFillColor(kWhite); 00340 // 00341 // //create one-dimensional TEfficiency object with fixed bin size 00342 // TEfficiency* pEff = new TEfficiency("eff","different confidence levels;x;#epsilon",20,0,10); 00343 // TRandom3 rand; 00344 // 00345 // bool bPassed; 00346 // double x; 00347 // for(int i=0; i<1000; ++i) 00348 // { 00349 // //simulate events with variable under investigation 00350 // x = rand.Uniform(10); 00351 // //check selection: bPassed = DoesEventPassSelection(x) 00352 // bPassed = rand.Rndm() < TMath::Gaus(x,5,4); 00353 // pEff->Fill(bPassed,x); 00354 // } 00355 // 00356 // //set style attributes 00357 // pEff->SetFillStyle(3004); 00358 // pEff->SetFillColor(kRed); 00359 // 00360 // //copy current TEfficiency object and set new confidence level 00361 // TEfficiency* pCopy = new TEfficiency(*pEff); 00362 // pCopy->SetConfidenceLevel(0.90); 00363 // 00364 // //set style attributes 00365 // pCopy->SetFillStyle(3005); 00366 // pCopy->SetFillColor(kBlue); 00367 // 00368 // c1->cd(1); 00369 // 00370 // //add legend 00371 // TLegend* leg1 = new TLegend(0.3,0.1,0.7,0.5); 00372 // leg1->AddEntry(pEff,"95%","F"); 00373 // leg1->AddEntry(pCopy,"68.3%","F"); 00374 // 00375 // pEff->Draw("A4"); 00376 // pCopy->Draw("same4"); 00377 // leg1->Draw("same"); 00378 // 00379 // //use same confidence level but different statistic methods 00380 // TEfficiency* pEff2 = new TEfficiency(*pEff); 00381 // TEfficiency* pCopy2 = new TEfficiency(*pEff); 00382 // 00383 // pEff2->SetStatisticOption(TEfficiency::kFNormal); 00384 // pCopy2->SetStatisticOption(TEfficiency::kFAC); 00385 // 00386 // pEff2->SetTitle("different statistic options;x;#epsilon"); 00387 // 00388 // //set style attributes 00389 // pCopy2->SetFillStyle(3005); 00390 // pCopy2->SetFillColor(kBlue); 00391 // 00392 // c1->cd(2); 00393 // 00394 // //add legend 00395 // TLegend* leg2 = new TLegend(0.3,0.1,0.7,0.5); 00396 // leg2->AddEntry(pEff2,"kFNormal","F"); 00397 // leg2->AddEntry(pCopy2,"kFAC","F"); 00398 // 00399 // pEff2->Draw("a4"); 00400 // pCopy2->Draw("same4"); 00401 // leg2->Draw("same"); 00402 // 00403 // //only for this documentation 00404 // c1->cd(0); 00405 // return c1; 00406 // } 00407 // End_Macro 00408 // Begin_Html 00409 // The prior probability of the efficiency in bayesian statistics can be given 00410 // in terms of a beta distribution. The beta distribution has to positive shape 00411 // parameters. The resulting priors for different combinations of these shape 00412 // parameters are shown in the plot below. 00413 // End_Html 00414 // Begin_Macro(source) 00415 // { 00416 // //canvas only needed for the documentation 00417 // TCanvas* c1 = new TCanvas("c1","",600,400); 00418 // c1->SetFillStyle(1001); 00419 // c1->SetFillColor(kWhite); 00420 // 00421 // //create different beta distributions 00422 // TF1* f1 = new TF1("f1","TMath::BetaDist(x,1,1)",0,1); 00423 // f1->SetLineColor(kBlue); 00424 // TF1* f2 = new TF1("f2","TMath::BetaDist(x,0.5,0.5)",0,1); 00425 // f2->SetLineColor(kRed); 00426 // TF1* f3 = new TF1("f3","TMath::BetaDist(x,1,5)",0,1); 00427 // f3->SetLineColor(kGreen+3); 00428 // f3->SetTitle("Beta distributions as priors;#epsilon;P(#epsilon)"); 00429 // TF1* f4 = new TF1("f4","TMath::BetaDist(x,4,3)",0,1); 00430 // f4->SetLineColor(kViolet); 00431 // 00432 // //add legend 00433 // TLegend* leg = new TLegend(0.25,0.5,0.85,0.89); 00434 // leg->SetFillColor(kWhite); 00435 // leg->SetFillStyle(1001); 00436 // leg->AddEntry(f1,"a=1, b=1","L"); 00437 // leg->AddEntry(f2,"a=0.5, b=0.5","L"); 00438 // leg->AddEntry(f3,"a=1, b=5","L"); 00439 // leg->AddEntry(f4,"a=4, b=3","L"); 00440 // 00441 // f3->Draw(); 00442 // f1->Draw("same"); 00443 // f2->Draw("Same"); 00444 // f4->Draw("same"); 00445 // leg->Draw("same"); 00446 // 00447 // //only for this documentation 00448 // return c1; 00449 // } 00450 // End_Macro 00451 // Begin_Html 00452 // 00453 // <h4><a name="compare">IV.1 Coverage probabilities for different methods</a></h4> 00454 // The following pictures illustrate the actual coverage probability for the 00455 // different values of the true efficiency and the total number of events when a 00456 // confidence level of 95% is desired. 00457 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/normal95.gif" alt="normal approximation" width="600" height="400" /></p> 00458 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/wilson95.gif" alt="wilson" width="600" height="400" /></p> 00459 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/ac95.gif" alt="agresti coull" width="600" height="400" /></p> 00460 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/cp95.gif" alt="clopper pearson" width="600" height="400" /></p> 00461 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/uni95.gif" alt="uniform prior" width="600" height="400" /></p> 00462 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/jeffrey95.gif" alt="jeffrey prior" width="600" height="400" /></p> 00463 // 00464 // The average (over all possible true efficiencies) coverage probability for 00465 // different number of total events is shown in the next picture. 00466 // <p><img src="http://root.cern.ch/drupal/sites/default/files/images/av_cov.png" alt="average coverage" width="600" height="400" /></p> 00467 // <h3><a name="cm">V. Merging and combining TEfficiency objects</a></h3> 00468 // In many applications the efficiency should be calculated for an inhomogenous 00469 // sample in the sense that it contains events with different weights. In order 00470 // to be able to determine the correct overall efficiency, it is necessary to 00471 // use for each subsample (= all events with the same weight) a different 00472 // TEfficiency object. After finsihing your analysis you can then construct the 00473 // overall efficiency with its uncertainty.<br /> 00474 // This procedure has the advantage that you can change the weight of one 00475 // subsample easily without rerunning the whole analysis. On the other hand more 00476 // efford is needed to handle several TEfficiency objects instead of one 00477 // histogram. In the case of many different or even continuously distributed 00478 // weights this approach becomes cumbersome. One possibility to overcome this 00479 // problem is the usage of binned weights.<br /><br /> 00480 // <b>Example</b><br /> 00481 // In high particle physics weights arises from the fact that you want to 00482 // normalise your results to a certain reference value. A very common formula for 00483 // calculating weights is End_Html 00484 // Begin_Latex(separator='-') 00485 // w = #frac{#sigma L}{N_{gen} #epsilon_{trig}} - #sigma ... cross section 00486 // - L ... luminosity 00487 // - N_{gen} ... number of generated events 00488 // - #epsilon_{trig} ... (known) trigger efficiency 00489 // End_Latex 00490 // Begin_Html 00491 // The reason for different weights can therefore be:<ul> 00492 // <li>different processes</li> 00493 // <li>other integrated luminosity</li> 00494 // <li>varying trigger efficiency</li> 00495 // <li>different sample sizes</li> 00496 // <li>...</li> 00497 // <li>or even combination of them</li> 00498 // </ul> 00499 // Depending on the actual meaning of different weights in your case, you 00500 // should either merge or combine them to get the overall efficiency. 00501 // 00502 // <h4><a name="merge">V.1 When should I use merging?</a></h4> 00503 // If the weights are artificial and do not represent real alternative hypotheses, 00504 // you should merge the different TEfficiency objects. That means especially for 00505 // the bayesian case that the prior probability should be the same for all merged 00506 // TEfficiency objects. The merging can be done by invoking one of the following 00507 // operations: 00508 // <ul> 00509 // <li> <b>eff1</b>.<a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Add">Add</a>(eff2)</li> 00510 // <li> <b>eff1</b> <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:operator+=">+=</a> eff2</li> 00511 // <li> <b>eff</b> = eff1 + eff2</li> 00512 // </ul> 00513 // The result of the merging is stored in the TEfficiency object which is marked 00514 // bold above. The contents of the internal histograms of both TEfficiency 00515 // objects are added and a new weight is assigned. The statistic options are not 00516 // changed. 00517 // End_Html 00518 // Begin_Latex #frac{1}{w_{new}} = #frac{1}{w_{1}} + #frac{1}{w_{2}}End_Latex 00519 // Begin_Html 00520 // <b>Example:</b><br /> 00521 // If you use two samples with different numbers of generated events for the same 00522 // process and you want to normalise both to the same integrated luminosity and 00523 // trigger efficiency, the different weights then arise just from the fact that 00524 // you have different numbers of events. The TEfficiency objects should be merged 00525 // because the samples do not represent true alternatives. You expect the same 00526 // result as if you would have a big sample with all events in it. 00527 // End_Html 00528 // Begin_Latex 00529 // 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}}} 00530 // End_Latex 00531 // Begin_Html 00532 // 00533 // <h4><a name="comb">V.2 When should I use combining?</a></h4> 00534 // You should combine TEfficiency objects whenever the weights represent 00535 // alternatives processes for the efficiency. As the combination of two TEfficiency 00536 // objects is not always consistent with the representation by two internal 00537 // histograms, the result is not stored in a TEfficiency object but a TGraphAsymmErrors 00538 // is returned which shows the estimated combined efficiency and its uncertainty 00539 // for each bin. At the moment the combination method <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Combine">Combine </a>only supports combination of 1-dimensional efficiencies in a bayesian approach.<br /> 00540 // For calculating the combined efficiency and its uncertainty for each bin only Bayesian statistics is used. No frequentists methods are presently 00541 // supported for computing the combined efficiency and its confidence interval. 00542 // 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 00543 // for all the TEfficiency objects. 00544 // End_Html 00545 // Begin_Latex 00546 // P_{comb}(#epsilon | {w_{i}}, {k_{i}} , {N_{i}}) = #frac{1}{norm} #prod_{i}{L(k_{i} | N_{i}, #epsilon)}^{w_{i}} #Pi( #epsilon ) 00547 // L(k_{i} | N_{i}, #epsilon) is the likelihood function for the sample i ( a Binomial distribution) 00548 // #Pi( #epsilon) is the prior, a beta distribution B(#epsilon, #alpha, #beta). 00549 // The resulting combined posterior is 00550 // 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) 00551 // #hat{#varepsilon} = #int_{0}^{1} #epsilon #times P_{comb}(#epsilon | {k_{i}} , {N_{i}}) d#epsilon 00552 // confidence level = 1 - #alpha 00553 // #frac{#alpha}{2} = #int_{0}^{#epsilon_{low}} P_{comb}(#epsilon | {k_{i}} , {N_{i}}) d#epsilon ... defines lower boundary 00554 // 1- #frac{#alpha}{2} = #int_{0}^{#epsilon_{up}} P_{comb}(#epsilon | {k_{i}} , {N_{i}}) d#epsilon ... defines upper boundary 00555 // End_Latex 00556 // Begin_Html 00557 // <b>Example:</b><br /> 00558 // If you use cuts to select electrons which can originate from two different 00559 // processes, you can determine the selection efficiency for each process. The 00560 // overall selection efficiency is then the combined efficiency. The weights to be used in the 00561 // combination should be the probability that an 00562 // electron comes from the corresponding process. 00563 // End_Html 00564 // Begin_Latex 00565 // p_{1} = #frac{#sigma_{1}}{#sigma_{1} + #sigma_{2}} = #frac{N_{1}w_{1}}{N_{1}w_{1} + N_{2}w_{2}} 00566 // p_{2} = #frac{#sigma_{2}}{#sigma_{1} + #sigma_{2}} = #frac{N_{2}w_{2}}{N_{1}w_{1} + N_{2}w_{2}} 00567 // End_Latex 00568 // Begin_Html 00569 // <h3><a name="other">VI. Further operations</a></h3> 00570 // 00571 // <h4><a name="histo">VI.1 Information about the internal histograms</a></h4> 00572 // The methods <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetPassedHistogram">GetPassedHistogram</a> and <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetTotalHistogram">GetTotalHistogram</a> 00573 // return a constant pointer to the internal histograms. They can be used to 00574 // obtain information about the internal histograms (e.g. the binning, number of passed / total events in a bin, mean values...). 00575 // One can obtain a clone of the internal histograms by calling <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetCopyPassedHisto">GetCopyPassedHisto</a> or <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetCopyTotalHisto">GetCopyTotalHisto</a>. 00576 // The returned histograms are completely independent from the current 00577 // TEfficiency object. By default, they are not attached to a directory to 00578 // avoid the duplication of data and the user is responsible for deleting them. 00579 // <div class="code"> 00580 // <pre> 00581 // <b>Example:</b> 00582 // //open a root file which contains a TEfficiency object 00583 // TFile* pFile = new TFile("myfile.root","update"); 00584 // 00585 // //get TEfficiency object with name "my_eff" 00586 // TEfficiency* pEff = (TEfficiency*)pFile->Get("my_eff"); 00587 // 00588 // //get clone of total histogram 00589 // TH1* clone = pEff->GetCopyTotalHisto(); 00590 // 00591 // //change clone... 00592 // //save changes of clone directly 00593 // clone->Write(); 00594 // //or append it to the current directoy and write the file 00595 // //clone->SetDirectory(gDirectory); 00596 // //pFile->Wrtie(); 00597 // 00598 // //delete histogram object 00599 // delete clone; 00600 // clone = 0; 00601 // </pre> 00602 // </div><div class="clear" /> 00603 // It is also possible to set the internal total or passed histogram by using the 00604 // methods <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetPassedHistogram">SetPassedHistogram</a> or 00605 // <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:SetTotalHistogram">SetTotalHistogram</a>. 00606 // In order to ensure the validity of the TEfficiency object, the consistency of the 00607 // new histogram and the stored histogram is checked. It sometimes might be 00608 // impossible to change the histograms in a consistent way. Therefore one can force 00609 // the replacement by passing the option "f". Then the user has to ensure that the 00610 // other internal histogram is replaced as well and that the TEfficiency object is 00611 // in a valid state. 00612 // 00613 // <h4><a name="fit">VI.2 Fitting</a></h4> 00614 // The efficiency can be fitted using the <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Fit">Fit</a> function which uses internally the <a href="http://root.cern.ch/root/html/TBinomialEfficiencyFitter.html#TBinomialEfficiencyFitter:Fit">TBinomialEfficiencyFitter::Fit</a> method. 00615 // As this method is using a maximum-likelihood-fit, it is necessary to initialise 00616 // the given fit function with reasonable start values. 00617 // The resulting fit function is attached to the list of associated functions and 00618 // will be drawn automatically during the next <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Draw">Draw</a> command. 00619 // The list of associated function can be modified by using the pointer returned 00620 // by <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:GetListOfFunctions">GetListOfFunctions</a>. 00621 // End_Html 00622 // Begin_Macro(source) 00623 // { 00624 // //canvas only needed for this documentation 00625 // TCanvas* c1 = new TCanvas("example","",600,400); 00626 // c1->SetFillStyle(1001); 00627 // c1->SetFillColor(kWhite); 00628 // 00629 // //create one-dimensional TEfficiency object with fixed bin size 00630 // TEfficiency* pEff = new TEfficiency("eff","my efficiency;x;#epsilon",20,0,10); 00631 // TRandom3 rand; 00632 // 00633 // bool bPassed; 00634 // double x; 00635 // for(int i=0; i<10000; ++i) 00636 // { 00637 // //simulate events with variable under investigation 00638 // x = rand.Uniform(10); 00639 // //check selection: bPassed = DoesEventPassSelection(x) 00640 // bPassed = rand.Rndm() < TMath::Gaus(x,5,4); 00641 // pEff->Fill(bPassed,x); 00642 // } 00643 // 00644 // //create a function for fitting and do the fit 00645 // TF1* f1 = new TF1("f1","gaus",0,10); 00646 // f1->SetParameters(1,5,2); 00647 // pEff->Fit(f1); 00648 // 00649 // //create a threshold function 00650 // TF1* f2 = new TF1("thres","0.8",0,10); 00651 // f2->SetLineColor(kRed); 00652 // //add it to the list of functions 00653 // //use add first because the parameters of the last function will be displayed 00654 // pEff->GetListOfFunctions()->AddFirst(f2); 00655 // 00656 // pEff->Draw("AP"); 00657 // 00658 // //only for this documentation 00659 // return c1; 00660 // } 00661 // End_Macro 00662 // Begin_Html 00663 // 00664 // <h4><a name="draw">VI.3 Draw a TEfficiency object</a></h4> 00665 // A TEfficiency object can be drawn by calling the usual <a href="http://root.cern.ch/root/html/TEfficiency.html#TEfficiency:Draw">Draw</a> method. 00666 // At the moment drawing is only supported for 1- and 2-dimensional TEfficiency 00667 // objects. In the 1-dimensional case you can use the same options as for the <br /> 00668 // <a href="http://root.cern.ch/root/html/TGraph.html#TGraph:Draw">TGraphAsymmErrors::Draw</a> 00669 // method. For 2-dimensional TEfficiency objects you can pass the same options as 00670 // for a <a href="http://root.cern.ch/root/html/TH1.html#TH1:Draw">TH2::Draw</a> object. 00671 // 00672 // <h3 style="margin-bottom:-3em;"><a name="class">VII. TEfficiency class</a></h3> 00673 // 00674 //End_Html 00675 //______________________________________________________________________________ 00676 00677 //______________________________________________________________________________ 00678 TEfficiency::TEfficiency(): 00679 fBeta_alpha(kDefBetaAlpha), 00680 fBeta_beta(kDefBetaBeta), 00681 fBoundary(0), 00682 fConfLevel(kDefConfLevel), 00683 fDirectory(0), 00684 fFunctions(0), 00685 fPaintGraph(0), 00686 fPaintHisto(0), 00687 fPassedHistogram(0), 00688 fTotalHistogram(0), 00689 fWeight(kDefWeight) 00690 { 00691 //default constructor 00692 // 00693 //should not be used explicitly 00694 00695 SetStatisticOption(kDefStatOpt); 00696 } 00697 00698 //______________________________________________________________________________ 00699 TEfficiency::TEfficiency(const TH1& passed,const TH1& total): 00700 fBeta_alpha(kDefBetaAlpha), 00701 fBeta_beta(kDefBetaBeta), 00702 fConfLevel(kDefConfLevel), 00703 fDirectory(0), 00704 fFunctions(new TList()), 00705 fPaintGraph(0), 00706 fPaintHisto(0), 00707 fWeight(kDefWeight) 00708 { 00709 //constructor using two existing histograms as input 00710 // 00711 //Input: passed - contains the events fullfilling some criteria 00712 // total - contains all investigated events 00713 // 00714 //Notes: - both histograms have to fullfill the conditions of CheckConsistency 00715 // - dimension of the resulating efficiency object depends 00716 // on the dimension of the given histograms 00717 // - Clones of both histograms are stored internally 00718 // - The function SetName(total.GetName() + "_clone") is called to set 00719 // the names of the new object and the internal histograms.. 00720 // - The created TEfficiency object is NOT appended to a directory. It 00721 // will not be written to disk during the next TFile::Write() command 00722 // in order to prevent duplication of data. If you want to save this 00723 // TEfficiency object anyway, you can either append it to a 00724 // directory by calling SetDirectory(TDirectory*) or write it 00725 // explicitly to disk by calling Write(). 00726 00727 //check consistency of histograms 00728 if(CheckConsistency(passed,total)) { 00729 Bool_t bStatus = TH1::AddDirectoryStatus(); 00730 TH1::AddDirectory(kFALSE); 00731 fTotalHistogram = (TH1*)total.Clone(); 00732 fPassedHistogram = (TH1*)passed.Clone(); 00733 TH1::AddDirectory(bStatus); 00734 00735 TString newName = total.GetName(); 00736 newName += TString("_clone"); 00737 SetName(newName); 00738 } 00739 else { 00740 Error("TEfficiency(const TH1&,const TH1&)","histograms are not consistent -> results are useless"); 00741 Warning("TEfficiency(const TH1&,const TH1&)","using two empty TH1D('h1','h1',10,0,10)"); 00742 00743 Bool_t bStatus = TH1::AddDirectoryStatus(); 00744 TH1::AddDirectory(kFALSE); 00745 fTotalHistogram = new TH1D("h1_total","h1 (total)",10,0,10); 00746 fPassedHistogram = new TH1D("h1_passed","h1 (passed)",10,0,10); 00747 TH1::AddDirectory(bStatus); 00748 } 00749 00750 SetStatisticOption(kDefStatOpt); 00751 SetDirectory(0); 00752 } 00753 00754 //______________________________________________________________________________ 00755 TEfficiency::TEfficiency(const char* name,const char* title,Int_t nbins, 00756 const Double_t* xbins): 00757 fBeta_alpha(kDefBetaAlpha), 00758 fBeta_beta(kDefBetaBeta), 00759 fConfLevel(kDefConfLevel), 00760 fDirectory(0), 00761 fFunctions(0), 00762 fPaintGraph(0), 00763 fPaintHisto(0), 00764 fWeight(kDefWeight) 00765 { 00766 //create 1-dimensional TEfficiency object with variable bin size 00767 // 00768 //constructor creates two new and empty histograms with a given binning 00769 // 00770 // Input: name - the common part of the name for both histograms (no blanks) 00771 // fTotalHistogram has name: name + "_total" 00772 // fPassedHistogram has name: name + "_passed" 00773 // title - the common part of the title for both histogram 00774 // fTotalHistogram has title: title + " (total)" 00775 // fPassedHistogram has title: title + " (passed)" 00776 // It is possible to label the axis by passing a title with 00777 // the following format: "title;xlabel;ylabel". 00778 // nbins - number of bins on the x-axis 00779 // xbins - array of length (nbins + 1) with low-edges for each bin 00780 // xbins[nbinsx] ... lower edge for overflow bin 00781 00782 Bool_t bStatus = TH1::AddDirectoryStatus(); 00783 TH1::AddDirectory(kFALSE); 00784 fTotalHistogram = new TH1D("total","total",nbins,xbins); 00785 fPassedHistogram = new TH1D("passed","passed",nbins,xbins); 00786 TH1::AddDirectory(bStatus); 00787 00788 Build(name,title); 00789 } 00790 00791 //______________________________________________________________________________ 00792 TEfficiency::TEfficiency(const char* name,const char* title,Int_t nbinsx, 00793 Double_t xlow,Double_t xup): 00794 fBeta_alpha(kDefBetaAlpha), 00795 fBeta_beta(kDefBetaBeta), 00796 fConfLevel(kDefConfLevel), 00797 fDirectory(0), 00798 fFunctions(0), 00799 fPaintGraph(0), 00800 fPaintHisto(0), 00801 fWeight(kDefWeight) 00802 { 00803 //create 1-dimensional TEfficiency object with fixed bins isze 00804 // 00805 //constructor creates two new and empty histograms with a fixed binning 00806 // 00807 // Input: name - the common part of the name for both histograms(no blanks) 00808 // fTotalHistogram has name: name + "_total" 00809 // fPassedHistogram has name: name + "_passed" 00810 // title - the common part of the title for both histogram 00811 // fTotalHistogram has title: title + " (total)" 00812 // fPassedHistogram has title: title + " (passed)" 00813 // It is possible to label the axis by passing a title with 00814 // the following format: "title;xlabel;ylabel". 00815 // nbinsx - number of bins on the x-axis 00816 // xlow - lower edge of first bin 00817 // xup - upper edge of last bin 00818 00819 Bool_t bStatus = TH1::AddDirectoryStatus(); 00820 TH1::AddDirectory(kFALSE); 00821 fTotalHistogram = new TH1D("total","total",nbinsx,xlow,xup); 00822 fPassedHistogram = new TH1D("passed","passed",nbinsx,xlow,xup); 00823 TH1::AddDirectory(bStatus); 00824 00825 Build(name,title); 00826 } 00827 00828 //______________________________________________________________________________ 00829 TEfficiency::TEfficiency(const char* name,const char* title,Int_t nbinsx, 00830 Double_t xlow,Double_t xup,Int_t nbinsy, 00831 Double_t ylow,Double_t yup): 00832 fBeta_alpha(kDefBetaAlpha), 00833 fBeta_beta(kDefBetaBeta), 00834 fConfLevel(kDefConfLevel), 00835 fDirectory(0), 00836 fFunctions(0), 00837 fPaintGraph(0), 00838 fPaintHisto(0), 00839 fWeight(kDefWeight) 00840 { 00841 //create 2-dimensional TEfficiency object with fixed bin size 00842 // 00843 //constructor creates two new and empty histograms with a fixed binning 00844 // 00845 // Input: name - the common part of the name for both histograms(no blanks) 00846 // fTotalHistogram has name: name + "_total" 00847 // fPassedHistogram has name: name + "_passed" 00848 // title - the common part of the title for both histogram 00849 // fTotalHistogram has title: title + " (total)" 00850 // fPassedHistogram has title: title + " (passed)" 00851 // It is possible to label the axis by passing a title with 00852 // the following format: "title;xlabel;ylabel;zlabel". 00853 // nbinsx - number of bins on the x-axis 00854 // xlow - lower edge of first x-bin 00855 // xup - upper edge of last x-bin 00856 // nbinsy - number of bins on the y-axis 00857 // ylow - lower edge of first y-bin 00858 // yup - upper edge of last y-bin 00859 00860 Bool_t bStatus = TH1::AddDirectoryStatus(); 00861 TH1::AddDirectory(kFALSE); 00862 fTotalHistogram = new TH2D("total","total",nbinsx,xlow,xup,nbinsy,ylow,yup); 00863 fPassedHistogram = new TH2D("passed","passed",nbinsx,xlow,xup,nbinsy,ylow,yup); 00864 TH1::AddDirectory(bStatus); 00865 00866 Build(name,title); 00867 } 00868 00869 //______________________________________________________________________________ 00870 TEfficiency::TEfficiency(const char* name,const char* title,Int_t nbinsx, 00871 const Double_t* xbins,Int_t nbinsy, 00872 const Double_t* ybins): 00873 fBeta_alpha(kDefBetaAlpha), 00874 fBeta_beta(kDefBetaBeta), 00875 fConfLevel(kDefConfLevel), 00876 fDirectory(0), 00877 fFunctions(0), 00878 fPaintGraph(0), 00879 fPaintHisto(0), 00880 fWeight(kDefWeight) 00881 { 00882 //create 2-dimensional TEfficiency object with variable bin size 00883 // 00884 //constructor creates two new and empty histograms with a given binning 00885 // 00886 // Input: name - the common part of the name for both histograms(no blanks) 00887 // fTotalHistogram has name: name + "_total" 00888 // fPassedHistogram has name: name + "_passed" 00889 // title - the common part of the title for both histogram 00890 // fTotalHistogram has title: title + " (total)" 00891 // fPassedHistogram has title: title + " (passed)" 00892 // It is possible to label the axis by passing a title with 00893 // the following format: "title;xlabel;ylabel;zlabel". 00894 // nbinsx - number of bins on the x-axis 00895 // xbins - array of length (nbins + 1) with low-edges for each bin 00896 // xbins[nbinsx] ... lower edge for overflow x-bin 00897 // nbinsy - number of bins on the y-axis 00898 // ybins - array of length (nbins + 1) with low-edges for each bin 00899 // ybins[nbinsy] ... lower edge for overflow y-bin 00900 00901 Bool_t bStatus = TH1::AddDirectoryStatus(); 00902 TH1::AddDirectory(kFALSE); 00903 fTotalHistogram = new TH2D("total","total",nbinsx,xbins,nbinsy,ybins); 00904 fPassedHistogram = new TH2D("passed","passed",nbinsx,xbins,nbinsy,ybins); 00905 TH1::AddDirectory(bStatus); 00906 00907 Build(name,title); 00908 } 00909 00910 //______________________________________________________________________________ 00911 TEfficiency::TEfficiency(const char* name,const char* title,Int_t nbinsx, 00912 Double_t xlow,Double_t xup,Int_t nbinsy, 00913 Double_t ylow,Double_t yup,Int_t nbinsz, 00914 Double_t zlow,Double_t zup): 00915 fBeta_alpha(kDefBetaAlpha), 00916 fBeta_beta(kDefBetaBeta), 00917 fConfLevel(kDefConfLevel), 00918 fDirectory(0), 00919 fFunctions(0), 00920 fPaintGraph(0), 00921 fPaintHisto(0), 00922 fWeight(kDefWeight) 00923 { 00924 //create 3-dimensional TEfficiency object with fixed bin size 00925 // 00926 //constructor creates two new and empty histograms with a fixed binning 00927 // 00928 // Input: name - the common part of the name for both histograms(no blanks) 00929 // fTotalHistogram has name: name + "_total" 00930 // fPassedHistogram has name: name + "_passed" 00931 // title - the common part of the title for both histogram 00932 // fTotalHistogram has title: title + " (total)" 00933 // fPassedHistogram has title: title + " (passed)" 00934 // It is possible to label the axis by passing a title with 00935 // the following format: "title;xlabel;ylabel;zlabel". 00936 // nbinsx - number of bins on the x-axis 00937 // xlow - lower edge of first x-bin 00938 // xup - upper edge of last x-bin 00939 // nbinsy - number of bins on the y-axis 00940 // ylow - lower edge of first y-bin 00941 // yup - upper edge of last y-bin 00942 // nbinsz - number of bins on the z-axis 00943 // zlow - lower edge of first z-bin 00944 // zup - upper edge of last z-bin 00945 00946 Bool_t bStatus = TH1::AddDirectoryStatus(); 00947 TH1::AddDirectory(kFALSE); 00948 fTotalHistogram = new TH3D("total","total",nbinsx,xlow,xup,nbinsy,ylow,yup,nbinsz,zlow,zup); 00949 fPassedHistogram = new TH3D("passed","passed",nbinsx,xlow,xup,nbinsy,ylow,yup,nbinsz,zlow,zup); 00950 TH1::AddDirectory(bStatus); 00951 00952 Build(name,title); 00953 } 00954 00955 //______________________________________________________________________________ 00956 TEfficiency::TEfficiency(const char* name,const char* title,Int_t nbinsx, 00957 const Double_t* xbins,Int_t nbinsy, 00958 const Double_t* ybins,Int_t nbinsz, 00959 const Double_t* zbins): 00960 fBeta_alpha(kDefBetaAlpha), 00961 fBeta_beta(kDefBetaBeta), 00962 fConfLevel(kDefConfLevel), 00963 fDirectory(0), 00964 fFunctions(0), 00965 fPaintGraph(0), 00966 fPaintHisto(0), 00967 fWeight(kDefWeight) 00968 { 00969 //create 3-dimensional TEfficiency object with variable bin size 00970 // 00971 //constructor creates two new and empty histograms with a given binning 00972 // 00973 // Input: name - the common part of the name for both histograms(no blanks) 00974 // fTotalHistogram has name: name + "_total" 00975 // fPassedHistogram has name: name + "_passed" 00976 // title - the common part of the title for both histogram 00977 // fTotalHistogram has title: title + " (total)" 00978 // fPassedHistogram has title: title + " (passed)" 00979 // It is possible to label the axis by passing a title with 00980 // the following format: "title;xlabel;ylabel;zlabel". 00981 // nbinsx - number of bins on the x-axis 00982 // xbins - array of length (nbins + 1) with low-edges for each bin 00983 // xbins[nbinsx] ... lower edge for overflow x-bin 00984 // nbinsy - number of bins on the y-axis 00985 // ybins - array of length (nbins + 1) with low-edges for each bin 00986 // xbins[nbinsx] ... lower edge for overflow y-bin 00987 // nbinsz - number of bins on the z-axis 00988 // zbins - array of length (nbins + 1) with low-edges for each bin 00989 // xbins[nbinsx] ... lower edge for overflow z-bin 00990 00991 Bool_t bStatus = TH1::AddDirectoryStatus(); 00992 TH1::AddDirectory(kFALSE); 00993 fTotalHistogram = new TH3D("total","total",nbinsx,xbins,nbinsy,ybins,nbinsz,zbins); 00994 fPassedHistogram = new TH3D("passed","passed",nbinsx,xbins,nbinsy,ybins,nbinsz,zbins); 00995 TH1::AddDirectory(bStatus); 00996 00997 Build(name,title); 00998 } 00999 01000 //______________________________________________________________________________ 01001 TEfficiency::TEfficiency(const TEfficiency& rEff): 01002 TNamed(), 01003 TAttLine(), 01004 TAttFill(), 01005 TAttMarker(), 01006 fBeta_alpha(rEff.fBeta_alpha), 01007 fBeta_beta(rEff.fBeta_beta), 01008 fConfLevel(rEff.fConfLevel), 01009 fDirectory(0), 01010 fFunctions(new TList()), 01011 fPaintGraph(0), 01012 fPaintHisto(0), 01013 fWeight(rEff.fWeight) 01014 { 01015 //copy constructor 01016 // 01017 //The list of associated objects (e.g. fitted functions) is not copied. 01018 // 01019 //Note: - SetName(rEff.GetName() + "_copy") is called to set the names of the 01020 // object and the histograms. 01021 // - The titles are set by calling SetTitle("[copy] " + rEff.GetTitle()). 01022 // - The copied TEfficiency object is NOT appended to a directory. It 01023 // will not be written to disk during the next TFile::Write() command 01024 // in order to prevent duplication of data. If you want to save this 01025 // TEfficiency object anyway, you can either append it to a directory 01026 // by calling SetDirectory(TDirectory*) or write it explicitly to disk 01027 // by calling Write(). 01028 01029 Bool_t bStatus = TH1::AddDirectoryStatus(); 01030 TH1::AddDirectory(kFALSE); 01031 fTotalHistogram = (TH1*)((rEff.fTotalHistogram)->Clone()); 01032 fPassedHistogram = (TH1*)((rEff.fPassedHistogram)->Clone()); 01033 TH1::AddDirectory(bStatus); 01034 01035 TString name = rEff.GetName(); 01036 name += "_copy"; 01037 SetName(name); 01038 TString title = "[copy] "; 01039 title += rEff.GetTitle(); 01040 SetTitle(title); 01041 01042 SetStatisticOption(rEff.GetStatisticOption()); 01043 01044 SetDirectory(0); 01045 01046 //copy style 01047 rEff.TAttLine::Copy(*this); 01048 rEff.TAttFill::Copy(*this); 01049 rEff.TAttMarker::Copy(*this); 01050 } 01051 01052 //______________________________________________________________________________ 01053 TEfficiency::~TEfficiency() 01054 { 01055 //default destructor 01056 01057 //delete all function in fFunctions 01058 if(fFunctions) { 01059 01060 TIter next(fFunctions); 01061 TObject* obj = 0; 01062 while((obj = next())) { 01063 delete obj; 01064 } 01065 01066 fFunctions->Delete(); 01067 } 01068 01069 if(fDirectory) 01070 fDirectory->Remove(this); 01071 01072 delete fFunctions; 01073 delete fTotalHistogram; 01074 delete fPassedHistogram; 01075 delete fPaintGraph; 01076 delete fPaintHisto; 01077 } 01078 01079 //______________________________________________________________________________ 01080 Double_t TEfficiency::AgrestiCoull(Int_t total,Int_t passed,Double_t level,Bool_t bUpper) 01081 { 01082 //calculates the boundaries for the frequentist Agresti-Coull interval 01083 // 01084 //Input: - total : number of total events 01085 // - passed: 0 <= number of passed events <= total 01086 // - level : confidence level 01087 // - bUpper: true - upper boundary is returned 01088 // false - lower boundary is returned 01089 // 01090 //calculation: 01091 //Begin_Latex(separator='=',align='rl') 01092 // #alpha = 1 - #frac{level}{2} 01093 // #kappa = #Phi^{-1}(1 - #alpha,1) ... normal quantile function 01094 // mode = #frac{passed + #frac{#kappa^{2}}{2}}{total + #kappa^{2}} 01095 // #Delta = #kappa * #sqrt{#frac{mode * (1 - mode)}{total + #kappa^{2}}} 01096 // return = max(0,mode - #Delta) or min(1,mode + #Delta) 01097 //End_Latex 01098 01099 Double_t alpha = (1.0 - level)/2; 01100 Double_t kappa = ROOT::Math::normal_quantile(1 - alpha,1); 01101 01102 Double_t mode = (passed + 0.5 * kappa * kappa) / (total + kappa * kappa); 01103 Double_t delta = kappa * std::sqrt(mode * (1 - mode) / (total + kappa * kappa)); 01104 01105 if(bUpper) 01106 return ((mode + delta) > 1) ? 1.0 : (mode + delta); 01107 else 01108 return ((mode - delta) < 0) ? 0.0 : (mode - delta); 01109 } 01110 01111 //______________________________________________________________________________ 01112 Double_t TEfficiency::FeldmanCousins(Int_t total,Int_t passed,Double_t level,Bool_t bUpper) 01113 { 01114 //calculates the boundaries for the frequentist Feldman-Cousins interval 01115 // 01116 //Input: - total : number of total events 01117 // - passed: 0 <= number of passed events <= total 01118 // - level : confidence level 01119 // - bUpper: true - upper boundary is returned 01120 // false - lower boundary is returned 01121 // 01122 // 01123 Double_t lower = 0; 01124 Double_t upper = 1; 01125 if (!FeldmanCousinsInterval(total,passed,level, lower, upper)) { 01126 ::Error("FeldmanCousins","Error running FC method - return 0 or 1"); 01127 } 01128 return (bUpper) ? upper : lower; 01129 } 01130 Bool_t TEfficiency::FeldmanCousinsInterval(Int_t total,Int_t passed,Double_t level,Double_t & lower, Double_t & upper) 01131 { 01132 //calculates the interval boundaries using the frequentist methods of Feldman-Cousins 01133 // 01134 //Input: - total : number of total events 01135 // - passed: 0 <= number of passed events <= total 01136 // - level : confidence level 01137 //Output: 01138 // - lower : lower boundary returned on exit 01139 // - upper : lower boundary returned on exit 01140 // 01141 //Return a flag with the status of the calculation 01142 // 01143 // Calculation: 01144 // The Feldman-Cousins is a frequentist method where the interval is estimated using a Neyman construction where the ordering 01145 // is based on the likelihood ratio: 01146 //Begin_Latex(separator='=',align='rl') 01147 // LR = #frac{Binomial(k | N, #epsilon)}{Binomial(k | N, #hat{#epsilon} ) } 01148 //End_Latex 01149 //See G. J. Feldman and R. D. Cousins, Phys. Rev. D57 (1998) 3873 01150 // and R. D. Cousins, K. E. Hymes, J. Tucker, Nuclear Instruments and Methods in Physics Research A 612 (2010) 388 01151 // 01152 // Implemented using classes developed by Jordan Tucker and Luca Lista 01153 // See File hist/hist/src/TEfficiencyHelper.h 01154 // 01155 FeldmanCousinsBinomialInterval fc; 01156 double alpha = 1.-level; 01157 fc.Init(alpha); 01158 fc.Calculate(passed, total); 01159 lower = fc.Lower(); 01160 upper = fc.Upper(); 01161 return true; 01162 } 01163 01164 //______________________________________________________________________________ 01165 Double_t TEfficiency::Bayesian(Int_t total,Int_t passed,Double_t level,Double_t alpha,Double_t beta,Bool_t bUpper, Bool_t bShortest) 01166 { 01167 //calculates the boundaries for a Bayesian confidence interval (shortest or central interval depending on the option) 01168 // 01169 //Input: - total : number of total events 01170 // - passed: 0 <= number of passed events <= total 01171 // - level : confidence level 01172 // - alpha : shape parameter > 0 for the prior distribution (fBeta_alpha) 01173 // - beta : shape parameter > 0 for the prior distribution (fBeta_beta) 01174 // - bUpper: true - upper boundary is returned 01175 // false - lower boundary is returned 01176 // 01177 //Note: In the case central confidence interval is calculated. 01178 // when passed = 0 (or passed = total) the lower (or upper) 01179 // interval values will be larger than 0 (or smaller than 1). 01180 // 01181 //Calculation: 01182 // 01183 //The posterior probability in bayesian statistics is given by: 01184 //Begin_Latex P(#varepsilon |k,N) #propto L(#varepsilon|k,N) #times Prior(#varepsilon)End_Latex 01185 //As an efficiency can be interpreted as probability of a positive outcome of 01186 //a Bernoullli trial the likelihood function is given by the binomial 01187 //distribution: 01188 //Begin_Latex L(#varepsilon|k,N) = Binomial(N,k) #varepsilon ^{k} (1 - #varepsilon)^{N-k}End_Latex 01189 //At the moment only beta distributions are supported as prior probabilities 01190 //of the efficiency (Begin_Latex #scale[0.8]{B(#alpha,#beta)}End_Latex is the beta function): 01191 //Begin_Latex Prior(#varepsilon) = #frac{1}{B(#alpha,#beta)} #varepsilon ^{#alpha - 1} (1 - #varepsilon)^{#beta - 1}End_Latex 01192 //The posterior probability is therefore again given by a beta distribution: 01193 //Begin_Latex P(#varepsilon |k,N) #propto #varepsilon ^{k + #alpha - 1} (1 - #varepsilon)^{N - k + #beta - 1} End_Latex 01194 //In case of central intervals 01195 //the lower boundary for the equal-tailed confidence interval is given by the 01196 //inverse cumulative (= quantile) function for the quantile Begin_Latex #frac{1 - level}{2} End_Latex. 01197 //The upper boundary for the equal-tailed confidence interval is given by the 01198 //inverse cumulative (= quantile) function for the quantile Begin_Latex #frac{1 + level}{2} End_Latex. 01199 //Hence it is the solution Begin_Latex #varepsilon End_Latex of the following equation: 01200 //Begin_Latex 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} End_Latex 01201 // In the case of shortest interval the minimum interval aorund the mode is found by minimizing the length of all intervals whith the 01202 // given probability content. See TEfficiency::BetaShortestInterval 01203 01204 Double_t a = double(passed)+alpha; 01205 Double_t b = double(total-passed)+beta; 01206 01207 if (bShortest) { 01208 double lower = 0; 01209 double upper = 1; 01210 BetaShortestInterval(level,a,b,lower,upper); 01211 return (bUpper) ? upper : lower; 01212 } 01213 else 01214 return BetaCentralInterval(level, a, b, bUpper); 01215 } 01216 //______________________________________________________________________________ 01217 Double_t TEfficiency::BetaCentralInterval(Double_t level,Double_t a,Double_t b,Bool_t bUpper) 01218 { 01219 //calculates the boundaries for a central confidence interval for a Beta distribution 01220 // 01221 //Input: - level : confidence level 01222 // - a : parameter > 0 for the beta distribution (for a posterior is passed + prior_alpha 01223 // - b : parameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta 01224 // - bUpper: true - upper boundary is returned 01225 // false - lower boundary is returned 01226 // 01227 01228 if(bUpper) { 01229 if((a > 0) && (b > 0)) 01230 return ROOT::Math::beta_quantile((1+level)/2,a,b); 01231 else { 01232 gROOT->Error("TEfficiency::BayesianCentral","Invalid input parameters - return 1"); 01233 return 1; 01234 } 01235 } 01236 else { 01237 if((a > 0) && (b > 0)) 01238 return ROOT::Math::beta_quantile((1-level)/2,a,b); 01239 else { 01240 gROOT->Error("TEfficiency::BayesianCentral","Invalid input parameters - return 0"); 01241 return 0; 01242 } 01243 } 01244 } 01245 01246 struct Beta_interval_length { 01247 Beta_interval_length(Double_t level,Double_t alpha,Double_t beta ) : 01248 fCL(level), fAlpha(alpha), fBeta(beta) 01249 {} 01250 01251 Double_t LowerMax() { 01252 // max allowed value of lower given the interval size 01253 return ROOT::Math::beta_quantile_c(fCL, fAlpha,fBeta); 01254 } 01255 01256 Double_t operator() (double lower) const { 01257 // return length of interval 01258 Double_t plow = ROOT::Math::beta_cdf(lower, fAlpha, fBeta); 01259 Double_t pup = plow + fCL; 01260 double upper = ROOT::Math::beta_quantile(pup, fAlpha,fBeta); 01261 return upper-lower; 01262 } 01263 Double_t fCL; // interval size (confidence level) 01264 Double_t fAlpha; // beta distribution alpha parameter 01265 Double_t fBeta; // beta distribution beta parameter 01266 01267 }; 01268 01269 //______________________________________________________________________________ 01270 Bool_t TEfficiency::BetaShortestInterval(Double_t level,Double_t a,Double_t b, Double_t & lower, Double_t & upper) 01271 { 01272 //calculates the boundaries for a shortest confidence interval for a Beta distribution 01273 // 01274 //Input: - level : confidence level 01275 // - a : parameter > 0 for the beta distribution (for a posterior is passed + prior_alpha 01276 // - b : parameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta 01277 // - bUpper: true - upper boundary is returned 01278 // false - lower boundary is returned 01279 // 01280 // 01281 //The lower/upper boundary are then obtained by finding the shortest interval of the beta distribbution 01282 // contained the desired probability level. 01283 // The length of all possible intervals is minimized in order to find the shortest one 01284 01285 if (a <= 0 || b <= 0) { 01286 lower = 0; upper = 1; 01287 gROOT->Error("TEfficiency::BayesianShortest","Invalid input parameters - return [0,1]"); 01288 return kFALSE; 01289 } 01290 01291 // treat here special cases when mode == 0 or 1 01292 double mode = BetaMode(a,b); 01293 if (mode == 0.0) { 01294 lower = 0; 01295 upper = ROOT::Math::beta_quantile(level, a, b); 01296 return kTRUE; 01297 } 01298 if (mode == 1.0) { 01299 lower = ROOT::Math::beta_quantile_c(level, a, b); 01300 upper = 1.0; 01301 return kTRUE; 01302 } 01303 // special case when the shortest interval is undefined return the central interval 01304 // can happen for a posterior when passed=total=0 01305 // 01306 if ( a==b && a<1.0) { 01307 lower = BetaCentralInterval(level,a,b,kFALSE); 01308 upper = BetaCentralInterval(level,a,b,kTRUE); 01309 return kTRUE; 01310 } 01311 01312 // for the other case perform a minimization 01313 // make a function of the length of the posterior interval as a function of lower bound 01314 Beta_interval_length intervalLength(level,a,b); 01315 // minimize the interval length 01316 ROOT::Math::WrappedFunction<const Beta_interval_length &> func(intervalLength); 01317 ROOT::Math::BrentMinimizer1D minim; 01318 minim.SetFunction(func, 0, intervalLength.LowerMax() ); 01319 minim.SetNpx(2); // no need to bracket with many iterations. Just do few times to estimate some better points 01320 bool ret = minim.Minimize(100, 1.E-10,1.E-10); 01321 if (!ret) { 01322 gROOT->Error("TEfficiency::BayesianShortes","Error finding the shortest interval"); 01323 return kFALSE; 01324 } 01325 lower = minim.XMinimum(); 01326 upper = lower + minim.FValMinimum(); 01327 return kTRUE; 01328 } 01329 01330 //______________________________________________________________________________ 01331 Double_t TEfficiency::BetaMean(Double_t a,Double_t b) 01332 { 01333 // compute the mean (average) of the beta distribution 01334 // 01335 //Input: a : parameter > 0 for the beta distribution (for a posterior is passed + prior_alpha 01336 // b : parameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta 01337 // 01338 01339 if (a <= 0 || b <= 0 ) { 01340 gROOT->Error("TEfficiency::BayesianMean","Invalid input parameters - return 0"); 01341 return 0; 01342 } 01343 01344 Double_t mean = a / (a + b); 01345 return mean; 01346 } 01347 01348 //______________________________________________________________________________ 01349 Double_t TEfficiency::BetaMode(Double_t a,Double_t b) 01350 { 01351 // compute the mode of the beta distribution 01352 // 01353 //Input: a : parameter > 0 for the beta distribution (for a posterior is passed + prior_alpha 01354 // b : parameter > 0 for the beta distribution (for a posterior is (total-passed) + prior_beta 01355 // 01356 // note the mode is defined for a Beta(a,b) only if (a,b)>1 (a = passed+alpha; b = total-passed+beta) 01357 // return then the following in case (a,b) < 1: 01358 // if (a==b) return 0.5 (it is really undefined) 01359 // if (a < b) return 0; 01360 // if (a > b) return 1; 01361 01362 if (a <= 0 || b <= 0 ) { 01363 gROOT->Error("TEfficiency::BayesianMode","Invalid input parameters - return 0"); 01364 return 0; 01365 } 01366 if ( a <= 1 || b <= 1) { 01367 if ( a < b) return 0; 01368 if ( a > b) return 1; 01369 if (a == b) return 0.5; // cannot do otherwise 01370 } 01371 01372 // since a and b are > 1 here denominator cannot be 0 or < 0 01373 Double_t mode = (a - 1.0) / (a + b -2.0); 01374 return mode; 01375 } 01376 //______________________________________________________________________________ 01377 void TEfficiency::Build(const char* name,const char* title) 01378 { 01379 //building standard data structure of a TEfficiency object 01380 // 01381 //Notes: - calls: SetName(name), SetTitle(title) 01382 // - set the statistic option to the default (kFCP) 01383 // - appends this object to the current directory 01384 // SetDirectory(gDirectory) 01385 01386 SetName(name); 01387 SetTitle(title); 01388 01389 SetStatisticOption(kDefStatOpt); 01390 SetDirectory(gDirectory); 01391 01392 //set normalisation factors to 0, otherwise the += may not work properly 01393 fPassedHistogram->SetNormFactor(0); 01394 fTotalHistogram->SetNormFactor(0); 01395 01396 fFunctions = new TList(); 01397 } 01398 01399 //______________________________________________________________________________ 01400 Bool_t TEfficiency::CheckBinning(const TH1& pass,const TH1& total) 01401 { 01402 //checks binning for each axis 01403 // 01404 //It is assumed that the passed histograms have the same dimension. 01405 01406 TAxis* ax1 = 0; 01407 TAxis* ax2 = 0; 01408 01409 //check binning along axis 01410 for(Int_t j = 0; j < pass.GetDimension(); ++j) { 01411 switch(j) { 01412 case 0: 01413 ax1 = pass.GetXaxis(); 01414 ax2 = total.GetXaxis(); 01415 break; 01416 case 1: 01417 ax1 = pass.GetYaxis(); 01418 ax2 = total.GetYaxis(); 01419 break; 01420 case 2: 01421 ax1 = pass.GetZaxis(); 01422 ax2 = total.GetZaxis(); 01423 break; 01424 } 01425 01426 if(ax1->GetNbins() != ax2->GetNbins()) 01427 return false; 01428 01429 for(Int_t i = 1; i <= ax1->GetNbins() + 1; ++i) 01430 if(!TMath::AreEqualRel(ax1->GetBinLowEdge(i), ax2->GetBinLowEdge(i), 1.E-15)) 01431 return false; 01432 01433 if(!TMath::AreEqualRel(ax1->GetXmax(), ax2->GetXmax(), 1.E-15)) 01434 return false; 01435 01436 01437 } 01438 01439 return true; 01440 } 01441 01442 //______________________________________________________________________________ 01443 Bool_t TEfficiency::CheckConsistency(const TH1& pass,const TH1& total,Option_t* opt) 01444 { 01445 //checks the consistence of the given histograms 01446 // 01447 //The histograms are considered as consistent if: 01448 //- both have the same dimension 01449 //- both have the same binning 01450 //- pass.GetBinContent(i) <= total.GetBinContent(i) for each bin i 01451 // 01452 //Option: - w: The check for unit weights is skipped and therefore histograms 01453 // filled with weights are accepted. 01454 01455 if(pass.GetDimension() != total.GetDimension()) { 01456 gROOT->Error("TEfficiency::CheckConsistency","passed TEfficiency objects have different dimensions"); 01457 return false; 01458 } 01459 01460 if(!CheckBinning(pass,total)) { 01461 gROOT->Error("TEfficiency::CheckConsistency","passed TEfficiency objects have different binning"); 01462 return false; 01463 } 01464 01465 if(!CheckEntries(pass,total,opt)) { 01466 gROOT->Error("TEfficiency::CheckConsistency","passed TEfficiency objects do not have consistent bin contents"); 01467 return false; 01468 } 01469 01470 return true; 01471 } 01472 01473 //______________________________________________________________________________ 01474 Bool_t TEfficiency::CheckEntries(const TH1& pass,const TH1& total,Option_t* opt) 01475 { 01476 //checks whether bin contents are compatible with binomial statistics 01477 // 01478 //The following inequality has to be valid for each bin i: 01479 // total.GetBinContent(i) >= pass.GetBinContent(i) 01480 // 01481 //and the histogram have to be filled with unit weights. 01482 // 01483 //Option: - w: Do not check for unit weights -> accept histograms filled with 01484 // weights 01485 // 01486 //Note: - It is assumed that both histograms have the same dimension and 01487 // binning. 01488 01489 TString option = opt; 01490 option.ToLower(); 01491 01492 //check for unit weights 01493 if(!option.Contains("w")) { 01494 Double_t statpass[10]; 01495 Double_t stattotal[10]; 01496 01497 pass.GetStats(statpass); 01498 total.GetStats(stattotal); 01499 01500 //require: sum of weights == sum of weights^2 01501 if((TMath::Abs(statpass[0]-statpass[1]) > 1e-5) || 01502 (TMath::Abs(stattotal[0]-stattotal[1]) > 1e-5)) 01503 return false; 01504 } 01505 01506 //check: pass <= total 01507 Int_t nbinsx, nbinsy, nbinsz, nbins; 01508 01509 nbinsx = pass.GetNbinsX(); 01510 nbinsy = pass.GetNbinsY(); 01511 nbinsz = pass.GetNbinsZ(); 01512 01513 switch(pass.GetDimension()) { 01514 case 1: nbins = nbinsx + 2; break; 01515 case 2: nbins = (nbinsx + 2) * (nbinsy + 2); break; 01516 case 3: nbins = (nbinsx + 2) * (nbinsy + 2) * (nbinsz + 2); break; 01517 default: nbins = 0; 01518 } 01519 01520 for(Int_t i = 0; i < nbins; ++i) { 01521 if(pass.GetBinContent(i) > total.GetBinContent(i)) 01522 return false; 01523 } 01524 01525 return true; 01526 } 01527 01528 //______________________________________________________________________________ 01529 Double_t TEfficiency::ClopperPearson(Int_t total,Int_t passed,Double_t level,Bool_t bUpper) 01530 { 01531 //calculates the boundaries for the frequentist Clopper-Pearson interval 01532 // 01533 //This interval is recommended by the PDG. 01534 // 01535 //Input: - total : number of total events 01536 // - passed: 0 <= number of passed events <= total 01537 // - level : confidence level 01538 // - bUpper: true - upper boundary is returned 01539 // false - lower boundary is returned 01540 // 01541 //calculation: 01542 // 01543 //The lower boundary of the Clopper-Pearson interval is the "exact" inversion 01544 //of the test: 01545 //Begin_Latex(separator='=',align='rl') 01546 //P(x #geq passed; total) = #frac{1 - level}{2} 01547 //P(x #geq passed; total) = 1 - P(x #leq passed - 1; total) 01548 // = 1 - #frac{1}{norm} * #int_{0}^{1 - #varepsilon} t^{total - passed} (1 - t)^{passed - 1} dt 01549 // = 1 - #frac{1}{norm} * #int_{#varepsilon}^{1} t^{passed - 1} (1 - t)^{total - passed} dt 01550 // = #frac{1}{norm} * #int_{0}^{#varepsilon} t^{passed - 1} (1 - t)^{total - passed} dt 01551 // = I_{#varepsilon}(passed,total - passed + 1) 01552 //End_Latex 01553 //The lower boundary is therfore given by the Begin_Latex #frac{1 - level}{2}End_Latex quantile 01554 //of the beta distribution. 01555 // 01556 //The upper boundary of the Clopper-Pearson interval is the "exact" inversion 01557 //of the test: 01558 //Begin_Latex(separator='=',align='rl') 01559 //P(x #leq passed; total) = #frac{1 - level}{2} 01560 //P(x #leq passed; total) = #frac{1}{norm} * #int_{0}^{1 - #varepsilon} t^{total - passed - 1} (1 - t)^{passed} dt 01561 // = #frac{1}{norm} * #int_{#varepsilon}^{1} t^{passed} (1 - t)^{total - passed - 1} dt 01562 // = 1 - #frac{1}{norm} * #int_{0}^{#varepsilon} t^{passed} (1 - t)^{total - passed - 1} dt 01563 // #Rightarrow 1 - #frac{1 - level}{2} = #frac{1}{norm} * #int_{0}^{#varepsilon} t^{passed} (1 - t)^{total - passed -1} dt 01564 // #frac{1 + level}{2} = I_{#varepsilon}(passed + 1,total - passed) 01565 //End_Latex 01566 //The upper boundary is therfore given by the Begin_Latex #frac{1 + level}{2}End_Latex quantile 01567 //of the beta distribution. 01568 // 01569 //Note: The connection between the binomial distribution and the regularized 01570 // incomplete beta function Begin_Latex I_{#varepsilon}(#alpha,#beta)End_Latex has been used. 01571 01572 Double_t alpha = (1.0 - level) / 2; 01573 if(bUpper) 01574 return ((passed == total) ? 1.0 : ROOT::Math::beta_quantile(1 - alpha,passed + 1,total-passed)); 01575 else 01576 return ((passed == 0) ? 0.0 : ROOT::Math::beta_quantile(alpha,passed,total-passed+1.0)); 01577 } 01578 //______________________________________________________________________________ 01579 Double_t TEfficiency::Combine(Double_t& up,Double_t& low,Int_t n, 01580 const Int_t* pass,const Int_t* total, 01581 Double_t alpha, Double_t beta, 01582 Double_t level,const Double_t* w,Option_t* opt) 01583 { 01584 //calculates the combined efficiency and its uncertainties 01585 // 01586 //This method does a bayesian combination of the given samples. 01587 // 01588 //Input: 01589 //- up : contains the upper limit of the confidence interval afterwards 01590 //- low : contains the lower limit of the confidence interval afterwards 01591 //- n : number of samples which are combined 01592 //- pass : array of length n containing the number of passed events 01593 //- total : array of length n containing the corresponding numbers of total 01594 // events 01595 //- alpha : shape parameters for the beta distribution as prior 01596 //- beta : shape parameters for the beta distribution as prior 01597 //- level : desired confidence level 01598 //- w : weights for each sample; if not given, all samples get the weight 1 01599 // The weights do not need to be normalized, since they are internally renormalized 01600 // to the number of effective entries. 01601 //- options: 01602 // 01603 // + mode : The mode is returned instead of the mean of the posterior as best value 01604 // When using the mode the shortest interval is also computed instead of the central one 01605 // + shortest: compute shortest interval (done by default if mode option is set) 01606 // + central: compute central interval (done by default if mode option is NOT set) 01607 // 01608 //Begin_Html 01609 //Calculation: 01610 //<ol> 01611 //<li>The combined posterior distributions is calculated from the Bayes theorem assuming a common prior Beta distribution. 01612 // It is easy to proof that the combined posterior is then:</li> 01613 //Begin_Latex(separator='=',align='rl') 01614 //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) 01615 //w_{i} = weight for each sample renormalized to the effective entries 01616 //w^{'}_{i} = w_{i} #frac{ #sum_{i} {w_{i} } } { #sum_{i} {w_{i}^{2} } } 01617 //End_Latex 01618 //Begin_Html 01619 //<li>The estimated efficiency is the mode (or the mean) of the obtained posterior distribution </li> 01620 //End_Html 01621 //Begin_Html 01622 //<li>The boundaries of the confidence interval for a confidence level (1 - a) 01623 //are given by the a/2 and 1-a/2 quantiles of the resulting cumulative 01624 //distribution.</li> 01625 //</ol> 01626 //End_Html 01627 //Example (uniform prior distribution): 01628 //Begin_Macro(source) 01629 //{ 01630 // TCanvas* c1 = new TCanvas("c1","",600,800); 01631 // c1->Divide(1,2); 01632 // c1->SetFillStyle(1001); 01633 // c1->SetFillColor(kWhite); 01634 // 01635 // TF1* p1 = new TF1("p1","TMath::BetaDist(x,19,9)",0,1); 01636 // TF1* p2 = new TF1("p2","TMath::BetaDist(x,4,8)",0,1); 01637 // TF1* comb = new TF1("comb2","TMath::BetaDist(x,[0],[1])",0,1); 01638 // double norm = 1./(0.6*0.6+0.4*0.4); // weight normalization 01639 // double a = 0.6*18.0 + 0.4*3.0 + 1.0; // new alpha parameter of combined beta dist. 01640 // double b = 0.6*10+0.4*7+1.0; // new beta parameter of combined beta dist. 01641 // comb->SetParameters(norm*a ,norm *b ); 01642 // TF1* const1 = new TF1("const1","0.05",0,1); 01643 // TF1* const2 = new TF1("const2","0.95",0,1); 01644 // 01645 // p1->SetLineColor(kRed); 01646 // p1->SetTitle("combined posteriors;#epsilon;P(#epsilon|k,N)"); 01647 // p2->SetLineColor(kBlue); 01648 // comb->SetLineColor(kGreen+2); 01649 // 01650 // TLegend* leg1 = new TLegend(0.12,0.65,0.5,0.85); 01651 // leg1->AddEntry(p1,"k1 = 18, N1 = 26","l"); 01652 // leg1->AddEntry(p2,"k2 = 3, N2 = 10","l"); 01653 // leg1->AddEntry(comb,"combined: p1 = 0.6, p2=0.4","l"); 01654 // 01655 // c1->cd(1); 01656 // comb->Draw(); 01657 // p1->Draw("same"); 01658 // p2->Draw("same"); 01659 // leg1->Draw("same"); 01660 // c1->cd(2); 01661 // const1->SetLineWidth(1); 01662 // const2->SetLineWidth(1); 01663 // TGraph* gr = (TGraph*)comb->DrawIntegral(); 01664 // gr->SetTitle("cumulative function of combined posterior with boundaries for cl = 95%;#epsilon;CDF"); 01665 // const1->Draw("same"); 01666 // const2->Draw("same"); 01667 // 01668 // c1->cd(0); 01669 // return c1; 01670 //} 01671 //End_Macro 01672 01673 TString option(opt); 01674 option.ToLower(); 01675 01676 //LM: new formula for combination 01677 // works only if alpha beta are the same always 01678 // the weights are normalized to w(i) -> N_eff w(i)/ Sum w(i) 01679 // i.e. w(i) -> Sum (w(i) / Sum (w(i)^2) * w(i) 01680 // norm = Sum (w(i) / Sum (w(i)^2) 01681 double ntot = 0; 01682 double ktot = 0; 01683 double sumw = 0; 01684 double sumw2 = 0; 01685 for (int i = 0; i < n ; ++i) { 01686 if(pass[i] > total[i]) { 01687 ::Error("TEfficiency::Combine","total events = %i < passed events %i",total[i],pass[i]); 01688 ::Info("TEfficiency::Combine","stop combining"); 01689 return -1; 01690 } 01691 01692 ntot += w[i] * total[i]; 01693 ktot += w[i] * pass[i]; 01694 sumw += w[i]; 01695 sumw2 += w[i]*w[i]; 01696 //mean += w[i] * (pass[i] + alpha[i])/(total[i] + alpha[i] + beta[i]); 01697 } 01698 double norm = sumw/sumw2; 01699 ntot *= norm; 01700 ktot *= norm; 01701 if(ktot > ntot) { 01702 ::Error("TEfficiency::Combine","total = %f < passed %f",ntot,ktot); 01703 ::Info("TEfficiency::Combine","stop combining"); 01704 return -1; 01705 } 01706 01707 double a = ktot + alpha; 01708 double b = ntot - ktot + beta; 01709 01710 double mean = a/(a+b); 01711 double mode = BetaMode(a,b); 01712 01713 01714 Bool_t shortestInterval = option.Contains("sh") || ( option.Contains("mode") && !option.Contains("cent") ); 01715 01716 if (shortestInterval) 01717 BetaShortestInterval(level, a, b, low, up); 01718 else { 01719 low = BetaCentralInterval(level, a, b, false); 01720 up = BetaCentralInterval(level, a, b, true); 01721 } 01722 01723 if (option.Contains("mode")) return mode; 01724 return mean; 01725 01726 } 01727 //______________________________________________________________________________ 01728 TGraphAsymmErrors* TEfficiency::Combine(TCollection* pList,Option_t* option, 01729 Int_t n,const Double_t* w) 01730 { 01731 //combines a list of 1-dimensional TEfficiency objects 01732 // 01733 //A TGraphAsymmErrors object is returned which contains the estimated 01734 //efficiency and its uncertainty for each bin. 01735 //If the combination fails, a zero pointer is returned. 01736 // 01737 //At the moment the combining is only implemented for bayesian statistics. 01738 // 01739 //Input: 01740 //- pList : list containing TEfficiency objects which should be combined 01741 // only one-dimensional efficiencies are taken into account 01742 //- options 01743 // + s : strict combining; only TEfficiency objects with the same beta 01744 // prior and the flag kIsBayesian == true are combined 01745 // If not specified the prior parameter of the first TEfficiency object is used 01746 // + v : verbose mode; print information about combining 01747 // + cl=x : set confidence level (0 < cl < 1). If not specified, the 01748 // confidence level of the first TEfficiency object is used. 01749 // + mode Use mode of combined posterior as estimated value for the efficiency 01750 // + shortest: compute shortest interval (done by default if mode option is set) 01751 // + central: compute central interval (done by default if mode option is NOT set) 01752 // 01753 //- n : number of weights (has to be the number of one-dimensional 01754 // TEfficiency objects in pList) 01755 // If no weights are passed, the internal weights GetWeight() of 01756 // the given TEfficiency objects are used. 01757 //- w : array of length n with weights for each TEfficiency object in 01758 // pList (w[0] correspond to pList->First ... w[n-1] -> pList->Last) 01759 // The weights do not have to be normalised. 01760 // 01761 //For each bin the calculation is done by the Combine(double&, double& ...) method. 01762 01763 TString opt = option; 01764 opt.ToLower(); 01765 01766 //parameter of prior distribution, confidence level and normalisation factor 01767 Double_t alpha = -1; 01768 Double_t beta = -1; 01769 Double_t level = 0; 01770 01771 //flags for combining 01772 Bool_t bStrict = false; 01773 Bool_t bOutput = false; 01774 Bool_t bWeights = false; 01775 //list of all information needed to weight and combine efficiencies 01776 std::vector<TH1*> vTotal; vTotal.reserve(n); 01777 std::vector<TH1*> vPassed; vPassed.reserve(n); 01778 std::vector<Double_t> vWeights; vWeights.reserve(n); 01779 // std::vector<Double_t> vAlpha; 01780 // std::vector<Double_t> vBeta; 01781 01782 if(opt.Contains("s")) { 01783 opt.ReplaceAll("s",""); 01784 bStrict = true; 01785 } 01786 01787 if(opt.Contains("v")) { 01788 opt.ReplaceAll("v",""); 01789 bOutput = true; 01790 } 01791 01792 if(opt.Contains("cl=")) { 01793 Ssiz_t pos = opt.Index("cl=") + 3; 01794 level = atof( opt(pos,opt.Length() ).Data() ); 01795 if((level <= 0) || (level >= 1)) 01796 level = 0; 01797 opt.ReplaceAll("cl=",""); 01798 } 01799 01800 //are weights explicitly given 01801 if(n && w) { 01802 bWeights = true; 01803 for(Int_t k = 0; k < n; ++k) { 01804 if(w[k] > 0) 01805 vWeights.push_back(w[k]); 01806 else { 01807 gROOT->Error("TEfficiency::Combine","invalid custom weight found w = %.2lf",w[k]); 01808 gROOT->Info("TEfficiency::Combine","stop combining"); 01809 return 0; 01810 } 01811 } 01812 } 01813 01814 TIter next(pList); 01815 TObject* obj = 0; 01816 TEfficiency* pEff = 0; 01817 while((obj = next())) { 01818 pEff = dynamic_cast<TEfficiency*>(obj); 01819 //is object a TEfficiency object? 01820 if(pEff) { 01821 if(pEff->GetDimension() > 1) 01822 continue; 01823 if(!level) level = pEff->GetConfidenceLevel(); 01824 01825 if(alpha<1) alpha = pEff->GetBetaAlpha(); 01826 if(beta<1) beta = pEff->GetBetaBeta(); 01827 01828 //if strict combining, check priors, confidence level and statistic 01829 if(bStrict) { 01830 if(alpha != pEff->GetBetaAlpha()) 01831 continue; 01832 if(beta != pEff->GetBetaBeta()) 01833 continue; 01834 if(!pEff->UsesBayesianStat()) 01835 continue; 01836 } 01837 01838 vTotal.push_back(pEff->fTotalHistogram); 01839 vPassed.push_back(pEff->fPassedHistogram); 01840 01841 //no weights given -> use weights of TEfficiency objects 01842 if(!bWeights) 01843 vWeights.push_back(pEff->fWeight); 01844 01845 //strict combining -> using global prior 01846 // if(bStrict) { 01847 // vAlpha.push_back(alpha); 01848 // vBeta.push_back(beta); 01849 // } 01850 // else { 01851 // vAlpha.push_back(pEff->GetBetaAlpha()); 01852 // vBeta.push_back(pEff->GetBetaBeta()); 01853 // } 01854 } 01855 } 01856 01857 //no TEfficiency objects found 01858 if(vTotal.empty()) { 01859 gROOT->Error("TEfficiency::Combine","no TEfficiency objects in given list"); 01860 gROOT->Info("TEfficiency::Combine","stop combining"); 01861 return 0; 01862 } 01863 01864 //invalid number of custom weights 01865 if(bWeights && (n != (Int_t)vTotal.size())) { 01866 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()); 01867 gROOT->Info("TEfficiency::Combine","stop combining"); 01868 return 0; 01869 } 01870 01871 Int_t nbins_max = vTotal.at(0)->GetNbinsX(); 01872 //check binning of all histograms 01873 for(UInt_t i=0; i<vTotal.size(); ++i) { 01874 try { 01875 TEfficiency::CheckBinning(*vTotal.at(0),*vTotal.at(i)); 01876 } 01877 catch(std::exception&) { 01878 gROOT->Warning("TEfficiency::Combine","histograms have not the same binning -> results may be useless"); 01879 } 01880 if(vTotal.at(i)->GetNbinsX() < nbins_max) nbins_max = vTotal.at(i)->GetNbinsX(); 01881 } 01882 01883 //display information about combining 01884 if(bOutput) { 01885 gROOT->Info("TEfficiency::Combine","combining %i TEfficiency objects",(Int_t)vTotal.size()); 01886 if(bWeights) 01887 gROOT->Info("TEfficiency::Combine","using custom weights"); 01888 if(bStrict) { 01889 gROOT->Info("TEfficiency::Combine","using the following prior probability for the efficiency: P(e) ~ Beta(e,%.3lf,%.3lf)",alpha,beta); 01890 } 01891 else 01892 gROOT->Info("TEfficiency::Combine","using individual priors of each TEfficiency object"); 01893 gROOT->Info("TEfficiency::Combine","confidence level = %.2lf",level); 01894 } 01895 01896 //create TGraphAsymmErrors with efficiency 01897 std::vector<Double_t> x(nbins_max); 01898 std::vector<Double_t> xlow(nbins_max); 01899 std::vector<Double_t> xhigh(nbins_max); 01900 std::vector<Double_t> eff(nbins_max); 01901 std::vector<Double_t> efflow(nbins_max); 01902 std::vector<Double_t> effhigh(nbins_max); 01903 01904 //parameters for combining: 01905 //number of objects 01906 Int_t num = vTotal.size(); 01907 std::vector<Int_t> pass(num); 01908 std::vector<Int_t> total(num); 01909 01910 //loop over all bins 01911 Double_t low, up; 01912 for(Int_t i=1; i <= nbins_max; ++i) { 01913 //the binning of the x-axis is taken from the first total histogram 01914 x[i-1] = vTotal.at(0)->GetBinCenter(i); 01915 xlow[i-1] = x[i-1] - vTotal.at(0)->GetBinLowEdge(i); 01916 xhigh[i-1] = vTotal.at(0)->GetBinWidth(i) - xlow[i-1]; 01917 01918 for(Int_t j = 0; j < num; ++j) { 01919 pass[j] = (Int_t)(vPassed.at(j)->GetBinContent(i) + 0.5); 01920 total[j] = (Int_t)(vTotal.at(j)->GetBinContent(i) + 0.5); 01921 } 01922 01923 //fill efficiency and errors 01924 eff[i-1] = Combine(up,low,num,&pass[0],&total[0],alpha,beta,level,&vWeights[0],opt.Data()); 01925 //did an error occured ? 01926 if(eff[i-1] == -1) { 01927 gROOT->Error("TEfficiency::Combine","error occured during combining"); 01928 gROOT->Info("TEfficiency::Combine","stop combining"); 01929 return 0; 01930 } 01931 efflow[i-1]= eff[i-1] - low; 01932 effhigh[i-1]= up - eff[i-1]; 01933 }//loop over all bins 01934 01935 TGraphAsymmErrors* gr = new TGraphAsymmErrors(nbins_max,&x[0],&eff[0],&xlow[0],&xhigh[0],&efflow[0],&effhigh[0]); 01936 01937 return gr; 01938 } 01939 01940 //______________________________________________________________________________ 01941 Int_t TEfficiency::DistancetoPrimitive(Int_t px, Int_t py) 01942 { 01943 // Compute distance from point px,py to a graph. 01944 // 01945 // Compute the closest distance of approach from point px,py to this line. 01946 // The distance is computed in pixels units. 01947 // 01948 // Forward the call to the painted graph 01949 01950 if (fPaintGraph) return fPaintGraph->DistancetoPrimitive(px,py); 01951 if (fPaintHisto) return fPaintHisto->DistancetoPrimitive(px,py); 01952 return 0; 01953 } 01954 01955 01956 //______________________________________________________________________________ 01957 void TEfficiency::Draw(Option_t* opt) 01958 { 01959 //draws the current TEfficiency object 01960 // 01961 //options: 01962 //- 1-dimensional case: same options as TGraphAsymmErrors::Draw() 01963 // but as default "AP" is used 01964 //- 2-dimensional case: same options as TH2::Draw() 01965 //- 3-dimensional case: not yet supported 01966 // 01967 // specific TEfficiency drawing options: 01968 // - E0 - plot bins where the total number of passed events is zero 01969 // (the error interval will be [0,1] ) 01970 01971 //check options 01972 TString option = opt; 01973 option.ToLower(); 01974 // use by default "AP" 01975 if (option.IsNull() ) option = "ap"; 01976 01977 if(gPad && !option.Contains("same")) 01978 gPad->Clear(); 01979 01980 AppendPad(option.Data()); 01981 } 01982 01983 //______________________________________________________________________________ 01984 void TEfficiency::ExecuteEvent(Int_t event, Int_t px, Int_t py) 01985 { 01986 // Execute action corresponding to one event. 01987 // 01988 // This member function is called when the drawn class is clicked with the locator 01989 // If Left button clicked on one of the line end points, this point 01990 // follows the cursor until button is released. 01991 // 01992 // if Middle button clicked, the line is moved parallel to itself 01993 // until the button is released. 01994 // Forward the call to the underlying graph 01995 if (fPaintGraph) fPaintGraph->ExecuteEvent(event,px,py); 01996 else if (fPaintHisto) fPaintHisto->ExecuteEvent(event,px,py); 01997 } 01998 01999 //______________________________________________________________________________ 02000 void TEfficiency::Fill(Bool_t bPassed,Double_t x,Double_t y,Double_t z) 02001 { 02002 //This function is used for filling the two histograms. 02003 // 02004 //Input: bPassed - flag whether the current event passed the selection 02005 // true: both histograms are filled 02006 // false: only the total histogram is filled 02007 // x - x value 02008 // y - y value (use default=0 for 1-D efficiencies) 02009 // z - z value (use default=0 for 2-D or 1-D efficiencies) 02010 02011 switch(GetDimension()) { 02012 case 1: 02013 fTotalHistogram->Fill(x); 02014 if(bPassed) 02015 fPassedHistogram->Fill(x); 02016 break; 02017 case 2: 02018 ((TH2*)(fTotalHistogram))->Fill(x,y); 02019 if(bPassed) 02020 ((TH2*)(fPassedHistogram))->Fill(x,y); 02021 break; 02022 case 3: 02023 ((TH3*)(fTotalHistogram))->Fill(x,y,z); 02024 if(bPassed) 02025 ((TH3*)(fPassedHistogram))->Fill(x,y,z); 02026 break; 02027 } 02028 } 02029 02030 //______________________________________________________________________________ 02031 Int_t TEfficiency::FindFixBin(Double_t x,Double_t y,Double_t z) const 02032 { 02033 //returns the global bin number containing the given values 02034 // 02035 //Note: - values which belong to dimensions higher than the current dimension 02036 // of the TEfficiency object are ignored (i.e. for 1-dimensional 02037 // efficiencies only the x-value is considered) 02038 02039 Int_t nx = fTotalHistogram->GetXaxis()->FindFixBin(x); 02040 Int_t ny = 0; 02041 Int_t nz = 0; 02042 02043 switch(GetDimension()) { 02044 case 3: nz = fTotalHistogram->GetZaxis()->FindFixBin(z); 02045 case 2: ny = fTotalHistogram->GetYaxis()->FindFixBin(y);break; 02046 } 02047 02048 return GetGlobalBin(nx,ny,nz); 02049 } 02050 02051 //______________________________________________________________________________ 02052 Int_t TEfficiency::Fit(TF1* f1,Option_t* opt) 02053 { 02054 //fits the efficiency using the TBinomialEfficiencyFitter class 02055 // 02056 //The resulting fit function is added to the list of associated functions. 02057 // 02058 //Options: - "+": previous fitted functions in the list are kept, by default 02059 // all functions in the list are deleted 02060 // - for more fitting options see TBinomialEfficiencyFitter::Fit 02061 02062 TString option = opt; 02063 option.ToLower(); 02064 02065 //replace existing functions in list with same name 02066 Bool_t bDeleteOld = true; 02067 if(option.Contains("+")) { 02068 option.ReplaceAll("+",""); 02069 bDeleteOld = false; 02070 } 02071 02072 TBinomialEfficiencyFitter Fitter(fPassedHistogram,fTotalHistogram); 02073 02074 Int_t result = Fitter.Fit(f1,option.Data()); 02075 02076 //create copy which is appended to the list 02077 TF1* pFunc = new TF1(*f1); 02078 02079 if(bDeleteOld) { 02080 TIter next(fFunctions); 02081 TObject* obj = 0; 02082 while((obj = next())) { 02083 if(obj->InheritsFrom(TF1::Class())) { 02084 fFunctions->Remove(obj); 02085 delete obj; 02086 } 02087 } 02088 } 02089 02090 fFunctions->Add(pFunc); 02091 02092 return result; 02093 } 02094 02095 //______________________________________________________________________________ 02096 TH1* TEfficiency::GetCopyPassedHisto() const 02097 { 02098 //returns a cloned version of fPassedHistogram 02099 // 02100 //Notes: - The histogram is filled with unit weights. You might want to scale 02101 // it with the global weight GetWeight(). 02102 // - The returned object is owned by the user who has to care about the 02103 // deletion of the new TH1 object. 02104 // - This histogram is by default NOT attached to the current directory 02105 // to avoid duplication of data. If you want to store it automatically 02106 // during the next TFile::Write() command, you have to attach it to 02107 // the corresponding directory. 02108 //Begin_html 02109 //<div class="code"><pre> 02110 // TFile* pFile = new TFile("passed.root","update"); 02111 // TEfficiency* pEff = (TEfficiency*)gDirectory->Get("my_eff"); 02112 // TH1* copy = pEff->GetCopyPassedHisto(); 02113 // copy->SetDirectory(gDirectory); 02114 // pFile->Write(); 02115 //</pre></div> 02116 //<div class="clear"></div> 02117 //End_Html 02118 02119 Bool_t bStatus = TH1::AddDirectoryStatus(); 02120 TH1::AddDirectory(kFALSE); 02121 TH1* tmp = (TH1*)(fPassedHistogram->Clone()); 02122 TH1::AddDirectory(bStatus); 02123 02124 return tmp; 02125 } 02126 02127 //______________________________________________________________________________ 02128 TH1* TEfficiency::GetCopyTotalHisto() const 02129 { 02130 //returns a cloned version of fTotalHistogram 02131 // 02132 //Notes: - The histogram is filled with unit weights. You might want to scale 02133 // it with the global weight GetWeight(). 02134 // - The returned object is owned by the user who has to care about the 02135 // deletion of the new TH1 object. 02136 // - This histogram is by default NOT attached to the current directory 02137 // to avoid duplication of data. If you want to store it automatically 02138 // during the next TFile::Write() command, you have to attach it to 02139 // the corresponding directory. 02140 //Begin_Html 02141 //<div class="code"><pre> 02142 // TFile* pFile = new TFile("total.root","update"); 02143 // TEfficiency* pEff = (TEfficiency*)gDirectory->Get("my_eff"); 02144 // TH1* copy = pEff->GetCopyTotalHisto(); 02145 // copy->SetDirectory(gDirectory); 02146 // pFile->Write(); 02147 //</pre></div> 02148 //<div class="clear"></div> 02149 //End_Html 02150 02151 Bool_t bStatus = TH1::AddDirectoryStatus(); 02152 TH1::AddDirectory(kFALSE); 02153 TH1* tmp = (TH1*)(fTotalHistogram->Clone()); 02154 TH1::AddDirectory(bStatus); 02155 02156 return tmp; 02157 } 02158 02159 //______________________________________________________________________________ 02160 Int_t TEfficiency::GetDimension() const 02161 { 02162 //returns the dimension of the current TEfficiency object 02163 02164 return fTotalHistogram->GetDimension(); 02165 } 02166 02167 //______________________________________________________________________________ 02168 Double_t TEfficiency::GetEfficiency(Int_t bin) const 02169 { 02170 //returns the efficiency in the given global bin 02171 // 02172 //Note: - The estimated efficiency depends on the chosen statistic option: 02173 // for frequentist ones: 02174 // Begin_Latex #hat{#varepsilon} = #frac{passed}{total} End_Latex 02175 // for bayesian ones the expectation value of the resulting posterior 02176 // distribution is returned: 02177 // Begin_Latex #hat{#varepsilon} = #frac{passed + #alpha}{total + #alpha + #beta} End_Latex 02178 // If the bit kPosteriorMode is set (or the method TEfficiency::UsePosteriorMode() has been called ) the 02179 // mode (most probable value) of the posterior is returned: 02180 // Begin_Latex #hat{#varepsilon} = #frac{passed + #alpha -1}{total + #alpha + #beta -2} End_Latex 02181 // 02182 // - If the denominator is equal to 0, an efficiency of 0 is returned. 02183 // - When Begin_Latex passed + #alpha < 1 End_Latex or Begin_Latex total - passed + #beta < 1 End_latex the above 02184 // formula for the mode is not valid. In these cases values the estimated efficiency is 0 or 1. 02185 02186 Int_t total = (Int_t)fTotalHistogram->GetBinContent(bin); 02187 Int_t passed = (Int_t)fPassedHistogram->GetBinContent(bin); 02188 02189 if(TestBit(kIsBayesian)) { 02190 02191 // parameters for the beta prior distribution 02192 Double_t alpha = TestBit(kUseBinPrior) ? GetBetaAlpha(bin) : GetBetaAlpha(); 02193 Double_t beta = TestBit(kUseBinPrior) ? GetBetaBeta(bin) : GetBetaBeta(); 02194 02195 if (!TestBit(kPosteriorMode) ) 02196 return BetaMean(double(passed) + alpha , double(total-passed) + beta); 02197 else 02198 return BetaMode(double(passed) + alpha , double(total-passed) + beta); 02199 02200 } 02201 else 02202 return (total)? ((Double_t)passed)/total : 0; 02203 } 02204 02205 //______________________________________________________________________________ 02206 Double_t TEfficiency::GetEfficiencyErrorLow(Int_t bin) const 02207 { 02208 //returns the lower error on the efficiency in the given global bin 02209 // 02210 //The result depends on the current confidence level fConfLevel and the 02211 //chosen statistic option fStatisticOption. See SetStatisticOption(Int_t) for 02212 //more details. 02213 02214 Int_t total = (Int_t)fTotalHistogram->GetBinContent(bin); 02215 Int_t passed = (Int_t)fPassedHistogram->GetBinContent(bin); 02216 02217 Double_t eff = GetEfficiency(bin); 02218 // parameters for the beta prior distribution 02219 Double_t alpha = TestBit(kUseBinPrior) ? GetBetaAlpha(bin) : GetBetaAlpha(); 02220 Double_t beta = TestBit(kUseBinPrior) ? GetBetaBeta(bin) : GetBetaBeta(); 02221 02222 if(TestBit(kIsBayesian)) 02223 return (eff - Bayesian(total,passed,fConfLevel,alpha,beta,false,TestBit(kShortestInterval))); 02224 else 02225 return (eff - fBoundary(total,passed,fConfLevel,false)); 02226 } 02227 02228 //______________________________________________________________________________ 02229 Double_t TEfficiency::GetEfficiencyErrorUp(Int_t bin) const 02230 { 02231 //returns the upper error on the efficiency in the given global bin 02232 // 02233 //The result depends on the current confidence level fConfLevel and the 02234 //chosen statistic option fStatisticOption. See SetStatisticOption(Int_t) for 02235 //more details. 02236 02237 Int_t total = (Int_t)fTotalHistogram->GetBinContent(bin); 02238 Int_t passed = (Int_t)fPassedHistogram->GetBinContent(bin); 02239 02240 Double_t eff = GetEfficiency(bin); 02241 // parameters for the beta prior distribution 02242 Double_t alpha = TestBit(kUseBinPrior) ? GetBetaAlpha(bin) : GetBetaAlpha(); 02243 Double_t beta = TestBit(kUseBinPrior) ? GetBetaBeta(bin) : GetBetaBeta(); 02244 02245 if(TestBit(kIsBayesian)) 02246 return (Bayesian(total,passed,fConfLevel,alpha,beta,true,TestBit(kShortestInterval)) - eff); 02247 else 02248 return fBoundary(total,passed,fConfLevel,true) - eff; 02249 } 02250 02251 //______________________________________________________________________________ 02252 Int_t TEfficiency::GetGlobalBin(Int_t binx,Int_t biny,Int_t binz) const 02253 { 02254 //returns the global bin number which can be used as argument for the 02255 //following functions: 02256 // 02257 // - GetEfficiency(bin), GetEfficiencyErrorLow(bin), GetEfficiencyErrorUp(bin) 02258 // - GetPassedEvents(bin), SetPassedEvents(bin), GetTotalEvents(bin), 02259 // SetTotalEvents(bin) 02260 // 02261 //see TH1::GetBin() for conventions on numbering bins 02262 02263 return fTotalHistogram->GetBin(binx,biny,binz); 02264 } 02265 02266 //______________________________________________________________________________ 02267 void TEfficiency::Merge(TCollection* pList) 02268 { 02269 //merges the TEfficiency objects in the given list to the given 02270 //TEfficiency object using the operator+=(TEfficiency&) 02271 // 02272 //The merged result is stored in the current object. The statistic options and 02273 //the confidence level are taken from the current object. 02274 // 02275 //This function should be used when all TEfficiency objects correspond to 02276 //the same process. 02277 // 02278 //The new weight is set according to: 02279 //Begin_Latex #frac{1}{w_{new}} = #sum_{i} \frac{1}{w_{i}}End_Latex 02280 02281 if(pList->IsEmpty()) 02282 return; 02283 02284 TIter next(pList); 02285 TObject* obj = 0; 02286 TEfficiency* pEff = 0; 02287 while((obj = next())) { 02288 pEff = dynamic_cast<TEfficiency*>(obj); 02289 if(pEff) 02290 *this += *pEff; 02291 } 02292 } 02293 02294 //______________________________________________________________________________ 02295 Double_t TEfficiency::Normal(Int_t total,Int_t passed,Double_t level,Bool_t bUpper) 02296 { 02297 //returns the confidence limits for the efficiency supposing that the 02298 //efficiency follows a normal distribution with the rms below 02299 // 02300 //Input: - total : number of total events 02301 // - passed: 0 <= number of passed events <= total 02302 // - level : confidence level 02303 // - bUpper: true - upper boundary is returned 02304 // false - lower boundary is returned 02305 // 02306 //calculation: 02307 //Begin_Latex(separator='=',align='rl') 02308 // #hat{#varepsilon} = #frac{passed}{total} 02309 // #sigma_{#varepsilon} = #sqrt{#frac{#hat{#varepsilon} (1 - #hat{#varepsilon})}{total}} 02310 // #varepsilon_{low} = #hat{#varepsilon} #pm #Phi^{-1}(#frac{level}{2},#sigma_{#varepsilon}) 02311 //End_Latex 02312 02313 Double_t alpha = (1.0 - level)/2; 02314 if (total == 0) return (bUpper) ? 1 : 0; 02315 Double_t average = ((Double_t)passed) / total; 02316 Double_t sigma = std::sqrt(average * (1 - average) / total); 02317 Double_t delta = ROOT::Math::normal_quantile(1 - alpha,sigma); 02318 02319 if(bUpper) 02320 return ((average + delta) > 1) ? 1.0 : (average + delta); 02321 else 02322 return ((average - delta) < 0) ? 0.0 : (average - delta); 02323 } 02324 02325 //______________________________________________________________________________ 02326 TEfficiency& TEfficiency::operator+=(const TEfficiency& rhs) 02327 { 02328 //adds the histograms of another TEfficiency object to current histograms 02329 // 02330 //The statistic options and the confidence level remain unchanged. 02331 // 02332 //fTotalHistogram += rhs.fTotalHistogram; 02333 //fPassedHistogram += rhs.fPassedHistogram; 02334 // 02335 //calculates a new weight: 02336 //current weight of this TEfficiency object = Begin_Latex w_{1} End_Latex 02337 //weight of rhs = Begin_Latex w_{2} End_Latex 02338 //Begin_Latex w_{new} = \frac{w_{1} \times w_{2}}{w_{1} + w_{2}}End_Latex 02339 02340 02341 if (fTotalHistogram == 0 && fPassedHistogram == 0) { 02342 // efficiency is empty just copy it over 02343 *this = rhs; 02344 return *this; 02345 } 02346 else if (fTotalHistogram == 0 || fPassedHistogram == 0) { 02347 Fatal("operator+=","Adding to a non consistent TEfficiency object which has not a total or a passed histogram "); 02348 return *this; 02349 } 02350 02351 if (rhs.fTotalHistogram == 0 && rhs.fTotalHistogram == 0 ) { 02352 Warning("operator+=","no operation: adding an empty object"); 02353 return *this; 02354 } 02355 else if (rhs.fTotalHistogram == 0 || rhs.fTotalHistogram == 0 ) { 02356 Fatal("operator+=","Adding a non consistent TEfficiency object which has not a total or a passed histogram "); 02357 return *this; 02358 } 02359 02360 fTotalHistogram->ResetBit(TH1::kIsAverage); 02361 fPassedHistogram->ResetBit(TH1::kIsAverage); 02362 02363 fTotalHistogram->Add(rhs.fTotalHistogram); 02364 fPassedHistogram->Add(rhs.fPassedHistogram); 02365 02366 SetWeight((fWeight * rhs.GetWeight())/(fWeight + rhs.GetWeight())); 02367 02368 return *this; 02369 } 02370 02371 //______________________________________________________________________________ 02372 TEfficiency& TEfficiency::operator=(const TEfficiency& rhs) 02373 { 02374 //assignment operator 02375 // 02376 //The histograms, statistic option, confidence level, weight and paint styles 02377 //of rhs are copied to the this TEfficiency object. 02378 // 02379 //Note: - The list of associated functions is not copied. After this 02380 // operation the list of associated functions is empty. 02381 02382 if(this != &rhs) 02383 { 02384 //statistic options 02385 SetStatisticOption(rhs.GetStatisticOption()); 02386 SetConfidenceLevel(rhs.GetConfidenceLevel()); 02387 SetBetaAlpha(rhs.GetBetaAlpha()); 02388 SetBetaBeta(rhs.GetBetaBeta()); 02389 SetWeight(rhs.GetWeight()); 02390 02391 //associated list of functions 02392 if(fFunctions) 02393 fFunctions->Delete(); 02394 02395 //copy histograms 02396 delete fTotalHistogram; 02397 delete fPassedHistogram; 02398 02399 Bool_t bStatus = TH1::AddDirectoryStatus(); 02400 TH1::AddDirectory(kFALSE); 02401 fTotalHistogram = (TH1*)(rhs.fTotalHistogram->Clone()); 02402 fPassedHistogram = (TH1*)(rhs.fPassedHistogram->Clone()); 02403 TH1::AddDirectory(bStatus); 02404 02405 //delete temporary paint objects 02406 delete fPaintHisto; 02407 delete fPaintGraph; 02408 fPaintHisto = 0; 02409 fPaintGraph = 0; 02410 02411 //copy style 02412 rhs.TAttLine::Copy(*this); 02413 rhs.TAttFill::Copy(*this); 02414 rhs.TAttMarker::Copy(*this); 02415 } 02416 02417 return *this; 02418 } 02419 02420 //______________________________________________________________________________ 02421 void TEfficiency::Paint(const Option_t* opt) 02422 { 02423 //paints this TEfficiency object 02424 // 02425 //For details on the possible option see Draw(Option_t*) 02426 // 02427 // Note for 1D classes 02428 // In 1D the TEfficiency uses a TGraphAsymmErrors for drawing 02429 // The TGraph is created only the first time Paint is used. The user can manipulate the 02430 // TGraph via the method TEfficiency::GetPaintedGraph() 02431 // The TGraph creates behing an histogram for the axis. The histogram is created also only the first time. 02432 // If the axis needs to be updated because in the meantime the class changed use this trick 02433 // which will trigger a re-calculation of the axis of the graph 02434 // TEfficiency::GetPaintedGraph()->Set(0) 02435 // 02436 02437 02438 02439 if(!gPad) 02440 return; 02441 02442 TString option = opt; 02443 option.ToLower(); 02444 02445 Bool_t plot0Bins = false; 02446 if (option.Contains("e0") ) plot0Bins = true; 02447 02448 //use TGraphAsymmErrors for painting 02449 if(GetDimension() == 1) { 02450 Int_t npoints = fTotalHistogram->GetNbinsX(); 02451 if(!fPaintGraph) { 02452 fPaintGraph = new TGraphAsymmErrors(npoints); 02453 fPaintGraph->SetName("eff_graph"); 02454 } 02455 02456 //errors for points 02457 Double_t x,y,xlow,xup,ylow,yup; 02458 //point i corresponds to bin i+1 in histogram 02459 // point j is point graph index 02460 // LM: cannot use TGraph::SetPoint because it deletes the underlying 02461 // histogram each time (see TGraph::SetPoint) 02462 // so use it only when extra points are added to the graph 02463 Int_t j = 0; 02464 double * px = fPaintGraph->GetX(); 02465 double * py = fPaintGraph->GetY(); 02466 double * exl = fPaintGraph->GetEXlow(); 02467 double * exh = fPaintGraph->GetEXhigh(); 02468 double * eyl = fPaintGraph->GetEYlow(); 02469 double * eyh = fPaintGraph->GetEYhigh(); 02470 for (Int_t i = 0; i < npoints; ++i) { 02471 if (!plot0Bins && fTotalHistogram->GetBinContent(i+1) == 0 ) continue; 02472 x = fTotalHistogram->GetBinCenter(i+1); 02473 y = GetEfficiency(i+1); 02474 xlow = fTotalHistogram->GetBinCenter(i+1) - fTotalHistogram->GetBinLowEdge(i+1); 02475 xup = fTotalHistogram->GetBinWidth(i+1) - xlow; 02476 ylow = GetEfficiencyErrorLow(i+1); 02477 yup = GetEfficiencyErrorUp(i+1); 02478 // in the case the graph already existed and extra points have been added 02479 if (j >= fPaintGraph->GetN() ) { 02480 fPaintGraph->SetPoint(j,x,y); 02481 fPaintGraph->SetPointError(j,xlow,xup,ylow,yup); 02482 } 02483 else { 02484 px[j] = x; 02485 py[j] = y; 02486 exl[j] = xlow; 02487 exh[j] = xup; 02488 eyl[j] = ylow; 02489 eyh[j] = yup; 02490 } 02491 j++; 02492 } 02493 02494 // tell the graph the effective number of points 02495 fPaintGraph->Set(j); 02496 //refresh title before painting if changed 02497 TString oldTitle = fPaintGraph->GetTitle(); 02498 TString newTitle = GetTitle(); 02499 if (oldTitle != newTitle ) 02500 fPaintGraph->SetTitle(newTitle); 02501 02502 //copying style information 02503 TAttLine::Copy(*fPaintGraph); 02504 TAttFill::Copy(*fPaintGraph); 02505 TAttMarker::Copy(*fPaintGraph); 02506 02507 // this method forces the graph to compute correctly the axis 02508 // according to the given points 02509 fPaintGraph->GetHistogram(); 02510 02511 //paint graph 02512 fPaintGraph->Paint(option.Data()); 02513 02514 //paint all associated functions 02515 if(fFunctions) { 02516 //paint box with fit parameters 02517 gStyle->SetOptFit(1); 02518 TIter next(fFunctions); 02519 TObject* obj = 0; 02520 while((obj = next())) { 02521 if(obj->InheritsFrom(TF1::Class())) { 02522 fPaintGraph->PaintStats((TF1*)obj); 02523 ((TF1*)obj)->Paint("sameC"); 02524 } 02525 } 02526 } 02527 02528 return; 02529 } 02530 02531 //use TH2 for painting 02532 if(GetDimension() == 2) { 02533 Int_t nbinsx = fTotalHistogram->GetNbinsX(); 02534 Int_t nbinsy = fTotalHistogram->GetNbinsY(); 02535 TAxis * xaxis = fTotalHistogram->GetXaxis(); 02536 TAxis * yaxis = fTotalHistogram->GetYaxis(); 02537 if(!fPaintHisto) { 02538 if (xaxis->IsVariableBinSize() && yaxis->IsVariableBinSize() ) 02539 fPaintHisto = new TH2F("eff_histo",GetTitle(),nbinsx,xaxis->GetXbins()->GetArray(), 02540 nbinsy,yaxis->GetXbins()->GetArray()); 02541 else if (xaxis->IsVariableBinSize() && ! yaxis->IsVariableBinSize() ) 02542 fPaintHisto = new TH2F("eff_histo",GetTitle(),nbinsx,xaxis->GetXbins()->GetArray(), 02543 nbinsy,yaxis->GetXmin(), yaxis->GetXmax()); 02544 else if (!xaxis->IsVariableBinSize() && yaxis->IsVariableBinSize() ) 02545 fPaintHisto = new TH2F("eff_histo",GetTitle(),nbinsx,xaxis->GetXmin(), xaxis->GetXmax(), 02546 nbinsy,yaxis->GetXbins()->GetArray()); 02547 else 02548 fPaintHisto = new TH2F("eff_histo",GetTitle(),nbinsx,xaxis->GetXmin(), xaxis->GetXmax(), 02549 nbinsy,yaxis->GetXmin(), yaxis->GetXmax()); 02550 02551 02552 02553 fPaintHisto->SetDirectory(0); 02554 } 02555 //refresh title before each painting 02556 fPaintHisto->SetTitle(GetTitle()); 02557 02558 Int_t bin; 02559 for(Int_t i = 0; i < nbinsx + 2; ++i) { 02560 for(Int_t j = 0; j < nbinsy + 2; ++j) { 02561 bin = GetGlobalBin(i,j); 02562 fPaintHisto->SetBinContent(bin,GetEfficiency(bin)); 02563 } 02564 } 02565 02566 //copying style information 02567 TAttLine::Copy(*fPaintHisto); 02568 TAttFill::Copy(*fPaintHisto); 02569 TAttMarker::Copy(*fPaintHisto); 02570 fPaintHisto->SetStats(0); 02571 02572 //paint histogram 02573 fPaintHisto->Paint(option.Data()); 02574 return; 02575 } 02576 } 02577 02578 //______________________________________________________________________________ 02579 void TEfficiency::SavePrimitive(ostream& out,Option_t* opt) 02580 { 02581 //have histograms fixed bins along each axis? 02582 Bool_t equi_bins = true; 02583 02584 //indentation 02585 TString indent = " "; 02586 //names for arrays containing the bin edges 02587 //static counter needed if more objects are saved 02588 static Int_t naxis = 0; 02589 TString sxaxis="xAxis",syaxis="yAxis",szaxis="zAxis"; 02590 02591 //note the missing break statements! 02592 switch(GetDimension()) { 02593 case 3: 02594 equi_bins = equi_bins && !fTotalHistogram->GetZaxis()->GetXbins()->fArray 02595 && !fTotalHistogram->GetZaxis()->GetXbins()->fN; 02596 case 2: 02597 equi_bins = equi_bins && !fTotalHistogram->GetYaxis()->GetXbins()->fArray 02598 && !fTotalHistogram->GetYaxis()->GetXbins()->fN; 02599 case 1: 02600 equi_bins = equi_bins && !fTotalHistogram->GetXaxis()->GetXbins()->fArray 02601 && !fTotalHistogram->GetXaxis()->GetXbins()->fN; 02602 } 02603 02604 //create arrays containing the variable binning 02605 if(!equi_bins) { 02606 Int_t i; 02607 ++naxis; 02608 sxaxis += naxis; 02609 syaxis += naxis; 02610 szaxis += naxis; 02611 //x axis 02612 out << indent << "Double_t " << sxaxis << "[" 02613 << fTotalHistogram->GetXaxis()->GetXbins()->fN << "] = {"; 02614 for (i = 0; i < fTotalHistogram->GetXaxis()->GetXbins()->fN; ++i) { 02615 if (i != 0) out << ", "; 02616 out << fTotalHistogram->GetXaxis()->GetXbins()->fArray[i]; 02617 } 02618 out << "}; " << std::endl; 02619 //y axis 02620 if(GetDimension() > 1) { 02621 out << indent << "Double_t " << syaxis << "[" 02622 << fTotalHistogram->GetYaxis()->GetXbins()->fN << "] = {"; 02623 for (i = 0; i < fTotalHistogram->GetYaxis()->GetXbins()->fN; ++i) { 02624 if (i != 0) out << ", "; 02625 out << fTotalHistogram->GetYaxis()->GetXbins()->fArray[i]; 02626 } 02627 out << "}; " << std::endl; 02628 } 02629 //z axis 02630 if(GetDimension() > 2) { 02631 out << indent << "Double_t " << szaxis << "[" 02632 << fTotalHistogram->GetZaxis()->GetXbins()->fN << "] = {"; 02633 for (i = 0; i < fTotalHistogram->GetZaxis()->GetXbins()->fN; ++i) { 02634 if (i != 0) out << ", "; 02635 out << fTotalHistogram->GetZaxis()->GetXbins()->fArray[i]; 02636 } 02637 out << "}; " << std::endl; 02638 } 02639 }//creating variable binning 02640 02641 //TEfficiency pointer has efficiency name + counter 02642 static Int_t eff_count = 0; 02643 ++eff_count; 02644 TString eff_name = GetName(); 02645 eff_name += eff_count; 02646 02647 const char* name = eff_name.Data(); 02648 02649 //construct TEfficiency object 02650 const char quote = '"'; 02651 out << indent << std::endl; 02652 out << indent << ClassName() << " * " << name << " = new " << ClassName() 02653 << "(" << quote << GetName() << quote << "," << quote 02654 << GetTitle() << quote <<","; 02655 //fixed bin size -> use n,min,max constructor 02656 if(equi_bins) { 02657 out << fTotalHistogram->GetXaxis()->GetNbins() << "," 02658 << fTotalHistogram->GetXaxis()->GetXmin() << "," 02659 << fTotalHistogram->GetXaxis()->GetXmax(); 02660 if(GetDimension() > 1) { 02661 out << "," << fTotalHistogram->GetYaxis()->GetNbins() << "," 02662 << fTotalHistogram->GetYaxis()->GetXmin() << "," 02663 << fTotalHistogram->GetYaxis()->GetXmax(); 02664 } 02665 if(GetDimension() > 2) { 02666 out << "," << fTotalHistogram->GetZaxis()->GetNbins() << "," 02667 << fTotalHistogram->GetZaxis()->GetXmin() << "," 02668 << fTotalHistogram->GetZaxis()->GetXmax(); 02669 } 02670 } 02671 //variable bin size -> use n,*bins constructor 02672 else { 02673 out << fTotalHistogram->GetXaxis()->GetNbins() << "," << sxaxis; 02674 if(GetDimension() > 1) 02675 out << "," << fTotalHistogram->GetYaxis()->GetNbins() << "," 02676 << syaxis; 02677 if(GetDimension() > 2) 02678 out << "," << fTotalHistogram->GetZaxis()->GetNbins() << "," 02679 << szaxis; 02680 } 02681 out << ");" << std::endl; 02682 out << indent << std::endl; 02683 02684 //set statistic options 02685 out << indent << name << "->SetConfidenceLevel(" << fConfLevel << ");" 02686 << std::endl; 02687 out << indent << name << "->SetBetaAlpha(" << fBeta_alpha << ");" 02688 << std::endl; 02689 out << indent << name << "->SetBetaBeta(" << fBeta_beta << ");" << std::endl; 02690 out << indent << name << "->SetWeight(" << fWeight << ");" << std::endl; 02691 out << indent << name << "->SetStatisticOption(" << fStatisticOption << ");" 02692 << std::endl; 02693 02694 //set bin contents 02695 Int_t nbins = fTotalHistogram->GetNbinsX() + 2; 02696 if(GetDimension() > 1) 02697 nbins *= fTotalHistogram->GetNbinsY() + 2; 02698 if(GetDimension() > 2) 02699 nbins *= fTotalHistogram->GetNbinsZ() + 2; 02700 02701 //important: set first total number than passed number 02702 for(Int_t i = 0; i < nbins; ++i) { 02703 out << indent << name <<"->SetTotalEvents(" << i << "," << 02704 fTotalHistogram->GetBinContent(i) << ");" << std::endl; 02705 out << indent << name <<"->SetPassedEvents(" << i << "," << 02706 fPassedHistogram->GetBinContent(i) << ");" << std::endl; 02707 } 02708 02709 //save list of functions 02710 TIter next(fFunctions); 02711 TObject* obj = 0; 02712 while((obj = next())) { 02713 obj->SavePrimitive(out,"nodraw"); 02714 if(obj->InheritsFrom(TF1::Class())) { 02715 out << indent << name << "->GetListOfFunctions()->Add(" 02716 << obj->GetName() << ");" << std::endl; 02717 } 02718 } 02719 02720 //set style 02721 SaveFillAttributes(out,name); 02722 SaveLineAttributes(out,name); 02723 SaveMarkerAttributes(out,name); 02724 02725 //draw TEfficiency object 02726 TString option = opt; 02727 option.ToLower(); 02728 if (!option.Contains("nodraw")) 02729 out<< indent << name<< "->Draw(" << quote << opt << quote << ");" 02730 << std::endl; 02731 } 02732 02733 //______________________________________________________________________________ 02734 void TEfficiency::SetBetaAlpha(Double_t alpha) 02735 { 02736 //sets the shape parameter Begin_Latex \alpha End_Latex 02737 // 02738 //The prior probability of the efficiency is given by the beta distribution: 02739 //Begin_Latex 02740 // f(\varepsilon;\alpha;\beta) = \frac{1}{B(\alpha,\beta)} \varepsilon^{\alpha-1} (1 - \varepsilon)^{\beta-1} 02741 //End_Latex 02742 // 02743 //Note: - both shape parameters have to be positive (i.e. > 0) 02744 02745 if(alpha > 0) 02746 fBeta_alpha = alpha; 02747 else 02748 Warning("SetBetaAlpha(Double_t)","invalid shape parameter %.2lf",alpha); 02749 } 02750 02751 //______________________________________________________________________________ 02752 void TEfficiency::SetBetaBeta(Double_t beta) 02753 { 02754 //sets the shape parameter Begin_Latex \beta End_Latex 02755 // 02756 //The prior probability of the efficiency is given by the beta distribution: 02757 //Begin_Latex 02758 // f(\varepsilon;\alpha,\beta) = \frac{1}{B(\alpha,\beta)} \varepsilon^{\alpha-1} (1 - \varepsilon)^{\beta-1} 02759 //End_Latex 02760 // 02761 //Note: - both shape parameters have to be positive (i.e. > 0) 02762 02763 if(beta > 0) 02764 fBeta_beta = beta; 02765 else 02766 Warning("SetBetaBeta(Double_t)","invalid shape parameter %.2lf",beta); 02767 } 02768 02769 //______________________________________________________________________________ 02770 void TEfficiency::SetBetaBinParameters(Int_t bin, Double_t alpha, Double_t beta) 02771 { 02772 //sets different shape parameter Begin_Latex \alpha and \beta End_Latex 02773 // for the prior distribution for each bin. By default the global parameter are used if they are not set 02774 // for the specific bin 02775 //The prior probability of the efficiency is given by the beta distribution: 02776 //Begin_Latex 02777 // f(\varepsilon;\alpha;\beta) = \frac{1}{B(\alpha,\beta)} \varepsilon^{\alpha-1} (1 - \varepsilon)^{\beta-1} 02778 //End_Latex 02779 // 02780 //Note: - both shape parameters have to be positive (i.e. > 0) 02781 02782 if (!fPassedHistogram || !fTotalHistogram) return; 02783 TH1 * h1 = fTotalHistogram; 02784 // doing this I get h1->fN which is available only for a TH1D 02785 UInt_t n = h1->GetBin(h1->GetNbinsX()+1, h1->GetNbinsY()+1, h1->GetNbinsZ()+1 ) + 1; 02786 02787 // in case vector is not created do with defult alpha, beta params 02788 if (fBeta_bin_params.size() != n ) 02789 fBeta_bin_params = std::vector<std::pair<Double_t, Double_t> >(n, std::make_pair(fBeta_alpha, fBeta_beta) ); 02790 02791 // vector contains also values for under/overflows 02792 fBeta_bin_params[bin] = std::make_pair(alpha,beta); 02793 SetBit(kUseBinPrior,true); 02794 02795 } 02796 02797 //______________________________________________________________________________ 02798 void TEfficiency::SetConfidenceLevel(Double_t level) 02799 { 02800 //sets the confidence level (0 < level < 1) 02801 // The default value is 1-sigma :~ 0.683 02802 02803 if((level > 0) && (level < 1)) 02804 fConfLevel = level; 02805 else 02806 Warning("SetConfidenceLevel(Double_t)","invalid confidence level %.2lf",level); 02807 } 02808 02809 //______________________________________________________________________________ 02810 void TEfficiency::SetDirectory(TDirectory* dir) 02811 { 02812 //sets the directory holding this TEfficiency object 02813 // 02814 //A reference to this TEfficiency object is removed from the current 02815 //directory (if it exists) and a new reference to this TEfficiency object is 02816 //added to the given directory. 02817 // 02818 //Notes: - If the given directory is 0, the TEfficiency object does not 02819 // belong to any directory and will not be written to file during the 02820 // next TFile::Write() command. 02821 02822 if(fDirectory == dir) 02823 return; 02824 if(fDirectory) 02825 fDirectory->Remove(this); 02826 fDirectory = dir; 02827 if(fDirectory) 02828 fDirectory->Append(this); 02829 } 02830 02831 //______________________________________________________________________________ 02832 void TEfficiency::SetName(const char* name) 02833 { 02834 //sets the name 02835 // 02836 //Note: The names of the internal histograms are set to "name + _total" or 02837 // "name + _passed" respectively. 02838 02839 TNamed::SetName(name); 02840 02841 //setting the names (appending the correct ending) 02842 TString name_total = name + TString("_total"); 02843 TString name_passed = name + TString("_passed"); 02844 fTotalHistogram->SetName(name_total); 02845 fPassedHistogram->SetName(name_passed); 02846 } 02847 02848 //______________________________________________________________________________ 02849 Bool_t TEfficiency::SetPassedEvents(Int_t bin,Int_t events) 02850 { 02851 //sets the number of passed events in the given global bin 02852 // 02853 //returns "true" if the number of passed events has been updated 02854 //otherwise "false" ist returned 02855 // 02856 //Note: - requires: 0 <= events <= fTotalHistogram->GetBinContent(bin) 02857 02858 if(events <= fTotalHistogram->GetBinContent(bin)) { 02859 fPassedHistogram->SetBinContent(bin,events); 02860 return true; 02861 } 02862 else { 02863 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); 02864 return false; 02865 } 02866 } 02867 02868 //______________________________________________________________________________ 02869 Bool_t TEfficiency::SetPassedHistogram(const TH1& rPassed,Option_t* opt) 02870 { 02871 //sets the histogram containing the passed events 02872 // 02873 //The given histogram is cloned and stored internally as histogram containing 02874 //the passed events. The given histogram has to be consistent with the current 02875 //fTotalHistogram (see CheckConsistency(const TH1&,const TH1&)). 02876 //The method returns whether the fPassedHistogram has been replaced (true) or 02877 //not (false). 02878 // 02879 //Note: The list of associated functions fFunctions is cleared. 02880 // 02881 //Option: - "f": force the replacement without checking the consistency 02882 // This can lead to inconsistent histograms and useless results 02883 // or unexpected behaviour. But sometimes it might be the only 02884 // way to change the histograms. If you use this option, you 02885 // should ensure that the fTotalHistogram is replaced by a 02886 // consistent one (with respect to rPassed) as well. 02887 02888 TString option = opt; 02889 option.ToLower(); 02890 02891 Bool_t bReplace = option.Contains("f"); 02892 02893 if(!bReplace) 02894 bReplace = CheckConsistency(rPassed,*fTotalHistogram); 02895 02896 if(bReplace) { 02897 delete fPassedHistogram; 02898 Bool_t bStatus = TH1::AddDirectoryStatus(); 02899 TH1::AddDirectory(kFALSE); 02900 fPassedHistogram = (TH1*)(rPassed.Clone()); 02901 fPassedHistogram->SetNormFactor(0); 02902 TH1::AddDirectory(bStatus); 02903 02904 if(fFunctions) 02905 fFunctions->Delete(); 02906 02907 return true; 02908 } 02909 else 02910 return false; 02911 } 02912 02913 //______________________________________________________________________________ 02914 void TEfficiency::SetStatisticOption(EStatOption option) 02915 { 02916 //sets the statistic option which affects the calculation of the confidence interval 02917 // 02918 //Options: 02919 //- kFCP (=0)(default): using the Clopper-Pearson interval (recommended by PDG) 02920 // sets kIsBayesian = false 02921 // see also ClopperPearson 02922 //- kFNormal (=1) : using the normal approximation 02923 // sets kIsBayesian = false 02924 // see also Normal 02925 //- kFWilson (=2) : using the Wilson interval 02926 // sets kIsBayesian = false 02927 // see also Wilson 02928 //- kFAC (=3) : using the Agresti-Coull interval 02929 // sets kIsBayesian = false 02930 // see also AgrestiCoull 02931 //- kFFC (=4) : using the Feldman-Cousins frequentist method 02932 // sets kIsBayesian = false 02933 // see also FeldmanCousins 02934 //- kBJeffrey (=5) : using the Jeffrey interval 02935 // sets kIsBayesian = true, fBeta_alpha = 0.5 and fBeta_beta = 0.5 02936 // see also Bayesian 02937 //- kBUniform (=6) : using a uniform prior 02938 // sets kIsBayesian = true, fBeta_alpha = 1 and fBeta_beta = 1 02939 // see also Bayesian 02940 //- kBBayesian (=7) : using a custom prior defined by fBeta_alpha and fBeta_beta 02941 // sets kIsBayesian = true 02942 // see also Bayesian 02943 02944 fStatisticOption = option; 02945 02946 switch(option) 02947 { 02948 case kFCP: 02949 fBoundary = &ClopperPearson; 02950 SetBit(kIsBayesian,false); 02951 break; 02952 case kFNormal: 02953 fBoundary = &Normal; 02954 SetBit(kIsBayesian,false); 02955 break; 02956 case kFWilson: 02957 fBoundary = &Wilson; 02958 SetBit(kIsBayesian,false); 02959 break; 02960 case kFAC: 02961 fBoundary = &AgrestiCoull; 02962 SetBit(kIsBayesian,false); 02963 break; 02964 case kFFC: 02965 fBoundary = &FeldmanCousins; 02966 SetBit(kIsBayesian,false); 02967 break; 02968 case kBJeffrey: 02969 fBeta_alpha = 0.5; 02970 fBeta_beta = 0.5; 02971 SetBit(kIsBayesian,true); 02972 SetBit(kUseBinPrior,false); 02973 break; 02974 case kBUniform: 02975 fBeta_alpha = 1; 02976 fBeta_beta = 1; 02977 SetBit(kIsBayesian,true); 02978 SetBit(kUseBinPrior,false); 02979 break; 02980 case kBBayesian: 02981 SetBit(kIsBayesian,true); 02982 break; 02983 default: 02984 fStatisticOption = kFCP; 02985 fBoundary = &ClopperPearson; 02986 SetBit(kIsBayesian,false); 02987 } 02988 } 02989 02990 //______________________________________________________________________________ 02991 void TEfficiency::SetTitle(const char* title) 02992 { 02993 //sets the title 02994 // 02995 //Notes: - The titles of the internal histograms are set to "title + (total)" 02996 // or "title + (passed)" respectively. 02997 // - It is possible to label the axis of the histograms as usual (see 02998 // TH1::SetTitle). 02999 // 03000 //Example: Setting the title to "My Efficiency" and label the axis 03001 //Begin_Html 03002 //<div class="code"><pre> 03003 //pEff->SetTitle("My Efficiency;x label;eff"); 03004 //</pre></div> 03005 //<div class="clear"></div> 03006 //End_Html 03007 03008 TNamed::SetTitle(title); 03009 03010 //setting the titles (looking for the first semicolon and insert the tokens there) 03011 TString title_passed = title; 03012 TString title_total = title; 03013 Ssiz_t pos = title_passed.First(";"); 03014 if (pos != kNPOS) { 03015 title_passed.Insert(pos," (passed)"); 03016 title_total.Insert(pos," (total)"); 03017 } 03018 else { 03019 title_passed.Append(" (passed)"); 03020 title_total.Append(" (total)"); 03021 } 03022 fPassedHistogram->SetTitle(title_passed); 03023 fTotalHistogram->SetTitle(title_total); 03024 } 03025 03026 //______________________________________________________________________________ 03027 Bool_t TEfficiency::SetTotalEvents(Int_t bin,Int_t events) 03028 { 03029 //sets the number of total events in the given global bin 03030 // 03031 //returns "true" if the number of total events has been updated 03032 //otherwise "false" ist returned 03033 // 03034 //Note: - requires: fPassedHistogram->GetBinContent(bin) <= events 03035 03036 if(events >= fPassedHistogram->GetBinContent(bin)) { 03037 fTotalHistogram->SetBinContent(bin,events); 03038 return true; 03039 } 03040 else { 03041 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); 03042 return false; 03043 } 03044 } 03045 03046 //______________________________________________________________________________ 03047 Bool_t TEfficiency::SetTotalHistogram(const TH1& rTotal,Option_t* opt) 03048 { 03049 //sets the histogram containing all events 03050 // 03051 //The given histogram is cloned and stored internally as histogram containing 03052 //all events. The given histogram has to be consistent with the current 03053 //fPassedHistogram (see CheckConsistency(const TH1&,const TH1&)). 03054 //The method returns whether the fTotalHistogram has been replaced (true) or 03055 //not (false). 03056 // 03057 //Note: The list of associated functions fFunctions is cleared. 03058 // 03059 //Option: - "f": force the replacement without checking the consistency 03060 // This can lead to inconsistent histograms and useless results 03061 // or unexpected behaviour. But sometimes it might be the only 03062 // way to change the histograms. If you use this option, you 03063 // should ensure that the fPassedHistogram is replaced by a 03064 // consistent one (with respect to rTotal) as well. 03065 03066 TString option = opt; 03067 option.ToLower(); 03068 03069 Bool_t bReplace = option.Contains("f"); 03070 03071 if(!bReplace) 03072 bReplace = CheckConsistency(*fPassedHistogram,rTotal); 03073 03074 if(bReplace) { 03075 delete fTotalHistogram; 03076 Bool_t bStatus = TH1::AddDirectoryStatus(); 03077 TH1::AddDirectory(kFALSE); 03078 fTotalHistogram = (TH1*)(rTotal.Clone()); 03079 fTotalHistogram->SetNormFactor(0); 03080 TH1::AddDirectory(bStatus); 03081 03082 if(fFunctions) 03083 fFunctions->Delete(); 03084 03085 return true; 03086 } 03087 else 03088 return false; 03089 } 03090 03091 //______________________________________________________________________________ 03092 void TEfficiency::SetWeight(Double_t weight) 03093 { 03094 //sets the global weight for this TEfficiency object 03095 // 03096 //Note: - weight has to be positive ( > 0) 03097 03098 if(weight > 0) 03099 fWeight = weight; 03100 else 03101 Warning("SetWeight","invalid weight %.2lf",weight); 03102 } 03103 03104 //______________________________________________________________________________ 03105 Double_t TEfficiency::Wilson(Int_t total,Int_t passed,Double_t level,Bool_t bUpper) 03106 { 03107 //calculates the boundaries for the frequentist Wilson interval 03108 // 03109 //Input: - total : number of total events 03110 // - passed: 0 <= number of passed events <= total 03111 // - level : confidence level 03112 // - bUpper: true - upper boundary is returned 03113 // false - lower boundary is returned 03114 // 03115 //calculation: 03116 //Begin_Latex(separator='=',align='rl') 03117 // #alpha = 1 - #frac{level}{2} 03118 // #kappa = #Phi^{-1}(1 - #alpha,1) ... normal quantile function 03119 // mode = #frac{passed + #frac{#kappa^{2}}{2}}{total + #kappa^{2}} 03120 // #Delta = #frac{#kappa}{total + #kappa^{2}} * #sqrt{passed (1 - #frac{passed}{total}) + #frac{#kappa^{2}}{4}} 03121 // return = max(0,mode - #Delta) or min(1,mode + #Delta) 03122 //End_Latex 03123 03124 Double_t alpha = (1.0 - level)/2; 03125 if (total == 0) return (bUpper) ? 1 : 0; 03126 Double_t average = ((Double_t)passed) / total; 03127 Double_t kappa = ROOT::Math::normal_quantile(1 - alpha,1); 03128 03129 Double_t mode = (passed + 0.5 * kappa * kappa) / (total + kappa * kappa); 03130 Double_t delta = kappa / (total + kappa*kappa) * std::sqrt(total * average 03131 * (1 - average) + kappa * kappa / 4); 03132 if(bUpper) 03133 return ((mode + delta) > 1) ? 1.0 : (mode + delta); 03134 else 03135 return ((mode - delta) < 0) ? 0.0 : (mode - delta); 03136 } 03137 03138 //______________________________________________________________________________ 03139 const TEfficiency operator+(const TEfficiency& lhs,const TEfficiency& rhs) 03140 { 03141 // addition operator 03142 // 03143 // adds the corresponding histograms: 03144 // lhs.GetTotalHistogram() + rhs.GetTotalHistogram() 03145 // lhs.GetPassedHistogram() + rhs.GetPassedHistogram() 03146 // 03147 // the statistic option and the confidence level are taken from lhs 03148 03149 TEfficiency tmp(lhs); 03150 tmp += rhs; 03151 return tmp; 03152 } 03153 03154 #endif