\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.13.1 - Planar Parameterization of Triangulated Surface Meshes
CGAL::Surface_mesh_parameterization::Orbifold_Tutte_parameterizer_3< SeamMesh, SolverTraits_ > Class Template Reference

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

Definition

template<typename SeamMesh, typename SolverTraits_ = Default>
class CGAL::Surface_mesh_parameterization::Orbifold_Tutte_parameterizer_3< SeamMesh, SolverTraits_ >

The class Orbifold_Tutte_parameterizer_3 implements Orbifold Tutte Planar Embeddings [1].

This is a borderless parameterization. A one-to-one mapping is guaranteed.

The main function of the class Orbifold_Tutte_parameterizer_3 is parameterize(), to which the user provides a Seam_mesh with marked edges (the seams) and a set of vertices of the mesh (the cones). The choice of cones influences the resulting parameterization, but not the choice of the seam path between these cones.

Some helper functions related to the class Orbifold_Tutte_parameterizer_3 (for example to read and compute paths between cones) can be found here .

The example orbifold.cpp shows how to select cones on the input mesh and automatically construct the seams and the cones on the Seam_mesh.

Is Model Of:
Parameterizer_3
Template Parameters
SeamMeshmust be a Seam_mesh, with underlying mesh any model of FaceListGraph and HalfedgeListGraph.
SolverTraits_must be a model of SparseLinearAlgebraTraits_d.
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::SparseLU<Eigen_sparse_matrix<double>::EigenType> >
Moreover, if SparseSuite solvers are available, which is greatly preferable for speed, then the default parameter is:
Eigen::UmfPackLU<Eigen_sparse_matrix<double>::EigenType> >
See also
Orbifold Helper Functions

Public Types

typedef SolverTraits_ Solver_traits
 

Public Member Functions

template<typename ConeMap , typename VertexIndexMap , typename VertexUVMap >
Error_code parameterize (SeamMesh &mesh, halfedge_descriptor bhd, const ConeMap &cmap, VertexUVMap uvmap, VertexIndexMap vimap) const
 Compute a one-to-one mapping from a triangular 3D surface mesh to a piece of the 2D space. More...
 
 Orbifold_Tutte_parameterizer_3 (const Orbifold_type orb_type=Square, const Weight_type weight_type=Cotangent)
 Constructor of the parameterizer. More...
 

Constructor & Destructor Documentation

◆ Orbifold_Tutte_parameterizer_3()

template<typename SeamMesh , typename SolverTraits_ = Default>
CGAL::Surface_mesh_parameterization::Orbifold_Tutte_parameterizer_3< SeamMesh, SolverTraits_ >::Orbifold_Tutte_parameterizer_3 ( const Orbifold_type  orb_type = Square,
const Weight_type  weight_type = Cotangent 
)

Constructor of the parameterizer.

The arguments allow to select the desired orbifold and weight types.

Member Function Documentation

◆ parameterize()

template<typename SeamMesh , typename SolverTraits_ = Default>
template<typename ConeMap , typename VertexIndexMap , typename VertexUVMap >
Error_code CGAL::Surface_mesh_parameterization::Orbifold_Tutte_parameterizer_3< SeamMesh, SolverTraits_ >::parameterize ( SeamMesh &  mesh,
halfedge_descriptor  bhd,
const ConeMap &  cmap,
VertexUVMap  uvmap,
VertexIndexMap  vimap 
) const

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

Template Parameters
ConeMapmust be a model of AssociativeContainer with key type boost::graph_traits<Seam_mesh>::vertex_descriptor and Cone_type as value type.
VertexUVmapmust be a model of ReadWritePropertyMap with boost::graph_traits<Seam_mesh>::vertex_descriptor as key type and Point_2 (type deduced from Seam_mesh using Kernel_traits) as value type.
VertexIndexMapmust be a model of ReadablePropertyMap with boost::graph_traits<Seam_mesh>::vertex_descriptor as key type and a unique integer as value type.
Parameters
mesha Seam_mesh parameterized by any model of a FaceListGraph and HalfedgeListGraph
bhda halfedge on the border of the seam mesh
cmapa mapping of the vertex_descriptors of mesh that are cones to their respective Cone_type classification.
uvmapan instanciation of the class VertexUVmap.
vimapan instanciation of the class VertexIndexMap.
Precondition
mesh must be a triangular mesh.
The underlying mesh of mesh is a topological ball.
The vertices must be indexed (vimap must be initialized).
The cones are vertices of mesh and their number is adapted to the orbifold type (4 for types I, II or III and 6 for type IV).
The seam edges form a set of segments that contains the different cones, and starts and ends at two different cones.
The seam edges form a set of segments that is homotopic to a line. Specifically, the paths between cones must not self intersect or intersect other paths (see Figure below).

orbifold_path.svg
Figure 66.1 Invalid (left) and valid (right) seam paths. Seam edges are drawn in teal. Cones are marked in yellow and blue.