PolyFEM
Loading...
Searching...
No Matches
MatrixUtils.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <Eigen/Dense>
6#include <Eigen/Sparse>
7
8namespace polyfem
9{
10 namespace utils
11 {
12 // Show some stats about the matrix M: det, singular values, condition number, etc
13 void show_matrix_stats(const Eigen::MatrixXd &M);
14
15 template <typename T>
16 T determinant(const Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> &mat)
17 {
18 assert(mat.rows() == mat.cols());
19
20 if (mat.rows() == 1)
21 return mat(0);
22 else if (mat.rows() == 2)
23 return mat(0, 0) * mat(1, 1) - mat(0, 1) * mat(1, 0);
24 else if (mat.rows() == 3)
25 return mat(0, 0) * (mat(1, 1) * mat(2, 2) - mat(1, 2) * mat(2, 1)) - mat(0, 1) * (mat(1, 0) * mat(2, 2) - mat(1, 2) * mat(2, 0)) + mat(0, 2) * (mat(1, 0) * mat(2, 1) - mat(1, 1) * mat(2, 0));
26
27 assert(false);
28 return T(0);
29 }
30
31 template <typename T>
32 Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> inverse(const Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> &mat)
33 {
34 assert(mat.rows() == mat.cols());
35
36 if (mat.rows() == 1)
37 {
38 Eigen::Matrix<T, 1, 1> inv;
39 inv(0, 0) = T(1.) / mat(0);
40
41 return inv;
42 }
43 else if (mat.rows() == 2)
44 {
45 Eigen::Matrix<T, 2, 2> inv;
46 T det = determinant(mat);
47 inv(0, 0) = mat(1, 1) / det;
48 inv(0, 1) = -mat(0, 1) / det;
49 inv(1, 0) = -mat(1, 0) / det;
50 inv(1, 1) = mat(0, 0) / det;
51
52 return inv;
53 }
54 else if (mat.rows() == 3)
55 {
56 Eigen::Matrix<T, 3, 3> inv;
57 T det = determinant(mat);
58 inv(0, 0) = (-mat(1, 2) * mat(2, 1) + mat(1, 1) * mat(2, 2)) / det;
59 inv(0, 1) = (mat(0, 2) * mat(2, 1) - mat(0, 1) * mat(2, 2)) / det;
60 inv(0, 2) = (-mat(0, 2) * mat(1, 1) + mat(0, 1) * mat(1, 2)) / det;
61 inv(1, 0) = (mat(1, 2) * mat(2, 0) - mat(1, 0) * mat(2, 2)) / det;
62 inv(1, 1) = (-mat(0, 2) * mat(2, 0) + mat(0, 0) * mat(2, 2)) / det;
63 inv(1, 2) = (mat(0, 2) * mat(1, 0) - mat(0, 0) * mat(1, 2)) / det;
64 inv(2, 0) = (-mat(1, 1) * mat(2, 0) + mat(1, 0) * mat(2, 1)) / det;
65 inv(2, 1) = (mat(0, 1) * mat(2, 0) - mat(0, 0) * mat(2, 1)) / det;
66 inv(2, 2) = (-mat(0, 1) * mat(1, 0) + mat(0, 0) * mat(1, 1)) / det;
67
68 return inv;
69 }
70
71 assert(false);
72 Eigen::Matrix<T, 1, 1> inv;
73 inv(0, 0) = T(0);
74 return inv;
75 }
76
77 inline Eigen::SparseMatrix<double> sparse_identity(int rows, int cols)
78 {
79 Eigen::SparseMatrix<double> I(rows, cols);
80 I.setIdentity();
81 return I;
82 }
83
85 Eigen::VectorXd flatten(const Eigen::MatrixXd &X);
86
88 Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim);
89
90 void vector2matrix(const Eigen::VectorXd &vec, Eigen::MatrixXd &mat);
91
95 Eigen::SparseMatrix<double> lump_matrix(const Eigen::SparseMatrix<double> &M);
96
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
116 Eigen::MatrixXd reorder_matrix(
117 const Eigen::MatrixXd &in,
118 const Eigen::VectorXi &in_to_out,
119 int out_blocks = -1,
120 const int block_size = 1);
121
128 Eigen::MatrixXd unreorder_matrix(
129 const Eigen::MatrixXd &out,
130 const Eigen::VectorXi &in_to_out,
131 int in_blocks = -1,
132 const int block_size = 1);
133
138 Eigen::MatrixXi map_index_matrix(
139 const Eigen::MatrixXi &in,
140 const Eigen::VectorXi &index_mapping);
141
142 template <typename DstMat, typename SrcMat>
143 void append_rows(DstMat &dst, const SrcMat &src)
144 {
145 if (src.rows() == 0)
146 return;
147 if (dst.cols() == 0)
148 dst.resize(dst.rows(), src.cols());
149 assert(dst.cols() == src.cols());
150 dst.conservativeResize(dst.rows() + src.rows(), dst.cols());
151 dst.bottomRows(src.rows()) = src;
152 }
153
154 template <typename DstMat>
155 void append_rows_of_zeros(DstMat &dst, const size_t n_zero_rows)
156 {
157 assert(dst.cols() > 0);
158 if (n_zero_rows == 0)
159 return;
160 dst.conservativeResize(dst.rows() + n_zero_rows, dst.cols());
161 dst.bottomRows(n_zero_rows).setZero();
162 }
163 } // namespace utils
164} // namespace polyfem
Eigen::MatrixXd vec
Definition Assembler.cpp:72
int x
Eigen::MatrixXd mat
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.
Eigen::SparseMatrix< double > sparse_identity(int rows, int cols)
void vector2matrix(const Eigen::VectorXd &vec, Eigen::MatrixXd &mat)
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > inverse(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &mat)
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 append_rows_of_zeros(DstMat &dst, const size_t n_zero_rows)
T determinant(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &mat)
void append_rows(DstMat &dst, const SrcMat &src)
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.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20