- muThe viscosity
C++ Type:MaterialPropertyName
Description:The viscosity
 - pressureThe pressure variable.
C++ Type:std::vector<VariableName>
Description:The pressure variable.
 - rhoThe value for the density
C++ Type:double
Description:The value for the density
 - uThe velocity in the x direction.
C++ Type:std::vector<VariableName>
Description:The velocity in the x direction.
 - variableThe name of the finite volume variable this kernel applies to
C++ Type:NonlinearVariableName
Description:The name of the finite volume variable this kernel applies to
 - veladvection velocity
C++ Type:MaterialPropertyName
Description:advection velocity
 
INSFVEnergyAdvection
This object adds a  term to a finite volume formulation of a heat transport equation. The user can control what (material) quantity is advected through the advected_quantity parameter. The default value is the name rho_cp_temp which corresponds to a material property name declared by INSFVMaterial.
Input Parameters
- advected_interp_methodupwindThe interpolation to use for the advected quantity. Options are 'upwind' and 'average', with the default being 'upwind'.
Default:upwind
C++ Type:MooseEnum
Options:average, upwind
Description:The interpolation to use for the advected quantity. Options are 'upwind' and 'average', with the default being 'upwind'.
 - advected_quantityrho_cp_tempAn optional parameter for specifying an advected quantity from a material property. If this is not specified, then the advected quantity will simply be the variable that this object is acting on
Default:rho_cp_temp
C++ Type:MaterialPropertyName
Options:
Description:An optional parameter for specifying an advected quantity from a material property. If this is not specified, then the advected quantity will simply be the variable that this object is acting on
 - 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
 - boundaries_to_forceThe set of boundaries to force execution of this FVFluxKernel on.
C++ Type:std::vector<BoundaryName>
Options:
Description:The set of boundaries to force execution of this FVFluxKernel on.
 - boundaries_to_not_forceThe set of boundaries to not force execution of this FVFluxKernel on.
C++ Type:std::vector<BoundaryName>
Options:
Description:The set of boundaries to not force execution of this FVFluxKernel on.
 - force_boundary_executionFalseWhether to force execution of this object on the boundary.
Default:False
C++ Type:bool
Options:
Description:Whether to force execution of this object on the boundary.
 - ghost_layers2The number of layers of elements to ghost.
Default:2
C++ Type:unsigned short
Options:
Description:The number of layers of elements to ghost.
 - use_point_neighborsFalseWhether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
Default:False
C++ Type:bool
Options:
Description:Whether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
 - vThe velocity in the y direction.
C++ Type:std::vector<VariableName>
Options:
Description:The velocity in the y direction.
 - velocity_interp_methodrcThe interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow.
Default:rc
C++ Type:MooseEnum
Options:average, rc
Description:The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow.
 - wThe velocity in the z direction.
C++ Type:std::vector<VariableName>
Options:
Description:The velocity in the z direction.
 
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.
 - 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
 - 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_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Description:The tag for the matrices this Kernel should fill
 - vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Description:The tag for the vectors this Kernel should fill
 
Tagging Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/cylindrical/2d-average-with-temp.i)
 - (modules/navier_stokes/test/tests/postprocessors/conservation_INSFV.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-ambient-convection.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/block_restriction/2d-rc.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/boussinesq/boussinesq.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/lid-driven-with-energy.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/transient-lid-driven-with-energy.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/2d-average-with-temp.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-scalar-transport.i)
 
(modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/cylindrical/2d-average-with-temp.i)
mu=1.1
rho=1.1
k=1.1
cp=1.1
advected_interp_method='average'
velocity_interp_method='average'
velocity='velocity'
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 1
    ymin = 0
    ymax = 1
    nx = 2
    ny = 2
  []
[]
[Problem]
  kernel_coverage_check = false
  fv_bcs_integrity_check = true
  coord_type = 'RZ'
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 1
    two_term_boundary_expansion = false
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
    two_term_boundary_expansion = false
  []
  [pressure]
    type = INSFVPressureVariable
    two_term_boundary_expansion = false
  []
  [temperature]
    type = INSFVEnergyVariable
    two_term_boundary_expansion = false
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = ${velocity}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [mass_forcing]
    type = FVBodyForce
    variable = pressure
    function = forcing_p
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = ${velocity}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [u_forcing]
    type = FVBodyForce
    variable = u
    function = forcing_u
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = ${velocity}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [v_forcing]
    type = FVBodyForce
    variable = v
    function = forcing_v
  []
  [temp_conduction]
    type = FVDiffusion
    coeff = 'k'
    variable = temperature
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    vel = ${velocity}
    variable = temperature
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [temp_forcing]
    type = FVBodyForce
    variable = temperature
    function = forcing_t
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'bottom'
    variable = u
    function = 'exact_u'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'bottom'
    variable = v
    function = 'exact_v'
  []
  [no-slip-wall-u]
    type = INSFVNoSlipWallBC
    boundary = 'right'
    variable = u
    function = 'exact_u'
  []
  [no-slip-wall-v]
    type = INSFVNoSlipWallBC
    boundary = 'right'
    variable = v
    function = 'exact_v'
  []
  [outlet-p]
    type = INSFVOutletPressureBC
    boundary = 'top'
    variable = pressure
    function = 'exact_p'
  []
  [axis-u]
    type = INSFVSymmetryVelocityBC
    boundary = 'left'
    variable = u
    u = u
    v = v
    mu = ${mu}
    momentum_component = x
  []
  [axis-v]
    type = INSFVSymmetryVelocityBC
    boundary = 'left'
    variable = v
    u = u
    v = v
    mu = ${mu}
    momentum_component = y
  []
  [axis-p]
    type = INSFVSymmetryPressureBC
    boundary = 'left'
    variable = pressure
  []
  [axis-inlet-wall-t]
    type = FVFunctionDirichletBC
    boundary = 'left bottom right'
    variable = temperature
    function = 'exact_t'
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'k cp'
    prop_values = '${k} ${cp}'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'temperature'
    rho = ${rho}
  []
[]
[Functions]
[exact_u]
  type = ParsedFunction
  value = 'sin(x*pi)*sin((1/2)*y*pi)'
[]
[exact_rhou]
  type = ParsedFunction
  value = 'rho*sin(x*pi)*sin((1/2)*y*pi)'
  vars = 'rho'
  vals = '${rho}'
[]
[forcing_u]
  type = ParsedFunction
  value = '(1/4)*pi^2*mu*sin(x*pi)*sin((1/2)*y*pi) - pi*sin(x*pi)*cos((1/2)*y*pi) - (-x*pi^2*mu*sin(x*pi)*sin((1/2)*y*pi) + pi*mu*sin((1/2)*y*pi)*cos(x*pi))/x + (2*x*pi*rho*sin(x*pi)*sin((1/2)*y*pi)^2*cos(x*pi) + rho*sin(x*pi)^2*sin((1/2)*y*pi)^2)/x + (-x*pi*rho*sin(x*pi)*sin((1/2)*y*pi)*sin(y*pi)*cos(x*pi) + (1/2)*x*pi*rho*sin(x*pi)*cos(x*pi)*cos((1/2)*y*pi)*cos(y*pi))/x'
  vars = 'mu rho'
  vals = '${mu} ${rho}'
[]
[exact_v]
  type = ParsedFunction
  value = 'cos(x*pi)*cos(y*pi)'
[]
[exact_rhov]
  type = ParsedFunction
  value = 'rho*cos(x*pi)*cos(y*pi)'
  vars = 'rho'
  vals = '${rho}'
[]
[forcing_v]
  type = ParsedFunction
  value = 'pi^2*mu*cos(x*pi)*cos(y*pi) - 2*pi*rho*sin(y*pi)*cos(x*pi)^2*cos(y*pi) - 1/2*pi*sin((1/2)*y*pi)*cos(x*pi) - (-x*pi^2*mu*cos(x*pi)*cos(y*pi) - pi*mu*sin(x*pi)*cos(y*pi))/x + (-x*pi*rho*sin(x*pi)^2*sin((1/2)*y*pi)*cos(y*pi) + x*pi*rho*sin((1/2)*y*pi)*cos(x*pi)^2*cos(y*pi) + rho*sin(x*pi)*sin((1/2)*y*pi)*cos(x*pi)*cos(y*pi))/x'
  vars = 'mu rho'
  vals = '${mu} ${rho}'
[]
[exact_p]
  type = ParsedFunction
  value = 'cos(x*pi)*cos((1/2)*y*pi)'
[]
[forcing_p]
  type = ParsedFunction
  value = '-pi*rho*sin(y*pi)*cos(x*pi) + (x*pi*rho*sin((1/2)*y*pi)*cos(x*pi) + rho*sin(x*pi)*sin((1/2)*y*pi))/x'
  vars = 'rho'
  vals = '${rho}'
[]
[exact_t]
  type = ParsedFunction
  value = 'sin(x*pi)*sin((1/2)*y*pi)'
[]
[forcing_t]
  type = ParsedFunction
  value = '(1/4)*pi^2*k*sin(x*pi)*sin((1/2)*y*pi) - (-x*pi^2*k*sin(x*pi)*sin((1/2)*y*pi) + pi*k*sin((1/2)*y*pi)*cos(x*pi))/x + (2*x*pi*cp*rho*sin(x*pi)*sin((1/2)*y*pi)^2*cos(x*pi) + cp*rho*sin(x*pi)^2*sin((1/2)*y*pi)^2)/x + (-x*pi*cp*rho*sin(x*pi)*sin((1/2)*y*pi)*sin(y*pi)*cos(x*pi) + (1/2)*x*pi*cp*rho*sin(x*pi)*cos(x*pi)*cos((1/2)*y*pi)*cos(y*pi))/x'
  vars = 'k rho cp'
  vals = '${k} ${rho} ${cp}'
[]
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      100                lu           NONZERO'
  line_search = 'none'
[]
[Outputs]
  exodus = true
  csv = true
  [dof]
    type = DOFMap
    execute_on = 'initial'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    outputs = 'console csv'
    execute_on = 'timestep_end'
  []
  [./L2u]
    type = ElementL2Error
    variable = u
    function = exact_u
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./L2v]
    type = ElementL2Error
    variable = v
    function = exact_v
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./L2p]
    variable = pressure
    function = exact_p
    type = ElementL2Error
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./L2t]
    variable = temperature
    function = exact_t
    type = ElementL2Error
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/postprocessors/conservation_INSFV.i)
mu=1.1
rho=1
advected_interp_method='average'
velocity_interp_method='average'
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  inactive = 'mesh internal_boundary_bot internal_boundary_top'
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1'
    dy = '1 1 1'
    ix = '5'
    iy = '5 5 5'
    subdomain_id = '1
                    2
                    3'
  []
  [internal_boundary_bot]
    type = SideSetsBetweenSubdomainsGenerator
    input = mesh
    new_boundary = 'internal_bot'
    primary_block = 1
    paired_block = 2
  []
  [internal_boundary_top]
    type = SideSetsBetweenSubdomainsGenerator
    input = internal_boundary_bot
    new_boundary = 'internal_top'
    primary_block = 2
    paired_block = 3
  []
  [diverging_mesh]
    type = FileMeshGenerator
    file = 'expansion_quad.e'
  []
[]
[Problem]
  fv_bcs_integrity_check = true
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [temperature]
    type = INSFVEnergyVariable
  []
[]
[AuxVariables]
  [advected_density]
    order = CONSTANT
    family = MONOMIAL
    fv = true
    initial_condition = ${rho}
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
    force_boundary_execution = true
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
    force_boundary_execution = true
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    vel = 'velocity'
    variable = temperature
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [temp_source]
    type = FVBodyForce
    variable = temperature
    function = 10
    block = 1
  []
[]
[FVBCs]
  inactive = 'noslip-u noslip-v'
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'bottom'
    variable = u
    function = 0
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'bottom'
    variable = v
    function = 1
  []
  [noslip-u]
    type = INSFVNoSlipWallBC
    boundary = 'right'
    variable = u
    function = 0
  []
  [noslip-v]
    type = INSFVNoSlipWallBC
    boundary = 'right'
    variable = v
    function = 0
  []
  [free-slip-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'right'
    variable = u
  []
  [free-slip-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'right'
    variable = v
  []
  [axis-u]
    type = INSFVSymmetryVelocityBC
    boundary = 'left'
    variable = u
    u = u
    v = v
    mu = ${mu}
    momentum_component = x
  []
  [axis-v]
    type = INSFVSymmetryVelocityBC
    boundary = 'left'
    variable = v
    u = u
    v = v
    mu = ${mu}
    momentum_component = y
  []
  [axis-p]
    type = INSFVSymmetryPressureBC
    boundary = 'left'
    variable = pressure
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'top'
    variable = pressure
    function = 0
  []
  [inlet_temp]
    type = FVNeumannBC
    boundary = 'bottom'
    variable = temperature
    value = 300
  []
[]
[Materials]
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'temperature'
    rho = ${rho}
  []
  [advected_material_property]
    type = ADGenericConstantMaterial
    prop_names = 'advected_rho cp'
    prop_values ='${rho} 1'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      200                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]
[Postprocessors]
  [inlet_mass_variable]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = advected_density
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_mass_constant]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = ${rho}
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_mass_matprop]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_mat_prop = 'advected_rho'
    fv = true
  []
  [mid1_mass]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_mass]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_mass]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_momentum_x]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid1_momentum_x]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_momentum_x]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_momentum_x]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    advected_variable = u
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_momentum_y]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid1_momentum_y]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_momentum_y]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_momentum_y]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    advected_variable = v
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [inlet_advected_energy]
    type = VolumetricFlowRate
    boundary = bottom
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
  []
  [mid1_advected_energy]
    type = InternalVolumetricFlowRate
    boundary = internal_bot
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [mid2_advected_energy]
    type = InternalVolumetricFlowRate
    boundary = internal_top
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
    advected_interp_method = ${advected_interp_method}
  []
  [outlet_advected_energy]
    type = VolumetricFlowRate
    boundary = top
    vel_x = u
    vel_y = v
    advected_mat_prop = 'rho_cp_temp'
    fv = true
  []
[]
[Outputs]
  exodus = false
  csv = true
  inactive = 'console_mass console_momentum_x console_momentum_y console_energy'
  [console_mass]
    type = Console
    start_step = 1
    show = 'inlet_mass_variable inlet_mass_constant inlet_mass_matprop mid1_mass mid2_mass outlet_mass'
  []
  [console_momentum_x]
    type = Console
    start_step = 1
    show = 'inlet_momentum_x mid1_momentum_x mid2_momentum_x outlet_momentum_x'
  []
  [console_momentum_y]
    type = Console
    start_step = 1
    show = 'inlet_momentum_y mid1_momentum_y mid2_momentum_y outlet_momentum_y'
  []
  [console_energy]
    type = Console
    start_step = 1
    show = 'inlet_advected_energy mid1_advected_energy mid2_advected_energy outlet_advected_energy'
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-ambient-convection.i)
mu=1
rho=1
k=1e-3
cp=1
advected_interp_method='average'
velocity_interp_method='rc'
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 5
    ymin = -1
    ymax = 1
    nx = 50
    ny = 16
  []
[]
[Problem]
  fv_bcs_integrity_check = true
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [temperature]
    type = INSFVEnergyVariable
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [energy_advection]
    type = INSFVEnergyAdvection
    variable = temperature
    vel = 'velocity'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [energy_diffusion]
    type = FVDiffusion
    coeff = ${k}
    variable = temperature
  []
  [ambient_convection]
    type = NSFVEnergyAmbientConvection
    variable = temperature
    T_ambient = 100
    alpha = 'alpha'
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 0
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = u
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = v
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0
  []
  [inlet_t]
    type = FVDirichletBC
    boundary = 'left'
    variable = temperature
    value = 1
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'cp alpha'
    prop_values = '${cp} 1'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    rho = ${rho}
    temperature = 'temperature'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      100                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
  csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/block_restriction/2d-rc.i)
mu=1.1
rho=1.1
advected_interp_method='average'
velocity_interp_method='rc'
restricted_blocks = '1'
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  parallel_type = 'replicated'
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1 1'
    dy = '1'
    ix = '7 7'
    iy = 10
    subdomain_id = '1 2'
  []
  [mid]
    type = SideSetsBetweenSubdomainsGenerator
    primary_block = 1
    paired_block = 2
    input = mesh
    new_boundary = 'middle'
  []
  [break_top]
    type = PatchSidesetGenerator
    boundary = 'top'
    n_patches = 2
    input = mid
  []
  [break_bottom]
    type = PatchSidesetGenerator
    boundary = 'bottom'
    n_patches = 2
    input = break_top
  []
[]
[Problem]
  kernel_coverage_check = false
  fv_bcs_integrity_check = true
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 1
    block = ${restricted_blocks}
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
    block = ${restricted_blocks}
  []
  [pressure]
    type = INSFVPressureVariable
    block = ${restricted_blocks}
  []
  [temperature]
    type = INSFVEnergyVariable
    block = ${restricted_blocks}
  []
  [scalar]
    type = INSFVScalarFieldVariable
    block = ${restricted_blocks}
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [energy_advection]
    type = INSFVEnergyAdvection
    variable = temperature
    vel = 'velocity'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [energy_diffusion]
    type = FVDiffusion
    coeff = 1.1
    variable = temperature
  []
  [energy_loss]
    type = FVBodyForce
    variable = temperature
    value = -0.1
  []
  [scalar_advection]
    type = INSFVScalarFieldAdvection
    variable = scalar
    vel = 'velocity'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [scalar_diffusion]
    type = FVDiffusion
    coeff = 1
    variable = scalar
  []
  [scalar_src]
    type = FVBodyForce
    variable = scalar
    value = 0.1
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 0
  []
  [top-wall-u]
    type = INSFVNoSlipWallBC
    boundary = 'top_0'
    variable = u
    function = 0
  []
  [top-wall-v]
    type = INSFVNoSlipWallBC
    boundary = 'top_0'
    variable = v
    function = 0
  []
  [bottom-wall-u]
    type = INSFVSymmetryVelocityBC
    boundary = 'bottom_0'
    variable = u
    mu = ${mu}
    u = u
    v = v
    momentum_component = 'x'
  []
  [bottom-wall-v]
    type = INSFVSymmetryVelocityBC
    boundary = 'bottom_0'
    variable = v
    mu = ${mu}
    u = u
    v = v
    momentum_component = 'y'
  []
  [bottom-wall-p]
    type = INSFVSymmetryPressureBC
    boundary = 'bottom_0'
    variable = pressure
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'middle'
    variable = pressure
    function = 0
  []
  [inlet_t]
    type = FVDirichletBC
    boundary = 'left'
    variable = temperature
    value = 1
  []
  [outlet_scalar]
    type = FVDirichletBC
    boundary = 'middle'
    variable = scalar
    value = 1
  []
[]
[Materials]
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'temperature'
    rho = ${rho}
    block = ${restricted_blocks}
  []
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'cp'
    prop_values = '2'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      100                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
  csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/boussinesq/boussinesq.i)
mu = 1
rho = 1
k = 1
cp = 1
alpha = 1
vel = 'velocity'
velocity_interp_method = 'rc'
advected_interp_method = 'upwind'
rayleigh=1e3
hot_temp=${rayleigh}
temp_ref=${fparse hot_temp / 2.}
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 1
    ymin = 0
    ymax = 1
    nx = 32
    ny = 32
  []
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
  []
  [v]
    type = INSFVVelocityVariable
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [T]
    type = INSFVEnergyVariable
    scaling = 1e-4
  []
  [lambda]
    family = SCALAR
    order = FIRST
  []
[]
[AuxVariables]
  [U]
    order = CONSTANT
    family = MONOMIAL
    fv = true
  []
  [vel_x]
    order = FIRST
    family = MONOMIAL
  []
  [vel_y]
    order = FIRST
    family = MONOMIAL
  []
  [viz_T]
    order = FIRST
    family = MONOMIAL
  []
[]
[AuxKernels]
  [mag]
    type = VectorMagnitudeAux
    variable = U
    x = u
    y = v
    execute_on = 'initial timestep_end'
  []
  [vel_x]
    type = ParsedAux
    variable = vel_x
    function = 'u'
    execute_on = 'initial timestep_end'
    args = 'u'
  []
  [vel_y]
    type = ParsedAux
    variable = vel_y
    function = 'v'
    execute_on = 'initial timestep_end'
    args = 'v'
  []
  [viz_T]
    type = ParsedAux
    variable = viz_T
    function = 'T'
    execute_on = 'initial timestep_end'
    args = 'T'
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    vel = ${vel}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    u = u
    v = v
    pressure = pressure
    mu = ${mu}
    rho = ${rho}
  []
  [mean_zero_pressure]
    type = FVScalarLagrangeMultiplier
    variable = pressure
    lambda = lambda
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [u_buoyancy]
    type = INSFVMomentumBoussinesq
    variable = u
    temperature = T
    gravity = '0 -1 0'
    rho = ${rho}
    ref_temperature = ${temp_ref}
    momentum_component = 'x'
  []
  [u_gravity]
    type = INSFVMomentumGravity
    variable = u
    gravity = '0 -1 0'
    rho = ${rho}
    momentum_component = 'x'
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [v_buoyancy]
    type = INSFVMomentumBoussinesq
    variable = v
    temperature = T
    gravity = '0 -1 0'
    rho = ${rho}
    ref_temperature = ${temp_ref}
    momentum_component = 'y'
  []
  [v_gravity]
    type = INSFVMomentumGravity
    variable = v
    gravity = '0 -1 0'
    rho = ${rho}
    momentum_component = 'y'
  []
  [temp_conduction]
    type = FVDiffusion
    coeff = 'k'
    variable = T
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    variable = T
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
[]
[FVBCs]
  [top_x]
    type = INSFVNoSlipWallBC
    variable = u
    boundary = 'top'
    function = 'lid_function'
  []
  [no_slip_x]
    type = INSFVNoSlipWallBC
    variable = u
    boundary = 'left right bottom'
    function = 0
  []
  [no_slip_y]
    type = INSFVNoSlipWallBC
    variable = v
    boundary = 'left right top bottom'
    function = 0
  []
  [T_hot]
    type = FVDirichletBC
    variable = T
    boundary = left
    value = ${hot_temp}
  []
  [T_cold]
    type = FVDirichletBC
    variable = T
    boundary = right
    value = 0
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'k cp alpha'
    prop_values = '${k} ${cp} ${alpha}'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'T'
    rho = ${rho}
  []
[]
[Functions]
  [lid_function]
    type = ParsedFunction
    value = '4*x*(1-x)'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      300                lu           NONZERO'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/lid-driven-with-energy.i)
mu = 1
rho = 1
k = .01
cp = 1
vel = 'velocity'
velocity_interp_method = 'rc'
advected_interp_method = 'average'
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 1
    ymin = 0
    ymax = 1
    nx = 32
    ny = 32
  []
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
  []
  [v]
    type = INSFVVelocityVariable
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [T]
    type = INSFVEnergyVariable
  []
  [lambda]
    family = SCALAR
    order = FIRST
  []
[]
[AuxVariables]
  [U]
    order = CONSTANT
    family = MONOMIAL
    fv = true
  []
[]
[AuxKernels]
  [mag]
    type = VectorMagnitudeAux
    variable = U
    x = u
    y = v
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    vel = ${vel}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    u = u
    v = v
    pressure = pressure
    mu = ${mu}
    rho = ${rho}
  []
  [mean_zero_pressure]
    type = FVScalarLagrangeMultiplier
    variable = pressure
    lambda = lambda
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [temp-condution]
    type = FVDiffusion
    coeff = 'k'
    variable = T
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    variable = T
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
[]
[FVBCs]
  [top_x]
    type = INSFVNoSlipWallBC
    variable = u
    boundary = 'top'
    function = 'lid_function'
  []
  [no_slip_x]
    type = INSFVNoSlipWallBC
    variable = u
    boundary = 'left right bottom'
    function = 0
  []
  [no_slip_y]
    type = INSFVNoSlipWallBC
    variable = v
    boundary = 'left right top bottom'
    function = 0
  []
  [T_hot]
    type = FVDirichletBC
    variable = T
    boundary = 'bottom'
    value = 1
  []
  [T_cold]
    type = FVDirichletBC
    variable = T
    boundary = 'top'
    value = 0
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'k cp'
    prop_values = '${k} ${cp}'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'T'
    rho = ${rho}
  []
[]
[Functions]
  [lid_function]
    type = ParsedFunction
    value = '4*x*(1-x)'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      300                lu           NONZERO'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/transient-lid-driven-with-energy.i)
mu = 1
rho = 1
k = .01
cp = 1
vel = 'velocity'
velocity_interp_method = 'rc'
advected_interp_method = 'average'
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 1
    ymin = 0
    ymax = 1
    nx = 32
    ny = 32
  []
  [pin]
    type = ExtraNodesetGenerator
    input = gen
    new_boundary = 'pin'
    nodes = '0'
  []
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
  []
  [v]
    type = INSFVVelocityVariable
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [T]
    type = INSFVEnergyVariable
  []
  [lambda]
    family = SCALAR
    order = FIRST
  []
[]
[ICs]
  [T]
    type = ConstantIC
    variable = T
    value = 1
  []
[]
[AuxVariables]
  [U]
    order = CONSTANT
    family = MONOMIAL
    fv = true
  []
[]
[AuxKernels]
  [mag]
    type = VectorMagnitudeAux
    variable = U
    x = u
    y = v
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    vel = ${vel}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    u = u
    v = v
    pressure = pressure
    mu = ${mu}
    rho = ${rho}
  []
  [mean_zero_pressure]
    type = FVScalarLagrangeMultiplier
    variable = pressure
    lambda = lambda
  []
  [u_time]
    type = INSFVMomentumTimeDerivative
    variable = 'u'
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [v_time]
    type = INSFVMomentumTimeDerivative
    variable = v
    rho = ${rho}
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [temp_time]
    type = INSFVEnergyTimeDerivative
    variable = T
    rho = ${rho}
  []
  [temp_conduction]
    type = FVDiffusion
    coeff = 'k'
    variable = T
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    variable = T
    vel = ${vel}
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
[]
[FVBCs]
  [top_x]
    type = INSFVNoSlipWallBC
    variable = u
    boundary = 'top'
    function = 'lid_function'
  []
  [no_slip_x]
    type = INSFVNoSlipWallBC
    variable = u
    boundary = 'left right bottom'
    function = 0
  []
  [no_slip_y]
    type = INSFVNoSlipWallBC
    variable = v
    boundary = 'left right top bottom'
    function = 0
  []
  [T_hot]
    type = FVDirichletBC
    variable = T
    boundary = 'bottom'
    value = 1
  []
  [T_cold]
    type = FVDirichletBC
    variable = T
    boundary = 'top'
    value = 0
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'k cp'
    prop_values = '${k} ${cp}'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'T'
    rho = ${rho}
  []
[]
[Functions]
  [lid_function]
    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_type -sub_pc_factor_shift_type -ksp_gmres_restart'
  petsc_options_value = 'asm      lu           NONZERO                   200'
  line_search = 'none'
  nl_rel_tol = 1e-12
  nl_max_its = 6
  l_max_its = 200
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/2d-average-with-temp.i)
mu=1.1
rho=1.1
k=1.1
cp=1.1
advected_interp_method='average'
velocity_interp_method='average'
velocity='velocity'
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 2
    ymin = -1
    ymax = 1
    nx = 2
    ny = 2
  []
[]
[Problem]
  fv_bcs_integrity_check = true
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 1
    two_term_boundary_expansion = false
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
    two_term_boundary_expansion = false
  []
  [pressure]
    type = INSFVPressureVariable
    two_term_boundary_expansion = false
  []
  [temperature]
    type = INSFVEnergyVariable
    two_term_boundary_expansion = false
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = ${velocity}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [mass_forcing]
    type = FVBodyForce
    variable = pressure
    function = forcing_p
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = ${velocity}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [u_forcing]
    type = FVBodyForce
    variable = u
    function = forcing_u
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = ${velocity}
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [v_forcing]
    type = FVBodyForce
    variable = v
    function = forcing_v
  []
  [temp_conduction]
    type = FVDiffusion
    coeff = 'k'
    variable = temperature
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    vel = ${velocity}
    variable = temperature
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [temp_forcing]
    type = FVBodyForce
    variable = temperature
    function = forcing_t
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = 'exact_u'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 'exact_v'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = u
    function = 'exact_u'
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = v
    function = 'exact_v'
  []
  [inlet-and-walls-t]
    type = FVFunctionDirichletBC
    boundary = 'left top bottom'
    variable = temperature
    function = 'exact_t'
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 'exact_p'
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'k cp'
    prop_values = '${k} ${cp}'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    temperature = 'temperature'
    rho = ${rho}
  []
[]
[Functions]
[exact_u]
  type = ParsedFunction
  value = 'sin((1/2)*y*pi)*cos((1/2)*x*pi)'
[]
[exact_rhou]
  type = ParsedFunction
  value = 'rho*sin((1/2)*y*pi)*cos((1/2)*x*pi)'
  vars = 'rho'
  vals = '${rho}'
[]
[forcing_u]
  type = ParsedFunction
  value = '(1/2)*pi^2*mu*sin((1/2)*y*pi)*cos((1/2)*x*pi) - 1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi)^2*cos((1/2)*x*pi) + (1/2)*pi*rho*sin((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi)^2 - pi*rho*sin((1/2)*x*pi)*sin((1/2)*y*pi)^2*cos((1/2)*x*pi) - 1/4*pi*sin((1/4)*x*pi)*sin((3/2)*y*pi)'
  vars = 'mu rho'
  vals = '${mu} ${rho}'
[]
[exact_v]
  type = ParsedFunction
  value = 'sin((1/4)*x*pi)*cos((1/2)*y*pi)'
[]
[exact_rhov]
  type = ParsedFunction
  value = 'rho*sin((1/4)*x*pi)*cos((1/2)*y*pi)'
  vars = 'rho'
  vals = '${rho}'
[]
[forcing_v]
  type = ParsedFunction
  value = '(5/16)*pi^2*mu*sin((1/4)*x*pi)*cos((1/2)*y*pi) - pi*rho*sin((1/4)*x*pi)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*y*pi) + (1/4)*pi*rho*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + (3/2)*pi*cos((1/4)*x*pi)*cos((3/2)*y*pi)'
  vars = 'mu rho'
  vals = '${mu} ${rho}'
[]
[exact_p]
  type = ParsedFunction
  value = 'sin((3/2)*y*pi)*cos((1/4)*x*pi)'
[]
[forcing_p]
  type = ParsedFunction
  value = '-1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi) - 1/2*pi*rho*sin((1/2)*x*pi)*sin((1/2)*y*pi)'
  vars = 'rho'
  vals = '${rho}'
[]
[exact_t]
  type = ParsedFunction
  value = 'sin((1/4)*x*pi)*cos((1/2)*y*pi)'
[]
[forcing_t]
  type = ParsedFunction
  value = '-pi*cp*rho*sin((1/4)*x*pi)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*cp*rho*sin((1/4)*x*pi)*sin((1/2)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*y*pi) + (1/4)*pi*cp*rho*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + (5/16)*pi^2*k*sin((1/4)*x*pi)*cos((1/2)*y*pi)'
  vars = 'k rho cp'
  vals = '${k} ${rho} ${cp}'
[]
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      100                lu           NONZERO'
  line_search = 'none'
[]
[Outputs]
  exodus = true
  csv = true
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    outputs = 'console csv'
    execute_on = 'timestep_end'
  []
  [./L2u]
    type = ElementL2Error
    variable = u
    function = exact_u
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./L2v]
    type = ElementL2Error
    variable = v
    function = exact_v
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./L2p]
    variable = pressure
    function = exact_p
    type = ElementL2Error
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
  [./L2t]
    variable = temperature
    function = exact_t
    type = ElementL2Error
    outputs = 'console csv'
    execute_on = 'timestep_end'
  [../]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-scalar-transport.i)
mu=1
rho=1
k=1e-3
diff=1e-3
cp=1
advected_interp_method='average'
velocity_interp_method='rc'
[GlobalParams]
  two_term_boundary_expansion = true
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 10
    ymin = -1
    ymax = 1
    nx = 100
    ny = 20
  []
[]
[Problem]
  fv_bcs_integrity_check = true
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [temperature]
    type = INSFVEnergyVariable
  []
  [scalar]
    type = INSFVScalarFieldVariable
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    vel = 'velocity'
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_quantity = 'rhou'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [u_viscosity]
    type = FVDiffusion
    variable = u
    coeff = ${mu}
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    p = pressure
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_quantity = 'rhov'
    vel = 'velocity'
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [v_viscosity]
    type = FVDiffusion
    variable = v
    coeff = ${mu}
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    p = pressure
  []
  [energy_advection]
    type = INSFVEnergyAdvection
    variable = temperature
    vel = 'velocity'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [energy_diffusion]
    type = FVDiffusion
    coeff = ${k}
    variable = temperature
  []
  [scalar_advection]
    type = INSFVScalarFieldAdvection
    variable = scalar
    vel = 'velocity'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
    pressure = pressure
    u = u
    v = v
    mu = ${mu}
    rho = ${rho}
  []
  [scalar_diffusion]
    type = FVDiffusion
    coeff = ${diff}
    variable = scalar
  []
  [scalar_src]
    type = FVBodyForce
    variable = scalar
    value = 0.1
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 0
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = u
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = v
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0
  []
  [inlet_t]
    type = FVDirichletBC
    boundary = 'left'
    variable = temperature
    value = 1
  []
  [inlet_scalar]
    type = FVDirichletBC
    boundary = 'left'
    variable = scalar
    value = 1
  []
[]
[Materials]
  [const]
    type = ADGenericConstantMaterial
    prop_names = 'cp'
    prop_values = '${cp}'
  []
  [ins_fv]
    type = INSFVMaterial
    u = 'u'
    v = 'v'
    pressure = 'pressure'
    rho = ${rho}
    temperature = 'temperature'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      100                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
  csv = true
[]