- variableThe variable this initial condition is supposed to provide values for.
C++ Type:VariableName
Unit:(no unit assumed)
Controllable:No
Description:The variable this initial condition is supposed to provide values for.
- x_valueThe x value to be set in IC
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The x value to be set in IC
VectorConstantIC
The VectorConstantIC
class is used to set initial values of components of a vector variable. The x component can be set through the x_value
parameter, the y component through y_value
, and the z component through z_value
. Note that x_value
is required. If y_value
or z_value
are not supplied, they are defaulted to zero.
Sets constant component values for a vector field variable.
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- stateCURRENTThis parameter is used to set old state solutions at the start of simulation. If specifying multiple states at the start of simulation, use one IC object for each state being specified. The states are CURRENT=0 OLD=1 OLDER=2. States older than 2 are not currently supported. When the user only specifies current state, the solution is copied to the old and older states, as expected. This functionality is mainly used for dynamic simulations with explicit time integration schemes, where old solution states are used in the velocity and acceleration approximations.
Default:CURRENT
C++ Type:MooseEnum
Options:CURRENT, OLD, OLDER
Controllable:No
Description:This parameter is used to set old state solutions at the start of simulation. If specifying multiple states at the start of simulation, use one IC object for each state being specified. The states are CURRENT=0 OLD=1 OLDER=2. States older than 2 are not currently supported. When the user only specifies current state, the solution is copied to the old and older states, as expected. This functionality is mainly used for dynamic simulations with explicit time integration schemes, where old solution states are used in the velocity and acceleration approximations.
- y_value0The y value to be set in IC
Default:0
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The y value to be set in IC
- z_value0The z value to be set in IC
Default:0
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The z value to be set in IC
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:No
Description:Set the enabled status of the MooseObject.
- ignore_uo_dependencyFalseWhen set to true, a UserObject retrieved by this IC will not be executed before the this IC
Default:False
C++ Type:bool
Controllable:No
Description:When set to true, a UserObject retrieved by this IC will not be executed before the this IC
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/energy_source/steady-var.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_stab_jac_test.i)
- (modules/fsi/test/tests/2d-finite-strain-steady/thermal-me.i)
- (modules/navier_stokes/test/tests/finite_element/ins/coupled-force/gravity-through-coupled-force.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_traction_steady_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/coupled-force/gravity-object.i)
- (modules/navier_stokes/test/tests/finite_element/ins/boussinesq/benchmark/benchmark.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad-rz-displacements.i)
- (test/tests/kernels/ad_transient_diffusion/ad_transient_vector_diffusion.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/mixed-transient-steady/mixed.i)
- (modules/navier_stokes/test/tests/finite_element/ins/coupled-force/steady.i)
- (test/tests/ics/vector_constant_ic/vector_constant_ic.i)
- (test/tests/kernels/transient_vector_diffusion/transient_vector_diffusion.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized_second_order.i)
- (modules/navier_stokes/test/tests/finite_element/ins/coupled-force/steady-function.i)
- (modules/navier_stokes/test/tests/finite_element/ins/energy_source/steady.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized_with_temp_transient.i)
- (modules/navier_stokes/test/tests/finite_element/ins/boussinesq/boussinesq_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized_with_temp.i)
- (modules/navier_stokes/test/tests/finite_element/ins/wall_convection/steady.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized_second_order.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/coupled-force/gravity-through-coupled-force-action.i)
(modules/navier_stokes/test/tests/finite_element/ins/energy_source/steady-var.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[AuxVariables]
[u]
initial_condition = 1
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[temperature][]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[./temperature_advection]
type = INSADEnergyAdvection
variable = temperature
[../]
[./temperature_conduction]
type = ADHeatConduction
variable = temperature
thermal_conductivity = 'k'
[../]
[temperature_source]
type = INSADEnergySource
variable = temperature
source_variable = u
[]
[temperature_supg]
type = INSADEnergySUPG
variable = temperature
velocity = velocity
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[./temperature_hot]
type = DirichletBC
variable = temperature
boundary = 'bottom'
value = 1
[../]
[./temperature_cold]
type = DirichletBC
variable = temperature
boundary = 'top'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temperature
[]
[]
[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 = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
[out]
type = Exodus
hide = 'u'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_stab_jac_test.i)
[GlobalParams]
order = SECOND
integrate_p_by_parts = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1.1
ymin = -1.1
ymax = 1.1
elem_type = QUAD9
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Problem]
coord_type = RZ
[]
[Preconditioning]
[./SMP_PJFNK]
type = SMP
full = true
solve_type = NEWTON
[../]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 1.1
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
order = FIRST
[../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_time]
type = INSADMomentumTimeDerivative
variable = velocity
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[../]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 1
[../]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[outlet]
type = INSADMomentumNoBCBC
variable = velocity
pressure = p
boundary = 'top'
[]
# When the NoBCBC is applied on the outlet boundary then there is nothing
# constraining the pressure. Thus we must pin the pressure somewhere to ensure
# that the problem is not singular. If the below BC is not applied then
# -pc_type svd -pc_svd_monitor reveals a singular value
[p_corner]
type = DirichletBC
boundary = pinned_node
value = 0
variable = p
[]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1.1 1.1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
(modules/fsi/test/tests/2d-finite-strain-steady/thermal-me.i)
# Units: specific_heat_capacity--cp--J/(kg.K); density--rho--kg/(cm^3);
# dynamic_viscosity--mu--kg/(cm.s); thermal_conductivity--k--W/(cm.K);
# pressure--kg/(cm.s^2); force--kg.cm/s^2
outlet_pressure = 0
inlet_velocity = 150 # cm/s
ini_temp = 593 # K
heat_transfer_coefficient = 9 # W/(cm2.K)
g = -981 # cm/s2
alpha_fluid = 2e-4 # thermal expansion coefficient of fluid used in INSADBoussinesqBodyForce
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
file = '2layers_2d_midline.msh'
[]
[Variables]
[velocity]
family = LAGRANGE_VEC
order = FIRST
block = 'fluid'
[]
[p]
family = LAGRANGE
order = FIRST
block = 'fluid'
[]
[Tf]
family = LAGRANGE
order = FIRST
block = 'fluid'
[]
[Ts]
family = LAGRANGE
order = FIRST
block = 'solid'
[]
[disp_x]
family = LAGRANGE
order = FIRST
block = 'solid fluid'
[]
[disp_y]
family = LAGRANGE
order = FIRST
block = 'solid fluid'
[]
[]
[AuxVariables]
[heat_source]
family = MONOMIAL
order = FIRST
block = 'solid'
[]
[]
[ICs]
[initial_velocity]
type = VectorConstantIC
variable = velocity
x_value = 0
y_value = ${inlet_velocity}
z_value = 0
[]
[initial_p]
type = FunctionIC
variable = p
function = ini_p
[]
[initial_Tf]
type = ConstantIC
variable = Tf
value = ${ini_temp}
[]
[initial_Ts]
type = ConstantIC
variable = Ts
value = ${ini_temp}
[]
[]
[Kernels]
[fluid_mass]
type = INSADMass
variable = p
use_displaced_mesh = true
[]
[fluid_mass_pspg]
type = INSADMassPSPG
variable = p
use_displaced_mesh = true
[]
[fluid_momentum_time]
type = INSADMomentumTimeDerivative
variable = velocity
use_displaced_mesh = true
[]
[fluid_momentum_convection]
type = INSADMomentumAdvection
variable = velocity
use_displaced_mesh = true
[]
[fluid_momentum_viscous]
type = INSADMomentumViscous
variable = velocity
use_displaced_mesh = true
[]
[fluid_momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
use_displaced_mesh = true
[]
[fluid_momentum_gravity]
type = INSADGravityForce
variable = velocity
gravity = '0 ${g} 0'
use_displaced_mesh = true
[]
[fluid_momentum_buoyancy]
type = INSADBoussinesqBodyForce
variable = velocity
gravity = '0 ${g} 0'
alpha_name = 'alpha_fluid'
ref_temp = 'T_ref'
temperature = Tf
use_displaced_mesh = true
[]
[fluid_momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
use_displaced_mesh = true
[]
[fluid_temperature_time]
type = INSADHeatConductionTimeDerivative
variable = Tf
use_displaced_mesh = true
[]
[fluid_temperature_conduction]
type = ADHeatConduction
variable = Tf
thermal_conductivity = 'k'
use_displaced_mesh = true
[]
[fluid_temperature_advection]
type = INSADEnergyAdvection
variable = Tf
use_displaced_mesh = true
[]
[fluid_temperature_supg]
type = INSADEnergySUPG
variable = Tf
velocity = velocity
use_displaced_mesh = true
[]
[solid_temperature_time]
type = ADHeatConductionTimeDerivative
variable = Ts
density_name = 'rho'
specific_heat = 'cp'
block = 'solid'
use_displaced_mesh = true
[]
[solid_temperature_conduction]
type = ADHeatConduction
variable = Ts
thermal_conductivity = 'k'
block = 'solid'
use_displaced_mesh = true
[]
[heat_source]
type = ADCoupledForce
variable = Ts
v = heat_source
block = 'solid'
use_displaced_mesh = true
[]
[disp_x_smooth]
type = Diffusion
variable = disp_x
block = fluid
[]
[disp_y_smooth]
type = Diffusion
variable = disp_y
block = fluid
[]
[]
[Physics/SolidMechanics/QuasiStatic]
strain = FINITE
material_output_order = FIRST
generate_output = 'vonmises_stress stress_xx stress_yy stress_zz strain_xx strain_yy strain_zz'
[solid]
block = 'solid'
temperature = Ts
automatic_eigenstrain_names = true
[]
[]
[InterfaceKernels]
[convection_heat_transfer]
type = ConjugateHeatTransfer
variable = Tf
T_fluid = Tf
neighbor_var = 'Ts'
boundary = 'solid_wall'
htc = 'htc'
use_displaced_mesh = true
[]
[]
[AuxKernels]
[heat_source_distribution_auxk]
type = FunctionAux
variable = heat_source
function = heat_source_distribution_function
block = 'solid'
use_displaced_mesh = true
execute_on = 'INITIAL TIMESTEP_BEGIN'
[]
[]
[BCs]
[no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'solid_wall'
use_displaced_mesh = true
[]
[inlet_velocity]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'fluid_bottom'
function_y = ${inlet_velocity}
use_displaced_mesh = true
[]
[symmetry]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'fluid_wall'
function_x = 0
set_x_comp = true
set_y_comp = false
set_z_comp = false
use_displaced_mesh = true
[]
[outlet_p]
type = DirichletBC
variable = p
boundary = 'fluid_top'
value = ${outlet_pressure}
use_displaced_mesh = true
[]
[inlet_T]
type = DirichletBC
variable = Tf
boundary = 'fluid_bottom'
value = ${ini_temp}
use_displaced_mesh = true
[]
[pin1_y]
type = DirichletBC
variable = disp_y
boundary = 'pin1'
value = 0
use_displaced_mesh = true
[]
[pin1_x]
type = DirichletBC
variable = disp_x
boundary = 'pin1'
value = 0
use_displaced_mesh = true
[]
[top_and_bottom_y]
type = DirichletBC
variable = disp_y
boundary = 'solid_bottom solid_top fluid_top fluid_bottom'
value = 0
use_displaced_mesh = true
[]
[left_and_right_x]
type = DirichletBC
variable = disp_x
boundary = 'fluid_wall fluid_bottom'
value = 0
use_displaced_mesh = true
[]
[]
[Materials]
[rho_solid]
type = ADParsedMaterial
property_name = rho
expression = '0.0110876 * pow(9.9672e-1 + 1.179e-5 * Ts - 2.429e-9 * pow(Ts,2) + 1.219e-12 * pow(Ts,3),-3)'
coupled_variables = 'Ts'
block = 'solid'
use_displaced_mesh = true
[]
[cp_solid]
type = ADParsedMaterial
property_name = cp
expression = '0.76 * ((302.27 * pow((548.68 / Ts),2) * exp(548.68 / Ts)) / pow((exp(548.68 / Ts) - 1),2) + 2 * 8.463e-3 * Ts + 8.741e7 * 18531.7 * exp(-18531.7 / Ts) / pow(Ts,2)) + 0.24 * ((322.49 * pow((587.41/Ts),2) * exp(587.41 / Ts)) / pow((exp(587.41 / Ts) - 1),2) + 2 * 1.4679e-2 * Ts)'
coupled_variables = 'Ts'
block = 'solid'
use_displaced_mesh = true
[]
[k_solid]
type = ADParsedMaterial
property_name = k
expression = '1.158/(7.5408 + 17.692 * (Ts / 1000) + 3.6142 * pow((Ts/1000),2)) + 74.105 * pow((Ts / 1000),-2.5) * exp(-16.35 / (Ts / 1000))'
coupled_variables = 'Ts'
block = 'solid'
use_displaced_mesh = true
[]
[rho_fluid]
type = ADParsedMaterial
property_name = rho
expression = '(11096 - 1.3236 * Tf) * 1e-6'
coupled_variables = 'Tf'
block = 'fluid'
use_displaced_mesh = true
[]
[cp_fluid]
type = ADParsedMaterial
property_name = cp
expression = '159 - 2.72e-2 * Tf + 7.12e-6 * pow(Tf,2)'
coupled_variables = 'Tf'
block = 'fluid'
use_displaced_mesh = true
[]
[k_fluid]
type = ADParsedMaterial
property_name = k
expression = '(3.61 + 1.517e-2 * Tf - 1.741e-6 * pow(Tf,2)) * 1e-2'
coupled_variables = 'Tf'
block = 'fluid'
use_displaced_mesh = true
[]
[mu_fluid]
type = ADParsedMaterial
property_name = mu
expression = '4.94e-6 * exp(754.1/Tf)'
coupled_variables = 'Tf'
block = 'fluid'
use_displaced_mesh = true
[]
[buoyancy_thermal_expansion_coefficient_fluid]
type = ADGenericConstantMaterial
prop_names = 'alpha_fluid'
prop_values = '${alpha_fluid}'
block = 'fluid'
use_displaced_mesh = true
[]
[buoyancy_reference_temperature_fluid]
type = GenericConstantMaterial
prop_names = 'T_ref'
prop_values = '${ini_temp}'
block = 'fluid'
use_displaced_mesh = true
[]
[ins_mat_fluid]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = Tf
block = 'fluid'
use_displaced_mesh = true
[]
[htc]
type = ADGenericFunctionMaterial
prop_names = htc
prop_values = htc_function
use_displaced_mesh = true
[]
[elasticity_solid]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 2e7
poissons_ratio = 0.32
block = 'solid'
use_displaced_mesh = true
[]
[thermal_expansion_solid]
type = ComputeThermalExpansionEigenstrain
temperature = Ts
thermal_expansion_coeff = 2e-4
stress_free_temperature = 593
eigenstrain_name = thermal_expansion
block = 'solid'
use_displaced_mesh = true
[]
[stress_solid]
type = ComputeFiniteStrainElasticStress
block = 'solid'
[]
[]
[Functions]
[htc_function]
type = ParsedFunction
expression = ${heat_transfer_coefficient}
[]
[ini_p]
type = ParsedFunction
expression = '0.010302 * 981 * (10 - y)'
[]
[heat_source_distribution_function]
type = ParsedFunction
expression = '300 * sin(pi * y / 10)'
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
solve_type = 'PJFNK'
[]
[]
[Executioner]
type = Transient
end_time = 1e4
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason -ksp_converged_reason -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
line_search = 'none'
nl_max_its = 30
l_max_its = 100
automatic_scaling = true
compute_scaling_once = true
off_diagonals_in_auto_scaling = true
dtmin = 1
nl_abs_tol = 1e-12
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 6
growth_factor = 1.5
dt = 1
[]
[]
[Outputs]
[csv]
type = CSV
file_base = 'thermal-me'
execute_on = 'final'
[]
[]
[Postprocessors]
[average_solid_Ts]
type = ElementAverageValue
variable = Ts
block = 'solid'
use_displaced_mesh = true
[]
[average_fluid_Tf]
type = ElementAverageValue
variable = Tf
block = 'fluid'
use_displaced_mesh = true
[]
[max_solid_Ts]
type = ElementExtremeValue
variable = Ts
value_type = max
block = 'solid'
use_displaced_mesh = true
[]
[max_fluid_Tf]
type = ElementExtremeValue
variable = Tf
value_type = max
block = 'fluid'
use_displaced_mesh = true
[]
[min_solid_Ts]
type = ElementExtremeValue
variable = Ts
value_type = min
block = 'solid'
use_displaced_mesh = true
[]
[min_fluid_Tf]
type = ElementExtremeValue
variable = Tf
value_type = min
block = 'fluid'
use_displaced_mesh = true
[]
[]
[Debug]
show_var_residual_norms = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/coupled-force/gravity-through-coupled-force.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[u]
family = LAGRANGE_VEC
[]
[]
[AuxVariables]
[gravity]
family = LAGRANGE_VEC
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[gravity]
type = VectorConstantIC
x_value = '0'
y_value = '-9.81'
variable = gravity
[]
[]
[Kernels]
inactive = 'momentum_coupled_forces_two_vars momentum_coupled_forces_two_funcs'
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[momentum_coupled_forces_var_and_func]
type = INSADMomentumCoupledForce
variable = velocity
coupled_vector_var = u
vector_function = 'vector_gravity_func'
[]
[momentum_coupled_forces_two_vars]
type = INSADMomentumCoupledForce
variable = velocity
coupled_vector_var = 'u gravity'
[]
[momentum_coupled_forces_two_funcs]
type = INSADMomentumCoupledForce
variable = velocity
vector_function = 'vector_func vector_gravity_func'
[]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[u_diff]
type = VectorDiffusion
variable = u
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left top'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[u_left]
type = VectorFunctionDirichletBC
variable = u
boundary = 'left'
function_x = 1
function_y = 1
[]
[u_right]
type = VectorFunctionDirichletBC
variable = u
boundary = 'right'
function_x = -1
function_y = -1
[]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
[out]
type = Exodus
hide = 'gravity'
[]
[]
[Functions]
[vector_func]
type = ParsedVectorFunction
expression_x = '-2*x + 1'
expression_y = '-2*x + 1'
[]
[vector_gravity_func]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '-9.81'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_traction_steady_stabilized.i)
[GlobalParams]
order = FIRST
integrate_p_by_parts = true
viscous_form = 'traction'
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[AuxVariables]
[vel_x]
[]
[vel_y]
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[../]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 'inlet_func'
[../]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
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]
console = true
[./out]
type = Exodus
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
execute_on = 'timestep_end'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/coupled-force/gravity-object.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[u]
family = LAGRANGE_VEC
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[momentum_coupled_force]
type = INSADMomentumCoupledForce
variable = velocity
coupled_vector_var = u
[]
[gravity]
type = INSADGravityForce
variable = velocity
gravity = '0 -9.81 0'
[]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[u_diff]
type = VectorDiffusion
variable = u
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left top'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[u_left]
type = VectorFunctionDirichletBC
variable = u
boundary = 'left'
function_x = 1
function_y = 1
[]
[u_right]
type = VectorFunctionDirichletBC
variable = u
boundary = 'right'
function_x = -1
function_y = -1
[]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/boussinesq/benchmark/benchmark.i)
rayleigh=1e3
hot_temp=${rayleigh}
temp_ref=${fparse hot_temp / 2.}
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 100
ny = 100
[]
[./bottom_left]
type = ExtraNodesetGenerator
new_boundary = corner
coord = '0 0'
input = gen
[../]
[]
[Preconditioning]
[./Newton_SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
petsc_options_value = 'bjacobi lu NONZERO 200'
[]
[Debug]
show_var_residual_norms = true
[]
[Outputs]
[out]
type = Exodus
[]
[]
[Variables]
[velocity]
family = LAGRANGE_VEC
[]
[p][]
[temp]
initial_condition = 340
scaling = 1e-4
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[BCs]
[./velocity_dirichlet]
type = VectorDirichletBC
boundary = 'left right bottom top'
variable = velocity
# The third entry is to satisfy RealVectorValue
values = '0 0 0'
[../]
# Even though we are integrating by parts, because there are no integrated
# boundary conditions on the velocity p doesn't appear in the system of
# equations. Thus we must pin the pressure somewhere in order to ensure a
# unique solution
[./p_zero]
type = DirichletBC
boundary = corner
variable = p
value = 0
[../]
[./hot]
type = DirichletBC
variable = temp
boundary = left
value = ${hot_temp}
[../]
[./cold]
type = DirichletBC
variable = temp
boundary = right
value = 0
[../]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[]
[./buoyancy]
type = INSADBoussinesqBodyForce
variable = velocity
temperature = temp
gravity = '0 -1 0'
[../]
[./gravity]
type = INSADGravityForce
variable = velocity
gravity = '0 -1 0'
[../]
[supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[temp_advection]
type = INSADEnergyAdvection
variable = temp
[]
[temp_conduction]
type = ADHeatConduction
variable = temp
thermal_conductivity = 'k'
[../]
[temp_supg]
type = INSADEnergySUPG
variable = temp
velocity = velocity
[]
[]
[Materials]
[./ad_const]
type = ADGenericConstantMaterial
# alpha = coefficient of thermal expansion where rho = rho0 -alpha * rho0 * delta T
prop_names = 'mu rho alpha k cp'
prop_values = '1 1 1 1 1'
[../]
[./const]
type = GenericConstantMaterial
prop_names = 'temp_ref'
prop_values = '${temp_ref}'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temp
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized.i)
[GlobalParams]
order = FIRST
integrate_p_by_parts = true
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[AuxVariables]
[vel_x]
[]
[vel_y]
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[../]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 'inlet_func'
[../]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
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]
console = true
[./out]
type = Exodus
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
execute_on = 'timestep_end'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad-rz-displacements.i)
[GlobalParams]
order = FIRST
integrate_p_by_parts = true
use_displaced_mesh = true
[]
[Mesh]
file = '2d_cone.msh'
displacements = 'disp_x disp_y'
[]
[Problem]
coord_type = RZ
[]
[AuxVariables]
[vel_x][]
[vel_y][]
[disp_x]
order = SECOND
[]
[disp_y]
order = SECOND
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[velocity]
family = LAGRANGE_VEC
[]
[p]
[]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[mass]
type = INSADMass
variable = p
[]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 'inlet_func'
[]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[]
[Functions]
[inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
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
[out]
type = Exodus
hide = 'disp_x disp_y'
[]
[]
[Postprocessors]
[flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
execute_on = 'timestep_end'
[]
[flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
execute_on = 'timestep_end'
[]
[]
(test/tests/kernels/ad_transient_diffusion/ad_transient_vector_diffusion.i)
[Mesh]
[./generator]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[../]
[./block1]
type = SubdomainBoundingBoxGenerator
input = generator
bottom_left = '0 0 -1'
top_right = '1 1 1'
block_id = 1
[../]
[./block2]
type = SubdomainBoundingBoxGenerator
input = block1
bottom_left = '0.33 0.33 -1'
top_right = '0.67 0.67 1'
block_id = 2
[../]
[]
[Variables]
[./u]
family = LAGRANGE_VEC
[../]
[]
[ICs]
[./u]
type = VectorConstantIC
variable = u
x_value = 1
y_value = 2
z_value = 3
block = 2
[../]
[]
[Kernels]
[./diff]
type = ADVectorDiffusion
variable = u
[../]
[./time]
type = ADVectorTimeDerivative
variable = u
[../]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.01
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/mixed-transient-steady/mixed.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[temperature]
[InitialCondition]
type = ConstantIC
value = 1.0
[]
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_time]
type = INSADMomentumTimeDerivative
variable = velocity
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[./temperature_advection]
type = INSADEnergyAdvection
variable = temperature
[../]
[./temperature_conduction]
type = ADHeatConduction
variable = temperature
thermal_conductivity = 'k'
[../]
[temperature_supg]
type = INSADEnergySUPG
variable = temperature
velocity = velocity
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[./temperature_hot]
type = DirichletBC
variable = temperature
boundary = 'bottom'
value = 1
[../]
[./temperature_cold]
type = DirichletBC
variable = temperature
boundary = 'top'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temperature
[]
[]
[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
solve_type = 'NEWTON'
# Run for 100+ timesteps to reach steady state.
num_steps = 5
dt = .5
dtmin = .5
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/coupled-force/steady.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[u]
family = LAGRANGE_VEC
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[momentum_coupled_force]
type = INSADMomentumCoupledForce
variable = velocity
coupled_vector_var = u
[]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[u_diff]
type = VectorDiffusion
variable = u
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left top'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[u_left]
type = VectorFunctionDirichletBC
variable = u
boundary = 'left'
function_x = 1
function_y = 1
[]
[u_right]
type = VectorFunctionDirichletBC
variable = u
boundary = 'right'
function_x = -1
function_y = -1
[]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(test/tests/ics/vector_constant_ic/vector_constant_ic.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 2
ny = 2
nz = 2
[]
[Problem]
solve = false
kernel_coverage_check = false
[]
[Variables]
[./A]
family = LAGRANGE_VEC
order = FIRST
[../]
[]
[ICs]
[./A]
type = VectorConstantIC
variable = A
x_value = 2
y_value = 3
z_value = 4
[../]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/transient_vector_diffusion/transient_vector_diffusion.i)
[Mesh]
[./generator]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
[../]
[./block1]
type = SubdomainBoundingBoxGenerator
input = generator
bottom_left = '0 0 -1'
top_right = '1 1 1'
block_id = 1
[../]
[./block2]
type = SubdomainBoundingBoxGenerator
input = block1
bottom_left = '0.33 0.33 -1'
top_right = '0.67 0.67 1'
block_id = 2
[../]
[]
[Variables]
[./u]
family = LAGRANGE_VEC
[../]
[]
[ICs]
[./u]
type = VectorConstantIC
variable = u
x_value = 1
y_value = 2
z_value = 3
block = 2
[../]
[]
[Kernels]
[./diff]
type = VectorDiffusion
variable = u
[../]
[./time]
type = VectorTimeDerivative
variable = u
[../]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.01
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized_second_order.i)
[GlobalParams]
order = SECOND
integrate_p_by_parts = true
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[AuxVariables]
[vel_x]
[]
[vel_y]
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
order = FIRST
[../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[../]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 'inlet_func'
[../]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
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]
console = true
[./out]
type = Exodus
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
execute_on = 'timestep_end'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/coupled-force/steady-function.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[u]
family = LAGRANGE_VEC
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[momentum_coupled_force]
type = INSADMomentumCoupledForce
variable = velocity
vector_function = 'vector_func'
[]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[u_diff]
type = VectorDiffusion
variable = u
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left top'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[u_left]
type = VectorFunctionDirichletBC
variable = u
boundary = 'left'
function_x = 1
function_y = 1
[]
[u_right]
type = VectorFunctionDirichletBC
variable = u
boundary = 'right'
function_x = -1
function_y = -1
[]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
[Functions]
[vector_func]
type = ParsedVectorFunction
expression_x = '-2*x + 1'
expression_y = '-2*x + 1'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/energy_source/steady.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[temperature][]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[./temperature_advection]
type = INSADEnergyAdvection
variable = temperature
[../]
[./temperature_conduction]
type = ADHeatConduction
variable = temperature
thermal_conductivity = 'k'
[../]
[temperature_source]
type = INSADEnergySource
variable = temperature
source_function = 1
[]
[temperature_supg]
type = INSADEnergySUPG
variable = temperature
velocity = velocity
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[./temperature_hot]
type = DirichletBC
variable = temperature
boundary = 'bottom'
value = 1
[../]
[./temperature_cold]
type = DirichletBC
variable = temperature
boundary = 'top'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temperature
[]
[]
[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 = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized_with_temp_transient.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[temperature]
[InitialCondition]
type = ConstantIC
value = 1.0
[]
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_time]
type = INSADMomentumTimeDerivative
variable = velocity
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[./temperature_advection]
type = INSADEnergyAdvection
variable = temperature
[../]
[temperature_time]
type = INSADHeatConductionTimeDerivative
variable = temperature
[]
[./temperature_conduction]
type = ADHeatConduction
variable = temperature
thermal_conductivity = 'k'
[../]
[temperature_supg]
type = INSADEnergySUPG
variable = temperature
velocity = velocity
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[./temperature_hot]
type = DirichletBC
variable = temperature
boundary = 'bottom'
value = 1
[../]
[./temperature_cold]
type = DirichletBC
variable = temperature
boundary = 'top'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temperature
[]
[]
[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
solve_type = 'NEWTON'
# Run for 100+ timesteps to reach steady state.
num_steps = 5
dt = .5
dtmin = .5
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/boussinesq/boussinesq_stabilized.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmax = .05
ymax = .05
nx = 20
ny = 20
elem_type = QUAD9
[]
[./bottom_left]
type = ExtraNodesetGenerator
new_boundary = corner
coord = '0 0'
input = gen
[../]
[]
[Preconditioning]
[./Newton_SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
petsc_options_value = 'bjacobi lu NONZERO 200'
[]
[Debug]
show_var_residual_norms = true
[]
[Outputs]
[out]
type = Exodus
execute_on = 'final'
[]
[]
[Variables]
[velocity]
family = LAGRANGE_VEC
[]
[p][]
[temp]
initial_condition = 340
scaling = 1e-4
[]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[BCs]
[./velocity_dirichlet]
type = VectorDirichletBC
boundary = 'left right bottom top'
variable = velocity
# The third entry is to satisfy RealVectorValue
values = '0 0 0'
[../]
# Even though we are integrating by parts, because there are no integrated
# boundary conditions on the velocity p doesn't appear in the system of
# equations. Thus we must pin the pressure somewhere in order to ensure a
# unique solution
[./p_zero]
type = DirichletBC
boundary = corner
variable = p
value = 0
[../]
[./cold]
type = DirichletBC
variable = temp
boundary = left
value = 300
[../]
[./hot]
type = DirichletBC
variable = temp
boundary = right
value = 400
[../]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[]
[./buoyancy]
type = INSADBoussinesqBodyForce
variable = velocity
temperature = temp
gravity = '0 -9.81 0'
[../]
[./gravity]
type = INSADGravityForce
variable = velocity
gravity = '0 -9.81 0'
[../]
[supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[temp_advection]
type = INSADEnergyAdvection
variable = temp
[]
[temp_conduction]
type = ADHeatConduction
variable = temp
thermal_conductivity = 'k'
[../]
[temp_supg]
type = INSADEnergySUPG
variable = temp
velocity = velocity
[]
[]
[Materials]
[./ad_const]
type = ADGenericConstantMaterial
# alpha = coefficient of thermal expansion where rho = rho0 -alpha * rho0 * delta T
prop_names = 'mu rho alpha k cp'
prop_values = '30.74e-6 .5757 2.9e-3 46.38e-3 1054'
[../]
[./const]
type = GenericConstantMaterial
prop_names = 'temp_ref'
prop_values = '900'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temp
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized_with_temp.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[temperature][]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[./temperature_advection]
type = INSADEnergyAdvection
variable = temperature
[../]
[./temperature_conduction]
type = ADHeatConduction
variable = temperature
thermal_conductivity = 'k'
[../]
[temperature_supg]
type = INSADEnergySUPG
variable = temperature
velocity = velocity
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[./temperature_hot]
type = DirichletBC
variable = temperature
boundary = 'bottom'
value = 1
[../]
[./temperature_cold]
type = DirichletBC
variable = temperature
boundary = 'top'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temperature
[]
[]
[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 = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/wall_convection/steady.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[temperature][]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[./temperature_advection]
type = INSADEnergyAdvection
variable = temperature
[../]
[./temperature_conduction]
type = ADHeatConduction
variable = temperature
thermal_conductivity = 'k'
[../]
[temperature_ambient_convection]
type = INSADEnergyAmbientConvection
variable = temperature
alpha = 1
T_ambient = 0.5
[]
[temperature_supg]
type = INSADEnergySUPG
variable = temperature
velocity = velocity
[]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[./temperature_hot]
type = DirichletBC
variable = temperature
boundary = 'bottom'
value = 1
[../]
[./temperature_cold]
type = DirichletBC
variable = temperature
boundary = 'top'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu cp k'
prop_values = '1 1 1 .01'
[../]
[ins_mat]
type = INSADStabilized3Eqn
velocity = velocity
pressure = p
temperature = temperature
[]
[]
[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 = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized_second_order.i)
[GlobalParams]
order = SECOND
integrate_p_by_parts = false
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[AuxVariables]
[vel_x]
[]
[vel_y]
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
order = FIRST
[../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[../]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[p_corner]
type = DirichletBC
boundary = top_right
value = 0
variable = p
[]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 'inlet_func'
[../]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
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]
console = true
[./out]
type = Exodus
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
execute_on = 'timestep_end'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 64
ny = 64
[]
[./corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[../]
[]
[AuxVariables]
[vel_x]
[]
[vel_y]
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[]
[ICs]
[velocity]
type = VectorConstantIC
x_value = 0
y_value = 0
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[./mass_pspg]
type = INSADMassPSPG
variable = p
[../]
[./momentum_convection]
type = INSADMomentumAdvection
variable = velocity
[../]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
integrate_p_by_parts = true
[../]
[./momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[../]
[]
[BCs]
[./no_slip]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom right left'
[../]
[./lid]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'top'
function_x = 'lid_function'
[../]
[./pressure_pin]
type = DirichletBC
variable = p
boundary = 'pinned_node'
value = 0
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
alpha = .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'
petsc_options_value = 'lu'
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
file_base = lid_driven_stabilized_out
[]
[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/RZ_cone/ad_rz_cone_no_parts_steady_stabilized.i)
[GlobalParams]
order = FIRST
integrate_p_by_parts = false
[]
[Mesh]
file = '2d_cone.msh'
[]
[Problem]
coord_type = RZ
[]
[AuxVariables]
[vel_x]
[]
[vel_y]
[]
[]
[AuxKernels]
[vel_x]
type = VectorVariableComponentAux
variable = vel_x
vector_variable = velocity
component = 'x'
[]
[vel_y]
type = VectorVariableComponentAux
variable = vel_y
vector_variable = velocity
component = 'y'
[]
[]
[Variables]
[./velocity]
family = LAGRANGE_VEC
[../]
[./p]
[../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
[velocity]
type = VectorConstantIC
x_value = 1e-15
y_value = 1e-15
variable = velocity
[]
[]
[Kernels]
[./mass]
type = INSADMass
variable = p
[../]
[mass_pspg]
type = INSADMassPSPG
variable = p
[]
[momentum_advection]
type = INSADMomentumAdvection
variable = velocity
[]
[./momentum_viscous]
type = INSADMomentumViscous
variable = velocity
[../]
[./momentum_pressure]
type = INSADMomentumPressure
variable = velocity
pressure = p
[../]
[momentum_supg]
type = INSADMomentumSUPG
variable = velocity
velocity = velocity
[]
[]
[BCs]
[p_corner]
type = DirichletBC
boundary = top_right
value = 0
variable = p
[]
[inlet]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'bottom'
function_x = 0
function_y = 'inlet_func'
[../]
[wall]
type = VectorFunctionDirichletBC
variable = velocity
boundary = 'right'
function_x = 0
function_y = 0
[]
[axis]
type = ADVectorFunctionDirichletBC
variable = velocity
boundary = 'left'
set_y_comp = false
function_x = 0
[]
[]
[Functions]
[./inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[../]
[]
[Materials]
[./const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[../]
[ins_mat]
type = INSADTauMaterial
velocity = velocity
pressure = p
[]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Steady
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]
console = true
[./out]
type = Exodus
[../]
[]
[Postprocessors]
[./flow_in]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'bottom'
execute_on = 'timestep_end'
[../]
[./flow_out]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
boundary = 'top'
execute_on = 'timestep_end'
[../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/coupled-force/gravity-through-coupled-force-action.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1.0
ymin = 0
ymax = 1.0
nx = 16
ny = 16
[]
[]
[Variables]
[u]
family = LAGRANGE_VEC
[]
[]
[AuxVariables]
[gravity]
family = LAGRANGE_VEC
[]
[]
[ICs]
[gravity]
type = VectorConstantIC
x_value = '0'
y_value = '-9.81'
variable = gravity
[]
[]
[Modules]
[IncompressibleNavierStokes]
equation_type = steady-state
velocity_boundary = 'bottom right top left'
velocity_function = '0 0 0 0 0 0 0 0'
add_standard_velocity_variables_for_ad = false
pressure_pinned_node = 0
density_name = rho
dynamic_viscosity_name = mu
use_ad = true
laplace = true
family = LAGRANGE
order = FIRST
supg = true
pspg = true
has_coupled_force = true
[]
[]
[Kernels]
[u_diff]
type = VectorDiffusion
variable = u
[]
[]
[BCs]
[u_left]
type = VectorFunctionDirichletBC
variable = u
boundary = 'left'
function_x = 1
function_y = 1
[]
[u_right]
type = VectorFunctionDirichletBC
variable = u
boundary = 'right'
function_x = -1
function_y = -1
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '1 1'
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_levels -ksp_gmres_restart'
petsc_options_value = 'asm 6 200'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
[out]
type = Exodus
hide = 'gravity'
[]
[]
[Functions]
[vector_func]
type = ParsedVectorFunction
expression_x = '-2*x + 1'
expression_y = '-2*x + 1'
[]
[vector_gravity_func]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '-9.81'
[]
[]