ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MultivariateGaussianTest.C
Go to the documentation of this file.
1 // comparison of MCMC and PLC in a multi-variate gaussian problem
2 
3 #include "RooGlobalFunc.h"
4 #include <stdlib.h>
5 #include "TMatrixDSym.h"
6 #include "RooMultiVarGaussian.h"
7 #include "RooArgList.h"
8 #include "RooRealVar.h"
9 #include "TH2F.h"
10 #include "TCanvas.h"
11 #include "RooAbsReal.h"
12 #include "RooFitResult.h"
13 #include "TStopwatch.h"
16 #include "RooStats/MarkovChain.h"
17 #include "RooStats/ConfInterval.h"
18 #include "RooStats/MCMCInterval.h"
20 #include "RooStats/ModelConfig.h"
23 #include "RooStats/PdfProposal.h"
27 
28 using namespace std;
29 using namespace RooFit;
30 using namespace RooStats;
31 
32 
33 void MultivariateGaussianTest(Int_t dim = 4, Int_t nPOI = 2)
34 {
35  /*
36  Authors: Kevin Belasco and Kyle Cranmer.
37 
38  This tutorial produces an N-dimensional multivariate Gaussian
39  with a non-trivial covariance matrix. By default N=4 (called "dim").
40 
41  A subset of these are considered parameters of interest.
42  This problem is tractable analytically.
43 
44  We use this mainly as a test of Markov Chain Monte Carlo
45  and we compare the result to the profile likelihood ratio.
46 
47  We use the proposal helper to create a customized
48  proposal function for this problem.
49 
50  For N=4 and 2 parameters of interest it takes about 10-20 seconds
51  and the acceptance rate is 37%
52  */
53 
54  // let's time this challenging example
55  TStopwatch t;
56  t.Start();
57 
58  RooArgList xVec;
59  RooArgList muVec;
60  RooArgSet poi;
61 
62  // make the observable and means
63  Int_t i,j;
64  RooRealVar* x;
65  RooRealVar* mu_x;
66  for (i = 0; i < dim; i++) {
67  char* name = Form("x%d", i);
68  x = new RooRealVar(name, name, 0, -3,3);
69  xVec.add(*x);
70 
71  char* mu_name = Form("mu_x%d",i);
72  mu_x = new RooRealVar(mu_name, mu_name, 0, -2,2);
73  muVec.add(*mu_x);
74  }
75 
76  // put them into the list of parameters of interest
77  for (i = 0; i < nPOI; i++) {
78  poi.add(*muVec.at(i));
79  }
80 
81  // make a covariance matrix that is all 1's
82  TMatrixDSym cov(dim);
83  for (i = 0; i < dim; i++) {
84  for (j = 0; j < dim; j++) {
85  if (i == j) cov(i,j) = 3.;
86  else cov(i,j) = 1.0;
87  }
88  }
89 
90  // now make the multivariate Gaussian
91  RooMultiVarGaussian mvg("mvg", "mvg", xVec, muVec, cov);
92 
93  ///////////////////////////////////////////
94  // make a toy dataset
95  RooDataSet* data = mvg.generate(xVec, 100);
96 
97  ///////////////////////////////////////////////
98  // now create the model config for this problem
99  RooWorkspace* w = new RooWorkspace("MVG");
100  ModelConfig modelConfig(w);
101  modelConfig.SetPdf(mvg);
102  modelConfig.SetParametersOfInterest(poi);
103 
104  ///////////////////////////////////////////
105  // Setup calculators
106 
107  // MCMC
108  // we want to setup an efficient proposal function
109  // using the covariance matrix from a fit to the data
110  RooFitResult* fit = mvg.fitTo(*data, Save(true));
111  ProposalHelper ph;
113  ph.SetCovMatrix(fit->covarianceMatrix());
115  ph.SetCacheSize(100);
116  ProposalFunction* pdfProp = ph.GetProposalFunction();
117 
118  // now create the calculator
119  MCMCCalculator mc(*data, modelConfig);
120  mc.SetConfidenceLevel(0.95);
121  mc.SetNumBurnInSteps(100);
122  mc.SetNumIters(10000);
123  mc.SetNumBins(50);
124  mc.SetProposalFunction(*pdfProp);
125 
126  MCMCInterval* mcInt = mc.GetInterval();
127  RooArgList* poiList = mcInt->GetAxes();
128 
129  // now setup the profile likelihood calculator
130  ProfileLikelihoodCalculator plc(*data, modelConfig);
131  plc.SetConfidenceLevel(0.95);
133 
134 
135  ///////////////////////////////////////////
136  // make some plots
137  MCMCIntervalPlot mcPlot(*mcInt);
138 
139  TCanvas* c1 = new TCanvas();
140  mcPlot.SetLineColor(kGreen);
141  mcPlot.SetLineWidth(2);
142  mcPlot.Draw();
143 
144  LikelihoodIntervalPlot plPlot(plInt);
145  plPlot.Draw("same");
146 
147  if (poiList->getSize() == 1) {
148  RooRealVar* p = (RooRealVar*)poiList->at(0);
149  Double_t ll = mcInt->LowerLimit(*p);
150  Double_t ul = mcInt->UpperLimit(*p);
151  cout << "MCMC interval: [" << ll << ", " << ul << "]" << endl;
152  }
153 
154  if (poiList->getSize() == 2) {
155  RooRealVar* p0 = (RooRealVar*)poiList->at(0);
156  RooRealVar* p1 = (RooRealVar*)poiList->at(1);
157  TCanvas* scatter = new TCanvas();
158  Double_t ll = mcInt->LowerLimit(*p0);
159  Double_t ul = mcInt->UpperLimit(*p0);
160  cout << "MCMC interval on p0: [" << ll << ", " << ul << "]" << endl;
161  ll = mcInt->LowerLimit(*p0);
162  ul = mcInt->UpperLimit(*p0);
163  cout << "MCMC interval on p1: [" << ll << ", " << ul << "]" << endl;
164 
165  //MCMC interval on p0: [-0.2, 0.6]
166  //MCMC interval on p1: [-0.2, 0.6]
167 
168  mcPlot.DrawChainScatter(*p0, *p1);
169  scatter->Update();
170  }
171 
172  t.Print();
173 
174 }
175 
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:217
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:52
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:56
void Draw(const Option_t *options=0)
draw the likelihood interval or contour for the 1D case a RooPlot is drawn by default of the profiled...
This class provides simple and straightforward utilities to plot a MCMCInterval object.
const RooArgList & floatParsFinal() const
Definition: RooFitResult.h:109
LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
TCanvas * c1
Definition: legend1.C:2
virtual ProposalFunction * GetProposalFunction()
virtual void SetUpdateProposalParameters(Bool_t updateParams)
This class provides simple and straightforward utilities to plot a LikelihoodInterval object...
ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface class f...
virtual void SetNumBurnInSteps(Int_t numBurnInSteps)
set the number of steps in the chain to discard as burn-in, starting from the first ...
virtual void SetNumBins(Int_t numBins)
set the number of bins to create for each axis when constructing the interval
int Int_t
Definition: RtypesCore.h:41
void Draw(const Option_t *options=NULL)
Definition: Rtypes.h:61
virtual LikelihoodInterval * GetInterval() const
Return a likelihood interval.
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval) ...
Double_t x[n]
Definition: legend1.C:17
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:97
void MultivariateGaussianTest(Int_t dim=4, Int_t nPOI=2)
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
TThread * t[5]
Definition: threadsh1.C:13
virtual void SetNumIters(Int_t numIters)
set the number of iterations to run the metropolis algorithm
virtual MCMCInterval * GetInterval() const
Main interface to get a ConfInterval.
virtual void SetCacheSize(Int_t size)
char * Form(const char *fmt,...)
void SetLineColor(Color_t color)
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
tuple w
Definition: qtexample.py:51
void SetLineWidth(Int_t width)
static double p1(double t, double a, double b)
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
Multivariate Gaussian p.d.f.
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
virtual void SetProposalFunction(ProposalFunction &proposalFunction)
set the proposal function for suggesting new points for the MCMC
double Double_t
Definition: RtypesCore.h:55
void DrawChainScatter(RooRealVar &xVar, RooRealVar &yVar)
virtual void SetCovMatrix(const TMatrixDSym &covMatrix)
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval) ...
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooCmdArg Save(Bool_t flag=kTRUE)
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
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
virtual void SetParametersOfInterest(const RooArgSet &set)
Definition: ModelConfig.h:115
Int_t getSize() const
virtual void SetVariables(RooArgList &vars)
Bayesian Calculator estimating an interval or a credible region using the Markov-Chain Monte Carlo me...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add element to non-owning set.
Definition: RooArgSet.cxx:448
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
Stopwatch class.
Definition: TStopwatch.h:30