hessian_energy.cpp 2.3 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 Alec Jacobson <alecjacobson@gmail.com>
  4. // and Oded Stein <oded.stein@columbia.edu>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla Public License
  7. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  8. // obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "hessian_energy.h"
  10. #include <vector>
  11. #include "hessian.h"
  12. #include "massmatrix.h"
  13. #include "boundary_loop.h"
  14. template <typename DerivedV, typename DerivedF, typename Scalar>
  15. IGL_INLINE void igl::hessian_energy(
  16. const Eigen::MatrixBase<DerivedV> & V,
  17. const Eigen::MatrixBase<DerivedF> & F,
  18. Eigen::SparseMatrix<Scalar>& Q)
  19. {
  20. typedef typename DerivedV::Scalar denseScalar;
  21. typedef typename Eigen::Matrix<denseScalar, Eigen::Dynamic, 1> VecXd;
  22. typedef typename Eigen::SparseMatrix<Scalar> SparseMat;
  23. typedef typename Eigen::DiagonalMatrix
  24. <Scalar, Eigen::Dynamic, Eigen::Dynamic> DiagMat;
  25. int dim = V.cols();
  26. assert((dim==2 || dim==3) &&
  27. "The dimension of the vertices should be 2 or 3");
  28. SparseMat M;
  29. igl::massmatrix(V,F,igl::MASSMATRIX_TYPE_VORONOI,M);
  30. //Kill non-interior DOFs
  31. VecXd Mint = M.diagonal();
  32. std::vector<std::vector<int> > bdryLoop;
  33. igl::boundary_loop(DerivedF(F),bdryLoop);
  34. for(const std::vector<int>& loop : bdryLoop)
  35. for(const int& bdryVert : loop)
  36. Mint(bdryVert) = 0.;
  37. //Invert Mint
  38. for(int i=0; i<Mint.rows(); ++i)
  39. if(Mint(i) > 0)
  40. Mint(i) = 1./Mint(i);
  41. //Repeat Mint to form diaginal matrix
  42. DiagMat stackedMinv = Mint.replicate(dim*dim,1).asDiagonal();
  43. //Compute squared Hessian
  44. SparseMat H;
  45. igl::hessian(V,F,H);
  46. Q = H.transpose()*stackedMinv*H;
  47. }
  48. #ifdef IGL_STATIC_LIBRARY
  49. // Explicit template instantiation
  50. template void igl::hessian_energy<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&, Eigen::SparseMatrix<double, 0, int>&);
  51. #endif