grad.cpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 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.h"
  9. #include <Eigen/Geometry>
  10. #include <vector>
  11. template <typename DerivedV, typename DerivedF>
  12. IGL_INLINE void igl::grad(const Eigen::PlainObjectBase<DerivedV>&V,
  13. const Eigen::PlainObjectBase<DerivedF>&F,
  14. Eigen::SparseMatrix<typename DerivedV::Scalar> &G)
  15. {
  16. Eigen::Matrix<typename DerivedV::Scalar,Eigen::Dynamic,3>
  17. eperp21(F.rows(),3), eperp13(F.rows(),3);
  18. for (int i=0;i<F.rows();++i)
  19. {
  20. // renaming indices of vertices of triangles for convenience
  21. int i1 = F(i,0);
  22. int i2 = F(i,1);
  23. int i3 = F(i,2);
  24. // #F x 3 matrices of triangle edge vectors, named after opposite vertices
  25. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> v32 = V.row(i3) - V.row(i2);
  26. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> v13 = V.row(i1) - V.row(i3);
  27. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> v21 = V.row(i2) - V.row(i1);
  28. // area of parallelogram is twice area of triangle
  29. // area of parallelogram is || v1 x v2 ||
  30. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> n = v32.cross(v13);
  31. // This does correct l2 norm of rows, so that it contains #F list of twice
  32. // triangle areas
  33. double dblA = std::sqrt(n.dot(n));
  34. // now normalize normals to get unit normals
  35. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> u = n / dblA;
  36. // rotate each vector 90 degrees around normal
  37. double norm21 = std::sqrt(v21.dot(v21));
  38. double norm13 = std::sqrt(v13.dot(v13));
  39. eperp21.row(i) = u.cross(v21);
  40. eperp21.row(i) = eperp21.row(i) / std::sqrt(eperp21.row(i).dot(eperp21.row(i)));
  41. eperp21.row(i) *= norm21 / dblA;
  42. eperp13.row(i) = u.cross(v13);
  43. eperp13.row(i) = eperp13.row(i) / std::sqrt(eperp13.row(i).dot(eperp13.row(i)));
  44. eperp13.row(i) *= norm13 / dblA;
  45. }
  46. std::vector<int> rs;
  47. rs.reserve(F.rows()*4*3);
  48. std::vector<int> cs;
  49. cs.reserve(F.rows()*4*3);
  50. std::vector<double> vs;
  51. vs.reserve(F.rows()*4*3);
  52. // row indices
  53. for(int r=0;r<3;r++)
  54. {
  55. for(int j=0;j<4;j++)
  56. {
  57. for(int i=r*F.rows();i<(r+1)*F.rows();i++) rs.push_back(i);
  58. }
  59. }
  60. // column indices
  61. for(int r=0;r<3;r++)
  62. {
  63. for(int i=0;i<F.rows();i++) cs.push_back(F(i,1));
  64. for(int i=0;i<F.rows();i++) cs.push_back(F(i,0));
  65. for(int i=0;i<F.rows();i++) cs.push_back(F(i,2));
  66. for(int i=0;i<F.rows();i++) cs.push_back(F(i,0));
  67. }
  68. // values
  69. for(int i=0;i<F.rows();i++) vs.push_back(eperp13(i,0));
  70. for(int i=0;i<F.rows();i++) vs.push_back(-eperp13(i,0));
  71. for(int i=0;i<F.rows();i++) vs.push_back(eperp21(i,0));
  72. for(int i=0;i<F.rows();i++) vs.push_back(-eperp21(i,0));
  73. for(int i=0;i<F.rows();i++) vs.push_back(eperp13(i,1));
  74. for(int i=0;i<F.rows();i++) vs.push_back(-eperp13(i,1));
  75. for(int i=0;i<F.rows();i++) vs.push_back(eperp21(i,1));
  76. for(int i=0;i<F.rows();i++) vs.push_back(-eperp21(i,1));
  77. for(int i=0;i<F.rows();i++) vs.push_back(eperp13(i,2));
  78. for(int i=0;i<F.rows();i++) vs.push_back(-eperp13(i,2));
  79. for(int i=0;i<F.rows();i++) vs.push_back(eperp21(i,2));
  80. for(int i=0;i<F.rows();i++) vs.push_back(-eperp21(i,2));
  81. // create sparse gradient operator matrix
  82. G.resize(3*F.rows(),V.rows());
  83. std::vector<Eigen::Triplet<typename DerivedV::Scalar> > triplets;
  84. for (int i=0;i<(int)vs.size();++i)
  85. {
  86. triplets.push_back(Eigen::Triplet<typename DerivedV::Scalar>(rs[i],cs[i],vs[i]));
  87. }
  88. G.setFromTriplets(triplets.begin(), triplets.end());
  89. }
  90. #ifdef IGL_STATIC_LIBRARY
  91. // Explicit template specialization
  92. // template void igl::grad<double, int>(Eigen::Matrix<double, -1, -1, 0, -1,-1> const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&,Eigen::SparseMatrix<double, 0, int>&);
  93. template void igl::grad<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::SparseMatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>::Scalar, 0, int>&);
  94. //template void igl::grad<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::SparseMatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>::Scalar, 0, int>&);
  95. template void igl::grad<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -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::SparseMatrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 0, int>&);
  96. #endif