v2.0.0
Loading...
Searching...
No Matches
fwd_eeg_sphere_model.h
Go to the documentation of this file.
1//=============================================================================================================
36
37#ifndef FWDEEGSPHEREMODEL_H
38#define FWDEEGSPHEREMODEL_H
39
40//=============================================================================================================
41// INCLUDES
42//=============================================================================================================
43
44#include "fwd_global.h"
46#include "fwd_coil_set.h"
47
48//=============================================================================================================
49// EIGEN INCLUDES
50//=============================================================================================================
51
52#include <Eigen/Core>
53
54//=============================================================================================================
55// QT INCLUDES
56//=============================================================================================================
57
58#include <QSharedPointer>
59#include <QList>
60#include <QDebug>
61
62//=============================================================================================================
63// DEFINE NAMESPACE FWDLIB
64//=============================================================================================================
65
66namespace FWDLIB
67{
68
69/*
70 * This is the beginning of the specific code
71 */
75typedef struct {
76 double *y;
77 double *resi;
78 double **M;
79 double **uu;
80 double **vv;
81 double *sing;
82 double *fn;
83 double *w;
84 int nfit;
85 int nterms;
86} *fitUser,fitUserRec;
87
88//=============================================================================================================
95{
96public:
97 typedef QSharedPointer<FwdEegSphereModel> SPtr;
98 typedef QSharedPointer<const FwdEegSphereModel> ConstSPtr;
99
100 //=========================================================================================================
106 explicit FwdEegSphereModel();
107
108 //=========================================================================================================
115 explicit FwdEegSphereModel(const FwdEegSphereModel& p_FwdEegSphereModel);
116
117 /*
118 * Produce a new sphere model structure
119 */
121 int nlayer,
122 const Eigen::VectorXf& rads,
123 const Eigen::VectorXf& sigmas);
124
125 //=========================================================================================================
130 virtual ~FwdEegSphereModel();
131
132 //=========================================================================================================
143 static FwdEegSphereModel* setup_eeg_sphere_model(const QString& eeg_model_file, QString eeg_model_name, float eeg_sphere_rad);
144
145 static fitUser new_fit_user(int nfit, int nterms);
146
147 //=========================================================================================================
158
159 static void next_legen (int n,
160 double x,
161 double *p0, /* Input: P0(n-1) Output: P0(n) */
162 double *p01, /* Input: P0(n-2) Output: P0(n-1) */
163 double *p1, /* Input: P1(n-1) Output: P1(n) */
164 double *p11);
165
166 static void calc_pot_components(double beta, /* rd/r */
167 double cgamma, /* Cosine of the angle between
168 * the source and field points */
169 double *Vrp, /* Potential component for the radial dipole */
170 double *Vtp, /* Potential component for the tangential dipole */
171 const Eigen::VectorXd& fn,
172 int nterms);
173
174 static int fwd_eeg_multi_spherepot(float *rd, /* Dipole position */
175 float *Q, /* Dipole moment */
176 float **el, /* Electrode positions */
177 int neeg, /* Number of electrodes */
178 float *Vval, /* The potential values */
179 void *client);
180
181 static int fwd_eeg_multi_spherepot_coil1(float *rd, /* Dipole position */
182 float *Q, /* Dipole moment */
183 FwdCoilSet* els, /* Electrode positions */
184 float *Vval, /* The potential values */
185 void *client);
186
187 //=========================================================================================================
212 static bool fwd_eeg_spherepot_vec (float *rd, float **el, int neeg, float **Vval_vec, void *client);
213
214 //=========================================================================================================
232 static int fwd_eeg_spherepot_coil_vec(float *rd, FwdCoilSet* els, float **Vval_vec, void *client);
233
234 static int fwd_eeg_spherepot_grad_coil( float *rd, /* The dipole location */
235 float Q[], /* The dipole components (xyz) */
236 FwdCoilSet* coils, /* The coil definitions */
237 float Vval[], /* Results */
238 float xgrad[], /* The derivatives with respect to */
239 float ygrad[], /* the dipole position coordinates */
240 float zgrad[],
241 void *client);
242
243 //=========================================================================================================
261 static int fwd_eeg_spherepot( float *rd, float *Q, float **el, int neeg, Eigen::VectorXf& Vval, void *client);
262
263 //=========================================================================================================
278 static int fwd_eeg_spherepot_coil(float *rd, float *Q, FwdCoilSet* els, float *Vval, void *client);
279
280 //=========================================================================================================
292 bool fwd_setup_eeg_sphere_model(float rad, bool fit_berg_scherg, int nfit);
293
294 // fwd_fit_berg_scherg.c
295
296 static void compose_linear_fitting_data(const Eigen::VectorXd& mu,fitUser u);
297
298 // fwd_fit_berg_scherg.c
299 /*
300 * Compute the best-fitting linear parameters
301 * Return the corresponding RV
302 */
303 static double compute_linear_parameters(const Eigen::VectorXd& mu, Eigen::VectorXd& lambda, fitUser u);
304
305 // fwd_fit_berg_scherg.c
306 /*
307 * Evaluate the residual sum of squares fit for one set of
308 * mu values
309 */
310 static double one_step (const Eigen::VectorXd& mu, const void *user_data);
311
312 /*
313 * This routine fits the Berg-Scherg equivalent spherical model
314 * dipole parameters by minimizing the difference between the
315 * actual and approximative series expansions
316 */
317 bool fwd_eeg_fit_berg_scherg(int nterms, /* Number of terms to use in the series expansion
318 * when fitting the parameters */
319 int nfit, /* Number of equivalent dipoles to fit */
320 float &rv);
321
323 int nlayer() const
324 {
325 return layers.size();
326 }
327
328public:
329 QString name;
330 QList<FwdEegSphereLayer> layers;
331 Eigen::Vector3f r0;
332
333 Eigen::VectorXd fn;
334 int nterms;
335
336 Eigen::VectorXf mu;
337 Eigen::VectorXf lambda;
338 int nfit;
340
341// ### OLD STRUCT ###
342// typedef struct {
343// char *name; /* Textual identifier */
344// int nlayer; /* Number of layers */
345// fwdEegSphereLayer layers; /* An array of layers */
346// float r0[3]; /* The origin */
347
348// double *fn; /* Coefficients saved to speed up the computations */
349// int nterms; /* How many? */
350
351// float *mu; /* The Berg-Scherg equivalence parameters */
352// float *lambda;
353// int nfit; /* How many? */
354// int scale_pos; /* Scale the positions to the surface of the sphere? */
355// } *fwdEegSphereModel,fwdEegSphereModelRec;
356};
357
358//=============================================================================================================
359// INLINE DEFINITIONS
360//=============================================================================================================
361} // NAMESPACE FWDLIB
362
363#endif // FWDEEGSPHEREMODEL_H
forward library export/import macros.
#define FWDSHARED_EXPORT
Definition fwd_global.h:53
FwdCoilSet class declaration.
FwdEegSphereLayer class declaration.
Forward modelling (BEM, MEG/EEG lead fields).
Definition compute_fwd.h:95
struct FWDLIB::fitUser fitUserRec
Collection of FwdCoil objects representing a full MEG or EEG sensor array.
Workspace for the linear least-squares fit of Berg-Scherg parameters in the EEG sphere model (SVD mat...
QList< FwdEegSphereLayer > layers
static double compute_linear_parameters(const Eigen::VectorXd &mu, Eigen::VectorXd &lambda, fitUser u)
QSharedPointer< FwdEegSphereModel > SPtr
static bool fwd_eeg_spherepot_vec(float *rd, float **el, int neeg, float **Vval_vec, void *client)
static int fwd_eeg_multi_spherepot_coil1(float *rd, float *Q, FwdCoilSet *els, float *Vval, void *client)
static int fwd_eeg_multi_spherepot(float *rd, float *Q, float **el, int neeg, float *Vval, void *client)
QSharedPointer< const FwdEegSphereModel > ConstSPtr
static FwdEegSphereModel * fwd_create_eeg_sphere_model(const QString &name, int nlayer, const Eigen::VectorXf &rads, const Eigen::VectorXf &sigmas)
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 int fwd_eeg_spherepot_coil_vec(float *rd, FwdCoilSet *els, float **Vval_vec, void *client)
bool fwd_eeg_fit_berg_scherg(int nterms, int nfit, float &rv)
static int fwd_eeg_spherepot_coil(float *rd, float *Q, FwdCoilSet *els, float *Vval, void *client)
static void compose_linear_fitting_data(const Eigen::VectorXd &mu, fitUser u)
static int fwd_eeg_spherepot(float *rd, float *Q, float **el, int neeg, Eigen::VectorXf &Vval, void *client)
static FwdEegSphereModel * setup_eeg_sphere_model(const QString &eeg_model_file, QString eeg_model_name, float eeg_sphere_rad)
static int fwd_eeg_spherepot_grad_coil(float *rd, float Q[], FwdCoilSet *coils, float Vval[], float xgrad[], float ygrad[], float zgrad[], void *client)
double fwd_eeg_get_multi_sphere_model_coeff(int n)
static void calc_pot_components(double beta, double cgamma, double *Vrp, double *Vtp, const Eigen::VectorXd &fn, int nterms)
bool fwd_setup_eeg_sphere_model(float rad, bool fit_berg_scherg, int nfit)
static double one_step(const Eigen::VectorXd &mu, const void *user_data)