smooth_corner_adjacency.cpp 3.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148
  1. #include "smooth_corner_adjacency.h"
  2. #include "vertex_triangle_adjacency.h"
  3. #include "matlab_format.h"
  4. #include "parallel_for.h"
  5. #include "unzip_corners.h"
  6. #include <iostream>
  7. void igl::smooth_corner_adjacency(
  8. const Eigen::MatrixXd & V,
  9. const Eigen::MatrixXi & F,
  10. const double corner_threshold_radians,
  11. Eigen::VectorXi & CI,
  12. Eigen::VectorXi & CC)
  13. {
  14. typedef double Scalar;
  15. typedef Eigen::Index Index;
  16. Eigen::Matrix<Index,Eigen::Dynamic,1> VF,NI;
  17. igl::vertex_triangle_adjacency(F,V.rows(),VF,NI);
  18. // unit normals
  19. Eigen::Matrix<Scalar,Eigen::Dynamic,3,Eigen::RowMajor> FN(F.rows(),3);
  20. igl::parallel_for(F.rows(),[&](const Index f)
  21. {
  22. const Eigen::Matrix<Scalar,1,3> v10 = V.row(F(f,1))-V.row(F(f,0));
  23. const Eigen::Matrix<Scalar,1,3> v20 = V.row(F(f,2))-V.row(F(f,0));
  24. const Eigen::Matrix<Scalar,1,3> n = v10.cross(v20);
  25. const Scalar a = n.norm();
  26. FN.row(f) = n/a;
  27. },10000);
  28. // number of faces
  29. const Index m = F.rows();
  30. // valence of faces
  31. const Index n = F.cols();
  32. assert(n == 3);
  33. CI.resize(m*n*8);
  34. CI.setConstant(-1);
  35. Index ncc = 0;
  36. Index ci = -1;
  37. // assumes that ci is strictly increasing and we're appending to CI
  38. const auto append_CI = [&](Index nf)
  39. {
  40. // make room
  41. if(ncc >= CI.size()) { CI.conservativeResize(CI.size()*2+1); }
  42. CI(ncc++) = nf;
  43. CC(ci+1)++;
  44. };
  45. CC.resize(m*3+1);
  46. CC.setConstant(-1);
  47. CC(0) = 0;
  48. const Scalar cos_thresh = cos(corner_threshold_radians);
  49. // parallelizing this probably requires map-reduce
  50. for(Index i = 0;i<m;i++)
  51. {
  52. // Normal of this face
  53. const auto & fnhat = FN.row(i);
  54. // loop over corners
  55. for(Index j = 0;j<n;j++)
  56. {
  57. // increment ci
  58. ci++;
  59. assert(ci == i*n+j);
  60. CC(ci+1) = CC(ci);
  61. const auto & v = F(i,j);
  62. for(int k = NI(v); k<NI(v+1); k++)
  63. {
  64. const int nf = VF(k);
  65. const auto & ifn = FN.row(nf);
  66. // dot product between face's normal and other face's normal
  67. const Scalar dp = fnhat.dot(ifn);
  68. // if difference in normal is slight then add to average
  69. if(dp > cos_thresh)
  70. {
  71. append_CI(nf);
  72. }
  73. }
  74. }
  75. }
  76. CI.conservativeResize(ncc);
  77. }
  78. void igl::smooth_corner_adjacency(
  79. const Eigen::MatrixXi & FV,
  80. const Eigen::MatrixXi & FN,
  81. Eigen::VectorXi & CI,
  82. Eigen::VectorXi & CC)
  83. {
  84. typedef double Scalar;
  85. typedef Eigen::Index Index;
  86. assert(FV.rows() == FN.rows());
  87. assert(FV.cols() == 3);
  88. assert(FN.cols() == 3);
  89. Eigen::VectorXi J;
  90. Index nu = -1;
  91. {
  92. Eigen::MatrixXi U;
  93. Eigen::MatrixXi _;
  94. igl::unzip_corners<const Eigen::MatrixXi>({FV,FN},U,_,J);
  95. nu = U.rows();
  96. assert(J.maxCoeff() == nu-1);
  97. }
  98. // could use linear arrays here if every becomes bottleneck
  99. std::vector<std::vector<Index>> U2F(nu);
  100. const Index m = FV.rows();
  101. for(Index j = 0;j<3;j++)
  102. {
  103. for(Index i = 0;i<m;i++)
  104. {
  105. // unzip_corners uses convention i+j*#F
  106. const Index ci = i+j*m;
  107. U2F[J(ci)].emplace_back(i);
  108. }
  109. }
  110. CC.resize(m*3+1);
  111. CC.setConstant(-1);
  112. CC(0) = 0;
  113. CI.resize(m*3*8);
  114. int ncc = 0;
  115. // assumes that ci is strictly increasing and we're appending to CI
  116. const auto append_CI = [&CI,&CC,&ncc](Index ci, Index nf)
  117. {
  118. // make room
  119. if(ncc >= CI.size()) { CI.conservativeResize(CI.size()*2+1); }
  120. CI(ncc++) = nf;
  121. CC(ci+1)++;
  122. };
  123. for(Index i = 0;i<m;i++)
  124. {
  125. for(Index j = 0;j<3;j++)
  126. {
  127. const Index J_ci = i+j*m;
  128. // CI,CC uses convention i*3 + j
  129. const Index C_ci = i*3+j;
  130. CC(C_ci+1) = CC(C_ci);
  131. for(const auto & nf : U2F[J(J_ci)])
  132. {
  133. append_CI(C_ci,nf);
  134. }
  135. }
  136. }
  137. CI.conservativeResize(ncc);
  138. }
  139. #ifdef IGL_STATIC_LIBRARY
  140. // Explicit template instantiations
  141. #endif