Time series analysis methods


Classes

class  SimpleLp
 A simple low pass. More...
class  TimeSeriesComponent
 TimeSeriesComponent is the base storage class for all types of time series data. More...
class  TsSpectrum
 The class CTsSpectrum is used to calculate spectra from (real) time series data. More...
struct  Hamming
 This functor returns the weighting factor for the Hamming window, given the relative position relpos [0..1] in the time series. More...
struct  Hanning
 This functor returns the weighting factor for the Hanning window, given the relative position (0..1) in the time series. More...
struct  Boxcar
 A functor for the simple Boxcar function. More...
struct  Steep
 This functor rises steeply at the edges and then has a wide range where it is unity. More...
struct  CosSq
 The cosine squared windows of fixed width. More...
class  TruncCosSq
 A vraible width cosine squared window that is zero outside. More...

Functions

double FreqToW (const double f)
 Transform from frequency domain to w-domain.
 SimpleLp::SimpleLp (const double cornerfreq)
 The constructor takes the dimensionless corner frequency, i.e. the corner frequency in Hz times the sampling rate.
double SimpleLp::operator() (const double currvalue)
 This version of the operator is suitable for use with std::transform, it returns a filtered value for each call with the current value.
template<typename _InputIterator, typename _OutputIterator>
void ShortCorr (_InputIterator masterbegin, _InputIterator masterend, _InputIterator shortbegin, _InputIterator shortend, _OutputIterator outbegin)
 Calculate the correlation between a short time series and a master time series.
template<typename InputIterator, typename OutputIterator, typename WindowFunctype>
void StackedSpectrum (InputIterator tsbegin, InputIterator tsend, OutputIterator freqbegin, const size_t seglength, WindowFunctype WFunc)
 This template is used to calculate stacked spectra for example for power spectrum estimation.
template<typename InputIterator, typename WindowFunctype>
gplib::cmat TimeFrequency (InputIterator tsbegin, InputIterator tsend, const size_t seglength, WindowFunctype WFunc)
 Calculate a sliding windowed fourier transform for a time series and store the results for each segment in a matrix.
template<typename InputIterator, typename OutputIterator, typename WindowFunctype>
void ApplyWindow (InputIterator inbegin, InputIterator inend, OutputIterator outbegin, WindowFunctype WFunc)
 Apply one of the above window functions to a range.
gplib::rvec TimeSeriesComponent::GetUblasData ()
 For some methods we prefer to get the data as a ublas vector, we return a copy.
template<typename _InputIterator, typename _OutputIterator>
void TsSpectrum::CalcSpectrum (_InputIterator tsbegin, _InputIterator tsend, _OutputIterator freqbegin, _OutputIterator freqend)
 The member function to calculate a spectrum from a time series.
template<typename _InputIterator, typename _OutputIterator>
void TsSpectrum::CalcTimeSeries (_InputIterator freqbegin, _InputIterator freqend, _OutputIterator tsbegin, _OutputIterator tsend)
 The member function to calculate a time series from (complex) spectra.

Detailed Description


Function Documentation

template<typename InputIterator, typename OutputIterator, typename WindowFunctype>
void ApplyWindow ( InputIterator  inbegin,
InputIterator  inend,
OutputIterator  outbegin,
WindowFunctype  WFunc 
) [inline]

Apply one of the above window functions to a range.

Apply one of the above window functions to a range given by the iterators inbegin and outbegin and write the result into a structure through the iterator outbegin, make sure the datastructure can hold enough elements. Input and output structures can be the same.

Definition at line 90 of file WFunc.h.

References length.

Referenced by main().

template<typename _InputIterator, typename _OutputIterator>
void TsSpectrum::CalcSpectrum ( _InputIterator  tsbegin,
_InputIterator  tsend,
_OutputIterator  freqbegin,
_OutputIterator  freqend 
) [inline, inherited]

The member function to calculate a spectrum from a time series.

Calc a spectrum from the input range given by tsbegin and tsend and write to output range given by freqbegin,freqend.

Takes three iterators for input: beginning of input, end of input and beginning of output Make sure the structure pointed to by freqbegin can hold at least (tsend -tsbegin)/2 +1 elements

Definition at line 67 of file TsSpectrum.h.

References size.

Referenced by main(), StackedSpectrum(), and TimeFrequency().

template<typename _InputIterator, typename _OutputIterator>
void TsSpectrum::CalcTimeSeries ( _InputIterator  freqbegin,
_InputIterator  freqend,
_OutputIterator  tsbegin,
_OutputIterator  tsend 
) [inline, inherited]

The member function to calculate a time series from (complex) spectra.

Calc a times series from the input range given by freqbegin and freqend and write to output range given by tsbegin,tsend.

Takes three iterators for input: beginning of input, end of input and beginning of output Make sure the structure pointed to by tsbegin can hold at least (freqend -freqbegin)*2 - 2 elements

Definition at line 83 of file TsSpectrum.h.

References size.

Referenced by MultiRecCalc::CalcRecData(), and main().

double FreqToW ( const double  f  ) 

Transform from frequency domain to w-domain.

Definition at line 17 of file FilterFunc.h.

Referenced by SimpleLp::SimpleLp().

gplib::rvec TimeSeriesComponent::GetUblasData (  )  [inline, inherited]

For some methods we prefer to get the data as a ublas vector, we return a copy.

We want to make this inline, so it appears in the header.

This function is quite expensive, as the ublas object is generated on the fly and we return a copy.

Definition at line 107 of file TimeSeriesComponent.h.

double SimpleLp::operator() ( const double  currvalue  )  [inline, inherited]

This version of the operator is suitable for use with std::transform, it returns a filtered value for each call with the current value.

Definition at line 48 of file FilterFunc.h.

template<typename _InputIterator, typename _OutputIterator>
void ShortCorr ( _InputIterator  masterbegin,
_InputIterator  masterend,
_InputIterator  shortbegin,
_InputIterator  shortend,
_OutputIterator  outbegin 
) [inline]

Calculate the correlation between a short time series and a master time series.

The short time series will be shifted along the master and the corresponding correlation will be stored in the structure referenced by outpebegin, make sure that structure can hold enough data

Definition at line 15 of file ShortCorr.h.

References SubMean().

Referenced by CalcCorrAndWrite().

SimpleLp::SimpleLp ( const double  cornerfreq  )  [inline, inherited]

The constructor takes the dimensionless corner frequency, i.e. the corner frequency in Hz times the sampling rate.

Definition at line 36 of file FilterFunc.h.

References FreqToW().

template<typename InputIterator, typename OutputIterator, typename WindowFunctype>
void StackedSpectrum ( InputIterator  tsbegin,
InputIterator  tsend,
OutputIterator  freqbegin,
const size_t  seglength,
WindowFunctype  WFunc 
) [inline]

This template is used to calculate stacked spectra for example for power spectrum estimation.

The template StackedSpectrum calculates stacked spectra from an input time series. Two input iterators define the input range, one output iterator the start of the output range, we have to provide the length of individual segments and a windowing function. Make sure the container the output is stored in can hold at least seglength/2+1. The Windowing function gets the relative position within the segment and should therefore be defined on the interval 0..1

Definition at line 21 of file StackedSpectrum.h.

References TsSpectrum::CalcSpectrum(), and SubMean().

Referenced by CalcPSpecAndWrite(), and main().

template<typename InputIterator, typename WindowFunctype>
gplib::cmat TimeFrequency ( InputIterator  tsbegin,
InputIterator  tsend,
const size_t  seglength,
WindowFunctype  WFunc 
) [inline]

Calculate a sliding windowed fourier transform for a time series and store the results for each segment in a matrix.

This function returns the time-frequency matrix for a real time series

Parameters:
tsbegin iterator to the start of the input time series
tsend iterator to the end of the input time series
length of segments used for fourier transform in points, if the last segment is longer than the time-series it is discarded
Function object for windowing the segments before fourier transform
See also:
WFunc.h for a selection of provided function objects

Definition at line 19 of file TimeFrequency.h.

References TsSpectrum::CalcSpectrum(), and SubMean().

Referenced by MultiRecCalc::CalcRecData(), CalcTfAndWrite(), and main().


Generated on Fri Jul 4 15:30:21 2008 for GPLIB++ by  doxygen 1.5.5