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