ROOT  6.06/09
Reference Guide
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 
13 
14 // include header file of this class
16 
17 #include "RooStats/HybridResult.h"
20 #include "RooMsgService.h"
21 #include "RooGlobalFunc.h"
22 
23 #include "TF1.h"
24 #include "TGraphErrors.h"
25 #include <cmath>
26 #include "Math/BrentRootFinder.h"
27 #include "Math/WrappedFunction.h"
28 
29 #include "TCanvas.h"
31 
32 #include <algorithm>
33 
35 
36 using namespace RooStats;
37 using namespace RooFit;
38 using namespace std;
39 
40 
41 // initialize static value
42 double HypoTestInverterResult::fgAsymptoticMaxSigma = 5;
43 
44 
45 HypoTestInverterResult::HypoTestInverterResult(const char * name ) :
46  SimpleInterval(name),
47  fUseCLs(false),
48  fIsTwoSided(false),
49  fInterpolateLowerLimit(true),
50  fInterpolateUpperLimit(true),
51  fFittedLowerLimit(false),
52  fFittedUpperLimit(false),
53  fInterpolOption(kLinear),
54  fLowerLimitError(-1),
55  fUpperLimitError(-1),
56  fCLsCleanupThreshold(0.005)
57 {
58  // default constructor
59  fLowerLimit = TMath::QuietNaN();
60  fUpperLimit = TMath::QuietNaN();
61 
62  fYObjects.SetOwner();
63  fExpPValues.SetOwner();
64 }
65 
66 
67 HypoTestInverterResult::HypoTestInverterResult( const HypoTestInverterResult& other, const char * name ) :
68  SimpleInterval(other,name),
69  fUseCLs(other.fUseCLs),
70  fIsTwoSided(other.fIsTwoSided),
71  fInterpolateLowerLimit(other.fInterpolateLowerLimit),
72  fInterpolateUpperLimit(other.fInterpolateUpperLimit),
73  fFittedLowerLimit(other.fFittedLowerLimit),
74  fFittedUpperLimit(other.fFittedUpperLimit),
75  fInterpolOption(other.fInterpolOption),
76  fLowerLimitError(other.fLowerLimitError),
77  fUpperLimitError(other.fUpperLimitError),
78  fCLsCleanupThreshold(other.fCLsCleanupThreshold)
79 {
80  // copy constructor
83  int nOther = other.ArraySize();
84 
85  fXValues = other.fXValues;
86  for (int i = 0; i < nOther; ++i)
87  fYObjects.Add( other.fYObjects.At(i)->Clone() );
88  for (int i = 0; i < fExpPValues.GetSize() ; ++i)
89  fExpPValues.Add( other.fExpPValues.At(i)->Clone() );
90 
93 }
94 
95 
98 {
99  if (&other==this) {
100  return *this ;
101  }
102 
104  fLowerLimit = other.fLowerLimit;
105  fUpperLimit = other.fUpperLimit;
106  fUseCLs = other.fUseCLs;
107  fIsTwoSided = other.fIsTwoSided;
113  fLowerLimitError = other.fLowerLimitError;
114  fUpperLimitError = other.fUpperLimitError;
115  fCLsCleanupThreshold = other.fCLsCleanupThreshold;
116 
117  int nOther = other.ArraySize();
118  fXValues = other.fXValues;
119 
121  for (int i=0; i < nOther; ++i) {
122  fYObjects.Add( other.fYObjects.At(i)->Clone() );
123  }
125  for (int i=0; i < fExpPValues.GetSize() ; ++i) {
126  fExpPValues.Add( other.fExpPValues.At(i)->Clone() );
127  }
128 
131 
132  return *this;
133 }
134 
135 
137  const RooRealVar& scannedVariable,
138  double cl ) :
139  SimpleInterval(name,scannedVariable,TMath::QuietNaN(),TMath::QuietNaN(),cl),
140  fUseCLs(false),
141  fIsTwoSided(false),
142  fInterpolateLowerLimit(true),
143  fInterpolateUpperLimit(true),
144  fFittedLowerLimit(false),
145  fFittedUpperLimit(false),
146  fInterpolOption(kLinear),
147  fLowerLimitError(-1),
148  fUpperLimitError(-1),
149  fCLsCleanupThreshold(0.005)
150 {
151  // constructor
154 
155  // put a cloned copy of scanned variable to set in the interval
156  // to avoid I/O problem of the Result class -
157  // make the set owning the cloned copy (use clone istead of Clone to not copying all links)
158  fParameters.removeAll();
159  fParameters.takeOwnership();
160  fParameters.addOwned(*((RooRealVar *) scannedVariable.clone(scannedVariable.GetName()) ));
161 }
162 
164 {
165  // destructor
166  // explicitly empty the TLists - these contain pointers, not objects
169 
170  fYObjects.Delete();
172 }
173 
174 int
176 {
177  const int nEntries = ArraySize();
178 
179  // initialization
180  double nsig1(1.0);
181  double nsig2(2.0);
182  double p[5];
183  double q[5];
184  std::vector<double> qv;
185  qv.resize(11,-1.0);
186 
187  p[0] = ROOT::Math::normal_cdf(-nsig2);
188  p[1] = ROOT::Math::normal_cdf(-nsig1);
189  p[2] = 0.5;
190  p[3] = ROOT::Math::normal_cdf(nsig1);
191  p[4] = ROOT::Math::normal_cdf(nsig2);
192 
193  bool resultIsAsymptotic(false);
194  if (nEntries>=1) {
195  HypoTestResult* r = dynamic_cast<HypoTestResult *> ( GetResult(0) );
196  assert(r!=0);
197  if ( !r->GetNullDistribution() && !r->GetAltDistribution() ) {
198  resultIsAsymptotic = true;
199  }
200  }
201 
202  int nPointsRemoved(0);
203 
204  double CLsobsprev(1.0);
205  std::vector<double>::iterator itr = fXValues.begin();
206 
207  for (; itr!=fXValues.end();) {
208 
209  double x = (*itr);
210  int i = FindIndex(x);
211  //HypoTestResult* oneresult = GetResult(i);
212 
214  if (!s) break;
215 
216  /////////////////////////////////////////////////////////////////////////////////////////
217 
218  const std::vector<double> & values = s->GetSamplingDistribution();
219 
220  /// expected p-values
221  // special case for asymptotic results (cannot use TMath::quantile in that case)
222  if (resultIsAsymptotic) {
223  double maxSigma = 5; // == HypoTestInverterResult::fgAsymptoticMaxSigma; // MB: HACK
224  double dsig = 2.*maxSigma / (values.size() -1) ;
225  int i0 = (int) TMath::Floor ( ( -nsig2 + maxSigma )/dsig + 0.5 );
226  int i1 = (int) TMath::Floor ( ( -nsig1 + maxSigma )/dsig + 0.5 );
227  int i2 = (int) TMath::Floor ( ( maxSigma )/dsig + 0.5 );
228  int i3 = (int) TMath::Floor ( ( nsig1 + maxSigma )/dsig + 0.5 );
229  int i4 = (int) TMath::Floor ( ( nsig2 + maxSigma )/dsig + 0.5 );
230  //
231  q[0] = values[i0];
232  q[1] = values[i1];
233  q[2] = values[i2];
234  q[3] = values[i3];
235  q[4] = values[i4];
236  } else {
237  double * z = const_cast<double *>( &values[0] ); // need to change TMath::Quantiles
238  TMath::Quantiles(values.size(), 5, z, q, p, false);
239  }
240 
241  delete s;
242 
243  /// store useful quantities for reuse later ...
244  /// http://root.cern.ch/root/html532/src/RooStats__HypoTestInverterPlot.cxx.html#197
245  for (int j=0; j<5; ++j) { qv[j]=q[j]; }
246  qv[5] = CLs(i) ; //
247  qv[6] = CLsError(i) ; //
248  qv[7] = CLb(i) ; //
249  qv[8] = CLbError(i) ; //
250  qv[9] = CLsplusb(i) ; //
251  qv[10] = CLsplusbError(i) ; //
252  double CLsobs = qv[5];
253 
254  /////////////////////////////////////////////////////////////////////////////////////////
255 
256  bool removeThisPoint(false);
257 
258  // 1. CLs should drop, else skip this point
259  if (!removeThisPoint && resultIsAsymptotic && i>=1 && CLsobs>CLsobsprev) {
260  //StatToolsLogger << kERROR << "Asymptotic. CLs not dropping: " << CLsobs << ". Remove this point." << GEndl;
261  removeThisPoint = true;
262  } else { CLsobsprev = CLsobs; }
263 
264  // 2. CLs should not spike, else skip this point
265  if (!removeThisPoint && i>=1 && CLsobs>=0.9999) {
266  //StatToolsLogger << kERROR << "CLs spiking at 1.0: " << CLsobs << ". Remove this point." << GEndl;
267  removeThisPoint = true;
268  }
269  // 3. Not interested in CLs values that become too low.
270  if (!removeThisPoint && i>=1 && qv[4]<fCLsCleanupThreshold) { removeThisPoint = true; }
271 
272  // to remove or not to remove
273  if (removeThisPoint) {
274  itr = fXValues.erase(itr); // returned itr has been updated.
275  fYObjects.RemoveAt(i);
277  nPointsRemoved++;
278  continue;
279  } else { // keep
280  CLsobsprev = CLsobs;
281  itr++;
282  }
283  }
284 
285  // after cleanup, reset existing limits
286  fFittedUpperLimit = false;
287  fFittedLowerLimit = false;
289 
290  return nPointsRemoved;
291 }
292 
294 {
295  // Merge this HypoTestInverterResult with another
296  // HypoTestInverterResult passed as argument
297  // The merge is done by combining the HypoTestResult when the same point value exist in both results.
298  // If results exist at different points these are added in the new result
299  // NOTE: Merging of the expected p-values obtained with pseudo-data.
300  // When expected p-values exist in the result (i.e. when rebuild option is used when getting the expected
301  // limit distribuiton in the HYpoTestInverter) then the expected p-values are also merged. This is equivalent
302  // at merging the pseudo-data. However there can be an inconsistency if the expected p-values have been
303  // obtained with different toys. In this case the merge is done but a warning message is printed.
304 
305 
306  int nThis = ArraySize();
307  int nOther = otherResult.ArraySize();
308  if (nOther == 0) return true;
309  if (nOther != otherResult.fYObjects.GetSize() ) return false;
310  if (nThis != fYObjects.GetSize() ) return false;
311 
312  // cannot merge in case of inconsistent memebrs
313  if (fExpPValues.GetSize() > 0 && fExpPValues.GetSize() != nThis ) return false;
314  if (otherResult.fExpPValues.GetSize() > 0 && otherResult.fExpPValues.GetSize() != nOther ) return false;
315 
316  oocoutI(this,Eval) << "HypoTestInverterResult::Add - merging result from " << otherResult.GetName()
317  << " in " << GetName() << std::endl;
318 
319  bool addExpPValues = (fExpPValues.GetSize() == 0 && otherResult.fExpPValues.GetSize() > 0);
320  bool mergeExpPValues = (fExpPValues.GetSize() > 0 && otherResult.fExpPValues.GetSize() > 0);
321 
322  if (addExpPValues || mergeExpPValues)
323  oocoutI(this,Eval) << "HypoTestInverterResult::Add - merging also the expected p-values from pseudo-data" << std::endl;
324 
325 
326  // case current result is empty
327  // just make a simple copy of the other result
328  if (nThis == 0) {
329  fXValues = otherResult.fXValues;
330  for (int i = 0; i < nOther; ++i)
331  fYObjects.Add( otherResult.fYObjects.At(i)->Clone() );
332  for (int i = 0; i < fExpPValues.GetSize() ; ++i)
333  fExpPValues.Add( otherResult.fExpPValues.At(i)->Clone() );
334  }
335  // now do teh real mergemerge combining point with same value or adding extra ones
336  else {
337  for (int i = 0; i < nOther; ++i) {
338  double otherVal = otherResult.fXValues[i];
339  HypoTestResult * otherHTR = (HypoTestResult*) otherResult.fYObjects.At(i);
340  if (otherHTR == 0) continue;
341  bool sameXFound = false;
342  for (int j = 0; j < nThis; ++j) {
343  double thisVal = fXValues[j];
344 
345  // if same value merge the result
346  if ( (std::abs(otherVal) < 1 && TMath::AreEqualAbs(otherVal, thisVal,1.E-12) ) ||
347  (std::abs(otherVal) >= 1 && TMath::AreEqualRel(otherVal, thisVal,1.E-12) ) ) {
348  HypoTestResult * thisHTR = (HypoTestResult*) fYObjects.At(j);
349  thisHTR->Append(otherHTR);
350  sameXFound = true;
351  if (mergeExpPValues) {
353  // check if same toys have been used for the test statistic distribution
354  int thisNToys = (thisHTR->GetNullDistribution() ) ? thisHTR->GetNullDistribution()->GetSize() : 0;
355  int otherNToys = (otherHTR->GetNullDistribution() ) ? otherHTR->GetNullDistribution()->GetSize() : 0;
356  if (thisNToys != otherNToys )
357  oocoutW(this,Eval) << "HypoTestInverterResult::Add expexted p values have been generated with different toys " << thisNToys << " , " << otherNToys << std::endl;
358  }
359  break;
360  }
361  }
362  if (!sameXFound) {
363  // add the new result
364  fYObjects.Add(otherHTR->Clone() );
365  fXValues.push_back( otherVal);
366  }
367  // add in any case also when same x found
368  if (addExpPValues)
369  fExpPValues.Add( otherResult.fExpPValues.At(i)->Clone() );
370 
371 
372  }
373  }
374 
375  if (ArraySize() > nThis)
376  oocoutI(this,Eval) << "HypoTestInverterResult::Add - new number of points is " << fXValues.size()
377  << std::endl;
378  else
379  oocoutI(this,Eval) << "HypoTestInverterResult::Add - new toys/point is "
380  << ((HypoTestResult*) fYObjects.At(0))->GetNullDistribution()->GetSize()
381  << std::endl;
382 
383  // reset cached limit values
386 
387  return true;
388 }
389 
391 {
392  // Add a single point result (an HypoTestResult)
393  int i= FindIndex(x);
394  if (i<0) {
395  fXValues.push_back(x);
396  fYObjects.Add(res.Clone());
397  } else {
398  HypoTestResult* r= GetResult(i);
399  if (!r) return false;
400  r->Append(&res);
401  }
402 
403  // reset cached limit values
406 
407  return true;
408 }
409 
410 
411 
412 double HypoTestInverterResult::GetXValue( int index ) const
413 {
414  // function to return the value of the parameter of interest for the i^th entry in the results
415  if ( index >= ArraySize() || index<0 ) {
416  oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
417  return -999;
418  }
419 
420  return fXValues[index];
421 }
422 
423 double HypoTestInverterResult::GetYValue( int index ) const
424 {
425  // function to return the value of the confidence level for the i^th entry in the results
426  auto result = GetResult(index);
427  if ( !result ) {
428  return -999;
429  }
430 
431  if (fUseCLs) {
432  return result->CLs();
433  } else {
434  return result->CLsplusb(); // CLs+b
435  }
436 }
437 
438 double HypoTestInverterResult::GetYError( int index ) const
439 {
440  // function to return the estimated error on the value of the confidence level for the i^th entry in the results
441  auto result = GetResult(index);
442  if ( !result ) {
443  return -999;
444  }
445 
446  if (fUseCLs) {
447  return result->CLsError();
448  } else {
449  return result->CLsplusbError();
450  }
451 }
452 
453 double HypoTestInverterResult::CLb( int index ) const
454 {
455  // function to return the observed CLb value for the i-th entry
456  auto result = GetResult(index);
457  if ( !result ) {
458  return -999;
459  }
460  return result->CLb();
461 }
462 
463 double HypoTestInverterResult::CLsplusb( int index ) const
464 {
465  // function to return the observed CLs+b value for the i-th entry
466  auto result = GetResult(index);
467  if ( !result) {
468  return -999;
469  }
470  return result->CLsplusb();
471 }
472 double HypoTestInverterResult::CLs( int index ) const
473 {
474  // function to return the observed CLs value for the i-th entry
475  auto result = GetResult(index);
476  if ( !result ) {
477  return -999;
478  }
479  return result->CLs();
480 }
481 
482 double HypoTestInverterResult::CLbError( int index ) const
483 {
484  // function to return the error on the observed CLb value for the i-th entry
485  auto result = GetResult(index);
486  if ( !result ) {
487  return -999;
488  }
489  return result->CLbError();
490 }
491 
492 double HypoTestInverterResult::CLsplusbError( int index ) const
493 {
494  // function to return the error on the observed CLs+b value for the i-th entry
495  auto result = GetResult(index);
496  if ( ! result ) {
497  return -999;
498  }
499  return result->CLsplusbError();
500 }
501 
502 double HypoTestInverterResult::CLsError( int index ) const
503 {
504  // function to return the error on the observed CLs value for the i-th entry
505  auto result = GetResult(index);
506  if ( ! result ){
507  return -999;
508  }
509  return result->CLsError();
510 }
511 
513 {
514  // get the HypoTestResult object at the given index point
515  if ( index >= ArraySize() || index<0 ) {
516  oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
517  return 0;
518  }
519 
520  return ((HypoTestResult*) fYObjects.At(index));
521 }
522 
523 int HypoTestInverterResult::FindIndex(double xvalue) const
524 {
525  // find the index corresponding at the poi value xvalue
526  // If no points is found return -1
527  // Note that a tolerance is used of 10^-12 to find the closest point
528  const double tol = 1.E-12;
529  for (int i=0; i<ArraySize(); i++) {
530  double xpoint = fXValues[i];
531  if ( (std::abs(xvalue) > 1 && TMath::AreEqualRel( xvalue, xpoint, tol) ) ||
532  (std::abs(xvalue) < 1 && TMath::AreEqualAbs( xvalue, xpoint, tol) ) )
533  return i;
534  }
535  return -1;
536 }
537 
538 
539 struct InterpolatedGraph {
540  InterpolatedGraph( const TGraph & g, double target, const char * interpOpt) :
541  fGraph(g), fTarget(target), fInterpOpt(interpOpt) {}
542 
543  // return interpolated value for x - target
544  double operator() (double x) const {
545  return fGraph.Eval(x, (TSpline*) 0, fInterpOpt) - fTarget;
546  }
547  const TGraph & fGraph;
548  double fTarget;
549  TString fInterpOpt;
550 };
551 
552 double HypoTestInverterResult::GetGraphX(const TGraph & graph, double y0, bool lowSearch, double &axmin, double &axmax) const {
553  // return the X value of the given graph for the target value y0
554  // the graph is evaluated using linear interpolation by default.
555  // if option = "S" a TSpline3 is used
556 
557 
558 
559  TString opt = "";
560  if (fInterpolOption == kSpline) opt = "S";
561 
562  InterpolatedGraph f(graph,y0,opt);
565 
566  // find reasanable xmin and xmax if not given
567  const double * y = graph.GetY();
568  int n = graph.GetN();
569  if (n < 2) {
570  ooccoutE(this,Eval) << "HypoTestInverterResult::GetGraphX - need at least 2 points for interpolation (n=" << n << ")\n";
571  return (n>0) ? y[0] : 0;
572  }
573 
574  double varmin = - TMath::Infinity();
575  double varmax = TMath::Infinity();
576  const RooRealVar* var = dynamic_cast<RooRealVar*>( fParameters.first() );
577  if (var) {
578  varmin = var->getMin();
579  varmax = var->getMax();
580  }
581 
582  double xmin = axmin;
583  double xmax = axmax;
584  // case no range is given check if need to extrapolate to lower/upper values
585  if (axmin >= axmax ) {
586  xmin = graph.GetX()[0];
587  xmax = graph.GetX()[n-1];
588 
589  // case we have lower/upper limits
590 
591  // find ymin and ymax and corresponding values
592  double ymin = TMath::MinElement(n,y);
593  double ymax = TMath::MaxElement(n,y);
594  // do lower extrapolation
595  if ( (ymax < y0 && !lowSearch) || ( ymin > y0 && lowSearch) ) {
596  xmin = varmin;
597  }
598  // do upper extrapolation
599  if ( (ymax < y0 && lowSearch) || ( ymin > y0 && !lowSearch) ) {
600  xmax = varmax;
601  }
602  }
603  else {
604 
605 #ifdef ISREALLYNEEDED //??
606  // in case of range given, check if all points are above or below the confidence level line
607  bool isCross = false;
608  bool first = true;
609  double prod = 1;
610  for (int i = 0; i< n; ++i) {
611  double xp, yp;
612  graph.GetPoint(i,xp,yp);
613  if (xp >= xmin && xp <= xmax) {
614  prod *= TMath::Sign(1., (yp - y0) );
615  if (prod < 0 && !first) {
616  isCross = true;
617  break;
618  }
619  first = false;
620  }
621  }
622  if (!isCross) {
623  return (lowSearch) ? varmin : varmax;
624  }
625 #endif
626  }
627 
628 
629 
630  brf.SetFunction(wf,xmin,xmax);
631  brf.SetNpx(20);
632  bool ret = brf.Solve(100, 1.E-16, 1.E-6);
633  if (!ret) {
634  ooccoutE(this,Eval) << "HypoTestInverterResult - interpolation failed - return inf" << std::endl;
635  return TMath::Infinity();
636  }
637  double limit = brf.Root();
638 
639 //#define DO_DEBUG
640 #ifdef DO_DEBUG
641  if (lowSearch) std::cout << "lower limit search : ";
642  else std::cout << "Upper limit search : ";
643  std::cout << "do interpolation between " << xmin << " and " << xmax << " limit is " << limit << std::endl;
644 #endif
645 
646  // look in case if a new interseption exists
647  // only when boundaries are not given
648  if (axmin >= axmax) {
649  int index = TMath::BinarySearch(n, graph.GetX(), limit);
650 #ifdef DO_DEBUG
651  std::cout << "do new interpolation dividing from " << index << " and " << y[index] << std::endl;
652 #endif
653 
654  if (lowSearch && index >= 1 && (y[0] - y0) * ( y[index]- y0) < 0) {
655  //search if another interseption exists at a lower value
656  limit = GetGraphX(graph, y0, lowSearch, graph.GetX()[0], graph.GetX()[index] );
657  }
658  else if (!lowSearch && index < n-2 && (y[n-1] - y0) * ( y[index+1]- y0) < 0) {
659  // another interseption exists at an higher value
660  limit = GetGraphX(graph, y0, lowSearch, graph.GetX()[index+1], graph.GetX()[n-1] );
661  }
662  }
663  // return also xmin, xmax values
664  axmin = xmin;
665  axmax = xmax;
666 
667  return limit;
668 }
669 
670 double HypoTestInverterResult::FindInterpolatedLimit(double target, bool lowSearch, double xmin, double xmax )
671 {
672  // interpolate to find a limit value
673  // Use a linear or a spline interpolation depending on the interpolation option
674 
675  ooccoutD(this,Eval) << "HypoTestInverterResult - "
676  << "Interpolate the upper limit between the 2 results closest to the target confidence level"
677  << std::endl;
678 
679  // variable minimum and maximum
680  double varmin = - TMath::Infinity();
681  double varmax = TMath::Infinity();
682  const RooRealVar* var = dynamic_cast<RooRealVar*>( fParameters.first() );
683  if (var) {
684  varmin = var->getMin();
685  varmax = var->getMax();
686  }
687 
688  if (ArraySize()<2) {
689  double val = (lowSearch) ? xmin : xmax;
690  oocoutW(this,Eval) << "HypoTestInverterResult::FindInterpolatedLimit"
691  << " - not enough points to get the inverted interval - return "
692  << val << std::endl;
693  fLowerLimit = varmin;
694  fUpperLimit = varmax;
695  return (lowSearch) ? fLowerLimit : fUpperLimit;
696  }
697 
698  // sort the values in x
699  int n = ArraySize();
700  std::vector<unsigned int> index(n );
701  TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);
702  // make a graph with the sorted point
703  TGraph graph(n);
704  for (int i = 0; i < n; ++i)
705  graph.SetPoint(i, GetXValue(index[i]), GetYValue(index[i] ) );
706 
707 
708  //std::cout << " search for " << lowSearch << std::endl;
709 
710 
711  // search first for min/max in the given range
712  if (xmin >= xmax) {
713 
714 
715  // search for maximum between the point
716  double * itrmax = std::max_element(graph.GetY() , graph.GetY() +n);
717  double ymax = *itrmax;
718  int iymax = itrmax - graph.GetY();
719  double xwithymax = graph.GetX()[iymax];
720 
721  //std::cout << " max of y " << iymax << " " << xwithymax << " " << ymax << " target is " << target << std::endl;
722 
723  // look if maximum is above/belove target
724  if (ymax > target) {
725  if (lowSearch) {
726  if ( iymax > 0) {
727  // low search (minimmum is first point or minimum range)
728  xmin = ( graph.GetY()[0] <= target ) ? graph.GetX()[0] : varmin;
729  xmax = xwithymax;
730  }
731  else {
732  // no room for lower limit
733  fLowerLimit = varmin;
734  lowSearch = false; // search now for upper limits
735  }
736  }
737  if (!lowSearch ) {
738  // up search
739  if ( iymax < n-1 ) {
740  xmin = xwithymax;
741  xmax = ( graph.GetY()[n-1] <= target ) ? graph.GetX()[n-1] : varmax;
742  }
743  else {
744  // no room for upper limit
745  fUpperLimit = varmax;
746  lowSearch = true; // search now for lower limits
747  xmin = varmin;
748  xmax = xwithymax;
749  }
750  }
751  }
752  else {
753  // in case is below the target
754  // find out if is a lower or upper search
755  if (iymax <= (n-1)/2 ) {
756  lowSearch = false;
757  fLowerLimit = varmin;
758  }
759  else {
760  lowSearch = true;
761  fUpperLimit = varmax;
762  }
763  }
764 
765 #ifdef DO_DEBUG
766  std::cout << " found xmin, xmax = " << xmin << " " << xmax << " for search " << lowSearch << std::endl;
767 #endif
768 
769  // now come here if I have already found a lower/upper limit
770  // i.e. I am calling routine for the second time
771 #ifdef ISNEEDED
772  // should not really come here
773  if (lowSearch && fUpperLimit < varmax) {
774  xmin = fXValues[ index.front() ];
775  // find xmax (is first point before upper limit)
776  int upI = FindClosestPointIndex(target, 2, fUpperLimit);
777  if (upI < 1) return xmin;
778  xmax = GetXValue(upI);
779  }
780  else if (!lowSearch && fLowerLimit > varmin ) {
781  // find xmin (is first point after lower limit)
782  int lowI = FindClosestPointIndex(target, 3, fLowerLimit);
783  if (lowI >= n-1) return xmax;
784  xmin = GetXValue(lowI);
785  xmax = fXValues[ index.back() ];
786  }
787 #endif
788  }
789 
790 #ifdef DO_DEBUG
791  std::cout << "finding " << lowSearch << " limit betweem " << xmin << " " << xmax << endl;
792 #endif
793 
794  double limit = GetGraphX(graph, target, lowSearch, xmin, xmax);
795  if (lowSearch) fLowerLimit = limit;
796  else fUpperLimit = limit;
797  CalculateEstimatedError( target, lowSearch, xmin, xmax);
798 
799 #ifdef DO_DEBUG
800  std::cout << "limit is " << limit << std::endl;
801 #endif
802 
803  if (lowSearch && !TMath::IsNaN(fUpperLimit)) return fLowerLimit;
804  if (!lowSearch && !TMath::IsNaN(fLowerLimit)) return fUpperLimit;
805 
806  // now perform the opposite search on the complement interval
807  if (lowSearch) {
808  xmin = xmax;
809  xmax = varmax;
810  } else {
811  xmax = xmin;
812  xmin = varmin;
813  }
814  double limit2 = GetGraphX(graph, target, !lowSearch, xmin, xmax);
815  if (!lowSearch) fLowerLimit = limit2;
816  else fUpperLimit = limit2;
817 
818  CalculateEstimatedError( target, !lowSearch, xmin, xmax);
819 
820 #ifdef DO_DEBUG
821  std::cout << "other limit is " << limit2 << std::endl;
822 #endif
823 
824  return (lowSearch) ? fLowerLimit : fUpperLimit;
825 
826 }
827 
828 
829 int HypoTestInverterResult::FindClosestPointIndex(double target, int mode, double xtarget)
830 {
831  // if mode = 0
832  // find closest point to target in Y, the object closest to the target which is 3 sigma from the target
833  // and has smaller error
834  // if mode = 1
835  // find 2 closest point to target in X and between these two take the one closer to the target
836  // if mode = 2 as in mode = 1 but return the lower point not the closest one
837  // if mode = 3 as in mode = 1 but return the upper point not the closest one
838 
839  int bestIndex = -1;
840  int closestIndex = -1;
841  if (mode == 0) {
842  double smallestError = 2; // error must be < 1
843  double bestValue = 2;
844  for (int i=0; i<ArraySize(); i++) {
845  double dist = fabs(GetYValue(i)-target);
846  if ( dist <3 *GetYError(i) ) { // less than 1 sigma from target CL
847  if (GetYError(i) < smallestError ) {
848  smallestError = GetYError(i);
849  bestIndex = i;
850  }
851  }
852  if ( dist < bestValue) {
853  bestValue = dist;
854  closestIndex = i;
855  }
856  }
857  if (bestIndex >=0) return bestIndex;
858  // if no points found just return the closest one to the target
859  return closestIndex;
860  }
861  // else mode = 1,2,3
862  // find the two closest points to limit value
863  // sort the array first
864  int n = fXValues.size();
865  std::vector<unsigned int> indx(n);
866  TMath::SortItr(fXValues.begin(), fXValues.end(), indx.begin(), false);
867  std::vector<double> xsorted( n);
868  for (int i = 0; i < n; ++i) xsorted[i] = fXValues[indx[i] ];
869  int index1 = TMath::BinarySearch( n, &xsorted[0], xtarget);
870 
871 #ifdef DO_DEBUG
872  std::cout << "finding closest point to " << xtarget << " is " << index1 << " " << indx[index1] << std::endl;
873 #endif
874 
875  // case xtarget is outside the range (bbefore or afterwards)
876  if (index1 < 0) return indx[0];
877  if (index1 >= n-1) return indx[n-1];
878  int index2 = index1 +1;
879 
880  if (mode == 2) return (GetXValue(indx[index1]) < GetXValue(indx[index2])) ? indx[index1] : indx[index2];
881  if (mode == 3) return (GetXValue(indx[index1]) > GetXValue(indx[index2])) ? indx[index1] : indx[index2];
882  // get smaller point of the two (mode == 1)
883  if (fabs(GetYValue(indx[index1])-target) <= fabs(GetYValue(indx[index2])-target) )
884  return indx[index1];
885  return indx[index2];
886 
887 }
888 
889 
891 {
892  if (fFittedLowerLimit) return fLowerLimit;
893  //std::cout << "finding point with cl = " << 1-(1-ConfidenceLevel())/2 << endl;
894  if ( fInterpolateLowerLimit ) {
895  // find both lower/upper limit
897  } else {
898  //LM: I think this is never called
900  }
901  return fLowerLimit;
902 }
903 
905 {
906  //std::cout << "finding point with cl = " << (1-ConfidenceLevel())/2 << endl;
907  if (fFittedUpperLimit) return fUpperLimit;
908  if ( fInterpolateUpperLimit ) {
910  } else {
911  //LM: I think this is never called
913  }
914  return fUpperLimit;
915 }
916 
917 Double_t HypoTestInverterResult::CalculateEstimatedError(double target, bool lower, double xmin, double xmax)
918 {
919  // Return an error estimate on the upper(lower) limit. This is the error on
920  // either CLs or CLsplusb divided by an estimate of the slope at this
921  // point.
922 
923  if (ArraySize()==0) {
924  oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimateError"
925  << "Empty result \n";
926  return 0;
927  }
928 
929  if (ArraySize()<2) {
930  oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimateError"
931  << " only points - return its error\n";
932  return GetYError(0);
933  }
934 
935  // it does not make sense in case of asymptotic which do not have point errors
936  if (!GetNullTestStatDist(0) ) return 0;
937 
938  TString type = (!lower) ? "upper" : "lower";
939 
940 #ifdef DO_DEBUG
941  std::cout << "calculate estimate error " << type << " between " << xmin << " and " << xmax << std::endl;
942  std::cout << "computed limit is " << ( (lower) ? fLowerLimit : fUpperLimit ) << std::endl;
943 #endif
944 
945  // make a TGraph Errors with the sorted points
946  std::vector<unsigned int> indx(fXValues.size());
947  TMath::SortItr(fXValues.begin(), fXValues.end(), indx.begin(), false);
948  // make a graph with the sorted point
949  TGraphErrors graph;
950  int ip = 0, np = 0;
951  for (int i = 0; i < ArraySize(); ++i) {
952  if ( (xmin < xmax) && ( GetXValue(indx[i]) >= xmin && GetXValue(indx[i]) <= xmax) ) {
953  np++;
954  // exclude points with zero or very small errors
955  if (GetYError(indx[i] ) > 1.E-6) {
956  graph.SetPoint(ip, GetXValue(indx[i]), GetYValue(indx[i] ) );
957  graph.SetPointError(ip, 0., GetYError(indx[i]) );
958  ip++;
959  }
960  }
961  }
962  if (graph.GetN() < 2) {
963  if (np >= 2) oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimatedError - no valid points - cannot estimate the " << type << " limit error " << std::endl;
964  return 0;
965  }
966 
967  double minX = xmin;
968  double maxX = xmax;
969  if (xmin >= xmax) {
970  minX = fXValues[ indx.front() ];
971  maxX = fXValues[ indx.back() ];
972  }
973 
974 
975 
976  TF1 fct("fct", "exp([0] * (x - [2] ) + [1] * (x-[2])**2)", minX, maxX);
977  double scale = maxX-minX;
978  if (lower) {
979  fct.SetParameters( 2./scale, 0.1/scale, graph.GetX()[0] );
980  fct.SetParLimits(0,0,100./scale);
981  fct.SetParLimits(1,0, 10./scale); }
982  else {
983  fct.SetParameters( -2./scale, -0.1/scale );
984  fct.SetParLimits(0,-100./scale, 0);
985  fct.SetParLimits(1,-100./scale, 0); }
986 
987  if (graph.GetN() < 3) fct.FixParameter(1,0.);
988 
989  // find the point closest to the limit
990  double limit = (!lower) ? fUpperLimit : fLowerLimit;
991  if (TMath::IsNaN(limit)) return 0; // cannot do if limit not computed
992 
993 
994 #ifdef DO_DEBUG
995  TCanvas * c1 = new TCanvas();
996  std::cout << "fitting for limit " << type << "between " << minX << " , " << maxX << " points considered " << graph.GetN() << std::endl;
997  int fitstat = graph.Fit(&fct," EX0");
998  graph.SetMarkerStyle(20);
999  graph.Draw("AP");
1000  graph.Print();
1001  c1->SaveAs(TString::Format("graphFit_%s.pdf",type.Data()) );
1002  delete c1;
1003 #else
1004  int fitstat = graph.Fit(&fct,"Q EX0");
1005 #endif
1006 
1007  int index = FindClosestPointIndex(target, 1, limit);
1008  double theError = 0;
1009  if (fitstat == 0) {
1010  double errY = GetYError(index);
1011  if (errY > 0) {
1012  double m = fct.Derivative( GetXValue(index) );
1013  theError = std::min(fabs( GetYError(index) / m), maxX-minX);
1014  }
1015  }
1016  else {
1017  oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimatedError - cannot estimate the " << type << " limit error " << std::endl;
1018  theError = 0;
1019  }
1020  if (lower)
1021  fLowerLimitError = theError;
1022  else
1023  fUpperLimitError = theError;
1024 
1025 #ifdef DO_DEBUG
1026  std::cout << "closes point to the limit is " << index << " " << GetXValue(index) << " and has error " << GetYError(index) << std::endl;
1027 #endif
1028 
1029  return theError;
1030 }
1031 
1032 
1034 {
1035  // need to have compute first lower limit
1037  if (fLowerLimitError >= 0) return fLowerLimitError;
1038  // try to recompute the error
1039  return CalculateEstimatedError(1-ConfidenceLevel(), true);
1040 }
1041 
1042 
1044 {
1046  if (fUpperLimitError >= 0) return fUpperLimitError;
1047  // try to recompute the error
1048  return CalculateEstimatedError(1-ConfidenceLevel(), false);
1049 }
1050 
1052  // get the background test statistic distribution
1053 
1054  HypoTestResult * firstResult = (HypoTestResult*) fYObjects.At(index);
1055  if (!firstResult) return 0;
1056  return firstResult->GetBackGroundIsAlt() ? firstResult->GetAltDistribution() : firstResult->GetNullDistribution();
1057 }
1058 
1060  // get the signal and background test statistic distribution
1062  if (!result) return 0;
1063  return !result->GetBackGroundIsAlt() ? result->GetAltDistribution() : result->GetNullDistribution();
1064 }
1065 
1067  // get the expected p-value distribution at the scanned point index
1068 
1069  if (index < 0 || index >= ArraySize() ) return 0;
1070 
1071  if (fExpPValues.GetSize() == ArraySize() ) {
1072  return (SamplingDistribution*) fExpPValues.At(index)->Clone();
1073  }
1074 
1075  static bool useFirstB = false;
1076  // get the S+B distribution
1077  int bIndex = (useFirstB) ? 0 : index;
1078 
1079  SamplingDistribution * bDistribution = GetBackgroundTestStatDist(bIndex);
1080  SamplingDistribution * sbDistribution = GetSignalAndBackgroundTestStatDist(index);
1081 
1083 
1084  if (bDistribution && sbDistribution) {
1085 
1086  // create a new HypoTestResult
1087  HypoTestResult tempResult;
1088  tempResult.SetPValueIsRightTail( result->GetPValueIsRightTail() );
1089  tempResult.SetBackgroundAsAlt( true);
1090  // ownership of SamplingDistribution is in HypoTestResult class now
1091  tempResult.SetNullDistribution( new SamplingDistribution(*sbDistribution) );
1092  tempResult.SetAltDistribution( new SamplingDistribution(*bDistribution ) );
1093 
1094  std::vector<double> values(bDistribution->GetSize());
1095  for (int i = 0; i < bDistribution->GetSize(); ++i) {
1096  tempResult.SetTestStatisticData( bDistribution->GetSamplingDistribution()[i] );
1097  values[i] = (fUseCLs) ? tempResult.CLs() : tempResult.CLsplusb();
1098  }
1099  return new SamplingDistribution("expected values","expected values",values);
1100  }
1101  // in case b abs sbDistribution are null assume is coming from the asymptotic calculator
1102  // hard -coded this value (no really needed to be used by user)
1103  fgAsymptoticMaxSigma = 5;
1104  const int npoints = 11;
1105  const double smax = fgAsymptoticMaxSigma;
1106  const double dsig = 2* smax/ (npoints-1) ;
1107  std::vector<double> values(npoints);
1108  for (int i = 0; i < npoints; ++i) {
1109  double nsig = -smax + dsig*i;
1110  double pval = AsymptoticCalculator::GetExpectedPValues( result->NullPValue(), result->AlternatePValue(), nsig, fUseCLs, !fIsTwoSided);
1111  if (pval < 0) { return 0;}
1112  values[i] = pval;
1113  }
1114  return new SamplingDistribution("Asymptotic expected values","Asymptotic expected values",values);
1115 
1116 }
1117 
1118 
1119 
1121  // get the limit distribution (lower/upper depending on the flag)
1122  // by interpolating the expected p values for each point
1123 
1124 
1125  if (ArraySize()<2) {
1126  oocoutE(this,Eval) << "HypoTestInverterResult::GetLimitDistribution"
1127  << " not enought points - return 0 " << std::endl;
1128  return 0;
1129  }
1130 
1131  ooccoutD(this,Eval) << "HypoTestInverterResult - computing limit distribution...." << std::endl;
1132 
1133 
1134  // find optimal size by looking at the PValue distribution obtained
1135  int npoints = ArraySize();
1136  std::vector<SamplingDistribution*> distVec( npoints );
1137  double sum = 0;
1138  for (unsigned int i = 0; i < distVec.size(); ++i) {
1139  distVec[i] = GetExpectedPValueDist(i);
1140  // sort the distributions
1141  // hack (by calling InverseCDF(0) will sort the sampling distribution
1142  if (distVec[i] ) {
1143  distVec[i]->InverseCDF(0);
1144  sum += distVec[i]->GetSize();
1145  }
1146  }
1147  int size = int( sum/ npoints);
1148 
1149  if (size < 10) {
1150  ooccoutW(this,InputArguments) << "HypoTestInverterResult - set a minimum size of 10 for limit distribution" << std::endl;
1151  size = 10;
1152  }
1153 
1154 
1155  double target = 1-fConfidenceLevel;
1156 
1157  // vector with the quantiles of the p-values for each scanned poi point
1158  std::vector< std::vector<double> > quantVec(npoints );
1159  for (int i = 0; i < npoints; ++i) {
1160 
1161  if (!distVec[i]) continue;
1162 
1163  // make quantiles from the sampling distributions of the expected p values
1164  std::vector<double> pvalues = distVec[i]->GetSamplingDistribution();
1165  delete distVec[i]; distVec[i] = 0;
1166  std::sort(pvalues.begin(), pvalues.end());
1167  // find the quantiles of the distribution
1168  double p[1] = {0};
1169  double q[1] = {0};
1170 
1171  quantVec[i] = std::vector<double>(size);
1172  for (int ibin = 0; ibin < size; ++ibin) {
1173  // exclude for a bug in TMath::Quantiles last item
1174  p[0] = std::min( (ibin+1) * 1./double(size), 1.0);
1175  // use the type 1 which give the point value
1176  TMath::Quantiles(pvalues.size(), 1, &pvalues[0], q, p, true, 0, 1 );
1177  (quantVec[i])[ibin] = q[0];
1178  }
1179 
1180  }
1181 
1182  // sort the values in x
1183  std::vector<unsigned int> index( npoints );
1184  TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);
1185 
1186  // SamplingDistribution * dist0 = distVec[index[0]];
1187  // SamplingDistribution * dist1 = distVec[index[1]];
1188 
1189 
1190  std::vector<double> limits(size);
1191  // loop on the p values and find the limit for each expected point in the quantiles vector
1192  for (int j = 0; j < size; ++j ) {
1193 
1194  TGraph g;
1195  for (int k = 0; k < npoints ; ++k) {
1196  if (quantVec[index[k]].size() > 0 )
1197  g.SetPoint(g.GetN(), GetXValue(index[k]), (quantVec[index[k]])[j] );
1198  }
1199 
1200  limits[j] = GetGraphX(g, target, lower);
1201 
1202  }
1203 
1204 
1205  if (lower)
1206  return new SamplingDistribution("Expected lower Limit","Expected lower limits",limits);
1207  else
1208  return new SamplingDistribution("Expected upper Limit","Expected upper limits",limits);
1209 
1210 }
1211 
1212 
1213 double HypoTestInverterResult::GetExpectedLowerLimit(double nsig, const char * opt ) const {
1214  // Get the expected lower limit
1215  // nsig is used to specify which expected value of the UpperLimitDistribution
1216  // For example
1217  // nsig = 0 (default value) returns the expected value
1218  // nsig = -1 returns the lower band value at -1 sigma
1219  // nsig + 1 return the upper value
1220  // opt = "" (default) : compute limit by interpolating all the p values, find the corresponding limit distribution
1221  // and then find the quantiles in the limit distribution
1222  // ioption = "P" is the method used for plotting. One Finds the corresponding nsig quantile in the p values and then
1223  // interpolates them
1224 
1225  return GetExpectedLimit(nsig, true, opt);
1226 }
1227 
1228 double HypoTestInverterResult::GetExpectedUpperLimit(double nsig, const char * opt ) const {
1229  // Get the expected upper limit
1230  // nsig is used to specify which expected value of the UpperLimitDistribution
1231  // For example
1232  // nsig = 0 (default value) returns the expected value
1233  // nsig = -1 returns the lower band value at -1 sigma
1234  // nsig + 1 return the upper value
1235  // opt is an option specifying the type of method used for computing the upper limit
1236  // opt = "" (default) : compute limit by interpolating all the p values, find the corresponding limit distribution
1237  // and then find the quantiles in the limit distribution
1238  // ioption = "P" is the method used for plotting. One Finds the corresponding nsig quantile in the p values and then
1239  // interpolates them
1240  //
1241 
1242  return GetExpectedLimit(nsig, false, opt);
1243 }
1244 
1245 
1246 
1247 double HypoTestInverterResult::GetExpectedLimit(double nsig, bool lower, const char * opt ) const {
1248  // get expected limit (lower/upper) depending on the flag
1249  // for asymptotic is a special case (the distribution is generated an step in sigma values)
1250  // distringuish asymptotic looking at the hypotest results
1251  // if option = "P" get expected limit using directly quantiles of p value distribution
1252  // else (default) find expected limit by obtaining first a full limit distributions
1253  // The last one is in general more correct
1254 
1255  const int nEntries = ArraySize();
1256  if (nEntries <= 0) return (lower) ? 1 : 0; // return 1 for lower, 0 for upper
1257 
1258  HypoTestResult * r = dynamic_cast<HypoTestResult *> (fYObjects.First() );
1259  assert(r != 0);
1260  if (!r->GetNullDistribution() && !r->GetAltDistribution() ) {
1261  // we are in the asymptotic case
1262  // get the limits obtained at the different sigma values
1263  SamplingDistribution * limitDist = GetLimitDistribution(lower);
1264  if (!limitDist) return 0;
1265  const std::vector<double> & values = limitDist->GetSamplingDistribution();
1266  if (values.size() <= 1) return 0;
1267  double dsig = 2* fgAsymptoticMaxSigma/ (values.size() -1) ;
1268  int i = (int) TMath::Floor ( (nsig + fgAsymptoticMaxSigma)/dsig + 0.5);
1269  return values[i];
1270  }
1271 
1272  double p[1] = {0};
1273  double q[1] = {0};
1274  p[0] = ROOT::Math::normal_cdf(nsig,1);
1275 
1276  // for CLs+b can get the quantiles of p-value distribution and
1277  // interpolate them
1278  // In the case of CLs (since it is not a real p-value anymore but a ratio)
1279  // then it is needed to get first all limit distribution values and then the quantiles
1280 
1281  TString option(opt);
1282  option.ToUpper();
1283  if (option.Contains("P")) {
1284 
1285  TGraph g;
1286 
1287  // sort the arrays based on the x values
1288  std::vector<unsigned int> index(nEntries);
1289  TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);
1290 
1291  for (int j=0; j<nEntries; ++j) {
1292  int i = index[j]; // i is the order index
1294  if (!s) {
1295  ooccoutI(this,Eval) << "HypoTestInverterResult - cannot compute expected p value distribution for point, x = "
1296  << GetXValue(i) << " skip it " << std::endl;
1297  continue;
1298  }
1299  const std::vector<double> & values = s->GetSamplingDistribution();
1300  double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
1301  TMath::Quantiles(values.size(), 1, x,q,p,false);
1302  g.SetPoint(g.GetN(), fXValues[i], q[0] );
1303  delete s;
1304  }
1305  if (g.GetN() < 2) {
1306  ooccoutE(this,Eval) << "HypoTestInverterResult - cannot compute limits , not enough points, n = " << g.GetN() << std::endl;
1307  return 0;
1308  }
1309 
1310  // interpolate the graph to obtain the limit
1311  double target = 1-fConfidenceLevel;
1312  return GetGraphX(g, target, lower);
1313 
1314  }
1315  // here need to use the limit distribution
1316  SamplingDistribution * limitDist = GetLimitDistribution(lower);
1317  if (!limitDist) return 0;
1318  const std::vector<double> & values = limitDist->GetSamplingDistribution();
1319  double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
1320  TMath::Quantiles(values.size(), 1, x,q,p,false);
1321  return q[0];
1322 
1323 }
1324 
SamplingDistribution * GetLimitDistribution(bool lower) const
void SetBackgroundAsAlt(Bool_t l=kTRUE)
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
void SetAltDistribution(SamplingDistribution *alt)
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:439
float xmin
Definition: THbookFile.cxx:93
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:404
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
#define ooccoutI(o, a)
Definition: RooMsgService.h:52
Double_t Floor(Double_t x)
Definition: TMath.h:473
virtual Double_t ConfidenceLevel() const
return confidence level
bool Solve(int maxIter=100, double absTol=1E-8, double relTol=1E-10)
Returns the X value corresponding to the function value fy for (xmin
T1 Sign(T1 a, T2 b)
Definition: TMathBase.h:155
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 ...
double CLbError(int index) const
return the observed CLb value for the i-th entry
HypoTestInverterResult(const char *name=0)
default constructor
int ArraySize() const
number of entries in the results array
virtual Double_t CLs() const
CLs is simply CLs+b/CLb (not a method, but a quantity)
TCanvas * c1
Definition: legend1.C:2
int FindClosestPointIndex(double target, int mode=0, double xtarget=0)
Double_t QuietNaN()
Definition: TMath.h:635
#define assert(cond)
Definition: unittest.h:542
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
Base class for spline implementation containing the Draw/Paint methods //.
Definition: TSpline.h:22
#define oocoutI(o, a)
Definition: RooMsgService.h:45
HypoTestResult is a base class for results from hypothesis tests.
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1101
SamplingDistribution * GetSignalAndBackgroundTestStatDist(int index) const
Basic string class.
Definition: TString.h:137
virtual Double_t getMin(const char *name=0) const
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...
SamplingDistribution * GetAltDistribution(void) const
Int_t GetSize() const
size of samples
virtual Double_t CLsplusb() const
Convert AlternatePValue into a "confidence level".
Int_t GetN() const
Definition: TGraph.h:132
SamplingDistribution * GetBackgroundTestStatDist(int index) const
#define ooccoutW(o, a)
Definition: RooMsgService.h:54
const std::vector< Double_t > & GetSamplingDistribution() const
Get test statistics values.
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:310
virtual TObject * RemoveAt(Int_t idx)
#define ooccoutE(o, a)
Definition: RooMsgService.h:55
Template class to wrap any C++ callable object which takes one argument i.e.
const char * Data() const
Definition: TString.h:349
Double_t * GetY() const
Definition: TGraph.h:140
SimpleInterval & operator=(const SimpleInterval &other)
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TObject.cxx:203
Double_t x[n]
Definition: legend1.C:17
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:2334
double FindInterpolatedLimit(double target, bool lowSearch=false, double xmin=1, double xmax=0)
double CLsError(int index) const
return the observed CLb value for the i-th entry
#define oocoutE(o, a)
Definition: RooMsgService.h:48
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=0, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities Parameters: x -the data sample n ...
Definition: TMath.cxx:1173
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
virtual void SaveAs(const char *filename="", Option_t *option="") const
Save Pad contents in a file in one of various formats.
Definition: TPad.cxx:5009
static Vc_ALWAYS_INLINE Vector< T > abs(const Vector< T > &x)
Definition: vector.h:450
virtual Double_t NullPValue() const
Return p-value for null hypothesis.
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:63
Double_t LowerLimit()
lower and upper bound of the confidence interval (to get upper/lower limits, multiply the size( = 1-c...
bool fInterpolateLowerLimit
two sided scan (look for lower/upper limit)
Double_t Infinity()
Definition: TMath.h:648
int ExclusionCleanup()
remove points that appear to have failed.
double GetGraphX(const TGraph &g, double y0, bool lowSearch, double &xmin, double &xmax) const
HypoTestResult * GetResult(int index) const
return a pointer to the i^th result object
double CLsplusb(int index) const
return the observed CLsplusb value for the i-th entry
double Root() const
Returns root value.
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Definition: TMath.h:192
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Definition: TMath.h:196
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set limits for parameter ipar.
Definition: TF1.cxx:3206
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
const double tol
Double_t * GetX() const
Definition: TGraph.h:139
SamplingDistribution * GetExpectedPValueDist(int index) const
return expected distribution of p-values (Cls or Clsplusb)
double GetExpectedLowerLimit(double nsig=0, const char *opt="") const
get Limit value correspnding at the desired nsigma level (0) is median -1 sigma is 1 sigma ...
TMarker * m
Definition: textangle.C:8
virtual void FixParameter(Int_t ipar, Double_t value)
Fix the value of a parameter The specified value will be used in a fit operation. ...
Definition: TF1.cxx:1276
Double_t E()
Definition: TMath.h:54
TList fExpPValues
list of HypoTestResult for each point
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
float xmax
Definition: THbookFile.cxx:93
Bool_t GetPValueIsRightTail(void) const
void SetTestStatisticData(const Double_t tsd)
This class simply holds a sampling distribution of some test statistic.
ClassImp(RooStats::HypoTestInverterResult) using namespace RooStats
virtual void Append(const HypoTestResult *other)
add values from another HypoTestResult
Class for finding the root of a one dimensional function using the Brent algorithm.
void SetPValueIsRightTail(Bool_t pr)
virtual Double_t Derivative(Double_t x, Double_t *params=0, Double_t epsilon=0.001) const
Returns the first derivative of the function at point x, computed by Richardson's extrapolation metho...
Definition: TF1.cxx:852
The Canvas class.
Definition: TCanvas.h:48
double GetXValue(int index) const
function to return the value of the parameter of interest for the i^th entry in the results ...
TRObject operator()(const T1 &t1) const
virtual Int_t GetSize() const
Definition: TCollection.h:95
double CLs(int index) const
return the observed CLb value for the i-th entry
double f(double x)
HypoTestInverterResult class: holds the array of hypothesis test results and compute a confidence int...
double Double_t
Definition: RtypesCore.h:55
Double_t LowerLimitEstimatedError()
rough estimation of the error on the computed bound of the confidence interval Estimate of lower limi...
int type
Definition: TGX11.cxx:120
double GetYValue(int index) const
function to return the value of the confidence level for the i^th entry in the results ...
T MaxElement(Long64_t n, const T *a)
Definition: TMath.h:688
SamplingDistribution * GetNullDistribution(void) const
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition: TGraph.cxx:1551
#define ooccoutD(o, a)
Definition: RooMsgService.h:51
double CLsplusbError(int index) const
return the observed CLsplusb value for the i-th entry
virtual Double_t AlternatePValue() const
Return p-value for alternate hypothesis.
#define oocoutW(o, a)
Definition: RooMsgService.h:47
virtual Double_t getMax(const char *name=0) const
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:556
Double_t UpperLimitEstimatedError()
Estimate of lower limit error function evaluates only a rought error on the lower limit...
virtual TObject * clone(const char *newname) const
Definition: RooRealVar.h:48
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition: TGraph.cxx:2127
double GetExpectedUpperLimit(double nsig=0, const char *opt="") const
get Limit value correspnding at the desired nsigma level (0) is median -1 sigma is 1 sigma ...
HypoTestInverterResult & operator=(const HypoTestInverterResult &other)
operator =
virtual void Add(TObject *obj)
Definition: TList.h:81
double fLowerLimitError
interpolatation option (linear or spline)
virtual void RemoveAll(TCollection *col)
Remove all objects in collection col from this collection.
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
1-Dim function class
Definition: TF1.h:149
Int_t IsNaN(Double_t x)
Definition: TMath.h:617
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
A TGraphErrors is a TGraph with error bars.
Definition: TGraphErrors.h:28
double GetExpectedLimit(double nsig, bool lower, const char *opt="") const
bool Add(const HypoTestInverterResult &otherResult)
merge with the content of another HypoTestInverterResult object
double result[121]
bool SetFunction(const ROOT::Math::IGenFunction &f, double xlow, double xup)
Sets the function for the rest of the algorithms.
void SetNpx(int npx)
Set the number of point used to bracket root using a grid.
double CalculateEstimatedError(double target, bool lower=true, double xmin=1, double xmax=0)
float * q
Definition: THbookFile.cxx:87
double CLb(int index) const
return the observed CLb value for the i-th entry
const Int_t n
Definition: legend1.C:16
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Definition: TMath.h:944
static const float lower
Definition: main.cpp:48
void SetNullDistribution(SamplingDistribution *null)
void SortItr(Iterator first, Iterator last, IndexIterator index, Bool_t down=kTRUE)
Definition: TMath.h:977
SamplingDistribution * GetNullTestStatDist(int index) const
same in terms of alt and null
std::vector< double > fXValues
max sigma value used to scan asymptotic expected p values
T MinElement(Long64_t n, const T *a)
Definition: TMath.h:681
Bool_t GetBackGroundIsAlt(void) const