Logo ROOT  
Reference Guide
SimpleLikelihoodRatioTestStat.h
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer and Sven Kreiss June 2010
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#ifndef ROOSTATS_SimpleLikelihoodRatioTestStat
12#define ROOSTATS_SimpleLikelihoodRatioTestStat
13
14#include "Rtypes.h"
15
16#include "RooAbsPdf.h"
17#include "RooRealVar.h"
18
20
21namespace RooStats {
22
24
25 public:
26
27 /// Constructor for proof. Do not use.
29 fNullPdf(nullptr), fAltPdf(nullptr)
30 {
31 fFirstEval = true;
33 fDetailedOutput = nullptr;
34 fNullParameters = nullptr;
35 fAltParameters = nullptr;
36 fReuseNll=false ;
37 fNllNull=nullptr ;
38 fNllAlt=nullptr ;
39 }
40
41 /// Takes null and alternate parameters from PDF. Can be overridden.
43 RooAbsPdf& nullPdf,
44 RooAbsPdf& altPdf
45 ) :
46 fFirstEval(true)
47 {
48 fNullPdf = &nullPdf;
49 fAltPdf = &altPdf;
50
51 RooArgSet * allNullVars = fNullPdf->getVariables();
52 fNullParameters = (RooArgSet*) allNullVars->snapshot();
53 delete allNullVars;
54
55 RooArgSet * allAltVars = fAltPdf->getVariables();
56 fAltParameters = (RooArgSet*) allAltVars->snapshot();
57 delete allAltVars;
58
60 fDetailedOutput = nullptr;
61
62 fReuseNll=false ;
63 fNllNull=nullptr ;
64 fNllAlt=nullptr ;
65 }
66
67 /// Takes null and alternate parameters from values in nullParameters
68 /// and altParameters. Can be overridden.
70 RooAbsPdf& nullPdf,
71 RooAbsPdf& altPdf,
72 const RooArgSet& nullParameters,
73 const RooArgSet& altParameters
74 ) :
75 fFirstEval(true)
76 {
77 fNullPdf = &nullPdf;
78 fAltPdf = &altPdf;
79
80 fNullParameters = (RooArgSet*) nullParameters.snapshot();
81 fAltParameters = (RooArgSet*) altParameters.snapshot();
82
84 fDetailedOutput = nullptr;
85
86 fReuseNll=false ;
87 fNllNull=nullptr ;
88 fNllAlt=nullptr ;
89 }
90
94 if (fNllNull) delete fNllNull ;
95 if (fNllAlt) delete fNllAlt ;
97 }
98
99 static void SetAlwaysReuseNLL(bool flag);
100
101 void SetReuseNLL(bool flag) { fReuseNll = flag ; }
102
103 void SetNullParameters(const RooArgSet& nullParameters) {
105 fFirstEval = true;
106 fNullParameters = (RooArgSet*) nullParameters.snapshot();
107 }
108
109 void SetAltParameters(const RooArgSet& altParameters) {
110 if (fAltParameters) delete fAltParameters;
111 fFirstEval = true;
112 fAltParameters = (RooArgSet*) altParameters.snapshot();
113 }
114
115 /// this should be possible with RooAbsCollection
117 if (!fNullParameters->equals(*fAltParameters)) return false;
118
119 bool ret = true;
120
121 for (auto nullIt = fNullParameters->begin(), altIt = fAltParameters->begin();
122 nullIt != fNullParameters->end() && altIt != fAltParameters->end(); ++nullIt, ++altIt) {
123 RooAbsReal *null = static_cast<RooAbsReal *>(*nullIt);
124 RooAbsReal *alt = static_cast<RooAbsReal *>(*altIt);
125 if (null->getVal() != alt->getVal())
126 ret = false;
127 }
128
129 return ret;
130 }
131
132
133 /// set the conditional observables which will be used when creating the NLL
134 /// so the pdf's will not be normalized on the conditional observables when computing the NLL
136
137 /// set the global observables which will be used when creating the NLL
138 /// so the constraint pdf's will be normalized correctly on the global observables when computing the NLL
140
141 double Evaluate(RooAbsData& data, RooArgSet& nullPOI) override;
142
143 virtual void EnableDetailedOutput( bool e=true ) { fDetailedOutputEnabled = e; fDetailedOutput = nullptr; }
144 const RooArgSet* GetDetailedOutput(void) const override { return fDetailedOutput; }
145
146 const TString GetVarName() const override {
147 return "log(L(#mu_{1}) / L(#mu_{0}))";
148 }
149
150 private:
151
159
162
163 RooAbsReal* fNllNull ; ///<! transient copy of the null NLL
164 RooAbsReal* fNllAlt ; ///<! transient copy of the alt NLL
165 static bool fgAlwaysReuseNll ;
167
168
169 protected:
171};
172
173}
174
175#endif
#define e(i)
Definition: RSha256.hxx:103
#define ClassDefOverride(name, id)
Definition: Rtypes.h:339
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
RooArgSet * getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2079
bool equals(const RooAbsCollection &otherColl) const
Check if this and other collection have identically-named contents.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
const_iterator end() const
const_iterator begin() const
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:61
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:61
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:105
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:56
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition: RooArgSet.h:179
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
RooAbsReal * fNllNull
! transient copy of the null NLL
double Evaluate(RooAbsData &data, RooArgSet &nullPOI) override
Main interface to evaluate the test statistic on a dataset given the values for the Null Parameters O...
SimpleLikelihoodRatioTestStat(RooAbsPdf &nullPdf, RooAbsPdf &altPdf, const RooArgSet &nullParameters, const RooArgSet &altParameters)
Takes null and alternate parameters from values in nullParameters and altParameters.
void SetConditionalObservables(const RooArgSet &set) override
set the conditional observables which will be used when creating the NLL so the pdf's will not be nor...
SimpleLikelihoodRatioTestStat()
Constructor for proof. Do not use.
void SetNullParameters(const RooArgSet &nullParameters)
SimpleLikelihoodRatioTestStat(RooAbsPdf &nullPdf, RooAbsPdf &altPdf)
Takes null and alternate parameters from PDF. Can be overridden.
RooAbsReal * fNllAlt
! transient copy of the alt NLL
bool ParamsAreEqual()
this should be possible with RooAbsCollection
void SetGlobalObservables(const RooArgSet &set) override
set the global observables which will be used when creating the NLL so the constraint pdf's will be n...
void SetAltParameters(const RooArgSet &altParameters)
const RooArgSet * GetDetailedOutput(void) const override
return detailed output: for fits this can be pulls, processing time, ... The returned pointer will no...
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
Definition: TestStatistic.h:31
Basic string class.
Definition: TString.h:136
Namespace for the RooStats classes.
Definition: Asimov.h:19
null_t< F > null()