GCpp general purpose C++ library  version 1.0
GRealSampleFFT Class Reference

#include <GRealSampleFFT.hh>

Inheritance diagram for GRealSampleFFT:
GBaseSampleFFT

Public Member Functions

Constructors, destructor, affectation
 GRealSampleFFT (u_int n=0, double dt=1.L)
 
 GRealSampleFFT (const GRealSampleFFT &original)
 
 GRealSampleFFT (u_int n, const double data[])
 
 GRealSampleFFT (u_int n, double t_lo, double t_hi, const double data[], const bool ctr=false)
 
 GRealSampleFFT (u_int n, double t_lo, double t_hi, double(*fct_ptr)(const double, const double[]), const double fct_par[], const bool ctr=false)
 
 GRealSampleFFT (u_int n, double t_lo, double t_hi, double(*fct_ptr)(const double[], const double[]), const double fct_par[], const bool ctr=false)
 
 GRealSampleFFT (const u_int n, const double data_re[], const double data_im[])
 
GRealSampleFFToperator= (const GRealSampleFFT &original)
 
virtual ~GRealSampleFFT ()
 
Functions to initialise and set data for function sample or FT
virtual void InitFunctionData (u_int n, const double data[])
 
virtual void InitFunctionData (u_int n, double dt, const double data[])
 
virtual void InitFunctionData (u_int n, double t_lo, double t_hi, double(*fct_ptr)(double, const double[]), const double fct_par[], bool ctr=false)
 
virtual void InitFunctionData (u_int n, double t_lo, double t_hi, double(*fct_ptr)(const double[], const double[]), const double fct_par[], bool ctr=false)
 
virtual void InitTransformData (u_int n, const double data[])
 
virtual void InitTransformData (u_int n, const double data_re[], const double data_im[])
 
virtual void SetFunctionData (double(*fct_ptr)(double, const double[]), const double fct_par[])
 
virtual void SetFunctionData (double(*fct_ptr)(const double[], const double[]), const double fct_par[])
 
virtual void SetFullBandFilter (double amp=1.L)
 
virtual void SetFunctionValue (u_int i, double val)
 
virtual void SetTransformValue (u_int i, double re, double im=0.L)
 
virtual void SetFFTValue (u_int i, double re, double im=0.L)
 
virtual int SetFunctionValues (u_int i1, u_int i2, double val)
 
virtual int SetFunctionValues (u_int i1, u_int i2, const double data[])
 
virtual int SetFunctionValues (u_int i1, u_int i2, double(*fct_ptr)(double, const double[]), const double fct_par[])
 
virtual int SetFunctionValues (const u_int i1, const u_int i2, double(*fct_ptr)(const double[], const double[]), const double fct_par[])
 
virtual int ComputeTransformData ()
 
virtual int ComputeFunctionData ()
 
virtual void Normalise (double s=1.L)
 
Function to extract samples information
virtual double GetSum ()
 
virtual double GetMean ()
 
virtual double GetRMS ()
 
virtual double GetIntegral ()
 
virtual double GetPower (bool f=false)
 
virtual double GetMinimum ()
 
virtual double GetMinimum (double &t)
 
virtual double GetMaximum ()
 
virtual double GetMaximum (double &t)
 
virtual u_int GetMinimumIndex ()
 
virtual u_int GetMaximumIndex ()
 
GRealSampleFFTDerivate ()
 
GRealSampleFFTDerivate (GRealSampleFFT &sample)
 
GRealSampleFFTIntegrate ()
 
GRealSampleFFTIntegrate (GRealSampleFFT &sample)
 
Function Functions for convolution / deconvolution
virtual int InterpolateFrom (GRealSampleFFT &fct, double dt, double t0, int mode=0)
 
virtual int InterpolateFrom (GRealSampleFFT &fct, double dt, int mode=0)
 
virtual int SetTimeResponse (GRealSampleFFT &signal, GRealSampleFFT &resp, int mode=0)
 
virtual int SetConvolution (GRealSampleFFT &signal, GRealSampleFFT &resp, bool inv=false, int mode=0)
 
virtual int SetProduct (GRealSampleFFT &s1, GRealSampleFFT &s2)
 
Modification functions
virtual GRealSampleFFTDoubleSymmetric (bool neg=false)
 
virtual GRealSampleFFTHalfDivide ()
 
Arithmetic operators
GRealSampleFFT operator+ (GRealSampleFFT &s)
 
GRealSampleFFT operator- (GRealSampleFFT &s)
 
GRealSampleFFToperator+= (GRealSampleFFT &s)
 
GRealSampleFFToperator-= (GRealSampleFFT &s)
 
GRealSampleFFT operator* (GRealSampleFFT &s)
 
GRealSampleFFToperator*= (GRealSampleFFT &s)
 
GRealSampleFFT operator* (double r)
 
GRealSampleFFT operator/ (double r)
 
GRealSampleFFToperator*= (double r)
 
GRealSampleFFToperator/= (double r)
 
GRealSampleFFT operator+ (double k)
 
GRealSampleFFT operator- (double k)
 
GRealSampleFFToperator+= (double k)
 
GRealSampleFFToperator-= (double k)
 
- Public Member Functions inherited from GBaseSampleFFT
 GBaseSampleFFT (u_int n=0, double dt=1.L)
 
 GBaseSampleFFT (const GBaseSampleFFT &original)
 
GBaseSampleFFToperator= (const GBaseSampleFFT &original)
 
virtual ~GBaseSampleFFT ()
 
virtual void ResetData ()
 
virtual void ClearData ()
 
virtual u_int InitData (u_int n, double dt=1.L, bool clear=true)
 
virtual int ResizeData (u_int n)
 
virtual void Scale (double s)
 
virtual double GetSumSquared (bool f=false)
 
bool IsDefined () const
 
bool IsFunctionSet () const
 
bool IsTransformSet () const
 
void FunctionSet (bool b=true)
 
void TransformSet (bool b=true)
 
int GetDataType () const
 
u_int GetSamplesNumber () const
 
u_int GetDimension () const
 
bool HasSameSampling (const GBaseSampleFFT &sample, double prec=0.0001) const
 
virtual int SetTimeRange (double t_lo, double t_up, bool ctr=false)
 
virtual double GetTimeMin () const
 
virtual double GetTimeMax () const
 
virtual double GetTimeEnd () const
 
virtual double GetTimeStep () const
 
virtual double GetCriticalFrequency () const
 
virtual double GetFrequencyMax () const
 
virtual double GetFrequencyStep () const
 
double * GetFunctionData (bool check=true)
 
double * GetTransformData (bool check=true)
 
double * CopyFunctionData () const
 
double * CopyTransformData () const
 
double * CreatePowerSpectralDensity ()
 
double * CreatePSD ()
 
virtual int Write (const string &fname)
 
virtual int Write (FILE *fp)
 
virtual int Read (const string &fname)
 
virtual int Read (FILE *fp)
 

Private Member Functions

 GObject (GRealSampleFFT)
 

Additional Inherited Members

- Protected Member Functions inherited from GBaseSampleFFT
virtual int WriteHeader (FILE *fp)
 
virtual int ReadHeader (FILE *fp)
 
void CopyData (const GBaseSampleFFT &original)
 
- Protected Attributes inherited from GBaseSampleFFT
int data_type
 0=undefined; 1=real; 2=complex
 
u_int samples_number
 Number of data samples.
 
u_int samples_dim
 Effective arrays dimension (power of 2)
 
u_int alloc_dim
 Effective allocated dim.
 
double * data_array
 Array of function samples.
 
double * fft_array
 Array of Fourier transform data.
 
bool data_set
 Indicate if the data samples are set.
 
bool fft_set
 Indicate if the FFT data are set.
 
double time_min
 Time of the first sample.
 
double time_max
 Time of the last sample (after padding with zeros)
 
double time_end
 Time of the last (defined) sample (before padding)
 
double time_step
 Time step.
 
double freq_step
 Frequency step.
 
double freq_crit
 Frequency maximum value (Nyquist critical frequency)
 
- Static Protected Attributes inherited from GBaseSampleFFT
static double epsilon = 1.e-30L
 Limit for 0 checking.
 
static string data_type_name [3] = { "undefined", "real", "complex" }
 Names of the data types.
 

Detailed Description

This is a class to manage of a real function samples and the corresponding Fourier transform.

Constructor & Destructor Documentation

GRealSampleFFT::GRealSampleFFT ( u_int  n = 0,
double  dt = 1.L 
)

Base constructor.

Parameters
nnumber of function samples
dtsampling time step

References GBaseSampleFFT::data_type, GDebugConst, and GGetString().

Referenced by SetConvolution().

GRealSampleFFT::GRealSampleFFT ( const GRealSampleFFT original)

Copy constructor.

Parameters
originalobject to be copied

References GBaseSampleFFT::CopyData(), GBaseSampleFFT::data_type, and GDebugConst.

GRealSampleFFT::GRealSampleFFT ( u_int  n,
const double  data[] 
)

Constructor with function data initialisation. The data array must have at least a length n.

Parameters
nnumber of real data
dataarray of input data for real function samples

References GBaseSampleFFT::data_type, GDebugConst, InitFunctionData(), and GBaseSampleFFT::SetTimeRange().

GRealSampleFFT::GRealSampleFFT ( u_int  n,
double  t_lo,
double  t_hi,
const double  data[],
const bool  ctr = false 
)

Constructor with function data initialisation. The data array must have at least a length n. The time limits are given with respect to the defined sample (The effective time range will be adjusted if n is not a power of 2).

Parameters
nnumber of real data
t_lolower time limit of the defined sample
t_hiupper time limit of the defined sample
dataarray of input data for real function samples
ctrindicates how the time limits are defined (see time scale)

References GBaseSampleFFT::data_type, GDebugConst, GLogWarning(), InitFunctionData(), and GBaseSampleFFT::SetTimeRange().

GRealSampleFFT::GRealSampleFFT ( u_int  n,
double  t_lo,
double  t_hi,
double(*)(const double, const double[])  fct_ptr,
const double  fct_par[],
const bool  ctr = false 
)

Constructor with function pointer. The function is real with 2 arguments: the first one is the time coordinate, the second one is an array of parameters. The time limits are used to set the proper values of the function argument. The time limits are given with respect to the defined sample (The effective time range will be adjusted if n is not a power of 2).

Parameters
nnumber of real data
t_lolower time limit of the defined sample
t_hiupper time limit of the defined sample
fct_ptrfunction pointer
fct_pararray of parameters for the function
ctrindicates how the time limits are defined (see time scale)

References GBaseSampleFFT::data_type, GDebugConst, and InitFunctionData().

GRealSampleFFT::GRealSampleFFT ( u_int  n,
double  t_lo,
double  t_hi,
double(*)(const double[], const double[])  fct_ptr,
const double  fct_par[],
const bool  ctr = false 
)

Constructor with function pointer. The function is real with 2 arguments: the first one is the function variables array (only the first one is defined), the second one is an array of parameters. The time limits are used to set the proper values of the function argument. The time limits are given with respect to the defined sample (The effective time range will be adjusted if n is not a power of 2).

Parameters
nnumber of real data
t_lolower time limit of the defined sample
t_hiupper time limit of the defined sample
fct_ptrfunction pointer
fct_pararray of parameters for the function
ctrindicates how the time limits are defined (see time scale)

References GBaseSampleFFT::data_type, GDebugConst, and InitFunctionData().

GRealSampleFFT::GRealSampleFFT ( const u_int  n,
const double  data_re[],
const double  data_im[] 
)

Constructor with transform data initialisation. The corresponding real function has 2*n samples. Since the sampled function is real, Im[0] (that should be 0) must contain Re[n]. The data_re and data_im arrays must have at least a length n.

Parameters
nnumber of complex data
data_rereal component input (n real data)
data_imimaginary component input (n real data)

References GBaseSampleFFT::data_type, GDebugConst, InitTransformData(), and GBaseSampleFFT::SetTimeRange().

GRealSampleFFT::~GRealSampleFFT ( )
virtual

Destructor: cleans-up the data.

References GDebugDest.

Member Function Documentation

int GRealSampleFFT::ComputeFunctionData ( )
virtual

Compute the inverse Fourier transform (function samples) from the Fourier transform data that are currently stored in the object. The function returns 0 if no error occured (Fourier transform data are defined, ...).

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GDebugClass, GFFTReal(), GGetString(), GLogWarning(), GBaseSampleFFT::IsTransformSet(), and GBaseSampleFFT::samples_dim.

Referenced by Derivate(), DoubleSymmetric(), GetMaximum(), GetMaximumIndex(), GetMinimum(), GetMinimumIndex(), GetSum(), HalfDivide(), Integrate(), and InterpolateFrom().

int GRealSampleFFT::ComputeTransformData ( )
virtual

Compute the Fourier transform from the function samples that are currently stored in the object. The function returns 0 if no error occured (function samples are defined, ...).

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_array, GBaseSampleFFT::fft_array, GBaseSampleFFT::fft_set, GDebugClass, GFFTReal(), GGetString(), GLogWarning(), GBaseSampleFFT::IsFunctionSet(), and GBaseSampleFFT::samples_dim.

GRealSampleFFT & GRealSampleFFT::DoubleSymmetric ( bool  neg = false)
virtual

Function that expands the sample to double its dimension and copy the first half as mirror symmetry in the second half.

Parameters
negif true, the second half contains the opposite of the first half

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsFunctionSet(), GBaseSampleFFT::IsTransformSet(), GBaseSampleFFT::samples_dim, GBaseSampleFFT::time_step, and GBaseSampleFFT::TransformSet().

double GRealSampleFFT::GetIntegral ( )
inlinevirtual

Compute the sample integral: sum of the real data multiplied by the time step.

References GetSum(), and GBaseSampleFFT::GetTimeStep().

double GRealSampleFFT::GetMaximum ( )
virtual

Return the maximum sample value. The function data is computed if needed.

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().

double GRealSampleFFT::GetMaximum ( double &  t)
virtual

Return the maximum sample value, and the corresponding time. The function data is computed if needed.

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::time_min, and GBaseSampleFFT::time_step.

u_int GRealSampleFFT::GetMaximumIndex ( )
virtual

Return the position (index) of the maximum sample value. The function data is computed if needed.

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().

double GRealSampleFFT::GetMean ( )
inlinevirtual

Compute the sample average value over the samples.

References GetSum(), and GBaseSampleFFT::samples_number.

Referenced by GetRMS().

double GRealSampleFFT::GetMinimum ( )
virtual

Return the minimum sample value. The function data is computed if needed.

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().

double GRealSampleFFT::GetMinimum ( double &  t)
virtual

Return the minimum sample value, and the corresponding time. The function data is computed if needed.

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::time_min, and GBaseSampleFFT::time_step.

u_int GRealSampleFFT::GetMinimumIndex ( )
virtual

Return the position (index) of the minimum sample value. The function data is computed if needed.

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().

double GRealSampleFFT::GetPower ( bool  f = false)
inlinevirtual

Return the total power of the sample (time-integral squared amplitude that does (almost) not depend on the sampling).

Parameters
findicates if calculation is done in frequency space (default is time), but the result is the same...

References GBaseSampleFFT::GetSumSquared(), and GBaseSampleFFT::GetTimeStep().

double GRealSampleFFT::GetRMS ( )
virtual

Compute the RMS of the real data.

References GBaseSampleFFT::data_array, GetMean(), and GBaseSampleFFT::GetSamplesNumber().

double GRealSampleFFT::GetSum ( )
virtual
GRealSampleFFT::GObject ( GRealSampleFFT  )
private

Macro from GCpp library that defines the following functions:

  • ClassName(): return the real class name of the object
  • StaticClassName(): return the used class name of the object (that may be a base class of the real object).
  • IsInstanceOf<T>(): return true if the current object is an instance of the template class name argument
  • Clone(): return an allocated copy of the object.
GRealSampleFFT & GRealSampleFFT::HalfDivide ( )
virtual
void GRealSampleFFT::InitFunctionData ( u_int  n,
const double  data[] 
)
inlinevirtual

Set the function samples to be n real data (double). The data array must have at least a length n. The function does not change the sampling time step.

After this function call, the functions samples are set and the Fourier transform data are not set.

Parameters
nnumber of real data
dataarray of input data

Implements GBaseSampleFFT.

References GBaseSampleFFT::time_step.

Referenced by GRealSampleFFT().

void GRealSampleFFT::InitFunctionData ( u_int  n,
double  dt,
const double  data[] 
)
virtual

Set the function samples to be n real data (double). The data array must have at least a length n.

After this function call, the functions samples are set and the Fourier transform data are not set.

Parameters
nnumber of real data
dttime step
dataarray of input data

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GDebugClass, GGetString(), GBaseSampleFFT::InitData(), and GBaseSampleFFT::samples_number.

void GRealSampleFFT::InitFunctionData ( u_int  n,
double  t_lo,
double  t_hi,
double(*)(double, const double[])  fct_ptr,
const double  fct_par[],
bool  ctr = false 
)
virtual

Set the function samples to be n real data (double) from a function defined by its pointer (see SetFunctionData). The function is real with 2 arguments: the first one is the time coordinate, the second one is an array of parameters. The time limits are used to set the proper values of the function argument.

After this function call, the functions samples are set and the Fourier transform data are not set.

Parameters
nnumber of real data
t_lolower time limit of the sample
t_hiupper time limit of the sample
fct_ptrfunction pointer
fct_pararray of parameters for the function
ctrindicates if the time limits are defined as the first and last points times (n-1 intervals, default case) or as a binning times if set to true (n intervals, points correspnding to the bin centers)

References GDebugClass, GGetString(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::samples_number, SetFunctionData(), and GBaseSampleFFT::SetTimeRange().

void GRealSampleFFT::InitFunctionData ( u_int  n,
double  t_lo,
double  t_hi,
double(*)(const double[], const double[])  fct_ptr,
const double  fct_par[],
bool  ctr = false 
)
virtual

Set the function samples to be n real data (double) from a function defined by its pointer (see SetFunctionData). The function is real with 2 arguments: the first one is the function variables array (only the first one is defined), the second one is an array of parameters. The time limits are used to set the proper values of the function argument.

After this function call, the functions samples are set and the Fourier transform data are not set.

Parameters
nnumber of real data
t_lolower time limit of the sample
t_hiupper time limit of the sample
fct_ptrfunction pointer
fct_pararray of parameters for the function
ctrindicates how the time limits are defined (see time scale)

References GDebugClass, GGetString(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::samples_number, SetFunctionData(), and GBaseSampleFFT::SetTimeRange().

void GRealSampleFFT::InitTransformData ( u_int  n,
const double  data[] 
)
virtual

Set the transform data to be n/2 complex data (2 double values for each data). The data array is in the following order: Re[0], Im[0], Re[1], Im[1], ... Re[n/2-1], Im[n/2-1]. Since the sampled function is real, Im[0] (that should be 0) must contain Re[n/2]. The data array must have at least a length n.

After this function call, the Fourier transform data are set and the functions samples are not set.

Parameters
nnumber of real data
dataarray of input transform data (n/2 complex data)

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GDebugClass, GGetString(), GBaseSampleFFT::InitData(), and GBaseSampleFFT::samples_number.

Referenced by GRealSampleFFT().

void GRealSampleFFT::InitTransformData ( u_int  n,
const double  data_re[],
const double  data_im[] 
)
virtual

Set the transform data to be n complex data (2 arrays of values for real and imaginary parts of each data). The corresponding real function has 2*n samples. Since the sampled function is real, Im[0] (that should be 0) must contain Re[n]. The data_re and data_im arrays must have at least a length n.

After this function call, the Fourier transform data are set and the functions samples are not set.

Parameters
nnumber of complex data
data_rereal component input (n real data)
data_imimaginary component input (n real data)

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GDebugClass, GGetString(), GBaseSampleFFT::InitData(), and GBaseSampleFFT::samples_number.

int GRealSampleFFT::InterpolateFrom ( GRealSampleFFT fct,
double  dt,
double  t0,
int  mode = 0 
)
virtual

Set the data from the interpolation of the function sample given as argument.

No point outside the input sample time range is computed. In case of a cubic interpolation mode, the derivative are kept constant at each point (but the interpolation is only quadratic for the first and last intervals). The function returns 0 if no error occured.

Parameters
fctinput function sample (to interpolate from)
dtnew time step
t0reference time (that will be exactly on a sample)
modeinterpolation mode : 0=linear, 1=cubic

References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GArrayInterpolationCubic(), GArrayInterpolationLinear(), GDebugClass, GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeEnd(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GGetString(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsDefined(), and GBaseSampleFFT::SetTimeRange().

Referenced by InterpolateFrom(), and SetConvolution().

int GRealSampleFFT::InterpolateFrom ( GRealSampleFFT fct,
double  dt,
int  mode = 0 
)
virtual

Set the data from the interpolation of the function sample given as argument. The result sample has the same first point than the argument sample.

No point outside the input sample time range is computed. In case of a cubic interpolation mode, the derivative are kept constant at each point (but the interpolation is only quadratic for the first and last intervals). The function returns 0 if no error occured.

Parameters
fctinput function sample (to interpolate from)
dtnew time step
modeinterpolation mode : 0=linear, 1=cubic

References GBaseSampleFFT::GetTimeMin(), and InterpolateFrom().

void GRealSampleFFT::Normalise ( double  s = 1.L)
virtual

Normalisation of the real samples sum.

Parameters
srequested normalisation

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_array, and GBaseSampleFFT::samples_dim.

GRealSampleFFT GRealSampleFFT::operator* ( GRealSampleFFT s)

Product of the sample with another sample (see SetProduct function).

Parameters
soperand sample

References SetProduct().

GRealSampleFFT & GRealSampleFFT::operator*= ( GRealSampleFFT s)

Product of the sample with another sample (see SetProduct function).

Parameters
soperand sample

References SetProduct().

GRealSampleFFT & GRealSampleFFT::operator*= ( double  r)

Scaling (multiplication factor) of the current sample.

Parameters
rscale factor

References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GLogWarning(), GBaseSampleFFT::IsDefined(), GBaseSampleFFT::samples_number, and GBaseSampleFFT::TransformSet().

GRealSampleFFT GRealSampleFFT::operator+ ( GRealSampleFFT s)

Sum of samples. The samples must have the same number of data. WARNING: there is no check of the sampling time is the same !

Parameters
ssample to be summed to the current one

References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::samples_number, GBaseSampleFFT::SetTimeRange(), GBaseSampleFFT::time_end, GBaseSampleFFT::time_min, and GBaseSampleFFT::TransformSet().

GRealSampleFFT & GRealSampleFFT::operator+= ( GRealSampleFFT s)

Add a sample to the current one. If the added sample has a different size than the current one, only the available data, up to maximum the current sample size, are added. WARNING: there is no check of the sampling time is the same !

Parameters
ssample to be summed to the current one

References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GLogWarning(), GBaseSampleFFT::samples_number, and GBaseSampleFFT::TransformSet().

GRealSampleFFT & GRealSampleFFT::operator+= ( double  k)
GRealSampleFFT GRealSampleFFT::operator- ( GRealSampleFFT s)

Difference of samples. The samples must have the same number of data. WARNING: there is no check of the sampling time is the same !

Parameters
ssample to be summed to the current one

References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::samples_number, GBaseSampleFFT::SetTimeRange(), GBaseSampleFFT::time_end, GBaseSampleFFT::time_min, and GBaseSampleFFT::TransformSet().

GRealSampleFFT GRealSampleFFT::operator- ( double  k)
inline

Subtracting a constant offset to current sample.

Parameters
kconstant offset
GRealSampleFFT & GRealSampleFFT::operator-= ( GRealSampleFFT s)

Subtract a sample to the current one. If the added sample has a different size than the current one, only the available data, up to maximum the current sample size, are added. WARNING: there is no check of the sampling time is the same !

Parameters
ssample to be summed to the current one

References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GLogWarning(), GBaseSampleFFT::samples_number, and GBaseSampleFFT::TransformSet().

GRealSampleFFT & GRealSampleFFT::operator-= ( double  k)
inline

Subtracting an offset to the current sample.

Parameters
kconstant offset
GRealSampleFFT GRealSampleFFT::operator/ ( double  r)
inline

Scaling operation the sample (division factor).

Parameters
rscale factor
GRealSampleFFT & GRealSampleFFT::operator/= ( double  r)
inline

Scaling operation the sample (division factor).

Parameters
rscale factor
GRealSampleFFT & GRealSampleFFT::operator= ( const GRealSampleFFT original)

Affectation operator.

Parameters
originalobject to be copied

References GBaseSampleFFT::CopyData(), G_AFFECT_MOTHER, and GDebugClass.

int GRealSampleFFT::SetConvolution ( GRealSampleFFT signal,
GRealSampleFFT resp,
bool  inv = false,
int  mode = 0 
)
virtual

Set data from the (de)convolution of the signal (S) and the response (R), using the Fourier transform. The result (0) verifies:

  • O=S*R for a convolution
  • S=0*R for a deconvolution

The response must have the same sampling than the signal. If it is not the case, the response function is first interpolated from the original one.

If the total number of samples is different for signal and response, the shortest one is padded with 0s to get the same dimension (It is set back to it's original dimension after processing).

The function returns 0 if no error. If the result is negative, this means there are 0's in the response function.

The output has its FFT set, but the function data is not set.

The following options can be choosen for the convolution (mode argument):

  • mode = 0 : proceed if signal and response have same sampling
  • mode = 1 : linear interpolation of response to signal sampling
  • mode = 2 : cubic interpolation of response to signal sampling
  • mode = 3 : ignore response sampling and consider it as the signal sampling (without check), with first point of the response data being at time 0

Note that the signal sample or the resp sample can be the current sample object, in which case it will be overwritten by the result.

Parameters
signalinput sampled signal
respsampled response function
invconvolution (false) or deconvolution (true)
modeconvolution mode option

References GBaseSampleFFT::data_set, GBaseSampleFFT::epsilon, GBaseSampleFFT::fft_array, GBaseSampleFFT::fft_set, GBaseSampleFFT::GetDimension(), GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GBaseSampleFFT::GetTransformData(), GLogWarning(), GRealSampleFFT(), GBaseSampleFFT::InitData(), InterpolateFrom(), GBaseSampleFFT::IsDefined(), GBaseSampleFFT::IsFunctionSet(), GBaseSampleFFT::IsTransformSet(), GBaseSampleFFT::ResizeData(), and GBaseSampleFFT::SetTimeRange().

void GRealSampleFFT::SetFFTValue ( u_int  i,
double  re,
double  im = 0.L 
)
inlinevirtual

Set a single (complex) value for the FFT. This function is similar to SetTransformValue, but it contains less index check: the index must be in the range from 0 to (n/2)-1, and the data for (n/2) that is pure real is set using the imaginary part for data index 0 (as effectively used in the FFT functions).

The use of this function implies that the full FFT is set using this function, and that the function samples have to be recomputed.

Parameters
iindex of the FFT value
rereal part of the data
imimaginary part of the data

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_dim.

void GRealSampleFFT::SetFullBandFilter ( double  amp = 1.L)
virtual

Function that sets the data for a full band filter, (or a dirac in time space), with optional amplification. Note that there is a division by the time staep, to get conservative integrals.

Parameters
ampamplification factor

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GBaseSampleFFT::samples_dim, GBaseSampleFFT::samples_number, and GBaseSampleFFT::time_step.

void GRealSampleFFT::SetFunctionData ( double(*)(double, const double[])  fct_ptr,
const double  fct_par[] 
)
virtual

Set the function full sample from a function defined by its pointer. The function is real with 2 arguments: the first one is the time coordinate, the second one is an array of parameters. The time limits are used to set the proper values of the function argument.

After this function call, the functions samples are set and the Fourier transform data are not set.

Parameters
fct_ptrfunction pointer
fct_pararray of parameters for the function

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GDebugClass, GBaseSampleFFT::samples_number, GBaseSampleFFT::time_min, and GBaseSampleFFT::time_step.

Referenced by InitFunctionData().

void GRealSampleFFT::SetFunctionData ( double(*)(const double[], const double[])  fct_ptr,
const double  fct_par[] 
)
virtual

Set the function full sample from a function defined by its pointer. The function is real with 2 arguments: the first one is the time coordinate, the second one is an array of parameters. The time limits are used to set the proper values of the function argument.

After this function call, the functions samples are set and the Fourier transform data are not set.

Parameters
fct_ptrfunction pointer
fct_pararray of parameters for the function

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GDebugClass, GGetString(), GBaseSampleFFT::samples_number, GBaseSampleFFT::time_min, and GBaseSampleFFT::time_step.

void GRealSampleFFT::SetFunctionValue ( u_int  i,
double  val 
)
inlinevirtual

Function to set a single value in the sample. The index must range between 0 and n-1.

The use of this function implies that the full sample is set using this function, and that the FFT data have to be recomputed.

Parameters
idata index
valdata (real) value

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_dim.

int GRealSampleFFT::SetFunctionValues ( u_int  i1,
u_int  i2,
double  val 
)
virtual

Function that sets a part of the sample array values, from index i1 to index i2 (included) with the same value. The function returns the number of values that have been set.

The use of this function implies that the sample is defined through the function data and that the FFT data have to be recomputed.

Parameters
i1lower index of the sample data
i2upper index of the sample data
valvalue to set in the function

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_number.

int GRealSampleFFT::SetFunctionValues ( u_int  i1,
u_int  i2,
const double  data[] 
)
virtual

Function that sets a part of the sample array values, from index i1 to index i2 (included) with the data array values. The function returns the number of values that have been set.

The use of this function implies that the sample is defined through the function data and that the FFT data have to be recomputed.

Parameters
i1lower index of the sample data
i2upper index of the sample data
dataarray of values (must be of size at least (i2-i1+1)

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_number.

int GRealSampleFFT::SetFunctionValues ( u_int  i1,
u_int  i2,
double(*)(double, const double[])  fct_ptr,
const double  fct_par[] 
)
virtual

Function that sets a part of the sample array values, from index i1 to index i2 (included) with the values computed from the function (time is from currently defined time sampling). The function returns the number of values that have been set.

The use of this function implies that the sample is defined through the function data and that the FFT data have to be recomputed.

Parameters
i1lower index of the sample data
i2upper index of the sample data
fct_ptrfunction pointer
fct_pararray of parameters for the function

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GBaseSampleFFT::samples_number, GBaseSampleFFT::time_min, and GBaseSampleFFT::time_step.

int GRealSampleFFT::SetFunctionValues ( const u_int  i1,
const u_int  i2,
double(*)(const double[], const double[])  fct_ptr,
const double  fct_par[] 
)
virtual

Function that sets a part of the sample array values, from index i1 to index i2 (included) with the values computed from the function (time is from currently defined time sampling). The function returns the number of values that have been set.

The use of this function implies that the sample is defined through the function data and that the FFT data have to be recomputed.

Parameters
i1lower index of the sample data
i2upper index of the sample data
fct_ptrfunction pointer
fct_pararray of parameters for the function

References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GBaseSampleFFT::samples_number, GBaseSampleFFT::time_min, and GBaseSampleFFT::time_step.

int GRealSampleFFT::SetProduct ( GRealSampleFFT s1,
GRealSampleFFT s2 
)
virtual

Compute the product of 2 samples.

The s1 and s2 must have the same sampling. The function returns 0 if no error.

The output has its function data set, but the FFT data is not set.

Note that the s1 sample or the s2 sample can be the current sample object, in which case it will be overwritten by the result.

Parameters
s1first signal
s2second signal

References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GBaseSampleFFT::GetDimension(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetTimeStep(), GLogWarning(), GBaseSampleFFT::HasSameSampling(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsDefined(), GBaseSampleFFT::IsFunctionSet(), and GBaseSampleFFT::IsTransformSet().

Referenced by operator*(), and operator*=().

int GRealSampleFFT::SetTimeResponse ( GRealSampleFFT signal,
GRealSampleFFT resp,
int  mode = 0 
)
virtual

Sets the current object data from the convolution of a sampled signal with a response function. If signal has Ns samples and response has Nr samples (with same time step), the result is significant for Ns+Nr-1 samples (the result sample is padded to a dimension that is a power of 2). The function returns 0 if no error occured.

The output has its function data set, but the FFT is not set.

The following options can be choosen for the convolution (mode argument):

  • mode = 0 : proceed if signal and response have same sampling
  • mode = 1 : linear interpolation of response to signal sampling
  • mode = 2 : cubic interpolation of response to signal sampling
  • mode = 3 : ignore response sampling and consider it as the signal sampling (without check), with first point of the response data being at time 0

Note that interpolations require to recompute the response function and slows the calculation. It is faster to use mode 0 or 3 depending if a compatibility test is required or not.

Parameters
signalinput sampled signal
respsampled response function
modeconvolution mode option

References GBaseSampleFFT::data_array, GBaseSampleFFT::fft_set, GArrayInterpolationCubic(), GArrayInterpolationLinear(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeEnd(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GGetString(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsDefined(), GBaseSampleFFT::IsFunctionSet(), GBaseSampleFFT::IsTransformSet(), and GBaseSampleFFT::SetTimeRange().

void GRealSampleFFT::SetTransformValue ( u_int  i,
double  re,
double  im = 0.L 
)
inlinevirtual

Set a single (complex) value for the FFT. The index must be in the range 0 to (n/2). For the index values from 0 and (n/2), only the real part of the value is effective.

The use of this function implies that the full FFT is set using this function, and that the function samples have to be recomputed.

Parameters
iindex of the FFT value
rereal part of the data
imimaginary part of the data

Implements GBaseSampleFFT.

References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_dim.


The documentation for this class was generated from the following files: