Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HybridInstructional.py
Go to the documentation of this file.
1# \file
2# \ingroup tutorial_roostats
3# \notebook -js
4# Example demonstrating usage of HybridCalcultor
5#
6# A hypothesis testing example based on number counting
7# with background uncertainty.
8#
9# NOTE: This example must be run with the ACLIC (the + option ) due to the
10# new class that is defined.
11#
12# This example:
13# - demonstrates the usage of the HybridCalcultor (Part 4-6)
14# - demonstrates the numerical integration of RooFit (Part 2)
15# - validates the RooStats against an example with a known analytic answer
16# - demonstrates usage of different test statistics
17# - explains subtle choices in the prior used for hybrid methods
18# - demonstrates usage of different priors for the nuisance parameters
19#
20# The basic setup here is that a main measurement has observed x events with an
21# expectation of s+b. One can choose an ad hoc prior for the uncertainty on b,
22# or try to base it on an auxiliary measurement. In this case, the auxiliary
23# measurement (aka control measurement, sideband) is another counting experiment
24# with measurement y and expectation tau*b. With an 'original prior' on b,
25# called \f$\eta(b)\f$ then one can obtain a posterior from the auxiliary measurement
26# \f$\pi(b) = \eta(b) * Pois(y|tau*b).\f$ This is a principled choice for a prior
27# on b in the main measurement of x, which can then be treated in a hybrid
28# Bayesian/Frequentist way. Additionally, one can try to treat the two
29# measurements simultaneously, which is detailed in Part 6 of the tutorial.
30#
31# This tutorial is related to the FourBin.C tutorial in the modeling, but
32# focuses on hypothesis testing instead of interval estimation.
33#
34# More background on this 'prototype problem' can be found in the
35# following papers:
36#
37# - Evaluation of three methods for calculating statistical significance
38# when incorporating a systematic uncertainty into a test of the
39# background-only hypothesis for a Poisson process
40# Authors: Robert D. Cousins, James T. Linnemann, Jordan Tucker
41# http://arxiv.org/abs/physics/0702156
42# NIM A 595 (2008) 480--501
43#
44# - Statistical Challenges for Searches for New Physics at the LHC
45# Author: Kyle Cranmer
46# http://arxiv.org/abs/physics/0511028
47#
48# - Measures of Significance in HEP and Astrophysics
49# Authors J. T. Linnemann
50# http://arxiv.org/abs/physics/0312059
51#
52# \macro_image
53# \macro_output
54# \macro_code
55#
56# \authors Kyle Cranmer, Wouter Verkerke, Sven Kreiss, and Jolly Chen (Python translation)
57
58import ROOT
59
60# User configuration parameters
61ntoys = 6000
62nToysRatio = 20 # ratio Ntoys Null/ntoys ALT
63
64
65# ----------------------------------
66# A New Test Statistic Class for this example.
67# It simply returns the sum of the values in a particular
68# column of a dataset.
69# You can ignore this class and focus on the macro below
71 """
72using namespace RooFit;
73using namespace RooStats;
74
75class BinCountTestStat : public TestStatistic {
76public:
77 BinCountTestStat(void) : fColumnName("tmp") {}
78 BinCountTestStat(string columnName) : fColumnName(columnName) {}
79
80 virtual Double_t Evaluate(RooAbsData &data, RooArgSet & /*nullPOI*/)
81 {
82 // This is the main method in the interface
83 Double_t value = 0.0;
84 for (int i = 0; i < data.numEntries(); i++) {
85 value += data.get(i)->getRealValue(fColumnName.c_str());
86 }
87 return value;
88 }
89 virtual const TString GetVarName() const { return fColumnName; }
90
91private:
92 string fColumnName;
93
94protected:
95 ClassDef(BinCountTestStat, 1)
96};
97"""
98)
99
100# ----------------------------------
101# The Actual Tutorial Macro
102
103# This tutorial has 6 parts
104# Table of Contents
105# Setup
106# 1. Make the model for the 'prototype problem'
107# Special cases
108# 2. Use RooFit's direct integration to get p-value & significance
109# 3. Use RooStats analytic solution for this problem
110# RooStats HybridCalculator -- can be generalized
111# 4. RooStats ToyMC version of 2. & 3.
112# 5. RooStats ToyMC with an equivalent test statistic
113# 6. RooStats ToyMC with simultaneous control & main measurement
114
115t = ROOT.TStopwatch()
116t.Start()
117c = ROOT.TCanvas()
118c.Divide(2, 2)
119
120# ----------------------------------------------------
121# P A R T 1 : D I R E C T I N T E G R A T I O N
122# ====================================================
123# Make model for prototype on/off problem
124# $Pois(x | s+b) * Pois(y | tau b )$
125# for Z_Gamma, use uniform prior on b.
126w = ROOT.RooWorkspace("w")
127w.factory("Poisson::px(x[150,0,500],sum::splusb(s[0,0,100],b[100,0.1,300]))")
128w.factory("Poisson::py(y[100,0.1,500],prod::taub(tau[1.],b))")
129w.factory("PROD::model(px,py)")
130w.factory("Uniform::prior_b(b)")
131
132# We will control the output level in a few places to avoid
133# verbose progress messages. We start by keeping track
134# of the current threshold on messages.
135msglevel = ROOT.RooMsgService.instance().globalKillBelow()
136
137# ----------------------------------------------------
138# P A R T 2 : D I R E C T I N T E G R A T I O N
139# ====================================================
140# This is not the 'RooStats' way, but in this case the distribution
141# of the test statistic is simply x and can be calculated directly
142# from the PDF using RooFit's built-in integration.
143# Note, this does not generalize to situations in which the test statistic
144# depends on many events (rows in a dataset).
145
146# construct the Bayesian-averaged model (eg. a projection pdf)
147# $p'(x|s) = \int db p(x|s+b) * [ p(y|b) * prior(b) ]$
148w.factory("PROJ::averagedModel(PROD::foo(px|b,py,prior_b),b)")
149
151# lower message level
152# plot it, red is averaged model, green is b known exactly, blue is s+b av model
153frame = w.var("x").frame(Range=(50, 230))
154w.pdf("averagedModel").plotOn(frame, LineColor="r")
155w.pdf("px").plotOn(frame, LineColor="g")
156w.var("s").setVal(50.0)
157w.pdf("averagedModel").plotOn(frame, LineColor="b")
158c.cd(1)
160w.var("s").setVal(0.0)
161
162# compare analytic calculation of Z_Bi
163# with the numerical RooFit implementation of Z_Gamma
164# for an example with x = 150, y = 100
165
166# numeric RooFit Z_Gamma
167w.var("y").setVal(100)
168w.var("x").setVal(150)
169cdf = w.pdf("averagedModel").createCdf(w.var("x"))
171# get ugly print messages out of the way
172print("-----------------------------------------")
173print("Part 2")
174print(f"Hybrid p-value from direct integration = {1 - cdf.getVal()}")
175print(f"Z_Gamma Significance = {ROOT.RooStats.PValueToSignificance(1 - cdf.getVal())}")
176ROOT.RooMsgService.instance().setGlobalKillBelow(msglevel)
177# set it back
178
179# ---------------------------------------------
180# P A R T 3 : A N A L Y T I C R E S U L T
181# =============================================
182# In this special case, the integrals are known analytically
183# and they are implemented in NumberCountingUtils
184
185# analytic Z_Bi
188print("-----------------------------------------")
189print("Part 3")
190print(f"Z_Bi p-value (analytic): {p_Bi}")
191print(f"Z_Bi significance (analytic): {Z_Bi}")
192t.Stop()
193t.Print()
194t.Reset()
195t.Start()
196
197# -------------------------------------------------------------
198# P A R T 4 : U S I N G H Y B R I D C A L C U L A T O R
199# =============================================================
200# Now we demonstrate the RooStats HybridCalculator.
201#
202# Like all RooStats calculators it needs the data and a ModelConfig
203# for the relevant hypotheses. Since we are doing hypothesis testing
204# we need a ModelConfig for the null (background only) and the alternate
205# (signal+background) hypotheses. We also need to specify the PDF,
206# the parameters of interest, and the observables. Furthermore, since
207# the parameter of interest is floating, we need to specify which values
208# of the parameter corresponds to the null and alternate (eg. s=0 and s=50)
209#
210# define some sets of variables obs={x} and poi={s}
211# note here, x is the only observable in the main measurement
212# and y is treated as a separate measurement, which is used
213# to produce the prior that will be used in this calculation
214# to randomize the nuisance parameters.
215w.defineSet("obs", "x")
216w.defineSet("poi", "s")
217
218# create a toy dataset with the x=150
219data = ROOT.RooDataSet("d", "d", w.set("obs"))
220data.add(w.set("obs"))
221
222# Part 3a : Setup ModelConfigs
223# -------------------------------------------------------
224# create the null (background-only) ModelConfig with s=0
225b_model = ROOT.RooStats.ModelConfig("B_model", w)
229w.var("s").setVal(0.0)
230# important!
232
233# create the alternate (signal+background) ModelConfig with s=50
234sb_model = ROOT.RooStats.ModelConfig("S+B_model", w)
238w.var("s").setVal(50.0)
239# important!
241
242# Part 3b : Choose Test Statistic
243# ----------------------------------
244# To make an equivalent calculation we need to use x as the test
245# statistic. This is not a built-in test statistic in RooStats
246# so we define it above. The new class inherits from the
247# TestStatistic interface, and simply returns the value
248# of x in the dataset.
249
250binCount = ROOT.BinCountTestStat("x")
251
252# Part 3c : Define Prior used to randomize nuisance parameters
253# -------------------------------------------------------------
254#
255# The prior used for the hybrid calculator is the posterior
256# from the auxiliary measurement y. The model for the aux.
257# measurement is Pois(y|tau*b), thus the likelihood function
258# is proportional to (has the form of) a Gamma distribution.
259# if the 'original prior' $\eta(b)$ is uniform, then from
260# Bayes's theorem we have the posterior:
261# $\pi(b) = Pois(y|tau*b) * \eta(b)$
262# If $\eta(b)$ is flat, then we arrive at a Gamma distribution.
263# Since RooFit will normalize the PDF we can actually supply
264# $py=Pois(y,tau*b)$ that will be equivalent to multiplying by a uniform.
265#
266# Alternatively, we could explicitly use a gamma distribution:
267# `w.factory("Gamma::gamma(b,sum::temp(y,1),1,0)");`
268#
269# or we can use some other ad hoc prior that do not naturally
270# follow from the known form of the auxiliary measurement.
271# The common choice is the equivalent Gaussian:
272w.factory("Gaussian::gauss_prior(b,y, expr::sqrty('sqrt(y)',y))")
273# this corresponds to the "Z_N" calculation.
274#
275# or one could use the analogous log-normal prior
276w.factory("Lognormal::lognorm_prior(b,y, expr::kappa('1+1./sqrt(y)',y))")
277#
278# Ideally, the HybridCalculator would be able to inspect the full
279# model $Pois(x | s+b) * Pois(y | tau b )$ and be given the original
280# prior $\eta(b)$ to form $\pi(b) = Pois(y|tau*b) * \eta(b)$.
281# This is not yet implemented because in the general case
282# it is not easy to identify the terms in the PDF that correspond
283# to the auxiliary measurement. So for now, it must be set
284# explicitly with:
285# - ForcePriorNuisanceNull()
286# - ForcePriorNuisanceAlt()
287# the name "ForcePriorNuisance" was chosen because we anticipate
288# this to be auto-detected, but will leave the option open
289# to force to a different prior for the nuisance parameters.
290
291# Part 3d : Construct and configure the HybridCalculator
292# -------------------------------------------------------
293
294hc1 = ROOT.RooStats.HybridCalculator(data, sb_model, b_model)
296toymcs1 = hc1.GetTestStatSampler()
298# because the model is in number counting form
300# set the test statistic
301hc1.SetToys(ntoys, ntoys // nToysRatio)
304# if you wanted to use the ad hoc Gaussian prior instead
305# ~~~
306# hc1.ForcePriorNuisanceAlt(w.pdf("gauss_prior"))
307# hc1.ForcePriorNuisanceNull(w.pdf("gauss_prior"))
308# ~~~
309# if you wanted to use the ad hoc log-normal prior instead
310# ~~~
311# hc1.ForcePriorNuisanceAlt(w.pdf("lognorm_prior"))
312# hc1.ForcePriorNuisanceNull(w.pdf("lognorm_prior"))
313# ~~~
314
315# these lines save current msg level and then kill any messages below ERROR
317# Get the result
318r1 = hc1.GetHypoTest()
319ROOT.RooMsgService.instance().setGlobalKillBelow(msglevel)
320# set it back
321print("-----------------------------------------")
322print("Part 4")
323r1.Print()
324t.Stop()
325t.Print()
326t.Reset()
327t.Start()
328
329c.cd(2)
331# 30 bins,TS is discrete
332p1.Draw()
333
334# -------------------------------------------------------------------------
335# # P A R T 5 : U S I N G H Y B R I D C A L C U L A T O R
336# # W I T H A N A L T E R N A T I V E T E S T S T A T I S T I C
337#
338# A likelihood ratio test statistics should be 1-to-1 with the count x
339# when the value of b is fixed in the likelihood. This is implemented
340# by the SimpleLikelihoodRatioTestStat
341
345
346# HYBRID CALCULATOR
347hc2 = ROOT.RooStats.HybridCalculator(data, sb_model, b_model)
349toymcs2 = hc2.GetTestStatSampler()
352hc2.SetToys(ntoys, ntoys // nToysRatio)
355# if you wanted to use the ad hoc Gaussian prior instead
356# ~~~
357# hc2.ForcePriorNuisanceAlt(w.pdf("gauss_prior"))
358# hc2.ForcePriorNuisanceNull(w.pdf("gauss_prior"))
359# ~~~
360# if you wanted to use the ad hoc log-normal prior instead
361# ~~~
362# hc2.ForcePriorNuisanceAlt(w.pdf("lognorm_prior"))
363# hc2.ForcePriorNuisanceNull(w.pdf("lognorm_prior"))
364# ~~~
365
366# these lines save current msg level and then kill any messages below ERROR
368# Get the result
369r2 = hc2.GetHypoTest()
370print("-----------------------------------------")
371print("Part 5")
372r2.Print()
373t.Stop()
374t.Print()
375t.Reset()
376t.Start()
377ROOT.RooMsgService.instance().setGlobalKillBelow(msglevel)
378
379c.cd(3)
381# 30 bins
382p2.Draw()
383
384# -----------------------------------------------------------------------------
385# # P A R T 6 : U S I N G H Y B R I D C A L C U L A T O R W I T H A N A L T E R N A T I V E T E S T
386# # S T A T I S T I C A N D S I M U L T A N E O U S M O D E L
387#
388# If one wants to use a test statistic in which the nuisance parameters
389# are profiled (in one way or another), then the PDF must constrain b.
390# Otherwise any observation x can always be explained with s=0 and b=x/tau.
391#
392# In this case, one is really thinking about the problem in a
393# different way. They are considering x,y simultaneously.
394# and the PDF should be $Pois(x | s+b) * Pois(y | tau b )$
395# and the set 'obs' should be {x,y}.
396
397w.defineSet("obsXY", "x,y")
398
399# create a toy dataset with the x=150, y=100
400w.var("x").setVal(150.0)
401w.var("y").setVal(100.0)
402dataXY = ROOT.RooDataSet("dXY", "dXY", w.set("obsXY"))
403dataXY.add(w.set("obsXY"))
404
405# now we need new model configs, with PDF="model"
406b_modelXY = ROOT.RooStats.ModelConfig("B_modelXY", w)
407b_modelXY.SetPdf(w.pdf("model"))
408# IMPORTANT
411w.var("s").setVal(0.0)
412# IMPORTANT
414
415# create the alternate (signal+background) ModelConfig with s=50
416sb_modelXY = ROOT.RooStats.ModelConfig("S+B_modelXY", w)
417sb_modelXY.SetPdf(w.pdf("model"))
418# IMPORTANT
421w.var("s").setVal(50.0)
422# IMPORTANT
424
425# Test statistics like the profile likelihood ratio
426# (or the ratio of profiled likelihoods (Tevatron) or the MLE for s)
427# will now work, since the nuisance parameter b is constrained by y.
428# ratio of alt and null likelihoods with background yield profiled.
429#
430# NOTE: These are slower because they have to run fits for each toy
431
432# Tevatron-style Ratio of profiled likelihoods
433# $Q_Tev = -log L(s=0,\hat\hat{b})/L(s=50,\hat\hat{b})$
436)
438
439# profile likelihood where alternate is best fit value of signal yield
440# $\lambda(0) = -log L(s=0,\hat\hat{b})/L(\hat{s},\hat{b})$
442
443# just use the maximum likelihood estimate of signal yield
444# $MLE = \hat{s}$
446
447# However, it is less clear how to justify the prior used in randomizing
448# the nuisance parameters (since that is a property of the ensemble,
449# and y is a property of each toy pseudo experiment. In that case,
450# one probably wants to consider a different y0 which will be held
451# constant and the prior $\pi(b) = Pois(y0 | tau b) * \eta(b)$.
452w.factory("y0[100]")
453w.factory("Gamma::gamma_y0(b,sum::temp0(y0,1),1,0)")
454w.factory("Gaussian::gauss_prior_y0(b,y0, expr::sqrty0('sqrt(y0)',y0))")
455
456# HYBRID CALCULATOR
457hc3 = ROOT.RooStats.HybridCalculator(dataXY, sb_modelXY, b_modelXY)
459toymcs3 = hc3.GetTestStatSampler()
462hc3.SetToys(ntoys, ntoys // nToysRatio)
465# if you wanted to use the ad hoc Gaussian prior instead
466# ~~~{.cpp}
467# hc3.ForcePriorNuisanceAlt(w.pdf("gauss_prior_y0"))
468# hc3.ForcePriorNuisanceNull(w.pdf("gauss_prior_y0"))
469# ~~~
470
471# choose fit-based test statistic
473# toymcs3.SetTestStatistic(ropl)
474# toymcs3.SetTestStatistic(mlets)
475
476# these lines save current msg level and then kill any messages below ERROR
478# Get the result
479r3 = hc3.GetHypoTest()
480print("-----------------------------------------")
481print("Part 6")
482r3.Print()
483t.Stop()
484t.Print()
485t.Reset()
486t.Start()
487ROOT.RooMsgService.instance().setGlobalKillBelow(msglevel)
488
489c.cd(4)
490c.GetPad(4).SetLogy()
492# 50 bins
493p3.Draw()
494
495c.SaveAs("zbi.pdf")
496
497# -----------------------------------------
498# OUTPUT W/O PROOF (2.66 GHz Intel Core i7)
499# =========================================
500
501# -----------------------------------------
502# Part 2
503# Hybrid p-value from direct integration = 0.00094165
504# Z_Gamma Significance = 3.10804
505# -----------------------------------------
506#
507# Part 3
508# Z_Bi p-value (analytic): 0.00094165
509# Z_Bi significance (analytic): 3.10804
510# Real time 0:00:00, CP time 0.610
511
512# -----------------------------------------
513# Part 4
514# Results HybridCalculator_result:
515# - Null p-value = 0.00115 +/- 0.000228984
516# - Significance = 3.04848 sigma
517# - Number of S+B toys: 1000
518# - Number of B toys: 20000
519# - Test statistic evaluated on data: 150
520# - CL_b: 0.99885 +/- 0.000239654
521# - CL_s+b: 0.476 +/- 0.0157932
522# - CL_s: 0.476548 +/- 0.0158118
523# Real time 0:00:07, CP time 7.620
524
525# -----------------------------------------
526# Part 5
527# Results HybridCalculator_result:
528# - Null p-value = 0.0009 +/- 0.000206057
529# - Significance = 3.12139 sigma
530# - Number of S+B toys: 1000
531# - Number of B toys: 20000
532# - Test statistic evaluated on data: 10.8198
533# - CL_b: 0.9991 +/- 0.000212037
534# - CL_s+b: 0.465 +/- 0.0157726
535# - CL_s: 0.465419 +/- 0.0157871
536# Real time 0:00:34, CP time 34.360
537
538# -----------------------------------------
539# Part 6
540# Results HybridCalculator_result:
541# - Null p-value = 0.000666667 +/- 0.000149021
542# - Significance = 3.20871 sigma
543# - Number of S+B toys: 1000
544# - Number of B toys: 30000
545# - Test statistic evaluated on data: 5.03388
546# - CL_b: 0.999333 +/- 0.000149021
547# - CL_s+b: 0.511 +/- 0.0158076
548# - CL_s: 0.511341 +/- 0.0158183
549# Real time 0:05:06, CP time 306.330
550
551# ---------------------------------------------------------
552# OUTPUT w/ PROOF (2.66 GHz Intel Core i7, 4 virtual cores)
553# =========================================================
554
555# -----------------------------------------
556# Part 5
557# Results HybridCalculator_result:
558# - Null p-value = 0.00075 +/- 0.000173124
559# - Significance = 3.17468 sigma
560# - Number of S+B toys: 1000
561# - Number of B toys: 20000
562# - Test statistic evaluated on data: 10.8198
563# - CL_b: 0.99925 +/- 0.000193577
564# - CL_s+b: 0.454 +/- 0.0157443
565# - CL_s: 0.454341 +/- 0.0157564
566# Real time 0:00:16, CP time 0.990
567
568# -----------------------------------------
569# Part 6
570# Results HybridCalculator_result:
571# - Null p-value = 0.0007 +/- 0.000152699
572# - Significance = 3.19465 sigma
573# - Number of S+B toys: 1000
574# - Number of B toys: 30000
575# - Test statistic evaluated on data: 5.03388
576# - CL_b: 0.9993 +/- 0.000152699
577# - CL_s+b: 0.518 +/- 0.0158011
578# - CL_s: 0.518363 +/- 0.0158124
579# Real time 0:01:25, CP time 0.580
580
581# ----------------------------------
582# Comparison
583# ----------------------------------
584# LEPStatToolsForLHC
585# https:#plone4.fnal.gov:4430/P0/phystat/packages/0703002
586# Uses Gaussian prior
587# CL_b = 6.218476e-04, Significance = 3.228665 sigma
588#
589# ----------------------------------
590# Comparison
591# ----------------------------------
592# Asymptotic
593# From the value of the profile likelihood ratio (5.0338)
594# The significance can be estimated using Wilks's theorem
595# significance = sqrt(2*profileLR) = 3.1729 sigma
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.