MFEM  v4.6.0
Finite element discretization library
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, bool provisional=false)
 Compute solution at the next time step t+dt. More...
 
ParGridFunctionGetProvisionalVelocity ()
 Return a pointer to the provisional velocity ParGridFunction. 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...
 
void UpdateTimestepHistory (double dt)
 Rotate entries in the time step and solution history arrays. More...
 
void SetMaxBDFOrder (int maxbdforder)
 Set the maximum order to use for the BDF method. More...
 
double ComputeCFL (ParGridFunction &u, double dt)
 Compute CFL. More...
 
void SetCutoffModes (int c)
 Set the number of modes to cut off in the interpolation filter. More...
 
void SetFilterAlpha (double a)
 Set the interpolation filter parameter a. More...
 

Protected Member Functions

void PrintInfo ()
 Print information 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...
 
IntegrationRules gll_rules
 
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 un_next
 
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 un_next_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 max_bdf_order = 3
 
int cur_step = 0
 
std::vector< double > dthist = {0.0, 0.0, 0.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
 
ParLORDiscretizationlor = nullptr
 
int filter_cutoff_modes = 1
 
double filter_alpha = 0.0
 
FiniteElementCollectionvfec_filter = nullptr
 
ParFiniteElementSpacevfes_filter = nullptr
 
ParGridFunction un_NM1_gf
 
ParGridFunction un_filtered_gf
 

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 incompressibility 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 141 of file navier_solver.hpp.

Constructor & Destructor Documentation

◆ NavierSolver()

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::~NavierSolver ( )

Definition at line 1165 of file navier_solver.cpp.

Member Function Documentation

◆ AddAccelTerm() [1/2]

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 1039 of file navier_solver.cpp.

◆ AddAccelTerm() [2/2]

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

Definition at line 1057 of file navier_solver.cpp.

◆ AddPresDirichletBC() [1/2]

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

Add a Dirichlet boundary condition to the pressure field.

Definition at line 1006 of file navier_solver.cpp.

◆ AddPresDirichletBC() [2/2]

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

Definition at line 1034 of file navier_solver.cpp.

◆ AddVelDirichletBC() [1/2]

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

Add a Dirichlet boundary condition to the velocity field.

Definition at line 973 of file navier_solver.cpp.

◆ AddVelDirichletBC() [2/2]

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

Definition at line 1001 of file navier_solver.cpp.

◆ ComputeCFL()

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

Compute CFL.

Definition at line 892 of file navier_solver.cpp.

◆ ComputeCurl2D()

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 791 of file navier_solver.cpp.

◆ ComputeCurl3D()

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

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

Definition at line 700 of file navier_solver.cpp.

◆ EliminateRHS()

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 667 of file navier_solver.cpp.

◆ EnableNI()

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

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

Definition at line 218 of file navier_solver.hpp.

◆ EnablePA()

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

Enable partial assembly for every operator.

Definition at line 214 of file navier_solver.hpp.

◆ GetCurrentPressure()

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

Return a pointer to the current pressure ParGridFunction.

Definition at line 192 of file navier_solver.hpp.

◆ GetCurrentVelocity()

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

Return a pointer to the current velocity ParGridFunction.

Definition at line 189 of file navier_solver.hpp.

◆ GetProvisionalVelocity()

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

Return a pointer to the provisional velocity ParGridFunction.

Definition at line 186 of file navier_solver.hpp.

◆ MeanZero()

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 638 of file navier_solver.cpp.

◆ Orthogonalize()

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 687 of file navier_solver.cpp.

◆ PrintInfo()

void NavierSolver::PrintInfo ( )
protected

Print information about the Navier version.

Definition at line 1150 of file navier_solver.cpp.

◆ PrintTimingData()

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 1115 of file navier_solver.cpp.

◆ SetCutoffModes()

void mfem::navier::NavierSolver::SetCutoffModes ( int  c)
inline

Set the number of modes to cut off in the interpolation filter.

Definition at line 274 of file navier_solver.hpp.

◆ SetFilterAlpha()

void mfem::navier::NavierSolver::SetFilterAlpha ( double  a)
inline

Set the interpolation filter parameter a.

If a is > 0, the filtering algorithm for the velocity field after every time step from [1] is used. The parameter should be 0 > >= 1.

[1] Paul Fischer, Julia Mullen (2001) Filter-based stabilization of spectral element methods

Definition at line 284 of file navier_solver.hpp.

◆ SetMaxBDFOrder()

void mfem::navier::NavierSolver::SetMaxBDFOrder ( int  maxbdforder)
inline

Set the maximum order to use for the BDF method.

Definition at line 268 of file navier_solver.hpp.

◆ SetTimeIntegrationCoefficients()

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 (or whichever order is set in SetMaxBDFOrder).

Definition at line 1062 of file navier_solver.cpp.

◆ Setup()

void NavierSolver::Setup ( double  dt)

Initialize forms, solvers and preconditioners.

Definition at line 100 of file navier_solver.cpp.

◆ Step()

void NavierSolver::Step ( double &  time,
double  dt,
int  cur_step,
bool  provisional = false 
)

Compute solution at the next time step t+dt.

This method can be called with the default value provisional which always accepts the computed time step by automatically calling UpdateTimestepHistory.

If provisional is set to true, the solution at t+dt is not accepted automatically and the application code has to call UpdateTimestepHistory and update the time variable accordingly.

The application code can check the provisional step by retrieving the GridFunction with the method GetProvisionalVelocity. If the check fails, it is possible to retry the step with a different time step by not calling UpdateTimestepHistory and calling this method with the previous time and cur_step.

The method and parameter choices are based on [1].

[1] D. Wang, S.J. Ruuth (2008) Variable step-size implicit-explicit linear multistep methods for time-dependent partial differential equations

Definition at line 367 of file navier_solver.cpp.

◆ UpdateTimestepHistory()

void NavierSolver::UpdateTimestepHistory ( double  dt)

Rotate entries in the time step and solution history arrays.

Definition at line 346 of file navier_solver.cpp.

Member Data Documentation

◆ ab1

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

Definition at line 428 of file navier_solver.hpp.

◆ ab2

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

Definition at line 429 of file navier_solver.hpp.

◆ ab3

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

Definition at line 430 of file navier_solver.hpp.

◆ accel_terms

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

Definition at line 417 of file navier_solver.hpp.

◆ bd0

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

Definition at line 424 of file navier_solver.hpp.

◆ bd1

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

Definition at line 425 of file navier_solver.hpp.

◆ bd2

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

Definition at line 426 of file navier_solver.hpp.

◆ bd3

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

Definition at line 427 of file navier_solver.hpp.

◆ cur_step

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

Definition at line 420 of file navier_solver.hpp.

◆ curlcurlu_gf

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

Definition at line 397 of file navier_solver.hpp.

◆ curlu_gf

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

Definition at line 397 of file navier_solver.hpp.

◆ D

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

Definition at line 377 of file navier_solver.hpp.

◆ D_form

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

Definition at line 351 of file navier_solver.hpp.

◆ debug

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

Enable/disable debug output.

Definition at line 311 of file navier_solver.hpp.

◆ dthist

std::vector<double> mfem::navier::NavierSolver::dthist = {0.0, 0.0, 0.0}
protected

Definition at line 421 of file navier_solver.hpp.

◆ f_form

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

Definition at line 361 of file navier_solver.hpp.

◆ Fext

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

Definition at line 391 of file navier_solver.hpp.

◆ filter_alpha

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

Definition at line 456 of file navier_solver.hpp.

◆ filter_cutoff_modes

int mfem::navier::NavierSolver::filter_cutoff_modes = 1
protected

Definition at line 455 of file navier_solver.hpp.

◆ fn

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

Definition at line 391 of file navier_solver.hpp.

◆ FText

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

Definition at line 391 of file navier_solver.hpp.

◆ FText_bdr

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

Definition at line 395 of file navier_solver.hpp.

◆ FText_bdr_form

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

Definition at line 359 of file navier_solver.hpp.

◆ FText_gf

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

Definition at line 397 of file navier_solver.hpp.

◆ FText_gfcoeff

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

Definition at line 357 of file navier_solver.hpp.

◆ G

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

Definition at line 378 of file navier_solver.hpp.

◆ g_bdr

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

Definition at line 395 of file navier_solver.hpp.

◆ g_bdr_form

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

Definition at line 363 of file navier_solver.hpp.

◆ G_form

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

Definition at line 353 of file navier_solver.hpp.

◆ gll_rules

IntegrationRules mfem::navier::NavierSolver::gll_rules
protected

Definition at line 331 of file navier_solver.hpp.

◆ H

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

Definition at line 379 of file navier_solver.hpp.

◆ H_bdfcoeff

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

Definition at line 373 of file navier_solver.hpp.

◆ H_form

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

Definition at line 355 of file navier_solver.hpp.

◆ H_lincoeff

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

Definition at line 372 of file navier_solver.hpp.

◆ HInv

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

Definition at line 389 of file navier_solver.hpp.

◆ HInvPC

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

Definition at line 388 of file navier_solver.hpp.

◆ iter_hsolve

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

Definition at line 446 of file navier_solver.hpp.

◆ iter_mvsolve

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

Definition at line 446 of file navier_solver.hpp.

◆ iter_spsolve

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

Definition at line 446 of file navier_solver.hpp.

◆ kin_vis

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

Kinematic viscosity (dimensionless).

Definition at line 329 of file navier_solver.hpp.

◆ Lext

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

Definition at line 391 of file navier_solver.hpp.

◆ Lext_gf

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

Definition at line 397 of file navier_solver.hpp.

◆ lor

ParLORDiscretization* mfem::navier::NavierSolver::lor = nullptr
protected

Definition at line 452 of file navier_solver.hpp.

◆ mass_lf

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

Linear form to compute the mass matrix in various subroutines.

Definition at line 366 of file navier_solver.hpp.

◆ max_bdf_order

int mfem::navier::NavierSolver::max_bdf_order = 3
protected

Definition at line 419 of file navier_solver.hpp.

◆ Mv

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

Definition at line 375 of file navier_solver.hpp.

◆ Mv_form

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

Definition at line 347 of file navier_solver.hpp.

◆ MvInv

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

Definition at line 382 of file navier_solver.hpp.

◆ MvInvPC

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

Definition at line 381 of file navier_solver.hpp.

◆ N

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

Definition at line 345 of file navier_solver.hpp.

◆ nlcoeff

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

Definition at line 370 of file navier_solver.hpp.

◆ numerical_integ

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

Enable/disable numerical integration rules of forms.

Definition at line 320 of file navier_solver.hpp.

◆ Nun

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

Definition at line 391 of file navier_solver.hpp.

◆ Nunm1

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

Definition at line 391 of file navier_solver.hpp.

◆ Nunm2

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

Definition at line 391 of file navier_solver.hpp.

◆ onecoeff

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

Definition at line 367 of file navier_solver.hpp.

◆ order

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

The order of the velocity and pressure space.

Definition at line 326 of file navier_solver.hpp.

◆ partial_assembly

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

Enable/disable partial assembly of forms.

Definition at line 317 of file navier_solver.hpp.

◆ pfec

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

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

Definition at line 337 of file navier_solver.hpp.

◆ pfes

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

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

Definition at line 343 of file navier_solver.hpp.

◆ pl_amg

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

Definition at line 439 of file navier_solver.hpp.

◆ pl_hsolve

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

Definition at line 438 of file navier_solver.hpp.

◆ pl_mvsolve

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

Definition at line 436 of file navier_solver.hpp.

◆ pl_spsolve

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

Definition at line 437 of file navier_solver.hpp.

◆ pmesh

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

The parallel mesh.

Definition at line 323 of file navier_solver.hpp.

◆ pn

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

Definition at line 395 of file navier_solver.hpp.

◆ pn_gf

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

Definition at line 400 of file navier_solver.hpp.

◆ pres_dbcs

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

Definition at line 414 of file navier_solver.hpp.

◆ pres_ess_attr

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

Definition at line 404 of file navier_solver.hpp.

◆ pres_ess_tdof

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

Definition at line 408 of file navier_solver.hpp.

◆ res_hsolve

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

Definition at line 449 of file navier_solver.hpp.

◆ res_mvsolve

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

Definition at line 449 of file navier_solver.hpp.

◆ res_spsolve

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

Definition at line 449 of file navier_solver.hpp.

◆ resp

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

Definition at line 395 of file navier_solver.hpp.

◆ resp_gf

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

Definition at line 400 of file navier_solver.hpp.

◆ resu

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

Definition at line 391 of file navier_solver.hpp.

◆ resu_gf

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

Definition at line 397 of file navier_solver.hpp.

◆ rtol_hsolve

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

Definition at line 443 of file navier_solver.hpp.

◆ rtol_spsolve

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

Definition at line 442 of file navier_solver.hpp.

◆ Sp

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

Definition at line 376 of file navier_solver.hpp.

◆ Sp_coeff

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

Definition at line 371 of file navier_solver.hpp.

◆ Sp_form

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

Definition at line 349 of file navier_solver.hpp.

◆ SpInv

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

Definition at line 386 of file navier_solver.hpp.

◆ SpInvOrthoPC

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

Definition at line 385 of file navier_solver.hpp.

◆ SpInvPC

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

Definition at line 384 of file navier_solver.hpp.

◆ sw_curlcurl

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

Definition at line 433 of file navier_solver.hpp.

◆ sw_extrap

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

Definition at line 433 of file navier_solver.hpp.

◆ sw_hsolve

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

Definition at line 433 of file navier_solver.hpp.

◆ sw_setup

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

Definition at line 433 of file navier_solver.hpp.

◆ sw_spsolve

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

Definition at line 433 of file navier_solver.hpp.

◆ sw_step

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

Definition at line 433 of file navier_solver.hpp.

◆ tmp1

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

Definition at line 393 of file navier_solver.hpp.

◆ un

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

Definition at line 391 of file navier_solver.hpp.

◆ un_filtered_gf

ParGridFunction mfem::navier::NavierSolver::un_filtered_gf
protected

Definition at line 460 of file navier_solver.hpp.

◆ un_gf

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

Definition at line 397 of file navier_solver.hpp.

◆ un_next

Vector mfem::navier::NavierSolver::un_next
protected

Definition at line 391 of file navier_solver.hpp.

◆ un_next_gf

ParGridFunction mfem::navier::NavierSolver::un_next_gf
protected

Definition at line 397 of file navier_solver.hpp.

◆ un_NM1_gf

ParGridFunction mfem::navier::NavierSolver::un_NM1_gf
protected

Definition at line 459 of file navier_solver.hpp.

◆ unm1

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

Definition at line 391 of file navier_solver.hpp.

◆ unm2

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

Definition at line 391 of file navier_solver.hpp.

◆ vel_dbcs

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

Definition at line 411 of file navier_solver.hpp.

◆ vel_ess_attr

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

Definition at line 403 of file navier_solver.hpp.

◆ vel_ess_tdof

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

Definition at line 407 of file navier_solver.hpp.

◆ verbose

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

Enable/disable verbose output.

Definition at line 314 of file navier_solver.hpp.

◆ vfec

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

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

Definition at line 334 of file navier_solver.hpp.

◆ vfec_filter

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

Definition at line 457 of file navier_solver.hpp.

◆ vfes

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

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

Definition at line 340 of file navier_solver.hpp.

◆ vfes_filter

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

Definition at line 458 of file navier_solver.hpp.

◆ volume

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

Definition at line 368 of file navier_solver.hpp.


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