v2.0.0
Loading...
Searching...
No Matches
INVLIB::InvBeamformerCompute Class Reference

Core beamformer math. More...

#include <inv_beamformer_compute.h>

Static Public Member Functions

static bool computeBeamformer (const Eigen::MatrixXd &G, const Eigen::MatrixXd &Cm, double reg, int nOrient, BeamformerWeightNorm weightNorm, BeamformerPickOri pickOri, bool reduceRank, BeamformerInversion invMethod, const Eigen::MatrixX3d &nn, Eigen::MatrixXd &W, Eigen::MatrixX3d &maxPowerOri)
static Eigen::VectorXd computePower (const Eigen::MatrixXd &Cm, const Eigen::MatrixXd &W, int nOrient)
static Eigen::MatrixXd symMatPow (const Eigen::MatrixXd &X, double p, bool reduceRank=false)

Detailed Description

Core beamformer math.

Shared computation routines used by both LCMV and DICS beamformers.

All methods are static — no state is kept between calls.

Definition at line 75 of file inv_beamformer_compute.h.

Member Function Documentation

◆ computeBeamformer()

bool InvBeamformerCompute::computeBeamformer ( const Eigen::MatrixXd & G,
const Eigen::MatrixXd & Cm,
double reg,
int nOrient,
BeamformerWeightNorm weightNorm,
BeamformerPickOri pickOri,
bool reduceRank,
BeamformerInversion invMethod,
const Eigen::MatrixX3d & nn,
Eigen::MatrixXd & W,
Eigen::MatrixX3d & maxPowerOri )
static

Compute beamformer spatial filter weights from a leadfield and data covariance (or CSD).

This implements the core LCMV/DICS algorithm:

  1. Regularize and invert covariance
  2. Reshape leadfield per source (n_channels x n_orient)
  3. Optionally reduce leadfield rank
  4. Select orientation (max-power, normal, or keep all)
  5. Invert denominator G^T Cm^{-1} G
  6. Apply weight normalization
Parameters
[in]GLeadfield matrix (n_channels, n_sources * n_orient).
[in]CmData covariance or CSD matrix (n_channels, n_channels). Real-valued.
[in]regRegularization parameter (fraction of trace to add, e.g. 0.05).
[in]nOrientOrientations per source: 1 (fixed) or 3 (free).
[in]weightNormWeight normalization strategy.
[in]pickOriOrientation selection mode.
[in]reduceRankIf true, reduce leadfield rank by 1 (for MEG sphere models).
[in]invMethodDenominator inversion strategy (matrix or scalar).
[in]nnSource normals (n_sources, 3). Used for max-power sign alignment.
[out]WOutput spatial filter weights (n_sources * n_orient_out, n_channels).
[out]maxPowerOriOutput max-power orientations (n_sources, 3). Empty if not max-power.
Returns
True on success.

Definition at line 193 of file inv_beamformer_compute.cpp.

◆ computePower()

VectorXd InvBeamformerCompute::computePower ( const Eigen::MatrixXd & Cm,
const Eigen::MatrixXd & W,
int nOrient )
static

Compute source power from spatial filter weights and a data covariance / CSD matrix.

power_i = trace(W_i Cm W_i^T)

Parameters
[in]CmData covariance or CSD (n_channels, n_channels).
[in]WSpatial filter weights (n_sources * n_orient, n_channels).
[in]nOrientOrientations per source.
Returns
Source power vector (n_sources).

Definition at line 382 of file inv_beamformer_compute.cpp.

◆ symMatPow()

MatrixXd InvBeamformerCompute::symMatPow ( const Eigen::MatrixXd & X,
double p,
bool reduceRank = false )
static

Symmetric matrix power: X^p via eigendecomposition.

Parameters
[in]XSymmetric matrix.
[in]pPower exponent (e.g. -1, -0.5).
[in]reduceRankDrop smallest eigenvalue before exponentiating.
Returns
X^p.

Definition at line 157 of file inv_beamformer_compute.cpp.


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