Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
StandardHypoTestInvDemo.py
Go to the documentation of this file.
1# \file
2# \ingroup tutorial_roostats
3# \notebook
4# Standard tutorial macro for performing an inverted hypothesis test for computing an interval
5#
6# This macro will perform a scan of the p-values for computing the interval or limit
7#
8# Usage:
9#
10# ~~~{.py}
11# ipython3> %run StandardHypoTestInvDemo.C
12# ipython3> StandardHypoTestInvDemo("fileName","workspace name","S+B modelconfig name","B model name","data set")
13# name",calculator type, test statistic type, use CLS,
14# number of points, xmin, xmax, number of toys, use number counting)
15#
16# type = 0 Freq calculator
17# type = 1 Hybrid calculator
18# type = 2 Asymptotic calculator
19# type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
20#
21# testStatType = 0 LEP
22# = 1 Tevatron
23# = 2 Profile Likelihood two sided
24# = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
25# = 4 Profile Likelihood signed ( pll = -pll if mu < mu_hat)
26# = 5 Max Likelihood Estimate as test statistic
27# = 6 Number of observed event as test statistic
28# ~~~
29#
30# \macro_image
31# \macro_output
32# \macro_code
33#
34# \author Lorenzo Moneta (C++ version), and P. P. (Python translation)
35
36import ROOT
37import os
38
39
40# structure defining the options
41class HypoTestInvOptions:
42
43 plotHypoTestResult = True # plot test statistic result at each point
44 writeResult = True # write HypoTestInverterResult in a file
45 resultFileName = (
46 ROOT.TString()
47 ) # file with results (by default is built automatically using the workspace input file name)
48 optimize = True # optimize evaluation of test statistic
49 useVectorStore = True # convert data to use roofit data store
50 generateBinned = False # generate binned data sets
51 noSystematics = False # force all systematics to be off (i.e. set all nuisance parameters as constat
52 # to their nominal values)
53 nToysRatio = 2 # ratio Ntoys S+b/ntoysB
54 maxPOI = -1 # max value used of POI (in case of auto scan)
55 enableDetailedOutput = (
56 False # enable detailed output with all fit information for each toys (output will be written in result file)
57 )
58 rebuild = False # re-do extra toys for computing expected limits and rebuild test stat
59 # distributions (N.B this requires much more CPU (factor is equivalent to nToyToRebuild)
60 nToyToRebuild = 100 # number of toys used to rebuild
61 rebuildParamValues = 0 # = 0 do a profile of all the parameters on the B (alt snapshot) before performing a
62 # rebuild operation (default)
63 # = 1 use initial workspace parameters with B snapshot values
64 # = 2 use all initial workspace parameters with B
65 # Otherwise the rebuild will be performed using
66 initialFit = -1 # do a first fit to the model (-1 : default, 0 skip fit, 1 do always fit)
67 randomSeed = -1 # random seed (if = -1: use default value, if = 0 always random )
68
69 nAsimovBins = 0 # number of bins in observables used for Asimov data sets (0 is the default and it is given by
70 # workspace, typically is 100)
71
72 reuseAltToys = False # reuse same toys for alternate hypothesis (if set one gets more stable bands)
73 confLevel = 0.95 # confidence level value
74
75 minimizerType = "" # minimizer type (default is what is in ROOT.Math.MinimizerOptions.DefaultMinimizerType()
76 massValue = "" # extra string to tag output file of result
77 printLevel = 0 # print level for debugging PL test statistics and calculators
78
79 useNLLOffset = False # use NLL offset when fitting (this increase stability of fits)
80
81
82optHTInv = HypoTestInvOptions()
83
84# internal class to run the inverter and more
85
86large_declaration = """ \
87namespace RooStats{
88
89class HypoTestInvTool {
90
91 public:
92 HypoTestInvTool();
93 ~HypoTestInvTool(){};
94
95 HypoTestInverterResult *RunInverter(RooWorkspace *w, const char *modelSBName, const char *modelBName,
96 const char *dataName, int type, int testStatType, bool useCLs, int npoints,
97 double poimin, double poimax, int ntoys, bool useNumberCounting = false,
98 const char *nuisPriorName = 0);
99
100 void AnalyzeResult(HypoTestInverterResult *r, int calculatorType, int testStatType, bool useCLs, int npoints,
101 const char *fileNameBase = 0);
102
103 void SetParameter(const char *name, const char *value);
104 void SetParameter(const char *name, bool value);
105 void SetParameter(const char *name, int value);
106 void SetParameter(const char *name, double value);
107
108 //private:
109 public:
110 bool mPlotHypoTestResult;
111 bool mWriteResult;
112 bool mOptimize;
113 bool mUseVectorStore;
114 bool mGenerateBinned;
115 bool mRebuild;
116 bool mReuseAltToys;
117 bool mEnableDetOutput;
118 int mNToyToRebuild;
119 int mRebuildParamValues;
120 int mPrintLevel;
121 int mInitialFit;
122 int mRandomSeed;
123 double mNToysRatio;
124 double mMaxPoi;
125 int mAsimovBins;
126 std::string mMassValue;
127 std::string
128 mMinimizerType; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType())
129 TString mResultFileName;
130};
131
132} // end using RooStats namespace
133
134RooStats::HypoTestInvTool::HypoTestInvTool()
135 : mPlotHypoTestResult(true), mWriteResult(false), mOptimize(true), mUseVectorStore(true), mGenerateBinned(false),
136 mEnableDetOutput(false), mRebuild(false), mReuseAltToys(false),
137 mNToyToRebuild(100), mRebuildParamValues(0), mPrintLevel(0), mInitialFit(-1), mRandomSeed(-1), mNToysRatio(2),
138 mMaxPoi(-1), mAsimovBins(0), mMassValue(""), mMinimizerType(""), mResultFileName()
139{
140}
141
142"""
143ROOT.gInterpreter.Declare(large_declaration)
144
145
146# Expanding definitions of class: HypoTestInvTool into HypoTestInvTool_plus
147class HypoTestInvTool_plus(ROOT.RooStats.HypoTestInvTool):
148 # RooStats.HypoTestInvTool.SetParameter(name, value):
149 def SetParameter(self, name, value): # SetParameter is polymorphic.
150 if (type(name) is str) and (type(value) is str):
151 #
152 # set boolean parameters
153 #
154
155 s_name = ROOT.std.string(name)
156
157 if s_name.find("PlotHypoTestResult") != ROOT.std.string.npos:
158 self.mPlotHypoTestResult = value
159 if s_name.find("WriteResult") != ROOT.std.string.npos:
160 self.mWriteResult = value
161 if s_name.find("Optimize") != ROOT.std.string.npos:
162 self.mOptimize = value
163 if s_name.find("UseVectorStore") != ROOT.std.string.npos:
164 self.mUseVectorStore = value
165 if s_name.find("GenerateBinned") != ROOT.std.string.npos:
166 self.mGenerateBinned = value
167 if s_name.find("EnableDetailedOutput") != ROOT.std.string.npos:
168 self.mEnableDetOutput = value
169 if s_name.find("Rebuild") != ROOT.std.string.npos:
170 self.mRebuild = value
171 if s_name.find("ReuseAltToys") != ROOT.std.string.npos:
172 self.mReuseAltToys = value
173
174 return
175 # RooStats.HypoTestInvTool.SetParameter = SetParameter
176
177 # RooStats.HypoTestInvTool.SetParameter(name, value):
178 # def SetParameter(self, name, value):
179 elif (type(name) is str) and (type(value) is bool):
180
181 #
182 # set integer parameters
183 #
184
185 s_name = ROOT.std.string(name)
186
187 if s_name.find("NToyToRebuild") != ROOT.std.string.npos:
188 self.mNToyToRebuild = value
189 if s_name.find("RebuildParamValues") != ROOT.std.string.npos:
190 self.mRebuildParamValues = value
191 if s_name.find("PrintLevel") != ROOT.std.string.npos:
192 self.mPrintLevel = value
193 if s_name.find("InitialFit") != ROOT.std.string.npos:
194 self.mInitialFit = value
195 if s_name.find("RandomSeed") != ROOT.std.string.npos:
196 self.mRandomSeed = value
197 if s_name.find("AsimovBins") != ROOT.std.string.npos:
198 self.mAsimovBins = value
199
200 return
201 # RooStats.HypoTestInvTool.SetParameter = SetParameter
202
203 # RooStats.HypoTestInvTool.SetParameter(name, value):
204 # def SetParameter(self, name, value):
205 elif (type(name) is str) and (type(value) is int):
206
207 #
208 # set double precision parameters
209 #
210
211 s_name = ROOT.std.string(name)
212
213 if s_name.find("NToysRatio") != ROOT.std.string.npos:
214 self.mNToysRatio = value
215 if s_name.find("MaxPOI") != ROOT.std.string.npos:
216 self.mMaxPoi = value
217
218 return
219 # RooStats.HypoTestInvTool.SetParameter = SetParameter
220
221 # RooStats.HypoTestInvTool.SetParameter(name, value):
222 # def SetParameter(self, name, value):
223 elif (type(name) is str) and (type(value) is (float or double)):
224
225 #
226 # set string parameters
227 #
228
229 s_name = ROOT.std.string(name)
230
231 if s_name.find("MassValue") != ROOT.std.string.npos:
232 global gselfmMassValue
233 gselfmMassValue = self.mMassValue
234 self.mMassValue = value # (self.mMassValue).assign(value)
235 if s_name.find("MinimizerType") != ROOT.std.string.npos:
236 self.mMassValue = value # self.mMinimizerType.assign(value)
237 if s_name.find("ResultFileName") != ROOT.std.string.npos:
238 self.mResultFileName = value
239
240 return
241
242 ROOT.RooStats.HypoTestInvTool.SetParameter = SetParameter
243
244 # --piece of code moved forwards...-
245 # ---
246
247 # RooStats.HypoTestInvTool.AnalyzeResult
248 def AnalyzeResult(self, r, calculatorType, testStatType, useCLs, npoints, fileNameBase):
249 # type(r) is HypoTestInverterResult
250
251 # analyze result produced by the inverter, optionally save it in a file
252
253 lowerLimit = 0
254 llError = 0
255 # if defined ROOT_SVN_VERSION && ROOT_SVN_VERSION >= 44126
256 if r.IsTwoSided():
257 lowerLimit = r.LowerLimit()
258 llError = r.LowerLimitEstimatedError()
259
260 # else
261 lowerLimit = r.LowerLimit()
262 llError = r.LowerLimitEstimatedError()
263 # endif
264
265 upperLimit = r.UpperLimit()
266 ulError = r.UpperLimitEstimatedError()
267
268 # ROOT.std::cout << "DEBUG : [ " << lowerLimit << " , " << upperLimit << " ] " << ROOT.std::endl;
269
270 if lowerLimit < upperLimit * (1.0 - 1.0e-4) and lowerLimit != 0:
271 print(f"The computed lower limit is: ", lowerLimit, " +/- ", llError)
272
273 print(f"The computed upper limit is: {upperLimit} +/- ", ulError)
274
275 # compute expected limit
276 print(f"Expected upper limits, using the B (alternate) model : ")
277 print(f" expected limit (median) ", r.GetExpectedUpperLimit(0))
278 print(f" expected limit (-1 sig) ", r.GetExpectedUpperLimit(-1))
279 print(f" expected limit (+1 sig) ", r.GetExpectedUpperLimit(1))
280 print(f" expected limit (-2 sig) ", r.GetExpectedUpperLimit(-2))
281 print(f" expected limit (+2 sig) ", r.GetExpectedUpperLimit(2))
282
283 # detailed output
284 if self.mEnableDetOutput:
285 self.mWriteResult = True
286 ROOT.Info("StandardHypoTestInvDemo", "detailed output will be written in output result file")
287
288 # write result in a file
289 if r != ROOT.kNone and self.mWriteResult:
290
291 # write to a file the results
292 calcType = "Freq" if (calculatorType == 0) else ("Hybr" if (calculatorType == 1) else "Asym")
293 limitType = "CLs" if (useCLs) else "Cls+b"
294 scanType = "auto" if (npoints < 0) else "grid"
295 if self.mResultFileName.IsNull():
296 self.mResultFileName = ROOT.TString.Format(
297 "%s_%s_%s_ts%d_", calcType, limitType, scanType, testStatType
298 )
299 # strip the / from the filename
300 if self.mMassValue.size() > 0:
301 self.mResultFileName += self.mMassValue
302 self.mResultFileName += "_"
303
304 name = fileNameBase
305 name.Replace(0, name.Last("/") + 1, "")
306 self.mResultFileName += name
307
308 # get (if existing) rebuilt UL distribution
309 uldistFile = "RULDist.root"
310 ulDist = 0
311 existULDist = not ROOT.gSystem.AccessPathName(uldistFile)
312 if existULDist:
313 fileULDist = TFile.Open(uldistFile)
314 if fileULDist:
315 ulDist = fileULDist.Get("RULDist")
316
317 fileOut = TFile(self.mResultFileName, "RECREATE")
318 r.Write()
319 if ulDist:
320 ulDist.Write()
321 ROOT.Info(
322 "StandardHypoTestInvDemo",
323 "HypoTestInverterResult has been written in the file %s",
324 self.mResultFileName.Data(),
325 )
326
327 fileOut.Close()
328
329 # plot the result ( p values vs scan points)
330 typeName = ""
331 if calculatorType == 0:
332 typeName = "Frequentist"
333 if calculatorType == 1:
334 typeName = "Hybrid"
335 elif calculatorType == 2 or calculatorType == 3:
336 typeName = "Asymptotic"
337 self.mPlotHypoTestResult = False
338
339 resultName = r.GetName()
340 plotTitle = "{} CL Scan for workspace {}".format(str(typeName), resultName)
341 global gr
342 gr = r
343 plot = ROOT.RooStats.HypoTestInverterPlot("HTI_Result_Plot", plotTitle, r)
344
345 # plot in a new canvas with style
346 c1Name = "{}_Scan".format(typeName)
347 c1 = ROOT.TCanvas(c1Name)
348 c1.SetLogy(False)
349
350 plot.Draw("CLb 2CL") # plot all and Clb
351 # if (useCLs):
352 # plot.Draw("CLb 2CL") # plot all and Clb ???
353 # else:
354 # plot.Draw("") # plot all and Clb ???
355 c1.Update()
356 c1.Draw()
357 c1.SaveAs("StandardHypoTestInvDemo.c1.png")
358
359 nEntries = r.ArraySize()
360
361 # plot test statistics distributions for the two hypothesis
362 if self.mPlotHypoTestResult:
363 c2 = ROOT.TCanvas("c2")
364 if nEntries > 1:
365 ny = ROOT.TMath.CeilNint(ROOT.TMath.Sqrt(nEntries))
366 nx = ROOT.TMath.CeilNint(float(nEntries) / ny)
367 c2.Divide(nx, ny)
368
369 for i in range(nEntries):
370 if nEntries > 1:
371 c2.cd(i + 1)
372 pl = plot.MakeTestStatPlot(i)
373 pl.SetLogYaxis(True)
374 pl.Draw()
375 c2.Update()
376 c2.Draw()
377 c2.SaveAs("StandardHypoTestInvDemo.c2.png")
378
379 ROOT.gPad = c1
380
381 ROOT.RooStats.HypoTestInvTool.AnalyzeResult = AnalyzeResult
382
383 # internal routine to run the inverter
384 # RooStats.HypoTestInvTool.RunInverter
385 def RunInverter(
386 self,
387 w,
388 modelSBName,
389 modelBName,
390 dataName,
391 Type,
392 testStatType,
393 useCLs,
394 npoints,
395 poimin,
396 poimax,
397 ntoys,
398 useNumberCounting,
399 nuisPriorName,
400 ):
401
402 print(f"Running HypoTestInverter on the workspace ", w.GetName())
403
404 w.Print()
405
406 data = w.data(dataName)
407 if not data:
408 Error("StandardHypoTestDemo", "Not existing data {}".format(dataName))
409 return 0
410 else:
411 print(f"Using data set ", dataName)
412
413 if self.mUseVectorStore:
414 ROOT.RooAbsData.setDefaultStorageType(ROOT.RooAbsData.Vector)
415 data.convertToVectorStore()
416
417 # get models from WS
418 # get the modelConfig out of the file
419 bModel = w.obj(modelBName)
420 sbModel = w.obj(modelSBName)
421
422 if not sbModel:
423 Error("StandardHypoTestDemo", "Not existing ModelConfig %s", modelSBName)
424 return 0
425
426 # check the model
427 if not sbModel.GetPdf():
428 Error("StandardHypoTestDemo", "Model %s has no pdf ", modelSBName)
429 return 0
430
431 if not sbModel.GetParametersOfInterest():
432 Error("StandardHypoTestDemo", "Model %s has no poi ", modelSBName)
433 return 0
434
435 if not sbModel.GetObservables():
436 Error("StandardHypoTestInvDemo", "Model %s has no observables ", modelSBName)
437 return 0
438
439 if not sbModel.GetSnapshot():
440 ROOT.Info(
441 "StandardHypoTestInvDemo", "Model {} has no snapshot - make one using model poi".format(modelSBName)
442 )
443 sbModel.SetSnapshot(sbModel.GetParametersOfInterest())
444
445 # case of no systematics
446 # remove nuisance parameters from model
447 if optHTInv.noSystematics:
448 nuisPar = sbModel.GetNuisanceParameters()
449 if nuisPar and nuisPar.getSize() > 0:
450 print(
451 f"StandardHypoTestInvDemo",
452 " - Switch off all systematics by setting them constant to their initial values",
453 )
454 ROOT.RooStats.SetAllConstant(nuisPar)
455
456 if bModel:
457 bnuisPar = bModel.GetNuisanceParameters()
458 if bnuisPar:
459 ROOT.RooStats.SetAllConstant(bnuisPar)
460
461 if (not bModel) or (bModel == sbModel):
462 ROOT.Info("StandardHypoTestInvDemo", "The background model {} does not exist".format(modelBName))
463 ROOT.Info("StandardHypoTestInvDemo", "Copy it from ModelConfig {} and set POI to zero".format(modelSBName))
464 bModel = sbModel.Clone()
465 bModel.SetName(modelSBName + "_with_poi_0")
466 var = bModel.GetParametersOfInterest().first()
467 if not var:
468 return 0
469 oldval = var.getVal()
470 var.setVal(0)
471 bModel.SetSnapshot(ROOT.RooArgSet(var))
472 var.setVal(oldval)
473 else:
474 if not bModel.GetSnapshot():
475 ROOT.Info(
476 "StandardHypoTestInvDemo",
477 "Model %s has no snapshot - make one using model poi and 0 values ",
478 modelBName,
479 )
480 var = bModel.GetParametersOfInterest().first()
481 if var:
482 oldval = var.getVal()
483 var.setVal(0)
484 bModel.SetSnapshot(ROOT.RooArgSet(var))
485 var.setVal(oldval)
486 else:
487 Error("StandardHypoTestInvDemo", "Model %s has no valid poi", modelBName)
488 return 0
489
490 # check model has global observables when there are nuisance pdf
491 # for the hybrid case the globals are not needed
492 if Type != 1:
493 hasNuisParam = sbModel.GetNuisanceParameters() and sbModel.GetNuisanceParameters().getSize() > 0
494 hasGlobalObs = sbModel.GetGlobalObservables() and sbModel.GetGlobalObservables().getSize() > 0
495 if hasNuisParam and not hasGlobalObs:
496 # try to see if model has nuisance parameters first
497 constrPdf = RooStats.MakeNuisancePdf(sbModel, "nuisanceConstraintPdf_sbmodel")
498 if constrPdf:
499 Warning(
500 "StandardHypoTestInvDemo",
501 "Model %s has nuisance parameters but no global observables associated",
502 sbModel.GetName(),
503 )
504 Warning(
505 "StandardHypoTestInvDemo",
506 "\tThe effect of the nuisance parameters will not be treated correctly ",
507 )
508
509 # save all initial parameters of the model including the global observables
510 initialParameters = ROOT.RooArgSet()
511 allParams = sbModel.GetPdf().getParameters(data)
512 allParams.snapshot(initialParameters)
513
514 # run first a data fit
515
516 poiSet = sbModel.GetParametersOfInterest()
517 poi = poiSet.first()
518
519 print("StandardHypoTestInvDemo : POI initial value: ", poi.GetName(), " = ", poi.getVal())
520
521 # fit the data first (need to use constraint )
522 tw = ROOT.TStopwatch()
523
524 doFit = self.mInitialFit
525 if testStatType == 0 and self.mInitialFit == -1:
526 doFit = False # case of LEP test statistic
527 if type == 3 and self.mInitialFit == -1:
528 doFit = False # case of Asymptoticcalculator with nominal Asimov
529 poihat = 0
530
531 if len(self.mMinimizerType) == 0:
533 else:
534 ROOT.Math.MinimizerOptions.SetDefaultMinimizer(str(self.mMinimizerType))
535
536 ROOT.Info(
537 "StandardHypoTestInvDemo",
538 "Using {} as minimizer for computing the test statistic".format(
540 ),
541 )
542
543 if doFit:
544
545 # do the fit : By doing a fit the POI snapshot (for S+B) is set to the fit value
546 # and the nuisance parameters nominal values will be set to the fit value.
547 # This is relevant when using LEP test statistics
548
549 ROOT.Info("StandardHypoTestInvDemo", " Doing a first fit to the observed data ")
550 constrainParams = RooArgSet()
551 if sbModel.GetNuisanceParameters():
552 constrainParams.add(sbModel.GetNuisanceParameters())
553 RooStats.RemoveConstantParameters(constrainParams)
554 tw.Start()
555 fitres = sbModel.GetPdf().fitTo(
556 data,
557 InitialHesse(False),
558 Hesse(False),
559 Minimizer(str(self.mMinimizerType), "Migrad"),
560 Strategy(0),
561 PrintLevel(self.mPrintLevel),
562 Constrain(constrainParams),
563 Save(True),
564 Offset(RooStats.IsNLLOffset()),
565 )
566
567 if fitres.status() != 0:
568 Warning(
569 "StandardHypoTestInvDemo",
570 "Fit to the model failed - try with strategy 1 and perform first an Hesse computation",
571 )
572 fitres = sbModel.GetPdf().fitTo(
573 data,
574 InitialHesse(True),
575 Hesse(False),
576 Minimizer(self.mMinimizerType, "Migrad"),
577 Strategy(1),
578 PrintLevel(self.mPrintLevel + 1),
579 Constrain(constrainParams),
580 Save(True),
581 Offset(RooStats.IsNLLOffset()),
582 )
583
584 if fitres.status() != 0:
585 Warning("StandardHypoTestInvDemo", " Fit still failed - continue anyway.....")
586
587 print("StandardHypoTestInvDemo - Best Fit value : ", poi.GetName(), " = ", poihat, " +/- ", poi.getError())
588 print(f"Time for fitting : ")
589 tw.Print()
590
591 # save best fit value in the poi snapshot
592 sbModel.SetSnapshot(sbModel.GetParametersOfInterest())
593 print(f"StandardHypoTestInvo: snapshot of S+B Model ", sbModel.GetName(), " is set to the best fit value")
594
595 # print a message in case of LEP test statistics because it affects result by doing or not doing a fit
596 if testStatType == 0:
597 if not doFit:
598 ROOT.Info(
599 "StandardHypoTestInvDemo",
600 "Using LEP test statistic - an initial fit is not done and the TS will use "
601 + "the nuisances at the model value",
602 )
603 else:
604 ROOT.Info(
605 "StandardHypoTestInvDemo",
606 "Using LEP test statistic - an initial fit has been done and the TS will use "
607 + "the nuisances at the best fit value",
608 )
609
610 # build test statistics and hypotest calculators for running the inverter
611
612 slrts = ROOT.RooStats.SimpleLikelihoodRatioTestStat(sbModel.GetPdf(), bModel.GetPdf())
613
614 # null parameters must includes snapshot of poi plus the nuisance values
615 nullParams = ROOT.RooArgSet(sbModel.GetSnapshot())
616 if sbModel.GetNuisanceParameters():
617 nullParams.add(sbModel.GetNuisanceParameters())
618 if sbModel.GetSnapshot():
619 slrts.SetNullParameters(nullParams)
620 altParams = ROOT.RooArgSet(bModel.GetSnapshot())
621 if bModel.GetNuisanceParameters():
622 altParams.add(bModel.GetNuisanceParameters())
623 if bModel.GetSnapshot():
624 slrts.SetAltParameters(altParams)
625 if self.mEnableDetOutput:
626 slrts.EnableDetailedOutput()
627
628 # ratio of profile likelihood - need to pass snapshot for the alt
629 ropl = ROOT.RooStats.RatioOfProfiledLikelihoodsTestStat(sbModel.GetPdf(), bModel.GetPdf(), bModel.GetSnapshot())
630 ropl.SetSubtractMLE(False)
631 if testStatType == 11:
632 ropl.SetSubtractMLE(True)
633 ropl.SetPrintLevel(self.mPrintLevel)
634 ropl.SetMinimizer(self.mMinimizerType.c_str())
635 if self.mEnableDetOutput:
636 ropl.EnableDetailedOutput()
637
638 profll = ROOT.RooStats.ProfileLikelihoodTestStat(sbModel.GetPdf())
639 if testStatType == 3:
640 profll.SetOneSided(True)
641 if testStatType == 4:
642 profll.SetSigned(True)
643 profll.SetMinimizer(self.mMinimizerType.c_str())
644 profll.SetPrintLevel(self.mPrintLevel)
645 if self.mEnableDetOutput:
646 profll.EnableDetailedOutput()
647
648 profll.SetReuseNLL(self.mOptimize)
649 slrts.SetReuseNLL(self.mOptimize)
650 ropl.SetReuseNLL(self.mOptimize)
651
652 if self.mOptimize:
653 profll.SetStrategy(0)
654 ropl.SetStrategy(0)
656
657 if self.mMaxPoi > 0:
658 poi.setMax(self.mMaxPoi) # increase limit
659
660 maxll = ROOT.RooStats.MaxLikelihoodEstimateTestStat(sbModel.GetPdf(), poi)
661 nevtts = ROOT.RooStats.NumEventsTestStat()
662
663 ROOT.RooStats.AsymptoticCalculator.SetPrintLevel(self.mPrintLevel)
664
665 # create the HypoTest calculator class
666 hc = ROOT.nullptr
667 if Type == 0:
668 hc = ROOT.RooStats.FrequentistCalculator(data, bModel, sbModel)
669 elif Type == 1:
670 hc = ROOT.RooStats.HybridCalculator(data, bModel, sbModel)
671 # elif (Type == 2 ):
672 # hc = AsymptoticCalculator(data, bModel, sbModel, false, self.mAsimovBins)
673 # elif (Type == 3 ):
674 # hc = AsymptoticCalculator(data, bModel, sbModel, True, self.mAsimovBins) # for using
675 # Asimov data generated with nominal values
676 elif Type == 2:
677 hc = ROOT.RooStats.AsymptoticCalculator(data, bModel, sbModel, False)
678 elif Type == 3:
679 hc = ROOT.RooStats.AsymptoticCalculator(
680 data, bModel, sbModel, True
681 ) # for using Asimov data generated with nominal values
682 else:
683 Error(
684 "StandardHypoTestInvDemo",
685 "Invalid - calculator type = {Type}supported values are only :\n\t\t\t 0 "
686 + "(Frequentist) , 1 (Hybrid) , 2 (Asymptotic) ",
687 )
688
689 return 0
690
691 # set the test statistic
692 testStat = 0
693 if testStatType == 0:
694 testStat = slrts
695 if testStatType == 1 or testStatType == 11:
696 testStat = ropl
697 if testStatType == 2 or testStatType == 3 or testStatType == 4:
698 testStat = profll
699 if testStatType == 5:
700 testStat = maxll
701 if testStatType == 6:
702 testStat = nevtts
703
704 if testStat == 0:
705 Error(
706 "StandardHypoTestInvDemo",
707 "Invalid - test statistic type = {testStatType} supported values are only :\n\t\t\t 0 (SLR) "
708 + ", 1 (Tevatron) , 2 (PLR), 3 (PLR1), 4(MLE)",
709 )
710
711 return 0
712
713 toymcs = hc.GetTestStatSampler()
714 if toymcs and (Type == 0 or Type == 1):
715 # look if pdf is number counting or extended
716 if sbModel.GetPdf().canBeExtended():
717 if useNumberCounting:
718 Warning("StandardHypoTestInvDemo", "Pdf is extended: but number counting flag is set: ignore it ")
719 else:
720 # for not extended pdf
721 if not useNumberCounting:
722 nEvents = data.numEntries()
723 ROOT.Info(
724 "StandardHypoTestInvDemo",
725 "Pdf is not extended: number of events to generate taken from observed data set is {nEvents}",
726 )
727 toymcs.SetNEventsPerToy(nEvents)
728 else:
729 ROOT.Info("StandardHypoTestInvDemo", "using a number counting pdf")
730 toymcs.SetNEventsPerToy(1)
731
732 toymcs.SetTestStatistic(testStat)
733
734 if data.isWeighted() and not self.mGenerateBinned:
735 ROOT.Info(
736 "StandardHypoTestInvDemo",
737 (
738 "Data set is weighted, nentries = {} and sum of weights = {:8.1f} but toy"
739 + "generation is unbinned - it would be faster to set self.mGenerateBinned to true\n"
740 ).format(data.numEntries(), data.sumEntries()),
741 )
742
743 toymcs.SetGenerateBinned(self.mGenerateBinned)
744
745 toymcs.SetUseMultiGen(self.mOptimize)
746
747 if self.mGenerateBinned and sbModel.GetObservables().getSize() > 2:
748 Warning(
749 "StandardHypoTestInvDemo",
750 (
751 "generate binned is activated but the number of observable is {}. Too much "
752 + "memory could be needed for allocating all the bins"
753 ).format(sbModel.GetObservables().getSize()),
754 )
755
756 # set the random seed if needed
757 if self.mRandomSeed >= 0:
758 RooRandom.randomGenerator().SetSeed(self.mRandomSeed)
759
760 # specify if need to re-use same toys
761 if self.mReuseAltToys:
762 hc.UseSameAltToys()
763
764 if Type == 1:
765 hhc = HybridCalculator(hc)
766 assert hhc
767
768 hhc.SetToys(ntoys, ntoys / self.mNToysRatio) # can use less ntoys for b hypothesis
769
770 # remove global observables from ModelConfig (this is probably not needed anymore in 5.32)
771 bModel.SetGlobalObservables(RooArgSet())
772 sbModel.SetGlobalObservables(RooArgSet())
773
774 # check for nuisance prior pdf in case of nuisance parameters
775 if bModel.GetNuisanceParameters() or sbModel.GetNuisanceParameters():
776
777 # fix for using multigen (does not work in this case)
778 toymcs.SetUseMultiGen(False)
779 ToyMCSampler.SetAlwaysUseMultiGen(False)
780
781 nuisPdf = 0
782 if nuisPriorName:
783 nuisPdf = w.pdf(nuisPriorName)
784 # use prior defined first in bModel (then in SbModel)
785 if not nuisPdf:
786 ROOT.Info(
787 "StandardHypoTestInvDemo",
788 "No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model",
789 )
790 if bModel.GetPdf() and bModel.GetObservables():
791 nuisPdf = RooStats.MakeNuisancePdf(bModel, "nuisancePdf_bmodel")
792 else:
793 nuisPdf = RooStats.MakeNuisancePdf(sbModel, "nuisancePdf_sbmodel")
794
795 if not nuisPdf:
796 if bModel.GetPriorPdf():
797 nuisPdf = bModel.GetPriorPdf()
798 ROOT.Info(
799 "StandardHypoTestInvDemo",
800 "No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",
801 nuisPdf.GetName(),
802 )
803 else:
804 Error(
805 "StandardHypoTestInvDemo",
806 "Cannot run Hybrid calculator because no prior on the nuisance "
807 "parameter is specified or can be derived",
808 )
809 return 0
810
811 assert nuisPdf
812 ROOT.Info("StandardHypoTestInvDemo", "Using as nuisance Pdf ... ")
813 nuisPdf.Print()
814
815 nuisParams = (
816 ibModel.GetNuisanceParameters()
817 if (bModel.GetNuisanceParameters())
818 else (sbModel.GetNuisanceParameters())
819 )
820 npnuisPdf.getObservables(nuisParams)
821 if np.getSize() == 0:
822 Warning(
823 "StandardHypoTestInvDemo",
824 "Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range",
825 )
826
827 hhc.ForcePriorNuisanceAlt(nuisPdf)
828 hhc.ForcePriorNuisanceNull(nuisPdf)
829
830 elif type == 2 or type == 3:
831 if testStatType == 3:
832 hc.SetOneSided(True)
833 if testStatType != 2 and testStatType != 3:
834 Warning(
835 "StandardHypoTestInvDemo",
836 "Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL",
837 )
838 elif type == 0:
839 hc.SetToys(ntoys, ntoys / self.mNToysRatio)
840 # store also the fit information for each poi point used by calculator based on toys
841 if self.mEnableDetOutput:
842 hc.StoreFitInfo(True)
843 elif type == 1:
844 hc.SetToys(ntoys, ntoys / self.mNToysRatio)
845 # store also the fit information for each poi point used by calculator based on toys
846 # if (self.mEnableDetOutput):
847 # hc.StoreFitInfo(True)
848
849 # Get the result
850 ROOT.RooMsgService.instance().getStream(1).removeTopic(ROOT.RooFit.NumIntegration)
851
852 calc = ROOT.RooStats.HypoTestInverter(hc)
853 calc.SetConfidenceLevel(optHTInv.confLevel)
854
855 calc.UseCLs(useCLs)
856 calc.SetVerbose(True)
857
858 if npoints > 0:
859 if poimin > poimax:
860 # if no min/max given scan between MLE and +4 sigma
861 poimin = int(poihat)
862 poimax = int(poihat + 4 * poi.getError())
863
864 print(f"Doing a fixed scan in interval : {poimin}, {poimax} ")
865 calc.SetFixedScan(npoints, poimin, poimax)
866 else:
867 # poi.setMax(10*int( (poihat+ 10 *poi.getError() )/10 ) )
868 print(f"Doing an automatic scan in interval : {poi.getMin()} , ", poi.getMax())
869
870 tw.Start()
871 r = calc.GetInterval()
872 print(f"Time to perform limit scan \n")
873 tw.Print()
874
875 if self.mRebuild:
876
877 print(f"\n\n")
878 print(
879 f"Rebuild the upper limit distribution by re-generating new set of pseudo-experiment and re-compute "
880 + "for each of them a new upper limit\n\n"
881 )
882
883 allParams = sbModel.GetPdf().getParameters(data)
884
885 # define on which value of nuisance parameters to do the rebuild
886 # default is best fit value for bmodel snapshot
887
888 if self.mRebuildParamValues != 0:
889 # set all parameters to their initial workspace values
890 allParams.assign(initialParameters)
891
892 if self.mRebuildParamValues == 0 or self.mRebuildParamValues == 1:
893 constrainParams = RooArgSet()
894 if sbModel.GetNuisanceParameters():
895 constrainParams.add(sbModel.GetNuisanceParameters())
896 RooStats.RemoveConstantParameters(constrainParams)
897
898 poiModel = sbModel.GetParametersOfInterest()
899 bModel.LoadSnapshot()
900
901 # do a profile using the B model snapshot
902 if self.mRebuildParamValues == 0:
903
904 RooStats.SetAllConstant(poiModel, True)
905
906 sbModel.GetPdf().fitTo(
907 data,
908 InitialHesse(False),
909 Hesse(False),
910 Minimizer(self.mMinimizerType, "Migrad"),
911 Strategy(0),
912 PrintLevel(self.mPrintLevel),
913 Constrain(constrainParams),
914 Offset(RooStats.IsNLLOffset()),
915 )
916
917 print(f"rebuild using fitted parameter value for B-model snapshot")
918 constrainParams.Print("v")
919
920 RooStats.SetAllConstant(poiModel, False)
921
922 print(f"StandardHypoTestInvDemo: Initial parameters used for rebuilding: ")
923 RooStats.PrintListContent(allParams, ROOT.std.cout)
924
925 tw.Start()
926 limDist = calc.GetUpperLimitDistribution(True, self.mNToyToRebuild)
927 print(f"Time to rebuild distributions :")
928 tw.Print()
929
930 if limDist:
931 print(f"Expected limits after rebuild distribution ")
932 print(f"expected upper limit (median of limit distribution) ", limDist.InverseCDF(0.5))
933 print(
934 f"expected -1 sig limit (0.16% quantile of limit dist) ",
935 limDist.InverseCDF(ROOT.Math.normal_cdf(-1)),
936 )
937 print(
938 f"expected +1 sig limit (0.84% quantile of limit dist) ",
939 limDist.InverseCDF(ROOT.Math.normal_cdf(1)),
940 )
941 print(
942 f"expected -2 sig limit (.025% quantile of limit dist) ",
943 limDist.InverseCDF(ROOT.Math.normal_cdf(-2)),
944 )
945 print(
946 f"expected +2 sig limit (.975% quantile of limit dist) ",
947 limDist.InverseCDF(ROOT.Math.normal_cdf(2)),
948 )
949
950 # Plot the upper limit distribution
951 limPlot = SamplingDistPlot(50 if (self.mNToyToRebuild < 200) else 100)
952 limPlot.AddSamplingDistribution(limDist)
953 limPlot.GetTH1F().SetStats(True) # display statistics
954 limPlot.SetLineColor(kBlue)
955 c1 = TCanvas("limPlot", "Upper Limit Distribution")
956 limPlot.Draw()
957 c1.Update()
958 c1.Draw()
959 c1.SaveAs("StandardHypoTestInvDemo.1.png")
960
961 # save result in a file
962 limDist.SetName("RULDist")
963 fileOut = TFile("RULDist.root", "RECREATE")
964 limDist.Write()
965 fileOut.Close()
966
967 # update r to a new updated result object containing the rebuilt expected p-values distributions
968 # (it will not recompute the expected limit)
969 if r:
970 del r # need to delete previous object since GetInterval will return a cloned copy
971 r = calc.GetInterval()
972
973 else:
974 print(f"ERROR : failed to re-build distributions ")
975
976 return r
977
978 ROOT.RooStats.HypoTestInvTool.RunInverter = RunInverter
979
980
981# end class HypoTestInvTool
982# Loading new definitions of HypoTestInvTool_plus into HypoTestInvTool ...
983ROOT.RooStats.HypoTestInvTool = HypoTestInvTool_plus
984# Could seem redundant but it would help you if you want expand new-new definitions.
985
986
987# --------------------------------------------------
989 infile=0,
990 wsName="combined",
991 modelSBName="ModelConfig",
992 modelBName="",
993 dataName="obsData",
994 calculatorType=0,
995 testStatType=0,
996 useCLs=True,
997 npoints=6,
998 poimin=0,
999 poimax=5,
1000 ntoys=1000,
1001 useNumberCounting=False,
1002 nuisPriorName=0,
1003):
1004 """
1005
1006 Other Parameter to pass in tutorial
1007 apart from standard for filename, ws, modelconfig and data
1008
1009 type = 0 Freq calculator
1010 type = 1 Hybrid calculator
1011 type = 2 Asymptotic calculator
1012 type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
1013
1014 testStatType = 0 LEP
1015 Tevatron
1016 Likelihood
1017 mu_hat)
1018 mu_hat)
1019 statistic
1020 statistic
1021
1022 useCLs scan for CLs (otherwise for CLs+b)
1023
1024 npoints = -1
1025
1026 poimin,poimax: min/max value to scan in case of fixed scans
1027 (if min > max, try to find automatically)
1028
1029 ntoys: number of toys to use
1030
1031 useNumberCounting: set to True when using number counting events
1032
1033 nuisPriorName: name of prior for the nuisance. This is often expressed as constraint term in the global model
1034 HybridCalculator (type=1)
1035 If not given by default the prior pdf from ModelConfig is used.
1036
1037 extra options are available as global parameters of the macro. They major ones are:
1038
1039 plotHypoTestResult plot result of tests at each point (TS distributions) (default is True)
1040 writeResult write result of scan (default is True)
1041 rebuild rebuild scan for expected limits (require extra toys) (default is False)
1042 generateBinned generate binned data sets for toys (default is False) - be careful not to activate with
1043 large (>=3) number of observables
1044 nToyRatio ratio of S+B/B toys (default is 2)
1045 """
1046
1047 filename = ROOT.TString(infile)
1048 if filename.IsNull():
1049 filename = "results/example_combined_GaussExample_model.root"
1050 fileExist = not ROOT.gSystem.AccessPathName(filename) # note opposite return code
1051 # if file does not exists generate with histfactory
1052 if not fileExist:
1053 if os.name == "nt":
1054 print("HistFactory file cannto be generated on Windows - exit")
1055 return
1056 # Normally this would be run on the command line
1057 print(f"will run standard hist2workspace example")
1058 ROOT.gROOT.ProcessLine(".! prepareHistFactory .")
1059 ROOT.gROOT.ProcessLine(".! hist2workspace config/example.xml")
1060 print(f"\n\n---------------------")
1061 print(f"Done creating example input")
1062 print(f"---------------------\n\n")
1063
1064 else:
1065 filename = infile
1066
1067 # Try to open the file
1068 file = ROOT.TFile.Open(filename)
1069
1070 # if input file was specified but not found, quit
1071 if not file:
1072 print(f"StandardRooStatsDemoMacro: Input file {filename} is not found")
1073 return
1074
1075 calc = ROOT.RooStats.HypoTestInvTool() # instance of HypoTestInvTool
1076 # calc = HypoTestInvTool_plus() # instance of HypoTestInvTool
1077
1078 # set parameters
1079 calc.SetParameter("PlotHypoTestResult", optHTInv.plotHypoTestResult)
1080 calc.SetParameter("WriteResult", optHTInv.writeResult)
1081 calc.SetParameter("Optimize", optHTInv.optimize)
1082 calc.SetParameter("UseVectorStore", optHTInv.useVectorStore)
1083 calc.SetParameter("GenerateBinned", optHTInv.generateBinned)
1084 calc.SetParameter("NToysRatio", optHTInv.nToysRatio)
1085 calc.SetParameter("MaxPOI", optHTInv.maxPOI)
1086 calc.SetParameter("EnableDetailedOutput", optHTInv.enableDetailedOutput)
1087 calc.SetParameter("Rebuild", optHTInv.rebuild)
1088 calc.SetParameter("ReuseAltToys", optHTInv.reuseAltToys)
1089 calc.SetParameter("NToyToRebuild", optHTInv.nToyToRebuild)
1090 calc.SetParameter("RebuildParamValues", optHTInv.rebuildParamValues)
1091 calc.SetParameter("MassValue", optHTInv.massValue)
1092 calc.SetParameter("MinimizerType", optHTInv.minimizerType)
1093 calc.SetParameter("PrintLevel", optHTInv.printLevel)
1094 calc.SetParameter("InitialFit", optHTInv.initialFit)
1095 calc.SetParameter("ResultFileName", optHTInv.resultFileName)
1096 calc.SetParameter("RandomSeed", optHTInv.randomSeed)
1097 calc.SetParameter("AsimovBins", optHTInv.nAsimovBins)
1098
1099 # enable offset for all roostats
1100 if optHTInv.useNLLOffset:
1102
1103 w = file.Get(wsName)
1104 r = 0
1105 print(w, "\t", filename)
1106 if w != ROOT.kNone:
1107 r = calc.RunInverter(
1108 w,
1109 modelSBName,
1110 modelBName,
1111 dataName,
1112 calculatorType,
1113 testStatType,
1114 useCLs,
1115 npoints,
1116 poimin,
1117 poimax,
1118 ntoys,
1119 useNumberCounting,
1120 nuisPriorName,
1121 )
1122 if not r:
1123 raise RuntimeError("Error running the HypoTestInverter - Exit ")
1124 return
1125
1126 else:
1127 # case workspace is not present look for the inverter result
1128 print(f"Reading an HypoTestInverterResult with name {wsName} from file " << filename)
1129 r = file.Get(wsName) # dynamic_cast
1130 if not r:
1131 raise RuntimeError("File {filename} does not contain a workspace or an HypoTestInverterResult - Exit ")
1132
1133 file.ls()
1134 return
1135
1136 calc.AnalyzeResult(r, calculatorType, testStatType, useCLs, npoints, infile)
1137
1138 return
1139
1140
1141# --------------------------------------------------
1142
1143
1144def ReadResult(fileName, resultName="", useCLs=True):
1145 # read a previous stored result from a file given the result name
1146 StandardHypoTestInvDemo(fileName, resultName, "", "", "", 0, 0, useCLs)
1147
1148
1149# ifdef USE_AS_MAIN
1150def main():
1152
1153
1154main()
1155# endif
int main()
Definition Prototype.cxx:12
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
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
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 type
static void SetDefaultMinimizer(const char *type, const char *algo=nullptr)
Set the default Minimizer type and corresponding algorithms.
static void SetDefaultStrategy(int strat)
Set the default strategy.
static const std::string & DefaultMinimizerType()
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:48
The Canvas class.
Definition TCanvas.h:23
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
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
utility function to set all variable constant in a collection (from G.
void RemoveConstantParameters(RooArgSet *set)
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
void UseNLLOffset(bool on)
function to set a global flag in RooStats to use NLL offset when performing nll computations Note tha...
bool IsNLLOffset()
function returning if the flag to check if the flag to use NLLOffset is set
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
useful function to print in one line the content of a set with their values