HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
xcurl.hpp
Go to the documentation of this file.
1#ifndef XCURL_HPP
2#define XCURL_HPP
3
4#include <globaldofspace.hpp>
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
7#include <xgrad.hpp>
8
9namespace HArDCore2D
10{
17 class XCurl : public GlobalDOFSpace
18 {
19 public:
20 typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> FunctionType;
21
24 {
26 const Eigen::MatrixXd & _curl,
27 const Eigen::MatrixXd & _potential
28 )
29 : curl(_curl),
30 potential(_potential)
31 {
32 // Do nothing
33 }
34
35 Eigen::MatrixXd curl;
36 Eigen::MatrixXd potential;
37 };
38
40 XCurl(const DDRCore & ddr_core, bool use_threads = true, std::ostream & output = std::cout);
41
43 const Mesh & mesh() const
44 {
45 return m_ddr_core.mesh();
46 }
47
49 const size_t & degree() const
50 {
51 return m_ddr_core.degree();
52 }
53
55 Eigen::VectorXd interpolate(
56 const FunctionType & v,
57 const int deg_quad = -1
58 ) const;
59
61 inline const LocalOperators & cellOperators(size_t iT) const
62 {
63 return *m_cell_operators[iT];
64 }
65
67 inline const LocalOperators & cellOperators(const Cell & T) const
68 {
69 return * m_cell_operators[T.global_index()];
70 }
71
73 inline const DDRCore::CellBases & cellBases(size_t iT) const
74 {
75 return m_ddr_core.cellBases(iT);
76 }
77
79 inline const DDRCore::CellBases & cellBases(const Cell & T) const
80 {
81 return m_ddr_core.cellBases(T.global_index());
82 }
83
85 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
86 {
87 return m_ddr_core.edgeBases(iE);
88 }
89
91 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
92 {
93 return m_ddr_core.edgeBases(E.global_index());
94 }
95
97 // The mass matrix of P^k(T)^2 is the most expensive mass matrix in the calculation of this norm, which
98 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
99 Eigen::MatrixXd computeL2Product(
100 const size_t iT,
101 const double & penalty_factor = 1.,
102 const Eigen::MatrixXd & mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
103 const IntegralWeight & weight = IntegralWeight(1.)
104 ) const;
105
107 Eigen::MatrixXd computeL2ProductGradient(
108 const size_t iT,
109 const XGrad & x_grad,
110 const std::string & side,
111 const double & penalty_factor = 1.,
112 const Eigen::MatrixXd & mass_Pk2_T = Eigen::MatrixXd::Zero(1,1),
113 const IntegralWeight & weight = IntegralWeight(1.)
114 ) const;
115
117 Eigen::MatrixXd computeL2Product_with_Ops(
118 const size_t iT,
119 const std::vector<Eigen::MatrixXd> & leftOp,
120 const std::vector<Eigen::MatrixXd> & rightOp,
121 const double & penalty_factor,
122 const Eigen::MatrixXd & w_mass_Pk2_T,
123 const IntegralWeight & weight
124 ) const;
125
126
127 private:
128 LocalOperators _compute_cell_curl_potential(size_t iT);
129
130 const DDRCore & m_ddr_core;
131 bool m_use_threads;
132 std::ostream & m_output;
133
134 // Containers for local operators
135 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
136 };
137
138} // end of namespace HArDCore2D
139#endif
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:63
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:18
Discrete H1 space: local operators, L2 product and global interpolator.
Definition xgrad.hpp:17
Definition Mesh2D.hpp:26
const Mesh & mesh() const
Return the mesh.
Definition xcurl.hpp:43
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:115
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:121
const size_t & degree() const
Return the polynomial degree.
Definition xcurl.hpp:49
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
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xcurl.hpp:79
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xcurl.hpp:67
LocalOperators(const Eigen::MatrixXd &_curl, const Eigen::MatrixXd &_potential)
Definition xcurl.hpp:25
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
Eigen::MatrixXd computeL2ProductGradient(const size_t iT, const XGrad &x_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 xcurl.cpp:233
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:127
std::function< Eigen::Vector2d(const Eigen::Vector2d &)> FunctionType
Definition xcurl.hpp:20
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xcurl.hpp:73
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xcurl.hpp:85
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:135
Eigen::VectorXd interpolate(const FunctionType &v, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition xcurl.cpp:48
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xcurl.hpp:91
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_Pk2_T, const IntegralWeight &weight) const
Compute the matrix of the L2 product, applying leftOp and rightOp to the variables....
Definition xcurl.cpp:300
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:102
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:33
A structure to store the local operators (curl and potential)
Definition xcurl.hpp:24