cr_vector_laplacian.cpp 5.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2020 Oded Stein <[email protected]>
  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 "cr_vector_laplacian.h"
  9. #include <vector>
  10. #include "orient_halfedges.h"
  11. #include "doublearea.h"
  12. #include "squared_edge_lengths.h"
  13. template <typename DerivedV, typename DerivedF, typename DerivedE,
  14. typename DerivedOE, typename ScalarL>
  15. IGL_INLINE void
  16. igl::cr_vector_laplacian(
  17. const Eigen::MatrixBase<DerivedV>& V,
  18. const Eigen::MatrixBase<DerivedF>& F,
  19. const Eigen::MatrixBase<DerivedE>& E,
  20. const Eigen::MatrixBase<DerivedOE>& oE,
  21. Eigen::SparseMatrix<ScalarL>& L)
  22. {
  23. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, Eigen::Dynamic>
  24. l_sq;
  25. squared_edge_lengths(V, F, l_sq);
  26. cr_vector_laplacian_intrinsic(F, l_sq, E, oE, L);
  27. }
  28. template <typename DerivedV, typename DerivedF, typename DerivedE,
  29. typename DerivedOE, typename ScalarL>
  30. IGL_INLINE void
  31. igl::cr_vector_laplacian(
  32. const Eigen::MatrixBase<DerivedV>& V,
  33. const Eigen::MatrixBase<DerivedF>& F,
  34. Eigen::PlainObjectBase<DerivedE>& E,
  35. Eigen::PlainObjectBase<DerivedOE>& oE,
  36. Eigen::SparseMatrix<ScalarL>& L)
  37. {
  38. if(E.rows()!=F.rows() || E.cols()!=F.cols() || oE.rows()!=F.rows() ||
  39. oE.cols()!=F.cols()) {
  40. orient_halfedges(F, E, oE);
  41. }
  42. const Eigen::PlainObjectBase<DerivedE>& cE = E;
  43. const Eigen::PlainObjectBase<DerivedOE>& coE = oE;
  44. cr_vector_laplacian(V, F, cE, coE, L);
  45. }
  46. template <typename DerivedF, typename DerivedL_sq, typename DerivedE,
  47. typename DerivedOE, typename ScalarL>
  48. IGL_INLINE void
  49. igl::cr_vector_laplacian_intrinsic(
  50. const Eigen::MatrixBase<DerivedF>& F,
  51. const Eigen::MatrixBase<DerivedL_sq>& l_sq,
  52. const Eigen::MatrixBase<DerivedE>& E,
  53. const Eigen::MatrixBase<DerivedOE>& oE,
  54. Eigen::SparseMatrix<ScalarL>& L)
  55. {
  56. Eigen::Matrix<typename DerivedL_sq::Scalar, Eigen::Dynamic, Eigen::Dynamic>
  57. dA;
  58. DerivedL_sq l_sqrt = l_sq.array().sqrt().matrix();
  59. doublearea(l_sqrt, dA);
  60. cr_vector_laplacian_intrinsic(F, l_sq, dA, E, oE, L);
  61. }
  62. template <typename DerivedF, typename DerivedL_sq, typename DeriveddA,
  63. typename DerivedE, typename DerivedOE, typename ScalarL>
  64. IGL_INLINE void
  65. igl::cr_vector_laplacian_intrinsic(
  66. const Eigen::MatrixBase<DerivedF>& F,
  67. const Eigen::MatrixBase<DerivedL_sq>& l_sq,
  68. const Eigen::MatrixBase<DeriveddA>& dA,
  69. const Eigen::MatrixBase<DerivedE>& E,
  70. const Eigen::MatrixBase<DerivedOE>& oE,
  71. Eigen::SparseMatrix<ScalarL>& L)
  72. {
  73. assert(F.cols()==3 && "Faces have three vertices");
  74. assert(E.rows()==F.rows() && E.cols()==F.cols() && oE.rows()==F.rows() &&
  75. oE.cols()==F.cols() && "Wrong dimension in edge vectors");
  76. assert(l_sq.rows()==F.rows() && l_sq.cols()==3 && "l_sq dimensions wrong");
  77. assert(dA.size()==F.rows() && "dA dimensions wrong");
  78. const Eigen::Index m = F.rows();
  79. const typename DerivedE::Scalar nE = E.maxCoeff() + 1;
  80. std::vector<Eigen::Triplet<ScalarL> > tripletList;
  81. tripletList.reserve(10*3*m);
  82. for(Eigen::Index f=0; f<m; ++f) {
  83. for(int e=0; e<3; ++e) {
  84. const ScalarL eij=l_sq(f,e), ejk=l_sq(f,(e+1)%3),
  85. eki=l_sq(f,(e+2)%3); //These are squared quantities.
  86. const ScalarL lens = sqrt(eij*eki);
  87. const ScalarL o = oE(f,e)*oE(f,(e+2)%3);
  88. tripletList.emplace_back(E(f,e), E(f,e), 2./dA(f) * eij);
  89. tripletList.emplace_back(E(f,e)+nE, E(f,e)+nE, 2./dA(f) * eij);
  90. const ScalarL Dijki = o * pow(eij-ejk+eki,2)/(2.*lens*dA(f));
  91. tripletList.emplace_back(E(f,e), E(f,(e+2)%3), Dijki);
  92. tripletList.emplace_back(E(f,(e+2)%3), E(f,e), Dijki);
  93. tripletList.emplace_back(E(f,e)+nE, E(f,(e+2)%3)+nE, Dijki);
  94. tripletList.emplace_back(E(f,(e+2)%3)+nE, E(f,e)+nE, Dijki);
  95. const ScalarL Dijkiperp = -o * (eij-ejk+eki)/lens;
  96. tripletList.emplace_back(E(f,e), E(f,(e+2)%3)+nE, Dijkiperp);
  97. tripletList.emplace_back(E(f,(e+2)%3)+nE, E(f,e), Dijkiperp);
  98. tripletList.emplace_back(E(f,e)+nE, E(f,(e+2)%3), -Dijkiperp);
  99. tripletList.emplace_back(E(f,(e+2)%3), E(f,e)+nE, -Dijkiperp);
  100. }
  101. }
  102. L.resize(2*nE, 2*nE);
  103. L.setFromTriplets(tripletList.begin(), tripletList.end());
  104. }
  105. #ifdef IGL_STATIC_LIBRARY
  106. // Explicit template instantiation
  107. template void igl::cr_vector_laplacian<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<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::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::SparseMatrix<double, 0, int>&);
  108. template void igl::cr_vector_laplacian_intrinsic<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::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  109. #endif