CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Protected Attributes

FastVerticalInterpHistPdf2D Class Reference

#include <VerticalInterpHistPdf.h>

Inheritance diagram for FastVerticalInterpHistPdf2D:
FastVerticalInterpHistPdfBase

List of all members.

Public Member Functions

virtual TObject * clone (const char *newname) const
Bool_t conditional () const
Double_t evaluate () const
 FastVerticalInterpHistPdf2D ()
 FastVerticalInterpHistPdf2D (const char *name, const char *title, const RooRealVar &x, const RooRealVar &y, bool conditional, const RooArgList &funcList, const RooArgList &coefList, Double_t smoothRegion=1., Int_t smoothAlgo=1)
 FastVerticalInterpHistPdf2D (const FastVerticalInterpHistPdf2D &other, const char *name=0)
Bool_t hasCache () const
Bool_t isCacheReady () const
Bool_t selfNormalized () const
virtual ~FastVerticalInterpHistPdf2D ()

Protected Member Functions

void setupCaches () const
 not to be serialized
void syncComponent (int which) const
void syncComponents (int dimension) const
void syncNominal () const
void syncTotal () const

Protected Attributes

FastHisto2D _cache
 Cache of the result.
FastHisto2D _cacheNominal
 Cache of nominal pdf (additive morphing) and its bin-by-bin logarithm (multiplicative)
FastHisto2D _cacheNominalLog
 not to be serialized
bool _conditional
RooRealProxy _x
RooRealProxy _y

Detailed Description

Definition at line 167 of file VerticalInterpHistPdf.h.


Constructor & Destructor Documentation

FastVerticalInterpHistPdf2D::FastVerticalInterpHistPdf2D ( ) [inline]

Definition at line 170 of file VerticalInterpHistPdf.h.

Referenced by clone().

FastVerticalInterpHistPdf2D::FastVerticalInterpHistPdf2D ( const char *  name,
const char *  title,
const RooRealVar &  x,
const RooRealVar &  y,
bool  conditional,
const RooArgList &  funcList,
const RooArgList &  coefList,
Double_t  smoothRegion = 1.,
Int_t  smoothAlgo = 1 
) [inline]

If conditional = false, the pdf is normalized integrating on (x,y). If conditional = true, the pdf is separately normalized integrating on (y) for each specific (x) bin

Definition at line 173 of file VerticalInterpHistPdf.h.

                                                                                                                                                                                                                                     :
    FastVerticalInterpHistPdfBase(name,title,RooArgSet(x,y),funcList,coefList,smoothRegion,smoothAlgo),
    _x("x","Independent variable",this,const_cast<RooRealVar&>(x)),
    _y("y","Independent variable",this,const_cast<RooRealVar&>(y)),
    _conditional(conditional),
    _cache(), _cacheNominal(), _cacheNominalLog()  {}
FastVerticalInterpHistPdf2D::FastVerticalInterpHistPdf2D ( const FastVerticalInterpHistPdf2D other,
const char *  name = 0 
) [inline]

Definition at line 179 of file VerticalInterpHistPdf.h.

                                                                                            :
    FastVerticalInterpHistPdfBase(other, name),
    _x("x",this,other._x), _y("y",this,other._y), _conditional(other._conditional),
#ifdef FastVerticalInterpHistPdf_CopyConstructor
    _cache(other._cache), _cacheNominal(other._cacheNominal), _cacheNominalLog(other._cacheNominalLog)  {}
virtual FastVerticalInterpHistPdf2D::~FastVerticalInterpHistPdf2D ( ) [inline, virtual]

Definition at line 188 of file VerticalInterpHistPdf.h.

{}

Member Function Documentation

virtual TObject* FastVerticalInterpHistPdf2D::clone ( const char *  newname) const [inline, virtual]

Implements FastVerticalInterpHistPdfBase.

Definition at line 187 of file VerticalInterpHistPdf.h.

References FastVerticalInterpHistPdf2D().

{ return new FastVerticalInterpHistPdf2D(*this,newname) ; }
Bool_t FastVerticalInterpHistPdf2D::conditional ( ) const [inline]

Definition at line 193 of file VerticalInterpHistPdf.h.

References _conditional.

{ return _conditional; }
Double_t FastVerticalInterpHistPdf2D::evaluate ( ) const [virtual]
Bool_t FastVerticalInterpHistPdf2D::hasCache ( ) const [inline]

Definition at line 195 of file VerticalInterpHistPdf.h.

References _cache, and FastTemplate::size().

{ return _cache.size() > 0; }
Bool_t FastVerticalInterpHistPdf2D::isCacheReady ( ) const [inline]

Definition at line 196 of file VerticalInterpHistPdf.h.

References _cache, FastVerticalInterpHistPdfBase::_init, and FastTemplate::size().

{ return _cache.size() > 0 && _init; }
Bool_t FastVerticalInterpHistPdf2D::selfNormalized ( ) const [inline]

Reimplemented from FastVerticalInterpHistPdfBase.

Definition at line 192 of file VerticalInterpHistPdf.h.

{ return kTRUE; }
void FastVerticalInterpHistPdf2D::setupCaches ( ) const [protected]

not to be serialized

Definition at line 510 of file VerticalInterpHistPdf.cc.

References _cache, _cacheNominal, FastVerticalInterpHistPdfBase::_coefIter, FastVerticalInterpHistPdfBase::_coefList, FastVerticalInterpHistPdfBase::_morphParams, FastVerticalInterpHistPdfBase::_morphs, FastVerticalInterpHistPdfBase::_sentry, SimpleCacheSentry::addVars(), SimpleCacheSentry::empty(), i, FastTemplate::size(), syncComponents(), syncNominal(), syncTotal(), and TRACEME.

Referenced by evaluate().

                                                    {
    TRACEME()
    int ndim = _coefList.getSize();

    _morphs.resize(ndim);
    _morphParams.resize(ndim);
    syncNominal();
    //printf("Nominal template has been set up: \n");  _cacheNominal.Dump();
    _coefIter->Reset();
    for (int i = 0; i < ndim; ++i) {
        _morphParams[i] = dynamic_cast<RooAbsReal *>(_coefIter->Next());
        _morphs[i].sum.Resize(_cacheNominal.size());
        _morphs[i].diff.Resize(_cacheNominal.size());
        syncComponents(i);
    } 
    _cache = FastHisto2D(_cacheNominal);

    if (_sentry.empty()) _sentry.addVars(_coefList); 
    syncTotal();
}
void FastVerticalInterpHistPdf2D::syncComponent ( int  which) const [protected]
void FastVerticalInterpHistPdf2D::syncComponents ( int  dimension) const [protected]

Definition at line 406 of file VerticalInterpHistPdf.cc.

References _cacheNominal, _conditional, FastVerticalInterpHistPdfBase::_funcList, FastVerticalInterpHistPdfBase::_morphs, _x, _y, none, FastHisto2D::Normalize(), FastHisto2D::NormalizeXSlices(), FastVerticalInterpHistPdfBase::syncMorph(), x, and detailsBasic3DVector::y.

Referenced by setupCaches().

                                                              {
    RooAbsPdf *pdfHi = dynamic_cast<RooAbsPdf *>(_funcList.at(2*dim+1));
    RooAbsPdf *pdfLo = dynamic_cast<RooAbsPdf *>(_funcList.at(2*dim+2));
    const RooRealVar &x = dynamic_cast<const RooRealVar &>(_x.arg());
    const RooRealVar &y = dynamic_cast<const RooRealVar &>(_y.arg());
    const RooCmdArg &cond = _conditional ? RooFit::ConditionalObservables(RooArgSet(x)) : RooCmdArg::none();
    std::auto_ptr<TH1> histHi(pdfHi->createHistogram("", x, RooFit::YVar(y), cond)); histHi->SetDirectory(0); 
    std::auto_ptr<TH1> histLo(pdfLo->createHistogram("", x, RooFit::YVar(y), cond)); histLo->SetDirectory(0);

    FastHisto2D hi(dynamic_cast<TH2&>(*histHi), _conditional), lo(dynamic_cast<TH2&>(*histLo), _conditional); 
    //printf("Un-normalized templates for dimension %d: \n", dim);  hi.Dump(); lo.Dump();
    if (_conditional) {
        hi.NormalizeXSlices(); lo.NormalizeXSlices();
    } else {
        hi.Normalize(); lo.Normalize();
    }
    //printf("Normalized templates for dimension %d: \n", dim);  hi.Dump(); lo.Dump();
    syncMorph(_morphs[dim], _cacheNominal, lo, hi);
}
void FastVerticalInterpHistPdf2D::syncNominal ( ) const [protected]

Definition at line 360 of file VerticalInterpHistPdf.cc.

References _cacheNominal, _cacheNominalLog, _conditional, FastVerticalInterpHistPdfBase::_funcList, FastVerticalInterpHistPdfBase::_smoothAlgo, _x, _y, estimatePileup::hist, FastTemplate::Log(), none, FastHisto2D::Normalize(), FastHisto2D::NormalizeXSlices(), TRACEME, x, and detailsBasic3DVector::y.

Referenced by setupCaches().

                                                    {
    TRACEME()
    RooAbsPdf *pdf = dynamic_cast<RooAbsPdf *>(_funcList.at(0));
    const RooRealVar &x = dynamic_cast<const RooRealVar &>(_x.arg());
    const RooRealVar &y = dynamic_cast<const RooRealVar &>(_y.arg());
    const RooCmdArg &cond = _conditional ? RooFit::ConditionalObservables(RooArgSet(x)) : RooCmdArg::none();
    std::auto_ptr<TH1> hist(pdf->createHistogram("", x, RooFit::YVar(y), cond));
    hist->SetDirectory(0); 
    _cacheNominal = FastHisto2D(dynamic_cast<TH2F&>(*hist), _conditional);
    if (_conditional) _cacheNominal.NormalizeXSlices(); 
    else              _cacheNominal.Normalize(); 

    if (_smoothAlgo < 0) {
        _cacheNominalLog = _cacheNominal;
        _cacheNominalLog.Log();
    }
}
void FastVerticalInterpHistPdf2D::syncTotal ( ) const [protected]

Definition at line 479 of file VerticalInterpHistPdf.cc.

References _cache, _cacheNominal, _cacheNominalLog, _conditional, FastHisto2D::Normalize(), and FastHisto2D::NormalizeXSlices().

Referenced by evaluate(), and setupCaches().

                                                  {
    FastVerticalInterpHistPdfBase::syncTotal(_cache, _cacheNominal, _cacheNominalLog);
    // normalize the result
    if (_conditional) _cache.NormalizeXSlices(); 
    else              _cache.Normalize(); 
    //printf("Normalized result\n");  _cache.Dump();
}

Member Data Documentation

Cache of the result.

Definition at line 202 of file VerticalInterpHistPdf.h.

Referenced by evaluate(), hasCache(), isCacheReady(), setupCaches(), and syncTotal().

Cache of nominal pdf (additive morphing) and its bin-by-bin logarithm (multiplicative)

not to be serialized

Definition at line 204 of file VerticalInterpHistPdf.h.

Referenced by setupCaches(), syncComponents(), syncNominal(), and syncTotal().

not to be serialized

Definition at line 205 of file VerticalInterpHistPdf.h.

Referenced by syncNominal(), and syncTotal().

Definition at line 199 of file VerticalInterpHistPdf.h.

Referenced by conditional(), syncComponents(), syncNominal(), and syncTotal().

RooRealProxy FastVerticalInterpHistPdf2D::_x [protected]

Reimplemented from FastVerticalInterpHistPdfBase.

Definition at line 198 of file VerticalInterpHistPdf.h.

Referenced by evaluate(), syncComponents(), and syncNominal().

RooRealProxy FastVerticalInterpHistPdf2D::_y [protected]

Definition at line 198 of file VerticalInterpHistPdf.h.

Referenced by evaluate(), syncComponents(), and syncNominal().