v2.0.0
Loading...
Searching...
No Matches
fwd_eeg_sphere_model.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
42
45
46#include <qmath.h>
47
48#include <Eigen/Core>
49
50#include <algorithm>
51#include <vector>
52#include <Eigen/Dense>
53
54//=============================================================================================================
55// USED NAMESPACES
56//=============================================================================================================
57
58using namespace Eigen;
59using namespace UTILSLIB;
60using namespace FWDLIB;
61
62//=============================================================================================================
63// Local constants and helpers
64//=============================================================================================================
65
66namespace {
67constexpr int FAIL = -1;
68constexpr int OK = 0;
69constexpr int MAXTERMS = 1000;
70constexpr double EPS = 1e-10;
71constexpr double SIN_EPS = 1e-3;
72} // anonymous namespace
73
74//=============================================================================================================
75// DEFINE MEMBER METHODS
76//=============================================================================================================
77
79: nterms (0)
80, nfit (0)
81, scale_pos (0)
82{
83 r0.setZero();
84}
85
86//=============================================================================================================
87
89{
90 int k;
91
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++)
96 this->layers.push_back(p_FwdEegSphereModel.layers[k]);
97 }
98 this->r0 = p_FwdEegSphereModel.r0;
99 if (p_FwdEegSphereModel.nterms > 0) {
100 this->fn = VectorXd(p_FwdEegSphereModel.nterms);
101 this->nterms = p_FwdEegSphereModel.nterms;
102 for (k = 0; k < p_FwdEegSphereModel.nterms; k++)
103 this->fn[k] = p_FwdEegSphereModel.fn[k];
104 }
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];
111 this->lambda[k] = p_FwdEegSphereModel.lambda[k];
112 }
113 }
114 this->scale_pos = p_FwdEegSphereModel.scale_pos;
115}
116
117//=============================================================================================================
118
120 int nlayer,
121 const VectorXf& rads,
122 const VectorXf& sigmas)
123/*
124 * Produce a new sphere model structure
125 */
126{
127 auto new_model = std::make_unique<FwdEegSphereModel>();
128
129 new_model->name = name;
130
131 for (int k = 0; k < nlayer; k++) {
132 FwdEegSphereLayer layer;
133 layer.rad = layer.rel_rad = rads[k];
134 layer.sigma = sigmas[k];
135 new_model->layers.push_back(layer);
136 }
137 /*
138 * Sort...
139 */
140 std::sort(new_model->layers.begin(), new_model->layers.end(), FwdEegSphereLayer::comp_layers);
141
142 /*
143 * Scale the radiuses
144 */
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;
150 }
151 return new_model;
152}
153
154//=============================================================================================================
155
159
160//=============================================================================================================
161
162FwdEegSphereModel::UPtr FwdEegSphereModel::setup_eeg_sphere_model(const QString& eeg_model_file, QString eeg_model_name, float eeg_sphere_rad)
163{
164 if (eeg_model_name.isEmpty())
165 eeg_model_name = QString("Default");
166
168 eeg_models->fwd_list_eeg_sphere_models();
169
170 FwdEegSphereModel::UPtr eeg_model(eeg_models->fwd_select_eeg_sphere_model(eeg_model_name));
171 if (!eeg_model) {
172 return nullptr;
173 }
174
175 if (!eeg_model->fwd_setup_eeg_sphere_model(eeg_sphere_rad,true,3)) {
176 return nullptr;
177 }
178
179 qInfo("Using EEG sphere model \"%s\" with scalp radius %7.1f mm",
180 eeg_model->name.toUtf8().constData(),1000*eeg_sphere_rad);
181 return eeg_model;
182}
183
184//=============================================================================================================
185
187
188{
189 fitUser u = new fitUserRec();
190 u->y.resize(nterms-1);
191 u->resi.resize(nterms-1);
192 u->M = MatrixXd::Zero(nterms-1,nfit-1);
193 u->uu = MatrixXd::Zero(nfit-1,nterms-1);
194 u->vv = MatrixXd::Zero(nfit-1,nfit-1);
195 u->sing.resize(nfit);
196 u->fn.resize(nterms);
197 u->w.resize(nterms);
198 u->nfit = nfit;
199 u->nterms = nterms;
200 return u;
201}
202
203//=============================================================================================================
204// fwd_multi_spherepot.c
206{
207 MatrixXd M,Mn,help,Mm;
208 static MatrixXd mat1;
209 static MatrixXd mat2;
210 static MatrixXd mat3;
211 static VectorXd c1;
212 static VectorXd c2;
213 static VectorXd cr;
214 static VectorXd cr_mult;
215 double div,div_mult;
216 double n1;
217#ifdef TEST
218 double rel1,rel2;
219 double b,c;
220#endif
221 int k;
222
223 if (this->nlayer() == 0 || this->nlayer() == 1)
224 return 1.0;
225 /*
226 * Now follows the tricky case
227 */
228#ifdef TEST
229 if (this->nlayer() == 2) {
230 rel1 = layers[0].sigma/layers[1].sigma;
231 n1 = n + 1.0;
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));
235 }
236 else if (this->nlayer() == 3) {
237 rel1 = this->layers[0].sigma/this->layers[1].sigma;
238 rel2 = this->layers[1].sigma/this->layers[2].sigma;
239 n1 = n + 1.0;
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));
246 return div_mult/div;
247 }
248#endif
249 if (n == 1) {
250 /*
251 * Initialize the arrays
252 */
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++) {
258 c1[k] = this->layers[k].sigma/this->layers[k+1].sigma;
259 c2[k] = c1[k] - 1.0;
260 cr_mult[k] = this->layers[k].rel_rad;
261 cr[k] = cr_mult[k];
262 cr_mult[k] = cr_mult[k]*cr_mult[k];
263 }
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);
270 }
271 /*
272 * Increment the radius coefficients
273 */
274 for (k = 0; k < this->nlayer()-1; k++)
275 cr[k] = cr[k]*cr_mult[k];
276 /*
277 * Multiply the matrices
278 */
279 M = mat1;
280 Mn = mat2;
281 Mm = mat3;
282 M(0,0) = M(1,1) = 1.0;
283 M(0,1) = M(1,0) = 0.0;
284 div = 1.0;
285 div_mult = 2.0*n + 1.0;
286 n1 = n + 1.0;
287
288 for (k = this->nlayer()-2; k >= 0; k--) {
289
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];
294
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);
299 help = M;
300 M = Mn;
301 Mn = help;
302 div = div*div_mult;
303
304 }
305 return n*div/(n*M(1,1) + n1*M(1,0));
306}
307
308//=============================================================================================================
309// fwd_multi_spherepot.c
310void FwdEegSphereModel::next_legen(int n, double x, double &p0, double &p01, double &p1, double &p11)
311/*
312 * Compute the next Legendre polynomials of the
313 * first and second kind using the recursion formulas.
314 *
315 * The routine initializes automatically with the known values
316 * when n = 1
317 */
318{
319 double help0,help1;
320
321 if (n > 1) {
322 help0 = p0;
323 help1 = p1;
324 p0 = ((2*n-1)*x*help0 - (n-1)*(p01))/n;
325 p1 = ((2*n-1)*x*help1 - n*(p11))/(n-1);
326 p01 = help0;
327 p11 = help1;
328 }
329 else if (n == 0) {
330 p0 = 1.0;
331 p1 = 0.0;
332 }
333 else if (n == 1) {
334 p01 = 1.0;
335 p0 = x;
336 p11 = 0.0;
337 p1 = sqrt(1.0-x*x);
338 }
339 return;
340}
341
342//=============================================================================================================
343
344void FwdEegSphereModel::calc_pot_components(double beta, double cgamma, double &Vrp, double &Vtp, const Eigen::VectorXd& fn, int nterms)
345{
346 double Vt = 0.0;
347 double Vr = 0.0;
348 double p0,p01,p1,p11;
349 double betan,multn;
350 int n;
351
352 betan = 1.0;
353 p0 = p01 = p1 = p11 = 0.0;
354 for (n = 1; n <= nterms; n++) {
355 if (betan < EPS) {
356 break;
357 }
358 next_legen (n,cgamma,p0,p01,p1,p11);
359 multn = betan*fn[n-1]; /* The 2*n + 1 factor is included in fn */
360 Vr = Vr + multn*p0;
361 Vt = Vt + multn*p1/n;
362 betan = beta*betan;
363 }
364 Vrp = Vr;
365 Vtp = Vt;
366 return;
367}
368
369//=============================================================================================================
370// fwd_multi_spherepot.c
371int FwdEegSphereModel::fwd_eeg_multi_spherepot(const Eigen::Vector3f& rd_in, const Eigen::Vector3f& Q_in, const Eigen::Matrix<float, Eigen::Dynamic, 3, Eigen::RowMajor>& el, int neeg, Eigen::VectorXf& Vval, void *client)
372/*
373 * Compute the electric potentials in a set of electrodes in spherically
374 * Symmetric head model.
375 *
376 * The code is based on the formulas presented in
377 *
378 * Z. Zhang, A fast method to compute surface potentials
379 * generated by dipoles within multilayer anisotropic spheres,
380 * Phys. Med. Biol., 40, 335 - 349, 1995.
381 *
382 * and
383 *
384 * J.C. Moscher, R.M. Leahy, and P.S. Lewis, Matrix Kernels for
385 * Modeling of EEG and MEG Data, Los Alamos Technical Report,
386 * LA-UR-96-1993, 1996.
387 *
388 * This version does not use the acceleration with help of equivalent sources
389 * in the homogeneous model
390 *
391 */
392{
393 auto* m = static_cast<FwdEegSphereModel*>(client);
394 Eigen::Vector3f rd = rd_in - m->r0;
395 Eigen::Vector3f Q = Q_in;
396 Eigen::Vector3f pos;
397 int k;
398 float pos2,rd_len,pos_len;
399 double beta,cos_gamma,Vr,Vt;
400 Eigen::Vector3f vec1, vec2;
401 float v1,v2;
402 float cos_beta,Qr,Qt,Q2,c;
403 float pi4_inv = 0.25/M_PI;
404 float sigmaM_inv;
405 /*
406 * Precompute the coefficients
407 */
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);
413 }
414 /*
415 * Move to the sphere coordinates
416 */
417 rd_len = rd.norm();
418 Q2 = Q.dot(Q);
419 /*
420 * Ignore dipoles outside the innermost sphere
421 */
422 if (rd_len >= m->layers[0].rad) {
423 for (k = 0; k < neeg; k++)
424 Vval[k] = 0.0;
425 return OK;
426 }
427 /*
428 * Special case: rd and Q are parallel
429 */
430 c = rd.dot(Q)/(rd_len*sqrt(Q2));
431 if ((1.0-c*c) < SIN_EPS) { /* Almost parallel:
432 * Q is purely radial */
433 Qr = sqrt(Q2);
434 Qt = 0.0;
435 cos_beta = 0.0;
436 v1 = 0.0;
437 vec1.setZero();
438 }
439 else {
440 vec1 = rd.cross(Q);
441 v1 = vec1.norm();
442 cos_beta = 0.0;
443 Qr = Qt = 0.0;
444 }
445 for (k = 0; k < neeg; k++) {
446 pos = el.row(k).transpose() - m->r0;
447 /*
448 * Should the position be scaled or not?
449 */
450 if (m->scale_pos) {
451 pos_len = m->layers[m->nlayer()-1].rad/pos.norm();
452#ifdef DEBUG
453 qInfo("%10.4f %10.4f %10.4f %10.4f",pos_len,1000*pos[0],1000*pos[1],1000*pos[2]);
454#endif
455 pos *= pos_len;
456 }
457 pos2 = pos.dot(pos);
458 pos_len = sqrt(pos2);
459 /*
460 * Calculate the two ingredients for the final result
461 */
462 cos_gamma = pos.dot(rd)/(rd_len*pos_len);
463 beta = rd_len/pos_len;
464 calc_pot_components(beta,cos_gamma,Vr,Vt,m->fn,m->nterms);
465 /*
466 * Then compute the combined result
467 */
468 if (v1 > 0.0) {
469 vec2 = rd.cross(pos);
470 v2 = vec2.norm();
471
472 if (v2 > 0.0)
473 cos_beta = vec1.dot(vec2)/(v1*v2);
474 else
475 cos_beta = 0.0;
476
477 Qr = Q.dot(rd)/rd_len;
478 Qt = sqrt(Q2 - Qr*Qr);
479 }
480 Vval[k] = pi4_inv*(Qr*Vr + Qt*cos_beta*Vt)/pos2;
481 }
482 /*
483 * Scale by the conductivity if we have the layers
484 * defined
485 */
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;
490 }
491 return OK;
492}
493
494//=============================================================================================================
495// fwd_multi_spherepot.c
496int FwdEegSphereModel::fwd_eeg_multi_spherepot_coil1(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, FwdCoilSet &els, Eigen::Ref<Eigen::VectorXf> Vval, void *client) /* Client data will be the sphere model definition */
497/*
498 * Calculate the EEG in the sphere model using the fwdCoilSet structure
499 *
500 * This version does not use the acceleration with help of equivalent sources
501 * in the homogeneous model
502 *
503 */
504{
505 VectorXf vval_one;
506 float val;
507 int nvval = 0;
508 int k,c;
509 FwdCoil* el;
510
511 Vval.resize(els.ncoil());
512 for (k = 0; k < els.ncoil(); k++, el++) {
513 el = els.coils[k].get();
514 if (el->coil_class == FWD_COILC_EEG) {
515 if (el->np > nvval) {
516 vval_one.resize(el->np);
517 nvval = el->np;
518 }
519 if (fwd_eeg_multi_spherepot(rd,Q,el->rmag,el->np,vval_one,client) != OK) {
520 return FAIL;
521 }
522 for (c = 0, val = 0.0; c < el->np; c++)
523 val += el->w[c]*vval_one[c];
524 Vval[k] = val;
525 }
526 }
527 return OK;
528}
529
530//=============================================================================================================
531// fwd_multi_spherepot.c
532bool FwdEegSphereModel::fwd_eeg_spherepot_vec(const Eigen::Vector3f& rd_in, const Eigen::Matrix<float, Eigen::Dynamic, 3, Eigen::RowMajor>& el, int neeg, Eigen::MatrixXf& Vval_vec, void *client)
533{
534 auto* m = static_cast<FwdEegSphereModel*>(client);
535 float fact = 0.25f / static_cast<float>(M_PI);
536 Eigen::Vector3f a_vec;
537 float a,a2,a3;
538 float rrd,rd2,rd2_inv,r,r2,ra,rda;
539 float F;
540 float c1,c2,m1,m2;
541 int k,eq;
542 Eigen::Vector3f orig_rd = rd_in - m->r0;
543 Eigen::Vector3f rd;
544 Eigen::Vector3f pos;
545 float pos_len;
546 /*
547 * Initialize the arrays
548 */
549 for (k = 0 ; k < neeg ; k++) {
550 Vval_vec(0,k) = 0.0;
551 Vval_vec(1,k) = 0.0;
552 Vval_vec(2,k) = 0.0;
553 }
554 /*
555 * Ignore dipoles outside the innermost sphere
556 */
557 if (orig_rd.norm() >= m->layers[0].rad)
558 return true;
559 /*
560 * Default to homogeneous model if no model was previously set
561 */
562#ifdef FOO
563 if (nequiv == 0) /* what to do */
564 eeg_set_homog_sphere_model();
565#endif
566 /*
567 * Make a weighted sum over the equivalence parameters
568 */
569 for (eq = 0; eq < m->nfit; eq++) {
570 /*
571 * Scale the dipole position
572 */
573 rd = m->mu[eq] * orig_rd;
574
575 rd2 = rd.dot(rd);
576 rd2_inv = 1.0/rd2;
577
578 /*
579 * Go over all electrodes
580 */
581 for (k = 0; k < neeg ; k++) {
582
583 pos = el.row(k).transpose() - m->r0;
584 /*
585 * Scale location onto the surface of the sphere
586 */
587 if (m->scale_pos) {
588 pos_len = m->layers[m->nlayer()-1].rad/pos.norm();
589 pos *= pos_len;
590 }
591
592 /* Vector from dipole to the field point */
593
594 a_vec = pos - rd;
595
596 /* Compute the dot products needed */
597
598 a2 = a_vec.dot(a_vec); a = sqrt(a2);
599 a3 = 2.0/(a2*a);
600 r2 = pos.dot(pos); r = sqrt(r2);
601 rrd = pos.dot(rd);
602 ra = r2 - rrd;
603 rda = rrd - rd2;
604
605 /* The main ingredients */
606
607 F = a*(r*a + ra);
608 c1 = a3*rda + 1.0/a - 1.0/r;
609 c2 = a3 + (a+r)/(r*F);
610
611 /* Mix them together and scale by lambda/(rd*rd) */
612
613 m1 = (c1 - c2*rrd);
614 m2 = c2*rd2;
615
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]);
619 } /* All electrodes done */
620 } /* All equivalent dipoles done */
621 /*
622 * Finish by scaling by 1/(4*M_PI);
623 */
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);
628 }
629 return true;
630}
631
632//=============================================================================================================
633// fwd_multi_spherepot.c
634int FwdEegSphereModel::fwd_eeg_spherepot_coil_vec(const Eigen::Vector3f& rd, FwdCoilSet& els, Eigen::Ref<Eigen::MatrixXf> Vval_vec, void *client)
635{
636 Eigen::MatrixXf vval_one;
637 float val;
638 int nvval = 0;
639 int k,c,p;
640 FwdCoil* el;
641
642 for (k = 0; k < els.ncoil(); k++, el++) {
643 el = els.coils[k].get();
644 if (el->coil_class == FWD_COILC_EEG) {
645 if (el->np > nvval) {
646 vval_one.resize(3, el->np);
647 nvval = el->np;
648 }
649 if (!fwd_eeg_spherepot_vec(rd,el->rmag,el->np,vval_one,client)) {
650 return FAIL;
651 }
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);
655 Vval_vec(p,k) = val;
656 }
657 }
658 }
659 return OK;
660}
661
662//=============================================================================================================
663
664int FwdEegSphereModel::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) /* Client data to be passed to some foward modelling routines */
665/*
666 * Quick and dirty solution: use differences
667 *
668 * This routine uses the acceleration with help of equivalent sources
669 * in the homogeneous sphere.
670 *
671 */
672{
673 Eigen::Vector3f my_rd;
674 float step = 0.0005;
675 float step2 = 2*step;
676 int p,q;
677
678 Eigen::Ref<Eigen::VectorXf>* grads[3] = { &xgrad, &ygrad, &zgrad };
679
680 for (p = 0; p < 3; p++) {
681 my_rd = rd;
682 my_rd[p] += step;
683 if (fwd_eeg_spherepot_coil(my_rd,Q,coils,*grads[p],client) == FAIL)
684 return FAIL;
685 my_rd = rd;
686 my_rd[p] -= step;
687 if (fwd_eeg_spherepot_coil(my_rd,Q,coils,Vval,client) == FAIL)
688 return FAIL;
689 for (q = 0; q < coils.ncoil(); q++)
690 (*grads[p])[q] = ((*grads[p])[q]-Vval[q])/step2;
691 }
692 if (fwd_eeg_spherepot_coil(rd,Q,coils,Vval,client) == FAIL)
693 return FAIL;
694 return OK;
695}
696
697//=============================================================================================================
698// fwd_multi_spherepot.c
699int FwdEegSphereModel::fwd_eeg_spherepot( const Eigen::Vector3f& rd_in,
700 const Eigen::Vector3f& Q_in,
701 const Eigen::Matrix<float, Eigen::Dynamic, 3, Eigen::RowMajor>& el,
702 int neeg,
703 VectorXf& Vval,
704 void *client)
705/*
706 * This routine calculates the potentials for a specific dipole direction
707 *
708 * This routine uses the acceleration with help of equivalent sources
709 * in the homogeneous sphere.
710 */
711{
712 auto* m = static_cast<FwdEegSphereModel*>(client);
713 float fact = 0.25f/M_PI;
714 Eigen::Vector3f a_vec;
715 float a,a2,a3;
716 float rrd,rd2,rd2_inv,r,r2,ra,rda;
717 float F;
718 float c1,c2,m1,m2,f1,f2;
719 int k,eq;
720 Eigen::Vector3f orig_rd = rd_in - m->r0;
721 Eigen::Vector3f rd;
722 Eigen::Vector3f Q = Q_in;
723 Eigen::Vector3f pos;
724 float pos_len;
725 /*
726 * Initialize the arrays
727 */
728 for (k = 0 ; k < neeg ; k++)
729 Vval[k] = 0.0;
730 /*
731 * Ignore dipoles outside the innermost sphere
732 */
733 if (orig_rd.norm() >= m->layers[0].rad)
734 return true;
735 /*
736 * Default to homogeneous model if no model was previously set
737 */
738#ifdef FOO
739 if (nequiv == 0) /* what to do */
740 eeg_set_homog_sphere_model();
741#endif
742 /*
743 * Make a weighted sum over the equivalence parameters
744 */
745 for (eq = 0; eq < m->nfit; eq++) {
746 /*
747 * Scale the dipole position
748 */
749 rd = m->mu[eq] * orig_rd;
750
751 rd2 = rd.dot(rd);
752 rd2_inv = 1.0/rd2;
753
754 f1 = rd.dot(Q);
755 /*
756 * Go over all electrodes
757 */
758 for (k = 0; k < neeg ; k++) {
759
760 pos = el.row(k).transpose() - m->r0;
761 /*
762 * Scale location onto the surface of the sphere
763 */
764 if (m->scale_pos) {
765 pos_len = m->layers[m->nlayer()-1].rad/pos.norm();
766 pos *= pos_len;
767 }
768
769 /* Vector from dipole to the field point */
770
771 a_vec = pos - rd;
772
773 /* Compute the dot products needed */
774
775 a2 = a_vec.dot(a_vec); a = sqrt(a2);
776 a3 = 2.0/(a2*a);
777 r2 = pos.dot(pos); r = sqrt(r2);
778 rrd = pos.dot(rd);
779 ra = r2 - rrd;
780 rda = rrd - rd2;
781
782 /* The main ingredients */
783
784 F = a*(r*a + ra);
785 c1 = a3*rda + 1.0/a - 1.0/r;
786 c2 = a3 + (a+r)/(r*F);
787
788 /* Mix them together and scale by lambda/(rd*rd) */
789
790 m1 = (c1 - c2*rrd);
791 m2 = c2*rd2;
792
793 f2 = pos.dot(Q);
794 Vval[k] = Vval[k] + m->lambda[eq]*rd2_inv*(m1*f1 + m2*f2);
795 } /* All electrodes done */
796 } /* All equivalent dipoles done */
797 /*
798 * Finish by scaling by 1/(4*M_PI);
799 */
800 for (k = 0; k < neeg; k++)
801 Vval[k] = fact*Vval[k];
802 return OK;
803}
804
805//=============================================================================================================
806// fwd_multi_spherepot.c
807int FwdEegSphereModel::fwd_eeg_spherepot_coil(const Eigen::Vector3f& rd, const Eigen::Vector3f& Q, FwdCoilSet& els, Eigen::Ref<Eigen::VectorXf> Vval, void *client)
808{
809 VectorXf vval_one;
810 float val;
811 int nvval = 0;
812 int k,c;
813 FwdCoil* el;
814
815 for (k = 0; k < els.ncoil(); k++, el++) {
816 el = els.coils[k].get();
817 if (el->coil_class == FWD_COILC_EEG) {
818 if (el->np > nvval) {
819 vval_one.resize(el->np);
820 nvval = el->np;
821 }
822 if (fwd_eeg_spherepot(rd,Q,el->rmag,el->np,vval_one,client) != OK) {
823 return FAIL;
824 }
825 for (c = 0, val = 0.0; c < el->np; c++)
826 val += el->w[c]*vval_one[c];
827 Vval[k] = val;
828 }
829 }
830 return OK;
831}
832
833//=============================================================================================================
834// fwd_eeg_sphere_models.c
835bool FwdEegSphereModel::fwd_setup_eeg_sphere_model(float rad, bool fit_berg_scherg, int nfit)
836{
837 int nterms = 200;
838 float rv;
839
840 /*
841 * Scale the relative radiuses
842 */
843 for (int k = 0; k < this->nlayer(); k++)
844 this->layers[k].rad = rad*this->layers[k].rel_rad;
845
846 if (fit_berg_scherg) {
847 if (this->fwd_eeg_fit_berg_scherg(nterms,nfit,rv)) {
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]);
851 }
852 else
853 return false;
854 }
855
856 qInfo("Defined EEG sphere model with rad = %7.2f mm", 1000.0*rad);
857 return true;
858}
859
860static void compute_svd(Eigen::MatrixXd& mat,
861 Eigen::VectorXd& sing,
862 Eigen::MatrixXd& uu,
863 Eigen::MatrixXd* vv)
864/*
865 * Compute the SVD of mat.
866 * Results are stored in sing, uu, and optionally vv.
867 * mat is not modified.
868 */
869{
870 int udim = std::min(static_cast<int>(mat.rows()), static_cast<int>(mat.cols()));
871
872 Eigen::JacobiSVD<Eigen::MatrixXd> svd(mat, Eigen::ComputeFullU | Eigen::ComputeFullV);
873
874 sing = svd.singularValues();
875 uu = svd.matrixU().transpose().topRows(udim);
876
877 if (vv != nullptr)
878 *vv = svd.matrixV().transpose();
879}
880
881/*
882 * Include the simplex and SVD code here.
883 * It is not too much of a problem
884 */
885
886
887namespace FWDLIB
888{
889
894 double lambda;
895 double mu;
896};
897
898} // namespace
899
900static void sort_parameters(VectorXd& mu,VectorXd& lambda,int nfit)
901{
902 std::vector<BergSchergPar> pars(nfit);
903
904 for (int k = 0; k < nfit; k++) {
905 pars[k].mu = mu[k];
906 pars[k].lambda = lambda[k];
907 }
908
909 std::sort(pars.begin(), pars.end(), [](const BergSchergPar& a, const BergSchergPar& b) {
910 return a.mu > b.mu;
911 });
912
913 for (int k = 0; k < nfit; k++) {
914 mu[k] = pars[k].mu;
915 lambda[k] = pars[k].lambda;
916 }
917}
918
919static bool report_fit(int loop,
920 const VectorXd &fitpar,
921 double Smin,
922 double /*fval_hi*/,
923 double /*par_diff*/)
924{
925#ifdef LOG_FIT
926 for (int k = 0; k < fitpar.size(); k++)
927 qInfo("%g ",mu[k]);
928 qInfo("%g",Smin);
929#endif
930 return true;
931}
932
933static MatrixXd get_initial_simplex(const VectorXd &pars,
934 double simplex_size)
935
936{
937 int npar = pars.size();
938
939 MatrixXd simplex = MatrixXd::Zero(npar+1,npar);
940
941 simplex.rowwise() += pars.transpose();
942
943 for (int k = 1; k < npar+1; k++)
944 simplex(k,k-1) += simplex_size;
945
946 return simplex;
947}
948
950{
951 double mu1n,k1;
952 int k,p;
953 /*
954 * y is the data to be fitted (nterms-1 x 1)
955 * M is the model matrix (nterms-1 x nfit-1)
956 */
957 for (k = 0; k < u->nterms-1; k++) {
958 k1 = k + 1;
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);
963 }
964}
965
966// fwd_fit_berg_scherg.c
968 VectorXd& lambda,
969 fitUser u)
970/*
971 * Compute the best-fitting linear parameters
972 * Return the corresponding RV
973 */
974{
975 int k,p,q;
976 VectorXd vec(u->nfit-1);
977 double sum;
978
980
981 compute_svd(u->M, u->sing, u->uu, &u->vv);
982 /*
983 * Compute the residuals
984 */
985 for (k = 0; k < u->nterms-1; k++)
986 u->resi[k] = u->y[k];
987
988 for (p = 0; p < u->nfit-1; p++) {
989 vec[p] = u->uu.row(p).head(u->nterms-1).dot(u->y.head(u->nterms-1));
990 for (k = 0; k < u->nterms-1; k++)
991 u->resi[k] = u->resi[k] - u->uu(p,k)*vec[p];
992 vec[p] = vec[p]/u->sing[p];
993 }
994
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];
998 lambda[p+1] = sum;
999 }
1000 for (p = 1, sum = 0.0; p < u->nfit; p++)
1001 sum += lambda[p];
1002 lambda[0] = u->fn[0] - sum;
1003 return u->resi.head(u->nterms-1).squaredNorm() / u->y.head(u->nterms-1).squaredNorm();
1004}
1005
1006// fwd_fit_berg_scherg.c
1007double FwdEegSphereModel::one_step (const VectorXd& mu, const void *user_data)
1008/*
1009 * Evaluate the residual sum of squares fit for one set of
1010 * mu values
1011 */
1012{
1013 int k,p;
1014 double dot;
1015 fitUser u = (fitUser)user_data;
1016
1017 for (k = 0; k < u->nfit; k++) {
1018 if (std::fabs(mu[k]) > 1.0)
1019 return 1.0;
1020 }
1021 /*
1022 * Compose the data for the linear fitting
1023 */
1025 /*
1026 * Compute SVD
1027 */
1028 compute_svd(u->M, u->sing, u->uu, nullptr);
1029 /*
1030 * Compute the residuals
1031 */
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++) {
1035 dot = u->uu.row(p).head(u->nterms-1).dot(u->y.head(u->nterms-1));
1036 for (k = 0; k < u->nterms-1; k++)
1037 u->resi[k] = u->resi[k] - u->uu(p,k)*dot;
1038 }
1039 /*
1040 * Return their sum of squares
1041 */
1042 return u->resi.head(u->nterms-1).squaredNorm();
1043}
1044
1045// fwd_fit_berg_scherg.c
1046bool FwdEegSphereModel::fwd_eeg_fit_berg_scherg(int nterms, /* Number of terms to use in the series expansion
1047 * when fitting the parameters */
1048 int nfit, /* Number of equivalent dipoles to fit */
1049 float &rv)
1050/*
1051 * This routine fits the Berg-Scherg equivalent spherical model
1052 * dipole parameters by minimizing the difference between the
1053 * actual and approximative series expansions
1054 */
1055{
1056 bool res = false;
1057 int k;
1058 double rd,R,f;
1059 double simplex_size = 0.01;
1060 MatrixXd simplex;
1061 VectorXd func_val;
1062 double ftol = 1e-9;
1063 VectorXd lambda;
1064 VectorXd mu;
1065 int neval;
1066 int max_eval = 1000;
1067 int report = 1;
1069
1070 if (nfit < 2) {
1071 qWarning("fwd_fit_berg_scherg does not work with less than two equivalent sources.");
1072 return false;
1073 }
1074
1075 /*
1076 * (1) Calculate the coefficients of the true expansion
1077 */
1078 for (k = 0; k < nterms; k++)
1079 u->fn[k] = this->fwd_eeg_get_multi_sphere_model_coeff(k+1);
1080
1081 /*
1082 * (2) Calculate the weighting
1083 */
1084 rd = R = this->layers[0].rad;
1085 for (k = 1; k < this->nlayer(); k++) {
1086 if (this->layers[k].rad > R)
1087 R = this->layers[k].rad;
1088 if (this->layers[k].rad < rd)
1089 rd = this->layers[k].rad;
1090 }
1091 f = rd/R;
1092
1093#ifdef ZHANG
1094 /*
1095 * This is the Zhang weighting
1096 */
1097 for (k = 1; k < nterms; k++)
1098 u->w[k-1] = pow(f,k);
1099#else
1100 /*
1101 * This is the correct weighting
1102 */
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));
1105#endif
1106
1107 /*
1108 * (3) Prepare for simplex minimization
1109 */
1110 func_val = VectorXd(nfit+1);
1111 lambda = VectorXd(nfit);
1112 mu = VectorXd(nfit);
1113 /*
1114 * (4) Rather arbitrary initial guess
1115 */
1116 for (k = 0; k < nfit; k++) {
1117 /*
1118 mu[k] = (k+1)*0.1*f;
1119 */
1120 mu[k] = (rand() / (RAND_MAX + 1.0))*f;//replacement for: mu[k] = drand48()*f;
1121 }
1122
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);
1126
1127 /*
1128 * (5) Do the nonlinear minimization
1129 */
1131 func_val,
1132 ftol,
1133 0.0,
1134 one_step,
1135 u,
1136 max_eval,
1137 neval,
1138 report,
1139 report_fit);
1140
1141 if (res) {
1142 for (k = 0; k < nfit; k++)
1143 mu[k] = simplex(0,k);
1144
1145 /*
1146 * (6) Do the final step: calculation of the linear parameters
1147 */
1149
1150 sort_parameters(mu,lambda,nfit);
1151#ifdef LOG_FIT
1152 qInfo("RV = %g %%",100*rv);
1153#endif
1154 this->mu.resize(nfit);
1155 this->lambda.resize(nfit);
1156 this->nfit = nfit;
1157 for (k = 0; k < nfit; k++) {
1158 this->mu[k] = mu[k];
1159 /*
1160 * This division takes into account the actual conductivities
1161 */
1162 this->lambda[k] = lambda[k]/this->layers[this->nlayer()-1].sigma;
1163#ifdef LOG_FIT
1164 qInfo("lambda%d = %g\tmu%d = %g",k+1,lambda[k],k+1,mu[k]);
1165#endif
1166 }
1167 }
1168
1169 delete u;
1170 return res;
1171}
#define OK
#define FAIL
#define M_PI
#define EPS
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).
Definition compute_fwd.h:91
constexpr int FWD_COILC_EEG
Definition fwd_coil.h:74
fitUserRec * fitUser
Single MEG or EEG sensor coil with integration points, weights, and coordinate frame.
Definition fwd_coil.h:93
Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > rmag
Definition fwd_coil.h:175
Eigen::VectorXf w
Definition fwd_coil.h:177
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
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))