134#ifndef ROOFIT_MATH_FFTW3
141 static void (*doFFT)(
int,
double *,
double *,
double *) =
nullptr;
146 bool success =
gInterpreter->Declare(
"#include \"fftw3.h\"");
149 std::stringstream ss;
150 ss <<
"RooFFTConvPdf evaluation Failed! The interpreter could not find the fftw3.h header.\n";
151 ss <<
"You have three options to fix this problem:\n";
152 ss <<
" 1) Install fftw3 on your system so that the interpreter can find it\n";
153 ss <<
" 2) In case fftw3.h is installed somewhere else,\n"
154 <<
" tell ROOT with gInterpreter->AddIncludePath() where to find it\n";
155 ss <<
" 3) Compile ROOT with the -Dfftw3=ON in the CMake configuration,\n"
156 <<
" such that ROOT comes with built-in fftw3 convolution routines\n";
157 oocoutE(
nullptr, Eval) << ss.str() << std::endl;
158 throw std::runtime_error(
"RooFFTConvPdf evaluation Failed! The interpreter could not find the fftw3.h header");
162void RooFFTConvPdf_doFFT(int n, double *input1, double *input2, double *output)
164 auto fftr2c1_Out = reinterpret_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex) * (n / 2 + 1)));
165 auto fftr2c2_Out = reinterpret_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex) * (n / 2 + 1)));
166 auto fftc2r_In = reinterpret_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex) * (n / 2 + 1)));
168 fftw_plan fftr2c1_plan = fftw_plan_dft_r2c(1, &n, input1, fftr2c1_Out, FFTW_ESTIMATE);
169 fftw_plan fftr2c2_plan = fftw_plan_dft_r2c(1, &n, input2, fftr2c2_Out, FFTW_ESTIMATE);
170 fftw_plan fftc2r_plan = fftw_plan_dft_c2r(1, &n, fftc2r_In, output, FFTW_ESTIMATE);
172 // Real->Complex FFT Transform on p.d.f. samplings
173 fftw_execute(fftr2c1_plan);
174 fftw_execute(fftr2c2_plan);
176 // Loop over first half +1 of complex output results, multiply
177 // and set as input of reverse transform
178 for (Int_t i = 0; i < n / 2 + 1; i++) {
179 double re1 = fftr2c1_Out[i][0];
180 double re2 = fftr2c2_Out[i][0];
181 double im1 = fftr2c1_Out[i][1];
182 double im2 = fftr2c2_Out[i][1];
183 fftc2r_In[i][0] = re1 * re2 - im1 * im2;
184 fftc2r_In[i][1] = re1 * im2 + re2 * im1;
187 // Reverse Complex->Real FFT transform product
188 fftw_execute(fftc2r_plan);
190 fftw_destroy_plan(fftr2c1_plan);
191 fftw_destroy_plan(fftr2c2_plan);
192 fftw_destroy_plan(fftc2r_plan);
194 fftw_free(fftr2c1_Out);
195 fftw_free(fftr2c2_Out);
196 fftw_free(fftc2r_In);
200 doFFT = reinterpret_cast<void(*)(
int,
double*,
double*,
double*)
>(
gInterpreter->ProcessLine(
"RooFFTConvPdf_doFFT;"));
225 _x(
"!x",
"Convolution Variable",this,convVar),
226 _xprime(
"!xprime",
"External Convolution Variable",this,false),
227 _pdf1(
"!pdf1",
"pdf1",this,pdf1,false),
228 _pdf2(
"!pdf2",
"pdf2",this,pdf2,false),
229 _params(
"!params",
"effective parameters",this),
234 _cacheObs(
"!cacheObs",
"Cached observables",this,false,false)
251 _x(
"!x",
"Convolution Variable",this,convVar,false,false),
252 _xprime(
"!xprime",
"External Convolution Variable",this,pdfConvVar),
253 _pdf1(
"!pdf1",
"pdf1",this,pdf1,false),
254 _pdf2(
"!pdf2",
"pdf2",this,pdf2,false),
255 _params(
"!params",
"effective parameters",this),
260 _cacheObs(
"!cacheObs",
"Cached observables",this,false,false)
276 _x(
"!x",this,other._x),
277 _xprime(
"!xprime",this,other._xprime),
278 _pdf1(
"!pdf1",this,other._pdf1),
279 _pdf2(
"!pdf2",this,other._pdf2),
280 _params(
"!params",this,other._params),
281 _bufFrac(other._bufFrac),
282 _bufStrat(other._bufStrat),
283 _shift1(other._shift1),
284 _shift2(other._shift2),
285 _cacheObs(
"!cacheObs",this,other._cacheObs)
311 coutI(Caching) <<
"Changing internal binning of variable '" << convVar.
GetName()
312 <<
"' in FFT '" <<
fName <<
"'"
313 <<
" from " << varBinning.
numBins()
314 <<
" to " << optimal <<
" to improve the precision of the numerical FFT."
315 <<
" This can be done manually by setting an additional binning named 'cache'." << std::endl;
318 coutE(Caching) <<
"The internal binning of variable " << convVar.
GetName()
319 <<
" is not uniform. The numerical FFT will likely yield wrong results." << std::endl;
333 name.Append(
"_CONV_") ;
367 string refName =
Form(
"refrange_fft_%s",self.
GetName()) ;
380 pdf1Clone->addOwnedComponents(*shiftObs1) ;
381 pdf1Clone->addOwnedComponents(*clonePdf1) ;
395 pdf1Clone->addOwnedComponents(*shiftObs2) ;
396 pdf1Clone->addOwnedComponents(*clonePdf2) ;
414 pdf1Clone->recursiveRedirectServers(fftParams) ;
415 pdf2Clone->recursiveRedirectServers(fftParams) ;
416 pdf1Clone->fixAddCoefRange(refName.c_str(),
true) ;
417 pdf2Clone->fixAddCoefRange(refName.c_str(),
true) ;
421 pdf1Clone->fixAddCoefNormalization(convSet,
true);
422 pdf2Clone->fixAddCoefNormalization(convSet,
true);
429 oocoutW(&self, Eval) <<
"The FFT convolution '" << self.
GetName() <<
"' will run with " <<
N
430 <<
" bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number"
431 " of bins of the observable '" << convObs->
GetName() <<
"'." << std::endl;
464 ret.
add(*pdf1Clone) ;
465 ret.
add(*pdf2Clone) ;
466 if (pdf1Clone->ownedComponents()) {
467 ret.
add(*pdf1Clone->ownedComponents()) ;
469 if (pdf2Clone->ownedComponents()) {
470 ret.
add(*pdf2Clone->ownedComponents()) ;
494 otherObs.
remove(*histArg,
true,
true) ;
500 if (otherObs.
empty()) {
510 std::vector<int> binCur(
n+1);
511 std::vector<int> binMax(
n+1);
514 std::vector<RooAbsLValue*> obsLV(
n);
516 for (
auto const& arg : otherObs) {
536 while(binCur[curObs]==binMax[curObs]) {
577 Int_t N,N2,binShift1,binShift2 ;
585#ifndef ROOFIT_MATH_FFTW3
588 std::vector<double>
output(N2);
590 auto doFFT = declareDoFFT();
591 doFFT(N2, input1.data(), input2.data(),
output.data());
604 coutF(Eval) <<
"RooFFTConvPdf::fillCacheSlice(" <<
GetName() <<
"Cannot get a handle to fftw. Maybe ROOT was built without it?" << std::endl;
605 throw std::runtime_error(
"Cannot get a handle to fftw.");
610 aux.
fftr2c1->SetPoints(input1.data());
614 aux.
fftr2c2->SetPoints(input2.data());
619 for (
Int_t i=0 ; i<N2/2+1 ; i++) {
620 double re1,re2,im1,im2 ;
621 aux.
fftr2c1->GetPointComplex(i,re1,im1) ;
622 aux.
fftr2c2->GetPointComplex(i,re2,im2) ;
623 double re = re1*re2 - im1*im2 ;
624 double im = re1*im2 + re2*im1 ;
626 aux.
fftc2r->SetPointComplex(i,t) ;
633 Int_t totalShift = binShift1 + (N2-
N)/2 ;
637 std::unique_ptr<TIterator> iter{
const_cast<RooDataHist&
>(cacheHist).sliceIterator(
const_cast<RooAbsReal&
>(
_x.
arg()),slicePos)};
638 for (
Int_t i =0 ; i<
N ; i++) {
641 Int_t j = i + totalShift ;
643 while (j>=N2) j-= N2 ;
646#ifndef ROOFIT_MATH_FFTW3
649 cacheHist.
set(aux.
fftc2r->GetPointReal(j));
674 std::vector<double> array(N2);
683 }
else if (histX->
getMin()>0) {
693 zeroBin += binShift ;
694 while(zeroBin>=N2) zeroBin-= N2 ;
695 while(zeroBin<0) zeroBin+= N2 ;
698 std::vector<double> tmp(N2);
704 for (k=0 ; k<N2 ; k++) {
716 for (k=0 ; k<Nbuf ; k++) {
719 for (k=0 ; k<
N ; k++) {
725 for (k=0 ; k<Nbuf ; k++) {
726 tmp[
N+Nbuf+k] = val ;
734 for (k=0 ; k<
N ; k++) {
738 for (k=1 ; k<=Nbuf ; k++) {
748 for (
Int_t i=0 ; i<N2 ; i++) {
750 Int_t j = i - (zeroBin) ;
781 obs1->add(obs2,
true) ;
788 for(
auto const& arg : *obs1) {
793 obs1->remove(killList) ;
796 obs1->add(
_x.
arg(),
true) ;
805 for(
auto const& arg : *obs1) {
810 obs1->remove(killList) ;
815 obs1->add(
_x.
arg(),
true) ;
848 return histObservable ;
861 const RooArgSet* auxProto,
bool verbose)
const
875 <<
" has internal generator that is safe to use in current context" << endl ;
879 <<
" has internal generator that is safe to use in current context" << endl ;
882 cxcoutI(Generation) <<
"RooFFTConvPdf::genContext() generation requested for observables other than the convolution observable " <<
_x.
arg().
GetName() << endl ;
886 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
889 cxcoutI(Generation) <<
"RooFFTConvPdf::genContext() selecting accept/reject generator context because one or both of the input "
890 <<
"p.d.f.s cannot use internal generator and/or "
891 <<
"observables other than the convolution variable are requested for generation" << endl ;
892 return new RooGenContext(*
this,vars,prototype,auxProto,verbose) ;
896 cxcoutI(Generation) <<
"RooFFTConvPdf::genContext() selecting specialized convolution generator context as both input "
897 <<
"p.d.fs are safe for internal generator and only "
898 <<
"the convolution observables is requested for generation" << endl ;
911 coutE(InputArguments) <<
"RooFFTConvPdf::setBufferFraction(" <<
GetName() <<
") fraction should be greater than or equal to zero" << endl ;
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Common abstract base class for objects that represent a value and a "shape" in RooFit.
void setOperMode(OperMode mode, bool recurseADirty=true)
Set the operation mode of this node.
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
virtual RooAbsArg * cloneTree(const char *newname=nullptr) const
Clone tree expression of objects.
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
int binNumber(double x) const
Returns the bin number corresponding to the value x.
Int_t numBins() const
Return number of bins.
virtual bool isUniform() const
virtual double highBound() const =0
virtual double lowBound() const =0
virtual RooAbsBinning * clone(const char *name=nullptr) const =0
Abstract base class for p.d.f.s that need or want to cache their evaluate() output in a RooHistPdf de...
virtual const char * binningName() const
RooObjCacheManager _cacheMgr
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
void setDirtyProp(bool flag)
Control propagation of dirty flags from observables in dataset.
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract base class for objects that are lvalues, i.e.
virtual Int_t numBins(const char *rangeName=nullptr) const =0
Abstract interface for all probability density functions.
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Int_t numBins(const char *rangeName=nullptr) const override
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
void setBin(Int_t ibin, const char *rangeName=nullptr) override
Set value to center of bin 'ibin' of binning 'rangeName' (or of default binning if no range is specif...
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
Abstract base class for objects that represent a real value and implements functionality common to al...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooAbsArg * absArg() const
Return pointer to contained argument.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
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...
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurrence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, bool verbose=false)
Build a clone of the prototype executing all registered 'replace' rules and 'split' rules for the mas...
The RooDataHist is a container class to hold N-dimensional binned data.
void set(std::size_t binNumber, double weight, double wgtErr)
Set bin content of bin that was last loaded with get(std::size_t).
const RooArgSet * get() const override
Get bin centre of current bin.
RooDataSet is a container class to hold unbinned data.
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
FFTCacheElem(const RooFFTConvPdf &self, const RooArgSet *nset)
Clone input pdf and attach to dataset.
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.
friend class RooConvGenContext
RooSetProxy _params
Effective parameters of this p.d.f.
void calcParams()
(Re)calculate effective parameters of this p.d.f.
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 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.
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.
~RooFFTConvPdf() override
Destructor.
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.
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
RooLinearVar is the most general form of a derived real-valued object that can be used by RooRealInte...
void sterilize() override
Clear the cache payload but retain slot mapping w.r.t to normalization and integration sets.
RooRealVar represents a variable that can be changed from the outside.
bool hasBinning(const char *name) const override
Returns true if variable has a binning named 'name'.
const RooAbsBinning & getBinning(const char *name=nullptr, bool verbose=true, bool createOnTheFly=false) const override
Return binning definition with name.
void setBinning(const RooAbsBinning &binning, const char *name=nullptr)
Add given binning under name 'name' with this variable.
void setRange(const char *name, double min, double max)
Set a fit or plotting range.
const T & arg() const
Return reference to object held in proxy.
const char * GetName() const override
Returns name of object.
static TVirtualFFT * FFT(Int_t ndim, Int_t *n, Option_t *option)
Returns a pointer to the FFT of requested size and type.
RooConstVar & RooConst(double val)
void(off) SmallVectorTemplateBase< T
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...