slice.cpp 5.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160
  1. #include "slice.h"
  2. #include "colon.h"
  3. #include <vector>
  4. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  5. #include <iostream>
  6. #include <unsupported/Eigen/SparseExtra>
  7. template <typename T>
  8. IGL_INLINE void igl::slice(
  9. const Eigen::SparseMatrix<T>& X,
  10. const Eigen::Matrix<int,Eigen::Dynamic,1> & R,
  11. const Eigen::Matrix<int,Eigen::Dynamic,1> & C,
  12. Eigen::SparseMatrix<T>& Y)
  13. {
  14. int xm = X.rows();
  15. int xn = X.cols();
  16. int ym = R.size();
  17. int yn = C.size();
  18. // special case when R or C is empty
  19. if(ym == 0 || yn == 0)
  20. {
  21. Y.resize(ym,yn);
  22. return;
  23. }
  24. assert(R.minCoeff() >= 0);
  25. assert(R.maxCoeff() < xm);
  26. assert(C.minCoeff() >= 0);
  27. assert(C.maxCoeff() < xn);
  28. // Build reindexing maps for columns and rows, -1 means not in map
  29. std::vector<std::vector<int> > RI;
  30. RI.resize(xm);
  31. for(int i = 0;i<ym;i++)
  32. {
  33. RI[R(i)].push_back(i);
  34. }
  35. std::vector<std::vector<int> > CI;
  36. CI.resize(xn);
  37. // initialize to -1
  38. for(int i = 0;i<yn;i++)
  39. {
  40. CI[C(i)].push_back(i);
  41. }
  42. // Resize output
  43. Eigen::DynamicSparseMatrix<T, Eigen::RowMajor> dyn_Y(ym,yn);
  44. // Take a guess at the number of nonzeros (this assumes uniform distribution
  45. // not banded or heavily diagonal)
  46. dyn_Y.reserve((X.nonZeros()/(X.rows()*X.cols())) * (ym*yn));
  47. // Iterate over outside
  48. for(int k=0; k<X.outerSize(); ++k)
  49. {
  50. // Iterate over inside
  51. for(typename Eigen::SparseMatrix<T>::InnerIterator it (X,k); it; ++it)
  52. {
  53. std::vector<int>::iterator rit, cit;
  54. for(rit = RI[it.row()].begin();rit != RI[it.row()].end(); rit++)
  55. {
  56. for(cit = CI[it.col()].begin();cit != CI[it.col()].end(); cit++)
  57. {
  58. dyn_Y.coeffRef(*rit,*cit) = it.value();
  59. }
  60. }
  61. }
  62. }
  63. Y = Eigen::SparseMatrix<T>(dyn_Y);
  64. }
  65. template <typename T>
  66. IGL_INLINE void igl::slice(
  67. const Eigen::SparseMatrix<T>& X,
  68. const Eigen::Matrix<int,Eigen::Dynamic,1> & R,
  69. const int dim,
  70. Eigen::SparseMatrix<T>& Y)
  71. {
  72. Eigen::VectorXi C;
  73. switch(dim)
  74. {
  75. case 1:
  76. igl::colon(0,X.cols()-1,C);
  77. return slice(X,R,C,Y);
  78. case 2:
  79. igl::colon(0,X.rows()-1,C);
  80. return slice(X,C,R,Y);
  81. default:
  82. assert(false && "Unsupported dimension");
  83. return;
  84. }
  85. }
  86. template <typename DerivedX>
  87. IGL_INLINE void igl::slice(
  88. const Eigen::PlainObjectBase<DerivedX> & X,
  89. const Eigen::Matrix<int,Eigen::Dynamic,1> & R,
  90. const Eigen::Matrix<int,Eigen::Dynamic,1> & C,
  91. Eigen::PlainObjectBase<DerivedX> & Y)
  92. {
  93. #ifndef NDEBUG
  94. int xm = X.rows();
  95. int xn = X.cols();
  96. #endif
  97. int ym = R.size();
  98. int yn = C.size();
  99. // special case when R or C is empty
  100. if(ym == 0 || yn == 0)
  101. {
  102. Y.resize(ym,yn);
  103. return;
  104. }
  105. assert(R.minCoeff() >= 0);
  106. assert(R.maxCoeff() < xm);
  107. assert(C.minCoeff() >= 0);
  108. assert(C.maxCoeff() < xn);
  109. // Resize output
  110. Y.resize(ym,yn);
  111. // loop over output rows, then columns
  112. for(int i = 0;i<ym;i++)
  113. {
  114. for(int j = 0;j<yn;j++)
  115. {
  116. Y(i,j) = X(R(i),C(j));
  117. }
  118. }
  119. }
  120. template <typename DerivedX>
  121. IGL_INLINE void igl::slice(
  122. const Eigen::PlainObjectBase<DerivedX> & X,
  123. const Eigen::Matrix<int,Eigen::Dynamic,1> & R,
  124. Eigen::PlainObjectBase<DerivedX> & Y)
  125. {
  126. // phony column indices
  127. Eigen::Matrix<int,Eigen::Dynamic,1> C;
  128. C.resize(1);
  129. C(0) = 0;
  130. return igl::slice(X,R,C,Y);
  131. }
  132. #ifndef IGL_HEADER_ONLY
  133. // Explicit template specialization
  134. // generated by autoexplicit.sh
  135. template void igl::slice<Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  136. // generated by autoexplicit.sh
  137. template void igl::slice<Eigen::Matrix<float, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> > const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> >&);
  138. // generated by autoexplicit.sh
  139. template void igl::slice<Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  140. // generated by autoexplicit.sh
  141. template void igl::slice<Eigen::Matrix<float, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&);
  142. // generated by autoexplicit.sh
  143. template void igl::slice<double>(Eigen::SparseMatrix<double, 0, int> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::SparseMatrix<double, 0, int>&);
  144. template void igl::slice<Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  145. template void igl::slice<double>(Eigen::SparseMatrix<double, 0, int> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, int, Eigen::SparseMatrix<double, 0, int>&);
  146. #endif