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// // Set up the inverse according to the parameters
131// //
132// qint32 nave = p_fiffEvoked.nave;
133
134// if(!m_inverseOperator.check_ch_names(p_fiffEvoked.info))
135// {
136// qWarning("Channel name check failed.");
137// return SourceEstimate();
138// }
139
140// //ToDo his could be heavily accelerated for real time calculation -> ToDo calculate inverse RT
141// MNEInverseOperator inv = m_inverseOperator.prepare_inverse_operator(nave, m_fLambda, m_bdSPM, m_bsLORETA);
142// //
143// // Pick the correct channels from the data
144// //
145// FiffEvoked t_fiffEvoked = p_fiffEvoked.pick_channels(inv.noise_cov->names);
146
147// printf("Picked %d channels from the data\n",t_fiffEvoked.info.nchan);
148// printf("Computing inverse...");
149
150// MatrixXd K;
151// SparseMatrix<double> noise_norm;
152// QList<VectorXi> vertno;
153// FsLabel label;
154// inv.assemble_kernel(label, m_sMethod, pick_normal, K, noise_norm, vertno);
155
156// MatrixXd sol = K * t_fiffEvoked.data; //apply imaging kernel
157
158// if (inv.source_ori == FIFFV_MNE_FREE_ORI)
159// {
160// printf("combining the current components...");
161// MatrixXd sol1(sol.rows()/3,sol.cols());
162// for(qint32 i = 0; i < sol.cols(); ++i)
163// {
164// VectorXd* tmp = Linalg::combine_xyz(sol.block(0,i,sol.rows(),1));
165// sol1.block(0,i,sol.rows()/3,1) = tmp->cwiseSqrt();
166// delete tmp;
167// }
168// sol.resize(sol1.rows(),sol1.cols());
169// sol = sol1;
170// }
171
172// if (m_bdSPM)
173// {
174// printf("(dSPM)...");
175// sol = inv.noisenorm*sol;
176// }
177// else if (m_bsLORETA)
178// {
179// printf("(sLORETA)...");
180// sol = inv.noisenorm*sol;
181// }
182// printf("[done]\n");
183
184// //Results
185// float tmin = ((float)t_fiffEvoked.first) / t_fiffEvoked.info.sfreq;
186// float tstep = 1/t_fiffEvoked.info.sfreq;
187
188// QList<VectorXi> t_qListVertices;
189// for(qint32 h = 0; h < inv.src.size(); ++h)
190// t_qListVertices.push_back(inv.src[h].vertno);
191
192// return SourceEstimate(sol, t_qListVertices, tmin, tstep);
193}
194
195//=============================================================================================================
196
197InvSourceEstimate InvMinimumNorm::calculateInverse(const MatrixXd &data, float tmin, float tstep, bool pick_normal) const
198{
199 if(!inverseSetup)
200 {
201 qWarning("InvMinimumNorm::calculateInverse - Inverse not setup -> call doInverseSetup first!");
202 return InvSourceEstimate();
203 }
204
205 if(K.cols() != data.rows()) {
206 qWarning() << "InvMinimumNorm::calculateInverse - Dimension mismatch between K.cols() and data.rows() -" << K.cols() << "and" << data.rows();
207 return InvSourceEstimate();
208 }
209
210 MatrixXd sol = K * data; //apply imaging kernel
211
212 if (inv.source_ori == FIFFV_MNE_FREE_ORI && pick_normal == false)
213 {
214 qInfo("combining the current components...");
215
216 MatrixXd sol1(sol.rows()/3,sol.cols());
217 for(qint32 i = 0; i < sol.cols(); ++i)
218 {
219 VectorXd tmp = Linalg::combine_xyz(sol.col(i));
220 sol1.block(0,i,sol.rows()/3,1) = tmp.cwiseSqrt();
221 }
222 sol.resize(sol1.rows(),sol1.cols());
223 sol = sol1;
224 }
225
226 if (m_bdSPM)
227 {
228 qInfo("(dSPM)...");
229 sol = inv.noisenorm*sol;
230 }
231 else if (m_bsLORETA)
232 {
233 qInfo("(sLORETA)...");
234 sol = inv.noisenorm*sol;
235 }
236 else if (m_beLoreta)
237 {
238 qInfo("(eLORETA)...");
239 sol = inv.noisenorm*sol;
240 }
241 qInfo("[done]");
242
243 //Results
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;
246
247// VectorXi p_vecVertices();
248// for(qint32 h = 0; h < inv.src.size(); ++h)
249// t_qListVertices.push_back(inv.src[h].vertno);
250
251 return InvSourceEstimate(sol, p_vecVertices, tmin, tstep);
252}
253
254//=============================================================================================================
255
256void InvMinimumNorm::doInverseSetup(qint32 nave, bool pick_normal)
257{
258 //
259 // Set up the inverse according to the parameters
260 //
261 inv = m_inverseOperator.prepare_inverse_operator(nave, m_fLambda, m_bdSPM || m_beLoreta, m_bsLORETA);
262
263 // For eLORETA: recompute source covariance weights before assembling kernel
264 if(m_beLoreta) {
265 computeELoreta();
266 }
267
268 qInfo("Computing inverse...");
269 inv.assemble_kernel(label, m_sMethod, pick_normal, K, noise_norm, vertno);
270
271 std::cout << "K " << K.rows() << " x " << K.cols() << std::endl;
272
273 inverseSetup = true;
274}
275
276//=============================================================================================================
277
278const char* InvMinimumNorm::getName() const
279{
280 return "Minimum Norm Estimate";
281}
282
283//=============================================================================================================
284
286{
287 return m_inverseOperator.src;
288}
289
290//=============================================================================================================
291
292void InvMinimumNorm::setMethod(QString method)
293{
294 if(method.compare("MNE") == 0)
295 setMethod(false, false);
296 else if(method.compare("dSPM") == 0)
297 setMethod(true, false);
298 else if(method.compare("sLORETA") == 0)
299 setMethod(false, true);
300 else if(method.compare("eLORETA") == 0)
301 setMethod(false, false, true);
302 else
303 {
304 qWarning("Method not recognized!");
305 method = "dSPM";
306 setMethod(true, false);
307 }
308
309 qInfo("\tSet minimum norm method to %s.", method.toUtf8().constData());
310}
311
312//=============================================================================================================
313
314void InvMinimumNorm::setMethod(bool dSPM, bool sLORETA, bool eLoreta)
315{
316 int nActive = (dSPM ? 1 : 0) + (sLORETA ? 1 : 0) + (eLoreta ? 1 : 0);
317 if(nActive > 1)
318 {
319 qWarning("Only one method can be active at a time! - Activating dSPM");
320 dSPM = true;
321 sLORETA = false;
322 eLoreta = false;
323 }
324
325 m_bdSPM = dSPM;
326 m_bsLORETA = sLORETA;
327 m_beLoreta = eLoreta;
328
329 if(dSPM)
330 m_sMethod = QString("dSPM");
331 else if(sLORETA)
332 m_sMethod = QString("sLORETA");
333 else if(eLoreta)
334 m_sMethod = QString("eLORETA");
335 else
336 m_sMethod = QString("MNE");
337}
338
339//=============================================================================================================
340
342{
343 m_fLambda = lambda;
344}
345
346//=============================================================================================================
347
348void InvMinimumNorm::setELoretaOptions(int maxIter, double eps, bool forceEqual)
349{
350 m_iELoretaMaxIter = maxIter;
351 m_dELoretaEps = eps;
352 m_bELoretaForceEqual = forceEqual;
353}
354
355//=============================================================================================================
356
357void InvMinimumNorm::computeELoreta()
358{
359 //
360 // eLORETA: Iteratively compute optimized source covariance weights R
361 // so that lambda2 acts consistently across depth.
362 //
363 // Reference: Pascual-Marqui (2007), Discrete, 3D distributed, linear imaging
364 // methods of electric neuronal activity.
365 // Adapted from MNE-Python: mne.minimum_norm._eloreta._compute_eloreta()
366 //
367
368 qInfo("Computing eLORETA source weights...");
369
370 // Reassemble the whitened gain matrix G from the SVD of the inverse operator:
371 // G = eigen_fields * diag(sing) * eigen_leads^T
372 // where eigen_fields is (n_channels, n_comp) and eigen_leads is (n_sources*n_orient, n_comp),
373 // giving G of shape (n_channels, n_sources*n_orient).
374 if(!inv.eigen_fields || !inv.eigen_leads) {
375 qWarning("InvMinimumNorm::computeELoreta - Inverse operator missing eigen structures!");
376 return;
377 }
378
379 const MatrixXd &eigenFields = inv.eigen_fields->data; // (n_channels, n_comp)
380 const MatrixXd &eigenLeads = inv.eigen_leads->data; // (n_sources*n_orient, n_comp)
381 const VectorXd &sing = inv.sing;
382
383 // G = eigenFields * diag(sing) * eigenLeads^T
384 // (n_channels, n_comp) * (n_comp, n_comp) * (n_comp, n_sources*n_orient) = (n_channels, n_sources*n_orient)
385 MatrixXd G = eigenFields * sing.asDiagonal() * eigenLeads.transpose();
386
387 const int nChan = static_cast<int>(G.rows());
388 const int nSrc = inv.nsource;
389 const int nOrient = static_cast<int>(G.cols()) / nSrc;
390
391 if(nOrient != 1 && nOrient != 3) {
392 qWarning("InvMinimumNorm::computeELoreta - Unexpected n_orient: %d", nOrient);
393 return;
394 }
395
396 // Divide by sqrt(source_cov) to undo source covariance weighting
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);
401 }
402 }
403
404 // Restore orientation prior
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);
410 }
411 }
412 for(int i = 0; i < G.cols(); ++i) {
413 G.col(i) *= sourceStd(i);
414 }
415
416 // Compute rank (number of non-zero singular values)
417 int nNonZero = 0;
418 for(int i = 0; i < sing.size(); ++i) {
419 if(std::abs(sing(i)) > 1e-10 * sing(0))
420 ++nNonZero;
421 }
422
423 double lambda2 = static_cast<double>(m_fLambda);
424
425 // Initialize weight matrix R
426 // For fixed orientation (or force_equal): R is a diagonal vector (nSrc * nOrient)
427 // For free orientation: R is block-diagonal (nSrc x 3 x 3)
428 const bool useScalar = (nOrient == 1 || m_bELoretaForceEqual);
429
430 VectorXd R_vec; // For scalar mode: (nSrc * nOrient)
431 std::vector<Matrix3d> R_mat; // For matrix mode: nSrc x (3x3)
432
433 if(useScalar) {
434 R_vec = VectorXd::Ones(nSrc * nOrient);
435 // Apply prior: R *= sourceStd^2
436 for(int i = 0; i < R_vec.size(); ++i) {
437 R_vec(i) *= sourceStd(i) * sourceStd(i);
438 }
439 } else {
440 R_mat.resize(nSrc);
441 for(int s = 0; s < nSrc; ++s) {
442 R_mat[s] = Matrix3d::Identity();
443 // Apply prior as outer product
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);
447 }
448 }
449 }
450 }
451
452 qInfo(" Fitting up to %d iterations (n_orient=%d, force_equal=%s)...",
453 m_iELoretaMaxIter, nOrient, m_bELoretaForceEqual ? "true" : "false");
454
455 // Lambda for computing G * R * G^T and normalizing
456 auto computeGRGt = [&]() -> MatrixXd {
457 MatrixXd GRGt;
458 if(useScalar) {
459 // G_R = G * diag(R)
460 MatrixXd GR = G;
461 for(int i = 0; i < G.cols(); ++i) {
462 GR.col(i) *= R_vec(i);
463 }
464 GRGt = GR * G.transpose();
465 } else {
466 // Block multiplication
467 MatrixXd RGt = MatrixXd::Zero(nSrc * 3, nChan);
468 for(int s = 0; s < nSrc; ++s) {
469 // G_s: (nChan, 3)
470 MatrixXd Gs = G.middleCols(s * 3, 3);
471 // R_s * G_s^T: (3, nChan)
472 RGt.middleRows(s * 3, 3) = R_mat[s] * Gs.transpose();
473 }
474 GRGt = G * RGt;
475 }
476 // Normalize so trace(GRGt) / nNonZero = 1
477 double trace = GRGt.trace();
478 double norm = trace / static_cast<double>(nNonZero);
479 if(norm > 1e-30) {
480 GRGt /= norm;
481 if(useScalar) R_vec /= norm;
482 else for(auto &Rm : R_mat) Rm /= norm;
483 }
484 return GRGt;
485 };
486
487 MatrixXd GRGt = computeGRGt();
488
489 for(int kk = 0; kk < m_iELoretaMaxIter; ++kk) {
490 // 1. Compute inverse of GRGt (stabilized eigendecomposition)
491 SelfAdjointEigenSolver<MatrixXd> eig(GRGt);
492 VectorXd s = eig.eigenvalues().cwiseAbs();
493 MatrixXd u = eig.eigenvectors();
494
495 // Keep top nNonZero eigenvalues
496 // Sort descending
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); });
500
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]);
506 }
507
508 // N = u * diag(1/(s + lambda2)) * u^T
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;
512 }
513 MatrixXd N = uKeep * sInv.asDiagonal() * uKeep.transpose();
514
515 // Save old R for convergence check
516 VectorXd R_old_vec;
517 std::vector<Matrix3d> R_old_mat;
518 if(useScalar) R_old_vec = R_vec;
519 else R_old_mat = R_mat;
520
521 // 2. Update R
522 if(nOrient == 1) {
523 // R_i = 1 / sqrt(sum_j(N * G)_ij * G_ij)
524 MatrixXd NG = N * G; // (nChan, nDipoles)
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;
528 }
529 } else if(m_bELoretaForceEqual) {
530 // For force_equal: compute M_s = G_s^T N G_s, then average eigenvalues of M_s^{-1/2}
531 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
532 MatrixXd Gs = G.middleCols(s_idx * 3, 3); // (nChan, 3)
533 Matrix3d M = Gs.transpose() * N * Gs;
534
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;
540 }
541 meanInvSqrt /= 3.0;
542 for(int d = 0; d < 3; ++d) {
543 R_vec(s_idx * 3 + d) = meanInvSqrt;
544 }
545 }
546 } else {
547 // Free orientation, independent: R_s = sqrtm(G_s^T N G_s)^{-1/2}
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;
551
552 // Symmetric matrix power: M^{-1/2}
553 SelfAdjointEigenSolver<Matrix3d> eigM(M);
554 Vector3d mEig = eigM.eigenvalues();
555 Matrix3d mVec = eigM.eigenvectors();
556 Vector3d mPow;
557 for(int d = 0; d < 3; ++d) {
558 mPow(d) = (mEig(d) > 1e-30) ? std::pow(mEig(d), -0.5) : 0.0;
559 }
560 R_mat[s_idx] = mVec * mPow.asDiagonal() * mVec.transpose();
561 }
562 }
563
564 // Reapply prior
565 if(useScalar) {
566 for(int i = 0; i < R_vec.size(); ++i) {
567 R_vec(i) *= sourceStd(i) * sourceStd(i);
568 }
569 } else {
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);
574 }
575 }
576 }
577 }
578
579 GRGt = computeGRGt();
580
581 // 3. Check convergence
582 double deltaNum = 0.0, deltaDen = 0.0;
583 if(useScalar) {
584 deltaNum = (R_vec - R_old_vec).norm();
585 deltaDen = R_old_vec.norm();
586 } else {
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();
591 }
592 deltaNum = std::sqrt(deltaNum);
593 deltaDen = std::sqrt(deltaDen);
594 }
595 double delta = (deltaDen > 1e-30) ? deltaNum / deltaDen : 0.0;
596
597 if(delta < m_dELoretaEps) {
598 qInfo(" eLORETA converged on iteration %d (delta=%.2e < eps=%.2e)", kk + 1, delta, m_dELoretaEps);
599 break;
600 }
601 if(kk == m_iELoretaMaxIter - 1) {
602 qWarning(" eLORETA weight fitting did not converge after %d iterations (delta=%.2e)", m_iELoretaMaxIter, delta);
603 }
604 }
605
606 // Undo source_std weighting on G
607 for(int i = 0; i < G.cols(); ++i) {
608 G.col(i) /= sourceStd(i);
609 }
610
611 // Compute R^{1/2}
612 VectorXd R_sqrt_vec;
613 std::vector<Matrix3d> R_sqrt_mat;
614 if(useScalar) {
615 R_sqrt_vec = R_vec.cwiseSqrt();
616 } else {
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();
622 Vector3d rSqrt;
623 for(int d = 0; d < 3; ++d) {
624 rSqrt(d) = (rEig(d) > 1e-30) ? std::sqrt(rEig(d)) : 0.0;
625 }
626 R_sqrt_mat[s_idx] = rVec * rSqrt.asDiagonal() * rVec.transpose();
627 }
628 }
629
630 // Compute weighted gain: A = G * R^{1/2}
631 MatrixXd A = G;
632 if(useScalar) {
633 for(int i = 0; i < A.cols(); ++i) {
634 A.col(i) *= R_sqrt_vec(i);
635 }
636 } else {
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];
640 }
641 }
642
643 // SVD of A = G * R^{1/2}, shape (nChan, nSrc*nOrient)
644 // A = U * Σ * V^T → U: (nChan, ncomp), V: (nSrc*nOrient, ncomp)
645 JacobiSVD<MatrixXd> svd(A, ComputeThinU | ComputeThinV);
646 const VectorXd newSing = svd.singularValues(); // (ncomp,)
647 const MatrixXd newU = svd.matrixU(); // (nChan, ncomp) — new eigen_fields
648 const MatrixXd newV = svd.matrixV(); // (nSrc*nOrient, ncomp)
649
650 // Build R^{1/2}-weighted eigen_leads: weightedLeads[i,:] = R_sqrt[i] * V[i,:]
651 MatrixXd weightedLeads = newV; // (nSrc*nOrient, ncomp)
652 if(useScalar) {
653 for(int i = 0; i < weightedLeads.rows(); ++i) {
654 weightedLeads.row(i) *= R_sqrt_vec(i);
655 }
656 } else {
657 // For each source s: rows [s*3, s*3+3) = R_sqrt_mat[s] * V[s*3, s*3+3)
658 for(int s_idx = 0; s_idx < nSrc; ++s_idx) {
659 MatrixXd Vs = newV.middleRows(s_idx * 3, 3); // (3, ncomp)
660 weightedLeads.middleRows(s_idx * 3, 3) = R_sqrt_mat[s_idx] * Vs;
661 }
662 }
663
664 // Update inverse operator:
665 // eigen_fields: (nChan, ncomp) = U
666 // eigen_leads: (nSrc*nOrient, ncomp) = R_sqrt * V (eLORETA-weighted)
667 // eigen_leads_weighted = true → prepare_inverse_operator uses K = eigenLeads * trans directly
668 inv.sing = newSing;
669 inv.eigen_fields->data = newU;
670 inv.eigen_leads->data = weightedLeads;
671 inv.eigen_leads_weighted = true;
672
673 // Recompute reginv: σ / (σ² + λ²)
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;
678 }
679 inv.reginv = reginv;
680
681 qInfo(" eLORETA inverse operator updated. [done]");
682}
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).
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.