v2.0.0
Loading...
Searching...
No Matches
multitaper_psd.cpp
Go to the documentation of this file.
1//=============================================================================================================
33
34//=============================================================================================================
35// INCLUDES
36//=============================================================================================================
37
38#include "multitaper_psd.h"
39#include "dpss.h"
40
41//=============================================================================================================
42// EIGEN INCLUDES
43//=============================================================================================================
44
45#include <Eigen/Core>
46//#ifndef EIGEN_FFTW_DEFAULT
47//#define EIGEN_FFTW_DEFAULT
48//#endif
49#include <unsupported/Eigen/FFT>
50
51//=============================================================================================================
52// STD INCLUDES
53//=============================================================================================================
54
55#include <cmath>
56
57//=============================================================================================================
58// USED NAMESPACES
59//=============================================================================================================
60
61using namespace UTILSLIB;
62using namespace Eigen;
63
64//=============================================================================================================
65// DEFINE MEMBER METHODS
66//=============================================================================================================
67
69 double sfreq,
70 double halfBandwidth,
71 int nTapers,
72 int nFft)
73{
74 const int nChannels = static_cast<int>(matData.rows());
75 const int nTimes = static_cast<int>(matData.cols());
76
77 if (nFft < 0)
78 nFft = nTimes;
79
80 const int nFreqs = nFft / 2 + 1;
81
82 // 1. Compute DPSS tapers
83 DpssResult dpss = Dpss::compute(nTimes, halfBandwidth, nTapers);
84 const int nTap = static_cast<int>(dpss.matTapers.rows());
85
86 // Compute eigenvalue weights for weighted averaging
87 double weightSum = 0.0;
88 for (int t = 0; t < nTap; ++t)
89 weightSum += dpss.vecEigenvalues[t];
90
92 result.matPsd = MatrixXd::Zero(nChannels, nFreqs);
93
94 // Build frequency axis
95 result.vecFreqs.resize(nFreqs);
96 for (int k = 0; k < nFreqs; ++k)
97 result.vecFreqs[k] = static_cast<double>(k) * sfreq / static_cast<double>(nFft);
98
99 Eigen::FFT<double> fft;
100
101 // 2. For each channel, apply each taper, FFT, accumulate weighted |X|^2
102 for (int ch = 0; ch < nChannels; ++ch) {
103 RowVectorXd psd = RowVectorXd::Zero(nFreqs);
104
105 for (int t = 0; t < nTap; ++t) {
106 // Apply taper element-wise
107 VectorXd tapered = matData.row(ch).transpose().array()
108 * dpss.matTapers.row(t).transpose().array();
109
110 // Zero-pad to nFft if needed
111 VectorXd padded = VectorXd::Zero(nFft);
112 const int copyLen = std::min(nTimes, nFft);
113 padded.head(copyLen) = tapered.head(copyLen);
114
115 // Forward FFT
116 VectorXcd spec;
117 fft.fwd(spec, padded);
118
119 // Accumulate weighted |X[k]|^2
120 const double w = dpss.vecEigenvalues[t];
121 for (int k = 0; k < nFreqs; ++k)
122 psd[k] += w * std::norm(spec[k]);
123 }
124
125 // Normalise: divide by (weightSum × sfreq)
126 const double dNorm = 1.0 / (weightSum * sfreq);
127 psd *= dNorm;
128
129 // One-sided: double all non-DC, non-Nyquist bins
130 for (int k = 1; k < nFreqs - 1; ++k)
131 psd[k] *= 2.0;
132
133 result.matPsd.row(ch) = psd;
134 }
135
136 return result;
137}
MultitaperPsd class declaration — multitaper power spectral density estimator.
Dpss class declaration — Discrete Prolate Spheroidal Sequences (Slepian tapers).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Result of a DPSS taper computation.
Definition dpss.h:62
Eigen::MatrixXd matTapers
nTapers × N, each row is a unit-norm taper
Definition dpss.h:63
Eigen::VectorXd vecEigenvalues
Concentration ratios, length nTapers.
Definition dpss.h:64
static DpssResult compute(int N, double halfBandwidth, int nTapers=-1)
Definition dpss.cpp:69
Result of a multitaper PSD computation.
Eigen::RowVectorXd vecFreqs
Frequency axis in Hz, length nFft/2+1.
Eigen::MatrixXd matPsd
n_channels × n_freqs; one-sided PSD in unit²/Hz
static MultitaperPsdResult compute(const Eigen::MatrixXd &matData, double sfreq, double halfBandwidth=4.0, int nTapers=-1, int nFft=-1)