ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HypoTestInverter.cxx
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
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 ////////////////////////////////////////////////////////////////////////////////
12 
13 
14 // include other header files
15 
16 #include "RooAbsData.h"
17 #
18 #include "TMath.h"
19 
20 #include "RooStats/HybridResult.h"
21 
23 #include <cassert>
24 #include <cmath>
25 
26 #include "TF1.h"
27 #include "TFile.h"
28 #include "TH1.h"
29 #include "TLine.h"
30 #include "TCanvas.h"
31 #include "TGraphErrors.h"
32 #include "RooRealVar.h"
33 #include "RooArgSet.h"
34 #include "RooAbsPdf.h"
35 #include "RooRandom.h"
36 #include "RooConstVar.h"
37 #include "RooMsgService.h"
38 #include "RooStats/ModelConfig.h"
45 #include "RooStats/ToyMCSampler.h"
46 #include "RooStats/HypoTestPlot.h"
48 
49 #include "RooStats/ProofConfig.h"
50 
52 
53 using namespace RooStats;
54 using namespace std;
55 
56 // static variable definitions
57 double HypoTestInverter::fgCLAccuracy = 0.005;
58 unsigned int HypoTestInverter::fgNToys = 500;
59 
60 double HypoTestInverter::fgAbsAccuracy = 0.05;
61 double HypoTestInverter::fgRelAccuracy = 0.05;
62 std::string HypoTestInverter::fgAlgo = "logSecant";
63 
64 bool HypoTestInverter::fgCloseProof = false;
65 
66 // helper class to wrap the functionality of the various HypoTestCalculators
67 
68 template<class HypoTestType>
69 struct HypoTestWrapper {
70 
71  static void SetToys(HypoTestType * h, int toyNull, int toyAlt) { h->SetToys(toyNull,toyAlt); }
72 
73 };
74 
75 
76 void HypoTestInverter::SetCloseProof(Bool_t flag) {
77 // set flag to close proof for every new run
78  fgCloseProof = flag;
79 }
80 
81 
82 RooRealVar * HypoTestInverter::GetVariableToScan(const HypoTestCalculatorGeneric &hc) {
83  // get the variable to scan
84  // try first with null model if not go to alternate model
85 
86  RooRealVar * varToScan = 0;
87  const ModelConfig * mc = hc.GetNullModel();
88  if (mc) {
89  const RooArgSet * poi = mc->GetParametersOfInterest();
90  if (poi) varToScan = dynamic_cast<RooRealVar*> (poi->first() );
91  }
92  if (!varToScan) {
93  mc = hc.GetAlternateModel();
94  if (mc) {
95  const RooArgSet * poi = mc->GetParametersOfInterest();
96  if (poi) varToScan = dynamic_cast<RooRealVar*> (poi->first() );
97  }
98  }
99  return varToScan;
100 }
101 
102 void HypoTestInverter::CheckInputModels(const HypoTestCalculatorGeneric &hc,const RooRealVar & scanVariable) {
103  // check the model given the given hypotestcalculator
104  const ModelConfig * modelSB = hc.GetNullModel();
105  const ModelConfig * modelB = hc.GetAlternateModel();
106  if (!modelSB || ! modelB)
107  oocoutF((TObject*)0,InputArguments) << "HypoTestInverter - model are not existing" << std::endl;
108  assert(modelSB && modelB);
109 
110  oocoutI((TObject*)0,InputArguments) << "HypoTestInverter ---- Input models: \n"
111  << "\t\t using as S+B (null) model : "
112  << modelSB->GetName() << "\n"
113  << "\t\t using as B (alternate) model : "
114  << modelB->GetName() << "\n" << std::endl;
115 
116  // check if scanVariable is included in B model pdf
117  RooAbsPdf * bPdf = modelB->GetPdf();
118  const RooArgSet * bObs = modelB->GetObservables();
119  if (!bPdf || !bObs) {
120  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - B model has no pdf or observables defined" << std::endl;
121  return;
122  }
123  RooArgSet * bParams = bPdf->getParameters(*bObs);
124  if (!bParams) {
125  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - pdf of B model has no parameters" << std::endl;
126  return;
127  }
128  if (bParams->find(scanVariable.GetName() ) ) {
129  const RooArgSet * poiB = modelB->GetSnapshot();
130  if (!poiB || !poiB->find(scanVariable.GetName()) ||
131  ( (RooRealVar*) poiB->find(scanVariable.GetName()) )->getVal() != 0 )
132  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter - using a B model with POI "
133  << scanVariable.GetName() << " not equal to zero "
134  << " user must check input model configurations " << endl;
135  if (poiB) delete poiB;
136  }
137  delete bParams;
138 }
139 
140 HypoTestInverter::HypoTestInverter( ) :
141  fTotalToysRun(0),
142  fMaxToys(0),
143  fCalculator0(0),
144  fScannedVariable(0),
145  fResults(0),
146  fUseCLs(false),
147  fScanLog(false),
148  fSize(0),
149  fVerbose(0),
150  fCalcType(kUndefined),
151  fNBins(0), fXmin(1), fXmax(1),
152  fNumErr(0)
153 {
154  // default constructor (doesn't do anything)
155 }
156 
158  RooRealVar* scannedVariable, double size ) :
159  fTotalToysRun(0),
160  fMaxToys(0),
161  fCalculator0(0),
162  fScannedVariable(scannedVariable),
163  fResults(0),
164  fUseCLs(false),
165  fScanLog(false),
166  fSize(size),
167  fVerbose(0),
168  fCalcType(kUndefined),
169  fNBins(0), fXmin(1), fXmax(1),
170  fNumErr(0)
171 {
172  // Constructor from a HypoTestCalculatorGeneric
173  // The HypoTest calculator must be a FrequentistCalculator or HybridCalculator type
174  // Other type of calculators are not supported.
175  // The calculator must be created before by using the S+B model for the null and
176  // the B model for the alt
177  // If no variable to scan are given they are assumed to be the first variable
178  // from the parameter of interests of the null model
179 
180  if (!fScannedVariable) {
182  }
183  if (!fScannedVariable)
184  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
185  else
187 
188  HybridCalculator * hybCalc = dynamic_cast<HybridCalculator*>(&hc);
189  if (hybCalc) {
190  fCalcType = kHybrid;
191  fCalculator0 = hybCalc;
192  return;
193  }
194  FrequentistCalculator * freqCalc = dynamic_cast<FrequentistCalculator*>(&hc);
195  if (freqCalc) {
197  fCalculator0 = freqCalc;
198  return;
199  }
200  AsymptoticCalculator * asymCalc = dynamic_cast<AsymptoticCalculator*>(&hc);
201  if (asymCalc) {
203  fCalculator0 = asymCalc;
204  return;
205  }
206  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - Type of hypotest calculator is not supported " <<std::endl;
207  fCalculator0 = &hc;
208 }
209 
210 
211 HypoTestInverter::HypoTestInverter( HybridCalculator& hc,
212  RooRealVar* scannedVariable, double size ) :
213  fTotalToysRun(0),
214  fMaxToys(0),
215  fCalculator0(&hc),
216  fScannedVariable(scannedVariable),
217  fResults(0),
218  fUseCLs(false),
219  fScanLog(false),
220  fSize(size),
221  fVerbose(0),
222  fCalcType(kHybrid),
223  fNBins(0), fXmin(1), fXmax(1),
224  fNumErr(0)
225 {
226  // Constructor from a reference to a HybridCalculator
227  // The calculator must be created before by using the S+B model for the null and
228  // the B model for the alt
229  // If no variable to scan are given they are assumed to be the first variable
230  // from the parameter of interests of the null model
231 
232  if (!fScannedVariable) {
234  }
235  if (!fScannedVariable)
236  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
237  else
239 
240 }
241 
242 
243 
244 
245 HypoTestInverter::HypoTestInverter( FrequentistCalculator& hc,
246  RooRealVar* scannedVariable, double size ) :
247  fTotalToysRun(0),
248  fMaxToys(0),
249  fCalculator0(&hc),
250  fScannedVariable(scannedVariable),
251  fResults(0),
252  fUseCLs(false),
253  fScanLog(false),
254  fSize(size),
255  fVerbose(0),
256  fCalcType(kFrequentist),
257  fNBins(0), fXmin(1), fXmax(1),
258  fNumErr(0)
259 {
260  // Constructor from a reference to a FrequentistCalculator
261  // The calculator must be created before by using the S+B model for the null and
262  // the B model for the alt
263  // If no variable to scan are given they are assumed to be the first variable
264  // from the parameter of interests of the null model
265 
266  if (!fScannedVariable) {
268  }
269  if (!fScannedVariable)
270  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
271  else
273 }
274 
275 HypoTestInverter::HypoTestInverter( AsymptoticCalculator& hc,
276  RooRealVar* scannedVariable, double size ) :
277  fTotalToysRun(0),
278  fMaxToys(0),
279  fCalculator0(&hc),
280  fScannedVariable(scannedVariable),
281  fResults(0),
282  fUseCLs(false),
283  fScanLog(false),
284  fSize(size),
285  fVerbose(0),
286  fCalcType(kAsymptotic),
287  fNBins(0), fXmin(1), fXmax(1),
288  fNumErr(0)
289 {
290  // Constructor from a reference to a AsymptoticCalculator
291  // The calculator must be created before by using the S+B model for the null and
292  // the B model for the alt
293  // If no variable to scan are given they are assumed to be the first variable
294  // from the parameter of interests of the null model
295 
296  if (!fScannedVariable) {
298  }
299  if (!fScannedVariable)
300  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
301  else
303 
304 }
305 
306 
307 ////////////////////////////////////////////////////////////////////////////////
308 /// Constructor from a model for B model and a model for S+B.
309 /// An HypoTestCalculator (Hybrid of Frequentis) will be created using the
310 /// S+B model as the null and the B model as the alternate
311 /// If no variable to scan are given they are assumed to be the first variable
312 /// from the parameter of interests of the null model
313 
314 HypoTestInverter::HypoTestInverter( RooAbsData& data, ModelConfig &sbModel, ModelConfig &bModel,
315  RooRealVar * scannedVariable, ECalculatorType type, double size) :
316  fTotalToysRun(0),
317  fMaxToys(0),
318  fCalculator0(0),
319  fScannedVariable(scannedVariable),
320  fResults(0),
321  fUseCLs(false),
322  fScanLog(false),
323  fSize(size),
324  fVerbose(0),
325  fCalcType(type),
326  fNBins(0), fXmin(1), fXmax(1),
327  fNumErr(0)
328 {
329  if(fCalcType==kFrequentist) fHC.reset(new FrequentistCalculator(data, bModel, sbModel));
330  if(fCalcType==kHybrid) fHC.reset( new HybridCalculator(data, bModel, sbModel)) ;
331  if(fCalcType==kAsymptotic) fHC.reset( new AsymptoticCalculator(data, bModel, sbModel));
332  fCalculator0 = fHC.get();
333  // get scanned variabke
334  if (!fScannedVariable) {
336  }
337  if (!fScannedVariable)
338  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
339  else
341 
342 }
343 
344 
347  fTotalToysRun(0),
348  fCalculator0(0), fScannedVariable(0), // add these for Coverity
349  fResults(0)
350 {
351  // copy-constructor
352  // NOTE: this class does not copy the contained result and
353  // the HypoTestCalculator, but only the pointers
354  // It requires the original HTI to be alive
355  (*this) = rhs;
356 }
357 
359  // assignment operator
360  // NOTE: this class does not copy the contained result and
361  // the HypoTestCalculator, but only the pointers
362  // It requires the original HTI to be alive
363  if (this == &rhs) return *this;
364  fTotalToysRun = 0;
365  fMaxToys = rhs.fMaxToys;
368  fUseCLs = rhs.fUseCLs;
369  fScanLog = rhs.fScanLog;
370  fSize = rhs.fSize;
371  fVerbose = rhs.fVerbose;
372  fCalcType = rhs.fCalcType;
373  fNBins = rhs.fNBins;
374  fXmin = rhs.fXmin;
375  fXmax = rhs.fXmax;
376  fNumErr = rhs.fNumErr;
377 
378  return *this;
379 }
380 
382 {
383  // destructor (delete the HypoTestInverterResult)
384 
385  if (fResults) delete fResults;
386  fCalculator0 = 0;
387 }
388 
389 
391 {
392  // return the test statistic which is or will be used by the class
395  }
396  else
397  return 0;
398 }
399 
401 {
402  // set the test statistic to use
405  return true;
406  }
407  else return false;
408 }
409 
411  // delete contained result and graph
412  if (fResults) delete fResults;
413  fResults = 0;
414  fLimitPlot.reset(nullptr);
415 }
416 
418  // create a new HypoTestInverterResult to hold all computed results
419  if (fResults == 0) {
420  TString results_name = "result_";
421  results_name += fScannedVariable->GetName();
423  TString title = "HypoTestInverter Result For ";
424  title += fScannedVariable->GetName();
425  fResults->SetTitle(title);
426  }
429  // check if one or two sided scan
430  if (fCalculator0) {
431  // if asymptotic calculator
433  if (ac)
435  else {
436  // in case of the other calculators
438  if (sampler) {
440  if (pl && pl->IsTwoSided() ) fResults->fIsTwoSided = true;
441  }
442  }
443  }
444 }
445 
447  // Run a fixed scan or the automatic scan depending on the configuration
448  // Return if needed a copy of the result object which will be managed by the user
449 
450  // if having a result with at least one point return it
451  if (fResults && fResults->ArraySize() >= 1) {
452  oocoutI((TObject*)0,Eval) << "HypoTestInverter::GetInterval - return an already existing interval " << std::endl;
453  return (HypoTestInverterResult*)(fResults->Clone());
454  }
455 
456  if (fNBins > 0) {
457  oocoutI((TObject*)0,Eval) << "HypoTestInverter::GetInterval - run a fixed scan" << std::endl;
458  bool ret = RunFixedScan(fNBins, fXmin, fXmax, fScanLog);
459  if (!ret)
460  oocoutE((TObject*)0,Eval) << "HypoTestInverter::GetInterval - error running a fixed scan " << std::endl;
461  }
462  else {
463  oocoutI((TObject*)0,Eval) << "HypoTestInverter::GetInterval - run an automatic scan" << std::endl;
464  double limit(0),err(0);
465  bool ret = RunLimit(limit,err);
466  if (!ret)
467  oocoutE((TObject*)0,Eval) << "HypoTestInverter::GetInterval - error running an auto scan " << std::endl;
468  }
469 
471 
472  return (HypoTestInverterResult*) (fResults->Clone());
473 }
474 
475 
476 
477 HypoTestResult * HypoTestInverter::Eval(HypoTestCalculatorGeneric &hc, bool adaptive, double clsTarget) const {
478  // Run the Hypothesis test at a previous configured point
479  //(internal function called by RunOnePoint)
480 
481 
482  //for debug
483  // std::cout << ">>>>>>>>>>> " << std::endl;
484  // std::cout << "alternate model " << std::endl;
485  // hc.GetAlternateModel()->GetNuisanceParameters()->Print("V");
486  // hc.GetAlternateModel()->GetParametersOfInterest()->Print("V");
487  // std::cout << "Null model " << std::endl;
488  // hc.GetNullModel()->GetNuisanceParameters()->Print("V");
489  // hc.GetNullModel()->GetParametersOfInterest()->Print("V");
490  // std::cout << "<<<<<<<<<<<<<<< " << std::endl;
491 
492  // run the hypothesis test
493  HypoTestResult * hcResult = hc.GetHypoTest();
494  if (hcResult == 0) {
495  oocoutE((TObject*)0,Eval) << "HypoTestInverter::Eval - HypoTest failed" << std::endl;
496  return hcResult;
497  }
498 
499  // since the b model is the alt need to set the flag
500  hcResult->SetBackgroundAsAlt(true);
501 
502 
503  // bool flipPvalue = false;
504  // if (flipPValues)
505  // hcResult->SetPValueIsRightTail(!hcResult->GetPValueIsRightTail());
506 
507  // adjust for some numerical error in discrete models and == is not anymore
508  if (hcResult->GetPValueIsRightTail() )
509  hcResult->SetTestStatisticData(hcResult->GetTestStatisticData()-fNumErr); // issue with < vs <= in discrete models
510  else
511  hcResult->SetTestStatisticData(hcResult->GetTestStatisticData()+fNumErr); // issue with < vs <= in discrete models
512 
513  double clsMid = (fUseCLs ? hcResult->CLs() : hcResult->CLsplusb());
514  double clsMidErr = (fUseCLs ? hcResult->CLsError() : hcResult->CLsplusbError());
515 
516  //if (fVerbose) std::cout << (fUseCLs ? "\tCLs = " : "\tCLsplusb = ") << clsMid << " +/- " << clsMidErr << std::endl;
517 
518  if (adaptive) {
519 
520  if (fCalcType == kHybrid) HypoTestWrapper<HybridCalculator>::SetToys((HybridCalculator*)&hc, fUseCLs ? fgNToys : 1, 4*fgNToys);
521  if (fCalcType == kFrequentist) HypoTestWrapper<FrequentistCalculator>::SetToys((FrequentistCalculator*)&hc, fUseCLs ? fgNToys : 1, 4*fgNToys);
522 
523  while (clsMidErr >= fgCLAccuracy && (clsTarget == -1 || fabs(clsMid-clsTarget) < 3*clsMidErr) ) {
524  std::unique_ptr<HypoTestResult> more(hc.GetHypoTest());
525 
526  // if (flipPValues)
527  // more->SetPValueIsRightTail(!more->GetPValueIsRightTail());
528 
529  hcResult->Append(more.get());
530  clsMid = (fUseCLs ? hcResult->CLs() : hcResult->CLsplusb());
531  clsMidErr = (fUseCLs ? hcResult->CLsError() : hcResult->CLsplusbError());
532  if (fVerbose) std::cout << (fUseCLs ? "\tCLs = " : "\tCLsplusb = ") << clsMid << " +/- " << clsMidErr << std::endl;
533  }
534 
535  }
536  if (fVerbose ) {
537  oocoutP((TObject*)0,Eval) << "P values for " << fScannedVariable->GetName() << " = " <<
538  fScannedVariable->getVal() << "\n" <<
539  "\tCLs = " << hcResult->CLs() << " +/- " << hcResult->CLsError() << "\n" <<
540  "\tCLb = " << hcResult->CLb() << " +/- " << hcResult->CLbError() << "\n" <<
541  "\tCLsplusb = " << hcResult->CLsplusb() << " +/- " << hcResult->CLsplusbError() << "\n" <<
542  std::endl;
543  }
544 
545  if (fCalcType == kFrequentist || fCalcType == kHybrid) {
546  fTotalToysRun += (hcResult->GetAltDistribution()->GetSize() + hcResult->GetNullDistribution()->GetSize());
547 
548  // set sampling distribution name
549  TString nullDistName = TString::Format("%s_%s_%4.2f",hcResult->GetNullDistribution()->GetName(),
551  TString altDistName = TString::Format("%s_%s_%4.2f",hcResult->GetAltDistribution()->GetName(),
553 
554  hcResult->GetNullDistribution()->SetName(nullDistName);
555  hcResult->GetAltDistribution()->SetName(altDistName);
556  }
557 
558 
559 
560 
561  return hcResult;
562 }
563 
564 
565 bool HypoTestInverter::RunFixedScan( int nBins, double xMin, double xMax, bool scanLog ) const
566 {
567  // Run a Fixed scan in npoints between min and max
568 
569  CreateResults();
570  // interpolate the limits
571  fResults->fFittedLowerLimit = false;
572  fResults->fFittedUpperLimit = false;
573 
574  // safety checks
575  if ( nBins<=0 ) {
576  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::RunFixedScan - Please provide nBins>0\n";
577  return false;
578  }
579  if ( nBins==1 && xMin!=xMax ) {
580  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter::RunFixedScan - nBins==1 -> I will run for xMin (" << xMin << ")\n";
581  }
582  if ( xMin==xMax && nBins>1 ) {
583  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter::RunFixedScan - xMin==xMax -> I will enforce nBins==1\n";
584  nBins = 1;
585  }
586  if ( xMin>xMax ) {
587  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::RunFixedScan - Please provide xMin ("
588  << xMin << ") smaller that xMax (" << xMax << ")\n";
589  return false;
590  }
591 
592  if (xMin < fScannedVariable->getMin()) {
593  xMin = fScannedVariable->getMin();
594  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter::RunFixedScan - xMin < lower bound, use xmin = "
595  << xMin << std::endl;
596  }
597  if (xMax > fScannedVariable->getMax()) {
598  xMax = fScannedVariable->getMax();
599  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter::RunFixedScan - xMax > upper bound, use xmax = "
600  << xMax << std::endl;
601  }
602 
603  double thisX = xMin;
604  for (int i=0; i<nBins; i++) {
605 
606  if (i > 0) { // avoids case of nBins = 1
607  if (scanLog)
608  thisX = exp( log(xMin) + i*(log(xMax)-log(xMin))/(nBins-1) ); // scan in log x
609  else
610  thisX = xMin + i*(xMax-xMin)/(nBins-1); // linear scan in x
611  }
612 
613  bool status = RunOnePoint(thisX);
614 
615  // check if failed status
616  if ( status==false ) {
617  std::cout << "\t\tLoop interrupted because of failed status\n";
618  return false;
619  }
620  }
621 
622 
623 
624  return true;
625 }
626 
627 
628 bool HypoTestInverter::RunOnePoint( double rVal, bool adaptive, double clTarget) const
629 {
630  // run only one point at the given POI value
631 
632  CreateResults();
633 
634  // check if rVal is in the range specified for fScannedVariable
635  if ( rVal < fScannedVariable->getMin() ) {
636  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::RunOnePoint - Out of range: using the lower bound "
637  << fScannedVariable->getMin()
638  << " on the scanned variable rather than " << rVal<< "\n";
639  rVal = fScannedVariable->getMin();
640  }
641  if ( rVal > fScannedVariable->getMax() ) {
642  // print a message when you have a significative difference since rval is computed
643  if ( rVal > fScannedVariable->getMax()*(1.+1.E-12) )
644  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::RunOnePoint - Out of range: using the upper bound "
645  << fScannedVariable->getMax()
646  << " on the scanned variable rather than " << rVal<< "\n";
647  rVal = fScannedVariable->getMax();
648  }
649 
650  // save old value
651  double oldValue = fScannedVariable->getVal();
652 
653  // evaluate hybrid calculator at a single point
654  fScannedVariable->setVal(rVal);
655  // need to set value of rval in hybridcalculator
656  // assume null model is S+B and alternate is B only
657  const ModelConfig * sbModel = fCalculator0->GetNullModel();
658  RooArgSet poi; poi.add(*sbModel->GetParametersOfInterest());
659  // set poi to right values
660  poi = RooArgSet(*fScannedVariable);
661  const_cast<ModelConfig*>(sbModel)->SetSnapshot(poi);
662 
663  if (fVerbose > 0)
664  oocoutP((TObject*)0,Eval) << "Running for " << fScannedVariable->GetName() << " = " << fScannedVariable->getVal() << endl;
665 
666  // compute the results
667  HypoTestResult* result = Eval(*fCalculator0,adaptive,clTarget);
668  if (!result) {
669  oocoutE((TObject*)0,Eval) << "HypoTestInverter - Error running point " << fScannedVariable->GetName() << " = " <<
670  fScannedVariable->getVal() << endl;
671  return false;
672  }
673  // in case of a dummy result
674  if (TMath::IsNaN(result->NullPValue() ) && TMath::IsNaN(result->AlternatePValue() ) ) {
675  oocoutW((TObject*)0,Eval) << "HypoTestInverter - Skip invalid result for point " << fScannedVariable->GetName() << " = " <<
676  fScannedVariable->getVal() << endl;
677  return true; // need to return true to avoid breaking the scan loop
678  }
679 
680  double lastXtested;
681  if ( fResults->ArraySize()!=0 ) lastXtested = fResults->GetXValue(fResults->ArraySize()-1);
682  else lastXtested = -999;
683 
684  if ( (std::abs(rVal) < 1 && TMath::AreEqualAbs(rVal, lastXtested,1.E-12) ) ||
685  (std::abs(rVal) >= 1 && TMath::AreEqualRel(rVal, lastXtested,1.E-12) ) ) {
686 
687  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunOnePoint - Merge with previous result for "
688  << fScannedVariable->GetName() << " = " << rVal << std::endl;
689  HypoTestResult* prevResult = fResults->GetResult(fResults->ArraySize()-1);
690  if (prevResult && prevResult->GetNullDistribution() && prevResult->GetAltDistribution()) {
691  prevResult->Append(result);
692  delete result; // we can delete the result
693  }
694  else {
695  // if it was empty we re-use it
696  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunOnePoint - replace previous empty result\n";
697  fResults->fYObjects.Remove( prevResult);
698  fResults->fYObjects.Add(result);
699  }
700 
701  } else {
702 
703  // fill the results in the HypoTestInverterResult array
704  fResults->fXValues.push_back(rVal);
705  fResults->fYObjects.Add(result);
706 
707  }
708 
709  // std::cout << "computed value for poi " << rVal << " : " << fResults->GetYValue(fResults->ArraySize()-1)
710  // << " +/- " << fResults->GetYError(fResults->ArraySize()-1) << endl;
711 
712  fScannedVariable->setVal(oldValue);
713 
714  return true;
715 }
716 
717 
718 
719 bool HypoTestInverter::RunLimit(double &limit, double &limitErr, double absAccuracy, double relAccuracy, const double*hint) const {
720  // run an automatic scan until the desired accurancy is reached
721  // Start by default from the full interval (min,max) of the POI and then via bisection find the line crossing
722  // the target line
723  // Optionally an hint can be provided and the scan will be done closer to that value
724  // If by bisection the desired accuracy will not be reached a fit to the points is performed
725 
726 
727  // routine from G. Petrucciani (from HiggsCombination CMS package)
728 // bool HybridNew::runLimit(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) {
729 
731  //w->loadSnapshot("clean");
732 
733  if ((hint != 0) && (*hint > r->getMin())) {
734  r->setMax(std::min<double>(3.0 * (*hint), r->getMax()));
735  r->setMin(std::max<double>(0.3 * (*hint), r->getMin()));
736  oocoutI((TObject*)0,InputArguments) << "HypoTestInverter::RunLimit - Use hint value " << *hint
737  << " search in interval " << r->getMin() << " , " << r->getMax() << std::endl;
738  }
739 
740  // if not specified use the default values for rel and absolute accuracy
741  if (absAccuracy <= 0) absAccuracy = fgAbsAccuracy;
742  if (relAccuracy <= 0) relAccuracy = fgRelAccuracy;
743 
744  typedef std::pair<double,double> CLs_t;
745  double clsTarget = fSize;
746  CLs_t clsMin(1,0), clsMax(0,0), clsMid(0,0);
747  double rMin = r->getMin(), rMax = r->getMax();
748  limit = 0.5*(rMax + rMin);
749  limitErr = 0.5*(rMax - rMin);
750  bool done = false;
751 
752  TF1 expoFit("expoFit","[0]*exp([1]*(x-[2]))", rMin, rMax);
753 
754  // if (readHybridResults_) {
755  // if (verbose > 0) std::cout << "Search for upper limit using pre-computed grid of p-values" << std::endl;
756 
757  // readAllToysFromFile();
758  // double minDist=1e3;
759  // for (int i = 0, n = limitPlot_->GetN(); i < n; ++i) {
760  // double x = limitPlot_->GetX()[i], y = limitPlot_->GetY()[i], ey = limitPlot_->GetErrorY(i);
761  // if (verbose > 0) std::cout << " r " << x << (CLs_ ? ", CLs = " : ", CLsplusb = ") << y << " +/- " << ey << std::endl;
762  // if (y-3*ey >= clsTarget) { rMin = x; clsMin = CLs_t(y,ey); }
763  // if (y+3*ey <= clsTarget) { rMax = x; clsMax = CLs_t(y,ey); }
764  // if (fabs(y-clsTarget) < minDist) { limit = x; minDist = fabs(y-clsTarget); }
765  // }
766  // if (verbose > 0) std::cout << " after scan x ~ " << limit << ", bounds [ " << rMin << ", " << rMax << "]" << std::endl;
767  // limitErr = std::max(limit-rMin, rMax-limit);
768  // expoFit.SetRange(rMin,rMax);
769 
770  // if (limitErr < std::max(rAbsAccuracy_, rRelAccuracy_ * limit)) {
771  // if (verbose > 1) std::cout << " reached accuracy " << limitErr << " below " << std::max(rAbsAccuracy_, rRelAccuracy_ * limit) << std::endl;
772  // done = true;
773  // }
774  // } else {
775 
776  fLimitPlot.reset(new TGraphErrors());
777 
778  if (fVerbose > 0) std::cout << "Search for upper limit to the limit" << std::endl;
779  for (int tries = 0; tries < 6; ++tries) {
780  //clsMax = eval(w, mc_s, mc_b, data, rMax);
781  if (! RunOnePoint(rMax) ) {
782  oocoutE((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Hypotest failed" << std::endl;
783  return false;
784  }
785  clsMax = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
786  if (clsMax.first == 0 || clsMax.first + 3 * fabs(clsMax.second) < clsTarget ) break;
787  rMax += rMax;
788  if (tries == 5) {
789  oocoutE((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Cannot set higher limit: at " << r->GetName()
790  << " = " << rMax << " still get "
791  << (fUseCLs ? "CLs" : "CLsplusb") << " = " << clsMax.first << std::endl;
792  return false;
793  }
794  }
795  if (fVerbose > 0) {
796  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Search for lower limit to the limit" << std::endl;
797  }
798  //clsMin = (fUseCLs && rMin == 0 ? CLs_t(1,0) : eval(w, mc_s, mc_b, data, rMin));
799  if ( fUseCLs && rMin == 0 ) {
800  clsMin = CLs_t(1,0);
801  }
802  else {
803  if (! RunOnePoint(rMin) ) return false;
804  clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
805  }
806  if (clsMin.first != 1 && clsMin.first - 3 * fabs(clsMin.second) < clsTarget) {
807  if (fUseCLs) {
808  rMin = 0;
809  clsMin = CLs_t(1,0); // this is always true for CLs
810  } else {
811  rMin = -rMax / 4;
812  for (int tries = 0; tries < 6; ++tries) {
813  if (! RunOnePoint(rMax) ) return false;
814  clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
815  if (clsMin.first == 1 || clsMin.first - 3 * fabs(clsMin.second) > clsTarget) break;
816  rMin += rMin;
817  if (tries == 5) {
818  oocoutE((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Cannot set lower limit: at " << r->GetName()
819  << " = " << rMin << " still get " << (fUseCLs ? "CLs" : "CLsplusb")
820  << " = " << clsMin.first << std::endl;
821  return false;
822  }
823  }
824  }
825  }
826 
827  if (fVerbose > 0)
828  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Now doing proper bracketing & bisection" << std::endl;
829  do {
830 
831  // break loop in case max toys is reached
832  if (fMaxToys > 0 && fTotalToysRun > fMaxToys ) {
833  oocoutW((TObject*)0,Eval) << "HypoTestInverter::RunLimit - maximum number of toys reached " << std::endl;
834  done = false; break;
835  }
836 
837 
838  // determine point by bisection or interpolation
839  limit = 0.5*(rMin+rMax); limitErr = 0.5*(rMax-rMin);
840  if (fgAlgo == "logSecant" && clsMax.first != 0) {
841  double logMin = log(clsMin.first), logMax = log(clsMax.first), logTarget = log(clsTarget);
842  limit = rMin + (rMax-rMin) * (logTarget - logMin)/(logMax - logMin);
843  if (clsMax.second != 0 && clsMin.second != 0) {
844  limitErr = hypot((logTarget-logMax) * (clsMin.second/clsMin.first), (logTarget-logMin) * (clsMax.second/clsMax.first));
845  limitErr *= (rMax-rMin)/((logMax-logMin)*(logMax-logMin));
846  }
847  }
848  r->setError(limitErr);
849 
850  // exit if reached accuracy on r
851  if (limitErr < std::max(absAccuracy, relAccuracy * limit)) {
852  if (fVerbose > 1)
853  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunLimit - reached accuracy " << limitErr
854  << " below " << std::max(absAccuracy, relAccuracy * limit) << std::endl;
855  done = true; break;
856  }
857 
858  // evaluate point
859 
860  //clsMid = eval(w, mc_s, mc_b, data, limit, true, clsTarget);
861  if (! RunOnePoint(limit, true, clsTarget) ) return false;
862  clsMid = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
863 
864  if (clsMid.second == -1) {
865  std::cerr << "Hypotest failed" << std::endl;
866  return false;
867  }
868 
869  // if sufficiently far away, drop one of the points
870  if (fabs(clsMid.first-clsTarget) >= 2*clsMid.second) {
871  if ((clsMid.first>clsTarget) == (clsMax.first>clsTarget)) {
872  rMax = limit; clsMax = clsMid;
873  } else {
874  rMin = limit; clsMin = clsMid;
875  }
876  } else {
877  if (fVerbose > 0) std::cout << "Trying to move the interval edges closer" << std::endl;
878  double rMinBound = rMin, rMaxBound = rMax;
879  // try to reduce the size of the interval
880  while (clsMin.second == 0 || fabs(rMin-limit) > std::max(absAccuracy, relAccuracy * limit)) {
881  rMin = 0.5*(rMin+limit);
882  if (!RunOnePoint(rMin,true, clsTarget) ) return false;
883  clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
884  //clsMin = eval(w, mc_s, mc_b, data, rMin, true, clsTarget);
885  if (fabs(clsMin.first-clsTarget) <= 2*clsMin.second) break;
886  rMinBound = rMin;
887  }
888  while (clsMax.second == 0 || fabs(rMax-limit) > std::max(absAccuracy, relAccuracy * limit)) {
889  rMax = 0.5*(rMax+limit);
890 // clsMax = eval(w, mc_s, mc_b, data, rMax, true, clsTarget);
891  if (!RunOnePoint(rMax,true,clsTarget) ) return false;
892  clsMax = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
893  if (fabs(clsMax.first-clsTarget) <= 2*clsMax.second) break;
894  rMaxBound = rMax;
895  }
896  expoFit.SetRange(rMinBound,rMaxBound);
897  break;
898  }
899  } while (true);
900 
901  if (!done) { // didn't reach accuracy with scan, now do fit
902  if (fVerbose) {
903  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Before fit --- \n";
904  std::cout << "Limit: " << r->GetName() << " < " << limit << " +/- " << limitErr << " [" << rMin << ", " << rMax << "]\n";
905  }
906 
907  expoFit.FixParameter(0,clsTarget);
908  expoFit.SetParameter(1,log(clsMax.first/clsMin.first)/(rMax-rMin));
909  expoFit.SetParameter(2,limit);
910  double rMinBound, rMaxBound; expoFit.GetRange(rMinBound, rMaxBound);
911  limitErr = std::max(fabs(rMinBound-limit), fabs(rMaxBound-limit));
912  int npoints = 0;
913 
914  HypoTestInverterPlot plot("plot","plot",fResults);
915  fLimitPlot.reset(plot.MakePlot() );
916 
917 
918  for (int j = 0; j < fLimitPlot->GetN(); ++j) {
919  if (fLimitPlot->GetX()[j] >= rMinBound && fLimitPlot->GetX()[j] <= rMaxBound) npoints++;
920  }
921  for (int i = 0, imax = /*(readHybridResults_ ? 0 : */ 8; i <= imax; ++i, ++npoints) {
922  fLimitPlot->Sort();
923  fLimitPlot->Fit(&expoFit,(fVerbose <= 1 ? "QNR EX0" : "NR EXO"));
924  if (fVerbose) {
925  oocoutI((TObject*)0,Eval) << "Fit to " << npoints << " points: " << expoFit.GetParameter(2) << " +/- " << expoFit.GetParError(2) << std::endl;
926  }
927  if ((rMin < expoFit.GetParameter(2)) && (expoFit.GetParameter(2) < rMax) && (expoFit.GetParError(2) < 0.5*(rMaxBound-rMinBound))) {
928  // sanity check fit result
929  limit = expoFit.GetParameter(2);
930  limitErr = expoFit.GetParError(2);
931  if (limitErr < std::max(absAccuracy, relAccuracy * limit)) break;
932  }
933  // add one point in the interval.
934  double rTry = RooRandom::uniform()*(rMaxBound-rMinBound)+rMinBound;
935  if (i != imax) {
936  if (!RunOnePoint(rTry,true,clsTarget) ) return false;
937  //eval(w, mc_s, mc_b, data, rTry, true, clsTarget);
938  }
939 
940  }
941  }
942 
943 //if (!plot_.empty() && fLimitPlot.get()) {
944  if (fLimitPlot.get() && fLimitPlot->GetN() > 0) {
945  //new TCanvas("c1","c1");
946  fLimitPlot->Sort();
947  fLimitPlot->SetLineWidth(2);
948  double xmin = r->getMin(), xmax = r->getMax();
949  for (int j = 0; j < fLimitPlot->GetN(); ++j) {
950  if (fLimitPlot->GetY()[j] > 1.4*clsTarget || fLimitPlot->GetY()[j] < 0.6*clsTarget) continue;
951  xmin = std::min(fLimitPlot->GetX()[j], xmin);
952  xmax = std::max(fLimitPlot->GetX()[j], xmax);
953  }
954  fLimitPlot->GetXaxis()->SetRangeUser(xmin,xmax);
955  fLimitPlot->GetYaxis()->SetRangeUser(0.5*clsTarget, 1.5*clsTarget);
956  fLimitPlot->Draw("AP");
957  expoFit.Draw("SAME");
958  TLine line(fLimitPlot->GetX()[0], clsTarget, fLimitPlot->GetX()[fLimitPlot->GetN()-1], clsTarget);
959  line.SetLineColor(kRed); line.SetLineWidth(2); line.Draw();
960  line.DrawLine(limit, 0, limit, fLimitPlot->GetY()[0]);
961  line.SetLineWidth(1); line.SetLineStyle(2);
962  line.DrawLine(limit-limitErr, 0, limit-limitErr, fLimitPlot->GetY()[0]);
963  line.DrawLine(limit+limitErr, 0, limit+limitErr, fLimitPlot->GetY()[0]);
964  //c1->Print(plot_.c_str());
965  }
966 
967  oocoutI((TObject*)0,Eval) << "HypoTestInverter::RunLimit - Result: \n"
968  << "\tLimit: " << r->GetName() << " < " << limit << " +/- " << limitErr << " @ " << (1-fSize) * 100 << "% CL\n";
969  if (fVerbose > 1) oocoutI((TObject*)0,Eval) << "Total toys: " << fTotalToysRun << std::endl;
970 
971  // set value in results
972  fResults->fUpperLimit = limit;
973  fResults->fUpperLimitError = limitErr;
974  fResults->fFittedUpperLimit = true;
975  // lower limit are always min of p value
978  fResults->fFittedLowerLimit = true;
979 
980  return true;
981 }
982 
984  // get the distribution of lower limit
985  // if rebuild = false (default) it will re-use the results of the scan done
986  // for obtained the observed limit and no extra toys will be generated
987  // if rebuild a new set of B toys will be done and the procedure will be repeted
988  // for each toy
989  if (!rebuild) {
990  if (!fResults) {
991  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::GetLowerLimitDistribution(false) - result not existing\n";
992  return 0;
993  }
995  }
996 
997  TList * clsDist = 0;
998  TList * clsbDist = 0;
999  if (fUseCLs) clsDist = &fResults->fExpPValues;
1000  else clsbDist = &fResults->fExpPValues;
1001 
1002  return RebuildDistributions(false, nToys,clsDist, clsbDist);
1003 
1004 }
1005 
1007  // get the distribution of lower limit
1008  // if rebuild = false (default) it will re-use the results of the scan done
1009  // for obtained the observed limit and no extra toys will be generated
1010  // if rebuild a new set of B toys will be done and the procedure will be repeted
1011  // for each toy
1012  // The nuisance parameter value used for rebuild is the current one in the model
1013  // so it is user responsability to set to the desired value (nomi
1014  if (!rebuild) {
1015  if (!fResults) {
1016  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::GetUpperLimitDistribution(false) - result not existing\n";
1017  return 0;
1018  }
1020  }
1021 
1022  TList * clsDist = 0;
1023  TList * clsbDist = 0;
1024  if (fUseCLs) clsDist = &fResults->fExpPValues;
1025  else clsbDist = &fResults->fExpPValues;
1026 
1027  return RebuildDistributions(true, nToys,clsDist, clsbDist);
1028 }
1029 
1031  if (fCalculator0) fCalculator0->SetData(data);
1032 }
1033 
1034 SamplingDistribution * HypoTestInverter::RebuildDistributions(bool isUpper, int nToys, TList * clsDist, TList * clsbDist, TList * clbDist, const char *outputfile) {
1035  // rebuild the sampling distributions by
1036  // generating some toys and find for each of theam a new upper limit
1037  // Return the upper limit distribution and optionally also the pValue distributions for Cls, Clsb and Clbxs
1038  // as a TList for each scanned point
1039  // The method uses the present parameter value. It is user responsability to give the current parameters to rebuild the distributions
1040  // It returns a upper or lower limit distribution depending on the isUpper flag, however it computes also the lower limit distribution and it is saved in the
1041  // output file as an histogram
1042 
1043  if (!fScannedVariable || !fCalculator0) return 0;
1044  // get first background snapshot
1045  const ModelConfig * bModel = fCalculator0->GetAlternateModel();
1046  const ModelConfig * sbModel = fCalculator0->GetNullModel();
1047  if (!bModel || ! sbModel) return 0;
1048  RooArgSet paramPoint;
1049  if (!sbModel->GetParametersOfInterest()) return 0;
1050  paramPoint.add(*sbModel->GetParametersOfInterest());
1051 
1052  const RooArgSet * poibkg = bModel->GetSnapshot();
1053  if (!poibkg) {
1054  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter::RebuildDistribution - bakground snapshot not existing"
1055  << " assume is for POI = 0" << std::endl;
1057  paramPoint = RooArgSet(*fScannedVariable);
1058  }
1059  else
1060  paramPoint = *poibkg;
1061  // generate data at bkg parameter point
1062 
1063  ToyMCSampler * toymcSampler = dynamic_cast<ToyMCSampler *>(fCalculator0->GetTestStatSampler() );
1064  if (!toymcSampler) {
1065  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter::RebuildDistribution - no toy MC sampler existing" << std::endl;
1066  return 0;
1067  }
1068  // set up test stat sampler in case of asymptotic calculator
1069  if (dynamic_cast<RooStats::AsymptoticCalculator*>(fCalculator0) ) {
1070  toymcSampler->SetObservables(*sbModel->GetObservables() );
1071  toymcSampler->SetParametersForTestStat(*sbModel->GetParametersOfInterest());
1072  toymcSampler->SetPdf(*sbModel->GetPdf());
1073  toymcSampler->SetNuisanceParameters(*sbModel->GetNuisanceParameters());
1074  if (sbModel->GetGlobalObservables() ) toymcSampler->SetGlobalObservables(*sbModel->GetGlobalObservables() );
1075  // set number of events
1076  if (!sbModel->GetPdf()->canBeExtended())
1077  toymcSampler->SetNEventsPerToy(1);
1078  }
1079 
1080  // loop on data to generate
1081  int nPoints = fNBins;
1082 
1083  bool storePValues = clsDist || clsbDist || clbDist;
1084  if (fNBins <=0 && storePValues) {
1085  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter::RebuildDistribution - cannot return p values distribution with the auto scan" << std::endl;
1086  storePValues = false;
1087  nPoints = 0;
1088  }
1089 
1090  if (storePValues) {
1091  if (fResults) nPoints = fResults->ArraySize();
1092  if (nPoints <=0) {
1093  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - result is not existing and number of point to scan is not set"
1094  << std::endl;
1095  return 0;
1096  }
1097  }
1098 
1099  if (nToys <= 0) nToys = 100; // default value
1100 
1101  std::vector<std::vector<double> > CLs_values(nPoints);
1102  std::vector<std::vector<double> > CLsb_values(nPoints);
1103  std::vector<std::vector<double> > CLb_values(nPoints);
1104 
1105  if (storePValues) {
1106  for (int i = 0; i < nPoints; ++i) {
1107  CLs_values[i].reserve(nToys);
1108  CLb_values[i].reserve(nToys);
1109  CLsb_values[i].reserve(nToys);
1110  }
1111  }
1112 
1113  std::vector<double> limit_values; limit_values.reserve(nToys);
1114 
1115  oocoutI((TObject*)0,InputArguments) << "HypoTestInverter - rebuilding the p value distributions by generating ntoys = "
1116  << nToys << std::endl;
1117 
1118 
1119  oocoutI((TObject*)0,InputArguments) << "Rebuilding using parameter of interest point: ";
1121  if (sbModel->GetNuisanceParameters() ) {
1122  oocoutI((TObject*)0,InputArguments) << "And using nuisance parameters: ";
1124  }
1125  // save all parameters to restore them later
1126  assert(bModel->GetPdf() );
1127  assert(bModel->GetObservables() );
1128  RooArgSet * allParams = bModel->GetPdf()->getParameters( *bModel->GetObservables() );
1129  RooArgSet saveParams;
1130  allParams->snapshot(saveParams);
1131 
1132  TFile * fileOut = TFile::Open(outputfile,"RECREATE");
1133  if (!fileOut) {
1134  oocoutE((TObject*)0,InputArguments) << "HypoTestInverter - RebuildDistributions - Error opening file " << outputfile
1135  << " - the resulting limits will not be stored" << std::endl;
1136  }
1137  // create temporary histograms to store the limit result
1138  TH1D * hL = new TH1D("lowerLimitDist","Rebuilt lower limit distribution",100,1.,0.);
1139  TH1D * hU = new TH1D("upperLimitDist","Rebuilt upper limit distribution",100,1.,0.);
1140  TH1D * hN = new TH1D("nObs","Observed events",100,1.,0.);
1141  hL->SetBuffer(2*nToys);
1142  hU->SetBuffer(2*nToys);
1143  std::vector<TH1*> hCLb;
1144  std::vector<TH1*> hCLsb;
1145  std::vector<TH1*> hCLs;
1146  if (storePValues) {
1147  for (int i = 0; i < nPoints; ++i) {
1148  hCLb.push_back(new TH1D(TString::Format("CLbDist_bin%d",i),"CLb distribution",100,1.,0.));
1149  hCLs.push_back(new TH1D(TString::Format("ClsDist_bin%d",i),"CLs distribution",100,1.,0.));
1150  hCLsb.push_back(new TH1D(TString::Format("CLsbDist_bin%d",i),"CLs+b distribution",100,1.,0.));
1151  }
1152  }
1153 
1154 
1155  // loop now on the toys
1156  for (int itoy = 0; itoy < nToys; ++itoy) {
1157 
1158  oocoutP((TObject*)0,Eval) << "\nHypoTestInverter - RebuildDistributions - running toy # " << itoy << " / "
1159  << nToys << std::endl;
1160 
1161 
1162  printf("\n\nshnapshot of s+b model \n");
1163  sbModel->GetSnapshot()->Print("v");
1164 
1165  // reset parameters to initial values to be sure in case they are not reset
1166  if (itoy> 0) *allParams = saveParams;
1167 
1168  // need to set th epdf to clear the cache in ToyMCSampler
1169  // pdf we must use is background pdf
1170  toymcSampler->SetPdf(*bModel->GetPdf() );
1171 
1172 
1173  RooAbsData * bkgdata = toymcSampler->GenerateToyData(paramPoint);
1174 
1175  double nObs = bkgdata->sumEntries();
1176  // for debugging in case of number counting models
1177  if (bkgdata->numEntries() ==1 && !bModel->GetPdf()->canBeExtended()) {
1178  oocoutP((TObject*)0,Generation) << "Generate observables are : ";
1179  RooArgList genObs(*bkgdata->get(0));
1181  nObs = 0;
1182  for (int i = 0; i < genObs.getSize(); ++i) {
1183  RooRealVar * x = dynamic_cast<RooRealVar*>(&genObs[i]);
1184  if (x) nObs += x->getVal();
1185  }
1186  }
1187  hN->Fill(nObs);
1188 
1189  // by copying I will have the same min/max as previous ones
1190  HypoTestInverter inverter = *this;
1191  inverter.SetData(*bkgdata);
1192 
1193  // print global observables
1194  auto gobs = bModel->GetPdf()->getVariables()->selectCommon(* sbModel->GetGlobalObservables() );
1195  gobs->Print("v");
1196 
1197  HypoTestInverterResult * r = inverter.GetInterval();
1198 
1199  if (r == 0) continue;
1200 
1201  double value = (isUpper) ? r->UpperLimit() : r->LowerLimit();
1202  limit_values.push_back( value );
1203  hU->Fill(r->UpperLimit() );
1204  hL->Fill(r->LowerLimit() );
1205 
1206 
1207  std::cout << "The computed upper limit for toy #" << itoy << " is " << value << std::endl;
1208 
1209  // write every 10 toys
1210  if (itoy%10 == 0 || itoy == nToys-1) {
1211  hU->Write("",TObject::kOverwrite);
1212  hL->Write("",TObject::kOverwrite);
1213  hN->Write("",TObject::kOverwrite);
1214  }
1215 
1216  if (!storePValues) continue;
1217 
1218  if (nPoints < r->ArraySize()) {
1219  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter: skip extra points" << std::endl;
1220  }
1221  else if (nPoints > r->ArraySize()) {
1222  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter: missing some points" << std::endl;
1223  }
1224 
1225 
1226  for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1227  HypoTestResult * hr = r->GetResult(ipoint);
1228  if (hr) {
1229  CLs_values[ipoint].push_back( hr->CLs() );
1230  CLsb_values[ipoint].push_back( hr->CLsplusb() );
1231  CLb_values[ipoint].push_back( hr->CLb() );
1232  hCLs[ipoint]->Fill( hr->CLs() );
1233  hCLb[ipoint]->Fill( hr->CLb() );
1234  hCLsb[ipoint]->Fill( hr->CLsplusb() );
1235  }
1236  else {
1237  oocoutW((TObject*)0,InputArguments) << "HypoTestInverter: missing result for point: x = "
1238  << fResults->GetXValue(ipoint) << std::endl;
1239  }
1240  }
1241  // write every 10 toys
1242  if (itoy%10 == 0 || itoy == nToys-1) {
1243  for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1244  hCLs[ipoint]->Write("",TObject::kOverwrite);
1245  hCLb[ipoint]->Write("",TObject::kOverwrite);
1246  hCLsb[ipoint]->Write("",TObject::kOverwrite);
1247  }
1248  }
1249 
1250 
1251  delete r;
1252  delete bkgdata;
1253  }
1254 
1255 
1256  if (storePValues) {
1257  if (clsDist) clsDist->SetOwner(true);
1258  if (clbDist) clbDist->SetOwner(true);
1259  if (clsbDist) clsbDist->SetOwner(true);
1260 
1261  oocoutI((TObject*)0,InputArguments) << "HypoTestInverter: storing rebuilt p values " << std::endl;
1262 
1263  for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1264  if (clsDist) {
1265  TString name = TString::Format("CLs_distrib_%d",ipoint);
1266  clsDist->Add( new SamplingDistribution(name,name,CLs_values[ipoint] ) );
1267  }
1268  if (clbDist) {
1269  TString name = TString::Format("CLb_distrib_%d",ipoint);
1270  clbDist->Add( new SamplingDistribution(name,name,CLb_values[ipoint] ) );
1271  }
1272  if (clsbDist) {
1273  TString name = TString::Format("CLsb_distrib_%d",ipoint);
1274  clsbDist->Add( new SamplingDistribution(name,name,CLsb_values[ipoint] ) );
1275  }
1276  }
1277  }
1278 
1279  if (fileOut) {
1280  fileOut->Close();
1281  }
1282  else {
1283  // delete all the histograms
1284  delete hL;
1285  delete hU;
1286  for (int i = 0; i < nPoints && storePValues; ++i) {
1287  delete hCLs[i];
1288  delete hCLb[i];
1289  delete hCLsb[i];
1290  }
1291  }
1292 
1293  const char * disName = (isUpper) ? "upperLimit_dist" : "lowerLimit_dist";
1294  return new SamplingDistribution(disName, disName, limit_values);
1295 }
virtual Double_t sumEntries() const =0
bool rebuild
virtual HypoTestInverterResult * GetInterval() const
Main interface to get a ConfInterval, pure virtual.
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:823
void SetBackgroundAsAlt(Bool_t l=kTRUE)
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
float xmin
Definition: THbookFile.cxx:93
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:52
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:244
Double_t CLbError() const
The error on the "confidence level" of the null hypothesis.
Small helper to implement the function to create,access and destroy iterators.
HypoTestInverter class for performing an hypothesis test inversion by scanning the hypothesis test re...
IntervalCalculator is an interface class for a tools which produce RooStats ConfIntervals.
virtual void SetParametersForTestStat(const RooArgSet &nullpoi)
Definition: ToyMCSampler.h:185
TLine * line
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
virtual TestStatistic * GetTestStatistic() const =0
virtual void SetGlobalObservables(const RooArgSet &o)
Definition: ToyMCSampler.h:205
int ArraySize() const
number of entries in the results array
virtual Double_t CLs() const
CLs is simply CLs+b/CLb (not a method, but a quantity)
Definition: Rtypes.h:61
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return NULL if not existing)
Definition: ModelConfig.h:265
bool RunOnePoint(double thisX, bool adaptive=false, double clTarget=-1) const
virtual Double_t ConfidenceLevel() const
Get the Confidence level for the test.
#define oocoutF(o, a)
Definition: RooMsgService.h:49
std::unique_ptr< TGraphErrors > fLimitPlot
#define assert(cond)
Definition: unittest.h:542
TestStatSampler is an interface class for a tools which produce RooStats SamplingDistributions.
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TH1 * h
Definition: legend2.C:5
TestStatistic * GetTestStatistic() const
#define oocoutI(o, a)
Definition: RooMsgService.h:45
HypoTestResult is a base class for results from hypothesis tests.
SamplingDistribution * GetUpperLimitDistribution(bool rebuild=false, int nToys=100)
virtual void SetNEventsPerToy(const Int_t nevents)
Definition: ToyMCSampler.h:177
virtual RooAbsData * GenerateToyData(RooArgSet &paramPoint, RooAbsPdf &pdf) const
Definition: ToyMCSampler.h:140
Basic string class.
Definition: TString.h:137
virtual Double_t getMin(const char *name=0) const
static void CheckInputModels(const HypoTestCalculatorGeneric &hc, const RooRealVar &scanVar)
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
bool Bool_t
Definition: RtypesCore.h:59
bool RunFixedScan(int nBins, double xMin, double xMax, bool scanLog=false) const
SamplingDistribution * GetAltDistribution(void) const
SamplingDistribution * GetUpperLimitDistribution() const
get expected upper limit distributions implemented using interpolation
Int_t GetSize() const
size of samples
virtual Double_t CLsplusb() const
Convert AlternatePValue into a "confidence level".
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
Definition: ModelConfig.h:250
HypoTestInverterResult * fResults
void setMax(const char *name, Double_t value)
Set maximum of name range to given value.
Definition: RooRealVar.cxx:416
virtual void SetBuffer(Int_t buffersize, Option_t *option="")
set the maximum number of entries to be kept in the buffer
Definition: TH1.cxx:7837
HypoTestInverter & operator=(const HypoTestInverter &rhs)
virtual void SetTestStatistic(TestStatistic *testStatistic)=0
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3851
RooAbsArg * first() const
void plot(TString fname="data.root", TString var0="var0", TString var1="var1")
Definition: createData.C:17
Common base class for the Hypothesis Test Calculators.
#define oocoutP(o, a)
Definition: RooMsgService.h:46
Double_t x[n]
Definition: legend1.C:17
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:2321
virtual HypoTestResult * GetHypoTest() const
inherited methods from HypoTestCalculator interface
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2082
TestStatSampler * GetTestStatSampler(void) const
#define oocoutE(o, a)
Definition: RooMsgService.h:48
static Vc_ALWAYS_INLINE Vector< T > abs(const Vector< T > &x)
Definition: vector.h:450
void UseCLs(bool on=true)
flag to switch between using CLsb (default) or CLs as confidence level
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual Double_t NullPValue() const
Return p-value for null hypothesis.
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:63
Double_t LowerLimit()
lower and upper bound of the confidence interval (to get upper/lower limits, multiply the size( = 1-c...
virtual void SetData(RooAbsData &data)
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:560
HypoTestResult * GetResult(int index) const
return a pointer to the i^th result object
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual const RooArgSet * get() const
Definition: RooAbsData.h:77
A doubly linked list.
Definition: TList.h:47
const ModelConfig * GetAlternateModel(void) const
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Definition: TMath.h:192
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Definition: TMath.h:196
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:203
HypoTestResult * Eval(HypoTestCalculatorGeneric &hc, bool adaptive, double clsTarget) const
virtual Double_t CLb() const
Convert NullPValue into a "confidence level".
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
const ModelConfig * GetNullModel(void) const
ROOT::R::TRInterface & r
Definition: Object.C:4
TGraphErrors * MakePlot(Option_t *opt="")
return a TGraphErrors with the obtained observed p-values resultinf from the scan By default (Option ...
RooAbsArg * find(const char *name) const
Find object with given name in list.
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
void setMin(const char *name, Double_t value)
Set minimum of name range to given value.
Definition: RooRealVar.cxx:386
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:675
int fTotalToysRun
plot of limits
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:293
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
Double_t E()
Definition: TMath.h:54
TList fExpPValues
list of HypoTestResult for each point
RooRealVar * fScannedVariable
pointer to the generic hypotest calculator used
HypoTestCalculatorGeneric * fCalculator0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
tuple pl
Definition: first.py:10
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
float xmax
Definition: THbookFile.cxx:93
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
ToyMCSampler is an implementation of the TestStatSampler interface.
Definition: ToyMCSampler.h:99
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
Bool_t GetPValueIsRightTail(void) const
void SetTestStatisticData(const Double_t tsd)
This class simply holds a sampling distribution of some test statistic.
virtual void Append(const HypoTestResult *other)
add values from another HypoTestResult
This class implements the Hypothesis test calculation using an hybrid (frequentist/bayesian) procedur...
ClassImp(RooStats::HypoTestInverter) using namespace RooStats
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
Definition: ModelConfig.h:259
double GetXValue(int index) const
function to return the value of the parameter of interest for the i^th entry in the results ...
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:84
HypoTestInverterResult class: holds the array of hypothesis test results and compute a confidence int...
virtual void SetObservables(const RooArgSet &o)
Definition: ToyMCSampler.h:203
Class to plot an HypoTestInverterResult, result of the HypoTestInverter calculator.
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:216
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
SamplingDistribution * GetNullDistribution(void) const
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return NULL if not existing) ...
virtual void SetPdf(RooAbsPdf &pdf)
Definition: ToyMCSampler.h:190
bool SetTestStatistic(TestStatistic &stat)
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
#define name(a, b)
Definition: linkTestLib0.cpp:5
virtual Double_t AlternatePValue() const
Return p-value for alternate hypothesis.
#define oocoutW(o, a)
Definition: RooMsgService.h:47
Mother of all ROOT objects.
Definition: TObject.h:58
virtual Double_t getMax(const char *name=0) const
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual void Add(TObject *obj)
Definition: TList.h:81
double fLowerLimitError
interpolatation option (linear or spline)
Int_t IsNaN(Double_t x)
Definition: TMath.h:617
std::unique_ptr< HypoTestCalculatorGeneric > fHC
Hypothesis Test Calculator using a full frequentist procedure for sampling the test statistic distrib...
SamplingDistribution * GetLowerLimitDistribution() const
get expected lower limit distributions implemented using interpolation The size for the sampling dist...
Double_t CLsplusbError() const
The error on the "confidence level" of the alternative hypothesis.
static RooRealVar * GetVariableToScan(const HypoTestCalculatorGeneric &hc)
bool RunLimit(double &limit, double &limitErr, double absTol=0, double relTol=0, const double *hint=0) const
Double_t GetTestStatisticData(void) const
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:247
double result[121]
Py_ssize_t ArraySize(const std::string &name)
Extract size from an array type, if available.
Definition: Utility.cxx:674
SamplingDistribution * RebuildDistributions(bool isUpper=true, int nToys=100, TList *clsDist=0, TList *clsbDist=0, TList *clbDist=0, const char *outputfile="HypoTestInverterRebuiltDist.root")
double exp(double)
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
virtual void SetNuisanceParameters(const RooArgSet &np)
Definition: ToyMCSampler.h:201
virtual void SetData(RooAbsData &)
Set the DataSet ( add to the the workspace if not already there ?)
float value
Definition: math.cpp:443
static void CloseProof(Option_t *option="s")
close all proof connections
Definition: ProofConfig.h:95
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
Definition: TestStatistic.h:33
SamplingDistribution * GetLowerLimitDistribution(bool rebuild=false, int nToys=100)
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio...
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval) ...
double log(double)
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add element to non-owning set.
Definition: RooArgSet.cxx:448
Double_t CLsError() const
The error on the ratio CLs+b/CLb.
void setError(Double_t value)
Definition: RooRealVar.h:56
std::vector< double > fXValues
max sigma value used to scan asymptotic expected p values