GCpp general purpose C++ library  version 1.0
GBaseSampleFFT.hh
Go to the documentation of this file.
1 //======================================================================
2 /*! \file GBaseSampleFFT.hh
3  *
4  * Include file for class GBaseSampleFFT.
5  */
6 //======================================================================
7 
8 #ifndef G_BASE_SAMPLE_FFT_HH
9 #define G_BASE_SAMPLE_FFT_HH
10 
11 //======================================================================
12 
13 #include "GFFTDianelsonLanczos.hh"
14 #include "GLogMessage.hh"
15 
16 //! Samples data type identifier
17 #define GFFT_DATA_UNDEFINED 0
18 //! Samples data type identifier
19 #define GFFT_DATA_REAL 1
20 //! Samples data type identifier
21 #define GFFT_DATA_COMPLEX 2
22 
23 //! File header identifier
24 #define GFFT_HEADER_ID "GFFT_SAMPLE "
25 
26 //! I/O error code
27 #define GFFT_IOERR_FOPEN 1
28 //! I/O error code
29 #define GFFT_IOERR_HEADER 2
30 //! I/O error code
31 #define GFFT_IOERR_HEADER_DATA 3
32 //! I/O error code
33 #define GFFT_IOERR_DATA_TYPE 4
34 //! I/O error code
35 #define GFFT_IOERR_DATA 5
36 
37 //----------------------------------------------------------------------
38 /*! \class GBaseSampleFFT
39  *
40  * This class defines a base class to manage data samples and the
41  * corresponding Fourier transform.
42  *
43  * The data are processed as real numbers with the FFT functions
44  * using the Danielson-Lanczos algorithm (see \ref FFT_Dianelson_Lanczos).
45  * Thus, the arrays (sample and Fourier transform) sizes depend whether
46  * the function is a real or a complex function.
47  *
48  * The FFT algorithm works with dimensions N=2^k (power of two), so if the
49  * data dimension is not a power of 2, it is padded with zeros.
50  * It is anyway suggested to use samples with a dimension that is a
51  * power of 2.
52  *
53  * In order to reduce computation time, the Fourier transform and/or the
54  * function samples are only computed when required.
55  * The status for the 2 data arrays (function samples and Fourier transform
56  * data) is defined by 2 boolean variables: \b data_set and \b fft_set.
57  *
58  * \ref GFFT_time_scale_mode
59  */
61 {
62  //------------------------------------------------------------
63  /*! \object_doc */
65  //------------------------------------------------------------
66 
67  protected:
68  static double epsilon; ///< Limit for 0 checking
69  static string data_type_name[3]; ///< Names of the data types
70 
71  int data_type; ///< 0=undefined; 1=real; 2=complex
72  u_int samples_number; ///< Number of data samples
73  u_int samples_dim; ///< Effective arrays dimension (power of 2)
74  u_int alloc_dim; ///< Effective allocated dim
75 
76  double * data_array; ///< Array of function samples
77  double * fft_array; ///< Array of Fourier transform data
78 
79  bool data_set; ///< Indicate if the data samples are set
80  bool fft_set; ///< Indicate if the FFT data are set
81 
82  // sampling rate data
83  double time_min; ///< Time of the first sample
84  double time_max; ///< Time of the last sample (after padding with zeros)
85  double time_end; ///< Time of the last (defined) sample (before padding)
86  double time_step; ///< Time step
87  double freq_step; ///< Frequency step
88  double freq_crit; ///< Frequency maximum value (Nyquist critical frequency)
89 
90  public:
91  // Constructors
92  GBaseSampleFFT ( u_int n = 0, double dt = 1.L );
93  GBaseSampleFFT ( const GBaseSampleFFT & original );
94 
95  // Affectation
96  GBaseSampleFFT & operator = ( const GBaseSampleFFT & original );
97 
98  // Destructors
99  virtual ~GBaseSampleFFT ( );
100 
101  //----------------------------------------------------------
102  // Initialisation and setting functions
103  virtual void ResetData ( );
104  virtual void ClearData ( );
105  virtual u_int InitData ( u_int n, double dt = 1.L, bool clear = true );
106  virtual int ResizeData ( u_int n );
107 
108  /*! Pure virtual function that must be defined in derived classes
109  * to fill the function samples data from an array.
110  * \param n number of data (expressed in real)
111  * \param data data array
112  */
113  virtual void InitFunctionData ( u_int n, const double data [] ) = 0;
114 
115  /*! Pure virtual function that must be defined in derived classes
116  * to fill the Fourier transform data from an array of real numbers.
117  * The array must be sorted in the following order:
118  * Re[0], Im[0], Re[1], Im[1], ... Re[n/2-1], Im[n/2-1]
119  * (see \ref GFFTComplex, \ref GFFTReal, ...).
120  * \param n number of data (expressed in real)
121  * \param data data array
122  */
123  virtual void InitTransformData ( u_int n, const double data [] ) = 0;
124 
125  /*! Pure virtual function that must be defined in derived classes
126  * to fill the Fourier transform data from 2 arrays of real numbers:
127  * 1 for the real terms of the Fourier transform, 1 for the imaginary
128  * terms.
129  * \param n number of complex numbers data
130  * \param data_re array for real part of data
131  * \param data_im array for imaginary part of data
132  */
133  virtual void InitTransformData ( u_int n, const double data_re [],
134  const double data_im [] ) = 0;
135 
136 
137  /*! Pure virtual function to set a single value in the sample.
138  * The index meaning is different for a real or a complex data sample.
139  *
140  * The use of this function implies that the full sample is set using this
141  * function, and that the FFT data have to be recomputed.
142  *
143  * \param i data index
144  * \param val data (real) value
145  */
146  virtual void SetFunctionValue ( u_int i, double val ) = 0;
147 
148  /*! Set a single (complex) value for the FFT
149  * (see derived class descriptions for SetFFTValue and SetTransformValue).
150  * \param i index of the FFT value
151  * \param re real part of the data
152  * \param im imaginary part of the data
153  */
154  virtual void SetTransformValue ( u_int i, double re, double im = 0.L ) = 0;
155 
156  /*! Set a single (complex) value for the FFT
157  * (see derived class descriptions for SetFFTValue and SetTransformValue).
158  * \param i index of the FFT value
159  * \param re real part of the data
160  * \param im imaginary part of the data
161  */
162  virtual void SetFFTValue ( u_int i, double re, double im = 0.L ) = 0;
163 
164  /*! Compute the Fourier transform from the function samples that are
165  * currently stored in the object.
166  * The function returns 0 if no error occured (dimensions and arrays
167  * are defined, ...).
168  */
169  virtual int ComputeTransformData ( ) = 0;
170 
171  /*! Compute the inverse Fourier transform (function samples) from the
172  * Fourier transform data that are currently stored in the object.
173  * The function returns 0 if no error occured (dimensions and arrays
174  * are defined, ...).
175  */
176  virtual int ComputeFunctionData ( ) = 0;
177 
178  /*! Normalisation of the samples sum.
179  * This function is virtual because the normalisation procedure
180  * is different for real or complex samples).
181  * \param s requested normalisation
182  */
183  virtual void Normalise ( double s = 1.L ) = 0;
184 
185  virtual void Scale ( double s );
186 
187  // total power
188  virtual double GetSumSquared ( bool f = false );
189 
190  //----------------------------------------------------------
191  // Retrieve sample information
192  bool IsDefined ( ) const; // inline
193  bool IsFunctionSet ( ) const; // inline
194  bool IsTransformSet ( ) const; // inline
195 
196  void FunctionSet ( bool b = true ); // inline
197  void TransformSet ( bool b = true ); // inline
198 
199  int GetDataType ( ) const; // inline
200  u_int GetSamplesNumber ( ) const; // inline
201  u_int GetDimension ( ) const; // inline
202 
203  bool HasSameSampling ( const GBaseSampleFFT & sample, double prec = 0.0001 ) const;
204 
205  // time sampling functions
206  virtual int SetTimeRange ( double t_lo, double t_up, bool ctr = false );
207 
208  virtual double GetTimeMin ( ) const;
209  virtual double GetTimeMax ( ) const;
210  virtual double GetTimeEnd ( ) const;
211  virtual double GetTimeStep ( ) const;
212 
213  virtual double GetCriticalFrequency ( ) const;
214  virtual double GetFrequencyMax ( ) const;
215  virtual double GetFrequencyStep ( ) const;
216 
217  // Get the data
218  double * GetFunctionData ( bool check = true );
219  double * GetTransformData ( bool check = true );
220 
221  // Dupplicate the arrays
222  double * CopyFunctionData ( ) const;
223  double * CopyTransformData ( ) const;
224 
225  double * CreatePowerSpectralDensity ( );
226  double * CreatePSD ( );
227 
228  //----------------------------------------------------------
229  // Write/read functions
230  virtual int Write ( const string & fname );
231  virtual int Write ( FILE * fp );
232 
233  virtual int Read ( const string & fname );
234  virtual int Read ( FILE * fp );
235 
236  //----------------------------------------------------------
237  protected:
238  // Functions for internal use
239  virtual int WriteHeader ( FILE * fp );
240  virtual int ReadHeader ( FILE * fp );
241 
242  void CopyData ( const GBaseSampleFFT & original );
243 
244  public:
245 };
246 
247 //======================================================================
248 // Inline functions
249 #include "icc/GBaseSampleFFT.icc"
250 
251 
252 //======================================================================
253 #endif
virtual void InitTransformData(u_int n, const double data[])=0
GBaseSampleFFT(u_int n=0, double dt=1.L)
Definition: GBaseSampleFFT.cpp:25
static string data_type_name[3]
Names of the data types.
Definition: GBaseSampleFFT.hh:69
virtual int Write(const string &fname)
Definition: GBaseSampleFFT.cpp:619
double freq_crit
Frequency maximum value (Nyquist critical frequency)
Definition: GBaseSampleFFT.hh:88
bool fft_set
Indicate if the FFT data are set.
Definition: GBaseSampleFFT.hh:80
virtual double GetTimeEnd() const
Definition: GBaseSampleFFT.icc:62
int data_type
0=undefined; 1=real; 2=complex
Definition: GBaseSampleFFT.hh:71
double time_step
Time step.
Definition: GBaseSampleFFT.hh:86
double * CreatePowerSpectralDensity()
Definition: GBaseSampleFFT.cpp:579
int GetDataType() const
Definition: GBaseSampleFFT.icc:34
virtual int ResizeData(u_int n)
Definition: GBaseSampleFFT.cpp:229
u_int GetSamplesNumber() const
Definition: GBaseSampleFFT.icc:40
double * fft_array
Array of Fourier transform data.
Definition: GBaseSampleFFT.hh:77
bool HasSameSampling(const GBaseSampleFFT &sample, double prec=0.0001) const
Definition: GBaseSampleFFT.cpp:513
double * CreatePSD()
Definition: GBaseSampleFFT.icc:119
u_int samples_dim
Effective arrays dimension (power of 2)
Definition: GBaseSampleFFT.hh:73
virtual void SetFFTValue(u_int i, double re, double im=0.L)=0
virtual double GetTimeStep() const
Definition: GBaseSampleFFT.icc:66
virtual ~GBaseSampleFFT()
Definition: GBaseSampleFFT.cpp:81
virtual double GetSumSquared(bool f=false)
Definition: GBaseSampleFFT.cpp:473
u_int alloc_dim
Effective allocated dim.
Definition: GBaseSampleFFT.hh:74
virtual void Scale(double s)
Definition: GBaseSampleFFT.cpp:458
double * data_array
Array of function samples.
Definition: GBaseSampleFFT.hh:76
virtual double GetTimeMin() const
Definition: GBaseSampleFFT.icc:52
void TransformSet(bool b=true)
Definition: GBaseSampleFFT.icc:30
void CopyData(const GBaseSampleFFT &original)
Definition: GBaseSampleFFT.cpp:353
virtual int SetTimeRange(double t_lo, double t_up, bool ctr=false)
Definition: GBaseSampleFFT.cpp:413
GObjectV(GBaseSampleFFT)
double time_min
Time of the first sample.
Definition: GBaseSampleFFT.hh:83
virtual double GetFrequencyStep() const
Definition: GBaseSampleFFT.icc:80
virtual void SetFunctionValue(u_int i, double val)=0
u_int GetDimension() const
Definition: GBaseSampleFFT.icc:46
virtual int ComputeFunctionData()=0
virtual void SetTransformValue(u_int i, double re, double im=0.L)=0
double time_end
Time of the last (defined) sample (before padding)
Definition: GBaseSampleFFT.hh:85
virtual double GetFrequencyMax() const
Definition: GBaseSampleFFT.icc:76
virtual void InitFunctionData(u_int n, const double data[])=0
GBaseSampleFFT & operator=(const GBaseSampleFFT &original)
Definition: GBaseSampleFFT.cpp:67
virtual int WriteHeader(FILE *fp)
Definition: GBaseSampleFFT.cpp:751
u_int samples_number
Number of data samples.
Definition: GBaseSampleFFT.hh:72
virtual double GetTimeMax() const
Definition: GBaseSampleFFT.icc:58
void FunctionSet(bool b=true)
Definition: GBaseSampleFFT.icc:24
double * CopyTransformData() const
Definition: GBaseSampleFFT.cpp:562
virtual int ComputeTransformData()=0
virtual double GetCriticalFrequency() const
Definition: GBaseSampleFFT.icc:72
double * GetTransformData(bool check=true)
Definition: GBaseSampleFFT.icc:110
double time_max
Time of the last sample (after padding with zeros)
Definition: GBaseSampleFFT.hh:84
virtual void ResetData()
Definition: GBaseSampleFFT.cpp:95
bool IsDefined() const
Definition: GBaseSampleFFT.icc:10
bool data_set
Indicate if the data samples are set.
Definition: GBaseSampleFFT.hh:79
virtual int Read(const string &fname)
Definition: GBaseSampleFFT.cpp:645
virtual int ReadHeader(FILE *fp)
Definition: GBaseSampleFFT.cpp:783
double freq_step
Frequency step.
Definition: GBaseSampleFFT.hh:87
unsigned int u_int
Definition: GTypes.hh:38
bool IsTransformSet() const
Definition: GBaseSampleFFT.icc:18
double * GetFunctionData(bool check=true)
Definition: GBaseSampleFFT.icc:96
virtual void ClearData()
Definition: GBaseSampleFFT.cpp:117
virtual u_int InitData(u_int n, double dt=1.L, bool clear=true)
Definition: GBaseSampleFFT.cpp:151
double * CopyFunctionData() const
Definition: GBaseSampleFFT.cpp:545
static double epsilon
Limit for 0 checking.
Definition: GBaseSampleFFT.hh:68
bool IsFunctionSet() const
Definition: GBaseSampleFFT.icc:14
Definition: GBaseSampleFFT.hh:60
virtual void Normalise(double s=1.L)=0