GPLIB++
TimeFrequency.h
Go to the documentation of this file.
1 #ifndef TIMEFREQUENCY_H_
2 #define TIMEFREQUENCY_H_
3 #include <complex>
4 #include "TsSpectrum.h"
5 #include "statutils.h"
6 #include "VecMat.h"
7 #include "FatalException.h"
8 
9 namespace gplib
10  {
11 
12  /** \addtogroup tstools Time series analysis methods */
13  /* @{ */
14 
15  //! Calculate a sliding windowed fourier transform for a time series and store the results for each segment in a matrix
16  /*! This function returns the time-frequency matrix for a real time series
17  * \param tsbegin iterator to the start of the input time series
18  * \param tsend iterator to the end of the input time series
19  * \param length of segments used for fourier transform in points, if the last segment is longer than the time-series it is discarded
20  * \param Function object for windowing the segments before fourier transform \see WFunc.h for a selection of provided function objects
21  */
22  template<typename InputIterator, typename WindowFunctype>
23  gplib::cmat TimeFrequency(InputIterator tsbegin, InputIterator tsend,
24  const size_t seglength, WindowFunctype WFunc)
25  {
26  const size_t insize = distance(tsbegin, tsend);
27  if (insize < seglength)
28  throw FatalException(
29  "Time Series shorter than requested segment length");
30  const unsigned int nsegments = insize / seglength; //we do integer division, so the last (not full) segment will be discarded
31  const unsigned int speclength = seglength / 2 + 1;
32  TsSpectrum SpecEst(false);
33  gplib::cmat output(nsegments, speclength);
34  typedef std::vector<
35  typename std::iterator_traits<InputIterator>::value_type> tinvector;
36 
37  InputIterator currstart = tsbegin; //the start of the current segment
38  InputIterator currend = tsbegin;
39  advance(currend, seglength);
40 
41  std::vector<std::complex<double> > segspec(speclength);
42  for (unsigned int i = 0; i < nsegments; ++i)
43  {
44  tinvector currseg(currstart, currend); //temporary copy of time series
45  SubMean(currseg.begin(), currseg.end());
46  ApplyWindow(currseg.begin(),currseg.end(),currseg.begin(),WFunc);
47  SpecEst.CalcSpectrum(currseg.begin(), currseg.end(),
48  segspec.begin(), segspec.end());
49  copy(segspec.begin(), segspec.end(), row(output, i).begin());
50  advance(currstart, seglength); //shift the boundaries of current segment
51  advance(currend, seglength);
52  }
53  return output;
54  }
55  /* @} */
56  }
57 #endif /*TIMEFREQUENCY_H_*/
void CalcSpectrum(_InputIterator tsbegin, _InputIterator tsend, _OutputIterator freqbegin, _OutputIterator freqend)
The member function to calculate a spectrum from a time series.
Definition: TsSpectrum.h:80
void SubMean(InputIterator begin, InputIterator end, typename std::iterator_traits< InputIterator >::value_type mean)
Substract the mean from each element in the container, mean is passed as a parameter.
Definition: statutils.h:85
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 segme...
Definition: TimeFrequency.h:23
void ApplyWindow(InputIterator inbegin, InputIterator inend, OutputIterator outbegin, WindowFunctype WFunc, double relshift=0.0)
Apply one of the above window functions to a range.
Definition: WFunc.h:104
The class CTsSpectrum is used to calculate spectra from (real) time series data.
Definition: TsSpectrum.h:21
The basic exception class for all errors that arise in gplib.