12 #ifndef MFEM_NAVIER_SOLVER_HPP
13 #define MFEM_NAVIER_SOLVER_HPP
15 #define NAVIER_VERSION 0.1
32 : attr(attr), coeff(coeff)
38 this->
attr = obj.attr;
41 this->
coeff = obj.coeff;
56 : attr(attr), coeff(coeff)
62 this->
attr = obj.attr;
65 this->
coeff = obj.coeff;
80 : attr(attr), coeff(coeff)
86 this->
attr = obj.attr;
89 this->
coeff = obj.coeff;
160 void Setup(
double dt);
222 bool assume_scalar =
false);
265 int copy_interior = 0);
346 Vector fn,
un,
unm1,
unm2,
Nun,
Nunm1,
Nunm2,
Fext,
FText,
Lext,
resu;
ParBilinearForm * Mv_form
VectorCoefficient * coeff
Container for a Dirichlet boundary condition of the velocity field.
Array< int > pres_ess_attr
void PrintTimingData()
Print timing summary of the solving routine.
void Orthogonalize(Vector &v)
Remove mean from a Vector.
Conjugate gradient method.
FiniteElementCollection * pfec
Pressure finite element collection.
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.
ParBilinearForm * Mv_form_lor
Container for a Dirichlet boundary condition of the pressure field.
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
void PrintInfo()
Print informations about the Navier version.
FiniteElementCollection * pfec_lor
ConstantCoefficient nlcoeff
Pointer to an Operator of a specified type.
void(const Vector &x, double t, Vector &u) VecFuncT
double(const Vector &x, double t) ScalarFuncT
ParLinearForm * FText_bdr_form
std::vector< PresDirichletBC_T > pres_dbcs
bool numerical_integ
Enable/disable numerical integration rules of forms.
Abstract parallel finite element space.
ConstantCoefficient H_lincoeff
ParLinearForm * g_bdr_form
Container for an acceleration term.
NavierSolver(ParMesh *mesh, int order, double kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
AccelTerm_T(Array< int > attr, VectorCoefficient *coeff)
PresDirichletBC_T(PresDirichletBC_T &&obj)
ParMixedBilinearForm * G_form
void Setup(double dt)
Initialize forms, solvers and preconditioners.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
VelDirichletBC_T(VelDirichletBC_T &&obj)
Array< int > vel_ess_attr
Array< int > vel_ess_tdof
bool partial_assembly
Enable/disable partial assembly of forms.
ConstantCoefficient Sp_coeff
void Step(double &time, double dt, int cur_step)
Compute solution at the next time step t+dt.
The BoomerAMG solver in hypre.
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
FiniteElementCollection * vfec
Velocity finite element collection.
ParFiniteElementSpace * pfes_lor
InterpolationGridTransfer * vgt
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
ParGridFunction * GetCurrentPressure()
Return a pointer to the current pressure ParGridFunction.
ConstantCoefficient onecoeff
ParBilinearForm * Sp_form
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
void EnablePA(bool pa)
Enable partial assembly for every operator.
ParBilinearForm * Sp_form_lor
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
ParBilinearForm * H_form_lor
ParFiniteElementSpace * vfes
Velocity finite element space.
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
ParLinearForm * mass_lf
Linear form to compute the mass matrix in various subroutines.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
std::vector< AccelTerm_T > accel_terms
OrthoSolver * SpInvOrthoPC
bool verbose
Enable/disable verbose output.
InterpolationGridTransfer * pgt
double kin_vis
Kinematic viscosity (dimensionless).
Array< int > pres_ess_tdof
std::vector< VelDirichletBC_T > vel_dbcs
VectorGridFunctionCoefficient * FText_gfcoeff
VectorCoefficient * coeff
ParFiniteElementSpace * pfes
Pressure finite element space.
PresDirichletBC_T(Array< int > attr, Coefficient *coeff)
AccelTerm_T(AccelTerm_T &&obj)
ParMixedBilinearForm * D_form
Vector coefficient defined by a vector GridFunction.
ParGridFunction curlcurlu_gf
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
Class for parallel grid function.
Square Operator for imposing essential boundary conditions using only the action, Mult()...
bool debug
Enable/disable debug output.
Class for parallel meshes.
Transient incompressible Navier Stokes solver in a split scheme formulation.
Solver wrapper which orthogonalizes the input and output vector.
ParMesh * pmesh
The parallel mesh.
void AddAccelTerm(VectorCoefficient *coeff, Array< int > &attr)
Add an acceleration term to the RHS of the equation.
int order
The order of the velocity and pressure space.
ConstantCoefficient H_bdfcoeff
double f(const Vector &p)
VelDirichletBC_T(Array< int > attr, VectorCoefficient *coeff)