ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rs101_limitexample.C
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'Limit Example' RooStats tutorial macro #101
4 // author: Kyle Cranmer
5 // date June. 2009
6 //
7 // This tutorial shows an example of creating a simple
8 // model for a number counting experiment with uncertainty
9 // on both the background rate and signal efficeincy. We then
10 // use a Confidence Interval Calculator to set a limit on the signal.
11 //
12 //
13 /////////////////////////////////////////////////////////////////////////
14 
15 #ifndef __CINT__
16 #include "RooGlobalFunc.h"
17 #endif
18 
19 #include "RooProfileLL.h"
20 #include "RooAbsPdf.h"
22 #include "RooRealVar.h"
23 #include "RooPlot.h"
24 #include "RooDataSet.h"
25 #include "RooTreeDataStore.h"
26 #include "TTree.h"
27 #include "TCanvas.h"
28 #include "TLine.h"
29 #include "TStopwatch.h"
30 
36 #include "RooStats/ConfInterval.h"
40 #include "RooStats/RooStatsUtils.h"
41 #include "RooStats/ModelConfig.h"
42 #include "RooStats/MCMCInterval.h"
46 #include "RooFitResult.h"
47 #include "TGraph2D.h"
48 
49 // use this order for safety on library loading
50 using namespace RooFit ;
51 using namespace RooStats ;
52 
53 
55 {
56  /////////////////////////////////////////
57  // An example of setting a limit in a number counting experiment with uncertainty on background and signal
58  /////////////////////////////////////////
59 
60  // to time the macro
61  TStopwatch t;
62  t.Start();
63 
64  /////////////////////////////////////////
65  // The Model building stage
66  /////////////////////////////////////////
67  RooWorkspace* wspace = new RooWorkspace();
68  wspace->factory("Poisson::countingModel(obs[150,0,300], sum(s[50,0,120]*ratioSigEff[1.,0,3.],b[100]*ratioBkgEff[1.,0.,3.]))"); // counting model
69  // wspace->factory("Gaussian::sigConstraint(ratioSigEff,1,0.05)"); // 5% signal efficiency uncertainty
70  // wspace->factory("Gaussian::bkgConstraint(ratioBkgEff,1,0.1)"); // 10% background efficiency uncertainty
71  wspace->factory("Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)"); // 5% signal efficiency uncertainty
72  wspace->factory("Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)"); // 10% background efficiency uncertainty
73  wspace->factory("PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)"); // product of terms
74  wspace->Print();
75 
76  RooAbsPdf* modelWithConstraints = wspace->pdf("modelWithConstraints"); // get the model
77  RooRealVar* obs = wspace->var("obs"); // get the observable
78  RooRealVar* s = wspace->var("s"); // get the signal we care about
79  RooRealVar* b = wspace->var("b"); // get the background and set it to a constant. Uncertainty included in ratioBkgEff
80  b->setConstant();
81 
82  RooRealVar* ratioSigEff = wspace->var("ratioSigEff"); // get uncertaint parameter to constrain
83  RooRealVar* ratioBkgEff = wspace->var("ratioBkgEff"); // get uncertaint parameter to constrain
84  RooArgSet constrainedParams(*ratioSigEff, *ratioBkgEff); // need to constrain these in the fit (should change default behavior)
85 
86  RooRealVar * gSigEff = wspace->var("gSigEff"); // global observables for signal efficiency
87  RooRealVar * gSigBkg = wspace->var("gSigBkg"); // global obs for background efficiency
88  gSigEff->setConstant();
89  gSigBkg->setConstant();
90 
91  // Create an example dataset with 160 observed events
92  obs->setVal(160.);
93  RooDataSet* data = new RooDataSet("exampleData", "exampleData", RooArgSet(*obs));
94  data->add(*obs);
95 
96  RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);
97 
98  // not necessary
99  modelWithConstraints->fitTo(*data, RooFit::Constrain(RooArgSet(*ratioSigEff, *ratioBkgEff)));
100 
101  // Now let's make some confidence intervals for s, our parameter of interest
102  RooArgSet paramOfInterest(*s);
103 
104  ModelConfig modelConfig(wspace);
105  modelConfig.SetPdf(*modelWithConstraints);
106  modelConfig.SetParametersOfInterest(paramOfInterest);
107  modelConfig.SetNuisanceParameters(constrainedParams);
108  modelConfig.SetObservables(*obs);
109  modelConfig.SetGlobalObservables( RooArgSet(*gSigEff,*gSigBkg));
110  modelConfig.SetName("ModelConfig");
111  wspace->import(modelConfig);
112  wspace->import(*data);
113  wspace->SetName("w");
114  wspace->writeToFile("rs101_ws.root");
115 
116 
117 
118  // First, let's use a Calculator based on the Profile Likelihood Ratio
119  //ProfileLikelihoodCalculator plc(*data, *modelWithConstraints, paramOfInterest);
120  ProfileLikelihoodCalculator plc(*data, modelConfig);
121  plc.SetTestSize(.05);
122  ConfInterval* lrint = plc.GetInterval(); // that was easy.
123 
124  // Let's make a plot
125  TCanvas* dataCanvas = new TCanvas("dataCanvas");
126  dataCanvas->Divide(2,1);
127 
128  dataCanvas->cd(1);
130  plotInt.SetTitle("Profile Likelihood Ratio and Posterior for S");
131  plotInt.Draw();
132 
133  // Second, use a Calculator based on the Feldman Cousins technique
134  FeldmanCousins fc(*data, modelConfig);
135  fc.UseAdaptiveSampling(true);
136  fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
137  fc.SetNBins(100); // number of points to test per parameter
138  fc.SetTestSize(.05);
139  // fc.SaveBeltToFile(true); // optional
140  ConfInterval* fcint = NULL;
141  fcint = fc.GetInterval(); // that was easy.
142 
143  RooFitResult* fit = modelWithConstraints->fitTo(*data, Save(true));
144 
145  // Third, use a Calculator based on Markov Chain monte carlo
146  // Before configuring the calculator, let's make a ProposalFunction
147  // that will achieve a high acceptance rate
148  ProposalHelper ph;
150  ph.SetCovMatrix(fit->covarianceMatrix());
152  ph.SetCacheSize(100);
153  ProposalFunction* pdfProp = ph.GetProposalFunction(); // that was easy
154 
155  MCMCCalculator mc(*data, modelConfig);
156  mc.SetNumIters(20000); // steps to propose in the chain
157  mc.SetTestSize(.05); // 95% CL
158  mc.SetNumBurnInSteps(40); // ignore first N steps in chain as "burn in"
159  mc.SetProposalFunction(*pdfProp);
160  mc.SetLeftSideTailFraction(0.5); // find a "central" interval
161  MCMCInterval* mcInt = (MCMCInterval*)mc.GetInterval(); // that was easy
162 
163 
164  // Get Lower and Upper limits from Profile Calculator
165  cout << "Profile lower limit on s = " << ((LikelihoodInterval*) lrint)->LowerLimit(*s) << endl;
166  cout << "Profile upper limit on s = " << ((LikelihoodInterval*) lrint)->UpperLimit(*s) << endl;
167 
168  // Get Lower and Upper limits from FeldmanCousins with profile construction
169  if (fcint != NULL) {
170  double fcul = ((PointSetInterval*) fcint)->UpperLimit(*s);
171  double fcll = ((PointSetInterval*) fcint)->LowerLimit(*s);
172  cout << "FC lower limit on s = " << fcll << endl;
173  cout << "FC upper limit on s = " << fcul << endl;
174  TLine* fcllLine = new TLine(fcll, 0, fcll, 1);
175  TLine* fculLine = new TLine(fcul, 0, fcul, 1);
176  fcllLine->SetLineColor(kRed);
177  fculLine->SetLineColor(kRed);
178  fcllLine->Draw("same");
179  fculLine->Draw("same");
180  dataCanvas->Update();
181  }
182 
183  // Plot MCMC interval and print some statistics
184  MCMCIntervalPlot mcPlot(*mcInt);
185  mcPlot.SetLineColor(kMagenta);
186  mcPlot.SetLineWidth(2);
187  mcPlot.Draw("same");
188 
189  double mcul = mcInt->UpperLimit(*s);
190  double mcll = mcInt->LowerLimit(*s);
191  cout << "MCMC lower limit on s = " << mcll << endl;
192  cout << "MCMC upper limit on s = " << mcul << endl;
193  cout << "MCMC Actual confidence level: "
194  << mcInt->GetActualConfidenceLevel() << endl;
195 
196  // 3-d plot of the parameter points
197  dataCanvas->cd(2);
198  // also plot the points in the markov chain
199  RooDataSet * chainData = mcInt->GetChainAsDataSet();
200 
201  assert(chainData);
202  std::cout << "plotting the chain data - nentries = " << chainData->numEntries() << std::endl;
203  TTree* chain = RooStats::GetAsTTree("chainTreeData","chainTreeData",*chainData);
204  assert(chain);
205  chain->SetMarkerStyle(6);
206  chain->SetMarkerColor(kRed);
207 
208  chain->Draw("s:ratioSigEff:ratioBkgEff","nll_MarkovChain_local_","box"); // 3-d box proporional to posterior
209 
210  // the points used in the profile construction
211  RooDataSet * parScanData = (RooDataSet*) fc.GetPointsToScan();
212  assert(parScanData);
213  std::cout << "plotting the scanned points used in the frequentist construction - npoints = " << parScanData->numEntries() << std::endl;
214  // getting the tree and drawing it -crashes (very strange....);
215  // TTree* parameterScan = RooStats::GetAsTTree("parScanTreeData","parScanTreeData",*parScanData);
216  // assert(parameterScan);
217  // parameterScan->Draw("s:ratioSigEff:ratioBkgEff","","candle goff");
218  TGraph2D *gr = new TGraph2D(parScanData->numEntries());
219  for (int ievt = 0; ievt < parScanData->numEntries(); ++ievt) {
220  const RooArgSet * evt = parScanData->get(ievt);
221  double x = evt->getRealValue("ratioBkgEff");
222  double y = evt->getRealValue("ratioSigEff");
223  double z = evt->getRealValue("s");
224  gr->SetPoint(ievt, x,y,z);
225  // std::cout << ievt << " " << x << " " << y << " " << z << std::endl;
226  }
227  gr->SetMarkerStyle(24);
228  gr->Draw("P SAME");
229 
230 
231  delete wspace;
232  delete lrint;
233  delete mcInt;
234  delete fcint;
235  delete data;
236 
237  /// print timing info
238  t.Stop();
239  t.Print();
240 }
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
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...
virtual void Draw(Option_t *option="")
Specific drawing options can be used to paint a TGraph2D:
Definition: TGraph2D.cxx:719
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.
Definition: Rtypes.h:61
virtual void SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval) ...
#define assert(cond)
Definition: unittest.h:542
virtual void SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval) ...
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
void rs101_limitexample()
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 ...
void Draw(const Option_t *options=NULL)
TTree * GetAsTTree(TString name, TString desc, const RooDataSet &data)
virtual void SetObservables(const RooArgSet &set)
specify the observables
Definition: ModelConfig.h:156
virtual LikelihoodInterval * GetInterval() const
Return a likelihood interval.
RooAbsData * GetPointsToScan()
TChain chain("h42")
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:839
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:75
Double_t x[n]
Definition: legend1.C:17
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
virtual void SetMarkerColor(Color_t mcolor=1)
Definition: TAttMarker.h:51
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:97
Float_t z[5]
Definition: Ifit.C:16
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:203
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 SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval) ...
virtual void SetNuisanceParameters(const RooArgSet &set)
specify the nuisance parameters (e.g. the rest of the parameters)
Definition: ModelConfig.h:130
void setConstant(Bool_t value=kTRUE)
void SetNBins(Int_t bins)
virtual void SetCacheSize(Int_t size)
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:293
void FluctuateNumDataEntries(bool flag=true)
void SetLineColor(Color_t color)
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
void SetLineWidth(Int_t width)
Bool_t writeToFile(const char *fileName, Bool_t recreate=kTRUE)
Save this current workspace into given file.
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
TGraphErrors * gr
Definition: legend1.C:25
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0)
Add a data point, with its coordinates specified in the 'data' argset, to the data set...
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
ConfInterval is an interface class for a generic interval in the RooStats framework.
Definition: ConfInterval.h:44
PointSetInterval is a concrete implementation of the ConfInterval interface.
virtual void SetProposalFunction(ProposalFunction &proposalFunction)
set the proposal function for suggesting new points for the MCMC
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Sets point number n.
Definition: TGraph2D.cxx:1707
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:360
Double_t y[n]
Definition: legend1.C:17
virtual PointSetInterval * GetInterval() const
Main interface to get a ConfInterval (will be a PointSetInterval)
virtual void SetCovMatrix(const TMatrixDSym &covMatrix)
RooFactoryWSTool & factory()
Return instance to factory tool.
RooCmdArg Save(Bool_t flag=kTRUE)
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.
void Print(Option_t *opts=0) const
Print contents of the workspace.
#define NULL
Definition: Rtypes.h:82
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 SetLeftSideTailFraction(Double_t a)
Set the left side tail fraction.
virtual void SetParametersOfInterest(const RooArgSet &set)
Definition: ModelConfig.h:115
A TTree object has a header with a name and a title.
Definition: TTree.h:98
virtual void SetVariables(RooArgList &vars)
virtual void SetGlobalObservables(const RooArgSet &set)
specify the global observables
Definition: ModelConfig.h:182
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:50
Bayesian Calculator estimating an interval or a credible region using the Markov-Chain Monte Carlo me...
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
tuple all
Definition: na49view.py:13
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event 'index'.
RooCmdArg Constrain(const RooArgSet &params)
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
void UseAdaptiveSampling(bool flag=true)
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Definition: RooArgSet.cxx:527
Stopwatch class.
Definition: TStopwatch.h:30