Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
30 /// Takes null and alternate parameters from PDF. Can be overridden.
32 : fNullPdf(&nullPdf),
33 fAltPdf(&altPdf)
34 {
35 std::unique_ptr<RooArgSet> allNullVars{fNullPdf->getVariables()};
36 fNullParameters = allNullVars->snapshot();
37
38 std::unique_ptr<RooArgSet> allAltVars{fAltPdf->getVariables()};
39 fAltParameters = allAltVars->snapshot();
40 }
41
42 /// Takes null and alternate parameters from values in nullParameters
43 /// and altParameters. Can be overridden.
44 SimpleLikelihoodRatioTestStat(RooAbsPdf &nullPdf, RooAbsPdf &altPdf, const RooArgSet &nullParameters,
45 const RooArgSet &altParameters)
46 : fNullPdf(&nullPdf),
47 fAltPdf(&altPdf),
48 fNullParameters(nullParameters.snapshot()),
49 fAltParameters(altParameters.snapshot())
50 {
51 }
52
56 }
57
58 static void SetAlwaysReuseNLL(bool flag);
59
60 void SetReuseNLL(bool flag) { fReuseNll = flag ; }
61
62 void SetNullParameters(const RooArgSet& nullParameters) {
64 fFirstEval = true;
65 fNullParameters = (RooArgSet*) nullParameters.snapshot();
66 }
67
68 void SetAltParameters(const RooArgSet& altParameters) {
70 fFirstEval = true;
71 fAltParameters = (RooArgSet*) altParameters.snapshot();
72 }
73
74 /// this should be possible with RooAbsCollection
76 if (!fNullParameters->equals(*fAltParameters)) return false;
77
78 bool ret = true;
79
80 for (auto nullIt = fNullParameters->begin(), altIt = fAltParameters->begin();
81 nullIt != fNullParameters->end() && altIt != fAltParameters->end(); ++nullIt, ++altIt) {
82 RooAbsReal *null = static_cast<RooAbsReal *>(*nullIt);
83 RooAbsReal *alt = static_cast<RooAbsReal *>(*altIt);
84 if (null->getVal() != alt->getVal())
85 ret = false;
86 }
87
88 return ret;
89 }
90
91
92 /// set the conditional observables which will be used when creating the NLL
93 /// so the pdf's will not be normalized on the conditional observables when computing the NLL
95
96 /// set the global observables which will be used when creating the NLL
97 /// so the constraint pdf's will be normalized correctly on the global observables when computing the NLL
98 void SetGlobalObservables(const RooArgSet& set) override {fGlobalObs.removeAll(); fGlobalObs.add(set);}
99
100 double Evaluate(RooAbsData& data, RooArgSet& nullPOI) override;
101
102 virtual void EnableDetailedOutput( bool e=true ) { fDetailedOutputEnabled = e; fDetailedOutput = nullptr; }
103 const RooArgSet* GetDetailedOutput(void) const override { return fDetailedOutput.get(); }
104
105 const TString GetVarName() const override {
106 return "log(L(#mu_{1}) / L(#mu_{0}))";
107 }
108
109 private:
110
111 RooAbsPdf* fNullPdf = nullptr;
112 RooAbsPdf* fAltPdf = nullptr;
117 bool fFirstEval = true;
118
120 std::unique_ptr<RooArgSet> fDetailedOutput; ///<!
121
122 std::unique_ptr<RooAbsReal> fNllNull; ///<! transient copy of the null NLL
123 std::unique_ptr<RooAbsReal> fNllAlt; ///<! transient copy of the alt NLL
124 static bool fgAlwaysReuseNll ;
125 bool fReuseNll = false;
126
127
129};
130
131}
132
133#endif
#define e(i)
Definition RSha256.hxx:103
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
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.
Storage_t const & get() const
Const access to the underlying stl container.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
TIterator begin()
TIterator end() and range-based for loops.")
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
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...
std::unique_ptr< RooAbsReal > fNllAlt
! transient copy of the alt NLL
std::unique_ptr< RooAbsReal > fNllNull
! transient copy of the null NLL
void SetNullParameters(const RooArgSet &nullParameters)
SimpleLikelihoodRatioTestStat(RooAbsPdf &nullPdf, RooAbsPdf &altPdf)
Takes null and alternate parameters from PDF. Can be overridden.
SimpleLikelihoodRatioTestStat()=default
Constructor for proof. Do not use.
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...
Basic string class.
Definition TString.h:139
Namespace for the RooStats classes.
Definition Asimov.h:19