MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::navier::NavierSolver Class Reference

Transient incompressible Navier Stokes solver in a split scheme formulation. More...

#include <navier_solver.hpp>

Collaboration diagram for mfem::navier::NavierSolver:
[legend]

Public Member Functions

 NavierSolver (ParMesh *mesh, int order, double kin_vis)
 Initialize data structures, set FE space order and kinematic viscosity. More...
 
void Setup (double dt)
 Initialize forms, solvers and preconditioners. More...
 
void Step (double &time, double dt, int cur_step)
 Compute solution at the next time step t+dt. More...
 
ParGridFunctionGetCurrentVelocity ()
 Return a pointer to the current velocity ParGridFunction. More...
 
ParGridFunctionGetCurrentPressure ()
 Return a pointer to the current pressure ParGridFunction. More...
 
void AddVelDirichletBC (VectorCoefficient *coeff, Array< int > &attr)
 Add a Dirichlet boundary condition to the velocity field. More...
 
void AddVelDirichletBC (VecFuncT *f, Array< int > &attr)
 
void AddPresDirichletBC (Coefficient *coeff, Array< int > &attr)
 Add a Dirichlet boundary condition to the pressure field. More...
 
void AddPresDirichletBC (ScalarFuncT *f, Array< int > &attr)
 
void AddAccelTerm (VectorCoefficient *coeff, Array< int > &attr)
 Add an acceleration term to the RHS of the equation. More...
 
void AddAccelTerm (VecFuncT *f, Array< int > &attr)
 
void EnablePA (bool pa)
 Enable partial assembly for every operator. More...
 
void EnableNI (bool ni)
 
void PrintTimingData ()
 Print timing summary of the solving routine. More...
 
 ~NavierSolver ()
 
void ComputeCurl2D (ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
 Compute \(\nabla \times \nabla \times u\) for \(u \in (H^1)^2\). More...
 
void ComputeCurl3D (ParGridFunction &u, ParGridFunction &cu)
 Compute \(\nabla \times \nabla \times u\) for \(u \in (H^1)^3\). More...
 
void Orthogonalize (Vector &v)
 Remove mean from a Vector. More...
 
void MeanZero (ParGridFunction &v)
 Remove the mean from a ParGridFunction. More...
 
double ComputeCFL (ParGridFunction &u, double dt)
 Compute CFL. More...
 

Protected Member Functions

void PrintInfo ()
 Print informations about the Navier version. More...
 
void SetTimeIntegrationCoefficients (int step)
 Update the EXTk/BDF time integration coefficient. More...
 
void EliminateRHS (Operator &A, ConstrainedOperator &constrainedA, const Array< int > &ess_tdof_list, Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0)
 Eliminate essential BCs in an Operator and apply to RHS. More...
 

Protected Attributes

bool debug = false
 Enable/disable debug output. More...
 
bool verbose = true
 Enable/disable verbose output. More...
 
bool partial_assembly = false
 Enable/disable partial assembly of forms. More...
 
bool numerical_integ = false
 Enable/disable numerical integration rules of forms. More...
 
ParMeshpmesh = nullptr
 The parallel mesh. More...
 
int order
 The order of the velocity and pressure space. More...
 
double kin_vis
 Kinematic viscosity (dimensionless). More...
 
FiniteElementCollectionvfec = nullptr
 Velocity \(H^1\) finite element collection. More...
 
FiniteElementCollectionpfec = nullptr
 Pressure \(H^1\) finite element collection. More...
 
ParFiniteElementSpacevfes = nullptr
 Velocity \((H^1)^d\) finite element space. More...
 
ParFiniteElementSpacepfes = nullptr
 Pressure \(H^1\) finite element space. More...
 
ParNonlinearFormN = nullptr
 
ParBilinearFormMv_form = nullptr
 
ParBilinearFormSp_form = nullptr
 
ParMixedBilinearFormD_form = nullptr
 
ParMixedBilinearFormG_form = nullptr
 
ParBilinearFormH_form = nullptr
 
VectorGridFunctionCoefficientFText_gfcoeff = nullptr
 
ParLinearFormFText_bdr_form = nullptr
 
ParLinearFormf_form = nullptr
 
ParLinearFormg_bdr_form = nullptr
 
ParLinearFormmass_lf = nullptr
 Linear form to compute the mass matrix in various subroutines. More...
 
ConstantCoefficient onecoeff
 
double volume = 0.0
 
ConstantCoefficient nlcoeff
 
ConstantCoefficient Sp_coeff
 
ConstantCoefficient H_lincoeff
 
ConstantCoefficient H_bdfcoeff
 
OperatorHandle Mv
 
OperatorHandle Sp
 
OperatorHandle D
 
OperatorHandle G
 
OperatorHandle H
 
SolverMvInvPC = nullptr
 
CGSolverMvInv = nullptr
 
HypreBoomerAMGSpInvPC = nullptr
 
OrthoSolverSpInvOrthoPC = nullptr
 
CGSolverSpInv = nullptr
 
SolverHInvPC = nullptr
 
CGSolverHInv = nullptr
 
Vector fn
 
Vector un
 
Vector unm1
 
Vector unm2
 
Vector Nun
 
Vector Nunm1
 
Vector Nunm2
 
Vector Fext
 
Vector FText
 
Vector Lext
 
Vector resu
 
Vector tmp1
 
Vector pn
 
Vector resp
 
Vector FText_bdr
 
Vector g_bdr
 
ParGridFunction un_gf
 
ParGridFunction curlu_gf
 
ParGridFunction curlcurlu_gf
 
ParGridFunction Lext_gf
 
ParGridFunction FText_gf
 
ParGridFunction resu_gf
 
ParGridFunction pn_gf
 
ParGridFunction resp_gf
 
Array< int > vel_ess_attr
 
Array< int > pres_ess_attr
 
Array< int > vel_ess_tdof
 
Array< int > pres_ess_tdof
 
std::vector< VelDirichletBC_Tvel_dbcs
 
std::vector< PresDirichletBC_Tpres_dbcs
 
std::vector< AccelTerm_Taccel_terms
 
int cur_step = 0
 
double bd0 = 0.0
 
double bd1 = 0.0
 
double bd2 = 0.0
 
double bd3 = 0.0
 
double ab1 = 0.0
 
double ab2 = 0.0
 
double ab3 = 0.0
 
StopWatch sw_setup
 
StopWatch sw_step
 
StopWatch sw_extrap
 
StopWatch sw_curlcurl
 
StopWatch sw_spsolve
 
StopWatch sw_hsolve
 
int pl_mvsolve = 0
 
int pl_spsolve = 0
 
int pl_hsolve = 0
 
int pl_amg = 0
 
double rtol_spsolve = 1e-6
 
double rtol_hsolve = 1e-8
 
int iter_mvsolve = 0
 
int iter_spsolve = 0
 
int iter_hsolve = 0
 
double res_mvsolve = 0.0
 
double res_spsolve = 0.0
 
double res_hsolve = 0.0
 
ParMeshpmesh_lor = nullptr
 
FiniteElementCollectionpfec_lor = nullptr
 
ParFiniteElementSpacepfes_lor = nullptr
 
InterpolationGridTransfervgt = nullptr
 
InterpolationGridTransferpgt = nullptr
 
ParBilinearFormMv_form_lor = nullptr
 
ParBilinearFormSp_form_lor = nullptr
 
ParBilinearFormH_form_lor = nullptr
 
OperatorHandle Mv_lor
 
OperatorHandle Sp_lor
 
OperatorHandle H_lor
 

Detailed Description

Transient incompressible Navier Stokes solver in a split scheme formulation.

This implementation of a transient incompressible Navier Stokes solver uses the non-dimensionalized formulation. The coupled momentum and incompressibilty equations are decoupled using the split scheme described in [1]. This leads to three solving steps.

  1. An extrapolation step for all nonlinear terms which are treated explicitly. This step avoids a fully coupled nonlinear solve and only requires a solve of the mass matrix in velocity space \(M_v^{-1}\). On the other hand this introduces a CFL stability condition on the maximum timestep.
  2. A Poisson solve \(S_p^{-1}\).
  3. A Helmholtz like solve \((M_v - \partial t K_v)^{-1}\).

The numerical solver setup for each step are as follows.

\(M_v^{-1}\) is solved using CG with Jacobi as preconditioner.

\(S_p^{-1}\) is solved using CG with AMG applied to the low order refined (LOR) assembled pressure Poisson matrix. To avoid assembling a matrix for preconditioning, one can use p-MG as an alternative (NYI).

\((M_v - \partial t K_v)^{-1}\) due to the CFL condition we expect the time step to be small. Therefore this is solved using CG with Jacobi as preconditioner. For large time steps a preconditioner like AMG or p-MG should be used (NYI).

Statements marked with NYI mean this feature is planned but Not Yet Implemented.

A detailed description is available in [1] section 4.2. The algorithm is originated from [2].

[1] Michael Franco, Jean-Sylvain Camier, Julian Andrej, Will Pazner (2020) High-order matrix-free incompressible flow solvers with GPU acceleration and low-order refined preconditioners (https://arxiv.org/abs/1910.03032)

[2] A. G. Tomboulides, J. C. Y. Lee & S. A. Orszag (1997) Numerical Simulation of Low Mach Number Reactive Flows

Definition at line 142 of file navier_solver.hpp.

Constructor & Destructor Documentation

NavierSolver::NavierSolver ( ParMesh mesh,
int  order,
double  kin_vis 
)

Initialize data structures, set FE space order and kinematic viscosity.

The ParMesh mesh can be a linear or curved parallel mesh. The order of the finite element spaces is this algorithm is of equal order \((P_N)^d P_N\) for velocity and pressure respectively. This means the pressure is in discretized in the same space (just scalar instead of a vector space) as the velocity.

Kinematic viscosity (dimensionless) is set using kin_vis and automatically converted to the Reynolds number. If you want to set the Reynolds number directly, you can provide the inverse.

Definition at line 29 of file navier_solver.cpp.

NavierSolver::~NavierSolver ( )

Definition at line 1050 of file navier_solver.cpp.

Member Function Documentation

void NavierSolver::AddAccelTerm ( VectorCoefficient coeff,
Array< int > &  attr 
)

Add an acceleration term to the RHS of the equation.

The VecFuncT f is evaluated at the current time t and extrapolated together with the nonlinear parts of the Navier Stokes equation.

Definition at line 943 of file navier_solver.cpp.

void NavierSolver::AddAccelTerm ( VecFuncT f,
Array< int > &  attr 
)

Definition at line 961 of file navier_solver.cpp.

void NavierSolver::AddPresDirichletBC ( Coefficient coeff,
Array< int > &  attr 
)

Add a Dirichlet boundary condition to the pressure field.

Definition at line 910 of file navier_solver.cpp.

void NavierSolver::AddPresDirichletBC ( ScalarFuncT f,
Array< int > &  attr 
)

Definition at line 938 of file navier_solver.cpp.

void NavierSolver::AddVelDirichletBC ( VectorCoefficient coeff,
Array< int > &  attr 
)

Add a Dirichlet boundary condition to the velocity field.

Definition at line 877 of file navier_solver.cpp.

void NavierSolver::AddVelDirichletBC ( VecFuncT f,
Array< int > &  attr 
)

Definition at line 905 of file navier_solver.cpp.

double NavierSolver::ComputeCFL ( ParGridFunction u,
double  dt 
)

Compute CFL.

Definition at line 797 of file navier_solver.cpp.

void NavierSolver::ComputeCurl2D ( ParGridFunction u,
ParGridFunction cu,
bool  assume_scalar = false 
)

Compute \(\nabla \times \nabla \times u\) for \(u \in (H^1)^2\).

Definition at line 696 of file navier_solver.cpp.

void NavierSolver::ComputeCurl3D ( ParGridFunction u,
ParGridFunction cu 
)

Compute \(\nabla \times \nabla \times u\) for \(u \in (H^1)^3\).

Definition at line 605 of file navier_solver.cpp.

void NavierSolver::EliminateRHS ( Operator A,
ConstrainedOperator constrainedA,
const Array< int > &  ess_tdof_list,
Vector x,
Vector b,
Vector X,
Vector B,
int  copy_interior = 0 
)
protected

Eliminate essential BCs in an Operator and apply to RHS.

Definition at line 572 of file navier_solver.cpp.

void mfem::navier::NavierSolver::EnableNI ( bool  ni)
inline

Enable numerical integration rules. This means collocated quadrature at the nodal points.

Definition at line 195 of file navier_solver.hpp.

void mfem::navier::NavierSolver::EnablePA ( bool  pa)
inline

Enable partial assembly for every operator.

Definition at line 191 of file navier_solver.hpp.

ParGridFunction* mfem::navier::NavierSolver::GetCurrentPressure ( )
inline

Return a pointer to the current pressure ParGridFunction.

Definition at line 169 of file navier_solver.hpp.

ParGridFunction* mfem::navier::NavierSolver::GetCurrentVelocity ( )
inline

Return a pointer to the current velocity ParGridFunction.

Definition at line 166 of file navier_solver.hpp.

void NavierSolver::MeanZero ( ParGridFunction v)

Remove the mean from a ParGridFunction.

Modify the ParGridFunction v by subtracting its mean using \( v = v - \int_\Omega \frac{v}{vol(\Omega)} dx \).

Definition at line 550 of file navier_solver.cpp.

void NavierSolver::Orthogonalize ( Vector v)

Remove mean from a Vector.

Modify the Vector v by subtracting its mean using \(v = v - \frac{\sum_i^N v_i}{N} \)

Definition at line 592 of file navier_solver.cpp.

void NavierSolver::PrintInfo ( )
protected

Print informations about the Navier version.

Definition at line 1035 of file navier_solver.cpp.

void NavierSolver::PrintTimingData ( )

Print timing summary of the solving routine.

The summary shows the timing in seconds in the first row of

  1. SETUP: Time spent for the setup of all forms, solvers and preconditioners.
  2. STEP: Time spent computing a full time step. It includes all three solves.
  3. EXTRAP: Time spent for extrapolation of all forcing and nonlinear terms.
  4. CURLCURL: Time spent for computing the curl curl term in the pressure Poisson equation (see references for detailed explanation).
  5. PSOLVE: Time spent in the pressure Poisson solve.
  6. HSOLVE: Time spent in the Helmholtz solve.

The second row shows a proportion of a column relative to the whole time step.

Definition at line 1000 of file navier_solver.cpp.

void NavierSolver::SetTimeIntegrationCoefficients ( int  step)
protected

Update the EXTk/BDF time integration coefficient.

Depending on which time step the computation is in, the EXTk/BDF time integration coefficients have to be set accordingly. This allows bootstrapping with a BDF scheme of order 1 and increasing the order each following time step, up to order 3.

Definition at line 966 of file navier_solver.cpp.

void NavierSolver::Setup ( double  dt)

Initialize forms, solvers and preconditioners.

Definition at line 95 of file navier_solver.cpp.

void NavierSolver::Step ( double &  time,
double  dt,
int  cur_step 
)

Compute solution at the next time step t+dt.

Definition at line 296 of file navier_solver.cpp.

Member Data Documentation

double mfem::navier::NavierSolver::ab1 = 0.0
protected

Definition at line 379 of file navier_solver.hpp.

double mfem::navier::NavierSolver::ab2 = 0.0
protected

Definition at line 380 of file navier_solver.hpp.

double mfem::navier::NavierSolver::ab3 = 0.0
protected

Definition at line 381 of file navier_solver.hpp.

std::vector<AccelTerm_T> mfem::navier::NavierSolver::accel_terms
protected

Definition at line 370 of file navier_solver.hpp.

double mfem::navier::NavierSolver::bd0 = 0.0
protected

Definition at line 375 of file navier_solver.hpp.

double mfem::navier::NavierSolver::bd1 = 0.0
protected

Definition at line 376 of file navier_solver.hpp.

double mfem::navier::NavierSolver::bd2 = 0.0
protected

Definition at line 377 of file navier_solver.hpp.

double mfem::navier::NavierSolver::bd3 = 0.0
protected

Definition at line 378 of file navier_solver.hpp.

int mfem::navier::NavierSolver::cur_step = 0
protected

Definition at line 372 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::curlcurlu_gf
protected

Definition at line 351 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::curlu_gf
protected

Definition at line 351 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::D
protected

Definition at line 332 of file navier_solver.hpp.

ParMixedBilinearForm* mfem::navier::NavierSolver::D_form = nullptr
protected

Definition at line 306 of file navier_solver.hpp.

bool mfem::navier::NavierSolver::debug = false
protected

Enable/disable debug output.

Definition at line 268 of file navier_solver.hpp.

ParLinearForm* mfem::navier::NavierSolver::f_form = nullptr
protected

Definition at line 316 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::Fext
protected

Definition at line 346 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::fn
protected

Definition at line 346 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::FText
protected

Definition at line 346 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::FText_bdr
protected

Definition at line 349 of file navier_solver.hpp.

ParLinearForm* mfem::navier::NavierSolver::FText_bdr_form = nullptr
protected

Definition at line 314 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::FText_gf
protected

Definition at line 351 of file navier_solver.hpp.

VectorGridFunctionCoefficient* mfem::navier::NavierSolver::FText_gfcoeff = nullptr
protected

Definition at line 312 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::G
protected

Definition at line 333 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::g_bdr
protected

Definition at line 349 of file navier_solver.hpp.

ParLinearForm* mfem::navier::NavierSolver::g_bdr_form = nullptr
protected

Definition at line 318 of file navier_solver.hpp.

ParMixedBilinearForm* mfem::navier::NavierSolver::G_form = nullptr
protected

Definition at line 308 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::H
protected

Definition at line 334 of file navier_solver.hpp.

ConstantCoefficient mfem::navier::NavierSolver::H_bdfcoeff
protected

Definition at line 328 of file navier_solver.hpp.

ParBilinearForm* mfem::navier::NavierSolver::H_form = nullptr
protected

Definition at line 310 of file navier_solver.hpp.

ParBilinearForm* mfem::navier::NavierSolver::H_form_lor = nullptr
protected

Definition at line 410 of file navier_solver.hpp.

ConstantCoefficient mfem::navier::NavierSolver::H_lincoeff
protected

Definition at line 327 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::H_lor
protected

Definition at line 414 of file navier_solver.hpp.

CGSolver* mfem::navier::NavierSolver::HInv = nullptr
protected

Definition at line 344 of file navier_solver.hpp.

Solver* mfem::navier::NavierSolver::HInvPC = nullptr
protected

Definition at line 343 of file navier_solver.hpp.

int mfem::navier::NavierSolver::iter_hsolve = 0
protected

Definition at line 397 of file navier_solver.hpp.

int mfem::navier::NavierSolver::iter_mvsolve = 0
protected

Definition at line 397 of file navier_solver.hpp.

int mfem::navier::NavierSolver::iter_spsolve = 0
protected

Definition at line 397 of file navier_solver.hpp.

double mfem::navier::NavierSolver::kin_vis
protected

Kinematic viscosity (dimensionless).

Definition at line 286 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::Lext
protected

Definition at line 346 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::Lext_gf
protected

Definition at line 351 of file navier_solver.hpp.

ParLinearForm* mfem::navier::NavierSolver::mass_lf = nullptr
protected

Linear form to compute the mass matrix in various subroutines.

Definition at line 321 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::Mv
protected

Definition at line 330 of file navier_solver.hpp.

ParBilinearForm* mfem::navier::NavierSolver::Mv_form = nullptr
protected

Definition at line 302 of file navier_solver.hpp.

ParBilinearForm* mfem::navier::NavierSolver::Mv_form_lor = nullptr
protected

Definition at line 408 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::Mv_lor
protected

Definition at line 412 of file navier_solver.hpp.

CGSolver* mfem::navier::NavierSolver::MvInv = nullptr
protected

Definition at line 337 of file navier_solver.hpp.

Solver* mfem::navier::NavierSolver::MvInvPC = nullptr
protected

Definition at line 336 of file navier_solver.hpp.

ParNonlinearForm* mfem::navier::NavierSolver::N = nullptr
protected

Definition at line 300 of file navier_solver.hpp.

ConstantCoefficient mfem::navier::NavierSolver::nlcoeff
protected

Definition at line 325 of file navier_solver.hpp.

bool mfem::navier::NavierSolver::numerical_integ = false
protected

Enable/disable numerical integration rules of forms.

Definition at line 277 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::Nun
protected

Definition at line 346 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::Nunm1
protected

Definition at line 346 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::Nunm2
protected

Definition at line 346 of file navier_solver.hpp.

ConstantCoefficient mfem::navier::NavierSolver::onecoeff
protected

Definition at line 322 of file navier_solver.hpp.

int mfem::navier::NavierSolver::order
protected

The order of the velocity and pressure space.

Definition at line 283 of file navier_solver.hpp.

bool mfem::navier::NavierSolver::partial_assembly = false
protected

Enable/disable partial assembly of forms.

Definition at line 274 of file navier_solver.hpp.

FiniteElementCollection* mfem::navier::NavierSolver::pfec = nullptr
protected

Pressure \(H^1\) finite element collection.

Definition at line 292 of file navier_solver.hpp.

FiniteElementCollection* mfem::navier::NavierSolver::pfec_lor = nullptr
protected

Definition at line 404 of file navier_solver.hpp.

ParFiniteElementSpace* mfem::navier::NavierSolver::pfes = nullptr
protected

Pressure \(H^1\) finite element space.

Definition at line 298 of file navier_solver.hpp.

ParFiniteElementSpace* mfem::navier::NavierSolver::pfes_lor = nullptr
protected

Definition at line 405 of file navier_solver.hpp.

InterpolationGridTransfer * mfem::navier::NavierSolver::pgt = nullptr
protected

Definition at line 406 of file navier_solver.hpp.

int mfem::navier::NavierSolver::pl_amg = 0
protected

Definition at line 390 of file navier_solver.hpp.

int mfem::navier::NavierSolver::pl_hsolve = 0
protected

Definition at line 389 of file navier_solver.hpp.

int mfem::navier::NavierSolver::pl_mvsolve = 0
protected

Definition at line 387 of file navier_solver.hpp.

int mfem::navier::NavierSolver::pl_spsolve = 0
protected

Definition at line 388 of file navier_solver.hpp.

ParMesh* mfem::navier::NavierSolver::pmesh = nullptr
protected

The parallel mesh.

Definition at line 280 of file navier_solver.hpp.

ParMesh* mfem::navier::NavierSolver::pmesh_lor = nullptr
protected

Definition at line 403 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::pn
protected

Definition at line 349 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::pn_gf
protected

Definition at line 353 of file navier_solver.hpp.

std::vector<PresDirichletBC_T> mfem::navier::NavierSolver::pres_dbcs
protected

Definition at line 367 of file navier_solver.hpp.

Array<int> mfem::navier::NavierSolver::pres_ess_attr
protected

Definition at line 357 of file navier_solver.hpp.

Array<int> mfem::navier::NavierSolver::pres_ess_tdof
protected

Definition at line 361 of file navier_solver.hpp.

double mfem::navier::NavierSolver::res_hsolve = 0.0
protected

Definition at line 400 of file navier_solver.hpp.

double mfem::navier::NavierSolver::res_mvsolve = 0.0
protected

Definition at line 400 of file navier_solver.hpp.

double mfem::navier::NavierSolver::res_spsolve = 0.0
protected

Definition at line 400 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::resp
protected

Definition at line 349 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::resp_gf
protected

Definition at line 353 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::resu
protected

Definition at line 346 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::resu_gf
protected

Definition at line 351 of file navier_solver.hpp.

double mfem::navier::NavierSolver::rtol_hsolve = 1e-8
protected

Definition at line 394 of file navier_solver.hpp.

double mfem::navier::NavierSolver::rtol_spsolve = 1e-6
protected

Definition at line 393 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::Sp
protected

Definition at line 331 of file navier_solver.hpp.

ConstantCoefficient mfem::navier::NavierSolver::Sp_coeff
protected

Definition at line 326 of file navier_solver.hpp.

ParBilinearForm* mfem::navier::NavierSolver::Sp_form = nullptr
protected

Definition at line 304 of file navier_solver.hpp.

ParBilinearForm* mfem::navier::NavierSolver::Sp_form_lor = nullptr
protected

Definition at line 409 of file navier_solver.hpp.

OperatorHandle mfem::navier::NavierSolver::Sp_lor
protected

Definition at line 413 of file navier_solver.hpp.

CGSolver* mfem::navier::NavierSolver::SpInv = nullptr
protected

Definition at line 341 of file navier_solver.hpp.

OrthoSolver* mfem::navier::NavierSolver::SpInvOrthoPC = nullptr
protected

Definition at line 340 of file navier_solver.hpp.

HypreBoomerAMG* mfem::navier::NavierSolver::SpInvPC = nullptr
protected

Definition at line 339 of file navier_solver.hpp.

StopWatch mfem::navier::NavierSolver::sw_curlcurl
protected

Definition at line 384 of file navier_solver.hpp.

StopWatch mfem::navier::NavierSolver::sw_extrap
protected

Definition at line 384 of file navier_solver.hpp.

StopWatch mfem::navier::NavierSolver::sw_hsolve
protected

Definition at line 384 of file navier_solver.hpp.

StopWatch mfem::navier::NavierSolver::sw_setup
protected

Definition at line 384 of file navier_solver.hpp.

StopWatch mfem::navier::NavierSolver::sw_spsolve
protected

Definition at line 384 of file navier_solver.hpp.

StopWatch mfem::navier::NavierSolver::sw_step
protected

Definition at line 384 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::tmp1
protected

Definition at line 347 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::un
protected

Definition at line 346 of file navier_solver.hpp.

ParGridFunction mfem::navier::NavierSolver::un_gf
protected

Definition at line 351 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::unm1
protected

Definition at line 346 of file navier_solver.hpp.

Vector mfem::navier::NavierSolver::unm2
protected

Definition at line 346 of file navier_solver.hpp.

std::vector<VelDirichletBC_T> mfem::navier::NavierSolver::vel_dbcs
protected

Definition at line 364 of file navier_solver.hpp.

Array<int> mfem::navier::NavierSolver::vel_ess_attr
protected

Definition at line 356 of file navier_solver.hpp.

Array<int> mfem::navier::NavierSolver::vel_ess_tdof
protected

Definition at line 360 of file navier_solver.hpp.

bool mfem::navier::NavierSolver::verbose = true
protected

Enable/disable verbose output.

Definition at line 271 of file navier_solver.hpp.

FiniteElementCollection* mfem::navier::NavierSolver::vfec = nullptr
protected

Velocity \(H^1\) finite element collection.

Definition at line 289 of file navier_solver.hpp.

ParFiniteElementSpace* mfem::navier::NavierSolver::vfes = nullptr
protected

Velocity \((H^1)^d\) finite element space.

Definition at line 295 of file navier_solver.hpp.

InterpolationGridTransfer* mfem::navier::NavierSolver::vgt = nullptr
protected

Definition at line 406 of file navier_solver.hpp.

double mfem::navier::NavierSolver::volume = 0.0
protected

Definition at line 323 of file navier_solver.hpp.


The documentation for this class was generated from the following files: