CGAL 5.6 - CGAL and the Boost Graph Library
User Manual

Authors
Andreas Fabri, Fernando Cacciola, Philipp Moeller, and Ron Wein

Many geometric data structures can be interpreted as graphs as they consist of vertices and edges. This is the case for the halfedge data structure, for the polyhedral surface, for the arrangement, and for the 2D triangulation classes. With means of duality one can also interpret faces as vertices and edges between adjacent faces as edges of the dual graph.

The scope of CGAL is geometry and not graph algorithms. Nevertheless, this package provides the necessary classes and functions that enable using the algorithms of the Boost Graph Library [4] (BGL for short) with CGAL data structures.

Furthermore, this package extends the BGL by introducing concepts such as HalfedgeGraph and FaceGraph allowing to handle halfedges and faces. These concepts reflect the design of the halfedge data structure described in Chapter Halfedge Data Structures, with opposite halfedges and circular sequences of halfedges around vertices and around faces.

This chapter is organized as follows:

A Short Introduction to the Boost Graph Library

The algorithms of the BGL operate on models of various graph concepts. The traits class boost::graph_traits enable algorithms to determine the types of vertices and edges (similar to std::iterator_traits for iterators). Free functions that operate on graphs enable algorithms to obtain, for example, the source vertex of an edge, or all edges incident to a vertex. The algorithms use property maps to associate information with vertices and edges. The algorithms enable visitors to register callbacks that are called later on during the execution of the algorithms. Finally, the graph algorithms use the named parameter mechanism, which enables passing the arguments in arbitrary order.

Graph Concepts

The BGL introduces several graph concepts, which have different sets of characteristics and requirements. For example, iterating through all vertices or all edges in a graph, obtaining the outgoing or in-going edges of a vertex, inserting vertices and edges into a graph, and removing vertices and edges from a graph.

The Graph Traits Class

An algorithm operating on a graph model determines types with the help of the traits class boost::graph_traits. Such types are the vertex_descriptor, which is similar to a vertex handle in CGAL data structures, or the edge_descriptor, which is similar to the halfedge handle in the halfedge data structure or to the type Edge in 2D triangulations. There are also iterators, such as the vertex_iterator, which is similar to a vertex iterator in CGAL data structures, and the out_edge_iterator, which is similar to the edge circulator; it enables to iterate through the edges incident to a vertex. The iterators are similar and not equivalent, because their value type is a vertex_descriptor, whereas in CGAL handles, iterators, and circulators all have the same value type, namely the vertex or edge types.

Given a graph type G, definitions of descriptors and iterators look as follows:

boost::graph_traits<Graph>::vertex_descriptor vd;
boost::graph_traits<Graph>::edge_iterator ei;
...

Free Functions for Exploring a Graph

Algorithms obtain incidence information in graphs with the help of global functions such as:

  • std::pair<vertex_iterator,vertex_iterator> vertices(const Graph& g); to obtain an iterator range providing access to all the vertices, or
  • int num_vertices(const Graph&); to obtain the number of vertices of a graph, or
  • vertex_descriptor source(edge_descriptor, const Graph&); to obtain the source vertex of an edge.

Note, that the way we have written the types is a simplification; in reality, the signature of the first of the above functions is:

typedef boost::graph_traits<Graph>::vertex_iterator vertex_iterator;
std::pair<vertex_iterator,vertex_iterator> vertices(const Graph& g);

Property Maps

Another feature extensively used in the BGL is the property map, which is offered by the Boost Property Map Library. Property maps are a general purpose interface for mapping key objects to corresponding value objects.

The BGL uses property maps to associate information with vertices and edges. This mechanism uses a traits class (boost::property_traits) and free functions to read (get) and write (put) information in vertices, edges, and also in halfedges and faces for models of the CGAL graph concepts. For example, the BGL Dijksta's shortest path algorithm writes the predecessor of each vertex, as well as the distance to the source in such a property map.

Some default property maps are associated with the graph types. They are called internal property maps and can be retrieved using an overload of the function get(). For example,

pm = get(boost::vertex_index, g)

returns a property map that associates an index in the range [0, num_vertices(g)) with each vertex of the graph. This reduces the number of parameters to pass. The data itself may be stored in the vertex or the edge, or it may be stored in an external data structure, or it may be computed on the fly. This is an implementation detail of a particular property map.

See also Chapter CGAL and Boost Property Maps.

Visitors

Visitors are objects that provide functions to be called at specified event points by the algorithm that they visit. The functions as well as the event points are algorithm-specific. Examples of such event points in graph algorithms are when a vertex is traversed the first time, or when all outgoing edges of a vertex have been traversed.

See also Section Visitor Concepts in the BGL manual.

Named Parameters

The notion of named parameters was introduced in the BGL, and allow the user to specify only those parameters which are really needed, by name, making the parameter ordering unimportant. See also this page in the manual of the BGL for more information.

Say there is a function f() that takes 3 parameters called name, age and gender, and you have variables n, a and g to pass as parameters to that function. Without named parameters, you would call it like this: f(n,a,g), whereas with named parameters, you call it like this: f(name(n).age(a).gender(g)). That is, you give each parameter a name by wrapping it into a function whose name matches that of the parameter. The entire list of named parameters is really a composition of function calls separated by a dot ("."). Thus, if the function takes a mix of mandatory and named parameters, you use a comma to separate the last non-named parameter from the first named parameters, like this:

f(non_named_par0, non_named_par1, name(n).age(a).gender(g))

When you use named parameters, the ordering is irrelevant, so f(name(n).age(a).gender(g)) is equivalent to f(age(a).gender(g).name(n)), and you can just omit any named parameter that has a default value.

The sequence of named parameters should start with CGAL::parameters::.

Example

Below is a sample call of a function that uses the optional BGL named parameters.

// pmesh : polygon mesh with patches to be refined
// faces : the range of faces defining the patches to refine
// faces_out : output iterator into which descriptors of new faces are put
// vertices_out : output iterator into which descriptors of new vertices are put
// vertex_point_map : the property map with the points associated to the vertices of `pmesh`
// density_control_factor : factor to control density of the output mesh
refine(pmesh,
faces,
faces_out,
vertices_out,
CGAL::parameters::vertex_point_map(vpmap)
.density_control_factor(d));

Header Files, Namespaces, and Naming Conventions

This package provides the necessary classes and functions that enable using CGAL data structures as models of the BGL graph concepts. To this end, we offer partial specializations of the boost::graph_traits<Graph> for various CGAL packages. For each such package, denoted PACKAGE, the partial specializations live in the namespace boost and are located in the header file CGAL/boost/graph/graph_traits_PACKAGE.h. Free functions are in the namespace CGAL, and the compiler uses argument-dependent lookup to find them. Euler operations, described in Section Euler Operations, are in the namespace CGAL::Euler, as the function remove_face() is at the same time a low-level and an Euler operation. Concerning the naming conventions, we have to use those of the BGL, as to fulfill the requirements of the concepts defined in the BGL.

Note that these partial specializations are often providing more than is required, making these classes not only models of the graph concepts of the BGL, but also models of the CGAL graph concepts, that will be described in detail in Section Extensions of the BGL. Correspondence tables between the types of a CGAL data structure and their BGL equivalents can be found in the Specializations of boost::graph_traits documentation page.

We present in the following sections some examples of utilization of some CGAL data structures as BGL graphs.

The Class Surface_mesh as Model of the Boost Graph Concept

The class Surface_mesh is a model of most of the graph concepts of the BGL as well as the concepts provided by CGAL. A complete list can be found in the documentation of boost::graph_traits . The examples show how to use some of the BGL algorithms with Surface_mesh and show how to use the concepts provided by CGAL to implement a simple algorithm.

Example: Minimum Spanning Tree of a Surface_mesh

The following example program computes the minimum spanning tree on a surface mesh. More examples can be found in Chapters Triangulated Surface Mesh Simplification, Triangulated Surface Mesh Segmentation, and Triangulated Surface Mesh Deformation.

The surface mesh class uses integer indices to address vertices and edges, and it comes with a built-in property mechanism that maps nicely on the BGL.


File BGL_surface_mesh/prim.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <iostream>
#include <fstream>
#include <CGAL/boost/graph/prim_minimum_spanning_tree.h>
typedef Kernel::Point_3 Point;
typedef boost::graph_traits<Mesh>::vertex_descriptor vertex_descriptor;
int main(int argc, char* argv[])
{
const std::string filename = (argc>1) ? argv[1] : CGAL::data_file_path("meshes/prim.off");
Mesh P;
if(!CGAL::IO::read_polygon_mesh(filename, P))
{
std::cerr << "Invalid input." << std::endl;
return 1;
}
Mesh::Property_map<vertex_descriptor,vertex_descriptor> predecessor;
predecessor = P.add_property_map<vertex_descriptor,vertex_descriptor>("v:predecessor").first;
boost::prim_minimum_spanning_tree(P, predecessor, boost::root_vertex(*vertices(P).first));
std::cout << "#VRML V2.0 utf8\n"
"DirectionalLight {\n"
"direction 0 -1 0\n"
"}\n"
"Shape {\n"
" appearance Appearance {\n"
" material Material { emissiveColor 1 0 0}}\n"
" geometry\n"
" IndexedLineSet {\n"
" coord Coordinate {\n"
" point [ \n";
for(vertex_descriptor vd : vertices(P)){
std::cout << " " << P.point(vd) << "\n";
}
std::cout << " ]\n"
" }\n"
" coordIndex [\n";
for(vertex_descriptor vd : vertices(P)){
if(predecessor[vd]!=vd){
std::cout << " " << std::size_t(vd) << ", " << std::size_t(predecessor[vd]) << ", -1\n";
}
}
std::cout << "]\n"
" }#IndexedLineSet\n"
"}# Shape\n";
P.remove_property_map(predecessor);
return 0;
}

The Class Polyhedron_3 as Model of the Boost Graph Concept

The class Polyhedron_3 is a model of most of the graph concepts of the BGL as well as the concepts provided by CGAL. A complete list can be found in the documentation of boost::graph_traits . The examples show how to use some of the BGL algorithms with Polyhedron_3 and show how to use the concepts provided by CGAL to implement a simple algorithm.

Example: Minimum Spanning Tree of a Polyhedral Surface

The following example program computes the minimum spanning tree on a polyhedral surface. More examples can be found in the Chapter Triangulated Surface Mesh Simplification.


File BGL_polyhedron_3/kruskal.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Polyhedron_3.h>
#include <iostream>
#include <list>
#include <CGAL/boost/graph/kruskal_min_spanning_tree.h>
typedef Kernel::Vector_3 Vector;
typedef Kernel::Point_3 Point;
typedef CGAL::Polyhedron_3<Kernel> Polyhedron;
typedef boost::graph_traits<Polyhedron>::vertex_descriptor vertex_descriptor;
typedef boost::graph_traits<Polyhedron>::vertex_iterator vertex_iterator;
typedef boost::graph_traits<Polyhedron>::edge_descriptor edge_descriptor;
// The BGL makes heavy use of indices associated to the vertices
// We use a std::map to store the index
typedef std::map<vertex_descriptor,int> Vertex_index_map;
Vertex_index_map vertex_index_map;
// A std::map is not a property map, because it is not lightweight
typedef boost::associative_property_map<Vertex_index_map> Vertex_index_pmap;
Vertex_index_pmap vertex_index_pmap(vertex_index_map);
void
kruskal(const Polyhedron& P)
{
// associate indices to the vertices
vertex_iterator vb, ve;
int index = 0;
// boost::tie assigns the first and second element of the std::pair
// returned by boost::vertices to the variables vb and ve
for(boost::tie(vb, ve)=vertices(P); vb!=ve; ++vb){
vertex_index_pmap[*vb]= index++;
}
// We use the default edge weight which is the length of the edge
// This property map is defined in graph_traits_Polyhedron_3.h
// In the function call you can see a named parameter: vertex_index_map
std::list<edge_descriptor> mst;
boost::kruskal_minimum_spanning_tree(P,
std::back_inserter(mst),
boost::vertex_index_map(vertex_index_pmap));
std::cout << "#VRML V2.0 utf8\n"
"Shape {\n"
" appearance Appearance {\n"
" material Material { emissiveColor 1 0 0}}\n"
" geometry\n"
" IndexedLineSet {\n"
" coord Coordinate {\n"
" point [ \n";
for(boost::tie(vb, ve) = vertices(P); vb!=ve; ++vb){
std::cout << " " << (*vb)->point() << "\n";
}
std::cout << " ]\n"
" }\n"
" coordIndex [\n";
for(std::list<edge_descriptor>::iterator it = mst.begin(); it != mst.end(); ++it)
{
edge_descriptor e = *it ;
vertex_descriptor s = source(e,P);
vertex_descriptor t = target(e,P);
std::cout << " " << vertex_index_pmap[s] << ", " << vertex_index_pmap[t] << ", -1\n";
}
std::cout << "]\n"
" }#IndexedLineSet\n"
"}# Shape\n";
}
int main() {
Polyhedron P;
Point a(1,0,0);
Point b(0,1,0);
Point c(0,0,1);
Point d(0,0,0);
P.make_tetrahedron(a,b,c,d);
kruskal(P);
return 0;
}

Example: Using Vertices, and Edges with an ID

The following example program shows a call to the BGL Kruskal's minimum spanning tree algorithm accessing the id() field stored in a polyhedron vertex.

The main function illustrates the access to the id() field.


File BGL_polyhedron_3/kruskal_with_stored_id.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/Polyhedron_items_with_id_3.h>
#include <iostream>
#include <list>
#include <CGAL/boost/graph/kruskal_min_spanning_tree.h>
typedef Kernel::Point_3 Point;
typedef boost::graph_traits<Polyhedron>::vertex_descriptor vertex_descriptor;
typedef boost::graph_traits<Polyhedron>::vertex_iterator vertex_iterator;
typedef boost::graph_traits<Polyhedron>::edge_descriptor edge_descriptor;
void
kruskal( const Polyhedron& P)
{
// We use the default edge weight which is the length of the edge
// This property map is defined in graph_traits_Polyhedron_3.h
// This function call requires a vertex_index_map named parameter which
// when omitted defaults to "get(vertex_index,graph)".
// That default works here because the vertex type has an "id()" method
// field which is used by the vertex_index internal property.
std::list<edge_descriptor> mst;
boost::kruskal_minimum_spanning_tree(P,std::back_inserter(mst));
std::cout << "#VRML V2.0 utf8\n"
"Shape {\n"
"appearance Appearance {\n"
"material Material { emissiveColor 1 0 0}}\n"
"geometry\n"
"IndexedLineSet {\n"
"coord Coordinate {\n"
"point [ \n";
vertex_iterator vb, ve;
for(boost::tie(vb,ve) = vertices(P); vb!=ve; ++vb){
std::cout << (*vb)->point() << "\n";
}
std::cout << "]\n"
"}\n"
"coordIndex [\n";
for(std::list<edge_descriptor>::iterator it = mst.begin(); it != mst.end(); ++it){
std::cout << source(*it,P)->id()
<< ", " << target(*it,P)->id() << ", -1\n";
}
std::cout << "]\n"
"}#IndexedLineSet\n"
"}# Shape\n";
}
int main() {
Polyhedron P;
Point a(1,0,0);
Point b(0,1,0);
Point c(0,0,1);
Point d(0,0,0);
P.make_tetrahedron(a,b,c,d);
// associate indices to the vertices using the "id()" field of the vertex.
vertex_iterator vb, ve;
int index = 0;
// boost::tie assigns the first and second element of the std::pair
// returned by boost::vertices to the variables vit and ve
for(boost::tie(vb,ve)=vertices(P); vb!=ve; ++vb ){
vertex_descriptor vd = *vb;
vd->id() = index++;
}
kruskal(P);
return 0;
}

Triangulations as Models of the Boost Graph Concept

Triangulations have vertices and faces, allowing for a direct translation as a graph. A halfedge is defined as a pair of a face handle and the index of the edge. A complete list can be found in the documentation of boost::graph_traits .

A classical example for an algorithm that is a combination of computational geometry and graph theory is the Euclidean Minimum Spanning Tree for a point set in the plane. It can be computed by running the minimum spanning tree algorithm on a Delaunay triangulation of the point set.

Example: Euclidean Minimum Spanning Tree

In the following example we create a Delaunay triangulation and run Kruskal's minimum spanning tree algorithm on it. Because the vertex handles of the triangulation are not indices in an array, we have to provide a property map that maps vertex handles to integers in the range [0, t.number_of_vertices()).


File BGL_triangulation_2/emst.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Delaunay_triangulation_2.h>
#include <CGAL/boost/graph/graph_traits_Delaunay_triangulation_2.h>
#include <CGAL/boost/graph/kruskal_min_spanning_tree.h>
#include <fstream>
#include <iostream>
#include <map>
typedef K::Point_2 Point;
typedef CGAL::Delaunay_triangulation_2<K> Triangulation;
typedef boost::graph_traits<Triangulation>::vertex_descriptor vertex_descriptor;
typedef boost::graph_traits<Triangulation>::vertex_iterator vertex_iterator;
typedef boost::graph_traits<Triangulation>::edge_descriptor edge_descriptor;
// The BGL makes use of indices associated to the vertices
// We use a std::map to store the index
typedef std::map<vertex_descriptor,int> VertexIndexMap;
// A std::map is not a property map, because it is not lightweight
typedef boost::associative_property_map<VertexIndexMap> VertexIdPropertyMap;
int main(int argc,char* argv[])
{
const char* filename = (argc > 1) ? argv[1] : "data/points.xy";
std::ifstream input(filename);
Triangulation tr;
Point p;
while(input >> p)
tr.insert(p);
// Associate indices to the vertices
VertexIndexMap vertex_id_map;
VertexIdPropertyMap vertex_index_pmap(vertex_id_map);
int index = 0;
for(vertex_descriptor vd : vertices(tr))
vertex_id_map[vd] = index++;
// We use the default edge weight which is the squared length of the edge
// This property map is defined in graph_traits_Triangulation_2.h
// In the function call you can see a named parameter: vertex_index_map
std::list<edge_descriptor> mst;
boost::kruskal_minimum_spanning_tree(tr, std::back_inserter(mst),
vertex_index_map(vertex_index_pmap));
std::cout << "The edges of the Euclidean minimum spanning tree:" << std::endl;
for(edge_descriptor ed : mst)
{
vertex_descriptor svd = source(ed, tr);
vertex_descriptor tvd = target(ed, tr);
Triangulation::Vertex_handle sv = svd;
Triangulation::Vertex_handle tv = tvd;
std::cout << "[ " << sv->point() << " | " << tv->point() << " ] " << std::endl;
}
return EXIT_SUCCESS;
}

Example: Storing the Vertex ID in the Vertex

The algorithms of the BGL extensively use of the indices of vertices. In the previous example we stored the indices in a std::map and turned that map in a property map. This property map was then passed as argument to the shortest path function.

If the user does not pass explicitly a property map, the graph algorithms use the property map returned by the call get(boost::vertex_index,ft). This property map assumes that the vertex has a member function id() that returns a reference to an int. Therefore CGAL offers a class Triangulation_vertex_base_with_id_2. It is in the user's responsibility to set the indices properly.

The example further illustrates that the graph traits also works for the Delaunay triangulation.


File BGL_triangulation_2/dijkstra_with_internal_properties.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Delaunay_triangulation_2.h>
#include <CGAL/Triangulation_vertex_base_with_id_2.h>
#include <CGAL/boost/graph/graph_traits_Delaunay_triangulation_2.h>
#include <CGAL/boost/graph/dijkstra_shortest_paths.h>
#include <fstream>
typedef K::Point_2 Point;
typedef CGAL::Triangulation_data_structure_2<Tvb, Tfb> Tds;
typedef boost::graph_traits<Triangulation>::vertex_descriptor vertex_descriptor;
typedef boost::graph_traits<Triangulation>::vertex_iterator vertex_iterator;
typedef boost::property_map<Triangulation, boost::vertex_index_t>::type VertexIdPropertyMap;
int main(int argc,char* argv[])
{
const char* filename = (argc > 1) ? argv[1] : "data/points.xy";
std::ifstream input(filename);
Triangulation tr;
Point p;
while(input >> p)
tr.insert(p);
// associate indices to the vertices
int index = 0;
for(vertex_descriptor vd : vertices(tr))
vd->id()= index++;
VertexIdPropertyMap vertex_index_pmap = get(boost::vertex_index, tr);
// Dijkstra's shortest path needs property maps for the predecessor and distance
std::vector<vertex_descriptor> predecessor(num_vertices(tr));
boost::iterator_property_map<std::vector<vertex_descriptor>::iterator, VertexIdPropertyMap>
predecessor_pmap(predecessor.begin(), vertex_index_pmap);
std::vector<double> distance(num_vertices(tr));
boost::iterator_property_map<std::vector<double>::iterator, VertexIdPropertyMap>
distance_pmap(distance.begin(), vertex_index_pmap);
vertex_descriptor source = *vertices(tr).first;
std::cout << "\nStart dijkstra_shortest_paths at " << source->point() << std::endl;
boost::dijkstra_shortest_paths(tr, source, distance_map(distance_pmap)
.predecessor_map(predecessor_pmap));
for(vertex_descriptor vd : vertices(tr))
{
std::cout << vd->point() << " [" << vd->id() << "] ";
std::cout << " has distance = " << get(distance_pmap,vd) << " and predecessor ";
vd = get(predecessor_pmap,vd);
std::cout << vd->point() << " [" << vd->id() << "]\n";
}
return EXIT_SUCCESS;
}

Arrangements as Models of the Boost Graph Concept

CGAL offers a partial specialization of the boost graph traits for its arrangement class as well as for its dual graph.

Example for the Arrangement as Graph

Arrangement instances are adapted to boost graphs by specializing the boost::graph_traits template for Arrangement_2 instances. The graph-traits states the graph concepts that the arrangement class models (see below) and defines the types required by these concepts.

In this specialization the Arrangement_2 vertices correspond to the graph vertices, where two vertices are adjacent if there is at least one halfedge connecting them. More precisely, Arrangement_2::Vertex_handle is the graph-vertex type, while Arrangement_2::Halfedge_handle is the graph-edge type. As halfedges are directed, we consider the graph to be directed as well. Moreover, as several interior-disjoint \( x\)-monotone curves (say circular arcs) may share two common endpoints, inducing an arrangement with two vertices that are connected with several edges, we allow parallel edges in our boost graph.

Given an Arrangement_2 instance, we can efficiently traverse its vertices and halfedges. Thus, the arrangement graph is a model of the concepts VertexListGraph and EdgeListGraph introduced by the BGL. At the same time, we use an iterator adapter of the circulator over the halfedges incident to a vertex (Halfedge_around_target_circulator - see Section Traversal Methods for an Arrangement Vertex of the chapter on arrangements), so it is possible to go over the ingoing and outgoing edges of a vertex in linear time. Thus, our arrangement graph is a model of the concept BidirectionalGraph (this concept refines IncidenceGraph, which requires only the traversal of outgoing edges).

It is important to notice that the vertex descriptors we use are Vertex_handle objects and not vertex indices. However, in order to gain more efficiency in most BGL algorithm, it is better to have them indexed \( 0, 1, \ldots, (n-1)\), where \( n\) is the number of vertices. We therefore introduce the Arr_vertex_index_map<Arrangement> class-template, which maintains a mapping of vertex handles to indices, as required by the BGL. An instance of this class must be attached to a valid arrangement vertex when it is created. It uses the notification mechanism (see Section The Notification Mechanism) to automatically maintain the mapping of vertices to indices, even when new vertices are inserted into the arrangement or existing vertices are removed.

A complete description of the types correspondences can be found in the documentation of boost::graph_traits .

In most algorithm provided by the BGL, the output is given by property maps, such that each map entry corresponds to a vertex. For example, when we compute the shortest paths from a given source vertex \( s\) to all other vertices we can obtain a map of distances and a map of predecessors - namely for each \( v\) vertex we have its distance from \( s\) and a descriptor of the vertex that precedes \( v\) in the shortest path from \( s\). If the vertex descriptors are simply indices, one can use vectors to efficiently represent the property maps. As this is not the case with the arrangement graph, we offer the Arr_vertex_property_map<Arrangement,Type> template allows for an efficient mapping of Vertex_handle objects to properties of type Type. Note however that unlike the Arr_vertex_index_map class, the vertex property-map class is not kept synchronized with the number of vertices in the arrangement, so it should not be reused in calls to the BGL functions in case the arrangement is modified in between these calls.

ex_bgl.png
Figure 103.1 An arrangement of 7 line segments, as constructed by bgl_primal_adapter.cpp and bgl_dual_adapter.cpp. The breadth-first visit times for the arrangement faces, starting from the unbounded face \( f_0\), are shown in brackets.

In the following example we construct an arrangement of 7 line segments, as shown in Figure 103.1, then use the BGL Dijkstra's shortest-paths algorithm to compute the graph distance of all vertices from the leftmost vertex in the arrangement \( v_0\). Note the usage of the Arr_vertex_index_map and the Arr_vertex_property_map classes. The latter one, instantiated by the type double is used to map vertices to their distances from \( v_0\).


File BGL_arrangement_2/primal.cpp

// Adapting an arrangement to a BGL graph.
#include "arr_rational_nt.h"
#include <CGAL/Cartesian.h>
#include <CGAL/Arr_segment_traits_2.h>
#include <CGAL/Arrangement_2.h>
#include <CGAL/graph_traits_Arrangement_2.h>
#include <CGAL/Arr_vertex_index_map.h>
#include <CGAL/boost/graph/dijkstra_shortest_paths.h>
#include <CGAL/property_map.h>
typedef Traits_2::Point_2 Point_2;
typedef Traits_2::X_monotone_curve_2 Segment_2;
typedef CGAL::Arrangement_2<Traits_2> Arrangement_2;
// A functor used to compute the length of an edge.
class Edge_length_func
{
public:
// Boost property type definitions:
typedef boost::readable_property_map_tag category;
typedef double value_type;
typedef value_type reference;
double operator()(Arrangement_2::Halfedge_handle e) const
{
const double x1 = CGAL::to_double (e->source()->point().x());
const double y1 = CGAL::to_double (e->source()->point().y());
const double x2 = CGAL::to_double (e->target()->point().x());
const double y2 = CGAL::to_double (e->target()->point().y());
const double diff_x = x2 - x1;
const double diff_y = y2 - y1;
return std::sqrt(diff_x*diff_x + diff_y*diff_y);
}
};
double get(const Edge_length_func& edge_length, Arrangement_2::Halfedge_handle e)
{
return edge_length(e);
}
int main()
{
Arrangement_2 arr;
// Construct an arrangement of seven intersecting line segments.
// We keep a handle for the vertex v_0 that corresponds to the point (1,1).
insert_non_intersecting_curve (arr, Segment_2 (Point_2 (1, 1),
Point_2 (7, 1)));
Arrangement_2::Vertex_handle v0 = e->source();
insert (arr, Segment_2 (Point_2 (1, 1), Point_2 (3, 7)));
insert (arr, Segment_2 (Point_2 (1, 4), Point_2 (7, 1)));
insert (arr, Segment_2 (Point_2 (2, 2), Point_2 (9, 3)));
insert (arr, Segment_2 (Point_2 (2, 2), Point_2 (4, 4)));
insert (arr, Segment_2 (Point_2 (7, 1), Point_2 (9, 3)));
insert (arr, Segment_2 (Point_2 (3, 7), Point_2 (9, 3)));
// Create a mapping of the arrangement vertices to indices.
// Perform Dijkstra's algorithm from the vertex v0.
Edge_length_func edge_length;
boost::vector_property_map<double, CGAL::Arr_vertex_index_map<Arrangement_2> > dist_map(static_cast<unsigned int>(arr.number_of_vertices()), index_map);
boost::dijkstra_shortest_paths(arr, v0,
boost::vertex_index_map(index_map).
weight_map(edge_length).
distance_map(dist_map));
// Print the results:
std::cout << "The distances of the arrangement vertices from ("
<< v0->point() << ") :" << std::endl;
for (vit = arr.vertices_begin(); vit != arr.vertices_end(); ++vit)
std::cout << "(" << vit->point() << ") at distance "
<< dist_map[vit] << std::endl;
return 0;
}

Example for the Dual of an Arrangement as Graph

It is possible to give a dual graph representation for an arrangement instance, such that each arrangement face corresponds to a graph vertex and two vertices are adjacent iff the corresponding faces share a common edge on their boundaries. This is done by specializing the boost::graph_traits template for Dual<Arrangement_2> instances, where Dual<Arrangement_2> is a template specialization that gives a dual interpretation to an arrangement instance.

In dual representation, Arrangement_2::Face_handle is the graph-vertex type, while Arrangement_2::Halfedge_handle is the graph-edge type. We treat the graph edges as directed, such that a halfedge e is directed from \( f_1\), which is its incident face, to \( f_2\), which is the incident face of its twin halfedge. As two arrangement faces may share more than a single edge on their boundary, we allow parallel edges in our boost graph. As is the case in the primal graph, the dual arrangement graph is also a model of the concepts VertexListGraph, EdgeListGraph and BidirectionalGraph (thus also of IncidenceGraph).

Since we use Face_handle objects as the vertex descriptors, we define the Arr_face_index_map<Arrangement> class-template, which maintains an efficient mapping of face handles to indices. We also provide the template Arr_face_property_map<Arrangement,Type> for associating arbitrary data with the arrangement faces.

In the following example we construct the same arrangement as in example ex_bgl_primal_adapter.cpp (see Figure 34.49), and perform breadth-first search on the graph faces, starting from the unbounded face. We extend the DCEL faces with an unsigned integer, marking the discover time of the face and use a breadth-first-search visitor to obtain these times and update the faces accordingly:


File BGL_arrangement_2/arrangement_dual.cpp

// Adapting the dual of an arrangement to a BGL graph.
#include "arr_rational_nt.h"
#include <CGAL/Cartesian.h>
#include <CGAL/Arr_segment_traits_2.h>
#include <CGAL/Arr_extended_dcel.h>
#include <CGAL/Arrangement_2.h>
#include <CGAL/graph_traits_dual_arrangement_2.h>
#include <CGAL/Arr_face_index_map.h>
#include <climits>
#include <CGAL/boost/graph/breadth_first_search.h>
#include <boost/graph/visitors.hpp>
#include "arr_print.h"
// A property map that reads/writes the information to/from the extended
// face.
template <typename Arrangement, class Type>
class Extended_face_property_map {
public:
typedef typename Arrangement::Face_handle Face_handle;
// Boost property type definitions.
typedef boost::read_write_property_map_tag category;
typedef Type value_type;
typedef value_type& reference;
typedef Face_handle key_type;
// The get function is required by the property map concept.
friend reference get(const Extended_face_property_map&, key_type key)
{ return key->data(); }
// The put function is required by the property map concept.
friend void put(const Extended_face_property_map&, key_type key, value_type val)
{ key->set_data(val); }
};
typedef CGAL::Arrangement_2<Traits_2, Dcel> Ex_arrangement;
typedef CGAL::Dual<Ex_arrangement> Dual_arrangement;
typedef Extended_face_property_map<Ex_arrangement,unsigned int>
Face_property_map;
typedef Kernel::Point_2 Point_2;
typedef Kernel::Segment_2 Segment_2;
int main()
{
// Construct an arrangement of seven intersecting line segments.
Point_2 p1(1, 1), p2(1, 4), p3(2, 2), p4(3, 7), p5(4, 4), p6(7, 1), p7(9, 3);
Ex_arrangement arr;
insert(arr, Segment_2(p1, p6));
insert(arr, Segment_2(p1, p4)); insert(arr, Segment_2(p2, p6));
insert(arr, Segment_2(p3, p7)); insert(arr, Segment_2(p3, p5));
insert(arr, Segment_2(p6, p7)); insert(arr, Segment_2(p4, p7));
// Create a mapping of the arrangement faces to indices.
Face_index_map index_map(arr);
// Perform breadth-first search from the unbounded face, using the event
// visitor to associate each arrangement face with its discover time.
unsigned int time = 0;
boost::breadth_first_search(Dual_arrangement(arr), arr.unbounded_face(),
boost::vertex_index_map(index_map).visitor
(boost::make_bfs_visitor
(stamp_times(Face_property_map(), time,
boost::on_discover_vertex()))));
// Print the discover time of each arrangement face.
Ex_arrangement::Face_iterator fit;
for (fit = arr.faces_begin(); fit != arr.faces_end(); ++fit) {
std::cout << "Discover time " << fit->data() << " for ";
if (fit != arr.unbounded_face()) {
std::cout << "face ";
print_ccb<Ex_arrangement>(fit->outer_ccb());
}
else std::cout << "the unbounded face." << std::endl;
}
return 0;
}

Extensions of the BGL

The previous sections introduced partial specializations and free functions so that several CGAL data structures are adapted as models of some of the BGL graph concepts. In this section, we introduce new concepts, iterators, and property maps inspired by the functionalities of the BGL.

Graph concepts

In order to match Halfedge Data Structures more closely and to enable writing generic algorithms which operate on data structures that have faces and halfedges, we define a set of new concepts to extend the graph concepts of the BGL:

  • HalfedgeGraph refines Graph with operations to accommodate halfedge data structures: given a halfedge, say h, the concept HalfedgeGraph requires the provision of the halfedge opposite to h, the halfedge that succeeds h, and the halfedge that precedes h.
  • HalfedgeListGraph adds the requirement for efficient traversal of the halfedges of the graph.
  • MutableHalfedgeGraph adds the requirement for operations to change next/previous relations and to adjust the target of a halfedge.
  • FaceGraph adds the requirements to explicitly handle faces in a graph, to provide quick access to the incident halfedges of a face, and to enable access from every halfedge to an adjacent face.
  • FaceListGraph adds the requirement for efficient traversal of the faces of a graph.
  • MutableFaceGraph adds requirements to change adjacency of faces and halfedges, and to remove and add faces.

A summary of the expressions and types associated with each of these concepts as well as a refinement relation graph can be found in the Concepts documentation page.

Iterators and Circulators

By combining basic operations on graphs, we create various useful iterators and circulators to traverse specific types of elements. For example:

  • Starting at a halfedge h of a halfedge graph g, applying several times next(h,g) brings us back to the halfedge where we started. All halfedges traversed on the way are incident to the same face.
  • Using the composition of the functions next(h,g) and opposite(h,g) results in another cycle, namely the cycle of halfedges which are incident to the same vertex.

A complete list of these traversal tools can be found in the reference manual.

For convenience, two iterator and circulator types enable the traversal of all halfedges incident to a given face, and all halfedges having a given vertex as target. These types are not part of the concept HalfedgeGraph, but they are class templates that work for any model of this concept.

Example: Finding Incident Vertices in a HalfedgeGraph

The following example shows several functions to navigate in a HalfedgeGraph. We have two implementations of the operation that finds the vertices adjacent to a vertex v.

Let us have a look at the first version. Given a vertex descriptor vd, we first call halfedge(vd,g) to obtain the halfedge with vd as target. Applying source() then gives us an adjacent vertex. We then get to the next halfedge with vd as target, by first going to the next halfedge around the face, and then going to the opposite halfedge.

The second version does the next() and opposite() jumping with an iterator. Note that when calling source() we have to dereference hi, as the function expects a halfedge descriptor and not a halfedge iterator. Also observe that halfedges_around_target() expects a halfedge, and not a vertex. This provides the user with the ability to start the traversal at a specific halfedge incident to the input vertex (and not the arbitrary incident halfedge stored in the vertex record.)


File BGL_polyhedron_3/incident_vertices.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/boost/graph/iterator.h>
#include <iostream>
#include <fstream>
typedef CGAL::Polyhedron_3<Kernel> Polyhedron;
typedef boost::graph_traits<Polyhedron> GraphTraits;
typedef GraphTraits::vertex_descriptor vertex_descriptor;
typedef CGAL::Halfedge_around_target_iterator<Polyhedron> halfedge_around_target_iterator;
template <typename OutputIterator>
adjacent_vertices_V1(const Polyhedron& g,
vertex_descriptor vd,
{
typename GraphTraits::halfedge_descriptor hb = halfedge(vd,g), done(hb);
do {
*out++ = source(hb,g);
hb = opposite(next(hb,g),g);
} while(hb!= done);
return out;
}
template <typename OutputIterator>
adjacent_vertices_V2(const Polyhedron& g,
vertex_descriptor vd,
{
halfedge_around_target_iterator hi, he;
for(boost::tie(hi, he) = halfedges_around_target(halfedge(vd,g),g); hi != he; ++hi)
{
*out++ = source(*hi,g);
}
return out;
}
int main(int argc, char** argv)
{
std::ifstream in((argc>1)?argv[1]:CGAL::data_file_path("meshes/cube_poly.off"));
Polyhedron P;
in >> P;
GraphTraits::vertex_iterator vi = vertices(P).first;
std::list<vertex_descriptor> V;
adjacent_vertices_V1(P, *vi, std::back_inserter(V));
++vi;
adjacent_vertices_V2(P, *vi, std::back_inserter(V));
std::cerr << "done\n";
return 0;
}

Example: Calculating Facet Normals using HalfedgeGraph

The following example program shows a simple algorithm for calculating facet normals for a polyhedron using the BGL API. A boost::vector_property_map is used to to store the calculated normals instead of changing the Polyhedron items class.


File BGL_polyhedron_3/normals.cpp

#include <fstream>
#include <boost/graph/graph_traits.hpp>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/Polyhedron_items_with_id_3.h>
#include <CGAL/property_map.h>
typedef Kernel::Point_3 Point;
typedef Kernel::Vector_3 Vector;
template<typename HalfedgeGraph,
typename PointMap,
typename NormalMap>
void calculate_face_normals(const HalfedgeGraph& g,
PointMap pm,
NormalMap nm)
{
typedef boost::graph_traits<HalfedgeGraph> GraphTraits;
typedef typename GraphTraits::face_iterator face_iterator;
typedef typename GraphTraits::halfedge_descriptor halfedge_descriptor;
typedef typename boost::property_traits<PointMap>::value_type point;
typedef typename boost::property_traits<NormalMap>::value_type normal;
face_iterator fb, fe;
for(boost::tie(fb, fe) = faces(g); fb != fe; ++fb)
{
halfedge_descriptor edg = halfedge(*fb, g);
halfedge_descriptor edgb = edg;
point p0 = pm[target(edg, g)];
edg = next(edg, g);
point p1 = pm[target(edg, g)];
edg = next(edg, g);
point p2 = pm[target(edg, g)];
edg = next(edg, g);
if(edg == edgb) {
// triangle
nm[*fb] = CGAL::unit_normal(p1, p2, p0);
} else {
// not a triangle
normal n(CGAL::NULL_VECTOR);
do {
n = n + CGAL::normal(p1, p2, p0);
p0 = p1;
p1 = p2;
edg = next(edg, g);
p2 = pm[target(edg, g)];
} while(edg != edgb);
nm[*fb] = n / CGAL::sqrt(n.squared_length());
}
}
}
int main(int argc, char** argv)
{
typedef boost::property_map<
Polyhedron,
>::const_type Face_index_map;
std::ifstream in((argc>1)?argv[1]:CGAL::data_file_path("meshes/cube_poly.off"));
Polyhedron P;
in >> P ;
// initialize facet indices
std::size_t i = 0;
for(Polyhedron::Facet_iterator it = P.facets_begin(); it != P.facets_end(); ++it, ++i)
{
it->id() = i;
}
// Ad hoc property_map to store normals. Face_index_map is used to
// map face_descriptors to a contiguous range of indices. See
// https://www.boost.org/libs/property_map/doc/vector_property_map.html
// for details.
boost::vector_property_map<Vector, Face_index_map>
normals(static_cast<unsigned>(num_faces(P)), get(CGAL::face_index, P));
calculate_face_normals(
P // Graph
, get(CGAL::vertex_point, P) // map from vertex_descriptor to point
, normals // map from face_descriptor to Vector_3
);
std::cout << "Normals" << std::endl;
for(Polyhedron::Facet_iterator it = P.facets_begin(); it != P.facets_end(); ++it) {
// Facet_iterator is a face_descriptor, so we can use it as the
// key here
std::cout << normals[it] << std::endl;
}
return 0;
}

Properties and Dynamic Properties

As the concepts HalfedgeGraph and FaceGraph add the notion of halfedges and faces, as well as a geometric embedding of the vertices, we have to add property tags such as face_index_t and vertex_point_t.

We further add dynamic properties that enable the user to add properties to vertices, halfedges, edges, and faces on the fly. The lifetime of a dynamic property is bound to the lifetime of the property map: reference counting is used to delete the property when no map refers to it.

Dynamic property tags, such as dynamic_vertex_property_t, are a generalization of boost::vertex_index_t, as they have a template parameter for the value type of the dynamic property map, and a default value. boost::property_map<G,T>::type is used to obtain the type of the dynamic property map for a graph of type G, for a dynamic property tag T. This type must be default constructible and assignable. As for ordinary properties, the function get() is overloaded and serves for retrieving a property map for a given graph and dynamic property tag, as well as for retrieving a value for a given key and property map.

The following example shows how to attach a string property to vertices and a double value to the halfedges of a graph.


File Property_map/dynamic_properties.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <string>
typedef K::Point_3 Point_3;
int main()
{
Mesh mesh;
CGAL::make_triangle(Point_3(0,0,0),Point_3(1,0,0),Point_3(1,1,0), mesh);
typedef boost::property_map<Mesh, CGAL::dynamic_vertex_property_t<std::string> >::type VertexNameMap;
VertexNameMap vnm = get(CGAL::dynamic_vertex_property_t<std::string>(), mesh);
put(vnm, *(vertices(mesh).first), "Paris");
std::cout << get(vnm, *(vertices(mesh).first)) << std::endl;
typedef boost::property_map<Mesh, CGAL::dynamic_halfedge_property_t<double> >::type TrafficDensityMap;
TrafficDensityMap tdm = get(CGAL::dynamic_halfedge_property_t<double>(), mesh);
put(tdm, *(halfedges(mesh).first), 0.7);
std::cout << get(tdm, *(halfedges(mesh).first)) << std::endl;
return 0;
}

Graph I/O

Using a common graph concept enables having common input/output functions for all the models of this concept. The following file formats are supported for models of FaceGraph:

See the page I/O Functions for an exhaustive description of the I/O functions of this package.

Euler Operations

There are two categories of mutating operations. The first category comprises low level operations that change incidences such as the target vertex of a halfedge. A halfedge graph might turn invalid by the application of inconsistent low lever operations. The second category of operations are called Euler Operations. These are high level operations such as adding a center vertex to a face, which means also adding halfedges and faces, and updating incidence information. The Euler operations enable manipulating models of MutableFaceGraph.

The complete list of Euler operations provided by this package can be found in the reference manual.

Graph Adaptors

Graph adaptors are classes that build an interface over an existing graph to provide new functionalities. By operating almost entirely on the input graph, adaptors can avoid potentially expensive operations, both in term of time and memory.

The Dual Graph

The dual graph of a FaceGraph G is a graph that has a vertex for each face of G. The dual graph has an edge whenever two faces of G are separated from each other by an edge. Thus, each edge e of G has a corresponding dual edge, the edge that connects the two faces on either side of e. Computing the dual graph of a graph has many uses, for example when one wishes to compute the connected components of a mesh.

The class template Dual is an adaptor that creates the dual view of a FaceGraph. Faces of the original graph correspond to vertices in the Dual and vice versa.

Note that border edges in a Dual have the null_face of the original graph as either source or target. This is unusual and might break other algorithms since edges are always assumed to have non-null vertices as a source and target. It is nevertheless possible to filter border edges using boost::filtered_graph, as shown in the following example.


File BGL_surface_mesh/surface_mesh_dual.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/boost/graph/Dual.h>
#include <CGAL/boost/graph/helpers.h>
#include <iostream>
#include <fstream>
#include <boost/graph/filtered_graph.hpp>
#include <boost/graph/connected_components.hpp>
typedef Kernel::Point_3 Point;
typedef CGAL::Dual<Mesh> Dual;
typedef boost::graph_traits<Dual>::edge_descriptor edge_descriptor;
template <typename G>
struct noborder {
noborder() : g(nullptr) {} // default-constructor required by filtered_graph
noborder(G& g) : g(&g) {}
bool operator()(const edge_descriptor& e) const
{ return !is_border(e,*g); }
G* g;
};
// A dual border edge has a null_face as the source or target "vertex"
// BGL algorithms won't like that, so we remove border edges through a
// boost::filtered_graph.
typedef boost::filtered_graph<Dual, noborder<Mesh> > FiniteDual;
typedef boost::graph_traits<Mesh>::vertex_descriptor vertex_descriptor;
typedef boost::graph_traits<Mesh>::face_descriptor face_descriptor;
typedef boost::graph_traits<Mesh>::edge_descriptor edge_descriptor;
int main(int argc, char* argv[])
{
const std::string filename = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/prim.off");
Mesh primal;
if(!CGAL::IO::read_polygon_mesh(filename, primal))
{
std::cerr << "Invalid input." << std::endl;
return 1;
}
Dual dual(primal);
FiniteDual finite_dual(dual,noborder<Mesh>(primal));
std::cout << "dual has " << num_vertices(dual) << " vertices" << std::endl;
std::cout << "The vertices of dual are faces in primal"<< std::endl;
for(boost::graph_traits<Dual>::vertex_descriptor dvd : vertices(dual)) {
std::cout << dvd << std::endl;
}
std::cout << "The edges in primal and dual with source and target" << std::endl;
for(edge_descriptor e : edges(dual)) {
std::cout << e << " in primal: " << source(e,primal) << " -- " << target(e,primal) << " "
<< " in dual : " << source(e,finite_dual) << " -- " << target(e,finite_dual) << std::endl;
}
std::cout << "edges of the finite dual graph" << std::endl;
for(boost::graph_traits<FiniteDual>::edge_descriptor e : CGAL::make_range(edges(finite_dual))) {
std::cout << e << " " << source(e,primal) << " " << source(e,finite_dual) << std::endl;
}
// the storage of a property map is in primal
Mesh::Property_map<face_descriptor,int> fccmap;
fccmap = primal.add_property_map<face_descriptor,int>("f:CC").first;
int num = connected_components(finite_dual, fccmap);
std::cout << "The graph has " << num << " connected components (face connectivity)" << std::endl;
for(face_descriptor f : faces(primal)) {
std::cout << f << " in connected component " << fccmap[f] << std::endl;
}
Mesh::Property_map<vertex_descriptor,int> vccmap;
vccmap = primal.add_property_map<vertex_descriptor,int>("v:CC").first;
num = connected_components(primal, vccmap);
std::cout << "The graph has " << num << " connected components (edge connectvity)" << std::endl;
for(vertex_descriptor v : vertices(primal)) {
std::cout << v << " in connected component " << vccmap[v] << std::endl;
}
return 0;
}

The Seam Mesh

The class Seam_mesh allows to mark edges of a mesh as seam edges so that they virtually become border edges when exploring a seam mesh with the BGL API. The input mesh is referred to as underlying mesh of the seam mesh. We denote tm and sm the underlying mesh and the seam mesh respectively.

Figure Figure 103.2 shows an example of mesh on which two edges, defined by the halfedge pairs h2-h3 and h6-h7, are marked as seams. The introduction of virtual borders modifies the elementary BGL graph traversal operations: when we circulate around the target of h7 in the underlying mesh, we traverse h7, h1, h3, h5, before arriving at h7 again. However, when we circulate in the seam mesh, we traverse h7, h1, h3*, before arriving at h7 again. Similarly, if we start at h3, we traverse h3, h5, h7*, and h3 again.

seam_mesh_1.svg
Figure 103.2 A seam mesh with two seam edges (h2, h3) and (h6, h7).

A vertex of the underlying mesh may correspond to multiple vertices in the seam mesh. For example in Figure Figure 103.2, the target of h7 corresponds to two vertices in the seam mesh, on either side of the virtual border created by the seam edges. For this reason, a vertex v of the seam mesh is internally represented as a halfedge h of the seam mesh. To obtain a canonical definition, the halfedge h is defined as the halfedge that has v as target, that lies on the seam, and that is not a border halfedge. The function target(hd, sm) will return this halfedge. For vertices v in the underlying mesh that are not on a seam edge, we choose halfedge(v, tm) as its canonical halfedge.

Seam Mesh Traversal

Using the function next(halfedge_descriptor, FaceGraph), we can walk around a face but also around a border of a mesh. For the seam mesh sm from Figure Figure 103.2, we have opposite(h2, sm) == h3*, and it holds that face(h3*, sm) == null_face(). We can walk along this virtual border: starting at h3* and repeatedly calling next(..,sm), we will traverse h6*, h7*, h2*, before reaching h3* again.

All other traversal functions, iterators, and circulators (see Iterators and Circulators) can be used on a seam mesh, but their behavior is similarly modified by the (virtual and real) border edges of the seam mesh.

Seams

A collection of seam edges, or simply a seam, is not necessarily a simple polyline as we can see in the next figure:

  • In (a), the seam forms a tree. Consequently, we pass at a vertex as often as there are incident seam edges.
  • In (b), the seam has a vertex v on the border of the underlying mesh. While walking along the border of the seam mesh, we leave the border of the underlying mesh when we reach v and walk on a virtual border until we reach v again, from the other side of the seam.
  • In (c), the seam forms a closed polyline. While the first two define a single border, a cycle defines two borders and splits the set of faces in two connected components. Something similar happens when the seam touches the same border more than once. A seam can also connect different borders, potentially changing the genus of the mesh. Finally, a seam may have more than one connected component.

seam_mesh_2.svg
Figure 103.3 Walking around a seam (a) with no seam vertex on the real border, (b) with a seam vertex on the real border, (c) with a closed polyline. Vertices of the seam mesh that are linked by a green dashed segment correspond to the same vertex in the underlying mesh.

Seam meshes are for example used in Chapter Triangulated Surface Mesh Parameterization to parameterize a topological sphere by first virtually cutting it into a topological disk.

Inheriting from a Model of a Face Graph Concept

In order for a type FG to be a model of one of the face graph concepts, a specialization of boost::graph_traits<FG> must exist. When you derive a class from FG this graph traits class specialization does not match. For such cases, a header file consuming some user defined macros is provided to make the derived class a valid model of the same concept. See CGAL/boost/graph/graph_traits_inheritance_macros.h for more details.

Graph Partitioning

For algorithms that operate locally, partitioning is often an easy way to parallelize the algorithm at little cost. The functions CGAL::METIS::partition_graph() and CGAL::METIS::partition_dual_graph() provide wrappers to the graph partitioning library METIS [3], allowing to split triangular meshes that are models of the concept FaceListGraph into a given number of subdomains.

The following example shows how to read, partition, and write a mesh using partition_dual_graph(). The class template CGAL::Face_filtered_graph and the free function copy_face_graph() are used to create an independent mesh from one of the subdomains of the partition. Note that the copy is optional as writing can be done directly using Face_filtered_graph.


File BGL_surface_mesh/surface_mesh_partition.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/boost/graph/Face_filtered_graph.h>
#include <CGAL/boost/graph/IO/polygon_mesh_io.h>
#include <fstream>
#include <iostream>
#include <cassert>
int main(int argc, char** argv)
{
const std::string filename = (argc>1) ? argv[1] : CGAL::data_file_path("meshes/blobby.off");
int number_of_parts = (argc>2) ? atoi(argv[2]) : 8;
SM sm;
if(!CGAL::IO::read_polygon_mesh(filename, sm))
{
std::cerr << "Invalid input." << std::endl;
return 1;
}
// The vertex <--> partition_id property map
typedef SM::Property_map<SM::Vertex_index, std::size_t> Vertex_id_map;
Vertex_id_map vertex_pid_map = sm.add_property_map<SM::Vertex_index, std::size_t>("v:pid").first;
// The face <--> partition_id property map
typedef SM::Property_map<SM::Face_index, std::size_t> Face_id_map;
Face_id_map face_pid_map = sm.add_property_map<SM::Face_index, std::size_t>("f:pid").first;
// Partition the mesh
CGAL::METIS::partition_dual_graph(sm, number_of_parts,
CGAL::parameters::vertex_partition_id_map(vertex_pid_map)
.face_partition_id_map(face_pid_map));
// Extract the part n°0 of the partition into a new, independent mesh
typedef CGAL::Face_filtered_graph<SM> Filtered_graph;
Filtered_graph filtered_sm(sm, 0 /*id of th part*/, face_pid_map);
assert(filtered_sm.is_selection_valid());
SM part_sm;
CGAL::copy_face_graph(filtered_sm, part_sm);
// Output the mesh extracted from subpart n°0
CGAL::IO::write_polygon_mesh("sm_part_0.off", part_sm, CGAL::parameters::stream_precision(17));
// Output all the vertices that are in the part n°0
std::ofstream outxyz("out.xyz");
outxyz.precision(17);
boost::graph_traits<SM>::vertex_iterator vit, ve;
boost::tie(vit, ve) = vertices(sm);
for(; vit!=ve; ++vit)
{
if(get(vertex_pid_map, *vit) == 0)
outxyz << sm.point(*vit) << std::endl;
}
return EXIT_SUCCESS;
}

Using Named Parameters some of the many options of METIS can be customized, as shown in this example.

Graph Cut

An optimal partition from a set of labels can be computed through a graph cut approach called alpha expansion [1]. CGAL provides CGAL::alpha_expansion_graphcut() which, for a graph \((V,E)\), computes the partition f that minimizes the following cost function:

\[ \mathrm{C}(f) = \sum_{\{v0,v1\} \in E} C_E(v0,v1) + \sum_{v \in V} C_V(f_v) \]

where \(C_E(v0,v1)\) is the edge cost of assigning a different label to \(v0\) and \(v1\), and \(C_V(f_v)\) is the vertex cost of assigning the label \(f\) to the vertex \(v\).

Three different implementations are provided and can be selected by using one of the following tags:

  • CGAL::Alpha_expansion_boost_adjacency_list_tag (default)
  • CGAL::Alpha_expansion_boost_compressed_sparse_raw_tag
  • CGAL::Alpha_expansion_MaxFlow_tag, released under GPL license and provided by the Triangulated Surface Mesh Segmentation Reference package

All these implementations produce the exact same result but behave differently in terms of timing and memory (see Figure 103.4). The MaxFlow implementation is the fastest, but it grows rapidly in memory when increasing the complexity of the input graph and labeling; the compressed sparse raw (CSR) is very efficient from a memory point of view but becomes very slow as the complexity of the input graph and labeling increases; the adjacency list version provides a good compromise and is therefore the default implementation.

alpha_expansion.png
Figure 103.4 Comparison of time and memory consumed by the different alpha expansion implementations.

Example

The following example shows how to apply the alpha expansion algorithm to a boost::adjacency_list describing a 2D array with 3 labels "X", " " and "O":


File BGL_graphcut/alpha_expansion_example.cpp

#include <CGAL/boost/graph/alpha_expansion_graphcut.h>
#include <boost/graph/adjacency_list.hpp>
struct Vertex_property
{
int label;
std::vector<double> cost;
};
struct Edge_property
{
double weight;
};
using Graph = boost::adjacency_list <boost::setS,
boost::vecS,
boost::undirectedS,
Vertex_property,
Edge_property>;
using GT = boost::graph_traits<Graph>;
using vertex_descriptor = GT::vertex_descriptor;
using edge_descriptor = GT::edge_descriptor;
int main()
{
std::array<char, 3> labels = { 'X', ' ', 'O' };
std::array<std::array<int, 6>, 5> input
= { { { 0, 2, 0, 1, 1, 1 },
{ 0, 0, 1, 0, 1, 2 },
{ 2, 0, 1, 1, 2, 2 },
{ 0, 1, 1, 2, 2, 0 },
{ 1, 1, 2, 0, 2, 2 } } };
std::array<std::array<vertex_descriptor, 6>, 5> vertices;
// Init vertices from values
Graph g;
for (std::size_t i = 0; i < input.size(); ++ i)
for (std::size_t j = 0; j < input[i].size(); ++ j)
{
vertices[i][j] = boost::add_vertex(g);
g[vertices[i][j]].label = input[i][j];
// Cost of assigning this vertex to any label is positive except
// for current label which is 0 (favor init solution)
g[vertices[i][j]].cost.resize(3, 1);
g[vertices[i][j]].cost[std::size_t(input[i][j])] = 0;
}
// Display input values
std::cerr << "Input:" << std::endl;
for (std::size_t i = 0; i < vertices.size(); ++ i)
{
for (std::size_t j = 0; j < vertices[i].size(); ++ j)
std::cerr << labels[std::size_t(g[vertices[i][j]].label)];
std::cerr << std::endl;
}
// Init adjacency
double weight = 0.5;
for (std::size_t i = 0; i < vertices.size(); ++ i)
for (std::size_t j = 0; j < vertices[i].size(); ++ j)
{
// Neighbor vertices are connected
if (i < vertices.size() - 1)
{
edge_descriptor ed = boost::add_edge (vertices[i][j], vertices[i+1][j], g).first;
g[ed].weight = weight;
}
if (j < vertices[i].size() - 1)
{
edge_descriptor ed = boost::add_edge (vertices[i][j], vertices[i][j+1], g).first;
g[ed].weight = weight;
}
}
std::cerr << std::endl << "Alpha expansion..." << std::endl << std::endl;
get (&Edge_property::weight, g),
get (&Vertex_property::cost, g),
get (&Vertex_property::label, g),
CGAL::parameters::vertex_index_map (get (boost::vertex_index, g)));
// Display output graph
std::cerr << "Output:" << std::endl;
for (std::size_t i = 0; i < vertices.size(); ++ i)
{
for (std::size_t j = 0; j < vertices[i].size(); ++ j)
std::cerr << labels[std::size_t(g[vertices[i][j]].label)];
std::cerr << std::endl;
}
return 0;
}

The output of this program shows how the initial 2D array is regularized spatially:

Input:
XOX
XX X O
OX OO
X OOX
OXOO
Alpha expansion...
Output:
XXX
XX O
XX OO
X OOO
OOOO

Application to Regularization of the Borders of a Face Selection

Manually selecting faces on a triangle mesh may create irregular borders (sawtooth) because of the shape of the triangles. Such borders can be regularized using the alpha expansion algorithm.

CGAL provides a function CGAL::regularize_face_selection_borders() to apply this algorithm to the borders of a face selection on a FaceGraph. Figure 103.5 shows how this function affects a selection depending on the parameters.

regularize_selection.png
Figure 103.5 Regularization of the borders of a face selection using alpha expansion. Different outputs are shown for different weight parameters, with and without preventing unselection.

The following example shows how to apply this alpha expansion regularization to the borders of a face selection of a CGAL::Surface_mesh object:


File BGL_graphcut/face_selection_borders_regularization_example.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/boost/graph/selection.h>
#include <CGAL/boost/graph/IO/OFF.h>
#include <fstream>
#include <iostream>
using Face_index = Mesh::Face_index;
int main(int argc, char** argv)
{
std::ifstream in((argc>1) ? argv[1] : CGAL::data_file_path("meshes/blobby.off"));
if(!in)
{
std::cerr << "Error: could not read input file" << std::endl;
return EXIT_FAILURE;
}
Mesh mesh;
CGAL::IO::read_OFF (in, mesh);
std::unordered_map<Face_index, bool> is_selected_map;
// randomly select 1/3 of faces
std::size_t nb_selected_before = 0;
CGAL::Random rand;
for (Face_index fi : faces(mesh))
{
bool selected = (rand.get_double() < 1. / 3.);
is_selected_map[fi] = selected;
if (selected)
nb_selected_before ++;
}
std::cerr << nb_selected_before << " selected before regularization" << std::endl;
boost::make_assoc_property_map(is_selected_map),
0.5); // using weight = 0.5
std::size_t nb_selected_after = 0;
for (const auto& sel : is_selected_map)
if (sel.second)
++ nb_selected_after;
std::cerr << nb_selected_after << " selected after regularization" << std::endl;
return EXIT_SUCCESS;
}