GCpp general purpose C++ library
version 1.0
|
#include <GRealSampleFFT.hh>
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[]) | |
GRealSampleFFT & | operator= (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 () |
GRealSampleFFT & | Derivate () |
GRealSampleFFT & | Derivate (GRealSampleFFT &sample) |
GRealSampleFFT & | Integrate () |
GRealSampleFFT & | Integrate (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 GRealSampleFFT & | DoubleSymmetric (bool neg=false) |
virtual GRealSampleFFT & | HalfDivide () |
Arithmetic operators | |
GRealSampleFFT | operator+ (GRealSampleFFT &s) |
GRealSampleFFT | operator- (GRealSampleFFT &s) |
GRealSampleFFT & | operator+= (GRealSampleFFT &s) |
GRealSampleFFT & | operator-= (GRealSampleFFT &s) |
GRealSampleFFT | operator* (GRealSampleFFT &s) |
GRealSampleFFT & | operator*= (GRealSampleFFT &s) |
GRealSampleFFT | operator* (double r) |
GRealSampleFFT | operator/ (double r) |
GRealSampleFFT & | operator*= (double r) |
GRealSampleFFT & | operator/= (double r) |
GRealSampleFFT | operator+ (double k) |
GRealSampleFFT | operator- (double k) |
GRealSampleFFT & | operator+= (double k) |
GRealSampleFFT & | operator-= (double k) |
![]() | |
GBaseSampleFFT (u_int n=0, double dt=1.L) | |
GBaseSampleFFT (const GBaseSampleFFT &original) | |
GBaseSampleFFT & | operator= (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 | |
![]() | |
virtual int | WriteHeader (FILE *fp) |
virtual int | ReadHeader (FILE *fp) |
void | CopyData (const GBaseSampleFFT &original) |
![]() | |
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 double | epsilon = 1.e-30L |
Limit for 0 checking. | |
static string | data_type_name [3] = { "undefined", "real", "complex" } |
Names of the data types. | |
This is a class to manage of a real function samples and the corresponding Fourier transform.
GRealSampleFFT::GRealSampleFFT | ( | u_int | n = 0 , |
double | dt = 1.L |
||
) |
Base constructor.
n | number of function samples |
dt | sampling time step |
References GBaseSampleFFT::data_type, GDebugConst, and GGetString().
Referenced by SetConvolution().
GRealSampleFFT::GRealSampleFFT | ( | const GRealSampleFFT & | original | ) |
Copy constructor.
original | object 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.
n | number of real data |
data | array 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).
n | number of real data |
t_lo | lower time limit of the defined sample |
t_hi | upper time limit of the defined sample |
data | array of input data for real function samples |
ctr | indicates 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).
n | number of real data |
t_lo | lower time limit of the defined sample |
t_hi | upper time limit of the defined sample |
fct_ptr | function pointer |
fct_par | array of parameters for the function |
ctr | indicates 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).
n | number of real data |
t_lo | lower time limit of the defined sample |
t_hi | upper time limit of the defined sample |
fct_ptr | function pointer |
fct_par | array of parameters for the function |
ctr | indicates 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.
n | number of complex data |
data_re | real component input (n real data) |
data_im | imaginary component input (n real data) |
References GBaseSampleFFT::data_type, GDebugConst, InitTransformData(), and GBaseSampleFFT::SetTimeRange().
|
virtual |
Destructor: cleans-up the data.
References GDebugDest.
|
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().
|
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::Derivate | ( | ) |
Set the current sample to its derivative.
References ComputeFunctionData(), GBaseSampleFFT::data_set, GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GLogWarning(), and GBaseSampleFFT::SetTimeRange().
Referenced by GExtremaFinder::Analyse().
GRealSampleFFT & GRealSampleFFT::Derivate | ( | GRealSampleFFT & | sample | ) |
Set the current sample as the derivative of the argument.
sample | initial sample to derivate |
References ComputeFunctionData(), GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsFunctionSet(), GBaseSampleFFT::SetTimeRange(), and GBaseSampleFFT::TransformSet().
|
virtual |
Function that expands the sample to double its dimension and copy the first half as mirror symmetry in the second half.
neg | if 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().
|
inlinevirtual |
Compute the sample integral: sum of the real data multiplied by the time step.
References GetSum(), and GBaseSampleFFT::GetTimeStep().
|
virtual |
Return the maximum sample value. The function data is computed if needed.
References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().
|
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.
|
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().
|
inlinevirtual |
Compute the sample average value over the samples.
References GetSum(), and GBaseSampleFFT::samples_number.
Referenced by GetRMS().
|
virtual |
Return the minimum sample value. The function data is computed if needed.
References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().
|
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.
|
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().
|
inlinevirtual |
Return the total power of the sample (time-integral squared amplitude that does (almost) not depend on the sampling).
f | indicates if calculation is done in frequency space (default is time), but the result is the same... |
References GBaseSampleFFT::GetSumSquared(), and GBaseSampleFFT::GetTimeStep().
|
virtual |
Compute the RMS of the real data.
References GBaseSampleFFT::data_array, GetMean(), and GBaseSampleFFT::GetSamplesNumber().
|
virtual |
Compute the sum of the real data.
References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, and GBaseSampleFFT::GetSamplesNumber().
Referenced by GetIntegral(), and GetMean().
|
private |
Macro from GCpp library that defines the following functions:
|
virtual |
Cut the sample in two, and keep only the first half. The function operates only if the dimension is larger than 2 !
References ComputeFunctionData(), GBaseSampleFFT::data_array, GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsFunctionSet(), GBaseSampleFFT::IsTransformSet(), GBaseSampleFFT::samples_dim, GBaseSampleFFT::time_step, and GBaseSampleFFT::TransformSet().
|
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.
n | number of real data |
data | array of input data |
Implements GBaseSampleFFT.
References GBaseSampleFFT::time_step.
Referenced by GRealSampleFFT().
|
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.
n | number of real data |
dt | time step |
data | array of input data |
References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GDebugClass, GGetString(), GBaseSampleFFT::InitData(), and GBaseSampleFFT::samples_number.
|
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.
n | number of real data |
t_lo | lower time limit of the sample |
t_hi | upper time limit of the sample |
fct_ptr | function pointer |
fct_par | array of parameters for the function |
ctr | indicates 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().
|
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.
n | number of real data |
t_lo | lower time limit of the sample |
t_hi | upper time limit of the sample |
fct_ptr | function pointer |
fct_par | array of parameters for the function |
ctr | indicates how the time limits are defined (see time scale) |
References GDebugClass, GGetString(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::samples_number, SetFunctionData(), and GBaseSampleFFT::SetTimeRange().
|
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.
n | number of real data |
data | array 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().
|
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.
n | number of complex data |
data_re | real component input (n real data) |
data_im | imaginary component input (n real data) |
Implements GBaseSampleFFT.
References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GDebugClass, GGetString(), GBaseSampleFFT::InitData(), and GBaseSampleFFT::samples_number.
GRealSampleFFT & GRealSampleFFT::Integrate | ( | ) |
Set the current sample to its integral.
References ComputeFunctionData(), GBaseSampleFFT::data_set, GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GLogWarning(), and GBaseSampleFFT::SetTimeRange().
GRealSampleFFT & GRealSampleFFT::Integrate | ( | GRealSampleFFT & | sample | ) |
Set the current sample to the integral of the argument.
sample | initial sample to integrate |
References ComputeFunctionData(), GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GBaseSampleFFT::GetTimeMin(), GBaseSampleFFT::GetTimeStep(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsFunctionSet(), GBaseSampleFFT::SetTimeRange(), and GBaseSampleFFT::TransformSet().
|
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.
fct | input function sample (to interpolate from) |
dt | new time step |
t0 | reference time (that will be exactly on a sample) |
mode | interpolation 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().
|
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.
fct | input function sample (to interpolate from) |
dt | new time step |
mode | interpolation mode : 0=linear, 1=cubic |
References GBaseSampleFFT::GetTimeMin(), and InterpolateFrom().
|
virtual |
Normalisation of the real samples sum.
s | requested 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).
s | operand sample |
References SetProduct().
GRealSampleFFT GRealSampleFFT::operator* | ( | double | r | ) |
Scaling (multiplication factor) of the current sample.
r | scale factor |
References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsDefined(), GBaseSampleFFT::samples_number, GBaseSampleFFT::SetTimeRange(), GBaseSampleFFT::time_end, GBaseSampleFFT::time_min, and GBaseSampleFFT::TransformSet().
GRealSampleFFT & GRealSampleFFT::operator*= | ( | GRealSampleFFT & | s | ) |
Product of the sample with another sample (see SetProduct function).
s | operand sample |
References SetProduct().
GRealSampleFFT & GRealSampleFFT::operator*= | ( | double | r | ) |
Scaling (multiplication factor) of the current sample.
r | scale 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 !
s | sample 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 | ) |
Summing current sample with a constant offset.
k | constant offset |
References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GLogWarning(), GBaseSampleFFT::InitData(), GBaseSampleFFT::IsDefined(), 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 !
s | sample 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 | ) |
Adding an offset to the current sample.
k | constant offset |
References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GLogWarning(), GBaseSampleFFT::IsDefined(), GBaseSampleFFT::samples_number, and GBaseSampleFFT::TransformSet().
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 !
s | sample 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().
|
inline |
Subtracting a constant offset to current sample.
k | constant 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 !
s | sample to be summed to the current one |
References GBaseSampleFFT::FunctionSet(), GBaseSampleFFT::GetFunctionData(), GBaseSampleFFT::GetSamplesNumber(), GLogWarning(), GBaseSampleFFT::samples_number, and GBaseSampleFFT::TransformSet().
|
inline |
Subtracting an offset to the current sample.
k | constant offset |
|
inline |
Scaling operation the sample (division factor).
r | scale factor |
|
inline |
Scaling operation the sample (division factor).
r | scale factor |
GRealSampleFFT & GRealSampleFFT::operator= | ( | const GRealSampleFFT & | original | ) |
Affectation operator.
original | object to be copied |
References GBaseSampleFFT::CopyData(), G_AFFECT_MOTHER, and GDebugClass.
|
virtual |
Set data from the (de)convolution of the signal (S) and the response (R), using the Fourier transform. The result (0) verifies:
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):
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.
signal | input sampled signal |
resp | sampled response function |
inv | convolution (false) or deconvolution (true) |
mode | convolution 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().
|
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.
i | index of the FFT value |
re | real part of the data |
im | imaginary part of the data |
Implements GBaseSampleFFT.
References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_dim.
|
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.
amp | amplification factor |
References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, GBaseSampleFFT::samples_dim, GBaseSampleFFT::samples_number, and GBaseSampleFFT::time_step.
|
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.
fct_ptr | function pointer |
fct_par | array 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().
|
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.
fct_ptr | function pointer |
fct_par | array 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.
|
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.
i | data index |
val | data (real) value |
Implements GBaseSampleFFT.
References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_dim.
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.
i1 | lower index of the sample data |
i2 | upper index of the sample data |
val | value to set in the function |
References GBaseSampleFFT::data_array, GBaseSampleFFT::data_set, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_number.
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.
i1 | lower index of the sample data |
i2 | upper index of the sample data |
data | array 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.
|
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.
i1 | lower index of the sample data |
i2 | upper index of the sample data |
fct_ptr | function pointer |
fct_par | array 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.
|
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.
i1 | lower index of the sample data |
i2 | upper index of the sample data |
fct_ptr | function pointer |
fct_par | array 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.
|
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.
s1 | first signal |
s2 | second 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*=().
|
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):
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.
signal | input sampled signal |
resp | sampled response function |
mode | convolution 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().
|
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.
i | index of the FFT value |
re | real part of the data |
im | imaginary part of the data |
Implements GBaseSampleFFT.
References GBaseSampleFFT::data_set, GBaseSampleFFT::fft_array, GBaseSampleFFT::fft_set, and GBaseSampleFFT::samples_dim.