point_areas.cpp 7.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181
  1. #include "point_areas.h"
  2. #include "delaunay_triangulation.h"
  3. #include "../../colon.h"
  4. #include "../../slice.h"
  5. #include "../../slice_mask.h"
  6. #include "../../parallel_for.h"
  7. #include "CGAL/Exact_predicates_inexact_constructions_kernel.h"
  8. #include "CGAL/Triangulation_vertex_base_with_info_2.h"
  9. #include "CGAL/Triangulation_data_structure_2.h"
  10. #include "CGAL/Delaunay_triangulation_2.h"
  11. typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
  12. typedef CGAL::Triangulation_vertex_base_with_info_2<unsigned int, Kernel> Vb;
  13. typedef CGAL::Triangulation_data_structure_2<Vb> Tds;
  14. typedef CGAL::Delaunay_triangulation_2<Kernel, Tds> Delaunay;
  15. typedef Kernel::Point_2 Point;
  16. namespace igl {
  17. namespace copyleft{
  18. namespace cgal{
  19. template <typename DerivedP, typename DerivedI, typename DerivedN,
  20. typename DerivedA>
  21. IGL_INLINE void point_areas(
  22. const Eigen::MatrixBase<DerivedP>& P,
  23. const Eigen::MatrixBase<DerivedI>& I,
  24. const Eigen::MatrixBase<DerivedN>& N,
  25. Eigen::PlainObjectBase<DerivedA> & A)
  26. {
  27. Eigen::MatrixXd T;
  28. point_areas(P,I,N,A,T);
  29. }
  30. template <typename DerivedP, typename DerivedI, typename DerivedN,
  31. typename DerivedA, typename DerivedT>
  32. IGL_INLINE void point_areas(
  33. const Eigen::MatrixBase<DerivedP>& P,
  34. const Eigen::MatrixBase<DerivedI>& I,
  35. const Eigen::MatrixBase<DerivedN>& N,
  36. Eigen::PlainObjectBase<DerivedA> & A,
  37. Eigen::PlainObjectBase<DerivedT> & T)
  38. {
  39. typedef typename DerivedP::Scalar real;
  40. typedef typename DerivedN::Scalar scalarN;
  41. typedef typename DerivedA::Scalar scalarA;
  42. typedef Eigen::Matrix<real,1,3> RowVec3;
  43. typedef Eigen::Matrix<real,1,2> RowVec2;
  44. typedef Eigen::Matrix<real, Eigen::Dynamic, Eigen::Dynamic> MatrixP;
  45. typedef Eigen::Matrix<scalarN, Eigen::Dynamic, Eigen::Dynamic> MatrixN;
  46. typedef Eigen::Matrix<typename DerivedN::Scalar,
  47. Eigen::Dynamic, Eigen::Dynamic> VecotorO;
  48. typedef Eigen::Matrix<typename DerivedI::Scalar,
  49. Eigen::Dynamic, Eigen::Dynamic> MatrixI;
  50. const int n = P.rows();
  51. assert(P.cols() == 3 && "P must have exactly 3 columns");
  52. assert(P.rows() == N.rows()
  53. && "P and N must have the same number of rows");
  54. assert(P.rows() == I.rows()
  55. && "P and I must have the same number of rows");
  56. A.setZero(n,1);
  57. T.setZero(n,3);
  58. igl::parallel_for(P.rows(),[&](int i)
  59. {
  60. MatrixI neighbor_index = I.row(i);
  61. MatrixP neighbors;
  62. igl::slice(P,neighbor_index,1,neighbors);
  63. if(N.rows() && neighbors.rows() > 1){
  64. MatrixN neighbor_normals;
  65. igl::slice(N,neighbor_index,1,neighbor_normals);
  66. Eigen::Matrix<scalarN,1,3> poi_normal = neighbor_normals.row(0);
  67. Eigen::Matrix<scalarN,Eigen::Dynamic,1> dotprod =
  68. poi_normal(0)*neighbor_normals.col(0)
  69. + poi_normal(1)*neighbor_normals.col(1)
  70. + poi_normal(2)*neighbor_normals.col(2);
  71. Eigen::Array<bool,Eigen::Dynamic,1> keep = dotprod.array() > 0;
  72. igl::slice_mask(Eigen::MatrixXd(neighbors),keep,1,neighbors);
  73. }
  74. if(neighbors.rows() <= 2){
  75. A(i) = 0;
  76. } else {
  77. //subtract the mean from neighbors, then take svd,
  78. //the scores will be U*S. This is our pca plane fitting
  79. RowVec3 mean = neighbors.colwise().mean();
  80. MatrixP mean_centered = neighbors.rowwise() - mean;
  81. Eigen::JacobiSVD<MatrixP> svd(mean_centered,
  82. Eigen::ComputeThinU | Eigen::ComputeThinV);
  83. MatrixP scores = svd.matrixU() * svd.singularValues().asDiagonal();
  84. T.row(i) = svd.matrixV().col(2).transpose();
  85. if(T.row(i).dot(N.row(i)) < 0){
  86. T.row(i) *= -1;
  87. }
  88. MatrixP plane;
  89. igl::slice(scores,igl::colon<int>(0,scores.rows()-1),
  90. igl::colon<int>(0,1),plane);
  91. std::vector< std::pair<Point,unsigned> > points;
  92. //This is where we obtain a delaunay triangulation of the points
  93. for(unsigned iter = 0; iter < plane.rows(); iter++){
  94. points.push_back( std::make_pair(
  95. Point(plane(iter,0),plane(iter,1)), iter ) );
  96. }
  97. Delaunay triangulation;
  98. triangulation.insert(points.begin(),points.end());
  99. Eigen::MatrixXi F(triangulation.number_of_faces(),3);
  100. int f_row = 0;
  101. for(Delaunay::Finite_faces_iterator fit =
  102. triangulation.finite_faces_begin();
  103. fit != triangulation.finite_faces_end(); ++fit) {
  104. Delaunay::Face_handle face = fit;
  105. F.row(f_row) = Eigen::RowVector3i((int)face->vertex(0)->info(),
  106. (int)face->vertex(1)->info(),
  107. (int)face->vertex(2)->info());
  108. f_row++;
  109. }
  110. //Here we calculate the voronoi area of the point
  111. scalarA area_accumulator = 0;
  112. for(int f = 0; f < F.rows(); f++){
  113. int X = -1;
  114. for(int face_iter = 0; face_iter < 3; face_iter++){
  115. if(F(f,face_iter) == 0){
  116. X = face_iter;
  117. }
  118. }
  119. if(X >= 0){
  120. //Triangle XYZ with X being the point we want the area of
  121. int Y = (X+1)%3;
  122. int Z = (X+2)%3;
  123. scalarA x = (plane.row(F(f,Y))-plane.row(F(f,Z))).norm();
  124. scalarA y = (plane.row(F(f,X))-plane.row(F(f,Z))).norm();
  125. scalarA z = (plane.row(F(f,Y))-plane.row(F(f,X))).norm();
  126. scalarA cosX = (z*z + y*y - x*x)/(2*y*z);
  127. scalarA cosY = (z*z + x*x - y*y)/(2*x*z);
  128. scalarA cosZ = (x*x + y*y - z*z)/(2*y*x);
  129. Eigen::Matrix<scalarA,1,3> barycentric;
  130. barycentric << x*cosX, y*cosY, z*cosZ;
  131. barycentric /= (barycentric(0)+barycentric(1)+barycentric(2));
  132. //TODO: to make numerically stable, reorder so that x≥y≥z:
  133. scalarA full_area = 0.25*std::sqrt(
  134. (x+(y+z))*(z-(x-y))*(z+(x-y))*(x+(y-z)));
  135. Eigen::Matrix<scalarA,1,3> partial_area =
  136. barycentric * full_area;
  137. if(cosX < 0){
  138. area_accumulator += 0.5*full_area;
  139. } else if (cosY < 0 || cosZ < 0){
  140. area_accumulator += 0.25*full_area;
  141. } else {
  142. area_accumulator += (partial_area(1) + partial_area(2))/2;
  143. }
  144. }
  145. }
  146. if(std::isfinite(area_accumulator)){
  147. A(i) = area_accumulator;
  148. } else {
  149. A(i) = 0;
  150. }
  151. }
  152. },1000);
  153. }
  154. }
  155. }
  156. }