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

Welch's averaged-periodogram power spectral density estimator. More...

#include <welch_psd.h>

Public Types

enum  WindowType { Hann =0 , Hamming =1 , Blackman =2 , FlatTop =3 }

Static Public Member Functions

static WelchPsdResult compute (const Eigen::MatrixXd &matData, double dSFreq, int iNfft=256, double dOverlap=0.5, WindowType window=Hann, const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi())
static Eigen::RowVectorXd computeVector (const Eigen::RowVectorXd &vecData, double dSFreq, int iNfft=256, double dOverlap=0.5, WindowType window=Hann)
static Eigen::RowVectorXd freqAxis (int iNfft, double dSFreq)

Detailed Description

Welch's averaged-periodogram power spectral density estimator.

Divides each channel into overlapping windowed segments, computes the FFT of each segment, and averages the squared magnitudes. The result is a one-sided PSD normalised so that integrating over frequency recovers the mean signal power.

// 600 Hz data, 512-sample FFT, 50 % overlap, Hann window
WelchPsdResult r = WelchPsd::compute(matData, 600.0, 512);
// r.matPsd → n_channels × 257
// r.vecFreqs → [0, 1.17, 2.34, …, 300] Hz
Result of a Welch PSD computation.
Definition welch_psd.h:64
static WelchPsdResult compute(const Eigen::MatrixXd &matData, double dSFreq, int iNfft=256, double dOverlap=0.5, WindowType window=Hann, const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi())

Definition at line 84 of file welch_psd.h.

Member Enumeration Documentation

◆ WindowType

Window function applied to each segment before the FFT.

Enumerator
Hann 
Hamming 
Blackman 
FlatTop 

Definition at line 88 of file welch_psd.h.

Member Function Documentation

◆ compute()

WelchPsdResult WelchPsd::compute ( const Eigen::MatrixXd & matData,
double dSFreq,
int iNfft = 256,
double dOverlap = 0.5,
WindowType window = Hann,
const Eigen::RowVectorXi & vecPicks = Eigen::RowVectorXi() )
static

Compute Welch PSD for every (selected) channel of a data matrix.

Parameters
[in]matDataData matrix (n_channels × n_samples).
[in]dSFreqSampling frequency in Hz.
[in]iNfftFFT length / segment length in samples (default 256).
[in]dOverlapFractional overlap between adjacent segments, in [0, 1) (default 0.5).
[in]windowWindow function (default Hann).
[in]vecPicksChannel row indices to process; empty = all channels.
Returns
WelchPsdResult with matPsd (n_picks × iNfft/2+1) and vecFreqs.

Definition at line 164 of file welch_psd.cpp.

◆ computeVector()

RowVectorXd WelchPsd::computeVector ( const Eigen::RowVectorXd & vecData,
double dSFreq,
int iNfft = 256,
double dOverlap = 0.5,
WindowType window = Hann )
static

Compute Welch PSD for a single channel row vector.

Parameters
[in]vecDataSingle-channel data (1 × n_samples row vector).
[in]dSFreqSampling frequency in Hz.
[in]iNfftFFT / segment length in samples.
[in]dOverlapFractional segment overlap in [0, 1).
[in]windowWindow function.
Returns
One-sided PSD row vector of length iNfft/2+1.

Definition at line 117 of file welch_psd.cpp.

◆ freqAxis()

RowVectorXd WelchPsd::freqAxis ( int iNfft,
double dSFreq )
static

Build the frequency axis for a given FFT length and sampling frequency.

Parameters
[in]iNfftFFT length.
[in]dSFreqSampling frequency in Hz.
Returns
Row vector of length iNfft/2+1 with frequencies in Hz.

Definition at line 106 of file welch_psd.cpp.


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