point_areas_and_normals.cpp 7.0 KB

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