Logo ROOT   6.12/07
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(),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();
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(), 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 
289 Double_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 
730  if (fOneSided || fOneSidedDiscovery) {
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 
789 struct 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 
803 double 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 
852 void 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 = false;
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  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  SetObsToExpected(*gaus, obs);
930  } else {
931  // should try to add also lognormal case ?
932  RooProdPdf *subprod = dynamic_cast<RooProdPdf *>(a);
933  if (subprod)
934  return SetObsToExpected(*subprod, obs);
935  else {
936  oocoutE((TObject*)0,InputArguments) << "Illegal term in counting model: depends on observables, but not Poisson or Gaussian or Product"
937  << endl;
938  return false;
939  }
940  }
941  ret = (pois != 0 || gaus != 0 );
942  }
943  return ret;
944 }
945 
946 ////////////////////////////////////////////////////////////////////////////////
947 /// set observed value to the expected one
948 /// works for Gaussian, Poisson or LogNormal
949 /// assumes mean parameter value is the argument not constant and not depending on observables
950 /// (if more than two arguments are not constant will use first one but print a warning !)
951 /// need to iterate on the components of the Poisson to get n and nu (nu can be a RooAbsReal)
952 /// (code from G. Petrucciani and extended by L.M.)
953 
955 {
956  RooRealVar *myobs = 0;
957  RooAbsReal *myexp = 0;
958  const char * pdfName = pdf.IsA()->GetName();
959  RooFIter iter(pdf.serverMIterator());
960  for (RooAbsArg *a = iter.next(); a != 0; a = iter.next()) {
961  if (obs.contains(*a)) {
962  if (myobs != 0) {
963  oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Has two observables ?? " << endl;
964  return false;
965  }
966  myobs = dynamic_cast<RooRealVar *>(a);
967  if (myobs == 0) {
968  oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Observable is not a RooRealVar??" << endl;
969  return false;
970  }
971  } else {
972  if (!a->isConstant() ) {
973  if (myexp != 0) {
974  oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Has two non-const arguments " << endl;
975  return false;
976  }
977  myexp = dynamic_cast<RooAbsReal *>(a);
978  if (myexp == 0) {
979  oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Expected is not a RooAbsReal??" << endl;
980  return false;
981  }
982  }
983  }
984  }
985  if (myobs == 0) {
986  oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : No observable?" << endl;
987  return false;
988  }
989  if (myexp == 0) {
990  oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : No observable?" << endl;
991  return false;
992  }
993 
994  myobs->setVal(myexp->getVal());
995 
996  if (fgPrintLevel > 2) {
997  std::cout << "SetObsToExpected : setting " << myobs->GetName() << " to expected value " << myexp->getVal() << " of " << myexp->GetName() << std::endl;
998  }
999 
1000  return true;
1001 }
1002 
1003 ////////////////////////////////////////////////////////////////////////////////
1004 /// generate counting Asimov data for the case when the pdf cannot be extended
1005 /// assume pdf is a RooPoisson or can be decomposed in a product of RooPoisson,
1006 /// otherwise we cannot know how to make the Asimov data sets in the other cases
1007 
1009  RooArgSet obs(observables);
1010  RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
1011  RooPoisson *pois = 0;
1012  RooGaussian *gaus = 0;
1013 
1014  if (fgPrintLevel > 1)
1015  std::cout << "generate counting Asimov data for pdf of type " << pdf.IsA()->GetName() << std::endl;
1016 
1017  bool r = false;
1018  if (prod != 0) {
1019  r = SetObsToExpected(*prod, observables);
1020  } else if ((pois = dynamic_cast<RooPoisson *>(&pdf)) != 0) {
1021  r = SetObsToExpected(*pois, observables);
1022  // we need in this case to set Poisson to real values
1023  pois->setNoRounding(true);
1024  } else if ((gaus = dynamic_cast<RooGaussian *>(&pdf)) != 0) {
1025  r = SetObsToExpected(*gaus, observables);
1026  } else {
1027  oocoutE((TObject*)0,InputArguments) << "A counting model pdf must be either a RooProdPdf or a RooPoisson or a RooGaussian" << endl;
1028  }
1029  if (!r) return 0;
1030  int icat = 0;
1031  if (channelCat) {
1032  icat = channelCat->getIndex();
1033  }
1034 
1035  RooDataSet *ret = new RooDataSet(TString::Format("CountingAsimovData%d",icat),TString::Format("CountingAsimovData%d",icat), obs);
1036  ret->add(obs);
1037  return ret;
1038 }
1039 
1040 ////////////////////////////////////////////////////////////////////////////////
1041 /// compute the asimov data set for an observable of a pdf
1042 /// use the number of bins sets in the observables
1043 /// to do : (possibility to change number of bins)
1044 /// implement integration over bin content
1045 
1046 RooAbsData * AsymptoticCalculator::GenerateAsimovDataSinglePdf(const RooAbsPdf & pdf, const RooArgSet & allobs, const RooRealVar & weightVar, RooCategory * channelCat) {
1047 
1048  int printLevel = fgPrintLevel;
1049 
1050  // Get observables defined by the pdf associated with this state
1051  std::unique_ptr<RooArgSet> obs(pdf.getObservables(allobs) );
1052 
1053 
1054  // if pdf cannot be extended assume is then a counting experiment
1055  if (!pdf.canBeExtended() ) return GenerateCountingAsimovData(const_cast<RooAbsPdf&>(pdf), *obs, weightVar, channelCat);
1056 
1057  RooArgSet obsAndWeight(*obs);
1058  obsAndWeight.add(weightVar);
1059 
1060  RooDataSet* asimovData = 0;
1061  if (channelCat) {
1062  int icat = channelCat->getIndex();
1063  asimovData = new RooDataSet(TString::Format("AsimovData%d",icat),TString::Format("combAsimovData%d",icat),
1064  RooArgSet(obsAndWeight,*channelCat),RooFit::WeightVar(weightVar));
1065  }
1066  else
1067  asimovData = new RooDataSet("AsimovData","AsimovData",RooArgSet(obsAndWeight),RooFit::WeightVar(weightVar));
1068 
1069  // This works only for 1D observables
1070  //RooRealVar* thisObs = ((RooRealVar*)obstmp->first());
1071 
1072  RooArgList obsList(*obs);
1073 
1074  // loop on observables and on the bins
1075  if (printLevel >= 2) {
1076  cout << "Generating Asimov data for pdf " << pdf.GetName() << endl;
1077  cout << "list of observables " << endl;
1078  obsList.Print();
1079  }
1080 
1081  int obsIndex = 0;
1082  double binVolume = 1;
1083  int nbins = 0;
1084  FillBins(pdf, obsList, *asimovData, obsIndex, binVolume, nbins);
1085  if (printLevel >= 2)
1086  cout << "filled from " << pdf.GetName() << " " << nbins << " nbins " << " volume is " << binVolume << endl;
1087 
1088  // for (int iobs = 0; iobs < obsList.getSize(); ++iobs) {
1089  // RooRealVar * thisObs = dynamic_cast<RooRealVar*> &obsList[i];
1090  // if (thisObs == 0) continue;
1091  // // loop on the bin contents
1092  // for(int ibin=0; ibin<thisObs->numBins(); ++ibin){
1093  // thisObs->setBin(ibin);
1094 
1095  // thisNorm=pdftmp->getVal(obstmp)*thisObs->getBinWidth(jj);
1096  // if (thisNorm*expectedEvents <= 0)
1097  // {
1098  // cout << "WARNING::Detected bin with zero expected events! Please check your inputs." << endl;
1099  // }
1100  // // have a cut off for overflows ??
1101  // obsDataUnbinned->add(*mc->GetObservables(), thisNorm*expectedEvents);
1102  // }
1103 
1104  if (printLevel >= 1)
1105  {
1106  asimovData->Print();
1107  //cout <<"sum entries "<< asimovData->sumEntries()<<endl;
1108  }
1109  if( TMath::IsNaN(asimovData->sumEntries()) ){
1110  cout << "sum entries is nan"<<endl;
1111  assert(0);
1112  delete asimovData;
1113  asimovData = 0;
1114  }
1115 
1116  return asimovData;
1117 
1118 }
1119 
1120 ////////////////////////////////////////////////////////////////////////////////
1121 /// generate the asimov data for the observables (not the global ones)
1122 /// need to deal with the case of a sim pdf
1123 
1125 
1126  int printLevel = fgPrintLevel;
1127 
1128  RooRealVar * weightVar = new RooRealVar("binWeightAsimov", "binWeightAsimov", 1, 0, 1.E30 );
1129 
1130  if (printLevel > 1) cout <<" Generate Asimov data for observables"<<endl;
1131  //RooDataSet* simData=NULL;
1132  const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf);
1133  if (!simPdf) {
1134  // generate data for non sim pdf
1135  return GenerateAsimovDataSinglePdf( pdf, observables, *weightVar, 0);
1136  }
1137 
1138  std::map<std::string, RooDataSet*> asimovDataMap;
1139 
1140  //look at category of simpdf
1141  RooCategory& channelCat = (RooCategory&)simPdf->indexCat();
1142  int nrIndices = channelCat.numTypes();
1143  if( nrIndices == 0 ) {
1144  oocoutW((TObject*)0,Generation) << "Simultaneous pdf does not contain any categories." << endl;
1145  }
1146  for (int i=0;i<nrIndices;i++){
1147  channelCat.setIndex(i);
1148  //iFrame++;
1149  // Get pdf associated with state from simpdf
1150  RooAbsPdf* pdftmp = simPdf->getPdf(channelCat.getLabel()) ;
1151  assert(pdftmp != 0);
1152 
1153  if (printLevel > 1)
1154  {
1155  cout << "on type " << channelCat.getLabel() << " " << channelCat.getIndex() << endl;
1156  }
1157 
1158  RooAbsData * dataSinglePdf = GenerateAsimovDataSinglePdf( *pdftmp, observables, *weightVar, &channelCat);
1159  //((RooRealVar*)obstmp->first())->Print();
1160  //cout << "expected events " << pdftmp->expectedEvents(*obstmp) << endl;
1161  if (!dataSinglePdf) {
1162  oocoutE((TObject*)0,Generation) << "Error generating an Asimov data set for pdf " << pdftmp->GetName() << endl;
1163  return 0;
1164  }
1165 
1166 
1167  asimovDataMap[string(channelCat.getLabel())] = (RooDataSet*) dataSinglePdf;
1168 
1169  if (printLevel > 1)
1170  {
1171  cout << "channel: " << channelCat.getLabel() << ", data: ";
1172  dataSinglePdf->Print();
1173  cout << endl;
1174  }
1175  }
1176 
1177  RooArgSet obsAndWeight(observables);
1178  obsAndWeight.add(*weightVar);
1179 
1180 
1181  RooDataSet* asimovData = new RooDataSet("asimovDataFullModel","asimovDataFullModel",RooArgSet(obsAndWeight,channelCat),
1182  RooFit::Index(channelCat),RooFit::Import(asimovDataMap),RooFit::WeightVar(*weightVar));
1183 
1184  delete weightVar;
1185  return asimovData;
1186 
1187 }
1188 
1189 ////////////////////////////////////////////////////////////////////////////////
1190 /// static function to the an Asimov data set
1191 /// given an observed dat set, a model and a snapshot of poi.
1192 /// Return the asimov data set + global observables set to values satisfying the constraints
1193 
1194 RooAbsData * AsymptoticCalculator::MakeAsimovData(RooAbsData & realData, const ModelConfig & model, const RooArgSet & paramValues, RooArgSet & asimovGlobObs, const RooArgSet * genPoiValues ) {
1195 
1196  int verbose = fgPrintLevel;
1197 
1198 
1199  RooArgSet poi(*model.GetParametersOfInterest());
1200  poi = paramValues;
1201  //RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
1202  // set poi constant for conditional MLE
1203  // need to fit nuisance parameters at their conditional MLE value
1204  RooLinkedListIter it = poi.iterator();
1205  RooRealVar* tmpPar = NULL;
1206  RooArgSet paramsSetConstant;
1207  while((tmpPar = (RooRealVar*)it.Next())){
1208  tmpPar->setConstant();
1209  if (verbose>0)
1210  std::cout << "MakeAsimov: Setting poi " << tmpPar->GetName() << " to a constant value = " << tmpPar->getVal() << std::endl;
1211  paramsSetConstant.add(*tmpPar);
1212  }
1213 
1214  // find conditional value of the nuisance parameters
1215  bool hasFloatParams = false;
1216  RooArgSet constrainParams;
1217  if (model.GetNuisanceParameters()) {
1218  constrainParams.add(*model.GetNuisanceParameters());
1219  RooStats::RemoveConstantParameters(&constrainParams);
1220  if (constrainParams.getSize() > 0) hasFloatParams = true;
1221 
1222  } else {
1223  // Do we have free parameters anyway that need fitting?
1224  std::unique_ptr<RooArgSet> params(model.GetPdf()->getParameters(realData));
1225  RooLinkedListIter iter(params->iterator());
1226  for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
1227  RooRealVar *rrv = dynamic_cast<RooRealVar *>(a);
1228  if ( rrv != 0 && rrv->isConstant() == false ) { hasFloatParams = true; break; }
1229  }
1230  }
1231  if (hasFloatParams) {
1232  // models need to be fitted to find best nuisance parameter values
1233 
1234  TStopwatch tw2; tw2.Start();
1235  int minimPrintLevel = ROOT::Math::MinimizerOptions::DefaultPrintLevel();
1236  if (verbose>0) {
1237  std::cout << "MakeAsimov: doing a conditional fit for finding best nuisance values " << std::endl;
1238  minimPrintLevel = verbose;
1239  if (verbose>1) {
1240  std::cout << "POI values:\n"; poi.Print("v");
1241  if (verbose > 2) {
1242  std::cout << "Nuis param values:\n";
1243  constrainParams.Print("v");
1244  }
1245  }
1246  }
1249 
1250  RooArgSet conditionalObs;
1251  if (model.GetConditionalObservables()) conditionalObs.add(*model.GetConditionalObservables());
1252  RooArgSet globalObs;
1253  if (model.GetGlobalObservables()) globalObs.add(*model.GetGlobalObservables());
1254 
1255  std::string minimizerType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
1256  std::string minimizerAlgo = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo();
1257  model.GetPdf()->fitTo(realData, RooFit::Minimizer(minimizerType.c_str(),minimizerAlgo.c_str()),
1259  RooFit::PrintLevel(minimPrintLevel-1), RooFit::Hesse(false),
1260  RooFit::Constrain(constrainParams),RooFit::GlobalObservables(globalObs),
1262  if (verbose>0) { std::cout << "fit time "; tw2.Print();}
1263  if (verbose > 1) {
1264  // after the fit the nuisance parameters will have their best fit value
1265  if (model.GetNuisanceParameters() ) {
1266  std::cout << "Nuisance parameters after fit for asimov dataset: " << std::endl;
1267  model.GetNuisanceParameters()->Print("V");
1268  }
1269  }
1270 
1271  if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
1272 
1273  }
1274 
1275  // restore the parameters which were set constant
1276  SetAllConstant(paramsSetConstant, false);
1277 
1278 
1279  RooArgSet * allParams = model.GetPdf()->getParameters(realData);
1281 
1282  // if a RooArgSet of poi is passed , different poi will be used for generating the Asimov data set
1283  if (genPoiValues) *allParams = *genPoiValues;
1284 
1285  // now do the actual generation of the AsimovData Set
1286  // no need to pass parameters values since we have set them before
1287  RooAbsData * asimovData = MakeAsimovData(model, *allParams, asimovGlobObs);
1288 
1289  delete allParams;
1290 
1291  return asimovData;
1292 
1293 }
1294 
1295 ////////////////////////////////////////////////////////////////////////////////
1296 /// static function to the an Asimov data set
1297 /// given the model and the values of all parameters including the nuisance
1298 /// Return the asimov data set + global observables set to values satisfying the constraints
1299 
1300 RooAbsData * AsymptoticCalculator::MakeAsimovData(const ModelConfig & model, const RooArgSet & allParamValues, RooArgSet & asimovGlobObs) {
1301 
1302  int verbose = fgPrintLevel;
1303 
1304  TStopwatch tw;
1305  tw.Start();
1306 
1307  // set the parameter values (do I need the poi to be constant ? )
1308  // the nuisance parameter values could be set at their fitted value (the MLE)
1309  if (allParamValues.getSize() > 0) {
1310  RooArgSet * allVars = model.GetPdf()->getVariables();
1311  *allVars = allParamValues;
1312  delete allVars;
1313  }
1314 
1315 
1316  // generate the Asimov data set for the observables
1317  RooAbsData * asimov = GenerateAsimovData(*model.GetPdf() , *model.GetObservables() );
1318 
1319  if (verbose>0) {
1320  std::cout << "Generated Asimov data for observables "; (model.GetObservables() )->Print();
1321  if (verbose > 1) {
1322  if (asimov->numEntries() == 1 ) {
1323  std::cout << "--- Asimov data values \n";
1324  asimov->get()->Print("v");
1325  }
1326  else {
1327  std::cout << "--- Asimov data numEntries = " << asimov->numEntries() << " sumOfEntries = " << asimov->sumEntries() << std::endl;
1328  }
1329  std::cout << "\ttime for generating : "; tw.Print();
1330  }
1331  }
1332 
1333 
1334  // Now need to have in ASIMOV the data sets also the global observables
1335  // Their values must be the one satisfying the constraint.
1336  // to do it make a nuisance pdf with all product of constraints and then
1337  // assign to each constraint a glob observable value = to the current fitted nuisance parameter value
1338  // IN general one should solve in general the system of equations f( gobs| nuispar ) = 0 where f are the
1339  // derivatives of the constraint with respect the nuisance parameter and they are evaluated at the best fit nuisance
1340  // parameter points
1341  // As simple solution assume that constrain has a direct dependence on the nuisance parameter, i.e.
1342  // Constraint (gobs, func( nuispar) ) and the condition is satisfied for
1343  // gobs = func( nuispar) where nunispar is at the MLE value
1344 
1345 
1346  if (model.GetGlobalObservables() && model.GetGlobalObservables()->getSize() > 0) {
1347 
1348  if (verbose>1) {
1349  std::cout << "Generating Asimov data for global observables " << std::endl;
1350  }
1351 
1352  RooArgSet gobs(*model.GetGlobalObservables());
1353 
1354  // snapshot data global observables
1355  RooArgSet snapGlobalObsData;
1356  SetAllConstant(gobs, true);
1357  gobs.snapshot(snapGlobalObsData);
1358 
1359 
1360  RooArgSet nuis;
1361  if (model.GetNuisanceParameters()) nuis.add(*model.GetNuisanceParameters());
1362  if (nuis.getSize() == 0) {
1363  oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: model does not have nuisance parameters but has global observables"
1364  << " set global observables to model values " << endl;
1365  asimovGlobObs = gobs;
1366  return asimov;
1367  }
1368 
1369  // part 1: create the nuisance pdf
1370  std::unique_ptr<RooAbsPdf> nuispdf(RooStats::MakeNuisancePdf(model,"TempNuisPdf") );
1371  if (nuispdf.get() == 0) {
1372  oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: model has nuisance parameters and global obs but no nuisance pdf "
1373  << std::endl;
1374  }
1375  // unfold the nuisance pdf if it is a prod pdf
1376  RooArgList pdfList;
1377  RooProdPdf *prod = dynamic_cast<RooProdPdf *>(nuispdf.get());
1378  if (prod ) {
1379  pdfList.add(prod->pdfList());
1380  }
1381  else
1382  // nothing to unfold - just use the pdf
1383  pdfList.add(*nuispdf.get());
1384 
1385  RooLinkedListIter iter(pdfList.iterator());
1386  for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
1387  RooAbsPdf *cterm = dynamic_cast<RooAbsPdf *>(a);
1388  assert(cterm && "AsimovUtils: a factor of the nuisance pdf is not a Pdf!");
1389  if (!cterm->dependsOn(nuis)) continue; // dummy constraints
1390  // skip also the case of uniform components
1391  if (typeid(*cterm) == typeid(RooUniform)) continue;
1392 
1393  std::unique_ptr<RooArgSet> cpars(cterm->getParameters(&gobs));
1394  std::unique_ptr<RooArgSet> cgobs(cterm->getObservables(&gobs));
1395  if (cgobs->getSize() > 1) {
1396  oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1397  << " has multiple global observables -cannot generate - skip it" << std::endl;
1398  continue;
1399  }
1400  else if (cgobs->getSize() == 0) {
1401  oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1402  << " has no global observables - skip it" << std::endl;
1403  continue;
1404  }
1405  // the variable representing the global observable
1406  RooRealVar &rrv = dynamic_cast<RooRealVar &>(*cgobs->first());
1407 
1408  // remove the constant parameters in cpars
1410  if (cpars->getSize() != 1) {
1411  oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1412  << cterm->GetName() << " has multiple floating params - cannot generate - skip it " << std::endl;
1413  continue;
1414  }
1415 
1416  bool foundServer = false;
1417  // note : this will work only for this type of constraints
1418  // expressed as RooPoisson, RooGaussian, RooLognormal, RooGamma
1419  TClass * cClass = cterm->IsA();
1420  if (verbose > 2) std::cout << "Constraint " << cterm->GetName() << " of type " << cClass->GetName() << std::endl;
1421  if ( cClass != RooGaussian::Class() && cClass != RooPoisson::Class() &&
1422  cClass != RooGamma::Class() && cClass != RooLognormal::Class() &&
1423  cClass != RooBifurGauss::Class() ) {
1424  TString className = (cClass) ? cClass->GetName() : "undefined";
1425  oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1426  << cterm->GetName() << " of type " << className
1427  << " is a non-supported type - result might be not correct " << std::endl;
1428  }
1429 
1430  // in case of a Poisson constraint make sure the rounding is not set
1431  if (cClass == RooPoisson::Class() ) {
1432  RooPoisson * pois = dynamic_cast<RooPoisson*>(cterm);
1433  assert(pois);
1434  pois->setNoRounding(true);
1435  }
1436 
1437  // look at server of the constraint term and check if the global observable is part of the server
1438  RooAbsArg * arg = cterm->findServer(rrv);
1439  if (!arg) {
1440  // special case is for the Gamma where one might define the global observable n and you have a Gamma(b, n+1, ...._
1441  // 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
1442  // so in case of the Gamma ignore this test
1443  if ( cClass != RooGamma::Class() ) {
1444  oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1445  << cterm->GetName() << " has no direct dependence on global observable- cannot generate it " << std::endl;
1446  continue;
1447  }
1448  }
1449 
1450  // loop on the server of the constraint term
1451  // need to treat the Gamma as a special case
1452  // the mode of the Gamma is (k-1)*theta where theta is the inverse of the rate parameter.
1453  // 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
1454  RooAbsReal * thetaGamma = 0;
1455  if ( cClass == RooGamma::Class() ) {
1456  RooFIter itc(cterm->serverMIterator() );
1457  for (RooAbsArg *a2 = itc.next(); a2 != 0; a2 = itc.next()) {
1458  if (TString(a2->GetName()).Contains("theta") ) {
1459  thetaGamma = dynamic_cast<RooAbsReal*>(a2);
1460  break;
1461  }
1462  }
1463  if (thetaGamma == 0) {
1464  oocoutI((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1465  << cterm->GetName() << " is a Gamma distribution and no server named theta is found. Assume that the Gamma scale is 1 " << std::endl;
1466  }
1467  else {
1468  if (verbose>2)
1469  std::cout << "Gamma constraint has a scale " << thetaGamma->GetName() << " = " << thetaGamma->getVal() << std::endl;
1470  }
1471  }
1472  RooFIter iter2(cterm->serverMIterator() );
1473  for (RooAbsArg *a2 = iter2.next(); a2 != 0; a2 = iter2.next()) {
1474  RooAbsReal * rrv2 = dynamic_cast<RooAbsReal *>(a2);
1475  if (verbose > 2) std::cout << "Loop on constraint server term " << a2->GetName() << std::endl;
1476  if (rrv2 && rrv2->dependsOn(nuis) ) {
1477 
1478 
1479  // found server depending on nuisance
1480  if (foundServer) {
1481  oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1482  << cterm->GetName() << " constraint term has more server depending on nuisance- cannot generate it " <<
1483  std::endl;
1484  foundServer = false;
1485  break;
1486  }
1487  if (thetaGamma && thetaGamma->getVal() > 0)
1488  rrv.setVal( rrv2->getVal() / thetaGamma->getVal() );
1489  else
1490  rrv.setVal( rrv2->getVal() );
1491  foundServer = true;
1492 
1493  if (verbose>2)
1494  std::cout << "setting global observable "<< rrv.GetName() << " to value " << rrv.getVal()
1495  << " which comes from " << rrv2->GetName() << std::endl;
1496  }
1497  }
1498 
1499  if (!foundServer) {
1500  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;
1501  std::cerr << "Parameters: " << std::endl;
1502  cpars->Print("V");
1503  std::cerr << "Observables: " << std::endl;
1504  cgobs->Print("V");
1505  }
1506 // rrv.setVal(match->getVal());
1507  }
1508 
1509  // make a snapshot of global observables
1510  // needed this ?? (LM)
1511 
1512  asimovGlobObs.removeAll();
1513  SetAllConstant(gobs, true);
1514  gobs.snapshot(asimovGlobObs);
1515 
1516  // revert global observables to the data value
1517  gobs = snapGlobalObsData;
1518 
1519  if (verbose>0) {
1520  std::cout << "Generated Asimov data for global observables ";
1521  if (verbose == 1) gobs.Print();
1522  }
1523 
1524  if (verbose > 1) {
1525  std::cout << "\nGlobal observables for data: " << std::endl;
1526  gobs.Print("V");
1527  std::cout << "\nGlobal observables for asimov: " << std::endl;
1528  asimovGlobObs.Print("V");
1529  }
1530 
1531 
1532  }
1533 
1534  return asimov;
1535 
1536 }
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Definition: RooAbsPdf.cxx:778
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:2078
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
static double EvaluateNLL(RooAbsPdf &pdf, RooAbsData &data, const RooArgSet *condObs, const RooArgSet *globObs, const RooArgSet *poiSet=0)
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:740
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:79
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:783
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)
Bool_t IsNaN(Double_t x)
Definition: TMath.h:777
static RooMsgService & instance()
Return reference to singleton instance.
STL namespace.
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:2365
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:2953
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
RooCmdArg GlobalObservables(const RooArgSet &globs)
double Root() const
Returns root value.
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
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
ROOT::R::TRInterface & r
Definition: Object.C:4
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)
auto * a
Definition: textangle.C:12
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:98
void setGlobalKillBelow(RooFit::MsgLevel level)
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:75
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:88
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:159
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:359
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:537
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
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
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:1079
virtual Int_t getBins(const char *name=0) const
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:87
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:285
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...