CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MarkovChainMC.cc
Go to the documentation of this file.
1 #include "../interface/MarkovChainMC.h"
2 #include <stdexcept>
3 #include <cmath>
4 #include "TKey.h"
5 #include "RooRealVar.h"
6 #include "RooArgSet.h"
7 #include "RooUniform.h"
8 #include "RooWorkspace.h"
9 #include "RooFitResult.h"
10 #include "RooRandom.h"
11 #ifndef ROOT_THnSparse
12 class THnSparse;
13 #define ROOT_THnSparse
14 #endif
15 #include "RooStats/MCMCCalculator.h"
16 #include "RooStats/MCMCInterval.h"
17 #include "RooStats/ModelConfig.h"
18 #include "RooStats/ProposalHelper.h"
19 #include "RooStats/ProposalFunction.h"
20 #include "RooStats/RooStatsUtils.h"
21 #include "../interface/Combine.h"
22 #include "../interface/TestProposal.h"
23 #include "../interface/DebugProposal.h"
24 #include "../interface/CloseCoutSentry.h"
25 #include "../interface/RooFitGlobalKillSentry.h"
26 #include "../interface/JacknifeQuantile.h"
27 
28 #include "../interface/ProfilingTools.h"
29 
30 using namespace RooStats;
31 
32 std::string MarkovChainMC::proposalTypeName_ = "ortho";
34 bool MarkovChainMC::runMinos_ = false;
35 bool MarkovChainMC::noReset_ = false;
37 bool MarkovChainMC::updateHint_ = false;
38 unsigned int MarkovChainMC::iterations_ = 10000;
39 unsigned int MarkovChainMC::burnInSteps_ = 200;
42 unsigned int MarkovChainMC::tries_ = 10;
46 bool MarkovChainMC::saveChain_ = false;
47 bool MarkovChainMC::noSlimChain_ = false;
48 bool MarkovChainMC::mergeChains_ = false;
49 bool MarkovChainMC::readChains_ = false;
55 
57  LimitAlgo("Markov Chain MC specific options")
58 {
59  options_.add_options()
60  ("iteration,i", boost::program_options::value<unsigned int>(&iterations_)->default_value(iterations_), "Number of iterations")
61  ("tries", boost::program_options::value<unsigned int>(&tries_)->default_value(tries_), "Number of times to run the MCMC on the same data")
62  ("burnInSteps,b", boost::program_options::value<unsigned int>(&burnInSteps_)->default_value(burnInSteps_), "Burn in steps (absolute number)")
63  ("burnInFraction", boost::program_options::value<float>(&burnInFraction_)->default_value(burnInFraction_), "Burn in steps (fraction of total accepted steps)")
64  ("adaptiveBurnIn", boost::program_options::value<bool>(&adaptiveBurnIn_)->default_value(adaptiveBurnIn_), "Adaptively determine burn in steps (experimental!).")
65  ("proposal", boost::program_options::value<std::string>(&proposalTypeName_)->default_value(proposalTypeName_),
66  "Proposal function to use: 'fit', 'uniform', 'gaus', 'ortho' (also known as 'test')")
67  ("runMinos", "Run MINOS when fitting the data")
68  ("noReset", "Don't reset variable state after fit")
69  ("updateHint", "Update hint with the results")
70  ("updateProposalParams",
71  boost::program_options::value<bool>(&updateProposalParams_)->default_value(updateProposalParams_),
72  "Control ProposalHelper::SetUpdateProposalParameters")
73  ("propHelperWidthRangeDivisor",
74  boost::program_options::value<float>(&proposalHelperWidthRangeDivisor_)->default_value(proposalHelperWidthRangeDivisor_),
75  "Sets the fractional size of the gaussians in the proposal")
76  ("alwaysStepPOI", boost::program_options::value<bool>(&alwaysStepPoi_)->default_value(alwaysStepPoi_),
77  "When using 'ortho' proposal, always step also the parameter of interest. On by default, as it improves convergence, but you can turn it off (e.g. if you turn off --optimizeSimPdf)")
78  ("propHelperUniformFraction",
79  boost::program_options::value<float>(&proposalHelperUniformFraction_)->default_value(proposalHelperUniformFraction_),
80  "Add a fraction of uniform proposals to the algorithm")
81  ("debugProposal", boost::program_options::value<int>(&debugProposal_)->default_value(debugProposal_), "Printout the first N proposals")
82  ("cropNSigmas",
83  boost::program_options::value<float>(&cropNSigmas_)->default_value(cropNSigmas_),
84  "crop range of all parameters to N times their uncertainty")
85  ("truncatedMeanFraction",
86  boost::program_options::value<float>(&truncatedMeanFraction_)->default_value(truncatedMeanFraction_),
87  "Discard this fraction of the results before computing the mean and rms")
88  ("adaptiveTruncation", boost::program_options::value<bool>(&adaptiveTruncation_)->default_value(adaptiveTruncation_),
89  "When averaging multiple runs, ignore results that are more far away from the median than the inter-quartile range")
90  ("hintSafetyFactor",
91  boost::program_options::value<float>(&hintSafetyFactor_)->default_value(hintSafetyFactor_),
92  "set range of integration equal to this number of times the hinted limit")
93  ("saveChain", "Save MarkovChain to output file")
94  ("noSlimChain", "Include also nuisance parameters in the chain that is saved to file")
95  ("mergeChains", "Merge MarkovChains instead of averaging limits")
96  ("readChains", "Just read MarkovChains from toysFile instead of running MCMC directly")
97  ;
98 }
99 
100 void MarkovChainMC::applyOptions(const boost::program_options::variables_map &vm) {
101  if (proposalTypeName_ == "fit") proposalType_ = FitP;
102  else if (proposalTypeName_ == "uniform") proposalType_ = UniformP;
103  else if (proposalTypeName_ == "gaus") proposalType_ = MultiGaussianP;
104  else if (proposalTypeName_ == "ortho") proposalType_ = TestP;
105  else if (proposalTypeName_ == "test") proposalType_ = TestP;
106  else {
107  std::cerr << "MarkovChainMC: proposal type " << proposalTypeName_ << " not known." << "\n" << options_ << std::endl;
108  throw std::invalid_argument("MarkovChainMC: unsupported proposal");
109  }
110 
111  runMinos_ = vm.count("runMinos");
112  noReset_ = vm.count("noReset");
113  updateHint_ = vm.count("updateHint");
114 
115  mass_ = vm["mass"].as<float>();
116  saveChain_ = vm.count("saveChain");
117  noSlimChain_ = vm.count("noSlimChain");
118  mergeChains_ = vm.count("mergeChains");
119  readChains_ = vm.count("readChains");
120 
121  if (mergeChains_ && !saveChain_ && !readChains_) chains_.SetOwner(true);
122 }
123 
124 bool MarkovChainMC::run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) {
126  std::cerr << "Sorry, the multi-gaussian proposal does not work without systematics.\n" <<
127  "Uniform proposal will be used instead.\n" << std::endl;
129  }
130 
132 
133  CloseCoutSentry coutSentry(verbose <= 0); // close standard output and error, so that we don't flood them with minuit messages
134 
135  // Get degrees of freedom
136  modelNDF_ = mc_s->GetParametersOfInterest()->getSize();
137  if (withSystematics) modelNDF_ += mc_s->GetNuisanceParameters()->getSize();
138 
139  double suma = 0; int num = 0;
140  double savhint = (hint ? *hint : -1); const double *thehint = hint;
141  std::vector<double> limits;
142  if (readChains_) {
143  readChains(*mc_s->GetParametersOfInterest(), limits);
144  } else {
145  for (unsigned int i = 0; i < tries_; ++i) {
146  if (int nacc = runOnce(w,mc_s,mc_b,data,limit,limitErr,thehint)) {
147  suma += nacc;
148  if (verbose > 1) std::cout << "Limit from this run: " << limit << std::endl;
149  limits.push_back(limit);
150  if (updateHint_ && tries_ > 1 && limit > savhint) {
151  if (verbose > 0) std::cout << "Updating hint from " << savhint << " to " << limit << std::endl;
152  savhint = limit; thehint = &savhint;
153  }
154  }
155  }
156  }
157  num = limits.size();
158  if (num == 0) return false;
159  // average acceptance
160  suma = suma / (num * double(iterations_));
161  limitAndError(limit, limitErr, limits);
162  if (mergeChains_) {
163  std::cout << "Limit from averaging: " << limit << " +/- " << limitErr << std::endl;
164  // copy constructors don't work, so we just have to leak memory :-(
165  RooStats::MarkovChain *merged = mergeChains(*mc_s->GetParametersOfInterest(), limits);
166  // set burn-in to zero, since steps have already been discarded when merging
167  limitFromChain(limit, limitErr, *mc_s->GetParametersOfInterest(), *merged, 0);
168  std::cout << "Limit from merged chain: " << limit << " +/- " << limitErr << std::endl;
169  }
170  coutSentry.clear();
171 
172  if (verbose >= 0) {
173  std::cout << "\n -- MarkovChainMC -- " << "\n";
174  RooRealVar *r = dynamic_cast<RooRealVar *>(mc_s->GetParametersOfInterest()->first());
175  if (num > 1) {
176  std::cout << "Limit: " << r->GetName() <<" < " << limit << " +/- " << limitErr << " @ " << cl * 100 << "% CL (" << num << " tries)" << std::endl;
177  if (verbose > 0 && !readChains_) std::cout << "Average chain acceptance: " << suma << std::endl;
178  } else {
179  std::cout << "Limit: " << r->GetName() <<" < " << limit << " @ " << cl * 100 << "% CL" << std::endl;
180  }
181  }
182  return true;
183 }
184 int MarkovChainMC::runOnce(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) const {
185  RooArgList poi(*mc_s->GetParametersOfInterest());
186  RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
187 
188  if ((hint != 0) && (*hint > r->getMin())) {
189  r->setMax(hintSafetyFactor_*(*hint));
190  }
191 
192  if (withSystematics && (mc_s->GetNuisanceParameters() == 0)) {
193  throw std::logic_error("MarkovChainMC: running with systematics enabled, but nuisances not defined.");
194  }
195 
196  w->loadSnapshot("clean");
197  std::auto_ptr<RooFitResult> fit(0);
198  if (proposalType_ == FitP || (cropNSigmas_ > 0)) {
199  CloseCoutSentry coutSentry(verbose <= 1); // close standard output and error, so that we don't flood them with minuit messages
200  fit.reset(mc_s->GetPdf()->fitTo(data, RooFit::Save(), RooFit::Minos(runMinos_)));
201  coutSentry.clear();
202  if (fit.get() == 0) { std::cerr << "Fit failed." << std::endl; return false; }
203  if (verbose > 1) fit->Print("V");
204  if (!noReset_) w->loadSnapshot("clean");
205  }
206 
207  if (cropNSigmas_ > 0) {
208  const RooArgList &fpf = fit->floatParsFinal();
209  for (int i = 0, n = fpf.getSize(); i < n; ++i) {
210  RooRealVar *fv = dynamic_cast<RooRealVar *>(fpf.at(i));
211  if (std::string(r->GetName()) == fv->GetName()) continue;
212  RooRealVar *v = w->var(fv->GetName());
213  double min = v->getMin(), max = v->getMax();
214  if (fv->hasAsymError(false)) {
215  min = (std::max(v->getMin(), fv->getVal() + cropNSigmas_ * fv->getAsymErrorLo()));
216  max = (std::min(v->getMax(), fv->getVal() + cropNSigmas_ * fv->getAsymErrorHi()));
217  } else if (fv->hasError(false)) {
218  min = (std::max(v->getMin(), fv->getVal() - cropNSigmas_ * fv->getError()));
219  max = (std::min(v->getMax(), fv->getVal() + cropNSigmas_ * fv->getError()));
220  }
221  if (verbose > 1) {
222  std::cout << " " << fv->GetName() << "[" << v->getMin() << ", " << v->getMax() << "] -> [" << min << ", " << max << "]" << std::endl;
223  }
224  v->setMin(min); v->setMax(max);
225  }
226  }
227 
228  std::auto_ptr<ProposalFunction> ownedPdfProp;
229  ProposalFunction* pdfProp = 0;
230  ProposalHelper ph;
231  switch (proposalType_) {
232  case UniformP:
233  if (verbose) std::cout << "Using uniform proposal" << std::endl;
234  ownedPdfProp.reset(new UniformProposal());
235  pdfProp = ownedPdfProp.get();
236  break;
237  case FitP:
238  if (verbose) std::cout << "Using fit proposal" << std::endl;
239  ph.SetVariables(fit->floatParsFinal());
240  ph.SetCovMatrix(fit->covarianceMatrix());
241  pdfProp = ph.GetProposalFunction();
242  break;
243  case MultiGaussianP:
244  if (verbose) std::cout << "Using multi-gaussian proposal" << std::endl;
245  ph.SetVariables(*mc_s->GetNuisanceParameters());
246  ph.SetWidthRangeDivisor(proposalHelperWidthRangeDivisor_);
247  pdfProp = ph.GetProposalFunction();
248  break;
249  case TestP:
250  ownedPdfProp.reset(new TestProposal(proposalHelperWidthRangeDivisor_, alwaysStepPoi_ ? poi : RooArgList()));
251  pdfProp = ownedPdfProp.get();
252  break;
253  }
254  if (proposalType_ != UniformP) {
255  ph.SetUpdateProposalParameters(updateProposalParams_);
256  if (proposalHelperUniformFraction_ > 0) ph.SetUniformFraction(proposalHelperUniformFraction_);
257  }
258 
259  std::auto_ptr<DebugProposal> pdfDebugProp(debugProposal_ > 0 ? new DebugProposal(pdfProp, mc_s->GetPdf(), &data, debugProposal_) : 0);
260 
261  MCMCCalculator mc(data, *mc_s);
262  mc.SetNumIters(iterations_);
263  mc.SetConfidenceLevel(cl);
264  mc.SetNumBurnInSteps(burnInSteps_);
265  mc.SetProposalFunction(debugProposal_ > 0 ? *pdfDebugProp : *pdfProp);
266  mc.SetLeftSideTailFraction(0);
267 
268  if (typeid(*mc_s->GetPriorPdf()) == typeid(RooUniform)) {
269  mc.SetPriorPdf(*((RooAbsPdf *)0));
270  }
271 
272  std::auto_ptr<MCMCInterval> mcInt;
273  try {
274  mcInt.reset((MCMCInterval*)mc.GetInterval());
275  } catch (std::length_error &ex) {
276  mcInt.reset(0);
277  }
278  if (mcInt.get() == 0) return false;
279  if (adaptiveBurnIn_) {
280  mcInt->SetNumBurnInSteps(guessBurnInSteps(*mcInt->GetChain()));
281  } else if (mcInt->GetChain()->Size() * burnInFraction_ > burnInSteps_) {
282  mcInt->SetNumBurnInSteps(mcInt->GetChain()->Size() * burnInFraction_);
283  }
284  limit = mcInt->UpperLimit(*r);
285 
286  if (saveChain_ || mergeChains_) {
287  // Copy-constructors don't work properly, so we just have to leak memory.
288  //RooStats::MarkovChain *chain = new RooStats::MarkovChain(*mcInt->GetChain());
289  RooStats::MarkovChain *chain = slimChain(*mc_s->GetParametersOfInterest(), *mcInt->GetChain());
290  if (mergeChains_) chains_.Add(chain);
291  if (saveChain_) writeToysHere->WriteTObject(chain, TString::Format("MarkovChain_mh%g_%u",mass_, RooRandom::integer(std::numeric_limits<UInt_t>::max() - 1)));
292  return chain->Size();
293  } else {
294  return mcInt->GetChain()->Size();
295  }
296 }
297 
298 void MarkovChainMC::limitAndError(double &limit, double &limitErr, const std::vector<double> &limitsIn) const {
299  std::vector<double> limits(limitsIn);
300  int num = limits.size();
301  // possibly remove outliers before computing mean
302  if (adaptiveTruncation_ && num >= 10) {
303  std::sort(limits.begin(), limits.end());
304  // determine location and size of the sample
305  double median = (num % 2 ? limits[num/2] : 0.5*(limits[num/2] + limits[num/2+1]));
306  double iqr = limits[3*num/4] - limits[num/4];
307  // determine range of plausible values
308  double min = median - iqr, max = median + iqr;
309  int start = 0, end = num-1;
310  while (start < end && limits[start] < min) ++start;
311  while (end > start && limits[end] > max) --end;
312  num = end-start+1;
313  // compute mean and rms of accepted part
314  limit = 0; limitErr = 0;
315  for (int k = start; k <= end; ++k) limit += limits[k];
316  limit /= num;
317  for (int k = start; k <= end; k++) limitErr += (limits[k]-limit)*(limits[k]-limit);
318  limitErr = (num > 1 ? sqrt(limitErr/(num*(num-1))) : 0);
319  std::cout << "Result from truncated mean: " << limit << " +/- " << limitErr << std::endl;
320 #if 0
321  QuantileCalculator qc(limits);
322  std::pair<double,double> qn = qc.quantileAndError(0.5, QuantileCalculator::Simple);
323  std::pair<double,double> qs = qc.quantileAndError(0.5, QuantileCalculator::Sectioning);
324  std::pair<double,double> qj = qc.quantileAndError(0.5, QuantileCalculator::Jacknife);
325  std::cout << "Median of limits (simple): " << qn.first << " +/- " << qn.second << std::endl;
326  std::cout << "Median of limits (sectioning): " << qs.first << " +/- " << qs.second << std::endl;
327  std::cout << "Median of limits (jacknife): " << qj.first << " +/- " << qj.second << std::endl;
328 #endif
329  } else {
330  int noutl = floor(truncatedMeanFraction_ * num);
331  if (noutl >= 1) {
332  std::sort(limits.begin(), limits.end());
333  double median = (num % 2 ? limits[num/2] : 0.5*(limits[num/2] + limits[num/2+1]));
334  for (int k = 0; k < noutl; ++k) {
335  // make sure outliers are all at the end
336  if (std::abs(limits[0]-median) > std::abs(limits[num-k-1]-median)) {
337  std::swap(limits[0], limits[num-k-1]);
338  }
339  }
340  num -= noutl;
341  }
342  // compute mean and rms
343  limit = 0; limitErr = 0;
344  for (int k = 0; k < num; k++) limit += limits[k];
345  limit /= num;
346  for (int k = 0; k < num; k++) limitErr += (limits[k]-limit)*(limits[k]-limit);
347  limitErr = (num > 1 ? sqrt(limitErr/(num*(num-1))) : 0);
348  }
349 }
350 
351 RooStats::MarkovChain *MarkovChainMC::mergeChains(const RooArgSet &poi, const std::vector<double> &limits) const
352 {
353  std::vector<double> limitsSorted(limits); std::sort(limitsSorted.begin(), limitsSorted.end());
354  double lmin = limitsSorted.front(), lmax = limitsSorted.back();
355  if (limitsSorted.size() > 5) {
356  int n = limitsSorted.size();
357  double lmedian = limitsSorted[n/2];
358  lmin = lmedian - 2*(+lmedian - limitsSorted[1*n/4]);
359  lmax = lmedian + 2*(-lmedian + limitsSorted[3*n/4]);
360  }
361  if (chains_.GetSize() == 0) throw std::runtime_error("No chains to merge");
362  if (verbose > 1) std::cout << "Will merge " << chains_.GetSize() << " chains." << std::endl;
363  RooArgSet pars(poi);
364  RooStats::MarkovChain *merged = new RooStats::MarkovChain("Merged","",pars);
365  TIter iter(&chains_);
366  int index = 0;
367  for (RooStats::MarkovChain *other = (RooStats::MarkovChain *) iter.Next();
368  other != 0;
369  other = (RooStats::MarkovChain *) iter.Next(), ++index) {
370  if (limits[index] < lmin || limits[index] > lmax) continue;
371  int burninSteps = adaptiveBurnIn_ ? guessBurnInSteps(*other) : max<int>(burnInSteps_, other->Size() * burnInFraction_);
372  if (verbose > 1) std::cout << "Adding chain of " << other->Size() << " entries, skipping the first " << burninSteps << "; individual limit " << limits[index] << std::endl;
373  for (int i = burninSteps, n = other->Size(); i < n; ++i) {
374  RooArgSet point(*other->Get(i));
375  double nllval = other->NLL();
376  double weight = other->Weight();
377  merged->Add(point,nllval,weight);
378  if (verbose > 2 && (i % 500 == 0)) std::cout << " added " << i << "/" << other->Size() << " entries." << std::endl;
379  }
380  }
381  return merged;
382 }
383 
384 void MarkovChainMC::readChains(const RooArgSet &poi, std::vector<double> &limits)
385 {
386  double mylim, myerr;
387  chains_.Clear();
388  chains_.SetOwner(false);
389  if (!readToysFromHere) throw std::logic_error("Cannot use readChains: option toysFile not specified, or input file empty");
390  TDirectory *toyDir = readToysFromHere->GetDirectory("toys");
391  if (!toyDir) throw std::logic_error("Cannot use readChains: empty toy dir in input file empty");
392  TString prefix = TString::Format("MarkovChain_mh%g_",mass_);
393  TIter next(toyDir->GetListOfKeys()); TKey *k;
394  while ((k = (TKey *) next()) != 0) {
395  if (TString(k->GetName()).Index(prefix) != 0) continue;
396  RooStats::MarkovChain *toy = dynamic_cast<RooStats::MarkovChain *>(toyDir->Get(k->GetName()));
397  if (toy == 0) continue;
398  limitFromChain(mylim, myerr, poi, *toy);
399  if (verbose > 1) std::cout << " limit " << mylim << " +/- " << myerr << std::endl;
400  // vvvvv ---- begin convergence test, still being developed, not recommended yet.
401  if (runtimedef::get("MCMC_STATIONARITY")) {
402  if (!stationarityTest(*toy, poi, 30)) {
403  if (verbose > 1) std::cout << " ---> rejecting chain!" << std::endl;
404  continue;
405  }
406  }
407  // ^^^^^ ---- end of convergence test
408  chains_.Add(toy);
409  limits.push_back(mylim);
410  }
411  if (verbose) { std::cout << "Read " << chains_.GetSize() << " Markov Chains from input file." << std::endl; }
412 }
413 
414 void
415 MarkovChainMC::limitFromChain(double &limit, double &limitErr, const RooArgSet &poi, RooStats::MarkovChain &chain, int burnInSteps)
416 {
417  if (burnInSteps < 0) {
418  if (adaptiveBurnIn_) burnInSteps = guessBurnInSteps(chain);
419  else burnInSteps = max<int>(burnInSteps_, chain.Size() * burnInFraction_);
420  }
421 #if 1 // This is much faster and gives the same result
422  QuantileCalculator qc(*chain.GetAsConstDataSet(), poi.first()->GetName(), burnInSteps);
424 #else
425  MCMCInterval interval("",poi,chain);
426  RooArgList axes(poi);
427  interval.SetConfidenceLevel(cl);
428  interval.SetIntervalType(MCMCInterval::kTailFraction);
429  interval.SetLeftSideTailFraction(0);
430  interval.SetNumBurnInSteps(burnInSteps);
431  interval.SetAxes(axes);
432  limit = interval.UpperLimit((RooRealVar&)*poi.first());
433  if (mergeChains_) {
434  // must avoid that MCMCInterval deletes the chain
435  interval.SetChain(*(RooStats::MarkovChain *)0);
436  }
437 #endif
438 }
439 
440 RooStats::MarkovChain *
441 MarkovChainMC::slimChain(const RooArgSet &poi, const RooStats::MarkovChain &chain) const
442 {
443  RooArgSet poilvalue(poi); // wtf they don't take a const & ??
444  if (noSlimChain_) poilvalue.add(*chain.Get());
445  RooStats::MarkovChain * ret = new RooStats::MarkovChain("","",poilvalue);
446  for (int i = 0, n = chain.Size(); i < n; ++i) {
447  RooArgSet entry(*chain.Get(i));
448  Double_t nll = chain.NLL();
449  Double_t weight = chain.Weight();
450  if (i) ret->AddFast(entry, nll, weight);
451  else ret->Add(entry, nll, weight);
452  }
453  return ret;
454 }
455 
456 int
457 MarkovChainMC::guessBurnInSteps(const RooStats::MarkovChain &chain) const
458 {
459  int n = chain.Size();
460  std::vector<double> nll(n);
461  for (int i = 0; i < n; ++i) {
462  nll[i] = chain.NLL(i);
463  }
464  // get minimum of nll
465  double minnll = nll[0];
466  for (int i = 0; i < n; ++i) { if (nll[i] < minnll) minnll = nll[i]; }
467  // subtract it from all
468  for (int i = 0; i < n; ++i) nll[i] -= minnll;
469  // the NLL looks like 0.5 * a chi2 with nparam degrees of freedom, plus an arbitrary constant
470  // so it should have a sigma of 0.5*sqrt(2*nparams)
471  // anything sensible should be between minimum and minimum + 10*sigma
472  double maxcut = 5*sqrt(2.0*modelNDF_);
473  // now loop backwards until you find something that's outside
474  int start = 0;
475  for (start = n-1; start >= 0; --start) {
476  if (nll[start] > maxcut) break;
477  }
478  return start;
479 }
480 
481 int
482 MarkovChainMC::stationarityTest(const RooStats::MarkovChain &chain, const RooArgSet &poi, int nchunks) const
483 {
484  std::vector<int> entries(nchunks, 0);
485  std::vector<double> mean(nchunks, .0);
486  const RooDataSet *data = chain.GetAsConstDataSet();
487  const RooRealVar *r = dynamic_cast<const RooRealVar *>(data->get()->find(poi.first()->GetName()));
488  int n = data->numEntries(), chunksize = ceil(n/double(nchunks));
489  for (int i = 0, chunk = 0; i < n; i++) {
490  data->get(i);
491  if (i > 0 && i % chunksize == 0) chunk++;
492  entries[chunk]++;
493  mean[chunk] += r->getVal();
494  }
495  for (int c = 0; c < nchunks; ++c) { mean[c] /= entries[c]; }
496 
497  std::vector<double> dists, dists25;
498  for (int c = 0; c < nchunks; ++c) {
499  for (int c2 = 0; c2 < nchunks; ++c2) {
500  if (c2 != c) dists.push_back(fabs(mean[c]-mean[c2]));
501  }
502  std::sort(dists.begin(), dists.end());
503  dists25.push_back(dists[ceil(0.25*nchunks)]/mean[c]);
504  dists.clear();
505  //printf("chunk %3d: mean %9.5f dist25 %9.5f abs, %9.5f real\n", c, mean[c], mean[c]*dists25.back(), dists25.back());
506  }
507  std::sort(dists25.begin(), dists25.end());
508  double tolerance = 10*dists25[ceil(0.25*nchunks)];
509  //printf("Q(25) is %9.5f\n", dists25[ceil(0.25*nchunks)]);
510  //printf("Q(50) is %9.5f\n", dists25[ceil(0.50*nchunks)]);
511  //printf("Tolerance set to %9.5f\n", tolerance);
512  bool converged = true;
513  std::vector<int> friends(nchunks, 0), foes(nchunks, 0);
514  for (int c = 0; c < nchunks; ++c) {
515  for (int c2 = c+1; c2 < nchunks; ++c2) {
516  if (c2 == c) continue;
517  if (fabs(mean[c] - mean[c2]) < tolerance*mean[c]) {
518  friends[c]++;
519  } else {
520  foes[c]++;
521  }
522  }
523  //printf("chunk %3d: mean %9.5f friends %3d foes %3d \n", c, mean[c], friends[c], foes[c]);
524  if (friends[c] >= 2 && foes[c] > 1) {
525  converged = false;
526  }
527  //fflush(stdout);
528  }
529  return converged;
530 }
static bool updateProposalParams_
Definition: MarkovChainMC.h:29
int i
Definition: DBlmapReader.cc:9
TDirectory * readToysFromHere
Definition: Combine.cc:66
static float cropNSigmas_
Definition: MarkovChainMC.h:63
static ProposalType proposalType_
Definition: MarkovChainMC.h:28
void limitFromChain(double &limit, double &limitErr, const RooArgSet &poi, RooStats::MarkovChain &chain, int burnInSteps=-1)
static float proposalHelperWidthRangeDivisor_
Definition: MarkovChainMC.h:62
static bool adaptiveTruncation_
do adaptive truncated mean
Definition: MarkovChainMC.h:43
static float proposalHelperUniformFraction_
Definition: MarkovChainMC.h:62
tuple interval
Definition: MergeJob_cfg.py:20
static float hintSafetyFactor_
Safety factor for hint (integrate up to this number of times the hinted limit)
Definition: MarkovChainMC.h:45
int get(const char *name)
static bool readChains_
Read chains from file instead of running them.
Definition: MarkovChainMC.h:53
long int integer
Definition: mlp_lapack.h:12
static const int WARNING
bool withSystematics
Definition: Combine.cc:68
#define abs(x)
Definition: mlp_lapack.h:159
static float truncatedMeanFraction_
Ignore up to this fraction of results if they&#39;re too far from the median.
Definition: MarkovChainMC.h:41
#define min(a, b)
Definition: mlp_lapack.h:161
static unsigned int iterations_
Propose this number of points for the chain.
Definition: MarkovChainMC.h:31
static bool alwaysStepPoi_
Definition: MarkovChainMC.h:61
TDirectory * writeToysHere
Definition: Combine.cc:65
RooStats::MarkovChain * mergeChains(const RooArgSet &poi, const std::vector< double > &limits) const
float mass_
Mass of the Higgs boson (goes into the name of the saved chains)
Definition: MarkovChainMC.h:55
static bool runMinos_
Definition: MarkovChainMC.h:29
const T & max(const T &a, const T &b)
static float burnInFraction_
Discard these fraction of points.
Definition: MarkovChainMC.h:35
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:46
static std::string proposalTypeName_
Definition: MarkovChainMC.h:27
std::pair< std::string, MonitorElement * > entry
Definition: ME_MAP.h:8
static bool adaptiveBurnIn_
Adaptive burn-in (experimental!)
Definition: MarkovChainMC.h:37
static bool saveChain_
Save Markov Chain in output file.
Definition: MarkovChainMC.h:47
static int debugProposal_
Definition: MarkovChainMC.h:64
std::pair< double, double > quantileAndError(double quantile, Method method)
bool run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
static bool updateHint_
Definition: MarkovChainMC.h:29
#define end
Definition: vmac.h:38
int stationarityTest(const RooStats::MarkovChain &chain, const RooArgSet &poi, int nchunks) const
note: this is still being developed
static bool noReset_
Definition: MarkovChainMC.h:29
void limitAndError(double &limit, double &limitErr, const std::vector< double > &limits) const
float cl
Definition: Combine.cc:71
int k[5][pyjets_maxn]
int runOnce(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint) const
static bool noSlimChain_
Leave all parameters in the markov chain, not just the POI.
Definition: MarkovChainMC.h:49
RooStats::MarkovChain * slimChain(const RooArgSet &poi, const RooStats::MarkovChain &chain) const
static unsigned int burnInSteps_
Discard these points.
Definition: MarkovChainMC.h:33
long long int num
Definition: procUtils.cc:71
virtual void applyOptions(const boost::program_options::variables_map &vm)
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
int modelNDF_
Number of degrees of freedom of the problem, approximately.
Definition: MarkovChainMC.h:57
static unsigned int tries_
compute the limit N times
Definition: MarkovChainMC.h:39
int guessBurnInSteps(const RooStats::MarkovChain &chain) const
tuple cout
Definition: gather_cfg.py:121
boost::program_options::options_description options_
Definition: LimitAlgo.h:31
static bool mergeChains_
Merge chains instead of averaging limits.
Definition: MarkovChainMC.h:51
mathSSE::Vec4< T > v
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
T w() const
void readChains(const RooArgSet &poi, std::vector< double > &limits)