Logo ROOT   6.08/07
Reference Guide
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 
22 #include <vector>
23 #include <string>
24 class RooChangeTracker ;
25 
26 class RooMomentMorph : public RooAbsPdf {
27 public:
28 
30 
31  RooMomentMorph() ;
32 
33  RooMomentMorph(const char *name, const char *title, RooAbsReal& _m, const RooArgList& varList,
34  const RooArgList& pdfList, const RooArgList& mrefList, Setting setting = NonLinearPosFractions);
35  RooMomentMorph(const char *name, const char *title, RooAbsReal& _m, const RooArgList& varList,
36  const RooArgList& pdfList, const TVectorD& mrefpoints, Setting setting = NonLinearPosFractions );
37  RooMomentMorph(const RooMomentMorph& other, const char* name=0) ;
38  virtual TObject* clone(const char* newname) const { return new RooMomentMorph(*this,newname); }
39  virtual ~RooMomentMorph();
40 
41  void setMode(const Setting& setting) { _setting = setting; }
42 
43  void useHorizontalMorphing(bool val) { _useHorizMorph = val; }
44 
45  virtual Bool_t selfNormalized() const {
46  // P.d.f is self normalized
47  return kTRUE ;
48  }
49 
50  virtual Double_t getVal(const RooArgSet* set=0) const ;
51  RooAbsPdf* sumPdf(const RooArgSet* nset) ;
52 
53 
54 protected:
55 
56  class CacheElem : public RooAbsCacheElement {
57  public:
58  CacheElem(RooAbsPdf& sumPdf, RooChangeTracker& tracker, const RooArgList& flist) : _sumPdf(&sumPdf), _tracker(&tracker) { _frac.add(flist) ; } ;
60  virtual ~CacheElem() ;
65 
66  RooRealVar* frac(Int_t i ) ;
67  const RooRealVar* frac(Int_t i ) const ;
68  void calculateFractions(const RooMomentMorph& self, Bool_t verbose=kTRUE) const;
69  } ;
70  mutable RooObjCacheManager _cacheMgr ; //! The cache manager
71  mutable RooArgSet* _curNormSet ; //! Current normalization set
72 
73  friend class CacheElem ; // Cache needs to be able to clear _norm pointer
74 
75  Double_t evaluate() const ;
76 
77  void initialize();
78  CacheElem* getCache(const RooArgSet* nset) const ;
79 
80  inline Int_t ij(const Int_t& i, const Int_t& j) const { return (i*_varList.getSize()+j); }
81  int idxmin(const double& m) const;
82  int idxmax(const double& m) const;
83 
87  mutable TVectorD* _mref;
88 
89  TIterator* _varItr ; //! do not persist
91  mutable TMatrixD* _M; //
92 
94 
96 
97  ClassDef(RooMomentMorph,3) // Your description goes here...
98 };
99 
100 #endif
101 
102 
virtual Double_t getVal(const RooArgSet *set=0) const
Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
CacheElem * getCache(const RooArgSet *nset) const
Double_t evaluate() const
void setMode(const Setting &setting)
void useHorizontalMorphing(bool val)
TMatrixD * _M
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooAbsPdf * sumPdf(const RooArgSet *nset)
Iterator abstract base class.
Definition: TIterator.h:32
TMatrixT.
Definition: TMatrixDfwd.h:24
virtual Bool_t selfNormalized() const
#define ClassDef(name, id)
Definition: Rtypes.h:254
TIterator * _pdfItr
do not persist
void operModeHook(RooAbsArg::OperMode)
Interface for operation mode change calls.
RooArgSet * _curNormSet
The cache manager.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooListProxy _pdfList
Int_t getSize() const
RooChangeTracker * _tracker
int idxmax(const double &m) const
RooRealVar * frac(Int_t i)
RooSetProxy _varList
TVectorD * _mref
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
bool verbose
int idxmin(const double &m) const
RooListProxy is the concrete proxy for RooArgList objects.
Definition: RooListProxy.h:26
RooObjCacheManager _cacheMgr
virtual TObject * clone(const char *newname) const
Int_t ij(const Int_t &i, const Int_t &j) const
virtual ~RooMomentMorph()
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
RooChangeTracker is a meta object that tracks value changes in a given set of RooAbsArgs by registeri...
Mother of all ROOT objects.
Definition: TObject.h:37
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooRealProxy m
RooSetProxy is the concrete proxy for RooArgSet objects.
Definition: RooSetProxy.h:25
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
const Bool_t kTRUE
Definition: Rtypes.h:91
TIterator * _varItr
void calculateFractions(const RooMomentMorph &self, Bool_t verbose=kTRUE) const
char name[80]
Definition: TGX11.cxx:109
CacheElem(RooAbsPdf &sumPdf, RooChangeTracker &tracker, const RooArgList &flist)
RooMomentMorph()
coverity[UNINIT_CTOR]
virtual RooArgList containedArgs(Action)