StandardBayesianMCMCDemo.C

Go to the documentation of this file.
00001 // Standard demo of the Bayesian MCMC calculator
00002  
00003 /*
00004 
00005 Author: Kyle Cranmer
00006 date: Dec. 2010
00007 
00008 This is a standard demo that can be used with any ROOT file 
00009 prepared in the standard way.  You specify:
00010  - name for input ROOT file
00011  - name of workspace inside ROOT file that holds model and data
00012  - name of ModelConfig that specifies details for calculator tools
00013  - name of dataset 
00014 
00015 With default parameters the macro will attempt to run the 
00016 standard hist2workspace example and read the ROOT file
00017 that it produces.
00018 
00019 The actual heart of the demo is only about 10 lines long.
00020 
00021 The MCMCCalculator is a Bayesian tool that uses
00022 the Metropolis-Hastings algorithm to efficiently integrate
00023 in many dimensions.  It is not as accurate as the BayesianCalculator
00024 for simple problems, but it scales to much more complicated cases.
00025 */
00026 
00027 #include "TFile.h"
00028 #include "TROOT.h"
00029 #include "RooWorkspace.h"
00030 #include "RooAbsData.h"
00031 
00032 #include "RooStats/ModelConfig.h"
00033 #include "RooStats/MCMCCalculator.h"
00034 #include "RooStats/MCMCInterval.h"
00035 #include "RooStats/MCMCIntervalPlot.h"
00036 
00037 using namespace RooFit;
00038 using namespace RooStats;
00039 
00040 void StandardBayesianMCMCDemo(const char* infile = "",
00041                       const char* workspaceName = "proto",
00042                       const char* modelConfigName = "ModelConfig",
00043                       const char* dataName = "expData"){
00044 
00045   /////////////////////////////////////////////////////////////
00046   // First part is just to access a user-defined file 
00047   // or create the standard example file if it doesn't exist
00048   ////////////////////////////////////////////////////////////
00049   const char* filename = "";
00050   if (!strcmp(infile,""))
00051     filename = "results/example_channel1_GammaExample_model.root";
00052   else
00053     filename = infile;
00054   // Check if example input file exists
00055   TFile *file = TFile::Open(filename);
00056 
00057   // if input file was specified byt not found, quit
00058   if(!file && strcmp(infile,"")){
00059     cout <<"file not found" << endl;
00060     return;
00061   } 
00062 
00063   // if default file not found, try to create it
00064   if(!file ){
00065     // Normally this would be run on the command line
00066     cout <<"will run standard hist2workspace example"<<endl;
00067     gROOT->ProcessLine(".! prepareHistFactory .");
00068     gROOT->ProcessLine(".! hist2workspace config/example.xml");
00069     cout <<"\n\n---------------------"<<endl;
00070     cout <<"Done creating example input"<<endl;
00071     cout <<"---------------------\n\n"<<endl;
00072   }
00073 
00074   // now try to access the file again
00075   file = TFile::Open(filename);
00076   if(!file){
00077     // if it is still not there, then we can't continue
00078     cout << "Not able to run hist2workspace to create example input" <<endl;
00079     return;
00080   }
00081 
00082   
00083   /////////////////////////////////////////////////////////////
00084   // Tutorial starts here
00085   ////////////////////////////////////////////////////////////
00086 
00087   // get the workspace out of the file
00088   RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
00089   if(!w){
00090     cout <<"workspace not found" << endl;
00091     return;
00092   }
00093 
00094   // get the modelConfig out of the file
00095   ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
00096 
00097   // get the modelConfig out of the file
00098   RooAbsData* data = w->data(dataName);
00099 
00100   // make sure ingredients are found
00101   if(!data || !mc){
00102     w->Print();
00103     cout << "data or ModelConfig was not found" <<endl;
00104     return;
00105   }
00106 
00107   /////////////////////////////////////////////
00108   // create and use the MCMCCalculator
00109   // to find and plot the 95% credible interval
00110   // on the parameter of interest as specified
00111   // in the model config
00112   MCMCCalculator mcmc(*data,*mc);
00113   mcmc.SetConfidenceLevel(0.95); // 95% interval
00114   mcmc.SetNumIters(1000000);         // Metropolis-Hastings algorithm iterations
00115   mcmc.SetNumBurnInSteps(500);       // first N steps to be ignored as burn-in
00116 
00117   // default is the shortest interval.  here use central
00118   mcmc.SetLeftSideTailFraction(0.5); // for central interval
00119 
00120   MCMCInterval* interval = mcmc.GetInterval();
00121 
00122   // make a plot
00123   MCMCIntervalPlot plot(*interval);
00124   plot.Draw();
00125   
00126   // print out the iterval on the first Parameter of Interest
00127   RooRealVar* firstPOI = (RooRealVar*) mc->GetParametersOfInterest()->first();
00128   cout << "\n95% interval on " <<firstPOI->GetName()<<" is : ["<<
00129     interval->LowerLimit(*firstPOI) << ", "<<
00130     interval->UpperLimit(*firstPOI) <<"] "<<endl;
00131 
00132 }

Generated on Tue Jul 5 15:45:10 2011 for ROOT_528-00b_version by  doxygen 1.5.1