grad.cpp 8.8 KB

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