PolyFEM
Loading...
Searching...
No Matches
svd.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <Eigen/Dense>
6
7#define USE_IQRSVD
8
9namespace polyfem
10{
11 namespace utils
12 {
13 template <typename MatrixType>
14 class AutoFlipSVD : Eigen::JacobiSVD<MatrixType> {
15 protected:
17
18 typename Eigen::JacobiSVD<MatrixType>::SingularValuesType singularValues_flipped;
20
21 public:
22 AutoFlipSVD(void) {}
23 AutoFlipSVD(const MatrixType& mtr, unsigned int computationOptions = 0)
24 {
25 compute(mtr, computationOptions);
26 }
27
28 public:
29 template <int dim = MatrixType::RowsAtCompileTime>
30 typename std::enable_if<dim == 3, AutoFlipSVD<MatrixType>>::type&
31 compute(const MatrixType& mtr, unsigned int computationOptions)
32 {
34 #ifdef USE_IQRSVD
39 #else
40 if ((computationOptions & Eigen::ComputeFullU) || (computationOptions & Eigen::ComputeFullV)) {
42 }
43 else {
45 }
46 #endif
47 return *this;
48 }
49 template <int dim = MatrixType::RowsAtCompileTime>
50 typename std::enable_if<dim == 2, AutoFlipSVD<MatrixType>>::type&
51 compute(const MatrixType& mtr, unsigned int computationOptions)
52 {
53 #ifdef USE_IQRSVD
59 #else
61 Eigen::JacobiSVD<MatrixType>::compute(mtr, computationOptions);
62 flip2d(mtr, computationOptions);
63 #endif
64 return *this;
65 }
66
67 void set(const Eigen::Matrix3d& U, const Eigen::Vector3d& Sigma, const Eigen::Matrix3d& V)
68 {
69 flipped_U = true, flipped_V = true, flipped_sigma = true;
73 }
74
75 protected:
76 void flip2d(const MatrixType& mtr, unsigned int computationOptions)
77 {
79 bool fullUComputed = (computationOptions & Eigen::ComputeFullU);
80 bool fullVComputed = (computationOptions & Eigen::ComputeFullV);
81 if (fullUComputed && fullVComputed) {
82 if (Eigen::JacobiSVD<MatrixType>::m_matrixU.determinant() < 0.0) {
83 matrixU_flipped = Eigen::JacobiSVD<MatrixType>::m_matrixU;
84 matrixU_flipped.col(1) *= -1.0;
85 flipped_U = true;
86
87 if (!flipped_sigma) {
88 singularValues_flipped = Eigen::JacobiSVD<MatrixType>::m_singularValues;
89 }
90 singularValues_flipped[1] *= -1.0;
91 flipped_sigma = true;
92 }
93 if (Eigen::JacobiSVD<MatrixType>::m_matrixV.determinant() < 0.0) {
94 matrixV_flipped = Eigen::JacobiSVD<MatrixType>::m_matrixV;
95 matrixV_flipped.col(1) *= -1.0;
96 flipped_V = true;
97
98 if (!flipped_sigma) {
99 singularValues_flipped = Eigen::JacobiSVD<MatrixType>::m_singularValues;
100 }
101 singularValues_flipped[1] *= -1.0;
102 flipped_sigma = true;
103 }
104 }
105 else if (mtr.determinant() < 0.0) {
106 singularValues_flipped = Eigen::JacobiSVD<MatrixType>::m_singularValues;
107 singularValues_flipped[1] *= -1.0;
108 flipped_sigma = true;
109 }
110
111 if (std::isnan(singularValues()[0]) || std::isnan(singularValues()[1])) {
112 // degenerated case
113 singularValues_flipped.setZero();
114 flipped_sigma = true;
115 if (fullUComputed && fullVComputed) {
116 matrixU_flipped.setIdentity();
117 matrixV_flipped.setIdentity();
118 flipped_U = flipped_V = true;
119 }
120 }
121 }
122
123 //TODO: merge with IglUtils::computeCofactorMtr
124 template <int dim>
125 void computeCofactorMtr(const Eigen::Matrix<double, dim, dim>& F,
126 Eigen::Matrix<double, dim, dim>& A)
127 {
128 switch (dim) {
129 case 2:
130 A(0, 0) = F(1, 1);
131 A(0, 1) = -F(1, 0);
132 A(1, 0) = -F(0, 1);
133 A(1, 1) = F(0, 0);
134 break;
135
136 case 3:
137 A(0, 0) = F(1, 1) * F(2, 2) - F(1, 2) * F(2, 1);
138 A(0, 1) = F(1, 2) * F(2, 0) - F(1, 0) * F(2, 2);
139 A(0, 2) = F(1, 0) * F(2, 1) - F(1, 1) * F(2, 0);
140 A(1, 0) = F(0, 2) * F(2, 1) - F(0, 1) * F(2, 2);
141 A(1, 1) = F(0, 0) * F(2, 2) - F(0, 2) * F(2, 0);
142 A(1, 2) = F(0, 1) * F(2, 0) - F(0, 0) * F(2, 1);
143 A(2, 0) = F(0, 1) * F(1, 2) - F(0, 2) * F(1, 1);
144 A(2, 1) = F(0, 2) * F(1, 0) - F(0, 0) * F(1, 2);
145 A(2, 2) = F(0, 0) * F(1, 1) - F(0, 1) * F(1, 0);
146 break;
147
148 default:
149 assert(0 && "dim not 2 or 3");
150 break;
151 }
152 }
153 void fastEigenvalues(const Eigen::Matrix3d& A_Sym,
154 Eigen::Vector3d& lambda)
155 // 24 mults, 20 adds, 1 atan2, 1 sincos, 2 sqrts
156 {
157 using T = double;
158 using std::max;
159 using std::swap;
160 T m = ((T)1 / 3) * (A_Sym(0, 0) + A_Sym(1, 1) + A_Sym(2, 2));
161 T a00 = A_Sym(0, 0) - m;
162 T a11 = A_Sym(1, 1) - m;
163 T a22 = A_Sym(2, 2) - m;
164 T a12_sqr = A_Sym(0, 1) * A_Sym(0, 1);
165 T a13_sqr = A_Sym(0, 2) * A_Sym(0, 2);
166 T a23_sqr = A_Sym(1, 2) * A_Sym(1, 2);
167 T p = ((T)1 / 6) * (a00 * a00 + a11 * a11 + a22 * a22 + 2 * (a12_sqr + a13_sqr + a23_sqr));
168 T q = (T).5 * (a00 * (a11 * a22 - a23_sqr) - a11 * a13_sqr - a22 * a12_sqr) + A_Sym(0, 1) * A_Sym(0, 2) * A_Sym(1, 2);
169 T sqrt_p = sqrt(p);
170 T disc = p * p * p - q * q;
171 T phi = ((T)1 / 3) * atan2(sqrt(max((T)0, disc)), q);
172 T c = cos(phi), s = sin(phi);
173 T sqrt_p_cos = sqrt_p * c;
174 T root_three_sqrt_p_sin = sqrt((T)3) * sqrt_p * s;
175 lambda(0) = m + 2 * sqrt_p_cos;
176 lambda(1) = m - sqrt_p_cos - root_three_sqrt_p_sin;
177 lambda(2) = m - sqrt_p_cos + root_three_sqrt_p_sin;
178 if (lambda(0) < lambda(1))
179 swap(lambda(0), lambda(1));
180 if (lambda(1) < lambda(2))
181 swap(lambda(1), lambda(2));
182 if (lambda(0) < lambda(1))
183 swap(lambda(0), lambda(1));
184 }
185 void fastEigenvectors(const Eigen::Matrix3d& A_Sym,
186 const Eigen::Vector3d& lambda,
187 Eigen::Matrix3d& V)
188 // 71 mults, 44 adds, 3 divs, 3 sqrts
189 {
190 // flip if necessary so that first eigenvalue is the most different
191 using T = double;
192 using std::sqrt;
193 using std::swap;
194 bool flipped = false;
195 Eigen::Vector3d lambda_flip(lambda);
196 if (lambda(0) - lambda(1) < lambda(1) - lambda(2)) { // 2a
197 swap(lambda_flip(0), lambda_flip(2));
198 flipped = true;
199 }
200
201 // get first eigenvector
202 Eigen::Matrix3d C1;
203 computeCofactorMtr<3>(A_Sym - lambda_flip(0) * Eigen::Matrix3d::Identity(), C1);
204 Eigen::Matrix3d::Index i;
205 T norm2 = C1.colwise().squaredNorm().maxCoeff(&i); // 3a + 12m+6a + 9m+6a+1d+1s = 21m+15a+1d+1s
206 Eigen::Vector3d v1;
207 if (norm2 != 0) {
208 T one_over_sqrt = (T)1 / sqrt(norm2);
209 v1 = C1.col(i) * one_over_sqrt;
210 }
211 else
212 v1 << 1, 0, 0;
213
214 // form basis for orthogonal complement to v1, and reduce A to this space
215 Eigen::Vector3d v1_orthogonal = v1.unitOrthogonal(); // 6m+2a+1d+1s (tweak: 5m+1a+1d+1s)
216 Eigen::Matrix<T, 3, 2> other_v;
217 other_v.col(0) = v1_orthogonal;
218 other_v.col(1) = v1.cross(v1_orthogonal); // 6m+3a (tweak: 4m+1a)
219 Eigen::Matrix2d A_reduced = other_v.transpose() * A_Sym * other_v; // 21m+12a (tweak: 18m+9a)
220
221 // find third eigenvector from A_reduced, and fill in second via cross product
222 Eigen::Matrix2d C3;
223 computeCofactorMtr<2>(A_reduced - lambda_flip(2) * Eigen::Matrix2d::Identity(), C3);
224 Eigen::Matrix2d::Index j;
225 norm2 = C3.colwise().squaredNorm().maxCoeff(&j); // 3a + 12m+6a + 9m+6a+1d+1s = 21m+15a+1d+1s
226 Eigen::Vector3d v3;
227 if (norm2 != 0) {
228 T one_over_sqrt = (T)1 / sqrt(norm2);
229 v3 = other_v * C3.col(j) * one_over_sqrt;
230 }
231 else
232 v3 = other_v.col(0);
233
234 Eigen::Vector3d v2 = v3.cross(v1); // 6m+3a
235
236 // finish
237 if (flipped) {
238 V.col(0) = v3;
239 V.col(1) = v2;
240 V.col(2) = -v1;
241 }
242 else {
243 V.col(0) = v1;
244 V.col(1) = v2;
245 V.col(2) = v3;
246 }
247 }
248 void fastSolveEigenproblem(const Eigen::Matrix3d& A_Sym,
249 Eigen::Vector3d& lambda,
250 Eigen::Matrix3d& V)
251 // 71 mults, 44 adds, 3 divs, 3 sqrts
252 {
253 fastEigenvalues(A_Sym, lambda);
254 fastEigenvectors(A_Sym, lambda, V);
255 }
256
257 void fastSVD3d(const Eigen::Matrix3d& A,
258 Eigen::Matrix3d& U,
259 Eigen::Vector3d& singular_values,
260 Eigen::Matrix3d& V)
261 // 182 mults, 112 adds, 6 divs, 11 sqrts, 1 atan2, 1 sincos
262 {
263 using T = double;
264 // decompose normal equations
265 Eigen::Vector3d lambda;
266 fastSolveEigenproblem(A.transpose() * A, lambda, V);
267
268 // compute singular values
269 if (lambda(2) < 0)
270 lambda = (lambda.array() >= (T)0).select(lambda, (T)0);
271 singular_values = lambda.array().sqrt();
272 if (A.determinant() < 0)
274
275 // compute singular vectors
276 U.col(0) = A * V.col(0);
277 T norm = U.col(0).norm();
278 if (norm != 0) {
279 T one_over_norm = (T)1 / norm;
280 U.col(0) = U.col(0) * one_over_norm;
281 }
282 else
283 U.col(0) << 1, 0, 0;
284 Eigen::Vector3d v1_orthogonal = U.col(0).unitOrthogonal();
285 Eigen::Matrix<T, 3, 2> other_v;
286 other_v.col(0) = v1_orthogonal;
287 other_v.col(1) = U.col(0).cross(v1_orthogonal);
288 Eigen::Vector2d w = other_v.transpose() * A * V.col(1);
289 norm = w.norm();
290 if (norm != 0) {
291 T one_over_norm = (T)1 / norm;
292 w = w * one_over_norm;
293 }
294 else
295 w << 1, 0;
296 U.col(1) = other_v * w;
297 U.col(2) = U.col(0).cross(U.col(1));
298 }
299
300 void fastComputeSingularValues3d(const Eigen::Matrix3d& A,
301 Eigen::Vector3d& singular_values)
302 {
303 using T = double;
304 // decompose normal equations
305 Eigen::Vector3d lambda;
306 fastEigenvalues(A.transpose() * A, lambda);
307
308 // compute singular values
309 if (lambda(2) < 0)
310 lambda = (lambda.array() >= (T)0).select(lambda, (T)0);
311 singular_values = lambda.array().sqrt();
312 if (A.determinant() < 0)
314 }
315
316 public:
317 const typename Eigen::JacobiSVD<MatrixType>::SingularValuesType& singularValues(void) const
318 {
319 if (flipped_sigma) {
321 }
322 else {
323 return Eigen::JacobiSVD<MatrixType>::singularValues();
324 }
325 }
326 const MatrixType& matrixU(void) const
327 {
328 if (flipped_U) {
329 return matrixU_flipped;
330 }
331 else {
332 return Eigen::JacobiSVD<MatrixType>::matrixU();
333 }
334 }
335 const MatrixType& matrixV(void) const
336 {
337 if (flipped_V) {
338 return matrixV_flipped;
339 }
340 else {
341 return Eigen::JacobiSVD<MatrixType>::matrixV();
342 }
343 }
344
345 void setIdentity(void)
346 {
347 flipped_sigma = true;
348 flipped_V = true;
349 flipped_U = true;
350
351 matrixU_flipped.setIdentity();
352 matrixV_flipped.setIdentity();
353 singularValues_flipped.setOnes();
354 }
355 };
356
357 template <int dim>
358 Eigen::Vector<double, dim> singular_values(const Eigen::Matrix<double, dim, dim> &A)
359 {
360 AutoFlipSVD<Eigen::Matrix<double, dim, dim>> svd(A, Eigen::ComputeFullU | Eigen::ComputeFullV);
361 return svd.singularValues();
362 }
363 }
364}
int V
const Eigen::JacobiSVD< MatrixType >::SingularValuesType & singularValues(void) const
Definition svd.hpp:317
std::enable_if< dim==2, AutoFlipSVD< MatrixType > >::type & compute(const MatrixType &mtr, unsigned int computationOptions)
Definition svd.hpp:51
void fastComputeSingularValues3d(const Eigen::Matrix3d &A, Eigen::Vector3d &singular_values)
Definition svd.hpp:300
void fastSolveEigenproblem(const Eigen::Matrix3d &A_Sym, Eigen::Vector3d &lambda, Eigen::Matrix3d &V)
Definition svd.hpp:248
MatrixType matrixU_flipped
Definition svd.hpp:19
void fastEigenvectors(const Eigen::Matrix3d &A_Sym, const Eigen::Vector3d &lambda, Eigen::Matrix3d &V)
Definition svd.hpp:185
const MatrixType & matrixV(void) const
Definition svd.hpp:335
AutoFlipSVD(const MatrixType &mtr, unsigned int computationOptions=0)
Definition svd.hpp:23
void flip2d(const MatrixType &mtr, unsigned int computationOptions)
Definition svd.hpp:76
void fastSVD3d(const Eigen::Matrix3d &A, Eigen::Matrix3d &U, Eigen::Vector3d &singular_values, Eigen::Matrix3d &V)
Definition svd.hpp:257
void computeCofactorMtr(const Eigen::Matrix< double, dim, dim > &F, Eigen::Matrix< double, dim, dim > &A)
Definition svd.hpp:125
const MatrixType & matrixU(void) const
Definition svd.hpp:326
Eigen::JacobiSVD< MatrixType >::SingularValuesType singularValues_flipped
Definition svd.hpp:18
void fastEigenvalues(const Eigen::Matrix3d &A_Sym, Eigen::Vector3d &lambda)
Definition svd.hpp:153
MatrixType matrixV_flipped
Definition svd.hpp:19
std::enable_if< dim==3, AutoFlipSVD< MatrixType > >::type & compute(const MatrixType &mtr, unsigned int computationOptions)
Definition svd.hpp:31
void set(const Eigen::Matrix3d &U, const Eigen::Vector3d &Sigma, const Eigen::Matrix3d &V)
Definition svd.hpp:67
enable_if_t< isSize< TA >(2, 2) &&isSize< Ts >(2, 1)> singularValueDecomposition(const Eigen::MatrixBase< TA > &A, GivensRotation< T > &U, const Eigen::MatrixBase< Ts > &Sigma, GivensRotation< T > &V, const ScalarType< TA > tol=64 *std::numeric_limits< ScalarType< TA > >::epsilon())
2x2 SVD (singular value decomposition) A=USV'
Eigen::Vector< double, dim > singular_values(const Eigen::Matrix< double, dim, dim > &A)
Definition svd.hpp:358
T determinant(const Eigen::Matrix< T, rows, cols, option, maxRow, maxCol > &mat)