00001 #include "../interface/FitterAlgoBase.h"
00002 #include <limits>
00003 #include <cmath>
00004
00005 #include "RooRealVar.h"
00006 #include "RooArgSet.h"
00007 #include "RooRandom.h"
00008 #include "RooDataSet.h"
00009 #include "RooFitResult.h"
00010 #include "RooSimultaneous.h"
00011 #include "RooAddPdf.h"
00012 #include "RooProdPdf.h"
00013 #include "RooConstVar.h"
00014 #include "RooPlot.h"
00015 #include "../interface/RooMinimizerOpt.h"
00016 #include "TCanvas.h"
00017 #include "TStyle.h"
00018 #include "TH2.h"
00019 #include "TFile.h"
00020 #include <RooStats/ModelConfig.h>
00021 #include "../interface/Combine.h"
00022 #include "../interface/ProfileLikelihood.h"
00023 #include "../interface/CascadeMinimizer.h"
00024 #include "../interface/CloseCoutSentry.h"
00025 #include "../interface/utils.h"
00026
00027 #include "../interface/ProfilingTools.h"
00028
00029
00030 #include <Math/MinimizerOptions.h>
00031 #include <Math/QuantFuncMathCore.h>
00032 #include <Math/ProbFunc.h>
00033
00034 using namespace RooStats;
00035
00036 std::string FitterAlgoBase::minimizerAlgo_ = "Minuit2";
00037 std::string FitterAlgoBase::minimizerAlgoForMinos_ = "Minuit2,simplex";
00038 #include <Math/QuantFuncMathCore.h>
00039 #include <Math/ProbFunc.h>
00040 float FitterAlgoBase::minimizerTolerance_ = 1e-2;
00041 float FitterAlgoBase::minimizerToleranceForMinos_ = 1e-4;
00042 int FitterAlgoBase::minimizerStrategy_ = 1;
00043 int FitterAlgoBase::minimizerStrategyForMinos_ = 0;
00044 float FitterAlgoBase::preFitValue_ = 1.0;
00045 float FitterAlgoBase::stepSize_ = 0.1;
00046 bool FitterAlgoBase::robustFit_ = false;
00047 int FitterAlgoBase::maxFailedSteps_ = 5;
00048 bool FitterAlgoBase::do95_ = false;
00049 bool FitterAlgoBase::saveNLL_ = false;
00050 bool FitterAlgoBase::keepFailures_ = false;
00051 float FitterAlgoBase::nllValue_ = std::numeric_limits<float>::quiet_NaN();
00052
00053 FitterAlgoBase::FitterAlgoBase(const char *title) :
00054 LimitAlgo(title)
00055 {
00056 options_.add_options()
00057 ("minimizerAlgo", boost::program_options::value<std::string>(&minimizerAlgo_)->default_value(minimizerAlgo_), "Choice of minimizer (Minuit vs Minuit2)")
00058 ("minimizerTolerance", boost::program_options::value<float>(&minimizerTolerance_)->default_value(minimizerTolerance_), "Tolerance for minimizer")
00059 ("minimizerStrategy", boost::program_options::value<int>(&minimizerStrategy_)->default_value(minimizerStrategy_), "Stragegy for minimizer")
00060 ("preFitValue", boost::program_options::value<float>(&preFitValue_)->default_value(preFitValue_), "Value of signal strength pre-fit")
00061 ("do95", boost::program_options::value<bool>(&do95_)->default_value(do95_), "Compute also 2-sigma interval from delta(nll) = 1.92 instead of 0.5")
00062 ("robustFit", boost::program_options::value<bool>(&robustFit_)->default_value(robustFit_), "Search manually for 1 and 2 sigma bands instead of using Minos")
00063 ("maxFailedSteps", boost::program_options::value<int>(&maxFailedSteps_)->default_value(maxFailedSteps_), "How many failed steps to retry before giving up")
00064 ("stepSize", boost::program_options::value<float>(&stepSize_)->default_value(stepSize_), "Step size for robust fits (multiplier of the range)")
00065 ("minimizerAlgoForMinos", boost::program_options::value<std::string>(&minimizerAlgoForMinos_)->default_value(minimizerAlgoForMinos_), "Choice of minimizer (Minuit vs Minuit2) for profiling in robust fits")
00066 ("minimizerStrategyForMinos", boost::program_options::value<int>(&minimizerStrategyForMinos_)->default_value(minimizerStrategyForMinos_), "Stragegy for minimizer for profiling in robust fits")
00067 ("minimizerToleranceForMinos", boost::program_options::value<float>(&minimizerToleranceForMinos_)->default_value(minimizerToleranceForMinos_), "Tolerance for minimizer for profiling in robust fits")
00068 ("saveNLL", "Save the negative log-likelihood at the minimum in the output tree (note: value is relative to the pre-fit state)")
00069 ("keepFailures", "Save the results even if the fit is declared as failed (for NLL studies)")
00070 ;
00071 }
00072
00073 void FitterAlgoBase::applyOptionsBase(const boost::program_options::variables_map &vm)
00074 {
00075 saveNLL_ = vm.count("saveNLL");
00076 keepFailures_ = vm.count("keepFailures");
00077 }
00078
00079 bool FitterAlgoBase::run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) {
00080 ProfileLikelihood::MinimizerSentry minimizerConfig(minimizerAlgo_, minimizerTolerance_);
00081 CloseCoutSentry sentry(verbose < 0);
00082
00083 static bool shouldCreateNLLBranch = saveNLL_;
00084 if (shouldCreateNLLBranch) { Combine::addBranch("nll", &nllValue_, "nll/F"); shouldCreateNLLBranch = false; }
00085
00086 return runSpecific(w, mc_s, mc_b, data, limit, limitErr, hint);
00087 }
00088
00089
00090 RooFitResult *FitterAlgoBase::doFit(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &r, const RooCmdArg &constrain, bool doHesse, int ndim, bool reuseNLL) {
00091 return doFit(pdf, data, RooArgList (r), constrain, doHesse, ndim, reuseNLL);
00092 }
00093
00094 RooFitResult *FitterAlgoBase::doFit(RooAbsPdf &pdf, RooAbsData &data, const RooArgList &rs, const RooCmdArg &constrain, bool doHesse, int ndim, bool reuseNLL) {
00095 RooFitResult *ret = 0;
00096 if (reuseNLL) nll->setData(data);
00097 else nll.reset(pdf.createNLL(data, constrain, RooFit::Extended(pdf.canBeExtended())));
00098
00099 double nll0 = nll->getVal();
00100 double delta68 = 0.5*ROOT::Math::chisquared_quantile_c(1-0.68,ndim);
00101 double delta95 = 0.5*ROOT::Math::chisquared_quantile_c(1-0.95,ndim);
00102 CascadeMinimizer minim(*nll, CascadeMinimizer::Unconstrained, rs.getSize() ? dynamic_cast<RooRealVar*>(rs.first()) : 0);
00103 minim.setStrategy(minimizerStrategy_);
00104 minim.setErrorLevel(delta68);
00105 CloseCoutSentry sentry(verbose < 3);
00106 bool ok = minim.minimize(verbose-1);
00107 nllValue_ = nll->getVal() - nll0;
00108 if (!ok && !keepFailures_) { std::cout << "Initial minimization failed. Aborting." << std::endl; return 0; }
00109 if (doHesse) minim.minimizer().hesse();
00110 sentry.clear();
00111 ret = minim.save();
00112 if (verbose > 1) { ret->Print("V"); }
00113
00114 std::auto_ptr<RooArgSet> allpars(pdf.getParameters(data));
00115
00116 for (int i = 0, n = rs.getSize(); i < n; ++i) {
00117
00118 if (i) {
00119 RooArgSet oldparams(ret->floatParsFinal());
00120 *allpars = oldparams;
00121 }
00122
00123
00124 RooRealVar &r = dynamic_cast<RooRealVar &>(*rs.at(i));
00125 RooRealVar &rf = dynamic_cast<RooRealVar &>(*ret->floatParsFinal().find(r.GetName()));
00126 double r0 = r.getVal(), rMin = r.getMin(), rMax = r.getMax();
00127
00128 if (!robustFit_) {
00129 if (do95_) {
00130 throw std::runtime_error("95% CL errors with Minos are not working at the moment.");
00131 minim.setErrorLevel(delta95);
00132 minim.improve(verbose-1);
00133 minim.setErrorLevel(delta95);
00134 if (minim.minimizer().minos(RooArgSet(r)) != -1) {
00135 rf.setRange("err95", r.getVal() + r.getAsymErrorLo(), r.getVal() + r.getAsymErrorHi());
00136 }
00137 minim.setErrorLevel(delta68);
00138 minim.improve(verbose-1);
00139 }
00140 if (minim.minimizer().minos(RooArgSet(r)) != -1) {
00141 rf.setRange("err68", r.getVal() + r.getAsymErrorLo(), r.getVal() + r.getAsymErrorHi());
00142 rf.setAsymError(r.getAsymErrorLo(), r.getAsymErrorHi());
00143 }
00144 } else {
00145 r.setVal(r0); r.setConstant(true);
00146
00147 CascadeMinimizer minim2(*nll, CascadeMinimizer::Constrained);
00148 minim2.setStrategy(minimizerStrategyForMinos_);
00149
00150 std::auto_ptr<RooArgSet> allpars(nll->getParameters((const RooArgSet *)0));
00151
00152 double nll0 = nll->getVal();
00153 double threshold68 = nll0 + delta68;
00154 double threshold95 = nll0 + delta95;
00155
00156
00157 assert(!std::isnan(r0));
00158
00159 double hi68 = findCrossing(minim2, *nll, r, threshold68, r0, rMax);
00160 double hi95 = do95_ ? findCrossing(minim2, *nll, r, threshold95, std::isnan(hi68) ? r0 : hi68, std::max(rMax, std::isnan(hi68*2-r0) ? r0 : hi68*2-r0)) : r0;
00161
00162 *allpars = RooArgSet(ret->floatParsFinal()); r.setVal(r0); r.setConstant(true);
00163 double lo68 = findCrossing(minim2, *nll, r, threshold68, r0, rMin);
00164 double lo95 = do95_ ? findCrossing(minim2, *nll, r, threshold95, std::isnan(lo68) ? r0 : lo68, rMin) : r0;
00165
00166 rf.setAsymError(!std::isnan(lo68) ? lo68 - r0 : 0, !std::isnan(hi68) ? hi68 - r0 : 0);
00167 rf.setRange("err68", !std::isnan(lo68) ? lo68 : r0, !std::isnan(hi68) ? hi68 : r0);
00168 if (do95_ && (!std::isnan(lo95) || !std::isnan(hi95))) {
00169 rf.setRange("err95", !std::isnan(lo95) ? lo95 : r0, !std::isnan(hi95) ? hi95 : r0);
00170 }
00171
00172 r.setVal(r0); r.setConstant(false);
00173 }
00174 }
00175
00176 return ret;
00177 }
00178
00179 double FitterAlgoBase::findCrossing(CascadeMinimizer &minim, RooAbsReal &nll, RooRealVar &r, double level, double rStart, double rBound) {
00180 ProfileLikelihood::MinimizerSentry minimizerConfig(minimizerAlgoForMinos_, minimizerToleranceForMinos_);
00181 if (verbose) std::cout << "Searching for crossing at nll = " << level << " in the interval " << rStart << ", " << rBound << std::endl;
00182 double rInc = stepSize_*(rBound - rStart);
00183 r.setVal(rStart);
00184 std::auto_ptr<RooFitResult> checkpoint;
00185 std::auto_ptr<RooArgSet> allpars;
00186 bool ok = false;
00187 {
00188 CloseCoutSentry sentry(verbose < 3);
00189 ok = minim.improve(verbose-1);
00190 checkpoint.reset(minim.save());
00191 }
00192 if (!ok) { std::cout << "Error: minimization failed at " << r.GetName() << " = " << rStart << std::endl; return NAN; }
00193 double here = nll.getVal();
00194 int nfail = 0;
00195 if (verbose > 0) { printf(" %s lvl-here lvl-there stepping\n", r.GetName()); fflush(stdout); }
00196 do {
00197 rStart += rInc;
00198 if (rInc*(rStart - rBound) > 0) {
00199 rStart -= rInc;
00200 rInc = 0.5*(rBound-rStart);
00201 }
00202 r.setVal(rStart);
00203 nll.clearEvalErrorLog(); nll.getVal();
00204 if (nll.numEvalErrors() > 0) {
00205 ok = false;
00206 } else {
00207 CloseCoutSentry sentry(verbose < 3);
00208 ok = minim.improve(verbose-1);
00209 }
00210 if (!ok) {
00211 nfail++;
00212 if (nfail >= maxFailedSteps_) { std::cout << "Error: minimization failed at " << r.GetName() << " = " << rStart << std::endl; return NAN; }
00213 RooArgSet oldparams(checkpoint->floatParsFinal());
00214 if (allpars.get() == 0) allpars.reset(nll.getParameters((const RooArgSet *)0));
00215 *allpars = oldparams;
00216 rStart -= rInc; rInc *= 0.5;
00217 continue;
00218 } else nfail = 0;
00219 double there = here;
00220 here = nll.getVal();
00221 if (verbose > 0) { printf("%f %+.5f %+.5f %f\n", rStart, level-here, level-there, rInc); fflush(stdout); }
00222 if ( fabs(here - level) < 4*minimizerToleranceForMinos_ ) {
00223
00224 r.setVal(rStart + (level-here)*(level-there)/(here-there));
00225 return r.getVal();
00226 } else if (here > level) {
00227
00228
00229 rStart -= rInc;
00230
00231 if (runtimedef::get("FITTER_DYN_STEP")) {
00232 if (fabs(there - level) > 0.05) {
00233 double rIncFactor = std::max(0.2, std::min(0.7, 0.75*(level-there)/(here-there)));
00234
00235 rInc *= rIncFactor;
00236 } else {
00237 double rIncFactor = std::max(0.05, std::min(0.95, 0.95*(level-there)/(here-there)));
00238
00239 rInc *= rIncFactor;
00240 }
00241 } else {
00242 rInc *= 0.3;
00243 }
00244 if (allpars.get() == 0) allpars.reset(nll.getParameters((const RooArgSet *)0));
00245 RooArgSet oldparams(checkpoint->floatParsFinal());
00246 *allpars = oldparams;
00247 } else if ((here-there)*(level-there) < 0 &&
00248 fabs(here-there) > 0.1) {
00249 if (allpars.get() == 0) allpars.reset(nll.getParameters((const RooArgSet *)0));
00250 RooArgSet oldparams(checkpoint->floatParsFinal());
00251 *allpars = oldparams;
00252 rStart -= rInc; rInc *= 0.5;
00253 } else {
00254
00255 if (runtimedef::get("FITTER_DYN_STEP")) {
00256 if (fabs(here - level) > 0.05) {
00257 if ((here-there)*(level-there) > 0) {
00258
00259 double rIncFactor = std::max(0.2, std::min(2.0, 0.75*(level-there)/(here-there)));
00260
00261 rInc *= rIncFactor;
00262 }
00263 } else {
00264 double rIncFactor = std::max(0.05, std::min(4.0, 0.95*(level-there)/(here-there)));
00265
00266 rInc *= rIncFactor;
00267 }
00268 } else {
00269
00270 }
00271 checkpoint.reset(minim.save());
00272 }
00273 } while (fabs(rInc) > minimizerToleranceForMinos_*stepSize_*std::max(1.0,rBound-rStart));
00274 if (fabs(here - level) > 0.01) {
00275 std::cout << "Error: closed range without finding crossing." << std::endl;
00276 return NAN;
00277 } else {
00278 return r.getVal();
00279 }
00280 }