Logo ROOT   6.16/01
Reference Guide
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/// and calculate the convolution by calculate a Real->Complex FFT of both input p.d.fs
23/// multiplying the complex coefficients and performing the reverse Complex->Real FFT
24/// to get the result in the input space. This class using 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). More instructions 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 & desirable behavior for cyclical observables such as angles,
33/// but it may not be for other observables. The effect that is observed is that if
34/// p.d.f is zero at xMin and non-zero at xMax some spillover occurs and
35/// a rising tail may appear at xMin.
36/// This is inevitable when using FFTs, because the FFT behaves as if an input with e.g. 3 bins was:
37/// ` ... 0 1 2 0 1 2 0 1 2 ... `
38///
39/// Therefore, if bins 0 and 2 are not equal, the FFT sees a step at the boundary, which causes
40/// problems in Fourier space.
41///
42/// The spillover or discontinuity can be reduced or eliminated by
43/// introducing a buffer zone in the FFT calculation. If this feature is activated,
44/// the sampling array for the FFT calculation is extended in both directions
45/// and padded with the lowest/highest bin.
46/// Example:
47/// ```
48/// original: -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5
49/// add buffer zones: U U -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5 O O
50/// rotate: 0 +1 +2 +3 +4 +5 O O U U -5 -4 -3 -2 -1
51/// ```
52/// The buffer bins are stripped away when the FFT output values
53/// are transferred to the p.d.f cache. The default buffer size is 10% of the
54/// observable domain size and can be changed with the `setBufferFraction()` member function.
55///
56/// This class is a caching p.d.f inheriting from RooAbsCachedPdf. If this p.d.f
57/// is evaluated for a particular value of x, the FFT calculates the values for the
58/// p.d.f at all points in observable space for the given choice of parameters,
59/// which are stored in the cache. Subsequent evaluations of RooFFTConvPdf with
60/// identical parameters will retrieve results from the cache. If one or more
61/// of the parameters change, the cache will be updated, i.e., a new FFT runs.
62///
63/// The sampling density of the cache is controlled by the binning of the
64/// the convolution observable, which can be changed from RooRealVar::setBins(N)
65/// For good results N should be large (>1000). Additional interpolation of
66/// cache values may improve the result if coarse binnings are chosen. These can be
67/// set in the constructor or through the `setInterpolationOrder()` member function.
68/// For N>1000 interpolation will not substantially improve the performance.
69///
70/// Additionial information on caching activities can be displayed by monitoring
71/// the message stream with topic "Caching" at the INFO level, i.e.
72/// do `RooMsgService::instance().addStream(RooMsgService::INFO,Topic("Caching"))`
73/// to see these message on stdout.
74///
75/// Multi-dimensional convolutions are not supported at the moment.
76///
77/// ---
78///
79/// Installing a copy of FFTW on Linux and compiling ROOT to use it
80/// -------
81///
82/// You have two options:
83/// * ROOT can automatically install FFTW for itself, see `builtin_fftw3` at https://root.cern.ch/building-root
84/// * Install FFTW and let ROOT discover it. `fftw3` is on by default (see https://root.cern.ch/building-root)
85///
86/// 1) Go to www.fftw.org and download the latest stable version (a .tar.gz file)
87///
88/// If you have root access to your machine and want to make a system installation of FFTW
89///
90/// 2) Untar fftw-XXX.tar.gz in /tmp, cd into the untarred directory
91/// and type './configure' followed by 'make install'.
92/// This will install fftw in /usr/local/bin,lib etc...
93///
94/// 3) Start from a source installation of ROOT. ROOT should discover it. See https://root.cern.ch/building-root
95///
96///
97/// If you do not have root access and want to make a private installation of FFTW
98///
99/// 2) Make a private install area for FFTW, e.g. /home/myself/fftw
100///
101/// 3) Untar fftw-XXX.tar.gz in /tmp, cd into the untarred directory
102/// and type './configure --prefix=/home/myself/fftw' followed by 'make install'.
103/// Substitute /home/myself/fftw with a directory of your choice. This
104/// procedure will install FFTW in the location designated by you
105///
106/// 4) Start from a source installation of ROOT.
107/// Look up and set the proper paths for ROOT to discover FFTW. See https://root.cern.ch/building-root
108///
109
110
111#include "Riostream.h"
112
113#include "RooFit.h"
114#include "RooFFTConvPdf.h"
115#include "RooAbsReal.h"
116#include "RooMsgService.h"
117#include "RooDataHist.h"
118#include "RooHistPdf.h"
119#include "RooRealVar.h"
120#include "TComplex.h"
121#include "TVirtualFFT.h"
122#include "RooGenContext.h"
123#include "RooConvGenContext.h"
124#include "RooBinning.h"
125#include "RooLinearVar.h"
126#include "RooCustomizer.h"
127#include "RooGlobalFunc.h"
128#include "RooLinearVar.h"
129#include "RooConstVar.h"
130#include "TClass.h"
131#include "TSystem.h"
132
133using namespace std ;
134
136
137
138
139////////////////////////////////////////////////////////////////////////////////
140/// Constructor for numerical (FFT) convolution of PDFs.
141/// \param[in] name Name of this PDF
142/// \param[in] title Title for plotting this PDF
143/// \param[in] convVar Observable to convolve the PDFs in \attention Use a high number of bins (>= 1000) for good accuracy.
144/// \param[in] pdf1 First PDF to be convolved
145/// \param[in] pdf2 Second PDF to be convolved
146/// \param[in] ipOrder Order for interpolation of the discrete FFT outputs
147/// The binning used for the FFT sampling is controlled by the binning named "cache" in the convolution observable `convVar`.
148/// If such a binning is not set, the same number of bins as for `convVar` will be used.
149
150RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder) :
151 RooAbsCachedPdf(name,title,ipOrder),
152 _x("!x","Convolution Variable",this,convVar),
153 _xprime("!xprime","External Convolution Variable",this,0),
154 _pdf1("!pdf1","pdf1",this,pdf1,kFALSE),
155 _pdf2("!pdf2","pdf2",this,pdf2,kFALSE),
156 _params("!params","effective parameters",this),
157 _bufFrac(0.1),
158 _bufStrat(Extend),
159 _shift1(0),
160 _shift2(0),
161 _cacheObs("!cacheObs","Cached observables",this,kFALSE,kFALSE)
162{
163 prepareFFTBinning(convVar);
164
165 _shift2 = (convVar.getMax("cache")+convVar.getMin("cache"))/2 ;
166
167 calcParams() ;
168
169}
170
171////////////////////////////////////////////////////////////////////////////////
172/// \copydoc RooFFTConvPdf(const char*, const char*, RooRealVar&, RooAbsPdf&, RooAbsPdf&, Int_t)
173/// \param[in] pdfConvVar If the variable used for convolution has a PDF, itself, this holds the PDF and
174/// convVar holds the variable. See also rf210_angularconv.C in the <a href="https://root.cern.ch/root/html/tutorials/roofit/index.html.">roofit tutorials</a>
175
176RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooAbsReal& pdfConvVar, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder) :
177 RooAbsCachedPdf(name,title,ipOrder),
178 _x("!x","Convolution Variable",this,convVar,kFALSE,kFALSE),
179 _xprime("!xprime","External Convolution Variable",this,pdfConvVar),
180 _pdf1("!pdf1","pdf1",this,pdf1,kFALSE),
181 _pdf2("!pdf2","pdf2",this,pdf2,kFALSE),
182 _params("!params","effective parameters",this),
183 _bufFrac(0.1),
184 _bufStrat(Extend),
185 _shift1(0),
186 _shift2(0),
187 _cacheObs("!cacheObs","Cached observables",this,kFALSE,kFALSE)
188{
189 prepareFFTBinning(convVar);
190
191 _shift2 = (convVar.getMax("cache")+convVar.getMin("cache"))/2 ;
192
193 calcParams() ;
194}
195
196
197
198////////////////////////////////////////////////////////////////////////////////
199/// Copy constructor
200
202 RooAbsCachedPdf(other,name),
203 _x("!x",this,other._x),
204 _xprime("!xprime",this,other._xprime),
205 _pdf1("!pdf1",this,other._pdf1),
206 _pdf2("!pdf2",this,other._pdf2),
207 _params("!params",this,other._params),
208 _bufFrac(other._bufFrac),
209 _bufStrat(other._bufStrat),
210 _shift1(other._shift1),
211 _shift2(other._shift2),
212 _cacheObs("!cacheObs",this,other._cacheObs)
213 {
214 }
215
216
217
218////////////////////////////////////////////////////////////////////////////////
219/// Destructor
220
222{
223}
224
225
226////////////////////////////////////////////////////////////////////////////////
227/// Try to improve the binning and inform user if possible.
228/// With a 10% buffer fraction, 930 raw bins yield 1024 FFT bins,
229/// a sweet spot for the speed of FFTW.
230
232 if (!convVar.hasBinning("cache")) {
233 const RooAbsBinning& varBinning = convVar.getBinning();
234 const int optimal = static_cast<Int_t>(1024/(1.+_bufFrac));
235
236 //Can improve precision if binning is uniform
237 if (varBinning.numBins() < optimal && varBinning.isUniform()) {
238 coutI(Caching) << "Changing internal binning of variable '" << convVar.GetName()
239 << "' in FFT '" << fName << "'"
240 << " from " << varBinning.numBins()
241 << " to " << optimal << " to improve the precision of the numerical FFT."
242 << " This can be done manually by setting an additional binning named 'cache'." << std::endl;
243 convVar.setBinning(RooUniformBinning(varBinning.lowBound(), varBinning.highBound(), optimal, "cache"), "cache");
244 } else {
245 coutW(Caching) << "The internal binning of variable " << convVar.GetName()
246 << " is not uniform. The numerical FFT will likely yield wrong results." << std::endl;
247 convVar.setBinning(varBinning, "cache");
248 }
249 }
250}
251
252
253////////////////////////////////////////////////////////////////////////////////
254/// Return base name component for cache components in this case 'PDF1_CONV_PDF2'
255
257{
258 static TString name ;
259 name = _pdf1.arg().GetName() ;
260 name.Append("_CONV_") ;
261 name.Append(_pdf2.arg().GetName()) ;
262 return name.Data() ;
263}
264
265
266
267
268////////////////////////////////////////////////////////////////////////////////
269/// Return specialized cache subclass for FFT calculations
270
272{
273 return new FFTCacheElem(*this,nset) ;
274}
275
276
277
278
279////////////////////////////////////////////////////////////////////////////////
280/// Clone input pdf and attach to dataset
281
283 PdfCacheElem(self,nsetIn),
284 fftr2c1(0),fftr2c2(0),fftc2r(0)
285{
286 RooAbsPdf* clonePdf1 = (RooAbsPdf*) self._pdf1.arg().cloneTree() ;
287 RooAbsPdf* clonePdf2 = (RooAbsPdf*) self._pdf2.arg().cloneTree() ;
288 clonePdf1->attachDataSet(*hist()) ;
289 clonePdf2->attachDataSet(*hist()) ;
290
291 // Shift observable
292 RooRealVar* convObs = (RooRealVar*) hist()->get()->find(self._x.arg().GetName()) ;
293
294 // Install FFT reference range
295 string refName = Form("refrange_fft_%s",self.GetName()) ;
296 convObs->setRange(refName.c_str(),convObs->getMin(),convObs->getMax()) ;
297
298 if (self._shift1!=0) {
299 RooLinearVar* shiftObs1 = new RooLinearVar(Form("%s_shifted_FFTBuffer1",convObs->GetName()),"shiftObs1",
300 *convObs,RooFit::RooConst(1),RooFit::RooConst(-1*self._shift1)) ;
301
302 RooArgSet clonedBranches1 ;
303 RooCustomizer cust(*clonePdf1,"fft") ;
304 cust.replaceArg(*convObs,*shiftObs1) ;
305
306 pdf1Clone = (RooAbsPdf*) cust.build() ;
307
308 pdf1Clone->addOwnedComponents(*shiftObs1) ;
309 pdf1Clone->addOwnedComponents(*clonePdf1) ;
310
311 } else {
312 pdf1Clone = clonePdf1 ;
313 }
314
315 if (self._shift2!=0) {
316 RooLinearVar* shiftObs2 = new RooLinearVar(Form("%s_shifted_FFTBuffer2",convObs->GetName()),"shiftObs2",
317 *convObs,RooFit::RooConst(1),RooFit::RooConst(-1*self._shift2)) ;
318
319 RooArgSet clonedBranches2 ;
320 RooCustomizer cust(*clonePdf2,"fft") ;
321 cust.replaceArg(*convObs,*shiftObs2) ;
322
323 pdf1Clone->addOwnedComponents(*shiftObs2) ;
324 pdf1Clone->addOwnedComponents(*clonePdf2) ;
325
326 pdf2Clone = (RooAbsPdf*) cust.build() ;
327
328 } else {
329 pdf2Clone = clonePdf2 ;
330 }
331
332
333 // Attach cloned pdf to all original parameters of self
334 RooArgSet* fftParams = self.getParameters(*convObs) ;
335
336 // Remove all cache histogram from fftParams as these
337 // observable need to remain attached to the histogram
338 fftParams->remove(*hist()->get(),kTRUE,kTRUE) ;
339
342 pdf1Clone->fixAddCoefRange(refName.c_str()) ;
343 pdf2Clone->fixAddCoefRange(refName.c_str()) ;
344
345 delete fftParams ;
346
347 // Save copy of original histX binning and make alternate binning
348 // for extended range scanning
349
350 const Int_t N = convObs->numBins();
351 if (N < 900) {
352 oocoutW(&self, Eval) << "The FFT convolution '" << self.GetName() << "' will run with " << N
353 << " bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number"
354 " of bins of the observable '" << convObs->GetName() << "'." << std::endl;
355 }
356 Int_t Nbuf = static_cast<Int_t>((N*self.bufferFraction())/2 + 0.5) ;
357 Double_t obw = (convObs->getMax() - convObs->getMin())/N ;
358 Int_t N2 = N+2*Nbuf ;
359
360 scanBinning = new RooUniformBinning (convObs->getMin()-Nbuf*obw,convObs->getMax()+Nbuf*obw,N2) ;
361 histBinning = convObs->getBinning().clone() ;
362
363 // Deactivate dirty state propagation on datahist observables
364 // and set all nodes on both pdfs to operMode AlwaysDirty
366 convObs->setOperMode(ADirty,kTRUE) ;
367}
368
369
370////////////////////////////////////////////////////////////////////////////////
371/// Suffix for cache histogram (added in addition to suffix for cache name)
372
374{
375 return TString(Form("_BufFrac%3.1f_BufStrat%d",_bufFrac,_bufStrat)) ;
376}
377
378
379
380////////////////////////////////////////////////////////////////////////////////
381/// Returns all RooAbsArg objects contained in the cache element
382
384{
386
387 ret.add(*pdf1Clone) ;
388 ret.add(*pdf2Clone) ;
389 if (pdf1Clone->ownedComponents()) {
390 ret.add(*pdf1Clone->ownedComponents()) ;
391 }
392 if (pdf2Clone->ownedComponents()) {
393 ret.add(*pdf2Clone->ownedComponents()) ;
394 }
395
396 return ret ;
397}
398
399
400////////////////////////////////////////////////////////////////////////////////
401
403{
404 delete fftr2c1 ;
405 delete fftr2c2 ;
406 delete fftc2r ;
407
408 delete pdf1Clone ;
409 delete pdf2Clone ;
410
411 delete histBinning ;
412 delete scanBinning ;
413
414}
415
416
417
418
419////////////////////////////////////////////////////////////////////////////////
420/// Fill the contents of the cache the FFT convolution output
421
423{
424 RooDataHist& cacheHist = *cache.hist() ;
425
426 ((FFTCacheElem&)cache).pdf1Clone->setOperMode(ADirty,kTRUE) ;
427 ((FFTCacheElem&)cache).pdf2Clone->setOperMode(ADirty,kTRUE) ;
428
429 // Determine if there other observables than the convolution observable in the cache
430 RooArgSet otherObs ;
431 RooArgSet(*cacheHist.get()).snapshot(otherObs) ;
432
433 RooAbsArg* histArg = otherObs.find(_x.arg().GetName()) ;
434 if (histArg) {
435 otherObs.remove(*histArg,kTRUE,kTRUE) ;
436 delete histArg ;
437 }
438
439 //cout << "RooFFTConvPdf::fillCacheObject() otherObs = " << otherObs << endl ;
440
441 // Handle trivial scenario -- no other observables
442 if (otherObs.getSize()==0) {
444 return ;
445 }
446
447 // Handle cases where there are other cache slices
448 // Iterator over available slice positions and fill each
449
450 // Determine number of bins for each slice position observable
451 Int_t n = otherObs.getSize() ;
452 Int_t* binCur = new Int_t[n+1] ;
453 Int_t* binMax = new Int_t[n+1] ;
454 Int_t curObs = 0 ;
455
456 RooAbsLValue** obsLV = new RooAbsLValue*[n] ;
457 TIterator* iter = otherObs.createIterator() ;
458 RooAbsArg* arg ;
459 Int_t i(0) ;
460 while((arg=(RooAbsArg*)iter->Next())) {
461 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(arg) ;
462 obsLV[i] = lvarg ;
463 binCur[i] = 0 ;
464 // coverity[FORWARD_NULL]
465 binMax[i] = lvarg->numBins(binningName())-1 ;
466 i++ ;
467 }
468 delete iter ;
469
470 Bool_t loop(kTRUE) ;
471 while(loop) {
472 // Set current slice position
473 for (Int_t j=0 ; j<n ; j++) { obsLV[j]->setBin(binCur[j],binningName()) ; }
474
475// cout << "filling slice: bin of obsLV[0] = " << obsLV[0]->getBin() << endl ;
476
477 // Fill current slice
478 fillCacheSlice((FFTCacheElem&)cache,otherObs) ;
479
480 // Determine which iterator to increment
481 while(binCur[curObs]==binMax[curObs]) {
482
483 // Reset current iterator and consider next iterator ;
484 binCur[curObs]=0 ;
485 curObs++ ;
486
487 // master termination condition
488 if (curObs==n) {
489 loop=kFALSE ;
490 break ;
491 }
492 }
493
494 // Increment current iterator
495 binCur[curObs]++ ;
496 curObs=0 ;
497
498 }
499
500 delete[] obsLV ;
501 delete[] binMax ;
502 delete[] binCur ;
503
504}
505
506
507////////////////////////////////////////////////////////////////////////////////
508/// Fill a slice of cachePdf with the output of the FFT convolution calculation
509
510void RooFFTConvPdf::fillCacheSlice(FFTCacheElem& aux, const RooArgSet& slicePos) const
511{
512 // Extract histogram that is the basis of the RooHistPdf
513 RooDataHist& cacheHist = *aux.hist() ;
514
515 // Sample array of input points from both pdfs
516 // Note that returned arrays have optional buffers zones below and above range ends
517 // to reduce cyclical effects and have been cyclically rotated so that bin containing
518 // zero value is at position zero. Example:
519 //
520 // original: -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5
521 // add buffer zones: U U -5 -4 -3 -2 -1 0 +1 +2 +3 +4 +5 O O
522 // rotate: 0 +1 +2 +3 +4 +5 O O U U -5 -4 -3 -2 -1
523 //
524 //
525
526 Int_t N,N2,binShift1,binShift2 ;
527
528 RooRealVar* histX = (RooRealVar*) cacheHist.get()->find(_x.arg().GetName()) ;
529 if (_bufStrat==Extend) histX->setBinning(*aux.scanBinning) ;
530 Double_t* input1 = scanPdf((RooRealVar&)_x.arg(),*aux.pdf1Clone,cacheHist,slicePos,N,N2,binShift1,_shift1) ;
531 Double_t* input2 = scanPdf((RooRealVar&)_x.arg(),*aux.pdf2Clone,cacheHist,slicePos,N,N2,binShift2,_shift2) ;
532 if (_bufStrat==Extend) histX->setBinning(*aux.histBinning) ;
533
534
535
536
537 // Retrieve previously defined FFT transformation plans
538 if (!aux.fftr2c1) {
539 aux.fftr2c1 = TVirtualFFT::FFT(1, &N2, "R2CK");
540 aux.fftr2c2 = TVirtualFFT::FFT(1, &N2, "R2CK");
541 aux.fftc2r = TVirtualFFT::FFT(1, &N2, "C2RK");
542 }
543
544 // Real->Complex FFT Transform on p.d.f. 1 sampling
545 aux.fftr2c1->SetPoints(input1);
546 aux.fftr2c1->Transform();
547
548 // Real->Complex FFT Transform on p.d.f 2 sampling
549 aux.fftr2c2->SetPoints(input2);
550 aux.fftr2c2->Transform();
551
552 // Loop over first half +1 of complex output results, multiply
553 // and set as input of reverse transform
554 for (Int_t i=0 ; i<N2/2+1 ; i++) {
555 Double_t re1,re2,im1,im2 ;
556 aux.fftr2c1->GetPointComplex(i,re1,im1) ;
557 aux.fftr2c2->GetPointComplex(i,re2,im2) ;
558 Double_t re = re1*re2 - im1*im2 ;
559 Double_t im = re1*im2 + re2*im1 ;
560 TComplex t(re,im) ;
561 aux.fftc2r->SetPointComplex(i,t) ;
562 }
563
564 // Reverse Complex->Real FFT transform product
565 aux.fftc2r->Transform() ;
566
567 Int_t totalShift = binShift1 + (N2-N)/2 ;
568
569 // Store FFT result in cache
570
571 TIterator* iter = const_cast<RooDataHist&>(cacheHist).sliceIterator(const_cast<RooAbsReal&>(_x.arg()),slicePos) ;
572 for (Int_t i =0 ; i<N ; i++) {
573
574 // Cyclically shift array back so that bin containing zero is back in zeroBin
575 Int_t j = i + totalShift ;
576 while (j<0) j+= N2 ;
577 while (j>=N2) j-= N2 ;
578
579 iter->Next() ;
580 cacheHist.set(aux.fftc2r->GetPointReal(j)) ;
581 }
582 delete iter ;
583
584 // cacheHist.dump2() ;
585
586 // Delete input arrays
587 delete[] input1 ;
588 delete[] input2 ;
589
590}
591
592
593////////////////////////////////////////////////////////////////////////////////
594/// Scan the values of 'pdf' in observable 'obs' using the bin values stored in 'hist' at slice position 'slicePos'
595/// N is filled with the number of bins defined in hist, N2 is filled with N plus the number of buffer bins
596/// The return value is an array of doubles of length N2 with the sampled values. The caller takes ownership
597/// of the array
598
599Double_t* RooFFTConvPdf::scanPdf(RooRealVar& obs, RooAbsPdf& pdf, const RooDataHist& hist, const RooArgSet& slicePos,
600 Int_t& N, Int_t& N2, Int_t& zeroBin, Double_t shift) const
601{
602
603 RooRealVar* histX = (RooRealVar*) hist.get()->find(obs.GetName()) ;
604
605 // Calculate number of buffer bins on each size to avoid cyclical flow
606 N = histX->numBins(binningName()) ;
607 Int_t Nbuf = static_cast<Int_t>((N*bufferFraction())/2 + 0.5) ;
608 N2 = N+2*Nbuf ;
609
610
611 // Allocate array of sampling size plus optional buffer zones
612 Double_t* array = new Double_t[N2] ;
613
614 // Set position of non-convolution observable to that of the cache slice that were are processing now
615 hist.get(slicePos) ;
616
617 // Find bin ID that contains zero value
618 zeroBin = 0 ;
619 if (histX->getMax()>=0 && histX->getMin()<=0) {
620 zeroBin = histX->getBinning().binNumber(0) ;
621 } else if (histX->getMin()>0) {
622 Double_t bw = (histX->getMax() - histX->getMin())/N2 ;
623 zeroBin = Int_t(-histX->getMin()/bw) ;
624 } else {
625 Double_t bw = (histX->getMax() - histX->getMin())/N2 ;
626 zeroBin = Int_t(-1*histX->getMax()/bw) ;
627 }
628
629 Int_t binShift = Int_t((N2* shift) / (histX->getMax()-histX->getMin())) ;
630
631 zeroBin += binShift ;
632 while(zeroBin>=N2) zeroBin-= N2 ;
633 while(zeroBin<0) zeroBin+= N2 ;
634
635 // First scan hist into temp array
636 Double_t *tmp = new Double_t[N2] ;
637 Int_t k(0) ;
638 switch(_bufStrat) {
639
640 case Extend:
641 // Sample entire extended range (N2 samples)
642 for (k=0 ; k<N2 ; k++) {
643 histX->setBin(k) ;
644 tmp[k] = pdf.getVal(hist.get()) ;
645 }
646 break ;
647
648 case Flat:
649 // Sample original range (N samples) and fill lower and upper buffer
650 // bins with p.d.f. value at respective boundary
651 {
652 histX->setBin(0) ;
653 Double_t val = pdf.getVal(hist.get()) ;
654 for (k=0 ; k<Nbuf ; k++) {
655 tmp[k] = val ;
656 }
657 for (k=0 ; k<N ; k++) {
658 histX->setBin(k) ;
659 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
660 }
661 histX->setBin(N-1) ;
662 val = pdf.getVal(hist.get()) ;
663 for (k=0 ; k<Nbuf ; k++) {
664 tmp[N+Nbuf+k] = val ;
665 }
666 }
667 break ;
668
669 case Mirror:
670 // Sample original range (N samples) and fill lower and upper buffer
671 // bins with mirror image of sampled range
672 for (k=0 ; k<N ; k++) {
673 histX->setBin(k) ;
674 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
675 }
676 for (k=1 ; k<=Nbuf ; k++) {
677 histX->setBin(k) ;
678 tmp[Nbuf-k] = pdf.getVal(hist.get()) ;
679 histX->setBin(N-k) ;
680 tmp[Nbuf+N+k-1] = pdf.getVal(hist.get()) ;
681 }
682 break ;
683 }
684
685 // Scan function and store values in array
686 for (Int_t i=0 ; i<N2 ; i++) {
687 // Cyclically shift writing location by zero bin position
688 Int_t j = i - (zeroBin) ;
689 if (j<0) j+= N2 ;
690 if (j>=N2) j-= N2 ;
691 array[i] = tmp[j] ;
692 }
693
694 // Cleanup
695 delete[] tmp ;
696 return array ;
697}
698
699
700
701////////////////////////////////////////////////////////////////////////////////
702/// Return the observables to be cached given the normalization set nset
703///
704/// If the cache observables is in nset then this is
705/// - the convolution observable plus
706/// - any member of nset that is either a RooCategory,
707/// - or was previously specified through setCacheObservables().
708///
709/// In case the cache observable is _not_ in nset, then it is
710/// - the convolution observable plus
711/// - all member of nset are observables of this p.d.f.
712///
713
715{
716 // Get complete list of observables
717 RooArgSet* obs1 = _pdf1.arg().getObservables(nset) ;
718 RooArgSet* obs2 = _pdf2.arg().getObservables(nset) ;
719 obs1->add(*obs2,kTRUE) ;
720
721 // Check if convolution observable is in nset
722 if (nset.contains(_x.arg())) {
723
724 // Now strip out all non-category observables
725 TIterator* iter = obs1->createIterator() ;
726 RooAbsArg* arg ;
727 RooArgSet killList ;
728 while((arg=(RooAbsArg*)iter->Next())) {
729 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
730 killList.add(*arg) ;
731 }
732 }
733 delete iter ;
734 obs1->remove(killList) ;
735
736 // And add back the convolution observables
737 obs1->add(_x.arg(),kTRUE) ;
738
739 obs1->add(_cacheObs) ;
740
741 delete obs2 ;
742
743 } else {
744
745 // If cacheObs was filled, cache only observables in there
746 if (_cacheObs.getSize()>0) {
747 TIterator* iter = obs1->createIterator() ;
748 RooAbsArg* arg ;
749 RooArgSet killList ;
750 while((arg=(RooAbsArg*)iter->Next())) {
751 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
752 killList.add(*arg) ;
753 }
754 }
755 delete iter ;
756 obs1->remove(killList) ;
757 }
758
759
760 // Make sure convolution observable is always in there
761 obs1->add(_x.arg(),kTRUE) ;
762 delete obs2 ;
763
764 }
765
766 return obs1 ;
767}
768
769
770
771////////////////////////////////////////////////////////////////////////////////
772/// Return the parameters on which the cache depends given normalization
773/// set nset. For this p.d.f these are the parameters of the input p.d.f.
774/// but never the convolution variable, it case it is not part of nset
775
777{
778 RooArgSet* vars = getVariables() ;
779 RooArgSet* obs = actualObservables(nset) ;
780 vars->remove(*obs) ;
781 delete obs ;
782
783 return vars ;
784}
785
786
787
788////////////////////////////////////////////////////////////////////////////////
789/// Return p.d.f. observable (which can be a function) to substitute given
790/// p.d.f. observable. Substitute x by xprime if xprime is set
791
793{
794 if (_xprime.absArg() && string(histObservable.GetName())==_x.absArg()->GetName()) {
795 return (*_xprime.absArg()) ;
796 }
797 return histObservable ;
798}
799
800
801
802////////////////////////////////////////////////////////////////////////////////
803/// Create appropriate generator context for this convolution. If both input p.d.f.s support
804/// internal generation, if it is safe to use them and if no observables other than the convolution
805/// observable are requested for generation, use the specialized convolution generator context
806/// which implements a smearing strategy in the convolution observable. If not return the
807/// regular accept/reject generator context
808
810 const RooArgSet* auxProto, Bool_t verbose) const
811{
812 RooArgSet vars2(vars) ;
813 vars2.remove(_x.arg(),kTRUE,kTRUE) ;
814 Int_t numAddDep = vars2.getSize() ;
815
817 Bool_t pdfCanDir = (((RooAbsPdf&)_pdf1.arg()).getGenerator(_x.arg(),dummy) != 0 && \
819 Bool_t resCanDir = (((RooAbsPdf&)_pdf2.arg()).getGenerator(_x.arg(),dummy) !=0 &&
821
822 if (pdfCanDir) {
823 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f " << _pdf1.arg().GetName()
824 << " has internal generator that is safe to use in current context" << endl ;
825 }
826 if (resCanDir) {
827 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f. " << _pdf2.arg().GetName()
828 << " has internal generator that is safe to use in current context" << endl ;
829 }
830 if (numAddDep>0) {
831 cxcoutI(Generation) << "RooFFTConvPdf::genContext() generation requested for observables other than the convolution observable " << _x.arg().GetName() << endl ;
832 }
833
834
835 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
836 // Any resolution model with more dependents than the convolution variable
837 // or pdf or resmodel do not support direct generation
838 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting accept/reject generator context because one or both of the input "
839 << "p.d.f.s cannot use internal generator and/or "
840 << "observables other than the convolution variable are requested for generation" << endl ;
841 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
842 }
843
844 // Any other resolution model: use specialized generator context
845 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting specialized convolution generator context as both input "
846 << "p.d.fs are safe for internal generator and only "
847 << "the convolution observables is requested for generation" << endl ;
848 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
849}
850
851
852
853////////////////////////////////////////////////////////////////////////////////
854/// Change the size of the buffer on either side of the observable range to frac times the
855/// size of the range of the convolution observable
856
858{
859 if (frac<0) {
860 coutE(InputArguments) << "RooFFTConvPdf::setBufferFraction(" << GetName() << ") fraction should be greater than or equal to zero" << endl ;
861 return ;
862 }
863 _bufFrac = frac ;
864
865 // Sterilize the cache as certain partial results depend on buffer fraction
867}
868
869
870////////////////////////////////////////////////////////////////////////////////
871/// Change strategy to fill the overflow buffer on either side of the convolution observable range.
872///
873/// 'Extend' means is that the input p.d.f convolution observable range is widened to include the buffer range
874/// 'Flat' means that the buffer is filled with the p.d.f. value at the boundary of the observable range
875/// 'Mirror' means that the buffer is filled with a ,irror image of the p.d.f. around the convolution observable boundary
876///
877/// The default strategy is extend. If one of the input p.d.f.s is a RooAddPdf, it is configured so that the interpretation
878/// range of the fraction coefficients is kept at the nominal convolutions observable range (instead of interpreting coefficients
879/// in the widened range including the buffer)
880
882{
883 _bufStrat = bs ;
884}
885
886
887
888////////////////////////////////////////////////////////////////////////////////
889/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
890/// product operator construction
891
892void RooFFTConvPdf::printMetaArgs(ostream& os) const
893{
894 os << _pdf1.arg().GetName() << "(" << _x.arg().GetName() << ") (*) " << _pdf2.arg().GetName() << "(" << _x.arg().GetName() << ") " ;
895}
896
897
898
899////////////////////////////////////////////////////////////////////////////////
900/// (Re)calculate effective parameters of this p.d.f.
901
903{
904 RooArgSet* params1 = _pdf1.arg().getParameters(_x.arg()) ;
905 RooArgSet* params2 = _pdf2.arg().getParameters(_x.arg()) ;
907 _params.add(*params1) ;
908 _params.add(*params2,kTRUE) ;
909 delete params1 ;
910 delete params2 ;
911}
912
913
914
915////////////////////////////////////////////////////////////////////////////////
916///calcParams() ;
917
918Bool_t RooFFTConvPdf::redirectServersHook(const RooAbsCollection& /*newServerList*/, Bool_t /*mustReplaceAll*/, Bool_t /*nameChange*/, Bool_t /*isRecursive*/)
919{
920 return kFALSE ;
921}
void Class()
Definition: Class.C:29
static RooMathCoreReg dummy
#define coutI(a)
Definition: RooMsgService.h:31
#define cxcoutI(a)
Definition: RooMsgService.h:83
#define oocoutW(o, a)
Definition: RooMsgService.h:46
#define coutW(a)
Definition: RooMsgService.h:33
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
#define N
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
virtual RooAbsArg * cloneTree(const char *newname=0) const
Clone tree expression of objects.
Definition: RooAbsArg.cxx:2288
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
friend class RooArgSet
Definition: RooAbsArg.h:471
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:533
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2074
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Change cache operation mode to given mode.
Definition: RooAbsArg.cxx:1746
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2274
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Definition: RooAbsArg.cxx:1080
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
Definition: RooAbsArg.cxx:1489
RooAbsBinning is the abstract base class for RooRealVar binning definitions This class defines the in...
Definition: RooAbsBinning.h:26
virtual RooAbsBinning * clone(const char *name=0) const =0
Int_t numBins() const
Definition: RooAbsBinning.h:37
virtual Bool_t isUniform() const
Definition: RooAbsBinning.h:48
virtual Double_t highBound() const =0
virtual Int_t binNumber(Double_t x) const =0
virtual Double_t lowBound() const =0
virtual RooArgList containedArgs(Action)
Returns all RooAbsArg objects contained in the cache element.
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
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
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.
Definition: RooAbsData.cxx:340
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
Abstract base class for objects that are lvalues, i.e.
Definition: RooAbsLValue.h:26
virtual Int_t numBins(const char *rangeName=0) const =0
virtual void setBin(Int_t ibin, const char *rangeName=0)=0
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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...
Definition: RooAbsPdf.cxx:2047
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
Definition: RooAbsPdf.cxx:2082
virtual Double_t getMax(const char *name=0) const
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
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
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 *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
RooAbsArg * absArg() const
Definition: RooArgProxy.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
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...
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
void set(Double_t weight, Double_t wgtErr=-1)
Increment the weight of the bin enclosing the coordinates given by 'row' by the specified amount.
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
virtual RooArgList containedArgs(Action)
Returns all RooAbsArg objects contained in the cache element.
RooAbsBinning * histBinning
Definition: RooFFTConvPdf.h:93
FFTCacheElem(const RooFFTConvPdf &self, const RooArgSet *nset)
Clone input pdf and attach to dataset.
RooAbsBinning * scanBinning
Definition: RooFFTConvPdf.h:94
PDF for the numerical (FFT) convolution of two PDFs.
Definition: RooFFTConvPdf.h:26
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
Definition: RooFFTConvPdf.h:72
BufStrat _bufStrat
void calcParams()
(Re)calculate effective parameters of this p.d.f.
Double_t _shift2
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
Definition: RooFFTConvPdf.h:69
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
Definition: RooFFTConvPdf.h:70
RooRealProxy _x
Definition: RooFFTConvPdf.h:68
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
Definition: RooFFTConvPdf.h:71
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.
Double_t _shift1
virtual RooArgSet * actualParameters(const RooArgSet &nset) const
Return the parameters on which the cache depends given normalization set nset.
Double_t bufferFraction() const
Definition: RooFFTConvPdf.h:42
friend class FFTCacheElem
Definition: RooFFTConvPdf.h:98
RooSetProxy _cacheObs
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
Definition: RooGenContext.h:30
RooLinearVar is the most general form of a derived real-valued object that can be used by RooRealInte...
Definition: RooLinearVar.h:29
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
Bool_t hasBinning(const char *name) const
Returns true if variable has a binning with 'name'.
Definition: RooRealVar.cxx:253
void setRange(const char *name, Double_t min, Double_t max)
Set range named 'name to [min,max].
Definition: RooRealVar.cxx:448
const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const
Return binning definition with name.
Definition: RooRealVar.cxx:266
void setBinning(const RooAbsBinning &binning, const char *name=0)
Add given binning under name 'name' with this variable.
Definition: RooRealVar.cxx:345
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
virtual void removeAll()
Remove all argument inset using remove(const RooAbsArg&).
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:443
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
const Int_t n
Definition: legend1.C:16
@ Generation
Definition: RooGlobalFunc.h:57
@ InputArguments
Definition: RooGlobalFunc.h:58
RooConstVar & RooConst(Double_t val)
STL namespace.
auto * a
Definition: textangle.C:12