orient_outward_ao.cpp 6.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217
  1. #include "orient_outward_ao.h"
  2. #include "../per_face_normals.h"
  3. #include "../doublearea.h"
  4. #include "../random_dir.h"
  5. #include "EmbreeIntersector.h"
  6. #include <iostream>
  7. #include <random>
  8. #include <ctime>
  9. #include <limits>
  10. template <
  11. typename DerivedV,
  12. typename DerivedF,
  13. typename DerivedC,
  14. typename DerivedFF,
  15. typename DerivedI>
  16. IGL_INLINE void igl::orient_outward_ao(
  17. const Eigen::PlainObjectBase<DerivedV> & V,
  18. const Eigen::PlainObjectBase<DerivedF> & F,
  19. const Eigen::PlainObjectBase<DerivedC> & C,
  20. const int min_num_rays_per_component,
  21. const int total_num_rays,
  22. Eigen::PlainObjectBase<DerivedFF> & FF,
  23. Eigen::PlainObjectBase<DerivedI> & I)
  24. {
  25. using namespace Eigen;
  26. using namespace std;
  27. assert(C.rows() == F.rows());
  28. assert(F.cols() == 3);
  29. assert(V.cols() == 3);
  30. EmbreeIntersector ei;
  31. ei.init(V.template cast<float>(),F);
  32. // number of faces
  33. const int m = F.rows();
  34. // number of patches
  35. const int num_cc = C.maxCoeff()+1;
  36. I.resize(num_cc);
  37. if(&FF != &F)
  38. {
  39. FF = F;
  40. }
  41. // face normal
  42. MatrixXd N;
  43. per_face_normals(V,F,N);
  44. // face area
  45. Matrix<typename DerivedV::Scalar,Dynamic,1> A;
  46. doublearea(V,F,A);
  47. double area_min = numeric_limits<double>::max();
  48. for (int f = 0; f < m; ++f)
  49. {
  50. area_min = A(f) != 0 && A(f) < area_min ? A(f) : area_min;
  51. }
  52. double area_total = A.sum();
  53. // determine number of rays per component according to its area
  54. VectorXd area_per_component;
  55. area_per_component.setZero(num_cc);
  56. for (int f = 0; f < m; ++f)
  57. {
  58. area_per_component(C(f)) += A(f);
  59. }
  60. VectorXi num_rays_per_component;
  61. num_rays_per_component.setZero(num_cc);
  62. for (int c = 0; c < num_cc; ++c)
  63. {
  64. num_rays_per_component(c) = max<int>(min_num_rays_per_component, static_cast<int>(total_num_rays * area_per_component(c) / area_total));
  65. }
  66. // generate all the rays
  67. //cout << "generating rays... ";
  68. uniform_real_distribution<float> rdist;
  69. mt19937 prng;
  70. prng.seed(time(nullptr));
  71. vector<int > ray_face;
  72. vector<Vector3f> ray_ori;
  73. vector<Vector3f> ray_dir;
  74. ray_face.reserve(total_num_rays);
  75. ray_ori .reserve(total_num_rays);
  76. ray_dir .reserve(total_num_rays);
  77. for (int c = 0; c < num_cc; ++c)
  78. {
  79. if (area_per_component[c] == 0)
  80. {
  81. continue;
  82. }
  83. vector<int> CF; // set of faces per component
  84. vector<int> CF_area;
  85. for (int f = 0; f < m; ++f)
  86. {
  87. if (C(f)==c)
  88. {
  89. CF.push_back(f);
  90. CF_area.push_back(static_cast<int>(100 * A(f) / area_min));
  91. }
  92. }
  93. // discrete distribution for random selection of faces with probability proportional to their areas
  94. auto ddist_func = [&] (double i) { return CF_area[static_cast<int>(i)]; };
  95. discrete_distribution<int> ddist(CF.size(), 0, CF.size(), ddist_func); // simple ctor of (Iter, Iter) not provided by the stupid VC11 impl...
  96. for (int i = 0; i < num_rays_per_component[c]; ++i)
  97. {
  98. int f = CF[ddist(prng)]; // select face with probability proportional to face area
  99. float t0 = rdist(prng); // random barycentric coordinate
  100. float t1 = rdist(prng);
  101. float t2 = rdist(prng);
  102. float t_sum = t0 + t1 + t2;
  103. t0 /= t_sum;
  104. t1 /= t_sum;
  105. t2 /= t_sum;
  106. Vector3f p = t0 * V.row(F(f,0)).template cast<float>().eval() // be careful with the index!!!
  107. + t1 * V.row(F(f,1)).template cast<float>().eval()
  108. + t2 * V.row(F(f,2)).template cast<float>().eval();
  109. Vector3f n = N.row(f).cast<float>();
  110. // random direction in hemisphere around n (avoid too grazing angle)
  111. Vector3f d;
  112. while (true) {
  113. d = random_dir().cast<float>();
  114. float ndotd = n.dot(d);
  115. if (fabsf(ndotd) < 0.1f)
  116. {
  117. continue;
  118. }
  119. if (ndotd < 0)
  120. {
  121. d *= -1.0f;
  122. }
  123. break;
  124. }
  125. ray_face.push_back(f);
  126. ray_ori .push_back(p);
  127. ray_dir .push_back(d);
  128. }
  129. }
  130. // per component voting: first=front, second=back
  131. vector<pair<float, float>> C_vote_distance(num_cc, make_pair(0, 0)); // sum of distance between ray origin and intersection
  132. vector<pair<int , int >> C_vote_infinity(num_cc, make_pair(0, 0)); // number of rays reaching infinity
  133. //cout << "shooting rays... ";
  134. #pragma omp parallel for
  135. for (int i = 0; i < (int)ray_face.size(); ++i)
  136. {
  137. int f = ray_face[i];
  138. Vector3f o = ray_ori [i];
  139. Vector3f d = ray_dir [i];
  140. int c = C(f);
  141. // shoot ray toward front & back
  142. vector<Hit> hits_front;
  143. vector<Hit> hits_back;
  144. int num_rays_front;
  145. int num_rays_back;
  146. ei.intersectRay(o, d, hits_front, num_rays_front);
  147. ei.intersectRay(o, -d, hits_back , num_rays_back );
  148. if (!hits_front.empty() && hits_front[0].id == f) hits_front.erase(hits_front.begin());
  149. if (!hits_back .empty() && hits_back [0].id == f) hits_back .erase(hits_back .begin());
  150. if (hits_front.empty())
  151. {
  152. #pragma omp atomic
  153. C_vote_infinity[c].first++;
  154. } else {
  155. #pragma omp atomic
  156. C_vote_distance[c].first += hits_front[0].t;
  157. }
  158. if (hits_back.empty())
  159. {
  160. #pragma omp atomic
  161. C_vote_infinity[c].second++;
  162. } else {
  163. #pragma omp atomic
  164. C_vote_distance[c].second += hits_back[0].t;
  165. }
  166. }
  167. for(int c = 0;c<num_cc;c++)
  168. {
  169. I(c) = C_vote_infinity[c].first == C_vote_infinity[c].second &&
  170. C_vote_distance[c].first < C_vote_distance[c].second ||
  171. C_vote_infinity[c].first < C_vote_infinity[c].second;
  172. }
  173. // flip according to I
  174. for(int f = 0;f<m;f++)
  175. {
  176. if(I(C(f)))
  177. {
  178. FF.row(f) = FF.row(f).reverse().eval();
  179. }
  180. }
  181. //cout << "done!\n";
  182. }
  183. // Call with default parameters
  184. template <
  185. typename DerivedV,
  186. typename DerivedF,
  187. typename DerivedC,
  188. typename DerivedFF,
  189. typename DerivedI>
  190. IGL_INLINE void igl::orient_outward_ao(
  191. const Eigen::PlainObjectBase<DerivedV> & V,
  192. const Eigen::PlainObjectBase<DerivedF> & F,
  193. const Eigen::PlainObjectBase<DerivedC> & C,
  194. Eigen::PlainObjectBase<DerivedFF> & FF,
  195. Eigen::PlainObjectBase<DerivedI> & I)
  196. {
  197. return orient_outward_ao(V, F, C, 100, F.rows() * 100, FF, I);
  198. }
  199. #ifndef IGL_HEADER_ONLY
  200. // Explicit template specialization
  201. template void igl::orient_outward_ao<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int, int, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  202. #endif