Cross-Spectral Density (CSD) estimator.
More...
#include <csd.h>
|
| static CsdResult | computeMultitaper (const Eigen::MatrixXd &matData, double sfreq, double fmin=0.0, double fmax=-1.0, double halfBandwidth=4.0, int nTapers=-1) |
| static CsdResult | computeFourier (const Eigen::MatrixXd &matData, double sfreq, double fmin=0.0, double fmax=-1.0, int nFft=256, double overlap=0.5) |
| static CsdResult | computeMorlet (const Eigen::MatrixXd &matData, double sfreq, const Eigen::RowVectorXd &frequencies, int nCycles=7) |
Cross-Spectral Density (CSD) estimator.
Provides three methods for computing the CSD matrix between channels: multitaper (DPSS tapers), Fourier (Welch-style segmented), and Morlet wavelet.
Result of a Cross-Spectral Density computation.
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 at line 89 of file csd.h.
◆ computeFourier()
| CsdResult Csd::computeFourier |
( |
const Eigen::MatrixXd & | matData, |
|
|
double | sfreq, |
|
|
double | fmin = 0.0, |
|
|
double | fmax = -1.0, |
|
|
int | nFft = 256, |
|
|
double | overlap = 0.5 ) |
|
static |
Compute CSD using Welch-style Fourier segmented estimation.
- Parameters
-
| [in] | matData | Data matrix (n_channels × n_samples). |
| [in] | sfreq | Sampling frequency in Hz. |
| [in] | fmin | Minimum frequency of interest in Hz (default 0). |
| [in] | fmax | Maximum frequency in Hz; -1 → Nyquist (default -1). |
| [in] | nFft | FFT / segment length in samples (default 256). |
| [in] | overlap | Fractional overlap between adjacent segments, in [0, 1) (default 0.5). |
- Returns
- CsdResult with per-frequency and mean CSD matrices.
Definition at line 162 of file csd.cpp.
◆ computeMorlet()
| CsdResult Csd::computeMorlet |
( |
const Eigen::MatrixXd & | matData, |
|
|
double | sfreq, |
|
|
const Eigen::RowVectorXd & | frequencies, |
|
|
int | nCycles = 7 ) |
|
static |
Compute CSD using Morlet wavelet convolution.
- Parameters
-
| [in] | matData | Data matrix (n_channels × n_samples). |
| [in] | sfreq | Sampling frequency in Hz. |
| [in] | frequencies | Target frequencies in Hz. |
| [in] | nCycles | Number of cycles in the Morlet wavelet (default 7). |
- Returns
- CsdResult with one CSD matrix per target frequency.
Definition at line 253 of file csd.cpp.
◆ computeMultitaper()
| CsdResult Csd::computeMultitaper |
( |
const Eigen::MatrixXd & | matData, |
|
|
double | sfreq, |
|
|
double | fmin = 0.0, |
|
|
double | fmax = -1.0, |
|
|
double | halfBandwidth = 4.0, |
|
|
int | nTapers = -1 ) |
|
static |
Compute CSD using multitaper (DPSS) spectral estimation.
- Parameters
-
| [in] | matData | Data matrix (n_channels × n_samples). |
| [in] | sfreq | Sampling frequency in Hz. |
| [in] | fmin | Minimum frequency of interest in Hz (default 0). |
| [in] | fmax | Maximum frequency in Hz; -1 → Nyquist (default -1). |
| [in] | halfBandwidth | Half-bandwidth parameter (NW) for DPSS tapers (default 4). |
| [in] | nTapers | Number of tapers; -1 → floor(2*halfBandwidth - 1) (default -1). |
- Returns
- CsdResult with per-frequency and mean CSD matrices.
Definition at line 77 of file csd.cpp.
The documentation for this class was generated from the following files: