ROOT  6.06/08
Reference Guide
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 class RooRealVar ;
22 
23 #include <map>
24 
26 public:
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=0) ;
34  virtual TObject* clone(const char* newname) const { return new RooFFTConvPdf(*this,newname); }
35  virtual ~RooFFTConvPdf() ;
36 
37  void setShift(Double_t val1, Double_t val2) { _shift1 = val1 ; _shift2 = val2 ; }
39  const RooArgSet& cacheObservables() const { return _cacheObs ; }
40 
42  // Return value of buffer fraction applied in FFT calculation array beyond either
43  // end of the observable domain to reduce cyclical effects
44  return _bufFrac ;
45  }
46 
47  enum BufStrat { Extend=0, Mirror=1, Flat=2 } ;
49  // Return the strategy currently used to fill the buffer:
50  // 'Extend' means is that the input p.d.f convolution observable range is widened to include the buffer range
51  // 'Flat' means that the buffer is filled with the p.d.f. value at the boundary of the observable range
52  // '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  }
55  void setBufferStrategy(BufStrat bs) ;
56  void setBufferFraction(Double_t frac) ;
57 
58  void printMetaArgs(std::ostream& os) const ;
59 
60  // Propagate maximum value estimate of pdf1 as convolution can only result in lower max values
61  virtual Int_t getMaxVal(const RooArgSet& vars) const { return _pdf1.arg().getMaxVal(vars) ; }
62  virtual Double_t maxVal(Int_t code) const { return _pdf1.arg().maxVal(code) ; }
63 
64 
65 protected:
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  Bool_t redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive) ;
75 
76  Double_t* scanPdf(RooRealVar& obs, RooAbsPdf& pdf, const RooDataHist& hist, const RooArgSet& slicePos, Int_t& N, Int_t& N2, Int_t& zeroBin, Double_t shift) const ;
77 
78  class FFTCacheElem : public PdfCacheElem {
79  public:
80  FFTCacheElem(const RooFFTConvPdf& self, const RooArgSet* nset) ;
81  ~FFTCacheElem() ;
82 
83  virtual RooArgList containedArgs(Action) ;
84 
88 
91 
94 
95  };
96 
97  friend class FFTCacheElem ;
98 
99  virtual Double_t evaluate() const { RooArgSet dummy(_x.arg()) ; return getVal(&dummy) ; } ; // dummy
100  virtual const char* inputBaseName() const ;
101  virtual RooArgSet* actualObservables(const RooArgSet& nset) const ;
102  virtual RooArgSet* actualParameters(const RooArgSet& nset) const ;
103  virtual RooAbsArg& pdfObservable(RooAbsArg& histObservable) const ;
104  virtual void fillCacheObject(PdfCacheElem& cache) const ;
105  void fillCacheSlice(FFTCacheElem& cache, const RooArgSet& slicePosition) const ;
106 
107  virtual PdfCacheElem* createCache(const RooArgSet* nset) const ;
108  virtual TString histNameSuffix() const ;
109 
110  // mutable std::map<const RooHistPdf*,CacheAuxInfo*> _cacheAuxInfo ; //! Auxilary Cache information (do not persist)
111  Double_t _bufFrac ; // Sampling buffer size as fraction of domain size
112  BufStrat _bufStrat ; // Strategy to fill the buffer
113 
116 
117  virtual RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=0,
118  const RooArgSet* auxProto=0, Bool_t verbose= kFALSE) const ;
119 
120  friend class RooConvGenContext ;
121  RooSetProxy _cacheObs ; // Non-convolution observables that are also cached
122 
123 private:
124 
125  ClassDef(RooFFTConvPdf,1) // Convolution operator p.d.f based on numeric Fourier transforms
126 };
127 
128 #endif
void calcParams()
(Re)calculate effective parameters of this p.d.f.
RooRealProxy _pdf1
Definition: RooFFTConvPdf.h:69
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
Definition: RooFFTConvPdf.h:62
RooSetProxy _cacheObs
RooAbsBinning * scanBinning
Definition: RooFFTConvPdf.h:93
virtual void fillCacheObject(PdfCacheElem &cache) const
Fill the contents of the cache the FFT convolution output.
Double_t _shift2
virtual RooArgSet * actualObservables(const RooArgSet &nset) const
Return the observables to be cached given the normalization set nset.
BufStrat _bufStrat
void setShift(Double_t val1, Double_t val2)
Definition: RooFFTConvPdf.h:37
virtual void removeAll()
Remove all argument inset using remove(const RooAbsArg&).
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
RooRealProxy _pdf2
Definition: RooFFTConvPdf.h:70
#define N
void setBufferStrategy(BufStrat bs)
Change strategy to fill the overflow buffer on either side of the convolution observable range...
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Create appropriate generator context for this convolution.
Basic string class.
Definition: TString.h:137
void printMetaArgs(std::ostream &os) const
Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the p...
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
FFTCacheElem(const RooFFTConvPdf &self, const RooArgSet *nset)
Clone input pdf and attach to dataset.
virtual PdfCacheElem * createCache(const RooArgSet *nset) const
Return specialized cache subclass for FFT calculations.
Double_t _bufFrac
#define ClassDef(name, id)
Definition: Rtypes.h:254
const RooArgSet & cacheObservables() const
Definition: RooFFTConvPdf.h:39
virtual TString histNameSuffix() const
Suffix for cache histogram (added in addition to suffix for cache name)
virtual RooArgList containedArgs(Action)
Returns all RooAbsArg objects contained in the cache element.
RooRealProxy _x
Definition: RooFFTConvPdf.h:67
Double_t bufferFraction() const
Definition: RooFFTConvPdf.h:41
BufStrat bufferStrategy() const
Definition: RooFFTConvPdf.h:48
Double_t * scanPdf(RooRealVar &obs, RooAbsPdf &pdf, const RooDataHist &hist, const RooArgSet &slicePos, Int_t &N, Int_t &N2, Int_t &zeroBin, Double_t shift) const
Scan the values of &#39;pdf&#39; in observable &#39;obs&#39; using the bin values stored in &#39;hist&#39; at slice position ...
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
void fillCacheSlice(FFTCacheElem &cache, const RooArgSet &slicePosition) const
Fill a slice of cachePdf with the output of the FFT convolution calculation.
bool verbose
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
Definition: RooFFTConvPdf.h:61
virtual RooAbsArg & pdfObservable(RooAbsArg &histObservable) const
Return p.d.f.
TVirtualFFT is an interface class for Fast Fourier Transforms.
Definition: TVirtualFFT.h:92
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
static RooMathCoreReg dummy
virtual RooArgSet * actualParameters(const RooArgSet &nset) const
Return the parameters on which the cache depends given normalization set nset.
RooSetProxy _params
Definition: RooFFTConvPdf.h:71
Double_t _shift1
#define name(a, b)
Definition: linkTestLib0.cpp:5
Mother of all ROOT objects.
Definition: TObject.h:58
RooAbsBinning * histBinning
Definition: RooFFTConvPdf.h:92
virtual TObject * clone(const char *newname) const
Definition: RooFFTConvPdf.h:34
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
void setBufferFraction(Double_t frac)
Change the size of the buffer on either side of the observable range to frac times the size of the ra...
virtual const char * inputBaseName() const
Return base name component for cache components in this case &#39;PDF1_CONV_PDF2&#39;.
virtual ~RooFFTConvPdf()
Destructor.
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
virtual Double_t evaluate() const
Definition: RooFFTConvPdf.h:99
Bool_t redirectServersHook(const RooAbsCollection &newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
calcParams() ;
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
RooRealProxy _xprime
Definition: RooFFTConvPdf.h:68
void setCacheObservables(const RooArgSet &obs)
Definition: RooFFTConvPdf.h:38
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts &#39;var&#39; into set and registers &#39;var&#39; as server to owner with...