Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
54
55using namespace RooFit;
56using namespace RooStats;
57using std::cout, std::endl;
58
59void StandardHistFactoryPlotsWithCategories(const char *infile = "", const char *workspaceName = "combined",
60 const char *modelConfigName = "ModelConfig",
61 const char *dataName = "obsData")
62{
63
64 double nSigmaToVary = 5.;
65 double muVal = 0;
66 bool doFit = false;
67
68 // -------------------------------------------------------
69 // First part is just to access a user-defined file
70 // or create the standard example file if it doesn't exist
71 const char *filename = "";
72 if (!strcmp(infile, "")) {
73 filename = "results/example_combined_GaussExample_model.root";
74 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
75 // if file does not exists generate with histfactory
76 if (!fileExist) {
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 } else
87 filename = infile;
88
89 // Try to open the file
90 TFile *file = TFile::Open(filename);
91
92 // if input file was specified but not found, quit
93 if (!file) {
94 cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
95 return;
96 }
97
98 // -------------------------------------------------------
99 // Tutorial starts here
100 // -------------------------------------------------------
101
102 // get the workspace out of the file
103 RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
104 if (!w) {
105 cout << "workspace not found" << endl;
106 return;
107 }
108
109 // get the modelConfig out of the file
110 ModelConfig *mc = (ModelConfig *)w->obj(modelConfigName);
111
112 // get the modelConfig out of the file
113 RooAbsData *data = w->data(dataName);
114
115 // make sure ingredients are found
116 if (!data || !mc) {
117 w->Print();
118 cout << "data or ModelConfig was not found" << endl;
119 return;
120 }
121
122 // -------------------------------------------------------
123 // now use the profile inspector
124
125 RooRealVar *obs = (RooRealVar *)mc->GetObservables()->first();
126 std::vector<RooPlot *> frameList;
127
128 RooRealVar *firstPOI = dynamic_cast<RooRealVar *>(mc->GetParametersOfInterest()->first());
129
130 firstPOI->setVal(muVal);
131 // firstPOI->setConstant();
132 if (doFit) {
133 mc->GetPdf()->fitTo(*data);
134 }
135
136 // -------------------------------------------------------
137
138 mc->GetNuisanceParameters()->Print("v");
139 int nPlotsMax = 1000;
140 cout << " check expectedData by category" << endl;
141 RooDataSet *simData = NULL;
142 RooSimultaneous *simPdf = NULL;
143 if (strcmp(mc->GetPdf()->ClassName(), "RooSimultaneous") == 0) {
144 cout << "Is a simultaneous PDF" << endl;
145 simPdf = (RooSimultaneous *)(mc->GetPdf());
146 } else {
147 cout << "Is not a simultaneous PDF" << endl;
148 }
149
150 if (doFit) {
151 RooCategory *channelCat = (RooCategory *)(&simPdf->indexCat());
152 auto const& catName = channelCat->begin()->first;
153 RooAbsPdf *pdftmp = ((RooSimultaneous *)mc->GetPdf())->getPdf(catName.c_str());
154 std::unique_ptr<RooArgSet> obstmp{pdftmp->getObservables(*mc->GetObservables())};
155 obs = ((RooRealVar *)obstmp->first());
156 RooPlot *frame = obs->frame();
157 cout << Form("%s==%s::%s", channelCat->GetName(), channelCat->GetName(), catName.c_str()) << endl;
158 cout << catName << " " << channelCat->getLabel() << endl;
159 data->plotOn(frame, MarkerSize(1),
160 Cut(Form("%s==%s::%s", channelCat->GetName(), channelCat->GetName(), catName.c_str())),
162
163 Double_t normCount =
164 data->sumEntries(Form("%s==%s::%s", channelCat->GetName(), channelCat->GetName(), catName.c_str()));
165
166 pdftmp->plotOn(frame, LineWidth(2.), Normalization(normCount, RooAbsReal::NumEvent));
167 frame->Draw();
168 cout << "expected events = " << mc->GetPdf()->expectedEvents(*data->get()) << endl;
169 return;
170 }
171
172 int nPlots = 0;
173 if (!simPdf) {
174
175 for (auto *var : static_range_cast<RooRealVar *>(*mc->GetNuisanceParameters())) {
176 RooPlot *frame = obs->frame();
177 frame->SetYTitle(var->GetName());
178 data->plotOn(frame, MarkerSize(1));
179 var->setVal(0);
180 mc->GetPdf()->plotOn(frame, LineWidth(1.));
181 var->setVal(1);
182 mc->GetPdf()->plotOn(frame, LineColor(kRed), LineStyle(kDashed), LineWidth(1));
183 var->setVal(-1);
185 frameList.push_back(frame);
186 var->setVal(0);
187 }
188
189 } else {
190 RooCategory *channelCat = (RooCategory *)(&simPdf->indexCat());
191 for (auto const& tt : *channelCat) {
192
193 if (nPlots == nPlotsMax) {
194 break;
195 }
196
197 auto const& catName = tt.first;
198
199 cout << "on type " << catName << " " << endl;
200 // Get pdf associated with state from simpdf
201 RooAbsPdf *pdftmp = simPdf->getPdf(catName.c_str());
202
203 // Generate observables defined by the pdf associated with this state
204 std::unique_ptr<RooArgSet> obstmp{pdftmp->getObservables(*mc->GetObservables())};
205 // obstmp->Print();
206
207 obs = ((RooRealVar *)obstmp->first());
208
209 for (auto *var : static_range_cast<RooRealVar*>(*mc->GetNuisanceParameters())) {
210 if (nPlots == nPlotsMax) break;
211
212 TCanvas *c2 = new TCanvas("c2");
213 RooPlot *frame = obs->frame();
214 frame->SetName(Form("frame%d", nPlots));
215 frame->SetYTitle(var->GetName());
216
217 cout << Form("%s==%s::%s", channelCat->GetName(), channelCat->GetName(), catName.c_str()) << endl;
218 cout << catName << " " << channelCat->getLabel() << endl;
219 data->plotOn(frame, MarkerSize(1),
220 Cut(Form("%s==%s::%s", channelCat->GetName(), channelCat->GetName(), catName.c_str())),
222
223 Double_t normCount =
224 data->sumEntries(Form("%s==%s::%s", channelCat->GetName(), channelCat->GetName(), catName.c_str()));
225
226 if (strcmp(var->GetName(), "Lumi") == 0) {
227 cout << "working on lumi" << endl;
228 var->setVal(w->var("nominalLumi")->getVal());
229 var->Print();
230 } else {
231 var->setVal(0);
232 }
233 // w->allVars().Print("v");
234 // mc->GetNuisanceParameters()->Print("v");
235 // pdftmp->plotOn(frame,LineWidth(2.));
236 // mc->GetPdf()->plotOn(frame,LineWidth(2.),Slice(*channelCat,catName.c_str()),ProjWData(*data));
237 // pdftmp->plotOn(frame,LineWidth(2.),Slice(*channelCat,catName.c_str()),ProjWData(*data));
238 normCount = pdftmp->expectedEvents(*obs);
239 pdftmp->plotOn(frame, LineWidth(2.), Normalization(normCount, RooAbsReal::NumEvent));
240
241 if (strcmp(var->GetName(), "Lumi") == 0) {
242 cout << "working on lumi" << endl;
243 var->setVal(w->var("nominalLumi")->getVal() + 0.05);
244 var->Print();
245 } else {
246 var->setVal(nSigmaToVary);
247 }
248 // pdftmp->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2));
249 // mc->GetPdf()->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(*channelCat,catName.c_str()),ProjWData(*data));
250 // pdftmp->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(*channelCat,catName.c_str()),ProjWData(*data));
251 normCount = pdftmp->expectedEvents(*obs);
252 pdftmp->plotOn(frame, LineWidth(2.), LineColor(kRed), LineStyle(kDashed),
254
255 if (strcmp(var->GetName(), "Lumi") == 0) {
256 cout << "working on lumi" << endl;
257 var->setVal(w->var("nominalLumi")->getVal() - 0.05);
258 var->Print();
259 } else {
260 var->setVal(-nSigmaToVary);
261 }
262 // pdftmp->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2));
263 // mc->GetPdf()->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(*channelCat,catName.c_str()),ProjWData(*data));
264 // pdftmp->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(*channelCat,catName.c_str()),ProjWData(*data));
265 normCount = pdftmp->expectedEvents(*obs);
266 pdftmp->plotOn(frame, LineWidth(2.), LineColor(kGreen), LineStyle(kDashed),
268
269 // set them back to normal
270 if (strcmp(var->GetName(), "Lumi") == 0) {
271 cout << "working on lumi" << endl;
272 var->setVal(w->var("nominalLumi")->getVal());
273 var->Print();
274 } else {
275 var->setVal(0);
276 }
277
278 frameList.push_back(frame);
279
280 // quit making plots
281 ++nPlots;
282
283 frame->Draw();
284 c2->SaveAs(Form("%s_%s_%s.pdf", catName.c_str(), obs->GetName(), var->GetName()));
285 delete c2;
286 }
287 }
288 }
289
290 // -------------------------------------------------------
291
292 // now make plots
293 TCanvas *c1 = new TCanvas("c1", "ProfileInspectorDemo", 800, 200);
294 int nFrames = frameList.size();
295 if (nFrames > 4) {
296 int nx = (int)sqrt(nFrames);
297 int ny = TMath::CeilNint(nFrames / nx);
298 nx = TMath::CeilNint(sqrt(nFrames));
299 c1->Divide(ny, nx);
300 } else
301 c1->Divide(nFrames);
302 for (int i = 0; i < nFrames; ++i) {
303 c1->cd(i + 1);
304 frameList[i]->Draw();
305 }
306}
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
#define gROOT
Definition TROOT.h:406
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
R__EXTERN TSystem * gSystem
Definition TSystem.h:561
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
const char * getLabel() const
Retrieve current label. Use getCurrentLabel() for more clarity.
std::map< std::string, value_type >::const_iterator begin() const
Iterator for category state names. Points to pairs of index and name.
RooAbsArg * first() const
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition RooAbsPdf.h:124
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:157
RooPlot * frame(const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create a new RooPlot on the heap with a drawing frame initialized for this object,...
Object to represent discrete states.
Definition RooCategory.h:28
Container class to hold unbinned data.
Definition RooDataSet.h:34
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
void SetName(const char *name) override
Set the name of the RooPlot to 'name'.
Definition RooPlot.cxx:1181
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Definition RooPlot.h:132
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:597
void SetYTitle(const char *title)
Definition RooPlot.cxx:1296
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
RooAbsPdf * getPdf(RooStringView catName) const
Return the p.d.f associated with the given index category name.
const RooAbsCategoryLValue & indexCat() const
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
Persistable container for RooFit projects.
The Canvas class.
Definition TCanvas.h:23
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4086
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
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:1296
RooCmdArg DataError(Int_t)
RooCmdArg LineWidth(Width_t width)
RooCmdArg Normalization(double scaleFactor)
RooCmdArg MarkerSize(Size_t size)
RooCmdArg Cut(const char *cutSpec)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:678
auto * tt
Definition textangle.C:16