grad.cpp 8.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216
  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. #include "per_face_normals.h"
  12. template <typename DerivedV, typename DerivedF>
  13. IGL_INLINE void igl::grad(const Eigen::PlainObjectBase<DerivedV>&V,
  14. const Eigen::PlainObjectBase<DerivedF>&F,
  15. Eigen::SparseMatrix<typename DerivedV::Scalar> &G)
  16. {
  17. Eigen::Matrix<typename DerivedV::Scalar,Eigen::Dynamic,3>
  18. eperp21(F.rows(),3), eperp13(F.rows(),3);
  19. for (int i=0;i<F.rows();++i)
  20. {
  21. // renaming indices of vertices of triangles for convenience
  22. int i1 = F(i,0);
  23. int i2 = F(i,1);
  24. int i3 = F(i,2);
  25. // #F x 3 matrices of triangle edge vectors, named after opposite vertices
  26. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> v32 = V.row(i3) - V.row(i2);
  27. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> v13 = V.row(i1) - V.row(i3);
  28. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> v21 = V.row(i2) - V.row(i1);
  29. // area of parallelogram is twice area of triangle
  30. // area of parallelogram is || v1 x v2 ||
  31. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> n = v32.cross(v13);
  32. // This does correct l2 norm of rows, so that it contains #F list of twice
  33. // triangle areas
  34. double dblA = std::sqrt(n.dot(n));
  35. // now normalize normals to get unit normals
  36. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> u = n / dblA;
  37. // rotate each vector 90 degrees around normal
  38. double norm21 = std::sqrt(v21.dot(v21));
  39. double norm13 = std::sqrt(v13.dot(v13));
  40. eperp21.row(i) = u.cross(v21);
  41. eperp21.row(i) = eperp21.row(i) / std::sqrt(eperp21.row(i).dot(eperp21.row(i)));
  42. eperp21.row(i) *= norm21 / dblA;
  43. eperp13.row(i) = u.cross(v13);
  44. eperp13.row(i) = eperp13.row(i) / std::sqrt(eperp13.row(i).dot(eperp13.row(i)));
  45. eperp13.row(i) *= norm13 / dblA;
  46. }
  47. std::vector<int> rs;
  48. rs.reserve(F.rows()*4*3);
  49. std::vector<int> cs;
  50. cs.reserve(F.rows()*4*3);
  51. std::vector<double> vs;
  52. vs.reserve(F.rows()*4*3);
  53. // row indices
  54. for(int r=0;r<3;r++)
  55. {
  56. for(int j=0;j<4;j++)
  57. {
  58. for(int i=r*F.rows();i<(r+1)*F.rows();i++) rs.push_back(i);
  59. }
  60. }
  61. // column indices
  62. for(int r=0;r<3;r++)
  63. {
  64. for(int i=0;i<F.rows();i++) cs.push_back(F(i,1));
  65. for(int i=0;i<F.rows();i++) cs.push_back(F(i,0));
  66. for(int i=0;i<F.rows();i++) cs.push_back(F(i,2));
  67. for(int i=0;i<F.rows();i++) cs.push_back(F(i,0));
  68. }
  69. // values
  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(-eperp13(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(-eperp21(i,0));
  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(-eperp13(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(-eperp21(i,1));
  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(-eperp13(i,2));
  80. for(int i=0;i<F.rows();i++) vs.push_back(eperp21(i,2));
  81. for(int i=0;i<F.rows();i++) vs.push_back(-eperp21(i,2));
  82. // create sparse gradient operator matrix
  83. G.resize(3*F.rows(),V.rows());
  84. std::vector<Eigen::Triplet<typename DerivedV::Scalar> > triplets;
  85. for (int i=0;i<(int)vs.size();++i)
  86. {
  87. triplets.push_back(Eigen::Triplet<typename DerivedV::Scalar>(rs[i],cs[i],vs[i]));
  88. }
  89. G.setFromTriplets(triplets.begin(), triplets.end());
  90. }
  91. /*
  92. IGL_INLINE void igl::grad_tet(const Eigen::PlainObjectBase<DerivedV>&V,
  93. const Eigen::PlainObjectBase<DerivedF>&T,
  94. Eigen::SparseMatrix<typename DerivedV::Scalar> &G)
  95. */
  96. template <typename DerivedV, typename DerivedF>
  97. IGL_INLINE void igl::grad_tet(const Eigen::PlainObjectBase<DerivedV>&V,
  98. const Eigen::PlainObjectBase<DerivedF>&T,
  99. Eigen::SparseMatrix<typename DerivedV::Scalar> &G,
  100. bool uniform) {
  101. using namespace Eigen;
  102. assert(T.cols() == 4);
  103. const int n = V.rows(); int m = T.rows();
  104. /*
  105. F = [ ...
  106. T(:,1) T(:,2) T(:,3); ...
  107. T(:,1) T(:,3) T(:,4); ...
  108. T(:,1) T(:,4) T(:,2); ...
  109. T(:,2) T(:,4) T(:,3)]; */
  110. MatrixXi F(4*m,3);
  111. for (int i = 0; i < m; i++) {
  112. F.row(0*m + i) << T(i,0), T(i,1), T(i,2);
  113. F.row(1*m + i) << T(i,0), T(i,2), T(i,3);
  114. F.row(2*m + i) << T(i,0), T(i,3), T(i,1);
  115. F.row(3*m + i) << T(i,1), T(i,3), T(i,2);
  116. }
  117. // compute volume of each tet
  118. VectorXd vol; igl::volume(V,T,vol);
  119. VectorXd A(F.rows());
  120. MatrixXd N(F.rows(),3);
  121. if (!uniform) {
  122. // compute tetrahedron face normals
  123. igl::per_face_normals(V,F,N); int norm_rows = N.rows();
  124. for (int i = 0; i < norm_rows; i++)
  125. N.row(i) /= N.row(i).norm();
  126. igl::doublearea(V,F,A); A/=2.;
  127. } else {
  128. // Use a uniform tetrahedra as a reference, with the same volume as the original one:
  129. //
  130. // Use normals of the uniform tet (V = h*[0,0,0;1,0,0;0.5,sqrt(3)/2.,0;0.5,sqrt(3)/6.,sqrt(2)/sqrt(3)])
  131. // 0 0 1.0000
  132. // 0.8165 -0.4714 -0.3333
  133. // 0 0.9428 -0.3333
  134. // -0.8165 -0.4714 -0.3333
  135. for (int i = 0; i < m; i++) {
  136. N.row(0*m+i) << 0,0,1;
  137. double a = sqrt(2)*std::cbrt(3*vol(i)); // area of a face in a uniform tet with volume = vol(i)
  138. A(0*m+i) = (pow(a,2)*sqrt(3))/4.;
  139. }
  140. for (int i = 0; i < m; i++) {
  141. N.row(1*m+i) << 0.8165,-0.4714,-0.3333;
  142. double a = sqrt(2)*std::cbrt(3*vol(i));
  143. A(1*m+i) = (pow(a,2)*sqrt(3))/4.;
  144. }
  145. for (int i = 0; i < m; i++) {
  146. N.row(2*m+i) << 0,0.9428,-0.3333;
  147. double a = sqrt(2)*std::cbrt(3*vol(i));
  148. A(2*m+i) = (pow(a,2)*sqrt(3))/4.;
  149. }
  150. for (int i = 0; i < m; i++) {
  151. N.row(3*m+i) << -0.8165,-0.4714,-0.3333;
  152. double a = sqrt(2)*std::cbrt(3*vol(i));
  153. A(3*m+i) = (pow(a,2)*sqrt(3))/4.;
  154. }
  155. }
  156. /* G = sparse( ...
  157. [0*m + repmat(1:m,1,4) ...
  158. 1*m + repmat(1:m,1,4) ...
  159. 2*m + repmat(1:m,1,4)], ...
  160. repmat([T(:,4);T(:,2);T(:,3);T(:,1)],3,1), ...
  161. repmat(A./(3*repmat(vol,4,1)),3,1).*N(:), ...
  162. 3*m,n);*/
  163. std::vector<Triplet<double> > G_t;
  164. for (int i = 0; i < 4*m; i++) {
  165. int T_j; // j indexes : repmat([T(:,4);T(:,2);T(:,3);T(:,1)],3,1)
  166. switch (i/m) {
  167. case 0:
  168. T_j = 3;
  169. break;
  170. case 1:
  171. T_j = 1;
  172. break;
  173. case 2:
  174. T_j = 2;
  175. break;
  176. case 3:
  177. T_j = 0;
  178. break;
  179. }
  180. int i_idx = i%m;
  181. int j_idx = T(i_idx,T_j);
  182. double val_before_n = A(i)/(3*vol(i_idx));
  183. G_t.push_back(Triplet<double>(0*m+i_idx, j_idx, val_before_n * N(i,0)));
  184. G_t.push_back(Triplet<double>(1*m+i_idx, j_idx, val_before_n * N(i,1)));
  185. G_t.push_back(Triplet<double>(2*m+i_idx, j_idx, val_before_n * N(i,2)));
  186. }
  187. G.resize(3*m,n);
  188. G.setFromTriplets(G_t.begin(), G_t.end());
  189. }
  190. #ifdef IGL_STATIC_LIBRARY
  191. // Explicit template specialization
  192. // 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>&);
  193. 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>&);
  194. //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>&);
  195. 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>&);
  196. #endif