slice_into.cpp 7.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157
  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_into.h"
  9. #include "IGL_ASSERT.h"
  10. #include "colon.h"
  11. #include "sort.h"
  12. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  13. #include <iostream>
  14. #include <unsupported/Eigen/SparseExtra>
  15. template <typename T, typename DerivedR, typename DerivedC>
  16. IGL_INLINE void igl::slice_into(
  17. const Eigen::SparseMatrix<T>& X,
  18. const Eigen::MatrixBase<DerivedR> & R,
  19. const Eigen::MatrixBase<DerivedC> & C,
  20. Eigen::SparseMatrix<T>& Y)
  21. {
  22. const int xm = X.rows();
  23. const int xn = X.cols();
  24. IGL_ASSERT(R.size() == xm);
  25. IGL_ASSERT(C.size() == xn);
  26. const int ym = Y.rows();
  27. const int yn = Y.cols();
  28. IGL_ASSERT(R.minCoeff() >= 0);
  29. IGL_ASSERT(R.maxCoeff() < ym);
  30. IGL_ASSERT(C.minCoeff() >= 0);
  31. IGL_ASSERT(C.maxCoeff() < yn);
  32. std::vector<bool> in_R(Y.rows());
  33. for(int r = 0;r<R.size();r++) { in_R[R(r)] = true; }
  34. // Rebuild each column in C
  35. for(int c = 0;c<C.size();c++)
  36. {
  37. int k = C(c);
  38. Eigen::SparseVector<T> Yk = Y.col(k);
  39. // implicit zeros
  40. for(typename Eigen::SparseMatrix<T>::InnerIterator yit(Y, k);yit;++yit)
  41. {
  42. if(in_R[yit.row()]){ Yk.coeffRef(yit.row()) = 0; }
  43. }
  44. // explicit values
  45. for(typename Eigen::SparseMatrix<T>::InnerIterator xit(X, c);xit;++xit)
  46. {
  47. // row in X
  48. int r = xit.row();
  49. // row in Y
  50. int s = R(r);
  51. Yk.coeffRef(s) = xit.value();
  52. }
  53. Y.col(k) = Yk;
  54. }
  55. }
  56. template <typename DerivedX, typename DerivedY, typename DerivedR, typename DerivedC>
  57. IGL_INLINE void igl::slice_into(
  58. const Eigen::MatrixBase<DerivedX> & X,
  59. const Eigen::MatrixBase<DerivedR> & R,
  60. const Eigen::MatrixBase<DerivedC> & C,
  61. Eigen::PlainObjectBase<DerivedY> & Y)
  62. {
  63. int xm = X.rows();
  64. int xn = X.cols();
  65. IGL_ASSERT(R.size() == xm);
  66. IGL_ASSERT(C.size() == xn);
  67. const int ym = Y.rows();
  68. const int yn = Y.cols();
  69. IGL_ASSERT(R.minCoeff() >= 0);
  70. IGL_ASSERT(R.maxCoeff() < ym);
  71. IGL_ASSERT(C.minCoeff() >= 0);
  72. IGL_ASSERT(C.maxCoeff() < yn);
  73. // Build reindexing maps for columns and rows, -1 means not in map
  74. Eigen::Matrix<typename DerivedR::Scalar,Eigen::Dynamic,1> RI;
  75. RI.resize(xm);
  76. for(int i = 0;i<xm;i++)
  77. {
  78. for(int j = 0;j<xn;j++)
  79. {
  80. Y(int(R(i)),int(C(j))) = X(i,j);
  81. }
  82. }
  83. }
  84. template <typename MatX, typename MatY, typename DerivedR>
  85. IGL_INLINE void igl::slice_into(
  86. const MatX & X,
  87. const Eigen::MatrixBase<DerivedR> & R,
  88. const int dim,
  89. MatY& Y)
  90. {
  91. Eigen::Matrix<typename MatX::Scalar, Eigen::Dynamic, 1> C;
  92. switch(dim)
  93. {
  94. case 1:
  95. IGL_ASSERT(R.size() == X.rows());
  96. // boring base case
  97. if(X.cols() == 0)
  98. {
  99. return;
  100. }
  101. igl::colon(0,X.cols()-1,C);
  102. return slice_into(X,R,C,Y);
  103. case 2:
  104. IGL_ASSERT(R.size() == X.cols());
  105. // boring base case
  106. if(X.rows() == 0)
  107. {
  108. return;
  109. }
  110. igl::colon(0,X.rows()-1,C);
  111. return slice_into(X,C,R,Y);
  112. default:
  113. IGL_ASSERT(false && "Unsupported dimension");
  114. return;
  115. }
  116. }
  117. template <typename DerivedX, typename DerivedR, typename DerivedY>
  118. IGL_INLINE void igl::slice_into(
  119. const Eigen::MatrixBase<DerivedX> & X,
  120. const Eigen::MatrixBase<DerivedR> & R,
  121. Eigen::PlainObjectBase<DerivedY> & Y)
  122. {
  123. // phony column indices
  124. Eigen::Matrix<typename DerivedR::Scalar, Eigen::Dynamic, 1> C;
  125. C.resize(1);
  126. C(0) = 0;
  127. return igl::slice_into(X,R,C,Y);
  128. }
  129. #ifdef IGL_STATIC_LIBRARY
  130. // Explicit template instantiation
  131. // generated by autoexplicit.sh
  132. template void igl::slice_into<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<int, -1, 1, 0, -1, 1> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int, Eigen::Matrix<int, -1, 1, 0, -1, 1>&);
  133. // generated by autoexplicit.sh
  134. template void igl::slice_into<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  135. template void igl::slice_into<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::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  136. template void igl::slice_into<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::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  137. template void igl::slice_into<double, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  138. template void igl::slice_into<Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  139. template void igl::slice_into<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  140. template void igl::slice_into<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<double, -1, -1, 0, -1, -1> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int, Eigen::Matrix<double, -1, -1, 0, -1, -1>&);
  141. template void igl::slice_into<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<double, -1, 1, 0, -1, 1> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int, Eigen::Matrix<double, -1, 1, 0, -1, 1>&);
  142. template void igl::slice_into<Eigen::SparseMatrix<double, 0, int>, Eigen::SparseMatrix<double, 0, int>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, int, Eigen::SparseMatrix<double, 0, int>&);
  143. template void igl::slice_into<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> >(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<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  144. template void igl::slice_into<Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::MatrixWrapper<Eigen::Array<int, -1, 1, 0, -1, 1> > >(Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, Eigen::MatrixBase<Eigen::MatrixWrapper<Eigen::Array<int, -1, 1, 0, -1, 1> > > const&, int, Eigen::Matrix<double, -1, 1, 0, -1, 1>&);
  145. #endif