Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HypoTestInverterResult.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11
12/** \class RooStats::HypoTestInverterResult
13 \ingroup Roostats
14
15HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence interval.
16Based on the RatioFinder code available in the RooStatsCms package developed by Gregory Schott and Danilo Piparo
17Ported and adapted to RooStats by Gregory Schott
18Some contributions to this class have been written by Matthias Wolf (error estimation)
19
20*/
21
22// include header file of this class
24
27#include "RooMsgService.h"
28#include "RooGlobalFunc.h"
29
30#include "TF1.h"
31#include "TGraphErrors.h"
32#include <cmath>
35#include "Math/Functor.h"
36
37#include "TCanvas.h"
38#include "TFile.h"
40
41#include <algorithm>
42
43
44using namespace RooStats;
45using namespace RooFit;
46
47// initialize static value
50
51////////////////////////////////////////////////////////////////////////////////
52/// default constructor
53
56 fUseCLs(false),
57 fIsTwoSided(false),
58 fInterpolateLowerLimit(true),
59 fInterpolateUpperLimit(true),
60 fFittedLowerLimit(false),
61 fFittedUpperLimit(false),
62 fInterpolOption(kLinear),
63 fLowerLimitError(-1),
64 fUpperLimitError(-1),
65 fCLsCleanupThreshold(0.005)
66{
69
72}
73
74////////////////////////////////////////////////////////////////////////////////
75/// copy constructor
76
79 fUseCLs(other.fUseCLs),
80 fIsTwoSided(other.fIsTwoSided),
81 fInterpolateLowerLimit(other.fInterpolateLowerLimit),
82 fInterpolateUpperLimit(other.fInterpolateUpperLimit),
83 fFittedLowerLimit(other.fFittedLowerLimit),
84 fFittedUpperLimit(other.fFittedUpperLimit),
85 fInterpolOption(other.fInterpolOption),
86 fLowerLimitError(other.fLowerLimitError),
87 fUpperLimitError(other.fUpperLimitError),
88 fCLsCleanupThreshold(other.fCLsCleanupThreshold),
89 fXValues(other.fXValues)
90{
93 int nOther = other.ArraySize();
94
95 for (int i = 0; i < nOther; ++i)
96 fYObjects.Add( other.fYObjects.At(i)->Clone() );
97 for (int i = 0; i < fExpPValues.GetSize() ; ++i)
98 fExpPValues.Add( other.fExpPValues.At(i)->Clone() );
99
102}
103
104////////////////////////////////////////////////////////////////////////////////
105
108{
109 if (&other==this) {
110 return *this ;
111 }
112
114 fLowerLimit = other.fLowerLimit;
115 fUpperLimit = other.fUpperLimit;
116 fUseCLs = other.fUseCLs;
117 fIsTwoSided = other.fIsTwoSided;
118 fInterpolateLowerLimit = other.fInterpolateLowerLimit;
119 fInterpolateUpperLimit = other.fInterpolateUpperLimit;
120 fFittedLowerLimit = other.fFittedLowerLimit;
121 fFittedUpperLimit = other.fFittedUpperLimit;
122 fInterpolOption = other.fInterpolOption;
123 fLowerLimitError = other.fLowerLimitError;
124 fUpperLimitError = other.fUpperLimitError;
125 fCLsCleanupThreshold = other.fCLsCleanupThreshold;
126
127 int nOther = other.ArraySize();
128 fXValues = other.fXValues;
129
131 for (int i=0; i < nOther; ++i) {
132 fYObjects.Add( other.fYObjects.At(i)->Clone() );
133 }
135 for (int i=0; i < fExpPValues.GetSize() ; ++i) {
136 fExpPValues.Add( other.fExpPValues.At(i)->Clone() );
137 }
138
141
142 return *this;
143}
144
145////////////////////////////////////////////////////////////////////////////////
146/// constructor
147
150 double cl ) :
151 SimpleInterval(name,scannedVariable,TMath::QuietNaN(),TMath::QuietNaN(),cl),
152 fUseCLs(false),
153 fIsTwoSided(false),
154 fInterpolateLowerLimit(true),
155 fInterpolateUpperLimit(true),
156 fFittedLowerLimit(false),
157 fFittedUpperLimit(false),
158 fInterpolOption(kLinear),
159 fLowerLimitError(-1),
160 fUpperLimitError(-1),
161 fCLsCleanupThreshold(0.005)
162{
165
166 // put a cloned copy of scanned variable to set in the interval
167 // to avoid I/O problem of the Result class -
168 // make the set owning the cloned copy (use clone instead of Clone to not copying all links)
170 fParameters.addOwned(std::unique_ptr<RooRealVar>{static_cast<RooRealVar *>(scannedVariable.clone(scannedVariable.GetName()))});
171}
172
173////////////////////////////////////////////////////////////////////////////////
174/// destructor
175
177{
178 // explicitly empty the TLists - these contain pointers, not objects
181
184}
185
186////////////////////////////////////////////////////////////////////////////////
187/// Remove problematic points from this result.
188///
189/// This function can be used to clean up a result that has failed fits, spiking CLs
190/// or similar problems. It removes
191/// - Points where CLs is not falling monotonously. These may result from a lack of numerical precision.
192/// - Points where CLs spikes to more than 0.999.
193/// - Points with very low CLs. These are not needed to run the inverter, which speeds up the process.
194/// - Points where CLs < 0. These occur when fits fail.
196{
197 const int nEntries = ArraySize();
198
199 // initialization
200 double nsig1(1.0);
201 double nsig2(2.0);
202 double p[5];
203 double q[5];
204
207 p[2] = 0.5;
210
211 bool resultIsAsymptotic(false);
212 if (nEntries>=1) {
213 HypoTestResult* r = dynamic_cast<HypoTestResult *> ( GetResult(0) );
214 assert(r!=nullptr);
215 if ( !r->GetNullDistribution() && !r->GetAltDistribution() ) {
216 resultIsAsymptotic = true;
217 }
218 }
219
220 int nPointsRemoved(0);
221
222 double CLsobsprev(1.0);
223
224 for (auto itr = fXValues.begin(); itr != fXValues.end(); ++itr) {
225 const double x = *itr;
226 const int i = FindIndex(x);
227
229 if (!s) break;
230
231 /////////////////////////////////////////////////////////////////////////////////////////
232
233 const std::vector<double> & values = s->GetSamplingDistribution();
234 if ((int) values.size() != fgAsymptoticNumPoints) {
235 coutE(Eval) << "HypoTestInverterResult::ExclusionCleanup - invalid size of sampling distribution" << std::endl;
236 delete s;
237 break;
238 }
239
240 /// expected p-values
241 // special case for asymptotic results (cannot use TMath::quantile in that case)
242 if (resultIsAsymptotic) {
244 double dsig = 2.*maxSigma / (values.size() -1) ;
245 int i0 = (int) TMath::Floor ( ( -nsig2 + maxSigma )/dsig + 0.5 );
246 int i1 = (int) TMath::Floor ( ( -nsig1 + maxSigma )/dsig + 0.5 );
247 int i2 = (int) TMath::Floor ( ( maxSigma )/dsig + 0.5 );
248 int i3 = (int) TMath::Floor ( ( nsig1 + maxSigma )/dsig + 0.5 );
249 int i4 = (int) TMath::Floor ( ( nsig2 + maxSigma )/dsig + 0.5 );
250 //
251 q[0] = values[i0];
252 q[1] = values[i1];
253 q[2] = values[i2];
254 q[3] = values[i3];
255 q[4] = values[i4];
256 } else {
257 double * z = const_cast<double *>(&values[0] ); // need to change TMath::Quantiles
258 TMath::Quantiles(values.size(), 5, z, q, p, false);
259 }
260
261 delete s;
262
263 const double CLsobs = CLs(i);
264
265 /////////////////////////////////////////////////////////////////////////////////////////
266
267 bool removeThisPoint(false);
268
269 // 1. CLs should drop, else skip this point
270 if (resultIsAsymptotic && i>=1 && CLsobs>CLsobsprev) {
271 removeThisPoint = true;
272 } else if (CLsobs >= 0.) {
274 }
275
276 // 2. CLs should not spike, else skip this point
277 removeThisPoint |= i>=1 && CLsobs >= 0.9999;
278
279 // 3. Not interested in CLs values that become too low.
280 removeThisPoint |= i>=1 && q[4] < fCLsCleanupThreshold;
281
282 // 4. Negative CLs indicate failed fits
283 removeThisPoint |= CLsobs < 0.;
284
285 // to remove or not to remove
286 if (removeThisPoint) {
287 itr = fXValues.erase(itr)--;
291 continue;
292 } else { // keep
294 }
295 }
296
297 // after cleanup, reset existing limits
298 fFittedUpperLimit = false;
299 fFittedLowerLimit = false;
301
302 return nPointsRemoved;
303}
304
305////////////////////////////////////////////////////////////////////////////////
306/// Merge this HypoTestInverterResult with another
307/// HypoTestInverterResult passed as argument
308/// The merge is done by combining the HypoTestResult when the same point value exist in both results.
309/// If results exist at different points these are added in the new result
310/// NOTE: Merging of the expected p-values obtained with pseudo-data.
311/// When expected p-values exist in the result (i.e. when rebuild option is used when getting the expected
312/// limit distribution in the HYpoTestInverter) then the expected p-values are also merged. This is equivalent
313/// at merging the pseudo-data. However there can be an inconsistency if the expected p-values have been
314/// obtained with different toys. In this case the merge is done but a warning message is printed.
315
317{
318 int nThis = ArraySize();
319 int nOther = otherResult.ArraySize();
320 if (nOther == 0) return true;
321 if (nOther != otherResult.fYObjects.GetSize() ) return false;
322 if (nThis != fYObjects.GetSize() ) return false;
323
324 // cannot merge in case of inconsistent members
325 if (fExpPValues.GetSize() > 0 && fExpPValues.GetSize() != nThis ) return false;
326 if (otherResult.fExpPValues.GetSize() > 0 && otherResult.fExpPValues.GetSize() != nOther ) return false;
327
328 coutI(Eval) << "HypoTestInverterResult::Add - merging result from " << otherResult.GetName()
329 << " in " << GetName() << std::endl;
330
331 bool addExpPValues = (fExpPValues.GetSize() == 0 && otherResult.fExpPValues.GetSize() > 0);
332 bool mergeExpPValues = (fExpPValues.GetSize() > 0 && otherResult.fExpPValues.GetSize() > 0);
333
335 coutI(Eval) << "HypoTestInverterResult::Add - merging also the expected p-values from pseudo-data" << std::endl;
336
337
338 // case current result is empty
339 // just make a simple copy of the other result
340 if (nThis == 0) {
341 fXValues = otherResult.fXValues;
342 for (int i = 0; i < nOther; ++i)
343 fYObjects.Add( otherResult.fYObjects.At(i)->Clone() );
344 for (int i = 0; i < fExpPValues.GetSize() ; ++i)
345 fExpPValues.Add( otherResult.fExpPValues.At(i)->Clone() );
346 }
347 // now do the real merge combining point with same value or adding extra ones
348 else {
349 for (int i = 0; i < nOther; ++i) {
350 double otherVal = otherResult.fXValues[i];
351 HypoTestResult * otherHTR = static_cast<HypoTestResult*>(otherResult.fYObjects.At(i));
352 if (otherHTR == nullptr) continue;
353 bool sameXFound = false;
354 for (int j = 0; j < nThis; ++j) {
355 double thisVal = fXValues[j];
356
357 // if same value merge the result
358 if ( (std::abs(otherVal) < 1 && TMath::AreEqualAbs(otherVal, thisVal,1.E-12) ) ||
359 (std::abs(otherVal) >= 1 && TMath::AreEqualRel(otherVal, thisVal,1.E-12) ) ) {
361 thisHTR->Append(otherHTR);
362 sameXFound = true;
363 if (mergeExpPValues) {
364 (static_cast<SamplingDistribution*>(fExpPValues.At(j)))->Add( static_cast<SamplingDistribution*>(otherResult.fExpPValues.At(i)) );
365 // check if same toys have been used for the test statistic distribution
366 int thisNToys = (thisHTR->GetNullDistribution() ) ? thisHTR->GetNullDistribution()->GetSize() : 0;
367 int otherNToys = (otherHTR->GetNullDistribution() ) ? otherHTR->GetNullDistribution()->GetSize() : 0;
368 if (thisNToys != otherNToys )
369 coutW(Eval) << "HypoTestInverterResult::Add expected p values have been generated with different toys " << thisNToys << " , " << otherNToys << std::endl;
370 }
371 break;
372 }
373 }
374 if (!sameXFound) {
375 // add the new result
376 fYObjects.Add(otherHTR->Clone() );
377 fXValues.push_back( otherVal);
378 }
379 // add in any case also when same x found
380 if (addExpPValues)
381 fExpPValues.Add( otherResult.fExpPValues.At(i)->Clone() );
382
383
384 }
385 }
386
387 if (ArraySize() > nThis) {
388 coutI(Eval) << "HypoTestInverterResult::Add - new number of points is " << fXValues.size()
389 << std::endl;
390 } else {
391 coutI(Eval) << "HypoTestInverterResult::Add - new toys/point is "
392 << (static_cast<HypoTestResult *>(fYObjects.At(0)))->GetNullDistribution()->GetSize() << std::endl;
393 }
394
395 // reset cached limit values
398
399 return true;
400}
401
402////////////////////////////////////////////////////////////////////////////////
403/// Add a single point result (an HypoTestResult)
404
406{
407 int i= FindIndex(x);
408 if (i<0) {
409 fXValues.push_back(x);
410 fYObjects.Add(res.Clone());
411 } else {
413 if (!r) return false;
414 r->Append(&res);
415 }
416
417 // reset cached limit values
420
421 return true;
422}
423
424////////////////////////////////////////////////////////////////////////////////
425/// function to return the value of the parameter of interest for the i^th entry in the results
426
428{
429 if ( index >= ArraySize() || index<0 ) {
430 coutE(InputArguments) << "Problem: You are asking for an impossible array index value\n";
431 return -999;
432 }
433
434 return fXValues[index];
435}
436
437////////////////////////////////////////////////////////////////////////////////
438/// function to return the value of the confidence level for the i^th entry in the results
439
441{
442 auto result = GetResult(index);
443 if ( !result ) {
444 return -999;
445 }
446
447 if (fUseCLs) {
448 return result->CLs();
449 } else {
450 return result->CLsplusb(); // CLs+b
451 }
452}
453
454////////////////////////////////////////////////////////////////////////////////
455/// function to return the estimated error on the value of the confidence level for the i^th entry in the results
456
458{
459 auto result = GetResult(index);
460 if ( !result ) {
461 return -999;
462 }
463
464 if (fUseCLs) {
465 return result->CLsError();
466 } else {
467 return result->CLsplusbError();
468 }
469}
470
471////////////////////////////////////////////////////////////////////////////////
472/// function to return the observed CLb value for the i-th entry
473
475{
476 auto result = GetResult(index);
477 if ( !result ) {
478 return -999;
479 }
480 return result->CLb();
481}
482
483////////////////////////////////////////////////////////////////////////////////
484/// function to return the observed CLs+b value for the i-th entry
485
487{
488 auto result = GetResult(index);
489 if ( !result) {
490 return -999;
491 }
492 return result->CLsplusb();
493}
494
495////////////////////////////////////////////////////////////////////////////////
496/// function to return the observed CLs value for the i-th entry
497
499{
500 auto result = GetResult(index);
501 if ( !result ) {
502 return -999;
503 }
504 return result->CLs();
505}
506
507////////////////////////////////////////////////////////////////////////////////
508/// function to return the error on the observed CLb value for the i-th entry
509
511{
512 auto result = GetResult(index);
513 if ( !result ) {
514 return -999;
515 }
516 return result->CLbError();
517}
518
519////////////////////////////////////////////////////////////////////////////////
520/// function to return the error on the observed CLs+b value for the i-th entry
521
523{
524 auto result = GetResult(index);
525 if ( ! result ) {
526 return -999;
527 }
528 return result->CLsplusbError();
529}
530
531////////////////////////////////////////////////////////////////////////////////
532/// function to return the error on the observed CLs value for the i-th entry
533
535{
536 auto result = GetResult(index);
537 if ( ! result ){
538 return -999;
539 }
540 return result->CLsError();
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// get the HypoTestResult object at the given index point
545
547{
548 if ( index >= ArraySize() || index<0 ) {
549 coutE(InputArguments) << "Problem: You are asking for an impossible array index value\n";
550 return nullptr;
551 }
552
553 return (static_cast<HypoTestResult*>(fYObjects.At(index)));
554}
555
556////////////////////////////////////////////////////////////////////////////////
557/// find the index corresponding at the poi value xvalue
558/// If no points is found return -1
559/// Note that a tolerance is used of 10^-12 to find the closest point
560
562{
563 const double tol = 1.E-12;
564 for (int i=0; i<ArraySize(); i++) {
565 double xpoint = fXValues[i];
566 if ( (std::abs(xvalue) > 1 && TMath::AreEqualRel( xvalue, xpoint, tol) ) ||
567 (std::abs(xvalue) <= 1 && TMath::AreEqualAbs( xvalue, xpoint, tol) ) )
568 return i;
569 }
570 return -1;
571}
572
573////////////////////////////////////////////////////////////////////////////////
574/// return the X value of the given graph for the target value y0
575/// the graph is evaluated using linear interpolation by default.
576/// if option = "S" a TSpline3 is used
577
578double HypoTestInverterResult::GetGraphX(const TGraph & graph, double y0, bool lowSearch, double &axmin, double &axmax) const {
579
580//#define DO_DEBUG
581#ifdef DO_DEBUG
582 std::cout << "using graph for search " << lowSearch << " min " << axmin << " max " << axmax << std::endl;
583#endif
584
585
586 // find reasonable xmin and xmax if not given
587 const double * y = graph.GetY();
588 int n = graph.GetN();
589 if (n < 2) {
590 ooccoutE(this,Eval) << "HypoTestInverterResult::GetGraphX - need at least 2 points for interpolation (n=" << n << ")\n";
591 return (n>0) ? y[0] : 0;
592 }
593
594 double varmin = - TMath::Infinity();
595 double varmax = TMath::Infinity();
596 const RooRealVar* var = dynamic_cast<RooRealVar*>(fParameters.first() );
597 if (var) {
598 varmin = var->getMin();
599 varmax = var->getMax();
600 }
601
602
603 // find ymin and ymax and corresponding values
604 double ymin = TMath::MinElement(n,y);
605 double ymax = TMath::MaxElement(n,y);
606 // cannot find intercept in the full range - return min or max valie
607 if (ymax < y0) {
608 return (lowSearch) ? varmax : varmin;
609 }
610 if (ymin > y0) {
611 return (lowSearch) ? varmin : varmax;
612 }
613
614 double xmin = axmin;
615 double xmax = axmax;
616
617 // case no range is given check if need to extrapolate to lower/upper values
618 if (axmin >= axmax ) {
619
620#ifdef DO_DEBUG
621 std::cout << "No range given - check if extrapolation is needed " << std::endl;
622#endif
623
624 xmin = graph.GetX()[0];
625 xmax = graph.GetX()[n-1];
626
627 double yfirst = graph.GetY()[0];
628 double ylast = graph.GetY()[n-1];
629
630 // distinguish the case we have lower /upper limits
631 // check if a possible crossing exists otherwise return variable min/max
632
633 // do lower extrapolation
634 if ( (ymax < y0 && !lowSearch) || ( yfirst > y0 && lowSearch) ) {
635 xmin = varmin;
636 }
637 // do upper extrapolation
638 if ( (ymax < y0 && lowSearch) || ( ylast > y0 && !lowSearch) ) {
639 xmax = varmax;
640 }
641 }
642
643 auto func = [&](double x) {
644 return (fInterpolOption == kSpline) ? graph.Eval(x, nullptr, "S") - y0 : graph.Eval(x) - y0;
645 };
647
649 brf.SetFunction(f1d,xmin,xmax);
650 brf.SetNpx(std::max(graph.GetN()*2,100) );
651#ifdef DO_DEBUG
652 std::cout << "findind root for " << xmin << " , "<< xmax << "f(x) : " << graph.Eval(xmin) << " , " << graph.Eval(0.5*(xmax+xmin))
653 << " , " << graph.Eval(xmax) << " target " << y0 << std::endl;
654#endif
655
656 bool ret = brf.Solve(100, 1.E-16, 1.E-6);
657 if (!ret) {
658 ooccoutE(this,Eval) << "HypoTestInverterResult - interpolation failed for interval [" << xmin << "," << xmax
659 << " ] g(xmin,xmax) =" << graph.Eval(xmin) << "," << graph.Eval(xmax)
660 << " target=" << y0 << " return inf" << std::endl
661 << "One may try to clean up invalid points using HypoTestInverterResult::ExclusionCleanup()." << std::endl;
662 return TMath::Infinity();
663 }
664 double limit = brf.Root();
665
666#ifdef DO_DEBUG
667 if (lowSearch) std::cout << "lower limit search : ";
668 else std::cout << "Upper limit search : ";
669 std::cout << "interpolation done between " << xmin << " and " << xmax
670 << "\n Found limit using RootFinder is " << limit << std::endl;
671
672 TString fname = "graph_upper.root";
673 if (lowSearch) fname = "graph_lower.root";
674 {
675 std::unique_ptr<TFile> file{TFile::Open(fname,"RECREATE")};
676 graph.Write("graph");
677 }
678#endif
679
680 // look in case if a new intersection exists
681 // only when boundaries are not given
682 if (axmin >= axmax) {
683 int index = TMath::BinarySearch(n, graph.GetX(), limit);
684#ifdef DO_DEBUG
685 std::cout << "do new interpolation dividing from " << index << " and " << y[index] << std::endl;
686#endif
687
688 if (lowSearch && index >= 1 && (y[0] - y0) * ( y[index]- y0) < 0) {
689 //search if another intersection exists at a lower value
690 limit = GetGraphX(graph, y0, lowSearch, graph.GetX()[0], graph.GetX()[index] );
691 }
692 else if (!lowSearch && index < n-2 && (y[n-1] - y0) * ( y[index+1]- y0) < 0) {
693 // another intersection exists at an higher value
694 limit = GetGraphX(graph, y0, lowSearch, graph.GetX()[index+1], graph.GetX()[n-1] );
695 }
696 }
697 // return also xmin, xmax values
698 axmin = xmin;
699 axmax = xmax;
700
701 return limit;
702}
703
704////////////////////////////////////////////////////////////////////////////////
705/// interpolate to find a limit value
706/// Use a linear or a spline interpolation depending on the interpolation option
707
709{
710
711 // variable minimum and maximum
712 double varmin = - TMath::Infinity();
713 double varmax = TMath::Infinity();
714 const RooRealVar* var = dynamic_cast<RooRealVar*>(fParameters.first() );
715 if (var) {
716 varmin = var->getMin();
717 varmax = var->getMax();
718 }
719
720 if (ArraySize()<2) {
721 double val = (lowSearch) ? xmin : xmax;
722 coutW(Eval) << "HypoTestInverterResult::FindInterpolatedLimit"
723 << " - not enough points to get the inverted interval - return "
724 << val << std::endl;
727 return (lowSearch) ? fLowerLimit : fUpperLimit;
728 }
729
730 // sort the values in x
731 int n = ArraySize();
732 std::vector<unsigned int> index(n );
733 TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);
734 // make a graph with the sorted point
735 TGraph graph(n);
736 for (int i = 0; i < n; ++i)
737 graph.SetPoint(i, GetXValue(index[i]), GetYValue(index[i] ) );
738
739
740 //std::cout << " search for " << lowSearch << " xmin = " << xmin << " xmax " << xmax << std::endl;
741
742
743 // search first for min/max in the given range
744 if (xmin >= xmax) {
745
746
747 // search for maximum between the point
748 double * itrmax = std::max_element(graph.GetY() , graph.GetY() +n);
749 double ymax = *itrmax;
750 int iymax = itrmax - graph.GetY();
751 double xwithymax = graph.GetX()[iymax];
752
753#ifdef DO_DEBUG
754 std::cout << " max of y " << iymax << " " << xwithymax << " " << ymax << " target is " << target << std::endl;
755#endif
756 // look if maximum is above/below target
757 if (ymax > target) {
758 if (lowSearch) {
759 if ( iymax > 0) {
760 // low search (minimum is first point or minimum range)
761 xmin = ( graph.GetY()[0] <= target ) ? graph.GetX()[0] : varmin;
762 xmax = xwithymax;
763 }
764 else {
765 // no room for lower limit
767 return fLowerLimit;
768 }
769 }
770 if (!lowSearch ) {
771 // up search
772 if ( iymax < n-1 ) {
773 xmin = xwithymax;
774 xmax = ( graph.GetY()[n-1] <= target ) ? graph.GetX()[n-1] : varmax;
775 }
776 else {
777 // no room for upper limit
779 return fUpperLimit;
780 }
781 }
782 }
783 else {
784 // in case is below the target
785 // find out if is a lower or upper search
786 if (iymax <= (n-1)/2 ) {
787 lowSearch = false;
789 }
790 else {
791 lowSearch = true;
793 }
794 }
795
796#ifdef DO_DEBUG
797 std::cout << " found xmin, xmax = " << xmin << " " << xmax << " for search " << lowSearch << std::endl;
798#endif
799
800 // now come here if I have already found a lower/upper limit
801 // i.e. I am calling routine for the second time
802#ifdef ISNEEDED
803 // should not really come here
804 if (lowSearch && fUpperLimit < varmax) {
805 xmin = fXValues[ index.front() ];
806 // find xmax (is first point before upper limit)
808 if (upI < 1) return xmin;
809 xmax = GetXValue(upI);
810 }
811 else if (!lowSearch && fLowerLimit > varmin ) {
812 // find xmin (is first point after lower limit)
814 if (lowI >= n-1) return xmax;
816 xmax = fXValues[ index.back() ];
817 }
818#endif
819 }
820
821#ifdef DO_DEBUG
822 std::cout << "finding " << lowSearch << " limit between " << xmin << " " << xmax << std::endl;
823#endif
824
825 // compute now the limit using the TGraph interpolations routine
826 double limit = GetGraphX(graph, target, lowSearch, xmin, xmax);
827 if (lowSearch) fLowerLimit = limit;
828 else fUpperLimit = limit;
829 // estimate the error
831
832 TString limitType = (lowSearch) ? "lower" : "upper";
833 ooccoutD(this,Eval) << "HypoTestInverterResult::FindInterpolateLimit "
834 << "the computed " << limitType << " limit is " << limit << " +/- " << error << std::endl;
835
836#ifdef DO_DEBUG
837 std::cout << "Found limit is " << limit << " +/- " << error << std::endl;
838#endif
839
840
841 if (lowSearch) return fLowerLimit;
842 return fUpperLimit;
843
844
845// if (lowSearch && !TMath::IsNaN(fUpperLimit)) return fLowerLimit;
846// if (!lowSearch && !TMath::IsNaN(fLowerLimit)) return fUpperLimit;
847// // is this needed ?
848// // we call again the function for the upper limits
849
850// // now perform the opposite search on the complement interval
851// if (lowSearch) {
852// xmin = xmax;
853// xmax = varmax;
854// } else {
855// xmax = xmin;
856// xmin = varmin;
857// }
858// double limit2 = GetGraphX(graph, target, !lowSearch, xmin, xmax);
859// if (!lowSearch) fLowerLimit = limit2;
860// else fUpperLimit = limit2;
861
862// CalculateEstimatedError( target, !lowSearch, xmin, xmax);
863
864// #ifdef DO_DEBUG
865// std::cout << "other limit is " << limit2 << std::endl;
866// #endif
867
868// return (lowSearch) ? fLowerLimit : fUpperLimit;
869
870}
871
872////////////////////////////////////////////////////////////////////////////////
873/// - if mode = 0
874/// find closest point to target in Y, the object closest to the target which is 3 sigma from the target
875/// and has smaller error
876/// - if mode = 1
877/// find 2 closest point to target in X and between these two take the one closer to the target
878/// - if mode = 2 as in mode = 1 but return the lower point not the closest one
879/// - if mode = 3 as in mode = 1 but return the upper point not the closest one
880
882{
883
884 int bestIndex = -1;
885 int closestIndex = -1;
886 if (mode == 0) {
887 double smallestError = 2; // error must be < 1
888 double bestValue = 2;
889 for (int i=0; i<ArraySize(); i++) {
890 double dist = std::abs(GetYValue(i)-target);
891 if ( dist <3 *GetYError(i) ) { // less than 1 sigma from target CL
892 if (GetYError(i) < smallestError ) {
894 bestIndex = i;
895 }
896 }
897 if ( dist < bestValue) {
898 bestValue = dist;
899 closestIndex = i;
900 }
901 }
902 if (bestIndex >=0) return bestIndex;
903 // if no points found just return the closest one to the target
904 return closestIndex;
905 }
906 // else mode = 1,2,3
907 // find the two closest points to limit value
908 // sort the array first
909 int n = fXValues.size();
910 std::vector<unsigned int> indx(n);
911 TMath::SortItr(fXValues.begin(), fXValues.end(), indx.begin(), false);
912 std::vector<double> xsorted( n);
913 for (int i = 0; i < n; ++i) xsorted[i] = fXValues[indx[i] ];
915
916#ifdef DO_DEBUG
917 std::cout << "finding closest point to " << xtarget << " is " << index1 << " " << indx[index1] << std::endl;
918#endif
919
920 // case xtarget is outside the range (before or afterwards)
921 if (index1 < 0) return indx[0];
922 if (index1 >= n-1) return indx[n-1];
923 int index2 = index1 +1;
924
925 if (mode == 2) return (GetXValue(indx[index1]) < GetXValue(indx[index2])) ? indx[index1] : indx[index2];
926 if (mode == 3) return (GetXValue(indx[index1]) > GetXValue(indx[index2])) ? indx[index1] : indx[index2];
927 // get smaller point of the two (mode == 1)
928 if (std::abs(GetYValue(indx[index1])-target) <= std::abs(GetYValue(indx[index2])-target) )
929 return indx[index1];
930 return indx[index2];
931
932}
933
934////////////////////////////////////////////////////////////////////////////////
935
937{
938 if (fFittedLowerLimit) return fLowerLimit;
939 //std::cout << "finding point with cl = " << 1-(1-ConfidenceLevel())/2 << std::endl;
941 // find both lower/upper limit
943 } else {
944 //LM: I think this is never called
946 }
947 return fLowerLimit;
948}
949
950////////////////////////////////////////////////////////////////////////////////
951
953{
954 //std::cout << "finding point with cl = " << (1-ConfidenceLevel())/2 << std::endl;
955 if (fFittedUpperLimit) return fUpperLimit;
958 } else {
959 //LM: I think this is never called
961 }
962 return fUpperLimit;
963}
964
965////////////////////////////////////////////////////////////////////////////////
966/// Return an error estimate on the upper(lower) limit. This is the error on
967/// either CLs or CLsplusb divided by an estimate of the slope at this
968/// point.
969
971{
972
973 if (ArraySize()==0) {
974 coutW(Eval) << "HypoTestInverterResult::CalculateEstimateError"
975 << "Empty result \n";
976 return 0;
977 }
978
979 if (ArraySize()<2) {
980 coutW(Eval) << "HypoTestInverterResult::CalculateEstimateError"
981 << " only points - return its error\n";
982 return GetYError(0);
983 }
984
985 // it does not make sense in case of asymptotic which do not have point errors
986 if (!GetNullTestStatDist(0) ) return 0;
987
988 TString type = (!lower) ? "upper" : "lower";
989
990#ifdef DO_DEBUG
991 std::cout << "calculate estimate error " << type << " between " << xmin << " and " << xmax << std::endl;
992 std::cout << "computed limit is " << ( (lower) ? fLowerLimit : fUpperLimit ) << std::endl;
993#endif
994
995 // make a TGraph Errors with the sorted points
996 std::vector<unsigned int> indx(fXValues.size());
997 TMath::SortItr(fXValues.begin(), fXValues.end(), indx.begin(), false);
998 // make a graph with the sorted point
999 TGraphErrors graph;
1000 int ip = 0;
1001 int np = 0;
1002 for (int i = 0; i < ArraySize(); ++i) {
1003 if ( (xmin < xmax) && ( GetXValue(indx[i]) >= xmin && GetXValue(indx[i]) <= xmax) ) {
1004 np++;
1005 // exclude points with zero or very small errors
1006 if (GetYError(indx[i] ) > 1.E-6) {
1007 graph.SetPoint(ip, GetXValue(indx[i]), GetYValue(indx[i] ) );
1008 graph.SetPointError(ip, 0., GetYError(indx[i]) );
1009 ip++;
1010 }
1011 }
1012 }
1013 if (graph.GetN() < 2) {
1014 if (np >= 2) coutW(Eval) << "HypoTestInverterResult::CalculateEstimatedError - no valid points - cannot estimate the " << type << " limit error " << std::endl;
1015 return 0;
1016 }
1017
1018 double minX = xmin;
1019 double maxX = xmax;
1020 if (xmin >= xmax) {
1021 minX = fXValues[ indx.front() ];
1022 maxX = fXValues[ indx.back() ];
1023 }
1024
1025 TF1 fct("fct", "exp([0] * (x - [2] ) + [1] * (x-[2])**2)", minX, maxX);
1026 double scale = maxX-minX;
1027 if (lower) {
1028 fct.SetParameters( 2./scale, 0.1/scale, graph.GetX()[0] );
1029 fct.SetParLimits(0,0,100./scale);
1030 fct.SetParLimits(1,0, 10./scale); }
1031 else {
1032 fct.SetParameters( -2./scale, -0.1/scale );
1033 fct.SetParLimits(0,-100./scale, 0);
1034 fct.SetParLimits(1,-100./scale, 0); }
1035
1036 if (graph.GetN() < 3) fct.FixParameter(1,0.);
1037
1038 // find the point closest to the limit
1039 double limit = (!lower) ? fUpperLimit : fLowerLimit;
1040 if (TMath::IsNaN(limit)) return 0; // cannot do if limit not computed
1041
1042
1043#ifdef DO_DEBUG
1044 TCanvas * c1 = new TCanvas();
1045 std::cout << "fitting for limit " << type << "between " << minX << " , " << maxX << " points considered " << graph.GetN() << std::endl;
1046 int fitstat = graph.Fit(&fct," EX0");
1047 graph.SetMarkerStyle(20);
1048 graph.Draw("AP");
1049 graph.Print();
1050 c1->SaveAs(TString::Format("graphFit_%s.pdf",type.Data()) );
1051 delete c1;
1052#else
1053 int fitstat = graph.Fit(&fct,"Q EX0");
1054#endif
1055
1056 int index = FindClosestPointIndex(target, 1, limit);
1057 double theError = 0;
1058 if (fitstat == 0) {
1059 double errY = GetYError(index);
1060 if (errY > 0) {
1061 double m = fct.Derivative( GetXValue(index) );
1062 theError = std::min(std::abs( GetYError(index) / m), maxX-minX);
1063 }
1064 }
1065 else {
1066 coutW(Eval) << "HypoTestInverterResult::CalculateEstimatedError - cannot estimate the " << type << " limit error " << std::endl;
1067 theError = 0;
1068 }
1069 if (lower) {
1071 } else {
1073 }
1074
1075#ifdef DO_DEBUG
1076 std::cout << "closes point to the limit is " << index << " " << GetXValue(index) << " and has error " << GetYError(index) << std::endl;
1077#endif
1078
1079 return theError;
1080}
1081
1082////////////////////////////////////////////////////////////////////////////////
1083/// need to have compute first lower limit
1084
1086{
1088 if (fLowerLimitError >= 0) return fLowerLimitError;
1089 // try to recompute the error
1090 return CalculateEstimatedError(1-ConfidenceLevel(), true);
1091}
1092
1093////////////////////////////////////////////////////////////////////////////////
1094
1096{
1098 if (fUpperLimitError >= 0) return fUpperLimitError;
1099 // try to recompute the error
1100 return CalculateEstimatedError(1-ConfidenceLevel(), false);
1101}
1102
1103////////////////////////////////////////////////////////////////////////////////
1104/// get the background test statistic distribution
1105
1107
1109 if (!firstResult) return nullptr;
1110 return firstResult->GetBackGroundIsAlt() ? firstResult->GetAltDistribution() : firstResult->GetNullDistribution();
1111}
1112
1113////////////////////////////////////////////////////////////////////////////////
1114/// get the signal and background test statistic distribution
1115
1118 if (!result) return nullptr;
1119 return !result->GetBackGroundIsAlt() ? result->GetAltDistribution() : result->GetNullDistribution();
1120}
1121
1122////////////////////////////////////////////////////////////////////////////////
1123/// get the expected p-value distribution at the scanned point index
1124
1126
1127 if (index < 0 || index >= ArraySize() ) return nullptr;
1128
1129 if (fExpPValues.GetSize() == ArraySize() ) {
1130 return static_cast<SamplingDistribution*>( fExpPValues.At(index)->Clone());
1131 }
1132
1133 static bool useFirstB = false;
1134 // get the S+B distribution
1135 int bIndex = (useFirstB) ? 0 : index;
1136
1139
1141
1143
1144 // create a new HypoTestResult
1146 tempResult.SetPValueIsRightTail( result->GetPValueIsRightTail() );
1147 tempResult.SetBackgroundAsAlt( true);
1148 // ownership of SamplingDistribution is in HypoTestResult class now
1149 tempResult.SetNullDistribution( new SamplingDistribution(*sbDistribution) );
1150 tempResult.SetAltDistribution( new SamplingDistribution(*bDistribution ) );
1151
1152 std::vector<double> values(bDistribution->GetSize());
1153 for (int i = 0; i < bDistribution->GetSize(); ++i) {
1154 tempResult.SetTestStatisticData( bDistribution->GetSamplingDistribution()[i] );
1155 values[i] = (fUseCLs) ? tempResult.CLs() : tempResult.CLsplusb();
1156 }
1157 return new SamplingDistribution("expected values","expected values",values);
1158 }
1159 // in case b abs sbDistribution are null assume is coming from the asymptotic calculator
1160 // hard -coded this value (no really needed to be used by user)
1163 const double smax = fgAsymptoticMaxSigma;
1164 const int npoints = fgAsymptoticNumPoints;
1165 const double dsig = 2* smax/ (npoints-1) ;
1166 std::vector<double> values(npoints);
1167 for (int i = 0; i < npoints; ++i) {
1168 double nsig = -smax + dsig*i;
1169 double pval = AsymptoticCalculator::GetExpectedPValues( result->NullPValue(), result->AlternatePValue(), nsig, fUseCLs, !fIsTwoSided);
1170 if (pval < 0) { return nullptr;}
1171 values[i] = pval;
1172 }
1173 return new SamplingDistribution("Asymptotic expected values","Asymptotic expected values",values);
1174
1175}
1176
1177////////////////////////////////////////////////////////////////////////////////
1178/// get the limit distribution (lower/upper depending on the flag)
1179/// by interpolating the expected p values for each point
1180
1182 if (ArraySize()<2) {
1183 coutE(Eval) << "HypoTestInverterResult::GetLimitDistribution"
1184 << " not enough points - return 0 " << std::endl;
1185 return nullptr;
1186 }
1187
1188 ooccoutD(this,Eval) << "HypoTestInverterResult - computing limit distribution...." << std::endl;
1189
1190
1191 // find optimal size by looking at the PValue distribution obtained
1192 int npoints = ArraySize();
1193 std::vector<SamplingDistribution*> distVec( npoints );
1194 double sum = 0;
1195 for (unsigned int i = 0; i < distVec.size(); ++i) {
1197 // sort the distributions
1198 // hack (by calling InverseCDF(0) will sort the sampling distribution
1199 if (distVec[i] ) {
1200 distVec[i]->InverseCDF(0);
1201 sum += distVec[i]->GetSize();
1202 }
1203 }
1204 int size = int( sum/ npoints);
1205
1206 if (size < 10) {
1207 ooccoutW(this,InputArguments) << "HypoTestInverterResult - set a minimum size of 10 for limit distribution" << std::endl;
1208 size = 10;
1209 }
1210
1211
1212 double target = 1-fConfidenceLevel;
1213
1214 // vector with the quantiles of the p-values for each scanned poi point
1215 std::vector< std::vector<double> > quantVec(npoints );
1216 for (int i = 0; i < npoints; ++i) {
1217
1218 if (!distVec[i]) continue;
1219
1220 // make quantiles from the sampling distributions of the expected p values
1221 std::vector<double> pvalues = distVec[i]->GetSamplingDistribution();
1222 delete distVec[i]; distVec[i] = nullptr;
1223 std::sort(pvalues.begin(), pvalues.end());
1224 // find the quantiles of the distribution
1225 double p[1] = {0};
1226 double q[1] = {0};
1227
1228 quantVec[i] = std::vector<double>(size);
1229 for (int ibin = 0; ibin < size; ++ibin) {
1230 // exclude for a bug in TMath::Quantiles last item
1231 p[0] = std::min( (ibin+1) * 1./double(size), 1.0);
1232 // use the type 1 which give the point value
1233 TMath::Quantiles(pvalues.size(), 1, &pvalues[0], q, p, true, nullptr, 1 );
1234 (quantVec[i])[ibin] = q[0];
1235 }
1236
1237 }
1238
1239 // sort the values in x
1240 std::vector<unsigned int> index( npoints );
1241 TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);
1242
1243 // SamplingDistribution * dist0 = distVec[index[0]];
1244 // SamplingDistribution * dist1 = distVec[index[1]];
1245
1246
1247 std::vector<double> limits(size);
1248 // loop on the p values and find the limit for each expected point in the quantiles vector
1249 for (int j = 0; j < size; ++j ) {
1250
1251 TGraph g;
1252 for (int k = 0; k < npoints ; ++k) {
1253 if (!quantVec[index[k]].empty() )
1254 g.SetPoint(g.GetN(), GetXValue(index[k]), (quantVec[index[k]])[j] );
1255 }
1256
1257 limits[j] = GetGraphX(g, target, lower);
1258
1259 }
1260
1261 if (lower) {
1262 return new SamplingDistribution("Expected lower Limit","Expected lower limits",limits);
1263 } else {
1264 return new SamplingDistribution("Expected upper Limit", "Expected upper limits", limits);
1265 }
1266}
1267
1268////////////////////////////////////////////////////////////////////////////////
1269/// Get the expected lower limit
1270/// nsig is used to specify which expected value of the UpperLimitDistribution
1271/// For example
1272/// - nsig = 0 (default value) returns the expected value
1273/// - nsig = -1 returns the lower band value at -1 sigma
1274/// - nsig + 1 return the upper value
1275/// - opt = "" (default) : compute limit by interpolating all the p values, find the corresponding limit distribution
1276/// and then find the quantiles in the limit distribution
1277/// ioption = "P" is the method used for plotting. One Finds the corresponding nsig quantile in the p values and then
1278/// interpolates them
1279
1280double HypoTestInverterResult::GetExpectedLowerLimit(double nsig, const char * opt ) const {
1281
1282 return GetExpectedLimit(nsig, true, opt);
1283}
1284
1285////////////////////////////////////////////////////////////////////////////////
1286/// Get the expected upper limit
1287/// nsig is used to specify which expected value of the UpperLimitDistribution
1288/// For example
1289/// - nsig = 0 (default value) returns the expected value
1290/// - nsig = -1 returns the lower band value at -1 sigma
1291/// - nsig + 1 return the upper value
1292/// - opt is an option specifying the type of method used for computing the upper limit
1293/// - opt = "" (default) : compute limit by interpolating all the p values, find the corresponding limit distribution
1294/// and then find the quantiles in the limit distribution
1295/// ioption = "P" is the method used for plotting. One Finds the corresponding nsig quantile in the p values and then
1296/// interpolates them
1297
1298double HypoTestInverterResult::GetExpectedUpperLimit(double nsig, const char * opt ) const {
1299
1300 return GetExpectedLimit(nsig, false, opt);
1301}
1302
1303////////////////////////////////////////////////////////////////////////////////
1304/// get expected limit (lower/upper) depending on the flag
1305/// for asymptotic is a special case (the distribution is generated an step in sigma values)
1306/// distinguish asymptotic looking at the hypotest results
1307/// if option = "P" get expected limit using directly quantiles of p value distribution
1308/// else (default) find expected limit by obtaining first a full limit distributions
1309/// The last one is in general more correct
1310
1311double HypoTestInverterResult::GetExpectedLimit(double nsig, bool lower, const char * opt ) const {
1312
1313 const int nEntries = ArraySize();
1314 if (nEntries <= 0) return (lower) ? 1 : 0; // return 1 for lower, 0 for upper
1315
1316 HypoTestResult * r = dynamic_cast<HypoTestResult *> (fYObjects.First() );
1317 assert(r != nullptr);
1318 if (!r->GetNullDistribution() && !r->GetAltDistribution() ) {
1319 // we are in the asymptotic case
1320 // get the limits obtained at the different sigma values
1322 if (!limitDist) return 0;
1323 const std::vector<double> & values = limitDist->GetSamplingDistribution();
1324 if (values.size() <= 1) return 0;
1325 double dsig = 2* fgAsymptoticMaxSigma/ (values.size() -1) ;
1326 int i = (int) TMath::Floor ( (nsig + fgAsymptoticMaxSigma)/dsig + 0.5);
1327 return values[i];
1328 }
1329
1330 double p[1] = {0};
1331 double q[1] = {0};
1333
1334 // for CLs+b can get the quantiles of p-value distribution and
1335 // interpolate them
1336 // In the case of CLs (since it is not a real p-value anymore but a ratio)
1337 // then it is needed to get first all limit distribution values and then the quantiles
1338
1339 TString option(opt);
1340 option.ToUpper();
1341 if (option.Contains("P")) {
1342
1343 TGraph g;
1344
1345 // sort the arrays based on the x values
1346 std::vector<unsigned int> index(nEntries);
1347 TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);
1348
1349 for (int j=0; j<nEntries; ++j) {
1350 int i = index[j]; // i is the order index
1352 if (!s) {
1353 ooccoutI(this,Eval) << "HypoTestInverterResult - cannot compute expected p value distribution for point, x = "
1354 << GetXValue(i) << " skip it " << std::endl;
1355 continue;
1356 }
1357 const std::vector<double> & values = s->GetSamplingDistribution();
1358 double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
1359 TMath::Quantiles(values.size(), 1, x,q,p,false);
1360 g.SetPoint(g.GetN(), fXValues[i], q[0] );
1361 delete s;
1362 }
1363 if (g.GetN() < 2) {
1364 ooccoutE(this,Eval) << "HypoTestInverterResult - cannot compute limits , not enough points, n = " << g.GetN() << std::endl;
1365 return 0;
1366 }
1367
1368 // interpolate the graph to obtain the limit
1369 double target = 1-fConfidenceLevel;
1370 return GetGraphX(g, target, lower);
1371
1372 }
1373 // here need to use the limit distribution
1375 if (!limitDist) return 0;
1376 const std::vector<double> & values = limitDist->GetSamplingDistribution();
1377 double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
1378 TMath::Quantiles(values.size(), 1, x,q,p,false);
1379 return q[0];
1380
1381}
#define g(i)
Definition RSha256.hxx:105
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutI(a)
#define coutW(a)
#define coutE(a)
#define ooccoutW(o, a)
#define ooccoutI(o, a)
#define ooccoutD(o, a)
#define ooccoutE(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
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 target
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 np
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 char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char mode
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:146
float xmin
float * q
float ymin
float xmax
float ymax
Class for finding the root of a one dimensional function using the Brent algorithm.
Functor1D class for one-dimensional functions.
Definition Functor.h:97
const_iterator begin() const
const_iterator end() const
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
RooAbsArg * first() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
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.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
static double GetExpectedPValues(double pnull, double palt, double nsigma, bool usecls, bool oneSided=true)
function given the null and the alt p value - return the expected one given the N - sigma value
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
double GetYValue(int index) const
function to return the value of the confidence level for the i^th entry in the results
double LowerLimitEstimatedError()
rough estimation of the error on the computed bound of the confidence interval Estimate of lower limi...
double FindInterpolatedLimit(double target, bool lowSearch=false, double xmin=1, double xmax=0.0)
interpolate to find a limit value Use a linear or a spline interpolation depending on the interpolati...
double UpperLimitEstimatedError()
Estimate of lower limit error function evaluates only a rough error on the lower limit.
SamplingDistribution * GetSignalAndBackgroundTestStatDist(int index) const
get the signal and background test statistic distribution
HypoTestResult * GetResult(int index) const
return a pointer to the i^th result object
HypoTestInverterResult & operator=(const HypoTestInverterResult &other)
operator =
InterpolOption_t fInterpolOption
interpolation option (linear or spline)
int ArraySize() const
number of entries in the results array
double CLsplusbError(int index) const
return the observed CLsplusb value for the i-th entry
double CLsError(int index) const
return the observed CLb value for the i-th entry
double GetGraphX(const TGraph &g, double y0, bool lowSearch, double &xmin, double &xmax) const
return the X value of the given graph for the target value y0 the graph is evaluated using linear int...
double CLbError(int index) const
return the observed CLb value for the i-th entry
double GetExpectedUpperLimit(double nsig=0, const char *opt="") const
get Limit value corresponding at the desired nsigma level (0) is median -1 sigma is 1 sigma
HypoTestInverterResult(const char *name=nullptr)
default constructor
double GetXValue(int index) const
function to return the value of the parameter of interest for the i^th entry in the results
double UpperLimit() override
return the interval upper limit
bool Add(const HypoTestInverterResult &otherResult)
merge with the content of another HypoTestInverterResult object
int FindClosestPointIndex(double target, int mode=0, double xtarget=0)
TList fExpPValues
list of expected sampling distribution for each point
double GetExpectedLowerLimit(double nsig=0, const char *opt="") const
get Limit value corresponding at the desired nsigma level (0) is median -1 sigma is 1 sigma
bool fIsTwoSided
two sided scan (look for lower/upper limit)
double CalculateEstimatedError(double target, bool lower=true, double xmin=1, double xmax=0.0)
Return an error estimate on the upper(lower) limit.
static int fgAsymptoticNumPoints
number of points used to build expected p-values
static double fgAsymptoticMaxSigma
max sigma value used to scan asymptotic expected p values
int ExclusionCleanup()
remove points that appear to have failed.
double LowerLimit() override
lower and upper bound of the confidence interval (to get upper/lower limits, multiply the size( = 1-c...
double CLs(int index) const
return the observed CLb value for the i-th entry
int FindIndex(double xvalue) const
find the index corresponding at the poi value xvalue If no points is found return -1 Note that a tole...
double GetYError(int index) const
function to return the estimated error on the value of the confidence level for the i^th entry in the...
double CLb(int index) const
return the observed CLb value for the i-th entry
TList fYObjects
list of HypoTestResult for each point
double GetExpectedLimit(double nsig, bool lower, const char *opt="") const
get expected limit (lower/upper) depending on the flag for asymptotic is a special case (the distribu...
SamplingDistribution * GetExpectedPValueDist(int index) const
return expected distribution of p-values (Cls or Clsplusb)
double CLsplusb(int index) const
return the observed CLsplusb value for the i-th entry
SamplingDistribution * GetLimitDistribution(bool lower) const
get the limit distribution (lower/upper depending on the flag) by interpolating the expected p values...
SamplingDistribution * GetNullTestStatDist(int index) const
same in terms of alt and null
SamplingDistribution * GetBackgroundTestStatDist(int index) const
get the background test statistic distribution
HypoTestResult is a base class for results from hypothesis tests.
TObject * Clone(const char *newname=nullptr) const override
clone method, required since some data members cannot rely on the streamers to copy them
This class simply holds a sampling distribution of some test statistic.
const std::vector< double > & GetSamplingDistribution() const
Get test statistics values.
SimpleInterval is a concrete implementation of the ConfInterval interface.
double fUpperLimit
upper interval limit
RooArgSet fParameters
set containing the parameter of interest
double ConfidenceLevel() const override
return the confidence interval
double fLowerLimit
lower interval limit
double fConfidenceLevel
confidence level
SimpleInterval & operator=(const SimpleInterval &other)
default constructor
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:43
The Canvas class.
Definition TCanvas.h:23
virtual void RemoveAll(TCollection *col)
Remove all objects in collection col from this collection.
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
1-Dim function class
Definition TF1.h:182
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:3787
A TGraphErrors is a TGraph with error bars.
void Print(Option_t *chopt="") const override
Print graph and errors values.
virtual void SetPointError(Double_t ex, Double_t ey)
Set ex and ey values for point pointed by the mouse.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2387
Double_t * GetY() const
Definition TGraph.h:139
Int_t GetN() const
Definition TGraph.h:131
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0)
Fit this graph with the global function named fname.
Definition TGraph.cxx:1285
Double_t * GetX() const
Definition TGraph.h:138
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
Interpolate points in this graph at x using a TSpline.
Definition TGraph.cxx:981
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:859
void Add(TObject *obj) override
Definition TList.h:81
TObject * First() const override
Return the first object in the list. Returns 0 when list is empty.
Definition TList.cxx:789
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
Definition TList.cxx:600
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:487
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition TObject.cxx:240
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:986
virtual TObject * RemoveAt(Int_t idx)
Basic string class.
Definition TString.h:138
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:2384
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:71
@ InputArguments
Namespace for the RooStats classes.
Definition CodegenImpl.h:65
TMath.
Definition TMathBase.h:35
Bool_t IsNaN(Double_t x)
Definition TMath.h:903
void SortItr(Iterator first, Iterator last, IndexIterator index, Bool_t down=kTRUE)
Sort the n1 elements of the Short_t array defined by its iterators.
Definition TMathBase.h:388
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:913
Double_t Floor(Double_t x)
Rounds x downward, returning the largest integral value that is not greater than x.
Definition TMath.h:691
T MinElement(Long64_t n, const T *a)
Returns minimum of array a of length n.
Definition TMath.h:971
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=nullptr, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities.
Definition TMath.cxx:1207
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Comparing floating points.
Definition TMath.h:429
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Comparing floating points.
Definition TMath.h:421
T MaxElement(Long64_t n, const T *a)
Returns maximum of array a of length n.
Definition TMath.h:979
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:329
Double_t Infinity()
Returns an infinity as defined by the IEEE standard.
Definition TMath.h:928
TMarker m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2338