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 1, while all other coordinates are 0. If query is on the face, the three corresponding coordinates are triangle coordinates, while all other coordinates are set to 0. If query is not on the boundary, all the coordinates are set to 0.
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
np
an optional sequence of Named Parameters among the ones listed below
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
Optional Named Parameters
a property map associating points to the vertices of tmesh
Type: a class model of ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key type and Point as value type
Default: boost::get(CGAL::vertex_point, tmesh)
Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in TriangleMesh.
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.
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 tmesh, 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.
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\).
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
Returns
an array std::array<GeomTraits::FT, 4> with the computed coordinates
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\).
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
Returns
an output iterator to the element in the destination range, one past the last coordinate stored
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 tmesh, 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.