crouzeix_raviart_cotmatrix.cpp 3.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 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_cotmatrix.h"
  9. #include "unique_simplices.h"
  10. #include "oriented_facets.h"
  11. #include "is_edge_manifold.h"
  12. #include "cotmatrix_entries.h"
  13. template <typename DerivedV, typename DerivedF, typename LT, typename DerivedE, typename DerivedEMAP>
  14. void igl::crouzeix_raviart_cotmatrix(
  15. const Eigen::MatrixBase<DerivedV> & V,
  16. const Eigen::MatrixBase<DerivedF> & F,
  17. Eigen::SparseMatrix<LT> & L,
  18. Eigen::PlainObjectBase<DerivedE> & E,
  19. Eigen::PlainObjectBase<DerivedEMAP> & EMAP)
  20. {
  21. // All occurrences of directed "facets"
  22. Eigen::MatrixXi allE;
  23. oriented_facets(F,allE);
  24. Eigen::VectorXi _1;
  25. unique_simplices(allE,E,_1,EMAP);
  26. return crouzeix_raviart_cotmatrix(V,F,E,EMAP,L);
  27. }
  28. template <typename DerivedV, typename DerivedF, typename DerivedE, typename DerivedEMAP, typename LT>
  29. void igl::crouzeix_raviart_cotmatrix(
  30. const Eigen::MatrixBase<DerivedV> & V,
  31. const Eigen::MatrixBase<DerivedF> & F,
  32. const Eigen::MatrixBase<DerivedE> & E,
  33. const Eigen::MatrixBase<DerivedEMAP> & EMAP,
  34. Eigen::SparseMatrix<LT> & L)
  35. {
  36. // number of rows
  37. const int m = F.rows();
  38. // Element simplex size
  39. const int ss = F.cols();
  40. // Mesh should be edge-manifold
  41. assert(F.cols() != 3 || is_edge_manifold(F));
  42. typedef Eigen::Matrix<LT,Eigen::Dynamic,Eigen::Dynamic> MatrixXS;
  43. MatrixXS C;
  44. cotmatrix_entries(V,F,C);
  45. Eigen::MatrixXi F2E(m,ss);
  46. {
  47. int k =0;
  48. for(int c = 0;c<ss;c++)
  49. {
  50. for(int f = 0;f<m;f++)
  51. {
  52. F2E(f,c) = k++;
  53. }
  54. }
  55. }
  56. // number of entries inserted per facet
  57. const int k = ss*(ss-1)*2;
  58. std::vector<Eigen::Triplet<LT> > LIJV;LIJV.reserve(k*m);
  59. Eigen::VectorXi LI(k),LJ(k),LV(k);
  60. // Compensation factor to match scales in matlab version
  61. double factor = 2.0;
  62. switch(ss)
  63. {
  64. default: assert(false && "unsupported simplex size");
  65. case 3:
  66. factor = 4.0;
  67. LI<<0,1,2,1,2,0,0,1,2,1,2,0;
  68. LJ<<1,2,0,0,1,2,0,1,2,1,2,0;
  69. LV<<2,0,1,2,0,1,2,0,1,2,0,1;
  70. break;
  71. case 4:
  72. factor *= -1.0;
  73. LI<<0,3,3,3,1,2,1,0,1,2,2,0,0,3,3,3,1,2,1,0,1,2,2,0;
  74. LJ<<1,0,1,2,2,0,0,3,3,3,1,2,0,3,3,3,1,2,1,0,1,2,2,0;
  75. LV<<2,3,4,5,0,1,2,3,4,5,0,1,2,3,4,5,0,1,2,3,4,5,0,1;
  76. break;
  77. }
  78. for(int f=0;f<m;f++)
  79. {
  80. for(int c = 0;c<k;c++)
  81. {
  82. LIJV.emplace_back(
  83. EMAP(F2E(f,LI(c))),
  84. EMAP(F2E(f,LJ(c))),
  85. (c<(k/2)?-1.:1.) * factor *C(f,LV(c)));
  86. }
  87. }
  88. L.resize(E.rows(),E.rows());
  89. L.setFromTriplets(LIJV.begin(),LIJV.end());
  90. }
  91. #ifdef IGL_STATIC_LIBRARY
  92. // Explicit template instantiation
  93. // generated by autoexplicit.sh
  94. template void igl::crouzeix_raviart_cotmatrix<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>, double>(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>&);
  95. #endif