CGAL 5.1 - 3D Generalized Barycentric Coordinates
CGAL::Barycentric_coordinates::Discrete_harmonic_coordinates_3< TriangleMesh, GeomTraits, VertexToPointMap > Class Template Reference

#include <CGAL/Barycentric_coordinates_3/Discrete_harmonic_coordinates_3.h>

Definition

template<typename TriangleMesh, typename GeomTraits, typename VertexToPointMap = typename property_map_selector<TriangleMesh, CGAL::vertex_point_t>::const_type>
class CGAL::Barycentric_coordinates::Discrete_harmonic_coordinates_3< TriangleMesh, GeomTraits, VertexToPointMap >

3D discrete harmonic coordinates.

This class implements 3D discrete harmonic coordinates [3], which can be computed at any point inside a convex polyhedron with triangular faces.

Discrete harmonic coordinates are well-defined in the closure of a convex polyhedron with triangular faces but they are not necessarily positive. The coordinates are computed analytically.

Template Parameters
TriangleMeshmust be a model of the concept FaceListGraph.
GeomTraitsa model of BarycentricTraits_3
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>.

Types

typedef GeomTraits::FT FT
 Number type.
 
typedef GeomTraits::Point_3 Point_3
 Point type.
 
typedef GeomTraits::Vector_3 Vector_3
 Vector type.
 

Initialization

 Discrete_harmonic_coordinates_3 (const TriangleMesh &triangle_mesh, const Computation_policy_3 policy, const VertexToPointMap vertex_to_point_map, const GeomTraits traits=GeomTraits())
 initializes all internal data structures. More...
 

Access

template<typename OutIterator >
OutIterator operator() (const Point_3 &query, OutIterator c_begin)
 computes 3D discrete harmonic coordinates. More...
 

Constructor & Destructor Documentation

◆ Discrete_harmonic_coordinates_3()

template<typename TriangleMesh , typename GeomTraits , typename VertexToPointMap = typename property_map_selector<TriangleMesh, CGAL::vertex_point_t>::const_type>
CGAL::Barycentric_coordinates::Discrete_harmonic_coordinates_3< TriangleMesh, GeomTraits, VertexToPointMap >::Discrete_harmonic_coordinates_3 ( const TriangleMesh &  triangle_mesh,
const Computation_policy_3  policy,
const VertexToPointMap  vertex_to_point_map,
const GeomTraits  traits = GeomTraits() 
)

initializes all internal data structures.

This class implements the behavior of discrete harmonic coordinates for 3D query points.

Parameters
triangle_meshan instance of TriangleMesh, which must be a convex simplicial polyhedron
policyone of the CGAL::Barycentric_coordinates::Computation_policy_3; the default is Computation_policy_3::FAST_WITH_EDGE_CASES
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
Precondition
num_vertices(triangle_mesh) >= 4.
triangle_mesh is strongly convex.
triangle_mesh is simplicial.

Member Function Documentation

◆ operator()()

template<typename TriangleMesh , typename GeomTraits , typename VertexToPointMap = typename property_map_selector<TriangleMesh, CGAL::vertex_point_t>::const_type>
template<typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::Discrete_harmonic_coordinates_3< TriangleMesh, GeomTraits, VertexToPointMap >::operator() ( const Point_3 query,
OutIterator  c_begin 
)

computes 3D discrete harmonic coordinates.

This function fills c_begin with 3D discrete harmonic coordinates computed at the query point with respect to the vertices of the input polyhedron.

The number of returned coordinates equals to the number of vertices.

After the coordinates \(b_i\) with \(i = 1\dots n\) are computed, where \(n\) is the number of vertices, the query point \(q\) can be obtained as \(q = \sum_{i = 1}^{n}b_ip_i\), where \(p_i\) are the polyhedron vertices.

Template Parameters
OutIteratora model of OutputIterator that accepts values of type FT
Parameters
querya query point
c_beginthe beginning of the destination range with the computed coordinates
Returns
an output iterator to the element in the destination range, one past the last coordinate stored