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