winding_number.cpp 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216
  1. #include "winding_number.h"
  2. #include "WindingNumberAABB.h"
  3. #include <igl/PI.h>
  4. #include <cmath>
  5. IGL_INLINE void igl::winding_number(
  6. const Eigen::MatrixXd & V,
  7. const Eigen::MatrixXi & F,
  8. const Eigen::MatrixXd & O,
  9. Eigen::VectorXd & W)
  10. {
  11. using namespace Eigen;
  12. // make room for output
  13. W.resize(O.rows(),1);
  14. switch(F.cols())
  15. {
  16. case 2:
  17. return winding_number_2(
  18. V.data(),
  19. V.rows(),
  20. F.data(),
  21. F.rows(),
  22. O.data(),
  23. O.rows(),
  24. W.data());
  25. case 3:
  26. {
  27. WindingNumberAABB<Vector3d> hier(V,F);
  28. hier.grow();
  29. // loop over origins
  30. const int no = O.rows();
  31. # pragma omp parallel for if (no>IGL_WINDING_NUMBER_OMP_MIN_VALUE)
  32. for(int o = 0;o<no;o++)
  33. {
  34. Vector3d p = O.row(o);
  35. W(o) = hier.winding_number(p);
  36. }
  37. break;
  38. }
  39. default: assert(false && "Bad simplex size"); break;
  40. }
  41. }
  42. template <typename DerivedF>
  43. IGL_INLINE void igl::winding_number_3(
  44. const double * V,
  45. const int n,
  46. const DerivedF * F,
  47. const int m,
  48. const double * O,
  49. const int no,
  50. double * S)
  51. {
  52. // Initialize output
  53. // loop over origins
  54. #pragma omp parallel for if (no>IGL_WINDING_NUMBER_OMP_MIN_VALUE)
  55. for(int o = 0;o<no;o++)
  56. {
  57. S[o] = 0;
  58. }
  59. // Only use parallel for if there are many facets and more than one origin.
  60. // Assumes that if there is exactly one origin then this is being called
  61. // within an outer for loop which may be parallel
  62. #pragma omp parallel for if (m>IGL_WINDING_NUMBER_OMP_MIN_VALUE && no>1)
  63. // loop over faces
  64. for(int f = 0;f<m;f++)
  65. {
  66. // Gather corners
  67. double C[3][3];
  68. // loop around triangle
  69. for(int t=0;t<3;t++)
  70. {
  71. // loop over dimensions
  72. for(int d = 0;d<3;d++)
  73. {
  74. // Indices are offset by 1
  75. int Ff = F[m*t + f];
  76. C[t][d] = V[d*n + Ff];
  77. }
  78. }
  79. // loop over origins
  80. for(int o = 0;o<no;o++)
  81. {
  82. // Gather vectors to corners
  83. double v[3][3];
  84. double vl[3];
  85. // loop around triangle
  86. for(int t=0;t<3;t++)
  87. {
  88. vl[t] = 0;
  89. // loop over dimensions
  90. for(int d = 0;d<3;d++)
  91. {
  92. v[t][d] = C[t][d] - O[d*no + o];
  93. // compute edge length contribution
  94. vl[t] += v[t][d]*v[t][d];
  95. }
  96. // finish edge length computation
  97. // Matlab crashes on NaN
  98. if(vl[t]!=0)
  99. {
  100. vl[t] /= sqrt(vl[t]);
  101. }
  102. }
  103. //printf("\n");
  104. // Compute determinant
  105. double detf =
  106. v[0][0]*v[1][1]*v[2][2]+
  107. v[1][0]*v[2][1]*v[0][2]+
  108. v[2][0]*v[0][1]*v[1][2]-
  109. v[2][0]*v[1][1]*v[0][2]-
  110. v[1][0]*v[0][1]*v[2][2]-
  111. v[0][0]*v[2][1]*v[1][2];
  112. // Compute pairwise dotproducts
  113. double dp[3];
  114. dp[0] = v[1][0]*v[2][0];
  115. dp[0] += v[1][1]*v[2][1];
  116. dp[0] += v[1][2]*v[2][2];
  117. dp[1] = v[2][0]*v[0][0];
  118. dp[1] += v[2][1]*v[0][1];
  119. dp[1] += v[2][2]*v[0][2];
  120. dp[2] = v[0][0]*v[1][0];
  121. dp[2] += v[0][1]*v[1][1];
  122. dp[2] += v[0][2]*v[1][2];
  123. // Compute winding number
  124. // Only divide by TWO_PI instead of 4*pi because there was a 2 out front
  125. double val = atan2(detf,
  126. vl[0]*vl[1]*vl[2] +
  127. dp[0]*vl[0] +
  128. dp[1]*vl[1] +
  129. dp[2]*vl[2]) / (2.*igl::PI);
  130. #pragma omp atomic
  131. S[o] += val;
  132. }
  133. }
  134. }
  135. template <typename DerivedF>
  136. IGL_INLINE void igl::winding_number_2(
  137. const double * V,
  138. const int n,
  139. const DerivedF * F,
  140. const int m,
  141. const double * O,
  142. const int no,
  143. double * S)
  144. {
  145. // Initialize output
  146. // loop over origins
  147. #pragma omp parallel for if (no>IGL_WINDING_NUMBER_OMP_MIN_VALUE)
  148. for(int o = 0;o<no;o++)
  149. {
  150. S[o] = 0;
  151. }
  152. #pragma omp parallel for if (m>IGL_WINDING_NUMBER_OMP_MIN_VALUE && no>1)
  153. // loop over faces
  154. for(int f = 0;f<m;f++)
  155. {
  156. // Index of source and destination
  157. int s = F[m*0 + f];
  158. int d = F[m*1 + f];
  159. // Positions of source and destination
  160. double vs[2];
  161. double vd[2];
  162. vs[0] = V[0*n + s];
  163. vs[1] = V[1*n + s];
  164. vd[0] = V[0*n + d];
  165. vd[1] = V[1*n + d];
  166. // loop over origins
  167. for(int o = 0;o<no;o++)
  168. {
  169. // Gather vectors to source and destination
  170. double o2vs[2];
  171. double o2vd[2];
  172. // and lengths
  173. double o2vsl = 0;
  174. double o2vdl = 0;
  175. for(int i = 0;i<2;i++)
  176. {
  177. o2vs[i] = O[i*no + o] - vs[i];
  178. o2vd[i] = O[i*no + o] - vd[i];
  179. o2vsl += o2vs[i]*o2vs[i];
  180. o2vdl += o2vd[i]*o2vd[i];
  181. }
  182. o2vsl = sqrt(o2vsl);
  183. o2vdl = sqrt(o2vdl);
  184. // Normalize
  185. for(int i = 0;i<2;i++)
  186. {
  187. // Matlab crashes on NaN
  188. if(o2vsl!=0)
  189. {
  190. o2vs[i] /= o2vsl;
  191. }
  192. if(o2vdl!=0)
  193. {
  194. o2vd[i] /= o2vdl;
  195. }
  196. }
  197. double val =
  198. atan2(o2vd[0]*o2vs[1]-o2vd[1]*o2vs[0],o2vd[0]*o2vs[0]+o2vd[1]*o2vs[1])/
  199. (2.*igl::PI);
  200. #pragma omp atomic
  201. S[o] += val;
  202. }
  203. }
  204. }
  205. #ifdef IGL_STATIC_LIBRARY
  206. // Explicit template specialization
  207. template void igl::winding_number_3<int>(double const*, int, int const*, int, double const*, int, double*);
  208. template void igl::winding_number_2<double>(double const*, int, double const*, int, double const*, int, double*);
  209. template void igl::winding_number_3<double>(double const*, int, double const*, int, double const*, int, double*);
  210. #endif