massmatrix.cpp 4.1 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394
  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 "massmatrix_intrinsic.h"
  10. #include "edge_lengths.h"
  11. #include "normalize_row_sums.h"
  12. #include "sparse.h"
  13. #include "doublearea.h"
  14. #include "repmat.h"
  15. #include <Eigen/Geometry>
  16. #include <iostream>
  17. template <typename DerivedV, typename DerivedF, typename Scalar>
  18. IGL_INLINE void igl::massmatrix(
  19. const Eigen::MatrixBase<DerivedV> & V,
  20. const Eigen::MatrixBase<DerivedF> & F,
  21. const MassMatrixType type,
  22. Eigen::SparseMatrix<Scalar>& M)
  23. {
  24. using namespace Eigen;
  25. using namespace std;
  26. const int n = V.rows();
  27. const int m = F.rows();
  28. const int simplex_size = F.cols();
  29. MassMatrixType eff_type = type;
  30. // Use voronoi of for triangles by default, otherwise barycentric
  31. if(type == MASSMATRIX_TYPE_DEFAULT)
  32. {
  33. eff_type = (simplex_size == 3?MASSMATRIX_TYPE_VORONOI:MASSMATRIX_TYPE_BARYCENTRIC);
  34. }
  35. // Not yet supported
  36. assert(type!=MASSMATRIX_TYPE_FULL);
  37. if(simplex_size == 3)
  38. {
  39. // Triangles
  40. // edge lengths numbered same as opposite vertices
  41. Matrix<Scalar,Dynamic,3> l;
  42. igl::edge_lengths(V,F,l);
  43. return massmatrix_intrinsic(l,F,type,M);
  44. }else if(simplex_size == 4)
  45. {
  46. Matrix<int,Dynamic,1> MI;
  47. Matrix<int,Dynamic,1> MJ;
  48. Matrix<Scalar,Dynamic,1> MV;
  49. assert(V.cols() == 3);
  50. assert(eff_type == MASSMATRIX_TYPE_BARYCENTRIC);
  51. MI.resize(m*4,1); MJ.resize(m*4,1); MV.resize(m*4,1);
  52. MI.block(0*m,0,m,1) = F.col(0);
  53. MI.block(1*m,0,m,1) = F.col(1);
  54. MI.block(2*m,0,m,1) = F.col(2);
  55. MI.block(3*m,0,m,1) = F.col(3);
  56. MJ = MI;
  57. // loop over tets
  58. for(int i = 0;i<m;i++)
  59. {
  60. // http://en.wikipedia.org/wiki/Tetrahedron#Volume
  61. Matrix<Scalar,3,1> v0m3,v1m3,v2m3;
  62. v0m3.head(V.cols()) = V.row(F(i,0)) - V.row(F(i,3));
  63. v1m3.head(V.cols()) = V.row(F(i,1)) - V.row(F(i,3));
  64. v2m3.head(V.cols()) = V.row(F(i,2)) - V.row(F(i,3));
  65. Scalar v = fabs(v0m3.dot(v1m3.cross(v2m3)))/6.0;
  66. MV(i+0*m) = v/4.0;
  67. MV(i+1*m) = v/4.0;
  68. MV(i+2*m) = v/4.0;
  69. MV(i+3*m) = v/4.0;
  70. }
  71. sparse(MI,MJ,MV,n,n,M);
  72. }else
  73. {
  74. // Unsupported simplex size
  75. assert(false && "Unsupported simplex size");
  76. }
  77. }
  78. #ifdef IGL_STATIC_LIBRARY
  79. // Explicit template instantiation
  80. // generated by autoexplicit.sh
  81. template void igl::massmatrix<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MassMatrixType, Eigen::SparseMatrix<double, 0, int>&);
  82. // generated by autoexplicit.sh
  83. 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>&);
  84. // generated by autoexplicit.sh
  85. 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>&);
  86. 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>&);
  87. 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>&);
  88. #endif