slice.cpp 8.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229
  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. #include "slice.h"
  9. #include "colon.h"
  10. #include <vector>
  11. template <
  12. typename TX,
  13. typename TY,
  14. typename DerivedR,
  15. typename DerivedC>
  16. IGL_INLINE void igl::slice(
  17. const Eigen::SparseMatrix<TX> &X,
  18. const Eigen::DenseBase<DerivedR> &R,
  19. const Eigen::DenseBase<DerivedC> &C,
  20. Eigen::SparseMatrix<TY> &Y)
  21. {
  22. int xm = X.rows();
  23. int xn = X.cols();
  24. int ym = R.size();
  25. int yn = C.size();
  26. // special case when R or C is empty
  27. if (ym == 0 || yn == 0)
  28. {
  29. Y.resize(ym, yn);
  30. return;
  31. }
  32. assert(R.minCoeff() >= 0);
  33. assert(R.maxCoeff() < xm);
  34. assert(C.minCoeff() >= 0);
  35. assert(C.maxCoeff() < xn);
  36. // Build reindexing maps for columns and rows
  37. std::vector<std::vector<typename DerivedR::Scalar>> RI;
  38. RI.resize(xm);
  39. for (int i = 0; i < ym; i++)
  40. {
  41. RI[R(i)].push_back(i);
  42. }
  43. std::vector<std::vector<typename DerivedC::Scalar>> CI;
  44. CI.resize(xn);
  45. for (int i = 0; i < yn; i++)
  46. {
  47. CI[C(i)].push_back(i);
  48. }
  49. // Take a guess at the number of nonzeros (this assumes uniform distribution
  50. // not banded or heavily diagonal)
  51. std::vector<Eigen::Triplet<TY>> entries;
  52. entries.reserve((X.nonZeros()/(X.rows()*X.cols())) * (ym*yn));
  53. // Iterate over outside
  54. for (int k = 0; k < X.outerSize(); ++k)
  55. {
  56. // Iterate over inside
  57. for (typename Eigen::SparseMatrix<TX>::InnerIterator it(X, k); it; ++it)
  58. {
  59. for (auto rit = RI[it.row()].begin(); rit != RI[it.row()].end(); rit++)
  60. {
  61. for (auto cit = CI[it.col()].begin(); cit != CI[it.col()].end(); cit++)
  62. {
  63. entries.emplace_back(*rit, *cit, it.value());
  64. }
  65. }
  66. }
  67. }
  68. Y.resize(ym, yn);
  69. Y.setFromTriplets(entries.begin(), entries.end());
  70. }
  71. template <typename MatX, typename DerivedR, typename MatY>
  72. IGL_INLINE void igl::slice(
  73. const MatX &X,
  74. const Eigen::DenseBase<DerivedR> &R,
  75. const int dim,
  76. MatY &Y)
  77. {
  78. Eigen::Matrix<typename DerivedR::Scalar, Eigen::Dynamic, 1> C;
  79. switch (dim)
  80. {
  81. case 1:
  82. // boring base case
  83. if (X.cols() == 0)
  84. {
  85. Y.resize(R.size(), 0);
  86. return;
  87. }
  88. igl::colon(0, X.cols() - 1, C);
  89. return slice(X, R, C, Y);
  90. case 2:
  91. // boring base case
  92. if (X.rows() == 0)
  93. {
  94. Y.resize(0, R.size());
  95. return;
  96. }
  97. igl::colon(0, X.rows() - 1, C);
  98. return slice(X, C, R, Y);
  99. default:
  100. assert(false && "Unsupported dimension");
  101. return;
  102. }
  103. }
  104. template <
  105. typename DerivedX,
  106. typename DerivedR,
  107. typename DerivedC,
  108. typename DerivedY>
  109. IGL_INLINE void igl::slice(
  110. const Eigen::DenseBase<DerivedX> &X,
  111. const Eigen::DenseBase<DerivedR> &R,
  112. const Eigen::DenseBase<DerivedC> &C,
  113. Eigen::PlainObjectBase<DerivedY> &Y)
  114. {
  115. #ifndef NDEBUG
  116. int xm = X.rows();
  117. int xn = X.cols();
  118. #endif
  119. int ym = R.size();
  120. int yn = C.size();
  121. // special case when R or C is empty
  122. if (ym == 0 || yn == 0)
  123. {
  124. Y.resize(ym, yn);
  125. return;
  126. }
  127. assert(R.minCoeff() >= 0);
  128. assert(R.maxCoeff() < xm);
  129. assert(C.minCoeff() >= 0);
  130. assert(C.maxCoeff() < xn);
  131. // Resize output
  132. Y.resize(ym, yn);
  133. // loop over output rows, then columns
  134. for (int i = 0; i < ym; i++)
  135. {
  136. for (int j = 0; j < yn; j++)
  137. {
  138. Y(i, j) = X(R(i), C(j));
  139. }
  140. }
  141. }
  142. template <typename DerivedX, typename DerivedY, typename DerivedR>
  143. IGL_INLINE void igl::slice(
  144. const Eigen::DenseBase<DerivedX> &X,
  145. const Eigen::DenseBase<DerivedR> &R,
  146. Eigen::PlainObjectBase<DerivedY> &Y)
  147. {
  148. // phony column indices
  149. Eigen::Matrix<typename DerivedR::Scalar, Eigen::Dynamic, 1> C;
  150. C.resize(1);
  151. C(0) = 0;
  152. return igl::slice(X, R, C, Y);
  153. }
  154. template <typename DerivedX, typename DerivedR>
  155. IGL_INLINE DerivedX igl::slice(
  156. const Eigen::DenseBase<DerivedX> &X,
  157. const Eigen::DenseBase<DerivedR> &R)
  158. {
  159. // This is not safe. See PlainMatrix
  160. DerivedX Y;
  161. igl::slice(X, R, Y);
  162. return Y;
  163. }
  164. template <typename DerivedX, typename DerivedR>
  165. IGL_INLINE DerivedX igl::slice(
  166. const Eigen::DenseBase<DerivedX> &X,
  167. const Eigen::DenseBase<DerivedR> &R,
  168. const int dim)
  169. {
  170. // This is not safe. See PlainMatrix
  171. DerivedX Y;
  172. igl::slice(X, R, dim, Y);
  173. return Y;
  174. }
  175. template< class T >
  176. IGL_INLINE void igl::slice(
  177. const std::vector<T> & unordered,
  178. std::vector<size_t> const & index_map,
  179. std::vector<T> & ordered)
  180. {
  181. // copy for the slice according to index_map, because unordered may also be
  182. // ordered
  183. std::vector<T> copy = unordered;
  184. ordered.resize(index_map.size());
  185. for(int i = 0; i<(int)index_map.size();i++)
  186. {
  187. ordered[i] = copy[index_map[i]];
  188. }
  189. }
  190. #ifdef IGL_STATIC_LIBRARY
  191. // Explicit template instantiation
  192. template void igl::slice<std::int64_t>(std::vector<std::int64_t, std::allocator<std::int64_t>> const&, std::vector<size_t, std::allocator<size_t>> const&, std::vector<std::int64_t, std::allocator<std::int64_t>>&);
  193. template void igl::slice<Eigen::SparseMatrix<bool, 0, int>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::SparseMatrix<bool, 0, int>>(Eigen::SparseMatrix<bool, 0, int> const &, Eigen::DenseBase<Eigen::Matrix<int, -1, 1, 0, -1, 1>> const &, int, Eigen::SparseMatrix<bool, 0, int> &);
  194. template void igl::slice<Eigen::SparseMatrix<double, 0, int>, Eigen::Array<int, -1, 1, 0, -1, 1>, Eigen::SparseMatrix<double, 0, int> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::DenseBase<Eigen::Array<int, -1, 1, 0, -1, 1> > const&, int, Eigen::SparseMatrix<double, 0, int>&);
  195. template void igl::slice<Eigen::SparseMatrix<double, 0, int>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::SparseMatrix<double, 0, int>>(Eigen::SparseMatrix<double, 0, int> const &, Eigen::DenseBase<Eigen::Matrix<double, -1, 1, 0, -1, 1>> const &, int, Eigen::SparseMatrix<double, 0, int> &);
  196. template void igl::slice<Eigen::SparseMatrix<double, 0, int>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::SparseMatrix<double, 0, int>>(Eigen::SparseMatrix<double, 0, int> const &, Eigen::DenseBase<Eigen::Matrix<int, -1, -1, 0, -1, -1>> const &, int, Eigen::SparseMatrix<double, 0, int> &);
  197. template void igl::slice<Eigen::SparseMatrix<double, 0, int>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::SparseMatrix<double, 0, int>>(Eigen::SparseMatrix<double, 0, int> const &, Eigen::DenseBase<Eigen::Matrix<int, -1, 1, 0, -1, 1>> const &, int, Eigen::SparseMatrix<double, 0, int> &);
  198. #ifdef WIN32
  199. template void igl::slice<unsigned int>(class std::vector<unsigned int,class std::allocator<unsigned int> > const &,class std::vector<unsigned __int64,class std::allocator<unsigned __int64> > const &,class std::vector<unsigned int,class std::allocator<unsigned int> > &);
  200. template void igl::slice<float>(class std::vector<float,class std::allocator<float> > const &,class std::vector<unsigned __int64,class std::allocator<unsigned __int64> > const &,class std::vector<float,class std::allocator<float> > &);
  201. template void igl::slice<__int64>(class std::vector<__int64,class std::allocator<__int64> > const &,class std::vector<unsigned __int64,class std::allocator<unsigned __int64> > const &,class std::vector<__int64,class std::allocator<__int64> > &);
  202. #endif
  203. template void igl::slice<unsigned int>(std::vector<unsigned int, std::allocator<unsigned int> > const&, std::vector<size_t, std::allocator<size_t> > const&, std::vector<unsigned int, std::allocator<unsigned int> >&);
  204. template void igl::slice<float>(std::vector<float, std::allocator<float> > const&, std::vector<size_t, std::allocator<size_t> > const&, std::vector<float, std::allocator<float> >&);
  205. template void igl::slice<double>(std::vector<double, std::allocator<double> > const&, std::vector<size_t, std::allocator<size_t> > const&, std::vector<double, std::allocator<double> >&);
  206. template void igl::slice<int>(std::vector<int, std::allocator<int> > const&, std::vector<size_t, std::allocator<size_t> > const&, std::vector<int, std::allocator<int> >&);
  207. #include "SortableRow.h"
  208. template void igl::slice<igl::SortableRow<Eigen::Matrix<int, -1, 1, 0, -1, 1> > >(std::vector<igl::SortableRow<Eigen::Matrix<int, -1, 1, 0, -1, 1> >, std::allocator<igl::SortableRow<Eigen::Matrix<int, -1, 1, 0, -1, 1> > > > const&, std::vector<size_t, std::allocator<size_t> > const&, std::vector<igl::SortableRow<Eigen::Matrix<int, -1, 1, 0, -1, 1> >, std::allocator<igl::SortableRow<Eigen::Matrix<int, -1, 1, 0, -1, 1> > > >&);
  209. template void igl::slice<igl::SortableRow<Eigen::Matrix<double, -1, 1, 0, -1, 1> > >(std::vector<igl::SortableRow<Eigen::Matrix<double, -1, 1, 0, -1, 1> >, std::allocator<igl::SortableRow<Eigen::Matrix<double, -1, 1, 0, -1, 1> > > > const&, std::vector<size_t, std::allocator<size_t> > const&, std::vector<igl::SortableRow<Eigen::Matrix<double, -1, 1, 0, -1, 1> >, std::allocator<igl::SortableRow<Eigen::Matrix<double, -1, 1, 0, -1, 1> > > >&);
  210. #endif