point_mesh_squared_distance.cpp 8.4 KB

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