November 2024
Main capabilities:
Fluid types:
Incompressible
Weakly-compressible
Compressible
Flow regimes:
Laminar
Limited turbulence handling (mixing length)
Flow types:
Free flows
Flows in porous media describing homogenized structures
Utilizes the following techniques:
Stabilized Continuous Finite Element Discretization (maintained, but not intensely developed) (Peterson et al. (2018))
Finite Volume Discretization (current development direction) (Lindsay et al. (2023))
Conservation of mass:
Conservation of (linear) momentum:
- porosity
- real (or interstitial) velocity
- superficial velocity
- density
- effective dynamic viscosity (laminar + turbulent)
- pressure
- friction tensor (using correlations)
If then and we get the original Navier-Stokes equations back.
These are supplemented by the conservation of energy:
Assumptions: Weakly-compressible fluid, kinetic energy of the fluid can be neglected
- specific heat
- fluid temperature
- effective thermal conductivity
- volumetric heat-transfer coefficient (using correlations)
- solid temperature
- volumetric source term
The equations are extended with initial and boundary conditions:
Few examples:
Velocity Inlet:
Pressure Outlet:
No-slip, insulated walls:
Heat flux on the walls:
Free slip on walls:
...
Good to know before navigating source code and documentation:
NS - Navier-Stokes
FV - Finite Volume
I - Incompressible
WC - Weakly-Compressible
C - Compressible
P - Porous-Medium
Few Examples:
PINSFV: Porous-Medium Incompressible Navier-Stokes Finite Volume
WCNSFV: Weakly-Compressible Navier-Stokes Finite Volume
The building blocks in MOOSE for terms in the PDEs are Kernels for FE or FVKernels for FV:
The building blocks in MOOSE for boundary conditions are BCs for FE or FVBCs for FV:
Inlet Velocity:
Inlet Temperature:
Outlet pressure:
Heatflux:
Freeslip:
...
Let us consider the following (fictional) material properties:
(we are using incompressible formulation in this case)
mu = 1.1
rho = 1.1
l = 2
U = 1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = ${fparse -l / 2}
ymax = ${fparse l / 2}
nx = 100
ny = 20
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 1
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1
[]
[pressure]
type = INSFVPressureVariable
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
function = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
function = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '0'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-12
[]
[Preconditioning]
active = FSP
[FSP]
type = FSP
# It is the starting point of splitting
topsplit = 'up' # 'up' should match the following block name
[up]
splitting = 'u p' # 'u' and 'p' are the names of subsolvers
splitting_type = schur
# Splitting type is set as schur, because the pressure part of Stokes-like systems
# is not diagonally dominant. CAN NOT use additive, multiplicative and etc.
#
# Original system:
#
# | Auu Aup | | u | = | f_u |
# | Apu 0 | | p | | f_p |
#
# is factorized into
#
# |I 0 | | Auu 0| | I Auu^{-1}*Aup | | u | = | f_u |
# |Apu*Auu^{-1} I | | 0 -S| | 0 I | | p | | f_p |
#
# where
#
# S = Apu*Auu^{-1}*Aup
#
# The preconditioning is accomplished via the following steps
#
# (1) p* = f_p - Apu*Auu^{-1}f_u,
# (2) p = (-S)^{-1} p*
# (3) u = Auu^{-1}(f_u-Aup*p)
petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition -ksp_gmres_restart -ksp_rtol -ksp_type'
petsc_options_value = 'full selfp 300 1e-4 fgmres'
[]
[u]
vars = 'vel_x vel_y'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_type -ksp_rtol -ksp_gmres_restart -ksp_pc_side'
petsc_options_value = 'hypre boomeramg gmres 5e-1 300 right'
[]
[p]
vars = 'pressure'
petsc_options_iname = '-ksp_type -ksp_gmres_restart -ksp_rtol -pc_type -ksp_pc_side'
petsc_options_value = 'gmres 300 5e-1 jacobi right'
[]
[]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
print_linear_residuals = true
print_nonlinear_residuals = true
[out]
type = Exodus
hide = 'Re lin cum_lin'
[]
[perf]
type = PerfGraphOutput
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${l} * ${U}'
[]
[lin]
type = NumLinearIterations
[]
[cum_lin]
type = CumulativeValuePostprocessor
postprocessor = lin
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-no-slip.i)When using linear interpolation for the advecting velocity and the pressure gradient, we encounter a checker-boarding effect for the pressure field. Neighboring cell center solution values are decoupled, which is undesirable.
Solution for this issue: Rhie-Chow interpolation for the advecting velocity. For more information on the issue and the derivation of the interpolation method, see Moukalled et al. (2016).
A simplified syntax has been developed, it relies on the Action
system in MOOSE. For the documentation of the corresponding action, click here!
mu = 1.1
rho = 1.1
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = -1
ymax = 1
nx = 100
ny = 20
[]
[]
[Modules]
[NavierStokesFV]
compressibility = 'incompressible'
density = ${rho}
dynamic_viscosity = ${mu}
initial_velocity = '1 1 0'
initial_pressure = 0.0
inlet_boundaries = 'left'
momentum_inlet_types = 'fixed-velocity'
momentum_inlet_function = '1 0'
wall_boundaries = 'top bottom'
momentum_wall_types = 'noslip noslip'
outlet_boundaries = 'right'
momentum_outlet_types = 'fixed-pressure'
pressure_function = '0'
mass_advection_interpolation = 'average'
momentum_advection_interpolation = 'average'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-no-slip-action.i)If possible, use the Navier Stokes Finite Volume Action syntax
Use Rhie-Chow interpolation for the advecting velocity
Other interpolation techniques may lead to checker-boarding/instability
Start with first-order advected-interpolation schemes (e.g. upwind)
Make sure that the pressure is pinned for incompressible/weakly-compressible simulations in closed systems
For monolithic solvers (the default at the moment) use a variant of LU preconditioner
For complex monolithic systems monitor the residuals of every variable
Try to keep the number of elements relatively low (limited scaling with LU preconditioner). Larger problems can be solved in a reasonable wall-time as segregated-solve techniques are introduced
Try to utilize porosity_smoothing_layers
or a higher value for the consistent_scaling
if you encounter oscillatory behavior in case of simulations using porous medium
Because the complexity of the LU preconditioner is in the worst case scenario, its performance in general largely depends on the matrix sparsity pattern, and it is much more expensive compared with an iterative approach.