massmatrix.cpp 6.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168
  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 "massmatrix.h"
  9. #include "normalize_row_sums.h"
  10. #include "sparse.h"
  11. #include "repmat.h"
  12. #include <Eigen/Geometry>
  13. #include <iostream>
  14. template <typename DerivedV, typename DerivedF, typename Scalar>
  15. IGL_INLINE void igl::massmatrix(
  16. const Eigen::MatrixBase<DerivedV> & V,
  17. const Eigen::MatrixBase<DerivedF> & F,
  18. const MassMatrixType type,
  19. Eigen::SparseMatrix<Scalar>& M)
  20. {
  21. using namespace Eigen;
  22. using namespace std;
  23. const int n = V.rows();
  24. const int m = F.rows();
  25. const int simplex_size = F.cols();
  26. MassMatrixType eff_type = type;
  27. // Use voronoi of for triangles by default, otherwise barycentric
  28. if(type == MASSMATRIX_TYPE_DEFAULT)
  29. {
  30. eff_type = (simplex_size == 3?MASSMATRIX_TYPE_VORONOI:MASSMATRIX_TYPE_BARYCENTRIC);
  31. }
  32. // Not yet supported
  33. assert(type!=MASSMATRIX_TYPE_FULL);
  34. Matrix<int,Dynamic,1> MI;
  35. Matrix<int,Dynamic,1> MJ;
  36. Matrix<Scalar,Dynamic,1> MV;
  37. if(simplex_size == 3)
  38. {
  39. // Triangles
  40. // edge lengths numbered same as opposite vertices
  41. Matrix<Scalar,Dynamic,3> l(m,3);
  42. // loop over faces
  43. for(int i = 0;i<m;i++)
  44. {
  45. l(i,0) = sqrt((V.row(F(i,1))-V.row(F(i,2))).array().pow(2).sum());
  46. l(i,1) = sqrt((V.row(F(i,2))-V.row(F(i,0))).array().pow(2).sum());
  47. l(i,2) = sqrt((V.row(F(i,0))-V.row(F(i,1))).array().pow(2).sum());
  48. }
  49. // semiperimeters
  50. Matrix<Scalar,Dynamic,1> s = l.rowwise().sum()*0.5;
  51. assert(s.rows() == m);
  52. // Heron's forumal for area
  53. Matrix<Scalar,Dynamic,1> dblA(m);
  54. for(int i = 0;i<m;i++)
  55. {
  56. dblA(i) = 2.0*sqrt(s(i)*(s(i)-l(i,0))*(s(i)-l(i,1))*(s(i)-l(i,2)));
  57. }
  58. switch(eff_type)
  59. {
  60. case MASSMATRIX_TYPE_BARYCENTRIC:
  61. // diagonal entries for each face corner
  62. MI.resize(m*3,1); MJ.resize(m*3,1); MV.resize(m*3,1);
  63. MI.block(0*m,0,m,1) = F.col(0);
  64. MI.block(1*m,0,m,1) = F.col(1);
  65. MI.block(2*m,0,m,1) = F.col(2);
  66. MJ = MI;
  67. repmat(dblA,3,1,MV);
  68. MV.array() /= 6.0;
  69. break;
  70. case MASSMATRIX_TYPE_VORONOI:
  71. {
  72. // diagonal entries for each face corner
  73. // http://www.alecjacobson.com/weblog/?p=874
  74. MI.resize(m*3,1); MJ.resize(m*3,1); MV.resize(m*3,1);
  75. MI.block(0*m,0,m,1) = F.col(0);
  76. MI.block(1*m,0,m,1) = F.col(1);
  77. MI.block(2*m,0,m,1) = F.col(2);
  78. MJ = MI;
  79. // Holy shit this needs to be cleaned up and optimized
  80. Matrix<Scalar,Dynamic,3> cosines(m,3);
  81. cosines.col(0) =
  82. (l.col(2).array().pow(2)+l.col(1).array().pow(2)-l.col(0).array().pow(2))/(l.col(1).array()*l.col(2).array()*2.0);
  83. cosines.col(1) =
  84. (l.col(0).array().pow(2)+l.col(2).array().pow(2)-l.col(1).array().pow(2))/(l.col(2).array()*l.col(0).array()*2.0);
  85. cosines.col(2) =
  86. (l.col(1).array().pow(2)+l.col(0).array().pow(2)-l.col(2).array().pow(2))/(l.col(0).array()*l.col(1).array()*2.0);
  87. Matrix<Scalar,Dynamic,3> barycentric = cosines.array() * l.array();
  88. normalize_row_sums(barycentric,barycentric);
  89. Matrix<Scalar,Dynamic,3> partial = barycentric;
  90. partial.col(0).array() *= dblA.array() * 0.5;
  91. partial.col(1).array() *= dblA.array() * 0.5;
  92. partial.col(2).array() *= dblA.array() * 0.5;
  93. Matrix<Scalar,Dynamic,3> quads(partial.rows(),partial.cols());
  94. quads.col(0) = (partial.col(1)+partial.col(2))*0.5;
  95. quads.col(1) = (partial.col(2)+partial.col(0))*0.5;
  96. quads.col(2) = (partial.col(0)+partial.col(1))*0.5;
  97. quads.col(0) = (cosines.col(0).array()<0).select( 0.25*dblA,quads.col(0));
  98. quads.col(1) = (cosines.col(0).array()<0).select(0.125*dblA,quads.col(1));
  99. quads.col(2) = (cosines.col(0).array()<0).select(0.125*dblA,quads.col(2));
  100. quads.col(0) = (cosines.col(1).array()<0).select(0.125*dblA,quads.col(0));
  101. quads.col(1) = (cosines.col(1).array()<0).select(0.25*dblA,quads.col(1));
  102. quads.col(2) = (cosines.col(1).array()<0).select(0.125*dblA,quads.col(2));
  103. quads.col(0) = (cosines.col(2).array()<0).select(0.125*dblA,quads.col(0));
  104. quads.col(1) = (cosines.col(2).array()<0).select(0.125*dblA,quads.col(1));
  105. quads.col(2) = (cosines.col(2).array()<0).select( 0.25*dblA,quads.col(2));
  106. MV.block(0*m,0,m,1) = quads.col(0);
  107. MV.block(1*m,0,m,1) = quads.col(1);
  108. MV.block(2*m,0,m,1) = quads.col(2);
  109. break;
  110. }
  111. case MASSMATRIX_TYPE_FULL:
  112. assert(false && "Implementation incomplete");
  113. break;
  114. default:
  115. assert(false && "Unknown Mass matrix eff_type");
  116. }
  117. }else if(simplex_size == 4)
  118. {
  119. assert(V.cols() == 3);
  120. assert(eff_type == MASSMATRIX_TYPE_BARYCENTRIC);
  121. MI.resize(m*4,1); MJ.resize(m*4,1); MV.resize(m*4,1);
  122. MI.block(0*m,0,m,1) = F.col(0);
  123. MI.block(1*m,0,m,1) = F.col(1);
  124. MI.block(2*m,0,m,1) = F.col(2);
  125. MI.block(3*m,0,m,1) = F.col(3);
  126. MJ = MI;
  127. // loop over tets
  128. for(int i = 0;i<m;i++)
  129. {
  130. // http://en.wikipedia.org/wiki/Tetrahedron#Volume
  131. Matrix<Scalar,3,1> v0m3 = V.row(F(i,0)) - V.row(F(i,3));
  132. Matrix<Scalar,3,1> v1m3 = V.row(F(i,1)) - V.row(F(i,3));
  133. Matrix<Scalar,3,1> v2m3 = V.row(F(i,2)) - V.row(F(i,3));
  134. Scalar v = fabs(v0m3.dot(v1m3.cross(v2m3)))/6.0;
  135. MV(i+0*m) = v/4.0;
  136. MV(i+1*m) = v/4.0;
  137. MV(i+2*m) = v/4.0;
  138. MV(i+3*m) = v/4.0;
  139. }
  140. }else
  141. {
  142. // Unsupported simplex size
  143. assert(false && "Unsupported simplex size");
  144. }
  145. sparse(MI,MJ,MV,n,n,M);
  146. }
  147. #ifdef IGL_STATIC_LIBRARY
  148. // Explicit template specialization
  149. // generated by autoexplicit.sh
  150. template void igl::massmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 4, 0, -1, 4>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 4, 0, -1, 4> > const&, igl::MassMatrixType, Eigen::SparseMatrix<double, 0, int>&);
  151. // generated by autoexplicit.sh
  152. template void igl::massmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, igl::MassMatrixType, Eigen::SparseMatrix<double, 0, int>&);
  153. template void igl::massmatrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MassMatrixType, Eigen::SparseMatrix<double, 0, int>&);
  154. #endif