massmatrix_intrinsic.cpp 6.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134
  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_intrinsic.h"
  9. #include "edge_lengths.h"
  10. #include "sparse.h"
  11. #include "doublearea.h"
  12. #include "repmat.h"
  13. #include <Eigen/Geometry>
  14. #include <iostream>
  15. #include <cassert>
  16. template <typename Derivedl, typename DerivedF, typename Scalar>
  17. IGL_INLINE void igl::massmatrix_intrinsic(
  18. const Eigen::MatrixBase<Derivedl> & l,
  19. const Eigen::MatrixBase<DerivedF> & F,
  20. const MassMatrixType type,
  21. Eigen::SparseMatrix<Scalar>& M)
  22. {
  23. const int n = F.maxCoeff()+1;
  24. return massmatrix_intrinsic(l,F,type,n,M);
  25. }
  26. template <typename Derivedl, typename DerivedF, typename Scalar>
  27. IGL_INLINE void igl::massmatrix_intrinsic(
  28. const Eigen::MatrixBase<Derivedl> & l,
  29. const Eigen::MatrixBase<DerivedF> & F,
  30. const MassMatrixType type,
  31. const int n,
  32. Eigen::SparseMatrix<Scalar>& M)
  33. {
  34. MassMatrixType eff_type = type;
  35. const int m = F.rows();
  36. const int simplex_size = F.cols();
  37. // Use voronoi of for triangles by default, otherwise barycentric
  38. if(type == MASSMATRIX_TYPE_DEFAULT)
  39. {
  40. eff_type = (simplex_size == 3?MASSMATRIX_TYPE_VORONOI:MASSMATRIX_TYPE_BARYCENTRIC);
  41. }
  42. assert(F.cols() == 3 && "only triangles supported");
  43. Eigen::Matrix<Scalar ,Eigen::Dynamic,1> dblA;
  44. doublearea(l,0.,dblA);
  45. Eigen::Matrix<typename DerivedF::Scalar ,Eigen::Dynamic,1> MI;
  46. Eigen::Matrix<typename DerivedF::Scalar ,Eigen::Dynamic,1> MJ;
  47. Eigen::Matrix<Scalar ,Eigen::Dynamic,1> MV;
  48. switch(eff_type)
  49. {
  50. case MASSMATRIX_TYPE_BARYCENTRIC:
  51. // diagonal entries for each face corner
  52. MI.resize(m*3,1); MJ.resize(m*3,1); MV.resize(m*3,1);
  53. MI.block(0*m,0,m,1) = F.col(0);
  54. MI.block(1*m,0,m,1) = F.col(1);
  55. MI.block(2*m,0,m,1) = F.col(2);
  56. MJ = MI;
  57. repmat(dblA,3,1,MV);
  58. MV.array() /= 6.0;
  59. break;
  60. case MASSMATRIX_TYPE_VORONOI:
  61. {
  62. // diagonal entries for each face corner
  63. // http://www.alecjacobson.com/weblog/?p=874
  64. MI.resize(m*3,1); MJ.resize(m*3,1); MV.resize(m*3,1);
  65. MI.block(0*m,0,m,1) = F.col(0);
  66. MI.block(1*m,0,m,1) = F.col(1);
  67. MI.block(2*m,0,m,1) = F.col(2);
  68. MJ = MI;
  69. // Holy shit this needs to be cleaned up and optimized
  70. Eigen::Matrix<Scalar ,Eigen::Dynamic,3> cosines(m,3);
  71. cosines.col(0) =
  72. (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);
  73. cosines.col(1) =
  74. (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);
  75. cosines.col(2) =
  76. (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);
  77. Eigen::Matrix<Scalar ,Eigen::Dynamic,3> barycentric = cosines.array() * l.array();
  78. // Replace this: normalize_row_sums(barycentric,barycentric);
  79. barycentric = (barycentric.array().colwise() / barycentric.array().rowwise().sum()).eval();
  80. Eigen::Matrix<Scalar ,Eigen::Dynamic,3> partial = barycentric;
  81. partial.col(0).array() *= dblA.array() * 0.5;
  82. partial.col(1).array() *= dblA.array() * 0.5;
  83. partial.col(2).array() *= dblA.array() * 0.5;
  84. Eigen::Matrix<Scalar ,Eigen::Dynamic,3> quads(partial.rows(),partial.cols());
  85. quads.col(0) = (partial.col(1)+partial.col(2))*0.5;
  86. quads.col(1) = (partial.col(2)+partial.col(0))*0.5;
  87. quads.col(2) = (partial.col(0)+partial.col(1))*0.5;
  88. quads.col(0) = (cosines.col(0).array()<0).select( 0.25*dblA,quads.col(0));
  89. quads.col(1) = (cosines.col(0).array()<0).select(0.125*dblA,quads.col(1));
  90. quads.col(2) = (cosines.col(0).array()<0).select(0.125*dblA,quads.col(2));
  91. quads.col(0) = (cosines.col(1).array()<0).select(0.125*dblA,quads.col(0));
  92. quads.col(1) = (cosines.col(1).array()<0).select(0.25*dblA,quads.col(1));
  93. quads.col(2) = (cosines.col(1).array()<0).select(0.125*dblA,quads.col(2));
  94. quads.col(0) = (cosines.col(2).array()<0).select(0.125*dblA,quads.col(0));
  95. quads.col(1) = (cosines.col(2).array()<0).select(0.125*dblA,quads.col(1));
  96. quads.col(2) = (cosines.col(2).array()<0).select( 0.25*dblA,quads.col(2));
  97. MV.block(0*m,0,m,1) = quads.col(0);
  98. MV.block(1*m,0,m,1) = quads.col(1);
  99. MV.block(2*m,0,m,1) = quads.col(2);
  100. break;
  101. }
  102. case MASSMATRIX_TYPE_FULL:
  103. MI.resize(m*9,1); MJ.resize(m*9,1); MV.resize(m*9,1);
  104. // indicies and values of the element mass matrix entries in the order
  105. // (0,1),(1,0),(1,2),(2,1),(2,0),(0,2),(0,0),(1,1),(2,2);
  106. MI<<F.col(0),F.col(1),F.col(1),F.col(2),F.col(2),F.col(0),F.col(0),F.col(1),F.col(2);
  107. MJ<<F.col(1),F.col(0),F.col(2),F.col(1),F.col(0),F.col(2),F.col(0),F.col(1),F.col(2);
  108. repmat(dblA,9,1,MV);
  109. MV.block(0*m,0,6*m,1) /= 24.0;
  110. MV.block(6*m,0,3*m,1) /= 12.0;
  111. break;
  112. default:
  113. assert(false && "Unknown Mass matrix eff_type");
  114. }
  115. sparse(MI,MJ,MV,n,n,M);
  116. }
  117. #ifdef IGL_STATIC_LIBRARY
  118. // Explicit template instantiation
  119. template void igl::massmatrix_intrinsic<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>&);
  120. template void igl::massmatrix_intrinsic<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>&);
  121. template void igl::massmatrix_intrinsic<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 1, -1, 3>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 1, -1, 3> > const&, igl::MassMatrixType, Eigen::SparseMatrix<double, 0, int>&);
  122. template void igl::massmatrix_intrinsic<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>&);
  123. #endif