|
OpenABF 2.1.0
|
Compute parameterized mesh using Hierarchical LSCM. More...
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 | set_level_ratio (std::size_t ratio) |
| Set the vertex ratio between consecutive hierarchy levels (default: 10) | |
| void | set_min_coarse_vertices (std::size_t count) |
| Set the minimum vertex count at the coarsest level (default: 100) | |
| void | setLevelRatio (std::size_t ratio) |
| void | setMinCoarseVertices (std::size_t count) |
| 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, const PinMap &pins) |
| Compute with caller-specified pin UVs. | |
| static void | Compute (typename Mesh::Pointer &mesh, const PinMap &pins, std::size_t levelRatio, std::size_t minCoarseVerts) |
| Compute with caller-specified pin UVs and hierarchy tuning. | |
| 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 | 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, const PinMap &pins, std::size_t levelRatio=10, std::size_t minCoarseVerts=100) |
| Core hierarchical LSCM solve given a resolved PinMap. | |
Private Attributes | |
| std::optional< PinMap > | pins_ |
| std::optional< std::pair< std::size_t, std::size_t > > | legacy_pin_indices_ |
| std::size_t | level_ratio_ {10} |
| std::size_t | min_coarse_vertices_ {100} |
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.
| T | Floating-point type |
| MeshType | HalfEdgeMesh type which implements the default mesh traits |
| Solver | An 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. |
Solver template parameter is valid but disables warm-starting — the initial guess is ignored.Solver template parameter, not AngleBasedLSCM's default (SparseLU). The result is numerically equivalent but may differ in convergence behavior from a plain AngleBasedLSCM call. | using OpenABF::HierarchicalLSCM< 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 pin survives every decimation level — DecimationMesh::try_collapse rejects any collapse that would remove a pinned vertex.
|
inlinestatic |
Compute with automatic pin selection.
Selects pins identically to AngleBasedLSCM::Compute() — first boundary vertex and its boundary-edge neighbor, placed via LSCM axis-snap.
| MeshException | if the mesh has no boundary vertices (pin selection fails) or the mesh is otherwise invalid |
| SolverException | if any hierarchy level fails to solve |
|
inline |
Compute parameterization using instance configuration.
If set_pins() was called, uses the supplied PinMap. Otherwise, if the deprecated setPinnedVertices() was called, builds a PinMap from those indices via the LSCM axis-snap convention. Otherwise auto-selects two boundary vertices using the same logic as Compute(mesh).
| MeshException | if pin selection fails (no boundary vertices) |
| SolverException | if any hierarchy level fails to solve |
|
inlinestatic |
Compute with caller-specified pin UVs.
| std::invalid_argument | If pins has fewer than two entries, a duplicate index, or an out-of-range index. |
| SolverException | If any hierarchy level fails to solve. |
|
inlinestatic |
Compute with caller-specified pin UVs and hierarchy tuning.
Static counterpart to configuring set_level_ratio() and set_min_coarse_vertices() on an instance.
| levelRatio | Target vertex ratio between consecutive hierarchy levels. Must be >= 2. |
| minCoarseVerts | Minimum vertex count at the coarsest level. Must be >= 3. |
| std::invalid_argument | If pins is invalid (see PinMap overload), or if levelRatio < 2, or if minCoarseVerts < 3. |
| SolverException | If any hierarchy level fails to solve. |
Compute(mesh, levelRatio, minCoarseVerts) is not provided because its signature would collide with the deprecated Compute(mesh, pin0Idx, pin1Idx) overload. It will be added when that overload is removed in 3.0; until then, use the instance API for auto-pin selection with custom tuning.
|
inlinestatic |
Deprecated: compute with an explicit pin pair by index.
Builds a 2-entry PinMap using the LSCM axis-snap convention and dispatches to the PinMap path.
Compute(mesh, PinMap). This overload will be removed in version 3.0.
|
inlinestaticprivate |
Core hierarchical LSCM solve given a resolved PinMap.
Builds the mesh hierarchy with each pinned vertex flagged as non-collapsible, solves LSCM at the coarsest level, then prolongates and refines at each finer level. Pin UVs come from the PinMap verbatim; non-pin vertex UVs are written by the solve.
|
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.
|
inline |
set_level_ratio; will be removed in 3.0.
|
inline |
set_min_coarse_vertices; will be removed in 3.0.
|
inline |
Deprecated: set a pin pair by index using the LSCM axis-snap convention at compute time.
set_pins(PinMap). This overload will be removed in version 3.0.
|
private |
Deprecated: legacy two-pin index pair set via setPinnedVertices.
|
private |
Ratio of vertices between consecutive hierarchy levels
|
private |
Minimum vertex count for the coarsest hierarchy level
|
private |
Optional explicit pin set configured via set_pins().