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
lasxcurl.hpp
Go to the documentation of this file.
1#ifndef LASXCURL_HPP
2#define LASXCURL_HPP
3
4#include <sxcurl.hpp>
5#include <liealgebra.hpp>
6#include <unsupported/Eigen/KroneckerProduct>
7
8namespace HArDCore3D
9{
16
18 {
19 public:
20 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
21 typedef std::vector<FunctionType> LAFunctionType;
22
24
26 {
28 const Eigen::MatrixXd & _serendipity,
29 const Eigen::MatrixXd & _extension,
30 const Eigen::MatrixXd & _reduction
31 )
35 {
36 // Do nothing
37 }
38
39 Eigen::MatrixXd serendipity;
40 Eigen::MatrixXd extension;
41 Eigen::MatrixXd reduction;
42 };
43
45 LASXCurl(const LieAlgebra & lie_algebra, const SXCurl & sx_curl, bool use_threads = true, std::ostream & output = std::cout);
46
48 const Mesh & mesh() const
49 {
50 return m_sxcurl.mesh();
51 }
52
54 const size_t & degree() const
55 {
56 return m_sxcurl.degree();
57 }
58
60 inline const LieAlgebra & lieAlg() const
61 {
62 return m_lie_algebra;
63 }
64
66 const SerendipityProblem & serPro() const
67 {
68 return m_sxcurl.serPro();
69 }
70
72 Eigen::VectorXd interpolate(
73 const LAFunctionType & v,
74 const int doe_cell = -1,
75 const int doe_face = -1,
76 const int doe_edge = -1
77 ) const;
78
79 //---------------------------------//
80 //---- Lie algebra operators-------//
81
82 inline const Eigen::VectorXd LaxcurlI(size_t iG, Eigen::VectorXd & v) const
83 {
84 return Eigen::Map<Eigen::VectorXd,0,Eigen::InnerStride<>>(v.data()+iG, m_sxcurl.dimension(),Eigen::InnerStride<>(m_lie_algebra.dimension()));
85 }
86
87 //---------------------------------//
88 //---- Face transfer operators-----//
89
91 inline const Eigen::MatrixXd & ScurlFace(size_t iF) const
92 {
93 return (*m_face_transfer_operators[iF]).serendipity;
94 }
95
97 inline const Eigen::MatrixXd & EcurlFace(size_t iF) const
98 {
99 return (*m_face_transfer_operators[iF]).extension;
100 }
101
103 inline const Eigen::MatrixXd & RcurlFace(size_t iF) const
104 {
105 return (*m_face_transfer_operators[iF]).reduction;
106 }
107
109 inline const Eigen::MatrixXd & ScurlFace(const Face & F) const
110 {
111 return ScurlFace(F.global_index());
112 }
113
115 inline const Eigen::MatrixXd & EcurlFace(const Face & F) const
116 {
117 return EcurlFace(F.global_index());
118 }
119
121 inline const Eigen::MatrixXd & RcurlFace(const Face & F) const
122 {
123 return RcurlFace(F.global_index());
124 }
125
126 //----------------------------------//
127 //---- Cell transfer operators -----//
129 inline const Eigen::MatrixXd & ScurlCell(size_t iT) const
130 {
131 return (*m_cell_transfer_operators[iT]).serendipity;
132 }
133
135 inline const Eigen::MatrixXd & EcurlCell(size_t iT) const
136 {
137 return (*m_cell_transfer_operators[iT]).extension;
138 }
139
141 inline const Eigen::MatrixXd & RcurlCell(size_t iT) const
142 {
143 return (*m_cell_transfer_operators[iT]).reduction;
144 }
145
147 inline const Eigen::MatrixXd & ScurlCell(const Cell & T) const
148 {
149 return ScurlCell(T.global_index());
150 }
151
153 inline const Eigen::MatrixXd & EcurlCell(const Cell & T) const
154 {
155 return EcurlCell(T.global_index());
156 }
157
159 inline const Eigen::MatrixXd & RcurlCell(const Cell & T) const
160 {
161 return RcurlCell(T.global_index());
162 }
163
164 //-----------------------------------------------------------------//
165 //---- Full curl and potential reconstructions, and L2 product ----//
166
168 inline const Eigen::MatrixXd faceCurl(size_t iF) const
169 {
170 Eigen::MatrixXd id = Eigen::MatrixXd::Identity(m_lie_algebra.dimension(), m_lie_algebra.dimension());
171 return Eigen::KroneckerProduct(m_sxcurl.faceCurl(iF), id);
172 }
173
175 inline const Eigen::MatrixXd faceCurl(const Face & F) const
176 {
177 return faceCurl(F.global_index());
178 }
179
181 inline const Eigen::MatrixXd facePotential(size_t iF) const
182 {
183 Eigen::MatrixXd id = Eigen::MatrixXd::Identity(m_lie_algebra.dimension(), m_lie_algebra.dimension());
184 return Eigen::KroneckerProduct(m_sxcurl.facePotential(iF), id);
185 }
186
188 inline const Eigen::MatrixXd facePotential(const Face & F) const
189 {
190 return facePotential(F.global_index());
191 }
192
194 inline const Eigen::MatrixXd cellCurl(size_t iT) const
195 {
196 Eigen::MatrixXd id = Eigen::MatrixXd::Identity(m_lie_algebra.dimension(), m_lie_algebra.dimension());
197 return Eigen::KroneckerProduct(m_sxcurl.cellCurl(iT), id);
198 }
199
201 inline const Eigen::MatrixXd cellCurl(const Cell & T) const
202 {
203 return cellCurl(T.global_index());
204 }
205
207 inline const Eigen::MatrixXd cellPotential(size_t iT) const
208 {
209 Eigen::MatrixXd id = Eigen::MatrixXd::Identity(m_lie_algebra.dimension(), m_lie_algebra.dimension());
210 return Eigen::KroneckerProduct(m_sxcurl.cellPotential(iT), id);
211 }
212
214 inline const Eigen::MatrixXd cellPotential(const Cell & T) const
215 {
216 return cellPotential(T.global_index());
217 }
218
220 Eigen::MatrixXd computeL2Product(
221 const size_t iT,
222 const double & penalty_factor = 1.,
223 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
225 ) const
226 {
227 return Eigen::KroneckerProduct(m_sxcurl.computeL2Product(iT, penalty_factor, mass_Pk3_T, weight), m_lie_algebra.massMatrix());
228 }
229
230
232 Eigen::MatrixXd computeL2ProductGradient(
233 const size_t iT,
234 const SXGrad & sx_grad,
235 const std::string & side,
236 const double & penalty_factor = 1.,
237 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
239 ) const;
240
241 //-----------------------//
242 //---- Getters ----------//
243
245 inline const DDRCore::CellBases & cellBases(size_t iT) const
246 {
247 return m_sxcurl.cellBases(iT);
248 }
249
251 inline const DDRCore::CellBases & cellBases(const Cell & T) const
252 {
253 return m_sxcurl.cellBases(T.global_index());
254 }
255
257 inline const DDRCore::FaceBases & faceBases(size_t iF) const
258 {
259 return m_sxcurl.faceBases(iF);
260 }
261
263 inline const DDRCore::FaceBases & faceBases(const Face & F) const
264 {
265 return m_sxcurl.faceBases(F.global_index());
266 }
267
269 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
270 {
271 return m_sxcurl.edgeBases(iE);
272 }
273
275 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
276 {
277 return m_sxcurl.edgeBases(E.global_index());
278 }
279
280 private:
281 TransferOperators _compute_face_transfer_operators(size_t iF);
282 TransferOperators _compute_cell_transfer_operators(size_t iT);
283
284 const LieAlgebra & m_lie_algebra;
285 const SXCurl & m_sxcurl;
286
287 // Containers for serendipity, extension and reduction operators
288 std::vector<std::unique_ptr<TransferOperators> > m_face_transfer_operators;
289 std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
290
291 bool m_use_threads;
292 std::ostream & m_output;
293
294 };
295
296} // end of namespace HArDCore3D
297#endif
Discrete Lie algebra valued Serendipity Hcurl space: local operators, L2 product and global interpola...
Definition lasxcurl.hpp:18
Lie algebra class: mass matrix, structure constants and Lie bracket.
Definition liealgebra.hpp:17
Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator.
Definition sxcurl.hpp:21
Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator.
Definition sxgrad.hpp:20
Construct all polynomial spaces for the DDR sequence.
Definition serendipity_problem.hpp:40
Base class for global DOF spaces.
Definition variabledofspace.hpp:17
Class to describe a mesh.
Definition MeshND.hpp:17
@ Matrix
Definition basis.hpp:67
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition variabledofspace.hpp:158
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition sxcurl.hpp:179
const Mesh & mesh() const
Return the mesh.
Definition sxcurl.hpp:50
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition sxcurl.hpp:191
const Eigen::MatrixXd faceCurl(size_t iF) const
Return the full curl operator on the face of index iF.
Definition sxcurl.hpp:155
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.
Definition sxcurl.hpp:203
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition sxcurl.hpp:235
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition sxcurl.hpp:259
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition sxcurl.hpp:247
const SerendipityProblem & serPro() const
Definition sxcurl.hpp:61
const size_t & degree() const
Return the polynomial degree.
Definition sxcurl.hpp:56
const Eigen::MatrixXd facePotential(size_t iF) const
Return the potential operator on the face of index iF.
Definition sxcurl.hpp:167
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition lasxcurl.hpp:269
Eigen::MatrixXd extension
Definition lasxcurl.hpp:40
const Eigen::MatrixXd & ScurlCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition lasxcurl.hpp:147
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition lasxcurl.hpp:194
const Eigen::MatrixXd facePotential(const Face &F) const
Return the potential operator on face F.
Definition lasxcurl.hpp:188
const Eigen::MatrixXd cellPotential(const Cell &T) const
Return the potential operator on cell T.
Definition lasxcurl.hpp:214
std::vector< FunctionType > LAFunctionType
Definition lasxcurl.hpp:21
const Mesh & mesh() const
Return the mesh.
Definition lasxcurl.hpp:48
const Eigen::MatrixXd & EcurlCell(size_t iT) const
Return the extension for the cell of index iT.
Definition lasxcurl.hpp:135
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.
Definition lasxcurl.hpp:220
const Eigen::MatrixXd & RcurlFace(const Face &F) const
Return cell reduction for cell T.
Definition lasxcurl.hpp:121
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition lasxcurl.hpp:275
const Eigen::MatrixXd & RcurlCell(size_t iT) const
Return the reduction for the cell of index iT.
Definition lasxcurl.hpp:141
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition lasxcurl.hpp:245
const Eigen::MatrixXd cellCurl(const Cell &T) const
Return the full curl operator on cell T.
Definition lasxcurl.hpp:201
const Eigen::MatrixXd & ScurlFace(const Face &F) const
Return the serendipity reconstruction for face F.
Definition lasxcurl.hpp:109
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition lasxcurl.hpp:263
const Eigen::MatrixXd facePotential(size_t iF) const
Return the potential operator on the face of index iF.
Definition lasxcurl.hpp:181
const Eigen::MatrixXd & ScurlFace(size_t iF) const
Return the serendipity reconstruction for the face of index iF.
Definition lasxcurl.hpp:91
const Eigen::MatrixXd & RcurlCell(const Cell &T) const
Return the reduction for cell T.
Definition lasxcurl.hpp:159
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition lasxcurl.hpp:251
const Eigen::MatrixXd & EcurlCell(const Cell &T) const
Return the extension for cell T.
Definition lasxcurl.hpp:153
Eigen::MatrixXd serendipity
Definition lasxcurl.hpp:39
const Eigen::VectorXd LaxcurlI(size_t iG, Eigen::VectorXd &v) const
Definition lasxcurl.hpp:82
const SerendipityProblem & serPro() const
Return the serendipity operators.
Definition lasxcurl.hpp:66
const Eigen::MatrixXd & massMatrix() const
Returns the Gram matrix of the Lie algebra.
Definition liealgebra.hpp:40
const Eigen::MatrixXd faceCurl(size_t iF) const
Return the full curl operator on the face of index iF.
Definition lasxcurl.hpp:168
const size_t dimension() const
Assuming orthonormal basis.
Definition liealgebra.hpp:28
TransferOperators(const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
Definition lasxcurl.hpp:27
const Eigen::MatrixXd & EcurlFace(const Face &F) const
Return the extension for face F.
Definition lasxcurl.hpp:115
Eigen::MatrixXd reduction
Definition lasxcurl.hpp:41
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition lasxcurl.hpp:257
const Eigen::MatrixXd & ScurlCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition lasxcurl.hpp:129
const LieAlgebra & lieAlg() const
Return the Lie algebra.
Definition lasxcurl.hpp:60
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition lasxcurl.hpp:207
const Eigen::MatrixXd & RcurlFace(size_t iF) const
Return the reduction for the face of index iF.
Definition lasxcurl.hpp:103
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const SXGrad &sx_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 lasxcurl.cpp:100
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition lasxcurl.hpp:20
const Eigen::MatrixXd faceCurl(const Face &F) const
Return the full curl operator on face F.
Definition lasxcurl.hpp:175
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 function.
Definition lasxcurl.cpp:61
const Eigen::MatrixXd & EcurlFace(size_t iF) const
Return the extension for the face of index iF.
Definition lasxcurl.hpp:97
const size_t & degree() const
Return the polynomial degree.
Definition lasxcurl.hpp:54
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 serendipity, extension and reduction operators.
Definition lasxcurl.hpp:26