NavierStokesFV System

warningwarning

This syntax is deprecated. Please refer to the section on how to transition to the new Physics syntax for guidance on how to use the current syntax.

Overview

The NavierStokesFV system is dedicated to decrease the effort required by the user to prepare simulations which need to solve the Navier Stokes equations. The action is capable of setting up:

  • Incompressible and weakly-compressible simulations for

  • Clean fluids and flows in porous media

The action handles boundary conditions, necessary materials for fluid properties, variables, variable initialization and various objects for turbulence modeling.

commentnote

This action only supports Rhie-Chow interpolation for the determination of face velocities in the advection terms. The face interpolation of the advected quantities (e.g. upwind, average) can be controlled through the *_advection_interpolation action parameters.

Automatically defined variables

The NavierStokesFV action automatically sets up the variables which are necessary for the solution of a given problem. These variables can then be used to couple fluid flow simulations with other physics. The list of variable names commonly used in the action syntax is presented below:

For the default names of other variables used in this action, visit this site.

Bernoulli pressure jump treatment

Please see the Bernoulli pressure variable documentation for more information.

Examples

Incompressible fluid flow in a lid-driven cavity

In the following examples we present how the NavierStokesFV action can be utilized to simplify input files. We start with the simple lid-driven cavity problem with an Incompressible Navier Stokes formulation. The original description of the problem is available under the following link. First, we present the input file where the simulation is set up manually, by defining every kernel, boundary condition and material explicitly.

mu = 1
rho = 1
k = .01
cp = 1
velocity_interp_method = 'rc'
advected_interp_method = 'average'

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

[UserObjects<<<{"href": "../../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": "../../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": "../../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": "../../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": "../../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": "../../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": "../../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": "../../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": "../../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": "../../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": "../../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)
commentnote

Careful! The utilization of central difference (average) advected interpolation may lead to oscillatory behavior in certain scenarios. Even though it is not the case for this example, if this phenomenon arises, we recommend using first order upwind or second order TVD schemes.

The same simulation can be set up using the action syntax which improves input file readability:

mu = 1
rho = 1
k = .01
cp = 1

[Mesh<<<{"href": "../../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
  []
[]

[Modules<<<{"href": "../index.html"}>>>]
  [NavierStokesFV<<<{"href": "index.html"}>>>]
    compressibility<<<{"description": "Compressibility constraint for the Navier-Stokes equations."}>>> = 'incompressible'
    add_energy_equation<<<{"description": "Whether to add the energy equation. This parameter is not necessary if using the Physics syntax"}>>> = true

    density<<<{"description": "The name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'rho'
    dynamic_viscosity<<<{"description": "The name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mu'
    thermal_conductivity<<<{"description": "The name of the fluid thermal conductivity for each block. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
    specific_heat<<<{"description": "The name of the specific heat. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'cp'

    initial_pressure<<<{"description": "The initial pressure, assumed constant everywhere"}>>> = 0.0
    initial_temperature<<<{"description": "The initial temperature, assumed constant everywhere"}>>> = 0.0

    inlet_boundaries<<<{"description": "Names of inlet boundaries"}>>> = 'top'
    momentum_inlet_types<<<{"description": "Types of inlet boundaries for the momentum equation."}>>> = 'fixed-velocity'
    momentum_inlet_function = 'lid_function 0'
    energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'fixed-temperature'
    energy_inlet_function = '0'

    wall_boundaries<<<{"description": "Names of wall boundaries"}>>> = 'left right bottom'
    momentum_wall_types<<<{"description": "Types of wall boundaries for the momentum equation"}>>> = 'noslip noslip noslip'
    energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux fixed-temperature'
    energy_wall_function = '0 0 1'

    pin_pressure<<<{"description": "Switch to enable pressure shifting for incompressible simulations."}>>> = true
    pinned_pressure_type<<<{"description": "Types for shifting (pinning) the pressure in case of incompressible simulations."}>>> = average
    pinned_pressure_value<<<{"description": "The value used for pinning the pressure (point value/domain average)."}>>> = 0

    mass_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating density, as an advected quantity, to the face."}>>> = 'average'
    momentum_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face."}>>> = 'average'
    energy_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face."}>>> = 'average'
  []
[]

[AuxVariables<<<{"href": "../../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": "../../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
  []
[]

[FunctorMaterials<<<{"href": "../../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 rho mu'
    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} ${rho} ${mu}'
  []
[]

[Functions<<<{"href": "../../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": "../../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": "../../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-action.i)

It is visible that in this case we defined the "compressibility" parameter to be incompressible. Furthermore, the energy (enthalpy) equation is solved as well. The user can request this by setting "add_energy_equation" to true. The boundary types are grouped into wall, inlet and outlet types. For more information on the available boundary types, see the list of parameters at the bottom of the page.

Incompressible fluid flow in porous medium

The following input file sets up a simulation of an incompressible fluid flow within a channel which contains a homogenized structure that is treated as a porous medium. The model accounts for the heat exchange between the fluid and the homogenized structure as well. For more description on the used model, visit the Porous medium Incompressible Navier Stokes page. First, the input file with the manually defined kernels and boundary conditions is presented:

mu = 1
rho = 1
k = 1e-3
cp = 1
u_inlet = 1
T_inlet = 200
advected_interp_method = 'average'
velocity_interp_method = 'rc'

[Mesh<<<{"href": "../../Mesh/index.html"}>>>]
  [mesh]
    type = CartesianMeshGenerator<<<{"description": "This CartesianMeshGenerator creates a non-uniform Cartesian mesh.", "href": "../../../source/meshgenerators/CartesianMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    dx<<<{"description": "Intervals in the X direction"}>>> = '5 5'
    dy<<<{"description": "Intervals in the Y direction (required when dim>1 otherwise ignored)"}>>> = '1.0'
    ix<<<{"description": "Number of grids in all intervals in the X direction (default to all one)"}>>> = '50 50'
    iy<<<{"description": "Number of grids in all intervals in the Y direction (default to all one)"}>>> = '20'
    subdomain_id<<<{"description": "Block IDs (default to all zero)"}>>> = '1 2'
  []
[]

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

[UserObjects<<<{"href": "../../UserObjects/index.html"}>>>]
  [rc]
    type = PINSFVRhieChowInterpolator<<<{"description": "Performs interpolations and reconstructions of porosity and computes the Rhie-Chow face velocities.", "href": "../../../source/userobjects/PINSFVRhieChowInterpolator.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = superficial_vel_x
    v<<<{"description": "The y-component of velocity"}>>> = superficial_vel_y
    pressure<<<{"description": "The pressure variable."}>>> = pressure
    porosity<<<{"description": "The porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
  []
[]

[Variables<<<{"href": "../../Variables/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'T_solid'
  [superficial_vel_x]
    type = PINSFVSuperficialVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/PINSFVSuperficialVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${u_inlet}
  []
  [superficial_vel_y]
    type = PINSFVSuperficialVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/PINSFVSuperficialVelocityVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e-6
  []
  [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"}>>>
  []
  [T_solid]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = 'MONOMIAL'
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = 'CONSTANT'
    fv = true
  []
[]

[AuxVariables<<<{"href": "../../AuxVariables/index.html"}>>>]
  [T_solid]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = 'MONOMIAL'
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = 'CONSTANT'
    fv = true
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 100
  []
  [porosity]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    fv = true
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.5
  []
[]

[FVKernels<<<{"href": "../../FVKernels/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'solid_energy_diffusion solid_energy_convection'
  [mass]
    type = PINSFVMassAdvection<<<{"description": "Object for advecting mass in porous media mass equation", "href": "../../../source/fvkernels/PINSFVMassAdvection.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 = PINSFVMomentumAdvection<<<{"description": "Object for advecting superficial momentum, e.g. rho*u_d, in the porous media momentum equation", "href": "../../../source/fvkernels/PINSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = superficial_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}
    porosity<<<{"description": "Porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion<<<{"description": "Viscous diffusion term, div(mu eps grad(u_d / eps)), in the porous media incompressible Navier-Stokes momentum equation.", "href": "../../../source/fvkernels/PINSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = superficial_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}
    porosity<<<{"description": "Porosity auxiliary variable. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_pressure]
    type = PINSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term $eps \nabla P$ into the Navier-Stokes porous media momentum equation.", "href": "../../../source/fvkernels/PINSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = superficial_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
    porosity<<<{"description": "Porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
  []

  [v_advection]
    type = PINSFVMomentumAdvection<<<{"description": "Object for advecting superficial momentum, e.g. rho*u_d, in the porous media momentum equation", "href": "../../../source/fvkernels/PINSFVMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = superficial_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}
    porosity<<<{"description": "Porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion<<<{"description": "Viscous diffusion term, div(mu eps grad(u_d / eps)), in the porous media incompressible Navier-Stokes momentum equation.", "href": "../../../source/fvkernels/PINSFVMomentumDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = superficial_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}
    porosity<<<{"description": "Porosity auxiliary variable. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_pressure]
    type = PINSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term $eps \nabla P$ into the Navier-Stokes porous media momentum equation.", "href": "../../../source/fvkernels/PINSFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = superficial_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
    porosity<<<{"description": "Porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
  []

  [energy_advection]
    type = PINSFVEnergyAdvection<<<{"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/PINSFVEnergyAdvection.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}
  []
  [energy_diffusion]
    type = PINSFVEnergyDiffusion<<<{"description": "Diffusion term in the porous media incompressible Navier-Stokes fluid energy equations :  $-div(eps * k * grad(T))$", "href": "../../../source/fvkernels/PINSFVEnergyDiffusion.html"}>>>
    k<<<{"description": "Thermal conductivity. 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
    porosity<<<{"description": "Porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
  []
  [energy_convection]
    type = PINSFVEnergyAmbientConvection<<<{"description": "Implements the solid-fluid ambient convection term in the porous media Navier Stokes energy equation.", "href": "../../../source/fvkernels/PINSFVEnergyAmbientConvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
    is_solid<<<{"description": "Whether this kernel acts on the solid temperature"}>>> = false
    T_fluid<<<{"description": "Fluid temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_fluid'
    T_solid<<<{"description": "Solid temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_solid'
    h_solid_fluid<<<{"description": "Name of the convective heat transfer coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'h_cv'
  []

  [solid_energy_diffusion]
    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_solid
  []
  [solid_energy_convection]
    type = PINSFVEnergyAmbientConvection<<<{"description": "Implements the solid-fluid ambient convection term in the porous media Navier Stokes energy equation.", "href": "../../../source/fvkernels/PINSFVEnergyAmbientConvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_solid
    is_solid<<<{"description": "Whether this kernel acts on the solid temperature"}>>> = true
    T_fluid<<<{"description": "Fluid temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_fluid'
    T_solid<<<{"description": "Solid temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_solid'
    h_solid_fluid<<<{"description": "Name of the convective heat transfer coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'h_cv'
  []
[]

[FVBCs<<<{"href": "../../FVBCs/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'heated-side'
  [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"}>>> = superficial_vel_x
    function = ${u_inlet}
  []
  [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"}>>> = superficial_vel_y
    function = 0
  []
  [inlet-T]
    type = FVNeumannBC<<<{"description": "Neumann boundary condition for finite volume method.", "href": "../../../source/fvbcs/FVNeumannBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = T_fluid
    value<<<{"description": "The value of the flux crossing the boundary."}>>> = '${fparse u_inlet * rho * cp * T_inlet}'
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
  []

  [no-slip-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'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = superficial_vel_x
    function<<<{"description": "The exact solution function."}>>> = 0
  []
  [no-slip-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'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = superficial_vel_y
    function<<<{"description": "The exact solution function."}>>> = 0
  []
  [heated-side]
    type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../../source/fvbcs/FVDirichletBC.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"}>>> = 'T_solid'
    value<<<{"description": "value to enforce at the boundary face"}>>> = 150
  []

  [symmetry-u]
    type = PINSFVSymmetryVelocityBC<<<{"description": "Implements a symmetry boundary condition for the velocity.", "href": "../../../source/fvbcs/PINSFVSymmetryVelocityBC.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"}>>> = superficial_vel_x
    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."}>>> = superficial_vel_x
    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."}>>> = superficial_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."}>>> = 'x'
  []
  [symmetry-v]
    type = PINSFVSymmetryVelocityBC<<<{"description": "Implements a symmetry boundary condition for the velocity.", "href": "../../../source/fvbcs/PINSFVSymmetryVelocityBC.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"}>>> = superficial_vel_y
    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."}>>> = superficial_vel_x
    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."}>>> = superficial_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'
  []
  [symmetry-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"}>>> = 'bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
  []

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

[FunctorMaterials<<<{"href": "../../FunctorMaterials/index.html"}>>>]
  [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"}>>> = 'h_cv'
    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."}>>> = '1'
  []
  [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'
    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}'
  []
  [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"}>>>
    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}
    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'
  []
[]

[Executioner<<<{"href": "../../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-14
[]

# Some basic Postprocessors to examine the solution
[Postprocessors<<<{"href": "../../Postprocessors/index.html"}>>>]
  [inlet-p]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../../../source/postprocessors/SideAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = pressure
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
  []
  [outlet-u]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../../../source/postprocessors/SideAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = superficial_vel_x
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
  []
  [outlet-temp]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../../../source/postprocessors/SideAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = T_fluid
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
  []
  [solid-temp]
    type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../../../source/postprocessors/ElementAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable that this object operates on"}>>> = T_solid
  []
[]

[Outputs<<<{"href": "../../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."}>>> = false
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/heated/2d-rc-heated.i)
commentnote

Careful! The utilization of central difference (average) advected interpolation may lead to oscillatory behavior in certain scenarios. Even though it is not the case for this example, if this phenomenon arises, we recommend using first order upwind or second order TVD schemes.

The same simulation can also be set up using the NavierStokesFV action syntax:

mu = 1
rho = 1
k = 1e-3
cp = 1
u_inlet = 1
T_inlet = 200
h_cv = 1.0

[Mesh<<<{"href": "../../Mesh/index.html"}>>>]
  [mesh]
    type = CartesianMeshGenerator<<<{"description": "This CartesianMeshGenerator creates a non-uniform Cartesian mesh.", "href": "../../../source/meshgenerators/CartesianMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    dx<<<{"description": "Intervals in the X direction"}>>> = '5 5'
    dy<<<{"description": "Intervals in the Y direction (required when dim>1 otherwise ignored)"}>>> = '1.0'
    ix<<<{"description": "Number of grids in all intervals in the X direction (default to all one)"}>>> = '50 50'
    iy<<<{"description": "Number of grids in all intervals in the Y direction (default to all one)"}>>> = '20'
    subdomain_id<<<{"description": "Block IDs (default to all zero)"}>>> = '1 2'
  []
[]

[Variables<<<{"href": "../../Variables/index.html"}>>>]
  [T_solid]
    type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/MooseVariableFV.html"}>>>
  []
[]

[AuxVariables<<<{"href": "../../AuxVariables/index.html"}>>>]
  [porosity]
    type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/MooseVariableFV.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.5
  []
[]

[Modules<<<{"href": "../index.html"}>>>]
  [NavierStokesFV<<<{"href": "index.html"}>>>]
    compressibility<<<{"description": "Compressibility constraint for the Navier-Stokes equations."}>>> = 'incompressible'
    porous_medium_treatment<<<{"description": "Whether to use porous medium kernels or not."}>>> = true
    add_energy_equation<<<{"description": "Whether to add the energy equation. This parameter is not necessary if using the Physics syntax"}>>> = true

    density<<<{"description": "The name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    dynamic_viscosity<<<{"description": "The name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    thermal_conductivity<<<{"description": "The name of the fluid thermal conductivity for each block. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${k}
    specific_heat<<<{"description": "The name of the specific heat. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${cp}
    porosity<<<{"description": "The name of the auxiliary variable for the porosity field. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'porosity'

    # Reference file sets effective_conductivity by default that way
    # so the conductivity is multiplied by the porosity in the kernel
    effective_conductivity<<<{"description": "Whether the conductivity should be multiplied by porosity, or whether the provided conductivity is an effective conductivity taking porosity effects into account"}>>> = false

    initial_velocity<<<{"description": "The initial velocity, assumed constant everywhere"}>>> = '${u_inlet} 1e-6 0'
    initial_pressure<<<{"description": "The initial pressure, assumed constant everywhere"}>>> = 0.0
    initial_temperature<<<{"description": "The initial temperature, assumed constant everywhere"}>>> = 0.0

    inlet_boundaries<<<{"description": "Names of inlet boundaries"}>>> = 'left'
    momentum_inlet_types<<<{"description": "Types of inlet boundaries for the momentum equation."}>>> = 'fixed-velocity'
    momentum_inlet_function = '${u_inlet} 0'
    energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'heatflux'
    energy_inlet_function = '${fparse u_inlet * rho * cp * T_inlet}'

    wall_boundaries<<<{"description": "Names of wall boundaries"}>>> = 'top bottom'
    momentum_wall_types<<<{"description": "Types of wall boundaries for the momentum equation"}>>> = 'noslip symmetry'
    energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux'
    energy_wall_function = '0 0'

    outlet_boundaries<<<{"description": "Names of outlet boundaries"}>>> = 'right'
    momentum_outlet_types<<<{"description": "Types of outlet boundaries for the momentum equation"}>>> = 'fixed-pressure'
    pressure_function = '0.1'

    ambient_convection_alpha<<<{"description": "The heat exchange coefficients for each block in 'ambient_convection_blocks'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${h_cv}
    ambient_temperature<<<{"description": "The ambient temperature for each block in 'ambient_convection_blocks'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_solid'

    mass_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating density, as an advected quantity, to the face."}>>> = 'average'
    momentum_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face."}>>> = 'average'
    energy_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face."}>>> = 'average'
  []
[]

[FVKernels<<<{"href": "../../FVKernels/index.html"}>>>]
  [solid_energy_diffusion]
    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_solid
  []
  [solid_energy_convection]
    type = PINSFVEnergyAmbientConvection<<<{"description": "Implements the solid-fluid ambient convection term in the porous media Navier Stokes energy equation.", "href": "../../../source/fvkernels/PINSFVEnergyAmbientConvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'T_solid'
    is_solid<<<{"description": "Whether this kernel acts on the solid temperature"}>>> = true
    T_fluid<<<{"description": "Fluid temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_fluid'
    T_solid<<<{"description": "Solid temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_solid'
    h_solid_fluid<<<{"description": "Name of the convective heat transfer coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${h_cv}
  []
[]

[FVBCs<<<{"href": "../../FVBCs/index.html"}>>>]
  [heated-side]
    type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../../source/fvbcs/FVDirichletBC.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"}>>> = 'T_solid'
    value<<<{"description": "value to enforce at the boundary face"}>>> = 150
  []
[]

[Executioner<<<{"href": "../../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-14
[]

# Some basic Postprocessors to examine the solution
[Postprocessors<<<{"href": "../../Postprocessors/index.html"}>>>]
  [inlet-p]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../../../source/postprocessors/SideAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = pressure
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
  []
  [outlet-u]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../../../source/postprocessors/SideAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = superficial_vel_x
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
  []
  [outlet-temp]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../../../source/postprocessors/SideAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = T_fluid
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
  []
  [solid-temp]
    type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../../../source/postprocessors/ElementAverageValue.html"}>>>
    variable<<<{"description": "The name of the variable that this object operates on"}>>> = T_solid
  []
[]

[Outputs<<<{"href": "../../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."}>>> = false
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/heated/2d-rc-heated-action.i)

Compared to the previous example, we see that in this case the porous medium treatment is enabled by setting "porous_medium_treatment" to true. The corresponding porosity can be supplied through the "porosity" parameter. Furthermore, the heat exchange between the fluid and the homogenized structure is enabled using the "ambient_temperature" and "ambient_convection_alpha" parameters.

Weakly-compressible fluid flow

The last example is dedicated to demonstrating a transient flow in a channel using a weakly-compressible approximation. The following examples shows how this simulation is set up by manually defining the kernels and boundary conditions. For more information on the weakly-compressible treatment, visit the Weakly-compressible Navier Stokes page.

rho = 'rho'
l = 10
velocity_interp_method = 'rc'
advected_interp_method = 'average'

# Artificial fluid properties
# For a real case, use a GeneralFluidFunctorProperties and a viscosity rampdown
# or initialize very well!
k = 1
cp = 1000
mu = 1e2

# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_v = 0.001

[Mesh<<<{"href": "../../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"}>>> = ${l}
    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"}>>> = 20
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 10
  []
[]

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

[UserObjects<<<{"href": "../../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": "../../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"}>>> = ${inlet_v}
  []
  [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"}>>> = 1e-15
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/INSFVPressureVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${outlet_pressure}
  []
  [T_fluid]
    type = INSFVEnergyVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/INSFVEnergyVariable.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${inlet_temp}
  []
[]

[AuxVariables<<<{"href": "../../AuxVariables/index.html"}>>>]
  [mixing_length]
    type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/MooseVariableFV.html"}>>>
  []
  [power_density]
    type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/MooseVariableFV.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e4
  []
[]

[FVKernels<<<{"href": "../../FVKernels/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'u_turb v_turb temp_turb'
  [mass_time]
    type = WCNSFVMassTimeDerivative<<<{"description": "Adds the time derivative term to the weakly-compressible Navier-Stokes continuity equation.", "href": "../../../source/fvkernels/WCNSFVMassTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
    drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
  []
  [mass]
    type = WCNSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../../source/fvkernels/WCNSFVMassAdvection.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_time]
    type = WCNSFVMomentumTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../../source/fvkernels/WCNSFVMomentumTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
    rho<<<{"description": "The density material property. 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_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
  []
  [u_turb]
    type = INSFVMixingLengthReynoldsStress<<<{"description": "Computes the force due to the Reynolds stress term in the incompressible Reynolds-averaged Navier-Stokes equations.", "href": "../../../source/fvkernels/INSFVMixingLengthReynoldsStress.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
    rho<<<{"description": "fluid density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    mixing_length<<<{"description": "Turbulent eddy mixing length. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mixing_length'
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    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."}>>> = vel_x
    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."}>>> = vel_y
  []

  [v_time]
    type = WCNSFVMomentumTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../../source/fvkernels/WCNSFVMomentumTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
    rho<<<{"description": "The density material property. 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_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
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = '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}
  []
  [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_turb]
    type = INSFVMixingLengthReynoldsStress<<<{"description": "Computes the force due to the Reynolds stress term in the incompressible Reynolds-averaged Navier-Stokes equations.", "href": "../../../source/fvkernels/INSFVMixingLengthReynoldsStress.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
    rho<<<{"description": "fluid density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    mixing_length<<<{"description": "Turbulent eddy mixing length. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mixing_length'
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    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."}>>> = vel_x
    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."}>>> = vel_y
  []

  [temp_time]
    type = WCNSFVEnergyTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../../source/fvkernels/WCNSFVEnergyTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
    rho<<<{"description": "Density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
    drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
    h<<<{"description": "The specific enthalpy. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = h
    dh_dt<<<{"description": "The time derivative of the specific enthalpy. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = dh_dt
  []
  [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}
  []
  [heat_source]
    type = FVCoupledForce<<<{"description": "Implements a source term proportional to the value of a coupled variable.", "href": "../../../source/fvkernels/FVCoupledForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
    v<<<{"description": "The coupled functor which provides the force. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = power_density
  []
  [temp_turb]
    type = WCNSFVMixingLengthEnergyDiffusion<<<{"description": "Computes the turbulent diffusive flux that appears in Reynolds-averaged fluid energy conservation equations.", "href": "../../../source/fvkernels/WCNSFVMixingLengthEnergyDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
    rho<<<{"description": "Density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
    cp<<<{"description": "Specific heat capacity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = cp
    mixing_length<<<{"description": "The turbulent mixing length. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mixing_length'
    schmidt_number<<<{"description": "The turbulent Schmidt number (or turbulent Prandtl number if the passive scalar is energy) that relates the turbulent scalar diffusivity to the turbulent momentum diffusivity."}>>> = 1
    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."}>>> = vel_x
    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."}>>> = vel_y
  []
[]

[FVBCs<<<{"href": "../../FVBCs/index.html"}>>>]
  [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"}>>> = 'top 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"}>>> = 'top bottom'
    function<<<{"description": "The exact solution function."}>>> = 0
  []

  # Inlet
  [inlet_u]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../../source/fvbcs/INSFVInletVelocityBC.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'
    function = ${inlet_v}
  []
  [inlet_v]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../../source/fvbcs/INSFVInletVelocityBC.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'
    function = 0
  []
  [inlet_T]
    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"}>>> = 'left'
    value<<<{"description": "value to enforce at the boundary face"}>>> = ${inlet_temp}
  []

  [outlet_p]
    type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
    function<<<{"description": "The boundary pressure as a regular function"}>>> = ${outlet_pressure}
  []
[]

[FluidProperties<<<{"href": "../../FluidProperties/index.html"}>>>]
  [fp]
    type = FlibeFluidProperties<<<{"description": "Fluid properties for flibe", "href": "../../../source/fluidproperties/FlibeFluidProperties.html"}>>>
  []
[]

[FunctorMaterials<<<{"href": "../../FunctorMaterials/index.html"}>>>]
  [const_functor]
    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}'
  []
  [rho]
    type = RhoFromPTFunctorMaterial<<<{"description": "Computes the density from coupled pressure and temperature functors (variables, functions, functor material properties", "href": "../../../source/functormaterials/RhoFromPTFunctorMaterial.html"}>>>
    fp<<<{"description": "fluid userobject"}>>> = fp
    temperature<<<{"description": "temperature functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = T_fluid
    pressure<<<{"description": "pressure functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
  [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}
  []
[]

[AuxKernels<<<{"href": "../../AuxKernels/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'mixing_len'
  [mixing_len]
    type = WallDistanceMixingLengthAux<<<{"description": "Computes the turbulent mixing length by assuming that it is proportional to the distance from the nearest wall. The mixinglength is capped at a distance proportional to inputted parameter delta.", "href": "../../../source/auxkernels/WallDistanceMixingLengthAux.html"}>>>
    walls<<<{"description": "Boundaries that correspond to solid walls."}>>> = 'top'
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = mixing_length
    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."}>>> = 'initial'
    delta<<<{"description": ". A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 0.5
  []
[]

[Executioner<<<{"href": "../../Executioner/index.html"}>>>]
  type = Transient
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu       NONZERO'

  [TimeStepper<<<{"href": "../../Executioner/TimeStepper/index.html"}>>>]
    type = IterationAdaptiveDT
    dt = 1e-3
    optimal_iterations = 6
  []
  end_time = 15

  nl_abs_tol = 1e-9
  nl_max_its = 50
  line_search = 'none'

  automatic_scaling = true
  off_diagonals_in_auto_scaling = true
  compute_scaling_once = false
[]

[Outputs<<<{"href": "../../Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/channel-flow/2d-transient.i)
commentnote

Careful! The utilization of central difference (average) advected interpolation may lead to oscillatory behavior in certain scenarios. Even though it is not the case for this example, if this phenomenon arises, we recommend using first order upwind or second order TVD schemes.

The same simulation can be set up using the action syntax as follows:

l = 10

# Artificial fluid properties
# For a real case, use a GeneralFluidFunctorProperties and a viscosity rampdown
# or initialize very well!
k = 1
cp = 1000
mu = 1e2

# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_v = 0.001

[Mesh<<<{"href": "../../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"}>>> = ${l}
    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"}>>> = 20
    ny<<<{"description": "Number of elements in the Y direction"}>>> = 10
  []
[]

[Modules<<<{"href": "../index.html"}>>>]
  [NavierStokesFV<<<{"href": "index.html"}>>>]
    compressibility<<<{"description": "Compressibility constraint for the Navier-Stokes equations."}>>> = 'weakly-compressible'
    add_energy_equation<<<{"description": "Whether to add the energy equation. This parameter is not necessary if using the Physics syntax"}>>> = true

    density<<<{"description": "The name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'rho'
    dynamic_viscosity<<<{"description": "The name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mu'
    thermal_conductivity<<<{"description": "The name of the fluid thermal conductivity for each block. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
    specific_heat<<<{"description": "The name of the specific heat. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'cp'

    initial_velocity<<<{"description": "The initial velocity, assumed constant everywhere"}>>> = '${inlet_v} 1e-15 0'
    initial_temperature<<<{"description": "The initial temperature, assumed constant everywhere"}>>> = '${inlet_temp}'
    initial_pressure<<<{"description": "The initial pressure, assumed constant everywhere"}>>> = '${outlet_pressure}'

    inlet_boundaries<<<{"description": "Names of inlet boundaries"}>>> = 'left'
    momentum_inlet_types<<<{"description": "Types of inlet boundaries for the momentum equation."}>>> = 'fixed-velocity'
    momentum_inlet_functors<<<{"description": "Functions for inlet boundary velocities or pressures (for fixed-pressure option). Provide a double vector where the leading dimension corresponds to the number of fixed-velocity and fixed-pressure entries in momentum_inlet_types and the second index runs either over dimensions for fixed-velocity boundaries or is a single function name for pressure inlets. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${inlet_v} 0'
    energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'fixed-temperature'
    energy_inlet_functors<<<{"description": "Functions for fixed-value boundaries in the energy equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${inlet_temp}'

    wall_boundaries<<<{"description": "Names of wall boundaries"}>>> = 'top bottom'
    momentum_wall_types<<<{"description": "Types of wall boundaries for the momentum equation"}>>> = 'noslip noslip'
    energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux'
    energy_wall_functors<<<{"description": "Functions for Dirichlet/Neumann boundaries in the energy equation. For wall types requiring multiple functions, the syntax is <function_1>:<function_2>:... So, 'convection' types are '<Tinf_function>:<htc_function>'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '0 0'

    outlet_boundaries<<<{"description": "Names of outlet boundaries"}>>> = 'right'
    momentum_outlet_types<<<{"description": "Types of outlet boundaries for the momentum equation"}>>> = 'fixed-pressure'
    pressure_functors<<<{"description": "Functions for boundary pressures at outlets. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${outlet_pressure}'

    external_heat_source<<<{"description": "The name of a functor which contains the external heat source for the energy equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'power_density'

    mass_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating density, as an advected quantity, to the face."}>>> = 'average'
    momentum_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face."}>>> = 'average'
    energy_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face."}>>> = 'average'
  []
[]

[AuxVariables<<<{"href": "../../AuxVariables/index.html"}>>>]
  [power_density]
    type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../../source/variables/MooseVariableFV.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e4
  []
[]

[FluidProperties<<<{"href": "../../FluidProperties/index.html"}>>>]
  [fp]
    type = FlibeFluidProperties<<<{"description": "Fluid properties for flibe", "href": "../../../source/fluidproperties/FlibeFluidProperties.html"}>>>
  []
[]

[FunctorMaterials<<<{"href": "../../FunctorMaterials/index.html"}>>>]
  [const_functor]
    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 mu'
    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} ${mu}'
  []
  [rho]
    type = RhoFromPTFunctorMaterial<<<{"description": "Computes the density from coupled pressure and temperature functors (variables, functions, functor material properties", "href": "../../../source/functormaterials/RhoFromPTFunctorMaterial.html"}>>>
    fp<<<{"description": "fluid userobject"}>>> = fp
    temperature<<<{"description": "temperature functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = T_fluid
    pressure<<<{"description": "pressure functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
  []
[]

[Executioner<<<{"href": "../../Executioner/index.html"}>>>]
  type = Transient
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu       NONZERO'

  [TimeStepper<<<{"href": "../../Executioner/TimeStepper/index.html"}>>>]
    type = IterationAdaptiveDT
    dt = 1e-3
    optimal_iterations = 6
  []
  end_time = 15

  nl_abs_tol = 1e-9
  nl_max_its = 50
  line_search = 'none'

  automatic_scaling = true
  off_diagonals_in_auto_scaling = true
  compute_scaling_once = false
[]

[Outputs<<<{"href": "../../Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/channel-flow/2d-transient-action.i)

We note that the weakly-compressible handling can be enabled by setting "compressibility" to weakly-compressible. As shown in the example, an arbitrary energy source function can also be supplied to the incorporated energy equation using the "external_heat_source" parameter.

How to transition to the Physics syntax

The /Modules/NavierStokes action syntax and (Physics/NavierStokes/Flow/.., Physics/NavierStokes/FluidHeatTransfer/.., Physics/NavierStokes/ScalarTransport/.., Physics/NavierStokes/Turbulence/..) syntax create the exact same objects in the background. We currently do not expect any difference in results, notably because /Modules/NavierStokes has been changed to create the relevant Physics under the hood!

To transition, you will have to split the /Modules/NavierStokes parameters into four groups below. Your simulation may only feature only one of these groups, in which case you will only need to create a single Physics:

  • mass and momentum (conservation) equations

  • heat transfer / energy (conservation) equation

  • scalar conservation equations

  • turbulence parameters

and copy paste the mass/momentum parameters into the [Physics/NavierStokes/Flow/<name>] syntax as shown in this example:

[Physics]
  [NavierStokes]
    [Flow]
      [flow]
        compressibility = 'incompressible'

        density = 'rho'
        dynamic_viscosity = 'mu'

        initial_velocity = '${u_inlet} 1e-12 0'
        initial_pressure = 0.0

        inlet_boundaries = 'left'
        momentum_inlet_types = 'fixed-velocity'
        momentum_inlet_function = '${u_inlet} 0'
        wall_boundaries = 'bottom top'
        momentum_wall_types = 'symmetry noslip'

        outlet_boundaries = 'right'
        momentum_outlet_types = 'fixed-pressure-zero-gradient'
        pressure_function = '${p_outlet}'

        mass_advection_interpolation = 'average'
        momentum_advection_interpolation = 'average'
      []
    []
  []
[]
# This separation is introduced for documentation purposes.
# Both Physics could be nested under Physics/NavierStokes
[Physics/NavierStokes]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-transient-physics.i)

the energy conservation parameters into this [Physics/NavierStokes/FluidHeatTransfer/<name>] syntax:

# Both Physics could be nested under Physics/NavierStokes
[Physics<<<{"href": "../../Physics/index.html"}>>>/NavierStokes<<<{"href": "../../Physics/NavierStokes/index.html"}>>>]
  [FluidHeatTransfer<<<{"href": "../../Physics/NavierStokes/FluidHeatTransfer/index.html"}>>>]
    [heat]
      thermal_conductivity<<<{"description": "The name of the fluid thermal conductivity for each block. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
      specific_heat<<<{"description": "The name of the specific heat. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'cp'

      fluid_temperature_variable<<<{"description": "Name of the fluid temperature variable"}>>> = 'T_fluid'
      initial_temperature<<<{"description": "The initial temperature, assumed constant everywhere"}>>> = '${T_inlet}'
      energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'heatflux'
      energy_inlet_functors<<<{"description": "Functions for fixed-value boundaries in the energy equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${fparse u_inlet * rho * cp * T_inlet}'

      energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux'
      energy_wall_functors<<<{"description": "Functions for Dirichlet/Neumann boundaries in the energy equation. For wall types requiring multiple functions, the syntax is <function_1>:<function_2>:... So, 'convection' types are '<Tinf_function>:<htc_function>'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '0 0'

      ambient_convection_alpha<<<{"description": "The heat exchange coefficients for each block in 'ambient_convection_blocks'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'h_cv'
      ambient_temperature<<<{"description": "The ambient temperature for each block in 'ambient_convection_blocks'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_solid'

      energy_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face."}>>> = 'average'
    []
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-transient-physics.i)

the scalar conservation equations parameters into this [Physics/NavierStokes/ScalarTransport/<name>] syntax:

# Both could be nested under Physics/NavierStokes
[Physics<<<{"href": "../../Physics/index.html"}>>>/NavierStokes<<<{"href": "../../Physics/NavierStokes/index.html"}>>>]
  [ScalarTransport<<<{"href": "../../Physics/NavierStokes/ScalarTransport/index.html"}>>>]
    [heat]
      passive_scalar_names<<<{"description": "Vector containing the names of the advected scalar variables."}>>> = 'scalar'

      passive_scalar_diffusivity<<<{"description": "Functor names for the diffusivities used for the passive scalar fields. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${diff}
      passive_scalar_source<<<{"description": "Passive scalar sources. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 0.1
      passive_scalar_coupled_source<<<{"description": "Coupled variable names for the sources used for the passive scalar fields. If multiple sources for each equation are specified, major (outer) ordering by equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = U
      passive_scalar_coupled_source_coeff<<<{"description": "Coupled variable multipliers for the sources used for the passive scalar fields. If multiple sources for each equation are specified, major (outer) ordering by equation."}>>> = 0.1

      passive_scalar_inlet_types<<<{"description": "Types for the inlet boundaries for the passive scalar equation."}>>> = 'fixed-value'
      passive_scalar_inlet_function = '1'

      passive_scalar_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating passive scalar field, as an advected quantity, to the face."}>>> = 'average'
    []
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-scalar-transport-physics.i)

and finally the turbulence parameters into this [Physics/NavierStokes/Turbulence/<name>] syntax:

# Both could be nested under Physics/NavierStokes
[Physics<<<{"href": "../../Physics/index.html"}>>>/NavierStokes<<<{"href": "../../Physics/NavierStokes/index.html"}>>>]
  [Turbulence<<<{"href": "../../Physics/NavierStokes/Turbulence/index.html"}>>>]
    [mixing-length]
      turbulence_handling<<<{"description": "The way turbulent diffusivities are determined in the turbulent regime."}>>> = 'mixing-length'
      coupled_flow_physics<<<{"description": "WCNSFVFlowPhysics generating the velocities"}>>> = flow
      scalar_transport_physics<<<{"description": "WCNSFVScalarTransportPhysics generating the scalar advection equations"}>>> = scalars

      passive_scalar_schmidt_number = 1.0

      von_karman_const<<<{"description": "Von Karman parameter for the mixing length model. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${von_karman_const}
      mixing_length_delta<<<{"description": "Tunable parameter related to the thickness of the boundary layer.When it is not specified, Prandtl's original unbounded wall distance mixing length model isretrieved. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 1e9
      mixing_length_walls = 'top bottom'
      mixing_length_aux_execute_on<<<{"description": "When the mixing length aux kernels should be executed."}>>> = 'initial'
    []
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-mixing-length-physics.i)
commentnote

All NavierStokes Physics may be nested under the same [NavierStokes] sub-block. The separation in the examples into different Physics blocks is only for the purpose of simplifying the examples.

If all goes well, the input will give exactly the same result. Some parameters have been renamed but the old names are currently still supported, even in the new syntax. If it does not go well, you can use the DumpObjectsProblem with the "dump_path" parameter set to:

  • Modules/NavierStokesFV for the initial input

  • Physics/NavierStokes/... for the new input

to compare all the objects and parameters of both syntaxes. They should match rigorously. Any difference should be reported on the MOOSE discussions forum. Please attach both inputs and a description of the differences found to facilitate our analysis.

Available Actions