- epsilonCoupled turbulent kinetic energy dissipation rate. 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:Coupled turbulent kinetic energy dissipation rate. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- muDynamic viscosity. 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:Dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- 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.
- tkeCoupled turbulent kinetic energy. 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:Coupled turbulent kinetic energy. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- uThe velocity in the x direction. 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:The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- variableThe name of the variable that this object applies to
C++ Type:AuxVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this object applies to
kEpsilonViscosityAux
This is the auxiliary kernel used to compute the dynamic turbulent viscosity
μt=ρCμkTe,
where:
ρ is the density,
Cμ=0.09 is a closure parameter,
k is the turbulent kinetic energy,
ϵ is the turbulent kinetic energy dissipation rate.
Te=max(ϵk,(ϵν)).
By setting parameter "bulk_wall_treatment" to true
, the kernel allows us to set the value of the cells on the boundaries specified in "walls" to the dynamic turbulent viscosity predicted from the law of the wall or non-equilibrium wall functions. See INSFVTurbulentViscosityWallFunction for more details about the near-wall implementation.
If the boundary conditions for the dynamic turbulent viscosity are already set via INSFVTurbulentViscosityWallFunction, there is no need to add bulk wall treatment, i.e., we should set "bulk_wall_treatment" to false
. This type of bulk wall treatment is mainly designed for porous media formulations with large computational cells.
Input Parameters
- C_muCoupled turbulent kinetic energy closure.
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Coupled turbulent kinetic energy closure.
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
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>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- bulk_wall_treatmentFalseActivate bulk wall treatment.
Default:False
C++ Type:bool
Controllable:No
Description:Activate bulk wall treatment.
- check_boundary_restrictedTrueWhether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
Default:True
C++ Type:bool
Controllable:No
Description:Whether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
- execute_onLINEAR TIMESTEP_ENDThe 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:LINEAR TIMESTEP_END
C++ Type:ExecFlagEnum
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, PRE_DISPLACE
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.
- mu_t_ratio_max100000Maximum allowable mu_t_ratio : mu/mu_t.
Default:100000
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Maximum allowable mu_t_ratio : mu/mu_t.
- scale_limiterstandardThe method used to limit the k-epsilon time scale.'none', 'standard'
Default:standard
C++ Type:MooseEnum
Options:none, standard
Controllable:No
Description:The method used to limit the k-epsilon time scale.'none', 'standard'
- vThe velocity in the y direction. 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:The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- wThe velocity in the z direction. 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:The velocity in the z direction. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- wall_treatmenteq_newtonThe method used for computing the wall functions.'eq_newton', 'eq_incremental', 'eq_linearized', 'neq'
Default:eq_newton
C++ Type:MooseEnum
Options:eq_newton, eq_incremental, eq_linearized, neq
Controllable:No
Description:The method used for computing the wall functions.'eq_newton', 'eq_incremental', 'eq_linearized', 'neq'
- wallsBoundaries that correspond to solid walls.
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Boundaries that correspond to solid walls.
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
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
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- newton_solveFalseWhether a Newton nonlinear solve is being used
Default:False
C++ Type:bool
Controllable:No
Description:Whether a Newton nonlinear solve is being used
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- 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.
- 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
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.
Material Property Retrieval Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-linear-wall.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/bfs/BFS_ERCOFTAC.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-inc-wall.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-std-wall.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-non-eq-bulk.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-energy.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-no-wall.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-capped.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-non-eq-wall.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/channel/channel_ERCOFTAC.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-energy-wall.i)
bulk_wall_treatment
Default:False
C++ Type:bool
Controllable:No
Description:Activate bulk wall treatment.
walls
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Boundaries that correspond to solid walls.
bulk_wall_treatment
Default:False
C++ Type:bool
Controllable:No
Description:Activate bulk wall treatment.
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-linear-wall.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# Linear wall function formulation (faster runs)
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'eq_linearized' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/bfs/BFS_ERCOFTAC.i)
##########################################################
# ERCOFTAC test case foe BFS
# Case Number: 031
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# Equilibrium + Newton wall treatement
# SIMPLE solve
##########################################################
Re = 5100
rho = 1.0
bulk_u = 1.0
H = 1.0
mu = '${fparse rho * bulk_u * H/ Re}'
advected_interp_method = 'upwind'
pressure_tag = "pressure_grad"
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Initial and Boundary Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * bulk_u)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / H}'
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'bottom wall-side top'
wall_treatment = 'eq_incremental' # Options: eq_newton, eq_incremental, eq_linearized, neq
[Mesh]
[gen]
type = CartesianMeshGenerator
dim = 2
dx = '${fparse 10.0*H} ${fparse 20.0*H}'
dy = '${H} ${fparse 5*H}'
ix = '8 16'
iy = '2 8'
subdomain_id = '
2 1
1 1
'
[]
[corner_walls]
type = SideSetsBetweenSubdomainsGenerator
input = gen
primary_block = '1'
paired_block = '2'
new_boundary = 'wall-side'
[]
[delete_bottom]
type = BlockDeletionGenerator
input = corner_walls
block = '2'
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = ${advected_interp_method}
velocity_interp_method = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = ${bulk_u}
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
initial_condition = 1e-8
solver_sys = pressure_system
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${bulk_u}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = 0
[]
[inlet_TKE]
type = INSFVInletIntensityTKEBC
boundary = 'left'
variable = TKE
u = vel_x
v = vel_y
intensity = ${intensity}
[]
[inlet_TKED]
type = INSFVMixingLengthTKEDBC
boundary = 'left'
variable = TKED
k = TKE
characteristic_length = '${fparse 2*H}'
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
functor = 0
[]
[walls-u]
type = FVDirichletBC
boundary = ${walls}
variable = vel_x
value = 0
[]
[walls-v]
type = FVDirichletBC
boundary = ${walls}
variable = vel_y
value = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = ${walls}
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.7
pressure_variable_relaxation = 0.3
turbulence_equation_relaxation = '0.3 0.3'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
[]
[Outputs]
exodus = true
[console]
type = Console
outlier_variable_norms = false
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-inc-wall.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# Incremental wall function formulation (similar to OpenFOAM)
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'eq_incremental' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-std-wall.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# Standard wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 1e-3
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
mu_interp_method = 'average'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
mu_interp_method = 'average'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
coeff_interp_method = 'average'
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
coeff_interp_method = 'average'
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
continue_on_max_its = true
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-non-eq-bulk.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# Standard wall functions with non-equilibrium bulk formaultion
# No wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-energy.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model with energy transport
# Standard wall functions without temperature wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
k = 0.01
cp = 10.0
Pr_t = 0.9
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system energy_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[T_fluid]
type = INSFVEnergyVariable
solver_sys = energy_system
initial_condition = 1.0
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T_fluid
[]
[temp_conduction]
type = FVDiffusion
coeff = ${k}
variable = T_fluid
[]
[temp_turb_conduction]
type = FVDiffusion
coeff = 'k_t'
variable = T_fluid
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[T_hot]
type = FVDirichletBC
variable = T_fluid
boundary = 'top'
value = 1
[]
[T_cold]
type = FVDirichletBC
variable = T_fluid
boundary = 'bottom'
value = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[FunctorMaterials]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T_fluid'
rho = ${rho}
cp = ${cp}
[]
[k_t]
type = ADParsedFunctorMaterial
expression = 'mu_t * cp / Pr_t'
functor_names = 'mu_t ${cp} ${Pr_t}'
functor_symbols = 'mu_t cp Pr_t'
property_name = 'k_t'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
energy_equation_relaxation = 0.9
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
energy_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
energy_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
energy_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
continue_on_max_its = true
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-no-wall.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# No wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = ''
wall_treatment = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[walls_TKED]
type = INSFVTKEDWallFunctionBC
boundary = 'left right top bottom'
variable = TKED
u = vel_x
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
[]
[walls_TKE]
type = FVDirichletBC
boundary = 'left right top bottom'
variable = TKE
value = ${k_init}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.7
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.9 0.9'
num_iterations = 1000
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-capped.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model with capped mixing length
# Standard wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
C_pl = 0.1
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C_pl = ${C_pl}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
C_pl = ${C_pl}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
continue_on_max_its = true
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-non-eq-wall.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# Standard wall functions with non-equilibrium wall formulation
# No wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment = 'neq' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/channel/channel_ERCOFTAC.i)
##########################################################
# ERCOFTAC test case foe turbulent channel flow
# Case Number: 032
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# Equilibrium + Newton wall treatement
# SIMPLE solve
##########################################################
H = 1 #halfwidth of the channel
L = 30
Re = 13700
rho = 1
bulk_u = 1
mu = '${fparse rho * bulk_u * 2 * H / Re}'
advected_interp_method = 'upwind'
pressure_tag = "pressure_grad"
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Initial and Boundary Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * bulk_u)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / H}'
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'top'
wall_treatment = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${L}
ymin = 0
ymax = ${H}
nx = 20
ny = 5
bias_y = 0.7
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = ${advected_interp_method}
velocity_interp_method = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = ${bulk_u}
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
initial_condition = 1e-8
solver_sys = pressure_system
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment}
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
function = '${bulk_u}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
function = 0
[]
[walls-u]
type = FVDirichletBC
boundary = 'top'
variable = vel_x
value = 0
[]
[walls-v]
type = FVDirichletBC
boundary = 'top'
variable = vel_y
value = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = 0
[]
[inlet_TKE]
type = INSFVInletIntensityTKEBC
boundary = 'left'
variable = TKE
u = vel_x
v = vel_y
intensity = ${intensity}
[]
[inlet_TKED]
type = INSFVMixingLengthTKEDBC
boundary = 'left'
variable = TKED
k = TKE
characteristic_length = '${fparse 2*H}'
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'top'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment}
[]
[sym-u]
type = INSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_x
u = vel_x
v = vel_y
mu = 'mu_t'
momentum_component = x
[]
[sym-v]
type = INSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_y
u = vel_x
v = vel_y
mu = 'mu_t'
momentum_component = y
[]
[symmetry_pressure]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[symmetry_TKE]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = TKE
[]
[symmetry_TKED]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = TKED
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[yplus]
type = MooseVariableFVReal
two_term_boundary_expansion = false
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[compute_y_plus]
type = RANSYPlusAux
variable = yplus
k = TKE
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
walls = ${walls}
wall_treatment = ${wall_treatment}
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.7
pressure_variable_relaxation = 0.3
turbulence_equation_relaxation = '0.25 0.25'
num_iterations = 1000
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
continue_on_max_its = true
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-energy-wall.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# Standard wall functions with temperature wall functions
# SIMPLE Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
k = 0.01
cp = 10.0
Pr_t = 0.9
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
bulk_wall_treatment = false
walls = 'left top right bottom'
wall_treatment_v = 'eq_newton' # Options: eq_newton, eq_incremental, eq_linearized, neq
wall_treatment_T = 'eq_linearized' # Options: eq_newton, eq_incremental, eq_linearized, neq
pressure_tag = "pressure_grad"
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 12
ny = 12
[]
[]
[Problem]
nl_sys_names = 'u_system v_system pressure_system energy_system TKE_system TKED_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolatorSegregated
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = u_system
two_term_boundary_expansion = false
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
solver_sys = v_system
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
solver_sys = pressure_system
initial_condition = 0.2
two_term_boundary_expansion = false
[]
[T_fluid]
type = INSFVEnergyVariable
solver_sys = energy_system
initial_condition = 1.0
two_term_boundary_expansion = false
[]
[TKE]
type = INSFVEnergyVariable
solver_sys = TKE_system
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
solver_sys = TKED_system
initial_condition = ${eps_init}
[]
[]
[FVKernels]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
extra_vector_tags = ${pressure_tag}
[]
[p_diffusion]
type = FVAnisotropicDiffusion
variable = pressure
coeff = "Ainv"
coeff_interp_method = 'average'
[]
[p_source]
type = FVDivergence
variable = pressure
vector_field = "HbyA"
force_boundary_execution = true
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T_fluid
[]
[temp_conduction]
type = FVDiffusion
coeff = ${k}
variable = T_fluid
[]
[temp_turb_conduction]
type = FVDiffusion
coeff = 'k_t'
variable = T_fluid
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
wall_treatment = ${wall_treatment_v}
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
k = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
wall_treatment = ${wall_treatment_v}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[T_hot]
type = INSFVTurbulentTemperatureWallFunction
variable = T_fluid
boundary = 'top'
T_w = 1
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
cp = ${cp}
kappa = ${k}
k = TKE
wall_treatment = ${wall_treatment_T}
[]
[T_cold]
type = INSFVTurbulentTemperatureWallFunction
variable = T_fluid
boundary = 'bottom'
T_w = 0
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
cp = ${cp}
kappa = ${k}
k = TKE
wall_treatment = ${wall_treatment_T}
[]
[walls_mu_t]
type = INSFVTurbulentViscosityWallFunction
boundary = 'left right top bottom'
variable = mu_t
u = vel_x
v = vel_y
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
k = TKE
wall_treatment = ${wall_treatment_v}
[]
[]
[AuxVariables]
[mu_t]
type = MooseVariableFVReal
initial_condition = '${fparse rho * C_mu * ${k_init}^2 / eps_init}'
two_term_boundary_expansion = false
[]
[k_t]
type = MooseVariableFVReal
initial_condition = 1.0
[]
[]
[AuxKernels]
[compute_mu_t]
type = kEpsilonViscosityAux
variable = mu_t
C_mu = ${C_mu}
k = TKE
epsilon = TKED
mu = ${mu}
rho = ${rho}
u = vel_x
v = vel_y
bulk_wall_treatment = ${bulk_wall_treatment}
walls = ${walls}
wall_treatment = ${wall_treatment_v}
execute_on = 'NONLINEAR'
[]
[compute_k_t]
type = TurbulentConductivityAux
variable = k_t
Pr_t = ${Pr_t}
cp = ${cp}
mu_t = 'mu_t'
execute_on = 'NONLINEAR'
[]
[]
[FunctorMaterials]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T_fluid'
rho = ${rho}
cp = ${cp}
[]
[]
[Executioner]
type = SIMPLENonlinearAssembly
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
turbulence_systems = 'TKED_system TKE_system'
pressure_gradient_tag = ${pressure_tag}
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.5
energy_equation_relaxation = 0.9
turbulence_equation_relaxation = '0.8 0.8'
num_iterations = 500
pressure_absolute_tolerance = 1e-12
momentum_absolute_tolerance = 1e-12
energy_absolute_tolerance = 1e-12
turbulence_absolute_tolerance = '1e-12 1e-12'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
momentum_l_abs_tol = 1e-14
energy_l_abs_tol = 1e-14
pressure_l_abs_tol = 1e-14
turbulence_l_abs_tol = 1e-14
momentum_l_max_its = 30
pressure_l_max_its = 30
momentum_l_tol = 0.0
energy_l_tol = 0.0
pressure_l_tol = 0.0
turbulence_l_tol = 0.0
print_fields = false
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '0.01 0.099 0.0'
continue_on_max_its = true
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = false
print_linear_residuals = true
[]