DGtal  1.5.beta
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. More...
 
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. More...
 
CountedConstPtrOrConstPtr< Eigen::MatrixXd > myNormals
 Const alias to the normals. More...
 
Eigen::VectorXd myPointAreas
 Const alias to point area measure. More...
 
std::vector< std::vector< int > > myO_PI
 libIGL octree for fast queries data structure More...
 
Eigen::MatrixXi myO_CH
 libIGL octree for fast queries data structure More...
 
Eigen::MatrixXd myO_CN
 libIGL octree for fast queries data structure More...
 
Eigen::VectorXd myO_W
 libIGL octree for fast queries data structure More...
 

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  {
101  myPoints = points;
102  myNormals = normals;
103  myPointAreas = Eigen::VectorXd::Ones(myPoints->rows());
104  // Build octree, from libIGL tutorials
105  igl::octree(*myPoints,myO_PI,myO_CH,myO_CN,myO_W);
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  {
111  Eigen::MatrixXi I;
112  igl::knn(*myPoints,(int)points->rows(),myO_PI,myO_CH,myO_CN,myO_W,I);
113  // CGAL is only used to help get point areas
114  igl::copyleft::cgal::point_areas(*myPoints,I,*myNormals,myPointAreas);
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  }
std::ostream & info()
std::ostream & warning()
Trace trace
Definition: Common.h:153
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 DGtal::Trace::info(), DGtal::WindingNumbersShape< TSpace >::myNormals, DGtal::WindingNumbersShape< TSpace >::myO_CH, DGtal::WindingNumbersShape< TSpace >::myO_CN, DGtal::WindingNumbersShape< TSpace >::myO_PI, DGtal::WindingNumbersShape< TSpace >::myO_W, DGtal::WindingNumbersShape< TSpace >::myPointAreas, DGtal::WindingNumbersShape< TSpace >::myPoints, DGtal::trace, and DGtal::Trace::warning().

◆ 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.

132  {
133  myPoints = points;
134  myNormals = normals;
135  myPointAreas = areas;
136  igl::octree(*myPoints,myO_PI,myO_CH,myO_CN,myO_W);
137  }

References DGtal::WindingNumbersShape< TSpace >::myNormals, DGtal::WindingNumbersShape< TSpace >::myO_CH, DGtal::WindingNumbersShape< TSpace >::myO_CN, DGtal::WindingNumbersShape< TSpace >::myO_PI, DGtal::WindingNumbersShape< TSpace >::myO_W, DGtal::WindingNumbersShape< TSpace >::myPointAreas, and DGtal::WindingNumbersShape< TSpace >::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  {
158  Eigen::MatrixXd queries(1,3);
159  queries << aPoint(0) , aPoint(1) , aPoint(2);
160  auto singlePoint = orientationBatch(queries, threshold);
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 DGtal::WindingNumbersShape< TSpace >::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  {
173  Eigen::VectorXd W;
174  std::vector<Orientation> results( queries.rows(), DGtal::OUTSIDE );
175 
176  rawWindingNumberBatch(queries, W);
177 
178  //Reformating the output
179  for(auto i=0u; i < queries.rows(); ++i)
180  {
181  if (std::abs(W(i)) < threshold )
182  results[i] = DGtal::OUTSIDE;
183  else
184  if (std::abs(W(i)) > threshold)
185  results[i] = DGtal::INSIDE;
186  else
187  results[i] = DGtal::ON;
188  }
189  return results;
190  }
@ INSIDE
Definition: Common.h:141
@ OUTSIDE
Definition: Common.h:141
@ ON
Definition: Common.h:141
void rawWindingNumberBatch(const Eigen::MatrixXd &queries, Eigen::VectorXd &W) const

References DGtal::INSIDE, DGtal::ON, DGtal::OUTSIDE, and DGtal::WindingNumbersShape< TSpace >::rawWindingNumberBatch().

Referenced by DGtal::WindingNumbersShape< TSpace >::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.

199  {
200  Eigen::MatrixXd O_CM;
201  Eigen::VectorXd O_R;
202  Eigen::MatrixXd O_EC;
203 
204  //Computing the WN values
205  igl::fast_winding_number(*myPoints,*myNormals,myPointAreas,myO_PI,myO_CH,2,O_CM,O_R,O_EC);
206  igl::fast_winding_number(*myPoints,*myNormals,myPointAreas,myO_PI,myO_CH,O_CM,O_R,O_EC,queries,2,W);
207  }

References DGtal::WindingNumbersShape< TSpace >::myNormals, DGtal::WindingNumbersShape< TSpace >::myO_CH, DGtal::WindingNumbersShape< TSpace >::myO_PI, DGtal::WindingNumbersShape< TSpace >::myPointAreas, and DGtal::WindingNumbersShape< TSpace >::myPoints.

Referenced by DGtal::WindingNumbersShape< TSpace >::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  {
144  myPointAreas = areas;
145  }

References DGtal::WindingNumbersShape< TSpace >::myPointAreas.

Field Documentation

◆ myNormals

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

◆ 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 DGtal::WindingNumbersShape< TSpace >::rawWindingNumberBatch(), and DGtal::WindingNumbersShape< TSpace >::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 DGtal::WindingNumbersShape< TSpace >::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 DGtal::WindingNumbersShape< TSpace >::rawWindingNumberBatch(), and DGtal::WindingNumbersShape< TSpace >::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 DGtal::WindingNumbersShape< TSpace >::WindingNumbersShape().

◆ myPointAreas

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

◆ myPoints

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

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