harmonic.cpp 6.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125
  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 "cotmatrix.h"
  10. #include "massmatrix.h"
  11. #include "invert_diag.h"
  12. #include "adjacency_matrix.h"
  13. #include "sum.h"
  14. #include "diag.h"
  15. #include "min_quad_with_fixed.h"
  16. #include <Eigen/Sparse>
  17. template <
  18. typename DerivedV,
  19. typename DerivedF,
  20. typename Derivedb,
  21. typename Derivedbc,
  22. typename DerivedW>
  23. IGL_INLINE bool igl::harmonic(
  24. const Eigen::PlainObjectBase<DerivedV> & V,
  25. const Eigen::PlainObjectBase<DerivedF> & F,
  26. const Eigen::PlainObjectBase<Derivedb> & b,
  27. const Eigen::PlainObjectBase<Derivedbc> & bc,
  28. const int k,
  29. Eigen::PlainObjectBase<DerivedW> & W)
  30. {
  31. using namespace Eigen;
  32. typedef typename DerivedV::Scalar Scalar;
  33. typedef Matrix<Scalar,Dynamic,1> VectorXS;
  34. SparseMatrix<Scalar> Q;
  35. SparseMatrix<Scalar> L,M,Mi;
  36. cotmatrix(V,F,L);
  37. switch(F.cols())
  38. {
  39. case 3:
  40. massmatrix(V,F,MASSMATRIX_TYPE_VORONOI,M);
  41. break;
  42. case 4:
  43. default:
  44. massmatrix(V,F,MASSMATRIX_TYPE_BARYCENTRIC,M);
  45. break;
  46. }
  47. invert_diag(M,Mi);
  48. Q = -L;
  49. for(int p = 1;p<k;p++)
  50. {
  51. Q = (Q*Mi*-L).eval();
  52. }
  53. min_quad_with_fixed_data<Scalar> data;
  54. min_quad_with_fixed_precompute(Q,b,SparseMatrix<Scalar>(),true,data);
  55. W.resize(V.rows(),bc.cols());
  56. const VectorXS B = VectorXS::Zero(V.rows(),1);
  57. for(int w = 0;w<bc.cols();w++)
  58. {
  59. const VectorXS bcw = bc.col(w);
  60. VectorXS Ww;
  61. if(!min_quad_with_fixed_solve(data,B,bcw,VectorXS(),Ww))
  62. {
  63. return false;
  64. }
  65. W.col(w) = Ww;
  66. }
  67. return true;
  68. }
  69. template <
  70. typename DerivedF,
  71. typename Derivedb,
  72. typename Derivedbc,
  73. typename DerivedW>
  74. IGL_INLINE bool igl::harmonic(
  75. const Eigen::PlainObjectBase<DerivedF> & F,
  76. const Eigen::PlainObjectBase<Derivedb> & b,
  77. const Eigen::PlainObjectBase<Derivedbc> & bc,
  78. const int k,
  79. Eigen::PlainObjectBase<DerivedW> & W)
  80. {
  81. using namespace Eigen;
  82. typedef typename Derivedbc::Scalar Scalar;
  83. typedef Matrix<Scalar,Dynamic,1> VectorXS;
  84. SparseMatrix<Scalar> Q;
  85. SparseMatrix<Scalar> A;
  86. adjacency_matrix(F,A);
  87. // sum each row
  88. SparseVector<Scalar> Asum;
  89. sum(A,1,Asum);
  90. // Convert row sums into diagonal of sparse matrix
  91. SparseMatrix<Scalar> Adiag;
  92. diag(Asum,Adiag);
  93. // Build uniform laplacian
  94. Q = -A+Adiag;
  95. min_quad_with_fixed_data<Scalar> data;
  96. min_quad_with_fixed_precompute(Q,b,SparseMatrix<Scalar>(),true,data);
  97. W.resize(A.rows(),bc.cols());
  98. const VectorXS B = VectorXS::Zero(A.rows(),1);
  99. for(int w = 0;w<bc.cols();w++)
  100. {
  101. const VectorXS bcw = bc.col(w);
  102. VectorXS Ww;
  103. if(!min_quad_with_fixed_solve(data,B,bcw,VectorXS(),Ww))
  104. {
  105. return false;
  106. }
  107. W.col(w) = Ww;
  108. }
  109. return true;
  110. }
  111. #ifdef IGL_STATIC_LIBRARY
  112. // Explicit template specialization
  113. 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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  114. 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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  115. 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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  116. 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::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  117. 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::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  118. 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::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  119. #endif