JG ROOT Tools libraries  version 5.0 - august 2014
RRealSampleFFT Class Reference

#include <RRealSampleFFT.hh>

Inheritance diagram for RRealSampleFFT:

Public Member Functions

 RRealSampleFFT (const u_int n=0, const double dt=1.L)
 
 RRealSampleFFT (const RRealSampleFFT &original)
 
 RRealSampleFFT (const GRealSampleFFT &original)
 
 RRealSampleFFT (const u_int n, const double data[])
 
 RRealSampleFFT (const u_int n, const double t_lo, const double t_hi, const double data[], bool ctr=false)
 
 RRealSampleFFT (const u_int n, const double t_lo, const double t_hi, double(*fct_ptr)(const double, const double[]), const double fct_par[], bool ctr=false)
 
 RRealSampleFFT (const u_int n, const double t_lo, const double t_hi, double(*fct_ptr)(const double[], const double[]), const double fct_par[], bool ctr=false)
 
 RRealSampleFFT (const u_int n, const double t_lo, const double t_hi, TF1 *fct, bool ctr=false)
 
 RRealSampleFFT (TH1 *hptr)
 
 RRealSampleFFT (const u_int n, const double data_re[], const double data_im[])
 
RRealSampleFFToperator= (const RRealSampleFFT &original)
 
virtual ~RRealSampleFFT ()
 
virtual void InitFunctionData (const u_int n, const double t_lo, const double t_hi, TF1 *fct, bool ctr=false)
 
virtual void InitFunctionData (TH1 *hptr)
 
virtual void SetFunctionData (TF1 *fct)
 
virtual void FillRandomGaus (double mean, double fwhm, TRandom *rnd=NULL)
 
virtual void FillRandomGaus (double fwhm, TRandom *rnd=NULL)
 
virtual int SetFunctionValues (u_int i1, u_int i2, TF1 *fct)
 
void SetDefaultValues ()
 
void SetColor (const Color_t col)
 
void SetFunctionTitles (const char tit[]="", const char tit_x[]="", const char tit_y[]="")
 
void SetTransformTitles (const char tit[]="", const char tit_f[]="", const char tit_amp[]="", const char tit_ph[]="", const char tit_pow[]="")
 
void SetTransformComplexTitles (const char tit_re[]="", const char tit_im[]="")
 
TGraph * CreateFunctionGraph (double scale=1.L)
 
TGraph * CreateTransformGraph (double scale=1.L)
 
TGraph * CreateTrfAmpliGraph (double scale=1.L)
 
TGraph * CreateTrfPhaseGraph (bool deg=false)
 
TGraph * CreatePowerGraph (double scale=1.L)
 
TH1D * CreateFunctionHisto (const string &name, double scale=1.L)
 
int FillFunctionHisto (TH1D *hptr, double scale=1.L, bool reset=true)
 
RRealSampleFFT operator+ (RRealSampleFFT &s)
 
RRealSampleFFT operator- (RRealSampleFFT &s)
 
RRealSampleFFT operator* (RRealSampleFFT &s)
 
RRealSampleFFT operator/ (RRealSampleFFT &s)
 
RRealSampleFFT operator* (double r)
 
RRealSampleFFT operator/ (double r)
 
RRealSampleFFT operator+ (double r)
 
RRealSampleFFT operator- (double r)
 
 ClassDef (RRealSampleFFT, 0)
 

Public Attributes

string fct_title
 Function title.
 
string fct_title_x
 Function title (X axis)
 
string fct_title_y
 Function title (Y axis)
 
string fft_title
 Fourier transform title.
 
string fft_title_freq
 Fourier transform frequency axis title.
 
string fft_title_amp
 Fourier transform amplitude axis title.
 
string fft_title_phase
 Fourier transform phase axis title.
 
string fft_title_power
 Fourier transform power axis title.
 
string fft_title_re
 Fourier transform real axis title.
 
string fft_title_im
 Fourier transform imaginary axis title.
 

Protected Member Functions

void CopyData (const RRealSampleFFT &original)
 

Private Member Functions

 GObject (RRealSampleFFT)
 

Detailed Description

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

This class is based on GRealSampleFFT general class, with features for ROOT. It defines some information for graphs and histograms to be created for ROOT.

Constructor & Destructor Documentation

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

Base constructor.

Parameters
nnumber of function samples
dtsampling time step
RRealSampleFFT::RRealSampleFFT ( const RRealSampleFFT original)

Copy constructor.

Parameters
originalobject to be copied

References CopyData().

RRealSampleFFT::RRealSampleFFT ( const GRealSampleFFT &  original)

Copy constructor.

Parameters
originalobject to be copied

References CopyData().

RRealSampleFFT::RRealSampleFFT ( const 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 SetDefaultValues().

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

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

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 SetDefaultValues().

RRealSampleFFT::RRealSampleFFT ( const u_int  n,
const double  t_lo,
const double  t_hi,
double(*)(const double, const double[])  fct_ptr,
const double  fct_par[],
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 SetDefaultValues().

RRealSampleFFT::RRealSampleFFT ( const u_int  n,
const double  t_lo,
const double  t_hi,
double(*)(const double[], const double[])  fct_ptr,
const double  fct_par[],
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 SetDefaultValues().

RRealSampleFFT::RRealSampleFFT ( const u_int  n,
const double  t_lo,
const double  t_hi,
TF1 *  fct,
bool  ctr = false 
)

Constructor with TF1 function pointer 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
fctROOT function pointer
ctrindicates how the time limits are defined (see time scale)

References InitFunctionData(), and SetDefaultValues().

RRealSampleFFT::RRealSampleFFT ( TH1 *  hptr)

Constructor from histogram (1D) pointer argument.

Parameters
hptrpointer to the histogram

References InitFunctionData(), and SetDefaultValues().

RRealSampleFFT::RRealSampleFFT ( 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 SetDefaultValues().

RRealSampleFFT::~RRealSampleFFT ( )
virtual

Destructor: cleans-up the data.

Member Function Documentation

RRealSampleFFT::GObject ( RRealSampleFFT  )
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.
RRealSampleFFT & RRealSampleFFT::operator= ( const RRealSampleFFT original)

Affectation operator.

Parameters
originalobject to be copied

References CopyData().

void RRealSampleFFT::InitFunctionData ( const u_int  n,
const double  t_lo,
const double  t_hi,
TF1 *  fct,
bool  ctr = false 
)
virtual

Set the function samples to be n real data (double) from a TF1 function defined by its pointer (see SetFunctionData). The time limits are given with respect to the defined sample

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

References SetFunctionData().

Referenced by RRealSampleFFT().

void RRealSampleFFT::InitFunctionData ( TH1 *  hptr)
virtual

Set the function samples to be n real data from a TH1 histogram defined by its pointer.

Parameters
hptrpointer to the histogram
void RRealSampleFFT::SetFunctionData ( TF1 *  fct)
virtual

Set the function samples from a function defined by its pointer.

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

Parameters
fctROOT function pointer

Referenced by InitFunctionData().

void RRealSampleFFT::FillRandomGaus ( double  mean,
double  fwhm,
TRandom *  rnd = NULL 
)
virtual

Fill the function sample with random noise with Gauss distribution. If a random generator is provided, it is used to generate the valuse (for samples reproductability), else, a local generator is used.

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

Parameters
meanmean samples value
fwhmFWHM of the samples
rndROOT random generator

Referenced by FillRandomGaus().

void RRealSampleFFT::FillRandomGaus ( double  fwhm,
TRandom *  rnd = NULL 
)
inlinevirtual

Fill the function sample with random noise with Gauss distribution with a mean value equal to 0. If a random generator is provided, it is used to generate the valuse (for samples reproductability), else, a local generator is used.

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

Parameters
fwhmFWHM of the samples
rndROOT random generator

References FillRandomGaus().

int RRealSampleFFT::SetFunctionValues ( u_int  i1,
u_int  i2,
TF1 *  fct 
)
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
fctROOT function pointer
void RRealSampleFFT::SetDefaultValues ( )
inline

Set the default titles for graphs and histograms.

References fct_title, fct_title_x, fct_title_y, fft_title, fft_title_amp, fft_title_freq, fft_title_im, fft_title_phase, fft_title_power, and fft_title_re.

Referenced by RRealSampleFFT().

void RRealSampleFFT::SetColor ( const Color_t  col)
inline

Set the default line and marker colors for graphs.

Parameters
colselected ROOT color
void RRealSampleFFT::SetFunctionTitles ( const char  tit[] = "",
const char  tit_x[] = "",
const char  tit_y[] = "" 
)
inline

Set titles for the function samples (for graphs and histograms).

Parameters
titfunction title
tit_xfunction time axis title
tit_yfunction amplitude axis title

References fct_title, fct_title_x, and fct_title_y.

void RRealSampleFFT::SetTransformTitles ( const char  tit[] = "",
const char  tit_f[] = "",
const char  tit_amp[] = "",
const char  tit_ph[] = "",
const char  tit_pow[] = "" 
)
inline

Set titles for the Fourier transform (for graphs and histograms).

Parameters
titFT title
tit_fFT frequency axis title
tit_ampFT amplitude axis title
tit_phFT phase axis title
tit_powFT power density axis title

References fft_title, fft_title_amp, fft_title_freq, fft_title_phase, and fft_title_power.

void RRealSampleFFT::SetTransformComplexTitles ( const char  tit_re[] = "",
const char  tit_im[] = "" 
)
inline

Set titles for the Fourier transform axis in complex plane (for graphs and histograms).

Parameters
tit_reFT real axis title
tit_imFT imaginary axis title

References fft_title_im, and fft_title_re.

TGraph * RRealSampleFFT::CreateFunctionGraph ( double  scale = 1.L)

Create a graph from function samples, with signal amplitude as a function of the time coordinate. The graph must be destroyed by the calling program.

Parameters
scalescaling factor

References fct_title, fct_title_x, and fct_title_y.

TGraph * RRealSampleFFT::CreateTransformGraph ( double  scale = 1.L)

Create a graph from the Fourier transform, in the complex plane. The graph dimension contains n/2-1 points, where n is the number of samples in the function. The first and last point have a null imaginary part. The graph must be destroyed by the calling program.

Parameters
scalescaling factor

References fft_title, fft_title_im, and fft_title_re.

TGraph * RRealSampleFFT::CreateTrfAmpliGraph ( double  scale = 1.L)

Create a graph from the amplitude of the Fourier transform. The graph dimension contains n/2-1 points, where n is the number of samples in the function. The first and last point have a null imaginary part. The graph must be destroyed by the calling program.

Parameters
scalescaling factor

References fft_title, fft_title_amp, and fft_title_freq.

TGraph * RRealSampleFFT::CreateTrfPhaseGraph ( bool  deg = false)

Create a graph from the phase of the Fourier transform. The graph dimension contains n/2-1 points, where n is the number of samples in the function. The first and last point have a null imaginary part. The graph must be destroyed by the calling program.

Parameters
degif the phase is expressed in deg. (default is rad.)

References fft_title, fft_title_freq, and fft_title_phase.

TGraph * RRealSampleFFT::CreatePowerGraph ( double  scale = 1.L)

Create a graph from the Power Spectral Density. The graph dimension contains n/2-1 points, where n is the number of samples in the function. The graph must be destroyed by the calling program.

Parameters
scalescaling factor

References fft_title, fft_title_freq, and fft_title_power.

TH1D * RRealSampleFFT::CreateFunctionHisto ( const string &  name,
double  scale = 1.L 
)

Create an histogram from function samples, with signal amplitude as a function of the time coordinate. The histogram must be destroyed by the calling program.

Parameters
namehistogram identifier
scalescaling factor

References fct_title_x, and fct_title_y.

int RRealSampleFFT::FillFunctionHisto ( TH1D *  hptr,
double  scale = 1.L,
bool  reset = true 
)

Fills an histogram from function samples, with signal amplitude as a function of the time coordinate. The function returns 0 if no error occured. The histogram is filled bin by bin, and no check is performed on the bound values, only on the dimension.

Parameters
hptrhistogram pointer
scalescaling factor
resetwhether the histogram should be resetted before filling
RRealSampleFFT RRealSampleFFT::operator+ ( RRealSampleFFT s)
inline

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
RRealSampleFFT RRealSampleFFT::operator- ( RRealSampleFFT s)
inline

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
RRealSampleFFT RRealSampleFFT::operator* ( RRealSampleFFT s)
inline

Product bin by bin 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
RRealSampleFFT RRealSampleFFT::operator/ ( RRealSampleFFT s)
inline

Division bin by bin 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
RRealSampleFFT RRealSampleFFT::operator* ( double  r)
inline

Scaling (multiplication factor) of the current sample.

Parameters
rscale factor
RRealSampleFFT RRealSampleFFT::operator/ ( double  r)
inline

Scaling (division factor) of the current sample.

Parameters
rscale factor
RRealSampleFFT RRealSampleFFT::operator+ ( double  r)
inline

Offset addition.

Parameters
roffset
RRealSampleFFT RRealSampleFFT::operator- ( double  r)
inline

Offset subtraction.

Parameters
roffset
void RRealSampleFFT::CopyData ( const RRealSampleFFT original)
protected

This function copies the data from original to the current object. It is not a virtual function, and it copies only the part corresponding to this stage of the heritance tree.

Parameters
originalobject to be copied

References fct_title, fct_title_x, fct_title_y, fft_title, fft_title_amp, fft_title_freq, fft_title_im, fft_title_phase, fft_title_power, and fft_title_re.

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

RRealSampleFFT::ClassDef ( RRealSampleFFT  ,
 
)

for use within ROOT.


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