Manicore
Library to implement schemes on n-dimensionnal manifolds.
pec.hpp
Go to the documentation of this file.
1 #ifndef PEC_HPP
2 #define PEC_HPP
3 
4 #include "dcell_integral.hpp"
5 
6 #include <iostream>
7 
18 namespace Manicore {
19 
22 
24 
34  template<size_t dimension>
35  class PEC {
36  public:
37  PEC(Mesh<dimension> const & mesh,
38  int r,
39  bool use_threads = true,
40  std::array<int,dimension> const * dqr_p = nullptr,
41  std::ostream & output = std::cout
42  );
43 
45 
46  double getTopScaling(size_t i_cell ) const {return _listScales[i_cell];}
48 
49  Eigen::MatrixXd get_mass(size_t k , size_t d ,size_t i_cell ) const;
51 
52  Eigen::MatrixXd get_trace(size_t k , size_t d ,size_t i_cell , size_t j_bd ) const;
54 
55  Eigen::MatrixXd get_starTrace(size_t k , size_t d ,size_t i_cell , size_t j_bd ) const;
56 
57  // Getter for the generic operators matrices
59 
60  const Eigen::MatrixXd & get_diff(size_t l , size_t d ) const {return _list_diff[_cmp_ind(l,d)];}
62 
63  const Eigen::MatrixXd & get_Koszul(size_t l , size_t d ) const {return _list_Koszul[_cmp_ind(l,d)];}
65 
66  const Eigen::MatrixXd & get_diff_as_degr(size_t l , size_t d ) const {return _list_diff_as_degr[_cmp_ind(l,d)];}
68 
69  const Eigen::MatrixXd & get_trimmed(size_t l , size_t d ) const {return _list_trimmed[_cmp_ind(l,d)];}
71 
72  const Eigen::MatrixXd & get_reduced_Koszul_m1(size_t l , size_t d ) const {return _list_reduced_Koszul_m1[_cmp_ind(l,d)];}
73 
74  private:
75  template<size_t d,size_t l> void _fill_lists();
76  inline int _cmp_ind(size_t l, size_t d) const {return _dim_table[d-1]+l;}
77 
78  int _r;
79  std::array<size_t,dimension+1> _dim_table;
80  std::vector<Eigen::MatrixXd> _list_diff, // Image of dPL(r,l,d) inside PL(r-1,l+1,d)
81  _list_Koszul, // Image of kPL(r,l,d) inside PL(r+1,l-1,d)
82  _list_diff_as_degr, // Image of dPL(r,l,d) inside PL(r,l+1,d)
83  _list_trimmed, // Image of PLtrimmed(r,l,d)
84  _list_reduced_Koszul_m1; // Image of kPL(r-1,l,d)
85  std::vector<double> _listScales;
86  template<size_t _dimension,size_t _d>
87  class dCellVariableList {
88  private:
89  std::vector<dCell_mass<_dimension,_d>> _mass;
90  std::vector<dCell_traces<_dimension,_d>> _traces;
91  dCellVariableList<_dimension,_d-1> _dCVL;
92  public:
93  template<size_t d> std::vector<dCell_mass<_dimension,d>> & mass() {
94  if constexpr(d==_d) return _mass;
95  else return _dCVL.template mass<d>();
96  }
97  template<size_t d> std::vector<dCell_mass<_dimension,d>> mass() const {
98  if constexpr(d==_d) return _mass;
99  else return _dCVL.template mass<d>();
100  }
101  template<size_t d> std::vector<dCell_traces<_dimension,d>> & traces() {
102  if constexpr(d==_d) return _traces;
103  else return _dCVL.template traces<d>();
104  }
105  template<size_t d> std::vector<dCell_traces<_dimension,d>> traces() const {
106  if constexpr(d==_d) return _traces;
107  else return _dCVL.template traces<d>();
108  }
109  };
110 
111  template<size_t _dimension>
112  class dCellVariableList<_dimension,1> {
113  private:
114  std::vector<dCell_mass<_dimension,1>> _mass;
115  std::vector<dCell_traces<_dimension,1>> _traces;
116  public:
117  template<size_t d> std::vector<dCell_mass<_dimension,d>> & mass() {
118  static_assert(d == 1);
119  return _mass;
120  }
121  template<size_t d> std::vector<dCell_mass<_dimension,d>> mass() const {
122  static_assert(d == 1);
123  return _mass;
124  }
125  template<size_t d> std::vector<dCell_traces<_dimension,d>> & traces() {
126  static_assert(d == 1);
127  return _traces;
128  }
129  template<size_t d> std::vector<dCell_traces<_dimension,d>> traces() const {
130  static_assert(d == 1);
131  return _traces;
132  }
133  };
134 
135  dCellVariableList<dimension,dimension> _dCellList;
136  };
138 }
139 #endif
140 
Main data structure for the mesh.
Definition: mesh.hpp:54
Implement the discrete spaces of DDR-PEC.
Definition: pec.hpp:35
const Eigen::MatrixXd & get_diff_as_degr(size_t l, size_t d) const
Return the image of the differential operator.
Definition: pec.hpp:66
Eigen::MatrixXd get_starTrace(size_t k, size_t d, size_t i_cell, size_t j_bd) const
Return the Hodge dual of the trace for the Hodge dual of -forms on the i-th d-cell onto its j-th (d-1...
Definition: pec.cpp:195
const Eigen::MatrixXd & get_reduced_Koszul_m1(size_t l, size_t d) const
Return the image of the Koszul operator.
Definition: pec.hpp:72
Eigen::MatrixXd get_trace(size_t k, size_t d, size_t i_cell, size_t j_bd) const
Return the trace for the k-forms on the i-th d-cell onto its j-th (d-1)-neighbour.
Definition: pec.cpp:176
const Eigen::MatrixXd & get_trimmed(size_t l, size_t d) const
Return the image of the trimmed polynomial basis.
Definition: pec.hpp:69
PEC(Mesh< dimension > const &mesh, int r, bool use_threads=true, std::array< int, dimension > const *dqr_p=nullptr, std::ostream &output=std::cout)
Definition: pec.cpp:30
Eigen::MatrixXd get_mass(size_t k, size_t d, size_t i_cell) const
Return the mass matrix for the k-forms on the i-th d-cell in the basis PL(r,k,d)
Definition: pec.cpp:158
const Eigen::MatrixXd & get_Koszul(size_t l, size_t d) const
Return the image of the Koszul operator.
Definition: pec.hpp:63
const Eigen::MatrixXd & get_diff(size_t l, size_t d) const
Return the image of the differential operator.
Definition: pec.hpp:60
double getTopScaling(size_t i_cell) const
Return the characteristic scale of a top dimensional cell .
Definition: pec.hpp:46
Compute the mass matrices and traces operator of a cell.
Definition: maxwell.hpp:21