00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include "RooAbsPdf.h"
00030 #include "RooAbsData.h"
00031 #include "RooRealVar.h"
00032 #include "TMath.h"
00033
00034 #include "RooStats/HybridCalculatorOriginal.h"
00035 #include "RooStats/HybridResult.h"
00036
00037
00038 #include "RooStats/HypoTestInverter.h"
00039
00040
00041 ClassImp(RooStats::HypoTestInverter)
00042
00043 using namespace RooStats;
00044
00045
00046 HypoTestInverter::HypoTestInverter( ) :
00047 fCalculator0(0),
00048 fScannedVariable(0),
00049 fResults(0),
00050 fUseCLs(false),
00051 fSize(0)
00052 {
00053
00054 }
00055
00056
00057 HypoTestInverter::HypoTestInverter( HypoTestCalculator& myhc0,
00058 RooRealVar& scannedVariable, double size ) :
00059 TNamed( ),
00060 fCalculator0(&myhc0),
00061 fScannedVariable(&scannedVariable),
00062 fResults(0),
00063 fUseCLs(false),
00064 fSize(size)
00065 {
00066
00067
00068 SetName("HypoTestInverter");
00069
00070
00071 HybridCalculatorOriginal * hc = dynamic_cast<HybridCalculatorOriginal *> (fCalculator0);
00072 if (hc == 0) {
00073 Fatal("HypoTestInverter","Using non HybridCalculatorOriginal class IS NOT SUPPORTED");
00074 }
00075
00076 }
00077
00078
00079 HypoTestInverter::~HypoTestInverter()
00080 {
00081
00082
00083
00084 if (fResults) delete fResults;
00085 }
00086
00087 void HypoTestInverter::CreateResults() {
00088
00089 if (fResults == 0) {
00090 TString results_name = this->GetName();
00091 results_name += "_results";
00092 fResults = new HypoTestInverterResult(results_name,*fScannedVariable,ConfidenceLevel());
00093 fResults->SetTitle("HypoTestInverter Result");
00094 }
00095 fResults->UseCLs(fUseCLs);
00096 }
00097
00098
00099 bool HypoTestInverter::RunAutoScan( double xMin, double xMax, double target, double epsilon, unsigned int numAlgorithm )
00100 {
00101
00102
00103
00104
00105
00106
00107
00108 if ( xMin>=xMax || xMin< fScannedVariable->getMin() || xMax>fScannedVariable->getMax() ) {
00109 std::cout << "Error: problem with the specified range\n";
00110 return false;
00111 }
00112 if ( target<=0 || target>=1 ) {
00113 std::cout << "Error: problem with target value\n";
00114 return false;
00115 }
00116 if ( epsilon>0.5-fabs(0.5-target) ) {
00117 std::cout << "Error: problem with error value\n";
00118 return false;
00119 }
00120 if ( numAlgorithm!=0 && numAlgorithm!=1 ) {
00121 std::cout << "Error: invalid interpolation algorithm\n";
00122 return false;
00123 }
00124
00125 CreateResults();
00126
00127
00128 if ( fabs(1-target/(1-Size()/2))<DBL_EPSILON ) {
00129 fResults->fInterpolateLowerLimit = false;
00130 std::cout << "Target matches lower limit: de-activate interpolation in HypoTestInverterResult\n";
00131 }
00132
00133 if ( fabs(1-target/((Size()/2)))<DBL_EPSILON ) {
00134 fResults->fInterpolateUpperLimit = false;
00135 std::cout << "Target matches upper limit: de-activate interpolation in HypoTestInverterResult\n";
00136 }
00137
00138
00139 const double nSigma = 1;
00140
00141
00142 const unsigned int nToys_backup = ((HybridCalculatorOriginal*)fCalculator0)->GetNumberOfToys();
00143
00144
00145 double leftX = xMin;
00146 if (!RunOnePoint(leftX)) return false;
00147 double leftCL = fResults->GetYValue(fResults->ArraySize()-1);
00148 double leftCLError = fResults->GetYError(fResults->ArraySize()-1);
00149
00150 double rightX = xMax;
00151 if (!RunOnePoint(rightX)) return false;
00152 double rightCL = fResults->GetYValue(fResults->ArraySize()-1);
00153 double rightCLError = fResults->GetYError(fResults->ArraySize()-1);
00154
00155 if ( rightCL>target && leftCL>target ) {
00156 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";
00157 return false;
00158 }
00159 if ( rightCL<target && leftCL<target ) {
00160 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";
00161 return false;
00162 }
00163
00164 unsigned int nIteration = 2;
00165 bool quitThisLoop = false;
00166
00167 double centerCL = 0;
00168 double centerCLError = 0;
00169
00170
00171
00172
00173 do {
00174 double x = 0;
00175
00176
00177 if (leftCL==rightCL) {
00178 std::cout << "This cannot (and should not) happen... quit\n";
00179 quitThisLoop = true;
00180 } else if (leftX==rightX) {
00181 std::cout << "This cannot (and should not) happen... quit\n";
00182 quitThisLoop = true;
00183 } else {
00184
00185
00186 if (numAlgorithm==0) {
00187
00188
00189
00190 if (!leftCL) leftCL = DBL_EPSILON;
00191 if (!rightCL) rightCL = DBL_EPSILON;
00192
00193 double a = (log(leftCL) - log(rightCL)) / (leftX - rightX);
00194 double b = leftCL / exp(a * leftX);
00195 x = (log(target) - log(b)) / a;
00196
00197
00198 if (x<xMin || x>xMax || isnan(x)) {
00199 std::cout << "Extrapolated value out of range or nan: exits\n";
00200 quitThisLoop = true;
00201 }
00202 } else if (numAlgorithm==1) {
00203
00204
00205 double a = (leftCL-rightCL)/(leftX-rightX);
00206 double b = leftCL-a*leftX;
00207 x = (target-b)/a;
00208
00209 if (x<xMin || x>xMax || isnan(x)) {
00210 std::cout << "Extrapolated value out of range or nan: exits\n";
00211 quitThisLoop = true;
00212 }
00213 }
00214 }
00215
00216 if ( x==leftX || x==rightX ) {
00217 std::cout << "Error: exit because interpolated value equals to a previous iteration\n";
00218 quitThisLoop = true;
00219 }
00220
00221
00222 bool success = false;
00223 if (!quitThisLoop) success = RunOnePoint(x);
00224
00225 if (success) {
00226
00227 nIteration++;
00228 centerCL = fResults->GetYValue(fResults->ArraySize()-1);
00229 centerCLError = fResults->GetYError(fResults->ArraySize()-1);
00230
00231
00232
00233
00234
00235 if ( (leftCL > target) == (rightCL < target) ) {
00236 if ( (centerCL > target) == (leftCL > target) ) {
00237 leftX = x;
00238 leftCL = centerCL;
00239 leftCLError = centerCLError;
00240 } else {
00241 rightX = x;
00242 rightCL = centerCL;
00243 rightCLError = centerCLError;
00244 }
00245
00246
00247 } else if ( (fabs(leftCL - target) / leftCLError) >
00248 (fabs(rightCL - target) / rightCLError) ) {
00249 leftX = x;
00250 leftCL = centerCL;
00251 leftCLError = centerCLError;
00252 } else {
00253 rightX = x;
00254 rightCL = centerCL;
00255 rightCLError = centerCLError;
00256 }
00257
00258
00259
00260 if ( fabs(centerCL-target) < nSigma*centerCLError && centerCLError > epsilon ) {
00261 do {
00262
00263 int nToys = ((HybridCalculatorOriginal*)fCalculator0)->GetNumberOfToys();
00264 int nToysTarget = (int) TMath::Max(nToys*1.5, 1.2*nToys*pow(centerCLError/epsilon,2));
00265
00266 std::cout << "Increasing the number of toys to: " << nToysTarget << std::endl;
00267
00268
00269 ((HybridCalculatorOriginal*)fCalculator0)->SetNumberOfToys(nToysTarget-nToys);
00270
00271 if (!RunOnePoint(x)) quitThisLoop=true;
00272 nIteration++;
00273 centerCL = fResults->GetYValue(fResults->ArraySize()-1);
00274 centerCLError = fResults->GetYError(fResults->ArraySize()-1);
00275
00276
00277 ((HybridCalculatorOriginal*)fCalculator0)->SetNumberOfToys(nToysTarget);
00278 } while ( fabs(centerCL-target) < nSigma*centerCLError && centerCLError > epsilon && quitThisLoop==false );
00279 }
00280
00281 if (leftCL==rightCL) {
00282 std::cout << "Algorithm failed: left and right CL are equal (no intrapolation possible or more toy-MC statistics needed)\n";
00283 quitThisLoop = true;
00284 }
00285 }
00286
00287 } while ( ( fabs(centerCL-target) > nSigma*centerCLError || centerCLError > epsilon ) && quitThisLoop==false );
00288
00289
00290 ((HybridCalculatorOriginal*)fCalculator0)->SetNumberOfToys(nToys_backup);
00291
00292 if ( quitThisLoop==true ) {
00293
00294 std::cout << "Aborted the search because something happened\n";
00295 return false;
00296 }
00297
00298 std::cout << "Converged in " << fResults->ArraySize() << " iterations\n";
00299
00300
00301 return true;
00302 }
00303
00304
00305 bool HypoTestInverter::RunFixedScan( int nBins, double xMin, double xMax )
00306 {
00307
00308
00309 CreateResults();
00310
00311 if ( nBins<=0 ) {
00312 std::cout << "Please provide nBins>0\n";
00313 return false;
00314 }
00315 if ( nBins==1 && xMin!=xMax ) {
00316 std::cout << "nBins==1 -> I will run for xMin (" << xMin << ")\n";
00317 }
00318 if ( xMin==xMax && nBins>1 ) {
00319 std::cout << "xMin==xMax -> I will enforce nBins==1\n";
00320 nBins = 1;
00321 }
00322 if ( xMin>xMax ) {
00323 std::cout << "Please provide xMin (" << xMin << ") smaller that xMax (" << xMax << ")\n";
00324 return false;
00325 }
00326
00327 for (int i=0; i<nBins; i++) {
00328 double thisX = xMin+i*(xMax-xMin)/(nBins-1);
00329 bool status = RunOnePoint(thisX);
00330
00331
00332 if ( status==false ) {
00333 std::cout << "Loop interupted because of failed status\n";
00334 return false;
00335 }
00336 }
00337
00338 return true;
00339 }
00340
00341
00342 bool HypoTestInverter::RunOnePoint( double thisX )
00343 {
00344
00345
00346 CreateResults();
00347
00348
00349 if ( thisX<fScannedVariable->getMin() ) {
00350 std::cout << "Out of range: using the lower bound on the scanned variable rather than " << thisX<< "\n";
00351 thisX = fScannedVariable->getMin();
00352 }
00353 if ( thisX>fScannedVariable->getMax() ) {
00354 std::cout << "Out of range: using the upper bound on the scanned variable rather than " << thisX<< "\n";
00355 thisX = fScannedVariable->getMax();
00356 }
00357
00358 double oldValue = fScannedVariable->getVal();
00359
00360 fScannedVariable->setVal(thisX);
00361 std::cout << "Running for " << fScannedVariable->GetName() << " = " << thisX << endl;
00362
00363
00364 HypoTestResult* myHybridResult = fCalculator0->GetHypoTest();
00365
00366 double lastXtested;
00367 if ( fResults->ArraySize()!=0 ) lastXtested = fResults->GetXValue(fResults->ArraySize()-1);
00368 else lastXtested = -999;
00369
00370 if ( lastXtested==thisX ) {
00371
00372 std::cout << "Merge with previous result\n";
00373 HybridResult* latestResult = (HybridResult*) fResults->GetResult(fResults->ArraySize()-1);
00374 latestResult->Add((HybridResult*)myHybridResult);
00375 delete myHybridResult;
00376
00377 } else {
00378
00379
00380 fResults->fXValues.push_back(thisX);
00381 fResults->fYObjects.Add(myHybridResult);
00382 }
00383
00384
00385 std::cout << "computed: " << fResults->GetYValue(fResults->ArraySize()-1) << endl;
00386
00387 fScannedVariable->setVal(oldValue);
00388
00389 return true;
00390 }