sparse.cpp 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123
  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 "sparse.h"
  9. #include <iostream>
  10. #include <vector>
  11. template <class IndexVector, class ValueVector, typename T>
  12. IGL_INLINE void igl::sparse(
  13. const IndexVector & I,
  14. const IndexVector & J,
  15. const ValueVector & V,
  16. Eigen::SparseMatrix<T>& X)
  17. {
  18. size_t m = (size_t)I.maxCoeff()+1;
  19. size_t n = (size_t)J.maxCoeff()+1;
  20. return igl::sparse(I,J,V,m,n,X);
  21. }
  22. #include "verbose.h"
  23. template <
  24. class IndexVectorI,
  25. class IndexVectorJ,
  26. class ValueVector,
  27. typename T>
  28. IGL_INLINE void igl::sparse(
  29. const IndexVectorI & I,
  30. const IndexVectorJ & J,
  31. const ValueVector & V,
  32. const size_t m,
  33. const size_t n,
  34. Eigen::SparseMatrix<T>& X)
  35. {
  36. using namespace std;
  37. using namespace Eigen;
  38. assert((int)I.maxCoeff() < (int)m);
  39. assert((int)I.minCoeff() >= 0);
  40. assert((int)J.maxCoeff() < (int)n);
  41. assert((int)J.minCoeff() >= 0);
  42. assert(I.size() == J.size());
  43. assert(J.size() == V.size());
  44. // Really we just need .size() to be the same, but this is safer
  45. assert(I.rows() == J.rows());
  46. assert(J.rows() == V.rows());
  47. assert(I.cols() == J.cols());
  48. assert(J.cols() == V.cols());
  49. //// number of values
  50. //int nv = V.size();
  51. //Eigen::DynamicSparseMatrix<T, Eigen::RowMajor> dyn_X(m,n);
  52. //// over estimate the number of entries
  53. //dyn_X.reserve(I.size());
  54. //for(int i = 0;i < nv;i++)
  55. //{
  56. // dyn_X.coeffRef((int)I(i),(int)J(i)) += (T)V(i);
  57. //}
  58. //X = Eigen::SparseMatrix<T>(dyn_X);
  59. vector<Triplet<T> > IJV;
  60. IJV.reserve(I.size());
  61. for(int x = 0;x<I.size();x++)
  62. {
  63. IJV.push_back(Triplet<T >(I(x),J(x),V(x)));
  64. }
  65. X.resize(m,n);
  66. X.setFromTriplets(IJV.begin(),IJV.end());
  67. }
  68. template <typename DerivedD, typename T>
  69. IGL_INLINE void igl::sparse(
  70. const Eigen::PlainObjectBase<DerivedD>& D,
  71. Eigen::SparseMatrix<T>& X)
  72. {
  73. assert(false && "Obsolete. Just call D.sparseView() directly");
  74. using namespace std;
  75. using namespace Eigen;
  76. vector<Triplet<T> > DIJV;
  77. const int m = D.rows();
  78. const int n = D.cols();
  79. for(int i = 0;i<m;i++)
  80. {
  81. for(int j = 0;j<n;j++)
  82. {
  83. if(D(i,j)!=0)
  84. {
  85. DIJV.push_back(Triplet<T>(i,j,D(i,j)));
  86. }
  87. }
  88. }
  89. X.resize(m,n);
  90. X.setFromTriplets(DIJV.begin(),DIJV.end());
  91. }
  92. template <typename DerivedD>
  93. IGL_INLINE Eigen::SparseMatrix<typename DerivedD::Scalar > igl::sparse(
  94. const Eigen::PlainObjectBase<DerivedD>& D)
  95. {
  96. assert(false && "Obsolete. Just call D.sparseView() directly");
  97. Eigen::SparseMatrix<typename DerivedD::Scalar > X;
  98. igl::sparse(D,X);
  99. return X;
  100. }
  101. #ifdef IGL_STATIC_LIBRARY
  102. // Explicit template instantiation
  103. // generated by autoexplicit.sh
  104. #ifndef WIN32
  105. //template void igl::sparse<Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<bool>, Eigen::Array<bool, -1, 2, 0, -1, 2> >, bool>(Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<bool>, Eigen::Array<bool, -1, 2, 0, -1, 2> > const&, unsigned long, unsigned long, Eigen::SparseMatrix<bool, 0, int>&);
  106. //template void igl::sparse<Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<bool>, Eigen::Array<bool, -1, 3, 0, -1, 3> >, bool>(Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<bool>, Eigen::Array<bool, -1, 3, 0, -1, 3> > const&, unsigned long, unsigned long, Eigen::SparseMatrix<bool, 0, int>&);
  107. #if EIGEN_VERSION_AT_LEAST(3,3,0)
  108. #else
  109. //template void igl::sparse<Eigen::CwiseNullaryOp<Eigen::internal::linspaced_op<int, true>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<bool>, Eigen::Array<bool, -1, 1, 0, -1, 1> >, bool>(Eigen::CwiseNullaryOp<Eigen::internal::linspaced_op<int, true>, Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<bool>, Eigen::Array<bool, -1, 1, 0, -1, 1> > const&, unsigned long, unsigned long, Eigen::SparseMatrix<bool, 0, int>&);
  110. #endif
  111. #endif
  112. template void igl::sparse<Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, std::complex<double> >(Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, size_t, size_t, Eigen::SparseMatrix<std::complex<double>, 0, int>&);
  113. template void igl::sparse<Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, double>(Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, Eigen::SparseMatrix<double, 0, int>&);
  114. template void igl::sparse<Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, double>(Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, size_t, size_t, Eigen::SparseMatrix<double, 0, int>&);
  115. #endif