CGAL 5.3.2 - Polygon Mesh Processing

Functions to locate points on a mesh, and manipulate such locations.

Typedefs

template<typename TriangleMesh >
using CGAL::Polygon_mesh_processing::descriptor_variant = boost::variant< typename boost::graph_traits< TriangleMesh >::vertex_descriptor, typename boost::graph_traits< TriangleMesh >::halfedge_descriptor, typename boost::graph_traits< TriangleMesh >::face_descriptor >
 A variant used in the function get_descriptor_from_location(). More...
 
template<typename FT >
using CGAL::Polygon_mesh_processing::Barycentric_coordinates = std::array< FT, 3 >
 A triplet of coordinates describing the barycentric coordinates of a point with respect to the vertices of a triangular face. More...
 
template<typename TriangleMesh , typename FT >
using CGAL::Polygon_mesh_processing::Face_location = std::pair< typename boost::graph_traits< TriangleMesh >::face_descriptor, Barycentric_coordinates< FT > >
 If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following: More...
 

Functions

template<typename GeomTraits , typename Point >
std::array< typename GeomTraits::FT, 3 > CGAL::Polygon_mesh_processing::barycentric_coordinates (const Point &p, const Point &q, const Point &r, const Point &query, const GeomTraits &gt)
 Given a set of three points and a query point, computes the barycentric coordinates of the query point with respect to the first three points. More...
 
template<typename FT , typename TriangleMesh >
descriptor_variant< TriangleMesh > CGAL::Polygon_mesh_processing::get_descriptor_from_location (const Face_location< TriangleMesh, FT > &loc, const TriangleMesh &tm)
 Given a location, returns a descriptor to the simplex of smallest dimension on which the point corresponding to the location lies. More...
 
template<typename FT , typename TriangleMesh , typename NamedParameters >
Point CGAL::Polygon_mesh_processing::construct_point (const Face_location< TriangleMesh, FT > &loc, const TriangleMesh &tm, const NamedParameters &np)
 Given a location in a face, returns the geometric position described by these coordinates, as a point. More...
 

Random Location Generation

template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::random_location_on_halfedge (typename boost::graph_traits< TriangleMesh >::halfedge_descriptor hd, const TriangleMesh &tm, CGAL::Random &rnd=get_default_random())
 returns a random point over the halfedge hd, as a location. More...
 
template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::random_location_on_face (typename boost::graph_traits< TriangleMesh >::face_descriptor fd, const TriangleMesh &tm, CGAL::Random &rnd=get_default_random())
 returns a random point over the face fd, as a location. More...
 
template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::random_location_on_mesh (const TriangleMesh &tm, CGAL::Random &rnd=get_default_random())
 returns a random point over the mesh tm. More...
 

Location Predicates

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_vertex (const Face_location< TriangleMesh, FT > &loc, const typename boost::graph_traits< TriangleMesh >::vertex_descriptor vd, const TriangleMesh &tm)
 Given a location, returns whether the location is on the vertex vd or not. More...
 
template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_halfedge (const Face_location< TriangleMesh, FT > &loc, const typename boost::graph_traits< TriangleMesh >::halfedge_descriptor hd, const TriangleMesh &tm)
 Given a location, returns whether this location is on the halfedge hd or not. More...
 
template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_in_face (const Barycentric_coordinates< FT > &bar, const TriangleMesh &tm)
 Given a set of barycentric coordinates, returns whether those barycentric coordinates correspond to a point within the face (boundary included), that is, if all the barycentric coordinates are positive. More...
 
template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_in_face (const Face_location< TriangleMesh, FT > &loc, const TriangleMesh &tm)
 Given a location, returns whether the location is in the face (boundary included) or not. More...
 
template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_face_border (const Face_location< TriangleMesh, FT > &loc, const TriangleMesh &tm)
 Given a location, returns whether the location is on the boundary of the face or not. More...
 
template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_mesh_border (const Face_location< TriangleMesh, FT > &loc, const TriangleMesh &tm)
 Given a location, returns whether the location is on the border of the mesh or not. More...
 

Point Location

template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_vertex (typename boost::graph_traits< TriangleMesh >::vertex_descriptor vd, const TriangleMesh &tm)
 returns the location of the given vertex vd as a location, that is an ordered pair specifying a face incident to vd and the barycentric coordinates of the vertex vd in that face. More...
 
template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_vertex (const typename boost::graph_traits< TriangleMesh >::vertex_descriptor vd, const typename boost::graph_traits< TriangleMesh >::face_descriptor fd, const TriangleMesh &tm)
 returns the location of a given vertex as a location in fd, that is an ordered pair composed of fd and of the barycentric coordinates of the vertex in fd. More...
 
template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_on_halfedge (const typename boost::graph_traits< TriangleMesh >::halfedge_descriptor hd, const FT t, const TriangleMesh &tm)
 Given a point described by a halfedge hd and a scalar t as p = (1 - t) * source(hd, tm) + t * target(hd, tm), returns this location along the given edge as a location, that is an ordered pair specifying a face containing the location and the barycentric coordinates of that location in that face. More...
 
template<typename TriangleMesh , typename NamedParameters >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_in_face (const Point &query, const typename boost::graph_traits< TriangleMesh >::face_descriptor fd, const TriangleMesh &tm, const NamedParameters &np)
 Given a point query and a face fd of a triangulated surface mesh, returns this location as a location, that is an ordered pair composed of fd and of the barycentric coordinates of query with respect to the vertices of fd. More...
 
template<typename FT , typename TriangleMesh >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_in_adjacent_face (const Face_location< TriangleMesh, FT > &loc, const typename boost::graph_traits< TriangleMesh >::face_descriptor fd, const TriangleMesh &tm)
 Given a location and a second face adjacent to the first, returns the location of the point in the second face. More...
 

Nearest Face Location Queries

The following functions can be used to find the closest point on a triangle mesh, given either a point or a ray.

This closest point is computed using a CGAL::AABB_tree. Users intending to call location functions on more than a single point (or ray) should first compute an AABB tree to store it (otherwise, it will be recomputed every time). Note that since the AABB tree class is a 3D structure, it might be required to wrap your point property map to convert your point type to the 3D point type (i.e., your traits' Point_3) if you are working with a 2D triangle structure.

template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
void CGAL::Polygon_mesh_processing::build_AABB_tree (const TriangleMesh &tm, AABB_tree< AABB_traits< Geom_traits, CGAL::AABB_face_graph_triangle_primitive< TriangleMesh, Point3VPM > > > &outTree, const NamedParameters &np)
 creates an AABB tree suitable for use with locate_with_AABB_tree(). More...
 
template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_with_AABB_tree (const Point &p, const AABB_tree< AABB_traits< Geom_traits, AABB_face_graph_triangle_primitive< TriangleMesh, Point3VPM > > > &tree, const TriangleMesh &tm, const NamedParameters &np)
 returns the face location nearest to the given point, as a location. More...
 
template<typename TriangleMesh , typename NamedParameters >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate (const Point &p, const TriangleMesh &tm, const NamedParameters &np)
 returns the nearest face location to the given point. More...
 
template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate_with_AABB_tree (const Ray &ray, const AABB_tree< AABB_traits< Geom_traits, AABB_face_graph_triangle_primitive< TriangleMesh, Point3VPM > > > &tree, const TriangleMesh &tm, const NamedParameters &np)
 returns the face location along ray nearest to its source point. More...
 
template<typename TriangleMesh , typename NamedParameters >
Face_location< TriangleMesh, FT > CGAL::Polygon_mesh_processing::locate (const Ray &ray, const TriangleMesh &tm, const NamedParameters &np)
 returns the face location along ray nearest to its source point. More...
 

Typedef Documentation

◆ Barycentric_coordinates

template<typename FT >
using CGAL::Polygon_mesh_processing::Barycentric_coordinates = typedef std::array<FT, 3>

#include <CGAL/Polygon_mesh_processing/locate.h>

A triplet of coordinates describing the barycentric coordinates of a point with respect to the vertices of a triangular face.

See also
Face_location
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ descriptor_variant

template<typename TriangleMesh >
using CGAL::Polygon_mesh_processing::descriptor_variant = typedef boost::variant<typename boost::graph_traits<TriangleMesh>::vertex_descriptor, typename boost::graph_traits<TriangleMesh>::halfedge_descriptor, typename boost::graph_traits<TriangleMesh>::face_descriptor>

#include <CGAL/Polygon_mesh_processing/locate.h>

A variant used in the function get_descriptor_from_location().

◆ Face_location

template<typename TriangleMesh , typename FT >
using CGAL::Polygon_mesh_processing::Face_location = typedef std::pair<typename boost::graph_traits<TriangleMesh>::face_descriptor, Barycentric_coordinates<FT> >

#include <CGAL/Polygon_mesh_processing/locate.h>

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Examples:
Polygon_mesh_processing/locate_example.cpp.

Function Documentation

◆ barycentric_coordinates()

template<typename GeomTraits , typename Point >
std::array<typename GeomTraits::FT, 3> CGAL::Polygon_mesh_processing::barycentric_coordinates ( const Point &  p,
const Point &  q,
const Point &  r,
const Point &  query,
const GeomTraits &  gt 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a set of three points and a query point, computes the barycentric coordinates of the query point with respect to the first three points.

Template Parameters
GeomTraitsthe type of a geometric traits. Must be a model of Kernel and be compatible with the template parameter Point.
Pointthe type of a geometric 2D or 3D point
Parameters
p,q,rthree points with respect to whom the barycentric coordinates of query will be computed
querythe query point whose barycentric coordinates will be computed
gtan instance of the geometric traits
Precondition
p, q, and r are not collinear.
query lies on the plane defined by p, q, and r.

◆ build_AABB_tree()

template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
void CGAL::Polygon_mesh_processing::build_AABB_tree ( const TriangleMesh &  tm,
AABB_tree< AABB_traits< Geom_traits, CGAL::AABB_face_graph_triangle_primitive< TriangleMesh, Point3VPM > > > &  outTree,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

creates an AABB tree suitable for use with locate_with_AABB_tree().

This function should first be called by users who intend to locate multiple points: in this case, it is better to first build an AABB tree, and use the function locate_with_AABB_tree() that takes as parameter an AABB tree, instead of calling locate() multiple times, which will build a new AABB tree on every call.

Template Parameters
TriangleMeshmust be a model of FaceListGraph
Point3VPMmust be a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and the CGAL 3D point type (your traits' Point_3) as value type.
NamedParametersa sequence of Named Parameters
Parameters
tma triangulated surface mesh
outTreeoutput parameter that stores the computed AABB_tree
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
  • Extra: Must be identical to the traits used in the template parameter of the AABB_traits.
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ construct_point()

template<typename FT , typename TriangleMesh , typename NamedParameters >
Point CGAL::Polygon_mesh_processing::construct_point ( const Face_location< TriangleMesh, FT > &  loc,
const TriangleMesh &  tm,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location in a face, returns the geometric position described by these coordinates, as a point.

Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
NamedParametersa sequence of Named Parameters
Parameters
locthe location from which a point is constructed
tma triangulated surface mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
  • Extra: If such traits class is provided, its type FT must be identical to the template parameter FT of this function.
Precondition
loc.first is a face descriptor corresponding to a face of tm.
Returns
a point whose type is the same as the value type of the vertex point property map provided by the user or via named parameters, or the internal point map of the mesh tm.
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ get_descriptor_from_location()

template<typename FT , typename TriangleMesh >
descriptor_variant<TriangleMesh> CGAL::Polygon_mesh_processing::get_descriptor_from_location ( const Face_location< TriangleMesh, FT > &  loc,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location, returns a descriptor to the simplex of smallest dimension on which the point corresponding to the location lies.

In other words:

  • if the point lies on a vertex, this function returns a boost::graph_traits<TriangleMesh>::vertex_descriptor v;
  • if the point lies on a halfedge, this function returns a boost::graph_traits<TriangleMesh>::halfedge_descriptor hd (note that in that case, loc.first == face(hd, tm) holds).
  • otherwise, this function returns a boost::graph_traits<TriangleMesh>::face_descriptor fd (equal to loc.first).
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
loca location with loc.first a face of tm
tma triangulated surface mesh
Precondition
loc.first is a face descriptor corresponding to a face of tm.
loc describes the barycentric coordinates of a point that lives within the face (boundary included), meaning the barycentric coordinates are all positive.

◆ is_in_face() [1/2]

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_in_face ( const Barycentric_coordinates< FT > &  bar,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a set of barycentric coordinates, returns whether those barycentric coordinates correspond to a point within the face (boundary included), that is, if all the barycentric coordinates are positive.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
baran array of barycentric coordinates
tma triangulated surface mesh

◆ is_in_face() [2/2]

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_in_face ( const Face_location< TriangleMesh, FT > &  loc,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location, returns whether the location is in the face (boundary included) or not.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
loca location with loc.first a face of tm
tma triangulated surface mesh
Precondition
loc.first is a face descriptor corresponding to a face of tm.

◆ is_on_face_border()

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_face_border ( const Face_location< TriangleMesh, FT > &  loc,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location, returns whether the location is on the boundary of the face or not.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
loca location with loc.first a face of tm
tma triangulated surface mesh
Precondition
loc.first is a face descriptor corresponding to a face of tm.
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ is_on_halfedge()

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_halfedge ( const Face_location< TriangleMesh, FT > &  loc,
const typename boost::graph_traits< TriangleMesh >::halfedge_descriptor  hd,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location, returns whether this location is on the halfedge hd or not.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
loca location with loc.first a face of tm
hda halfedge of tm
tma triangulated surface mesh
Precondition
loc.first is a face descriptor corresponding to a face of tm.

◆ is_on_mesh_border()

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_mesh_border ( const Face_location< TriangleMesh, FT > &  loc,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location, returns whether the location is on the border of the mesh or not.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
loca location with loc.first a face of tm
tma triangulated surface mesh
Precondition
loc.first is a face descriptor corresponding to a face of tm.
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ is_on_vertex()

template<typename FT , typename TriangleMesh >
bool CGAL::Polygon_mesh_processing::is_on_vertex ( const Face_location< TriangleMesh, FT > &  loc,
const typename boost::graph_traits< TriangleMesh >::vertex_descriptor  vd,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location, returns whether the location is on the vertex vd or not.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
loca location with loc.first a face of tm
vda vertex of tm
tma triangulated surface mesh
Precondition
loc.first is a face descriptor corresponding to a face of tm.

◆ locate() [1/2]

template<typename TriangleMesh , typename NamedParameters >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate ( const Point &  p,
const TriangleMesh &  tm,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns the nearest face location to the given point.

Note that this function will build an AABB tree on each call. If you need to call this function more than once, first use build_AABB_tree() to create a an AABB tree that you can store and use the function locate_with_AABB_tree().

Template Parameters
TriangleMeshmust be a model of FaceListGraph.
NamedParametersa sequence of Named Parameters
Parameters
pthe point to locate on the input triangulated surface mesh
tma triangulated surface mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.

  • a tolerance value used to snap barycentric coordinates
  • Type: double
  • Default: 0
  • Extra: Depending on the geometric traits used, the computation of the barycentric coordinates might be an inexact construction, thus leading to sometimes surprising values (e.g. a triplet [0.5, 0.5, -1-e17] for a point at the middle of an edge). The coordinates will be snapped towards 0 and 1 if the difference is smaller than the tolerance value, while still ensuring that the total sum of the coordinates is 1.
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ locate() [2/2]

template<typename TriangleMesh , typename NamedParameters >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate ( const Ray &  ray,
const TriangleMesh &  tm,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns the face location along ray nearest to its source point.

If the ray does not intersect the mesh, a default constructed location is returned.

Note that this function will 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.

Template Parameters
TriangleMeshmust be a model of FaceListGraph.
NamedParametersa sequence of Named Parameters
Parameters
raya ray to intersect with the input triangulated surface mesh
tmthe input triangulated surface mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.

  • a tolerance value used to snap barycentric coordinates
  • Type: double
  • Default: 0
  • Extra: Depending on the geometric traits used, the computation of the barycentric coordinates might be an inexact construction, thus leading to sometimes surprising values (e.g. a triplet [0.5, 0.5, -1-e17] for a point at the middle of an edge). The coordinates will be snapped towards 0 and 1 if the difference is smaller than the tolerance value, while still ensuring that the total sum of the coordinates is 1.
Precondition
ray is an object with the same ambient dimension as the point type (the value type of the vertex point map).

◆ locate_in_adjacent_face()

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_in_adjacent_face ( const Face_location< TriangleMesh, FT > &  loc,
const typename boost::graph_traits< TriangleMesh >::face_descriptor  fd,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a location and a second face adjacent to the first, returns the location of the point in the second face.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
locthe first location, with loc.first being a face of tm
fdthe second face, adjacent to loc.first
tmthe triangle mesh to which fd belongs
Precondition
loc corresponds to a point that lies on a face incident to both loc.first and fd.

◆ locate_in_face()

template<typename TriangleMesh , typename NamedParameters >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_in_face ( const Point &  query,
const typename boost::graph_traits< TriangleMesh >::face_descriptor  fd,
const TriangleMesh &  tm,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a point query and a face fd of a triangulated surface mesh, returns this location as a location, that is an ordered pair composed of fd and of the barycentric coordinates of query with respect to the vertices of fd.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
TriangleMeshmust be a model of FaceGraph
NamedParametersa sequence of Named Parameters
Parameters
querya point, whose type is equal to the value type of the vertex point property map (either user-provided via named parameters or the internal point map of the mesh tm)
fda face of tm
tma triangulated surface mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
  • Extra: If such traits class is provided, its type FT must be identical to the template parameter FT of this function.

  • a tolerance value used to snap barycentric coordinates
  • Type: double
  • Default: 0
  • Extra: Depending on the geometric traits used, the computation of the barycentric coordinates might be an inexact construction, thus leading to sometimes surprising values (e.g. a triplet [0.5, 0.5, -1-e17] for a point at the middle of an edge). The coordinates will be snapped towards 0 and 1 if the difference is smaller than the tolerance value, while still ensuring that the total sum of the coordinates is 1.
Precondition
fd is not the null face
Returns
a face location. The type FT is deduced from the geometric traits, either provided by the user via named parameters (with geom_traits) or using CGAL::Kernel_traits and the point type of the vertex point property map in use.

◆ locate_on_halfedge()

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_on_halfedge ( const typename boost::graph_traits< TriangleMesh >::halfedge_descriptor  hd,
const FT  t,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

Given a point described by a halfedge hd and a scalar t as p = (1 - t) * source(hd, tm) + t * target(hd, tm), returns this location along the given edge as a location, that is an ordered pair specifying a face containing the location and the barycentric coordinates of that location in that face.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
hda halfedge of tm
tthe parametric distance of the desired point along hd
tma triangulated surface mesh

◆ locate_vertex() [1/2]

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_vertex ( typename boost::graph_traits< TriangleMesh >::vertex_descriptor  vd,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns the location of the given vertex vd as a location, that is an ordered pair specifying a face incident to vd and the barycentric coordinates of the vertex vd in that face.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
vda vertex of tm
tma triangulated surface mesh
Precondition
vd is not an isolated vertex

◆ locate_vertex() [2/2]

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_vertex ( const typename boost::graph_traits< TriangleMesh >::vertex_descriptor  vd,
const typename boost::graph_traits< TriangleMesh >::face_descriptor  fd,
const TriangleMesh &  tm 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns the location of a given vertex as a location in fd, that is an ordered pair composed of fd and of the barycentric coordinates of the vertex in fd.

If tm is the input triangulated surface mesh and given the pair (f, bc) such that bc is the triplet of barycentric coordinates (w0, w1, w2), the correspondance between the coordinates in bc and the vertices of the face f is the following:

  • w0 corresponds to source(halfedge(f, tm), tm)
  • w1 corresponds to target(halfedge(f, tm), tm)
  • w2 corresponds to target(next(halfedge(f, tm), tm), tm)
Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
vda vertex of tm and a vertex of the face fd
fda face of tm
tma triangulated surface mesh
Precondition
fd is not the null face

◆ locate_with_AABB_tree() [1/2]

template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_with_AABB_tree ( const Point &  p,
const AABB_tree< AABB_traits< Geom_traits, AABB_face_graph_triangle_primitive< TriangleMesh, Point3VPM > > > &  tree,
const TriangleMesh &  tm,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns the face location nearest to the given point, as a location.

Note that it is possible for the triangle mesh to have ambiant dimension 2 (e.g. the mesh is a 2D triangulation, or a CGAL::Surface_mesh<CGAL::Point_2<Kernel> >), as long as an appropriate vertex point property map is passed in the AABB tree, which will convert from 2D to 3D.

Template Parameters
TriangleMeshmust be a model of FaceListGraph
Point3VPMmust be a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and the CGAL 3D point type (your traits' Point_3) as value type.
NamedParametersa sequence of Named Parameters
Parameters
pthe point to locate on the input triangulated surface mesh
treean AABB tree containing the triangular faces of the input surface mesh to perform the point location with
tma triangulated surface mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
  • Extra: Must be identical to the traits used in the template parameter of the AABB_traits.

  • a tolerance value used to snap barycentric coordinates
  • Type: double
  • Default: 0
  • Extra: Depending on the geometric traits used, the computation of the barycentric coordinates might be an inexact construction, thus leading to sometimes surprising values (e.g. a triplet [0.5, 0.5, -1-e17] for a point at the middle of an edge). The coordinates will be snapped towards 0 and 1 if the difference is smaller than the tolerance value, while still ensuring that the total sum of the coordinates is 1.
Returns
a face location. The type FT is deduced from the geometric traits, either provided by the user via named parameters (with geom_traits) or using CGAL::Kernel_traits and the point type of the vertex point property map in use.
Examples:
Polygon_mesh_processing/locate_example.cpp.

◆ locate_with_AABB_tree() [2/2]

template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::locate_with_AABB_tree ( const Ray &  ray,
const AABB_tree< AABB_traits< Geom_traits, AABB_face_graph_triangle_primitive< TriangleMesh, Point3VPM > > > &  tree,
const TriangleMesh &  tm,
const NamedParameters &  np 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns the face location along ray nearest to its source point.

If the ray does not intersect the mesh, a default constructed location is returned.

Template Parameters
TriangleMeshmust be a model of FaceListGraph.
Point3VPMmust be a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and the CGAL 3D point type (your traits' Point_3) as value type.
NamedParametersa sequence of Named Parameters
Parameters
raya ray to intersect with the input triangulated surface mesh
treean AABB tree containing the triangular faces of the input surface mesh to perform the point location with
tma triangulated surface mesh
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of tm
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tm)

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.
  • Extra: Must be identical to the traits used in the template parameter of the AABB_traits.

  • a tolerance value used to snap barycentric coordinates
  • Type: double
  • Default: 0
  • Extra: Depending on the geometric traits used, the computation of the barycentric coordinates might be an inexact construction, thus leading to sometimes surprising values (e.g. a triplet [0.5, 0.5, -1-e17] for a point at the middle of an edge). The coordinates will be snapped towards 0 and 1 if the difference is smaller than the tolerance value, while still ensuring that the total sum of the coordinates is 1.
Precondition
ray is an object with the same ambient dimension as the point type (the value type of the vertex point map).

◆ random_location_on_face()

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::random_location_on_face ( typename boost::graph_traits< TriangleMesh >::face_descriptor  fd,
const TriangleMesh &  tm,
CGAL::Random rnd = get_default_random() 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns a random point over the face fd, as a location.

The random point is on the face, meaning that all its barycentric coordinates are positive. It is constructed by uniformly picking a value u between 0 and 1, a value v between 1-u, and setting the barycentric coordinates to u, v, and 1-u-v for respectively the source and target of halfedge(fd, tm), and the third point.

Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
fda face of tm
tma triangulated surface mesh
rndoptional random number generator

◆ random_location_on_halfedge()

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::random_location_on_halfedge ( typename boost::graph_traits< TriangleMesh >::halfedge_descriptor  hd,
const TriangleMesh &  tm,
CGAL::Random rnd = get_default_random() 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns a random point over the halfedge hd, as a location.

The random point is chosen on the halfedge, meaning that all its barycentric coordinates are positive. It is constructed by uniformly generating a value t between 0 and 1 and setting the barycentric coordinates to t, 1-t, and 0 for respetively the source and target of hd, and the third vertex.

Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
hda halfedge of tm
tma triangulated surface mesh
rndoptional random number generator

◆ random_location_on_mesh()

template<typename FT , typename TriangleMesh >
Face_location<TriangleMesh, FT> CGAL::Polygon_mesh_processing::random_location_on_mesh ( const TriangleMesh &  tm,
CGAL::Random rnd = get_default_random() 
)

#include <CGAL/Polygon_mesh_processing/locate.h>

returns a random point over the mesh tm.

The returned location is obtained by choosing a random face of the mesh and a random point on that face. The barycentric coordinates of the point in the face are thus all positive. Note that all faces have the same probability to be chosen.

Template Parameters
FTmust be a model of FieldNumberType
TriangleMeshmust be a model of FaceGraph
Parameters
tma triangulated surface mesh
rndoptional random number generator
See also
random_location_on_face()