v2.0.0
Loading...
Searching...
No Matches
fiff_sparse_matrix.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "fiff_sparse_matrix.h"
42#include <fiff/fiff_file.h>
43#include <fiff/fiff_types.h>
44
45#include <vector>
46#include <QDebug>
47
48//=============================================================================================================
49// USED NAMESPACES
50//=============================================================================================================
51
52using namespace Eigen;
53using namespace FIFFLIB;
54
55//============================= fiff_type_spec.h =============================
56
57/*
58 * These return information about a fiff type.
59 */
60
62{
63 return type & FIFFTS_BASE_MASK;
64}
65
70
75
76//============================= fiff_matrix.c =============================
77
78std::vector<int> fiff_get_matrix_dims(const FiffTag::UPtr& tag)
79/*
80 * Interpret dimensions from matrix data (dense and sparse)
81 */
82{
83 int ndim;
84 int *dims;
85 unsigned int tsize = tag->size();
86 /*
87 * Initial checks
88 */
89 if (tag->data() == nullptr) {
90 qCritical("fiff_get_matrix_dims: no data available!");
91 return {};
92 }
93 if (fiff_type_fundamental(tag->type) != FIFFTS_FS_MATRIX) {
94 qCritical("fiff_get_matrix_dims: tag does not contain a matrix!");
95 return {};
96 }
97 if (tsize < sizeof(fiff_int_t)) {
98 qCritical("fiff_get_matrix_dims: too small matrix data!");
99 return {};
100 }
101 /*
102 * Get the number of dimensions and check
103 */
104 ndim = *((fiff_int_t *)((fiff_byte_t *)(tag->data())+tag->size()-sizeof(fiff_int_t)));
105 if (ndim <= 0 || ndim > FIFFC_MATRIX_MAX_DIM) {
106 qCritical("fiff_get_matrix_dims: unreasonable # of dimensions!");
107 return {};
108 }
109 if (fiff_type_matrix_coding(tag->type) == FIFFTS_MC_DENSE) {
110 if (tsize < (ndim+1)*sizeof(fiff_int_t)) {
111 qCritical("fiff_get_matrix_dims: too small matrix data!");
112 return {};
113 }
114 std::vector<int> res(ndim + 1);
115 res[0] = ndim;
116 dims = ((fiff_int_t *)((fiff_byte_t *)(tag->data())+tag->size())) - ndim - 1;
117 for (int k = 0; k < ndim; k++)
118 res[k+1] = dims[k];
119 return res;
120 }
121 else if (fiff_type_matrix_coding(tag->type) == FIFFTS_MC_CCS ||
123 if (tsize < (ndim+2)*sizeof(fiff_int_t)) {
124 qCritical("fiff_get_matrix_sparse_dims: too small matrix data!");
125 return {};
126 }
127 std::vector<int> res(ndim + 2);
128 res[0] = ndim;
129 dims = ((fiff_int_t *)((fiff_byte_t *)(tag->data())+tag->size())) - ndim - 1;
130 for (int k = 0; k < ndim; k++)
131 res[k+1] = dims[k];
132 res[ndim+1] = dims[-1];
133 return res;
134 }
135 else {
136 qCritical("fiff_get_matrix_dims: unknown matrix coding.");
137 return {};
138 }
139}
140
141//=============================================================================================================
142// DEFINE MEMBER METHODS
143//=============================================================================================================
144
149
150//=============================================================================================================
151
152FiffSparseMatrix::FiffSparseMatrix(Eigen::SparseMatrix<float>&& mat,
154: coding(coding)
155, m_eigen(std::move(mat))
156{
157}
158
159//=============================================================================================================
160
162{
163 return fiff_get_matrix_dims(tag);
164}
165
166//=============================================================================================================
167
169{
170 int m,n,nz;
171 int cod,correct_size;
172
173 if ( fiff_type_fundamental(tag->type) != FIFFT_MATRIX ||
174 fiff_type_base(tag->type) != FIFFT_FLOAT ||
175 (fiff_type_matrix_coding(tag->type) != FIFFTS_MC_CCS &&
176 fiff_type_matrix_coding(tag->type) != FIFFTS_MC_RCS) ) {
177 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] wrong data type!");
178 return nullptr;
179 }
180
181 auto dims = fiff_get_matrix_sparse_dims(tag);
182 if (dims.empty())
183 return nullptr;
184
185 if (dims[0] != 2) {
186 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] wrong # of dimensions!");
187 return nullptr;
188 }
189
190 m = dims[1];
191 n = dims[2];
192 nz = dims[3];
193
194 cod = fiff_type_matrix_coding(tag->type);
195 if (cod == FIFFTS_MC_CCS)
196 correct_size = nz*(sizeof(fiff_float_t) + sizeof(fiff_int_t)) +
197 (n+1+dims[0]+2)*(sizeof(fiff_int_t));
198 else if (cod == FIFFTS_MC_RCS)
199 correct_size = nz*(sizeof(fiff_float_t) + sizeof(fiff_int_t)) +
200 (m+1+dims[0]+2)*(sizeof(fiff_int_t));
201 else {
202 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] Incomprehensible sparse matrix coding");
203 return nullptr;
204 }
205 if (tag->size() != correct_size) {
206 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] wrong data size!");
207 return nullptr;
208 }
209 /*
210 * Parse tag data into triplets and build Eigen sparse matrix directly
211 */
212 const float* src_data = reinterpret_cast<const float*>(tag->data());
213 const int* src_inds = reinterpret_cast<const int*>(src_data + nz);
214 const int* src_ptrs = src_inds + nz;
215
216 using T = Eigen::Triplet<float>;
217 std::vector<T> triplets;
218 triplets.reserve(nz);
219
220 if (cod == FIFFTS_MC_RCS) {
221 for (int row = 0; row < m; ++row) {
222 for (int j = src_ptrs[row]; j < src_ptrs[row + 1]; ++j) {
223 triplets.push_back(T(row, src_inds[j], src_data[j]));
224 }
225 }
226 } else { // CCS
227 for (int col = 0; col < n; ++col) {
228 for (int j = src_ptrs[col]; j < src_ptrs[col + 1]; ++j) {
229 triplets.push_back(T(src_inds[j], col, src_data[j]));
230 }
231 }
232 }
233
234 Eigen::SparseMatrix<float> eigenMat(m, n);
235 eigenMat.setFromTriplets(triplets.begin(), triplets.end());
236 eigenMat.makeCompressed();
237
238 auto res = std::make_unique<FiffSparseMatrix>(std::move(eigenMat), cod);
239 return res;
240}
241
242//=============================================================================================================
243
244FiffSparseMatrix::UPtr FiffSparseMatrix::create_sparse_rcs(int nrow, int ncol, int *nnz, int **colindex, float **vals)
245{
246 int j,k,totalNz;
247
248 for (j = 0, totalNz = 0; j < nrow; j++)
249 totalNz = totalNz + nnz[j];
250
251 if (totalNz <= 0) {
252 qWarning("[FiffSparseMatrix::create_sparse_rcs] No nonzero elements specified.");
253 return nullptr;
254 }
255
256 using T = Eigen::Triplet<float>;
257 std::vector<T> triplets;
258 triplets.reserve(totalNz);
259
260 for (j = 0; j < nrow; j++) {
261 for (k = 0; k < nnz[j]; k++) {
262 int col = colindex[j][k];
263 if (col < 0 || col >= ncol) {
264 qWarning("[FiffSparseMatrix::create_sparse_rcs] Column index out of range");
265 return nullptr;
266 }
267 float val = vals ? vals[j][k] : 0.0f;
268 triplets.push_back(T(j, col, val));
269 }
270 }
271
272 Eigen::SparseMatrix<float> eigenMat(nrow, ncol);
273 eigenMat.setFromTriplets(triplets.begin(), triplets.end());
274 eigenMat.makeCompressed();
275
276 return std::make_unique<FiffSparseMatrix>(std::move(eigenMat), FIFFTS_MC_RCS);
277}
278
279//=============================================================================================================
280
282/*
283 * Fill in upper triangle with the lower triangle values
284 */
285{
286 int nRows = rows();
287 int nCols = cols();
288
289 if (nRows != nCols) {
290 qWarning("[FiffSparseMatrix::mne_add_upper_triangle_rcs] input must be square");
291 return nullptr;
292 }
293
294 // Build full (lower + upper) by adding transpose
295 Eigen::SparseMatrix<float> full = m_eigen + Eigen::SparseMatrix<float>(m_eigen.transpose());
296
297 // The diagonal was counted twice — fix by subtracting the diagonal once
298 for (int k = 0; k < full.outerSize(); ++k) {
299 for (Eigen::SparseMatrix<float>::InnerIterator it(full, k); it; ++it) {
300 if (it.row() == it.col()) {
301 // Original diagonal value is in m_eigen; full has 2x, so set back to 1x
302 it.valueRef() = 0.5f * it.value();
303 }
304 }
305 }
306
307 full.makeCompressed();
308 return std::make_unique<FiffSparseMatrix>(std::move(full), FIFFTS_MC_RCS);
309}
310
311//=============================================================================================================
312
313FiffSparseMatrix FiffSparseMatrix::fromEigenSparse(const Eigen::SparseMatrix<double>& mat)
314{
315 FiffSparseMatrix result;
316 if (mat.nonZeros() == 0)
317 return result;
318
319 result.coding = FIFFTS_MC_RCS;
320 result.m_eigen = mat.cast<float>();
321 result.m_eigen.makeCompressed();
322 return result;
323}
324
325//=============================================================================================================
326
327FiffSparseMatrix FiffSparseMatrix::fromEigenSparse(const Eigen::SparseMatrix<float>& mat)
328{
329 FiffSparseMatrix result;
330 if (mat.nonZeros() == 0)
331 return result;
332
333 result.coding = FIFFTS_MC_RCS;
334 result.m_eigen = mat;
335 result.m_eigen.makeCompressed();
336 return result;
337}
338
339//=============================================================================================================
340
342{
343 int nRows = rows();
344 int nCols = cols();
345
346 if (nRows != nCols) {
347 qWarning("[FiffSparseMatrix::pickLowerTriangleRcs] input must be square");
348 return nullptr;
349 }
350
351 using T = Eigen::Triplet<float>;
352 std::vector<T> triplets;
353 triplets.reserve(m_eigen.nonZeros());
354
355 for (int k = 0; k < m_eigen.outerSize(); ++k) {
356 for (Eigen::SparseMatrix<float>::InnerIterator it(m_eigen, k); it; ++it) {
357 if (it.row() >= it.col()) { // lower triangle including diagonal
358 triplets.push_back(T(it.row(), it.col(), it.value()));
359 }
360 }
361 }
362
363 Eigen::SparseMatrix<float> lower(nRows, nCols);
364 lower.setFromTriplets(triplets.begin(), triplets.end());
365 lower.makeCompressed();
366
367 return std::make_unique<FiffSparseMatrix>(std::move(lower), FIFFTS_MC_RCS);
368}
Old fiff_type declarations - replace them.
fiff_int_t fiff_type_base(fiff_int_t type)
std::vector< int > fiff_get_matrix_dims(const FiffTag::UPtr &tag)
fiff_int_t fiff_type_matrix_coding(fiff_int_t type)
fiff_int_t fiff_type_fundamental(fiff_int_t type)
Header file describing the numerical values used in fif files.
#define FIFFTS_FS_MATRIX
Definition fiff_file.h:266
#define FIFFTS_MC_RCS
Definition fiff_file.h:270
#define FIFFTS_MC_MASK
Definition fiff_file.h:263
#define FIFFT_MATRIX
Definition fiff_file.h:272
#define FIFFC_MATRIX_MAX_DIM
Definition fiff_file.h:259
#define FIFFTS_FS_MASK
Definition fiff_file.h:261
#define FIFFTS_MC_CCS
Definition fiff_file.h:269
#define FIFFT_FLOAT
Definition fiff_file.h:232
#define FIFFTS_MC_DENSE
Definition fiff_file.h:268
#define FIFFTS_BASE_MASK
Definition fiff_file.h:262
FiffSparseMatrix class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
float fiff_float_t
Definition fiff_types.h:93
unsigned char fiff_byte_t
Definition fiff_types.h:85
static FiffSparseMatrix fromEigenSparse(const Eigen::SparseMatrix< double > &mat)
static std::vector< int > fiff_get_matrix_sparse_dims(const FIFFLIB::FiffTag::UPtr &tag)
FiffSparseMatrix::UPtr pickLowerTriangleRcs() const
static FiffSparseMatrix::UPtr create_sparse_rcs(int nrow, int ncol, int *nnz, int **colindex, float **vals)
static FiffSparseMatrix::UPtr fiff_get_float_sparse_matrix(const FIFFLIB::FiffTag::UPtr &tag)
std::unique_ptr< FiffSparseMatrix > UPtr
FiffSparseMatrix::UPtr mne_add_upper_triangle_rcs()
std::unique_ptr< FiffTag > UPtr
Definition fiff_tag.h:158