ROOT  6.06/09
Reference Guide
HypoTestInverterOriginal.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 
15 // include other header files
16 
17 #include "RooAbsPdf.h"
18 #include "RooAbsData.h"
19 #include "RooRealVar.h"
20 #include "TMath.h"
21 
23 #include "RooStats/HybridResult.h"
24 
25 // include header file of this class
27 
28 
30 
31 using namespace RooStats;
32 using namespace std;
33 
34 
36  fCalculator0(0),
37  fScannedVariable(0),
38  fResults(0),
39  fUseCLs(false),
40  fSize(0)
41 {
42  // default constructor (doesn't do anything)
43 }
44 
45 
46 HypoTestInverterOriginal::HypoTestInverterOriginal( HypoTestCalculator& myhc0,
47  RooRealVar& scannedVariable, double size ) :
48  TNamed( ),
49  fCalculator0(&myhc0),
50  fScannedVariable(&scannedVariable),
51  fResults(0),
52  fUseCLs(false),
53  fSize(size)
54 {
55  // constructor from a reference to an HypoTestCalculator
56  // (it must be an HybridCalculator type) and a RooRealVar for the variable
57  SetName("HypoTestInverterOriginal");
58 
59 
61  if (hc == 0) {
62  Fatal("HypoTestInverterOriginal","Using non HybridCalculatorOriginal class IS NOT SUPPORTED");
63  }
64 
65 }
66 
67 
69 {
70  // destructor
71 
72  // delete the HypoTestInverterResult
73  if (fResults) delete fResults;
74 }
75 
77  // create a new HypoTestInverterResult to hold all computed results
78  if (fResults == 0) {
79  TString results_name = this->GetName();
80  results_name += "_results";
82  fResults->SetTitle("HypoTestInverterOriginal Result");
83  }
85 }
86 
87 
88 bool HypoTestInverterOriginal::RunAutoScan( double xMin, double xMax, double target, double epsilon, unsigned int numAlgorithm )
89 {
90  /// Search for the value of the parameter of interest (vary the
91  /// hypothesis being tested) in the specified range [xMin,xMax]
92  /// until the confidence level is compatible with the target value
93  /// within one time the estimated error (and the estimated error
94  /// should also become smaller than the specified parameter epsilon)
95 
96  // various sanity checks on the input parameters
97  if ( xMin>=xMax || xMin< fScannedVariable->getMin() || xMax>fScannedVariable->getMax() ) {
98  std::cout << "Error: problem with the specified range\n";
99  return false;
100  }
101  if ( target<=0 || target>=1 ) {
102  std::cout << "Error: problem with target value\n";
103  return false;
104  }
105  if ( epsilon>0.5-fabs(0.5-target) ) {
106  std::cout << "Error: problem with error value\n";
107  return false;
108  }
109  if ( numAlgorithm!=0 && numAlgorithm!=1 ) {
110  std::cout << "Error: invalid interpolation algorithm\n";
111  return false;
112  }
113 
114  CreateResults();
115 
116  // if ( TMath::AreEqualRel(target,1-Size()/2,DBL_EPSILON) ) { // to uncomment for ROOT 5.26
117  if ( fabs(1-target/(1-Size()/2))<DBL_EPSILON ) {
119  std::cout << "Target matches lower limit: de-activate interpolation in HypoTestInverterResult\n";
120  }
121  // if ( TMath::AreEqualRel(target,Size()/2,DBL_EPSILON) ) { // to uncomment for ROOT 5.26
122  if ( fabs(1-target/((Size()/2)))<DBL_EPSILON ) {
124  std::cout << "Target matches upper limit: de-activate interpolation in HypoTestInverterResult\n";
125  }
126 
127  // parameters of the algorithm that are hard-coded
128  const double nSigma = 1; // number of times the estimated error the final p-value should be from the target
129 
130  // backup some values to be restored at the end
131  const unsigned int nToys_backup = ((HybridCalculatorOriginal*)fCalculator0)->GetNumberOfToys();
132 
133  // check the 2 hypothesis tests specified as extrema in the constructor
134  double leftX = xMin;
135  if (!RunOnePoint(leftX)) return false;
136  double leftCL = fResults->GetYValue(fResults->ArraySize()-1);
137  double leftCLError = fResults->GetYError(fResults->ArraySize()-1);
138 
139  double rightX = xMax;
140  if (!RunOnePoint(rightX)) return false;
141  double rightCL = fResults->GetYValue(fResults->ArraySize()-1);
142  double rightCLError = fResults->GetYError(fResults->ArraySize()-1);
143 
144  if ( rightCL>target && leftCL>target ) {
145  std::cout << "The confidence level at both boundaries are both too large ( " << leftCL << " and " << rightCL << std::endl << "Run again with other boundaries or larger toy-MC statistics\n";
146  return false;
147  }
148  if ( rightCL<target && leftCL<target ) {
149  std::cout << "The confidence level at both boundaries are both too small ( " << leftCL << " and " << rightCL << std::endl << "Run again with other boundaries or larger toy-MC statistics\n";
150  return false;
151  }
152 
153  unsigned int nIteration = 2; // number of iteration performed by the algorithm
154  bool quitThisLoop = false; // flag to interrupt the search and quit cleanly
155 
156  double centerCL = 0;
157  double centerCLError = 0;
158 
159  // search for the value of the searched variable where the CL is
160  // within 1 sigma of the desired level and sigma smaller than
161  // epsilon.
162  do {
163  double x = 0;
164 
165  // safety checks
166  if (leftCL==rightCL) {
167  std::cout << "This cannot (and should not) happen... quit\n";
168  quitThisLoop = true;
169  } else if (leftX==rightX) {
170  std::cout << "This cannot (and should not) happen... quit\n";
171  quitThisLoop = true;
172  } else {
173 
174  // apply chosen type of interpolation algorithm
175  if (numAlgorithm==0) {
176  // exponential interpolation
177 
178  // add safety checks
179  if (!leftCL) leftCL = DBL_EPSILON;
180  if (!rightCL) rightCL = DBL_EPSILON;
181 
182  double a = (log(leftCL) - log(rightCL)) / (leftX - rightX);
183  double b = leftCL / exp(a * leftX);
184  x = (log(target) - log(b)) / a;
185 
186  // to do: do not allow next iteration outside the xMin,xMax interval
187  if (x<xMin || x>xMax || TMath::IsNaN(x)) {
188  std::cout << "Extrapolated value out of range or nan: exits\n";
189  quitThisLoop = true;
190  }
191  } else if (numAlgorithm==1) {
192  // linear interpolation
193 
194  double a = (leftCL-rightCL)/(leftX-rightX);
195  double b = leftCL-a*leftX;
196  x = (target-b)/a;
197 
198  if (x<xMin || x>xMax || TMath::IsNaN(x)) {
199  std::cout << "Extrapolated value out of range or nan: exits\n";
200  quitThisLoop = true;
201  }
202  } // end of interpolation algorithms
203  }
204 
205  if ( x==leftX || x==rightX ) {
206  std::cout << "Error: exit because interpolated value equals to a previous iteration\n";
207  quitThisLoop = true;
208  }
209 
210  // perform another hypothesis-test for value x
211  bool success = false;
212  if (!quitThisLoop) success = RunOnePoint(x);
213 
214  if (success) {
215 
216  nIteration++; // succeeded, increase the iteration counter
217  centerCL = fResults->GetYValue(fResults->ArraySize()-1);
218  centerCLError = fResults->GetYError(fResults->ArraySize()-1);
219 
220  // replace either the left or right point by this new point
221 
222  // test if the interval points are on different sides, then
223  // replace the one on the correct side with the center
224  if ( (leftCL > target) == (rightCL < target) ) {
225  if ( (centerCL > target) == (leftCL > target) ) {
226  leftX = x;
227  leftCL = centerCL;
228  leftCLError = centerCLError;
229  } else {
230  rightX = x;
231  rightCL = centerCL;
232  rightCLError = centerCLError;
233  }
234  // Otherwise replace the point farest away from target (measured in
235  // sigmas)
236  } else if ( (fabs(leftCL - target) / leftCLError) >
237  (fabs(rightCL - target) / rightCLError) ) {
238  leftX = x;
239  leftCL = centerCL;
240  leftCLError = centerCLError;
241  } else {
242  rightX = x;
243  rightCL = centerCL;
244  rightCLError = centerCLError;
245  }
246 
247  // if a point is found compatible with the target CL but with too
248  // large error, increase the number of toyMC
249  if ( fabs(centerCL-target) < nSigma*centerCLError && centerCLError > epsilon ) {
250  do {
251 
252  int nToys = ((HybridCalculatorOriginal*)fCalculator0)->GetNumberOfToys(); // current number of toys
253  int nToysTarget = (int) TMath::Max(nToys*1.5, 1.2*nToys*pow(centerCLError/epsilon,2)); // estimated number of toys until the target precision is reached
254 
255  std::cout << "Increasing the number of toys to: " << nToysTarget << std::endl;
256 
257  // run again the same point with more toyMC (run the complement number of toys)
258  ((HybridCalculatorOriginal*)fCalculator0)->SetNumberOfToys(nToysTarget-nToys);
259 
260  if (!RunOnePoint(x)) quitThisLoop=true;
261  nIteration++; // succeeded, increase the iteration counter
262  centerCL = fResults->GetYValue(fResults->ArraySize()-1);
263  centerCLError = fResults->GetYError(fResults->ArraySize()-1);
264 
265  // set the number of toys to reach the target
266  ((HybridCalculatorOriginal*)fCalculator0)->SetNumberOfToys(nToysTarget);
267  } while ( fabs(centerCL-target) < nSigma*centerCLError && centerCLError > epsilon && quitThisLoop==false ); // run this block again if it's still compatible with the target and the error still too large
268  }
269 
270  if (leftCL==rightCL) {
271  std::cout << "Algorithm failed: left and right CL are equal (no intrapolation possible or more toy-MC statistics needed)\n";
272  quitThisLoop = true;
273  }
274  } // end running one more iteration
275 
276  } while ( ( fabs(centerCL-target) > nSigma*centerCLError || centerCLError > epsilon ) && quitThisLoop==false ); // end of the main 'do' loop
277 
278  // restore some parameters that might have been changed by the algorithm
279  ((HybridCalculatorOriginal*)fCalculator0)->SetNumberOfToys(nToys_backup);
280 
281  if ( quitThisLoop==true ) {
282  // abort and return 'false' to indicate fail status
283  std::cout << "Aborted the search because something happened\n";
284  return false;
285  }
286 
287  std::cout << "Converged in " << fResults->ArraySize() << " iterations\n";
288 
289  // finished: return 'true' for success status
290  return true;
291 }
292 
293 
294 bool HypoTestInverterOriginal::RunFixedScan( int nBins, double xMin, double xMax )
295 {
296  // Run a Fixed scan in npoints between min and max
297 
298  CreateResults();
299  // safety checks
300  if ( nBins<=0 ) {
301  std::cout << "Please provide nBins>0\n";
302  return false;
303  }
304  if ( nBins==1 && xMin!=xMax ) {
305  std::cout << "nBins==1 -> I will run for xMin (" << xMin << ")\n";
306  }
307  if ( xMin==xMax && nBins>1 ) {
308  std::cout << "xMin==xMax -> I will enforce nBins==1\n";
309  nBins = 1;
310  }
311  if ( xMin>xMax ) {
312  std::cout << "Please provide xMin (" << xMin << ") smaller that xMax (" << xMax << ")\n";
313  return false;
314  }
315 
316  for (int i=0; i<nBins; i++) {
317  double thisX = xMin+i*(xMax-xMin)/(nBins-1);
318  bool status = RunOnePoint(thisX);
319 
320  // check if failed status
321  if ( status==false ) {
322  std::cout << "Loop interrupted because of failed status\n";
323  return false;
324  }
325  }
326 
327  return true;
328 }
329 
330 
332 {
333  // run only one point
334 
335  CreateResults();
336 
337  // check if thisX is in the range specified for fScannedVariable
338  if ( thisX<fScannedVariable->getMin() ) {
339  std::cout << "Out of range: using the lower bound on the scanned variable rather than " << thisX<< "\n";
340  thisX = fScannedVariable->getMin();
341  }
342  if ( thisX>fScannedVariable->getMax() ) {
343  std::cout << "Out of range: using the upper bound on the scanned variable rather than " << thisX<< "\n";
344  thisX = fScannedVariable->getMax();
345  }
346 
347  double oldValue = fScannedVariable->getVal();
348 
349  fScannedVariable->setVal(thisX);
350  std::cout << "Running for " << fScannedVariable->GetName() << " = " << thisX << endl;
351 
352  // compute the results
353  HypoTestResult* myHybridResult = fCalculator0->GetHypoTest();
354 
355  double lastXtested;
356  if ( fResults->ArraySize()!=0 ) lastXtested = fResults->GetXValue(fResults->ArraySize()-1);
357  else lastXtested = -999;
358 
359  if ( lastXtested==thisX ) {
360 
361  std::cout << "Merge with previous result\n";
362  HybridResult* latestResult = (HybridResult*) fResults->GetResult(fResults->ArraySize()-1);
363  latestResult->Add((HybridResult*)myHybridResult);
364  delete myHybridResult;
365 
366  } else {
367 
368  // fill the results in the HypoTestInverterResult array
369  fResults->fXValues.push_back(thisX);
370  fResults->fYObjects.Add(myHybridResult);
371  }
372 
373 
374  std::cout << "computed: " << fResults->GetYValue(fResults->ArraySize()-1) << endl;
375 
376  fScannedVariable->setVal(oldValue);
377 
378  return true;
379 }
bool RunAutoScan(double xMin, double xMax, double target, double epsilon=0.005, unsigned int numAlgorithm=0)
int ArraySize() const
number of entries in the results array
void Add(HybridResult *other)
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
HypoTestResult is a base class for results from hypothesis tests.
Basic string class.
Definition: TString.h:137
virtual Double_t getMin(const char *name=0) const
HybridCalculatorOriginal class.
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
TArc * a
Definition: textangle.C:12
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.
ClassImp(TIterator) Bool_t TIterator return false
Compare two iterator objects.
Definition: TIterator.cxx:20
ClassImp(RooStats::HypoTestInverterOriginal) using namespace RooStats
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:946
Double_t x[n]
Definition: legend1.C:17
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
double pow(double, double)
void UseCLs(bool on=true)
flag to switch between using CLsb (default) or CLs as confidence level
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
bool fInterpolateLowerLimit
two sided scan (look for lower/upper limit)
HypoTestCalculator is an interface class for a tools which produce RooStats HypoTestResults.
HypoTestResult * GetResult(int index) const
return a pointer to the i^th result object
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:202
virtual Double_t Size() const
Get the size of the test (eg. rate of Type I error)
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
REAL epsilon
Definition: triangle.c:617
virtual Double_t ConfidenceLevel() const
Get the Confidence level for the test.
virtual HypoTestResult * GetHypoTest() const =0
double GetXValue(int index) const
function to return the value of the parameter of interest for the i^th entry in the results ...
Namespace for the RooStats classes.
Definition: Asimov.h:20
This class is now depratcated and to be replaced by the HypoTestInverter.
HypoTestInverterResult class: holds the array of hypothesis test results and compute a confidence int...
double GetYValue(int index) const
function to return the value of the confidence level for the i^th entry in the results ...
virtual Double_t getMax(const char *name=0) const
virtual void Add(TObject *obj)
Definition: TList.h:81
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Int_t IsNaN(Double_t x)
Definition: TMath.h:617
bool RunFixedScan(int nBins, double xMin, double xMax)
double exp(double)
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
double log(double)
std::vector< double > fXValues
max sigma value used to scan asymptotic expected p values