massmatrix_intrinsic.cpp 5.9 KB

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