v2.0.0
Loading...
Searching...
No Matches
spectral.h
Go to the documentation of this file.
1//=============================================================================================================
39
40#ifndef SPECTRAL_H
41#define SPECTRAL_H
42
43//=============================================================================================================
44// INCLUDES
45//=============================================================================================================
46
47#include "utils_global.h"
48
49#include <vector>
50#include <utility>
51
52//=============================================================================================================
53// EIGEN INCLUDES
54//=============================================================================================================
55
56#include <Eigen/Core>
57
58//=============================================================================================================
59// QT INCLUDES
60//=============================================================================================================
61
62#include <QString>
63#include <QPair>
64#include <QSharedPointer>
65
66//=============================================================================================================
67// DEFINE NAMESPACE UTILSLIB
68//=============================================================================================================
69
70namespace UTILSLIB
71{
72
77 Eigen::RowVectorXd vecData;
78 Eigen::MatrixXd matTaper;
79 int iNfft;
80};
81
82//=============================================================================================================
89{
90
91public:
92 //=========================================================================================================
96 Spectral() = delete;
97
98 //=========================================================================================================
108 static Eigen::MatrixXcd computeTaperedSpectraRow(const Eigen::RowVectorXd &vecData,
109 const Eigen::MatrixXd &matTaper,
110 int iNfft);
111
112 //=========================================================================================================
123 static QVector<Eigen::MatrixXcd> computeTaperedSpectraMatrix(const Eigen::MatrixXd &matData,
124 const Eigen::MatrixXd &matTaper,
125 int iNfft,
126 bool bUseThreads = true);
127
128// //=========================================================================================================
129// /**
130// * Calculates the full tapered spectra of a given input matrix data. This function calculates each row in parallel.
131// *
132// * @param[in] matData input matrix data (time domain), for which the spectrum is computed.
133// * @param[in] matTaper tapers used to compute the spectra.
134// * @param[in] iNfft FFT length.
135// * @param[in] bUseThreads Whether to use multiple threads.
136// *
137// * @return tapered spectra of the input data.
138// */
139// static std::vector<Eigen::MatrixXcd> computeTaperedSpectraMatrix(const Eigen::MatrixXd &matData,
140// const Eigen::MatrixXd &matTaper,
141// int iNfft,
142// bool bUseThreads = true);
143
144 //=========================================================================================================
152 static Eigen::MatrixXcd compute(const TaperedSpectraInputData& inputData);
153
154 //=========================================================================================================
161 static void reduce(QVector<Eigen::MatrixXcd>& finalData,
162 const Eigen::MatrixXcd& resultData);
163
164// //=========================================================================================================
165// /**
166// * Reduces the taperedSpectra results to a final result. This function gets called in parallel.
167// *
168// * @param[out] finalData The final data data.
169// * @param[in] resultData The resulting data from the computation step.
170// */
171// static void reduce(std::vector<Eigen::MatrixXcd>& finalData,
172// const Eigen::MatrixXcd& resultData);
173
174 //=========================================================================================================
185 static Eigen::RowVectorXd psdFromTaperedSpectra(const Eigen::MatrixXcd &matTapSpectrum,
186 const Eigen::VectorXd &vecTapWeights,
187 int iNfft,
188 double dSampFreq=1.0);
189
190 //=========================================================================================================
203 static Eigen::RowVectorXcd csdFromTaperedSpectra(const Eigen::MatrixXcd &vecTapSpectrumSeed,
204 const Eigen::MatrixXcd &vecTapSpectrumTarget,
205 const Eigen::VectorXd &vecTapWeightsSeed,
206 const Eigen::VectorXd &vecTapWeightsTarget,
207 int iNfft,
208 double dSampFreq = 1.0);
209
210 //=========================================================================================================
219 static Eigen::VectorXd calculateFFTFreqs(int iNfft, double dSampFreq);
220
221 //=========================================================================================================
230 static QPair<Eigen::MatrixXd, Eigen::VectorXd> generateTapers(int iSignalLength,
231 const QString &sWindowType = "hanning");
232
233 //=========================================================================================================
242 static std::pair<Eigen::MatrixXd, Eigen::VectorXd> generateTapers(int iSignalLength,
243 const std::string &sWindowType = "hanning");
244
245private:
246 //=========================================================================================================
254 static Eigen::MatrixXd hanningWindow(int iSignalLength);
255};
256
257//=============================================================================================================
258// INLINE DEFINITIONS
259//=============================================================================================================
260}//namespace
261
262#endif // SPECTRAL_H
utils library export/import macros.
#define UTILSSHARED_EXPORT
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Definition buildinfo.h:45
Input parameters for multi-taper spectral estimation (data matrix, taper weights, frequencies,...
Definition spectral.h:76
Eigen::RowVectorXd vecData
Definition spectral.h:77
static Eigen::RowVectorXcd csdFromTaperedSpectra(const Eigen::MatrixXcd &vecTapSpectrumSeed, const Eigen::MatrixXcd &vecTapSpectrumTarget, const Eigen::VectorXd &vecTapWeightsSeed, const Eigen::VectorXd &vecTapWeightsTarget, int iNfft, double dSampFreq=1.0)
Definition spectral.cpp:224
static QVector< Eigen::MatrixXcd > computeTaperedSpectraMatrix(const Eigen::MatrixXd &matData, const Eigen::MatrixXd &matTaper, int iNfft, bool bUseThreads=true)
Definition spectral.cpp:100
static Eigen::RowVectorXd psdFromTaperedSpectra(const Eigen::MatrixXcd &matTapSpectrum, const Eigen::VectorXd &vecTapWeights, int iNfft, double dSampFreq=1.0)
Definition spectral.cpp:198
static QPair< Eigen::MatrixXd, Eigen::VectorXd > generateTapers(int iSignalLength, const QString &sWindowType="hanning")
Definition spectral.cpp:292
static Eigen::MatrixXcd compute(const TaperedSpectraInputData &inputData)
Definition spectral.cpp:171
static Eigen::MatrixXcd computeTaperedSpectraRow(const Eigen::RowVectorXd &vecData, const Eigen::MatrixXd &matTaper, int iNfft)
Definition spectral.cpp:72
static Eigen::VectorXd calculateFFTFreqs(int iNfft, double dSampFreq)
Definition spectral.cpp:278
static void reduce(QVector< Eigen::MatrixXcd > &finalData, const Eigen::MatrixXcd &resultData)
Definition spectral.cpp:181