orient_outward_ao.cpp 5.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154
  1. #include "orient_outward_ao.h"
  2. #include "../per_face_normals.h"
  3. #include "../barycenter.h"
  4. #include "../doublearea.h"
  5. #include "../matlab_format.h"
  6. #include "ambient_occlusion.h"
  7. #include "EmbreeIntersector.h"
  8. #include <iostream>
  9. #include <random>
  10. #include <omp.h>
  11. template <
  12. typename DerivedV,
  13. typename DerivedF,
  14. typename DerivedC,
  15. typename PointMatrixType,
  16. typename FaceMatrixType,
  17. typename RowVector3,
  18. typename DerivedFF,
  19. typename DerivedI>
  20. IGL_INLINE void igl::orient_outward_ao(
  21. const Eigen::PlainObjectBase<DerivedV> & V,
  22. const Eigen::PlainObjectBase<DerivedF> & F,
  23. const Eigen::PlainObjectBase<DerivedC> & C,
  24. const igl::EmbreeIntersector<PointMatrixType,FaceMatrixType,RowVector3> & ei,
  25. const int num_samples,
  26. Eigen::PlainObjectBase<DerivedFF> & FF,
  27. Eigen::PlainObjectBase<DerivedI> & I)
  28. {
  29. using namespace Eigen;
  30. using namespace std;
  31. assert(C.rows() == F.rows());
  32. assert(F.cols() == 3);
  33. assert(V.cols() == 3);
  34. // number of faces
  35. const int m = F.rows();
  36. // number of patches
  37. const int num_cc = C.maxCoeff()+1;
  38. I.resize(num_cc);
  39. if(&FF != &F)
  40. {
  41. FF = F;
  42. }
  43. // face normal
  44. PlainObjectBase<DerivedV> N;
  45. per_face_normals(V,F,N);
  46. // random number generator/distribution for each thread
  47. int max_threads = omp_get_max_threads();
  48. // prng
  49. vector<mt19937> engine(max_threads);
  50. for (int i = 0; i < max_threads; ++i)
  51. engine[i].seed(time(0) * (i + 1));
  52. // discrete distribution for random selection of faces with probability proportional to their areas
  53. Matrix<typename DerivedV::Scalar,Dynamic,1> A;
  54. doublearea(V,F,A);
  55. double minarea = A.minCoeff();
  56. Matrix<int, Dynamic, 1> A_int = (A * 100.0 / minarea).template cast<int>(); // only integer is allowed for weight
  57. auto ddist_func = [&] (double i) { return A_int(static_cast<int>(i)); };
  58. vector<discrete_distribution<int>> ddist(max_threads, discrete_distribution<int>(m, 0, m, ddist_func)); // simple ctor of (Iter, Iter) not provided by the stupid VC11 impl...
  59. // uniform real between in [0, 1]
  60. vector<uniform_real_distribution<double>> rdist(max_threads);
  61. //// occlusion count per component: +1 when front ray is occluded, -1 when back ray is occluded
  62. // occlussion count per component, per back/front: C(c,0) --> number of
  63. // front-side rays occluded, C(c,1) --> number of back-side rays occluded
  64. Matrix<int, Dynamic, 2> C_occlude_count;
  65. C_occlude_count.setZero(num_cc, 2);
  66. #pragma omp parallel for
  67. for (int i = 0; i < num_samples; ++i)
  68. {
  69. int thread_num = omp_get_thread_num();
  70. int f = ddist[thread_num](engine[thread_num]); // select face with probability proportional to face area
  71. double t0 = rdist[thread_num](engine[thread_num]);
  72. double t1 = rdist[thread_num](engine[thread_num]);
  73. double t2 = rdist[thread_num](engine[thread_num]);
  74. double t_sum = t0 + t1 + t2;
  75. t0 /= t_sum;
  76. t1 /= t_sum;
  77. t2 /= t_sum;
  78. RowVector3d p = t0 * V.row(F(f,0)) + t1 * V.row(F(f,1)) + t1 * V.row(F(f,2));
  79. RowVector3d n = N.row(f);
  80. //bool is_backside = rdist[thread_num](engine[thread_num]) < 0.5;
  81. // Loop over front or back side
  82. for(int s = 0;s<2;s++)
  83. {
  84. if(s==1)
  85. {
  86. n *= -1;
  87. }
  88. Matrix<typename DerivedV::Scalar,Dynamic,1> S;
  89. ambient_occlusion(ei, p, n, 1, S);
  90. if (S(0) > 0)
  91. {
  92. #pragma omp atomic
  93. C_occlude_count(C(f),s)++;
  94. }
  95. }
  96. }
  97. for(int c = 0;c<num_cc;c++)
  98. {
  99. //I(c) = C_occlude_count(c) > 0;
  100. I(c) = C_occlude_count(c,0) > C_occlude_count(c,1);
  101. }
  102. // flip according to I
  103. for(int f = 0;f<m;f++)
  104. {
  105. if(I(C(f)))
  106. {
  107. FF.row(f) = FF.row(f).reverse().eval();
  108. }
  109. }
  110. }
  111. // EmbreeIntersector generated on the fly
  112. template <
  113. typename DerivedV,
  114. typename DerivedF,
  115. typename DerivedC,
  116. typename DerivedFF,
  117. typename DerivedI>
  118. IGL_INLINE void igl::orient_outward_ao(
  119. const Eigen::PlainObjectBase<DerivedV> & V,
  120. const Eigen::PlainObjectBase<DerivedF> & F,
  121. const Eigen::PlainObjectBase<DerivedC> & C,
  122. const int num_samples,
  123. Eigen::PlainObjectBase<DerivedFF> & FF,
  124. Eigen::PlainObjectBase<DerivedI> & I)
  125. {
  126. using namespace igl;
  127. using namespace Eigen;
  128. // Both sides
  129. MatrixXi F2;
  130. F2.resize(F.rows()*2,F.cols());
  131. F2 << F, F.rowwise().reverse().eval();
  132. EmbreeIntersector<
  133. PlainObjectBase<DerivedV>,
  134. PlainObjectBase<DerivedF>,
  135. Matrix<typename DerivedV::Scalar,3,1> > ei(V,F2);
  136. return orient_outward_ao(V, F, C, ei, num_samples, FF, I);
  137. }
  138. #ifndef IGL_HEADER_ONLY
  139. // Explicit template specialization
  140. 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, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  141. #endif