harmonic.cpp 7.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175
  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. if(k>1)
  38. {
  39. massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M);
  40. }
  41. return harmonic(L,M,b,bc,k,W);
  42. }
  43. template <
  44. typename DerivedF,
  45. typename Derivedb,
  46. typename Derivedbc,
  47. typename DerivedW>
  48. IGL_INLINE bool igl::harmonic(
  49. const Eigen::MatrixBase<DerivedF> & F,
  50. const Eigen::MatrixBase<Derivedb> & b,
  51. const Eigen::MatrixBase<Derivedbc> & bc,
  52. const int k,
  53. Eigen::PlainObjectBase<DerivedW> & W)
  54. {
  55. using namespace Eigen;
  56. typedef typename Derivedbc::Scalar Scalar;
  57. SparseMatrix<Scalar> A;
  58. adjacency_matrix(F,A);
  59. // sum each row
  60. SparseVector<Scalar> Asum;
  61. sum(A,1,Asum);
  62. // Convert row sums into diagonal of sparse matrix
  63. SparseMatrix<Scalar> Adiag;
  64. diag(Asum,Adiag);
  65. SparseMatrix<Scalar> L = A-Adiag;
  66. SparseMatrix<Scalar> M;
  67. speye(L.rows(),M);
  68. return harmonic(L,M,b,bc,k,W);
  69. }
  70. template <
  71. typename DerivedL,
  72. typename DerivedM,
  73. typename Derivedb,
  74. typename Derivedbc,
  75. typename DerivedW>
  76. IGL_INLINE bool igl::harmonic(
  77. const Eigen::SparseMatrix<DerivedL> & L,
  78. const Eigen::SparseMatrix<DerivedM> & M,
  79. const Eigen::MatrixBase<Derivedb> & b,
  80. const Eigen::MatrixBase<Derivedbc> & bc,
  81. const int k,
  82. Eigen::PlainObjectBase<DerivedW> & W)
  83. {
  84. const int n = L.rows();
  85. assert(n == L.cols() && "L must be square");
  86. assert((k==1 || n == M.cols() ) && "M must be same size as L");
  87. assert((k==1 || n == M.rows() ) && "M must be square");
  88. assert((k==1 || igl::isdiag(M)) && "Mass matrix should be diagonal");
  89. Eigen::SparseMatrix<DerivedL> Q;
  90. igl::harmonic(L,M,k,Q);
  91. typedef DerivedL Scalar;
  92. min_quad_with_fixed_data<Scalar> data;
  93. min_quad_with_fixed_precompute(Q,b,Eigen::SparseMatrix<Scalar>(),true,data);
  94. W.resize(n,bc.cols());
  95. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> VectorXS;
  96. const VectorXS B = VectorXS::Zero(n,1);
  97. for(int w = 0;w<bc.cols();w++)
  98. {
  99. const VectorXS bcw = bc.col(w);
  100. VectorXS Ww;
  101. if(!min_quad_with_fixed_solve(data,B,bcw,VectorXS(),Ww))
  102. {
  103. return false;
  104. }
  105. W.col(w) = Ww;
  106. }
  107. return true;
  108. }
  109. template <
  110. typename DerivedL,
  111. typename DerivedM,
  112. typename DerivedQ>
  113. IGL_INLINE void igl::harmonic(
  114. const Eigen::SparseMatrix<DerivedL> & L,
  115. const Eigen::SparseMatrix<DerivedM> & M,
  116. const int k,
  117. Eigen::SparseMatrix<DerivedQ> & Q)
  118. {
  119. assert(L.rows() == L.cols()&&"L should be square");
  120. Q = -L;
  121. if(k == 1) return;
  122. assert(L.rows() == M.rows()&&"L should match M's dimensions");
  123. assert(M.rows() == M.cols()&&"M should be square");
  124. Eigen::SparseMatrix<DerivedM> Mi;
  125. invert_diag(M,Mi);
  126. // This is **not** robust for k>2. See KKT system in [Jacobson et al. 2010]
  127. // of the kharmonic function in gptoolbox
  128. for(int p = 1;p<k;p++)
  129. {
  130. Q = (Q*Mi*-L).eval();
  131. }
  132. }
  133. #include "find.h"
  134. template <
  135. typename DerivedV,
  136. typename DerivedF,
  137. typename DerivedQ>
  138. IGL_INLINE void igl::harmonic(
  139. const Eigen::MatrixBase<DerivedV> & V,
  140. const Eigen::MatrixBase<DerivedF> & F,
  141. const int k,
  142. Eigen::SparseMatrix<DerivedQ> & Q)
  143. {
  144. Eigen::SparseMatrix<DerivedQ> L,M;
  145. cotmatrix(V,F,L);
  146. if(k>1)
  147. {
  148. massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M);
  149. }
  150. return harmonic(L,M,k,Q);
  151. }
  152. #ifdef IGL_STATIC_LIBRARY
  153. // Explicit template instantiation
  154. // generated by autoexplicit.sh
  155. 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, 0, -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, 0, -1, -1> >&);
  156. // generated by autoexplicit.sh
  157. 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> >&);
  158. // generated by autoexplicit.sh
  159. 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> >&);
  160. // generated by autoexplicit.sh
  161. 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>&);
  162. // generated by autoexplicit.sh
  163. 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> >&);
  164. // generated by autoexplicit.sh
  165. 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> >&);
  166. 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> >&);
  167. #endif