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
43
45
47{
51
52public:
54 fFunc1(&f1),
55 fFunc2(&f2),
56 fT0(t)
57 {}
59 {
60 // use EvalPar that is faster
61 Double_t xx[2];
62 xx[0] = x;
63 xx[1] = fT0-x;
64 return fFunc1->EvalPar(xx,nullptr) * fFunc2->EvalPar(xx+1,nullptr);
65 }
66};
67
68////////////////////////////////////////////////////////////////////////////////
69/// Internal function to initialize data members.
70/// Use TF1::Copy instead of Clone.
71
73{
74 if (function1) {
75 // functions must be 1d- if not flag an error
76 if (function1->GetNdim() != 1)
77 Error("InitializeDataMembers","function1 %s is not of dimension 1 ",function1->GetName());
78 //TF1 * fnew1 = (TF1*) function1->IsA()->New();
79 // since function1 is a TF1 (cannot be a derived class) we can instantiate it directly
80 fFunction1 = std::make_unique<TF1> ();
81 function1->Copy(*fFunction1);
82 }
83 if (function2) {
84 if (function2->GetNdim() != 1)
85 Error("InitializeDataMembers","function2 %s is not of dimension 1 ",function2->GetName());
86 //TF1 * fnew2 = (TF1*) function2->IsA()->New();
87 fFunction2 = std::make_unique<TF1>();
88 function2->Copy(*fFunction2);
89 }
90 if (fFunction1 == nullptr|| fFunction2 == nullptr)
91 Fatal("InitializeDataMembers","Invalid functions - Abort");
92
93 // Set kNotGlobal bit
96
97 // use by default range of first function
98 fFunction1->GetRange(fXmin, fXmax);
99 // when using FFT add by default an extra 10% on each side
100 if (useFFT) {
102 }
103 fNofParams1 = fFunction1->GetNpar();
104 fNofParams2 = fFunction2->GetNpar();
105 fParams1 = std::vector<Double_t>(fNofParams1);
106 fParams2 = std::vector<Double_t>(fNofParams2);
107 fCstIndex = (fFunction1->GetParNumber("Constant") == -1)
108 ? -1
109 : fFunction2->GetParNumber("Constant"); // TODO: add dropConstantParam flag?
111 fFlagGraph = false;
112 fNofPoints = 10000;
113
115 for (int i=0; i<fNofParams1; i++)
116 {
117 fParams1[i] = fFunction1 -> GetParameter(i);
118 fParNames.push_back(fFunction1 -> GetParName(i) );
119 }
120 for (int i=0; i<fNofParams2; i++)
121 {
122 fParams2[i] = fFunction2 -> GetParameter(i);
123 if (i != fCstIndex) fParNames.push_back(fFunction2 -> GetParName(i) );
124 }
125 if (fCstIndex!=-1)
126 {
127 fFunction2 -> FixParameter(fCstIndex,1.);
129 fParams2.erase(fParams2.begin()+fCstIndex);
130 }
131}
132
133////////////////////////////////////////////////////////////////////////////////
134/// constructor without arguments.
135
137{
138 // Nothing to do
139}
140
141////////////////////////////////////////////////////////////////////////////////
142/// constructor from the two function pointer and a flag is using FFT.
143
148
149////////////////////////////////////////////////////////////////////////////////
150/// Constructor from the two function pointer and the convolution range.
151
153{
155 if (xmin < xmax) {
156 fXmin = xmin;
157 fXmax = xmax;
159 } else {
160 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
162 }
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// Constructor from a formula expression as f1 * f2 where f1 and f2 are two functions known to ROOT.
167
169{
171
172 TObjArray *objarray = formula.Tokenize("*");
173 std::vector < TString > stringarray(2);
174 std::vector < TF1* > funcarray(2);
175 for (int i=0; i<2; i++)
176 {
177 stringarray[i] = ((TObjString*)((*objarray)[i])) -> GetString();
178 stringarray[i].ReplaceAll(" ","");
179 funcarray[i] = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(stringarray[i]));
180 // case function is not found try to use as a TFormula
181 if (funcarray[i] == nullptr) {
182 TF1 * f = new TF1(TString::Format("f_conv_%d",i+1),stringarray[i]);
183 if (!f->GetFormula()->IsValid() )
184 Error("TF1Convolution","Invalid formula : %s",stringarray[i].Data() );
185 if (i == 0)
186 fFunction1 = std::unique_ptr<TF1>(f);
187 else
188 fFunction2 = std::unique_ptr<TF1>(f);
189 }
190 }
192 if (xmin < xmax) {
193 fXmin = xmin;
194 fXmax = xmax;
196 } else {
197 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
199 }
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Constructor from 2 function names where f1 and f2 are two functions known to
204/// ROOT.
205///
206/// If the function names are not known to ROOT, tries to interpret them as
207/// TFormula.
209{
211 (TString)formula1.ReplaceAll(" ","");
212 (TString)formula2.ReplaceAll(" ","");
213 TF1* f1 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula1));
214 TF1* f2 = (TF1*)(gROOT -> GetListOfFunctions() -> FindObject(formula2));
215 // if function do not exists try using TFormula
216 if (!f1) {
217 fFunction1 = std::make_unique<TF1>("f_conv_1", formula1);
218 if (!fFunction1->GetFormula()->IsValid() )
219 Error("TF1Convolution","Invalid formula for : %s",formula1.Data() );
220 }
221 if (!f2) {
222 fFunction2 = std::make_unique<TF1>("f_conv_1", formula2);
223 if (!fFunction2->GetFormula()->IsValid() )
224 Error("TF1Convolution","Invalid formula for : %s",formula2.Data() );
225 }
226 // if f1 or f2 are null ptr are not used in InitializeDataMembers
228 if (xmin < xmax) {
229 fXmin = xmin;
230 fXmax = xmax;
231 } else {
232 Info("TF1Convolution", "Using default range [-inf, inf] for TF1Convolution");
234 }
235}
236
237////////////////////////////////////////////////////////////////////////////////
238/// Copy constructor (necessary to hold unique_ptr as member variable).
239
241{
242 conv.TF1Convolution::Copy(*this);
243}
244
245////////////////////////////////////////////////////////////////////////////////
246/// Operator =
247
249{
250 if (this != &rhs)
251 rhs.TF1Convolution::Copy(*this);
252 return *this;
253}
254
255////////////////////////////////////////////////////////////////////////////////
256/// Perform the FFT of the two functions.
257
259{
260 if (gDebug)
261 Info("MakeFFTConv","Making FFT convolution using %d points in range [%g,%g]",fNofPoints,fXmin,fXmax);
262
263 std::vector < Double_t > x (fNofPoints);
264 std::vector < Double_t > in1(fNofPoints);
265 std::vector < Double_t > in2(fNofPoints);
266
269 if (fft1 == nullptr || fft2 == nullptr) {
270 Warning("MakeFFTConv","Cannot use FFT, probably FFTW package is not available. Switch to numerical convolution");
271 fFlagFFT = false;
272 return;
273 }
274
275 // apply a shift in order to have the second function centered around middle of the range of the convolution
276 Double_t shift2 = 0.5*(fXmin+fXmax);
277 Double_t x2;
278 for (int i=0; i<fNofPoints; i++)
279 {
280 x[i] = fXmin + (fXmax-fXmin)/(fNofPoints-1)*i;
281 x2 = x[i] - shift2;
282 in1[i] = fFunction1 -> EvalPar( &x[i], nullptr);
283 in2[i] = fFunction2 -> EvalPar( &x2, nullptr);
284 fft1 -> SetPoint(i, in1[i]);
285 fft2 -> SetPoint(i, in2[i]);
286 }
287 fft1 -> Transform();
288 fft2 -> Transform();
289
290 //inverse transformation of the product
291
294
295 for (int i=0;i<=fNofPoints/2.;i++)
296 {
297 fft1 -> GetPointComplex(i,re1,im1);
298 fft2 -> GetPointComplex(i,re2,im2);
299 out_re = re1*re2 - im1*im2;
300 out_im = re1*im2 + re2*im1;
301 fftinverse -> SetPoint(i, out_re, out_im);
302 }
303 fftinverse -> Transform();
304
305 // fill a graph with the result of the convolution
306 if (!fGraphConv)
307 fGraphConv = std::make_unique<TGraph>(fNofPoints);
308
309 for (int i=0;i<fNofPoints;i++)
310 {
311 // we need this since we have applied a shift in the middle of f2
312 int j = i + fNofPoints/2;
313 if (j >= fNofPoints) j -= fNofPoints;
314 // need to normalize by dividing by the number of points and multiply by the bin width = Range/Number of points
315 fGraphConv->SetPoint(i, x[i], fftinverse->GetPointReal(j)*(fXmax-fXmin)/(fNofPoints*fNofPoints) );
316 }
317 fGraphConv->SetBit(TGraph::kIsSortedX); // indicate that points are sorted in X to speed up TGraph::Eval
318 fFlagGraph = true; // we can use the graph
319
320 // delete the fft objects
321 delete fft1;
322 delete fft2;
323 delete fftinverse;
324}
325
326////////////////////////////////////////////////////////////////////////////////
327/// Perform FFT convolution.
328
330{
331 if (!fFlagGraph) {
332 MakeFFTConv();
333 }
334 // if cannot make FFT use numconv
335 if (fGraphConv)
336 return fGraphConv -> Eval(t);
337 else
338
339 return EvalNumConv(t);
340}
341
342////////////////////////////////////////////////////////////////////////////////
343/// Perform numerical convolution.
344///
345
347{
348 /// Could in principle cache the integral in a Graph as it is done for the FFTW
350 Double_t result = 0;
351
353 if (fXmin != - TMath::Infinity() && fXmax != TMath::Infinity() )
354 result = integrator.Integral(fXmin, fXmax);
355 else if (fXmin == - TMath::Infinity() && fXmax != TMath::Infinity() )
356 result = integrator.IntegralLow(fXmax);
357 else if (fXmin != - TMath::Infinity() && fXmax == TMath::Infinity() )
358 result = integrator.IntegralUp(fXmin);
359 else if (fXmin == - TMath::Infinity() && fXmax == TMath::Infinity() )
360 result = integrator.Integral();
361
362 return result;
363}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Used in TF1 when doing the fit, will be evaluated at each point.
367
369{
370 if (p!=nullptr) TF1Convolution::SetParameters(p); // first refresh the parameters
371
372 Double_t result = 0.;
373 if (fFlagFFT)
374 result = EvalFFTConv(x[0]);
375 else
376 result = EvalNumConv(x[0]);
377 return result;
378}
379
380////////////////////////////////////////////////////////////////////////////////
381/// Set the number of points used for the FFT convolution.
382
384{
385 if (n<0) return;
386 fNofPoints = n;
387 if (fGraphConv) fGraphConv -> Set(fNofPoints);
388 fFlagGraph = false; // to indicate we need to re-do the graph
389}
390
391////////////////////////////////////////////////////////////////////////////////
392/// Set the vector of parameters p for the convolution function g(x,p) = f1 * f2.
393
395{
396 bool equalParams = true;
397 for (int i=0; i<fNofParams1; i++) {
398 fFunction1->SetParameter(i, params[i]);
399 equalParams &= (fParams1[i] == params[i]);
400 fParams1[i] = params[i];
401 }
402 Int_t k = 0;
403 Int_t offset = 0;
404 Int_t offset2 = 0;
405 if (fCstIndex!=-1) offset = 1;
407 for (int i=fNofParams1; i<totalnofparams; i++) {
408 if (k==fCstIndex)
409 {
410 k++;
411 offset2=1;
412 continue;
413 }
414 fFunction2->SetParameter(k, params[i - offset2]);
415 equalParams &= (fParams2[k - offset2] == params[i - offset2]);
416 fParams2[k - offset2] = params[i - offset2];
417 k++;
418 }
419
420 if (!equalParams) fFlagGraph = false; // to indicate we need to re-do the convolution
421}
422
423////////////////////////////////////////////////////////////////////////////////
424/// Set the parameter values for the convolution function.
425
432
433////////////////////////////////////////////////////////////////////////////////
434/// Set the fraction of extra range used when doing an convolution.
435/// The extra range is often needed to avoid mirroring effect of the resulting convolution
436/// function at the borders.
437/// By default an extra range of 0.1 is used when doing FFT and it is not use for numerical convolution
438
440{
441 if (percentage<0) return;
442 double range = fXmax - fXmin;
445 fFlagGraph = false; // to indicate we need to re-do the convolution
446}
447
448////////////////////////////////////////////////////////////////////////////////
449/// Set the actual range used for the convolution.
450/// In case a or b are -inf or +inf and FFT convolution is used, then the
451/// range of the first function will be used and extended by the default extra range fraction.
452
454{
455 if (a >= b) {
456 Warning("SetRange", "Invalid range: %f >= %f", a, b);
457 return;
458 }
459
460 fXmin = a;
461 fXmax = b;
462 if (fFlagFFT && ( a==-TMath::Infinity() || b==TMath::Infinity() ) )
463 {
464 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.");
465 if (a ==-TMath::Infinity()) fXmin = fFunction1 -> GetXmin();
466 if ( b== TMath::Infinity()) fXmax = fFunction1 -> GetXmax();
467 // add a spill over of 10% in this case
469 }
470 fFlagGraph = false; // to indicate we need to re-do the convolution
471}
472
473////////////////////////////////////////////////////////////////////////////////
474/// Set the default extra range fraction used when doing a FFT convolution.
475/// By default the value is 0.1 (10%).
476/// The function return the previous default defined value.
477
484
485////////////////////////////////////////////////////////////////////////////////
486/// Get the range used for the convolution.
487
489{
490 a = fXmin;
491 b = fXmax;
492}
493
494////////////////////////////////////////////////////////////////////////////////
495/// Update the two component functions of the convolution.
496
498{
499 fFunction1->Update();
500 fFunction2->Update();
501}
502
503////////////////////////////////////////////////////////////////////////////////
504
506{
507 // copy numbers
508 ((TF1Convolution &)obj).fXmin = fXmin;
509 ((TF1Convolution &)obj).fXmax = fXmax;
510 ((TF1Convolution &)obj).fNofParams1 = fNofParams1;
511 ((TF1Convolution &)obj).fNofParams2 = fNofParams2;
512 ((TF1Convolution &)obj).fCstIndex = fCstIndex;
513 ((TF1Convolution &)obj).fNofPoints = fNofPoints;
514 ((TF1Convolution &)obj).fFlagFFT = fFlagFFT;
515 ((TF1Convolution &)obj).fFlagGraph = false; // since we're not copying the graph
516
517 // copy vectors
518 ((TF1Convolution &)obj).fParams1 = fParams1;
519 ((TF1Convolution &)obj).fParams2 = fParams2;
520 ((TF1Convolution &)obj).fParNames = fParNames;
521
522 // we need to copy the content of the unique_ptr's
523 ((TF1Convolution &)obj).fFunction1 = std::make_unique<TF1>();
524 ((TF1Convolution &)obj).fFunction2 = std::make_unique<TF1>();
525 fFunction1->Copy(*(((TF1Convolution &)obj).fFunction1));
526 fFunction2->Copy(*(((TF1Convolution &)obj).fFunction2));
527 // fGraphConv is transient anyway, so we don't bother to copy it
528}
#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
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
Global variable setting the debug level. Set to 0 to disable, increase it in steps of 1 to increase t...
Definition TROOT.cxx:627
#define gROOT
Definition TROOT.h:411
static IntegrationOneDim::Type DefaultIntegratorType()
User Class for performing numerical integration of a function in one dimension.
Definition Integrator.h:94
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
double Integral(Function &f, double a, double b)
evaluate the Integral of a function f over the defined interval (a,b)
Definition Integrator.h:495
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:292
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 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:182
static void InitStandardFunctions()
Create the basic function objects.
Definition TF1.cxx:2522
@ kNotGlobal
Definition TF1.h:296
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1475
@ kIsSortedX
Graph is sorted in X points.
Definition TGraph.h:78
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:1057
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:421
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1099
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1045
Basic string class.
Definition TString.h:138
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition TString.cxx:2270
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:2384
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.
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:928