CGAL 6.0 - Triangulated Surface Mesh Shortest Paths
|
#include <CGAL/Surface_mesh_shortest_path/Surface_mesh_shortest_path.h>
Computes shortest surface paths from one or more source points on a surface mesh.
Uses an optimized variation of Chen and Han's \(O(n^2)\) algorithm by Xin and Wang. Refer to those respective papers for the details of the implementation.
Traits | a model of SurfaceMeshShortestPathTraits . |
VIM | a model of ReadablePropertyMap with vertex_descriptor as key and unsigned int as value type. The default is boost::property_map<HG, boost::vertex_index_t>::const_type . |
HIM | a model of ReadablePropertyMap with halfedge_descriptor as key and unsigned int as value type. The default is boost::property_map<HG, boost::halfedge_index_t>::const_type . |
FIM | a model of ReadablePropertyMap with face_descriptor as key and unsigned int as value type. The default is boost::property_map<HG, boost::face_index_t>::const_type . |
VPM | a model of ReadablePropertyMap with vertex_descriptor as key and Traits::Point_3 as value type. The default is boost::property_map<HG, CGAL::vertex_point_t>::const_type . |
If index property maps are not provided through the constructor of the class, internal property maps must be available and initialized.
CGAL::set_halfedgeds_items_id()
Classes | |
class | Source_point_iterator |
A model of BidirectionalIterator to access the source points. More... | |
Types | |
typedef Traits::Triangle_mesh | Triangle_mesh |
The triangle mesh type which this algorithm acts on. | |
typedef boost::graph_traits< Triangle_mesh > | Graph_traits |
The BGL graph traits for this triangle mesh. | |
typedef Graph_traits::vertex_descriptor | vertex_descriptor |
Descriptors for the vertices of Triangle_mesh | |
typedef Graph_traits::halfedge_descriptor | halfedge_descriptor |
Descriptors for the halfedges of Triangle_mesh | |
typedef Graph_traits::face_descriptor | face_descriptor |
Descriptors of the faces of Triangle_mesh | |
typedef VIM | Vertex_index_map |
The vertex index property map class. | |
typedef HIM | Halfedge_index_map |
The halfedge index property map class. | |
typedef FIM | Face_index_map |
The face index property map class. | |
typedef VPM | Vertex_point_map |
The vertex point property map class. | |
typedef Traits::FT | FT |
The numeric type used by this algorithm. | |
typedef Traits::Point_3 | Point_3 |
The 3-dimensional point type, which must coincide with the value type of Vertex_point_map . | |
typedef Traits::Barycentric_coordinates | Barycentric_coordinates |
An ordered triple which specifies a location inside a triangle as a convex combination of its three vertices. | |
typedef Barycentric_coordinates | Barycentric_coordinate |
typedef std::pair< face_descriptor, Barycentric_coordinates > | Face_location |
An ordered pair specifying a location on the surface of the Triangle_mesh . | |
typedef std::pair< FT, Source_point_iterator > | Shortest_path_result |
The return type from shortest path distance queries. | |
Constructors | |
Surface_mesh_shortest_path (const Triangle_mesh &tm, const Traits &traits=Traits()) | |
creates a shortest paths object using tm as input. | |
Surface_mesh_shortest_path (const Triangle_mesh &tm, Vertex_index_map vertexIndexMap, Halfedge_index_map halfedgeIndexMap, Face_index_map faceIndexMap, Vertex_point_map vertexPointMap, const Traits &traits=Traits()) | |
creates a shortest paths object using tm as input. | |
Addition and Removal of Source Points | |
Source_point_iterator | add_source_point (vertex_descriptor v) |
adds v as a source for the shortest path queries. | |
Source_point_iterator | add_source_point (const face_descriptor f, const Barycentric_coordinates &location) |
adds a point inside the face f as a source for the shortest path queries. | |
Source_point_iterator | add_source_point (const Face_location &location) |
adds a point inside a face as a source for the shortest path queries, equivalent to Surface_mesh_shortest_path::add_source_point(location.first, location.second); | |
template<class InputIterator > | |
Source_point_iterator | add_source_points (InputIterator begin, InputIterator end) |
adds a range of points as sources for the shortest path queries. | |
void | remove_source_point (Source_point_iterator it) |
removes a source point for the shortest path queries. | |
void | remove_all_source_points () |
removes all source points for the shortest path queries. | |
Creation and Destruction of the Shortest Paths Sequence Tree | |
void | build_sequence_tree () |
Computes all pending changes to the internal sequence tree. | |
void | clear () |
removes all data, the class is as if it was constructed. | |
Accessors | |
Source_point_iterator | source_points_begin () const |
returns an iterator to the first source point location | |
Source_point_iterator | source_points_end () const |
returns an iterator to one past the last source point location | |
std::size_t | number_of_source_points () const |
returns the total number of source points used for the shortest path queries. | |
bool | changed_since_last_build () const |
determines if the internal sequence tree is valid (already built and no new source point has been added). | |
Shortest Distance Queries | |
Shortest_path_result | shortest_distance_to_source_points (const vertex_descriptor v) |
Computes the shortest surface distance from a vertex to any source point. | |
Shortest_path_result | shortest_distance_to_source_points (const face_descriptor f, const Barycentric_coordinates &location) |
Computes the shortest surface distance from any surface location to any source point. | |
Shortest Path Sequence Queries | |
template<class Visitor > | |
Shortest_path_result | shortest_path_sequence_to_source_points (const vertex_descriptor v, Visitor &visitor) |
visits the sequence of edges, vertices and faces traversed by the shortest path from a vertex to any source point. | |
template<class Visitor > | |
Shortest_path_result | shortest_path_sequence_to_source_points (const face_descriptor f, const Barycentric_coordinates &location, Visitor &visitor) |
visits the sequence of edges, vertices and faces traversed by the shortest path from any surface location to any source point. | |
Shortest Path Point Queries | |
template<class OutputIterator > | |
Shortest_path_result | shortest_path_points_to_source_points (const vertex_descriptor v, OutputIterator output) |
Computes the sequence of points in the shortest path along the surface of the input face graph from the given vertex to the closest source point. | |
template<class OutputIterator > | |
Shortest_path_result | shortest_path_points_to_source_points (const face_descriptor f, const Barycentric_coordinates &location, OutputIterator output) |
Computes the sequence of points in the shortest path along the surface of the input face graph from the given query location to the closest source point. | |
Surface Point Constructions | |
Point_3 | point (const face_descriptor f, const Barycentric_coordinates &location) const |
returns the 3-dimensional coordinates at the barycentric coordinates of the given face. | |
Point_3 | point (const halfedge_descriptor edge, const FT t) const |
returns the 3-dimensional coordinates at the parametric location along the given edge. | |
decltype(auto) | point (const vertex_descriptor v) const |
returns the 3-dimensional coordinates of the given vertex. | |
Surface Face Location Constructions | |
Face_location | face_location (const vertex_descriptor vertex) const |
returns the location of the given vertex as a Face_location | |
Face_location | face_location (const halfedge_descriptor he, const FT t) const |
returns a location along the given edge as a Face_location . | |
Nearest Face Location Queries | |
template<class AABBTraits > | |
Face_location | locate (const Point_3 &p) const |
returns the nearest face location to the given point. | |
template<class AABBTraits > | |
Face_location | locate (const Point_3 &p, const AABB_tree< AABBTraits > &tree) const |
returns the face location nearest to the given point. | |
template<class AABBTraits > | |
Face_location | locate (const Ray_3 &ray) const |
returns the face location along ray nearest to its source point. | |
template<class AABBTraits > | |
Face_location | locate (const Ray_3 &ray, const AABB_tree< AABBTraits > &tree) const |
returns the face location along ray nearest to its source point. | |
template<class AABBTraits > | |
void | build_aabb_tree (AABB_tree< AABBTraits > &outTree) const |
creates an AABB_tree suitable for use with locate . | |
typedef Barycentric_coordinates CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Barycentric_coordinate |
typedef std::pair<face_descriptor, Barycentric_coordinates> CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Face_location |
An ordered pair specifying a location on the surface of the Triangle_mesh
.
If tm
is the input graph and given the pair (f
, bc
) such that bc
is (w0, w1, w2)
, the correspondence with the weights in bc
and the vertices of the face f
is the following:
w0 = source(halfedge(f,tm),tm)
w1 = target(halfedge(f,tm),tm)
w2 = target(next(halfedge(f,tm),tm),tm)
typedef std::pair<FT, Source_point_iterator> CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Shortest_path_result |
The return type from shortest path distance queries.
Stores the distance to the nearest source point, and a Source_point_iterator
to the source point itself.
CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Surface_mesh_shortest_path | ( | const Triangle_mesh & | tm, |
const Traits & | traits = Traits() |
||
) |
creates a shortest paths object using tm
as input.
Equivalent to Surface_mesh_shortest_path(tm, get(boost::vertex_index, tm), get(boost::halfedge_index, tm), get(boost::face_index, tm), get(CGAL::vertex_point, tm), traits)
.
Internal property maps must be available and initialized.
CGAL::set_halfedgeds_items_id()
CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Surface_mesh_shortest_path | ( | const Triangle_mesh & | tm, |
Vertex_index_map | vertexIndexMap, | ||
Halfedge_index_map | halfedgeIndexMap, | ||
Face_index_map | faceIndexMap, | ||
Vertex_point_map | vertexPointMap, | ||
const Traits & | traits = Traits() |
||
) |
creates a shortest paths object using tm
as input.
No copy of the Triangle_mesh
is made, only a reference to the tm
is held.
tm | The surface mesh to compute shortest paths on. Note that it must be triangulated. |
vertexIndexMap | Property map associating an id to each vertex, from 0 to num_vertices(tm) - 1 . |
halfedgeIndexMap | Property map associating an id to each halfedge, from 0 to num_halfedges(tm) - 1 . |
faceIndexMap | Property map associating an id to each face, from 0 to num_faces(tm) - 1 . |
vertexPointMap | Property map used to access the points associated to each vertex of the graph. |
traits | Optional instance of the traits class to use. |
Source_point_iterator CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::add_source_point | ( | const face_descriptor | f, |
const Barycentric_coordinates & | location | ||
) |
adds a point inside the face f
as a source for the shortest path queries.
No change to the internal shortest paths data structure occurs until either Surface_mesh_shortest_path::build_sequence_tree()
or the first shortest path query is done.
f | A face of the input face graph |
location | Barycentric coordinates in face f specifying the source point. |
Source_point_iterator CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::add_source_point | ( | vertex_descriptor | v | ) |
adds v
as a source for the shortest path queries.
No change to the internal shortest paths data structure occurs until either Surface_mesh_shortest_path::build_sequence_tree()
or the first shortest path query is done.
Source_point_iterator CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::add_source_points | ( | InputIterator | begin, |
InputIterator | end | ||
) |
adds a range of points as sources for the shortest path queries.
No change to the internal shortest paths data structure occurs until either Surface_mesh_shortest_path::build_sequence_tree()
or the first shortest path query is done.
InputIterator | A ForwardIterator which dereferences to either Surface_mesh_shortest_path::Face_location , or Surface_mesh_shortest_path::vertex_descriptor . |
begin | iterator to the first in the list of source point locations. |
end | iterator to one past the end of the list of source point locations. |
void CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::build_aabb_tree | ( | AABB_tree< AABBTraits > & | outTree | ) | const |
creates an AABB_tree
suitable for use with locate
.
The following static overload is also available:
static void build_aabb_tree(const Triangle_mesh& tm, AABB_tree<AABBTraits>& outTree)
AABBTraits | A model of AABBTraits used to define a CGAL AABB_tree . |
outTree | Output parameter to store the computed AABB_tree |
void CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::build_sequence_tree | ( | ) |
Computes all pending changes to the internal sequence tree.
A call to this method will only trigger a computation only if some change to the set of source points occurred since the last time the sequence tree was computed.
bool CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::changed_since_last_build | ( | ) | const |
determines if the internal sequence tree is valid (already built and no new source point has been added).
void CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::clear | ( | ) |
removes all data, the class is as if it was constructed.
All internal containers are cleared and the internal sequence tree is also cleared. For a version which defers deletion until it is necessary, use Surface_mesh_shortest_path::remove_all_source_points()
.
Face_location CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::face_location | ( | const halfedge_descriptor | he, |
const FT | t | ||
) | const |
returns a location along the given edge as a Face_location
.
The following static overload is also available:
static Face_location face_location(halfedge_descriptor he, FT t, const Triangle_mesh& tm, const Traits& traits = Traits())
he | A halfedge of the input face graph |
t | Parametric distance of the desired point along he |
Face_location CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::face_location | ( | const vertex_descriptor | vertex | ) | const |
returns the location of the given vertex as a Face_location
The following static overload is also available:
static Face_location face_location(vertex_descriptor vertex, const Triangle_mesh& tm, const Traits& traits = Traits())
vertex | A vertex of the input face graph |
Face_location CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::locate | ( | const Point_3 & | p | ) | const |
returns the nearest face location to the given point.
Note that this will (re-)build an AABB_tree
on each call. If you need to call this function more than once, use build_aabb_tree()
to cache a copy of the AABB_tree
, and use the overloads of this function that accept a reference to an AABB_tree
as input.
The following static overload is also available:
static Face_location locate(const Point_3& p, const Triangle_mesh& tm, Vertex_point_map vertexPointMap, const Traits& traits = Traits())
AABBTraits | A model of AABBTraits used to define a CGAL AABB_tree . |
p | Point to locate on the input face graph |
Face_location CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::locate | ( | const Point_3 & | p, |
const AABB_tree< AABBTraits > & | tree | ||
) | const |
returns the face location nearest to the given point.
The following static overload is also available:
AABBTraits | A model of AABBTraits used to define a CGAL AABB_tree . |
p | Point to locate on the input face graph |
tree | A AABB_tree containing the triangular faces of the input surface mesh to perform the point location with |
Face_location CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::locate | ( | const Ray_3 & | ray | ) | const |
returns the face location along ray
nearest to its source point.
Note that this will (re-)build an AABB_tree
on each call. If you need to call this function more than once, use build_aabb_tree()
to cache a copy of the AABB_tree
, and use the overloads of this function that accept a reference to an AABB_tree
as input.
The following static overload is also available:
static Face_location locate(const Ray_3& ray, const Triangle_mesh& tm, Vertex_point_map vertexPointMap, const Traits& traits = Traits())
AABBTraits | A model of AABBTraits used to define an AABB_tree . |
ray | Ray to intersect with the input face graph |
Face_location CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::locate | ( | const Ray_3 & | ray, |
const AABB_tree< AABBTraits > & | tree | ||
) | const |
returns the face location along ray
nearest to its source point.
The following static overload is also available:
AABBTraits | A model of AABBTraits used to define a CGAL AABB_tree . |
ray | Ray to intersect with the input face graph |
tree | A AABB_tree containing the triangular faces of the input surface mesh to perform the point location with |
Point_3 CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::point | ( | const face_descriptor | f, |
const Barycentric_coordinates & | location | ||
) | const |
returns the 3-dimensional coordinates at the barycentric coordinates of the given face.
The following static overloads are also available:
static Point_3 point(face_descriptor f, Barycentric_coordinates location, const Triangle_mesh& tm, const Traits& traits = Traits())
static Point_3 point(face_descriptor f, Barycentric_coordinates location, const Triangle_mesh& tm, Vertex_point_map vertexPointMap, const Traits& traits = Traits())
f | A face of on the input face graph |
location | The barycentric coordinates of the query point on face f |
Point_3 CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::point | ( | const halfedge_descriptor | edge, |
const FT | t | ||
) | const |
returns the 3-dimensional coordinates at the parametric location along the given edge.
The following static overloads are also available:
static Point_3 point(halfedge_descriptor edge, FT t, const Triangle_mesh& tm, const Traits& traits = Traits())
static Point_3 point(halfedge_descriptor edge, FT t, const Triangle_mesh& tm, Vertex_point_map vertexPointMap, const Traits& traits = Traits())
edge | An edge of the input face graph |
t | The parametric distance along edge of the desired point |
decltype(auto) CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::point | ( | const vertex_descriptor | v | ) | const |
returns the 3-dimensional coordinates of the given vertex.
v | A vertex of the input face graph |
void CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::remove_all_source_points | ( | ) |
removes all source points for the shortest path queries.
No change to the internal shortest paths data structure occurs until either Surface_mesh_shortest_path::build_sequence_tree()
or the first shortest path query is done. For a version which deletes all data immediately, use clear()
instead.
void CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::remove_source_point | ( | Source_point_iterator | it | ) |
removes a source point for the shortest path queries.
No change to the internal shortest paths data structure occurs until either Surface_mesh_shortest_path::build_sequence_tree()
or the first shortest path query is done. Behaviour is undefined if the source point it
was already removed.
it | iterator to the source point to be removed |
Shortest_path_result CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::shortest_distance_to_source_points | ( | const face_descriptor | f, |
const Barycentric_coordinates & | location | ||
) |
Computes the shortest surface distance from any surface location to any source point.
f | A face of the input face graph |
location | Barycentric coordinates of the query point on face f |
source_points_end()
. Shortest_path_result CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::shortest_distance_to_source_points | ( | const vertex_descriptor | v | ) |
Computes the shortest surface distance from a vertex to any source point.
v | A vertex of the input face graph |
source_points_end()
. Shortest_path_result CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::shortest_path_points_to_source_points | ( | const face_descriptor | f, |
const Barycentric_coordinates & | location, | ||
OutputIterator | output | ||
) |
Computes the sequence of points in the shortest path along the surface of the input face graph from the given query location to the closest source point.
f | A face of on the input face graph |
location | The barycentric coordinates of the query point on face f |
output | An OutputIterator to receive the shortest path points as Point_3 objects |
source_points_end()
. Shortest_path_result CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::shortest_path_points_to_source_points | ( | const vertex_descriptor | v, |
OutputIterator | output | ||
) |
Computes the sequence of points in the shortest path along the surface of the input face graph from the given vertex to the closest source point.
v | A vertex of the input face graph |
output | An OutputIterator to receive the shortest path points as Point_3 objects |
source_points_end()
. Shortest_path_result CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::shortest_path_sequence_to_source_points | ( | const face_descriptor | f, |
const Barycentric_coordinates & | location, | ||
Visitor & | visitor | ||
) |
visits the sequence of edges, vertices and faces traversed by the shortest path from any surface location to any source point.
Visits simplices, starting from the query point, back to the nearest source point. If no shortest path could be found (for example, the surface is disconnected), then no calls to the visitor will be made (not even for the query point).
f | A face of the input face graph |
location | Barycentric coordinates of the query point on face f |
visitor | A model of SurfaceMeshShortestPathVisitor to receive the shortest path |
source_points_end()
. Shortest_path_result CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::shortest_path_sequence_to_source_points | ( | const vertex_descriptor | v, |
Visitor & | visitor | ||
) |
visits the sequence of edges, vertices and faces traversed by the shortest path from a vertex to any source point.
Visits simplices, starting from the query vertex, back to the nearest source point. If no shortest path could be found (for example, the surface is disconnected), then no calls to the visitor will be made (not even for the query vertex).
v | A vertex of the input face graph |
visitor | A model of SurfaceMeshShortestPathVisitor to receive the shortest path |
source_points_end()
. Source_point_iterator CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::source_points_begin | ( | ) | const |
returns an iterator to the first source point location
The elements will appear in the order they were inserted to the structure by calls to add_source_point()
or add_source_points()
. Deleted points will not appear in the sequence.
Source_point_iterator CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::source_points_end | ( | ) | const |
returns an iterator to one past the last source point location