orient_outward_ao.cpp 6.5 KB

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