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
vcurl.hpp
Go to the documentation of this file.
1#ifndef VCURL_HPP
2#define VCURL_HPP
3
4#include <globaldofspace.hpp>
5#include <vemcore.hpp>
6#include <integralweight.hpp>
7
8namespace HArDCore3D
9{
16
27 class VCurl : public GlobalDOFSpace
28 {
29 public:
30 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
31
34 {
36 const Eigen::MatrixXd & _proj_curl,
37 const Eigen::MatrixXd & _dofs_curl,
38 const Eigen::MatrixXd & _proj_function
39 )
43 {
44 // Do nothing
45 }
46
47 Eigen::MatrixXd proj_curl;
48 Eigen::MatrixXd dofs_curl;
49 Eigen::MatrixXd proj_function;
50 };
51
53 VCurl(const VEMCore & vem_core, bool use_threads = true, std::ostream & output = std::cout);
54
56 const Mesh & mesh() const
57 {
58 return m_vem_core.mesh();
59 }
60
62 const size_t & degree() const
63 {
64 return m_vem_core.degree();
65 }
66
68 Eigen::VectorXd interpolate(
69 const FunctionType & v,
70 const FunctionType & curl_v,
71 const int doe_cell = -1,
72 const int doe_face = -1,
73 const int doe_edge = -1
74 ) const;
75
77 inline const LocalOperators & cellOperators(size_t iT) const
78 {
79 return *m_cell_operators[iT];
80 }
81
83 inline const LocalOperators & cellOperators(const Cell & T) const
84 {
85 return * m_cell_operators[T.global_index()];
86 }
87
89 inline const LocalOperators & faceOperators(size_t iF) const
90 {
91 return *m_face_operators[iF];
92 }
93
95 inline const LocalOperators & faceOperators(const Face & F) const
96 {
97 return * m_face_operators[F.global_index()];
98 }
99
101 inline const VEMCore::CellBases & cellBases(size_t iT) const
102 {
103 return m_vem_core.cellBases(iT);
104 }
105
107 inline const VEMCore::CellBases & cellBases(const Cell & T) const
108 {
109 return m_vem_core.cellBases(T.global_index());
110 }
111
113 inline const VEMCore::FaceBases & faceBases(size_t iF) const
114 {
115 return m_vem_core.faceBases(iF);
116 }
117
119 inline const VEMCore::FaceBases & faceBases(const Face & F) const
120 {
121 return m_vem_core.faceBases(F.global_index());
122 }
123
125 inline const VEMCore::EdgeBases & edgeBases(size_t iE) const
126 {
127 return m_vem_core.edgeBases(iE);
128 }
129
131 inline const VEMCore::EdgeBases & edgeBases(const Edge & E) const
132 {
133 return m_vem_core.edgeBases(E.global_index());
134 }
135
137 // The mass matrix of P^k(T)^3 is the most expensive mass matrix in the calculation of this norm, which
138 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
139 Eigen::MatrixXd computeL2Product(
140 const size_t iT,
141 const double & penalty_factor = 1.,
142 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
144 ) const;
145
146 private:
147 LocalOperators _compute_face_operators(size_t iF);
148 LocalOperators _compute_cell_operators(size_t iT);
149
150 const VEMCore & m_vem_core;
151 bool m_use_threads;
152 std::ostream & m_output;
153
154 // Containers for local operators
155 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
156 std::vector<std::unique_ptr<LocalOperators> > m_face_operators;
157
158 };
159
160} // end of namespace HArDCore3D
161#endif
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Virtual Hcurl space: local operators, L2 product and global interpolator.
Definition vcurl.hpp:28
Construct all polynomial spaces for the VEM sequence.
Definition vemcore.hpp:39
Class to describe a mesh.
Definition MeshND.hpp:17
@ Matrix
Definition basis.hpp:67
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
LocalOperators(const Eigen::MatrixXd &_proj_curl, const Eigen::MatrixXd &_dofs_curl, const Eigen::MatrixXd &_proj_function)
Definition vcurl.hpp:35
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition vemcore.hpp:145
const size_t & degree() const
Return the polynomial degree.
Definition vemcore.hpp:123
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition vemcore.hpp:129
const VEMCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition vcurl.hpp:113
const Mesh & mesh() const
Return the mesh.
Definition vcurl.hpp:56
Eigen::VectorXd interpolate(const FunctionType &v, const FunctionType &curl_v, const int doe_cell=-1, const int doe_face=-1, const int doe_edge=-1) const
Interpolator of a continuous function.
Definition vcurl.cpp:64
const LocalOperators & faceOperators(size_t iF) const
Return face operators for the face of index iF.
Definition vcurl.hpp:89
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition vcurl.hpp:77
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition vcurl.hpp:83
const Mesh & mesh() const
Return a const reference to the mesh.
Definition vemcore.hpp:117
const size_t & degree() const
Return the polynomial degree.
Definition vcurl.hpp:62
Eigen::MatrixXd dofs_curl
Definition vcurl.hpp:48
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition vcurl.hpp:30
const VEMCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition vcurl.hpp:131
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition vemcore.hpp:137
const VEMCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition vcurl.hpp:107
const LocalOperators & faceOperators(const Face &F) const
Return face operators for face F.
Definition vcurl.hpp:95
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.
Definition vcurl.cpp:359
Eigen::MatrixXd proj_function
Definition vcurl.hpp:49
const VEMCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition vcurl.hpp:125
Eigen::MatrixXd proj_curl
Definition vcurl.hpp:47
const VEMCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition vcurl.hpp:119
const VEMCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition vcurl.hpp:101
Definition ddr-magnetostatics.hpp:41
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:36
A structure to store the local operators (projections of curl and function, dof of curl in Vdiv)
Definition vcurl.hpp:34
Structure to store element bases.
Definition vemcore.hpp:63
Structure to store edge bases.
Definition vemcore.hpp:104
Structure to store face bases.
Definition vemcore.hpp:84