Functions to locate points on a mesh, and manipulate such locations.
|
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 >) |
| 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...
|
|
|
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...
|
|
|
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...
|
|
|
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...
|
|
|
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...
|
|
template<typename TriangleMesh , typename 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.
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
-
GeomTraits | the type of a geometric traits. Must be a model of Kernel and be compatible with the template parameter Point . |
Point | the type of a geometric 2D or 3D point |
- Parameters
-
p,q,r | three points with respect to whom the barycentric coordinates of query will be computed |
query | the query point whose barycentric coordinates will be computed |
gt | an instance of the geometric traits |
- Precondition
p
, q
, and r
are not collinear.
-
query
lies on the plane defined by p
, q
, and r
.
template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
#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
-
TriangleMesh | must be a model of FaceListGraph |
Point3VPM | must 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. |
NamedParameters | a sequence of Named Parameters |
- Parameters
-
tm | a triangulated surface mesh |
outTree | output parameter that stores the computed AABB_tree |
np | an optional sequence of Named Parameters among the ones listed below: |
- Named Parameters
vertex_point_map | the property map with the points associated to the vertices of tm . If this parameter is omitted, an internal property map for boost::vertex_point_t must be available in TriangleMesh . |
geom_traits | a geometric traits class instance, model of Kernel compatible with the point type held in the vertex point property map (either user-provided or internal to the mesh). Must be identical to the traits used in the template parameter of the AABB_traits . |
- Examples
- Polygon_mesh_processing/locate_example.cpp.
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
-
- Parameters
-
bar | an array of barycentric coordinates |
tm | a triangulated surface mesh |
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
-
- Parameters
-
loc | a location with loc.first a face of tm |
tm | a triangulated surface mesh |
- Precondition
loc.first
is a face descriptor corresponding to a face of tm
.
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
-
- Parameters
-
loc | a location with loc.first a face of tm |
tm | a triangulated surface mesh |
- Precondition
loc.first
is a face descriptor corresponding to a face of tm
.
- Examples
- Polygon_mesh_processing/locate_example.cpp.
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
-
- Parameters
-
loc | a location with loc.first a face of tm |
hd | a halfedge of tm |
tm | a triangulated surface mesh |
- Precondition
loc.first
is a face descriptor corresponding to a face of tm
.
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
-
- Parameters
-
loc | a location with loc.first a face of tm |
tm | a triangulated surface mesh |
- Precondition
loc.first
is a face descriptor corresponding to a face of tm
.
- Examples
- Polygon_mesh_processing/locate_example.cpp.
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
-
- Parameters
-
loc | a location with loc.first a face of tm |
vd | a vertex of tm |
tm | a triangulated surface mesh |
- Precondition
loc.first
is a face descriptor corresponding to a face of tm
.
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
-
- Parameters
-
p | the point to locate on the input triangulated surface mesh |
tm | a triangulated surface mesh |
np | an optional sequence of Named Parameters among the ones listed below: |
- Named Parameters
vertex_point_map | the property map with the points associated to the vertices of tm . If this parameter is omitted, an internal property map for boost::vertex_point_t must be available in TriangleMesh . |
geom_traits | a geometric traits class instance, model of Kernel compatible with the point type held in the vertex point property map (either user-provided or internal to the mesh). |
snapping_tolerance | a tolerance value used to snap barycentric coordinates. 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 . By default, the tolerance is 0 . |
- Examples
- Polygon_mesh_processing/locate_example.cpp.
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
-
- Parameters
-
ray | a ray to intersect with the input triangulated surface mesh |
tm | the input triangulated surface mesh |
np | an optional sequence of Named Parameters among the ones listed below: |
- Named Parameters
vertex_point_map | the property map with the points associated to the vertices of tm . If this parameter is omitted, an internal property map for boost::vertex_point_t must be available in TriangleMesh . |
geom_traits | a geometric traits class instance, model of Kernel compatible with the point type held in the vertex point property map (either user-provided or internal to the mesh). |
snapping_tolerance | a tolerance value used to snap barycentric coordinates. 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 . By default, the tolerance is 0 . |
- Precondition
ray
is an object with the same ambient dimension as the point type (the value type of the vertex point map).
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
-
- Parameters
-
loc | the first location, with loc.first being a face of tm |
fd | the second face, adjacent to loc.first |
tm | the triangle mesh to which fd belongs |
- Precondition
loc
corresponds to a point that lies on a face incident to both loc.first
and fd
.
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
-
- Parameters
-
query | a 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 ) |
fd | a face of tm |
tm | a triangulated surface mesh |
np | an optional sequence of Named Parameters among the ones listed below: |
- Named Parameters
vertex_point_map | the property map with the points associated to the vertices of tm . If this parameter is omitted, an internal property map for boost::vertex_point_t must be available in TriangleMesh . |
geom_traits | a geometric traits class instance, model of Kernel . If provided, the types FT and Kernel::FT must be identical and the traits class must be compatible with the value type of the vertex point property map. |
snapping_tolerance | a tolerance value used to snap barycentric coordinates. 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 . By default, the tolerance is 0 . |
- 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.
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
-
- Parameters
-
hd | a halfedge of tm |
t | the parametric distance of the desired point along hd |
tm | a triangulated surface mesh |
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
-
- Parameters
-
vd | a vertex of tm and a vertex of the face fd |
fd | a face of tm |
tm | a triangulated surface mesh |
- Precondition
fd
is not the null face
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
-
- Parameters
-
vd | a vertex of tm |
tm | a triangulated surface mesh |
- Precondition
vd
is not an isolated vertex
template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
#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
-
TriangleMesh | must be a model of FaceListGraph |
Point3VPM | must 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. |
NamedParameters | a sequence of Named Parameters |
- Parameters
-
p | the point to locate on the input triangulated surface mesh |
tree | an AABB tree containing the triangular faces of the input surface mesh to perform the point location with |
tm | a triangulated surface mesh |
np | an optional sequence of Named Parameters among the ones listed below: |
- Named Parameters
vertex_point_map | the property map with the points associated to the vertices of tm . If this parameter is omitted, an internal property map for boost::vertex_point_t must be available in TriangleMesh . |
geom_traits | a geometric traits class instance, model of Kernel compatible with the point type held in the vertex point property map (either user-provided or internal to the mesh). Must be identical to the traits used in the template parameter of the AABB_traits . |
snapping_tolerance | a tolerance value used to snap barycentric coordinates. 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 . By default, the tolerance is 0 . |
- 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.
template<typename TriangleMesh , typename Point3VPM , typename NamedParameters >
#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
-
TriangleMesh | must be a model of FaceListGraph . |
Point3VPM | must 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. |
NamedParameters | a sequence of Named Parameters |
- Parameters
-
ray | a ray to intersect with the input triangulated surface mesh |
tree | an AABB tree containing the triangular faces of the input surface mesh to perform the point location with |
tm | a triangulated surface mesh |
np | an optional sequence of Named Parameters among the ones listed below: |
- Named Parameters
vertex_point_map | the property map with the points associated to the vertices of tm . If this parameter is omitted, an internal property map for boost::vertex_point_t must be available in TriangleMesh . |
geom_traits | a geometric traits class instance, model of Kernel compatible with the point type held in the vertex point property map (either user-provided or internal to the mesh). Must be identical to the traits used in the template parameter of the AABB_traits . |
snapping_tolerance | a tolerance value used to snap barycentric coordinates. 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 . By default, the tolerance is 0 . |
- Precondition
ray
is an object with the same ambient dimension as the point type (the value type of the vertex point map).
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
-
- Parameters
-
fd | a face of tm |
tm | a triangulated surface mesh |
rnd | optional random number generator |
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
-
- Parameters
-
hd | a halfedge of tm |
tm | a triangulated surface mesh |
rnd | optional random number generator |