Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
StandardHistFactoryPlotsWithCategories.py
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#
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 (C++ version), and P. P. (Python translation)
38
39
40import ROOT
41
42
43def StandardHistFactoryPlotsWithCategories(
44 infile="", workspaceName="combined", modelConfigName="ModelConfig", dataName="obsData"
45):
46
47 nSigmaToVary = 5.0
48 muVal = 0
49 doFit = False
50
51 # -------------------------------------------------------
52 # First part is just to access a user-defined file
53 # or create the standard example file if it doesn't exist
54 filename = ""
55 if infile == "":
56 filename = "results/example_combined_GaussExample_model.root"
57 fileExist = not ROOT.gSystem.AccessPathName(filename) # note opposite return code
58 # if file does not exists generate with histfactory
59 if not fileExist:
60 # Normally this would be run on the command line
61 print(f"will run standard hist2workspace example")
62 ROOT.gROOT.ProcessLine(".! prepareHistFactory .")
63 ROOT.gROOT.ProcessLine(".! hist2workspace config/example.xml")
64 print(f"\n\n---------------------")
65 print(f"Done creating example input")
66 print(f"---------------------\n\n")
67
68 else:
69 filename = infile
70
71 # Try to open the file
72 file = ROOT.TFile.Open(filename)
73
74 # if input file was specified but not found, quit
75 if not file:
76 print(f"StandardRooStatsDemoMacro: Input file {filename} is not found")
77 return
78
79 # -------------------------------------------------------
80 # Tutorial starts here
81 # -------------------------------------------------------
82
83 # get the workspace out of the file
84 w = file.Get(workspaceName)
85 if not w:
86 print(f"workspace not found")
87 return
88
89 # get the modelConfig out of the file
90 mc = w.obj(modelConfigName)
91
92 # get the modelConfig out of the file
93 data = w.data(dataName)
94
95 # make sure ingredients are found
96 if not data or not mc:
97 w.Print()
98 print(f"data or ModelConfig was not found")
99 return
100
101 # -------------------------------------------------------
102 # now use the profile inspector
103
104 obs = mc.GetObservables().first()
105 frameList = []
106
107 firstPOI = mc.GetParametersOfInterest().first()
108
109 firstPOI.setVal(muVal)
110 # firstPOI.setConstant()
111 if doFit:
112 mc.GetPdf().fitTo(data)
113
114 # -------------------------------------------------------
115
117 nPlotsMax = 1000
118 print(f" check expectedData by category")
119 simData = ROOT.kNone
120 simPdf = ROOT.kNone
121 if str(mc.GetPdf().ClassName()) == "RooSimultaneous":
122 print(f"Is a simultaneous PDF")
123 simPdf = mc.GetPdf()
124 else:
125 print(f"Is not a simultaneous PDF")
126
127 if doFit:
128 channelCat = simPdf.indexCat()
129 catName = channelCat.begin().first
130 pdftmp = (mc.GetPdf()).getPdf(str(catName))
132 obs = obstmp.first()
133 frame = obs.frame()
134 print("{0}=={0}::{1}".format(channelCat.GetName(), catName))
135 print(catName, " ", channelCat.getLabel())
137 frame,
138 MarkerSize=1,
139 Cut="{0}=={0}::{1}".format(channelCat.GetName(), catName),
140 DataError="None",
141 )
142
143 normCount = data.sumEntries("{0}=={0}::{1}".format(channelCat.GetName(), catName))
144
145 pdftmp.plotOn(frame, Normalization(normCount, ROOT.RooAbsReal.NumEvent), LineWidth=2)
146 frame.Draw()
147 print("expected events = ", mc.GetPdf().expectedEvents(data.get()))
148 return
149
150 nPlots = 0
151 if not simPdf:
152
153 for var in mc.GetNuisanceParameters():
154 frame = obs.frame()
156 data.plotOn(frame, MarkerSize(1))
157 value = var.getVal()
158 mc.GetPdf().plotOn(frame, LineWidth(1))
159 var.setVal(value + var.getError());
160 mc.GetPdf().plotOn(frame, LineColor(kRed), LineStyle(kDashed), LineWidth(1))
161 var.setVal(value - var.getError());
162 mc.GetPdf().plotOn(frame, LineColor(kGreen), LineStyle(kDashed), LineWidth(1))
163 frameList.append(frame)
164 var.setVal(value)
165
166 else:
167 channelCat = simPdf.indexCat()
168 for tt in channelCat:
169
170 if nPlots == nPlotsMax:
171 break
172
173 catName = tt.first
174
175 print("on type ", catName, " ")
176 # Get pdf associated with state from simpdf
177 pdftmp = simPdf.getPdf(str(catName))
178
179 # Generate observables defined by the pdf associated with this state
181 # obstmp.Print()
182
183 obs = obstmp.first()
184
185 for var in mc.GetNuisanceParameters():
186
187 if nPlots >= nPlotsMax:
188 break
189
190 c2 = ROOT.TCanvas("c2")
191 frame = obs.frame()
192 frame.SetName(f"frame{nPlots}")
194
195 cut = "{0}=={0}::{1}".format(channelCat.GetName(), catName)
196 print(cut)
197 print(catName, " ", channelCat.getLabel())
199 frame,
200 MarkerSize=1,
201 Cut=cut,
202 DataError="None",
203 )
204
205 normCount = data.sumEntries(cut)
206
207 # remember the nominal value
208 value = var.getVal()
209
210 # w.allVars().Print("v")
211 # mc.GetNuisanceParameters().Print("v")
212 # pdftmp.plotOn(frame,LineWidth(2))
213 # mc.GetPdf().plotOn(frame,LineWidth(2.),Slice(channelCat,catName.c_str()),ProjWData(data))
214 # pdftmp.plotOn(frame,LineWidth(2.),Slice(channelCat,catName.c_str()),ProjWData(data))
215 global gobs
216 gobs = obs
217 global gpdftmp
218 gpdftmp = pdftmp
219
220 # notes: obs is RooRealVar / obstmp is RooArgSet
221 # pdftmp.expectedEvents receives RooArgSet as an argument
222 # in C++ automatic conversion is possible.
223 # in python the conversion is not possible.
224 # C-code : normCount = pdftmp->expectedEvents(*obs);
225 # Python : normCount = pdftmp.expectedEvents(obs) #doesn't work properly
226 # RooArgSet(obs) doesn´t reproduce well the results
227 # instead, we have to use obstmp
228 # normCount = pdftmp.expectedEvents(RooArgSet(obstmp)) #doesn´t work properly
229 normCount = pdftmp.expectedEvents(obstmp)
230 pdftmp.plotOn(frame, ROOT.RooFit.Normalization(normCount, ROOT.RooAbsReal.NumEvent), LineWidth=2) # nominal
231
232 var.setVal(value + nSigmaToVary * var.getError())
233
234 # pdftmp.plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2))
235 # mc.GetPdf().plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(channelCat,catName.c_str()),ProjWData(data))
236 # pdftmp.plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(channelCat,catName.c_str()),ProjWData(data))
237 normCount = pdftmp.expectedEvents(obstmp)
239 frame,
241 LineWidth=2,
242 LineColor="r",
243 LineStyle="--",
244 ) # +n sigma
245
246 var.setVal(value - nSigmaToVary * var.getError())
247
248 # pdftmp.plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2))
249 # mc.GetPdf().plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(channelCat,catName.c_str()),ProjWData(data))
250 # pdftmp.plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(channelCat,catName.c_str()),ProjWData(data))
251 normCount = pdftmp.expectedEvents(obstmp)
253 frame,
255 LineWidth=2,
256 LineColor="g",
257 LineStyle="--",
258 ) # -n sigma
259
260 # set them back to normal
261 var.setVal(value)
262
263 frameList.append(frame)
264
265 # quit making plots
266 nPlots += 1
267
268 frame.Draw()
269 c2.Update()
270 c2.Draw()
271 c2.SaveAs(f"StandardHistFactoryPlotsWithCategories.1.{catName}_{obs.GetName()}_{var.GetName()}.png")
272 del c2
273
274 # -------------------------------------------------------
275
276 # now make plots
277 c1 = ROOT.TCanvas("c1", "ProfileInspectorDemo", 800, 200)
278 nFrames = len(frameList)
279 if nFrames > 4:
280 nx = int(sqrt(nFrames))
281 ny = ROOT.TMath.CeilNint(nFrames / nx)
282 nx = ROOT.TMath.CeilNint(sqrt(nFrames))
283 c1.Divide(ny, nx)
284 else:
285 c1.Divide(nFrames)
286 for i in range(nFrames):
287 c1.cd(i + 1)
288 frameList[i].Draw()
289 c1.Update()
290 c1.Draw()
291 c1.SaveAs("StandardHistFactoryPlotsWithCategories.2.pdf")
292
293 file.Close()
294
295
296StandardHistFactoryPlotsWithCategories(
297 infile="", workspaceName="combined", modelConfigName="ModelConfig", dataName="obsData"
298)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format
void Print(Option_t *option="") const override
th1 Draw()