Logo ROOT   6.12/06
Reference Guide
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 "Riostream.h"
13 #include "TROOT.h"
14 #include "TObject.h"
15 #include "TObjString.h"
16 #include "TMath.h"
17 #include "Math/Integrator.h"
19 #include "Math/IntegratorOptions.h"
20 #include "Math/GaussIntegrator.h"
23 #include "Math/Functor.h"
24 #include "TVirtualFFT.h"
25 #include "TClass.h"
26 
27 /** \class TF1Convolution
28  \ingroup Hist
29  \brief Class wrapping convolution of two functions
30 
31 Class wrapping convolution of two functions: evaluation of \f$\int f(x)g(x-t)dx\f$
32 
33 The convolution is performed by default using FFTW if it is available .
34 One can pass optionally the range of the convolution (by default the first function range is used).
35 Note that when using Discrete Fourier Transform (as FFTW), it is a circular transform, so the functions should be
36 approximately zero at the end of the range. If they are significantly different than zero on one side (e.g. the left side)
37 a spill over will occur on the other side (e.g right side).
38 If no function range is given by default the function1 range + 10% is used
39 One should use also a not too small number of points for the DFT (a minimum of 1000). By default 10000 points are used.
40 */
41 
43 
44 class TF1Convolution_EvalWrapper
45 {
46  std::unique_ptr<TF1> fFunction1;
47  std::unique_ptr<TF1> fFunction2;
48  Double_t fT0;
49 
50 public:
51  TF1Convolution_EvalWrapper(std::unique_ptr<TF1> &f1, std::unique_ptr<TF1> &f2, Double_t t) : fT0(t)
52  {
53  fFunction1 = std::unique_ptr<TF1>((TF1 *)f1->Clone());
54  fFunction2 = std::unique_ptr<TF1>((TF1 *)f2->Clone());
55  }
57  {
58  // use EvalPar that is faster
59  Double_t xx[2];
60  xx[0] = x;
61  xx[1] = fT0-x;
62  return fFunction1->EvalPar(xx,nullptr) * fFunction2->EvalPar(xx+1,nullptr);
63  }
64 };
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 /// Use copy instead of Clone
68 
69 void TF1Convolution::InitializeDataMembers(TF1* function1, TF1* function2, Bool_t useFFT)
70 {
71  if (function1) {
72  TF1 * fnew1 = (TF1*) function1->IsA()->New();
73  function1->Copy(*fnew1);
74  fFunction1 = std::unique_ptr<TF1>(fnew1);
75  }
76  if (function2) {
77  TF1 * fnew2 = (TF1*) function2->IsA()->New();
78  function2->Copy(*fnew2);
79  fFunction2 = std::unique_ptr<TF1>(fnew2);
80  }
81  if (fFunction1.get() == nullptr|| fFunction2.get() == nullptr)
82  Fatal("InitializeDataMembers","Invalid functions - Abort");
83 
84  // Set kNotGlobal bit
87 
88  // add by default an extra 10% on each side
89  fFunction1->GetRange(fXmin, fXmax);
90  Double_t range = fXmax - fXmin;
91  fXmin -= 0.1*range;
92  fXmax += 0.1*range;
93  fNofParams1 = fFunction1->GetNpar();
94  fNofParams2 = fFunction2->GetNpar();
95  fParams1 = std::vector<Double_t>(fNofParams1);
96  fParams2 = std::vector<Double_t>(fNofParams2);
97  fCstIndex = (fFunction1->GetParNumber("Constant") == -1)
98  ? -1
99  : fFunction2->GetParNumber("Constant"); // TODO: add dropConstantParam flag?
100  fFlagFFT = useFFT;
101  fFlagGraph = false;
102  fNofPoints = 10000;
103 
104  fParNames.reserve( fNofParams1 + fNofParams2);
105  for (int i=0; i<fNofParams1; i++)
106  {
107  fParams1[i] = fFunction1 -> GetParameter(i);
108  fParNames.push_back(fFunction1 -> GetParName(i) );
109  }
110  for (int i=0; i<fNofParams2; i++)
111  {
112  fParams2[i] = fFunction2 -> GetParameter(i);
113  if (i != fCstIndex) fParNames.push_back(fFunction2 -> GetParName(i) );
114  }
115  if (fCstIndex!=-1)
116  {
117  fFunction2 -> FixParameter(fCstIndex,1.);
119  fParams2.erase(fParams2.begin()+fCstIndex);
120  }
121 }
122 
123 ////////////////////////////////////////////////////////////////////////////////
124 /// constructor without arguments
125 
127 {
128  // Nothing to do
129 }
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// constructor from the two function pointer and a flag is using FFT
133 
134 TF1Convolution::TF1Convolution(TF1* function1, TF1* function2, Bool_t useFFT)
135 {
136  InitializeDataMembers(function1,function2, useFFT);
137 }
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// Constructor from the two function pointer and the convolution range
141 
143 {
144  InitializeDataMembers(function1, function2,useFFT);
145  if (xmin < xmax) {
146  fXmin = xmin;
147  fXmax = xmax;
148  } else {
149  Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
151  }
152 }
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 /// Constructor from a formula expression as f1 * f2 where f1 and f2 are two functions known to ROOT
156 
158 {
160 
161  TObjArray *objarray = formula.Tokenize("*");
162  std::vector < TString > stringarray(2);
163  std::vector < TF1* > funcarray(2);
164  for (int i=0; i<2; i++)
165  {
166  stringarray[i] = ((TObjString*)((*objarray)[i])) -> GetString();
167  stringarray[i].ReplaceAll(" ","");
168  funcarray[i] = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(stringarray[i]));
169  // case function is not found try to use as a TFormula
170  if (funcarray[i] == nullptr) {
171  TF1 * f = new TF1(TString::Format("f_conv_%d",i+1),stringarray[i]);
172  if (!f->GetFormula()->IsValid() )
173  Error("TF1Convolution","Invalid formula : %s",stringarray[i].Data() );
174  if (i == 0)
175  fFunction1 = std::unique_ptr<TF1>(f);
176  else
177  fFunction2 = std::unique_ptr<TF1>(f);
178  }
179  }
180  InitializeDataMembers(funcarray[0], funcarray[1],useFFT);
181  if (xmin < xmax) {
182  fXmin = xmin;
183  fXmax = xmax;
184  } else {
185  Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
187  }
188 }
189 
190 ////////////////////////////////////////////////////////////////////////////////
191 /// constructor from 2 function names where f1 and f2 are two functions known to
192 /// ROOT
193 ///
194 /// if the function names are not known to ROOT, tries to interpret them as
195 /// TFormula
197 {
199  (TString)formula1.ReplaceAll(" ","");
200  (TString)formula2.ReplaceAll(" ","");
201  TF1* f1 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula1));
202  TF1* f2 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula2));
203  // if function do not exists try using TFormula
204  if (!f1) {
205  fFunction1 = std::unique_ptr<TF1>(new TF1("f_conv_1", formula1));
206  if (!fFunction1->GetFormula()->IsValid() )
207  Error("TF1Convolution","Invalid formula for : %s",formula1.Data() );
208  }
209  if (!f2) {
210  fFunction2 = std::unique_ptr<TF1>(new TF1("f_conv_1", formula2));
211  if (!fFunction2->GetFormula()->IsValid() )
212  Error("TF1Convolution","Invalid formula for : %s",formula2.Data() );
213  }
214  // if f1 or f2 are null ptr are not used in InitializeDataMembers
215  InitializeDataMembers(f1, f2,useFFT);
216  if (xmin < xmax) {
217  fXmin = xmin;
218  fXmax = xmax;
219  } else {
220  Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
222  }
223 }
224 
225 ////////////////////////////////////////////////////////////////////////////////
226 /// Copy constructor (necessary to hold unique_ptr as member variable)
227 
229 {
230  conv.Copy((TObject &)*this);
231 }
232 
233 ////////////////////////////////////////////////////////////////////////////////
234 /// Operator =
235 
237 {
238  if (this != &rhs)
239  rhs.Copy(*this);
240  return *this;
241 }
242 
243 ////////////////////////////////////////////////////////////////////////////////
244 /// Perform the FFT of the two functions
245 
247 {
248  if (gDebug)
249  Info("MakeFFTConv","Making FFT convolution using %d points in range [%g,%g]",fNofPoints,fXmin,fXmax);
250 
251  std::vector < Double_t > x (fNofPoints);
252  std::vector < Double_t > in1(fNofPoints);
253  std::vector < Double_t > in2(fNofPoints);
254 
255  TVirtualFFT *fft1 = TVirtualFFT::FFT(1, &fNofPoints, "R2C K");
256  TVirtualFFT *fft2 = TVirtualFFT::FFT(1, &fNofPoints, "R2C K");
257  if (fft1 == nullptr || fft2 == nullptr) {
258  Warning("MakeFFTConv","Cannot use FFT, probably FFTW package is not available. Switch to numerical convolution");
259  fFlagFFT = false;
260  return;
261  }
262 
263  // apply a shift in order to have the second function centered around middle of the range of the convolution
264  Double_t shift2 = 0.5*(fXmin+fXmax);
265  Double_t x2;
266  for (int i=0; i<fNofPoints; i++)
267  {
268  x[i] = fXmin + (fXmax-fXmin)/(fNofPoints-1)*i;
269  x2 = x[i] - shift2;
270  in1[i] = fFunction1 -> EvalPar( &x[i], nullptr);
271  in2[i] = fFunction2 -> EvalPar( &x2, nullptr);
272  fft1 -> SetPoint(i, in1[i]);
273  fft2 -> SetPoint(i, in2[i]);
274  }
275  fft1 -> Transform();
276  fft2 -> Transform();
277 
278  //inverse transformation of the product
279 
280  TVirtualFFT *fftinverse = TVirtualFFT::FFT(1, &fNofPoints, "C2R K");
281  Double_t re1, re2, im1, im2, out_re, out_im;
282 
283  for (int i=0;i<=fNofPoints/2.;i++)
284  {
285  fft1 -> GetPointComplex(i,re1,im1);
286  fft2 -> GetPointComplex(i,re2,im2);
287  out_re = re1*re2 - im1*im2;
288  out_im = re1*im2 + re2*im1;
289  fftinverse -> SetPoint(i, out_re, out_im);
290  }
291  fftinverse -> Transform();
292 
293  // fill a graph with the result of the convolution
294  if (!fGraphConv)
295  fGraphConv = std::unique_ptr<TGraph>(new TGraph(fNofPoints));
296 
297  for (int i=0;i<fNofPoints;i++)
298  {
299  // we need this since we have applied a shift in the middle of f2
300  int j = i + fNofPoints/2;
301  if (j >= fNofPoints) j -= fNofPoints;
302  // need to normalize by dividing by the number of points and multiply by the bin width = Range/Number of points
303  fGraphConv->SetPoint(i, x[i], fftinverse->GetPointReal(j)*(fXmax-fXmin)/(fNofPoints*fNofPoints) );
304  }
305  fGraphConv->SetBit(TGraph::kIsSortedX); // indicate that points are sorted in X to speed up TGraph::Eval
306  fFlagGraph = true; // we can use the graph
307 }
308 
309 ////////////////////////////////////////////////////////////////////////////////
310 
312 {
313  if (!fFlagGraph) MakeFFTConv();
314  // if cannot make FFT use numconv
315  if (fGraphConv)
316  return fGraphConv -> Eval(t);
317  else
318  return EvalNumConv(t);
319 }
320 
321 ////////////////////////////////////////////////////////////////////////////////
322 /// Perform numerical convolution.
323 /// Could in principle cache the integral in a Graph as it is done for the FFTW
324 
326 {
327  TF1Convolution_EvalWrapper fconv( fFunction1, fFunction2, t);
328  Double_t result = 0;
329 
331  if (fXmin != - TMath::Infinity() && fXmax != TMath::Infinity() )
332  result = integrator.Integral(fXmin, fXmax);
333  else if (fXmin == - TMath::Infinity() && fXmax != TMath::Infinity() )
334  result = integrator.IntegralLow(fXmax);
335  else if (fXmin != - TMath::Infinity() && fXmax == TMath::Infinity() )
336  result = integrator.IntegralUp(fXmin);
337  else if (fXmin == - TMath::Infinity() && fXmax == TMath::Infinity() )
338  result = integrator.Integral();
339 
340  return result;
341 }
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 /// Used in TF1 when doing the fit, will be evaluated at each point.
345 
347 {
348  if (p!=0) TF1Convolution::SetParameters(p); // first refresh the parameters
349 
350  Double_t result = 0.;
351  if (fFlagFFT)
352  result = EvalFFTConv(x[0]);
353  else
354  result = EvalNumConv(x[0]);
355  return result;
356 }
357 
358 ////////////////////////////////////////////////////////////////////////////////
359 
361 {
362  if (n<0) return;
363  fNofPoints = n;
364  if (fGraphConv) fGraphConv -> Set(fNofPoints);
365  fFlagGraph = false; // to indicate we need to re-do the graph
366 }
367 
368 ////////////////////////////////////////////////////////////////////////////////
369 
371 {
372  bool equalParams = true;
373  for (int i=0; i<fNofParams1; i++) {
374  fFunction1->SetParameter(i, params[i]);
375  equalParams &= (fParams1[i] == params[i]);
376  fParams1[i] = params[i];
377  }
378  Int_t k = 0;
379  Int_t offset = 0;
380  Int_t offset2 = 0;
381  if (fCstIndex!=-1) offset = 1;
382  Int_t totalnofparams = fNofParams1+fNofParams2+offset;
383  for (int i=fNofParams1; i<totalnofparams; i++) {
384  if (k==fCstIndex)
385  {
386  k++;
387  offset2=1;
388  continue;
389  }
390  fFunction2->SetParameter(k, params[i - offset2]);
391  equalParams &= (fParams2[k - offset2] == params[i - offset2]);
392  fParams2[k - offset2] = params[i - offset2];
393  k++;
394  }
395 
396  if (!equalParams) fFlagGraph = false; // to indicate we need to re-do the convolution
397 }
398 
399 ////////////////////////////////////////////////////////////////////////////////
400 
402  Double_t p4, Double_t p5, Double_t p6, Double_t p7)
403 {
404  Double_t params[]={p0,p1,p2,p3,p4,p5,p6,p7};
406 }
407 
408 ////////////////////////////////////////////////////////////////////////////////
409 
411 {
412  if (percentage<0) return;
413  double range = fXmax = fXmin;
414  fXmin -= percentage * range;
415  fXmax += percentage * range;
416  fFlagGraph = false; // to indicate we need to re-do the convolution
417 }
418 
419 ////////////////////////////////////////////////////////////////////////////////
420 
422 {
423  if (a >= b) {
424  Warning("SetRange", "Invalid range: %f >= %f", a, b);
425  return;
426  }
427 
428  fXmin = a;
429  fXmax = b;
430  if (fFlagFFT && ( a==-TMath::Infinity() || b==TMath::Infinity() ) )
431  {
432  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.");
433  if (a ==-TMath::Infinity()) fXmin = fFunction1 -> GetXmin();
434  if ( b== TMath::Infinity()) fXmax = fFunction1 -> GetXmax();
435  // add a spill over of 10% in this case
436  SetExtraRange(0.1);
437  }
438  fFlagGraph = false; // to indicate we need to re-do the convolution
439 }
440 
441 ////////////////////////////////////////////////////////////////////////////////
442 
444 {
445  a = fXmin;
446  b = fXmax;
447 }
448 
449 ////////////////////////////////////////////////////////////////////////////////
450 /// Update the two component functions of the convolution
451 
453 {
454  fFunction1->Update();
455  fFunction2->Update();
456 }
457 
458 ////////////////////////////////////////////////////////////////////////////////
459 
461 {
462  // copy numbers
463  ((TF1Convolution &)obj).fXmin = fXmin;
464  ((TF1Convolution &)obj).fXmax = fXmax;
465  ((TF1Convolution &)obj).fNofParams1 = fNofParams1;
466  ((TF1Convolution &)obj).fNofParams2 = fNofParams2;
467  ((TF1Convolution &)obj).fCstIndex = fCstIndex;
468  ((TF1Convolution &)obj).fNofPoints = fNofPoints;
469  ((TF1Convolution &)obj).fFlagFFT = fFlagFFT;
470  ((TF1Convolution &)obj).fFlagGraph = false; // since we're not copying the graph
471 
472  // copy vectors
473  ((TF1Convolution &)obj).fParams1 = fParams1;
474  ((TF1Convolution &)obj).fParams2 = fParams2;
475  ((TF1Convolution &)obj).fParNames = fParNames;
476 
477  // Clone unique_ptr's
478  ((TF1Convolution &)obj).fFunction1 = std::unique_ptr<TF1>((TF1 *)fFunction1->Clone());
479  ((TF1Convolution &)obj).fFunction2 = std::unique_ptr<TF1>((TF1 *)fFunction2->Clone());
480  // fGraphConv is transient anyway, so we don't bother to copy it
481 }
An array of TObjects.
Definition: TObjArray.h:37
std::unique_ptr< TF1 > fFunction2
Second function to be convolved.
float xmin
Definition: THbookFile.cxx:93
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:854
virtual TFormula * GetFormula()
Definition: TF1.h:437
static double p3(double t, double a, double b, double c, double d)
virtual void Copy(TObject &f1) const
Copy this F1 to a new F1.
Definition: TF1.cxx:904
const char * GetParName(Int_t ipar) const
void Update()
Update the two component functions of the convolution.
Collectable string class.
Definition: TObjString.h:28
std::unique_ptr< TF1 > fFunction1
First function to be convolved.
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:638
TF1Convolution()
constructor without arguments
#define gROOT
Definition: TROOT.h:393
Basic string class.
Definition: TString.h:125
std::vector< Double_t > fParams2
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t operator()(const Double_t *x, const Double_t *p)
Used in TF1 when doing the fit, will be evaluated at each point.
TF1Convolution & operator=(const TF1Convolution &rhs)
Operator =.
TRObject operator()(const T1 &t1) const
Int_t fNofPoints
Number of point for FFT array.
Class wrapping convolution of two functions.
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
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:2365
static double p2(double t, double a, double b, double c)
void MakeFFTConv()
Perform the FFT of the two functions.
void InitializeDataMembers(TF1 *function1, TF1 *function2, Bool_t useFFT)
Use copy instead of Clone.
static IntegrationOneDim::Type DefaultIntegratorType()
Double_t Infinity()
Definition: TMath.h:796
double IntegralLow(const IGenFunction &f, double b)
evaluate the Integral of a function f over the over the semi-infinite interval (-inf,b)
Definition: Integrator.h:292
Double_t fXmax
Maximal bound of the range of the convolution.
double Integral(Function &f, double a, double b)
evaluate the Integral of a function f over the defined interval (a,b)
Definition: Integrator.h:496
void SetRange(Double_t a, Double_t b)
double IntegralUp(const IGenFunction &f, double a)
evaluate the Integral of a function f over the semi-infinite interval (a,+inf)
Definition: Integrator.h:274
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition: TObject.cxx:321
virtual Double_t GetPointReal(Int_t ipoint, Bool_t fromInput=kFALSE) const =0
auto * a
Definition: textangle.C:12
static TVirtualFFT * FFT(Int_t ndim, Int_t *n, Option_t *option)
Returns a pointer to the FFT of requested size and type.
Double_t EvalNumConv(Double_t t)
Perform numerical convolution.
void SetNofPointsFFT(Int_t n)
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
User Class for performing numerical integration of a function in one dimension.
Definition: Integrator.h:94
static void InitStandardFunctions()
Create the basic function objects.
Definition: TF1.cxx:2369
Bool_t fFlagGraph
! Tells if the graph is already done or not
std::unique_ptr< TGraph > fGraphConv
! Graph of the convolution
static double p1(double t, double a, double b)
float xmax
Definition: THbookFile.cxx:93
TVirtualFFT is an interface class for Fast Fourier Transforms.
Definition: TVirtualFFT.h:88
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition: TString.cxx:2251
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
graph is sorted in X points
Definition: TGraph.h:73
Double_t fXmin
Minimal bound of the range of the convolution.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
Double_t EvalFFTConv(Double_t t)
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:74
Mother of all ROOT objects.
Definition: TObject.h:37
void Copy(TObject &obj) const
Copy this to obj.
1-Dim function class
Definition: TF1.h:211
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
TF1 * f1
Definition: legend1.C:11
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
Int_t fCstIndex
Index of the constant parameter f the first function.
R__EXTERN Int_t gDebug
Definition: Rtypes.h:86
Double_t GetXmin() const
std::vector< TString > fParNames
Parameters&#39; names.
Bool_t fFlagFFT
Choose FFT or numerical convolution.
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:908
Double_t GetXmax() const
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=0)
Evaluate function with given coordinates and parameters.
Definition: TF1.cxx:1363
Bool_t IsValid() const
Definition: TFormula.h:201
std::vector< Double_t > fParams1
const Bool_t kTRUE
Definition: RtypesCore.h:87
const Int_t n
Definition: legend1.C:16
void GetRange(Double_t &a, Double_t &b) const
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
void SetExtraRange(Double_t percentage)
const char * Data() const
Definition: TString.h:345
void SetParameters(const Double_t *params)