v2.0.0
Loading...
Searching...
No Matches
inv_minimum_norm.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "inv_minimum_norm.h"
42
44#include <fiff/fiff_evoked.h>
45#include <math/linalg.h>
46
47#include <iostream>
48#include <cmath>
49#include <algorithm>
50
51//=============================================================================================================
52// EIGEN INCLUDES
53//=============================================================================================================
54
55#include <Eigen/Core>
56#include <Eigen/Dense>
57#include <Eigen/Eigenvalues>
58#include <Eigen/SVD>
59
60//=============================================================================================================
61// USED NAMESPACES
62//=============================================================================================================
63
64using namespace Eigen;
65using namespace MNELIB;
66using namespace INVLIB;
67using namespace UTILSLIB;
68using namespace FIFFLIB;
69
70//=============================================================================================================
71// DEFINE MEMBER METHODS
72//=============================================================================================================
73
74InvMinimumNorm::InvMinimumNorm(const MNEInverseOperator &p_inverseOperator, float lambda, const QString method)
75: m_inverseOperator(p_inverseOperator)
76, m_beLoreta(false)
77, m_iELoretaMaxIter(20)
78, m_dELoretaEps(1e-6)
79, m_bELoretaForceEqual(false)
80, inverseSetup(false)
81{
82 this->setRegularization(lambda);
83 this->setMethod(method);
84}
85
86//=============================================================================================================
87
88InvMinimumNorm::InvMinimumNorm(const MNEInverseOperator &p_inverseOperator, float lambda, bool dSPM, bool sLORETA)
89: m_inverseOperator(p_inverseOperator)
90, m_beLoreta(false)
91, m_iELoretaMaxIter(20)
92, m_dELoretaEps(1e-6)
93, m_bELoretaForceEqual(false)
94, inverseSetup(false)
95{
96 this->setRegularization(lambda);
97 this->setMethod(dSPM, sLORETA);
98}
99
100//=============================================================================================================
101
103{
104 //
105 // Set up the inverse according to the parameters
106 //
107 qint32 nave = p_fiffEvoked.nave;
108
109 if(!m_inverseOperator.check_ch_names(p_fiffEvoked.info)) {
110 qWarning("Channel name check failed.");
111 return InvSourceEstimate();
112 }
113
114 doInverseSetup(nave,pick_normal);
115
116 //
117 // Pick the correct channels from the data
118 //
119 FiffEvoked t_fiffEvoked = p_fiffEvoked.pick_channels(inv.noise_cov->names);
120
121 qInfo("Picked %d channels from the data", t_fiffEvoked.info.nchan);
122
123 //Results
124 float tmin = p_fiffEvoked.times[0];
125 float tstep = 1/t_fiffEvoked.info.sfreq;
126
127 return calculateInverse(t_fiffEvoked.data, tmin, tstep, pick_normal);
128}
129
130//=============================================================================================================
131
132InvSourceEstimate InvMinimumNorm::calculateInverse(const MatrixXd &data, float tmin, float tstep, bool pick_normal) const
133{
134 if(!inverseSetup)
135 {
136 qWarning("InvMinimumNorm::calculateInverse - Inverse not setup -> call doInverseSetup first!");
137 return InvSourceEstimate();
138 }
139
140 if(K.cols() != data.rows()) {
141 qWarning() << "InvMinimumNorm::calculateInverse - Dimension mismatch between K.cols() and data.rows() -" << K.cols() << "and" << data.rows();
142 return InvSourceEstimate();
143 }
144
145 MatrixXd sol = K * data; //apply imaging kernel
146
147 if (inv.source_ori == FIFFV_MNE_FREE_ORI && pick_normal == false)
148 {
149 qInfo("combining the current components...");
150
151 MatrixXd sol1(sol.rows()/3,sol.cols());
152 for(qint32 i = 0; i < sol.cols(); ++i)
153 {
154 VectorXd tmp = Linalg::combine_xyz(sol.col(i));
155 sol1.block(0,i,sol.rows()/3,1) = tmp.cwiseSqrt();
156 }
157 sol.resize(sol1.rows(),sol1.cols());
158 sol = sol1;
159 }
160
161 if (m_bdSPM)
162 {
163 qInfo("(dSPM)...");
164 sol = inv.noisenorm*sol;
165 }
166 else if (m_bsLORETA)
167 {
168 qInfo("(sLORETA)...");
169 sol = inv.noisenorm*sol;
170 }
171 else if (m_beLoreta)
172 {
173 qInfo("(eLORETA)...");
174 sol = inv.noisenorm*sol;
175 }
176 qInfo("[done]");
177
178 //Results
179 VectorXi p_vecVertices(inv.src[0].vertno.size() + inv.src[1].vertno.size());
180 p_vecVertices << inv.src[0].vertno, inv.src[1].vertno;
181
182// VectorXi p_vecVertices();
183// for(qint32 h = 0; h < inv.src.size(); ++h)
184// t_qListVertices.push_back(inv.src[h].vertno);
185
186 return InvSourceEstimate(sol, p_vecVertices, tmin, tstep);
187}
188
189//=============================================================================================================
190
191void InvMinimumNorm::doInverseSetup(qint32 nave, bool pick_normal)
192{
193 //
194 // Set up the inverse according to the parameters
195 //
196 inv = m_inverseOperator.prepare_inverse_operator(nave, m_fLambda, m_bdSPM || m_beLoreta, m_bsLORETA);
197
198 // For eLORETA: recompute source covariance weights before assembling kernel
199 if(m_beLoreta) {
200 computeELoreta();
201 }
202
203 qInfo("Computing inverse...");
204 inv.assemble_kernel(label, m_sMethod, pick_normal, K, noise_norm, vertno);
205
206 std::cout << "K " << K.rows() << " x " << K.cols() << std::endl;
207
208 inverseSetup = true;
209}
210
211//=============================================================================================================
212
213const char* InvMinimumNorm::getName() const
214{
215 return "Minimum Norm Estimate";
216}
217
218//=============================================================================================================
219
221{
222 return m_inverseOperator.src;
223}
224
225//=============================================================================================================
226
227void InvMinimumNorm::setMethod(QString method)
228{
229 if(method.compare("MNE") == 0)
230 setMethod(false, false);
231 else if(method.compare("dSPM") == 0)
232 setMethod(true, false);
233 else if(method.compare("sLORETA") == 0)
234 setMethod(false, true);
235 else if(method.compare("eLORETA") == 0)
236 setMethod(false, false, true);
237 else
238 {
239 qWarning("Method not recognized!");
240 method = "dSPM";
241 setMethod(true, false);
242 }
243
244 qInfo("\tSet minimum norm method to %s.", method.toUtf8().constData());
245}
246
247//=============================================================================================================
248
249void InvMinimumNorm::setMethod(bool dSPM, bool sLORETA, bool eLoreta)
250{
251 int nActive = (dSPM ? 1 : 0) + (sLORETA ? 1 : 0) + (eLoreta ? 1 : 0);
252 if(nActive > 1)
253 {
254 qWarning("Only one method can be active at a time! - Activating dSPM");
255 dSPM = true;
256 sLORETA = false;
257 eLoreta = false;
258 }
259
260 m_bdSPM = dSPM;
261 m_bsLORETA = sLORETA;
262 m_beLoreta = eLoreta;
263
264 if(dSPM)
265 m_sMethod = QString("dSPM");
266 else if(sLORETA)
267 m_sMethod = QString("sLORETA");
268 else if(eLoreta)
269 m_sMethod = QString("eLORETA");
270 else
271 m_sMethod = QString("MNE");
272}
273
274//=============================================================================================================
275
277{
278 m_fLambda = lambda;
279}
280
281//=============================================================================================================
282
283void InvMinimumNorm::setELoretaOptions(int maxIter, double eps, bool forceEqual)
284{
285 m_iELoretaMaxIter = maxIter;
286 m_dELoretaEps = eps;
287 m_bELoretaForceEqual = forceEqual;
288}
289
290//=============================================================================================================
291
292void InvMinimumNorm::computeELoreta()
293{
294 //
295 // eLORETA: Iteratively compute optimized source covariance weights R
296 // so that lambda2 acts consistently across depth.
297 //
298 // Reference: Pascual-Marqui (2007), Discrete, 3D distributed, linear imaging
299 // methods of electric neuronal activity.
300 // Adapted from MNE-Python: mne.minimum_norm._eloreta._compute_eloreta()
301 //
302
303 qInfo("Computing eLORETA source weights...");
304
305 // Reassemble the whitened gain matrix G from the SVD of the inverse operator:
306 // G = eigen_fields * diag(sing) * eigen_leads^T
307 // where eigen_fields is (n_channels, n_comp) and eigen_leads is (n_sources*n_orient, n_comp),
308 // giving G of shape (n_channels, n_sources*n_orient).
309 if(!inv.eigen_fields || !inv.eigen_leads) {
310 qWarning("InvMinimumNorm::computeELoreta - Inverse operator missing eigen structures!");
311 return;
312 }
313
314 const MatrixXd &eigenFields = inv.eigen_fields->data; // (n_channels, n_comp)
315 const MatrixXd &eigenLeads = inv.eigen_leads->data; // (n_sources*n_orient, n_comp)
316 const VectorXd &sing = inv.sing;
317
318 // G = eigenFields * diag(sing) * eigenLeads^T
319 // (n_channels, n_comp) * (n_comp, n_comp) * (n_comp, n_sources*n_orient) = (n_channels, n_sources*n_orient)
320 MatrixXd G = eigenFields * sing.asDiagonal() * eigenLeads.transpose();
321
322 const int nChan = static_cast<int>(G.rows());
323 const int nSrc = inv.nsource;
324 const int nOrient = static_cast<int>(G.cols()) / nSrc;
325
326 if(nOrient != 1 && nOrient != 3) {
327 qWarning("InvMinimumNorm::computeELoreta - Unexpected n_orient: %d", nOrient);
328 return;
329 }
330
331 // Divide by sqrt(source_cov) to undo source covariance weighting
332 if(inv.source_cov && inv.source_cov->data.size() > 0) {
333 for(int i = 0; i < G.cols(); ++i) {
334 double sc = inv.source_cov->data(i, 0);
335 if(sc > 0) G.col(i) /= std::sqrt(sc);
336 }
337 }
338
339 // Restore orientation prior
340 VectorXd sourceStd = VectorXd::Ones(G.cols());
341 if(inv.orient_prior && inv.orient_prior->data.size() > 0) {
342 for(int i = 0; i < G.cols() && i < inv.orient_prior->data.rows(); ++i) {
343 double op = inv.orient_prior->data(i, 0);
344 if(op > 0) sourceStd(i) *= std::sqrt(op);
345 }
346 }
347 for(int i = 0; i < G.cols(); ++i) {
348 G.col(i) *= sourceStd(i);
349 }
350
351 // Compute rank (number of non-zero singular values)
352 int nNonZero = 0;
353 for(int i = 0; i < sing.size(); ++i) {
354 if(std::abs(sing(i)) > 1e-10 * sing(0))
355 ++nNonZero;
356 }
357
358 double lambda2 = static_cast<double>(m_fLambda);
359
360 // Initialize weight matrix R
361 // For fixed orientation (or force_equal): R is a diagonal vector (nSrc * nOrient)
362 // For free orientation: R is block-diagonal (nSrc x 3 x 3)
363 const bool useScalar = (nOrient == 1 || m_bELoretaForceEqual);
364
365 VectorXd R_vec; // For scalar mode: (nSrc * nOrient)
366 std::vector<Matrix3d> R_mat; // For matrix mode: nSrc x (3x3)
367
368 if(useScalar) {
369 R_vec = VectorXd::Ones(nSrc * nOrient);
370 // Apply prior: R *= sourceStd^2
371 for(int i = 0; i < R_vec.size(); ++i) {
372 R_vec(i) *= sourceStd(i) * sourceStd(i);
373 }
374 } else {
375 R_mat.resize(nSrc);
376 for(int s = 0; s < nSrc; ++s) {
377 R_mat[s] = Matrix3d::Identity();
378 // Apply prior as outer product
379 for(int a = 0; a < 3; ++a) {
380 for(int b = 0; b < 3; ++b) {
381 R_mat[s](a, b) *= sourceStd(s * 3 + a) * sourceStd(s * 3 + b);
382 }
383 }
384 }
385 }
386
387 qInfo(" Fitting up to %d iterations (n_orient=%d, force_equal=%s)...",
388 m_iELoretaMaxIter, nOrient, m_bELoretaForceEqual ? "true" : "false");
389
390 // Lambda for computing G * R * G^T and normalizing
391 auto computeGRGt = [&]() -> MatrixXd {
392 MatrixXd GRGt;
393 if(useScalar) {
394 // G_R = G * diag(R)
395 MatrixXd GR = G;
396 for(int i = 0; i < G.cols(); ++i) {
397 GR.col(i) *= R_vec(i);
398 }
399 GRGt = GR * G.transpose();
400 } else {
401 // Block multiplication
402 MatrixXd RGt = MatrixXd::Zero(nSrc * 3, nChan);
403 for(int s = 0; s < nSrc; ++s) {
404 // G_s: (nChan, 3)
405 MatrixXd Gs = G.middleCols(s * 3, 3);
406 // R_s * G_s^T: (3, nChan)
407 RGt.middleRows(s * 3, 3) = R_mat[s] * Gs.transpose();
408 }
409 GRGt = G * RGt;
410 }
411 // Normalize so trace(GRGt) / nNonZero = 1
412 double trace = GRGt.trace();
413 double norm = trace / static_cast<double>(nNonZero);
414 if(norm > 1e-30) {
415 GRGt /= norm;
416 if(useScalar) R_vec /= norm;
417 else for(auto &Rm : R_mat) Rm /= norm;
418 }
419 return GRGt;
420 };
421
422 MatrixXd GRGt = computeGRGt();
423
424 for(int kk = 0; kk < m_iELoretaMaxIter; ++kk) {
425 // 1. Compute inverse of GRGt (stabilized eigendecomposition)
426 SelfAdjointEigenSolver<MatrixXd> eig(GRGt);
427 VectorXd s = eig.eigenvalues().cwiseAbs();
428 MatrixXd u = eig.eigenvectors();
429
430 // Keep top nNonZero eigenvalues
431 // Sort descending
432 std::vector<int> idx(s.size());
433 std::iota(idx.begin(), idx.end(), 0);
434 std::sort(idx.begin(), idx.end(), [&s](int a, int b) { return s(a) > s(b); });
435
436 MatrixXd uKeep(nChan, nNonZero);
437 VectorXd sKeep(nNonZero);
438 for(int i = 0; i < nNonZero && i < static_cast<int>(idx.size()); ++i) {
439 uKeep.col(i) = u.col(idx[i]);
440 sKeep(i) = s(idx[i]);
441 }
442
443 // N = u * diag(1/(s + lambda2)) * u^T
444 VectorXd sInv(nNonZero);
445 for(int i = 0; i < nNonZero; ++i) {
446 sInv(i) = (sKeep(i) > 0) ? 1.0 / (sKeep(i) + lambda2) : 0.0;
447 }
448 MatrixXd N = uKeep * sInv.asDiagonal() * uKeep.transpose();
449
450 // Save old R for convergence check
451 VectorXd R_old_vec;
452 std::vector<Matrix3d> R_old_mat;
453 if(useScalar) R_old_vec = R_vec;
454 else R_old_mat = R_mat;
455
456 // 2. Update R
457 if(nOrient == 1) {
458 // R_i = 1 / sqrt(sum_j(N * G)_ij * G_ij)
459 MatrixXd NG = N * G; // (nChan, nDipoles)
460 for(int i = 0; i < nSrc; ++i) {
461 double val = (NG.col(i).array() * G.col(i).array()).sum();
462 R_vec(i) = (val > 1e-30) ? 1.0 / std::sqrt(val) : 1.0;
463 }
464 } else if(m_bELoretaForceEqual) {
465 // For force_equal: compute M_s = G_s^T N G_s, then average eigenvalues of M_s^{-1/2}
466 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
467 MatrixXd Gs = G.middleCols(s_idx * 3, 3); // (nChan, 3)
468 Matrix3d M = Gs.transpose() * N * Gs;
469
470 SelfAdjointEigenSolver<Matrix3d> eigM(M);
471 Vector3d mEig = eigM.eigenvalues();
472 double meanInvSqrt = 0;
473 for(int d = 0; d < 3; ++d) {
474 meanInvSqrt += (mEig(d) > 1e-30) ? 1.0 / std::sqrt(mEig(d)) : 0.0;
475 }
476 meanInvSqrt /= 3.0;
477 for(int d = 0; d < 3; ++d) {
478 R_vec(s_idx * 3 + d) = meanInvSqrt;
479 }
480 }
481 } else {
482 // Free orientation, independent: R_s = sqrtm(G_s^T N G_s)^{-1/2}
483 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
484 MatrixXd Gs = G.middleCols(s_idx * 3, 3);
485 Matrix3d M = Gs.transpose() * N * Gs;
486
487 // Symmetric matrix power: M^{-1/2}
488 SelfAdjointEigenSolver<Matrix3d> eigM(M);
489 Vector3d mEig = eigM.eigenvalues();
490 Matrix3d mVec = eigM.eigenvectors();
491 Vector3d mPow;
492 for(int d = 0; d < 3; ++d) {
493 mPow(d) = (mEig(d) > 1e-30) ? std::pow(mEig(d), -0.5) : 0.0;
494 }
495 R_mat[s_idx] = mVec * mPow.asDiagonal() * mVec.transpose();
496 }
497 }
498
499 // Reapply prior
500 if(useScalar) {
501 for(int i = 0; i < R_vec.size(); ++i) {
502 R_vec(i) *= sourceStd(i) * sourceStd(i);
503 }
504 } else {
505 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
506 for(int a = 0; a < 3; ++a) {
507 for(int b = 0; b < 3; ++b) {
508 R_mat[s_idx](a, b) *= sourceStd(s_idx * 3 + a) * sourceStd(s_idx * 3 + b);
509 }
510 }
511 }
512 }
513
514 GRGt = computeGRGt();
515
516 // 3. Check convergence
517 double deltaNum = 0.0, deltaDen = 0.0;
518 if(useScalar) {
519 deltaNum = (R_vec - R_old_vec).norm();
520 deltaDen = R_old_vec.norm();
521 } else {
522 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
523 Matrix3d diff = R_mat[s_idx] - R_old_mat[s_idx];
524 deltaNum += diff.squaredNorm();
525 deltaDen += R_old_mat[s_idx].squaredNorm();
526 }
527 deltaNum = std::sqrt(deltaNum);
528 deltaDen = std::sqrt(deltaDen);
529 }
530 double delta = (deltaDen > 1e-30) ? deltaNum / deltaDen : 0.0;
531
532 if(delta < m_dELoretaEps) {
533 qInfo(" eLORETA converged on iteration %d (delta=%.2e < eps=%.2e)", kk + 1, delta, m_dELoretaEps);
534 break;
535 }
536 if(kk == m_iELoretaMaxIter - 1) {
537 qWarning(" eLORETA weight fitting did not converge after %d iterations (delta=%.2e)", m_iELoretaMaxIter, delta);
538 }
539 }
540
541 // Undo source_std weighting on G
542 for(int i = 0; i < G.cols(); ++i) {
543 G.col(i) /= sourceStd(i);
544 }
545
546 // Compute R^{1/2}
547 VectorXd R_sqrt_vec;
548 std::vector<Matrix3d> R_sqrt_mat;
549 if(useScalar) {
550 R_sqrt_vec = R_vec.cwiseSqrt();
551 } else {
552 R_sqrt_mat.resize(nSrc);
553 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
554 SelfAdjointEigenSolver<Matrix3d> eigR(R_mat[s_idx]);
555 Vector3d rEig = eigR.eigenvalues();
556 Matrix3d rVec = eigR.eigenvectors();
557 Vector3d rSqrt;
558 for(int d = 0; d < 3; ++d) {
559 rSqrt(d) = (rEig(d) > 1e-30) ? std::sqrt(rEig(d)) : 0.0;
560 }
561 R_sqrt_mat[s_idx] = rVec * rSqrt.asDiagonal() * rVec.transpose();
562 }
563 }
564
565 // Compute weighted gain: A = G * R^{1/2}
566 MatrixXd A = G;
567 if(useScalar) {
568 for(int i = 0; i < A.cols(); ++i) {
569 A.col(i) *= R_sqrt_vec(i);
570 }
571 } else {
572 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
573 MatrixXd Gs = G.middleCols(s_idx * 3, 3);
574 A.middleCols(s_idx * 3, 3) = Gs * R_sqrt_mat[s_idx];
575 }
576 }
577
578 // SVD of A = G * R^{1/2}, shape (nChan, nSrc*nOrient)
579 // A = U * Σ * V^T → U: (nChan, ncomp), V: (nSrc*nOrient, ncomp)
580 JacobiSVD<MatrixXd> svd(A, ComputeThinU | ComputeThinV);
581 const VectorXd newSing = svd.singularValues(); // (ncomp,)
582 const MatrixXd newU = svd.matrixU(); // (nChan, ncomp) — new eigen_fields
583 const MatrixXd newV = svd.matrixV(); // (nSrc*nOrient, ncomp)
584
585 // Build R^{1/2}-weighted eigen_leads: weightedLeads[i,:] = R_sqrt[i] * V[i,:]
586 MatrixXd weightedLeads = newV; // (nSrc*nOrient, ncomp)
587 if(useScalar) {
588 for(int i = 0; i < weightedLeads.rows(); ++i) {
589 weightedLeads.row(i) *= R_sqrt_vec(i);
590 }
591 } else {
592 // For each source s: rows [s*3, s*3+3) = R_sqrt_mat[s] * V[s*3, s*3+3)
593 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
594 MatrixXd Vs = newV.middleRows(s_idx * 3, 3); // (3, ncomp)
595 weightedLeads.middleRows(s_idx * 3, 3) = R_sqrt_mat[s_idx] * Vs;
596 }
597 }
598
599 // Update inverse operator:
600 // eigen_fields: (nChan, ncomp) = U
601 // eigen_leads: (nSrc*nOrient, ncomp) = R_sqrt * V (eLORETA-weighted)
602 // eigen_leads_weighted = true → prepare_inverse_operator uses K = eigenLeads * trans directly
603 inv.sing = newSing;
604 inv.eigen_fields->data = newU;
605 inv.eigen_leads->data = weightedLeads;
606 inv.eigen_leads_weighted = true;
607
608 // Recompute reginv: σ / (σ² + λ²)
609 VectorXd reginv(newSing.size());
610 for(int i = 0; i < newSing.size(); ++i) {
611 const double s2 = newSing(i) * newSing(i);
612 reginv(i) = (s2 > 1e-30) ? newSing(i) / (s2 + lambda2) : 0.0;
613 }
614 inv.reginv = reginv;
615
616 qInfo(" eLORETA inverse operator updated. [done]");
617}
Linalg class declaration.
InvSourceEstimate class declaration.
Minimum norm class declaration.
FiffEvoked class declaration.
#define FIFFV_MNE_FREE_ORI
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Inverse source estimation (MNE, dSPM, sLORETA, dipole fitting).
Eigen::RowVectorXf times
Eigen::MatrixXd data
FiffEvoked pick_channels(const QStringList &include=defaultQStringList, const QStringList &exclude=defaultQStringList) const
virtual const MNELIB::MNESourceSpaces & getSourceSpace() const
void setELoretaOptions(int maxIter=20, double eps=1e-6, bool forceEqual=false)
virtual InvSourceEstimate calculateInverse(const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal=false)
virtual void doInverseSetup(qint32 nave, bool pick_normal=false)
void setMethod(QString method)
void setRegularization(float lambda)
virtual const char * getName() const
InvMinimumNorm(const MNELIB::MNEInverseOperator &p_inverseOperator, float lambda, const QString method)
static Eigen::VectorXd combine_xyz(const Eigen::VectorXd &vec)
Definition linalg.cpp:64
MNE-style inverse operator.
FIFFLIB::FiffNamedMatrix::SDPtr eigen_leads
FIFFLIB::FiffNamedMatrix::SDPtr eigen_fields
Source Space descritpion.