Logo ROOT   6.14/05
Reference Guide
rs101_limitexample.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roostats
3 /// \notebook
4 /// 'Limit Example' RooStats tutorial macro #101
5 /// This tutorial shows an example of creating a simple
6 /// model for a number counting experiment with uncertainty
7 /// on both the background rate and signal efficiency. We then
8 /// use a Confidence Interval Calculator to set a limit on the signal.
9 ///
10 /// \macro_image
11 /// \macro_output
12 /// \macro_code
13 ///
14 /// \author Kyle Cranmer
15 
16 #include "RooProfileLL.h"
17 #include "RooAbsPdf.h"
19 #include "RooRealVar.h"
20 #include "RooPlot.h"
21 #include "RooDataSet.h"
22 #include "RooTreeDataStore.h"
23 #include "TTree.h"
24 #include "TCanvas.h"
25 #include "TLine.h"
26 #include "TStopwatch.h"
27 
33 #include "RooStats/ConfInterval.h"
37 #include "RooStats/RooStatsUtils.h"
38 #include "RooStats/ModelConfig.h"
39 #include "RooStats/MCMCInterval.h"
43 #include "RooFitResult.h"
44 #include "TGraph2D.h"
45 
46 #include <cassert>
47 
48 // use this order for safety on library loading
49 using namespace RooFit;
50 using namespace RooStats;
51 
52 
53 void rs101_limitexample()
54 {
55  // --------------------------------------
56  // An example of setting a limit in a number counting experiment with uncertainty on background and signal
57 
58  // to time the macro
59  TStopwatch t;
60  t.Start();
61 
62  // --------------------------------------
63  // The Model building stage
64  // --------------------------------------
65  RooWorkspace* wspace = new RooWorkspace();
66  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
67  // wspace->factory("Gaussian::sigConstraint(ratioSigEff,1,0.05)"); // 5% signal efficiency uncertainty
68  // wspace->factory("Gaussian::bkgConstraint(ratioBkgEff,1,0.1)"); // 10% background efficiency uncertainty
69  wspace->factory("Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)"); // 5% signal efficiency uncertainty
70  wspace->factory("Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)"); // 10% background efficiency uncertainty
71  wspace->factory("PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)"); // product of terms
72  wspace->Print();
73 
74  RooAbsPdf* modelWithConstraints = wspace->pdf("modelWithConstraints"); // get the model
75  RooRealVar* obs = wspace->var("obs"); // get the observable
76  RooRealVar* s = wspace->var("s"); // get the signal we care about
77  RooRealVar* b = wspace->var("b"); // get the background and set it to a constant. Uncertainty included in ratioBkgEff
78  b->setConstant();
79 
80  RooRealVar* ratioSigEff = wspace->var("ratioSigEff"); // get uncertain parameter to constrain
81  RooRealVar* ratioBkgEff = wspace->var("ratioBkgEff"); // get uncertain parameter to constrain
82  RooArgSet constrainedParams(*ratioSigEff, *ratioBkgEff); // need to constrain these in the fit (should change default behavior)
83 
84  RooRealVar * gSigEff = wspace->var("gSigEff"); // global observables for signal efficiency
85  RooRealVar * gSigBkg = wspace->var("gSigBkg"); // global obs for background efficiency
86  gSigEff->setConstant();
87  gSigBkg->setConstant();
88 
89  // Create an example dataset with 160 observed events
90  obs->setVal(160.);
91  RooDataSet* data = new RooDataSet("exampleData", "exampleData", RooArgSet(*obs));
92  data->add(*obs);
93 
94  RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);
95 
96  // not necessary
97  modelWithConstraints->fitTo(*data, RooFit::Constrain(RooArgSet(*ratioSigEff, *ratioBkgEff)));
98 
99  // Now let's make some confidence intervals for s, our parameter of interest
100  RooArgSet paramOfInterest(*s);
101 
102  ModelConfig modelConfig(wspace);
103  modelConfig.SetPdf(*modelWithConstraints);
104  modelConfig.SetParametersOfInterest(paramOfInterest);
105  modelConfig.SetNuisanceParameters(constrainedParams);
106  modelConfig.SetObservables(*obs);
107  modelConfig.SetGlobalObservables( RooArgSet(*gSigEff,*gSigBkg));
108  modelConfig.SetName("ModelConfig");
109  wspace->import(modelConfig);
110  wspace->import(*data);
111  wspace->SetName("w");
112  wspace->writeToFile("rs101_ws.root");
113 
114 
115 
116  // First, let's use a Calculator based on the Profile Likelihood Ratio
117  //ProfileLikelihoodCalculator plc(*data, *modelWithConstraints, paramOfInterest);
118  ProfileLikelihoodCalculator plc(*data, modelConfig);
119  plc.SetTestSize(.05);
120  ConfInterval* lrinterval = plc.GetInterval(); // that was easy.
121 
122  // Let's make a plot
123  TCanvas* dataCanvas = new TCanvas("dataCanvas");
124  dataCanvas->Divide(2,1);
125 
126  dataCanvas->cd(1);
127  LikelihoodIntervalPlot plotInt((LikelihoodInterval*)lrinterval);
128  plotInt.SetTitle("Profile Likelihood Ratio and Posterior for S");
129  plotInt.Draw();
130 
131  // Second, use a Calculator based on the Feldman Cousins technique
132  FeldmanCousins fc(*data, modelConfig);
133  fc.UseAdaptiveSampling(true);
134  fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
135  fc.SetNBins(100); // number of points to test per parameter
136  fc.SetTestSize(.05);
137  // fc.SaveBeltToFile(true); // optional
138  ConfInterval* fcint = NULL;
139  fcint = fc.GetInterval(); // that was easy.
140 
141  RooFitResult* fit = modelWithConstraints->fitTo(*data, Save(true));
142 
143  // Third, use a Calculator based on Markov Chain monte carlo
144  // Before configuring the calculator, let's make a ProposalFunction
145  // that will achieve a high acceptance rate
146  ProposalHelper ph;
148  ph.SetCovMatrix(fit->covarianceMatrix());
150  ph.SetCacheSize(100);
151  ProposalFunction* pdfProp = ph.GetProposalFunction(); // that was easy
152 
153  MCMCCalculator mc(*data, modelConfig);
154  mc.SetNumIters(20000); // steps to propose in the chain
155  mc.SetTestSize(.05); // 95% CL
156  mc.SetNumBurnInSteps(40); // ignore first N steps in chain as "burn in"
157  mc.SetProposalFunction(*pdfProp);
158  mc.SetLeftSideTailFraction(0.5); // find a "central" interval
159  MCMCInterval* mcInt = (MCMCInterval*)mc.GetInterval(); // that was easy
160 
161 
162  // Get Lower and Upper limits from Profile Calculator
163  cout << "Profile lower limit on s = " << ((LikelihoodInterval*) lrinterval)->LowerLimit(*s) << endl;
164  cout << "Profile upper limit on s = " << ((LikelihoodInterval*) lrinterval)->UpperLimit(*s) << endl;
165 
166  // Get Lower and Upper limits from FeldmanCousins with profile construction
167  if (fcint != NULL) {
168  double fcul = ((PointSetInterval*) fcint)->UpperLimit(*s);
169  double fcll = ((PointSetInterval*) fcint)->LowerLimit(*s);
170  cout << "FC lower limit on s = " << fcll << endl;
171  cout << "FC upper limit on s = " << fcul << endl;
172  TLine* fcllLine = new TLine(fcll, 0, fcll, 1);
173  TLine* fculLine = new TLine(fcul, 0, fcul, 1);
174  fcllLine->SetLineColor(kRed);
175  fculLine->SetLineColor(kRed);
176  fcllLine->Draw("same");
177  fculLine->Draw("same");
178  dataCanvas->Update();
179  }
180 
181  // Plot MCMC interval and print some statistics
182  MCMCIntervalPlot mcPlot(*mcInt);
183  mcPlot.SetLineColor(kMagenta);
184  mcPlot.SetLineWidth(2);
185  mcPlot.Draw("same");
186 
187  double mcul = mcInt->UpperLimit(*s);
188  double mcll = mcInt->LowerLimit(*s);
189  cout << "MCMC lower limit on s = " << mcll << endl;
190  cout << "MCMC upper limit on s = " << mcul << endl;
191  cout << "MCMC Actual confidence level: "
192  << mcInt->GetActualConfidenceLevel() << endl;
193 
194  // 3-d plot of the parameter points
195  dataCanvas->cd(2);
196  // also plot the points in the markov chain
197  RooDataSet * chainData = mcInt->GetChainAsDataSet();
198 
199  assert(chainData);
200  std::cout << "plotting the chain data - nentries = " << chainData->numEntries() << std::endl;
201  TTree* chain = RooStats::GetAsTTree("chainTreeData","chainTreeData",*chainData);
202  assert(chain);
203  chain->SetMarkerStyle(6);
204  chain->SetMarkerColor(kRed);
205 
206  chain->Draw("s:ratioSigEff:ratioBkgEff","nll_MarkovChain_local_","box"); // 3-d box proportional to posterior
207 
208  // the points used in the profile construction
209  RooDataSet * parScanData = (RooDataSet*) fc.GetPointsToScan();
210  assert(parScanData);
211  std::cout << "plotting the scanned points used in the frequentist construction - npoints = " << parScanData->numEntries() << std::endl;
212  // getting the tree and drawing it -crashes (very strange....);
213  // TTree* parameterScan = RooStats::GetAsTTree("parScanTreeData","parScanTreeData",*parScanData);
214  // assert(parameterScan);
215  // parameterScan->Draw("s:ratioSigEff:ratioBkgEff","","goff");
216  TGraph2D *gr = new TGraph2D(parScanData->numEntries());
217  for (int ievt = 0; ievt < parScanData->numEntries(); ++ievt) {
218  const RooArgSet * evt = parScanData->get(ievt);
219  double x = evt->getRealValue("ratioBkgEff");
220  double y = evt->getRealValue("ratioSigEff");
221  double z = evt->getRealValue("s");
222  gr->SetPoint(ievt, x,y,z);
223  // std::cout << ievt << " " << x << " " << y << " " << z << std::endl;
224  }
225  gr->SetMarkerStyle(24);
226  gr->Draw("P SAME");
227 
228 
229  delete wspace;
230  delete lrinterval;
231  delete mcInt;
232  delete fcint;
233  delete data;
234 
235  // print timing info
236  t.Stop();
237  t.Print();
238 }
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
const RooArgList & floatParsFinal() const
Definition: RooFitResult.h:110
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
This class provides simple and straightforward utilities to plot a MCMCInterval object.
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
Definition: Rtypes.h:59
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:140
virtual ProposalFunction * GetProposalFunction()
virtual void SetUpdateProposalParameters(Bool_t updateParams)
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:688
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 Draw(Option_t *option="P0")
Specific drawing options can be used to paint a TGraph2D:
Definition: TGraph2D.cxx:704
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:195
TTree * GetAsTTree(TString name, TString desc, const RooDataSet &data)
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:3352
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
Double_t x[n]
Definition: legend1.C:17
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void setVal(Double_t value)
Set value of variable to &#39;value&#39;.
Definition: RooRealVar.cxx:204
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
void setConstant(Bool_t value=kTRUE)
virtual void SetCacheSize(Int_t size)
A simple line.
Definition: TLine.h:23
Bool_t writeToFile(const char *fileName, Bool_t recreate=kTRUE)
Save this current workspace into given file.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
TGraphErrors * gr
Definition: legend1.C:25
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event &#39;index&#39;.
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 &#39;data&#39; argset, to the data set...
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
The Canvas class.
Definition: TCanvas.h:31
ConfInterval is an interface class for a generic interval in the RooStats framework.
Definition: ConfInterval.h:35
Namespace for the RooStats classes.
Definition: Asimov.h:20
PointSetInterval is a concrete implementation of the ConfInterval interface.
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Sets point number n.
Definition: TGraph2D.cxx:1706
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
virtual Double_t LowerLimit(RooRealVar &param)
get the lowest value of param that is within the confidence interval
Double_t y[n]
Definition: legend1.C:17
static constexpr double s
virtual void SetCovMatrix(const TMatrixDSym &covMatrix)
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
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
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
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.
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:526
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
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:1079
MCMCInterval is a concrete implementation of the RooStats::ConfInterval interface.
Definition: MCMCInterval.h:30
virtual void SetVariables(RooArgList &vars)
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2248
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:40
Bayesian Calculator estimating an interval or a credible region using the Markov-Chain Monte Carlo me...
void Print(Option_t *opts=0) const
Print contents of the workspace.
RooCmdArg Constrain(const RooArgSet &params)
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:285
Stopwatch class.
Definition: TStopwatch.h:28