harmonic.cpp 3.0 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071
  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 Eigen;
  29. typedef typename DerivedV::Scalar Scalar;
  30. typedef Matrix<Scalar,Dynamic,1> VectorXS;
  31. SparseMatrix<Scalar> L,M,Mi;
  32. cotmatrix(V,F,L);
  33. switch(F.cols())
  34. {
  35. case 3:
  36. massmatrix(V,F,MASSMATRIX_TYPE_VORONOI,M);
  37. break;
  38. case 4:
  39. default:
  40. massmatrix(V,F,MASSMATRIX_TYPE_BARYCENTRIC,M);
  41. break;
  42. }
  43. invert_diag(M,Mi);
  44. SparseMatrix<Scalar> Q = -L;
  45. for(int p = 1;p<k;p++)
  46. {
  47. Q = (Q*Mi*-L).eval();
  48. }
  49. const VectorXS B = VectorXS::Zero(V.rows(),1);
  50. min_quad_with_fixed_data<Scalar> data;
  51. min_quad_with_fixed_precompute(Q,b,SparseMatrix<Scalar>(),true,data);
  52. W.resize(V.rows(),bc.cols());
  53. for(int w = 0;w<bc.cols();w++)
  54. {
  55. const VectorXS bcw = bc.col(w);
  56. VectorXS Ww;
  57. if(!min_quad_with_fixed_solve(data,B,bcw,VectorXS(),Ww))
  58. {
  59. return false;
  60. }
  61. W.col(w) = Ww;
  62. }
  63. return true;
  64. }
  65. #ifdef IGL_STATIC_LIBRARY
  66. // Explicit template specialization
  67. 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> >&);
  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. #endif