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
43
45
46////////////////////////////////////////////////////////////////////////////////
47/// internal function to evaluate test statistics
48/// can do depending on type:
49/// - type = 0 standard evaluation,
50/// - type = 1 find only unconditional NLL minimum,
51/// - type = 2 conditional MLL
52
54
55 if (fDetailedOutputEnabled) {
56 fDetailedOutput = std::make_unique<RooArgSet>();
57 }
58
59 //data.Print("V");
60
62 tsw.Start();
63
64 double initial_mu_value = 0;
65 RooRealVar* firstPOI = dynamic_cast<RooRealVar*>(paramsOfInterest.first());
66 if (firstPOI) initial_mu_value = firstPOI->getVal();
67 //paramsOfInterest.getRealValue(firstPOI->GetName());
68 if (fPrintLevel > 1) {
69 std::cout << "POIs: " << std::endl;
70 paramsOfInterest.Print("v");
71 }
72
74 if (fPrintLevel < 3) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
75
76 // simple
77 bool reuse=(fReuseNll || fgAlwaysReuseNll) ;
78
79 bool created(false) ;
80 if (!reuse || fNll==nullptr) {
81 std::unique_ptr<RooArgSet> allParams{fPdf->getParameters(data)};
83
84 // need to call constrain for RooSimultaneous until stripDisconnected problem fixed
85 fNll = std::unique_ptr<RooAbsReal>{fPdf->createNLL(data, RooFit::CloneData(false),RooFit::Constrain(*allParams),
86 RooFit::GlobalObservables(fGlobalObs), RooFit::ConditionalObservables(fConditionalObs), RooFit::Offset(fLOffset))};
87
88 if (fPrintLevel > 0) {
89 std::cout << "ProfileLikelihoodTestStat::Evaluate - Use Offset mode \""
90 << fLOffset << "\" in creating NLL" << std::endl;
91 }
92
93 created = true ;
94 if (fPrintLevel > 1) std::cout << "creating NLL " << &*fNll << " with data = " << &data << std::endl ;
95 }
96 if (reuse && !created) {
97 if (fPrintLevel > 1) std::cout << "reusing NLL " << &*fNll << " new data = " << &data << std::endl ;
98 fNll->setData(data,false) ;
99 }
100 // print data in case of number counting (simple data sets)
101 if (fPrintLevel > 1 && data.numEntries() == 1) {
102 std::cout << "Data set used is: ";
103 RooStats::PrintListContent(*data.get(0), std::cout);
104 }
105
106
107 // make sure we set the variables attached to this nll
108 std::unique_ptr<RooArgSet> attachedSet{fNll->getVariables()};
109
112
113 ///////////////////////////////////////////////////////////////////////
114 // New profiling based on RooMinimizer (allows for Minuit2)
115 // based on major speed increases seen by CMS for complex problems
116
117
118 // other order
119 // get the numerator
121
122 tsw.Stop();
123 double createTime = tsw.CpuTime();
124 tsw.Start();
125
126 // get the denominator
127 double uncondML = 0;
128 double fit_favored_mu = 0;
129 int statusD = 0;
130 RooArgSet * detOutput = nullptr;
131 if (type != 2) {
132 // minimize and count eval errors
133 fNll->clearEvalErrorLog();
134 if (fPrintLevel>1) std::cout << "Do unconditional fit" << std::endl;
135 std::unique_ptr<RooFitResult> result{GetMinNLL()};
136 if (result) {
137 uncondML = result->minNll();
138 statusD = result->status();
139
140 // get best fit value for one-sided interval
141 if (firstPOI) fit_favored_mu = attachedSet->getRealValue(firstPOI->GetName()) ;
142
143 // save this snapshot
144 if( fDetailedOutputEnabled ) {
145 detOutput = DetailedOutputAggregator::GetAsArgSet(result.get(), "fitUncond_", fDetailedOutputWithErrorsAndPulls);
146 fDetailedOutput->addOwned(*detOutput);
147 delete detOutput;
148 }
149 }
150 else {
151 return TMath::SignalingNaN(); // this should not really happen
152 }
153 }
154 tsw.Stop();
155 double fitTime1 = tsw.CpuTime();
156
157 //double ret = 0;
158 int statusN = 0;
159 tsw.Start();
160
161 double condML = 0;
162
163 bool doConditionalFit = (type != 1);
164
165 // skip the conditional ML (the numerator) only when fit value is smaller than test value
166 if (!fSigned && type==0 &&
167 ((fLimitType==oneSided && fit_favored_mu >= initial_mu_value) ||
168 (fLimitType==oneSidedDiscovery && fit_favored_mu <= initial_mu_value))) {
169 doConditionalFit = false;
171 }
172
173 if (doConditionalFit) {
174
175 if (fPrintLevel>1) std::cout << "Do conditional fit " << std::endl;
176
177
178 // std::cout <<" reestablish snapshot"<< std::endl;
179 attachedSet->assign(*snap);
180
181
182 // set the POI to constant
183 for (auto *tmpPar : paramsOfInterest) {
184 RooRealVar *tmpParA = dynamic_cast<RooRealVar *>(attachedSet->find(tmpPar->GetName()));
185 if (tmpParA) tmpParA->setConstant();
186 }
187
188
189 // check if there are non-const parameters so it is worth to do the minimization
190 RooArgSet allParams(*attachedSet);
192
193 // in case no nuisance parameters are present
194 // no need to minimize just evaluate the nll
195 if (allParams.empty() ) {
196 // be sure to evaluate with offsets
197 if (fLOffset == "initial") RooAbsReal::setHideOffset(false);
198 condML = fNll->getVal();
199 if (fLOffset == "initial") RooAbsReal::setHideOffset(true);
200 }
201 else {
202 fNll->clearEvalErrorLog();
203 std::unique_ptr<RooFitResult> result{GetMinNLL()};
204 if (result) {
205 condML = result->minNll();
206 statusN = result->status();
207 if( fDetailedOutputEnabled ) {
208 detOutput = DetailedOutputAggregator::GetAsArgSet(result.get(), "fitCond_", fDetailedOutputWithErrorsAndPulls);
209 fDetailedOutput->addOwned(*detOutput);
210 delete detOutput;
211 }
212 }
213 else {
214 return TMath::SignalingNaN(); // this should not really happen
215 }
216 }
217
218 }
219
220 tsw.Stop();
221 double fitTime2 = tsw.CpuTime();
222
223 double pll = 0;
224 if (type != 0) {
225 // for conditional only or unconditional fits
226 // need to compute nll value without the offset
227 if (fLOffset == "initial") {
229 pll = fNll->getVal();
230 }
231 else {
232 if (type == 1) {
233 pll = uncondML;
234 } else if (type == 2) {
235 pll = condML;
236 }
237 }
238 }
239 else { // type == 0
240 // for standard profile likelihood evaluations
241 pll = condML-uncondML;
242
243 if (fSigned) {
244 if (pll<0.0) {
245 if (fPrintLevel > 0) std::cout << "pll is negative - setting it to zero " << std::endl;
246 pll = 0.0; // bad fit
247 }
248 if (fLimitType==oneSidedDiscovery ? (fit_favored_mu < initial_mu_value)
250 pll = -pll;
251 }
252 }
253
254 if (fPrintLevel > 0) {
255 std::cout << "EvaluateProfileLikelihood - ";
256 if (type <= 1)
257 std::cout << "mu hat = " << fit_favored_mu << ", uncond ML = " << uncondML;
258 if (type != 1)
259 std::cout << ", cond ML = " << condML;
260 if (type == 0)
261 std::cout << " pll = " << pll;
262 std::cout << " time (create/fit1/2) " << createTime << " , " << fitTime1 << " , " << fitTime2
263 << std::endl;
264 }
265
266
267 // need to restore the values ?
269
270 delete origAttachedSet;
271 delete snap;
272
273 if (!reuse) {
274 fNll.reset();
275 }
276
277 RooMsgService::instance().setGlobalKillBelow(msglevel);
278
279 if(statusN!=0 || statusD!=0) {
280 return -1; // indicate failed fit (WVE is not used anywhere yet)
281 }
282
283 return pll;
284
285 }
286
287////////////////////////////////////////////////////////////////////////////////
288/// find minimum of NLL using RooMinimizer
289
290std::unique_ptr<RooFitResult> RooStats::ProfileLikelihoodTestStat::GetMinNLL() {
291
292 const auto& config = GetGlobalRooStatsConfig();
293 RooMinimizer minim(*fNll);
294 minim.setStrategy(fStrategy);
295 minim.setEvalErrorWall(config.useEvalErrorWall);
296 //LM: RooMinimizer.setPrintLevel has +1 offset - so subtract here -1 + an extra -1
297 int level = (fPrintLevel == 0) ? -1 : fPrintLevel -2;
298 minim.setPrintLevel(level);
299 minim.setEps(fTolerance);
300 // this causes a memory leak
301 minim.optimizeConst(2);
302 TString minimizer = fMinimizer;
304 if (algorithm == "Migrad") algorithm = "Minimize"; // prefer to use Minimize instead of Migrad
305 int status;
306 for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
307 status = minim.minimize(minimizer,algorithm);
308 if (status%1000 == 0) { // ignore errors from Improve
309 break;
310 } else if (tries < maxtries) {
311 std::cout << " ----> Doing a re-scan first" << std::endl;
312 minim.minimize(minimizer,"Scan");
313 if (tries == 2) {
314 if (fStrategy == 0 ) {
315 std::cout << " ----> trying with strategy = 1" << std::endl;
316 minim.setStrategy(1);
317 }
318 else
319 tries++; // skip this trial if strategy is already 1
320 }
321 if (tries == 3) {
322 std::cout << " ----> trying with improve" << std::endl;
323 minimizer = "Minuit";
324 algorithm = "migradimproved";
325 }
326 }
327 }
328
329 //how to get cov quality faster?
330 return std::unique_ptr<RooFitResult>{minim.save()};
331 //minim.optimizeConst(false);
332}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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()
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
static void setHideOffset(bool flag)
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Wrapper class around ROOT::Math::Minimizer that provides a seamless interface between the minimizer f...
static RooMsgService & instance()
Return reference to singleton instance.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
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
Stopwatch class.
Definition TStopwatch.h:28
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:914