MFEM v4.7.0
Finite element discretization library
|
Transient incompressible Navier Stokes solver in a split scheme formulation. More...
#include <navier_solver.hpp>
Public Member Functions | |
NavierSolver (ParMesh *mesh, int order, real_t kin_vis) | |
Initialize data structures, set FE space order and kinematic viscosity. | |
void | Setup (real_t dt) |
Initialize forms, solvers and preconditioners. | |
void | Step (real_t &time, real_t dt, int cur_step, bool provisional=false) |
Compute solution at the next time step t+dt. | |
ParGridFunction * | GetProvisionalVelocity () |
Return a pointer to the provisional velocity ParGridFunction. | |
ParGridFunction * | GetCurrentVelocity () |
Return a pointer to the current velocity ParGridFunction. | |
ParGridFunction * | GetCurrentPressure () |
Return a pointer to the current pressure ParGridFunction. | |
void | AddVelDirichletBC (VectorCoefficient *coeff, Array< int > &attr) |
Add a Dirichlet boundary condition to the velocity field. | |
void | AddVelDirichletBC (VecFuncT *f, Array< int > &attr) |
void | AddPresDirichletBC (Coefficient *coeff, Array< int > &attr) |
Add a Dirichlet boundary condition to the pressure field. | |
void | AddPresDirichletBC (ScalarFuncT *f, Array< int > &attr) |
void | AddAccelTerm (VectorCoefficient *coeff, Array< int > &attr) |
Add an acceleration term to the RHS of the equation. | |
void | AddAccelTerm (VecFuncT *f, Array< int > &attr) |
void | EnablePA (bool pa) |
Enable partial assembly for every operator. | |
void | EnableNI (bool ni) |
void | PrintTimingData () |
Print timing summary of the solving routine. | |
~NavierSolver () | |
void | ComputeCurl2D (ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false) |
Compute \(\nabla \times \nabla \times u\) for \(u \in (H^1)^2\). | |
void | ComputeCurl3D (ParGridFunction &u, ParGridFunction &cu) |
Compute \(\nabla \times \nabla \times u\) for \(u \in (H^1)^3\). | |
void | Orthogonalize (Vector &v) |
Remove mean from a Vector. | |
void | MeanZero (ParGridFunction &v) |
Remove the mean from a ParGridFunction. | |
void | UpdateTimestepHistory (real_t dt) |
Rotate entries in the time step and solution history arrays. | |
void | SetMaxBDFOrder (int maxbdforder) |
Set the maximum order to use for the BDF method. | |
real_t | ComputeCFL (ParGridFunction &u, real_t dt) |
Compute CFL. | |
void | SetCutoffModes (int c) |
Set the number of modes to cut off in the interpolation filter. | |
void | SetFilterAlpha (real_t a) |
Set the interpolation filter parameter a. | |
Protected Member Functions | |
void | PrintInfo () |
Print information about the Navier version. | |
void | SetTimeIntegrationCoefficients (int step) |
Update the EXTk/BDF time integration coefficient. | |
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. | |
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.
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.
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 1167 of file navier_solver.cpp.
Definition at line 1058 of file navier_solver.cpp.
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 1040 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 1007 of file navier_solver.cpp.
void NavierSolver::AddPresDirichletBC | ( | ScalarFuncT * | f, |
Array< int > & | attr ) |
Definition at line 1035 of file navier_solver.cpp.
Definition at line 1002 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 974 of file navier_solver.cpp.
real_t NavierSolver::ComputeCFL | ( | ParGridFunction & | u, |
real_t | dt ) |
Compute CFL.
Definition at line 893 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 792 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 701 of file navier_solver.cpp.
|
protected |
Eliminate essential BCs in an Operator and apply to RHS.
Definition at line 667 of file navier_solver.cpp.
|
inline |
Enable numerical integration rules. This means collocated quadrature at the nodal points.
Definition at line 218 of file navier_solver.hpp.
|
inline |
Enable partial assembly for every operator.
Definition at line 214 of file navier_solver.hpp.
|
inline |
Return a pointer to the current pressure ParGridFunction.
Definition at line 192 of file navier_solver.hpp.
|
inline |
Return a pointer to the current velocity ParGridFunction.
Definition at line 189 of file navier_solver.hpp.
|
inline |
Return a pointer to the provisional velocity ParGridFunction.
Definition at line 186 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 638 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 687 of file navier_solver.cpp.
|
protected |
Print information about the Navier version.
Definition at line 1152 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
The second row shows a proportion of a column relative to the whole time step.
Definition at line 1116 of file navier_solver.cpp.
|
inline |
Set the number of modes to cut off in the interpolation filter.
Definition at line 274 of file navier_solver.hpp.
|
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.
|
inline |
Set the maximum order to use for the BDF method.
Definition at line 268 of file navier_solver.hpp.
|
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 1063 of file navier_solver.cpp.
void NavierSolver::Setup | ( | real_t | dt | ) |
Initialize forms, solvers and preconditioners.
Definition at line 100 of file navier_solver.cpp.
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.
void NavierSolver::UpdateTimestepHistory | ( | real_t | dt | ) |
Rotate entries in the time step and solution history arrays.
Definition at line 346 of file navier_solver.cpp.
|
protected |
Definition at line 428 of file navier_solver.hpp.
|
protected |
Definition at line 429 of file navier_solver.hpp.
|
protected |
Definition at line 430 of file navier_solver.hpp.
|
protected |
Definition at line 417 of file navier_solver.hpp.
|
protected |
Definition at line 424 of file navier_solver.hpp.
|
protected |
Definition at line 425 of file navier_solver.hpp.
|
protected |
Definition at line 426 of file navier_solver.hpp.
|
protected |
Definition at line 427 of file navier_solver.hpp.
|
protected |
Definition at line 420 of file navier_solver.hpp.
|
protected |
Definition at line 397 of file navier_solver.hpp.
|
protected |
Definition at line 397 of file navier_solver.hpp.
|
protected |
Definition at line 377 of file navier_solver.hpp.
|
protected |
Definition at line 351 of file navier_solver.hpp.
|
protected |
Enable/disable debug output.
Definition at line 311 of file navier_solver.hpp.
|
protected |
Definition at line 421 of file navier_solver.hpp.
|
protected |
Definition at line 361 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 468 of file navier_solver.hpp.
|
protected |
Definition at line 467 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 395 of file navier_solver.hpp.
|
protected |
Definition at line 359 of file navier_solver.hpp.
|
protected |
Definition at line 397 of file navier_solver.hpp.
|
protected |
Definition at line 357 of file navier_solver.hpp.
|
protected |
Definition at line 378 of file navier_solver.hpp.
|
protected |
Definition at line 395 of file navier_solver.hpp.
|
protected |
Definition at line 363 of file navier_solver.hpp.
|
protected |
Definition at line 353 of file navier_solver.hpp.
|
protected |
Definition at line 331 of file navier_solver.hpp.
|
protected |
Definition at line 379 of file navier_solver.hpp.
|
protected |
Definition at line 373 of file navier_solver.hpp.
|
protected |
Definition at line 355 of file navier_solver.hpp.
|
protected |
Definition at line 372 of file navier_solver.hpp.
|
protected |
Definition at line 389 of file navier_solver.hpp.
|
protected |
Definition at line 388 of file navier_solver.hpp.
|
protected |
Definition at line 458 of file navier_solver.hpp.
|
protected |
Definition at line 458 of file navier_solver.hpp.
|
protected |
Definition at line 458 of file navier_solver.hpp.
|
protected |
Kinematic viscosity (dimensionless).
Definition at line 329 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 397 of file navier_solver.hpp.
|
protected |
Definition at line 464 of file navier_solver.hpp.
|
protected |
Linear form to compute the mass matrix in various subroutines.
Definition at line 366 of file navier_solver.hpp.
|
protected |
Definition at line 419 of file navier_solver.hpp.
|
protected |
Definition at line 375 of file navier_solver.hpp.
|
protected |
Definition at line 347 of file navier_solver.hpp.
|
protected |
Definition at line 382 of file navier_solver.hpp.
|
protected |
Definition at line 381 of file navier_solver.hpp.
|
protected |
Definition at line 345 of file navier_solver.hpp.
|
protected |
Definition at line 370 of file navier_solver.hpp.
|
protected |
Enable/disable numerical integration rules of forms.
Definition at line 320 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 367 of file navier_solver.hpp.
|
protected |
The order of the velocity and pressure space.
Definition at line 326 of file navier_solver.hpp.
|
protected |
Enable/disable partial assembly of forms.
Definition at line 317 of file navier_solver.hpp.
|
protected |
Pressure \(H^1\) finite element collection.
Definition at line 337 of file navier_solver.hpp.
|
protected |
Pressure \(H^1\) finite element space.
Definition at line 343 of file navier_solver.hpp.
|
protected |
Definition at line 439 of file navier_solver.hpp.
|
protected |
Definition at line 438 of file navier_solver.hpp.
|
protected |
Definition at line 436 of file navier_solver.hpp.
|
protected |
Definition at line 437 of file navier_solver.hpp.
|
protected |
The parallel mesh.
Definition at line 323 of file navier_solver.hpp.
|
protected |
Definition at line 395 of file navier_solver.hpp.
|
protected |
Definition at line 400 of file navier_solver.hpp.
|
protected |
Definition at line 414 of file navier_solver.hpp.
|
protected |
Definition at line 404 of file navier_solver.hpp.
|
protected |
Definition at line 408 of file navier_solver.hpp.
|
protected |
Definition at line 461 of file navier_solver.hpp.
|
protected |
Definition at line 461 of file navier_solver.hpp.
|
protected |
Definition at line 461 of file navier_solver.hpp.
|
protected |
Definition at line 395 of file navier_solver.hpp.
|
protected |
Definition at line 400 of file navier_solver.hpp.
|
protected |
Definition at line 392 of file navier_solver.hpp.
|
protected |
Definition at line 398 of file navier_solver.hpp.
|
protected |
Definition at line 445 of file navier_solver.hpp.
|
protected |
Definition at line 443 of file navier_solver.hpp.
|
protected |
Definition at line 444 of file navier_solver.hpp.
|
protected |
Definition at line 376 of file navier_solver.hpp.
|
protected |
Definition at line 371 of file navier_solver.hpp.
|
protected |
Definition at line 349 of file navier_solver.hpp.
|
protected |
Definition at line 386 of file navier_solver.hpp.
|
protected |
Definition at line 385 of file navier_solver.hpp.
|
protected |
Definition at line 384 of file navier_solver.hpp.
|
protected |
Definition at line 433 of file navier_solver.hpp.
|
protected |
Definition at line 433 of file navier_solver.hpp.
|
protected |
Definition at line 433 of file navier_solver.hpp.
|
protected |
Definition at line 433 of file navier_solver.hpp.
|
protected |
Definition at line 433 of file navier_solver.hpp.
|
protected |
Definition at line 433 of file navier_solver.hpp.
|
protected |
Definition at line 393 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 472 of file navier_solver.hpp.
|
protected |
Definition at line 397 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 397 of file navier_solver.hpp.
|
protected |
Definition at line 471 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 391 of file navier_solver.hpp.
|
protected |
Definition at line 411 of file navier_solver.hpp.
|
protected |
Definition at line 403 of file navier_solver.hpp.
|
protected |
Definition at line 407 of file navier_solver.hpp.
|
protected |
Enable/disable verbose output.
Definition at line 314 of file navier_solver.hpp.
|
protected |
Velocity \(H^1\) finite element collection.
Definition at line 334 of file navier_solver.hpp.
|
protected |
Definition at line 469 of file navier_solver.hpp.
|
protected |
Velocity \((H^1)^d\) finite element space.
Definition at line 340 of file navier_solver.hpp.
|
protected |
Definition at line 470 of file navier_solver.hpp.
|
protected |
Definition at line 368 of file navier_solver.hpp.