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