Logo ROOT   6.10/09
Reference Guide
ProfileLikelihoodTestStat.cxx
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
3 // Additional Contributions: Giovanni Petrucciani
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::ProfileLikelihoodTestStat
13  \ingroup Roostats
14 
15 ProfileLikelihoodTestStat is an implementation of the TestStatistic interface
16 that calculates the profile likelihood ratio at a particular parameter point
17 given a dataset. It does not constitute a statistical test, for that one may
18 either use:
19 
20  - the ProfileLikelihoodCalculator that relies on asymptotic properties of the
21  Profile Likelihood Ratio
22  - the Neyman Construction classes with this class as a test statistic
23  - the Hybrid Calculator class with this class as a test statistic
24 
25 
26 */
27 
29 #include "RooFitResult.h"
30 #include "RooPullVar.h"
32 
33 #include "RooProfileLL.h"
34 #include "RooNLLVar.h"
35 #include "RooMsgService.h"
36 #include "RooMinimizer.h"
37 #include "RooArgSet.h"
38 #include "RooDataSet.h"
39 #include "TStopwatch.h"
40 
41 #include "RooStats/RooStatsUtils.h"
42 
43 using namespace std;
44 
46 
47 void RooStats::ProfileLikelihoodTestStat::SetAlwaysReuseNLL(Bool_t flag) { fgAlwaysReuseNll = flag ; }
48 
49 ////////////////////////////////////////////////////////////////////////////////
50 /// internal function to evaluate test statistics
51 /// can do depending on type:
52 /// - type = 0 standard evaluation,
53 /// - type = 1 find only unconditional NLL minimum,
54 /// - type = 2 conditional MLL
55 
57 
58  if( fDetailedOutputEnabled && fDetailedOutput ) {
59  delete fDetailedOutput;
60  fDetailedOutput = 0;
61  }
62  if( fDetailedOutputEnabled && !fDetailedOutput ) {
63  fDetailedOutput = new RooArgSet();
64  }
65 
66  //data.Print("V");
67 
68  TStopwatch tsw;
69  tsw.Start();
70 
71  double initial_mu_value = 0;
72  RooRealVar* firstPOI = dynamic_cast<RooRealVar*>( paramsOfInterest.first());
73  if (firstPOI) initial_mu_value = firstPOI->getVal();
74  //paramsOfInterest.getRealValue(firstPOI->GetName());
75  if (fPrintLevel > 1) {
76  cout << "POIs: " << endl;
77  paramsOfInterest.Print("v");
78  }
79 
82 
83  // simple
84  Bool_t reuse=(fReuseNll || fgAlwaysReuseNll) ;
85 
86  Bool_t created(kFALSE) ;
87  if (!reuse || fNll==0) {
88  RooArgSet* allParams = fPdf->getParameters(data);
90 
91  // need to call constrain for RooSimultaneous until stripDisconnected problem fixed
92  fNll = fPdf->createNLL(data, RooFit::CloneData(kFALSE),RooFit::Constrain(*allParams),RooFit::ConditionalObservables(fConditionalObs), RooFit::Offset(fLOffset));
93 
94  if (fPrintLevel > 0 && fLOffset) cout << "ProfileLikelihoodTestStat::Evaluate - Use Offset in creating NLL " << endl ;
95 
96  created = kTRUE ;
97  delete allParams;
98  if (fPrintLevel > 1) cout << "creating NLL " << fNll << " with data = " << &data << endl ;
99  }
100  if (reuse && !created) {
101  if (fPrintLevel > 1) cout << "reusing NLL " << fNll << " new data = " << &data << endl ;
102  fNll->setData(data,kFALSE) ;
103  }
104  // print data in case of number counting (simple data sets)
105  if (fPrintLevel > 1 && data.numEntries() == 1) {
106  std::cout << "Data set used is: ";
107  RooStats::PrintListContent(*data.get(0), std::cout);
108  }
109 
110 
111  // make sure we set the variables attached to this nll
112  RooArgSet* attachedSet = fNll->getVariables();
113 
114  *attachedSet = paramsOfInterest;
115  RooArgSet* origAttachedSet = (RooArgSet*) attachedSet->snapshot();
116 
117  ///////////////////////////////////////////////////////////////////////
118  // New profiling based on RooMinimizer (allows for Minuit2)
119  // based on major speed increases seen by CMS for complex problems
120 
121 
122  // other order
123  // get the numerator
124  RooArgSet* snap = (RooArgSet*)paramsOfInterest.snapshot();
125 
126  tsw.Stop();
127  double createTime = tsw.CpuTime();
128  tsw.Start();
129 
130  // get the denominator
131  double uncondML = 0;
132  double fit_favored_mu = 0;
133  int statusD = 0;
134  RooArgSet * detOutput = 0;
135  if (type != 2) {
136  // minimize and count eval errors
137  fNll->clearEvalErrorLog();
138  if (fPrintLevel>1) std::cout << "Do unconditional fit" << std::endl;
139  RooFitResult* result = GetMinNLL();
140  if (result) {
141  uncondML = result->minNll();
142  statusD = result->status();
143 
144  // get best fit value for one-sided interval
145  if (firstPOI) fit_favored_mu = attachedSet->getRealValue(firstPOI->GetName()) ;
146 
147  // save this snapshot
148  if( fDetailedOutputEnabled ) {
149  detOutput = DetailedOutputAggregator::GetAsArgSet(result, "fitUncond_", fDetailedOutputWithErrorsAndPulls);
150  fDetailedOutput->addOwned(*detOutput);
151  delete detOutput;
152  }
153  delete result;
154  }
155  else {
156  return TMath::SignalingNaN(); // this should not really happen
157  }
158  }
159  tsw.Stop();
160  double fitTime1 = tsw.CpuTime();
161 
162  //double ret = 0;
163  int statusN = 0;
164  tsw.Start();
165 
166  double condML = 0;
167 
168  bool doConditionalFit = (type != 1);
169 
170  // skip the conditional ML (the numerator) only when fit value is smaller than test value
171  if (!fSigned && type==0 &&
172  ((fLimitType==oneSided && fit_favored_mu >= initial_mu_value) ||
173  (fLimitType==oneSidedDiscovery && fit_favored_mu <= initial_mu_value))) {
174  doConditionalFit = false;
175  condML = uncondML;
176  }
177 
178  if (doConditionalFit) {
179 
180  if (fPrintLevel>1) std::cout << "Do conditional fit " << std::endl;
181 
182 
183  // cout <<" reestablish snapshot"<<endl;
184  *attachedSet = *snap;
185 
186 
187  // set the POI to constant
188  RooLinkedListIter it = paramsOfInterest.iterator();
189  RooRealVar* tmpPar = NULL, *tmpParA=NULL;
190  while((tmpPar = (RooRealVar*)it.Next())){
191  tmpParA = dynamic_cast<RooRealVar*>( attachedSet->find(tmpPar->GetName()));
192  if (tmpParA) tmpParA->setConstant();
193  }
194 
195 
196  // check if there are non-const parameters so it is worth to do the minimization
197  RooArgSet allParams(*attachedSet);
199 
200  // in case no nuisance parameters are present
201  // no need to minimize just evaluate the nll
202  if (allParams.getSize() == 0 ) {
203  // be sure to evaluate with offsets
204  if (fLOffset) RooAbsReal::setHideOffset(false);
205  condML = fNll->getVal();
206  if (fLOffset) RooAbsReal::setHideOffset(true);
207  }
208  else {
209  fNll->clearEvalErrorLog();
210  RooFitResult* result = GetMinNLL();
211  if (result) {
212  condML = result->minNll();
213  statusN = result->status();
214  if( fDetailedOutputEnabled ) {
215  detOutput = DetailedOutputAggregator::GetAsArgSet(result, "fitCond_", fDetailedOutputWithErrorsAndPulls);
216  fDetailedOutput->addOwned(*detOutput);
217  delete detOutput;
218  }
219  delete result;
220  }
221  else {
222  return TMath::SignalingNaN(); // this should not really happen
223  }
224  }
225 
226  }
227 
228  tsw.Stop();
229  double fitTime2 = tsw.CpuTime();
230 
231  double pll = 0;
232  if (type != 0) {
233  // for conditional only or unconditional fits
234  // need to compute nll value without the offset
235  if (fLOffset) {
237  pll = fNll->getVal();
238  }
239  else {
240  if (type == 1)
241  pll = uncondML;
242  else if (type == 2)
243  pll = condML;
244  }
245  }
246  else { // type == 0
247  // for standard profile likelihood evaluations
248  pll = condML-uncondML;
249 
250  if (fSigned) {
251  if (pll<0.0) {
252  if (fPrintLevel > 0) std::cout << "pll is negative - setting it to zero " << std::endl;
253  pll = 0.0; // bad fit
254  }
255  if (fLimitType==oneSidedDiscovery ? (fit_favored_mu < initial_mu_value)
256  : (fit_favored_mu > initial_mu_value))
257  pll = -pll;
258  }
259  }
260 
261  if (fPrintLevel > 0) {
262  std::cout << "EvaluateProfileLikelihood - ";
263  if (type <= 1)
264  std::cout << "mu hat = " << fit_favored_mu << ", uncond ML = " << uncondML;
265  if (type != 1)
266  std::cout << ", cond ML = " << condML;
267  if (type == 0)
268  std::cout << " pll = " << pll;
269  std::cout << " time (create/fit1/2) " << createTime << " , " << fitTime1 << " , " << fitTime2
270  << std::endl;
271  }
272 
273 
274  // need to restore the values ?
275  *attachedSet = *origAttachedSet;
276 
277  delete attachedSet;
278  delete origAttachedSet;
279  delete snap;
280 
281  if (!reuse) {
282  delete fNll;
283  fNll = 0;
284  }
285 
287 
288  if(statusN!=0 || statusD!=0) {
289  return -1; // indicate failed fit (WVE is not used anywhere yet)
290  }
291 
292  return pll;
293 
294  }
295 
296 ////////////////////////////////////////////////////////////////////////////////
297 /// find minimum of NLL using RooMinimizer
298 
300 
301  RooMinimizer minim(*fNll);
302  minim.setStrategy(fStrategy);
303  //LM: RooMinimizer.setPrintLevel has +1 offset - so subtract here -1 + an extra -1
304  int level = (fPrintLevel == 0) ? -1 : fPrintLevel -2;
305  minim.setPrintLevel(level);
306  minim.setEps(fTolerance);
307  // this causes a memory leak
308  minim.optimizeConst(2);
309  TString minimizer = fMinimizer;
311  if (algorithm == "Migrad") algorithm = "Minimize"; // prefer to use Minimize instead of Migrad
312  int status;
313  for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
314  status = minim.minimize(minimizer,algorithm);
315  if (status%1000 == 0) { // ignore erros from Improve
316  break;
317  } else if (tries < maxtries) {
318  cout << " ----> Doing a re-scan first" << endl;
319  minim.minimize(minimizer,"Scan");
320  if (tries == 2) {
321  if (fStrategy == 0 ) {
322  cout << " ----> trying with strategy = 1" << endl;;
323  minim.setStrategy(1);
324  }
325  else
326  tries++; // skip this trial if strategy is already 1
327  }
328  if (tries == 3) {
329  cout << " ----> trying with improve" << endl;;
330  minimizer = "Minuit";
331  algorithm = "migradimproved";
332  }
333  }
334  }
335 
336  //how to get cov quality faster?
337  return minim.save();
338  //minim.optimizeConst(false);
339 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
RooCmdArg Offset(Bool_t flag=kTRUE)
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void optimizeConst(Int_t flag)
If flag is true, perform constant term optimization on function being minimized.
virtual const RooArgSet * get() const
Definition: RooAbsData.h:77
RooCmdArg CloneData(Bool_t flag)
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
RooFit::MsgLevel globalKillBelow() const
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
Definition: TStopwatch.cxx:125
bool Bool_t
Definition: RtypesCore.h:59
virtual Double_t EvaluateProfileLikelihood(int type, RooAbsData &data, RooArgSet &paramsOfInterest)
internal function to evaluate test statistics can do depending on type:
static RooMsgService & instance()
Return reference to singleton instance.
STL namespace.
#define NULL
Definition: RtypesCore.h:88
static void setHideOffset(Bool_t flag)
Definition: RooAbsReal.cxx:116
void setEps(Double_t eps)
Change MINUIT epsilon.
void setStrategy(Int_t strat)
Change MINUIT strategy to istrat.
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
virtual TObject * Next()
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
Int_t getSize() const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
RooAbsArg * first() const
void setConstant(Bool_t value=kTRUE)
RooFitResult * save(const char *name=0, const char *title=0)
Save and return a RooFitResult snaphot of current minimizer status.
Double_t minNll() const
Definition: RooFitResult.h:96
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
void setGlobalKillBelow(RooFit::MsgLevel level)
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
Double_t SignalingNaN()
Definition: TMath.h:791
const Bool_t kFALSE
Definition: RtypesCore.h:92
static const std::string & DefaultMinimizerAlgo()
double Double_t
Definition: RtypesCore.h:55
int type
Definition: TGX11.cxx:120
Int_t minimize(const char *type, const char *alg=0)
Int_t setPrintLevel(Int_t newLevel)
Change the MINUIT internal printing level.
void RemoveConstantParameters(RooArgSet *set)
Definition: RooStatsUtils.h:62
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
Definition: RooMinimizer.h:38
RooLinkedListIter iterator(Bool_t dir=kIterForward) const
RooFitResult * GetMinNLL()
find minimum of NLL using RooMinimizer
double result[121]
RooCmdArg ConditionalObservables(const RooArgSet &set)
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooLinkedListIter is the TIterator implementation for RooLinkedList.
RooCmdArg Constrain(const RooArgSet &params)
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:269
Stopwatch class.
Definition: TStopwatch.h:28
Int_t status() const
Definition: RooFitResult.h:75