HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
vhhospace.hpp
Go to the documentation of this file.
1 // Core data structures and methods required to implement the Hybrid High-Order in 3D, for vector-valued functions
2 //
3 // Provides:
4 // - Polynomial spaces on the element and faces
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 VHHOSPACE_HPP
27 #define VHHOSPACE_HPP
28 
29 #include <iostream>
30 #include <globaldofspace.hpp>
31 #include <basis.hpp>
33 
34 
35 namespace HArDCore3D
36 {
37 
44  //------------------------------------------------------------------------------
45 
46  // Type for selecting some cells (return 0 or 1 for each cell T)
47  typedef std::function<bool(const Cell &)> CellSelection;
48  static const CellSelection allcells = [](const Cell &)->bool {return true;};
49 
51  class VHHOSpace : public GlobalDOFSpace
52  {
53  public:
54  // Types for element bases
58 
59  // Types for face bases
62 
63  // Types for functions to interpolate
64  typedef std::function<VectorRd(const VectorRd &)> FunctionType;
65 
67 
71  struct CellBases
72  {
75 
76  std::unique_ptr<PolyBasisCellType> Polykpo;
77  std::unique_ptr<PolyBasisCellType> Polyk;
78  std::unique_ptr<PolydBasisCellType> Polykpod;
79  std::unique_ptr<PolydBasisCellType> Polykd;
80  std::unique_ptr<PolydxdBasisCellType> Polykdxd;
81  };
82 
84 
85  struct FaceBases
86  {
89 
90  std::unique_ptr<PolyBasisFaceType> Polyk;
91  std::unique_ptr<PolydBasisFaceType> Polykd;
92  };
93 
96  {
98  const Eigen::MatrixXd & _gradient,
99  const Eigen::MatrixXd & _potential,
100  const Eigen::MatrixXd & _stabilisation,
101  const Eigen::MatrixXd & _potential_div,
102  const Eigen::MatrixXd & _stabilisation_div
103  )
104  : gradient(_gradient),
105  potential(_potential),
106  stabilisation(_stabilisation),
107  potential_div(_potential_div),
108  stabilisation_div(_stabilisation_div)
109  {
110  // Do nothing
111  }
112 
113  Eigen::MatrixXd gradient;
114  Eigen::MatrixXd potential;
115  Eigen::MatrixXd stabilisation;
116  Eigen::MatrixXd potential_div;
117  Eigen::MatrixXd stabilisation_div;
118  };
119 
121  VHHOSpace(const Mesh & mesh, size_t K, const CellSelection & BoundaryStab, bool use_threads = true, std::ostream & output = std::cout);
122 
124  VHHOSpace(const Mesh & mesh, size_t K, bool use_threads = true, std::ostream & output = std::cout)
125  : VHHOSpace(mesh, K, allcells, use_threads, output) {};
126 
128  const Mesh & mesh() const
129  {
130  return m_mesh;
131  }
132 
134  const size_t & degree() const
135  {
136  return m_K;
137  }
138 
140  const CellSelection & boundaryStab() const
141  {
142  return m_boundary_stab;
143  }
144 
146  Eigen::VectorXd interpolate(
147  const FunctionType & q,
148  const int doe_cell = -1,
149  const int doe_face = -1
150  ) const;
151 
153  inline const CellBases & cellBases(size_t iT) const
154  {
155  // Make sure that the basis has been created
156  assert( m_cell_bases[iT] );
157  return *m_cell_bases[iT].get();
158  }
159 
161  inline const CellBases & cellBases(const Cell & T) const
162  {
163  return cellBases(T.global_index());
164  }
165 
167  inline const FaceBases & faceBases(size_t iF) const
168  {
169  // Make sure that the basis has been created
170  assert( m_face_bases[iF] );
171  return *m_face_bases[iF].get();
172  }
173 
175  inline const FaceBases & faceBases(const Face & F) const
176  {
177  return faceBases(F.global_index());
178  }
179 
181  inline const LocalOperators & operators(size_t iT) const
182  {
183  assert( m_operators[iT] );
184  return *m_operators[iT];
185  }
186 
188  inline const LocalOperators & operators(const Cell & T) const
189  {
190  return operators(T.global_index());
191  }
192 
194  std::vector<std::pair<double,double>> computeNorms(
195  const std::vector<Eigen::VectorXd> & list_dofs
196  ) const;
197 
199  std::vector<VectorRd> computeVertexValues(
200  const Eigen::VectorXd & u
201  ) const;
202 
203  private:
205  CellBases _construct_cell_bases(size_t iT);
206 
208  FaceBases _construct_face_bases(size_t iF);
209 
211  LocalOperators _compute_operators(size_t iT);
212 
213  // Pointer to the mesh
214  const Mesh & m_mesh;
215  // Degrees
216  const size_t m_K;
217  // Choice of boundary stabilisation
218  CellSelection m_boundary_stab;
219  // Parallel or not
220  bool m_use_threads;
221  // Output stream
222  std::ostream & m_output;
223 
224  // Cell bases
225  std::vector<std::unique_ptr<CellBases> > m_cell_bases;
226  // Face bases
227  std::vector<std::unique_ptr<FaceBases> > m_face_bases;
228 
229  // Local operators
230  std::vector<std::unique_ptr<LocalOperators> > m_operators;
231 
232  };
233 
234 } // end of namespace HArDCore3D
235 
236 #endif // VHHOSPACE_HPP
Family of functions expressed as linear combination of the functions of a given basis.
Definition: basis.hpp:388
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition: globaldofspace.hpp:16
Matrix family obtained from a scalar family.
Definition: basis.hpp:776
Vector family obtained by tensorization of a scalar family.
Definition: basis.hpp:609
Class definition: polynomial bases and operators.
Definition: vhhospace.hpp:52
Class to describe a mesh.
Definition: MeshND.hpp:17
std::unique_ptr< PolyBasisCellType > Polyk
Definition: vhhospace.hpp:77
VHHOSpace(const Mesh &mesh, size_t K, bool use_threads=true, std::ostream &output=std::cout)
Overloaded constructor when the selection of boundary stabilisation is not entered (all boundary face...
Definition: vhhospace.hpp:124
std::unique_ptr< PolyBasisFaceType > Polyk
Definition: vhhospace.hpp:90
TensorizedVectorFamily< PolyBasisCellType, dimspace > PolydBasisCellType
Definition: vhhospace.hpp:56
MatrixFamily< PolyBasisCellType, dimspace > PolydxdBasisCellType
Definition: vhhospace.hpp:57
Family< MonomialScalarBasisCell > PolyBasisCellType
Definition: vhhospace.hpp:55
Eigen::VectorXd interpolate(const FunctionType &q, const int doe_cell=-1, const int doe_face=-1) const
Interpolator of a continuous function.
Definition: vhhospace.cpp:157
Eigen::MatrixXd stabilisation_div
Definition: vhhospace.hpp:117
const FaceBases & faceBases(size_t iF) const
Return face bases for face iF.
Definition: vhhospace.hpp:167
std::unique_ptr< PolydBasisCellType > Polykd
Definition: vhhospace.hpp:79
Face GeometricSupport
Geometric support.
Definition: vhhospace.hpp:88
std::function< VectorRd(const VectorRd &)> FunctionType
Definition: vhhospace.hpp:64
const FaceBases & faceBases(const Face &F) const
Return cell bases for face F.
Definition: vhhospace.hpp:175
std::vector< VectorRd > computeVertexValues(const Eigen::VectorXd &u) const
Computes the values of the potential reconstruction at the mesh vertices.
Definition: vhhospace.cpp:526
const LocalOperators & operators(size_t iT) const
Return operators for the cell of index iT.
Definition: vhhospace.hpp:181
std::unique_ptr< PolyBasisCellType > Polykpo
Definition: vhhospace.hpp:76
Eigen::MatrixXd potential
Definition: vhhospace.hpp:114
std::unique_ptr< PolydBasisFaceType > Polykd
Definition: vhhospace.hpp:91
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: vhhospace.hpp:128
static const CellSelection allcells
Definition: vhhospace.hpp:48
std::unique_ptr< PolydxdBasisCellType > Polykdxd
Definition: vhhospace.hpp:80
const CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition: vhhospace.hpp:161
const size_t & degree() const
Return the polynomial degree (common face and elements)
Definition: vhhospace.hpp:134
Eigen::MatrixXd gradient
Definition: vhhospace.hpp:113
Cell GeometricSupport
Geometric support.
Definition: vhhospace.hpp:74
const CellSelection & boundaryStab() const
Return the function to select the cells with boundary stabilisation.
Definition: vhhospace.hpp:140
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: vhhospace.cpp:445
std::function< bool(const Cell &)> CellSelection
Definition: vhhospace.hpp:47
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition: vhhospace.hpp:153
LocalOperators(const Eigen::MatrixXd &_gradient, const Eigen::MatrixXd &_potential, const Eigen::MatrixXd &_stabilisation, const Eigen::MatrixXd &_potential_div, const Eigen::MatrixXd &_stabilisation_div)
Definition: vhhospace.hpp:97
const LocalOperators & operators(const Cell &T) const
Return cell operators for cell T.
Definition: vhhospace.hpp:188
Eigen::MatrixXd potential_div
Definition: vhhospace.hpp:116
Eigen::MatrixXd stabilisation
Definition: vhhospace.hpp:115
VHHOSpace(const Mesh &mesh, size_t K, const CellSelection &BoundaryStab, bool use_threads=true, std::ostream &output=std::cout)
Constructor (with function to select cells in which boundary faces are used in the stabilisation)
Definition: vhhospace.cpp:12
Family< TensorizedVectorFamily< PolyBasisFaceType, dimspace > > PolydBasisFaceType
Definition: vhhospace.hpp:61
Family< MonomialScalarBasisFace > PolyBasisFaceType
Definition: vhhospace.hpp:60
std::unique_ptr< PolydBasisCellType > Polykpod
Definition: vhhospace.hpp:78
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
size_t K
Definition: HHO_DiffAdvecReac.hpp:46
Definition: ddr-magnetostatics.hpp:40
MeshND::VectorRd< 2 > VectorRd
Definition: Mesh2D.hpp:14
MeshND::Face< 2 > Face
Definition: Mesh2D.hpp:12
MeshND::Cell< 2 > Cell
Definition: Mesh2D.hpp:13
Structure to store element bases.
Definition: vhhospace.hpp:72
Structure to store face bases.
Definition: vhhospace.hpp:86
A structure to store local operators (gradient, potential, stabilisation)
Definition: vhhospace.hpp:96