grad.cpp 8.2 KB

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