CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiDimFit.cc
Go to the documentation of this file.
1 #include "../interface/MultiDimFit.h"
2 #include <stdexcept>
3 #include <cmath>
4 
5 #include "TMath.h"
6 #include "RooArgSet.h"
7 #include "RooArgList.h"
8 #include "RooRandom.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"
17 
18 #include <Math/Minimizer.h>
19 #include <Math/MinimizerOptions.h>
20 #include <Math/QuantFuncMathCore.h>
21 #include <Math/ProbFunc.h>
22 
23 using namespace RooStats;
24 
26 std::vector<std::string> MultiDimFit::poi_;
27 std::vector<RooRealVar *> MultiDimFit::poiVars_;
28 std::vector<float> MultiDimFit::poiVals_;
29 RooArgList MultiDimFit::poiList_;
30 float MultiDimFit::deltaNLL_ = 0;
31 unsigned int MultiDimFit::points_ = 50;
32 unsigned int MultiDimFit::firstPoint_ = 0;
34 bool MultiDimFit::floatOtherPOIs_ = false;
35 unsigned int MultiDimFit::nOtherFloatingPoi_ = 0;
36 bool MultiDimFit::fastScan_ = false;
37 
38 
40  FitterAlgoBase("MultiDimFit specific options")
41 {
42  options_.add_options()
43  ("algo", boost::program_options::value<std::string>()->default_value("none"), "Algorithm to compute uncertainties")
44  ("poi,P", boost::program_options::value<std::vector<std::string> >(&poi_), "Parameters of interest to fit (default = all)")
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.")
50  ;
51 }
52 
53 void MultiDimFit::applyOptions(const boost::program_options::variables_map &vm)
54 {
55  applyOptionsBase(vm);
56  std::string algo = vm["algo"].as<std::string>();
57  if (algo == "none") {
58  algo_ = None;
59  } else if (algo == "singles") {
60  algo_ = Singles;
61  } else if (algo == "cross") {
62  algo_ = Cross;
63  } else if (algo == "grid") {
64  algo_ = Grid;
65  } else if (algo == "random") {
67  } else if (algo == "contour2d") {
68  algo_ = Contour2D;
69  } else throw std::invalid_argument(std::string("Unknown algorithm: "+algo));
70  fastScan_ = (vm.count("fastScan") > 0);
71 }
72 
73 bool MultiDimFit::runSpecific(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) {
74  // one-time initialization of POI variables, TTree branches, ...
75  static int isInit = false;
76  if (!isInit) { initOnce(w, mc_s); isInit = true; }
77 
78  // Get PDF
79  RooAbsPdf &pdf = *mc_s->GetPdf();
80 
81  // Process POI not in list
83  RooLinkedListIter iterP = mc_s->GetParametersOfInterest()->iterator();
84  for (RooAbsArg *a = (RooAbsArg*) iterP.Next(); a != 0; a = (RooAbsArg*) iterP.Next()) {
85  if (poiList_.contains(*a)) continue;
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; }
88  rrv->setConstant(!floatOtherPOIs_);
90  }
91 
92  // start with a best fit
93  const RooCmdArg &constrainCmdArg = withSystematics ? RooFit::Constrain(*mc_s->GetNuisanceParameters()) : RooFit::NumCPU(1);
94  std::auto_ptr<RooFitResult> res(doFit(pdf, data, (algo_ == Singles ? poiList_ : RooArgList()), constrainCmdArg));
95  if (res.get() || keepFailures_) {
96  for (int i = 0, n = poi_.size(); i < n; ++i) {
97  poiVals_[i] = poiVars_[i]->getVal();
98  }
99  Combine::commitPoint(/*expected=*/false, /*quantile=*/1.);
100  }
101  std::auto_ptr<RooAbsReal> nll;
102  if (algo_ != None && algo_ != Singles) {
103  nll.reset(pdf.createNLL(data, constrainCmdArg, RooFit::Extended(pdf.canBeExtended())));
104  }
105  switch(algo_) {
106  case None:
107  if (verbose > 0) {
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());
113  }
114  for (int i = 0, n = poi_.size(); i < n; ++i) {
115  printf(" %*s : %+8.3f\n", len, poi_[i].c_str(), poiVals_[i]);
116  }
117  }
118  break;
119  case Singles: if (res.get()) doSingles(*res); break;
120  case Cross: doBox(*nll, cl, "box", true); break;
121  case Grid: doGrid(*nll); break;
122  case RandomPoints: doRandomPoints(*nll); break;
123  case Contour2D: doContour2D(*nll); break;
124  }
125 
126  return true;
127 }
128 
129 void MultiDimFit::initOnce(RooWorkspace *w, RooStats::ModelConfig *mc_s) {
130  RooArgSet mcPoi(*mc_s->GetParametersOfInterest());
131  if (poi_.empty()) {
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());
135  }
136  }
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.");
142  poiVars_.push_back(rrv);
143  poiVals_.push_back(rrv->getVal());
144  poiList_.add(*rrv);
145  }
146  // then add the branches to the tree (at the end, so there are no resizes)
147  for (int i = 0, n = poi_.size(); i < n; ++i) {
148  Combine::addBranch(poi_[i].c_str(), &poiVals_[i], (poi_[i]+"/F").c_str());
149  }
150  Combine::addBranch("deltaNLL", &deltaNLL_, "deltaNLL/F");
151 }
152 
153 void MultiDimFit::doSingles(RooFitResult &res)
154 {
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());
160  }
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();
164 
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());
168 
169  if (fabs(hiErr) < 0.001*maxError) hiErr = -bestFitVal + rf->getMax();
170  if (fabs(loErr) < 0.001*maxError) loErr = +bestFitVal - rf->getMin();
171 
172  double hiErr95 = +(do95_ && rf->hasRange("err95") ? rf->getMax("err95") - bestFitVal : 0);
173  double loErr95 = -(do95_ && rf->hasRange("err95") ? rf->getMin("err95") - bestFitVal : 0);
174 
175  poiVals_[i] = bestFitVal - loErr; Combine::commitPoint(true, /*quantile=*/0.32);
176  poiVals_[i] = bestFitVal + hiErr; Combine::commitPoint(true, /*quantile=*/0.32);
177  if (do95_ && rf->hasRange("err95")) {
178  poiVals_[i] = rf->getMax("err95"); Combine::commitPoint(true, /*quantile=*/0.05);
179  poiVals_[i] = rf->getMin("err95"); Combine::commitPoint(true, /*quantile=*/0.05);
180  poiVals_[i] = bestFitVal;
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);
183  } else {
184  poiVals_[i] = bestFitVal;
185  printf(" %*s : %+8.3f %+6.3f/%+6.3f (68%%)\n", len, poi_[i].c_str(),
186  poiVals_[i], -loErr, hiErr);
187  }
188  }
189 }
190 
191 void MultiDimFit::doGrid(RooAbsReal &nll)
192 {
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();
196 
197  std::vector<double> p0(n), pmin(n), pmax(n);
198  for (unsigned int i = 0; i < n; ++i) {
199  p0[i] = poiVars_[i]->getVal();
200  pmin[i] = poiVars_[i]->getMin();
201  pmax[i] = poiVars_[i]->getMax();
202  poiVars_[i]->setConstant(true);
203  }
204 
207  std::auto_ptr<RooArgSet> params(nll.getParameters((const RooArgSet *)0));
208  RooArgSet snap; params->snapshot(snap);
209  //snap.Print("V");
210  if (n == 1) {
211  for (unsigned int i = 0; i < points_; ++i) {
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;
214  *params = snap;
215  poiVals_[0] = x;
216  poiVars_[0]->setVal(x);
217  // now we minimize
218  {
219  //CloseCoutSentry sentry(verbose < 3);
220  bool ok = fastScan_ ? true : minim.minimize(verbose-1);
221  if (ok) {
222  deltaNLL_ = nll.getVal() - nll0;
223  double qN = 2*(deltaNLL_);
224  double prob = ROOT::Math::chisquared_cdf_c(qN, n+nOtherFloatingPoi_);
225  Combine::commitPoint(true, /*quantile=*/prob);
226  }
227  }
228  }
229  } else if (n == 2) {
230  unsigned int sqrn = ceil(sqrt(double(points_)));
231  unsigned int ipoint = 0, nprint = ceil(0.005*sqrn*sqrn);
232  RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CountErrors);
233  CloseCoutSentry sentry(verbose < 2);
234  for (unsigned int i = 0; i < sqrn; ++i) {
235  for (unsigned int j = 0; j < sqrn; ++j, ++ipoint) {
236  if (ipoint < firstPoint_) continue;
237  if (ipoint > lastPoint_) break;
238  *params = snap;
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",
243  ipoint,sqrn*sqrn, i,j, poiVars_[0]->GetName(), x, poiVars_[1]->GetName(), y);
244  }
245  poiVals_[0] = x;
246  poiVals_[1] = y;
247  poiVars_[0]->setVal(x);
248  poiVars_[1]->setVal(y);
249  nll.clearEvalErrorLog(); nll.getVal();
250  if (nll.numEvalErrors() > 0) {
251  deltaNLL_ = 9999; Combine::commitPoint(true, /*quantile=*/0);
252  continue;
253  }
254  // now we minimize
255  {
256  bool ok = fastScan_ ? true : minim.minimize(verbose-1);
257  if (ok) {
258  deltaNLL_ = nll.getVal() - nll0;
259  double qN = 2*(deltaNLL_);
260  double prob = ROOT::Math::chisquared_cdf_c(qN, n+nOtherFloatingPoi_);
261  Combine::commitPoint(true, /*quantile=*/prob);
262  }
263  }
264  }
265  }
266  }
267 }
268 
269 void MultiDimFit::doRandomPoints(RooAbsReal &nll)
270 {
271  double nll0 = nll.getVal();
272  for (unsigned int i = 0, n = poi_.size(); i < n; ++i) {
273  poiVars_[i]->setConstant(true);
274  }
275 
278  unsigned int n = poi_.size();
279  for (unsigned int j = 0; j < points_; ++j) {
280  for (unsigned int i = 0; i < n; ++i) {
281  poiVars_[i]->randomize();
282  poiVals_[i] = poiVars_[i]->getVal();
283  }
284  // now we minimize
285  {
286  CloseCoutSentry sentry(verbose < 3);
287  bool ok = minim.minimize(verbose-1);
288  if (ok) {
289  double qN = 2*(nll.getVal() - nll0);
290  double prob = ROOT::Math::chisquared_cdf_c(qN, n+nOtherFloatingPoi_);
291  Combine::commitPoint(true, /*quantile=*/prob);
292  }
293  }
294  }
295 }
296 
297 void MultiDimFit::doContour2D(RooAbsReal &nll)
298 {
299  if (poi_.size() != 2) throw std::logic_error("Contour2D works only in 2 dimensions");
300  RooRealVar *xv = poiVars_[0]; double x0 = poiVals_[0]; float &x = poiVals_[0];
301  RooRealVar *yv = poiVars_[1]; double y0 = poiVals_[1]; float &y = poiVals_[1];
302 
303  double threshold = nll.getVal() + 0.5*ROOT::Math::chisquared_quantile_c(1-cl,2+nOtherFloatingPoi_);
304  if (verbose>0) std::cout << "Best fit point is for " << xv->GetName() << ", " << yv->GetName() << " = " << x0 << ", " << y0 << std::endl;
305 
306  // make a box
307  doBox(nll, cl, "box");
308  double xMin = xv->getMin("box"), xMax = xv->getMax("box");
309  double yMin = yv->getMin("box"), yMax = yv->getMax("box");
310 
311  verbose--; // reduce verbosity to avoid messages from findCrossing
312  // ===== Get relative min/max of x for several fixed y values =====
313  yv->setConstant(true);
314  for (unsigned int j = 0; j <= points_; ++j) {
315  // take points uniformly spaced in polar angle in the case of a perfect circle
316  double yc = 0.5*(yMax + yMin), yr = 0.5*(yMax - yMin);
317  yv->setVal( yc + yr * std::cos(j*M_PI/double(points_)) );
318  // ===== Get the best fit x (could also do without profiling??) =====
319  xv->setConstant(false); xv->setVal(x0);
322  {
323  CloseCoutSentry sentry(verbose < 3);
324  minimXI.minimize(verbose-1);
325  }
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;
328  // ===== Then get the range =====
330  double xup = findCrossing(minim, nll, *xv, threshold, xc, xMax);
331  if (!std::isnan(xup)) {
332  x = xup; y = yv->getVal(); Combine::commitPoint(true, /*quantile=*/1-cl);
333  if (verbose>-1) std::cout << "Minimum of " << xv->GetName() << " at " << cl << " CL for " << yv->GetName() << " = " << y << " is " << x << std::endl;
334  }
335 
336  double xdn = findCrossing(minim, nll, *xv, threshold, xc, xMin);
337  if (!std::isnan(xdn)) {
338  x = xdn; y = yv->getVal(); Combine::commitPoint(true, /*quantile=*/1-cl);
339  if (verbose>-1) std::cout << "Maximum of " << xv->GetName() << " at " << cl << " CL for " << yv->GetName() << " = " << y << " is " << x << std::endl;
340  }
341  }
342 
343  verbose++; // restore verbosity
344 }
345 void MultiDimFit::doBox(RooAbsReal &nll, double cl, const char *name, bool commitPoints) {
346  unsigned int n = poi_.size();
347  double nll0 = nll.getVal(), threshold = nll0 + 0.5*ROOT::Math::chisquared_quantile_c(1-cl,n+nOtherFloatingPoi_);
348 
349  std::vector<double> p0(n);
350  for (unsigned int i = 0; i < n; ++i) {
351  p0[i] = poiVars_[i]->getVal();
352  poiVars_[i]->setConstant(false);
353  }
354 
355  verbose--; // reduce verbosity due to findCrossing
356  for (unsigned int i = 0; i < n; ++i) {
357  RooRealVar *xv = poiVars_[i];
358  xv->setConstant(true);
361 
362  for (unsigned int j = 0; j < n; ++j) poiVars_[j]->setVal(p0[j]);
363  double xMin = findCrossing(minimX, nll, *xv, threshold, p0[i], xv->getMin());
364  if (!std::isnan(xMin)) {
365  if (verbose > -1) std::cout << "Minimum of " << xv->GetName() << " at " << cl << " CL for all others floating is " << xMin << std::endl;
366  for (unsigned int j = 0; j < n; ++j) poiVals_[j] = poiVars_[j]->getVal();
367  if (commitPoints) Combine::commitPoint(true, /*quantile=*/1-cl);
368  } else {
369  xMin = xv->getMin();
370  for (unsigned int j = 0; j < n; ++j) poiVals_[j] = poiVars_[j]->getVal();
371  double prob = ROOT::Math::chisquared_cdf_c(2*(nll.getVal() - nll0), n+nOtherFloatingPoi_);
372  if (commitPoints) Combine::commitPoint(true, /*quantile=*/prob);
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;
374  }
375 
376  for (unsigned int j = 0; j < n; ++j) poiVars_[j]->setVal(p0[j]);
377  double xMax = findCrossing(minimX, nll, *xv, threshold, p0[i], xv->getMax());
378  if (!std::isnan(xMax)) {
379  if (verbose > -1) std::cout << "Maximum of " << xv->GetName() << " at " << cl << " CL for all others floating is " << xMax << std::endl;
380  for (unsigned int j = 0; j < n; ++j) poiVals_[j] = poiVars_[j]->getVal();
381  if (commitPoints) Combine::commitPoint(true, /*quantile=*/1-cl);
382  } else {
383  xMax = xv->getMax();
384  double prob = ROOT::Math::chisquared_cdf_c(2*(nll.getVal() - nll0), n+nOtherFloatingPoi_);
385  for (unsigned int j = 0; j < n; ++j) poiVals_[j] = poiVars_[j]->getVal();
386  if (commitPoints) Combine::commitPoint(true, /*quantile=*/prob);
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;
388  }
389 
390  xv->setRange(name, xMin, xMax);
391  xv->setConstant(false);
392  }
393  verbose++; // restore verbosity
394 }
static std::vector< std::string > poi_
Definition: MultiDimFit.h:30
static unsigned int points_
Definition: MultiDimFit.h:38
static void addBranch(const char *name, void *address, const char *leaflist)
Add a branch to the output tree (for advanced use or debugging only)
Definition: Combine.cc:564
void applyOptionsBase(const boost::program_options::variables_map &vm)
int i
Definition: DBlmapReader.cc:9
static bool floatOtherPOIs_
Definition: MultiDimFit.h:39
static float deltaNLL_
Definition: MultiDimFit.h:35
void doGrid(RooAbsReal &nll)
Definition: MultiDimFit.cc:191
bool minimize(int verbose=0, bool cascade=true)
static unsigned int lastPoint_
Definition: MultiDimFit.h:38
bool withSystematics
Definition: Combine.cc:68
static std::vector< float > poiVals_
Definition: MultiDimFit.h:32
static bool do95_
void initOnce(RooWorkspace *w, RooStats::ModelConfig *mc_s)
Definition: MultiDimFit.cc:129
void doBox(RooAbsReal &nll, double cl, const char *name="box", bool commitPoints=true)
for each RooRealVar, set a range &#39;box&#39; from the PL profiling all other parameters ...
Definition: MultiDimFit.cc:345
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 ...
Definition: Combine.cc:557
static Algo algo_
Definition: MultiDimFit.h:28
double findCrossing(CascadeMinimizer &minim, RooAbsReal &nll, RooRealVar &r, double level, double rStart, double rBound)
const T & max(const T &a, const T &b)
bool isnan(float x)
Definition: math.h:13
T sqrt(T t)
Definition: SSEVec.h:46
virtual void applyOptions(const boost::program_options::variables_map &vm)
Definition: MultiDimFit.cc:53
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::auto_ptr< RooAbsReal > nll
int j
Definition: DBlmapReader.cc:9
static unsigned int nOtherFloatingPoi_
Definition: MultiDimFit.h:34
virtual bool runSpecific(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
Definition: MultiDimFit.cc:73
void doSingles(RooFitResult &res)
Definition: MultiDimFit.cc:153
static int minimizerStrategy_
float cl
Definition: Combine.cc:71
#define M_PI
Definition: BFit3D.cc:3
static std::vector< RooRealVar * > poiVars_
Definition: MultiDimFit.h:31
void doRandomPoints(RooAbsReal &nll)
Definition: MultiDimFit.cc:269
static RooArgList poiList_
Definition: MultiDimFit.h:33
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double a
Definition: hdecay.h:121
void doContour2D(RooAbsReal &nll)
Definition: MultiDimFit.cc:297
static unsigned int firstPoint_
Definition: MultiDimFit.h:38
tuple cout
Definition: gather_cfg.py:121
LimitAlgo * algo
Definition: Combine.cc:60
Definition: DDAxes.h:10
boost::program_options::options_description options_
Definition: LimitAlgo.h:31
static bool fastScan_
Definition: MultiDimFit.h:40
void setStrategy(int strategy)
T w() const