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