- variableThe name of the variable that this Kernel operates on
C++ Type:NonlinearVariableName
Description:The name of the variable that this Kernel operates on
This class computes the time derivative for the incompressible Navier-Stokes momentum equation.
Input Parameters
- blockThe list of block ids (SubdomainID) that this object will be applied
C++ Type:std::vector
Options:
Description:The list of block ids (SubdomainID) that this object will be applied
- displacementsThe displacements
C++ Type:std::vector
Options:
Description:The displacements
- lumpingFalseTrue for mass matrix lumping, false otherwise
Default:False
C++ Type:bool
Options:
Description:True for mass matrix lumping, false otherwise
- rho_namerhodensity name
Default:rho
C++ Type:MaterialPropertyName
Options:
Description:density name
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector
Options:
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
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
Options:
Description:Determines whether this object is calculated using an implicit or explicit form
- save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector
Options:
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Options:
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
Options:
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
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector
Options:
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector
Options:
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystem timeThe tag for the matrices this Kernel should fill
Default:system time
C++ Type:MultiMooseEnum
Options:nontime system time
Description:The tag for the matrices this Kernel should fill
- vector_tagstimeThe tag for the vectors this Kernel should fill
Default:time
C++ Type:MultiMooseEnum
Options:nontime time
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_high_reynolds.i
- modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_stab_jac_test.i
- modules/navier_stokes/test/tests/ins/mms/supg/supg_pspg_adv_dominated_mms.i
- modules/navier_stokes/test/tests/ins/mms/supg/supg_adv_dominated_mms.i
- modules/navier_stokes/test/tests/ins/jacobian_test/jacobian_test.i
- modules/navier_stokes/test/tests/ins/stagnation/stagnation.i
- modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_no_parts.i
- modules/navier_stokes/test/tests/ins/jeffery_hamel/wedge_dirichlet.i
- modules/navier_stokes/test/tests/ins/jeffery_hamel/wedge_natural.i
- modules/navier_stokes/test/tests/ins/lid_driven/lid_driven.i
- modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_by_parts.i
modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_high_reynolds.i
[GlobalParams]
gravity = '0 0 0'
laplace = true
transient_term = false
supg = true
pspg = true
family = LAGRANGE
order = FIRST
[]
[Mesh]
file = 'cone_linear_alltri.e'
[]
[Problem]
coord_type = RZ
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
# type = Transient
# dt = 0.005
# dtmin = 0.005
# num_steps = 5
# l_max_its = 100
# Block Jacobi works well for this problem, as does "-pc_type asm
# -pc_asm_overlap 2", but an overlap of 1 does not work for some
# reason?
# petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
# petsc_options_value = 'bjacobi ilu 4'
# Note: The Steady executioner can be used for this problem, if you
# drop the INSMomentumTimeDerivative kernels and use the following
# direct solver options.
type = Steady
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
nl_max_its = 20
[]
[Outputs]
csv = true
console = true
[./out]
type = Exodus
[../]
[]
[Variables]
[./vel_x]
# Velocity in radial (r) direction
[../]
[./vel_y]
# Velocity in axial (z) direction
[../]
[./p]
order = FIRST
[../]
[]
[BCs]
[./u_in]
type = DirichletBC
boundary = bottom
variable = vel_x
value = 0
[../]
[./v_in]
type = FunctionDirichletBC
boundary = bottom
variable = vel_y
function = 'inlet_func'
[../]
[./u_axis_and_walls]
type = DirichletBC
boundary = 'left right'
variable = vel_x
value = 0
[../]
[./v_no_slip]
type = DirichletBC
boundary = 'right'
variable = vel_y
value = 0
[../]
[]
[Kernels]
# [./x_momentum_time]
# type = INSMomentumTimeDerivative
# variable = vel_x
# [../]
# [./y_momentum_time]
# type = INSMomentumTimeDerivative
# variable = vel_y
# [../]
[./mass]
type = INSMassRZ
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceFormRZ
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceFormRZ
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 1
prop_names = 'rho mu'
prop_values = '1 1e-3'
[../]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
value = '-4 * x^2 + 1'
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]
modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_stab_jac_test.i
[GlobalParams]
gravity = '0 0 0'
laplace = true
transient_term = true
supg = true
pspg = true
family = LAGRANGE
order = SECOND
[]
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 1
xmin = 0
xmax = 1.1
ymin = -1.1
ymax = 1.1
elem_type = QUAD9
[]
[Problem]
coord_type = RZ
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 1.1
# petsc_options = '-snes_test_display'
petsc_options_iname = '-snes_type'
petsc_options_value = 'test'
[]
[Variables]
[./vel_x]
# Velocity in radial (r) direction
[../]
[./vel_y]
# Velocity in axial (z) direction
[../]
[./p]
order = FIRST
[../]
[]
[Kernels]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./mass]
type = INSMassRZ
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceFormRZ
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceFormRZ
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1.1 1.1'
[../]
[]
[ICs]
[./vel_x]
type = RandomIC
variable = vel_x
min = 0.1
max = 0.9
[../]
[./vel_y]
type = RandomIC
variable = vel_y
min = 0.1
max = 0.9
[../]
[./p]
type = RandomIC
variable = p
min = 0.1
max = 0.9
[../]
[]
[Outputs]
dofmap = true
[]
modules/navier_stokes/test/tests/ins/mms/supg/supg_pspg_adv_dominated_mms.i
mu=1.5e-4
rho=2.5
[GlobalParams]
gravity = '0 0 0'
supg = true
pspg = true
convective_term = true
integrate_p_by_parts = false
transient_term = true
laplace = true
u = vel_x
v = vel_y
p = p
alpha = 1e0
order = FIRST
family = LAGRANGE
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
elem_type = QUAD9
nx = 4
ny = 4
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./vel_x]
[../]
[./vel_y]
[../]
[./p]
order = FIRST
[../]
[]
[Kernels]
# mass
[./mass]
type = INSMass
variable = p
x_vel_forcing_func = vel_x_source_func
y_vel_forcing_func = vel_y_source_func
[../]
[./x_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
# x-momentum, space
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
component = 0
forcing_func = vel_x_source_func
[../]
# y-momentum, space
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
component = 1
forcing_func = vel_y_source_func
[../]
[./p_source]
type = BodyForce
function = p_source_func
variable = p
[../]
[]
[BCs]
[./vel_x]
type = FunctionDirichletBC
boundary = 'left right top bottom'
function = vel_x_func
variable = vel_x
[../]
[./vel_y]
type = FunctionDirichletBC
boundary = 'left right top bottom'
function = vel_y_func
variable = vel_y
[../]
[./p]
type = FunctionDirichletBC
boundary = 'left right top bottom'
function = p_func
variable = p
[../]
[]
[Functions]
[./vel_x_source_func]
type = ParsedFunction
value = '-${mu}*(-0.028*pi^2*x^2*sin(0.2*pi*x*y) - 0.028*pi^2*y^2*sin(0.2*pi*x*y) - 0.1*pi^2*sin(0.5*pi*x) - 0.4*pi^2*sin(pi*y)) + ${rho}*(0.14*pi*x*cos(0.2*pi*x*y) + 0.4*pi*cos(pi*y))*(0.6*sin(0.8*pi*x) + 0.3*sin(0.3*pi*y) + 0.2*sin(0.3*pi*x*y) + 0.3) + ${rho}*(0.14*pi*y*cos(0.2*pi*x*y) + 0.2*pi*cos(0.5*pi*x))*(0.4*sin(0.5*pi*x) + 0.4*sin(pi*y) + 0.7*sin(0.2*pi*x*y) + 0.5) + 0.1*pi*y*cos(0.2*pi*x*y) + 0.25*pi*cos(0.5*pi*x)'
[../]
[./vel_y_source_func]
type = ParsedFunction
value = '-${mu}*(-0.018*pi^2*x^2*sin(0.3*pi*x*y) - 0.018*pi^2*y^2*sin(0.3*pi*x*y) - 0.384*pi^2*sin(0.8*pi*x) - 0.027*pi^2*sin(0.3*pi*y)) + ${rho}*(0.06*pi*x*cos(0.3*pi*x*y) + 0.09*pi*cos(0.3*pi*y))*(0.6*sin(0.8*pi*x) + 0.3*sin(0.3*pi*y) + 0.2*sin(0.3*pi*x*y) + 0.3) + ${rho}*(0.06*pi*y*cos(0.3*pi*x*y) + 0.48*pi*cos(0.8*pi*x))*(0.4*sin(0.5*pi*x) + 0.4*sin(pi*y) + 0.7*sin(0.2*pi*x*y) + 0.5) + 0.1*pi*x*cos(0.2*pi*x*y) + 0.3*pi*cos(0.3*pi*y)'
[../]
[./p_source_func]
type = ParsedFunction
value = '-0.06*pi*x*cos(0.3*pi*x*y) - 0.14*pi*y*cos(0.2*pi*x*y) - 0.2*pi*cos(0.5*pi*x) - 0.09*pi*cos(0.3*pi*y)'
[../]
[./vel_x_func]
type = ParsedFunction
value = '0.4*sin(0.5*pi*x) + 0.4*sin(pi*y) + 0.7*sin(0.2*pi*x*y) + 0.5'
[../]
[./vel_y_func]
type = ParsedFunction
value = '0.6*sin(0.8*pi*x) + 0.3*sin(0.3*pi*y) + 0.2*sin(0.3*pi*x*y) + 0.3'
[../]
[./p_func]
type = ParsedFunction
value = '0.5*sin(0.5*pi*x) + 1.0*sin(0.3*pi*y) + 0.5*sin(0.2*pi*x*y) + 0.5'
[../]
[./vxx_func]
type = ParsedFunction
value = '0.14*pi*y*cos(0.2*pi*x*y) + 0.2*pi*cos(0.5*pi*x)'
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_view'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu NONZERO superlu_dist'
line_search = 'none'
nl_rel_tol = 1e-8
nl_abs_tol = 1e-12
nl_max_its = 10
l_tol = 1e-6
l_max_its = 10
# To run to steady-state, set num-steps to some large number (1000000 for example)
type = Transient
num_steps = 10
steady_state_detection = true
steady_state_tolerance = 1e-10
[./TimeStepper]
dt = .1
type = IterationAdaptiveDT
cutback_factor = 0.4
growth_factor = 1.2
optimal_iterations = 20
[../]
[]
[Outputs]
execute_on = 'final'
[./exodus]
type = Exodus
[../]
[./csv]
type = CSV
[../]
[]
[Postprocessors]
[./L2vel_x]
type = ElementL2Error
variable = vel_x
function = vel_x_func
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2vel_y]
variable = vel_y
function = vel_y_func
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2p]
variable = p
function = p_func
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2vxx]
variable = vxx
function = vxx_func
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]
[AuxVariables]
[./vxx]
family = MONOMIAL
order = FIRST
[../]
[]
[AuxKernels]
[./vxx]
type = VariableGradientComponent
component = x
variable = vxx
gradient_variable = vel_x
[../]
[]
modules/navier_stokes/test/tests/ins/mms/supg/supg_adv_dominated_mms.i
mu=1.5e-2
rho=2.5
[GlobalParams]
gravity = '0 0 0'
supg = true
convective_term = true
integrate_p_by_parts = false
transient_term = true
laplace = true
u = vel_x
v = vel_y
p = p
alpha = 1e0
order = SECOND
family = LAGRANGE
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
elem_type = QUAD9
nx = 4
ny = 4
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./vel_x]
[../]
[./vel_y]
[../]
[./p]
order = FIRST
[../]
[]
[Kernels]
# mass
[./mass]
type = INSMass
variable = p
[../]
[./x_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
# x-momentum, space
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
component = 0
forcing_func = vel_x_source_func
[../]
# y-momentum, space
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
component = 1
forcing_func = vel_y_source_func
[../]
[./p_source]
type = BodyForce
function = p_source_func
variable = p
[../]
[]
[BCs]
[./vel_x]
type = FunctionDirichletBC
boundary = 'left right top bottom'
function = vel_x_func
variable = vel_x
[../]
[./vel_y]
type = FunctionDirichletBC
boundary = 'left right top bottom'
function = vel_y_func
variable = vel_y
[../]
[./p]
type = FunctionDirichletBC
boundary = 'left right top bottom'
function = p_func
variable = p
[../]
[]
[Functions]
[./vel_x_source_func]
type = ParsedFunction
value = '-${mu}*(-0.028*pi^2*x^2*sin(0.2*pi*x*y) - 0.028*pi^2*y^2*sin(0.2*pi*x*y) - 0.1*pi^2*sin(0.5*pi*x) - 0.4*pi^2*sin(pi*y)) + ${rho}*(0.14*pi*x*cos(0.2*pi*x*y) + 0.4*pi*cos(pi*y))*(0.6*sin(0.8*pi*x) + 0.3*sin(0.3*pi*y) + 0.2*sin(0.3*pi*x*y) + 0.3) + ${rho}*(0.14*pi*y*cos(0.2*pi*x*y) + 0.2*pi*cos(0.5*pi*x))*(0.4*sin(0.5*pi*x) + 0.4*sin(pi*y) + 0.7*sin(0.2*pi*x*y) + 0.5) + 0.1*pi*y*cos(0.2*pi*x*y) + 0.25*pi*cos(0.5*pi*x)'
[../]
[./vel_y_source_func]
type = ParsedFunction
value = '-${mu}*(-0.018*pi^2*x^2*sin(0.3*pi*x*y) - 0.018*pi^2*y^2*sin(0.3*pi*x*y) - 0.384*pi^2*sin(0.8*pi*x) - 0.027*pi^2*sin(0.3*pi*y)) + ${rho}*(0.06*pi*x*cos(0.3*pi*x*y) + 0.09*pi*cos(0.3*pi*y))*(0.6*sin(0.8*pi*x) + 0.3*sin(0.3*pi*y) + 0.2*sin(0.3*pi*x*y) + 0.3) + ${rho}*(0.06*pi*y*cos(0.3*pi*x*y) + 0.48*pi*cos(0.8*pi*x))*(0.4*sin(0.5*pi*x) + 0.4*sin(pi*y) + 0.7*sin(0.2*pi*x*y) + 0.5) + 0.1*pi*x*cos(0.2*pi*x*y) + 0.3*pi*cos(0.3*pi*y)'
[../]
[./p_source_func]
type = ParsedFunction
value = '-0.06*pi*x*cos(0.3*pi*x*y) - 0.14*pi*y*cos(0.2*pi*x*y) - 0.2*pi*cos(0.5*pi*x) - 0.09*pi*cos(0.3*pi*y)'
[../]
[./vel_x_func]
type = ParsedFunction
value = '0.4*sin(0.5*pi*x) + 0.4*sin(pi*y) + 0.7*sin(0.2*pi*x*y) + 0.5'
[../]
[./vel_y_func]
type = ParsedFunction
value = '0.6*sin(0.8*pi*x) + 0.3*sin(0.3*pi*y) + 0.2*sin(0.3*pi*x*y) + 0.3'
[../]
[./p_func]
type = ParsedFunction
value = '0.5*sin(0.5*pi*x) + 1.0*sin(0.3*pi*y) + 0.5*sin(0.2*pi*x*y) + 0.5'
[../]
[./vxx_func]
type = ParsedFunction
value = '0.14*pi*y*cos(0.2*pi*x*y) + 0.2*pi*cos(0.5*pi*x)'
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Transient
num_steps = 10
petsc_options = '-snes_converged_reason -ksp_converged_reason'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
line_search = 'none'
nl_rel_tol = 1e-8
nl_abs_tol = 1e-14
nl_max_its = 10
l_tol = 1e-6
l_max_its = 10
[./TimeStepper]
dt = .05
type = IterationAdaptiveDT
cutback_factor = 0.4
growth_factor = 1.2
optimal_iterations = 20
[../]
[]
[Outputs]
execute_on = 'final'
[./exodus]
type = Exodus
[../]
[./csv]
type = CSV
[../]
[]
[Postprocessors]
[./L2vel_x]
type = ElementL2Error
variable = vel_x
function = vel_x_func
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2vel_y]
variable = vel_y
function = vel_y_func
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2p]
variable = p
function = p_func
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2vxx]
variable = vxx
function = vxx_func
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]
[AuxVariables]
[./vxx]
family = MONOMIAL
order = FIRST
[../]
[]
[AuxKernels]
[./vxx]
type = VariableGradientComponent
component = x
variable = vxx
gradient_variable = vel_x
[../]
[]
modules/navier_stokes/test/tests/ins/jacobian_test/jacobian_test.i
# This input file tests the jacobians of many of the INS kernels
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.5
nx = 1
ny = 1
elem_type = QUAD9
[]
[Variables]
[./vel_x]
order = SECOND
family = LAGRANGE
[../]
[./vel_y]
order = SECOND
family = LAGRANGE
[../]
[./p]
order = FIRST
family = LAGRANGE
[../]
[./temp]
order = SECOND
family = LAGRANGE
[../]
[]
[Kernels]
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
integrate_p_by_parts = false
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
integrate_p_by_parts = false
[../]
[./x_mom_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_mom_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./temp]
type = INSTemperature
variable = temp
u = vel_x
v = vel_y
[../]
[./temp_time_deriv]
type = INSTemperatureTimeDerivative
variable = temp
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu k cp'
prop_values = '0.5 1.5 0.7 1.3'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
[../]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
[]
[ICs]
[./p]
type = RandomIC
variable = p
min = 0.5
max = 1.5
[../]
[./vel_x]
type = RandomIC
variable = vel_x
min = 0.5
max = 1.5
[../]
[./vel_y]
type = RandomIC
variable = vel_y
min = 0.5
max = 1.5
[../]
[./temp]
type = RandomIC
variable = temp
min = 0.5
max = 1.5
[../]
[]
modules/navier_stokes/test/tests/ins/stagnation/stagnation.i
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 2.0
ymin = 0
ymax = 2.0
nx = 20
ny = 20
elem_type = QUAD9
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = Newton
[../]
[]
[Executioner]
type = Transient
dt = 1.0
dtmin = 1.e-6
num_steps = 5
l_max_its = 100
nl_max_its = 10
nl_rel_tol = 1.e-9
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = 'asm 2 ilu 4'
[]
[Variables]
[./vel_x]
family = LAGRANGE
order = SECOND
[../]
[./vel_y]
family = LAGRANGE
order = SECOND
[../]
[./p]
family = LAGRANGE
order = FIRST
[../]
[]
[BCs]
[./u_in]
type = FunctionDirichletBC
boundary = 'top'
variable = vel_x
function = vel_x_inlet
[../]
[./v_in]
type = FunctionDirichletBC
boundary = 'top'
variable = vel_y
function = vel_y_inlet
[../]
[./vel_x_no_slip]
type = DirichletBC
boundary = 'left bottom'
variable = vel_x
value = 0
[../]
[./vel_y_no_slip]
type = DirichletBC
boundary = 'bottom'
variable = vel_y
value = 0
[../]
# Note: setting INSMomentumNoBCBC on the outlet boundary causes the
# matrix to be singular. The natural BC, on the other hand, is
# sufficient to specify the value of the pressure without requiring
# a pressure pin.
[]
[Functions]
[./vel_x_inlet]
type = ParsedFunction
value = 'k*x'
vars = 'k'
vals = '1'
[../]
[./vel_y_inlet]
type = ParsedFunction
value = '-k*y'
vars = 'k'
vals = '1'
[../]
[]
[Kernels]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 .01389' # 2/144
[../]
[]
[Outputs]
exodus = true
[./out]
type = CSV
execute_on = 'final'
[../]
[]
[VectorPostprocessors]
[./nodal_sample]
# Pick off the wall pressure values.
type = NodalValueSampler
variable = p
boundary = 'bottom'
sort_by = x
[../]
[]
modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_no_parts.i
# This input file tests several different things:
# .) The axisymmetric (RZ) form of the governing equations.
# .) An open boundary.
# .) Not integrating the pressure by parts, thereby requiring a pressure pin.
# .) Natural boundary condition at the outlet.
[GlobalParams]
integrate_p_by_parts = false
laplace = false
gravity = '0 0 0'
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = Newton
[../]
[]
[Executioner]
type = Transient
dt = 0.005
dtmin = 0.005
num_steps = 5
l_max_its = 100
# Note: The Steady executioner can be used for this problem, if you
# drop the INSMomentumTimeDerivative kernels and use the following
# direct solver options.
# petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount -ksp_type'
# petsc_options_value = 'lu NONZERO 1.e-10 preonly'
# Block Jacobi works well for this problem, as does "-pc_type asm
# -pc_asm_overlap 2", but an overlap of 1 does not work for some
# reason?
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = 'bjacobi ilu 4'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
csv = true
console = true
[./out]
type = Exodus
[../]
[]
[Variables]
[./vel_x]
# Velocity in radial (r) direction
family = LAGRANGE
order = SECOND
[../]
[./vel_y]
# Velocity in axial (z) direction
family = LAGRANGE
order = SECOND
[../]
[./p]
family = LAGRANGE
order = FIRST
[../]
[]
[BCs]
[./p_corner]
# This is required, because pressure term is *not* integrated by parts.
type = DirichletBC
boundary = top_right
value = 0
variable = p
[../]
[./u_out]
type = INSMomentumNoBCBCTractionForm
boundary = top
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./v_out]
type = INSMomentumNoBCBCTractionForm
boundary = top
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[./u_in]
type = DirichletBC
boundary = bottom
variable = vel_x
value = 0
[../]
[./v_in]
type = FunctionDirichletBC
boundary = bottom
variable = vel_y
function = 'inlet_func'
[../]
[./u_axis_and_walls]
type = DirichletBC
boundary = 'left right'
variable = vel_x
value = 0
[../]
[./v_no_slip]
type = DirichletBC
boundary = 'right'
variable = vel_y
value = 0
[../]
[]
[Kernels]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./mass]
type = INSMassRZ
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_space]
type = INSMomentumTractionFormRZ
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumTractionFormRZ
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 'volume'
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
value = '-4 * x^2 + 1'
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]
modules/navier_stokes/test/tests/ins/jeffery_hamel/wedge_dirichlet.i
# This input file tests whether we can converge to the semi-analytical
# solution for flow in a 2D wedge.
[GlobalParams]
gravity = '0 0 0'
# Params used by the WedgeFunction for computing the exact solution.
# The value of K is only required for comparing the pressure to the
# exact solution, and is computed by the associated jeffery_hamel.py
# script.
alpha_degrees = 15
Re = 30
K = -9.78221333616
f = f_theta
[]
[Mesh]
[file]
type = FileMeshGenerator
# file = wedge_4x6.e
file = wedge_8x12.e
# file = wedge_16x24.e
# file = wedge_32x48.e
# file = wedge_64x96.e
[]
[./corner_node]
# Pin is on the centerline of the channel on the left-hand side of
# the domain at r=1. If you change the domain, you will need to
# update this pin location for the pressure exact solution to
# work.
type = ExtraNodesetGenerator
new_boundary = pinned_node
coord = '1 0'
input = file
[../]
[]
[Variables]
[./vel_x]
order = SECOND
family = LAGRANGE
[../]
[./vel_y]
order = SECOND
family = LAGRANGE
[../]
[./p]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[BCs]
[./vel_x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top_wall bottom_wall'
value = 0.0
[../]
[./vel_y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'top_wall bottom_wall'
value = 0.0
[../]
[./vel_x_inlet]
type = FunctionDirichletBC
variable = vel_x
boundary = 'inlet outlet'
function = 'vel_x_exact'
[../]
[./vel_y_inlet]
type = FunctionDirichletBC
variable = vel_y
boundary = 'inlet outlet'
function = 'vel_y_exact'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 1
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Transient
dt = 1.e-2
dtmin = 1.e-2
num_steps = 5
petsc_options_iname = '-ksp_gmres_restart -pc_type -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = '300 bjacobi ilu 4'
line_search = none
nl_rel_tol = 1e-13
nl_abs_tol = 1e-11
nl_max_its = 10
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
exodus = true
[]
[Functions]
[./f_theta]
# Non-dimensional solution values f(eta), 0 <= eta <= 1 for
# alpha=15 deg, Re=30. Note: this introduces an input file
# ordering dependency: this Function must appear *before* the two
# functions below which use it since apparently proper dependency
# resolution is not done in this scenario.
type = PiecewiseLinear
data_file = 'f.csv'
format = 'columns'
[../]
[./vel_x_exact]
type = WedgeFunction
var_num = 0
mu = 1
rho = 1
[../]
[./vel_y_exact]
type = WedgeFunction
var_num = 1
mu = 1
rho = 1
[../]
[./p_exact]
type = WedgeFunction
var_num = 2
mu = 1
rho = 1
[../]
[]
[Postprocessors]
[./vel_x_L2_error]
type = ElementL2Error
variable = vel_x
function = vel_x_exact
execute_on = 'initial timestep_end'
[../]
[./vel_y_L2_error]
type = ElementL2Error
variable = vel_y
function = vel_y_exact
execute_on = 'initial timestep_end'
[../]
[./p_L2_error]
type = ElementL2Error
variable = p
function = p_exact
execute_on = 'initial timestep_end'
[../]
[]
modules/navier_stokes/test/tests/ins/jeffery_hamel/wedge_natural.i
# This input file solves the Jeffery-Hamel problem with the exact
# solution's outlet BC replaced by a natural BC. This problem does
# not converge to the analytical solution, although the flow at the
# outlet still "looks" reasonable.
[GlobalParams]
gravity = '0 0 0'
# Params used by the WedgeFunction for computing the exact solution.
# The value of K is only required for comparing the pressure to the
# exact solution, and is computed by the associated jeffery_hamel.py
# script.
alpha_degrees = 15
Re = 30
K = -9.78221333616
f = f_theta
[]
[Mesh]
file = wedge_8x12.e
[]
[Variables]
[./vel_x]
order = SECOND
family = LAGRANGE
[../]
[./vel_y]
order = SECOND
family = LAGRANGE
[../]
[./p]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[BCs]
[./vel_x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top_wall bottom_wall'
value = 0.0
[../]
[./vel_y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'top_wall bottom_wall'
value = 0.0
[../]
[./vel_x_inlet]
type = FunctionDirichletBC
variable = vel_x
boundary = 'inlet'
function = 'vel_x_exact'
[../]
[./vel_y_inlet]
type = FunctionDirichletBC
variable = vel_y
boundary = 'inlet'
function = 'vel_y_exact'
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 1
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Preconditioning]
[./SMP_NEWTON]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Transient
dt = 1.e-2
dtmin = 1.e-2
num_steps = 5
petsc_options_iname = '-ksp_gmres_restart -pc_type -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = '300 bjacobi ilu 4'
line_search = none
nl_rel_tol = 1e-13
nl_abs_tol = 1e-11
nl_max_its = 10
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
exodus = true
[]
[Functions]
[./f_theta]
# Non-dimensional solution values f(eta), 0 <= eta <= 1 for
# alpha=15deg, Re=30. Note: this introduces an input file
# ordering dependency: this Function must appear *before* the two
# function below which use it since apparently proper dependency
# resolution is not done in this scenario.
type = PiecewiseLinear
data_file = 'f.csv'
format = 'columns'
[../]
[./vel_x_exact]
type = WedgeFunction
var_num = 0
mu = 1
rho = 1
[../]
[./vel_y_exact]
type = WedgeFunction
var_num = 1
mu = 1
rho = 1
[../]
[]
modules/navier_stokes/test/tests/ins/lid_driven/lid_driven.i
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
elem_type = QUAD9
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./vel_x]
order = SECOND
family = LAGRANGE
[../]
[./vel_y]
order = SECOND
family = LAGRANGE
[../]
[./T]
order = SECOND
family = LAGRANGE
[./InitialCondition]
type = ConstantIC
value = 1.0
[../]
[../]
[./p]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
# mass
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
p = p
[../]
# x-momentum, time
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
# x-momentum, space
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
# y-momentum, time
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
# y-momentum, space
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
# temperature
[./temperature_time]
type = INSTemperatureTimeDerivative
variable = T
[../]
[./temperature_space]
type = INSTemperature
variable = T
u = vel_x
v = vel_y
[../]
[]
[BCs]
[./x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'bottom right left'
value = 0.0
[../]
[./lid]
type = FunctionDirichletBC
variable = vel_x
boundary = 'top'
function = 'lid_function'
[../]
[./y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'bottom right top left'
value = 0.0
[../]
[./T_hot]
type = DirichletBC
variable = T
boundary = 'bottom'
value = 1
[../]
[./T_cold]
type = DirichletBC
variable = T
boundary = 'top'
value = 0
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[]
[Functions]
[./lid_function]
# We pick a function that is exactly represented in the velocity
# space so that the Dirichlet conditions are the same regardless
# of the mesh spacing.
type = ParsedFunction
value = '4*x*(1-x)'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Transient
# Run for 100+ timesteps to reach steady state.
num_steps = 5
dt = .5
dtmin = .5
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = 'asm 2 ilu 4'
line_search = 'none'
nl_rel_tol = 1e-12
nl_abs_tol = 1e-13
nl_max_its = 6
l_tol = 1e-6
l_max_its = 500
[]
[Outputs]
file_base = lid_driven_out
exodus = true
perf_graph = true
[]
modules/navier_stokes/test/tests/ins/RZ_cone/RZ_cone_by_parts.i
# This input file tests several different things:
# .) The axisymmetric (RZ) form of the governing equations.
# .) An open boundary.
# .) Integrating the pressure by parts.
# .) Natural boundary condition at the outlet.
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = Newton
[../]
[]
[Executioner]
type = Transient
dt = 0.005
dtmin = 0.005
num_steps = 5
l_max_its = 100
# Note: The Steady executioner can be used for this problem, if you
# drop the INSMomentumTimeDerivative kernels and use the following
# direct solver options.
# petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount -ksp_type'
# petsc_options_value = 'lu NONZERO 1.e-10 preonly'
# Block Jacobi works well for this problem, as does "-pc_type asm
# -pc_asm_overlap 2", but an overlap of 1 does not work for some
# reason?
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = 'bjacobi ilu 4'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
csv = true
console = true
[./out]
type = Exodus
[../]
[]
[Variables]
[./vel_x]
# Velocity in radial (r) direction
family = LAGRANGE
order = SECOND
[../]
[./vel_y]
# Velocity in axial (z) direction
family = LAGRANGE
order = SECOND
[../]
[./p]
family = LAGRANGE
order = FIRST
[../]
[]
[BCs]
[./u_in]
type = DirichletBC
boundary = bottom
variable = vel_x
value = 0
[../]
[./v_in]
type = FunctionDirichletBC
boundary = bottom
variable = vel_y
function = 'inlet_func'
[../]
[./u_axis_and_walls]
type = DirichletBC
boundary = 'left right'
variable = vel_x
value = 0
[../]
[./v_no_slip]
type = DirichletBC
boundary = 'right'
variable = vel_y
value = 0
[../]
[]
[Kernels]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./y_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_y
[../]
[./mass]
type = INSMassRZ
variable = p
u = vel_x
v = vel_y
p = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceFormRZ
variable = vel_x
u = vel_x
v = vel_y
p = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceFormRZ
variable = vel_y
u = vel_x
v = vel_y
p = p
component = 1
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 'volume'
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
value = '-4 * x^2 + 1'
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]