|
CGAL 6.2 - 3D Generalized Barycentric Coordinates
|
Free functions to compute barycentric coordinates.
Functions | |
| template<typename TriangleMesh , typename CoordinateRange , typename VertexPointMap , typename GeomTraits = typename CGAL::Kernel_traits<typename boost::property_traits<VertexPointMap>::value_type>::type> | |
| boost::property_traits< VertexPointMap >::value_type | CGAL::Barycentric_coordinates::apply_barycentric_coordinates (const TriangleMesh &tmesh, const CoordinateRange &coordinates, VertexPointMap vpm) |
| computes a point location from barycentric coordinates with respect to a triangle mesh. | |
| template<typename PointRange , typename CoordinateRange , typename GeomTraits = typename boost::range_value<PointRange>::type> | |
| boost::range_value< PointRange >::type | CGAL::Barycentric_coordinates::apply_barycentric_coordinates (const PointRange &points, const CoordinateRange &coordinates) |
| computes a point location from barycentric coordinates with respect to a set of points. | |
| template<typename TriangleMesh , typename OutputIterator , typename GeomTraits , typename VertexPointMap > | |
| std::pair< OutputIterator, bool > | CGAL::Barycentric_coordinates::boundary_coordinates_3 (const TriangleMesh &tmesh, const typename GeomTraits::Point_3 &query, OutputIterator oi, const GeomTraits &traits, const VertexPointMap vertex_point_map) |
| computes boundary barycentric coordinates with respect to a closed convex triangle mesh. | |
| template<typename TriangleMesh , typename Point_3 , typename OutputIterator , typename VertexPointMap = typename boost::property_map<TriangleMesh, CGAL::vertex_point_t>::const_type> | |
| std::pair< OutputIterator, bool > | CGAL::Barycentric_coordinates::boundary_coordinates_3 (const TriangleMesh &tmesh, const Point_3 &query, OutputIterator oi, const VertexPointMap vertex_point_map) |
| computes boundary barycentric coordinates with respect to a closed convex triangle mesh. | |
| template<typename Point_3 , typename TriangleMesh , typename OutputIterator > | |
| OutputIterator | CGAL::Barycentric_coordinates::discrete_harmonic_coordinates_3 (const TriangleMesh &tmesh, const Point_3 &query, OutputIterator oi, const Computation_policy_3 policy=Computation_policy_3::FAST) |
| computes 3D discrete harmonic coordinates with respect to a closed convex triangle mesh. | |
| template<typename Point_3 , typename TriangleMesh , typename OutputIterator > | |
| OutputIterator | CGAL::Barycentric_coordinates::mean_value_coordinates_3 (const TriangleMesh &tmesh, const Point_3 &query, OutputIterator oi, const Computation_policy_3 policy=Computation_policy_3::FAST_WITH_EDGE_CASES) |
| computes 3D mean value barycentric coordinates with respect to a closed triangle mesh. | |
| template<typename OutputIterator , typename GeomTraits > | |
| OutputIterator | 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, OutputIterator oi, const GeomTraits &traits) |
| computes barycentric coordinates with respect to a tetrahedron. | |
| template<typename GeomTraits > | |
| std::array< typename GeomTraits::FT, 4 > | 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, const GeomTraits &traits) |
| computes barycentric coordinates with respect to a tetrahedron. | |
| template<typename Point_3 , typename TriangleMesh , typename OutputIterator > | |
| OutputIterator | CGAL::Barycentric_coordinates::wachspress_coordinates_3 (const TriangleMesh &tmesh, const Point_3 &query, OutputIterator oi, const Computation_policy_3 policy=Computation_policy_3::FAST) |
| computes 3D Wachspress coordinates with respect to a closed convex triangle mesh. | |
| boost::range_value< PointRange >::type CGAL::Barycentric_coordinates::apply_barycentric_coordinates | ( | const PointRange & | points, |
| const CoordinateRange & | coordinates | ||
| ) |
#include <CGAL/Barycentric_coordinates_3.h>
computes a point location from barycentric coordinates with respect to a set of points.
This function a computes point location from barycentric coordinates with respect to a set of points.
| PointRange | is a model of ConstRange and RandomAccessRange with value type GeomTraits::Point_3 |
| CoordinateRange | a range whose iterator is a model of ForwardIterator with value type GeomTraits::FT |
| GeomTraits | a model of BarycentricTraits_3 |
| points | a range of input points |
| coordinates | barycentric coordinates of the point |
GeomTraits::Point_3pts.size() == coords.size() | boost::property_traits< VertexPointMap >::value_type CGAL::Barycentric_coordinates::apply_barycentric_coordinates | ( | const TriangleMesh & | tmesh, |
| const CoordinateRange & | coordinates, | ||
| VertexPointMap | vpm | ||
| ) |
#include <CGAL/Barycentric_coordinates_3.h>
computes a point location from barycentric coordinates with respect to a triangle mesh.
This function computes a point location from barycentric coordinates with respect to the vertices of tmesh
| TriangleMesh | must be a model of the concept FaceListGraph |
| CoordinateRange | a range whose iterator is a model of ForwardIterator with value type GeomTraits::FT |
| VertexPointMap | a property map with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and GeomTraits::Point_3 as value type. |
| GeomTraits | a model of BarycentricTraits_3 |
| tmesh | an instance of TriangleMesh |
| coordinates | barycentric coordinates of the point |
| vpm | an instance of VertexPointMap that maps a vertex from tmesh to GeomTraits::Point_3 |
boost::property_traits<VertexPointMap>::value_typevertices(tmesh).size() == coords.size() | std::pair< OutputIterator, bool > CGAL::Barycentric_coordinates::boundary_coordinates_3 | ( | const TriangleMesh & | tmesh, |
| const Point_3 & | query, | ||
| OutputIterator | oi, | ||
| const VertexPointMap | vertex_point_map | ||
| ) |
#include <CGAL/Barycentric_coordinates_3/boundary_coordinates_3.h>
computes boundary barycentric coordinates with respect to a closed convex triangle mesh.
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 oi.
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 |
| OutputIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
| VertexPointMap | 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>. |
| tmesh | an instance of TriangleMesh, which must be a convex simplicial polyhedron |
| query | a query point |
| oi | the beginning of the destination range with the computed coordinates |
| vertex_point_map | an instance of VertexPointMap that maps a vertex from tmesh to Point_3; the default initialization is provided |
tmesh) >= 4. tmesh) tmesh). tmesh). | std::pair< OutputIterator, bool > CGAL::Barycentric_coordinates::boundary_coordinates_3 | ( | const TriangleMesh & | tmesh, |
| const typename GeomTraits::Point_3 & | query, | ||
| OutputIterator | oi, | ||
| const GeomTraits & | traits, | ||
| const VertexPointMap | vertex_point_map | ||
| ) |
#include <CGAL/Barycentric_coordinates_3/boundary_coordinates_3.h>
computes boundary barycentric coordinates with respect to a closed convex triangle mesh.
This function computes boundary barycentric coordinates at a given query point with respect to the vertices of a convex simplicial polyhedron, that is one coordinate per vertex. The coordinates are stored in a destination range beginning at oi.
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 |
| OutputIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
| VertexPointMap | a property map with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and GeomTraits::Point_3 as value type |
| tmesh | an instance of TriangleMesh |
| query | a query point |
| oi | 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_point_map | an instance of VertexPointMap that maps a vertex from tmesh to Point_3; the default initialization is provided |
tmesh) >= 4. tmesh). tmesh). tmesh). | OutputIterator CGAL::Barycentric_coordinates::discrete_harmonic_coordinates_3 | ( | const TriangleMesh & | tmesh, |
| const Point_3 & | query, | ||
| OutputIterator | oi, | ||
| const Computation_policy_3 | policy = Computation_policy_3::FAST |
||
| ) |
#include <CGAL/Barycentric_coordinates_3/Discrete_harmonic_coordinates_3.h>
computes 3D discrete harmonic coordinates with respect to a closed convex triangle mesh.
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 oi.
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 |
| OutputIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
| tmesh | an instance of TriangleMesh |
| query | a query point |
| oi | 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 |
tmesh) >= 4. tmesh). tmesh). tmesh). | OutputIterator CGAL::Barycentric_coordinates::mean_value_coordinates_3 | ( | const TriangleMesh & | tmesh, |
| const Point_3 & | query, | ||
| OutputIterator | oi, | ||
| 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 barycentric coordinates with respect to a closed triangle mesh.
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 oi.
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 |
| OutputIterator | a model of OutputIterator that accepts values of type Kernel_traits<Point_3>::Kernel::FT |
| tmesh | an instance of TriangleMesh |
| query | a query point |
| oi | 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 |
tmesh) >= 4. tmesh) tmesh). | std::array< typename GeomTraits::FT, 4 > 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, | ||
| const GeomTraits & | traits | ||
| ) |
#include <CGAL/Barycentric_coordinates_3/tetrahedron_coordinates.h>
computes barycentric coordinates with respect to a tetrahedron.
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 an array.
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 coordinatestraits.compute_volume_3_object()(p0, p1, p2, p3) != 0 | OutputIterator 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, | ||
| OutputIterator | oi, | ||
| const GeomTraits & | traits | ||
| ) |
#include <CGAL/Barycentric_coordinates_3/tetrahedron_coordinates.h>
computes barycentric coordinates with respect to a tetrahedron.
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 oi.
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\).
| OutputIterator | 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 |
| oi | 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 |
traits.compute_volume_3_object()(p0, p1, p2, p3) != 0 | OutputIterator CGAL::Barycentric_coordinates::wachspress_coordinates_3 | ( | const TriangleMesh & | tmesh, |
| const Point_3 & | query, | ||
| OutputIterator | oi, | ||
| const Computation_policy_3 | policy = Computation_policy_3::FAST |
||
| ) |
#include <CGAL/Barycentric_coordinates_3/Wachspress_coordinates_3.h>
computes 3D Wachspress coordinates with respect to a closed convex triangle mesh.
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 oi.
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 |
| TriangleMesh | must be a model of the concept FaceListGraph |
| OutputIterator | a model of OutputIterator that accepts values of type GeomTraits::FT |
| tmesh | an instance of TriangleMesh |
| query | a query point |
| oi | 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 |
tmesh) >= 4. tmesh). tmesh). tmesh).