Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorph.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * This code was autogenerated by RooClassFactory *
5 *****************************************************************************/
6
7#ifndef ROOMOMENTMORPH
8#define ROOMOMENTMORPH
9
10#include "RooAbsPdf.h"
11#include "RooRealProxy.h"
12#include "RooCategoryProxy.h"
13#include "RooAbsReal.h"
14#include "RooAbsCategory.h"
15#include "RooSetProxy.h"
16#include "RooListProxy.h"
17#include "RooArgList.h"
18
19#include "TMatrixD.h"
20#include "TVectorD.h"
21
23
24class RooMomentMorph : public RooAbsPdf {
25public:
26
28
30
31 RooMomentMorph(const char *name, const char *title, RooAbsReal& _m, const RooArgList& varList,
32 const RooArgList& pdfList, const RooArgList& mrefList, Setting setting = NonLinearPosFractions);
33 RooMomentMorph(const char *name, const char *title, RooAbsReal& _m, const RooArgList& varList,
34 const RooArgList& pdfList, const TVectorD& mrefpoints, Setting setting = NonLinearPosFractions );
35 RooMomentMorph(const RooMomentMorph& other, const char* name=nullptr) ;
36 TObject* clone(const char* newname) const override { return new RooMomentMorph(*this,newname); }
37 ~RooMomentMorph() override;
38
39 void setMode(const Setting& setting) { _setting = setting; }
40
41 void useHorizontalMorphing(bool val) { _useHorizMorph = val; }
42
43 bool selfNormalized() const override {
44 // P.d.f is self normalized
45 return true ;
46 }
47
48 virtual double getVal(const RooArgSet* set=nullptr) const ;
49 RooAbsPdf* sumPdf(const RooArgSet* nset) ;
50
51
52protected:
53
55 public:
56 CacheElem(std::unique_ptr<RooAbsPdf> && sumPdf,
57 std::unique_ptr<RooChangeTracker> && tracker,
58 const RooArgList& flist);
59 ~CacheElem() override ;
61 std::unique_ptr<RooAbsPdf> _sumPdf ;
62 std::unique_ptr<RooChangeTracker> _tracker ;
64
65 RooRealVar* frac(Int_t i ) ;
66 const RooRealVar* frac(Int_t i ) const ;
67 void calculateFractions(const RooMomentMorph& self, bool verbose=true) const;
68 } ;
69 mutable RooObjCacheManager _cacheMgr ; //! The cache manager
70 mutable RooArgSet* _curNormSet ; //! Current normalization set
71
72 friend class CacheElem ; // Cache needs to be able to clear _norm pointer
73
74 double evaluate() const override ;
75
76 void initialize();
77 CacheElem* getCache(const RooArgSet* nset) const ;
78
79 inline Int_t ij(const Int_t& i, const Int_t& j) const { return (i*_varList.getSize()+j); }
80 int idxmin(const double& m) const;
81 int idxmax(const double& m) const;
82
86 mutable TVectorD* _mref;
87
88 mutable TMatrixD* _M; //
89
91
93
94 ClassDefOverride(RooMomentMorph,3) // Your description goes here...
95};
96
97#endif
98
99
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
char name[80]
Definition TGX11.cxx:110
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
Int_t getSize() const
Return the number of elements in the collection.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Meta object that tracks value changes in a given set of RooAbsArgs by registering itself as value cli...
RooArgList containedArgs(Action) override
RooRealVar * frac(Int_t i)
void calculateFractions(const RooMomentMorph &self, bool verbose=true) const
std::unique_ptr< RooAbsPdf > _sumPdf
std::unique_ptr< RooChangeTracker > _tracker
RooObjCacheManager _cacheMgr
virtual double getVal(const RooArgSet *set=nullptr) const
Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set...
~RooMomentMorph() override
RooSetProxy _varList
void setMode(const Setting &setting)
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooMomentMorph()
coverity[UNINIT_CTOR]
RooRealProxy m
RooAbsPdf * sumPdf(const RooArgSet *nset)
Int_t ij(const Int_t &i, const Int_t &j) const
RooArgSet * _curNormSet
The cache manager.
bool selfNormalized() const override
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
RooListProxy _pdfList
TObject * clone(const char *newname) const override
int idxmin(const double &m) const
int idxmax(const double &m) const
void useHorizontalMorphing(bool val)
CacheElem * getCache(const RooArgSet *nset) const
TVectorD * _mref
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41