Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooFFTConvPdf.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * Copyright (c) 2000-2005, Regents of the University of California *
5 * and Stanford University. All rights reserved. *
6 * *
7 * Redistribution and use in source and binary forms, *
8 * with or without modification, are permitted according to the terms *
9 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
10 *****************************************************************************/
11
12#ifndef ROOFFTCONVPDF
13#define ROOFFTCONVPDF
14
15#include "RooAbsCachedPdf.h"
16#include "RooRealProxy.h"
17#include "RooSetProxy.h"
18#include "RooAbsReal.h"
19#include "RooHistPdf.h"
20#include "TVirtualFFT.h"
21
22class RooRealVar;
23
24///PDF for the numerical (FFT) convolution of two PDFs.
26public:
27
29 // coverity[UNINIT_CTOR]
30 } ;
31 RooFFTConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder=2);
32 RooFFTConvPdf(const char *name, const char *title, RooAbsReal& pdfConvVar, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder=2);
33 RooFFTConvPdf(const RooFFTConvPdf& other, const char* name=nullptr) ;
34 TObject* clone(const char* newname) const override { return new RooFFTConvPdf(*this,newname); }
35 ~RooFFTConvPdf() override ;
36
37 void setShift(double val1, double val2) { _shift1 = val1 ; _shift2 = val2 ; }
39 const RooArgSet& cacheObservables() const { return _cacheObs ; }
40
41 /// Return value of buffer fraction applied in FFT calculation array beyond either
42 /// end of the observable domain to reduce cyclical effects
43 double bufferFraction() const {
44 return _bufFrac ;
45 }
46
47 enum BufStrat { Extend=0, Mirror=1, Flat=2 } ;
48 /// Return the strategy currently used to fill the buffer:
49 /// 'Extend' means is that the input p.d.f convolution observable range is widened to include the buffer range
50 /// 'Flat' means that the buffer is filled with the p.d.f. value at the boundary of the observable range
51 /// 'Mirror' means that the buffer is filled with a mirror image of the p.d.f. around the convolution observable boundary
53 return _bufStrat ;
54 }
56 void setBufferFraction(double frac) ;
57
58 void printMetaArgs(std::ostream& os) const override ;
59
60 // Propagate maximum value estimate of pdf1 as convolution can only result in lower max values
61 Int_t getMaxVal(const RooArgSet& vars) const override { return _pdf1.arg().getMaxVal(vars) ; }
62 double maxVal(Int_t code) const override { return _pdf1.arg().maxVal(code) ; }
63
64
65protected:
66
67 RooRealProxy _x ; ///< Convolution observable
68 RooRealProxy _xprime ; ///< Input function representing value of convolution observable
69 RooRealProxy _pdf1 ; ///< First input p.d.f
70 RooRealProxy _pdf2 ; ///< Second input p.d.f
71 RooSetProxy _params ; ///< Effective parameters of this p.d.f.
72
73 void calcParams() ;
74
75 std::vector<double> scanPdf(RooRealVar& obs, RooAbsPdf& pdf, const RooDataHist& hist, const RooArgSet& slicePos, Int_t& N, Int_t& N2, Int_t& zeroBin, double shift) const ;
76
77 class FFTCacheElem : public PdfCacheElem {
78 public:
79 FFTCacheElem(const RooFFTConvPdf& self, const RooArgSet* nset) ;
80
82
83 std::unique_ptr<TVirtualFFT> fftr2c1;
84 std::unique_ptr<TVirtualFFT> fftr2c2;
85 std::unique_ptr<TVirtualFFT> fftc2r;
86
87 std::unique_ptr<RooAbsPdf> pdf1Clone;
88 std::unique_ptr<RooAbsPdf> pdf2Clone;
89
90 std::unique_ptr<RooAbsBinning> histBinning;
91 std::unique_ptr<RooAbsBinning> scanBinning;
92 };
93
94 friend class FFTCacheElem ;
95
96 double evaluate() const override { RooArgSet dummy(_x.arg()) ; return getVal(&dummy) ; } ; // dummy
97 const char* inputBaseName() const override ;
98 RooFit::OwningPtr<RooArgSet> actualObservables(const RooArgSet& nset) const override ;
99 RooFit::OwningPtr<RooArgSet> actualParameters(const RooArgSet& nset) const override ;
100 RooAbsArg& pdfObservable(RooAbsArg& histObservable) const override ;
101 void fillCacheObject(PdfCacheElem& cache) const override ;
102 void fillCacheSlice(FFTCacheElem& cache, const RooArgSet& slicePosition) const ;
103
104 PdfCacheElem* createCache(const RooArgSet* nset) const override ;
105 TString histNameSuffix() const override ;
106
107 // mutable std::map<const RooHistPdf*,CacheAuxInfo*> _cacheAuxInfo ; //! Auxiliary Cache information (do not persist)
108 double _bufFrac ; // Sampling buffer size as fraction of domain size
109 BufStrat _bufStrat ; // Strategy to fill the buffer
110
111 double _shift1 ;
112 double _shift2 ;
113
114 RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr,
115 const RooArgSet* auxProto=nullptr, bool verbose= false) const override ;
116
117 friend class RooConvGenContext ;
118 RooSetProxy _cacheObs ; ///< Non-convolution observables that are also cached
119
120private:
121
122 void prepareFFTBinning(RooRealVar& convVar) const;
123
124 ClassDefOverride(RooFFTConvPdf,1) // Convolution operator p.d.f based on numeric Fourier transforms
125};
126
127#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
#define N
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
Abstract base class for p.d.f.s that need or want to cache their evaluate() output in a RooHistPdf de...
Abstract base class for generator contexts of RooAbsPdf objects.
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
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual double maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
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
void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
RooConvGenContext is an efficient implementation of the generator context specific for RooAbsAnaConvP...
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
std::unique_ptr< TVirtualFFT > fftr2c2
std::unique_ptr< TVirtualFFT > fftc2r
RooArgList containedArgs(Action) override
Returns all RooAbsArg objects contained in the cache element.
std::unique_ptr< RooAbsPdf > pdf2Clone
std::unique_ptr< RooAbsBinning > histBinning
std::unique_ptr< RooAbsBinning > scanBinning
std::unique_ptr< RooAbsPdf > pdf1Clone
std::unique_ptr< TVirtualFFT > fftr2c1
PDF for the numerical (FFT) convolution of two PDFs.
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
RooSetProxy _params
Effective parameters of this p.d.f.
TObject * clone(const char *newname) const override
BufStrat _bufStrat
void calcParams()
(Re)calculate effective parameters of this p.d.f.
BufStrat bufferStrategy() const
Return the strategy currently used to fill the buffer: 'Extend' means is that the input p....
void setBufferFraction(double frac)
Change the size of the buffer on either side of the observable range to frac times the size of the ra...
double bufferFraction() const
Return value of buffer fraction applied in FFT calculation array beyond either end of the observable ...
TString histNameSuffix() const override
Suffix for cache histogram (added in addition to suffix for cache name)
void setCacheObservables(const RooArgSet &obs)
void prepareFFTBinning(RooRealVar &convVar) const
Try to improve the binning and inform user if possible.
void fillCacheSlice(FFTCacheElem &cache, const RooArgSet &slicePosition) const
Fill a slice of cachePdf with the output of the FFT convolution calculation.
RooRealProxy _xprime
Input function representing value of convolution observable.
std::vector< double > scanPdf(RooRealVar &obs, RooAbsPdf &pdf, const RooDataHist &hist, const RooArgSet &slicePos, Int_t &N, Int_t &N2, Int_t &zeroBin, double shift) const
Scan the values of 'pdf' in observable 'obs' using the bin values stored in 'hist' at slice position ...
RooRealProxy _pdf1
First input p.d.f.
RooRealProxy _x
Convolution observable.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
const RooArgSet & cacheObservables() const
const char * inputBaseName() const override
Return base name component for cache components in this case 'PDF1_CONV_PDF2'.
RooAbsArg & pdfObservable(RooAbsArg &histObservable) const override
Return p.d.f.
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise capability to determine maximum value of function for given set of observables.
~RooFFTConvPdf() override
Destructor.
void setShift(double val1, double val2)
RooFit::OwningPtr< RooArgSet > actualParameters(const RooArgSet &nset) const override
Return the parameters on which the cache depends given normalization set nset.
RooRealProxy _pdf2
Second input p.d.f.
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the p...
void setBufferStrategy(BufStrat bs)
Change strategy to fill the overflow buffer on either side of the convolution observable range.
PdfCacheElem * createCache(const RooArgSet *nset) const override
Return specialized cache subclass for FFT calculations.
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Create appropriate generator context for this convolution.
RooFit::OwningPtr< RooArgSet > actualObservables(const RooArgSet &nset) const override
Return the observables to be cached given the normalization set nset.
void fillCacheObject(PdfCacheElem &cache) const override
Fill the contents of the cache the FFT convolution output.
RooSetProxy _cacheObs
Non-convolution observables that are also cached.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
const T & arg() const
Return reference to object held in proxy.
Mother of all ROOT objects.
Definition TObject.h:41
Basic string class.
Definition TString.h:139
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:43