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 "TF1Convolution.h"
12#include "TROOT.h"
13#include "TObject.h"
14#include "TObjString.h"
15#include "TObjArray.h"
16#include "TMath.h"
17#include "Math/Integrator.h"
23#include "Math/Functor.h"
24#include "TVirtualFFT.h"
25
26/** \class TF1Convolution
27 \ingroup Functions
28 \brief Class wrapping convolution of two functions
29
30Class wrapping convolution of two functions: evaluation of \f$\int f(x)g(x-t)dx\f$
31
32The convolution is performed by default using FFTW if it is available .
33One can pass optionally the range of the convolution (by default the first function range is used).
34Note that when using Discrete Fourier Transform (as FFTW), it is a circular transform, so the functions should be
35approximately zero at the end of the range. If they are significantly different than zero on one side (e.g. the left side)
36a spill over will occur on the other side (e.g right side).
37If no function range is given by default the function1 range + 10% is used
38One should use also a not too small number of points for the DFT (a minimum of 1000). By default 10000 points are used.
39*/
40
42
44
46{
50
51public:
53 fFunc1(&f1),
54 fFunc2(&f2),
55 fT0(t)
56 {}
58 {
59 // use EvalPar that is faster
60 Double_t xx[2];
61 xx[0] = x;
62 xx[1] = fT0-x;
63 return fFunc1->EvalPar(xx,nullptr) * fFunc2->EvalPar(xx+1,nullptr);
64 }
65};
66
67////////////////////////////////////////////////////////////////////////////////
68/// Internal function to initialize data members.
69/// Use TF1::Copy instead of Clone.
70
71void TF1Convolution::InitializeDataMembers(TF1* function1, TF1* function2, Bool_t useFFT)
72{
73 if (function1) {
74 // functions must be 1d- if not flag an error
75 if (function1->GetNdim() != 1)
76 Error("InitializeDataMembers","function1 %s is not of dimension 1 ",function1->GetName());
77 //TF1 * fnew1 = (TF1*) function1->IsA()->New();
78 // since function1 is a TF1 (cannot be a derived class) we can instantiate it directly
79 fFunction1 = std::unique_ptr<TF1> (new TF1());
80 function1->Copy(*fFunction1);
81 }
82 if (function2) {
83 if (function2->GetNdim() != 1)
84 Error("InitializeDataMembers","function2 %s is not of dimension 1 ",function2->GetName());
85 //TF1 * fnew2 = (TF1*) function2->IsA()->New();
86 fFunction2 = std::unique_ptr<TF1>(new TF1());
87 function2->Copy(*fFunction2);
88 }
89 if (fFunction1.get() == nullptr|| fFunction2.get() == nullptr)
90 Fatal("InitializeDataMembers","Invalid functions - Abort");
91
92 // Set kNotGlobal bit
95
96 // use by default range of first function
97 fFunction1->GetRange(fXmin, fXmax);
98 // when using FFT add by default an extra 10% on each side
99 if (useFFT) {
101 }
102 fNofParams1 = fFunction1->GetNpar();
103 fNofParams2 = fFunction2->GetNpar();
104 fParams1 = std::vector<Double_t>(fNofParams1);
105 fParams2 = std::vector<Double_t>(fNofParams2);
106 fCstIndex = (fFunction1->GetParNumber("Constant") == -1)
107 ? -1
108 : fFunction2->GetParNumber("Constant"); // TODO: add dropConstantParam flag?
109 fFlagFFT = useFFT;
110 fFlagGraph = false;
111 fNofPoints = 10000;
112
114 for (int i=0; i<fNofParams1; i++)
115 {
116 fParams1[i] = fFunction1 -> GetParameter(i);
117 fParNames.push_back(fFunction1 -> GetParName(i) );
118 }
119 for (int i=0; i<fNofParams2; i++)
120 {
121 fParams2[i] = fFunction2 -> GetParameter(i);
122 if (i != fCstIndex) fParNames.push_back(fFunction2 -> GetParName(i) );
123 }
124 if (fCstIndex!=-1)
125 {
126 fFunction2 -> FixParameter(fCstIndex,1.);
128 fParams2.erase(fParams2.begin()+fCstIndex);
129 }
130}
131
132////////////////////////////////////////////////////////////////////////////////
133/// constructor without arguments.
134
136{
137 // Nothing to do
138}
139
140////////////////////////////////////////////////////////////////////////////////
141/// constructor from the two function pointer and a flag is using FFT.
142
143TF1Convolution::TF1Convolution(TF1* function1, TF1* function2, Bool_t useFFT)
144{
145 InitializeDataMembers(function1,function2, useFFT);
146}
147
148////////////////////////////////////////////////////////////////////////////////
149/// Constructor from the two function pointer and the convolution range.
150
152{
153 InitializeDataMembers(function1, function2,useFFT);
154 if (xmin < xmax) {
155 fXmin = xmin;
156 fXmax = xmax;
158 } else {
159 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
161 }
162}
163
164////////////////////////////////////////////////////////////////////////////////
165/// Constructor from a formula expression as f1 * f2 where f1 and f2 are two functions known to ROOT.
166
168{
170
171 TObjArray *objarray = formula.Tokenize("*");
172 std::vector < TString > stringarray(2);
173 std::vector < TF1* > funcarray(2);
174 for (int i=0; i<2; i++)
175 {
176 stringarray[i] = ((TObjString*)((*objarray)[i])) -> GetString();
177 stringarray[i].ReplaceAll(" ","");
178 funcarray[i] = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(stringarray[i]));
179 // case function is not found try to use as a TFormula
180 if (funcarray[i] == nullptr) {
181 TF1 * f = new TF1(TString::Format("f_conv_%d",i+1),stringarray[i]);
182 if (!f->GetFormula()->IsValid() )
183 Error("TF1Convolution","Invalid formula : %s",stringarray[i].Data() );
184 if (i == 0)
185 fFunction1 = std::unique_ptr<TF1>(f);
186 else
187 fFunction2 = std::unique_ptr<TF1>(f);
188 }
189 }
190 InitializeDataMembers(funcarray[0], funcarray[1],useFFT);
191 if (xmin < xmax) {
192 fXmin = xmin;
193 fXmax = xmax;
195 } else {
196 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
198 }
199}
200
201////////////////////////////////////////////////////////////////////////////////
202/// Constructor from 2 function names where f1 and f2 are two functions known to
203/// ROOT.
204///
205/// If the function names are not known to ROOT, tries to interpret them as
206/// TFormula.
208{
210 (TString)formula1.ReplaceAll(" ","");
211 (TString)formula2.ReplaceAll(" ","");
212 TF1* f1 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula1));
213 TF1* f2 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula2));
214 // if function do not exists try using TFormula
215 if (!f1) {
216 fFunction1 = std::unique_ptr<TF1>(new TF1("f_conv_1", formula1));
217 if (!fFunction1->GetFormula()->IsValid() )
218 Error("TF1Convolution","Invalid formula for : %s",formula1.Data() );
219 }
220 if (!f2) {
221 fFunction2 = std::unique_ptr<TF1>(new TF1("f_conv_1", formula2));
222 if (!fFunction2->GetFormula()->IsValid() )
223 Error("TF1Convolution","Invalid formula for : %s",formula2.Data() );
224 }
225 // if f1 or f2 are null ptr are not used in InitializeDataMembers
226 InitializeDataMembers(f1, f2,useFFT);
227 if (xmin < xmax) {
228 fXmin = xmin;
229 fXmax = xmax;
230 } else {
231 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
233 }
234}
235
236////////////////////////////////////////////////////////////////////////////////
237/// Copy constructor (necessary to hold unique_ptr as member variable).
238
240{
241 conv.TF1Convolution::Copy(*this);
242}
243
244////////////////////////////////////////////////////////////////////////////////
245/// Operator =
246
248{
249 if (this != &rhs)
250 rhs.TF1Convolution::Copy(*this);
251 return *this;
252}
253
254////////////////////////////////////////////////////////////////////////////////
255/// Perform the FFT of the two functions.
256
258{
259 if (gDebug)
260 Info("MakeFFTConv","Making FFT convolution using %d points in range [%g,%g]",fNofPoints,fXmin,fXmax);
261
262 std::vector < Double_t > x (fNofPoints);
263 std::vector < Double_t > in1(fNofPoints);
264 std::vector < Double_t > in2(fNofPoints);
265
266 TVirtualFFT *fft1 = TVirtualFFT::FFT(1, &fNofPoints, "R2C K");
267 TVirtualFFT *fft2 = TVirtualFFT::FFT(1, &fNofPoints, "R2C K");
268 if (fft1 == nullptr || fft2 == nullptr) {
269 Warning("MakeFFTConv","Cannot use FFT, probably FFTW package is not available. Switch to numerical convolution");
270 fFlagFFT = false;
271 return;
272 }
273
274 // apply a shift in order to have the second function centered around middle of the range of the convolution
275 Double_t shift2 = 0.5*(fXmin+fXmax);
276 Double_t x2;
277 for (int i=0; i<fNofPoints; i++)
278 {
279 x[i] = fXmin + (fXmax-fXmin)/(fNofPoints-1)*i;
280 x2 = x[i] - shift2;
281 in1[i] = fFunction1 -> EvalPar( &x[i], nullptr);
282 in2[i] = fFunction2 -> EvalPar( &x2, nullptr);
283 fft1 -> SetPoint(i, in1[i]);
284 fft2 -> SetPoint(i, in2[i]);
285 }
286 fft1 -> Transform();
287 fft2 -> Transform();
288
289 //inverse transformation of the product
290
291 TVirtualFFT *fftinverse = TVirtualFFT::FFT(1, &fNofPoints, "C2R K");
292 Double_t re1, re2, im1, im2, out_re, out_im;
293
294 for (int i=0;i<=fNofPoints/2.;i++)
295 {
296 fft1 -> GetPointComplex(i,re1,im1);
297 fft2 -> GetPointComplex(i,re2,im2);
298 out_re = re1*re2 - im1*im2;
299 out_im = re1*im2 + re2*im1;
300 fftinverse -> SetPoint(i, out_re, out_im);
301 }
302 fftinverse -> Transform();
303
304 // fill a graph with the result of the convolution
305 if (!fGraphConv)
306 fGraphConv = std::unique_ptr<TGraph>(new TGraph(fNofPoints));
307
308 for (int i=0;i<fNofPoints;i++)
309 {
310 // we need this since we have applied a shift in the middle of f2
311 int j = i + fNofPoints/2;
312 if (j >= fNofPoints) j -= fNofPoints;
313 // need to normalize by dividing by the number of points and multiply by the bin width = Range/Number of points
314 fGraphConv->SetPoint(i, x[i], fftinverse->GetPointReal(j)*(fXmax-fXmin)/(fNofPoints*fNofPoints) );
315 }
316 fGraphConv->SetBit(TGraph::kIsSortedX); // indicate that points are sorted in X to speed up TGraph::Eval
317 fFlagGraph = true; // we can use the graph
318
319 // delete the fft objects
320 delete fft1;
321 delete fft2;
322 delete fftinverse;
323}
324
325////////////////////////////////////////////////////////////////////////////////
326/// Perform FFT convolution.
327
329{
330 if (!fFlagGraph) MakeFFTConv();
331 // if cannot make FFT use numconv
332 if (fGraphConv)
333 return fGraphConv -> Eval(t);
334 else
335
336 return EvalNumConv(t);
337}
338
339////////////////////////////////////////////////////////////////////////////////
340/// Perform numerical convolution.
341///
342
344{
345 /// Could in principle cache the integral in a Graph as it is done for the FFTW
347 Double_t result = 0;
348
350 if (fXmin != - TMath::Infinity() && fXmax != TMath::Infinity() )
351 result = integrator.Integral(fXmin, fXmax);
352 else if (fXmin == - TMath::Infinity() && fXmax != TMath::Infinity() )
353 result = integrator.IntegralLow(fXmax);
354 else if (fXmin != - TMath::Infinity() && fXmax == TMath::Infinity() )
355 result = integrator.IntegralUp(fXmin);
356 else if (fXmin == - TMath::Infinity() && fXmax == TMath::Infinity() )
357 result = integrator.Integral();
358
359 return result;
360}
361
362////////////////////////////////////////////////////////////////////////////////
363/// Used in TF1 when doing the fit, will be evaluated at each point.
364
366{
367 if (p!=0) TF1Convolution::SetParameters(p); // first refresh the parameters
368
369 Double_t result = 0.;
370 if (fFlagFFT)
371 result = EvalFFTConv(x[0]);
372 else
373 result = EvalNumConv(x[0]);
374 return result;
375}
376
377////////////////////////////////////////////////////////////////////////////////
378/// Set the number of points used for the FFT convolution.
379
381{
382 if (n<0) return;
383 fNofPoints = n;
384 if (fGraphConv) fGraphConv -> Set(fNofPoints);
385 fFlagGraph = false; // to indicate we need to re-do the graph
386}
387
388////////////////////////////////////////////////////////////////////////////////
389/// Set the vector of parameters p for the convolution function g(x,p) = f1 * f2.
390
392{
393 bool equalParams = true;
394 for (int i=0; i<fNofParams1; i++) {
395 fFunction1->SetParameter(i, params[i]);
396 equalParams &= (fParams1[i] == params[i]);
397 fParams1[i] = params[i];
398 }
399 Int_t k = 0;
400 Int_t offset = 0;
401 Int_t offset2 = 0;
402 if (fCstIndex!=-1) offset = 1;
403 Int_t totalnofparams = fNofParams1+fNofParams2+offset;
404 for (int i=fNofParams1; i<totalnofparams; i++) {
405 if (k==fCstIndex)
406 {
407 k++;
408 offset2=1;
409 continue;
410 }
411 fFunction2->SetParameter(k, params[i - offset2]);
412 equalParams &= (fParams2[k - offset2] == params[i - offset2]);
413 fParams2[k - offset2] = params[i - offset2];
414 k++;
415 }
416
417 if (!equalParams) fFlagGraph = false; // to indicate we need to re-do the convolution
418}
419
420////////////////////////////////////////////////////////////////////////////////
421/// Set the parameter values for the convolution function.
422
424 Double_t p4, Double_t p5, Double_t p6, Double_t p7)
425{
426 Double_t params[]={p0,p1,p2,p3,p4,p5,p6,p7};
428}
429
430////////////////////////////////////////////////////////////////////////////////
431/// Set the fraction of extra range used when doing an FFT convolution.
432/// The extra range is often needed to avoid mirroring effect of the resulting convolution
433/// function at the borders.
434/// By default an extra range of 0.1 is used.
435
437{
438 if (percentage<0) return;
439 double range = fXmax - fXmin;
440 fXmin -= percentage * range;
441 fXmax += percentage * range;
442 fFlagGraph = false; // to indicate we need to re-do the convolution
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// Set the actual range used for the convolution.
447/// In case a or b are -inf or +inf and FFT convolution is used, then the
448/// range of the first function will be used and extended by the default extra range fraction.
449
451{
452 if (a >= b) {
453 Warning("SetRange", "Invalid range: %f >= %f", a, b);
454 return;
455 }
456
457 fXmin = a;
458 fXmax = b;
459 if (fFlagFFT && ( a==-TMath::Infinity() || b==TMath::Infinity() ) )
460 {
461 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.");
462 if (a ==-TMath::Infinity()) fXmin = fFunction1 -> GetXmin();
463 if ( b== TMath::Infinity()) fXmax = fFunction1 -> GetXmax();
464 // add a spill over of 10% in this case
466 }
467 fFlagGraph = false; // to indicate we need to re-do the convolution
468}
469
470////////////////////////////////////////////////////////////////////////////////
471/// Set the default extra range fraction used when doing a FFT convolution.
472/// By default the value is 0.1 (10%).
473/// The function return the previous default defined value.
474
476{
477 Double_t prevValue = fgExtraRangeFraction;
478 fgExtraRangeFraction = fraction;
479 return prevValue;
480}
481
482////////////////////////////////////////////////////////////////////////////////
483/// Get the range used for the convolution.
484
486{
487 a = fXmin;
488 b = fXmax;
489}
490
491////////////////////////////////////////////////////////////////////////////////
492/// Update the two component functions of the convolution.
493
495{
496 fFunction1->Update();
497 fFunction2->Update();
498}
499
500////////////////////////////////////////////////////////////////////////////////
501
503{
504 // copy numbers
505 ((TF1Convolution &)obj).fXmin = fXmin;
506 ((TF1Convolution &)obj).fXmax = fXmax;
507 ((TF1Convolution &)obj).fNofParams1 = fNofParams1;
508 ((TF1Convolution &)obj).fNofParams2 = fNofParams2;
509 ((TF1Convolution &)obj).fCstIndex = fCstIndex;
510 ((TF1Convolution &)obj).fNofPoints = fNofPoints;
511 ((TF1Convolution &)obj).fFlagFFT = fFlagFFT;
512 ((TF1Convolution &)obj).fFlagGraph = false; // since we're not copying the graph
513
514 // copy vectors
515 ((TF1Convolution &)obj).fParams1 = fParams1;
516 ((TF1Convolution &)obj).fParams2 = fParams2;
517 ((TF1Convolution &)obj).fParNames = fParNames;
518
519 // we need to copy the content of the unique_ptr's
520 ((TF1Convolution &)obj).fFunction1 = std::make_unique<TF1>();
521 ((TF1Convolution &)obj).fFunction2 = std::make_unique<TF1>();
522 fFunction1->Copy(*(((TF1Convolution &)obj).fFunction1));
523 fFunction2->Copy(*(((TF1Convolution &)obj).fFunction2));
524 // fGraphConv is transient anyway, so we don't bother to copy it
525}
#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:585
#define gROOT
Definition TROOT.h:405
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:213
static void InitStandardFunctions()
Create the basic function objects.
Definition TF1.cxx:2496
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1006
@ kNotGlobal
Definition TF1.h:325
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1469
virtual Int_t GetNdim() const
Definition TF1.h:490
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
@ kIsSortedX
Graph is sorted in X points.
Definition TGraph.h:77
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:956
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:970
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:998
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:944
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:915