cat.cpp 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263
  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 "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. // Special case if B or A is empty
  25. if(A.size() == 0)
  26. {
  27. C = B;
  28. return;
  29. }
  30. if(B.size() == 0)
  31. {
  32. C = A;
  33. return;
  34. }
  35. // This is faster than using DynamicSparseMatrix or setFromTriplets
  36. C = Eigen::SparseMatrix<Scalar>(
  37. dim == 1 ? A.rows()+B.rows() : A.rows(),
  38. dim == 1 ? A.cols() : A.cols()+B.cols());
  39. Eigen::VectorXi per_col = Eigen::VectorXi::Zero(C.cols());
  40. if(dim == 1)
  41. {
  42. assert(A.outerSize() == B.outerSize());
  43. for(int k = 0;k<A.outerSize();++k)
  44. {
  45. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  46. {
  47. per_col(k)++;
  48. }
  49. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  50. {
  51. per_col(k)++;
  52. }
  53. }
  54. }else
  55. {
  56. for(int k = 0;k<A.outerSize();++k)
  57. {
  58. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  59. {
  60. per_col(k)++;
  61. }
  62. }
  63. for(int k = 0;k<B.outerSize();++k)
  64. {
  65. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  66. {
  67. per_col(A.cols() + k)++;
  68. }
  69. }
  70. }
  71. C.reserve(per_col);
  72. if(dim == 1)
  73. {
  74. for(int k = 0;k<A.outerSize();++k)
  75. {
  76. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  77. {
  78. C.insert(it.row(),k) = it.value();
  79. }
  80. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  81. {
  82. C.insert(A.rows()+it.row(),k) = it.value();
  83. }
  84. }
  85. }else
  86. {
  87. for(int k = 0;k<A.outerSize();++k)
  88. {
  89. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  90. {
  91. C.insert(it.row(),k) = it.value();
  92. }
  93. }
  94. for(int k = 0;k<B.outerSize();++k)
  95. {
  96. for(typename Eigen::SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  97. {
  98. C.insert(it.row(),A.cols()+k) = it.value();
  99. }
  100. }
  101. }
  102. C.makeCompressed();
  103. }
  104. template <typename Derived, class MatC>
  105. IGL_INLINE void igl::cat(
  106. const int dim,
  107. const Eigen::MatrixBase<Derived> & A,
  108. const Eigen::MatrixBase<Derived> & B,
  109. MatC & C)
  110. {
  111. assert(dim == 1 || dim == 2);
  112. // Special case if B or A is empty
  113. if(A.size() == 0)
  114. {
  115. C = B;
  116. return;
  117. }
  118. if(B.size() == 0)
  119. {
  120. C = A;
  121. return;
  122. }
  123. if(dim == 1)
  124. {
  125. assert(A.cols() == B.cols());
  126. C.resize(A.rows()+B.rows(),A.cols());
  127. C << A,B;
  128. }else if(dim == 2)
  129. {
  130. assert(A.rows() == B.rows());
  131. C.resize(A.rows(),A.cols()+B.cols());
  132. C << A,B;
  133. }else
  134. {
  135. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  136. }
  137. }
  138. template <class Mat>
  139. IGL_INLINE Mat igl::cat(const int dim, const Mat & A, const Mat & B)
  140. {
  141. assert(dim == 1 || dim == 2);
  142. Mat C;
  143. igl::cat(dim,A,B,C);
  144. return C;
  145. }
  146. template <class Mat>
  147. IGL_INLINE void igl::cat(const std::vector<std::vector< Mat > > & A, Mat & C)
  148. {
  149. // Start with empty matrix
  150. C.resize(0,0);
  151. for(const auto & row_vec : A)
  152. {
  153. // Concatenate each row horizontally
  154. // Start with empty matrix
  155. Mat row(0,0);
  156. for(const auto & element : row_vec)
  157. {
  158. row = cat(2,row,element);
  159. }
  160. // Concatenate rows vertically
  161. C = cat(1,C,row);
  162. }
  163. }
  164. template <typename T, typename DerivedC>
  165. IGL_INLINE void igl::cat(const int dim, const std::vector<T> & A, Eigen::PlainObjectBase<DerivedC> & C)
  166. {
  167. assert(dim == 1 || dim == 2);
  168. const int num_mat = A.size();
  169. if(num_mat == 0)
  170. {
  171. C.resize(0,0);
  172. return;
  173. }
  174. if(dim == 1)
  175. {
  176. const int A_cols = A[0].cols();
  177. int tot_rows = 0;
  178. for(const auto & m : A)
  179. {
  180. tot_rows += m.rows();
  181. }
  182. C.resize(tot_rows, A_cols);
  183. int cur_row = 0;
  184. for(int i = 0; i < num_mat; i++)
  185. {
  186. assert(A_cols == A[i].cols());
  187. C.block(cur_row,0,A[i].rows(),A_cols) = A[i];
  188. cur_row += A[i].rows();
  189. }
  190. }
  191. else if(dim == 2)
  192. {
  193. const int A_rows = A[0].rows();
  194. int tot_cols = 0;
  195. for(const auto & m : A)
  196. {
  197. tot_cols += m.cols();
  198. }
  199. C.resize(A_rows,tot_cols);
  200. int cur_col = 0;
  201. for(int i = 0; i < num_mat; i++)
  202. {
  203. assert(A_rows == A[i].rows());
  204. C.block(0,cur_col,A_rows,A[i].cols()) = A[i];
  205. cur_col += A[i].cols();
  206. }
  207. }
  208. else
  209. {
  210. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  211. }
  212. }
  213. #ifdef IGL_STATIC_LIBRARY
  214. // Explicit template instantiation
  215. // generated by autoexplicit.sh
  216. template void igl::cat<Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(int, std::vector<Eigen::Matrix<int, -1, 1, 0, -1, 1>, std::allocator<Eigen::Matrix<int, -1, 1, 0, -1, 1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  217. // generated by autoexplicit.sh
  218. template void igl::cat<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<double, -1, -1, 0, -1, -1>, std::allocator<Eigen::Matrix<double, -1, -1, 0, -1, -1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  219. // generated by autoexplicit.sh
  220. template void igl::cat<Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<float, -1, -1, 0, -1, -1>, std::allocator<Eigen::Matrix<float, -1, -1, 0, -1, -1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&);
  221. // generated by autoexplicit.sh
  222. 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&);
  223. // generated by autoexplicit.sh
  224. 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&);
  225. // generated by autoexplicit.sh
  226. 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&);
  227. 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>&);
  228. 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&);
  229. 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&);
  230. 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>&);
  231. 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>&);
  232. 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>&);
  233. template void igl::cat<Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>, std::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  234. template void igl::cat<Eigen::Matrix<int, 1, 4, 1, 1, 4>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 1, 4, 1, 1, 4>, std::allocator<Eigen::Matrix<int, 1, 4, 1, 1, 4> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  235. template void igl::cat<Eigen::Matrix<int, 1, 15, 1, 1, 15>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 1, 15, 1, 1, 15>, std::allocator<Eigen::Matrix<int, 1, 15, 1, 1, 15> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  236. template void igl::cat<Eigen::Matrix<int, 1, 2, 1, 1, 2>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 1, 2, 1, 1, 2>, std::allocator<Eigen::Matrix<int, 1, 2, 1, 1, 2> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  237. template void igl::cat<Eigen::Matrix<int, 1, 27, 1, 1, 27>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 1, 27, 1, 1, 27>, std::allocator<Eigen::Matrix<int, 1, 27, 1, 1, 27> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  238. template void igl::cat<Eigen::Matrix<int, 1, 3, 1, 1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 1, 3, 1, 1, 3>, std::allocator<Eigen::Matrix<int, 1, 3, 1, 1, 3> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  239. template void igl::cat<Eigen::Matrix<int, 3, 1, 0, 3, 1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 3, 1, 0, 3, 1>, std::allocator<Eigen::Matrix<int, 3, 1, 0, 3, 1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  240. template void igl::cat<Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<double, 1, 3, 1, 1, 3>, std::allocator<Eigen::Matrix<double, 1, 3, 1, 1, 3> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  241. template void igl::cat<Eigen::Matrix<double, 3, 1, 0, 3, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<double, 3, 1, 0, 3, 1>, std::allocator<Eigen::Matrix<double, 3, 1, 0, 3, 1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  242. template void igl::cat<Eigen::Matrix<int, 1, -1, 1, 1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<int, 1, -1, 1, 1, -1>, std::allocator<Eigen::Matrix<int, 1, -1, 1, 1, -1> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  243. template void igl::cat<Eigen::Matrix<double, 1, 2, 1, 1, 2>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, std::vector<Eigen::Matrix<double, 1, 2, 1, 1, 2>, std::allocator<Eigen::Matrix<double, 1, 2, 1, 1, 2> > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  244. #endif