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#include "vsxgrad.hpp"
11
12namespace HArDCore2D
13{
20 class XHess : public GlobalDOFSpace
21 {
22 public:
23 typedef std::function<double(const Eigen::Vector2d &)> FunctionType;
24 typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> GradFunctionType;
25
26
29 {
31 const Eigen::MatrixXd & _gradient,
32 const Eigen::MatrixXd & _potential
33 )
34 : gradient(_gradient),
35 potential(_potential)
36 {
37 // Do nothing
38 }
39
40 Eigen::MatrixXd gradient;
41 Eigen::MatrixXd potential;
42 };
44 {
46 const Eigen::MatrixXd & _serendipity
47 )
48 : serendipity(_serendipity)
49 {
50 // Do nothing
51 }
52
53 Eigen::MatrixXd serendipity;
54 };
56 XHess(const DDRCore & ddr_core, const SerendipityProblem & ser_pro, bool use_threads = true, std::ostream & output = std::cout);//ajout senredipity ser_pro
57
59 const Mesh & mesh() const
60 {
61 return m_ddr_core.mesh();
62 }
63
65 const size_t & degree() const
66 {
67 return m_ddr_core.degree();
68 }
69
71 Eigen::VectorXd interpolate(
72 const FunctionType & q,
73 const GradFunctionType & Dq,
74 const int deg_quad = -1
75 ) const;
76
77 //---------------------------------//
78 //---- Cell transfer operators -----//
80 inline const Eigen::MatrixXd & SgradCell(size_t iT) const
81 {
82 return (*m_cell_transfer_operators[iT]).serendipity;
83 }
84
86 inline const Eigen::MatrixXd & SgradCell(const Cell & T) const
87 {
88 return SgradCell(T.global_index());
89 }
90
92 inline const TransferOperators & TcellOperators(size_t iT) const
93 {
94 return *m_cell_transfer_operators[iT];
95 }
96
98 inline const TransferOperators & TcellOperators(const Cell & T) const
99 {
100 return *m_cell_transfer_operators[T.global_index()];
101 }
102
104 inline const LocalOperators & edgeOperators(size_t iE) const
105 {
106 return *m_edge_operators[iE];
107 }
108
110 inline const LocalOperators & edgeOperators(const Edge & E) const
111 {
112 return *m_edge_operators[E.global_index()];
113 }
114
116 inline const LocalOperators & cellOperators(size_t iT) const
117 {
118 return *m_cell_operators[iT];
119 }
120
122 inline const LocalOperators & cellOperators(const Cell & T) const
123 {
124 return *m_cell_operators[T.global_index()];
125 }
126
128 inline const DDRCore::CellBases & cellBases(size_t iT) const
129 {
130 return m_ddr_core.cellBases(iT);
131 }
132
134 inline const DDRCore::CellBases & cellBases(const Cell & T) const
135 {
136 return m_ddr_core.cellBases(T.global_index());
137 }
138
140 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
141 {
142 return m_ddr_core.edgeBases(iE);
143 }
144
146 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
147 {
148 return m_ddr_core.edgeBases(E.global_index());
149 }
150
152 // The mass matrix of P^{k+1}(T) is the most expensive mass matrix in the calculation of this norm, which
153 // is why there's the option of passing it as parameter if it's been already pre-computed when the norm is called.
154 Eigen::MatrixXd computeL2Product(
155 const size_t iT,
156 const double & penalty_factor = 1.,
157 const Eigen::MatrixXd & mass_Pkpo_T = Eigen::MatrixXd::Zero(1,1),
158 const IntegralWeight & weight = IntegralWeight(1.)
159 ) const;
160
161
162 Eigen::MatrixXd computeGradientFull(
163 const size_t iT,
164 const VSXGrad & vsx_grad
165 ) const;
166
168 double evaluatePotential(
169 const size_t iT,
170 const Eigen::VectorXd & vT,
171 const VectorRd & x
172 ) const;
173
175 double computeL2Norm(const Eigen::VectorXd & v) const;
176
177 private:
178 Eigen::MatrixXd _compute_cell_serendipity_operator(size_t iT); //ajout serendipity
179
180 LocalOperators _compute_edge_gradient_potential(size_t iE);
181 LocalOperators _compute_cell_gradient_potential(size_t iT);
182 double _compute_squared_l2_norm(size_t iT, const Eigen::VectorXd & vT) const;
183
184 const DDRCore & m_ddr_core;
185 const SerendipityProblem & m_ser_pro;
186 const SXGrad m_sxgrad;
187
188 bool m_use_threads;
189 std::ostream & m_output;
190
191 // Containers for local operators
192 std::vector<std::unique_ptr<LocalOperators> > m_edge_operators;
193 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
194
195 // Containers for serendipity
196 std::vector<std::shared_ptr<TransferOperators> > m_cell_transfer_operators; //ajout serendipity
197 };
198
199} // end of namespace HArDCore2D
200
201#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
Vector version of sXgrad, the arbitrary order space with nodal primal unknowns.
Definition vsxgrad.hpp:32
Discrete H2 space: local operators, L2 product and global interpolator.
Definition xhess.hpp:21
Definition Mesh2D.hpp:26
Create grid points x
Definition generate_cartesian_mesh.m:22
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:116
const TransferOperators & TcellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xhess.hpp:98
const size_t & degree() const
Return the polynomial degree.
Definition xhess.hpp:65
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:346
std::function< Eigen::Vector2d(const Eigen::Vector2d &)> GradFunctionType
Definition xhess.hpp:24
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:146
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:30
Eigen::MatrixXd potential
Definition xhess.hpp:41
Eigen::MatrixXd gradient
Definition xhess.hpp:40
const LocalOperators & edgeOperators(const Edge &E) const
Return edge operators for edge E.
Definition xhess.hpp:110
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xhess.hpp:140
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xhess.hpp:134
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:362
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:128
Eigen::MatrixXd computeGradientFull(const size_t iT, const VSXGrad &vsx_grad) const
Definition xhess.cpp:511
const Eigen::MatrixXd & SgradCell(size_t iT) const
Return the serendipity reconstruction for the cell of index iT.
Definition xhess.hpp:80
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xhess.hpp:128
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xhess.hpp:122
const TransferOperators & TcellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xhess.hpp:92
std::function< double(const Eigen::Vector2d &)> FunctionType
Definition xhess.hpp:23
const Mesh & mesh() const
Return the mesh.
Definition xhess.hpp:59
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:45
const LocalOperators & edgeOperators(size_t iE) const
Return edge operators for the edge of index iE.
Definition xhess.hpp:104
Eigen::MatrixXd serendipity
Definition xhess.hpp:53
double computeL2Norm(const Eigen::VectorXd &v) const
Compute the L2-norm of a vector of the space.
Definition xhess.cpp:552
bool use_threads
Definition HHO_DiffAdvecReac.hpp:45
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:29