- componentThe component to computeC++ Type:MooseEnum Controllable:No Description:The component to compute 
- variableThe name of the variable that this object applies toC++ Type:AuxVariableName Unit:(no unit assumed) Controllable:No Description:The name of the variable that this object applies to 
- vector_variableThe variable from which to compute the componentC++ Type:std::vector<VariableName> Unit:(no unit assumed) Controllable:No Description:The variable from which to compute the component 
VectorVariableComponentAux
The VectorVariableComponentAux class takes a vector variable, specified through the vector_variable parameter, and generates an auxiliary variable corresponding to one of the vector variable's components; the component is specified through the component parameter. This object is only meant to be used with LAGRANGE_VEC vector variables, and hence the auxiliary variable should be of type LAGRANGE.
Creates a field consisting of one component of a coupled vector variable.
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 
- boundaryThe list of boundaries (ids or names) from the mesh where this object appliesC++ Type:std::vector<BoundaryName> Controllable:No Description:The list of boundaries (ids or names) from the mesh where this object applies 
- check_boundary_restrictedTrueWhether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a meshDefault:True C++ Type:bool Controllable:No Description:Whether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh 
- execute_onLINEAR TIMESTEP_ENDThe list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html.Default:LINEAR TIMESTEP_END C++ Type:ExecFlagEnum Options:XFEM_MARK, FORWARD, ADJOINT, HOMOGENEOUS_FORWARD, ADJOINT_TIMESTEP_BEGIN, ADJOINT_TIMESTEP_END, NONE, INITIAL, LINEAR, LINEAR_CONVERGENCE, NONLINEAR, NONLINEAR_CONVERGENCE, POSTCHECK, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, MULTIAPP_FIXED_POINT_CONVERGENCE, FINAL, CUSTOM, PRE_DISPLACE Controllable:No Description:The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html. 
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.C++ Type:std::vector<std::string> Controllable:No Description:Adds user-defined labels for accessing object parameters via control logic. 
- enableTrueSet the enabled status of the MooseObject.Default:True C++ Type:bool Controllable:Yes Description:Set the enabled status of the MooseObject. 
- 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/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady.i)
- (modules/combined/test/tests/electromagnetic_joule_heating/fusing_current_through_copper_wire.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_mean_zero_pressure.i)
- (test/tests/auxkernels/functor_elemental_gradient/functor_gradient.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_traction_steady_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized_second_order.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts.i)
- (test/tests/transfers/multiapp_copy_transfer/vector-variable-transfer/sub_L2_LagrangeVec.i)
- (modules/electromagnetics/test/tests/auxkernels/azimuthal_Faradays_law/scalar_azim_magnetic_time_deriv.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_nobcbc.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized_second_order.i)
- (test/tests/vectorpostprocessors/element_value_sampler/mixed_fe_fv_sampler.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized.i)
- (modules/navier_stokes/examples/laser-welding/2d.i)
- (modules/electromagnetics/test/tests/auxkernels/azimuthal_Faradays_law/error_azim_magnetic_time_deriv.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad-rz-displacements.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_nobcbc.i)
- (modules/fsi/test/tests/newmark-beta/test_ALE.i)
- (modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven.i)
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady.i)
[GlobalParams]
  integrate_p_by_parts = true
[]
[Mesh]
  file = '2d_cone.msh'
  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_convection]
    type = INSADMomentumAdvection
    variable = velocity
  [../]
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/combined/test/tests/electromagnetic_joule_heating/fusing_current_through_copper_wire.i)
# This test is a simpified coupled case between the electromagnetic and
# heat transfer modules. While the file microwave_heating.i is a test
# utilizing the method of manufactured solutions, where both real and
# complex components of the electromagnetic properties are provided
# (such that no term is zeroed out), this test involves only the
# real components of the electromagnetic properties. In particular,
# this test supplies the fusing current to a copper wire and simulations
# the spatial and temporal heating profile until the wire reaches its
# melting point. The PDE's of this test file are as follows:
#
#   curl(curl(A)) + j*mu*omega*(sigma*A) = J
#   mag(E) = mag(-j*omega*A) + mag(J/sigma)
#   rho*C*dT/dt - div(k*grad(T)) = Q
#   Q = 0.5*sigma*mag(E)^2
#
# Where:
#   - A is the magnetic vector potential
#   - j is the sqrt(-1)
#   - mu is the permeability of free space
#   - omega is the angular frequency of the system
#   - sigma is the electric conductivity of the wire
#   - J is the supplied DC current
#   - E is the electric field
#   - rho is the density of copper
#   - C is the heat capacity of copper
#   - T is the temperature
#   - k is the thermal conductivity of the wire
#   - Q is the Joule heating
#
# The BCs are as follows:
#
#   curl(n) x curl(A) = 0,  where n is the normal vector
#   q * n = h (T - T_infty), where q is the heat flux,
#                            h is the convective heat transfer coefficient,
#                            and T_infty is the far-field temperature.
[Mesh]
  # Mesh of the copper wire
  [fmg]
    type = FileMeshGenerator
    file = copper_wire.msh
  []
[]
[Variables]
  # The real and complex components of the magnetic vector
  # potential in the frequency domain
  [A_real]
    family = NEDELEC_ONE
    order = FIRST
  []
  [A_imag]
    family = NEDELEC_ONE
    order = FIRST
  []
  # The temperature of the air in the copper wire
  [T]
    initial_condition = 293.0 #in K
  []
[]
[Kernels]
  ### Physics to determine the magnetic vector potential propagation ###
  # The propagation of the real component
  [curl_curl_real]
    type = CurlCurlField
    variable = A_real
  []
  # Current induced by the electrical conductivity
  # of the copper wire
  [conduction_real]
    type = ADConductionCurrent
    variable = A_real
    field_imag =  A_imag
    field_real =  A_real
    conductivity_real = electrical_conductivity
    conductivity_imag = 0.0
    ang_freq_real = omega_real
    ang_freq_imag = 0.0
    permeability_real = mu_real
    permeability_imag = 0.0
    component = real
  []
  # Current supplied to the wire
  [source_real]
    type = VectorBodyForce
    variable = A_real
    function = mu_curr_real
  []
  # The propagation of the complex component
  [curl_curl_imag]
    type = CurlCurlField
    variable = A_imag
  []
  # Current induced by the electrical conductivity
  # of the copper wire
  [conduction_imag]
    type = ADConductionCurrent
    variable = A_imag
    field_imag =  A_imag
    field_real =  A_real
    conductivity_real = electrical_conductivity
    conductivity_imag = 0.0
    ang_freq_real = omega_real
    ang_freq_imag = 0.0
    permeability_real = mu_real
    permeability_imag = 0.0
    component = imaginary
  []
  ### Physics to determine the heat transfer ###
  # Heat transfer in the copper wire
  [HeatTdot_in_copper]
    type = ADHeatConductionTimeDerivative
    variable = T
    specific_heat = specific_heat_copper
    density_name = density_copper
  []
  [HeatDiff_in_copper]
    type = ADHeatConduction
    variable = T
    thermal_conductivity = thermal_conductivity_copper
  []
  # Heating due the total current
  [HeatSrc]
    type = ADJouleHeatingSource
    variable = T
    heating_term = 'electric_field_heating'
  []
[]
[AuxVariables]
  # Decomposing the magnetic vector potential
  # for the electric field calculations
  [A_x_real]
    family = MONOMIAL
    order = FIRST
  []
  [A_y_real]
    family = MONOMIAL
    order = FIRST
  []
  [A_x_imag]
    family = MONOMIAL
    order = FIRST
  []
  [A_y_imag]
    family = MONOMIAL
    order = FIRST
  []
  # The electrical conductivity for the electric
  # field calculations
  [elec_cond]
    family = MONOMIAL
    order = FIRST
  []
  # The electric field profile determined from
  # the magnetic vector potential
  [E_real]
    family = NEDELEC_ONE
    order = FIRST
  []
  [E_imag]
    family = NEDELEC_ONE
    order = FIRST
  []
[]
[AuxKernels]
  # Decomposing the magnetic vector potential
  # for the electric field calculations
  [A_x_real]
    type = VectorVariableComponentAux
    variable = A_x_real
    vector_variable = A_real
    component = X
  []
  [A_y_real]
    type = VectorVariableComponentAux
    variable = A_y_real
    vector_variable = A_real
    component = Y
  []
  [A_x_imag]
    type = VectorVariableComponentAux
    variable = A_x_imag
    vector_variable = A_imag
    component = X
  []
  [A_y_imag]
    type = VectorVariableComponentAux
    variable = A_y_imag
    vector_variable = A_imag
    component = Y
  []
  # The electrical conductivity for the electric
  # field calculations
  [cond]
    type = ADMaterialRealAux
    property = electrical_conductivity
    variable = elec_cond
    execute_on = 'INITIAL LINEAR TIMESTEP_END'
  []
  # The magnitude of electric field profile determined
  # from the magnetic vector potential using:
  # abs(E) = abs(-j*omega*A) + abs(supplied current / elec_cond)
  # NOTE: The reason for calculating the magnitude of the electric
  #       field is the heating term is defined as:
  #       Q = 1/2 abs(E)^2 for frequency domain field formulations
  [E_real]
    type = ParsedVectorAux
    coupled_variables = 'A_x_imag A_y_imag elec_cond'
    expression_x = 'abs(2*3.14*60*A_x_imag) + abs(60e6/elec_cond)'
    expression_y = 'abs(2*3.14*60*A_y_imag)'
    variable = E_real
  []
  [E_imag]
    type = ParsedVectorAux
    coupled_variables = 'A_x_real A_y_real'
    expression_x = 'abs(-2*3.14*60*A_x_real)'
    expression_y = 'abs(-2*3.14*60*A_y_real)'
    variable = E_imag
  []
[]
[Functions]
  # The supplied current density to the wire
  # where only the real x-component is considered
  [curr_real_x]
    type = ParsedFunction
    expression = '60e6' # Units in A/m^2, equivalent to 1178 A in a 5mm diameter wire
  []
  # Permeability of free space
  [mu_real_func]
    type = ParsedFunction
    expression = '4*pi*1e-7' # Units in N/A^2
  []
  # The angular drive frequency of the system
  [omega_real_func]
    type = ParsedFunction
    expression = '2*pi*60' # Units in rad/s
  []
  # The angular frequency time permeability of free space
  [omegaMu]
    type = ParsedFunction
    symbol_names = 'omega mu'
    symbol_values = 'omega_real_func mu_real_func'
    expression = 'omega*mu'
  []
  # The supplied current density time permeability of free space
  [mu_curr_real]
    type = ParsedVectorFunction
    symbol_names = 'current_mag mu'
    symbol_values = 'curr_real_x mu_real_func'
    expression_x = 'mu * current_mag'
  []
[]
[BCs]
  ### Temperature boundary conditions ###
  # Convective heat flux BC with copper wire
  # exposed to air
  [surface]
    type = ADConvectiveHeatFluxBC
    variable = T
    boundary = walls
    T_infinity = 293
    heat_transfer_coefficient = 10
  []
  ### Magnetic vector potential boundary conditions ###
  # No defined boundary conditions represents
  # zero curl conditions at the boundaries, such that:
  # A x n = 0
[]
[Materials]
  [k]
    type = ADGenericConstantMaterial
    prop_names = 'thermal_conductivity_copper'
    prop_values = '397.48' #in W/(m K)
  []
  [cp]
    type = ADGenericConstantMaterial
    prop_names = 'specific_heat_copper'
    prop_values = '385.0' #in J/(kg K)
  []
  [rho]
    type = ADGenericConstantMaterial
    prop_names = 'density_copper'
    prop_values = '8920.0' #in kg/(m^3)
  []
  # Electrical conductivity (copper is default material)
  [sigma]
    type = ADElectricalConductivity
    temperature = T
    block = copper
  []
  # Material that supplies the correct Joule heating formulation
  [ElectromagneticMaterial]
    type = ElectromagneticHeatingMaterial
    electric_field = E_real
    complex_electric_field = E_imag
    electric_field_heating_name = electric_field_heating
    electrical_conductivity = electrical_conductivity
    formulation = FREQUENCY
    solver = ELECTROMAGNETIC
    block = copper
  []
  # Coefficient for wave propagation
  [mu_real]
    type = ADGenericFunctionMaterial
    prop_names = mu_real
    prop_values = mu_real_func
  []
  [omega_real]
    type = ADGenericFunctionMaterial
    prop_names = omega_real
    prop_values = omega_real_func
  []
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  scheme = bdf2
  solve_type = NEWTON
  line_search = NONE
  petsc_options_iname = '-pc_type'
  petsc_options_value = 'lu'
  dt = 1.0
  # NOTE: Change 'end_time' to 10s to accurately simulate the fusing current
  # end_time = 10
  end_time = 5
  automatic_scaling = true
[]
[Outputs]
  exodus = true
  perf_graph = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized.i)
[GlobalParams]
  order = FIRST
  integrate_p_by_parts = true
[]
[Mesh]
  file = '2d_cone.msh'
  coord_type = RZ
[]
[AuxVariables]
  [vel_x]
  []
  [vel_y]
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
  [../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts.i)
[GlobalParams]
  integrate_p_by_parts = false
  viscous_form = 'traction'
[]
[Mesh]
  file = '2d_cone.msh'
  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]
  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
    pressure = 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
    pressure = 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
    expression = '-4 * x^2 + 1'
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    outputs = 'console'    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    outputs = 'console'    execute_on = 'timestep_end'
  [../]
[]
(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
[]
(test/tests/auxkernels/functor_elemental_gradient/functor_gradient.i)
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 10
  ny = 10
  xmin = 0
  xmax = 3.141
  ymin = 0
  ymax = 3.141
[]
[Variables]
  [u]
  []
  [v]
  []
  [w]
    type = MooseVariableFVReal
  []
[]
[ICs]
  [u_ic]
    type = FunctionIC
    variable = 'u'
    function = parsed_function
  []
  [v_ic]
    type = FunctionIC
    variable = 'v'
    function = 'x'
  []
  [w_ic]
    type = FunctionIC
    variable = 'w'
    function = 'x + y'
  []
[]
[Functions]
  [parsed_function]
    type = ParsedFunction
    value = 'sin(x)-cos(y/2)'
  []
  [parsed_grad_function]
    type = ParsedVectorFunction
    expression_x = 'cos(x)'
    expression_y = 'sin(y/2)/2'
  []
  [parsed_gradx_function]
    type = ParsedFunction
    value = 'cos(x)'
  []
[]
[AuxVariables]
  [funcGrad_u]
    order = CONSTANT
    family = MONOMIAL_VEC
  []
  [auxGrad_u]
    order = CONSTANT
    family = MONOMIAL_VEC
  []
  [auxGrad_v]
    order = CONSTANT
    family = MONOMIAL_VEC
  []
  [auxGrad_fv]
    order = CONSTANT
    family = MONOMIAL_VEC
  []
  [auxGrad_function]
    order = CONSTANT
    family = MONOMIAL_VEC
  []
  [funcGrad_u_x]
    order = CONSTANT
    family = MONOMIAL
  []
  [auxGrad_u_x]
    order = CONSTANT
    family = MONOMIAL
  []
  [auxGrad_v_x]
    order = CONSTANT
    family = MONOMIAL
  []
  [auxGrad_fv_x]
    order = CONSTANT
    family = MONOMIAL
  []
  [auxGrad_function_x]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  # Verification
  [vec]
    type = VectorFunctionAux
    variable = funcGrad_u
    function = parsed_grad_function
  []
  # Finite element variables with and without scaling by material
  [grad_u]
    type = ADFunctorElementalGradientAux
    variable = auxGrad_u
    functor = u
  []
  [grad_v]
    type = ADFunctorElementalGradientAux
    variable = auxGrad_v
    functor = v
    factor_matprop = 'trig_material'
  []
  # Finite volume variable
  [grad_w]
    type = ADFunctorElementalGradientAux
    variable = auxGrad_fv
    functor = w
    factor = w
  []
  # Functions
  [grad_function]
    type = FunctorElementalGradientAux
    variable = auxGrad_function
    functor = parsed_gradx_function
  []
  # Output a component, line sampler does not do vector variables
  [funcGrad_u_x]
    type = VectorVariableComponentAux
    variable = funcGrad_u_x
    vector_variable = funcGrad_u
    component = 'x'
  []
  [auxGrad_u_x]
    type = VectorVariableComponentAux
    variable = auxGrad_u_x
    vector_variable = auxGrad_u
    component = 'x'
  []
  [auxGrad_v_x]
    type = VectorVariableComponentAux
    variable = auxGrad_v_x
    vector_variable = auxGrad_v
    component = 'x'
  []
  [funcGrad_fv_x]
    type = VectorVariableComponentAux
    variable = auxGrad_fv_x
    vector_variable = auxGrad_fv
    component = 'x'
  []
  [auxGrad_function_x]
    type = VectorVariableComponentAux
    variable = auxGrad_function_x
    vector_variable = auxGrad_function
    component = 'x'
  []
[]
[Materials]
  [steel]
    type = ADGenericFunctionMaterial
    prop_names = 'trig_material'
    prop_values = 'parsed_gradx_function'
  []
[]
[VectorPostprocessors]
  [results]
    type = LineValueSampler
    start_point = '0 1 0'
    end_point = '3.141 1 0'
    variable = 'funcGrad_u_x auxGrad_u_x auxGrad_v_x auxGrad_fv_x auxGrad_function_x'
    num_points = 20
    sort_by = x
  []
[]
[Problem]
  solve = false
[]
[Executioner]
  type = Steady
[]
[Outputs]
  csv = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_traction_steady_stabilized.i)
[GlobalParams]
  order = FIRST
  integrate_p_by_parts = true
  viscous_form = 'traction'
[]
[Mesh]
  file = '2d_cone.msh'
  coord_type = RZ
[]
[AuxVariables]
  [vel_x]
  []
  [vel_y]
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
  [../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady.i)
[GlobalParams]
  integrate_p_by_parts = false
[]
[Mesh]
  file = '2d_cone.msh'
  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_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
[]
[BCs]
  [p_corner]
    type = DirichletBC
    boundary = top_right
    value = 0
    variable = p
  []
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_stabilized_second_order.i)
[GlobalParams]
  order = SECOND
  integrate_p_by_parts = true
[]
[Mesh]
  file = '2d_cone.msh'
  coord_type = RZ
[]
[AuxVariables]
  [vel_x]
  []
  [vel_y]
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
    order = FIRST
  [../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts.i)
[Mesh]
  file = '2d_cone.msh'
  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
    pressure = 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
    expression = '-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]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    outputs = 'console'    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    outputs = 'console'    execute_on = 'timestep_end'
  [../]
[]
(test/tests/transfers/multiapp_copy_transfer/vector-variable-transfer/sub_L2_LagrangeVec.i)
[Mesh]
  type = GeneratedMesh
  dim = 3
  nx = 5
  ny = 5
  nz = 5
[]
[AuxVariables]
  [received_vector]
    family = LAGRANGE_VEC
    order = FIRST
  []
  [expected_vector_x]
    family = LAGRANGE
    order = FIRST
  []
  [expected_vector_y]
    family = LAGRANGE
    order = FIRST
  []
  [expected_vector_z]
    family = LAGRANGE
    order = FIRST
  []
  [received_vector_x]
    family = LAGRANGE
    order = FIRST
  []
  [received_vector_y]
    family = LAGRANGE
    order = FIRST
  []
  [received_vector_z]
    family = LAGRANGE
    order = FIRST
  []
[]
[ICs]
  # Set the expected components. If everything works, the received vector components should match.
  [set_expected_vector_x]
    type = FunctionIC
    variable = expected_vector_x
    function = "100*x*x"
  []
  [set_expected_vector_y]
    type = FunctionIC
    variable = expected_vector_y
    function = "100*y*y"
  []
  [set_expected_vector_z]
    type = FunctionIC
    variable = expected_vector_z
    function = "100*z*z"
  []
[]
[AuxKernels]
  # Set the components from the received vector.
  [set_received_vector_x]
    type = VectorVariableComponentAux
    vector_variable = received_vector
    variable = received_vector_x
    component = 'x'
    execute_on = timestep_begin
  []
  [set_received_vector_y]
    type = VectorVariableComponentAux
    vector_variable = received_vector
    variable = received_vector_y
    component = 'y'
    execute_on = timestep_begin
  []
  [set_received_vector_z]
    type = VectorVariableComponentAux
    vector_variable = received_vector
    variable = received_vector_z
    component = 'z'
    execute_on = timestep_begin
  []
[]
[Postprocessors]
  [ensure_something_happened]
    type = ElementAverageValue
    variable = received_vector_x
  []
  # Compare the received vector against the expected components.
  [l2_difference_x]
    type = ElementL2Difference
    variable = received_vector_x
    other_variable = expected_vector_x
  []
  [l2_difference_y]
    type = ElementL2Difference
    variable = received_vector_y
    other_variable = expected_vector_y
  []
  [l2_difference_z]
    type = ElementL2Difference
    variable = received_vector_z
    other_variable = expected_vector_z
  []
[]
[Executioner]
  type = Transient
  dt = 1.0
  start_time = 0.0
  end_time = 1.0
[]
[Outputs]
  csv = true
[]
[Problem]
  solve = false
[]
(modules/electromagnetics/test/tests/auxkernels/azimuthal_Faradays_law/scalar_azim_magnetic_time_deriv.i)
# Test for AzimuthMagneticTimeDerivRZ with scalar inputs
# Manufactured solution: E_real = y^2 * x_hat - x^2 * y_hat
#                        E_imag = y^2 * x_hat - x^2 * y_hat
#                        dB_theta_real / dt = -(2*y + 2*x)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 5
    ny = 5
    xmax = 1
    xmin = 0
    ymax = 1
    ymin = -1
    elem_type = QUAD9
  []
  coord_type = RZ
  rz_coord_axis = Y
[]
[Functions]
  #The exact solution for the heated species and electric field real and imag. component
  [exact_real]
    type = ParsedVectorFunction
    expression_x = 'y^2'
    expression_y = '-x^2'
  []
  [exact_imag]
    type = ParsedVectorFunction
    expression_x = 'y^2'
    expression_y = '-x^2'
  []
  #The forcing terms for the heated species and electric field real and imag. component
  [source_real]
    type = ParsedVectorFunction
    expression_x = '-y^2 - cos(pi*y) - 2'
    expression_y = 'x^2 + cos(pi*x) + 4 + 2*y/x'
  []
  [source_imag]
    type = ParsedVectorFunction
    expression_x = '-y^2 + sin(pi*y) - 2'
    expression_y = 'x^2 - sin(pi*x) + 4 + 2*y/x'
  []
  [current_real]
    type = ParsedVectorFunction
    expression_x = 'sin(pi*y)'
    expression_y = '-sin(pi*x)'
  []
  [current_imag]
    type = ParsedVectorFunction
    expression_x = 'cos(pi*y)'
    expression_y = '-cos(pi*x)'
  []
  [azim_dB_dt_func]
    type = ParsedFunction
    expression = '-(2*y + 2*x)'
  []
[]
[Materials]
  [WaveCoeff]
    type = WaveEquationCoefficient
    eps_rel_real = 1.0
    eps_rel_imag = 0.0
    k_real = 1.0
    k_imag = 0.0
    mu_rel_real = 1.0
    mu_rel_imag = 0.0
  []
[]
[Variables]
  [E_real]
    family = NEDELEC_ONE
    order = FIRST
  []
  [E_imag]
    family = NEDELEC_ONE
    order = FIRST
  []
[]
[Kernels]
  [curl_curl_real]
    type = CurlCurlField
    variable = E_real
  []
  [coeff_real]
    type = ADMatWaveReaction
    variable = E_real
    field_real =  E_real
    field_imag =  E_imag
    wave_coef_real = wave_equation_coefficient_real
    wave_coef_imag = wave_equation_coefficient_imaginary
    component = real
  []
  [current_real]
    type = VectorCurrentSource
    variable = E_real
    source_real = current_real
    source_imag = current_imag
    component = real
  []
  [body_force_real]
    type = VectorBodyForce
    variable = E_real
    function = source_real
  []
  [curl_curl_imag]
    type = CurlCurlField
    variable = E_imag
  []
  [coeff_imag]
    type = ADMatWaveReaction
    variable = E_imag
    field_real =  E_real
    field_imag =  E_imag
    wave_coef_real = wave_equation_coefficient_real
    wave_coef_imag = wave_equation_coefficient_imaginary
    component = imaginary
  []
  [current_imag]
    type = VectorCurrentSource
    variable = E_imag
    source_real = current_real
    source_imag = current_imag
    component = imaginary
  []
  [body_force_imag]
    type = VectorBodyForce
    variable = E_imag
    function = source_imag
  []
[]
[AuxVariables]
  [aux_E_real_x]
    family = MONOMIAL
    order = FIRST
  []
  [aux_E_real_y]
    family = MONOMIAL
    order = FIRST
  []
  [azim_dB_dt_term_scalar]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [aux_E_real_x]
    type = VectorVariableComponentAux
    variable = aux_E_real_x
    vector_variable = E_real
    component = X
  []
  [aux_E_real_y]
    type = VectorVariableComponentAux
    variable = aux_E_real_y
    vector_variable = E_real
    component = Y
  []
  [aux_azim_dB_dt_scalar]
    type = AzimuthMagneticTimeDerivRZ
    Efield_X = aux_E_real_x
    Efield_Y = aux_E_real_y
    variable = azim_dB_dt_term_scalar
  []
[]
[BCs]
  [sides_real]
    type = VectorCurlPenaltyDirichletBC
    variable = E_real
    function = exact_real
    penalty = 1e8
    boundary = 'left right top bottom'
  []
  [sides_imag]
    type = VectorCurlPenaltyDirichletBC
    variable = E_imag
    function = exact_imag
    penalty = 1e8
    boundary = 'left right top bottom'
  []
[]
[Postprocessors]
  [error_real]
    type = ElementVectorL2Error
    variable = E_real
    function = exact_real
  []
  [error_imag]
    type = ElementVectorL2Error
    variable = E_imag
    function = exact_imag
  []
  [error_azim_dB_dt_scalar]
    type = ElementL2Error
    variable = azim_dB_dt_term_scalar
    function = azim_dB_dt_func
  []
  [h]
    type = AverageElementSize
  []
  [h_squared]
    type = ParsedPostprocessor
    pp_names = 'h'
    expression = 'h * h'
  []
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type'
  petsc_options_value = 'lu'
  nl_rel_tol = 1e-16
[]
[Outputs]
  exodus = true
  csv = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/lid_driven/ad_lid_driven_stabilized.i)
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 1.0
    ymin = 0
    ymax = 1.0
    nx = 64
    ny = 64
  []
  [./corner_node]
    type = ExtraNodesetGenerator
    new_boundary = 'pinned_node'
    nodes = '0'
    input = gen
  [../]
[]
[AuxVariables]
  [vel_x]
  []
  [vel_y]
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
  [../]
[]
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 0
    y_value = 0
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [./mass_pspg]
    type = INSADMassPSPG
    variable = p
  [../]
  [./momentum_convection]
    type = INSADMomentumAdvection
    variable = velocity
  [../]
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
    integrate_p_by_parts = true
  [../]
  [./momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  [../]
[]
[BCs]
  [./no_slip]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom right left'
  [../]
  [./lid]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'top'
    function_x = 'lid_function'
  [../]
  [./pressure_pin]
    type = DirichletBC
    variable = p
    boundary = 'pinned_node'
    value = 0
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
    alpha = .1
  []
[]
[Functions]
  [./lid_function]
    # We pick a function that is exactly represented in the velocity
    # space so that the Dirichlet conditions are the same regardless
    # of the mesh spacing.
    type = ParsedFunction
    expression = '4*x*(1-x)'
  [../]
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type'
  petsc_options_value = 'lu'
  line_search = 'none'
  nl_rel_tol = 1e-12
  nl_abs_tol = 1e-13
  nl_max_its = 6
  l_tol = 1e-6
  l_max_its = 500
[]
[Outputs]
  exodus = true
  file_base = lid_driven_stabilized_out
[]
[Postprocessors]
  [lin]
    type = NumLinearIterations
  []
  [nl]
    type = NumNonlinearIterations
  []
  [lin_tot]
    type = CumulativeValuePostprocessor
    postprocessor = 'lin'
  []
  [nl_tot]
    type = CumulativeValuePostprocessor
    postprocessor = 'nl'
  []
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_nobcbc.i)
[GlobalParams]
  integrate_p_by_parts = false
[]
[Mesh]
  file = '2d_cone.msh'
  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_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
[]
[BCs]
  [p_corner]
    type = DirichletBC
    boundary = top_right
    value = 0
    variable = p
  []
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
  [outlet]
    type = INSADMomentumNoBCBC
    variable = velocity
    pressure = p
    boundary = 'top'
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized_second_order.i)
[GlobalParams]
  order = SECOND
  integrate_p_by_parts = false
[]
[Mesh]
  file = '2d_cone.msh'
  coord_type = RZ
[]
[AuxVariables]
  [vel_x]
  []
  [vel_y]
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
    order = FIRST
  [../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [p_corner]
    type = DirichletBC
    boundary = top_right
    value = 0
    variable = p
  []
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(test/tests/vectorpostprocessors/element_value_sampler/mixed_fe_fv_sampler.i)
[Mesh]
  [cmg]
    type = CartesianMeshGenerator
    dim = 1
    dx = 3
    ix = 10
  []
[]
[AuxVariables]
  [T]
    type = MooseVariableFVReal
    [FVInitialCondition]
      type = FVFunctionIC
      function = '10 * x*x'
    []
  []
  [T_linear]
    type = MooseLinearVariableFVReal
    initial_condition = 0.15
  []
  [grad_T]
    order = CONSTANT
    family = MONOMIAL_VEC
  []
  [auxGrad_T_x]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  [grad_T_aux]
    type = FunctorElementalGradientAux
    variable = grad_T
    functor = T
  []
  [grad_T_x_aux]
    type = VectorVariableComponentAux
    variable = auxGrad_T_x
    vector_variable = grad_T
    component = 'x'
  []
[]
[VectorPostprocessors]
  [element_value_sampler]
    type = ElementValueSampler
    variable = 'T auxGrad_T_x T_linear'
    sort_by = id
  []
[]
[Outputs]
  csv = true
[]
[Executioner]
  type = Steady
[]
[Problem]
  solve = false
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_no_parts_steady_stabilized.i)
[GlobalParams]
  order = FIRST
  integrate_p_by_parts = false
[]
[Mesh]
  file = '2d_cone.msh'
  coord_type = RZ
[]
[AuxVariables]
  [vel_x]
  []
  [vel_y]
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [./velocity]
    family = LAGRANGE_VEC
  [../]
  [./p]
  [../]
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [./mass]
    type = INSADMass
    variable = p
  [../]
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [p_corner]
    type = DirichletBC
    boundary = top_right
    value = 0
    variable = p
  []
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  [../]
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/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/electromagnetics/test/tests/auxkernels/azimuthal_Faradays_law/error_azim_magnetic_time_deriv.i)
# Test for AzimuthMagneticTimeDerivRZ with scalar inputs
# Manufactured solution: E_real = y^2 * x_hat - x^2 * y_hat
#                        E_imag = y^2 * x_hat - x^2 * y_hat
#                        dB_theta_real / dt = -(2*y + 2*x)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 10
    ny = 10
    xmax = 1
    xmin = 0
    ymax = 1
    ymin = -1
    elem_type = QUAD9
  []
  coord_type = RZ
  rz_coord_axis = Y
[]
[Functions]
  #The exact solution for the heated species and electric field real and imag. component
  [exact_real]
    type = ParsedVectorFunction
    expression_x = 'y^2'
    expression_y = '-x^2'
  []
  [exact_imag]
    type = ParsedVectorFunction
    expression_x = 'y^2'
    expression_y = '-x^2'
  []
  #The forcing terms for the heated species and electric field real and imag. component
  [source_real]
    type = ParsedVectorFunction
    expression_x = '-y^2 - cos(pi*y) - 2'
    expression_y = 'x^2 + cos(pi*x) + 4 + 2*y/x'
  []
  [source_imag]
    type = ParsedVectorFunction
    expression_x = '-y^2 + sin(pi*y) - 2'
    expression_y = 'x^2 - sin(pi*x) + 4 + 2*y/x'
  []
  [current_real]
    type = ParsedVectorFunction
    expression_x = 'sin(pi*y)'
    expression_y = '-sin(pi*x)'
  []
  [current_imag]
    type = ParsedVectorFunction
    expression_x = 'cos(pi*y)'
    expression_y = '-cos(pi*x)'
  []
  [azim_dB_dt_func]
    type = ParsedFunction
    expression = '-(2*y + 2*x)'
  []
[]
[Materials]
  [WaveCoeff]
    type = WaveEquationCoefficient
    eps_rel_real = 1.0
    eps_rel_imag = 0.0
    k_real = 1.0
    k_imag = 0.0
    mu_rel_real = 1.0
    mu_rel_imag = 0.0
  []
[]
[Variables]
  [E_real]
    family = NEDELEC_ONE
    order = FIRST
  []
  [E_imag]
    family = NEDELEC_ONE
    order = FIRST
  []
[]
[Kernels]
  [curl_curl_real]
    type = CurlCurlField
    variable = E_real
  []
  [coeff_real]
    type = ADMatWaveReaction
    variable = E_real
    field_real =  E_real
    field_imag =  E_imag
    wave_coef_real = wave_equation_coefficient_real
    wave_coef_imag = wave_equation_coefficient_imaginary
    component = real
  []
  [current_real]
    type = VectorCurrentSource
    variable = E_real
    source_real = current_real
    source_imag = current_imag
    component = real
  []
  [body_force_real]
    type = VectorBodyForce
    variable = E_real
    function = source_real
  []
  [curl_curl_imag]
    type = CurlCurlField
    variable = E_imag
  []
  [coeff_imag]
    type = ADMatWaveReaction
    variable = E_imag
    field_real =  E_real
    field_imag =  E_imag
    wave_coef_real = wave_equation_coefficient_real
    wave_coef_imag = wave_equation_coefficient_imaginary
    component = imaginary
  []
  [current_imag]
    type = VectorCurrentSource
    variable = E_imag
    source_real = current_real
    source_imag = current_imag
    component = imaginary
  []
  [body_force_imag]
    type = VectorBodyForce
    variable = E_imag
    function = source_imag
  []
[]
[AuxVariables]
  [aux_E_real_x]
    family = MONOMIAL
    order = FIRST
  []
  [aux_E_real_y]
    family = MONOMIAL
    order = FIRST
  []
  [azim_dB_dt_term]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [aux_E_real_x]
    type = VectorVariableComponentAux
    variable = aux_E_real_x
    vector_variable = E_real
    component = X
  []
  [aux_E_real_y]
    type = VectorVariableComponentAux
    variable = aux_E_real_y
    vector_variable = E_real
    component = Y
  []
  [aux_azim_dB_dt]
    type = AzimuthMagneticTimeDerivRZ
    # Efield = E_real
    # Efield_X = aux_E_real_x
    # Efield_Y = aux_E_real_y
    variable = azim_dB_dt_term
  []
[]
[BCs]
  [sides_real]
    type = VectorCurlPenaltyDirichletBC
    variable = E_real
    function = exact_real
    penalty = 1e8
    boundary = 'left right top bottom'
  []
  [sides_imag]
    type = VectorCurlPenaltyDirichletBC
    variable = E_imag
    function = exact_imag
    penalty = 1e8
    boundary = 'left right top bottom'
  []
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type'
  petsc_options_value = 'lu'
  nl_rel_tol = 1e-16
[]
[Outputs]
  exodus = true
  csv = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad-rz-displacements.i)
[GlobalParams]
  order = FIRST
  integrate_p_by_parts = true
  use_displaced_mesh = true
[]
[Mesh]
  file = '2d_cone.msh'
  displacements = 'disp_x disp_y'
  coord_type = RZ
[]
[AuxVariables]
  [vel_x][]
  [vel_y][]
  [disp_x]
    order = SECOND
  []
  [disp_y]
    order = SECOND
  []
[]
[AuxKernels]
  [vel_x]
    type = VectorVariableComponentAux
    variable = vel_x
    vector_variable = velocity
    component = 'x'
  []
  [vel_y]
    type = VectorVariableComponentAux
    variable = vel_y
    vector_variable = velocity
    component = 'y'
  []
[]
[Variables]
  [velocity]
    family = LAGRANGE_VEC
  []
  [p]
  []
[]
# Need to set a non-zero initial condition because we have a velocity norm in
# the denominator for the tau coefficient of the stabilization term
[ICs]
  [velocity]
    type = VectorConstantIC
    x_value = 1e-15
    y_value = 1e-15
    variable = velocity
  []
[]
[Kernels]
  [mass]
    type = INSADMass
    variable = p
  []
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
  []
  [momentum_advection]
    type = INSADMomentumAdvection
    variable = velocity
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = velocity
    velocity = velocity
  []
[]
[BCs]
  [inlet]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'bottom'
    function_x = 0
    function_y = 'inlet_func'
  []
  [wall]
    type = VectorFunctionDirichletBC
    variable = velocity
    boundary = 'right'
    function_x = 0
    function_y = 0
  []
  [axis]
    type = ADVectorFunctionDirichletBC
    variable = velocity
    boundary = 'left'
    set_y_comp = false
    function_x = 0
  []
[]
[Functions]
  [inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  []
  [ins_mat]
    type = INSADTauMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  []
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  csv = true
  [out]
    type = Exodus
    hide = 'disp_x disp_y'
  []
[]
[Postprocessors]
  [flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  []
  [flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  []
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/ad_rz_cone_by_parts_steady_nobcbc.i)
[GlobalParams]
  integrate_p_by_parts = true
[]
[Mesh]
  file = '2d_cone.msh'
  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_convection]
    type = INSADMomentumAdvection
    variable = velocity
  [../]
  [./momentum_viscous]
    type = INSADMomentumViscous
    variable = velocity
  [../]
  [./momentum_pressure]
    type = INSADMomentumPressure
    variable = velocity
    pressure = p
  [../]
[]
[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
  []
  [outlet]
    type = INSADMomentumNoBCBC
    variable = velocity
    pressure = p
    boundary = 'top'
  []
  # When the NoBCBC is applied on the outlet boundary then there is nothing
  # constraining the pressure. Thus we must pin the pressure somewhere to ensure
  # that the problem is not singular. If the below BC is not applied then
  # -pc_type svd -pc_svd_monitor reveals a singular value
  [p_corner]
    type = DirichletBC
    boundary = top_right
    value = 0
    variable = p
  []
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    expression = '-4 * x^2 + 1'
  [../]
[]
[Materials]
  [./const]
    type = ADGenericConstantMaterial
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
  [ins_mat]
    type = INSADMaterial
    velocity = velocity
    pressure = p
  []
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
    solve_type = 'NEWTON'
  [../]
[]
[Executioner]
  type = Steady
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
  petsc_options_value = 'bjacobi  ilu          4'
  nl_rel_tol = 1e-12
  nl_max_its = 6
[]
[Outputs]
  console = true
  [./out]
    type = Exodus
  [../]
[]
[Postprocessors]
  [./flow_in]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'bottom'
    execute_on = 'timestep_end'
  [../]
  [./flow_out]
    type = VolumetricFlowRate
    vel_x = vel_x
    vel_y = vel_y
    boundary = 'top'
    execute_on = 'timestep_end'
  [../]
[]
(modules/fsi/test/tests/newmark-beta/test_ALE.i)
beta = 0.25
gamma = 0.5
eta = 19.63
zeta = 0.000025
youngs_modulus = 1e8
[GlobalParams]
  displacements = 'disp_x disp_y'
[]
[Mesh]
  [file]
    type = FileMeshGenerator
    file = tmesh_HR.msh
  []
  [convert]
    type = ElementOrderConversionGenerator
    input = file
    conversion_type = FIRST_ORDER
  []
  [matrix_side_interface]
    type = SideSetsBetweenSubdomainsGenerator
    input = convert
    new_boundary = interface_matrix_side
    paired_block = 'inclusion'
    primary_block = 'matrix'
  []
[]
[Variables]
  [disp_x]
    scaling = '${fparse 1/youngs_modulus}'
  []
  [disp_y]
    scaling = '${fparse 1/youngs_modulus}'
  []
  [vel]
    family = LAGRANGE_VEC
    block = 'matrix'
  []
  [p]
    block = 'matrix'
  []
  [lambda]
    family = SCALAR
    block = 'matrix'
  []
[]
[AuxVariables]
  [accel_x]
    block = 'inclusion'
  []
  [accel_y]
    block = 'inclusion'
  []
  [vel_x_solid]
    block = 'inclusion'
  []
  [vel_y_solid]
    block = 'inclusion'
  []
  [vel_x_fluid]
    block = 'matrix'
  []
  [vel_y_fluid]
    block = 'matrix'
  []
[]
[AuxKernels]
  [accel_x] # Calculates and stores acceleration at the end of time step
    type = NewmarkAccelAux
    variable = accel_x
    displacement = disp_x
    velocity = vel_x_solid
    beta = ${beta}
    execute_on = timestep_end
    block = 'inclusion'
  []
  [accel_y] # Calculates and stores acceleration at the end of time step
    type = NewmarkAccelAux
    variable = accel_y
    displacement = disp_y
    velocity = vel_y_solid
    beta = ${beta}
    execute_on = timestep_end
    block = 'inclusion'
  []
  [vel_x_solid]
    type = NewmarkVelAux
    variable = vel_x_solid
    acceleration = accel_x
    gamma = ${gamma}
    execute_on = timestep_end
    block = 'inclusion'
  []
  [vel_y_solid]
    type = NewmarkVelAux
    variable = vel_y_solid
    acceleration = accel_y
    gamma = ${gamma}
    execute_on = timestep_end
    block = 'inclusion'
  []
  [vel_x_fluid]
    type = VectorVariableComponentAux
    variable = vel_x_fluid
    vector_variable = vel
    execute_on = timestep_end
    component = 'x'
  []
  [vel_y_fluid]
    type = VectorVariableComponentAux
    variable = vel_y_fluid
    vector_variable = vel
    execute_on = timestep_end
    component = 'y'
  []
[]
[ScalarKernels]
  [mean_zero_pressure_lm]
    type = AverageValueConstraint
    variable = lambda
    pp_name = pressure_integral
    value = 0
  []
[]
[Kernels]
  [mat_disp_x]
    type = MatDiffusion
    variable = disp_x
    block = 'matrix'
    use_displaced_mesh = false
    diffusivity = ${youngs_modulus}
  []
  [mat_disp_y]
    type = MatDiffusion
    variable = disp_y
    block = 'matrix'
    use_displaced_mesh = false
    diffusivity = ${youngs_modulus}
  []
  [mass]
    type = INSADMass
    variable = p
    use_displaced_mesh = true
    block = 'matrix'
  []
  [mass_pspg]
    type = INSADMassPSPG
    variable = p
    use_displaced_mesh = true
    block = 'matrix'
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative
    variable = vel
    block = 'matrix'
  []
  [momentum_convection]
    type = INSADMomentumAdvection
    variable = vel
    block = 'matrix'
  []
  [momentum_viscous]
    type = INSADMomentumViscous
    variable = vel
    use_displaced_mesh = true
    block = 'matrix'
  []
  [momentum_pressure]
    type = INSADMomentumPressure
    variable = vel
    pressure = p
    integrate_p_by_parts = true
    use_displaced_mesh = true
    block = 'matrix'
  []
  [momentum_supg]
    type = INSADMomentumSUPG
    variable = vel
    material_velocity = relative_velocity
    block = 'matrix'
    use_displaced_mesh = true
  []
  [momentum_mesh_advection]
    type = INSADMomentumMeshAdvection
    variable = vel
    disp_x = 'disp_x'
    disp_y = 'disp_y'
    use_displaced_mesh = true
    block = 'matrix'
  []
  [mean_zero_pressure]
    type = ScalarLagrangeMultiplier
    variable = p
    lambda = lambda
    block = 'matrix'
  []
  # zeta*K*vel + K * disp
  [dynamic_stress_x]
    type = DynamicStressDivergenceTensors
    block = inclusion
    component = 0
    variable = disp_x
    zeta = ${zeta}
  []
  [dynamic_stress_y]
    type = DynamicStressDivergenceTensors
    block = inclusion
    component = 1
    variable = disp_y
    zeta = ${zeta}
  []
  # M*accel + eta*M*vel
  [inertia_x]
    type = InertialForce
    variable = disp_x
    velocity = vel_x_solid
    acceleration = accel_x
    beta = ${beta} # Newmark time integration
    gamma = ${gamma} # Newmark time integration
    eta = ${eta}
    block = 'inclusion'
  []
  [inertia_y]
    type = InertialForce
    variable = disp_y
    velocity = vel_y_solid
    acceleration = accel_y
    beta = ${beta}
    gamma = ${gamma}
    eta = ${eta}
    block = 'inclusion'
  []
[]
[InterfaceKernels]
  [penalty]
    type = ADPenaltyVelocityContinuityNewmarkBeta
    variable = vel
    fluid_velocity = vel
    displacements = 'disp_x disp_y'
    solid_velocities = 'vel_x_solid vel_y_solid'
    solid_accelerations = 'accel_x accel_y'
    boundary = 'interface_matrix_side'
    penalty = ${youngs_modulus}
    beta = ${beta}
    gamma = ${gamma}
  []
[]
[Materials]
  [viscous_mat]
    type = ADGenericConstantMaterial
    block = 'matrix'
    prop_names = 'rho mu'
    prop_values = '1  1'
  []
  [ins_mat]
    type = INSADTauMaterial
    velocity = vel
    pressure = p
    block = 'matrix'
  []
  [elasticity_tensor]
    type = ComputeIsotropicElasticityTensor
    youngs_modulus = ${youngs_modulus}
    poissons_ratio = 0.3
    block = 'inclusion'
  []
  [strain]
    type = ComputeFiniteStrain
    displacements = 'disp_x disp_y'
    block = 'inclusion'
  []
  [small_stress]
    type = ComputeFiniteStrainElasticStress
    block = 'inclusion'
  []
  [density]
    type = GenericConstantMaterial
    block = 'inclusion'
    prop_names = density
    prop_values = 3 # kg/m3
  []
[]
[BCs] # mesh boundaries remain still so I dont think we need to use deformed mesh for vel
  [no_disp_x]
    type = DirichletBC
    variable = disp_x
    boundary = 'bottom top left right'
    value = 0
  []
  [no_disp_y]
    type = DirichletBC
    variable = disp_y
    boundary = 'bottom top left right'
    value = 0
  []
  [shear_top_x]
    type = ADVectorFunctionDirichletBC
    boundary = top
    variable = vel
    function_x = '-0.001'
  []
  [shear_bottom_x]
    type = ADVectorFunctionDirichletBC
    boundary = 'bottom'
    variable = vel
    function_x = '0.001'
  []
  [Periodic]
    [vel]
      variable = vel
      primary = 'left'
      secondary = 'right'
      translation = '1 0 0'
    []
    [x_p]
      variable = p
      primary = 'left'
      secondary = 'right'
      translation = '1 0 0'
    []
  []
[]
[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
  solve_type = 'NEWTON'
  end_time = 100.0
  nl_abs_tol = 1e-12
  [TimeStepper]
    type = IterationAdaptiveDT
    optimal_iterations = 5
    dt = 0.005
    growth_factor = 1.5
    cutback_factor = 0.9
  []
[]
[Postprocessors]
  [pressure_integral]
    type = ElementIntegralVariablePostprocessor
    variable = p
    execute_on = linear
    block = 'matrix'
  []
  [max_vel_y]
    type = ElementExtremeValue
    variable = vel_y_fluid
    block = 'matrix'
    value_type = max
  []
  [min_vel_y]
    type = ElementExtremeValue
    variable = vel_y_fluid
    block = 'matrix'
    value_type = min
  []
[]
[Outputs]
  hide = 'pressure_integral lambda'
  [csv]
    type = CSV
    execute_on = 'final'
  []
[]
(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
[]