tetrahedralize.h 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 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. #ifndef IGL_COPYLEFT_TETGEN_TETRAHEDRALIZE_H
  9. #define IGL_COPYLEFT_TETGEN_TETRAHEDRALIZE_H
  10. #include "../../igl_inline.h"
  11. #include <vector>
  12. #include <string>
  13. #include <Eigen/Core>
  14. #ifndef TETLIBRARY
  15. #define TETLIBRARY
  16. #endif
  17. #include <tetgen.h> // Defined REAL
  18. namespace igl
  19. {
  20. namespace copyleft
  21. {
  22. namespace tetgen
  23. {
  24. /// Mesh the interior of a surface mesh (V,F) using tetgen
  25. ///
  26. /// @param[in] V #V by 3 vertex position list
  27. /// @param[in] F #F list of polygon face indices into V (0-indexed)
  28. /// @param[in] H #H by 3 list of seed points inside holes
  29. /// @param[in] VM #VM list of vertex markers
  30. /// @param[in] FM #FM list of face markers
  31. /// @param[in] R #R by 5 list of region attributes
  32. /// @param[in] switches string of tetgen options (See tetgen documentation) e.g.
  33. /// "pq1.414a0.01" tries to mesh the interior of a given surface with
  34. /// quality and area constraints
  35. /// "" will mesh the convex hull constrained to pass through V (ignores F)
  36. /// @param[out] TV #TV by 3 vertex position list
  37. /// @param[out] TT #TT by 4 list of tet face indices
  38. /// @param[out] TF #TF by 3 list of triangle face indices ('f', else
  39. /// `boundary_facets` is called on TT)
  40. /// @param[out] TR #TT list of region ID for each tetrahedron
  41. /// @param[out] TN #TT by 4 list of indices neighbors for each tetrahedron ('n')
  42. /// @param[out] PT #TV list of incident tetrahedron for a vertex ('m')
  43. /// @param[out] FT #TF by 2 list of tetrahedrons sharing a triface ('nn')
  44. /// @param[out] num_regions Number of regions in output mesh
  45. /// @return status:
  46. /// 0 success
  47. /// 1 tetgen threw exception
  48. /// 2 tetgen did not crash but could not create any tets (probably there are
  49. /// holes, duplicate faces etc.)
  50. /// -1 other error
  51. ///
  52. /// \note The polygons F can contain polygons with different number of vertices.
  53. /// Trailing unused columns are filled with -1. For example, triangles and
  54. /// segments can be specified using a #F x 3 matrix: for segments the third
  55. /// column contains -1.
  56. ///
  57. /// \note Tetgen mixes integer region ids in with other region data `attr
  58. /// = (int) in->regionlist[i + 3];`. So it's declared safe to use integer
  59. /// types for `TR` since this also assumes that there's a single tet
  60. /// attribute and that it's the region id.
  61. ///
  62. /// #### Example
  63. ///
  64. /// ```cpp
  65. /// Eigen::MatrixXd V;
  66. /// Eigen::MatrixXi F;
  67. /// …
  68. /// Eigen::VectorXi VM,FM;
  69. /// Eigen::MatrixXd H,R;
  70. /// Eigen::VectorXi TM,TR,PT;
  71. /// Eigen::MatrixXi FT,TN;
  72. /// int numRegions;
  73. /// tetrahedralize(V,F,H,VM,FM,R,switches,TV,TT,TF,TM,TR,TN,PT,FT,numRegions);
  74. /// ```
  75. template <
  76. typename DerivedV,
  77. typename DerivedF,
  78. typename DerivedH,
  79. typename DerivedVM,
  80. typename DerivedFM,
  81. typename DerivedR,
  82. typename DerivedTV,
  83. typename DerivedTT,
  84. typename DerivedTF,
  85. typename DerivedTM,
  86. typename DerivedTR,
  87. typename DerivedTN,
  88. typename DerivedPT,
  89. typename DerivedFT>
  90. IGL_INLINE int tetrahedralize(
  91. const Eigen::MatrixBase<DerivedV>& V,
  92. const Eigen::MatrixBase<DerivedF>& F,
  93. const Eigen::MatrixBase<DerivedH>& H,
  94. const Eigen::MatrixBase<DerivedVM>& VM,
  95. const Eigen::MatrixBase<DerivedFM>& FM,
  96. const Eigen::MatrixBase<DerivedR>& R,
  97. const std::string switches,
  98. Eigen::PlainObjectBase<DerivedTV>& TV,
  99. Eigen::PlainObjectBase<DerivedTT>& TT,
  100. Eigen::PlainObjectBase<DerivedTF>& TF,
  101. Eigen::PlainObjectBase<DerivedTM>& TM,
  102. Eigen::PlainObjectBase<DerivedTR>& TR,
  103. Eigen::PlainObjectBase<DerivedTN>& TN,
  104. Eigen::PlainObjectBase<DerivedPT>& PT,
  105. Eigen::PlainObjectBase<DerivedFT>& FT,
  106. int & num_regions);
  107. /// \overload
  108. template <
  109. typename DerivedV,
  110. typename DerivedF,
  111. typename DerivedTV,
  112. typename DerivedTT,
  113. typename DerivedTF>
  114. IGL_INLINE int tetrahedralize(
  115. const Eigen::MatrixBase<DerivedV>& V,
  116. const Eigen::MatrixBase<DerivedF>& F,
  117. const std::string switches,
  118. Eigen::PlainObjectBase<DerivedTV>& TV,
  119. Eigen::PlainObjectBase<DerivedTT>& TT,
  120. Eigen::PlainObjectBase<DerivedTF>& TF);
  121. }
  122. }
  123. }
  124. #ifndef IGL_STATIC_LIBRARY
  125. # include "tetrahedralize.cpp"
  126. #endif
  127. #endif