Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions BGL/doc/BGL/Doxyfile.in
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@ INPUT += ${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/IO/polygon_mesh_io.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/Euler_operations.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/iterator.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/helpers.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/border.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/generators.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/selection.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/boost/graph/split_graph_into_polylines.h \
Expand Down
3 changes: 2 additions & 1 deletion BGL/doc/BGL/PackageDescription.txt
Original file line number Diff line number Diff line change
Expand Up @@ -673,13 +673,14 @@ user might encounter.
\cgalCRPSection{Helper Functions}
- `CGAL::is_border()`
- `CGAL::is_border_edge()`
- `CGAL::border_halfedges()`
- `CGAL::extract_boundary_cycles()`
- `CGAL::is_bivalent()`
- `CGAL::is_bivalent_mesh()`
- `CGAL::is_trivalent()`
- `CGAL::is_trivalent_mesh()`
- `CGAL::is_isolated_triangle()`
- `CGAL::is_closed()`

- `CGAL::is_triangle()`
- `CGAL::is_triangle_mesh()`
- `CGAL::is_quad()`
Expand Down
267 changes: 267 additions & 0 deletions BGL/include/CGAL/boost/graph/border.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,267 @@
// Copyright (c) 2015 GeometryFactory (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL$
// $Id$
// SPDX-License-Identifier: LGPL-3.0-or-later OR LicenseRef-Commercial
//
//
// Author(s) : Jane Tournois

#ifndef CGAL_BGL_BORDER_H
#define CGAL_BGL_BORDER_H

#include <CGAL/algorithm.h>
#include <CGAL/boost/graph/iterator.h>
#include <CGAL/boost/graph/helpers.h>
#include <CGAL/Named_function_parameters.h>
#include <CGAL/boost/graph/named_params_helper.h>

#include <boost/graph/graph_traits.hpp>

#include <map>
#include <unordered_set>
#include <utility>

namespace CGAL {
namespace internal {

template <typename Graph>
std::size_t border_size(typename boost::graph_traits<Graph>::halfedge_descriptor h,
const Graph& g)
{
// if you want to use it on a non-border halfedge, just use `degree(face, mesh)`
CGAL_precondition(is_border(h, g));

std::size_t res = 0;
typename boost::graph_traits<Graph>::halfedge_descriptor done = h;
do
{
++res;
h = next(h, g);
}
while(h != done);

return res;
}

template<typename Graph,
typename FaceRange,
typename HalfedgeOutputIterator>
HalfedgeOutputIterator border_halfedges_impl(const FaceRange& face_range,
HalfedgeOutputIterator out,
const Graph& g)
{
typedef typename boost::graph_traits<Graph>::halfedge_descriptor halfedge_descriptor;
typedef typename boost::graph_traits<Graph>::face_descriptor face_descriptor;

// collect halfedges that appear only once
// the bool is true if the halfedge stored is the one of the face,
// false if it is its opposite
std::map<halfedge_descriptor, bool> border;
for(face_descriptor f : face_range)
{
for(halfedge_descriptor h :
halfedges_around_face(halfedge(f, g), g))
{
// halfedge_descriptor is model of `LessThanComparable`
bool from_face = (h < opposite(h, g));
halfedge_descriptor he = from_face ? h : opposite(h, g);
if (border.find(he) != border.end())
border.erase(he); //even number of appearances
else
border.insert(std::make_pair(he, from_face));//odd number of appearances
}
}
//copy them in out
typedef typename std::map<halfedge_descriptor, bool>::value_type HD_bool;
for(const HD_bool& hd : border)
{
if (!hd.second) // to get the border halfedge (which is not on the face)
*out++ = hd.first;
else
*out++ = opposite(hd.first, g);
}
return out;
}

template <typename Graph,
typename FaceRange,
typename HalfedgeOutputIterator,
typename NamedParameters>
HalfedgeOutputIterator border_halfedges_impl(const FaceRange& face_range,
typename boost::cgal_no_property::type,
HalfedgeOutputIterator out,
const Graph& g,
const NamedParameters& /* np */)
{
return border_halfedges_impl(face_range, out, g);
}

template <typename Graph,
typename FaceRange,
typename FaceIndexMap,
typename HalfedgeOutputIterator,
typename NamedParameters>
HalfedgeOutputIterator border_halfedges_impl(const FaceRange& face_range,
const FaceIndexMap& fmap,
HalfedgeOutputIterator out,
const Graph& g,
const NamedParameters& /* np */)
{
typedef typename boost::graph_traits<Graph>::halfedge_descriptor halfedge_descriptor;
typedef typename boost::graph_traits<Graph>::face_descriptor face_descriptor;

CGAL_assertion(BGL::internal::is_index_map_valid(fmap, num_faces(g), faces(g)));

std::vector<bool> present(num_faces(g), false);
for(face_descriptor fd : face_range)
present[get(fmap, fd)] = true;

for(face_descriptor fd : face_range)
{
for(halfedge_descriptor hd : halfedges_around_face(halfedge(fd, g), g))
{
halfedge_descriptor opp = opposite(hd, g);
if (is_border(opp, g) || !present[get(fmap, face(opp,g))])
*out++ = opp;
}
}

return out;
}

struct Dummy_PM
{
public:
typedef bool vertex_property_type;
};

} // namespace internal

/*!
* \ingroup PkgBGLHelperFct
*
* \brief collects the border halfedges of a surface patch defined as a face range.
*
* For each returned halfedge `h`, `opposite(h, g)` belongs to a face of the patch,
* but `face(h, g)` does not belong to the patch.
*
* @tparam Graph model of `HalfedgeGraph`
* @tparam FaceRange a model of `Range` with value type `boost::graph_traits<Graph>::%face_descriptor`
* @tparam HalfedgeOutputIterator model of `OutputIterator` holding `boost::graph_traits<Graph>::%halfedge_descriptor` for patch border
* @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
*
* @param g the polygon mesh to which the faces in `face_range` belong
* @param face_range the range of faces defining the patch whose border halfedges are collected
* @param out the output iterator that collects the border halfedges of the patch, seen from outside
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
* \cgalNamedParamsBegin
* \cgalParamNBegin{face_index_map}
* \cgalParamDescription{a property map associating to each face of `g` a unique index between `0` and `num_faces(g) - 1`}
* \cgalParamType{a class model of `ReadablePropertyMap` with `boost::graph_traits<Graph>::%face_descriptor`
* as key type and `std::size_t` as value type}
* \cgalParamDefault{an automatically indexed internal map}
* \cgalParamNEnd
* \cgalNamedParamsEnd
*
* @returns `out`
*
* @see `extract_boundary_cycles()`
*/
template <typename Graph,
typename FaceRange,
typename HalfedgeOutputIterator,
typename NamedParameters = parameters::Default_named_parameters>
HalfedgeOutputIterator border_halfedges(const FaceRange& face_range,
const Graph& g,
HalfedgeOutputIterator out,
const NamedParameters& np = parameters::default_values())
{
if (face_range.empty())
return out;

typedef typename CGAL::GetInitializedFaceIndexMap<Graph, NamedParameters>::const_type FIMap;
FIMap fim = CGAL::get_initialized_face_index_map(g, np);

return internal::border_halfedges_impl(face_range, fim, out, g, np);
}

template <typename Graph,
typename HalfedgeOutputIterator>
HalfedgeOutputIterator border_halfedges(const Graph& g,
HalfedgeOutputIterator out)
{
typedef typename boost::graph_traits<Graph>::halfedge_descriptor halfedge_descriptor;
for(halfedge_descriptor hd : halfedges(g))
if(is_border(hd, g))
*out++ = hd;
return out;
}

// counts the number of connected components of the boundary of the mesh.
//
// @tparam Graph model of `HalfedgeGraph`.
//
// @param g the polygon mesh to which `face_range` belong
//
template<typename Graph>
unsigned int number_of_borders(const Graph& g)
{
typedef typename boost::graph_traits<Graph>::halfedge_descriptor halfedge_descriptor;

unsigned int border_counter = 0;
std::unordered_set<halfedge_descriptor> visited;
for(halfedge_descriptor h : halfedges(g)){
if(visited.find(h)== visited.end()){
if(is_border(h,g)){
++border_counter;
for(halfedge_descriptor haf : halfedges_around_face(h, g)){
visited.insert(haf);
}
}
}
}

return border_counter;
}

/*!
* @ingroup PkgBGLHelperFct
*
* extracts boundary cycles as a list of halfedges, with one halfedge per border.
*
* @tparam Graph a model of `HalfedgeListGraph`
* @tparam OutputIterator a model of `OutputIterator` holding objects of type
* `boost::graph_traits<Graph>::%halfedge_descriptor`
*
* @param g a polygon mesh
* @param out an output iterator where the border halfedges will be put
*
* @see `border_halfedges()`
*/
template <typename Graph, typename OutputIterator>
OutputIterator extract_boundary_cycles(const Graph& g,
OutputIterator out)
{
typedef typename boost::graph_traits<Graph>::halfedge_descriptor halfedge_descriptor;

std::unordered_set<halfedge_descriptor> hedge_handled;
for(halfedge_descriptor h : halfedges(g))
{
if(is_border(h, g) && hedge_handled.insert(h).second)
{
*out++ = h;
for(halfedge_descriptor h2 : halfedges_around_face(h, g))
hedge_handled.insert(h2);
}
}
return out;
}

} // namespace CGAL

#endif // CGAL_BGL_BORDER_H
4 changes: 2 additions & 2 deletions BGL/include/CGAL/boost/graph/helpers.h
Original file line number Diff line number Diff line change
Expand Up @@ -945,7 +945,7 @@ void collect_garbage(Graph&)
/**
* \ingroup PkgBGLHelperFct
*
* removes all vertices, faces and halfedges from a graph. Calls
* removes all vertices, edges, and faces from a graph. Calls
* \link MutableHalfedgeGraph `remove_vertex()`\endlink,
* \link MutableHalfedgeGraph `remove_edge()`\endlink, and
* \link MutableFaceGraph `remove_face()`\endlink, for each vertex, edge, and face.
Expand Down Expand Up @@ -1000,7 +1000,7 @@ clear_impl(FaceGraph& g)
/**
* \ingroup PkgBGLHelperFct
*
* removes all vertices, faces and halfedges from a graph. Calls
* removes all vertices, edges, and faces from a graph. Calls
* \link MutableHalfedgeGraph `remove_vertex()`\endlink,
* \link MutableHalfedgeGraph `remove_edge()`\endlink, and
* \link MutableFaceGraph `remove_face()`\endlink, for each vertex, edge, and face.
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -229,7 +229,7 @@ triangulation.
The following example demonstrates how to detect planar surface patches and remesh them as coarsely
as possible using
\link CGAL::Polygon_mesh_processing::remesh_planar_patches(const TriangleMeshIn&,PolygonMeshOut&,const NamedParametersIn&,const NamedParametersOut&) `CGAL::Polygon_mesh_processing::remesh_planar_patches()` \endlink
from the \ref PkgPolygonMeshProcessing package.
from the \ref PkgPMPRemeshing package.
The resulting patches and segmentation are then used to build a conforming constrained Delaunay triangulation.

When the named parameter `plc_face_id` is specified, each constrained facet in the 3D triangulation
Expand Down Expand Up @@ -270,7 +270,7 @@ is used to associate each facet of the triangulation with the corresponding inpu
Given a PLC, the algorithms in this package can construct a conforming constrained Delaunay triangulation, provided
the input surface can be represented as a valid surface mesh or a collection of surface meshes,
and does not contain self-intersections.
Several preprocessing functions are available in the \ref PkgPolygonMeshProcessing package to help ensure these
Several preprocessing functions are available in the \ref PkgPMPMeshRepair package to help ensure these
preconditions are met.

The following example demonstrates how to construct a conforming constrained Delaunay triangulation from
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,8 @@ BGL
Kernel_23
Manual
Polygon_mesh_processing
PMP_Boolean_operations
PMP_Remeshing
Property_map
SMDS_3
STL_Extension
Expand Down
2 changes: 1 addition & 1 deletion Documentation/doc/Documentation/Third_party.txt
Original file line number Diff line number Diff line change
Expand Up @@ -278,7 +278,7 @@ for instance, on Windows, you can download it from <A HREF="https://www.zlib.net

\ceres is an open source C++ library for modeling and solving large, complicated optimization problems.

In \cgal, \ceres is used by the \ref PkgPolygonMeshProcessing package for mesh smoothing, which
In \cgal, \ceres is used by the \ref PkgPMPRemeshing package for mesh smoothing, which
requires solving complex non-linear least squares problems.

Visit the official website of the library at <A HREF="http://ceres-solver.org/index.html">`ceres-solver.org`</A>
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -256,7 +256,7 @@ through the points.
Each of these methods produce a triangle mesh stored in different
ways. If this output mesh is hampered by defects such as holes or
self-intersections, \cgal provide several algorithms to postprocess
it (hole filling, remeshing, etc.) in the package \ref PkgPolygonMeshProcessing "Polygon Mesh Processing".
it (hole filling, remeshing, etc.) in the package \ref PkgPMPMeshRepair "Polygon Mesh Processing".

We do not discuss these functions here as there are many
postprocessing possibilities whose relevance strongly depends on the
Expand Down
2 changes: 1 addition & 1 deletion Documentation/doc/Documentation/main.txt
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,7 @@ in the form of a C++ library.
\cgal offers data structures and algorithms like \ref PartTriangulationsAndDelaunayTriangulations "triangulations",
\ref PartVoronoiDiagrams "Voronoi diagrams", \ref PartPolygons, \ref PartPolyhedra,
\ref PartArrangements "arrangements of curves", \ref PartMeshing "mesh generation",
\ref PartGeometryProcessing "geometry processing", \ref PartConvexHullAlgorithms "convex hull algorithms",
\ref PartPolygonMeshProcessing "geometry processing (polygon meshes)", \ref PartPointSetProcessing "geometry processing (point sets)", \ref PartConvexHullAlgorithms "convex hull algorithms",
to name just a few.

All these data structures and algorithms operate on geometric objects
Expand Down
12 changes: 9 additions & 3 deletions Documentation/doc/Documentation/packages.txt
Original file line number Diff line number Diff line change
Expand Up @@ -111,9 +111,12 @@
\package_listing{Kinetic_surface_reconstruction}
\package_listing{Optimal_transportation_reconstruction_2}

\cgalPackageSection{PartGeometryProcessing,Geometry Processing}
\cgalPackageSection{PartPolygonMeshProcessing,Polygon Mesh Processing}

\package_listing{Polygon_mesh_processing}
\package_listing{PMP_Boolean_operations}
\package_listing{PMP_Remeshing}
\package_listing{PMP_Mesh_repair}
\package_listing{Subdivision_method_3}
\package_listing{Surface_mesh_segmentation}
\package_listing{Surface_mesh_simplification}
Expand All @@ -122,6 +125,11 @@
\package_listing{Surface_mesh_shortest_path}
\package_listing{Surface_mesh_skeletonization}
\package_listing{Surface_mesh_approximation}
\package_listing{Heat_method_3}
\package_listing{Surface_mesh_topology}

\cgalPackageSection{PartPointSetProcessing,Point Set Processing}

\package_listing{Ridges_3}
\package_listing{Jet_fitting_3}
\package_listing{Point_set_3}
Expand All @@ -130,8 +138,6 @@
\package_listing{Shape_regularization}
\package_listing{Stream_lines_2}
\package_listing{Classification}
\package_listing{Heat_method_3}
\package_listing{Surface_mesh_topology}

\cgalPackageSection{PartSearchStructures,Spatial Searching and Sorting}

Expand Down
Loading
Loading