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