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
laxcurl.hpp
Go to the documentation of this file.
1#ifndef LAXCURL_HPP
2#define LAXCURL_HPP
3
4#include <globaldofspace.hpp>
5#include <xcurl.hpp>
6#include <liealgebra.hpp>
7
8namespace HArDCore3D
9{
16
17 class LAXCurl : public GlobalDOFSpace
18 {
19 public:
20 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
21 typedef std::vector<FunctionType> LAFunctionType;
22
24 inline const LieAlgebra & lieAlgebra() const
25 {
26 return m_lie_algebra;
27 }
29 inline const XCurl & xCurl() const
30 {
31 return m_xcurl;
32 }
33
36 {
38 const Eigen::MatrixXd & _curl,
39 const Eigen::MatrixXd & _potential
40 )
41 : curl(_curl),
43 {
44 // Do nothing
45 }
46 Eigen::MatrixXd curl;
47 Eigen::MatrixXd potential;
48 };
49
51 LAXCurl(const LieAlgebra & lie_algebra, const XCurl & xcurl, bool use_threads = true, std::ostream & output = std::cout);
52
54 const Mesh & mesh() const
55 {
56 return m_xcurl.mesh();
57 }
58
60 const size_t & degree() const
61 {
62 return m_xcurl.degree();
63 }
64
66 Eigen::VectorXd interpolate(
67 const LAFunctionType & v,
68 const int doe_cell = -1,
69 const int doe_face = -1,
70 const int doe_edge = -1
71 ) const;
72
74 inline const LocalOperators & cellOperators(size_t iT) const
75 {
76 return *m_cell_operators[iT];
77 }
78
80 inline const LocalOperators & cellOperators(const Cell & T) const
81 {
82 return * m_cell_operators[T.global_index()];
83 }
84
86 inline const LocalOperators & faceOperators(size_t iF) const
87 {
88 return *m_face_operators[iF];
89 }
90
92 inline const LocalOperators & faceOperators(const Face & F) const
93 {
94 return * m_face_operators[F.global_index()];
95 }
96
98 inline const DDRCore::CellBases & cellBases(size_t iT) const
99 {
100 return m_xcurl.cellBases(iT);
101 }
102
104 inline const DDRCore::CellBases & cellBases(const Cell & T) const
105 {
106 return m_xcurl.cellBases(T.global_index());
107 }
108
110 inline const DDRCore::FaceBases & faceBases(size_t iF) const
111 {
112 return m_xcurl.faceBases(iF);
113 }
114
116 inline const DDRCore::FaceBases & faceBases(const Face & F) const
117 {
118 return m_xcurl.faceBases(F.global_index());
119 }
120
122 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
123 {
124 return m_xcurl.edgeBases(iE);
125 }
126
128 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
129 {
130 return m_xcurl.edgeBases(E.global_index());
131 }
132
136 Eigen::MatrixXd computeL2Product(
137 const size_t iT,
138 const double & penalty_factor = 1.,
139 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
141 ) const;
142
144 Eigen::MatrixXd computeL2ProductGradient(
145 const size_t iT,
146 const XGrad & x_grad,
147 const std::string & side,
148 const double & penalty_factor = 1.,
149 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
151 ) const;
152
153 private:
154 LocalOperators _compute_face_curl_potential(size_t iF);
155 LocalOperators _compute_cell_curl_potential(size_t iT);
156
157 const XCurl & m_xcurl;
158 const LieAlgebra & m_lie_algebra;
159 bool m_use_threads;
160 std::ostream & m_output;
161
162 // Containers for local operators
163 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
164 std::vector<std::unique_ptr<LocalOperators> > m_face_operators;
165
166 };
167
168}// end of namespace HArDCore3D
169
170#endif
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Discrete Lie algebra valued Hcurl space.
Definition laxcurl.hpp:18
Lie algebra class: mass matrix, structure constants and Lie bracket.
Definition liealgebra.hpp:17
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(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
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xcurl.hpp:112
const size_t & degree() const
Return the polynomial degree.
Definition xcurl.hpp:50
const Mesh & mesh() const
Return the mesh.
Definition xcurl.hpp:44
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition laxcurl.hpp:80
LocalOperators(const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential)
Definition laxcurl.hpp:37
Eigen::MatrixXd curl
Definition laxcurl.hpp:46
const size_t & degree() const
Return the polynomial degree.
Definition laxcurl.hpp:60
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 laxcurl.cpp:124
const LieAlgebra & lieAlgebra() const
Returns the Lie algebra.
Definition laxcurl.hpp:24
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition laxcurl.hpp:98
std::vector< FunctionType > LAFunctionType
Definition laxcurl.hpp:21
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition laxcurl.hpp:86
Eigen::MatrixXd potential
Definition laxcurl.hpp:47
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition laxcurl.hpp:104
const XCurl & xCurl() const
Returns the space XCurl.
Definition laxcurl.hpp:29
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition laxcurl.hpp:74
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition laxcurl.hpp:122
const Mesh & mesh() const
Return the mesh.
Definition laxcurl.hpp:54
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
Definition laxcurl.cpp:110
const LocalOperators & faceOperators(const Face &F) const
Return face operators for face F.
Definition laxcurl.hpp:92
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition laxcurl.hpp:128
Eigen::VectorXd interpolate(const LAFunctionType &v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
Interpolator of a continuous Lie algebra valued function decomposed on the basis of the LieAlgebra,...
Definition laxcurl.cpp:63
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition laxcurl.hpp:20
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition laxcurl.hpp:110
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition laxcurl.hpp:116
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 laxcurl.hpp:36