slice_tets.cpp 8.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322
  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 "edges.h"
  12. #include "slice.h"
  13. #include "cat.h"
  14. #include "ismember.h"
  15. #include "unique_rows.h"
  16. #include <cassert>
  17. #include <algorithm>
  18. #include <vector>
  19. template <
  20. typename DerivedV,
  21. typename DerivedT,
  22. typename DerivedS,
  23. typename DerivedSV,
  24. typename DerivedSF,
  25. typename DerivedJ,
  26. typename BCType>
  27. IGL_INLINE void igl::slice_tets(
  28. const Eigen::MatrixBase<DerivedV>& V,
  29. const Eigen::MatrixBase<DerivedT>& T,
  30. const Eigen::MatrixBase<DerivedS> & S,
  31. Eigen::PlainObjectBase<DerivedSV>& SV,
  32. Eigen::PlainObjectBase<DerivedSF>& SF,
  33. Eigen::PlainObjectBase<DerivedJ>& J,
  34. Eigen::SparseMatrix<BCType> & BC)
  35. {
  36. Eigen::MatrixXi sE;
  37. Eigen::Matrix<typename DerivedSV::Scalar,Eigen::Dynamic,1> lambda;
  38. igl::slice_tets(V,T,S,SV,SF,J,sE,lambda);
  39. const int ns = SV.rows();
  40. std::vector<Eigen::Triplet<BCType> > BCIJV(ns*2);
  41. for(int i = 0;i<ns;i++)
  42. {
  43. BCIJV[2*i+0] = Eigen::Triplet<BCType>(i,sE(i,0), lambda(i));
  44. BCIJV[2*i+1] = Eigen::Triplet<BCType>(i,sE(i,1),1.0-lambda(i));
  45. }
  46. BC.resize(SV.rows(),V.rows());
  47. BC.setFromTriplets(BCIJV.begin(),BCIJV.end());
  48. }
  49. template <
  50. typename DerivedV,
  51. typename DerivedT,
  52. typename DerivedS,
  53. typename DerivedSV,
  54. typename DerivedSF,
  55. typename DerivedJ,
  56. typename DerivedsE,
  57. typename Derivedlambda>
  58. IGL_INLINE void igl::slice_tets(
  59. const Eigen::MatrixBase<DerivedV>& V,
  60. const Eigen::MatrixBase<DerivedT>& T,
  61. const Eigen::MatrixBase<DerivedS> & S,
  62. Eigen::PlainObjectBase<DerivedSV>& SV,
  63. Eigen::PlainObjectBase<DerivedSF>& SF,
  64. Eigen::PlainObjectBase<DerivedJ>& J,
  65. Eigen::PlainObjectBase<DerivedsE>& sE,
  66. Eigen::PlainObjectBase<Derivedlambda>& lambda)
  67. {
  68. using namespace Eigen;
  69. using namespace std;
  70. assert(V.cols() == 3 && "V should be #V by 3");
  71. assert(T.cols() == 4 && "T should be #T by 4");
  72. static const Eigen::Matrix<int,12,4> flipped_order =
  73. (Eigen::Matrix<int,12,4>(12,4)<<
  74. 3,2,0,1,
  75. 3,1,2,0,
  76. 3,0,1,2,
  77. 2,3,1,0,
  78. 2,1,0,3,
  79. 2,0,3,1,
  80. 1,3,0,2,
  81. 1,2,3,0,
  82. 1,0,2,3,
  83. 0,3,2,1,
  84. 0,2,1,3,
  85. 0,1,3,2
  86. ).finished();
  87. // number of tets
  88. const size_t m = T.rows();
  89. typedef typename DerivedV::Scalar Scalar;
  90. typedef typename DerivedT::Scalar Index;
  91. typedef Matrix<Scalar,Dynamic,1> VectorXS;
  92. typedef Matrix<Scalar,Dynamic,4> MatrixX4S;
  93. typedef Matrix<Scalar,Dynamic,3> MatrixX3S;
  94. typedef Matrix<Scalar,Dynamic,2> MatrixX2S;
  95. typedef Matrix<Index,Dynamic,4> MatrixX4I;
  96. typedef Matrix<Index,Dynamic,3> MatrixX3I;
  97. typedef Matrix<Index,Dynamic,2> MatrixX2I;
  98. typedef Matrix<Index,Dynamic,1> VectorXI;
  99. typedef Array<bool,Dynamic,1> ArrayXb;
  100. MatrixX4S IT(m,4);
  101. for(size_t t = 0;t<m;t++)
  102. {
  103. for(size_t c = 0;c<4;c++)
  104. {
  105. IT(t,c) = S(T(t,c));
  106. }
  107. }
  108. // Essentially, just a glorified slice(X,1)
  109. //
  110. // Inputs:
  111. // T #T by 4 list of tet indices into V
  112. // IT #IT by 4 list of isosurface values at each tet
  113. // I #I list of bools whether to grab data corresponding to each tet
  114. const auto & extract_rows = [](
  115. const MatrixBase<DerivedT> & T,
  116. const MatrixX4S & IT,
  117. const ArrayXb & I,
  118. MatrixX4I & TI,
  119. MatrixX4S & ITI,
  120. VectorXI & JI)
  121. {
  122. const Index num_I = std::count(I.data(),I.data()+I.size(),true);
  123. TI.resize(num_I,4);
  124. ITI.resize(num_I,4);
  125. JI.resize(num_I,1);
  126. {
  127. size_t k = 0;
  128. for(size_t t = 0;t<(size_t)T.rows();t++)
  129. {
  130. if(I(t))
  131. {
  132. TI.row(k) = T.row(t);
  133. ITI.row(k) = IT.row(t);
  134. JI(k) = t;
  135. k++;
  136. }
  137. }
  138. assert(k == num_I);
  139. }
  140. };
  141. ArrayXb I13 = (IT.array()<0).rowwise().count()==1;
  142. ArrayXb I31 = (IT.array()>0).rowwise().count()==1;
  143. ArrayXb I22 = (IT.array()<0).rowwise().count()==2;
  144. MatrixX4I T13,T31,T22;
  145. MatrixX4S IT13,IT31,IT22;
  146. VectorXI J13,J31,J22;
  147. extract_rows(T,IT,I13,T13,IT13,J13);
  148. extract_rows(T,IT,I31,T31,IT31,J31);
  149. extract_rows(T,IT,I22,T22,IT22,J22);
  150. const auto & apply_sort4 = [] (
  151. const MatrixX4I & T,
  152. const MatrixX4I & sJ,
  153. MatrixX4I & sT)
  154. {
  155. sT.resize(T.rows(),4);
  156. for(size_t t = 0;t<(size_t)T.rows();t++)
  157. {
  158. for(size_t c = 0;c<4;c++)
  159. {
  160. sT(t,c) = T(t,sJ(t,c));
  161. }
  162. }
  163. };
  164. const auto & apply_sort2 = [] (
  165. const MatrixX2I & E,
  166. const MatrixX2I & sJ,
  167. Eigen::PlainObjectBase<DerivedsE>& sE)
  168. {
  169. sE.resize(E.rows(),2);
  170. for(size_t t = 0;t<(size_t)E.rows();t++)
  171. {
  172. for(size_t c = 0;c<2;c++)
  173. {
  174. sE(t,c) = E(t,sJ(t,c));
  175. }
  176. }
  177. };
  178. const auto & one_below = [&apply_sort4](
  179. const MatrixX4I & T,
  180. const MatrixX4S & IT,
  181. MatrixX2I & U,
  182. MatrixX3I & SF)
  183. {
  184. // Number of tets
  185. const size_t m = T.rows();
  186. MatrixX4S sIT;
  187. MatrixX4I sJ;
  188. sort(IT,2,true,sIT,sJ);
  189. MatrixX4I sT;
  190. apply_sort4(T,sJ,sT);
  191. U.resize(3*m,2);
  192. U<<
  193. sT.col(0),sT.col(1),
  194. sT.col(0),sT.col(2),
  195. sT.col(0),sT.col(3);
  196. SF.resize(m,3);
  197. for(size_t c = 0;c<3;c++)
  198. {
  199. SF.col(c) =
  200. igl::LinSpaced<
  201. Eigen::Matrix<typename DerivedSF::Scalar,Eigen::Dynamic,1> >
  202. (m,0+c*m,(m-1)+c*m);
  203. }
  204. ArrayXb flip;
  205. {
  206. VectorXi _;
  207. ismember_rows(sJ,flipped_order,flip,_);
  208. }
  209. for(int i = 0;i<m;i++)
  210. {
  211. if(flip(i))
  212. {
  213. SF.row(i) = SF.row(i).reverse().eval();
  214. }
  215. }
  216. };
  217. const auto & two_below = [&apply_sort4](
  218. const MatrixX4I & T,
  219. const MatrixX4S & IT,
  220. MatrixX2I & U,
  221. MatrixX3I & SF)
  222. {
  223. // Number of tets
  224. const size_t m = T.rows();
  225. MatrixX4S sIT;
  226. MatrixX4I sJ;
  227. sort(IT,2,true,sIT,sJ);
  228. MatrixX4I sT;
  229. apply_sort4(T,sJ,sT);
  230. U.resize(4*m,2);
  231. U<<
  232. sT.col(0),sT.col(2),
  233. sT.col(0),sT.col(3),
  234. sT.col(1),sT.col(2),
  235. sT.col(1),sT.col(3);
  236. SF.resize(2*m,3);
  237. SF.block(0,0,m,1) = igl::LinSpaced<VectorXI >(m,0+0*m,(m-1)+0*m);
  238. SF.block(0,1,m,1) = igl::LinSpaced<VectorXI >(m,0+1*m,(m-1)+1*m);
  239. SF.block(0,2,m,1) = igl::LinSpaced<VectorXI >(m,0+3*m,(m-1)+3*m);
  240. SF.block(m,0,m,1) = igl::LinSpaced<VectorXI >(m,0+0*m,(m-1)+0*m);
  241. SF.block(m,1,m,1) = igl::LinSpaced<VectorXI >(m,0+3*m,(m-1)+3*m);
  242. SF.block(m,2,m,1) = igl::LinSpaced<VectorXI >(m,0+2*m,(m-1)+2*m);
  243. ArrayXb flip;
  244. {
  245. VectorXi _;
  246. ismember_rows(sJ,flipped_order,flip,_);
  247. }
  248. for(int i = 0;i<m;i++)
  249. {
  250. if(flip(i))
  251. {
  252. SF.row(i ) = SF.row(i ).reverse().eval();
  253. SF.row(i+m) = SF.row(i+m).reverse().eval();
  254. }
  255. }
  256. };
  257. MatrixX3I SF13,SF31,SF22;
  258. MatrixX2I U13,U31,U22;
  259. one_below(T13, IT13,U13,SF13);
  260. one_below(T31,-IT31,U31,SF31);
  261. two_below(T22, IT22,U22,SF22);
  262. const MatrixX2I U =
  263. (MatrixX2I(U13.rows()+ U31.rows()+ U22.rows(),2)<<U13,U31,U22).finished();
  264. MatrixX2I sU;
  265. {
  266. MatrixX2I _;
  267. sort(U,2,true,sU,_);
  268. }
  269. MatrixX2I E;
  270. VectorXI uI,uJ;
  271. unique_rows(sU,E,uI,uJ);
  272. MatrixX2S IE(E.rows(),2);
  273. for(size_t t = 0;t<E.rows();t++)
  274. {
  275. for(size_t c = 0;c<2;c++)
  276. {
  277. IE(t,c) = S(E(t,c));
  278. }
  279. }
  280. MatrixX2S sIE;
  281. MatrixX2I sJ;
  282. sort(IE,2,true,sIE,sJ);
  283. apply_sort2(E,sJ,sE);
  284. lambda = sIE.col(1).array() / (sIE.col(1)-sIE.col(0)).array();
  285. SV.resize(sE.rows(),V.cols());
  286. for(int e = 0;e<sE.rows();e++)
  287. {
  288. SV.row(e) = V.row(sE(e,0))*lambda(e) + V.row(sE(e,1))*(1.0-lambda(e));
  289. }
  290. SF.resize( SF13.rows()+SF31.rows()+SF22.rows(),3);
  291. SF<<
  292. SF13,
  293. U13.rows()+ SF31.rowwise().reverse().array(),
  294. U13.rows()+U31.rows()+SF22.array();
  295. std::for_each(
  296. SF.data(),
  297. SF.data()+SF.size(),
  298. [&uJ](typename DerivedSF::Scalar & i){i=uJ(i);});
  299. J.resize(SF.rows());
  300. J<<J13,J31,J22,J22;
  301. }
  302. #ifdef IGL_STATIC_LIBRARY
  303. // Explicit template instantiation
  304. 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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<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>&);
  305. #endif