hessian_energy.cpp 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 Alec Jacobson <alecjacobson@gmail.com> and Oded Stein <oded.stein@columbia.edu>
  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 "hessian_energy.h"
  9. #include <vector>
  10. #include <igl/grad.h>
  11. #include <igl/doublearea.h>
  12. #include <igl/massmatrix.h>
  13. #include <igl/boundary_loop.h>
  14. #include <igl/repdiag.h>
  15. template <typename DerivedV, typename DerivedF, typename Scalar>
  16. IGL_INLINE void igl::hessian_energy(
  17. const Eigen::PlainObjectBase<DerivedV> & V,
  18. const Eigen::PlainObjectBase<DerivedF> & F,
  19. Eigen::SparseMatrix<Scalar>& Q)
  20. {
  21. typedef typename DerivedV::Scalar denseScalar;
  22. typedef typename Eigen::Matrix<denseScalar, Eigen::Dynamic, 1> VecXd;
  23. typedef typename Eigen::SparseMatrix<Scalar> SparseMat;
  24. typedef typename Eigen::DiagonalMatrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> DiagMat;
  25. int dim = V.cols();
  26. assert((dim==2 || dim==3) && "The dimension of the vertices should be 2 or 3");
  27. SparseMat M;
  28. igl::massmatrix(V,F,igl::MASSMATRIX_TYPE_VORONOI,M);
  29. //Kill non-interior DOFs
  30. VecXd Mint = M.diagonal();
  31. std::vector<std::vector<int> > bdryLoop;
  32. igl::boundary_loop(F,bdryLoop);
  33. for(const std::vector<int>& loop : bdryLoop)
  34. for(const int& bdryVert : loop)
  35. Mint(bdryVert) = 0.;
  36. //Invert Mint
  37. for(int i=0; i<Mint.rows(); ++i)
  38. if(Mint(i) > 0)
  39. Mint(i) = 1./Mint(i);
  40. //Repeat Mint to form diaginal matrix
  41. DiagMat stackedMinv = Mint.replicate(dim*dim,1).asDiagonal();
  42. //Compute squared Hessian
  43. SparseMat H;
  44. igl::fem_hessian(V,F,H);
  45. Q = H.transpose()*stackedMinv*H;
  46. }
  47. template <typename DerivedV, typename DerivedF, typename Scalar>
  48. IGL_INLINE void igl::fem_hessian(
  49. const Eigen::PlainObjectBase<DerivedV> & V,
  50. const Eigen::PlainObjectBase<DerivedF> & F,
  51. Eigen::SparseMatrix<Scalar>& H)
  52. {
  53. typedef typename DerivedV::Scalar denseScalar;
  54. typedef typename Eigen::Matrix<denseScalar, Eigen::Dynamic, 1> VecXd;
  55. typedef typename Eigen::SparseMatrix<Scalar> SparseMat;
  56. typedef typename Eigen::DiagonalMatrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> DiagMat;
  57. int dim = V.cols();
  58. assert((dim==2 || dim==3) && "The dimension of the vertices should be 2 or 3");
  59. //Construct the combined gradient matric
  60. SparseMat G;
  61. igl::grad(V, F, G, false);
  62. SparseMat GG(F.rows(), dim*V.rows());
  63. GG.reserve(G.nonZeros());
  64. for(int i=0; i<dim; ++i)
  65. GG.middleCols(i*G.cols(), G.cols()) = G.middleRows(i*F.rows(), F.rows());
  66. SparseMat D;
  67. igl::repdiag(GG,dim,D);
  68. //Compute area matrix
  69. VecXd areas;
  70. igl::doublearea(V, F, areas);
  71. DiagMat A = (0.5*areas).replicate(dim,1).asDiagonal();
  72. //Compute FEM Hessian
  73. H = D.transpose()*A*G;
  74. }
  75. #ifdef IGL_STATIC_LIBRARY
  76. // Explicit template instantiation
  77. template void igl::hessian_energy<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  78. template void igl::fem_hessian<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  79. #endif