Incompressible Finite Volume Navier Stokes
MOOSE's Incompressible Navier Stokes Finite Volume (INSFV) implementation uses a colocated grid. To suppress the checkerboard pattern in the pressure field, INSFV
objects support a Rhie-Chow interpolation for the velocity. Users can get a feel for INSFV by looking at some tests. In addition, to ease the burden of preparing long input files, the NavierStokesFV action syntax can also be used to set up INSFV simulations.
Lid Driven Cavity Flow
This example solves the INS equations for mass, momentum, and energy in a closed cavity. Because there are no inlet or outlet boundaries, one pressure degree of freedom must be constrained in order to eliminate the nullspace. This is done using the FVScalarLagrangeMultiplier
object which implements the mean-zero pressure approach. The finite element theory of the mean-zero approach is described here. The finite volume implementation is completed by simply substituting unity for the test functions.
mu = 1
rho = 1
k = .01
cp = 1
velocity_interp_method = 'rc'
advected_interp_method = 'average'
[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 = 1
ymin = 0
ymax = 1
nx = 32
ny = 32
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
[]
[vel_y]
type = INSFVVelocityVariable
[]
[pressure]
type = INSFVPressureVariable
[]
[T_fluid]
type = INSFVEnergyVariable
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[AuxVariables]
[U]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[]
[AuxKernels]
[mag]
type = VectorMagnitudeAux
variable = U
x = vel_x
y = vel_y
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[mean_zero_pressure]
type = FVIntegralValueConstraint
variable = pressure
lambda = lambda
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T_fluid
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T_fluid
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = 'lid_function'
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[T_hot]
type = FVDirichletBC
variable = T_fluid
boundary = 'bottom'
value = 1
[]
[T_cold]
type = FVDirichletBC
variable = T_fluid
boundary = 'top'
value = 0
[]
[]
[Materials]
[functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
[]
[ins_fv]
type = INSFVEnthalpyMaterial
temperature = 'T_fluid'
rho = ${rho}
[]
[]
[Functions]
[lid_function]
type = ParsedFunction
expression = '4*x*(1-x)'
[]
[]
[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
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/lid-driven-with-energy.i)Channel Flow
There are examples of both no-slip and free-slip channel flow. The only difference between the two inputs is in the boundary condition block. For the no-slip case, a Dirichlet 0 condition is applied to the tangential velocity (in this case u
); this condition is not applied for free-slip. A Dirichlet 0 condition is applied on the pressure at the outlet boundary for both inputs. This is necessary, in particular for the no-slip case, because if a Dirichlet condition is not applied for a finite volume variable at a boundary, then a zero-gradient condition is implicitly applied. This is an invalid boundary condition for no-slip; in fully developed flow we know that pressure, while constant across the channel cross-section, decreases in the direction of flow. (Otherwise without that pressure driving-force, how are you going to drive flow with the walls slowing you down?)
One may reasonably ask why we implicitly apply a zero normal-gradient condition when Dirichlet conditions are not applied. This is so that FVFluxKernels
executed along a boundary have a value for the field in the neighboring ghost cell. FVFluxKernels
are always executed along a boundary if FVDirichletBCs
are active; their execution ensures that that the Dirichlet condition is weakly enforced. When FVDirichletBCs
are not active, FVFluxKernels
may still be forced to execute along a boundary by specifying force_boundary_execution = true
in the respective block. Forcing execution of a FVFluxKernel
on a boundary when no Dirichlet BC is present, e.g. with an implicit application of the zero normal-gradient condition, is typically identical to applying a corresponding FVFluxBC
because the latter only has access to the cell center value adjacent to the boundary. By directly using a cell-center value in a FVFluxBC
(see for example FVConstantScalarOutflowBC
), you are implicitly applying the same zero normal-gradient condition. In the future (see MOOSE issue #16169), we hope to be able to apply more appropriate outflow conditions.
no slip
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
function = '${rho} * ${l} * ${U}'
pp_names = ''
[]
[lin]
type = NumLinearIterations
[]
[cum_lin]
type = CumulativeValuePostprocessor
postprocessor = lin
[]
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-no-slip.i)free slip
mu = 1.1
rho = 1.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = -1
ymax = 1
nx = 100
ny = 20
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[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 = '1'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
function = 0
[]
[walls-u]
type = INSFVNaturalFreeSlipBC
boundary = 'top bottom'
variable = vel_x
momentum_component = 'x'
[]
[walls-v]
type = INSFVNaturalFreeSlipBC
boundary = 'top bottom'
variable = vel_y
momentum_component = 'y'
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = 0
[]
[]
[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
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc.i)Axisymmetric Channel Flow
Channel flow in axisymmetric coordinates is also implemented. Below is an example of solving the no-slip problem using an average
interpolation for the velocity:
mu = 1
rho = 1
advected_interp_method = 'average'
velocity_interp_method = 'average'
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 4
nx = 10
ny = 40
[]
[]
[Problem]
coord_type = 'RZ'
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = 1
[]
[v]
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 = u
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = u
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = u
momentum_component = 'x'
pressure = pressure
[]
[v_advection]
type = INSFVMomentumAdvection
variable = v
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = v
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = v
momentum_component = 'y'
pressure = pressure
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'bottom'
variable = u
function = 0
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'bottom'
variable = v
function = 1
[]
[no-slip-wall-u]
type = INSFVNoSlipWallBC
boundary = 'right'
variable = u
function = 0
[]
[no-slip-wall-v]
type = INSFVNoSlipWallBC
boundary = 'right'
variable = v
function = 0
[]
[outlet-p]
type = INSFVOutletPressureBC
boundary = 'top'
variable = pressure
function = 0
[]
[axis-u]
type = INSFVSymmetryVelocityBC
boundary = 'left'
variable = u
u = u
v = v
mu = ${mu}
momentum_component = x
[]
[axis-v]
type = INSFVSymmetryVelocityBC
boundary = 'left'
variable = v
u = u
v = v
mu = ${mu}
momentum_component = y
[]
[axis-p]
type = INSFVSymmetryPressureBC
boundary = 'left'
variable = pressure
[]
[]
[Postprocessors]
[in]
type = SideIntegralVariablePostprocessor
variable = v
boundary = 'bottom'
[]
[out]
type = SideIntegralVariablePostprocessor
variable = v
boundary = 'top'
[]
[]
[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
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/cylindrical/2d-average-no-slip.i)The Rhie-Chow interpolation can be used simply by specifying velocity_interp_method='rc'
either in the input file or from the command line. An axisymmetric example with free slip conditions, using the Rhie-Chow interpolation is shown below:
mu = 1.1
rho = 1.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 2
ymin = 0
ymax = 10
nx = 10
ny = 50
[]
[]
[Problem]
coord_type = 'RZ'
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = 1
[]
[v]
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 = u
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = u
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = u
momentum_component = 'x'
pressure = pressure
[]
[v_advection]
type = INSFVMomentumAdvection
variable = v
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = v
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = v
momentum_component = 'y'
pressure = pressure
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'bottom'
variable = u
function = 0
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'bottom'
variable = v
function = 1
[]
[free-slip-wall-u]
type = INSFVNaturalFreeSlipBC
boundary = 'right'
variable = u
momentum_component = 'x'
[]
[free-slip-wall-v]
type = INSFVNaturalFreeSlipBC
boundary = 'right'
variable = v
momentum_component = 'y'
[]
[outlet-p]
type = INSFVOutletPressureBC
boundary = 'top'
variable = pressure
function = 0
[]
[axis-u]
type = INSFVSymmetryVelocityBC
boundary = 'left'
variable = u
u = u
v = v
mu = ${mu}
momentum_component = x
[]
[axis-v]
type = INSFVSymmetryVelocityBC
boundary = 'left'
variable = v
u = u
v = v
mu = ${mu}
momentum_component = y
[]
[axis-p]
type = INSFVSymmetryPressureBC
boundary = 'left'
variable = pressure
[]
[]
[Postprocessors]
[in]
type = SideIntegralVariablePostprocessor
variable = v
boundary = 'bottom'
outputs = 'csv'
[]
[out]
type = SideIntegralVariablePostprocessor
variable = v
boundary = 'top'
outputs = 'csv'
[]
[]
[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
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/cylindrical/2d-rc-slip.i)Skewness-correction
The skewness-correction of different variables can be enabled by defining the face_interp_method=skewness-corrected
parameter for the INSFVVariables and selecting it as an option in the advection kernels. It has proven to increase accuracy on unstructured grids. In case of a skewed 2D triangulation, it increases the order of the error from to for velocity, and from to for pressure. For an example see:
mu = 1.0
rho = 1.0
[Problem]
error_on_jacobian_nonzero_reallocation = true
[]
[Mesh]
[gen_mesh]
type = FileMeshGenerator
file = skewed.msh
[]
coord_type = 'XYZ'
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 1
face_interp_method = 'skewness-corrected'
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1
face_interp_method = 'skewness-corrected'
[]
[pressure]
type = INSFVPressureVariable
face_interp_method = 'skewness-corrected'
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = 'skewness-corrected'
velocity_interp_method = 'rc'
rho = ${rho}
[]
[mean_zero_pressure]
type = FVIntegralValueConstraint
variable = pressure
lambda = lambda
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = 'skewness-corrected'
velocity_interp_method = 'rc'
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
[]
[u_forcing]
type = INSFVBodyForce
variable = vel_x
functor = forcing_u
momentum_component = 'x'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = 'skewness-corrected'
velocity_interp_method = 'rc'
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
[]
[v_forcing]
type = INSFVBodyForce
variable = vel_y
functor = forcing_v
momentum_component = 'y'
[]
[]
[FVBCs]
[no-slip-wall-u]
type = INSFVNoSlipWallBC
boundary = 'left right top bottom'
variable = vel_x
function = '0'
[]
[no-slip-wall-v]
type = INSFVNoSlipWallBC
boundary = 'left right top bottom'
variable = vel_y
function = '0'
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
expression = 'x^2*(1-x)^2*(2*y-6*y^2+4*y^3)'
[]
[exact_v]
type = ParsedFunction
expression = '-y^2*(1-y)^2*(2*x-6*x^2+4*x^3)'
[]
[exact_p]
type = ParsedFunction
expression = 'x*(1-x)-2/12'
[]
[forcing_u]
type = ParsedFunction
expression = '-4*mu/rho*(-1+2*y)*(y^2-6*x*y^2+6*x^2*y^2-y+6*x*y-6*x^2*y+3*x^2-6*x^3+3*x^4)+1-2*x+4*x^3'
'*y^2*(2*y^2-2*y+1)*(y-1)^2*(-1+2*x)*(x-1)^3'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[forcing_v]
type = ParsedFunction
expression = '4*mu/rho*(-1+2*x)*(x^2-6*y*x^2+6*x^2*y^2-x+6*x*y-6*x*y^2+3*y^2-6*y^3+3*y^4)+4*y^3*x^2*(2'
'*x^2-2*x+1)*(x-1)^2*(-1+2*y)*(y-1)^3'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-8
[]
[Outputs]
[out]
type = Exodus
hide = lambda
[]
csv = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2u]
type = ElementL2FunctorError
approximate = vel_x
exact = exact_u
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2v]
type = ElementL2FunctorError
approximate = vel_y
exact = exact_v
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2p]
approximate = pressure
exact = exact_p
type = ElementL2FunctorError
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/mms/skew-correction/skewed-vortex.i)Cell-centered vector field reconstruction using face fluxes
Weller's method
For a detailed description on the origins and properties of this method, see Weller (2014) and Aguerre et al. (2018). In short, this reconstruction can be used to obtain cell-center velocities and pressure gradients based on face fluxes and normal pressure face gradients, respectively. It exhibits a second order convergence () spatially. The expression for the vector value at the cell centers is the following:
(1)where denotes the surface normal, the surface area and the face flux. Latter can be computed by the face vector values by .
References
- Horacio J Aguerre, Cesar I Pairetti, Cesar M Venier, Santiago Márquez Damián, and Norberto M Nigro.
An oscillation-free flow solver based on flux reconstruction.
Journal of Computational Physics, 365:135–148, 2018.[BibTeX]
- Hilary Weller.
Non-orthogonal version of the arbitrary polygonal c-grid and a new diamond grid.
Geoscientific Model Development, 7(3):779–797, 2014.[BibTeX]