tetrahedralize.cpp 6.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 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 "tetrahedralize.h"
  9. #include "mesh_to_tetgenio.h"
  10. #include "tetgenio_to_tetmesh.h"
  11. // IGL includes
  12. #include "../../matrix_to_list.h"
  13. #include "../../list_to_matrix.h"
  14. #include "../../boundary_facets.h"
  15. // STL includes
  16. #include <cassert>
  17. #include <iostream>
  18. IGL_INLINE int igl::copyleft::tetgen::tetrahedralize(
  19. const std::vector<std::vector<REAL > > & V,
  20. const std::vector<std::vector<int> > & F,
  21. const std::string switches,
  22. std::vector<std::vector<REAL > > & TV,
  23. std::vector<std::vector<int > > & TT,
  24. std::vector<std::vector<int> > & TF)
  25. {
  26. using namespace std;
  27. tetgenio in,out;
  28. bool success;
  29. success = mesh_to_tetgenio(V,F,in);
  30. if(!success)
  31. {
  32. return -1;
  33. }
  34. try
  35. {
  36. char * cswitches = new char[switches.size() + 1];
  37. std::strcpy(cswitches,switches.c_str());
  38. ::tetrahedralize(cswitches,&in, &out);
  39. delete[] cswitches;
  40. }catch(int e)
  41. {
  42. cerr<<"^"<<__FUNCTION__<<": TETGEN CRASHED... KABOOOM!!!"<<endl;
  43. return 1;
  44. }
  45. if(out.numberoftetrahedra == 0)
  46. {
  47. cerr<<"^"<<__FUNCTION__<<": Tetgen failed to create tets"<<endl;
  48. return 2;
  49. }
  50. success = tetgenio_to_tetmesh(out,TV,TT,TF);
  51. if(!success)
  52. {
  53. return -1;
  54. }
  55. //boundary_facets(TT,TF);
  56. return 0;
  57. }
  58. template <
  59. typename DerivedV,
  60. typename DerivedF,
  61. typename DerivedTV,
  62. typename DerivedTT,
  63. typename DerivedTF>
  64. IGL_INLINE int igl::copyleft::tetgen::tetrahedralize(
  65. const Eigen::PlainObjectBase<DerivedV>& V,
  66. const Eigen::PlainObjectBase<DerivedF>& F,
  67. const std::string switches,
  68. Eigen::PlainObjectBase<DerivedTV>& TV,
  69. Eigen::PlainObjectBase<DerivedTT>& TT,
  70. Eigen::PlainObjectBase<DerivedTF>& TF)
  71. {
  72. using namespace std;
  73. vector<vector<REAL> > vV,vTV;
  74. vector<vector<int> > vF,vTT,vTF;
  75. matrix_to_list(V,vV);
  76. matrix_to_list(F,vF);
  77. int e = tetrahedralize(vV,vF,switches,vTV,vTT,vTF);
  78. if(e == 0)
  79. {
  80. bool TV_rect = list_to_matrix(vTV,TV);
  81. if(!TV_rect)
  82. {
  83. return 3;
  84. }
  85. bool TT_rect = list_to_matrix(vTT,TT);
  86. if(!TT_rect)
  87. {
  88. return 3;
  89. }
  90. bool TF_rect = list_to_matrix(vTF,TF);
  91. if(!TF_rect)
  92. {
  93. return 3;
  94. }
  95. }
  96. return e;
  97. }
  98. template <
  99. typename DerivedV,
  100. typename DerivedF,
  101. typename DerivedVM,
  102. typename DerivedFM,
  103. typename DerivedTV,
  104. typename DerivedTT,
  105. typename DerivedTF,
  106. typename DerivedTM>
  107. IGL_INLINE int igl::copyleft::tetgen::tetrahedralize(
  108. const Eigen::PlainObjectBase<DerivedV>& V,
  109. const Eigen::PlainObjectBase<DerivedF>& F,
  110. const Eigen::PlainObjectBase<DerivedVM>& VM,
  111. const Eigen::PlainObjectBase<DerivedFM>& FM,
  112. const std::string switches,
  113. Eigen::PlainObjectBase<DerivedTV>& TV,
  114. Eigen::PlainObjectBase<DerivedTT>& TT,
  115. Eigen::PlainObjectBase<DerivedTF>& TF,
  116. Eigen::PlainObjectBase<DerivedTM>& TM)
  117. {
  118. using namespace std;
  119. vector<vector<REAL> > vV,vTV;
  120. vector<vector<int> > vF,vTT,vTF;
  121. vector<int> vTM;
  122. matrix_to_list(V,vV);
  123. matrix_to_list(F,vF);
  124. vector<int> vVM = matrix_to_list(VM);
  125. vector<int> vFM = matrix_to_list(FM);
  126. int e = tetrahedralize(vV,vF,vVM,vFM,switches,vTV,vTT,vTF,vTM);
  127. if(e == 0)
  128. {
  129. bool TV_rect = list_to_matrix(vTV,TV);
  130. if(!TV_rect)
  131. {
  132. return false;
  133. }
  134. bool TT_rect = list_to_matrix(vTT,TT);
  135. if(!TT_rect)
  136. {
  137. return false;
  138. }
  139. bool TF_rect = list_to_matrix(vTF,TF);
  140. if(!TF_rect)
  141. {
  142. return false;
  143. }
  144. bool TM_rect = list_to_matrix(vTM,TM);
  145. if(!TM_rect)
  146. {
  147. return false;
  148. }
  149. }
  150. return e;
  151. }
  152. IGL_INLINE int igl::copyleft::tetgen::tetrahedralize(
  153. const std::vector<std::vector<REAL > > & V,
  154. const std::vector<std::vector<int> > & F,
  155. const std::vector<int> & VM,
  156. const std::vector<int> & FM,
  157. const std::string switches,
  158. std::vector<std::vector<REAL > > & TV,
  159. std::vector<std::vector<int > > & TT,
  160. std::vector<std::vector<int> > & TF,
  161. std::vector<int> & TM)
  162. {
  163. using namespace std;
  164. tetgenio in,out;
  165. bool success;
  166. success = mesh_to_tetgenio(V,F,in);
  167. if(!success)
  168. {
  169. return -1;
  170. }
  171. in.pointmarkerlist = new int[VM.size()];
  172. for (int i = 0; i < VM.size(); ++i) {
  173. in.pointmarkerlist[i] = VM[i];
  174. }
  175. // These have already been created in mesh_to_tetgenio.
  176. // Reset them here.
  177. for (int i = 0; i < FM.size(); ++i) {
  178. in.facetmarkerlist[i] = FM[i];
  179. }
  180. try
  181. {
  182. char * cswitches = new char[switches.size() + 1];
  183. std::strcpy(cswitches,switches.c_str());
  184. ::tetrahedralize(cswitches,&in, &out);
  185. delete[] cswitches;
  186. }catch(int e)
  187. {
  188. cerr<<"^"<<__FUNCTION__<<": TETGEN CRASHED... KABOOOM!!!"<<endl;
  189. return 1;
  190. }
  191. if(out.numberoftetrahedra == 0)
  192. {
  193. cerr<<"^"<<__FUNCTION__<<": Tetgen failed to create tets"<<endl;
  194. return 2;
  195. }
  196. success = tetgenio_to_tetmesh(out,TV,TT,TF);
  197. if(!success)
  198. {
  199. return -1;
  200. }
  201. TM.resize(out.numberofpoints);
  202. for (int i = 0; i < out.numberofpoints; ++i) {
  203. TM[i] = out.pointmarkerlist[i];
  204. }
  205. //boundary_facets(TT,TF);
  206. return 0;
  207. }
  208. #ifdef IGL_STATIC_LIBRARY
  209. // Explicit template specialization
  210. template int igl::copyleft::tetgen::tetrahedralize<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<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&, std::basic_string<char, std::char_traits<char>, std::allocator<char> >, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  211. template int igl::copyleft::tetgen::tetrahedralize<Eigen::Matrix<double, -1, -1, 0, -1, -1>,Eigen::Matrix<int, -1, -1, 0, -1, -1>,Eigen::Matrix<int, -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<int, -1, -1, 0, -1, -1>,Eigen::Matrix<int, -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<int, -1, 1, 0, -1, 1> > &,const Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > &,const std::basic_string<char, std::char_traits<char>, std::allocator<char> >,Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > &,Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &,Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > &);
  212. #endif