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
xdiv.hpp
Go to the documentation of this file.
1#ifndef XDIV_HPP
2#define XDIV_HPP
3
4#include <globaldofspace.hpp>
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
7#include <xcurl.hpp>
8
9namespace HArDCore3D
10{
17
19 class XDiv : public GlobalDOFSpace
20 {
21 public:
22 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
23
26 {
28 const Eigen::MatrixXd & _divergence,
29 const Eigen::MatrixXd & _divergence_rhs,
30 const Eigen::MatrixXd & _potential
31 )
35 {
36 // Do nothing
37 }
38
39 Eigen::MatrixXd divergence;
40 Eigen::MatrixXd divergence_rhs;
41 Eigen::MatrixXd potential;
42 };
43
45 XDiv(const DDRCore & ddr_core, bool use_threads = true, std::ostream & output = std::cout);
46
48 const Mesh & mesh() const
49 {
50 return m_ddr_core.mesh();
51 }
52
54 const size_t & degree() const
55 {
56 return m_ddr_core.degree();
57 }
58
60 Eigen::VectorXd interpolate(
61 const FunctionType & v,
62 const int doe_cell = -1,
63 const int doe_face = -1
64 ) const;
65
66
68 inline const LocalOperators & cellOperators(size_t iT) const
69 {
70 return *m_cell_operators[iT];
71 }
72
74 inline const LocalOperators & cellOperators(const Cell & T) const
75 {
76 return * m_cell_operators[T.global_index()];
77 }
78
80 inline const DDRCore::CellBases & cellBases(size_t iT) const
81 {
82 return m_ddr_core.cellBases(iT);
83 }
84
86 inline const DDRCore::CellBases & cellBases(const Cell & T) const
87 {
88 return m_ddr_core.cellBases(T.global_index());
89 }
90
92 inline const DDRCore::FaceBases & faceBases(size_t iF) const
93 {
94 return m_ddr_core.faceBases(iF);
95 }
96
98 inline const DDRCore::FaceBases & faceBases(const Face & F) const
99 {
100 return m_ddr_core.faceBases(F.global_index());
101 }
102
104 // The mass matrix of P^k(T)^3 is the most expensive mass matrix in the calculation of this norm, which
105 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
106 Eigen::MatrixXd computeL2Product(
107 const size_t iT,
108 const double & penalty_factor = 1.,
109 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
111 ) const;
112
114 Eigen::MatrixXd computeL2ProductCurl(
115 const size_t iT,
116 const XCurl & x_curl,
117 const std::string & side,
118 const double & penalty_factor = 1.,
119 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
121 ) const;
122
124 Eigen::MatrixXd computeL2Product_with_Ops(
125 const size_t iT,
126 const std::vector<Eigen::MatrixXd> & leftOp,
127 const std::vector<Eigen::MatrixXd> & rightOp,
128 const double & penalty_factor,
129 const Eigen::MatrixXd & w_mass_Pk3_T,
130 const IntegralWeight & weight
131 ) const;
132
133
134 protected:
135 LocalOperators _compute_cell_divergence_potential(size_t iT);
136
139 std::ostream & m_output;
140
141 // Containers for local operators
142 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
143 };
144} // end of namespace HArDCore3D
145
146#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 Hdiv space: local operators, L2 product and global interpolator.
Definition xdiv.hpp:20
Class to describe a mesh.
Definition MeshND.hpp:17
@ Matrix
Definition basis.hpp:67
LocalOperators _compute_cell_divergence_potential(size_t iT)
Definition xdiv.cpp:108
const DDRCore & m_ddr_core
Definition xdiv.hpp:137
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition xdiv.hpp:92
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 xdiv.cpp:326
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition xdiv.hpp:22
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xdiv.hpp:86
Eigen::MatrixXd computeL2ProductCurl(const size_t iT, const XCurl &x_curl, 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 xdiv.cpp:259
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition xdiv.hpp:80
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:134
Eigen::MatrixXd divergence_rhs
Definition xdiv.hpp:40
const Mesh & mesh() const
Return the mesh.
Definition xdiv.hpp:48
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xdiv.hpp:68
Eigen::MatrixXd divergence
Definition xdiv.hpp:39
bool m_use_threads
Definition xdiv.hpp:138
std::ostream & m_output
Definition xdiv.hpp:139
Eigen::MatrixXd potential
Definition xdiv.hpp:41
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 xdiv.cpp:216
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:140
std::vector< std::unique_ptr< LocalOperators > > m_cell_operators
Definition xdiv.hpp:142
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition xdiv.hpp:98
Eigen::VectorXd interpolate(const FunctionType &v, const int doe_cell=-1, const int doe_face=-1) const
Interpolator of a continuous function.
Definition xdiv.cpp:48
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xdiv.hpp:74
LocalOperators(const Eigen::MatrixXd &_divergence, const Eigen::MatrixXd &_divergence_rhs, const Eigen::MatrixXd &_potential)
Definition xdiv.hpp:27
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
const size_t & degree() const
Return the polynomial degree.
Definition xdiv.hpp:54
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 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 (divergence and potential)
Definition xdiv.hpp:26