ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rs801_HypoTestInverterOriginal.C
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'Hypothesis Test Inversion' RooStats tutorial macro #801
4 // author: Gregory Schott
5 // date Sep 2009
6 //
7 // This tutorial shows an example of using the HypoTestInverterOriginal class
8 //
9 /////////////////////////////////////////////////////////////////////////
10 
11 #include "RooRealVar.h"
12 #include "RooConstVar.h"
13 #include "RooProdPdf.h"
14 #include "RooWorkspace.h"
15 #include "RooDataSet.h"
16 #include "RooPolynomial.h"
17 #include "RooAddPdf.h"
18 #include "RooExtendPdf.h"
19 
24 
25 #include "TGraphErrors.h"
26 
27 using namespace RooFit;
28 using namespace RooStats;
29 
30 
32 {
33  // prepare the model
34  RooRealVar lumi("lumi","luminosity",1);
35  RooRealVar r("r","cross-section ratio",3.74,0,50);
36  RooFormulaVar ns("ns","1*r*lumi",RooArgList(lumi,r));
37  RooRealVar nb("nb","background yield",1);
38  RooRealVar x("x","dummy observable",0,1);
40  RooPolynomial flatPdf("flatPdf","flat PDF",x,p0);
41  RooAddPdf totPdf("totPdf","S+B model",RooArgList(flatPdf,flatPdf),RooArgList(ns,nb));
42  RooExtendPdf bkgPdf("bkgPdf","B-only model",flatPdf,nb);
43  RooDataSet* data = totPdf.generate(x,1);
44 
45  // prepare the calculator
46  HybridCalculatorOriginal myhc(*data, totPdf, bkgPdf,0,0);
47  myhc.SetTestStatistic(2);
48  myhc.SetNumberOfToys(1000);
49  myhc.UseNuisance(false);
50 
51  // run the hypothesis-test invertion
52  HypoTestInverterOriginal myInverter(myhc,r);
53  myInverter.SetTestSize(0.10);
54  myInverter.UseCLs(true);
55  // myInverter.RunFixedScan(5,1,6);
56  // scan for a 95% UL
57  myInverter.RunAutoScan(3.,5,myInverter.Size()/2,0.005);
58  // run an alternative autoscan algorithm
59  // myInverter.RunAutoScan(1,6,myInverter.Size()/2,0.005,1);
60  //myInverter.RunOnePoint(3.9);
61 
62 
63  HypoTestInverterResult* results = myInverter.GetInterval();
64 
65  HypoTestInverterPlot myInverterPlot("myInverterPlot","",results);
66  TGraphErrors* gr1 = myInverterPlot.MakePlot();
67  gr1->Draw("ALP");
68 
69  double ulError = results->UpperLimitEstimatedError();
70 
71  double upperLimit = results->UpperLimit();
72  std::cout << "The computed upper limit is: " << upperLimit << std::endl;
73  std::cout << "an estimated error on this upper limit is: " << ulError << std::endl;
74  // expected result: 4.10
75 }
76 int main() {
78 }
void rs801_HypoTestInverterOriginal()
virtual HypoTestInverterResult * GetInterval() const
Main interface to get a ConfInterval, pure virtual.
bool RunAutoScan(double xMin, double xMax, double target, double epsilon=0.005, unsigned int numAlgorithm=0)
HybridCalculatorOriginal class.
Double_t x[n]
Definition: legend1.C:17
RooConstVar represent a constant real-valued object.
Definition: RooConstVar.h:25
virtual void SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( e.g. 0.05 for a 95% Confidence Interval) ...
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual Double_t Size() const
Get the size of the test (eg. rate of Type I error)
ROOT::R::TRInterface & r
Definition: Object.C:4
TGraphErrors * MakePlot(Option_t *opt="")
return a TGraphErrors with the obtained observed p-values resultinf from the scan By default (Option ...
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
This class is now depratcated and to be replaced by the HypoTestInverter.
HypoTestInverterResult class: holds the array of hypothesis test results and compute a confidence int...
Class to plot an HypoTestInverterResult, result of the HypoTestInverter calculator.
Double_t UpperLimitEstimatedError()
Estimate of lower limit error function evaluates only a rought error on the lower limit...
RooConstVar & RooConst(Double_t val)
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
Generate a new dataset containing the specified variables with events sampled from our distribution...
Definition: RooAbsPdf.cxx:1702
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
void SetTestStatistic(int index)
set the desired test statistics: index=1 : 2 * log( L_sb / L_b ) (DEFAULT) index=2 : number of genera...