-
a property map associating points to the vertices of
pmesh -
Type: a class model of
ReadablePropertyMapwithboost::graph_traits<PolygonMesh>::vertex_descriptoras key type andPoint_3as value type -
Default:
boost::get(CGAL::vertex_point, pmesh)
|
CGAL 6.1 - Dynamic Skeletonization Via Variational Medial Axis Sampling
|
#include <CGAL/fast_winding_number.h>
Fast evaluation of the (normalized) winding number of a triangle mesh.
This class implements a hierarchical (multipole / cluster–far field) approximation of the generalized winding number for a triangle mesh.
| TriangleMesh | a model of FaceListGraph |
| FaceNormalMap | a readable property map contains the normal of each face: face_descriptor -> Vector_3 |
| FaceAreaMap | a readable property map contains the area of each face: face_descriptor -> FT |
| FaceCentroidMap | a readable property map contains the centroid of each face: face_descriptor -> Point_3 |
| Tree | a built AABB tree over the faces of TriangleMesh |
| ORDER | the order of the taylor expansion, must be 1, 2, or 3. Default: 3 |
| GeomTraits | a model of KernelDefault: boost::property_traits<
boost::property_map<TriangleMesh, CGAL::vertex_point_t>::type
>::value_type
ReadWritePropertyMap with boost::graph_traits<TriangleMesh_>::vertex_descriptor as key and GeomTraits_::Point_3 as value type.Default: boost::property_map<TriangleMesh_, CGAL::vertex_point_t>::const_type.
|
Public Member Functions | |
| template<class NamedParameters = parameters::Default_named_parameters> | |
| FT | evaluate_fast_winding_number (const Point_3 &p, const NamedParameters &np=parameters::default_values()) |
| computes the fast winding number of a given query point. | |
| FT | exact_winding_number (const Point_3 &p) const |
computes the exact winding number of a given query point p and the input mesh. | |
| bool | is_inside (const Point_3 &p) |
| tests whether a given query point is inside the input mesh. | |
Constructor | |||||||||||||||||||
The constructor of a vmas object.
| |||||||||||||||||||
| template<class NamedParameters = parameters::Default_named_parameters> | |||||||||||||||||||
| Fast_winding_number (const TriangleMesh &tmesh, const FaceNormalMap &fnm, const FaceAreaMap &fam, const FaceCentroidMap &fcm, Tree &tree, const NamedParameters &np=parameters::default_values()) | |||||||||||||||||||
| FT CGAL::Fast_winding_number< TriangleMesh, FaceNormalMap, FaceAreaMap, FaceCentroidMap, Tree, ORDER, GeomTraits, VertexPointMap >::evaluate_fast_winding_number | ( | const Point_3 & | p, |
| const NamedParameters & | np = parameters::default_values() |
||
| ) |
computes the fast winding number of a given query point.
This function computes the fast winding number of a given query point p recursively. For a given query point p, we first compute its distance r to the center of root node in the AABB tree. If r is larger than beta times the radius of the root node, we consider the root node as a single dipole and use the precomputed Taylor expansion coefficients to evaluate the winding number. Otherwise, we traverse down the tree recursively and repeat the same process for each child node.
| p | the query point | |
| np | an optional sequence of Named Parameters, listed below:
|
p | FT CGAL::Fast_winding_number< TriangleMesh, FaceNormalMap, FaceAreaMap, FaceCentroidMap, Tree, ORDER, GeomTraits, VertexPointMap >::exact_winding_number | ( | const Point_3 & | p | ) | const |
computes the exact winding number of a given query point p and the input mesh.
| p | the query point |
p | bool CGAL::Fast_winding_number< TriangleMesh, FaceNormalMap, FaceAreaMap, FaceCentroidMap, Tree, ORDER, GeomTraits, VertexPointMap >::is_inside | ( | const Point_3 & | p | ) |
tests whether a given query point is inside the input mesh.
| p | the query point |
true if p is inside the input mesh, false otherwise.