CGAL 5.1.2 - Classification
|
#include <CGAL/Classification/Local_eigen_analysis.h>
Class that precomputes and stores the eigenvectors and eigenvalues of the covariance matrices of all points of a point set using a local neighborhood.
This class can be used to compute eigen features (see Predefined Features) and to estimate local normal vectors and tangent planes.
Public Types | |
typedef std::array< float, 3 > | Eigenvalues |
Eigenvalues (sorted in ascending order) | |
Named Constructors | |
template<typename PointRange , typename PointMap , typename NeighborQuery , typename ConcurrencyTag , typename DiagonalizeTraits > | |
static Local_eigen_analysis | create_from_point_set (const PointRange &input, PointMap point_map, const NeighborQuery &neighbor_query, const ConcurrencyTag &=ConcurrencyTag(), const DiagonalizeTraits &=DiagonalizeTraits()) |
Computes the local eigen analysis of an input point set based on a local neighborhood. More... | |
template<typename FaceListGraph , typename NeighborQuery , typename ConcurrencyTag , typename DiagonalizeTraits > | |
static Local_eigen_analysis | create_from_face_graph (const FaceListGraph &input, const NeighborQuery &neighbor_query, const ConcurrencyTag &=ConcurrencyTag(), const DiagonalizeTraits &=DiagonalizeTraits()) |
Computes the local eigen analysis of an input face graph based on a local neighborhood. More... | |
template<typename ClusterRange , typename ConcurrencyTag , typename DiagonalizeTraits > | |
static Local_eigen_analysis | create_from_point_clusters (const ClusterRange &input, const ConcurrencyTag &=ConcurrencyTag(), const DiagonalizeTraits &=DiagonalizeTraits()) |
Computes the local eigen analysis of an input set of point clusters based on a local neighborhood. More... | |
Analysis | |
template<typename GeomTraits > | |
GeomTraits::Vector_3 | normal_vector (std::size_t index) const |
Returns the estimated unoriented normal vector of the point at position index . More... | |
template<typename GeomTraits > | |
GeomTraits::Plane_3 | plane (std::size_t index) const |
Returns the estimated local tangent plane of the point at position index . More... | |
Eigenvalues | eigenvalue (std::size_t index) const |
Returns the normalized eigenvalues of the point at position index . | |
|
static |
Computes the local eigen analysis of an input face graph based on a local neighborhood.
FaceListGraph | model of FaceListGraph . |
NeighborQuery | model of NeighborQuery |
ConcurrencyTag | enables sequential versus parallel algorithm. Possible values are Parallel_tag (default value if CGAL is linked with TBB) or Sequential_tag (default value otherwise). |
DiagonalizeTraits | model of DiagonalizeTraits used for matrix diagonalization. It can be omitted: if Eigen 3 (or greater) is available and CGAL_EIGEN3_ENABLED is defined then an overload using Eigen_diagonalize_traits is provided. Otherwise, the internal implementation Diagonalize_traits is used. |
input | face graph. |
neighbor_query | object used to access neighborhoods of points. |
|
static |
Computes the local eigen analysis of an input set of point clusters based on a local neighborhood.
ClusterRange | model of ConstRange . Its iterator type is RandomAccessIterator and its value type is the key type of PointMap . |
ConcurrencyTag | enables sequential versus parallel algorithm. Possible values are Parallel_tag (default value if CGAL is linked with TBB) or Sequential_tag (default value otherwise). |
DiagonalizeTraits | model of DiagonalizeTraits used for matrix diagonalization. It can be omitted: if Eigen 3 (or greater) is available and CGAL_EIGEN3_ENABLED is defined then an overload using Eigen_diagonalize_traits is provided. Otherwise, the internal implementation Diagonalize_traits is used. |
input | cluster range. |
|
static |
Computes the local eigen analysis of an input point set based on a local neighborhood.
PointRange | model of ConstRange . Its iterator type is RandomAccessIterator and its value type is the key type of PointMap . |
PointMap | model of ReadablePropertyMap whose key type is the value type of the iterator of PointRange and value type is CGAL::Point_3 . |
NeighborQuery | model of NeighborQuery |
ConcurrencyTag | enables sequential versus parallel algorithm. Possible values are Parallel_tag (default value if CGAL is linked with TBB) or Sequential_tag (default value otherwise). |
DiagonalizeTraits | model of DiagonalizeTraits used for matrix diagonalization. It can be omitted if Eigen 3 (or greater) is available and CGAL_EIGEN3_ENABLED is defined. In that case, an overload using Eigen_diagonalize_traits is provided. |
input | point range. |
point_map | property map to access the input points. |
neighbor_query | object used to access neighborhoods of points. |
GeomTraits::Vector_3 CGAL::Classification::Local_eigen_analysis::normal_vector | ( | std::size_t | index | ) | const |
Returns the estimated unoriented normal vector of the point at position index
.
GeomTraits | model of CGAL Kernel. |
GeomTraits::Plane_3 CGAL::Classification::Local_eigen_analysis::plane | ( | std::size_t | index | ) | const |
Returns the estimated local tangent plane of the point at position index
.
GeomTraits | model of CGAL Kernel. |