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