Go to the documentation of this file.00001 #include <stdexcept>
00002 #include "../interface/FeldmanCousins.h"
00003 #include "../interface/Combine.h"
00004 #include <RooRealVar.h>
00005 #include <RooArgSet.h>
00006 #include <RooArgSet.h>
00007 #include <RooWorkspace.h>
00008 #include <RooDataHist.h>
00009 #include <RooStats/ModelConfig.h>
00010 #include <RooStats/FeldmanCousins.h>
00011 #include <RooStats/PointSetInterval.h>
00012
00013 float FeldmanCousins::toysFactor_ = 1;
00014 float FeldmanCousins::rAbsAccuracy_ = 0.1;
00015 float FeldmanCousins::rRelAccuracy_ = 0.02;
00016
00017 FeldmanCousins::FeldmanCousins() :
00018 LimitAlgo("FeldmanCousins specific options") {
00019 options_.add_options()
00020 ("rAbsAcc", boost::program_options::value<float>(&rAbsAccuracy_)->default_value(rAbsAccuracy_), "Absolute accuracy on r to reach to terminate the scan")
00021 ("rRelAcc", boost::program_options::value<float>(&rRelAccuracy_)->default_value(rRelAccuracy_), "Relative accuracy on r to reach to terminate the scan")
00022 ("toysFactor", boost::program_options::value<float>(&toysFactor_)->default_value(toysFactor_), "Increase the toys per point by this factor w.r.t. the minimum from adaptive sampling")
00023 ;
00024 }
00025
00026 void FeldmanCousins::applyOptions(const boost::program_options::variables_map &vm)
00027 {
00028 }
00029
00030 bool FeldmanCousins::run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) {
00031 RooArgSet poi(*mc_s->GetParametersOfInterest());
00032 RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
00033
00034 if ((hint != 0) && (*hint > r->getMin())) {
00035 r->setMax(std::min<double>(3*(*hint), r->getMax()));
00036 }
00037
00038 RooStats::ModelConfig modelConfig(*mc_s);
00039 modelConfig.SetSnapshot(poi);
00040
00041 RooStats::FeldmanCousins fc(data, modelConfig);
00042 fc.FluctuateNumDataEntries(mc_s->GetPdf()->canBeExtended());
00043 fc.UseAdaptiveSampling(true);
00044 fc.SetConfidenceLevel(cl);
00045 fc.AdditionalNToysFactor(toysFactor_);
00046
00047 fc.SetNBins(10);
00048 do {
00049 if (verbose > 1) std::cout << "scan in range [" << r->getMin() << ", " << r->getMax() << "]" << std::endl;
00050 std::auto_ptr<RooStats::PointSetInterval> fcInterval((RooStats::PointSetInterval *)fc.GetInterval());
00051 if (fcInterval.get() == 0) return false;
00052 if (verbose > 1) fcInterval->GetParameterPoints()->Print("V");
00053 RooDataHist* parameterScan = (RooDataHist*) fc.GetPointsToScan();
00054 int found = -1;
00055 if (lowerLimit_) {
00056 for(int i=0; i<parameterScan->numEntries(); ++i){
00057 const RooArgSet * tmpPoint = parameterScan->get(i);
00058 if (!fcInterval->IsInInterval(*tmpPoint)) found = i;
00059 else break;
00060 }
00061 } else {
00062 for(int i=0; i<parameterScan->numEntries(); ++i){
00063 const RooArgSet * tmpPoint = parameterScan->get(i);
00064 bool inside = fcInterval->IsInInterval(*tmpPoint);
00065 if (inside) found = i;
00066 else if (found != -1) break;
00067 }
00068 if (found == -1) {
00069 if (verbose) std::cout << "Points are either all inside or all outside the bound." << std::endl;
00070 return false;
00071 }
00072 }
00073 double fcBefore = (found > -1 ? parameterScan->get(found)->getRealValue(r->GetName()) : r->getMin());
00074 double fcAfter = (found < parameterScan->numEntries()-1 ?
00075 parameterScan->get(found+1)->getRealValue(r->GetName()) : r->getMax());
00076 limit = 0.5*(fcAfter+fcBefore);
00077 limitErr = 0.5*(fcAfter-fcBefore);
00078 if (verbose > 0) std::cout << " would be " << r->GetName() << " < " << limit << " +/- "<<limitErr << std::endl;
00079 r->setMin(std::max(r->getMin(), limit-3*limitErr));
00080 r->setMax(std::min(r->getMax(), limit+3*limitErr));
00081 if (limitErr < 4*std::max<float>(rAbsAccuracy_, rRelAccuracy_ * limit)) {
00082 fc.AdditionalNToysFactor(4*toysFactor_);
00083 }
00084 } while (limitErr > std::max<float>(rAbsAccuracy_, rRelAccuracy_ * limit));
00085
00086 if (verbose > -1) {
00087 std::cout << "\n -- FeldmanCousins++ -- \n";
00088 std::cout << "Limit: " << r->GetName() << (lowerLimit_ ? "> " : "< ") << limit << " +/- " << limitErr << " @ " << cl * 100 << "% CL" << std::endl;
00089 }
00090 return true;
00091 }