CGAL 5.4.5 - Planar Parameterization of Triangulated Surface Meshes

Functions

template<typename TriangleMesh , typename VertexIndexMap , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones (const TriangleMesh &tm, std::ifstream &in, VertexIndexMap vpmap, ConeOutputIterator out)
 reads a serie of cones from an input stream. More...
 
template<typename TriangleMesh , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones (const TriangleMesh &tm, std::ifstream &in, ConeOutputIterator out)
 Same as above, using the default indexation of the vertices of tm: vertices are numbered from 0 to num_vertices(tm)-1, in the order that they appear while calling vertices(tm). More...
 
template<typename TriangleMesh , typename VertexIndexMap , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones (const TriangleMesh &tm, const char *filename, VertexIndexMap vpmap, ConeOutputIterator out)
 Same as above, but from a file instead of a stream. More...
 
template<typename TriangleMesh , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones (const TriangleMesh &tm, const char *filename, ConeOutputIterator out)
 Same as above, but from a file instead of a stream. More...
 
template<typename SeamMesh , typename ConeInputBidirectionalIterator , typename ConeMap >
bool CGAL::Surface_mesh_parameterization::locate_cones (const SeamMesh &mesh, ConeInputBidirectionalIterator first, ConeInputBidirectionalIterator beyond, ConeMap &cones)
 locates the cones on the seam mesh (that is, find the corresponding seam mesh vertex_descriptor) and mark them with a tag to indicate whether the cone is a simple cone or a duplicated cone (see Cone_type ). More...
 
template<typename SeamMesh , typename ConeInputBidirectionalIterator , typename ConeMap >
bool CGAL::Surface_mesh_parameterization::locate_unordered_cones (const SeamMesh &mesh, ConeInputBidirectionalIterator first, ConeInputBidirectionalIterator beyond, ConeMap &cones)
 Same as above, but the cones are not ordered and we thus use seam mesh information to determine which cones are extremities of the seam (so-called unique cones) or not (so-called duplicate cones). More...
 
template<typename TriangleMesh , typename EdgeOutputIterator >
void CGAL::Surface_mesh_parameterization::compute_shortest_paths_between_two_cones (const TriangleMesh &mesh, typename boost::graph_traits< TriangleMesh >::vertex_descriptor source, typename boost::graph_traits< TriangleMesh >::vertex_descriptor target, EdgeOutputIterator oi)
 computes the shortest path between source and target over mesh, using boost::dijkstra_shortest_paths(). More...
 
template<typename TriangleMesh , typename InputConesForwardIterator , typename SeamContainer >
void CGAL::Surface_mesh_parameterization::compute_shortest_paths_between_cones (const TriangleMesh &mesh, InputConesForwardIterator first, InputConesForwardIterator beyond, SeamContainer &seams)
 Given a range [first; beyond[ of cones (described as vertex descriptors), compute the shortest path for all pairs of consecutive entries in the range and add them to the container seams. More...
 

Function Documentation

◆ compute_shortest_paths_between_cones()

template<typename TriangleMesh , typename InputConesForwardIterator , typename SeamContainer >
void CGAL::Surface_mesh_parameterization::compute_shortest_paths_between_cones ( const TriangleMesh &  mesh,
InputConesForwardIterator  first,
InputConesForwardIterator  beyond,
SeamContainer &  seams 
)

#include <CGAL/Surface_mesh_parameterization/orbifold_shortest_path.h>

Given a range [first; beyond[ of cones (described as vertex descriptors), compute the shortest path for all pairs of consecutive entries in the range and add them to the container seams.

Template Parameters
TriangleMeshA triangle mesh, model of FaceListGraph and HalfedgeListGraph.
InputConesForwardIteratorA model of ForwardIterator with value type boost::graph_traits<TriangleMesh>::vertex_descriptor.
SeamContainerA model of SequenceContainer with value type boost::graph_traits<TriangleMesh>::edge_descriptor.
Parameters
meshthe triangular mesh on which paths are computed
first,beyonda range of cones
seamsa container that will store the paths, as a sequence of edges of the mesh.
Precondition
std::distance(first,beyond) > 1
Examples:
Surface_mesh_parameterization/orbifold.cpp.

◆ compute_shortest_paths_between_two_cones()

template<typename TriangleMesh , typename EdgeOutputIterator >
void CGAL::Surface_mesh_parameterization::compute_shortest_paths_between_two_cones ( const TriangleMesh &  mesh,
typename boost::graph_traits< TriangleMesh >::vertex_descriptor  source,
typename boost::graph_traits< TriangleMesh >::vertex_descriptor  target,
EdgeOutputIterator  oi 
)

#include <CGAL/Surface_mesh_parameterization/orbifold_shortest_path.h>

computes the shortest path between source and target over mesh, using boost::dijkstra_shortest_paths().

Template Parameters
TriangleMeshA triangle mesh, model of FaceListGraph and HalfedgeListGraph.
EdgeOutputIteratorA model of OutputIterator with value type boost::graph_traits<TriangleMesh>::edge_descriptor.
Parameters
meshthe triangular mesh to be parameterized
source,targetthe extremities of the path to be computed
oithe output iterator
Precondition
source and target are vertices of mesh.
source != target

◆ locate_cones()

template<typename SeamMesh , typename ConeInputBidirectionalIterator , typename ConeMap >
bool CGAL::Surface_mesh_parameterization::locate_cones ( const SeamMesh &  mesh,
ConeInputBidirectionalIterator  first,
ConeInputBidirectionalIterator  beyond,
ConeMap &  cones 
)

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

locates the cones on the seam mesh (that is, find the corresponding seam mesh vertex_descriptor) and mark them with a tag to indicate whether the cone is a simple cone or a duplicated cone (see Cone_type ).

Attention
The cones must be ordered: the first and last cones are the extremities of the seam.
Template Parameters
SeamMeshis the same mesh that is passed to the parameterizer. It is an object of typeCGAL::Seam_mesh, but is passed here as a template parameter for convenience, to avoid having to pass the multiple template parameters of the class CGAL::Seam_mesh.
ConeInputBidirectionalIteratormust be a model of BidirectionalIterator with value type boost::graph_traits<SeamMesh::Triangle_mesh>::vertex_descriptor.
ConeMapmust be a model of AssociativeContainer with boost::graph_traits<SeamMesh>::vertex_descriptor as key type and Cone_type as value type.
Parameters
meshthe seam mesh
first,beyondthe range of cones, as vertex descriptors of the base mesh.
conesan object of type ConeMap. Cones will be stored in this container as vertex descriptors of the seam mesh, along with their associated cone types.
Examples:
Surface_mesh_parameterization/orbifold.cpp.

◆ locate_unordered_cones()

template<typename SeamMesh , typename ConeInputBidirectionalIterator , typename ConeMap >
bool CGAL::Surface_mesh_parameterization::locate_unordered_cones ( const SeamMesh &  mesh,
ConeInputBidirectionalIterator  first,
ConeInputBidirectionalIterator  beyond,
ConeMap &  cones 
)

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

Same as above, but the cones are not ordered and we thus use seam mesh information to determine which cones are extremities of the seam (so-called unique cones) or not (so-called duplicate cones).

◆ read_cones() [1/4]

template<typename TriangleMesh , typename VertexIndexMap , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones ( const TriangleMesh &  tm,
std::ifstream &  in,
VertexIndexMap  vpmap,
ConeOutputIterator  out 
)

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

reads a serie of cones from an input stream.

Cones are passed as an integer value that is the index of a vertex handle in the mesh tm, using the vertex index property mapvpmap` for correspondency.

Attention
The mesh is here tm, it is the base mesh of the CGAL::Seam_mesh that is passed in input, not the seam mesh itself.
Template Parameters
TriangleMeshA triangle mesh, model of FaceListGraph and HalfedgeListGraph.
VertexIndexMapmust be a model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and a unique integer as value type.
ConeOutputIteratora model of OutputIterator with value type boost::graph_traits<TriangleMesh>::vertex_descriptor.
Parameters
tmthe triangular mesh to be parameterized
inthe input stream
vpmapan initialized vertex index map
outthe output iterator
Precondition
The number of cones must match the chosen Orbifold_type .
No two cones correspond to the same vertex (all cones have different index).
Returns
The corresponding vertex descriptors are output, in the same order as the input integers, in out. The function checks if the input is valid (no duplicate, correct number of cones) and returns an Error_code.

◆ read_cones() [2/4]

template<typename TriangleMesh , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones ( const TriangleMesh &  tm,
std::ifstream &  in,
ConeOutputIterator  out 
)

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

Same as above, using the default indexation of the vertices of tm: vertices are numbered from 0 to num_vertices(tm)-1, in the order that they appear while calling vertices(tm).

◆ read_cones() [3/4]

template<typename TriangleMesh , typename VertexIndexMap , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones ( const TriangleMesh &  tm,
const char *  filename,
VertexIndexMap  vpmap,
ConeOutputIterator  out 
)

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

Same as above, but from a file instead of a stream.

◆ read_cones() [4/4]

template<typename TriangleMesh , typename ConeOutputIterator >
Error_code CGAL::Surface_mesh_parameterization::read_cones ( const TriangleMesh &  tm,
const char *  filename,
ConeOutputIterator  out 
)

#include <CGAL/Surface_mesh_parameterization/Orbifold_Tutte_parameterizer_3.h>

Same as above, but from a file instead of a stream.

The default indexation of the vertices of tm is used: vertices are numbered from 0 to num_vertices(tm)-1, in the order that they appear while calling vertices(tm).