v2.0.0
Loading...
Searching...
No Matches
fwd_eeg_sphere_model.h
Go to the documentation of this file.
1//=============================================================================================================
36
37#ifndef FWD_EEG_SPHERE_MODEL_H
38#define FWD_EEG_SPHERE_MODEL_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 <QDebug>
59
60#include <memory>
61#include <vector>
62
63//=============================================================================================================
64// DEFINE NAMESPACE FWDLIB
65//=============================================================================================================
66
67namespace FWDLIB
68{
69
70/*
71 * This is the beginning of the specific code
72 */
76struct fitUserRec {
77 Eigen::VectorXd y;
78 Eigen::VectorXd resi;
79 Eigen::MatrixXd M;
80 Eigen::MatrixXd uu;
81 Eigen::MatrixXd vv;
82 Eigen::VectorXd sing;
83 Eigen::VectorXd fn;
84 Eigen::VectorXd w;
85 int nfit;
86 int nterms;
87};
89
90//=============================================================================================================
97{
98public:
99 typedef std::unique_ptr<FwdEegSphereModel> UPtr;
100
101 //=========================================================================================================
106 explicit FwdEegSphereModel();
107
108 //=========================================================================================================
114 explicit FwdEegSphereModel(const FwdEegSphereModel& p_FwdEegSphereModel);
115
116 //=========================================================================================================
130 int nlayer,
131 const Eigen::VectorXf& rads,
132 const Eigen::VectorXf& sigmas);
133
134 //=========================================================================================================
138 virtual ~FwdEegSphereModel();
139
140 //=========================================================================================================
150 static FwdEegSphereModel::UPtr setup_eeg_sphere_model(const QString& eeg_model_file, QString eeg_model_name, float eeg_sphere_rad);
151
152 //=========================================================================================================
161 static fitUser new_fit_user(int nfit, int nterms);
162
163 //=========================================================================================================
173
174 //=========================================================================================================
186 static void next_legen (int n,
187 double x,
188 double &p0,
189 double &p01,
190 double &p1,
191 double &p11);
192
193 //=========================================================================================================
205 static void calc_pot_components(double beta,
206 double cgamma,
207 double &Vrp,
208 double &Vtp,
209 const Eigen::VectorXd& fn,
210 int nterms);
211
212 //=========================================================================================================
228 static int fwd_eeg_multi_spherepot(const Eigen::Vector3f& rd,
229 const Eigen::Vector3f& Q,
230 const Eigen::Matrix<float, Eigen::Dynamic, 3, Eigen::RowMajor>& el,
231 int neeg,
232 Eigen::VectorXf& Vval,
233 void *client);
234
235 //=========================================================================================================
248 static int fwd_eeg_multi_spherepot_coil1(const Eigen::Vector3f& rd,
249 const Eigen::Vector3f& Q,
250 FwdCoilSet& els,
251 Eigen::Ref<Eigen::VectorXf> Vval,
252 void *client);
253
254 //=========================================================================================================
278 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);
279
280 //=========================================================================================================
298 static int fwd_eeg_spherepot_coil_vec(const Eigen::Vector3f& rd, FwdCoilSet& els, Eigen::Ref<Eigen::MatrixXf> Vval_vec, void *client);
299
300 //=========================================================================================================
316 static int fwd_eeg_spherepot_grad_coil( const Eigen::Vector3f& rd,
317 const Eigen::Vector3f& Q,
318 FwdCoilSet& coils,
319 Eigen::Ref<Eigen::VectorXf> Vval,
320 Eigen::Ref<Eigen::VectorXf> xgrad,
321 Eigen::Ref<Eigen::VectorXf> ygrad,
322 Eigen::Ref<Eigen::VectorXf> zgrad,
323 void *client);
324
325 //=========================================================================================================
343 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);
344
345 //=========================================================================================================
360 static int fwd_eeg_spherepot_coil(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, FwdCoilSet& els, Eigen::Ref<Eigen::VectorXf> Vval, void *client);
361
362 //=========================================================================================================
374 bool fwd_setup_eeg_sphere_model(float rad, bool fit_berg_scherg, int nfit);
375
376 //=========================================================================================================
384 static void compose_linear_fitting_data(const Eigen::VectorXd& mu,fitUser u);
385
386 //=========================================================================================================
397 static double compute_linear_parameters(const Eigen::VectorXd& mu, Eigen::VectorXd& lambda, fitUser u);
398
399 //=========================================================================================================
409 static double one_step (const Eigen::VectorXd& mu, const void *user_data);
410
411 //=========================================================================================================
426 int nfit,
427 float &rv);
428
430 int nlayer() const
431 {
432 return layers.size();
433 }
434
435public:
436 QString name;
437 std::vector<FwdEegSphereLayer> layers;
438 Eigen::Vector3f r0;
439
440 Eigen::VectorXd fn;
441 int nterms;
442
443 Eigen::VectorXf mu;
444 Eigen::VectorXf lambda;
445 int nfit;
447};
448
449//=============================================================================================================
450// INLINE DEFINITIONS
451//=============================================================================================================
452} // NAMESPACE FWDLIB
453
454#endif // FWD_EEG_SPHERE_MODEL_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:91
fitUserRec * fitUser
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...
static double compute_linear_parameters(const Eigen::VectorXd &mu, Eigen::VectorXd &lambda, fitUser u)
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 FwdEegSphereModel::UPtr setup_eeg_sphere_model(const QString &eeg_model_file, QString eeg_model_name, float eeg_sphere_rad)
static void calc_pot_components(double beta, double cgamma, double &Vrp, double &Vtp, const Eigen::VectorXd &fn, int nterms)
static fitUser new_fit_user(int nfit, int nterms)
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 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)
bool fwd_eeg_fit_berg_scherg(int nterms, int nfit, float &rv)
static void compose_linear_fitting_data(const Eigen::VectorXd &mu, fitUser u)
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_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_spherepot_coil_vec(const Eigen::Vector3f &rd, FwdCoilSet &els, Eigen::Ref< Eigen::MatrixXf > Vval_vec, void *client)
static FwdEegSphereModel::UPtr fwd_create_eeg_sphere_model(const QString &name, int nlayer, const Eigen::VectorXf &rads, const Eigen::VectorXf &sigmas)
std::vector< FwdEegSphereLayer > layers
static void next_legen(int n, double x, double &p0, double &p01, double &p1, double &p11)
std::unique_ptr< FwdEegSphereModel > UPtr
static int fwd_eeg_spherepot_coil(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &els, Eigen::Ref< Eigen::VectorXf > Vval, void *client)
double fwd_eeg_get_multi_sphere_model_coeff(int n)
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)