v2.0.0
Loading...
Searching...
No Matches
fiff_proj.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "fiff_proj.h"
42#include "fiff_raw_data.h"
43#include "fiff_constants.h"
44#include "fiff_file.h"
45#include <stdio.h>
46#include <cmath>
47#include <utils/mnemath.h>
48
49//=============================================================================================================
50// EIGEN INCLUDES
51//=============================================================================================================
52
53#include <Eigen/SVD>
54
55//=============================================================================================================
56// USED NAMESPACES
57//=============================================================================================================
58
59using namespace UTILSLIB;
60using namespace FIFFLIB;
61using namespace Eigen;
62
63//=============================================================================================================
64// DEFINE MEMBER METHODS
65//=============================================================================================================
66
68: kind(-1)
69, active(false)
70, desc("")
72{
73
74}
75
76//=============================================================================================================
77
78FiffProj::FiffProj(const FiffProj& p_FiffProj)
79: kind(p_FiffProj.kind)
80, active(p_FiffProj.active)
81, desc(p_FiffProj.desc)
82, data(p_FiffProj.data)
83{
84
85}
86
87//=============================================================================================================
88
89FiffProj::FiffProj( fiff_int_t p_kind, bool p_active, QString p_desc, FiffNamedMatrix& p_data)
90: kind(p_kind)
91, active(p_active)
92, desc(p_desc)
93, data(new FiffNamedMatrix(p_data))
94{
95
96}
97
98//=============================================================================================================
99
101{
102
103}
104
105//=============================================================================================================
106
107void FiffProj::activate_projs(QList<FiffProj> &p_qListFiffProj)
108{
109 // Activate the projection items
110 QList<FiffProj>::Iterator it;
111 for(it = p_qListFiffProj.begin(); it != p_qListFiffProj.end(); ++it)
112 it->active = true;
113
114 printf("\t%lld projection items activated.\n", p_qListFiffProj.size());
115}
116
117//=============================================================================================================
118
119fiff_int_t FiffProj::make_projector(const QList<FiffProj>& projs, const QStringList& ch_names, MatrixXd& proj, const QStringList& bads, MatrixXd& U)
120{
121 fiff_int_t nchan = ch_names.size();
122 if (nchan == 0)
123 {
124 printf("No channel names specified\n");//ToDo throw here
125 return 0;
126 }
127
128// if(proj)
129// delete proj;
130 proj = MatrixXd::Identity(nchan,nchan);
131 fiff_int_t nproj = 0;
132 U = MatrixXd();
133
134 //
135 // Check trivial cases first
136 //
137 if (projs.size() == 0)
138 return 0;
139
140 fiff_int_t nvec = 0;
141 fiff_int_t k, l;
142 for (k = 0; k < projs.size(); ++k)
143 {
144 if (projs[k].active)
145 {
146 ++nproj;
147 nvec += projs[k].data->nrow;
148 }
149 }
150
151 if (nproj == 0) {
152 printf("FiffProj::make_projector - No projectors nproj=0\n");
153 return 0;
154 }
155
156 if (nvec <= 0) {
157 printf("FiffProj::make_projector - No rows in projector matrices found nvec<=0\n");
158 return 0;
159 }
160
161 //
162 // Pick the appropriate entries
163 //
164 MatrixXd vecs = MatrixXd::Zero(nchan,nvec);
165 nvec = 0;
166 fiff_int_t nonzero = 0;
167 qint32 p, c, i, j, v;
168 double onesize;
169 bool isBad = false;
170 RowVectorXi sel(nchan);
171 RowVectorXi vecSel(nchan);
172 sel.setConstant(-1);
173 vecSel.setConstant(-1);
174 for (k = 0; k < projs.size(); ++k)
175 {
176 if (projs[k].active)
177 {
178 FiffProj one = projs[k];
179
180 QMap<QString, int> uniqueMap;
181 for(l = 0; l < one.data->col_names.size(); ++l)
182 uniqueMap[one.data->col_names[l] ] = 0;
183
184 if (one.data->col_names.size() != uniqueMap.keys().size())
185 {
186 printf("Channel name list in projection item %d contains duplicate items",k);
187 return 0;
188 }
189
190 //
191 // Get the two selection vectors to pick correct elements from
192 // the projection vectors omitting bad channels
193 //
194 sel.resize(nchan);
195 vecSel.resize(nchan);
196 sel.setConstant(-1);
197 vecSel.setConstant(-1);
198 p = 0;
199 for (c = 0; c < nchan; ++c)
200 {
201 for (i = 0; i < one.data->col_names.size(); ++i)
202 {
203 if (QString::compare(ch_names.at(c),one.data->col_names[i]) == 0)
204 {
205 isBad = false;
206 for (j = 0; j < bads.size(); ++j)
207 {
208 if (QString::compare(ch_names.at(c),bads.at(j)) == 0)
209 {
210 isBad = true;
211 }
212 }
213
214 if (!isBad && sel[p] != c)
215 {
216 sel[p] = c;
217 vecSel[p] = i;
218 ++p;
219 }
220
221 }
222 }
223 }
224 sel.conservativeResize(p);
225 vecSel.conservativeResize(p);
226 //
227 // If there is something to pick, pickit
228 //
229 if (sel.cols() > 0)
230 for (v = 0; v < one.data->nrow; ++v)
231 for (i = 0; i < p; ++i)
232 vecs(sel[i],nvec+v) = one.data->data(v,vecSel[i]);
233
234 //
235 // Rescale for more straightforward detection of small singular values
236 //
237 for (v = 0; v < one.data->nrow; ++v)
238 {
239 onesize = sqrt((vecs.col(nvec+v).transpose()*vecs.col(nvec+v))(0,0));
240 if (onesize > 0.0)
241 {
242 vecs.col(nvec+v) = vecs.col(nvec+v)/onesize;
243 ++nonzero;
244 }
245 }
246 nvec += one.data->nrow;
247 }
248 }
249 //
250 // Check whether all of the vectors are exactly zero
251 //
252 if (nonzero == 0)
253 return 0;
254
255 //
256 // Reorthogonalize the vectors
257 //
258 JacobiSVD<MatrixXd> svd(vecs.block(0,0,vecs.rows(),nvec), ComputeFullU);
259 //Sort singular values and singular vectors
260 VectorXd S = svd.singularValues();
261 MatrixXd t_U = svd.matrixU();
262 MNEMath::sort<double>(S, t_U);
263
264 //
265 // Throw away the linearly dependent guys
266 //
267 nproj = 0;
268 for(k = 0; k < S.size(); ++k)
269 if (S[k]/S[0] > 1e-2)
270 ++nproj;
271
272 U = t_U.block(0, 0, t_U.rows(), nproj);
273
274 //
275 // Here is the celebrated result
276 //
277 proj -= U*U.transpose();
278
279 return nproj;
280}
281
282//=============================================================================================================
283
284QList<FiffProj> FiffProj::compute_from_raw(const FiffRawData &raw,
285 const MatrixXi &events,
286 int eventCode,
287 float tmin,
288 float tmax,
289 int nGrad,
290 int nMag,
291 int nEeg,
292 const QMap<QString,double> &mapReject)
293{
294 QList<FiffProj> projs;
295 float sfreq = raw.info.sfreq;
296 int nchan = raw.info.nchan;
297
298 int minSamp = static_cast<int>(std::round(tmin * sfreq));
299 int maxSamp = static_cast<int>(std::round(tmax * sfreq));
300 int ns = maxSamp - minSamp + 1;
301
302 if (ns <= 0) {
303 qWarning() << "[FiffProj::compute_from_raw] Invalid time window.";
304 return projs;
305 }
306
307 // Classify channels
308 QList<int> gradIdx, magIdx, eegIdx;
309 for (int k = 0; k < nchan; ++k) {
310 if (raw.info.bads.contains(raw.info.ch_names[k]))
311 continue;
312 if (raw.info.chs[k].kind == FIFFV_MEG_CH) {
313 if (raw.info.chs[k].unit == FIFF_UNIT_T)
314 magIdx.append(k);
315 else
316 gradIdx.append(k);
317 } else if (raw.info.chs[k].kind == FIFFV_EEG_CH) {
318 eegIdx.append(k);
319 }
320 }
321
322 // Collect matching epochs
323 QList<MatrixXd> epochs;
324 double gradReject = mapReject.value("grad", 0.0);
325 double magReject = mapReject.value("mag", 0.0);
326 double eegReject = mapReject.value("eeg", 0.0);
327
328 for (int k = 0; k < events.rows(); ++k) {
329 if (events(k, 1) != 0 || events(k, 2) != eventCode)
330 continue;
331
332 int evSample = events(k, 0);
333 int epochStart = evSample + minSamp;
334 int epochEnd = evSample + maxSamp;
335
336 if (epochStart < raw.first_samp || epochEnd > raw.last_samp)
337 continue;
338
339 MatrixXd epochData, epochTimes;
340 if (!raw.read_raw_segment(epochData, epochTimes, epochStart, epochEnd))
341 continue;
342
343 // Simple peak-to-peak rejection
344 bool ok = true;
345 for (int c = 0; c < nchan && ok; ++c) {
346 if (raw.info.bads.contains(raw.info.ch_names[c]))
347 continue;
348 double pp = epochData.row(c).maxCoeff() - epochData.row(c).minCoeff();
349 if (raw.info.chs[c].kind == FIFFV_MEG_CH) {
350 if (raw.info.chs[c].unit == FIFF_UNIT_T && magReject > 0 && pp > magReject)
351 ok = false;
352 else if (raw.info.chs[c].unit != FIFF_UNIT_T && gradReject > 0 && pp > gradReject)
353 ok = false;
354 } else if (raw.info.chs[c].kind == FIFFV_EEG_CH && eegReject > 0 && pp > eegReject) {
355 ok = false;
356 }
357 }
358 if (!ok) continue;
359
360 epochs.append(epochData);
361 }
362
363 if (epochs.isEmpty()) {
364 qWarning() << "[FiffProj::compute_from_raw] No valid epochs found for event" << eventCode;
365 return projs;
366 }
367
368 qInfo() << "[FiffProj::compute_from_raw]" << epochs.size() << "epochs collected for event" << eventCode;
369
370 // Lambda: compute SVD-based projectors for a channel subset
371 auto computeProjForChannels = [&](const QList<int> &chIdx, int nVec, const QString &desc) {
372 if (nVec <= 0 || chIdx.isEmpty())
373 return;
374
375 int nRows = epochs.size() * ns;
376 MatrixXd dataMat(nRows, chIdx.size());
377
378 for (int e = 0; e < epochs.size(); ++e) {
379 for (int c = 0; c < chIdx.size(); ++c) {
380 dataMat.block(e * ns, c, ns, 1) = epochs[e].row(chIdx[c]).transpose();
381 }
382 }
383
384 // Remove column mean
385 VectorXd colMean = dataMat.colwise().mean();
386 dataMat.rowwise() -= colMean.transpose();
387
388 // SVD
389 Eigen::JacobiSVD<MatrixXd> svd(dataMat, Eigen::ComputeThinV);
390 MatrixXd V = svd.matrixV();
391
392 int nComp = qMin(nVec, static_cast<int>(V.cols()));
393
394 for (int v = 0; v < nComp; ++v) {
395 FiffProj proj;
397 proj.active = false;
398 proj.desc = QString("%1-v%2").arg(desc).arg(v + 1);
399
400 FiffNamedMatrix::SDPtr namedMatrix(new FiffNamedMatrix());
401 namedMatrix->nrow = 1;
402 namedMatrix->ncol = nchan;
403 namedMatrix->row_names.clear();
404 namedMatrix->col_names = raw.info.ch_names;
405 namedMatrix->data = MatrixXd::Zero(1, nchan);
406
407 for (int c = 0; c < chIdx.size(); ++c) {
408 namedMatrix->data(0, chIdx[c]) = V(c, v);
409 }
410
411 proj.data = namedMatrix;
412 projs.append(proj);
413 }
414
415 qInfo() << "[FiffProj::compute_from_raw] Created" << nComp << desc << "projection vector(s)";
416 };
417
418 computeProjForChannels(gradIdx, nGrad, "PCA-grad");
419 computeProjForChannels(magIdx, nMag, "PCA-mag");
420 computeProjForChannels(eegIdx, nEeg, "PCA-eeg");
421
422 return projs;
423}
Fiff constants.
#define FIFFV_EEG_CH
#define FIFFV_MEG_CH
#define FIFF_UNIT_T
FiffRawData class declaration.
FiffProj class declaration.
Header file describing the numerical values used in fif files.
#define FIFFV_PROJ_ITEM_FIELD
Definition fiff_file.h:823
MNEMath class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Definition buildinfo.h:45
QList< FiffChInfo > chs
QSharedDataPointer< FiffNamedMatrix > SDPtr
FiffNamedMatrix::SDPtr data
Definition fiff_proj.h:192
fiff_int_t kind
Definition fiff_proj.h:188
static QList< FiffProj > compute_from_raw(const FiffRawData &raw, const Eigen::MatrixXi &events, int eventCode, float tmin, float tmax, int nGrad, int nMag, int nEeg, const QMap< QString, double > &mapReject=QMap< QString, double >())
static fiff_int_t make_projector(const QList< FiffProj > &projs, const QStringList &ch_names, Eigen::MatrixXd &proj, const QStringList &bads=defaultQStringList, Eigen::MatrixXd &U=defaultMatrixXd)
static void activate_projs(QList< FiffProj > &p_qListFiffProj)
FIFF raw measurement data.
bool read_raw_segment(Eigen::MatrixXd &data, Eigen::MatrixXd &times, fiff_int_t from=-1, fiff_int_t to=-1, const Eigen::RowVectorXi &sel=defaultRowVectorXi, bool do_debug=false) const
static Eigen::VectorXi sort(Eigen::Matrix< T, Eigen::Dynamic, 1 > &v, bool desc=true)
Definition mnemath.h:499