ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
StandardBayesianNumericalDemo.C
Go to the documentation of this file.
1 // Standard demo of the numerical Bayesian calculator
2 
3 /*
4 
5 
6 Author: Kyle Cranmer
7 date: Dec. 2010
8 
9 This is a standard demo that can be used with any ROOT file
10 prepared in the standard way. You specify:
11  - name for input ROOT file
12  - name of workspace inside ROOT file that holds model and data
13  - name of ModelConfig that specifies details for calculator tools
14  - name of dataset
15 
16 With default parameters the macro will attempt to run the
17 standard hist2workspace example and read the ROOT file
18 that it produces.
19 
20 The actual heart of the demo is only about 10 lines long.
21 
22 The BayesianCalculator is based on Bayes's theorem
23 and performs the integration using ROOT's numeric integration utilities
24 */
25 
26 #include "TFile.h"
27 #include "TROOT.h"
28 #include "RooWorkspace.h"
29 #include "RooAbsData.h"
30 #include "RooRealVar.h"
31 
32 #include "RooUniform.h"
33 #include "RooStats/ModelConfig.h"
36 #include "RooStats/RooStatsUtils.h"
37 #include "RooPlot.h"
38 #include "TSystem.h"
39 
40 using namespace RooFit;
41 using namespace RooStats;
42 
43 
44 TString integrationType = ""; // integration Type (default is adaptive (numerical integration)
45  // possible values are "TOYMC" (toy MC integration, work when nuisances have a constraints pdf)
46  // "VEGAS" , "MISER", or "PLAIN" (these are all possible MC integration)
47 int nToys = 10000; // number of toys used for the MC integrations - for Vegas should be probably set to an higher value
48 bool scanPosterior = false; // flag to compute interval by scanning posterior (it is more robust but maybe less precise)
49 int nScanPoints = 20; // 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
50 int intervalType = 1; // type of interval (0 is shortest, 1 central, 2 upper limit)
51 double maxPOI = -999; // force a different value of POI for doing the scan (default is given value)
52 double nSigmaNuisance = -1; // force integration of nuisance parameters to be withing nSigma of their error (do first a model fit to find nuisance error)
53 
54 void StandardBayesianNumericalDemo(const char* infile = "",
55  const char* workspaceName = "combined",
56  const char* modelConfigName = "ModelConfig",
57  const char* dataName = "obsData") {
58 
59  /////////////////////////////////////////////////////////////
60  // First part is just to access a user-defined file
61  // or create the standard example file if it doesn't exist
62  ////////////////////////////////////////////////////////////
63 
64  const char* filename = "";
65  if (!strcmp(infile,"")) {
66  filename = "results/example_combined_GaussExample_model.root";
67  bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
68  // if file does not exists generate with histfactory
69  if (!fileExist) {
70 #ifdef _WIN32
71  cout << "HistFactory file cannot be generated on Windows - exit" << endl;
72  return;
73 #endif
74  // Normally this would be run on the command line
75  cout <<"will run standard hist2workspace example"<<endl;
76  gROOT->ProcessLine(".! prepareHistFactory .");
77  gROOT->ProcessLine(".! hist2workspace config/example.xml");
78  cout <<"\n\n---------------------"<<endl;
79  cout <<"Done creating example input"<<endl;
80  cout <<"---------------------\n\n"<<endl;
81  }
82 
83  }
84  else
85  filename = infile;
86 
87  // Try to open the file
88  TFile *file = TFile::Open(filename);
89 
90  // if input file was specified byt not found, quit
91  if(!file ){
92  cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
93  return;
94  }
95 
96 
97  /////////////////////////////////////////////////////////////
98  // Tutorial starts here
99  ////////////////////////////////////////////////////////////
100 
101  // get the workspace out of the file
102  RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
103  if(!w){
104  cout <<"workspace not found" << endl;
105  return;
106  }
107 
108  // get the modelConfig out of the file
109  ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
110 
111  // get the modelConfig out of the file
112  RooAbsData* data = w->data(dataName);
113 
114  // make sure ingredients are found
115  if(!data || !mc){
116  w->Print();
117  cout << "data or ModelConfig was not found" <<endl;
118  return;
119  }
120 
121  /////////////////////////////////////////////
122  // create and use the BayesianCalculator
123  // to find and plot the 95% credible interval
124  // on the parameter of interest as specified
125  // in the model config
126 
127  // before we do that, we must specify our prior
128  // it belongs in the model config, but it may not have
129  // been specified
130  RooUniform prior("prior","",*mc->GetParametersOfInterest());
131  w->import(prior);
132  mc->SetPriorPdf(*w->pdf("prior"));
133 
134  // do without systematics
135  //mc->SetNuisanceParameters(RooArgSet() );
136  if (nSigmaNuisance > 0) {
137  RooAbsPdf * pdf = mc->GetPdf();
138  assert(pdf);
139  RooFitResult * res = pdf->fitTo(*data, Save(true), Minimizer(ROOT::Math::MinimizerOptions::DefaultMinimizerType().c_str()), Hesse(true),
141 
142  res->Print();
143  RooArgList nuisPar(*mc->GetNuisanceParameters());
144  for (int i = 0; i < nuisPar.getSize(); ++i) {
145  RooRealVar * v = dynamic_cast<RooRealVar*> (&nuisPar[i] );
146  assert( v);
147  v->setMin( TMath::Max( v->getMin(), v->getVal() - nSigmaNuisance * v->getError() ) );
148  v->setMax( TMath::Min( v->getMax(), v->getVal() + nSigmaNuisance * v->getError() ) );
149  std::cout << "setting interval for nuisance " << v->GetName() << " : [ " << v->getMin() << " , " << v->getMax() << " ]" << std::endl;
150  }
151  }
152 
153 
154  BayesianCalculator bayesianCalc(*data,*mc);
155  bayesianCalc.SetConfidenceLevel(0.95); // 95% interval
156 
157  // default of the calculator is central interval. here use shortest , central or upper limit depending on input
158  // doing a shortest interval might require a longer time since it requires a scan of the posterior function
159  if (intervalType == 0) bayesianCalc.SetShortestInterval(); // for shortest interval
160  if (intervalType == 1) bayesianCalc.SetLeftSideTailFraction(0.5); // for central interval
161  if (intervalType == 2) bayesianCalc.SetLeftSideTailFraction(0.); // for upper limit
162 
163  if (!integrationType.IsNull() ) {
164  bayesianCalc.SetIntegrationType(integrationType); // set integrationType
165  bayesianCalc.SetNumIters(nToys); // set number of ietrations (i.e. number of toys for MC integrations)
166  }
167 
168  // in case of toyMC make a nnuisance pdf
169  if (integrationType.Contains("TOYMC") ) {
170  RooAbsPdf * nuisPdf = RooStats::MakeNuisancePdf(*mc, "nuisance_pdf");
171  cout << "using TOYMC integration: make nuisance pdf from the model " << std::endl;
172  nuisPdf->Print();
173  bayesianCalc.ForceNuisancePdf(*nuisPdf);
174  scanPosterior = true; // for ToyMC the posterior is scanned anyway so used given points
175  }
176 
177  // compute interval by scanning the posterior function
178  if (scanPosterior)
179  bayesianCalc.SetScanOfPosterior(nScanPoints);
180 
182  if (maxPOI != -999 && maxPOI > poi->getMin())
183  poi->setMax(maxPOI);
184 
185 
186  SimpleInterval* interval = bayesianCalc.GetInterval();
187 
188  // print out the iterval on the first Parameter of Interest
189  cout << "\n95% interval on " << poi->GetName()<<" is : ["<<
190  interval->LowerLimit() << ", "<<
191  interval->UpperLimit() <<"] "<<endl;
192 
193 
194  // make a plot
195  // since plotting may take a long time (it requires evaluating
196  // the posterior in many points) this command will speed up
197  // by reducing the number of points to plot - do 50
198 
199  cout << "\nDrawing plot of posterior function....." << endl;
200 
201  bayesianCalc.SetScanOfPosterior(nScanPoints);
202 
203  RooPlot * plot = bayesianCalc.GetPosteriorPlot();
204  plot->Draw();
205 
206 }
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:65
RooPlot * GetPosteriorPlot(bool norm=false, double precision=0.01) const
get the plot with option to get it normalized
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:1213
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:52
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:244
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found...
void SetScanOfPosterior(int nbin=100)
use directly the approximate posterior function obtained by binning it in nbins by default the cdf is...
RooCmdArg PrintLevel(Int_t code)
#define assert(cond)
Definition: unittest.h:542
void SetIntegrationType(const char *type)
set the integration type (possible type are) :
static const char * filename()
#define gROOT
Definition: TROOT.h:344
Basic string class.
Definition: TString.h:137
virtual Double_t getMin(const char *name=0) const
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
Definition: ModelConfig.h:250
void setMax(const char *name, Double_t value)
Set maximum of name range to given value.
Definition: RooRealVar.cxx:416
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3851
virtual SimpleInterval * GetInterval() const
compute the interval.
RooAbsArg * first() const
void plot(TString fname="data.root", TString var0="var0", TString var1="var1")
Definition: createData.C:17
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
virtual void SetNumIters(Int_t numIters)
set the number of iterations when running a MC integration algorithm If not set use default algorithm...
virtual Double_t LowerLimit()
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval) ...
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
tuple infile
Definition: mrt.py:15
void SetShortestInterval()
set the Bayesian calculator to compute the shorest interval (default is central interval) to switch o...
void SetLeftSideTailFraction(Double_t leftSideFraction)
set the fraction of probability content on the left tail Central limits use 0.5 (default case) for up...
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
static const std::string & DefaultMinimizerType()
SVector< double, 2 > v
Definition: Dict.h:5
RooCmdArg Minimizer(const char *type, const char *alg=0)
void setMin(const char *name, Double_t value)
Set minimum of name range to given value.
Definition: RooRealVar.cxx:386
tuple w
Definition: qtexample.py:51
double nSigmaNuisance
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
Bool_t IsNull() const
Definition: TString.h:387
void StandardBayesianNumericalDemo(const char *infile="", const char *workspaceName="combined", const char *modelConfigName="ModelConfig", const char *dataName="obsData")
Flat p.d.f.
Definition: RooUniform.h:24
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
TString integrationType
RooCmdArg Hesse(Bool_t flag=kTRUE)
tuple file
Definition: fildir.py:20
virtual void SetPriorPdf(const RooAbsPdf &pdf)
Set the Prior Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:103
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name) ...
virtual Double_t getMax(const char *name=0) const
RooCmdArg Save(Bool_t flag=kTRUE)
virtual Double_t UpperLimit()
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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.
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
void Print(Option_t *opts=0) const
Print contents of the workspace.
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
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:1056
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:247
void ForceNuisancePdf(RooAbsPdf &pdf)
force the nuisance pdf when using the toy mc sampling
BayesianCalculator is a concrete implementation of IntervalCalculator, providing the computation of a...
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559
Double_t getError() const
Definition: RooRealVar.h:54