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