Logo ROOT   6.08/07
Reference Guide
HybridCalculator.cxx
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer, Sven Kreiss 23/05/10
3 /*************************************************************************
4  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5  * All rights reserved. *
6  * *
7  * For the licensing terms see $ROOTSYS/LICENSE. *
8  * For the list of contributors see $ROOTSYS/README/CREDITS. *
9  *************************************************************************/
10 
11 /**
12 Same purpose as HybridCalculatorOriginal, but different implementation.
13 */
14 
16 #include "RooStats/ToyMCSampler.h"
17 
18 
20 
21 using namespace RooStats;
22 using namespace std;
23 
24 
25 
26 int HybridCalculator::CheckHook(void) const {
27 
28  if( fPriorNuisanceNull && (!fNullModel->GetNuisanceParameters() || fNullModel->GetNuisanceParameters()->getSize() == 0) ) {
29  oocoutE((TObject*)0,InputArguments) << "HybridCalculator - Nuisance PDF has been specified, but is unaware of which parameters are the nuisance parameters. Must set nuisance parameters in the Null ModelConfig." << endl;
30  return -1; // error
31  }
32  if( fPriorNuisanceAlt && (!fAltModel->GetNuisanceParameters() || fAltModel->GetNuisanceParameters()->getSize() == 0) ) {
33  oocoutE((TObject*)0,InputArguments) << "HybridCalculator - Nuisance PDF has been specified, but is unaware of which parameters are the nuisance parameters. Must set nuisance parameters in the Alt ModelConfig" << endl;
34  return -1; // error
35  }
36 
37  return 0; // ok
38 }
39 
40 
41 int HybridCalculator::PreNullHook(RooArgSet* /*parameterPoint*/, double obsTestStat) const {
42 
43 
44  // ****** any TestStatSampler ********
45 
46  if(fPriorNuisanceNull) {
47  // Setup Priors for ad hoc Hybrid
48  fTestStatSampler->SetPriorNuisance(fPriorNuisanceNull);
49  } else if(
50  fNullModel->GetNuisanceParameters() == NULL ||
51  fNullModel->GetNuisanceParameters()->getSize() == 0
52  ) {
54  << "HybridCalculator - No nuisance parameters specified for Null model and no prior forced. "
55  << "Case is reduced to simple hypothesis testing with no uncertainty." << endl;
56  } else {
57  oocoutI((TObject*)0,InputArguments) << "HybridCalculator - Using uniform prior on nuisance parameters (Null model)." << endl;
58  }
59 
60 
61  // ***** ToyMCSampler specific *******
62 
63  // check whether TestStatSampler is a ToyMCSampler
64  ToyMCSampler *toymcs = dynamic_cast<ToyMCSampler*>(GetTestStatSampler());
65  if(toymcs) {
66  oocoutI((TObject*)0,InputArguments) << "Using a ToyMCSampler. Now configuring for Null." << endl;
67 
68  // variable number of toys
69  if(fNToysNull >= 0) toymcs->SetNToys(fNToysNull);
70 
71  // adaptive sampling
72  if(fNToysNullTail) {
73  oocoutI((TObject*)0,InputArguments) << "Adaptive Sampling" << endl;
74  if(GetTestStatSampler()->GetTestStatistic()->PValueIsRightTail()) {
75  toymcs->SetToysRightTail(fNToysNullTail, obsTestStat);
76  }else{
77  toymcs->SetToysLeftTail(fNToysNullTail, obsTestStat);
78  }
79  }else{
80  toymcs->SetToysBothTails(0, 0, obsTestStat); // disable adaptive sampling
81  }
82 
83  GetNullModel()->LoadSnapshot();
84  }
85 
86  return 0;
87 }
88 
89 
90 int HybridCalculator::PreAltHook(RooArgSet* /*parameterPoint*/, double obsTestStat) const {
91 
92  // ****** any TestStatSampler ********
93 
94  if(fPriorNuisanceAlt){
95  // Setup Priors for ad hoc Hybrid
96  fTestStatSampler->SetPriorNuisance(fPriorNuisanceAlt);
97  } else if (
98  fAltModel->GetNuisanceParameters()==NULL ||
99  fAltModel->GetNuisanceParameters()->getSize()==0
100  ) {
102  << "HybridCalculator - No nuisance parameters specified for Alt model and no prior forced. "
103  << "Case is reduced to simple hypothesis testing with no uncertainty." << endl;
104  } else {
105  oocoutI((TObject*)0,InputArguments) << "HybridCalculator - Using uniform prior on nuisance parameters (Alt model)." << endl;
106  }
107 
108 
109  // ***** ToyMCSampler specific *******
110 
111  // check whether TestStatSampler is a ToyMCSampler
112  ToyMCSampler *toymcs = dynamic_cast<ToyMCSampler*>(GetTestStatSampler());
113  if(toymcs) {
114  oocoutI((TObject*)0,InputArguments) << "Using a ToyMCSampler. Now configuring for Alt." << endl;
115 
116  // variable number of toys
117  if(fNToysAlt >= 0) toymcs->SetNToys(fNToysAlt);
118 
119  // adaptive sampling
120  if(fNToysAltTail) {
121  oocoutI((TObject*)0,InputArguments) << "Adaptive Sampling" << endl;
122  if(GetTestStatSampler()->GetTestStatistic()->PValueIsRightTail()) {
123  toymcs->SetToysLeftTail(fNToysAltTail, obsTestStat);
124  }else{
125  toymcs->SetToysRightTail(fNToysAltTail, obsTestStat);
126  }
127  }else{
128  toymcs->SetToysBothTails(0, 0, obsTestStat); // disable adaptive sampling
129  }
130  }
131 
132  return 0;
133 }
134 
135 
136 
137 
int CheckHook(void) const
check whether all input is consistent
void SetToysLeftTail(Double_t toys, Double_t threshold)
Definition: ToyMCSampler.h:246
#define oocoutI(o, a)
Definition: RooMsgService.h:45
int PreNullHook(RooArgSet *, double obsTestStat) const
configure TestStatSampler for the Null run
STL namespace.
void SetToysRightTail(Double_t toys, Double_t threshold)
Definition: ToyMCSampler.h:251
void SetToysBothTails(Double_t toys, Double_t low_threshold, Double_t high_threshold)
Definition: ToyMCSampler.h:256
#define oocoutE(o, a)
Definition: RooMsgService.h:48
int PreAltHook(RooArgSet *, double obsTestStat) const
configure TestStatSampler for the Alt run
virtual void SetNToys(const Int_t ntoy)
Definition: ToyMCSampler.h:176
ToyMCSampler is an implementation of the TestStatSampler interface.
Definition: ToyMCSampler.h:99
This class implements the Hypothesis test calculation using an hybrid (frequentist/bayesian) procedur...
Namespace for the RooStats classes.
Definition: Asimov.h:20
#define ClassImp(name)
Definition: Rtypes.h:279
Mother of all ROOT objects.
Definition: TObject.h:37
#define NULL
Definition: Rtypes.h:82