Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TF1Convolution.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Authors: Lorenzo Moneta, Aurélie Flandi 27/08/14
3//
4/**********************************************************************
5 * *
6 * Copyright (c) 2015 ROOT Team, CERN/PH-SFT *
7 * *
8 * *
9 **********************************************************************/
10
11#include <memory>
12
13#include "TF1Convolution.h"
14#include "TROOT.h"
15#include "TObject.h"
16#include "TObjString.h"
17#include "TObjArray.h"
18#include "TMath.h"
19#include "Math/Integrator.h"
25#include "Math/Functor.h"
26#include "TVirtualFFT.h"
27
28/** \class TF1Convolution
29 \ingroup Functions
30 \brief Class wrapping convolution of two functions
31
32Class wrapping convolution of two functions: evaluation of \f$\int f(x)g(x-t)dx\f$
33
34The convolution is performed by default using FFTW if it is available .
35One can pass optionally the range of the convolution (by default the first function range is used).
36Note that when using Discrete Fourier Transform (as FFTW), it is a circular transform, so the functions should be
37approximately zero at the end of the range. If they are significantly different than zero on one side (e.g. the left side)
38a spill over will occur on the other side (e.g right side).
39If no function range is given by default the function1 range + 10% is used
40One should use also a not too small number of points for the DFT (a minimum of 1000). By default 10000 points are used.
41*/
42
44
46
48{
52
53public:
55 fFunc1(&f1),
56 fFunc2(&f2),
57 fT0(t)
58 {}
60 {
61 // use EvalPar that is faster
62 Double_t xx[2];
63 xx[0] = x;
64 xx[1] = fT0-x;
65 return fFunc1->EvalPar(xx,nullptr) * fFunc2->EvalPar(xx+1,nullptr);
66 }
67};
68
69////////////////////////////////////////////////////////////////////////////////
70/// Internal function to initialize data members.
71/// Use TF1::Copy instead of Clone.
72
73void TF1Convolution::InitializeDataMembers(TF1* function1, TF1* function2, Bool_t useFFT)
74{
75 if (function1) {
76 // functions must be 1d- if not flag an error
77 if (function1->GetNdim() != 1)
78 Error("InitializeDataMembers","function1 %s is not of dimension 1 ",function1->GetName());
79 //TF1 * fnew1 = (TF1*) function1->IsA()->New();
80 // since function1 is a TF1 (cannot be a derived class) we can instantiate it directly
81 fFunction1 = std::make_unique<TF1> ();
82 function1->Copy(*fFunction1);
83 }
84 if (function2) {
85 if (function2->GetNdim() != 1)
86 Error("InitializeDataMembers","function2 %s is not of dimension 1 ",function2->GetName());
87 //TF1 * fnew2 = (TF1*) function2->IsA()->New();
88 fFunction2 = std::make_unique<TF1>();
89 function2->Copy(*fFunction2);
90 }
91 if (fFunction1.get() == nullptr|| fFunction2.get() == nullptr)
92 Fatal("InitializeDataMembers","Invalid functions - Abort");
93
94 // Set kNotGlobal bit
97
98 // use by default range of first function
99 fFunction1->GetRange(fXmin, fXmax);
100 // when using FFT add by default an extra 10% on each side
101 if (useFFT) {
103 }
104 fNofParams1 = fFunction1->GetNpar();
105 fNofParams2 = fFunction2->GetNpar();
106 fParams1 = std::vector<Double_t>(fNofParams1);
107 fParams2 = std::vector<Double_t>(fNofParams2);
108 fCstIndex = (fFunction1->GetParNumber("Constant") == -1)
109 ? -1
110 : fFunction2->GetParNumber("Constant"); // TODO: add dropConstantParam flag?
111 fFlagFFT = useFFT;
112 fFlagGraph = false;
113 fNofPoints = 10000;
114
116 for (int i=0; i<fNofParams1; i++)
117 {
118 fParams1[i] = fFunction1 -> GetParameter(i);
119 fParNames.push_back(fFunction1 -> GetParName(i) );
120 }
121 for (int i=0; i<fNofParams2; i++)
122 {
123 fParams2[i] = fFunction2 -> GetParameter(i);
124 if (i != fCstIndex) fParNames.push_back(fFunction2 -> GetParName(i) );
125 }
126 if (fCstIndex!=-1)
127 {
128 fFunction2 -> FixParameter(fCstIndex,1.);
130 fParams2.erase(fParams2.begin()+fCstIndex);
131 }
132}
133
134////////////////////////////////////////////////////////////////////////////////
135/// constructor without arguments.
136
138{
139 // Nothing to do
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// constructor from the two function pointer and a flag is using FFT.
144
145TF1Convolution::TF1Convolution(TF1* function1, TF1* function2, Bool_t useFFT)
146{
147 InitializeDataMembers(function1,function2, useFFT);
148}
149
150////////////////////////////////////////////////////////////////////////////////
151/// Constructor from the two function pointer and the convolution range.
152
154{
155 InitializeDataMembers(function1, function2,useFFT);
156 if (xmin < xmax) {
157 fXmin = xmin;
158 fXmax = xmax;
160 } else {
161 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
163 }
164}
165
166////////////////////////////////////////////////////////////////////////////////
167/// Constructor from a formula expression as f1 * f2 where f1 and f2 are two functions known to ROOT.
168
170{
172
173 TObjArray *objarray = formula.Tokenize("*");
174 std::vector < TString > stringarray(2);
175 std::vector < TF1* > funcarray(2);
176 for (int i=0; i<2; i++)
177 {
178 stringarray[i] = ((TObjString*)((*objarray)[i])) -> GetString();
179 stringarray[i].ReplaceAll(" ","");
180 funcarray[i] = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(stringarray[i]));
181 // case function is not found try to use as a TFormula
182 if (funcarray[i] == nullptr) {
183 TF1 * f = new TF1(TString::Format("f_conv_%d",i+1),stringarray[i]);
184 if (!f->GetFormula()->IsValid() )
185 Error("TF1Convolution","Invalid formula : %s",stringarray[i].Data() );
186 if (i == 0)
187 fFunction1 = std::unique_ptr<TF1>(f);
188 else
189 fFunction2 = std::unique_ptr<TF1>(f);
190 }
191 }
192 InitializeDataMembers(funcarray[0], funcarray[1],useFFT);
193 if (xmin < xmax) {
194 fXmin = xmin;
195 fXmax = xmax;
197 } else {
198 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
200 }
201}
202
203////////////////////////////////////////////////////////////////////////////////
204/// Constructor from 2 function names where f1 and f2 are two functions known to
205/// ROOT.
206///
207/// If the function names are not known to ROOT, tries to interpret them as
208/// TFormula.
210{
212 (TString)formula1.ReplaceAll(" ","");
213 (TString)formula2.ReplaceAll(" ","");
214 TF1* f1 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula1));
215 TF1* f2 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula2));
216 // if function do not exists try using TFormula
217 if (!f1) {
218 fFunction1 = std::make_unique<TF1>("f_conv_1", formula1);
219 if (!fFunction1->GetFormula()->IsValid() )
220 Error("TF1Convolution","Invalid formula for : %s",formula1.Data() );
221 }
222 if (!f2) {
223 fFunction2 = std::make_unique<TF1>("f_conv_1", formula2);
224 if (!fFunction2->GetFormula()->IsValid() )
225 Error("TF1Convolution","Invalid formula for : %s",formula2.Data() );
226 }
227 // if f1 or f2 are null ptr are not used in InitializeDataMembers
228 InitializeDataMembers(f1, f2,useFFT);
229 if (xmin < xmax) {
230 fXmin = xmin;
231 fXmax = xmax;
232 } else {
233 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
235 }
236}
237
238////////////////////////////////////////////////////////////////////////////////
239/// Copy constructor (necessary to hold unique_ptr as member variable).
240
242{
243 conv.TF1Convolution::Copy(*this);
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Operator =
248
250{
251 if (this != &rhs)
252 rhs.TF1Convolution::Copy(*this);
253 return *this;
254}
255
256////////////////////////////////////////////////////////////////////////////////
257/// Perform the FFT of the two functions.
258
260{
261 if (gDebug)
262 Info("MakeFFTConv","Making FFT convolution using %d points in range [%g,%g]",fNofPoints,fXmin,fXmax);
263
264 std::vector < Double_t > x (fNofPoints);
265 std::vector < Double_t > in1(fNofPoints);
266 std::vector < Double_t > in2(fNofPoints);
267
268 TVirtualFFT *fft1 = TVirtualFFT::FFT(1, &fNofPoints, "R2C K");
269 TVirtualFFT *fft2 = TVirtualFFT::FFT(1, &fNofPoints, "R2C K");
270 if (fft1 == nullptr || fft2 == nullptr) {
271 Warning("MakeFFTConv","Cannot use FFT, probably FFTW package is not available. Switch to numerical convolution");
272 fFlagFFT = false;
273 return;
274 }
275
276 // apply a shift in order to have the second function centered around middle of the range of the convolution
277 Double_t shift2 = 0.5*(fXmin+fXmax);
278 Double_t x2;
279 for (int i=0; i<fNofPoints; i++)
280 {
281 x[i] = fXmin + (fXmax-fXmin)/(fNofPoints-1)*i;
282 x2 = x[i] - shift2;
283 in1[i] = fFunction1 -> EvalPar( &x[i], nullptr);
284 in2[i] = fFunction2 -> EvalPar( &x2, nullptr);
285 fft1 -> SetPoint(i, in1[i]);
286 fft2 -> SetPoint(i, in2[i]);
287 }
288 fft1 -> Transform();
289 fft2 -> Transform();
290
291 //inverse transformation of the product
292
293 TVirtualFFT *fftinverse = TVirtualFFT::FFT(1, &fNofPoints, "C2R K");
294 Double_t re1, re2, im1, im2, out_re, out_im;
295
296 for (int i=0;i<=fNofPoints/2.;i++)
297 {
298 fft1 -> GetPointComplex(i,re1,im1);
299 fft2 -> GetPointComplex(i,re2,im2);
300 out_re = re1*re2 - im1*im2;
301 out_im = re1*im2 + re2*im1;
302 fftinverse -> SetPoint(i, out_re, out_im);
303 }
304 fftinverse -> Transform();
305
306 // fill a graph with the result of the convolution
307 if (!fGraphConv)
308 fGraphConv = std::make_unique<TGraph>(fNofPoints);
309
310 for (int i=0;i<fNofPoints;i++)
311 {
312 // we need this since we have applied a shift in the middle of f2
313 int j = i + fNofPoints/2;
314 if (j >= fNofPoints) j -= fNofPoints;
315 // need to normalize by dividing by the number of points and multiply by the bin width = Range/Number of points
316 fGraphConv->SetPoint(i, x[i], fftinverse->GetPointReal(j)*(fXmax-fXmin)/(fNofPoints*fNofPoints) );
317 }
318 fGraphConv->SetBit(TGraph::kIsSortedX); // indicate that points are sorted in X to speed up TGraph::Eval
319 fFlagGraph = true; // we can use the graph
320
321 // delete the fft objects
322 delete fft1;
323 delete fft2;
324 delete fftinverse;
325}
326
327////////////////////////////////////////////////////////////////////////////////
328/// Perform FFT convolution.
329
331{
332 if (!fFlagGraph) MakeFFTConv();
333 // if cannot make FFT use numconv
334 if (fGraphConv)
335 return fGraphConv -> Eval(t);
336 else
337
338 return EvalNumConv(t);
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// Perform numerical convolution.
343///
344
346{
347 /// Could in principle cache the integral in a Graph as it is done for the FFTW
349 Double_t result = 0;
350
352 if (fXmin != - TMath::Infinity() && fXmax != TMath::Infinity() )
353 result = integrator.Integral(fXmin, fXmax);
354 else if (fXmin == - TMath::Infinity() && fXmax != TMath::Infinity() )
355 result = integrator.IntegralLow(fXmax);
356 else if (fXmin != - TMath::Infinity() && fXmax == TMath::Infinity() )
357 result = integrator.IntegralUp(fXmin);
358 else if (fXmin == - TMath::Infinity() && fXmax == TMath::Infinity() )
359 result = integrator.Integral();
360
361 return result;
362}
363
364////////////////////////////////////////////////////////////////////////////////
365/// Used in TF1 when doing the fit, will be evaluated at each point.
366
368{
369 if (p!=nullptr) TF1Convolution::SetParameters(p); // first refresh the parameters
370
371 Double_t result = 0.;
372 if (fFlagFFT)
373 result = EvalFFTConv(x[0]);
374 else
375 result = EvalNumConv(x[0]);
376 return result;
377}
378
379////////////////////////////////////////////////////////////////////////////////
380/// Set the number of points used for the FFT convolution.
381
383{
384 if (n<0) return;
385 fNofPoints = n;
386 if (fGraphConv) fGraphConv -> Set(fNofPoints);
387 fFlagGraph = false; // to indicate we need to re-do the graph
388}
389
390////////////////////////////////////////////////////////////////////////////////
391/// Set the vector of parameters p for the convolution function g(x,p) = f1 * f2.
392
394{
395 bool equalParams = true;
396 for (int i=0; i<fNofParams1; i++) {
397 fFunction1->SetParameter(i, params[i]);
398 equalParams &= (fParams1[i] == params[i]);
399 fParams1[i] = params[i];
400 }
401 Int_t k = 0;
402 Int_t offset = 0;
403 Int_t offset2 = 0;
404 if (fCstIndex!=-1) offset = 1;
405 Int_t totalnofparams = fNofParams1+fNofParams2+offset;
406 for (int i=fNofParams1; i<totalnofparams; i++) {
407 if (k==fCstIndex)
408 {
409 k++;
410 offset2=1;
411 continue;
412 }
413 fFunction2->SetParameter(k, params[i - offset2]);
414 equalParams &= (fParams2[k - offset2] == params[i - offset2]);
415 fParams2[k - offset2] = params[i - offset2];
416 k++;
417 }
418
419 if (!equalParams) fFlagGraph = false; // to indicate we need to re-do the convolution
420}
421
422////////////////////////////////////////////////////////////////////////////////
423/// Set the parameter values for the convolution function.
424
426 Double_t p4, Double_t p5, Double_t p6, Double_t p7)
427{
428 Double_t params[]={p0,p1,p2,p3,p4,p5,p6,p7};
430}
431
432////////////////////////////////////////////////////////////////////////////////
433/// Set the fraction of extra range used when doing an FFT convolution.
434/// The extra range is often needed to avoid mirroring effect of the resulting convolution
435/// function at the borders.
436/// By default an extra range of 0.1 is used.
437
439{
440 if (percentage<0) return;
441 double range = fXmax - fXmin;
442 fXmin -= percentage * range;
443 fXmax += percentage * range;
444 fFlagGraph = false; // to indicate we need to re-do the convolution
445}
446
447////////////////////////////////////////////////////////////////////////////////
448/// Set the actual range used for the convolution.
449/// In case a or b are -inf or +inf and FFT convolution is used, then the
450/// range of the first function will be used and extended by the default extra range fraction.
451
453{
454 if (a >= b) {
455 Warning("SetRange", "Invalid range: %f >= %f", a, b);
456 return;
457 }
458
459 fXmin = a;
460 fXmax = b;
461 if (fFlagFFT && ( a==-TMath::Infinity() || b==TMath::Infinity() ) )
462 {
463 Warning("TF1Convolution::SetRange()","In FFT mode, range can not be infinite. Infinity has been replaced by range of first function plus a bufferzone to avoid spillover.");
464 if (a ==-TMath::Infinity()) fXmin = fFunction1 -> GetXmin();
465 if ( b== TMath::Infinity()) fXmax = fFunction1 -> GetXmax();
466 // add a spill over of 10% in this case
468 }
469 fFlagGraph = false; // to indicate we need to re-do the convolution
470}
471
472////////////////////////////////////////////////////////////////////////////////
473/// Set the default extra range fraction used when doing a FFT convolution.
474/// By default the value is 0.1 (10%).
475/// The function return the previous default defined value.
476
478{
479 Double_t prevValue = fgExtraRangeFraction;
480 fgExtraRangeFraction = fraction;
481 return prevValue;
482}
483
484////////////////////////////////////////////////////////////////////////////////
485/// Get the range used for the convolution.
486
488{
489 a = fXmin;
490 b = fXmax;
491}
492
493////////////////////////////////////////////////////////////////////////////////
494/// Update the two component functions of the convolution.
495
497{
498 fFunction1->Update();
499 fFunction2->Update();
500}
501
502////////////////////////////////////////////////////////////////////////////////
503
505{
506 // copy numbers
507 ((TF1Convolution &)obj).fXmin = fXmin;
508 ((TF1Convolution &)obj).fXmax = fXmax;
509 ((TF1Convolution &)obj).fNofParams1 = fNofParams1;
510 ((TF1Convolution &)obj).fNofParams2 = fNofParams2;
511 ((TF1Convolution &)obj).fCstIndex = fCstIndex;
512 ((TF1Convolution &)obj).fNofPoints = fNofPoints;
513 ((TF1Convolution &)obj).fFlagFFT = fFlagFFT;
514 ((TF1Convolution &)obj).fFlagGraph = false; // since we're not copying the graph
515
516 // copy vectors
517 ((TF1Convolution &)obj).fParams1 = fParams1;
518 ((TF1Convolution &)obj).fParams2 = fParams2;
519 ((TF1Convolution &)obj).fParNames = fParNames;
520
521 // we need to copy the content of the unique_ptr's
522 ((TF1Convolution &)obj).fFunction1 = std::make_unique<TF1>();
523 ((TF1Convolution &)obj).fFunction2 = std::make_unique<TF1>();
524 fFunction1->Copy(*(((TF1Convolution &)obj).fFunction1));
525 fFunction2->Copy(*(((TF1Convolution &)obj).fFunction2));
526 // fGraphConv is transient anyway, so we don't bother to copy it
527}
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:377
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char x2
float xmin
float xmax
Int_t gDebug
Definition TROOT.cxx:597
#define gROOT
Definition TROOT.h:407
static IntegrationOneDim::Type DefaultIntegratorType()
User Class for performing numerical integration of a function in one dimension.
Definition Integrator.h:98
double IntegralUp(const IGenFunction &f, double a)
evaluate the Integral of a function f over the semi-infinite interval (a,+inf)
Definition Integrator.h:278
double Integral(Function &f, double a, double b)
evaluate the Integral of a function f over the defined interval (a,b)
Definition Integrator.h:500
double IntegralLow(const IGenFunction &f, double b)
evaluate the Integral of a function f over the over the semi-infinite interval (-inf,...
Definition Integrator.h:296
TF1Convolution_EvalWrapper(TF1 &f1, TF1 &f2, Double_t t)
Double_t operator()(Double_t x) const
Class wrapping convolution of two functions.
std::vector< Double_t > fParams1
Double_t operator()(const Double_t *x, const Double_t *p) override
Used in TF1 when doing the fit, will be evaluated at each point.
void SetParameters(const Double_t *params) override
Set the vector of parameters p for the convolution function g(x,p) = f1 * f2.
Int_t fCstIndex
Index of the constant parameter f the first function.
void GetRange(Double_t &a, Double_t &b) const
Get the range used for the convolution.
std::vector< TString > fParNames
Parameters' names.
std::unique_ptr< TF1 > fFunction1
First function to be convolved.
Int_t fNofPoints
Number of point for FFT array.
static Double_t SetDefaultExtraRange(Double_t percentage)
Set the default extra range fraction used when doing a FFT convolution.
Double_t GetXmin() const
std::vector< Double_t > fParams2
Double_t fXmin
Minimal bound of the range of the convolution.
void SetExtraRange(Double_t percentage)
Set the fraction of extra range used when doing an FFT convolution.
void Copy(TObject &obj) const override
Copy this to obj.
Double_t EvalNumConv(Double_t t)
Perform numerical convolution.
void Update() override
Update the two component functions of the convolution.
void SetRange(Double_t a, Double_t b) override
Set the actual range used for the convolution.
Bool_t fFlagFFT
Choose FFT or numerical convolution.
Double_t GetXmax() const
TF1Convolution & operator=(const TF1Convolution &rhs)
Operator =.
void MakeFFTConv()
Perform the FFT of the two functions.
static Double_t fgExtraRangeFraction
! Additional default fraction of the range used for FFT convolution
void SetNofPointsFFT(Int_t n)
Set the number of points used for the FFT convolution.
TF1Convolution()
constructor without arguments.
const char * GetParName(Int_t ipar) const
Double_t EvalFFTConv(Double_t t)
Perform FFT convolution.
Double_t fXmax
Maximal bound of the range of the convolution.
std::unique_ptr< TF1 > fFunction2
Second function to be convolved.
Bool_t fFlagGraph
! Tells if the graph is already done or not
void InitializeDataMembers(TF1 *function1, TF1 *function2, Bool_t useFFT)
Internal function to initialize data members.
std::unique_ptr< TGraph > fGraphConv
! Graph of the convolution
1-Dim function class
Definition TF1.h:214
static void InitStandardFunctions()
Create the basic function objects.
Definition TF1.cxx:2497
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1007
@ kNotGlobal
Definition TF1.h:326
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1470
virtual Int_t GetNdim() const
Definition TF1.h:491
@ kIsSortedX
Graph is sorted in X points.
Definition TGraph.h:78
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
An array of TObjects.
Definition TObjArray.h:31
Collectable string class.
Definition TObjString.h:28
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:973
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:403
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1015
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:961
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition TString.cxx:2242
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
TVirtualFFT is an interface class for Fast Fourier Transforms.
Definition TVirtualFFT.h:88
static TVirtualFFT * FFT(Int_t ndim, Int_t *n, Option_t *option)
Returns a pointer to the FFT of requested size and type.
virtual Double_t GetPointReal(Int_t ipoint, Bool_t fromInput=kFALSE) const =0
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TF1 * f1
Definition legend1.C:11
Double_t Infinity()
Returns an infinity as defined by the IEEE standard.
Definition TMath.h:917