sparse.h 2.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990
  1. #ifndef IGL_SPARSE_H
  2. #define IGL_SPARSE_H
  3. #include <Eigen/Dense>
  4. #include <Eigen/Sparse>
  5. namespace igl
  6. {
  7. // Build a sparse matrix from list of indices and values (I,J,V), functions
  8. // like the sparse function in matlab
  9. //
  10. // Templates:
  11. // IndexVector list of indices, value should be non-negative and should
  12. // expect to be cast to an index. Must implement operator(i) to retrieve
  13. // ith element
  14. // ValueVector list of values, value should be expect to be cast to type
  15. // T. Must implement operator(i) to retrieve ith element
  16. // T should be a eigen sparse matrix primitive type like int or double
  17. // Input:
  18. // I nnz vector of row indices of non zeros entries in X
  19. // J nnz vector of column indices of non zeros entries in X
  20. // V nnz vector of non-zeros entries in X
  21. // Optional:
  22. // m number of rows
  23. // n number of cols
  24. // Outputs:
  25. // X m by n matrix of type T whose entries are to be found
  26. //
  27. template <class IndexVector, class ValueVector, typename T>
  28. inline void sparse(
  29. const IndexVector & I,
  30. const IndexVector & J,
  31. const ValueVector & V,
  32. Eigen::SparseMatrix<T>& X);
  33. template <class IndexVector, class ValueVector, typename T>
  34. inline void sparse(
  35. const IndexVector & I,
  36. const IndexVector & J,
  37. const ValueVector & V,
  38. const size_t m,
  39. const size_t n,
  40. Eigen::SparseMatrix<T>& X);
  41. }
  42. // Implementation
  43. template <class IndexVector, class ValueVector, typename T>
  44. inline void igl::sparse(
  45. const IndexVector & I,
  46. const IndexVector & J,
  47. const ValueVector & V,
  48. Eigen::SparseMatrix<T>& X)
  49. {
  50. size_t m = (size_t)I.maxCoeff()+1;
  51. size_t n = (size_t)J.maxCoeff()+1;
  52. return igl::sparse(I,J,V,m,n,X);
  53. }
  54. #include "verbose.h"
  55. template <class IndexVector, class ValueVector, typename T>
  56. inline void igl::sparse(
  57. const IndexVector & I,
  58. const IndexVector & J,
  59. const ValueVector & V,
  60. const size_t m,
  61. const size_t n,
  62. Eigen::SparseMatrix<T>& X)
  63. {
  64. assert((int)I.maxCoeff() < (int)m);
  65. assert((int)I.minCoeff() >= 0);
  66. assert((int)J.maxCoeff() < (int)n);
  67. assert((int)J.minCoeff() >= 0);
  68. assert(I.size() == J.size());
  69. assert(J.size() == V.size());
  70. // Really we just need .size() to be the same, but this is safer
  71. assert(I.rows() == J.rows());
  72. assert(J.rows() == V.rows());
  73. assert(I.cols() == J.cols());
  74. assert(J.cols() == V.cols());
  75. // number of values
  76. int nv = V.size();
  77. Eigen::DynamicSparseMatrix<T, Eigen::RowMajor>
  78. dyn_X(m,n);
  79. for(int i = 0;i < nv;i++)
  80. {
  81. dyn_X.coeffRef((int)I(i),(int)J(i)) += (T)V(i);
  82. }
  83. X = Eigen::SparseMatrix<T>(dyn_X);
  84. }
  85. #endif