- T_fluidFluid temperature. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:Fluid temperature. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- characteristic_lengthcharacteristic length for Reynolds number calculation. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:characteristic length for Reynolds number calculation. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- fpFluid properties functor userobject
C++ Type:UserObjectName
Unit:(no unit assumed)
Controllable:No
Description:Fluid properties functor userobject
- porosityporosity. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:porosity. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- pressurePressure. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:Pressure. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- speedVelocity norm. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:Velocity norm. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
GeneralFunctorFluidProps
Creates functor fluid properties using a (P, T) formulation
Overview
This object uses a SinglePhaseFluidProperties
derived-object to compute the following properties:
specific heat at constant volume,
specific heat at constant pressure,
dynamic viscosity,
thermal conductivity,
Prandtl number,
pore/particle Reynolds number
hydraulic Reynolds number
interstitial Reynolds number
the time derivatives of the:
specific heat at constant pressure,
density
and the pressure and temperature derivatives of the:
specific heat at constant pressure,
density
dynamic viscosity,
thermal conductivity,
Prandtl number,
pore Reynolds number
In order to use this with some fluid properties that do not compute the AD version of the density derivatives, such as the Spline Base Table Lookup fluid properties, you can use the "neglect_derivatives_of_density_time_derivative" to neglect the derivatives with regards to the nonlinear variables (usually pressure, temperature) of the time derivative of the density.
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Unit:(no unit assumed)
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Unit:(no unit assumed)
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
Default:NONE
C++ Type:MooseEnum
Unit:(no unit assumed)
Options:NONE, ELEMENT, SUBDOMAIN
Controllable:No
Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
- declare_suffixAn optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
- execute_onALWAYSThe list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html.
Default:ALWAYS
C++ Type:ExecFlagEnum
Unit:(no unit assumed)
Options:XFEM_MARK, FORWARD, ADJOINT, HOMOGENEOUS_FORWARD, ADJOINT_TIMESTEP_BEGIN, ADJOINT_TIMESTEP_END, NONE, INITIAL, LINEAR, NONLINEAR_CONVERGENCE, NONLINEAR, POSTCHECK, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, ALWAYS
Controllable:No
Description:The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html.
- force_define_densityFalseWhether to force the definition of a density functor from the fluid properties
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether to force the definition of a density functor from the fluid properties
- mu_rampdown1A function describing a ramp down of viscosity over time
Default:1
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:A function describing a ramp down of viscosity over time
- neglect_derivatives_of_density_time_derivativeFalseWhether to neglect the derivatives with regards to nonlinear variables of the density time derivatives
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether to neglect the derivatives with regards to nonlinear variables of the density time derivatives
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- rhoDensity. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:Density. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:The seed for the master random number generator
Advanced Parameters
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:List of material properties, from this material, to output (outputs must also be defined to an output type)
- outputsnone Vector of output names where you would like to restrict the output of variables(s) associated with this object
Default:none
C++ Type:std::vector<OutputName>
Unit:(no unit assumed)
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/turbulent_driven_growth.i)
- (modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-physics.i)
- (modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-gas.i)
- (modules/navier_stokes/test/tests/finite_volume/materials/ergun/ergun.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/materials/2d-transient.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/natural_convection/natural_circulation_pipe.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/materials/functorfluidprops.i)
- (modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient.i)
- (modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-action.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/pressure_driven_growth_transient.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/pressure_driven_growth.i)
neglect_derivatives_of_density_time_derivative
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether to neglect the derivatives with regards to nonlinear variables of the density time derivatives
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/turbulent_driven_growth.i)
###############################################################################
# Validation test based on Hibiki and Ishii experiment [1] reported in Figure 5
# [1] Hibiki, T., & Ishii, M. (2000). One-group interfacial area transport of
# bubbly flows in vertical round tubes.
# International Journal of Heat and Mass Transfer, 43(15), 2711-2726.
###############################################################################
mu = 1.0
rho = 1000.0
mu_d = 1.0
rho_d = 1.0
l = ${fparse 50.8/1000.0}
U = 5.031429
dp = 0.005
inlet_phase_2 = 0.442
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
mass_exchange_coeff = 0.0
inlet_interface_area = ${fparse 6.0*inlet_phase_2/dp}
outlet_pressure = 1e5
[GlobalParams]
rhie_chow_user_object = 'rc'
density_interp_method = 'average'
mu_interp_method = 'average'
[]
[Problem]
identify_variable_groups_in_nl = false
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
coord_type = 'RZ'
rz_coord_axis = 'X'
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = '${fparse l * 60}'
ymin = 0
ymax = '${fparse l / 2}'
nx = 20
ny = 5
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_phase_2}
[]
[interface_area]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_interface_area}
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_x
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_y
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_x'
v_slip = 'vel_y'
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1.0
[]
[phase_2_src]
type = NSFVMixturePhaseInterface
variable = phase_2
phase_coupled = phase_1
alpha = ${mass_exchange_coeff}
[]
[interface_area_advection]
type = INSFVScalarFieldAdvection
variable = interface_area
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[interface_area_diffusion]
type = FVDiffusion
variable = interface_area
coeff = 0.1
[]
[interface_area_source_sink]
type = WCNSFV2PInterfaceAreaSourceSink
variable = interface_area
u = 'vel_x'
v = 'vel_y'
L = ${fparse l/2}
rho = 'rho_mixture'
rho_d = 'rho'
pressure = 'pressure'
k_c = '${fparse mass_exchange_coeff}'
fd = 'phase_2'
sigma = 1e-3
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '${outlet_pressure}'
[]
[inlet_phase_2]
type = FVDirichletBC
boundary = 'left'
variable = phase_2
value = ${inlet_phase_2}
[]
[inlet_interface_area]
type = FVDirichletBC
boundary = 'left'
variable = interface_area
value = ${inlet_interface_area}
[]
[symmetry-u]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_x
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'x'
[]
[symmetry-v]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_y
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[symmetry-phase-2]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = phase_2
[]
[symmetry-interface-area]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = interface_area
[]
[]
[AuxVariables]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[FunctorMaterials]
[bubble_properties]
type = GeneralFunctorFluidProps
fp = 'fp'
pressure = 'pressure'
T_fluid = 300.0
speed = 1.0
characteristic_length = 1.0
porosity = 1.0
output_properties = 'rho'
outputs = 'out'
[]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${rho} ${mu}'
phase_1_names = 'rho ${mu_d}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
line_search = 'none'
[]
[Debug]
show_var_residual_norms = true
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
[out]
type = Exodus
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${l} * ${U}'
pp_names = ''
[]
[rho_outlet]
type = SideAverageValue
boundary = 'right'
variable = 'rho_mixture_var'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-physics.i)
# Solid properties
cp_s = 2
rho_s = 4
k_s = '${fparse 1e-2 / 0.5}'
h_fs = 10
# thermal diffusivity is divided by 0.5 to match the reference using the action
# Operating conditions
u_inlet = 1
T_inlet = 200
p_outlet = 10
top_side_temperature = 150
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 1
nx = 20
ny = 5
[]
[]
[AuxVariables]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.5
[]
[velocity_norm]
type = MooseVariableFVReal
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[Physics]
[NavierStokes]
[Flow]
[flow]
compressibility = 'weakly-compressible'
porous_medium_treatment = true
define_variables = true
pressure_variable = 'pressure'
density = 'rho'
dynamic_viscosity = 'mu'
initial_velocity = '${u_inlet} 1e-6 0'
initial_pressure = '${p_outlet}'
inlet_boundaries = 'left'
momentum_inlet_types = 'fixed-velocity'
momentum_inlet_function = '${u_inlet} 0'
wall_boundaries = 'top bottom'
momentum_wall_types = 'noslip symmetry'
outlet_boundaries = 'right'
momentum_outlet_types = 'fixed-pressure'
pressure_function = '${p_outlet}'
mass_advection_interpolation = 'average'
momentum_advection_interpolation = 'average'
[]
[]
[FluidHeatTransfer]
[fluid]
thermal_conductivity = 'k'
effective_conductivity = true
specific_heat = 'cp'
initial_temperature = '${T_inlet}'
# See 'flow' for inlet boundaries
energy_inlet_types = 'fixed-temperature'
energy_inlet_function = '${T_inlet}'
# See 'flow' for wall boundaries
energy_wall_types = 'heatflux heatflux'
energy_wall_function = '0 0'
ambient_convection_alpha = 'h_cv'
ambient_temperature = 'T_solid'
energy_advection_interpolation = 'average'
[]
[]
[SolidHeatTransfer]
[solid]
block = 0
initial_temperature = 100
transient = true
# To match the previous test results
solid_temperature_two_term_bc_expansion = true
thermal_conductivity_solid = '${k_s}'
cp_solid = ${cp_s}
rho_solid = ${rho_s}
fixed_temperature_boundaries = 'top'
boundary_temperatures = '${top_side_temperature}'
ambient_convection_alpha = 'h_cv'
ambient_convection_temperature = 'T_fluid'
verbose = true
[]
[]
[]
[]
[FunctorMaterials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'h_cv'
prop_values = '${h_fs}'
[]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T_fluid'
speed = 'velocity_norm'
# To initialize with a high viscosity
mu_rampdown = 'mu_rampdown'
# For porous flow
characteristic_length = 1
porosity = 'porosity'
[]
[]
[Functions]
[mu_rampdown]
type = PiecewiseLinear
x = '1 2 3 4'
y = '1e3 1e2 1e1 1'
[]
[]
[AuxKernels]
[speed]
type = ParsedAux
variable = 'velocity_norm'
coupled_variables = 'superficial_vel_x superficial_vel_y porosity'
expression = 'sqrt(superficial_vel_x*superficial_vel_x + superficial_vel_y*superficial_vel_y) / porosity'
[]
[]
[Executioner]
type = Transient
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
end_time = 3.0
[]
# 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
[]
(modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-gas.i)
# Fluid properties
mu = 'mu'
rho = 'rho'
k = 'k'
# Solid properties
cp_s = 2
rho_s = 4
k_s = 1e-2
h_fs = 10
# Operating conditions
u_inlet = 1
T_inlet = 200
p_outlet = 10
top_side_temperature = 150
# Numerical scheme
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 1
nx = 20
ny = 5
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = superficial_vel_x
v = superficial_vel_y
pressure = pressure
porosity = porosity
[]
[]
[Variables]
[superficial_vel_x]
type = PINSFVSuperficialVelocityVariable
initial_condition = ${u_inlet}
[]
[superficial_vel_y]
type = PINSFVSuperficialVelocityVariable
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${p_outlet}
[]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = ${T_inlet}
[]
[T_solid]
type = MooseVariableFVReal
initial_condition = 100
[]
[]
[AuxVariables]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.5
[]
[]
[FVKernels]
[mass_time]
type = PWCNSFVMassTimeDerivative
variable = pressure
porosity = 'porosity'
drho_dt = 'drho_dt'
[]
[mass]
type = PWCNSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = superficial_vel_x
rho = ${rho}
drho_dt = 'drho_dt'
momentum_component = 'x'
[]
[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_time]
type = WCNSFVMomentumTimeDerivative
variable = superficial_vel_y
rho = ${rho}
drho_dt = 'drho_dt'
momentum_component = 'y'
[]
[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_time]
type = PINSFVEnergyTimeDerivative
variable = T_fluid
h = 'h'
dh_dt = 'dh_dt'
rho = ${rho}
drho_dt = 'drho_dt'
is_solid = false
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
variable = T_fluid
k = ${k}
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_time]
type = PINSFVEnergyTimeDerivative
variable = T_solid
cp = ${cp_s}
rho = ${rho_s}
is_solid = true
porosity = porosity
[]
[solid_energy_diffusion]
type = FVDiffusion
variable = T_solid
coeff = ${k_s}
[]
[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]
[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 = FVDirichletBC
variable = T_fluid
value = ${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 = ${top_side_temperature}
[]
[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 = ${p_outlet}
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
[]
[]
[FunctorMaterials]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T_fluid'
speed = 'speed'
# To initialize with a high viscosity
mu_rampdown = 'mu_rampdown'
# For porous flow
characteristic_length = 1
porosity = 'porosity'
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
rho = ${rho}
temperature = 'T_fluid'
[]
[constants]
type = ADGenericFunctorMaterial
prop_names = 'h_cv'
prop_values = '${h_fs}'
[]
[speed]
type = PINSFVSpeedFunctorMaterial
porosity = 'porosity'
superficial_vel_x = 'superficial_vel_x'
superficial_vel_y = 'superficial_vel_y'
[]
[]
[Functions]
[mu_rampdown]
type = PiecewiseLinear
x = '1 2 3 4'
y = '1e3 1e2 1e1 1'
[]
[]
[Executioner]
type = Transient
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
nl_abs_tol = 1e-10
automatic_scaling = true
end_time = 3.0
[]
# Some basic Postprocessors to examine the solution
[Postprocessors]
[inlet-p]
type = SideAverageValue
variable = pressure
boundary = 'left'
[]
[outlet-u]
type = VolumetricFlowRate
boundary = 'right'
advected_quantity = '1'
advected_interp_method = ${advected_interp_method}
vel_x = 'superficial_vel_x'
vel_y = 'superficial_vel_y'
[]
[outlet-temp]
type = SideAverageValue
variable = T_fluid
boundary = 'right'
[]
[solid-temp]
type = ElementAverageValue
variable = T_solid
[]
[]
[Outputs]
exodus = true
csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/materials/ergun/ergun.i)
# This file simulates flow of fluid in a porous elbow for the purpose of verifying
# correct implementation of the various different solution variable sets. This input
# tests correct implementation of the primitive superficial variable set. Flow enters on the top
# and exits on the right. Because the purpose is only to test the equivalence of
# different equation sets, no solid energy equation is included.
porosity_left = 0.4
porosity_right = 0.6
pebble_diameter = 0.06
mu = 1.81e-5 # This has been increased to avoid refining the mesh
M = 28.97e-3
R = 8.3144598
# inlet mass flowrate, kg/s
mdot = -10.0
# inlet mass flux (superficial)
mflux_in_superficial = ${fparse mdot / (pi * 0.5 * 0.5)}
# inlet mass flux (interstitial)
mflux_in_interstitial = ${fparse mflux_in_superficial / porosity_left}
p_initial = 201325.0
T_initial = 300.0
rho_initial = ${fparse p_initial / T_initial * M / R}
vel_y_initial = ${fparse mflux_in_interstitial / rho_initial}
vel_x_initial = 0.0
superficial_vel_y_initial = ${fparse mflux_in_superficial / rho_initial}
superficial_vel_x_initial = 1e-12
# Computation parameters
velocity_interp_method = 'rc'
advected_interp_method = 'upwind'
# ==============================================================================
# GEOMETRY AND MESH
# ==============================================================================
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'ergun_in.e'
[]
coord_type = RZ
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = superficial_vel_x
v = superficial_vel_y
pressure = pressure
porosity = porosity
[]
[]
[GlobalParams]
porosity = porosity
pebble_diameter = ${pebble_diameter}
fp = fp
# rho for the kernels. Must match fluid property!
rho = ${rho_initial}
fv = true
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
# behavior at time of test creation
two_term_boundary_expansion = false
rhie_chow_user_object = 'rc'
[]
# ==============================================================================
# VARIABLES AND KERNELS
# ==============================================================================
[Variables]
[pressure]
type = INSFVPressureVariable
initial_condition = ${p_initial}
[]
[superficial_vel_x]
type = PINSFVSuperficialVelocityVariable
initial_condition = ${superficial_vel_x_initial}
[]
[superficial_vel_y]
type = PINSFVSuperficialVelocityVariable
initial_condition = ${superficial_vel_y_initial}
[]
[]
[FVKernels]
# Mass Equation.
[mass]
type = PINSFVMassAdvection
variable = 'pressure'
[]
# Momentum x component equation.
[vel_x_time]
type = PINSFVMomentumTimeDerivative
variable = 'superficial_vel_x'
momentum_component = 'x'
[]
[vel_x_advection]
type = PINSFVMomentumAdvection
variable = 'superficial_vel_x'
momentum_component = 'x'
[]
[vel_x_viscosity]
type = PINSFVMomentumDiffusion
variable = 'superficial_vel_x'
momentum_component = 'x'
mu = 'mu'
[]
[u_pressure]
type = PINSFVMomentumPressure
variable = 'superficial_vel_x'
pressure = pressure
momentum_component = 'x'
[]
[u_friction]
type = PINSFVMomentumFriction
variable = 'superficial_vel_x'
Darcy_name = 'Darcy_coefficient'
Forchheimer_name = 'Forchheimer_coefficient'
momentum_component = 'x'
speed = speed
mu = 'mu'
[]
# Momentum y component equation.
[vel_y_time]
type = PINSFVMomentumTimeDerivative
variable = 'superficial_vel_y'
momentum_component = 'y'
[]
[vel_y_advection]
type = PINSFVMomentumAdvection
variable = 'superficial_vel_y'
momentum_component = 'y'
[]
[vel_y_viscosity]
type = PINSFVMomentumDiffusion
variable = 'superficial_vel_y'
momentum_component = 'y'
mu = 'mu'
[]
[v_pressure]
type = PINSFVMomentumPressure
variable = 'superficial_vel_y'
pressure = pressure
momentum_component = 'y'
[]
[v_friction]
type = PINSFVMomentumFriction
variable = 'superficial_vel_y'
Darcy_name = 'Darcy_coefficient'
Forchheimer_name = 'Forchheimer_coefficient'
momentum_component = 'y'
mu = 'mu'
speed = speed
[]
[gravity]
type = PINSFVMomentumGravity
variable = 'superficial_vel_y'
gravity = '0 -9.81 0'
momentum_component = 'y'
[]
[]
# ==============================================================================
# AUXVARIABLES AND AUXKERNELS
# ==============================================================================
[AuxVariables]
[T_fluid]
initial_condition = ${T_initial}
order = CONSTANT
family = MONOMIAL
[]
[vel_x]
initial_condition = ${fparse vel_x_initial}
order = CONSTANT
family = MONOMIAL
[]
[vel_y]
initial_condition = ${fparse vel_y_initial}
order = CONSTANT
family = MONOMIAL
[]
[porosity_out]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[vel_x]
type = FunctorAux
variable = vel_x
functor = vel_x_mat
[]
[vel_y]
type = FunctorAux
variable = vel_y
functor = vel_y_mat
[]
[porosity_out]
type = FunctorAux
variable = porosity_out
functor = porosity
[]
[]
# ==============================================================================
# FLUID PROPERTIES, MATERIALS AND USER OBJECTS
# ==============================================================================
[FluidProperties]
[fp]
type = IdealGasFluidProperties
k = 0.0
mu = ${mu}
gamma = 1.4
molar_mass = ${M}
[]
[]
[FunctorMaterials]
[enthalpy]
type = INSFVEnthalpyMaterial
temperature = 'T_fluid'
[]
[speed]
type = PINSFVSpeedFunctorMaterial
superficial_vel_x = 'superficial_vel_x'
superficial_vel_y = 'superficial_vel_y'
porosity = porosity
vel_x = vel_x_mat
vel_y = vel_y_mat
[]
[kappa]
type = FunctorKappaFluid
[]
[const_Fdrags_mat]
type = FunctorErgunDragCoefficients
porosity = porosity
[]
[fluidprops]
type = GeneralFunctorFluidProps
mu_rampdown = mu_func
porosity = porosity
characteristic_length = ${pebble_diameter}
T_fluid = 'T_fluid'
pressure = 'pressure'
speed = 'speed'
[]
[]
d = 0.05
[Functions]
[mu_func]
type = PiecewiseLinear
x = '1 3 5 10 15 20'
y = '1e5 1e4 1e3 1e2 1e1 1'
[]
[real_porosity_function]
type = ParsedFunction
expression = 'if (x < 0.6 - ${d}, ${porosity_left}, if (x > 0.6 + ${d}, ${porosity_right},
(x-(0.6-${d}))/(2*${d})*(${porosity_right}-${porosity_left}) + ${porosity_left}))'
[]
[porosity]
type = ParsedFunction
expression = 'if (x < 0.6 - ${d}, ${porosity_left}, if (x > 0.6 + ${d}, ${porosity_right},
(x-(0.6-${d}))/(2*${d})*(${porosity_right}-${porosity_left}) + ${porosity_left}))'
[]
[]
# ==============================================================================
# BOUNDARY CONDITIONS
# ==============================================================================
[FVBCs]
[outlet_p]
type = INSFVOutletPressureBC
variable = 'pressure'
function = ${p_initial}
boundary = 'right'
[]
## No or Free slip BC
[free-slip-wall-x]
type = INSFVNaturalFreeSlipBC
boundary = 'bottom wall_1 wall_2 left'
variable = superficial_vel_x
momentum_component = 'x'
[]
[free-slip-wall-y]
type = INSFVNaturalFreeSlipBC
boundary = 'bottom wall_1 wall_2 left'
variable = superficial_vel_y
momentum_component = 'y'
[]
## Symmetry
[symmetry-x]
type = PINSFVSymmetryVelocityBC
boundary = 'left'
variable = superficial_vel_x
u = superficial_vel_x
v = superficial_vel_y
mu = 'mu'
momentum_component = 'x'
[]
[symmetry-y]
type = PINSFVSymmetryVelocityBC
boundary = 'left'
variable = superficial_vel_y
u = superficial_vel_x
v = superficial_vel_y
mu = 'mu'
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'left'
variable = 'pressure'
[]
## inlet
[inlet_vel_x]
type = INSFVInletVelocityBC
variable = 'superficial_vel_x'
function = ${superficial_vel_x_initial}
boundary = 'top'
[]
[inlet_vel_y]
type = INSFVInletVelocityBC
variable = 'superficial_vel_y'
function = ${superficial_vel_y_initial}
boundary = 'top'
[]
[]
# ==============================================================================
# EXECUTION PARAMETERS
# ==============================================================================
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
petsc_options_value = 'asm lu NONZERO 200'
line_search = 'none'
# Problem time parameters
dtmin = 0.01
dtmax = 2000
end_time = 3000
# must be the same as the fluid
# Iterations parameters
l_max_its = 50
l_tol = 1e-8
nl_max_its = 25
# nl_rel_tol = 5e-7
nl_abs_tol = 2e-7
# Automatic scaling
automatic_scaling = true
verbose = true
[TimeStepper]
type = IterationAdaptiveDT
dt = 0.025
cutback_factor = 0.5
growth_factor = 2.0
[]
# Steady state detection.
steady_state_detection = true
steady_state_tolerance = 1e-7
steady_state_start_time = 400
[]
# ==============================================================================
# POSTPROCESSORS DEBUG AND OUTPUTS
# ==============================================================================
[Postprocessors]
[mass_flow_in]
type = VolumetricFlowRate
boundary = 'top'
vel_x = 'superficial_vel_x'
vel_y = 'superficial_vel_y'
advected_quantity = ${rho_initial}
execute_on = 'INITIAL TIMESTEP_END'
[]
[mass_flow_out]
type = VolumetricFlowRate
boundary = 'right'
vel_x = 'superficial_vel_x'
vel_y = 'superficial_vel_y'
advected_quantity = ${rho_initial}
execute_on = 'INITIAL TIMESTEP_END'
[]
[p_in]
type = SideAverageValue
variable = pressure
boundary = 'top'
[]
[dP]
type = LinearCombinationPostprocessor
pp_names = 'p_in'
pp_coefs = '1.0'
b = ${fparse -p_initial}
[]
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/materials/2d-transient.i)
l = 10
velocity_interp_method = 'rc'
advected_interp_method = 'average'
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_v = 0.001
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = 1
nx = 20
ny = 10
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
rho = 'rho'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = ${inlet_v}
[]
[v]
type = INSFVVelocityVariable
initial_condition = 1e-15
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[]
[AuxVariables]
[velocity_norm]
type = MooseVariableFVReal
[]
[power_density]
type = MooseVariableFVReal
initial_condition = 1e4
[]
[]
[FVKernels]
[mass_time]
type = WCNSFVMassTimeDerivative
variable = pressure
drho_dt = drho_dt
[]
[mass]
type = WCNSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho'
[]
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = u
drho_dt = drho_dt
rho = rho
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = u
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
rho = 'rho'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = u
mu = 'mu'
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = u
momentum_component = 'x'
pressure = pressure
[]
[v_time]
type = WCNSFVMomentumTimeDerivative
variable = v
drho_dt = drho_dt
rho = rho
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = v
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
rho = 'rho'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = v
mu = 'mu'
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = v
momentum_component = 'y'
pressure = pressure
[]
[temp_time]
type = WCNSFVEnergyTimeDerivative
variable = T
rho = rho
drho_dt = drho_dt
h = h
dh_dt = dh_dt
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[heat_source]
type = FVCoupledForce
variable = T
v = power_density
[]
[]
[FVBCs]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = u
boundary = 'top bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = v
boundary = 'top bottom'
function = 0
[]
# Inlet
[inlet_u]
type = INSFVInletVelocityBC
variable = u
boundary = 'left'
function = ${inlet_v}
[]
[inlet_v]
type = INSFVInletVelocityBC
variable = v
boundary = 'left'
function = 0
[]
[inlet_T]
type = FVDirichletBC
variable = T
boundary = 'left'
value = ${inlet_temp}
[]
[outlet_p]
type = INSFVOutletPressureBC
variable = pressure
boundary = 'right'
function = ${outlet_pressure}
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T'
rho = 'rho'
[]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T'
speed = 'velocity_norm'
# even though we provide rho from the parameters, we
# want to get rho from the fluid properties
force_define_density = true
# To initialize with a high viscosity
mu_rampdown = 'mu_rampdown'
# For porous flow
characteristic_length = 1
porosity = 1
[]
[]
[AuxKernels]
[speed]
type = VectorMagnitudeAux
variable = 'velocity_norm'
x = u
y = v
[]
[]
[Functions]
[mu_rampdown]
type = PiecewiseLinear
x = '1 2 3 4'
y = '1e3 1e2 1e1 1'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-3
optimal_iterations = 6
[]
end_time = 15
nl_abs_tol = 1e-12
nl_max_its = 50
line_search = 'none'
automatic_scaling = true
off_diagonals_in_auto_scaling = true
compute_scaling_once = false
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/natural_convection/natural_circulation_pipe.i)
# natural convection through a pipe
# Reference solution in "reference_pipe_natural_convection.py"
# Reference mdot: 0.0792 kg/s
# this input
# iy mdot
# 10 8.302364e-02
# 20 8.111192e-02
# 40 8.007924e-02
# 80 7.954403e-02
# 160 7.927201e-02
# Convergence to the analytical result is observed
height = 10.0
gravity = 9.81
p0 = 1e5
molar_mass = 29.0e-3
T0 = 328
Ru = 8.3145
Ri = '${fparse Ru / molar_mass}'
density = '${fparse p0 / (Ri * T0)}'
head = '${fparse height * density * gravity}'
k = 25.68e-3
gamma = 1.4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '0.1'
ix = '2'
dy = '${height}'
iy = '5'
[]
[]
[GlobalParams]
rhie_chow_user_object = pins_rhie_chow_interpolator
[]
[FluidProperties]
[air]
type = IdealGasFluidProperties
molar_mass = ${molar_mass}
k = ${k}
gamma = ${gamma}
[]
[]
[Modules]
[NavierStokesFV]
compressibility = 'weakly-compressible'
add_energy_equation = true
gravity = '0 -${gravity} 0'
density = rho
dynamic_viscosity = mu
specific_heat = cp
thermal_conductivity = k
initial_velocity = '0 1e-6 0'
initial_pressure = ${p0}
initial_temperature = ${T0}
inlet_boundaries = 'bottom'
momentum_inlet_types = 'fixed-pressure'
momentum_inlet_function = '${fparse p0 + head}'
energy_inlet_types = 'fixed-temperature'
energy_inlet_function = '${T0}'
energy_scaling = 1e-5
wall_boundaries = 'left right'
momentum_wall_types = 'slip slip'
energy_wall_types = 'heatflux heatflux'
energy_wall_function = '300 300'
outlet_boundaries = 'top'
momentum_outlet_types = 'fixed-pressure'
pressure_function = '${fparse p0}'
momentum_advection_interpolation = 'upwind'
mass_advection_interpolation = 'upwind'
porous_medium_treatment = true
porosity = porosity
energy_advection_interpolation = 'average'
[]
[]
[FVKernels]
[u_friction]
type = PINSFVMomentumFriction
variable = superficial_vel_x
Darcy_name = linear_friction_coeff
momentum_component = 'x'
standard_friction_formulation = false
rho = rho
[]
[v_friction]
type = PINSFVMomentumFriction
variable = superficial_vel_y
Darcy_name = linear_friction_coeff
momentum_component = 'y'
standard_friction_formulation = false
rho = rho
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-8
nl_abs_tol = 1e-6
end_time = 1e4
[TimeStepper]
type = IterationAdaptiveDT
dt = 0.1
growth_factor = 2
iteration_window = 2
optimal_iterations = 6
[]
[]
[Functions]
[mu_rampdown_fn]
type = PiecewiseLinear
x = '0 0.5 1 5 10 100 1000 2000'
y = '1000 1000 100 10 1 1 1 0'
[]
[]
[FunctorMaterials]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = air
pressure = pressure
T_fluid = T_fluid
speed = speed
force_define_density = true
mu_rampdown = 'mu_rampdown_fn'
characteristic_length = 1
porosity = porosity
[]
[scalar_props]
type = ADGenericFunctorMaterial
prop_names = 'porosity loss_coeff'
prop_values = '1 1.3'
[]
[linear_friction]
type = ADParsedFunctorMaterial
property_name = 'linear_friction'
expression = 'loss_coeff * rho'
functor_names = 'loss_coeff rho'
[]
[linear_friction_coeff]
type = ADGenericVectorFunctorMaterial
prop_names = 'linear_friction_coeff'
prop_values = 'linear_friction linear_friction linear_friction'
[]
[]
[AuxVariables]
[rho_var]
type = MooseVariableFVReal
[]
[cp_var]
type = MooseVariableFVReal
[]
[rho_cp_T_fluid_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[rho_var_aux]
type = FunctorAux
variable = rho_var
functor = rho
[]
[cp_var_aux]
type = FunctorAux
variable = cp_var
functor = cp
[]
[rho_cp_T_fluid_var_aux]
type = ParsedAux
variable = rho_cp_T_fluid_var
coupled_variables = 'rho_var cp_var T_fluid'
expression = 'rho_var * cp_var * T_fluid'
[]
[]
[Postprocessors]
[inlet_mfr]
type = VolumetricFlowRate
vel_x = superficial_vel_x
vel_y = superficial_vel_y
advected_quantity = rho
boundary = bottom
advected_interp_method = average
[]
[outlet_mfr]
type = VolumetricFlowRate
vel_x = superficial_vel_x
vel_y = superficial_vel_y
advected_quantity = rho
boundary = top
advected_interp_method = average
[]
[inlet_energy]
type = VolumetricFlowRate
vel_x = superficial_vel_x
vel_y = superficial_vel_y
advected_quantity = rho_cp_T_fluid_var
boundary = bottom
advected_interp_method = average
[]
[outlet_energy]
type = VolumetricFlowRate
vel_x = superficial_vel_x
vel_y = superficial_vel_y
advected_quantity = rho_cp_T_fluid_var
boundary = top
advected_interp_method = average
[]
[]
[Debug]
show_var_residual_norms = true
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/materials/functorfluidprops.i)
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_v = 4
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 2
ymin = 0
ymax = 1
nx = 5
ny = 5
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = ${inlet_v}
[]
[v]
type = INSFVVelocityVariable
initial_condition = 2
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[]
[FVKernels]
[u_time]
type = FVFunctorTimeKernel
variable = u
[]
[v_time]
type = FVFunctorTimeKernel
variable = v
[]
[p_time]
type = FVFunctorTimeKernel
variable = pressure
[]
[T_time]
type = FVFunctorTimeKernel
variable = T
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T'
speed = 'velocity_norm'
# For porous flow
characteristic_length = 2
porosity = 'porosity'
[]
[]
[AuxVariables]
[velocity_norm]
type = MooseVariableFVReal
[]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.4
[]
[rho_var]
type = MooseVariableFVReal
[]
[drho_dp_var]
type = MooseVariableFVReal
[]
[drho_dT_var]
type = MooseVariableFVReal
[]
[rho_dot_var]
type = MooseVariableFVReal
[]
[cp_var]
type = MooseVariableFVReal
[]
[dcp_dp_var]
type = MooseVariableFVReal
[]
[dcp_dT_var]
type = MooseVariableFVReal
[]
[cp_dot_var]
type = MooseVariableFVReal
[]
[cv_var]
type = MooseVariableFVReal
[]
[mu_var]
type = MooseVariableFVReal
[]
[dmu_dp_var]
type = MooseVariableFVReal
[]
[dmu_dT_var]
type = MooseVariableFVReal
[]
[k_var]
type = MooseVariableFVReal
[]
[dk_dp_var]
type = MooseVariableFVReal
[]
[dk_dT_var]
type = MooseVariableFVReal
[]
[Pr_var]
type = MooseVariableFVReal
[]
[dPr_dp_var]
type = MooseVariableFVReal
[]
[dPr_dT_var]
type = MooseVariableFVReal
[]
[Re_var]
type = MooseVariableFVReal
[]
[dRe_dp_var]
type = MooseVariableFVReal
[]
[dRe_dT_var]
type = MooseVariableFVReal
[]
[Re_h_var]
type = MooseVariableFVReal
[]
[Re_i_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[speed]
type = VectorMagnitudeAux
variable = 'velocity_norm'
x = u
y = v
[]
# To output the functor material properties
[rho_out]
type = FunctorAux
functor = 'rho'
variable = 'rho_var'
execute_on = 'timestep_begin'
[]
[drho_dp_out]
type = FunctorAux
functor = 'drho/dpressure'
variable = 'drho_dp_var'
execute_on = 'timestep_begin'
[]
[drho_dT_out]
type = FunctorAux
functor = 'drho/dT_fluid'
variable = 'drho_dT_var'
execute_on = 'timestep_begin'
[]
[drho_dt_out]
type = FunctorAux
functor = 'drho_dt'
variable = 'rho_dot_var'
execute_on = 'timestep_begin'
[]
[cp_out]
type = FunctorAux
functor = 'cp'
variable = 'cp_var'
execute_on = 'timestep_begin'
[]
[dcp_dp_out]
type = FunctorAux
functor = 'dcp/dpressure'
variable = 'dcp_dp_var'
execute_on = 'timestep_begin'
[]
[dcp_dT_out]
type = FunctorAux
functor = 'dcp/dT_fluid'
variable = 'dcp_dT_var'
execute_on = 'timestep_begin'
[]
[dcp_dt_out]
type = FunctorAux
functor = 'dcp_dt'
variable = 'cp_dot_var'
execute_on = 'timestep_begin'
[]
[cv_out]
type = FunctorAux
functor = 'cv'
variable = 'cv_var'
execute_on = 'timestep_begin'
[]
[mu_out]
type = FunctorAux
functor = 'mu'
variable = 'mu_var'
execute_on = 'timestep_begin'
[]
[dmu_dp_out]
type = FunctorAux
functor = 'dmu/dpressure'
variable = 'dmu_dp_var'
execute_on = 'timestep_begin'
[]
[dmu_dT_out]
type = FunctorAux
functor = 'dmu/dT_fluid'
variable = 'dmu_dT_var'
execute_on = 'timestep_begin'
[]
[k_out]
type = FunctorAux
functor = 'k'
variable = 'k_var'
execute_on = 'timestep_begin'
[]
[dk_dp_out]
type = FunctorAux
functor = 'dk/dpressure'
variable = 'dk_dp_var'
execute_on = 'timestep_begin'
[]
[dk_dT_out]
type = FunctorAux
functor = 'dk/dT_fluid'
variable = 'dk_dT_var'
execute_on = 'timestep_begin'
[]
[Pr_out]
type = FunctorAux
functor = 'Pr'
variable = 'Pr_var'
execute_on = 'timestep_begin'
[]
[dPr_dp_out]
type = FunctorAux
functor = 'dPr/dpressure'
variable = 'dPr_dp_var'
execute_on = 'timestep_begin'
[]
[dPr_dT_out]
type = FunctorAux
functor = 'dPr/dT_fluid'
variable = 'dPr_dT_var'
execute_on = 'timestep_begin'
[]
[Re_out]
type = FunctorAux
functor = 'Re'
variable = 'Re_var'
execute_on = 'timestep_begin'
[]
[dRe_dp_out]
type = FunctorAux
functor = 'dRe/dpressure'
variable = 'dRe_dp_var'
execute_on = 'timestep_begin'
[]
[dRe_dT_out]
type = FunctorAux
functor = 'dRe/dT_fluid'
variable = 'dRe_dT_var'
execute_on = 'timestep_begin'
[]
[Re_h_out]
type = FunctorAux
functor = 'Re_h'
variable = 'Re_h_var'
execute_on = 'timestep_begin'
[]
[Re_i_out]
type = FunctorAux
functor = 'Re_i'
variable = 'Re_i_var'
execute_on = 'timestep_begin'
[]
[]
[Executioner]
type = Transient
end_time = 0.1
dt = 0.1
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient.i)
# Fluid properties
mu = 'mu'
rho = 'rho'
cp = 'cp'
k = 'k'
# Solid properties
cp_s = 2
rho_s = 4
k_s = 1e-2
h_fs = 10
# Operating conditions
u_inlet = 1
T_inlet = 200
p_outlet = 10
top_side_temperature = 150
# Numerical scheme
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 1
nx = 20
ny = 5
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = superficial_vel_x
v = superficial_vel_y
pressure = pressure
porosity = porosity
[]
[]
[Variables]
[superficial_vel_x]
type = PINSFVSuperficialVelocityVariable
initial_condition = ${u_inlet}
[]
[superficial_vel_y]
type = PINSFVSuperficialVelocityVariable
initial_condition = 1e-6
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${p_outlet}
[]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = ${T_inlet}
[]
[T_solid]
type = MooseVariableFVReal
initial_condition = 100
[]
[]
[AuxVariables]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.5
[]
[velocity_norm]
type = MooseVariableFVReal
[]
[]
[FVKernels]
[mass_time]
type = PWCNSFVMassTimeDerivative
variable = pressure
porosity = 'porosity'
drho_dt = 'drho_dt'
[]
[mass]
type = PWCNSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = superficial_vel_x
rho = ${rho}
drho_dt = 'drho_dt'
momentum_component = 'x'
[]
[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_time]
type = WCNSFVMomentumTimeDerivative
variable = superficial_vel_y
rho = ${rho}
drho_dt = 'drho_dt'
momentum_component = 'y'
[]
[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_time]
type = PINSFVEnergyTimeDerivative
variable = T_fluid
cp = ${cp}
rho = ${rho}
drho_dt = 'drho_dt'
is_solid = false
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
variable = T_fluid
k = ${k}
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_time]
type = PINSFVEnergyTimeDerivative
variable = T_solid
cp = ${cp_s}
rho = ${rho_s}
is_solid = true
porosity = porosity
[]
[solid_energy_diffusion]
type = FVDiffusion
variable = T_solid
coeff = ${k_s}
[]
[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]
[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 = FVDirichletBC
variable = T_fluid
value = ${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 = ${top_side_temperature}
[]
[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 = ${p_outlet}
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T_fluid'
speed = 'velocity_norm'
# To initialize with a high viscosity
mu_rampdown = 'mu_rampdown'
# For porous flow
characteristic_length = 1
porosity = 'porosity'
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
rho = ${rho}
temperature = 'T_fluid'
[]
[constants]
type = ADGenericFunctorMaterial
prop_names = 'h_cv'
prop_values = '${h_fs}'
[]
[]
[Functions]
[mu_rampdown]
type = PiecewiseLinear
x = '1 2 3 4'
y = '1e3 1e2 1e1 1'
[]
[]
[AuxKernels]
[speed]
type = ParsedAux
variable = 'velocity_norm'
coupled_variables = 'superficial_vel_x superficial_vel_y porosity'
expression = 'sqrt(superficial_vel_x*superficial_vel_x + superficial_vel_y*superficial_vel_y) / '
'porosity'
[]
[]
[Executioner]
type = Transient
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
end_time = 3.0
[]
# 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
[]
(modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-action.i)
# Solid properties
cp_s = 2
rho_s = 4
k_s = 1e-2
h_fs = 10
# Operating conditions
u_inlet = 1
T_inlet = 200
p_outlet = 10
top_side_temperature = 150
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 1
nx = 20
ny = 5
[]
[]
[Variables]
[T_solid]
type = INSFVEnergyVariable
initial_condition = 100
[]
[]
[AuxVariables]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.5
[]
[velocity_norm]
type = MooseVariableFVReal
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[Modules]
[NavierStokesFV]
compressibility = 'weakly-compressible'
add_energy_equation = true
porous_medium_treatment = true
density = 'rho'
dynamic_viscosity = 'mu'
thermal_conductivity = 'k'
specific_heat = 'cp'
initial_velocity = '${u_inlet} 1e-6 0'
initial_pressure = '${p_outlet}'
initial_temperature = '${T_inlet}'
inlet_boundaries = 'left'
momentum_inlet_types = 'fixed-velocity'
momentum_inlet_function = '${u_inlet} 0'
energy_inlet_types = 'fixed-temperature'
energy_inlet_function = '${T_inlet}'
wall_boundaries = 'top bottom'
momentum_wall_types = 'noslip symmetry'
energy_wall_types = 'heatflux heatflux'
energy_wall_function = '0 0'
outlet_boundaries = 'right'
momentum_outlet_types = 'fixed-pressure'
pressure_function = '${p_outlet}'
ambient_convection_alpha = 'h_cv'
ambient_temperature = 'T_solid'
mass_advection_interpolation = 'average'
momentum_advection_interpolation = 'average'
energy_advection_interpolation = 'average'
[]
[]
[FVKernels]
[solid_energy_time]
type = PINSFVEnergyTimeDerivative
variable = T_solid
cp = ${cp_s}
rho = ${rho_s}
is_solid = true
porosity = 'porosity'
[]
[solid_energy_diffusion]
type = FVDiffusion
variable = T_solid
# this should use eps * k instead of k
coeff = ${k_s}
[]
[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]
[heated-side]
type = FVDirichletBC
boundary = 'top'
variable = 'T_solid'
value = ${top_side_temperature}
[]
[]
[FunctorMaterials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'h_cv'
prop_values = '${h_fs}'
[]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T_fluid'
speed = 'velocity_norm'
# To initialize with a high viscosity
mu_rampdown = 'mu_rampdown'
# For porous flow
characteristic_length = 1
porosity = 'porosity'
[]
[]
[Functions]
[mu_rampdown]
type = PiecewiseLinear
x = '1 2 3 4'
y = '1e3 1e2 1e1 1'
[]
[]
[AuxKernels]
[speed]
type = ParsedAux
variable = 'velocity_norm'
coupled_variables = 'superficial_vel_x superficial_vel_y porosity'
expression = 'sqrt(superficial_vel_x*superficial_vel_x + superficial_vel_y*superficial_vel_y) / porosity'
[]
[]
[Executioner]
type = Transient
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
end_time = 3.0
[]
# 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
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/pressure_driven_growth_transient.i)
###############################################################################
# Validation test based on Hibiki and Ishii experiment [1] reported in Figure 3
# [1] Hibiki, T., & Ishii, M. (2000). One-group interfacial area transport of bubbly flows in vertical round tubes.
# International Journal of Heat and Mass Transfer, 43(15), 2711-2726.
###############################################################################
mu = 1.0
rho = 1000.0
mu_d = 1.0
rho_d = 1.0
l = ${fparse 50.8/1000.0}
U = 0.491230114
dp = 0.001
inlet_phase_2 = 0.049
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
mass_exchange_coeff = 0.0
inlet_interface_area = ${fparse 6.0*inlet_phase_2/dp}
outlet_pressure = 1e6
[GlobalParams]
rhie_chow_user_object = 'rc'
density_interp_method = 'average'
mu_interp_method = 'average'
[]
[Problem]
identify_variable_groups_in_nl = false
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
coord_type = 'RZ'
rz_coord_axis = 'X'
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = '${fparse l * 60}'
ymin = 0
ymax = '${fparse l / 2}'
nx = 20
ny = 5
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_phase_2}
[]
[interface_area]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_interface_area}
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_x
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_y
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[phase_2_time]
type = FVFunctorTimeKernel
variable = phase_2
functor = phase_2
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_x'
v_slip = 'vel_y'
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1.0
[]
[phase_2_src]
type = NSFVMixturePhaseInterface
variable = phase_2
phase_coupled = phase_1
alpha = ${mass_exchange_coeff}
[]
[interface_area_time]
type = FVFunctorTimeKernel
variable = interface_area
functor = interface_area
[]
[interface_area_advection]
type = INSFVScalarFieldAdvection
variable = interface_area
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[interface_area_diffusion]
type = FVDiffusion
variable = interface_area
coeff = 0.1
[]
[interface_area_source_sink]
type = WCNSFV2PInterfaceAreaSourceSink
variable = interface_area
u = 'vel_x'
v = 'vel_y'
L = ${fparse l/2}
rho = 'rho_mixture'
rho_d = 'rho'
pressure = 'pressure'
k_c = '${fparse mass_exchange_coeff}'
fd = 'phase_2'
sigma = 1e-3
cutoff_fraction = 0.0
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '${outlet_pressure}'
[]
[inlet_phase_2]
type = FVDirichletBC
boundary = 'left'
variable = phase_2
value = ${inlet_phase_2}
[]
[inlet_interface_area]
type = FVDirichletBC
boundary = 'left'
variable = interface_area
value = ${inlet_interface_area}
[]
[symmetry-u]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_x
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'x'
[]
[symmetry-v]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_y
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[symmetry-phase-2]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = phase_2
[]
[symmetry-interface-area]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = interface_area
[]
[]
[AuxVariables]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[FunctorMaterials]
[bubble_properties]
type = GeneralFunctorFluidProps
fp = 'fp'
pressure = 'pressure'
T_fluid = 300.0
speed = 1.0
characteristic_length = 1.0
porosity = 1.0
output_properties = 'rho'
outputs = 'out'
[]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${rho} ${mu}'
phase_1_names = 'rho ${mu_d}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
nl_abs_tol = 1e-7
dt = 0.1
end_time = 1.0
nl_max_its = 10
line_search = 'none'
[]
[Debug]
show_var_residual_norms = true
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
[out]
type = Exodus
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${l} * ${U}'
pp_names = ''
[]
[rho_outlet]
type = SideAverageValue
boundary = 'right'
variable = 'rho_mixture_var'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/pressure_driven_growth.i)
###############################################################################
# Validation test based on Hibiki and Ishii experiment [1] reported in Figure 3
# [1] Hibiki, T., & Ishii, M. (2000). One-group interfacial area transport of bubbly flows in vertical round tubes.
# International Journal of Heat and Mass Transfer, 43(15), 2711-2726.
###############################################################################
mu = 1.0
rho = 1000.0
mu_d = 1.0
rho_d = 1.0
l = ${fparse 50.8/1000.0}
U = 0.491230114
dp = 0.001
inlet_phase_2 = 0.049
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
mass_exchange_coeff = 0.0
inlet_interface_area = ${fparse 6.0*inlet_phase_2/dp}
outlet_pressure = 1e5
[GlobalParams]
rhie_chow_user_object = 'rc'
density_interp_method = 'average'
mu_interp_method = 'average'
[]
[Problem]
identify_variable_groups_in_nl = false
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
coord_type = 'RZ'
rz_coord_axis = 'X'
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = '${fparse l * 60}'
ymin = 0
ymax = '${fparse l / 2}'
nx = 20
ny = 5
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_phase_2}
[]
[interface_area]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_interface_area}
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_x
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_y
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_x'
v_slip = 'vel_y'
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1.0
[]
[phase_2_src]
type = NSFVMixturePhaseInterface
variable = phase_2
phase_coupled = phase_1
alpha = ${mass_exchange_coeff}
[]
[interface_area_advection]
type = INSFVScalarFieldAdvection
variable = interface_area
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[interface_area_diffusion]
type = FVDiffusion
variable = interface_area
coeff = 0.1
[]
[interface_area_source_sink]
type = WCNSFV2PInterfaceAreaSourceSink
variable = interface_area
u = 'vel_x'
v = 'vel_y'
L = ${fparse l/2}
rho = 'rho_mixture'
rho_d = 'rho'
pressure = 'pressure'
k_c = '${fparse mass_exchange_coeff}'
fd = 'phase_2'
sigma = 1e-3
cutoff_fraction = 0.0
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '${outlet_pressure}'
[]
[inlet_phase_2]
type = FVDirichletBC
boundary = 'left'
variable = phase_2
value = ${inlet_phase_2}
[]
[inlet_interface_area]
type = FVDirichletBC
boundary = 'left'
variable = interface_area
value = ${inlet_interface_area}
[]
[symmetry-u]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_x
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'x'
[]
[symmetry-v]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_y
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[symmetry-phase-2]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = phase_2
[]
[symmetry-interface-area]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = interface_area
[]
[]
[AuxVariables]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[FunctorMaterials]
[bubble_properties]
type = GeneralFunctorFluidProps
fp = 'fp'
pressure = 'pressure'
T_fluid = 300.0
speed = 1.0
characteristic_length = 1.0
porosity = 1.0
output_properties = 'rho'
outputs = 'out'
[]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${rho} ${mu}'
phase_1_names = 'rho ${mu_d}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
line_search = 'none'
[]
[Debug]
show_var_residual_norms = true
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
[out]
type = Exodus
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${l} * ${U}'
pp_names = ''
[]
[rho_outlet]
type = SideAverageValue
boundary = 'right'
variable = 'rho_mixture_var'
[]
[]