Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RatioOfProfiledLikelihoodsTestStat.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
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::RatioOfProfiledLikelihoodsTestStat
12 \ingroup Roostats
13
14TestStatistic that returns the ratio of profiled likelihoods.
15
16By default the calculation is:
17
18\f[
19 \log{ \frac{ \lambda(\mu_{alt} , {conditional \: MLE \: for \: alt \: nuisance}) }
20 { \lambda(\mu_{null} , {conditional \: MLE \: for \: null \: nuisance}) } }
21\f]
22
23where \f$ \lambda \f$ is the profile likelihood ratio, so the
24MLE for the null and alternate are subtracted off.
25
26If `SetSubtractMLE(false)` then it calculates:
27
28\f[
29 \log{ \frac{ L(\mu_{alt} , {conditional \: MLE \: for \: alt \: nuisance}) }
30 { L(\mu_{null} , {conditional \: MLE \: for \: null \: nuisance}) } }
31\f]
32
33where \f$ L \f$ is the Likelihood function.
34
35The values of the parameters of interest for the alternative
36hypothesis are taken at the time of the construction.
37If empty, it treats all free parameters as nuisance parameters.
38
39The value of the parameters of interest for the null hypotheses
40are given at each call of Evaluate.
41
42This test statistic is often called the Tevatron test statistic, because it has
43been used by the Tevatron experiments.
44*/
45
47
49
50#include "RooArgSet.h"
51#include "RooAbsData.h"
52#include "TMath.h"
53#include "RooMsgService.h"
54#include "RooGlobalFunc.h"
55
56
58
60
61////////////////////////////////////////////////////////////////////////////////
62/// returns -logL(poi, conditional MLE of nuisance params)
63/// subtract off the global MLE or not depending on the option
64/// It is the numerator or the denominator of the ratio (depending on the pdf)
65///
66/// L.M. : not sure why this method is needed now
67
69 int type = (fSubtractMLE) ? 0 : 2;
70
71 // null
72 if (&pdf == fNullProfile.GetPdf()) {
73 return fNullProfile.EvaluateProfileLikelihood(type, data, poi);
74 } else if (&pdf == fAltProfile.GetPdf()) {
75 return fAltProfile.EvaluateProfileLikelihood(type, data, poi);
76 }
77
78 oocoutE(nullptr,InputArguments) << "RatioOfProfiledLikelihoods::ProfileLikelihood - invalid pdf used for computing the profiled likelihood - return NaN"
79 << std::endl;
80
81 return TMath::QuietNaN();
82
83}
84
85////////////////////////////////////////////////////////////////////////////////
86/// evaluate the ratio of profile likelihood
87
89
90 int type = (fSubtractMLE) ? 0 : 2;
91
92 // null
93 double nullNLL = fNullProfile.EvaluateProfileLikelihood(type, data, nullParamsOfInterest);
94 const RooArgSet *nullset = fNullProfile.GetDetailedOutput();
95
96 // alt
97 double altNLL = fAltProfile.EvaluateProfileLikelihood(type, data, *fAltPOI);
98 const RooArgSet *altset = fAltProfile.GetDetailedOutput();
99
100 if (fDetailedOutput != nullptr) {
101 delete fDetailedOutput;
102 fDetailedOutput = nullptr;
103 }
106 for (auto const *var : static_range_cast<RooRealVar *>(*nullset)) {
107 auto cloneVar = std::make_unique<RooRealVar>(TString::Format("nullprof_%s", var->GetName()),
108 TString::Format("%s for null", var->GetTitle()), var->getVal());
109 fDetailedOutput->addOwned(std::move(cloneVar));
110 }
111 for (auto const *var : static_range_cast<RooRealVar *>(*altset)) {
112 auto cloneVar = std::make_unique<RooRealVar>(TString::Format("altprof_%s", var->GetName()),
113 TString::Format("%s for null", var->GetTitle()), var->getVal());
114 fDetailedOutput->addOwned(std::move(cloneVar));
115 }
116 }
117
118/*
119// set variables back to where they were
120nullParamsOfInterest = *saveNullPOI;
121*allVars = *saveAll;
122delete saveAll;
123delete allVars;
124*/
125
126 return nullNLL -altNLL;
127}
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
#define oocoutE(o, a)
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:56
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
double Evaluate(RooAbsData &data, RooArgSet &nullParamsOfInterest) override
evaluate the ratio of profile likelihood
double ProfiledLikelihood(RooAbsData &data, RooArgSet &poi, RooAbsPdf &pdf)
returns -logL(poi, conditional MLE of nuisance params) it does not subtract off the global MLE becaus...
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2385
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:913