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(nullptr,InputArguments) << "HypoTestInverter - model are not existing" << std::endl;
132 assert(modelSB && modelB);
133
134 oocoutI(nullptr,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(nullptr,InputArguments) << "HypoTestInverter - B model has no pdf or observables defined" << std::endl;
145 return;
146 }
147 std::unique_ptr<RooArgSet> bParams{bPdf->getParameters(*bObs)};
148 if (!bParams) {
149 oocoutE(nullptr,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(nullptr,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}
162
163////////////////////////////////////////////////////////////////////////////////
164/// default constructor (doesn't do anything)
165
167 fTotalToysRun(0),
168 fMaxToys(0),
169 fCalculator0(0),
170 fScannedVariable(0),
171 fResults(0),
172 fUseCLs(false),
173 fScanLog(false),
174 fSize(0),
175 fVerbose(0),
176 fCalcType(kUndefined),
177 fNBins(0), fXmin(1), fXmax(1),
178 fNumErr(0)
179{
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Constructor from a HypoTestCalculatorGeneric
184/// The HypoTest calculator must be a FrequentistCalculator or HybridCalculator type
185/// Other type of calculators are not supported.
186/// The calculator must be created before by using the S+B model for the null and
187/// the B model for the alt
188/// If no variable to scan are given they are assumed to be the first variable
189/// from the parameter of interests of the null model
190
192 RooRealVar* scannedVariable, double size ) :
193 fTotalToysRun(0),
194 fMaxToys(0),
195 fCalculator0(0),
196 fScannedVariable(scannedVariable),
197 fResults(0),
198 fUseCLs(false),
199 fScanLog(false),
200 fSize(size),
201 fVerbose(0),
202 fCalcType(kUndefined),
203 fNBins(0), fXmin(1), fXmax(1),
204 fNumErr(0)
205{
206
207 if (!fScannedVariable) {
209 }
210 if (!fScannedVariable)
211 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
212 else
214
215 HybridCalculator * hybCalc = dynamic_cast<HybridCalculator*>(&hc);
216 if (hybCalc) {
218 fCalculator0 = hybCalc;
219 return;
220 }
221 FrequentistCalculator * freqCalc = dynamic_cast<FrequentistCalculator*>(&hc);
222 if (freqCalc) {
224 fCalculator0 = freqCalc;
225 return;
226 }
227 AsymptoticCalculator * asymCalc = dynamic_cast<AsymptoticCalculator*>(&hc);
228 if (asymCalc) {
230 fCalculator0 = asymCalc;
231 return;
232 }
233 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Type of hypotest calculator is not supported " <<std::endl;
234 fCalculator0 = &hc;
235}
236
237////////////////////////////////////////////////////////////////////////////////
238/// Constructor from a reference to a HybridCalculator
239/// The calculator must be created before by using the S+B model for the null and
240/// the B model for the alt
241/// If no variable to scan are given they are assumed to be the first variable
242/// from the parameter of interests of the null model
243
245 RooRealVar* scannedVariable, double size ) :
246 fTotalToysRun(0),
247 fMaxToys(0),
248 fCalculator0(&hc),
249 fScannedVariable(scannedVariable),
250 fResults(0),
251 fUseCLs(false),
252 fScanLog(false),
253 fSize(size),
254 fVerbose(0),
255 fCalcType(kHybrid),
256 fNBins(0), fXmin(1), fXmax(1),
257 fNumErr(0)
258{
259
260 if (!fScannedVariable) {
262 }
263 if (!fScannedVariable)
264 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
265 else
267
268}
269
270////////////////////////////////////////////////////////////////////////////////
271/// Constructor from a reference to a FrequentistCalculator
272/// The calculator must be created before by using the S+B model for the null and
273/// the B model for the alt
274/// If no variable to scan are given they are assumed to be the first variable
275/// from the parameter of interests of the null model
276
278 RooRealVar* scannedVariable, double size ) :
279 fTotalToysRun(0),
280 fMaxToys(0),
281 fCalculator0(&hc),
282 fScannedVariable(scannedVariable),
283 fResults(0),
284 fUseCLs(false),
285 fScanLog(false),
286 fSize(size),
287 fVerbose(0),
288 fCalcType(kFrequentist),
289 fNBins(0), fXmin(1), fXmax(1),
290 fNumErr(0)
291{
292
293 if (!fScannedVariable) {
295 }
296 if (!fScannedVariable)
297 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
298 else
300}
301
302////////////////////////////////////////////////////////////////////////////////
303/// Constructor from a reference to a AsymptoticCalculator
304/// The calculator must be created before by using the S+B model for the null and
305/// the B model for the alt
306/// If no variable to scan are given they are assumed to be the first variable
307/// from the parameter of interests of the null model
308
310 RooRealVar* scannedVariable, double size ) :
311 fTotalToysRun(0),
312 fMaxToys(0),
313 fCalculator0(&hc),
314 fScannedVariable(scannedVariable),
315 fResults(0),
316 fUseCLs(false),
317 fScanLog(false),
318 fSize(size),
319 fVerbose(0),
320 fCalcType(kAsymptotic),
321 fNBins(0), fXmin(1), fXmax(1),
322 fNumErr(0)
323{
324
325 if (!fScannedVariable) {
327 }
328 if (!fScannedVariable)
329 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
330 else
332
333}
334
335////////////////////////////////////////////////////////////////////////////////
336/// Constructor from a model for B model and a model for S+B.
337/// An HypoTestCalculator (Hybrid of Frequentis) will be created using the
338/// S+B model as the null and the B model as the alternate
339/// If no variable to scan are given they are assumed to be the first variable
340/// from the parameter of interests of the null model
341
343 RooRealVar * scannedVariable, ECalculatorType type, double size) :
344 fTotalToysRun(0),
345 fMaxToys(0),
346 fCalculator0(0),
347 fScannedVariable(scannedVariable),
348 fResults(0),
349 fUseCLs(false),
350 fScanLog(false),
351 fSize(size),
352 fVerbose(0),
353 fCalcType(type),
354 fNBins(0), fXmin(1), fXmax(1),
355 fNumErr(0)
356{
357 if(fCalcType==kFrequentist) fHC = std::make_unique<FrequentistCalculator>(data, bModel, sbModel);
358 if(fCalcType==kHybrid) fHC = std::make_unique<HybridCalculator>(data, bModel, sbModel);
359 if(fCalcType==kAsymptotic) fHC = std::make_unique<AsymptoticCalculator>(data, bModel, sbModel);
360 fCalculator0 = fHC.get();
361 // get scanned variable
362 if (!fScannedVariable) {
364 }
365 if (!fScannedVariable)
366 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
367 else
369
370}
371
372////////////////////////////////////////////////////////////////////////////////
373/// copy-constructor
374/// NOTE: this class does not copy the contained result and
375/// the HypoTestCalculator, but only the pointers
376/// It requires the original HTI to be alive
377
380 fTotalToysRun(0),
381 fCalculator0(0), fScannedVariable(0), // add these for Coverity
382 fResults(0)
383{
384 (*this) = rhs;
385}
386
387////////////////////////////////////////////////////////////////////////////////
388/// assignment operator
389/// NOTE: this class does not copy the contained result and
390/// the HypoTestCalculator, but only the pointers
391/// It requires the original HTI to be alive
392
394 if (this == &rhs) return *this;
395 fTotalToysRun = 0;
396 fMaxToys = rhs.fMaxToys;
399 fUseCLs = rhs.fUseCLs;
400 fScanLog = rhs.fScanLog;
401 fSize = rhs.fSize;
402 fVerbose = rhs.fVerbose;
403 fCalcType = rhs.fCalcType;
404 fNBins = rhs.fNBins;
405 fXmin = rhs.fXmin;
406 fXmax = rhs.fXmax;
407 fNumErr = rhs.fNumErr;
408
409 return *this;
410}
411
412////////////////////////////////////////////////////////////////////////////////
413/// destructor (delete the HypoTestInverterResult)
414
416{
417 if (fResults) delete fResults;
418 fCalculator0 = 0;
419}
420
421////////////////////////////////////////////////////////////////////////////////
422/// return the test statistic which is or will be used by the class
423
425{
428 }
429 else
430 return 0;
431}
432
433////////////////////////////////////////////////////////////////////////////////
434/// set the test statistic to use
435
437{
440 return true;
441 }
442 else return false;
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// delete contained result and graph
447
449 if (fResults) delete fResults;
450 fResults = 0;
451 fLimitPlot.reset(nullptr);
452}
453
454////////////////////////////////////////////////////////////////////////////////
455/// create a new HypoTestInverterResult to hold all computed results
456
458 if (fResults == 0) {
459 TString results_name = "result_";
460 results_name += fScannedVariable->GetName();
462 TString title = "HypoTestInverter Result For ";
463 title += fScannedVariable->GetName();
464 fResults->SetTitle(title);
465 }
468 // check if one or two sided scan
469 if (fCalculator0) {
470 // if asymptotic calculator
472 if (ac)
474 else {
475 // in case of the other calculators
477 if (sampler) {
479 if (pl && pl->IsTwoSided() ) fResults->fIsTwoSided = true;
480 }
481 }
482 }
483}
484
485////////////////////////////////////////////////////////////////////////////////
486/// Run a fixed scan or the automatic scan depending on the configuration.
487/// Return if needed a copy of the result object which will be managed by the user.
488
490
491 // if having a result with at least one point return it
492 if (fResults && fResults->ArraySize() >= 1) {
493 oocoutI(nullptr,Eval) << "HypoTestInverter::GetInterval - return an already existing interval " << std::endl;
495 }
496
497 if (fNBins > 0) {
498 oocoutI(nullptr,Eval) << "HypoTestInverter::GetInterval - run a fixed scan" << std::endl;
499 bool ret = RunFixedScan(fNBins, fXmin, fXmax, fScanLog);
500 if (!ret)
501 oocoutE(nullptr,Eval) << "HypoTestInverter::GetInterval - error running a fixed scan " << std::endl;
502 }
503 else {
504 oocoutI(nullptr,Eval) << "HypoTestInverter::GetInterval - run an automatic scan" << std::endl;
505 double limit(0),err(0);
506 bool ret = RunLimit(limit,err);
507 if (!ret)
508 oocoutE(nullptr,Eval) << "HypoTestInverter::GetInterval - error running an auto scan " << std::endl;
509 }
510
512
514}
515
516////////////////////////////////////////////////////////////////////////////////
517/// Run the Hypothesis test at a previous configured point
518/// (internal function called by RunOnePoint)
519
520HypoTestResult * HypoTestInverter::Eval(HypoTestCalculatorGeneric &hc, bool adaptive, double clsTarget) const {
521 //for debug
522 // std::cout << ">>>>>>>>>>> " << std::endl;
523 // std::cout << "alternate model " << std::endl;
524 // hc.GetAlternateModel()->GetNuisanceParameters()->Print("V");
525 // hc.GetAlternateModel()->GetParametersOfInterest()->Print("V");
526 // std::cout << "Null model " << std::endl;
527 // hc.GetNullModel()->GetNuisanceParameters()->Print("V");
528 // hc.GetNullModel()->GetParametersOfInterest()->Print("V");
529 // std::cout << "<<<<<<<<<<<<<<< " << std::endl;
530
531 // run the hypothesis test
532 HypoTestResult * hcResult = hc.GetHypoTest();
533 if (hcResult == 0) {
534 oocoutE(nullptr,Eval) << "HypoTestInverter::Eval - HypoTest failed" << std::endl;
535 return hcResult;
536 }
537
538 // since the b model is the alt need to set the flag
539 hcResult->SetBackgroundAsAlt(true);
540
541
542 // bool flipPvalue = false;
543 // if (flipPValues)
544 // hcResult->SetPValueIsRightTail(!hcResult->GetPValueIsRightTail());
545
546 // adjust for some numerical error in discrete models and == is not anymore
547 if (hcResult->GetPValueIsRightTail() )
548 hcResult->SetTestStatisticData(hcResult->GetTestStatisticData()-fNumErr); // issue with < vs <= in discrete models
549 else
550 hcResult->SetTestStatisticData(hcResult->GetTestStatisticData()+fNumErr); // issue with < vs <= in discrete models
551
552 double clsMid = (fUseCLs ? hcResult->CLs() : hcResult->CLsplusb());
553 double clsMidErr = (fUseCLs ? hcResult->CLsError() : hcResult->CLsplusbError());
554
555 //if (fVerbose) std::cout << (fUseCLs ? "\tCLs = " : "\tCLsplusb = ") << clsMid << " +/- " << clsMidErr << std::endl;
556
557 if (adaptive) {
558
561
562 while (clsMidErr >= fgCLAccuracy && (clsTarget == -1 || std::abs(clsMid-clsTarget) < 3*clsMidErr) ) {
563 std::unique_ptr<HypoTestResult> more(hc.GetHypoTest());
564
565 // if (flipPValues)
566 // more->SetPValueIsRightTail(!more->GetPValueIsRightTail());
567
568 hcResult->Append(more.get());
569 clsMid = (fUseCLs ? hcResult->CLs() : hcResult->CLsplusb());
570 clsMidErr = (fUseCLs ? hcResult->CLsError() : hcResult->CLsplusbError());
571 if (fVerbose) std::cout << (fUseCLs ? "\tCLs = " : "\tCLsplusb = ") << clsMid << " +/- " << clsMidErr << std::endl;
572 }
573
574 }
575 if (fVerbose ) {
576 oocoutP(nullptr,Eval) << "P values for " << fScannedVariable->GetName() << " = " <<
577 fScannedVariable->getVal() << "\n" <<
578 "\tCLs = " << hcResult->CLs() << " +/- " << hcResult->CLsError() << "\n" <<
579 "\tCLb = " << hcResult->CLb() << " +/- " << hcResult->CLbError() << "\n" <<
580 "\tCLsplusb = " << hcResult->CLsplusb() << " +/- " << hcResult->CLsplusbError() << "\n" <<
581 std::endl;
582 }
583
585 fTotalToysRun += (hcResult->GetAltDistribution()->GetSize() + hcResult->GetNullDistribution()->GetSize());
586
587 // set sampling distribution name
588 TString nullDistName = TString::Format("%s_%s_%4.2f",hcResult->GetNullDistribution()->GetName(),
590 TString altDistName = TString::Format("%s_%s_%4.2f",hcResult->GetAltDistribution()->GetName(),
592
593 hcResult->GetNullDistribution()->SetName(nullDistName);
594 hcResult->GetAltDistribution()->SetName(altDistName);
595 }
596
597 return hcResult;
598}
599
600////////////////////////////////////////////////////////////////////////////////
601/// Run a Fixed scan in npoints between min and max
602
603bool HypoTestInverter::RunFixedScan( int nBins, double xMin, double xMax, bool scanLog ) const
604{
605
607 // interpolate the limits
610
611 // safety checks
612 if ( nBins<=0 ) {
613 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - Please provide nBins>0\n";
614 return false;
615 }
616 if ( nBins==1 && xMin!=xMax ) {
617 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - nBins==1 -> I will run for xMin (" << xMin << ")\n";
618 }
619 if ( xMin==xMax && nBins>1 ) {
620 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - xMin==xMax -> I will enforce nBins==1\n";
621 nBins = 1;
622 }
623 if ( xMin>xMax ) {
624 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - Please provide xMin ("
625 << xMin << ") smaller than xMax (" << xMax << ")\n";
626 return false;
627 }
628
629 if (xMin < fScannedVariable->getMin()) {
630 xMin = fScannedVariable->getMin();
631 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - xMin < lower bound, using xmin = "
632 << xMin << std::endl;
633 }
634 if (xMax > fScannedVariable->getMax()) {
635 xMax = fScannedVariable->getMax();
636 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - xMax > upper bound, using xmax = "
637 << xMax << std::endl;
638 }
639
640 if (xMin <= 0. && scanLog) {
641 oocoutE(nullptr, InputArguments) << "HypoTestInverter::RunFixedScan - cannot go in log steps if xMin <= 0" << std::endl;
642 return false;
643 }
644
645 double thisX = xMin;
646 for (int i=0; i<nBins; i++) {
647
648 if (i > 0) { // avoids case of nBins = 1
649 if (scanLog)
650 thisX = exp( log(xMin) + i*(log(xMax)-log(xMin))/(nBins-1) ); // scan in log x
651 else
652 thisX = xMin + i*(xMax-xMin)/(nBins-1); // linear scan in x
653 }
654
655 const bool status = RunOnePoint(thisX);
656
657 // check if failed status
658 if ( status==false ) {
659 oocoutW(nullptr,Eval) << "HypoTestInverter::RunFixedScan - The hypo test for point " << thisX << " failed. Skipping." << std::endl;
660 }
661 }
662
663 return true;
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// run only one point at the given POI value
668
669bool HypoTestInverter::RunOnePoint( double rVal, bool adaptive, double clTarget) const
670{
671
673
674 // check if rVal is in the range specified for fScannedVariable
675 if ( rVal < fScannedVariable->getMin() ) {
676 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunOnePoint - Out of range: using the lower bound "
678 << " on the scanned variable rather than " << rVal<< "\n";
679 rVal = fScannedVariable->getMin();
680 }
681 if ( rVal > fScannedVariable->getMax() ) {
682 // print a message when you have a significative difference since rval is computed
683 if ( rVal > fScannedVariable->getMax()*(1.+1.E-12) )
684 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunOnePoint - Out of range: using the upper bound "
686 << " on the scanned variable rather than " << rVal<< "\n";
687 rVal = fScannedVariable->getMax();
688 }
689
690 // save old value
691 double oldValue = fScannedVariable->getVal();
692
693 // evaluate hybrid calculator at a single point
695 // need to set value of rval in hybridcalculator
696 // assume null model is S+B and alternate is B only
697 const ModelConfig * sbModel = fCalculator0->GetNullModel();
698 RooArgSet poi; poi.add(*sbModel->GetParametersOfInterest());
699 // set poi to right values
701 const_cast<ModelConfig*>(sbModel)->SetSnapshot(poi);
702
703 if (fVerbose > 0)
704 oocoutP(nullptr,Eval) << "Running for " << fScannedVariable->GetName() << " = " << fScannedVariable->getVal() << endl;
705
706 // compute the results
707 std::unique_ptr<HypoTestResult> result( Eval(*fCalculator0,adaptive,clTarget) );
708 if (!result) {
709 oocoutE(nullptr,Eval) << "HypoTestInverter - Error running point " << fScannedVariable->GetName() << " = " <<
710 fScannedVariable->getVal() << endl;
711 return false;
712 }
713 // in case of a dummy result
714 const double nullPV = result->NullPValue();
715 const double altPV = result->AlternatePValue();
716 if (!std::isfinite(nullPV) || nullPV < 0. || nullPV > 1. || !std::isfinite(altPV) || altPV < 0. || altPV > 1.) {
717 oocoutW(nullptr,Eval) << "HypoTestInverter - Skipping invalid result for point " << fScannedVariable->GetName() << " = " <<
718 fScannedVariable->getVal() << ". null p-value=" << nullPV << ", alternate p-value=" << altPV << endl;
719 return false;
720 }
721
722 double lastXtested;
723 if ( fResults->ArraySize()!=0 ) lastXtested = fResults->GetXValue(fResults->ArraySize()-1);
724 else lastXtested = -999;
725
726 if ( (std::abs(rVal) < 1 && TMath::AreEqualAbs(rVal, lastXtested,1.E-12) ) ||
727 (std::abs(rVal) >= 1 && TMath::AreEqualRel(rVal, lastXtested,1.E-12) ) ) {
728
729 oocoutI(nullptr,Eval) << "HypoTestInverter::RunOnePoint - Merge with previous result for "
730 << fScannedVariable->GetName() << " = " << rVal << std::endl;
732 if (prevResult && prevResult->GetNullDistribution() && prevResult->GetAltDistribution()) {
733 prevResult->Append(result.get());
734 }
735 else {
736 // if it was empty we re-use it
737 oocoutI(nullptr,Eval) << "HypoTestInverter::RunOnePoint - replace previous empty result\n";
738 auto oldObj = fResults->fYObjects.Remove(prevResult);
739 delete oldObj;
740
741 fResults->fYObjects.Add(result.release());
742 }
743
744 } else {
745
746 // fill the results in the HypoTestInverterResult array
747 fResults->fXValues.push_back(rVal);
748 fResults->fYObjects.Add(result.release());
749
750 }
751
752 fScannedVariable->setVal(oldValue);
753
754 return true;
755}
756
757////////////////////////////////////////////////////////////////////////////////
758/// Run an automatic scan until the desired accuracy is reached.
759/// Start by default from the full interval (min,max) of the POI and then via bisection find the line crossing
760/// the target line.
761/// Optionally, a hint can be provided and the scan will be done closer to that value.
762/// If by bisection the desired accuracy will not be reached, a fit to the points is performed.
763/// \param[out] limit The limit.
764/// \param[out] limitErr The error of the limit.
765/// \param[in] absAccuracy Desired absolute accuracy.
766/// \param[in] relAccuracy Desired relative accuracy.
767/// \param[in] hint Hint to start from or nullptr for no hint.
768
769bool HypoTestInverter::RunLimit(double &limit, double &limitErr, double absAccuracy, double relAccuracy, const double*hint) const {
770
771
772 // routine from G. Petrucciani (from HiggsCombination CMS package)
773
775
776 if ((hint != 0) && (*hint > r->getMin())) {
777 r->setMax(std::min<double>(3.0 * (*hint), r->getMax()));
778 r->setMin(std::max<double>(0.3 * (*hint), r->getMin()));
779 oocoutI(nullptr,InputArguments) << "HypoTestInverter::RunLimit - Use hint value " << *hint
780 << " search in interval " << r->getMin() << " , " << r->getMax() << std::endl;
781 }
782
783 // if not specified use the default values for rel and absolute accuracy
784 if (absAccuracy <= 0) absAccuracy = fgAbsAccuracy;
785 if (relAccuracy <= 0) relAccuracy = fgRelAccuracy;
786
787 typedef std::pair<double,double> CLs_t;
788 double clsTarget = fSize;
789 CLs_t clsMin(1,0), clsMax(0,0), clsMid(0,0);
790 double rMin = r->getMin(), rMax = r->getMax();
791 limit = 0.5*(rMax + rMin);
792 limitErr = 0.5*(rMax - rMin);
793 bool done = false;
794
795 TF1 expoFit("expoFit","[0]*exp([1]*(x-[2]))", rMin, rMax);
796
797 fLimitPlot = std::make_unique<TGraphErrors>();
798
799 if (fVerbose > 0) std::cout << "Search for upper limit to the limit" << std::endl;
800 for (int tries = 0; tries < 6; ++tries) {
801 if (! RunOnePoint(rMax) ) {
802 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypotest failed at upper limit of scan range: " << rMax << std::endl;
803 rMax *= 0.95;
804 continue;
805 }
806 clsMax = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
807 if (clsMax.first == 0 || clsMax.first + 3 * std::abs(clsMax.second) < clsTarget ) break;
808 rMax += rMax;
809 if (tries == 5) {
810 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Cannot determine upper limit of scan range. At " << r->GetName()
811 << " = " << rMax << " still getting "
812 << (fUseCLs ? "CLs" : "CLsplusb") << " = " << clsMax.first << std::endl;
813 return false;
814 }
815 }
816 if (fVerbose > 0) {
817 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Search for lower limit to the limit" << std::endl;
818 }
819
820 if ( fUseCLs && rMin == 0 ) {
821 clsMin = CLs_t(1,0);
822 }
823 else {
824 if (! RunOnePoint(rMin) ) {
825 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypotest failed at lower limit of scan range: " << rMin << std::endl;
826 return false;
827 }
828 clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
829 }
830 if (clsMin.first != 1 && clsMin.first - 3 * std::abs(clsMin.second) < clsTarget) {
831 if (fUseCLs) {
832 rMin = 0;
833 clsMin = CLs_t(1,0); // this is always true for CLs
834 } else {
835 rMin = -rMax / 4;
836 for (int tries = 0; tries < 6; ++tries) {
837 if (! RunOnePoint(rMin) ) {
838 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypotest failed at lower limit of scan range: " << rMin << std::endl;
839 rMin = rMin == 0. ? 0.1 : rMin * 1.1;
840 continue;
841 }
842 clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
843 if (clsMin.first == 1 || clsMin.first - 3 * std::abs(clsMin.second) > clsTarget) break;
844 rMin += rMin;
845 if (tries == 5) {
846 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Cannot determine lower limit of scan range. At " << r->GetName()
847 << " = " << rMin << " still get " << (fUseCLs ? "CLs" : "CLsplusb")
848 << " = " << clsMin.first << std::endl;
849 return false;
850 }
851 }
852 }
853 }
854
855 if (fVerbose > 0)
856 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Now doing proper bracketing & bisection" << std::endl;
857 do {
858
859 // break loop in case max toys is reached
860 if (fMaxToys > 0 && fTotalToysRun > fMaxToys ) {
861 oocoutW(nullptr,Eval) << "HypoTestInverter::RunLimit - maximum number of toys reached " << std::endl;
862 done = false; break;
863 }
864
865
866 // determine point by bisection or interpolation
867 limit = 0.5*(rMin+rMax); limitErr = 0.5*(rMax-rMin);
868 if (fgAlgo == "logSecant" && clsMax.first != 0) {
869 double logMin = log(clsMin.first), logMax = log(clsMax.first), logTarget = log(clsTarget);
870 limit = rMin + (rMax-rMin) * (logTarget - logMin)/(logMax - logMin);
871 if (clsMax.second != 0 && clsMin.second != 0) {
872 limitErr = hypot((logTarget-logMax) * (clsMin.second/clsMin.first), (logTarget-logMin) * (clsMax.second/clsMax.first));
873 limitErr *= (rMax-rMin)/((logMax-logMin)*(logMax-logMin));
874 }
875 }
876 r->setError(limitErr);
877
878 // exit if reached accuracy on r
879 if (limitErr < std::max(absAccuracy, relAccuracy * limit)) {
880 if (fVerbose > 1)
881 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - reached accuracy " << limitErr
882 << " below " << std::max(absAccuracy, relAccuracy * limit) << std::endl;
883 done = true; break;
884 }
885
886 // evaluate point
887 if (! RunOnePoint(limit, true, clsTarget) ) {
888 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypo test failed at x=" << limit << " when trying to find limit." << std::endl;
889 return false;
890 }
891 clsMid = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
892
893 if (clsMid.second == -1) {
894 std::cerr << "Hypotest failed" << std::endl;
895 return false;
896 }
897
898 // if sufficiently far away, drop one of the points
899 if (std::abs(clsMid.first-clsTarget) >= 2*clsMid.second) {
900 if ((clsMid.first>clsTarget) == (clsMax.first>clsTarget)) {
901 rMax = limit; clsMax = clsMid;
902 } else {
903 rMin = limit; clsMin = clsMid;
904 }
905 } else {
906 if (fVerbose > 0) std::cout << "Trying to move the interval edges closer" << std::endl;
907 double rMinBound = rMin, rMaxBound = rMax;
908 // try to reduce the size of the interval
909 while (clsMin.second == 0 || std::abs(rMin-limit) > std::max(absAccuracy, relAccuracy * limit)) {
910 rMin = 0.5*(rMin+limit);
911 if (!RunOnePoint(rMin,true, clsTarget) ) {
912 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypo test failed at x=" << rMin << " when trying to find limit from below." << std::endl;
913 return false;
914 }
915 clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
916 if (std::abs(clsMin.first-clsTarget) <= 2*clsMin.second) break;
917 rMinBound = rMin;
918 }
919 while (clsMax.second == 0 || std::abs(rMax-limit) > std::max(absAccuracy, relAccuracy * limit)) {
920 rMax = 0.5*(rMax+limit);
921 if (!RunOnePoint(rMax,true,clsTarget) ) {
922 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypo test failed at x=" << rMin << " when trying to find limit from above." << std::endl;
923 return false;
924 }
925 clsMax = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
926 if (std::abs(clsMax.first-clsTarget) <= 2*clsMax.second) break;
927 rMaxBound = rMax;
928 }
929 expoFit.SetRange(rMinBound,rMaxBound);
930 break;
931 }
932 } while (true);
933
934 if (!done) { // didn't reach accuracy with scan, now do fit
935 if (fVerbose) {
936 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Before fit --- \n";
937 std::cout << "Limit: " << r->GetName() << " < " << limit << " +/- " << limitErr << " [" << rMin << ", " << rMax << "]\n";
938 }
939
940 expoFit.FixParameter(0,clsTarget);
941 expoFit.SetParameter(1,log(clsMax.first/clsMin.first)/(rMax-rMin));
942 expoFit.SetParameter(2,limit);
943 double rMinBound, rMaxBound; expoFit.GetRange(rMinBound, rMaxBound);
944 limitErr = std::max(std::abs(rMinBound-limit), std::abs(rMaxBound-limit));
945 int npoints = 0;
946
947 HypoTestInverterPlot plot("plot","plot",fResults);
948 fLimitPlot.reset(plot.MakePlot() );
949
950
951 for (int j = 0; j < fLimitPlot->GetN(); ++j) {
952 if (fLimitPlot->GetX()[j] >= rMinBound && fLimitPlot->GetX()[j] <= rMaxBound) npoints++;
953 }
954 for (int i = 0, imax = /*(readHybridResults_ ? 0 : */ 8; i <= imax; ++i, ++npoints) {
955 fLimitPlot->Sort();
956 fLimitPlot->Fit(&expoFit,(fVerbose <= 1 ? "QNR EX0" : "NR EXO"));
957 if (fVerbose) {
958 oocoutI(nullptr,Eval) << "Fit to " << npoints << " points: " << expoFit.GetParameter(2) << " +/- " << expoFit.GetParError(2) << std::endl;
959 }
960 if ((rMin < expoFit.GetParameter(2)) && (expoFit.GetParameter(2) < rMax) && (expoFit.GetParError(2) < 0.5*(rMaxBound-rMinBound))) {
961 // sanity check fit result
962 limit = expoFit.GetParameter(2);
963 limitErr = expoFit.GetParError(2);
964 if (limitErr < std::max(absAccuracy, relAccuracy * limit)) break;
965 }
966 // add one point in the interval.
967 double rTry = RooRandom::uniform()*(rMaxBound-rMinBound)+rMinBound;
968 if (i != imax) {
969 if (!RunOnePoint(rTry,true,clsTarget) ) return false;
970 //eval(w, mc_s, mc_b, data, rTry, true, clsTarget);
971 }
972
973 }
974 }
975
976//if (!plot_.empty() && fLimitPlot.get()) {
977 if (fLimitPlot.get() && fLimitPlot->GetN() > 0) {
978 //new TCanvas("c1","c1");
979 fLimitPlot->Sort();
980 fLimitPlot->SetLineWidth(2);
981 double xmin = r->getMin(), xmax = r->getMax();
982 for (int j = 0; j < fLimitPlot->GetN(); ++j) {
983 if (fLimitPlot->GetY()[j] > 1.4*clsTarget || fLimitPlot->GetY()[j] < 0.6*clsTarget) continue;
984 xmin = std::min(fLimitPlot->GetX()[j], xmin);
985 xmax = std::max(fLimitPlot->GetX()[j], xmax);
986 }
987 fLimitPlot->GetXaxis()->SetRangeUser(xmin,xmax);
988 fLimitPlot->GetYaxis()->SetRangeUser(0.5*clsTarget, 1.5*clsTarget);
989 fLimitPlot->Draw("AP");
990 expoFit.Draw("SAME");
991 TLine line(fLimitPlot->GetX()[0], clsTarget, fLimitPlot->GetX()[fLimitPlot->GetN()-1], clsTarget);
993 line.DrawLine(limit, 0, limit, fLimitPlot->GetY()[0]);
995 line.DrawLine(limit-limitErr, 0, limit-limitErr, fLimitPlot->GetY()[0]);
996 line.DrawLine(limit+limitErr, 0, limit+limitErr, fLimitPlot->GetY()[0]);
997 //c1->Print(plot_.c_str());
998 }
999
1000 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Result: \n"
1001 << "\tLimit: " << r->GetName() << " < " << limit << " +/- " << limitErr << " @ " << (1-fSize) * 100 << "% CL\n";
1002 if (fVerbose > 1) oocoutI(nullptr,Eval) << "Total toys: " << fTotalToysRun << std::endl;
1003
1004 // set value in results
1005 fResults->fUpperLimit = limit;
1006 fResults->fUpperLimitError = limitErr;
1008 // lower limit are always min of p value
1012
1013 return true;
1014}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// get the distribution of lower limit
1018/// if rebuild = false (default) it will re-use the results of the scan done
1019/// for obtained the observed limit and no extra toys will be generated
1020/// if rebuild a new set of B toys will be done and the procedure will be repeated
1021/// for each toy
1022
1024 if (!rebuild) {
1025 if (!fResults) {
1026 oocoutE(nullptr,InputArguments) << "HypoTestInverter::GetLowerLimitDistribution(false) - result not existing\n";
1027 return 0;
1028 }
1030 }
1031
1032 TList * clsDist = 0;
1033 TList * clsbDist = 0;
1034 if (fUseCLs) clsDist = &fResults->fExpPValues;
1035 else clsbDist = &fResults->fExpPValues;
1036
1037 return RebuildDistributions(false, nToys,clsDist, clsbDist);
1038
1039}
1040
1041////////////////////////////////////////////////////////////////////////////////
1042/// get the distribution of lower limit
1043/// if rebuild = false (default) it will re-use the results of the scan done
1044/// for obtained the observed limit and no extra toys will be generated
1045/// if rebuild a new set of B toys will be done and the procedure will be repeated
1046/// for each toy
1047/// The nuisance parameter value used for rebuild is the current one in the model
1048/// so it is user responsibility to set to the desired value (nomi
1049
1051 if (!rebuild) {
1052 if (!fResults) {
1053 oocoutE(nullptr,InputArguments) << "HypoTestInverter::GetUpperLimitDistribution(false) - result not existing\n";
1054 return 0;
1055 }
1057 }
1058
1059 TList * clsDist = 0;
1060 TList * clsbDist = 0;
1061 if (fUseCLs) clsDist = &fResults->fExpPValues;
1062 else clsbDist = &fResults->fExpPValues;
1063
1064 return RebuildDistributions(true, nToys,clsDist, clsbDist);
1065}
1066
1067////////////////////////////////////////////////////////////////////////////////
1068
1071}
1072
1073////////////////////////////////////////////////////////////////////////////////
1074/// rebuild the sampling distributions by
1075/// generating some toys and find for each of them a new upper limit
1076/// Return the upper limit distribution and optionally also the pValue distributions for Cls, Clsb and Clbxs
1077/// as a TList for each scanned point
1078/// The method uses the present parameter value. It is user responsibility to give the current parameters to rebuild the distributions
1079/// 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
1080/// output file as an histogram
1081
1082SamplingDistribution * HypoTestInverter::RebuildDistributions(bool isUpper, int nToys, TList * clsDist, TList * clsbDist, TList * clbDist, const char *outputfile) {
1083
1084 if (!fScannedVariable || !fCalculator0) return 0;
1085 // get first background snapshot
1086 const ModelConfig * bModel = fCalculator0->GetAlternateModel();
1087 const ModelConfig * sbModel = fCalculator0->GetNullModel();
1088 if (!bModel || ! sbModel) return 0;
1089 RooArgSet paramPoint;
1090 if (!sbModel->GetParametersOfInterest()) return 0;
1091 paramPoint.add(*sbModel->GetParametersOfInterest());
1092
1093 const RooArgSet * poibkg = bModel->GetSnapshot();
1094 if (!poibkg) {
1095 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RebuildDistribution - background snapshot not existing"
1096 << " assume is for POI = 0" << std::endl;
1098 paramPoint.assign(RooArgSet(*fScannedVariable));
1099 }
1100 else
1101 paramPoint.assign(*poibkg);
1102 // generate data at bkg parameter point
1103
1104 ToyMCSampler * toymcSampler = dynamic_cast<ToyMCSampler *>(fCalculator0->GetTestStatSampler() );
1105 if (!toymcSampler) {
1106 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RebuildDistribution - no toy MC sampler existing" << std::endl;
1107 return 0;
1108 }
1109 // set up test stat sampler in case of asymptotic calculator
1110 if (dynamic_cast<RooStats::AsymptoticCalculator*>(fCalculator0) ) {
1111 toymcSampler->SetObservables(*sbModel->GetObservables() );
1112 toymcSampler->SetParametersForTestStat(*sbModel->GetParametersOfInterest());
1113 toymcSampler->SetPdf(*sbModel->GetPdf());
1114 toymcSampler->SetNuisanceParameters(*sbModel->GetNuisanceParameters());
1115 if (sbModel->GetGlobalObservables() ) toymcSampler->SetGlobalObservables(*sbModel->GetGlobalObservables() );
1116 // set number of events
1117 if (!sbModel->GetPdf()->canBeExtended())
1118 toymcSampler->SetNEventsPerToy(1);
1119 }
1120
1121 // loop on data to generate
1122 int nPoints = fNBins;
1123
1124 bool storePValues = clsDist || clsbDist || clbDist;
1125 if (fNBins <=0 && storePValues) {
1126 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RebuildDistribution - cannot return p values distribution with the auto scan" << std::endl;
1127 storePValues = false;
1128 nPoints = 0;
1129 }
1130
1131 if (storePValues) {
1132 if (fResults) nPoints = fResults->ArraySize();
1133 if (nPoints <=0) {
1134 oocoutE(nullptr,InputArguments) << "HypoTestInverter - result is not existing and number of point to scan is not set"
1135 << std::endl;
1136 return 0;
1137 }
1138 }
1139
1140 if (nToys <= 0) nToys = 100; // default value
1141
1142 std::vector<std::vector<double> > CLs_values(nPoints);
1143 std::vector<std::vector<double> > CLsb_values(nPoints);
1144 std::vector<std::vector<double> > CLb_values(nPoints);
1145
1146 if (storePValues) {
1147 for (int i = 0; i < nPoints; ++i) {
1148 CLs_values[i].reserve(nToys);
1149 CLb_values[i].reserve(nToys);
1150 CLsb_values[i].reserve(nToys);
1151 }
1152 }
1153
1154 std::vector<double> limit_values; limit_values.reserve(nToys);
1155
1156 oocoutI(nullptr,InputArguments) << "HypoTestInverter - rebuilding the p value distributions by generating ntoys = "
1157 << nToys << std::endl;
1158
1159
1160 oocoutI(nullptr,InputArguments) << "Rebuilding using parameter of interest point: ";
1161 RooStats::PrintListContent(paramPoint, oocoutI(nullptr,InputArguments) );
1162 if (sbModel->GetNuisanceParameters() ) {
1163 oocoutI(nullptr,InputArguments) << "And using nuisance parameters: ";
1164 RooStats::PrintListContent(*sbModel->GetNuisanceParameters(), oocoutI(nullptr,InputArguments) );
1165 }
1166 // save all parameters to restore them later
1167 assert(bModel->GetPdf() );
1168 assert(bModel->GetObservables() );
1169 std::unique_ptr<RooArgSet> allParams{bModel->GetPdf()->getParameters( *bModel->GetObservables() )};
1170 RooArgSet saveParams;
1171 allParams->snapshot(saveParams);
1172
1173 std::unique_ptr<TFile> fileOut{TFile::Open(outputfile,"RECREATE")};
1174 if (!fileOut) {
1175 oocoutE(nullptr,InputArguments) << "HypoTestInverter - RebuildDistributions - Error opening file " << outputfile
1176 << " - the resulting limits will not be stored" << std::endl;
1177 }
1178 // create temporary histograms to store the limit result
1179 TH1D * hL = new TH1D("lowerLimitDist","Rebuilt lower limit distribution",100,1.,0.);
1180 TH1D * hU = new TH1D("upperLimitDist","Rebuilt upper limit distribution",100,1.,0.);
1181 TH1D * hN = new TH1D("nObs","Observed events",100,1.,0.);
1182 hL->SetBuffer(2*nToys);
1183 hU->SetBuffer(2*nToys);
1184 std::vector<TH1*> hCLb;
1185 std::vector<TH1*> hCLsb;
1186 std::vector<TH1*> hCLs;
1187 if (storePValues) {
1188 for (int i = 0; i < nPoints; ++i) {
1189 hCLb.push_back(new TH1D(TString::Format("CLbDist_bin%d",i),"CLb distribution",100,1.,0.));
1190 hCLs.push_back(new TH1D(TString::Format("ClsDist_bin%d",i),"CLs distribution",100,1.,0.));
1191 hCLsb.push_back(new TH1D(TString::Format("CLsbDist_bin%d",i),"CLs+b distribution",100,1.,0.));
1192 }
1193 }
1194
1195
1196 // loop now on the toys
1197 for (int itoy = 0; itoy < nToys; ++itoy) {
1198
1199 oocoutP(nullptr,Eval) << "\nHypoTestInverter - RebuildDistributions - running toy # " << itoy << " / "
1200 << nToys << std::endl;
1201
1202
1203 printf("\n\nshnapshot of s+b model \n");
1204 sbModel->GetSnapshot()->Print("v");
1205
1206 // reset parameters to initial values to be sure in case they are not reset
1207 if (itoy> 0) {
1208 allParams->assign(saveParams);
1209 }
1210
1211 // need to set th epdf to clear the cache in ToyMCSampler
1212 // pdf we must use is background pdf
1213 toymcSampler->SetPdf(*bModel->GetPdf() );
1214
1215
1216 RooAbsData * bkgdata = toymcSampler->GenerateToyData(paramPoint);
1217
1218 double nObs = bkgdata->sumEntries();
1219 // for debugging in case of number counting models
1220 if (bkgdata->numEntries() ==1 && !bModel->GetPdf()->canBeExtended()) {
1221 oocoutP(nullptr,Generation) << "Generate observables are : ";
1222 RooArgList genObs(*bkgdata->get(0));
1223 RooStats::PrintListContent(genObs, oocoutP(nullptr,Generation) );
1224 nObs = 0;
1225 for (int i = 0; i < genObs.getSize(); ++i) {
1226 RooRealVar * x = dynamic_cast<RooRealVar*>(&genObs[i]);
1227 if (x) nObs += x->getVal();
1228 }
1229 }
1230 hN->Fill(nObs);
1231
1232 // by copying I will have the same min/max as previous ones
1233 HypoTestInverter inverter = *this;
1234 inverter.SetData(*bkgdata);
1235
1236 // print global observables
1237 auto gobs = bModel->GetPdf()->getVariables()->selectCommon(* sbModel->GetGlobalObservables() );
1238 gobs->Print("v");
1239
1240 HypoTestInverterResult * r = inverter.GetInterval();
1241
1242 if (r == 0) continue;
1243
1244 double value = (isUpper) ? r->UpperLimit() : r->LowerLimit();
1245 limit_values.push_back( value );
1246 hU->Fill(r->UpperLimit() );
1247 hL->Fill(r->LowerLimit() );
1248
1249
1250 std::cout << "The computed upper limit for toy #" << itoy << " is " << value << std::endl;
1251
1252 // write every 10 toys
1253 if (itoy%10 == 0 || itoy == nToys-1) {
1254 hU->Write("",TObject::kOverwrite);
1255 hL->Write("",TObject::kOverwrite);
1256 hN->Write("",TObject::kOverwrite);
1257 }
1258
1259 if (!storePValues) continue;
1260
1261 if (nPoints < r->ArraySize()) {
1262 oocoutW(nullptr,InputArguments) << "HypoTestInverter: skip extra points" << std::endl;
1263 }
1264 else if (nPoints > r->ArraySize()) {
1265 oocoutW(nullptr,InputArguments) << "HypoTestInverter: missing some points" << std::endl;
1266 }
1267
1268
1269 for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1270 HypoTestResult * hr = r->GetResult(ipoint);
1271 if (hr) {
1272 CLs_values[ipoint].push_back( hr->CLs() );
1273 CLsb_values[ipoint].push_back( hr->CLsplusb() );
1274 CLb_values[ipoint].push_back( hr->CLb() );
1275 hCLs[ipoint]->Fill( hr->CLs() );
1276 hCLb[ipoint]->Fill( hr->CLb() );
1277 hCLsb[ipoint]->Fill( hr->CLsplusb() );
1278 }
1279 else {
1280 oocoutW(nullptr,InputArguments) << "HypoTestInverter: missing result for point: x = "
1281 << fResults->GetXValue(ipoint) << std::endl;
1282 }
1283 }
1284 // write every 10 toys
1285 if (itoy%10 == 0 || itoy == nToys-1) {
1286 for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1287 hCLs[ipoint]->Write("",TObject::kOverwrite);
1288 hCLb[ipoint]->Write("",TObject::kOverwrite);
1289 hCLsb[ipoint]->Write("",TObject::kOverwrite);
1290 }
1291 }
1292
1293
1294 delete r;
1295 delete bkgdata;
1296 }
1297
1298
1299 if (storePValues) {
1300 if (clsDist) clsDist->SetOwner(true);
1301 if (clbDist) clbDist->SetOwner(true);
1302 if (clsbDist) clsbDist->SetOwner(true);
1303
1304 oocoutI(nullptr,InputArguments) << "HypoTestInverter: storing rebuilt p values " << std::endl;
1305
1306 for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1307 if (clsDist) {
1308 TString name = TString::Format("CLs_distrib_%d",ipoint);
1309 clsDist->Add( new SamplingDistribution(name,name,CLs_values[ipoint] ) );
1310 }
1311 if (clbDist) {
1312 TString name = TString::Format("CLb_distrib_%d",ipoint);
1313 clbDist->Add( new SamplingDistribution(name,name,CLb_values[ipoint] ) );
1314 }
1315 if (clsbDist) {
1316 TString name = TString::Format("CLsb_distrib_%d",ipoint);
1317 clsbDist->Add( new SamplingDistribution(name,name,CLsb_values[ipoint] ) );
1318 }
1319 }
1320 }
1321
1322 else {
1323 // delete all the histograms
1324 delete hL;
1325 delete hU;
1326 for (int i = 0; i < nPoints && storePValues; ++i) {
1327 delete hCLs[i];
1328 delete hCLb[i];
1329 delete hCLsb[i];
1330 }
1331 }
1332
1333 const char * disName = (isUpper) ? "upperLimit_dist" : "lowerLimit_dist";
1334 return new SamplingDistribution(disName, disName, limit_values);
1335}
size_t fSize
#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:377
@ kRed
Definition Rtypes.h:66
winID h TVirtualViewer3D TVirtualGLPainter char TVirtualGLPainter plot
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
RooFit::OwningPtr< 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...
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
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
RooAbsArg * find(const char *name) const
Find object with given name in list.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:59
virtual double sumEntries() const =0
Return effective number of entries in dataset, i.e., sum all weights.
virtual const RooArgSet * get() const
Definition RooAbsData.h:103
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:278
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
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:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:81
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setVal(double value) override
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
void SetData(RooAbsData &data) override
Set the DataSet.
HypoTestResult * GetHypoTest() const override
inherited methods from HypoTestCalculator interface
const ModelConfig * GetAlternateModel(void) const
TestStatSampler * GetTestStatSampler(void) const
Returns instance of TestStatSampler.
Class to plot a HypoTestInverterResult, the output of the HypoTestInverter calculator.
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
int ArraySize() const
number of entries in the results array
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 expected sampling distribution for each point
bool fIsTwoSided
two sided scan (look for lower/upper limit)
TList fYObjects
list of HypoTestResult for each point
void SetConfidenceLevel(double cl) override
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval)
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...
static void CheckInputModels(const HypoTestCalculatorGeneric &hc, const RooRealVar &scanVar)
check the model given the given hypotestcalculator
RooRealVar * fScannedVariable
pointer to the constrained variable
std::unique_ptr< HypoTestCalculatorGeneric > fHC
! pointer to the generic hypotest calculator used
bool RunLimit(double &limit, double &limitErr, double absTol=0, double relTol=0, const double *hint=nullptr) const
Run an automatic scan until the desired accuracy is reached.
void SetData(RooAbsData &) override
Set the DataSet ( add to the workspace if not already there ?)
HypoTestInverterResult * fResults
pointer to the result
SamplingDistribution * RebuildDistributions(bool isUpper=true, int nToys=100, TList *clsDist=nullptr, TList *clsbDist=nullptr, TList *clbDist=nullptr, const char *outputfile="HypoTestInverterRebuiltDist.root")
function to rebuild the distributions
int fMaxToys
maximum number of toys to run
HypoTestInverter()
default constructor (used only for I/O)
~HypoTestInverter() override
destructor
double ConfidenceLevel() const override
Get the Confidence level for the test.
bool SetTestStatistic(TestStatistic &stat)
set the test statistic
HypoTestInverterResult * GetInterval() const override
Run a fixed scan or the automatic scan depending on the configuration.
static void SetCloseProof(bool flag)
set flag to close proof for every new run
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
get the test statistic
std::unique_ptr< TGraphErrors > fLimitPlot
! plot of limits
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 upper/lower limit distribution
HypoTestInverter & operator=(const HypoTestInverter &rhs)
assignment
HypoTestCalculatorGeneric * fCalculator0
pointer to the calculator passed in the constructor
void CreateResults() const
create a new HypoTestInverterResult to hold all computed results
static RooRealVar * GetVariableToScan(const HypoTestCalculatorGeneric &hc)
helper functions
HypoTestResult * Eval(HypoTestCalculatorGeneric &hc, bool adaptive, double clsTarget) const
run the hybrid at a single point
HypoTestResult is a base class for results from hypothesis tests.
virtual double CLsplusb() const
Convert AlternatePValue into a "confidence level".
double GetTestStatisticData(void) const
virtual void Append(const HypoTestResult *other)
add values from another HypoTestResult
double CLsError() const
The error on the ratio .
void SetBackgroundAsAlt(bool l=true)
bool GetPValueIsRightTail(void) const
double CLbError() const
The error on the "confidence level" of the null hypothesis.
void SetTestStatisticData(const double tsd)
double CLsplusbError() const
The error on the "confidence level" of the alternative hypothesis.
SamplingDistribution * GetNullDistribution(void) const
virtual double CLs() const
is simply (not a method, but a quantity)
virtual double CLb() const
Convert NullPValue into a "confidence level".
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:35
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr 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:91
This class simply holds a sampling distribution of some test statistic.
Int_t GetSize() const
size of samples
double fUpperLimit
upper interval limit
double fLowerLimit
lower interval limit
TestStatSampler is an interface class for a tools which produce RooStats SamplingDistributions.
virtual TestStatistic * GetTestStatistic() const =0
Get the TestStatistic.
virtual void SetTestStatistic(TestStatistic *testStatistic)=0
Set the TestStatistic (want the argument to be a function of the data & parameter points.
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
ToyMCSampler is an implementation of the TestStatSampler interface.
void SetParametersForTestStat(const RooArgSet &nullpoi) override
Set the Pdf, add to the workspace if not already there.
void SetObservables(const RooArgSet &o) override
specify the observables in the dataset (needed to evaluate the test statistic)
virtual RooAbsData * GenerateToyData(RooArgSet &paramPoint, RooAbsPdf &pdf) const
generates toy data without weight
void SetNuisanceParameters(const RooArgSet &np) override
specify the nuisance parameters (eg. the rest of the parameters)
void SetPdf(RooAbsPdf &pdf) override
Set the Pdf, add to the workspace if not already there.
void SetGlobalObservables(const RooArgSet &o) override
specify the conditional observables
virtual void SetNEventsPerToy(const Int_t nevents)
Forces the generation of exactly n events even for extended PDFs.
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:1931
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF1.cxx:3520
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1334
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition TF1.cxx:2280
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:639
virtual void FixParameter(Int_t ipar, Double_t value)
Fix the value of a parameter for a fit operation The specified value will be used in the fit and the ...
Definition TF1.cxx:1558
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:517
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:4053
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:620
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3338
virtual void SetBuffer(Int_t buffersize, Option_t *option="")
Set the maximum number of entries to be kept in the buffer.
Definition TH1.cxx:8324
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
void Add(TObject *obj) override
Definition TList.h:81
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:822
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
Definition TNamed.cxx:74
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
@ kOverwrite
overwrite existing object with same name
Definition TObject.h:92
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:874
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
Basic string class.
Definition TString.h:139
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:2356
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)
useful function to print in one line the content of a set with their values
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Comparing floating points.
Definition TMath.h:425
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Comparing floating points.
Definition TMath.h:418
static void SetToys(HypoTestType *h, int toyNull, int toyAlt)