HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
hypre.hpp
Go to the documentation of this file.
1#ifndef HYPRE_HPP
2#define HYPRE_HPP
3
4#include <ostream>
5
6#include <Eigen/Sparse>
7
9
11
12#include "discrete-space.hpp"
14#include "hho-interpolate.hpp"
15#include "ns-solutions.hpp"
16
17namespace HArDCore2D {
18
19 namespace DSL {
20
23 bool use_threads = true;
25 double viscosity = 1.;
26 double flow_index = 2.;
27 double degen_param = 0.;
28 std::ostream & output = std::cout;
29 };
30
31 class HYPRE {
32 public:
33 typedef Eigen::SparseMatrix<double> SystemMatrixType;
34
35 typedef std::function<VectorRd(const VectorRd &)> MomentumForcingTermType;
36 typedef std::function<double(const VectorRd &)> CompressibilityForcingTermType;
37 typedef std::function<VectorRd(const VectorRd &)> VelocityType;
38 typedef std::function<MatrixRd(const VectorRd &)> VelocityGradientType;
39 typedef std::function<double(const VectorRd &)> PressureType;
40 typedef std::function<VectorRd(const VectorRd &)> PressureGradientType;
41
42 HYPRE(
43 const Mesh & mesh,
44 size_t degree,
45 const BoundaryConditions & bc,
46 const HYPREParameters & parameters = {}
47 );
48 virtual ~HYPRE() {}
49
51 virtual Eigen::VectorXd interpolate(
52 const VelocityType & u,
53 const PressureType & p
54 ) const = 0;
55
57 template<typename F>
58 Eigen::VectorXd pkpoPkPkPkInterpolate(
59 const F & v,
60 const InterpolateParameters & parameters = {}
61 ) const;
62
64 size_t degree() const
65 {
66 return m_degree;
67 }
68
71 {
73 }
74
76 const std::vector<Eigen::MatrixXd> potential() const
77 {
78 return m_potential;
79 }
80
82 const std::vector<Eigen::MatrixXd> stabilization() const
83 {
84 return m_stabilization;
85 }
86
88 const Eigen::MatrixXd & potential(size_t iT) const
89 {
90 return m_potential[iT];
91 }
92
95 {
96 return m_nloc_sc_u;
97 }
98
101 {
102 return m_nloc_sc_p;
103 }
104
107 {
108 return m_velocity_space->mesh().n_cells() * m_nloc_sc_u;
109 }
110
113 {
114 return m_velocity_space->mesh().n_cells() * m_nloc_sc_p;
115 }
116
122
124 size_t numDirichletDofs() const
125 {
126 return m_bc.n_dir_edges() * m_velocity_space->numberOfLocalEdgeDofs();
127 }
128
130 size_t dimVelocity() const
131 {
132 return m_velocity_space->dimension();
133 }
134
136 size_t dimPressure() const
137 {
138 return m_pressure_space->dimension();
139 }
140
143 {
145 }
146
148 size_t sizeSystem() const
149 {
151 }
152
154 inline const DiscreteSpace & velocitySpace() const
155 {
156 return *m_velocity_space;
157 }
158
160 inline const DiscreteSpace & pressureSpace() const
161 {
162 return *m_pressure_space;
163 }
164
166 inline const Mesh & mesh() const
167 {
168 return m_velocity_space->mesh();
169 }
170
172 inline const SystemMatrixType & systemMatrix() const {
173 return m_A;
174 }
175
178 return m_A;
179 }
180
182 inline const Eigen::VectorXd & systemVector() const {
183 return m_b;
184 }
185
187 inline Eigen::VectorXd & systemVector() {
188 return m_b;
189 }
190
193 return m_sc_A;
194 }
195
197 inline const Eigen::VectorXd & staticallyCondensedVector() const {
198 return m_sc_b;
199 }
200
202 inline const BoundaryConditions & boundaryConditions() const {
203 return m_bc;
204 }
205
209 const Eigen::VectorXd & uph_old,
210 const Eigen::VectorXd & uph_prev,
211 const double time_step,
212 const double current_time);
213
216 const Eigen::VectorXd & uph_old) {
217 Eigen::VectorXd uph_prev = Eigen::VectorXd::Zero(uph_old.size());
218 assembleResidualSystem(isolution,
219 uph_old,
220 uph_prev,
221 std::numeric_limits<double>::infinity(),
222 0.);
223 };
224
226 Eigen::VectorXd reconstructSolution(const Eigen::VectorXd & uhp_solsystem) const;
227
228 // double convective_stabilization_parameter(size_t iT,
229 // const HArDCore2D::NavierStokesSolutions::IExactSolution * isolution,
230 // const size_t cell_degree,
231 // const double viscosity,
232 // const Eigen::VectorXd uT_old,
233 // const double current_time = 0.) const {
234 // // Construct quadrature scheme on element
235 // const Cell & T = *m_velocity_space->mesh().cell(iT);
236
237 // // Local convective stabilization bilinear form assumed to be given by
238 // // j_T(w_T, v_T) = Beta_T Sum_{F in F_T} Int_{F} (w_F - w_T).(v_F - v_T)
239 // // where Beta_T = max(c_s, ||u_T||_{L^infty(T)}). We therefore need an
240 // // estimate of ||u_T||_{L^infty(T)} which we obtain by taking the largest
241 // // l^infty norm of u_T on the quadrature points on T.
242 // QuadratureRule quad = generate_quadrature_rule(T, 3 * cell_degree);
243 // const auto & cell_dofs = m_velocity_space->get<PolynCellType>("CellDOFs")[iT];
244 // auto cell_dofs_quad = evaluate_quad<Function>::compute(*cell_dofs, quad);
245 // double beta_T = 1.e-4; // Safeguard parameter
246 // for (size_t iqn = 0; iqn < quad.size(); iqn++) {
247 // // Stabilisation based on max-norm of discrete solution
248 // VectorRd uT_old_iqn = VectorRd::Zero();
249 // for (size_t i = 0; i < cell_dofs->dimension(); i++) {
250 // double uT_old_i = uT_old(m_velocity_space->localOffset(T) + i);
251 // uT_old_iqn += uT_old_i * cell_dofs_quad[i][iqn];
252 // } // for i
253 // // Stabilisation based on max-norm of analytic solution
254 // // Eigen::VectorXd uT_old_iqn = isolution->velocity(quad[iqn].vector(), current_time);
255 // beta_T = std::max(beta_T, uT_old_iqn.lpNorm<Eigen::Infinity>());
256 // } // for iqn
257
258 // // ALTERNATIVE APPROACH (based on maximum absolute value of u.n on the boundary of T)
259 // // double beta_T = 1.e-4; // Safeguard parameter
260 // // for (size_t iE = 0; iE < T.n_edges(); iE++) {
261 // // const Edge & E = *T.edge(iE);
262 // // auto nTE = T.edge_normal(iE);
263 // // QuadratureRule quad_E = generate_quadrature_rule(E, 3 * cell_degree);
264 // // for (size_t iqn = 0; iqn < quad_E.size(); iqn++) {
265 // // Eigen::VectorXd uE_old_iqn = isolution->velocity(quad_E[iqn].vector());
266 // // beta_T = std::max(beta_T, std::abs(uE_old_iqn.dot(nTE)));
267 // // }
268 // // }
269
270 // return beta_T;
271 // };
272
274 virtual double convective_stabilization_parameter(size_t iT,
276 const size_t cell_degree,
277 const double viscosity,
278 const Eigen::VectorXd uT_old,
279 const double current_time = 0.) const;
280
281 protected:
282 void _construct_local_operators(size_t iT);
283
284 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _time_deriv_term(size_t iT,
285 const Eigen::VectorXd & uT_old,
286 const Eigen::VectorXd & uT_prev,
287 const double time_step); // = std::numeric_limits<double>::infinity()
288 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _viscous_term(size_t iT, const Eigen::VectorXd & uT_old);
289 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd> _nonlinear_viscous_term(size_t iT, const Eigen::VectorXd & uT_old);
290 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd> _velocity_pressure_coupling(
291 size_t iT,
293 const Eigen::VectorXd & uT_old,
294 const Eigen::VectorXd & pT_old,
295 const double & current_time = 0.
296 );
297 virtual Eigen::MatrixXd _forcing_terms(
298 size_t iT,
300 const double current_time = 0.
301 );
302 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
303 _linearized_convective_term(size_t iT, const Eigen::VectorXd & uT_old);
304 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
306 const Eigen::VectorXd & uT_old,
308 const double current_time = 0.);
309
310 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
312 size_t iT,
314 const Eigen::VectorXd & uph_old,
315 const Eigen::VectorXd & uph_prev,
316 const double time_step,
317 const double current_time);
318 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
320 size_t iT,
322 const Eigen::VectorXd & uph_old) {
323 Eigen::VectorXd uph_prev = uph_old;
325 isolution,
326 uph_old,
327 uph_prev,
328 std::numeric_limits<double>::infinity(),
329 0.);
330
331 }
332
334
336 size_t iT,
337 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
338 std::list<Eigen::Triplet<double> > & A1,
339 Eigen::VectorXd & b1,
340 std::list<Eigen::Triplet<double> > & A2,
341 Eigen::VectorXd & b2
342 );
343
344
345 size_t m_degree;
353 std::ostream & m_output;
354
357
358 std::vector<Eigen::MatrixXd> m_discrete_l2;
359 std::vector<Eigen::MatrixXd> m_gradient;
360 std::vector<Eigen::MatrixXd> m_gradient_rhs;
361 std::vector<Eigen::MatrixXd> m_symgradient;
362 std::vector<Eigen::MatrixXd> m_symgradient_rhs;
363 std::vector<Eigen::MatrixXd> m_potential;
364 std::vector<Eigen::MatrixXd> m_stabilization;
365
366 std::unique_ptr<DiscreteSpace> m_velocity_space;
367 std::unique_ptr<DiscreteSpace> m_pressure_space;
368
370 Eigen::VectorXd m_b;
371
373 Eigen::VectorXd m_sc_b;
374
375 std::vector<int> m_nsc_to_dof_map;
376 std::vector<int> m_sc_to_dof_map;
377 };
378
379 } // namespace DSL
380
381} // namespace HArDCore2D
382
383#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
Definition discrete-space.hpp:20
Definition hypre.hpp:31
const Eigen::MatrixXd & potential(size_t iT) const
Returns the potential reconstruction.
Definition hypre.hpp:88
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _compute_local_contribution(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uph_old, const Eigen::VectorXd &uph_prev, const double time_step, const double current_time)
Definition hypre.cpp:834
Eigen::VectorXd m_b
Definition hypre.hpp:370
const std::vector< Eigen::MatrixXd > potential() const
Returns the vector of potential reconstructions.
Definition hypre.hpp:76
virtual Eigen::MatrixXd _forcing_terms(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const double current_time=0.)
Definition hypre.cpp:616
bool m_static_condensation
Definition hypre.hpp:347
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _linearized_convective_term(size_t iT, const Eigen::VectorXd &uT_old)
Definition hypre.cpp:657
LocalStaticCondensation _compute_static_condensation(const size_t &iT)
Definition hypre.cpp:921
std::function< VectorRd(const VectorRd &)> VelocityType
Definition hypre.hpp:37
const SystemMatrixType & staticallyCondensedMatrix() const
Returns the static condensation recovery operator.
Definition hypre.hpp:192
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _velocity_pressure_coupling(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uT_old, const Eigen::VectorXd &pT_old, const double &current_time=0.)
Definition hypre.cpp:555
std::vector< Eigen::MatrixXd > m_gradient_rhs
Definition hypre.hpp:360
size_t degree() const
Returns the degree.
Definition hypre.hpp:64
virtual double convective_stabilization_parameter(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const size_t cell_degree, const double viscosity, const Eigen::VectorXd uT_old, const double current_time=0.) const
Evaluates convective stabilisation parameter.
Definition hypre.cpp:1155
std::function< VectorRd(const VectorRd &)> PressureGradientType
Definition hypre.hpp:40
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition hypre.hpp:172
const std::vector< Eigen::MatrixXd > stabilization() const
Returns the vector of stabilization matrices.
Definition hypre.hpp:82
const Mesh & mesh() const
Returns the mesh.
Definition hypre.hpp:166
std::ostream & m_output
Definition hypre.hpp:353
std::function< MatrixRd(const VectorRd &)> VelocityGradientType
Definition hypre.hpp:38
const DiscreteSpace & pressureSpace() const
Returns the pressure space.
Definition hypre.hpp:160
double stabilizationParameter() const
Returns the stabilization parameter.
Definition hypre.hpp:70
std::vector< Eigen::MatrixXd > m_symgradient
Definition hypre.hpp:361
size_t numStaticallyCondensedDofsVelocity() const
Returns the number of statically condensed DOFs for the velocity.
Definition hypre.hpp:106
double m_stabilization_parameter
Definition hypre.hpp:349
Eigen::VectorXd m_sc_b
Definition hypre.hpp:373
size_t m_nloc_sc_u
Definition hypre.hpp:355
size_t numStaticallyCondensedDofs() const
Returns the number of statically condensed DOFs.
Definition hypre.hpp:118
size_t numLocalStaticallyCondensedDofsVelocity() const
Returns the local number of velocity statically condensed DOFs.
Definition hypre.hpp:94
std::unique_ptr< DiscreteSpace > m_pressure_space
Definition hypre.hpp:367
void _construct_local_operators(size_t iT)
Definition hypre.cpp:42
void assembleResidualSystem(const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uph_old, const Eigen::VectorXd &uph_prev, const double time_step, const double current_time)
Assemble the global system.
Definition hypre.cpp:1042
const DiscreteSpace & velocitySpace() const
Returns the velocity space.
Definition hypre.hpp:154
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _time_deriv_term(size_t iT, const Eigen::VectorXd &uT_old, const Eigen::VectorXd &uT_prev, const double time_step)
Definition hypre.cpp:253
Eigen::SparseMatrix< double > SystemMatrixType
Definition hypre.hpp:33
std::vector< Eigen::MatrixXd > m_discrete_l2
Definition hypre.hpp:358
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition hypre.hpp:182
const BoundaryConditions & boundaryConditions() const
Returns the boundary condition handler.
Definition hypre.hpp:202
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _convective_stabilization_term(size_t iT, const Eigen::VectorXd &uT_old, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const double current_time=0.)
Definition hypre.cpp:775
std::function< VectorRd(const VectorRd &)> MomentumForcingTermType
Definition hypre.hpp:35
std::vector< Eigen::MatrixXd > m_symgradient_rhs
Definition hypre.hpp:362
void _assemble_local_contribution(size_t iT, const std::pair< Eigen::MatrixXd, Eigen::VectorXd > &lsT, std::list< Eigen::Triplet< double > > &A1, Eigen::VectorXd &b1, std::list< Eigen::Triplet< double > > &A2, Eigen::VectorXd &b2)
Definition hypre.cpp:1095
double m_degen_param
Definition hypre.hpp:352
std::vector< int > m_nsc_to_dof_map
Definition hypre.hpp:375
const Eigen::VectorXd & staticallyCondensedVector() const
Returns the static condensation rhs.
Definition hypre.hpp:197
size_t dimVelocity() const
Returns the dimension of velocity space.
Definition hypre.hpp:130
std::unique_ptr< DiscreteSpace > m_velocity_space
Definition hypre.hpp:366
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition hypre.hpp:177
SystemMatrixType m_A
Definition hypre.hpp:369
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _compute_local_contribution(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uph_old)
Definition hypre.hpp:319
size_t m_nloc_sc_p
Definition hypre.hpp:356
std::function< double(const VectorRd &)> CompressibilityForcingTermType
Definition hypre.hpp:36
std::vector< Eigen::MatrixXd > m_potential
Definition hypre.hpp:363
size_t sizeSystem() const
Returns the size of the final system with Lagrange multiplier, after application of SC and removal of...
Definition hypre.hpp:148
std::vector< Eigen::MatrixXd > m_stabilization
Definition hypre.hpp:364
size_t dimPressure() const
Returns the dimension of pressure space.
Definition hypre.hpp:136
size_t numLocalStaticallyCondensedDofsPressure() const
Returns the local number of pressure statically condensed DOFs.
Definition hypre.hpp:100
size_t m_degree
Definition hypre.hpp:345
size_t numNonStaticallyCondensedDofs() const
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DO...
Definition hypre.hpp:142
Eigen::VectorXd reconstructSolution(const Eigen::VectorXd &uhp_solsystem) const
Reconstruct the solution from the vector of non-statically-condensed DOFs.
Definition hypre.cpp:1133
const BoundaryConditions & m_bc
Definition hypre.hpp:346
Eigen::VectorXd pkpoPkPkPkInterpolate(const F &v, const InterpolateParameters &parameters={}) const
Compute the interpolate for the Pk+1-Pk-Pk-Pk method.
SystemMatrixType m_sc_A
Definition hypre.hpp:372
double m_viscosity
Definition hypre.hpp:350
double m_flow_index
Definition hypre.hpp:351
virtual ~HYPRE()
Definition hypre.hpp:48
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _nonlinear_viscous_term(size_t iT, const Eigen::VectorXd &uT_old)
Definition hypre.cpp:274
std::function< double(const VectorRd &)> PressureType
Definition hypre.hpp:39
bool m_use_threads
Definition hypre.hpp:348
std::vector< int > m_sc_to_dof_map
Definition hypre.hpp:376
std::vector< Eigen::MatrixXd > m_gradient
Definition hypre.hpp:359
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _viscous_term(size_t iT, const Eigen::VectorXd &uT_old)
Definition hypre.cpp:264
size_t numStaticallyCondensedDofsPressure() const
Returns the number of statically condensed DOFs for the pressure.
Definition hypre.hpp:112
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition hypre.hpp:187
virtual Eigen::VectorXd interpolate(const VelocityType &u, const PressureType &p) const =0
Interpolates velocity and pressure.
size_t numDirichletDofs() const
Returns the number of Dirichlet DOFs.
Definition hypre.hpp:124
void assembleResidualSystem(const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uph_old)
Definition hypre.hpp:214
Definition Mesh2D.hpp:26
Eigen::Vector2d VectorRd
Definition basis.hpp:55
Eigen::Matrix2d MatrixRd
Definition basis.hpp:54
Definition ddr-klplate.hpp:27
static auto v
Definition ddrcore-test.hpp:32
Definition hypre.hpp:21
double flow_index
Definition hypre.hpp:26
bool static_condensation
Definition hypre.hpp:22
double degen_param
Definition hypre.hpp:27
double stabilization_parameter
Definition hypre.hpp:24
std::ostream & output
Definition hypre.hpp:28
bool use_threads
Definition hypre.hpp:23
double viscosity
Definition hypre.hpp:25
Definition hho-interpolate.hpp:15
Structure to store information for, and perform, local static condensation.
Definition local_static_condensation.hpp:25