Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ProfileLikelihoodCalculator.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer 28/07/2008
3
4/*************************************************************************
5 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class RooStats::ProfileLikelihoodCalculator
13 \ingroup Roostats
14
15The ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator
16(the interface class for tools which can produce both a RooStats HypoTestResult
17and ConfInterval). The tool uses the profile likelihood ratio as a test statistic,
18and assumes that Wilks' theorem is valid. Wilks' theorem states that \f$ -2 \cdot \ln(\lambda) \f$
19(profile likelihood ratio) is asymptotically distributed as a \f$ \chi^2 \f$ distribution
20with \f$ N \f$ degrees of freedom. Thus, \f$p\f$-values can be
21constructed, and the profile likelihood ratio can be used to construct a
22LikelihoodInterval. (In the future, this class could be extended to use toy
23Monte Carlo to calibrate the distribution of the test statistic).
24
25Usage: It uses the interface of the CombinedCalculator, so that it can be
26configured by specifying:
27
28 - A model common model (*e.g.* a family of specific models, which includes both
29 the null and alternate)
30 - A data set
31 - A set of parameters of interest. The nuisance parameters will be all other
32 parameters of the model.
33 - A set of parameters which specify the null hypothesis (including values
34 and const/non-const status).
35
36The interface allows one to pass the model, data, and parameters either directly
37or via a ModelConfig class. The alternate hypothesis leaves the parameter free
38to take any value other than those specified by the null hypothesis. There is
39therefore no need to specify the alternate parameters.
40
41After configuring the calculator, one only needs to call GetHypoTest() (which
42will return a HypoTestResult pointer) or GetInterval() (which will return a
43ConfInterval pointer).
44
45This calculator can work with both one-dimensional intervals or multi-
46dimensional ones (contours).
47
48Note that for hypothesis tests, it is often better to use the
49AsymptoticCalculator, which can compute in addition the expected
50\f$p\f$-value using an Asimov data set.
51
52*/
53
55
57
60
61#include "RooFitResult.h"
62#include "RooRealVar.h"
63#include "RooProfileLL.h"
64#include "RooGlobalFunc.h"
65#include "RooMsgService.h"
66
68#include "RooMinimizer.h"
69//#include "RooProdPdf.h"
70
71using namespace std;
72
74
75using namespace RooFit;
76using namespace RooStats;
77
78
79////////////////////////////////////////////////////////////////////////////////
80/// default constructor
81
83 CombinedCalculator(), fGlobalFitDone(false)
84{
85}
86
88 double size, const RooArgSet* nullParams ) :
89 CombinedCalculator(data,pdf, paramsOfInterest, size, nullParams ),
90 fGlobalFitDone(false)
91{
92 // constructor from pdf and parameters
93 // the pdf must contain eventually the nuisance parameters
94}
95
98 fGlobalFitDone(false)
99{
100 // construct from a ModelConfig. Assume data model.GetPdf() will provide full description of model including
101 // constraint term on the nuisances parameters
102 assert(model.GetPdf() );
103}
104
105
106////////////////////////////////////////////////////////////////////////////////
107/// destructor
108/// cannot delete prod pdf because it will delete all the composing pdf's
109/// if (fOwnPdf) delete fPdf;
110/// fPdf = 0;
111
113
115 // reset and clear fit result
116 // to be called when a new model or data are set in the calculator
117 fFitResult.reset();
118}
119
121 // perform a global fit of the likelihood letting with all parameter of interest and
122 // nuisance parameters
123 // keep the list of fitted parameters
124
125 DoReset();
126 RooAbsPdf * pdf = GetPdf();
128 if (!data || !pdf ) return 0;
129
130 // get all non-const parameters
131 std::unique_ptr<RooArgSet> constrainedParams{pdf->getParameters(*data)};
132 if (!constrainedParams) return 0;
133 RemoveConstantParameters(&*constrainedParams);
134
135 const auto& config = GetGlobalRooStatsConfig();
137 RooFit::Offset(config.useLikelihoodOffset) );
138
139 // check if global fit has been already done
140 if (fFitResult && fGlobalFitDone) {
141 return RooFit::OwningPtr<RooAbsReal>{std::move(nll)};
142 }
143
144 // calculate MLE
145 oocoutP(nullptr,Minimization) << "ProfileLikelihoodCalcultor::DoGLobalFit - find MLE " << std::endl;
146
147 fFitResult = std::unique_ptr<RooFitResult>{DoMinimizeNLL(&*nll)};
148
149 // print fit result
150 if (fFitResult) {
151 fFitResult->printStream( oocoutI(nullptr,Minimization), fFitResult->defaultPrintContents(0), fFitResult->defaultPrintStyle(0) );
152
153 if (fFitResult->status() != 0)
154 oocoutW(nullptr,Minimization) << "ProfileLikelihoodCalcultor::DoGlobalFit - Global fit failed - status = " << fFitResult->status() << std::endl;
155 else
156 fGlobalFitDone = true;
157 }
158
159 return RooFit::OwningPtr<RooAbsReal>{std::move(nll)};
160}
161
163 // Minimizer the given NLL using the default options
164
165 const char * minimType = ""; // empty string to select RooMinimizer default
166 const char * minimAlgo = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo().c_str();
168 int level = ROOT::Math::MinimizerOptions::DefaultPrintLevel() -1;// RooFit level starts from -1
170 // do global fit and store fit result for further use
171
172 const auto& config = GetGlobalRooStatsConfig();
173
174 RooMinimizer minim(*nll);
175 minim.setStrategy(strategy);
176 minim.setEps(tolerance);
177 minim.setPrintLevel(level);
178 minim.optimizeConst(2); // to optimize likelihood calculations
179 minim.setEvalErrorWall(config.useEvalErrorWall);
180
181 oocoutP(nullptr,Minimization) << "ProfileLikelihoodCalcultor::DoMinimizeNLL - using " << minim.minimizerType()
182 << " / " << minimAlgo << " with strategy " << strategy << std::endl;
183
184 int status = -1;
185 for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
186 status = minim.minimize(minimType,minimAlgo);
187 if (status%1000 == 0) { // ignore erros from Improve
188 break;
189 } else if (tries < maxtries) {
190 cout << " ----> Doing a re-scan first" << endl;
191 minim.minimize(minimType,"Scan");
192 if (tries == 2) {
193 if (strategy == 0 ) {
194 cout << " ----> trying with strategy = 1" << endl;;
195 minim.setStrategy(1);
196 }
197 else
198 tries++; // skip this trial if strategy is already 1
199 }
200 if (tries == 3) {
201 cout << " ----> trying with improve" << endl;;
202 minimType = "Minuit";
203 minimAlgo = "migradimproved";
204 }
205 }
206 }
207
208 return minim.save();
209}
210
211////////////////////////////////////////////////////////////////////////////////
212/// Main interface to get a RooStats::ConfInterval.
213/// It constructs a profile likelihood ratio, and uses that to construct a RooStats::LikelihoodInterval.
214
216// RooAbsPdf* pdf = fWS->pdf(fPdfName);
217// RooAbsData* data = fWS->data(fDataName);
218 RooAbsPdf * pdf = GetPdf();
220 if (!data || !pdf || fPOI.empty()) return nullptr;
221
222 std::unique_ptr<RooArgSet> constrainedParams{pdf->getParameters(*data)};
223 RemoveConstantParameters(&*constrainedParams);
224
225
226 /*
227 RooNLLVar* nll = new RooNLLVar("nll","",*pdf,*data, Extended(),Constrain(*constrainedParams));
228 RooProfileLL* profile = new RooProfileLL("pll","",*nll, *fPOI);
229 profile->addOwnedComponents(*nll) ; // to avoid memory leak
230 */
231
232 // do a global fit cloning the data
233 std::unique_ptr<RooAbsReal> nll{DoGlobalFit()};
234 if (!nll) return nullptr;
235
236 if (!fFitResult) {
237 return nullptr;
238 }
239
240 RooAbsReal* profile = nll->createProfile(fPOI);
241 profile->addOwnedComponents(std::move(nll)) ; // to avoid memory leak
242
243 // t.b.f. " RooProfileLL should keep and provide possibility to query on global minimum
244 // set POI to fit value (this will speed up profileLL calculation of global minimum)
245 const RooArgList & fitParams = fFitResult->floatParsFinal();
246 for (int i = 0; i < fitParams.getSize(); ++i) {
247 RooRealVar & fitPar = (RooRealVar &) fitParams[i];
248 RooRealVar * par = (RooRealVar*) fPOI.find( fitPar.GetName() );
249 if (par) {
250 par->setVal( fitPar.getVal() );
251 par->setError( fitPar.getError() );
252 }
253 }
254
255 // do this so profile will cache inside the absolute minimum and
256 // minimum values of nuisance parameters
257 // (no need to this here)
258 // profile->getVal();
259 //RooMsgService::instance().setGlobalKillBelow(RooFit::DEBUG) ;
260 // profile->Print();
261
262 TString name = TString("LikelihoodInterval_");// + TString(GetName() );
263
264 // make a list of fPOI with fit result values and pass to LikelihoodInterval class
265 // bestPOI is a cloned list of POI only with their best fit values
266 RooArgSet fitParSet(fitParams);
267 RooArgSet * bestPOI = new RooArgSet();
268 for (auto const *arg : fPOI){
269 RooAbsArg * p = fitParSet.find( arg->GetName() );
270 if (p) bestPOI->addClone(*p);
271 else bestPOI->addClone(*arg);
272 }
273 // fPOI contains the parameter of interest of the PL object
274 // and bestPOI contains a snapshot with the best fit values
275 LikelihoodInterval* interval = new LikelihoodInterval(name, profile, &fPOI, bestPOI);
276 interval->SetConfidenceLevel(1.-fSize);
277 return interval;
278}
279
280////////////////////////////////////////////////////////////////////////////////
281/// Main interface to get a HypoTestResult.
282/// It does two fits:
283/// 1. The first lets the null parameters float, so it's a maximum likelihood estimate.
284/// 2. The second is to the null model (fixing null parameters to their specified values): *e.g.* a conditional maximum likelihood.
285/// Since not all parameters are floating, this likelihood will be lower than the unconditional model.
286///
287/// The ratio of the likelihood obtained from the conditional MLE to the MLE is the profile likelihood ratio.
288/// Wilks' theorem is used to get \f$p\f$-values.
289
291// RooAbsPdf* pdf = fWS->pdf(fPdfName);
292// RooAbsData* data = fWS->data(fDataName);
293 RooAbsPdf * pdf = GetPdf();
295
296
297 if (!data || !pdf) return 0;
298
299 if (fNullParams.empty()) return 0;
300
301 // make a clone and ordered list since a vector will be associated to keep parameter values
302 // clone the list since first fit will changes the fNullParams values
303 RooArgList poiList;
304 poiList.addClone(fNullParams); // make a clone list
305
306
307 // do a global fit
308 std::unique_ptr<RooAbsReal> nll{DoGlobalFit()};
309 if (!nll) return nullptr;
310
311 if (!fFitResult) {
312 return nullptr;
313 }
314
315 std::unique_ptr<RooArgSet> constrainedParams{pdf->getParameters(*data)};
316 RemoveConstantParameters(&*constrainedParams);
317
318 double nLLatMLE = fFitResult->minNll();
319 // in case of using offset need to save offset value
320 double nlloffset = (RooStats::IsNLLOffset() ) ? nll->getVal() - nLLatMLE : 0;
321
322 // set POI to given values, set constant, calculate conditional MLE
323 std::vector<double> oldValues(poiList.getSize() );
324 for (unsigned int i = 0; i < oldValues.size(); ++i) {
325 RooRealVar * mytarget = (RooRealVar*) constrainedParams->find(poiList[i].GetName());
326 if (mytarget) {
327 oldValues[i] = mytarget->getVal();
328 mytarget->setVal( ( (RooRealVar&) poiList[i] ).getVal() );
329 mytarget->setConstant(true);
330 }
331 }
332
333
334
335 // perform the fit only if nuisance parameters are available
336 // get nuisance parameters
337 // nuisance parameters are the non const parameters from the likelihood parameters
338 RooArgSet nuisParams(*constrainedParams);
339
340 // need to remove the parameter of interest
341 RemoveConstantParameters(&nuisParams);
342
343 // check there are variable parameter in order to do a fit
344 bool existVarParams = false;
345 for (auto const *myarg : static_range_cast<RooRealVar *> (nuisParams)) {
346 if ( !myarg->isConstant() ) {
347 existVarParams = true;
348 break;
349 }
350 }
351
352 double nLLatCondMLE = nLLatMLE;
353 if (existVarParams) {
354 oocoutP(nullptr,Minimization) << "ProfileLikelihoodCalcultor::GetHypoTest - do conditional fit " << std::endl;
355
356 std::unique_ptr<RooFitResult> fit2{DoMinimizeNLL(&*nll)};
357
358 // print fit result
359 if (fit2) {
360 nLLatCondMLE = fit2->minNll();
361 fit2->printStream( oocoutI(nullptr,Minimization), fit2->defaultPrintContents(0), fit2->defaultPrintStyle(0) );
362
363 if (fit2->status() != 0)
364 oocoutW(nullptr,Minimization) << "ProfileLikelihoodCalcultor::GetHypotest - Conditional fit failed - status = " << fit2->status() << std::endl;
365 }
366
367 }
368 else {
369 // get just the likelihood value (no need to do a fit since the likelihood is a constant function)
370 nLLatCondMLE = nll->getVal();
371 // this value contains the offset
372 if (RooStats::IsNLLOffset() ) nLLatCondMLE -= nlloffset;
373 }
374
375 // Use Wilks' theorem to translate -2 log lambda into a significance/p-value
376 double deltaNLL = std::max( nLLatCondMLE-nLLatMLE, 0.);
377
378 // get number of free parameter of interest
380 int ndf = poiList.getSize();
381
382 double pvalue = ROOT::Math::chisquared_cdf_c( 2* deltaNLL, ndf);
383
384 // in case of one dimension (1 poi) do the one-sided p-value (need to divide by 2)
385 if (ndf == 1) pvalue = 0.5 * pvalue;
386
387 TString name = TString("ProfileLRHypoTestResult_");// + TString(GetName() );
388 HypoTestResult* htr = new HypoTestResult(name, pvalue, 0 );
389
390 // restore previous value of poi
391 for (unsigned int i = 0; i < oldValues.size(); ++i) {
392 RooRealVar * mytarget = (RooRealVar*) constrainedParams->find(poiList[i].GetName());
393 if (mytarget) {
394 mytarget->setVal(oldValues[i] );
395 mytarget->setConstant(false);
396 }
397 }
398
399 return htr;
400
401}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define oocoutW(o, a)
#define oocoutI(o, a)
#define oocoutP(o, a)
#define ClassImp(name)
Definition Rtypes.h:377
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
static const std::string & DefaultMinimizerAlgo()
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
Int_t getSize() const
Return the number of elements in the collection.
const char * GetName() const override
Returns name of object.
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:59
virtual RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, const RooLinkedList &cmdList={})
Construct representation of -log(L) of PDF with given dataset.
void setConstant(bool value=true)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
void optimizeConst(int flag)
If flag is true, perform constant term optimization on function being minimized.
RooFit::OwningPtr< RooFitResult > save(const char *name=nullptr, const char *title=nullptr)
Save and return a RooFitResult snapshot of current minimizer status.
void setEvalErrorWall(bool flag)
void setEps(double eps)
Change MINUIT epsilon.
void setPrintLevel(int newLevel)
Change the MINUIT internal printing level.
std::string const & minimizerType() const
int minimize(const char *type, const char *alg=nullptr)
Minimise the function passed in the constructor.
void setStrategy(int strat)
Change MINUIT strategy to istrat.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setVal(double value) override
Set value of variable to 'value'.
void setError(double value)
Definition RooRealVar.h:64
double getError() const
Definition RooRealVar.h:62
CombinedCalculator is an interface class for a tools which can produce both RooStats HypoTestResults ...
RooArgSet fNullParams
RooArgSet specifying null parameters for hypothesis test.
RooArgSet fPOI
RooArgSet specifying parameters of interest for interval.
RooArgSet fGlobalObs
RooArgSet specifying the global observables.
double fSize
size of the test (eg. specified rate of Type I error)
RooArgSet fConditionalObs
RooArgSet specifying the conditional observables.
HypoTestResult is a base class for results from hypothesis tests.
LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
void SetConfidenceLevel(double cl) override
set the confidence level for the interval (e.g 0.682 for a 1-sigma interval)
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
The ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface cla...
static RooFit::OwningPtr< RooFitResult > DoMinimizeNLL(RooAbsReal *nll)
minimize likelihood
HypoTestResult * GetHypoTest() const override
Return the hypothesis test result obtained from the likelihood ratio of the maximum likelihood value ...
RooFit::OwningPtr< RooAbsReal > DoGlobalFit() const
perform a global fit
bool fGlobalFitDone
flag to control if a global fit has been done
~ProfileLikelihoodCalculator() override
destructor cannot delete prod pdf because it will delete all the composing pdf's if (fOwnPdf) delete ...
void DoReset() const
clear internal fit result
LikelihoodInterval * GetInterval() const override
Return a likelihood interval.
ProfileLikelihoodCalculator()
Default constructor (needed for I/O)
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:139
RooCmdArg Offset(std::string const &mode)
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg GlobalObservables(Args_t &&... argsOrArgSet)
RooCmdArg CloneData(bool flag)
RooCmdArg ConditionalObservables(Args_t &&... argsOrArgSet)
Create a RooCmdArg to declare conditional observables.
double chisquared_cdf_c(double x, double r, double x0=0)
Complement of the cumulative distribution function of the distribution with degrees of freedom (upp...
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition Common.h:18
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:43
Namespace for the RooStats classes.
Definition Asimov.h:19
void RemoveConstantParameters(RooArgSet *set)
RooStatsConfig & GetGlobalRooStatsConfig()
Retrieve the config object which can be used to set flags for things like offsetting the likelihood o...
bool IsNLLOffset()
function returning if the flag to check if the flag to use NLLOffset is set