CGAL 5.2.1 - Triangulated Surface Mesh Shortest Paths
CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM > Class Template Reference

#include <CGAL/Surface_mesh_shortest_path/Surface_mesh_shortest_path.h>

Definition

template<class Traits, class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
class CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >

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.

Template Parameters
Traitsa model of SurfaceMeshShortestPathTraits.
VIMa 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.
HIMa 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.
FIMa 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.
VPMa 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.

See also
CGAL::set_halfedgeds_items_id()
Examples:
Surface_mesh_shortest_path/shortest_path_sequence.cpp, Surface_mesh_shortest_path/shortest_path_with_locate.cpp, Surface_mesh_shortest_path/shortest_paths.cpp, Surface_mesh_shortest_path/shortest_paths_multiple_sources.cpp, Surface_mesh_shortest_path/shortest_paths_no_id.cpp, Surface_mesh_shortest_path/shortest_paths_OpenMesh.cpp, and Surface_mesh_shortest_path/shortest_paths_with_id.cpp.

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_meshGraph_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. More...
 
typedef Barycentric_coordinates Barycentric_coordinate
 
typedef std::pair< face_descriptor, Barycentric_coordinatesFace_location
 An ordered pair specifying a location on the surface of the Triangle_mesh. More...
 
typedef std::pair< FT, Source_point_iteratorShortest_path_result
 The return type from shortest path distance queries. More...
 

Constructors

 Surface_mesh_shortest_path (const Triangle_mesh &tm, const Traits &traits=Traits())
 Creates a shortest paths object using tm as input. More...
 
 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. More...
 

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. More...
 
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. More...
 
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. More...
 
void remove_source_point (Source_point_iterator it)
 Removes a source point for the shortest path queries. More...
 
void remove_all_source_points ()
 Removes all source points for the shortest path queries. More...
 

Creation and Destruction of the Shortest Paths Sequence Tree

void build_sequence_tree ()
 Computes all pending changes to the internal sequence tree. More...
 
void clear ()
 Removes all data, the class is as if it was constructed. More...
 

Accessors

Source_point_iterator source_points_begin () const
 Returns an iterator to the first source point location. More...
 
Source_point_iterator source_points_end () const
 Returns an iterator to one past the last source point location. More...
 
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). More...
 

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. More...
 
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. More...
 

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. More...
 
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. More...
 

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. More...
 
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. More...
 

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. More...
 
Point_3 point (const halfedge_descriptor edge, const FT t) const
 Returns the 3-dimensional coordinates at the parametric location along the given edge. More...
 
Point_3 point (const vertex_descriptor vertex) const
 Returns the 3-dimensional coordinates of the given vertex. More...
 

Surface Face Location Constructions

Face_location face_location (const vertex_descriptor vertex) const
 Returns the location of the given vertex as a Face_location More...
 
Face_location face_location (const halfedge_descriptor he, const FT t) const
 Returns a location along the given edge as a Face_location. More...
 

Nearest Face Location Queries

template<class AABBTraits >
Face_location locate (const Point_3 &p) const
 Returns the nearest face location to the given point. More...
 
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. More...
 
template<class AABBTraits >
Face_location locate (const Ray_3 &ray) const
 Returns the face location along ray nearest to its source point. More...
 
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. More...
 
template<class AABBTraits >
void build_aabb_tree (AABB_tree< AABBTraits > &outTree) const
 Creates an AABB_tree suitable for use with locate. More...
 

Member Typedef Documentation

◆ Barycentric_coordinate

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
typedef Barycentric_coordinates CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Barycentric_coordinate

◆ Barycentric_coordinates

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
typedef Traits::Barycentric_coordinates CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::Barycentric_coordinates

An ordered triple which specifies a location inside a triangle as a convex combination of its three vertices.

◆ Face_location

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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 correspondance 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)

◆ Shortest_path_result

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Constructor & Destructor Documentation

◆ Surface_mesh_shortest_path() [1/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

See also
CGAL::set_halfedgeds_items_id()

◆ Surface_mesh_shortest_path() [2/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Parameters
tmThe surface mesh to compute shortest paths on. Note that it must be triangulated.
vertexIndexMapProperty map associating an id to each vertex, from 0 to num_vertices(tm) - 1.
halfedgeIndexMapProperty map associating an id to each halfedge, from 0 to num_halfedges(tm) - 1.
faceIndexMapProperty map associating an id to each face, from 0 to num_faces(tm) - 1.
vertexPointMapProperty map used to access the points associated to each vertex of the graph.
traitsOptional instance of the traits class to use.

Member Function Documentation

◆ add_source_point() [1/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Returns
An iterator to the source point added

◆ add_source_point() [2/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Parameters
fA face of the input face graph
locationBarycentric coordinates in face f specifying the source point.
Returns
An iterator to the source point added

◆ add_source_points()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class InputIterator >
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.

Template Parameters
InputIteratorA ForwardIterator which dereferences to either Surface_mesh_shortest_path::Face_location, or Surface_mesh_shortest_path::vertex_descriptor.
Parameters
beginiterator to the first in the list of source point locations.
enditerator to one past the end of the list of source point locations.
Returns
An iterator to the first source point added.

◆ build_aabb_tree()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class AABBTraits >
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:

Template Parameters
AABBTraitsA model of AABBTraits used to define a CGAL AABB_tree.
Parameters
outTreeOutput parameter to store the computed AABB_tree

◆ build_sequence_tree()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

◆ changed_since_last_build()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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).

Returns
true if the structure needs to be rebuilt, false otherwise

◆ clear()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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() [1/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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())
Parameters
vertexA vertex of the input face graph

◆ face_location() [2/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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())
Parameters
heA halfedge of the input face graph
tParametric distance of the desired point along he

◆ locate() [1/4]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class AABBTraits >
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())
Template Parameters
AABBTraitsA model of AABBTraits used to define a CGAL AABB_tree.
Parameters
pPoint to locate on the input face graph

◆ locate() [2/4]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class AABBTraits >
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:

  • static Face_location locate(const Point_3& p, const AABB_tree<AABBTraits>& tree, const Triangle_mesh& tm, Vertex_point_map vertexPointMap, const Traits& traits = Traits())
Template Parameters
AABBTraitsA model of AABBTraits used to define a CGAL AABB_tree.
Parameters
pPoint to locate on the input face graph
treeA AABB_tree containing the triangular faces of the input surface mesh to perform the point location with

◆ locate() [3/4]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class AABBTraits >
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())
Template Parameters
AABBTraitsA model of AABBTraits used to define an AABB_tree.
Parameters
rayRay to intersect with the input face graph

◆ locate() [4/4]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class AABBTraits >
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:

  • static Face_location locate(const Ray_3& ray, const AABB_tree<AABBTraits>& tree, const Triangle_mesh& tm, Vertex_point_map vertexPointMap, const Traits& traits = Traits())
Template Parameters
AABBTraitsA model of AABBTraits used to define a CGAL AABB_tree.
Parameters
rayRay to intersect with the input face graph
treeA AABB_tree containing the triangular faces of the input surface mesh to perform the point location with

◆ point() [1/3]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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())
Parameters
fA face of on the input face graph
locationThe barycentric coordinates of the query point on face f

◆ point() [2/3]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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())
Parameters
edgeAn edge of the input face graph
tThe parametric distance along edge of the desired point

◆ point() [3/3]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
Point_3 CGAL::Surface_mesh_shortest_path< Traits, VIM, HIM, FIM, VPM >::point ( const vertex_descriptor  vertex) const

Returns the 3-dimensional coordinates of the given vertex.

Parameters
vertexA vertex of the input face graph

◆ remove_all_source_points()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

◆ remove_source_point()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Parameters
ititerator to the source point to be removed

◆ shortest_distance_to_source_points() [1/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Parameters
vA vertex of the input face graph
Returns
A pair, containing the distance to the source point, and an iterator to the source point. If no source point was reachable (can occur when the graph is disconnected), the distance will be a negative value and the source point iterator will be equal to source_points_end().

◆ shortest_distance_to_source_points() [2/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Parameters
fA face of the input face graph
locationBarycentric coordinates of the query point on face f
Returns
A pair, containing the distance to the source point, and an iterator to the source point. If no source point was reachable (can occur when the graph is disconnected), the distance will be a negative value and the source point iterator will be equal to source_points_end().

◆ shortest_path_points_to_source_points() [1/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class OutputIterator >
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.

Parameters
vA vertex of the input face graph
outputAn OutputIterator to receive the shortest path points as Point_3 objects
Returns
A pair, containing the distance to the source point, and an iterator to the source point. If no source point was reachable (can occur when the graph is disconnected), the distance will be a negative value and the source point iterator will be equal to source_points_end().

◆ shortest_path_points_to_source_points() [2/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class OutputIterator >
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.

Parameters
fA face of on the input face graph
locationThe barycentric coordinates of the query point on face f
outputAn OutputIterator to receive the shortest path points as Point_3 objects
Returns
A pair, containing the distance to the source point, and an iterator to the source point. If no source point was reachable (can occur when the graph is disconnected), the distance will be a negative value and the source point iterator will be equal to source_points_end().

◆ shortest_path_sequence_to_source_points() [1/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class Visitor >
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).

Parameters
vA vertex of the input face graph
visitorA model of SurfaceMeshShortestPathVisitor to receive the shortest path
Returns
A pair, containing the distance to the source point, and an iterator to the source point. If no source point was reachable (can occur when the graph is disconnected), the distance will be a negative value and the source point iterator will be equal to source_points_end().

◆ shortest_path_sequence_to_source_points() [2/2]

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
template<class Visitor >
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).

Parameters
fA face of the input face graph
locationBarycentric coordinates of the query point on face f
visitorA model of SurfaceMeshShortestPathVisitor to receive the shortest path
Returns
A pair, containing the distance to the source point, and an iterator to the source point. If no source point was reachable (can occur when the graph is disconnected), the distance will be a negative value and the source point iterator will be equal to source_points_end().

◆ source_points_begin()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Returns
An iterator to the first of the stored source points.

◆ source_points_end()

template<class Traits , class VIM = Default, class HIM = Default, class FIM = Default, class VPM = Default>
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.

Returns
An iterator to one past-the-end in the list of stored source points.