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#include "RooFitImplHelpers.h"
127
128#include "TClass.h"
129#include "TComplex.h"
130#include "TVirtualFFT.h"
131
132#include <iostream>
133#include <stdexcept>
134
135#ifndef ROOFIT_MATH_FFTW3
136#include "TInterpreter.h"
137
138namespace {
139
140auto declareDoFFT()
141{
142 static void (*doFFT)(int, double *, double *, double *) = nullptr;
143
144 if (doFFT)
145 return doFFT;
146
147 bool success = gInterpreter->Declare("#include \"fftw3.h\"");
148
149 if (!success) {
150 std::stringstream ss;
151 ss << "RooFFTConvPdf evaluation Failed! The interpreter could not find the fftw3.h header.\n";
152 ss << "You have three options to fix this problem:\n";
153 ss << " 1) Install fftw3 on your system so that the interpreter can find it\n";
154 ss << " 2) In case fftw3.h is installed somewhere else,\n"
155 << " tell ROOT with gInterpreter->AddIncludePath() where to find it\n";
156 ss << " 3) Compile ROOT with the -Dfftw3=ON in the CMake configuration,\n"
157 << " such that ROOT comes with built-in fftw3 convolution routines\n";
158 oocoutE(nullptr, Eval) << ss.str() << std::endl;
159 throw std::runtime_error("RooFFTConvPdf evaluation Failed! The interpreter could not find the fftw3.h header");
160 }
161
162 gInterpreter->Declare(R"(
163void RooFFTConvPdf_doFFT(int n, double *input1, double *input2, double *output)
164{
165 auto fftr2c1_Out = reinterpret_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex) * (n / 2 + 1)));
166 auto fftr2c2_Out = reinterpret_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex) * (n / 2 + 1)));
167 auto fftc2r_In = reinterpret_cast<fftw_complex *>(fftw_malloc(sizeof(fftw_complex) * (n / 2 + 1)));
168
169 fftw_plan fftr2c1_plan = fftw_plan_dft_r2c(1, &n, input1, fftr2c1_Out, FFTW_ESTIMATE);
170 fftw_plan fftr2c2_plan = fftw_plan_dft_r2c(1, &n, input2, fftr2c2_Out, FFTW_ESTIMATE);
171 fftw_plan fftc2r_plan = fftw_plan_dft_c2r(1, &n, fftc2r_In, output, FFTW_ESTIMATE);
172
173 // Real->Complex FFT Transform on p.d.f. samplings
174 fftw_execute(fftr2c1_plan);
175 fftw_execute(fftr2c2_plan);
176
177 // Loop over first half +1 of complex output results, multiply
178 // and set as input of reverse transform
179 for (Int_t i = 0; i < n / 2 + 1; i++) {
180 double re1 = fftr2c1_Out[i][0];
181 double re2 = fftr2c2_Out[i][0];
182 double im1 = fftr2c1_Out[i][1];
183 double im2 = fftr2c2_Out[i][1];
184 fftc2r_In[i][0] = re1 * re2 - im1 * im2;
185 fftc2r_In[i][1] = re1 * im2 + re2 * im1;
186 }
187
188 // Reverse Complex->Real FFT transform product
189 fftw_execute(fftc2r_plan);
190
191 fftw_destroy_plan(fftr2c1_plan);
192 fftw_destroy_plan(fftr2c2_plan);
193 fftw_destroy_plan(fftc2r_plan);
194
195 fftw_free(fftr2c1_Out);
196 fftw_free(fftr2c2_Out);
197 fftw_free(fftc2r_In);
198}
199)");
200
201 doFFT = reinterpret_cast<void(*)(int, double*, double*, double*)>(gInterpreter->ProcessLine("RooFFTConvPdf_doFFT;"));
202 return doFFT;
203}
204
205} // namespace
206
207#endif
208
209using std::endl, std::string, std::ostream;
210
211
212////////////////////////////////////////////////////////////////////////////////
213/// Constructor for numerical (FFT) convolution of PDFs.
214/// \param[in] name Name of this PDF
215/// \param[in] title Title for plotting this PDF
216/// \param[in] convVar Observable to convolve the PDFs in \attention Use a high number of bins (>= 1000) for good accuracy.
217/// \param[in] pdf1 First PDF to be convolved
218/// \param[in] pdf2 Second PDF to be convolved
219/// \param[in] ipOrder Order for interpolation between bins (since FFT is discrete)
220/// The binning used for the FFT sampling is controlled by the binning named "cache" in the convolution observable `convVar`.
221/// If such a binning is not set, the same number of bins as for `convVar` will be used.
222
223RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooRealVar &convVar, RooAbsPdf &pdf1, RooAbsPdf &pdf2,
225 : RooAbsCachedPdf(name, title, ipOrder),
226 _x("!x", "Convolution Variable", this, convVar),
227 _xprime("!xprime", "External Convolution Variable", this, false),
228 _pdf1("!pdf1", "pdf1", this, pdf1, false),
229 _pdf2("!pdf2", "pdf2", this, pdf2, false),
230 _params("!params", "effective parameters", this),
231 _bufFrac(0.1),
232 _bufStrat(Extend),
233 _shift1(0),
234 _shift2((convVar.getMax("cache") + convVar.getMin("cache")) / 2),
235 _cacheObs("!cacheObs", "Cached observables", this, false, false)
236{
237 prepareFFTBinning(convVar);
238
239 calcParams() ;
240
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// \copydoc RooFFTConvPdf(const char*, const char*, RooRealVar&, RooAbsPdf&, RooAbsPdf&, Int_t)
245/// \param[in] pdfConvVar If the variable used for convolution is a function, itself, pass the function here, and pass
246/// the convolution variable to `convVar`. See also rf210_angularconv.C in the <a
247/// href="https://root.cern/root/html/tutorials/roofit/index.html.">roofit tutorials</a>
248
249RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooAbsReal &pdfConvVar, RooRealVar &convVar,
250 RooAbsPdf &pdf1, RooAbsPdf &pdf2, Int_t ipOrder)
251 : RooAbsCachedPdf(name, title, ipOrder),
252 _x("!x", "Convolution Variable", this, convVar, false, false),
253 _xprime("!xprime", "External Convolution Variable", this, pdfConvVar),
254 _pdf1("!pdf1", "pdf1", this, pdf1, false),
255 _pdf2("!pdf2", "pdf2", this, pdf2, false),
256 _params("!params", "effective parameters", this),
257 _bufFrac(0.1),
258 _bufStrat(Extend),
259 _shift1(0),
260 _shift2((convVar.getMax("cache") + convVar.getMin("cache")) / 2),
261 _cacheObs("!cacheObs", "Cached observables", this, false, false)
262{
263 prepareFFTBinning(convVar);
264
265 calcParams() ;
266}
267
268
269
270////////////////////////////////////////////////////////////////////////////////
271/// Copy constructor
272
275 _x("!x",this,other._x),
276 _xprime("!xprime",this,other._xprime),
277 _pdf1("!pdf1",this,other._pdf1),
278 _pdf2("!pdf2",this,other._pdf2),
279 _params("!params",this,other._params),
280 _bufFrac(other._bufFrac),
281 _bufStrat(other._bufStrat),
282 _shift1(other._shift1),
283 _shift2(other._shift2),
284 _cacheObs("!cacheObs",this,other._cacheObs)
285 {
286 }
287
288
289
290////////////////////////////////////////////////////////////////////////////////
291/// Destructor
292
296
297
298////////////////////////////////////////////////////////////////////////////////
299/// Try to improve the binning and inform user if possible.
300/// With a 10% buffer fraction, 930 raw bins yield 1024 FFT bins,
301/// a sweet spot for the speed of FFTW.
302
304 if (!convVar.hasBinning("cache")) {
305 const RooAbsBinning& varBinning = convVar.getBinning();
306 const int optimal = static_cast<Int_t>(1024/(1.+_bufFrac));
307
308 //Can improve precision if binning is uniform
309 if (varBinning.numBins() < optimal && varBinning.isUniform()) {
310 coutI(Caching) << "Changing internal binning of variable '" << convVar.GetName()
311 << "' in FFT '" << fName << "'"
312 << " from " << varBinning.numBins()
313 << " to " << optimal << " to improve the precision of the numerical FFT."
314 << " This can be done manually by setting an additional binning named 'cache'." << std::endl;
315 convVar.setBinning(RooUniformBinning(varBinning.lowBound(), varBinning.highBound(), optimal, "cache"), "cache");
316 } else {
317 coutE(Caching) << "The internal binning of variable " << convVar.GetName()
318 << " is not uniform. The numerical FFT will likely yield wrong results." << std::endl;
319 convVar.setBinning(varBinning, "cache");
320 }
321 }
322}
323
324
325////////////////////////////////////////////////////////////////////////////////
326/// Return base name component for cache components in this case 'PDF1_CONV_PDF2'
327
329{
330 static TString name ;
331 name = _pdf1.arg().GetName() ;
332 name.Append("_CONV_") ;
333 name.Append(_pdf2.arg().GetName()) ;
334 return name.Data() ;
335}
336
337
338
339
340////////////////////////////////////////////////////////////////////////////////
341/// Return specialized cache subclass for FFT calculations
342
344{
345 return new FFTCacheElem(*this,nset) ;
346}
347
348
349
350
351////////////////////////////////////////////////////////////////////////////////
352/// Clone input pdf and attach to dataset
353
356{
357 RooAbsPdf* clonePdf1 = static_cast<RooAbsPdf*>(self._pdf1.arg().cloneTree()) ;
358 RooAbsPdf* clonePdf2 = static_cast<RooAbsPdf*>(self._pdf2.arg().cloneTree()) ;
359 clonePdf1->attachDataSet(*hist()) ;
360 clonePdf2->attachDataSet(*hist()) ;
361
362 // Shift observable
363 RooRealVar* convObs = static_cast<RooRealVar*>(hist()->get()->find(self._x.arg().GetName())) ;
364
365 // Install FFT reference range
366 string refName = Form("refrange_fft_%s",self.GetName()) ;
367 convObs->setRange(refName.c_str(),convObs->getMin(),convObs->getMax()) ;
368
369 if (self._shift1!=0) {
370 RooLinearVar* shiftObs1 = new RooLinearVar(Form("%s_shifted_FFTBuffer1",convObs->GetName()),"shiftObs1",
372
375 cust.replaceArg(*convObs,*shiftObs1) ;
376
377 pdf1Clone.reset(static_cast<RooAbsPdf*>(cust.build()));
378
379 pdf1Clone->addOwnedComponents(*shiftObs1) ;
380 pdf1Clone->addOwnedComponents(*clonePdf1) ;
381
382 } else {
383 pdf1Clone.reset(clonePdf1) ;
384 }
385
386 if (self._shift2!=0) {
387 RooLinearVar* shiftObs2 = new RooLinearVar(Form("%s_shifted_FFTBuffer2",convObs->GetName()),"shiftObs2",
389
392 cust.replaceArg(*convObs,*shiftObs2) ;
393
394 pdf1Clone->addOwnedComponents(*shiftObs2) ;
395 pdf1Clone->addOwnedComponents(*clonePdf2) ;
396
397 pdf2Clone.reset(static_cast<RooAbsPdf*>(cust.build()));
398
399 } else {
400 pdf2Clone.reset(clonePdf2) ;
401 }
402
403 // Cache normalization integral values, since we know they don't change for
404 // the given normalization set in this cache object. When using this cache,
405 // we evaluate the pdfs without normalization set and then do the
406 // normalization manually using these cached values. This has less overhead
407 // compared to letting RooAbsPdf::getVal(normSet) figure out if the normSet
408 // has changed and get the caching right.
409 normVal1 = pdf1Clone->getNorm(hist()->get());
410 normVal2 = pdf2Clone->getNorm(hist()->get());
411
412 // Attach cloned pdf to all original parameters of self
415 self.getParameters(&convObsSet, fftParams) ;
416
417 // Remove all cache histogram from fftParams as these
418 // observable need to remain attached to the histogram
419 fftParams.remove(*hist()->get(),true,true) ;
420
421 pdf1Clone->recursiveRedirectServers(fftParams) ;
422 pdf2Clone->recursiveRedirectServers(fftParams) ;
423 pdf1Clone->fixAddCoefRange(refName.c_str(), true) ;
424 pdf2Clone->fixAddCoefRange(refName.c_str(), true) ;
425
426 // Ensure that coefficients for Add PDFs are only interpreted with respect to the convolution observable
427 RooArgSet convSet(self._x.arg());
428 pdf1Clone->fixAddCoefNormalization(convSet, true);
429 pdf2Clone->fixAddCoefNormalization(convSet, true);
430
431 // Save copy of original histX binning and make alternate binning
432 // for extended range scanning
433
434 const Int_t N = convObs->numBins();
435 if (N < 900) {
436 oocoutW(&self, Eval) << "The FFT convolution '" << self.GetName() << "' will run with " << N
437 << " bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number"
438 " of bins of the observable '" << convObs->GetName() << "'." << std::endl;
439 }
440 Int_t Nbuf = static_cast<Int_t>((N*self.bufferFraction())/2 + 0.5) ;
441 double obw = (convObs->getMax() - convObs->getMin())/N ;
442 Int_t N2 = N+2*Nbuf ;
443
444 scanBinning = std::make_unique<RooUniformBinning>(convObs->getMin()-Nbuf*obw,convObs->getMax()+Nbuf*obw,N2);
445 histBinning.reset(convObs->getBinning().clone());
446
447 // Deactivate dirty state propagation on datahist observables
448 // and set all nodes on both pdfs to operMode AlwaysDirty
449 hist()->setDirtyProp(false) ;
450 convObs->setOperMode(ADirty,true) ;
451}
452
453
454////////////////////////////////////////////////////////////////////////////////
455/// Suffix for cache histogram (added in addition to suffix for cache name)
456
458{
459 return TString(Form("_BufFrac%3.1f_BufStrat%d",_bufFrac,_bufStrat)) ;
460}
461
462
463
464////////////////////////////////////////////////////////////////////////////////
465/// Returns all RooAbsArg objects contained in the cache element
466
468{
469 RooArgList ret(PdfCacheElem::containedArgs(a)) ;
470
471 ret.add(*pdf1Clone) ;
472 ret.add(*pdf2Clone) ;
473 if (pdf1Clone->ownedComponents()) {
474 ret.add(*pdf1Clone->ownedComponents()) ;
475 }
476 if (pdf2Clone->ownedComponents()) {
477 ret.add(*pdf2Clone->ownedComponents()) ;
478 }
479
480 return ret ;
481}
482
483
484
485////////////////////////////////////////////////////////////////////////////////
486/// Fill the contents of the cache the FFT convolution output
487
489{
490 RooDataHist& cacheHist = *cache.hist() ;
491
492 (static_cast<FFTCacheElem&>(cache)).pdf1Clone->setOperMode(ADirty,true) ;
493 (static_cast<FFTCacheElem&>(cache)).pdf2Clone->setOperMode(ADirty,true) ;
494
495 // Determine if there other observables than the convolution observable in the cache
498
499 RooAbsArg* histArg = otherObs.find(_x.arg().GetName()) ;
500 if (histArg) {
501 otherObs.remove(*histArg,true,true) ;
502 }
503
504 //cout << "RooFFTConvPdf::fillCacheObject() otherObs = " << otherObs << std::endl ;
505
506 // Handle trivial scenario -- no other observables
507 if (otherObs.empty()) {
508 fillCacheSlice(static_cast<FFTCacheElem&>(cache),RooArgSet()) ;
509 return ;
510 }
511
512 // Handle cases where there are other cache slices
513 // Iterator over available slice positions and fill each
514
515 // Determine number of bins for each slice position observable
516 Int_t n = otherObs.size() ;
517 std::vector<int> binCur(n+1);
518 std::vector<int> binMax(n+1);
519 Int_t curObs = 0 ;
520
521 std::vector<RooAbsLValue*> obsLV(n);
522 Int_t i(0) ;
523 for (auto const& arg : otherObs) {
524 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(arg) ;
525 obsLV[i] = lvarg ;
526 binCur[i] = 0 ;
527 // coverity[FORWARD_NULL]
528 binMax[i] = lvarg->numBins(binningName())-1 ;
529 i++ ;
530 }
531
532 bool loop(true) ;
533 while(loop) {
534 // Set current slice position
535 for (Int_t j=0 ; j<n ; j++) { obsLV[j]->setBin(binCur[j],binningName()) ; }
536
537// std::cout << "filling slice: bin of obsLV[0] = " << obsLV[0]->getBin() << std::endl ;
538
539 // Fill current slice
540 fillCacheSlice(static_cast<FFTCacheElem&>(cache),otherObs) ;
541
542 // Determine which iterator to increment
543 while(binCur[curObs]==binMax[curObs]) {
544
545 // Reset current iterator and consider next iterator ;
546 binCur[curObs]=0 ;
547 curObs++ ;
548
549 // master termination condition
550 if (curObs==n) {
551 loop=false ;
552 break ;
553 }
554 }
555
556 // Increment current iterator
557 binCur[curObs]++ ;
558 curObs=0 ;
559
560 }
561
562}
563
564
565////////////////////////////////////////////////////////////////////////////////
566/// Fill a slice of cachePdf with the output of the FFT convolution calculation
567
569{
570 // Extract histogram that is the basis of the RooHistPdf
571 RooDataHist& cacheHist = *aux.hist() ;
572
573 // Sample array of input points from both pdfs
574 // Note that returned arrays have optional buffers zones below and above range ends
575 // to reduce cyclical effects and have been cyclically rotated so that bin containing
576 // zero value is at position zero. Example:
577 //
578 // original: -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5
579 // add buffer zones: U U -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5 O O
580 // rotate: 0 +1 +2 +3 +4 +5 O O U U -5 -4 -3 -2 -1
581 //
582 //
583
584 Int_t N;
585 Int_t N2;
588
589 RooRealVar* histX = static_cast<RooRealVar*>(cacheHist.get()->find(_x.arg().GetName())) ;
590 if (_bufStrat==Extend) histX->setBinning(*aux.scanBinning) ;
591 RooRealVar &xVar = const_cast<RooRealVar &>(static_cast<RooRealVar const&>(_x.arg()));
592 std::vector<double> input1 = scanPdf(xVar,*aux.pdf1Clone,aux.normVal1,cacheHist,slicePos,N,N2,binShift1,_shift1) ;
593 std::vector<double> input2 = scanPdf(xVar,*aux.pdf2Clone,aux.normVal2,cacheHist,slicePos,N,N2,binShift2,_shift2) ;
594 if (_bufStrat==Extend) histX->setBinning(*aux.histBinning) ;
595
596#ifndef ROOFIT_MATH_FFTW3
597 // If ROOT was NOT built with the fftw3 interface, we try to include fftw3.h
598 // with the interpreter and run the concolution in the interpreter.
599 std::vector<double> output(N2);
600
601 auto doFFT = declareDoFFT();
602 doFFT(N2, input1.data(), input2.data(), output.data());
603#else
604 // If ROOT was built with the fftw3 interface, we can use it as a TVirtualFFT
605 // plugin. The advantage here is that nothing can go wrong if fftw3.h wahs
606 // not istalled by the user separately.
607
608 // Retrieve previously defined FFT transformation plans
609 if (!aux.fftr2c1) {
610 aux.fftr2c1.reset(TVirtualFFT::FFT(1, &N2, "R2CK"));
611 aux.fftr2c2.reset(TVirtualFFT::FFT(1, &N2, "R2CK"));
612 aux.fftc2r.reset(TVirtualFFT::FFT(1, &N2, "C2RK"));
613
614 if (aux.fftr2c1 == nullptr || aux.fftr2c2 == nullptr || aux.fftc2r == nullptr) {
615 coutF(Eval) << "RooFFTConvPdf::fillCacheSlice(" << GetName() << "Cannot get a handle to fftw. Maybe ROOT was built without it?" << std::endl;
616 throw std::runtime_error("Cannot get a handle to fftw.");
617 }
618 }
619
620 // Real->Complex FFT Transform on p.d.f. 1 sampling
621 aux.fftr2c1->SetPoints(input1.data());
622 aux.fftr2c1->Transform();
623
624 // Real->Complex FFT Transform on p.d.f 2 sampling
625 aux.fftr2c2->SetPoints(input2.data());
626 aux.fftr2c2->Transform();
627
628 // Loop over first half +1 of complex output results, multiply
629 // and set as input of reverse transform
630 for (Int_t i=0 ; i<N2/2+1 ; i++) {
631 double re1;
632 double re2;
633 double im1;
634 double im2;
635 aux.fftr2c1->GetPointComplex(i,re1,im1) ;
636 aux.fftr2c2->GetPointComplex(i,re2,im2) ;
637 double re = re1*re2 - im1*im2 ;
638 double im = re1*im2 + re2*im1 ;
639 TComplex t(re,im) ;
640 aux.fftc2r->SetPointComplex(i,t) ;
641 }
642
643 // Reverse Complex->Real FFT transform product
644 aux.fftc2r->Transform() ;
645#endif
646
647 Int_t totalShift = binShift1 + (N2-N)/2 ;
648
649 // Store FFT result in cache
650
651 std::unique_ptr<TIterator> iter{const_cast<RooDataHist&>(cacheHist).sliceIterator(const_cast<RooAbsReal&>(_x.arg()),slicePos)};
652 for (Int_t i =0 ; i<N ; i++) {
653
654 // Cyclically shift array back so that bin containing zero is back in zeroBin
655 Int_t j = i + totalShift ;
656 while (j<0) j+= N2 ;
657 while (j>=N2) j-= N2 ;
658
659 iter->Next() ;
660#ifndef ROOFIT_MATH_FFTW3
661 cacheHist.set(output[j]);
662#else
663 cacheHist.set(aux.fftc2r->GetPointReal(j));
664#endif
665 }
666}
667
668
669////////////////////////////////////////////////////////////////////////////////
670/// Scan the values of 'pdf' in observable 'obs' using the bin values stored in 'hist' at slice position 'slicePos'
671/// N is filled with the number of bins defined in hist, N2 is filled with N plus the number of buffer bins
672/// The return value is an array of doubles of length N2 with the sampled values. The caller takes ownership
673/// of the array
674
675std::vector<double> RooFFTConvPdf::scanPdf(RooRealVar &obs, RooAbsPdf &pdf, double normVal, const RooDataHist &hist,
676 const RooArgSet &slicePos, Int_t &N, Int_t &N2, Int_t &zeroBin,
677 double shift) const
678{
679
680 RooRealVar* histX = static_cast<RooRealVar*>(hist.get()->find(obs.GetName())) ;
681
682 // Calculate number of buffer bins on each size to avoid cyclical flow
683 N = histX->numBins(binningName()) ;
684 Int_t Nbuf = static_cast<Int_t>((N*bufferFraction())/2 + 0.5) ;
685 N2 = N+2*Nbuf ;
686
687
688 // Allocate array of sampling size plus optional buffer zones
689 std::vector<double> array(N2);
690
691 // Set position of non-convolution observable to that of the cache slice that were are processing now
692 hist.get(slicePos) ;
693
694 // Find bin ID that contains zero value
695 zeroBin = 0 ;
696 if (histX->getMax()>=0 && histX->getMin()<=0) {
697 zeroBin = histX->getBinning().binNumber(0) ;
698 } else if (histX->getMin()>0) {
699 double bw = (histX->getMax() - histX->getMin())/N2 ;
700 zeroBin = Int_t(-histX->getMin()/bw) ;
701 } else {
702 double bw = (histX->getMax() - histX->getMin())/N2 ;
703 zeroBin = Int_t(-1*histX->getMax()/bw) ;
704 }
705
706 Int_t binShift = Int_t((N2* shift) / (histX->getMax()-histX->getMin())) ;
707
708 zeroBin += binShift ;
709 while(zeroBin>=N2) zeroBin-= N2 ;
710 while(zeroBin<0) zeroBin+= N2 ;
711
712 // To mimic exactly the normalization code in RooAbsPdf::getValV()
713 auto getPdfVal = [&]() {
714 double rawVal = pdf.getVal();
716 };
717
718 // First scan hist into temp array
719 std::vector<double> tmp(N2);
720 Int_t k(0) ;
721 switch(_bufStrat) {
722
723 case Extend:
724 // Sample entire extended range (N2 samples)
725 for (k=0 ; k<N2 ; k++) {
726 histX->setBin(k) ;
727 tmp[k] = getPdfVal();
728 }
729 break ;
730
731 case Flat:
732 // Sample original range (N samples) and fill lower and upper buffer
733 // bins with p.d.f. value at respective boundary
734 {
735 histX->setBin(0) ;
736 double val = getPdfVal();
737 for (k=0 ; k<Nbuf ; k++) {
738 tmp[k] = val ;
739 }
740 for (k=0 ; k<N ; k++) {
741 histX->setBin(k) ;
742 tmp[k+Nbuf] = getPdfVal();
743 }
744 histX->setBin(N-1) ;
745 val = getPdfVal();
746 for (k=0 ; k<Nbuf ; k++) {
747 tmp[N+Nbuf+k] = val ;
748 }
749 }
750 break ;
751
752 case Mirror:
753 // Sample original range (N samples) and fill lower and upper buffer
754 // bins with mirror image of sampled range
755 for (k=0 ; k<N ; k++) {
756 histX->setBin(k) ;
757 tmp[k+Nbuf] = getPdfVal();
758 }
759 for (k=1 ; k<=Nbuf ; k++) {
760 histX->setBin(k) ;
761 tmp[Nbuf-k] = getPdfVal();
762 histX->setBin(N-k) ;
763 tmp[Nbuf+N+k-1] = getPdfVal();
764 }
765 break ;
766 }
767
768 // Scan function and store values in array
769 for (Int_t i=0 ; i<N2 ; i++) {
770 // Cyclically shift writing location by zero bin position
771 Int_t j = i - (zeroBin) ;
772 if (j<0) j+= N2 ;
773 if (j>=N2) j-= N2 ;
774 array[i] = tmp[j] ;
775 }
776
777 return array ;
778}
779
780
781
782////////////////////////////////////////////////////////////////////////////////
783/// Return the observables to be cached given the normalization set nset.
784///
785/// If the cache observable is in nset then this is
786/// - the convolution observable plus
787/// - any member of nset that is either a RooCategory,
788/// - or was previously specified through setCacheObservables().
789///
790/// In case the cache observable is *not* in nset, then it is
791/// - the convolution observable plus
792/// - all member of nset that are observables of this p.d.f.
793///
794
796{
797 // Get complete list of observables
798 auto obs1 = new RooArgSet{};
800 _pdf1.arg().getObservables(&nset, *obs1) ;
801 _pdf2.arg().getObservables(&nset, obs2) ;
802 obs1->add(obs2, true) ;
803
804 // Check if convolution observable is in nset
805 if (nset.contains(_x.arg())) {
806
807 // Now strip out all non-category observables
809 for(auto const& arg : *obs1) {
810 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
811 killList.add(*arg) ;
812 }
813 }
814 obs1->remove(killList) ;
815
816 // And add back the convolution observables
817 obs1->add(_x.arg(),true) ;
818
819 obs1->add(_cacheObs) ;
820
821 } else {
822
823 // If cacheObs was filled, cache only observables in there
824 if (!_cacheObs.empty()) {
826 for(auto const& arg : *obs1) {
827 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
828 killList.add(*arg) ;
829 }
830 }
831 obs1->remove(killList) ;
832 }
833
834
835 // Make sure convolution observable is always in there
836 obs1->add(_x.arg(),true) ;
837
838 }
839
841}
842
843
844
845////////////////////////////////////////////////////////////////////////////////
846/// Return the parameters on which the cache depends given normalization
847/// set nset. For this p.d.f these are the parameters of the input p.d.f.
848/// but never the convolution variable, in case it is not part of nset.
849
851{
852 std::unique_ptr<RooArgSet> vars{getVariables()};
853 vars->remove(*std::unique_ptr<RooArgSet>{actualObservables(nset)});
854
855 return RooFit::makeOwningPtr(std::move(vars));
856}
857
858
859
860////////////////////////////////////////////////////////////////////////////////
861/// Return p.d.f. observable (which can be a function) to substitute given
862/// p.d.f. observable. Substitutes x by xprime if xprime is set.
863
865{
866 if (_xprime.absArg() && string(histObservable.GetName())==_x.absArg()->GetName()) {
867 return (*_xprime.absArg()) ;
868 }
869 return histObservable ;
870}
871
872
873
874////////////////////////////////////////////////////////////////////////////////
875/// Create appropriate generator context for this convolution. If both input p.d.f.s support
876/// internal generation, if it is safe to use them and if no observables other than the convolution
877/// observable are requested for generation, use the specialized convolution generator context
878/// which implements a smearing strategy in the convolution observable. If not return the
879/// regular accept/reject generator context
880
882 const RooArgSet* auxProto, bool verbose) const
883{
884 RooArgSet vars2(vars) ;
885 vars2.remove(_x.arg(),true,true) ;
886 Int_t numAddDep = vars2.size() ;
887
888 RooArgSet dummy ;
889 bool pdfCanDir = ((static_cast<RooAbsPdf const &>(_pdf1.arg())).getGenerator(_x.arg(), dummy) != 0 &&
890 (static_cast<RooAbsPdf const &>(_pdf1.arg())).isDirectGenSafe(_x.arg()));
891 bool resCanDir = ((static_cast<RooAbsPdf const &>(_pdf2.arg())).getGenerator(_x.arg(), dummy) != 0 &&
892 (static_cast<RooAbsPdf const &>(_pdf2.arg())).isDirectGenSafe(_x.arg()));
893
894 if (pdfCanDir) {
895 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f " << _pdf1.arg().GetName()
896 << " has internal generator that is safe to use in current context" << std::endl ;
897 }
898 if (resCanDir) {
899 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f. " << _pdf2.arg().GetName()
900 << " has internal generator that is safe to use in current context" << std::endl ;
901 }
902 if (numAddDep>0) {
903 cxcoutI(Generation) << "RooFFTConvPdf::genContext() generation requested for observables other than the convolution observable " << _x.arg().GetName() << std::endl ;
904 }
905
906
907 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
908 // Any resolution model with more dependents than the convolution variable
909 // or pdf or resmodel do not support direct generation
910 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting accept/reject generator context because one or both of the input "
911 << "p.d.f.s cannot use internal generator and/or "
912 << "observables other than the convolution variable are requested for generation" << std::endl ;
913 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
914 }
915
916 // Any other resolution model: use specialized generator context
917 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting specialized convolution generator context as both input "
918 << "p.d.fs are safe for internal generator and only "
919 << "the convolution observables is requested for generation" << std::endl ;
920 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
921}
922
923
924
925////////////////////////////////////////////////////////////////////////////////
926/// Change the size of the buffer on either side of the observable range to `frac` times the
927/// size of the range of the convolution observable.
928
930{
931 if (frac<0) {
932 coutE(InputArguments) << "RooFFTConvPdf::setBufferFraction(" << GetName() << ") fraction should be greater than or equal to zero" << std::endl ;
933 return ;
934 }
935 _bufFrac = frac ;
936
937 // Sterilize the cache as certain partial results depend on buffer fraction
939}
940
941
942////////////////////////////////////////////////////////////////////////////////
943/// Change strategy to fill the overflow buffer on either side of the convolution observable range.
944///
945/// - `Extend` means is that the input p.d.f convolution observable range is widened to include the buffer range
946/// - `Flat` means that the buffer is filled with the p.d.f. value at the boundary of the observable range
947/// - `Mirror` means that the buffer is filled with a mirror image of the p.d.f. around the convolution observable boundary
948///
949/// The default strategy is extend. If one of the input p.d.f.s is a RooAddPdf, it is configured so that the interpretation
950/// range of the fraction coefficients is kept at the nominal convolutions observable range (instead of interpreting coefficients
951/// in the widened range including the buffer).
952
957
958
959
960////////////////////////////////////////////////////////////////////////////////
961/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
962/// product operator construction
963
964void RooFFTConvPdf::printMetaArgs(ostream& os) const
965{
966 os << _pdf1.arg().GetName() << "(" << _x.arg().GetName() << ") (*) " << _pdf2.arg().GetName() << "(" << _x.arg().GetName() << ") " ;
967}
968
969
970
971////////////////////////////////////////////////////////////////////////////////
972/// (Re)calculate effective parameters of this p.d.f.
973
975{
976 RooArgSet argSet{_x.arg()};
979 _pdf1.arg().getParameters(&argSet, params1);
980 _pdf2.arg().getParameters(&argSet, params2);
983 _params.add(params2,true) ;
984}
#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
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
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:2495
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
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
! The cache manager
bool contains(const char *name) const
Check if collection contains an argument with a specific name.
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:32
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
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:32
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.
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.
std::vector< double > scanPdf(RooRealVar &obs, RooAbsPdf &pdf, double normVal, 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 ...
~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:138
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
double normalizeWithNaNPacking(RooAbsPdf const &pdf, double rawVal, double normVal)
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()