CGAL 6.0 - Planar Parameterization of Triangulated Surface Meshes
|
#include <CGAL/Surface_mesh_parameterization/Iterative_authalic_parameterizer_3.h>
The class Iterative_authalic_parameterizer_3
implements the Iterative Parameterization algorithm, as described by Jain et al.
[6].
This parameterization is a fixed border parameterization and is part of the authalic parameterization family, meaning that it aims to minimize area distortion between the input surface mesh and the parameterized output. More precisely, the approach used by this parameterizer is to iteratively redistribute the \( L_2\) stretch - as defined by Sander et al. [11] - over the mesh.
TriangleMesh_ | must be a model of FaceGraph . |
BorderParameterizer_ | is a strategy to parameterize the surface border and must be a model of Parameterizer_3 .Default: This class parameterizes the border of a 3D surface onto a circle, with an arc-length parameterizatio... Definition: Circular_border_parameterizer_3.h:244 |
SolverTraits_ | must be a model of SparseLinearAlgebraTraits_d .Note that the system is not symmetric because border vertices are not removed from the system. Default: If Eigen 3.1 (or greater) is available and CGAL_EIGEN3_ENABLED is defined, then an overload of Eigen_solver_traits is provided as default parameter: Eigen::BiCGSTAB<Eigen_sparse_matrix<double>::EigenType,
Eigen::IncompleteLUT< double > > >
|
CGAL::Surface_mesh_parameterization::Discrete_authalic_parameterizer_3<TriangleMesh, BorderParameterizer, SolverTraits>
CGAL::Surface_mesh_parameterization::ARAP_parameterizer_3<TriangleMesh, BorderParameterizer, SolverTraits>
Public Member Functions | |
Iterative_authalic_parameterizer_3 (Border_parameterizer border_parameterizer=Border_parameterizer(), Solver_traits sparse_la=Solver_traits()) | |
Constructor. | |
template<typename VertexUVmap , typename VertexIndexMap > | |
void | initialize_system_from_mesh_border (Matrix &A, Vector &Bu, Vector &Bv, const Triangle_mesh &tmesh, halfedge_descriptor bhd, VertexUVmap uvmap, VertexIndexMap vimap) const |
initializes A , Bu , and Bv after border parameterization. | |
template<typename VertexUVmap , typename VertexIndexMap , typename VertexParameterizedMap > | |
Error_code | parameterize (Triangle_mesh &tmesh, halfedge_descriptor bhd, VertexUVmap uvmap, VertexIndexMap vimap, VertexParameterizedMap vpmap, const unsigned int &iterations=15) |
computes a one-to-one mapping from a triangular 3D surface mesh to a piece of the 2D space. | |
Protected Member Functions | |
Border_parameterizer & | get_border_parameterizer () |
Get the object that maps the surface's border onto a 2D space. | |
Solver_traits & | get_linear_algebra_traits () |
Get the sparse linear algebra (traits object to access the linear system). | |
CGAL::Surface_mesh_parameterization::Iterative_authalic_parameterizer_3< TriangleMesh_, BorderParameterizer_, SolverTraits_ >::Iterative_authalic_parameterizer_3 | ( | Border_parameterizer | border_parameterizer = Border_parameterizer() , |
Solver_traits | sparse_la = Solver_traits() |
||
) |
Constructor.
border_parameterizer | Object that maps the surface's border to 2D space |
sparse_la | Traits object to access a sparse linear system |
void CGAL::Surface_mesh_parameterization::Iterative_authalic_parameterizer_3< TriangleMesh_, BorderParameterizer_, SolverTraits_ >::initialize_system_from_mesh_border | ( | Matrix & | A, |
Vector & | Bu, | ||
Vector & | Bv, | ||
const Triangle_mesh & | tmesh, | ||
halfedge_descriptor | bhd, | ||
VertexUVmap | uvmap, | ||
VertexIndexMap | vimap | ||
) | const |
initializes A
, Bu
, and Bv
after border parameterization.
Fill the border vertices' lines in both linear systems: "u = constant" and "v = constant".
VertexUVmap | must be a model of ReadWritePropertyMap with boost::graph_traits<Triangle_mesh>::vertex_descriptor as key type and Point_2 (type deduced from Triangle_mesh using Kernel_traits ) as value type. |
VertexIndexMap | must be a model of ReadablePropertyMap with boost::graph_traits<Triangle_mesh>::vertex_descriptor as key type and a unique integer as value type. |
A | the matrix in both linear system |
Bu | the right hand side vector in the linear system of x coordinates |
Bv | the right hand side vector in the linear system of y coordinates |
tmesh | a triangulated surface |
bhd | a halfedge descriptor on the boundary of mesh |
uvmap | an instantiation of the class VertexUVmap |
vimap | an instantiation of the class VertexIndexMap |
vimap
must be initialized). A
, Bu
, and Bv
must be allocated. Error_code CGAL::Surface_mesh_parameterization::Iterative_authalic_parameterizer_3< TriangleMesh_, BorderParameterizer_, SolverTraits_ >::parameterize | ( | Triangle_mesh & | tmesh, |
halfedge_descriptor | bhd, | ||
VertexUVmap | uvmap, | ||
VertexIndexMap | vimap, | ||
VertexParameterizedMap | vpmap, | ||
const unsigned int & | iterations = 15 |
||
) |
computes a one-to-one mapping from a triangular 3D surface mesh to a piece of the 2D space.
The mapping is piecewise linear (linear in each triangle). The result is the (u,v)
pair image of each vertex of the 3D surface.
VertexUVmap | must be a model of ReadWritePropertyMap with boost::graph_traits<Triangle_mesh>::vertex_descriptor as key type and Point_2 (type deduced from Triangle_mesh using Kernel_traits ) as value type. |
VertexIndexMap | must be a model of ReadablePropertyMap with boost::graph_traits<Triangle_mesh>::vertex_descriptor as key type and a unique integer as value type. |
VertexParameterizedMap | must be a model of ReadWritePropertyMap with boost::graph_traits<Triangle_mesh>::vertex_descriptor as key type and a Boolean as value type. |
tmesh | a triangulated surface |
bhd | a halfedge descriptor on the boundary of mesh |
uvmap | an instantiation of the class VertexUVmap |
vimap | an instantiation of the class VertexIndexMap |
vpmap | an instantiation of the class VertexParameterizedMap |
iterations | an integer number of iterations to run the parameterization |
tmesh
must be a triangular mesh. vimap
must be initialized).