Finite Volume Incompressible Porous media Navier Stokes

commentnote

The weakly compressible formulation is also implemented for porous flow finite volume.

Equations

This module implements the porous media Navier Stokes equations. They are expressed in terms of the superficial velocity where is the porosity and the interstitial velocity. The superficial velocity is also known as the extrinsic or Darcy velocity. The other non-linear variables used are pressure and temperature. This is known as the primitive superficial set of variables.

Mass equation:

Momentum equation, with friction and gravity force as example forces: (1)

Fluid phase energy equation, with a convective heat transfer term:

Solid phase energy equation, with convective heat transfer and an energy source :

where is the fluid density, the dynamic viscosity, the specific heat capacity the convective heat transfer coefficient. The effective diffusivities include the effect of porosity, which is often approximated as with the thermal diffusivity.

Implementation for a straight channel

We define a straight channel using a simple Cartesian mesh.

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [gen]
    type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
    dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
    xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
    xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = 5
    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/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

The non linear variables are specified among the INSFV and PINSFV sets of variables as their evaluation and the computation of their gradients on faces is specific to the treatment of the Navier Stokes equations.

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'lambda'
  [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"}>>> = 1
  []
  [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"}>>>
  []
  [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
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

The porous media Navier Stokes equations are implemented by using FVKernels which correspond to each term of the equations. In the following example, the first kernel is the mass advection kernel. This kernel corresponds to the conservation of mass. It acts on the pressure non-linear variable which appears in the mass equation through the Rhie Chow interpolation of the mass fluxes on the element faces, and through its action on the velocity in the momentum equation.

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [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
    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}
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

Then we add the kernels acting on the 'x' component of the momentum, so essentially transcribing the first component of the momentum equation into MOOSE.

The momentum advection term is defined by a PINSFVMomentumAdvection kernel. The momentum equation requires a velocity and an advected quantity interpolation method. This is because this kernel is a FVFlux kernel, it uses the divergence theorem to compute the divergence based on face fluxes rather than on the volumetric divergence value. The velocity interpolation method may be average or rc (Rhie Chow). The latter is preferred to avoid oscillations in pressure due to the grid collocation.

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [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
    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'
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

The next term is the momentum diffusion. This term arises from viscous stress in the fluid. The incompressible approximation simplifies its treatment. This term is also evaluated using the divergence theorem.

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [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'
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

The pressure gradient term in the momentum equation is expressed using a PINSFVMomentumPressure term.

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [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
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

The same set of kernels are defined for the 'y' component of the superficial velocity to transcribe the 'y' component of the momentum equation.

This test features a variety of boundary conditions listed in the FVBCs block. They may not all be used at the same time, the user has to selectively activate the desired one. We list the following example boundary conditions, for the first component of the superficial velocity:

  • no slip walls

    [FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
      [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 bottom'
        variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = superficial_vel_x
        function<<<{"description": "The exact solution function."}>>> = 0
      []
    []
    (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

  • free slip walls

    [FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
      [free-slip-u]
        type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
        boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
        variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = superficial_vel_x
        momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
      []
    []
    (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

  • symmetry axis. This symmetry condition should also be indicated for the pressure variable.

    [FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
      [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'
      []
    []
    (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)
    [FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
      [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
      []
    []
    (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

  • inlet velocity, to specify mass flux given that density is constant

    [FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
      [inlet-u]
        type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
        boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
        variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = superficial_vel_x
        function = '1'
      []
    []
    (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

  • momentum advection outflow (only for a mean-pressure approach, equivalent to executing the momentum advection kernel on the boundary)

    [FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
      [outlet-u]
        type = PINSFVMomentumAdvectionOutflowBC<<<{"description": "Outflow boundary condition for advecting momentum in the porous media momentum equation. This will impose a zero normal gradient on the boundary velocity.", "href": "../../source/fvbcs/PINSFVMomentumAdvectionOutflowBC.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"}>>> = superficial_vel_x
        u<<<{"description": "The superficial 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 superficial 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
        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
        momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
        rho<<<{"description": "The density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
      []
    []
    (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)
commentnote

If the PINSFV version of a boundary condition does not exist, it may be because the INSFV version is valid to use, replacing velocity by superficial velocity.

The pressure boundary condition is usually only set at the outlet:

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [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
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

For a mean-pressure approach, usually for cavity problems, the user may specify a mass advection boundary condition. This is equivalent to executing the mass advection kernel on boundaries.

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [outlet-p-novalue]
    type = INSFVMassAdvectionOutflowBC<<<{"description": "Outflow boundary condition for advecting mass.", "href": "../../source/fvbcs/INSFVMassAdvectionOutflowBC.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
    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
    rho<<<{"description": "The functor for the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

Example inputs : heated straight channel

This channel features a porous media with temperature exchanged by convection between the fluid and the solid phase. The fluid temperature obeys the porous media Navier Stokes energy equation with both advection and diffusion of the fluid energy. The inlet boundary condition is set to a FVNeumannBC to specify an incoming flux. An alternative would be to use a FVDirichletBC to specify an inlet temperature, but depending on the numerical scheme for the computation of advected quantities on faces, this is not as exact.

The solid temperature is by default constant in this input file, however by activating the relevant non-linear Variable, FVKernels and FVBCs, this input also allows solving for the solid phase temperature.

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

[UserObjects<<<{"href": "../../syntax/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": "../../syntax/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": "../../syntax/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": "../../syntax/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": "../../syntax/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": "../../syntax/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": "../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  nl_rel_tol = 1e-14
[]

# Some basic Postprocessors to examine the solution
[Postprocessors<<<{"href": "../../syntax/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": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = false
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/heated/2d-rc-heated.i)

Example inputs : straight channel with a porosity change

This channel features a change in porosity between 1 and 0.5 with either a sharp discontinuity or a continuous variation between the two regions, depending on the porosity initial condition (ICs block).

Both Rhie Chow interpolation and the momentum equation feature a pressure gradient, and the superficial momentum advection and diffusion terms feature a porosity or inverse porosity gradient. These are ill-defined near the discontinuity. To deal with this issue, the discontinuous case uses a flux-based formulation of the pressure gradient, PINSFVMomentumPressureFlux.

To study cases with continuous porosity gradients, the smooth_porosity boolean parameter may be provided to the kernels, which will enable the full treatment of pressure and porosity gradients.

mu=1.1
rho=1.1
advected_interp_method='upwind'
velocity_interp_method='rc'

[Mesh<<<{"href": "../../syntax/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"}>>> = '1 1'
    dy<<<{"description": "Intervals in the Y direction (required when dim>1 otherwise ignored)"}>>> = '0.5'
    ix<<<{"description": "Number of grids in all intervals in the X direction (default to all one)"}>>> = '30 30'
    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": "../../syntax/GlobalParams/index.html"}>>>]
  rhie_chow_user_object = 'rc'
  porosity = porosity
[]

[UserObjects<<<{"href": "../../syntax/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"}>>> = u
    v<<<{"description": "The y-component of velocity"}>>> = v
    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
    pressure<<<{"description": "The pressure variable."}>>> = pressure
  []
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [u]
    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"}>>> = 1
  []
  [v]
    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"}>>>
  []
[]

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

[ICs<<<{"href": "../../syntax/ICs/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'porosity_continuous'
  [porosity_1]
    type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../source/ics/ConstantIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = porosity
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
    value<<<{"description": "The value to be set in IC"}>>> = 1
  []
  [porosity_2]
    type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../source/ics/ConstantIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = porosity
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
    value<<<{"description": "The value to be set in IC"}>>> = 0.5
  []
  [porosity_continuous]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = porosity
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1 2'
    function<<<{"description": "The initial condition function."}>>> = smooth_jump
  []
[]

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

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [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"}>>> = u
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_viscosity]
    type = 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"}>>> = u
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [u_pressure]
    type = 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"}>>> = u
    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
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []

  [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"}>>> = v
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_viscosity]
    type = 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"}>>> = v
    mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [v_pressure]
    type = 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"}>>> = v
    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
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
[]

[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
  [inlet-u]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    function = 0
  []

  [walls-u]
    type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [walls-v]
    type = INSFVNaturalFreeSlipBC<<<{"description": "Implements a free slip boundary condition naturally.", "href": "../../source/fvbcs/INSFVNaturalFreeSlipBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
    variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []

  [outlet_p]
    type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '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.4
  []
[]

[FunctorMaterials<<<{"href": "../../syntax/FunctorMaterials/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'smooth'
  [jump]
    type = ADPiecewiseByBlockFunctorMaterial<<<{"description": "Computes a property value on a per-subdomain basis", "href": "../../source/functormaterials/PiecewiseByBlockFunctorMaterial.html"}>>>
    prop_name<<<{"description": "The name of the property to declare. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'porosity'
    subdomain_to_prop_value<<<{"description": "Map from subdomain to property value. The value may be a constant or any kind of functor (functions, variables, functor material properties)"}>>> = '1 1
                               2 0.5'
  []
  [smooth]
    type = ADGenericFunctionFunctorMaterial<<<{"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"}>>> = 'porosity'
    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."}>>> = 'smooth_jump'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu       NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-10
[]

[Postprocessors<<<{"href": "../../syntax/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 = SideIntegralVariablePostprocessor<<<{"description": "Computes a surface integral of the specified variable", "href": "../../source/postprocessors/SideIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = u
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
  []
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = false
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/porosity_jump/2d-rc-epsjump.i)

Other inputs

One of the tests to verify the conservation of mass, momentum and energy, features both a straight and a diverging channel. The test also features postprocessors to measure the flow of the conserved quantities on both internal and external boundaries. The geometry can be switched using the inactive parameter of the Mesh block.

mu=1
rho=1
advected_interp_method='average'
velocity_interp_method='rc'

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'mesh internal_boundary_bot internal_boundary_top'
  [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"}>>> = '1'
    dy<<<{"description": "Intervals in the Y direction (required when dim>1 otherwise ignored)"}>>> = '1 1 1'
    ix<<<{"description": "Number of grids in all intervals in the X direction (default to all one)"}>>> = '5'
    iy<<<{"description": "Number of grids in all intervals in the Y direction (default to all one)"}>>> = '5 5 5'
    subdomain_id<<<{"description": "Block IDs (default to all zero)"}>>> = '1
                    2
                    3'
  []
  [internal_boundary_bot]
    type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
    input<<<{"description": "The mesh we want to modify"}>>> = mesh
    new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = 'internal_bot'
    primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = 1
    paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = 2
  []
  [internal_boundary_top]
    type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
    input<<<{"description": "The mesh we want to modify"}>>> = internal_boundary_bot
    new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = 'internal_top'
    primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = 2
    paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = 3
  []
  [diverging_mesh]
    type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../../source/meshgenerators/FileMeshGenerator.html"}>>>
    file<<<{"description": "The filename to read."}>>> = 'expansion_quad.e'
  []
[]

[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
  fv_bcs_integrity_check = true
[]

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

[UserObjects<<<{"href": "../../syntax/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"}>>> = u
    v<<<{"description": "The y-component of velocity"}>>> = v
    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": "../../syntax/Variables/index.html"}>>>]
  [u]
    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"}>>> = 0
  []
  [v]
    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"}>>> = 1
  []
  [pressure]
    type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
  []
  [temperature]
    type = INSFVEnergyVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVEnergyVariable.html"}>>>
  []
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  [advected_density]
    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
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${rho}
  []
  [porosity]
    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
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.5
  []
[]

[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
  [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
    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"}>>> = u
    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"}>>> = u
    force_boundary_execution<<<{"description": "Whether to force execution of this object on all external boundaries."}>>> = true
    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
    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 = 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"}>>> = u
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
    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"}>>> = v
    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"}>>> = v
    force_boundary_execution<<<{"description": "Whether to force execution of this object on all external boundaries."}>>> = true
    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
    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 = 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"}>>> = v
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
    porosity<<<{"description": "Porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = porosity
  []

  [temp_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"}>>> = temperature
    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'."}>>> = 'upwind'
  []
  [temp_source]
    type = FVBodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into finite volume kernels, e.g. (controllable) constants, functions, and postprocessors.", "href": "../../source/fvkernels/FVBodyForce.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
    function<<<{"description": "A function that describes the body force"}>>> = 10
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
  []
[]

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

[FunctorMaterials<<<{"href": "../../syntax/FunctorMaterials/index.html"}>>>]
  [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."}>>> = 'temperature'
    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}
  []
  [advected_material_property]
    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"}>>> = 'advected_rho 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."}>>> ='${rho} 1'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      200                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [inlet_mass_variable]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = bottom
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = advected_density
  []
  [inlet_mass_constant]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = bottom
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
  [inlet_mass_matprop]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = bottom
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'advected_rho'
  []
  [mid1_mass]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = internal_bot
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
  [mid2_mass]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = internal_top
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []
  [outlet_mass]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = top
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
  []

  [inlet_momentum_x]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = bottom
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = u
  []

  [inlet_momentum_y]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = bottom
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = v
  []

  [mid1_advected_energy]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = internal_bot
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'rho_cp_temp'
    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'."}>>> = 'upwind'
  []
  [mid2_advected_energy]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = internal_top
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'rho_cp_temp'
    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'."}>>> = 'upwind'
  []
  [outlet_advected_energy]
    type = VolumetricFlowRate<<<{"description": "Computes the volumetric flow rate of an advected quantity through a sideset.", "href": "../../source/postprocessors/VolumetricFlowRate.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = top
    vel_x<<<{"description": "The x-axis velocity"}>>> = u
    vel_y<<<{"description": "The y-axis velocity"}>>> = v
    advected_quantity<<<{"description": "The quantity to advect. This is the canonical parameter to set the advected quantity when finite volume is being used. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'rho_cp_temp'
    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'."}>>> = 'upwind'
  []
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(modules/navier_stokes/test/tests/postprocessors/flow_rates/conservation_PINSFV.i)

Action syntax

To ease the burden of preparing long input files, the NavierStokesFV action syntax can also be used to set up porous medium INSFV simulations.

Pronghorn

Predefined temperature and pressure dependent material properties, closure relations for the effective viscosity, thermal conductivity, drag models and other models required for the modeling of nuclear reactors can be found in the Pronghorn application. Access to Pronghorn may be requested here by creating a software request.