1 #include "../interface/FitterAlgoBase.h"
5 #include "RooRealVar.h"
8 #include "RooDataSet.h"
9 #include "RooFitResult.h"
10 #include "RooSimultaneous.h"
11 #include "RooAddPdf.h"
12 #include "RooProdPdf.h"
13 #include "RooConstVar.h"
15 #include "../interface/RooMinimizerOpt.h"
20 #include <RooStats/ModelConfig.h>
21 #include "../interface/Combine.h"
22 #include "../interface/ProfileLikelihood.h"
23 #include "../interface/CascadeMinimizer.h"
24 #include "../interface/CloseCoutSentry.h"
25 #include "../interface/utils.h"
27 #include "../interface/ProfilingTools.h"
30 #include <Math/MinimizerOptions.h>
31 #include <Math/QuantFuncMathCore.h>
32 #include <Math/ProbFunc.h>
34 using namespace RooStats;
38 #include <Math/QuantFuncMathCore.h>
39 #include <Math/ProbFunc.h>
57 (
"minimizerAlgo", boost::program_options::value<std::string>(&
minimizerAlgo_)->default_value(
minimizerAlgo_),
"Choice of minimizer (Minuit vs Minuit2)")
60 (
"preFitValue", boost::program_options::value<float>(&
preFitValue_)->default_value(
preFitValue_),
"Value of signal strength pre-fit")
61 (
"do95", boost::program_options::value<bool>(&
do95_)->default_value(
do95_),
"Compute also 2-sigma interval from delta(nll) = 1.92 instead of 0.5")
62 (
"robustFit", boost::program_options::value<bool>(&
robustFit_)->default_value(
robustFit_),
"Search manually for 1 and 2 sigma bands instead of using Minos")
63 (
"maxFailedSteps", boost::program_options::value<int>(&
maxFailedSteps_)->default_value(
maxFailedSteps_),
"How many failed steps to retry before giving up")
64 (
"stepSize", boost::program_options::value<float>(&
stepSize_)->default_value(
stepSize_),
"Step size for robust fits (multiplier of the range)")
68 (
"saveNLL",
"Save the negative log-likelihood at the minimum in the output tree (note: value is relative to the pre-fit state)")
69 (
"keepFailures",
"Save the results even if the fit is declared as failed (for NLL studies)")
79 bool FitterAlgoBase::run(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
83 static bool shouldCreateNLLBranch =
saveNLL_;
86 return runSpecific(w, mc_s, mc_b, data, limit, limitErr, hint);
90 RooFitResult *
FitterAlgoBase::doFit(RooAbsPdf &pdf, RooAbsData &
data, RooRealVar &
r,
const RooCmdArg &constrain,
bool doHesse,
int ndim,
bool reuseNLL) {
91 return doFit(pdf, data, RooArgList (r), constrain, doHesse, ndim, reuseNLL);
94 RooFitResult *
FitterAlgoBase::doFit(RooAbsPdf &pdf, RooAbsData &
data,
const RooArgList &rs,
const RooCmdArg &constrain,
bool doHesse,
int ndim,
bool reuseNLL) {
95 RooFitResult *
ret = 0;
96 if (reuseNLL)
nll->setData(data);
97 else nll.reset(pdf.createNLL(data, constrain, RooFit::Extended(pdf.canBeExtended())));
99 double nll0 =
nll->getVal();
100 double delta68 = 0.5*ROOT::Math::chisquared_quantile_c(1-0.68,ndim);
101 double delta95 = 0.5*ROOT::Math::chisquared_quantile_c(1-0.95,ndim);
108 if (!ok && !
keepFailures_) {
std::cout <<
"Initial minimization failed. Aborting." << std::endl;
return 0; }
112 if (
verbose > 1) { ret->Print(
"V"); }
114 std::auto_ptr<RooArgSet> allpars(pdf.getParameters(data));
116 for (
int i = 0,
n = rs.getSize();
i <
n; ++
i) {
119 RooArgSet oldparams(ret->floatParsFinal());
120 *allpars = oldparams;
124 RooRealVar &
r =
dynamic_cast<RooRealVar &
>(*rs.at(
i));
125 RooRealVar &rf =
dynamic_cast<RooRealVar &
>(*ret->floatParsFinal().find(r.GetName()));
126 double r0 = r.getVal(), rMin = r.getMin(), rMax = r.getMax();
130 throw std::runtime_error(
"95% CL errors with Minos are not working at the moment.");
134 if (minim.
minimizer().minos(RooArgSet(r)) != -1) {
135 rf.setRange(
"err95", r.getVal() + r.getAsymErrorLo(), r.getVal() + r.getAsymErrorHi());
140 if (minim.
minimizer().minos(RooArgSet(r)) != -1) {
141 rf.setRange(
"err68", r.getVal() + r.getAsymErrorLo(), r.getVal() + r.getAsymErrorHi());
142 rf.setAsymError(r.getAsymErrorLo(), r.getAsymErrorHi());
145 r.setVal(r0); r.setConstant(
true);
150 std::auto_ptr<RooArgSet> allpars(
nll->getParameters((
const RooArgSet *)0));
152 double nll0 =
nll->getVal();
153 double threshold68 = nll0 + delta68;
154 double threshold95 = nll0 + delta95;
162 *allpars = RooArgSet(ret->floatParsFinal()); r.setVal(r0); r.setConstant(
true);
172 r.setVal(r0); r.setConstant(
false);
181 if (
verbose)
std::cout <<
"Searching for crossing at nll = " << level <<
" in the interval " << rStart <<
", " << rBound << std::endl;
182 double rInc =
stepSize_*(rBound - rStart);
184 std::auto_ptr<RooFitResult> checkpoint;
185 std::auto_ptr<RooArgSet> allpars;
190 checkpoint.reset(minim.
save());
192 if (!ok) {
std::cout <<
"Error: minimization failed at " << r.GetName() <<
" = " << rStart << std::endl;
return NAN; }
193 double here = nll.getVal();
195 if (
verbose > 0) { printf(
" %s lvl-here lvl-there stepping\n", r.GetName()); fflush(stdout); }
198 if (rInc*(rStart - rBound) > 0) {
200 rInc = 0.5*(rBound-rStart);
203 nll.clearEvalErrorLog(); nll.getVal();
204 if (nll.numEvalErrors() > 0) {
212 if (nfail >=
maxFailedSteps_) {
std::cout <<
"Error: minimization failed at " << r.GetName() <<
" = " << rStart << std::endl;
return NAN; }
213 RooArgSet oldparams(checkpoint->floatParsFinal());
214 if (allpars.get() == 0) allpars.reset(nll.getParameters((
const RooArgSet *)0));
215 *allpars = oldparams;
216 rStart -= rInc; rInc *= 0.5;
221 if (
verbose > 0) { printf(
"%f %+.5f %+.5f %f\n", rStart, level-here, level-there, rInc); fflush(stdout); }
224 r.setVal(rStart + (level-here)*(level-there)/(here-there));
226 }
else if (here > level) {
232 if (fabs(there - level) > 0.05) {
233 double rIncFactor =
std::max(0.2,
std::min(0.7, 0.75*(level-there)/(here-there)));
237 double rIncFactor =
std::max(0.05,
std::min(0.95, 0.95*(level-there)/(here-there)));
244 if (allpars.get() == 0) allpars.reset(nll.getParameters((
const RooArgSet *)0));
245 RooArgSet oldparams(checkpoint->floatParsFinal());
246 *allpars = oldparams;
247 }
else if ((here-there)*(level-there) < 0 &&
248 fabs(here-there) > 0.1) {
249 if (allpars.get() == 0) allpars.reset(nll.getParameters((
const RooArgSet *)0));
250 RooArgSet oldparams(checkpoint->floatParsFinal());
251 *allpars = oldparams;
252 rStart -= rInc; rInc *= 0.5;
256 if (fabs(here - level) > 0.05) {
257 if ((here-there)*(level-there) > 0) {
259 double rIncFactor =
std::max(0.2,
std::min(2.0, 0.75*(level-there)/(here-there)));
264 double rIncFactor =
std::max(0.05,
std::min(4.0, 0.95*(level-there)/(here-there)));
271 checkpoint.reset(minim.
save());
274 if (fabs(here - level) > 0.01) {
275 std::cout <<
"Error: closed range without finding crossing." << std::endl;
static void addBranch(const char *name, void *address, const char *leaflist)
Add a branch to the output tree (for advanced use or debugging only)
void applyOptionsBase(const boost::program_options::variables_map &vm)
int get(const char *name)
bool minimize(int verbose=0, bool cascade=true)
static std::string minimizerAlgoForMinos_
bool improve(int verbose=0, bool cascade=true)
static float minimizerTolerance_
RooFitResult * doFit(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &r, const RooCmdArg &constrain, bool doHesse=true, int ndim=1, bool reuseNLL=false)
static bool keepFailures_
double findCrossing(CascadeMinimizer &minim, RooAbsReal &nll, RooRealVar &r, double level, double rStart, double rBound)
const T & max(const T &a, const T &b)
static float preFitValue_
std::auto_ptr< RooAbsReal > nll
static std::string minimizerAlgo_
static int minimizerStrategy_
FitterAlgoBase(const char *title="<FillMe> specific options")
virtual bool runSpecific(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)=0
virtual bool run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
Setup Minimizer configuration on creation, reset the previous one on destruction. ...
static float minimizerToleranceForMinos_
RooMinimizerOpt & minimizer()
char data[epos_bytes_allocation]
static int maxFailedSteps_
static int minimizerStrategyForMinos_
boost::program_options::options_description options_
void setErrorLevel(float errorLevel)
void setStrategy(int strategy)