Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorph.h
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2023, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11#ifndef ROOMOMENTMORPH
12#define ROOMOMENTMORPH
13
14#include "RooAbsPdf.h"
15#include "RooRealProxy.h"
16#include "RooCategoryProxy.h"
17#include "RooAbsReal.h"
18#include "RooAbsCategory.h"
19#include "RooSetProxy.h"
20#include "RooListProxy.h"
21#include "RooArgList.h"
22
23#include "TMatrixD.h"
24#include "TVectorD.h"
25
27
28class RooMomentMorph : public RooAbsPdf {
29public:
30
32
34
35 RooMomentMorph(const char *name, const char *title, RooAbsReal& _m, const RooArgList& varList,
36 const RooArgList& pdfList, const RooArgList& mrefList, Setting setting = NonLinearPosFractions);
37 RooMomentMorph(const char *name, const char *title, RooAbsReal& _m, const RooArgList& varList,
38 const RooArgList& pdfList, const TVectorD& mrefpoints, Setting setting = NonLinearPosFractions );
39 RooMomentMorph(const RooMomentMorph& other, const char* name=nullptr) ;
40 TObject* clone(const char* newname) const override { return new RooMomentMorph(*this,newname); }
41 ~RooMomentMorph() override;
42
43 void setMode(const Setting& setting) { _setting = setting; }
44
45 void useHorizontalMorphing(bool val) { _useHorizMorph = val; }
46
47 bool selfNormalized() const override {
48 // P.d.f is self normalized
49 return true ;
50 }
51
52 double getValV(const RooArgSet* set=nullptr) const override;
53 RooAbsPdf* sumPdf(const RooArgSet* nset) ;
54
55
56protected:
57
59 public:
60 CacheElem(std::unique_ptr<RooAbsPdf> && sumPdf,
61 std::unique_ptr<RooChangeTracker> && tracker,
62 const RooArgList& flist);
63 ~CacheElem() override ;
65 std::unique_ptr<RooAbsPdf> _sumPdf ;
66 std::unique_ptr<RooChangeTracker> _tracker ;
68
69 RooRealVar* frac(Int_t i ) ;
70 const RooRealVar* frac(Int_t i ) const ;
71 void calculateFractions(const RooMomentMorph& self, bool verbose=true) const;
72 } ;
73 mutable RooObjCacheManager _cacheMgr ; //! The cache manager
74 mutable RooArgSet* _curNormSet = nullptr; //! Current normalization set
75
76 friend class CacheElem ; // Cache needs to be able to clear _norm pointer
77
78 double evaluate() const override ;
79
80 void initialize();
81 CacheElem* getCache(const RooArgSet* nset) const ;
82
83 inline Int_t ij(const Int_t& i, const Int_t& j) const { return (i*_varList.size()+j); }
84 int idxmin(const double& m) const;
85 int idxmax(const double& m) const;
86
90 mutable TVectorD* _mref = nullptr;
91
92 mutable TMatrixD* _M = nullptr;
93
95
96 bool _useHorizMorph = true;
97
99};
100
101#endif
102
103
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
char name[80]
Definition TGX11.cxx:110
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
Storage_t::size_type size() const
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:24
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
~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)
double getValV(const RooArgSet *set=nullptr) const override
Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set...
CacheElem * getCache(const RooArgSet *nset) const
TVectorD * _mref
Implementation of a RooCacheManager<RooAbsCacheElement> that specializes in the storage of cache elem...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41