grad_intrinsic.cpp 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2018 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 "grad_intrinsic.h"
  9. #include "grad.h"
  10. template <typename Derivedl, typename DerivedF, typename Gtype>
  11. IGL_INLINE void igl::grad_intrinsic(
  12. const Eigen::MatrixBase<Derivedl>&l,
  13. const Eigen::MatrixBase<DerivedF>&F,
  14. Eigen::SparseMatrix<Gtype> &G)
  15. {
  16. assert(F.cols() ==3 && "Only triangles supported");
  17. // number of vertices
  18. const int n = F.maxCoeff()+1;
  19. // number of faces
  20. const int m = F.rows();
  21. // JD: There is a pretty subtle bug when using a fixed column size for this matrix.
  22. // When calling igl::grad(V, ...), the two code paths `grad_tet` and `grad_tri`
  23. // will be compiled. It turns out that `igl::grad_tet` calls `igl::volume`, which
  24. // reads the coordinates of the `V` matrix into `RowVector3d`. If the matrix `V`
  25. // has a known compile-time size of 2, this produces a compilation error when
  26. // libigl is compiled in header-only mode. In static mode this doesn't happen
  27. // because the matrix `V` is probably implicitly copied into a `Eigen::MatrixXd`.
  28. // This is a situation that could be solved using `if constexpr` in C++17.
  29. // In C++11, the alternative is to use SFINAE and `std::enable_if` (ugh).
  30. typedef Eigen::Matrix<Gtype,Eigen::Dynamic,Eigen::Dynamic> MatrixX2S;
  31. MatrixX2S V2 = MatrixX2S::Zero(3*m,2);
  32. // 1=[x,y]
  33. // /\
  34. // l3 / \ l2
  35. // / \
  36. // / \
  37. // 2-----------3
  38. // l1
  39. //
  40. // x = (l2²-l1²-l3²)/(-2*l1)
  41. // y = sqrt(l3² - x²)
  42. //
  43. //
  44. // Place 3rd vertex at [l(:,1) 0]
  45. V2.block(2*m,0,m,1) = l.col(0);
  46. // Place second vertex at [0 0]
  47. // Place third vertex at [x y]
  48. V2.block(0,0,m,1) =
  49. (l.col(1).cwiseAbs2()-l.col(0).cwiseAbs2()-l.col(2).cwiseAbs2()).array()/
  50. (-2.*l.col(0)).array();
  51. V2.block(0,1,m,1) =
  52. (l.col(2).cwiseAbs2() - V2.block(0,0,m,1).cwiseAbs2()).array().sqrt();
  53. DerivedF F2(F.rows(),F.cols());
  54. std::vector<Eigen::Triplet<Gtype> > Pijv;
  55. Pijv.reserve(F.size());
  56. for(int f = 0;f<m;f++)
  57. {
  58. for(int c = 0;c<F.cols();c++)
  59. {
  60. F2(f,c) = f+c*m;
  61. Pijv.emplace_back(F2(f,c),F(f,c),1);
  62. }
  63. }
  64. Eigen::SparseMatrix<Gtype> P(m*3,n);
  65. P.setFromTriplets(Pijv.begin(),Pijv.end());
  66. Eigen::SparseMatrix<Gtype> G2;
  67. grad(V2,F2,G2);
  68. G = G2*P;
  69. }
  70. #ifdef IGL_STATIC_LIBRARY
  71. // Explicit template instantiation
  72. // generated by autoexplicit.sh
  73. template void igl::grad_intrinsic<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  74. // generated by autoexplicit.sh
  75. template void igl::grad_intrinsic<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>&);
  76. #endif