Multi-layer spherical head model for EEG forward computation. More...
#include <fwd_eeg_sphere_model.h>
Public Types | |
| typedef std::unique_ptr< FwdEegSphereModel > | UPtr |
Public Member Functions | |
| FwdEegSphereModel () | |
| FwdEegSphereModel (const FwdEegSphereModel &p_FwdEegSphereModel) | |
| virtual | ~FwdEegSphereModel () |
| double | fwd_eeg_get_multi_sphere_model_coeff (int n) |
| bool | fwd_setup_eeg_sphere_model (float rad, bool fit_berg_scherg, int nfit) |
| bool | fwd_eeg_fit_berg_scherg (int nterms, int nfit, float &rv) |
| int | nlayer () const |
Static Public Member Functions | |
| static FwdEegSphereModel::UPtr | fwd_create_eeg_sphere_model (const QString &name, int nlayer, const Eigen::VectorXf &rads, const Eigen::VectorXf &sigmas) |
| static FwdEegSphereModel::UPtr | setup_eeg_sphere_model (const QString &eeg_model_file, QString eeg_model_name, float eeg_sphere_rad) |
| static fitUser | new_fit_user (int nfit, int nterms) |
| static void | next_legen (int n, double x, double &p0, double &p01, double &p1, double &p11) |
| static void | calc_pot_components (double beta, double cgamma, double &Vrp, double &Vtp, const Eigen::VectorXd &fn, int nterms) |
| static int | fwd_eeg_multi_spherepot (const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, const Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > &el, int neeg, Eigen::VectorXf &Vval, void *client) |
| static int | fwd_eeg_multi_spherepot_coil1 (const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &els, Eigen::Ref< Eigen::VectorXf > Vval, void *client) |
| static bool | fwd_eeg_spherepot_vec (const Eigen::Vector3f &rd, const Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > &el, int neeg, Eigen::MatrixXf &Vval_vec, void *client) |
| static int | fwd_eeg_spherepot_coil_vec (const Eigen::Vector3f &rd, FwdCoilSet &els, Eigen::Ref< Eigen::MatrixXf > Vval_vec, void *client) |
| static int | fwd_eeg_spherepot_grad_coil (const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > Vval, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad, void *client) |
| static int | fwd_eeg_spherepot (const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, const Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > &el, int neeg, Eigen::VectorXf &Vval, void *client) |
| static int | fwd_eeg_spherepot_coil (const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &els, Eigen::Ref< Eigen::VectorXf > Vval, void *client) |
| static void | compose_linear_fitting_data (const Eigen::VectorXd &mu, fitUser u) |
| static double | compute_linear_parameters (const Eigen::VectorXd &mu, Eigen::VectorXd &lambda, fitUser u) |
| static double | one_step (const Eigen::VectorXd &mu, const void *user_data) |
Public Attributes | |
| QString | name |
| std::vector< FwdEegSphereLayer > | layers |
| Eigen::Vector3f | r0 |
| Eigen::VectorXd | fn |
| int | nterms |
| Eigen::VectorXf | mu |
| Eigen::VectorXf | lambda |
| int | nfit |
| int | scale_pos |
Multi-layer spherical head model for EEG forward computation.
Implements FwdEegSphereModel (Replaces *fwdEegSphereModel,fwdEegSphereModelRec struct of MNE-C fwd_types.h).
Definition at line 96 of file fwd_eeg_sphere_model.h.
| typedef std::unique_ptr<FwdEegSphereModel> FWDLIB::FwdEegSphereModel::UPtr |
Unique pointer type for FwdEegSphereModel.
Definition at line 99 of file fwd_eeg_sphere_model.h.
|
explicit |
Constructs the Forward EEG Sphere Model
Definition at line 78 of file fwd_eeg_sphere_model.cpp.
|
explicit |
Copy constructor.
| [in] | p_FwdEegSphereModel | Forward EEG Sphere Model which should be copied. |
Definition at line 88 of file fwd_eeg_sphere_model.cpp.
|
virtual |
Destroys the Electric Current Dipole description
Definition at line 156 of file fwd_eeg_sphere_model.cpp.
|
static |
Calculate the radial and tangential potential components for a dipole in a layered sphere model using a Legendre polynomial series expansion.
| [in] | beta | Ratio of dipole distance to field point distance. |
| [in] | cgamma | Cosine of the angle between dipole and field point vectors. |
| [out] | Vrp | Radial component of the potential. |
| [out] | Vtp | Tangential component of the potential. |
| [in] | fn | Pre-computed expansion coefficients. |
| [in] | nterms | Maximum number of terms in the series. |
Definition at line 344 of file fwd_eeg_sphere_model.cpp.
|
static |
Build the data matrix and target vector for the linear part of the Berg-Scherg parameter fit.
| [in] | mu | Distance multipliers for the equivalent dipoles. |
| [in] | u | Fitting workspace (M and y are populated on output). |
Definition at line 949 of file fwd_eeg_sphere_model.cpp.
|
static |
Solve for the optimal lambda (dipole magnitude) parameters given fixed mu values using SVD. Returns the relative variance goodness-of-fit metric.
| [in] | mu | Distance multipliers for the equivalent dipoles. |
| [out] | lambda | Computed dipole magnitudes. |
| [in] | u | Fitting workspace. |
Definition at line 967 of file fwd_eeg_sphere_model.cpp.
|
static |
Create a new multi-layer EEG sphere model structure.
Sorts layers by radius and scales radii relative to the outermost layer.
| [in] | name | Textual identifier for the model. |
| [in] | nlayer | Number of concentric spherical layers. |
| [in] | rads | Radius values for each layer. |
| [in] | sigmas | Conductivity values for each layer. |
Definition at line 119 of file fwd_eeg_sphere_model.cpp.
| bool FwdEegSphereModel::fwd_eeg_fit_berg_scherg | ( | int | nterms, |
| int | nfit, | ||
| float & | rv ) |
Fit Berg-Scherg equivalent spherical model parameters (mu and lambda) by minimizing the difference between actual and approximated series expansions using simplex optimization.
On success, updates the internal mu and lambda member vectors.
| [in] | nterms | Number of terms in the series expansion. |
| [in] | nfit | Number of equivalent dipoles to fit (must be >= 2). |
| [out] | rv | Relative variance of the final fit. |
Definition at line 1046 of file fwd_eeg_sphere_model.cpp.
| double FwdEegSphereModel::fwd_eeg_get_multi_sphere_model_coeff | ( | int | n | ) |
fwd_multi_spherepot.c Get the model depended weighting factor for n
| [in] | n | coefficient to which the expansion shopuld be calculated. |
Definition at line 205 of file fwd_eeg_sphere_model.cpp.
|
static |
Compute EEG potentials at multiple electrodes for a current dipole in a multi-layer spherical head model.
Based on the formulas in Zhang (1995) and Mosher et al. (1996).
| [in] | rd | Dipole position. |
| [in] | Q | Dipole moment. |
| [in] | el | Electrode positions (neeg x 3). |
| [in] | neeg | Number of electrodes. |
| [out] | Vval | Computed potential values at each electrode. |
| [in] | client | Pointer to the FwdEegSphereModel. |
Definition at line 371 of file fwd_eeg_sphere_model.cpp.
|
static |
Calculate EEG potentials at compound electrodes (coils) by summing weighted contributions from individual integration points.
| [in] | rd | Dipole position. |
| [in] | Q | Dipole moment. |
| [in] | els | Set of coils/electrodes. |
| [out] | Vval | Computed potentials (one per coil). |
| [in] | client | Pointer to the FwdEegSphereModel. |
Definition at line 496 of file fwd_eeg_sphere_model.cpp.
|
static |
fwd_multi_spherepot.c
This routine calculates the potentials for a specific dipole direction
This routine uses the acceleration with help of equivalent sources in the homogeneous sphere.
| [in] | rd | Dipole position. |
| [in] | Q | Dipole moment. |
| [in] | el | Electrode positions. |
| [in] | neeg | Number of electrodes. |
| [in] | Vval | The potential values. |
| [in] | client. |
Definition at line 699 of file fwd_eeg_sphere_model.cpp.
|
static |
fwd_multi_spherepot.c
Calculate the EEG in the sphere model using the megCoil structure MEG channels are skipped
| [in] | rd | Dipole position. |
| [in] | Q | Dipole moment. |
| [in] | els | Electrode positions. |
| [in] | Vval | The potential values. |
| [in] | client. |
Definition at line 807 of file fwd_eeg_sphere_model.cpp.
|
static |
fwd_multi_spherepot.c
Calculate the EEG in the sphere model using the fwdCoilSet structure MEG channels are skipped
This routine uses the acceleration with help of equivalent sources in the homogeneous sphere.
| [in] | rd | Dipole position. |
| [in] | els | Electrode positions. |
| [in] | Vval_vec | The potential values; Vval_vec[0][k] potentials given by Q = (1.0,0.0,0.0) at electrode k; Vval_vec[1][k] potentials given by Q = (0.0,1.0,0.0) at electrode k; Vval_vec[2][k] potentials given by Q = (0.0,0.0,1.0) at electrode k. |
| [in] | client. |
Definition at line 634 of file fwd_eeg_sphere_model.cpp.
|
static |
Compute potentials and spatial gradients of the potential field at coils using finite differences.
| [in] | rd | Dipole position. |
| [in] | Q | Dipole moment. |
| [in] | coils | Set of coils/electrodes. |
| [out] | Vval | Computed potentials (one per coil). |
| [out] | xgrad | Gradient with respect to x dipole coordinate. |
| [out] | ygrad | Gradient with respect to y dipole coordinate. |
| [out] | zgrad | Gradient with respect to z dipole coordinate. |
| [in] | client | Pointer to the FwdEegSphereModel. |
Definition at line 664 of file fwd_eeg_sphere_model.cpp.
|
static |
Compute the electric potentials in a set of electrodes in spherically Symmetric head model. This routine calculates the fields for all dipole directions.
The code is based on the formulas presented in
J.C. Moscher, R.M. Leahy, and P.S. Lewis, Matrix Kernels for Modeling of EEG and MEG Data, Los Alamos Technical Report, LA-UR-96-1993, 1996.
This routine uses the acceleration with help of equivalent sources in the homogeneous sphere.
| [in] | rd | Dipole position. |
| [in] | el | Electrode positions. |
| [in] | neeg | Number of electrodes. |
| [in] | Vval_vec | The potential values Vval_vec[0][k] potentials given by Q = (1.0,0.0,0.0) at electrode k; Vval_vec[1][k] potentials given by Q = (0.0,1.0,0.0) at electrode k; Vval_vec[2][k] potentials given by Q = (0.0,0.0,1.0) at electrode k. |
Definition at line 532 of file fwd_eeg_sphere_model.cpp.
| bool FwdEegSphereModel::fwd_setup_eeg_sphere_model | ( | float | rad, |
| bool | fit_berg_scherg, | ||
| int | nfit ) |
fwd_eeg_sphere_models.c
Setup the EEG sphere model calculations
| [in] | rad. | |
| [in] | fit_berg_scherg | If Fit Berg Scherg should be performed. |
| [in] | nfit. |
Definition at line 835 of file fwd_eeg_sphere_model.cpp.
|
static |
Allocate and initialize a fitting workspace for Berg-Scherg parameter estimation.
| [in] | nfit | Number of equivalent dipoles to fit. |
| [in] | nterms | Number of terms in the series expansion. |
Definition at line 186 of file fwd_eeg_sphere_model.cpp.
|
static |
Compute the next Legendre polynomials of the first and second kind using a recursion formula. Self-initializes for n = 0 and n = 1.
| [in] | n | Order of the Legendre polynomial. |
| [in] | x | Evaluation point (cosine of the angle). |
| [out] | p0 | Legendre polynomial of the first kind. |
| [out] | p01 | Previous value of p0. |
| [out] | p1 | Legendre polynomial of the second kind. |
| [out] | p11 | Previous value of p1. |
Definition at line 310 of file fwd_eeg_sphere_model.cpp.
|
inline |
Definition at line 430 of file fwd_eeg_sphere_model.h.
|
static |
Objective function for the simplex optimizer: evaluate the residual sum of squares for a given set of mu values.
Definition at line 1007 of file fwd_eeg_sphere_model.cpp.
|
static |
Set up the desired sphere model for EEG
| [in] | eeg_model_file | Contains the model specifications. |
| [in] | eeg_model_name | Name of the model to use. |
| [in] | eeg_sphere_rad | Outer surface radius. |
Definition at line 162 of file fwd_eeg_sphere_model.cpp.
| Eigen::VectorXd FWDLIB::FwdEegSphereModel::fn |
Coefficients saved to speed up the computations.
Definition at line 440 of file fwd_eeg_sphere_model.h.
| Eigen::VectorXf FWDLIB::FwdEegSphereModel::lambda |
Definition at line 444 of file fwd_eeg_sphere_model.h.
| std::vector<FwdEegSphereLayer> FWDLIB::FwdEegSphereModel::layers |
An array of layers.
Definition at line 437 of file fwd_eeg_sphere_model.h.
| Eigen::VectorXf FWDLIB::FwdEegSphereModel::mu |
The Berg-Scherg equivalence parameters.
Definition at line 443 of file fwd_eeg_sphere_model.h.
| QString FWDLIB::FwdEegSphereModel::name |
Textual identifier.
Definition at line 436 of file fwd_eeg_sphere_model.h.
| int FWDLIB::FwdEegSphereModel::nfit |
How many?.
Definition at line 445 of file fwd_eeg_sphere_model.h.
| int FWDLIB::FwdEegSphereModel::nterms |
How many?.
Definition at line 441 of file fwd_eeg_sphere_model.h.
| Eigen::Vector3f FWDLIB::FwdEegSphereModel::r0 |
The origin.
Definition at line 438 of file fwd_eeg_sphere_model.h.
| int FWDLIB::FwdEegSphereModel::scale_pos |
Scale the positions to the surface of the sphere?.
Definition at line 446 of file fwd_eeg_sphere_model.h.