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