74 ClassImp(RooStats::ProfileLikelihoodCalculator); ;
76 using namespace RooFit;
77 using namespace RooStats;
83 ProfileLikelihoodCalculator::ProfileLikelihoodCalculator() :
84 CombinedCalculator(), fFitResult(0), fGlobalFitDone(false)
88 ProfileLikelihoodCalculator::ProfileLikelihoodCalculator(RooAbsData& data, RooAbsPdf& pdf,
const RooArgSet& paramsOfInterest,
89 Double_t size,
const RooArgSet* nullParams ) :
90 CombinedCalculator(data,pdf, paramsOfInterest, size, nullParams ),
91 fFitResult(0), fGlobalFitDone(false)
97 ProfileLikelihoodCalculator::ProfileLikelihoodCalculator(RooAbsData& data, ModelConfig& model, Double_t size) :
98 CombinedCalculator(data, model, size),
99 fFitResult(0), fGlobalFitDone(false)
103 assert(model.GetPdf() );
113 ProfileLikelihoodCalculator::~ProfileLikelihoodCalculator(){
114 if (fFitResult)
delete fFitResult;
117 void ProfileLikelihoodCalculator::DoReset()
const {
120 if (fFitResult)
delete fFitResult;
124 RooAbsReal * ProfileLikelihoodCalculator::DoGlobalFit()
const {
130 RooAbsPdf * pdf = GetPdf();
131 RooAbsData* data = GetData();
132 if (!data || !pdf )
return 0;
135 RooArgSet* constrainedParams = pdf->getParameters(*data);
136 if (!constrainedParams)
return 0;
137 RemoveConstantParameters(constrainedParams);
139 const auto& config = GetGlobalRooStatsConfig();
140 RooAbsReal * nll = pdf->createNLL(*data, CloneData(
true), Constrain(*constrainedParams),ConditionalObservables(fConditionalObs), GlobalObservables(fGlobalObs),
141 RooFit::Offset(config.useLikelihoodOffset) );
144 if (fFitResult && fGlobalFitDone) {
145 delete constrainedParams;
150 oocoutP((TObject*)0,Minimization) <<
"ProfileLikelihoodCalcultor::DoGLobalFit - find MLE " << std::endl;
152 if (fFitResult)
delete fFitResult;
153 fFitResult = DoMinimizeNLL(nll);
157 fFitResult->printStream( oocoutI((TObject*)0,Minimization), fFitResult->defaultPrintContents(0), fFitResult->defaultPrintStyle(0) );
159 if (fFitResult->status() != 0)
160 oocoutW((TObject*)0,Minimization) <<
"ProfileLikelihoodCalcultor::DoGlobalFit - Global fit failed - status = " << fFitResult->status() << std::endl;
162 fGlobalFitDone =
true;
165 delete constrainedParams;
169 RooFitResult * ProfileLikelihoodCalculator::DoMinimizeNLL(RooAbsReal * nll) {
172 const char * minimType = ROOT::Math::MinimizerOptions::DefaultMinimizerType().c_str();
173 const char * minimAlgo = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo().c_str();
174 int strategy = ROOT::Math::MinimizerOptions::DefaultStrategy();
175 int level = ROOT::Math::MinimizerOptions::DefaultPrintLevel() -1;
176 int tolerance = ROOT::Math::MinimizerOptions::DefaultTolerance();
177 oocoutP((TObject*)0,Minimization) <<
"ProfileLikelihoodCalcultor::DoMinimizeNLL - using " << minimType <<
" / " << minimAlgo <<
" with strategy " << strategy << std::endl;
180 const auto& config = GetGlobalRooStatsConfig();
182 RooMinimizer minim(*nll);
183 minim.setStrategy(strategy);
184 minim.setEps(tolerance);
185 minim.setPrintLevel(level);
186 minim.optimizeConst(2);
187 minim.setEvalErrorWall(config.useEvalErrorWall);
190 for (
int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
191 status = minim.minimize(minimType,minimAlgo);
192 if (status%1000 == 0) {
194 }
else if (tries < maxtries) {
195 cout <<
" ----> Doing a re-scan first" << endl;
196 minim.minimize(minimType,
"Scan");
198 if (strategy == 0 ) {
199 cout <<
" ----> trying with strategy = 1" << endl;;
200 minim.setStrategy(1);
206 cout <<
" ----> trying with improve" << endl;;
207 minimType =
"Minuit";
208 minimAlgo =
"migradimproved";
213 RooFitResult * result = minim.save();
223 LikelihoodInterval* ProfileLikelihoodCalculator::GetInterval()
const {
226 RooAbsPdf * pdf = GetPdf();
227 RooAbsData* data = GetData();
228 if (!data || !pdf || fPOI.getSize() == 0)
return 0;
230 RooArgSet* constrainedParams = pdf->getParameters(*data);
231 RemoveConstantParameters(constrainedParams);
241 RooAbsReal * nll = DoGlobalFit();
249 RooAbsReal* profile = nll->createProfile(fPOI);
250 profile->addOwnedComponents(*nll) ;
254 const RooArgList & fitParams = fFitResult->floatParsFinal();
255 for (
int i = 0; i < fitParams.getSize(); ++i) {
256 RooRealVar & fitPar = (RooRealVar &) fitParams[i];
257 RooRealVar * par = (RooRealVar*) fPOI.find( fitPar.GetName() );
259 par->setVal( fitPar.getVal() );
260 par->setError( fitPar.getError() );
271 TString name = TString(
"LikelihoodInterval_");
275 TIter iter = fPOI.createIterator();
276 RooArgSet fitParSet(fitParams);
277 RooArgSet * bestPOI =
new RooArgSet();
278 while (RooAbsArg * arg = (RooAbsArg*) iter.Next() ) {
279 RooAbsArg * p = fitParSet.find( arg->GetName() );
280 if (p) bestPOI->addClone(*p);
281 else bestPOI->addClone(*arg);
285 LikelihoodInterval* interval =
new LikelihoodInterval(name, profile, &fPOI, bestPOI);
286 interval->SetConfidenceLevel(1.-fSize);
287 delete constrainedParams;
301 HypoTestResult* ProfileLikelihoodCalculator::GetHypoTest()
const {
304 RooAbsPdf * pdf = GetPdf();
305 RooAbsData* data = GetData();
308 if (!data || !pdf)
return 0;
310 if (fNullParams.getSize() == 0)
return 0;
315 poiList.addClone(fNullParams);
319 RooAbsReal * nll = DoGlobalFit();
327 RooArgSet* constrainedParams = pdf->getParameters(*data);
328 RemoveConstantParameters(constrainedParams);
330 Double_t nLLatMLE = fFitResult->minNll();
332 Double_t nlloffset = (RooStats::IsNLLOffset() ) ? nll->getVal() - nLLatMLE : 0;
335 std::vector<double> oldValues(poiList.getSize() );
336 for (
unsigned int i = 0; i < oldValues.size(); ++i) {
337 RooRealVar * mytarget = (RooRealVar*) constrainedParams->find(poiList[i].GetName());
339 oldValues[i] = mytarget->getVal();
340 mytarget->setVal( ( (RooRealVar&) poiList[i] ).getVal() );
341 mytarget->setConstant(kTRUE);
350 RooArgSet nuisParams(*constrainedParams);
353 RemoveConstantParameters(&nuisParams);
356 bool existVarParams =
false;
357 TIter it = nuisParams.createIterator();
358 RooRealVar * myarg = 0;
359 while ((myarg = (RooRealVar *)it.Next())) {
360 if ( !myarg->isConstant() ) {
361 existVarParams =
true;
366 Double_t nLLatCondMLE = nLLatMLE;
367 if (existVarParams) {
368 oocoutP((TObject*)0,Minimization) <<
"ProfileLikelihoodCalcultor::GetHypoTest - do conditional fit " << std::endl;
370 RooFitResult * fit2 = DoMinimizeNLL(nll);
374 nLLatCondMLE = fit2->minNll();
375 fit2->printStream( oocoutI((TObject*)0,Minimization), fit2->defaultPrintContents(0), fit2->defaultPrintStyle(0) );
377 if (fit2->status() != 0)
378 oocoutW((TObject*)0,Minimization) <<
"ProfileLikelihoodCalcultor::GetHypotest - Conditional fit failed - status = " << fit2->status() << std::endl;
384 nLLatCondMLE = nll->getVal();
386 if (RooStats::IsNLLOffset() ) nLLatCondMLE -= nlloffset;
390 Double_t deltaNLL = std::max( nLLatCondMLE-nLLatMLE, 0.);
393 RemoveConstantParameters(poiList);
394 int ndf = poiList.getSize();
396 Double_t pvalue = ROOT::Math::chisquared_cdf_c( 2* deltaNLL, ndf);
399 if (ndf == 1) pvalue = 0.5 * pvalue;
401 TString name = TString(
"ProfileLRHypoTestResult_");
402 HypoTestResult* htr =
new HypoTestResult(name, pvalue, 0 );
405 for (
unsigned int i = 0; i < oldValues.size(); ++i) {
406 RooRealVar * mytarget = (RooRealVar*) constrainedParams->find(poiList[i].GetName());
408 mytarget->setVal(oldValues[i] );
409 mytarget->setConstant(
false);
413 delete constrainedParams;