v2.0.0
Loading...
Searching...
No Matches
csd.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "csd.h"
40#include "dpss.h"
41
42//=============================================================================================================
43// EIGEN INCLUDES
44//=============================================================================================================
45
46#include <Eigen/Core>
47#include <Eigen/Dense>
48//#ifndef EIGEN_FFTW_DEFAULT
49//#define EIGEN_FFTW_DEFAULT
50//#endif
51#include <unsupported/Eigen/FFT>
52
53//=============================================================================================================
54// STD INCLUDES
55//=============================================================================================================
56
57#include <cmath>
58#include <complex>
59#include <vector>
60
61//=============================================================================================================
62// USED NAMESPACES
63//=============================================================================================================
64
65using namespace UTILSLIB;
66using namespace Eigen;
67
68namespace
69{
70constexpr double CSD_PI = 3.14159265358979323846;
71}
72
73//=============================================================================================================
74// DEFINE MEMBER METHODS
75//=============================================================================================================
76
77CsdResult Csd::computeMultitaper(const MatrixXd& matData,
78 double sfreq,
79 double fmin,
80 double fmax,
81 double halfBandwidth,
82 int nTapers)
83{
84 const int nChannels = static_cast<int>(matData.rows());
85 const int nTimes = static_cast<int>(matData.cols());
86 const int nFreqsFull = nTimes / 2 + 1;
87
88 if (fmax < 0.0)
89 fmax = sfreq / 2.0;
90
91 // 1. Compute DPSS tapers
92 DpssResult dpss = Dpss::compute(nTimes, halfBandwidth, nTapers);
93 const int nTap = static_cast<int>(dpss.matTapers.rows());
94
95 // Compute eigenvalue weight sum
96 double weightSum = 0.0;
97 for (int t = 0; t < nTap; ++t)
98 weightSum += dpss.vecEigenvalues[t];
99
100 // Build full frequency axis
101 const double freqRes = sfreq / static_cast<double>(nTimes);
102
103 // Determine frequency bin range
104 const int iBinMin = std::max(0, static_cast<int>(std::ceil(fmin / freqRes)));
105 const int iBinMax = std::min(nFreqsFull - 1, static_cast<int>(std::floor(fmax / freqRes)));
106 const int nFreqsSel = iBinMax - iBinMin + 1;
107
108 // Initialize per-frequency CSD accumulators
109 std::vector<MatrixXcd> csdAccum(static_cast<std::size_t>(nFreqsSel),
110 MatrixXcd::Zero(nChannels, nChannels));
111
112 Eigen::FFT<double> fft;
113
114 // 2. For each taper, compute FFT for all channels, then accumulate CSD
115 for (int t = 0; t < nTap; ++t) {
116 const double w = dpss.vecEigenvalues[t];
117
118 // Compute FFT spectra for all channels with this taper
119 MatrixXcd spectra(nChannels, nFreqsFull);
120
121 for (int ch = 0; ch < nChannels; ++ch) {
122 VectorXd tapered = matData.row(ch).transpose().array()
123 * dpss.matTapers.row(t).transpose().array();
124
125 VectorXcd spec;
126 fft.fwd(spec, tapered);
127
128 spectra.row(ch) = spec.head(nFreqsFull).transpose();
129 }
130
131 // Accumulate weighted CSD: X(:,f) * X(:,f)^H for selected frequency bins
132 for (int fi = 0; fi < nFreqsSel; ++fi) {
133 const int fBin = iBinMin + fi;
134 VectorXcd col = spectra.col(fBin);
135 csdAccum[static_cast<std::size_t>(fi)] += w * (col * col.adjoint());
136 }
137 }
138
139 // 3. Assemble result
140 CsdResult result;
141 result.vecFreqs.resize(nFreqsSel);
142 result.csdByFreq.resize(nFreqsSel);
143
144 MatrixXcd meanCsd = MatrixXcd::Zero(nChannels, nChannels);
145
146 for (int fi = 0; fi < nFreqsSel; ++fi) {
147 csdAccum[static_cast<std::size_t>(fi)] /= weightSum;
148 result.csdByFreq[fi] = csdAccum[static_cast<std::size_t>(fi)];
149 result.vecFreqs[fi] = (iBinMin + fi) * freqRes;
150 meanCsd += result.csdByFreq[fi];
151 }
152
153 if (nFreqsSel > 0)
154 meanCsd /= static_cast<double>(nFreqsSel);
155
156 result.matCsd = meanCsd;
157 return result;
158}
159
160//=============================================================================================================
161
162CsdResult Csd::computeFourier(const MatrixXd& matData,
163 double sfreq,
164 double fmin,
165 double fmax,
166 int nFft,
167 double overlap)
168{
169 const int nChannels = static_cast<int>(matData.rows());
170 const int nTimes = static_cast<int>(matData.cols());
171
172 if (nFft > nTimes)
173 nFft = nTimes;
174
175 const int nFreqsFull = nFft / 2 + 1;
176
177 if (fmax < 0.0)
178 fmax = sfreq / 2.0;
179
180 const double freqRes = sfreq / static_cast<double>(nFft);
181 const int iBinMin = std::max(0, static_cast<int>(std::ceil(fmin / freqRes)));
182 const int iBinMax = std::min(nFreqsFull - 1, static_cast<int>(std::floor(fmax / freqRes)));
183 const int nFreqsSel = iBinMax - iBinMin + 1;
184
185 // Build Hann window
186 VectorXd hannWin(nFft);
187 for (int n = 0; n < nFft; ++n) {
188 const double t = static_cast<double>(n) / static_cast<double>(nFft - 1);
189 hannWin[n] = 0.5 * (1.0 - std::cos(2.0 * CSD_PI * t));
190 }
191 const double winPow = hannWin.squaredNorm();
192
193 const int iStep = std::max(1, static_cast<int>(std::round(
194 static_cast<double>(nFft) * (1.0 - overlap))));
195
196 // Initialize per-frequency CSD accumulators
197 std::vector<MatrixXcd> csdAccum(static_cast<std::size_t>(nFreqsSel),
198 MatrixXcd::Zero(nChannels, nChannels));
199
200 Eigen::FFT<double> fft;
201 int nSeg = 0;
202
203 // Process each segment
204 for (int start = 0; start + nFft <= nTimes; start += iStep) {
205 // FFT all channels for this segment
206 MatrixXcd spectra(nChannels, nFreqsFull);
207
208 for (int ch = 0; ch < nChannels; ++ch) {
209 VectorXd seg = matData.row(ch).segment(start, nFft).transpose().array()
210 * hannWin.array();
211
212 VectorXcd spec;
213 fft.fwd(spec, seg);
214
215 spectra.row(ch) = spec.head(nFreqsFull).transpose();
216 }
217
218 // Accumulate CSD for selected frequency bins
219 for (int fi = 0; fi < nFreqsSel; ++fi) {
220 const int fBin = iBinMin + fi;
221 VectorXcd col = spectra.col(fBin);
222 csdAccum[static_cast<std::size_t>(fi)] += col * col.adjoint();
223 }
224
225 ++nSeg;
226 }
227
228 // 3. Assemble result
229 CsdResult result;
230 result.vecFreqs.resize(nFreqsSel);
231 result.csdByFreq.resize(nFreqsSel);
232
233 MatrixXcd meanCsd = MatrixXcd::Zero(nChannels, nChannels);
234
235 const double dNorm = (nSeg > 0) ? 1.0 / (static_cast<double>(nSeg) * winPow) : 0.0;
236
237 for (int fi = 0; fi < nFreqsSel; ++fi) {
238 csdAccum[static_cast<std::size_t>(fi)] *= dNorm;
239 result.csdByFreq[fi] = csdAccum[static_cast<std::size_t>(fi)];
240 result.vecFreqs[fi] = (iBinMin + fi) * freqRes;
241 meanCsd += result.csdByFreq[fi];
242 }
243
244 if (nFreqsSel > 0)
245 meanCsd /= static_cast<double>(nFreqsSel);
246
247 result.matCsd = meanCsd;
248 return result;
249}
250
251//=============================================================================================================
252
253CsdResult Csd::computeMorlet(const MatrixXd& matData,
254 double sfreq,
255 const RowVectorXd& frequencies,
256 int nCycles)
257{
258 const int nChannels = static_cast<int>(matData.rows());
259 const int nTimes = static_cast<int>(matData.cols());
260 const int nFreqs = static_cast<int>(frequencies.size());
261
262 CsdResult result;
263 result.vecFreqs = frequencies;
264 result.csdByFreq.resize(nFreqs);
265
266 MatrixXcd meanCsd = MatrixXcd::Zero(nChannels, nChannels);
267
268 for (int fi = 0; fi < nFreqs; ++fi) {
269 const double f = frequencies[fi];
270 const double sigma = static_cast<double>(nCycles) / (2.0 * CSD_PI * f);
271
272 // Determine wavelet length: ±3σ in samples
273 const int halfLen = static_cast<int>(std::ceil(3.0 * sigma * sfreq));
274 const int wavLen = 2 * halfLen + 1;
275
276 // Build Morlet wavelet
277 VectorXcd wavelet(wavLen);
278 double normFactor = 0.0;
279 for (int n = 0; n < wavLen; ++n) {
280 const double t = (static_cast<double>(n) - static_cast<double>(halfLen)) / sfreq;
281 const double gauss = std::exp(-t * t / (2.0 * sigma * sigma));
282 wavelet[n] = std::complex<double>(gauss * std::cos(2.0 * CSD_PI * f * t),
283 gauss * std::sin(2.0 * CSD_PI * f * t));
284 normFactor += gauss * gauss;
285 }
286 // Normalize wavelet for unit energy
287 normFactor = std::sqrt(normFactor);
288 if (normFactor > 0.0)
289 wavelet /= normFactor;
290
291 // Convolve each channel with the wavelet → complex analytic signal
292 MatrixXcd analytic(nChannels, nTimes);
293
294 for (int ch = 0; ch < nChannels; ++ch) {
295 for (int t = 0; t < nTimes; ++t) {
296 std::complex<double> sum(0.0, 0.0);
297 for (int k = 0; k < wavLen; ++k) {
298 const int idx = t - halfLen + k;
299 if (idx >= 0 && idx < nTimes) {
300 sum += matData(ch, idx) * std::conj(wavelet[k]);
301 }
302 }
303 analytic(ch, t) = sum;
304 }
305 }
306
307 // CSD at this frequency = mean over time of C(:,t) * C(:,t)^H
308 MatrixXcd csdAtFreq = MatrixXcd::Zero(nChannels, nChannels);
309 for (int t = 0; t < nTimes; ++t) {
310 VectorXcd col = analytic.col(t);
311 csdAtFreq += col * col.adjoint();
312 }
313 csdAtFreq /= static_cast<double>(nTimes);
314
315 result.csdByFreq[fi] = csdAtFreq;
316 meanCsd += csdAtFreq;
317 }
318
319 if (nFreqs > 0)
320 meanCsd /= static_cast<double>(nFreqs);
321
322 result.matCsd = meanCsd;
323 return result;
324}
Csd class declaration — Cross-Spectral Density computation.
Dpss class declaration — Discrete Prolate Spheroidal Sequences (Slepian tapers).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Result of a Cross-Spectral Density computation.
Definition csd.h:69
Eigen::MatrixXcd matCsd
n_channels × n_channels mean CSD across selected frequencies
Definition csd.h:70
Eigen::RowVectorXd vecFreqs
Frequency axis (Hz) for csdByFreq entries.
Definition csd.h:72
QVector< Eigen::MatrixXcd > csdByFreq
One n_ch × n_ch CSD matrix per selected frequency bin.
Definition csd.h:71
static CsdResult computeMultitaper(const Eigen::MatrixXd &matData, double sfreq, double fmin=0.0, double fmax=-1.0, double halfBandwidth=4.0, int nTapers=-1)
Definition csd.cpp:77
static CsdResult computeMorlet(const Eigen::MatrixXd &matData, double sfreq, const Eigen::RowVectorXd &frequencies, int nCycles=7)
Definition csd.cpp:253
static CsdResult computeFourier(const Eigen::MatrixXd &matData, double sfreq, double fmin=0.0, double fmax=-1.0, int nFft=256, double overlap=0.5)
Definition csd.cpp:162
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