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