Logo ROOT   6.10/09
Reference Guide
StandardHistFactoryPlotsWithCategories.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roostats
3 /// \notebook -js
4 /// StandardHistFactoryPlotsWithCategories
5 ///
6 /// This is a standard demo that can be used with any ROOT file
7 /// prepared in the standard way. You specify:
8 /// - name for input ROOT file
9 /// - name of workspace inside ROOT file that holds model and data
10 /// - name of ModelConfig that specifies details for calculator tools
11 /// - name of dataset
12 ///
13 /// With default parameters the macro will attempt to run the
14 /// standard hist2workspace example and read the ROOT file
15 /// that it produces.
16 ///
17 /// The macro will scan through all the categories in a simPdf find the corresponding
18 /// observable. For each category, it will loop through each of the nuisance parameters
19 /// and plot
20 /// - the data
21 /// - the nominal model (blue)
22 /// - the +Nsigma (red)
23 /// - the -Nsigma (green)
24 ///
25 /// You can specify how many sigma to vary by changing nSigmaToVary.
26 /// You can also change the signal rate by changing muVal.
27 ///
28 /// The script produces a lot plots, you can merge them by doing:
29 /// ~~~{.cpp}
30 /// gs -q -dNOPAUSE -dBATCH -sDEVICE=pdfwrite -sOutputFile=merged.pdf `ls *pdf`
31 /// ~~~
32 ///
33 /// \macro_image
34 /// \macro_output
35 /// \macro_code
36 ///
37 /// \author Kyle Cranmer
38 
39 #include "TFile.h"
40 #include "TROOT.h"
41 #include "TCanvas.h"
42 #include "TList.h"
43 #include "TMath.h"
44 #include "TSystem.h"
45 #include "RooWorkspace.h"
46 #include "RooAbsData.h"
47 #include "RooRealVar.h"
48 #include "RooPlot.h"
49 #include "RooSimultaneous.h"
50 #include "RooCategory.h"
51 
52 #include "RooStats/ModelConfig.h"
54 
55 using namespace RooFit;
56 using namespace RooStats;
57 using namespace std;
58 
59 void StandardHistFactoryPlotsWithCategories(const char* infile = "",
60  const char* workspaceName = "combined",
61  const char* modelConfigName = "ModelConfig",
62  const char* dataName = "obsData"){
63 
64 
65  double nSigmaToVary=5.;
66  double muVal=0;
67  bool doFit=false;
68 
69  // -------------------------------------------------------
70  // First part is just to access a user-defined file
71  // or create the standard example file if it doesn't exist
72  const char* filename = "";
73  if (!strcmp(infile,"")) {
74  filename = "results/example_combined_GaussExample_model.root";
75  bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
76  // if file does not exists generate with histfactory
77  if (!fileExist) {
78 #ifdef _WIN32
79  cout << "HistFactory file cannot be generated on Windows - exit" << endl;
80  return;
81 #endif
82  // Normally this would be run on the command line
83  cout <<"will run standard hist2workspace example"<<endl;
84  gROOT->ProcessLine(".! prepareHistFactory .");
85  gROOT->ProcessLine(".! hist2workspace config/example.xml");
86  cout <<"\n\n---------------------"<<endl;
87  cout <<"Done creating example input"<<endl;
88  cout <<"---------------------\n\n"<<endl;
89  }
90 
91  }
92  else
93  filename = infile;
94 
95  // Try to open the file
96  TFile *file = TFile::Open(filename);
97 
98  // if input file was specified byt not found, quit
99  if(!file ){
100  cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
101  return;
102  }
103 
104  // -------------------------------------------------------
105  // Tutorial starts here
106  // -------------------------------------------------------
107 
108  // get the workspace out of the file
109  RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
110  if(!w){
111  cout <<"workspace not found" << endl;
112  return;
113  }
114 
115  // get the modelConfig out of the file
116  ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
117 
118  // get the modelConfig out of the file
119  RooAbsData* data = w->data(dataName);
120 
121  // make sure ingredients are found
122  if(!data || !mc){
123  w->Print();
124  cout << "data or ModelConfig was not found" <<endl;
125  return;
126  }
127 
128  // -------------------------------------------------------
129  // now use the profile inspector
130 
131  RooRealVar* obs = (RooRealVar*)mc->GetObservables()->first();
132  TList* list = new TList();
133 
134 
135  RooRealVar * firstPOI = dynamic_cast<RooRealVar*>(mc->GetParametersOfInterest()->first());
136 
137  firstPOI->setVal(muVal);
138  // firstPOI->setConstant();
139  if(doFit){
140  mc->GetPdf()->fitTo(*data);
141  }
142 
143  // -------------------------------------------------------
144 
145 
146  mc->GetNuisanceParameters()->Print("v");
147  int nPlotsMax = 1000;
148  cout <<" check expectedData by category"<<endl;
149  RooDataSet* simData=NULL;
150  RooSimultaneous* simPdf = NULL;
151  if(strcmp(mc->GetPdf()->ClassName(),"RooSimultaneous")==0){
152  cout <<"Is a simultaneous PDF"<<endl;
153  simPdf = (RooSimultaneous *)(mc->GetPdf());
154  } else {
155  cout <<"Is not a simultaneous PDF"<<endl;
156  }
157 
158 
159 
160  if(doFit) {
161  RooCategory* channelCat = (RooCategory*) (&simPdf->indexCat());
162  TIterator* iter = channelCat->typeIterator() ;
163  RooCatType* tt = NULL;
164  tt=(RooCatType*) iter->Next();
165  RooAbsPdf* pdftmp = ((RooSimultaneous*)mc->GetPdf())->getPdf(tt->GetName()) ;
166  RooArgSet* obstmp = pdftmp->getObservables(*mc->GetObservables()) ;
167  obs = ((RooRealVar*)obstmp->first());
168  RooPlot* frame = obs->frame();
169  cout <<Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())<<endl;
170  cout << tt->GetName() << " " << channelCat->getLabel() <<endl;
171  data->plotOn(frame,MarkerSize(1),Cut(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())),DataError(RooAbsData::None));
172 
173  Double_t normCount = data->sumEntries(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())) ;
174 
175  pdftmp->plotOn(frame,LineWidth(2.),Normalization(normCount,RooAbsReal::NumEvent)) ;
176  frame->Draw();
177  cout <<"expected events = " << mc->GetPdf()->expectedEvents(*data->get()) <<endl;
178  return;
179  }
180 
181 
182 
183  int nPlots=0;
184  if(!simPdf){
185 
187  RooRealVar* var = NULL;
188  while( (var = (RooRealVar*) it->Next()) != NULL){
189  RooPlot* frame = obs->frame();
190  frame->SetYTitle(var->GetName());
191  data->plotOn(frame,MarkerSize(1));
192  var->setVal(0);
193  mc->GetPdf()->plotOn(frame,LineWidth(1.));
194  var->setVal(1);
196  var->setVal(-1);
198  list->Add(frame);
199  var->setVal(0);
200  }
201 
202 
203  } else {
204  RooCategory* channelCat = (RooCategory*) (&simPdf->indexCat());
205  // TIterator* iter = simPdf->indexCat().typeIterator() ;
206  TIterator* iter = channelCat->typeIterator() ;
207  RooCatType* tt = NULL;
208  while(nPlots<nPlotsMax && (tt=(RooCatType*) iter->Next())) {
209 
210  cout << "on type " << tt->GetName() << " " << endl;
211  // Get pdf associated with state from simpdf
212  RooAbsPdf* pdftmp = simPdf->getPdf(tt->GetName()) ;
213 
214  // Generate observables defined by the pdf associated with this state
215  RooArgSet* obstmp = pdftmp->getObservables(*mc->GetObservables()) ;
216  // obstmp->Print();
217 
218 
219  obs = ((RooRealVar*)obstmp->first());
220 
222  RooRealVar* var = NULL;
223  while(nPlots<nPlotsMax && (var = (RooRealVar*) it->Next())){
224  TCanvas* c2 = new TCanvas("c2");
225  RooPlot* frame = obs->frame();
226  frame->SetName(Form("frame%d",nPlots));
227  frame->SetYTitle(var->GetName());
228 
229  cout <<Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())<<endl;
230  cout << tt->GetName() << " " << channelCat->getLabel() <<endl;
231  data->plotOn(frame,MarkerSize(1),Cut(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())),DataError(RooAbsData::None));
232 
233  Double_t normCount = data->sumEntries(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())) ;
234 
235  if(strcmp(var->GetName(),"Lumi")==0){
236  cout <<"working on lumi"<<endl;
237  var->setVal(w->var("nominalLumi")->getVal());
238  var->Print();
239  } else{
240  var->setVal(0);
241  }
242  // w->allVars().Print("v");
243  // mc->GetNuisanceParameters()->Print("v");
244  // pdftmp->plotOn(frame,LineWidth(2.));
245  // mc->GetPdf()->plotOn(frame,LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data));
246  //pdftmp->plotOn(frame,LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data));
247  normCount = pdftmp->expectedEvents(*obs);
248  pdftmp->plotOn(frame,LineWidth(2.),Normalization(normCount,RooAbsReal::NumEvent)) ;
249 
250  if(strcmp(var->GetName(),"Lumi")==0){
251  cout <<"working on lumi"<<endl;
252  var->setVal(w->var("nominalLumi")->getVal()+0.05);
253  var->Print();
254  } else{
255  var->setVal(nSigmaToVary);
256  }
257  // pdftmp->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2));
258  // mc->GetPdf()->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data));
259  //pdftmp->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data));
260  normCount = pdftmp->expectedEvents(*obs);
262 
263  if(strcmp(var->GetName(),"Lumi")==0){
264  cout <<"working on lumi"<<endl;
265  var->setVal(w->var("nominalLumi")->getVal()-0.05);
266  var->Print();
267  } else{
268  var->setVal(-nSigmaToVary);
269  }
270  // pdftmp->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2));
271  // mc->GetPdf()->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(*channelCat,tt->GetName()),ProjWData(*data));
272  //pdftmp->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(*channelCat,tt->GetName()),ProjWData(*data));
273  normCount = pdftmp->expectedEvents(*obs);
275 
276 
277 
278  // set them back to normal
279  if(strcmp(var->GetName(),"Lumi")==0){
280  cout <<"working on lumi"<<endl;
281  var->setVal(w->var("nominalLumi")->getVal());
282  var->Print();
283  } else{
284  var->setVal(0);
285  }
286 
287  list->Add(frame);
288 
289  // quit making plots
290  ++nPlots;
291 
292  frame->Draw();
293  c2->SaveAs(Form("%s_%s_%s.pdf",tt->GetName(),obs->GetName(),var->GetName()));
294  delete c2;
295  }
296  }
297  }
298 
299 
300 
301  // -------------------------------------------------------
302 
303 
304  // now make plots
305  TCanvas* c1 = new TCanvas("c1","ProfileInspectorDemo",800,200);
306  if(list->GetSize()>4){
307  double n = list->GetSize();
308  int nx = (int)sqrt(n) ;
309  int ny = TMath::CeilNint(n/nx);
310  nx = TMath::CeilNint( sqrt(n) );
311  c1->Divide(ny,nx);
312  } else
313  c1->Divide(list->GetSize());
314  for(int i=0; i<list->GetSize(); ++i){
315  c1->cd(i+1);
316  list->At(i)->Draw();
317  }
318 
319 
320 
321 
322 
323 }
const int nx
Definition: kalman.C:16
virtual Double_t sumEntries() const =0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1272
TIterator * createIterator(Bool_t dir=kIterForward) const
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
Definition: ModelConfig.h:237
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:552
void doFit(int n, const char *fitter)
RooCmdArg Cut(const char *cutSpec)
void SetName(const char *name)
Set the name of the RooPlot to &#39;name&#39;.
Definition: RooPlot.cxx:1077
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
RooCmdArg LineColor(Color_t color)
virtual const RooArgSet * get() const
Definition: RooAbsData.h:77
return c1
Definition: legend1.C:41
Definition: Rtypes.h:56
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
void SetYTitle(const char *title)
Definition: RooPlot.cxx:1155
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:679
#define gROOT
Definition: TROOT.h:375
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:202
Definition: Rtypes.h:56
STL namespace.
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
#define NULL
Definition: RtypesCore.h:88
Iterator abstract base class.
Definition: TIterator.h:30
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3909
double sqrt(double)
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:135
TText * tt
Definition: textangle.C:16
RooCmdArg MarkerSize(Size_t size)
const int ny
Definition: kalman.C:17
RooCmdArg DataError(Int_t)
RooCmdArg LineStyle(Style_t style)
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
RooCatType is an auxilary class for RooAbsCategory and defines a a single category state...
Definition: RooCatType.h:22
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
Definition: RooAbsPdf.cxx:2930
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
A doubly linked list.
Definition: TList.h:43
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
R__EXTERN TSystem * gSystem
Definition: TSystem.h:539
TIterator * typeIterator() const
Return iterator over all defined states.
virtual const Text_t * GetName() const
Returns name of object.
Definition: RooCatType.h:44
const RooAbsCategoryLValue & indexCat() const
RooAbsArg * first() const
char * Form(const char *fmt,...)
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:315
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
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
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name) ...
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
The Canvas class.
Definition: TCanvas.h:31
Namespace for the RooStats classes.
Definition: Asimov.h:20
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:222
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:225
return c2
Definition: legend2.C:14
double Double_t
Definition: RtypesCore.h:55
RooCmdArg Normalization(Double_t scaleFactor)
RooCmdArg LineWidth(Width_t width)
virtual const char * getLabel() const
Return label string of current state.
Definition: RooCategory.h:39
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
RooAbsPdf * getPdf(const char *catName) const
Return the p.d.f associated with the given index category name.
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:1135
virtual void Add(TObject *obj)
Definition: TList.h:77
Definition: file.py:1
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
Definition: ModelConfig.h:228
virtual TObject * Next()=0
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
virtual void SaveAs(const char *filename="", Option_t *option="") const
Save Pad contents in a file in one of various formats.
Definition: TPad.cxx:5468
virtual Int_t GetSize() const
Definition: TCollection.h:89
void Print(Option_t *opts=0) const
Print contents of the workspace.
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset...
const Int_t n
Definition: legend1.C:16
Int_t CeilNint(Double_t x)
Definition: TMath.h:597
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