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
146: coding(0)
147, m(0)
148, n(0)
149, nz(0)
150{
151}
152
153//=============================================================================================================
154
156{
157 return fiff_get_matrix_dims(tag);
158}
159
160//=============================================================================================================
161
163{
164 int m,n,nz;
165 int coding,correct_size;
166
167 if ( fiff_type_fundamental(tag->type) != FIFFT_MATRIX ||
168 fiff_type_base(tag->type) != FIFFT_FLOAT ||
169 (fiff_type_matrix_coding(tag->type) != FIFFTS_MC_CCS &&
170 fiff_type_matrix_coding(tag->type) != FIFFTS_MC_RCS) ) {
171 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] wrong data type!");
172 return nullptr;
173 }
174
175 auto dims = fiff_get_matrix_sparse_dims(tag);
176 if (dims.empty())
177 return nullptr;
178
179 if (dims[0] != 2) {
180 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] wrong # of dimensions!");
181 return nullptr;
182 }
183
184 m = dims[1];
185 n = dims[2];
186 nz = dims[3];
187
188 coding = fiff_type_matrix_coding(tag->type);
189 if (coding == FIFFTS_MC_CCS)
190 correct_size = nz*(sizeof(fiff_float_t) + sizeof(fiff_int_t)) +
191 (n+1+dims[0]+2)*(sizeof(fiff_int_t));
192 else if (coding == FIFFTS_MC_RCS)
193 correct_size = nz*(sizeof(fiff_float_t) + sizeof(fiff_int_t)) +
194 (m+1+dims[0]+2)*(sizeof(fiff_int_t));
195 else {
196 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] Incomprehensible sparse matrix coding");
197 return nullptr;
198 }
199 if (tag->size() != correct_size) {
200 qWarning("[FiffSparseMatrix::fiff_get_float_sparse_matrix] wrong data size!");
201 return nullptr;
202 }
203 /*
204 * Parse tag data into Eigen vectors via Map (zero-copy wrap + assignment copy)
205 */
206 auto res = std::make_unique<FiffSparseMatrix>();
207 res->m = m;
208 res->n = n;
209 res->nz = nz;
210 res->coding = coding;
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 int ptrs_count = (coding == FIFFTS_MC_CCS) ? (n + 1) : (m + 1);
216
217 res->data = Eigen::Map<const Eigen::VectorXf>(src_data, nz);
218 res->inds = Eigen::Map<const Eigen::VectorXi>(src_inds, nz);
219 res->ptrs = Eigen::Map<const Eigen::VectorXi>(src_ptrs, ptrs_count);
220
221 return res;
222}
223
224//=============================================================================================================
225
226FiffSparseMatrix::UPtr FiffSparseMatrix::create_sparse_rcs(int nrow, int ncol, int *nnz, int **colindex, float **vals)
227{
228 int j,k,nz,ptr,ind;
229
230 for (j = 0, nz = 0; j < nrow; j++)
231 nz = nz + nnz[j];
232
233 if (nz <= 0) {
234 qWarning("[FiffSparseMatrix::create_sparse_rcs] No nonzero elements specified.");
235 return nullptr;
236 }
237
238 auto sparse = std::make_unique<FiffSparseMatrix>();
239 sparse->coding = FIFFTS_MC_RCS;
240 sparse->m = nrow;
241 sparse->n = ncol;
242 sparse->nz = nz;
243 sparse->data = Eigen::VectorXf::Zero(nz);
244 sparse->inds = Eigen::VectorXi::Zero(nz);
245 sparse->ptrs = Eigen::VectorXi::Zero(nrow + 1);
246
247 for (j = 0, nz = 0; j < nrow; j++) {
248 ptr = -1;
249 for (k = 0; k < nnz[j]; k++) {
250 if (ptr < 0)
251 ptr = nz;
252 ind = sparse->inds[nz] = colindex[j][k];
253 if (ind < 0 || ind >= ncol) {
254 qWarning("[FiffSparseMatrix::create_sparse_rcs] Column index out of range");
255 return nullptr;
256 }
257 if (vals)
258 sparse->data[nz] = vals[j][k];
259 else
260 sparse->data[nz] = 0.0;
261 nz++;
262 }
263 sparse->ptrs[j] = ptr;
264 }
265 sparse->ptrs[nrow] = nz;
266 for (j = nrow-1; j >= 0; j--) /* Take care of the empty rows */
267 if (sparse->ptrs[j] < 0)
268 sparse->ptrs[j] = sparse->ptrs[j+1];
269 return sparse;
270}
271
272//=============================================================================================================
273
275/*
276 * Fill in upper triangle with the lower triangle values
277 */
278{
279 if (this->coding != FIFFTS_MC_RCS) {
280 qWarning("[FiffSparseMatrix::mne_add_upper_triangle_rcs] input must be in RCS format");
281 return nullptr;
282 }
283 if (this->m != this->n) {
284 qWarning("[FiffSparseMatrix::mne_add_upper_triangle_rcs] input must be square");
285 return nullptr;
286 }
287
288 // Build per-row data from existing lower triangle
289 std::vector<int> nnz_vec(this->m);
290 std::vector<std::vector<int>> colindex(this->m);
291 std::vector<std::vector<float>> vals(this->m);
292
293 for (int i = 0; i < this->m; i++) {
294 nnz_vec[i] = this->ptrs[i+1] - this->ptrs[i];
295 if (nnz_vec[i] > 0) {
296 colindex[i].resize(nnz_vec[i]);
297 vals[i].resize(nnz_vec[i]);
298 for (int j = this->ptrs[i], k = 0; j < this->ptrs[i+1]; j++, k++) {
299 vals[i][k] = this->data[j];
300 colindex[i][k] = this->inds[j];
301 }
302 }
303 }
304
305 // Count additional upper-triangle entries per row
306 std::vector<int> nadd(this->m, 0);
307 for (int i = 0; i < this->m; i++)
308 for (int j = this->ptrs[i]; j < this->ptrs[i+1]; j++)
309 nadd[this->inds[j]]++;
310
311 // Expand per-row storage and add upper-triangle entries
312 for (int i = 0; i < this->m; i++) {
313 colindex[i].resize(nnz_vec[i] + nadd[i]);
314 vals[i].resize(nnz_vec[i] + nadd[i]);
315 }
316 for (int i = 0; i < this->m; i++)
317 for (int j = this->ptrs[i]; j < this->ptrs[i+1]; j++) {
318 int row = this->inds[j];
319 colindex[row][nnz_vec[row]] = i;
320 vals[row][nnz_vec[row]] = this->data[j];
321 nnz_vec[row]++;
322 }
323
324 // Build raw pointer arrays for create_sparse_rcs
325 std::vector<int*> ci_ptrs(this->m);
326 std::vector<float*> val_ptrs(this->m);
327 for (int i = 0; i < this->m; i++) {
328 ci_ptrs[i] = colindex[i].data();
329 val_ptrs[i] = vals[i].data();
330 }
331
332 return create_sparse_rcs(this->m, this->n, nnz_vec.data(), ci_ptrs.data(), val_ptrs.data());
333}
334
335//=============================================================================================================
336
337Eigen::SparseMatrix<double> FiffSparseMatrix::toEigenSparse() const
338{
339 if (is_empty())
340 return Eigen::SparseMatrix<double>();
341
342 using T = Eigen::Triplet<double>;
343 std::vector<T> tripletList;
344 tripletList.reserve(nz);
345
346 if (coding == FIFFTS_MC_RCS) {
347 for (int row = 0; row < m; ++row) {
348 for (int j = ptrs[row]; j < ptrs[row + 1]; ++j) {
349 tripletList.push_back(T(row, inds[j], static_cast<double>(data[j])));
350 }
351 }
352 } else if (coding == FIFFTS_MC_CCS) {
353 for (int col = 0; col < n; ++col) {
354 for (int j = ptrs[col]; j < ptrs[col + 1]; ++j) {
355 tripletList.push_back(T(inds[j], col, static_cast<double>(data[j])));
356 }
357 }
358 } else {
359 qWarning("[FiffSparseMatrix::toEigenSparse] Unknown coding type: %d", coding);
360 return Eigen::SparseMatrix<double>();
361 }
362
363 Eigen::SparseMatrix<double> result(m, n);
364 result.setFromTriplets(tripletList.begin(), tripletList.end());
365 return result;
366}
367
368//=============================================================================================================
369
370FiffSparseMatrix FiffSparseMatrix::fromEigenSparse(const Eigen::SparseMatrix<double>& mat)
371{
372 FiffSparseMatrix result;
373 if (mat.nonZeros() == 0) {
374 return result;
375 }
376
377 // Store in RCS (row-compressed) format
378 // Eigen's default SparseMatrix is column-major, so we iterate by row
379 result.coding = FIFFTS_MC_RCS;
380 result.m = static_cast<fiff_int_t>(mat.rows());
381 result.n = static_cast<fiff_int_t>(mat.cols());
382 result.nz = static_cast<fiff_int_t>(mat.nonZeros());
383 result.data = Eigen::VectorXf::Zero(result.nz);
384 result.inds = Eigen::VectorXi::Zero(result.nz);
385 result.ptrs = Eigen::VectorXi::Zero(result.m + 1);
386
387 // Collect triplets sorted by row
388 using T = Eigen::Triplet<double>;
389 std::vector<T> triplets;
390 triplets.reserve(mat.nonZeros());
391 for (int k = 0; k < mat.outerSize(); ++k) {
392 for (Eigen::SparseMatrix<double>::InnerIterator it(mat, k); it; ++it) {
393 triplets.push_back(T(it.row(), it.col(), it.value()));
394 }
395 }
396 std::sort(triplets.begin(), triplets.end(),
397 [](const T& a, const T& b) {
398 return a.row() < b.row() || (a.row() == b.row() && a.col() < b.col());
399 });
400
401 int idx = 0;
402 int row = 0;
403 result.ptrs[0] = 0;
404 for (const auto& t : triplets) {
405 while (row < t.row()) {
406 result.ptrs[++row] = idx;
407 }
408 result.data[idx] = static_cast<float>(t.value());
409 result.inds[idx] = static_cast<int>(t.col());
410 ++idx;
411 }
412 while (row < result.m) {
413 result.ptrs[++row] = idx;
414 }
415
416 return result;
417}
418
419//=============================================================================================================
420
422{
423 if (coding != FIFFTS_MC_RCS) {
424 qWarning("[FiffSparseMatrix::pickLowerTriangleRcs] input must be in RCS format");
425 return nullptr;
426 }
427 if (m != n) {
428 qWarning("[FiffSparseMatrix::pickLowerTriangleRcs] input must be square");
429 return nullptr;
430 }
431
432 std::vector<int> nnz_vec(m);
433 std::vector<std::vector<int>> colindex(m);
434 std::vector<std::vector<float>> vals(m);
435
436 for (int i = 0; i < m; i++) {
437 int count = ptrs[i+1] - ptrs[i];
438 if (count > 0) {
439 colindex[i].resize(count);
440 vals[i].resize(count);
441 int k = 0;
442 for (int j = ptrs[i]; j < ptrs[i+1]; j++) {
443 if (inds[j] <= i) {
444 vals[i][k] = data[j];
445 colindex[i][k] = inds[j];
446 k++;
447 }
448 }
449 nnz_vec[i] = k;
450 } else {
451 nnz_vec[i] = 0;
452 }
453 }
454
455 std::vector<int*> ci_ptrs(m);
456 std::vector<float*> val_ptrs(m);
457 for (int i = 0; i < m; i++) {
458 ci_ptrs[i] = colindex[i].empty() ? nullptr : colindex[i].data();
459 val_ptrs[i] = vals[i].empty() ? nullptr : vals[i].data();
460 }
461
462 return create_sparse_rcs(m, n, nnz_vec.data(), ci_ptrs.data(), val_ptrs.data());
463}
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
Old fiff_type declarations - replace them.
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()
Eigen::SparseMatrix< double > toEigenSparse() const
std::unique_ptr< FiffTag > UPtr
Definition fiff_tag.h:158