grad.cpp 8.8 KB

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