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<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = vel_x
    v<<<{"description": "The y-component of velocity"}>>> = vel_y
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen]
    type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
    xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = 1
    ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = 0
    ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = 1
    nx<<<{"description": "Number of elements in the X direction"}>>> = 32
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 32
  []
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [vel_x]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
  []
  [vel_y]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
  []
  [T_fluid]
    type = INSFVEnergyVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVEnergyVariable.html"}>>>
  []
  [lambda]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = SCALAR
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
  []
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  [U]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    fv = true
  []
[]

[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
  [mag]
    type = VectorMagnitudeAux<<<{"description": "Creates a field representing the magnitude of three coupled variables using an Euclidean norm.", "href": "../../source/auxkernels/VectorMagnitudeAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = U
    x<<<{"description": "x-component of the vector"}>>> = vel_x
    y<<<{"description": "y-component of the vector"}>>> = vel_y
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [mass]
    type = INSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/INSFVMassAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
  [mean_zero_pressure]
    type = FVIntegralValueConstraint<<<{"description": "This class is used to enforce integral of phi = volume * phi_0 with a Lagrange multiplier approach.", "href": "../../source/fvkernels/FVIntegralValueConstraint.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    lambda<<<{"description": "Lagrange multiplier variable"}>>> = lambda
  []

  [u_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []

  [u_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []

  [u_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []

  [v_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []

  [v_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []

  [v_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []

  [temp_conduction]
    type = FVDiffusion<<<{"description": "Computes residual for diffusion operator for finite volume method.", "href": "../../source/fvkernels/FVDiffusion.html"}>>>
    coeff<<<{"description": "diffusion coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
  []

  [temp_advection]
    type = INSFVEnergyAdvection<<<{"description": "Advects energy, e.g. rho*cp*T. A user may still override what quantity is advected, but the default is rho*cp*T", "href": "../../source/fvkernels/INSFVEnergyAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [top_x]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
    function<<<{"description": "The exact solution function."}>>> = 'lid_function'
  []

  [no_slip_x]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right bottom'
    function<<<{"description": "The exact solution function."}>>> = 0
  []

  [no_slip_y]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
    function<<<{"description": "The exact solution function."}>>> = 0
  []

  [T_hot]
    type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/FVDirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = T_fluid
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
    value<<<{"description": "value to enforce at the boundary face"}>>> = 1
  []

  [T_cold]
    type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/FVDirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = T_fluid
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
    value<<<{"description": "value to enforce at the boundary face"}>>> = 0
  []
[]

[FunctorMaterials<<<{"href": "../../syntax/FunctorMaterials/index.html"}>>>]
  [functor_constants]
    type = ADGenericFunctorMaterial<<<{"description": "FunctorMaterial object for declaring properties that are populated by evaluation of a Functor (a constant, variable, function or functor material property) objects.", "href": "../../source/functormaterials/GenericFunctorMaterial.html"}>>>
    prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'cp k'
    prop_values<<<{"description": "The corresponding names of the functors that are going to provide the values for the variables. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${cp} ${k}'
  []
  [ins_fv]
    type = INSFVEnthalpyFunctorMaterial<<<{"description": "This is the material class used to compute enthalpy for the incompressible/weakly-compressible finite-volume implementation of the Navier-Stokes equations.", "href": "../../source/functormaterials/INSFVEnthalpyFunctorMaterial.html"}>>>
    temperature<<<{"description": "the temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_fluid'
    rho<<<{"description": "The value for the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
[]

[Functions<<<{"href": "../../syntax/Functions/index.html"}>>>]
  [lid_function]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = '4*x*(1-x)'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(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 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<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = vel_x
    v<<<{"description": "The y-component of velocity"}>>> = vel_y
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen]
    type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
    xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = 10
    ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = ${fparse -l / 2}
    ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = ${fparse l / 2}
    nx<<<{"description": "Number of elements in the X direction"}>>> = 100
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 20
  []
  uniform_refine<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 0
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [vel_x]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [vel_y]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [mass]
    type = INSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/INSFVMassAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []

  [u_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []

  [v_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [inlet-u]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    function = '${U}'
  []
  [inlet-v]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
    function = '0'
  []
  [walls-u]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    function<<<{"description": "The exact solution function."}>>> = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
    function<<<{"description": "The exact solution function."}>>> = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
    function<<<{"description": "The boundary pressure as a regular function"}>>> = '0'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = 'NEWTON'
  nl_rel_tol = 1e-12
[]

[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
  active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = FSP
  [FSP]
    type = FSP<<<{"description": "Preconditioner designed to map onto PETSc's PCFieldSplit.", "href": "../../source/preconditioners/FieldSplitPreconditioner.html"}>>>
    # It is the starting point of splitting
    topsplit<<<{"description": "Entrance to splits, the top split will specify how splits will go."}>>> = '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<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -pc_factor_shift_type'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'lu       NONZERO'
  []
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  print_linear_residuals<<<{"description": "Enable printing of linear residuals to the screen (Console)"}>>> = true
  print_nonlinear_residuals<<<{"description": "Enable printing of nonlinear residuals to the screen (Console)"}>>> = true
  [out]
    type = Exodus<<<{"description": "Object for output data in the Exodus format", "href": "../../source/outputs/Exodus.html"}>>>
    hide<<<{"description": "A list of the variables and postprocessors that should NOT be output to the Exodus file (may include Variables, ScalarVariables, and Postprocessor names)."}>>> = 'Re lin cum_lin'
  []
  [perf]
    type = PerfGraphOutput<<<{"description": "Controls output of the PerfGraph: the performance log for MOOSE", "href": "../../source/outputs/PerfGraphOutput.html"}>>>
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [Re]
    type = ParsedPostprocessor<<<{"description": "Computes a parsed expression with post-processors", "href": "../../source/postprocessors/ParsedPostprocessor.html"}>>>
    expression<<<{"description": "function expression"}>>> = '${rho} * ${l} * ${U}'
  []
  [lin]
    type = NumLinearIterations<<<{"description": "Compute the number of linear iterations.", "href": "../../source/postprocessors/NumLinearIterations.html"}>>>
  []
  [cum_lin]
    type = CumulativeValuePostprocessor<<<{"description": "Creates a cumulative sum of a Postprocessor value with time.", "href": "../../source/postprocessors/CumulativeValuePostprocessor.html"}>>>
    postprocessor<<<{"description": "The name of the postprocessor"}>>> = lin
  []
[]
(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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen]
    type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
    xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = 10
    ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = -1
    ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = 1
    nx<<<{"description": "Number of elements in the X direction"}>>> = 100
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 20
  []
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = vel_x
    v<<<{"description": "The y-component of velocity"}>>> = vel_y
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [vel_x]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [vel_y]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [mass]
    type = INSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/INSFVMassAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []

  [u_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []

  [v_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [inlet-u]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
    function = 0
  []
  [walls-u]
    type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [walls-v]
    type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [outlet_p]
    type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
    function<<<{"description": "The boundary pressure as a regular function"}>>> = 0
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(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<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = u
    v<<<{"description": "The y-component of velocity"}>>> = v
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen]
    type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
    xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = 1
    ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = 0
    ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = 4
    nx<<<{"description": "Number of elements in the X direction"}>>> = 10
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 40
  []
[]

[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
  coord_type = 'RZ'
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [u]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [v]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [mass]
    type = INSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/INSFVMassAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []

  [u_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []

  [v_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [inlet-u]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    function = 0
  []
  [inlet-v]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    function = 1
  []
  [no-slip-wall-u]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    function<<<{"description": "The exact solution function."}>>> = 0
  []
  [no-slip-wall-v]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    function<<<{"description": "The exact solution function."}>>> = 0
  []
  [outlet-p]
    type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
    function<<<{"description": "The boundary pressure as a regular function"}>>> = 0
  []
  [axis-u]
    type = INSFVSymmetryVelocityBC<<<{"description": "Implements a symmetry boundary condition for the velocity.", "href": "../../source/fvbcs/INSFVSymmetryVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = u
    v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = x
  []
  [axis-v]
    type = INSFVSymmetryVelocityBC<<<{"description": "Implements a symmetry boundary condition for the velocity.", "href": "../../source/fvbcs/INSFVSymmetryVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = u
    v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = y
  []
  [axis-p]
    type = INSFVSymmetryPressureBC<<<{"description": "Though not applied to velocity, this object ensures that the flux (velocity times the advected quantity) into a symmetry boundary is zero. When applied to pressure for the mass equation, this makes the normal velocity zero since density is constant", "href": "../../source/fvbcs/INSFVSymmetryPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [in]
    type = SideIntegralVariablePostprocessor<<<{"description": "Computes a surface integral of the specified variable", "href": "../../source/postprocessors/SideIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = v
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
  []
  [out]
    type = SideIntegralVariablePostprocessor<<<{"description": "Computes a surface integral of the specified variable", "href": "../../source/postprocessors/SideIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = v
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen]
    type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
    xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = 2
    ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = 0
    ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = 10
    nx<<<{"description": "Number of elements in the X direction"}>>> = 10
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 50
  []
[]

[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
  coord_type = 'RZ'
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = u
    v<<<{"description": "The y-component of velocity"}>>> = v
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [u]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [v]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [mass]
    type = INSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/INSFVMassAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []

  [u_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []

  [v_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [inlet-u]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    function = 0
  []
  [inlet-v]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    function = 1
  []
  [free-slip-wall-u]
    type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [free-slip-wall-v]
    type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [outlet-p]
    type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
    function<<<{"description": "The boundary pressure as a regular function"}>>> = 0
  []
  [axis-u]
    type = INSFVSymmetryVelocityBC<<<{"description": "Implements a symmetry boundary condition for the velocity.", "href": "../../source/fvbcs/INSFVSymmetryVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = u
    v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = x
  []
  [axis-v]
    type = INSFVSymmetryVelocityBC<<<{"description": "Implements a symmetry boundary condition for the velocity.", "href": "../../source/fvbcs/INSFVSymmetryVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = u
    v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = y
  []
  [axis-p]
    type = INSFVSymmetryPressureBC<<<{"description": "Though not applied to velocity, this object ensures that the flux (velocity times the advected quantity) into a symmetry boundary is zero. When applied to pressure for the mass equation, this makes the normal velocity zero since density is constant", "href": "../../source/fvbcs/INSFVSymmetryPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [in]
    type = SideIntegralVariablePostprocessor<<<{"description": "Computes a surface integral of the specified variable", "href": "../../source/postprocessors/SideIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = v
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'csv'
  []
  [out]
    type = SideIntegralVariablePostprocessor<<<{"description": "Computes a surface integral of the specified variable", "href": "../../source/postprocessors/SideIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = v
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top'
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'csv'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(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<<<{"href": "../../syntax/Problem/index.html"}>>>]
  error_on_jacobian_nonzero_reallocation = true
[]

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen_mesh]
    type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../../source/meshgenerators/FileMeshGenerator.html"}>>>
    file<<<{"description": "The filename to read."}>>> = skewed.msh
  []
  coord_type = 'XYZ'
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = vel_x
    v<<<{"description": "The y-component of velocity"}>>> = vel_y
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [vel_x]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
    face_interp_method<<<{"description": "Switch that can select between face interpolation methods."}>>> = 'skewness-corrected'
  []
  [vel_y]
    type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1
    face_interp_method<<<{"description": "Switch that can select between face interpolation methods."}>>> = 'skewness-corrected'
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
    face_interp_method<<<{"description": "Switch that can select between face interpolation methods."}>>> = 'skewness-corrected'
  []
  [lambda]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = SCALAR
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [mass]
    type = INSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/INSFVMassAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = 'skewness-corrected'
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = 'rc'
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
  [mean_zero_pressure]
    type = FVIntegralValueConstraint<<<{"description": "This class is used to enforce integral of phi = volume * phi_0 with a Lagrange multiplier approach.", "href": "../../source/fvkernels/FVIntegralValueConstraint.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    lambda<<<{"description": "Lagrange multiplier variable"}>>> = lambda
  []

  [u_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = 'skewness-corrected'
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = 'rc'
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    variable_interp_method<<<{"description": "Switch that can select between face interpolation methods for the variable."}>>> = skewness-corrected
  []
  [u_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
  [u_forcing]
    type = INSFVBodyForce<<<{"description": "Body force that contributes to the Rhie-Chow interpolation", "href": "../../source/fvkernels/INSFVBodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    functor<<<{"description": "A functor that describes the body force. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = forcing_u
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []

  [v_advection]
    type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = 'skewness-corrected'
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = 'rc'
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    variable_interp_method<<<{"description": "Switch that can select between face interpolation methods for the variable."}>>> = skewness-corrected
  []
  [v_pressure]
    type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
  [v_forcing]
    type = INSFVBodyForce<<<{"description": "Body force that contributes to the Rhie-Chow interpolation", "href": "../../source/fvkernels/INSFVBodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    functor<<<{"description": "A functor that describes the body force. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = forcing_v
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [no-slip-wall-u]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
    function<<<{"description": "The exact solution function."}>>> = '0'
  []
  [no-slip-wall-v]
    type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
    function<<<{"description": "The exact solution function."}>>> = '0'
  []
[]

[Functions<<<{"href": "../../syntax/Functions/index.html"}>>>]
  [exact_u]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = 'x^2*(1-x)^2*(2*y-6*y^2+4*y^3)'
  []
  [exact_v]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = '-y^2*(1-y)^2*(2*x-6*x^2+4*x^3)'
  []
  [exact_p]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = 'x*(1-x)-2/12'
  []
  [forcing_u]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = '-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<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = 'mu rho'
    symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = '${mu} ${rho}'
  []
  [forcing_v]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = '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<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = 'mu rho'
    symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = '${mu} ${rho}'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  [out]
    type = Exodus<<<{"description": "Object for output data in the Exodus format", "href": "../../source/outputs/Exodus.html"}>>>
    hide<<<{"description": "A list of the variables and postprocessors that should NOT be output to the Exodus file (may include Variables, ScalarVariables, and Postprocessor names)."}>>> = lambda
  []
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [h]
    type = AverageElementSize<<<{"description": "Computes the average element size.", "href": "../../source/postprocessors/AverageElementSize.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  []
  [L2u]
    type = ElementL2FunctorError<<<{"description": "Computes L2 error between an 'approximate' functor and an 'exact' functor", "href": "../../source/postprocessors/ElementL2FunctorError.html"}>>>
    approximate<<<{"description": "The approximate functor. This functor has to be an ADFunctor, like a variable or an ADFunction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_x
    exact<<<{"description": "The analytic solution to compare against. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = exact_u
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  []
  [L2v]
    type = ElementL2FunctorError<<<{"description": "Computes L2 error between an 'approximate' functor and an 'exact' functor", "href": "../../source/postprocessors/ElementL2FunctorError.html"}>>>
    approximate<<<{"description": "The approximate functor. This functor has to be an ADFunctor, like a variable or an ADFunction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_y
    exact<<<{"description": "The analytic solution to compare against. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = exact_v
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  []
  [L2p]
    approximate<<<{"description": "The approximate functor. This functor has to be an ADFunctor, like a variable or an ADFunction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
    exact<<<{"description": "The analytic solution to compare against. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = exact_p
    type = ElementL2FunctorError<<<{"description": "Computes L2 error between an 'approximate' functor and an 'exact' functor", "href": "../../source/postprocessors/ElementL2FunctorError.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = 'console csv'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  []
[]
(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

  1. 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]
  2. 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]