Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HybridInstructional.C
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, and Sven Kreiss
58
59#include "RooGlobalFunc.h"
60#include "RooRealVar.h"
61#include "RooProdPdf.h"
62#include "RooWorkspace.h"
63#include "RooDataSet.h"
64#include "TCanvas.h"
65#include "TStopwatch.h"
66#include "TH1.h"
67#include "RooPlot.h"
68#include "RooMsgService.h"
69
71
75
80
81using namespace RooFit;
82using namespace RooStats;
83
84// ----------------------------------
85// A New Test Statistic Class for this example.
86// It simply returns the sum of the values in a particular
87// column of a dataset.
88// You can ignore this class and focus on the macro below
89class BinCountTestStat : public TestStatistic {
90public:
91 BinCountTestStat(void) : fColumnName("tmp") {}
92 BinCountTestStat(string columnName) : fColumnName(columnName) {}
93
94 virtual Double_t Evaluate(RooAbsData &data, RooArgSet & /*nullPOI*/)
95 {
96 // This is the main method in the interface
97 Double_t value = 0.0;
98 for (int i = 0; i < data.numEntries(); i++) {
99 value += data.get(i)->getRealValue(fColumnName.c_str());
100 }
101 return value;
102 }
103 virtual const TString GetVarName() const { return fColumnName; }
104
105private:
106 string fColumnName;
107
108protected:
109 ClassDef(BinCountTestStat, 1)
110};
111
112ClassImp(BinCountTestStat)
113
114// ----------------------------------
115// The Actual Tutorial Macro
116
117void HybridInstructional(int ntoys = 6000)
118{
119 double nToysRatio = 20; // ratio Ntoys Null/ntoys ALT
120
121 // This tutorial has 6 parts
122 // Table of Contents
123 // Setup
124 // 1. Make the model for the 'prototype problem'
125 // Special cases
126 // 2. Use RooFit's direct integration to get p-value & significance
127 // 3. Use RooStats analytic solution for this problem
128 // RooStats HybridCalculator -- can be generalized
129 // 4. RooStats ToyMC version of 2. & 3.
130 // 5. RooStats ToyMC with an equivalent test statistic
131 // 6. RooStats ToyMC with simultaneous control & main measurement
132
133 // It takes ~4 min without PROOF and ~2 min with PROOF on 4 cores.
134 // Of course, everything looks nicer with more toys, which takes longer.
135
136 TStopwatch t;
137 t.Start();
138 TCanvas *c = new TCanvas;
139 c->Divide(2, 2);
140
141 // ----------------------------------------------------
142 // P A R T 1 : D I R E C T I N T E G R A T I O N
143 // ====================================================
144 // Make model for prototype on/off problem
145 // $Pois(x | s+b) * Pois(y | tau b )$
146 // for Z_Gamma, use uniform prior on b.
147 RooWorkspace *w = new RooWorkspace("w");
148 w->factory("Poisson::px(x[150,0,500],sum::splusb(s[0,0,100],b[100,0.1,300]))");
149 w->factory("Poisson::py(y[100,0.1,500],prod::taub(tau[1.],b))");
150 w->factory("PROD::model(px,py)");
151 w->factory("Uniform::prior_b(b)");
152
153 // We will control the output level in a few places to avoid
154 // verbose progress messages. We start by keeping track
155 // of the current threshold on messages.
157
158 // Use PROOF-lite on multi-core machines
159 ProofConfig *pc = NULL;
160 // uncomment below if you want to use PROOF
161 // ~~~
162 // pc = new ProofConfig(*w, 4, "workers=4", kFALSE); // machine with 4 cores
163 // pc = new ProofConfig(*w, 2, "workers=2", kFALSE); // machine with 2 cores
164 // ~~~
165
166 // ----------------------------------------------------
167 // P A R T 2 : D I R E C T I N T E G R A T I O N
168 // ====================================================
169 // This is not the 'RooStats' way, but in this case the distribution
170 // of the test statistic is simply x and can be calculated directly
171 // from the PDF using RooFit's built-in integration.
172 // Note, this does not generalize to situations in which the test statistic
173 // depends on many events (rows in a dataset).
174
175 // construct the Bayesian-averaged model (eg. a projection pdf)
176 // $p'(x|s) = \int db p(x|s+b) * [ p(y|b) * prior(b) ]$
177 w->factory("PROJ::averagedModel(PROD::foo(px|b,py,prior_b),b)");
178
180 // plot it, red is averaged model, green is b known exactly, blue is s+b av model
181 RooPlot *frame = w->var("x")->frame(Range(50, 230));
182 w->pdf("averagedModel")->plotOn(frame, LineColor(kRed));
183 w->pdf("px")->plotOn(frame, LineColor(kGreen));
184 w->var("s")->setVal(50.);
185 w->pdf("averagedModel")->plotOn(frame, LineColor(kBlue));
186 c->cd(1);
187 frame->Draw();
188 w->var("s")->setVal(0.);
189
190 // compare analytic calculation of Z_Bi
191 // with the numerical RooFit implementation of Z_Gamma
192 // for an example with x = 150, y = 100
193
194 // numeric RooFit Z_Gamma
195 w->var("y")->setVal(100);
196 w->var("x")->setVal(150);
197 std::unique_ptr<RooAbsReal> cdf{w->pdf("averagedModel")->createCdf(*w->var("x"))};
198 cdf->getVal(); // get ugly print messages out of the way
199 cout << "-----------------------------------------" << endl;
200 cout << "Part 2" << endl;
201 cout << "Hybrid p-value from direct integration = " << 1 - cdf->getVal() << endl;
202 cout << "Z_Gamma Significance = " << PValueToSignificance(1 - cdf->getVal()) << endl;
203 RooMsgService::instance().setGlobalKillBelow(msglevel); // set it back
204
205 // ---------------------------------------------
206 // P A R T 3 : A N A L Y T I C R E S U L T
207 // =============================================
208 // In this special case, the integrals are known analytically
209 // and they are implemented in RooStats::NumberCountingUtils
210
211 // analytic Z_Bi
212 double p_Bi = NumberCountingUtils::BinomialWithTauObsP(150, 100, 1);
213 double Z_Bi = NumberCountingUtils::BinomialWithTauObsZ(150, 100, 1);
214 cout << "-----------------------------------------" << endl;
215 cout << "Part 3" << endl;
216 std::cout << "Z_Bi p-value (analytic): " << p_Bi << std::endl;
217 std::cout << "Z_Bi significance (analytic): " << Z_Bi << std::endl;
218 t.Stop();
219 t.Print();
220 t.Reset();
221 t.Start();
222
223 // -------------------------------------------------------------
224 // 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
225 // =============================================================
226 // Now we demonstrate the RooStats HybridCalculator.
227 //
228 // Like all RooStats calculators it needs the data and a ModelConfig
229 // for the relevant hypotheses. Since we are doing hypothesis testing
230 // we need a ModelConfig for the null (background only) and the alternate
231 // (signal+background) hypotheses. We also need to specify the PDF,
232 // the parameters of interest, and the observables. Furthermore, since
233 // the parameter of interest is floating, we need to specify which values
234 // of the parameter corresponds to the null and alternate (eg. s=0 and s=50)
235 //
236 // define some sets of variables obs={x} and poi={s}
237 // note here, x is the only observable in the main measurement
238 // and y is treated as a separate measurement, which is used
239 // to produce the prior that will be used in this calculation
240 // to randomize the nuisance parameters.
241 w->defineSet("obs", "x");
242 w->defineSet("poi", "s");
243
244 // create a toy dataset with the x=150
245 RooDataSet *data = new RooDataSet("d", "d", *w->set("obs"));
246 data->add(*w->set("obs"));
247
248 // Part 3a : Setup ModelConfigs
249 // -------------------------------------------------------
250 // create the null (background-only) ModelConfig with s=0
251 ModelConfig b_model("B_model", w);
252 b_model.SetPdf(*w->pdf("px"));
253 b_model.SetObservables(*w->set("obs"));
254 b_model.SetParametersOfInterest(*w->set("poi"));
255 w->var("s")->setVal(0.0); // important!
256 b_model.SetSnapshot(*w->set("poi"));
257
258 // create the alternate (signal+background) ModelConfig with s=50
259 ModelConfig sb_model("S+B_model", w);
260 sb_model.SetPdf(*w->pdf("px"));
261 sb_model.SetObservables(*w->set("obs"));
262 sb_model.SetParametersOfInterest(*w->set("poi"));
263 w->var("s")->setVal(50.0); // important!
264 sb_model.SetSnapshot(*w->set("poi"));
265
266 // Part 3b : Choose Test Statistic
267 // ----------------------------------
268 // To make an equivalent calculation we need to use x as the test
269 // statistic. This is not a built-in test statistic in RooStats
270 // so we define it above. The new class inherits from the
271 // RooStats::TestStatistic interface, and simply returns the value
272 // of x in the dataset.
273
274 BinCountTestStat binCount("x");
275
276 // Part 3c : Define Prior used to randomize nuisance parameters
277 // -------------------------------------------------------------
278 //
279 // The prior used for the hybrid calculator is the posterior
280 // from the auxiliary measurement y. The model for the aux.
281 // measurement is Pois(y|tau*b), thus the likelihood function
282 // is proportional to (has the form of) a Gamma distribution.
283 // if the 'original prior' $\eta(b)$ is uniform, then from
284 // Bayes's theorem we have the posterior:
285 // $\pi(b) = Pois(y|tau*b) * \eta(b)$
286 // If $\eta(b)$ is flat, then we arrive at a Gamma distribution.
287 // Since RooFit will normalize the PDF we can actually supply
288 // $py=Pois(y,tau*b)$ that will be equivalent to multiplying by a uniform.
289 //
290 // Alternatively, we could explicitly use a gamma distribution:
291 // `w->factory("Gamma::gamma(b,sum::temp(y,1),1,0)");`
292 //
293 // or we can use some other ad hoc prior that do not naturally
294 // follow from the known form of the auxiliary measurement.
295 // The common choice is the equivalent Gaussian:
296 w->factory("Gaussian::gauss_prior(b,y, expr::sqrty('sqrt(y)',y))");
297 // this corresponds to the "Z_N" calculation.
298 //
299 // or one could use the analogous log-normal prior
300 w->factory("Lognormal::lognorm_prior(b,y, expr::kappa('1+1./sqrt(y)',y))");
301 //
302 // Ideally, the HybridCalculator would be able to inspect the full
303 // model $Pois(x | s+b) * Pois(y | tau b )$ and be given the original
304 // prior $\eta(b)$ to form $\pi(b) = Pois(y|tau*b) * \eta(b)$.
305 // This is not yet implemented because in the general case
306 // it is not easy to identify the terms in the PDF that correspond
307 // to the auxiliary measurement. So for now, it must be set
308 // explicitly with:
309 // - ForcePriorNuisanceNull()
310 // - ForcePriorNuisanceAlt()
311 // the name "ForcePriorNuisance" was chosen because we anticipate
312 // this to be auto-detected, but will leave the option open
313 // to force to a different prior for the nuisance parameters.
314
315 // Part 3d : Construct and configure the HybridCalculator
316 // -------------------------------------------------------
317
318 HybridCalculator hc1(*data, sb_model, b_model);
319 ToyMCSampler *toymcs1 = (ToyMCSampler *)hc1.GetTestStatSampler();
320 toymcs1->SetNEventsPerToy(1); // because the model is in number counting form
321 toymcs1->SetTestStatistic(&binCount); // set the test statistic
322 hc1.SetToys(ntoys, ntoys / nToysRatio);
323 hc1.ForcePriorNuisanceAlt(*w->pdf("py"));
324 hc1.ForcePriorNuisanceNull(*w->pdf("py"));
325 // if you wanted to use the ad hoc Gaussian prior instead
326 // ~~~
327 // hc1.ForcePriorNuisanceAlt(*w->pdf("gauss_prior"));
328 // hc1.ForcePriorNuisanceNull(*w->pdf("gauss_prior"));
329 // ~~~
330 // if you wanted to use the ad hoc log-normal prior instead
331 // ~~~
332 // hc1.ForcePriorNuisanceAlt(*w->pdf("lognorm_prior"));
333 // hc1.ForcePriorNuisanceNull(*w->pdf("lognorm_prior"));
334 // ~~~
335
336 // enable proof
337 // NOTE: This test statistic is defined in this macro, and is not
338 // working with PROOF currently. Luckily test stat is fast to evaluate.
339 // `if(pc) toymcs1->SetProofConfig(pc);`
340
341 // these lines save current msg level and then kill any messages below ERROR
343 // Get the result
344 HypoTestResult *r1 = hc1.GetHypoTest();
345 RooMsgService::instance().setGlobalKillBelow(msglevel); // set it back
346 cout << "-----------------------------------------" << endl;
347 cout << "Part 4" << endl;
348 r1->Print();
349 t.Stop();
350 t.Print();
351 t.Reset();
352 t.Start();
353
354 c->cd(2);
355 HypoTestPlot *p1 = new HypoTestPlot(*r1, 30); // 30 bins, TS is discrete
356 p1->Draw();
357
358 // -------------------------------------------------------------------------
359 // # 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
360 // # 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
361 //
362 // A likelihood ratio test statistics should be 1-to-1 with the count x
363 // when the value of b is fixed in the likelihood. This is implemented
364 // by the SimpleLikelihoodRatioTestStat
365
366 SimpleLikelihoodRatioTestStat slrts(*b_model.GetPdf(), *sb_model.GetPdf());
367 slrts.SetNullParameters(*b_model.GetSnapshot());
368 slrts.SetAltParameters(*sb_model.GetSnapshot());
369
370 // HYBRID CALCULATOR
371 HybridCalculator hc2(*data, sb_model, b_model);
372 ToyMCSampler *toymcs2 = (ToyMCSampler *)hc2.GetTestStatSampler();
373 toymcs2->SetNEventsPerToy(1);
374 toymcs2->SetTestStatistic(&slrts);
375 hc2.SetToys(ntoys, ntoys / nToysRatio);
376 hc2.ForcePriorNuisanceAlt(*w->pdf("py"));
377 hc2.ForcePriorNuisanceNull(*w->pdf("py"));
378 // if you wanted to use the ad hoc Gaussian prior instead
379 // ~~~
380 // hc2.ForcePriorNuisanceAlt(*w->pdf("gauss_prior"));
381 // hc2.ForcePriorNuisanceNull(*w->pdf("gauss_prior"));
382 // ~~~
383 // if you wanted to use the ad hoc log-normal prior instead
384 // ~~~
385 // hc2.ForcePriorNuisanceAlt(*w->pdf("lognorm_prior"));
386 // hc2.ForcePriorNuisanceNull(*w->pdf("lognorm_prior"));
387 // ~~~
388
389 // enable proof
390 if (pc)
391 toymcs2->SetProofConfig(pc);
392
393 // these lines save current msg level and then kill any messages below ERROR
395 // Get the result
396 HypoTestResult *r2 = hc2.GetHypoTest();
397 cout << "-----------------------------------------" << endl;
398 cout << "Part 5" << endl;
399 r2->Print();
400 t.Stop();
401 t.Print();
402 t.Reset();
403 t.Start();
405
406 c->cd(3);
407 HypoTestPlot *p2 = new HypoTestPlot(*r2, 30); // 30 bins
408 p2->Draw();
409
410 // -----------------------------------------------------------------------------
411 // # 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
412 // # 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
413 //
414 // If one wants to use a test statistic in which the nuisance parameters
415 // are profiled (in one way or another), then the PDF must constrain b.
416 // Otherwise any observation x can always be explained with s=0 and b=x/tau.
417 //
418 // In this case, one is really thinking about the problem in a
419 // different way. They are considering x,y simultaneously.
420 // and the PDF should be $Pois(x | s+b) * Pois(y | tau b )$
421 // and the set 'obs' should be {x,y}.
422
423 w->defineSet("obsXY", "x,y");
424
425 // create a toy dataset with the x=150, y=100
426 w->var("x")->setVal(150.);
427 w->var("y")->setVal(100.);
428 RooDataSet *dataXY = new RooDataSet("dXY", "dXY", *w->set("obsXY"));
429 dataXY->add(*w->set("obsXY"));
430
431 // now we need new model configs, with PDF="model"
432 ModelConfig b_modelXY("B_modelXY", w);
433 b_modelXY.SetPdf(*w->pdf("model")); // IMPORTANT
434 b_modelXY.SetObservables(*w->set("obsXY"));
435 b_modelXY.SetParametersOfInterest(*w->set("poi"));
436 w->var("s")->setVal(0.0); // IMPORTANT
437 b_modelXY.SetSnapshot(*w->set("poi"));
438
439 // create the alternate (signal+background) ModelConfig with s=50
440 ModelConfig sb_modelXY("S+B_modelXY", w);
441 sb_modelXY.SetPdf(*w->pdf("model")); // IMPORTANT
442 sb_modelXY.SetObservables(*w->set("obsXY"));
443 sb_modelXY.SetParametersOfInterest(*w->set("poi"));
444 w->var("s")->setVal(50.0); // IMPORTANT
445 sb_modelXY.SetSnapshot(*w->set("poi"));
446
447 // without this print, their can be a crash when using PROOF. Strange.
448 // w->Print();
449
450 // Test statistics like the profile likelihood ratio
451 // (or the ratio of profiled likelihoods (Tevatron) or the MLE for s)
452 // will now work, since the nuisance parameter b is constrained by y.
453 // ratio of alt and null likelihoods with background yield profiled.
454 //
455 // NOTE: These are slower because they have to run fits for each toy
456
457 // Tevatron-style Ratio of profiled likelihoods
458 // $Q_Tev = -log L(s=0,\hat\hat{b})/L(s=50,\hat\hat{b})$
459 RatioOfProfiledLikelihoodsTestStat ropl(*b_modelXY.GetPdf(), *sb_modelXY.GetPdf(), sb_modelXY.GetSnapshot());
460 ropl.SetSubtractMLE(false);
461
462 // profile likelihood where alternate is best fit value of signal yield
463 // $\lambda(0) = -log L(s=0,\hat\hat{b})/L(\hat{s},\hat{b})$
464 ProfileLikelihoodTestStat profll(*b_modelXY.GetPdf());
465
466 // just use the maximum likelihood estimate of signal yield
467 // $MLE = \hat{s}$
468 MaxLikelihoodEstimateTestStat mlets(*sb_modelXY.GetPdf(), *w->var("s"));
469
470 // However, it is less clear how to justify the prior used in randomizing
471 // the nuisance parameters (since that is a property of the ensemble,
472 // and y is a property of each toy pseudo experiment. In that case,
473 // one probably wants to consider a different y0 which will be held
474 // constant and the prior $\pi(b) = Pois(y0 | tau b) * \eta(b)$.
475 w->factory("y0[100]");
476 w->factory("Gamma::gamma_y0(b,sum::temp0(y0,1),1,0)");
477 w->factory("Gaussian::gauss_prior_y0(b,y0, expr::sqrty0('sqrt(y0)',y0))");
478
479 // HYBRID CALCULATOR
480 HybridCalculator hc3(*dataXY, sb_modelXY, b_modelXY);
481 ToyMCSampler *toymcs3 = (ToyMCSampler *)hc3.GetTestStatSampler();
482 toymcs3->SetNEventsPerToy(1);
483 toymcs3->SetTestStatistic(&slrts);
484 hc3.SetToys(ntoys, ntoys / nToysRatio);
485 hc3.ForcePriorNuisanceAlt(*w->pdf("gamma_y0"));
486 hc3.ForcePriorNuisanceNull(*w->pdf("gamma_y0"));
487 // if you wanted to use the ad hoc Gaussian prior instead
488 // ~~~{.cpp}
489 // hc3.ForcePriorNuisanceAlt(*w->pdf("gauss_prior_y0"));
490 // hc3.ForcePriorNuisanceNull(*w->pdf("gauss_prior_y0"));
491 // ~~~
492
493 // choose fit-based test statistic
494 toymcs3->SetTestStatistic(&profll);
495 // toymcs3->SetTestStatistic(&ropl);
496 // toymcs3->SetTestStatistic(&mlets);
497
498 // enable proof
499 if (pc)
500 toymcs3->SetProofConfig(pc);
501
502 // these lines save current msg level and then kill any messages below ERROR
504 // Get the result
505 HypoTestResult *r3 = hc3.GetHypoTest();
506 cout << "-----------------------------------------" << endl;
507 cout << "Part 6" << endl;
508 r3->Print();
509 t.Stop();
510 t.Print();
511 t.Reset();
512 t.Start();
514
515 c->cd(4);
516 c->GetPad(4)->SetLogy();
517 HypoTestPlot *p3 = new HypoTestPlot(*r3, 50); // 50 bins
518 p3->Draw();
519
520 c->SaveAs("zbi.pdf");
521
522 // -----------------------------------------
523 // OUTPUT W/O PROOF (2.66 GHz Intel Core i7)
524 // =========================================
525
526 // -----------------------------------------
527 // Part 2
528 // Hybrid p-value from direct integration = 0.00094165
529 // Z_Gamma Significance = 3.10804
530 // -----------------------------------------
531 //
532 // Part 3
533 // Z_Bi p-value (analytic): 0.00094165
534 // Z_Bi significance (analytic): 3.10804
535 // Real time 0:00:00, CP time 0.610
536
537 // -----------------------------------------
538 // Part 4
539 // Results HybridCalculator_result:
540 // - Null p-value = 0.00115 +/- 0.000228984
541 // - Significance = 3.04848 sigma
542 // - Number of S+B toys: 1000
543 // - Number of B toys: 20000
544 // - Test statistic evaluated on data: 150
545 // - CL_b: 0.99885 +/- 0.000239654
546 // - CL_s+b: 0.476 +/- 0.0157932
547 // - CL_s: 0.476548 +/- 0.0158118
548 // Real time 0:00:07, CP time 7.620
549
550 // -----------------------------------------
551 // Part 5
552 // Results HybridCalculator_result:
553 // - Null p-value = 0.0009 +/- 0.000206057
554 // - Significance = 3.12139 sigma
555 // - Number of S+B toys: 1000
556 // - Number of B toys: 20000
557 // - Test statistic evaluated on data: 10.8198
558 // - CL_b: 0.9991 +/- 0.000212037
559 // - CL_s+b: 0.465 +/- 0.0157726
560 // - CL_s: 0.465419 +/- 0.0157871
561 // Real time 0:00:34, CP time 34.360
562
563 // -----------------------------------------
564 // Part 6
565 // Results HybridCalculator_result:
566 // - Null p-value = 0.000666667 +/- 0.000149021
567 // - Significance = 3.20871 sigma
568 // - Number of S+B toys: 1000
569 // - Number of B toys: 30000
570 // - Test statistic evaluated on data: 5.03388
571 // - CL_b: 0.999333 +/- 0.000149021
572 // - CL_s+b: 0.511 +/- 0.0158076
573 // - CL_s: 0.511341 +/- 0.0158183
574 // Real time 0:05:06, CP time 306.330
575
576 // ---------------------------------------------------------
577 // OUTPUT w/ PROOF (2.66 GHz Intel Core i7, 4 virtual cores)
578 // =========================================================
579
580 // -----------------------------------------
581 // Part 5
582 // Results HybridCalculator_result:
583 // - Null p-value = 0.00075 +/- 0.000173124
584 // - Significance = 3.17468 sigma
585 // - Number of S+B toys: 1000
586 // - Number of B toys: 20000
587 // - Test statistic evaluated on data: 10.8198
588 // - CL_b: 0.99925 +/- 0.000193577
589 // - CL_s+b: 0.454 +/- 0.0157443
590 // - CL_s: 0.454341 +/- 0.0157564
591 // Real time 0:00:16, CP time 0.990
592
593 // -----------------------------------------
594 // Part 6
595 // Results HybridCalculator_result:
596 // - Null p-value = 0.0007 +/- 0.000152699
597 // - Significance = 3.19465 sigma
598 // - Number of S+B toys: 1000
599 // - Number of B toys: 30000
600 // - Test statistic evaluated on data: 5.03388
601 // - CL_b: 0.9993 +/- 0.000152699
602 // - CL_s+b: 0.518 +/- 0.0158011
603 // - CL_s: 0.518363 +/- 0.0158124
604 // Real time 0:01:25, CP time 0.580
605
606 // ----------------------------------
607 // Comparison
608 // ----------------------------------
609 // LEPStatToolsForLHC
610 // https://plone4.fnal.gov:4430/P0/phystat/packages/0703002
611 // Uses Gaussian prior
612 // CL_b = 6.218476e-04, Significance = 3.228665 sigma
613 //
614 // ----------------------------------
615 // Comparison
616 // ----------------------------------
617 // Asymptotic
618 // From the value of the profile likelihood ratio (5.0338)
619 // The significance can be estimated using Wilks's theorem
620 // significance = sqrt(2*profileLR) = 3.1729 sigma
621}
#define c(i)
Definition RSha256.hxx:101
double Double_t
Definition RtypesCore.h:59
#define ClassDef(name, id)
Definition Rtypes.h:337
#define ClassImp(name)
Definition Rtypes.h:377
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Container class to hold unbinned data.
Definition RooDataSet.h:57
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
static RooMsgService & instance()
Return reference to singleton instance.
void setGlobalKillBelow(RooFit::MsgLevel level)
RooFit::MsgLevel globalKillBelow() const
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:237
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
Same purpose as HybridCalculatorOriginal, but different implementation.
This class provides the plots for the result of a study performed with any of the HypoTestCalculatorG...
HypoTestResult is a base class for results from hypothesis tests.
void Print(const Option_t *="") const override
Print out some information about the results Note: use Alt/Null labels for the hypotheses here as the...
MaxLikelihoodEstimateTestStat: TestStatistic that returns maximum likelihood estimate of a specified ...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
Holds configuration options for proof and proof-lite.
Definition ProofConfig.h:45
TestStatistic that returns the ratio of profiled likelihoods.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
ToyMCSampler is an implementation of the TestStatSampler interface.
void SetProofConfig(ProofConfig *pc=nullptr)
calling with argument or nullptr deactivates proof
virtual void SetTestStatistic(TestStatistic *testStatistic, unsigned int i)
Set the TestStatistic (want the argument to be a function of the data & parameter points.
virtual void SetNEventsPerToy(const Int_t nevents)
Forces the generation of exactly n events even for extended PDFs.
Persistable container for RooFit projects.
The Canvas class.
Definition TCanvas.h:23
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
Definition TPad.cxx:1196
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Stop()
Stop the stopwatch.
void Reset()
Definition TStopwatch.h:52
void Print(Option_t *option="") const override
Print the real and cpu time passed between the start and stop events.
Basic string class.
Definition TString.h:139
RooCmdArg LineColor(Color_t color)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
MsgLevel
Verbosity level for RooMsgService::StreamConfig in RooMsgService.
Namespace for the RooStats classes.
Definition Asimov.h:19
double PValueToSignificance(double pvalue)
returns one-sided significance corresponding to a p-value
Ta Range(0, 0, 1, 1)