12#ifndef MFEM_NAVIER_SOLVER_HPP
13#define MFEM_NAVIER_SOLVER_HPP
15#define NAVIER_VERSION 0.1
37 this->
attr = obj.attr;
40 this->
coeff = obj.coeff;
61 this->
attr = obj.attr;
64 this->
coeff = obj.coeff;
85 this->
attr = obj.attr;
88 this->
coeff = obj.coeff;
245 bool assume_scalar =
false);
308 int copy_interior = 0);
391 Vector fn,
un,
un_next,
unm1,
unm2,
Nun,
Nunm1,
Nunm2,
Fext,
FText,
Lext,
421 std::vector<real_t>
dthist = {0.0, 0.0, 0.0};
442#if defined(MFEM_USE_DOUBLE)
446#elif defined(MFEM_USE_SINGLE)
451#error "Only single and double precision are supported!"
Conjugate gradient method.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
A coefficient that is constant across space and time.
Square Operator for imposing essential boundary conditions using only the action, Mult(),...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
The BoomerAMG solver in hypre.
Container class for integration rules.
Pointer to an Operator of a specified type.
Solver wrapper which orthogonalizes the input and output vector.
Abstract parallel finite element space.
Class for parallel grid function.
Create and assemble a low-order refined version of a ParBilinearForm.
Class for parallel meshes.
Base class for vector Coefficients that optionally depend on time and space.
Vector coefficient defined by a vector GridFunction.
Container for an acceleration term.
AccelTerm_T(AccelTerm_T &&obj)
AccelTerm_T(Array< int > attr, VectorCoefficient *coeff)
VectorCoefficient * coeff
Transient incompressible Navier Stokes solver in a split scheme formulation.
ParGridFunction * GetProvisionalVelocity()
Return a pointer to the provisional velocity ParGridFunction.
void Setup(real_t dt)
Initialize forms, solvers and preconditioners.
ParBilinearForm * Sp_form
ParFiniteElementSpace * vfes_filter
Array< int > pres_ess_attr
real_t kin_vis
Kinematic viscosity (dimensionless).
real_t ComputeCFL(ParGridFunction &u, real_t dt)
Compute CFL.
ParGridFunction un_filtered_gf
ParLORDiscretization * lor
bool verbose
Enable/disable verbose output.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
ConstantCoefficient Sp_coeff
OrthoSolver * SpInvOrthoPC
FiniteElementCollection * vfec
Velocity finite element collection.
std::vector< PresDirichletBC_T > pres_dbcs
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.
VectorGridFunctionCoefficient * FText_gfcoeff
ParGridFunction curlcurlu_gf
void SetMaxBDFOrder(int maxbdforder)
Set the maximum order to use for the BDF method.
int order
The order of the velocity and pressure space.
bool partial_assembly
Enable/disable partial assembly of forms.
IntegrationRules gll_rules
ParGridFunction * GetCurrentPressure()
Return a pointer to the current pressure ParGridFunction.
void PrintTimingData()
Print timing summary of the solving routine.
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
FiniteElementCollection * vfec_filter
Array< int > vel_ess_attr
ConstantCoefficient onecoeff
ParGridFunction un_next_gf
void EnablePA(bool pa)
Enable partial assembly for every operator.
ParLinearForm * FText_bdr_form
ParFiniteElementSpace * pfes
Pressure finite element space.
ConstantCoefficient nlcoeff
ParBilinearForm * Mv_form
void PrintInfo()
Print information about the Navier version.
void SetFilterAlpha(real_t a)
Set the interpolation filter parameter a.
std::vector< AccelTerm_T > accel_terms
void SetCutoffModes(int c)
Set the number of modes to cut off in the interpolation filter.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
Array< int > pres_ess_tdof
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
ParLinearForm * g_bdr_form
Array< int > vel_ess_tdof
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
bool numerical_integ
Enable/disable numerical integration rules of forms.
ParMixedBilinearForm * D_form
void UpdateTimestepHistory(real_t dt)
Rotate entries in the time step and solution history arrays.
ParLinearForm * mass_lf
Linear form to compute the mass matrix in various subroutines.
bool debug
Enable/disable debug output.
ConstantCoefficient H_lincoeff
ParMesh * pmesh
The parallel mesh.
ConstantCoefficient H_bdfcoeff
void AddAccelTerm(VectorCoefficient *coeff, Array< int > &attr)
Add an acceleration term to the RHS of the equation.
std::vector< VelDirichletBC_T > vel_dbcs
void Step(real_t &time, real_t dt, int cur_step, bool provisional=false)
Compute solution at the next time step t+dt.
ParMixedBilinearForm * G_form
FiniteElementCollection * pfec
Pressure finite element collection.
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
void Orthogonalize(Vector &v)
Remove mean from a Vector.
ParFiniteElementSpace * vfes
Velocity finite element space.
ParGridFunction un_NM1_gf
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
NavierSolver(ParMesh *mesh, int order, real_t kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
std::vector< real_t > dthist
Container for a Dirichlet boundary condition of the pressure field.
PresDirichletBC_T(Array< int > attr, Coefficient *coeff)
PresDirichletBC_T(PresDirichletBC_T &&obj)
Container for a Dirichlet boundary condition of the velocity field.
VelDirichletBC_T(VelDirichletBC_T &&obj)
VelDirichletBC_T(Array< int > attr, VectorCoefficient *coeff)
VectorCoefficient * coeff
real_t(const Vector &x, real_t t) ScalarFuncT
void(const Vector &x, real_t t, Vector &u) VecFuncT
real_t u(const Vector &xvec)
std::function< real_t(const Vector &)> f(real_t mass_coeff)