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;
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(
"----------------------------------------");
27 inline bool nanproof_equals(
const double x,
const double y)
29 return x ==
y || (std::isnan(
x) == std::isnan(
y));
37 return Eigen::VectorXd();
39 Eigen::VectorXd
x(X.size());
40 for (
int i = 0; i < X.rows(); ++i)
42 for (
int j = 0; j < X.cols(); ++j)
44 x(i * X.cols() + j) = X(i, j);
47 assert(nanproof_equals(X(0, 0),
x(0)));
48 assert(X.cols() <= 1 || nanproof_equals(X(0, 1),
x(1)));
56 return Eigen::MatrixXd(0, dim);
58 assert(
x.size() % dim == 0);
59 Eigen::MatrixXd X(
x.size() / dim, dim);
60 for (
int i = 0; i <
x.size(); ++i)
62 X(i / dim, i % dim) =
x(i);
64 assert(nanproof_equals(X(0, 0),
x(0)));
65 assert(X.cols() <= 1 || nanproof_equals(X(0, 1),
x(1)));
74 else if (
vec.size() == 4)
77 throw std::runtime_error(
"Invalid size in vector2matrix!");
79 assert(size * size ==
vec.size());
86 std::vector<Eigen::Triplet<double>> triplets;
88 for (
int k = 0; k < M.outerSize(); ++k)
90 for (Eigen::SparseMatrix<double>::InnerIterator it(M, k); it; ++it)
92 triplets.emplace_back(it.row(), it.row(), it.value());
96 Eigen::SparseMatrix<double> lumped(M.rows(), M.rows());
97 lumped.setFromTriplets(triplets.begin(), triplets.end());
98 lumped.makeCompressed();
105 const int reduced_size,
106 const std::vector<int> &removed_vars,
112 if (reduced_size == full_size || reduced_size == full.rows())
114 assert(reduced_size == full.rows() && reduced_size == full.cols());
118 assert(full.rows() == full_size && full.cols() == full_size);
120 Eigen::VectorXi indices(full_size);
123 for (
int i = 0; i < full_size; ++i)
125 if (kk < removed_vars.size() && removed_vars[kk] == i)
132 indices(i) = index++;
135 assert(index == reduced_size);
137 std::vector<Eigen::Triplet<double>>
entries;
138 entries.reserve(full.nonZeros());
139 for (
int k = 0; k < full.outerSize(); ++k)
144 for (StiffnessMatrix::InnerIterator it(full, k); it; ++it)
146 assert(it.col() == k);
147 if (indices(it.row()) < 0 || indices(it.col()) < 0)
150 assert(indices(it.row()) >= 0);
151 assert(indices(it.col()) >= 0);
153 entries.emplace_back(indices(it.row()), indices(it.col()), it.value());
157 reduced.resize(reduced_size, reduced_size);
159 reduced.makeCompressed();
163 const Eigen::MatrixXd &in,
164 const Eigen::VectorXi &in_to_out,
166 const int block_size)
168 constexpr double NaN = std::numeric_limits<double>::quiet_NaN();
170 assert(in.rows() % block_size == 0);
171 assert(in_to_out.size() == in.rows() / block_size);
174 out_blocks = in.rows() / block_size;
176 Eigen::MatrixXd out = Eigen::MatrixXd::Constant(
177 out_blocks * block_size, in.cols(), NaN);
179 const int in_blocks = in.rows() / block_size;
180 for (
int i = 0; i < in_blocks; ++i)
182 const int j = in_to_out[i];
186 out.middleRows(block_size * j, block_size) =
187 in.middleRows(block_size * i, block_size);
194 const Eigen::MatrixXd &out,
195 const Eigen::VectorXi &in_to_out,
197 const int block_size)
199 constexpr double NaN = std::numeric_limits<double>::quiet_NaN();
201 assert(out.rows() % block_size == 0);
204 in_blocks = out.rows() / block_size;
205 assert(in_to_out.size() == in_blocks);
207 Eigen::MatrixXd in = Eigen::MatrixXd::Constant(
208 in_blocks * block_size, out.cols(), NaN);
210 for (
int i = 0; i < in_blocks; i++)
212 const int j = in_to_out[i];
216 in.middleRows(block_size * i, block_size) =
217 out.middleRows(block_size * j, block_size);
224 const Eigen::MatrixXi &in,
225 const Eigen::VectorXi &index_mapping)
227 Eigen::MatrixXi out(in.rows(), in.cols());
228 for (
int i = 0; i < in.rows(); i++)
230 for (
int j = 0; j < in.cols(); j++)
232 out(i, j) = index_mapping[in(i, j)];
240 const Eigen::MatrixXd &A,
241 const Eigen::MatrixXd &b,
242 const std::vector<int> &local_to_global,
244 Eigen::MatrixXd &bout)
246 assert(A.rows() == b.rows());
247 std::vector<Eigen::Triplet<double>> Ae;
251 for (
int i = 0; i < A.rows(); ++i)
253 for (
int j = 0; j < A.cols(); ++j)
255 const auto global_j = (local_to_global.empty() ? j : local_to_global[j]) * dim;
259 for (
int d = 0; d < dim; ++d)
261 Ae.push_back(Eigen::Triplet<double>(
270 bout.resize(b.rows() * dim, 1);
271 for (
int i = 0; i < b.rows(); ++i)
273 for (
int d = 0; d < dim; ++d)
275 bout(i * dim + d) = b(i, d);
281 assert(b.cols() == 1);
282 assert(b.size() % dim == 0);
284 for (
int i = 0; i < A.rows(); ++i)
286 for (
int j = 0; j < A.cols(); ++j)
290 const auto nid = j / dim;
291 const auto noffset = j % dim;
293 const auto global_j = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
295 Ae.push_back(Eigen::Triplet<double>(
306 Aout.resize(bout.size(), n_dofs);
307 Aout.setFromTriplets(Ae.begin(), Ae.end());
308 Aout.makeCompressed();
313 const Eigen::MatrixXd &A,
314 const Eigen::MatrixXd &b,
315 const std::vector<int> &local_to_global,
317 Eigen::MatrixXd &bout)
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,
331 Eigen::MatrixXd &bout)
333 assert(shape.size() == 2);
334 assert(rows.size() == cols.size());
335 assert(rows.size() ==
vals.size());
337 std::vector<Eigen::Triplet<double>> Ae;
340 for (
int k = 0; k < rows.size(); ++k)
342 const auto i = rows[k];
343 const auto j = cols[k];
346 const auto global_j = (local_to_global.empty() ? j : local_to_global[j]) * dim;
350 for (
int d = 0; d < dim; ++d)
352 Ae.push_back(Eigen::Triplet<double>(
360 bout.resize(b.rows() * dim, 1);
361 for (
int i = 0; i < b.rows(); ++i)
363 for (
int d = 0; d < dim; ++d)
365 bout(i * dim + d) = b(i, d);
371 assert(b.cols() == 1);
372 assert(b.size() % dim == 0);
374 for (
int k = 0; k < rows.size(); ++k)
376 const auto i = rows[k];
377 const auto j = cols[k];
382 const auto nid = j / dim;
383 const auto noffset = j % dim;
385 const auto global_j = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
387 Ae.push_back(Eigen::Triplet<double>(
397 Aout.resize(bout.size(), n_dofs);
398 Aout.setFromTriplets(Ae.begin(), Ae.end());
399 Aout.makeCompressed();
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,
411 Eigen::MatrixXd &bout)
413 assert(shape.size() == 2);
414 assert(rows.size() == cols.size());
415 assert(rows.size() ==
vals.size());
417 std::vector<Eigen::Triplet<double>> Ae;
423 assert(n_dofs == b.rows() * dim);
425 for (
int k = 0; k < rows.size(); ++k)
427 const auto i = rows[k];
428 const auto j = cols[k];
431 const auto global_i = (local_to_global.empty() ? i : local_to_global[i]) * dim;
435 for (
int d = 0; d < dim; ++d)
437 Ae.push_back(Eigen::Triplet<double>(
445 bout.resize(b.rows() * dim, 1);
446 for (
int i = 0; i < b.rows(); ++i)
448 const auto global_i = (local_to_global.empty() ? i : local_to_global[i]) * dim;
450 for (
int d = 0; d < dim; ++d)
452 bout(global_i + d) = b(i, d);
456 A_cols = shape[1] * dim;
457 assert(shape[0] * dim == n_dofs);
461 assert(b.cols() == 1);
462 assert(b.size() % dim == 0);
463 assert(n_dofs == b.size());
465 for (
int k = 0; k < rows.size(); ++k)
467 const auto i = rows[k];
468 const auto j = cols[k];
473 const auto nid = i / dim;
474 const auto noffset = i % dim;
476 const auto global_i = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
478 Ae.push_back(Eigen::Triplet<double>(
484 bout.resize(b.size(), 1);
485 for (
int i = 0; i < b.size(); ++i)
487 const auto nid = i / dim;
488 const auto noffset = i % dim;
490 const auto global_i = (local_to_global.empty() ? nid : local_to_global[nid]) * dim;
492 bout(global_i + noffset) = b(i);
495 assert(shape[0] == n_dofs);
499 Aout.resize(n_dofs, A_cols);
500 Aout.setFromTriplets(Ae.begin(), Ae.end());
501 Aout.makeCompressed();
ElementAssemblyValues vals
std::vector< Eigen::Triplet< double > > entries
#define POLYFEM_SCOPED_TIMER(...)
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.
void log_and_throw_error(const std::string &msg)
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix