ROOT logo

From $ROOTSYS/tutorials/roostats/StandardBayesianMCMCDemo.C

// Standard demo of the Bayesian MCMC calculator
 
/*

Author: Kyle Cranmer
date: Dec. 2010
updated: July 2011 for 1-sided upper limit and SequentialProposalFunction

This is a standard demo that can be used with any ROOT file 
prepared in the standard way.  You specify:
 - name for input ROOT file
 - name of workspace inside ROOT file that holds model and data
 - name of ModelConfig that specifies details for calculator tools
 - name of dataset 

With default parameters the macro will attempt to run the 
standard hist2workspace example and read the ROOT file
that it produces.

The actual heart of the demo is only about 10 lines long.

The MCMCCalculator is a Bayesian tool that uses
the Metropolis-Hastings algorithm to efficiently integrate
in many dimensions.  It is not as accurate as the BayesianCalculator
for simple problems, but it scales to much more complicated cases.
*/

#include "TFile.h"
#include "TROOT.h"
#include "TCanvas.h"
#include "TMath.h"
#include "TSystem.h"
#include "RooWorkspace.h"
#include "RooAbsData.h"

#include "RooStats/ModelConfig.h"
#include "RooStats/MCMCCalculator.h"
#include "RooStats/MCMCInterval.h"
#include "RooStats/MCMCIntervalPlot.h"
#include "RooStats/SequentialProposal.h"
#include "RooStats/ProposalHelper.h"
#include "RooStats/ProposalHelper.h"
#include "RooFitResult.h"


using namespace RooFit;
using namespace RooStats;

void StandardBayesianMCMCDemo(const char* infile = "",
		      const char* workspaceName = "combined",
		      const char* modelConfigName = "ModelConfig",
		      const char* dataName = "obsData"){

  /////////////////////////////////////////////////////////////
  // First part is just to access a user-defined file 
  // or create the standard example file if it doesn't exist
  ////////////////////////////////////////////////////////////
   
   const char* filename = "";
   if (!strcmp(infile,"")) {
      filename = "results/example_combined_GaussExample_model.root";
      bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
      // if file does not exists generate with histfactory
      if (!fileExist) {
#ifdef _WIN32
         cout << "HistFactory file cannot be generated on Windows - exit" << endl;
         return;
#endif
         // Normally this would be run on the command line
         cout <<"will run standard hist2workspace example"<<endl;
         gROOT->ProcessLine(".! prepareHistFactory .");
         gROOT->ProcessLine(".! hist2workspace config/example.xml");
         cout <<"\n\n---------------------"<<endl;
         cout <<"Done creating example input"<<endl;
         cout <<"---------------------\n\n"<<endl;
      }
      
   }
   else
      filename = infile;
   
   // Try to open the file
   TFile *file = TFile::Open(filename);
   
   // if input file was specified byt not found, quit
   if(!file ){
      cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
      return;
   } 
 
 
  
  /////////////////////////////////////////////////////////////
  // Tutorial starts here
  ////////////////////////////////////////////////////////////

  // get the workspace out of the file
  RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
  if(!w){
    cout <<"workspace not found" << endl;
    return;
  }

  // get the modelConfig out of the file
  ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);

  // get the modelConfig out of the file
  RooAbsData* data = w->data(dataName);

  // make sure ingredients are found
  if(!data || !mc){
    w->Print();
    cout << "data or ModelConfig was not found" <<endl;
    return;
  }

  // Want an efficient proposal function
  // default is uniform.

  /*
  // this one is based on the covariance matrix of fit
  RooFitResult* fit = mc->GetPdf()->fitTo(*data,Save());
  ProposalHelper ph;
  ph.SetVariables((RooArgSet&)fit->floatParsFinal());
  ph.SetCovMatrix(fit->covarianceMatrix());
  ph.SetUpdateProposalParameters(kTRUE); // auto-create mean vars and add mappings
  ph.SetCacheSize(100);
  ProposalFunction* pf = ph.GetProposalFunction();
  */
  
  // this proposal function seems fairly robust
  SequentialProposal sp(0.1);
  /////////////////////////////////////////////
  // create and use the MCMCCalculator
  // to find and plot the 95% credible interval
  // on the parameter of interest as specified
  // in the model config
  MCMCCalculator mcmc(*data,*mc);
  mcmc.SetConfidenceLevel(0.95); // 95% interval
  //  mcmc.SetProposalFunction(*pf);
  mcmc.SetProposalFunction(sp);
  mcmc.SetNumIters(1000000);         // Metropolis-Hastings algorithm iterations
  mcmc.SetNumBurnInSteps(50);       // first N steps to be ignored as burn-in

  // default is the shortest interval.  here use central
  mcmc.SetLeftSideTailFraction(0); // for one-sided Bayesian interval

  RooRealVar* firstPOI = (RooRealVar*) mc->GetParametersOfInterest()->first();
  firstPOI->setMax(10.);

  MCMCInterval* interval = mcmc.GetInterval();

  // make a plot
  //TCanvas* c1 = 
  new TCanvas("IntervalPlot");
  MCMCIntervalPlot plot(*interval);
  plot.Draw();

  TCanvas* c2 = new TCanvas("extraPlots");
  const RooArgSet* list = mc->GetNuisanceParameters();
  if(list->getSize()>1){
    double n = list->getSize();
    int ny = TMath::CeilNint( sqrt(n) );
    int nx = TMath::CeilNint(double(n)/ny);
    c2->Divide( nx,ny);
  }

  // draw a scatter plot of chain results for poi vs each nuisance parameters
  TIterator* it = mc->GetNuisanceParameters()->createIterator();
  RooRealVar* nuis = NULL;
  int iPad=1; // iPad, that's funny
  while( (nuis = (RooRealVar*) it->Next() )){
    c2->cd(iPad++);
    plot.DrawChainScatter(*firstPOI,*nuis);
  }
  
  // print out the iterval on the first Parameter of Interest
  cout << "\n95% interval on " <<firstPOI->GetName()<<" is : ["<<
    interval->LowerLimit(*firstPOI) << ", "<<
    interval->UpperLimit(*firstPOI) <<"] "<<endl;

}
 StandardBayesianMCMCDemo.C:1
 StandardBayesianMCMCDemo.C:2
 StandardBayesianMCMCDemo.C:3
 StandardBayesianMCMCDemo.C:4
 StandardBayesianMCMCDemo.C:5
 StandardBayesianMCMCDemo.C:6
 StandardBayesianMCMCDemo.C:7
 StandardBayesianMCMCDemo.C:8
 StandardBayesianMCMCDemo.C:9
 StandardBayesianMCMCDemo.C:10
 StandardBayesianMCMCDemo.C:11
 StandardBayesianMCMCDemo.C:12
 StandardBayesianMCMCDemo.C:13
 StandardBayesianMCMCDemo.C:14
 StandardBayesianMCMCDemo.C:15
 StandardBayesianMCMCDemo.C:16
 StandardBayesianMCMCDemo.C:17
 StandardBayesianMCMCDemo.C:18
 StandardBayesianMCMCDemo.C:19
 StandardBayesianMCMCDemo.C:20
 StandardBayesianMCMCDemo.C:21
 StandardBayesianMCMCDemo.C:22
 StandardBayesianMCMCDemo.C:23
 StandardBayesianMCMCDemo.C:24
 StandardBayesianMCMCDemo.C:25
 StandardBayesianMCMCDemo.C:26
 StandardBayesianMCMCDemo.C:27
 StandardBayesianMCMCDemo.C:28
 StandardBayesianMCMCDemo.C:29
 StandardBayesianMCMCDemo.C:30
 StandardBayesianMCMCDemo.C:31
 StandardBayesianMCMCDemo.C:32
 StandardBayesianMCMCDemo.C:33
 StandardBayesianMCMCDemo.C:34
 StandardBayesianMCMCDemo.C:35
 StandardBayesianMCMCDemo.C:36
 StandardBayesianMCMCDemo.C:37
 StandardBayesianMCMCDemo.C:38
 StandardBayesianMCMCDemo.C:39
 StandardBayesianMCMCDemo.C:40
 StandardBayesianMCMCDemo.C:41
 StandardBayesianMCMCDemo.C:42
 StandardBayesianMCMCDemo.C:43
 StandardBayesianMCMCDemo.C:44
 StandardBayesianMCMCDemo.C:45
 StandardBayesianMCMCDemo.C:46
 StandardBayesianMCMCDemo.C:47
 StandardBayesianMCMCDemo.C:48
 StandardBayesianMCMCDemo.C:49
 StandardBayesianMCMCDemo.C:50
 StandardBayesianMCMCDemo.C:51
 StandardBayesianMCMCDemo.C:52
 StandardBayesianMCMCDemo.C:53
 StandardBayesianMCMCDemo.C:54
 StandardBayesianMCMCDemo.C:55
 StandardBayesianMCMCDemo.C:56
 StandardBayesianMCMCDemo.C:57
 StandardBayesianMCMCDemo.C:58
 StandardBayesianMCMCDemo.C:59
 StandardBayesianMCMCDemo.C:60
 StandardBayesianMCMCDemo.C:61
 StandardBayesianMCMCDemo.C:62
 StandardBayesianMCMCDemo.C:63
 StandardBayesianMCMCDemo.C:64
 StandardBayesianMCMCDemo.C:65
 StandardBayesianMCMCDemo.C:66
 StandardBayesianMCMCDemo.C:67
 StandardBayesianMCMCDemo.C:68
 StandardBayesianMCMCDemo.C:69
 StandardBayesianMCMCDemo.C:70
 StandardBayesianMCMCDemo.C:71
 StandardBayesianMCMCDemo.C:72
 StandardBayesianMCMCDemo.C:73
 StandardBayesianMCMCDemo.C:74
 StandardBayesianMCMCDemo.C:75
 StandardBayesianMCMCDemo.C:76
 StandardBayesianMCMCDemo.C:77
 StandardBayesianMCMCDemo.C:78
 StandardBayesianMCMCDemo.C:79
 StandardBayesianMCMCDemo.C:80
 StandardBayesianMCMCDemo.C:81
 StandardBayesianMCMCDemo.C:82
 StandardBayesianMCMCDemo.C:83
 StandardBayesianMCMCDemo.C:84
 StandardBayesianMCMCDemo.C:85
 StandardBayesianMCMCDemo.C:86
 StandardBayesianMCMCDemo.C:87
 StandardBayesianMCMCDemo.C:88
 StandardBayesianMCMCDemo.C:89
 StandardBayesianMCMCDemo.C:90
 StandardBayesianMCMCDemo.C:91
 StandardBayesianMCMCDemo.C:92
 StandardBayesianMCMCDemo.C:93
 StandardBayesianMCMCDemo.C:94
 StandardBayesianMCMCDemo.C:95
 StandardBayesianMCMCDemo.C:96
 StandardBayesianMCMCDemo.C:97
 StandardBayesianMCMCDemo.C:98
 StandardBayesianMCMCDemo.C:99
 StandardBayesianMCMCDemo.C:100
 StandardBayesianMCMCDemo.C:101
 StandardBayesianMCMCDemo.C:102
 StandardBayesianMCMCDemo.C:103
 StandardBayesianMCMCDemo.C:104
 StandardBayesianMCMCDemo.C:105
 StandardBayesianMCMCDemo.C:106
 StandardBayesianMCMCDemo.C:107
 StandardBayesianMCMCDemo.C:108
 StandardBayesianMCMCDemo.C:109
 StandardBayesianMCMCDemo.C:110
 StandardBayesianMCMCDemo.C:111
 StandardBayesianMCMCDemo.C:112
 StandardBayesianMCMCDemo.C:113
 StandardBayesianMCMCDemo.C:114
 StandardBayesianMCMCDemo.C:115
 StandardBayesianMCMCDemo.C:116
 StandardBayesianMCMCDemo.C:117
 StandardBayesianMCMCDemo.C:118
 StandardBayesianMCMCDemo.C:119
 StandardBayesianMCMCDemo.C:120
 StandardBayesianMCMCDemo.C:121
 StandardBayesianMCMCDemo.C:122
 StandardBayesianMCMCDemo.C:123
 StandardBayesianMCMCDemo.C:124
 StandardBayesianMCMCDemo.C:125
 StandardBayesianMCMCDemo.C:126
 StandardBayesianMCMCDemo.C:127
 StandardBayesianMCMCDemo.C:128
 StandardBayesianMCMCDemo.C:129
 StandardBayesianMCMCDemo.C:130
 StandardBayesianMCMCDemo.C:131
 StandardBayesianMCMCDemo.C:132
 StandardBayesianMCMCDemo.C:133
 StandardBayesianMCMCDemo.C:134
 StandardBayesianMCMCDemo.C:135
 StandardBayesianMCMCDemo.C:136
 StandardBayesianMCMCDemo.C:137
 StandardBayesianMCMCDemo.C:138
 StandardBayesianMCMCDemo.C:139
 StandardBayesianMCMCDemo.C:140
 StandardBayesianMCMCDemo.C:141
 StandardBayesianMCMCDemo.C:142
 StandardBayesianMCMCDemo.C:143
 StandardBayesianMCMCDemo.C:144
 StandardBayesianMCMCDemo.C:145
 StandardBayesianMCMCDemo.C:146
 StandardBayesianMCMCDemo.C:147
 StandardBayesianMCMCDemo.C:148
 StandardBayesianMCMCDemo.C:149
 StandardBayesianMCMCDemo.C:150
 StandardBayesianMCMCDemo.C:151
 StandardBayesianMCMCDemo.C:152
 StandardBayesianMCMCDemo.C:153
 StandardBayesianMCMCDemo.C:154
 StandardBayesianMCMCDemo.C:155
 StandardBayesianMCMCDemo.C:156
 StandardBayesianMCMCDemo.C:157
 StandardBayesianMCMCDemo.C:158
 StandardBayesianMCMCDemo.C:159
 StandardBayesianMCMCDemo.C:160
 StandardBayesianMCMCDemo.C:161
 StandardBayesianMCMCDemo.C:162
 StandardBayesianMCMCDemo.C:163
 StandardBayesianMCMCDemo.C:164
 StandardBayesianMCMCDemo.C:165
 StandardBayesianMCMCDemo.C:166
 StandardBayesianMCMCDemo.C:167
 StandardBayesianMCMCDemo.C:168
 StandardBayesianMCMCDemo.C:169
 StandardBayesianMCMCDemo.C:170
 StandardBayesianMCMCDemo.C:171
 StandardBayesianMCMCDemo.C:172
 StandardBayesianMCMCDemo.C:173
 StandardBayesianMCMCDemo.C:174
 StandardBayesianMCMCDemo.C:175
 StandardBayesianMCMCDemo.C:176
 StandardBayesianMCMCDemo.C:177
 StandardBayesianMCMCDemo.C:178
 StandardBayesianMCMCDemo.C:179
 StandardBayesianMCMCDemo.C:180
 StandardBayesianMCMCDemo.C:181
 StandardBayesianMCMCDemo.C:182
 StandardBayesianMCMCDemo.C:183