Logo ROOT   6.16/01
Reference Guide
AsymptoticCalculator.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, Sven Kreiss 23/05/10
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11/** \class RooStats::AsymptoticCalculator
12 \ingroup Roostats
13
14Hypothesis Test Calculator based on the asymptotic formulae for the profile
15likelihood ratio.
16
17 Performs hypothesis tests using asymptotic formula for the profile likelihood and
18Asimov data set
19
20See G. Cowan, K. Cranmer, E. Gross and O. Vitells: Asymptotic formulae for
21likelihood- based tests of new physics. Eur. Phys. J., C71:1–19, 2011.
22It provides method to perform an hypothesis tests using the likelihood function
23and computes the p values for the null and the alternate using the asymptotic
24formulae for the profile likelihood ratio described in the given paper.
25
26The calculator provides methods to produce the Asimov dataset, i.e a dataset
27generated where the observed values are equal to the expected ones.
28The Asimov data set is then used to compute the observed asymptotic p-value for
29the alternate hypothesis and the asymptotic expected p-values.
30
31The asymptotic formulae are valid only for one POI (parameter of interest). So
32the calculator works only for one-dimensional (one POI) model.
33If more than one POI exists consider as POI only the first one is used.
34
35*/
36
37
43
44#include "RooArgSet.h"
45#include "RooArgList.h"
46#include "RooProdPdf.h"
47#include "RooSimultaneous.h"
48#include "RooDataSet.h"
49#include "RooCategory.h"
50#include "RooRealVar.h"
51#include "RooMinimizer.h"
52#include "RooFitResult.h"
53#include "RooNLLVar.h"
55#include "RooPoisson.h"
56#include "RooUniform.h"
57#include "RooGamma.h"
58#include "RooGaussian.h"
59#include "RooBifurGauss.h"
60#include "RooLognormal.h"
61#include "RooDataHist.h"
62#include <cmath>
63#include <typeinfo>
64
67
68#include "TStopwatch.h"
69
70using namespace RooStats;
71using namespace std;
72
73
75
76int AsymptoticCalculator::fgPrintLevel = 1;
77
78////////////////////////////////////////////////////////////////////////////////
79/// set print level (static function)
80///
81/// - 0 minimal,
82/// - 1 normal,
83/// - 2 debug
84
85void AsymptoticCalculator::SetPrintLevel(int level) {
86 fgPrintLevel = level;
87}
88
89////////////////////////////////////////////////////////////////////////////////
90/// constructor for asymptotic calculator from Data set and ModelConfig
91
94 const ModelConfig &altModel,
95 const ModelConfig &nullModel, bool nominalAsimov) :
96 HypoTestCalculatorGeneric(data, altModel, nullModel, 0),
97 fOneSided(false), fOneSidedDiscovery(false), fNominalAsimov(nominalAsimov),
98 fUseQTilde(-1),
99 fNLLObs(0), fNLLAsimov(0),
100 fAsimovData(0)
101{
102 if (!Initialize()) return;
103
104 int verbose = fgPrintLevel;
105 // try to guess default configuration
106 // (this part should be only in constructor because the null snapshot might change during HypoTestInversion
107 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
108 assert(nullSnapshot);
109 RooRealVar * muNull = dynamic_cast<RooRealVar*>( nullSnapshot->first() );
110 assert(muNull);
111 if (muNull->getVal() == muNull->getMin()) {
112 fOneSidedDiscovery = true;
113 if (verbose > 0)
114 oocoutI((TObject*)0,InputArguments) << "AsymptotiCalculator: Minimum of POI is " << muNull->getMin() << " corresponds to null snapshot - default configuration is one-sided discovery formulae " << std::endl;
115 }
116
117}
118
119////////////////////////////////////////////////////////////////////////////////
120/// Initialize the calculator
121/// The initialization will perform a global fit of the model to the data
122/// and build an Asimov data set.
123/// It will then also fit the model to the Asimov data set to find the likelihood value
124/// of the Asimov data set
125/// nominalAsimov is an option for using Asimov data set obtained using nominal nuisance parameter values
126/// By default the nuisance parameters are fitted to the data
127/// NOTE: If a fit has been done before, one for speeding up could set all the initial parameters
128/// to the fit value and in addition set the null snapshot to the best fit
129
131
132 int verbose = fgPrintLevel;
133 if (verbose >= 0)
134 oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::Initialize...." << std::endl;
135
136
137 RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
138 if (!nullPdf) {
139 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has not a pdf defined" << std::endl;
140 return false;
141 }
142 RooAbsData * obsData = const_cast<RooAbsData *>(GetData() );
143 if (!obsData ) {
144 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - data set has not been defined" << std::endl;
145 return false;
146 }
147 RooAbsData & data = *obsData;
148
149
150
152 if (!poi || poi->getSize() == 0) {
153 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has not POI defined." << endl;
154 return false;
155 }
156 if (poi->getSize() > 1) {
157 oocoutW((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has more than one POI defined \n\t"
158 << "The asymptotic calculator works for only one POI - consider as POI only the first parameter"
159 << std::endl;
160 }
161
162
163 // This will set the poi value to the null snapshot value in the ModelConfig
164 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
165 if(nullSnapshot == NULL || nullSnapshot->getSize() == 0) {
166 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - Null model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << endl;
167 return false;
168 }
169
170 // GetNullModel()->Print();
171 // printf("ASymptotic calc: null snapshot\n");
172 // nullSnapshot->Print("v");
173 // printf("PDF variables " );
174 // nullPdf->getVariables()->Print("v");
175
176 // keep snapshot for the initial parameter values (need for nominal Asimov)
177 RooArgSet nominalParams;
178 RooArgSet * allParams = nullPdf->getParameters(data);
179 RemoveConstantParameters(allParams);
180 if (fNominalAsimov) {
181 allParams->snapshot(nominalParams);
182 }
186
187 // evaluate the unconditional nll for the full model on the observed data
188 if (verbose >= 0)
189 oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::Initialize - Find best unconditional NLL on observed data" << endl;
190 fNLLObs = EvaluateNLL( *nullPdf, data, GetNullModel()->GetConditionalObservables(),GetNullModel()->GetGlobalObservables());
191 // fill also snapshot of best poi
192 poi->snapshot(fBestFitPoi);
193 RooRealVar * muBest = dynamic_cast<RooRealVar*>(fBestFitPoi.first());
194 assert(muBest);
195 if (verbose >= 0)
196 oocoutP((TObject*)0,Eval) << "Best fitted POI value = " << muBest->getVal() << " +/- " << muBest->getError() << std::endl;
197 // keep snapshot of all best fit parameters
198 allParams->snapshot(fBestFitParams);
199 delete allParams;
200
201 // compute Asimov data set for the background (alt poi ) value
202 const RooArgSet * altSnapshot = GetAlternateModel()->GetSnapshot();
203 if(altSnapshot == NULL || altSnapshot->getSize() == 0) {
204 oocoutE((TObject*)0,InputArguments) << "Alt (Background) model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << endl;
205 return false;
206 }
207
208 RooArgSet poiAlt(*altSnapshot); // this is the poi snapshot of B (i.e. for mu=0)
209
210 oocoutP((TObject*)0,Eval) << "AsymptoticCalculator: Building Asimov data Set" << endl;
211
212 // check that in case of binned models the n number of bins of the observables are consistent
213 // with the number of bins in the observed data
214 // This number will be used for making the Asimov data set so it will be more consistent with the
215 // observed data
216 int prevBins = 0;
217 RooRealVar * xobs = 0;
218 if (GetNullModel()->GetObservables() && GetNullModel()->GetObservables()->getSize() == 1 ) {
219 xobs = (RooRealVar*) (GetNullModel()->GetObservables())->first();
220 if (data.IsA() == RooDataHist::Class() ) {
221 if (data.numEntries() != xobs->getBins() ) {
222 prevBins = xobs->getBins();
223 oocoutW((TObject*)0,InputArguments) << "AsymptoticCalculator: number of bins in " << xobs->GetName() << " are different than data bins "
224 << " set the same data bins " << data.numEntries() << " in range "
225 << " [ " << xobs->getMin() << " , " << xobs->getMax() << " ]" << std::endl;
226 xobs->setBins(data.numEntries());
227 }
228 }
229 }
230
231 if (!fNominalAsimov) {
232 if (verbose >= 0)
233 oocoutI((TObject*)0,InputArguments) << "AsymptoticCalculator: Asimov data will be generated using fitted nuisance parameter values" << endl;
234 RooArgSet * tmp = (RooArgSet*) poiAlt.snapshot();
236 }
237
238 else {
239 // assume use current value of nuisance as nominal ones
240 if (verbose >= 0)
241 oocoutI((TObject*)0,InputArguments) << "AsymptoticCalculator: Asimovdata set will be generated using nominal (current) nuisance parameter values" << endl;
242 nominalParams = poiAlt; // set poi to alt value but keep nuisance at the nominal one
244 }
245
246 if (!fAsimovData) {
247 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator: Error : Asimov data set could not be generated " << endl;
248 return false;
249 }
250
251 // set global observables to their Asimov values
252 RooArgSet globObs;
253 RooArgSet globObsSnapshot;
254 if (GetNullModel()->GetGlobalObservables() ) {
255 globObs.add(*GetNullModel()->GetGlobalObservables());
256 assert(globObs.getSize() == fAsimovGlobObs.getSize() );
257 // store previous snapshot value
258 globObs.snapshot(globObsSnapshot);
259 globObs = fAsimovGlobObs;
260 }
261
262
263 // evaluate the likelihood. Since we use on Asimov data , conditional and unconditional values should be the same
264 // do conditional fit since is faster
265
266 RooRealVar * muAlt = (RooRealVar*) poiAlt.first();
267 assert(muAlt);
268 if (verbose>=0)
269 oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::Initialize Find best conditional NLL on ASIMOV data set for given alt POI ( " <<
270 muAlt->GetName() << " ) = " << muAlt->getVal() << std::endl;
271
272 fNLLAsimov = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables(), &poiAlt );
273 // for unconditional fit
274 //fNLLAsimov = EvaluateNLL( *nullPdf, *fAsimovData);
275 //poi->Print("v");
276
277 // restore previous value
278 globObs = globObsSnapshot;
279
280 // restore number of bins
281 if (prevBins > 0 && xobs) xobs->setBins(prevBins);
282
283 fIsInitialized = true;
284 return true;
285}
286
287////////////////////////////////////////////////////////////////////////////////
288
289Double_t AsymptoticCalculator::EvaluateNLL(RooAbsPdf & pdf, RooAbsData& data, const RooArgSet * condObs, const RooArgSet * globObs, const RooArgSet *poiSet) {
290 int verbose = fgPrintLevel;
291
292
295
296
297 RooArgSet* allParams = pdf.getParameters(data);
299 // add constraint terms for all non-constant parameters
300
301 RooArgSet conditionalObs;
302 if (condObs) conditionalObs.add(*condObs);
303 RooArgSet globalObs;
304 if (globObs) globalObs.add(*globObs);
305
306 // need to call constrain for RooSimultaneous until stripDisconnected problem fixed
308
309 RooArgSet* attachedSet = nll->getVariables();
310
311 // if poi are specified - do a conditional fit
312 RooArgSet paramsSetConstant;
313 // support now only one POI
314 if (poiSet && poiSet->getSize() > 0) {
315 RooRealVar * muTest = (RooRealVar*) (poiSet->first());
316 RooRealVar * poiVar = dynamic_cast<RooRealVar*>( attachedSet->find( muTest->GetName() ) );
317 if (poiVar && !poiVar->isConstant() ) {
318 poiVar->setVal( muTest->getVal() );
319 poiVar->setConstant();
320 paramsSetConstant.add(*poiVar);
321 }
322 if (poiSet->getSize() > 1)
323 std::cout << "Model with more than one POI are not supported - ignore extra parameters, consider only first one" << std::endl;
324
325
326
327 // This for more than one POI (not yet supported)
328 //
329 // RooLinkedListIter it = poiSet->iterator();
330 // RooRealVar* tmpPar = NULL, *tmpParA=NULL;
331 // while((tmpPar = (RooRealVar*)it.Next())){
332 // tmpParA = ((RooRealVar*)attachedSet->find(tmpPar->GetName()));
333 // tmpParA->setVal( tmpPar->getVal() );
334 // if (!tmpParA->isConstant() ) {
335 // tmpParA->setConstant();
336 // paramsSetConstant.add(*tmpParA);
337 // }
338 // }
339
340 // check if there are non-const parameters so it is worth to do the minimization
341
342 }
343
344 TStopwatch tw;
345 tw.Start();
346 double val = -1;
347
348 //check if needed to skip the fit
349 RooArgSet nllParams(*attachedSet);
351 delete attachedSet;
352 bool skipFit = (nllParams.getSize() == 0);
353
354 if (skipFit)
355 val = nll->getVal(); // just evaluate nll in conditional fits with model without nuisance params
356 else {
357
358
359 int minimPrintLevel = verbose;
360
361 RooMinimizer minim(*nll);
363 minim.setStrategy( strategy);
364 // use tolerance - but never smaller than 1 (default in RooMinimizer)
366 tol = std::max(tol,1.0); // 1.0 is the minimum value used in RooMinimizer
367 minim.setEps( tol );
368 //LM: RooMinimizer.setPrintLevel has +1 offset - so subtract here -1
369 minim.setPrintLevel(minimPrintLevel-1);
370 int status = -1;
371 minim.optimizeConst(2);
374
375 if (verbose > 0 )
376 std::cout << "AsymptoticCalculator::EvaluateNLL ........ using " << minimizer << " / " << algorithm
377 << " with strategy " << strategy << " and tolerance " << tol << std::endl;
378
379
380 for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
381 // status = minim.minimize(fMinimizer, ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo().c_str());
382 status = minim.minimize(minimizer, algorithm);
383 // RooMinimizer::minimize returns -1 when the fit fails
384 if (status >= 0) {
385 break;
386 } else {
387 if (tries == 1) {
388 printf(" ----> Doing a re-scan first\n");
389 minim.minimize(minimizer,"Scan");
390 }
391 if (tries == 2) {
393 printf(" ----> trying with strategy = 1\n");
394 minim.setStrategy(1);
395 }
396 else
397 tries++; // skip this trial if strategy is already 1
398 }
399 if (tries == 3) {
400 printf(" ----> trying with improve\n");
401 minimizer = "Minuit";
402 algorithm = "migradimproved";
403 }
404 }
405 }
406
407 RooFitResult * result = 0;
408
409 // ignore errors in Hesse or in Improve and also when matrix was made pos def (status returned = 1)
410 if (status >= 0) {
411 result = minim.save();
412 }
413 if (result){
414 if (!RooStats::IsNLLOffset() )
415 val = result->minNll();
416 else {
417 bool previous = RooAbsReal::hideOffset();
419 val = nll->getVal();
420 if (!previous) RooAbsReal::setHideOffset(kFALSE) ;
421 }
422
423 }
424 else {
425 oocoutE((TObject*)0,Fitting) << "FIT FAILED !- return a NaN NLL " << std::endl;
426 val = TMath::QuietNaN();
427 }
428
429 minim.optimizeConst(false);
430 if (result) delete result;
431
432
433 }
434
435 double muTest = 0;
436 if (verbose > 0) {
437 std::cout << "AsymptoticCalculator::EvaluateNLL - value = " << val;
438 if (poiSet) {
439 muTest = ( (RooRealVar*) poiSet->first() )->getVal();
440 std::cout << " for poi fixed at = " << muTest;
441 }
442 if (!skipFit) {
443 std::cout << "\tfit time : ";
444 tw.Print();
445 }
446 else
447 std::cout << std::endl;
448 }
449
450 // reset the parameter free which where set as constant
451 if (poiSet && paramsSetConstant.getSize() > 0) SetAllConstant(paramsSetConstant,false);
452
453
454 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
455
456 delete allParams;
457 delete nll;
458
459 return val;
460}
461
462////////////////////////////////////////////////////////////////////////////////
463/// It performs an hypothesis tests using the likelihood function
464/// and computes the p values for the null and the alternate using the asymptotic
465/// formulae for the profile likelihood ratio.
466/// See G. Cowan, K. Cranmer, E. Gross and O. Vitells.
467/// Asymptotic formulae for likelihood- based tests of new physics. Eur. Phys. J., C71:1–19, 2011.
468/// The formulae are valid only for one POI. If more than one POI exists consider as POI only the
469/// first one
470
472 int verbose = fgPrintLevel;
473
474 // re-initialized the calculator in case it is needed (pdf or data modified)
475 if (!fIsInitialized) {
476 if (!Initialize() ) {
477 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::GetHypoTest - Error initializing Asymptotic calculator - return NULL result " << endl;
478 return 0;
479 }
480 }
481
482 if (!fAsimovData) {
483 oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::GetHypoTest - Asimov data set has not been generated - return NULL result " << endl;
484 return 0;
485 }
486
487 assert(GetNullModel() );
488 assert(GetData() );
489
490 RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
491 assert(nullPdf);
492
493 // make conditional fit on null snapshot of poi
494
495 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
496 assert(nullSnapshot && nullSnapshot->getSize() > 0);
497
498 // use as POI the nullSnapshot
499 // if more than one POI exists, consider only the first one
500 RooArgSet poiTest(*nullSnapshot);
501
502 if (poiTest.getSize() > 1) {
503 oocoutW((TObject*)0,InputArguments) << "AsymptoticCalculator::GetHypoTest: snapshot has more than one POI - assume as POI first parameter " << std::endl;
504 }
505
506 RooArgSet * allParams = nullPdf->getParameters(*GetData() );
507 *allParams = fBestFitParams;
508 delete allParams;
509
510 // set the one-side condition
511 // (this works when we have only one params of interest
512 RooRealVar * muHat = dynamic_cast<RooRealVar*> ( fBestFitPoi.first() );
513 assert(muHat && "no best fit parameter defined");
514 RooRealVar * muTest = dynamic_cast<RooRealVar*> ( nullSnapshot->find(muHat->GetName() ) );
515 assert(muTest && "poi snapshot is not existing");
516
517
518
519 if (verbose> 0) {
520 std::cout << std::endl;
521 oocoutI((TObject*)0,Eval) << "AsymptoticCalculator::GetHypoTest: - perform an hypothesis test for POI ( " << muTest->GetName() << " ) = " << muTest->getVal() << std::endl;
522 oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::GetHypoTest - Find best conditional NLL on OBSERVED data set ..... " << std::endl;
523 }
524
525 // evaluate the conditional NLL on the observed data for the snapshot value
526 double condNLL = EvaluateNLL( *nullPdf, const_cast<RooAbsData&>(*GetData()), GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables(), &poiTest);
527
528 double qmu = 2.*(condNLL - fNLLObs);
529
530
531
532 if (verbose > 0)
533 oocoutP((TObject*)0,Eval) << "\t OBSERVED DATA : qmu = " << qmu << " condNLL = " << condNLL << " uncond " << fNLLObs << std::endl;
534
535
536 // this tolerance is used to avoid having negative qmu due to numerical errors
537 double tol = 2.E-3 * std::max(1.,ROOT::Math::MinimizerOptions::DefaultTolerance());
538 if (qmu < -tol || TMath::IsNaN(fNLLObs) ) {
539
540 if (qmu < 0)
541 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: Found a negative value of the qmu - retry to do the unconditional fit "
542 << std::endl;
543 else
544 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: unconditional fit failed before - retry to do it now "
545 << std::endl;
546
547
548 double nll = EvaluateNLL( *nullPdf, const_cast<RooAbsData&>(*GetData()),GetNullModel()->GetConditionalObservables(),GetNullModel()->GetGlobalObservables());
549
550 if (nll < fNLLObs || (TMath::IsNaN(fNLLObs) && !TMath::IsNaN(nll) ) ) {
551 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: Found a better unconditional minimum "
552 << " old NLL = " << fNLLObs << " old muHat " << muHat->getVal() << std::endl;
553
554 // update values
555 fNLLObs = nll;
557 assert(poi);
559 poi->snapshot(fBestFitPoi);
560 // restore also muHad since previous pointer has been deleted
561 muHat = dynamic_cast<RooRealVar*> ( fBestFitPoi.first() );
562 assert(muHat);
563
564 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: New minimum found for "
565 << " NLL = " << fNLLObs << " muHat " << muHat->getVal() << std::endl;
566
567
568 qmu = 2.*(condNLL - fNLLObs);
569
570 if (verbose > 0)
571 oocoutP((TObject*)0,Eval) << "After unconditional refit, new qmu value is " << qmu << std::endl;
572
573 }
574 }
575
576 if (qmu < -tol ) {
577 oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator: qmu is still < 0 for mu = "
578 << muTest->getVal() << " return a dummy result "
579 << std::endl;
580 return new HypoTestResult();
581 }
582 if (TMath::IsNaN(qmu) ) {
583 oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator: failure in fitting for qmu or qmuA "
584 << muTest->getVal() << " return a dummy result "
585 << std::endl;
586 return new HypoTestResult();
587 }
588
589
590
591
592
593 // compute conditional ML on Asimov data set
594 // (need to const cast because it uses fitTo which is a non const method
595 // RooArgSet asimovGlobObs;
596 // RooAbsData * asimovData = (const_cast<AsymptoticCalculator*>(this))->MakeAsimovData( poi, asimovGlobObs);
597 // set global observables to their Asimov values
598 RooArgSet globObs;
599 RooArgSet globObsSnapshot;
600 if (GetNullModel()->GetGlobalObservables() ) {
601 globObs.add(*GetNullModel()->GetGlobalObservables());
602 // store previous snapshot value
603 globObs.snapshot(globObsSnapshot);
604 globObs = fAsimovGlobObs;
605 }
606
607
608 if (verbose > 0) oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::GetHypoTest -- Find best conditional NLL on ASIMOV data set .... " << std::endl;
609
610 double condNLL_A = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables(), &poiTest);
611
612
613 double qmu_A = 2.*(condNLL_A - fNLLAsimov );
614
615 if (verbose > 0)
616 oocoutP((TObject*)0,Eval) << "\t ASIMOV data qmu_A = " << qmu_A << " condNLL = " << condNLL_A << " uncond " << fNLLAsimov << std::endl;
617
618 if (qmu_A < -tol || TMath::IsNaN(fNLLAsimov) ) {
619
620 if (qmu_A < 0)
621 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: Found a negative value of the qmu Asimov- retry to do the unconditional fit "
622 << std::endl;
623 else
624 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: Fit failed for unconditional the qmu Asimov- retry unconditional fit "
625 << std::endl;
626
627
628 double nll = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables() );
629
630 if (nll < fNLLAsimov || (TMath::IsNaN(fNLLAsimov) && !TMath::IsNaN(nll) )) {
631 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: Found a better unconditional minimum for Asimov data set"
632 << " old NLL = " << fNLLAsimov << std::endl;
633
634 // update values
635 fNLLAsimov = nll;
636
637 oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator: New minimum found for "
638 << " NLL = " << fNLLAsimov << std::endl;
639 qmu_A = 2.*(condNLL_A - fNLLAsimov);
640
641 if (verbose > 0)
642 oocoutP((TObject*)0,Eval) << "After unconditional Asimov refit, new qmu_A value is " << qmu_A << std::endl;
643
644 }
645 }
646
647 if (qmu_A < - tol) {
648 oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator: qmu_A is still < 0 for mu = "
649 << muTest->getVal() << " return a dummy result "
650 << std::endl;
651 return new HypoTestResult();
652 }
653 if (TMath::IsNaN(qmu) ) {
654 oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator: failure in fitting for qmu or qmuA "
655 << muTest->getVal() << " return a dummy result "
656 << std::endl;
657 return new HypoTestResult();
658 }
659
660
661 // restore previous value of global observables
662 globObs = globObsSnapshot;
663
664 // now we compute p-values using the asymptotic formulae
665 // described in the paper
666 // Cowan et al, Eur.Phys.J. C (2011) 71:1554
667
668 // first try to guess automatically if needed to use qtilde (or ttilde in case of two sided)
669 // if explicitly fUseQTilde this was not set
670 // qtilde is in this case used if poi is bounded at the value of the alt hypothesis
671 // for Qtilde (need to distinguish case when qmu > qmuA = mu^2/ sigma^2)
672 // (see Cowan et al, Eur.Phys.J. C(2011) 71:1554 paper equations 64 and 65
673 // (remember qmu_A = mu^2/sigma^2 )
674 bool useQTilde = false;
675 // default case (check if poi is limited or not to a zero value)
676 if (!fOneSidedDiscovery) { // qtilde is not a discovery test
677 if (fUseQTilde == -1 && !fOneSidedDiscovery) {
678 // alternate snapshot is value for which background is zero (for limits)
679 RooRealVar * muAlt = dynamic_cast<RooRealVar*>( GetAlternateModel()->GetSnapshot()->first() );
680 // null snapshot is value for which background is zero (for discovery)
681 //RooRealVar * muNull = dynamic_cast<RooRealVar*>( GetNullModel()->GetSnapshot()->first() );
682 assert(muAlt != 0 );
683 if (muTest->getMin() == muAlt->getVal() ) {
684 fUseQTilde = 1;
685 oocoutI((TObject*)0,InputArguments) << "Minimum of POI is " << muTest->getMin() << " corresponds to alt snapshot - using qtilde asymptotic formulae " << std::endl;
686 } else {
687 fUseQTilde = 0;
688 oocoutI((TObject*)0,InputArguments) << "Minimum of POI is " << muTest->getMin() << " is different to alt snapshot " << muAlt->getVal()
689 << " - using standard q asymptotic formulae " << std::endl;
690 }
691 }
692 useQTilde = fUseQTilde;
693 }
694
695
696 //check for one side condition (remember this is valid only for one poi)
697 if (fOneSided ) {
698 if ( muHat->getVal() > muTest->getVal() ) {
699 oocoutI((TObject*)0,Eval) << "Using one-sided qmu - setting qmu to zero muHat = " << muHat->getVal()
700 << " muTest = " << muTest->getVal() << std::endl;
701 qmu = 0;
702 }
703 }
704 if (fOneSidedDiscovery ) {
705 if ( muHat->getVal() < muTest->getVal() ) {
706 oocoutI((TObject*)0,Eval) << "Using one-sided discovery qmu - setting qmu to zero muHat = " << muHat->getVal()
707 << " muTest = " << muTest->getVal() << std::endl;
708 qmu = 0;
709 }
710 }
711
712 // fix for negative qmu values due to numerical errors
713 if (qmu < 0 && qmu > -tol) qmu = 0;
714 if (qmu_A < 0 && qmu_A > -tol) qmu_A = 0;
715
716 // asymptotic formula for pnull and from paper Eur.Phys.J C 2011 71:1554
717 // we have 4 different cases:
718 // t(mu), t_tilde(mu) for the 2-sided
719 // q(mu) and q_tilde(mu) for the one -sided test statistics
720
721 double pnull = -1, palt = -1;
722
723 // asymptotic formula for pnull (for only one POI)
724 // From fact that qmu is a chi2 with ndf=1
725
726 double sqrtqmu = (qmu > 0) ? std::sqrt(qmu) : 0;
727 double sqrtqmu_A = (qmu_A > 0) ? std::sqrt(qmu_A) : 0;
728
729
731 // for one-sided PL (q_mu : equations 56,57)
732 if (verbose>2) {
733 if (fOneSided)
734 oocoutI((TObject*)0,Eval) << "Using one-sided limit asymptotic formula (qmu)" << endl;
735 else
736 oocoutI((TObject*)0,Eval) << "Using one-sided discovery asymptotic formula (q0)" << endl;
737 }
738 pnull = ROOT::Math::normal_cdf_c( sqrtqmu, 1.);
739 palt = ROOT::Math::normal_cdf( sqrtqmu_A - sqrtqmu, 1.);
740 }
741 else {
742 // for 2-sided PL (t_mu : equations 35,36 in asymptotic paper)
743 if (verbose > 2) oocoutI((TObject*)0,Eval) << "Using two-sided asymptotic formula (tmu)" << endl;
744 pnull = 2.*ROOT::Math::normal_cdf_c( sqrtqmu, 1.);
745 palt = ROOT::Math::normal_cdf_c( sqrtqmu + sqrtqmu_A, 1.) +
746 ROOT::Math::normal_cdf_c( sqrtqmu - sqrtqmu_A, 1.);
747
748 }
749
750 if (useQTilde ) {
751 if (fOneSided) {
752 // for bounded one-sided (q_mu_tilde: equations 64,65)
753 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) { // to avoid case 0/0
754 if (verbose > 2) oocoutI((TObject*)0,Eval) << "Using qmu_tilde (qmu is greater than qmu_A)" << endl;
755 pnull = ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
756 palt = ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
757 }
758 }
759 else {
760 // for 2 sided bounded test statistic (N.B there is no one sided discovery qtilde)
761 // t_mu_tilde: equations 43,44 in asymptotic paper
762 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) {
763 if (verbose > 2) oocoutI((TObject*)0,Eval) << "Using tmu_tilde (qmu is greater than qmu_A)" << endl;
764 pnull = ROOT::Math::normal_cdf_c(sqrtqmu,1.) +
765 ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
766 palt = ROOT::Math::normal_cdf_c( sqrtqmu_A + sqrtqmu, 1.) +
767 ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
768 }
769 }
770 }
771
772
773
774 // create an HypoTest result but where the sampling distributions are set to zero
775 string resultname = "HypoTestAsymptotic_result";
776 HypoTestResult* res = new HypoTestResult(resultname.c_str(), pnull, palt);
777
778 if (verbose > 0)
779 //std::cout
780 oocoutP((TObject*)0,Eval)
781 << "poi = " << muTest->getVal() << " qmu = " << qmu << " qmu_A = " << qmu_A
782 << " sigma = " << muTest->getVal()/sqrtqmu_A
783 << " CLsplusb = " << pnull << " CLb = " << palt << " CLs = " << res->CLs() << std::endl;
784
785 return res;
786
787}
788
789struct PaltFunction {
790 PaltFunction( double offset, double pval, int icase) :
791 fOffset(offset), fPval(pval), fCase(icase) {}
792 double operator() (double x) const {
793 return ROOT::Math::normal_cdf_c(x + fOffset) + ROOT::Math::normal_cdf_c(fCase*(x - fOffset)) - fPval;
794 }
795 double fOffset;
796 double fPval;
797 int fCase;
798};
799
800////////////////////////////////////////////////////////////////////////////////
801/// function given the null and the alt p value - return the expected one given the N - sigma value
802
803double AsymptoticCalculator::GetExpectedPValues(double pnull, double palt, double nsigma, bool useCls, bool oneSided ) {
804 if (oneSided) {
805 double sqrtqmu = ROOT::Math::normal_quantile_c( pnull,1.);
806 double sqrtqmu_A = ROOT::Math::normal_quantile( palt,1.) + sqrtqmu;
807 double clsplusb = ROOT::Math::normal_cdf_c( sqrtqmu_A - nsigma, 1.);
808 if (!useCls) return clsplusb;
809 double clb = ROOT::Math::normal_cdf( nsigma, 1.);
810 return (clb == 0) ? -1 : clsplusb / clb;
811 }
812
813 // case of 2 sided test statistic
814 // need to compute numerically
815 double sqrttmu = ROOT::Math::normal_quantile_c( 0.5*pnull,1.);
816 if (sqrttmu == 0) {
817 // here cannot invert the function - skip the point
818 return -1;
819 }
820 // invert formula for palt to get sqrttmu_A
821 PaltFunction f( sqrttmu, palt, -1);
824 brf.SetFunction( wf, 0, 20);
825 bool ret = brf.Solve();
826 if (!ret) {
827 oocoutE((TObject*)0,Eval) << "Error finding expected p-values - return -1" << std::endl;
828 return -1;
829 }
830 double sqrttmu_A = brf.Root();
831
832 // now invert for expected value
833 PaltFunction f2( sqrttmu_A, ROOT::Math::normal_cdf( nsigma, 1.), 1);
835 brf.SetFunction(wf2,0,20);
836 ret = brf.Solve();
837 if (!ret) {
838 oocoutE((TObject*)0,Eval) << "Error finding expected p-values - return -1" << std::endl;
839 return -1;
840 }
841 return 2*ROOT::Math::normal_cdf_c( brf.Root(),1.);
842}
843
844// void GetExpectedLimit(double nsigma, double alpha, double &clsblimit, double &clslimit) {
845// // get expected limit
846// double
847// }
848
849////////////////////////////////////////////////////////////////////////////////
850/// fill bins by looping recursively on observables
851
852void AsymptoticCalculator::FillBins(const RooAbsPdf & pdf, const RooArgList &obs, RooAbsData & data, int &index, double &binVolume, int &ibin) {
853
854 bool debug = (fgPrintLevel >= 2);
855
856 RooRealVar * v = dynamic_cast<RooRealVar*>( &(obs[index]) );
857 if (!v) return;
858
859 RooArgSet obstmp(obs);
860 double expectedEvents = pdf.expectedEvents(obstmp);
861 // if (debug) {
862 // std::cout << "expected events = " << expectedEvents << std::endl;
863 // }
864
865 if (debug) cout << "looping on observable " << v->GetName() << endl;
866 for (int i = 0; i < v->getBins(); ++i) {
867 v->setBin(i);
868 if (index < obs.getSize() -1) {
869 index++; // increase index
870 double prevBinVolume = binVolume;
871 binVolume *= v->getBinWidth(i); // increase bin volume
872 FillBins(pdf, obs, data, index, binVolume, ibin);
873 index--; // decrease index
874 binVolume = prevBinVolume; // decrease also bin volume
875 }
876 else {
877
878 // this is now a new bin - compute the pdf in this bin
879 double totBinVolume = binVolume * v->getBinWidth(i);
880 double fval = pdf.getVal(&obstmp)*totBinVolume;
881
882 //if (debug) std::cout << "pdf value in the bin " << fval << " bin volume = " << totBinVolume << " " << fval*expectedEvents << std::endl;
883 if (fval*expectedEvents <= 0)
884 {
885 if (fval*expectedEvents < 0)
886 cout << "WARNING::Detected a bin with negative expected events! Please check your inputs." << endl;
887 else
888 cout << "WARNING::Detected a bin with zero expected events- skip it" << endl;
889 }
890 // have a cut off for overflows ??
891 else
892 data.add(obs, fval*expectedEvents);
893
894 if (debug) {
895 cout << "bin " << ibin << "\t";
896 for (int j=0; j < obs.getSize(); ++j) { cout << " " << ((RooRealVar&) obs[j]).getVal(); }
897 cout << " w = " << fval*expectedEvents;
898 cout << endl;
899 }
900 // RooArgSet xxx(obs);
901 // h3->Fill(((RooRealVar&) obs[0]).getVal(), ((RooRealVar&) obs[1]).getVal(), ((RooRealVar&) obs[2]).getVal() ,
902 // pdf->getVal(&xxx) );
903 ibin++;
904 }
905 }
906 //reset bin values
907 if (debug)
908 cout << "ending loop on .. " << v->GetName() << endl;
909
910 v->setBin(0);
911
912}
913
914////////////////////////////////////////////////////////////////////////////////
915/// iterate a Prod pdf to find all the Poisson or Gaussian part to set the observed value to expected one
916
918{
919 RooLinkedListIter iter(prod.pdfList().iterator());
920 bool ret = true;
921 for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
922 if (!a->dependsOn(obs)) continue;
923 RooPoisson *pois = 0;
924 RooGaussian * gaus = 0;
925 if ((pois = dynamic_cast<RooPoisson *>(a)) != 0) {
926 ret &= SetObsToExpected(*pois, obs);
927 pois->setNoRounding(true); //needed since expected value is not an integer
928 } else if ((gaus = dynamic_cast<RooGaussian *>(a)) != 0) {
929 ret &= SetObsToExpected(*gaus, obs);
930 } else {
931 // should try to add also lognormal case ?
932 RooProdPdf *subprod = dynamic_cast<RooProdPdf *>(a);
933 if (subprod)
934 ret &= SetObsToExpected(*subprod, obs);
935 else {
936 oocoutE((TObject*)0,InputArguments) << "Illegal term in counting model: "
937 << "the PDF " << a->GetName()
938 << " depends on the observables, but is not a Poisson, Gaussian or Product"
939 << endl;
940 return false;
941 }
942 }
943 }
944
945 return ret;
946}
947
948////////////////////////////////////////////////////////////////////////////////
949/// set observed value to the expected one
950/// works for Gaussian, Poisson or LogNormal
951/// assumes mean parameter value is the argument not constant and not depending on observables
952/// (if more than two arguments are not constant will use first one but print a warning !)
953/// need to iterate on the components of the Poisson to get n and nu (nu can be a RooAbsReal)
954/// (code from G. Petrucciani and extended by L.M.)
955
957{
958 RooRealVar *myobs = 0;
959 RooAbsReal *myexp = 0;
960 const char * pdfName = pdf.IsA()->GetName();
961 RooFIter iter(pdf.serverMIterator());
962 for (RooAbsArg *a = iter.next(); a != 0; a = iter.next()) {
963 if (obs.contains(*a)) {
964 if (myobs != 0) {
965 oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Has two observables ?? " << endl;
966 return false;
967 }
968 myobs = dynamic_cast<RooRealVar *>(a);
969 if (myobs == 0) {
970 oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Observable is not a RooRealVar??" << endl;
971 return false;
972 }
973 } else {
974 if (!a->isConstant() ) {
975 if (myexp != 0) {
976 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Has two non-const arguments " << endl;
977 return false;
978 }
979 myexp = dynamic_cast<RooAbsReal *>(a);
980 if (myexp == 0) {
981 oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Expected is not a RooAbsReal??" << endl;
982 return false;
983 }
984 }
985 }
986 }
987 if (myobs == 0) {
988 oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : No observable?" << endl;
989 return false;
990 }
991 if (myexp == 0) {
992 oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : No observable?" << endl;
993 return false;
994 }
995
996 myobs->setVal(myexp->getVal());
997
998 if (fgPrintLevel > 2) {
999 std::cout << "SetObsToExpected : setting " << myobs->GetName() << " to expected value " << myexp->getVal() << " of " << myexp->GetName() << std::endl;
1000 }
1001
1002 return true;
1003}
1004
1005////////////////////////////////////////////////////////////////////////////////
1006/// generate counting Asimov data for the case when the pdf cannot be extended
1007/// assume pdf is a RooPoisson or can be decomposed in a product of RooPoisson,
1008/// otherwise we cannot know how to make the Asimov data sets in the other cases
1009
1011 RooArgSet obs(observables);
1012 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
1013 RooPoisson *pois = 0;
1014 RooGaussian *gaus = 0;
1015
1016 if (fgPrintLevel > 1)
1017 std::cout << "generate counting Asimov data for pdf of type " << pdf.IsA()->GetName() << std::endl;
1018
1019 bool r = false;
1020 if (prod != 0) {
1021 r = SetObsToExpected(*prod, observables);
1022 } else if ((pois = dynamic_cast<RooPoisson *>(&pdf)) != 0) {
1023 r = SetObsToExpected(*pois, observables);
1024 // we need in this case to set Poisson to real values
1025 pois->setNoRounding(true);
1026 } else if ((gaus = dynamic_cast<RooGaussian *>(&pdf)) != 0) {
1027 r = SetObsToExpected(*gaus, observables);
1028 } else {
1029 oocoutE((TObject*)0,InputArguments) << "A counting model pdf must be either a RooProdPdf or a RooPoisson or a RooGaussian" << endl;
1030 }
1031 if (!r) return 0;
1032 int icat = 0;
1033 if (channelCat) {
1034 icat = channelCat->getIndex();
1035 }
1036
1037 RooDataSet *ret = new RooDataSet(TString::Format("CountingAsimovData%d",icat),TString::Format("CountingAsimovData%d",icat), obs);
1038 ret->add(obs);
1039 return ret;
1040}
1041
1042////////////////////////////////////////////////////////////////////////////////
1043/// compute the asimov data set for an observable of a pdf
1044/// use the number of bins sets in the observables
1045/// to do : (possibility to change number of bins)
1046/// implement integration over bin content
1047
1048RooAbsData * AsymptoticCalculator::GenerateAsimovDataSinglePdf(const RooAbsPdf & pdf, const RooArgSet & allobs, const RooRealVar & weightVar, RooCategory * channelCat) {
1049
1050 int printLevel = fgPrintLevel;
1051
1052 // Get observables defined by the pdf associated with this state
1053 std::unique_ptr<RooArgSet> obs(pdf.getObservables(allobs) );
1054
1055
1056 // if pdf cannot be extended assume is then a counting experiment
1057 if (!pdf.canBeExtended() ) return GenerateCountingAsimovData(const_cast<RooAbsPdf&>(pdf), *obs, weightVar, channelCat);
1058
1059 RooArgSet obsAndWeight(*obs);
1060 obsAndWeight.add(weightVar);
1061
1062 RooDataSet* asimovData = 0;
1063 if (channelCat) {
1064 int icat = channelCat->getIndex();
1065 asimovData = new RooDataSet(TString::Format("AsimovData%d",icat),TString::Format("combAsimovData%d",icat),
1066 RooArgSet(obsAndWeight,*channelCat),RooFit::WeightVar(weightVar));
1067 }
1068 else
1069 asimovData = new RooDataSet("AsimovData","AsimovData",RooArgSet(obsAndWeight),RooFit::WeightVar(weightVar));
1070
1071 // This works only for 1D observables
1072 //RooRealVar* thisObs = ((RooRealVar*)obstmp->first());
1073
1074 RooArgList obsList(*obs);
1075
1076 // loop on observables and on the bins
1077 if (printLevel >= 2) {
1078 cout << "Generating Asimov data for pdf " << pdf.GetName() << endl;
1079 cout << "list of observables " << endl;
1080 obsList.Print();
1081 }
1082
1083 int obsIndex = 0;
1084 double binVolume = 1;
1085 int nbins = 0;
1086 FillBins(pdf, obsList, *asimovData, obsIndex, binVolume, nbins);
1087 if (printLevel >= 2)
1088 cout << "filled from " << pdf.GetName() << " " << nbins << " nbins " << " volume is " << binVolume << endl;
1089
1090 // for (int iobs = 0; iobs < obsList.getSize(); ++iobs) {
1091 // RooRealVar * thisObs = dynamic_cast<RooRealVar*> &obsList[i];
1092 // if (thisObs == 0) continue;
1093 // // loop on the bin contents
1094 // for(int ibin=0; ibin<thisObs->numBins(); ++ibin){
1095 // thisObs->setBin(ibin);
1096
1097 // thisNorm=pdftmp->getVal(obstmp)*thisObs->getBinWidth(jj);
1098 // if (thisNorm*expectedEvents <= 0)
1099 // {
1100 // cout << "WARNING::Detected bin with zero expected events! Please check your inputs." << endl;
1101 // }
1102 // // have a cut off for overflows ??
1103 // obsDataUnbinned->add(*mc->GetObservables(), thisNorm*expectedEvents);
1104 // }
1105
1106 if (printLevel >= 1)
1107 {
1108 asimovData->Print();
1109 //cout <<"sum entries "<< asimovData->sumEntries()<<endl;
1110 }
1111 if( TMath::IsNaN(asimovData->sumEntries()) ){
1112 cout << "sum entries is nan"<<endl;
1113 assert(0);
1114 delete asimovData;
1115 asimovData = 0;
1116 }
1117
1118 return asimovData;
1119
1120}
1121
1122////////////////////////////////////////////////////////////////////////////////
1123/// generate the asimov data for the observables (not the global ones)
1124/// need to deal with the case of a sim pdf
1125
1127
1128 int printLevel = fgPrintLevel;
1129
1130 unique_ptr<RooRealVar> weightVar (new RooRealVar("binWeightAsimov", "binWeightAsimov", 1, 0, 1.E30 ));
1131
1132 if (printLevel > 1) cout <<" Generate Asimov data for observables"<<endl;
1133 //RooDataSet* simData=NULL;
1134 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf);
1135 if (!simPdf) {
1136 // generate data for non sim pdf
1137 return GenerateAsimovDataSinglePdf( pdf, observables, *weightVar, 0);
1138 }
1139
1140 std::map<std::string, RooDataSet*> asimovDataMap;
1141
1142 //look at category of simpdf
1143 RooCategory& channelCat = (RooCategory&)simPdf->indexCat();
1144 int nrIndices = channelCat.numTypes();
1145 if( nrIndices == 0 ) {
1146 oocoutW((TObject*)0,Generation) << "Simultaneous pdf does not contain any categories." << endl;
1147 }
1148 for (int i=0;i<nrIndices;i++){
1149 channelCat.setIndex(i);
1150 //iFrame++;
1151 // Get pdf associated with state from simpdf
1152 RooAbsPdf* pdftmp = simPdf->getPdf(channelCat.getLabel()) ;
1153 assert(pdftmp != 0);
1154
1155 if (printLevel > 1)
1156 {
1157 cout << "on type " << channelCat.getLabel() << " " << channelCat.getIndex() << endl;
1158 }
1159
1160 RooAbsData * dataSinglePdf = GenerateAsimovDataSinglePdf( *pdftmp, observables, *weightVar, &channelCat);
1161 //((RooRealVar*)obstmp->first())->Print();
1162 //cout << "expected events " << pdftmp->expectedEvents(*obstmp) << endl;
1163 if (!dataSinglePdf) {
1164 oocoutE((TObject*)0,Generation) << "Error generating an Asimov data set for pdf " << pdftmp->GetName() << endl;
1165 return 0;
1166 }
1167
1168 if (asimovDataMap.count(string(channelCat.getLabel())) != 0) {
1169 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::GenerateAsimovData(): The PDF for " << channelCat.getLabel()
1170 << " was already defined. It will be overridden. The faulty category definitions follow:" << endl;
1171 channelCat.Print("V");
1172 }
1173
1174 asimovDataMap[string(channelCat.getLabel())] = (RooDataSet*) dataSinglePdf;
1175
1176 if (printLevel > 1)
1177 {
1178 cout << "channel: " << channelCat.getLabel() << ", data: ";
1179 dataSinglePdf->Print();
1180 cout << endl;
1181 }
1182 }
1183
1184 RooArgSet obsAndWeight(observables);
1185 obsAndWeight.add(*weightVar);
1186
1187
1188 RooDataSet* asimovData = new RooDataSet("asimovDataFullModel","asimovDataFullModel",RooArgSet(obsAndWeight,channelCat),
1189 RooFit::Index(channelCat),RooFit::Import(asimovDataMap),RooFit::WeightVar(*weightVar));
1190
1191 for (auto element : asimovDataMap) {
1192 delete element.second;
1193 }
1194
1195 return asimovData;
1196
1197}
1198
1199////////////////////////////////////////////////////////////////////////////////
1200/// static function to the an Asimov data set
1201/// given an observed dat set, a model and a snapshot of poi.
1202/// Return the asimov data set + global observables set to values satisfying the constraints
1203
1204RooAbsData * AsymptoticCalculator::MakeAsimovData(RooAbsData & realData, const ModelConfig & model, const RooArgSet & paramValues, RooArgSet & asimovGlobObs, const RooArgSet * genPoiValues ) {
1205
1206 int verbose = fgPrintLevel;
1207
1208
1209 RooArgSet poi(*model.GetParametersOfInterest());
1210 poi = paramValues;
1211 //RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
1212 // set poi constant for conditional MLE
1213 // need to fit nuisance parameters at their conditional MLE value
1214 RooLinkedListIter it = poi.iterator();
1215 RooRealVar* tmpPar = NULL;
1216 RooArgSet paramsSetConstant;
1217 while((tmpPar = (RooRealVar*)it.Next())){
1218 tmpPar->setConstant();
1219 if (verbose>0)
1220 std::cout << "MakeAsimov: Setting poi " << tmpPar->GetName() << " to a constant value = " << tmpPar->getVal() << std::endl;
1221 paramsSetConstant.add(*tmpPar);
1222 }
1223
1224 // find conditional value of the nuisance parameters
1225 bool hasFloatParams = false;
1226 RooArgSet constrainParams;
1227 if (model.GetNuisanceParameters()) {
1228 constrainParams.add(*model.GetNuisanceParameters());
1229 RooStats::RemoveConstantParameters(&constrainParams);
1230 if (constrainParams.getSize() > 0) hasFloatParams = true;
1231
1232 } else {
1233 // Do we have free parameters anyway that need fitting?
1234 std::unique_ptr<RooArgSet> params(model.GetPdf()->getParameters(realData));
1235 RooLinkedListIter iter(params->iterator());
1236 for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
1237 RooRealVar *rrv = dynamic_cast<RooRealVar *>(a);
1238 if ( rrv != 0 && rrv->isConstant() == false ) { hasFloatParams = true; break; }
1239 }
1240 }
1241 if (hasFloatParams) {
1242 // models need to be fitted to find best nuisance parameter values
1243
1244 TStopwatch tw2; tw2.Start();
1246 if (verbose>0) {
1247 std::cout << "MakeAsimov: doing a conditional fit for finding best nuisance values " << std::endl;
1248 minimPrintLevel = verbose;
1249 if (verbose>1) {
1250 std::cout << "POI values:\n"; poi.Print("v");
1251 if (verbose > 2) {
1252 std::cout << "Nuis param values:\n";
1253 constrainParams.Print("v");
1254 }
1255 }
1256 }
1259
1260 RooArgSet conditionalObs;
1261 if (model.GetConditionalObservables()) conditionalObs.add(*model.GetConditionalObservables());
1262 RooArgSet globalObs;
1263 if (model.GetGlobalObservables()) globalObs.add(*model.GetGlobalObservables());
1264
1265 std::string minimizerType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
1266 std::string minimizerAlgo = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo();
1267 model.GetPdf()->fitTo(realData, RooFit::Minimizer(minimizerType.c_str(),minimizerAlgo.c_str()),
1269 RooFit::PrintLevel(minimPrintLevel-1), RooFit::Hesse(false),
1270 RooFit::Constrain(constrainParams),RooFit::GlobalObservables(globalObs),
1272 if (verbose>0) { std::cout << "fit time "; tw2.Print();}
1273 if (verbose > 1) {
1274 // after the fit the nuisance parameters will have their best fit value
1275 if (model.GetNuisanceParameters() ) {
1276 std::cout << "Nuisance parameters after fit for asimov dataset: " << std::endl;
1277 model.GetNuisanceParameters()->Print("V");
1278 }
1279 }
1280
1281 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
1282
1283 }
1284
1285 // restore the parameters which were set constant
1286 SetAllConstant(paramsSetConstant, false);
1287
1288
1289 RooArgSet * allParams = model.GetPdf()->getParameters(realData);
1291
1292 // if a RooArgSet of poi is passed , different poi will be used for generating the Asimov data set
1293 if (genPoiValues) *allParams = *genPoiValues;
1294
1295 // now do the actual generation of the AsimovData Set
1296 // no need to pass parameters values since we have set them before
1297 RooAbsData * asimovData = MakeAsimovData(model, *allParams, asimovGlobObs);
1298
1299 delete allParams;
1300
1301 return asimovData;
1302
1303}
1304
1305////////////////////////////////////////////////////////////////////////////////
1306/// static function to the an Asimov data set
1307/// given the model and the values of all parameters including the nuisance
1308/// Return the asimov data set + global observables set to values satisfying the constraints
1309
1310RooAbsData * AsymptoticCalculator::MakeAsimovData(const ModelConfig & model, const RooArgSet & allParamValues, RooArgSet & asimovGlobObs) {
1311
1312 int verbose = fgPrintLevel;
1313
1314 TStopwatch tw;
1315 tw.Start();
1316
1317 // set the parameter values (do I need the poi to be constant ? )
1318 // the nuisance parameter values could be set at their fitted value (the MLE)
1319 if (allParamValues.getSize() > 0) {
1320 RooArgSet * allVars = model.GetPdf()->getVariables();
1321 *allVars = allParamValues;
1322 delete allVars;
1323 }
1324
1325
1326 // generate the Asimov data set for the observables
1327 RooAbsData * asimov = GenerateAsimovData(*model.GetPdf() , *model.GetObservables() );
1328
1329 if (verbose>0) {
1330 std::cout << "Generated Asimov data for observables "; (model.GetObservables() )->Print();
1331 if (verbose > 1) {
1332 if (asimov->numEntries() == 1 ) {
1333 std::cout << "--- Asimov data values \n";
1334 asimov->get()->Print("v");
1335 }
1336 else {
1337 std::cout << "--- Asimov data numEntries = " << asimov->numEntries() << " sumOfEntries = " << asimov->sumEntries() << std::endl;
1338 }
1339 std::cout << "\ttime for generating : "; tw.Print();
1340 }
1341 }
1342
1343
1344 // Now need to have in ASIMOV the data sets also the global observables
1345 // Their values must be the one satisfying the constraint.
1346 // to do it make a nuisance pdf with all product of constraints and then
1347 // assign to each constraint a glob observable value = to the current fitted nuisance parameter value
1348 // IN general one should solve in general the system of equations f( gobs| nuispar ) = 0 where f are the
1349 // derivatives of the constraint with respect the nuisance parameter and they are evaluated at the best fit nuisance
1350 // parameter points
1351 // As simple solution assume that constrain has a direct dependence on the nuisance parameter, i.e.
1352 // Constraint (gobs, func( nuispar) ) and the condition is satisfied for
1353 // gobs = func( nuispar) where nunispar is at the MLE value
1354
1355
1356 if (model.GetGlobalObservables() && model.GetGlobalObservables()->getSize() > 0) {
1357
1358 if (verbose>1) {
1359 std::cout << "Generating Asimov data for global observables " << std::endl;
1360 }
1361
1362 RooArgSet gobs(*model.GetGlobalObservables());
1363
1364 // snapshot data global observables
1365 RooArgSet snapGlobalObsData;
1366 SetAllConstant(gobs, true);
1367 gobs.snapshot(snapGlobalObsData);
1368
1369
1370 RooArgSet nuis;
1371 if (model.GetNuisanceParameters()) nuis.add(*model.GetNuisanceParameters());
1372 if (nuis.getSize() == 0) {
1373 oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: model does not have nuisance parameters but has global observables"
1374 << " set global observables to model values " << endl;
1375 asimovGlobObs = gobs;
1376 return asimov;
1377 }
1378
1379 // part 1: create the nuisance pdf
1380 std::unique_ptr<RooAbsPdf> nuispdf(RooStats::MakeNuisancePdf(model,"TempNuisPdf") );
1381 if (nuispdf.get() == 0) {
1382 oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: model has nuisance parameters and global obs but no nuisance pdf "
1383 << std::endl;
1384 }
1385 // unfold the nuisance pdf if it is a prod pdf
1386 RooArgList pdfList;
1387 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(nuispdf.get());
1388 if (prod ) {
1389 pdfList.add(prod->pdfList());
1390 }
1391 else
1392 // nothing to unfold - just use the pdf
1393 pdfList.add(*nuispdf.get());
1394
1395 RooLinkedListIter iter(pdfList.iterator());
1396 for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
1397 RooAbsPdf *cterm = dynamic_cast<RooAbsPdf *>(a);
1398 assert(cterm && "AsimovUtils: a factor of the nuisance pdf is not a Pdf!");
1399 if (!cterm->dependsOn(nuis)) continue; // dummy constraints
1400 // skip also the case of uniform components
1401 if (typeid(*cterm) == typeid(RooUniform)) continue;
1402
1403 std::unique_ptr<RooArgSet> cpars(cterm->getParameters(&gobs));
1404 std::unique_ptr<RooArgSet> cgobs(cterm->getObservables(&gobs));
1405 if (cgobs->getSize() > 1) {
1406 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1407 << " has multiple global observables -cannot generate - skip it" << std::endl;
1408 continue;
1409 }
1410 else if (cgobs->getSize() == 0) {
1411 oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1412 << " has no global observables - skip it" << std::endl;
1413 continue;
1414 }
1415 // the variable representing the global observable
1416 RooRealVar &rrv = dynamic_cast<RooRealVar &>(*cgobs->first());
1417
1418 // remove the constant parameters in cpars
1420 if (cpars->getSize() != 1) {
1421 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1422 << cterm->GetName() << " has multiple floating params - cannot generate - skip it " << std::endl;
1423 continue;
1424 }
1425
1426 bool foundServer = false;
1427 // note : this will work only for this type of constraints
1428 // expressed as RooPoisson, RooGaussian, RooLognormal, RooGamma
1429 TClass * cClass = cterm->IsA();
1430 if (verbose > 2) std::cout << "Constraint " << cterm->GetName() << " of type " << cClass->GetName() << std::endl;
1431 if ( cClass != RooGaussian::Class() && cClass != RooPoisson::Class() &&
1432 cClass != RooGamma::Class() && cClass != RooLognormal::Class() &&
1433 cClass != RooBifurGauss::Class() ) {
1434 TString className = (cClass) ? cClass->GetName() : "undefined";
1435 oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1436 << cterm->GetName() << " of type " << className
1437 << " is a non-supported type - result might be not correct " << std::endl;
1438 }
1439
1440 // in case of a Poisson constraint make sure the rounding is not set
1441 if (cClass == RooPoisson::Class() ) {
1442 RooPoisson * pois = dynamic_cast<RooPoisson*>(cterm);
1443 assert(pois);
1444 pois->setNoRounding(true);
1445 }
1446
1447 // look at server of the constraint term and check if the global observable is part of the server
1448 RooAbsArg * arg = cterm->findServer(rrv);
1449 if (!arg) {
1450 // special case is for the Gamma where one might define the global observable n and you have a Gamma(b, n+1, ...._
1451 // in this case n+1 is the server and we don;t have a direct dependency, but we want to set n to the b value
1452 // so in case of the Gamma ignore this test
1453 if ( cClass != RooGamma::Class() ) {
1454 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1455 << cterm->GetName() << " has no direct dependence on global observable- cannot generate it " << std::endl;
1456 continue;
1457 }
1458 }
1459
1460 // loop on the server of the constraint term
1461 // need to treat the Gamma as a special case
1462 // the mode of the Gamma is (k-1)*theta where theta is the inverse of the rate parameter.
1463 // we assume that the global observable is defined as ngobs = k-1 and the theta parameter has the name theta otherwise we use other procedure which might be wrong
1464 RooAbsReal * thetaGamma = 0;
1465 if ( cClass == RooGamma::Class() ) {
1466 RooFIter itc(cterm->serverMIterator() );
1467 for (RooAbsArg *a2 = itc.next(); a2 != 0; a2 = itc.next()) {
1468 if (TString(a2->GetName()).Contains("theta") ) {
1469 thetaGamma = dynamic_cast<RooAbsReal*>(a2);
1470 break;
1471 }
1472 }
1473 if (thetaGamma == 0) {
1474 oocoutI((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1475 << cterm->GetName() << " is a Gamma distribution and no server named theta is found. Assume that the Gamma scale is 1 " << std::endl;
1476 }
1477 else {
1478 if (verbose>2)
1479 std::cout << "Gamma constraint has a scale " << thetaGamma->GetName() << " = " << thetaGamma->getVal() << std::endl;
1480 }
1481 }
1482 RooFIter iter2(cterm->serverMIterator() );
1483 for (RooAbsArg *a2 = iter2.next(); a2 != 0; a2 = iter2.next()) {
1484 RooAbsReal * rrv2 = dynamic_cast<RooAbsReal *>(a2);
1485 if (verbose > 2) std::cout << "Loop on constraint server term " << a2->GetName() << std::endl;
1486 if (rrv2 && rrv2->dependsOn(nuis) ) {
1487
1488
1489 // found server depending on nuisance
1490 if (foundServer) {
1491 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1492 << cterm->GetName() << " constraint term has more server depending on nuisance- cannot generate it " <<
1493 std::endl;
1494 foundServer = false;
1495 break;
1496 }
1497 if (thetaGamma && thetaGamma->getVal() > 0)
1498 rrv.setVal( rrv2->getVal() / thetaGamma->getVal() );
1499 else
1500 rrv.setVal( rrv2->getVal() );
1501 foundServer = true;
1502
1503 if (verbose>2)
1504 std::cout << "setting global observable "<< rrv.GetName() << " to value " << rrv.getVal()
1505 << " which comes from " << rrv2->GetName() << std::endl;
1506 }
1507 }
1508
1509 if (!foundServer) {
1510 oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData - can't find nuisance for constraint term - global observables will not be set to Asimov value " << cterm->GetName() << std::endl;
1511 std::cerr << "Parameters: " << std::endl;
1512 cpars->Print("V");
1513 std::cerr << "Observables: " << std::endl;
1514 cgobs->Print("V");
1515 }
1516// rrv.setVal(match->getVal());
1517 }
1518
1519 // make a snapshot of global observables
1520 // needed this ?? (LM)
1521
1522 asimovGlobObs.removeAll();
1523 SetAllConstant(gobs, true);
1524 gobs.snapshot(asimovGlobObs);
1525
1526 // revert global observables to the data value
1527 gobs = snapGlobalObsData;
1528
1529 if (verbose>0) {
1530 std::cout << "Generated Asimov data for global observables ";
1531 if (verbose == 1) gobs.Print();
1532 }
1533
1534 if (verbose > 1) {
1535 std::cout << "\nGlobal observables for data: " << std::endl;
1536 gobs.Print("V");
1537 std::cout << "\nGlobal observables for asimov: " << std::endl;
1538 asimovGlobObs.Print("V");
1539 }
1540
1541
1542 }
1543
1544 return asimov;
1545
1546}
void Class()
Definition: Class.C:29
SVector< double, 2 > v
Definition: Dict.h:5
ROOT::R::TRInterface & r
Definition: Object.C:4
#define f(i)
Definition: RSha256.hxx:104
#define oocoutW(o, a)
Definition: RooMsgService.h:46
#define oocoutE(o, a)
Definition: RooMsgService.h:47
#define oocoutF(o, a)
Definition: RooMsgService.h:48
#define oocoutI(o, a)
Definition: RooMsgService.h:44
#define oocoutP(o, a)
Definition: RooMsgService.h:45
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
double sqrt(double)
TRObject operator()(const T1 &t1) const
Class for finding the root of a one dimensional function using the Brent algorithm.
double Root() const
Returns root value.
bool Solve(int maxIter=100, double absTol=1E-8, double relTol=1E-10)
Returns the X value corresponding to the function value fy for (xmin<x<xmax).
bool SetFunction(const ROOT::Math::IGenFunction &f, double xlow, double xup)
Sets the function for the rest of the algorithms.
static const std::string & DefaultMinimizerType()
static const std::string & DefaultMinimizerAlgo()
Template class to wrap any C++ callable object which takes one argument i.e.
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
Definition: RooAbsArg.cxx:736
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:533
RooFIter serverMIterator() const
Definition: RooAbsArg.h:119
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2074
Bool_t isConstant() const
Definition: RooAbsArg.h:266
RooAbsArg * findServer(const char *name) const
Definition: RooAbsArg.h:122
Int_t numTypes(const char *=0) const
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
RooAbsArg * first() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
RooLinkedListIter iterator(Bool_t dir=kIterForward) const
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:37
virtual const RooArgSet * get() const
Definition: RooAbsData.h:79
virtual Double_t sumEntries() const =0
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsData.h:161
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:285
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Definition: RooAbsPdf.cxx:781
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:230
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
Definition: RooAbsPdf.cxx:2855
virtual Double_t getMax(const char *name=0) const
virtual Int_t getBins(const char *name=0) const
void setConstant(Bool_t value=kTRUE)
virtual Double_t getMin(const char *name=0) const
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
static void setHideOffset(Bool_t flag)
Definition: RooAbsReal.cxx:117
static Bool_t hideOffset()
Definition: RooAbsReal.cxx:118
Double_t getVal(const RooArgSet *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
virtual const char * getLabel() const
Return label string of current state.
Definition: RooCategory.h:39
virtual Int_t getIndex() const
Return index number of current state.
Definition: RooCategory.h:34
virtual Bool_t setIndex(Int_t index, Bool_t printError=kTRUE)
Set value by specifying the index code of the desired state.
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
virtual Double_t sumEntries() const
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0)
Add a data point, with its coordinates specified in the 'data' argset, to the data set.
RooAbsArg * next()
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
Double_t minNll() const
Definition: RooFitResult.h:98
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooLinkedListIter is the TIterator implementation for RooLinkedList.
virtual TObject * Next()
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
Definition: RooMinimizer.h:38
RooFitResult * save(const char *name=0, const char *title=0)
Save and return a RooFitResult snaphot of current minimizer status.
Int_t minimize(const char *type, const char *alg=0)
void setEps(Double_t eps)
Change MINUIT epsilon.
Int_t setPrintLevel(Int_t newLevel)
Change the MINUIT internal printing level.
void optimizeConst(Int_t flag)
If flag is true, perform constant term optimization on function being minimized.
void setStrategy(Int_t strat)
Change MINUIT strategy to istrat.
static RooMsgService & instance()
Return reference to singleton instance.
void setGlobalKillBelow(RooFit::MsgLevel level)
RooFit::MsgLevel globalKillBelow() const
Poisson pdf.
Definition: RooPoisson.h:19
void setNoRounding(bool flag=kTRUE)
Definition: RooPoisson.h:33
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
const RooArgList & pdfList() const
Definition: RooProdPdf.h:69
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
void setBins(Int_t nBins, const char *name=0)
Definition: RooRealVar.h:77
Double_t getError() const
Definition: RooRealVar.h:53
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:204
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const RooAbsCategoryLValue & indexCat() const
RooAbsPdf * getPdf(const char *catName) const
Return the p.d.f associated with the given index category name.
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio.
static RooAbsData * GenerateAsimovDataSinglePdf(const RooAbsPdf &pdf, const RooArgSet &obs, const RooRealVar &weightVar, RooCategory *channelCat=0)
compute the asimov data set for an observable of a pdf use the number of bins sets in the observables...
virtual HypoTestResult * GetHypoTest() const
re-implement HypoTest computation using the asymptotic
static double GetExpectedPValues(double pnull, double palt, double nsigma, bool usecls, bool oneSided=true)
function given the null and the alt p value - return the expected one given the N - sigma value
static RooAbsData * MakeAsimovData(RooAbsData &data, const ModelConfig &model, const RooArgSet &poiValues, RooArgSet &globObs, const RooArgSet *genPoiValues=0)
make the asimov data from the ModelConfig and list of poi - return data set and snapshot of global ob...
static bool SetObsToExpected(RooAbsPdf &pdf, const RooArgSet &obs)
set observed value to the expected one works for Gaussian, Poisson or LogNormal assumes mean paramete...
static RooAbsData * GenerateAsimovData(const RooAbsPdf &pdf, const RooArgSet &observables)
generate the asimov data for the observables (not the global ones) need to deal with the case of a si...
static RooAbsData * GenerateCountingAsimovData(RooAbsPdf &pdf, const RooArgSet &obs, const RooRealVar &weightVar, RooCategory *channelCat=0)
generate counting Asimov data for the case when the pdf cannot be extended assume pdf is a RooPoisson...
int fUseQTilde
flag to check if calculator is initialized
AsymptoticCalculator(RooAbsData &data, const ModelConfig &altModel, const ModelConfig &nullModel, bool nominalAsimov=false)
constructor for asymptotic calculator from Data set and ModelConfig
static void FillBins(const RooAbsPdf &pdf, const RooArgList &obs, RooAbsData &data, int &index, double &binVolume, int &ibin)
fill bins by looping recursively on observables
static double EvaluateNLL(RooAbsPdf &pdf, RooAbsData &data, const RooArgSet *condObs, const RooArgSet *globObs, const RooArgSet *poiSet=0)
bool Initialize() const
initialize the calculator by performing a global fit and make the Asimov data set
Common base class for the Hypothesis Test Calculators.
const ModelConfig * GetNullModel(void) const
const ModelConfig * GetAlternateModel(void) const
const RooAbsData * GetData(void) const
HypoTestResult is a base class for results from hypothesis tests.
virtual Double_t CLs() const
is simply (not a method, but a quantity)
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:231
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
Definition: ModelConfig.h:243
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return NULL if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:228
Flat p.d.f.
Definition: RooUniform.h:24
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:75
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
Stopwatch class.
Definition: TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition: TString.cxx:2286
double normal_cdf_c(double x, double sigma=1, double x0=0)
Complement of the cumulative distribution function of the normal (Gaussian) distribution (upper tail)...
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
double normal_quantile(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the lower tail of the normal (Gaussian) distri...
double normal_quantile_c(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the upper tail of the normal (Gaussian) distri...
Double_t x[n]
Definition: legend1.C:17
void Print(std::ostream &os, const OptionType &opt)
RooCmdArg WeightVar(const char *name, Bool_t reinterpretAsWeight=kFALSE)
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg Strategy(Int_t code)
RooCmdArg Index(RooCategory &icat)
RooCmdArg GlobalObservables(const RooArgSet &globs)
RooCmdArg Hesse(Bool_t flag=kTRUE)
@ Minimization
Definition: RooGlobalFunc.h:57
@ Generation
Definition: RooGlobalFunc.h:57
@ InputArguments
Definition: RooGlobalFunc.h:58
RooCmdArg CloneData(Bool_t flag)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg ConditionalObservables(const RooArgSet &set)
RooCmdArg Offset(Bool_t flag=kTRUE)
RooCmdArg Minimizer(const char *type, const char *alg=0)
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
Definition: RooStatsUtils.h:82
void RemoveConstantParameters(RooArgSet *set)
Definition: RooStatsUtils.h:62
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
bool IsNLLOffset()
Bool_t IsNaN(Double_t x)
Definition: TMath.h:880
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754
Definition: TMath.h:889
Definition: first.py:1
STL namespace.
auto * a
Definition: textangle.C:12