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