fast_winding_number.h 8.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213
  1. #ifndef IGL_FAST_WINDING_NUMBER
  2. #define IGL_FAST_WINDING_NUMBER
  3. #include "igl_inline.h"
  4. #include "FastWindingNumberForSoups.h"
  5. #include <Eigen/Core>
  6. #include <vector>
  7. namespace igl
  8. {
  9. /// Generate the precomputation for the fast winding number for point data
  10. /// [Barill et. al 2018].
  11. ///
  12. /// Given a set of 3D points P, with normals N, areas A, along with octree
  13. /// data, and an expansion order, we define a taylor series expansion at each
  14. /// octree cell.
  15. ///
  16. /// The octree data is designed to come from igl::octree, and the areas (if not
  17. /// obtained at scan time), may be calculated using
  18. /// igl::copyleft::cgal::point_areas.
  19. ///
  20. /// @param[in] P #P by 3 list of point locations
  21. /// @param[in] N #P by 3 list of point normals
  22. /// @param[in] A #P by 1 list of point areas
  23. /// @param[in] point_indices a vector of vectors, where the ith entry is a vector of
  24. /// the indices into P that are the ith octree cell's points
  25. /// @param[in] CH #OctreeCells by 8, where the ith row is the indices of
  26. /// the ith octree cell's children
  27. /// @param[in] expansion_order the order of the taylor expansion. We support 0,1,2.
  28. /// @param[out] CM #OctreeCells by 3 list of each cell's center of mass
  29. /// @param[out] R #OctreeCells by 1 list of each cell's maximum distance of any point
  30. /// to the center of mass
  31. /// @param[out] EC #OctreeCells by #TaylorCoefficients list of expansion coefficients.
  32. /// (Note that #TaylorCoefficients = ∑_{i=1}^{expansion_order} 3^i)
  33. ///
  34. /// \see copyleft::cgal::point_areas, knn
  35. template <
  36. typename DerivedP,
  37. typename DerivedA,
  38. typename DerivedN,
  39. typename Index,
  40. typename DerivedCH,
  41. typename DerivedCM,
  42. typename DerivedR,
  43. typename DerivedEC>
  44. IGL_INLINE void fast_winding_number(
  45. const Eigen::MatrixBase<DerivedP>& P,
  46. const Eigen::MatrixBase<DerivedN>& N,
  47. const Eigen::MatrixBase<DerivedA>& A,
  48. const std::vector<std::vector<Index> > & point_indices,
  49. const Eigen::MatrixBase<DerivedCH>& CH,
  50. const int expansion_order,
  51. Eigen::PlainObjectBase<DerivedCM>& CM,
  52. Eigen::PlainObjectBase<DerivedR>& R,
  53. Eigen::PlainObjectBase<DerivedEC>& EC);
  54. /// Evaluate the fast winding number for point data, having already done the
  55. /// the precomputation
  56. ///
  57. /// @param[in] P #P by 3 list of point locations
  58. /// @param[in] N #P by 3 list of point normals
  59. /// @param[in] A #P by 1 list of point areas
  60. /// @param[in] point_indices a vector of vectors, where the ith entry is a vector of
  61. /// the indices into P that are the ith octree cell's points
  62. /// @param[in] CH #OctreeCells by 8, where the ith row is the indices of
  63. /// the ith octree cell's children
  64. /// @param[in] CM #OctreeCells by 3 list of each cell's center of mass
  65. /// @param[in] R #OctreeCells by 1 list of each cell's maximum distance of any point
  66. /// to the center of mass
  67. /// @param[in] EC #OctreeCells by #TaylorCoefficients list of expansion coefficients.
  68. /// (Note that #TaylorCoefficients = ∑_{i=1}^{expansion_order} 3^i)
  69. /// @param[in] Q #Q by 3 list of query points for the winding number
  70. /// @param[in] beta This is a Barnes-Hut style accuracy term that separates near feild
  71. /// from far field. The higher the beta, the more accurate and slower
  72. /// the evaluation. We reccommend using a beta value of 2. Note that
  73. /// for a beta value ≤ 0, we use the direct evaluation, rather than
  74. /// the fast approximation
  75. /// @param[out] WN #Q by 1 list of windinng number values at each query point
  76. template <
  77. typename DerivedP,
  78. typename DerivedA,
  79. typename DerivedN,
  80. typename Index,
  81. typename DerivedCH,
  82. typename DerivedCM,
  83. typename DerivedR,
  84. typename DerivedEC,
  85. typename DerivedQ,
  86. typename BetaType,
  87. typename DerivedWN>
  88. IGL_INLINE void fast_winding_number(
  89. const Eigen::MatrixBase<DerivedP>& P,
  90. const Eigen::MatrixBase<DerivedN>& N,
  91. const Eigen::MatrixBase<DerivedA>& A,
  92. const std::vector<std::vector<Index> > & point_indices,
  93. const Eigen::MatrixBase<DerivedCH>& CH,
  94. const Eigen::MatrixBase<DerivedCM>& CM,
  95. const Eigen::MatrixBase<DerivedR>& R,
  96. const Eigen::MatrixBase<DerivedEC>& EC,
  97. const Eigen::MatrixBase<DerivedQ>& Q,
  98. const BetaType beta,
  99. Eigen::PlainObjectBase<DerivedWN>& WN);
  100. /// \overload
  101. ///
  102. /// \brief Evaluate the fast winding number for point data without caching the
  103. /// precomputation.
  104. template <
  105. typename DerivedP,
  106. typename DerivedA,
  107. typename DerivedN,
  108. typename DerivedQ,
  109. typename BetaType,
  110. typename DerivedWN>
  111. IGL_INLINE void fast_winding_number(
  112. const Eigen::MatrixBase<DerivedP>& P,
  113. const Eigen::MatrixBase<DerivedN>& N,
  114. const Eigen::MatrixBase<DerivedA>& A,
  115. const Eigen::MatrixBase<DerivedQ>& Q,
  116. const int expansion_order,
  117. const BetaType beta,
  118. Eigen::PlainObjectBase<DerivedWN>& WN);
  119. /// \overload
  120. template <
  121. typename DerivedP,
  122. typename DerivedA,
  123. typename DerivedN,
  124. typename DerivedQ,
  125. typename DerivedWN>
  126. IGL_INLINE void fast_winding_number(
  127. const Eigen::MatrixBase<DerivedP>& P,
  128. const Eigen::MatrixBase<DerivedN>& N,
  129. const Eigen::MatrixBase<DerivedA>& A,
  130. const Eigen::MatrixBase<DerivedQ>& Q,
  131. Eigen::PlainObjectBase<DerivedWN>& WN);
  132. /// @private
  133. namespace FastWindingNumber {
  134. /// @private
  135. namespace HDK_Sample{
  136. /// @private
  137. template <typename T1, typename T2> class UT_SolidAngle;} }
  138. /// Structure for caching precomputation for fast winding number for triangle
  139. /// soups
  140. struct FastWindingNumberBVH {
  141. /// @private
  142. FastWindingNumber::HDK_Sample::UT_SolidAngle<float,float> ut_solid_angle;
  143. // Need copies of these so they stay alive between calls.
  144. /// @private
  145. std::vector<FastWindingNumber::HDK_Sample::UT_Vector3T<float> > U;
  146. std::vector<int> F;
  147. };
  148. /// Compute approximate winding number of a triangle soup mesh according to
  149. /// "Fast Winding Numbers for Soups and Clouds" [Barill et al. 2018].
  150. ///
  151. /// @param[in] V #V by 3 list of mesh vertex positions
  152. /// @param[in] F #F by 3 list of triangle mesh indices into rows of V
  153. /// @param[in] Q #Q by 3 list of query positions
  154. /// @param[out] W #Q list of winding number values
  155. template <
  156. typename DerivedV,
  157. typename DerivedF,
  158. typename DerivedQ,
  159. typename DerivedW>
  160. IGL_INLINE void fast_winding_number(
  161. const Eigen::MatrixBase<DerivedV> & V,
  162. const Eigen::MatrixBase<DerivedF> & F,
  163. const Eigen::MatrixBase<DerivedQ> & Q,
  164. Eigen::PlainObjectBase<DerivedW> & W);
  165. /// Precomputation for computing approximate winding numbers of a triangle
  166. /// soup.
  167. ///
  168. /// @param[in] V #V by 3 list of mesh vertex positions
  169. /// @param[in] F #F by 3 list of triangle mesh indices into rows of V
  170. /// @param[in] order Taylor series expansion order to use (e.g., 2)
  171. /// @param[out] fwn_bvh Precomputed bounding volume hierarchy
  172. ///
  173. template <
  174. typename DerivedV,
  175. typename DerivedF>
  176. IGL_INLINE void fast_winding_number(
  177. const Eigen::MatrixBase<DerivedV> & V,
  178. const Eigen::MatrixBase<DerivedF> & F,
  179. const int order,
  180. FastWindingNumberBVH & fwn_bvh);
  181. /// After precomputation, compute winding number at a each of many points in a
  182. /// list.
  183. ///
  184. /// @param[in] fwn_bvh Precomputed bounding volume hierarchy
  185. /// @param[in] accuracy_scale parameter controlling accuracy (e.g., 2)
  186. /// @param[in] Q #Q by 3 list of query positions
  187. /// @param[out] W #Q list of winding number values
  188. template <
  189. typename DerivedQ,
  190. typename DerivedW>
  191. IGL_INLINE void fast_winding_number(
  192. const FastWindingNumberBVH & fwn_bvh,
  193. const float accuracy_scale,
  194. const Eigen::MatrixBase<DerivedQ> & Q,
  195. Eigen::PlainObjectBase<DerivedW> & W);
  196. /// After precomputation, compute winding number at a single point
  197. ///
  198. /// @param[in] fwn_bvh Precomputed bounding volume hierarchy
  199. /// @param[in] accuracy_scale parameter controlling accuracy (e.g., 2)
  200. /// @param[in] p single position
  201. /// @return w winding number of this point
  202. template <typename Derivedp>
  203. IGL_INLINE typename Derivedp::Scalar fast_winding_number(
  204. const FastWindingNumberBVH & fwn_bvh,
  205. const float accuracy_scale,
  206. const Eigen::MatrixBase<Derivedp> & p);
  207. }
  208. #ifndef IGL_STATIC_LIBRARY
  209. # include "fast_winding_number.cpp"
  210. #endif
  211. #endif