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 "math_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 //=========================================================================================================
136 static Eigen::MatrixXcd compute(const TaperedSpectraInputData& inputData);
137
138 //=========================================================================================================
145 static void reduce(QVector<Eigen::MatrixXcd>& finalData,
146 const Eigen::MatrixXcd& resultData);
147
148 //=========================================================================================================
159 static Eigen::RowVectorXd psdFromTaperedSpectra(const Eigen::MatrixXcd &matTapSpectrum,
160 const Eigen::VectorXd &vecTapWeights,
161 int iNfft,
162 double dSampFreq=1.0);
163
164 //=========================================================================================================
177 static Eigen::RowVectorXcd csdFromTaperedSpectra(const Eigen::MatrixXcd &vecTapSpectrumSeed,
178 const Eigen::MatrixXcd &vecTapSpectrumTarget,
179 const Eigen::VectorXd &vecTapWeightsSeed,
180 const Eigen::VectorXd &vecTapWeightsTarget,
181 int iNfft,
182 double dSampFreq = 1.0);
183
184 //=========================================================================================================
193 static Eigen::VectorXd calculateFFTFreqs(int iNfft, double dSampFreq);
194
195 //=========================================================================================================
204 static QPair<Eigen::MatrixXd, Eigen::VectorXd> generateTapers(int iSignalLength,
205 const QString &sWindowType = "hanning");
206
207 //=========================================================================================================
216 static std::pair<Eigen::MatrixXd, Eigen::VectorXd> generateTapers(int iSignalLength,
217 const std::string &sWindowType = "hanning");
218
219private:
220 //=========================================================================================================
228 static Eigen::MatrixXd hanningWindow(int iSignalLength);
229};
230
231//=============================================================================================================
232// INLINE DEFINITIONS
233//=============================================================================================================
234}//namespace
235
236#endif // SPECTRAL_H
math library export/import macros.
#define MATHSHARED_EXPORT
Definition math_global.h:55
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
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:216
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:190
static QPair< Eigen::MatrixXd, Eigen::VectorXd > generateTapers(int iSignalLength, const QString &sWindowType="hanning")
Definition spectral.cpp:284
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:270
static void reduce(QVector< Eigen::MatrixXcd > &finalData, const Eigen::MatrixXcd &resultData)
Definition spectral.cpp:181