|
OpenABF 2.1.0
|
Compute parameterized interior angles using Angle-based flattening. More...
#include <OpenABF/OpenABF.hpp>
Public Types | |
| using | Mesh = MeshType |
| Mesh type alias. | |
Public Member Functions | |
| void | setMaxIterations (const std::size_t it) |
| Set the maximum number of iterations. | |
| void | setGradientThreshold (T t) |
| Set the gradient convergence threshold. | |
| auto | gradient () const -> T |
| Get the mesh gradient. | |
| auto | iterations () const -> std::size_t |
| Get the number of iterations of the last computation. | |
| void | compute (typename Mesh::Pointer &mesh) |
| Compute parameterized interior angles. | |
Static Public Member Functions | |
| static void | Compute (typename Mesh::Pointer &mesh, std::size_t &iters, T &gradient, const std::size_t maxIters=10, T gradThreshold=T(0.001)) |
| Compute parameterized interior angles. | |
| static void | Compute (typename Mesh::Pointer &mesh) |
| Compute parameterized interior angles. | |
Protected Attributes | |
| T | grad_ {0} |
| std::size_t | iters_ {0} |
| std::size_t | maxIters_ {10} |
| T | gradThreshold_ {0.001} |
Compute parameterized interior angles using Angle-based flattening.
Iteratively computes a new set of interior angles which minimize the total angular error of the parameterized mesh. This follows the standard angled-based flattening formulation, which directly solves for the objective functions and constraints. ABFPlusPlus is generally preferred as it dramatically simplifies the size of the solved problem without introducing more error.
This class does not compute a parameterized mesh. Rather, it calculates the optimal interior angles for such a mesh. To convert this information into a full parameterization, pass the processed HalfEdgeMesh to AngleBasedLSCM.
Implements "Parameterization of faceted surfaces for meshing using angle-based flattening" by Sheffer and de Sturler (2001) [4].
| T | Floating-point type |
| MeshType | HalfEdgeMesh type which implements the ABF traits |
| Solver | A solver implementing the Eigen Sparse solver concept and templated on Eigen::SparseMatrix<T> |
|
inline |
Compute parameterized interior angles.
| SolverException | If matrix cannot be decomposed or if solver fails to find a solution. |
| MeshException | If mesh gradient cannot be calculated. |
|
inlinestatic |
Compute parameterized interior angles.
| SolverException | If matrix cannot be decomposed or if solver fails to find a solution. |
| MeshException | If mesh gradient cannot be calculated. |
|
inlinestatic |
Compute parameterized interior angles.
| SolverException | If matrix cannot be decomposed or if solver fails to find a solution. |
| MeshException | If mesh gradient cannot be calculated. |
|
inline |
Get the mesh gradient.
Note: Result is only valid after running compute().
|
inline |
Get the number of iterations of the last computation.
Note: Result is only valid after running compute().
|
protected |
Gradient
|
protected |
Gradient convergence threshold
|
protected |
Number of executed iterations