AtA_cached.cpp 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 Daniele Panozzo <daniele.panozzo@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 "AtA_cached.h"
  9. #include <iostream>
  10. #include <vector>
  11. #include <utility>
  12. template <typename Scalar>
  13. IGL_INLINE void igl::AtA_cached_precompute(
  14. const Eigen::SparseMatrix<Scalar>& A,
  15. igl::AtA_cached_data& data,
  16. Eigen::SparseMatrix<Scalar>& AtA)
  17. {
  18. // 1 Compute At (this could be avoided, but performance-wise it will not make a difference)
  19. std::vector<std::vector<int> > Col_RowPtr;
  20. std::vector<std::vector<int> > Col_IndexPtr;
  21. Col_RowPtr.resize(A.cols());
  22. Col_IndexPtr.resize(A.cols());
  23. for (unsigned k=0; k<A.outerSize(); ++k)
  24. {
  25. unsigned outer_index = *(A.outerIndexPtr()+k);
  26. unsigned next_outer_index = (k+1 == A.outerSize()) ? A.nonZeros() : *(A.outerIndexPtr()+k+1);
  27. for (unsigned l=outer_index; l<next_outer_index; ++l)
  28. {
  29. int col = k;
  30. int row = *(A.innerIndexPtr()+l);
  31. int value_index = l;
  32. assert(col < A.cols());
  33. assert(col >= 0);
  34. assert(row < A.rows());
  35. assert(row >= 0);
  36. assert(value_index >= 0);
  37. assert(value_index < A.nonZeros());
  38. Col_RowPtr[col].push_back(row);
  39. Col_IndexPtr[col].push_back(value_index);
  40. }
  41. }
  42. Eigen::SparseMatrix<Scalar> At = A.transpose();
  43. At.makeCompressed();
  44. AtA = At * A;
  45. AtA.makeCompressed();
  46. assert(AtA.isCompressed());
  47. // If weights are not provided, use 1
  48. if (data.W.size() == 0)
  49. data.W = Eigen::VectorXd::Ones(A.rows());
  50. assert(data.W.size() == A.rows());
  51. data.I_outer.reserve(AtA.outerSize());
  52. data.I_row.reserve(2*AtA.nonZeros());
  53. data.I_col.reserve(2*AtA.nonZeros());
  54. data.I_w.reserve(2*AtA.nonZeros());
  55. // 2 Construct the rules
  56. for (unsigned k=0; k<AtA.outerSize(); ++k)
  57. {
  58. unsigned outer_index = *(AtA.outerIndexPtr()+k);
  59. unsigned next_outer_index = (k+1 == AtA.outerSize()) ? AtA.nonZeros() : *(AtA.outerIndexPtr()+k+1);
  60. for (unsigned l=outer_index; l<next_outer_index; ++l)
  61. {
  62. int col = k;
  63. int row = *(AtA.innerIndexPtr()+l);
  64. int value_index = l;
  65. assert(col < AtA.cols());
  66. assert(col >= 0);
  67. assert(row < AtA.rows());
  68. assert(row >= 0);
  69. assert(value_index >= 0);
  70. assert(value_index < AtA.nonZeros());
  71. data.I_outer.push_back(data.I_row.size());
  72. // Find correspondences
  73. unsigned i=0;
  74. unsigned j=0;
  75. while (i<Col_RowPtr[row].size() && j<Col_RowPtr[col].size())
  76. {
  77. if (Col_RowPtr[row][i] == Col_RowPtr[col][j])
  78. {
  79. data.I_row.push_back(Col_IndexPtr[row][i]);
  80. data.I_col.push_back(Col_IndexPtr[col][j]);
  81. data.I_w.push_back(Col_RowPtr[col][j]);
  82. ++i;
  83. ++j;
  84. } else
  85. if (Col_RowPtr[row][i] > Col_RowPtr[col][j])
  86. ++j;
  87. else
  88. ++i;
  89. }
  90. }
  91. }
  92. data.I_outer.push_back(data.I_row.size()); // makes it more efficient to iterate later on
  93. igl::AtA_cached(A,data,AtA);
  94. }
  95. template <typename Scalar>
  96. IGL_INLINE void igl::AtA_cached(
  97. const Eigen::SparseMatrix<Scalar>& A,
  98. const igl::AtA_cached_data& data,
  99. Eigen::SparseMatrix<Scalar>& AtA)
  100. {
  101. for (unsigned i=0; i<data.I_outer.size()-1; ++i)
  102. {
  103. *(AtA.valuePtr() + i) = 0;
  104. for (unsigned j=data.I_outer[i]; j<data.I_outer[i+1]; ++j)
  105. *(AtA.valuePtr() + i) += *(A.valuePtr() + data.I_row[j]) * data.W[data.I_w[j]] * *(A.valuePtr() + data.I_col[j]);
  106. }
  107. }
  108. #ifdef IGL_STATIC_LIBRARY
  109. template void igl::AtA_cached<double>(Eigen::SparseMatrix<double, 0, int> const&, igl::AtA_cached_data const&, Eigen::SparseMatrix<double, 0, int>&);
  110. template void igl::AtA_cached_precompute<double>(Eigen::SparseMatrix<double, 0, int> const&, igl::AtA_cached_data&, Eigen::SparseMatrix<double, 0, int>&);
  111. #endif