crouzeix_raviart_massmatrix.cpp 3.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869
  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 "all_edges.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::PlainObjectBase<DerivedV> & V,
  18. const Eigen::PlainObjectBase<DerivedF> & F,
  19. Eigen::SparseMatrix<MT> & M,
  20. Eigen::PlainObjectBase<DerivedE> & E,
  21. Eigen::PlainObjectBase<DerivedEMAP> & EMAP)
  22. {
  23. // All occurances of directed edges
  24. Eigen::MatrixXi allE;
  25. all_edges(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::PlainObjectBase<DerivedV> & V,
  33. const Eigen::PlainObjectBase<DerivedF> & F,
  34. const Eigen::PlainObjectBase<DerivedE> & E,
  35. const Eigen::PlainObjectBase<DerivedEMAP> & EMAP,
  36. Eigen::SparseMatrix<MT> & M)
  37. {
  38. using namespace igl;
  39. using namespace Eigen;
  40. using namespace std;
  41. assert(F.cols() == 3);
  42. // Mesh should be edge-manifold
  43. assert(is_edge_manifold(V,F));
  44. // number of elements (triangles)
  45. int m = F.rows();
  46. // Get triangle areas
  47. VectorXd TA;
  48. doublearea(V,F,TA);
  49. vector<Triplet<MT> > MIJV(3*m);
  50. for(int f = 0;f<m;f++)
  51. {
  52. for(int c = 0;c<3;c++)
  53. {
  54. MIJV[f+m*c] = Triplet<MT>(EMAP(f+m*c),EMAP(f+m*c),TA(f)*0.5);
  55. }
  56. }
  57. M.resize(E.rows(),E.rows());
  58. M.setFromTriplets(MIJV.begin(),MIJV.end());
  59. }
  60. #ifdef IGL_STATIC_LIBRARY
  61. // Explicit instanciation
  62. 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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  63. 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::PlainObjectBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<float, 0, int>&);
  64. #endif