HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
sxcurl.hpp
Go to the documentation of this file.
1#ifndef SXCURL_HPP
2#define SXCURL_HPP
3
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
7#include <xcurl.hpp>
8#include <sxgrad.hpp>
10
11namespace HArDCore2D
12{
14
15 class SXCurl : public VariableDOFSpace
16 {
17 public:
18 typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> FunctionType;
19
21
23 {
25 const Eigen::MatrixXd & _serendipity,
26 const Eigen::MatrixXd & _extension,
27 const Eigen::MatrixXd & _reduction
28 )
29 : serendipity(_serendipity),
30 extension(_extension),
31 reduction(_reduction)
32 {
33 // Do nothing
34 }
35
36 Eigen::MatrixXd serendipity;
37 Eigen::MatrixXd extension;
38 Eigen::MatrixXd reduction;
39 };
40
42 SXCurl(const DDRCore & ddr_core, const SerendipityProblem & ser_pro, bool use_threads = true, std::ostream & output = std::cout);
43
45 const Mesh & mesh() const
46 {
47 return m_ddr_core.mesh();
48 }
49
51 const size_t & degree() const
52 {
53 return m_ddr_core.degree();
54 }
56 Eigen::VectorXd interpolate(
57 const FunctionType & v,
58 const int deg_quad = -1
59 ) const;
60 //----------------------------------//
61 //---- Cell transfer operators -----//
63 inline const Eigen::MatrixXd & ScurlCell(size_t iT) const
64 {
65 return (*m_cell_transfer_operators[iT]).serendipity;
66 }
67
69 inline const Eigen::MatrixXd & EcurlCell(size_t iT) const
70 {
71 return (*m_cell_transfer_operators[iT]).extension;
72 }
73
75 inline const Eigen::MatrixXd & RcurlCell(size_t iT) const
76 {
77 return (*m_cell_transfer_operators[iT]).reduction;
78 }
79
81 inline const Eigen::MatrixXd & ScurlCell(const Cell & T) const
82 {
83 return ScurlCell(T.global_index());
84 }
85
87 inline const Eigen::MatrixXd & EcurlCell(const Cell & T) const
88 {
89 return EcurlCell(T.global_index());
90 }
91
93 inline const Eigen::MatrixXd & RcurlCell(const Cell & T) const
94 {
95 return RcurlCell(T.global_index());
96 }
97
98 //-----------------------------------------------------------------//
99 //---- Full curl and potential reconstructions, and L2 product ----//
100
102 inline const Eigen::MatrixXd cellCurl(size_t iT) const
103 {
104 return m_xcurl.cellOperators(iT).curl * EcurlCell(iT);
105 }
106
108 inline const Eigen::MatrixXd cellCurl(const Cell & T) const
109 {
110 return cellCurl(T.global_index());
111 }
112
114 inline const Eigen::MatrixXd cellPotential(size_t iT) const
115 {
116 return m_xcurl.cellOperators(iT).potential * EcurlCell(iT);
117 }
118
120 inline const Eigen::MatrixXd cellPotential(const Cell & T) const
121 {
122 return cellPotential(T.global_index());
123 }
124
126 Eigen::MatrixXd computeL2Product(
127 const size_t iT,
128 const double & penalty_factor = 1.,
129 const Eigen::MatrixXd & mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
130 const IntegralWeight & weight = IntegralWeight(1.)
131 ) const
132 {
133 return EcurlCell(iT).transpose()
134 * m_xcurl.computeL2Product(iT, penalty_factor, mass_Pk2_T, weight)
135 *EcurlCell(iT);
136 }
137
138
140 Eigen::MatrixXd computeL2ProductGradient(
141 const size_t iT,
142 const SXGrad & sx_grad,
143 const std::string & side,
144 const double & penalty_factor = 1.,
145 const Eigen::MatrixXd & mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
146 const IntegralWeight & weight = IntegralWeight(1.)
147 ) const;
148
149 //-----------------------//
150 //---- Getters ----------//
151
153 inline const DDRCore::CellBases & cellBases(size_t iT) const
154 {
155 return m_ddr_core.cellBases(iT);
156 }
157
159 inline const DDRCore::CellBases & cellBases(const Cell & T) const
160 {
161 return m_ddr_core.cellBases(T.global_index());
162 }
163
165 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
166 {
167 return m_ddr_core.edgeBases(iE);
168 }
169
171 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
172 {
173 return m_ddr_core.edgeBases(E.global_index());
174 }
175
176 private:
177 TransferOperators _compute_cell_transfer_operators(size_t iT);
178
179 const DDRCore & m_ddr_core;
180 const SerendipityProblem & m_ser_pro;
181 const XCurl m_xcurl;
182
183 // Containers for serendipity, extension and reduction operators
184 std::vector<std::unique_ptr<TransferOperators> > m_cell_transfer_operators;
185
186 bool m_use_threads;
187 std::ostream & m_output;
188
189 };
190
191} // end of namespace HArDCore2D
192#endif
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:63
Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator.
Definition sxcurl.hpp:16
const Eigen::MatrixXd & RcurlCell(size_t iT) const
Return the reduction for the cell of index iT.
Definition sxcurl.hpp:75
const Eigen::MatrixXd & ScurlCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition sxcurl.hpp:63
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the face of index iT.
Definition sxcurl.hpp:153
const Eigen::MatrixXd cellPotential(const Cell &T) const
Return the potential operator on cell T.
Definition sxcurl.hpp:120
Eigen::VectorXd interpolate(const FunctionType &v, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition sxcurl.cpp:51
const Eigen::MatrixXd cellCurl(size_t iT) const
Return the full curl operator on the cell of index iT.
Definition sxcurl.hpp:102
const Eigen::MatrixXd cellPotential(size_t iT) const
Return the potential operator on the cell of index iT.
Definition sxcurl.hpp:114
const Eigen::MatrixXd cellCurl(const Cell &T) const
Return the full curl operator on cell T.
Definition sxcurl.hpp:108
const Eigen::MatrixXd & EcurlCell(size_t iT) const
Return the extension for the cell of index iT.
Definition sxcurl.hpp:69
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition sxcurl.hpp:159
const Eigen::MatrixXd & EcurlCell(const Cell &T) const
Return the extension for cell T.
Definition sxcurl.hpp:87
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const SXGrad &sx_grad, const std::string &side, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_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 sxcurl.cpp:186
const size_t & degree() const
Return the polynomial degree.
Definition sxcurl.hpp:51
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition sxcurl.hpp:171
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition sxcurl.hpp:165
const Mesh & mesh() const
Return the mesh.
Definition sxcurl.hpp:45
std::function< Eigen::Vector2d(const Eigen::Vector2d &)> FunctionType
Definition sxcurl.hpp:18
const Eigen::MatrixXd & ScurlCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition sxcurl.hpp:81
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_T=Eigen::MatrixXd::Zero(1, 1), const IntegralWeight &weight=IntegralWeight(1.)) const
Compute the matrix of the (weighted) L2-product.
Definition sxcurl.hpp:126
const Eigen::MatrixXd & RcurlCell(const Cell &T) const
Return the reduction for cell T.
Definition sxcurl.hpp:93
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:20
Base class for global DOF spaces.
Definition variabledofspace.hpp:17
Discrete Hcurl space: local operators, L2 product and global interpolator.
Definition xcurl.hpp:18
Definition Mesh2D.hpp:26
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:116
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:122
Eigen::MatrixXd potential
Definition xcurl.hpp:36
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xcurl.hpp:61
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pk2_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 xcurl.cpp:189
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:128
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:136
Eigen::MatrixXd curl
Definition xcurl.hpp:35
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
Definition ddr-klplate.hpp:27
static auto v
Definition ddrcore-test.hpp:32
Structure to store element bases.
Definition ddrcore.hpp:81
Structure to store edge bases.
Definition ddrcore.hpp:103
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:33
A structure to store the serendipity, extension and reduction operators.
Definition sxcurl.hpp:23
TransferOperators(const Eigen::MatrixXd &_serendipity, const Eigen::MatrixXd &_extension, const Eigen::MatrixXd &_reduction)
Definition sxcurl.hpp:24
Eigen::MatrixXd extension
Definition sxcurl.hpp:37
Eigen::MatrixXd reduction
Definition sxcurl.hpp:38
Eigen::MatrixXd serendipity
Definition sxcurl.hpp:36