template<
typename T,
class MeshType = HalfEdgeMesh<T>,
class Solver = Eigen::SparseLU<Eigen::SparseMatrix<T>, Eigen::COLAMDOrdering<int>>, std::enable_if_t< std::is_floating_point_v< T >,
bool > = true>
class OpenABF::AngleBasedLSCM< T, MeshType, Solver, >
Compute parameterized mesh using Angle-based LSCM.
Computes a least-squares conformal parameterization of a mesh. Unlike the original LSCM algorithm, this class ignores the 3D vertex positions and instead uses the angle associated with the mesh's edge trait (MeshType::EdgeTraits::alpha) to calculate the initial per-triangle edge lengths. Without previously modifying the angles of the provided mesh, this class produces the same result as a vertex-based LSCM implementation. However, by first processing the mesh with a parameterized angle optimizer, such as ABFPlusPlus, the parameterization can be improved, sometimes significantly.
Implements the angle-based variant of "Least squares conformal maps for
automatic texture atlas generation" by Lévy et al. (2002) [1].
- Template Parameters
-
| T | Floating-point type |
| MeshType | HalfEdgeMesh type which implements the default mesh traits |
| Solver | A solver implementing the Eigen Sparse solver concept and templated on Eigen::SparseMatrix<T>. The default SparseLU is robust but slow for large meshes. For iterative solving, prefer Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower|EigenUpper> over the default Lower-only variant: the Lower|Upper template argument enables Eigen's full-matrix SpMV code path, which is faster and — when compiled with OpenMP — multi-threaded. Using only Lower (the Eigen default) routes through selfadjointView<Lower>, which is a different internal code path that is never OpenMP-parallelized regardless of Eigen::setNbThreads(). |
- Examples
- AlternativeSolvers.cpp, BasicFlattening.cpp, CloneMesh.cpp, FlattenTool.cpp, and SplitPath.cpp.