HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
ddr_pec.hpp
Go to the documentation of this file.
1#ifndef DDR_PEC_HPP
2#define DDR_PEC_HPP
3
4#include "exterior_cell.hpp"
5#include "mesh.hpp"
6
7namespace HArDCore3D {
8
14 class DDR_PEC
15 {
16 public:
17
19 DDR_PEC(Mesh const & mesh,int r,bool use_threads = true, std::ostream & output = std::cout);
20
21 typedef Eigen::KroneckerProduct<Eigen::MatrixXd,Eigen::MatrixXd> Kronecker;
23 Kronecker get_mass(size_t k, size_t d, size_t i) const;
25 Kronecker get_trace(size_t k, size_t d, size_t i, size_t j) const;
26
28 const Eigen::MatrixXd & get_diff(size_t l, size_t d) const {return _list_diff[_cmp_ind(l,d)];}
29 const Eigen::MatrixXd & get_Koszul(size_t l, size_t d) const {return _list_Koszul[_cmp_ind(l,d)];}
30 const Eigen::MatrixXd & get_diff_as_degr(size_t l, size_t d) const {return _list_diff_as_degr[_cmp_ind(l,d)];}
31 const Eigen::MatrixXd & get_trimmed(size_t l, size_t d) const {return _list_trimmed[_cmp_ind(l,d)];}
32 const Eigen::MatrixXd & get_reduced_Koszul_m1(size_t l, size_t d) const {return _list_reduced_Koszul_m1[_cmp_ind(l,d)];}
33
35 Eigen::Matrix<double,-1,-1,0,3,3> get_exterior_trace(size_t k, size_t d, size_t i_cell) const;
36
38 Eigen::Matrix<double,-1,-1,0,3,3> get_hodge_star(size_t k, size_t d, size_t i_cell) const;
39
41 Eigen::MatrixXd discrete_hodge_star(size_t k, size_t d, size_t i_cell) const;
42 // Takes the hodge star of a d-k form, i.e. a k form on the k basis to the original d-k form on the d-k basis (deals with scaling and sign change) on the mesh element of dimension d with index i_cell
44 Eigen::MatrixXd discrete_inv_hodge_star(size_t k, size_t d, size_t i_cell) const;
45
46 // Gives the evaluation into the global basis (before the scaling)
47 double evaluate_scalar_basis(Eigen::Vector3d const &x, size_t d, size_t i_cell, int i_basis) const;
49 Eigen::Matrix<double,-1,1,0,3,1> evaluate_basis(Eigen::Vector3d const &x, size_t k, size_t d, size_t i_cell, int i_basis) const;
51 Eigen::Matrix<double,-1,1,0,3,1> evaluate_basis(Eigen::Vector3d const &x, size_t k, size_t d, size_t i, Eigen::VectorXd const &b) const;
53 double get_scaling(size_t d, size_t i) const;
54
55 private:
56 template<size_t d,size_t l> void _fill_lists();
57 inline int _cmp_ind(size_t l, size_t d) const {return _dim_table[d-1]+l;}
58
59 int _r;
60 std::array<size_t,4> _nbelem;
61 std::vector<Cell_basis<3,1>> _1cells;
62 std::vector<Cell_basis<3,2>> _2cells;
63 std::vector<Cell_basis<3,3>> _3cells;
64 std::vector<size_t> _dim_table;
65 std::vector<Eigen::MatrixXd> _list_diff, // Image of dPL(r,l,d) inside PL(r-1,l+1,d)
66 _list_Koszul, // Image of kPL(r,l,d) inside PL(r+1,l-1,d)
67 _list_diff_as_degr, // Image of dPL(r,l,d) inside PL(r,l+1,d)
68 _list_trimmed, // Image of PLtrimmed(r,l,d)
69 _list_reduced_Koszul_m1; // Image of kPL(r-1,l,d)
70
71 };
72
73}
74#endif
75
Definition ddr_pec.hpp:15
Class to describe a mesh.
Definition MeshND.hpp:17
@ Matrix
Definition basis.hpp:67
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
Kronecker get_trace(size_t k, size_t d, size_t i, size_t j) const
Return the trace for the k-forms on the i-th d-cell onto its j-th (d-1)-neighbour.
Definition ddr_pec.cpp:244
const Eigen::MatrixXd & get_trimmed(size_t l, size_t d) const
Definition ddr_pec.hpp:31
Eigen::Matrix< double,-1,-1, 0, 3, 3 > get_hodge_star(size_t k, size_t d, size_t i_cell) const
Gets the continuous Hodge star on k forms on the mesh element of dimension d with index i.
Definition ddr_pec.cpp:300
Eigen::Matrix< double,-1,-1, 0, 3, 3 > get_exterior_trace(size_t k, size_t d, size_t i_cell) const
Trace of a k form from the global space to the mesh element of dimension d with index i_cell.
Definition ddr_pec.cpp:274
Eigen::MatrixXd discrete_hodge_star(size_t k, size_t d, size_t i_cell) const
Takes a k form on the polynomial k basis to a d-k form on the d-k basis (deals with scaling) on the m...
Definition ddr_pec.cpp:316
Eigen::MatrixXd discrete_inv_hodge_star(size_t k, size_t d, size_t i_cell) const
e.g. discrete_inv_hodge_star(d-k,d,i) * potential(k,d,i) is the potential P^k in PL(r,...
Definition ddr_pec.cpp:333
double get_scaling(size_t d, size_t i) const
Gets scaling of the i-th d-cell.
Definition ddr_pec.cpp:370
Eigen::KroneckerProduct< Eigen::MatrixXd, Eigen::MatrixXd > Kronecker
Definition ddr_pec.hpp:21
const Eigen::MatrixXd & get_Koszul(size_t l, size_t d) const
Definition ddr_pec.hpp:29
const Eigen::MatrixXd & get_reduced_Koszul_m1(size_t l, size_t d) const
Definition ddr_pec.hpp:32
double evaluate_scalar_basis(Eigen::Vector3d const &x, size_t d, size_t i_cell, int i_basis) const
Definition ddr_pec.cpp:257
Kronecker get_mass(size_t k, size_t d, size_t i) const
Return the mass matrix for the k-forms on the i-th d-cell.
Definition ddr_pec.cpp:231
const Eigen::MatrixXd & get_diff(size_t l, size_t d) const
Getters for the generic operators matrices.
Definition ddr_pec.hpp:28
const Eigen::MatrixXd & get_diff_as_degr(size_t l, size_t d) const
Definition ddr_pec.hpp:30
Eigen::Matrix< double,-1, 1, 0, 3, 1 > evaluate_basis(Eigen::Vector3d const &x, size_t k, size_t d, size_t i_cell, int i_basis) const
Evaluate at point x, the i th basis function of PL(r,k,d) in the cell with index i_cell.
Definition ddr_pec.cpp:337
Definition ddr-magnetostatics.hpp:41