67constexpr int FAIL = -1;
69constexpr int MAXTERMS = 1000;
70constexpr double EPS = 1e-10;
71constexpr double SIN_EPS = 1e-3;
92 if (!p_FwdEegSphereModel.
name.isEmpty())
93 this->
name = p_FwdEegSphereModel.
name;
94 if (p_FwdEegSphereModel.
nlayer() > 0) {
95 for (k = 0; k < p_FwdEegSphereModel.
nlayer(); k++)
98 this->
r0 = p_FwdEegSphereModel.
r0;
99 if (p_FwdEegSphereModel.
nterms > 0) {
100 this->
fn = VectorXd(p_FwdEegSphereModel.
nterms);
102 for (k = 0; k < p_FwdEegSphereModel.
nterms; k++)
103 this->
fn[k] = p_FwdEegSphereModel.
fn[k];
105 if (p_FwdEegSphereModel.
nfit > 0) {
106 this->
mu = VectorXf(p_FwdEegSphereModel.
nfit);
107 this->
lambda = VectorXf(p_FwdEegSphereModel.
nfit);
108 this->
nfit = p_FwdEegSphereModel.
nfit;
109 for (k = 0; k < p_FwdEegSphereModel.
nfit; k++) {
110 this->
mu[k] = p_FwdEegSphereModel.
mu[k];
121 const VectorXf& rads,
122 const VectorXf& sigmas)
127 auto new_model = std::make_unique<FwdEegSphereModel>();
129 new_model->name =
name;
131 for (
int k = 0; k <
nlayer; k++) {
134 layer.
sigma = sigmas[k];
135 new_model->layers.push_back(layer);
145 float R = new_model->layers[
nlayer-1].rad;
146 float rR = new_model->layers[
nlayer-1].rel_rad;
147 for (
int k = 0; k <
nlayer; k++) {
148 new_model->layers[k].rad = new_model->layers[k].rad/R;
149 new_model->layers[k].rel_rad = new_model->layers[k].rel_rad/rR;
164 if (eeg_model_name.isEmpty())
165 eeg_model_name = QString(
"Default");
168 eeg_models->fwd_list_eeg_sphere_models();
175 if (!eeg_model->fwd_setup_eeg_sphere_model(eeg_sphere_rad,
true,3)) {
179 qInfo(
"Using EEG sphere model \"%s\" with scalp radius %7.1f mm",
180 eeg_model->name.toUtf8().constData(),1000*eeg_sphere_rad);
207 MatrixXd M,Mn,help,Mm;
208 static MatrixXd mat1;
209 static MatrixXd mat2;
210 static MatrixXd mat3;
214 static VectorXd cr_mult;
229 if (this->
nlayer() == 2) {
232 div_mult = 2.0*n + 1;
233 b = pow(this->
layers[0].rel_rad,div_mult);
234 return div_mult/((n1 + n*rel1) + b*n1*(rel1-1.0));
236 else if (this->
nlayer() == 3) {
240 div_mult = 2.0*n + 1.0;
241 b = pow(this->
layers[0].rel_rad,div_mult);
242 c = pow(this->
layers[1].rel_rad,div_mult);
243 div_mult = div_mult*div_mult;
244 div = (b*n*n1*(rel1-1.0)*(rel2-1.0) + c*(rel1*n + n1)*(rel2*n + n1))/c +
245 n1*(b*(rel1-1.0)*(rel2*n1 + n) + c*(rel1*n + n1)*(rel2-1.0));
253 c1.resize(this->
nlayer()-1);
254 c2.resize(this->
nlayer()-1);
255 cr.resize(this->
nlayer()-1);
256 cr_mult.resize(this->
nlayer()-1);
257 for (k = 0; k < this->
nlayer()-1; k++) {
260 cr_mult[k] = this->
layers[k].rel_rad;
262 cr_mult[k] = cr_mult[k]*cr_mult[k];
264 if (mat1.cols() == 0)
265 mat1 = MatrixXd(2,2);
266 if (mat2.cols() == 0)
267 mat2 = MatrixXd(2,2);
268 if (mat3.cols() == 0)
269 mat3 = MatrixXd(2,2);
274 for (k = 0; k < this->
nlayer()-1; k++)
275 cr[k] = cr[k]*cr_mult[k];
282 M(0,0) = M(1,1) = 1.0;
283 M(0,1) = M(1,0) = 0.0;
285 div_mult = 2.0*n + 1.0;
288 for (k = this->
nlayer()-2; k >= 0; k--) {
290 Mm(0,0) = (n + n1*c1[k]);
291 Mm(0,1) = n1*c2[k]/cr[k];
292 Mm(1,0) = n*c2[k]*cr[k];
293 Mm(1,1) = n1 + n*c1[k];
295 Mn(0,0) = Mm(0,0)*M(0,0) + Mm(0,1)*M(1,0);
296 Mn(0,1) = Mm(0,0)*M(0,1) + Mm(0,1)*M(1,1);
297 Mn(1,0) = Mm(1,0)*M(0,0) + Mm(1,1)*M(1,0);
298 Mn(1,1) = Mm(1,0)*M(0,1) + Mm(1,1)*M(1,1);
305 return n*div/(n*M(1,1) + n1*M(1,0));
324 p0 = ((2*n-1)*x*help0 - (n-1)*(p01))/n;
325 p1 = ((2*n-1)*x*help1 - n*(p11))/(n-1);
348 double p0,p01,p1,p11;
353 p0 = p01 = p1 = p11 = 0.0;
354 for (n = 1; n <=
nterms; n++) {
359 multn = betan*
fn[n-1];
361 Vt = Vt + multn*p1/n;
394 Eigen::Vector3f rd = rd_in - m->r0;
395 Eigen::Vector3f Q = Q_in;
398 float pos2,rd_len,pos_len;
399 double beta,cos_gamma,Vr,Vt;
400 Eigen::Vector3f vec1, vec2;
402 float cos_beta,Qr,Qt,Q2,c;
403 float pi4_inv = 0.25/
M_PI;
408 if (m->fn.size() == 0 || m->nterms != MAXTERMS) {
409 m->fn.resize(MAXTERMS);
410 m->nterms = MAXTERMS;
411 for (k = 0; k < MAXTERMS; k++)
412 m->fn[k] = (2*k+3)*m->fwd_eeg_get_multi_sphere_model_coeff(k+1);
422 if (rd_len >= m->layers[0].rad) {
423 for (k = 0; k < neeg; k++)
430 c = rd.dot(Q)/(rd_len*sqrt(Q2));
431 if ((1.0-c*c) < SIN_EPS) {
445 for (k = 0; k < neeg; k++) {
446 pos = el.row(k).transpose() - m->r0;
451 pos_len = m->layers[m->nlayer()-1].rad/pos.norm();
453 qInfo(
"%10.4f %10.4f %10.4f %10.4f",pos_len,1000*pos[0],1000*pos[1],1000*pos[2]);
458 pos_len = sqrt(pos2);
462 cos_gamma = pos.dot(rd)/(rd_len*pos_len);
463 beta = rd_len/pos_len;
469 vec2 = rd.cross(pos);
473 cos_beta = vec1.dot(vec2)/(v1*v2);
477 Qr = Q.dot(rd)/rd_len;
478 Qt = sqrt(Q2 - Qr*Qr);
480 Vval[k] = pi4_inv*(Qr*Vr + Qt*cos_beta*Vt)/pos2;
486 if (m->nlayer() > 0) {
487 sigmaM_inv = 1.0/m->layers[m->nlayer()-1].sigma;
488 for (k = 0; k < neeg; k++)
489 Vval[k] = Vval[k]*sigmaM_inv;
511 Vval.resize(els.
ncoil());
512 for (k = 0; k < els.
ncoil(); k++, el++) {
513 el = els.
coils[k].get();
515 if (el->
np > nvval) {
516 vval_one.resize(el->
np);
522 for (c = 0, val = 0.0; c < el->
np; c++)
523 val += el->
w[c]*vval_one[c];
535 float fact = 0.25f /
static_cast<float>(
M_PI);
536 Eigen::Vector3f a_vec;
538 float rrd,rd2,rd2_inv,r,r2,ra,rda;
542 Eigen::Vector3f orig_rd = rd_in - m->r0;
549 for (k = 0 ; k < neeg ; k++) {
557 if (orig_rd.norm() >= m->layers[0].rad)
564 eeg_set_homog_sphere_model();
569 for (eq = 0; eq < m->nfit; eq++) {
573 rd = m->mu[eq] * orig_rd;
581 for (k = 0; k < neeg ; k++) {
583 pos = el.row(k).transpose() - m->r0;
588 pos_len = m->layers[m->nlayer()-1].rad/pos.norm();
598 a2 = a_vec.dot(a_vec); a = sqrt(a2);
600 r2 = pos.dot(pos); r = sqrt(r2);
608 c1 = a3*rda + 1.0/a - 1.0/r;
609 c2 = a3 + (a+r)/(r*F);
616 Vval_vec(0,k) = Vval_vec(0,k) + m->lambda[eq]*rd2_inv*(m1*rd[0] + m2*pos[0]);
617 Vval_vec(1,k) = Vval_vec(1,k) + m->lambda[eq]*rd2_inv*(m1*rd[1] + m2*pos[1]);
618 Vval_vec(2,k) = Vval_vec(2,k) + m->lambda[eq]*rd2_inv*(m1*rd[2] + m2*pos[2]);
624 for (k = 0; k < neeg; k++) {
625 Vval_vec(0,k) = fact*Vval_vec(0,k);
626 Vval_vec(1,k) = fact*Vval_vec(1,k);
627 Vval_vec(2,k) = fact*Vval_vec(2,k);
636 Eigen::MatrixXf vval_one;
642 for (k = 0; k < els.
ncoil(); k++, el++) {
643 el = els.
coils[k].get();
645 if (el->
np > nvval) {
646 vval_one.resize(3, el->
np);
652 for (p = 0; p < 3; p++) {
653 for (c = 0, val = 0.0; c < el->
np; c++)
654 val += el->
w[c]*vval_one(p,c);
673 Eigen::Vector3f my_rd;
675 float step2 = 2*step;
678 Eigen::Ref<Eigen::VectorXf>* grads[3] = { &xgrad, &ygrad, &zgrad };
680 for (p = 0; p < 3; p++) {
689 for (q = 0; q < coils.
ncoil(); q++)
690 (*grads[p])[q] = ((*grads[p])[q]-Vval[q])/step2;
700 const Eigen::Vector3f& Q_in,
701 const Eigen::Matrix<float, Eigen::Dynamic, 3, Eigen::RowMajor>& el,
713 float fact = 0.25f/
M_PI;
714 Eigen::Vector3f a_vec;
716 float rrd,rd2,rd2_inv,r,r2,ra,rda;
718 float c1,c2,m1,m2,f1,f2;
720 Eigen::Vector3f orig_rd = rd_in - m->r0;
722 Eigen::Vector3f Q = Q_in;
728 for (k = 0 ; k < neeg ; k++)
733 if (orig_rd.norm() >= m->layers[0].rad)
740 eeg_set_homog_sphere_model();
745 for (eq = 0; eq < m->nfit; eq++) {
749 rd = m->mu[eq] * orig_rd;
758 for (k = 0; k < neeg ; k++) {
760 pos = el.row(k).transpose() - m->r0;
765 pos_len = m->layers[m->nlayer()-1].rad/pos.norm();
775 a2 = a_vec.dot(a_vec); a = sqrt(a2);
777 r2 = pos.dot(pos); r = sqrt(r2);
785 c1 = a3*rda + 1.0/a - 1.0/r;
786 c2 = a3 + (a+r)/(r*F);
794 Vval[k] = Vval[k] + m->lambda[eq]*rd2_inv*(m1*f1 + m2*f2);
800 for (k = 0; k < neeg; k++)
801 Vval[k] = fact*Vval[k];
815 for (k = 0; k < els.
ncoil(); k++, el++) {
816 el = els.
coils[k].get();
818 if (el->
np > nvval) {
819 vval_one.resize(el->
np);
825 for (c = 0, val = 0.0; c < el->
np; c++)
826 val += el->
w[c]*vval_one[c];
843 for (
int k = 0; k < this->
nlayer(); k++)
846 if (fit_berg_scherg) {
848 qInfo(
"Equiv. model fitting -> RV = %g %%",100*rv);
849 for (
int k = 0; k <
nfit; k++)
850 qInfo(
"mu%d = %g\tlambda%d = %g", k+1,this->
mu[k],k+1,this->
layers[this->
nlayer()-1].sigma*this->
lambda[k]);
856 qInfo(
"Defined EEG sphere model with rad = %7.2f mm", 1000.0*rad);
860static void compute_svd(Eigen::MatrixXd& mat,
861 Eigen::VectorXd& sing,
870 int udim = std::min(
static_cast<int>(mat.rows()),
static_cast<int>(mat.cols()));
872 Eigen::JacobiSVD<Eigen::MatrixXd> svd(mat, Eigen::ComputeFullU | Eigen::ComputeFullV);
874 sing = svd.singularValues();
875 uu = svd.matrixU().transpose().topRows(udim);
878 *vv = svd.matrixV().transpose();
900static void sort_parameters(VectorXd& mu,VectorXd& lambda,
int nfit)
902 std::vector<BergSchergPar> pars(nfit);
904 for (
int k = 0; k < nfit; k++) {
906 pars[k].lambda = lambda[k];
913 for (
int k = 0; k < nfit; k++) {
915 lambda[k] = pars[k].lambda;
919static bool report_fit(
int loop,
920 const VectorXd &fitpar,
926 for (
int k = 0; k < fitpar.size(); k++)
933static MatrixXd get_initial_simplex(
const VectorXd &pars,
937 int npar = pars.size();
939 MatrixXd simplex = MatrixXd::Zero(npar+1,npar);
941 simplex.rowwise() += pars.transpose();
943 for (
int k = 1; k < npar+1; k++)
944 simplex(k,k-1) += simplex_size;
957 for (k = 0; k < u->
nterms-1; k++) {
959 mu1n = pow(
mu[0],k1);
960 u->
y[k] = u->
w[k]*(u->
fn[k+1] - mu1n*u->
fn[0]);
961 for (p = 0; p < u->
nfit-1; p++)
962 u->
M(k,p) = u->
w[k]*(pow(
mu[p+1],k1)-mu1n);
976 VectorXd vec(u->
nfit-1);
981 compute_svd(u->
M, u->
sing, u->
uu, &u->
vv);
985 for (k = 0; k < u->
nterms-1; k++)
986 u->
resi[k] = u->
y[k];
988 for (p = 0; p < u->
nfit-1; p++) {
990 for (k = 0; k < u->
nterms-1; k++)
992 vec[p] = vec[p]/u->
sing[p];
995 for (p = 0; p < u->
nfit-1; p++) {
996 for (q = 0, sum = 0.0; q < u->
nfit-1; q++)
997 sum += u->
vv(q,p)*vec[q];
1000 for (p = 1, sum = 0.0; p < u->
nfit; p++)
1003 return u->
resi.head(u->
nterms-1).squaredNorm() / u->
y.head(u->
nterms-1).squaredNorm();
1017 for (k = 0; k < u->
nfit; k++) {
1018 if (std::fabs(
mu[k]) > 1.0)
1028 compute_svd(u->
M, u->
sing, u->
uu,
nullptr);
1032 for (k = 0; k < u->
nterms-1; k++)
1033 u->
resi[k] = u->
y[k];
1034 for (p = 0; p < u->
nfit-1; p++) {
1036 for (k = 0; k < u->
nterms-1; k++)
1042 return u->
resi.head(u->
nterms-1).squaredNorm();
1059 double simplex_size = 0.01;
1066 int max_eval = 1000;
1071 qWarning(
"fwd_fit_berg_scherg does not work with less than two equivalent sources.");
1078 for (k = 0; k <
nterms; k++)
1079 u->
fn[k] = this->fwd_eeg_get_multi_sphere_model_coeff(k+1);
1084 rd = R = this->
layers[0].rad;
1085 for (k = 1; k < this->
nlayer(); k++) {
1086 if (this->
layers[k].rad > R)
1088 if (this->
layers[k].rad < rd)
1089 rd = this->
layers[k].rad;
1097 for (k = 1; k <
nterms; k++)
1098 u->
w[k-1] = pow(f,k);
1103 for (k = 1; k <
nterms; k++)
1104 u->
w[k-1] = sqrt((2.0*k+1)*(3.0*k+1.0)/k)*pow(f,(k-1.0));
1110 func_val = VectorXd(
nfit+1);
1116 for (k = 0; k <
nfit; k++) {
1120 mu[k] = (rand() / (RAND_MAX + 1.0))*f;
1123 simplex = get_initial_simplex(
mu,simplex_size);
1124 for (k = 0; k <
nfit+1; k++)
1125 func_val[k] =
one_step(
static_cast<VectorXd
>(simplex.row(k)),u);
1142 for (k = 0; k <
nfit; k++)
1143 mu[k] = simplex(0,k);
1152 qInfo(
"RV = %g %%",100*rv);
1154 this->mu.resize(
nfit);
1155 this->lambda.resize(
nfit);
1157 for (k = 0; k <
nfit; k++) {
1158 this->mu[k] =
mu[k];
1164 qInfo(
"lambda%d = %g\tmu%d = %g",k+1,
lambda[k],k+1,
mu[k]);
SimplexAlgorithm Template Implementation.
FwdEegSphereModelSet class declaration.
FwdEegSphereModel class declaration.
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
struct UTILSLIB::fitUser fitUserRec
Forward modelling (BEM, MEG/EEG lead fields).
constexpr int FWD_COILC_EEG
Single MEG or EEG sensor coil with integration points, weights, and coordinate frame.
Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > rmag
Collection of FwdCoil objects representing a full MEG or EEG sensor array.
std::vector< FwdCoil::UPtr > coils
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...
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
virtual ~FwdEegSphereModel()
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)
static FwdEegSphereModelSet * fwd_load_eeg_sphere_models(const QString &p_sFileName, FwdEegSphereModelSet *now)
std::unique_ptr< FwdEegSphereModelSet > UPtr
static bool simplex_minimize(Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &p, Eigen::Matrix< T, Eigen::Dynamic, 1 > &y, T ftol, T stol, 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_lo, double fval_hi, double par_diff))