point_mesh_squared_distance.cpp 9.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <[email protected]>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "point_mesh_squared_distance.h"
  9. #include "mesh_to_cgal_triangle_list.h"
  10. #include "assign_scalar.h"
  11. #include "../../parallel_for.h"
  12. template <
  13. typename Kernel,
  14. typename DerivedP,
  15. typename DerivedV,
  16. typename DerivedF,
  17. typename DerivedsqrD,
  18. typename DerivedI,
  19. typename DerivedC>
  20. IGL_INLINE void igl::copyleft::cgal::point_mesh_squared_distance(
  21. const Eigen::MatrixBase<DerivedP> & P,
  22. const Eigen::MatrixBase<DerivedV> & V,
  23. const Eigen::MatrixBase<DerivedF> & F,
  24. Eigen::PlainObjectBase<DerivedsqrD> & sqrD,
  25. Eigen::PlainObjectBase<DerivedI> & I,
  26. Eigen::PlainObjectBase<DerivedC> & C)
  27. {
  28. using namespace std;
  29. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  30. typedef typename std::vector<Triangle_3>::iterator Iterator;
  31. typedef CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  32. typedef CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  33. typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
  34. Tree tree;
  35. vector<Triangle_3> T;
  36. point_mesh_squared_distance_precompute(V,F,tree,T);
  37. return point_mesh_squared_distance(P,tree,T,sqrD,I,C);
  38. }
  39. template <typename Kernel, typename DerivedV, typename DerivedF>
  40. IGL_INLINE void igl::copyleft::cgal::point_mesh_squared_distance_precompute(
  41. const Eigen::MatrixBase<DerivedV> & V,
  42. const Eigen::MatrixBase<DerivedF> & F,
  43. CGAL::AABB_tree<
  44. CGAL::AABB_traits<Kernel,
  45. CGAL::AABB_triangle_primitive<Kernel,
  46. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  47. >
  48. >
  49. > & tree,
  50. std::vector<CGAL::Triangle_3<Kernel> > & T)
  51. {
  52. using namespace std;
  53. typedef CGAL::Point_3<Kernel> Point_3;
  54. // Must be 3D
  55. assert(V.cols() == 3);
  56. // Must be triangles
  57. assert(F.cols() == 3);
  58. // WTF ALERT!!!!
  59. //
  60. // There's a bug in clang probably or at least in cgal. Without calling this
  61. // line (I guess invoking some compilation from <vector>), clang will vomit
  62. // errors inside CGAL.
  63. //
  64. // http://stackoverflow.com/questions/27748442/is-clangs-c11-support-reliable
  65. T.reserve(0);
  66. // Make list of cgal triangles
  67. mesh_to_cgal_triangle_list(V,F,T);
  68. tree.clear();
  69. tree.insert(T.begin(),T.end());
  70. tree.accelerate_distance_queries();
  71. // accelerate_distance_queries doesn't seem actually to do _all_ of the
  72. // precomputation. the tree (despite being const) will still do more
  73. // precomputation and reorganizing on the first call of `closest_point` or
  74. // `closest_point_and_primitive`. Therefore, call it once here.
  75. tree.closest_point_and_primitive(Point_3(0,0,0));
  76. }
  77. template <
  78. typename Kernel,
  79. typename DerivedP,
  80. typename DerivedsqrD,
  81. typename DerivedI,
  82. typename DerivedC>
  83. IGL_INLINE void igl::copyleft::cgal::point_mesh_squared_distance(
  84. const Eigen::MatrixBase<DerivedP> & P,
  85. const CGAL::AABB_tree<
  86. CGAL::AABB_traits<Kernel,
  87. CGAL::AABB_triangle_primitive<Kernel,
  88. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  89. >
  90. >
  91. > & tree,
  92. const std::vector<CGAL::Triangle_3<Kernel> > & T,
  93. Eigen::PlainObjectBase<DerivedsqrD> & sqrD,
  94. Eigen::PlainObjectBase<DerivedI> & I,
  95. Eigen::PlainObjectBase<DerivedC> & C)
  96. {
  97. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  98. typedef typename std::vector<Triangle_3>::iterator Iterator;
  99. typedef CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  100. typedef CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  101. typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
  102. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  103. typedef CGAL::Point_3<Kernel> Point_3;
  104. assert(P.cols() == 3);
  105. const int n = P.rows();
  106. sqrD.resize(n,1);
  107. I.resize(n,1);
  108. C.resize(n,P.cols());
  109. //for(int p = 0;p < n;p++)
  110. igl::parallel_for(n,[&](const int p)
  111. {
  112. Point_3 query(P(p,0),P(p,1),P(p,2));
  113. // Find closest point and primitive id
  114. Point_and_primitive_id pp = tree.closest_point_and_primitive(query);
  115. Point_3 closest_point = pp.first;
  116. assign_scalar(closest_point[0],C(p,0));
  117. assign_scalar(closest_point[1],C(p,1));
  118. assign_scalar(closest_point[2],C(p,2));
  119. assign_scalar((closest_point-query).squared_length(),sqrD(p));
  120. I(p) = pp.second - T.begin();
  121. },1000);
  122. }
  123. #ifdef IGL_STATIC_LIBRARY
  124. // Explicit template instantiation
  125. // generated by autoexplicit.sh
  126. template void igl::copyleft::cgal::point_mesh_squared_distance<CGAL::Epeck, Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<CGAL::Epeck::FT, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>>(Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3>> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1>> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1>> const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 1, 0, -1, 1>>&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1>>&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1>>&);
  127. template void igl::copyleft::cgal::point_mesh_squared_distance<CGAL::Epeck, Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<CGAL::Epeck::FT, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  128. template void igl::copyleft::cgal::point_mesh_squared_distance<CGAL::Epeck, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  129. template void igl::copyleft::cgal::point_mesh_squared_distance<CGAL::Simple_cartesian<double>, Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 1, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> >, CGAL::Default> > const&, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> >&);
  130. template void igl::copyleft::cgal::point_mesh_squared_distance_precompute<CGAL::Simple_cartesian<double>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> >, CGAL::Default> >&, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >&);
  131. template void igl::copyleft::cgal::point_mesh_squared_distance_precompute<CGAL::Simple_cartesian<double>, Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> >, CGAL::Default> >&, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >&);
  132. #endif