HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
hhospace.hpp
Go to the documentation of this file.
1// Core data structures and methods required to implement the Hybrid High-Order in 2D
2//
3// Provides:
4// - Polynomial spaces on the element and edges
5// - Interpolator of smooth functions
6// - Full gradient, potential and stabilisation bilinear form in the elements
7//
8// Author: Jerome Droniou (jerome.droniou@monash.edu)
9//
10
11/*
12 *
13 * This library was developed around HHO methods, although some parts of it have a more
14 * general purpose. If you use this code or part of it in a scientific publication,
15 * please mention the following book as a reference for the underlying principles
16 * of HHO schemes:
17 *
18 * The Hybrid High-Order Method for Polytopal Meshes: Design, Analysis, and Applications.
19 * D. A. Di Pietro and J. Droniou. Modeling, Simulation and Applications, vol. 19.
20 * Springer International Publishing, 2020, xxxi + 525p. doi: 10.1007/978-3-030-37203-3.
21 * url: https://hal.archives-ouvertes.fr/hal-02151813.
22 *
23 */
24
25
26#ifndef HHOSPACE_HPP
27#define HHOSPACE_HPP
28
29#include <iostream>
30#include <globaldofspace.hpp>
31#include <basis.hpp>
33
40namespace HArDCore2D
41{
42
49 //------------------------------------------------------------------------------
50
52 class HHOSpace : public GlobalDOFSpace
53 {
54 public:
55 // Types for element bases
58
59 // Types for edge bases
61
62 // Types for functions to interpolate
63 typedef std::function<double(const VectorRd &)> FunctionType;
64
66
70 struct CellBases
71 {
73 typedef Cell GeometricSupport;
74
75 std::unique_ptr<PolyBasisCellType> Polykpo;
76 std::unique_ptr<PolyBasisCellType> Polyk;
77 std::unique_ptr<PolydBasisCellType> Polykd;
78 };
79
81
82 struct EdgeBases
83 {
85 typedef Edge GeometricSupport;
86
87 std::unique_ptr<PolyBasisEdgeType> Polyk;
88 };
89
92 {
94 const Eigen::MatrixXd & _gradient,
95 const Eigen::MatrixXd & _potential,
96 const Eigen::MatrixXd & _stabilisation
97 )
98 : gradient(_gradient),
99 potential(_potential),
100 stabilisation(_stabilisation)
101 {
102 // Do nothing
103 }
104
105 Eigen::MatrixXd gradient;
106 Eigen::MatrixXd potential;
107 Eigen::MatrixXd stabilisation;
108 };
109
111 HHOSpace(const Mesh & mesh, size_t K, bool use_threads = true, std::ostream & output = std::cout);
112
114 const Mesh & mesh() const
115 {
116 return m_mesh;
117 }
118
120 const size_t & degree() const
121 {
122 return m_K;
123 }
124
126 Eigen::VectorXd interpolate(
127 const FunctionType & q,
128 const int doe_cell = -1,
129 const int doe_edge = -1
130 ) const;
131
133 inline const CellBases & cellBases(size_t iT) const
134 {
135 // Make sure that the basis has been created
136 assert( m_cell_bases[iT] );
137 return *m_cell_bases[iT].get();
138 }
139
141 inline const CellBases & cellBases(const Cell & T) const
142 {
143 return cellBases(T.global_index());
144 }
145
147 inline const EdgeBases & edgeBases(size_t iE) const
148 {
149 // Make sure that the basis has been created
150 assert( m_edge_bases[iE] );
151 return *m_edge_bases[iE].get();
152 }
153
155 inline const EdgeBases & edgeBases(const Edge & E) const
156 {
157 return edgeBases(E.global_index());
158 }
159
161 inline const LocalOperators & operators(size_t iT) const
162 {
163 assert( m_operators[iT] );
164 return *m_operators[iT];
165 }
166
168 inline const LocalOperators & operators(const Cell & T) const
169 {
170 return operators(T.global_index());
171 }
172
174 std::vector<std::pair<double,double> > computeNorms(
175 const std::vector<Eigen::VectorXd> & list_dofs
176 ) const;
177
179 Eigen::VectorXd computeVertexValues(
180 const Eigen::VectorXd & u
181 ) const;
182
183 private:
185 CellBases _construct_cell_bases(size_t iT);
186
188 EdgeBases _construct_edge_bases(size_t iE);
189
191 LocalOperators _compute_operators(size_t iT);
192
193 // Pointer to the mesh
194 const Mesh & m_mesh;
195 // Degrees
196 const size_t m_K;
197 // Parallel or not
198 bool m_use_threads;
199 // Output stream
200 std::ostream & m_output;
201
202 // Cell bases
203 std::vector<std::unique_ptr<CellBases> > m_cell_bases;
204 // Edge bases
205 std::vector<std::unique_ptr<EdgeBases> > m_edge_bases;
206
207 // Local operators
208 std::vector<std::unique_ptr<LocalOperators> > m_operators;
209
210 };
211
212} // end of namespace HArDCore2D
213
214#endif // HHOSPACE_HPP
Definition basis.hpp:311
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Class definition: polynomial bases and operators.
Definition hhospace.hpp:53
Vector family obtained by tensorization of a scalar family.
Definition basis.hpp:559
Definition Mesh2D.hpp:26
Eigen::Vector2d VectorRd
Definition basis.hpp:55
const size_t & degree() const
Return the polynomial degree (common edge and elements)
Definition hhospace.hpp:120
const CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition hhospace.hpp:141
Eigen::MatrixXd gradient
Definition hhospace.hpp:105
const Mesh & mesh() const
Return a const reference to the mesh.
Definition hhospace.hpp:114
std::function< double(const VectorRd &)> FunctionType
Definition hhospace.hpp:63
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition hhospace.hpp:133
Edge GeometricSupport
Geometric support.
Definition hhospace.hpp:85
Cell GeometricSupport
Geometric support.
Definition hhospace.hpp:73
std::vector< std::pair< double, double > > computeNorms(const std::vector< Eigen::VectorXd > &list_dofs) const
Computes the discrete L2 (cell unknowns only) and H1 norms of a list of vectors.
Definition hhospace.cpp:323
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition hhospace.hpp:147
Family< MonomialScalarBasisEdge > PolyBasisEdgeType
Definition hhospace.hpp:60
std::unique_ptr< PolyBasisCellType > Polyk
Definition hhospace.hpp:76
std::unique_ptr< PolydBasisCellType > Polykd
Definition hhospace.hpp:77
Eigen::VectorXd interpolate(const FunctionType &q, const int doe_cell=-1, const int doe_edge=-1) const
Interpolator of a continuous function.
Definition hhospace.cpp:128
std::unique_ptr< PolyBasisCellType > Polykpo
Definition hhospace.hpp:75
const LocalOperators & operators(const Cell &T) const
Return cell operators for cell T.
Definition hhospace.hpp:168
TensorizedVectorFamily< PolyBasisCellType, dimspace > PolydBasisCellType
Definition hhospace.hpp:57
const EdgeBases & edgeBases(const Edge &E) const
Return cell bases for edge E.
Definition hhospace.hpp:155
Eigen::VectorXd computeVertexValues(const Eigen::VectorXd &u) const
Computes the values of the potential reconstruction at the mesh vertices.
Definition hhospace.cpp:394
Eigen::MatrixXd potential
Definition hhospace.hpp:106
Family< MonomialScalarBasisCell > PolyBasisCellType
Definition hhospace.hpp:56
LocalOperators(const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation)
Definition hhospace.hpp:93
const LocalOperators & operators(size_t iT) const
Return operators for the cell of index iT.
Definition hhospace.hpp:161
Eigen::MatrixXd stabilisation
Definition hhospace.hpp:107
std::unique_ptr< PolyBasisEdgeType > Polyk
Definition hhospace.hpp:87
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
size_t K
Definition HHO_DiffAdvecReac.hpp:46
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
Definition ddr-klplate.hpp:27
static auto q
Definition ddrcore-test.hpp:14
Structure to store element bases.
Definition hhospace.hpp:71
Structure to store edge bases.
Definition hhospace.hpp:83
A structure to store local operators (gradient, potential, stabilisation)
Definition hhospace.hpp:92