Logo ROOT   6.16/01
Reference Guide
rs102_hypotestwithshapes.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook -js
4/// rs102_hypotestwithshapes for RooStats project
5///
6/// This tutorial macro shows a typical search for a new particle
7/// by studying an invariant mass distribution.
8/// The macro creates a simple signal model and two background models,
9/// which are added to a RooWorkspace.
10/// The macro creates a toy dataset, and then uses a RooStats
11/// ProfileLikleihoodCalculator to do a hypothesis test of the
12/// background-only and signal+background hypotheses.
13/// In this example, shape uncertainties are not taken into account, but
14/// normalization uncertainties are.
15///
16/// \macro_image
17/// \macro_output
18/// \macro_code
19///
20/// \author Kyle Cranmer
21
22#include "RooDataSet.h"
23#include "RooRealVar.h"
24#include "RooGaussian.h"
25#include "RooAddPdf.h"
26#include "RooProdPdf.h"
27#include "RooAddition.h"
28#include "RooProduct.h"
29#include "TCanvas.h"
30#include "RooChebychev.h"
31#include "RooAbsPdf.h"
32#include "RooFit.h"
33#include "RooFitResult.h"
34#include "RooPlot.h"
35#include "RooAbsArg.h"
36#include "RooWorkspace.h"
39#include <string>
40
41
42// use this order for safety on library loading
43using namespace RooFit;
44using namespace RooStats;
45
46// see below for implementation
47void AddModel(RooWorkspace*);
48void AddData(RooWorkspace*);
49void DoHypothesisTest(RooWorkspace*);
50void MakePlots(RooWorkspace*);
51
52//____________________________________
53void rs102_hypotestwithshapes() {
54
55 // The main macro.
56
57 // Create a workspace to manage the project.
58 RooWorkspace* wspace = new RooWorkspace("myWS");
59
60 // add the signal and background models to the workspace
61 AddModel(wspace);
62
63 // add some toy data to the workspace
64 AddData(wspace);
65
66 // inspect the workspace if you wish
67 // wspace->Print();
68
69 // do the hypothesis test
70 DoHypothesisTest(wspace);
71
72 // make some plots
73 MakePlots(wspace);
74
75 // cleanup
76 delete wspace;
77}
78
79//____________________________________
80void AddModel(RooWorkspace* wks){
81
82 // Make models for signal (Higgs) and background (Z+jets and QCD)
83 // In real life, this part requires an intelligent modeling
84 // of signal and background -- this is only an example.
85
86 // set range of observable
87 Double_t lowRange = 60, highRange = 200;
88
89 // make a RooRealVar for the observable
90 RooRealVar invMass("invMass", "M_{inv}", lowRange, highRange,"GeV");
91
92
93 // --------------------------------------
94 // make a simple signal model.
95 RooRealVar mH("mH","Higgs Mass",130,90,160) ;
96 RooRealVar sigma1("sigma1","Width of Gaussian",12.,2,100) ;
97 RooGaussian sigModel("sigModel", "Signal Model", invMass, mH, sigma1);
98 // we will test this specific mass point for the signal
99 mH.setConstant();
100 // and we assume we know the mass resolution
101 sigma1.setConstant();
102
103 // --------------------------------------
104 // make zjj model. Just like signal model
105 RooRealVar mZ("mZ", "Z Mass", 91.2, 0, 100);
106 RooRealVar sigma1_z("sigma1_z","Width of Gaussian",10.,6,100) ;
107 RooGaussian zjjModel("zjjModel", "Z+jets Model", invMass, mZ, sigma1_z);
108 // we know Z mass
109 mZ.setConstant();
110 // assume we know resolution too
111 sigma1_z.setConstant();
112
113
114 // --------------------------------------
115 // make QCD model
116 RooRealVar a0("a0","a0",0.26,-1,1) ;
117 RooRealVar a1("a1","a1",-0.17596,-1,1) ;
118 RooRealVar a2("a2","a2",0.018437,-1,1) ;
119 RooRealVar a3("a3","a3",0.02,-1,1) ;
120 RooChebychev qcdModel("qcdModel","A Polynomial for QCD",invMass,RooArgList(a0,a1,a2)) ;
121
122 // let's assume this shape is known, but the normalization is not
123 a0.setConstant();
124 a1.setConstant();
125 a2.setConstant();
126
127 // --------------------------------------
128 // combined model
129
130 // Setting the fraction of Zjj to be 40% for initial guess.
131 RooRealVar fzjj("fzjj","fraction of zjj background events",.4,0.,1) ;
132
133 // Set the expected fraction of signal to 20%.
134 RooRealVar fsigExpected("fsigExpected","expected fraction of signal events",.2,0.,1) ;
135 fsigExpected.setConstant(); // use mu as main parameter, so fix this.
136
137 // Introduce mu: the signal strength in units of the expectation.
138 // eg. mu = 1 is the SM, mu = 0 is no signal, mu=2 is 2x the SM
139 RooRealVar mu("mu","signal strength in units of SM expectation",1,0.,2) ;
140
141 // Introduce ratio of signal efficiency to nominal signal efficiency.
142 // This is useful if you want to do limits on cross section.
143 RooRealVar ratioSigEff("ratioSigEff","ratio of signal efficiency to nominal signal efficiency",1. ,0.,2) ;
144 ratioSigEff.setConstant(kTRUE);
145
146 // finally the signal fraction is the product of the terms above.
147 RooProduct fsig("fsig","fraction of signal events",RooArgSet(mu,ratioSigEff,fsigExpected)) ;
148
149 // full model
150 RooAddPdf model("model","sig+zjj+qcd background shapes",RooArgList(sigModel,zjjModel, qcdModel),RooArgList(fsig,fzjj)) ;
151
152 // interesting for debugging and visualizing the model
153 // model.printCompactTree("","fullModel.txt");
154 // model.graphVizTree("fullModel.dot");
155
156 wks->import(model);
157}
158
159//____________________________________
160void AddData(RooWorkspace* wks){
161 // Add a toy dataset
162
163 Int_t nEvents = 150;
164 RooAbsPdf* model = wks->pdf("model");
165 RooRealVar* invMass = wks->var("invMass");
166
167 RooDataSet* data = model->generate(*invMass,nEvents);
168
169 wks->import(*data, Rename("data"));
170
171}
172
173//____________________________________
174void DoHypothesisTest(RooWorkspace* wks){
175
176 // Use a RooStats ProfileLikleihoodCalculator to do the hypothesis test.
178 model.SetWorkspace(*wks);
179 model.SetPdf("model");
180
181 //plc.SetData("data");
182
184 plc.SetData( *(wks->data("data") ));
185
186 // here we explicitly set the value of the parameters for the null.
187 // We want no signal contribution, eg. mu = 0
188 RooRealVar* mu = wks->var("mu");
189 // RooArgSet* nullParams = new RooArgSet("nullParams");
190 // nullParams->addClone(*mu);
191 RooArgSet poi(*mu);
192 RooArgSet * nullParams = (RooArgSet*) poi.snapshot();
193 nullParams->setRealValue("mu",0);
194
195
196 //plc.SetNullParameters(*nullParams);
197 plc.SetModel(model);
198 // NOTE: using snapshot will import nullparams
199 // in the WS and merge with existing "mu"
200 // model.SetSnapshot(*nullParams);
201
202 //use instead setNuisanceParameters
203 plc.SetNullParameters( *nullParams);
204
205
206
207 // We get a HypoTestResult out of the calculator, and we can query it.
208 HypoTestResult* htr = plc.GetHypoTest();
209 cout << "-------------------------------------------------" << endl;
210 cout << "The p-value for the null is " << htr->NullPValue() << endl;
211 cout << "Corresponding to a significance of " << htr->Significance() << endl;
212 cout << "-------------------------------------------------\n\n" << endl;
213
214
215}
216
217//____________________________________
218void MakePlots(RooWorkspace* wks) {
219
220 // Make plots of the data and the best fit model in two cases:
221 // first the signal+background case
222 // second the background-only case.
223
224 // get some things out of workspace
225 RooAbsPdf* model = wks->pdf("model");
226 RooAbsPdf* sigModel = wks->pdf("sigModel");
227 RooAbsPdf* zjjModel = wks->pdf("zjjModel");
228 RooAbsPdf* qcdModel = wks->pdf("qcdModel");
229
230 RooRealVar* mu = wks->var("mu");
231 RooRealVar* invMass = wks->var("invMass");
232 RooAbsData* data = wks->data("data");
233
234
235 // --------------------------------------
236 // Make plots for the Alternate hypothesis, eg. let mu float
237
238 mu->setConstant(kFALSE);
239
241
242 //plot sig candidates, full model, and individual components
243 new TCanvas();
244 RooPlot* frame = invMass->frame() ;
245 data->plotOn(frame ) ;
246 model->plotOn(frame) ;
247 model->plotOn(frame,Components(*sigModel),LineStyle(kDashed), LineColor(kRed)) ;
248 model->plotOn(frame,Components(*zjjModel),LineStyle(kDashed),LineColor(kBlack)) ;
249 model->plotOn(frame,Components(*qcdModel),LineStyle(kDashed),LineColor(kGreen)) ;
250
251 frame->SetTitle("An example fit to the signal + background model");
252 frame->Draw() ;
253 // cdata->SaveAs("alternateFit.gif");
254
255 // --------------------------------------
256 // Do Fit to the Null hypothesis. Eg. fix mu=0
257
258 mu->setVal(0); // set signal fraction to 0
259 mu->setConstant(kTRUE); // set constant
260
262
263 // plot signal candidates with background model and components
264 new TCanvas();
265 RooPlot* xframe2 = invMass->frame() ;
266 data->plotOn(xframe2, DataError(RooAbsData::SumW2)) ;
267 model->plotOn(xframe2) ;
268 model->plotOn(xframe2, Components(*zjjModel),LineStyle(kDashed),LineColor(kBlack)) ;
269 model->plotOn(xframe2, Components(*qcdModel),LineStyle(kDashed),LineColor(kGreen)) ;
270
271 xframe2->SetTitle("An example fit to the background-only model");
272 xframe2->Draw() ;
273 // cbkgonly->SaveAs("nullFit.gif");
274
275}
276
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kRed
Definition: Rtypes.h:63
@ kBlack
Definition: Rtypes.h:62
@ kGreen
Definition: Rtypes.h:63
@ kDashed
Definition: TAttLine.h:48
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooPlot * frame(const RooCmdArg &arg1, 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()) const
Create a new RooPlot on the heap with a drawing frame initialized for this object,...
void setConstant(Bool_t value=kTRUE)
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Bool_t setRealValue(const char *name, Double_t newVal=0, Bool_t verbose=kFALSE)
Set value of a RooAbsRealLValye stored in set with given name to newVal No error messages are printed...
Definition: RooArgSet.cxx:493
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition: RooPlot.cxx:1104
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooProduct a RooAbsReal implementation that represent the product of a given set of other RooAbsReal ...
Definition: RooProduct.h:32
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:204
virtual void SetNullParameters(const RooArgSet &set)
set parameter values for the null if using a common PDF
virtual void SetModel(const ModelConfig &model)
set the model (in this case can set only the parameters for the null hypothesis)
virtual void SetData(RooAbsData &data)
Set the DataSet, add to the the workspace if not already there.
HypoTestResult is a base class for results from hypothesis tests.
virtual Double_t Significance() const
familiar name for the Null p-value in terms of 1-sided Gaussian significance
virtual Double_t NullPValue() const
Return p-value for null hypothesis.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface class f...
virtual HypoTestResult * GetHypoTest() const
Return the hypothesis test result obtained from the likelihood ratio of the maximum likelihood value ...
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
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.
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
The Canvas class.
Definition: TCanvas.h:31
RooCmdArg Hesse(Bool_t flag=kTRUE)
RooCmdArg Rename(const char *suffix)
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg DataError(Int_t)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
RooCmdArg Minos(Bool_t flag=kTRUE)
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20