signed_distance.cpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <alecjacobson@gmail.com>
  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 "signed_distance.h"
  9. template <typename Kernel>
  10. IGL_INLINE typename Kernel::FT igl::signed_distance_pseudonormal(
  11. const CGAL::AABB_tree<
  12. CGAL::AABB_traits<Kernel,
  13. CGAL::AABB_triangle_primitive<Kernel,
  14. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  15. >
  16. >
  17. > & tree,
  18. const std::vector<CGAL::Triangle_3<Kernel> > & T,
  19. const Eigen::MatrixXi & F,
  20. const Eigen::MatrixXd & FN,
  21. const Eigen::MatrixXd & VN,
  22. const Eigen::MatrixXd & EN,
  23. const Eigen::VectorXi & EMAP,
  24. const typename Kernel::Point_3 & q)
  25. {
  26. using namespace Eigen;
  27. using namespace std;
  28. typedef typename Kernel::FT FT;
  29. typedef typename Kernel::Point_3 Point_3;
  30. typedef typename CGAL::Triangle_3<Kernel> Triangle_3;
  31. typedef typename std::vector<Triangle_3>::iterator Iterator;
  32. typedef typename CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  33. typedef typename CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  34. typedef typename CGAL::AABB_tree<AABB_triangle_traits> Tree;
  35. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  36. Point_and_primitive_id pp = tree.closest_point_and_primitive(q);
  37. Point_3 & p = pp.first;
  38. const auto & qp = q-p;
  39. const FT sqrd = qp.squared_length();
  40. Vector3d v(qp.x(),qp.y(),qp.z());
  41. const int f = pp.second - T.begin();
  42. const Triangle_3 & t = *pp.second;
  43. // barycentric coordinates
  44. const auto & area = [&p,&t](const int i, const int j)->FT
  45. {
  46. return sqrt(Triangle_3(p,t.vertex(i),t.vertex(j)).squared_area());
  47. };
  48. Vector3d b(area(1,2),area(2,0),area(0,1));
  49. b /= b.sum();
  50. // Determine which normal to use
  51. Vector3d n;
  52. const double epsilon = 1e-12;
  53. const int type = (b.array()<=epsilon).cast<int>().sum();
  54. switch(type)
  55. {
  56. case 2:
  57. // Find vertex
  58. for(int c = 0;c<3;c++)
  59. {
  60. if(b(c)>epsilon)
  61. {
  62. n = VN.row(F(f,c));
  63. break;
  64. }
  65. }
  66. break;
  67. case 1:
  68. // Find edge
  69. for(int c = 0;c<3;c++)
  70. {
  71. if(b(c)<=epsilon)
  72. {
  73. n = EN.row(EMAP(F.rows()*c+f));
  74. break;
  75. }
  76. }
  77. break;
  78. default:
  79. assert(false && "all barycentric coords zero.");
  80. case 0:
  81. n = FN.row(f);
  82. break;
  83. }
  84. const double s = (v.dot(n) >= 0 ? 1. : -1.);
  85. return s*sqrt(sqrd);
  86. }
  87. template <typename Kernel>
  88. IGL_INLINE typename Kernel::FT igl::signed_distance_winding_number(
  89. const CGAL::AABB_tree<
  90. CGAL::AABB_traits<Kernel,
  91. CGAL::AABB_triangle_primitive<Kernel,
  92. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  93. >
  94. >
  95. > & tree,
  96. const igl::WindingNumberAABB<Eigen::Vector3d> & hier,
  97. const typename Kernel::Point_3 & q)
  98. {
  99. typedef typename Kernel::FT FT;
  100. typedef typename Kernel::Point_3 Point_3;
  101. typedef typename CGAL::Triangle_3<Kernel> Triangle_3;
  102. typedef typename std::vector<Triangle_3>::iterator Iterator;
  103. typedef typename CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  104. typedef typename CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  105. typedef typename CGAL::AABB_tree<AABB_triangle_traits> Tree;
  106. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  107. using namespace Eigen;
  108. using namespace std;
  109. using namespace igl;
  110. Point_and_primitive_id pp = tree.closest_point_and_primitive(q);
  111. Point_3 & p = pp.first;
  112. const auto & qp = q-p;
  113. const Vector3d eq(q.x(),q.y(),q.z());
  114. const FT sqrd = qp.squared_length();
  115. const double w = hier.winding_number(eq);
  116. const FT s = 1.-2.*w;
  117. return s*sqrt(sqrd);
  118. }
  119. #ifdef IGL_STATIC_LIBRARY
  120. template CGAL::Epick::FT igl::signed_distance_winding_number<CGAL::Epick>(CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Epick, CGAL::AABB_triangle_primitive<CGAL::Epick, std::vector<CGAL::Triangle_3<CGAL::Epick>, std::allocator<CGAL::Triangle_3<CGAL::Epick> > >::iterator, CGAL::Boolean_tag<false> > > > const&, igl::WindingNumberAABB<Eigen::Matrix<double, 3, 1, 0, 3, 1> > const&, CGAL::Epick::Point_3 const&);
  121. template CGAL::Epick::FT igl::signed_distance_pseudonormal<CGAL::Epick>(CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Epick, CGAL::AABB_triangle_primitive<CGAL::Epick, std::vector<CGAL::Triangle_3<CGAL::Epick>, std::allocator<CGAL::Triangle_3<CGAL::Epick> > >::iterator, CGAL::Boolean_tag<false> > > > const&, std::vector<CGAL::Triangle_3<CGAL::Epick>, std::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, CGAL::Epick::Point_3 const&);
  122. #endif