File failed to load: https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.3/config/TeX-MML-AM_CHTML/MathJax.js
DGtal 2.0.0
DGtal::WindingNumbersShape< TSpace > Struct Template Reference

Aim: model of a CEuclideanOrientedShape from an implicit function from an oriented point cloud. The implicit function is given by the generalized winding number of the oriented point cloud [10] . We use the libIGL implementation. More...

#include <DGtal/shapes/WindingNumbersShape.h>

Public Types

using Space = TSpace
 Types.
using RealPoint = typename Space::RealPoint
using RealVector = typename Space::RealVector
using Orientation = DGtal::Orientation

Public Member Functions

 WindingNumbersShape ()=delete
 WindingNumbersShape (ConstAlias< Eigen::MatrixXd > points, ConstAlias< Eigen::MatrixXd > normals, bool skipPointAreas=false)
 WindingNumbersShape (ConstAlias< Eigen::MatrixXd > points, ConstAlias< Eigen::MatrixXd > normals, const Eigen::VectorXd &areas)
void setPointAreas (ConstAlias< Eigen::VectorXd > areas)
Orientation orientation (const RealPoint aPoint, const double threshold=0.3) const
std::vector< OrientationorientationBatch (const Eigen::MatrixXd &queries, const double threshold=0.3) const
void rawWindingNumberBatch (const Eigen::MatrixXd &queries, Eigen::VectorXd &W) const

Data Fields

CountedConstPtrOrConstPtr< Eigen::MatrixXd > myPoints
 Const alias to the points.
CountedConstPtrOrConstPtr< Eigen::MatrixXd > myNormals
 Const alias to the normals.
Eigen::VectorXd myPointAreas
 Const alias to point area measure.
std::vector< std::vector< int > > myO_PI
 libIGL octree for fast queries data structure
Eigen::MatrixXi myO_CH
 libIGL octree for fast queries data structure
Eigen::MatrixXd myO_CN
 libIGL octree for fast queries data structure
Eigen::VectorXd myO_W
 libIGL octree for fast queries data structure

Detailed Description

template<typename TSpace>
struct DGtal::WindingNumbersShape< TSpace >

Aim: model of a CEuclideanOrientedShape from an implicit function from an oriented point cloud. The implicit function is given by the generalized winding number of the oriented point cloud [10] . We use the libIGL implementation.

Description of template class 'WindingNumbersShape'

See also
testWindingNumberShape, windingNumberShape
Template Parameters
TSPacethe digital space type (a model of CSpace)

Definition at line 77 of file WindingNumbersShape.h.

Member Typedef Documentation

◆ Orientation

template<typename TSpace>
using DGtal::WindingNumbersShape< TSpace >::Orientation = DGtal::Orientation

Definition at line 83 of file WindingNumbersShape.h.

◆ RealPoint

template<typename TSpace>
using DGtal::WindingNumbersShape< TSpace >::RealPoint = typename Space::RealPoint

Definition at line 81 of file WindingNumbersShape.h.

◆ RealVector

template<typename TSpace>
using DGtal::WindingNumbersShape< TSpace >::RealVector = typename Space::RealVector

Definition at line 82 of file WindingNumbersShape.h.

◆ Space

template<typename TSpace>
using DGtal::WindingNumbersShape< TSpace >::Space = TSpace

Types.

Definition at line 80 of file WindingNumbersShape.h.

Constructor & Destructor Documentation

◆ WindingNumbersShape() [1/3]

template<typename TSpace>
DGtal::WindingNumbersShape< TSpace >::WindingNumbersShape ( )
delete

◆ WindingNumbersShape() [2/3]

template<typename TSpace>
DGtal::WindingNumbersShape< TSpace >::WindingNumbersShape ( ConstAlias< Eigen::MatrixXd > points,
ConstAlias< Eigen::MatrixXd > normals,
bool skipPointAreas = false )
inline

Construct a WindingNumberShape Euclidean shape from an oriented point cloud. This constructor estimates the area of each point using CGAL.

If the number of points is greater than 20, CGAL is used to estimate the per-sample area (if skipPointAreas is set to false).

Parameters
pointsa "nx3" matrix with the sample coordinates.
normalsa "nx3" matrix for the normal vectors.
skipPointAreasif true, we do not estimate the point areas using CGAL (default: false)

Definition at line 97 of file WindingNumbersShape.h.

100 {
104 // Build octree, from libIGL tutorials
106 if (skipPointAreas)
107 trace.warning()<<"[WindingNumberShape] Skipping the CGAL point area estimation. By default, point areas are set to 1.0"<<std::endl;
108 else
109 if (points->rows()> 20)
110 {
113 // CGAL is only used to help get point areas
115 trace.info()<<"[WindingNumberShape] Min/max point area : "<<myPointAreas.minCoeff()<<" -- "<<myPointAreas.maxCoeff()<<std::endl;
116 }
117 else
118 {
119 trace.warning()<<"[WindingNumberShape] Too few points to use CGAL point_areas. Using the constant area setting."<<std::endl;
120 }
121 }
Aim: model of a CEuclideanOrientedShape from an implicit function from an oriented point cloud....
Eigen::VectorXd myPointAreas
Const alias to point area measure.
Eigen::MatrixXd myO_CN
libIGL octree for fast queries data structure
Eigen::VectorXd myO_W
libIGL octree for fast queries data structure
CountedConstPtrOrConstPtr< Eigen::MatrixXd > myPoints
Const alias to the points.
Eigen::MatrixXi myO_CH
libIGL octree for fast queries data structure
std::vector< std::vector< int > > myO_PI
libIGL octree for fast queries data structure
CountedConstPtrOrConstPtr< Eigen::MatrixXd > myNormals
Const alias to the normals.

References myNormals, myO_CH, myO_CN, myO_PI, myO_W, myPointAreas, myPoints, and DGtal::trace.

◆ WindingNumbersShape() [3/3]

template<typename TSpace>
DGtal::WindingNumbersShape< TSpace >::WindingNumbersShape ( ConstAlias< Eigen::MatrixXd > points,
ConstAlias< Eigen::MatrixXd > normals,
const Eigen::VectorXd & areas )
inline

Construct a WindingNumberShape Euclidean shape from an oriented point cloud. For this constructor, the area of each point is given by the user.

Parameters
pointsa "nx3" matrix with the sample coordinates.
normalsa "nx3" matrix for the normal vectors.
areasa "n" vector with the area of each point.

Definition at line 129 of file WindingNumbersShape.h.

References myNormals, myO_CH, myO_CN, myO_PI, myO_W, myPointAreas, and myPoints.

Member Function Documentation

◆ orientation()

template<typename TSpace>
Orientation DGtal::WindingNumbersShape< TSpace >::orientation ( const RealPoint aPoint,
const double threshold = 0.3 ) const
inline

Orientation of a point using the winding number value from an oriented pointcloud.

Note
For multiple queries, orientationBatch() should be used.
Parameters
aPoint[in] a point in space
threshold[in] the iso-value of the surface of the winding number implicit map (default = 0.3).
Returns
a DGtal::Orientation value

Definition at line 155 of file WindingNumbersShape.h.

157 {
159 queries << aPoint(0) , aPoint(1) , aPoint(2);
161 return singlePoint[0];
162 }
std::vector< Orientation > orientationBatch(const Eigen::MatrixXd &queries, const double threshold=0.3) const
const Point aPoint(3, 4)

References aPoint, and orientationBatch().

◆ orientationBatch()

template<typename TSpace>
std::vector< Orientation > DGtal::WindingNumbersShape< TSpace >::orientationBatch ( const Eigen::MatrixXd & queries,
const double threshold = 0.3 ) const
inline

Orientation of a set of points (queries) using the winding number value from an oriented pointcloud.

Parameters
queries[in] a "nx3" matrix with the query points in space.
threshold[in] the iso-value of the surface of the winding number implicit map (default = 0.3).
Returns
a DGtal::Orientation value vector for each query point.

Definition at line 170 of file WindingNumbersShape.h.

172 {
175
177
178 //Reformating the output
179 for(auto i=0u; i < queries.rows(); ++i)
180 {
181 if (std::abs(W(i)) < threshold )
183 else
184 if (std::abs(W(i)) > threshold)
186 else
188 }
189 return results;
190 }
void rawWindingNumberBatch(const Eigen::MatrixXd &queries, Eigen::VectorXd &W) const

References DGtal::INSIDE, DGtal::ON, DGtal::OUTSIDE, and rawWindingNumberBatch().

Referenced by orientation().

◆ rawWindingNumberBatch()

template<typename TSpace>
void DGtal::WindingNumbersShape< TSpace >::rawWindingNumberBatch ( const Eigen::MatrixXd & queries,
Eigen::VectorXd & W ) const
inline

Returns the raw value of the Winding Number funciton at a set of points (queries).

Parameters
queries[in] a "nx3" matrix with the query points in space.
W[out] a vector with all windung number values.

Definition at line 197 of file WindingNumbersShape.h.

References myNormals, myO_CH, myO_PI, myPointAreas, and myPoints.

Referenced by orientationBatch().

◆ setPointAreas()

template<typename TSpace>
void DGtal::WindingNumbersShape< TSpace >::setPointAreas ( ConstAlias< Eigen::VectorXd > areas)
inline

Set the area map for each point.

Parameters
areasa Eigen vector of estimated area for each input point.

Definition at line 142 of file WindingNumbersShape.h.

143 {
145 }

References myPointAreas.

Field Documentation

◆ myNormals

template<typename TSpace>
CountedConstPtrOrConstPtr<Eigen::MatrixXd> DGtal::WindingNumbersShape< TSpace >::myNormals

Const alias to the normals.

Definition at line 212 of file WindingNumbersShape.h.

Referenced by rawWindingNumberBatch(), WindingNumbersShape(), and WindingNumbersShape().

◆ myO_CH

template<typename TSpace>
Eigen::MatrixXi DGtal::WindingNumbersShape< TSpace >::myO_CH

libIGL octree for fast queries data structure

Definition at line 219 of file WindingNumbersShape.h.

Referenced by rawWindingNumberBatch(), WindingNumbersShape(), and WindingNumbersShape().

◆ myO_CN

template<typename TSpace>
Eigen::MatrixXd DGtal::WindingNumbersShape< TSpace >::myO_CN

libIGL octree for fast queries data structure

Definition at line 221 of file WindingNumbersShape.h.

Referenced by WindingNumbersShape(), and WindingNumbersShape().

◆ myO_PI

template<typename TSpace>
std::vector<std::vector<int > > DGtal::WindingNumbersShape< TSpace >::myO_PI

libIGL octree for fast queries data structure

Definition at line 217 of file WindingNumbersShape.h.

Referenced by rawWindingNumberBatch(), WindingNumbersShape(), and WindingNumbersShape().

◆ myO_W

template<typename TSpace>
Eigen::VectorXd DGtal::WindingNumbersShape< TSpace >::myO_W

libIGL octree for fast queries data structure

Definition at line 223 of file WindingNumbersShape.h.

Referenced by WindingNumbersShape(), and WindingNumbersShape().

◆ myPointAreas

template<typename TSpace>
Eigen::VectorXd DGtal::WindingNumbersShape< TSpace >::myPointAreas

Const alias to point area measure.

Definition at line 214 of file WindingNumbersShape.h.

Referenced by rawWindingNumberBatch(), setPointAreas(), WindingNumbersShape(), and WindingNumbersShape().

◆ myPoints

template<typename TSpace>
CountedConstPtrOrConstPtr<Eigen::MatrixXd> DGtal::WindingNumbersShape< TSpace >::myPoints

Const alias to the points.

Definition at line 210 of file WindingNumbersShape.h.

Referenced by rawWindingNumberBatch(), WindingNumbersShape(), and WindingNumbersShape().


The documentation for this struct was generated from the following file: