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