signed_distance_isosurface.cpp 5.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149
  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_isosurface.h"
  9. #include "point_mesh_squared_distance.h"
  10. #include "complex_to_mesh.h"
  11. #include "signed_distance.h"
  12. #include <igl/per_face_normals.h>
  13. #include <igl/per_edge_normals.h>
  14. #include <igl/per_vertex_normals.h>
  15. #include <igl/centroid.h>
  16. #include <igl/WindingNumberAABB.h>
  17. #include <igl/matlab_format.h>
  18. #include <igl/remove_unreferenced.h>
  19. #include <CGAL/Surface_mesh_default_triangulation_3.h>
  20. #include <CGAL/Complex_2_in_triangulation_3.h>
  21. #include <CGAL/make_surface_mesh.h>
  22. #include <CGAL/Implicit_surface_3.h>
  23. #include <CGAL/Polyhedron_3.h>
  24. #include <CGAL/IO/output_surface_facets_to_polyhedron.h>
  25. // Axis-aligned bounding box tree for tet tri intersection
  26. #include <CGAL/AABB_tree.h>
  27. #include <CGAL/AABB_traits.h>
  28. #include <CGAL/AABB_triangle_primitive.h>
  29. #include <vector>
  30. IGL_INLINE bool igl::signed_distance_isosurface(
  31. const Eigen::MatrixXd & IV,
  32. const Eigen::MatrixXi & IF,
  33. const double level,
  34. const double angle_bound,
  35. const double radius_bound,
  36. const double distance_bound,
  37. const SignedDistanceType sign_type,
  38. Eigen::MatrixXd & V,
  39. Eigen::MatrixXi & F)
  40. {
  41. using namespace std;
  42. // default triangulation for Surface_mesher
  43. typedef CGAL::Surface_mesh_default_triangulation_3 Tr;
  44. // c2t3
  45. typedef CGAL::Complex_2_in_triangulation_3<Tr> C2t3;
  46. typedef Tr::Geom_traits GT;//Kernel
  47. typedef GT::Sphere_3 Sphere_3;
  48. typedef GT::Point_3 Point_3;
  49. typedef GT::FT FT;
  50. typedef std::function<FT (Point_3)> Function;
  51. typedef CGAL::Implicit_surface_3<GT, Function> Surface_3;
  52. typedef CGAL::Polyhedron_3<GT> Polyhedron;
  53. typedef GT::Kernel Kernel;
  54. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  55. typedef typename std::vector<Triangle_3>::iterator Iterator;
  56. typedef CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  57. typedef CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  58. typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
  59. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  60. Tree tree;
  61. vector<Triangle_3 > T;
  62. point_mesh_squared_distance_precompute(IV,IF,tree,T);
  63. Eigen::MatrixXd FN,VN,EN;
  64. Eigen::MatrixXi E;
  65. Eigen::VectorXi EMAP;
  66. WindingNumberAABB<Eigen::Vector3d> hier;
  67. switch(sign_type)
  68. {
  69. default:
  70. assert(false && "Unknown SignedDistanceType");
  71. case SIGNED_DISTANCE_TYPE_DEFAULT:
  72. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  73. hier.set_mesh(IV,IF);
  74. hier.grow();
  75. break;
  76. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  77. // "Signed Distance Computation Using the Angle Weighted Pseudonormal"
  78. // [Bærentzen & Aanæs 2005]
  79. per_face_normals(IV,IF,FN);
  80. per_vertex_normals(IV,IF,PER_VERTEX_NORMALS_WEIGHTING_TYPE_ANGLE,FN,VN);
  81. per_edge_normals(
  82. IV,IF,PER_EDGE_NORMALS_WEIGHTING_TYPE_UNIFORM,FN,EN,E,EMAP);
  83. break;
  84. }
  85. Tr tr; // 3D-Delaunay triangulation
  86. C2t3 c2t3 (tr); // 2D-complex in 3D-Delaunay triangulation
  87. // defining the surface
  88. const auto & IVmax = IV.colwise().maxCoeff();
  89. const auto & IVmin = IV.colwise().minCoeff();
  90. const double bbd = (IVmax-IVmin).norm();
  91. const double r = bbd/2.;
  92. const auto & IVmid = 0.5*(IVmax + IVmin);
  93. // Supposedly the implict needs to evaluate to <0 at cmid...
  94. // http://doc.cgal.org/latest/Surface_mesher/classCGAL_1_1Implicit__surface__3.html
  95. Point_3 cmid(IVmid(0),IVmid(1),IVmid(2));
  96. Function fun;
  97. switch(sign_type)
  98. {
  99. default:
  100. assert(false && "Unknown SignedDistanceType");
  101. case SIGNED_DISTANCE_TYPE_DEFAULT:
  102. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  103. fun =
  104. [&tree,&hier,&level](const Point_3 & q) -> FT
  105. {
  106. return signed_distance_winding_number(tree,hier,q)-level;
  107. };
  108. break;
  109. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  110. fun = [&tree,&T,&IF,&FN,&VN,&EN,&EMAP,&level](const Point_3 & q) -> FT
  111. {
  112. return
  113. igl::signed_distance_pseudonormal(tree,T,IF,FN,VN,EN,EMAP,q) -
  114. level;
  115. };
  116. break;
  117. }
  118. //[&tree,&hier,&T,&IF,&FN,&VN,&EN,&EMAP,&level](const Point_3 & q) ->FT
  119. //{
  120. // const FT p = signed_distance_pseudonormal(tree,T,IF,FN,VN,EN,EMAP,q);
  121. // const FT w = signed_distance_winding_number(tree,hier,q);
  122. // if(w*p < 0 && (fabs(w) > 0.1 || fabs(p) > 0.1))
  123. // {
  124. // cout<<"q=["<<q.x()<<","<<q.y()<<","<<q.z()<<"];"<<endl;
  125. // cout<<matlab_format(n.transpose().eval(),"n")<<endl;
  126. // cout<<matlab_format(b.transpose().eval(),"b")<<endl;
  127. // cout<<"Sig difference: "<<type<<endl;
  128. // cout<<"w: "<<w<<endl;
  129. // cout<<"p: "<<p<<endl;
  130. // exit(1);
  131. // }
  132. // return w;
  133. //},
  134. Sphere_3 bounding_sphere(cmid, (r+level)*(r+level));
  135. Surface_3 surface(fun,bounding_sphere);
  136. CGAL::Surface_mesh_default_criteria_3<Tr>
  137. criteria(angle_bound,radius_bound,distance_bound);
  138. // meshing surface
  139. CGAL::make_surface_mesh(c2t3, surface, criteria, CGAL::Manifold_tag());
  140. // complex to (V,F)
  141. return igl::complex_to_mesh(c2t3,V,F);
  142. }