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];
201 qWarning(
"InvMinimumNorm::calculateInverse - Inverse not setup -> call doInverseSetup first!");
205 if(K.cols() != data.rows()) {
206 qWarning() <<
"InvMinimumNorm::calculateInverse - Dimension mismatch between K.cols() and data.rows() -" << K.cols() <<
"and" << data.rows();
210 MatrixXd sol = K * data;
214 qInfo(
"combining the current components...");
216 MatrixXd sol1(sol.rows()/3,sol.cols());
217 for(qint32 i = 0; i < sol.cols(); ++i)
220 sol1.block(0,i,sol.rows()/3,1) = tmp.cwiseSqrt();
222 sol.resize(sol1.rows(),sol1.cols());
229 sol = inv.noisenorm*sol;
233 qInfo(
"(sLORETA)...");
234 sol = inv.noisenorm*sol;
238 qInfo(
"(eLORETA)...");
239 sol = inv.noisenorm*sol;
244 VectorXi p_vecVertices(inv.src[0].vertno.size() + inv.src[1].vertno.size());
245 p_vecVertices << inv.src[0].vertno, inv.src[1].vertno;
251 return InvSourceEstimate(sol, p_vecVertices, tmin, tstep);
261 inv = m_inverseOperator.prepare_inverse_operator(nave, m_fLambda, m_bdSPM || m_beLoreta, m_bsLORETA);
268 qInfo(
"Computing inverse...");
269 inv.assemble_kernel(label, m_sMethod, pick_normal, K, noise_norm, vertno);
271 std::cout <<
"K " << K.rows() <<
" x " << K.cols() << std::endl;
280 return "Minimum Norm Estimate";
287 return m_inverseOperator.src;
294 if(method.compare(
"MNE") == 0)
296 else if(method.compare(
"dSPM") == 0)
298 else if(method.compare(
"sLORETA") == 0)
300 else if(method.compare(
"eLORETA") == 0)
304 qWarning(
"Method not recognized!");
309 qInfo(
"\tSet minimum norm method to %s.", method.toUtf8().constData());
316 int nActive = (
dSPM ? 1 : 0) + (
sLORETA ? 1 : 0) + (eLoreta ? 1 : 0);
319 qWarning(
"Only one method can be active at a time! - Activating dSPM");
327 m_beLoreta = eLoreta;
330 m_sMethod = QString(
"dSPM");
332 m_sMethod = QString(
"sLORETA");
334 m_sMethod = QString(
"eLORETA");
336 m_sMethod = QString(
"MNE");
350 m_iELoretaMaxIter = maxIter;
352 m_bELoretaForceEqual = forceEqual;
357void InvMinimumNorm::computeELoreta()
368 qInfo(
"Computing eLORETA source weights...");
375 qWarning(
"InvMinimumNorm::computeELoreta - Inverse operator missing eigen structures!");
380 const MatrixXd &eigenLeads = inv.
eigen_leads->data;
381 const VectorXd &sing = inv.
sing;
385 MatrixXd G = eigenFields * sing.asDiagonal() * eigenLeads.transpose();
387 const int nChan =
static_cast<int>(G.rows());
389 const int nOrient =
static_cast<int>(G.cols()) / nSrc;
391 if(nOrient != 1 && nOrient != 3) {
392 qWarning(
"InvMinimumNorm::computeELoreta - Unexpected n_orient: %d", nOrient);
397 if(inv.source_cov && inv.source_cov->data.size() > 0) {
398 for(
int i = 0; i < G.cols(); ++i) {
399 double sc = inv.source_cov->data(i, 0);
400 if(sc > 0) G.col(i) /= std::sqrt(sc);
405 VectorXd sourceStd = VectorXd::Ones(G.cols());
406 if(inv.orient_prior && inv.orient_prior->data.size() > 0) {
407 for(
int i = 0; i < G.cols() && i < inv.orient_prior->data.rows(); ++i) {
408 double op = inv.orient_prior->data(i, 0);
409 if(op > 0) sourceStd(i) *= std::sqrt(op);
412 for(
int i = 0; i < G.cols(); ++i) {
413 G.col(i) *= sourceStd(i);
418 for(
int i = 0; i < sing.size(); ++i) {
419 if(std::abs(sing(i)) > 1e-10 * sing(0))
423 double lambda2 =
static_cast<double>(m_fLambda);
428 const bool useScalar = (nOrient == 1 || m_bELoretaForceEqual);
431 std::vector<Matrix3d> R_mat;
434 R_vec = VectorXd::Ones(nSrc * nOrient);
436 for(
int i = 0; i < R_vec.size(); ++i) {
437 R_vec(i) *= sourceStd(i) * sourceStd(i);
441 for(
int s = 0; s < nSrc; ++s) {
442 R_mat[s] = Matrix3d::Identity();
444 for(
int a = 0; a < 3; ++a) {
445 for(
int b = 0; b < 3; ++b) {
446 R_mat[s](a, b) *= sourceStd(s * 3 + a) * sourceStd(s * 3 + b);
452 qInfo(
" Fitting up to %d iterations (n_orient=%d, force_equal=%s)...",
453 m_iELoretaMaxIter, nOrient, m_bELoretaForceEqual ?
"true" :
"false");
456 auto computeGRGt = [&]() -> MatrixXd {
461 for(
int i = 0; i < G.cols(); ++i) {
462 GR.col(i) *= R_vec(i);
464 GRGt = GR * G.transpose();
467 MatrixXd RGt = MatrixXd::Zero(nSrc * 3, nChan);
468 for(
int s = 0; s < nSrc; ++s) {
470 MatrixXd Gs = G.middleCols(s * 3, 3);
472 RGt.middleRows(s * 3, 3) = R_mat[s] * Gs.transpose();
477 double trace = GRGt.trace();
478 double norm = trace /
static_cast<double>(nNonZero);
481 if(useScalar) R_vec /= norm;
482 else for(
auto &Rm : R_mat) Rm /= norm;
487 MatrixXd GRGt = computeGRGt();
489 for(
int kk = 0; kk < m_iELoretaMaxIter; ++kk) {
491 SelfAdjointEigenSolver<MatrixXd> eig(GRGt);
492 VectorXd s = eig.eigenvalues().cwiseAbs();
493 MatrixXd u = eig.eigenvectors();
497 std::vector<int> idx(s.size());
498 std::iota(idx.begin(), idx.end(), 0);
499 std::sort(idx.begin(), idx.end(), [&s](
int a,
int b) { return s(a) > s(b); });
501 MatrixXd uKeep(nChan, nNonZero);
502 VectorXd sKeep(nNonZero);
503 for(
int i = 0; i < nNonZero && i < static_cast<int>(idx.size()); ++i) {
504 uKeep.col(i) = u.col(idx[i]);
505 sKeep(i) = s(idx[i]);
509 VectorXd sInv(nNonZero);
510 for(
int i = 0; i < nNonZero; ++i) {
511 sInv(i) = (sKeep(i) > 0) ? 1.0 / (sKeep(i) + lambda2) : 0.0;
513 MatrixXd N = uKeep * sInv.asDiagonal() * uKeep.transpose();
517 std::vector<Matrix3d> R_old_mat;
518 if(useScalar) R_old_vec = R_vec;
519 else R_old_mat = R_mat;
525 for(
int i = 0; i < nSrc; ++i) {
526 double val = (NG.col(i).array() * G.col(i).array()).sum();
527 R_vec(i) = (val > 1e-30) ? 1.0 / std::sqrt(val) : 1.0;
529 }
else if(m_bELoretaForceEqual) {
531 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
532 MatrixXd Gs = G.middleCols(s_idx * 3, 3);
533 Matrix3d M = Gs.transpose() * N * Gs;
535 SelfAdjointEigenSolver<Matrix3d> eigM(M);
536 Vector3d mEig = eigM.eigenvalues();
537 double meanInvSqrt = 0;
538 for(
int d = 0; d < 3; ++d) {
539 meanInvSqrt += (mEig(d) > 1e-30) ? 1.0 / std::sqrt(mEig(d)) : 0.0;
542 for(
int d = 0; d < 3; ++d) {
543 R_vec(s_idx * 3 + d) = meanInvSqrt;
548 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
549 MatrixXd Gs = G.middleCols(s_idx * 3, 3);
550 Matrix3d M = Gs.transpose() * N * Gs;
553 SelfAdjointEigenSolver<Matrix3d> eigM(M);
554 Vector3d mEig = eigM.eigenvalues();
555 Matrix3d mVec = eigM.eigenvectors();
557 for(
int d = 0; d < 3; ++d) {
558 mPow(d) = (mEig(d) > 1e-30) ? std::pow(mEig(d), -0.5) : 0.0;
560 R_mat[s_idx] = mVec * mPow.asDiagonal() * mVec.transpose();
566 for(
int i = 0; i < R_vec.size(); ++i) {
567 R_vec(i) *= sourceStd(i) * sourceStd(i);
570 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
571 for(
int a = 0; a < 3; ++a) {
572 for(
int b = 0; b < 3; ++b) {
573 R_mat[s_idx](a, b) *= sourceStd(s_idx * 3 + a) * sourceStd(s_idx * 3 + b);
579 GRGt = computeGRGt();
582 double deltaNum = 0.0, deltaDen = 0.0;
584 deltaNum = (R_vec - R_old_vec).norm();
585 deltaDen = R_old_vec.norm();
587 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
588 Matrix3d diff = R_mat[s_idx] - R_old_mat[s_idx];
589 deltaNum += diff.squaredNorm();
590 deltaDen += R_old_mat[s_idx].squaredNorm();
592 deltaNum = std::sqrt(deltaNum);
593 deltaDen = std::sqrt(deltaDen);
595 double delta = (deltaDen > 1e-30) ? deltaNum / deltaDen : 0.0;
597 if(delta < m_dELoretaEps) {
598 qInfo(
" eLORETA converged on iteration %d (delta=%.2e < eps=%.2e)", kk + 1, delta, m_dELoretaEps);
601 if(kk == m_iELoretaMaxIter - 1) {
602 qWarning(
" eLORETA weight fitting did not converge after %d iterations (delta=%.2e)", m_iELoretaMaxIter, delta);
607 for(
int i = 0; i < G.cols(); ++i) {
608 G.col(i) /= sourceStd(i);
613 std::vector<Matrix3d> R_sqrt_mat;
615 R_sqrt_vec = R_vec.cwiseSqrt();
617 R_sqrt_mat.resize(nSrc);
618 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
619 SelfAdjointEigenSolver<Matrix3d> eigR(R_mat[s_idx]);
620 Vector3d rEig = eigR.eigenvalues();
621 Matrix3d rVec = eigR.eigenvectors();
623 for(
int d = 0; d < 3; ++d) {
624 rSqrt(d) = (rEig(d) > 1e-30) ? std::sqrt(rEig(d)) : 0.0;
626 R_sqrt_mat[s_idx] = rVec * rSqrt.asDiagonal() * rVec.transpose();
633 for(
int i = 0; i < A.cols(); ++i) {
634 A.col(i) *= R_sqrt_vec(i);
637 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
638 MatrixXd Gs = G.middleCols(s_idx * 3, 3);
639 A.middleCols(s_idx * 3, 3) = Gs * R_sqrt_mat[s_idx];
645 JacobiSVD<MatrixXd> svd(A, ComputeThinU | ComputeThinV);
646 const VectorXd newSing = svd.singularValues();
647 const MatrixXd newU = svd.matrixU();
648 const MatrixXd newV = svd.matrixV();
651 MatrixXd weightedLeads = newV;
653 for(
int i = 0; i < weightedLeads.rows(); ++i) {
654 weightedLeads.row(i) *= R_sqrt_vec(i);
658 for(
int s_idx = 0; s_idx < nSrc; ++s_idx) {
659 MatrixXd Vs = newV.middleRows(s_idx * 3, 3);
660 weightedLeads.middleRows(s_idx * 3, 3) = R_sqrt_mat[s_idx] * Vs;
669 inv.eigen_fields->data = newU;
670 inv.eigen_leads->data = weightedLeads;
671 inv.eigen_leads_weighted =
true;
674 VectorXd reginv(newSing.size());
675 for(
int i = 0; i < newSing.size(); ++i) {
676 const double s2 = newSing(i) * newSing(i);
677 reginv(i) = (s2 > 1e-30) ? newSing(i) / (s2 + lambda2) : 0.0;
681 qInfo(
" eLORETA inverse operator updated. [done]");
Minimum norm class declaration.
InvSourceEstimate class declaration.
#define FIFFV_MNE_FREE_ORI
FiffEvoked class declaration.
Linalg class declaration.
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.