Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
15ProfileLikelihoodTestStat is an implementation of the TestStatistic interface
16that calculates the profile likelihood ratio at a particular parameter point
17given a dataset. It does not constitute a statistical test, for that one may
18either use:
19
20 - the ProfileLikelihoodCalculator that relies on asymptotic properties of the
21 Profile Likelihood Ratio
22 - the NeymanConstruction class with this class as a test statistic
23 - the HybridCalculator 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 "RooMsgService.h"
35#include "RooMinimizer.h"
36#include "RooArgSet.h"
37#include "RooDataSet.h"
38#include "TStopwatch.h"
39
41
42using namespace std;
43
45
47
48////////////////////////////////////////////////////////////////////////////////
49/// internal function to evaluate test statistics
50/// can do depending on type:
51/// - type = 0 standard evaluation,
52/// - type = 1 find only unconditional NLL minimum,
53/// - type = 2 conditional MLL
54
56
57 if( fDetailedOutputEnabled && fDetailedOutput ) {
58 delete fDetailedOutput;
59 fDetailedOutput = 0;
60 }
61 if( fDetailedOutputEnabled && !fDetailedOutput ) {
62 fDetailedOutput = new RooArgSet();
63 }
64
65 //data.Print("V");
66
67 TStopwatch tsw;
68 tsw.Start();
69
70 double initial_mu_value = 0;
71 RooRealVar* firstPOI = dynamic_cast<RooRealVar*>( paramsOfInterest.first());
72 if (firstPOI) initial_mu_value = firstPOI->getVal();
73 //paramsOfInterest.getRealValue(firstPOI->GetName());
74 if (fPrintLevel > 1) {
75 cout << "POIs: " << endl;
76 paramsOfInterest.Print("v");
77 }
78
81
82 // simple
83 bool reuse=(fReuseNll || fgAlwaysReuseNll) ;
84
85 bool created(false) ;
86 if (!reuse || fNll==nullptr) {
87 std::unique_ptr<RooArgSet> allParams{fPdf->getParameters(data)};
89
90 // need to call constrain for RooSimultaneous until stripDisconnected problem fixed
91 fNll = std::unique_ptr<RooAbsReal>{fPdf->createNLL(data, RooFit::CloneData(false),RooFit::Constrain(*allParams),
92 RooFit::GlobalObservables(fGlobalObs), RooFit::ConditionalObservables(fConditionalObs), RooFit::Offset(fLOffset))};
93
94 if (fPrintLevel > 0 && fLOffset) cout << "ProfileLikelihoodTestStat::Evaluate - Use Offset in creating NLL " << endl ;
95
96 created = true ;
97 if (fPrintLevel > 1) cout << "creating NLL " << &*fNll << " with data = " << &data << endl ;
98 }
99 if (reuse && !created) {
100 if (fPrintLevel > 1) cout << "reusing NLL " << &*fNll << " new data = " << &data << endl ;
101 fNll->setData(data,false) ;
102 }
103 // print data in case of number counting (simple data sets)
104 if (fPrintLevel > 1 && data.numEntries() == 1) {
105 std::cout << "Data set used is: ";
106 RooStats::PrintListContent(*data.get(0), std::cout);
107 }
108
109
110 // make sure we set the variables attached to this nll
111 std::unique_ptr<RooArgSet> attachedSet{fNll->getVariables()};
112
113 attachedSet->assign(paramsOfInterest);
114 RooArgSet* origAttachedSet = (RooArgSet*) attachedSet->snapshot();
115
116 ///////////////////////////////////////////////////////////////////////
117 // New profiling based on RooMinimizer (allows for Minuit2)
118 // based on major speed increases seen by CMS for complex problems
119
120
121 // other order
122 // get the numerator
123 RooArgSet* snap = (RooArgSet*)paramsOfInterest.snapshot();
124
125 tsw.Stop();
126 double createTime = tsw.CpuTime();
127 tsw.Start();
128
129 // get the denominator
130 double uncondML = 0;
131 double fit_favored_mu = 0;
132 int statusD = 0;
133 RooArgSet * detOutput = 0;
134 if (type != 2) {
135 // minimize and count eval errors
136 fNll->clearEvalErrorLog();
137 if (fPrintLevel>1) std::cout << "Do unconditional fit" << std::endl;
138 std::unique_ptr<RooFitResult> result{GetMinNLL()};
139 if (result) {
140 uncondML = result->minNll();
141 statusD = result->status();
142
143 // get best fit value for one-sided interval
144 if (firstPOI) fit_favored_mu = attachedSet->getRealValue(firstPOI->GetName()) ;
145
146 // save this snapshot
147 if( fDetailedOutputEnabled ) {
148 detOutput = DetailedOutputAggregator::GetAsArgSet(result.get(), "fitUncond_", fDetailedOutputWithErrorsAndPulls);
149 fDetailedOutput->addOwned(*detOutput);
150 delete detOutput;
151 }
152 }
153 else {
154 return TMath::SignalingNaN(); // this should not really happen
155 }
156 }
157 tsw.Stop();
158 double fitTime1 = tsw.CpuTime();
159
160 //double ret = 0;
161 int statusN = 0;
162 tsw.Start();
163
164 double condML = 0;
165
166 bool doConditionalFit = (type != 1);
167
168 // skip the conditional ML (the numerator) only when fit value is smaller than test value
169 if (!fSigned && type==0 &&
170 ((fLimitType==oneSided && fit_favored_mu >= initial_mu_value) ||
171 (fLimitType==oneSidedDiscovery && fit_favored_mu <= initial_mu_value))) {
172 doConditionalFit = false;
173 condML = uncondML;
174 }
175
176 if (doConditionalFit) {
177
178 if (fPrintLevel>1) std::cout << "Do conditional fit " << std::endl;
179
180
181 // cout <<" reestablish snapshot"<<endl;
182 attachedSet->assign(*snap);
183
184
185 // set the POI to constant
186 for (auto *tmpPar : paramsOfInterest) {
187 RooRealVar *tmpParA = dynamic_cast<RooRealVar *>(attachedSet->find(tmpPar->GetName()));
188 if (tmpParA) tmpParA->setConstant();
189 }
190
191
192 // check if there are non-const parameters so it is worth to do the minimization
193 RooArgSet allParams(*attachedSet);
195
196 // in case no nuisance parameters are present
197 // no need to minimize just evaluate the nll
198 if (allParams.empty() ) {
199 // be sure to evaluate with offsets
200 if (fLOffset) RooAbsReal::setHideOffset(false);
201 condML = fNll->getVal();
202 if (fLOffset) RooAbsReal::setHideOffset(true);
203 }
204 else {
205 fNll->clearEvalErrorLog();
206 std::unique_ptr<RooFitResult> result{GetMinNLL()};
207 if (result) {
208 condML = result->minNll();
209 statusN = result->status();
210 if( fDetailedOutputEnabled ) {
211 detOutput = DetailedOutputAggregator::GetAsArgSet(result.get(), "fitCond_", fDetailedOutputWithErrorsAndPulls);
212 fDetailedOutput->addOwned(*detOutput);
213 delete detOutput;
214 }
215 }
216 else {
217 return TMath::SignalingNaN(); // this should not really happen
218 }
219 }
220
221 }
222
223 tsw.Stop();
224 double fitTime2 = tsw.CpuTime();
225
226 double pll = 0;
227 if (type != 0) {
228 // for conditional only or unconditional fits
229 // need to compute nll value without the offset
230 if (fLOffset) {
232 pll = fNll->getVal();
233 }
234 else {
235 if (type == 1)
236 pll = uncondML;
237 else if (type == 2)
238 pll = condML;
239 }
240 }
241 else { // type == 0
242 // for standard profile likelihood evaluations
243 pll = condML-uncondML;
244
245 if (fSigned) {
246 if (pll<0.0) {
247 if (fPrintLevel > 0) std::cout << "pll is negative - setting it to zero " << std::endl;
248 pll = 0.0; // bad fit
249 }
250 if (fLimitType==oneSidedDiscovery ? (fit_favored_mu < initial_mu_value)
251 : (fit_favored_mu > initial_mu_value))
252 pll = -pll;
253 }
254 }
255
256 if (fPrintLevel > 0) {
257 std::cout << "EvaluateProfileLikelihood - ";
258 if (type <= 1)
259 std::cout << "mu hat = " << fit_favored_mu << ", uncond ML = " << uncondML;
260 if (type != 1)
261 std::cout << ", cond ML = " << condML;
262 if (type == 0)
263 std::cout << " pll = " << pll;
264 std::cout << " time (create/fit1/2) " << createTime << " , " << fitTime1 << " , " << fitTime2
265 << std::endl;
266 }
267
268
269 // need to restore the values ?
270 attachedSet->assign(*origAttachedSet);
271
272 delete origAttachedSet;
273 delete snap;
274
275 if (!reuse) {
276 fNll.reset();
277 }
278
280
281 if(statusN!=0 || statusD!=0) {
282 return -1; // indicate failed fit (WVE is not used anywhere yet)
283 }
284
285 return pll;
286
287 }
288
289////////////////////////////////////////////////////////////////////////////////
290/// find minimum of NLL using RooMinimizer
291
292std::unique_ptr<RooFitResult> RooStats::ProfileLikelihoodTestStat::GetMinNLL() {
293
294 const auto& config = GetGlobalRooStatsConfig();
295 RooMinimizer minim(*fNll);
296 minim.setStrategy(fStrategy);
297 minim.setEvalErrorWall(config.useEvalErrorWall);
298 //LM: RooMinimizer.setPrintLevel has +1 offset - so subtract here -1 + an extra -1
299 int level = (fPrintLevel == 0) ? -1 : fPrintLevel -2;
300 minim.setPrintLevel(level);
301 minim.setEps(fTolerance);
302 // this causes a memory leak
303 minim.optimizeConst(2);
304 TString minimizer = fMinimizer;
306 if (algorithm == "Migrad") algorithm = "Minimize"; // prefer to use Minimize instead of Migrad
307 int status;
308 for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
309 status = minim.minimize(minimizer,algorithm);
310 if (status%1000 == 0) { // ignore erros from Improve
311 break;
312 } else if (tries < maxtries) {
313 cout << " ----> Doing a re-scan first" << endl;
314 minim.minimize(minimizer,"Scan");
315 if (tries == 2) {
316 if (fStrategy == 0 ) {
317 cout << " ----> trying with strategy = 1" << endl;;
318 minim.setStrategy(1);
319 }
320 else
321 tries++; // skip this trial if strategy is already 1
322 }
323 if (tries == 3) {
324 cout << " ----> trying with improve" << endl;;
325 minimizer = "Minuit";
326 algorithm = "migradimproved";
327 }
328 }
329 }
330
331 //how to get cov quality faster?
332 return std::unique_ptr<RooFitResult>{minim.save()};
333 //minim.optimizeConst(false);
334}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
static const std::string & DefaultMinimizerAlgo()
RooAbsArg * first() const
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:59
void setConstant(bool value=true)
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
static void setHideOffset(bool flag)
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
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.
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.
static RooMsgService & instance()
Return reference to singleton instance.
void setGlobalKillBelow(RooFit::MsgLevel level)
RooFit::MsgLevel globalKillBelow() const
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
static RooArgSet * GetAsArgSet(RooFitResult *result, TString prefix="", bool withErrorsAndPulls=false)
Translate the given fit result to a RooArgSet in a generic way.
virtual double EvaluateProfileLikelihood(int type, RooAbsData &data, RooArgSet &paramsOfInterest)
evaluate the profile likelihood ratio (type = 0) or the minimum of likelihood (type=1) or the conditi...
std::unique_ptr< RooFitResult > GetMinNLL()
find minimum of NLL using RooMinimizer
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
void Stop()
Stop the stopwatch.
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.
MsgLevel
Verbosity level for RooMsgService::StreamConfig in RooMsgService.
void RemoveConstantParameters(RooArgSet *set)
RooStatsConfig & GetGlobalRooStatsConfig()
Retrieve the config object which can be used to set flags for things like offsetting the likelihood o...
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
useful function to print in one line the content of a set with their values
Double_t SignalingNaN()
Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).
Definition TMath.h:908