crouzeix_raviart_massmatrix.cpp 3.3 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 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 "crouzeix_raviart_massmatrix.h"
  9. #include "unique_simplices.h"
  10. #include "oriented_facets.h"
  11. #include "is_edge_manifold.h"
  12. #include "doublearea.h"
  13. #include <cassert>
  14. #include <vector>
  15. template <typename MT, typename DerivedV, typename DerivedF, typename DerivedE, typename DerivedEMAP>
  16. void igl::crouzeix_raviart_massmatrix(
  17. const Eigen::MatrixBase<DerivedV> & V,
  18. const Eigen::MatrixBase<DerivedF> & F,
  19. Eigen::SparseMatrix<MT> & M,
  20. Eigen::PlainObjectBase<DerivedE> & E,
  21. Eigen::PlainObjectBase<DerivedEMAP> & EMAP)
  22. {
  23. // All occurances of directed "facets"
  24. Eigen::MatrixXi allE;
  25. oriented_facets(F,allE);
  26. Eigen::VectorXi _1;
  27. unique_simplices(allE,E,_1,EMAP);
  28. return crouzeix_raviart_massmatrix(V,F,E,EMAP,M);
  29. }
  30. template <typename MT, typename DerivedV, typename DerivedF, typename DerivedE, typename DerivedEMAP>
  31. void igl::crouzeix_raviart_massmatrix(
  32. const Eigen::MatrixBase<DerivedV> & V,
  33. const Eigen::MatrixBase<DerivedF> & F,
  34. const Eigen::MatrixBase<DerivedE> & E,
  35. const Eigen::MatrixBase<DerivedEMAP> & EMAP,
  36. Eigen::SparseMatrix<MT> & M)
  37. {
  38. using namespace Eigen;
  39. using namespace std;
  40. // Mesh should be edge-manifold (TODO: replace `is_edge_manifold` with
  41. // `is_facet_manifold`)
  42. assert(F.cols() != 3 || is_edge_manifold(F));
  43. // number of elements (triangles)
  44. const int m = F.rows();
  45. // Get triangle areas/volumes
  46. VectorXd TA;
  47. // Element simplex size
  48. const int ss = F.cols();
  49. switch(ss)
  50. {
  51. default:
  52. assert(false && "Unsupported simplex size");
  53. case 3:
  54. doublearea(V,F,TA);
  55. TA *= 0.5;
  56. break;
  57. case 4:
  58. volume(V,F,TA);
  59. break;
  60. }
  61. vector<Triplet<MT> > MIJV(ss*m);
  62. assert(EMAP.size() == m*ss);
  63. for(int f = 0;f<m;f++)
  64. {
  65. for(int c = 0;c<ss;c++)
  66. {
  67. MIJV[f+m*c] = Triplet<MT>(EMAP(f+m*c),EMAP(f+m*c),TA(f)/(double)(ss));
  68. }
  69. }
  70. M.resize(E.rows(),E.rows());
  71. M.setFromTriplets(MIJV.begin(),MIJV.end());
  72. }
  73. #ifdef IGL_STATIC_LIBRARY
  74. // Explicit template instantiation
  75. template void igl::crouzeix_raviart_massmatrix<double, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  76. template void igl::crouzeix_raviart_massmatrix<float, Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<float, 0, int>&);
  77. #endif