harmonic.cpp 3.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172
  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 "min_quad_with_fixed.h"
  13. #include <Eigen/Sparse>
  14. template <
  15. typename DerivedV,
  16. typename DerivedF,
  17. typename Derivedb,
  18. typename Derivedbc,
  19. typename DerivedW>
  20. IGL_INLINE bool igl::harmonic(
  21. const Eigen::PlainObjectBase<DerivedV> & V,
  22. const Eigen::PlainObjectBase<DerivedF> & F,
  23. const Eigen::PlainObjectBase<Derivedb> & b,
  24. const Eigen::PlainObjectBase<Derivedbc> & bc,
  25. const int k,
  26. Eigen::PlainObjectBase<DerivedW> & W)
  27. {
  28. using namespace igl;
  29. using namespace Eigen;
  30. typedef typename DerivedV::Scalar Scalar;
  31. typedef Matrix<Scalar,Dynamic,1> VectorXS;
  32. SparseMatrix<Scalar> L,M,Mi;
  33. cotmatrix(V,F,L);
  34. switch(F.cols())
  35. {
  36. case 3:
  37. massmatrix(V,F,MASSMATRIX_TYPE_VORONOI,M);
  38. break;
  39. case 4:
  40. default:
  41. massmatrix(V,F,MASSMATRIX_TYPE_BARYCENTRIC,M);
  42. break;
  43. }
  44. invert_diag(M,Mi);
  45. SparseMatrix<Scalar> Q = -L;
  46. for(int p = 1;p<k;p++)
  47. {
  48. Q = (Q*Mi*-L).eval();
  49. }
  50. const VectorXS B = VectorXS::Zero(V.rows(),1);
  51. min_quad_with_fixed_data<Scalar> data;
  52. min_quad_with_fixed_precompute(Q,b,SparseMatrix<Scalar>(),true,data);
  53. W.resize(V.rows(),bc.cols());
  54. for(int w = 0;w<bc.cols();w++)
  55. {
  56. const VectorXS bcw = bc.col(w);
  57. VectorXS Ww;
  58. if(!min_quad_with_fixed_solve(data,B,bcw,VectorXS(),Ww))
  59. {
  60. return false;
  61. }
  62. W.col(w) = Ww;
  63. }
  64. return true;
  65. }
  66. #ifdef IGL_STATIC_LIBRARY
  67. // Explicit template instanciation
  68. 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> >&);
  69. 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> >&);
  70. #endif