46#include <QtAlgorithms>
61#define MIN_1(a,b) ((a) < (b) ? (a) : (b))
70#define VEC_DOT_1(x,y) ((x)[X_1]*(y)[X_1] + (x)[Y_1]*(y)[Y_1] + (x)[Z_1]*(y)[Z_1])
71#define VEC_LEN_1(x) sqrt(VEC_DOT_1(x,x))
75#define VEC_DIFF_1(from,to,diff) {\
76 (diff)[X_1] = (to)[X_1] - (from)[X_1];\
77 (diff)[Y_1] = (to)[Y_1] - (from)[Y_1];\
78 (diff)[Z_1] = (to)[Z_1] - (from)[Z_1];\
81#define VEC_COPY_1(to,from) {\
82 (to)[X_1] = (from)[X_1];\
83 (to)[Y_1] = (from)[Y_1];\
84 (to)[Z_1] = (from)[Z_1];\
87#define CROSS_PRODUCT_1(x,y,xy) {\
88 (xy)[X_1] = (x)[Y_1]*(y)[Z_1]-(y)[Y_1]*(x)[Z_1];\
89 (xy)[Y_1] = -((x)[X_1]*(y)[Z_1]-(y)[X_1]*(x)[Z_1]);\
90 (xy)[Z_1] = (x)[X_1]*(y)[Y_1]-(y)[X_1]*(x)[Y_1];\
93#define MALLOC_1(x,t) (t *)malloc((x)*sizeof(t))
95#define REALLOC_1(x,y,t) (t *)((x == NULL) ? malloc((y)*sizeof(t)) : realloc((x),(y)*sizeof(t)))
97#define FREE(x) if ((char *)(x) != NULL) free((char *)(x))
99#define ALLOC_DCMATRIX_1(x,y) mne_dmatrix_1((x),(y))
100#define FREE_DCMATRIX_1(m) mne_free_dcmatrix_1((m))
102#define ALLOC_CMATRIX_1(x,y) mne_cmatrix_1((x),(y))
107#define FREE_CMATRIX_1(m) mne_free_cmatrix_1((m))
125static void matrix_error_1(
int kind,
int nr,
int nc)
129 printf(
"Failed to allocate memory pointers for a %d x %d matrix\n",nr,nc);
131 printf(
"Failed to allocate memory for a %d x %d matrix\n",nr,nc);
133 printf(
"Allocation error for a %d x %d matrix\n",nr,nc);
134 if (
sizeof(
void *) == 4) {
135 printf(
"This is probably because you seem to be using a computer with 32-bit architecture.\n");
136 printf(
"Please consider moving to a 64-bit platform.");
138 printf(
"Cannot continue. Sorry.\n");
150 if (!m) matrix_error_1(1,nr,nc);
152 if (!whole) matrix_error_1(2,nr,nc);
167 if (!m) matrix_error_1(1,nr,nc);
169 if (!whole) matrix_error_1(2,nr,nc);
223 if (!p_FwdEegSphereModel.
name.isEmpty())
224 this->
name = p_FwdEegSphereModel.
name;
225 if (p_FwdEegSphereModel.
nlayer() > 0) {
226 for (k = 0; k < p_FwdEegSphereModel.
nlayer(); k++)
230 if (p_FwdEegSphereModel.
nterms > 0) {
231 this->
fn = VectorXd(p_FwdEegSphereModel.
nterms);
233 for (k = 0; k < p_FwdEegSphereModel.
nterms; k++)
234 this->
fn[k] = p_FwdEegSphereModel.
fn[k];
236 if (p_FwdEegSphereModel.
nfit > 0) {
237 this->
mu = VectorXf(p_FwdEegSphereModel.
nfit);
238 this->
lambda = VectorXf(p_FwdEegSphereModel.
nfit);
239 this->
nfit = p_FwdEegSphereModel.
nfit;
240 for (k = 0; k < p_FwdEegSphereModel.
nfit; k++) {
241 this->
mu[k] = p_FwdEegSphereModel.
mu[k];
252 const VectorXf& rads,
253 const VectorXf& sigmas)
262 for (
int k = 0; k <
nlayer; k++) {
265 layer.
sigma = sigmas[k];
266 new_model->
layers.append(layer);
278 for (
int k = 0; k <
nlayer; k++) {
280 new_model->
layers[k].rel_rad = new_model->
layers[k].rel_rad/rR;
298 if (eeg_model_name.isEmpty())
299 eeg_model_name = QString(
"Default");
309 printf(
"Using EEG sphere model \"%s\" with scalp radius %7.1f mm\n",
310 eeg_model->
name.toUtf8().constData(),1000*eeg_sphere_rad);
345 MatrixXd M,Mn,help,Mm;
346 static MatrixXd mat1;
347 static MatrixXd mat2;
348 static MatrixXd mat3;
352 static VectorXd cr_mult;
367 if (this->
nlayer() == 2) {
370 div_mult = 2.0*n + 1;
371 b = pow(this->
layers[0].rel_rad,div_mult);
372 return div_mult/((n1 + n*rel1) + b*n1*(rel1-1.0));
374 else if (this->
nlayer() == 3) {
378 div_mult = 2.0*n + 1.0;
379 b = pow(this->
layers[0].rel_rad,div_mult);
380 c = pow(this->
layers[1].rel_rad,div_mult);
381 div_mult = div_mult*div_mult;
382 div = (b*n*n1*(rel1-1.0)*(rel2-1.0) + c*(rel1*n + n1)*(rel2*n + n1))/c +
383 n1*(b*(rel1-1.0)*(rel2*n1 + n) + c*(rel1*n + n1)*(rel2-1.0));
391 c1.resize(this->
nlayer()-1);
392 c2.resize(this->
nlayer()-1);
393 cr.resize(this->
nlayer()-1);
394 cr_mult.resize(this->
nlayer()-1);
395 for (k = 0; k < this->
nlayer()-1; k++) {
398 cr_mult[k] = this->
layers[k].rel_rad;
400 cr_mult[k] = cr_mult[k]*cr_mult[k];
402 if (mat1.cols() == 0)
403 mat1 = MatrixXd(2,2);
404 if (mat2.cols() == 0)
405 mat2 = MatrixXd(2,2);
406 if (mat3.cols() == 0)
407 mat3 = MatrixXd(2,2);
412 for (k = 0; k < this->
nlayer()-1; k++)
413 cr[k] = cr[k]*cr_mult[k];
420 M(0,0) = M(1,1) = 1.0;
421 M(0,1) = M(1,0) = 0.0;
423 div_mult = 2.0*n + 1.0;
426 for (k = this->
nlayer()-2; k >= 0; k--) {
428 Mm(0,0) = (n + n1*c1[k]);
429 Mm(0,1) = n1*c2[k]/cr[k];
430 Mm(1,0) = n*c2[k]*cr[k];
431 Mm(1,1) = n1 + n*c1[k];
433 Mn(0,0) = Mm(0,0)*M(0,0) + Mm(0,1)*M(1,0);
434 Mn(0,1) = Mm(0,0)*M(0,1) + Mm(0,1)*M(1,1);
435 Mn(1,0) = Mm(1,0)*M(0,0) + Mm(1,1)*M(1,0);
436 Mn(1,1) = Mm(1,0)*M(0,1) + Mm(1,1)*M(1,1);
443 return n*div/(n*M(1,1) + n1*M(1,0));
462 *p0 = ((2*n-1)*x*help0 - (n-1)*(*p01))/n;
463 *p1 = ((2*n-1)*x*help1 - n*(*p11))/(n-1);
486 double p0,p01,p1,p11;
491 p0 = p01 = p1 = p11 = 0.0;
492 for (n = 1; n <=
nterms; n++) {
499 multn = betan*
fn[n-1];
501 Vt = Vt + multn*p1/n;
534 float my_rd[3],pos[3];
536 float pos2,rd_len,pos_len;
537 double beta,cos_gamma,Vr,Vt;
538 float vec1[3],vec2[3],v1,v2;
539 float cos_beta,Qr,Qt,Q2,c;
540 float pi4_inv = 0.25/
M_PI;
554 for (p = 0; p < 3; p++)
555 my_rd[p] = rd[p] - m->
r0[p];
562 if (rd_len >= m->
layers[0].rad) {
563 for (k = 0; k < neeg; k++)
577 vec1[0] = vec1[1] = vec1[2] = 0.0;
585 for (k = 0; k < neeg; k++) {
586 for (p = 0; p < 3; p++)
587 pos[p] = el[k][p] - m->
r0[p];
594 printf(
"%10.4f %10.4f %10.4f %10.4f\n",pos_len,1000*pos[0],1000*pos[1],1000*pos[2]);
596 for (p = 0; p < 3; p++)
597 pos[p] = pos_len*pos[p];
600 pos_len = sqrt(pos2);
604 cos_gamma =
VEC_DOT_1(pos,rd)/(rd_len*pos_len);
605 beta = rd_len/pos_len;
620 Qt = sqrt(Q2 - Qr*Qr);
622 Vval[k] = pi4_inv*(Qr*Vr + Qt*cos_beta*Vt)/pos2;
630 for (k = 0; k < neeg; k++)
631 Vval[k] = Vval[k]*sigmaM_inv;
647 float *vval_one = NULL,val;
652 for (k = 0; k < els->
ncoil; k++, el++) {
655 if (el->
np > nvval) {
663 for (c = 0, val = 0.0; c < el->
np; c++)
664 val += el->
w[c]*vval_one[c];
678 float fact = 0.25f/(float)
M_PI;
681 float rrd,rd2,rd2_inv,r,r2,ra,rda;
686 float orig_rd[3],scaled_rd[3];
687 float pos[3],pos_len;
691 for (p = 0; p < 3; p++)
692 orig_rd[p] = rd[p] - m->
r0[p];
697 for (k = 0 ; k < neeg ; k++) {
698 Vval_vec[
X_1][k] = 0.0;
699 Vval_vec[
Y_1][k] = 0.0;
700 Vval_vec[
Z_1][k] = 0.0;
712 eeg_set_homog_sphere_model();
717 for (eq = 0; eq < m->
nfit; eq++) {
721 for (p = 0; p < 3; p++)
722 rd[p] = m->
mu[eq]*orig_rd[p];
730 for (k = 0; k < neeg ; k++) {
733 for (p = 0; p < 3; p++)
734 pos[p] = this_pos[p] - m->
r0[p];
740 for (p = 0; p < 3; p++)
741 pos[p] = pos_len*pos[p];
751 a2 =
VEC_DOT_1(a_vec,a_vec); a = sqrt(a2);
753 r2 =
VEC_DOT_1(this_pos,this_pos); r = sqrt(r2);
761 c1 = a3*rda + 1.0/a - 1.0/r;
762 c2 = a3 + (a+r)/(r*F);
769 Vval_vec[
X_1][k] = Vval_vec[
X_1][k] + m->
lambda[eq]*rd2_inv*(m1*rd[
X_1] + m2*this_pos[
X_1]);
770 Vval_vec[
Y_1][k] = Vval_vec[
Y_1][k] + m->
lambda[eq]*rd2_inv*(m1*rd[
Y_1] + m2*this_pos[
Y_1]);
771 Vval_vec[
Z_1][k] = Vval_vec[
Z_1][k] + m->
lambda[eq]*rd2_inv*(m1*rd[
Z_1] + m2*this_pos[
Z_1]);
777 for (k = 0; k < neeg; k++) {
778 Vval_vec[
X_1][k] = fact*Vval_vec[
X_1][k];
779 Vval_vec[
Y_1][k] = fact*Vval_vec[
Y_1][k];
780 Vval_vec[
Z_1][k] = fact*Vval_vec[
Z_1][k];
789 float **vval_one = NULL;
795 for (k = 0; k < els->
ncoil; k++, el++) {
798 if (el->
np > nvval) {
807 for (p = 0; p < 3; p++) {
808 for (c = 0, val = 0.0; c < el->
np; c++)
809 val += el->
w[c]*vval_one[p][c];
810 Vval_vec[p][k] = val;
831 float step2 = 2*step;
839 for (p = 0; p < 3; p++) {
841 my_rd[p] = my_rd[p] + step;
845 my_rd[p] = my_rd[p] - step;
848 for (q = 0; q < coils->
ncoil; q++)
849 grads[p][q] = (grads[p][q]-Vval[q])/step2;
874 float fact = 0.25f/
M_PI;
877 float rrd,rd2,rd2_inv,r,r2,ra,rda;
879 float c1,c2,m1,m2,f1,f2;
882 float orig_rd[3],scaled_rd[3];
883 float pos[3],pos_len;
887 for (p = 0; p < 3; p++)
888 orig_rd[p] = rd[p] - m->
r0[p];
893 for (k = 0 ; k < neeg ; k++)
905 eeg_set_homog_sphere_model();
910 for (eq = 0; eq < m->
nfit; eq++) {
914 for (p = 0; p < 3; p++)
915 rd[p] = m->
mu[eq]*orig_rd[p];
924 for (k = 0; k < neeg ; k++) {
927 for (p = 0; p < 3; p++)
928 pos[p] = this_pos[p] - m->
r0[p];
934 for (p = 0; p < 3; p++)
935 pos[p] = pos_len*pos[p];
945 a2 =
VEC_DOT_1(a_vec,a_vec); a = sqrt(a2);
947 r2 =
VEC_DOT_1(this_pos,this_pos); r = sqrt(r2);
955 c1 = a3*rda + 1.0/a - 1.0/r;
956 c2 = a3 + (a+r)/(r*F);
964 Vval[k] = Vval[k] + m->
lambda[eq]*rd2_inv*(m1*f1 + m2*f2);
970 for (k = 0; k < neeg; k++)
971 Vval[k] = fact*Vval[k];
985 for (k = 0; k < els->
ncoil; k++, el++) {
988 if (el->
np > nvval) {
989 vval_one.resize(el->
np);
995 for (c = 0, val = 0.0; c < el->
np; c++)
996 val += el->
w[c]*vval_one[c];
1014 for (
int k = 0; k < this->
nlayer(); k++)
1017 if (fit_berg_scherg) {
1019 printf(
"Equiv. model fitting -> ");
1020 printf(
"RV = %g %%\n",100*rv);
1021 for (
int k = 0; k <
nfit; k++)
1022 printf(
"mu%d = %g\tlambda%d = %g\n", k+1,this->
mu[k],k+1,this->
layers[this->
nlayer()-1].sigma*this->
lambda[k]);
1028 printf(
"Defined EEG sphere model with rad = %7.2f mm\n", 1000.0*rad);
1034 Eigen::MatrixXd eigen_mat(m,n);
1036 for (
int i = 0; i < m; ++i)
1037 for (
int j = 0; j < n; ++j)
1038 eigen_mat(i,j) = mat[i][j];
1045 for (
int i = 0; i < m; ++i)
1046 for (
int j = 0; j < n; ++j)
1047 to_mat[i][j] = from_mat(i,j);
1057 for (
int i = 0; i < n; ++i)
1058 to_vec[i] = from_vec[i];
1067static double dot_dvectors (
double *v1,
1071 double result = 0.0;
1074 for (k = 0; k < nn; k++)
1075 result = result + v1[k]*v2[k];
1079static int c_dsvd(
double **mat,
1102 int udim =
MIN_1(m,n);
1107 Eigen::JacobiSVD< Eigen::MatrixXd > svd(eigen_mat,Eigen::ComputeFullU | Eigen::ComputeFullV);
1143static int comp_pars(
const void *p1,
const void *p2)
1151 if (v1->
mu > v2->
mu)
1153 else if (v1->
mu < v2->
mu)
1160static void sort_parameters(VectorXd& mu,VectorXd& lambda,
int nfit)
1167 for (
int k = 0; k < nfit; k++) {
1169 pars[k].
lambda = lambda[k];
1174 for (
int k = 0; k < nfit; k++) {
1176 lambda[k] = pars[k].
lambda;
1183static bool report_fit(
int loop,
1184 const VectorXd &fitpar,
1192 for (
int k = 0; k < fitpar.size(); k++)
1193 printf(
"%g ",mu[k]);
1194 printf(
"%g\n",Smin);
1200static MatrixXd get_initial_simplex(
const VectorXd &pars,
1201 double simplex_size)
1204 int npar = pars.size();
1206 MatrixXd simplex = MatrixXd::Zero(npar+1,npar);
1208 simplex.rowwise() += pars.transpose();
1210 for (
int k = 1; k < npar+1; k++)
1211 simplex(k,k-1) += simplex_size;
1224 for (k = 0; k < u->
nterms-1; k++) {
1226 mu1n = pow(
mu[0],k1);
1227 u->
y[k] = u->
w[k]*(u->
fn[k+1] - mu1n*u->
fn[0]);
1228 for (p = 0; p < u->
nfit-1; p++)
1229 u->
M[k][p] = u->
w[k]*(pow(
mu[p+1],k1)-mu1n);
1243 VectorXd vec(u->
nfit-1);
1252 for (k = 0; k < u->
nterms-1; k++)
1253 u->
resi[k] = u->
y[k];
1255 for (p = 0; p < u->
nfit-1; p++) {
1256 vec[p] = dot_dvectors(u->
uu[p],u->
y,u->
nterms-1);
1257 for (k = 0; k < u->
nterms-1; k++)
1258 u->
resi[k] = u->
resi[k] - u->
uu[p][k]*vec[p];
1259 vec[p] = vec[p]/u->
sing[p];
1262 for (p = 0; p < u->
nfit-1; p++) {
1263 for (q = 0, sum = 0.0; q < u->
nfit-1; q++)
1264 sum += u->
vv[q][p]*vec[q];
1267 for (p = 1, sum = 0.0; p < u->
nfit; p++)
1284 for (k = 0; k < u->
nfit; k++) {
1285 if (std::fabs(
mu[k]) > 1.0)
1299 for (k = 0; k < u->
nterms-1; k++)
1300 u->
resi[k] = u->
y[k];
1301 for (p = 0; p < u->
nfit-1; p++) {
1302 dot = dot_dvectors(u->
uu[p],u->
y,u->
nterms-1);
1303 for (k = 0; k < u->
nterms-1; k++)
1326 double simplex_size = 0.01;
1333 int max_eval = 1000;
1338 printf(
"fwd_fit_berg_scherg does not work with less than two equivalent sources.");
1345 for (k = 0; k <
nterms; k++)
1346 u->
fn[k] = this->fwd_eeg_get_multi_sphere_model_coeff(k+1);
1351 rd = R = this->
layers[0].rad;
1352 for (k = 1; k < this->
nlayer(); k++) {
1353 if (this->
layers[k].rad > R)
1355 if (this->
layers[k].rad < rd)
1356 rd = this->
layers[k].rad;
1364 for (k = 1; k <
nterms; k++)
1365 u->
w[k-1] = pow(f,k);
1370 for (k = 1; k <
nterms; k++)
1371 u->
w[k-1] = sqrt((2.0*k+1)*(3.0*k+1.0)/k)*pow(f,(k-1.0));
1377 func_val = VectorXd(
nfit+1);
1383 for (k = 0; k <
nfit; k++) {
1387 mu[k] = (rand() / (RAND_MAX + 1.0))*f;
1390 simplex = get_initial_simplex(
mu,simplex_size);
1391 for (k = 0; k <
nfit+1; k++)
1392 func_val[k] =
one_step(
static_cast<VectorXd
>(simplex.row(k)),u);
1408 for (k = 0; k <
nfit; k++)
1409 mu[k] = simplex(0,k);
1418 printf(
"RV = %g %%\n",100*rv);
1420 this->mu.resize(
nfit);
1421 this->lambda.resize(
nfit);
1423 for (k = 0; k <
nfit; k++) {
1424 this->mu[k] =
mu[k];
1430 printf(
"lambda%d = %g\tmu%d = %g\n",k+1,
lambda[k],k+1,
mu[k]);
SimplexAlgorithm Template Implementation.
FwdEegSphereModelSet class declaration.
void fromDoubleEigenMatrix(const Eigen::MatrixXd &from_mat, double **to_mat, const int m, const int n)
#define ALLOC_CMATRIX_1(x, y)
void mne_free_cmatrix_1(float **m)
#define VEC_COPY_1(to, from)
#define CROSS_PRODUCT_1(x, y, xy)
#define FREE_DCMATRIX_1(m)
void mne_free_dcmatrix_1(double **m)
#define REALLOC_1(x, y, t)
Eigen::MatrixXd toDoubleEigenMatrix(double **mat, const int m, const int n)
double ** mne_dmatrix_1(int nr, int nc)
void fromDoubleEigenVector(const Eigen::VectorXd &from_vec, double *to_vec, const int n)
#define ALLOC_DCMATRIX_1(x, y)
#define FREE_CMATRIX_1(m)
float ** mne_cmatrix_1(int nr, int nc)
#define VEC_DIFF_1(from, to, diff)
FwdEegSphereModel class declaration.
Forward modelling (BEM, MEG/EEG lead fields).
struct FWDLIB::bergSchergPar bergSchergParRec
struct FWDLIB::fitUser fitUserRec
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Single MEG or EEG sensor coil with integration points, weights, and coordinate frame.
Collection of FwdCoil objects representing a full MEG or EEG sensor array.
Single concentric sphere layer with conductivity ratio for the EEG forward model.
static bool comp_layers(const FwdEegSphereLayer &v1, const FwdEegSphereLayer &v2)
Berg-Scherg parameter pair (magnitude and distance multiplier) for an equivalent dipole in the EEG sp...
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)
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)
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)
virtual ~FwdEegSphereModel()
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)
Collection of FwdEegSphereModel objects for multi-model EEG forward solutions.
FwdEegSphereModel * fwd_select_eeg_sphere_model(const QString &p_sName)
void fwd_list_eeg_sphere_models(FILE *f)
static FwdEegSphereModelSet * fwd_load_eeg_sphere_models(const QString &p_sFileName, FwdEegSphereModelSet *now)
static bool simplex_minimize(Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &p, Eigen::Matrix< T, Eigen::Dynamic, 1 > &y, T ftol, T(*func)(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, const void *user_data), const void *user_data, int max_eval, int &neval, int report, bool(*report_func)(int loop, const Eigen::Matrix< T, Eigen::Dynamic, 1 > &fitpar, double fval))