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::HierarchicalLSCM< T, MeshType, Solver, > Class Template Reference

Compute parameterized mesh using Hierarchical 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 setLevelRatio (std::size_t ratio)
 Set the vertex ratio between consecutive hierarchy levels (default: 10)
 
void setMinCoarseVertices (std::size_t count)
 Set the minimum vertex count at the coarsest level (default: 100)
 
void compute (typename Mesh::Pointer &mesh) const
 Compute parameterization using instance configuration.
 

Static Public Member Functions

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

Static Private Member Functions

static void AutoSelectPins (const typename Mesh::Pointer &mesh, std::size_t &p0, std::size_t &p1)
 
static void CopyAnglesFromOriginal (const typename Mesh::Pointer &original, const typename HalfEdgeMesh< T >::Pointer &levelMesh)
 Copy edge angles from the original mesh to a level mesh.
 
static void ComputeImpl (typename Mesh::Pointer &mesh, std::size_t pin0Idx, std::size_t pin1Idx, std::size_t levelRatio=10, std::size_t minCoarseVerts=100)
 

Private Attributes

std::optional< std::pair< std::size_t, std::size_t > > pinnedVertices_
 
std::size_t levelRatio_ {10}
 
std::size_t minCoarseVertices_ {100}
 

Detailed Description

template<typename T, class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
class OpenABF::HierarchicalLSCM< T, MeshType, Solver, >

Compute parameterized mesh using Hierarchical LSCM.

Implements the HLSCM algorithm from Ray & Lévy, "Hierarchical Least Squares Conformal Map" (2003) [3]. Uses cascadic multigrid to accelerate LSCM convergence: the mesh is decimated into a hierarchy, LSCM is solved on the coarsest level, and the solution is prolongated and refined at each finer level using conjugate gradient with the prolongated UVs as initial guess.

For small meshes the hierarchy has a single level and HLSCM degrades gracefully to a standard LSCM solve.

Template Parameters
TFloating-point type
MeshTypeHalfEdgeMesh type which implements the default mesh traits
SolverAn Eigen iterative or direct solver. Iterative solvers (ConjugateGradient, LeastSquaresConjugateGradient) support warm- starting from the coarser-level solution; direct solvers ignore the initial guess. Defaults to ConjugateGradient<SparseMatrix<T>, Lower|Upper>, which solves the normal equations (AᵀA x = Aᵀb) of the overdetermined LSCM system. The Lower|Upper flag enables OpenMP-parallelized SpMV on the symmetric AᵀA matrix, giving the best multi-thread performance. LeastSquaresConjugateGradient is a valid alternative; it operates on the same normal equations internally but without the OpenMP benefit.
Note
Solver default differs from AngleBasedLSCM. AngleBasedLSCM defaults to SparseLU (a direct solver); HierarchicalLSCM defaults to ConjugateGradient so it can warm-start from the coarser-level solution. Using a direct solver via the Solver template parameter is valid but disables warm-starting — the initial guess is ignored.
Single-level fallback. When the mesh is too small to decimate (all vertices are boundary, or minCoarseVertices is already reached), HLSCM falls back to a standard LSCM solve using this class's Solver template parameter, not AngleBasedLSCM's default (SparseLU). The result is numerically equivalent but may differ in convergence behavior from a plain AngleBasedLSCM call.

Member Function Documentation

◆ AutoSelectPins()

template<typename T , class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
static void OpenABF::HierarchicalLSCM< T, MeshType, Solver, >::AutoSelectPins ( const typename Mesh::Pointer &  mesh,
std::size_t &  p0,
std::size_t &  p1 
)
inlinestaticprivate

Select two pinned boundary vertices (same logic as AngleBasedLSCM)

◆ Compute() [1/2]

template<typename T , class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
static void OpenABF::HierarchicalLSCM< T, MeshType, Solver, >::Compute ( typename Mesh::Pointer &  mesh)
inlinestatic

Compute with automatic pin selection.

Selects pins identically to AngleBasedLSCM::Compute().

Exceptions
MeshExceptionif the mesh has no boundary vertices (pin selection fails) or the mesh is otherwise invalid
SolverExceptionif any hierarchy level fails to solve

◆ compute()

template<typename T , class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
void OpenABF::HierarchicalLSCM< T, MeshType, Solver, >::compute ( typename Mesh::Pointer &  mesh) const
inline

Compute parameterization using instance configuration.

Uses the pinned vertices, level ratio, and minimum coarse vertices configured via setPinnedVertices(), setLevelRatio(), and setMinCoarseVertices(). If no pins are set, selects them automatically using the same boundary-walk logic as Compute(mesh).

Exceptions
MeshExceptionif pin selection fails (no boundary vertices)
SolverExceptionif any hierarchy level fails to solve

◆ Compute() [2/2]

template<typename T , class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
static void OpenABF::HierarchicalLSCM< T, MeshType, Solver, >::Compute ( typename Mesh::Pointer &  mesh,
std::size_t  pin0Idx,
std::size_t  pin1Idx 
)
inlinestatic

Compute with explicit pinned vertex indices.

Exceptions
SolverExceptionif any hierarchy level fails to solve

◆ CopyAnglesFromOriginal()

template<typename T , class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
static void OpenABF::HierarchicalLSCM< T, MeshType, Solver, >::CopyAnglesFromOriginal ( const typename Mesh::Pointer &  original,
const typename HalfEdgeMesh< T >::Pointer &  levelMesh 
)
inlinestaticprivate

Copy edge angles from the original mesh to a level mesh.

At the finest hierarchy level (k=0), the level mesh has the same face/edge structure as the original mesh. If ABF was run beforehand, we must use the ABF-optimized angles rather than recomputing from geometry. Coarser levels always use geometry angles.

Member Data Documentation

◆ pinnedVertices_

template<typename T , class MeshType = HalfEdgeMesh<T>, class Solver = Eigen::ConjugateGradient<Eigen::SparseMatrix<T>, Eigen::Lower | Eigen::Upper>, std::enable_if_t< std::is_floating_point_v< T >, bool > = true>
std::optional<std::pair<std::size_t, std::size_t> > OpenABF::HierarchicalLSCM< T, MeshType, Solver, >::pinnedVertices_
private

Optional explicit pin pair