massmatrix.cpp 5.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <[email protected]>
  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 "sparse.h"
  12. #include "doublearea.h"
  13. #include "volume.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. if(simplex_size == 3)
  36. {
  37. // Triangles
  38. // edge lengths numbered same as opposite vertices
  39. Matrix<Scalar,Dynamic,3> l;
  40. igl::edge_lengths(V,F,l);
  41. return massmatrix_intrinsic(l,F,type,M);
  42. }else if(simplex_size == 4)
  43. {
  44. // tetrahedra
  45. assert(V.cols() == 3 && "vertices must be defined in 3D for tetrahedra");
  46. Matrix<Scalar,Dynamic,1> vol;
  47. volume(V,F,vol);
  48. vol = vol.array().abs();
  49. Matrix<typename DerivedF::Scalar,Dynamic,1> MI;
  50. Matrix<typename DerivedF::Scalar,Dynamic,1> MJ;
  51. Matrix<Scalar,Dynamic,1> MV;
  52. switch (eff_type)
  53. {
  54. case MASSMATRIX_TYPE_BARYCENTRIC:
  55. MI.resize(m*4,1); MJ.resize(m*4,1); MV.resize(m*4,1);
  56. MI.block(0*m,0,m,1) = F.col(0);
  57. MI.block(1*m,0,m,1) = F.col(1);
  58. MI.block(2*m,0,m,1) = F.col(2);
  59. MI.block(3*m,0,m,1) = F.col(3);
  60. MJ = MI;
  61. repmat(vol,4,1,MV);
  62. assert(MV.rows()==m*4&&MV.cols()==1);
  63. MV.array() /= 4.;
  64. break;
  65. case MASSMATRIX_TYPE_VORONOI:
  66. {
  67. assert(false && "Implementation incomplete");
  68. break;
  69. }
  70. case MASSMATRIX_TYPE_FULL:
  71. MI.resize(m*16,1); MJ.resize(m*16,1); MV.resize(m*16,1);
  72. // indicies and values of the element mass matrix entries in the order
  73. // (1,0),(2,0),(3,0),(2,1),(3,1),(0,1),(3,2),(0,2),(1,2),(0,3),(1,3),(2,3),(0,0),(1,1),(2,2),(3,3);
  74. MI<<F.col(1),F.col(2),F.col(3),F.col(2),F.col(3),F.col(0),F.col(3),F.col(0),F.col(1),F.col(0),F.col(1),F.col(2),F.col(0),F.col(1),F.col(2),F.col(3);
  75. MJ<<F.col(0),F.col(0),F.col(0),F.col(1),F.col(1),F.col(1),F.col(2),F.col(2),F.col(2),F.col(3),F.col(3),F.col(3),F.col(0),F.col(1),F.col(2),F.col(3);
  76. repmat(vol,16,1,MV);
  77. assert(MV.rows()==m*16&&MV.cols()==1);
  78. MV.block(0*m,0,12*m,1) /= 20.;
  79. MV.block(12*m,0,4*m,1) /= 10.;
  80. break;
  81. default:
  82. assert(false && "Unknown Mass matrix eff_type");
  83. }
  84. sparse(MI,MJ,MV,n,n,M);
  85. }else
  86. {
  87. // Unsupported simplex size
  88. assert(false && "Unsupported simplex size");
  89. }
  90. }
  91. #ifdef IGL_STATIC_LIBRARY
  92. // Explicit template instantiation
  93. // generated by autoexplicit.sh
  94. template void igl::massmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MassMatrixType, Eigen::SparseMatrix<double, 0, int>&);
  95. // generated by autoexplicit.sh
  96. 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>&);
  97. // generated by autoexplicit.sh
  98. 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>&);
  99. // generated by autoexplicit.sh
  100. 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>&);
  101. 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>&);
  102. 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>&);
  103. #endif