ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rs301_splot.C
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // SPlot tutorial
4 // author: Kyle Cranmer
5 // date Dec. 2008
6 //
7 // This tutorial shows an example of using SPlot to unfold two distributions.
8 // The physics context for the example is that we want to know
9 // the isolation distribution for real electrons from Z events
10 // and fake electrons from QCD. Isolation is our 'control' variable
11 // To unfold them, we need a model for an uncorrelated variable that
12 // discriminates between Z and QCD. To do this, we use the invariant
13 // mass of two electrons. We model the Z with a Gaussian and the QCD
14 // with a falling exponential.
15 //
16 // Note, since we don't have real data in this tutorial, we need to generate
17 // toy data. To do that we need a model for the isolation variable for
18 // both Z and QCD. This is only used to generate the toy data, and would
19 // not be needed if we had real data.
20 /////////////////////////////////////////////////////////////////////////
21 
22 #ifndef __CINT__
23 #include "RooGlobalFunc.h"
24 #endif
25 #include "RooRealVar.h"
26 #include "RooStats/SPlot.h"
27 #include "RooDataSet.h"
28 #include "RooRealVar.h"
29 #include "RooGaussian.h"
30 #include "RooExponential.h"
31 #include "RooChebychev.h"
32 #include "RooAddPdf.h"
33 #include "RooProdPdf.h"
34 #include "RooAddition.h"
35 #include "RooProduct.h"
36 #include "TCanvas.h"
37 #include "RooAbsPdf.h"
38 #include "RooFit.h"
39 #include "RooFitResult.h"
40 #include "RooWorkspace.h"
41 #include "RooConstVar.h"
42 
43 // use this order for safety on library loading
44 using namespace RooFit ;
45 using namespace RooStats ;
46 
47 
48 // see below for implementation
49 void AddModel(RooWorkspace*);
50 void AddData(RooWorkspace*);
51 void DoSPlot(RooWorkspace*);
52 void MakePlots(RooWorkspace*);
53 
55 {
56 
57  // Create a new workspace to manage the project.
58  RooWorkspace* wspace = new RooWorkspace("myWS");
59 
60  // add the signal and background models to the workspace.
61  // Inside this function you will find a discription our model.
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 sPlot.
71  //This wil make a new dataset with sWeights added for every event.
72  DoSPlot(wspace);
73 
74  // Make some plots showing the discriminating variable and
75  // the control variable after unfolding.
76  MakePlots(wspace);
77 
78  // cleanup
79  delete wspace;
80 
81 }
82 
83 
84 //____________________________________
86 
87  // Make models for signal (Higgs) and background (Z+jets and QCD)
88  // In real life, this part requires an intellegent modeling
89  // of signal and background -- this is only an example.
90 
91  // set range of observable
92  Double_t lowRange = 00, highRange = 200;
93 
94  // make a RooRealVar for the observables
95  RooRealVar invMass("invMass", "M_{inv}", lowRange, highRange,"GeV");
96  RooRealVar isolation("isolation", "isolation", 0., 20., "GeV");
97 
98 
99  /////////////////////////////////////////////
100  // make 2-d model for Z including the invariant mass
101  // distribution and an isolation distribution which we want to
102  // unfold from QCD.
103  std::cout << "make z model" << std::endl;
104  // mass model for Z
105  RooRealVar mZ("mZ", "Z Mass", 91.2, lowRange, highRange);
106  RooRealVar sigmaZ("sigmaZ", "Width of Gaussian", 2,0,10,"GeV");
107  RooGaussian mZModel("mZModel", "Z+jets Model", invMass, mZ, sigmaZ);
108  // we know Z mass
109  mZ.setConstant();
110  // we leave the width of the Z free during the fit in this example.
111 
112  // isolation model for Z. Only used to generate toy MC.
113  // the exponential is of the form exp(c*x). If we want
114  // the isolation to decay an e-fold every R GeV, we use
115  // c = -1/R.
116  RooConstVar zIsolDecayConst("zIsolDecayConst",
117  "z isolation decay constant", -1);
118  RooExponential zIsolationModel("zIsolationModel", "z isolation model",
119  isolation, zIsolDecayConst);
120 
121  // make the combined Z model
122  RooProdPdf zModel("zModel", "4-d model for Z",
123  RooArgSet(mZModel, zIsolationModel));
124 
125  //////////////////////////////////////////////
126  // make QCD model
127 
128  std::cout << "make qcd model" << std::endl;
129  // mass model for QCD.
130  // the exponential is of the form exp(c*x). If we want
131  // the mass to decay an e-fold every R GeV, we use
132  // c = -1/R.
133  // We can leave this parameter free during the fit.
134  RooRealVar qcdMassDecayConst("qcdMassDecayConst",
135  "Decay const for QCD mass spectrum",
136  -0.01, -100, 100,"1/GeV");
137  RooExponential qcdMassModel("qcdMassModel", "qcd Mass Model",
138  invMass, qcdMassDecayConst);
139 
140 
141  // isolation model for QCD. Only used to generate toy MC
142  // the exponential is of the form exp(c*x). If we want
143  // the isolation to decay an e-fold every R GeV, we use
144  // c = -1/R.
145  RooConstVar qcdIsolDecayConst("qcdIsolDecayConst",
146  "Et resolution constant", -.1);
147  RooExponential qcdIsolationModel("qcdIsolationModel", "QCD isolation model",
148  isolation, qcdIsolDecayConst);
149 
150  // make the 2-d model
151  RooProdPdf qcdModel("qcdModel", "2-d model for QCD",
152  RooArgSet(qcdMassModel, qcdIsolationModel));
153 
154  //////////////////////////////////////////////
155  // combined model
156 
157  // These variables represent the number of Z or QCD events
158  // They will be fitted.
159  RooRealVar zYield("zYield","fitted yield for Z",50 ,0.,1000) ;
160  RooRealVar qcdYield("qcdYield","fitted yield for QCD", 100 ,0.,1000) ;
161 
162  // now make the combined model
163  std::cout << "make full model" << std::endl;
164  RooAddPdf model("model","z+qcd background models",
165  RooArgList(zModel, qcdModel),
166  RooArgList(zYield,qcdYield));
167 
168 
169  // interesting for debugging and visualizing the model
170  model.graphVizTree("fullModel.dot");
171 
172  std::cout << "import model" << std::endl;
173 
174  ws->import(model);
175 }
176 
177 //____________________________________
179  // Add a toy dataset
180 
181  // how many events do we want?
182  Int_t nEvents = 1000;
183 
184  // get what we need out of the workspace to make toy data
185  RooAbsPdf* model = ws->pdf("model");
186  RooRealVar* invMass = ws->var("invMass");
187  RooRealVar* isolation = ws->var("isolation");
188 
189  // make the toy data
190  std::cout << "make data set and import to workspace" << std::endl;
191  RooDataSet* data = model->generate(RooArgSet(*invMass, *isolation),nEvents);
192 
193  // import data into workspace
194  ws->import(*data, Rename("data"));
195 
196 }
197 
198 //____________________________________
200  std::cout << "Calculate sWeights" << std::endl;
201 
202  // get what we need out of the workspace to do the fit
203  RooAbsPdf* model = ws->pdf("model");
204  RooRealVar* zYield = ws->var("zYield");
205  RooRealVar* qcdYield = ws->var("qcdYield");
206  RooDataSet* data = (RooDataSet*) ws->data("data");
207 
208  // fit the model to the data.
209  model->fitTo(*data, Extended() );
210 
211  // The sPlot technique requires that we fix the parameters
212  // of the model that are not yields after doing the fit.
213  RooRealVar* sigmaZ = ws->var("sigmaZ");
214  RooRealVar* qcdMassDecayConst = ws->var("qcdMassDecayConst");
215  sigmaZ->setConstant();
216  qcdMassDecayConst->setConstant();
217 
218 
220 
221 
222  // Now we use the SPlot class to add SWeights to our data set
223  // based on our model and our yield variables
224  RooStats::SPlot* sData = new RooStats::SPlot("sData","An SPlot",
225  *data, model, RooArgList(*zYield,*qcdYield) );
226 
227 
228  // Check that our weights have the desired properties
229 
230  std::cout << "Check SWeights:" << std::endl;
231 
232 
233  std::cout << std::endl << "Yield of Z is "
234  << zYield->getVal() << ". From sWeights it is "
235  << sData->GetYieldFromSWeight("zYield") << std::endl;
236 
237 
238  std::cout << "Yield of QCD is "
239  << qcdYield->getVal() << ". From sWeights it is "
240  << sData->GetYieldFromSWeight("qcdYield") << std::endl
241  << std::endl;
242 
243  for(Int_t i=0; i < 10; i++)
244  {
245  std::cout << "z Weight " << sData->GetSWeight(i,"zYield")
246  << " qcd Weight " << sData->GetSWeight(i,"qcdYield")
247  << " Total Weight " << sData->GetSumOfEventSWeight(i)
248  << std::endl;
249  }
250 
251  std::cout << std::endl;
252 
253  // import this new dataset with sWeights
254  std::cout << "import new dataset with sWeights" << std::endl;
255  ws->import(*data, Rename("dataWithSWeights"));
256 
257 
258 }
259 
261 
262  // Here we make plots of the discriminating variable (invMass) after the fit
263  // and of the control variable (isolation) after unfolding with sPlot.
264  std::cout << "make plots" << std::endl;
265 
266  // make our canvas
267  TCanvas* cdata = new TCanvas("sPlot","sPlot demo", 400, 600);
268  cdata->Divide(1,3);
269 
270  // get what we need out of the workspace
271  RooAbsPdf* model = ws->pdf("model");
272  RooAbsPdf* zModel = ws->pdf("zModel");
273  RooAbsPdf* qcdModel = ws->pdf("qcdModel");
274 
275  RooRealVar* isolation = ws->var("isolation");
276  RooRealVar* invMass = ws->var("invMass");
277 
278  // note, we get the dataset with sWeights
279  RooDataSet* data = (RooDataSet*) ws->data("dataWithSWeights");
280 
281  // this shouldn't be necessary, need to fix something with workspace
282  // do this to set parameters back to their fitted values.
283  model->fitTo(*data, Extended() );
284 
285  //plot invMass for data with full model and individual componenets overlayed
286  // TCanvas* cdata = new TCanvas();
287  cdata->cd(1);
288  RooPlot* frame = invMass->frame() ;
289  data->plotOn(frame ) ;
290  model->plotOn(frame) ;
291  model->plotOn(frame,Components(*zModel),LineStyle(kDashed), LineColor(kRed)) ;
292  model->plotOn(frame,Components(*qcdModel),LineStyle(kDashed),LineColor(kGreen)) ;
293 
294  frame->SetTitle("Fit of model to discriminating variable");
295  frame->Draw() ;
296 
297 
298  // Now use the sWeights to show isolation distribution for Z and QCD.
299  // The SPlot class can make this easier, but here we demonstrait in more
300  // detail how the sWeights are used. The SPlot class should make this
301  // very easy and needs some more development.
302 
303  // Plot isolation for Z component.
304  // Do this by plotting all events weighted by the sWeight for the Z component.
305  // The SPlot class adds a new variable that has the name of the corresponding
306  // yield + "_sw".
307  cdata->cd(2);
308 
309  // create weightfed data set
310  RooDataSet * dataw_z = new RooDataSet(data->GetName(),data->GetTitle(),data,*data->get(),0,"zYield_sw") ;
311 
312  RooPlot* frame2 = isolation->frame() ;
313  dataw_z->plotOn(frame2, DataError(RooAbsData::SumW2) ) ;
314 
315  frame2->SetTitle("isolation distribution for Z");
316  frame2->Draw() ;
317 
318  // Plot isolation for QCD component.
319  // Eg. plot all events weighted by the sWeight for the QCD component.
320  // The SPlot class adds a new variable that has the name of the corresponding
321  // yield + "_sw".
322  cdata->cd(3);
323  RooDataSet * dataw_qcd = new RooDataSet(data->GetName(),data->GetTitle(),data,*data->get(),0,"qcdYield_sw") ;
324  RooPlot* frame3 = isolation->frame() ;
325  dataw_qcd->plotOn(frame3,DataError(RooAbsData::SumW2) ) ;
326 
327  frame3->SetTitle("isolation distribution for QCD");
328  frame3->Draw() ;
329 
330  // cdata->SaveAs("SPlot.gif");
331 
332 }
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
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...
void AddData(RooWorkspace *)
Definition: rs301_splot.C:178
Definition: Rtypes.h:61
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
Double_t GetSumOfEventSWeight(Int_t numEvent) const
Definition: SPlot.cxx:216
void MakePlots(RooWorkspace *)
Definition: rs301_splot.C:260
void AddModel(RooWorkspace *)
Definition: rs301_splot.C:85
int Int_t
Definition: RtypesCore.h:41
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition: RooPlot.cxx:1099
Definition: Rtypes.h:61
static RooMsgService & instance()
Return reference to singleton instance.
RooCmdArg Extended(Bool_t flag=kTRUE)
Double_t GetYieldFromSWeight(const char *sVariable) const
Definition: SPlot.cxx:250
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
RooConstVar represent a constant real-valued object.
Definition: RooConstVar.h:25
Exponential p.d.f.
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
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
void graphVizTree(const char *fileName, const char *delimiter="\n", bool useTitle=false, bool useLatex=false)
Create a GraphViz .dot file visualizing the expression tree headed by this RooAbsArg object...
Definition: RooAbsArg.cxx:2119
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
const int nEvents
Definition: testRooFit.cxx:42
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)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
RooCmdArg Rename(const char *suffix)
void setSilentMode(Bool_t flag)
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
Double_t GetSWeight(Int_t numEvent, const char *sVariable) const
Definition: SPlot.cxx:170
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
double Double_t
Definition: RtypesCore.h:55
RooCmdArg Components(const RooArgSet &compSet)
void rs301_splot()
Definition: rs301_splot.C:54
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
This class calculates sWeights used to create an sPlot.
Definition: SPlot.h:58
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event 'index'.
void DoSPlot(RooWorkspace *)
Definition: rs301_splot.C:199
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