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