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