CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

cacheutils::CachingAddNLL Class Reference

#include <CachingNLL.h>

List of all members.

Public Member Functions

 CachingAddNLL (const char *name, const char *title, RooAbsPdf *pdf, RooAbsData *data)
 CachingAddNLL (const CachingAddNLL &other, const char *name=0)
void clearZeroPoint ()
virtual CachingAddNLLclone (const char *name=0) const
virtual Double_t defaultErrorLevel () const
virtual Double_t evaluate () const
virtual RooArgSet * getObservables (const RooArgSet *depList, Bool_t valueOnly=kTRUE) const
virtual RooArgSet * getParameters (const RooArgSet *depList, Bool_t stripDisconnected=kTRUE) const
virtual Bool_t isDerived () const
RooSetProxy & params ()
const RooAbsPdf * pdf () const
void setData (const RooAbsData &data)
void setZeroPoint ()
double sumWeights () const
virtual ~CachingAddNLL ()

Private Member Functions

void setup_ ()

Private Attributes

std::vector< RooAbsReal * > coeffs_
const RooAbsData * data_
std::vector< RooAbsReal * > integrals_
bool isRooRealSum_
RooSetProxy params_
std::vector< Double_t > partialSum_
RooAbsPdf * pdf_
std::vector< CachingPdfpdfs_
double sumWeights_
std::vector< Double_t > weights_
double zeroPoint_

Detailed Description

Definition at line 72 of file CachingNLL.h.


Constructor & Destructor Documentation

cacheutils::CachingAddNLL::CachingAddNLL ( const char *  name,
const char *  title,
RooAbsPdf *  pdf,
RooAbsData *  data 
)

Definition at line 238 of file CachingNLL.cc.

References setData(), and setup_().

                                                                                                            :
    RooAbsReal(name, title),
    pdf_(pdf),
    params_("params","parameters",this),
    zeroPoint_(0)
{
    if (pdf == 0) throw std::invalid_argument(std::string("Pdf passed to ")+name+" is null");
    setData(*data);
    setup_();
}
cacheutils::CachingAddNLL::CachingAddNLL ( const CachingAddNLL other,
const char *  name = 0 
)

Definition at line 249 of file CachingNLL.cc.

References data_, setData(), and setup_().

                                                                                   :
    RooAbsReal(name ? name : (TString("nll_")+other.pdf_->GetName()).Data(), ""),
    pdf_(other.pdf_),
    params_("params","parameters",this),
    zeroPoint_(0)
{
    setData(*other.data_);
    setup_();
}
cacheutils::CachingAddNLL::~CachingAddNLL ( ) [virtual]

Definition at line 259 of file CachingNLL.cc.

References i, and n.

{
    for (int i = 0, n = integrals_.size(); i < n; ++i) delete integrals_[i];
    integrals_.clear();
}

Member Function Documentation

void cacheutils::CachingAddNLL::clearZeroPoint ( ) [inline]

Definition at line 87 of file CachingNLL.h.

References zeroPoint_.

{ zeroPoint_ = 0.0; setValueDirty();  }
cacheutils::CachingAddNLL * cacheutils::CachingAddNLL::clone ( const char *  name = 0) const [virtual]

Definition at line 266 of file CachingNLL.cc.

{
    return new cacheutils::CachingAddNLL(*this, name);
}
virtual Double_t cacheutils::CachingAddNLL::defaultErrorLevel ( ) const [inline, virtual]

Definition at line 80 of file CachingNLL.h.

{ return 0.5; }
Double_t cacheutils::CachingAddNLL::evaluate ( ) const [virtual]

Definition at line 343 of file CachingNLL.cc.

References Clusterizer1DCommons::add(), dtNoiseDBValidation_cfg::cerr, alignCSCRings::e, lumiContext::fill, reco::get(), cacheutils::CachingSimNLL::hasError_, create_public_lumi_plots::log, cacheutils::CachingSimNLL::noDeepLEE_, run_regression::ret, and TRACE_NLL.

{
#ifdef DEBUG_CACHE
    PerfCounter::add("CachingAddNLL::evaluate called");
#endif
    std::fill( partialSum_.begin(), partialSum_.end(), 0.0 );

    std::vector<RooAbsReal*>::iterator  itc = coeffs_.begin(), edc = coeffs_.end();
    std::vector<CachingPdf>::iterator   itp = pdfs_.begin();//,   edp = pdfs_.end();
    std::vector<Double_t>::const_iterator itw, bgw = weights_.begin();//,    edw = weights_.end();
    std::vector<Double_t>::iterator       its, bgs = partialSum_.begin(), eds = partialSum_.end();
    double sumCoeff = 0;
    //std::cout << "Performing evaluation of " << GetName() << std::endl;
    for ( ; itc != edc; ++itp, ++itc ) {
        // get coefficient
        Double_t coeff = (*itc)->getVal();
        if (isRooRealSum_) {
            sumCoeff += coeff * integrals_[itc - coeffs_.begin()]->getVal();
            //std::cout << "  coefficient = " << coeff << ", integral = " << integrals_[itc - coeffs_.begin()]->getVal() << std::endl;
        } else {
            sumCoeff += coeff;
        }
        // get vals
        const std::vector<Double_t> &pdfvals = itp->eval(*data_);
        // update running sum
        std::vector<Double_t>::const_iterator itv = pdfvals.begin();
        for (its = bgs; its != eds; ++its, ++itv) {
             *its += coeff * (*itv); // sum (n_i * pdf_i)
        }
    }
    // then get the final nll
    double ret = 0;
    for ( its = bgs, itw = bgw ; its != eds ; ++its, ++itw ) {
        if (*itw == 0) continue;
        if (!isnormal(*its) || *its <= 0) {
            std::cerr << "WARNING: underflow to " << *its << " in " << GetName() << std::endl; 
            if (!CachingSimNLL::noDeepLEE_) logEvalError("Number of events is negative or error"); else CachingSimNLL::hasError_ = true;
        }
        double thispiece = (*itw) * (*its <= 0 ? -9e9 : log( ((*its) / sumCoeff) ));
        #ifdef ADDNLL_KAHAN_SUM
        static bool do_kahan = runtimedef::get("ADDNLL_KAHAN_SUM");
        if (do_kahan) {
            double kahan_y = thispiece  - compensation;
            double kahan_t = ret + kahan_y;
            double kahan_d = (kahan_t - ret);
            compensation = kahan_d - kahan_y;
            ret  = kahan_t;
        } else {
            ret += thispiece;
        }
        #else
        ret += thispiece;
        #endif
    }
    // then flip sign
    ret = -ret;
    // std::cout << "AddNLL for " << pdf_->GetName() << ": " << ret << std::endl;
    // and add extended term: expected - observed*log(expected);
    double expectedEvents = (isRooRealSum_ ? pdf_->getNorm(data_->get()) : sumCoeff);
    if (expectedEvents <= 0) {
        if (!CachingSimNLL::noDeepLEE_) logEvalError("Expected number of events is negative"); else CachingSimNLL::hasError_ = true;
        expectedEvents = 1e-6;
    }
    //ret += expectedEvents - UInt_t(sumWeights_) * log(expectedEvents); // no, doesn't work with Asimov dataset
    ret += expectedEvents - sumWeights_ * log(expectedEvents);
    ret += zeroPoint_;
    // std::cout << "     plus extended term: " << ret << std::endl;
    TRACE_NLL("AddNLL for " << pdf_->GetName() << ": " << ret)
    return ret;
}
RooArgSet * cacheutils::CachingAddNLL::getObservables ( const RooArgSet *  depList,
Bool_t  valueOnly = kTRUE 
) const [virtual]

Definition at line 452 of file CachingNLL.cc.

{
    return new RooArgSet();
}
RooArgSet * cacheutils::CachingAddNLL::getParameters ( const RooArgSet *  depList,
Bool_t  stripDisconnected = kTRUE 
) const [virtual]

Definition at line 458 of file CachingNLL.cc.

{
    return new RooArgSet(params_); 
}
virtual Bool_t cacheutils::CachingAddNLL::isDerived ( ) const [inline, virtual]

Definition at line 79 of file CachingNLL.h.

{ return kTRUE; }
RooSetProxy& cacheutils::CachingAddNLL::params ( ) [inline]

Definition at line 88 of file CachingNLL.h.

References params_.

{ return params_; }
const RooAbsPdf* cacheutils::CachingAddNLL::pdf ( ) const [inline]

Definition at line 85 of file CachingNLL.h.

References pdf_.

{ return pdf_; }
void cacheutils::CachingAddNLL::setData ( const RooAbsData &  data)

Definition at line 415 of file CachingNLL.cc.

References data, reco::get(), i, and n.

Referenced by CachingAddNLL(), and cacheutils::CachingSimNLL::setData().

{
    //std::cout << "Setting data for pdf " << pdf_->GetName() << std::endl;
    //utils::printRAD(&data);
    data_ = &data;
    setValueDirty();
    sumWeights_ = 0.0;
    weights_.resize(data.numEntries());
    partialSum_.resize(data.numEntries());
    std::vector<Double_t>::iterator itw = weights_.begin();
    #ifdef ADDNLL_KAHAN_SUM
    double compensation = 0;
    #endif
    for (int i = 0, n = data.numEntries(); i < n; ++i, ++itw) {
        data.get(i);
        *itw = data.weight();
        #ifdef ADDNLL_KAHAN_SUM
        static bool do_kahan = runtimedef::get("ADDNLL_KAHAN_SUM");
        if (do_kahan) {
            double kahan_y = *itw - compensation;
            double kahan_t = sumWeights_ + kahan_y;
            double kahan_d = (kahan_t - sumWeights_);
            compensation = kahan_d - kahan_y;
            sumWeights_  = kahan_t;
        } else {
            sumWeights_ += *itw;
        }
        #else
        sumWeights_ += *itw;
        #endif
    }
    for (std::vector<CachingPdf>::iterator itp = pdfs_.begin(), edp = pdfs_.end(); itp != edp; ++itp) {
        itp->setDataDirty();
    }
}
void cacheutils::CachingAddNLL::setup_ ( ) [private]

Temporarily switch this off, it doesn't work. Don't know why, however.

Definition at line 272 of file CachingNLL.cc.

References a, gather_cfg::cout, utils::factorizeFunc(), i, and n.

Referenced by CachingAddNLL().

{
    const RooArgSet *obs = data_->get();
    for (int i = 0, n = integrals_.size(); i < n; ++i) delete integrals_[i];
    integrals_.clear();
    RooAddPdf *addpdf = 0;
    RooRealSumPdf *sumpdf = 0;
    if ((addpdf = dynamic_cast<RooAddPdf *>(pdf_)) != 0) {
        isRooRealSum_ = false;
        int npdf = addpdf->coefList().getSize();
        coeffs_.reserve(npdf);
        pdfs_.reserve(npdf);
        for (int i = 0; i < npdf; ++i) {
            RooAbsReal * coeff = dynamic_cast<RooAbsReal*>(addpdf->coefList().at(i));
            RooAbsPdf  * pdfi  = dynamic_cast<RooAbsPdf *>(addpdf->pdfList().at(i));
            coeffs_.push_back(coeff);
            pdfs_.push_back(CachingPdf(pdfi, obs));
        }
    } else if ((sumpdf = dynamic_cast<RooRealSumPdf *>(pdf_)) != 0) {
        isRooRealSum_ = true;
        int npdf = sumpdf->coefList().getSize();
        coeffs_.reserve(npdf);
        pdfs_.reserve(npdf);
        integrals_.reserve(npdf);
        for (int i = 0; i < npdf; ++i) {
            RooAbsReal * coeff = dynamic_cast<RooAbsReal*>(sumpdf->coefList().at(i));
            RooAbsReal * funci = dynamic_cast<RooAbsReal*>(sumpdf->funcList().at(i));
            if (0 && typeid(*funci) == typeid(RooProduct)) {
                RooArgList obsDep, obsInd;
                obsInd.add(*coeff);
                utils::factorizeFunc(*obs, *funci, obsDep, obsInd);
                std::cout << "Entry " << i << ": coef name " << (coeff ? coeff->GetName()   : "null") << 
                                              "  type " << (coeff ? coeff->ClassName() :  "n/a") << std::endl;
                std::cout << "       " <<     "; func name " << (funci ? funci->GetName()   : "null") << 
                                              "  type " << (funci ? funci->ClassName() :  "n/a") << std::endl;
                std::cout << "Terms depending on observables: " << std::endl; obsDep.Print("V");
                std::cout << "Terms not depending on observables: " << std::endl; obsInd.Print("V");
                if (obsInd.getSize() > 1) {
                    coeff = new RooProduct(TString::Format("%s_x_%s_obsIndep", coeff->GetName(), funci->GetName()), "", RooArgSet(obsInd));
                    addOwnedComponents(RooArgSet(*coeff));
                }
                if (obsDep.getSize() > 1) {
                    funci = new RooProduct(TString::Format("%s_obsDep", funci->GetName()), "", RooArgSet(obsInd));
                    addOwnedComponents(RooArgSet(*funci));
                } else if (obsDep.getSize() == 1) {
                    funci = (RooAbsReal *) obsDep.first();
                } else throw std::logic_error("No part of pdf depends on observables?");
            }
            coeffs_.push_back(coeff);
            pdfs_.push_back(CachingPdf(funci, obs));
            integrals_.push_back(funci->createIntegral(*obs));
        }
    } else {
        std::string errmsg = "ERROR: CachingAddNLL: Pdf ";
        errmsg += pdf_->GetName();
        errmsg += " is neither a RooAddPdf nor a RooRealSumPdf, but a ";
        errmsg += pdf_->ClassName();
        throw std::invalid_argument(errmsg);
    }

    std::auto_ptr<RooArgSet> params(pdf_->getParameters(*data_));
    std::auto_ptr<TIterator> iter(params->createIterator());
    for (RooAbsArg *a = (RooAbsArg *) iter->Next(); a != 0; a = (RooAbsArg *) iter->Next()) {
        RooRealVar *rrv = dynamic_cast<RooRealVar *>(a);
        //if (rrv != 0 && !rrv->isConstant()) params_.add(*rrv);
        if (rrv != 0) params_.add(*rrv);
    }
}
void cacheutils::CachingAddNLL::setZeroPoint ( ) [inline]

Definition at line 86 of file CachingNLL.h.

References zeroPoint_.

{ zeroPoint_ = -this->getVal(); setValueDirty(); }
double cacheutils::CachingAddNLL::sumWeights ( ) const [inline]

Definition at line 84 of file CachingNLL.h.

References sumWeights_.

{ return sumWeights_; }

Member Data Documentation

std::vector<RooAbsReal*> cacheutils::CachingAddNLL::coeffs_ [mutable, private]

Definition at line 96 of file CachingNLL.h.

const RooAbsData* cacheutils::CachingAddNLL::data_ [private]

Definition at line 93 of file CachingNLL.h.

Referenced by CachingAddNLL().

std::vector<RooAbsReal*> cacheutils::CachingAddNLL::integrals_ [mutable, private]

Definition at line 98 of file CachingNLL.h.

Definition at line 100 of file CachingNLL.h.

RooSetProxy cacheutils::CachingAddNLL::params_ [private]

Definition at line 92 of file CachingNLL.h.

Referenced by params().

std::vector<Double_t> cacheutils::CachingAddNLL::partialSum_ [mutable, private]

Definition at line 99 of file CachingNLL.h.

RooAbsPdf* cacheutils::CachingAddNLL::pdf_ [private]

Definition at line 91 of file CachingNLL.h.

Referenced by pdf().

std::vector<CachingPdf> cacheutils::CachingAddNLL::pdfs_ [mutable, private]

Definition at line 97 of file CachingNLL.h.

Definition at line 95 of file CachingNLL.h.

Referenced by sumWeights().

std::vector<Double_t> cacheutils::CachingAddNLL::weights_ [private]

Definition at line 94 of file CachingNLL.h.

Definition at line 101 of file CachingNLL.h.

Referenced by clearZeroPoint(), and setZeroPoint().