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
laxdiv.hpp
Go to the documentation of this file.
1#ifndef LAXDIV_HPP
2#define LAXDIV_HPP
3
4#include <globaldofspace.hpp>
5#include <xdiv.hpp>
6#include <liealgebra.hpp>
7
8namespace HArDCore3D
9{
16
19 class LAXDiv : public GlobalDOFSpace
20 {
21 public:
22
23 typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType;
24 typedef std::vector<FunctionType> LAFunctionType;
25
28 {
30 const Eigen::MatrixXd & _divergence,
31 const Eigen::MatrixXd & _divergence_rhs,
32 const Eigen::MatrixXd & _potential
33 )
37 {
38 // Do nothing
39 }
40
41 Eigen::MatrixXd divergence;
42 Eigen::MatrixXd divergence_rhs;
43 Eigen::MatrixXd potential;
44 };
45
47 LAXDiv(const LieAlgebra & lie_algebra, const XDiv & xdiv, bool use_threads = true, std::ostream & output = std::cout);
48
50 const Mesh & mesh() const
51 {
52 return m_xdiv.mesh();
53 }
54
56 const size_t & degree() const
57 {
58 return m_xdiv.degree();
59 }
60
62 Eigen::VectorXd interpolate(
63 const LAFunctionType & v,
64 const int doe_cell = -1,
65 const int doe_face = -1
66 ) const;
67
69 inline const LocalOperators & cellOperators(size_t iT) const
70 {
71 return *m_cell_operators[iT];
72 }
73
75 inline const LocalOperators & cellOperators(const Cell & T) const
76 {
77 return * m_cell_operators[T.global_index()];
78 }
79
81 inline const DDRCore::CellBases & cellBases(size_t iT) const
82 {
83 return m_xdiv.cellBases(iT);
84 }
85
87 inline const DDRCore::CellBases & cellBases(const Cell & T) const
88 {
89 return m_xdiv.cellBases(T.global_index());
90 }
91
93 inline const DDRCore::FaceBases & faceBases(size_t iF) const
94 {
95 return m_xdiv.faceBases(iF);
96 }
97
99 inline const DDRCore::FaceBases & faceBases(const Face & F) const
100 {
101 return m_xdiv.faceBases(F.global_index());
102 }
103
105 // The mass matrix of P^k(T)^3 is the most expensive mass matrix in the calculation of this norm, which
106 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
107 Eigen::MatrixXd computeL2Product(
108 const size_t iT,
109 const double & penalty_factor = 1.,
110 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
112 ) const;
113
115 Eigen::MatrixXd computeL2ProductCurl(
116 const size_t iT,
117 const XCurl & x_curl,
118 const std::string & side,
119 const double & penalty_factor = 1.,
120 const Eigen::MatrixXd & mass_Pk3_T = Eigen::MatrixXd::Zero(1,1),
122 ) const;
123
124 protected:
125 LocalOperators _compute_cell_divergence_potential(size_t iT);
126
128 const XDiv & m_xdiv;
130 std::ostream & m_output;
131
132 // Containers for local operators
133 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
134 };
135} // end of namespace HArDCore3D
136
137#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 Hdiv space: local operators, L2 product and global interpolator.
Definition laxdiv.hpp:20
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 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
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition xdiv.hpp:92
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 the mesh.
Definition xdiv.hpp:48
const size_t & degree() const
Return the polynomial degree.
Definition xdiv.hpp:54
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
const LieAlgebra & m_lie_algebra
Definition laxdiv.hpp:127
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition laxdiv.hpp:69
std::ostream & m_output
Definition laxdiv.hpp:130
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition laxdiv.hpp:81
const size_t & degree() const
Return the polynomial degree.
Definition laxdiv.hpp:56
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition laxdiv.hpp:87
const Mesh & mesh() const
Return the mesh.
Definition laxdiv.hpp:50
Eigen::MatrixXd divergence
Definition laxdiv.hpp:41
Eigen::MatrixXd divergence_rhs
Definition laxdiv.hpp:42
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> FunctionType
Definition laxdiv.hpp:23
Eigen::VectorXd interpolate(const LAFunctionType &v, const int doe_cell=-1, const int doe_face=-1) const
Interpolator of a continuous function.
Definition laxdiv.cpp:51
const DDRCore::FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition laxdiv.hpp:99
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 laxdiv.cpp:101
bool m_use_threads
Definition laxdiv.hpp:129
std::vector< FunctionType > LAFunctionType
Definition laxdiv.hpp:24
Eigen::MatrixXd potential
Definition laxdiv.hpp:43
LocalOperators(const Eigen::MatrixXd &_divergence, const Eigen::MatrixXd &_divergence_rhs, const Eigen::MatrixXd &_potential)
Definition laxdiv.hpp:29
const XDiv & m_xdiv
Definition laxdiv.hpp:128
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 laxdiv.cpp:87
LocalOperators _compute_cell_divergence_potential(size_t iT)
Definition laxdiv.cpp:70
std::vector< std::unique_ptr< LocalOperators > > m_cell_operators
Definition laxdiv.hpp:133
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition laxdiv.hpp:75
const DDRCore::FaceBases & faceBases(size_t iF) const
Return face bases for the face of index iF.
Definition laxdiv.hpp:93
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 laxdiv.hpp:28