v2.0.0
Loading...
Searching...
No Matches
stcloadingworker.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "stcloadingworker.h"
41
44
45#include <QFile>
46#include <QDebug>
47
48//=============================================================================================================
49// DEFINE MEMBER METHODS
50//=============================================================================================================
51
53 const QString &rhPath,
54 BrainSurface *lhSurface,
55 BrainSurface *rhSurface,
56 double cancelDist,
57 QObject *parent)
58 : QObject(parent)
59 , m_lhPath(lhPath)
60 , m_rhPath(rhPath)
61 , m_lhSurface(lhSurface)
62 , m_rhSurface(rhSurface)
63 , m_cancelDist(cancelDist)
64{
65}
66
67//=============================================================================================================
68
72
73//=============================================================================================================
74
75QSharedPointer<Eigen::SparseMatrix<float>> StcLoadingWorker::computeInterpolationMatrix(
76 const Eigen::MatrixX3f &matVertices,
77 Eigen::VectorXi &vecSourceVertices,
78 double cancelDist,
79 const QString &hemiLabel,
80 int progressStart,
81 int progressEnd)
82{
83 emit progress(progressStart, QString("Computing %1 geodesic interpolation...").arg(hemiLabel));
84
85 // Compute neighbors from the appropriate surface
86 std::vector<Eigen::VectorXi> vecNeighbors;
87 if (hemiLabel == "LH" && m_lhSurface) {
88 vecNeighbors = m_lhSurface->computeNeighbors();
89 } else if (hemiLabel == "RH" && m_rhSurface) {
90 vecNeighbors = m_rhSurface->computeNeighbors();
91 }
92
93 // Use sparse SCDC+interpolation: runs Dijkstra per source vertex and directly
94 // builds the sparse interpolation matrix without a dense intermediate distance table.
95 // This produces identical results to the old scdc() + createInterpolationMat() pipeline
96 // but uses ~18 MB instead of ~8.7 GB for a typical fsaverage source space.
97 const int progressRange = progressEnd - progressStart;
98 auto progressCallback = [this, &hemiLabel, progressStart, progressRange](int current, int total) {
99 int pct = progressStart + (progressRange * current) / total;
100 emit progress(pct, QString("%1 geodesic interpolation %2/%3").arg(hemiLabel).arg(current).arg(total));
101 };
102
104 matVertices,
105 vecNeighbors,
106 vecSourceVertices,
108 cancelDist,
109 progressCallback
110 );
111
112 emit progress(progressEnd, QString("%1 interpolation matrix ready").arg(hemiLabel));
113 return interpMat;
114}
115
116//=============================================================================================================
117
119{
120 emit progress(0, "Loading STC files...");
121
122 // Load LH STC file
123 if (!m_lhPath.isEmpty()) {
124 QFile file(m_lhPath);
126 if (INVLIB::InvSourceEstimate::read(file, stc)) {
127 m_stcLh = stc;
128 m_hasLh = true;
129 } else {
130 emit error(QString("Failed to read LH STC file: %1").arg(m_lhPath));
131 }
132 }
133
134 emit progress(5, "Loading RH STC file...");
135
136 // Load RH STC file
137 if (!m_rhPath.isEmpty()) {
138 QFile file(m_rhPath);
140 if (INVLIB::InvSourceEstimate::read(file, stc)) {
141 m_stcRh = stc;
142 m_hasRh = true;
143 } else {
144 emit error(QString("Failed to read RH STC file: %1").arg(m_rhPath));
145 }
146 }
147
148 if (!m_hasLh && !m_hasRh) {
149 emit error("No STC data could be loaded.");
150 emit finished(false);
151 return;
152 }
153
154 // Compute LH interpolation matrix
155 if (m_hasLh && m_lhSurface) {
156 Eigen::MatrixX3f matVertices = m_lhSurface->verticesAsMatrix();
157 Eigen::VectorXi vecSourceVertices = m_stcLh.vertices;
158
159 if (vecSourceVertices.size() > 0) {
160 m_interpMatLh = computeInterpolationMatrix(matVertices, vecSourceVertices, m_cancelDist,
161 "LH", 10, 45);
162 }
163 }
164
165 // Compute RH interpolation matrix
166 if (m_hasRh && m_rhSurface) {
167 Eigen::MatrixX3f matVertices = m_rhSurface->verticesAsMatrix();
168 Eigen::VectorXi vecSourceVertices = m_stcRh.vertices;
169
170 if (vecSourceVertices.size() > 0) {
171 m_interpMatRh = computeInterpolationMatrix(matVertices, vecSourceVertices, m_cancelDist,
172 "RH", 50, 90);
173 }
174 }
175
176 emit progress(95, "Finalizing...");
177 emit progress(100, "Complete");
178 emit finished(true);
179}
StcLoadingWorker class declaration.
GeometryInfo class declaration.
Interpolation class declaration.
BrainSurface class declaration.
static QSharedPointer< Eigen::SparseMatrix< float > > scdcInterpolationMat(const Eigen::MatrixX3f &matVertices, const std::vector< Eigen::VectorXi > &vecNeighborVertices, const Eigen::VectorXi &vecVertSubset, double(*interpolationFunction)(double), double dCancelDist, std::function< void(int, int)> progressCallback=nullptr)
scdcInterpolationMat Computes geodesic distances (SCDC) and builds the sparse interpolation matrix in...
static double cubic(const double dIn)
cubic Cubic hyperbola interpolation function.
Renderable cortical surface mesh with per-vertex color, curvature data, and GPU buffer management.
std::vector< Eigen::VectorXi > computeNeighbors() const
void error(const QString &message)
void progress(int percent, const QString &message)
void finished(bool success)
StcLoadingWorker(const QString &lhPath, const QString &rhPath, BrainSurface *lhSurface, BrainSurface *rhSurface, double cancelDist=0.05, QObject *parent=nullptr)
static bool read(QIODevice &p_IODevice, InvSourceEstimate &p_stc)