| 
    libigl v2.5.0
    
   | 
 
Functions | |
| template<typename EigsScalar , typename DerivedU , typename DerivedS , typename Solver = Eigen::SparseLU<Eigen::SparseMatrix<EigsScalar>>> | |
| bool | eigs (const Eigen::SparseMatrix< EigsScalar > &A, const Eigen::SparseMatrix< EigsScalar > &B, const int k, const igl::EigsType type, Eigen::PlainObjectBase< DerivedU > &U, Eigen::PlainObjectBase< DerivedS > &S) | 
| Act like MATLAB's eigs function.   | |
| template<typename EigsScalar , typename DerivedU , typename DerivedS , typename Solver = Eigen::SparseLU<Eigen::SparseMatrix<EigsScalar>>> | |
| bool | eigs (const Eigen::SparseMatrix< EigsScalar > &A, const Eigen::SparseMatrix< EigsScalar > &B, const int k, const EigsScalar sigma, Eigen::PlainObjectBase< DerivedU > &U, Eigen::PlainObjectBase< DerivedS > &S) | 
| This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.   | |
| template<typename DerivedV , typename DerivedF , typename DerivedV_uv > | |
| bool | lscm (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, Eigen::PlainObjectBase< DerivedV_uv > &V_uv) | 
| Compute a free-boundary least-squares conformal map parametrization.   | |
| bool igl::spectra::eigs | ( | const Eigen::SparseMatrix< EigsScalar > & | A, | 
| const Eigen::SparseMatrix< EigsScalar > & | B, | ||
| const int | k, | ||
| const igl::EigsType | type, | ||
| Eigen::PlainObjectBase< DerivedU > & | U, | ||
| Eigen::PlainObjectBase< DerivedS > & | S | ||
| ) | 
Act like MATLAB's eigs function.
Compute the first/last k eigen pairs of the generalized eigen value problem:
A u = s B u
Solutions are approximate and sorted.
Ideally one should use ARPACK and the Eigen unsupported ARPACK module. This implementation does simple, naive power iterations.
| [in] | A | #A by #A symmetric matrix | 
| [in] | B | #A by #A symmetric positive-definite matrix | 
| [in] | k | number of eigen pairs to compute | 
| [in] | type | whether to extract from the high or low end | 
| [out] | sU | #A by k list of sorted eigen vectors (descending) | 
| [out] | sS | k list of sorted eigen values (descending) | 
| bool igl::spectra::eigs | ( | const Eigen::SparseMatrix< EigsScalar > & | A, | 
| const Eigen::SparseMatrix< EigsScalar > & | B, | ||
| const int | k, | ||
| const EigsScalar | sigma, | ||
| Eigen::PlainObjectBase< DerivedU > & | U, | ||
| Eigen::PlainObjectBase< DerivedS > & | S | ||
| ) | 
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
| [in] | sigma | shift to apply to A, as in A ← A + sigma B | 
| bool igl::spectra::lscm | ( | const Eigen::MatrixBase< DerivedV > & | V, | 
| const Eigen::MatrixBase< DerivedF > & | F, | ||
| Eigen::PlainObjectBase< DerivedV_uv > & | V_uv | ||
| ) | 
Compute a free-boundary least-squares conformal map parametrization.
Equivalently derived in "Intrinsic Parameterizations of Surface Meshes" [Desbrun et al. 2002] and "Least Squares Conformal Maps for Automatic Texture Atlas Generation" [Lévy et al. 2002], though this implementation follows the derivation in: "Spectral Conformal Parameterization" [Mullen et al. 2008] Free boundary version. "Spectral Conformal Parameterization" using Eigen decomposition. Assumes mesh is a single connected component topologically equivalent to a chunk of the plane.
| [in] | V | #V by 3 list of mesh vertex positions | 
| [in] | F | #F by 3 list of mesh faces (must be triangles) | 
| [out] | UV | #V by 2 list of 2D mesh vertex positions in UV space |