Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
65#include <cassert>
66#include <cmath>
67#include <memory>
68
69
70using namespace RooStats;
71using std::endl;
72
73// static variable definitions
75unsigned int HypoTestInverter::fgNToys = 500;
76
79std::string HypoTestInverter::fgAlgo = "logSecant";
80
81// helper class to wrap the functionality of the various HypoTestCalculators
82
83template<class HypoTestType>
85
86 static void SetToys(HypoTestType * h, int toyNull, int toyAlt) { h->SetToys(toyNull,toyAlt); }
87
88};
89
90////////////////////////////////////////////////////////////////////////////////
91/// get the variable to scan
92/// try first with null model if not go to alternate model
93
95
96 RooRealVar * varToScan = nullptr;
97 const ModelConfig * mc = hc.GetNullModel();
98 if (mc) {
99 const RooArgSet * poi = mc->GetParametersOfInterest();
100 if (poi) varToScan = dynamic_cast<RooRealVar*> (poi->first() );
101 }
102 if (!varToScan) {
103 mc = hc.GetAlternateModel();
104 if (mc) {
105 const RooArgSet * poi = mc->GetParametersOfInterest();
106 if (poi) varToScan = dynamic_cast<RooRealVar*> (poi->first() );
107 }
108 }
109 return varToScan;
110}
111
112////////////////////////////////////////////////////////////////////////////////
113/// check the model given the given hypotestcalculator
114
116 const ModelConfig * modelSB = hc.GetNullModel();
117 const ModelConfig * modelB = hc.GetAlternateModel();
118 if (!modelSB || ! modelB)
119 oocoutF(nullptr,InputArguments) << "HypoTestInverter - model are not existing" << std::endl;
121
122 oocoutI(nullptr,InputArguments) << "HypoTestInverter ---- Input models: \n"
123 << "\t\t using as S+B (null) model : "
124 << modelSB->GetName() << "\n"
125 << "\t\t using as B (alternate) model : "
126 << modelB->GetName() << "\n" << std::endl;
127
128 // check if scanVariable is included in B model pdf
129 RooAbsPdf * bPdf = modelB->GetPdf();
130 const RooArgSet * bObs = modelB->GetObservables();
131 if (!bPdf || !bObs) {
132 oocoutE(nullptr,InputArguments) << "HypoTestInverter - B model has no pdf or observables defined" << std::endl;
133 return;
134 }
135 std::unique_ptr<RooArgSet> bParams{bPdf->getParameters(*bObs)};
136 if (!bParams) {
137 oocoutE(nullptr,InputArguments) << "HypoTestInverter - pdf of B model has no parameters" << std::endl;
138 return;
139 }
140 if (bParams->find(scanVariable.GetName() ) ) {
141 const RooArgSet * poiB = modelB->GetSnapshot();
142 if (!poiB || !poiB->find(scanVariable.GetName()) ||
143 (static_cast<RooRealVar *>(poiB->find(scanVariable.GetName())))->getVal() != 0) {
144 oocoutW(nullptr, InputArguments)
145 << "HypoTestInverter - using a B model with POI " << scanVariable.GetName() << " not equal to zero "
146 << " user must check input model configurations " << std::endl;
147 }
148 if (poiB) delete poiB;
149 }
150}
151
152////////////////////////////////////////////////////////////////////////////////
153/// default constructor (doesn't do anything)
154
156 fTotalToysRun(0),
157 fMaxToys(0),
158 fCalculator0(nullptr),
159 fScannedVariable(nullptr),
160 fResults(nullptr),
161 fUseCLs(false),
162 fScanLog(false),
163 fSize(0),
164 fVerbose(0),
165 fCalcType(kUndefined),
166 fNBins(0), fXmin(1), fXmax(1),
167 fNumErr(0)
168{
169}
170
171////////////////////////////////////////////////////////////////////////////////
172/// Constructor from a HypoTestCalculatorGeneric
173/// The HypoTest calculator must be a FrequentistCalculator or HybridCalculator type
174/// Other type of calculators are not supported.
175/// The calculator must be created before by using the S+B model for the null and
176/// the B model for the alt
177/// If no variable to scan are given they are assumed to be the first variable
178/// from the parameter of interests of the null model
179
181 RooRealVar* scannedVariable, double size ) :
182 fTotalToysRun(0),
183 fMaxToys(0),
184 fCalculator0(nullptr),
185 fScannedVariable(scannedVariable),
186 fResults(nullptr),
187 fUseCLs(false),
188 fScanLog(false),
189 fSize(size),
190 fVerbose(0),
191 fCalcType(kUndefined),
192 fNBins(0), fXmin(1), fXmax(1),
193 fNumErr(0)
194{
195
196 if (!fScannedVariable) {
198 }
199 if (!fScannedVariable) {
200 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
201 } else {
203 }
204
205 HybridCalculator * hybCalc = dynamic_cast<HybridCalculator*>(&hc);
206 if (hybCalc) {
209 return;
210 }
212 if (freqCalc) {
215 return;
216 }
218 if (asymCalc) {
221 return;
222 }
223 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Type of hypotest calculator is not supported " <<std::endl;
224 fCalculator0 = &hc;
225}
226
227////////////////////////////////////////////////////////////////////////////////
228/// Constructor from a reference to a HybridCalculator
229/// The calculator must be created before by using the S+B model for the null and
230/// the B model for the alt
231/// If no variable to scan are given they are assumed to be the first variable
232/// from the parameter of interests of the null model
233
235 RooRealVar* scannedVariable, double size ) :
236 fTotalToysRun(0),
237 fMaxToys(0),
238 fCalculator0(&hc),
239 fScannedVariable(scannedVariable),
240 fResults(nullptr),
241 fUseCLs(false),
242 fScanLog(false),
243 fSize(size),
244 fVerbose(0),
245 fCalcType(kHybrid),
246 fNBins(0), fXmin(1), fXmax(1),
247 fNumErr(0)
248{
249
250 if (!fScannedVariable) {
252 }
253 if (!fScannedVariable) {
254 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
255 } else {
257 }
258}
259
260////////////////////////////////////////////////////////////////////////////////
261/// Constructor from a reference to a FrequentistCalculator
262/// The calculator must be created before by using the S+B model for the null and
263/// the B model for the alt
264/// If no variable to scan are given they are assumed to be the first variable
265/// from the parameter of interests of the null model
266
268 RooRealVar* scannedVariable, double size ) :
269 fTotalToysRun(0),
270 fMaxToys(0),
271 fCalculator0(&hc),
272 fScannedVariable(scannedVariable),
273 fResults(nullptr),
274 fUseCLs(false),
275 fScanLog(false),
276 fSize(size),
277 fVerbose(0),
278 fCalcType(kFrequentist),
279 fNBins(0), fXmin(1), fXmax(1),
280 fNumErr(0)
281{
282
283 if (!fScannedVariable) {
285 }
286 if (!fScannedVariable) {
287 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
288 } else {
290 }
291}
292
293////////////////////////////////////////////////////////////////////////////////
294/// Constructor from a reference to a AsymptoticCalculator
295/// The calculator must be created before by using the S+B model for the null and
296/// the B model for the alt
297/// If no variable to scan are given they are assumed to be the first variable
298/// from the parameter of interests of the null model
299
301 RooRealVar* scannedVariable, double size ) :
302 fTotalToysRun(0),
303 fMaxToys(0),
304 fCalculator0(&hc),
305 fScannedVariable(scannedVariable),
306 fResults(nullptr),
307 fUseCLs(false),
308 fScanLog(false),
309 fSize(size),
310 fVerbose(0),
311 fCalcType(kAsymptotic),
312 fNBins(0), fXmin(1), fXmax(1),
313 fNumErr(0)
314{
315
316 if (!fScannedVariable) {
318 }
319 if (!fScannedVariable) {
320 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
321 } else {
323 }
324}
325
326////////////////////////////////////////////////////////////////////////////////
327/// Constructor from a model for B model and a model for S+B.
328/// An HypoTestCalculator (Hybrid of Frequentis) will be created using the
329/// S+B model as the null and the B model as the alternate
330/// If no variable to scan are given they are assumed to be the first variable
331/// from the parameter of interests of the null model
332
335 fTotalToysRun(0),
336 fMaxToys(0),
337 fCalculator0(nullptr),
338 fScannedVariable(scannedVariable),
339 fResults(nullptr),
340 fUseCLs(false),
341 fScanLog(false),
342 fSize(size),
343 fVerbose(0),
344 fCalcType(type),
345 fNBins(0), fXmin(1), fXmax(1),
346 fNumErr(0)
347{
348 if(fCalcType==kFrequentist) fHC = std::make_unique<FrequentistCalculator>(data, bModel, sbModel);
349 if(fCalcType==kHybrid) fHC = std::make_unique<HybridCalculator>(data, bModel, sbModel);
350 if(fCalcType==kAsymptotic) fHC = std::make_unique<AsymptoticCalculator>(data, bModel, sbModel);
351 fCalculator0 = fHC.get();
352 // get scanned variable
353 if (!fScannedVariable) {
355 }
356 if (!fScannedVariable) {
357 oocoutE(nullptr,InputArguments) << "HypoTestInverter - Cannot guess the variable to scan " << std::endl;
358 } else {
360 }
361}
362
363////////////////////////////////////////////////////////////////////////////////
364/// copy-constructor
365/// NOTE: this class does not copy the contained result and
366/// the HypoTestCalculator, but only the pointers
367/// It requires the original HTI to be alive
368
371 fTotalToysRun(0),
372 fCalculator0(nullptr), fScannedVariable(nullptr), // add these for Coverity
373 fResults(nullptr)
374{
375 (*this) = rhs;
376}
377
378////////////////////////////////////////////////////////////////////////////////
379/// assignment operator
380/// NOTE: this class does not copy the contained result and
381/// the HypoTestCalculator, but only the pointers
382/// It requires the original HTI to be alive
383
385 if (this == &rhs) return *this;
386 fTotalToysRun = 0;
387 fMaxToys = rhs.fMaxToys;
388 fCalculator0 = rhs.fCalculator0;
389 fScannedVariable = rhs.fScannedVariable;
390 fUseCLs = rhs.fUseCLs;
391 fScanLog = rhs.fScanLog;
392 fSize = rhs.fSize;
393 fVerbose = rhs.fVerbose;
394 fCalcType = rhs.fCalcType;
395 fNBins = rhs.fNBins;
396 fXmin = rhs.fXmin;
397 fXmax = rhs.fXmax;
398 fNumErr = rhs.fNumErr;
399
400 return *this;
401}
402
403////////////////////////////////////////////////////////////////////////////////
404/// destructor (delete the HypoTestInverterResult)
405
407{
408 if (fResults) delete fResults;
409 fCalculator0 = nullptr;
410}
411
412////////////////////////////////////////////////////////////////////////////////
413/// return the test statistic which is or will be used by the class
414
416{
419 }
420 else
421 return nullptr;
422}
423
424////////////////////////////////////////////////////////////////////////////////
425/// set the test statistic to use
426
428{
431 return true;
432 }
433 else return false;
434}
435
436////////////////////////////////////////////////////////////////////////////////
437/// delete contained result and graph
438
440 if (fResults) delete fResults;
441 fResults = nullptr;
442 fLimitPlot.reset();
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// create a new HypoTestInverterResult to hold all computed results
447
449 if (fResults == nullptr) {
450 TString results_name = "result_";
453 TString title = "HypoTestInverter Result For ";
454 title += fScannedVariable->GetName();
455 fResults->SetTitle(title);
456 }
459 // check if one or two sided scan
460 if (fCalculator0) {
461 // if asymptotic calculator
463 if (ac) {
464 fResults->fIsTwoSided = ac->IsTwoSided();
465 } else {
466 // in case of the other calculators
468 if (sampler) {
469 ProfileLikelihoodTestStat * pl = dynamic_cast<ProfileLikelihoodTestStat*>(sampler->GetTestStatistic());
470 if (pl && pl->IsTwoSided() ) fResults->fIsTwoSided = true;
471 }
472 }
473 }
474}
475
476////////////////////////////////////////////////////////////////////////////////
477/// Run a fixed scan or the automatic scan depending on the configuration.
478/// Return if needed a copy of the result object which will be managed by the user.
479
481
482 // if having a result with at least one point return it
483 if (fResults && fResults->ArraySize() >= 1) {
484 oocoutI(nullptr,Eval) << "HypoTestInverter::GetInterval - return an already existing interval " << std::endl;
485 return static_cast<HypoTestInverterResult*>(fResults->Clone());
486 }
487
488 if (fNBins > 0) {
489 oocoutI(nullptr,Eval) << "HypoTestInverter::GetInterval - run a fixed scan" << std::endl;
491 if (!ret)
492 oocoutE(nullptr,Eval) << "HypoTestInverter::GetInterval - error running a fixed scan " << std::endl;
493 }
494 else {
495 oocoutI(nullptr,Eval) << "HypoTestInverter::GetInterval - run an automatic scan" << std::endl;
496 double limit(0);
497 double err(0);
498 bool ret = RunLimit(limit,err);
499 if (!ret)
500 oocoutE(nullptr,Eval) << "HypoTestInverter::GetInterval - error running an auto scan " << std::endl;
501 }
502
503 return static_cast<HypoTestInverterResult*> (fResults->Clone());
504}
505
506////////////////////////////////////////////////////////////////////////////////
507/// Run the Hypothesis test at a previous configured point
508/// (internal function called by RunOnePoint)
509
511 //for debug
512 // std::cout << ">>>>>>>>>>> " << std::endl;
513 // std::cout << "alternate model " << std::endl;
514 // hc.GetAlternateModel()->GetNuisanceParameters()->Print("V");
515 // hc.GetAlternateModel()->GetParametersOfInterest()->Print("V");
516 // std::cout << "Null model " << std::endl;
517 // hc.GetNullModel()->GetNuisanceParameters()->Print("V");
518 // hc.GetNullModel()->GetParametersOfInterest()->Print("V");
519 // std::cout << "<<<<<<<<<<<<<<< " << std::endl;
520
521 // run the hypothesis test
522 HypoTestResult * hcResult = hc.GetHypoTest();
523 if (hcResult == nullptr) {
524 oocoutE(nullptr,Eval) << "HypoTestInverter::Eval - HypoTest failed" << std::endl;
525 return hcResult;
526 }
527
528 // since the b model is the alt need to set the flag
529 hcResult->SetBackgroundAsAlt(true);
530
531
532 // bool flipPvalue = false;
533 // if (flipPValues)
534 // hcResult->SetPValueIsRightTail(!hcResult->GetPValueIsRightTail());
535
536 // adjust for some numerical error in discrete models and == is not anymore
537 if (hcResult->GetPValueIsRightTail()) {
538 hcResult->SetTestStatisticData(hcResult->GetTestStatisticData()-fNumErr); // issue with < vs <= in discrete models
539 } else {
540 hcResult->SetTestStatisticData(hcResult->GetTestStatisticData() +
541 fNumErr); // issue with < vs <= in discrete models
542 }
543
544 double clsMid = (fUseCLs ? hcResult->CLs() : hcResult->CLsplusb());
545 double clsMidErr = (fUseCLs ? hcResult->CLsError() : hcResult->CLsplusbError());
546
547 //if (fVerbose) std::cout << (fUseCLs ? "\tCLs = " : "\tCLsplusb = ") << clsMid << " +/- " << clsMidErr << std::endl;
548
549 if (adaptive) {
550
553
554 while (clsMidErr >= fgCLAccuracy && (clsTarget == -1 || std::abs(clsMid-clsTarget) < 3*clsMidErr) ) {
555 std::unique_ptr<HypoTestResult> more(hc.GetHypoTest());
556
557 // if (flipPValues)
558 // more->SetPValueIsRightTail(!more->GetPValueIsRightTail());
559
560 hcResult->Append(more.get());
561 clsMid = (fUseCLs ? hcResult->CLs() : hcResult->CLsplusb());
562 clsMidErr = (fUseCLs ? hcResult->CLsError() : hcResult->CLsplusbError());
563 if (fVerbose) std::cout << (fUseCLs ? "\tCLs = " : "\tCLsplusb = ") << clsMid << " +/- " << clsMidErr << std::endl;
564 }
565
566 }
567 if (fVerbose ) {
568 oocoutP(nullptr,Eval) << "P values for " << fScannedVariable->GetName() << " = " <<
569 fScannedVariable->getVal() << "\n" <<
570 "\tCLs = " << hcResult->CLs() << " +/- " << hcResult->CLsError() << "\n" <<
571 "\tCLb = " << hcResult->CLb() << " +/- " << hcResult->CLbError() << "\n" <<
572 "\tCLsplusb = " << hcResult->CLsplusb() << " +/- " << hcResult->CLsplusbError() << "\n" <<
573 std::endl;
574 }
575
577 fTotalToysRun += (hcResult->GetAltDistribution()->GetSize() + hcResult->GetNullDistribution()->GetSize());
578
579 // set sampling distribution name
580 TString nullDistName = TString::Format("%s_%s_%4.2f",hcResult->GetNullDistribution()->GetName(),
582 TString altDistName = TString::Format("%s_%s_%4.2f",hcResult->GetAltDistribution()->GetName(),
584
585 hcResult->GetNullDistribution()->SetName(nullDistName);
586 hcResult->GetAltDistribution()->SetName(altDistName);
587 }
588
589 return hcResult;
590}
591
592////////////////////////////////////////////////////////////////////////////////
593/// Run a Fixed scan in npoints between min and max
594
595bool HypoTestInverter::RunFixedScan( int nBins, double xMin, double xMax, bool scanLog ) const
596{
597
599 // interpolate the limits
602
603 // safety checks
604 if ( nBins<=0 ) {
605 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - Please provide nBins>0\n";
606 return false;
607 }
608 if ( nBins==1 && xMin!=xMax ) {
609 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - nBins==1 -> I will run for xMin (" << xMin << ")\n";
610 }
611 if ( xMin==xMax && nBins>1 ) {
612 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - xMin==xMax -> I will enforce nBins==1\n";
613 nBins = 1;
614 }
615 if ( xMin>xMax ) {
616 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - Please provide xMin ("
617 << xMin << ") smaller than xMax (" << xMax << ")\n";
618 return false;
619 }
620
621 if (xMin < fScannedVariable->getMin()) {
623 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - xMin < lower bound, using xmin = "
624 << xMin << std::endl;
625 }
626 if (xMax > fScannedVariable->getMax()) {
628 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RunFixedScan - xMax > upper bound, using xmax = "
629 << xMax << std::endl;
630 }
631
632 if (xMin <= 0. && scanLog) {
633 oocoutE(nullptr, InputArguments) << "HypoTestInverter::RunFixedScan - cannot go in log steps if xMin <= 0" << std::endl;
634 return false;
635 }
636
637 double thisX = xMin;
638 for (int i=0; i<nBins; i++) {
639
640 if (i > 0) { // avoids case of nBins = 1
641 if (scanLog) {
642 thisX = exp( log(xMin) + i*(log(xMax)-log(xMin))/(nBins-1) ); // scan in log x
643 } else {
644 thisX = xMin + i * (xMax - xMin) / (nBins - 1); // linear scan in x
645 }
646 }
647
648 const bool status = RunOnePoint(thisX);
649
650 // check if failed status
651 if ( status==false ) {
652 oocoutW(nullptr,Eval) << "HypoTestInverter::RunFixedScan - The hypo test for point " << thisX << " failed. Skipping." << std::endl;
653 }
654 }
655
656 return true;
657}
658
659////////////////////////////////////////////////////////////////////////////////
660/// run only one point at the given POI value
661
662bool HypoTestInverter::RunOnePoint( double rVal, bool adaptive, double clTarget) const
663{
664
666
667 // check if rVal is in the range specified for fScannedVariable
668 if ( rVal < fScannedVariable->getMin() ) {
669 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RunOnePoint - Out of range: using the lower bound "
671 << " on the scanned variable rather than " << rVal<< "\n";
673 }
674 if ( rVal > fScannedVariable->getMax() ) {
675 // print a message when you have a significative difference since rval is computed
676 if (rVal > fScannedVariable->getMax() * (1. + 1.E-12)) {
677 oocoutE(nullptr, InputArguments) << "HypoTestInverter::RunOnePoint - Out of range: using the upper bound "
678 << fScannedVariable->getMax() << " on the scanned variable rather than "
679 << rVal << "\n";
680 }
682 }
683
684 // save old value
685 double oldValue = fScannedVariable->getVal();
686
687 // evaluate hybrid calculator at a single point
689 // need to set value of rval in hybridcalculator
690 // assume null model is S+B and alternate is B only
692 RooArgSet poi; poi.add(*sbModel->GetParametersOfInterest());
693 // set poi to right values
695 const_cast<ModelConfig*>(sbModel)->SetSnapshot(poi);
696
697 if (fVerbose > 0)
698 oocoutP(nullptr,Eval) << "Running for " << fScannedVariable->GetName() << " = " << fScannedVariable->getVal() << std::endl;
699
700 // compute the results
701 std::unique_ptr<HypoTestResult> result( Eval(*fCalculator0,adaptive,clTarget) );
702 if (!result) {
703 oocoutE(nullptr,Eval) << "HypoTestInverter - Error running point " << fScannedVariable->GetName() << " = " <<
704 fScannedVariable->getVal() << std::endl;
705 return false;
706 }
707 // in case of a dummy result
708 const double nullPV = result->NullPValue();
709 const double altPV = result->AlternatePValue();
710 if (!std::isfinite(nullPV) || nullPV < 0. || nullPV > 1. || !std::isfinite(altPV) || altPV < 0. || altPV > 1.) {
711 oocoutW(nullptr,Eval) << "HypoTestInverter - Skipping invalid result for point " << fScannedVariable->GetName() << " = " <<
712 fScannedVariable->getVal() << ". null p-value=" << nullPV << ", alternate p-value=" << altPV << std::endl;
713 return false;
714 }
715
716 double lastXtested;
718 else lastXtested = -999;
719
720 if ( (std::abs(rVal) < 1 && TMath::AreEqualAbs(rVal, lastXtested,1.E-12) ) ||
721 (std::abs(rVal) >= 1 && TMath::AreEqualRel(rVal, lastXtested,1.E-12) ) ) {
722
723 oocoutI(nullptr,Eval) << "HypoTestInverter::RunOnePoint - Merge with previous result for "
724 << fScannedVariable->GetName() << " = " << rVal << std::endl;
726 if (prevResult && prevResult->GetNullDistribution() && prevResult->GetAltDistribution()) {
727 prevResult->Append(result.get());
728 }
729 else {
730 // if it was empty we re-use it
731 oocoutI(nullptr,Eval) << "HypoTestInverter::RunOnePoint - replace previous empty result\n";
733 delete oldObj;
734
735 fResults->fYObjects.Add(result.release());
736 }
737
738 } else {
739
740 // fill the results in the HypoTestInverterResult array
741 fResults->fXValues.push_back(rVal);
742 fResults->fYObjects.Add(result.release());
743
744 }
745
747
748 return true;
749}
750
751////////////////////////////////////////////////////////////////////////////////
752/// Run an automatic scan until the desired accuracy is reached.
753/// Start by default from the full interval (min,max) of the POI and then via bisection find the line crossing
754/// the target line.
755/// Optionally, a hint can be provided and the scan will be done closer to that value.
756/// If by bisection the desired accuracy will not be reached, a fit to the points is performed.
757/// \param[out] limit The limit.
758/// \param[out] limitErr The error of the limit.
759/// \param[in] absAccuracy Desired absolute accuracy.
760/// \param[in] relAccuracy Desired relative accuracy.
761/// \param[in] hint Hint to start from or nullptr for no hint.
762
763bool HypoTestInverter::RunLimit(double &limit, double &limitErr, double absAccuracy, double relAccuracy, const double*hint) const {
764
765
766 // routine from G. Petrucciani (from HiggsCombination CMS package)
767
769
770 if ((hint != nullptr) && (*hint > r->getMin())) {
771 r->setMax(std::min<double>(3.0 * (*hint), r->getMax()));
772 r->setMin(std::max<double>(0.3 * (*hint), r->getMin()));
773 oocoutI(nullptr,InputArguments) << "HypoTestInverter::RunLimit - Use hint value " << *hint
774 << " search in interval " << r->getMin() << " , " << r->getMax() << std::endl;
775 }
776
777 // if not specified use the default values for rel and absolute accuracy
780
781 typedef std::pair<double,double> CLs_t;
782 double clsTarget = fSize;
783 CLs_t clsMin(1, 0);
784 CLs_t clsMax(0, 0);
785 CLs_t clsMid(0, 0);
786 double rMin = r->getMin();
787 double rMax = r->getMax();
788 limit = 0.5*(rMax + rMin);
789 limitErr = 0.5*(rMax - rMin);
790 bool done = false;
791
792 TF1 expoFit("expoFit","[0]*exp([1]*(x-[2]))", rMin, rMax);
793
794 fLimitPlot = std::make_unique<TGraphErrors>();
795
796 if (fVerbose > 0) std::cout << "Search for upper limit to the limit" << std::endl;
797 for (int tries = 0; tries < 6; ++tries) {
798 if (! RunOnePoint(rMax) ) {
799 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypotest failed at upper limit of scan range: " << rMax << std::endl;
800 rMax *= 0.95;
801 continue;
802 }
803 clsMax = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
804 if (clsMax.first == 0 || clsMax.first + 3 * std::abs(clsMax.second) < clsTarget ) break;
805 rMax += rMax;
806 if (tries == 5) {
807 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Cannot determine upper limit of scan range. At " << r->GetName()
808 << " = " << rMax << " still getting "
809 << (fUseCLs ? "CLs" : "CLsplusb") << " = " << clsMax.first << std::endl;
810 return false;
811 }
812 }
813 if (fVerbose > 0) {
814 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Search for lower limit to the limit" << std::endl;
815 }
816
817 if ( fUseCLs && rMin == 0 ) {
818 clsMin = CLs_t(1,0);
819 }
820 else {
821 if (! RunOnePoint(rMin) ) {
822 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypotest failed at lower limit of scan range: " << rMin << std::endl;
823 return false;
824 }
825 clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
826 }
827 if (clsMin.first != 1 && clsMin.first - 3 * std::abs(clsMin.second) < clsTarget) {
828 if (fUseCLs) {
829 rMin = 0;
830 clsMin = CLs_t(1,0); // this is always true for CLs
831 } else {
832 rMin = -rMax / 4;
833 for (int tries = 0; tries < 6; ++tries) {
834 if (! RunOnePoint(rMin) ) {
835 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Hypotest failed at lower limit of scan range: " << rMin << std::endl;
836 rMin = rMin == 0. ? 0.1 : rMin * 1.1;
837 continue;
838 }
839 clsMin = std::make_pair( fResults->GetLastYValue(), fResults->GetLastYError() );
840 if (clsMin.first == 1 || clsMin.first - 3 * std::abs(clsMin.second) > clsTarget) break;
841 rMin += rMin;
842 if (tries == 5) {
843 oocoutE(nullptr,Eval) << "HypoTestInverter::RunLimit - Cannot determine lower limit of scan range. At " << r->GetName()
844 << " = " << rMin << " still get " << (fUseCLs ? "CLs" : "CLsplusb")
845 << " = " << clsMin.first << std::endl;
846 return false;
847 }
848 }
849 }
850 }
851
852 if (fVerbose > 0)
853 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Now doing proper bracketing & bisection" << std::endl;
854 do {
855
856 // break loop in case max toys is reached
857 if (fMaxToys > 0 && fTotalToysRun > fMaxToys ) {
858 oocoutW(nullptr,Eval) << "HypoTestInverter::RunLimit - maximum number of toys reached " << std::endl;
859 done = false; break;
860 }
861
862
863 // determine point by bisection or interpolation
864 limit = 0.5*(rMin+rMax); limitErr = 0.5*(rMax-rMin);
865 if (fgAlgo == "logSecant" && clsMax.first != 0) {
866 double logMin = log(clsMin.first);
867 double logMax = log(clsMax.first);
868 double logTarget = log(clsTarget);
869 limit = rMin + (rMax-rMin) * (logTarget - logMin)/(logMax - logMin);
870 if (clsMax.second != 0 && clsMin.second != 0) {
871 limitErr = hypot((logTarget-logMax) * (clsMin.second/clsMin.first), (logTarget-logMin) * (clsMax.second/clsMax.first));
873 }
874 }
875 r->setError(limitErr);
876
877 // exit if reached accuracy on r
878 if (limitErr < std::max(absAccuracy, relAccuracy * limit)) {
879 if (fVerbose > 1) {
880 oocoutI(nullptr, Eval) << "HypoTestInverter::RunLimit - reached accuracy " << limitErr << " below "
881 << std::max(absAccuracy, relAccuracy * limit) << std::endl;
882 }
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;
908 double rMaxBound = rMax;
909 // try to reduce the size of the interval
910 while (clsMin.second == 0 || std::abs(rMin-limit) > std::max(absAccuracy, relAccuracy * limit)) {
911 rMin = 0.5*(rMin+limit);
912 if (!RunOnePoint(rMin,true, clsTarget) ) {
913 oocoutE(nullptr,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 (std::abs(clsMin.first-clsTarget) <= 2*clsMin.second) break;
918 rMinBound = rMin;
919 }
920 while (clsMax.second == 0 || std::abs(rMax-limit) > std::max(absAccuracy, relAccuracy * limit)) {
921 rMax = 0.5*(rMax+limit);
922 if (!RunOnePoint(rMax,true,clsTarget) ) {
923 oocoutE(nullptr,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 (std::abs(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(nullptr,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;
945 double rMaxBound;
946 expoFit.GetRange(rMinBound, rMaxBound);
947 limitErr = std::max(std::abs(rMinBound-limit), std::abs(rMaxBound-limit));
948 int npoints = 0;
949
950 HypoTestInverterPlot plot("plot","plot",fResults);
951 fLimitPlot.reset(plot.MakePlot() );
952
953
954 for (int j = 0; j < fLimitPlot->GetN(); ++j) {
955 if (fLimitPlot->GetX()[j] >= rMinBound && fLimitPlot->GetX()[j] <= rMaxBound) npoints++;
956 }
957 for (int i = 0, imax = /*(readHybridResults_ ? 0 : */ 8; i <= imax; ++i, ++npoints) {
958 fLimitPlot->Sort();
959 fLimitPlot->Fit(&expoFit,(fVerbose <= 1 ? "QNR EX0" : "NR EXO"));
960 if (fVerbose) {
961 oocoutI(nullptr,Eval) << "Fit to " << npoints << " points: " << expoFit.GetParameter(2) << " +/- " << expoFit.GetParError(2) << std::endl;
962 }
963 if ((rMin < expoFit.GetParameter(2)) && (expoFit.GetParameter(2) < rMax) && (expoFit.GetParError(2) < 0.5*(rMaxBound-rMinBound))) {
964 // sanity check fit result
965 limit = expoFit.GetParameter(2);
966 limitErr = expoFit.GetParError(2);
967 if (limitErr < std::max(absAccuracy, relAccuracy * limit)) break;
968 }
969 // add one point in the interval.
971 if (i != imax) {
972 if (!RunOnePoint(rTry,true,clsTarget) ) return false;
973 //eval(w, mc_s, mc_b, data, rTry, true, clsTarget);
974 }
975
976 }
977 }
978
979//if (!plot_.empty() && fLimitPlot.get()) {
980 if (fLimitPlot.get() && fLimitPlot->GetN() > 0) {
981 //new TCanvas("c1","c1");
982 fLimitPlot->Sort();
983 fLimitPlot->SetLineWidth(2);
984 double xmin = r->getMin();
985 double xmax = r->getMax();
986 for (int j = 0; j < fLimitPlot->GetN(); ++j) {
987 if (fLimitPlot->GetY()[j] > 1.4*clsTarget || fLimitPlot->GetY()[j] < 0.6*clsTarget) continue;
988 xmin = std::min(fLimitPlot->GetX()[j], xmin);
989 xmax = std::max(fLimitPlot->GetX()[j], xmax);
990 }
991 fLimitPlot->GetXaxis()->SetRangeUser(xmin,xmax);
992 fLimitPlot->GetYaxis()->SetRangeUser(0.5*clsTarget, 1.5*clsTarget);
993 fLimitPlot->Draw("AP");
994 expoFit.Draw("SAME");
995 TLine line(fLimitPlot->GetX()[0], clsTarget, fLimitPlot->GetX()[fLimitPlot->GetN()-1], clsTarget);
997 line.DrawLine(limit, 0, limit, fLimitPlot->GetY()[0]);
999 line.DrawLine(limit-limitErr, 0, limit-limitErr, fLimitPlot->GetY()[0]);
1000 line.DrawLine(limit+limitErr, 0, limit+limitErr, fLimitPlot->GetY()[0]);
1001 //c1->Print(plot_.c_str());
1002 }
1003
1004 oocoutI(nullptr,Eval) << "HypoTestInverter::RunLimit - Result: \n"
1005 << "\tLimit: " << r->GetName() << " < " << limit << " +/- " << limitErr << " @ " << (1-fSize) * 100 << "% CL\n";
1006 if (fVerbose > 1) oocoutI(nullptr,Eval) << "Total toys: " << fTotalToysRun << std::endl;
1007
1008 // set value in results
1009 fResults->fUpperLimit = limit;
1012 // lower limit are always min of p value
1016
1017 return true;
1018}
1019
1020////////////////////////////////////////////////////////////////////////////////
1021/// get the distribution of lower limit
1022/// if rebuild = false (default) it will re-use the results of the scan done
1023/// for obtained the observed limit and no extra toys will be generated
1024/// if rebuild a new set of B toys will be done and the procedure will be repeated
1025/// for each toy
1026
1028 if (!rebuild) {
1029 if (!fResults) {
1030 oocoutE(nullptr,InputArguments) << "HypoTestInverter::GetLowerLimitDistribution(false) - result not existing\n";
1031 return nullptr;
1032 }
1034 }
1035
1036 TList * clsDist = nullptr;
1037 TList * clsbDist = nullptr;
1039 else clsbDist = &fResults->fExpPValues;
1040
1042
1043}
1044
1045////////////////////////////////////////////////////////////////////////////////
1046/// get the distribution of lower limit
1047/// if rebuild = false (default) it will re-use the results of the scan done
1048/// for obtained the observed limit and no extra toys will be generated
1049/// if rebuild a new set of B toys will be done and the procedure will be repeated
1050/// for each toy
1051/// The nuisance parameter value used for rebuild is the current one in the model
1052/// so it is user responsibility to set to the desired value (nomi
1053
1055 if (!rebuild) {
1056 if (!fResults) {
1057 oocoutE(nullptr,InputArguments) << "HypoTestInverter::GetUpperLimitDistribution(false) - result not existing\n";
1058 return nullptr;
1059 }
1061 }
1062
1063 TList * clsDist = nullptr;
1064 TList * clsbDist = nullptr;
1066 else clsbDist = &fResults->fExpPValues;
1067
1069}
1070
1071////////////////////////////////////////////////////////////////////////////////
1072
1076
1077////////////////////////////////////////////////////////////////////////////////
1078/// rebuild the sampling distributions by
1079/// generating some toys and find for each of them a new upper limit
1080/// Return the upper limit distribution and optionally also the pValue distributions for Cls, Clsb and Clbxs
1081/// as a TList for each scanned point
1082/// The method uses the present parameter value. It is user responsibility to give the current parameters to rebuild the distributions
1083/// 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
1084/// output file as an histogram
1085
1087
1088 if (!fScannedVariable || !fCalculator0) return nullptr;
1089 // get first background snapshot
1092 if (!bModel || ! sbModel) return nullptr;
1094 if (!sbModel->GetParametersOfInterest()) return nullptr;
1095 paramPoint.add(*sbModel->GetParametersOfInterest());
1096
1097 const RooArgSet * poibkg = bModel->GetSnapshot();
1098 if (!poibkg) {
1099 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RebuildDistribution - background snapshot not existing"
1100 << " assume is for POI = 0" << std::endl;
1103 }
1104 else
1105 paramPoint.assign(*poibkg);
1106 // generate data at bkg parameter point
1107
1109 if (!toymcSampler) {
1110 oocoutE(nullptr,InputArguments) << "HypoTestInverter::RebuildDistribution - no toy MC sampler existing" << std::endl;
1111 return nullptr;
1112 }
1113 // set up test stat sampler in case of asymptotic calculator
1114 if (dynamic_cast<RooStats::AsymptoticCalculator*>(fCalculator0) ) {
1115 toymcSampler->SetObservables(*sbModel->GetObservables() );
1116 toymcSampler->SetParametersForTestStat(*sbModel->GetParametersOfInterest());
1117 toymcSampler->SetPdf(*sbModel->GetPdf());
1118 toymcSampler->SetNuisanceParameters(*sbModel->GetNuisanceParameters());
1119 if (sbModel->GetGlobalObservables() ) toymcSampler->SetGlobalObservables(*sbModel->GetGlobalObservables() );
1120 // set number of events
1121 if (!sbModel->GetPdf()->canBeExtended())
1122 toymcSampler->SetNEventsPerToy(1);
1123 }
1124
1125 // loop on data to generate
1126 int nPoints = fNBins;
1127
1128 bool storePValues = clsDist || clsbDist || clbDist;
1129 if (fNBins <=0 && storePValues) {
1130 oocoutW(nullptr,InputArguments) << "HypoTestInverter::RebuildDistribution - cannot return p values distribution with the auto scan" << std::endl;
1131 storePValues = false;
1132 nPoints = 0;
1133 }
1134
1135 if (storePValues) {
1137 if (nPoints <=0) {
1138 oocoutE(nullptr,InputArguments) << "HypoTestInverter - result is not existing and number of point to scan is not set"
1139 << std::endl;
1140 return nullptr;
1141 }
1142 }
1143
1144 if (nToys <= 0) nToys = 100; // default value
1145
1146 std::vector<std::vector<double> > CLs_values(nPoints);
1147 std::vector<std::vector<double> > CLsb_values(nPoints);
1148 std::vector<std::vector<double> > CLb_values(nPoints);
1149
1150 if (storePValues) {
1151 for (int i = 0; i < nPoints; ++i) {
1152 CLs_values[i].reserve(nToys);
1153 CLb_values[i].reserve(nToys);
1154 CLsb_values[i].reserve(nToys);
1155 }
1156 }
1157
1158 std::vector<double> limit_values; limit_values.reserve(nToys);
1159
1160 oocoutI(nullptr,InputArguments) << "HypoTestInverter - rebuilding the p value distributions by generating ntoys = "
1161 << nToys << std::endl;
1162
1163
1164 oocoutI(nullptr,InputArguments) << "Rebuilding using parameter of interest point: ";
1165 RooStats::PrintListContent(paramPoint, oocoutI(nullptr,InputArguments) );
1166 if (sbModel->GetNuisanceParameters() ) {
1167 oocoutI(nullptr,InputArguments) << "And using nuisance parameters: ";
1168 RooStats::PrintListContent(*sbModel->GetNuisanceParameters(), oocoutI(nullptr,InputArguments) );
1169 }
1170 // save all parameters to restore them later
1171 assert(bModel->GetPdf() );
1172 assert(bModel->GetObservables() );
1173 std::unique_ptr<RooArgSet> allParams{bModel->GetPdf()->getParameters( *bModel->GetObservables() )};
1175 allParams->snapshot(saveParams);
1176
1177 std::unique_ptr<TFile> fileOut{TFile::Open(outputfile,"RECREATE")};
1178 if (!fileOut) {
1179 oocoutE(nullptr,InputArguments) << "HypoTestInverter - RebuildDistributions - Error opening file " << outputfile
1180 << " - the resulting limits will not be stored" << std::endl;
1181 }
1182 // create temporary histograms to store the limit result
1183 TH1D * hL = new TH1D("lowerLimitDist","Rebuilt lower limit distribution",100,1.,0.);
1184 TH1D * hU = new TH1D("upperLimitDist","Rebuilt upper limit distribution",100,1.,0.);
1185 TH1D * hN = new TH1D("nObs","Observed events",100,1.,0.);
1186 hL->SetBuffer(2*nToys);
1187 hU->SetBuffer(2*nToys);
1188 std::vector<TH1*> hCLb;
1189 std::vector<TH1*> hCLsb;
1190 std::vector<TH1*> hCLs;
1191 if (storePValues) {
1192 for (int i = 0; i < nPoints; ++i) {
1193 hCLb.push_back(new TH1D(TString::Format("CLbDist_bin%d",i),"CLb distribution",100,1.,0.));
1194 hCLs.push_back(new TH1D(TString::Format("ClsDist_bin%d",i),"CLs distribution",100,1.,0.));
1195 hCLsb.push_back(new TH1D(TString::Format("CLsbDist_bin%d",i),"CLs+b distribution",100,1.,0.));
1196 }
1197 }
1198
1199
1200 // loop now on the toys
1201 for (int itoy = 0; itoy < nToys; ++itoy) {
1202
1203 oocoutP(nullptr,Eval) << "\nHypoTestInverter - RebuildDistributions - running toy # " << itoy << " / "
1204 << nToys << std::endl;
1205
1206
1207 std::cout << "\n\nshnapshot of s+b model \n";
1208 sbModel->GetSnapshot()->Print("v");
1209
1210 // reset parameters to initial values to be sure in case they are not reset
1211 if (itoy> 0) {
1212 allParams->assign(saveParams);
1213 }
1214
1215 // need to set th epdf to clear the cache in ToyMCSampler
1216 // pdf we must use is background pdf
1217 toymcSampler->SetPdf(*bModel->GetPdf() );
1218
1219
1220 RooAbsData * bkgdata = toymcSampler->GenerateToyData(paramPoint);
1221
1222 double nObs = bkgdata->sumEntries();
1223 // for debugging in case of number counting models
1224 if (bkgdata->numEntries() ==1 && !bModel->GetPdf()->canBeExtended()) {
1225 oocoutP(nullptr,Generation) << "Generate observables are : ";
1226 RooArgList genObs(*bkgdata->get(0));
1227 RooStats::PrintListContent(genObs, oocoutP(nullptr,Generation) );
1228 nObs = 0;
1229 for (std::size_t i = 0; i < genObs.size(); ++i) {
1230 RooRealVar * x = dynamic_cast<RooRealVar*>(&genObs[i]);
1231 if (x) nObs += x->getVal();
1232 }
1233 }
1234 hN->Fill(nObs);
1235
1236 // by copying I will have the same min/max as previous ones
1237 HypoTestInverter inverter = *this;
1238 inverter.SetData(*bkgdata);
1239
1240 // print global observables
1241 auto gobs = bModel->GetPdf()->getVariables()->selectCommon(* sbModel->GetGlobalObservables() );
1242 gobs->Print("v");
1243
1244 HypoTestInverterResult * r = inverter.GetInterval();
1245
1246 if (r == nullptr) continue;
1247
1248 double value = (isUpper) ? r->UpperLimit() : r->LowerLimit();
1249 limit_values.push_back( value );
1250 hU->Fill(r->UpperLimit() );
1251 hL->Fill(r->LowerLimit() );
1252
1253
1254 std::cout << "The computed upper limit for toy #" << itoy << " is " << value << std::endl;
1255
1256 // write every 10 toys
1257 if (itoy%10 == 0 || itoy == nToys-1) {
1258 hU->Write("",TObject::kOverwrite);
1259 hL->Write("",TObject::kOverwrite);
1260 hN->Write("",TObject::kOverwrite);
1261 }
1262
1263 if (!storePValues) continue;
1264
1265 if (nPoints < r->ArraySize()) {
1266 oocoutW(nullptr,InputArguments) << "HypoTestInverter: skip extra points" << std::endl;
1267 }
1268 else if (nPoints > r->ArraySize()) {
1269 oocoutW(nullptr,InputArguments) << "HypoTestInverter: missing some points" << std::endl;
1270 }
1271
1272
1273 for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1274 HypoTestResult * hr = r->GetResult(ipoint);
1275 if (hr) {
1276 CLs_values[ipoint].push_back( hr->CLs() );
1277 CLsb_values[ipoint].push_back( hr->CLsplusb() );
1278 CLb_values[ipoint].push_back( hr->CLb() );
1279 hCLs[ipoint]->Fill( hr->CLs() );
1280 hCLb[ipoint]->Fill( hr->CLb() );
1281 hCLsb[ipoint]->Fill( hr->CLsplusb() );
1282 }
1283 else {
1284 oocoutW(nullptr,InputArguments) << "HypoTestInverter: missing result for point: x = "
1285 << fResults->GetXValue(ipoint) << std::endl;
1286 }
1287 }
1288 // write every 10 toys
1289 if (itoy%10 == 0 || itoy == nToys-1) {
1290 for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1291 hCLs[ipoint]->Write("",TObject::kOverwrite);
1292 hCLb[ipoint]->Write("",TObject::kOverwrite);
1293 hCLsb[ipoint]->Write("",TObject::kOverwrite);
1294 }
1295 }
1296
1297
1298 delete r;
1299 delete bkgdata;
1300 }
1301
1302
1303 if (storePValues) {
1304 if (clsDist) clsDist->SetOwner(true);
1305 if (clbDist) clbDist->SetOwner(true);
1306 if (clsbDist) clsbDist->SetOwner(true);
1307
1308 oocoutI(nullptr,InputArguments) << "HypoTestInverter: storing rebuilt p values " << std::endl;
1309
1310 for (int ipoint = 0; ipoint < nPoints; ++ipoint) {
1311 if (clsDist) {
1312 TString name = TString::Format("CLs_distrib_%d",ipoint);
1314 }
1315 if (clbDist) {
1316 TString name = TString::Format("CLb_distrib_%d",ipoint);
1318 }
1319 if (clsbDist) {
1320 TString name = TString::Format("CLsb_distrib_%d",ipoint);
1322 }
1323 }
1324 }
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";
1339}
dim_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)
@ kRed
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
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
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
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:103
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:24
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:77
Variable that can be changed from the outside.
Definition RooRealVar.h:37
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.
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.
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.
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
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
This class simply holds a sampling distribution of some test statistic.
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.
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
1-Dim function class
Definition TF1.h:234
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:4131
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:698
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:103
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:820
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:174
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
@ kOverwrite
overwrite existing object with same name
Definition TObject.h:98
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
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:2378
TLine * line
Double_t x[n]
Definition legend1.C:17
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
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:426
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)