// Copyright (c) 2013 GeometryFactory (France). All rights reserved. // // This file is part of CGAL (www.cgal.org); you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public License as // published by the Free Software Foundation; either version 3 of the License, // or (at your option) any later version. // // Licensees holding a valid commercial license may use this file in // accordance with the commercial license agreement provided with the software. // // This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE // WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE. // // $URL$ // $Id$ // // Author(s) : Laurent Rineau, Xiang Gao // #ifndef CGAL_SPLIT_GRAPH_INTO_POLYLINES #define CGAL_SPLIT_GRAPH_INTO_POLYLINES #include #include #include #include #include #include namespace CGAL { namespace internal{ struct IsTerminalDefault { template bool operator ()(VertexDescriptor& , const Graph& ) { return false; } }; template struct Dummy_visitor_for_split_graph_into_polylines { void start_new_polyline(){} void add_node(typename boost::graph_traits::vertex_descriptor){} void end_polyline(){} }; } //end of namespace internal namespace internal { template class Less_on_G_copy_vertex_descriptors { const G_copy& g_copy; const Less_on_orig_vertex_descriptors& less; public: Less_on_G_copy_vertex_descriptors ( const G_copy& g_copy, const Less_on_orig_vertex_descriptors& less) : g_copy(g_copy), less(less) {} typedef typename boost::graph_traits::vertex_descriptor g_copy_vertex_descriptor; bool operator()(g_copy_vertex_descriptor v1, g_copy_vertex_descriptor v2) const { return less(g_copy[v1], g_copy[v2]); } }; // end class Less_on_G_copy_vertex_descriptors /// Splits a graph at vertices with degree higher than two and at vertices where `is_terminal` returns `true` /// The vertices are duplicated, and new incident edges created. /// `OrigGraph` must be undirected template void duplicate_terminal_vertices(Graph& graph, const OrigGraph& orig, IsTerminal is_terminal) { typedef typename boost::graph_traits::vertex_descriptor vertex_descriptor; typedef typename boost::graph_traits::edge_descriptor edge_descriptor; typedef typename boost::graph_traits::vertex_iterator vertex_iterator; typedef typename boost::graph_traits::out_edge_iterator out_edge_iterator; vertex_iterator b,e; boost::tie(b,e) = vertices(graph); std::vector V(b,e); BOOST_FOREACH(vertex_descriptor v, V) { if (degree(v, graph) > 2 || is_terminal(graph[v], orig)) { out_edge_iterator b, e; boost::tie(b, e) = out_edges(v, graph); std::vector E(b, e); for (unsigned int i = 1; i < E.size(); ++i) { edge_descriptor e = E[i]; vertex_descriptor w = target(e, graph); remove_edge(e, graph); vertex_descriptor vc = add_vertex(graph); graph[vc] = graph[v]; add_edge(vc, w, graph); } } } // check all vertices are of degree 1 or 2 and that the source // and target of each edge are different vertices with different ids CGAL_assertion_code( BOOST_FOREACH(vertex_descriptor v, vertices(graph)){ typename boost::graph_traits::degree_size_type n = degree(v, graph); CGAL_assertion( n == 0 || n == 1 || n == 2); } BOOST_FOREACH(edge_descriptor e, edges(graph)){ vertex_descriptor v = target(e, graph); vertex_descriptor w = source(e, graph); CGAL_assertion(v != w); } ) // end of CGAL_assertion_code } // end of duplicate_terminal_vertices } // namespace internal /*! \ingroup PkgBGL splits into polylines the graph `g` at vertices of degree greater than 2 and at vertices for which `is_terminal(v,graph)==true`. The polylines are reported using a visitor. \tparam Graph a model of the `boost` concepts `VertexListGraph` and `EdgeListGraph`. \tparam Visitor a class that provides: - void start_new_polyline() called when starting the description of a polyline. - void add_node(typename boost::graph_traits::%vertex_descriptor v) called for each vertex `v` of the polyline currently described. - void end_polyline() called when the description of a polyline is finished. \tparam IsTerminal A functor providing `bool operator()(boost::graph_traits::%vertex_descriptor v, const Graph& g) const` returning true if the vertex `v` of degree 2 is a polyline endpoint and false otherwise. An overload without `is_terminal` is provided if no vertices but those of degree different from 2 are polyline endpoints. @todo Document the version with four parameters */ template void split_graph_into_polylines(const Graph& graph, Visitor& polyline_visitor, IsTerminal is_terminal) { typedef typename boost::graph_traits::vertex_descriptor Graph_vertex_descriptor; std::less less; split_graph_into_polylines(graph, polyline_visitor, is_terminal, less); } template void split_graph_into_polylines(const Graph& graph, Visitor& polyline_visitor, IsTerminal is_terminal, LessForVertexDescriptors less) { typedef typename boost::graph_traits::vertex_descriptor Graph_vertex_descriptor; typedef typename boost::graph_traits::edge_descriptor Graph_edge_descriptor; typedef boost::adjacency_list G_copy; typedef typename boost::graph_traits::vertex_descriptor vertex_descriptor; typedef typename boost::graph_traits::edge_descriptor edge_descriptor; typedef typename boost::graph_traits::out_edge_iterator out_edge_iterator; // we make a copy of the input graph G_copy g_copy; { typedef std::map::vertex_descriptor, typename boost::graph_traits::vertex_descriptor> V2vmap; V2vmap v2vmap; BOOST_FOREACH(Graph_vertex_descriptor v, vertices(graph)){ vertex_descriptor vc = add_vertex(g_copy); g_copy[vc] = v; v2vmap[v] = vc; } BOOST_FOREACH(Graph_edge_descriptor e, edges(graph)){ Graph_vertex_descriptor vs = source(e,graph); Graph_vertex_descriptor vt = target(e,graph); vertex_descriptor vsc, vtc; if(vs == vt){ std::cerr << "ignore self loop\n"; }else{ typename V2vmap::iterator it; if((it = v2vmap.find(vs)) == v2vmap.end()){ vsc = add_vertex(g_copy); g_copy[vsc] = vs; v2vmap[vs] = vsc; }else{ vsc = it->second; } if((it = v2vmap.find(vt)) == v2vmap.end()){ vtc = add_vertex(g_copy); g_copy[vtc] = vt; v2vmap[vt] = vtc; }else{ vtc = it->second; } add_edge(vsc,vtc,g_copy); } } } // duplicate terminal vertices and vertices of degree more than 2 internal::duplicate_terminal_vertices(g_copy, graph, is_terminal); // put polylines endpoint in a set typedef internal::Less_on_G_copy_vertex_descriptors< G_copy, LessForVertexDescriptors> G_copy_less; G_copy_less g_copy_less(g_copy, less); std::set terminal(g_copy_less); BOOST_FOREACH(vertex_descriptor v, vertices(g_copy)){ typename boost::graph_traits::degree_size_type n = degree(v, g_copy); if ( n == 1 ) terminal.insert(v); if ( n ==0 ){ //isolated vertex polyline_visitor.start_new_polyline(); polyline_visitor.add_node(g_copy[v]); polyline_visitor.end_polyline(); } } // go over all polylines and provide the description to the visitor while (!terminal.empty()) { typename std::set::iterator it = terminal.begin(); vertex_descriptor u = *it; terminal.erase(it); polyline_visitor.start_new_polyline(); polyline_visitor.add_node(g_copy[u]); while (degree(u,g_copy) != 0) { CGAL_assertion(degree(u,g_copy) == 1); out_edge_iterator b = out_edges(u, g_copy).first; vertex_descriptor v = target(*b, g_copy); CGAL_assertion(u!=v); polyline_visitor.add_node(g_copy[v]); remove_edge(b, g_copy); u = v; } terminal.erase(u); polyline_visitor.end_polyline(); } // do the same but for cycles while (num_edges(g_copy) != 0) { edge_descriptor first_edge = *edges(g_copy).first; vertex_descriptor u = source(first_edge, g_copy); polyline_visitor.start_new_polyline(); polyline_visitor.add_node(g_copy[u]); u = target(first_edge, g_copy); remove_edge(first_edge, g_copy); polyline_visitor.add_node(g_copy[u]); while (degree(u,g_copy) != 0) { CGAL_assertion(degree(u,g_copy) == 1); out_edge_iterator b = out_edges(u, g_copy).first; vertex_descriptor v = target(*b, g_copy); CGAL_assertion(u!=v); polyline_visitor.add_node(g_copy[v]); remove_edge(b, g_copy); u = v; } polyline_visitor.end_polyline(); } } /// \cond SKIP_IN_MANUAL /// Split graph into polylines delimited by vertices of degree different from 2. /// Then the graph is visited and Visitor is called to describe the polylines /// Graph must be undirected template void split_graph_into_polylines(const Graph& graph, Visitor& polyline_visitor) { split_graph_into_polylines(graph, polyline_visitor, internal::IsTerminalDefault()); } /// \endcond } //end of namespace CGAL #endif //CGAL_SPLIT_GRAPH_INTO_POLYLINES