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
xgrad.hpp
Go to the documentation of this file.
1#ifndef XGRAD_HPP
2#define XGRAD_HPP
3
4#include <globaldofspace.hpp>
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
7
8namespace HArDCore3D
9{
16
17 class XGrad : public GlobalDOFSpace
18 {
19 public:
20 typedef std::function<double(const Eigen::Vector3d &)> FunctionType;
21
24 {
26 const Eigen::MatrixXd & _gradient,
27 const Eigen::MatrixXd & _potential
28 )
31 {
32 // Do nothing
33 }
34
35 Eigen::MatrixXd gradient;
36 Eigen::MatrixXd potential;
37 };
38
40 XGrad(const DDRCore & ddr_core, bool use_threads = true, std::ostream & output = std::cout);
41
43 const Mesh & mesh() const
44 {
45 return m_ddr_core.mesh();
46 }
47
49 const size_t & degree() const
50 {
51 return m_ddr_core.degree();
52 }
53
55 Eigen::VectorXd interpolate(
56 const FunctionType & q,
57 const int doe_cell = -1,
58 const int doe_face = -1,
59 const int doe_edge = -1
60 ) const;
61
63 inline const LocalOperators & edgeOperators(size_t iE) const
64 {
65 return *m_edge_operators[iE];
66 }
67
69 inline const LocalOperators & edgeOperators(const Edge & E) const
70 {
71 return *m_edge_operators[E.global_index()];
72 }
73
75 inline const LocalOperators & faceOperators(size_t iF) const
76 {
77 return *m_face_operators[iF];
78 }
79
81 inline const LocalOperators & faceOperators(const Face & F) const
82 {
83 return *m_face_operators[F.global_index()];
84 }
85
87 inline const LocalOperators & cellOperators(size_t iT) const
88 {
89 return *m_cell_operators[iT];
90 }
91
93 inline const LocalOperators & cellOperators(const Cell & T) const
94 {
95 return *m_cell_operators[T.global_index()];
96 }
97
99 inline const DDRCore::CellBases & cellBases(size_t iT) const
100 {
101 return m_ddr_core.cellBases(iT);
102 }
103
105 inline const DDRCore::CellBases & cellBases(const Cell & T) const
106 {
107 return m_ddr_core.cellBases(T.global_index());
108 }
109
111 inline const DDRCore::FaceBases & faceBases(size_t iF) const
112 {
113 return m_ddr_core.faceBases(iF);
114 }
115
117 inline const DDRCore::FaceBases & faceBases(const Face & F) const
118 {
119 return m_ddr_core.faceBases(F.global_index());
120 }
121
123 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
124 {
125 return m_ddr_core.edgeBases(iE);
126 }
127
129 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
130 {
131 return m_ddr_core.edgeBases(E.global_index());
132 }
133
135 // This essentially calls computeStabilisation and adds the consistent part (weighted mass matrix).
136 Eigen::MatrixXd computeL2Product(
137 const size_t iT,
138 const double & penalty_factor = 1.,
139 const Eigen::MatrixXd & mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1),
141 ) const;
142
144 Eigen::MatrixXd computeStabilisation(
145 const size_t iT,
147 ) const;
148
150 Eigen::MatrixXd buildGradientComponentsCell(size_t iT) const;
151
152 private:
153 LocalOperators _compute_edge_gradient_potential(size_t iE);
154 LocalOperators _compute_face_gradient_potential(size_t iF);
155 LocalOperators _compute_cell_gradient_potential(size_t iT);
156
157 const DDRCore & m_ddr_core;
158 bool m_use_threads;
159 std::ostream & m_output;
160
161 // Containers for local operators
162 std::vector<std::unique_ptr<LocalOperators> > m_edge_operators;
163 std::vector<std::unique_ptr<LocalOperators> > m_face_operators;
164 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
165 };
166
167} // end of namespace HArDCore3D
168
169#endif
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:62
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Discrete H1 space: local operators, L2 product and global interpolator.
Definition xgrad.hpp:18
Class to describe a mesh.
Definition MeshND.hpp:17
@ Matrix
Definition basis.hpp:67
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition xgrad.hpp:117
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xgrad.hpp:123
LocalOperators(const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential)
Definition xgrad.hpp:25
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xgrad.hpp:99
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xgrad.hpp:105
std::function< double(const Eigen::Vector3d &)> FunctionType
Definition xgrad.hpp:20
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:134
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:162
const Mesh & mesh() const
Return the mesh.
Definition xgrad.hpp:43
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xgrad.hpp:129
const LocalOperators & edgeOperators(size_t iE) const
Return edge operators for the edge of index iE.
Definition xgrad.hpp:63
Eigen::MatrixXd potential
Definition xgrad.hpp:36
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pkpo_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
Compute the matrix of the (weighted) L2-product for the cell of index iT. The stabilisation here is b...
Definition xgrad.cpp:401
Eigen::MatrixXd buildGradientComponentsCell(size_t iT) const
Build the components of the gradient operator (probably not useful in practice to implement schemes)
Definition xgrad.cpp:540
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:140
Eigen::MatrixXd gradient
Definition xgrad.hpp:35
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition xgrad.hpp:75
const LocalOperators & faceOperators(const Face &F) const
Return face operators for face F.
Definition xgrad.hpp:81
Eigen::MatrixXd computeStabilisation(const size_t iT, const IntegralWeight &weight=IntegralWeight(1.)) const
Computes only the stabilisation matrix of the (weighted) L2-product for the cell of index iT....
Definition xgrad.cpp:441
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition xgrad.hpp:111
const size_t & degree() const
Return the polynomial degree.
Definition xgrad.hpp:49
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xgrad.hpp:87
const LocalOperators & edgeOperators(const Edge &E) const
Return edge operators for edge E.
Definition xgrad.hpp:69
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xgrad.hpp:93
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition ddrcore.hpp:154
Eigen::VectorXd interpolate(const FunctionType &q, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
Interpolator of a continuous function.
Definition xgrad.cpp:77
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:146
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
Definition ddr-magnetostatics.hpp:41
Structure to store element bases.
Definition ddrcore.hpp:86
Structure to store edge bases.
Definition ddrcore.hpp:121
Structure to store face bases.
Definition ddrcore.hpp:105
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:36
A structure to store local operators (gradient and potential)
Definition xgrad.hpp:24