harmonic.cpp 6.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 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 "harmonic.h"
  9. #include "adjacency_matrix.h"
  10. #include "cotmatrix.h"
  11. #include "diag.h"
  12. #include "invert_diag.h"
  13. #include "isdiag.h"
  14. #include "massmatrix.h"
  15. #include "min_quad_with_fixed.h"
  16. #include "speye.h"
  17. #include "sum.h"
  18. #include <Eigen/Sparse>
  19. template <
  20. typename DerivedV,
  21. typename DerivedF,
  22. typename Derivedb,
  23. typename Derivedbc,
  24. typename DerivedW>
  25. IGL_INLINE bool igl::harmonic(
  26. const Eigen::MatrixBase<DerivedV> & V,
  27. const Eigen::MatrixBase<DerivedF> & F,
  28. const Eigen::MatrixBase<Derivedb> & b,
  29. const Eigen::MatrixBase<Derivedbc> & bc,
  30. const int k,
  31. Eigen::PlainObjectBase<DerivedW> & W)
  32. {
  33. using namespace Eigen;
  34. typedef typename DerivedV::Scalar Scalar;
  35. SparseMatrix<Scalar> L,M;
  36. cotmatrix(V,F,L);
  37. massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M);
  38. return harmonic(L,M,b,bc,k,W);
  39. }
  40. template <
  41. typename DerivedF,
  42. typename Derivedb,
  43. typename Derivedbc,
  44. typename DerivedW>
  45. IGL_INLINE bool igl::harmonic(
  46. const Eigen::MatrixBase<DerivedF> & F,
  47. const Eigen::MatrixBase<Derivedb> & b,
  48. const Eigen::MatrixBase<Derivedbc> & bc,
  49. const int k,
  50. Eigen::PlainObjectBase<DerivedW> & W)
  51. {
  52. using namespace Eigen;
  53. typedef typename Derivedbc::Scalar Scalar;
  54. SparseMatrix<Scalar> A;
  55. adjacency_matrix(F,A);
  56. // sum each row
  57. SparseVector<Scalar> Asum;
  58. sum(A,1,Asum);
  59. // Convert row sums into diagonal of sparse matrix
  60. SparseMatrix<Scalar> Adiag;
  61. diag(Asum,Adiag);
  62. SparseMatrix<Scalar> L = A-Adiag;
  63. SparseMatrix<Scalar> M;
  64. speye(L.rows(),M);
  65. return harmonic(L,M,b,bc,k,W);
  66. }
  67. template <
  68. typename DerivedL,
  69. typename DerivedM,
  70. typename Derivedb,
  71. typename Derivedbc,
  72. typename DerivedW>
  73. IGL_INLINE bool igl::harmonic(
  74. const Eigen::SparseMatrix<DerivedL> & L,
  75. const Eigen::SparseMatrix<DerivedM> & M,
  76. const Eigen::MatrixBase<Derivedb> & b,
  77. const Eigen::MatrixBase<Derivedbc> & bc,
  78. const int k,
  79. Eigen::PlainObjectBase<DerivedW> & W)
  80. {
  81. const int n = L.rows();
  82. assert(n == L.cols() && "L must be square");
  83. assert(n == M.cols() && "M must be same size as L");
  84. assert(n == M.rows() && "M must be square");
  85. assert(igl::isdiag(M) && "Mass matrix should be diagonal");
  86. Eigen::SparseMatrix<DerivedL> Q;
  87. igl::harmonic(L,M,k,Q);
  88. typedef DerivedL Scalar;
  89. min_quad_with_fixed_data<Scalar> data;
  90. min_quad_with_fixed_precompute(Q,b,Eigen::SparseMatrix<Scalar>(),true,data);
  91. W.resize(n,bc.cols());
  92. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> VectorXS;
  93. const VectorXS B = VectorXS::Zero(n,1);
  94. for(int w = 0;w<bc.cols();w++)
  95. {
  96. const VectorXS bcw = bc.col(w);
  97. VectorXS Ww;
  98. if(!min_quad_with_fixed_solve(data,B,bcw,VectorXS(),Ww))
  99. {
  100. return false;
  101. }
  102. W.col(w) = Ww;
  103. }
  104. return true;
  105. }
  106. template <
  107. typename DerivedL,
  108. typename DerivedM,
  109. typename DerivedQ>
  110. IGL_INLINE void igl::harmonic(
  111. const Eigen::SparseMatrix<DerivedL> & L,
  112. const Eigen::SparseMatrix<DerivedM> & M,
  113. const int k,
  114. Eigen::SparseMatrix<DerivedQ> & Q)
  115. {
  116. assert(L.rows() == L.cols()&&"L should be square");
  117. assert(M.rows() == M.cols()&&"M should be square");
  118. assert(L.rows() == M.rows()&&"L should match M's dimensions");
  119. Q = -L;
  120. if(k == 1) return;
  121. Eigen::SparseMatrix<DerivedM> Mi;
  122. invert_diag(M,Mi);
  123. for(int p = 1;p<k;p++)
  124. {
  125. Q = (Q*Mi*-L).eval();
  126. }
  127. }
  128. #include "find.h"
  129. template <
  130. typename DerivedV,
  131. typename DerivedF,
  132. typename DerivedQ>
  133. IGL_INLINE void igl::harmonic(
  134. const Eigen::MatrixBase<DerivedV> & V,
  135. const Eigen::MatrixBase<DerivedF> & F,
  136. const int k,
  137. Eigen::SparseMatrix<DerivedQ> & Q)
  138. {
  139. Eigen::SparseMatrix<DerivedQ> L,M;
  140. cotmatrix(V,F,L);
  141. massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M);
  142. return harmonic(L,M,k,Q);
  143. }
  144. #ifdef IGL_STATIC_LIBRARY
  145. // Explicit template specialization
  146. // generated by autoexplicit.sh
  147. template bool igl::harmonic<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 1, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> >&);
  148. // generated by autoexplicit.sh
  149. template bool igl::harmonic<Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  150. // generated by autoexplicit.sh
  151. template void igl::harmonic<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, int, Eigen::SparseMatrix<double, 0, int>&);
  152. // generated by autoexplicit.sh
  153. template bool igl::harmonic<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  154. // generated by autoexplicit.sh
  155. template bool igl::harmonic<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  156. #endif