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
xcurl.hpp
Go to the documentation of this file.
1#ifndef XCURL_HPP
2#define XCURL_HPP
3
4#include <globaldofspace.hpp>
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
7#include <xgrad.hpp>
8
9namespace HArDCore3D
10{
17
18 class XCurl : public GlobalDOFSpace
19 {
20 public:
21 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
22
25 {
27 const Eigen::MatrixXd & _curl,
28 const Eigen::MatrixXd & _potential
29 )
30 : curl(_curl),
32 {
33 // Do nothing
34 }
35
36 Eigen::MatrixXd curl;
37 Eigen::MatrixXd potential;
38 };
39
41 XCurl(const DDRCore & ddr_core, bool use_threads = true, std::ostream & output = std::cout);
42
44 const Mesh & mesh() const
45 {
46 return m_ddr_core.mesh();
47 }
48
50 const size_t & degree() const
51 {
52 return m_ddr_core.degree();
53 }
54
56 Eigen::VectorXd interpolate(
57 const FunctionType & v,
58 const int doe_cell = -1,
59 const int doe_face = -1,
60 const int doe_edge = -1
61 ) const;
62
64 inline const LocalOperators & cellOperators(size_t iT) const
65 {
66 return *m_cell_operators[iT];
67 }
68
70 inline const LocalOperators & cellOperators(const Cell & T) const
71 {
72 return * m_cell_operators[T.global_index()];
73 }
74
76 inline const LocalOperators & faceOperators(size_t iF) const
77 {
78 return *m_face_operators[iF];
79 }
80
82 inline const LocalOperators & faceOperators(const Face & F) const
83 {
84 return * m_face_operators[F.global_index()];
85 }
86
88 inline const DDRCore::CellBases & cellBases(size_t iT) const
89 {
90 return m_ddr_core.cellBases(iT);
91 }
92
94 inline const DDRCore::CellBases & cellBases(const Cell & T) const
95 {
96 return m_ddr_core.cellBases(T.global_index());
97 }
98
100 inline const DDRCore::FaceBases & faceBases(size_t iF) const
101 {
102 return m_ddr_core.faceBases(iF);
103 }
104
106 inline const DDRCore::FaceBases & faceBases(const Face & F) const
107 {
108 return m_ddr_core.faceBases(F.global_index());
109 }
110
112 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
113 {
114 return m_ddr_core.edgeBases(iE);
115 }
116
118 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
119 {
120 return m_ddr_core.edgeBases(E.global_index());
121 }
122
124 // The mass matrix of P^k(T)^3 is the most expensive mass matrix in the calculation of this norm, which
125 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
126 Eigen::MatrixXd computeL2Product(
127 const size_t iT,
128 const double & penalty_factor = 1.,
129 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
131 ) const;
132
134 Eigen::MatrixXd computeL2ProductGradient(
135 const size_t iT,
136 const XGrad & x_grad,
137 const std::string & side,
138 const double & penalty_factor = 1.,
139 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
141 ) const;
142
144 Eigen::MatrixXd computeL2Product_with_Ops(
145 const size_t iT,
146 const std::vector<Eigen::MatrixXd> & leftOp,
147 const std::vector<Eigen::MatrixXd> & rightOp,
148 const double & penalty_factor,
149 const Eigen::MatrixXd & w_mass_Pk3_T,
150 const IntegralWeight & weight
151 ) const;
152
153
155 Eigen::MatrixXd computeL2ProductDOFs(
156 size_t iT,
157 const double & penalty_factor = 1.,
158 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
160 );
161
162 private:
163 LocalOperators _compute_face_curl_potential(size_t iF);
164 LocalOperators _compute_cell_curl_potential(size_t iT);
165
166 const DDRCore & m_ddr_core;
167 bool m_use_threads;
168 std::ostream & m_output;
169
170 // Containers for local operators
171 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
172 std::vector<std::unique_ptr<LocalOperators> > m_face_operators;
173
174 };
175
176} // end of namespace HArDCore3D
177#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 Hcurl space: local operators, L2 product and global interpolator.
Definition xcurl.hpp:19
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::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xcurl.hpp:94
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition xcurl.hpp:76
Eigen::MatrixXd computeL2Product_with_Ops(const size_t iT, const std::vector< Eigen::MatrixXd > &leftOp, const std::vector< Eigen::MatrixXd > &rightOp, const double &penalty_factor, const Eigen::MatrixXd &w_mass_Pk3_T, const IntegralWeight &weight) const
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables....
Definition xcurl.cpp:470
Eigen::MatrixXd computeL2ProductDOFs(size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.))
Alternate version of L2 product, using only cell potentials and 'DOFs' (as in Di-Pietro & Droniou,...
Definition xcurl.cpp:589
Eigen::VectorXd interpolate(const FunctionType &v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
Interpolator of a continuous function.
Definition xcurl.cpp:63
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition xcurl.hpp:88
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition xcurl.hpp:100
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const XGrad &x_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
Compute the matrix of the (weighted) L2-product as 'computeL2Product', with application of the discre...
Definition xcurl.cpp:395
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xcurl.hpp:112
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xcurl.hpp:70
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:134
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition xcurl.hpp:106
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:162
const LocalOperators & faceOperators(const Face &F) const
Return face operators for face F.
Definition xcurl.hpp:82
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition xcurl.hpp:21
const size_t & degree() const
Return the polynomial degree.
Definition xcurl.hpp:50
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xcurl.hpp:64
Eigen::MatrixXd potential
Definition xcurl.hpp:37
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:140
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xcurl.hpp:118
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk3_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 xcurl.cpp:346
Eigen::MatrixXd curl
Definition xcurl.hpp:36
const Mesh & mesh() const
Return the mesh.
Definition xcurl.hpp:44
LocalOperators(const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential)
Definition xcurl.hpp:26
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition ddrcore.hpp:154
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 the local operators (curl and potential)
Definition xcurl.hpp:25