Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooFFTConvPdf.cxx
Go to the documentation of this file.
1
2 /*****************************************************************************
3 * Project: RooFit *
4 * *
5 * Copyright (c) 2000-2005, Regents of the University of California *
6 * and Stanford University. All rights reserved. *
7 * *
8 * Redistribution and use in source and binary forms, *
9 * with or without modification, are permitted according to the terms *
10 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
11 *****************************************************************************/
12
13//////////////////////////////////////////////////////////////////////////////
14/// \class RooFFTConvPdf
15/// \ingroup Roofitcore
16///
17/// This class implements a generic one-dimensional numeric convolution of two PDFs,
18/// and can convolve any two RooAbsPdfs. The class exploits the convolution theorem
19/// \f[
20/// f(x) * g(x) \rightarrow F(k_i) \cdot G(k_i)
21/// \f]
22/// to calculate the convolution by calculating a Real->Complex FFT of both input PDFs,
23/// multiplying the complex coefficients and performing the reverse Complex->Real FFT
24/// to get the result in the input space. This class uses the ROOT FFT interface to
25/// the (free) FFTW3 package (www.fftw.org), and requires that your ROOT installation is
26/// compiled with the `fftw3=ON` (default). Instructions for manually installing fftw below.
27///
28/// Note that the performance in terms of speed and stability of RooFFTConvPdf is
29/// vastly superior to that of RooNumConvPdf.
30///
31/// An important feature of FFT convolutions is that the observable is assumed to be
32/// cyclical. This is correct for cyclical observables such as angles,
33/// but does not hold in general. For non-cyclical variables, wrap-around artifacts may be
34/// encountered, *e.g.* if the PDF is zero at xMin and non-zero at xMax. A rising tail may appear at xMin.
35/// This is inevitable when using FFTs. A distribution with 3 bins therefore looks like:
36/// ```
37/// ... 0 1 2 0 1 2 0 1 2 ...
38/// ```
39///
40/// Therefore, if bins 0 and 2 are not equal, the FFT sees a cyclical function with a step at the 2|0 boundary, which causes
41/// artifacts in Fourier space.
42///
43/// The spillover or discontinuity can be reduced or eliminated by
44/// introducing a buffer zone in the FFT calculation. If this feature is activated (on by default),
45/// the sampling array for the FFT calculation is extended in both directions,
46/// and padded with the lowest/highest bin.
47/// Example:
48/// ```
49/// original: -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5
50/// add buffer zones: U U -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5 O O
51/// rotate: 0 +1 +2 +3 +4 +5 O O U U -5 -4 -3 -2 -1
52/// ```
53/// The buffer bins are stripped away when the FFT output values
54/// are transferred back to the p.d.f cache. The default buffer size is 10% of the
55/// observable domain size, and can be changed with the `setBufferFraction()` member function.
56///
57/// The RooFFTConvPdf uses caching inherited from a RooAbsCachedPdf. If it is
58/// evaluated for a particular value of x, the FFT and convolution is calculated
59/// for all bins in the observable space for the given choice of parameters,
60/// which are also stored in the cache. Subsequent evaluations for different values of the convolution observable and
61/// identical parameters will be retrieved from the cache. If one or more
62/// of the parameters change, the cache will be updated, *i.e.*, a new FFT runs.
63///
64/// The sampling density of the FFT is controlled by the binning of the
65/// the convolution observable, which can be changed using RooRealVar::setBins(N).
66/// For good results, N should be large (>=1000). Additional interpolation
67/// between the bins may improve the result if coarse binnings are chosen. These can be
68/// activated in the constructor or by calling `setInterpolationOrder()`.
69/// For N >> 1000, interpolation will not substantially improve the accuracy.
70///
71/// Additionial information on caching can be displayed by monitoring
72/// the message stream with topic "Caching" at the INFO level, *i.e.*
73/// by calling `RooMsgService::instance().addStream(RooMsgService::INFO,Topic("Caching"))`
74/// to see these message on stdout.
75///
76/// Multi-dimensional convolutions are not supported at the moment.
77///
78/// ---
79///
80/// Installing an external version of FFTW on Linux and compiling ROOT to use it
81/// -------
82///
83/// You have two options:
84/// * **Recommended**: ROOT can automatically install FFTW for itself, see `builtin_fftw3` at https://root.cern/building-root
85/// * Install FFTW and let ROOT discover it. `fftw3` is on by default (see https://root.cern/building-root)
86///
87/// 1) Go to www.fftw.org and download the latest stable version (a .tar.gz file)
88///
89/// If you have root access to your machine and want to make a system installation of FFTW
90///
91/// 2) Untar fftw-XXX.tar.gz in /tmp, cd into the untarred directory
92/// and type './configure' followed by 'make install'.
93/// This will install fftw in /usr/local/bin,lib etc...
94///
95/// 3) Start from a source installation of ROOT. ROOT should discover it. See https://root.cern/building-root
96///
97///
98/// If you do not have root access and want to make a private installation of FFTW
99///
100/// 2) Make a private install area for FFTW, e.g. /home/myself/fftw
101///
102/// 3) Untar fftw-XXX.tar.gz in /tmp, cd into the untarred directory
103/// and type './configure --prefix=/home/myself/fftw' followed by 'make install'.
104/// Substitute /home/myself/fftw with a directory of your choice. This
105/// procedure will install FFTW in the location designated by you
106///
107/// 4) Start from a source installation of ROOT.
108/// Look up and set the proper paths for ROOT to discover FFTW. See https://root.cern/building-root
109///
110
111#include "RooFFTConvPdf.h"
112
113#include "RooAbsReal.h"
114#include "RooMsgService.h"
115#include "RooDataHist.h"
116#include "RooHistPdf.h"
117#include "RooRealVar.h"
118#include "RooGenContext.h"
119#include "RooConvGenContext.h"
120#include "RooBinning.h"
121#include "RooLinearVar.h"
122#include "RooCustomizer.h"
123#include "RooGlobalFunc.h"
124#include "RooConstVar.h"
125#include "RooUniformBinning.h"
126
127#include "TClass.h"
128#include "TComplex.h"
129#include "TVirtualFFT.h"
130
131#include <iostream>
132#include <stdexcept>
133
134#ifndef ROOFIT_MATH_FFTW3
135#include "TInterpreter.h"
136
137namespace {
138
139auto declareDoFFT()
140{
141 static void (*doFFT)(int, double *, double *, double *) = nullptr;
142
143 if (doFFT)
144 return doFFT;
145
146 bool success = gInterpreter->Declare("#include \"fftw3.h\"");
147
148 if (!success) {
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");
159 }
160
161 gInterpreter->Declare(R"(
162void RooFFTConvPdf_doFFT(int n, double *input1, double *input2, double *output)
163{
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)));
167
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);
171
172 // Real->Complex FFT Transform on p.d.f. samplings
173 fftw_execute(fftr2c1_plan);
174 fftw_execute(fftr2c2_plan);
175
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;
185 }
186
187 // Reverse Complex->Real FFT transform product
188 fftw_execute(fftc2r_plan);
189
190 fftw_destroy_plan(fftr2c1_plan);
191 fftw_destroy_plan(fftr2c2_plan);
192 fftw_destroy_plan(fftc2r_plan);
193
194 fftw_free(fftr2c1_Out);
195 fftw_free(fftr2c2_Out);
196 fftw_free(fftc2r_In);
197}
198)");
199
200 doFFT = reinterpret_cast<void(*)(int, double*, double*, double*)>(gInterpreter->ProcessLine("RooFFTConvPdf_doFFT;"));
201 return doFFT;
202}
203
204} // namespace
205
206#endif
207
208using std::endl, std::string, std::ostream;
209
210
211////////////////////////////////////////////////////////////////////////////////
212/// Constructor for numerical (FFT) convolution of PDFs.
213/// \param[in] name Name of this PDF
214/// \param[in] title Title for plotting this PDF
215/// \param[in] convVar Observable to convolve the PDFs in \attention Use a high number of bins (>= 1000) for good accuracy.
216/// \param[in] pdf1 First PDF to be convolved
217/// \param[in] pdf2 Second PDF to be convolved
218/// \param[in] ipOrder Order for interpolation between bins (since FFT is discrete)
219/// The binning used for the FFT sampling is controlled by the binning named "cache" in the convolution observable `convVar`.
220/// If such a binning is not set, the same number of bins as for `convVar` will be used.
221
222RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooRealVar &convVar, RooAbsPdf &pdf1, RooAbsPdf &pdf2,
224 : RooAbsCachedPdf(name, title, ipOrder),
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),
230 _bufFrac(0.1),
231 _bufStrat(Extend),
232 _shift1(0),
233 _shift2((convVar.getMax("cache") + convVar.getMin("cache")) / 2),
234 _cacheObs("!cacheObs", "Cached observables", this, false, false)
235{
236 prepareFFTBinning(convVar);
237
238 calcParams() ;
239
240}
241
242////////////////////////////////////////////////////////////////////////////////
243/// \copydoc RooFFTConvPdf(const char*, const char*, RooRealVar&, RooAbsPdf&, RooAbsPdf&, Int_t)
244/// \param[in] pdfConvVar If the variable used for convolution is a PDF, itself, pass the PDF here, and pass the convolution variable to
245/// `convVar`. See also rf210_angularconv.C in the <a href="https://root.cern/root/html/tutorials/roofit/index.html.">roofit tutorials</a>
246
247RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooAbsReal &pdfConvVar, RooRealVar &convVar,
248 RooAbsPdf &pdf1, RooAbsPdf &pdf2, Int_t ipOrder)
249 : RooAbsCachedPdf(name, title, ipOrder),
250 _x("!x", "Convolution Variable", this, convVar, false, false),
251 _xprime("!xprime", "External Convolution Variable", this, pdfConvVar),
252 _pdf1("!pdf1", "pdf1", this, pdf1, false),
253 _pdf2("!pdf2", "pdf2", this, pdf2, false),
254 _params("!params", "effective parameters", this),
255 _bufFrac(0.1),
256 _bufStrat(Extend),
257 _shift1(0),
258 _shift2((convVar.getMax("cache") + convVar.getMin("cache")) / 2),
259 _cacheObs("!cacheObs", "Cached observables", this, false, false)
260{
261 prepareFFTBinning(convVar);
262
263 calcParams() ;
264}
265
266
267
268////////////////////////////////////////////////////////////////////////////////
269/// Copy constructor
270
273 _x("!x",this,other._x),
274 _xprime("!xprime",this,other._xprime),
275 _pdf1("!pdf1",this,other._pdf1),
276 _pdf2("!pdf2",this,other._pdf2),
277 _params("!params",this,other._params),
278 _bufFrac(other._bufFrac),
279 _bufStrat(other._bufStrat),
280 _shift1(other._shift1),
281 _shift2(other._shift2),
282 _cacheObs("!cacheObs",this,other._cacheObs)
283 {
284 }
285
286
287
288////////////////////////////////////////////////////////////////////////////////
289/// Destructor
290
294
295
296////////////////////////////////////////////////////////////////////////////////
297/// Try to improve the binning and inform user if possible.
298/// With a 10% buffer fraction, 930 raw bins yield 1024 FFT bins,
299/// a sweet spot for the speed of FFTW.
300
302 if (!convVar.hasBinning("cache")) {
303 const RooAbsBinning& varBinning = convVar.getBinning();
304 const int optimal = static_cast<Int_t>(1024/(1.+_bufFrac));
305
306 //Can improve precision if binning is uniform
307 if (varBinning.numBins() < optimal && varBinning.isUniform()) {
308 coutI(Caching) << "Changing internal binning of variable '" << convVar.GetName()
309 << "' in FFT '" << fName << "'"
310 << " from " << varBinning.numBins()
311 << " to " << optimal << " to improve the precision of the numerical FFT."
312 << " This can be done manually by setting an additional binning named 'cache'." << std::endl;
313 convVar.setBinning(RooUniformBinning(varBinning.lowBound(), varBinning.highBound(), optimal, "cache"), "cache");
314 } else {
315 coutE(Caching) << "The internal binning of variable " << convVar.GetName()
316 << " is not uniform. The numerical FFT will likely yield wrong results." << std::endl;
317 convVar.setBinning(varBinning, "cache");
318 }
319 }
320}
321
322
323////////////////////////////////////////////////////////////////////////////////
324/// Return base name component for cache components in this case 'PDF1_CONV_PDF2'
325
327{
328 static TString name ;
329 name = _pdf1.arg().GetName() ;
330 name.Append("_CONV_") ;
331 name.Append(_pdf2.arg().GetName()) ;
332 return name.Data() ;
333}
334
335
336
337
338////////////////////////////////////////////////////////////////////////////////
339/// Return specialized cache subclass for FFT calculations
340
342{
343 return new FFTCacheElem(*this,nset) ;
344}
345
346
347
348
349////////////////////////////////////////////////////////////////////////////////
350/// Clone input pdf and attach to dataset
351
354{
355 RooAbsPdf* clonePdf1 = static_cast<RooAbsPdf*>(self._pdf1.arg().cloneTree()) ;
356 RooAbsPdf* clonePdf2 = static_cast<RooAbsPdf*>(self._pdf2.arg().cloneTree()) ;
357 clonePdf1->attachDataSet(*hist()) ;
358 clonePdf2->attachDataSet(*hist()) ;
359
360 // Shift observable
361 RooRealVar* convObs = static_cast<RooRealVar*>(hist()->get()->find(self._x.arg().GetName())) ;
362
363 // Install FFT reference range
364 string refName = Form("refrange_fft_%s",self.GetName()) ;
365 convObs->setRange(refName.c_str(),convObs->getMin(),convObs->getMax()) ;
366
367 if (self._shift1!=0) {
368 RooLinearVar* shiftObs1 = new RooLinearVar(Form("%s_shifted_FFTBuffer1",convObs->GetName()),"shiftObs1",
370
373 cust.replaceArg(*convObs,*shiftObs1) ;
374
375 pdf1Clone.reset(static_cast<RooAbsPdf*>(cust.build()));
376
377 pdf1Clone->addOwnedComponents(*shiftObs1) ;
378 pdf1Clone->addOwnedComponents(*clonePdf1) ;
379
380 } else {
381 pdf1Clone.reset(clonePdf1) ;
382 }
383
384 if (self._shift2!=0) {
385 RooLinearVar* shiftObs2 = new RooLinearVar(Form("%s_shifted_FFTBuffer2",convObs->GetName()),"shiftObs2",
387
390 cust.replaceArg(*convObs,*shiftObs2) ;
391
392 pdf1Clone->addOwnedComponents(*shiftObs2) ;
393 pdf1Clone->addOwnedComponents(*clonePdf2) ;
394
395 pdf2Clone.reset(static_cast<RooAbsPdf*>(cust.build()));
396
397 } else {
398 pdf2Clone.reset(clonePdf2) ;
399 }
400
401
402 // Attach cloned pdf to all original parameters of self
405 self.getParameters(&convObsSet, fftParams) ;
406
407 // Remove all cache histogram from fftParams as these
408 // observable need to remain attached to the histogram
409 fftParams.remove(*hist()->get(),true,true) ;
410
411 pdf1Clone->recursiveRedirectServers(fftParams) ;
412 pdf2Clone->recursiveRedirectServers(fftParams) ;
413 pdf1Clone->fixAddCoefRange(refName.c_str(), true) ;
414 pdf2Clone->fixAddCoefRange(refName.c_str(), true) ;
415
416 // Ensure that coefficients for Add PDFs are only interpreted with respect to the convolution observable
417 RooArgSet convSet(self._x.arg());
418 pdf1Clone->fixAddCoefNormalization(convSet, true);
419 pdf2Clone->fixAddCoefNormalization(convSet, true);
420
421 // Save copy of original histX binning and make alternate binning
422 // for extended range scanning
423
424 const Int_t N = convObs->numBins();
425 if (N < 900) {
426 oocoutW(&self, Eval) << "The FFT convolution '" << self.GetName() << "' will run with " << N
427 << " bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number"
428 " of bins of the observable '" << convObs->GetName() << "'." << std::endl;
429 }
430 Int_t Nbuf = static_cast<Int_t>((N*self.bufferFraction())/2 + 0.5) ;
431 double obw = (convObs->getMax() - convObs->getMin())/N ;
432 Int_t N2 = N+2*Nbuf ;
433
434 scanBinning = std::make_unique<RooUniformBinning>(convObs->getMin()-Nbuf*obw,convObs->getMax()+Nbuf*obw,N2);
435 histBinning.reset(convObs->getBinning().clone());
436
437 // Deactivate dirty state propagation on datahist observables
438 // and set all nodes on both pdfs to operMode AlwaysDirty
439 hist()->setDirtyProp(false) ;
440 convObs->setOperMode(ADirty,true) ;
441}
442
443
444////////////////////////////////////////////////////////////////////////////////
445/// Suffix for cache histogram (added in addition to suffix for cache name)
446
448{
449 return TString(Form("_BufFrac%3.1f_BufStrat%d",_bufFrac,_bufStrat)) ;
450}
451
452
453
454////////////////////////////////////////////////////////////////////////////////
455/// Returns all RooAbsArg objects contained in the cache element
456
458{
459 RooArgList ret(PdfCacheElem::containedArgs(a)) ;
460
461 ret.add(*pdf1Clone) ;
462 ret.add(*pdf2Clone) ;
463 if (pdf1Clone->ownedComponents()) {
464 ret.add(*pdf1Clone->ownedComponents()) ;
465 }
466 if (pdf2Clone->ownedComponents()) {
467 ret.add(*pdf2Clone->ownedComponents()) ;
468 }
469
470 return ret ;
471}
472
473
474
475////////////////////////////////////////////////////////////////////////////////
476/// Fill the contents of the cache the FFT convolution output
477
479{
480 RooDataHist& cacheHist = *cache.hist() ;
481
482 (static_cast<FFTCacheElem&>(cache)).pdf1Clone->setOperMode(ADirty,true) ;
483 (static_cast<FFTCacheElem&>(cache)).pdf2Clone->setOperMode(ADirty,true) ;
484
485 // Determine if there other observables than the convolution observable in the cache
488
489 RooAbsArg* histArg = otherObs.find(_x.arg().GetName()) ;
490 if (histArg) {
491 otherObs.remove(*histArg,true,true) ;
492 }
493
494 //cout << "RooFFTConvPdf::fillCacheObject() otherObs = " << otherObs << std::endl ;
495
496 // Handle trivial scenario -- no other observables
497 if (otherObs.empty()) {
498 fillCacheSlice(static_cast<FFTCacheElem&>(cache),RooArgSet()) ;
499 return ;
500 }
501
502 // Handle cases where there are other cache slices
503 // Iterator over available slice positions and fill each
504
505 // Determine number of bins for each slice position observable
506 Int_t n = otherObs.size() ;
507 std::vector<int> binCur(n+1);
508 std::vector<int> binMax(n+1);
509 Int_t curObs = 0 ;
510
511 std::vector<RooAbsLValue*> obsLV(n);
512 Int_t i(0) ;
513 for (auto const& arg : otherObs) {
514 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(arg) ;
515 obsLV[i] = lvarg ;
516 binCur[i] = 0 ;
517 // coverity[FORWARD_NULL]
518 binMax[i] = lvarg->numBins(binningName())-1 ;
519 i++ ;
520 }
521
522 bool loop(true) ;
523 while(loop) {
524 // Set current slice position
525 for (Int_t j=0 ; j<n ; j++) { obsLV[j]->setBin(binCur[j],binningName()) ; }
526
527// std::cout << "filling slice: bin of obsLV[0] = " << obsLV[0]->getBin() << std::endl ;
528
529 // Fill current slice
530 fillCacheSlice(static_cast<FFTCacheElem&>(cache),otherObs) ;
531
532 // Determine which iterator to increment
533 while(binCur[curObs]==binMax[curObs]) {
534
535 // Reset current iterator and consider next iterator ;
536 binCur[curObs]=0 ;
537 curObs++ ;
538
539 // master termination condition
540 if (curObs==n) {
541 loop=false ;
542 break ;
543 }
544 }
545
546 // Increment current iterator
547 binCur[curObs]++ ;
548 curObs=0 ;
549
550 }
551
552}
553
554
555////////////////////////////////////////////////////////////////////////////////
556/// Fill a slice of cachePdf with the output of the FFT convolution calculation
557
559{
560 // Extract histogram that is the basis of the RooHistPdf
561 RooDataHist& cacheHist = *aux.hist() ;
562
563 // Sample array of input points from both pdfs
564 // Note that returned arrays have optional buffers zones below and above range ends
565 // to reduce cyclical effects and have been cyclically rotated so that bin containing
566 // zero value is at position zero. Example:
567 //
568 // original: -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5
569 // add buffer zones: U U -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5 O O
570 // rotate: 0 +1 +2 +3 +4 +5 O O U U -5 -4 -3 -2 -1
571 //
572 //
573
574 Int_t N;
575 Int_t N2;
578
579 RooRealVar* histX = static_cast<RooRealVar*>(cacheHist.get()->find(_x.arg().GetName())) ;
580 if (_bufStrat==Extend) histX->setBinning(*aux.scanBinning) ;
581 std::vector<double> input1 = scanPdf(const_cast<RooRealVar &>(static_cast<RooRealVar const&>(_x.arg())),*aux.pdf1Clone,cacheHist,slicePos,N,N2,binShift1,_shift1) ;
582 std::vector<double> input2 = scanPdf(const_cast<RooRealVar &>(static_cast<RooRealVar const&>(_x.arg())),*aux.pdf2Clone,cacheHist,slicePos,N,N2,binShift2,_shift2) ;
583 if (_bufStrat==Extend) histX->setBinning(*aux.histBinning) ;
584
585#ifndef ROOFIT_MATH_FFTW3
586 // If ROOT was NOT built with the fftw3 interface, we try to include fftw3.h
587 // with the interpreter and run the concolution in the interpreter.
588 std::vector<double> output(N2);
589
590 auto doFFT = declareDoFFT();
591 doFFT(N2, input1.data(), input2.data(), output.data());
592#else
593 // If ROOT was built with the fftw3 interface, we can use it as a TVirtualFFT
594 // plugin. The advantage here is that nothing can go wrong if fftw3.h wahs
595 // not istalled by the user separately.
596
597 // Retrieve previously defined FFT transformation plans
598 if (!aux.fftr2c1) {
599 aux.fftr2c1.reset(TVirtualFFT::FFT(1, &N2, "R2CK"));
600 aux.fftr2c2.reset(TVirtualFFT::FFT(1, &N2, "R2CK"));
601 aux.fftc2r.reset(TVirtualFFT::FFT(1, &N2, "C2RK"));
602
603 if (aux.fftr2c1 == nullptr || aux.fftr2c2 == nullptr || aux.fftc2r == nullptr) {
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.");
606 }
607 }
608
609 // Real->Complex FFT Transform on p.d.f. 1 sampling
610 aux.fftr2c1->SetPoints(input1.data());
611 aux.fftr2c1->Transform();
612
613 // Real->Complex FFT Transform on p.d.f 2 sampling
614 aux.fftr2c2->SetPoints(input2.data());
615 aux.fftr2c2->Transform();
616
617 // Loop over first half +1 of complex output results, multiply
618 // and set as input of reverse transform
619 for (Int_t i=0 ; i<N2/2+1 ; i++) {
620 double re1;
621 double re2;
622 double im1;
623 double im2;
624 aux.fftr2c1->GetPointComplex(i,re1,im1) ;
625 aux.fftr2c2->GetPointComplex(i,re2,im2) ;
626 double re = re1*re2 - im1*im2 ;
627 double im = re1*im2 + re2*im1 ;
628 TComplex t(re,im) ;
629 aux.fftc2r->SetPointComplex(i,t) ;
630 }
631
632 // Reverse Complex->Real FFT transform product
633 aux.fftc2r->Transform() ;
634#endif
635
636 Int_t totalShift = binShift1 + (N2-N)/2 ;
637
638 // Store FFT result in cache
639
640 std::unique_ptr<TIterator> iter{const_cast<RooDataHist&>(cacheHist).sliceIterator(const_cast<RooAbsReal&>(_x.arg()),slicePos)};
641 for (Int_t i =0 ; i<N ; i++) {
642
643 // Cyclically shift array back so that bin containing zero is back in zeroBin
644 Int_t j = i + totalShift ;
645 while (j<0) j+= N2 ;
646 while (j>=N2) j-= N2 ;
647
648 iter->Next() ;
649#ifndef ROOFIT_MATH_FFTW3
650 cacheHist.set(output[j]);
651#else
652 cacheHist.set(aux.fftc2r->GetPointReal(j));
653#endif
654 }
655}
656
657
658////////////////////////////////////////////////////////////////////////////////
659/// Scan the values of 'pdf' in observable 'obs' using the bin values stored in 'hist' at slice position 'slicePos'
660/// N is filled with the number of bins defined in hist, N2 is filled with N plus the number of buffer bins
661/// The return value is an array of doubles of length N2 with the sampled values. The caller takes ownership
662/// of the array
663
664std::vector<double> RooFFTConvPdf::scanPdf(RooRealVar& obs, RooAbsPdf& pdf, const RooDataHist& hist, const RooArgSet& slicePos,
665 Int_t& N, Int_t& N2, Int_t& zeroBin, double shift) const
666{
667
668 RooRealVar* histX = static_cast<RooRealVar*>(hist.get()->find(obs.GetName())) ;
669
670 // Calculate number of buffer bins on each size to avoid cyclical flow
671 N = histX->numBins(binningName()) ;
672 Int_t Nbuf = static_cast<Int_t>((N*bufferFraction())/2 + 0.5) ;
673 N2 = N+2*Nbuf ;
674
675
676 // Allocate array of sampling size plus optional buffer zones
677 std::vector<double> array(N2);
678
679 // Set position of non-convolution observable to that of the cache slice that were are processing now
680 hist.get(slicePos) ;
681
682 // Find bin ID that contains zero value
683 zeroBin = 0 ;
684 if (histX->getMax()>=0 && histX->getMin()<=0) {
685 zeroBin = histX->getBinning().binNumber(0) ;
686 } else if (histX->getMin()>0) {
687 double bw = (histX->getMax() - histX->getMin())/N2 ;
688 zeroBin = Int_t(-histX->getMin()/bw) ;
689 } else {
690 double bw = (histX->getMax() - histX->getMin())/N2 ;
691 zeroBin = Int_t(-1*histX->getMax()/bw) ;
692 }
693
694 Int_t binShift = Int_t((N2* shift) / (histX->getMax()-histX->getMin())) ;
695
696 zeroBin += binShift ;
697 while(zeroBin>=N2) zeroBin-= N2 ;
698 while(zeroBin<0) zeroBin+= N2 ;
699
700 // First scan hist into temp array
701 std::vector<double> tmp(N2);
702 Int_t k(0) ;
703 switch(_bufStrat) {
704
705 case Extend:
706 // Sample entire extended range (N2 samples)
707 for (k=0 ; k<N2 ; k++) {
708 histX->setBin(k) ;
709 tmp[k] = pdf.getVal(hist.get()) ;
710 }
711 break ;
712
713 case Flat:
714 // Sample original range (N samples) and fill lower and upper buffer
715 // bins with p.d.f. value at respective boundary
716 {
717 histX->setBin(0) ;
718 double val = pdf.getVal(hist.get()) ;
719 for (k=0 ; k<Nbuf ; k++) {
720 tmp[k] = val ;
721 }
722 for (k=0 ; k<N ; k++) {
723 histX->setBin(k) ;
724 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
725 }
726 histX->setBin(N-1) ;
727 val = pdf.getVal(hist.get()) ;
728 for (k=0 ; k<Nbuf ; k++) {
729 tmp[N+Nbuf+k] = val ;
730 }
731 }
732 break ;
733
734 case Mirror:
735 // Sample original range (N samples) and fill lower and upper buffer
736 // bins with mirror image of sampled range
737 for (k=0 ; k<N ; k++) {
738 histX->setBin(k) ;
739 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
740 }
741 for (k=1 ; k<=Nbuf ; k++) {
742 histX->setBin(k) ;
743 tmp[Nbuf-k] = pdf.getVal(hist.get()) ;
744 histX->setBin(N-k) ;
745 tmp[Nbuf+N+k-1] = pdf.getVal(hist.get()) ;
746 }
747 break ;
748 }
749
750 // Scan function and store values in array
751 for (Int_t i=0 ; i<N2 ; i++) {
752 // Cyclically shift writing location by zero bin position
753 Int_t j = i - (zeroBin) ;
754 if (j<0) j+= N2 ;
755 if (j>=N2) j-= N2 ;
756 array[i] = tmp[j] ;
757 }
758
759 return array ;
760}
761
762
763
764////////////////////////////////////////////////////////////////////////////////
765/// Return the observables to be cached given the normalization set nset.
766///
767/// If the cache observable is in nset then this is
768/// - the convolution observable plus
769/// - any member of nset that is either a RooCategory,
770/// - or was previously specified through setCacheObservables().
771///
772/// In case the cache observable is *not* in nset, then it is
773/// - the convolution observable plus
774/// - all member of nset that are observables of this p.d.f.
775///
776
778{
779 // Get complete list of observables
780 auto obs1 = new RooArgSet{};
782 _pdf1.arg().getObservables(&nset, *obs1) ;
783 _pdf2.arg().getObservables(&nset, obs2) ;
784 obs1->add(obs2, true) ;
785
786 // Check if convolution observable is in nset
787 if (nset.contains(_x.arg())) {
788
789 // Now strip out all non-category observables
791 for(auto const& arg : *obs1) {
792 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
793 killList.add(*arg) ;
794 }
795 }
796 obs1->remove(killList) ;
797
798 // And add back the convolution observables
799 obs1->add(_x.arg(),true) ;
800
801 obs1->add(_cacheObs) ;
802
803 } else {
804
805 // If cacheObs was filled, cache only observables in there
806 if (!_cacheObs.empty()) {
808 for(auto const& arg : *obs1) {
809 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
810 killList.add(*arg) ;
811 }
812 }
813 obs1->remove(killList) ;
814 }
815
816
817 // Make sure convolution observable is always in there
818 obs1->add(_x.arg(),true) ;
819
820 }
821
823}
824
825
826
827////////////////////////////////////////////////////////////////////////////////
828/// Return the parameters on which the cache depends given normalization
829/// set nset. For this p.d.f these are the parameters of the input p.d.f.
830/// but never the convolution variable, in case it is not part of nset.
831
833{
834 std::unique_ptr<RooArgSet> vars{getVariables()};
835 vars->remove(*std::unique_ptr<RooArgSet>{actualObservables(nset)});
836
837 return RooFit::makeOwningPtr(std::move(vars));
838}
839
840
841
842////////////////////////////////////////////////////////////////////////////////
843/// Return p.d.f. observable (which can be a function) to substitute given
844/// p.d.f. observable. Substitutes x by xprime if xprime is set.
845
847{
848 if (_xprime.absArg() && string(histObservable.GetName())==_x.absArg()->GetName()) {
849 return (*_xprime.absArg()) ;
850 }
851 return histObservable ;
852}
853
854
855
856////////////////////////////////////////////////////////////////////////////////
857/// Create appropriate generator context for this convolution. If both input p.d.f.s support
858/// internal generation, if it is safe to use them and if no observables other than the convolution
859/// observable are requested for generation, use the specialized convolution generator context
860/// which implements a smearing strategy in the convolution observable. If not return the
861/// regular accept/reject generator context
862
864 const RooArgSet* auxProto, bool verbose) const
865{
866 RooArgSet vars2(vars) ;
867 vars2.remove(_x.arg(),true,true) ;
868 Int_t numAddDep = vars2.size() ;
869
870 RooArgSet dummy ;
871 bool pdfCanDir = ((static_cast<RooAbsPdf const &>(_pdf1.arg())).getGenerator(_x.arg(), dummy) != 0 &&
872 (static_cast<RooAbsPdf const &>(_pdf1.arg())).isDirectGenSafe(_x.arg()));
873 bool resCanDir = ((static_cast<RooAbsPdf const &>(_pdf2.arg())).getGenerator(_x.arg(), dummy) != 0 &&
874 (static_cast<RooAbsPdf const &>(_pdf2.arg())).isDirectGenSafe(_x.arg()));
875
876 if (pdfCanDir) {
877 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f " << _pdf1.arg().GetName()
878 << " has internal generator that is safe to use in current context" << std::endl ;
879 }
880 if (resCanDir) {
881 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f. " << _pdf2.arg().GetName()
882 << " has internal generator that is safe to use in current context" << std::endl ;
883 }
884 if (numAddDep>0) {
885 cxcoutI(Generation) << "RooFFTConvPdf::genContext() generation requested for observables other than the convolution observable " << _x.arg().GetName() << std::endl ;
886 }
887
888
889 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
890 // Any resolution model with more dependents than the convolution variable
891 // or pdf or resmodel do not support direct generation
892 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting accept/reject generator context because one or both of the input "
893 << "p.d.f.s cannot use internal generator and/or "
894 << "observables other than the convolution variable are requested for generation" << std::endl ;
895 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
896 }
897
898 // Any other resolution model: use specialized generator context
899 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting specialized convolution generator context as both input "
900 << "p.d.fs are safe for internal generator and only "
901 << "the convolution observables is requested for generation" << std::endl ;
902 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
903}
904
905
906
907////////////////////////////////////////////////////////////////////////////////
908/// Change the size of the buffer on either side of the observable range to `frac` times the
909/// size of the range of the convolution observable.
910
912{
913 if (frac<0) {
914 coutE(InputArguments) << "RooFFTConvPdf::setBufferFraction(" << GetName() << ") fraction should be greater than or equal to zero" << std::endl ;
915 return ;
916 }
917 _bufFrac = frac ;
918
919 // Sterilize the cache as certain partial results depend on buffer fraction
921}
922
923
924////////////////////////////////////////////////////////////////////////////////
925/// Change strategy to fill the overflow buffer on either side of the convolution observable range.
926///
927/// - `Extend` means is that the input p.d.f convolution observable range is widened to include the buffer range
928/// - `Flat` means that the buffer is filled with the p.d.f. value at the boundary of the observable range
929/// - `Mirror` means that the buffer is filled with a mirror image of the p.d.f. around the convolution observable boundary
930///
931/// The default strategy is extend. If one of the input p.d.f.s is a RooAddPdf, it is configured so that the interpretation
932/// range of the fraction coefficients is kept at the nominal convolutions observable range (instead of interpreting coefficients
933/// in the widened range including the buffer).
934
939
940
941
942////////////////////////////////////////////////////////////////////////////////
943/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
944/// product operator construction
945
946void RooFFTConvPdf::printMetaArgs(ostream& os) const
947{
948 os << _pdf1.arg().GetName() << "(" << _x.arg().GetName() << ") (*) " << _pdf2.arg().GetName() << "(" << _x.arg().GetName() << ") " ;
949}
950
951
952
953////////////////////////////////////////////////////////////////////////////////
954/// (Re)calculate effective parameters of this p.d.f.
955
957{
958 RooArgSet argSet{_x.arg()};
961 _pdf1.arg().getParameters(&argSet, params1);
962 _pdf2.arg().getParameters(&argSet, params2);
965 _params.add(params2,true) ;
966}
#define a(i)
Definition RSha256.hxx:99
#define coutI(a)
#define cxcoutI(a)
#define oocoutW(o, a)
#define coutF(a)
#define oocoutE(o, a)
#define coutE(a)
int Int_t
Definition RtypesCore.h:45
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define N
char name[80]
Definition TGX11.cxx:110
#define gInterpreter
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
Abstract base class for RooRealVar binning definitions.
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
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
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.
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
static TClass * Class()
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * absArg() const
Return pointer to contained argument.
Definition RooArgProxy.h:46
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
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 ...
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
Container class to hold unbinned data.
Definition RooDataSet.h:34
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
PDF for the numerical (FFT) convolution of two PDFs.
friend class RooConvGenContext
RooSetProxy _params
Effective parameters of this p.d.f.
BufStrat _bufStrat
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.
friend class FFTCacheElem
RooSetProxy _cacheObs
Non-convolution observables that are also cached.
Implements a universal generator context for all RooAbsPdf classes that do not have or need a special...
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.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
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.
const T & arg() const
Return reference to object held in proxy.
Implementation of RooAbsBinning that provides a uniform binning in 'n' bins between the range end poi...
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
TString fName
Definition TNamed.h:32
Basic string class.
Definition TString.h:139
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)
const Int_t n
Definition legend1.C:16
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:35
OwningPtr< T > makeOwningPtr(std::unique_ptr< T > &&ptr)
Internal helper to turn a std::unique_ptr<T> into an OwningPtr.
Definition Config.h:40
static void output()