Logo ROOT   6.16/01
Reference Guide
StandardBayesianNumericalDemo.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook -js
4/// Standard demo of the numerical Bayesian calculator
5///
6/// This is a standard demo that can be used with any ROOT file
7/// prepared in the standard way. You specify:
8/// - name for input ROOT file
9/// - name of workspace inside ROOT file that holds model and data
10/// - name of ModelConfig that specifies details for calculator tools
11/// - name of dataset
12///
13/// With default parameters the macro will attempt to run the
14/// standard hist2workspace example and read the ROOT file
15/// that it produces.
16///
17/// The actual heart of the demo is only about 10 lines long.
18///
19/// The BayesianCalculator is based on Bayes's theorem
20/// and performs the integration using ROOT's numeric integration utilities
21///
22/// \macro_image
23/// \macro_output
24/// \macro_code
25///
26/// \author Kyle Cranmer
27
28#include "TFile.h"
29#include "TROOT.h"
30#include "RooWorkspace.h"
31#include "RooAbsData.h"
32#include "RooRealVar.h"
33
34#include "RooUniform.h"
39#include "RooPlot.h"
40#include "TSystem.h"
41
42#include <cassert>
43
44using namespace RooFit;
45using namespace RooStats;
46
47
48struct BayesianNumericalOptions {
49
50 double confLevel = 0.95 ; // interval CL
51 TString integrationType = ""; // integration Type (default is adaptive (numerical integration)
52 // possible values are "TOYMC" (toy MC integration, work when nuisances have a constraints pdf)
53 // "VEGAS" , "MISER", or "PLAIN" (these are all possible MC integration)
54 int nToys = 10000; // number of toys used for the MC integrations - for Vegas should be probably set to an higher value
55 bool scanPosterior = false; // flag to compute interval by scanning posterior (it is more robust but maybe less precise)
56 bool plotPosterior = false; // plot posterior function after having computed the interval
57 int nScanPoints = 50; // number of points for scanning the posterior (if scanPosterior = false it is used only for plotting). Use by default a low value to speed-up tutorial
58 int intervalType = 1; // type of interval (0 is shortest, 1 central, 2 upper limit)
59 double maxPOI = -999; // force a different value of POI for doing the scan (default is given value)
60 double nSigmaNuisance = -1; // force integration of nuisance parameters to be within nSigma of their error (do first a model fit to find nuisance error)
61
62};
63
64BayesianNumericalOptions optBayes;
65
66void StandardBayesianNumericalDemo(const char* infile = "",
67 const char* workspaceName = "combined",
68 const char* modelConfigName = "ModelConfig",
69 const char* dataName = "obsData") {
70
71 // option definitions
72 double confLevel = optBayes.confLevel;
73 TString integrationType = optBayes.integrationType;
74 int nToys = optBayes.nToys;
75 bool scanPosterior = optBayes.scanPosterior;
76 bool plotPosterior = optBayes.plotPosterior;
77 int nScanPoints = optBayes.nScanPoints;
78 int intervalType = optBayes.intervalType;
79 int maxPOI = optBayes.maxPOI;
80 double nSigmaNuisance = optBayes.nSigmaNuisance;
81
82
83
84 // -------------------------------------------------------
85 // First part is just to access a user-defined file
86 // or create the standard example file if it doesn't exist
87
88 const char* filename = "";
89 if (!strcmp(infile,"")) {
90 filename = "results/example_combined_GaussExample_model.root";
91 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
92 // if file does not exists generate with histfactory
93 if (!fileExist) {
94#ifdef _WIN32
95 cout << "HistFactory file cannot be generated on Windows - exit" << endl;
96 return;
97#endif
98 // Normally this would be run on the command line
99 cout <<"will run standard hist2workspace example"<<endl;
100 gROOT->ProcessLine(".! prepareHistFactory .");
101 gROOT->ProcessLine(".! hist2workspace config/example.xml");
102 cout <<"\n\n---------------------"<<endl;
103 cout <<"Done creating example input"<<endl;
104 cout <<"---------------------\n\n"<<endl;
105 }
106
107 }
108 else
109 filename = infile;
110
111 // Try to open the file
112 TFile *file = TFile::Open(filename);
113
114 // if input file was specified byt not found, quit
115 if(!file ){
116 cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
117 return;
118 }
119
120
121 // -------------------------------------------------------
122 // Tutorial starts here
123 // -------------------------------------------------------
124
125 // get the workspace out of the file
126 RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
127 if(!w){
128 cout <<"workspace not found" << endl;
129 return;
130 }
131
132 // get the modelConfig out of the file
133 ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
134
135 // get the modelConfig out of the file
136 RooAbsData* data = w->data(dataName);
137
138 // make sure ingredients are found
139 if(!data || !mc){
140 w->Print();
141 cout << "data or ModelConfig was not found" <<endl;
142 return;
143 }
144
145 // ------------------------------------------
146 // create and use the BayesianCalculator
147 // to find and plot the 95% credible interval
148 // on the parameter of interest as specified
149 // in the model config
150
151 // before we do that, we must specify our prior
152 // it belongs in the model config, but it may not have
153 // been specified
154 RooUniform prior("prior","",*mc->GetParametersOfInterest());
155 w->import(prior);
156 mc->SetPriorPdf(*w->pdf("prior"));
157
158 // do without systematics
159 //mc->SetNuisanceParameters(RooArgSet() );
160 if (nSigmaNuisance > 0) {
161 RooAbsPdf * pdf = mc->GetPdf();
162 assert(pdf);
165
166 res->Print();
167 RooArgList nuisPar(*mc->GetNuisanceParameters());
168 for (int i = 0; i < nuisPar.getSize(); ++i) {
169 RooRealVar * v = dynamic_cast<RooRealVar*> (&nuisPar[i] );
170 assert( v);
171 v->setMin( TMath::Max( v->getMin(), v->getVal() - nSigmaNuisance * v->getError() ) );
172 v->setMax( TMath::Min( v->getMax(), v->getVal() + nSigmaNuisance * v->getError() ) );
173 std::cout << "setting interval for nuisance " << v->GetName() << " : [ " << v->getMin() << " , " << v->getMax() << " ]" << std::endl;
174 }
175 }
176
177
178 BayesianCalculator bayesianCalc(*data,*mc);
179 bayesianCalc.SetConfidenceLevel(confLevel); // 95% interval
180
181 // default of the calculator is central interval. here use shortest , central or upper limit depending on input
182 // doing a shortest interval might require a longer time since it requires a scan of the posterior function
183 if (intervalType == 0) bayesianCalc.SetShortestInterval(); // for shortest interval
184 if (intervalType == 1) bayesianCalc.SetLeftSideTailFraction(0.5); // for central interval
185 if (intervalType == 2) bayesianCalc.SetLeftSideTailFraction(0.); // for upper limit
186
187 if (!integrationType.IsNull() ) {
188 bayesianCalc.SetIntegrationType(integrationType); // set integrationType
189 bayesianCalc.SetNumIters(nToys); // set number of iterations (i.e. number of toys for MC integrations)
190 }
191
192 // in case of toyMC make a nuisance pdf
193 if (integrationType.Contains("TOYMC") ) {
194 RooAbsPdf * nuisPdf = RooStats::MakeNuisancePdf(*mc, "nuisance_pdf");
195 cout << "using TOYMC integration: make nuisance pdf from the model " << std::endl;
196 nuisPdf->Print();
197 bayesianCalc.ForceNuisancePdf(*nuisPdf);
198 scanPosterior = true; // for ToyMC the posterior is scanned anyway so used given points
199 }
200
201 // compute interval by scanning the posterior function
202 if (scanPosterior)
203 bayesianCalc.SetScanOfPosterior(nScanPoints);
204
206 if (maxPOI != -999 && maxPOI > poi->getMin())
207 poi->setMax(maxPOI);
208
209
210 SimpleInterval* interval = bayesianCalc.GetInterval();
211
212 // print out the interval on the first Parameter of Interest
213 cout << "\n>>>> RESULT : " << confLevel*100 << "% interval on " << poi->GetName()<<" is : ["<<
214 interval->LowerLimit() << ", "<<
215 interval->UpperLimit() <<"] "<<endl;
216
217 // end in case plotting is not requested
218 if (!plotPosterior) return;
219
220 // make a plot
221 // since plotting may take a long time (it requires evaluating
222 // the posterior in many points) this command will speed up
223 // by reducing the number of points to plot - do 50
224
225 // ignore errors of PDF if is zero
227
228
229 cout << "\nDrawing plot of posterior function....." << endl;
230
231 // always plot using numer of scan points
232 bayesianCalc.SetScanOfPosterior(nScanPoints);
233
234 RooPlot * plot = bayesianCalc.GetPosteriorPlot();
235 plot->Draw();
236
237}
SVector< double, 2 > v
Definition: Dict.h:5
#define gROOT
Definition: TROOT.h:410
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
static const std::string & DefaultMinimizerType()
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
RooAbsArg * first() const
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1081
virtual Double_t getMin(const char *name=0) const
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:66
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
void setMax(const char *name, Double_t value)
Set maximum of name range to given value.
Definition: RooRealVar.cxx:417
BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation of a...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
virtual void SetPriorPdf(const RooAbsPdf &pdf)
Set the Prior Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:87
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:231
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
Definition: ModelConfig.h:234
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:228
SimpleInterval is a concrete implementation of the ConfInterval interface.
virtual Double_t UpperLimit()
virtual Double_t LowerLimit()
Flat p.d.f.
Definition: RooUniform.h:24
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
void Print(Option_t *opts=0) const
Print contents of the workspace.
Bool_t import(const RooAbsArg &arg, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg())
Import a RooAbsArg object, e.g.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseGeneralPurpose, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3975
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1286
RooCmdArg Hesse(Bool_t flag=kTRUE)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Minimizer(const char *type, const char *alg=0)
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
Definition: file.py:1