ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HybridOriginalDemo.C
Go to the documentation of this file.
1 // Example on how to use the HybridCalculatorOriginal class
2 //
3 // Author: Gregory Schott
4 #include "RooRandom.h"
5 #include "RooRealVar.h"
6 #include "RooGaussian.h"
7 #include "RooPolynomial.h"
8 #include "RooArgSet.h"
9 #include "RooAddPdf.h"
10 #include "RooDataSet.h"
11 #include "RooExtendPdf.h"
12 #include "RooConstVar.h"
13 
14 #ifndef __CINT__ // problem including this file with CINT
15 #include "RooGlobalFunc.h"
16 #endif
17 
19 #include "RooStats/HybridResult.h"
20 #include "RooStats/HybridPlot.h"
21 
22 void HybridOriginalDemo(int ntoys = 1000)
23 {
24  //***********************************************************************//
25  // This macro show an example on how to use RooStats/HybridCalculatorOriginal //
26  //***********************************************************************//
27  //
28  // With this example, you should get: CL_sb = 0.130 and CL_b = 0.946
29  // (if data had -2lnQ = -3.0742). You can compare to the expected plot:
30  // http://www-ekp.physik.uni-karlsruhe.de/~schott/roostats/hybridplot_example.png
31 
32  using namespace RooFit;
33  using namespace RooStats;
34 
35  /// set RooFit random seed
37 
38  /// build the models for background and signal+background
39  RooRealVar x("x","",-3,3);
40  RooArgList observables(x); // variables to be generated
41 
42  // gaussian signal
43 // RooRealVar sig_mean("sig_mean","",0);
44 // RooRealVar sig_sigma("sig_sigma","",0.8);
45 // RooGaussian sig_pdf("sig_pdf","",x,sig_mean,sig_sigma);
46  RooGaussian sig_pdf("sig_pdf","",x, RooConst(0.0),RooConst(0.8));
47  RooRealVar sig_yield("sig_yield","",20,0,300);
48 
49  // flat background (extended PDF)
50 // RooRealVar bkg_slope("bkg_slope","",0);
51 // RooPolynomial bkg_pdf("bkg_pdf","",x,bkg_slope);
52  RooPolynomial bkg_pdf("bkg_pdf","", x, RooConst(0));
53  RooRealVar bkg_yield("bkg_yield","",40,0,300);
54  RooExtendPdf bkg_ext_pdf("bkg_ext_pdf","",bkg_pdf,bkg_yield);
55 
56 // bkg_yield.setConstant(kTRUE);
57  sig_yield.setConstant(kTRUE);
58 
59  // total sig+bkg (extended PDF)
60  RooAddPdf tot_pdf("tot_pdf","",RooArgList(sig_pdf,bkg_pdf),RooArgList(sig_yield,bkg_yield));
61 
62  /// build the prior PDF on the parameters to be integrated
63  // gaussian contraint on the background yield ( N_B = 40 +/- 10 ie. 25% )
64  RooGaussian bkg_yield_prior("bkg_yield_prior","",bkg_yield,RooConst(bkg_yield.getVal()),RooConst(10.));
65 
66  RooArgSet nuisance_parameters(bkg_yield); // variables to be integrated
67 
68  /// generate a data sample
69  RooDataSet* data = tot_pdf.generate(observables,RooFit::Extended());
70 
71  //***********************************************************************//
72 
73  /// run HybridCalculator on those inputs
74 
75  // use interface from HypoTest calculator by default
76 
77  HybridCalculatorOriginal myHybridCalc(*data, tot_pdf , bkg_ext_pdf ,
78  &nuisance_parameters, &bkg_yield_prior);
79 
80  // here I use the default test statistics: 2*lnQ (optional)
81  myHybridCalc.SetTestStatistic(1);
82  //myHybridCalc.SetTestStatistic(3); // profile likelihood ratio
83 
84  myHybridCalc.SetNumberOfToys(ntoys);
85  myHybridCalc.UseNuisance(true);
86 
87  // for speed up generation (do binned data)
88  myHybridCalc.SetGenerateBinned(false);
89 
90  // calculate by running ntoys for the S+B and B hypothesis and retrieve the result
91  HybridResult* myHybridResult = myHybridCalc.GetHypoTest();
92 
93  if (! myHybridResult) {
94  std::cerr << "\nError returned from Hypothesis test" << std::endl;
95  return;
96  }
97 
98  /// run 1000 toys without gaussian prior on the background yield
99  //HybridResult* myHybridResult = myHybridCalc.Calculate(*data,1000,false);
100 
101  /// save the toy-MC results to file, this way splitting into sub-batch jobs is possible
102  //TFile fileOut("some_hybridresult.root","RECREATE");
103  //fileOut.cd();
104  //myHybridResult.Write();
105  //fileOut.Close();
106 
107  /// read the results from a file
108  //TFile fileIn("some_hybridresult.root");
109  //HybridResult* myOtherHybridResult = (HybridResult*) fileIn.Get("myHybridCalc");
110 
111  /// example on how to merge with toy-MC results obtained in another job
112  //HybridResult* mergedHybridResult = new HybridResult("mergedHybridResult","this object holds merged results");
113  //mergedHybridResult->Add(myHybridResult);
114  //mergedHybridResult->Add(myOtherHybridResult);
115  /// or
116  //myHybridResult->Add(myOtherHybridResult);
117 
118  /// nice plot of the results
119  HybridPlot* myHybridPlot = myHybridResult->GetPlot("myHybridPlot","Plot of results with HybridCalculatorOriginal",100);
120  myHybridPlot->Draw();
121 
122  /// recover and display the results
123  double clsb_data = myHybridResult->CLsplusb();
124  double clb_data = myHybridResult->CLb();
125  double cls_data = myHybridResult->CLs();
126  double data_significance = myHybridResult->Significance();
127  double min2lnQ_data = myHybridResult->GetTestStat_data();
128 
129  /// compute the mean expected significance from toys
130  double mean_sb_toys_test_stat = myHybridPlot->GetSBmean();
131  myHybridResult->SetDataTestStatistics(mean_sb_toys_test_stat);
132  double toys_significance = myHybridResult->Significance();
133 
134  std::cout << "Completed HybridCalculatorOriginal example:\n";
135  std::cout << " - -2lnQ = " << min2lnQ_data << endl;
136  std::cout << " - CL_sb = " << clsb_data << std::endl;
137  std::cout << " - CL_b = " << clb_data << std::endl;
138  std::cout << " - CL_s = " << cls_data << std::endl;
139  std::cout << " - significance of data = " << data_significance << std::endl;
140  std::cout << " - mean significance of toys = " << toys_significance << std::endl;
141 }
142 
143 
HybridCalculatorOriginal class.
void Draw(const char *options="")
Draw on current pad.
Definition: HybridPlot.cxx:131
virtual void SetSeed(UInt_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:568
RooCmdArg Extended(Bool_t flag=kTRUE)
void HybridOriginalDemo(int ntoys=1000)
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:54
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
This class provides the plots for the result of a study performed with the HybridCalculatorOriginal c...
Definition: HybridPlot.h:53
virtual HybridResult * GetHypoTest() const
inherited methods from HypoTestCalculator interface
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
void setConstant(Bool_t value=kTRUE)
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
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
const Bool_t kTRUE
Definition: Rtypes.h:91
void SetTestStatistic(int index)
set the desired test statistics: index=1 : 2 * log( L_sb / L_b ) (DEFAULT) index=2 : number of genera...