OpenABF 2.1.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Static Private Member Functions | Private Attributes | List of all members
OpenABF::AngleBasedLSCM< T, MeshType, Solver, > Class Template Reference

Compute parameterized mesh using Angle-based LSCM. More...

Public Types

using Mesh = MeshType
 Mesh type alias.
 

Public Member Functions

void setPinnedVertices (std::size_t pin0Idx, std::size_t pin1Idx)
 Set the pinned vertex indices used by compute()
 
void compute (typename Mesh::Pointer &mesh) const
 Compute the parameterized mesh using automatic pin selection.
 

Static Public Member Functions

static void Compute (typename Mesh::Pointer &mesh)
 Compute the parameterized mesh using automatic pin selection.
 
static void Compute (typename Mesh::Pointer &mesh, std::size_t pin0Idx, std::size_t pin1Idx)
 Compute the parameterized mesh with explicit pinned vertex indices.
 

Static Private Member Functions

static void ComputeImpl (typename Mesh::Pointer &mesh, const typename Mesh::VertPtr &p0, const typename Mesh::VertPtr &p1)
 Core solver: place p0/p1 on the UV axes then solve for free vertices.
 

Private Attributes

std::optional< std::pair< std::size_t, std::size_t > > pinnedVertices_
 

Detailed Description

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
TFloating-point type
MeshTypeHalfEdgeMesh type which implements the default mesh traits
SolverA 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.

Member Function Documentation

◆ Compute() [1/2]

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>
static void OpenABF::AngleBasedLSCM< T, MeshType, Solver, >::Compute ( typename Mesh::Pointer &  mesh)
inlinestatic

Compute the parameterized mesh using automatic pin selection.

Selects the first boundary vertex and its boundary-edge neighbor as pinned vertices.

Exceptions
MeshExceptionIf pinned vertex is not on boundary.
SolverExceptionIf matrix cannot be decomposed or if solver fails to find a solution.
Examples
AlternativeSolvers.cpp, BasicFlattening.cpp, CloneMesh.cpp, FlattenTool.cpp, and SplitPath.cpp.

◆ compute()

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>
void OpenABF::AngleBasedLSCM< T, MeshType, Solver, >::compute ( typename Mesh::Pointer &  mesh) const
inline

Compute the parameterized mesh using automatic pin selection.

Selects the first boundary vertex and its boundary-edge neighbor as pinned vertices.

Exceptions
MeshExceptionIf pinned vertex is not on boundary.
SolverExceptionIf matrix cannot be decomposed or if solver fails to find a solution.

◆ Compute() [2/2]

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>
static void OpenABF::AngleBasedLSCM< T, MeshType, Solver, >::Compute ( typename Mesh::Pointer &  mesh,
std::size_t  pin0Idx,
std::size_t  pin1Idx 
)
inlinestatic

Compute the parameterized mesh with explicit pinned vertex indices.

Parameters
pin0IdxIndex of the first pinned vertex (placed at the UV origin)
pin1IdxIndex of the second pinned vertex (placed on the nearest axis)
Exceptions
SolverExceptionIf matrix cannot be decomposed or if solver fails to find a solution.

Member Data Documentation

◆ pinnedVertices_

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>
std::optional<std::pair<std::size_t, std::size_t> > OpenABF::AngleBasedLSCM< T, MeshType, Solver, >::pinnedVertices_
private

Optional explicit pin pair set via setPinnedVertices()