slice_tets.cpp 7.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 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 "slice_tets.h"
  9. #include "LinSpaced.h"
  10. #include "sort.h"
  11. #include "cat.h"
  12. #include "per_face_normals.h"
  13. #include <cassert>
  14. #include <algorithm>
  15. #include <vector>
  16. template <
  17. typename DerivedV,
  18. typename DerivedT,
  19. typename Derivedplane,
  20. typename DerivedU,
  21. typename DerivedG,
  22. typename DerivedJ,
  23. typename BCType>
  24. IGL_INLINE void igl::slice_tets(
  25. const Eigen::PlainObjectBase<DerivedV>& V,
  26. const Eigen::PlainObjectBase<DerivedT>& T,
  27. const Eigen::PlainObjectBase<Derivedplane> & plane,
  28. Eigen::PlainObjectBase<DerivedU>& U,
  29. Eigen::PlainObjectBase<DerivedG>& G,
  30. Eigen::PlainObjectBase<DerivedJ>& J,
  31. Eigen::SparseMatrix<BCType> & BC)
  32. {
  33. using namespace Eigen;
  34. using namespace std;
  35. assert(V.cols() == 3 && "V should be #V by 3");
  36. assert(T.cols() == 4 && "T should be #T by 4");
  37. assert(plane.size() == 4 && "Plane equation should be 4 coefficients");
  38. // number of tets
  39. const size_t m = T.rows();
  40. typedef typename DerivedV::Scalar Scalar;
  41. typedef typename DerivedT::Scalar Index;
  42. typedef Matrix<Scalar,Dynamic,1> VectorXS;
  43. typedef Matrix<Scalar,Dynamic,4> MatrixX4S;
  44. typedef Matrix<Scalar,Dynamic,3> MatrixX3S;
  45. typedef Matrix<Scalar,Dynamic,2> MatrixX2S;
  46. typedef Matrix<Index,Dynamic,4> MatrixX4I;
  47. typedef Matrix<Index,Dynamic,3> MatrixX3I;
  48. typedef Matrix<Index,Dynamic,1> VectorXI;
  49. typedef Matrix<bool,Dynamic,1> VectorXb;
  50. // Value of plane's implicit function at all vertices
  51. VectorXS IV =
  52. (V.col(0)*plane(0) +
  53. V.col(1)*plane(1) +
  54. V.col(2)*plane(2)).array()
  55. + plane(3);
  56. MatrixX4S IT(m,4);
  57. for(size_t t = 0;t<m;t++)
  58. {
  59. for(size_t c = 0;c<4;c++)
  60. {
  61. IT(t,c) = IV(T(t,c));
  62. }
  63. }
  64. const auto & extract_rows = [](
  65. const PlainObjectBase<DerivedT> & T,
  66. const MatrixX4S & IT,
  67. const VectorXb & I,
  68. MatrixX4I & TI,
  69. MatrixX4S & ITI,
  70. VectorXI & JI)
  71. {
  72. const Index num_I = std::count(I.data(),I.data()+I.size(),true);
  73. TI.resize(num_I,4);
  74. ITI.resize(num_I,4);
  75. JI.resize(num_I,1);
  76. {
  77. size_t k = 0;
  78. for(size_t t = 0;t<(size_t)T.rows();t++)
  79. {
  80. if(I(t))
  81. {
  82. TI.row(k) = T.row(t);
  83. ITI.row(k) = IT.row(t);
  84. JI(k) = t;
  85. k++;
  86. }
  87. }
  88. assert(k == num_I);
  89. }
  90. };
  91. VectorXb I13 = (IT.array()<0).rowwise().count()==1;
  92. VectorXb I31 = (IT.array()>0).rowwise().count()==1;
  93. VectorXb I22 = (IT.array()<0).rowwise().count()==2;
  94. MatrixX4I T13,T31,T22;
  95. MatrixX4S IT13,IT31,IT22;
  96. VectorXI J13,J31,J22;
  97. extract_rows(T,IT,I13,T13,IT13,J13);
  98. extract_rows(T,IT,I31,T31,IT31,J31);
  99. extract_rows(T,IT,I22,T22,IT22,J22);
  100. const auto & apply_sort = [] (
  101. const MatrixX4I & T,
  102. const MatrixX4I & sJ,
  103. MatrixX4I & sT)
  104. {
  105. sT.resize(T.rows(),4);
  106. for(size_t t = 0;t<(size_t)T.rows();t++)
  107. {
  108. for(size_t c = 0;c<4;c++)
  109. {
  110. sT(t,c) = T(t,sJ(t,c));
  111. }
  112. }
  113. };
  114. const auto & one_below = [&V,&apply_sort](
  115. const MatrixX4I & T,
  116. const MatrixX4S & IT,
  117. MatrixX3I & G,
  118. SparseMatrix<BCType> & BC)
  119. {
  120. // Number of tets
  121. const size_t m = T.rows();
  122. MatrixX4S sIT;
  123. MatrixX4I sJ;
  124. sort(IT,2,true,sIT,sJ);
  125. MatrixX4I sT;
  126. apply_sort(T,sJ,sT);
  127. MatrixX3S lambda =
  128. sIT.rightCols(3).array() /
  129. (sIT.rightCols(3).colwise()-sIT.col(0)).array();
  130. vector<Triplet<BCType> > IJV;
  131. IJV.reserve(m*3*2);
  132. for(size_t c = 0;c<3;c++)
  133. {
  134. for(size_t t = 0;t<(size_t)m;t++)
  135. {
  136. IJV.push_back(Triplet<BCType>(c*m+t, sT(t,0), lambda(t,c)));
  137. IJV.push_back(Triplet<BCType>(c*m+t,sT(t,c+1),1-lambda(t,c)));
  138. }
  139. }
  140. BC.resize(m*3,V.rows());
  141. BC.reserve(m*3*2);
  142. BC.setFromTriplets(IJV.begin(),IJV.end());
  143. G.resize(m,3);
  144. for(size_t c = 0;c<3;c++)
  145. {
  146. G.col(c) =
  147. igl::LinSpaced<
  148. Eigen::Matrix<typename DerivedG::Scalar,Eigen::Dynamic,1> >
  149. (m,0+c*m,(m-1)+c*m);
  150. }
  151. };
  152. const auto & two_below = [&V,&apply_sort](
  153. const MatrixX4I & T,
  154. const MatrixX4S & IT,
  155. MatrixX3I & G,
  156. SparseMatrix<BCType> & BC)
  157. {
  158. // Number of tets
  159. const size_t m = T.rows();
  160. MatrixX4S sIT;
  161. MatrixX4I sJ;
  162. sort(IT,2,true,sIT,sJ);
  163. MatrixX4I sT;
  164. apply_sort(T,sJ,sT);
  165. MatrixX2S lambda =
  166. sIT.rightCols(2).array() /
  167. (sIT.rightCols(2).colwise()-sIT.col(0)).array();
  168. MatrixX2S gamma =
  169. sIT.rightCols(2).array() /
  170. (sIT.rightCols(2).colwise()-sIT.col(1)).array();
  171. vector<Triplet<BCType> > IJV;
  172. IJV.reserve(m*4*2);
  173. for(size_t c = 0;c<2;c++)
  174. {
  175. for(size_t t = 0;t<(size_t)m;t++)
  176. {
  177. IJV.push_back(Triplet<BCType>(0*2*m+c*m+t, sT(t,0), lambda(t,c)));
  178. IJV.push_back(Triplet<BCType>(0*2*m+c*m+t,sT(t,c+2),1-lambda(t,c)));
  179. IJV.push_back(Triplet<BCType>(1*2*m+c*m+t, sT(t,1), gamma(t,c)));
  180. IJV.push_back(Triplet<BCType>(1*2*m+c*m+t,sT(t,c+2),1- gamma(t,c)));
  181. }
  182. }
  183. BC.resize(m*4,V.rows());
  184. BC.reserve(m*4*2);
  185. BC.setFromTriplets(IJV.begin(),IJV.end());
  186. G.resize(2*m,3);
  187. G.block(0,0,m,1) = igl::LinSpaced<VectorXI >(m,0+0*m,(m-1)+0*m);
  188. G.block(0,1,m,1) = igl::LinSpaced<VectorXI >(m,0+1*m,(m-1)+1*m);
  189. G.block(0,2,m,1) = igl::LinSpaced<VectorXI >(m,0+3*m,(m-1)+3*m);
  190. G.block(m,0,m,1) = igl::LinSpaced<VectorXI >(m,0+0*m,(m-1)+0*m);
  191. G.block(m,1,m,1) = igl::LinSpaced<VectorXI >(m,0+3*m,(m-1)+3*m);
  192. G.block(m,2,m,1) = igl::LinSpaced<VectorXI >(m,0+2*m,(m-1)+2*m);
  193. };
  194. MatrixX3I G13,G31,G22;
  195. SparseMatrix<BCType> BC13,BC31,BC22;
  196. one_below(T13,IT13,G13,BC13);
  197. one_below(T31,-IT31,G31,BC31);
  198. two_below(T22,IT22,G22,BC22);
  199. BC = cat(1,cat(1,BC13,BC31),BC22);
  200. U = BC*V;
  201. G.resize(G13.rows()+G31.rows()+G22.rows(),3);
  202. G<<G13,(G31.array()+BC13.rows()),(G22.array()+BC13.rows()+BC31.rows());
  203. MatrixX3S N;
  204. per_face_normals(U,G,N);
  205. Matrix<Scalar,1,3> planeN(plane(0),plane(1),plane(2));
  206. VectorXb flip = (N.array().rowwise() * planeN.array()).rowwise().sum()<0;
  207. for(size_t g = 0;g<(size_t)G.rows();g++)
  208. {
  209. if(flip(g))
  210. {
  211. G.row(g) = G.row(g).reverse().eval();
  212. }
  213. }
  214. J.resize(G.rows());
  215. J<<J13,J31,J22,J22;
  216. }
  217. #ifdef IGL_STATIC_LIBRARY
  218. // Explicit template instantiation
  219. template void igl::slice_tets<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, 4, 1, 0, 4, 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>, double>(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, 4, 1, 0, 4, 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<int, -1, 1, 0, -1, 1> >&, Eigen::SparseMatrix<double, 0, int>&);
  220. template void igl::slice_tets<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<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, double>(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> > 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<int, -1, 1, 0, -1, 1> >&, Eigen::SparseMatrix<double, 0, int>&);
  221. #endif