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