1 #include "../interface/CachingNLL.h"
2 #include "../interface/utils.h"
4 #include <RooCategory.h>
5 #include <RooDataSet.h>
6 #include <RooProduct.h>
7 #include "../interface/ProfilingTools.h"
18 #include "../interface/ProfilingTools.h"
25 #ifdef DEBUG_TRACE_POINTS
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());
42 RooLinkedListIter iter = point.iterator();
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();
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;
57 #define TRACE_POINT2(x,i)
58 #define TRACE_POINT(x)
64 std::auto_ptr<TIterator> iter(set.createIterator());
65 for (RooAbsArg *
a = dynamic_cast<RooAbsArg *>(iter->Next());
67 a =
dynamic_cast<RooAbsArg *
>(iter->Next())) {
68 RooRealVar *rrv =
dynamic_cast<RooRealVar *
>(
a);
71 vals_.push_back(rrv->getVal());
79 std::vector<RooRealVar *>::const_iterator it = vars_.begin(), ed = vars_.end();
80 std::vector<double>::iterator itv = vals_.begin();
82 for ( ; it != ed; ++it, ++itv) {
83 double val = (*it)->getVal();
87 if (updateIfChanged) { *itv = val; }
106 std::auto_ptr<RooArgSet> params(pdf.getParameters(obs));
113 for (
int i = 0; i <
size_; ++
i)
delete items[i];
118 for (
int i = 0; i <
size_; ++
i) items[i]->good =
false;
123 int found = -1;
bool good =
false;
124 for (
int i = 0; i <
size_; ++
i) {
125 if (items[i]->good) {
127 if (!items[i]->checker.changed()) {
129 PerfCounter::add(i == 0 ?
"ValuesCache::get hit first" :
"ValuesCache::get hit other");
136 }
else if (found == -1) {
149 if (size_ < maxSize_) {
151 items[
size_] =
new Item(items[0]->checker);
164 while (found > 0) { items[
found] = items[found-1]; --
found; }
168 if (!good) items[
found]->checker.changed(
true);
169 items[
found]->good =
true;
170 return std::pair<std::vector<Double_t> *,
bool>(&items[
found]->values, good);
185 pdfOriginal_(other.pdfOriginal_),
197 const std::vector<Double_t> &
203 bool newdata = (lastData_ != &
data);
206 pdf_->optimizeCacheMode(*data.get());
207 pdf_->attachDataSet(data);
208 const_cast<RooAbsData*
>(lastData_)->setDirtyProp(
false);
211 std::pair<std::vector<Double_t> *,
bool>
hit = cache_.get();
213 realFill_(data, *hit.first);
224 int n = data.numEntries();
227 std::vector<Double_t>::iterator itv = vals.begin();
228 for (
int i = 0; i <
n; ++
i, ++itv) {
231 *itv = pdf_->getVal(obs_);
233 TRACE_NLL(
"PDF value for " << pdf_->GetName() <<
" is " << *itv <<
" at this point.")
239 RooAbsReal(name, title),
241 params_(
"params",
"parameters",this),
244 if (pdf == 0)
throw std::invalid_argument(std::string(
"Pdf passed to ")+name+
" is null");
250 RooAbsReal(name ? name : (TString(
"nll_")+other.pdf_->GetName()).
Data(),
""),
252 params_(
"params",
"parameters",this),
261 for (
int i = 0,
n = integrals_.size(); i <
n; ++
i)
delete integrals_[i];
274 const RooArgSet *obs = data_->get();
275 for (
int i = 0,
n = integrals_.size(); i <
n; ++
i)
delete integrals_[i];
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);
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);
290 }
else if ((sumpdf = dynamic_cast<RooRealSumPdf *>(pdf_)) != 0) {
291 isRooRealSum_ =
true;
292 int npdf = sumpdf->coefList().getSize();
293 coeffs_.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;
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));
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?");
321 coeffs_.push_back(coeff);
323 integrals_.push_back(funci->createIntegral(*obs));
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);
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);
338 if (rrv != 0) params_.add(*rrv);
348 std::fill( partialSum_.begin(), partialSum_.end(), 0.0 );
350 std::vector<RooAbsReal*>::iterator itc = coeffs_.begin(), edc = coeffs_.end();
351 std::vector<CachingPdf>::iterator itp = pdfs_.begin();
352 std::vector<Double_t>::const_iterator itw, bgw = weights_.begin();
353 std::vector<Double_t>::iterator its, bgs = partialSum_.begin(), eds = partialSum_.end();
356 for ( ; itc != edc; ++itp, ++itc ) {
358 Double_t coeff = (*itc)->getVal();
360 sumCoeff += coeff * integrals_[itc - coeffs_.begin()]->getVal();
366 const std::vector<Double_t> &pdfvals = itp->eval(*data_);
368 std::vector<Double_t>::const_iterator itv = pdfvals.begin();
369 for (its = bgs; its != eds; ++its, ++itv) {
370 *its += coeff * (*itv);
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;
381 double thispiece = (*itw) * (*its <= 0 ? -9e9 :
log( ((*its) / sumCoeff) ));
382 #ifdef ADDNLL_KAHAN_SUM
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;
401 double expectedEvents = (isRooRealSum_ ? pdf_->getNorm(data_->get()) : sumCoeff);
402 if (expectedEvents <= 0) {
404 expectedEvents = 1
e-6;
407 ret += expectedEvents - sumWeights_ *
log(expectedEvents);
410 TRACE_NLL(
"AddNLL for " << pdf_->GetName() <<
": " <<
ret)
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;
428 for (
int i = 0,
n = data.numEntries(); i <
n; ++
i, ++itw) {
430 *itw = data.weight();
431 #ifdef ADDNLL_KAHAN_SUM
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;
446 for (std::vector<CachingPdf>::iterator itp = pdfs_.begin(), edp = pdfs_.end(); itp != edp; ++itp) {
454 return new RooArgSet();
460 return new RooArgSet(params_);
468 params_(
"params",
"parameters",this)
474 pdfOriginal_(other.pdfOriginal_),
475 dataOriginal_(other.dataOriginal_),
477 params_(
"params",
"parameters",this)
490 for (std::vector<SimpleGaussianConstraint*>::iterator it = constrainPdfsFast_.begin(), ed = constrainPdfsFast_.end(); it != ed; ++it) {
502 RooAbsPdf *pdfclone = pdfOriginal_;
512 RooSimultaneous *simpdf = factorizedPdf_.get();
513 constrainPdfs_.clear();
514 if (constraints.getSize()) {
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)) {
521 constrainZeroPointsFast_.push_back(0);
523 constrainPdfs_.push_back(pdfi);
524 constrainZeroPoints_.push_back(0);
527 std::auto_ptr<RooArgSet> params(pdfi->getParameters(*dataOriginal_));
528 params_.add(*params,
false);
531 std::cerr <<
"PDF didn't factorize!" << std::endl;
532 std::cout <<
"Parameters: " << std::endl;
533 std::auto_ptr<RooArgSet> params(pdfclone->getParameters(*dataOriginal_));
536 dataOriginal_->get()->Print(
"V");
537 factorizedPdf_.release();
538 simpdf =
dynamic_cast<RooSimultaneous *
>(pdfclone);
542 std::auto_ptr<RooAbsCategoryLValue> catClone((RooAbsCategoryLValue*) simpdf->indexCat().Clone());
543 pdfs_.resize(catClone->numBins(
NULL), 0);
545 datasets_.resize(pdfs_.size(), 0);
546 splitWithWeights(*dataOriginal_, simpdf->indexCat(),
true);
548 for (
int ib = 0, nb = pdfs_.size(); ib < nb; ++ib) {
549 catClone->setBin(ib);
550 RooAbsPdf *pdf = simpdf->getPdf(catClone->getLabel());
552 RooAbsData *
data = (RooAbsData *) datasets_[ib];
555 if (data == 0) {
throw std::logic_error(
"Error: no data"); }
557 params_.add(pdfs_[ib]->params(),
true);
575 for (std::vector<CachingAddNLL*>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
577 double nllval = (*it)->getVal();
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());
591 ret -= (
log(pdfval) + *itz);
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);
600 #ifdef TRACE_NLL_EVALS
601 static unsigned long _trace_ = 0; _trace_++;
602 if (_trace_ % 10 == 0) { putchar(
'.'); fflush(stdout); }
605 TRACE_NLL(
"SimNLL for " << GetName() <<
": " << ret)
612 dataOriginal_ = &
data;
616 splitWithWeights(*dataOriginal_, pdfOriginal_->indexCat(),
true);
617 for (
int ib = 0, nb = pdfs_.size(); ib < nb; ++ib) {
619 if (canll == 0)
continue;
620 RooAbsData *data = datasets_[ib];
622 if (data == 0) {
throw std::logic_error(
"Error: no data"); }
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");
637 for (
int ib = 0; ib < nb; ++ib) {
638 if (datasets_[ib] == 0) datasets_[ib] =
new RooDataSet(
"",
"", obsplus,
"_weight_");
639 else datasets_[ib]->reset();
642 for (
int i = 0; i <
ne; ++
i) {
643 data.get(i);
if (data.weight() == 0)
continue;
644 int ib = cat->getBin();
646 if (data.weight() > 0) datasets_[ib]->
add(obs, data.weight());
651 for (std::vector<CachingAddNLL*>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
652 if (*it != 0) (*it)->setZeroPoint();
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);
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();
668 for (std::vector<CachingAddNLL*>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
669 if (*it != 0) (*it)->clearZeroPoint();
671 std::fill(constrainZeroPoints_.begin(), constrainZeroPoints_.end(), 0.0);
672 std::fill(constrainZeroPointsFast_.begin(), constrainZeroPointsFast_.end(), 0.0);
679 return new RooArgSet();
685 return new RooArgSet(params_);
virtual CachingSimNLL * clone(const char *name=0) const
void add(double increment=1.0)
virtual CachingAddNLL * clone(const char *name=0) const
virtual RooArgSet * getObservables(const RooArgSet *depList, Bool_t valueOnly=kTRUE) const
int get(const char *name)
void setData(const RooAbsData &data)
void splitWithWeights(const RooAbsData &data, const RooAbsCategory &splitCat, Bool_t createEmptyDataSets)
void factorizeFunc(const RooArgSet &observables, RooAbsReal &pdf, RooArgList &obsTerms, RooArgList &otherTerms, bool debug=false)
factorize a RooAbsReal
void add(const std::vector< const T * > &source, std::vector< const T * > &dest)
virtual RooArgSet * getParameters(const RooArgSet *depList, Bool_t stripDisconnected=kTRUE) const
CachingAddNLL(const char *name, const char *title, RooAbsPdf *pdf, RooAbsData *data)
virtual Double_t evaluate() const
ValuesCache(const RooAbsReal &pdf, const RooArgSet &obs, int size=MaxItems_)
CachingSimNLL(RooSimultaneous *pdf, RooAbsData *data, const RooArgSet *nuis=0)
void setData(const RooAbsData &data)
RooAbsPdf * factorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &constraints)
bool changed(bool updateIfChanged=false)
std::pair< std::vector< Double_t > *, bool > get()
CachingPdf(RooAbsReal *pdf, const RooArgSet *obs)
#define TRACE_POINT2(x, i)
const std::vector< Double_t > & eval(const RooAbsData &data)
RooAbsReal * fullCloneFunc(const RooAbsReal *pdf, RooArgSet &holder, bool cloneLeafNodes=false)
std::vector< double > vals_
std::vector< RooRealVar * > vars_
void realFill_(const RooAbsData &data, std::vector< Double_t > &values)
char data[epos_bytes_allocation]
virtual RooArgSet * getObservables(const RooArgSet *depList, Bool_t valueOnly=kTRUE) const
virtual Double_t evaluate() const
tuple size
Write out results.
virtual RooArgSet * getParameters(const RooArgSet *depList, Bool_t stripDisconnected=kTRUE) const
*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
void set(const std::string &name, int value)
set the flag, with a run-time name