v2.0.0
Loading...
Searching...
No Matches
mne_forwardsolution.cpp
Go to the documentation of this file.
1//=============================================================================================================
37
38//=============================================================================================================
39// INCLUDES
40//=============================================================================================================
41
42#include "mne_forwardsolution.h"
43
44#include <utils/ioutils.h>
45
46//=============================================================================================================
47// EIGEN INCLUDES
48//=============================================================================================================
49
50#include <Eigen/SVD>
51#include <Eigen/Sparse>
52#include <unsupported/Eigen/KroneckerProduct>
53
54//=============================================================================================================
55// FIFF INCLUDES
56//=============================================================================================================
57
58#include <fs/colortable.h>
59#include <fs/label.h>
60#include <fs/surfaceset.h>
61#include <utils/mnemath.h>
62#include <utils/kmeans.h>
63
64#include <iostream>
65#include <QtConcurrent>
66#include <QFuture>
67
68//=============================================================================================================
69// USED NAMESPACES
70//=============================================================================================================
71
72using namespace MNELIB;
73using namespace UTILSLIB;
74using namespace FSLIB;
75using namespace Eigen;
76using namespace FIFFLIB;
77
78//=============================================================================================================
79// DEFINE MEMBER METHODS
80//=============================================================================================================
81
83: source_ori(-1)
84, surf_ori(false)
85, coord_frame(-1)
86, nsource(-1)
87, nchan(-1)
90//, mri_head_t(NULL)
91//, src(NULL)
92, source_rr(MatrixX3f::Zero(0,3))
93, source_nn(MatrixX3f::Zero(0,3))
94{
95}
96
97//=============================================================================================================
98
99MNEForwardSolution::MNEForwardSolution(QIODevice &p_IODevice, bool force_fixed, bool surf_ori, const QStringList& include, const QStringList& exclude, bool bExcludeBads)
100: source_ori(-1)
102, coord_frame(-1)
103, nsource(-1)
104, nchan(-1)
105, sol(new FiffNamedMatrix)
107//, mri_head_t(NULL)
108//, src(NULL)
109, source_rr(MatrixX3f::Zero(0,3))
110, source_nn(MatrixX3f::Zero(0,3))
111{
112 if(!read(p_IODevice, *this, force_fixed, surf_ori, include, exclude, bExcludeBads))
113 {
114 printf("\tForward solution not found.\n");//ToDo Throw here
115 return;
116 }
117}
118
119//=============================================================================================================
120
122: info(p_MNEForwardSolution.info)
123, source_ori(p_MNEForwardSolution.source_ori)
124, surf_ori(p_MNEForwardSolution.surf_ori)
125, coord_frame(p_MNEForwardSolution.coord_frame)
126, nsource(p_MNEForwardSolution.nsource)
127, nchan(p_MNEForwardSolution.nchan)
128, sol(p_MNEForwardSolution.sol)
129, sol_grad(p_MNEForwardSolution.sol_grad)
130, mri_head_t(p_MNEForwardSolution.mri_head_t)
131, src(p_MNEForwardSolution.src)
132, source_rr(p_MNEForwardSolution.source_rr)
133, source_nn(p_MNEForwardSolution.source_nn)
134{
135}
136
137//=============================================================================================================
138
142
143//=============================================================================================================
144
146{
147 info.clear();
148 source_ori = -1;
149 surf_ori = false;
150 coord_frame = -1;
151 nsource = -1;
152 nchan = -1;
155 mri_head_t.clear();
156 src.clear();
157 source_rr = MatrixX3f(0,3);
158 source_nn = MatrixX3f(0,3);
159}
160
161//=============================================================================================================
162
164 qint32 p_iClusterSize,
165 MatrixXd& p_D,
166 const FiffCov &p_pNoise_cov,
167 const FiffInfo &p_pInfo,
168 QString p_sMethod) const
169{
170 printf("Cluster forward solution using %s.\n", p_sMethod.toUtf8().constData());
171
172 MNEForwardSolution p_fwdOut = MNEForwardSolution(*this);
173
174 //Check if cov naming conventions are matching
175 if(!IOUtils::check_matching_chnames_conventions(p_pNoise_cov.names, p_pInfo.ch_names) && !p_pNoise_cov.names.isEmpty() && !p_pInfo.ch_names.isEmpty()) {
176 if(IOUtils::check_matching_chnames_conventions(p_pNoise_cov.names, p_pInfo.ch_names, true)) {
177 qWarning("MNEForwardSolution::cluster_forward_solution - Cov names do match with info channel names but have a different naming convention.");
178 //return p_fwdOut;
179 } else {
180 qWarning("MNEForwardSolution::cluster_forward_solution - Cov channel names do not match with info channel names.");
181 //return p_fwdOut;
182 }
183 }
184
185// qDebug() << "this->sol->data" << this->sol->data.rows() << "x" << this->sol->data.cols();
186
187 //
188 // Check consisty
189 //
190 if(this->isFixedOrient())
191 {
192 printf("Error: Fixed orientation not implemented jet!\n");
193 return p_fwdOut;
194 }
195
196// for(qint32 h = 0; h < this->src.hemispheres.size(); ++h )//obj.sizeForwardSolution)
197// {
198// if(this->src[h]->vertno.rows() != t_listAnnotation[h]->getLabel()->rows())
199// {
200// printf("Error: Annotation doesn't fit to Forward Solution: Vertice number is different!");
201// return false;
202// }
203// }
204
205 MatrixXd t_G_Whitened(0,0);
206 bool t_bUseWhitened = false;
207 //
208 //Whiten gain matrix before clustering -> cause diffenerent units Magnetometer, Gradiometer and EEG
209 //
210 if(!p_pNoise_cov.isEmpty() && !p_pInfo.isEmpty())
211 {
212 FiffInfo p_outFwdInfo;
213 FiffCov p_outNoiseCov;
214 MatrixXd p_outWhitener;
215 qint32 p_outNumNonZero;
216 //do whitening with noise cov
217 this->prepare_forward(p_pInfo, p_pNoise_cov, false, p_outFwdInfo, t_G_Whitened, p_outNoiseCov, p_outWhitener, p_outNumNonZero);
218 printf("\tWhitening the forward solution.\n");
219
220 t_G_Whitened = p_outWhitener*t_G_Whitened;
221 t_bUseWhitened = true;
222 }
223
224 //
225 // Assemble input data
226 //
227 qint32 count;
228 qint32 offset;
229
230 MatrixXd t_G_new;
231
232 for(qint32 h = 0; h < this->src.size(); ++h )
233 {
234
235 count = 0;
236 offset = 0;
237
238 // Offset for continuous indexing;
239 if(h > 0)
240 for(qint32 j = 0; j < h; ++j)
241 offset += this->src[j].nuse;
242
243 if(h == 0)
244 printf("Cluster Left Hemisphere\n");
245 else
246 printf("Cluster Right Hemisphere\n");
247
248 const Annotation annotation = p_AnnotationSet[h];
249 Colortable t_CurrentColorTable = annotation.getColortable();
250 VectorXi label_ids = t_CurrentColorTable.getLabelIds();
251
252 // Get label ids for every vertex
253 VectorXi vertno_labeled = VectorXi::Zero(this->src[h].vertno.rows());
254
255 //ToDo make this more universal -> using Label instead of annotations - obsolete when using Labels
256 for(qint32 i = 0; i < vertno_labeled.rows(); ++i)
257 vertno_labeled[i] = p_AnnotationSet[h].getLabelIds()[this->src[h].vertno[i]];
258
259 //Qt Concurrent List
260 QList<RegionData> m_qListRegionDataIn;
261
262 //
263 // Generate cluster input data
264 //
265 for (qint32 i = 0; i < label_ids.rows(); ++i)
266 {
267 if (label_ids[i] != 0)
268 {
269 QString curr_name = t_CurrentColorTable.struct_names[i];//obj.label2AtlasName(label(i));
270 printf("\tCluster %d / %ld %s...", i+1, label_ids.rows(), curr_name.toUtf8().constData());
271
272 //
273 // Get source space indeces
274 //
275 VectorXi idcs = VectorXi::Zero(vertno_labeled.rows());
276 qint32 c = 0;
277
278 //Select ROIs //change this use label info with a hash tabel
279 for(qint32 j = 0; j < vertno_labeled.rows(); ++j)
280 {
281 if(vertno_labeled[j] == label_ids[i])
282 {
283 idcs[c] = j;
284 ++c;
285 }
286 }
287 idcs.conservativeResize(c);
288
289 //get selected G
290 MatrixXd t_G(this->sol->data.rows(), idcs.rows()*3);
291 MatrixXd t_G_Whitened_Roi(t_G_Whitened.rows(), idcs.rows()*3);
292
293 for(qint32 j = 0; j < idcs.rows(); ++j)
294 {
295 t_G.block(0, j*3, t_G.rows(), 3) = this->sol->data.block(0, (idcs[j]+offset)*3, t_G.rows(), 3);
296 if(t_bUseWhitened)
297 t_G_Whitened_Roi.block(0, j*3, t_G_Whitened_Roi.rows(), 3) = t_G_Whitened.block(0, (idcs[j]+offset)*3, t_G_Whitened_Roi.rows(), 3);
298 }
299
300 qint32 nSens = t_G.rows();
301 qint32 nSources = t_G.cols()/3;
302
303 if (nSources > 0)
304 {
305 RegionData t_sensG;
306
307 t_sensG.idcs = idcs;
308 t_sensG.iLabelIdxIn = i;
309 t_sensG.nClusters = ceil((double)nSources/(double)p_iClusterSize);
310
311 t_sensG.matRoiGOrig = t_G;
312// if(t_bUseWhitened)
313// t_sensG.matRoiGOrigWhitened = t_G_Whitened_Roi;
314
315 printf("%d Cluster(s)... ", t_sensG.nClusters);
316
317 // Reshape Input data -> sources rows; sensors columns
318 t_sensG.matRoiG = MatrixXd(t_G.cols()/3, 3*nSens);
319 if(t_bUseWhitened)
320 t_sensG.matRoiGWhitened = MatrixXd(t_G_Whitened_Roi.cols()/3, 3*nSens);
321
322 for(qint32 j = 0; j < nSens; ++j)
323 {
324 for(qint32 k = 0; k < t_sensG.matRoiG.rows(); ++k)
325 t_sensG.matRoiG.block(k,j*3,1,3) = t_G.block(j,k*3,1,3);
326 if(t_bUseWhitened)
327 for(qint32 k = 0; k < t_sensG.matRoiGWhitened.rows(); ++k)
328 t_sensG.matRoiGWhitened.block(k,j*3,1,3) = t_G_Whitened_Roi.block(j,k*3,1,3);
329 }
330
331 t_sensG.bUseWhitened = t_bUseWhitened;
332
333 t_sensG.sDistMeasure = p_sMethod;
334
335 m_qListRegionDataIn.append(t_sensG);
336
337 printf("[added]\n");
338 }
339 else
340 {
341 printf("failed! Label contains no sources.\n");
342 }
343 }
344 }
345
346 //
347 // Calculate clusters
348 //
349 printf("Clustering... ");
350 QFuture< RegionDataOut > res;
351 res = QtConcurrent::mapped(m_qListRegionDataIn, &RegionData::cluster);
352 res.waitForFinished();
353
354 //
355 // Assign results
356 //
357 MatrixXd t_G_partial;
358
359 qint32 nClusters;
360 qint32 nSens;
361 QList<RegionData>::const_iterator itIn;
362 itIn = m_qListRegionDataIn.begin();
363 QFuture<RegionDataOut>::const_iterator itOut;
364 for (itOut = res.constBegin(); itOut != res.constEnd(); ++itOut)
365 {
366 nClusters = itOut->ctrs.rows();
367 nSens = itOut->ctrs.cols()/3;
368 t_G_partial = MatrixXd::Zero(nSens, nClusters*3);
369
370// std::cout << "Number of Clusters: " << nClusters << " x " << nSens << std::endl;//itOut->iLabelIdcsOut << std::endl;
371
372 //
373 // Assign the centroid for each cluster to the partial G
374 //
375 //ToDo change this use indeces found with whitened data
376 for(qint32 j = 0; j < nSens; ++j)
377 for(qint32 k = 0; k < nClusters; ++k)
378 t_G_partial.block(j, k*3, 1, 3) = itOut->ctrs.block(k,j*3,1,3);
379
380 //
381 // Get cluster indizes and its distances to the centroid
382 //
383 for(qint32 j = 0; j < nClusters; ++j)
384 {
385 VectorXi clusterIdcs = VectorXi::Zero(itOut->roiIdx.rows());
386 VectorXd clusterDistance = VectorXd::Zero(itOut->roiIdx.rows());
387 MatrixX3f clusterSource_rr = MatrixX3f::Zero(itOut->roiIdx.rows(), 3);
388 qint32 nClusterIdcs = 0;
389 for(qint32 k = 0; k < itOut->roiIdx.rows(); ++k)
390 {
391 if(itOut->roiIdx[k] == j)
392 {
393 clusterIdcs[nClusterIdcs] = itIn->idcs[k];
394
395 qint32 offset = h == 0 ? 0 : this->src[0].nuse;
396 clusterSource_rr.row(nClusterIdcs) = this->source_rr.row(offset + itIn->idcs[k]);
397 clusterDistance[nClusterIdcs] = itOut->D(k,j);
398 ++nClusterIdcs;
399 }
400 }
401 clusterIdcs.conservativeResize(nClusterIdcs);
402 clusterSource_rr.conservativeResize(nClusterIdcs,3);
403 clusterDistance.conservativeResize(nClusterIdcs);
404
405 VectorXi clusterVertnos = VectorXi::Zero(clusterIdcs.size());
406 for(qint32 k = 0; k < clusterVertnos.size(); ++k)
407 clusterVertnos(k) = this->src[h].vertno[clusterIdcs(k)];
408
409 p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterVertnos.append(clusterVertnos);
410 p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterSource_rr.append(clusterSource_rr);
411 p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterDistances.append(clusterDistance);
412 p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterLabelIds.append(label_ids[itOut->iLabelIdxOut]);
413 p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterLabelNames.append(t_CurrentColorTable.getNames()[itOut->iLabelIdxOut]);
414 }
415
416 //
417 // Assign partial G to new LeadField
418 //
419 if(t_G_partial.rows() > 0 && t_G_partial.cols() > 0)
420 {
421 t_G_new.conservativeResize(t_G_partial.rows(), t_G_new.cols() + t_G_partial.cols());
422 t_G_new.block(0, t_G_new.cols() - t_G_partial.cols(), t_G_new.rows(), t_G_partial.cols()) = t_G_partial;
423
424 // Map the centroids to the closest rr
425 for(qint32 k = 0; k < nClusters; ++k)
426 {
427 qint32 j = 0;
428
429 double sqec = sqrt((itIn->matRoiGOrig.block(0, j*3, itIn->matRoiGOrig.rows(), 3) - t_G_partial.block(0, k*3, t_G_partial.rows(), 3)).array().pow(2).sum());
430 double sqec_min = sqec;
431 qint32 j_min = 0;
432// MatrixXd matGainDiff;
433 for(qint32 j = 1; j < itIn->idcs.rows(); ++j)
434 {
435 sqec = sqrt((itIn->matRoiGOrig.block(0, j*3, itIn->matRoiGOrig.rows(), 3) - t_G_partial.block(0, k*3, t_G_partial.rows(), 3)).array().pow(2).sum());
436
437 if(sqec < sqec_min)
438 {
439 sqec_min = sqec;
440 j_min = j;
441// matGainDiff = itIn->matRoiGOrig.block(0, j*3, itIn->matRoiGOrig.rows(), 3) - t_G_partial.block(0, k*3, t_G_partial.rows(), 3);
442 }
443 }
444
445// qListGainDist.append(matGainDiff);
446
447 // Take the closest coordinates
448 qint32 sel_idx = itIn->idcs[j_min];
449
450 p_fwdOut.src.hemisphereAt(h)->cluster_info.centroidVertno.append(this->src[h].vertno[sel_idx]);
451 p_fwdOut.src.hemisphereAt(h)->cluster_info.centroidSource_rr.append(this->src[h].rr.row(this->src[h].vertno[sel_idx]));
452// p_fwdOut.src[h].nn.row(count) = MatrixXd::Zero(1,3);
453
454// // Option 1 closest vertno
455// p_fwdOut.src[h].vertno[count] = this->src[h].vertno[sel_idx]; //ToDo resizing necessary?
456 // Option 2 label ID
457 p_fwdOut.src[h].vertno[count] = p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterLabelIds[count];
458
459// //vertices
460// std::cout << this->src[h].vertno[sel_idx] << ", ";
461
462 ++count;
463 }
464 }
465
466 ++itIn;
467 }
468
469 //
470 // Assemble new hemisphere information
471 //
472// p_fwdOut.src[h].rr.conservativeResize(count, 3);
473// p_fwdOut.src[h].nn.conservativeResize(count, 3);
474 p_fwdOut.src[h].vertno.conservativeResize(count);
475
476 printf("[done]\n");
477 }
478
479 //
480 // Cluster operator D (sources x clusters)
481 //
482 qint32 totalNumOfClust = 0;
483 for (qint32 h = 0; h < 2; ++h)
484 totalNumOfClust += p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterVertnos.size();
485
486 if(this->isFixedOrient())
487 p_D = MatrixXd::Zero(this->sol->data.cols(), totalNumOfClust);
488 else
489 p_D = MatrixXd::Zero(this->sol->data.cols(), totalNumOfClust*3);
490
491 QList<VectorXi> t_vertnos = this->src.get_vertno();
492
493// qDebug() << "Size: " << t_vertnos[0].size() << t_vertnos[1].size();
494// qDebug() << "this->sol->data.cols(): " << this->sol->data.cols();
495
496 qint32 currentCluster = 0;
497 for (qint32 h = 0; h < 2; ++h)
498 {
499 int hemiOffset = h == 0 ? 0 : t_vertnos[0].size();
500 for(qint32 i = 0; i < p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterVertnos.size(); ++i)
501 {
502 VectorXi idx_sel;
503 MNEMath::intersect(t_vertnos[h], p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterVertnos[i], idx_sel);
504
505// std::cout << "\nVertnos:\n" << t_vertnos[h] << std::endl;
506
507// std::cout << "clusterVertnos[i]:\n" << p_fwdOut.src.hemisphereAt(h)->cluster_info.clusterVertnos[i] << std::endl;
508
509 idx_sel.array() += hemiOffset;
510
511// std::cout << "idx_sel]:\n" << idx_sel << std::endl;
512
513 double selectWeight = 1.0/idx_sel.size();
514 if(this->isFixedOrient())
515 {
516 for(qint32 j = 0; j < idx_sel.size(); ++j)
517 p_D.col(currentCluster)[idx_sel(j)] = selectWeight;
518 }
519 else
520 {
521 qint32 clustOffset = currentCluster*3;
522 for(qint32 j = 0; j < idx_sel.size(); ++j)
523 {
524 qint32 idx_sel_Offset = idx_sel(j)*3;
525 //x
526 p_D(idx_sel_Offset,clustOffset) = selectWeight;
527 //y
528 p_D(idx_sel_Offset+1, clustOffset+1) = selectWeight;
529 //z
530 p_D(idx_sel_Offset+2, clustOffset+2) = selectWeight;
531 }
532 }
533 ++currentCluster;
534 }
535 }
536
537// std::cout << "D:\n" << D.row(0) << std::endl << D.row(1) << std::endl << D.row(2) << std::endl << D.row(3) << std::endl << D.row(4) << std::endl << D.row(5) << std::endl;
538
539// //
540// // get mean and std of original gain matrix
541// //
542// double mean = 0;
543// qint32 c = 0;
544// for(qint32 i = 0; i < this->sol->data.rows(); ++i)
545// {
546// if(i % 3 == 1 || i%3 == 2)
547// {
548// for(qint32 j = 0; j < this->sol->data.cols(); ++j)
549// {
550// mean += this->sol->data(i,j);
551// ++c;
552// }
553// }
554// }
555
556// mean /= c;
557// double var = 0;
558// c = 0;
559// for(qint32 i = 0; i < this->sol->data.rows(); ++i)
560// {
561// if(i % 3 == 1 || i%3 == 2)
562// {
563// for(qint32 j = 0; j < this->sol->data.cols(); ++j)
564// {
565// var += pow(this->sol->data(i,j) - mean,2);
566// ++c;
567// }
568// }
569// }
570// var /= c - 1;
571// var = sqrt(var);
572// std::cout << "Original gain matrix (gradiometer):\n mean: " << mean << "\n var: " << var << std::endl;
573
574// mean = 0;
575// c = 0;
576// for(qint32 i = 0; i < this->sol->data.rows(); ++i)
577// {
578// if(i % 3 == 0)
579// {
580// for(qint32 j = 0; j < this->sol->data.cols(); ++j)
581// {
582// mean += this->sol->data(i,j);
583// ++c;
584// }
585// }
586// }
587// mean /= c;
588
589// var = 0;
590// c = 0;
591// for(qint32 i = 0; i < this->sol->data.rows(); ++i)
592// {
593// if(i % 3 == 0)
594// {
595// for(qint32 j = 0; j < this->sol->data.cols(); ++j)
596// {
597// var += pow(this->sol->data(i,j) - mean,2);
598// ++c;
599// }
600// }
601// }
602// var /= c - 1;
603// var = sqrt(var);
604// std::cout << "Original gain matrix (magnetometer):\n mean: " << mean << "\n var: " << var << std::endl;
605
606// //
607// // get mean and std of original gain matrix mapping
608// //
609// mean = 0;
610// c = 0;
611// for(qint32 h = 0; h < qListGainDist.size(); ++h)
612// {
613// for(qint32 i = 0; i < qListGainDist[h].rows(); ++i)
614// {
615// if(i % 3 == 1 || i%3 == 2)
616// {
617// for(qint32 j = 0; j < qListGainDist[h].cols(); ++j)
618// {
619// mean += qListGainDist[h](i,j);
620// ++c;
621// }
622// }
623// }
624// }
625// mean /= c;
626
627// var = 0;
628// c = 0;
629// for(qint32 h = 0; h < qListGainDist.size(); ++h)
630// {
631// for(qint32 i = 0; i < qListGainDist[h].rows(); ++i)
632// {
633// if(i % 3 == 1 || i%3 == 2)
634// {
635// for(qint32 j = 0; j < qListGainDist[h].cols(); ++j)
636// {
637// var += pow(qListGainDist[i](i,j) - mean,2);
638// ++c;
639// }
640// }
641// }
642// }
643
644// var /= c - 1;
645// var = sqrt(var);
646
647// std::cout << "Gain matrix offset mapping (gradiometer):\n mean: " << mean << "\n var: " << var << std::endl;
648
649// mean = 0;
650// c = 0;
651// for(qint32 h = 0; h < qListGainDist.size(); ++h)
652// {
653// for(qint32 i = 0; i < qListGainDist[h].rows(); ++i)
654// {
655// if(i % 3 == 0)
656// {
657// for(qint32 j = 0; j < qListGainDist[h].cols(); ++j)
658// {
659// mean += qListGainDist[h](i,j);
660// ++c;
661// }
662// }
663// }
664// }
665// mean /= c;
666
667// var = 0;
668// c = 0;
669// for(qint32 h = 0; h < qListGainDist.size(); ++h)
670// {
671// for(qint32 i = 0; i < qListGainDist[h].rows(); ++i)
672// {
673// if(i % 3 == 0)
674// {
675// for(qint32 j = 0; j < qListGainDist[h].cols(); ++j)
676// {
677// var += pow(qListGainDist[i](i,j) - mean,2);
678// ++c;
679// }
680// }
681// }
682// }
683
684// var /= c - 1;
685// var = sqrt(var);
686
687// std::cout << "Gain matrix offset mapping (magnetometer):\n mean: " << mean << "\n var: " << var << std::endl;
688
689 //
690 // Put it all together
691 //
692 p_fwdOut.sol->data = t_G_new;
693 p_fwdOut.sol->ncol = t_G_new.cols();
694
695 p_fwdOut.nsource = p_fwdOut.sol->ncol/3;
696
697 return p_fwdOut;
698}
699
700//=============================================================================================================
701
702MNEForwardSolution MNEForwardSolution::reduce_forward_solution(qint32 p_iNumDipoles, MatrixXd& p_D) const
703{
704 MNEForwardSolution p_fwdOut = MNEForwardSolution(*this);
705
706 bool isFixed = p_fwdOut.isFixedOrient();
707 qint32 np = isFixed ? p_fwdOut.sol->data.cols() : p_fwdOut.sol->data.cols()/3;
708
709 if(p_iNumDipoles > np)
710 return p_fwdOut;
711
712 VectorXi sel(p_iNumDipoles);
713
714 float t_fStep = (float)np/(float)p_iNumDipoles;
715
716 for(qint32 i = 0; i < p_iNumDipoles; ++i)
717 {
718 float t_fCurrent = ((float)i)*t_fStep;
719 sel[i] = (quint32)floor(t_fCurrent);
720 }
721
722 if(isFixed)
723 {
724 p_D = MatrixXd::Zero(p_fwdOut.sol->data.cols(), p_iNumDipoles);
725 for(qint32 i = 0; i < p_iNumDipoles; ++i)
726 p_D(sel[i], i) = 1;
727 }
728 else
729 {
730 p_D = MatrixXd::Zero(p_fwdOut.sol->data.cols(), p_iNumDipoles*3);
731 for(qint32 i = 0; i < p_iNumDipoles; ++i)
732 for(qint32 j = 0; j < 3; ++j)
733 p_D((sel[i]*3)+j, (i*3)+j) = 1;
734 }
735
736// //find idx of hemi switch
737// qint32 vertno_size = this->src[0].nuse;
738// qint32 hIdx = 0;
739
740// //LH
741// VectorXi vertnosLH(this->src[0].nuse);
742// for(qint32 i = 0; i < sel.size(); ++i)
743// {
744// if(sel[i] >= vertno_size)
745// {
746// hIdx = i;
747// break;
748// }
749// vertnosLH[i] = this->src[0].vertno(sel[i]);
750// }
751// vertnosLH.conservativeResize(hIdx);
752
753// QFile file_centroids_LH("./centroids_LH_sel.txt");
754// file_centroids_LH.open(QIODevice::WriteOnly | QIODevice::Text);
755// QTextStream out_centroids_LH(&file_centroids_LH);
756// for(qint32 i = 0; i < vertnosLH.size(); ++i)
757// out_centroids_LH << vertnosLH[i] << ", ";
758// file_centroids_LH.close();
759
760// QFile file_centroids_LH_full("./centroids_LH_full.txt");
761// file_centroids_LH_full.open(QIODevice::WriteOnly | QIODevice::Text);
762// QTextStream out_centroids_LH_full(&file_centroids_LH_full);
763// for(qint32 i = 0; i < this->src[0].vertno.size(); ++i)
764// out_centroids_LH_full << this->src[0].vertno[i] << ", ";
765// file_centroids_LH_full.close();
766
767// //RH
768// VectorXi vertnosRH(sel.size() - hIdx);
769// for(qint32 i = hIdx; i < sel.size(); ++i)
770// {
771// vertnosRH[i - hIdx] = this->src[1].vertno(sel[i] - hIdx);
772// }
773
774// QFile file_centroids_RH("./centroids_RH_sel.txt");
775// file_centroids_RH.open(QIODevice::WriteOnly | QIODevice::Text);
776// QTextStream out_centroids_RH(&file_centroids_RH);
777// for(qint32 i = 0; i < vertnosRH.size(); ++i)
778// out_centroids_RH << vertnosRH[i] << ", ";
779// file_centroids_RH.close();
780// //vertno end
781
782 // New gain matrix
783 p_fwdOut.sol->data = this->sol->data * p_D;
784
785 MatrixX3f rr(p_iNumDipoles,3);
786
787 MatrixX3f nn(p_iNumDipoles,3);
788
789 for(qint32 i = 0; i < p_iNumDipoles; ++i)
790 {
791 rr.row(i) = this->source_rr.row(sel(i));
792 nn.row(i) = this->source_nn.row(sel(i));
793 }
794
795 p_fwdOut.source_rr = rr;
796 p_fwdOut.source_nn = nn;
797
798 p_fwdOut.sol->ncol = p_fwdOut.sol->data.cols();
799
800 p_fwdOut.nsource = p_iNumDipoles;
801
802 return p_fwdOut;
803}
804
805//=============================================================================================================
806
807FiffCov MNEForwardSolution::compute_depth_prior(const MatrixXd &Gain, const FiffInfo &gain_info, bool is_fixed_ori, double exp, double limit, const MatrixXd &patch_areas, bool limit_depth_chs)
808{
809 printf("\tCreating the depth weighting matrix...\n");
810
811 MatrixXd G(Gain);
812 // If possible, pick best depth-weighting channels
813 if(limit_depth_chs)
815
816 VectorXd d;
817 // Compute the gain matrix
818 if(is_fixed_ori)
819 {
820 d = (G.array().square()).rowwise().sum(); //ToDo: is this correct - is G squared?
821// d = np.sum(G ** 2, axis=0)
822 }
823 else
824 {
825 qint32 n_pos = G.cols() / 3;
826 d = VectorXd::Zero(n_pos);
827 MatrixXd Gk;
828 for (qint32 k = 0; k < n_pos; ++k)
829 {
830 Gk = G.block(0,3*k, G.rows(), 3);
831 JacobiSVD<MatrixXd> svd(Gk.transpose()*Gk);
832 d[k] = svd.singularValues().maxCoeff();
833 }
834 }
835
836 // ToDo Currently the fwd solns never have "patch_areas" defined
837 if(patch_areas.cols() > 0)
838 {
839// d /= patch_areas ** 2
840 printf("\tToDo!!!!! >>> Patch areas taken into account in the depth weighting\n");
841 }
842
843 qint32 n_limit;
844 VectorXd w = d.cwiseInverse();
845 VectorXd ws = w;
846 VectorXd wpp;
847 MNEMath::sort<double>(ws, false);
848 double weight_limit = pow(limit, 2);
849 if (!limit_depth_chs)
850 {
851 // match old mne-python behavor
852 qint32 ind = 0;
853 ws.minCoeff(&ind);
854 n_limit = ind;
855 limit = ws[ind] * weight_limit;
856 }
857 else
858 {
859 // match C code behavior
860 limit = ws[ws.size()-1];
861 qint32 ind = 0;
862 n_limit = d.size();
863 if (ws[ws.size()-1] > weight_limit * ws[0])
864 {
865 double th = weight_limit * ws[0];
866 for(qint32 i = 0; i < ws.size(); ++i)
867 {
868 if(ws[i] > th)
869 {
870 ind = i;
871 break;
872 }
873 }
874 limit = ws[ind];
875 n_limit = ind;
876 }
877 }
878
879 printf("\tlimit = %d/%ld = %f", n_limit + 1, d.size(), sqrt(limit / ws[0]));
880 double scale = 1.0 / limit;
881 printf("\tscale = %g exp = %g", scale, exp);
882
883 VectorXd t_w = w.array() / limit;
884 for(qint32 i = 0; i < t_w.size(); ++i)
885 t_w[i] = t_w[i] > 1 ? 1 : t_w[i];
886 wpp = t_w.array().pow(exp);
887
888 FiffCov depth_prior;
889 if(is_fixed_ori)
890 depth_prior.data = wpp;
891 else
892 {
893 depth_prior.data.resize(wpp.rows()*3, 1);
894 qint32 idx = 0;
895 double v;
896 for(qint32 i = 0; i < wpp.rows(); ++i)
897 {
898 idx = i*3;
899 v = wpp[i];
900 depth_prior.data(idx, 0) = v;
901 depth_prior.data(idx+1, 0) = v;
902 depth_prior.data(idx+2, 0) = v;
903 }
904 }
905
906 depth_prior.kind = FIFFV_MNE_DEPTH_PRIOR_COV;
907 depth_prior.diag = true;
908 depth_prior.dim = depth_prior.data.rows();
909 depth_prior.nfree = 1;
910
911 return depth_prior;
912}
913
914//=============================================================================================================
915
917{
918 bool is_fixed_ori = this->isFixedOrient();
919 qint32 n_sources = this->sol->data.cols();
920
921 if (0 <= loose && loose <= 1)
922 {
923 qDebug() << "this->surf_ori" << this->surf_ori;
924 if(loose < 1 && !this->surf_ori)
925 {
926 printf("\tForward operator is not oriented in surface coordinates. loose parameter should be None not %f.", loose);//ToDo Throw here
927 loose = 1;
928 printf("\tSetting loose to %f.\n", loose);
929 }
930
931 if(is_fixed_ori)
932 {
933 printf("\tIgnoring loose parameter with forward operator with fixed orientation.\n");
934 loose = 0.0;
935 }
936 }
937 else
938 {
939 if(loose < 0 || loose > 1)
940 {
941 qWarning("Warning: Loose value should be in interval [0,1] not %f.\n", loose);
942 loose = loose > 1 ? 1 : 0;
943 printf("Setting loose to %f.\n", loose);
944 }
945 }
946
947 FiffCov orient_prior;
948 orient_prior.data = VectorXd::Ones(n_sources);
949 if(!is_fixed_ori && (0 <= loose && loose <= 1))
950 {
951 printf("\tApplying loose dipole orientations. Loose value of %f.\n", loose);
952 for(qint32 i = 0; i < n_sources; i+=3)
953 orient_prior.data.block(i,0,2,1).array() *= loose;
954
955 orient_prior.kind = FIFFV_MNE_ORIENT_PRIOR_COV;
956 orient_prior.diag = true;
957 orient_prior.dim = orient_prior.data.size();
958 orient_prior.nfree = 1;
959 }
960 return orient_prior;
961}
962
963//=============================================================================================================
964
966 const QStringList& exclude) const
967{
968 MNEForwardSolution fwd(*this);
969
970 if(include.size() == 0 && exclude.size() == 0)
971 return fwd;
972
973 RowVectorXi sel = FiffInfo::pick_channels(fwd.sol->row_names, include, exclude);
974
975 // Do we have something?
976 quint32 nuse = sel.size();
977
978 if (nuse == 0)
979 {
980 printf("Nothing remains after picking. Returning original forward solution.\n");
981 return fwd;
982 }
983 printf("\t%d out of %d channels remain after picking\n", nuse, fwd.nchan);
984
985 // Pick the correct rows of the forward operator
986 MatrixXd newData(nuse, fwd.sol->data.cols());
987 for(quint32 i = 0; i < nuse; ++i)
988 newData.row(i) = fwd.sol->data.row(sel[i]);
989
990 fwd.sol->data = newData;
991 fwd.sol->nrow = nuse;
992
993 QStringList ch_names;
994 for(qint32 i = 0; i < sel.cols(); ++i)
995 ch_names << fwd.sol->row_names[sel(i)];
996 fwd.nchan = nuse;
997 fwd.sol->row_names = ch_names;
998
999 QList<FiffChInfo> chs;
1000 for(qint32 i = 0; i < sel.cols(); ++i)
1001 chs.append(fwd.info.chs[sel(i)]);
1002 fwd.info.chs = chs;
1003 fwd.info.nchan = nuse;
1004
1005 QStringList bads;
1006 for(qint32 i = 0; i < fwd.info.bads.size(); ++i)
1007 if(ch_names.contains(fwd.info.bads[i]))
1008 bads.append(fwd.info.bads[i]);
1009 fwd.info.bads = bads;
1010
1011 if(!fwd.sol_grad->isEmpty())
1012 {
1013 newData.resize(nuse, fwd.sol_grad->data.cols());
1014 for(quint32 i = 0; i < nuse; ++i)
1015 newData.row(i) = fwd.sol_grad->data.row(sel[i]);
1016 fwd.sol_grad->data = newData;
1017 fwd.sol_grad->nrow = nuse;
1018 QStringList row_names;
1019 for(qint32 i = 0; i < sel.cols(); ++i)
1020 row_names << fwd.sol_grad->row_names[sel(i)];
1021 fwd.sol_grad->row_names = row_names;
1022 }
1023
1024 return fwd;
1025}
1026
1027//=============================================================================================================
1028
1029MNEForwardSolution MNEForwardSolution::pick_regions(const QList<Label> &p_qListLabels) const
1030{
1031 VectorXi selVertices;
1032
1033 qint32 iSize = 0;
1034 for(qint32 i = 0; i < p_qListLabels.size(); ++i)
1035 {
1036 VectorXi currentSelection;
1037 this->src.label_src_vertno_sel(p_qListLabels[i], currentSelection);
1038
1039 selVertices.conservativeResize(iSize+currentSelection.size());
1040 selVertices.block(iSize,0,currentSelection.size(),1) = currentSelection;
1041 iSize = selVertices.size();
1042 }
1043
1044 MNEMath::sort(selVertices, false);
1045
1046 MNEForwardSolution selectedFwd(*this);
1047
1048 MatrixX3f rr(selVertices.size(),3);
1049 MatrixX3f nn(selVertices.size(),3);
1050
1051 for(qint32 i = 0; i < selVertices.size(); ++i)
1052 {
1053 rr.block(i, 0, 1, 3) = selectedFwd.source_rr.row(selVertices[i]);
1054 nn.block(i, 0, 1, 3) = selectedFwd.source_nn.row(selVertices[i]);
1055 }
1056
1057 selectedFwd.source_rr = rr;
1058 selectedFwd.source_nn = nn;
1059
1060 VectorXi selSolIdcs = tripletSelection(selVertices);
1061 MatrixXd G(selectedFwd.sol->data.rows(),selSolIdcs.size());
1062// selectedFwd.sol_grad; //ToDo
1063 qint32 rows = G.rows();
1064
1065 for(qint32 i = 0; i < selSolIdcs.size(); ++i)
1066 G.block(0, i, rows, 1) = selectedFwd.sol->data.col(selSolIdcs[i]);
1067
1068 selectedFwd.sol->data = G;
1069 selectedFwd.sol->nrow = selectedFwd.sol->data.rows();
1070 selectedFwd.sol->ncol = selectedFwd.sol->data.cols();
1071 selectedFwd.nsource = selectedFwd.sol->ncol / 3;
1072
1073 selectedFwd.src = selectedFwd.src.pick_regions(p_qListLabels);
1074
1075 return selectedFwd;
1076}
1077
1078//=============================================================================================================
1079
1080MNEForwardSolution MNEForwardSolution::pick_types(bool meg, bool eeg, const QStringList& include, const QStringList& exclude) const
1081{
1082 RowVectorXi sel = info.pick_types(meg, eeg, false, include, exclude);
1083
1084 QStringList include_ch_names;
1085 for(qint32 i = 0; i < sel.cols(); ++i)
1086 include_ch_names << info.ch_names[sel[i]];
1087
1088 return this->pick_channels(include_ch_names);
1089}
1090
1091//=============================================================================================================
1092
1094 const FiffCov &p_noise_cov,
1095 bool p_pca,
1096 FiffInfo &p_outFwdInfo,
1097 MatrixXd &gain,
1098 FiffCov &p_outNoiseCov,
1099 MatrixXd &p_outWhitener,
1100 qint32 &p_outNumNonZero) const
1101{
1102 QStringList fwd_ch_names, ch_names;
1103 for(qint32 i = 0; i < this->info.chs.size(); ++i)
1104 fwd_ch_names << this->info.chs[i].ch_name;
1105
1106 ch_names.clear();
1107 for(qint32 i = 0; i < p_info.chs.size(); ++i)
1108 if(!p_info.bads.contains(p_info.chs[i].ch_name)
1109 && !p_noise_cov.bads.contains(p_info.chs[i].ch_name)
1110 && p_noise_cov.names.contains(p_info.chs[i].ch_name)
1111 && fwd_ch_names.contains(p_info.chs[i].ch_name))
1112 ch_names << p_info.chs[i].ch_name;
1113
1114 qint32 n_chan = ch_names.size();
1115 printf("Computing inverse operator with %d channels.\n", n_chan);
1116
1117 //
1118 // Handle noise cov
1119 //
1120 p_outNoiseCov = p_noise_cov.prepare_noise_cov(p_info, ch_names);
1121
1122 // Omit the zeroes due to projection
1123 p_outNumNonZero = 0;
1124 VectorXi t_vecNonZero = VectorXi::Zero(n_chan);
1125 for(qint32 i = 0; i < p_outNoiseCov.eig.rows(); ++i)
1126 {
1127 if(p_outNoiseCov.eig[i] > 0)
1128 {
1129 t_vecNonZero[p_outNumNonZero] = i;
1130 ++p_outNumNonZero;
1131 }
1132 }
1133 if(p_outNumNonZero > 0)
1134 t_vecNonZero.conservativeResize(p_outNumNonZero);
1135
1136 if(p_outNumNonZero > 0)
1137 {
1138 if (p_pca)
1139 {
1140 qWarning("Warning in MNEForwardSolution::prepare_forward: if (p_pca) havent been debugged.");
1141 p_outWhitener = MatrixXd::Zero(n_chan, p_outNumNonZero);
1142 // Rows of eigvec are the eigenvectors
1143 for(qint32 i = 0; i < p_outNumNonZero; ++i)
1144 p_outWhitener.col(t_vecNonZero[i]) = p_outNoiseCov.eigvec.col(t_vecNonZero[i]).array() / sqrt(p_outNoiseCov.eig(t_vecNonZero[i]));
1145 printf("\tReducing data rank to %d.\n", p_outNumNonZero);
1146 }
1147 else
1148 {
1149 printf("Creating non pca whitener.\n");
1150 p_outWhitener = MatrixXd::Zero(n_chan, n_chan);
1151 for(qint32 i = 0; i < p_outNumNonZero; ++i)
1152 p_outWhitener(t_vecNonZero[i],t_vecNonZero[i]) = 1.0 / sqrt(p_outNoiseCov.eig(t_vecNonZero[i]));
1153 // Cols of eigvec are the eigenvectors
1154 p_outWhitener *= p_outNoiseCov.eigvec;
1155 }
1156 }
1157
1158 VectorXi fwd_idx = VectorXi::Zero(ch_names.size());
1159 VectorXi info_idx = VectorXi::Zero(ch_names.size());
1160 qint32 idx;
1161 qint32 count_fwd_idx = 0;
1162 qint32 count_info_idx = 0;
1163 for(qint32 i = 0; i < ch_names.size(); ++i)
1164 {
1165 idx = fwd_ch_names.indexOf(ch_names[i]);
1166 if(idx > -1)
1167 {
1168 fwd_idx[count_fwd_idx] = idx;
1169 ++count_fwd_idx;
1170 }
1171 idx = p_info.ch_names.indexOf(ch_names[i]);
1172 if(idx > -1)
1173 {
1174 info_idx[count_info_idx] = idx;
1175 ++count_info_idx;
1176 }
1177 }
1178 fwd_idx.conservativeResize(count_fwd_idx);
1179 info_idx.conservativeResize(count_info_idx);
1180
1181 gain.resize(count_fwd_idx, this->sol->data.cols());
1182 for(qint32 i = 0; i < count_fwd_idx; ++i)
1183 gain.row(i) = this->sol->data.row(fwd_idx[i]);
1184
1185 p_outFwdInfo = p_info.pick_info(info_idx);
1186
1187 printf("\tTotal rank is %d\n", p_outNumNonZero);
1188}
1189
1190//=============================================================================================================
1191
1192bool MNEForwardSolution::read(QIODevice& p_IODevice,
1193 MNEForwardSolution& fwd,
1194 bool force_fixed,
1195 bool surf_ori,
1196 const QStringList& include,
1197 const QStringList& exclude,
1198 bool bExcludeBads)
1199{
1200 FiffStream::SPtr t_pStream(new FiffStream(&p_IODevice));
1201
1202 printf("Reading forward solution from %s...\n", t_pStream->streamName().toUtf8().constData());
1203 if(!t_pStream->open())
1204 return false;
1205 //
1206 // Find all forward solutions
1207 //
1208 QList<FiffDirNode::SPtr> fwds = t_pStream->dirtree()->dir_tree_find(FIFFB_MNE_FORWARD_SOLUTION);
1209
1210 if (fwds.size() == 0)
1211 {
1212 t_pStream->close();
1213 std::cout << "No forward solutions in " << t_pStream->streamName().toUtf8().constData(); // ToDo throw error
1214 return false;
1215 }
1216 //
1217 // Parent MRI data
1218 //
1219 QList<FiffDirNode::SPtr> parent_mri = t_pStream->dirtree()->dir_tree_find(FIFFB_MNE_PARENT_MRI_FILE);
1220 if (parent_mri.size() == 0)
1221 {
1222 t_pStream->close();
1223 std::cout << "No parent MRI information in " << t_pStream->streamName().toUtf8().constData(); // ToDo throw error
1224 return false;
1225 }
1226
1227 MNESourceSpaces t_SourceSpace;// = NULL;
1228 if(!MNESourceSpaces::readFromStream(t_pStream, true, t_SourceSpace))
1229 {
1230 t_pStream->close();
1231 std::cout << "Could not read the source spaces\n"; // ToDo throw error
1232 //ToDo error(me,'Could not read the source spaces (%s)',mne_omit_first_line(lasterr));
1233 return false;
1234 }
1235
1236 for(qint32 k = 0; k < t_SourceSpace.size(); ++k)
1237 t_SourceSpace[k].id = t_SourceSpace[k].find_source_space_hemi();
1238
1239 //
1240 // Bad channel list
1241 //
1242 QStringList bads;
1243 if(bExcludeBads)
1244 {
1245 bads = t_pStream->read_bad_channels(t_pStream->dirtree());
1246 if(bads.size() > 0)
1247 {
1248 printf("\t%lld bad channels ( ",bads.size());
1249 for(qint32 i = 0; i < bads.size(); ++i)
1250 printf("\"%s\" ", bads[i].toUtf8().constData());
1251 printf(") read\n");
1252 }
1253 }
1254
1255 //
1256 // Locate and read the forward solutions
1257 //
1258 FiffTag::SPtr t_pTag;
1259 FiffDirNode::SPtr megnode;
1260 FiffDirNode::SPtr eegnode;
1261 for(qint32 k = 0; k < fwds.size(); ++k)
1262 {
1263 if(!fwds[k]->find_tag(t_pStream, FIFF_MNE_INCLUDED_METHODS, t_pTag))
1264 {
1265 t_pStream->close();
1266 std::cout << "Methods not listed for one of the forward solutions\n"; // ToDo throw error
1267 return false;
1268 }
1269 if (*t_pTag->toInt() == FIFFV_MNE_MEG)
1270 {
1271 printf("MEG solution found\n");
1272 megnode = fwds[k];
1273 }
1274 else if(*t_pTag->toInt() == FIFFV_MNE_EEG)
1275 {
1276 printf("EEG solution found\n");
1277 eegnode = fwds[k];
1278 }
1279 }
1280
1281 MNEForwardSolution megfwd;
1282 QString ori;
1283 if (read_one(t_pStream, megnode, megfwd))
1284 {
1285 if (megfwd.source_ori == FIFFV_MNE_FIXED_ORI)
1286 ori = QString("fixed");
1287 else
1288 ori = QString("free");
1289 printf("\tRead MEG forward solution (%d sources, %d channels, %s orientations)\n", megfwd.nsource,megfwd.nchan,ori.toUtf8().constData());
1290 }
1291 MNEForwardSolution eegfwd;
1292 if (read_one(t_pStream, eegnode, eegfwd))
1293 {
1294 if (eegfwd.source_ori == FIFFV_MNE_FIXED_ORI)
1295 ori = QString("fixed");
1296 else
1297 ori = QString("free");
1298 printf("\tRead EEG forward solution (%d sources, %d channels, %s orientations)\n", eegfwd.nsource,eegfwd.nchan,ori.toUtf8().constData());
1299 }
1300
1301 //
1302 // Merge the MEG and EEG solutions together
1303 //
1304 fwd.clear();
1305
1306 if (!megfwd.isEmpty() && !eegfwd.isEmpty())
1307 {
1308 if (megfwd.sol->data.cols() != eegfwd.sol->data.cols() ||
1309 megfwd.source_ori != eegfwd.source_ori ||
1310 megfwd.nsource != eegfwd.nsource ||
1311 megfwd.coord_frame != eegfwd.coord_frame)
1312 {
1313 t_pStream->close();
1314 std::cout << "The MEG and EEG forward solutions do not match\n"; // ToDo throw error
1315 return false;
1316 }
1317
1318 fwd = MNEForwardSolution(megfwd);
1319 fwd.sol->data = MatrixXd(megfwd.sol->nrow + eegfwd.sol->nrow, megfwd.sol->ncol);
1320
1321 fwd.sol->data.block(0,0,megfwd.sol->nrow,megfwd.sol->ncol) = megfwd.sol->data;
1322 fwd.sol->data.block(megfwd.sol->nrow,0,eegfwd.sol->nrow,eegfwd.sol->ncol) = eegfwd.sol->data;
1323 fwd.sol->nrow = megfwd.sol->nrow + eegfwd.sol->nrow;
1324 fwd.sol->row_names.append(eegfwd.sol->row_names);
1325
1326 if (!fwd.sol_grad->isEmpty())
1327 {
1328 fwd.sol_grad->data.resize(megfwd.sol_grad->data.rows() + eegfwd.sol_grad->data.rows(), megfwd.sol_grad->data.cols());
1329
1330 fwd.sol->data.block(0,0,megfwd.sol_grad->data.rows(),megfwd.sol_grad->data.cols()) = megfwd.sol_grad->data;
1331 fwd.sol->data.block(megfwd.sol_grad->data.rows(),0,eegfwd.sol_grad->data.rows(),eegfwd.sol_grad->data.cols()) = eegfwd.sol_grad->data;
1332
1333 fwd.sol_grad->nrow = megfwd.sol_grad->nrow + eegfwd.sol_grad->nrow;
1334 fwd.sol_grad->row_names.append(eegfwd.sol_grad->row_names);
1335 }
1336 fwd.nchan = megfwd.nchan + eegfwd.nchan;
1337 printf("\tMEG and EEG forward solutions combined\n");
1338 }
1339 else if (!megfwd.isEmpty())
1340 fwd = megfwd; //new MNEForwardSolution(megfwd);//not copied for the sake of speed
1341 else
1342 fwd = eegfwd; //new MNEForwardSolution(eegfwd);//not copied for the sake of speed
1343
1344 //
1345 // Get the MRI <-> head coordinate transformation
1346 //
1347 if(!parent_mri[0]->find_tag(t_pStream, FIFF_COORD_TRANS, t_pTag))
1348 {
1349 t_pStream->close();
1350 std::cout << "MRI/head coordinate transformation not found\n"; // ToDo throw error
1351 return false;
1352 }
1353 else
1354 {
1355 fwd.mri_head_t = t_pTag->toCoordTrans();
1356
1358 {
1361 {
1362 t_pStream->close();
1363 std::cout << "MRI/head coordinate transformation not found\n"; // ToDo throw error
1364 return false;
1365 }
1366 }
1367 }
1368
1369 //
1370 // get parent MEG info -> from python package
1371 //
1372 t_pStream->read_meas_info_base(t_pStream->dirtree(), fwd.info);
1373
1374 t_pStream->close();
1375
1376 //
1377 // Transform the source spaces to the correct coordinate frame
1378 // if necessary
1379 //
1381 {
1382 std::cout << "Only forward solutions computed in MRI or head coordinates are acceptable";
1383 return false;
1384 }
1385
1386 //
1387 qint32 nuse = 0;
1388 t_SourceSpace.transform_source_space_to(fwd.coord_frame,fwd.mri_head_t);
1389 for(qint32 k = 0; k < t_SourceSpace.size(); ++k)
1390 nuse += t_SourceSpace[k].nuse;
1391
1392 if (nuse != fwd.nsource){
1393 qDebug() << "Source spaces do not match the forward solution.\n";
1394 return false;
1395 }
1396
1397 printf("\tSource spaces transformed to the forward solution coordinate frame\n");
1398 fwd.src = t_SourceSpace; //not new MNESourceSpaces(t_SourceSpace); for sake of speed
1399 //
1400 // Handle the source locations and orientations
1401 //
1402 if (fwd.isFixedOrient() || force_fixed == true)
1403 {
1404 nuse = 0;
1405 fwd.source_rr = MatrixXf::Zero(fwd.nsource,3);
1406 fwd.source_nn = MatrixXf::Zero(fwd.nsource,3);
1407 for(qint32 k = 0; k < t_SourceSpace.size();++k)
1408 {
1409 for(qint32 q = 0; q < t_SourceSpace[k].nuse; ++q)
1410 {
1411 fwd.source_rr.block(q,0,1,3) = t_SourceSpace[k].rr.block(t_SourceSpace[k].vertno(q),0,1,3);
1412 fwd.source_nn.block(q,0,1,3) = t_SourceSpace[k].nn.block(t_SourceSpace[k].vertno(q),0,1,3);
1413 }
1414 nuse += t_SourceSpace[k].nuse;
1415 }
1416 //
1417 // Modify the forward solution for fixed source orientations
1418 //
1420 {
1421 printf("\tChanging to fixed-orientation forward solution...");
1422
1423 MatrixXd tmp = fwd.source_nn.transpose().cast<double>();
1424 SparseMatrix<double>* fix_rot = MNEMath::make_block_diag(tmp,1);
1425 fwd.sol->data *= (*fix_rot);
1426 fwd.sol->ncol = fwd.nsource;
1428
1429 if (!fwd.sol_grad->isEmpty())
1430 {
1431 SparseMatrix<double> t_matKron;
1432 SparseMatrix<double> t_eye(3,3);
1433 for (qint32 i = 0; i < 3; ++i)
1434 t_eye.insert(i,i) = 1.0f;
1435 t_matKron = kroneckerProduct(*fix_rot,t_eye);//kron(fix_rot,eye(3));
1436 fwd.sol_grad->data *= t_matKron;
1437 fwd.sol_grad->ncol = 3*fwd.nsource;
1438 }
1439 delete fix_rot;
1440 printf("[done]\n");
1441 }
1442 }
1443 else if (surf_ori)
1444 {
1445 //
1446 // Rotate the local source coordinate systems
1447 //
1448 printf("\tConverting to surface-based source orientations...");
1449
1450 bool use_ave_nn = false;
1451 auto* hemi0 = t_SourceSpace.hemisphereAt(0);
1452 if(hemi0 && hemi0->patch_inds.size() > 0)
1453 {
1454 use_ave_nn = true;
1455 printf("\tAverage patch normals will be employed in the rotation to the local surface coordinates...\n");
1456 }
1457
1458 nuse = 0;
1459 qint32 pp = 0;
1460 fwd.source_rr = MatrixXf::Zero(fwd.nsource,3);
1461 fwd.source_nn = MatrixXf::Zero(fwd.nsource*3,3);
1462
1463 qWarning("Warning source_ori: Rotating the source coordinate system haven't been verified --> Singular Vectors U are different from MATLAB!");
1464
1465 for(qint32 k = 0; k < t_SourceSpace.size();++k)
1466 {
1467
1468 for (qint32 q = 0; q < t_SourceSpace[k].nuse; ++q)
1469 fwd.source_rr.block(q+nuse,0,1,3) = t_SourceSpace[k].rr.block(t_SourceSpace[k].vertno(q),0,1,3);
1470
1471 for (qint32 p = 0; p < t_SourceSpace[k].nuse; ++p)
1472 {
1473 //
1474 // Project out the surface normal and compute SVD
1475 //
1476 Vector3f nn;
1477 if(use_ave_nn)
1478 {
1479 auto* hemiK = t_SourceSpace.hemisphereAt(k);
1480 VectorXi t_vIdx = hemiK->pinfo[hemiK->patch_inds[p]];
1481 Matrix3Xf t_nn(3, t_vIdx.size());
1482 for(qint32 i = 0; i < t_vIdx.size(); ++i)
1483 t_nn.col(i) = t_SourceSpace[k].nn.block(t_vIdx[i],0,1,3).transpose();
1484 nn = t_nn.rowwise().sum();
1485 nn.array() /= nn.norm();
1486 }
1487 else
1488 nn = t_SourceSpace[k].nn.block(t_SourceSpace[k].vertno(p),0,1,3).transpose();
1489
1490 Matrix3f tmp = Matrix3f::Identity(nn.rows(), nn.rows()) - nn*nn.transpose();
1491
1492 JacobiSVD<MatrixXf> t_svd(tmp, Eigen::ComputeThinU);
1493 //Sort singular values and singular vectors
1494 VectorXf t_s = t_svd.singularValues();
1495 MatrixXf U = t_svd.matrixU();
1496 MNEMath::sort<float>(t_s, U);
1497
1498 //
1499 // Make sure that ez is in the direction of nn
1500 //
1501 if ((nn.transpose() * U.block(0,2,3,1))(0,0) < 0)
1502 U *= -1;
1503 fwd.source_nn.block(pp, 0, 3, 3) = U.transpose();
1504 pp += 3;
1505 }
1506 nuse += t_SourceSpace[k].nuse;
1507 }
1508 MatrixXd tmp = fwd.source_nn.transpose().cast<double>();
1509 SparseMatrix<double>* surf_rot = MNEMath::make_block_diag(tmp,3);
1510
1511 fwd.sol->data *= *surf_rot;
1512
1513 if (!fwd.sol_grad->isEmpty())
1514 {
1515 SparseMatrix<double> t_matKron;
1516 SparseMatrix<double> t_eye(3,3);
1517 for (qint32 i = 0; i < 3; ++i)
1518 t_eye.insert(i,i) = 1.0f;
1519 t_matKron = kroneckerProduct(*surf_rot,t_eye);//kron(surf_rot,eye(3));
1520 fwd.sol_grad->data *= t_matKron;
1521 }
1522 delete surf_rot;
1523 printf("[done]\n");
1524 }
1525 else
1526 {
1527 printf("\tCartesian source orientations...");
1528 nuse = 0;
1529 fwd.source_rr = MatrixXf::Zero(fwd.nsource,3);
1530 for(qint32 k = 0; k < t_SourceSpace.size(); ++k)
1531 {
1532 for (qint32 q = 0; q < t_SourceSpace[k].nuse; ++q)
1533 fwd.source_rr.block(q+nuse,0,1,3) = t_SourceSpace[k].rr.block(t_SourceSpace[k].vertno(q),0,1,3);
1534
1535 nuse += t_SourceSpace[k].nuse;
1536 }
1537
1538 MatrixXf t_ones = MatrixXf::Ones(fwd.nsource,1);
1539 Matrix3f t_eye = Matrix3f::Identity();
1540 fwd.source_nn = kroneckerProduct(t_ones,t_eye);
1541
1542 printf("[done]\n");
1543 }
1544
1545 //
1546 // Do the channel selection
1547 //
1548 QStringList exclude_bads = exclude;
1549 if (bads.size() > 0)
1550 {
1551 for(qint32 k = 0; k < bads.size(); ++k)
1552 if(!exclude_bads.contains(bads[k],Qt::CaseInsensitive))
1553 exclude_bads << bads[k];
1554 }
1555
1556 fwd.surf_ori = surf_ori;
1557 fwd = fwd.pick_channels(include, exclude_bads);
1558
1559// //
1560// // Do the channel selection - OLD VERSION
1561// //
1562// if (include.size() > 0 || exclude.size() > 0 || bads.size() > 0)
1563// {
1564// //
1565// // First do the channels to be included
1566// //
1567// RowVectorXi pick;
1568// if (include.size() == 0)
1569// pick = RowVectorXi::Ones(fwd.nchan);
1570// else
1571// {
1572// pick = RowVectorXi::Zero(fwd.nchan);
1573// for(qint32 k = 0; k < include.size(); ++k)
1574// {
1575// QList<int> c;
1576// for(qint32 l = 0; l < fwd.sol->row_names.size(); ++l)
1577// if(fwd.sol->row_names.at(l).contains(include.at(k),Qt::CaseInsensitive))
1578// pick(l) = 1;
1579
1584// }
1585// }
1586// //
1587// // Then exclude what needs to be excluded
1588// //
1589// if (exclude.size() > 0)
1590// {
1591// for(qint32 k = 0; k < exclude.size(); ++k)
1592// {
1593// QList<int> c;
1594// for(qint32 l = 0; l < fwd.sol->row_names.size(); ++l)
1595// if(fwd.sol->row_names.at(l).contains(exclude.at(k),Qt::CaseInsensitive))
1596// pick(l) = 0;
1597
1602// }
1603// }
1604// if (bads.size() > 0)
1605// {
1606// for(qint32 k = 0; k < bads.size(); ++k)
1607// {
1608// QList<int> c;
1609// for(qint32 l = 0; l < fwd.sol->row_names.size(); ++l)
1610// if(fwd.sol->row_names.at(l).contains(bads.at(k),Qt::CaseInsensitive))
1611// pick(l) = 0;
1612
1617// }
1618// }
1619// //
1620// // Do we have something?
1621// //
1622// nuse = pick.sum();
1623// if (nuse == 0)
1624// throw("Nothing remains after picking");
1625// //
1626// // Create the selector
1627// //
1628// qint32 p = 0;
1629// MatrixXd t_solData(nuse,fwd.sol->data.cols());
1630// QStringList t_solRowNames;
1631
1632// MatrixXd t_solGradData;// = NULL;
1633// QStringList t_solGradRowNames;
1634
1635// QList<FiffChInfo> chs;
1636
1637// if (!fwd.sol_grad->isEmpty())
1638// t_solGradData.resize(nuse, fwd.sol_grad->data.cols());
1639
1640// for(qint32 k = 0; k < fwd.nchan; ++k)
1641// {
1642// if(pick(k) > 0)
1643// {
1644// t_solData.block(p, 0, 1, fwd.sol->data.cols()) = fwd.sol->data.block(k, 0, 1, fwd.sol->data.cols());
1645// t_solRowNames.append(fwd.sol->row_names[k]);
1646
1647// if (!fwd.sol_grad->isEmpty())
1648// {
1649// t_solGradData.block(p, 0, 1, fwd.sol_grad->data.cols()) = fwd.sol_grad->data.block(k, 0, 1, fwd.sol_grad->data.cols());
1650// t_solGradRowNames.append(fwd.sol_grad->row_names[k]);
1651// }
1652
1653// chs.append(fwd.info.chs[k]);
1654
1655// ++p;
1656// }
1657// }
1658// printf("\t%d out of %d channels remain after picking\n",nuse,fwd.nchan);
1659// //
1660// // Pick the correct rows of the forward operator
1661// //
1662// fwd.nchan = nuse;
1665// fwd.sol->data = t_solData;
1666// fwd.sol->nrow = nuse;
1667// fwd.sol->row_names = t_solRowNames;
1668
1669// if (!fwd.sol_grad->isEmpty())
1670// {
1671// fwd.sol_grad->data = t_solGradData;
1672// fwd.sol_grad->nrow = nuse;
1673// fwd.sol_grad->row_names = t_solGradRowNames;
1674// }
1675
1676// fwd.info.chs = chs;
1677// }
1678
1679 //garbage collecting
1680 t_pStream->close();
1681
1682 return true;
1683}
1684
1685//=============================================================================================================
1686
1687bool MNEForwardSolution::read_one(FiffStream::SPtr& p_pStream,
1688 const FiffDirNode::SPtr& p_Node,
1689 MNEForwardSolution& one)
1690{
1691 //
1692 // Read all interesting stuff for one forward solution
1693 //
1694 if(!p_Node)
1695 return false;
1696
1697 one.clear();
1698 FiffTag::SPtr t_pTag;
1699
1700 if(!p_Node->find_tag(p_pStream, FIFF_MNE_SOURCE_ORIENTATION, t_pTag))
1701 {
1702 p_pStream->close();
1703 std::cout << "Source orientation tag not found."; //ToDo: throw error.
1704 return false;
1705 }
1706
1707 one.source_ori = *t_pTag->toInt();
1708
1709 if(!p_Node->find_tag(p_pStream, FIFF_MNE_COORD_FRAME, t_pTag))
1710 {
1711 p_pStream->close();
1712 std::cout << "Coordinate frame tag not found."; //ToDo: throw error.
1713 return false;
1714 }
1715
1716 one.coord_frame = *t_pTag->toInt();
1717
1718 if(!p_Node->find_tag(p_pStream, FIFF_MNE_SOURCE_SPACE_NPOINTS, t_pTag))
1719 {
1720 p_pStream->close();
1721 std::cout << "Number of sources not found."; //ToDo: throw error.
1722 return false;
1723 }
1724
1725 one.nsource = *t_pTag->toInt();
1726
1727 if(!p_Node->find_tag(p_pStream, FIFF_NCHAN, t_pTag))
1728 {
1729 p_pStream->close();
1730 printf("Number of channels not found."); //ToDo: throw error.
1731 return false;
1732 }
1733
1734 one.nchan = *t_pTag->toInt();
1735
1736 if(p_pStream->read_named_matrix(p_Node, FIFF_MNE_FORWARD_SOLUTION, *one.sol.data()))
1737 one.sol->transpose_named_matrix();
1738 else
1739 {
1740 p_pStream->close();
1741 printf("Forward solution data not found ."); //ToDo: throw error.
1742 //error(me,'Forward solution data not found (%s)',mne_omit_first_line(lasterr));
1743 return false;
1744 }
1745
1746 if(p_pStream->read_named_matrix(p_Node, FIFF_MNE_FORWARD_SOLUTION_GRAD, *one.sol_grad.data()))
1747 one.sol_grad->transpose_named_matrix();
1748 else
1749 one.sol_grad->clear();
1750
1751 if (one.sol->data.rows() != one.nchan ||
1752 (one.sol->data.cols() != one.nsource && one.sol->data.cols() != 3*one.nsource))
1753 {
1754 p_pStream->close();
1755 printf("Forward solution matrix has wrong dimensions.\n"); //ToDo: throw error.
1756 //error(me,'Forward solution matrix has wrong dimensions');
1757 return false;
1758 }
1759 if (!one.sol_grad->isEmpty())
1760 {
1761 if (one.sol_grad->data.rows() != one.nchan ||
1762 (one.sol_grad->data.cols() != 3*one.nsource && one.sol_grad->data.cols() != 3*3*one.nsource))
1763 {
1764 p_pStream->close();
1765 printf("Forward solution gradient matrix has wrong dimensions.\n"); //ToDo: throw error.
1766 //error(me,'Forward solution gradient matrix has wrong dimensions');
1767 }
1768 }
1769 return true;
1770}
1771
1772//=============================================================================================================
1773
1775{
1776 // Figure out which ones have been used
1777 if(info.chs.size() != G.rows())
1778 {
1779 printf("Error G.rows() and length of info.chs do not match: %ld != %lli", G.rows(), info.chs.size()); //ToDo throw
1780 return;
1781 }
1782
1783 RowVectorXi sel = info.pick_types(QString("grad"));
1784 if(sel.size() > 0)
1785 {
1786 for(qint32 i = 0; i < sel.size(); ++i)
1787 G.row(i) = G.row(sel[i]);
1788 G.conservativeResize(sel.size(), G.cols());
1789 printf("\t%ld planar channels", sel.size());
1790 }
1791 else
1792 {
1793 sel = info.pick_types(QString("mag"));
1794 if (sel.size() > 0)
1795 {
1796 for(qint32 i = 0; i < sel.size(); ++i)
1797 G.row(i) = G.row(sel[i]);
1798 G.conservativeResize(sel.size(), G.cols());
1799 printf("\t%ld magnetometer or axial gradiometer channels", sel.size());
1800 }
1801 else
1802 {
1803 sel = info.pick_types(false, true);
1804 if(sel.size() > 0)
1805 {
1806 for(qint32 i = 0; i < sel.size(); ++i)
1807 G.row(i) = G.row(sel[i]);
1808 G.conservativeResize(sel.size(), G.cols());
1809 printf("\t%ld EEG channels\n", sel.size());
1810 }
1811 else
1812 printf("Could not find MEG or EEG channels\n");
1813 }
1814 }
1815}
1816
1817//=============================================================================================================
1818
1820{
1821 if(!this->surf_ori || this->isFixedOrient())
1822 {
1823 qWarning("Warning: Only surface-oriented, free-orientation forward solutions can be converted to fixed orientaton.\n");//ToDo: Throw here//qCritical//qFatal
1824 return;
1825 }
1826 qint32 count = 0;
1827 for(qint32 i = 2; i < this->sol->data.cols(); i += 3)
1828 this->sol->data.col(count) = this->sol->data.col(i);//ToDo: is this right? - just take z?
1829 this->sol->data.conservativeResize(this->sol->data.rows(), count);
1830 this->sol->ncol = this->sol->ncol / 3;
1832 printf("\tConverted the forward solution into the fixed-orientation mode.\n");
1833}
1834
1835//=============================================================================================================
1836
1837MatrixX3f MNEForwardSolution::getSourcePositionsByLabel(const QList<Label> &lPickedLabels, const SurfaceSet& tSurfSetInflated)
1838{
1839 MatrixX3f matSourceVertLeft, matSourceVertRight, matSourcePositions;
1840
1841 if(lPickedLabels.isEmpty()) {
1842 qWarning() << "MNEForwardSolution::getSourcePositionsByLabel - picked label list is empty. Returning.";
1843 return matSourcePositions;
1844 }
1845
1846 if(tSurfSetInflated.isEmpty()) {
1847 qWarning() << "MNEForwardSolution::getSourcePositionsByLabel - tSurfSetInflated is empty. Returning.";
1848 return matSourcePositions;
1849 }
1850
1851 if(isClustered()) {
1852 for(int j = 0; j < this->src[0].vertno.rows(); ++j) {
1853 for(int k = 0; k < lPickedLabels.size(); k++) {
1854 if(this->src[0].vertno(j) == lPickedLabels.at(k).label_id) {
1855 matSourceVertLeft.conservativeResize(matSourceVertLeft.rows()+1,3);
1856 matSourceVertLeft.row(matSourceVertLeft.rows()-1) = tSurfSetInflated[0].rr().row(this->src.hemisphereAt(0)->cluster_info.centroidVertno.at(j)) - tSurfSetInflated[0].offset().transpose();
1857 break;
1858 }
1859 }
1860 }
1861
1862 for(int j = 0; j < this->src[1].vertno.rows(); ++j) {
1863 for(int k = 0; k < lPickedLabels.size(); k++) {
1864 if(this->src[1].vertno(j) == lPickedLabels.at(k).label_id) {
1865 matSourceVertRight.conservativeResize(matSourceVertRight.rows()+1,3);
1866 matSourceVertRight.row(matSourceVertRight.rows()-1) = tSurfSetInflated[1].rr().row(this->src.hemisphereAt(1)->cluster_info.centroidVertno.at(j)) - tSurfSetInflated[1].offset().transpose();
1867 break;
1868 }
1869 }
1870 }
1871 } else {
1872 for(int j = 0; j < this->src[0].vertno.rows(); ++j) {
1873 for(int k = 0; k < lPickedLabels.size(); k++) {
1874 for(int l = 0; l < lPickedLabels.at(k).vertices.rows(); l++) {
1875 if(this->src[0].vertno(j) == lPickedLabels.at(k).vertices(l) && lPickedLabels.at(k).hemi == 0) {
1876 matSourceVertLeft.conservativeResize(matSourceVertLeft.rows()+1,3);
1877 matSourceVertLeft.row(matSourceVertLeft.rows()-1) = tSurfSetInflated[0].rr().row(this->src[0].vertno(j)) - tSurfSetInflated[0].offset().transpose();
1878 break;
1879 }
1880 }
1881 }
1882 }
1883
1884 for(int j = 0; j < this->src[1].vertno.rows(); ++j) {
1885 for(int k = 0; k < lPickedLabels.size(); k++) {
1886 for(int l = 0; l < lPickedLabels.at(k).vertices.rows(); l++) {
1887 if(this->src[1].vertno(j) == lPickedLabels.at(k).vertices(l) && lPickedLabels.at(k).hemi == 1) {
1888 matSourceVertRight.conservativeResize(matSourceVertRight.rows()+1,3);
1889 matSourceVertRight.row(matSourceVertRight.rows()-1) = tSurfSetInflated[1].rr().row(this->src[1].vertno(j)) - tSurfSetInflated[1].offset().transpose();
1890 break;
1891 }
1892 }
1893 }
1894 }
1895 }
1896
1897 matSourcePositions.resize(matSourceVertLeft.rows()+matSourceVertRight.rows(),3);
1898 matSourcePositions << matSourceVertLeft, matSourceVertRight;
1899
1900 return matSourcePositions;
1901}
#define FIFF_MNE_COORD_FRAME
#define FIFF_MNE_FORWARD_SOLUTION_GRAD
#define FIFF_MNE_SOURCE_ORIENTATION
#define FIFF_MNE_FORWARD_SOLUTION
#define FIFF_MNE_INCLUDED_METHODS
#define FIFF_MNE_SOURCE_SPACE_NPOINTS
#define FIFFV_COORD_HEAD
#define FIFFV_MNE_FIXED_ORI
#define FIFFV_COORD_MRI
#define FIFFB_MNE_FORWARD_SOLUTION
#define FIFFV_MNE_MEG
#define FIFFV_MNE_ORIENT_PRIOR_COV
#define FIFFV_MNE_EEG
#define FIFFB_MNE_PARENT_MRI_FILE
#define FIFFV_MNE_DEPTH_PRIOR_COV
#define FIFF_NCHAN
Definition fiff_file.h:453
#define FIFF_COORD_TRANS
Definition fiff_file.h:475
IOUtils class declaration.
MNEMath class declaration.
KMeans class declaration.
MNEForwardSolution class declaration, which provides the forward solution including the source space ...
Colortable class declaration.
Label class declaration.
SurfaceSet class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FreeSurfer surface and annotation I/O.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Definition buildinfo.h:45
covariance data
Definition fiff_cov.h:84
fiff_int_t nfree
Definition fiff_cov.h:254
fiff_int_t dim
Definition fiff_cov.h:249
Eigen::MatrixXd eigvec
Definition fiff_cov.h:256
bool isEmpty() const
Definition fiff_cov.h:264
fiff_int_t kind
Definition fiff_cov.h:246
QStringList bads
Definition fiff_cov.h:253
QStringList names
Definition fiff_cov.h:250
Eigen::VectorXd eig
Definition fiff_cov.h:255
Eigen::MatrixXd data
Definition fiff_cov.h:251
FiffCov prepare_noise_cov(const FiffInfo &p_info, const QStringList &p_chNames) const
Definition fiff_cov.cpp:179
QSharedPointer< FiffDirNode > SPtr
FIFF measurement file information.
Definition fiff_info.h:85
FiffInfo pick_info(const Eigen::RowVectorXi &sel=defaultVectorXi) const
static Eigen::RowVectorXi pick_channels(const QStringList &ch_names, const QStringList &include=defaultQStringList, const QStringList &exclude=defaultQStringList)
QList< FiffChInfo > chs
QSharedDataPointer< FiffNamedMatrix > SDPtr
FIFF File I/O routines.
QSharedPointer< FiffStream > SPtr
QSharedPointer< FiffTag > SPtr
Definition fiff_tag.h:155
Free surfer annotation.
Definition annotation.h:81
Colortable & getColortable()
Definition annotation.h:329
Annotation set.
Vertices label based lookup table.
Definition colortable.h:68
QStringList getNames() const
Definition colortable.h:131
Eigen::VectorXi getLabelIds() const
Definition colortable.h:120
QStringList struct_names
Definition colortable.h:112
A hemisphere set of surfaces.
Definition surfaceset.h:72
bool isEmpty() const
Definition surfaceset.h:247
QList< Eigen::VectorXd > clusterDistances
QList< QString > clusterLabelNames
QList< qint32 > clusterLabelIds
QList< Eigen::MatrixX3f > clusterSource_rr
QList< Eigen::VectorXi > clusterVertnos
QList< qint32 > centroidVertno
QList< Eigen::Vector3f > centroidSource_rr
Input parameters for cluster-based forward solution computation on a single cortical region.
Eigen::MatrixXd matRoiGWhitened
Eigen::MatrixXd matRoiGOrig
RegionDataOut cluster() const
static void restrict_gain_matrix(Eigen::MatrixXd &G, const FIFFLIB::FiffInfo &info)
MNEForwardSolution reduce_forward_solution(qint32 p_iNumDipoles, Eigen::MatrixXd &p_D) const
static FIFFLIB::FiffCov compute_depth_prior(const Eigen::MatrixXd &Gain, const FIFFLIB::FiffInfo &gain_info, bool is_fixed_ori, double exp=0.8, double limit=10.0, const Eigen::MatrixXd &patch_areas=FIFFLIB::defaultConstMatrixXd, bool limit_depth_chs=false)
Eigen::MatrixX3f getSourcePositionsByLabel(const QList< FSLIB::Label > &lPickedLabels, const FSLIB::SurfaceSet &tSurfSetInflated)
void prepare_forward(const FIFFLIB::FiffInfo &p_info, const FIFFLIB::FiffCov &p_noise_cov, bool p_pca, FIFFLIB::FiffInfo &p_outFwdInfo, Eigen::MatrixXd &gain, FIFFLIB::FiffCov &p_outNoiseCov, Eigen::MatrixXd &p_outWhitener, qint32 &p_outNumNonZero) const
FIFFLIB::FiffCoordTrans mri_head_t
MNEForwardSolution pick_regions(const QList< FSLIB::Label > &p_qListLabels) const
MNEForwardSolution pick_channels(const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList) const
FIFFLIB::FiffNamedMatrix::SDPtr sol_grad
MNEForwardSolution cluster_forward_solution(const FSLIB::AnnotationSet &p_AnnotationSet, qint32 p_iClusterSize, Eigen::MatrixXd &p_D=defaultD, const FIFFLIB::FiffCov &p_pNoise_cov=defaultCov, const FIFFLIB::FiffInfo &p_pInfo=defaultInfo, QString p_sMethod="cityblock") const
static bool read(QIODevice &p_IODevice, MNEForwardSolution &fwd, bool force_fixed=false, bool surf_ori=false, const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList, bool bExcludeBads=true)
Eigen::VectorXi tripletSelection(const Eigen::VectorXi &p_vecIdxSelection) const
FIFFLIB::FiffCov compute_orient_prior(float loose=0.2)
MNEForwardSolution pick_types(bool meg, bool eeg, const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList) const
FIFFLIB::FiffNamedMatrix::SDPtr sol
MNEClusterInfo cluster_info
QList< Eigen::VectorXi > pinfo
Source Space descritpion.
MNESourceSpaces pick_regions(const QList< FSLIB::Label > &p_qListLabels) const
bool transform_source_space_to(FIFFLIB::fiff_int_t dest, FIFFLIB::FiffCoordTrans &trans)
MNEHemisphere * hemisphereAt(qint32 idx)
static bool readFromStream(FIFFLIB::FiffStream::SPtr &p_pStream, bool add_geom, MNESourceSpaces &p_SourceSpace)
static bool check_matching_chnames_conventions(const QStringList &chNamesA, const QStringList &chNamesB, bool bCheckForNewNamingConvention=false)
Definition ioutils.cpp:334
static Eigen::VectorXi sort(Eigen::Matrix< T, Eigen::Dynamic, 1 > &v, bool desc=true)
Definition mnemath.h:499
static Eigen::SparseMatrix< double > * make_block_diag(const Eigen::MatrixXd &A, qint32 n)
Definition mnemath.cpp:347
static Eigen::VectorXi intersect(const Eigen::VectorXi &v1, const Eigen::VectorXi &v2, Eigen::VectorXi &idx_sel)
Definition mnemath.cpp:194