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#include <Eigen/Sparse>
6
8
10
11#include <unordered_set>
12
13
14#include "discrete-space.hpp"
16#include "hho-interpolate.hpp"
17#include "ns-solutions.hpp"
18
19namespace HArDCore2D {
20
21 namespace DSL {
22
23 struct HYPREParameters {
24 bool static_condensation = false; // true;
25 bool use_threads = true;
26 bool upwinding = true;
27 double stabilization_parameter= 1.0;
28 double viscosity = 1.;
30 std::ostream & output = std::cout;
31 };
32
33 class HYPRE {
34 public:
35 typedef Eigen::SparseMatrix<double> SystemMatrixType;
36
37 typedef std::function<VectorRd(const VectorRd &)> MomentumForcingTermType;
38 typedef std::function<VectorRd(const VectorRd &)> MagneticForcingTermType;
39 typedef std::function<double(const VectorRd &)> CompressibilityForcingTermType;
40
41 typedef std::function<VectorRd(const VectorRd &)> VelocityType;
42 typedef std::function<MatrixRd(const VectorRd &)> VelocityGradientType;
43
44 typedef std::function<double(const VectorRd &)> PressureType;
45 typedef std::function<VectorRd(const VectorRd &)> PressureGradientType;
46
47
49 const Mesh & mesh,
50 size_t degree,
51 const BoundaryConditions & bc,
52 const HYPREParameters & parameters = {}
53 );
54 virtual ~HYPRE() {}
55
57 virtual Eigen::VectorXd interpolate(
58 const VelocityType & u,
59 const PressureType & p
60 ) const = 0;
61
63 template<typename F>
64 Eigen::VectorXd pkpoPkPkPkInterpolate(
65 const F & v,
66 const InterpolateParameters & parameters = {}
67 ) const;
68
70 size_t degree() const
71 {
72 return m_degree;
73 }
74
77 {
79 }
80
82 const std::vector<Eigen::MatrixXd> potential() const
83 {
84 return m_potential;
85 }
86
88 const std::vector<Eigen::MatrixXd> stabilization() const
89 {
90 return m_stabilization;
91 }
92
94 const Eigen::MatrixXd & potential(size_t iT) const
95 {
96 return m_potential[iT];
97 }
98
101 {
102 return m_nloc_sc_u;
103 }
104
107 {
108 return m_nloc_sc_p;
109 }
110
113 {
114 return m_velocity_space->mesh().n_cells() * m_nloc_sc_u;
115 }
116
119 {
120 return m_velocity_space->mesh().n_cells() * m_nloc_sc_p;
121 }
122
128
130 size_t numDirichletDofs() const
131 {
132 return 2*(m_bc.n_dir_edges() * m_velocity_space->numberOfLocalEdgeDofs());
133 }
134
136 size_t dimVelocity() const
137 {
138 return m_velocity_space->dimension();
139 }
140
142 size_t dimPressure() const
143 {
144 return m_pressure_space->dimension();
145 }
146
149 {
150 return 2*dimVelocity() + 2*dimPressure() + 2 - numStaticallyCondensedDofs();
151 }
152
154 size_t sizeSystem() const
155 {
157 }
158
164
170
171
173 inline const DiscreteSpace & velocitySpace() const
174 {
175 return *m_velocity_space;
176 }
177
179 inline const DiscreteSpace & pressureSpace() const
180 {
181 return *m_pressure_space;
182 }
183
185 inline const Mesh & mesh() const
186 {
187 return m_velocity_space->mesh();
188 }
189
191 inline const SystemMatrixType & systemMatrix() const {
192 return m_A;
193 }
194
197 return m_A;
198 }
199
201 inline const Eigen::VectorXd & systemVector() const {
202 return m_b;
203 }
204
206 inline Eigen::VectorXd & systemVector() {
207 return m_b;
208 }
209
212 return m_sc_A;
213 }
214
216 inline const Eigen::VectorXd & staticallyCondensedVector() const {
217 return m_sc_b;
218 }
219
221 inline const BoundaryConditions & boundaryConditions() const {
222 return m_bc;
223 }
224
228 const Eigen::VectorXd & uph_old,
229 const Eigen::VectorXd & uph_prev,
230 const Eigen::VectorXd & brh_old,
231 const Eigen::VectorXd & brh_prev,
232 const double time_step,
233 const double current_time);
234
235
236 //This should be the stationary case assembleResidual case!
237
240 const Eigen::VectorXd & uph_old,
241 const Eigen::VectorXd & brh_old) {
242 Eigen::VectorXd uph_prev = Eigen::VectorXd::Zero(uph_old.size());
243 Eigen::VectorXd brh_prev = Eigen::VectorXd::Zero(brh_old.size());
244 assembleResidualSystem(isolution,
245 uph_old,
246 uph_prev,
247 brh_old,
248 brh_prev,
249 std::numeric_limits<double>::infinity(),
250 0.);
251 };
252
254 Eigen::VectorXd reconstructSolution(const Eigen::VectorXd & uhp_solsystem) const;
255
256
258 virtual double convective_stabilization_parameter(size_t iT,
260 const Eigen::VectorXd & uT_old,
261 const Eigen::VectorXd & bT_old,
262 const size_t cell_degree,
263 const double viscosity,
264 const double magnetic_diffusivity,
265 const double current_time) const;
266
267 virtual double magnetic_convective_stabilization_parameter(size_t iT,
269 const Eigen::VectorXd & uT_old,
270 const Eigen::VectorXd & bT_old,
271 const size_t cell_degree,
272 const double viscosity,
273 const double magnetic_diffusivity,
274 const double current_time) const;
275
276
277 protected:
279
280 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _time_deriv_term(size_t iT,
281 const Eigen::VectorXd & uT_old,
282 const Eigen::VectorXd & uT_prev,
283 const double time_step); // = std::numeric_limits<double>::infinity()
284
285
286 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _viscous_term(size_t iT, const Eigen::VectorXd & uT_old);
287
288 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _magnetic_diffusive_term(size_t iT, const Eigen::VectorXd & bT_old);
289
290
291 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd> _velocity_pressure_coupling(
292 size_t iT,
294 const Eigen::VectorXd & uT_old,
295 const Eigen::VectorXd & pT_old,
296 const double & current_time = 0.
297 );
298
299 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd> _magnetic_field_magnetic_pressure_coupling(
300 size_t iT,
302 const Eigen::VectorXd & bT_old,
303 const Eigen::VectorXd & rT_old,
304 const double & current_time = 0.
305 );
306
307
308 virtual Eigen::MatrixXd _forcing_terms(
309 size_t iT,
311 const double current_time = 0.
312 );
313
314 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
315 _linearized_convective_term(size_t iT, const Eigen::VectorXd & uT_old);
316
317 virtual std::tuple<Eigen::MatrixXd, Eigen::MatrixXd, Eigen::VectorXd>
318 _linearized_convective_mixed_term(size_t iT, const Eigen::VectorXd & wT_old, const Eigen::VectorXd & vT_old);
319
320 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
322 size_t iT,
323 const Eigen::VectorXd & uT_old,
324 const Eigen::VectorXd & bT_old,
326 const double current_time);
327
328 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
330 size_t iT,
331 const Eigen::VectorXd & uT_old,
332 const Eigen::VectorXd & bT_old,
334 const double current_time);
335
336 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
338 size_t iT,
340 const Eigen::VectorXd & uph_old,
341 const Eigen::VectorXd & uph_prev,
342 const Eigen::VectorXd & brh_old,
343 const Eigen::VectorXd & brh_prev,
344 const double time_step,
345 const double current_time);
346
347 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
349 size_t iT,
351 const Eigen::VectorXd & uph_old,
352 const Eigen::VectorXd & brh_old) {
353 Eigen::VectorXd uph_prev = uph_old;
354 Eigen::VectorXd brh_prev = brh_old;
355
357 isolution,
358 uph_old,
359 uph_prev,
360 brh_old,
361 brh_prev,
362 std::numeric_limits<double>::infinity(),
363 0.);
364
365 }
366
368
370 size_t iT,
371 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
372 std::list<Eigen::Triplet<double> > & A1,
373 Eigen::VectorXd & b1,
374 std::list<Eigen::Triplet<double> > & A2,
375 Eigen::VectorXd & b2
376 );
377
378 size_t m_degree;
379 const BoundaryConditions & m_bc;
381 bool m_use_threads;
383
385 double m_viscosity;
387
388 std::ostream & m_output;
389
390 size_t m_nloc_sc_u;
391 size_t m_nloc_sc_p;
392
393 std::vector<Eigen::MatrixXd> m_discrete_l2;
394 std::vector<Eigen::MatrixXd> m_gradient;
395 std::vector<Eigen::MatrixXd> m_gradient_rhs;
396 std::vector<Eigen::MatrixXd> m_potential;
397 std::vector<Eigen::MatrixXd> m_stabilization;
398
399 std::unique_ptr<DiscreteSpace> m_velocity_space;
400 std::unique_ptr<DiscreteSpace> m_pressure_space;
401
402
404 Eigen::VectorXd m_b;
405
407 Eigen::VectorXd m_sc_b;
408
409 std::vector<int> m_nsc_to_dof_map;
410 std::vector<int> m_sc_to_dof_map;
411
412 };
413
414 } // namespace DSL
415
416} // namespace HArDCore2D
417
418#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:94
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
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _magnetic_field_magnetic_pressure_coupling(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &bT_old, const Eigen::VectorXd &rT_old, const double &current_time=0.)
Definition hypre.cpp:251
double m_magnetic_diffusivity
Definition hypre.hpp:386
const std::vector< Eigen::MatrixXd > potential() const
Returns the vector of potential reconstructions.
Definition hypre.hpp:82
bool m_static_condensation
Definition hypre.hpp:347
std::function< VectorRd(const VectorRd &)> VelocityType
Definition hypre.hpp:41
const SystemMatrixType & staticallyCondensedMatrix() const
Returns the static condensation recovery operator.
Definition hypre.hpp:211
virtual double magnetic_convective_stabilization_parameter(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uT_old, const Eigen::VectorXd &bT_old, const size_t cell_degree, const double viscosity, const double magnetic_diffusivity, const double current_time) const
Definition hypre.cpp:1342
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 &brh_old)
Definition hypre.hpp:348
std::vector< Eigen::MatrixXd > m_gradient_rhs
Definition hypre.hpp:360
size_t degree() const
Returns the degree.
Definition hypre.hpp:64
Eigen::VectorXd reconstructSolution(const Eigen::VectorXd &uhp_solsystem) const
Reconstruct the solution from the vector of non-statically-condensed DOFs.
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:45
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition hypre.hpp:191
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _magnetic_diffusive_term(size_t iT, const Eigen::VectorXd &bT_old)
Definition hypre.cpp:177
const std::vector< Eigen::MatrixXd > stabilization() const
Returns the vector of stabilization matrices.
Definition hypre.hpp:88
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:42
const DiscreteSpace & pressureSpace() const
Returns the pressure space.
Definition hypre.hpp:179
size_t numNonStaticallyCondensedDofsPressure() const
Returns the number of non statically condensed DOFs for the pressure.
Definition hypre.hpp:166
size_t numNonStaticallyCondensedDofsVelocity() const
Returns the number of non statically condensed DOFs for the velocity.
Definition hypre.hpp:160
virtual Eigen::MatrixXd _forcing_terms(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const double current_time=0.)
double stabilizationParameter() const
Returns the stabilization parameter.
Definition hypre.hpp:76
size_t numStaticallyCondensedDofsVelocity() const
Returns the number of statically condensed DOFs for the velocity.
Definition hypre.hpp:112
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _viscous_term(size_t iT, const Eigen::VectorXd &uT_old)
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:124
size_t numLocalStaticallyCondensedDofsVelocity() const
Returns the local number of velocity statically condensed DOFs.
Definition hypre.hpp:100
std::unique_ptr< DiscreteSpace > m_pressure_space
Definition hypre.hpp:367
bool m_upwinding
Definition hypre.hpp:382
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:173
Eigen::SparseMatrix< double > SystemMatrixType
Definition hypre.hpp:35
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:201
const BoundaryConditions & boundaryConditions() const
Returns the boundary condition handler.
Definition hypre.hpp:221
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _linearized_convective_term(size_t iT, const Eigen::VectorXd &uT_old)
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:37
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:216
size_t dimVelocity() const
Returns the dimension of velocity space.
Definition hypre.hpp:136
std::unique_ptr< DiscreteSpace > m_velocity_space
Definition hypre.hpp:366
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition hypre.hpp:196
SystemMatrixType m_A
Definition hypre.hpp:369
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)
size_t m_nloc_sc_p
Definition hypre.hpp:356
std::function< double(const VectorRd &)> CompressibilityForcingTermType
Definition hypre.hpp:39
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:154
std::vector< Eigen::MatrixXd > m_stabilization
Definition hypre.hpp:364
size_t dimPressure() const
Returns the dimension of pressure space.
Definition hypre.hpp:142
size_t numLocalStaticallyCondensedDofsPressure() const
Returns the local number of pressure statically condensed DOFs.
Definition hypre.hpp:106
std::function< VectorRd(const VectorRd &)> MagneticForcingTermType
Definition hypre.hpp:38
size_t m_degree
Definition hypre.hpp:345
virtual std::tuple< Eigen::MatrixXd, Eigen::MatrixXd, Eigen::VectorXd > _linearized_convective_mixed_term(size_t iT, const Eigen::VectorXd &wT_old, const Eigen::VectorXd &vT_old)
Definition hypre.cpp:483
size_t numNonStaticallyCondensedDofs() const
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DO...
Definition hypre.hpp:148
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
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)
double m_viscosity
Definition hypre.hpp:350
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.)
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _magnetic_convective_stabilization_term(size_t iT, const Eigen::VectorXd &uT_old, const Eigen::VectorXd &bT_old, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const double current_time)
Definition hypre.cpp:693
virtual ~HYPRE()
Definition hypre.hpp:54
std::function< double(const VectorRd &)> PressureType
Definition hypre.hpp:44
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
void _construct_local_operators(size_t iT)
size_t numStaticallyCondensedDofsPressure() const
Returns the number of statically condensed DOFs for the pressure.
Definition hypre.hpp:118
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition hypre.hpp:206
virtual Eigen::VectorXd interpolate(const VelocityType &u, const PressureType &p) const =0
Interpolates velocity and pressure.
HYPRE(const Mesh &mesh, size_t degree, const BoundaryConditions &bc, const HYPREParameters &parameters={})
void assembleResidualSystem(const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const Eigen::VectorXd &uph_old, const Eigen::VectorXd &brh_old)
Definition hypre.hpp:238
LocalStaticCondensation _compute_static_condensation(const size_t &iT)
size_t numDirichletDofs() const
Returns the number of Dirichlet DOFs.
Definition hypre.hpp:130
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
bool static_condensation
Definition hypre.hpp:22
double stabilization_parameter
Definition hypre.hpp:24
std::ostream & output
Definition hypre.hpp:28
bool use_threads
Definition hypre.hpp:23
bool upwinding
Definition hypre.hpp:26
double viscosity
Definition hypre.hpp:25
double magnetic_diffusivity
Definition hypre.hpp:29
Definition hho-interpolate.hpp:15
Structure to store information for, and perform, local static condensation.
Definition local_static_condensation.hpp:25