1 #include "../interface/MultiDimFit.h"
7 #include "RooArgList.h"
9 #include "RooAbsData.h"
10 #include "RooFitResult.h"
11 #include "../interface/RooMinimizerOpt.h"
12 #include <RooStats/ModelConfig.h>
13 #include "../interface/Combine.h"
14 #include "../interface/CascadeMinimizer.h"
15 #include "../interface/CloseCoutSentry.h"
16 #include "../interface/utils.h"
18 #include <Math/Minimizer.h>
19 #include <Math/MinimizerOptions.h>
20 #include <Math/QuantFuncMathCore.h>
21 #include <Math/ProbFunc.h>
23 using namespace RooStats;
43 (
"algo", boost::program_options::value<std::string>()->default_value(
"none"),
"Algorithm to compute uncertainties")
45 (
"floatOtherPOIs", boost::program_options::value<bool>(&
floatOtherPOIs_)->default_value(
floatOtherPOIs_),
"POIs other than the selected ones will be kept freely floating (1) or fixed (0, default)")
46 (
"points", boost::program_options::value<unsigned int>(&
points_)->default_value(
points_),
"Points to use for grid or contour scans")
47 (
"firstPoint", boost::program_options::value<unsigned int>(&
firstPoint_)->default_value(
firstPoint_),
"First point to use")
48 (
"lastPoint", boost::program_options::value<unsigned int>(&
lastPoint_)->default_value(
lastPoint_),
"Last point to use")
49 (
"fastScan",
"Do a fast scan, evaluating the likelihood without profiling it.")
56 std::string
algo = vm[
"algo"].as<std::string>();
59 }
else if (algo ==
"singles") {
61 }
else if (algo ==
"cross") {
63 }
else if (algo ==
"grid") {
65 }
else if (algo ==
"random") {
67 }
else if (algo ==
"contour2d") {
69 }
else throw std::invalid_argument(std::string(
"Unknown algorithm: "+algo));
73 bool MultiDimFit::runSpecific(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
75 static int isInit =
false;
76 if (!isInit) {
initOnce(w, mc_s); isInit =
true; }
79 RooAbsPdf &pdf = *mc_s->GetPdf();
83 RooLinkedListIter iterP = mc_s->GetParametersOfInterest()->iterator();
84 for (RooAbsArg *
a = (RooAbsArg*) iterP.Next();
a != 0;
a = (RooAbsArg*) iterP.Next()) {
86 RooRealVar *rrv =
dynamic_cast<RooRealVar *
>(
a);
87 if (rrv == 0) {
std::cerr <<
"MultiDimFit: Parameter of interest " <<
a->GetName() <<
" which is not a RooRealVar will be ignored" << std::endl;
continue; }
93 const RooCmdArg &constrainCmdArg =
withSystematics ? RooFit::Constrain(*mc_s->GetNuisanceParameters()) : RooFit::NumCPU(1);
96 for (
int i = 0,
n =
poi_.size();
i <
n; ++
i) {
101 std::auto_ptr<RooAbsReal>
nll;
103 nll.reset(pdf.createNLL(data, constrainCmdArg, RooFit::Extended(pdf.canBeExtended())));
108 std::cout <<
"\n --- MultiDimFit ---" << std::endl;
109 std::cout <<
"best fit parameter values: " << std::endl;
110 int len =
poi_[0].length();
111 for (
int i = 0,
n =
poi_.size();
i <
n; ++
i) {
112 len = std::max<int>(len,
poi_[
i].length());
114 for (
int i = 0, n =
poi_.size();
i <
n; ++
i) {
130 RooArgSet mcPoi(*mc_s->GetParametersOfInterest());
132 RooLinkedListIter iterP = mc_s->GetParametersOfInterest()->iterator();
133 for (RooAbsArg *
a = (RooAbsArg*) iterP.Next();
a != 0;
a = (RooAbsArg*) iterP.Next()) {
134 poi_.push_back(
a->GetName());
137 for (std::vector<std::string>::const_iterator it =
poi_.begin(), ed =
poi_.end(); it != ed; ++it) {
138 RooAbsArg *
a = mcPoi.find(it->c_str());
139 if (a == 0)
throw std::invalid_argument(std::string(
"Parameter of interest ")+*it+
" not in model.");
140 RooRealVar *rrv =
dynamic_cast<RooRealVar *
>(
a);
141 if (rrv == 0)
throw std::invalid_argument(std::string(
"Parameter of interest ")+*it+
" not a RooRealVar.");
147 for (
int i = 0,
n =
poi_.size();
i <
n; ++
i) {
155 std::cout <<
"\n --- MultiDimFit ---" << std::endl;
156 std::cout <<
"best fit parameter values and profile-likelihood uncertainties: " << std::endl;
157 int len =
poi_[0].length();
158 for (
int i = 0,
n =
poi_.size();
i <
n; ++
i) {
159 len = std::max<int>(len,
poi_[
i].length());
161 for (
int i = 0, n =
poi_.size();
i <
n; ++
i) {
162 RooRealVar *rf =
dynamic_cast<RooRealVar*
>(res.floatParsFinal().find(
poi_[
i].c_str()));
163 double bestFitVal = rf->getVal();
165 double hiErr = +(rf->hasRange(
"err68") ? rf->getMax(
"err68") - bestFitVal : rf->getAsymErrorHi());
166 double loErr = -(rf->hasRange(
"err68") ? rf->getMin(
"err68") - bestFitVal : rf->getAsymErrorLo());
167 double maxError = std::max<double>(std::max<double>(hiErr, loErr), rf->getError());
169 if (fabs(hiErr) < 0.001*maxError) hiErr = -bestFitVal + rf->getMax();
170 if (fabs(loErr) < 0.001*maxError) loErr = +bestFitVal - rf->getMin();
172 double hiErr95 = +(
do95_ && rf->hasRange(
"err95") ? rf->getMax(
"err95") - bestFitVal : 0);
173 double loErr95 = -(
do95_ && rf->hasRange(
"err95") ? rf->getMin(
"err95") - bestFitVal : 0);
177 if (
do95_ && rf->hasRange(
"err95")) {
181 printf(
" %*s : %+8.3f %+6.3f/%+6.3f (68%%) %+6.3f/%+6.3f (95%%) \n", len,
poi_[
i].c_str(),
182 poiVals_[
i], -loErr, hiErr, loErr95, -hiErr95);
185 printf(
" %*s : %+8.3f %+6.3f/%+6.3f (68%%)\n", len,
poi_[
i].c_str(),
193 unsigned int n =
poi_.size();
194 if (
poi_.size() > 2)
throw std::logic_error(
"Don't know how to do a grid with more than 2 POIs.");
195 double nll0 = nll.getVal();
197 std::vector<double> p0(n), pmin(n), pmax(n);
198 for (
unsigned int i = 0;
i <
n; ++
i) {
207 std::auto_ptr<RooArgSet> params(nll.getParameters((
const RooArgSet *)0));
208 RooArgSet snap; params->snapshot(snap);
212 double x = pmin[0] + (
i+0.5)*(pmax[0]-pmin[0])/points_;
213 if (
verbose > 1)
std::cout <<
"Point " <<
i <<
"/" << points_ <<
" " <<
poiVars_[0]->GetName() <<
" = " << x << std::endl;
231 unsigned int ipoint = 0, nprint = ceil(0.005*sqrn*sqrn);
232 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CountErrors);
234 for (
unsigned int i = 0;
i < sqrn; ++
i) {
235 for (
unsigned int j = 0;
j < sqrn; ++
j, ++
ipoint) {
239 double x = pmin[0] + (
i+0.5)*(pmax[0]-pmin[0])/sqrn;
240 double y = pmin[1] + (
j+0.5)*(pmax[1]-pmin[1])/sqrn;
241 if (
verbose && (ipoint % nprint == 0)) {
242 fprintf(sentry.
trueStdOut(),
"Point %d/%d, (i,j) = (%d,%d), %s = %f, %s = %f\n",
249 nll.clearEvalErrorLog(); nll.getVal();
250 if (nll.numEvalErrors() > 0) {
271 double nll0 = nll.getVal();
272 for (
unsigned int i = 0,
n =
poi_.size();
i <
n; ++
i) {
278 unsigned int n =
poi_.size();
280 for (
unsigned int i = 0;
i <
n; ++
i) {
289 double qN = 2*(nll.getVal() - nll0);
299 if (
poi_.size() != 2)
throw std::logic_error(
"Contour2D works only in 2 dimensions");
304 if (
verbose>0)
std::cout <<
"Best fit point is for " << xv->GetName() <<
", " << yv->GetName() <<
" = " << x0 <<
", " << y0 << std::endl;
308 double xMin = xv->getMin(
"box"), xMax = xv->getMax(
"box");
309 double yMin = yv->getMin(
"box"), yMax = yv->getMax(
"box");
313 yv->setConstant(
true);
316 double yc = 0.5*(yMax + yMin), yr = 0.5*(yMax - yMin);
319 xv->setConstant(
false); xv->setVal(x0);
326 double xc = xv->getVal(); xv->setConstant(
true);
327 if (
verbose>-1)
std::cout <<
"Best fit " << xv->GetName() <<
" for " << yv->GetName() <<
" = " << yv->getVal() <<
" is at " << xc << std::endl;
330 double xup =
findCrossing(minim, nll, *xv, threshold, xc, xMax);
333 if (
verbose>-1)
std::cout <<
"Minimum of " << xv->GetName() <<
" at " <<
cl <<
" CL for " << yv->GetName() <<
" = " << y <<
" is " << x << std::endl;
336 double xdn =
findCrossing(minim, nll, *xv, threshold, xc, xMin);
339 if (
verbose>-1)
std::cout <<
"Maximum of " << xv->GetName() <<
" at " <<
cl <<
" CL for " << yv->GetName() <<
" = " << y <<
" is " << x << std::endl;
346 unsigned int n =
poi_.size();
349 std::vector<double> p0(n);
350 for (
unsigned int i = 0;
i <
n; ++
i) {
356 for (
unsigned int i = 0;
i <
n; ++
i) {
358 xv->setConstant(
true);
365 if (
verbose > -1)
std::cout <<
"Minimum of " << xv->GetName() <<
" at " << cl <<
" CL for all others floating is " << xMin << std::endl;
371 double prob = ROOT::Math::chisquared_cdf_c(2*(nll.getVal() - nll0), n+
nOtherFloatingPoi_);
373 if (
verbose > -1)
std::cout <<
"Minimum of " << xv->GetName() <<
" at " << cl <<
" CL for all others floating is " << xMin <<
" (on the boundary, p-val " << prob <<
")" << std::endl;
376 for (
unsigned int j = 0; j <
n; ++
j)
poiVars_[j]->setVal(p0[j]);
379 if (
verbose > -1)
std::cout <<
"Maximum of " << xv->GetName() <<
" at " << cl <<
" CL for all others floating is " << xMax << std::endl;
384 double prob = ROOT::Math::chisquared_cdf_c(2*(nll.getVal() - nll0), n+
nOtherFloatingPoi_);
387 if (
verbose > -1)
std::cout <<
"Maximum of " << xv->GetName() <<
" at " << cl <<
" CL for all others floating is " << xMax <<
" (on the boundary, p-val " << prob <<
")" << std::endl;
390 xv->setRange(name, xMin, xMax);
391 xv->setConstant(
false);
static std::vector< std::string > poi_
static unsigned int points_
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)
static bool floatOtherPOIs_
void doGrid(RooAbsReal &nll)
bool minimize(int verbose=0, bool cascade=true)
static unsigned int lastPoint_
static std::vector< float > poiVals_
void initOnce(RooWorkspace *w, RooStats::ModelConfig *mc_s)
void doBox(RooAbsReal &nll, double cl, const char *name="box", bool commitPoints=true)
for each RooRealVar, set a range 'box' from the PL profiling all other parameters ...
RooFitResult * doFit(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &r, const RooCmdArg &constrain, bool doHesse=true, int ndim=1, bool reuseNLL=false)
static bool keepFailures_
static void commitPoint(bool expected, float quantile)
Save a point into the output tree. Usually if expected = false, quantile should be set to -1 (except ...
double findCrossing(CascadeMinimizer &minim, RooAbsReal &nll, RooRealVar &r, double level, double rStart, double rBound)
const T & max(const T &a, const T &b)
virtual void applyOptions(const boost::program_options::variables_map &vm)
Cos< T >::type cos(const T &t)
std::auto_ptr< RooAbsReal > nll
static unsigned int nOtherFloatingPoi_
virtual bool runSpecific(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
void doSingles(RooFitResult &res)
static int minimizerStrategy_
static std::vector< RooRealVar * > poiVars_
void doRandomPoints(RooAbsReal &nll)
static RooArgList poiList_
char data[epos_bytes_allocation]
void doContour2D(RooAbsReal &nll)
static unsigned int firstPoint_
boost::program_options::options_description options_
void setStrategy(int strategy)