- mu_tildeSpalart-Allmaras turbulence viscosity variable
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:Spalart-Allmaras turbulence viscosity variable
- pressureThe pressure
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The pressure
- velocityThe velocity
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The velocity
- wall_distance_varWall distance aux variable name
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:Wall distance aux variable name
SATauMaterial
This is the material class used to compute the stabilization parameter tau_viscosity for the Spalart-Allmaras turbulent viscosity equation.
Overview
This material class computes strong residuals for the Spalart-Allmaras turbulence model equation and adjusts the stabilization parameter for use in pressure-stabilized (PSPG) and streamline-upwind (SUPG) kernels. SATauMaterial
includes all functionality from INSADTauMaterial for incompressible Navier-Stokes flow modeling with INSAD kernels from the MOOSE Navier-Stokes
physics module.
Example Input File Syntax
Input Parameters
- alpha1Multiplicative factor on the stabilization parameter tau.
Default:1
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Multiplicative factor on the stabilization parameter tau.
- 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
- computeTrueWhen false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:When false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
- 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.
- mu_namemuThe name of the dynamic viscosity
Default:mu
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The name of the dynamic viscosity
- 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.
- rho_namerhoThe name of the density
Default:rho
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The name of the density
- use_ft2_termFalseWhether to apply the f_t2 term in the equation
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether to apply the f_t2 term in the equation
- 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
- 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
Unit:(no unit assumed)
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
- 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
- (problems/2023-basic-turbulence-cases/pipe/pipe_flow.i)
- (tests/sa-model/pipe_flow.i)
- (problems/2023-basic-turbulence-cases/channel/channel_flow.i)
- (problems/2023-basic-turbulence-cases/backward-facing-step/bfs_flow_upstream.i)
- (tests/sa-model/channel_flow_with_precursors.i)
- (problems/2023-basic-turbulence-cases/backward-facing-step/bfs_flow_step.i)
- (tests/sa-model/channel_flow.i)
- (tests/sa-model/precursor_turbulent_diffusion.i)
(problems/2023-basic-turbulence-cases/pipe/pipe_flow.i)
Re = 4e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = 0.33333 # INS SUPG and PSPG stabilization parameter
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
[]
[Mesh]
coord_type = 'RZ'
[pipe]
type = CartesianMeshGenerator
dim = 2
dx = '0.3 0.1 0.05 0.034 0.016'
dy = '0.025 0.475 149.5'
ix = '20 10 10 17 40'
iy = '20 1 299'
[]
[corner_node]
type = ExtraNodesetGenerator
input = pipe
new_boundary = 'pinned_node'
coord = '0.5 150'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 5}'
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'right'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = right
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[Functions]
[yfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[mu_func]
type = ParsedFunction
expression = 'if(x=0.5, 0, ${fparse viscosity * 5})'
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = 0
function_y = yfunc
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'right'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom'
function_x = 0
function_y = yfunc
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'top'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'top'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'top'
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'bottom'
function = mu_func
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'right'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'top'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'top'
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
[]
[]
[VectorPostprocessors]
[outlet]
type = NodalValueSampler
variable = vely
boundary = 'top'
sort_by = x
execute_on = final
[]
[stress]
type = SideValueSampler
variable = turbulent_stress
boundary = 'top'
sort_by = x
execute_on = final
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = 0.1
scaling_group_variables = 'vel; p; mu_tilde'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-12
nl_max_its = 10
steady_state_detection = true
steady_state_tolerance = 1e-8
dtmin = 1e-2
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-2
cutback_factor = 0.9
growth_factor = 1.1
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
[]
[csv]
type = CSV
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(tests/sa-model/pipe_flow.i)
Re = 4e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = 0.33333 # INS SUPG and PSPG stabilization parameter
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
[]
[Mesh]
coord_type = 'RZ'
[pipe]
type = CartesianMeshGenerator
dim = 2
dx = '0.3 0.1 0.05 0.034 0.016'
dy = '0.025 0.475 19.5'
ix = '10 5 5 4 10'
iy = '20 1 39'
[]
[corner_node]
type = ExtraNodesetGenerator
input = pipe
new_boundary = 'pinned_node'
coord = '0.5 20'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 5}'
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'right'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = right
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[Functions]
[yfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[mu_func]
type = ParsedFunction
expression = 'if(x=0.5, 0, ${fparse viscosity * 5})'
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = 0
function_y = yfunc
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'right'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom'
function_x = 0
function_y = yfunc
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'top'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'top'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'top'
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'bottom'
function = mu_func
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'right'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'top'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'top'
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = 0.1
scaling_group_variables = 'vel; p; mu_tilde'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-14
nl_max_its = 10
steady_state_detection = true
steady_state_tolerance = 1e-10
dtmin = 1e-0
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-0
cutback_factor = 0.9
growth_factor = 1.1
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(problems/2023-basic-turbulence-cases/channel/channel_flow.i)
Re = 1.375e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = 0.33333 # INS SUPG and PSPG stabilization parameter
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
[]
[Mesh]
[channel]
type = CartesianMeshGenerator
dim = 2
dx = '0.05 0.2 139.75'
dy = '0.3 0.1 0.05 0.05'
ix = '20 1 559'
iy = '20 10 10 40'
[]
[corner_node]
type = ExtraNodesetGenerator
input = channel
new_boundary = 'pinned_node'
coord = '140 0.5'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 5}'
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'top'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = top
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[Functions]
[xfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[mu_func]
type = ParsedFunction
expression = 'if(y=0.5, 0, ${fparse viscosity * 5})'
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = xfunc
function_y = 0
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'top'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'left'
function_x = xfunc
function_y = 0
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'right'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'right'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'right'
[]
[symmetry]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom'
set_x_comp = false
function_y = 0
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'left'
function = mu_func
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'top'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'right'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'right'
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
[]
[]
[VectorPostprocessors]
[outlet]
type = NodalValueSampler
variable = velx
boundary = 'right'
sort_by = y
execute_on = final
[]
[prior]
type = NodalValueSampler
variable = velx
boundary = 'right'
sort_by = y
execute_on = final
[]
[stress]
type = SideValueSampler
variable = turbulent_stress
boundary = 'right'
sort_by = y
execute_on = final
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = 0.1
scaling_group_variables = 'vel; p; mu_tilde'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-12
nl_max_its = 10
steady_state_detection = true
steady_state_tolerance = 1e-8
dtmin = 1e-2
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-1
cutback_factor = 0.9
growth_factor = 1.1
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
[]
[csv]
type = CSV
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(problems/2023-basic-turbulence-cases/backward-facing-step/bfs_flow_upstream.i)
Re = 3.6e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = 0.33333 # INS SUPG and PSPG stabilization parameter
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
[]
[Mesh]
[block]
type = CartesianMeshGenerator
dim = 2
dx = '0.25 0.75 103'
ix = '40 1 103'
dy = '0.025 0.025 0.9 0.025 0.025 0.025 0.025 0.2 1.25 1 2.5 1.5 1.25 0.2 0.025 0.025'
iy = '16 8 36 5 10 20 5 10 25 20 5 30 25 10 5 20'
[]
[rename]
type = SubdomainBoundingBoxGenerator
input = block
bottom_left = '0 0 0'
top_right = '104 1 0'
block_id = 1
[]
[delete]
type = BlockDeletionGenerator
input = rename
block = 1
new_boundary = 'bottom'
[]
[corner_node]
type = ExtraNodesetGenerator
input = delete
new_boundary = 'pinned_node'
coord = '104 9'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 3}'
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'top bottom'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = 'top bottom'
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[Functions]
[xfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[mu_func]
type = ParsedFunction
expression = 'if(y=1, 0, if(y=9, 0, ${fparse viscosity * 3}))'
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = xfunc
function_y = 0
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'top bottom'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'left'
function_x = xfunc
function_y = 0
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'right'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'right'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'right'
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'left'
function = mu_func
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'top bottom'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'right'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'right'
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = 0.1
scaling_group_variables = 'vel; p; mu_tilde'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-14
nl_max_its = 10
l_max_its = 200
steady_state_detection = true
steady_state_tolerance = 1e-6
dtmin = 1e-2
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-1
cutback_factor = 0.9
growth_factor = 1.1
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
[]
[csv]
type = CSV
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(tests/sa-model/channel_flow_with_precursors.i)
Re = 1.375e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = .33333 # INS SUPG and PSPG stabilization parameter
diff = '${fparse viscosity / density / 0.85}'
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
use_exp_form = false
[]
[Mesh]
[channel]
type = CartesianMeshGenerator
dim = 2
dx = '0.05 0.2 19.75'
dy = '0.3 0.1 0.05 0.05'
ix = '20 1 79'
iy = '10 5 5 10'
[]
[corner_node]
type = ExtraNodesetGenerator
input = channel
new_boundary = 'pinned_node'
coord = '20 0.5'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 5}'
[]
[prec]
family = MONOMIAL
order = CONSTANT
initial_condition = 300
scaling = 1e-3
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[prec_source]
type = BodyForce
variable = prec
[]
[prec_time]
type = ScalarTransportTimeDerivative
variable = prec
[]
[]
[DGKernels]
[diffusion]
type = DGDiffusion
variable = prec
epsilon = -1
sigma = 1e-3
diff = '${diff}'
[]
[turbulent_diffusion]
type = DGTurbulentDiffusion
variable = prec
mu_tilde = mu_tilde
epsilon = -1
sigma = 1e-3
[]
[advection]
type = DGCoupledAdvection
variable = prec
uvel = velx
vvel = vely
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'top'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = top
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[Functions]
[xfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[mu_func]
type = ParsedFunction
expression = 'if(y=0.5, 0, ${fparse viscosity * 5})'
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = xfunc
function_y = 0
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'top'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'left'
function_x = xfunc
function_y = 0
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'right'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'right'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'right'
[]
[symmetry]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom'
set_x_comp = false
function_y = 0
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'left'
function = mu_func
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'top'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'right'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'right'
[]
[prec_inlet]
type = InflowBC
variable = prec
boundary = 'left'
uu = '${inlet}'
inlet_conc = 300
[]
[prec_outlet]
type = CoupledScalarAdvectionNoBCBC
variable = prec
boundary = 'right'
u = velx
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
use_ft2_term = true
[]
[]
[Executioner]
type = Transient
end_time = 30
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = .1
ignore_variables_for_autoscaling = 'prec'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-13
nl_max_its = 10
dtmin = 1
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1
cutback_factor = 0.9
growth_factor = 1.2
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(problems/2023-basic-turbulence-cases/backward-facing-step/bfs_flow_step.i)
Re = 3.6e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = 0.33333 # INS SUPG and PSPG stabilization parameter
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
[]
[Mesh]
[block]
type = CartesianMeshGenerator
dim = 2
dx = '4 0.8 0.4 1.6 0.4 6.8 4 2 2 2.5 3 3.5 25'
ix = '64 16 10 50 10 136 64 25 20 20 20 20 125'
dy = '0.025 0.025 0.9 0.025 0.025 0.025 0.025 0.2 1.25 1 2.5 1.5 1.25 0.2 0.025 0.025'
iy = '16 8 36 5 10 20 5 10 25 20 5 30 25 10 5 20'
[]
[rename]
type = SubdomainBoundingBoxGenerator
input = block
bottom_left = '0 0 0'
top_right = '6 1 0'
block_id = 1
[]
[delete]
type = BlockDeletionGenerator
input = rename
block = 1
new_boundary = 'step'
[]
[bottom]
type = SideSetsFromNormalsGenerator
input = delete
normals = '0 -1 0'
new_boundary = 'bottom'
fixed_normal = true
replace = true
[]
[transform]
type = TransformGenerator
input = bottom
transform = TRANSLATE
vector_value = '104 0 0'
[]
[corner_node]
type = ExtraNodesetGenerator
input = transform
new_boundary = 'pinned_node'
coord = '160 0 0'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 3}'
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'top bottom step'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = 'top bottom'
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[UserObjects]
[inlet]
type = SolutionUserObject
mesh = bfs_flow_upstream_exodus.e
system_variables = 'velx vely mu_tilde'
timestep = LATEST
execute_on = initial
[]
[]
[Functions]
[xfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[velx_bc]
type = SolutionFunction
solution = inlet
from_variable = velx
[]
[vely_bc]
type = SolutionFunction
solution = inlet
from_variable = vely
[]
[mu_tilde_bc]
type = SolutionFunction
solution = inlet
from_variable = mu_tilde
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = xfunc
function_y = 0
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'top bottom step'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'left'
function_x = velx_bc
function_y = vely_bc
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'right'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'right'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'right'
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'left'
function = mu_tilde_bc
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'top bottom step'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'right'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'right'
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
[]
[]
[VectorPostprocessors]
[vel1]
type = LineValueSampler
variable = 'velx mu_turbulence'
start_point = '106 1 0'
end_point = '106 9 0'
num_points = 25601
sort_by = y
execute_on = final
[]
[vel2]
type = LineValueSampler
variable = 'velx mu_turbulence'
start_point = '111 0 0'
end_point = '111 9 0'
num_points = 28801
sort_by = y
execute_on = final
[]
[vel3]
type = LineValueSampler
variable = 'velx mu_turbulence'
start_point = '114 0 0'
end_point = '114 9 0'
num_points = 28801
sort_by = y
execute_on = final
[]
[vel4]
type = LineValueSampler
variable = 'velx mu_turbulence'
start_point = '116 0 0'
end_point = '116 9 0'
num_points = 28801
sort_by = y
execute_on = final
[]
[vel5]
type = LineValueSampler
variable = 'velx mu_turbulence'
start_point = '120 0 0'
end_point = '120 9 0'
num_points = 28801
sort_by = y
execute_on = final
[]
[upstream]
type = LineValueSampler
variable = 'velx p'
start_point = '104 1.00125 0'
end_point = '110 1.00125 0'
num_points = 7201
sort_by = x
execute_on = final
[]
[downstream]
type = LineValueSampler
variable = 'velx p'
start_point = '110 0.0015625 0'
end_point = '160 0.0015625 0'
num_points = 60001
sort_by = x
execute_on = final
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = .1
scaling_group_variables = 'vel; p; mu_tilde'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-12
nl_max_its = 10
l_max_its = 200
steady_state_detection = true
steady_state_tolerance = 1e-6
dtmin = 1e-2
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-1
cutback_factor = 0.9
growth_factor = 1.1
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
[]
[csv]
type = CSV
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(tests/sa-model/channel_flow.i)
Re = 1.375e4
density = 1 # kg cm-3
D = 1
inlet = 1
viscosity = '${fparse density * inlet * D / Re}' # dynamic viscosity, mu = nu * rho, kg cm-1 s-1
alpha = 0.33333 # INS SUPG and PSPG stabilization parameter
[GlobalParams]
integrate_p_by_parts = false
viscous_form = 'traction'
pressure = p
velocity = vel
[]
[Mesh]
[channel]
type = CartesianMeshGenerator
dim = 2
dx = '0.05 0.2 19.75'
dy = '0.3 0.1 0.05 0.05'
ix = '20 1 79'
iy = '10 5 5 10'
[]
[corner_node]
type = ExtraNodesetGenerator
input = channel
new_boundary = 'pinned_node'
coord = '20 0.5'
[]
[]
[Problem]
type = FEProblem
[]
[Variables]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[mu_tilde]
initial_condition = '${fparse viscosity * 5}'
[]
[]
[AuxVariables]
[mu_turbulence]
[]
[wall_dist]
[]
[velx]
[]
[vely]
[]
[wall_shear_stress]
family = MONOMIAL
order = CONSTANT
[]
[turbulent_stress]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = vel
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = vel
[]
[momentum_turbulent_viscous]
type = INSADMomentumTurbulentViscous
variable = vel
mu_tilde = mu_tilde
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = vel
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = vel
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = vel
[]
[mu_tilde_time]
type = SATimeDerivative
variable = mu_tilde
[]
[mu_tilde_space]
type = SATurbulentViscosity
variable = mu_tilde
[]
[mu_tilde_supg]
type = SATurbulentViscositySUPG
variable = mu_tilde
[]
[]
[AuxKernels]
[mu_space]
type = SATurbulentViscosityAux
variable = mu_turbulence
mu_tilde = mu_tilde
[]
[wall_distance]
type = WallDistanceAux
variable = wall_dist
walls = 'top'
execute_on = initial
[]
[velx]
type = VectorVariableComponentAux
variable = velx
vector_variable = vel
component = x
[]
[vely]
type = VectorVariableComponentAux
variable = vely
vector_variable = vel
component = y
[]
[wall_shear_stress]
type = WallShearStressAux
variable = wall_shear_stress
boundary = top
[]
[turbulent_stress]
type = TurbulentStressAux
variable = turbulent_stress
mu_tilde = mu_tilde
[]
[]
[Functions]
[xfunc]
type = ParsedFunction
expression = '${fparse inlet}'
[]
[mu_func]
type = ParsedFunction
expression = 'if(y=0.5, 0, ${fparse viscosity * 5})'
[]
[]
[ICs]
[velocity]
type = VectorFunctionIC
function_x = xfunc
function_y = 0
variable = vel
[]
[]
[BCs]
[no_slip]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'top'
function_x = 0
function_y = 0
[]
[inlet]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'left'
function_x = xfunc
function_y = 0
[]
[outlet]
type = INSADMomentumTurbulentNoBCBC
variable = vel
boundary = 'right'
mu_tilde = mu_tilde
[]
[outlet_supg]
type = INSADMomentumSUPGBC
variable = vel
boundary = 'right'
[]
[outlet_pspg]
type = INSADMassPSPGBC
variable = p
boundary = 'right'
[]
[symmetry]
type = ADVectorFunctionDirichletBC
variable = vel
boundary = 'bottom'
set_x_comp = false
function_y = 0
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[mu_inlet]
type = FunctionDirichletBC
variable = mu_tilde
boundary = 'left'
function = mu_func
[]
[mu_wall]
type = DirichletBC
variable = mu_tilde
boundary = 'top'
value = 0
[]
[mu_outlet]
type = SATurbulentViscosityNoBCBC
variable = mu_tilde
boundary = 'right'
[]
[mu_outlet_supg]
type = SATurbulentViscositySUPGBC
variable = mu_tilde
boundary = 'right'
[]
[]
[Materials]
[ad_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${density} ${viscosity}'
[]
[sa_mat]
type = SATauMaterial
alpha = ${alpha}
mu_tilde = mu_tilde
wall_distance_var = wall_dist
use_ft2_term = true
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = none
automatic_scaling = true
compute_scaling_once = false
resid_vs_jac_scaling_param = 0.1
scaling_group_variables = 'vel; p; mu_tilde'
off_diagonals_in_auto_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-14
nl_max_its = 10
steady_state_detection = true
steady_state_tolerance = 1e-10
dtmin = 1e-0
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-0
cutback_factor = 0.9
growth_factor = 1.1
optimal_iterations = 6
iteration_window = 1
linear_iteration_ratio = 1000
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Outputs]
[exodus]
type = Exodus
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]
(tests/sa-model/precursor_turbulent_diffusion.i)
[GlobalParams]
pressure = p
velocity = vel
[]
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[]
[]
[Variables]
[prec]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxVariables]
[mu_tilde]
initial_condition = 1
[]
[vel]
family = LAGRANGE_VEC
order = FIRST
[]
[p]
[]
[wall_dist]
[]
[]
[DGKernels]
[turbulent_diffusion]
type = DGTurbulentDiffusion
variable = prec
mu_tilde = mu_tilde
[]
[]
[BCs]
[left]
type = PenaltyDirichletBC
variable = prec
boundary = 'left'
value = 1
penalty = 1e5
[]
[right]
type = PenaltyDirichletBC
variable = prec
boundary = 'right'
value = 0
penalty = 1e6
[]
[]
[Materials]
[const_mat]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[]
[sa_mat]
type = SATauMaterial
alpha = 0.33333
mu_tilde = mu_tilde
wall_distance_var = wall_dist
use_ft2_term = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
[exodus]
type = Exodus
execute_on = final
[]
[]
[Debug]
show_var_residual_norms = true
[]