cat.h 4.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192
  1. #ifndef IGL_CAT_H
  2. #define IGL_CAT_H
  3. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  4. #include <Eigen/Sparse>
  5. #include <Eigen/Dense>
  6. namespace igl
  7. {
  8. // If you're using Dense matrices you might be better off using the << operator
  9. // This is an attempt to act like matlab's cat function.
  10. // Perform concatenation of a two matrices along a single dimension
  11. // If dim == 1, then C = [A;B]. If dim == 2 then C = [A B]
  12. //
  13. // Template:
  14. // Scalar scalar data type for sparse matrices like double or int
  15. // Mat matrix type for all matrices (e.g. MatrixXd, SparseMatrix)
  16. // MatC matrix type for ouput matrix (e.g. MatrixXd) needs to support
  17. // resize
  18. // Inputs:
  19. // A first input matrix
  20. // B second input matrix
  21. // dim dimension along which to concatenate, 0 or 1
  22. // Outputs:
  23. // C output matrix
  24. //
  25. template <typename Scalar>
  26. inline void cat(
  27. const int dim,
  28. const Eigen::SparseMatrix<Scalar> & A,
  29. const Eigen::SparseMatrix<Scalar> & B,
  30. Eigen::SparseMatrix<Scalar> & C);
  31. template <typename Derived, class MatC>
  32. inline void cat(
  33. const int dim,
  34. const Eigen::MatrixBase<Derived> & A,
  35. const Eigen::MatrixBase<Derived> & B,
  36. MatC & C);
  37. // Wrapper that returns C
  38. template <class Mat>
  39. inline Mat cat(const int dim, const Mat & A, const Mat & B);
  40. // Note: Maybe we can autogenerate a bunch of overloads D = cat(int,A,B,C),
  41. // E = cat(int,A,B,C,D), etc.
  42. // Concatenate a "matrix" of blocks
  43. // C = [A0;A1;A2;...;An] where Ai = [A[i][0] A[i][1] ... A[i][m]];
  44. //
  45. // Inputs:
  46. // A a matrix (vector of row vectors)
  47. // Output:
  48. // C
  49. template <class Mat>
  50. inline void cat(const std::vector<std::vector< Mat > > & A, Mat & C);
  51. }
  52. // Implementation
  53. // Sparse matrices need to be handled carefully. Because C++ does not
  54. // Template:
  55. // Scalar sparse matrix scalar type, e.g. double
  56. template <typename Scalar>
  57. inline void igl::cat(
  58. const int dim,
  59. const Eigen::SparseMatrix<Scalar> & A,
  60. const Eigen::SparseMatrix<Scalar> & B,
  61. Eigen::SparseMatrix<Scalar> & C)
  62. {
  63. assert(dim == 1 || dim == 2);
  64. using namespace Eigen;
  65. // Special case if B or A is empty
  66. if(A.size() == 0)
  67. {
  68. C = B;
  69. return;
  70. }
  71. if(B.size() == 0)
  72. {
  73. C = A;
  74. return;
  75. }
  76. DynamicSparseMatrix<Scalar, RowMajor> dyn_C;
  77. if(dim == 1)
  78. {
  79. assert(A.cols() == B.cols());
  80. dyn_C.resize(A.rows()+B.rows(),A.cols());
  81. }else if(dim == 2)
  82. {
  83. assert(A.rows() == B.rows());
  84. dyn_C.resize(A.rows(),A.cols()+B.cols());
  85. }else
  86. {
  87. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  88. }
  89. dyn_C.reserve(A.nonZeros()+B.nonZeros());
  90. // Iterate over outside of A
  91. for(int k=0; k<A.outerSize(); ++k)
  92. {
  93. // Iterate over inside
  94. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  95. {
  96. dyn_C.coeffRef(it.row(),it.col()) += it.value();
  97. }
  98. }
  99. // Iterate over outside of B
  100. for(int k=0; k<B.outerSize(); ++k)
  101. {
  102. // Iterate over inside
  103. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  104. {
  105. int r = (dim == 1 ? A.rows()+it.row() : it.row());
  106. int c = (dim == 2 ? A.cols()+it.col() : it.col());
  107. dyn_C.coeffRef(r,c) += it.value();
  108. }
  109. }
  110. C = SparseMatrix<Scalar>(dyn_C);
  111. }
  112. template <typename Derived, class MatC>
  113. inline void igl::cat(
  114. const int dim,
  115. const Eigen::MatrixBase<Derived> & A,
  116. const Eigen::MatrixBase<Derived> & B,
  117. MatC & C)
  118. {
  119. assert(dim == 1 || dim == 2);
  120. // Special case if B or A is empty
  121. if(A.size() == 0)
  122. {
  123. C = B;
  124. return;
  125. }
  126. if(B.size() == 0)
  127. {
  128. C = A;
  129. return;
  130. }
  131. if(dim == 1)
  132. {
  133. assert(A.cols() == B.cols());
  134. C.resize(A.rows()+B.rows(),A.cols());
  135. C << A,B;
  136. }else if(dim == 2)
  137. {
  138. assert(A.rows() == B.rows());
  139. C.resize(A.rows(),A.cols()+B.cols());
  140. C << A,B;
  141. }else
  142. {
  143. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  144. }
  145. }
  146. template <class Mat>
  147. inline Mat igl::cat(const int dim, const Mat & A, const Mat & B)
  148. {
  149. assert(dim == 1 || dim == 2);
  150. Mat C;
  151. igl::cat(dim,A,B,C);
  152. return C;
  153. }
  154. template <class Mat>
  155. inline void cat(const std::vector<std::vector< Mat > > & A, Mat & C)
  156. {
  157. using namespace igl;
  158. using namespace std;
  159. // Start with empty matrix
  160. C.resize(0,0);
  161. for(typename vector<vector< Mat > >::const_iterator rit = A.begin(); rit != A.end(); rit++)
  162. {
  163. // Concatenate each row horizontally
  164. // Start with empty matrix
  165. Mat row(0,0);
  166. for(typename vector<vector< Mat > >::iterator cit = A.begin(); rit != A.end(); rit++)
  167. {
  168. row = cat(2,row,*cit);
  169. }
  170. // Concatenate rows vertically
  171. C = cat(1,C,row);
  172. }
  173. }
  174. #endif