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 double stabilization_parameter= true;
27 double viscosity = 1.;
29 std::ostream & output = std::cout;
30 };
31
32 class HYPRE {
33 public:
34 typedef Eigen::SparseMatrix<double> SystemMatrixType;
35
36 typedef std::function<VectorRd(const VectorRd &)> MomentumForcingTermType;
37 typedef std::function<VectorRd(const VectorRd &)> MagneticForcingTermType;
38 typedef std::function<double(const VectorRd &)> CompressibilityForcingTermType;
39
40 typedef std::function<VectorRd(const VectorRd &)> VelocityType;
41 typedef std::function<MatrixRd(const VectorRd &)> VelocityGradientType;
42
43 typedef std::function<double(const VectorRd &)> PressureType;
44 typedef std::function<VectorRd(const VectorRd &)> PressureGradientType;
45
46
48 const Mesh & mesh,
49 size_t degree,
50 const BoundaryConditions & bc,
51 const HYPREParameters & parameters = {}
52 );
53 virtual ~HYPRE() {}
54
56 virtual Eigen::VectorXd interpolate(
57 const VelocityType & u,
58 const PressureType & p
59 ) const = 0;
60
62 template<typename F>
63 Eigen::VectorXd pkpoPkPkPkInterpolate(
64 const F & v,
65 const InterpolateParameters & parameters = {}
66 ) const;
67
69 size_t degree() const
70 {
71 return m_degree;
72 }
73
76 {
78 }
79
81 const std::vector<Eigen::MatrixXd> potential() const
82 {
83 return m_potential;
84 }
85
87 const std::vector<Eigen::MatrixXd> stabilization() const
88 {
89 return m_stabilization;
90 }
91
93 const Eigen::MatrixXd & potential(size_t iT) const
94 {
95 return m_potential[iT];
96 }
97
100 {
101 return m_nloc_sc_u;
102 }
103
106 {
107 return m_nloc_sc_p;
108 }
109
112 {
113 return m_velocity_space->mesh().n_cells() * m_nloc_sc_u;
114 }
115
118 {
119 return m_velocity_space->mesh().n_cells() * m_nloc_sc_p;
120 }
121
127
129 size_t numDirichletDofs() const
130 {
131 return 2*(m_bc.n_dir_edges() * m_velocity_space->numberOfLocalEdgeDofs());
132 }
133
135 size_t dimVelocity() const
136 {
137 return m_velocity_space->dimension();
138 }
139
141 size_t dimPressure() const
142 {
143 return m_pressure_space->dimension();
144 }
145
148 {
149 return 2*dimVelocity() + 2*dimPressure() + 2 - numStaticallyCondensedDofs();
150 }
151
153 size_t sizeSystem() const
154 {
156 }
157
163
169
170
172 inline const DiscreteSpace & velocitySpace() const
173 {
174 return *m_velocity_space;
175 }
176
178 inline const DiscreteSpace & pressureSpace() const
179 {
180 return *m_pressure_space;
181 }
182
184 inline const Mesh & mesh() const
185 {
186 return m_velocity_space->mesh();
187 }
188
190 inline const SystemMatrixType & systemMatrix() const {
191 return m_A;
192 }
193
196 return m_A;
197 }
198
200 inline const Eigen::VectorXd & systemVector() const {
201 return m_b;
202 }
203
205 inline Eigen::VectorXd & systemVector() {
206 return m_b;
207 }
208
211 return m_sc_A;
212 }
213
215 inline const Eigen::VectorXd & staticallyCondensedVector() const {
216 return m_sc_b;
217 }
218
220 inline const BoundaryConditions & boundaryConditions() const {
221 return m_bc;
222 }
223
227 const Eigen::VectorXd & uph_old,
228 const Eigen::VectorXd & uph_prev,
229 const Eigen::VectorXd & brh_old,
230 const Eigen::VectorXd & brh_prev,
231 const double time_step,
232 const double current_time);
233
234
235 //[Vito] This should be the stationary case assembleResidual case!
236
239 const Eigen::VectorXd & uph_old,
240 const Eigen::VectorXd & brh_old) {
241 Eigen::VectorXd uph_prev = Eigen::VectorXd::Zero(uph_old.size());
242 Eigen::VectorXd brh_prev = Eigen::VectorXd::Zero(brh_old.size());
243 assembleResidualSystem(isolution,
244 uph_old,
245 uph_prev,
246 brh_old,
247 brh_prev,
248 std::numeric_limits<double>::infinity(),
249 0.);
250 };
251
253 Eigen::VectorXd reconstructSolution(const Eigen::VectorXd & uhp_solsystem) const;
254
257 const size_t cell_degree,
258 const double viscosity,
259 const double current_time = 0.) const {
260 // Construct quadrature scheme on element
261 const Cell & T = *m_velocity_space->mesh().cell(iT);
262 QuadratureRule quad = generate_quadrature_rule(T, 3 * cell_degree);
263
264 // Local convective stabilization bilinear form assumed to be given by
265 // j_T(w_T, v_T) = Beta_T Sum_{F in F_T} Int_{F} (w_F - w_T).(v_F - v_T)
266 // where Beta_T = max(c_s, ||u_T||_{L^infty(T)}). We therefore need an
267 // estimate of ||u_T||_{L^infty(T)} which we obtain by taking the largest
268 // l^infty norm of u_T on the quadrature points on T.
269 double beta_T = 1.e-4; // Safeguard parameter
270 for (size_t iqn = 0; iqn < quad.size(); iqn++) {
271 // Stabilisation based on max-norm of discrete solution
272 // VectorRd uT_old_iqn = VectorRd::Zero();
273 // for (size_t i = 0; i < cell_dofs->dimension(); i++) {
274 // double uT_old_i = uT_old(m_velocity_space->localOffset(T) + i);
275 // uT_old_iqn += uT_old_i * cell_dofs_quad[i][iqn];
276 // } // for i
277 // Stabilisation based on max-norm of analytic solution
278 Eigen::VectorXd uT_old_iqn = isolution->velocity(quad[iqn].vector(), current_time);
279 beta_T = std::max(beta_T, uT_old_iqn.lpNorm<Eigen::Infinity>());
280 } // for iqn
281
282 // ALTERNATIVE APPROACH (based on maximum absolute value of u.n on the boundary of T)
283 // double beta_T = 1.e-4; // Safeguard parameter
284 // for (size_t iE = 0; iE < T.n_edges(); iE++) {
285 // const Edge & E = *T.edge(iE);
286 // auto nTE = T.edge_normal(iE);
287 // QuadratureRule quad_E = generate_quadrature_rule(E, 3 * cell_degree);
288 // for (size_t iqn = 0; iqn < quad_E.size(); iqn++) {
289 // Eigen::VectorXd uE_old_iqn = isolution->velocity(quad_E[iqn].vector());
290 // beta_T = std::max(beta_T, std::abs(uE_old_iqn.dot(nTE)));
291 // }
292 // }
293
294 return beta_T;
295 };
296
299 const size_t cell_degree,
300 const double magnetic_diffusivity,
301 const double current_time = 0.) const {
302 // Construct quadrature scheme on element
303 const Cell & T = *m_velocity_space->mesh().cell(iT);
304 QuadratureRule quad = generate_quadrature_rule(T, 3 * cell_degree);
305
306 // Local convective stabilization bilinear form assumed to be given by
307 // j_T(w_T, v_T) = Beta_T Sum_{F in F_T} Int_{F} (w_F - w_T).(v_F - v_T)
308 // where Beta_T = max(c_s, ||u_T||_{L^infty(T)}). We therefore need an
309 // estimate of ||u_T||_{L^infty(T)} which we obtain by taking the largest
310 // l^infty norm of u_T on the quadrature points on T.
311 double gamma_T = 1.e-4; // Safeguard parameter
312 for (size_t iqn = 0; iqn < quad.size(); iqn++) {
313 // Stabilisation based on max-norm of analytic solution
314 Eigen::VectorXd bT_old_iqn = isolution->magnetic_field(quad[iqn].vector(), current_time);
315 gamma_T = std::max(gamma_T, bT_old_iqn.lpNorm<Eigen::Infinity>());
316 } // for iqn
317 return gamma_T;
318 };
319
320
321
322
323 protected:
325
326 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _time_deriv_term(size_t iT,
327 const Eigen::VectorXd & uT_old,
328 const Eigen::VectorXd & uT_prev,
329 const double time_step); // = std::numeric_limits<double>::infinity()
330
331
332 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _viscous_term(size_t iT, const Eigen::VectorXd & uT_old);
333
334 std::pair<Eigen::MatrixXd, Eigen::VectorXd> _magnetic_diffusive_term(size_t iT, const Eigen::VectorXd & bT_old);
335
336
337 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd> _velocity_pressure_coupling(
338 size_t iT,
340 const Eigen::VectorXd & uT_old,
341 const Eigen::VectorXd & pT_old,
342 const double & current_time = 0.
343 );
344
345 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd> _magnetic_field_magnetic_pressure_coupling(
346 size_t iT,
348 const Eigen::VectorXd & bT_old,
349 const Eigen::VectorXd & rT_old,
350 const double & current_time = 0.
351 );
352
353
354 virtual Eigen::MatrixXd _forcing_terms(
355 size_t iT,
357 const double current_time = 0.
358 );
359
360 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
361 _linearized_convective_term(size_t iT, const Eigen::VectorXd & uT_old);
362
363 virtual std::tuple<Eigen::MatrixXd, Eigen::MatrixXd, Eigen::VectorXd>
364 _linearized_convective_mixed_term(size_t iT, const Eigen::VectorXd & wT_old, const Eigen::VectorXd & vT_old);
365
366 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
368 const Eigen::VectorXd & uT_old,
370 const double current_time = 0.);
371
372 virtual std::pair<Eigen::MatrixXd, Eigen::VectorXd>
374 const Eigen::VectorXd & bT_old,
376 const double current_time = 0.);
377
378 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
380 size_t iT,
382 const Eigen::VectorXd & uph_old,
383 const Eigen::VectorXd & uph_prev,
384 const Eigen::VectorXd & brh_old,
385 const Eigen::VectorXd & brh_prev,
386 const double time_step,
387 const double current_time);
388
389 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
391 size_t iT,
393 const Eigen::VectorXd & uph_old,
394 const Eigen::VectorXd & brh_old) {
395 Eigen::VectorXd uph_prev = uph_old;
396 Eigen::VectorXd brh_prev = brh_old;
397
399 isolution,
400 uph_old,
401 uph_prev,
402 brh_old,
403 brh_prev,
404 std::numeric_limits<double>::infinity(),
405 0.);
406
407 }
408
410
412 size_t iT,
413 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
414 std::list<Eigen::Triplet<double> > & A1,
415 Eigen::VectorXd & b1,
416 std::list<Eigen::Triplet<double> > & A2,
417 Eigen::VectorXd & b2
418 );
419
420 size_t m_degree;
421 const BoundaryConditions & m_bc;
423 bool m_use_threads;
425
427 double m_viscosity;
429
430 std::ostream & m_output;
431
432 size_t m_nloc_sc_u;
433 size_t m_nloc_sc_p;
434
435 std::vector<Eigen::MatrixXd> m_discrete_l2;
436 std::vector<Eigen::MatrixXd> m_gradient;
437 std::vector<Eigen::MatrixXd> m_gradient_rhs;
438 std::vector<Eigen::MatrixXd> m_potential;
439 std::vector<Eigen::MatrixXd> m_stabilization;
440
441 std::unique_ptr<DiscreteSpace> m_velocity_space;
442 std::unique_ptr<DiscreteSpace> m_pressure_space;
443
444
446 Eigen::VectorXd m_b;
447
449 Eigen::VectorXd m_sc_b;
450
451 std::vector<int> m_nsc_to_dof_map;
452 std::vector<int> m_sc_to_dof_map;
453
454 };
455
456 } // namespace DSL
457
458} // namespace HArDCore2D
459
460#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:93
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:250
double m_magnetic_diffusivity
Definition hypre.hpp:428
const std::vector< Eigen::MatrixXd > potential() const
Returns the vector of potential reconstructions.
Definition hypre.hpp:81
double magnetic_convective_stabilization_parameter(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const size_t cell_degree, const double magnetic_diffusivity, const double current_time=0.) const
Definition hypre.hpp:297
bool m_static_condensation
Definition hypre.hpp:347
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _magnetic_convective_stabilization_term(size_t iT, const Eigen::VectorXd &bT_old, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const double current_time=0.)
Definition hypre.cpp:692
std::function< VectorRd(const VectorRd &)> VelocityType
Definition hypre.hpp:40
const SystemMatrixType & staticallyCondensedMatrix() const
Returns the static condensation recovery operator.
Definition hypre.hpp:210
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:390
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.
std::function< VectorRd(const VectorRd &)> PressureGradientType
Definition hypre.hpp:44
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition hypre.hpp:190
std::pair< Eigen::MatrixXd, Eigen::VectorXd > _magnetic_diffusive_term(size_t iT, const Eigen::VectorXd &bT_old)
Definition hypre.cpp:176
const std::vector< Eigen::MatrixXd > stabilization() const
Returns the vector of stabilization matrices.
Definition hypre.hpp:87
const Mesh & mesh() const
Returns the mesh.
Definition hypre.hpp:166
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.)
std::ostream & m_output
Definition hypre.hpp:353
std::function< MatrixRd(const VectorRd &)> VelocityGradientType
Definition hypre.hpp:41
const DiscreteSpace & pressureSpace() const
Returns the pressure space.
Definition hypre.hpp:178
size_t numNonStaticallyCondensedDofsPressure() const
Returns the number of non statically condensed DOFs for the pressure.
Definition hypre.hpp:165
size_t numNonStaticallyCondensedDofsVelocity() const
Returns the number of non statically condensed DOFs for the velocity.
Definition hypre.hpp:159
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:75
size_t numStaticallyCondensedDofsVelocity() const
Returns the number of statically condensed DOFs for the velocity.
Definition hypre.hpp:111
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:123
size_t numLocalStaticallyCondensedDofsVelocity() const
Returns the local number of velocity statically condensed DOFs.
Definition hypre.hpp:99
std::unique_ptr< DiscreteSpace > m_pressure_space
Definition hypre.hpp:367
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:172
Eigen::SparseMatrix< double > SystemMatrixType
Definition hypre.hpp:34
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:200
const BoundaryConditions & boundaryConditions() const
Returns the boundary condition handler.
Definition hypre.hpp:220
virtual std::pair< Eigen::MatrixXd, Eigen::VectorXd > _linearized_convective_term(size_t iT, const Eigen::VectorXd &uT_old)
std::function< VectorRd(const VectorRd &)> MomentumForcingTermType
Definition hypre.hpp:36
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:215
size_t dimVelocity() const
Returns the dimension of velocity space.
Definition hypre.hpp:135
std::unique_ptr< DiscreteSpace > m_velocity_space
Definition hypre.hpp:366
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition hypre.hpp:195
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:38
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:153
std::vector< Eigen::MatrixXd > m_stabilization
Definition hypre.hpp:364
size_t dimPressure() const
Returns the dimension of pressure space.
Definition hypre.hpp:141
size_t numLocalStaticallyCondensedDofsPressure() const
Returns the local number of pressure statically condensed DOFs.
Definition hypre.hpp:105
std::function< VectorRd(const VectorRd &)> MagneticForcingTermType
Definition hypre.hpp:37
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:482
size_t numNonStaticallyCondensedDofs() const
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DO...
Definition hypre.hpp:147
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 ~HYPRE()
Definition hypre.hpp:53
std::function< double(const VectorRd &)> PressureType
Definition hypre.hpp:43
bool m_VitoParameter
Definition hypre.hpp:424
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
double convective_stabilization_parameter(size_t iT, const HArDCore2D::NavierStokesSolutions::IExactSolution *isolution, const size_t cell_degree, const double viscosity, const double current_time=0.) const
Definition hypre.hpp:255
void _construct_local_operators(size_t iT)
size_t numStaticallyCondensedDofsPressure() const
Returns the number of statically condensed DOFs for the pressure.
Definition hypre.hpp:117
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition hypre.hpp:205
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:237
LocalStaticCondensation _compute_static_condensation(const size_t &iT)
size_t numDirichletDofs() const
Returns the number of Dirichlet DOFs.
Definition hypre.hpp:129
Definition Mesh2D.hpp:26
Eigen::Vector2d VectorRd
Definition basis.hpp:55
Eigen::Matrix2d MatrixRd
Definition basis.hpp:54
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:55
QuadratureRule generate_quadrature_rule(const Cell &T, const int doe, const bool force_split)
Generate quadrature rule on mesh element.
Definition quadraturerule.cpp:10
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
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
double viscosity
Definition hypre.hpp:25
double magnetic_diffusivity
Definition hypre.hpp:28
Definition hho-interpolate.hpp:15
Structure to store information for, and perform, local static condensation.
Definition local_static_condensation.hpp:25
virtual Eigen::Vector2d velocity(const Eigen::Vector2d &x, const double &t=0) const =0
virtual Eigen::Vector2d magnetic_field(const Eigen::Vector2d &x, const double &t=0) const =0