CGAL 5.1 - 3D Generalized Barycentric Coordinates

Free functions to compute barycentric coordinates.

Functions

template<typename TriangleMesh , typename OutIterator , typename GeomTraits , typename VertexToPointMap >
std::pair< OutIterator, bool > CGAL::Barycentric_coordinates::boundary_coordinates_3 (const TriangleMesh &triangle_mesh, const typename GeomTraits::Point_3 &query, OutIterator c_begin, const GeomTraits &traits, const VertexToPointMap vertex_to_point_map)
 computes boundary barycentric coordinates. More...
 
template<typename TriangleMesh , typename Point_3 , typename OutIterator , typename VertexToPointMap = typename property_map_selector<TriangleMesh, CGAL::vertex_point_t>::const_type>
std::pair< OutIterator, bool > CGAL::Barycentric_coordinates::boundary_coordinates_3 (const TriangleMesh &triangle_mesh, const Point_3 &query, OutIterator c_begin, const VertexToPointMap vertex_to_point_map)
 computes boundary barycentric coordinates. More...
 
template<typename Point_3 , typename TriangleMesh , typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::discrete_harmonic_coordinates_3 (const TriangleMesh &triangle_mesh, const Point_3 &query, OutIterator c_begin, const Computation_policy_3 policy=Computation_policy_3::FAST)
 computes 3D discrete harmonic coordinates. More...
 
template<typename Point_3 , typename TriangleMesh , typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::mean_value_coordinates_3 (const TriangleMesh &triangle_mesh, const Point_3 &query, OutIterator c_begin, const Computation_policy_3 policy=Computation_policy_3::FAST_WITH_EDGE_CASES)
 computes 3D mean value coordinates. More...
 
template<typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::tetrahedron_coordinates (const typename GeomTraits::Point_3 &p0, const typename GeomTraits::Point_3 &p1, const typename GeomTraits::Point_3 &p2, const typename GeomTraits::Point_3 &p3, const typename GeomTraits::Point_3 &query, OutIterator c_begin, const GeomTraits &traits)
 computes tetrahedron coordinates. More...
 
template<typename GeomTraits >
std::array< typename GeomTraits::FT, 4 > CGAL::Barycentric_coordinates::tetrahedron_coordinates_in_array (const typename GeomTraits::Point_3 &p0, const typename GeomTraits::Point_3 &p1, const typename GeomTraits::Point_3 &p2, const typename GeomTraits::Point_3 &p3, const typename GeomTraits::Point_3 &query, const GeomTraits &traits)
 computes tetrahedron coordinates. More...
 
template<typename Point_3 , typename PolygonMesh , typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::wachspress_coordinates_3 (const PolygonMesh &polygon_mesh, const Point_3 &query, OutIterator c_begin, const Computation_policy_3 policy=Computation_policy_3::FAST)
 computes 3D Wachspress coordinates. More...
 

Function Documentation

◆ boundary_coordinates_3() [1/2]

template<typename TriangleMesh , typename Point_3 , typename OutIterator , typename VertexToPointMap = typename property_map_selector<TriangleMesh, CGAL::vertex_point_t>::const_type>
std::pair<OutIterator, bool> CGAL::Barycentric_coordinates::boundary_coordinates_3 ( const TriangleMesh &  triangle_mesh,
const Point_3 query,
OutIterator  c_begin,
const VertexToPointMap  vertex_to_point_map 
)

#include <CGAL/Barycentric_coordinates_3/boundary_coordinates_3.h>

computes boundary barycentric coordinates.

This function computes boundary barycentric coordinates at a given query point with respect to the vertices of a simple polyhedron, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

If query is at the vertex, the corresponding coordinate is set to one, while all other coordinates are zero. If query is on the face, the three corresponding coordinates are triangle coordinates, while all other coordinates are set to zero. If query is not on the boundary, all the coordinates are set to zero.

Template Parameters
TriangleMeshmust be a model of the concept FaceListGraph.
Point_3a model of Kernel::Point_3
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
VertexToPointMapa property map with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type. The default is property_map_selector<TriangleMesh, CGAL::vertex_point_t>.
Parameters
triangle_meshan instance of TriangleMesh, which must be a convex simplicial polyhedron
querya query point
c_beginthe beginning of the destination range with the computed coordinates
vertex_to_point_mapan instance of VertexToPointMap that maps a vertex from triangle_mesh to Point_3; the default initialization is provided
Returns
an output iterator to the element in the destination range, one past the last coordinate stored + the flag indicating whether the query point belongs to the polyhedron boundary
Precondition
num_vertices(triangle_mesh) >= 4.
triangle_mesh is simplicial.

◆ boundary_coordinates_3() [2/2]

template<typename TriangleMesh , typename OutIterator , typename GeomTraits , typename VertexToPointMap >
std::pair<OutIterator, bool> CGAL::Barycentric_coordinates::boundary_coordinates_3 ( const TriangleMesh &  triangle_mesh,
const typename GeomTraits::Point_3 &  query,
OutIterator  c_begin,
const GeomTraits &  traits,
const VertexToPointMap  vertex_to_point_map 
)

#include <CGAL/Barycentric_coordinates_3/boundary_coordinates_3.h>

computes boundary barycentric coordinates.

This function computes boundary barycentric coordinates at a given query point with respect to the vertices of a simple polyhedron, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

If query is at the vertex, the corresponding coordinate is set to one, while all other coordinates are zero. If query is on the face, the three corresponding coordinates are triangle coordinates, while all other coordinates are set to zero. If query is not on the boundary, all the coordinates are set to zero.

Template Parameters
TriangleMeshmust be a model of the concept FaceListGraph.
GeomTraitsa model of BarycentricTraits_3
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
VertexToPointMapa property map with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point_3 as value type. The default is property_map_selector<TriangleMesh, CGAL::vertex_point_t>.
Parameters
triangle_meshan instance of TriangleMesh, which must be a convex simplicial polyhedron
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; the default initialization is provided
vertex_to_point_mapan instance of VertexToPointMap that maps a vertex from triangle_mesh to Point_3; the default initialization is provided
Returns
an output iterator to the element in the destination range, one past the last coordinate stored + the flag indicating whether the query point belongs to the polyhedron boundary
Precondition
num_vertices(triangle_mesh) >= 4.
triangle_mesh is simplicial.

◆ discrete_harmonic_coordinates_3()

template<typename Point_3 , typename TriangleMesh , typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::discrete_harmonic_coordinates_3 ( const TriangleMesh &  triangle_mesh,
const Point_3 query,
OutIterator  c_begin,
const Computation_policy_3  policy = Computation_policy_3::FAST 
)

#include <CGAL/Barycentric_coordinates_3/Discrete_harmonic_coordinates_3.h>

computes 3D discrete harmonic coordinates.

This function computes 3D discrete harmonic coordinates at a given query point with respect to the vertices of a convex polyhedron with triangular faces, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

Internally, the class Discrete_harmonic_coordinates_3 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
Point_3A model of Kernel::Point_3.
TriangleMeshmust be a model of the concept FaceListGraph.
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
Parameters
triangle_meshan instance of TriangleMesh, which must be a convex simplicial polyhedron
querya query point
c_beginthe beginning of the destination range with the computed coordinates
policyone of the CGAL::Barycentric_coordinates::Computation_policy_3; the default is Computation_policy_3::FAST_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
num_vertices(triangle_mesh) >= 4.
triangle_mesh is strongly convex.
triangle_mesh is simplicial.
Examples
Barycentric_coordinates_3/discrete_harmonic_coordinates.cpp.

◆ mean_value_coordinates_3()

template<typename Point_3 , typename TriangleMesh , typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::mean_value_coordinates_3 ( const TriangleMesh &  triangle_mesh,
const Point_3 query,
OutIterator  c_begin,
const Computation_policy_3  policy = Computation_policy_3::FAST_WITH_EDGE_CASES 
)

#include <CGAL/Barycentric_coordinates_3/Mean_value_coordinates_3.h>

computes 3D mean value coordinates.

This function computes 3D mean value coordinates at a given query point with respect to the vertices of a polyhedron with triangular faces, that is one weight per vertex. The coordinates are stored in a destination range beginning at c_begin.

Internally, the class Mean_value_coordinates_3 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
Point_3A model of Kernel::Point_3.
TriangleMeshmust be a model of the concept FaceListGraph.
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
Parameters
triangle_meshan instance of TriangleMesh, which must be a simplicial polyhedron
querya query point
c_beginthe beginning of the destination range with the computed coordinates
policyone of the CGAL::Barycentric_coordinates::Computation_policy_3; the default is Computation_policy_3::FAST_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last coordinates stored
Precondition
num_vertices(triangle_mesh) >= 4.
triangle_mesh is simplicial.
Examples
Barycentric_coordinates_3/mean_value_coordinates.cpp.

◆ tetrahedron_coordinates()

template<typename OutIterator , typename GeomTraits >
OutIterator CGAL::Barycentric_coordinates::tetrahedron_coordinates ( const typename GeomTraits::Point_3 &  p0,
const typename GeomTraits::Point_3 &  p1,
const typename GeomTraits::Point_3 &  p2,
const typename GeomTraits::Point_3 &  p3,
const typename GeomTraits::Point_3 &  query,
OutIterator  c_begin,
const GeomTraits &  traits 
)

#include <CGAL/Barycentric_coordinates_3/tetrahedron_coordinates.h>

computes tetrahedron coordinates.

This function computes barycentric coordinates at a given query point with respect to the points p0, p1, p2, and p3, which form a tetrahedron, that is one coordinate per point. The coordinates are stored in a destination range beginning at c_begin.

After the coordinates \(b_0\), \(b_1\), \(b_2\), and \(b_2\) are computed, the query point \(q\) can be obtained as \(q = b_0p_0 + b_1p_1 + b_2p_2 + b_3p_3\).

Template Parameters
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
GeomTraitsa model of BarycentricTraits_3
Parameters
p0the first vertex of a tetrahedron
p1the second vertex of a tetrahedron
p2the third vertex of a tetrahedron
p3the fourth vertex of a tetrahedron
querya query point
c_beginthe beginning of the destination range with the computed coordinates
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
Compute_volume_3(p0, p1, p2, p3) != 0

◆ tetrahedron_coordinates_in_array()

template<typename GeomTraits >
std::array<typename GeomTraits::FT, 4> CGAL::Barycentric_coordinates::tetrahedron_coordinates_in_array ( const typename GeomTraits::Point_3 &  p0,
const typename GeomTraits::Point_3 &  p1,
const typename GeomTraits::Point_3 &  p2,
const typename GeomTraits::Point_3 &  p3,
const typename GeomTraits::Point_3 &  query,
const GeomTraits &  traits 
)

#include <CGAL/Barycentric_coordinates_3/tetrahedron_coordinates.h>

computes tetrahedron coordinates.

This function computes barycentric coordinates at a given query point with respect to the points p0, p1, p2, and p3, which form a tetrahedron, that is one coordinate per point. The coordinates are returned in a tuple.

After the coordinates \(b_0\), \(b_1\), \(b_2\), and \(b_3\) are computed, the query point \(q\) can be obtained as \(q = b_0p_0 + b_1p_1 + b_2p_2 + b_3p_3\).

Template Parameters
GeomTraitsa model of BarycentricTraits_3
Parameters
p0the first vertex of a tetrahedron
p1the second vertex of a tetrahedron
p2the third vertex of a tetrahedron
p3the fourth vertex of a tetrahedron
querya query point
traitsa traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type
Returns
a array std::array<GeomTraits::FT, 4> with the computed coordinates
Precondition
Compute_volume_3(p0, p1, p2, p3) != 0
Examples
Barycentric_coordinates_3/tetrahedron_coordinates.cpp.

◆ wachspress_coordinates_3()

template<typename Point_3 , typename PolygonMesh , typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::wachspress_coordinates_3 ( const PolygonMesh &  polygon_mesh,
const Point_3 query,
OutIterator  c_begin,
const Computation_policy_3  policy = Computation_policy_3::FAST 
)

#include <CGAL/Barycentric_coordinates_3/Wachspress_coordinates_3.h>

computes 3D Wachspress coordinates.

This function computes 3D Wachspress coordinates at a given query point with respect to the vertices of a convex polyhedron with triangular faces, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at c_begin.

Internally, the class Wachspress_coordinates_3 is used. If one wants to process multiple query points, it is better to use that class. When using the free function, internal memory is allocated for each query point, while when using the class, it is allocated only once, which is much more efficient. However, for a few query points, it is easier to use this function. It can also be used when the processing time is not a concern.

Template Parameters
Point_3A model of Kernel::Point_3.
PolygonMeshmust be a model of the concept FaceListGraph.
OutIteratora model of OutputIterator that accepts values of type GeomTraits::FT
Parameters
polygon_meshan instance of PolygonMesh, which must be a convex simplicial polyhedron
querya query point
c_beginthe beginning of the destination range with the computed coordinates
policyone of the CGAL::Barycentric_coordinates::Computation_policy_3; the default is Computation_policy_3::FAST_WITH_EDGE_CASES
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
Precondition
num_vertices(polygon_mesh) >= 4.
polygon_mesh is strongly convex.
polygon_mesh is simplicial.