Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
LikelihoodInterval.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 /*****************************************************************************
12  * Project: RooStats
13  * Package: RooFit/RooStats
14  * @(#)root/roofit/roostats:$Id$
15  * Authors:
16  * Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
17  *
18  *****************************************************************************/
19 
20 
21 /** \class RooStats::LikelihoodInterval
22  \ingroup Roostats
23 
24  LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
25  It implements a connected N-dimensional intervals based on the contour of a likelihood ratio.
26  The boundary of the interval is equivalent to a MINUIT/MINOS contour about the maximum likelihood estimator
27 
28  The interval does not need to be an ellipse (eg. it is not the HESSE error matrix).
29  The level used to make the contour is the same as that used in MINOS, eg. it uses Wilks' theorem,
30  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
31  N is the number of parameters of interest.
32 
33 
34  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'
35  theorem to be true.
36 
37  Also note, one can use any RooAbsReal as the function that will be used in the contour; however, the level of the contour
38  is based on Wilks' theorem as stated above.
39 
40 
41 #### References
42 
43 * 1. F. James., Minuit.Long writeup D506, CERN, 1998.
44 
45 */
46 
47 
49 #include "RooStats/RooStatsUtils.h"
50 
51 #include "RooAbsReal.h"
52 #include "RooMsgService.h"
53 
54 #include "Math/WrappedFunction.h"
55 #include "Math/Minimizer.h"
56 #include "Math/Factory.h"
57 #include "Math/MinimizerOptions.h"
58 #include "RooFunctor.h"
59 #include "RooProfileLL.h"
60 
61 #include "TMinuitMinimizer.h"
62 
63 #include <string>
64 #include <algorithm>
65 #include <functional>
66 #include <ctype.h> // need to use c version of toupper defined here
67 
68 /*
69 // for debugging
70 #include "RooNLLVar.h"
71 #include "RooProfileLL.h"
72 #include "RooDataSet.h"
73 #include "RooAbsData.h"
74 */
75 
76 ClassImp(RooStats::LikelihoodInterval); ;
77 
78 using namespace RooStats;
79 using namespace std;
80 
81 
82 ////////////////////////////////////////////////////////////////////////////////
83 /// Default constructor with name and title
84 
85 LikelihoodInterval::LikelihoodInterval(const char* name) :
86  ConfInterval(name), fBestFitParams(0), fLikelihoodRatio(0), fConfidenceLevel(0.95)
87 {
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// Alternate constructor taking a pointer to the profile likelihood ratio, parameter of interest and
92 /// optionally a snapshot of best parameter of interest for interval
93 
94 LikelihoodInterval::LikelihoodInterval(const char* name, RooAbsReal* lr, const RooArgSet* params, RooArgSet * bestParams) :
95  ConfInterval(name),
96  fParameters(*params),
97  fBestFitParams(bestParams),
98  fLikelihoodRatio(lr),
99  fConfidenceLevel(0.95)
100 {
101 }
102 
103 
104 ////////////////////////////////////////////////////////////////////////////////
105 /// Destructor
106 
107 LikelihoodInterval::~LikelihoodInterval()
108 {
109  if (fBestFitParams) delete fBestFitParams;
110  if (fLikelihoodRatio) delete fLikelihoodRatio;
111 }
112 
113 
114 ////////////////////////////////////////////////////////////////////////////////
115 /// This is the main method to satisfy the RooStats::ConfInterval interface.
116 /// It returns true if the parameter point is in the interval.
117 
118 Bool_t LikelihoodInterval::IsInInterval(const RooArgSet &parameterPoint) const
119 {
120  RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
121  RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
122  // Method to determine if a parameter point is in the interval
123  if( !this->CheckParameters(parameterPoint) ) {
124  std::cout << "parameters don't match" << std::endl;
125  RooMsgService::instance().setGlobalKillBelow(msglevel);
126  return false;
127  }
128 
129  // make sure likelihood ratio is set
130  if(!fLikelihoodRatio) {
131  std::cout << "likelihood ratio not set" << std::endl;
132  RooMsgService::instance().setGlobalKillBelow(msglevel);
133  return false;
134  }
135 
136 
137 
138  // set parameters
139  SetParameters(&parameterPoint, fLikelihoodRatio->getVariables() );
140 
141 
142  // evaluate likelihood ratio, see if it's bigger than threshold
143  if (fLikelihoodRatio->getVal()<0){
144  std::cout << "The likelihood ratio is < 0, indicates a bad minimum or numerical precision problems. Will return true" << std::endl;
145  RooMsgService::instance().setGlobalKillBelow(msglevel);
146  return true;
147  }
148 
149 
150  // here we use Wilks' theorem.
151  if ( TMath::Prob( 2* fLikelihoodRatio->getVal(), parameterPoint.getSize()) < (1.-fConfidenceLevel) ){
152  RooMsgService::instance().setGlobalKillBelow(msglevel);
153  return false;
154  }
155 
156 
157  RooMsgService::instance().setGlobalKillBelow(msglevel);
158 
159  return true;
160 
161 }
162 
163 ////////////////////////////////////////////////////////////////////////////////
164 /// returns list of parameters
165 
166 RooArgSet* LikelihoodInterval::GetParameters() const
167 {
168  return new RooArgSet(fParameters);
169 }
170 
171 ////////////////////////////////////////////////////////////////////////////////
172 /// check that the parameters are correct
173 
174 Bool_t LikelihoodInterval::CheckParameters(const RooArgSet &parameterPoint) const
175 {
176  if (parameterPoint.getSize() != fParameters.getSize() ) {
177  std::cout << "size is wrong, parameters don't match" << std::endl;
178  return false;
179  }
180  if ( ! parameterPoint.equals( fParameters ) ) {
181  std::cout << "size is ok, but parameters don't match" << std::endl;
182  return false;
183  }
184  return true;
185 }
186 
187 
188 
189 ////////////////////////////////////////////////////////////////////////////////
190 /// Compute lower limit, check first if limit has been computed
191 /// status is a boolean flag which will b set to false in case of error
192 /// and is true if calculation is successful
193 /// in case of error return also a lower limit value of zero
194 
195 Double_t LikelihoodInterval::LowerLimit(const RooRealVar& param, bool & status)
196 {
197  double lower = 0;
198  double upper = 0;
199  status = FindLimits(param, lower, upper);
200  return lower;
201 }
202 
203 ////////////////////////////////////////////////////////////////////////////////
204 /// Compute upper limit, check first if limit has been computed
205 /// status is a boolean flag which will b set to false in case of error
206 /// and is true if calculation is successful
207 /// in case of error return also a lower limit value of zero
208 
209 Double_t LikelihoodInterval::UpperLimit(const RooRealVar& param, bool & status)
210 {
211  double lower = 0;
212  double upper = 0;
213  status = FindLimits(param, lower, upper);
214  return upper;
215 }
216 
217 
218 void LikelihoodInterval::ResetLimits() {
219  // reset map with cached limits - called every time the test size or CL has been changed
220  fLowerLimits.clear();
221  fUpperLimits.clear();
222 }
223 
224 
225 bool LikelihoodInterval::CreateMinimizer() {
226  // internal function to create minimizer object needed to find contours or interval limits
227  // (running MINOS).
228  // Minimizer must be Minuit or Minuit2
229 
230  RooProfileLL * profilell = dynamic_cast<RooProfileLL*>(fLikelihoodRatio);
231  if (!profilell) return false;
232 
233  RooAbsReal & nll = profilell->nll();
234  // bind the nll function in the right interface for the Minimizer class
235  // as a function of only the parameters (poi + nuisance parameters)
236 
237  RooArgSet * partmp = profilell->getVariables();
238  // need to remove constant parameters
239  RemoveConstantParameters(partmp);
240 
241  RooArgList params(*partmp);
242  delete partmp;
243 
244  // need to restore values and errors for POI
245  if (fBestFitParams) {
246  for (int i = 0; i < params.getSize(); ++i) {
247  RooRealVar & par = (RooRealVar &) params[i];
248  RooRealVar * fitPar = (RooRealVar *) (fBestFitParams->find(par.GetName() ) );
249  if (fitPar) {
250  par.setVal( fitPar->getVal() );
251  par.setError( fitPar->getError() );
252  }
253  }
254  }
255 
256  const auto& config = GetGlobalRooStatsConfig();
257 
258  // now do binding of NLL with a functor for Minimizer
259  if (config.useLikelihoodOffset) {
260  ccoutI(InputArguments) << "LikelihoodInterval: using nll offset - set all RooAbsReal to hide the offset " << std::endl;
261  RooAbsReal::setHideOffset(kFALSE); // need to keep this false
262  }
263  fFunctor = std::make_shared<RooFunctor>(nll, RooArgSet(), params);
264 
265  std::string minimType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
266  std::transform(minimType.begin(), minimType.end(), minimType.begin(), (int(*)(int)) tolower );
267  *minimType.begin() = toupper( *minimType.begin());
268 
269  if (minimType != "Minuit" && minimType != "Minuit2") {
270  ccoutE(InputArguments) << minimType << " is wrong type of minimizer for getting interval limits or contours - must use Minuit or Minuit2" << std::endl;
271  return false;
272  }
273  // do not use static instance of TMInuit which could interfere with RooFit
274  if (minimType == "Minuit") TMinuitMinimizer::UseStaticMinuit(false);
275  // create minimizer class
276  fMinimizer = std::shared_ptr<ROOT::Math::Minimizer>(ROOT::Math::Factory::CreateMinimizer(minimType, "Migrad"));
277 
278  if (!fMinimizer.get()) return false;
279 
280  fMinFunc = std::static_pointer_cast<ROOT::Math::IMultiGenFunction>(
281  std::make_shared<ROOT::Math::WrappedMultiFunction<RooFunctor &>>(*fFunctor, fFunctor->nPar()) );
282  fMinimizer->SetFunction(*fMinFunc);
283 
284  // set minimizer parameters
285  assert( params.getSize() == int(fMinFunc->NDim()) );
286 
287  for (unsigned int i = 0; i < fMinFunc->NDim(); ++i) {
288  RooRealVar & v = (RooRealVar &) params[i];
289  fMinimizer->SetLimitedVariable( i, v.GetName(), v.getVal(), v.getError(), v.getMin(), v.getMax() );
290  }
291  // for finding the contour need to find first global minimum
292  bool iret = fMinimizer->Minimize();
293  if (!iret || fMinimizer->X() == 0) {
294  ccoutE(Minimization) << "Error: Minimization failed " << std::endl;
295  return false;
296  }
297 
298  //std::cout << "print minimizer result..........." << std::endl;
299  //fMinimizer->PrintResults();
300 
301  return true;
302 }
303 
304 bool LikelihoodInterval::FindLimits(const RooRealVar & param, double &lower, double & upper)
305 {
306  // Method to find both lower and upper limits using MINOS
307  // If cached values exist (limits have been already found) return them in that case
308  // check first if limit has been computed
309  // otherwise compute limit using MINOS
310  // in case of failure lower and upper will maintain previous value (will not be modified)
311 
312  std::map<std::string, double>::const_iterator itrl = fLowerLimits.find(param.GetName());
313  std::map<std::string, double>::const_iterator itru = fUpperLimits.find(param.GetName());
314  if ( itrl != fLowerLimits.end() && itru != fUpperLimits.end() ) {
315  lower = itrl->second;
316  upper = itru->second;
317  return true;
318  }
319 
320 
321  RooArgSet * partmp = fLikelihoodRatio->getVariables();
322  RemoveConstantParameters(partmp);
323  RooArgList params(*partmp);
324  delete partmp;
325  int ix = params.index(&param);
326  if (ix < 0 ) {
327  ccoutE(InputArguments) << "Error - invalid parameter " << param.GetName() << " specified for finding the interval limits " << std::endl;
328  return false;
329  }
330 
331  bool ret = true;
332  if (!fMinimizer.get()) ret = CreateMinimizer();
333  if (!ret) {
334  ccoutE(Eval) << "Error returned from minimization of likelihood function - cannot find interval limits " << std::endl;
335  return false;
336  }
337 
338  assert(fMinimizer.get());
339 
340  // getting a 1D interval so ndf = 1
341  double err_level = TMath::ChisquareQuantile(ConfidenceLevel(),1); // level for -2log LR
342  err_level = err_level/2; // since we are using -log LR
343  fMinimizer->SetErrorDef(err_level);
344 
345  unsigned int ivarX = ix;
346 
347  double elow = 0;
348  double eup = 0;
349  ret = fMinimizer->GetMinosError(ivarX, elow, eup );
350  if (!ret) {
351  ccoutE(Minimization) << "Error running Minos for parameter " << param.GetName() << std::endl;
352  return false;
353  }
354 
355  // WHEN error is zero normally is at limit
356  if (elow == 0) {
357  lower = param.getMin();
358  ccoutW(Minimization) << "Warning: lower value for " << param.GetName() << " is at limit " << lower << std::endl;
359  }
360  else
361  lower = fMinimizer->X()[ivarX] + elow; // elow is negative
362 
363  if (eup == 0) {
364  ccoutW(Minimization) << "Warning: upper value for " << param.GetName() << " is at limit " << upper << std::endl;
365  upper = param.getMax();
366  }
367  else
368  upper = fMinimizer->X()[ivarX] + eup;
369 
370  // store limits in the map
371  // minos return error limit = minValue +/- error
372  fLowerLimits[param.GetName()] = lower;
373  fUpperLimits[param.GetName()] = upper;
374 
375  return true;
376 }
377 
378 
379 Int_t LikelihoodInterval::GetContourPoints(const RooRealVar & paramX, const RooRealVar & paramY, Double_t * x, Double_t *y, Int_t npoints ) {
380  // use Minuit to find the contour of the likelihood function at the desired CL
381 
382  // check the parameters
383  // variable index in minimizer
384  // is index in the RooArgList obtained from the profileLL variables
385  RooArgSet * partmp = fLikelihoodRatio->getVariables();
386  RemoveConstantParameters(partmp);
387  RooArgList params(*partmp);
388  delete partmp;
389  int ix = params.index(&paramX);
390  int iy = params.index(&paramY);
391  if (ix < 0 || iy < 0) {
392  coutE(InputArguments) << "LikelihoodInterval - Error - invalid parameters specified for finding the contours; parX = " << paramX.GetName()
393  << " parY = " << paramY.GetName() << std::endl;
394  return 0;
395  }
396 
397  bool ret = true;
398  if (!fMinimizer.get()) ret = CreateMinimizer();
399  if (!ret) {
400  coutE(Eval) << "LikelihoodInterval - Error returned creating minimizer for likelihood function - cannot find contour points " << std::endl;
401  return 0;
402  }
403 
404  assert(fMinimizer.get());
405 
406  // getting a 2D contour so ndf = 2
407  double cont_level = TMath::ChisquareQuantile(ConfidenceLevel(),2); // level for -2log LR
408  cont_level = cont_level/2; // since we are using -log LR
409  fMinimizer->SetErrorDef(cont_level);
410 
411  unsigned int ncp = npoints;
412  unsigned int ivarX = ix;
413  unsigned int ivarY = iy;
414  coutI(Minimization) << "LikelihoodInterval - Finding the contour of " << paramX.GetName() << " ( " << ivarX << " ) and " << paramY.GetName() << " ( " << ivarY << " ) " << std::endl;
415  ret = fMinimizer->Contour(ivarX, ivarY, ncp, x, y );
416  if (!ret) {
417  coutE(Minimization) << "LikelihoodInterval - Error finding contour for parameters " << paramX.GetName() << " and " << paramY.GetName() << std::endl;
418  return 0;
419  }
420  if (int(ncp) < npoints) {
421  coutW(Minimization) << "LikelihoodInterval -Warning - Less points calculated in contours np = " << ncp << " / " << npoints << std::endl;
422  }
423 
424  return ncp;
425  }