- componentThe velocity component that this is applied to.
C++ Type:unsigned int
Controllable:No
Description:The velocity component that this is applied to.
- pressurepressure
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:pressure
- ux-velocity
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:x-velocity
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
INSMomentumLaplaceForm
This class computes momentum equation residual and Jacobian viscous contributions for the 'Laplacian' form of the governing equations.
This class implements the same terms as described in INSMomentumTractionForm but with a simplification applied to the viscous stress. For constant density and viscosity, the following simplification can be applied, as described in (Peterson et al., 2018):
\begin{aligned} \nabla\cdot\sigma &= \nabla\cdot\left(-p\bm{I} + \mu\left(\nabla\vec{u} + \left(\nabla\vec{u}\right)^T\right) &= -\nabla p + \mu\left(\nabla\cdot\nabla\vec{u} + \nabla\cdot\left(\nabla\vec{u}\right)^T\right) &= -\nabla p + \mu[left(\nabla\left(\nabla\cdot\vec{u}\right) + \nabla\cdot\left(\nabla\vec{u}\right)^T\right) &= -\nabla p + \mu\nabla\cdot\left(\nabla\vec{u}\right)^T &= -\nabla p + \mu\nabla^2\vec{u} \end{aligned}
Moving from lines 2 to 3 in the above equation we assumed sufficient smoothness in u to interchange the divergence and gradient operations. In moving from lines 3 to 4 we applied the incompressibility constraint ∇⋅u=0. The final form is a vector Laplacian (which can be different in non-Cartesian coordinate systems from a componentwise scalar Laplacian), hence the designation the "Laplace" form.
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>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- convective_termTrueWhether to include the convective term.
Default:True
C++ Type:bool
Controllable:No
Description:Whether to include the convective term.
- disp_xThe x displacement
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The x displacement
- disp_yThe y displacement
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The y displacement
- disp_zThe z displacement
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The z displacement
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
- forcing_func0The mms forcing function.
Default:0
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:The mms forcing function.
- gravity0 0 0Direction of the gravity vector
Default:0 0 0
C++ Type:libMesh::VectorValue<double>
Unit:(no unit assumed)
Controllable:No
Description:Direction of the gravity vector
- integrate_p_by_partsTrueWhether to integrate the pressure term by parts.
Default:True
C++ Type:bool
Controllable:No
Description:Whether to integrate the pressure term by parts.
- laplaceTrueWhether the viscous term of the momentum equations is in laplace form.
Default:True
C++ Type:bool
Controllable:No
Description:Whether the viscous term of the momentum equations is in laplace form.
- 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
- picardFalseWhether we are applying a Picard strategy in which case we will linearize the nonlinear convective term.
Default:False
C++ Type:bool
Controllable:No
Description:Whether we are applying a Picard strategy in which case we will linearize the nonlinear convective term.
- rho_namerhoThe name of the density
Default:rho
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The name of the density
- supgFalseWhether to perform SUPG stabilization of the momentum residuals
Default:False
C++ Type:bool
Controllable:No
Description:Whether to perform SUPG stabilization of the momentum residuals
- transient_termFalseWhether there should be a transient term in the momentum residuals.
Default:False
C++ Type:bool
Controllable:No
Description:Whether there should be a transient term in the momentum residuals.
- vy-velocity
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:y-velocity
- wz-velocity
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:z-velocity
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
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<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Contribution To Tagged Field Data 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.
- 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<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
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
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
Controllable:No
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<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
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
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_element/ins/mms/supg/supg_mms_test.i)
- (modules/navier_stokes/test/tests/finite_element/ins/velocity_channel/velocity_inletBC_by_parts.i)
- (modules/navier_stokes/test/tests/finite_element/ins/pressure_channel/open_bc_pressure_BC_fieldSplit.i)
- (modules/fsi/test/tests/2d-small-strain-transient/fsi_flat_channel.i)
- (modules/navier_stokes/test/tests/finite_element/ins/stagnation/stagnation.i)
- (modules/navier_stokes/test/tests/finite_element/ins/mms/supg/supg_pspg_adv_dominated_mms.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/transient_fsp.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/lid_driven_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/jacobian_test/jacobian_test.i)
- (modules/navier_stokes/test/tests/auxkernels/peclet-number-functor-aux/fe-thermal.i)
- (modules/navier_stokes/test/tests/finite_element/ins/velocity_channel/velocity_inletBC_no_parts.i)
- (modules/navier_stokes/test/tests/finite_element/ins/mms/supg/supg_adv_dominated_mms.i)
- (modules/navier_stokes/test/tests/auxkernels/reynolds-number-functor-aux/fe.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/steady_fsp.i)
- (modules/navier_stokes/test/tests/finite_element/ins/nonzero-malloc/test.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hydrostatic/gravity.i)
- (modules/navier_stokes/test/tests/finite_element/ins/jacobian_test/jacobian_stabilized_test.i)
- (modules/navier_stokes/test/tests/finite_element/ins/jeffery_hamel/wedge_dirichlet.i)
- (modules/navier_stokes/test/tests/finite_element/ins/pressure_channel/open_bc_pressure_BC.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/steady_fsp_diagonal_of_a_for_scaling.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/lid_driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/jeffery_hamel/wedge_natural.i)
- (modules/navier_stokes/test/tests/finite_element/ins/mms/pspg/pspg_mms_test.i)
Child Objects
References
- John W. Peterson, Alexander D. Lindsay, and Fande Kong.
Overview of the incompressible navier–stokes simulation capabilities in the moose framework.
Advances in Engineering Software, 119:68–92, 2018.[BibTeX]
@article{peterson2018overview, author = "Peterson, John W. and Lindsay, Alexander D. and Kong, Fande", title = "Overview of the incompressible Navier--Stokes simulation capabilities in the MOOSE framework", year = "2018", journal = "Advances in Engineering Software", volume = "119", pages = "68--92", publisher = "Elsevier" }
(modules/navier_stokes/test/tests/finite_element/ins/mms/supg/supg_mms_test.i)
mu=1.5
rho=2.5
[GlobalParams]
gravity = '0 0 0'
supg = true
convective_term = true
integrate_p_by_parts = false
laplace = true
u = vel_x
v = vel_y
pressure = p
alpha = 1
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-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
preset = false
boundary = 'left right top bottom'
function = vel_x_func
variable = vel_x
[../]
[./vel_y]
type = FunctionDirichletBC
preset = false
boundary = 'left right top bottom'
function = vel_y_func
variable = vel_y
[../]
[./p]
type = FunctionDirichletBC
preset = false
boundary = 'left right top bottom'
function = p_func
variable = p
[../]
[]
[Functions]
[./vel_x_source_func]
type = ParsedFunction
expression = '-${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
expression = '-${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
expression = '-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
expression = '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
expression = '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
expression = '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
expression = '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 = Steady
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-12
nl_abs_tol = 1e-13
nl_max_its = 6
l_tol = 1e-6
l_max_its = 500
[]
[Outputs]
[./exodus]
type = Exodus
[../]
[./csv]
type = CSV
[../]
[]
[Postprocessors]
[./L2vel_x]
type = ElementL2Error
variable = vel_x
function = vel_x_func
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2vel_y]
variable = vel_y
function = vel_y_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2p]
variable = p
function = p_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2vxx]
variable = vxx
function = vxx_func
type = ElementL2Error
outputs = 'console' 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/finite_element/ins/velocity_channel/velocity_inletBC_by_parts.i)
# This input file tests outflow boundary conditions for the incompressible NS equations.
[GlobalParams]
gravity = '0 0 0'
integrate_p_by_parts = true
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.0
nx = 30
ny = 10
elem_type = QUAD9
[]
[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
pressure = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[../]
[]
[BCs]
[./x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top bottom'
value = 0.0
[../]
[./y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'left top bottom'
value = 0.0
[../]
[./x_inlet]
type = FunctionDirichletBC
variable = vel_x
boundary = 'left'
function = 'inlet_func'
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Steady
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-12
nl_max_its = 6
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
[./out]
type = Exodus
[../]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * (y - 0.5)^2 + 1'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/pressure_channel/open_bc_pressure_BC_fieldSplit.i)
# This input file tests Dirichlet pressure in/outflow boundary conditions for the incompressible NS equations.
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.0
nx = 30
ny = 10
elem_type = QUAD9
[]
[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
pressure = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
integrate_p_by_parts = false
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
integrate_p_by_parts = false
[../]
[]
[BCs]
[./x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top bottom'
value = 0.0
[../]
[./y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'left top bottom'
value = 0.0
[../]
[./inlet_p]
type = DirichletBC
variable = p
boundary = left
value = 1.0
[../]
[./outlet_p]
type = DirichletBC
variable = p
boundary = right
value = 0.0
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Preconditioning]
active = FSP
[./FSP]
type = FSP
# It is the starting point of splitting
topsplit = 'up' # 'up' should match the following block name
[./up]
splitting = 'u p' # 'u' and 'p' are the names of subsolvers
splitting_type = schur
# Splitting type is set as schur, because the pressure part of Stokes-like systems
# is not diagonally dominant. CAN NOT use additive, multiplicative and etc.
# Original system:
# | A B | | u | = | f_u |
# | C 0 | | p | | f_v |
# is factorized into
# |I 0 | | A 0| | I A^{-1}B | | u | = | f_u |
# |CA^{-1} I | | 0 -S| | 0 I | | p | | f_v |
# S = CA^{-1}B
# The preconditioning is accomplished via the following steps
# (1) p^{(0)} = f_v - CA^{-1}f_u,
# (2) pressure = (-S)^{-1} p^{(0)}
# (3) u = A^{-1}(f_u-Bp)
petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition'
petsc_options_value = 'full selfp'
# Factorization type here is full, which means we approximate the original system
# exactly. There are three other options:
# diag:
# | A 0|
# | 0 -S|
# lower:
# |I 0 |
# |CA^{-1} -S |
# upper:
# | I A^{-1}B |
# | 0 -S |
# The preconditioning matrix is set as selfp, which means we explicitly form a
# matrix \hat{S} = C(diag(A))^{-1}B. We do not compute the inverse of A, but instead, we compute
# the inverse of diag(A).
[../]
[./u]
vars = 'vel_x vel_y'
# PETSc options for this subsolver
# A prefix will be applied, so just put the options for this subsolver only
petsc_options_iname = '-pc_type -ksp_type -ksp_rtol'
petsc_options_value = ' hypre gmres 1e-4'
# Specify options to solve A^{-1} in the steps (1), (2) and (3).
# Solvers for A^{-1} could be different in different steps. We could
# choose in the following pressure block.
[../]
[./p]
vars = 'p'
# PETSc options for this subsolver in the step (2)
petsc_options_iname = '-pc_type -ksp_type -ksp_rtol'
petsc_options_value = ' jacobi gmres 1e-4'
# Use -inner_ksp_type and -inner_pc_type to override A^{-1} in the step (2)
# Use -lower_ksp_type and -lower_pc_type to override A^{-1} in the step (1)
[../]
[../]
[]
[Executioner]
type = Steady
solve_type = PJFNK
nl_rel_tol = 1e-12
nl_max_its = 6
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
file_base = open_bc_out_pressure_BC_fieldSplit
exodus = true
[]
(modules/fsi/test/tests/2d-small-strain-transient/fsi_flat_channel.i)
[GlobalParams]
gravity = '0 0 0'
integrate_p_by_parts = true
laplace = true
convective_term = true
transient_term = true
pspg = true
supg = true
displacements = 'disp_x disp_y'
preset = false
order = FIRST
use_displaced_mesh = true
[]
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.0
nx = 10
ny = 15
elem_type = QUAD4
[]
[subdomain1]
type = SubdomainBoundingBoxGenerator
bottom_left = '0.0 0.5 0'
block_id = 1
top_right = '3.0 1.0 0'
input = gmg
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
primary_block = '0'
paired_block = '1'
new_boundary = 'master0_interface'
input = subdomain1
[]
[break_boundary]
type = BreakBoundaryOnSubdomainGenerator
input = interface
[]
[]
[Variables]
[./vel_x]
block = 0
[../]
[./vel_y]
block = 0
[../]
[./p]
block = 0
[../]
[./disp_x]
[../]
[./disp_y]
[../]
[./vel_x_solid]
block = 1
[../]
[./vel_y_solid]
block = 1
[../]
[]
[Kernels]
[./vel_x_time]
type = INSMomentumTimeDerivative
variable = vel_x
block = 0
[../]
[./vel_y_time]
type = INSMomentumTimeDerivative
variable = vel_y
block = 0
[../]
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
pressure = p
block = 0
disp_x = disp_x
disp_y = disp_y
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
block = 0
disp_x = disp_x
disp_y = disp_y
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
block = 0
disp_x = disp_x
disp_y = disp_y
[../]
[./vel_x_mesh]
type = ConvectedMesh
disp_x = disp_x
disp_y = disp_y
variable = vel_x
u = vel_x
v = vel_y
pressure = p
block = 0
[../]
[./vel_y_mesh]
type = ConvectedMesh
disp_x = disp_x
disp_y = disp_y
variable = vel_y
u = vel_x
v = vel_y
pressure = p
block = 0
[../]
[./p_mesh]
type = ConvectedMeshPSPG
disp_x = disp_x
disp_y = disp_y
variable = p
u = vel_x
v = vel_y
pressure = p
block = 0
[../]
[./disp_x_fluid]
type = Diffusion
variable = disp_x
block = 0
use_displaced_mesh = false
[../]
[./disp_y_fluid]
type = Diffusion
variable = disp_y
block = 0
use_displaced_mesh = false
[../]
[./accel_tensor_x]
type = CoupledTimeDerivative
variable = disp_x
v = vel_x_solid
block = 1
use_displaced_mesh = false
[../]
[./accel_tensor_y]
type = CoupledTimeDerivative
variable = disp_y
v = vel_y_solid
block = 1
use_displaced_mesh = false
[../]
[./vxs_time_derivative_term]
type = CoupledTimeDerivative
variable = vel_x_solid
v = disp_x
block = 1
use_displaced_mesh = false
[../]
[./vys_time_derivative_term]
type = CoupledTimeDerivative
variable = vel_y_solid
v = disp_y
block = 1
use_displaced_mesh = false
[../]
[./source_vxs]
type = MatReaction
variable = vel_x_solid
block = 1
reaction_rate = 1
use_displaced_mesh = false
[../]
[./source_vys]
type = MatReaction
variable = vel_y_solid
block = 1
reaction_rate = 1
use_displaced_mesh = false
[../]
[]
[InterfaceKernels]
[./penalty_interface_x]
type = CoupledPenaltyInterfaceDiffusion
variable = vel_x
neighbor_var = disp_x
secondary_coupled_var = vel_x_solid
boundary = master0_interface
penalty = 1e6
[../]
[./penalty_interface_y]
type = CoupledPenaltyInterfaceDiffusion
variable = vel_y
neighbor_var = disp_y
secondary_coupled_var = vel_y_solid
boundary = master0_interface
penalty = 1e6
[../]
[]
[Physics/SolidMechanics/QuasiStatic]
[./solid_domain]
strain = SMALL
incremental = false
# generate_output = 'strain_xx strain_yy strain_zz' ## Not at all necessary, but nice
block = '1'
[../]
[]
[Materials]
[./elasticity_tensor]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 1e2
poissons_ratio = 0.3
block = '1'
use_displaced_mesh = false
[../]
[./small_stress]
type = ComputeLinearElasticStress
block = 1
[../]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 1'
use_displaced_mesh = false
[../]
[]
[BCs]
[./fluid_x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'bottom'
value = 0.0
[../]
[./fluid_y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'bottom left_to_0'
value = 0.0
[../]
[./x_inlet]
type = FunctionDirichletBC
variable = vel_x
boundary = 'left_to_0'
function = 'inlet_func'
[../]
[./no_disp_x]
type = DirichletBC
variable = disp_x
boundary = 'bottom top left_to_1 right_to_1 left_to_0 right_to_0'
value = 0
[../]
[./no_disp_y]
type = DirichletBC
variable = disp_y
boundary = 'bottom top left_to_1 right_to_1 left_to_0 right_to_0'
value = 0
[../]
[./solid_x_no_slip]
type = DirichletBC
variable = vel_x_solid
boundary = 'top left_to_1 right_to_1'
value = 0.0
[../]
[./solid_y_no_slip]
type = DirichletBC
variable = vel_y_solid
boundary = 'top left_to_1 right_to_1'
value = 0.0
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
num_steps = 5
# num_steps = 60
dt = 0.1
dtmin = 0.1
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
line_search = none
nl_rel_tol = 1e-50
nl_abs_tol = 1e-10
[]
[Outputs]
[./out]
type = Exodus
[../]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '(-16 * (y - 0.25)^2 + 1) * (1 + cos(t))'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/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 = 15
nl_rel_tol = 1.e-9
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
petsc_options_value = 'asm 2 lu NONZERO 1000'
line_search = none
[]
[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
expression = 'k*x'
symbol_names = 'k'
symbol_values = '1'
[../]
[./vel_y_inlet]
type = ParsedFunction
expression = '-k*y'
symbol_names = 'k'
symbol_values = '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
pressure = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = 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/finite_element/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
pressure = 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
expression = '-${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
expression = '-${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
expression = '-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
expression = '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
expression = '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
expression = '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
expression = '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' execute_on = 'timestep_end'
[../]
[./L2vel_y]
variable = vel_y
function = vel_y_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2p]
variable = p
function = p_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2vxx]
variable = vxx
function = vxx_func
type = ElementL2Error
outputs = 'console' 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/finite_element/ins/lid_driven/transient_fsp.i)
n=64
mu=2e-3
[GlobalParams]
gravity = '0 0 0'
preset = true
supg = false
[]
[Problem]
extra_tag_matrices = 'mass'
previous_nl_solution_required = true
type = NavierStokesProblem
mass_matrix = 'mass'
schur_fs_index = '1'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = ${n}
ny = ${n}
elem_type = QUAD9
[]
[]
[Variables]
[vel_x]
order = SECOND
family = LAGRANGE
[]
[vel_y]
order = SECOND
family = LAGRANGE
[]
[p]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
# mass
[mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
pressure = p
[]
[x_time]
type = INSMomentumTimeDerivative
variable = vel_x
[]
[x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[]
[x_mass]
type = MassMatrix
variable = vel_x
matrix_tags = 'mass'
[]
[y_time]
type = INSMomentumTimeDerivative
variable = vel_y
[]
[y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[]
[y_mass]
type = MassMatrix
variable = vel_y
matrix_tags = 'mass'
[]
[]
[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
[]
[]
[Materials]
[const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 ${mu}'
[]
[]
[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
expression = '4*x*(1-x)'
[]
[]
[Preconditioning]
[FSP]
type = FSP
topsplit = 'by_diri_others'
[by_diri_others]
splitting = 'diri others'
splitting_type = additive
petsc_options_iname = '-ksp_type'
petsc_options_value = 'preonly'
[]
[diri]
sides = 'left right top bottom'
vars = 'vel_x vel_y'
petsc_options_iname = '-pc_type'
petsc_options_value = 'jacobi'
[]
[others]
splitting = 'u p'
splitting_type = schur
petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition -ksp_gmres_restart -ksp_rtol -ksp_type -ksp_atol'
petsc_options_value = 'full self 300 1e-5 fgmres 1e-9'
unside_by_var_boundary_name = 'left top right bottom left top right bottom'
unside_by_var_var_name = 'vel_x vel_x vel_x vel_x vel_y vel_y vel_y vel_y'
[]
[u]
vars = 'vel_x vel_y'
unside_by_var_boundary_name = 'left top right bottom left top right bottom'
unside_by_var_var_name = 'vel_x vel_x vel_x vel_x vel_y vel_y vel_y vel_y'
# petsc_options = '-ksp_converged_reason'
petsc_options_iname = '-pc_type -ksp_pc_side -ksp_type -ksp_rtol -pc_hypre_type -ksp_gmres_restart'
petsc_options_value = 'hypre right gmres 1e-2 boomeramg 300'
[]
[p]
vars = 'p'
petsc_options = '-pc_lsc_scale_diag -ksp_converged_reason'# -lsc_ksp_converged_reason -lsc_ksp_monitor_true_residual
petsc_options_iname = '-ksp_type -ksp_gmres_restart -ksp_rtol -pc_type -ksp_pc_side -lsc_pc_type -lsc_pc_hypre_type -lsc_ksp_type -lsc_ksp_rtol -lsc_ksp_pc_side -lsc_ksp_gmres_restart'
petsc_options_value = 'fgmres 300 1e-2 lsc right hypre boomeramg gmres 1e-1 right 300'
[]
[]
[]
[Postprocessors]
[pavg]
type = ElementAverageValue
variable = p
[]
[]
[UserObjects]
[set_pressure]
type = NSPressurePin
pin_type = 'average'
variable = p
pressure_average = 'pavg'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
petsc_options_iname = '-snes_max_it'
petsc_options_value = '100'
line_search = 'none'
nl_rel_tol = 1e-8
nl_abs_tol = 1e-9
abort_on_solve_fail = true
normalize_solution_diff_norm_by_dt = false
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 6
dt = 1e-2
[]
steady_state_detection = true
[]
[Outputs]
[exo]
type = Exodus
execute_on = 'final'
hide = 'pavg'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/lid_driven_stabilized.i)
[GlobalParams]
gravity = '0 0 0'
laplace = true
integrate_p_by_parts = true
family = LAGRANGE
order = FIRST
# There are multiple types of stabilization possible in incompressible
# Navier Stokes. The user can specify supg = true to apply streamline
# upwind petrov-galerkin stabilization to the momentum equations. This
# is most useful for high Reynolds numbers, e.g. when inertial effects
# dominate over viscous effects. The user can also specify pspg = true
# to apply pressure stabilized petrov-galerkin stabilization to the mass
# equation. PSPG is a form of Galerkin Least Squares. This stabilization
# allows equal order interpolations to be used for pressure and velocity.
# Finally, the alpha parameter controls the amount of stabilization.
# For PSPG, decreasing alpha leads to increased accuracy but may induce
# spurious oscillations in the pressure field. Some numerical experiments
# suggest that alpha between .1 and 1 may be optimal for accuracy and
# robustness.
supg = true
pspg = true
alpha = 1e-1
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 64
ny = 64
elem_type = QUAD4
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./vel_x]
[../]
[./vel_y]
[../]
[./p]
[../]
[]
[Kernels]
# mass
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
pressure = p
[../]
# x-momentum, space
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[../]
# y-momentum, space
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[../]
[]
[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
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[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
expression = '4*x*(1-x)'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
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]
exodus = true
print_linear_converged_reason = false
print_nonlinear_converged_reason = false
[]
[Postprocessors]
[lin]
type = NumLinearIterations
[]
[nl]
type = NumNonlinearIterations
[]
[lin_tot]
type = CumulativeValuePostprocessor
postprocessor = 'lin'
[]
[nl_tot]
type = CumulativeValuePostprocessor
postprocessor = 'nl'
[]
[]
(modules/navier_stokes/test/tests/finite_element/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
pressure = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
integrate_p_by_parts = false
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = 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/auxkernels/peclet-number-functor-aux/fe-thermal.i)
rho = 1
mu = 1
k = 1
cp = 1
[GlobalParams]
gravity = '0 0 0'
pspg = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[]
[corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[]
[]
[AuxVariables]
[Pe]
family = MONOMIAL
order = FIRST
[]
[]
[AuxKernels]
[Pe]
type = PecletNumberFunctorAux
variable = Pe
speed = speed
thermal_diffusivity = 'thermal_diffusivity'
[]
[]
[Variables]
[vel_x][]
[vel_y][]
[p][]
[T][]
[]
[Kernels]
# mass
[mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
pressure = p
[]
# x-momentum, space
[x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[]
# y-momentum, space
[y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[]
[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
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[T_hot]
type = DirichletBC
variable = T
boundary = 'bottom'
value = 1
[]
[T_cold]
type = DirichletBC
variable = T
boundary = 'top'
value = 0
[]
[]
[Materials]
[const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu k cp'
prop_values = '${rho} ${mu} ${k} ${cp}'
[]
[speed]
type = ADVectorMagnitudeFunctorMaterial
x_functor = vel_x
y_functor = vel_y
vector_magnitude_name = speed
[]
[thermal_diffusivity]
type = ThermalDiffusivityFunctorMaterial
k = ${k}
rho = ${rho}
cp = ${cp}
[]
[]
[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
expression = '4*x*(1-x)'
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type'
petsc_options_value = 'asm 2 lu'
line_search = 'none'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/velocity_channel/velocity_inletBC_no_parts.i)
# This input file tests outflow boundary conditions for the incompressible NS equations.
[GlobalParams]
gravity = '0 0 0'
integrate_p_by_parts = false
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.0
nx = 30
ny = 10
elem_type = QUAD9
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = top_right
coord = '3 1'
input = gen
[../]
[]
[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
pressure = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[../]
[]
[BCs]
[./x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top bottom'
value = 0.0
[../]
[./y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'left top bottom'
value = 0.0
[../]
[./x_inlet]
type = FunctionDirichletBC
variable = vel_x
boundary = 'left'
function = 'inlet_func'
[../]
[./p_corner]
# Since the pressure is not integrated by parts in this example,
# it is only specified up to a constant by the natural outflow BC.
# Therefore, we need to pin its value at a single location.
type = DirichletBC
boundary = top_right
value = 0
variable = p
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Steady
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-12
nl_max_its = 6
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
[./out]
type = Exodus
[../]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * (y - 0.5)^2 + 1'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/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
pressure = 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
expression = '-${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
expression = '-${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
expression = '-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
expression = '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
expression = '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
expression = '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
expression = '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' execute_on = 'timestep_end'
[../]
[./L2vel_y]
variable = vel_y
function = vel_y_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2p]
variable = p
function = p_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2vxx]
variable = vxx
function = vxx_func
type = ElementL2Error
outputs = 'console' 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/auxkernels/reynolds-number-functor-aux/fe.i)
rho=1
mu=1
[GlobalParams]
gravity = '0 0 0'
pspg = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[]
[corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[]
[]
[AuxVariables]
[Reynolds]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[Reynolds]
type = ReynoldsNumberFunctorAux
variable = Reynolds
speed = speed
rho = ${rho}
mu = ${mu}
[]
[]
[Variables]
[vel_x]
[]
[vel_y]
[]
[p]
[]
[]
[Kernels]
# mass
[mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
pressure = p
[]
# x-momentum, space
[x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[]
# y-momentum, space
[y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[]
[]
[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
[]
[pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[]
[]
[Materials]
[const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[speed]
type = ADVectorMagnitudeFunctorMaterial
x_functor = vel_x
y_functor = vel_y
vector_magnitude_name = speed
[]
[]
[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
expression = '4*x*(1-x)'
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type'
petsc_options_value = 'asm 2 lu'
line_search = 'none'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/steady_fsp.i)
rho=1
mu=2e-3
U=1
l=1
prefactor=${fparse 1/(l/2)^2}
n=64
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
[gen]
type = DistributedRectilinearMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = QUAD4
[]
second_order = true
parallel_type = distributed
[]
[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
pressure = p
[]
[x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[]
[momentum_x_mass]
type = MassMatrix
variable = vel_x
density = ${rho}
matrix_tags = 'mass'
[]
[y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[]
[momentum_y_mass]
type = MassMatrix
variable = vel_y
density = ${rho}
matrix_tags = 'mass'
[]
[]
[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
[]
[]
[Materials]
[const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[]
[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
expression = '${prefactor}*${U}*x*(${l}-x)'
[]
[]
[Problem]
type = NavierStokesProblem
mass_matrix = 'mass'
extra_tag_matrices = 'mass'
[]
[Preconditioning]
[FSP]
type = FSP
topsplit = 'up'
[up]
splitting = 'u p'
splitting_type = schur
petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition -ksp_gmres_restart -ksp_type -ksp_pc_side -ksp_rtol'
petsc_options_value = 'full self 300 fgmres right 1e-4'
[]
[u]
vars = 'vel_x vel_y'
# petsc_options = '-ksp_converged_reason'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_type -ksp_rtol -ksp_gmres_restart -ksp_pc_side'
petsc_options_value = 'hypre boomeramg gmres 1e-2 300 right'
[]
[p]
vars = 'p'
petsc_options = '-pc_lsc_scale_diag -ksp_converged_reason'# -lsc_ksp_converged_reason -lsc_ksp_monitor_true_residual
petsc_options_iname = '-ksp_type -ksp_gmres_restart -ksp_rtol -pc_type -ksp_pc_side -lsc_pc_type -lsc_pc_hypre_type -lsc_ksp_type -lsc_ksp_rtol -lsc_ksp_pc_side -lsc_ksp_gmres_restart'
petsc_options_value = 'fgmres 300 1e-2 lsc right hypre boomeramg gmres 1e-1 right 300'
[]
[]
[]
[Postprocessors]
[pavg]
type = ElementAverageValue
variable = p
[]
[]
[UserObjects]
[set_pressure]
type = NSPressurePin
pin_type = 'average'
variable = p
pressure_average = 'pavg'
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
[exo]
type = Exodus
execute_on = 'final'
hide = 'pavg'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/nonzero-malloc/test.i)
[GlobalParams]
gravity = '0 0 0'
pspg = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 5
ny = 5
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./vel_x]
[../]
[./vel_y]
[../]
[./T]
[./InitialCondition]
type = ConstantIC
value = 1.0
[../]
[../]
[./p]
[../]
[]
[Kernels]
# mass
[./mass]
type = INSMass
variable = p
u = vel_x
v = vel_y
pressure = 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
pressure = 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
pressure = p
component = 1
[../]
# temperature
[./temperature_time]
type = INSTemperatureTimeDerivative
variable = T
[../]
[./temperature_space]
type = INSTemperature
variable = T
u = vel_x
v = vel_y
[../]
[malloc]
type = MallocKernel
# Variable choice doesn't matter
variable = vel_x
[]
[]
[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
expression = '4*x*(1-x)'
[../]
[]
[Executioner]
type = Transient
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
perf_graph = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/hydrostatic/gravity.i)
[GlobalParams]
gravity = '0 -0.001 0'
convective_term = false
integrate_p_by_parts = false
u = vel_x
v = vel_y
pressure = p
[]
[Mesh]
second_order = true
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 1
ny = 5
ymax = 5
[../]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = top_right
coord = '0 5'
input = gen
[../]
[]
[Variables]
[./vel_x]
order = SECOND
[../]
[./vel_y]
order = SECOND
[../]
[./p]
[../]
[]
[Kernels]
[./mass]
type = INSMass
variable = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
component = 1
[../]
[]
[BCs]
[./x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top bottom left right'
value = 0.0
[../]
[./y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'top bottom left right'
value = 0.0
[../]
[./p_corner]
type = DirichletBC
boundary = top_right
value = 0
variable = p
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
prop_names = 'rho mu'
prop_values = '100 1'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Steady
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-12
nl_max_its = 6
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
exodus = true
execute_on = TIMESTEP_END
[]
(modules/navier_stokes/test/tests/finite_element/ins/jacobian_test/jacobian_stabilized_test.i)
# This input file tests the jacobians of many of the INS kernels
[GlobalParams]
gravity = '1.1 1.1 1.1'
u = vel_x
v = vel_y
w = vel_z
pressure = p
integrate_p_by_parts = true
laplace = true
pspg = true
supg = true
alpha = 1.1
[]
[Mesh]
type = GeneratedMesh
dim = 3
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.5
zmax = 1.1
nx = 1
ny = 1
nz = 1
elem_type = HEX27
[]
[Variables]
[./vel_x]
order = SECOND
family = LAGRANGE
[../]
[./vel_y]
order = SECOND
family = LAGRANGE
[../]
[./vel_z]
order = SECOND
family = LAGRANGE
[../]
[./p]
order = SECOND
family = LAGRANGE
[../]
[]
[Kernels]
[./mass]
type = INSMass
variable = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
component = 0
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
component = 1
[../]
[./z_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_z
component = 2
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '0.5 1.5'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
[../]
[]
[Executioner]
solve_type = NEWTON
type = Steady
petsc_options_iname = '-snes_type'
petsc_options_value = 'test'
[]
[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
[../]
[./vel_z]
type = RandomIC
variable = vel_z
min = 0.5
max = 1.5
[../]
[]
(modules/navier_stokes/test/tests/finite_element/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
pressure = p
[../]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = 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
pressure = 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/finite_element/ins/pressure_channel/open_bc_pressure_BC.i)
# This input file tests Dirichlet pressure in/outflow boundary conditions for the incompressible NS equations.
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 3.0
ymin = 0
ymax = 1.0
nx = 30
ny = 10
elem_type = QUAD9
[]
[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
pressure = p
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
integrate_p_by_parts = false
[../]
[./y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
integrate_p_by_parts = false
[../]
[]
[BCs]
[./x_no_slip]
type = DirichletBC
variable = vel_x
boundary = 'top bottom'
value = 0.0
[../]
[./y_no_slip]
type = DirichletBC
variable = vel_y
boundary = 'left top bottom'
value = 0.0
[../]
[./inlet_p]
type = DirichletBC
variable = p
boundary = left
value = 1.0
[../]
[./outlet_p]
type = DirichletBC
variable = p
boundary = right
value = 0.0
[../]
[]
[Materials]
[./const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = PJFNK
[../]
[]
[Executioner]
type = Steady
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-12
nl_max_its = 6
l_tol = 1e-6
l_max_its = 300
[]
[Outputs]
file_base = open_bc_out_pressure_BC
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/steady_fsp_diagonal_of_a_for_scaling.i)
rho=1
mu=2e-3
U=1
l=1
prefactor=${fparse 1/(l/2)^2}
n=64
[GlobalParams]
gravity = '0 0 0'
[]
[Mesh]
[gen]
type = DistributedRectilinearMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = QUAD4
[]
second_order = true
parallel_type = distributed
[]
[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
pressure = p
[]
[x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = p
component = 0
[]
[y_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_y
u = vel_x
v = vel_y
pressure = p
component = 1
[]
[]
[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
[]
[]
[Materials]
[const]
type = GenericConstantMaterial
block = 0
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[]
[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
expression = '${prefactor}*${U}*x*(${l}-x)'
[]
[]
[Problem]
type = NavierStokesProblem
[]
[Preconditioning]
[FSP]
type = FSP
topsplit = 'up'
[up]
splitting = 'u p'
splitting_type = schur
petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition -ksp_gmres_restart -ksp_type -ksp_pc_side -ksp_rtol'
petsc_options_value = 'full self 300 fgmres right 1e-4'
[]
[u]
vars = 'vel_x vel_y'
# petsc_options = '-ksp_converged_reason'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_type -ksp_rtol -ksp_gmres_restart -ksp_pc_side'
petsc_options_value = 'hypre boomeramg gmres 1e-2 300 right'
[]
[p]
vars = 'p'
petsc_options = '-pc_lsc_scale_diag -ksp_converged_reason'# -lsc_ksp_converged_reason -lsc_ksp_monitor_true_residual
petsc_options_iname = '-ksp_type -ksp_gmres_restart -ksp_rtol -pc_type -ksp_pc_side -lsc_pc_type -lsc_pc_hypre_type -lsc_ksp_type -lsc_ksp_rtol -lsc_ksp_pc_side -lsc_ksp_gmres_restart'
petsc_options_value = 'fgmres 300 1e-2 lsc right hypre boomeramg gmres 1e-1 right 300'
[]
[]
[]
[Postprocessors]
[pavg]
type = ElementAverageValue
variable = p
[]
[]
[UserObjects]
[set_pressure]
type = NSPressurePin
pin_type = 'average'
variable = p
pressure_average = 'pavg'
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
[exo]
type = Exodus
execute_on = 'final'
hide = 'pavg'
[]
[]
(modules/navier_stokes/test/tests/finite_element/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
pressure = 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
pressure = 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
pressure = 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
expression = '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/finite_element/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
pressure = p
[../]
[./x_momentum_time]
type = INSMomentumTimeDerivative
variable = vel_x
[../]
[./x_momentum_space]
type = INSMomentumLaplaceForm
variable = vel_x
u = vel_x
v = vel_y
pressure = 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
pressure = 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/finite_element/ins/mms/pspg/pspg_mms_test.i)
mu=1.5
rho=2.5
[GlobalParams]
gravity = '0 0 0'
pspg = true
convective_term = true
integrate_p_by_parts = true
laplace = true
u = vel_x
v = vel_y
pressure = p
alpha = 1e-6
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]
[../]
[]
[Kernels]
# mass
[./mass]
type = INSMass
variable = p
x_vel_forcing_func = vel_x_source_func
y_vel_forcing_func = vel_y_source_func
[../]
# 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
expression = '-${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
expression = '-${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
expression = '-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
expression = '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
expression = '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
expression = '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
expression = '0.14*pi*y*cos(0.2*pi*x*y) + 0.2*pi*cos(0.5*pi*x)'
[../]
[./px_func]
type = ParsedFunction
expression = '0.1*pi*y*cos(0.2*pi*x*y) + 0.25*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 = Steady
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-12
nl_abs_tol = 1e-13
nl_max_its = 6
l_tol = 1e-6
l_max_its = 500
[]
[Outputs]
[./exodus]
type = Exodus
[../]
[./csv]
type = CSV
[../]
[]
[Postprocessors]
[./L2vel_x]
type = ElementL2Error
variable = vel_x
function = vel_x_func
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2vel_y]
variable = vel_y
function = vel_y_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2p]
variable = p
function = p_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2vxx]
variable = vxx
function = vxx_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[./L2px]
variable = px
function = px_func
type = ElementL2Error
outputs = 'console' execute_on = 'timestep_end'
[../]
[]
[AuxVariables]
[./vxx]
family = MONOMIAL
order = FIRST
[../]
[./px]
family = MONOMIAL
order = FIRST
[../]
[]
[AuxKernels]
[./vxx]
type = VariableGradientComponent
component = x
variable = vxx
gradient_variable = vel_x
[../]
[./px]
type = VariableGradientComponent
component = x
variable = px
gradient_variable = p
[../]
[]
(modules/navier_stokes/include/kernels/INSMomentumLaplaceFormRZ.h)
// This file is part of the MOOSE framework
// https://mooseframework.inl.gov
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "INSMomentumLaplaceForm.h"
// Forward Declarations
/**
* This class computes additional momentum equation residual and
* Jacobian contributions for the incompressible Navier-Stokes
* momentum equation in RZ (axisymmetric cylindrical) coordinates,
* using the "Laplace" form of the governing equations.
*/
class INSMomentumLaplaceFormRZ : public INSMomentumLaplaceForm
{
public:
static InputParameters validParams();
INSMomentumLaplaceFormRZ(const InputParameters & parameters);
virtual ~INSMomentumLaplaceFormRZ() {}
protected:
virtual RealVectorValue strongViscousTermLaplace() override;
virtual RealVectorValue dStrongViscDUCompLaplace(unsigned comp) override;
virtual Real computeQpResidual() override;
virtual Real computeQpJacobian() override;
virtual Real computeQpOffDiagJacobian(unsigned jvar) override;
};