ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
LikelihoodInterval.h
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 #ifndef RooStats_LikelihoodInterval
12 #define RooStats_LikelihoodInterval
13 
14 #ifndef RooStats_ConfInterval
15 #include "RooStats/ConfInterval.h"
16 #endif
17 
18 #ifndef ROO_ARG_SET
19 #include "RooArgSet.h"
20 #endif
21 
22 #ifndef ROO_ABS_REAL
23 #include "RooAbsReal.h"
24 #endif
25 
26 #ifndef ROOT_Math_IFunctionfwd
27 #include "Math/IFunctionfwd.h"
28 #endif
29 
30 #include <map>
31 #include <memory>
32 
33 namespace ROOT {
34  namespace Math {
35  class Minimizer;
36  }
37 }
38 
39 namespace RooStats {
40 
41 
42 /**
43 
44  \ingroup Roostats
45 
46  LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
47  It implements a connected N-dimensional intervals based on the contour of a likelihood ratio.
48  The boundary of the inteval is equivalent to a MINUIT/MINOS contour about the maximum likelihood estimator
49 
50  The interval does not need to be an ellipse (eg. it is not the HESSE error matrix).
51  The level used to make the contour is the same as that used in MINOS, eg. it uses Wilks' theorem,
52  which states that under certain regularity conditions the function -2* log (profile likelihood ratio) is asymptotically distributed as a chi^2 with N-dof, where
53  N is the number of parameters of interest.
54 
55 
56  Note, a boundary on the parameter space (eg. s>= 0) or a degeneracy (eg. mass of signal if Nsig = 0) can lead to violations of the conditions necessary for Wilks'
57  theorem to be true.
58 
59  Also note, one can use any RooAbsReal as the function that will be used in the contour; however, the level of the contour
60  is based on Wilks' theorem as stated above.
61 
62 
63 #### References
64 
65 * 1. F. James., Minuit.Long writeup D506, CERN, 1998.
66 
67 
68 */
69 
70 
72 
73  public:
74 
75  /// default constructor
76  explicit LikelihoodInterval(const char* name = 0);
77 
78  //// construct the interval from a Profile Likelihood object, parameter of interest and optionally a snapshot of
79  //// POI with their best fit values
80  LikelihoodInterval(const char* name, RooAbsReal*, const RooArgSet*, RooArgSet * = 0);
81 
82  /// destructor
83  virtual ~LikelihoodInterval();
84 
85  /// check if given point is in the interval
86  virtual Bool_t IsInInterval(const RooArgSet&) const;
87 
88  /// set the confidence level for the interval (e.g 0.682 for a 1-sigma interval)
90 
91  /// return confidence level
92  virtual Double_t ConfidenceLevel() const {return fConfidenceLevel;}
93 
94  /// return a cloned list of parameters of interest. User manages the return object
95  virtual RooArgSet* GetParameters() const;
96 
97  /// check if parameters are correct (i.e. they are the POI of this interval)
98  Bool_t CheckParameters(const RooArgSet&) const ;
99 
100 
101  /// return the lower bound of the interval on a given parameter
102  Double_t LowerLimit(const RooRealVar& param) { bool ok; return LowerLimit(param,ok); }
103  Double_t LowerLimit(const RooRealVar& param, bool & status) ;
104 
105  /// return the upper bound of the interval on a given parameter
106  Double_t UpperLimit(const RooRealVar& param) { bool ok; return UpperLimit(param,ok); }
107  Double_t UpperLimit(const RooRealVar& param, bool & status) ;
108 
109  /// find both lower and upper interval boundaries for a given parameter
110  /// retun false if the bounds have not been found
111  Bool_t FindLimits(const RooRealVar & param, double & lower, double &upper);
112 
113  /**
114  return the 2D-contour points for the given subset of parameters
115  by default make the contour using 30 points. The User has to preallocate the x and y array which will return
116  the set of x and y points defining the contour.
117  The return value of the funciton specify the number of contour point found.
118  In case of error a zero is returned
119  */
120  Int_t GetContourPoints(const RooRealVar & paramX, const RooRealVar & paramY, Double_t * x, Double_t *y, Int_t npoints = 30);
121 
122  /// return the profile log-likelihood ratio function
124 
125  /// return a pointer to a snapshot with best fit parameter of interest
126  const RooArgSet * GetBestFitParameters() const { return fBestFitParams; }
127 
128  protected:
129 
130  /// reset the cached limit values
131  void ResetLimits();
132 
133  /// internal function to create the minimizer for finding the contours
134  bool CreateMinimizer();
135 
136  private:
137 
138  RooArgSet fParameters; /// parameters of interest for this interval
139  RooArgSet * fBestFitParams; /// snapshot of the model parameters with best fit value (managed internally)
140  RooAbsReal* fLikelihoodRatio; /// likelihood ratio function used to make contours (managed internally)
141  Double_t fConfidenceLevel; /// Requested confidence level (eg. 0.95 for 95% CL)
142  std::map<std::string, double> fLowerLimits; /// map with cached lower bound values
143  std::map<std::string, double> fUpperLimits; /// map with cached upper bound values
144  std::shared_ptr<ROOT::Math::Minimizer > fMinimizer; //! transient pointer to minimizer class used to find limits and contour
145  std::shared_ptr<RooFunctor> fFunctor; //! transient pointer to functor class used by the minimizer
146  std::shared_ptr<ROOT::Math::IMultiGenFunction> fMinFunc; //! transient pointer to the minimization function
147 
148  ClassDef(LikelihoodInterval,1) /// Concrete implementation of a ConfInterval based on a likelihood ratio
149 
150  };
151 }
152 
153 #endif
virtual Double_t ConfidenceLevel() const
return confidence level
RooAbsReal * fLikelihoodRatio
snapshot of the model parameters with best fit value (managed internally)
std::map< std::string, double > fLowerLimits
Requested confidence level (eg. 0.95 for 95% CL)
std::shared_ptr< ROOT::Math::Minimizer > fMinimizer
map with cached upper bound values
std::map< std::string, double > fUpperLimits
map with cached lower bound values
LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
Bool_t CheckParameters(const RooArgSet &) const
check if parameters are correct (i.e. they are the POI of this interval)
std::shared_ptr< ROOT::Math::IMultiGenFunction > fMinFunc
transient pointer to functor class used by the minimizer
LikelihoodInterval(const char *name=0)
default constructor
virtual RooArgSet * GetParameters() const
return a cloned list of parameters of interest. User manages the return object
int Int_t
Definition: RtypesCore.h:41
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
bool Bool_t
Definition: RtypesCore.h:59
Bool_t FindLimits(const RooRealVar &param, double &lower, double &upper)
find both lower and upper interval boundaries for a given parameter retun false if the bounds have no...
static const float upper
Definition: main.cpp:49
RooArgSet * fBestFitParams
parameters of interest for this interval
Double_t LowerLimit(const RooRealVar &param)
return the lower bound of the interval on a given parameter
RooAbsReal * GetLikelihoodRatio()
return the profile log-likelihood ratio function
Double_t x[n]
Definition: legend1.C:17
#define ClassDef(name, id)
Definition: Rtypes.h:254
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
Double_t UpperLimit(const RooRealVar &param)
return the upper bound of the interval on a given parameter
bool CreateMinimizer()
internal function to create the minimizer for finding the contours
void ResetLimits()
reset the cached limit values
RooCmdArg Minimizer(const char *type, const char *alg=0)
const RooArgSet * GetBestFitParameters() const
return a pointer to a snapshot with best fit parameter of interest
std::shared_ptr< RooFunctor > fFunctor
transient pointer to minimizer class used to find limits and contour
ConfInterval is an interface class for a generic interval in the RooStats framework.
Definition: ConfInterval.h:44
virtual ~LikelihoodInterval()
destructor
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
virtual Bool_t IsInInterval(const RooArgSet &) const
check if given point is in the interval
Double_t y[n]
Definition: legend1.C:17
Int_t GetContourPoints(const RooRealVar &paramX, const RooRealVar &paramY, Double_t *x, Double_t *y, Int_t npoints=30)
return the 2D-contour points for the given subset of parameters by default make the contour using 30 ...
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (e.g 0.682 for a 1-sigma interval)
#define name(a, b)
Definition: linkTestLib0.cpp:5
Double_t fConfidenceLevel
likelihood ratio function used to make contours (managed internally)
static const float lower
Definition: main.cpp:48