v2.0.0
Loading...
Searching...
No Matches
UTILSLIB::MorletTfr Class Reference

Complex Morlet wavelet time-frequency representation. More...

#include <morlet_tfr.h>

Static Public Member Functions

static MorletTfrResult compute (const Eigen::RowVectorXd &vecData, double dSFreq, const Eigen::RowVectorXd &vecFreqs, double dNCycles=7.0)
static QVector< MorletTfrResultcomputeMultiChannel (const Eigen::MatrixXd &matData, double dSFreq, const Eigen::RowVectorXd &vecFreqs, double dNCycles=7.0, const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi())

Detailed Description

Complex Morlet wavelet time-frequency representation.

For each requested centre frequency f the signal is convolved (via FFT) with a complex Morlet wavelet whose time-domain standard deviation is σ_t = nCycles / (2π·f) The instantaneous power at every time sample is |convolution|².

// 30 log-spaced frequencies from 4 to 80 Hz, 7 cycles per wavelet
RowVectorXd freqs = RowVectorXd::LinSpaced(30, 4.0, 80.0);
MorletTfrResult r = MorletTfr::compute(vecSignal, 600.0, freqs);
// r.matPower → 30 × n_samples instantaneous-power map
Result of a Morlet TFR computation for one channel.
Definition morlet_tfr.h:67
static MorletTfrResult compute(const Eigen::RowVectorXd &vecData, double dSFreq, const Eigen::RowVectorXd &vecFreqs, double dNCycles=7.0)

Definition at line 88 of file morlet_tfr.h.

Member Function Documentation

◆ compute()

MorletTfrResult MorletTfr::compute ( const Eigen::RowVectorXd & vecData,
double dSFreq,
const Eigen::RowVectorXd & vecFreqs,
double dNCycles = 7.0 )
static

Compute Morlet TFR for a single channel time series.

Parameters
[in]vecDataSingle-channel data row vector (1 × n_samples).
[in]dSFreqSampling frequency in Hz.
[in]vecFreqsCentre frequencies in Hz (row vector).
[in]dNCyclesNumber of wavelet cycles (controls time–frequency trade-off; default 7).
Returns
MorletTfrResult with matPower (n_freqs × n_samples).

Definition at line 109 of file morlet_tfr.cpp.

◆ computeMultiChannel()

QVector< MorletTfrResult > MorletTfr::computeMultiChannel ( const Eigen::MatrixXd & matData,
double dSFreq,
const Eigen::RowVectorXd & vecFreqs,
double dNCycles = 7.0,
const Eigen::RowVectorXi & vecPicks = Eigen::RowVectorXi() )
static

Compute Morlet TFR for every (selected) channel of a data matrix.

Parameters
[in]matDataData matrix (n_channels × n_samples).
[in]dSFreqSampling frequency in Hz.
[in]vecFreqsCentre frequencies in Hz.
[in]dNCyclesNumber of wavelet cycles (default 7).
[in]vecPicksChannel row indices; empty = all channels.
Returns
One MorletTfrResult per selected channel.

Definition at line 166 of file morlet_tfr.cpp.


The documentation for this class was generated from the following files: