HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
hho-uns.hpp
Go to the documentation of this file.
1// Implementation of the HHO scheme, with full gradient, for the Stokes equations -div(mu grad u)+grad p=f, div u = g.
2//
3// Author: Jerome Droniou (jerome.droniou@umontpellier.fr), Daniele Di Pietro (daniele.di-pietro@umontpellier.fr)
4//
5
6/*
7 *
8 * This implementation of HHO was developped following the principles described in
9 * Appendix B of the book
10 *
11 * The Hybrid High-Order Method for Polytopal Meshes: Design, Analysis, and Applications.
12 * D. A. Di Pietro and J. Droniou. Modeling, Simulation and Applications, vol. 19.
13 * Springer International Publishing, 2020, xxxi + 525p. doi: 10.1007/978-3-030-37203-3.
14 * url: https://hal.archives-ouvertes.fr/hal-02151813.
15 *
16 * If you use this code or part of it for a scientific publication, please cite the book
17 * above as a reference for the implementation.
18 *
19 */
20
21
22#ifndef HHOUNS_HPP
23#define HHOUNS_HPP
24
25#include <iostream>
26
27#include <boost/math/constants/constants.hpp>
28
29#include <Eigen/Sparse>
30#include <unsupported/Eigen/SparseExtra>
31
32#include <mesh.hpp>
33// #include <mesh_builder.hpp>
35#include <linearsolver.hpp>
36
37#include <vhhospace.hpp>
39#include <integralweight.hpp>
40
41
42namespace HArDCore2D
43{
44
45 //------------------------------------------------------------------------------
46 // Navier-Stokes model
47 //------------------------------------------------------------------------------
48
50
58 struct NavierStokes
59 {
60 typedef Eigen::SparseMatrix<double> SystemMatrixType;
61
62 typedef std::function<VectorRd(const VectorRd &)> MomentumForcingTermType;
63 typedef std::function<double(const VectorRd &)> CompressibilityForcingTermType;
64 typedef std::function<VectorRd(const VectorRd &)> VelocityType;
65 typedef std::function<MatrixRd(const VectorRd &)> VelocityGradientType;
66 typedef std::function<double(const VectorRd &)> PressureType;
67 typedef std::function<VectorRd(const VectorRd &)> PressureGradientType;
69
72 const VHHOSpace & vhho_space,
73 const GlobalDOFSpace & p_space,
74 const BoundaryConditions & BC,
75 bool use_threads,
76 bool stokes = false,
77 std::ostream & output = std::cout
78 );
79
82 const MomentumForcingTermType & f,
84 const ViscosityType & mu,
85 const Eigen::VectorXd & u0,
86 const Eigen::VectorXd & uold,
87 const double & time_step
88 );
89
91 inline size_t nloc_sc_u() const
92 {
93 return m_nloc_sc_u;
94 }
95
97 inline size_t nloc_sc_p() const
98 {
99 return m_nloc_sc_p;
100 }
101
103 inline size_t numSCDOFs_u() const
104 {
105 return m_vhhospace.mesh().n_cells() * m_nloc_sc_u;
106 }
107
109 inline size_t numSCDOFs_p() const
110 {
111 return m_vhhospace.mesh().n_cells() * m_nloc_sc_p;
112 }
113
115 inline size_t numSCDOFs() const
116 {
117 return numSCDOFs_u() + numSCDOFs_p();
118 }
119
121 inline size_t numDirDOFs() const
122 {
123 return m_BC.n_dir_edges()*m_vhhospace.numLocalDofsEdge();
124 }
125
127 inline size_t dimVelocity() const
128 {
129 return m_vhhospace.dimension();
130 }
131
133 inline size_t dimPressure() const
134 {
135 return m_pspace.dimension();
136 }
137
139 inline size_t numNonSCDOFs() const
140 {
141 return dimVelocity() + dimPressure() - numSCDOFs() + 1;
142 }
143
145 inline size_t sizeSystem() const
146 {
147 return numNonSCDOFs() - numDirDOFs();
148 }
149
151 inline const VHHOSpace & vhhospace() const
152 {
153 return m_vhhospace;
154 }
155
157 inline const GlobalDOFSpace & pspace() const
158 {
159 return m_pspace;
160 }
161
163 inline const Mesh & mesh() const
164 {
165 return m_vhhospace.mesh();
166 }
167
169 inline const SystemMatrixType & systemMatrix() const {
170 return m_A;
171 }
172
175 return m_A;
176 }
177
179 inline const Eigen::VectorXd & systemVector() const {
180 return m_b;
181 }
182
184 inline Eigen::VectorXd & systemVector() {
185 return m_b;
186 }
187
189 inline const double & stabilizationParameter() const {
190 return m_stab_par;
191 }
192
194 inline double & stabilizationParameter() {
195 return m_stab_par;
196 }
197
199 inline const SystemMatrixType & scMatrix() const {
200 return m_sc_A;
201 }
202
204 inline Eigen::VectorXd & scVector() {
205 return m_sc_b;
206 }
207
209 inline bool & isStokes() {
210 return m_stokes;
211 }
212
214 Eigen::VectorXd interpolate(
215 const VelocityType & u,
216 const PressureType & p,
217 const int doe_cell = -1,
218 const int doe_face = -1
219 ) const;
220
222 std::vector<double> computeEnergyNorms(
223 const std::vector<Eigen::VectorXd> & list_dofs
224 ) const;
225
227 Eigen::VectorXd pressureVertexValues(
228 const Eigen::VectorXd & p
229 ) const;
230
232 void writeToVtuFile(
233 const Eigen::VectorXd & uph,
234 const Eigen::VectorXd & upI,
235 const std::unique_ptr<Mesh> & mesh_ptr,
236 const std::string & solution_name,
237 const std::string & mesh_name,
238 const size_t & degree,
239 const int & time_it
240 );
241
243 Eigen::VectorXd imposingBC(
244 const VelocityType & u,
245 const Eigen::VectorXd & uph
246 ) const;
248 Eigen::VectorXd newtonRaphson(
250 const std::unique_ptr<Mesh> & mesh_ptr,
251 const size_t & newton_maxit,
252 const double & newton_tol,
253 const Eigen::VectorXd & u_newt,
254 const Eigen::VectorXd & un,
255 const MomentumForcingTermType & f,
257 const ViscosityType & mu,
258 const double & time_step
259 );
260
261 private:
263 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
264 _compute_local_contribution(
265 size_t iT,
266 const MomentumForcingTermType & f,
268 const ViscosityType & mu,
269 const Eigen::VectorXd & u0,
270 const Eigen::VectorXd & uold,
271 const double & time_step
272 );
273
275 LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
276
278 void _assemble_local_contribution(
279 size_t iT,
280 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
281 std::list<Eigen::Triplet<double> > & A1,
282 Eigen::VectorXd & b1,
283 std::list<Eigen::Triplet<double> > & A2,
284 Eigen::VectorXd & b2
285 );
286
288 const VHHOSpace & m_vhhospace;
289 const GlobalDOFSpace & m_pspace;
290
291 const BoundaryConditions & m_BC;
292 bool m_use_threads;
293 bool m_stokes;
294 std::ostream & m_output;
295 const size_t m_nloc_sc_u; // Number of statically condensed velocity DOFs in each cell
296 const size_t m_nloc_sc_p; // Number of statically condensed pressure DOFs in each cell
297 SystemMatrixType m_A; // Matrix and RHS of statically condensed system
298 Eigen::VectorXd m_b;
299 SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
300 Eigen::VectorXd m_sc_b;
301 double m_stab_par;
302 };
303
304} // end of namespace HArDCore3D
305
306#endif
The BoundaryConditions class provides definition of boundary conditions.
Definition BoundaryConditions.hpp:45
const size_t n_dir_edges() const
Returns the number of Dirichlet edges.
Definition BoundaryConditions.hpp:70
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
This structure is a wrapper to allow a given code to select various linear solvers.
Definition linearsolver.hpp:146
Class definition: polynomial bases and operators.
Definition vhhospace.hpp:52
Definition Mesh2D.hpp:26
Eigen::Vector2d VectorRd
Definition basis.hpp:55
Eigen::Matrix2d MatrixRd
Definition basis.hpp:54
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition localdofspace.hpp:61
size_t numLocalDofsEdge() const
Returns the number of local edge DOFs.
Definition localdofspace.hpp:45
const Mesh & mesh() const
Return a const reference to the mesh.
Definition vhhospace.hpp:140
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
std::string mesh_name
Definition HHO_DiffAdvecReac.hpp:44
std::size_t n_cells() const
number of cells in the mesh.
Definition Mesh2D.hpp:63
Definition ddr-klplate.hpp:27
Structure for weights (scalar, at the moment) in integral.
Definition integralweight.hpp:33
Structure to store information for, and perform, local static condensation.
Definition local_static_condensation.hpp:25
Class for the Navier-Stokes model.
Definition hho-ns.hpp:59
bool & isStokes()
Returns the Stokes status.
Definition hho-uns.hpp:209
const GlobalDOFSpace & pspace() const
Returns the pressure space.
Definition hho-uns.hpp:157
size_t nloc_sc_p() const
Returns the local number of pressure statically condensed DOFs.
Definition hho-uns.hpp:97
size_t numNonSCDOFs() const
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DO...
Definition hho-uns.hpp:139
const Mesh & mesh() const
Returns the mesh.
Definition hho-uns.hpp:163
size_t numSCDOFs_p() const
Returns the number of pressure statically condensed DOFs.
Definition hho-uns.hpp:109
NavierStokes(const VHHOSpace &vhho_space, const GlobalDOFSpace &p_space, const BoundaryConditions &BC, bool use_threads, bool stokes=false, std::ostream &output=std::cout)
Constructor.
size_t numSCDOFs_u() const
Returns the number of velocity statically condensed DOFs.
Definition hho-uns.hpp:103
std::function< VectorRd(const VectorRd &)> MomentumForcingTermType
Definition hho-uns.hpp:62
size_t nloc_sc_u() const
Returns the local number of velocity statically condensed DOFs.
Definition hho-uns.hpp:91
void writeToVtuFile(const Eigen::VectorXd &uph, const Eigen::VectorXd &upI, const std::unique_ptr< Mesh > &mesh_ptr, const std::string &solution_name, const std::string &mesh_name, const size_t &degree, const int &time_it)
Write solution to vut file.
Definition hho-uns.cpp:812
IntegralWeight ViscosityType
Definition hho-uns.hpp:68
std::function< VectorRd(const VectorRd &)> PressureGradientType
Definition hho-uns.hpp:67
double & stabilizationParameter()
Returns the stabilization parameter.
Definition hho-uns.hpp:194
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition hho-uns.hpp:169
const double & stabilizationParameter() const
Returns the stabilization parameter (scaling)
Definition hho-uns.hpp:189
size_t numDirDOFs() const
Returns the number of Dirichlet DOFs.
Definition hho-uns.hpp:121
size_t dimVelocity() const
Returns the dimension of velocity space.
Definition hho-uns.hpp:127
std::function< VectorRd(const VectorRd &)> VelocityType
Definition hho-uns.hpp:64
Eigen::SparseMatrix< double > SystemMatrixType
Definition hho-uns.hpp:60
void assembleLinearSystem(const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const ViscosityType &mu, const Eigen::VectorXd &u0, const Eigen::VectorXd &uold, const double &time_step)
Assemble the global system
Definition hho-uns.cpp:361
std::vector< double > computeEnergyNorms(const std::vector< Eigen::VectorXd > &list_dofs) const
Compute the discrete energy norm of a family of vectors representing the dofs.
size_t dimPressure() const
Returns the dimension of pressure space.
Definition hho-uns.hpp:133
Eigen::VectorXd interpolate(const VelocityType &u, const PressureType &p, const int doe_cell=-1, const int doe_face=-1) const
Interpolates velocity and pressure.
size_t sizeSystem() const
Returns the size of the final system with Lagrange multiplier, after application of SC and removal of...
Definition hho-uns.hpp:145
std::function< MatrixRd(const VectorRd &)> VelocityGradientType
Definition hho-uns.hpp:65
size_t numSCDOFs() const
Returns the number of statically condensed DOFs.
Definition hho-uns.hpp:115
Eigen::VectorXd newtonRaphson(LinearSolver< NavierStokes::SystemMatrixType > &solver, const std::unique_ptr< Mesh > &mesh_ptr, const size_t &newton_maxit, const double &newton_tol, const Eigen::VectorXd &u_newt, const Eigen::VectorXd &un, const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const ViscosityType &mu, const double &time_step)
newton method
Definition hho-uns.cpp:885
Eigen::VectorXd imposingBC(const VelocityType &u, const Eigen::VectorXd &uph) const
impose boundary conditions
Definition hho-uns.cpp:858
std::function< double(const VectorRd &)> PressureType
Definition hho-uns.hpp:66
const VHHOSpace & vhhospace() const
Returns the velocity space.
Definition hho-uns.hpp:151
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition hho-uns.hpp:179
std::function< double(const VectorRd &)> CompressibilityForcingTermType
Definition hho-uns.hpp:63
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition hho-uns.hpp:199
Eigen::VectorXd pressureVertexValues(const Eigen::VectorXd &p) const
Create vertex values for the pressure (from the element values), for plotting.
Definition hho-uns.cpp:786
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition hho-uns.hpp:204
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition hho-uns.hpp:174
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition hho-uns.hpp:184