cat.cpp 7.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 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 "cat.h"
  9. #include <cstdio>
  10. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  11. #include <iostream>
  12. #include <unsupported/Eigen/SparseExtra>
  13. // Sparse matrices need to be handled carefully. Because C++ does not
  14. // Template:
  15. // Scalar sparse matrix scalar type, e.g. double
  16. template <typename Scalar>
  17. IGL_INLINE void igl::cat(
  18. const int dim,
  19. const Eigen::SparseMatrix<Scalar> & A,
  20. const Eigen::SparseMatrix<Scalar> & B,
  21. Eigen::SparseMatrix<Scalar> & C)
  22. {
  23. assert(dim == 1 || dim == 2);
  24. using namespace Eigen;
  25. // Special case if B or A is empty
  26. if(A.size() == 0)
  27. {
  28. C = B;
  29. return;
  30. }
  31. if(B.size() == 0)
  32. {
  33. C = A;
  34. return;
  35. }
  36. #if false
  37. // This **must** be DynamicSparseMatrix, otherwise this implementation is
  38. // insanely slow
  39. DynamicSparseMatrix<Scalar, RowMajor> dyn_C;
  40. if(dim == 1)
  41. {
  42. assert(A.cols() == B.cols());
  43. dyn_C.resize(A.rows()+B.rows(),A.cols());
  44. }else if(dim == 2)
  45. {
  46. assert(A.rows() == B.rows());
  47. dyn_C.resize(A.rows(),A.cols()+B.cols());
  48. }else
  49. {
  50. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  51. }
  52. dyn_C.reserve(A.nonZeros()+B.nonZeros());
  53. // Iterate over outside of A
  54. for(int k=0; k<A.outerSize(); ++k)
  55. {
  56. // Iterate over inside
  57. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  58. {
  59. dyn_C.coeffRef(it.row(),it.col()) += it.value();
  60. }
  61. }
  62. // Iterate over outside of B
  63. for(int k=0; k<B.outerSize(); ++k)
  64. {
  65. // Iterate over inside
  66. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  67. {
  68. int r = (dim == 1 ? A.rows()+it.row() : it.row());
  69. int c = (dim == 2 ? A.cols()+it.col() : it.col());
  70. dyn_C.coeffRef(r,c) += it.value();
  71. }
  72. }
  73. C = SparseMatrix<Scalar>(dyn_C);
  74. #elif false
  75. std::vector<Triplet<Scalar> > CIJV;
  76. CIJV.reserve(A.nonZeros() + B.nonZeros());
  77. {
  78. // Iterate over outside of A
  79. for(int k=0; k<A.outerSize(); ++k)
  80. {
  81. // Iterate over inside
  82. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  83. {
  84. CIJV.emplace_back(it.row(),it.col(),it.value());
  85. }
  86. }
  87. // Iterate over outside of B
  88. for(int k=0; k<B.outerSize(); ++k)
  89. {
  90. // Iterate over inside
  91. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  92. {
  93. int r = (dim == 1 ? A.rows()+it.row() : it.row());
  94. int c = (dim == 2 ? A.cols()+it.col() : it.col());
  95. CIJV.emplace_back(r,c,it.value());
  96. }
  97. }
  98. }
  99. C = SparseMatrix<Scalar>(
  100. dim == 1 ? A.rows()+B.rows() : A.rows(),
  101. dim == 1 ? A.cols() : A.cols()+B.cols());
  102. C.reserve(A.nonZeros() + B.nonZeros());
  103. C.setFromTriplets(CIJV.begin(),CIJV.end());
  104. #else
  105. C = SparseMatrix<Scalar>(
  106. dim == 1 ? A.rows()+B.rows() : A.rows(),
  107. dim == 1 ? A.cols() : A.cols()+B.cols());
  108. Eigen::VectorXi per_col = Eigen::VectorXi::Zero(C.cols());
  109. if(dim == 1)
  110. {
  111. assert(A.outerSize() == B.outerSize());
  112. for(int k = 0;k<A.outerSize();++k)
  113. {
  114. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  115. {
  116. per_col(k)++;
  117. }
  118. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  119. {
  120. per_col(k)++;
  121. }
  122. }
  123. }else
  124. {
  125. for(int k = 0;k<A.outerSize();++k)
  126. {
  127. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  128. {
  129. per_col(k)++;
  130. }
  131. }
  132. for(int k = 0;k<B.outerSize();++k)
  133. {
  134. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  135. {
  136. per_col(A.cols() + k)++;
  137. }
  138. }
  139. }
  140. C.reserve(per_col);
  141. if(dim == 1)
  142. {
  143. for(int k = 0;k<A.outerSize();++k)
  144. {
  145. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  146. {
  147. C.insert(it.row(),k) = it.value();
  148. }
  149. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  150. {
  151. C.insert(A.rows()+it.row(),k) = it.value();
  152. }
  153. }
  154. }else
  155. {
  156. for(int k = 0;k<A.outerSize();++k)
  157. {
  158. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  159. {
  160. C.insert(it.row(),k) = it.value();
  161. }
  162. }
  163. for(int k = 0;k<B.outerSize();++k)
  164. {
  165. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  166. {
  167. C.insert(it.row(),A.cols()+k) = it.value();
  168. }
  169. }
  170. }
  171. C.makeCompressed();
  172. #endif
  173. }
  174. template <typename Derived, class MatC>
  175. IGL_INLINE void igl::cat(
  176. const int dim,
  177. const Eigen::MatrixBase<Derived> & A,
  178. const Eigen::MatrixBase<Derived> & B,
  179. MatC & C)
  180. {
  181. assert(dim == 1 || dim == 2);
  182. // Special case if B or A is empty
  183. if(A.size() == 0)
  184. {
  185. C = B;
  186. return;
  187. }
  188. if(B.size() == 0)
  189. {
  190. C = A;
  191. return;
  192. }
  193. if(dim == 1)
  194. {
  195. assert(A.cols() == B.cols());
  196. C.resize(A.rows()+B.rows(),A.cols());
  197. C << A,B;
  198. }else if(dim == 2)
  199. {
  200. assert(A.rows() == B.rows());
  201. C.resize(A.rows(),A.cols()+B.cols());
  202. C << A,B;
  203. }else
  204. {
  205. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  206. }
  207. }
  208. template <class Mat>
  209. IGL_INLINE Mat igl::cat(const int dim, const Mat & A, const Mat & B)
  210. {
  211. assert(dim == 1 || dim == 2);
  212. Mat C;
  213. igl::cat(dim,A,B,C);
  214. return C;
  215. }
  216. template <class Mat>
  217. IGL_INLINE void igl::cat(const std::vector<std::vector< Mat > > & A, Mat & C)
  218. {
  219. using namespace std;
  220. // Start with empty matrix
  221. C.resize(0,0);
  222. for(const auto & row_vec : A)
  223. {
  224. // Concatenate each row horizontally
  225. // Start with empty matrix
  226. Mat row(0,0);
  227. for(const auto & element : row_vec)
  228. {
  229. row = cat(2,row,element);
  230. }
  231. // Concatenate rows vertically
  232. C = cat(1,C,row);
  233. }
  234. }
  235. #ifdef IGL_STATIC_LIBRARY
  236. // Explicit template instantiation
  237. // generated by autoexplicit.sh
  238. template Eigen::Matrix<double, -1, -1, 0, -1, -1> igl::cat<Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&);
  239. // generated by autoexplicit.sh
  240. template Eigen::SparseMatrix<double, 0, int> igl::cat<Eigen::SparseMatrix<double, 0, int> >(int, Eigen::SparseMatrix<double, 0, int> const&, Eigen::SparseMatrix<double, 0, int> const&);
  241. // generated by autoexplicit.sh
  242. template Eigen::Matrix<int, -1, -1, 0, -1, -1> igl::cat<Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&);
  243. template void igl::cat<Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(int, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::Matrix<double, -1, 1, 0, -1, 1>&);
  244. template Eigen::Matrix<int, -1, 1, 0, -1, 1> igl::cat<Eigen::Matrix<int, -1, 1, 0, -1, 1> >(int, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&);
  245. template Eigen::Matrix<double, -1, 1, 0, -1, 1> igl::cat<Eigen::Matrix<double, -1, 1, 0, -1, 1> >(int, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&);
  246. template void igl::cat<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::Matrix<double, -1, -1, 0, -1, -1>&);
  247. template void igl::cat<Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::Matrix<int, -1, -1, 0, -1, -1>&);
  248. #endif