massmatrix.cpp 6.9 KB

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