- variableThe name of the variable that this residual object operates onC++ Type:NonlinearVariableName Unit:(no unit assumed) Controllable:No Description:The name of the variable that this residual object operates on 
INSADEnergyAdvection
This object adds a  term, likely to a heat conduction-convection equation, where  is the density,  is the specific heat capacity,  is the velocity (represented by _U in the code), and  is the temperature gradient (represented by _grad_u in the code). The divergence form of this term is given by
which assuming constant and (a big assumption) can be split into the two terms
Applying the incompressibility constraint  yields the final form used in the INSADEnergyAdvection object
stated above.
This class computes the residual and Jacobian contributions for temperature advection for a divergence free velocity field.
Input Parameters
- blockThe list of blocks (ids or names) that this object will be appliedC++ Type:std::vector<SubdomainName> Controllable:No Description:The list of blocks (ids or names) that this object will be applied 
- displacementsThe displacementsC++ Type:std::vector<VariableName> Unit:(no unit assumed) Controllable:No Description:The displacements 
- matrix_onlyFalseWhether this object is only doing assembly to matrices (no vectors)Default:False C++ Type:bool Controllable:No Description:Whether this object is only doing assembly to matrices (no vectors) 
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contributionC++ 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 fillC++ 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 fillC++ 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 fillDefault: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 fillDefault: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 formDefault: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.) 
- search_methodnearest_node_connected_sidesChoice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes).Default:nearest_node_connected_sides C++ Type:MooseEnum Options:nearest_node_connected_sides, all_proximate_sides Controllable:No Description:Choice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes). 
- seed0The seed for the master random number generatorDefault: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/fsi/test/tests/2d-finite-strain-steady/thermal-me.i)
- (modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q1q1.i)
- (modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q2q1.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_mean_zero_pressure.i)
- (modules/navier_stokes/test/tests/finite_element/ins/block-restriction/two-mats-one-eqn-set.i)
- (modules/navier_stokes/test/tests/finite_element/ins/boussinesq/benchmark/benchmark.i)
- (modules/navier_stokes/test/tests/finite_element/ins/boussinesq/boussinesq_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/block-restriction/two-mats-two-eqn-sets.i)
- (modules/navier_stokes/test/tests/finite_element/ins/wall_convection/steady.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/mixed-transient-steady/mixed.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/lid_driven/ad_lid_driven_stabilized_with_temp_transient.i)
- (modules/navier_stokes/test/tests/finite_element/ins/boussinesq/boussinesq_square_constant_names.i)
- (modules/combined/examples/stochastic/laser_welding_dimred/physics_objects.i)
- (modules/navier_stokes/examples/laser-welding/3d.i)
- (modules/navier_stokes/examples/laser-welding/2d.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/boussinesq/boussinesq_square.i)
- (modules/navier_stokes/test/tests/finite_element/ins/energy_source/steady.i)
- (modules/navier_stokes/test/tests/finite_element/ins/block-restriction/one-mat-two-eqn-sets.i)
- (modules/navier_stokes/test/tests/finite_element/ins/energy_source/steady-var.i)
(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/energy-conservation/q1q1.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    nx = 10
    ny = 10
    dim = 2
  []
  [subdomain]
    type = SubdomainBoundingBoxGenerator
    bottom_left = '0.5 0 0'
    top_right = '1 1 0'
    block_id = 1
    input = gen
  []
  [break_boundary]
    input = subdomain
    type = BreakBoundaryOnSubdomainGenerator
    boundaries = 'bottom top'
  []
  [sideset]
    type = SideSetsBetweenSubdomainsGenerator
    input = break_boundary
    primary_block = '1'
    paired_block = '0'
    new_boundary = 'fluid_left'
  []
  coord_type = RZ
[]
[Variables]
  [T][]
  [velocity]
    family = LAGRANGE_VEC
    block = 1
  []
  [pressure]
    block = 1
  []
[]
[Kernels]
  [mass]
    type = INSADMass
    variable = pressure
    block = 1
  []
  [pspg]
    type = INSADMassPSPG
    variable = pressure
    block = 1
  []
  [momentum_convection]
    type = INSADMomentumAdvection
    variable = velocity
    block = 1
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
    block = 1
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = pressure
    integrate_p_by_parts = true
    block = 1
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
    block = 1
  []
  [temperature_advection]
    type = INSADEnergyAdvection
    variable = T
     block = 1
  []
  [temperature_supg]
    type = INSADEnergySUPG
    variable = T
    velocity = velocity
    block = 1
  []
  [temperature_conduction]
    type = ADHeatConduction
    variable = T
    thermal_conductivity = 'k'
  []
  [heat_source]
    type = BodyForce
    variable = T
    block = 0
    function = 'x + y'
  []
[]
[BCs]
  [velocity_inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    function_y = 1
    boundary = 'bottom_to_1'
  []
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'fluid_left right'
  []
  [convective_heat_transfer]
    type = ConvectiveHeatFluxBC
    variable = T
    T_infinity = 0
    heat_transfer_coefficient = 1
    boundary = 'right'
  []
[]
[Materials]
  [constant]
    type = ADGenericConstantMaterial
    prop_names = 'cp rho k mu'
    prop_values = '1 1   1 1'
  []
  [ins]
    type = INSADStabilized3Eqn
    pressure = pressure
    velocity = velocity
    temperature = T
    block = 1
  []
[]
[Executioner]
  type = Steady
  solve_type = NEWTON
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu       NONZERO'
[]
[Outputs]
  csv = true
[]
[Postprocessors]
  [convective_heat_transfer]
    type = ConvectiveHeatTransferSideIntegral
    T_solid = T
    T_fluid = 0
    htc = 1
    boundary = 'right'
  []
  [advection]
    type = INSADElementIntegralEnergyAdvection
    temperature = T
    velocity = velocity
    cp = cp
    rho = rho
    block = 1
  []
  [source]
    type = FunctionElementIntegral
    function = 'x + y'
    block = 0
  []
  [energy_balance]
    type = ParsedPostprocessor
    expression = 'convective_heat_transfer + advection - source'
    pp_names = 'convective_heat_transfer advection source'
  []
[]
(modules/navier_stokes/test/tests/finite_element/ins/energy-conservation/q2q1.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    nx = 10
    ny = 10
    dim = 2
  []
  [subdomain]
    type = SubdomainBoundingBoxGenerator
    bottom_left = '0.5 0 0'
    top_right = '1 1 0'
    block_id = 1
    input = gen
  []
  [break_boundary]
    input = subdomain
    type = BreakBoundaryOnSubdomainGenerator
    boundaries = 'bottom top'
  []
  [sideset]
    type = SideSetsBetweenSubdomainsGenerator
    input = break_boundary
    primary_block = '1'
    paired_block = '0'
    new_boundary = 'fluid_left'
  []
  coord_type = RZ
  second_order = true
[]
[Variables]
  [T]
    order = SECOND
  []
  [velocity]
    family = LAGRANGE_VEC
    order = SECOND
    block = 1
  []
  [pressure]
    block = 1
  []
[]
[Kernels]
  [mass]
    type = INSADMass
    variable = pressure
    block = 1
  []
  [momentum_convection]
    type = INSADMomentumAdvection
    variable = velocity
    block = 1
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
    block = 1
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = pressure
    integrate_p_by_parts = true
    block = 1
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
    block = 1
  []
  [temperature_advection]
    type = INSADEnergyAdvection
    variable = T
     block = 1
  []
  [temperature_supg]
    type = INSADEnergySUPG
    variable = T
    velocity = velocity
    block = 1
  []
  [temperature_conduction]
    type = ADHeatConduction
    variable = T
    thermal_conductivity = 'k'
  []
  [heat_source]
    type = BodyForce
    variable = T
    block = 0
    function = 'x + y'
  []
[]
[BCs]
  [velocity_inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    function_y = 1
    boundary = 'bottom_to_1'
  []
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'fluid_left right'
  []
  [convective_heat_transfer]
    type = ConvectiveHeatFluxBC
    variable = T
    T_infinity = 0
    heat_transfer_coefficient = 1
    boundary = 'right'
  []
[]
[Materials]
  [constant]
    type = ADGenericConstantMaterial
    prop_names = 'cp rho k mu'
    prop_values = '1 1   1 1'
  []
  [ins]
    type = INSADStabilized3Eqn
    pressure = pressure
    velocity = velocity
    temperature = T
    block = 1
  []
[]
[Executioner]
  type = Steady
  solve_type = NEWTON
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu       NONZERO'
[]
[Outputs]
  csv = true
[]
[Postprocessors]
  [convective_heat_transfer]
    type = ConvectiveHeatTransferSideIntegral
    T_solid = T
    T_fluid = 0
    htc = 1
    boundary = 'right'
  []
  [advection]
    type = INSADElementIntegralEnergyAdvection
    temperature = T
    velocity = velocity
    cp = cp
    rho = rho
    block = 1
  []
  [source]
    type = FunctionElementIntegral
    function = 'x + y'
    block = 0
  []
  [energy_balance]
    type = ParsedPostprocessor
    expression = 'convective_heat_transfer + advection - source'
    pp_names = 'convective_heat_transfer advection source'
  []
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_mean_zero_pressure.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 1.0
    ymin = 0
    ymax = 1.0
    nx = 16
    ny = 16
    elem_type = QUAD9
  []
[]
[AuxVariables]
  [vel_x]
    order = SECOND
  []
  [vel_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]
    order = SECOND
    family = LAGRANGE_VEC
  [../]
  [./T]
    order = SECOND
    [./InitialCondition]
      type = ConstantIC
      value = 1.0
    [../]
  [../]
  [./p]
  [../]
  [./lambda]
    family = SCALAR
    order = FIRST
  [../]
[]
[Kernels]
  [./mass]
    type = INSADMass
    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
  [../]
 [./temperature_time]
   type = INSADHeatConductionTimeDerivative
   variable = T
 [../]
 [./temperature_advection]
   type = INSADEnergyAdvection
   variable = T
 [../]
 [./temperature_conduction]
   type = ADHeatConduction
   variable = T
   thermal_conductivity = 'k'
 [../]
 [./mean_zero_pressure]
    type = ScalarLagrangeMultiplier
    variable = p
    lambda = lambda
  [../]
[]
[ScalarKernels]
  [./mean_zero_pressure_lm]
    type = AverageValueConstraint
    variable = lambda
    pp_name = pressure_integral
    value = 0
  [../]
[]
[BCs]
  [./no_slip]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom right left'
  [../]
  [./lid]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'top'
    function_x = 'lid_function'
  [../]
  [./T_hot]
    type = DirichletBC
    variable = T
    boundary = 'bottom'
    value = 1
  [../]
  [./T_cold]
    type = DirichletBC
    variable = T
    boundary = 'top'
    value = 0
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu cp k'
    prop_values = '1  1  1  .01'
  [../]
  [ins_mat]
    type = INSAD3Eqn
    velocity = velocity
    pressure = p
    temperature = T
  []
[]
[Postprocessors]
  [./pressure_integral]
    type = ElementIntegralVariablePostprocessor
    variable = p
    execute_on = linear
  [../]
[]
[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 -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      2               ilu          4                     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 = true
  perf_graph = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/block-restriction/two-mats-one-eqn-set.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 2
    ymin = 0
    ymax = 1
    nx = 16
    ny = 8
    elem_type = QUAD9
  []
  [./corner_node_0]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node_0'
    coord = '0 0 0'
    input = gen
  [../]
  [./corner_node_1]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node_1'
    coord = '1 0 0'
    input = corner_node_0
  [../]
  [./subdomain1]
    input = corner_node_1
    type = SubdomainBoundingBoxGenerator
    bottom_left = '1 0 0'
    top_right = '2 1 0'
    block_id = 1
  [../]
  [./break_boundary]
    input = subdomain1
    type = BreakBoundaryOnSubdomainGenerator
  [../]
  [./interface0]
    type = SideSetsBetweenSubdomainsGenerator
    input = break_boundary
    primary_block = '0'
    paired_block = '1'
    new_boundary = 'interface0'
  [../]
  [./interface1]
    type = SideSetsBetweenSubdomainsGenerator
    input = interface0
    primary_block = '1'
    paired_block = '0'
    new_boundary = 'interface1'
  [../]
[]
[Variables]
  [velocity0]
    order = SECOND
    family = LAGRANGE_VEC
  []
  [T0]
    order = SECOND
    [InitialCondition]
      type = ConstantIC
      value = 1.0
    []
  []
  [p0]
  []
[]
[Kernels]
  [./mass0]
    type = INSADMass
    variable = p0
  [../]
  [./momentum_time0]
    type = INSADMomentumTimeDerivative
    variable = velocity0
  [../]
  [./momentum_convection0]
    type = INSADMomentumAdvection
    variable = velocity0
  [../]
  [./momentum_viscous0]
    type = INSADMomentumViscous
    variable = velocity0
  [../]
  [./momentum_pressure0]
    type = INSADMomentumPressure
    variable = velocity0
    pressure = p0
    integrate_p_by_parts = true
  [../]
  [./temperature_time0]
    type = INSADHeatConductionTimeDerivative
    variable = T0
  [../]
  [./temperature_advection0]
    type = INSADEnergyAdvection
    variable = T0
  [../]
  [./temperature_conduction0]
    type = ADHeatConduction
    variable = T0
    thermal_conductivity = 'k'
  [../]
[]
[BCs]
  [./no_slip0]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'bottom_to_0 interface0 left'
  [../]
  [./lid0]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'top_to_0'
    function_x = 'lid_function0'
  [../]
  [./T_hot0]
    type = DirichletBC
    variable = T0
    boundary = 'bottom_to_0'
    value = 1
  [../]
  [./T_cold0]
    type = DirichletBC
    variable = T0
    boundary = 'top_to_0'
    value = 0
  [../]
  [./pressure_pin0]
    type = DirichletBC
    variable = p0
    boundary = 'pinned_node_0'
    value = 0
  [../]
  [./no_slip1]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'bottom_to_1 interface1 right'
  [../]
  [./lid1]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'top_to_1'
    function_x = 'lid_function1'
  [../]
  [./T_hot1]
    type = DirichletBC
    variable = T0
    boundary = 'bottom_to_1'
    value = 1
  [../]
  [./T_cold1]
    type = DirichletBC
    variable = T0
    boundary = 'top_to_1'
    value = 0
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu cp k'
    prop_values = '1  1  1  .01'
  [../]
  [ins_mat0]
    type = INSAD3Eqn
    velocity = velocity0
    pressure = p0
    temperature = T0
    block = '0'
  []
  [ins_mat1]
    type = INSAD3Eqn
    velocity = velocity0
    pressure = p0
    temperature = T0
    block = '1'
  []
[]
[Functions]
    # 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.
  [./lid_function0]
    type = ParsedFunction
    expression = '4*x*(1-x)'
  [../]
  [./lid_function1]
    type = ParsedFunction
    expression = '4*(x-1)*(2-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 -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      2               ilu          4                     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 = 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/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/block-restriction/two-mats-two-eqn-sets.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 2
    ymin = 0
    ymax = 1
    nx = 16
    ny = 8
    elem_type = QUAD9
  []
  [./corner_node_0]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node_0'
    coord = '0 0 0'
    input = gen
  [../]
  [./corner_node_1]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node_1'
    coord = '1 0 0'
    input = corner_node_0
  [../]
  [./subdomain1]
    input = corner_node_1
    type = SubdomainBoundingBoxGenerator
    bottom_left = '1 0 0'
    top_right = '2 1 0'
    block_id = 1
  [../]
  [./break_boundary]
    input = subdomain1
    type = BreakBoundaryOnSubdomainGenerator
  [../]
  [./interface0]
    type = SideSetsBetweenSubdomainsGenerator
    input = break_boundary
    primary_block = '0'
    paired_block = '1'
    new_boundary = 'interface0'
  [../]
  [./interface1]
    type = SideSetsBetweenSubdomainsGenerator
    input = interface0
    primary_block = '1'
    paired_block = '0'
    new_boundary = 'interface1'
  [../]
[]
[Variables]
  [velocity0]
    order = SECOND
    family = LAGRANGE_VEC
    block = 0
  []
  [T0]
    order = SECOND
    [InitialCondition]
      type = ConstantIC
      value = 1.0
    []
    block = 0
  []
  [p0]
    block = 0
  []
  [velocity1]
    order = SECOND
    family = LAGRANGE_VEC
    block = 1
  []
  [T1]
    order = SECOND
    [InitialCondition]
      type = ConstantIC
      value = 1.0
    []
    block = 1
  []
  [p1]
    block = 1
  []
[]
[Kernels]
  [./mass0]
    type = INSADMass
    variable = p0
    block = 0
  [../]
  [./momentum_time0]
    type = INSADMomentumTimeDerivative
    variable = velocity0
    block = 0
  [../]
  [./momentum_convection0]
    type = INSADMomentumAdvection
    variable = velocity0
    block = 0
  [../]
  [./momentum_viscous0]
    type = INSADMomentumViscous
    variable = velocity0
    block = 0
  [../]
  [./momentum_pressure0]
    type = INSADMomentumPressure
    variable = velocity0
    pressure = p0
    integrate_p_by_parts = true
    block = 0
  [../]
  [./temperature_time0]
    type = INSADHeatConductionTimeDerivative
    variable = T0
    block = 0
  [../]
  [./temperature_advection0]
    type = INSADEnergyAdvection
    variable = T0
    block = 0
  [../]
  [./temperature_conduction0]
    type = ADHeatConduction
    variable = T0
    thermal_conductivity = 'k'
    block = 0
  [../]
  [./mass1]
    type = INSADMass
    variable = p1
    block = 1
  [../]
  [./momentum_time1]
    type = INSADMomentumTimeDerivative
    variable = velocity1
    block = 1
  [../]
  [./momentum_convection1]
    type = INSADMomentumAdvection
    variable = velocity1
    block = 1
  [../]
  [./momentum_viscous1]
    type = INSADMomentumViscous
    variable = velocity1
    block = 1
  [../]
  [./momentum_pressure1]
    type = INSADMomentumPressure
    variable = velocity1
    pressure = p1
    integrate_p_by_parts = true
    block = 1
  [../]
  [./temperature_time1]
    type = INSADHeatConductionTimeDerivative
    variable = T1
    block = 1
  [../]
  [./temperature_advection1]
    type = INSADEnergyAdvection
    variable = T1
    block = 1
  [../]
  [./temperature_conduction1]
    type = ADHeatConduction
    variable = T1
    thermal_conductivity = 'k'
    block = 1
  [../]
[]
[BCs]
  [./no_slip0]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'bottom_to_0 interface0 left'
  [../]
  [./lid0]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'top_to_0'
    function_x = 'lid_function0'
  [../]
  [./T_hot0]
    type = DirichletBC
    variable = T0
    boundary = 'bottom_to_0'
    value = 1
  [../]
  [./T_cold0]
    type = DirichletBC
    variable = T0
    boundary = 'top_to_0'
    value = 0
  [../]
  [./pressure_pin0]
    type = DirichletBC
    variable = p0
    boundary = 'pinned_node_0'
    value = 0
  [../]
  [./no_slip1]
    type = VectorFunctionDirichletBC
    variable = velocity1
    boundary = 'bottom_to_1 interface1 right'
  [../]
  [./lid1]
    type = VectorFunctionDirichletBC
    variable = velocity1
    boundary = 'top_to_1'
    function_x = 'lid_function1'
  [../]
  [./T_hot1]
    type = DirichletBC
    variable = T1
    boundary = 'bottom_to_1'
    value = 1
  [../]
  [./T_cold1]
    type = DirichletBC
    variable = T1
    boundary = 'top_to_1'
    value = 0
  [../]
  [./pressure_pin1]
    type = DirichletBC
    variable = p1
    boundary = 'pinned_node_1'
    value = 0
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu cp k'
    prop_values = '1  1  1  .01'
  [../]
  [ins_mat0]
    type = INSAD3Eqn
    velocity = velocity0
    pressure = p0
    temperature = T0
    block = 0
  []
  [ins_mat1]
    type = INSAD3Eqn
    velocity = velocity1
    pressure = p1
    temperature = T1
    block = 1
  []
[]
[Functions]
    # 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.
  [./lid_function0]
    type = ParsedFunction
    expression = '4*x*(1-x)'
  [../]
  [./lid_function1]
    type = ParsedFunction
    expression = '4*(x-1)*(2-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 -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      2               ilu          4                     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 = 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/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/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/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_square_constant_names.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
  []
[]
[Executioner]
  type = Steady
  solve_type = NEWTON
  nl_rel_tol = 1e-12
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu       NONZERO'
[]
[Debug]
  show_var_residual_norms = true
[]
[Outputs]
  [out]
    type = Exodus
    execute_on = 'final'
  []
[]
[Variables]
  [velocity]
    family = LAGRANGE_VEC
    order = SECOND
  []
  [p][]
  [temp]
    order = SECOND
    initial_condition = 340
    scaling = 1e-4
  []
[]
[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
  []
  [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
  []
  [temp_advection]
    type = INSADEnergyAdvection
    variable = temp
  []
  [temp_conduction]
    type = ADHeatConduction
    variable = temp
    thermal_conductivity = 'k'
  []
  [buoyancy]
    type = INSADBoussinesqBodyForce
    variable = velocity
    temperature = temp
    gravity = '0 -9.81 0'
    ref_temp = 900
    alpha_name = 2.9e-3
  []
  [gravity]
    type = INSADGravityForce
    variable = velocity
    gravity = '0 -9.81 0'
  []
[]
[Materials]
  [ad_const]
    type = ADGenericConstantMaterial
    prop_names =  'mu        rho    k        cp'
    prop_values = '30.74e-6  .5757  46.38e-3 1054'
  []
  [ins_mat]
    type = INSAD3Eqn
    velocity = velocity
    pressure = p
    temperature = temp
  []
[]
(modules/combined/examples/stochastic/laser_welding_dimred/physics_objects.i)
[ICs]
  [T]
    type = FunctionIC
    variable = T
    function = '(${surfacetemp} - 300) / ${thickness} * y + ${surfacetemp}'
  []
[]
[Kernels]
  [disp_x]
    type = Diffusion
    variable = disp_x
  []
  [disp_y]
    type = Diffusion
    variable = disp_y
  []
  [mass]
    type = INSADMass
    variable = p
    use_displaced_mesh = true
  []
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
    use_displaced_mesh = true
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_mesh_advection]
    type = INSADMomentumMeshAdvection
    variable = vel
    disp_x = disp_x
    disp_y = disp_y
    use_displaced_mesh = true
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = vel
    pressure = p
    integrate_p_by_parts = true
    use_displaced_mesh = true
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = vel
    material_velocity = relative_velocity
    use_displaced_mesh = true
  []
  [temperature_time]
    type = INSADHeatConductionTimeDerivative
    variable = T
    use_displaced_mesh = true
  []
  [temperature_advection]
    type = INSADEnergyAdvection
    variable = T
    use_displaced_mesh = true
  []
  [temperature_mesh_advection]
    type = INSADEnergyMeshAdvection
    variable = T
    disp_x = disp_x
    disp_y = disp_y
    use_displaced_mesh = true
  []
  [temperature_conduction]
    type = ADHeatConduction
    variable = T
    thermal_conductivity = 'k'
    use_displaced_mesh = true
  []
  [temperature_supg]
    type = INSADEnergySUPG
    variable = T
    velocity = vel
    use_displaced_mesh = true
  []
[]
[BCs]
  [x_no_disp]
    type = DirichletBC
    variable = disp_x
    boundary = 'bottom'
    value = 0
  []
  [y_no_disp]
    type = DirichletBC
    variable = disp_y
    boundary = 'bottom'
    value = 0
  []
  [no_slip]
    type = ADVectorFunctionDirichletBC
    variable = vel
    boundary = 'bottom right left'
  []
  [T_cold]
    type = DirichletBC
    variable = T
    boundary = 'bottom'
    value = 300
  []
  [radiation_flux]
    type = FunctionRadiativeBC
    variable = T
    boundary = 'top'
    emissivity_function = '1'
    Tinfinity = 300
    stefan_boltzmann_constant = 5.67e-8
    use_displaced_mesh = true
  []
  [weld_flux]
    type = GaussianEnergyFluxBC
    variable = T
    boundary = 'top'
    P0 = ${power}
    R = ${R}
    x_beam_coord = '${scanning_speed}*t'
    y_beam_coord = '0'
    use_displaced_mesh = true
  []
  [vapor_recoil]
    type = INSADVaporRecoilPressureMomentumFluxBC
    variable = vel
    boundary = 'top'
    use_displaced_mesh = true
  []
  [surface_tension]
    type = INSADSurfaceTensionBC
    variable = vel
    boundary = 'top'
    use_displaced_mesh = true
  []
  [displace_x_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'top'
    variable = 'disp_x'
    velocity = 'vel'
    component = 0
    associated_subdomain = 0
  []
  [displace_y_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'top'
    variable = 'disp_y'
    velocity = 'vel'
    component = 1
    associated_subdomain = 0
  []
  [displace_x_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'top'
    variable = 'disp_x'
    velocity = 'vel'
    component = 0
  []
  [displace_y_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'top'
    variable = 'disp_y'
    velocity = 'vel'
    component = 1
  []
[]
[Materials]
  [ins_mat]
    type = INSADStabilized3Eqn
    velocity = vel
    pressure = p
    temperature = T
    use_displaced_mesh = true
  []
  [steel]
    type = LaserWeld316LStainlessSteel
    temperature = T
    use_displaced_mesh = true
  []
  [steel_boundary]
    type = LaserWeld316LStainlessSteelBoundary
    boundary = 'top'
    temperature = T
    use_displaced_mesh = true
  []
  [const]
    type = GenericConstantMaterial
    prop_names = 'abs sb_constant'
    prop_values = '1 5.67e-8'
    use_displaced_mesh = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)
period=1.25e-3
endtime=${period}
timestep=1.25e-5
surfacetemp=300
sb=5.67e-8
[GlobalParams]
  temperature = T
[]
[Mesh]
  type = GeneratedMesh
  dim = 3
  xmin = -.35e-3
  xmax = 0.35e-3
  ymin = -.35e-3
  ymax = .35e-3
  zmin = -.7e-3
  zmax = 0
  nx = 2
  ny = 2
  nz = 2
  displacements = 'disp_x disp_y disp_z'
  uniform_refine = 2
[]
[Variables]
  [vel]
    family = LAGRANGE_VEC
  []
  [T]
  []
  [p]
  []
  [disp_x]
  []
  [disp_y]
  []
  [disp_z]
  []
[]
[ICs]
  [T]
    type = FunctionIC
    variable = T
    function = '(${surfacetemp} - 300) / .7e-3 * z + ${surfacetemp}'
  []
[]
[Kernels]
  [disp_x]
    type = Diffusion
    variable = disp_x
  []
  [disp_y]
    type = Diffusion
    variable = disp_y
  []
  [disp_z]
    type = Diffusion
    variable = disp_z
  []
  [mass]
    type = INSADMass
    variable = p
    use_displaced_mesh = true
  []
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
    use_displaced_mesh = true
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_mesh_advection]
    type = INSADMomentumMeshAdvection
    variable = vel
    disp_x = disp_x
    disp_y = disp_y
    disp_z = disp_z
    use_displaced_mesh = true
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = vel
    pressure = p
    integrate_p_by_parts = true
    use_displaced_mesh = true
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = vel
    material_velocity = relative_velocity
    use_displaced_mesh = true
  []
  [temperature_time]
    type = INSADHeatConductionTimeDerivative
    variable = T
    use_displaced_mesh = true
  []
  [temperature_advection]
    type = INSADEnergyAdvection
    variable = T
    use_displaced_mesh = true
  []
  [temperature_mesh_advection]
    type = INSADEnergyMeshAdvection
    variable = T
    disp_x = disp_x
    disp_y = disp_y
    disp_z = disp_z
    use_displaced_mesh = true
  []
  [temperature_conduction]
    type = ADHeatConduction
    variable = T
    thermal_conductivity = 'k'
    use_displaced_mesh = true
  []
  [temperature_supg]
    type = INSADEnergySUPG
    variable = T
    velocity = vel
    use_displaced_mesh = true
  []
[]
[BCs]
  [x_no_disp]
    type = DirichletBC
    variable = disp_x
    boundary = 'back'
    value = 0
  []
  [y_no_disp]
    type = DirichletBC
    variable = disp_y
    boundary = 'back'
    value = 0
  []
  [z_no_disp]
    type = DirichletBC
    variable = disp_z
    boundary = 'back'
    value = 0
  []
  [no_slip]
    type = ADVectorFunctionDirichletBC
    variable = vel
    boundary = 'bottom right left top back'
  []
  [T_cold]
    type = DirichletBC
    variable = T
    boundary = 'back'
    value = 300
  []
  [radiation_flux]
    type = FunctionRadiativeBC
    variable = T
    boundary = 'front'
    emissivity_function = '1'
    Tinfinity = 300
    stefan_boltzmann_constant = ${sb}
    use_displaced_mesh = true
  []
  [weld_flux]
    type = GaussianEnergyFluxBC
    variable = T
    boundary = 'front'
    P0 = 159.96989792079225
    R = 1.8257418583505537e-4
    x_beam_coord = '2e-4 * cos(t * 2 * pi / ${period})'
    y_beam_coord = '2e-4 * sin(t * 2 * pi / ${period})'
    z_beam_coord = 0
    use_displaced_mesh = true
  []
  [vapor_recoil]
    type = INSADVaporRecoilPressureMomentumFluxBC
    variable = vel
    boundary = 'front'
    use_displaced_mesh = true
  []
  [surface_tension]
    type = INSADSurfaceTensionBC
    variable = vel
    boundary = 'front'
    use_displaced_mesh = true
  []
  [displace_x_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'front'
    variable = 'disp_x'
    velocity = 'vel'
    component = 0
    associated_subdomain = 0
  []
  [displace_y_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'front'
    variable = 'disp_y'
    velocity = 'vel'
    component = 1
    associated_subdomain = 0
  []
  [displace_z_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'front'
    variable = 'disp_z'
    velocity = 'vel'
    component = 2
    associated_subdomain = 0
  []
  [displace_x_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'front'
    variable = 'disp_x'
    velocity = 'vel'
    component = 0
  []
  [displace_y_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'front'
    variable = 'disp_y'
    velocity = 'vel'
    component = 1
  []
  [displace_z_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'front'
    variable = 'disp_z'
    velocity = 'vel'
    component = 2
  []
[]
[Materials]
  [ins_mat]
    type = INSADStabilized3Eqn
    velocity = vel
    pressure = p
    temperature = T
    use_displaced_mesh = true
  []
  [steel]
    type = AriaLaserWeld304LStainlessSteel
    temperature = T
    beta = 1e7
    use_displaced_mesh = true
  []
  [steel_boundary]
    type = AriaLaserWeld304LStainlessSteelBoundary
    boundary = 'front'
    temperature = T
    use_displaced_mesh = true
  []
  [const]
    type = GenericConstantMaterial
    prop_names = 'abs sb_constant'
    prop_values = '1 ${sb}'
    use_displaced_mesh = true
  []
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type'
    petsc_options_value = 'lu       NONZERO               strumpack'
  []
[]
[Executioner]
  type = Transient
  end_time = ${endtime}
  dtmin = 1e-8
  dtmax = ${timestep}
  petsc_options = '-snes_converged_reason -ksp_converged_reason -options_left'
  solve_type = 'NEWTON'
  line_search = 'none'
  nl_max_its = 12
  l_max_its = 100
  [TimeStepper]
    type = IterationAdaptiveDT
    optimal_iterations = 7
    dt = ${timestep}
    linear_iteration_ratio = 1e6
    growth_factor = 1.5
  []
[]
[Outputs]
  [exodus]
    type = Exodus
    output_material_properties = true
    show_material_properties = 'mu'
  []
  checkpoint = true
  perf_graph = true
[]
[Debug]
  show_var_residual_norms = true
[]
[Adaptivity]
  marker = combo
  max_h_level = 4
  [Indicators]
    [error_T]
      type = GradientJumpIndicator
      variable = T
    []
    [error_dispz]
      type = GradientJumpIndicator
      variable = disp_z
    []
  []
  [Markers]
    [errorfrac_T]
      type = ErrorFractionMarker
      refine = 0.4
      coarsen = 0.2
      indicator = error_T
    []
    [errorfrac_dispz]
      type = ErrorFractionMarker
      refine = 0.4
      coarsen = 0.2
      indicator = error_dispz
    []
    [combo]
      type = ComboMarker
      markers = 'errorfrac_T errorfrac_dispz'
    []
  []
[]
[Postprocessors]
  [num_dofs]
    type = NumDOFs
    system = 'NL'
  []
  [nl]
    type = NumNonlinearIterations
  []
  [tot_nl]
    type = CumulativeValuePostprocessor
    postprocessor = 'nl'
  []
[]
(modules/navier_stokes/examples/laser-welding/2d.i)
endtime=5e-4 # s
timestep=${fparse endtime/100} # s
surfacetemp=300 # K
power=190 # W
R=1.8257418583505537e-4 # m
[Mesh]
  type = GeneratedMesh
  dim = 2
  xmin = -.45e-3 # m
  xmax = 0.45e-3 # m
  ymin = -.9e-4 # m
  ymax = 0
  nx = 25
  ny = 5
  displacements = 'disp_x disp_y'
[]
[GlobalParams]
  temperature = T
[]
[Variables]
  [vel]
    family = LAGRANGE_VEC
  []
  [T]
  []
  [p]
  []
  [disp_x]
  []
  [disp_y]
  []
[]
[AuxVariables]
  [vel_x_aux]
    [InitialCondition]
      type = ConstantIC
      value = 1e-15
    []
  []
  [vel_y_aux]
    [InitialCondition]
      type = ConstantIC
      value = 1e-15
    []
  []
[]
[AuxKernels]
  [vel_x_value]
    type = VectorVariableComponentAux
    variable = vel_x_aux
    vector_variable = vel
    component = x
  []
  [vel_y_value]
    type = VectorVariableComponentAux
    variable = vel_y_aux
    vector_variable = vel
    component = y
  []
[]
[ICs]
  [T]
    type = FunctionIC
    variable = T
    function = '(${surfacetemp} - 300) / .7e-3 * y + ${surfacetemp}'
  []
[]
[Kernels]
  [disp_x]
    type = Diffusion
    variable = disp_x
  []
  [disp_y]
    type = Diffusion
    variable = disp_y
  []
  [mass]
    type = INSADMass
    variable = p
    use_displaced_mesh = true
  []
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
    use_displaced_mesh = true
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_mesh_advection]
    type = INSADMomentumMeshAdvection
    variable = vel
    disp_x = disp_x
    disp_y = disp_y
    use_displaced_mesh = true
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = vel
    use_displaced_mesh = true
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = vel
    pressure = p
    integrate_p_by_parts = true
    use_displaced_mesh = true
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = vel
    material_velocity = relative_velocity
    use_displaced_mesh = true
  []
  [temperature_time]
    type = INSADHeatConductionTimeDerivative
    variable = T
    use_displaced_mesh = true
  []
  [temperature_advection]
    type = INSADEnergyAdvection
    variable = T
    use_displaced_mesh = true
  []
  [temperature_mesh_advection]
    type = INSADEnergyMeshAdvection
    variable = T
    disp_x = disp_x
    disp_y = disp_y
    use_displaced_mesh = true
  []
  [temperature_conduction]
    type = ADHeatConduction
    variable = T
    thermal_conductivity = 'k'
    use_displaced_mesh = true
  []
  [temperature_supg]
    type = INSADEnergySUPG
    variable = T
    velocity = vel
    use_displaced_mesh = true
  []
[]
[BCs]
  [x_no_disp]
    type = DirichletBC
    variable = disp_x
    boundary = 'bottom'
    value = 0
  []
  [y_no_disp]
    type = DirichletBC
    variable = disp_y
    boundary = 'bottom'
    value = 0
  []
  [no_slip]
    type = ADVectorFunctionDirichletBC
    variable = vel
    boundary = 'bottom right left'
  []
  [T_cold]
    type = DirichletBC
    variable = T
    boundary = 'bottom'
    value = 300
  []
  [radiation_flux]
    type = FunctionRadiativeBC
    variable = T
    boundary = 'top'
    emissivity_function = '1'
    Tinfinity = 300
    stefan_boltzmann_constant = 5.67e-8
    use_displaced_mesh = true
  []
  [weld_flux]
    type = GaussianEnergyFluxBC
    variable = T
    boundary = 'top'
    P0 = ${power}
    R = ${R}
    x_beam_coord = '-0.35e-3 +0.7e-3*t/${endtime}'
    y_beam_coord = '0'
    use_displaced_mesh = true
  []
  [vapor_recoil]
    type = INSADVaporRecoilPressureMomentumFluxBC
    variable = vel
    boundary = 'top'
    use_displaced_mesh = true
  []
  [surface_tension]
    type = INSADSurfaceTensionBC
    variable = vel
    boundary = 'top'
    use_displaced_mesh = true
    include_gradient_terms = true
  []
  [displace_x_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'top'
    variable = 'disp_x'
    velocity = 'vel'
    component = 0
    associated_subdomain = 0
  []
  [displace_y_top]
    type = INSADDisplaceBoundaryBC
    boundary = 'top'
    variable = 'disp_y'
    velocity = 'vel'
    component = 1
    associated_subdomain = 0
  []
  [displace_x_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'top'
    variable = 'disp_x'
    velocity = 'vel'
    component = 0
  []
  [displace_y_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC
    boundary = 'top'
    variable = 'disp_y'
    velocity = 'vel'
    component = 1
  []
[]
[Materials]
  [ins_mat]
    type = INSADStabilized3Eqn
    velocity = vel
    pressure = p
    temperature = T
    use_displaced_mesh = true
  []
  [steel]
    type = AriaLaserWeld304LStainlessSteel
    temperature = T
    beta = 1e7
    use_displaced_mesh = true
  []
  [steel_boundary]
    type = AriaLaserWeld304LStainlessSteelBoundary
    boundary = 'top'
    temperature = T
    use_displaced_mesh = true
  []
  [const]
    type = GenericConstantMaterial
    prop_names = 'abs sb_constant'
    prop_values = '1 5.67e-8'
    use_displaced_mesh = true
  []
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type'
    petsc_options_value = 'lu       NONZERO               strumpack'
  []
[]
[Executioner]
  type = Transient
  end_time = ${endtime}
  dtmin = 1e-8
  dtmax = ${timestep}
  petsc_options = '-snes_converged_reason -ksp_converged_reason -options_left'
  solve_type = 'NEWTON'
  line_search = 'none'
  nl_max_its = 12
  l_max_its = 100
  [TimeStepper]
    type = IterationAdaptiveDT
    optimal_iterations = 7
    dt = ${timestep}
    linear_iteration_ratio = 1e6
    growth_factor = 1.5
  []
[]
[Outputs]
  [exodus]
    type = Exodus
    output_material_properties = true
    show_material_properties = 'mu'
  []
  checkpoint = true
  perf_graph = true
[]
[Debug]
  show_var_residual_norms = true
[]
[Adaptivity]
  marker = combo
  max_h_level = 4
  [Indicators]
    [error_T]
      type = GradientJumpIndicator
      variable = T
    []
    [error_dispz]
      type = GradientJumpIndicator
      variable = disp_y
    []
  []
  [Markers]
    [errorfrac_T]
      type = ErrorFractionMarker
      refine = 0.4
      coarsen = 0.2
      indicator = error_T
    []
    [errorfrac_dispz]
      type = ErrorFractionMarker
      refine = 0.4
      coarsen = 0.2
      indicator = error_dispz
    []
    [combo]
      type = ComboMarker
      markers = 'errorfrac_T errorfrac_dispz'
    []
  []
[]
[Postprocessors]
  [num_dofs]
    type = NumDOFs
    system = 'NL'
  []
  [nl]
    type = NumNonlinearIterations
  []
  [tot_nl]
    type = CumulativeValuePostprocessor
    postprocessor = 'nl'
  []
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven.i)
[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
  [../]
[]
[AuxVariables]
  [vel_x]
    order = SECOND
  []
  [vel_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]
    order = SECOND
    family = LAGRANGE_VEC
  [../]
  [./T]
    order = SECOND
    [./InitialCondition]
      type = ConstantIC
      value = 1.0
    [../]
  [../]
  [./p]
  [../]
[]
[Kernels]
  [./mass]
    type = INSADMass
    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
  [../]
 [./temperature_time]
   type = INSADHeatConductionTimeDerivative
   variable = T
 [../]
 [./temperature_advection]
   type = INSADEnergyAdvection
   variable = T
 [../]
 [./temperature_conduction]
   type = ADHeatConduction
   variable = T
   thermal_conductivity = 'k'
 [../]
[]
[BCs]
  [./no_slip]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom right left'
  [../]
  [./lid]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'top'
    function_x = 'lid_function'
  [../]
  [./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 = ADGenericConstantMaterial
    prop_names = 'rho mu cp k'
    prop_values = '1  1  1  .01'
  [../]
  [ins_mat]
    type = INSAD3Eqn
    velocity = velocity
    pressure = p
    temperature = T
  []
[]
[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/boussinesq/boussinesq_square.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
    order = SECOND
  []
  [p][]
  [./temp]
    order = SECOND
    initial_condition = 340
    scaling = 1e-4
  [../]
[]
[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
  [../]
  [./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
  []
  [temp_advection]
    type = INSADEnergyAdvection
    variable = temp
  []
  [temp_conduction]
    type = ADHeatConduction
    variable = temp
    thermal_conductivity = 'k'
  [../]
  [./buoyancy]
    type = INSADBoussinesqBodyForce
    variable = velocity
    temperature = temp
    gravity = '0 -9.81 0'
  [../]
  [./gravity]
    type = INSADGravityForce
    variable = velocity
    gravity = '0 -9.81 0'
  [../]
[]
[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 = INSAD3Eqn
    velocity = velocity
    pressure = p
    temperature = temp
  []
[]
(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/block-restriction/one-mat-two-eqn-sets.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 2
    ymin = 0
    ymax = 1
    nx = 16
    ny = 8
    elem_type = QUAD9
  []
  [./corner_node_0]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node_0'
    coord = '0 0 0'
    input = gen
  [../]
  [./corner_node_1]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node_1'
    coord = '1 0 0'
    input = corner_node_0
  [../]
  [./subdomain1]
    input = corner_node_1
    type = SubdomainBoundingBoxGenerator
    bottom_left = '1 0 0'
    top_right = '2 1 0'
    block_id = 1
  [../]
  [./break_boundary]
    input = subdomain1
    type = BreakBoundaryOnSubdomainGenerator
  [../]
  [./interface0]
    type = SideSetsBetweenSubdomainsGenerator
    input = break_boundary
    primary_block = '0'
    paired_block = '1'
    new_boundary = 'interface0'
  [../]
  [./interface1]
    type = SideSetsBetweenSubdomainsGenerator
    input = interface0
    primary_block = '1'
    paired_block = '0'
    new_boundary = 'interface1'
  [../]
[]
[Variables]
  [velocity0]
    order = SECOND
    family = LAGRANGE_VEC
  []
  [T0]
    order = SECOND
    [InitialCondition]
      type = ConstantIC
      value = 1.0
    []
  []
  [p0]
  []
[]
[Kernels]
  [./mass0]
    type = INSADMass
    variable = p0
    block = 0
  [../]
  [./momentum_time0]
    type = INSADMomentumTimeDerivative
    variable = velocity0
    block = 0
  [../]
  [./momentum_convection0]
    type = INSADMomentumAdvection
    variable = velocity0
    block = 0
  [../]
  [./momentum_viscous0]
    type = INSADMomentumViscous
    variable = velocity0
    block = 0
  [../]
  [./momentum_pressure0]
    type = INSADMomentumPressure
    variable = velocity0
    pressure = p0
    integrate_p_by_parts = true
    block = 0
  [../]
  [./temperature_time0]
    type = INSADHeatConductionTimeDerivative
    variable = T0
    block = 0
  [../]
  [./temperature_advection0]
    type = INSADEnergyAdvection
    variable = T0
    block = 0
  [../]
  [./temperature_conduction0]
    type = ADHeatConduction
    variable = T0
    thermal_conductivity = 'k'
    block = 0
  [../]
  [./mass1]
    type = INSADMass
    variable = p0
    block = 1
  [../]
  [./momentum_time1]
    type = INSADMomentumTimeDerivative
    variable = velocity0
    block = 1
  [../]
  [./momentum_convection1]
    type = INSADMomentumAdvection
    variable = velocity0
    block = 1
  [../]
  [./momentum_viscous1]
    type = INSADMomentumViscous
    variable = velocity0
    block = 1
  [../]
  [./momentum_pressure1]
    type = INSADMomentumPressure
    variable = velocity0
    pressure = p0
    integrate_p_by_parts = true
    block = 1
  [../]
  [./temperature_time1]
    type = INSADHeatConductionTimeDerivative
    variable = T0
    block = 1
  [../]
  [./temperature_advection1]
    type = INSADEnergyAdvection
    variable = T0
    block = 1
  [../]
  [./temperature_conduction1]
    type = ADHeatConduction
    variable = T0
    thermal_conductivity = 'k'
    block = 1
  [../]
[]
[BCs]
  [./no_slip0]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'bottom_to_0 interface0 left'
  [../]
  [./lid0]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'top_to_0'
    function_x = 'lid_function0'
  [../]
  [./T_hot0]
    type = DirichletBC
    variable = T0
    boundary = 'bottom_to_0'
    value = 1
  [../]
  [./T_cold0]
    type = DirichletBC
    variable = T0
    boundary = 'top_to_0'
    value = 0
  [../]
  [./pressure_pin0]
    type = DirichletBC
    variable = p0
    boundary = 'pinned_node_0'
    value = 0
  [../]
  [./no_slip1]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'bottom_to_1 interface1 right'
  [../]
  [./lid1]
    type = VectorFunctionDirichletBC
    variable = velocity0
    boundary = 'top_to_1'
    function_x = 'lid_function1'
  [../]
  [./T_hot1]
    type = DirichletBC
    variable = T0
    boundary = 'bottom_to_1'
    value = 1
  [../]
  [./T_cold1]
    type = DirichletBC
    variable = T0
    boundary = 'top_to_1'
    value = 0
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu cp k'
    prop_values = '1  1  1  .01'
  [../]
  [ins_mat0]
    type = INSAD3Eqn
    velocity = velocity0
    pressure = p0
    temperature = T0
    block = '0 1'
  []
[]
[Functions]
    # 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.
  [./lid_function0]
    type = ParsedFunction
    expression = '4*x*(1-x)'
  [../]
  [./lid_function1]
    type = ParsedFunction
    expression = '4*(x-1)*(2-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 -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      2               ilu          4                     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 = true
[]
(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'
  []
[]