v2.0.0
Loading...
Searching...
No Matches
fwd_bem_model.h
Go to the documentation of this file.
1//=============================================================================================================
36
37#ifndef FWD_BEM_MODEL_H
38#define FWD_BEM_MODEL_H
39
40//=============================================================================================================
41// INCLUDES
42//=============================================================================================================
43
44#include "fwd_global.h"
45#include "fwd_coil_set.h"
46
48#include <fiff/fiff_dir_node.h>
49#include <fiff/fiff_tag.h>
51
52#include <memory>
53#include <vector>
54
55//=============================================================================================================
56// EIGEN INCLUDES
57//=============================================================================================================
58
59#include <Eigen/Core>
60
61//=============================================================================================================
62// QT INCLUDES
63//=============================================================================================================
64
65#include <QString>
66
67//=============================================================================================================
68// BEM METHOD CONSTANTS
69//=============================================================================================================
70
71namespace FWDLIB
72{
73
74constexpr int FWD_BEM_UNKNOWN = -1;
75constexpr int FWD_BEM_CONSTANT_COLL = 1;
76constexpr int FWD_BEM_LINEAR_COLL = 2;
77
78constexpr double FWD_BEM_IP_APPROACH_LIMIT = 0.1;
79
80constexpr int FWD_BEM_LIN_FIELD_SIMPLE = 1;
82constexpr int FWD_BEM_LIN_FIELD_URANKAR = 3;
83
84} // namespace FWDLIB
85
86//=============================================================================================================
87// FORWARD DECLARATIONS
88//=============================================================================================================
89
90namespace MNELIB
91{
92 class MNETriangle;
93 class MNESurface;
94 class MNESourceSpace;
95 class MNECTFCompDataSet;
96 class MNENamedMatrix;
97}
98
99namespace FIFFLIB {
100 class FiffNamedMatrix;
101}
102//=============================================================================================================
103// DEFINE NAMESPACE FWDLIB
104//=============================================================================================================
105
106namespace FWDLIB
107{
108
109//=============================================================================================================
110// FWDLIB FORWARD DECLARATIONS
111//=============================================================================================================
112
114class FwdThreadArg;
115
116//=============================================================================================================
129{
130public:
131 typedef std::unique_ptr<FwdBemModel> UPtr;
132
133 //=========================================================================================================
140 explicit FwdBemModel();
141
142 //=========================================================================================================
149 virtual ~FwdBemModel();
150
151 //=========================================================================================================
159
166 static QString fwd_bem_make_bem_sol_name(const QString& name);
167
168 //============================= fwd_bem_model.c =============================
169
176 static const QString& fwd_bem_explain_surface(int kind);
177
184 static const QString& fwd_bem_explain_method(int method);
185
189 static int get_int( FIFFLIB::FiffStream::SPtr& stream, const FIFFLIB::FiffDirNode::SPtr& node,int what,int *res);
190
198
199 //=========================================================================================================
207 static FwdBemModel::UPtr fwd_bem_load_surfaces(const QString& name,
208 const std::vector<int>& kinds);
209
210 //=========================================================================================================
219 static FwdBemModel::UPtr fwd_bem_load_homog_surface(const QString& name);
220
221 //=========================================================================================================
228 static FwdBemModel::UPtr fwd_bem_load_three_layer_surfaces(const QString& name);
229
230 //=========================================================================================================
241 int fwd_bem_load_solution(const QString& name, int bem_method);
242
243 //=========================================================================================================
251
252 //============================= dipole_fit_guesses.c =============================
253
254 //=========================================================================================================
269 static std::unique_ptr<MNELIB::MNESurface> make_guesses(MNELIB::MNESurface* guess_surf,
270 float guessrad,
271 const Eigen::Vector3f& guess_r0,
272 float grid,
273 float exclude,
274 float mindist);
275
276 //============================= fwd_bem_linear_collocation.c =============================
277
278 /*
279 * The following approach is based on:
280 *
281 * de Munck JC: "A linear discretization of the volume conductor boundary integral equation
282 * using analytically integrated elements",
283 * IEEE Trans Biomed Eng. 1992 39(9) : 986 - 990
284 */
285
286 //=========================================================================================================
294 static double calc_beta(const Eigen::Vector3d& rk, const Eigen::Vector3d& rk1);
295
296 //=========================================================================================================
304 static void lin_pot_coeff(const Eigen::Vector3f& from,
306 Eigen::Vector3d& omega);
307
308 //=========================================================================================================
316 Eigen::MatrixXf& mat);
317
318 //=========================================================================================================
325 static Eigen::MatrixXf fwd_bem_lin_pot_coeff(const std::vector<MNELIB::MNESurface*>& surfs);
326
327 //=========================================================================================================
334
335 //============================= fwd_bem_solution.c =============================
336
337 //=========================================================================================================
350 static Eigen::MatrixXf fwd_bem_multi_solution(Eigen::MatrixXf& solids,
351 const Eigen::MatrixXf *gamma,
352 int nsurf,
353 const Eigen::VectorXi& ntri);
354
355 //=========================================================================================================
363 static Eigen::MatrixXf fwd_bem_homog_solution(Eigen::MatrixXf& solids, int ntri);
364
365 //=========================================================================================================
378 static void fwd_bem_ip_modify_solution(Eigen::MatrixXf &solution,
379 Eigen::MatrixXf& ip_solution,
380 float ip_mult,
381 int nsurf,
382 const Eigen::VectorXi &ntri);
383
384 //============================= fwd_bem_constant_collocation.c =============================
385
386 //=========================================================================================================
396 static int fwd_bem_check_solids(const Eigen::MatrixXf& angles, int ntri1, int ntri2, float desired);
397
398 //=========================================================================================================
405 static Eigen::MatrixXf fwd_bem_solid_angles(const std::vector<MNELIB::MNESurface*>& surfs);
406
407 //=========================================================================================================
414
415 //============================= fwd_bem_model.c =============================
416
417 //=========================================================================================================
425
426 //=========================================================================================================
438 int fwd_bem_load_recompute_solution(const QString& name,
439 int bem_method,
440 int force_recompute);
441
442 //============================= fwd_bem_pot.c =============================
443
444 //=========================================================================================================
457 static float fwd_bem_inf_field(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, const Eigen::Vector3f& rp, const Eigen::Vector3f& dir);
458
459 //=========================================================================================================
470 static float fwd_bem_inf_pot(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, const Eigen::Vector3f& rp);
471
472 //=========================================================================================================
483
484 //=========================================================================================================
496 void fwd_bem_pot_grad_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
497 FwdCoilSet* els, int all_surfs,
498 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad);
499
500 //=========================================================================================================
510 void fwd_bem_lin_pot_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
511 FwdCoilSet* els, int all_surfs,
512 Eigen::Ref<Eigen::VectorXf> pot);
513
514 //=========================================================================================================
526 void fwd_bem_lin_pot_grad_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
527 FwdCoilSet* els, int all_surfs,
528 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad);
529
530 //=========================================================================================================
540 void fwd_bem_pot_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
541 FwdCoilSet* els, int all_surfs,
542 Eigen::Ref<Eigen::VectorXf> pot);
543
544 //=========================================================================================================
557 static int fwd_bem_pot_els(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
558 FwdCoilSet& els, Eigen::Ref<Eigen::VectorXf> pot,
559 void *client);
560
561 //=========================================================================================================
577 static int fwd_bem_pot_grad_els(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
578 FwdCoilSet& els, Eigen::Ref<Eigen::VectorXf> pot,
579 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad,
580 void *client);
581
582 //============================= fwd_bem_field.c =============================
583
584 /*
585 * Integration formulas from:
586 * L. Urankar, "Common compact analytical formulas for computation of
587 * geometry integrals on a basic Cartesian sub-domain in boundary and
588 * volume integral methods", Engineering Analysis with Boundary Elements,
589 * 7(3), 1990, 124-129.
590 */
591
592 //=========================================================================================================
602 static void calc_f(const Eigen::Vector3d& xx, const Eigen::Vector3d& yy,
603 Eigen::Vector3d& f0, Eigen::Vector3d& fx, Eigen::Vector3d& fy);
604
605 //=========================================================================================================
616 static void calc_magic(double u, double z,
617 double A, double B,
618 Eigen::Vector3d& beta, double& D);
619
620 //=========================================================================================================
634 static void field_integrals(const Eigen::Vector3f& from,
636 double& I1p,
637 Eigen::Vector2d& T, Eigen::Vector2d& S1, Eigen::Vector2d& S2,
638 Eigen::Vector3d& f0, Eigen::Vector3d& fx, Eigen::Vector3d& fy);
639
640 //=========================================================================================================
649 static double one_field_coeff(const Eigen::Vector3f& dest, const Eigen::Vector3f& normal,
651
652 //=========================================================================================================
659 Eigen::MatrixXf fwd_bem_field_coeff(FwdCoilSet* coils);
660
661 /*
662 * Linear field formulas from:
663 * Ferguson et al., "A Complete Linear Discretization for Calculating
664 * the Magnetic Field Using the Boundary-Element Method",
665 * IEEE Trans. Biomed. Eng., submitted.
666 */
667
668 //=========================================================================================================
676 static double calc_gamma(const Eigen::Vector3d& rk, const Eigen::Vector3d& rk1);
677
678 //=========================================================================================================
687 static void fwd_bem_one_lin_field_coeff_ferg(const Eigen::Vector3f& dest, const Eigen::Vector3f& dir,
689 Eigen::Vector3d& res);
690
691 //=========================================================================================================
700 static void fwd_bem_one_lin_field_coeff_uran(const Eigen::Vector3f& dest, const Eigen::Vector3f& dir,
702 Eigen::Vector3d& res);
703
704 //=========================================================================================================
713 static void fwd_bem_one_lin_field_coeff_simple(const Eigen::Vector3f& dest, const Eigen::Vector3f& normal,
714 MNELIB::MNETriangle& source,
715 Eigen::Vector3d& res);
716
717 //=========================================================================================================
721 typedef void (*linFieldIntFunc)(const Eigen::Vector3f& dest, const Eigen::Vector3f& dir,
722 MNELIB::MNETriangle& tri, Eigen::Vector3d& res);
723
724 //=========================================================================================================
732 Eigen::MatrixXf fwd_bem_lin_field_coeff(FwdCoilSet* coils, int method);
733
734 //=========================================================================================================
745
746 static constexpr double MAG_FACTOR = 1e-7;
747
748 //=========================================================================================================
757 void fwd_bem_lin_field_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
758 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> B);
759
760 //=========================================================================================================
769 void fwd_bem_field_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
770 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> B);
771
772 //=========================================================================================================
783 void fwd_bem_field_grad_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
784 FwdCoilSet& coils,
785 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad);
786
787 //=========================================================================================================
800 static float fwd_bem_inf_field_der(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, const Eigen::Vector3f& rp,
801 const Eigen::Vector3f& dir, const Eigen::Vector3f& comp);
802
803 //=========================================================================================================
813 static float fwd_bem_inf_pot_der(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, const Eigen::Vector3f& rp,
814 const Eigen::Vector3f& comp);
815
816 //=========================================================================================================
827 void fwd_bem_lin_field_grad_calc(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
828 FwdCoilSet& coils,
829 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad);
830
831 //=========================================================================================================
845 static int fwd_bem_field(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
846 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> B,
847 void *client);
848
849 //=========================================================================================================
865 static int fwd_bem_field_grad(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
866 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> Bval,
867 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad,
868 void *client);
869
870 //============================= compute_forward.c =============================
871
872 //=========================================================================================================
882
883 //=========================================================================================================
899 int compute_forward_meg(std::vector<std::unique_ptr<MNELIB::MNESourceSpace>>& spaces,
900 FwdCoilSet* coils,
901 FwdCoilSet* comp_coils,
902 MNELIB::MNECTFCompDataSet* comp_data,
903 bool fixed_ori,
904 const Eigen::Vector3f& r0,
905 bool use_threads,
907 FIFFLIB::FiffNamedMatrix& resp_grad,
908 bool bDoGRad);
909
910 //=========================================================================================================
924 int compute_forward_eeg(std::vector<std::unique_ptr<MNELIB::MNESourceSpace>>& spaces,
925 FwdCoilSet* els,
926 bool fixed_ori,
927 FwdEegSphereModel* eeg_model,
928 bool use_threads,
930 FIFFLIB::FiffNamedMatrix& resp_grad,
931 bool bDoGrad);
932
933 //============================= fwd_spherefield.c =============================
934
935 //=========================================================================================================
948 static int fwd_sphere_field(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
949 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> Bval,
950 void *client);
951
952 //=========================================================================================================
965 static int fwd_sphere_field_vec(const Eigen::Vector3f& rd,
966 FwdCoilSet& coils, Eigen::Ref<Eigen::MatrixXf> Bval,
967 void *client);
968
969 //=========================================================================================================
985 static int fwd_sphere_field_grad(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q,
986 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> Bval,
987 Eigen::Ref<Eigen::VectorXf> xgrad, Eigen::Ref<Eigen::VectorXf> ygrad, Eigen::Ref<Eigen::VectorXf> zgrad,
988 void *client);
989
990 //============================= fwd_mag_dipole_field.c =============================
991
992 //=========================================================================================================
1003 static int fwd_mag_dipole_field(const Eigen::Vector3f& rm, const Eigen::Vector3f& M,
1004 FwdCoilSet& coils, Eigen::Ref<Eigen::VectorXf> Bval,
1005 void *client);
1006
1007 //=========================================================================================================
1019 static int fwd_mag_dipole_field_vec(const Eigen::Vector3f& rm,
1020 FwdCoilSet& coils, Eigen::Ref<Eigen::MatrixXf> Bval,
1021 void *client);
1022
1023public:
1024 QString surf_name;
1025
1026 std::vector<std::shared_ptr<MNELIB::MNESurface>> surfs;
1027
1028 Eigen::VectorXi ntri;
1029 Eigen::VectorXi np;
1030 int nsurf;
1031
1032 Eigen::VectorXf sigma;
1033 Eigen::MatrixXf gamma;
1034 Eigen::VectorXf source_mult;
1035 Eigen::VectorXf field_mult;
1036
1038 QString sol_name;
1039
1040 Eigen::MatrixXf solution;
1041 Eigen::VectorXf v0;
1042 int nsol;
1043
1045
1048};
1049
1050//=============================================================================================================
1051// INLINE DEFINITIONS
1052//=============================================================================================================
1053} // NAMESPACE FWDLIB
1054
1055#endif // FWD_BEM_MODEL_H
FiffTag class declaration, which provides fiff tag I/O and processing methods.
FiffNamedMatrix class declaration.
FiffDirNode class declaration, which provides fiff dir tree processing methods.
FiffCoordTrans class declaration.
Forward library export/import macros.
#define FWDSHARED_EXPORT
Definition fwd_global.h:53
FwdCoilSet class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Forward modelling (BEM, MEG/EEG lead fields).
Definition compute_fwd.h:91
constexpr int FWD_BEM_CONSTANT_COLL
constexpr double FWD_BEM_IP_APPROACH_LIMIT
constexpr int FWD_BEM_LIN_FIELD_URANKAR
constexpr int FWD_BEM_LINEAR_COLL
constexpr int FWD_BEM_LIN_FIELD_FERGUSON
constexpr int FWD_BEM_UNKNOWN
constexpr int FWD_BEM_LIN_FIELD_SIMPLE
Coordinate transformation description.
QSharedPointer< FiffDirNode > SPtr
QSharedPointer< FiffStream > SPtr
static FwdBemModel::UPtr fwd_bem_load_three_layer_surfaces(const QString &name)
Load a three-layer BEM model (scalp, outer skull, inner skull) from a FIFF file.
static void field_integrals(const Eigen::Vector3f &from, MNELIB::MNETriangle &to, double &I1p, Eigen::Vector2d &T, Eigen::Vector2d &S1, Eigen::Vector2d &S2, Eigen::Vector3d &f0, Eigen::Vector3d &fx, Eigen::Vector3d &fy)
Compute the geometry integrals for the magnetic field from a triangle.
void fwd_bem_field_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > B)
Compute BEM magnetic fields at coils using constant collocation.
static double calc_gamma(const Eigen::Vector3d &rk, const Eigen::Vector3d &rk1)
Compute the gamma angle for the linear field integration (Ferguson).
Eigen::VectorXi np
void fwd_bem_pot_grad_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet *els, int all_surfs, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad)
Compute the gradient of BEM potentials with respect to dipole position (constant collocation).
Eigen::VectorXf source_mult
static double calc_beta(const Eigen::Vector3d &rk, const Eigen::Vector3d &rk1)
Compute the beta angle used in the linear collocation integration.
static int get_int(FIFFLIB::FiffStream::SPtr &stream, const FIFFLIB::FiffDirNode::SPtr &node, int what, int *res)
Read an integer tag from a FIFF node.
Eigen::MatrixXf gamma
static void fwd_bem_one_lin_field_coeff_uran(const Eigen::Vector3f &dest, const Eigen::Vector3f &dir, MNELIB::MNETriangle &tri, Eigen::Vector3d &res)
Compute linear field coefficients using the Urankar method.
static Eigen::MatrixXf fwd_bem_multi_solution(Eigen::MatrixXf &solids, const Eigen::MatrixXf *gamma, int nsurf, const Eigen::VectorXi &ntri)
Compute the multi-surface BEM solution from solid-angle coefficients.
static FwdBemModel::UPtr fwd_bem_load_surfaces(const QString &name, const std::vector< int > &kinds)
Load BEM surfaces of specified kinds from a FIFF file.
static void correct_auto_elements(MNELIB::MNESurface &surf, Eigen::MatrixXf &mat)
Correct the auto (self-coupling) elements of the linear collocation matrix.
static int fwd_sphere_field_grad(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > Bval, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad, void *client)
Callback: compute the spherical-model magnetic field and its position gradient at coils.
static int fwd_bem_check_solids(const Eigen::MatrixXf &angles, int ntri1, int ntri2, float desired)
Verify that solid-angle sums match the expected value.
int fwd_bem_load_solution(const QString &name, int bem_method)
Load a pre-computed BEM solution from a FIFF file.
static int fwd_mag_dipole_field_vec(const Eigen::Vector3f &rm, FwdCoilSet &coils, Eigen::Ref< Eigen::MatrixXf > Bval, void *client)
Callback: compute the vector magnetic field of a magnetic dipole at coils.
int fwd_bem_compute_solution(int bem_method)
Compute the BEM solution matrix using the specified method.
static constexpr double MAG_FACTOR
int fwd_bem_specify_coils(FwdCoilSet *coils)
Precompute the coil-specific BEM solution for MEG.
static Eigen::MatrixXf fwd_bem_solid_angles(const std::vector< MNELIB::MNESurface * > &surfs)
Compute the solid-angle matrix for all BEM surfaces.
int fwd_bem_specify_els(FwdCoilSet *els)
Precompute the electrode-specific BEM solution.
static void fwd_bem_one_lin_field_coeff_simple(const Eigen::Vector3f &dest, const Eigen::Vector3f &normal, MNELIB::MNETriangle &source, Eigen::Vector3d &res)
Compute linear field coefficients using the simple (direct) method.
std::vector< std::shared_ptr< MNELIB::MNESurface > > surfs
Eigen::VectorXf field_mult
std::unique_ptr< FwdBemModel > UPtr
static void fwd_bem_one_lin_field_coeff_ferg(const Eigen::Vector3f &dest, const Eigen::Vector3f &dir, MNELIB::MNETriangle &tri, Eigen::Vector3d &res)
Compute linear field coefficients using the Ferguson method.
static QString fwd_bem_make_bem_sol_name(const QString &name)
Build a standard BEM solution file name from a model name.
static Eigen::MatrixXf fwd_bem_homog_solution(Eigen::MatrixXf &solids, int ntri)
Compute the homogeneous (single-layer) BEM solution.
void fwd_bem_lin_field_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > B)
Compute BEM magnetic fields at coils using linear collocation.
static int fwd_mag_dipole_field(const Eigen::Vector3f &rm, const Eigen::Vector3f &M, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > Bval, void *client)
Callback: compute the magnetic field of a magnetic dipole at coils.
MNELIB::MNESurface * fwd_bem_find_surface(int kind)
Find a surface of the given kind in this BEM model.
void fwd_bem_free_solution()
Release the potential solution matrix and associated workspace.
Eigen::MatrixXf fwd_bem_field_coeff(FwdCoilSet *coils)
Assemble the constant-collocation magnetic field coefficient matrix.
static int fwd_sphere_field(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > Bval, void *client)
Callback: compute the spherical-model magnetic field at coils.
static void lin_pot_coeff(const Eigen::Vector3f &from, MNELIB::MNETriangle &to, Eigen::Vector3d &omega)
Compute the linear potential coefficients for one source-destination pair.
Eigen::VectorXi ntri
void fwd_bem_field_grad_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad)
Compute the gradient of BEM magnetic fields with respect to dipole position (constant collocation).
static const QString & fwd_bem_explain_method(int method)
Return a human-readable label for a BEM method.
void fwd_bem_pot_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet *els, int all_surfs, Eigen::Ref< Eigen::VectorXf > pot)
Compute BEM potentials at electrodes using constant collocation.
void fwd_bem_lin_field_grad_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad)
Compute the gradient of BEM magnetic fields with respect to dipole position (linear collocation).
static void calc_f(const Eigen::Vector3d &xx, const Eigen::Vector3d &yy, Eigen::Vector3d &f0, Eigen::Vector3d &fx, Eigen::Vector3d &fy)
Compute the f0, fx, fy integration helper values from corner coordinates.
int fwd_bem_linear_collocation_solution()
Compute the linear-collocation BEM solution for this model.
int fwd_bem_load_recompute_solution(const QString &name, int bem_method, int force_recompute)
Load a BEM solution from file, recomputing if necessary.
static int fwd_bem_field(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > B, void *client)
Callback: compute BEM magnetic fields at coils for a dipole.
static int fwd_bem_field_grad(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &coils, Eigen::Ref< Eigen::VectorXf > Bval, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad, void *client)
Callback: compute BEM magnetic fields and position gradients at coils.
int compute_forward_eeg(std::vector< std::unique_ptr< MNELIB::MNESourceSpace > > &spaces, FwdCoilSet *els, bool fixed_ori, FwdEegSphereModel *eeg_model, bool use_threads, FIFFLIB::FiffNamedMatrix &resp, FIFFLIB::FiffNamedMatrix &resp_grad, bool bDoGrad)
Compute the EEG forward solution for one or more source spaces.
static float fwd_bem_inf_pot(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, const Eigen::Vector3f &rp)
Compute the infinite-medium electric potential at a single point.
static int fwd_bem_pot_els(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &els, Eigen::Ref< Eigen::VectorXf > pot, void *client)
Callback: compute BEM potentials at electrodes for a dipole.
static Eigen::MatrixXf fwd_bem_lin_pot_coeff(const std::vector< MNELIB::MNESurface * > &surfs)
Assemble the full linear-collocation potential coefficient matrix.
int fwd_bem_set_head_mri_t(const FIFFLIB::FiffCoordTrans &t)
Set the Head-to-MRI coordinate transform for this BEM model.
void(* linFieldIntFunc)(const Eigen::Vector3f &dest, const Eigen::Vector3f &dir, MNELIB::MNETriangle &tri, Eigen::Vector3d &res)
Function pointer type for linear field coefficient integration methods.
static float fwd_bem_inf_field_der(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, const Eigen::Vector3f &rp, const Eigen::Vector3f &dir, const Eigen::Vector3f &comp)
Compute the derivative of the infinite-medium magnetic field with respect to dipole position.
static FwdBemModel::UPtr fwd_bem_load_homog_surface(const QString &name)
Load a single-layer (homogeneous) BEM model from a FIFF file.
static float fwd_bem_inf_field(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, const Eigen::Vector3f &rp, const Eigen::Vector3f &dir)
Compute the infinite-medium magnetic field at a single point.
static float fwd_bem_inf_pot_der(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, const Eigen::Vector3f &rp, const Eigen::Vector3f &comp)
Compute the derivative of the infinite-medium electric potential with respect to dipole position.
static int fwd_sphere_field_vec(const Eigen::Vector3f &rd, FwdCoilSet &coils, Eigen::Ref< Eigen::MatrixXf > Bval, void *client)
Callback: compute the spherical-model vector magnetic field at coils.
void fwd_bem_lin_pot_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet *els, int all_surfs, Eigen::Ref< Eigen::VectorXf > pot)
Compute BEM potentials at electrodes using linear collocation.
FwdBemModel()
Constructs an empty BEM model.
static double one_field_coeff(const Eigen::Vector3f &dest, const Eigen::Vector3f &normal, MNELIB::MNETriangle &tri)
Compute the constant-collocation magnetic field coefficient for one triangle.
static const QString & fwd_bem_explain_surface(int kind)
Return a human-readable label for a BEM surface kind.
static void meg_eeg_fwd_one_source_space(FwdThreadArg *arg)
Thread worker: compute the forward solution for one source space.
int compute_forward_meg(std::vector< std::unique_ptr< MNELIB::MNESourceSpace > > &spaces, FwdCoilSet *coils, FwdCoilSet *comp_coils, MNELIB::MNECTFCompDataSet *comp_data, bool fixed_ori, const Eigen::Vector3f &r0, bool use_threads, FIFFLIB::FiffNamedMatrix &resp, FIFFLIB::FiffNamedMatrix &resp_grad, bool bDoGRad)
Compute the MEG forward solution for one or more source spaces.
static void fwd_bem_ip_modify_solution(Eigen::MatrixXf &solution, Eigen::MatrixXf &ip_solution, float ip_mult, int nsurf, const Eigen::VectorXi &ntri)
Modify the BEM solution with the isolated-problem (IP) approach.
static std::unique_ptr< MNELIB::MNESurface > make_guesses(MNELIB::MNESurface *guess_surf, float guessrad, const Eigen::Vector3f &guess_r0, float grid, float exclude, float mindist)
Generate a set of dipole guess locations inside a boundary surface.
int fwd_bem_constant_collocation_solution()
Compute the constant-collocation BEM solution for this model.
Eigen::MatrixXf solution
void fwd_bem_lin_pot_grad_calc(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet *els, int all_surfs, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad)
Compute the gradient of BEM potentials with respect to dipole position (linear collocation).
static int fwd_bem_pot_grad_els(const Eigen::Vector3f &rd, const Eigen::Vector3f &Q, FwdCoilSet &els, Eigen::Ref< Eigen::VectorXf > pot, Eigen::Ref< Eigen::VectorXf > xgrad, Eigen::Ref< Eigen::VectorXf > ygrad, Eigen::Ref< Eigen::VectorXf > zgrad, void *client)
Callback: compute BEM potentials and position gradients at electrodes.
Eigen::MatrixXf fwd_bem_lin_field_coeff(FwdCoilSet *coils, int method)
Assemble the linear-collocation magnetic field coefficient matrix.
Eigen::VectorXf sigma
static void calc_magic(double u, double z, double A, double B, Eigen::Vector3d &beta, double &D)
Compute the "magic" beta and D factors for the Urankar field integration.
Eigen::VectorXf v0
FIFFLIB::FiffCoordTrans head_mri_t
Collection of FwdCoil objects representing a full MEG or EEG sensor array.
Multi-layer spherical head model for EEG forward computation.
Thread-local arguments for parallel forward field computation (source range, coils,...
Collection of CTF third-order gradient compensation operators.
This defines a surface.
Definition mne_surface.h:79
Per-triangle geometric data for a cortical or BEM surface.