cat.cpp 4.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158
  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. // This **must** be DynamicSparseMatrix, otherwise this implementation is
  37. // insanely slow
  38. DynamicSparseMatrix<Scalar, RowMajor> dyn_C;
  39. if(dim == 1)
  40. {
  41. assert(A.cols() == B.cols());
  42. dyn_C.resize(A.rows()+B.rows(),A.cols());
  43. }else if(dim == 2)
  44. {
  45. assert(A.rows() == B.rows());
  46. dyn_C.resize(A.rows(),A.cols()+B.cols());
  47. }else
  48. {
  49. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  50. }
  51. dyn_C.reserve(A.nonZeros()+B.nonZeros());
  52. // Iterate over outside of A
  53. for(int k=0; k<A.outerSize(); ++k)
  54. {
  55. // Iterate over inside
  56. for(typename SparseMatrix<Scalar>::InnerIterator it (A,k); it; ++it)
  57. {
  58. dyn_C.coeffRef(it.row(),it.col()) += it.value();
  59. }
  60. }
  61. // Iterate over outside of B
  62. for(int k=0; k<B.outerSize(); ++k)
  63. {
  64. // Iterate over inside
  65. for(typename SparseMatrix<Scalar>::InnerIterator it (B,k); it; ++it)
  66. {
  67. int r = (dim == 1 ? A.rows()+it.row() : it.row());
  68. int c = (dim == 2 ? A.cols()+it.col() : it.col());
  69. dyn_C.coeffRef(r,c) += it.value();
  70. }
  71. }
  72. C = SparseMatrix<Scalar>(dyn_C);
  73. }
  74. template <typename Derived, class MatC>
  75. IGL_INLINE void igl::cat(
  76. const int dim,
  77. const Eigen::MatrixBase<Derived> & A,
  78. const Eigen::MatrixBase<Derived> & B,
  79. MatC & C)
  80. {
  81. assert(dim == 1 || dim == 2);
  82. // Special case if B or A is empty
  83. if(A.size() == 0)
  84. {
  85. C = B;
  86. return;
  87. }
  88. if(B.size() == 0)
  89. {
  90. C = A;
  91. return;
  92. }
  93. if(dim == 1)
  94. {
  95. assert(A.cols() == B.cols());
  96. C.resize(A.rows()+B.rows(),A.cols());
  97. C << A,B;
  98. }else if(dim == 2)
  99. {
  100. assert(A.rows() == B.rows());
  101. C.resize(A.rows(),A.cols()+B.cols());
  102. C << A,B;
  103. }else
  104. {
  105. fprintf(stderr,"cat.h: Error: Unsupported dimension %d\n",dim);
  106. }
  107. }
  108. template <class Mat>
  109. IGL_INLINE Mat igl::cat(const int dim, const Mat & A, const Mat & B)
  110. {
  111. assert(dim == 1 || dim == 2);
  112. Mat C;
  113. igl::cat(dim,A,B,C);
  114. return C;
  115. }
  116. template <class Mat>
  117. IGL_INLINE void igl::cat(const std::vector<std::vector< Mat > > & A, Mat & C)
  118. {
  119. using namespace std;
  120. // Start with empty matrix
  121. C.resize(0,0);
  122. for(typename vector<vector< Mat > >::const_iterator rit = A.begin(); rit != A.end(); rit++)
  123. {
  124. // Concatenate each row horizontally
  125. // Start with empty matrix
  126. Mat row(0,0);
  127. for(typename vector<vector< Mat > >::iterator cit = A.begin(); rit != A.end(); rit++)
  128. {
  129. row = cat(2,row,*cit);
  130. }
  131. // Concatenate rows vertically
  132. C = cat(1,C,row);
  133. }
  134. }
  135. #ifdef IGL_STATIC_LIBRARY
  136. // Explicit template specialization
  137. // generated by autoexplicit.sh
  138. 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&);
  139. // generated by autoexplicit.sh
  140. 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&);
  141. // generated by autoexplicit.sh
  142. 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&);
  143. 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>&);
  144. 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&);
  145. 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&);
  146. #endif