DGtal  1.4.beta
WindingNumbersShape.h
1 
17 #pragma once
18 
31 #if !defined WindingNumbersShape_h
33 #define WindingNumbersShape_h
34 
35 #ifndef WITH_LIBIGL
36 #error You need to have activated LIBIGL (WITH_LIBIGL flag) to include this file.
37 #endif
38 
39 #include <vector>
40 #include <algorithm>
41 #include <DGtal/base/Common.h>
42 #include <DGtal/base/CountedConstPtrOrConstPtr.h>
43 #include <DGtal/base/ConstAlias.h>
44 
45 #include <DGtal/shapes/CEuclideanOrientedShape.h>
46 #include <igl/fast_winding_number.h>
47 #include <igl/octree.h>
48 #include <igl/knn.h>
49 #include <igl/copyleft/cgal/point_areas.h>
50 namespace DGtal
51 {
53 
65  template<typename TSpace>
67  {
69  using Space = TSpace;
70  using RealPoint = typename Space::RealPoint;
71  using RealVector = typename Space::RealVector;
73 
74  //Removing Default constructor
75  WindingNumbersShape() = delete;
76 
88  bool skipPointAreas = false)
89  {
90  myPoints = points;
91  myNormals = normals;
92  myPointAreas = Eigen::VectorXd::Ones(myPoints->rows());
93  // Build octree, from libIGL tutorials
94  igl::octree(*myPoints,myO_PI,myO_CH,myO_CN,myO_W);
95  if (skipPointAreas)
96  trace.warning()<<"[WindingNumberShape] Skipping the CGAL point area estimation. By default, point areas are set to 1.0"<<std::endl;
97  else
98  if (points->rows()> 20)
99  {
100  Eigen::MatrixXi I;
101  igl::knn(*myPoints,(int)points->rows(),myO_PI,myO_CH,myO_CN,myO_W,I);
102  // CGAL is only used to help get point areas
103  igl::copyleft::cgal::point_areas(*myPoints,I,*myNormals,myPointAreas);
104  trace.info()<<"[WindingNumberShape] Min/max point area : "<<myPointAreas.minCoeff()<<" -- "<<myPointAreas.maxCoeff()<<std::endl;
105  }
106  else
107  {
108  trace.warning()<<"[WindingNumberShape] Too few points to use CGAL point_areas. Using the constant area setting."<<std::endl;
109  }
110  }
111 
120  const Eigen::VectorXd &areas)
121  {
122  myPoints = points;
123  myNormals = normals;
124  myPointAreas = areas;
125  igl::octree(*myPoints,myO_PI,myO_CH,myO_CN,myO_W);
126  }
127 
128 
132  {
133  myPointAreas = areas;
134  }
135 
145  const double threshold = 0.3) const
146  {
147  Eigen::MatrixXd queries(1,3);
148  queries << aPoint(0) , aPoint(1) , aPoint(2);
149  auto singlePoint = orientationBatch(queries, threshold);
150  return singlePoint[0];
151  }
152 
159  std::vector<Orientation> orientationBatch(const Eigen::MatrixXd & queries,
160  const double threshold = 0.3) const
161  {
162  Eigen::VectorXd W;
163  std::vector<Orientation> results( queries.rows(), DGtal::OUTSIDE );
164  Eigen::MatrixXd O_CM;
165  Eigen::VectorXd O_R;
166  Eigen::MatrixXd O_EC;
167 
168  //Checking if the areas
169  igl::fast_winding_number(*myPoints,*myNormals,myPointAreas,myO_PI,myO_CH,2,O_CM,O_R,O_EC);
170  igl::fast_winding_number(*myPoints,*myNormals,myPointAreas,myO_PI,myO_CH,O_CM,O_R,O_EC,queries,2,W);
171 
172  //Reformating the output
173  for(auto i=0u; i < queries.rows(); ++i)
174  {
175  if (std::abs(W(i)) < threshold )
176  results[i] = DGtal::OUTSIDE;
177  else
178  if (std::abs(W(i)) > threshold)
179  results[i] = DGtal::INSIDE;
180  else
181  results[i] = DGtal::ON;
182  }
183  return results;
184  }
185 
191  Eigen::VectorXd myPointAreas;
192 
194  std::vector<std::vector<int > > myO_PI;
196  Eigen::MatrixXi myO_CH;
198  Eigen::MatrixXd myO_CN;
200  Eigen::VectorXd myO_W;
201 
202 
203  };
204 }
205 
206 #endif // !defined WindingNumbersShape_h
Aim: This class encapsulates its parameter class so that to indicate to the user that the object/poin...
Definition: ConstAlias.h:187
std::ostream & info()
std::ostream & warning()
Space::RealVector RealVector
DGtal is the top-level namespace which contains all DGtal functions and types.
Trace trace
Definition: Common.h:153
Orientation
Definition: Common.h:141
@ INSIDE
Definition: Common.h:141
@ OUTSIDE
Definition: Common.h:141
@ ON
Definition: Common.h:141
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
Orientation orientation(const RealPoint aPoint, const double threshold=0.3) const
typename Space::RealPoint RealPoint
std::vector< Orientation > orientationBatch(const Eigen::MatrixXd &queries, const double threshold=0.3) const
void setPointAreas(ConstAlias< Eigen::VectorXd > areas)
CountedConstPtrOrConstPtr< Eigen::MatrixXd > myPoints
Const alias to the points.
Eigen::MatrixXi myO_CH
libIGL octree for fast queries data structure
typename Space::RealVector RealVector
WindingNumbersShape(ConstAlias< Eigen::MatrixXd > points, ConstAlias< Eigen::MatrixXd > normals, const Eigen::VectorXd &areas)
WindingNumbersShape(ConstAlias< Eigen::MatrixXd > points, ConstAlias< Eigen::MatrixXd > normals, bool skipPointAreas=false)
std::vector< std::vector< int > > myO_PI
libIGL octree for fast queries data structure
CountedConstPtrOrConstPtr< Eigen::MatrixXd > myNormals
Const alias to the normals.
const Point aPoint(3, 4)
PointVector< 3, double > RealPoint