ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HybridInstructional.C
Go to the documentation of this file.
1 //+ example demostrating usage of HybridCalcultor
2 /*
3 HybridInstructional
4 
5 Authors: Kyle Cranmer, Wouter Verkerke, and Sven Kreiss
6 date May 2010 Part 1-3
7 date Dec 2010 Part 4-6
8 
9 A hypothesis testing example based on number counting
10 with background uncertainty.
11 
12 NOTE: This example must be run with the ACLIC (the + option ) due to the
13 new class that is defined.
14 
15 This example:
16  - demonstrates the usage of the HybridCalcultor (Part 4-6)
17  - demonstrates the numerical integration of RooFit (Part 2)
18  - validates the RooStats against an example with a known analytic answer
19  - demonstrates usage of different test statistics
20  - explains subtle choices in the prior used for hybrid methods
21  - demonstrates usage of different priors for the nuisance parameters
22  - demonstrates usage of PROOF
23 
24 The basic setup here is that a main measurement has observed x events with an
25 expectation of s+b. One can choose an ad hoc prior for the uncertainty on b,
26 or try to base it on an auxiliary measurement. In this case, the auxiliary
27 measurement (aka control measurement, sideband) is another counting experiment
28 with measurement y and expectation tau*b. With an 'original prior' on b,
29 called \eta(b) then one can obtain a posterior from the auxiliary measurement
30 \pi(b) = \eta(b) * Pois(y|tau*b). This is a principled choice for a prior
31 on b in the main measurement of x, which can then be treated in a hybrid
32 Bayesian/Frequentist way. Additionally, one can try to treat the two
33 measurements simultaneously, which is detailed in Part 6 of the tutorial.
34 
35 This tutorial is related to the FourBin.C tutorial in the modeling, but
36 focuses on hypothesis testing instead of interval estimation.
37 
38 More background on this 'prototype problem' can be found in the
39 following papers:
40 
41 Evaluation of three methods for calculating statistical significance
42 when incorporating a systematic uncertainty into a test of the
43 background-only hypothesis for a Poisson process
44 Authors: Robert D. Cousins, James T. Linnemann, Jordan Tucker
45 http://arxiv.org/abs/physics/0702156
46 NIM A 595 (2008) 480--501
47 
48 Statistical Challenges for Searches for New Physics at the LHC
49 Authors: Kyle Cranmer
50 http://arxiv.org/abs/physics/0511028
51 
52  Measures of Significance in HEP and Astrophysics
53  Authors: J. T. Linnemann
54  http://arxiv.org/abs/physics/0312059
55 */
56 
57 #include "RooGlobalFunc.h"
58 #include "RooRealVar.h"
59 #include "RooProdPdf.h"
60 #include "RooWorkspace.h"
61 #include "RooDataSet.h"
62 #include "TCanvas.h"
63 #include "TStopwatch.h"
64 #include "TH1.h"
65 #include "RooPlot.h"
66 #include "RooMsgService.h"
67 
69 
71 #include "RooStats/ToyMCSampler.h"
72 #include "RooStats/HypoTestPlot.h"
73 
78 
79 using namespace RooFit;
80 using namespace RooStats;
81 
82 //////////////////////////////////////////////////
83 // A New Test Statistic Class for this example.
84 // It simply returns the sum of the values in a particular
85 // column of a dataset.
86 // You can ignore this class and focus on the macro below
87 //////////////////////////////////////////////////
88 class BinCountTestStat : public TestStatistic {
89 public:
90  BinCountTestStat(void) : fColumnName("tmp") {}
91  BinCountTestStat(string columnName) : fColumnName(columnName) {}
92 
93  virtual Double_t Evaluate(RooAbsData& data, RooArgSet& /*nullPOI*/) {
94  // This is the main method in the interface
95  Double_t value = 0.0;
96  for(int i=0; i < data.numEntries(); i++) {
97  value += data.get(i)->getRealValue(fColumnName.c_str());
98  }
99  return value;
100  }
101  virtual const TString GetVarName() const { return fColumnName; }
102 
103 private:
104  string fColumnName;
105 
106 protected:
107  ClassDef(BinCountTestStat,1)
108 };
109 
110 ClassImp(BinCountTestStat)
111 
112 //////////////////////////////////////////////////
113 // The Actual Tutorial Macro
114 //////////////////////////////////////////////////
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 #ifdef __CINT__
134  cout << "DO NOT RUN WITH CINT: we are using a custom test statistic ";
135  cout << "which requires that this tutorial must be compiled ";
136  cout << "with ACLIC" << endl;
137  return;
138 #endif
139 
140 
141  TStopwatch t;
142  t.Start();
143  TCanvas *c = new TCanvas;
144  c->Divide(2,2);
145 
146  ///////////////////////////////////////////////////////
147  // P A R T 1 : D I R E C T I N T E G R A T I O N
148  //////////////////////////////////////////////////////
149  // Make model for prototype on/off problem
150  // Pois(x | s+b) * Pois(y | tau b )
151  // for Z_Gamma, use uniform prior on b.
152  RooWorkspace* w = new RooWorkspace("w");
153  w->factory("Poisson::px(x[150,0,500],sum::splusb(s[0,0,100],b[100,0,300]))");
154  w->factory("Poisson::py(y[100,0,500],prod::taub(tau[1.],b))");
155  w->factory("PROD::model(px,py)");
156  w->factory("Uniform::prior_b(b)");
157 
158  // We will control the output level in a few places to avoid
159  // verbose progress messages. We start by keeping track
160  // of the current threshold on messages.
162 
163  // Use PROOF-lite on multi-core machines
164  ProofConfig* pc = NULL;
165  // uncomment below if you want to use PROOF
166  // pc = new ProofConfig(*w, 4, "workers=4", kFALSE); // machine with 4 cores
167  // pc = new ProofConfig(*w, 2, "workers=2", kFALSE); // machine with 2 cores
168 
169  ///////////////////////////////////////////////////////
170  // P A R T 2 : D I R E C T I N T E G R A T I O N
171  //////////////////////////////////////////////////////
172  // This is not the 'RooStats' way, but in this case the distribution
173  // of the test statistic is simply x and can be calculated directly
174  // from the PDF using RooFit's built-in integration.
175  // Note, this does not generalize to situations in which the test statistic
176  // depends on many events (rows in a dataset).
177 
178  // construct the Bayesian-averaged model (eg. a projection pdf)
179  // p'(x|s) = \int db p(x|s+b) * [ p(y|b) * prior(b) ]
180  w->factory("PROJ::averagedModel(PROD::foo(px|b,py,prior_b),b)") ;
181 
182  RooMsgService::instance().setGlobalKillBelow(RooFit::ERROR); // lower message level
183  // plot it, red is averaged model, green is b known exactly, blue is s+b av model
184  RooPlot* frame = w->var("x")->frame(Range(50,230)) ;
185  w->pdf("averagedModel")->plotOn(frame,LineColor(kRed)) ;
186  w->pdf("px")->plotOn(frame,LineColor(kGreen)) ;
187  w->var("s")->setVal(50.);
188  w->pdf("averagedModel")->plotOn(frame,LineColor(kBlue)) ;
189  c->cd(1);
190  frame->Draw() ;
191  w->var("s")->setVal(0.);
192 
193  // compare analytic calculation of Z_Bi
194  // with the numerical RooFit implementation of Z_Gamma
195  // for an example with x = 150, y = 100
196 
197  // numeric RooFit Z_Gamma
198  w->var("y")->setVal(100);
199  w->var("x")->setVal(150);
200  RooAbsReal* cdf = w->pdf("averagedModel")->createCdf(*w->var("x"));
201  cdf->getVal(); // get ugly print messages out of the way
202  cout << "-----------------------------------------"<<endl;
203  cout << "Part 2" << endl;
204  cout << "Hybrid p-value from direct integration = " << 1-cdf->getVal() << endl;
205  cout << "Z_Gamma Significance = " <<
206  PValueToSignificance(1-cdf->getVal()) << endl;
207  RooMsgService::instance().setGlobalKillBelow(msglevel); // set it back
208 
209  /////////////////////////////////////////////////
210  // P A R T 3 : A N A L Y T I C R E S U L T
211  /////////////////////////////////////////////////
212  // In this special case, the integrals are known analytically
213  // and they are implemented in RooStats::NumberCountingUtils
214 
215  // analytic Z_Bi
216  double p_Bi = NumberCountingUtils::BinomialWithTauObsP(150, 100, 1);
217  double Z_Bi = NumberCountingUtils::BinomialWithTauObsZ(150, 100, 1);
218  cout << "-----------------------------------------"<<endl;
219  cout << "Part 3" << endl;
220  std::cout << "Z_Bi p-value (analytic): " << p_Bi << std::endl;
221  std::cout << "Z_Bi significance (analytic): " << Z_Bi << std::endl;
222  t.Stop(); t.Print(); t.Reset(); t.Start();
223 
224  ////////////////////////////////////////////////////////////////
225  // 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
226  ////////////////////////////////////////////////////////////////
227  // Now we demonstrate the RooStats HybridCalculator.
228  //
229  // Like all RooStats calculators it needs the data and a ModelConfig
230  // for the relevant hypotheses. Since we are doing hypothesis testing
231  // we need a ModelConfig for the null (background only) and the alternate
232  // (signal+background) hypotheses. We also need to specify the PDF,
233  // the parameters of interest, and the observables. Furthermore, since
234  // the parameter of interest is floating, we need to specify which values
235  // of the parameter corresponds to the null and alternate (eg. s=0 and s=50)
236  //
237  // define some sets of variables obs={x} and poi={s}
238  // note here, x is the only observable in the main measurement
239  // and y is treated as a separate measurement, which is used
240  // to produce the prior that will be used in this calculation
241  // to randomize the nuisance parameters.
242  w->defineSet("obs","x");
243  w->defineSet("poi","s");
244 
245  // create a toy dataset with the x=150
246  RooDataSet *data = new RooDataSet("d", "d", *w->set("obs"));
247  data->add(*w->set("obs"));
248 
249  //////////////////////////////////////////////////////////
250  // Part 3a : Setup ModelConfigs
251  // create the null (background-only) ModelConfig with s=0
252  ModelConfig b_model("B_model", w);
253  b_model.SetPdf(*w->pdf("px"));
254  b_model.SetObservables(*w->set("obs"));
255  b_model.SetParametersOfInterest(*w->set("poi"));
256  w->var("s")->setVal(0.0); // important!
257  b_model.SetSnapshot(*w->set("poi"));
258 
259  // create the alternate (signal+background) ModelConfig with s=50
260  ModelConfig sb_model("S+B_model", w);
261  sb_model.SetPdf(*w->pdf("px"));
262  sb_model.SetObservables(*w->set("obs"));
263  sb_model.SetParametersOfInterest(*w->set("poi"));
264  w->var("s")->setVal(50.0); // important!
265  sb_model.SetSnapshot(*w->set("poi"));
266 
267 
268  //////////////////////////////////////////////////////////
269  // Part 3b : Choose Test Statistic
270  // To make an equivalent calculation we need to use x as the test
271  // statistic. This is not a built-in test statistic in RooStats
272  // so we define it above. The new class inherits from the
273  // RooStats::TestStatistic interface, and simply returns the value
274  // of x in the dataset.
275 
276  BinCountTestStat binCount("x");
277 
278  //////////////////////////////////////////////////////////
279  // Part 3c : Define Prior used to randomize nuisance parameters
280  //
281  // The prior used for the hybrid calculator is the posterior
282  // from the auxiliary measurement y. The model for the aux.
283  // measurement is Pois(y|tau*b), thus the likleihood function
284  // is proportional to (has the form of) a Gamma distribution.
285  // if the 'original prior' \eta(b) is uniform, then from
286  // Bayes's theorem we have the posterior:
287  // \pi(b) = Pois(y|tau*b) * \eta(b)
288  // If \eta(b) is flat, then we arrive at a Gamma distribution.
289  // Since RooFit will normalize the PDF we can actually supply
290  // py=Pois(y,tau*b) that will be equivalent to multiplying by a uniform.
291  //
292  // Alternatively, we could explicitly use a gamma distribution:
293  // w->factory("Gamma::gamma(b,sum::temp(y,1),1,0)");
294  //
295  // or we can use some other ad hoc prior that do not naturally
296  // follow from the known form of the auxiliary measurement.
297  // The common choice is the equivlaent Gaussian:
298  w->factory("Gaussian::gauss_prior(b,y, expr::sqrty('sqrt(y)',y))");
299  // this corresponds to the "Z_N" calculation.
300  //
301  // or one could use the analogous log-normal prior
302  w->factory("Lognormal::lognorm_prior(b,y, expr::kappa('1+1./sqrt(y)',y))");
303  //
304  // Ideally, the HybridCalculator would be able to inspect the full
305  // model Pois(x | s+b) * Pois(y | tau b ) and be given the original
306  // prior \eta(b) to form \pi(b) = Pois(y|tau*b) * \eta(b).
307  // This is not yet implemented because in the general case
308  // it is not easy to identify the terms in the PDF that correspond
309  // to the auxiliary measurement. So for now, it must be set
310  // explicitly with:
311  // - ForcePriorNuisanceNull()
312  // - ForcePriorNuisanceAlt()
313  // the name "ForcePriorNuisance" was chosen because we anticipate
314  // this to be auto-detected, but will leave the option open
315  // to force to a different prior for the nuisance parameters.
316 
317  //////////////////////////////////////////////////////////
318  // Part 3d : Construct and configure the HybridCalculator
319 
320  HybridCalculator hc1(*data, sb_model, b_model);
321  ToyMCSampler *toymcs1 = (ToyMCSampler*)hc1.GetTestStatSampler();
322  toymcs1->SetNEventsPerToy(1); // because the model is in number counting form
323  toymcs1->SetTestStatistic(&binCount); // set the test statistic
324  hc1.SetToys(20000,1000);
325  hc1.ForcePriorNuisanceAlt(*w->pdf("py"));
326  hc1.ForcePriorNuisanceNull(*w->pdf("py"));
327  // if you wanted to use the ad hoc Gaussian prior instead
328  // hc1.ForcePriorNuisanceAlt(*w->pdf("gauss_prior"));
329  // hc1.ForcePriorNuisanceNull(*w->pdf("gauss_prior"));
330  // if you wanted to use the ad hoc log-normal prior instead
331  // hc1.ForcePriorNuisanceAlt(*w->pdf("lognorm_prior"));
332  // hc1.ForcePriorNuisanceNull(*w->pdf("lognorm_prior"));
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 W I T H
355  // 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  //
358  // A likelihood ratio test statistics should be 1-to-1 with the count x
359  // when the value of b is fixed in the likelihood. This is implemented
360  // by the SimpleLikelihoodRatioTestStat
361 
362  SimpleLikelihoodRatioTestStat slrts(*b_model.GetPdf(),*sb_model.GetPdf());
363  slrts.SetNullParameters(*b_model.GetSnapshot());
364  slrts.SetAltParameters(*sb_model.GetSnapshot());
365 
366  // HYBRID CALCULATOR
367  HybridCalculator hc2(*data, sb_model, b_model);
368  ToyMCSampler *toymcs2 = (ToyMCSampler*)hc2.GetTestStatSampler();
369  toymcs2->SetNEventsPerToy(1);
370  toymcs2->SetTestStatistic(&slrts);
371  hc2.SetToys(20000,1000);
372  hc2.ForcePriorNuisanceAlt(*w->pdf("py"));
373  hc2.ForcePriorNuisanceNull(*w->pdf("py"));
374  // if you wanted to use the ad hoc Gaussian prior instead
375  // hc2.ForcePriorNuisanceAlt(*w->pdf("gauss_prior"));
376  // hc2.ForcePriorNuisanceNull(*w->pdf("gauss_prior"));
377  // if you wanted to use the ad hoc log-normal prior instead
378  // hc2.ForcePriorNuisanceAlt(*w->pdf("lognorm_prior"));
379  // hc2.ForcePriorNuisanceNull(*w->pdf("lognorm_prior"));
380 
381  // enable proof
382  if(pc) toymcs2->SetProofConfig(pc);
383 
384  // these lines save current msg level and then kill any messages below ERROR
386  // Get the result
387  HypoTestResult *r2 = hc2.GetHypoTest();
388  cout << "-----------------------------------------"<<endl;
389  cout << "Part 5" << endl;
390  r2->Print();
391  t.Stop(); t.Print(); t.Reset(); t.Start();
393 
394  c->cd(3);
395  HypoTestPlot *p2 = new HypoTestPlot(*r2,30); // 30 bins
396  p2->Draw();
397 
398  ////////////////////////////////////////////////////////////////////////////
399  // 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
400  // 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
401  // A N D S I M U L T A N E O U S M O D E L
402  /////////////////////////////////////////////////////////////////////////////
403  //
404  // If one wants to use a test statistic in which the nuisance parameters
405  // are profiled (in one way or another), then the PDF must constrain b.
406  // Otherwise any observation x can always be explained with s=0 and b=x/tau.
407  //
408  // In this case, one is really thinking about the problem in a
409  // different way. They are considering x,y simultaneously.
410  // and the PDF should be Pois(x | s+b) * Pois(y | tau b )
411  // and the set 'obs' should be {x,y}.
412 
413  w->defineSet("obsXY","x,y");
414 
415  // create a toy dataset with the x=150, y=100
416  w->var("x")->setVal(150.);
417  w->var("y")->setVal(100.);
418  RooDataSet *dataXY = new RooDataSet("dXY", "dXY", *w->set("obsXY"));
419  dataXY->add(*w->set("obsXY"));
420 
421  // now we need new model configs, with PDF="model"
422  ModelConfig b_modelXY("B_modelXY", w);
423  b_modelXY.SetPdf(*w->pdf("model")); // IMPORTANT
424  b_modelXY.SetObservables(*w->set("obsXY"));
425  b_modelXY.SetParametersOfInterest(*w->set("poi"));
426  w->var("s")->setVal(0.0); // IMPORTANT
427  b_modelXY.SetSnapshot(*w->set("poi"));
428 
429  // create the alternate (signal+background) ModelConfig with s=50
430  ModelConfig sb_modelXY("S+B_modelXY", w);
431  sb_modelXY.SetPdf(*w->pdf("model")); // IMPORTANT
432  sb_modelXY.SetObservables(*w->set("obsXY"));
433  sb_modelXY.SetParametersOfInterest(*w->set("poi"));
434  w->var("s")->setVal(50.0); // IMPORTANT
435  sb_modelXY.SetSnapshot(*w->set("poi"));
436 
437  // without this print, their can be a crash when using PROOF. Strange.
438  // w->Print();
439 
440  // Test statistics like the profile likelihood ratio
441  // (or the ratio of profiled likelihoods (Tevatron) or the MLE for s)
442  // will now work, since the nuisance parameter b is constrained by y.
443  // ratio of alt and null likelihoods with background yield profiled.
444  //
445  // NOTE: These are slower because they have to run fits for each toy
446 
447  // Tevatron-style Ratio of profiled likelihoods
448  // Q_Tev = -log L(s=0,\hat\hat{b})/L(s=50,\hat\hat{b})
450  ropl(*b_modelXY.GetPdf(), *sb_modelXY.GetPdf(), sb_modelXY.GetSnapshot());
451  ropl.SetSubtractMLE(false);
452 
453  // profile likelihood where alternate is best fit value of signal yield
454  // \lambda(0) = -log L(s=0,\hat\hat{b})/L(\hat{s},\hat{b})
455  ProfileLikelihoodTestStat profll(*b_modelXY.GetPdf());
456 
457  // just use the maximum likelihood estimate of signal yield
458  // MLE = \hat{s}
459  MaxLikelihoodEstimateTestStat mlets(*sb_modelXY.GetPdf(), *w->var("s"));
460 
461  // However, it is less clear how to justify the prior used in randomizing
462  // the nuisance parameters (since that is a property of the ensemble,
463  // and y is a property of each toy pseudo experiment. In that case,
464  // one probably wants to consider a different y0 which will be held
465  // constant and the prior \pi(b) = Pois(y0 | tau b) * \eta(b).
466  w->factory("y0[100]");
467  w->factory("Gamma::gamma_y0(b,sum::temp0(y0,1),1,0)");
468  w->factory("Gaussian::gauss_prior_y0(b,y0, expr::sqrty0('sqrt(y0)',y0))");
469 
470 
471  // HYBRID CALCULATOR
472  HybridCalculator hc3(*dataXY, sb_modelXY, b_modelXY);
473  ToyMCSampler *toymcs3 = (ToyMCSampler*)hc3.GetTestStatSampler();
474  toymcs3->SetNEventsPerToy(1);
475  toymcs3->SetTestStatistic(&slrts);
476  hc3.SetToys(30000,1000);
477  hc3.ForcePriorNuisanceAlt(*w->pdf("gamma_y0"));
478  hc3.ForcePriorNuisanceNull(*w->pdf("gamma_y0"));
479  // if you wanted to use the ad hoc Gaussian prior instead
480  // hc3.ForcePriorNuisanceAlt(*w->pdf("gauss_prior_y0"));
481  // hc3.ForcePriorNuisanceNull(*w->pdf("gauss_prior_y0"));
482 
483  // choose fit-based test statistic
484  toymcs3->SetTestStatistic(&profll);
485  //toymcs3->SetTestStatistic(&ropl);
486  //toymcs3->SetTestStatistic(&mlets);
487 
488  // enable proof
489  if(pc) toymcs3->SetProofConfig(pc);
490 
491  // these lines save current msg level and then kill any messages below ERROR
493  // Get the result
494  HypoTestResult *r3 = hc3.GetHypoTest();
495  cout << "-----------------------------------------"<<endl;
496  cout << "Part 6" << endl;
497  r3->Print();
498  t.Stop(); t.Print(); t.Reset(); t.Start();
500 
501  c->cd(4);
502  c->GetPad(4)->SetLogy();
503  HypoTestPlot *p3 = new HypoTestPlot(*r3,50); // 50 bins
504  p3->Draw();
505 
506  c->SaveAs("zbi.pdf");
507 
508 
509  ///////////////////////////////////////////////////////////
510  // OUTPUT W/O PROOF (2.66 GHz Intel Core i7)
511  ///////////////////////////////////////////////////////////
512 
513  /*
514 -----------------------------------------
515 Part 2
516 Hybrid p-value from direct integration = 0.00094165
517 Z_Gamma Significance = 3.10804
518 -----------------------------------------
519 Part 3
520 Z_Bi p-value (analytic): 0.00094165
521 Z_Bi significance (analytic): 3.10804
522 Real time 0:00:00, CP time 0.610
523 
524 -----------------------------------------
525 Part 4
526 Results HybridCalculator_result:
527  - Null p-value = 0.00115 +/- 0.000228984
528  - Significance = 3.04848 sigma
529  - Number of S+B toys: 1000
530  - Number of B toys: 20000
531  - Test statistic evaluated on data: 150
532  - CL_b: 0.99885 +/- 0.000239654
533  - CL_s+b: 0.476 +/- 0.0157932
534  - CL_s: 0.476548 +/- 0.0158118
535 Real time 0:00:07, CP time 7.620
536 
537 -----------------------------------------
538 Part 5
539 Results HybridCalculator_result:
540  - Null p-value = 0.0009 +/- 0.000206057
541  - Significance = 3.12139 sigma
542  - Number of S+B toys: 1000
543  - Number of B toys: 20000
544  - Test statistic evaluated on data: 10.8198
545  - CL_b: 0.9991 +/- 0.000212037
546  - CL_s+b: 0.465 +/- 0.0157726
547  - CL_s: 0.465419 +/- 0.0157871
548 Real time 0:00:34, CP time 34.360
549 
550 -----------------------------------------
551 Part 6
552 Results HybridCalculator_result:
553  - Null p-value = 0.000666667 +/- 0.000149021
554  - Significance = 3.20871 sigma
555  - Number of S+B toys: 1000
556  - Number of B toys: 30000
557  - Test statistic evaluated on data: 5.03388
558  - CL_b: 0.999333 +/- 0.000149021
559  - CL_s+b: 0.511 +/- 0.0158076
560  - CL_s: 0.511341 +/- 0.0158183
561 Real time 0:05:06, CP time 306.330
562 
563  */
564 
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  //////////////////////////////////////////
600  // Comparison
601  ///////////////////////////////////////////
602  // LEPStatToolsForLHC
603  // https://plone4.fnal.gov:4430/P0/phystat/packages/0703002
604  // Uses Gaussian prior
605  // CL_b = 6.218476e-04, Significance = 3.228665 sigma
606  //
607  //////////////////////////////////////////
608  // Comparison
609  ///////////////////////////////////////////
610  // Asymptotics
611  // From the value of the profile likelihood ratio (5.0338)
612  // The significance can be estimated using Wilks's theorem
613  // significance = sqrt(2*profileLR) = 3.1729 sigma
614 
615 
616 }
Holds configuration options for proof and proof-lite.
Definition: ProofConfig.h:49
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:217
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:52
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:244
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:56
static double p3(double t, double a, double b, double c, double d)
RooCmdArg LineColor(Color_t color)
return c
Definition: Rtypes.h:61
Double_t PValueToSignificance(Double_t pvalue)
Definition: RooStatsUtils.h:59
HypoTestResult is a base class for results from hypothesis tests.
virtual void SetNEventsPerToy(const Int_t nevents)
Definition: ToyMCSampler.h:177
Basic string class.
Definition: TString.h:137
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 void SetObservables(const RooArgSet &set)
specify the observables
Definition: ModelConfig.h:156
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:75
#define ClassDef(name, id)
Definition: Rtypes.h:254
virtual HypoTestResult * GetHypoTest() const
inherited methods from HypoTestCalculator interface
TestStatSampler * GetTestStatSampler(void) const
double cdf(double *x, double *p)
Definition: unuranDistr.cxx:44
void SetNullParameters(const RooArgSet &nullParameters)
unsigned int r3[N_CITIES]
Definition: simanTSP.cxx:323
void Print(const Option_t *="") const
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
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.
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:97
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual const RooArgSet * get() const
Definition: RooAbsData.h:77
Double_t BinomialWithTauObsZ(Double_t nObs, Double_t bExp, Double_t tau)
RooFit::MsgLevel globalKillBelow() const
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:203
TThread * t[5]
Definition: threadsh1.C:13
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
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.
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:293
void SetProofConfig(ProofConfig *pc=NULL)
Definition: ToyMCSampler.h:263
tuple w
Definition: qtexample.py:51
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
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
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 'data' argset, to the data set...
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
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
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return NULL if not existing) ...
void SetToys(int toysNull, int toysAlt)
set number of toys
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
MaxLikelihoodEstimateTestStat: TestStatistic that returns maximum likelihood estimate of a specified ...
RooFactoryWSTool & factory()
Return instance to factory tool.
typedef void((*Func_t)())
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
virtual void SetParametersOfInterest(const RooArgSet &set)
Definition: ModelConfig.h:115
Definition: Rtypes.h:61
ClassImp(BinCountTestStat)
float value
Definition: math.cpp:443
virtual void SetSnapshot(const RooArgSet &set)
set parameter values for a particular hypothesis if using a common PDF by saving a snapshot in the wo...
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
virtual void ForcePriorNuisanceAlt(RooAbsPdf &priorNuisance)
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
virtual void ForcePriorNuisanceNull(RooAbsPdf &priorNuisance)
Override the distribution used for marginalizing nuisance parameters that is inferred from ModelConfi...
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559
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
Stopwatch class.
Definition: TStopwatch.h:30