CGAL 5.5 - Polygon Mesh Processing
Combinatorial Repairing

Functions to repair polygon soups and polygon meshes.

Functions

template<typename PolygonMesh >
bool CGAL::Polygon_mesh_processing::is_non_manifold_vertex (typename boost::graph_traits< PolygonMesh >::vertex_descriptor v, const PolygonMesh &pm)
 returns whether a vertex of a polygon mesh is non-manifold. More...
 
template<typename PolygonMesh , typename OutputIterator >
OutputIterator CGAL::Polygon_mesh_processing::non_manifold_vertices (const PolygonMesh &pm, OutputIterator out)
 collects the non-manifold vertices (if any) present in the mesh. More...
 
template<typename PolygonMesh , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::duplicate_non_manifold_vertices (PolygonMesh &pm, const NamedParameters &np=parameters::default_values())
 duplicates all the non-manifold vertices of the input mesh. More...
 
template<class PolygonMesh , class NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::merge_duplicated_vertices_in_boundary_cycle (typename boost::graph_traits< PolygonMesh >::halfedge_descriptor h, PolygonMesh &pm, const NamedParameters &np=parameters::default_values())
 merges identical vertices around a cycle of boundary edges. More...
 
template<class PolygonMesh , class NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::merge_duplicated_vertices_in_boundary_cycles (PolygonMesh &pm, const NamedParameters &np=parameters::default_values())
 extracts boundary cycles and merges the duplicated vertices of each cycle. More...
 
template<typename PolygonMesh , typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::polygon_mesh_to_polygon_soup (const PolygonMesh &mesh, PointRange &points, PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
 adds the vertices and faces of a mesh into a (possibly non-empty) polygon soup. More...
 
template<typename PolygonRange >
bool CGAL::Polygon_mesh_processing::is_polygon_soup_a_polygon_mesh (const PolygonRange &polygons)
 returns true if the soup of polygons defines a valid polygon mesh that can be handled by CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh(). More...
 
template<typename PolygonMesh , typename PointRange , typename PolygonRange , typename NamedParameters_PS = parameters::Default_named_parameters, typename NamedParameters_PM = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh (const PointRange &points, const PolygonRange &polygons, PolygonMesh &out, const NamedParameters_PS &np_ps=parameters::default_values(), const NamedParameters_PM &np_pm=parameters::default_values())
 builds a polygon mesh from a soup of polygons. More...
 
template<class PolygonMesh >
std::size_t CGAL::Polygon_mesh_processing::remove_isolated_vertices (PolygonMesh &pmesh)
 removes the isolated vertices from any polygon mesh. More...
 
template<typename TriangleMesh , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::remove_connected_components_of_negligible_size (TriangleMesh &tmesh, const NamedParameters &np=parameters::default_values())
 removes connected components whose area or volume is under a certain threshold value. More...
 
template<typename PointRange , typename PolygonRange >
std::size_t CGAL::Polygon_mesh_processing::remove_isolated_points_in_polygon_soup (PointRange &points, PolygonRange &polygons)
 removes the isolated points from a polygon soup. More...
 
template<typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::merge_duplicate_points_in_polygon_soup (PointRange &points, PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
 merges the duplicate points in a polygon soup. More...
 
template<typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::merge_duplicate_polygons_in_polygon_soup (const PointRange &points, PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
 merges the duplicate polygons in a polygon soup. More...
 
template<typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::repair_polygon_soup (PointRange &points, PolygonRange &polygons, const NamedParameters &np=parameters::default_values())
 cleans a given polygon soup through various repairing operations. More...
 
template<typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_boundary_cycle (const typename boost::graph_traits< PolygonMesh >::halfedge_descriptor h, PolygonMesh &pmesh, const NamedParameters &np=parameters::default_values())
 stitches together, whenever possible, two halfedges belonging to the boundary cycle described by the halfedge h. More...
 
template<typename BorderHalfedgeRange , typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_boundary_cycles (const BorderHalfedgeRange &boundary_cycle_representatives, PolygonMesh &pmesh, const NamedParameters &np=parameters::default_values())
 stitches together, whenever possible, two halfedges belonging to the same boundary cycle. More...
 
template<typename PolygonMesh , typename HalfedgePairsRange , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_borders (PolygonMesh &pmesh, const HalfedgePairsRange &hedge_pairs_to_stitch, const NamedParameters &np=parameters::default_values(), typename boost::enable_if< typename boost::has_range_iterator< HalfedgePairsRange > >::type *=0)
 stitches together border halfedges in a polygon mesh. More...
 
template<typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_borders (PolygonMesh &pmesh, const NamedParameters &np=parameters::default_values())
 Same as the other overload, but the pairs of halfedges to be stitched are automatically found amongst all border halfedges. More...
 
template<typename BorderHalfedgeRange , typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_borders (const BorderHalfedgeRange &boundary_cycle_representatives, PolygonMesh &pmesh, const NamedParameters &np=parameters::default_values())
 Same as the other overload, but the pairs of halfedges to be stitched are automatically found amongst halfedges in cycles described by boundary_cycle_representatives. More...
 

Function Documentation

◆ duplicate_non_manifold_vertices()

template<typename PolygonMesh , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::duplicate_non_manifold_vertices ( PolygonMesh &  pm,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/manifoldness.h>

duplicates all the non-manifold vertices of the input mesh.

Template Parameters
PolygonMesha model of HalfedgeListGraph and MutableHalfedgeGraph
NamedParametersa sequence of Named Parameters
Parameters
pmthe surface mesh to be repaired
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pm
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pm)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.

  • a property map containing the constrained-or-not status of each vertex of pm.
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and bool as value type. It must be default constructible.
  • Default: a default property map where no vertex is constrained
  • Extra: put(vcm, v, true) will be called for each duplicated vertices, as well as the original non-manifold vertex in the input mesh.

  • an output iterator to collect the duplicated vertices
  • Type: a model of OutputIterator with value type std::vector<vertex_descriptor>
  • Default: unused
  • Extra: The first vertex of each vector is a non-manifold vertex of the input mesh, followed by the new vertices that were created to fix the given non-manifold configuration.
Returns
the number of vertices created
See also
non_manifold_vertices()
Examples:
Polygon_mesh_processing/manifoldness_repair_example.cpp.

◆ is_non_manifold_vertex()

template<typename PolygonMesh >
bool CGAL::Polygon_mesh_processing::is_non_manifold_vertex ( typename boost::graph_traits< PolygonMesh >::vertex_descriptor  v,
const PolygonMesh &  pm 
)

#include <CGAL/Polygon_mesh_processing/manifoldness.h>

returns whether a vertex of a polygon mesh is non-manifold.

Warning
This function has linear runtime with respect to the size of the mesh. The function non_manifold_vertices() should be used when gathering all non manifold vertices.
Template Parameters
PolygonMesha model of HalfedgeListGraph
Parameters
va vertex of pm
pma triangle mesh containing v
Returns
true if the vertex is non-manifold, false otherwise
See also
duplicate_non_manifold_vertices()
Examples:
Polygon_mesh_processing/manifoldness_repair_example.cpp.

◆ is_polygon_soup_a_polygon_mesh()

template<typename PolygonRange >
bool CGAL::Polygon_mesh_processing::is_polygon_soup_a_polygon_mesh ( const PolygonRange &  polygons)

#include <CGAL/Polygon_mesh_processing/polygon_soup_to_polygon_mesh.h>

returns true if the soup of polygons defines a valid polygon mesh that can be handled by CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh().

It checks that each edge has at most two incident faces and such an edge is visited in opposite direction along the two face boundaries, no polygon has twice the same vertex, and the polygon soup describes a manifold surface. This function does not require a range of points as an argument since the check is purely topological. To each vertex of the mesh is associated an index that is used in the description of the boundaries of the polygons provided in polygons.

Template Parameters
PolygonRangea model of the concept RandomAccessContainer whose value_type is a model of the concept RandomAccessContainer whose value_type is std::size_t.
Parameters
polygonseach element in the range describes a polygon using the indices of the vertices.
See also
CGAL::Polygon_mesh_processing::orient_polygon_soup()
Examples:
Polygon_mesh_processing/orientation_pipeline_example.cpp.

◆ merge_duplicate_points_in_polygon_soup()

template<typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::merge_duplicate_points_in_polygon_soup ( PointRange &  points,
PolygonRange &  polygons,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/repair_polygon_soup.h>

merges the duplicate points in a polygon soup.

Note that the index of a point that is merged with another point will thus change in all the polygons that the point appears in.

Template Parameters
PointRangea model of the concepts SequenceContainer and Swappable whose value type is the point type.
PolygonRangea model of the concept RandomAccessContainer whose value_type is itself a model of the concept RandomAccessContainer whose value_type is std::size_t.
NamedParametersa sequence of Named Parameters
Parameters
pointspoints of the soup of polygons
polygonsa vector of polygons. Each element in the vector describes a polygon using the indices of the points in points.
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • an instance of a geometric traits class
  • Type: The traits class must provide the nested functor Less_xyz_3 to compare lexicographically two points a function Less_xyz_3 less_xyz_3_object().
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
Returns
the number of removed points
See also
repair_polygon_soup()

◆ merge_duplicate_polygons_in_polygon_soup()

template<typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::merge_duplicate_polygons_in_polygon_soup ( const PointRange &  points,
PolygonRange &  polygons,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/repair_polygon_soup.h>

merges the duplicate polygons in a polygon soup.

Two polygons are duplicate if they share the same vertices in the same order. Note that the first vertex of the polygon does not matter, that is the triangle 0,1,2 is a duplicate of the triangle 2,0,1.

Template Parameters
PointRangea model of the concept RandomAccessContainer whose value type is the point type.
PolygonRangea model of the concept SequenceContainer whose value_type is itself a model of the concepts RandomAccessContainer and ReversibleContainer whose value_type is std::size_t.
NamedParametersa sequence of Named Parameters
Parameters
pointspoints of the soup of polygons
polygonsa vector of polygons. Each element in the vector describes a polygon using the indices of the points in points.
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • an instance of a geometric traits class
  • Type: The traits class must provide the nested functor Less_xyz_3 to compare lexicographically two points a function Less_xyz_3 less_xyz_3_object().
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.

  • Parameter to indicate, when multiple polygons are duplicates, whether all the duplicate polygons should be removed or if one (arbitrarily chosen) face should be kept.
  • Type: Boolean
  • Default: false

  • Parameter to indicate if polygon orientation should be taken into account when determining whether two polygons are duplicates, that is, whether e.g. the triangles 0,1,2 and 0,2,1 are duplicates.
  • Type: Boolean
  • Default: false
Returns
the number of removed polygons
See also
repair_polygon_soup()

◆ merge_duplicated_vertices_in_boundary_cycle()

template<class PolygonMesh , class NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::merge_duplicated_vertices_in_boundary_cycle ( typename boost::graph_traits< PolygonMesh >::halfedge_descriptor  h,
PolygonMesh &  pm,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/merge_border_vertices.h>

merges identical vertices around a cycle of boundary edges.

Template Parameters
PolygonMesha model of FaceListGraph and MutableFaceGraph.
NamedParametersa sequence of Named Parameters.
Parameters
ha halfedge that belongs to a boundary cycle.
pmthe polygon mesh which contains the boundary cycle.
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pm
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pm)
See also
merge_duplicated_vertices_in_boundary_cycles()

◆ merge_duplicated_vertices_in_boundary_cycles()

template<class PolygonMesh , class NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::merge_duplicated_vertices_in_boundary_cycles ( PolygonMesh &  pm,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/merge_border_vertices.h>

extracts boundary cycles and merges the duplicated vertices of each cycle.

Template Parameters
PolygonMesha model of FaceListGraph and MutableFaceGraph.
NamedParametersa sequence of Named Parameters.
Parameters
pmthe polygon mesh which contains the cycles.
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pm
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pm)
See also
merge_duplicated_vertices_in_boundary_cycle()

◆ non_manifold_vertices()

template<typename PolygonMesh , typename OutputIterator >
OutputIterator CGAL::Polygon_mesh_processing::non_manifold_vertices ( const PolygonMesh &  pm,
OutputIterator  out 
)

#include <CGAL/Polygon_mesh_processing/manifoldness.h>

collects the non-manifold vertices (if any) present in the mesh.

A non-manifold vertex v is returned via one incident halfedge h such that target(h, pm) = v for all the umbrellas that v appears in (an umbrella being the set of faces incident to all the halfedges reachable by walking around v using hnext = prev(opposite(h, pm), pm), starting from h).

Template Parameters
PolygonMesha model of HalfedgeListGraph
OutputIteratora model of OutputIterator holding objects of type boost::graph_traits<PolygonMesh>::halfedge_descriptor
Parameters
pma triangle mesh
outthe output iterator that collects halfedges incident to v
Returns
the output iterator
See also
is_non_manifold_vertex()
duplicate_non_manifold_vertices()

◆ polygon_mesh_to_polygon_soup()

template<typename PolygonMesh , typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::polygon_mesh_to_polygon_soup ( const PolygonMesh &  mesh,
PointRange &  points,
PolygonRange &  polygons,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/polygon_mesh_to_polygon_soup.h>

adds the vertices and faces of a mesh into a (possibly non-empty) polygon soup.

Template Parameters
PolygonMesha model of FaceListGraph
PointRangea model of the concepts RandomAccessContainer and BackInsertionSequence whose value type can be constructed from the point type of the polygon mesh
PolygonRangea model of the concepts RandomAccessContainer and BackInsertionSequence whose value type is itself a model of the concepts RandomAccessContainer and BackInsertionSequence whose value type is std::size_t
NamedParametersa sequence of Named Parameters
Parameters
meshthe mesh whose faces are being put in the polygon soup
pointspoints making the polygons of the soup
polygonseach element in the vector describes a polygon using the indices of the points in points
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of mesh
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, mesh)

Advanced

PolygonRange can also be a model of the concepts RandomAccessContainer and BackInsertionSequence whose value type is an array, but it is the user's responsability to ensure that all faces have the same number of vertices, and that this number is equal to the size of the array.

See also
CGAL::Polygon_mesh_processing::orient_polygon_soup()
CGAL::Polygon_mesh_processing::is_polygon_soup_a_polygon_mesh()
CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh()

◆ polygon_soup_to_polygon_mesh()

template<typename PolygonMesh , typename PointRange , typename PolygonRange , typename NamedParameters_PS = parameters::Default_named_parameters, typename NamedParameters_PM = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh ( const PointRange &  points,
const PolygonRange &  polygons,
PolygonMesh &  out,
const NamedParameters_PS &  np_ps = parameters::default_values(),
const NamedParameters_PM &  np_pm = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/polygon_soup_to_polygon_mesh.h>

builds a polygon mesh from a soup of polygons.

Precondition
the input polygon soup describes a consistently oriented polygon mesh. This can be checked using the function CGAL::Polygon_mesh_processing::is_polygon_soup_a_polygon_mesh(polygons) .
Template Parameters
PolygonMesha model of MutableFaceGraph
PointRangea model of the concept RandomAccessContainer whose value type is the point type
PolygonRangea model of the concept RandomAccessContainer whose value type is a model of the concept RandomAccessContainer whose value type is std::size_t
NamedParameters_PSa sequence of Named Parameters
NamedParameters_PMa sequence of Named Parameters
Parameters
pointspoints of the soup of polygons
polygonseach element in the range describes a polygon using the indices of the points in points
outthe polygon mesh to be built
np_psan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the elements of the range points
  • Type: a model of ReadablePropertyMap whose value type is a point type convertible to the point type of the vertex point map associated to the polygon mesh
  • Default: CGAL::Identity_property_map
Parameters
np_pman optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of out
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, out)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.
See also
CGAL::Polygon_mesh_processing::orient_polygon_soup()
CGAL::Polygon_mesh_processing::is_polygon_soup_a_polygon_mesh()
CGAL::Polygon_mesh_processing::polygon_mesh_to_polygon_soup()
Examples:
Polygon_mesh_processing/cc_compatible_orientations.cpp, Polygon_mesh_processing/orient_polygon_soup_example.cpp, Polygon_mesh_processing/orientation_pipeline_example.cpp, and Polygon_mesh_processing/repair_polygon_soup_example.cpp.

◆ remove_connected_components_of_negligible_size()

template<typename TriangleMesh , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::remove_connected_components_of_negligible_size ( TriangleMesh &  tmesh,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/repair.h>

removes connected components whose area or volume is under a certain threshold value.

Thresholds are provided via Named Parameters. (see below). If thresholds are not provided by the user, default values are computed as follows:

  • the area threshold is taken as the square of one percent of the length of the diagonal of the bounding box of the mesh.
  • the volume threshold is taken as the third power of one percent of the length of the diagonal of the bounding box of the mesh.

The area and volume of a connected component will always be positive values (regardless of the orientation of the mesh).

As a consequence of the last sentence, the area or volume criteria can be disabled by passing zero (0) as threshold value.

Template Parameters
TriangleMesha model of FaceListGraph and MutableFaceGraph
NamedParametersa sequence of Named Parameters
Parameters
tmeshthe triangulated polygon mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tmesh
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in TriangleMesh.

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
  • Extra: Exact constructions kernels are not supported by this function.

  • a property map associating to each face of tmesh a unique index between 0 and num_faces(tmesh) - 1
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::face_descriptor as key type and std::size_t as value type
  • Default: an automatically indexed internal map

  • a fixed value such that only connected components whose area is larger than this value are kept
  • Type: geom_traits::FT
  • Default: 1% of the length of the diagonal of the axis-aligned bounding box of the mesh, squared

  • a fixed value such that only connected components whose volume is larger than this value are kept (only applies to closed connected components)
  • Type: geom_traits::FT
  • Default: 1% of the length of the diagonal of the axis-aligned bounding box of the mesh, cubed
  • Extra: The mesh must be closed.

  • a property map containing the constrained-or-not status of each edge of tmesh
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::edge_descriptor as key type and bool as value type. It must be default constructible.
  • Default: a default property map where no edge is constrained
  • Extra: A constrained edge can be split or collapsed, but not flipped, nor its endpoints moved by smoothing.

  • If true, the mesh will not be altered, but the number of components that would be removed is returned.
  • Type: Boolean
  • Default: false

  • An output iterator to collect the faces that would be removed by the algorithm, when using the "dry run" mode (see parameter dry_run)
  • Type: a model of OutputIterator with value type face_descriptor
  • Default: unused
Returns
the number of connected components removed (ignoring isolated vertices).
See also
keep_connected_components()
remove_connected_components()

◆ remove_isolated_points_in_polygon_soup()

template<typename PointRange , typename PolygonRange >
std::size_t CGAL::Polygon_mesh_processing::remove_isolated_points_in_polygon_soup ( PointRange &  points,
PolygonRange &  polygons 
)

#include <CGAL/Polygon_mesh_processing/repair_polygon_soup.h>

removes the isolated points from a polygon soup.

A point is considered isolated if it does not appear in any polygon of the soup.

Template Parameters
PointRangea model of the concept SequenceContainer whose value type is the point type.
PolygonRangea model of the concept RandomAccessContainer whose value_type is itself a model of the concept RandomAccessContainer whose value_type is std::size_t.
Parameters
pointspoints of the soup of polygons
polygonsa vector of polygons. Each element in the vector describes a polygon using the indices of the points in points.
Returns
the number of removed isolated points
See also
repair_polygon_soup()

◆ remove_isolated_vertices()

template<class PolygonMesh >
std::size_t CGAL::Polygon_mesh_processing::remove_isolated_vertices ( PolygonMesh &  pmesh)

#include <CGAL/Polygon_mesh_processing/repair.h>

removes the isolated vertices from any polygon mesh.

A vertex is considered isolated if it is not incident to a simplex of higher dimension.

Template Parameters
PolygonMesha model of FaceListGraph and MutableFaceGraph
Parameters
pmeshthe polygon mesh to be repaired
Returns
the number of removed isolated vertices

◆ repair_polygon_soup()

template<typename PointRange , typename PolygonRange , typename NamedParameters = parameters::Default_named_parameters>
void CGAL::Polygon_mesh_processing::repair_polygon_soup ( PointRange &  points,
PolygonRange &  polygons,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/repair_polygon_soup.h>

cleans a given polygon soup through various repairing operations.

More precisely, this function carries out the following tasks, in the same order as they are listed:

Note that the point and polygon containers will be modified by the repairing operations, and thus the indexing of the polygons will also be changed.

Template Parameters
PointRangea model of the concepts SequenceContainer and Swappable and whose value type is the point type.
PolygonRangea model of the concept SequenceContainer. whose value_type is itself a model of the concepts SequenceContainer, Swappable, and ReversibleContainer whose value_type is std::size_t.
NamedParametersa sequence of Named Parameters
Parameters
pointspoints of the soup of polygons
polygonsa vector of polygons. Each element in the vector describes a polygon using the indices of the points in points.
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • an instance of a geometric traits class
  • Type: The traits class must provide the nested functors Less_xyz_3 and Equal_3 to respectivelycompare lexicographically two points and to check if two points are identical. For each functor Foo, a function Foo foo_object() must be provided.
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.

  • Parameter to indicate, when multiple polygons are duplicates, whether all the duplicate polygons should be removed or if one (arbitrarily chosen) face should be kept.
  • Type: Boolean
  • Default: false

  • Parameter to indicate if polygon orientation should be taken into account when determining whether two polygons are duplicates, that is, whether e.g. the triangles 0,1,2 and 0,2,1 are duplicates.
  • Type: Boolean
  • Default: false
Examples:
Polygon_mesh_processing/repair_polygon_soup_example.cpp.

◆ stitch_borders() [1/3]

template<typename PolygonMesh , typename HalfedgePairsRange , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_borders ( PolygonMesh &  pmesh,
const HalfedgePairsRange &  hedge_pairs_to_stitch,
const NamedParameters &  np = parameters::default_values(),
typename boost::enable_if< typename boost::has_range_iterator< HalfedgePairsRange > >::type *  = 0 
)

#include <CGAL/Polygon_mesh_processing/stitch_borders.h>

stitches together border halfedges in a polygon mesh.

The halfedges to be stitched are provided in hedge_pairs_to_stitch. For each pair p in this vector, p.second and its opposite will be removed from pmesh.

Template Parameters
PolygonMesha model of MutableFaceGraph
HalfedgePairsRangea range of std::pair<boost::graph_traits<PolygonMesh>::halfedge_descriptor, boost::graph_traits<PolygonMesh>::halfedge_descriptor>, model of Range. Its iterator type is InputIterator.
Parameters
pmeshthe polygon mesh to be modified by stitching
hedge_pairs_to_stitcha range of std::pair of halfedges to be stitched together
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pm
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.
Returns
the number of pairs of halfedges that were stitched.
See also
stitch_boundary_cycle()
stitch_boundary_cycles()
Examples:
Polygon_mesh_processing/cc_compatible_orientations.cpp, and Polygon_mesh_processing/stitch_borders_example.cpp.

◆ stitch_borders() [2/3]

template<typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_borders ( PolygonMesh &  pmesh,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/stitch_borders.h>

Same as the other overload, but the pairs of halfedges to be stitched are automatically found amongst all border halfedges.

Two border halfedges h1 and h2 are set to be stitched if the points associated to the source and target vertices of h1 are the same as those of the target and source vertices of h2, respectively.

Template Parameters
BorderHalfedgeRangea model of Range with value type boost::graph_traits<PolygonMesh>::halfedge_descriptor
PolygonMesha model of FaceListGraph and MutableFaceGraph
NamedParametersa sequence of Named Parameters
Parameters
pmeshthe polygon mesh to be modified by the stitching procedure
npoptional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pmesh
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.

  • specifies if the borders should only be stitched only within their own connected component.
  • Type: Boolean
  • Default: false

  • a property map associating to each face of pmesh a unique index between 0 and num_faces(pmesh) - 1
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<PolygonMesh>::face_descriptor as key type and std::size_t as value type
  • Default: an automatically indexed internal map
Returns
the number of pairs of halfedges that were stitched.
See also
stitch_boundary_cycle()
stitch_boundary_cycles()

◆ stitch_borders() [3/3]

template<typename BorderHalfedgeRange , typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_borders ( const BorderHalfedgeRange &  boundary_cycle_representatives,
PolygonMesh &  pmesh,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/stitch_borders.h>

Same as the other overload, but the pairs of halfedges to be stitched are automatically found amongst halfedges in cycles described by boundary_cycle_representatives.

Two border halfedges h1 and h2 are set to be stitched if the points associated to the source and target vertices of h1 are the same as those of the target and source vertices of h2, respectively.

Template Parameters
BorderHalfedgeRangea model of Range with value type boost::graph_traits<PolygonMesh>::halfedge_descriptor
PolygonMesha model of FaceListGraph and MutableFaceGraph
NamedParametersa sequence of Named Parameters
Parameters
boundary_cycle_representativesa range of border halfedges, each describing a boundary cycle whose halfedges will be considered for stitching
pmeshthe polygon mesh to be modified by the stitching procedure
npoptional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • specifies if the borders should only be stitched only within their own connected component.
  • Type: Boolean
  • Default: false

  • a property map associating to each face of pmesh a unique index between 0 and num_faces(pmesh) - 1
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<PolygonMesh>::face_descriptor as key type and std::size_t as value type
  • Default: an automatically indexed internal map

  • a property map associating points to the vertices of pmesh
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.

  • an instance of a geometric traits class
  • Type: The traits class must provide the nested type Point_3, and the nested functors:
    • Less_xyz_3 to compare lexicographically two points
    • Equal_3 to check whether two points are identical. For each functor Foo, a function Foo foo_object() must be provided.
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
Returns
the number of pairs of halfedges that were stitched.
See also
stitch_boundary_cycle()
stitch_boundary_cycles()

◆ stitch_boundary_cycle()

template<typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_boundary_cycle ( const typename boost::graph_traits< PolygonMesh >::halfedge_descriptor  h,
PolygonMesh &  pmesh,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/stitch_borders.h>

stitches together, whenever possible, two halfedges belonging to the boundary cycle described by the halfedge h.

Two border halfedges h1 and h2 can be stitched if the points associated to the source and target vertices of h1 are the same as those of the target and source vertices of h2, respectively.

Template Parameters
PolygonMesha model of MutableFaceGraph
NamedParametersa sequence of Named Parameters
Parameters
ha border halfedge of the polygon mesh pmesh
pmeshthe polygon mesh to be stitched
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pm
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.

  • an instance of a geometric traits class
  • Type: The traits class must provide the nested type Point_3, and the nested functors:
    • Less_xyz_3 to compare lexicographically two points
    • Equal_3 to check whether two points are identical. For each functor Foo, a function Foo foo_object() must be provided.
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
Returns
the number of pairs of halfedges that were stitched.
See also
stitch_boundary_cycles()
stitch_borders()

◆ stitch_boundary_cycles()

template<typename BorderHalfedgeRange , typename PolygonMesh , typename NamedParameters = CGAL::parameters::Default_named_parameters>
std::size_t CGAL::Polygon_mesh_processing::stitch_boundary_cycles ( const BorderHalfedgeRange &  boundary_cycle_representatives,
PolygonMesh &  pmesh,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Polygon_mesh_processing/stitch_borders.h>

stitches together, whenever possible, two halfedges belonging to the same boundary cycle.

Two border halfedges h1 and h2 can be stitched if the points associated to the source and target vertices of h1 are the same as those of the target and source vertices of h2, respectively.

Template Parameters
BorderHalfedgeRangea model of Range with value type boost::graph_traits<PolygonMesh>::halfedge_descriptor
PolygonMesha model of MutableFaceGraph
NamedParametersa sequence of Named Parameters
Parameters
boundary_cycle_representativesa range of border halfedges, each describing a boundary cycle of the mesh pmesh
pmeshthe polygon mesh to be modified by stitching
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of pm
  • Type: a class model of ReadWritePropertyMap with boost::graph_traits<PolygonMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, pmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in PolygonMesh.

  • an instance of a geometric traits class
  • Type: The traits class must provide the nested type Point_3, and the nested functors:
    • Less_xyz_3 to compare lexicographically two points
    • Equal_3 to check whether two points are identical. For each functor Foo, a function Foo foo_object() must be provided.
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
Returns
the number of pairs of halfedges that were stitched.
See also
stitch_boundary_cycle()
stitch_borders()