Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MatrixUtils.cpp
Go to the documentation of this file.
1#include "MatrixUtils.hpp"
2
5
6#include <vector>
7
8void polyfem::utils::show_matrix_stats(const Eigen::MatrixXd &M)
9{
10 Eigen::FullPivLU<Eigen::MatrixXd> lu(M);
11 Eigen::JacobiSVD<Eigen::MatrixXd> svd(M);
12 double s1 = svd.singularValues()(0);
13 double s2 = svd.singularValues()(svd.singularValues().size() - 1);
14 double cond = s1 / s2;
15
16 logger().trace("----------------------------------------");
17 logger().trace("-- Determinant: {}", M.determinant());
18 logger().trace("-- Singular values: {} {}", s1, s2);
19 logger().trace("-- Cond: {}", cond);
20 logger().trace("-- Invertible: {}", lu.isInvertible());
21 logger().trace("----------------------------------------");
22 // logger().trace("{}", lu.solve(M) );
23}
24
25namespace
26{
27 inline bool nanproof_equals(const double x, const double y)
28 {
29 return x == y || (std::isnan(x) == std::isnan(y));
30 }
31} // namespace
32
33// Flatten rowwises
34Eigen::VectorXd polyfem::utils::flatten(const Eigen::MatrixXd &X)
35{
36 if (X.size() == 0)
37 return Eigen::VectorXd();
38
39 Eigen::VectorXd x(X.size());
40 for (int i = 0; i < X.rows(); ++i)
41 {
42 for (int j = 0; j < X.cols(); ++j)
43 {
44 x(i * X.cols() + j) = X(i, j);
45 }
46 }
47 assert(nanproof_equals(X(0, 0), x(0)));
48 assert(X.cols() <= 1 || nanproof_equals(X(0, 1), x(1)));
49 return x;
50}
51
52// Unflatten rowwises, so every dim elements in x become a row.
53Eigen::MatrixXd polyfem::utils::unflatten(const Eigen::VectorXd &x, int dim)
54{
55 if (x.size() == 0)
56 return Eigen::MatrixXd(0, dim);
57
58 assert(x.size() % dim == 0);
59 Eigen::MatrixXd X(x.size() / dim, dim);
60 for (int i = 0; i < x.size(); ++i)
61 {
62 X(i / dim, i % dim) = x(i);
63 }
64 assert(nanproof_equals(X(0, 0), x(0)));
65 assert(X.cols() <= 1 || nanproof_equals(X(0, 1), x(1)));
66 return X;
67}
68
69void polyfem::utils::vector2matrix(const Eigen::VectorXd &vec, Eigen::MatrixXd &mat)
70{
71 int size = 1;
72 if (vec.size() == 9)
73 size = 3;
74 else if (vec.size() == 4)
75 size = 2;
76 else
77 throw std::runtime_error("Invalid size in vector2matrix!");
78
79 assert(size * size == vec.size());
80
81 mat = unflatten(vec, size);
82}
83
84Eigen::SparseMatrix<double> polyfem::utils::lump_matrix(const Eigen::SparseMatrix<double> &M)
85{
86 std::vector<Eigen::Triplet<double>> triplets;
87
88 for (int k = 0; k < M.outerSize(); ++k)
89 {
90 for (Eigen::SparseMatrix<double>::InnerIterator it(M, k); it; ++it)
91 {
92 triplets.emplace_back(it.row(), it.row(), it.value());
93 }
94 }
95
96 Eigen::SparseMatrix<double> lumped(M.rows(), M.rows());
97 lumped.setFromTriplets(triplets.begin(), triplets.end());
98 lumped.makeCompressed();
99
100 return lumped;
101}
102
104 const int full_size,
105 const int reduced_size,
106 const std::vector<int> &removed_vars,
107 const StiffnessMatrix &full,
108 StiffnessMatrix &reduced)
109{
110 POLYFEM_SCOPED_TIMER("full to reduced matrix");
111
112 if (reduced_size == full_size || reduced_size == full.rows())
113 {
114 assert(reduced_size == full.rows() && reduced_size == full.cols());
115 reduced = full;
116 return;
117 }
118 assert(full.rows() == full_size && full.cols() == full_size);
119
120 Eigen::VectorXi indices(full_size);
121 int index = 0;
122 size_t kk = 0;
123 for (int i = 0; i < full_size; ++i)
124 {
125 if (kk < removed_vars.size() && removed_vars[kk] == i)
126 {
127 ++kk;
128 indices(i) = -1;
129 }
130 else
131 {
132 indices(i) = index++;
133 }
134 }
135 assert(index == reduced_size);
136
137 std::vector<Eigen::Triplet<double>> entries;
138 entries.reserve(full.nonZeros()); // Conservative estimate
139 for (int k = 0; k < full.outerSize(); ++k)
140 {
141 if (indices(k) < 0)
142 continue;
143
144 for (StiffnessMatrix::InnerIterator it(full, k); it; ++it)
145 {
146 assert(it.col() == k);
147 if (indices(it.row()) < 0 || indices(it.col()) < 0)
148 continue;
149
150 assert(indices(it.row()) >= 0);
151 assert(indices(it.col()) >= 0);
152
153 entries.emplace_back(indices(it.row()), indices(it.col()), it.value());
154 }
155 }
156
157 reduced.resize(reduced_size, reduced_size);
158 reduced.setFromTriplets(entries.begin(), entries.end());
159 reduced.makeCompressed();
160}
161
163 const Eigen::MatrixXd &in,
164 const Eigen::VectorXi &in_to_out,
165 int out_blocks,
166 const int block_size)
167{
168 constexpr double NaN = std::numeric_limits<double>::quiet_NaN();
169
170 assert(in.rows() % block_size == 0);
171 assert(in_to_out.size() == in.rows() / block_size);
172
173 if (out_blocks < 0)
174 out_blocks = in.rows() / block_size;
175
176 Eigen::MatrixXd out = Eigen::MatrixXd::Constant(
177 out_blocks * block_size, in.cols(), NaN);
178
179 const int in_blocks = in.rows() / block_size;
180 for (int i = 0; i < in_blocks; ++i)
181 {
182 const int j = in_to_out[i];
183 if (j < 0)
184 continue;
185
186 out.middleRows(block_size * j, block_size) =
187 in.middleRows(block_size * i, block_size);
188 }
189
190 return out;
191}
192
194 const Eigen::MatrixXd &out,
195 const Eigen::VectorXi &in_to_out,
196 int in_blocks,
197 const int block_size)
198{
199 constexpr double NaN = std::numeric_limits<double>::quiet_NaN();
200
201 assert(out.rows() % block_size == 0);
202
203 if (in_blocks < 0)
204 in_blocks = out.rows() / block_size;
205 assert(in_to_out.size() == in_blocks);
206
207 Eigen::MatrixXd in = Eigen::MatrixXd::Constant(
208 in_blocks * block_size, out.cols(), NaN);
209
210 for (int i = 0; i < in_blocks; i++)
211 {
212 const int j = in_to_out[i];
213 if (j < 0)
214 continue;
215
216 in.middleRows(block_size * i, block_size) =
217 out.middleRows(block_size * j, block_size);
218 }
219
220 return in;
221}
222
224 const Eigen::MatrixXi &in,
225 const Eigen::VectorXi &index_mapping)
226{
227 Eigen::MatrixXi out(in.rows(), in.cols());
228 for (int i = 0; i < in.rows(); i++)
229 {
230 for (int j = 0; j < in.cols(); j++)
231 {
232 out(i, j) = index_mapping[in(i, j)];
233 }
234 }
235 return out;
236}
237
238void polyfem::utils::scatter_matrix(const int n_dofs,
239 const int dim,
240 const Eigen::MatrixXd &A,
241 const Eigen::MatrixXd &b,
242 const std::vector<int> &local_to_global,
243 StiffnessMatrix &Aout,
244 Eigen::MatrixXd &bout)
245{
246 assert(A.rows() == b.rows());
247 std::vector<Eigen::Triplet<double>> Ae;
248
249 if (b.cols() == dim)
250 {
251 for (int i = 0; i < A.rows(); ++i)
252 {
253 for (int j = 0; j < A.cols(); ++j)
254 {
255 const auto global_j = (local_to_global.empty() ? j : local_to_global[j]) * dim;
256
257 if (A(i, j) != 0)
258 {
259 for (int d = 0; d < dim; ++d)
260 {
261 Ae.push_back(Eigen::Triplet<double>(
262 i * dim + d,
263 global_j + d,
264 A(i, j)));
265 }
266 }
267 }
268 }
269
270 bout.resize(b.rows() * dim, 1);
271 for (int i = 0; i < b.rows(); ++i)
272 {
273 for (int d = 0; d < dim; ++d)
274 {
275 bout(i * dim + d) = b(i, d);
276 }
277 }
278 }
279 else
280 {
281 assert(b.cols() == 1);
282 assert(b.size() % dim == 0);
283
284 for (int i = 0; i < A.rows(); ++i)
285 {
286 for (int j = 0; j < A.cols(); ++j)
287 {
288 if (A(i, j) != 0)
289 {
290 const auto nid = j / dim;
291 const auto noffset = j % dim;
292
293 const auto global_j = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
294
295 Ae.push_back(Eigen::Triplet<double>(
296 i,
297 global_j + noffset,
298 A(i, j)));
299 }
300 }
301 }
302
303 bout = b;
304 }
305
306 Aout.resize(bout.size(), n_dofs);
307 Aout.setFromTriplets(Ae.begin(), Ae.end());
308 Aout.makeCompressed();
309}
310
312 const int dim,
313 const Eigen::MatrixXd &A,
314 const Eigen::MatrixXd &b,
315 const std::vector<int> &local_to_global,
316 StiffnessMatrix &Aout,
317 Eigen::MatrixXd &bout)
318{
319 log_and_throw_error("Not implemented");
320}
321
322void polyfem::utils::scatter_matrix(const int n_dofs,
323 const int dim,
324 const std::vector<long> &shape,
325 const std::vector<int> &rows,
326 const std::vector<int> &cols,
327 const std::vector<double> &vals,
328 const Eigen::MatrixXd &b,
329 const std::vector<int> &local_to_global,
330 StiffnessMatrix &Aout,
331 Eigen::MatrixXd &bout)
332{
333 assert(shape.size() == 2);
334 assert(rows.size() == cols.size());
335 assert(rows.size() == vals.size());
336
337 std::vector<Eigen::Triplet<double>> Ae;
338 if (b.cols() == dim)
339 {
340 for (int k = 0; k < rows.size(); ++k)
341 {
342 const auto i = rows[k];
343 const auto j = cols[k];
344 const auto val = vals[k];
345
346 const auto global_j = (local_to_global.empty() ? j : local_to_global[j]) * dim;
347
348 if (val != 0)
349 {
350 for (int d = 0; d < dim; ++d)
351 {
352 Ae.push_back(Eigen::Triplet<double>(
353 i * dim + d,
354 global_j + d,
355 val));
356 }
357 }
358 }
359
360 bout.resize(b.rows() * dim, 1);
361 for (int i = 0; i < b.rows(); ++i)
362 {
363 for (int d = 0; d < dim; ++d)
364 {
365 bout(i * dim + d) = b(i, d);
366 }
367 }
368 }
369 else
370 {
371 assert(b.cols() == 1);
372 assert(b.size() % dim == 0);
373
374 for (int k = 0; k < rows.size(); ++k)
375 {
376 const auto i = rows[k];
377 const auto j = cols[k];
378 const auto val = vals[k];
379
380 if (val != 0)
381 {
382 const auto nid = j / dim;
383 const auto noffset = j % dim;
384
385 const auto global_j = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
386
387 Ae.push_back(Eigen::Triplet<double>(
388 i,
389 global_j + noffset,
390 val));
391 }
392 }
393
394 bout = b;
395 }
396
397 Aout.resize(bout.size(), n_dofs);
398 Aout.setFromTriplets(Ae.begin(), Ae.end());
399 Aout.makeCompressed();
400}
401
403 const int dim,
404 const std::vector<long> &shape,
405 const std::vector<int> &rows,
406 const std::vector<int> &cols,
407 const std::vector<double> &vals,
408 const Eigen::MatrixXd &b,
409 const std::vector<int> &local_to_global,
410 StiffnessMatrix &Aout,
411 Eigen::MatrixXd &bout)
412{
413 assert(shape.size() == 2);
414 assert(rows.size() == cols.size());
415 assert(rows.size() == vals.size());
416
417 std::vector<Eigen::Triplet<double>> Ae;
418
419 int A_cols = -1;
420
421 if (b.cols() == dim)
422 {
423 assert(n_dofs == b.rows() * dim);
424
425 for (int k = 0; k < rows.size(); ++k)
426 {
427 const auto i = rows[k];
428 const auto j = cols[k];
429 const auto val = vals[k];
430
431 const auto global_i = (local_to_global.empty() ? i : local_to_global[i]) * dim;
432
433 if (val != 0)
434 {
435 for (int d = 0; d < dim; ++d)
436 {
437 Ae.push_back(Eigen::Triplet<double>(
438 global_i + d,
439 j * dim + d,
440 val));
441 }
442 }
443 }
444
445 bout.resize(b.rows() * dim, 1);
446 for (int i = 0; i < b.rows(); ++i)
447 {
448 const auto global_i = (local_to_global.empty() ? i : local_to_global[i]) * dim;
449
450 for (int d = 0; d < dim; ++d)
451 {
452 bout(global_i + d) = b(i, d);
453 }
454 }
455
456 A_cols = shape[1] * dim;
457 assert(shape[0] * dim == n_dofs);
458 }
459 else
460 {
461 assert(b.cols() == 1);
462 assert(b.size() % dim == 0);
463 assert(n_dofs == b.size());
464
465 for (int k = 0; k < rows.size(); ++k)
466 {
467 const auto i = rows[k];
468 const auto j = cols[k];
469 const auto val = vals[k];
470
471 if (val != 0)
472 {
473 const auto nid = i / dim;
474 const auto noffset = i % dim;
475
476 const auto global_i = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
477
478 Ae.push_back(Eigen::Triplet<double>(
479 global_i + noffset,
480 j,
481 val));
482 }
483 }
484 bout.resize(b.size(), 1);
485 for (int i = 0; i < b.size(); ++i)
486 {
487 const auto nid = i / dim;
488 const auto noffset = i % dim;
489
490 const auto global_i = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
491
492 bout(global_i + noffset) = b(i);
493 }
494
495 assert(shape[0] == n_dofs);
496 A_cols = shape[1];
497 }
498
499 Aout.resize(n_dofs, A_cols);
500 Aout.setFromTriplets(Ae.begin(), Ae.end());
501 Aout.makeCompressed();
502}
Eigen::MatrixXd vec
Definition Assembler.cpp:72
double val
Definition Assembler.cpp:86
ElementAssemblyValues vals
Definition Assembler.cpp:22
std::vector< Eigen::Triplet< double > > entries
int y
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
Eigen::SparseMatrix< double > lump_matrix(const Eigen::SparseMatrix< double > &M)
Lump each row of a matrix into the diagonal.
void show_matrix_stats(const Eigen::MatrixXd &M)
Eigen::MatrixXd reorder_matrix(const Eigen::MatrixXd &in, const Eigen::VectorXi &in_to_out, int out_blocks=-1, const int block_size=1)
Reorder row blocks in a matrix.
void vector2matrix(const Eigen::VectorXd &vec, Eigen::MatrixXd &mat)
void scatter_matrix(const int n_dofs, const int dim, const Eigen::MatrixXd &A, const Eigen::MatrixXd &b, const std::vector< int > &local_to_global, StiffnessMatrix &Aout, Eigen::MatrixXd &bout)
Eigen::MatrixXd unreorder_matrix(const Eigen::MatrixXd &out, const Eigen::VectorXi &in_to_out, int in_blocks=-1, const int block_size=1)
Undo the reordering of row blocks in a matrix.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
void scatter_matrix_col(const int n_dofs, const int dim, const Eigen::MatrixXd &A, const Eigen::MatrixXd &b, const std::vector< int > &local_to_global, StiffnessMatrix &Aout, Eigen::MatrixXd &bout)
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
Eigen::MatrixXi map_index_matrix(const Eigen::MatrixXi &in, const Eigen::VectorXi &index_mapping)
Map the entrys of an index matrix to new indices.
void full_to_reduced_matrix(const int full_size, const int reduced_size, const std::vector< int > &removed_vars, const StiffnessMatrix &full, StiffnessMatrix &reduced)
Map a full size matrix to a reduced one by dropping rows and columns.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22