Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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/** \class RooStats::HybridCalculator
12 \ingroup Roostats
13
14Same purpose as HybridCalculatorOriginal, but different implementation.
15
16This class implements the Hypothesis test calculation using an hybrid
17(frequentist/bayesian) procedure.A frequentist sampling of the test statistic
18distribution is obtained but with marginalization of the nuisance parameters.
19The toys are generated by sampling the nuisance parameters according to their
20prior distribution.
21
22The use of the of ToyMCSampler as the TestStatSampler is assumed.
23
24*/
25
28
29
31
32using namespace RooStats;
33using namespace std;
34
35////////////////////////////////////////////////////////////////////////////////
36
38
40 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;
41 return -1; // error
42 }
44 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;
45 return -1; // error
46 }
47
48 return 0; // ok
49}
50
51////////////////////////////////////////////////////////////////////////////////
52
53int HybridCalculator::PreNullHook(RooArgSet* /*parameterPoint*/, double obsTestStat) const {
54
55
56 // ****** any TestStatSampler ********
57
59 // Setup Priors for ad hoc Hybrid
61 } else if(
64 ) {
65 oocoutI((TObject*)0,InputArguments)
66 << "HybridCalculator - No nuisance parameters specified for Null model and no prior forced. "
67 << "Case is reduced to simple hypothesis testing with no uncertainty." << endl;
68 } else {
69 oocoutI((TObject*)0,InputArguments) << "HybridCalculator - Using uniform prior on nuisance parameters (Null model)." << endl;
70 }
71
72
73 // ***** ToyMCSampler specific *******
74
75 // check whether TestStatSampler is a ToyMCSampler
76 ToyMCSampler *toymcs = dynamic_cast<ToyMCSampler*>(GetTestStatSampler());
77 if(toymcs) {
78 oocoutI((TObject*)0,InputArguments) << "Using a ToyMCSampler. Now configuring for Null." << endl;
79
80 // variable number of toys
81 if(fNToysNull >= 0) toymcs->SetNToys(fNToysNull);
82
83 // adaptive sampling
84 if(fNToysNullTail) {
85 oocoutI((TObject*)0,InputArguments) << "Adaptive Sampling" << endl;
86 if(GetTestStatSampler()->GetTestStatistic()->PValueIsRightTail()) {
87 toymcs->SetToysRightTail(fNToysNullTail, obsTestStat);
88 }else{
89 toymcs->SetToysLeftTail(fNToysNullTail, obsTestStat);
90 }
91 }else{
92 toymcs->SetToysBothTails(0, 0, obsTestStat); // disable adaptive sampling
93 }
94
96 }
97
98 return 0;
99}
100
101////////////////////////////////////////////////////////////////////////////////
102
103int HybridCalculator::PreAltHook(RooArgSet* /*parameterPoint*/, double obsTestStat) const {
104
105 // ****** any TestStatSampler ********
106
108 // Setup Priors for ad hoc Hybrid
110 } else if (
113 ) {
114 oocoutI((TObject*)0,InputArguments)
115 << "HybridCalculator - No nuisance parameters specified for Alt model and no prior forced. "
116 << "Case is reduced to simple hypothesis testing with no uncertainty." << endl;
117 } else {
118 oocoutI((TObject*)0,InputArguments) << "HybridCalculator - Using uniform prior on nuisance parameters (Alt model)." << endl;
119 }
120
121
122 // ***** ToyMCSampler specific *******
123
124 // check whether TestStatSampler is a ToyMCSampler
125 ToyMCSampler *toymcs = dynamic_cast<ToyMCSampler*>(GetTestStatSampler());
126 if(toymcs) {
127 oocoutI((TObject*)0,InputArguments) << "Using a ToyMCSampler. Now configuring for Alt." << endl;
128
129 // variable number of toys
130 if(fNToysAlt >= 0) toymcs->SetNToys(fNToysAlt);
131
132 // adaptive sampling
133 if(fNToysAltTail) {
134 oocoutI((TObject*)0,InputArguments) << "Adaptive Sampling" << endl;
135 if(GetTestStatSampler()->GetTestStatistic()->PValueIsRightTail()) {
136 toymcs->SetToysLeftTail(fNToysAltTail, obsTestStat);
137 }else{
138 toymcs->SetToysRightTail(fNToysAltTail, obsTestStat);
139 }
140 }else{
141 toymcs->SetToysBothTails(0, 0, obsTestStat); // disable adaptive sampling
142 }
143 }
144
145 return 0;
146}
#define oocoutE(o, a)
#define oocoutI(o, a)
#define ClassImp(name)
Definition Rtypes.h:364
Int_t getSize() const
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
Same purpose as HybridCalculatorOriginal, but different implementation.
int PreAltHook(RooArgSet *, double obsTestStat) const
configure TestStatSampler for the Alt run
int PreNullHook(RooArgSet *, double obsTestStat) const
configure TestStatSampler for the Null run
int CheckHook(void) const
check whether all input is consistent
const ModelConfig * GetNullModel(void) const
TestStatSampler * GetTestStatSampler(void) const
Returns instance of TestStatSampler.
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
void LoadSnapshot() const
load the snapshot from ws if it exists
virtual void SetPriorNuisance(RooAbsPdf *)=0
ToyMCSampler is an implementation of the TestStatSampler interface.
virtual void SetNToys(const Int_t ntoy)
void SetToysBothTails(Double_t toys, Double_t low_threshold, Double_t high_threshold)
void SetToysLeftTail(Double_t toys, Double_t threshold)
void SetToysRightTail(Double_t toys, Double_t threshold)
Mother of all ROOT objects.
Definition TObject.h:41
Namespace for the RooStats classes.
Definition Asimov.h:19