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_spaces.hpp
Go to the documentation of this file.
1#ifndef DDR_SPACES_HPP
2#define DDR_SPACES_HPP
3
4#include "mesh.hpp"
5#include "globaldofspace.hpp"
6#include "ddr_pec.hpp"
7
8#include <memory>
9#include <variant>
10
11namespace HArDCore3D {
12
19 {
20 public:
21
23 DDR_Spaces(Mesh const & mesh, int r, bool use_threads = true, std::ostream & output = std::cout);
24
26 typedef std::function<double(const Eigen::Vector3d &)> scalar;
27 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> vector;
28 std::variant<scalar,vector> func;
29 int dqr[4] = {-1,-1,-1,-1}; // Quadrature degree on {vertex (unused),edge,face,cell}. Default to 2*r + 3 if negative
30 };
31
33 const Mesh & mesh() const
34 {
35 return _mesh;
36 }
37
39 inline int degree() const {return _r;}
40
42 Eigen::VectorXd interpolate(const DDR_function_type &, size_t k) const;
43
45 const Eigen::MatrixXd & full_diff(size_t k, size_t d, size_t i_cell) const;
47 Eigen::MatrixXd compose_diff(size_t k, size_t d, size_t i_cell) const;
49 const Eigen::MatrixXd & potential(size_t k, size_t d, size_t i_cell) const;
50
52 inline GlobalDOFSpace const & dofspace(size_t k) const {
53 return _dofspace[k];
54 }
55
57 Eigen::Matrix<double,-1,1,0,3,1> evaluate_basis(Eigen::Vector3d const &x, size_t k, size_t d, size_t i_cell, Eigen::VectorXd const &b) const;
58
60 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;
61
63 DDR_PEC::Kronecker get_mass(size_t k, size_t d, size_t i_cell) const;
64
66 Eigen::Matrix<double,-1,-1,0,3,3> get_hodge_star(size_t k, size_t d, size_t i_cell) const;
67
69 Eigen::MatrixXd discrete_hodge_star(size_t k, size_t d, size_t i_cell) const;
70
71 // 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
73 Eigen::MatrixXd discrete_inv_hodge_star(size_t k, size_t d, size_t i_cell) const;
74
76 Eigen::Matrix<double,-1,-1,0,3,3> get_exterior_trace(size_t k, size_t d, size_t i_cell) const;
77
79 Eigen::MatrixXd computeL2Product(size_t k, size_t i_cell, const double penalty_factor = 0.1) const;
80
82 std::vector<VectorRd> computeVertexValues(size_t k, const Eigen::VectorXd & u) const;
83
84 private:
85 // More consistent notation to access, but generate many useless uninitialized entries
86 struct DDR_Operators { // one for each form degree k
87 std::array<Eigen::MatrixXd,4> full_diff; // equal to \star d, PL(r,d-k-1,d) valued
88 std::array<Eigen::MatrixXd,4> diff; // equal to \star \ul{d}, PLtrimmed(r,d-k-1,d) valued
89 std::array<Eigen::MatrixXd,4> P; // equal to \star P, PL(r,d-k,d) valued
90 };
91
92 // Pointer to the mesh
93 const Mesh & _mesh;
94
95 int _r;
96 bool _use_threads;
97 std::unique_ptr<DDR_PEC> _ddr, _ddr_po;
98 std::array<GlobalDOFSpace,3+1> _dofspace; // one for all form degree 0 <= k <= n=3 dimensions
99 std::array<std::vector<DDR_Operators>,4> _ops; // one for each cell dimension 0 <= d <= n=3 dimensions
100 };
101}
102#endif
103
Definition ddr_spaces.hpp:19
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Class to describe a mesh.
Definition MeshND.hpp:17
@ Matrix
Definition basis.hpp:67
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
std::variant< scalar, vector > func
Definition ddr_spaces.hpp:28
int dqr[4]
Definition ddr_spaces.hpp:29
Eigen::MatrixXd computeL2Product(size_t k, size_t i_cell, const double penalty_factor=0.1) const
Computes the L2 product for discrete k-forms on the highest dimensional cell with index i_cell.
Definition ddr_spaces.cpp:328
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_spaces.cpp:318
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_spaces.cpp:303
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_cell) * potential(k, d, i_cell) is the potential P^k in PL(r,...
Definition ddr_spaces.cpp:313
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_spaces.cpp:308
std::array< Eigen::MatrixXd, 4 > diff
Definition ddr_spaces.hpp:88
const Eigen::MatrixXd & full_diff(size_t k, size_t d, size_t i_cell) const
Returns \star d in PL(r,d-k-1,d), the Hodge star of the full discrete exterior derivative of a discre...
Definition ddr_spaces.cpp:285
std::array< Eigen::MatrixXd, 4 > full_diff
Definition ddr_spaces.hpp:87
DDR_PEC::Kronecker get_mass(size_t k, size_t d, size_t i_cell) const
Returns the mass matrix for the k-forms on the mesh element of dimension d with index i_cell.
Definition ddr_spaces.cpp:323
GlobalDOFSpace const & dofspace(size_t k) const
Returns the dofspace associated to discrete k-forms.
Definition ddr_spaces.hpp:52
const Eigen::MatrixXd & potential(size_t k, size_t d, size_t i_cell) const
Returns \star P^k in PL(r,d-k,d), the Hodge star of the potential of a discrete k-form on the mesh el...
Definition ddr_spaces.cpp:290
Eigen::MatrixXd compose_diff(size_t k, size_t d, size_t i_cell) const
Returns \star \ul{d}^k in PLtrimmed(r,d-k-1,d) on the cell (including its boundary),...
Definition ddr_spaces.cpp:238
Eigen::KroneckerProduct< Eigen::MatrixXd, Eigen::MatrixXd > Kronecker
Definition ddr_pec.hpp:21
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> vector
Definition ddr_spaces.hpp:27
std::function< double(const Eigen::Vector3d &)> scalar
Definition ddr_spaces.hpp:26
std::array< Eigen::MatrixXd, 4 > P
Definition ddr_spaces.hpp:89
Eigen::Matrix< double,-1, 1, 0, 3, 1 > evaluate_basis(Eigen::Vector3d const &x, size_t k, size_t d, size_t i_cell, Eigen::VectorXd const &b) const
Evaluate at point x, the sum of the coefficients b times the basis of PL(r,k,d) in the cell with inde...
Definition ddr_spaces.cpp:295
std::vector< VectorRd > computeVertexValues(size_t k, const Eigen::VectorXd &u) const
Computes the averaged vertex values of a discrete k-form u using the potential.
Definition ddr_spaces.cpp:415
const Mesh & mesh() const
Return the mesh.
Definition ddr_spaces.hpp:33
int degree() const
Return the polynomial degree.
Definition ddr_spaces.hpp:39
Eigen::VectorXd interpolate(const DDR_function_type &, size_t k) const
Returns the interpolate of the given function as a k-form.
Definition ddr_spaces.cpp:177
Definition ddr-magnetostatics.hpp:41
Definition ddr_spaces.hpp:25