Finite Volume Incompressible Porous media Navier Stokes

Equations

This module implements the porous media Navier Stokes equations. They are expressed in terms of the superficial viscosity where is the porosity and the interstitial velocity. The superficial viscosity 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:

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 viscosity, the specific heat capacity the convective heat transfer coefficient.

Implementation for a straight channel

We define a straight channel using a simple Cartesian mesh.

[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 5
    ymin = 0
    ymax = 1
    nx = 20
    ny = 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]
  inactive = 'lambda'
  [u]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [v]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1e-6
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [lambda]
    family = SCALAR
    order = FIRST
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)

The porous media Navier Stokes equation 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.

[mass]
  type = PINSFVMassAdvection
  variable = pressure
  advected_interp_method = ${advected_interp_method}
  velocity_interp_method = ${velocity_interp_method}
  vel = 'velocity'
  pressure = pressure
  u = u
  v = v
  mu = ${mu}
  rho = ${rho}
  porosity = porosity
[]
(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. Here the porosity is constant, so the porosity gradient term, defined by a PINSFVMomentumAdvectionPorosityGradient kernel is not included. 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.

[u_advection]
  type = PINSFVMomentumAdvection
  variable = u
  advected_quantity = 'rhou'
  vel = 'velocity'
  advected_interp_method = ${advected_interp_method}
  velocity_interp_method = ${velocity_interp_method}
  pressure = pressure
  u = u
  v = v
  mu = ${mu}
  rho = ${rho}
  porosity = porosity
[]
(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.

[u_viscosity]
  type = PINSFVMomentumDiffusion
  variable = u
  mu = ${mu}
  porosity = porosity
[]
(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.

[u_pressure]
  type = PINSFVMomentumPressure
  variable = u
  momentum_component = 'x'
  pressure = pressure
  porosity = 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 !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/no-slip-u

  • free slip walls !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/free-slip-u

  • symmetry axis. This symmetry condition should also be indicated for the pressure variable. !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/symmetry-u !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/symmetry-p

  • inlet velocity, to specify mass flux given that density is constant !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/inlet-u

  • momentum advection outflow (only for a mean-pressure approach, equivalent to executing the momentum advection kernel on the boundary) !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/outlet-u

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: !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/outlet-p

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. !listing modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i block=FVBCs/outlet-p-novalue

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]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 10
    ymin = 0
    ymax = 1
    nx = 100
    ny = 20
  []
[]

[GlobalParams]
  two_term_boundary_expansion = true
[]

[Variables]
  inactive = 'temp_solid'
  [u]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = ${u_inlet}
  []
  [v]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1e-6
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [temperature]
    type = INSFVEnergyVariable
  []
  [temp_solid]
    family = 'MONOMIAL'
    order = 'CONSTANT'
    fv = true
  []
[]

[AuxVariables]
  [temp_solid]
    family = 'MONOMIAL'
    order = 'CONSTANT'
    fv = true
    initial_condition = 100
  []
  [porosity]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_condition = 0.5
  []
[]

[FVKernels]
  inactive = 'solid_energy_diffusion solid_energy_convection'
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []

  [u_advection]
    type = PINSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = u
    mu = ${mu}
    porosity = porosity
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    pressure = pressure
    porosity = porosity
  []

  [v_advection]
    type = PINSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = v
    mu = ${mu}
    porosity = porosity
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    pressure = pressure
    porosity = porosity
  []

  [energy_advection]
    type = PINSFVEnergyAdvection
    variable = temperature
    vel = 'velocity'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [energy_diffusion]
    type = PINSFVEnergyDiffusion
    k = ${k}
    variable = temperature
    porosity = porosity
  []
  [energy_convection]
    type = PINSFVEnergyAmbientConvection
    variable = temperature
    is_solid = false
    temp_fluid = temperature
    temp_solid = temp_solid
    h_solid_fluid = 'h_cv'
  []

  [solid_energy_diffusion]
    type = FVDiffusion
    coeff = ${k}
    variable = temp_solid
  []
  [solid_energy_convection]
    type = PINSFVEnergyAmbientConvection
    variable = temp_solid
    is_solid = true
    temp_fluid = temperature
    temp_solid = temp_solid
    h_solid_fluid = 'h_cv'
  []
[]

[FVBCs]
  inactive = 'heated-side'
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = ${u_inlet}
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 0
  []
  [inlet-T]
    type = FVNeumannBC
    variable = temperature
    value = ${fparse u_inlet * rho * cp * T_inlet}
    boundary = 'left'
  []

  [no-slip-u]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = u
    function = 0
  []
  [no-slip-v]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = v
    function = 0
  []
  [heated-side]
    type = FVDirichletBC
    boundary = 'top'
    variable = 'temp_solid'
    value = 150
  []

  [symmetry-u]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = u
    u = u
    v = v
    mu = ${mu}
    momentum_component = 'x'
    porosity = porosity
  []
  [symmetry-v]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = v
    u = u
    v = v
    mu = ${mu}
    momentum_component = 'y'
    porosity = porosity
  []
  [symmetry-p]
    type = INSFVSymmetryPressureBC
    boundary = 'bottom'
    variable = pressure
  []

  [outlet-p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0.1
  []
[]

[Materials]
  [constants]
    type = ADGenericConstantMaterial
    prop_names = 'cp h_cv'
    prop_values = '${cp} 1'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    rho = ${rho}
    temperature = 'temperature'
  []
[]

[Executioner]
  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      100                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]

# Some basic Postprocessors to examine the solution
[Postprocessors]
  [inlet-p]
    type = SideAverageValue
    variable = pressure
    boundary = 'left'
  []
  [outlet-u]
    type = SideAverageValue
    variable = u
    boundary = 'right'
  []
  [outlet-temp]
    type = SideAverageValue
    variable = temperature
    boundary = 'right'
  []
  [solid-temp]
    type = ElementAverageValue
    variable = temp_solid
  []
[]

[Outputs]
  exodus = true
  csv = 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]
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1 1'
    dy = '0.5'
    ix = '30 30'
    iy = '20'
    subdomain_id = '1 2'
  []
[]

[GlobalParams]
  two_term_boundary_expansion = true
[]

[Variables]
  [u]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [v]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1e-6
  []
  [pressure]
    type = INSFVPressureVariable
  []
[]

[AuxVariables]
  [porosity]
    family = MONOMIAL
    order = CONSTANT
    fv = true
  []
[]

[ICs]
  inactive = 'porosity_continuous'
  [porosity_1]
    type = ConstantIC
    variable = porosity
    block = 1
    value = 1
  []
  [porosity_2]
    type = ConstantIC
    variable = porosity
    block = 2
    value = 0.5
  []
  [porosity_continuous]
    type = FunctionIC
    variable = porosity
    block = '1 2'
    function = smooth_jump
  []
[]

[Functions]
  [smooth_jump]
    type = ParsedFunction
    value = '1 - 0.5 * 1 / (1 + exp(-30*(x-1)))'
  []
[]

[FVKernels]
  inactive = 'u_advection_porosity_gradient v_advection_porosity_gradient'
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []

  [u_advection]
    type = PINSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = u
    mu = ${mu}
    porosity = porosity

    momentum_component = 'x'
    vel = 'velocity'
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = u
    pressure = pressure
    porosity = porosity
    momentum_component = 'x'
  []
  [u_advection_porosity_gradient]
    type = PINSFVMomentumAdvectionPorosityGradient
    variable = u
    u = u
    v = v
    rho = ${rho}
    porosity = porosity
    momentum_component = 'x'
  []

  [v_advection]
    type = PINSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = v
    mu = ${mu}
    porosity = porosity

    momentum_component = 'y'
    vel = 'velocity'
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = v
    pressure = pressure
    porosity = porosity
    momentum_component = 'y'
  []
  [v_advection_porosity_gradient]
    type = PINSFVMomentumAdvectionPorosityGradient
    variable = v
    u = u
    v = v
    rho = ${rho}
    porosity = porosity
    momentum_component = 'y'
  []
[]

[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 0
  []

  [walls-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top bottom'
    variable = u
  []
  [walls-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top bottom'
    variable = v
  []

  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0.4
  []
[]

[Materials]
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    rho = ${rho}
  []
[]

[Executioner]
  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      100                lu           NONZERO'
  line_search = 'none'
[]

[Postprocessors]
  [inlet_p]
    type = SideAverageValue
    variable = 'pressure'
    boundary = 'left'
  []
  [outlet-u]
    type = SideIntegralVariablePostprocessor
    variable = u
    boundary = 'right'
  []
[]

[Outputs]
  exodus = true
  csv = 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.1
rho=1
advected_interp_method='average'
velocity_interp_method='average'

[Mesh]
  inactive = 'mesh internal_boundary_bot internal_boundary_top'
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1'
    dy = '1 1 1'
    ix = '5'
    iy = '5 5 5'
    subdomain_id = '1
                    2
                    3'
  []
  [internal_boundary_bot]
    type = SideSetsBetweenSubdomainsGenerator
    input = mesh
    new_boundary = 'internal_bot'
    primary_block = 1
    paired_block = 2
  []
  [internal_boundary_top]
    type = SideSetsBetweenSubdomainsGenerator
    input = internal_boundary_bot
    new_boundary = 'internal_top'
    primary_block = 2
    paired_block = 3
  []
  [diverging_mesh]
    type = FileMeshGenerator
    file = 'expansion_quad.e'
  []
[]

[Problem]
  kernel_coverage_check = false
  fv_bcs_integrity_check = true
[]

[GlobalParams]
  two_term_boundary_expansion = true
[]

[Variables]
  [u]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 0
  []
  [v]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [temperature]
    type = INSFVEnergyVariable
  []
[]

[AuxVariables]
  [advected_density]
    order = CONSTANT
    family = MONOMIAL
    fv = true
    initial_condition = ${rho}
  []
  [porosity]
    order = CONSTANT
    family = MONOMIAL
    fv = true
    initial_condition = 0.5
  []
[]

[FVKernels]
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []

  [u_advection]
    type = PINSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = u
    force_boundary_execution = true
    porosity = porosity
    mu = ${mu}
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    pressure = pressure
    porosity = porosity
  []

  [v_advection]
    type = PINSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = v
    force_boundary_execution = true
    porosity = porosity
    mu = ${mu}
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    pressure = pressure
    porosity = porosity
  []

  [temp_advection]
    type = PINSFVEnergyAdvection
    vel = 'velocity'
    variable = temperature
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
    porosity = porosity
  []
  [temp_source]
    type = FVBodyForce
    variable = temperature
    function = 10
    block = 1
  []
[]

[FVBCs]
  inactive = 'noslip-u noslip-v'
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'bottom'
    variable = u
    function = 0
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'bottom'
    variable = v
    function = 1
  []
  [noslip-u]
    type = INSFVNoSlipWallBC
    boundary = 'right'
    variable = u
    function = 0
  []
  [noslip-v]
    type = INSFVNoSlipWallBC
    boundary = 'right'
    variable = v
    function = 0
  []
  [free-slip-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'right'
    variable = u
  []
  [free-slip-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'right'
    variable = v
  []
  [axis-u]
    type = PINSFVSymmetryVelocityBC
    boundary = 'left'
    variable = u
    u = u
    v = v
    mu = ${mu}
    momentum_component = x
    porosity = porosity
  []
  [axis-v]
    type = PINSFVSymmetryVelocityBC
    boundary = 'left'
    variable = v
    u = u
    v = v
    mu = ${mu}
    momentum_component = y
    porosity = porosity
  []
  [axis-p]
    type = INSFVSymmetryPressureBC
    boundary = 'left'
    variable = pressure
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'top'
    variable = pressure
    function = 0
  []
  [inlet_temp]
    type = FVNeumannBC
    boundary = 'bottom'
    variable = temperature
    value = 300
  []
[]

[Materials]
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'temperature'
    rho = ${rho}
  []
  [advected_material_property]
    type = ADGenericConstantMaterial
    prop_names = 'advected_rho cp'
    prop_values ='${rho} 1'
  []
[]

[Executioner]
  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]
  [inlet_mass_variable]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = advected_density
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_mass_constant]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = ${rho}
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_mass_matprop]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_mat_prop = 'advected_rho'
    fv = true
  []
  [mid1_mass]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_mass]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_mass]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []

  [inlet_momentum_x]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid1_momentum_x]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_momentum_x]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_momentum_x]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []

  [inlet_momentum_y]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid1_momentum_y]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_momentum_y]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_momentum_y]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []

  [inlet_advected_energy]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
  []
  [mid1_advected_energy]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
  []
  [mid2_advected_energy]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
  []
  [outlet_advected_energy]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
  []
[]

[Outputs]
  exodus = false
  csv = true
  inactive = 'console_mass console_momentum_x console_momentum_y console_energy'
  [console_mass]
    type = Console
    start_step = 1
    show = 'inlet_mass_variable inlet_mass_constant inlet_mass_matprop mid1_mass mid2_mass outlet_mass'
  []
  [console_momentum_x]
    type = Console
    start_step = 1
    show = 'inlet_momentum_x mid1_momentum_x mid2_momentum_x outlet_momentum_x'
  []
  [console_momentum_y]
    type = Console
    start_step = 1
    show = 'inlet_momentum_y mid1_momentum_y mid2_momentum_y outlet_momentum_y'
  []
  [console_energy]
    type = Console
    start_step = 1
    show = 'inlet_advected_energy mid1_advected_energy mid2_advected_energy outlet_advected_energy'
  []
[]
(modules/navier_stokes/test/tests/postprocessors/conservation_PINSFV.i)

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.