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...

#include <OpenABF/OpenABF.hpp>

Public Types

using Mesh = MeshType
 Mesh type alias.
 
using PinMap = detail::lscm::PinMap< T >
 Per-pin entry: (mesh vertex index, target UV).
 

Public Member Functions

void set_pins (PinMap pins)
 Set the explicit pin set used by compute()
 
void setPinnedVertices (std::size_t pin0Idx, std::size_t pin1Idx)
 Deprecated: set a pin pair by index using the LSCM axis-snap convention at compute time.
 
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, const PinMap &pins)
 Compute the parameterized mesh with caller-specified pin UVs.
 
static void Compute (typename Mesh::Pointer &mesh, std::size_t pin0Idx, std::size_t pin1Idx)
 Deprecated: compute with an explicit pin pair by index.
 

Static Private Member Functions

static void ComputeImpl (typename Mesh::Pointer &mesh, const PinMap &pins)
 Core solver: build the LSCM system from the PinMap, solve for free vertices, and write UVs back to the mesh.
 

Private Attributes

std::optional< PinMappins_
 
std::optional< std::pair< std::size_t, std::size_t > > legacy_pin_indices_
 

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 and the fastest option we have benchmarked on the LSCM normal equations across 50k–200k-face meshes (see examples/src/BenchmarkFlattening.cpp); for very large meshes where SparseLU exhausts memory, switch to an iterative solver via this template parameter.

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().

IncompleteCholesky is also available as a preconditioner via Eigen::ConjugateGradient<..., Eigen::IncompleteCholesky<T>>, but in our benchmarks the per-iteration overhead of applying IC on the LSCM normal equations dwarfs the iteration-count savings it provides over the default DiagonalPreconditioner (Jacobi), and both are much slower than SparseLU at the mesh sizes we target. Use IC only if you have profiled it favorably against Diagonal for your specific mesh class.

Examples
AlternativeSolvers.cpp, BasicFlattening.cpp, BenchmarkFlattening.cpp, CloneMesh.cpp, FlattenTool.cpp, MultiChartFlatten.cpp, and SplitPath.cpp.

Member Typedef Documentation

◆ PinMap

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>
using OpenABF::AngleBasedLSCM< T, MeshType, Solver, >::PinMap = detail::lscm::PinMap<T>

Per-pin entry: (mesh vertex index, target UV).

A PinMap of size N ≥ 2 specifies an explicit LSCM pin set; each pinned vertex's final UV equals its entry in the map.

Member Function Documentation

◆ Compute() [1/3]

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, placing pin0 at the UV origin and pin1 at distance |p1 - p0| along the dominant world-axis of (p1 - p0).

Exceptions
MeshExceptionIf pin selection fails (no boundary vertices).
SolverExceptionIf matrix cannot be decomposed or if solver fails to find a solution.
Examples
AlternativeSolvers.cpp, BasicFlattening.cpp, CloneMesh.cpp, FlattenTool.cpp, MultiChartFlatten.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, placing pin0 at the UV origin and pin1 at distance |p1 - p0| along the dominant world-axis of (p1 - p0).

Exceptions
MeshExceptionIf pin selection fails (no boundary vertices).
SolverExceptionIf matrix cannot be decomposed or if solver fails to find a solution.

◆ Compute() [2/3]

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,
const PinMap pins 
)
inlinestatic

Compute the parameterized mesh with caller-specified pin UVs.

Parameters
meshTriangle mesh whose vertex positions will be overwritten with computed 2D UV coordinates (z component set to 0).
pinsSequence of (vertex index, target UV) pairs. Must contain at least two unique, in-range vertex indices. Each pinned vertex's final UV equals the supplied target verbatim.
Exceptions
std::invalid_argumentIf pins has fewer than two entries, a duplicate index, or an out-of-range index.
SolverExceptionIf matrix cannot be decomposed or if solver fails to find a solution.

◆ Compute() [3/3]

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

Deprecated: compute with an explicit pin pair by index.

Builds a 2-entry PinMap using the LSCM axis-snap convention (pin0 at the UV origin, pin1 on the dominant world-axis of (p1 - p0)) and dispatches to the PinMap path.

Deprecated:
Prefer Compute(mesh, PinMap). This overload will be removed in version 3.0.

◆ set_pins()

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, >::set_pins ( PinMap  pins)
inline

Set the explicit pin set used by compute()

The PinMap must contain at least two unique, in-range vertex indices; compute() rejects malformed inputs at solve time.

◆ setPinnedVertices()

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, >::setPinnedVertices ( std::size_t  pin0Idx,
std::size_t  pin1Idx 
)
inline

Deprecated: set a pin pair by index using the LSCM axis-snap convention at compute time.

Deprecated:
Prefer set_pins(PinMap). This overload will be removed in version 3.0.

Member Data Documentation

◆ legacy_pin_indices_

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, >::legacy_pin_indices_
private

Deprecated: legacy two-pin index pair set via setPinnedVertices.

◆ pins_

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<PinMap> OpenABF::AngleBasedLSCM< T, MeshType, Solver, >::pins_
private

Optional explicit pin set configured via set_pins().