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... | |
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.
TriangleMesh | must be a model of the concept FaceListGraph . |
Point_3 | a model of Kernel::Point_3 |
OutIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
VertexToPointMap | a 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> . |
triangle_mesh | an instance of TriangleMesh , which must be a convex simplicial polyhedron |
query | a query point |
c_begin | the beginning of the destination range with the computed coordinates |
vertex_to_point_map | an instance of VertexToPointMap that maps a vertex from triangle_mesh to Point_3 ; the default initialization is provided |
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.
TriangleMesh | must be a model of the concept FaceListGraph . |
GeomTraits | a model of BarycentricTraits_3 |
OutIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
VertexToPointMap | a 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> . |
triangle_mesh | an instance of TriangleMesh , which must be a convex simplicial polyhedron |
query | a query point |
c_begin | the beginning of the destination range with the computed coordinates |
traits | a traits class with geometric objects, predicates, and constructions; the default initialization is provided |
vertex_to_point_map | an instance of VertexToPointMap that maps a vertex from triangle_mesh to Point_3 ; the default initialization is provided |
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.
Point_3 | A model of Kernel::Point_3 . |
TriangleMesh | must be a model of the concept FaceListGraph . |
OutIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
triangle_mesh | an instance of TriangleMesh , which must be a convex simplicial polyhedron |
query | a query point |
c_begin | the beginning of the destination range with the computed coordinates |
policy | one of the CGAL::Barycentric_coordinates::Computation_policy_3 ; the default is Computation_policy_3::FAST_WITH_EDGE_CASES |
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.
Point_3 | A model of Kernel::Point_3 . |
TriangleMesh | must be a model of the concept FaceListGraph . |
OutIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
triangle_mesh | an instance of TriangleMesh , which must be a simplicial polyhedron |
query | a query point |
c_begin | the beginning of the destination range with the computed coordinates |
policy | one of the CGAL::Barycentric_coordinates::Computation_policy_3 ; the default is Computation_policy_3::FAST_WITH_EDGE_CASES |
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\).
OutIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
GeomTraits | a model of BarycentricTraits_3 |
p0 | the first vertex of a tetrahedron |
p1 | the second vertex of a tetrahedron |
p2 | the third vertex of a tetrahedron |
p3 | the fourth vertex of a tetrahedron |
query | a query point |
c_begin | the beginning of the destination range with the computed coordinates |
traits | a traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type |
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\).
GeomTraits | a model of BarycentricTraits_3 |
p0 | the first vertex of a tetrahedron |
p1 | the second vertex of a tetrahedron |
p2 | the third vertex of a tetrahedron |
p3 | the fourth vertex of a tetrahedron |
query | a query point |
traits | a traits class with geometric objects, predicates, and constructions; this parameter can be omitted if the traits class can be deduced from the point type |
std::array<GeomTraits::FT, 4>
with the computed coordinatesOutIterator 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.
Point_3 | A model of Kernel::Point_3 . |
PolygonMesh | must be a model of the concept FaceListGraph . |
OutIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
polygon_mesh | an instance of PolygonMesh , which must be a convex simplicial polyhedron |
query | a query point |
c_begin | the beginning of the destination range with the computed coordinates |
policy | one of the CGAL::Barycentric_coordinates::Computation_policy_3 ; the default is Computation_policy_3::FAST_WITH_EDGE_CASES |