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.ch/building-root
85/// * Install FFTW and let ROOT discover it. `fftw3` is on by default (see https://root.cern.ch/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.ch/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.ch/building-root
109///
110
111
112#include <iostream>
113
114#include "RooFit.h"
115#include "RooFFTConvPdf.h"
116#include "RooAbsReal.h"
117#include "RooMsgService.h"
118#include "RooDataHist.h"
119#include "RooHistPdf.h"
120#include "RooRealVar.h"
121#include "TComplex.h"
122#include "TVirtualFFT.h"
123#include "RooGenContext.h"
124#include "RooConvGenContext.h"
125#include "RooBinning.h"
126#include "RooLinearVar.h"
127#include "RooCustomizer.h"
128#include "RooGlobalFunc.h"
129#include "RooConstVar.h"
130#include "TClass.h"
131#include "RooUniformBinning.h"
132
133using namespace std;
134
136
137
138////////////////////////////////////////////////////////////////////////////////
139/// Constructor for numerical (FFT) convolution of PDFs.
140/// \param[in] name Name of this PDF
141/// \param[in] title Title for plotting this PDF
142/// \param[in] convVar Observable to convolve the PDFs in \attention Use a high number of bins (>= 1000) for good accuracy.
143/// \param[in] pdf1 First PDF to be convolved
144/// \param[in] pdf2 Second PDF to be convolved
145/// \param[in] ipOrder Order for interpolation between bins (since FFT is discrete)
146/// The binning used for the FFT sampling is controlled by the binning named "cache" in the convolution observable `convVar`.
147/// If such a binning is not set, the same number of bins as for `convVar` will be used.
148
149RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder) :
150 RooAbsCachedPdf(name,title,ipOrder),
151 _x("!x","Convolution Variable",this,convVar),
152 _xprime("!xprime","External Convolution Variable",this,0),
153 _pdf1("!pdf1","pdf1",this,pdf1,kFALSE),
154 _pdf2("!pdf2","pdf2",this,pdf2,kFALSE),
155 _params("!params","effective parameters",this),
156 _bufFrac(0.1),
157 _bufStrat(Extend),
158 _shift1(0),
159 _shift2(0),
160 _cacheObs("!cacheObs","Cached observables",this,kFALSE,kFALSE)
161{
162 prepareFFTBinning(convVar);
163
164 _shift2 = (convVar.getMax("cache")+convVar.getMin("cache"))/2 ;
165
166 calcParams() ;
167
168}
169
170////////////////////////////////////////////////////////////////////////////////
171/// \copydoc RooFFTConvPdf(const char*, const char*, RooRealVar&, RooAbsPdf&, RooAbsPdf&, Int_t)
172/// \param[in] pdfConvVar If the variable used for convolution is a PDF, itself, pass the PDF here, and pass the convolution variable to
173/// `convVar`. See also rf210_angularconv.C in the <a href="https://root.cern.ch/root/html/tutorials/roofit/index.html.">roofit tutorials</a>
174
175RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooAbsReal& pdfConvVar, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder) :
176 RooAbsCachedPdf(name,title,ipOrder),
177 _x("!x","Convolution Variable",this,convVar,kFALSE,kFALSE),
178 _xprime("!xprime","External Convolution Variable",this,pdfConvVar),
179 _pdf1("!pdf1","pdf1",this,pdf1,kFALSE),
180 _pdf2("!pdf2","pdf2",this,pdf2,kFALSE),
181 _params("!params","effective parameters",this),
182 _bufFrac(0.1),
183 _bufStrat(Extend),
184 _shift1(0),
185 _shift2(0),
186 _cacheObs("!cacheObs","Cached observables",this,kFALSE,kFALSE)
187{
188 prepareFFTBinning(convVar);
189
190 _shift2 = (convVar.getMax("cache")+convVar.getMin("cache"))/2 ;
191
192 calcParams() ;
193}
194
195
196
197////////////////////////////////////////////////////////////////////////////////
198/// Copy constructor
199
201 RooAbsCachedPdf(other,name),
202 _x("!x",this,other._x),
203 _xprime("!xprime",this,other._xprime),
204 _pdf1("!pdf1",this,other._pdf1),
205 _pdf2("!pdf2",this,other._pdf2),
206 _params("!params",this,other._params),
207 _bufFrac(other._bufFrac),
208 _bufStrat(other._bufStrat),
209 _shift1(other._shift1),
210 _shift2(other._shift2),
211 _cacheObs("!cacheObs",this,other._cacheObs)
212 {
213 }
214
215
216
217////////////////////////////////////////////////////////////////////////////////
218/// Destructor
219
221{
222}
223
224
225////////////////////////////////////////////////////////////////////////////////
226/// Try to improve the binning and inform user if possible.
227/// With a 10% buffer fraction, 930 raw bins yield 1024 FFT bins,
228/// a sweet spot for the speed of FFTW.
229
231 if (!convVar.hasBinning("cache")) {
232 const RooAbsBinning& varBinning = convVar.getBinning();
233 const int optimal = static_cast<Int_t>(1024/(1.+_bufFrac));
234
235 //Can improve precision if binning is uniform
236 if (varBinning.numBins() < optimal && varBinning.isUniform()) {
237 coutI(Caching) << "Changing internal binning of variable '" << convVar.GetName()
238 << "' in FFT '" << fName << "'"
239 << " from " << varBinning.numBins()
240 << " to " << optimal << " to improve the precision of the numerical FFT."
241 << " This can be done manually by setting an additional binning named 'cache'." << std::endl;
242 convVar.setBinning(RooUniformBinning(varBinning.lowBound(), varBinning.highBound(), optimal, "cache"), "cache");
243 } else {
244 coutE(Caching) << "The internal binning of variable " << convVar.GetName()
245 << " is not uniform. The numerical FFT will likely yield wrong results." << std::endl;
246 convVar.setBinning(varBinning, "cache");
247 }
248 }
249}
250
251
252////////////////////////////////////////////////////////////////////////////////
253/// Return base name component for cache components in this case 'PDF1_CONV_PDF2'
254
256{
257 static TString name ;
258 name = _pdf1.arg().GetName() ;
259 name.Append("_CONV_") ;
260 name.Append(_pdf2.arg().GetName()) ;
261 return name.Data() ;
262}
263
264
265
266
267////////////////////////////////////////////////////////////////////////////////
268/// Return specialized cache subclass for FFT calculations
269
271{
272 return new FFTCacheElem(*this,nset) ;
273}
274
275
276
277
278////////////////////////////////////////////////////////////////////////////////
279/// Clone input pdf and attach to dataset
280
282 PdfCacheElem(self,nsetIn),
283 fftr2c1(0),fftr2c2(0),fftc2r(0)
284{
285 RooAbsPdf* clonePdf1 = (RooAbsPdf*) self._pdf1.arg().cloneTree() ;
286 RooAbsPdf* clonePdf2 = (RooAbsPdf*) self._pdf2.arg().cloneTree() ;
287 clonePdf1->attachDataSet(*hist()) ;
288 clonePdf2->attachDataSet(*hist()) ;
289
290 // Shift observable
291 RooRealVar* convObs = (RooRealVar*) hist()->get()->find(self._x.arg().GetName()) ;
292
293 // Install FFT reference range
294 string refName = Form("refrange_fft_%s",self.GetName()) ;
295 convObs->setRange(refName.c_str(),convObs->getMin(),convObs->getMax()) ;
296
297 if (self._shift1!=0) {
298 RooLinearVar* shiftObs1 = new RooLinearVar(Form("%s_shifted_FFTBuffer1",convObs->GetName()),"shiftObs1",
299 *convObs,RooFit::RooConst(1),RooFit::RooConst(-1*self._shift1)) ;
300
301 RooArgSet clonedBranches1 ;
302 RooCustomizer cust(*clonePdf1,"fft") ;
303 cust.replaceArg(*convObs,*shiftObs1) ;
304
305 pdf1Clone = (RooAbsPdf*) cust.build() ;
306
307 pdf1Clone->addOwnedComponents(*shiftObs1) ;
308 pdf1Clone->addOwnedComponents(*clonePdf1) ;
309
310 } else {
311 pdf1Clone = clonePdf1 ;
312 }
313
314 if (self._shift2!=0) {
315 RooLinearVar* shiftObs2 = new RooLinearVar(Form("%s_shifted_FFTBuffer2",convObs->GetName()),"shiftObs2",
316 *convObs,RooFit::RooConst(1),RooFit::RooConst(-1*self._shift2)) ;
317
318 RooArgSet clonedBranches2 ;
319 RooCustomizer cust(*clonePdf2,"fft") ;
320 cust.replaceArg(*convObs,*shiftObs2) ;
321
322 pdf1Clone->addOwnedComponents(*shiftObs2) ;
323 pdf1Clone->addOwnedComponents(*clonePdf2) ;
324
325 pdf2Clone = (RooAbsPdf*) cust.build() ;
326
327 } else {
328 pdf2Clone = clonePdf2 ;
329 }
330
331
332 // Attach cloned pdf to all original parameters of self
333 RooArgSet* fftParams = self.getParameters(*convObs) ;
334
335 // Remove all cache histogram from fftParams as these
336 // observable need to remain attached to the histogram
337 fftParams->remove(*hist()->get(),kTRUE,kTRUE) ;
338
341 pdf1Clone->fixAddCoefRange(refName.c_str(), true) ;
342 pdf2Clone->fixAddCoefRange(refName.c_str(), true) ;
343
344 // Ensure that coefficients for Add PDFs are only interpreted with respect to the convolution observable
345 RooArgSet convSet(self._x.arg());
346 pdf1Clone->fixAddCoefNormalization(convSet, true);
347 pdf2Clone->fixAddCoefNormalization(convSet, true);
348
349 delete fftParams ;
350
351 // Save copy of original histX binning and make alternate binning
352 // for extended range scanning
353
354 const Int_t N = convObs->numBins();
355 if (N < 900) {
356 oocoutW(&self, Eval) << "The FFT convolution '" << self.GetName() << "' will run with " << N
357 << " bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number"
358 " of bins of the observable '" << convObs->GetName() << "'." << std::endl;
359 }
360 Int_t Nbuf = static_cast<Int_t>((N*self.bufferFraction())/2 + 0.5) ;
361 Double_t obw = (convObs->getMax() - convObs->getMin())/N ;
362 Int_t N2 = N+2*Nbuf ;
363
364 scanBinning = new RooUniformBinning (convObs->getMin()-Nbuf*obw,convObs->getMax()+Nbuf*obw,N2) ;
365 histBinning = convObs->getBinning().clone() ;
366
367 // Deactivate dirty state propagation on datahist observables
368 // and set all nodes on both pdfs to operMode AlwaysDirty
370 convObs->setOperMode(ADirty,kTRUE) ;
371}
372
373
374////////////////////////////////////////////////////////////////////////////////
375/// Suffix for cache histogram (added in addition to suffix for cache name)
376
378{
379 return TString(Form("_BufFrac%3.1f_BufStrat%d",_bufFrac,_bufStrat)) ;
380}
381
382
383
384////////////////////////////////////////////////////////////////////////////////
385/// Returns all RooAbsArg objects contained in the cache element
386
388{
389 RooArgList ret(PdfCacheElem::containedArgs(a)) ;
390
391 ret.add(*pdf1Clone) ;
392 ret.add(*pdf2Clone) ;
393 if (pdf1Clone->ownedComponents()) {
394 ret.add(*pdf1Clone->ownedComponents()) ;
395 }
396 if (pdf2Clone->ownedComponents()) {
397 ret.add(*pdf2Clone->ownedComponents()) ;
398 }
399
400 return ret ;
401}
402
403
404////////////////////////////////////////////////////////////////////////////////
405
407{
408 delete fftr2c1 ;
409 delete fftr2c2 ;
410 delete fftc2r ;
411
412 delete pdf1Clone ;
413 delete pdf2Clone ;
414
415 delete histBinning ;
416 delete scanBinning ;
417
418}
419
420
421
422
423////////////////////////////////////////////////////////////////////////////////
424/// Fill the contents of the cache the FFT convolution output
425
427{
428 RooDataHist& cacheHist = *cache.hist() ;
429
430 ((FFTCacheElem&)cache).pdf1Clone->setOperMode(ADirty,kTRUE) ;
431 ((FFTCacheElem&)cache).pdf2Clone->setOperMode(ADirty,kTRUE) ;
432
433 // Determine if there other observables than the convolution observable in the cache
434 RooArgSet otherObs ;
435 RooArgSet(*cacheHist.get()).snapshot(otherObs) ;
436
437 RooAbsArg* histArg = otherObs.find(_x.arg().GetName()) ;
438 if (histArg) {
439 otherObs.remove(*histArg,kTRUE,kTRUE) ;
440 }
441
442 //cout << "RooFFTConvPdf::fillCacheObject() otherObs = " << otherObs << endl ;
443
444 // Handle trivial scenario -- no other observables
445 if (otherObs.getSize()==0) {
447 return ;
448 }
449
450 // Handle cases where there are other cache slices
451 // Iterator over available slice positions and fill each
452
453 // Determine number of bins for each slice position observable
454 Int_t n = otherObs.getSize() ;
455 Int_t* binCur = new Int_t[n+1] ;
456 Int_t* binMax = new Int_t[n+1] ;
457 Int_t curObs = 0 ;
458
459 RooAbsLValue** obsLV = new RooAbsLValue*[n] ;
460 TIterator* iter = otherObs.createIterator() ;
461 RooAbsArg* arg ;
462 Int_t i(0) ;
463 while((arg=(RooAbsArg*)iter->Next())) {
464 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(arg) ;
465 obsLV[i] = lvarg ;
466 binCur[i] = 0 ;
467 // coverity[FORWARD_NULL]
468 binMax[i] = lvarg->numBins(binningName())-1 ;
469 i++ ;
470 }
471 delete iter ;
472
473 Bool_t loop(kTRUE) ;
474 while(loop) {
475 // Set current slice position
476 for (Int_t j=0 ; j<n ; j++) { obsLV[j]->setBin(binCur[j],binningName()) ; }
477
478// cout << "filling slice: bin of obsLV[0] = " << obsLV[0]->getBin() << endl ;
479
480 // Fill current slice
481 fillCacheSlice((FFTCacheElem&)cache,otherObs) ;
482
483 // Determine which iterator to increment
484 while(binCur[curObs]==binMax[curObs]) {
485
486 // Reset current iterator and consider next iterator ;
487 binCur[curObs]=0 ;
488 curObs++ ;
489
490 // master termination condition
491 if (curObs==n) {
492 loop=kFALSE ;
493 break ;
494 }
495 }
496
497 // Increment current iterator
498 binCur[curObs]++ ;
499 curObs=0 ;
500
501 }
502
503 delete[] obsLV ;
504 delete[] binMax ;
505 delete[] binCur ;
506
507}
508
509
510////////////////////////////////////////////////////////////////////////////////
511/// Fill a slice of cachePdf with the output of the FFT convolution calculation
512
513void RooFFTConvPdf::fillCacheSlice(FFTCacheElem& aux, const RooArgSet& slicePos) const
514{
515 // Extract histogram that is the basis of the RooHistPdf
516 RooDataHist& cacheHist = *aux.hist() ;
517
518 // Sample array of input points from both pdfs
519 // Note that returned arrays have optional buffers zones below and above range ends
520 // to reduce cyclical effects and have been cyclically rotated so that bin containing
521 // zero value is at position zero. Example:
522 //
523 // original: -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5
524 // add buffer zones: U U -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5 O O
525 // rotate: 0 +1 +2 +3 +4 +5 O O U U -5 -4 -3 -2 -1
526 //
527 //
528
529 Int_t N,N2,binShift1,binShift2 ;
530
531 RooRealVar* histX = (RooRealVar*) cacheHist.get()->find(_x.arg().GetName()) ;
532 if (_bufStrat==Extend) histX->setBinning(*aux.scanBinning) ;
533 Double_t* input1 = scanPdf((RooRealVar&)_x.arg(),*aux.pdf1Clone,cacheHist,slicePos,N,N2,binShift1,_shift1) ;
534 Double_t* input2 = scanPdf((RooRealVar&)_x.arg(),*aux.pdf2Clone,cacheHist,slicePos,N,N2,binShift2,_shift2) ;
535 if (_bufStrat==Extend) histX->setBinning(*aux.histBinning) ;
536
537
538
539
540 // Retrieve previously defined FFT transformation plans
541 if (!aux.fftr2c1) {
542 aux.fftr2c1 = TVirtualFFT::FFT(1, &N2, "R2CK");
543 aux.fftr2c2 = TVirtualFFT::FFT(1, &N2, "R2CK");
544 aux.fftc2r = TVirtualFFT::FFT(1, &N2, "C2RK");
545 }
546
547 // Real->Complex FFT Transform on p.d.f. 1 sampling
548 aux.fftr2c1->SetPoints(input1);
549 aux.fftr2c1->Transform();
550
551 // Real->Complex FFT Transform on p.d.f 2 sampling
552 aux.fftr2c2->SetPoints(input2);
553 aux.fftr2c2->Transform();
554
555 // Loop over first half +1 of complex output results, multiply
556 // and set as input of reverse transform
557 for (Int_t i=0 ; i<N2/2+1 ; i++) {
558 Double_t re1,re2,im1,im2 ;
559 aux.fftr2c1->GetPointComplex(i,re1,im1) ;
560 aux.fftr2c2->GetPointComplex(i,re2,im2) ;
561 Double_t re = re1*re2 - im1*im2 ;
562 Double_t im = re1*im2 + re2*im1 ;
563 TComplex t(re,im) ;
564 aux.fftc2r->SetPointComplex(i,t) ;
565 }
566
567 // Reverse Complex->Real FFT transform product
568 aux.fftc2r->Transform() ;
569
570 Int_t totalShift = binShift1 + (N2-N)/2 ;
571
572 // Store FFT result in cache
573
574 TIterator* iter = const_cast<RooDataHist&>(cacheHist).sliceIterator(const_cast<RooAbsReal&>(_x.arg()),slicePos) ;
575 for (Int_t i =0 ; i<N ; i++) {
576
577 // Cyclically shift array back so that bin containing zero is back in zeroBin
578 Int_t j = i + totalShift ;
579 while (j<0) j+= N2 ;
580 while (j>=N2) j-= N2 ;
581
582 iter->Next() ;
583 cacheHist.set(aux.fftc2r->GetPointReal(j)) ;
584 }
585 delete iter ;
586
587 // cacheHist.dump2() ;
588
589 // Delete input arrays
590 delete[] input1 ;
591 delete[] input2 ;
592
593}
594
595
596////////////////////////////////////////////////////////////////////////////////
597/// Scan the values of 'pdf' in observable 'obs' using the bin values stored in 'hist' at slice position 'slicePos'
598/// N is filled with the number of bins defined in hist, N2 is filled with N plus the number of buffer bins
599/// The return value is an array of doubles of length N2 with the sampled values. The caller takes ownership
600/// of the array
601
602Double_t* RooFFTConvPdf::scanPdf(RooRealVar& obs, RooAbsPdf& pdf, const RooDataHist& hist, const RooArgSet& slicePos,
603 Int_t& N, Int_t& N2, Int_t& zeroBin, Double_t shift) const
604{
605
606 RooRealVar* histX = (RooRealVar*) hist.get()->find(obs.GetName()) ;
607
608 // Calculate number of buffer bins on each size to avoid cyclical flow
609 N = histX->numBins(binningName()) ;
610 Int_t Nbuf = static_cast<Int_t>((N*bufferFraction())/2 + 0.5) ;
611 N2 = N+2*Nbuf ;
612
613
614 // Allocate array of sampling size plus optional buffer zones
615 Double_t* array = new Double_t[N2] ;
616
617 // Set position of non-convolution observable to that of the cache slice that were are processing now
618 hist.get(slicePos) ;
619
620 // Find bin ID that contains zero value
621 zeroBin = 0 ;
622 if (histX->getMax()>=0 && histX->getMin()<=0) {
623 zeroBin = histX->getBinning().binNumber(0) ;
624 } else if (histX->getMin()>0) {
625 Double_t bw = (histX->getMax() - histX->getMin())/N2 ;
626 zeroBin = Int_t(-histX->getMin()/bw) ;
627 } else {
628 Double_t bw = (histX->getMax() - histX->getMin())/N2 ;
629 zeroBin = Int_t(-1*histX->getMax()/bw) ;
630 }
631
632 Int_t binShift = Int_t((N2* shift) / (histX->getMax()-histX->getMin())) ;
633
634 zeroBin += binShift ;
635 while(zeroBin>=N2) zeroBin-= N2 ;
636 while(zeroBin<0) zeroBin+= N2 ;
637
638 // First scan hist into temp array
639 Double_t *tmp = new Double_t[N2] ;
640 Int_t k(0) ;
641 switch(_bufStrat) {
642
643 case Extend:
644 // Sample entire extended range (N2 samples)
645 for (k=0 ; k<N2 ; k++) {
646 histX->setBin(k) ;
647 tmp[k] = pdf.getVal(hist.get()) ;
648 }
649 break ;
650
651 case Flat:
652 // Sample original range (N samples) and fill lower and upper buffer
653 // bins with p.d.f. value at respective boundary
654 {
655 histX->setBin(0) ;
656 Double_t val = pdf.getVal(hist.get()) ;
657 for (k=0 ; k<Nbuf ; k++) {
658 tmp[k] = val ;
659 }
660 for (k=0 ; k<N ; k++) {
661 histX->setBin(k) ;
662 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
663 }
664 histX->setBin(N-1) ;
665 val = pdf.getVal(hist.get()) ;
666 for (k=0 ; k<Nbuf ; k++) {
667 tmp[N+Nbuf+k] = val ;
668 }
669 }
670 break ;
671
672 case Mirror:
673 // Sample original range (N samples) and fill lower and upper buffer
674 // bins with mirror image of sampled range
675 for (k=0 ; k<N ; k++) {
676 histX->setBin(k) ;
677 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
678 }
679 for (k=1 ; k<=Nbuf ; k++) {
680 histX->setBin(k) ;
681 tmp[Nbuf-k] = pdf.getVal(hist.get()) ;
682 histX->setBin(N-k) ;
683 tmp[Nbuf+N+k-1] = pdf.getVal(hist.get()) ;
684 }
685 break ;
686 }
687
688 // Scan function and store values in array
689 for (Int_t i=0 ; i<N2 ; i++) {
690 // Cyclically shift writing location by zero bin position
691 Int_t j = i - (zeroBin) ;
692 if (j<0) j+= N2 ;
693 if (j>=N2) j-= N2 ;
694 array[i] = tmp[j] ;
695 }
696
697 // Cleanup
698 delete[] tmp ;
699 return array ;
700}
701
702
703
704////////////////////////////////////////////////////////////////////////////////
705/// Return the observables to be cached given the normalization set nset.
706///
707/// If the cache observable is in nset then this is
708/// - the convolution observable plus
709/// - any member of nset that is either a RooCategory,
710/// - or was previously specified through setCacheObservables().
711///
712/// In case the cache observable is *not* in nset, then it is
713/// - the convolution observable plus
714/// - all member of nset that are observables of this p.d.f.
715///
716
718{
719 // Get complete list of observables
720 RooArgSet* obs1 = _pdf1.arg().getObservables(nset) ;
721 RooArgSet* obs2 = _pdf2.arg().getObservables(nset) ;
722 obs1->add(*obs2,kTRUE) ;
723
724 // Check if convolution observable is in nset
725 if (nset.contains(_x.arg())) {
726
727 // Now strip out all non-category observables
728 TIterator* iter = obs1->createIterator() ;
729 RooAbsArg* arg ;
730 RooArgSet killList ;
731 while((arg=(RooAbsArg*)iter->Next())) {
732 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
733 killList.add(*arg) ;
734 }
735 }
736 delete iter ;
737 obs1->remove(killList) ;
738
739 // And add back the convolution observables
740 obs1->add(_x.arg(),kTRUE) ;
741
742 obs1->add(_cacheObs) ;
743
744 delete obs2 ;
745
746 } else {
747
748 // If cacheObs was filled, cache only observables in there
749 if (_cacheObs.getSize()>0) {
750 TIterator* iter = obs1->createIterator() ;
751 RooAbsArg* arg ;
752 RooArgSet killList ;
753 while((arg=(RooAbsArg*)iter->Next())) {
754 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
755 killList.add(*arg) ;
756 }
757 }
758 delete iter ;
759 obs1->remove(killList) ;
760 }
761
762
763 // Make sure convolution observable is always in there
764 obs1->add(_x.arg(),kTRUE) ;
765 delete obs2 ;
766
767 }
768
769 return obs1 ;
770}
771
772
773
774////////////////////////////////////////////////////////////////////////////////
775/// Return the parameters on which the cache depends given normalization
776/// set nset. For this p.d.f these are the parameters of the input p.d.f.
777/// but never the convolution variable, in case it is not part of nset.
778
780{
781 RooArgSet* vars = getVariables() ;
782 RooArgSet* obs = actualObservables(nset) ;
783 vars->remove(*obs) ;
784 delete obs ;
785
786 return vars ;
787}
788
789
790
791////////////////////////////////////////////////////////////////////////////////
792/// Return p.d.f. observable (which can be a function) to substitute given
793/// p.d.f. observable. Substitutes x by xprime if xprime is set.
794
796{
797 if (_xprime.absArg() && string(histObservable.GetName())==_x.absArg()->GetName()) {
798 return (*_xprime.absArg()) ;
799 }
800 return histObservable ;
801}
802
803
804
805////////////////////////////////////////////////////////////////////////////////
806/// Create appropriate generator context for this convolution. If both input p.d.f.s support
807/// internal generation, if it is safe to use them and if no observables other than the convolution
808/// observable are requested for generation, use the specialized convolution generator context
809/// which implements a smearing strategy in the convolution observable. If not return the
810/// regular accept/reject generator context
811
813 const RooArgSet* auxProto, Bool_t verbose) const
814{
815 RooArgSet vars2(vars) ;
816 vars2.remove(_x.arg(),kTRUE,kTRUE) ;
817 Int_t numAddDep = vars2.getSize() ;
818
819 RooArgSet dummy ;
820 Bool_t pdfCanDir = (((RooAbsPdf&)_pdf1.arg()).getGenerator(_x.arg(),dummy) != 0 && \
822 Bool_t resCanDir = (((RooAbsPdf&)_pdf2.arg()).getGenerator(_x.arg(),dummy) !=0 &&
824
825 if (pdfCanDir) {
826 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f " << _pdf1.arg().GetName()
827 << " has internal generator that is safe to use in current context" << endl ;
828 }
829 if (resCanDir) {
830 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f. " << _pdf2.arg().GetName()
831 << " has internal generator that is safe to use in current context" << endl ;
832 }
833 if (numAddDep>0) {
834 cxcoutI(Generation) << "RooFFTConvPdf::genContext() generation requested for observables other than the convolution observable " << _x.arg().GetName() << endl ;
835 }
836
837
838 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
839 // Any resolution model with more dependents than the convolution variable
840 // or pdf or resmodel do not support direct generation
841 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting accept/reject generator context because one or both of the input "
842 << "p.d.f.s cannot use internal generator and/or "
843 << "observables other than the convolution variable are requested for generation" << endl ;
844 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
845 }
846
847 // Any other resolution model: use specialized generator context
848 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting specialized convolution generator context as both input "
849 << "p.d.fs are safe for internal generator and only "
850 << "the convolution observables is requested for generation" << endl ;
851 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
852}
853
854
855
856////////////////////////////////////////////////////////////////////////////////
857/// Change the size of the buffer on either side of the observable range to `frac` times the
858/// size of the range of the convolution observable.
859
861{
862 if (frac<0) {
863 coutE(InputArguments) << "RooFFTConvPdf::setBufferFraction(" << GetName() << ") fraction should be greater than or equal to zero" << endl ;
864 return ;
865 }
866 _bufFrac = frac ;
867
868 // Sterilize the cache as certain partial results depend on buffer fraction
870}
871
872
873////////////////////////////////////////////////////////////////////////////////
874/// Change strategy to fill the overflow buffer on either side of the convolution observable range.
875///
876/// - `Extend` means is that the input p.d.f convolution observable range is widened to include the buffer range
877/// - `Flat` means that the buffer is filled with the p.d.f. value at the boundary of the observable range
878/// - `Mirror` means that the buffer is filled with a mirror image of the p.d.f. around the convolution observable boundary
879///
880/// The default strategy is extend. If one of the input p.d.f.s is a RooAddPdf, it is configured so that the interpretation
881/// range of the fraction coefficients is kept at the nominal convolutions observable range (instead of interpreting coefficients
882/// in the widened range including the buffer).
883
885{
886 _bufStrat = bs ;
887}
888
889
890
891////////////////////////////////////////////////////////////////////////////////
892/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
893/// product operator construction
894
895void RooFFTConvPdf::printMetaArgs(ostream& os) const
896{
897 os << _pdf1.arg().GetName() << "(" << _x.arg().GetName() << ") (*) " << _pdf2.arg().GetName() << "(" << _x.arg().GetName() << ") " ;
898}
899
900
901
902////////////////////////////////////////////////////////////////////////////////
903/// (Re)calculate effective parameters of this p.d.f.
904
906{
907 RooArgSet* params1 = _pdf1.arg().getParameters(_x.arg()) ;
908 RooArgSet* params2 = _pdf2.arg().getParameters(_x.arg()) ;
910 _params.add(*params1) ;
911 _params.add(*params2,kTRUE) ;
912 delete params1 ;
913 delete params2 ;
914}
915
916
917
918////////////////////////////////////////////////////////////////////////////////
919///calcParams() ;
920
921Bool_t RooFFTConvPdf::redirectServersHook(const RooAbsCollection& /*newServerList*/, Bool_t /*mustReplaceAll*/, Bool_t /*nameChange*/, Bool_t /*isRecursive*/)
922{
923 return kFALSE ;
924}
#define a(i)
Definition RSha256.hxx:99
#define coutI(a)
#define cxcoutI(a)
#define oocoutW(o, a)
#define coutE(a)
int Int_t
Definition RtypesCore.h:45
const Bool_t kFALSE
Definition RtypesCore.h:92
const Bool_t kTRUE
Definition RtypesCore.h:91
#define ClassImp(name)
Definition Rtypes.h:364
#define N
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:72
virtual RooAbsArg * cloneTree(const char *newname=0) const
Clone tree expression of objects.
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:313
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Set the operation mode of this node.
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
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...
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Recursively replace all servers with the new servers in newSet.
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.
virtual RooAbsBinning * clone(const char *name=0) const =0
Int_t numBins() const
Return number of bins.
virtual Bool_t isUniform() const
virtual Double_t highBound() const =0
virtual Int_t binNumber(Double_t x) const =0
virtual Double_t lowBound() const =0
RooAbsCachedPdf is the abstract base class for p.d.f.s that need or want to cache their evaluate() ou...
virtual const char * binningName() const
RooObjCacheManager _cacheMgr
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsArg * find(const char *name) const
Find object with given name in list.
void setDirtyProp(Bool_t flag)
Control propagation of dirty flags from observables in dataset.
RooAbsGenContext is the 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=0) const =0
virtual void setBin(Int_t ibin, const char *rangeName=0)=0
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual void setBin(Int_t ibin, const char *rangeName=0)
Set value to center of bin 'ibin' of binning 'rangeName' (or of default binning if no range is specif...
virtual Int_t numBins(const char *rangeName=0) const
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual void fixAddCoefRange(const char *rangeName=0, Bool_t force=kTRUE)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), Bool_t force=kTRUE)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooAbsArg * absArg() const
Definition RooArgProxy.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:118
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
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 occurence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, Bool_t verbose=kFALSE)
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.
Definition RooDataHist.h:37
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.
Definition RooDataHist.h:74
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
virtual RooArgList containedArgs(Action)
Returns all RooAbsArg objects contained in the cache element.
FFTCacheElem(const RooFFTConvPdf &self, const RooArgSet *nset)
Clone input pdf and attach to dataset.
PDF for the numerical (FFT) convolution of two PDFs.
Bool_t redirectServersHook(const RooAbsCollection &newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
calcParams() ;
void printMetaArgs(std::ostream &os) const
Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the p...
virtual PdfCacheElem * createCache(const RooArgSet *nset) const
Return specialized cache subclass for FFT calculations.
friend class RooConvGenContext
virtual ~RooFFTConvPdf()
Destructor.
RooSetProxy _params
BufStrat _bufStrat
void calcParams()
(Re)calculate effective parameters of this p.d.f.
virtual void fillCacheObject(PdfCacheElem &cache) const
Fill the contents of the cache the FFT convolution output.
void prepareFFTBinning(RooRealVar &convVar) const
Try to improve the binning and inform user if possible.
virtual TString histNameSuffix() const
Suffix for cache histogram (added in addition to suffix for cache name)
Double_t _bufFrac
void fillCacheSlice(FFTCacheElem &cache, const RooArgSet &slicePosition) const
Fill a slice of cachePdf with the output of the FFT convolution calculation.
RooRealProxy _xprime
Double_t * scanPdf(RooRealVar &obs, RooAbsPdf &pdf, const RooDataHist &hist, const RooArgSet &slicePos, Int_t &N, Int_t &N2, Int_t &zeroBin, Double_t shift) const
Scan the values of 'pdf' in observable 'obs' using the bin values stored in 'hist' at slice position ...
RooRealProxy _pdf1
RooRealProxy _x
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Create appropriate generator context for this convolution.
virtual const char * inputBaseName() const
Return base name component for cache components in this case 'PDF1_CONV_PDF2'.
void setBufferFraction(Double_t frac)
Change the size of the buffer on either side of the observable range to frac times the size of the ra...
virtual RooArgSet * actualObservables(const RooArgSet &nset) const
Return the observables to be cached given the normalization set nset.
RooRealProxy _pdf2
virtual RooAbsArg & pdfObservable(RooAbsArg &histObservable) const
Return p.d.f.
void setBufferStrategy(BufStrat bs)
Change strategy to fill the overflow buffer on either side of the convolution observable range.
virtual RooArgSet * actualParameters(const RooArgSet &nset) const
Return the parameters on which the cache depends given normalization set nset.
Double_t bufferFraction() const
RooSetProxy _cacheObs
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()
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.
Definition RooRealVar.h:39
Bool_t hasBinning(const char *name) const
Returns true if variable has a binning named 'name'.
void setRange(const char *name, Double_t min, Double_t max)
Set a fit or plotting range.
const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const
Return binning definition with name.
void setBinning(const RooAbsBinning &binning, const char *name=0)
Add given binning under name 'name' with this variable.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
virtual void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
const T & arg() const
Return reference to object held in proxy.
RooUniformBinning is an implementation of RooAbsBinning that provides a uniform binning in 'n' bins b...
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
TString fName
Definition TNamed.h:32
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:445
Basic string class.
Definition TString.h:136
static TVirtualFFT * FFT(Int_t ndim, Int_t *n, Option_t *option)
Returns a pointer to the FFT of requested size and type.
virtual void SetPoints(const Double_t *data)=0
virtual void Transform()=0
virtual void GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const =0
virtual void SetPointComplex(Int_t ipoint, TComplex &c)=0
virtual Double_t GetPointReal(Int_t ipoint, Bool_t fromInput=kFALSE) const =0
RooConstVar & RooConst(Double_t val)
const Int_t n
Definition legend1.C:16