PolyFEM
Loading...
Searching...
No Matches
JIXIE Namespace Reference

Copyright (c) 2016 Theodore Gast, Chuyuan Fu, Chenfanfu Jiang, Joseph Teran. More...

Namespaces

namespace  INTERNAL
 
namespace  MATH_TOOLS
 

Classes

class  GivensRotation
 Class for givens rotation. More...
 

Typedefs

template<bool B, class T = void>
using enable_if_t = typename std::enable_if< B, T >::type
 
template<class T >
using ScalarType = typename INTERNAL::ScalarTypeHelper< T >::type
 

Functions

template<class T >
void zeroChase (Eigen::Matrix< T, 3, 3 > &H, Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 3 > &V)
 zero chasing the 3X3 matrix to bidiagonal form original form of H: x x 0 x x x 0 0 x after zero chase: x x 0 0 x x 0 0 x
 
template<class T >
void makeUpperBidiag (Eigen::Matrix< T, 3, 3 > &H, Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 3 > &V)
 make a 3X3 matrix to upper bidiagonal form original form of H: x x x x x x x x x after zero chase: x x 0 0 x x 0 0 x
 
template<class T >
void makeLambdaShape (Eigen::Matrix< T, 3, 3 > &H, Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 3 > &V)
 make a 3X3 matrix to lambda shape original form of H: x x x x x x x x x after : x 0 0 x x 0 x 0 x
 
template<class TA , class T , class TS >
enable_if_t< isSize< TA >(2, 2) &&isSize< TS >(2, 2)> polarDecomposition (const Eigen::MatrixBase< TA > &A, GivensRotation< T > &R, const Eigen::MatrixBase< TS > &S_Sym)
 2x2 polar decomposition.
 
template<class TA , class TR , class TS >
enable_if_t< isSize< TA >(2, 2) &&isSize< TR >(2, 2) &&isSize< TS >(2, 2)> polarDecomposition (const Eigen::MatrixBase< TA > &A, const Eigen::MatrixBase< TR > &R, const Eigen::MatrixBase< TS > &S_Sym)
 2x2 polar decomposition.
 
template<class TA , class T , class Ts >
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'
 
template<class TA , class TU , class Ts , class TV >
enable_if_t< isSize< TA >(2, 2) &&isSize< TU >(2, 2) &&isSize< TV >(2, 2) &&isSize< Ts >(2, 1)> singularValueDecomposition (const Eigen::MatrixBase< TA > &A, const Eigen::MatrixBase< TU > &U, const Eigen::MatrixBase< Ts > &Sigma, const Eigen::MatrixBase< TV > &V, const ScalarType< TA > tol=64 *std::numeric_limits< ScalarType< TA > >::epsilon())
 2x2 SVD (singular value decomposition) A=USV'
 
template<class T >
wilkinsonShift (const T a1, const T b1, const T a2)
 compute wilkinsonShift of the block a1 b1 b1 a2 based on the wilkinsonShift formula mu = c + d - sign (d) \ sqrt (d*d + b*b), where d = (a-c)/2
 
template<int t, class T >
void process (Eigen::Matrix< T, 3, 3 > &B, Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 1 > &sigma, Eigen::Matrix< T, 3, 3 > &V)
 Helper function of 3X3 SVD for processing 2X2 SVD.
 
template<class T >
void flipSign (int i, Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 1 > &sigma)
 Helper function of 3X3 SVD for flipping signs due to flipping signs of sigma.
 
template<int t, class T >
enable_if_t< t==0 > sort (Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 1 > &sigma, Eigen::Matrix< T, 3, 3 > &V)
 Helper function of 3X3 SVD for sorting singular values.
 
template<int t, class T >
enable_if_t< t==1 > sort (Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 1 > &sigma, Eigen::Matrix< T, 3, 3 > &V)
 Helper function of 3X3 SVD for sorting singular values.
 
template<class T >
int singularValueDecomposition (const Eigen::Matrix< T, 3, 3 > &A, Eigen::Matrix< T, 3, 3 > &U, Eigen::Matrix< T, 3, 1 > &sigma, Eigen::Matrix< T, 3, 3 > &V, T tol=1024 *std::numeric_limits< T >::epsilon())
 3X3 SVD (singular value decomposition) A=USV'
 
template<class T >
void polarDecomposition (const Eigen::Matrix< T, 3, 3 > &A, Eigen::Matrix< T, 3, 3 > &R, Eigen::Matrix< T, 3, 3 > &S_Sym)
 3X3 polar decomposition.
 
template<class MatrixType >
constexpr bool isSize (int m, int n)
 

Detailed Description

Copyright (c) 2016 Theodore Gast, Chuyuan Fu, Chenfanfu Jiang, Joseph Teran.

Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

If the code is used in an article, the following paper shall be cited: @techreport{qrsvd:2016, title={Implicit-shifted Symmetric QR Singular Value Decomposition of 3x3 Matrices}, author={Gast, Theodore and Fu, Chuyuan and Jiang, Chenfanfu and Teran, Joseph}, year={2016}, institution={University of California Los Angeles} }

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

This file implements 2D and 3D polar decompositions and SVDs.

T may be float or double.

2D Polar: Eigen::Matrix<T, 2, 2> A,R,S; A<<1,2,3,4; JIXIE::polarDecomposition(A, R, S); R will be the closest rotation to A S will be symmetric

2D SVD: Eigen::Matrix<T, 2, 2> A; A<<1,2,3,4; Eigen::Matrix<T, 2, 1> S; Eigen::Matrix<T, 2, 2> U; Eigen::Matrix<T, 2, 2> V; JIXIE::singularValueDecomposition(A,U,S,V); A = U S V' U and V will be rotations S will be singular values sorted by decreasing magnitude. Only the last one may be negative.

3D Polar: Eigen::Matrix<T, 3, 3> A,R,S; A<<1,2,3,4,5,6; JIXIE::polarDecomposition(A, R, S); R will be the closest rotation to A S will be symmetric

3D SVD: Eigen::Matrix<T, 3, 3> A; A<<1,2,3,4,5,6; Eigen::Matrix<T, 3, 1> S; Eigen::Matrix<T, 3, 3> U; Eigen::Matrix<T, 3, 3> V; JIXIE::singularValueDecomposition(A,U,S,V); A = U S V' U and V will be rotations S will be singular values sorted by decreasing magnitude. Only the last one may be negative.

SVD based on implicit QR with Wilkinson Shift

Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

If the code is used in an article, the following paper shall be cited: @techreport{qrsvd:2016, title={Implicit-shifted Symmetric QR Singular Value Decomposition of 3x3 Matrices}, author={Gast, Theodore and Fu, Chuyuan and Jiang, Chenfanfu and Teran, Joseph}, year={2016}, institution={University of California Los Angeles} }

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

This file provides a random number generator and a timer. Sample usage: RandomNumber<float> rand; float x = randReal(-0.5, 0.8);

Timer timer; timer.start(); SOME CODE A std::cout<<"CODE A took "<<timer.click()<<" seconds"<<std::endl; SOME CODE B std::cout<<"CODE B took "<<timer.click()<<" seconds"<<std::endl;

Typedef Documentation

◆ enable_if_t

template<bool B, class T = void>
using JIXIE::enable_if_t = typedef typename std::enable_if<B, T>::type

Definition at line 68 of file Tools.h.

◆ ScalarType

template<class T >
using JIXIE::ScalarType = typedef typename INTERNAL::ScalarTypeHelper<T>::type

Definition at line 126 of file Tools.h.

Function Documentation

◆ flipSign()

template<class T >
void JIXIE::flipSign ( int  i,
Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 1 > &  sigma 
)
inline

Helper function of 3X3 SVD for flipping signs due to flipping signs of sigma.

Definition at line 590 of file ImplicitQRSVD.h.

Referenced by sort().

Here is the caller graph for this function:

◆ isSize()

template<class MatrixType >
constexpr bool JIXIE::isSize ( int  m,
int  n 
)
constexpr

Definition at line 129 of file Tools.h.

◆ makeLambdaShape()

template<class T >
void JIXIE::makeLambdaShape ( Eigen::Matrix< T, 3, 3 > &  H,
Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 3 > &  V 
)
inline

make a 3X3 matrix to lambda shape original form of H: x x x x x x x x x after : x 0 0 x x 0 x 0 x

Reduce H to of form x x 0 x x x x x x

Reduce H to of form x x 0 x x 0 x x x

Reduce H to of form x x 0 x x 0 x 0 x

Reduce H to of form x 0 0 x x 0 x 0 x

Definition at line 341 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::columnRotation(), JIXIE::GivensRotation< T >::computeUnconventional(), JIXIE::GivensRotation< T >::rowRotation(), and V.

Here is the call graph for this function:

◆ makeUpperBidiag()

template<class T >
void JIXIE::makeUpperBidiag ( Eigen::Matrix< T, 3, 3 > &  H,
Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 3 > &  V 
)
inline

make a 3X3 matrix to upper bidiagonal form original form of H: x x x x x x x x x after zero chase: x x 0 0 x x 0 0 x

Reduce H to of form x x x x x x 0 x x

Definition at line 310 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::columnRotation(), JIXIE::GivensRotation< T >::rowRotation(), V, and zeroChase().

Referenced by singularValueDecomposition().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ polarDecomposition() [1/3]

template<class T >
void JIXIE::polarDecomposition ( const Eigen::Matrix< T, 3, 3 > &  A,
Eigen::Matrix< T, 3, 3 > &  R,
Eigen::Matrix< T, 3, 3 > &  S_Sym 
)
inline

3X3 polar decomposition.

Parameters
[in]Amatrix.
[out]RRobustly a rotation matrix.
[out]S_SymSymmetric. Whole matrix is stored

Whole matrix S is stored Polar guarantees negative sign is on the small magnitude singular value. S is guaranteed to be the closest one to identity. R is guaranteed to be the closest rotation to A.

Definition at line 864 of file ImplicitQRSVD.h.

References singularValueDecomposition(), and V.

Here is the call graph for this function:

◆ polarDecomposition() [2/3]

template<class TA , class TR , class TS >
enable_if_t< isSize< TA >(2, 2) &&isSize< TR >(2, 2) &&isSize< TS >(2, 2)> JIXIE::polarDecomposition ( const Eigen::MatrixBase< TA > &  A,
const Eigen::MatrixBase< TR > &  R,
const Eigen::MatrixBase< TS > &  S_Sym 
)
inline

2x2 polar decomposition.

Parameters
[in]Amatrix.
[out]RRobustly a rotation matrix.
[out]S_SymSymmetric. Whole matrix is stored

Whole matrix S is stored since its faster to calculate due to simd vectorization Polar guarantees negative sign is on the small magnitude singular value. S is guaranteed to be the closest one to identity. R is guaranteed to be the closest rotation to A.

Definition at line 437 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::fill(), and polarDecomposition().

Here is the call graph for this function:

◆ polarDecomposition() [3/3]

template<class TA , class T , class TS >
enable_if_t< isSize< TA >(2, 2) &&isSize< TS >(2, 2)> JIXIE::polarDecomposition ( const Eigen::MatrixBase< TA > &  A,
GivensRotation< T > &  R,
const Eigen::MatrixBase< TS > &  S_Sym 
)
inline

2x2 polar decomposition.

Parameters
[in]Amatrix.
[out]RRobustly a rotation matrix in givens form
[out]S_SymSymmetric. Whole matrix is stored

Whole matrix S is stored since its faster to calculate due to simd vectorization Polar guarantees negative sign is on the small magnitude singular value. S is guaranteed to be the closest one to identity. R is guaranteed to be the closest rotation to A.

Definition at line 403 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::c, JIXIE::GivensRotation< T >::rowRotation(), JIXIE::GivensRotation< T >::s, and x.

Referenced by polarDecomposition(), and singularValueDecomposition().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ process()

template<int t, class T >
void JIXIE::process ( Eigen::Matrix< T, 3, 3 > &  B,
Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 1 > &  sigma,
Eigen::Matrix< T, 3, 3 > &  V 
)
inline

Helper function of 3X3 SVD for processing 2X2 SVD.

Definition at line 571 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::columnRotation(), JIXIE::GivensRotation< T >::rowi, JIXIE::GivensRotation< T >::rowk, singularValueDecomposition(), and V.

Here is the call graph for this function:

◆ singularValueDecomposition() [1/3]

template<class T >
int JIXIE::singularValueDecomposition ( const Eigen::Matrix< T, 3, 3 > &  A,
Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 1 > &  sigma,
Eigen::Matrix< T, 3, 3 > &  V,
tol = 1024 * std::numeric_limits<T>::epsilon() 
)
inline

3X3 SVD (singular value decomposition) A=USV'

Parameters
[in]AInput matrix.
[out]Uis a rotation matrix.
[out]sigmaDiagonal matrix, sorted with decreasing magnitude. The third one can be negative.
[out]Vis a rotation matrix.

Do implicit shift QR until A^T A is block diagonal

Handle the cases of one of the alphas and betas being 0 Sorted by ease of handling and then frequency of occurrence

If B is of form x x 0 0 x 0 0 0 x

If B is of form x 0 0 0 x x 0 0 x

If B is of form x x 0 0 0 x 0 0 x

Reduce B to x x 0 0 0 0 0 0 x

If B is of form x x 0 0 x x 0 0 0

Reduce B to x x + 0 x 0 0 0 0

Reduce B to x x 0

  • x 0 0 0 0

If B is of form 0 x 0 0 x x 0 0 x

Reduce B to 0 0 + 0 x x 0 0 x

Reduce B to 0 0 0 0 x x 0 + x

Definition at line 688 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::columnRotation(), JIXIE::GivensRotation< T >::compute(), JIXIE::GivensRotation< T >::computeUnconventional(), makeUpperBidiag(), JIXIE::GivensRotation< T >::rowRotation(), V, wilkinsonShift(), and zeroChase().

Here is the call graph for this function:

◆ singularValueDecomposition() [2/3]

template<class TA , class TU , class Ts , class TV >
enable_if_t< isSize< TA >(2, 2) &&isSize< TU >(2, 2) &&isSize< TV >(2, 2) &&isSize< Ts >(2, 1)> JIXIE::singularValueDecomposition ( const Eigen::MatrixBase< TA > &  A,
const Eigen::MatrixBase< TU > &  U,
const Eigen::MatrixBase< Ts > &  Sigma,
const Eigen::MatrixBase< TV > &  V,
const ScalarType< TA >  tol = 64 * std::numeric_limits<ScalarType<TA>>::epsilon() 
)
inline

2x2 SVD (singular value decomposition) A=USV'

Parameters
[in]AInput matrix.
[out]URobustly a rotation matrix.
[out]SigmaVector of singular values sorted with decreasing magnitude. The second one can be negative.
[out]VRobustly a rotation matrix.

Definition at line 528 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::fill(), singularValueDecomposition(), and V.

Here is the call graph for this function:

◆ singularValueDecomposition() [3/3]

template<class TA , class T , class Ts >
enable_if_t< isSize< TA >(2, 2) &&isSize< Ts >(2, 1)> JIXIE::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() 
)
inline

2x2 SVD (singular value decomposition) A=USV'

Parameters
[in]AInput matrix.
[out]URobustly a rotation matrix in Givens form
[out]SigmaVector of singular values sorted with decreasing magnitude. The second one can be negative.
[out]VRobustly a rotation matrix in Givens form

Definition at line 456 of file ImplicitQRSVD.h.

References polarDecomposition(), V, x, y, and z.

Referenced by polyfem::utils::AutoFlipSVD< MatrixType >::compute(), polarDecomposition(), process(), and singularValueDecomposition().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sort() [1/2]

template<int t, class T >
enable_if_t< t==0 > JIXIE::sort ( Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 1 > &  sigma,
Eigen::Matrix< T, 3, 3 > &  V 
)

Helper function of 3X3 SVD for sorting singular values.

Definition at line 600 of file ImplicitQRSVD.h.

References flipSign(), and V.

Here is the call graph for this function:

◆ sort() [2/2]

template<int t, class T >
enable_if_t< t==1 > JIXIE::sort ( Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 1 > &  sigma,
Eigen::Matrix< T, 3, 3 > &  V 
)

Helper function of 3X3 SVD for sorting singular values.

Definition at line 642 of file ImplicitQRSVD.h.

References flipSign(), and V.

Here is the call graph for this function:

◆ wilkinsonShift()

template<class T >
T JIXIE::wilkinsonShift ( const T  a1,
const T  b1,
const T  a2 
)

compute wilkinsonShift of the block a1 b1 b1 a2 based on the wilkinsonShift formula mu = c + d - sign (d) \ sqrt (d*d + b*b), where d = (a-c)/2

Definition at line 553 of file ImplicitQRSVD.h.

Referenced by singularValueDecomposition().

Here is the caller graph for this function:

◆ zeroChase()

template<class T >
void JIXIE::zeroChase ( Eigen::Matrix< T, 3, 3 > &  H,
Eigen::Matrix< T, 3, 3 > &  U,
Eigen::Matrix< T, 3, 3 > &  V 
)
inline

zero chasing the 3X3 matrix to bidiagonal form original form of H: x x 0 x x x 0 0 x after zero chase: x x 0 0 x x 0 0 x

Reduce H to of form x x + 0 x x 0 0 x

Reduce H to of form x x 0 0 x x 0 + x Can calculate r2 without multiplying by r1 since both entries are in first two rows thus no need to divide by sqrt(a^2+b^2)

Reduce H to of form x x 0 0 x x 0 0 x

Definition at line 253 of file ImplicitQRSVD.h.

References JIXIE::GivensRotation< T >::columnRotation(), JIXIE::GivensRotation< T >::compute(), JIXIE::GivensRotation< T >::rowRotation(), and V.

Referenced by makeUpperBidiag(), and singularValueDecomposition().

Here is the call graph for this function:
Here is the caller graph for this function: