57#include <Eigen/Eigenvalues>
75: m_inverseOperator(p_inverseOperator)
77, m_iELoretaMaxIter(20)
79, m_bELoretaForceEqual(false)
89: m_inverseOperator(p_inverseOperator)
91, m_iELoretaMaxIter(20)
93, m_bELoretaForceEqual(false)
107 qint32 nave = p_fiffEvoked.
nave;
109 if(!m_inverseOperator.check_ch_names(p_fiffEvoked.
info)) {
110 qWarning(
"Channel name check failed.");
121 qInfo(
"Picked %d channels from the data", t_fiffEvoked.
info.
nchan);
124 float tmin = p_fiffEvoked.
times[0];
136 qWarning(
"InvMinimumNorm::calculateInverse - Inverse not setup -> call doInverseSetup first!");
140 if(K.cols() != data.rows()) {
141 qWarning() <<
"InvMinimumNorm::calculateInverse - Dimension mismatch between K.cols() and data.rows() -" << K.cols() <<
"and" << data.rows();
145 MatrixXd sol = K * data;
149 qInfo(
"combining the current components...");
151 MatrixXd sol1(sol.rows()/3,sol.cols());
152 for(qint32 i = 0; i < sol.cols(); ++i)
155 sol1.block(0,i,sol.rows()/3,1) = tmp.cwiseSqrt();
157 sol.resize(sol1.rows(),sol1.cols());
164 sol = inv.noisenorm*sol;
168 qInfo(
"(sLORETA)...");
169 sol = inv.noisenorm*sol;
173 qInfo(
"(eLORETA)...");
174 sol = inv.noisenorm*sol;
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;
186 return InvSourceEstimate(sol, p_vecVertices, tmin, tstep);
196 inv = m_inverseOperator.prepare_inverse_operator(nave, m_fLambda, m_bdSPM || m_beLoreta, m_bsLORETA);
203 qInfo(
"Computing inverse...");
204 inv.assemble_kernel(label, m_sMethod, pick_normal, K, noise_norm, vertno);
206 std::cout <<
"K " << K.rows() <<
" x " << K.cols() << std::endl;
215 return "Minimum Norm Estimate";
222 return m_inverseOperator.src;
229 if(method.compare(
"MNE") == 0)
231 else if(method.compare(
"dSPM") == 0)
233 else if(method.compare(
"sLORETA") == 0)
235 else if(method.compare(
"eLORETA") == 0)
239 qWarning(
"Method not recognized!");
244 qInfo(
"\tSet minimum norm method to %s.", method.toUtf8().constData());
251 int nActive = (
dSPM ? 1 : 0) + (
sLORETA ? 1 : 0) + (eLoreta ? 1 : 0);
254 qWarning(
"Only one method can be active at a time! - Activating dSPM");
262 m_beLoreta = eLoreta;
265 m_sMethod = QString(
"dSPM");
267 m_sMethod = QString(
"sLORETA");
269 m_sMethod = QString(
"eLORETA");
271 m_sMethod = QString(
"MNE");
285 m_iELoretaMaxIter = maxIter;
287 m_bELoretaForceEqual = forceEqual;
292void InvMinimumNorm::computeELoreta()
303 qInfo(
"Computing eLORETA source weights...");
310 qWarning(
"InvMinimumNorm::computeELoreta - Inverse operator missing eigen structures!");
315 const MatrixXd &eigenLeads = inv.
eigen_leads->data;
316 const VectorXd &sing = inv.
sing;
320 MatrixXd G = eigenFields * sing.asDiagonal() * eigenLeads.transpose();
322 const int nChan =
static_cast<int>(G.rows());
324 const int nOrient =
static_cast<int>(G.cols()) / nSrc;
326 if(nOrient != 1 && nOrient != 3) {
327 qWarning(
"InvMinimumNorm::computeELoreta - Unexpected n_orient: %d", nOrient);
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);
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);
347 for(
int i = 0; i < G.cols(); ++i) {
348 G.col(i) *= sourceStd(i);
353 for(
int i = 0; i < sing.size(); ++i) {
354 if(std::abs(sing(i)) > 1e-10 * sing(0))
358 double lambda2 =
static_cast<double>(m_fLambda);
363 const bool useScalar = (nOrient == 1 || m_bELoretaForceEqual);
366 std::vector<Matrix3d> R_mat;
369 R_vec = VectorXd::Ones(nSrc * nOrient);
371 for(
int i = 0; i < R_vec.size(); ++i) {
372 R_vec(i) *= sourceStd(i) * sourceStd(i);
376 for(
int s = 0; s < nSrc; ++s) {
377 R_mat[s] = Matrix3d::Identity();
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);
387 qInfo(
" Fitting up to %d iterations (n_orient=%d, force_equal=%s)...",
388 m_iELoretaMaxIter, nOrient, m_bELoretaForceEqual ?
"true" :
"false");
391 auto computeGRGt = [&]() -> MatrixXd {
396 for(
int i = 0; i < G.cols(); ++i) {
397 GR.col(i) *= R_vec(i);
399 GRGt = GR * G.transpose();
402 MatrixXd RGt = MatrixXd::Zero(nSrc * 3, nChan);
403 for(
int s = 0; s < nSrc; ++s) {
405 MatrixXd Gs = G.middleCols(s * 3, 3);
407 RGt.middleRows(s * 3, 3) = R_mat[s] * Gs.transpose();
412 double trace = GRGt.trace();
413 double norm = trace /
static_cast<double>(nNonZero);
416 if(useScalar) R_vec /= norm;
417 else for(
auto &Rm : R_mat) Rm /= norm;
422 MatrixXd GRGt = computeGRGt();
424 for(
int kk = 0; kk < m_iELoretaMaxIter; ++kk) {
426 SelfAdjointEigenSolver<MatrixXd> eig(GRGt);
427 VectorXd s = eig.eigenvalues().cwiseAbs();
428 MatrixXd u = eig.eigenvectors();
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); });
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]);
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;
448 MatrixXd N = uKeep * sInv.asDiagonal() * uKeep.transpose();
452 std::vector<Matrix3d> R_old_mat;
453 if(useScalar) R_old_vec = R_vec;
454 else R_old_mat = R_mat;
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;
464 }
else if(m_bELoretaForceEqual) {
466 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
467 MatrixXd Gs = G.middleCols(s_idx * 3, 3);
468 Matrix3d M = Gs.transpose() * N * Gs;
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;
477 for(
int d = 0; d < 3; ++d) {
478 R_vec(s_idx * 3 + d) = meanInvSqrt;
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;
488 SelfAdjointEigenSolver<Matrix3d> eigM(M);
489 Vector3d mEig = eigM.eigenvalues();
490 Matrix3d mVec = eigM.eigenvectors();
492 for(
int d = 0; d < 3; ++d) {
493 mPow(d) = (mEig(d) > 1e-30) ? std::pow(mEig(d), -0.5) : 0.0;
495 R_mat[s_idx] = mVec * mPow.asDiagonal() * mVec.transpose();
501 for(
int i = 0; i < R_vec.size(); ++i) {
502 R_vec(i) *= sourceStd(i) * sourceStd(i);
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);
514 GRGt = computeGRGt();
517 double deltaNum = 0.0, deltaDen = 0.0;
519 deltaNum = (R_vec - R_old_vec).norm();
520 deltaDen = R_old_vec.norm();
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();
527 deltaNum = std::sqrt(deltaNum);
528 deltaDen = std::sqrt(deltaDen);
530 double delta = (deltaDen > 1e-30) ? deltaNum / deltaDen : 0.0;
532 if(delta < m_dELoretaEps) {
533 qInfo(
" eLORETA converged on iteration %d (delta=%.2e < eps=%.2e)", kk + 1, delta, m_dELoretaEps);
536 if(kk == m_iELoretaMaxIter - 1) {
537 qWarning(
" eLORETA weight fitting did not converge after %d iterations (delta=%.2e)", m_iELoretaMaxIter, delta);
542 for(
int i = 0; i < G.cols(); ++i) {
543 G.col(i) /= sourceStd(i);
548 std::vector<Matrix3d> R_sqrt_mat;
550 R_sqrt_vec = R_vec.cwiseSqrt();
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();
558 for(
int d = 0; d < 3; ++d) {
559 rSqrt(d) = (rEig(d) > 1e-30) ? std::sqrt(rEig(d)) : 0.0;
561 R_sqrt_mat[s_idx] = rVec * rSqrt.asDiagonal() * rVec.transpose();
568 for(
int i = 0; i < A.cols(); ++i) {
569 A.col(i) *= R_sqrt_vec(i);
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];
580 JacobiSVD<MatrixXd> svd(A, ComputeThinU | ComputeThinV);
581 const VectorXd newSing = svd.singularValues();
582 const MatrixXd newU = svd.matrixU();
583 const MatrixXd newV = svd.matrixV();
586 MatrixXd weightedLeads = newV;
588 for(
int i = 0; i < weightedLeads.rows(); ++i) {
589 weightedLeads.row(i) *= R_sqrt_vec(i);
593 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
594 MatrixXd Vs = newV.middleRows(s_idx * 3, 3);
595 weightedLeads.middleRows(s_idx * 3, 3) = R_sqrt_mat[s_idx] * Vs;
604 inv.eigen_fields->data = newU;
605 inv.eigen_leads->data = weightedLeads;
606 inv.eigen_leads_weighted =
true;
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;
616 qInfo(
" eLORETA inverse operator updated. [done]");
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).
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)
MNE-style inverse operator.
FIFFLIB::FiffNamedMatrix::SDPtr eigen_leads
FIFFLIB::fiff_int_t nsource
FIFFLIB::FiffNamedMatrix::SDPtr eigen_fields
Source Space descritpion.