v2.0.0
Loading...
Searching...
No Matches
inv_lcmv.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "inv_lcmv.h"
41
43#include <fiff/fiff_cov.h>
44#include <fiff/fiff_evoked.h>
45#include <fiff/fiff_info.h>
46#include <math/linalg.h>
47
48#include <QDebug>
49
50//=============================================================================================================
51// EIGEN INCLUDES
52//=============================================================================================================
53
54#include <Eigen/Dense>
55
56//=============================================================================================================
57// USED NAMESPACES
58//=============================================================================================================
59
60using namespace Eigen;
61using namespace INVLIB;
62using namespace MNELIB;
63using namespace FIFFLIB;
64using namespace UTILSLIB;
65
66//=============================================================================================================
67// DEFINE MEMBER METHODS
68//=============================================================================================================
69
71 const MNEForwardSolution &forward,
72 const FiffCov &dataCov,
73 double reg,
74 const FiffCov &noiseCov,
75 BeamformerPickOri pickOri,
76 BeamformerWeightNorm weightNorm,
77 bool reduceRank,
78 BeamformerInversion invMethod)
79{
80 InvBeamformer result;
81 result.kind = "LCMV";
82
83 // -----------------------------------------------------------------------
84 // Extract leadfield G from forward solution
85 // -----------------------------------------------------------------------
86 if(!forward.sol || forward.sol->data.size() == 0) {
87 qWarning("InvLCMV::makeLCMV - Forward solution has no gain matrix!");
88 return result;
89 }
90
91 MatrixXd G = forward.sol->data; // (n_channels, n_sources * n_orient)
92 const int nChannels = static_cast<int>(G.rows());
93 const int nOrient = (forward.source_ori == FIFFV_MNE_FREE_ORI) ? 3 : 1;
94 const int nSources = static_cast<int>(G.cols()) / nOrient;
95
96 qInfo("InvLCMV::makeLCMV - Leadfield: %d channels x %d sources (n_orient=%d)",
97 nChannels, nSources, nOrient);
98
99 // -----------------------------------------------------------------------
100 // Build whitening matrix from noise covariance
101 // -----------------------------------------------------------------------
102 MatrixXd whitener;
103 if(noiseCov.data.size() > 0) {
104 // Compute whitener from noise covariance eigendecomposition
105 // whitener = diag(1/sqrt(eig)) @ eigvec^T
106 if(noiseCov.eig.size() > 0 && noiseCov.eigvec.size() > 0) {
107 VectorXd invSqrtEig(noiseCov.eig.size());
108 for(int i = 0; i < noiseCov.eig.size(); ++i) {
109 invSqrtEig(i) = (noiseCov.eig(i) > 1e-30)
110 ? 1.0 / std::sqrt(noiseCov.eig(i))
111 : 0.0;
112 }
113 whitener = invSqrtEig.asDiagonal() * noiseCov.eigvec.transpose();
114 } else {
115 // Fallback: identity whitening
116 whitener = MatrixXd::Identity(nChannels, nChannels);
117 }
118 } else {
119 whitener = MatrixXd::Identity(nChannels, nChannels);
120 }
121
122 // -----------------------------------------------------------------------
123 // Build SSP projection matrix
124 // -----------------------------------------------------------------------
125 MatrixXd projMat = MatrixXd::Identity(nChannels, nChannels);
126 // Note: SSP projections from info.projs are typically pre-applied to
127 // the forward solution. If not, they should be applied here.
128
129 // -----------------------------------------------------------------------
130 // Whiten leadfield and data covariance
131 // G_w = whitener @ G
132 // Cm_w = whitener @ Cm @ whitener^T
133 // -----------------------------------------------------------------------
134 MatrixXd Gw = whitener * G;
135
136 MatrixXd CmData = dataCov.data;
137 if(CmData.rows() != nChannels || CmData.cols() != nChannels) {
138 qWarning("InvLCMV::makeLCMV - Data covariance dimension (%d x %d) "
139 "does not match leadfield channels (%d)!",
140 static_cast<int>(CmData.rows()), static_cast<int>(CmData.cols()), nChannels);
141 return result;
142 }
143
144 MatrixXd CmW = whitener * CmData * whitener.transpose();
145 // Ensure Hermitian (for numerical stability)
146 CmW = (CmW + CmW.transpose()) * 0.5;
147
148 // -----------------------------------------------------------------------
149 // Source normals for orientation picking
150 // -----------------------------------------------------------------------
151 MatrixX3d nn = forward.source_nn.cast<double>();
152
153 // -----------------------------------------------------------------------
154 // Compute spatial filter
155 // -----------------------------------------------------------------------
156 MatrixXd W;
157 MatrixX3d mpOri;
158
160 Gw, CmW, reg, nOrient,
161 weightNorm, pickOri, reduceRank, invMethod,
162 nn, W, mpOri);
163
164 if(!ok) {
165 qWarning("InvLCMV::makeLCMV - Beamformer computation failed!");
166 return result;
167 }
168
169 // -----------------------------------------------------------------------
170 // Populate result
171 // -----------------------------------------------------------------------
172 result.weights.push_back(W);
173 result.whitener = whitener;
174 result.proj = projMat;
175 result.chNames = forward.sol->row_names;
176 result.isFreOri = (nOrient == 3 && pickOri != BeamformerPickOri::Normal
177 && pickOri != BeamformerPickOri::MaxPower);
178 result.nSourcesTotal = nSources;
179 result.srcType = "surface";
180 result.weightNorm = weightNorm;
181 result.pickOri = pickOri;
182 result.inversion = invMethod;
183 result.reg = reg;
184 result.rank = static_cast<int>(CmW.rows());
185 result.maxPowerOri = mpOri;
186 result.sourceNn = forward.source_nn;
187
188 // Vertex indices
189 VectorXi verts(0);
190 if(forward.src.size() >= 2) {
191 verts.resize(forward.src[0].vertno.size() + forward.src[1].vertno.size());
192 verts << forward.src[0].vertno, forward.src[1].vertno;
193 } else if(forward.src.size() == 1) {
194 verts = forward.src[0].vertno;
195 }
196 result.vertices = verts;
197
198 qInfo("InvLCMV::makeLCMV - Done. Filter: %d x %d (sources=%d, orient=%d)",
199 static_cast<int>(W.rows()), static_cast<int>(W.cols()), nSources, result.nOrient());
200
201 return result;
202}
203
204//=============================================================================================================
205
206MatrixXd InvLCMV::applyFilter(const MatrixXd &data, const InvBeamformer &filters)
207{
208 // Apply projection + whitening + spatial filter
209 MatrixXd processed = data;
210
211 // Project
212 if(filters.proj.size() > 0 && filters.proj.rows() == data.rows()) {
213 processed = filters.proj * processed;
214 }
215
216 // Whiten
217 if(filters.whitener.size() > 0 && filters.whitener.rows() == processed.rows()) {
218 processed = filters.whitener * processed;
219 }
220
221 // Apply spatial filter: sol = W @ processed
222 return filters.weights[0] * processed;
223}
224
225//=============================================================================================================
226
228{
229 if(!filters.isValid() || filters.kind != "LCMV") {
230 qWarning("InvLCMV::applyLCMV - Invalid or non-LCMV filters!");
231 return InvSourceEstimate();
232 }
233
234 MatrixXd sol = applyFilter(evoked.data, filters);
235
236 // Combine XYZ for free orientation if needed
237 const int nOrient = filters.nOrient();
238 if(nOrient == 3 && filters.pickOri != BeamformerPickOri::Vector) {
239 // Combine: sqrt(x^2 + y^2 + z^2) per source per time
240 const int nSources = static_cast<int>(sol.rows()) / 3;
241 const int nTimes = static_cast<int>(sol.cols());
242 MatrixXd combined(nSources, nTimes);
243 for(int s = 0; s < nSources; ++s) {
244 combined.row(s) = sol.middleRows(s * 3, 3).colwise().norm();
245 }
246 sol = combined;
247 }
248
249 float tmin = evoked.times.size() > 0 ? evoked.times[0] : 0.0f;
250 float tstep = (evoked.info.sfreq > 0) ? 1.0f / evoked.info.sfreq : 1.0f;
251
252 InvSourceEstimate stc(sol, filters.vertices, tmin, tstep);
256
257 return stc;
258}
259
260//=============================================================================================================
261
263 float tmin,
264 float tstep,
265 const InvBeamformer &filters)
266{
267 if(!filters.isValid() || filters.kind != "LCMV") {
268 qWarning("InvLCMV::applyLCMVRaw - Invalid or non-LCMV filters!");
269 return InvSourceEstimate();
270 }
271
272 MatrixXd sol = applyFilter(data, filters);
273
274 const int nOrient = filters.nOrient();
275 if(nOrient == 3 && filters.pickOri != BeamformerPickOri::Vector) {
276 const int nSources = static_cast<int>(sol.rows()) / 3;
277 const int nTimes = static_cast<int>(sol.cols());
278 MatrixXd combined(nSources, nTimes);
279 for(int s = 0; s < nSources; ++s) {
280 combined.row(s) = sol.middleRows(s * 3, 3).colwise().norm();
281 }
282 sol = combined;
283 }
284
285 InvSourceEstimate stc(sol, filters.vertices, tmin, tstep);
289
290 return stc;
291}
292
293//=============================================================================================================
294
296 const InvBeamformer &filters)
297{
298 if(!filters.isValid() || filters.kind != "LCMV") {
299 qWarning("InvLCMV::applyLCMVCov - Invalid or non-LCMV filters!");
300 return InvSourceEstimate();
301 }
302
303 // Whiten data covariance
304 MatrixXd CmW = dataCov.data;
305 if(filters.whitener.size() > 0) {
306 CmW = filters.whitener * CmW * filters.whitener.transpose();
307 }
308
309 const int nOrient = filters.nOrient();
310 VectorXd power = InvBeamformerCompute::computePower(CmW, filters.weights[0], nOrient);
311
312 // Return as 1-column source estimate
313 MatrixXd powerMat = power; // (nSources, 1) implicitly via VectorXd
314
315 InvSourceEstimate stc(powerMat, filters.vertices, 0.0f, 1.0f);
319
320 return stc;
321}
Shared beamformer computation routines for LCMV and DICS.
InvLCMV class declaration — Linearly Constrained Minimum Variance beamformer.
FiffInfo class declaration.
#define FIFFV_MNE_FREE_ORI
FiffEvoked class declaration.
FiffCov class declaration.
MNEForwardSolution 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).
covariance data
Definition fiff_cov.h:85
Eigen::MatrixXd eigvec
Definition fiff_cov.h:259
Eigen::VectorXd eig
Definition fiff_cov.h:258
Eigen::MatrixXd data
Definition fiff_cov.h:254
Eigen::RowVectorXf times
Eigen::MatrixXd data
FIFF measurement file information.
Definition fiff_info.h:86
Computed beamformer spatial filter container.
BeamformerPickOri pickOri
BeamformerWeightNorm weightNorm
Eigen::MatrixX3f sourceNn
BeamformerInversion inversion
Eigen::MatrixX3d maxPowerOri
Eigen::VectorXi vertices
std::vector< Eigen::MatrixXd > weights
Eigen::MatrixXd whitener
static Eigen::VectorXd computePower(const Eigen::MatrixXd &Cm, const Eigen::MatrixXd &W, int nOrient)
static bool computeBeamformer(const Eigen::MatrixXd &G, const Eigen::MatrixXd &Cm, double reg, int nOrient, BeamformerWeightNorm weightNorm, BeamformerPickOri pickOri, bool reduceRank, BeamformerInversion invMethod, const Eigen::MatrixX3d &nn, Eigen::MatrixXd &W, Eigen::MatrixX3d &maxPowerOri)
static InvSourceEstimate applyLCMVCov(const FIFFLIB::FiffCov &dataCov, const InvBeamformer &filters)
Definition inv_lcmv.cpp:295
static InvBeamformer makeLCMV(const FIFFLIB::FiffInfo &info, const MNELIB::MNEForwardSolution &forward, const FIFFLIB::FiffCov &dataCov, double reg=0.05, const FIFFLIB::FiffCov &noiseCov=FIFFLIB::FiffCov(), BeamformerPickOri pickOri=BeamformerPickOri::None, BeamformerWeightNorm weightNorm=BeamformerWeightNorm::UnitNoiseGain, bool reduceRank=false, BeamformerInversion invMethod=BeamformerInversion::Matrix)
Definition inv_lcmv.cpp:70
static InvSourceEstimate applyLCMV(const FIFFLIB::FiffEvoked &evoked, const InvBeamformer &filters)
Definition inv_lcmv.cpp:227
static InvSourceEstimate applyLCMVRaw(const Eigen::MatrixXd &data, float tmin, float tstep, const InvBeamformer &filters)
Definition inv_lcmv.cpp:262
InvSourceSpaceType sourceSpaceType
InvOrientationType orientationType
MNELIB::MNESourceSpaces src
FIFFLIB::FiffNamedMatrix::SDPtr sol