- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Description:The name of the variable that this residual object operates on
 
INSADMomentumTimeDerivative
This object adds the term of the incompressible Navier Stokes momentum equation where is the density, is the velocity, and is time.
This class computes the time derivative for the incompressible Navier-Stokes momentum equation.
Input Parameters
- blockThe list of block ids (SubdomainID) that this object will be applied
C++ Type:std::vector<SubdomainName>
Options:
Description:The list of block ids (SubdomainID) that this object will be applied
 - displacementsThe displacements
C++ Type:std::vector<VariableName>
Options:
Description:The displacements
 - temperatureThe temperature on which material properties may depend. If properties do depend on temperature, this variable must be coupled in in order to correctly resize the element matrix
C++ Type:std::vector<VariableName>
Options:
Description:The temperature on which material properties may depend. If properties do depend on temperature, this variable must be coupled in in order to correctly resize the element matrix
 
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
 - diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Options:
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
 - enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
 - implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Options:
Description:Determines whether this object is calculated using an implicit or explicit form
 - save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Options:
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
 - seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Options:
Description:The seed for the master random number generator
 - use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Options:
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
 
Advanced Parameters
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Options:
Description:The extra tags for the matrices this Kernel should fill
 - extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Options:
Description:The extra tags for the vectors this Kernel should fill
 - matrix_tagssystem timeThe tag for the matrices this Kernel should fill
Default:system time
C++ Type:MultiMooseEnum
Options:nontime, system, time
Description:The tag for the matrices this Kernel should fill
 - vector_tagstimeThe tag for the vectors this Kernel should fill
Default:time
C++ Type:MultiMooseEnum
Options:nontime, time
Description:The tag for the vectors this Kernel should fill
 
Tagging Parameters
Input Files
- (modules/navier_stokes/test/tests/ins/block-restriction/two-mats-two-eqn-sets.i)
 - (modules/navier_stokes/test/tests/ins/block-restriction/one-mat-two-eqn-sets.i)
 - (modules/navier_stokes/test/tests/ins/lid_driven/ad_lid_driven.i)
 - (modules/navier_stokes/test/tests/ins/block-restriction/two-mats-one-eqn-set.i)
 - (modules/navier_stokes/test/tests/ins/lid_driven/mixed-transient-steady/mixed.i)
 - (modules/navier_stokes/test/tests/ins/RZ_cone/ad_rz_cone_no_parts.i)
 - (modules/navier_stokes/test/tests/ins/RZ_cone/ad_rz_cone_by_parts.i)
 - (modules/navier_stokes/test/tests/ins/RZ_cone/ad_rz_cone_stab_jac_test.i)
 - (modules/navier_stokes/test/tests/ins/lid_driven/ad_lid_driven_mean_zero_pressure.i)
 - (modules/navier_stokes/test/tests/ins/lid_driven/ad_lid_driven_stabilized_with_temp_transient.i)
 
(modules/navier_stokes/test/tests/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
    p = 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
    p = 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
    value = '4*x*(1-x)'
  [../]
  [./lid_function1]
    type = ParsedFunction
    value = '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/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
    p = 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
    p = 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
    value = '4*x*(1-x)'
  [../]
  [./lid_function1]
    type = ParsedFunction
    value = '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/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
    p = 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
    value = '4*x*(1-x)'
  [../]
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Transient
  # Run for 100+ timesteps to reach steady state.
  num_steps = 5
  dt = .5
  dtmin = .5
  petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'asm      2               ilu          4'
  line_search = 'none'
  nl_rel_tol = 1e-12
  nl_abs_tol = 1e-13
  nl_max_its = 6
  l_tol = 1e-6
  l_max_its = 500
[]
[Outputs]
  file_base = lid_driven_out
  exodus = true
  perf_graph = true
[]
(modules/navier_stokes/test/tests/ins/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
    p = 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
    value = '4*x*(1-x)'
  [../]
  [./lid_function1]
    type = ParsedFunction
    value = '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/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
    p = 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
    value = '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/ins/RZ_cone/ad_rz_cone_no_parts.i)
[GlobalParams]
  integrate_p_by_parts = false
  viscous_form = 'traction'
[]
[Mesh]
  file = '2d_cone.msh'
[]
[Problem]
  coord_type = RZ
[]
[Preconditioning]
  [./SMP_PJFNK]
    type = SMP
    full = true
    solve_type = Newton
  [../]
[]
[Executioner]
  type = Transient
  dt = 0.005
  dtmin = 0.005
  num_steps = 5
  l_max_its = 100
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  csv = true
  console = true
  [./out]
    type = Exodus
  [../]
[]
[AuxVariables]
  [./vel_x]
    # Velocity in radial (r) direction
    family = LAGRANGE
    order = SECOND
  [../]
  [./vel_y]
    # Velocity in axial (z) direction
    family = LAGRANGE
    order = SECOND
  [../]
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [velocity]
    family = LAGRANGE_VEC
    order = SECOND
  []
  [./p]
    family = LAGRANGE
    order = FIRST
  [../]
[]
[BCs]
  [./p_corner]
    # This is required because of the no bcs
    type = DirichletBC
    boundary = top_right
    value = 0
    variable = p
  [../]
  [./velocity_out]
    type = INSADMomentumNoBCBC
    boundary = top
    variable = velocity
    p = p
  [../]
  [./velocity_in]
    type = VectorFunctionDirichletBC
    boundary = bottom
    variable = velocity
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [./wall]
    type = VectorFunctionDirichletBC
    boundary = 'right'
    variable = velocity
    function_x = 0
    function_y = 0
  [../]
  [./axis]
    type = ADVectorFunctionDirichletBC
    boundary = 'left'
    variable = velocity
    set_y_comp = false
    function_x = 0
  [../]
[]
[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
    p = p
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADMaterial
    velocity = velocity
    pressure = p
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    value = '-4 * x^2 + 1'
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/ins/RZ_cone/ad_rz_cone_by_parts.i)
[Mesh]
  file = '2d_cone.msh'
[]
[Problem]
  coord_type = RZ
[]
[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
  [../]
  [./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
    p = p
    integrate_p_by_parts = true
  [../]
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    value = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Transient
  dt = 0.005
  dtmin = 0.005
  num_steps = 5
  l_max_its = 100
  # Note: The Steady executioner can be used for this problem, if you
  # drop the INSMomentumTimeDerivative kernels and use the following
  # direct solver options.
  # petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount -ksp_type'
  # petsc_options_value = 'lu NONZERO 1.e-10 preonly'
  # Block Jacobi works well for this problem, as does "-pc_type asm
  # -pc_asm_overlap 2", but an overlap of 1 does not work for some
  # reason?
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  csv = true
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/ins/RZ_cone/ad_rz_cone_stab_jac_test.i)
[GlobalParams]
  order = SECOND
  integrate_p_by_parts = true
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 2
    xmin = 0
    xmax = 1.1
    ymin = -1.1
    ymax = 1.1
    elem_type = QUAD9
  []
  [./corner_node]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node'
    nodes = '0'
    input = gen
  [../]
[]
[Problem]
  coord_type = RZ
[]
[Preconditioning]
  [./SMP_PJFNK]
    type = SMP
    full = true
    solve_type = NEWTON
  [../]
[]
[Executioner]
  type = Transient
  num_steps = 1
  dt = 1.1
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
    order = FIRST
  [../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative
    variable = velocity
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    p = p
  [../]
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 1
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
  [outlet]
    type = INSADMomentumNoBCBC
    variable = velocity
    p = p
    boundary = 'top'
  []
  # When the NoBCBC is applied on the outlet boundary then there is nothing
  # constraining the pressure. Thus we must pin the pressure somewhere to ensure
  # that the problem is not singular. If the below BC is not applied then
  # -pc_type svd -pc_svd_monitor reveals a singular value
  [p_corner]
    type = DirichletBC
    boundary = pinned_node
    value = 0
    variable = p
  []
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1.1 1.1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
(modules/navier_stokes/test/tests/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
    p = 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
    value = '4*x*(1-x)'
  [../]
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Transient
  # Run for 100+ timesteps to reach steady state.
  num_steps = 5
  dt = .5
  dtmin = .5
  petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -sub_pc_factor_levels -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/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
    p = 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
    value = '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
[]