HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
xhess.hpp
Go to the documentation of this file.
1#ifndef XHESS_HPP
2#define XHESS_HPP
3
4#include <globaldofspace.hpp>
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
8#include "xgrad.hpp"
9#include "sxgrad.hpp"
10
11namespace HArDCore2D
12{
19 class XHess : public GlobalDOFSpace
20 {
21 public:
22 typedef std::function<double(const Eigen::Vector2d &)> FunctionType;
23 typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> GradFunctionType;
24
25
28 {
30 const Eigen::MatrixXd & _gradient,
31 const Eigen::MatrixXd & _potential
32 )
33 : gradient(_gradient),
34 potential(_potential)
35 {
36 // Do nothing
37 }
38
39 Eigen::MatrixXd gradient;
40 Eigen::MatrixXd potential;
41 };
43 {
45 const Eigen::MatrixXd & _serendipity
46 )
47 : serendipity(_serendipity)
48 {
49 // Do nothing
50 }
51
52 Eigen::MatrixXd serendipity;
53 };
55 XHess(const DDRCore & ddr_core, const SerendipityProblem & ser_pro, bool use_threads = true, std::ostream & output = std::cout);//ajout senredipity ser_pro
56
58 const Mesh & mesh() const
59 {
60 return m_ddr_core.mesh();
61 }
62
64 const size_t & degree() const
65 {
66 return m_ddr_core.degree();
67 }
68
70 Eigen::VectorXd interpolate(
71 const FunctionType & q,
72 const GradFunctionType & Dq,
73 const int deg_quad = -1
74 ) const;
75
76 //---------------------------------//
77 //---- Cell transfer operators -----//
79 inline const Eigen::MatrixXd & SgradCell(size_t iT) const
80 {
81 return (*m_cell_transfer_operators[iT]).serendipity;
82 }
83
84
86 inline const Eigen::MatrixXd & SgradCell(const Cell & T) const
87 {
88 return SgradCell(T.global_index());
89 }
90
91
93 inline const TransferOperators & TcellOperators(size_t iT) const
94 {
95 return *m_cell_transfer_operators[iT];
96 }
97
99 inline const TransferOperators & TcellOperators(const Cell & T) const
100 {
101 return *m_cell_transfer_operators[T.global_index()];
102 }
103
104
106 inline const LocalOperators & edgeOperators(size_t iE) const
107 {
108 return *m_edge_operators[iE];
109 }
110
112 inline const LocalOperators & edgeOperators(const Edge & E) const
113 {
114 return *m_edge_operators[E.global_index()];
115 }
116
118 inline const LocalOperators & cellOperators(size_t iT) const
119 {
120 return *m_cell_operators[iT];
121 }
122
124 inline const LocalOperators & cellOperators(const Cell & T) const
125 {
126 return *m_cell_operators[T.global_index()];
127 }
128
130 inline const DDRCore::CellBases & cellBases(size_t iT) const
131 {
132 return m_ddr_core.cellBases(iT);
133 }
134
136 inline const DDRCore::CellBases & cellBases(const Cell & T) const
137 {
138 return m_ddr_core.cellBases(T.global_index());
139 }
140
142 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
143 {
144 return m_ddr_core.edgeBases(iE);
145 }
146
148 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
149 {
150 return m_ddr_core.edgeBases(E.global_index());
151 }
152
154 // The mass matrix of P^{k+1}(T) is the most expensive mass matrix in the calculation of this norm, which
155 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
156 Eigen::MatrixXd computeL2Product(
157 const size_t iT,
158 const double & penalty_factor = 1.,
159 const Eigen::MatrixXd & mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1),
160 const IntegralWeight & weight = IntegralWeight(1.)
161 ) const;
162
164 double evaluatePotential(
165 const size_t iT,
166 const Eigen::VectorXd & vT,
167 const VectorRd & x
168 ) const;
169
171 double computeL2Norm(const Eigen::VectorXd & v) const;
172
173 private:
174 Eigen::MatrixXd _compute_cell_serendipity_operator(size_t iT); //ajout serendipity
175
176 LocalOperators _compute_edge_gradient_potential(size_t iE);
177 LocalOperators _compute_cell_gradient_potential(size_t iT);
178 double _compute_squared_l2_norm(size_t iT, const Eigen::VectorXd & vT) const;
179
180
181
182 const DDRCore & m_ddr_core;
183
184 const SerendipityProblem & m_ser_pro;
185 const SXGrad m_sxgrad;
186
187 bool m_use_threads;
188 std::ostream & m_output;
189
190 // Containers for local operators
191 std::vector<std::unique_ptr<LocalOperators> > m_edge_operators;
192 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
193
194 // Containers for serendipity
195 std::vector<std::shared_ptr<TransferOperators> > m_cell_transfer_operators; //ajout serendipity
196 };
197
198} // end of namespace HArDCore2D
199
200#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 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
Discrete H2 space: local operators, L2 product and global interpolator.
Definition xhess.hpp:20
Definition Mesh2D.hpp:26
Eigen::Vector2d VectorRd
Definition basis.hpp:55
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xhess.hpp:118
const TransferOperators & TcellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xhess.hpp:99
const size_t & degree() const
Return the polynomial degree.
Definition xhess.hpp:64
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:116
double evaluatePotential(const size_t iT, const Eigen::VectorXd &vT, const VectorRd &x) const
Evaluate the value of the potential at a point x.
Definition xhess.cpp:355
std::function< Eigen::Vector2d(const Eigen::Vector2d &)> GradFunctionType
Definition xhess.hpp:23
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:122
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xhess.hpp:148
Eigen::VectorXd interpolate(const FunctionType &q, const GradFunctionType &Dq, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition xhess.cpp:65
const Eigen::MatrixXd & SgradCell(const Cell &T) const
Return the serendipity reconstruction for cell T.
Definition xhess.hpp:86
LocalOperators(const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential)
Definition xhess.hpp:29
Eigen::MatrixXd potential
Definition xhess.hpp:40
Eigen::MatrixXd gradient
Definition xhess.hpp:39
const LocalOperators & edgeOperators(const Edge &E) const
Return edge operators for edge E.
Definition xhess.hpp:112
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xhess.hpp:142
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xhess.hpp:136
Eigen::MatrixXd computeL2Product(const size_t iT, const double &penalty_factor=1., const Eigen::MatrixXd &mass_Pkpo_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 xhess.cpp:372
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:128
const Eigen::MatrixXd & SgradCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition xhess.hpp:79
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xhess.hpp:130
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xhess.hpp:124
const TransferOperators & TcellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xhess.hpp:93
std::function< double(const Eigen::Vector2d &)> FunctionType
Definition xhess.hpp:22
const Mesh & mesh() const
Return the mesh.
Definition xhess.hpp:58
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:136
TransferOperators(const Eigen::MatrixXd &_serendipity)
Definition xhess.hpp:44
const LocalOperators & edgeOperators(size_t iE) const
Return edge operators for the edge of index iE.
Definition xhess.hpp:106
Eigen::MatrixXd serendipity
Definition xhess.hpp:52
double computeL2Norm(const Eigen::VectorXd &v) const
Compute the L2-norm of a vector of the space.
Definition xhess.cpp:531
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
static auto q
Definition ddrcore-test.hpp:14
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 local operators (gradient and potential)
Definition xhess.hpp:28