CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CachingNLL.cc
Go to the documentation of this file.
1 #include "../interface/CachingNLL.h"
2 #include "../interface/utils.h"
3 #include <stdexcept>
4 #include <RooCategory.h>
5 #include <RooDataSet.h>
6 #include <RooProduct.h>
7 #include "../interface/ProfilingTools.h"
8 
9 //---- Uncomment this to get a '.' printed every some evals
10 //#define TRACE_NLL_EVALS
11 
12 //---- Uncomment this and run with --perfCounters to get cache statistics
13 //#define DEBUG_CACHE
14 
15 //---- Uncomment to enable Kahan's summation (if enabled at runtime with --X-rtd = ...
16 // http://en.wikipedia.org/wiki/Kahan_summation_algorithm
17 //#define ADDNLL_KAHAN_SUM
18 #include "../interface/ProfilingTools.h"
19 
20 //std::map<std::string,double> cacheutils::CachingAddNLL::offsets_;
23 
24 //#define DEBUG_TRACE_POINTS
25 #ifdef DEBUG_TRACE_POINTS
26 namespace {
27  template<unsigned int>
28  void tracePoint(const RooAbsCollection &point) {
29  static const RooAbsCollection *lastPoint = 0;
30  static std::vector<double> values;
31  if (&point != lastPoint) {
32  std::cout << "Arrived in a completely new point. " << std::endl;
33  values.resize(point.getSize());
34  RooLinkedListIter iter = point.iterator();
35  for (RooAbsArg *a = (RooAbsArg*)(iter.Next()); a != 0; a = (RooAbsArg*)(iter.Next())) {
36  RooRealVar *rrv = dynamic_cast<RooRealVar *>(a); if (!rrv) continue;
37  values.push_back(rrv->getVal());
38  }
39  lastPoint = &point;
40  } else {
41  std::cout << "Moved: ";
42  RooLinkedListIter iter = point.iterator();
43  int i = 0;
44  for (RooAbsArg *a = (RooAbsArg*)(iter.Next()); a != 0; a = (RooAbsArg*)(iter.Next())) {
45  RooRealVar *rrv = dynamic_cast<RooRealVar *>(a); if (!rrv) continue;
46  if (values[i] != rrv->getVal()) std::cout << a->GetName() << " " << values[i] << " => " << rrv->getVal() << " ";
47  values[i++] = rrv->getVal();
48  }
49  std::cout << std::endl;
50  }
51  }
52 }
53 #define TRACE_POINT2(x,i) ::tracePoint<i>(x);
54 #define TRACE_POINT(x) ::tracePoint<0>(x);
55 #define TRACE_NLL(x) std::cout << x << std::endl;
56 #else
57 #define TRACE_POINT2(x,i)
58 #define TRACE_POINT(x)
59 #define TRACE_NLL(x)
60 #endif
61 
63 {
64  std::auto_ptr<TIterator> iter(set.createIterator());
65  for (RooAbsArg *a = dynamic_cast<RooAbsArg *>(iter->Next());
66  a != 0;
67  a = dynamic_cast<RooAbsArg *>(iter->Next())) {
68  RooRealVar *rrv = dynamic_cast<RooRealVar *>(a);
69  if (rrv) { // && !rrv->isConstant()) {
70  vars_.push_back(rrv);
71  vals_.push_back(rrv->getVal());
72  }
73  }
74 }
75 
76 bool
77 cacheutils::ArgSetChecker::changed(bool updateIfChanged)
78 {
79  std::vector<RooRealVar *>::const_iterator it = vars_.begin(), ed = vars_.end();
80  std::vector<double>::iterator itv = vals_.begin();
81  bool changed = false;
82  for ( ; it != ed; ++it, ++itv) {
83  double val = (*it)->getVal();
84  if (val != *itv) {
85  //std::cerr << "var::CachingPdfable " << (*it)->GetName() << " changed: " << *itv << " -> " << val << std::endl;
86  changed = true;
87  if (updateIfChanged) { *itv = val; }
88  else break;
89  }
90  }
91  return changed;
92 }
93 
94 cacheutils::ValuesCache::ValuesCache(const RooAbsCollection &params, int size) :
95  size_(1),
96  maxSize_(size)
97 {
98  assert(size <= MaxItems_);
99  items[0] = new Item(params);
100 }
101 cacheutils::ValuesCache::ValuesCache(const RooAbsReal &pdf, const RooArgSet &obs, int size) :
102  size_(1),
103  maxSize_(size)
104 {
105  assert(size <= MaxItems_);
106  std::auto_ptr<RooArgSet> params(pdf.getParameters(obs));
107  items[0] = new Item(*params);
108 }
109 
110 
112 {
113  for (int i = 0; i < size_; ++i) delete items[i];
114 }
115 
117 {
118  for (int i = 0; i < size_; ++i) items[i]->good = false;
119 }
120 
121 std::pair<std::vector<Double_t> *, bool> cacheutils::ValuesCache::get()
122 {
123  int found = -1; bool good = false;
124  for (int i = 0; i < size_; ++i) {
125  if (items[i]->good) {
126  // valid entry, check if fresh
127  if (!items[i]->checker.changed()) {
128 #ifdef DEBUG_CACHE
129  PerfCounter::add(i == 0 ? "ValuesCache::get hit first" : "ValuesCache::get hit other");
130 #endif
131  // fresh: done!
132  found = i;
133  good = true;
134  break;
135  }
136  } else if (found == -1) {
137  // invalid entry, can be replaced
138  found = i;
139 #ifdef DEBUG_CACHE
140  PerfCounter::add("ValuesCache::get hit invalid");
141 #endif
142  }
143  }
144  if (found == -1) {
145  // all entries are valid but old
146 #ifdef DEBUG_CACHE
147  PerfCounter::add("ValuesCache::get miss");
148 #endif
149  if (size_ < maxSize_) {
150  // if I can, make a new entry
151  items[size_] = new Item(items[0]->checker); // create a new item, copying the ArgSetChecker from the first one
152  found = size_;
153  size_++;
154  } else {
155  // otherwise, pick the last one
156  found = size_-1;
157  }
158  }
159  // make sure new entry is the first one
160  if (found != 0) {
161  // remember what found is pointing to
162  Item *f = items[found];
163  // shift the other items down one place
164  while (found > 0) { items[found] = items[found-1]; --found; }
165  // and put found on top
166  items[found] = f;
167  }
168  if (!good) items[found]->checker.changed(true); // store new values in cache sentry
169  items[found]->good = true; // mark this as valid entry
170  return std::pair<std::vector<Double_t> *, bool>(&items[found]->values, good);
171 }
172 
173 cacheutils::CachingPdf::CachingPdf(RooAbsReal *pdf, const RooArgSet *obs) :
174  obs_(obs),
175  pdfOriginal_(pdf),
176  pdfPieces_(),
177  pdf_(utils::fullCloneFunc(pdf, pdfPieces_)),
178  lastData_(0),
179  cache_(*pdf_,*obs_)
180 {
181 }
182 
184  obs_(other.obs_),
185  pdfOriginal_(other.pdfOriginal_),
186  pdfPieces_(),
187  pdf_(utils::fullCloneFunc(pdfOriginal_, pdfPieces_)),
188  lastData_(0),
189  cache_(*pdf_,*obs_)
190 {
191 }
192 
194 {
195 }
196 
197 const std::vector<Double_t> &
199 {
200 #ifdef DEBUG_CACHE
201  PerfCounter::add("CachingPdf::eval called");
202 #endif
203  bool newdata = (lastData_ != &data);
204  if (newdata) {
205  lastData_ = &data;
206  pdf_->optimizeCacheMode(*data.get());
207  pdf_->attachDataSet(data);
208  const_cast<RooAbsData*>(lastData_)->setDirtyProp(false);
209  cache_.clear();
210  }
211  std::pair<std::vector<Double_t> *, bool> hit = cache_.get();
212  if (!hit.second) {
213  realFill_(data, *hit.first);
214  }
215  return *hit.first;
216 }
217 
218 void
219 cacheutils::CachingPdf::realFill_(const RooAbsData &data, std::vector<Double_t> &vals)
220 {
221 #ifdef DEBUG_CACHE
222  PerfCounter::add("CachingPdf::realFill_ called");
223 #endif
224  int n = data.numEntries();
225  vals.resize(n); // should be a no-op if size is already >= n.
226  //std::auto_ptr<RooArgSet> params(pdf_->getObservables(*obs)); // for non-smart handling of pointers
227  std::vector<Double_t>::iterator itv = vals.begin();
228  for (int i = 0; i < n; ++i, ++itv) {
229  data.get(i);
230  //*params = *data.get(i); // for non-smart handling of pointers
231  *itv = pdf_->getVal(obs_);
232  //std::cout << " at i = " << i << " pdf = " << *itv << std::endl;
233  TRACE_NLL("PDF value for " << pdf_->GetName() << " is " << *itv << " at this point.")
234  TRACE_POINT2(*obs_,1)
235  }
236 }
237 
238 cacheutils::CachingAddNLL::CachingAddNLL(const char *name, const char *title, RooAbsPdf *pdf, RooAbsData *data) :
239  RooAbsReal(name, title),
240  pdf_(pdf),
241  params_("params","parameters",this),
242  zeroPoint_(0)
243 {
244  if (pdf == 0) throw std::invalid_argument(std::string("Pdf passed to ")+name+" is null");
245  setData(*data);
246  setup_();
247 }
248 
250  RooAbsReal(name ? name : (TString("nll_")+other.pdf_->GetName()).Data(), ""),
251  pdf_(other.pdf_),
252  params_("params","parameters",this),
253  zeroPoint_(0)
254 {
255  setData(*other.data_);
256  setup_();
257 }
258 
260 {
261  for (int i = 0, n = integrals_.size(); i < n; ++i) delete integrals_[i];
262  integrals_.clear();
263 }
264 
267 {
268  return new cacheutils::CachingAddNLL(*this, name);
269 }
270 
271 void
273 {
274  const RooArgSet *obs = data_->get();
275  for (int i = 0, n = integrals_.size(); i < n; ++i) delete integrals_[i];
276  integrals_.clear();
277  RooAddPdf *addpdf = 0;
278  RooRealSumPdf *sumpdf = 0;
279  if ((addpdf = dynamic_cast<RooAddPdf *>(pdf_)) != 0) {
280  isRooRealSum_ = false;
281  int npdf = addpdf->coefList().getSize();
282  coeffs_.reserve(npdf);
283  pdfs_.reserve(npdf);
284  for (int i = 0; i < npdf; ++i) {
285  RooAbsReal * coeff = dynamic_cast<RooAbsReal*>(addpdf->coefList().at(i));
286  RooAbsPdf * pdfi = dynamic_cast<RooAbsPdf *>(addpdf->pdfList().at(i));
287  coeffs_.push_back(coeff);
288  pdfs_.push_back(CachingPdf(pdfi, obs));
289  }
290  } else if ((sumpdf = dynamic_cast<RooRealSumPdf *>(pdf_)) != 0) {
291  isRooRealSum_ = true;
292  int npdf = sumpdf->coefList().getSize();
293  coeffs_.reserve(npdf);
294  pdfs_.reserve(npdf);
295  integrals_.reserve(npdf);
296  for (int i = 0; i < npdf; ++i) {
297  RooAbsReal * coeff = dynamic_cast<RooAbsReal*>(sumpdf->coefList().at(i));
298  RooAbsReal * funci = dynamic_cast<RooAbsReal*>(sumpdf->funcList().at(i));
300  if (0 && typeid(*funci) == typeid(RooProduct)) {
301  RooArgList obsDep, obsInd;
302  obsInd.add(*coeff);
303  utils::factorizeFunc(*obs, *funci, obsDep, obsInd);
304  std::cout << "Entry " << i << ": coef name " << (coeff ? coeff->GetName() : "null") <<
305  " type " << (coeff ? coeff->ClassName() : "n/a") << std::endl;
306  std::cout << " " << "; func name " << (funci ? funci->GetName() : "null") <<
307  " type " << (funci ? funci->ClassName() : "n/a") << std::endl;
308  std::cout << "Terms depending on observables: " << std::endl; obsDep.Print("V");
309  std::cout << "Terms not depending on observables: " << std::endl; obsInd.Print("V");
310  if (obsInd.getSize() > 1) {
311  coeff = new RooProduct(TString::Format("%s_x_%s_obsIndep", coeff->GetName(), funci->GetName()), "", RooArgSet(obsInd));
312  addOwnedComponents(RooArgSet(*coeff));
313  }
314  if (obsDep.getSize() > 1) {
315  funci = new RooProduct(TString::Format("%s_obsDep", funci->GetName()), "", RooArgSet(obsInd));
316  addOwnedComponents(RooArgSet(*funci));
317  } else if (obsDep.getSize() == 1) {
318  funci = (RooAbsReal *) obsDep.first();
319  } else throw std::logic_error("No part of pdf depends on observables?");
320  }
321  coeffs_.push_back(coeff);
322  pdfs_.push_back(CachingPdf(funci, obs));
323  integrals_.push_back(funci->createIntegral(*obs));
324  }
325  } else {
326  std::string errmsg = "ERROR: CachingAddNLL: Pdf ";
327  errmsg += pdf_->GetName();
328  errmsg += " is neither a RooAddPdf nor a RooRealSumPdf, but a ";
329  errmsg += pdf_->ClassName();
330  throw std::invalid_argument(errmsg);
331  }
332 
333  std::auto_ptr<RooArgSet> params(pdf_->getParameters(*data_));
334  std::auto_ptr<TIterator> iter(params->createIterator());
335  for (RooAbsArg *a = (RooAbsArg *) iter->Next(); a != 0; a = (RooAbsArg *) iter->Next()) {
336  RooRealVar *rrv = dynamic_cast<RooRealVar *>(a);
337  //if (rrv != 0 && !rrv->isConstant()) params_.add(*rrv);
338  if (rrv != 0) params_.add(*rrv);
339  }
340 }
341 
342 Double_t
344 {
345 #ifdef DEBUG_CACHE
346  PerfCounter::add("CachingAddNLL::evaluate called");
347 #endif
348  std::fill( partialSum_.begin(), partialSum_.end(), 0.0 );
349 
350  std::vector<RooAbsReal*>::iterator itc = coeffs_.begin(), edc = coeffs_.end();
351  std::vector<CachingPdf>::iterator itp = pdfs_.begin();//, edp = pdfs_.end();
352  std::vector<Double_t>::const_iterator itw, bgw = weights_.begin();//, edw = weights_.end();
353  std::vector<Double_t>::iterator its, bgs = partialSum_.begin(), eds = partialSum_.end();
354  double sumCoeff = 0;
355  //std::cout << "Performing evaluation of " << GetName() << std::endl;
356  for ( ; itc != edc; ++itp, ++itc ) {
357  // get coefficient
358  Double_t coeff = (*itc)->getVal();
359  if (isRooRealSum_) {
360  sumCoeff += coeff * integrals_[itc - coeffs_.begin()]->getVal();
361  //std::cout << " coefficient = " << coeff << ", integral = " << integrals_[itc - coeffs_.begin()]->getVal() << std::endl;
362  } else {
363  sumCoeff += coeff;
364  }
365  // get vals
366  const std::vector<Double_t> &pdfvals = itp->eval(*data_);
367  // update running sum
368  std::vector<Double_t>::const_iterator itv = pdfvals.begin();
369  for (its = bgs; its != eds; ++its, ++itv) {
370  *its += coeff * (*itv); // sum (n_i * pdf_i)
371  }
372  }
373  // then get the final nll
374  double ret = 0;
375  for ( its = bgs, itw = bgw ; its != eds ; ++its, ++itw ) {
376  if (*itw == 0) continue;
377  if (!isnormal(*its) || *its <= 0) {
378  std::cerr << "WARNING: underflow to " << *its << " in " << GetName() << std::endl;
379  if (!CachingSimNLL::noDeepLEE_) logEvalError("Number of events is negative or error"); else CachingSimNLL::hasError_ = true;
380  }
381  double thispiece = (*itw) * (*its <= 0 ? -9e9 : log( ((*its) / sumCoeff) ));
382  #ifdef ADDNLL_KAHAN_SUM
383  static bool do_kahan = runtimedef::get("ADDNLL_KAHAN_SUM");
384  if (do_kahan) {
385  double kahan_y = thispiece - compensation;
386  double kahan_t = ret + kahan_y;
387  double kahan_d = (kahan_t - ret);
388  compensation = kahan_d - kahan_y;
389  ret = kahan_t;
390  } else {
391  ret += thispiece;
392  }
393  #else
394  ret += thispiece;
395  #endif
396  }
397  // then flip sign
398  ret = -ret;
399  // std::cout << "AddNLL for " << pdf_->GetName() << ": " << ret << std::endl;
400  // and add extended term: expected - observed*log(expected);
401  double expectedEvents = (isRooRealSum_ ? pdf_->getNorm(data_->get()) : sumCoeff);
402  if (expectedEvents <= 0) {
403  if (!CachingSimNLL::noDeepLEE_) logEvalError("Expected number of events is negative"); else CachingSimNLL::hasError_ = true;
404  expectedEvents = 1e-6;
405  }
406  //ret += expectedEvents - UInt_t(sumWeights_) * log(expectedEvents); // no, doesn't work with Asimov dataset
407  ret += expectedEvents - sumWeights_ * log(expectedEvents);
408  ret += zeroPoint_;
409  // std::cout << " plus extended term: " << ret << std::endl;
410  TRACE_NLL("AddNLL for " << pdf_->GetName() << ": " << ret)
411  return ret;
412 }
413 
414 void
416 {
417  //std::cout << "Setting data for pdf " << pdf_->GetName() << std::endl;
418  //utils::printRAD(&data);
419  data_ = &data;
420  setValueDirty();
421  sumWeights_ = 0.0;
422  weights_.resize(data.numEntries());
423  partialSum_.resize(data.numEntries());
424  std::vector<Double_t>::iterator itw = weights_.begin();
425  #ifdef ADDNLL_KAHAN_SUM
426  double compensation = 0;
427  #endif
428  for (int i = 0, n = data.numEntries(); i < n; ++i, ++itw) {
429  data.get(i);
430  *itw = data.weight();
431  #ifdef ADDNLL_KAHAN_SUM
432  static bool do_kahan = runtimedef::get("ADDNLL_KAHAN_SUM");
433  if (do_kahan) {
434  double kahan_y = *itw - compensation;
435  double kahan_t = sumWeights_ + kahan_y;
436  double kahan_d = (kahan_t - sumWeights_);
437  compensation = kahan_d - kahan_y;
438  sumWeights_ = kahan_t;
439  } else {
440  sumWeights_ += *itw;
441  }
442  #else
443  sumWeights_ += *itw;
444  #endif
445  }
446  for (std::vector<CachingPdf>::iterator itp = pdfs_.begin(), edp = pdfs_.end(); itp != edp; ++itp) {
447  itp->setDataDirty();
448  }
449 }
450 
451 RooArgSet*
452 cacheutils::CachingAddNLL::getObservables(const RooArgSet* depList, Bool_t valueOnly) const
453 {
454  return new RooArgSet();
455 }
456 
457 RooArgSet*
458 cacheutils::CachingAddNLL::getParameters(const RooArgSet* depList, Bool_t stripDisconnected) const
459 {
460  return new RooArgSet(params_);
461 }
462 
463 
464 cacheutils::CachingSimNLL::CachingSimNLL(RooSimultaneous *pdf, RooAbsData *data, const RooArgSet *nuis) :
465  pdfOriginal_(pdf),
466  dataOriginal_(data),
467  nuis_(nuis),
468  params_("params","parameters",this)
469 {
470  setup_();
471 }
472 
474  pdfOriginal_(other.pdfOriginal_),
475  dataOriginal_(other.dataOriginal_),
476  nuis_(other.nuis_),
477  params_("params","parameters",this)
478 {
479  setup_();
480 }
481 
484 {
485  return new cacheutils::CachingSimNLL(*this, name);
486 }
487 
489 {
490  for (std::vector<SimpleGaussianConstraint*>::iterator it = constrainPdfsFast_.begin(), ed = constrainPdfsFast_.end(); it != ed; ++it) {
491  delete *it;
492  }
493 }
494 
495 void
497 {
498  // Allow runtime-flag to switch off logEvalErrors
499  noDeepLEE_ = runtimedef::get("SIMNLL_NO_LEE");
500 
501  //RooAbsPdf *pdfclone = runtimedef::get("SIMNLL_CLONE") ? pdfOriginal_ : utils::fullClonePdf(pdfOriginal_, piecesForCloning_);
502  RooAbsPdf *pdfclone = pdfOriginal_; // never clone
503 
504  //---- Instead of getting the parameters here, we get them from the individual constraint terms and single pdfs ----
505  //---- This seems to save memory.
506  //std::auto_ptr<RooArgSet> params(pdfclone->getParameters(*dataOriginal_));
507  //params_.add(*params);
508 
509  RooArgList constraints;
510  factorizedPdf_.reset(dynamic_cast<RooSimultaneous *>(utils::factorizePdf(*dataOriginal_->get(), *pdfclone, constraints)));
511 
512  RooSimultaneous *simpdf = factorizedPdf_.get();
513  constrainPdfs_.clear();
514  if (constraints.getSize()) {
515  bool FastConstraints = runtimedef::get("SIMNLL_FASTGAUSS");
516  //constrainPdfs_.push_back(new RooProdPdf("constraints","constraints", constraints));
517  for (int i = 0, n = constraints.getSize(); i < n; ++i) {
518  RooAbsPdf *pdfi = dynamic_cast<RooAbsPdf*>(constraints.at(i));
519  if (FastConstraints && typeid(*pdfi) == typeid(RooGaussian)) {
520  constrainPdfsFast_.push_back(new SimpleGaussianConstraint(dynamic_cast<const RooGaussian&>(*pdfi)));
521  constrainZeroPointsFast_.push_back(0);
522  } else {
523  constrainPdfs_.push_back(pdfi);
524  constrainZeroPoints_.push_back(0);
525  }
526  //std::cout << "Constraint pdf: " << constraints.at(i)->GetName() << std::endl;
527  std::auto_ptr<RooArgSet> params(pdfi->getParameters(*dataOriginal_));
528  params_.add(*params, false);
529  }
530  } else {
531  std::cerr << "PDF didn't factorize!" << std::endl;
532  std::cout << "Parameters: " << std::endl;
533  std::auto_ptr<RooArgSet> params(pdfclone->getParameters(*dataOriginal_));
534  params->Print("V");
535  std::cout << "Obs: " << std::endl;
536  dataOriginal_->get()->Print("V");
537  factorizedPdf_.release();
538  simpdf = dynamic_cast<RooSimultaneous *>(pdfclone);
539  }
540 
541 
542  std::auto_ptr<RooAbsCategoryLValue> catClone((RooAbsCategoryLValue*) simpdf->indexCat().Clone());
543  pdfs_.resize(catClone->numBins(NULL), 0);
544  //dataSets_.reset(dataOriginal_->split(pdfOriginal_->indexCat(), true));
545  datasets_.resize(pdfs_.size(), 0);
546  splitWithWeights(*dataOriginal_, simpdf->indexCat(), true);
547  //std::cout << "Pdf " << simpdf->GetName() <<" is a SimPdf over category " << catClone->GetName() << ", with " << pdfs_.size() << " bins" << std::endl;
548  for (int ib = 0, nb = pdfs_.size(); ib < nb; ++ib) {
549  catClone->setBin(ib);
550  RooAbsPdf *pdf = simpdf->getPdf(catClone->getLabel());
551  if (pdf != 0) {
552  RooAbsData *data = (RooAbsData *) datasets_[ib]; //dataSets_->FindObject(catClone->getLabel());
553  //RooAbsData *data = (RooAbsData *) dataSets_->FindObject(catClone->getLabel());
554  //std::cout << " bin " << ib << " (label " << catClone->getLabel() << ") has pdf " << pdf->GetName() << " of type " << pdf->ClassName() << " and " << (data ? data->numEntries() : -1) << " dataset entries" << std::endl;
555  if (data == 0) { throw std::logic_error("Error: no data"); }
556  pdfs_[ib] = new CachingAddNLL(catClone->getLabel(), "", pdf, data);
557  params_.add(pdfs_[ib]->params(), /*silent=*/true);
558  } else {
559  pdfs_[ib] = 0;
560  //std::cout << " bin " << ib << " (label " << catClone->getLabel() << ") has no pdf" << std::endl;
561  }
562  }
563 
564  setValueDirty();
565 }
566 
567 Double_t
569 {
570  TRACE_POINT(params_)
571 #ifdef DEBUG_CACHE
572  PerfCounter::add("CachingSimNLL::evaluate called");
573 #endif
574  double ret = 0;
575  for (std::vector<CachingAddNLL*>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
576  if (*it != 0) {
577  double nllval = (*it)->getVal();
578  // what sanity check could I put here?
579  ret += nllval;
580  }
581  }
582  if (!constrainPdfs_.empty()) {
584  std::vector<double>::const_iterator itz = constrainZeroPoints_.begin();
585  for (std::vector<RooAbsPdf *>::const_iterator it = constrainPdfs_.begin(), ed = constrainPdfs_.end(); it != ed; ++it, ++itz) {
586  double pdfval = (*it)->getVal(nuis_);
587  if (!isnormal(pdfval) || pdfval <= 0) {
588  if (!noDeepLEE_) logEvalError((std::string("Constraint pdf ")+(*it)->GetName()+" evaluated to zero, negative or error").c_str());
589  pdfval = 1e-9;
590  }
591  ret -= (log(pdfval) + *itz);
592  }
594  itz = constrainZeroPointsFast_.begin();
595  for (std::vector<SimpleGaussianConstraint*>::const_iterator it = constrainPdfsFast_.begin(), ed = constrainPdfsFast_.end(); it != ed; ++it, ++itz) {
596  double logpdfval = (*it)->getLogValFast();
597  ret -= (logpdfval + *itz);
598  }
599  }
600 #ifdef TRACE_NLL_EVALS
601  static unsigned long _trace_ = 0; _trace_++;
602  if (_trace_ % 10 == 0) { putchar('.'); fflush(stdout); }
603  //if (_trace_ % 250 == 0) { printf(" NLL % 10.4f after %10lu evals.\n", ret, _trace_); fflush(stdout); }
604 #endif
605  TRACE_NLL("SimNLL for " << GetName() << ": " << ret)
606  return ret;
607 }
608 
609 void
611 {
612  dataOriginal_ = &data;
613  //std::cout << "combined data has " << data.numEntries() << " dataset entries (sumw " << data.sumEntries() << ", weighted " << data.isWeighted() << ")" << std::endl;
614  //utils::printRAD(&data);
615  //dataSets_.reset(dataOriginal_->split(pdfOriginal_->indexCat(), true));
616  splitWithWeights(*dataOriginal_, pdfOriginal_->indexCat(), true);
617  for (int ib = 0, nb = pdfs_.size(); ib < nb; ++ib) {
618  CachingAddNLL *canll = pdfs_[ib];
619  if (canll == 0) continue;
620  RooAbsData *data = datasets_[ib];
621  //RooAbsData *data = (RooAbsData *) dataSets_->FindObject(canll->GetName());
622  if (data == 0) { throw std::logic_error("Error: no data"); }
623  //std::cout << " bin " << ib << " (label " << canll->GetName() << ") has pdf " << canll->pdf()->GetName() << " of type " << canll->pdf()->ClassName() <<
624  // " and " << (data ? data->numEntries() : -1) << " dataset entries (sumw " << data->sumEntries() << ", weighted " << data->isWeighted() << ")" << std::endl;
625  canll->setData(*data);
626  }
627 }
628 
629 void cacheutils::CachingSimNLL::splitWithWeights(const RooAbsData &data, const RooAbsCategory& splitCat, Bool_t createEmptyDataSets) {
630  RooCategory *cat = dynamic_cast<RooCategory *>(data.get()->find(splitCat.GetName()));
631  if (cat == 0) throw std::logic_error("Error: no category");
632  int nb = cat->numBins((const char *)0), ne = data.numEntries();
633  RooArgSet obs(*data.get()); obs.remove(*cat, true, true);
634  RooRealVar weight("_weight_","",1);
635  RooArgSet obsplus(obs); obsplus.add(weight);
636  if (nb != int(datasets_.size())) throw std::logic_error("Number of categories changed"); // this can happen due to bugs in RooDataSet
637  for (int ib = 0; ib < nb; ++ib) {
638  if (datasets_[ib] == 0) datasets_[ib] = new RooDataSet("", "", obsplus, "_weight_");
639  else datasets_[ib]->reset();
640  }
641  //utils::printRDH((RooAbsData*)&data);
642  for (int i = 0; i < ne; ++i) {
643  data.get(i); if (data.weight() == 0) continue;
644  int ib = cat->getBin();
645  //std::cout << "Event " << i << " of weight " << data.weight() << " is in bin " << ib << " label " << cat->getLabel() << std::endl;
646  if (data.weight() > 0) datasets_[ib]->add(obs, data.weight());
647  }
648 }
649 
651  for (std::vector<CachingAddNLL*>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
652  if (*it != 0) (*it)->setZeroPoint();
653  }
654  std::vector<double>::iterator itz = constrainZeroPoints_.begin();
655  for (std::vector<RooAbsPdf *>::const_iterator it = constrainPdfs_.begin(), ed = constrainPdfs_.end(); it != ed; ++it, ++itz) {
656  double pdfval = (*it)->getVal(nuis_);
657  if (isnormal(pdfval) || pdfval > 0) *itz = -log(pdfval);
658  }
659  itz = constrainZeroPointsFast_.begin();
660  for (std::vector<SimpleGaussianConstraint*>::const_iterator it = constrainPdfsFast_.begin(), ed = constrainPdfsFast_.end(); it != ed; ++it, ++itz) {
661  double logpdfval = (*it)->getLogValFast();
662  *itz = -logpdfval;
663  }
664  setValueDirty();
665 }
666 
668  for (std::vector<CachingAddNLL*>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
669  if (*it != 0) (*it)->clearZeroPoint();
670  }
671  std::fill(constrainZeroPoints_.begin(), constrainZeroPoints_.end(), 0.0);
672  std::fill(constrainZeroPointsFast_.begin(), constrainZeroPointsFast_.end(), 0.0);
673  setValueDirty();
674 }
675 
676 RooArgSet*
677 cacheutils::CachingSimNLL::getObservables(const RooArgSet* depList, Bool_t valueOnly) const
678 {
679  return new RooArgSet();
680 }
681 
682 RooArgSet*
683 cacheutils::CachingSimNLL::getParameters(const RooArgSet* depList, Bool_t stripDisconnected) const
684 {
685  return new RooArgSet(params_);
686 }
virtual CachingSimNLL * clone(const char *name=0) const
Definition: CachingNLL.cc:483
#define TRACE_NLL(x)
Definition: CachingNLL.cc:59
int i
Definition: DBlmapReader.cc:9
string fill
Definition: lumiContext.py:319
void add(double increment=1.0)
virtual CachingAddNLL * clone(const char *name=0) const
Definition: CachingNLL.cc:266
virtual RooArgSet * getObservables(const RooArgSet *depList, Bool_t valueOnly=kTRUE) const
Definition: CachingNLL.cc:677
int get(const char *name)
void setData(const RooAbsData &data)
Definition: CachingNLL.cc:415
void splitWithWeights(const RooAbsData &data, const RooAbsCategory &splitCat, Bool_t createEmptyDataSets)
Definition: CachingNLL.cc:629
void factorizeFunc(const RooArgSet &observables, RooAbsReal &pdf, RooArgList &obsTerms, RooArgList &otherTerms, bool debug=false)
factorize a RooAbsReal
Definition: utils.cc:177
void add(const std::vector< const T * > &source, std::vector< const T * > &dest)
#define NULL
Definition: scimark2.h:8
const RooAbsData * data_
Definition: CachingNLL.h:93
virtual RooArgSet * getParameters(const RooArgSet *depList, Bool_t stripDisconnected=kTRUE) const
Definition: CachingNLL.cc:458
CachingAddNLL(const char *name, const char *title, RooAbsPdf *pdf, RooAbsData *data)
Definition: CachingNLL.cc:238
virtual Double_t evaluate() const
Definition: CachingNLL.cc:568
Item * items[MaxItems_]
Definition: CachingNLL.h:51
ValuesCache(const RooAbsReal &pdf, const RooArgSet &obs, int size=MaxItems_)
Definition: CachingNLL.cc:101
CachingSimNLL(RooSimultaneous *pdf, RooAbsData *data, const RooArgSet *nuis=0)
Definition: CachingNLL.cc:464
void setData(const RooAbsData &data)
Definition: CachingNLL.cc:610
RooAbsPdf * factorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &constraints)
Definition: utils.cc:80
bool changed(bool updateIfChanged=false)
Definition: CachingNLL.cc:77
std::pair< std::vector< Double_t > *, bool > get()
Definition: CachingNLL.cc:121
double f[11][100]
CachingPdf(RooAbsReal *pdf, const RooArgSet *obs)
Definition: CachingNLL.cc:173
#define TRACE_POINT2(x, i)
Definition: CachingNLL.cc:57
const std::vector< Double_t > & eval(const RooAbsData &data)
Definition: CachingNLL.cc:198
#define TRACE_POINT(x)
Definition: CachingNLL.cc:58
RooAbsReal * fullCloneFunc(const RooAbsReal *pdf, RooArgSet &holder, bool cloneLeafNodes=false)
Definition: utils.cc:226
std::vector< double > vals_
Definition: CachingNLL.h:26
std::vector< RooRealVar * > vars_
Definition: CachingNLL.h:25
void realFill_(const RooAbsData &data, std::vector< Double_t > &values)
Definition: CachingNLL.cc:219
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
virtual RooArgSet * getObservables(const RooArgSet *depList, Bool_t valueOnly=kTRUE) const
Definition: CachingNLL.cc:452
virtual Double_t evaluate() const
Definition: CachingNLL.cc:343
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
size_(0)
Definition: OwnArray.h:181
tuple size
Write out results.
virtual RooArgSet * getParameters(const RooArgSet *depList, Bool_t stripDisconnected=kTRUE) const
Definition: CachingNLL.cc:683
*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
void set(const std::string &name, int value)
set the flag, with a run-time name