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 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:
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
[]
[]
(moose/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'
[superficial_vel_x]
type = PINSFVSuperficialVelocityVariable
initial_condition = 1
[]
[superficial_vel_y]
type = PINSFVSuperficialVelocityVariable
initial_condition = 1e-6
[]
[pressure]
type = INSFVPressureVariable
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
(moose/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
rho = ${rho}
[]
(moose/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.
[u_advection]
type = PINSFVMomentumAdvection
variable = superficial_vel_x
rho = ${rho}
porosity = porosity
momentum_component = 'x'
[]
(moose/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 = superficial_vel_x
mu = ${mu}
porosity = porosity
momentum_component = 'x'
[]
(moose/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 = superficial_vel_x
momentum_component = 'x'
pressure = pressure
porosity = porosity
[]
(moose/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
(moose/modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)[no-slip-u] type = INSFVNoSlipWallBC boundary = 'top bottom' variable = superficial_vel_x function = 0 []
free slip walls
(moose/modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)[free-slip-u] type = INSFVNaturalFreeSlipBC boundary = 'top bottom' variable = superficial_vel_x momentum_component = 'x' []
symmetry axis. This symmetry condition should also be indicated for the pressure variable.
(moose/modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)[symmetry-u] type = PINSFVSymmetryVelocityBC boundary = 'bottom' variable = superficial_vel_x u = superficial_vel_x v = superficial_vel_y mu = ${mu} momentum_component = 'x' []
(moose/modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)[symmetry-p] type = INSFVSymmetryPressureBC boundary = 'bottom' variable = pressure []
inlet velocity, to specify mass flux given that density is constant
(moose/modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)[inlet-u] type = INSFVInletVelocityBC boundary = 'left' variable = superficial_vel_x function = '1' []
momentum advection outflow (only for a mean-pressure approach, equivalent to executing the momentum advection kernel on the boundary)
(moose/modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc.i)[outlet-u] type = PINSFVMomentumAdvectionOutflowBC boundary = 'right' variable = superficial_vel_x u = superficial_vel_x v = superficial_vel_y porosity = porosity momentum_component = 'x' rho = ${rho} []
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:
[outlet-p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = 0
[]
(moose/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.
[outlet-p-novalue]
type = INSFVMassAdvectionOutflowBC
boundary = 'right'
variable = pressure
u = superficial_vel_x
v = superficial_vel_y
rho = ${rho}
[]
(moose/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]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '5 5'
dy = '1.0'
ix = '50 50'
iy = '20'
subdomain_id = '1 2'
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = superficial_vel_x
v = superficial_vel_y
pressure = pressure
porosity = porosity
[]
[]
[Variables]
inactive = 'T_solid'
[superficial_vel_x]
type = PINSFVSuperficialVelocityVariable
initial_condition = ${u_inlet}
[]
[superficial_vel_y]
type = PINSFVSuperficialVelocityVariable
initial_condition = 1e-6
[]
[pressure]
type = INSFVPressureVariable
[]
[T_fluid]
type = INSFVEnergyVariable
[]
[T_solid]
family = 'MONOMIAL'
order = 'CONSTANT'
fv = true
[]
[]
[AuxVariables]
[T_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}
rho = ${rho}
[]
[u_advection]
type = PINSFVMomentumAdvection
variable = superficial_vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
porosity = porosity
momentum_component = 'x'
[]
[u_viscosity]
type = PINSFVMomentumDiffusion
variable = superficial_vel_x
mu = ${mu}
porosity = porosity
momentum_component = 'x'
[]
[u_pressure]
type = PINSFVMomentumPressure
variable = superficial_vel_x
momentum_component = 'x'
pressure = pressure
porosity = porosity
[]
[v_advection]
type = PINSFVMomentumAdvection
variable = superficial_vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
porosity = porosity
momentum_component = 'y'
[]
[v_viscosity]
type = PINSFVMomentumDiffusion
variable = superficial_vel_y
mu = ${mu}
porosity = porosity
momentum_component = 'y'
[]
[v_pressure]
type = PINSFVMomentumPressure
variable = superficial_vel_y
momentum_component = 'y'
pressure = pressure
porosity = porosity
[]
[energy_advection]
type = PINSFVEnergyAdvection
variable = T_fluid
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[energy_diffusion]
type = PINSFVEnergyDiffusion
k = ${k}
variable = T_fluid
porosity = porosity
[]
[energy_convection]
type = PINSFVEnergyAmbientConvection
variable = T_fluid
is_solid = false
T_fluid = 'T_fluid'
T_solid = 'T_solid'
h_solid_fluid = 'h_cv'
[]
[solid_energy_diffusion]
type = FVDiffusion
coeff = ${k}
variable = T_solid
[]
[solid_energy_convection]
type = PINSFVEnergyAmbientConvection
variable = T_solid
is_solid = true
T_fluid = 'T_fluid'
T_solid = 'T_solid'
h_solid_fluid = 'h_cv'
[]
[]
[FVBCs]
inactive = 'heated-side'
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = superficial_vel_x
function = ${u_inlet}
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = superficial_vel_y
function = 0
[]
[inlet-T]
type = FVNeumannBC
variable = T_fluid
value = '${fparse u_inlet * rho * cp * T_inlet}'
boundary = 'left'
[]
[no-slip-u]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = superficial_vel_x
function = 0
[]
[no-slip-v]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = superficial_vel_y
function = 0
[]
[heated-side]
type = FVDirichletBC
boundary = 'top'
variable = 'T_solid'
value = 150
[]
[symmetry-u]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = superficial_vel_x
u = superficial_vel_x
v = superficial_vel_y
mu = ${mu}
momentum_component = 'x'
[]
[symmetry-v]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = superficial_vel_y
u = superficial_vel_x
v = superficial_vel_y
mu = ${mu}
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[outlet-p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = 0.1
[]
[]
[Materials]
[constants]
type = ADGenericFunctorMaterial
prop_names = 'h_cv'
prop_values = '1'
[]
[functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp'
prop_values = '${cp}'
[]
[ins_fv]
type = INSFVEnthalpyMaterial
rho = ${rho}
temperature = 'T_fluid'
[]
[]
[Executioner]
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]
[inlet-p]
type = SideAverageValue
variable = pressure
boundary = 'left'
[]
[outlet-u]
type = SideAverageValue
variable = superficial_vel_x
boundary = 'right'
[]
[outlet-temp]
type = SideAverageValue
variable = T_fluid
boundary = 'right'
[]
[solid-temp]
type = ElementAverageValue
variable = T_solid
[]
[]
[Outputs]
exodus = true
csv = false
[]
(moose/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]
rhie_chow_user_object = 'rc'
porosity = porosity
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = u
v = v
porosity = porosity
pressure = pressure
[]
[]
[Variables]
[u]
type = PINSFVSuperficialVelocityVariable
initial_condition = 1
[]
[v]
type = PINSFVSuperficialVelocityVariable
initial_condition = 1e-6
[]
[pressure]
type = INSFVPressureVariable
[]
[]
[AuxVariables]
[porosity]
type = MooseVariableFVReal
[]
[]
[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
expression = '1 - 0.5 * 1 / (1 + exp(-30*(x-1)))'
[]
[]
[FVKernels]
[mass]
type = PINSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_advection]
type = PINSFVMomentumAdvection
variable = u
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = PINSFVMomentumDiffusion
variable = u
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = PINSFVMomentumPressure
variable = u
pressure = pressure
momentum_component = 'x'
[]
[v_advection]
type = PINSFVMomentumAdvection
variable = v
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = PINSFVMomentumDiffusion
variable = v
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = PINSFVMomentumPressure
variable = v
pressure = pressure
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
momentum_component = 'x'
[]
[walls-v]
type = INSFVNaturalFreeSlipBC
boundary = 'top bottom'
variable = v
momentum_component = 'y'
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = 0.4
[]
[]
[Materials]
inactive = 'smooth'
[jump]
type = ADPiecewiseByBlockFunctorMaterial
prop_name = 'porosity'
subdomain_to_prop_value = '1 1
2 0.5'
[]
[smooth]
type = ADGenericFunctionFunctorMaterial
prop_names = 'porosity'
prop_values = 'smooth_jump'
[]
[]
[Executioner]
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]
[inlet_p]
type = SideAverageValue
variable = 'pressure'
boundary = 'left'
[]
[outlet-u]
type = SideIntegralVariablePostprocessor
variable = u
boundary = 'right'
[]
[]
[Outputs]
exodus = true
csv = false
[]
(moose/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]
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]
fv_bcs_integrity_check = true
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
porosity = porosity
[]
[]
[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
rho = ${rho}
[]
[u_advection]
type = PINSFVMomentumAdvection
variable = u
rho = ${rho}
porosity = porosity
momentum_component = 'x'
[]
[u_viscosity]
type = PINSFVMomentumDiffusion
variable = u
force_boundary_execution = true
porosity = porosity
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = PINSFVMomentumPressure
variable = u
momentum_component = 'x'
pressure = pressure
porosity = porosity
[]
[v_advection]
type = PINSFVMomentumAdvection
variable = v
rho = ${rho}
porosity = porosity
momentum_component = 'y'
[]
[v_viscosity]
type = PINSFVMomentumDiffusion
variable = v
force_boundary_execution = true
porosity = porosity
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = PINSFVMomentumPressure
variable = v
momentum_component = 'y'
pressure = pressure
porosity = porosity
[]
[temp_advection]
type = PINSFVEnergyAdvection
variable = temperature
advected_interp_method = 'upwind'
[]
[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
momentum_component = 'x'
[]
[free-slip-v]
type = INSFVNaturalFreeSlipBC
boundary = 'right'
variable = v
momentum_component = 'y'
[]
[axis-u]
type = PINSFVSymmetryVelocityBC
boundary = 'left'
variable = u
u = u
v = v
mu = ${mu}
momentum_component = x
[]
[axis-v]
type = PINSFVSymmetryVelocityBC
boundary = 'left'
variable = v
u = u
v = v
mu = ${mu}
momentum_component = y
[]
[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 = INSFVEnthalpyMaterial
temperature = 'temperature'
rho = ${rho}
[]
[advected_material_property]
type = ADGenericFunctorMaterial
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_quantity = advected_density
[]
[inlet_mass_constant]
type = VolumetricFlowRate
boundary = bottom
vel_x = u
vel_y = v
advected_quantity = ${rho}
[]
[inlet_mass_matprop]
type = VolumetricFlowRate
boundary = bottom
vel_x = u
vel_y = v
advected_quantity = 'advected_rho'
[]
[mid1_mass]
type = VolumetricFlowRate
boundary = internal_bot
vel_x = u
vel_y = v
advected_quantity = ${rho}
[]
[mid2_mass]
type = VolumetricFlowRate
boundary = internal_top
vel_x = u
vel_y = v
advected_quantity = ${rho}
[]
[outlet_mass]
type = VolumetricFlowRate
boundary = top
vel_x = u
vel_y = v
advected_quantity = ${rho}
[]
[inlet_momentum_x]
type = VolumetricFlowRate
boundary = bottom
vel_x = u
vel_y = v
advected_quantity = u
[]
[inlet_momentum_y]
type = VolumetricFlowRate
boundary = bottom
vel_x = u
vel_y = v
advected_quantity = v
[]
[mid1_advected_energy]
type = VolumetricFlowRate
boundary = internal_bot
vel_x = u
vel_y = v
advected_quantity = 'rho_cp_temp'
advected_interp_method = 'upwind'
[]
[mid2_advected_energy]
type = VolumetricFlowRate
boundary = internal_top
vel_x = u
vel_y = v
advected_quantity = 'rho_cp_temp'
advected_interp_method = 'upwind'
[]
[outlet_advected_energy]
type = VolumetricFlowRate
boundary = top
vel_x = u
vel_y = v
advected_quantity = 'rho_cp_temp'
advected_interp_method = 'upwind'
[]
[]
[Outputs]
csv = true
[]
(moose/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.