- variableThe name of the variable whose linear system this object contributes to
C++ Type:LinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable whose linear system this object contributes to
 
LinearFVSource
Description
This kernel contributes to the right hand side of a system which is solved for a linear finite volume variable MooseVariableLinearFVReal. The contribution for each cell is the numerical integral of the source density function () in the following form:
where and denote the source density at the cell centroid and the cell volume, respectively. This integral is added to the corresponding entry of the right hand side of the linear system. The source density parameter ("source_density") accepts anything that supports functor-based evaluations. For more information on functors in MOOSE, see Functor system.
Example Syntax
The case below demonstrates the use of LinearFVSource where the force term is supplied based upon a function form:
[LinearFVKernels<<<{"href": "../../syntax/LinearFVKernels/index.html"}>>>]
  [reaction]
    type = LinearFVReaction<<<{"description": "Represents the matrix and right hand side contributions of a reaction term ($c u$) in a partial differential equation.", "href": "LinearFVReaction.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = u
    coeff<<<{"description": "The reaction coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = coeff_func
  []
  [source]
    type = LinearFVSource<<<{"description": "Represents the matrix and right hand side contributions of a solution-independent source term in a partial differential equation.", "href": "LinearFVSource.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = u
    source_density<<<{"description": "The source density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = source_func
  []
[](test/tests/linearfvkernels/reaction/reaction-1d.i)Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
 - matrix_onlyFalseWhether this object is only doing assembly to matrices (no vectors)
Default:False
C++ Type:bool
Controllable:No
Description:Whether this object is only doing assembly to matrices (no vectors)
 - scaling_factor1Coefficient to multiply the body force term with. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:1
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:Coefficient to multiply the body force term with. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
 - source_density1The source density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:1
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The source density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
 
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
 - extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
 - extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
 - matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
 - vector_tagsrhsThe tag for the vectors this Kernel should fill
Default:rhs
C++ Type:MultiMooseEnum
Options:rhs, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
 
Contribution To Tagged Field Data Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
 - enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
 - implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
 - seed0The seed for the master random number generator
Default: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
- ghost_layers1The number of layers of elements to ghost.
Default:1
C++ Type:unsigned short
Controllable:No
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
Controllable:No
Description:Whether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
 
Parallel Ghosting Parameters
Input Files
- (test/tests/linearfvkernels/advection/advection-1d.i)
 - (test/tests/linearfvkernels/diffusion/diffusion-2d.i)
 - (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/segregated/channel-drift-flux.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/cht/bulk_heat_transfer/flow-around-square-linear.i)
 - (modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h.i)
 - (test/tests/variables/linearfv/diffusion-1d-pp.i)
 - (test/tests/linearfvbcs/robin/diffusion-1d-robin.i)
 - (test/tests/linearfvkernels/diffusion/diffusion-2d-rz.i)
 - (test/tests/multiapps/linearfv_nonlinearfv/linearfv.i)
 - (test/tests/variables/linearfv/diffusion-1d-aux.i)
 - (test/tests/linearfvkernels/diffusion/diffusion-2d_neumann.i)
 - (test/tests/linearfvkernels/anisotropic-diffusion/anisotropic-diffusion-2d.i)
 - (test/tests/linearfvbcs/robin/advection-1d-robin.i)
 - (test/tests/linearfvkernels/block-restriction/block-restricted-adr.i)
 - (test/tests/linearfvkernels/advection/advection-2d.i)
 - (test/tests/multisystem/restore_multiapp/sub.i)
 - (test/tests/time_steppers/iteration_adaptive/adapt_linear_systems.i)
 - (test/tests/linearfvkernels/block-restriction/block-restricted-diffusion.i)
 - (test/tests/tag/tag-linearfv.i)
 - (test/tests/time_integrators/implicit-euler/ie-linearfv.i)
 - (test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-1d.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/mms/linear-segregated/2d-vortex/2d-vortex.i)
 - (test/tests/outputs/debug/show_execution_linear_fv_elemental.i)
 - (test/tests/linearfvkernels/diffusion/diffusion-1d_neumann.i)
 - (test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-2d.i)
 - (test/tests/linearfvkernels/diffusion/diffusion-1d.i)
 - (modules/navier_stokes/test/tests/finite_volume/ins/mms/linear-segregated/2d-vortex/spacedependent_mu/sl.i)
 - (test/tests/linearfvbcs/robin/advection-2d-robin.i)
 - (test/tests/multisystem/picard/linearfv_nonlinearfv/same_input.i)
 - (test/tests/linearfvbcs/robin/diffusion-2d-robin.i)
 - (test/tests/multisystem/restore_multiapp/parent.i)
 - (test/tests/linearfvkernels/reaction/reaction-1d.i)
 - (test/tests/linearfvkernels/advection/advection-2d-rz.i)
 - (test/tests/outputs/debug/show_execution_linear_fv_flux.i)
 - (test/tests/linearfvkernels/block-restriction/block-restricted-diffusion-react.i)
 - (modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h_fp.i)
 - (test/tests/transfers/multiapp_copy_transfer/linear_sys_to_aux/linear_sub.i)
 
source_density
Default:1
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The source density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
(test/tests/linearfvkernels/reaction/reaction-1d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 10
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '1+sin(x)'
  []
  [source_func]
    type = ParsedFunction
    expression = '(1+sin(x))*(1+cos(x))'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1+cos(x)'
  []
[]
[Postprocessors]
  [l2error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = TIMESTEP_END
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
[]
[Outputs]
  [exodus]
    type = Exodus
    execute_on = TIMESTEP_END
  []
[]
(test/tests/linearfvkernels/advection/advection-1d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 2
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "0.5 0 0"
    advected_interp_method = upwind
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [inflow]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left"
    functor = analytic_solution
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = false
  []
[]
[Functions]
  [source_func]
    type = ParsedFunction
    expression = '0.5*x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '0.5+0.5*x*x'
  []
[]
[Postprocessors]
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvkernels/diffusion/diffusion-2d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 1
    ymax = 0.5
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
    use_nonorthogonal_correction = false
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right top bottom"
    functor = analytic_solution
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '1+0.5*x*y'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*(1.5-y*y)+2*x*y*(1.5-y*y)+2*(1.5-x*x)+2*x*y*(1.5-x*x)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '(1.5-x*x)*(1.5-y*y)'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/segregated/channel-drift-flux.i)
mu = 1.0
rho = 10.0
mu_d = 0.1
rho_d = 1.0
l = 2
# 'average' leads to slight oscillations, upwind may be preferred
# This method is selected for consistency with the original nonlinear input
advected_interp_method = 'average'
# TODO remove need for those
cp = 1
k = 1
cp_d = 1
k_d = 1
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = '${fparse l * 5}'
    ymin = '${fparse -l / 2}'
    ymax = '${fparse l / 2}'
    nx = 10
    ny = 4
  []
  uniform_refine = 0
[]
[Problem]
  linear_sys_names = 'u_system v_system pressure_system phi_system'
  previous_nl_solution_required = true
[]
[Variables]
  [vel_x]
    type = MooseLinearVariableFVReal
    solver_sys = u_system
    initial_condition = 1
  []
  [vel_y]
    type = MooseLinearVariableFVReal
    solver_sys = v_system
  []
  [pressure]
    type = MooseLinearVariableFVReal
    solver_sys = pressure_system
  []
  [phase_2]
    type = MooseLinearVariableFVReal
    solver_sys = phi_system
  []
[]
[LinearFVKernels]
  [flow_p_diffusion]
    type = LinearFVAnisotropicDiffusion
    diffusion_tensor = Ainv
    use_nonorthogonal_correction = false
    variable = pressure
  []
  [flow_HbyA_divergence]
    type = LinearFVDivergence
    face_flux = HbyA
    force_boundary_execution = true
    variable = pressure
  []
  [flow_ins_momentum_flux_x]
    type = LinearWCNSFVMomentumFlux
    advected_interp_method = ${advected_interp_method}
    momentum_component = x
    mu = mu_mixture
    rhie_chow_user_object = ins_rhie_chow_interpolator
    u = vel_x
    use_deviatoric_terms = false
    use_nonorthogonal_correction = false
    v = vel_y
    variable = vel_x
  []
  [flow_ins_momentum_flux_y]
    type = LinearWCNSFVMomentumFlux
    advected_interp_method = ${advected_interp_method}
    momentum_component = y
    mu = mu_mixture
    rhie_chow_user_object = ins_rhie_chow_interpolator
    u = vel_x
    use_deviatoric_terms = false
    use_nonorthogonal_correction = false
    v = vel_y
    variable = vel_y
  []
  [mixture_drift_flux_x]
    type = LinearWCNSFV2PMomentumDriftFlux
    density_interp_method = average
    fraction_dispersed = phase_2
    momentum_component = x
    rhie_chow_user_object = ins_rhie_chow_interpolator
    rho_d = ${rho_d}
    u_slip = vel_slip_x
    v_slip = vel_slip_y
    variable = vel_x
  []
  [mixture_drift_flux_y]
    type = LinearWCNSFV2PMomentumDriftFlux
    density_interp_method = average
    fraction_dispersed = phase_2
    momentum_component = y
    rhie_chow_user_object = ins_rhie_chow_interpolator
    rho_d = ${rho_d}
    u_slip = vel_slip_x
    v_slip = vel_slip_y
    variable = vel_y
  []
  [flow_ins_momentum_pressure_x]
    type = LinearFVMomentumPressure
    momentum_component = x
    pressure = pressure
    variable = vel_x
  []
  [flow_ins_momentum_pressure_y]
    type = LinearFVMomentumPressure
    momentum_component = y
    pressure = pressure
    variable = vel_y
  []
  [flow_momentum_friction_0_x]
    type = LinearFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    momentum_component = x
    mu = mu_mixture
    variable = vel_x
  []
  [flow_momentum_friction_0_y]
    type = LinearFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    momentum_component = y
    mu = mu_mixture
    variable = vel_y
  []
  # Mixture phase equation
  [mixture_ins_phase_2_advection]
    type = LinearFVScalarAdvection
    advected_interp_method = upwind
    rhie_chow_user_object = ins_rhie_chow_interpolator
    u_slip = vel_slip_x
    v_slip = vel_slip_y
    variable = phase_2
  []
  [mixture_phase_interface_reaction]
    type = LinearFVReaction
    coeff = 0.1
    variable = phase_2
  []
  [mixture_phase_interface_source]
    type = LinearFVSource
    scaling_factor = 0.1
    source_density = phase_1
    variable = phase_2
  []
[]
[LinearFVBCs]
  [vel_x_left]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = left
    functor = 1
    variable = vel_x
  []
  [vel_y_left]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = left
    functor = 0
    variable = vel_y
  []
  [pressure_extrapolation_inlet_left]
    type = LinearFVExtrapolatedPressureBC
    boundary = left
    use_two_term_expansion = true
    variable = pressure
  []
  [vel_x_right]
    type = LinearFVAdvectionDiffusionOutflowBC
    boundary = right
    use_two_term_expansion = true
    variable = vel_x
  []
  [vel_y_right]
    type = LinearFVAdvectionDiffusionOutflowBC
    boundary = right
    use_two_term_expansion = true
    variable = vel_y
  []
  [pressure_right]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = right
    functor = 0
    variable = pressure
  []
  [vel_x_bottom]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = bottom
    functor = 0
    variable = vel_x
  []
  [vel_y_bottom]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = bottom
    functor = 0
    variable = vel_y
  []
  [vel_x_top]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = top
    functor = 0
    variable = vel_x
  []
  [vel_y_top]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = top
    functor = 0
    variable = vel_y
  []
  [pressure_extrapolation_top_bottom]
    type = LinearFVExtrapolatedPressureBC
    boundary = 'top bottom'
    use_two_term_expansion = true
    variable = pressure
  []
  [phase_2_left]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = left
    functor = 0.1
    variable = phase_2
  []
  [phase_2_right]
    type = LinearFVAdvectionDiffusionOutflowBC
    boundary = right
    use_two_term_expansion = true
    variable = phase_2
  []
[]
[FunctorMaterials]
  [flow_ins_speed_material]
    type = ADVectorMagnitudeFunctorMaterial
    execute_on = ALWAYS
    outputs = none
    vector_magnitude_name = speed
    x_functor = vel_x
    y_functor = vel_y
  []
  [mixture_phase_1_fraction]
    type = ParsedFunctorMaterial
    execute_on = ALWAYS
    expression = '1 - phase_2'
    functor_names = phase_2
    output_properties = phase_1
    outputs = all
    property_name = phase_1
  []
  [mixture_mixture_material]
    type = WCNSLinearFVMixtureFunctorMaterial
    execute_on = ALWAYS
    limit_phase_fraction = true
    outputs = all
    phase_1_fraction = phase_2
    phase_1_names = '${rho_d} ${mu_d} ${cp_d} ${k_d}'
    phase_2_names = '${rho}   ${mu}   ${cp}   ${k}'
    prop_names = 'rho_mixture mu_mixture cp_mixture k_mixture'
  []
  [mixture_slip_x]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    execute_on = ALWAYS
    gravity = '0 0 0'
    linear_coef_name = Darcy_coefficient
    momentum_component = x
    mu = mu_mixture
    outputs = all
    particle_diameter = 0.01
    rho = ${rho}
    rho_d = ${rho_d}
    slip_velocity_name = vel_slip_x
    u = vel_x
    v = vel_y
  []
  [mixture_slip_y]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    execute_on = ALWAYS
    gravity = '0 0 0'
    linear_coef_name = Darcy_coefficient
    momentum_component = y
    mu = mu_mixture
    outputs = all
    particle_diameter = 0.01
    rho = ${rho}
    rho_d = ${rho_d}
    slip_velocity_name = vel_slip_y
    u = vel_x
    v = vel_y
  []
  [mixture_dispersed_drag]
    type = NSFVDispersePhaseDragFunctorMaterial
    drag_coef_name = Darcy_coefficient
    execute_on = ALWAYS
    mu = mu_mixture
    outputs = all
    particle_diameter = 0.01
    rho = rho_mixture
    u = vel_x
    v = vel_y
  []
[]
[UserObjects]
  [ins_rhie_chow_interpolator]
    type = RhieChowMassFlux
    p_diffusion_kernel = flow_p_diffusion
    pressure = pressure
    rho = rho_mixture
    u = vel_x
    v = vel_y
  []
[]
[Executioner]
  type = SIMPLE
  rhie_chow_user_object = 'ins_rhie_chow_interpolator'
  # Systems
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  active_scalar_systems = 'phi_system'
  momentum_equation_relaxation = 0.8
  active_scalar_equation_relaxation = '0.7'
  pressure_variable_relaxation = 0.3
  # We need to converge the problem to show conservation
  num_iterations = 200
  pressure_absolute_tolerance = 1e-10
  momentum_absolute_tolerance = 1e-10
  active_scalar_absolute_tolerance = '1e-10'
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  active_scalar_petsc_options_iname = '-pc_type -pc_hypre_type'
  active_scalar_petsc_options_value = 'hypre boomeramg'
  momentum_l_abs_tol = 1e-13
  pressure_l_abs_tol = 1e-13
  active_scalar_l_abs_tol = 1e-13
  momentum_l_tol = 0
  pressure_l_tol = 0
  active_scalar_l_tol = 0
  # print_fields = true
  continue_on_max_its = true
[]
[Outputs]
  csv = true
[]
[Postprocessors]
  [Re]
    type = ParsedPostprocessor
    expression = '10.0 * 2 * 1'
  []
  [average_phase2]
    type = ElementAverageValue
    variable = phase_2
  []
  [dp]
    type = PressureDrop
    boundary = 'left right'
    downstream_boundary = right
    pressure = pressure
    upstream_boundary = left
  []
  [max_phase2]
    type = ElementExtremeValue
    variable = phase_2
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/bulk_heat_transfer/flow-around-square-linear.i)
mu = 0.01
rho = 1.1
k = 0.0005
cp = 10
k_s = 3.0
h_conv = 5
power_density = 10000
advected_interp_method = 'upwind'
[Mesh]
  [generated_mesh]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 10
    ny = 10
    xmin = 0
    ymin = 0
    ymax = 0.1
    xmax = 0.1
  []
  [subdomain1]
    type = SubdomainBoundingBoxGenerator
    input = generated_mesh
    block_name = subdomain1
    bottom_left = '0.04 0.04 0'
    block_id = 1
    top_right = '0.06 0.06 0'
  []
  [interface]
    type = SideSetsBetweenSubdomainsGenerator
    input = subdomain1
    primary_block = 0
    paired_block = 1
    new_boundary = interface
  []
[]
[Problem]
  linear_sys_names = 'u_system v_system pressure_system energy_system solid_energy_system'
  previous_nl_solution_required = true
[]
[UserObjects]
  [rc]
    type = RhieChowMassFlux
    u = vel_x
    v = vel_y
    pressure = pressure
    rho = ${rho}
    p_diffusion_kernel = p_diffusion
    block = 0
  []
[]
[Variables]
  [vel_x]
    type = MooseLinearVariableFVReal
    initial_condition = 0.1
    solver_sys = u_system
    block = 0
  []
  [vel_y]
    type = MooseLinearVariableFVReal
    solver_sys = v_system
    initial_condition = 0.01
    block = 0
  []
  [pressure]
    type = MooseLinearVariableFVReal
    solver_sys = pressure_system
    initial_condition = 0.2
    block = 0
  []
  [T_fluid]
    type = MooseLinearVariableFVReal
    solver_sys = energy_system
    initial_condition = 300
    block = 0
  []
  [T_solid]
    type = MooseLinearVariableFVReal
    solver_sys = solid_energy_system
    initial_condition = 500
    block = 1
  []
[]
[LinearFVKernels]
  [u_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    mu = ${mu}
    u = vel_x
    v = vel_y
    momentum_component = 'x'
    rhie_chow_user_object = 'rc'
    use_nonorthogonal_correction = true
  []
  [v_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    mu = ${mu}
    u = vel_x
    v = vel_y
    momentum_component = 'y'
    rhie_chow_user_object = 'rc'
    use_nonorthogonal_correction = true
  []
  [u_pressure]
    type = LinearFVMomentumPressure
    variable = vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [v_pressure]
    type = LinearFVMomentumPressure
    variable = vel_y
    pressure = pressure
    momentum_component = 'y'
  []
  [p_diffusion]
    type = LinearFVAnisotropicDiffusion
    variable = pressure
    diffusion_tensor = Ainv
    use_nonorthogonal_correction = true
  []
  [HbyA_divergence]
    type = LinearFVDivergence
    variable = pressure
    face_flux = HbyA
    force_boundary_execution = true
  []
  [h_advection]
    type = LinearFVEnergyAdvection
    variable = T_fluid
    advected_quantity = temperature
    cp = ${cp}
    advected_interp_method = ${advected_interp_method}
    rhie_chow_user_object = 'rc'
  []
  [conduction]
    type = LinearFVDiffusion
    variable = T_fluid
    diffusion_coeff = ${k}
    use_nonorthogonal_correction = true
  []
  [solid-conduction]
    type = LinearFVDiffusion
    variable = T_solid
    diffusion_coeff = ${k_s}
    use_nonorthogonal_correction = true
  []
  [solid-source]
    type = LinearFVSource
    variable = T_solid
    source_density = ${power_density}
  []
[]
[LinearFVBCs]
  [inlet-u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left'
    variable = vel_x
    functor = '0.1'
  []
  [inlet-v]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left'
    variable = vel_y
    functor = '0.0'
  []
  [walls-u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'top bottom interface'
    variable = vel_x
    functor = 0.0
  []
  [walls-v]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'top bottom interface'
    variable = vel_y
    functor = 0.0
  []
  [outlet_p]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'right'
    variable = pressure
    functor = 1.4
  []
  [outlet_u]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = vel_x
    use_two_term_expansion = false
    boundary = right
  []
  [outlet_v]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = vel_y
    use_two_term_expansion = false
    boundary = right
  []
  [inlet_T]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = T_fluid
    functor = 300
    boundary = 'left'
  []
  [walls_T]
    type = LinearFVAdvectionDiffusionFunctorNeumannBC
    variable = T_fluid
    functor = 0.0
    boundary = 'top bottom'
  []
  [outlet_T]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = T_fluid
    use_two_term_expansion = false
    boundary = right
  []
  [fluid_solid]
    type = LinearFVConvectiveHeatTransferBC
    variable = T_fluid
    T_solid = T_solid
    T_fluid = T_fluid
    boundary = interface
    h = ${h_conv}
  []
  [solid_fluid]
    type = LinearFVConvectiveHeatTransferBC
    variable = T_solid
    T_solid = T_solid
    T_fluid = T_fluid
    boundary = interface
    h = ${h_conv}
  []
[]
[FunctorMaterials]
  [rhocpT]
    property_name = 'rhocpT'
    type = ParsedFunctorMaterial
    functor_names = 'T_fluid'
    expression = '${rho}*${cp}*T_fluid'
  []
[]
[Postprocessors]
  [h_in]
    type = VolumetricFlowRate
    boundary = left
    vel_x = vel_x
    vel_y = vel_y
    rhie_chow_user_object = rc
    advected_quantity = 'rhocpT'
    subtract_mesh_velocity = false
  []
  [h_out]
    type = VolumetricFlowRate
    boundary = right
    vel_x = vel_x
    vel_y = vel_y
    rhie_chow_user_object = rc
    advected_quantity = 'rhocpT'
    advected_interp_method = upwind
    subtract_mesh_velocity = false
  []
  [power]
    type = ElementIntegralFunctorPostprocessor
    functor = ${power_density}
    block = 1
  []
[]
[Executioner]
  type = SIMPLE
  momentum_l_abs_tol = 1e-13
  pressure_l_abs_tol = 1e-13
  energy_l_abs_tol = 1e-13
  solid_energy_l_abs_tol = 1e-13
  momentum_l_tol = 0
  pressure_l_tol = 0
  energy_l_tol = 0
  solid_energy_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  energy_system = 'energy_system'
  solid_energy_system = 'solid_energy_system'
  momentum_equation_relaxation = 0.8
  energy_equation_relaxation = 1.0
  pressure_variable_relaxation = 0.3
  num_iterations = 1000
  pressure_absolute_tolerance = 1e-10
  momentum_absolute_tolerance = 1e-10
  energy_absolute_tolerance = 1e-10
  solid_energy_absolute_tolerance = 1e-10
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  energy_petsc_options_iname = '-pc_type -pc_hypre_type'
  energy_petsc_options_value = 'hypre boomeramg'
  solid_energy_petsc_options_iname = '-pc_type -pc_hypre_type'
  solid_energy_petsc_options_value = 'hypre boomeramg'
  print_fields = false
  continue_on_max_its = true
[]
[Outputs]
  exodus = true
  execute_on = timestep_end
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h.i)
L = 30
nx = 600
bulk_u = 0.01
q_source = 50000.
A_cp = 976.78
B_cp = 1.0634
T_in = 860.
p_ref = 101325.0
rho = 2000.
advected_interp_method = 'upwind'
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    xmin = 0
    xmax = ${L}
    nx = ${nx}
  []
  allow_renumbering = false
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
  advected_interp_method = ${advected_interp_method}
  u = vel_x
[]
[Problem]
  linear_sys_names = 'u_system pressure_system energy_system'
  previous_nl_solution_required = true
[]
[UserObjects]
  [rc]
    type = RhieChowMassFlux
    u = vel_x
    pressure = pressure
    rho = 'rho'
    p_diffusion_kernel = p_diffusion
  []
[]
[Variables]
  [vel_x]
    type = MooseLinearVariableFVReal
    solver_sys = u_system
    initial_condition = ${bulk_u}
  []
  [pressure]
    type = MooseLinearVariableFVReal
    solver_sys = pressure_system
    initial_condition = ${p_ref}
  []
  [h]
    type = MooseLinearVariableFVReal
    solver_sys = energy_system
    initial_condition = ${fparse 860.*1900.}
  []
[]
[AuxVariables]
  [rho_var]
    type = MooseLinearVariableFVReal
  []
  [cp_var]
    type = MooseLinearVariableFVReal
  []
  [mu_var]
    type = MooseLinearVariableFVReal
  []
  [k_var]
    type = MooseLinearVariableFVReal
  []
  [T]
    type = MooseLinearVariableFVReal
    initial_condition = 860.
  []
  [h_aux]
    type = MooseLinearVariableFVReal
  []
[]
[LinearFVKernels]
  [u_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_x
    mu = 'mu'
    momentum_component = 'x'
    use_nonorthogonal_correction = false
  []
  [u_pressure]
    type = LinearFVMomentumPressure
    variable = vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [p_diffusion]
    type = LinearFVAnisotropicDiffusion
    variable = pressure
    diffusion_tensor = Ainv
    use_nonorthogonal_correction = false
  []
  [HbyA_divergence]
    type = LinearFVDivergence
    variable = pressure
    face_flux = HbyA
    force_boundary_execution = true
  []
  [temp_advection]
    type = LinearFVEnergyAdvection
    variable = h
  []
  [source]
    type = LinearFVSource
    variable = h
    source_density = source_func
  []
[]
[LinearFVBCs]
  [inlet_u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left'
    variable = vel_x
    functor = ${bulk_u} #${bulk_u} #'fully_developed_velocity'
  []
  [inlet_h]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = h
    boundary = 'left'
    functor = 'h_from_p_T'
  []
  [inlet_T]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = T
    boundary = 'left'
    functor = ${T_in}
  []
  [outlet_p]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'right'
    variable = pressure
    functor = ${p_ref}
  []
  [outlet_h]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = h
    use_two_term_expansion = false
    boundary = 'right'
  []
  [outlet_u]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = vel_x
    use_two_term_expansion = false
    boundary = 'right'
  []
[]
[Functions]
  [source_func]
    type = ParsedFunction
    expression = ${q_source}
  []
  [T_analytical]
    type = ParsedFunction
    expression = ${fparse (-A_cp+sqrt(A_cp^2-2*B_cp*(-q_source/rho/bulk_u*L-A_cp*T_in-B_cp/2*T_in*T_in)))/B_cp}
  []
[]
[FunctorMaterials]
  [enthalpy_material]
    type = LinearFVEnthalpyFunctorMaterial
    pressure = ${p_ref}
    T_fluid = T
    h = h
    h_from_p_T_functor = h_from_p_T_functor
    T_from_p_h_functor = T_from_p_h_functor
  []
  [h_from_p_T_functor]
    type = ParsedFunctorMaterial
    property_name = 'h_from_p_T_functor'
    functor_names = 'T'
    expression = '${A_cp}*T+${B_cp}/2*(T^2)'
  []
  [T_from_p_h_functor]
    type = ParsedFunctorMaterial
    property_name = 'T_from_p_h_functor'
    functor_names = 'h'
    expression = '(-${A_cp}+sqrt(${A_cp}^2+2*h*${B_cp}))/${B_cp}'
  []
  [rho]
    type = ADParsedFunctorMaterial
    property_name = 'rho'
    functor_names = 'T'
    expression = ${rho}
  []
  [cp]
    type = ADParsedFunctorMaterial
    property_name = 'cp'
    functor_names = 'T'
    expression = '${A_cp}+${B_cp}*T'
  []
  [mu]
    type = ADParsedFunctorMaterial
    property_name = 'mu'
    functor_names = 'T'
    expression = '4.5e-3'
  []
  [k]
    type = ADParsedFunctorMaterial
    property_name = 'k'
    functor_names = 'T'
    expression = 0.7
  []
[]
[AuxKernels]
  [rho_out]
    type = FunctorAux
    functor = 'rho'
    variable = 'rho_var'
    execute_on = 'NONLINEAR'
  []
  [cp_out]
    type = FunctorAux
    functor = 'cp'
    variable = 'cp_var'
    execute_on = 'NONLINEAR'
  []
  [mu_out]
    type = FunctorAux
    functor = 'mu'
    variable = 'mu_var'
    execute_on = 'NONLINEAR'
  []
  [k_out]
    type = FunctorAux
    functor = 'k'
    variable = 'k_var'
    execute_on = 'NONLINEAR'
  []
  [T_from_h_functor_aux]
    type = FunctorAux
    functor = 'T_from_p_h'
    variable = 'T'
    execute_on = 'NONLINEAR'
  []
  [h_from_T_functor_aux]
    type = FunctorAux
    functor = 'h_from_p_T'
    variable = 'h_aux'
    execute_on = 'NONLINEAR'
  []
[]
[Postprocessors]
  [T_out_sim]
    type = ElementalVariableValue
    variable = T
    elementid = ${fparse nx-1}
  []
  [T_out_analytic]
    type = FunctionValuePostprocessor
    function = T_analytical
  []
[]
[Executioner]
  type = SIMPLE
  momentum_l_abs_tol = 1e-12
  pressure_l_abs_tol = 1e-12
  energy_l_abs_tol = 1e-12
  momentum_l_tol = 0
  pressure_l_tol = 0
  energy_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system'
  pressure_system = 'pressure_system'
  energy_system = 'energy_system'
  momentum_equation_relaxation = 0.7
  pressure_variable_relaxation = 0.3
  energy_equation_relaxation = 0.95
  num_iterations = 100
  pressure_absolute_tolerance = 1e-8
  momentum_absolute_tolerance = 1e-8
  energy_absolute_tolerance = 1e-6
  print_fields = false
  momentum_l_max_its = 200
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  energy_petsc_options_iname = '-pc_type -pc_hypre_type'
  energy_petsc_options_value = 'hypre boomeramg'
  continue_on_max_its = true
[]
[Outputs]
  [out]
    type = CSV
  []
[]
(test/tests/variables/linearfv/diffusion-1d-pp.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 50
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = analytic_solution
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '0.5*x'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1-x*x'
  []
[]
[Postprocessors]
  [average]
    type = ElementAverageValue
    variable = u
    execute_on = TIMESTEP_END
    outputs = csv
  []
  [min]
    type = ElementExtremeValue
    variable = u
    value_type = min
    execute_on = TIMESTEP_END
    outputs = csv
  []
  [max]
    type = ElementExtremeValue
    variable = u
    value_type = max
    execute_on = TIMESTEP_END
    outputs = csv
  []
  [num_dofs]
    type = NumDOFs
    execute_on = TIMESTEP_END
    outputs = csv
  []
  [elem_value]
    type = ElementalVariableValue
    variable = u
    elementid = 10
    execute_on = TIMESTEP_END
    outputs = csv
  []
  [point_value]
    type = PointValue
    variable = u
    point = '0.33333 0 0'
    execute_on = TIMESTEP_END
    outputs = csv
  []
[]
[VectorPostprocessors]
  [line-sample]
    type = LineValueSampler
    variable = u
    start_point = '0.13333 0 0'
    end_point = '0.766666 0 0'
    num_points = 9
    sort_by = x
    execute_on = TIMESTEP_END
    outputs = vpp_csv
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = TIMESTEP_END
  []
  [vpp_csv]
    type = CSV
    execute_on = TIMESTEP_END
  []
[]
(test/tests/linearfvbcs/robin/diffusion-1d-robin.i)
##################################################################
k     = 7.0 # diffusion coeff.
amp   = 3.6 # sinusoid amplitude, for u_exact
x_l   = ${fparse 0.0*pi} # domain bound (left)
x_r   = ${fparse 0.9*pi}     # domain bound (right)
alpha = 5.000 # robin BC coeff for gradient term
beta  = 2.000 # robin BC coeff for variable term
gamma = ${fparse (alpha*amp*cos(x_r) ) + (beta*amp*sin(x_r))} # RHS of Robin BC, applied at right boundary
npts = 2
##################################################################
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim  = 1
    nx   = ${npts}
    xmin = ${x_l}
    xmax = ${x_r}
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 0.0
  []
[]
[Functions]
  [u_exact]
    type = ParsedFunction
    expression = '${amp}*sin(x)'
  []
  [source_fn]
    type = ParsedFunction
    expression = '${fparse k*amp}*sin(x)'
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = ${k}
    use_nonorthogonal_correction = False
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_fn
  []
[]
[LinearFVBCs]
  [dir_r]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left"
    functor = 0
  []
  [rob_l]
    type = LinearFVAdvectionDiffusionFunctorRobinBC
    variable = u
    boundary = "right"
    alpha = ${alpha}
    beta  = ${beta}
    gamma = ${gamma}
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = u_exact
    execute_on = FINAL
  []
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 4
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-7
  linear_convergence = linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
(test/tests/linearfvkernels/diffusion/diffusion-2d-rz.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 1
    ymax = 0.5
  []
  coord_type = RZ
  rz_coord_axis = Y
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
    use_nonorthogonal_correction = true
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right top bottom"
    functor = analytic_solution
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '1+0.5*x*y'
  []
  [source_func]
    type = ParsedFunction
    expression = '-(-1.0*x^2*y*(1.5 - x^2) + x*(1.5 - x^2)*(-1.0*x*y - 2))/x - (-1.0*x^2*y*(1.5 - y^2) - 4*x*(1.5 - y^2)*(0.5*x*y + 1))/x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '(1.5-x*x)*(1.5-y*y)'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
  [exo]
    type = Exodus
    execute_on = FINAL
  []
[]
(test/tests/multiapps/linearfv_nonlinearfv/linearfv.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 6
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[AuxVariables]
  [diff_var]
    type = MooseVariableFVReal
    initial_condition = 2.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = diff_var
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = 1
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = 1
  []
[]
[MultiApps]
  inactive = 'nonlinear'
  [nonlinear]
    type = FullSolveMultiApp
    input_files = nonlinearfv.i
    execute_on = timestep_begin
    no_restore = true
  []
[]
[Transfers]
  inactive = 'from_nonlinear to_nonlinear'
  [from_nonlinear]
    type = MultiAppCopyTransfer
    from_multi_app = nonlinear
    source_variable = 'v'
    variable = 'diff_var'
    execute_on = timestep_begin
  []
  [to_nonlinear]
    type = MultiAppCopyTransfer
    to_multi_app = nonlinear
    source_variable = 'u'
    variable = 'diff_var'
    execute_on = timestep_begin
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_factor_mat_solver_type'
  petsc_options_value = 'lu mumps'
  fixed_point_rel_tol = 1e-10
[]
[Outputs]
  exodus = true
  execute_on = timestep_end
[]
(test/tests/variables/linearfv/diffusion-1d-aux.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 10
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[AuxVariables]
  [v_volume]
    type = MooseLinearVariableFVReal
    initial_condition = 50
  []
  [v_functor]
    type = MooseLinearVariableFVReal
    initial_condition = 25
  []
  [v_parsed]
    type = MooseLinearVariableFVReal
    initial_condition = 12.5
  []
[]
[AuxKernels]
  [volume]
    type = VolumeAux
    variable = v_volume
  []
  [functor]
    type = FunctorAux
    variable = v_functor
    functor = u
  []
  [parsed]
    type = ParsedAux
    variable = v_parsed
    coupled_variables = 'v_volume v_functor'
    expression = '0.5*v_volume+0.5*v_functor'
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = analytic_solution
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '0.5*x'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1-x*x'
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  exodus = true
  execute_on = TIMESTEP_END
[]
(test/tests/linearfvkernels/diffusion/diffusion-2d_neumann.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 1
    ymax = 0.5
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
    use_nonorthogonal_correction = false
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = analytic_solution
  []
  [neu_bottom]
    type = LinearFVAdvectionDiffusionFunctorNeumannBC
    variable = u
    boundary = "bottom"
    functor = analytic_solution_neumann_bottom
    diffusion_coeff = coeff_func
  []
  [neu_top]
    type = LinearFVAdvectionDiffusionFunctorNeumannBC
    variable = u
    boundary = "top"
    functor = analytic_solution_neumann_top
    diffusion_coeff = coeff_func
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '1+0.5*x*y'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*(1.5-y*y)+2*x*y*(1.5-y*y)+2*(1.5-x*x)+2*x*y*(1.5-x*x)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '(1.5-x*x)*(1.5-y*y)'
  []
  [analytic_solution_neumann_bottom]
    type = ParsedFunction
    expression = '-(1+0.5*x*y)*(1.5-x*x)*(-2*y)'
  []
  [analytic_solution_neumann_top]
    type = ParsedFunction
    expression = '(1+0.5*x*y)*(1.5-x*x)*(-2*y)'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvkernels/anisotropic-diffusion/anisotropic-diffusion-2d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 1
    ymax = 0.5
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVAnisotropicDiffusion
    variable = u
    diffusion_tensor = diffusivity_tensor
    use_nonorthogonal_correction = false
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right top bottom"
    functor = analytic_solution
  []
[]
[FunctorMaterials]
  [diff_tensor]
    type = GenericVectorFunctorMaterial
    prop_names = diffusivity_tensor
    prop_values = 'coeff_func_x coeff_func_y 0.0'
  []
[]
[Functions]
  [coeff_func_x]
    type = ParsedFunction
    expression = '1+0.5*x*y'
  []
  [coeff_func_y]
    type = ParsedFunction
    expression = '1+x*y'
  []
  [source_func]
    type = ParsedFunction
    expression = '(1.5-y*y)*(2+2*x*y)+(1.5-x*x)*(2+4*x*y)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '(1.5-x*x)*(1.5-y*y)'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvbcs/robin/advection-1d-robin.i)
##################################################################
c     = 0.1 # advection velocity (+x direction)
amp   = 7.0 # sinusoid amplitude, for u_exact
x_l   = ${fparse 0.0*pi} # domain bound (left)
x_r   = ${fparse pi}     # domain bound (right)
alpha = 5.000 # robin BC coeff for gradient term
beta  = 2.000 # robin BC coeff for variable term
u0 = 3 # some positive constant for the solution, > 1
gamma = ${fparse (-alpha*amp*sin(x_l)) + beta*amp*(u0 - cos(x_l))} # RHS of Robin BC, applied at left boundary
npts = 2
##################################################################
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim  = 1
    nx   = ${npts}
    xmin = ${x_l}
    xmax = ${x_r}
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 0.0
  []
[]
[Functions]
  [u_exact]
    type = ParsedFunction
    expression = '${amp}*(${u0} - cos(x))'
  []
  [source_fn]
    type = ParsedFunction
    expression = '${fparse c*amp}*sin(x)'
  []
[]
[LinearFVKernels]
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "${c} 0 0"
    advected_interp_method = average
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_fn
  []
[]
[LinearFVBCs]
  [rob_l]
    type = LinearFVAdvectionDiffusionFunctorRobinBC
    variable = u
    boundary = "left"
    alpha = ${alpha}
    beta  = ${beta}
    gamma = ${gamma}
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = true
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = u_exact
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 4
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-7
  petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = 'lu       NONZERO               1e-10'
  linear_convergence = linear
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvkernels/block-restriction/block-restricted-adr.i)
[Mesh]
  [cmg]
    type = CartesianMeshGenerator
    dim = 2
    dx = '0.1 1 0.1'
    dy = '0.1 0.5 0.1'
    ix = '1 2 1'
    iy = '1 1 1'
    subdomain_id = '1 1 1 1 2 3 1 1 1'
  []
  [transform]
    type = TransformGenerator
    input = cmg
    transform = TRANSLATE
    vector_value = '-0.1 -0.1 0.0'
  []
  [create_sides]
    type = SideSetsBetweenSubdomainsGenerator
    input = transform
    new_boundary = sides
    primary_block = 2
    paired_block = 1
  []
  [create_outlet]
    type = SideSetsBetweenSubdomainsGenerator
    input = create_sides
    new_boundary = outlet
    primary_block = 2
    paired_block = 3
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
    block = 2
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = diff_coeff_func
    use_nonorthogonal_correction = false
  []
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "0.5 0 0"
    advected_interp_method = average
  []
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  inactive = "outflow"
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "sides outlet"
    functor = analytic_solution
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = true
  []
[]
[Functions]
  [diff_coeff_func]
    type = ParsedFunction
    expression = '1.0+0.5*x*y'
  []
  [coeff_func]
    type = ParsedFunction
    expression = '1.0+1.0/(1+x*y)'
  []
  [source_func]
    type = ParsedFunction
    expression = '-1.0*x*pi*sin((1/2)*x*pi)*cos(2*y*pi) - 0.25*y*pi*sin(2*y*pi)*cos((1/2)*x*pi) + (1.0 + 1.0/(x*y + 1))*(sin((1/2)*x*pi)*sin(2*y*pi) + 1.5) + (17/4)*pi^2*(0.5*x*y + 1.0)*sin((1/2)*x*pi)*sin(2*y*pi) + 0.25*pi*sin(2*y*pi)*cos((1/2)*x*pi)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = 'sin((1/2)*x*pi)*sin(2*y*pi) + 1.5'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
    block = 2
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
    block = 2
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvkernels/advection/advection-2d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny= 1
    ymax = 0.5
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "0.5 0 0"
    advected_interp_method = upwind
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [inflow]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left top bottom"
    functor = analytic_solution
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = false
  []
[]
[Functions]
  [source_func]
    type = ParsedFunction
    expression = '0.5*pi*sin(2*y*pi)*cos(x*pi)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = 'sin(x*pi)*sin(2*y*pi) + 1.5'
  []
[]
[Postprocessors]
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = 'lu       NONZERO               1e-10'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/multisystem/restore_multiapp/sub.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 6
  []
[]
[Problem]
  nl_sys_names = 'v_sys'
  linear_sys_names = 'u_sys'
[]
[Variables]
  [v]
    type = MooseVariableFVReal
    initial_condition = 2.0
    solver_sys = v_sys
  []
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[FVKernels]
  [diffusion]
    type = FVDiffusion
    variable = v
    coeff = u
  []
  [source]
    type = FVBodyForce
    variable = v
    function = 3
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = v
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = 1
  []
[]
[FVBCs]
  [dir]
    type = FVFunctorDirichletBC
    variable = v
    boundary = "left right"
    functor = 2
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = 1
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 6
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  #type = Transient
  #steady_state_detection = true
  system_names = 'v_sys u_sys'
  l_abs_tol = 1e-12
  l_tol = 1e-10
  nl_abs_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  exodus = true
  #execute_on = timestep_end
[]
(test/tests/time_steppers/iteration_adaptive/adapt_linear_systems.i)
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 50
  ny = 2
  xmax = 5
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [time]
    type = LinearFVTimeDerivative
    variable = 'u'
  []
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = 5
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = 2
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = 12
  []
[]
[Executioner]
  type = Transient
  system_names = u_sys
  start_time = 0.0
  end_time = 19
  n_startup_steps = 2
  dtmax = 6.0
  [TimeStepper]
    type = IterationAdaptiveDT
    optimal_iterations = 8
    dt = 1.0
  []
  verbose = true
  linear_convergence = much_logic
[]
[Convergence]
  [force_grow]
    type = IterationCountConvergence
    min_iterations = 0
    max_iterations = 4
    converge_at_max_iterations = true
  []
  [force_shrink]
    type = IterationCountConvergence
    min_iterations = 12
    max_iterations = 13
    converge_at_max_iterations = true
  []
  [much_logic]
    type = ParsedConvergence
    convergence_expression = 'if(time < 5, force_grow, force_shrink)'
    symbol_names = 'time force_grow force_shrink'
    symbol_values = 'time force_grow force_shrink'
  []
[]
[Postprocessors]
  [_dt]
    type = TimestepSize
  []
  [time]
    type = TimePostprocessor
  []
[]
[Outputs]
  csv = true
[]
(test/tests/linearfvkernels/block-restriction/block-restricted-diffusion.i)
[Mesh]
  [cmg]
    type = CartesianMeshGenerator
    dim = 2
    dx = '0.1 1 0.1'
    dy = '0.1 0.5 0.1'
    ix = '1 2 1'
    iy = '1 1 1'
    subdomain_id = '1 1 1 1 2 3 1 1 1'
  []
  [transform]
    type = TransformGenerator
    input = cmg
    transform = TRANSLATE
    vector_value = '-0.1 -0.1 0.0'
  []
  [create_sides]
    type = SideSetsBetweenSubdomainsGenerator
    input = transform
    new_boundary = sides
    primary_block = 2
    paired_block = 1
  []
  [create_outlet]
    type = SideSetsBetweenSubdomainsGenerator
    input = create_sides
    new_boundary = outlet
    primary_block = 2
    paired_block = 3
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
    block = 2
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = diff_coeff_func
    use_nonorthogonal_correction = false
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "sides outlet"
    functor = analytic_solution
  []
[]
[Functions]
  [diff_coeff_func]
    type = ParsedFunction
    expression = '1.0+0.5*x*y'
  []
  [source_func]
    type = ParsedFunction
    expression = '-1.0*x*pi*sin(x*pi)*cos(2*y*pi) - 0.5*y*pi*sin(2*y*pi)*cos(x*pi) + 5*pi^2*(0.5*x*y + 1.0)*sin(x*pi)*sin(2*y*pi)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = 'sin(x*pi)*sin(2*y*pi) + 1.5'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
    block = 2
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
    block = 2
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/tag/tag-linearfv.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 3
  []
[]
[Problem]
  linear_sys_names = 'u_sys v_sys'
  extra_tag_matrices = 'mat_tag_u; mat_tag_v'
  extra_tag_vectors = 'vec_tag_u; vec_tag_v'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    initial_condition = 1.0
    solver_sys = u_sys
  []
  [v]
    type = MooseLinearVariableFVReal
    initial_condition = 0.5
    solver_sys = v_sys
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = 2.0
  []
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = 3.0
    matrix_tags = 'system mat_tag_u'
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = 60.0
    vector_tags = 'rhs vec_tag_u'
  []
  [diffusion_v]
    type = LinearFVDiffusion
    variable = v
    diffusion_coeff = 1.0
  []
  [reaction_v]
    type = LinearFVReaction
    variable = v
    coeff = 1.5
    matrix_tags = 'system mat_tag_v'
  []
  [source_v]
    type = LinearFVSource
    variable = v
    source_density = 20.0
    vector_tags = 'rhs vec_tag_v'
  []
[]
[LinearFVBCs]
  [left_u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = left
    functor = 1.0
  []
  [right_u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = right
    functor = 3.0
  []
  [left_v]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = v
    boundary = left
    functor = 1.0
  []
  [right_v]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = v
    boundary = right
    functor = 3.0
  []
[]
[AuxVariables]
  [soln_u_dof]
    type = MooseLinearVariableFVReal
  []
  [soln_u]
    type = MooseLinearVariableFVReal
  []
  [rhs_u_dof]
    type = MooseLinearVariableFVReal
  []
  [rhs_u]
    type = MooseLinearVariableFVReal
  []
  [vector_tag_u]
    type = MooseLinearVariableFVReal
  []
  [matrix_u_diag]
    type = MooseLinearVariableFVReal
  []
  [soln_v_dof]
    type = MooseLinearVariableFVReal
  []
  [soln_v]
    type = MooseLinearVariableFVReal
  []
  [rhs_v_dof]
    type = MooseLinearVariableFVReal
  []
  [rhs_v]
    type = MooseLinearVariableFVReal
  []
  [vector_tag_v]
    type = MooseLinearVariableFVReal
  []
  [matrix_v_diag]
    type = MooseLinearVariableFVReal
  []
[]
[AuxKernels]
  [soln_u_dof]
    type = TagVectorDofValueAux
    variable = soln_u_dof
    v = u
    vector_tag = 'solution'
  []
  [soln_u]
    type = TagVectorAux
    variable = soln_u
    v = u
    vector_tag = 'solution'
  []
  [rhs_u_dof]
    type = TagVectorDofValueAux
    variable = rhs_u_dof
    v = u
    vector_tag = 'rhs'
  []
  [rhs_u]
    type = TagVectorAux
    variable = rhs_u
    v = u
    vector_tag = 'rhs'
  []
  [extra_vector_u_dof]
    type = TagVectorDofValueAux
    variable = vector_tag_u
    v = u
    vector_tag = 'vec_tag_u'
  []
  [extra_vector_u]
    type = TagVectorAux
    variable = vector_tag_u
    v = u
    vector_tag = 'vec_tag_u'
  []
  [extra_matrix_u]
    type = TagMatrixAux
    variable = matrix_u_diag
    v = u
    matrix_tag = 'mat_tag_u'
  []
  [soln_v_dof]
    type = TagVectorDofValueAux
    variable = soln_v_dof
    v = v
    vector_tag = 'solution'
  []
  [soln_v]
    type = TagVectorAux
    variable = soln_v
    v = v
    vector_tag = 'solution'
  []
  [rhs_v_dof]
    type = TagVectorDofValueAux
    variable = rhs_v_dof
    v = v
    vector_tag = 'rhs'
  []
  [rhs_v]
    type = TagVectorAux
    variable = rhs_v
    v = v
    vector_tag = 'rhs'
  []
  [extra_vector_v_dof]
    type = TagVectorDofValueAux
    variable = vector_tag_v
    v = v
    vector_tag = 'vec_tag_v'
  []
  [extra_vector_v]
    type = TagVectorAux
    variable = vector_tag_v
    v = v
    vector_tag = 'vec_tag_v'
  []
  [extra_matrix_v]
    type = TagMatrixAux
    variable = matrix_v_diag
    v = v
    matrix_tag = 'mat_tag_v'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'LINEAR'
  system_names = 'u_sys v_sys'
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  exodus = true
  execute_on = timestep_end
[]
(test/tests/time_integrators/implicit-euler/ie-linearfv.i)
###########################################################
# This is a simple test with a time-dependent problem
# demonstrating the use of the TimeIntegrator system.
#
# Testing a solution that is second order in space
# and first order in time
#
# @Requirement F1.30
###########################################################
[Mesh]
  type = GeneratedMesh
  dim = 2
  xmin = -1
  xmax = 1
  ymin = -1
  ymax = 1
  nx = 10
  ny = 10
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 0.0
  []
[]
[Functions]
  [forcing_fn]
    type = ParsedFunction
    expression = ((x*x)+(y*y))-(4*t)
  []
  [exact_fn]
    type = ParsedFunction
    expression = t*((x*x)+(y*y))
  []
[]
[LinearFVKernels]
  [ie]
    type = LinearFVTimeDerivative
    variable = u
  []
  [diff]
    type = LinearFVDiffusion
    variable = u
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = forcing_fn
  []
[]
[LinearFVBCs]
  [all]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = '0 1 2 3'
    functor = exact_fn
  []
[]
[Postprocessors]
  [l2_err]
    type = ElementL2Error
    variable = u
    function = exact_fn
  []
[]
[Executioner]
  type = Transient
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
  # Test of the TimeIntegrator System
  scheme = 'implicit-euler'
  start_time = 0.0
  num_steps = 5
  dt = 0.25
[]
[Outputs]
  exodus = true
[]
(test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-1d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 2
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = diff_coeff_func
    use_nonorthogonal_correction = false
  []
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "0.5 0 0"
    advected_interp_method = average
  []
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  inactive = "outflow neumann"
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = analytic_solution
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = true
  []
  [neumann]
    type = LinearFVAdvectionDiffusionFunctorNeumannBC
    variable = u
    boundary = "left"
    functor = analytic_solution_neumann_left
    diffusion_coeff = diff_coeff_func
  []
[]
[Functions]
  [diff_coeff_func]
    type = ParsedFunction
    expression = '1+0.5*x'
  []
  [coeff_func]
    type = ParsedFunction
    expression = '1+1/(1+x)'
  []
  [source_func]
    type = ParsedFunction
    expression = '(1+1/(x+1))*(sin(pi/2*x)+1.5)+0.25*pi*pi*(0.5*x+1)*sin(pi/2*x)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = 'sin(pi/2*x)+1.5'
  []
  [analytic_solution_neumann_left]
    type = ParsedFunction
    expression = '-(1+0.5*x)*cos(pi/2*x)*pi/2'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type -ksp_rtol'
  petsc_options_value = 'hypre boomeramg 1e-10'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/mms/linear-segregated/2d-vortex/2d-vortex.i)
mu = 1
rho = 1
advected_interp_method = 'average'
[Problem]
  linear_sys_names = 'u_system v_system pressure_system'
  previous_nl_solution_required = true
[]
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 2
  []
[]
[UserObjects]
  [rc]
    type = RhieChowMassFlux
    u = vel_x
    v = vel_y
    pressure = pressure
    rho = ${rho}
    p_diffusion_kernel = p_diffusion
  []
[]
[Variables]
  [vel_x]
    type = MooseLinearVariableFVReal
    initial_condition = 0.0
    solver_sys = u_system
  []
  [vel_y]
    type = MooseLinearVariableFVReal
    solver_sys = v_system
    initial_condition = 0.0
  []
  [pressure]
    type = MooseLinearVariableFVReal
    solver_sys = pressure_system
    initial_condition = 0
  []
[]
[LinearFVKernels]
  [u_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    mu = ${mu}
    u = vel_x
    v = vel_y
    momentum_component = 'x'
    rhie_chow_user_object = 'rc'
    use_nonorthogonal_correction = false
  []
  [v_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    mu = ${mu}
    u = vel_x
    v = vel_y
    momentum_component = 'y'
    rhie_chow_user_object = 'rc'
    use_nonorthogonal_correction = false
  []
  [u_pressure]
    type = LinearFVMomentumPressure
    variable = vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [v_pressure]
    type = LinearFVMomentumPressure
    variable = vel_y
    pressure = pressure
    momentum_component = 'y'
  []
  [u_forcing]
    type = LinearFVSource
    variable = vel_x
    source_density = forcing_u
  []
  [v_forcing]
    type = LinearFVSource
    variable = vel_y
    source_density = forcing_v
  []
  [p_diffusion]
    type = LinearFVAnisotropicDiffusion
    variable = pressure
    diffusion_tensor = Ainv
    use_nonorthogonal_correction = false
    use_nonorthogonal_correction_on_boundary = false
  []
  [HbyA_divergence]
    type = LinearFVDivergence
    variable = pressure
    face_flux = HbyA
    force_boundary_execution = true
  []
[]
[LinearFVBCs]
  [no-slip-wall-u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left right top bottom'
    variable = vel_x
    functor = '0'
  []
  [no-slip-wall-v]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left right top bottom'
    variable = vel_y
    functor = '0'
  []
  [pressure-extrapolation]
    type = LinearFVExtrapolatedPressureBC
    boundary = 'left right top bottom'
    variable = pressure
    use_two_term_expansion = true
  []
[]
[Functions]
  [exact_u]
    type = ParsedFunction
    expression = 'x^2*(1-x)^2*(2*y-6*y^2+4*y^3)'
  []
  [exact_v]
    type = ParsedFunction
    expression = '-y^2*(1-y)^2*(2*x-6*x^2+4*x^3)'
  []
  [exact_p]
    type = ParsedFunction
    expression = 'x*(1-x)'
  []
  [forcing_u]
    type = ParsedFunction
    expression = '-4*mu*(-1+2*y)*(y^2-6*x*y^2+6*x^2*y^2-y+6*x*y-6*x^2*y+3*x^2-6*x^3+3*x^4)+1-2*x+rho*4*x^3'
            '*y^2*(2*y^2-2*y+1)*(y-1)^2*(-1+2*x)*(x-1)^3'
    symbol_names = 'mu rho'
    symbol_values = '${mu} ${rho}'
  []
  [forcing_v]
    type = ParsedFunction
    expression = '4*mu*(-1+2*x)*(x^2-6*y*x^2+6*x^2*y^2-x+6*x*y-6*x*y^2+3*y^2-6*y^3+3*y^4)+rho*4*y^3*x^2*(2'
            '*x^2-2*x+1)*(x-1)^2*(-1+2*y)*(y-1)^3'
    symbol_names = 'mu rho'
    symbol_values = '${mu} ${rho}'
  []
[]
[Executioner]
  type = SIMPLE
  momentum_l_abs_tol = 1e-8
  pressure_l_abs_tol = 1e-8
  momentum_l_tol = 0
  pressure_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  momentum_equation_relaxation = 0.8
  pressure_variable_relaxation = 0.3
  num_iterations = 2000
  pressure_absolute_tolerance = 1e-8
  momentum_absolute_tolerance = 1e-8
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  print_fields = false
  pin_pressure = true
  pressure_pin_value = 0.25
  pressure_pin_point = '0.5 0.5 0.0'
[]
[Outputs]
  exodus = true
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    outputs = 'csv'
    execute_on = FINAL
  []
  [L2u]
    type = ElementL2FunctorError
    approximate = vel_x
    exact = exact_u
    outputs = 'csv'
    execute_on = FINAL
  []
  [L2v]
    type = ElementL2FunctorError
    approximate = vel_y
    exact = exact_v
    outputs = 'csv'
    execute_on = FINAL
  []
  [L2p]
    approximate = pressure
    exact = exact_p
    type = ElementL2FunctorError
    outputs = 'csv'
    execute_on = FINAL
  []
[]
(test/tests/outputs/debug/show_execution_linear_fv_elemental.i)
[Mesh]
  [gen_mesh]
    type = GeneratedMeshGenerator
    dim = 1
    xmin = 0
    xmax = 10
    nx = 50
  []
  [left]
    type = ParsedSubdomainMeshGenerator
    input = 'gen_mesh'
    combinatorial_geometry = 'x < 0.5'
    block_id = '1'
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
  []
[]
[LinearFVKernels]
  [reaction_1]
    type = LinearFVReaction
    variable = u
    coeff = 1.5
    block = 0
  []
  [reaction_2]
    type = LinearFVReaction
    variable = u
    coeff = 2.5
    block = 1
  []
  [source_1]
    type = LinearFVSource
    variable = u
    source_density = 3.5
    block = 0
  []
  [source_2]
    type = LinearFVSource
    variable = u
    source_density = 4.5
    block = 1
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
[]
[Debug]
  show_execution_order = 'NONLINEAR'
[]
(test/tests/linearfvkernels/diffusion/diffusion-1d_neumann.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 2
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "right"
    functor = analytic_solution
  []
  [neu]
    type = LinearFVAdvectionDiffusionFunctorNeumannBC
    variable = u
    boundary = "left"
    functor = analytic_solution_neumann
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '0.5*x'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1-x*x'
  []
  [analytic_solution_neumann]
    type = ParsedFunction
    expression = '-(0.5*x)*(-2*x)'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-2d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 1
    ymax = 0.5
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = diff_coeff_func
    use_nonorthogonal_correction = false
  []
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "0.5 0 0"
    advected_interp_method = average
  []
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  inactive = "outflow neumann"
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right top bottom"
    functor = analytic_solution
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = true
  []
  [neumann]
    type = LinearFVAdvectionDiffusionFunctorNeumannBC
    variable = u
    boundary = "top"
    functor = analytic_solution_neumann_top
    diffusion_coeff = diff_coeff_func
  []
[]
[Functions]
  [diff_coeff_func]
    type = ParsedFunction
    expression = '1.0+0.5*x*y'
  []
  [coeff_func]
    type = ParsedFunction
    expression = '1.0+1.0/(1+x*y)'
  []
  [source_func]
    type = ParsedFunction
    expression = '-1.0*x*pi*sin((1/2)*x*pi)*cos(2*y*pi) - 0.25*y*pi*sin(2*y*pi)*cos((1/2)*x*pi) + (1.0 + 1.0/(x*y + 1))*(sin((1/2)*x*pi)*sin(2*y*pi) + 1.5) + (17/4)*pi^2*(0.5*x*y + 1.0)*sin((1/2)*x*pi)*sin(2*y*pi) + 0.25*pi*sin(2*y*pi)*cos((1/2)*x*pi)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = 'sin((1/2)*x*pi)*sin(2*y*pi) + 1.5'
  []
  [analytic_solution_neumann_top]
    type = ParsedFunction
    expression = '(1.0+0.5*x*y)*sin((1/2)*x*pi)*cos(2*y*pi)*2*pi'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/linearfvkernels/diffusion/diffusion-1d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 2
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = analytic_solution
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '0.5*x'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1-x*x'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/mms/linear-segregated/2d-vortex/spacedependent_mu/sl.i)
rho = 1
advected_interp_method = 'average'
[Problem]
  linear_sys_names = 'u_system v_system pressure_system'
  previous_nl_solution_required = true
[]
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny = 2
  []
[]
[UserObjects]
  [rc]
    type = RhieChowMassFlux
    u = vel_x
    v = vel_y
    pressure = pressure
    rho = ${rho}
    p_diffusion_kernel = p_diffusion
  []
[]
[Variables]
  [vel_x]
    type = MooseLinearVariableFVReal
    initial_condition = 0.0
    solver_sys = u_system
  []
  [vel_y]
    type = MooseLinearVariableFVReal
    solver_sys = v_system
    initial_condition = 0.0
  []
  [pressure]
    type = MooseLinearVariableFVReal
    solver_sys = pressure_system
    initial_condition = 0
  []
[]
[LinearFVKernels]
  [u_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    mu = 'mu'
    u = vel_x
    v = vel_y
    momentum_component = 'x'
    rhie_chow_user_object = 'rc'
    use_nonorthogonal_correction = false
    use_deviatoric_terms = false
  []
  [v_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    mu = 'mu'
    u = vel_x
    v = vel_y
    momentum_component = 'y'
    rhie_chow_user_object = 'rc'
    use_nonorthogonal_correction = false
    use_deviatoric_terms = false
  []
  [u_pressure]
    type = LinearFVMomentumPressure
    variable = vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [v_pressure]
    type = LinearFVMomentumPressure
    variable = vel_y
    pressure = pressure
    momentum_component = 'y'
  []
  [u_forcing]
    type = LinearFVSource
    variable = vel_x
    source_density = forcing_u
  []
  [v_forcing]
    type = LinearFVSource
    variable = vel_y
    source_density = forcing_v
  []
  [p_diffusion]
    type = LinearFVAnisotropicDiffusion
    variable = pressure
    diffusion_tensor = Ainv
    use_nonorthogonal_correction = false
    use_nonorthogonal_correction_on_boundary = false
  []
  [HbyA_divergence]
    type = LinearFVDivergence
    variable = pressure
    face_flux = HbyA
    force_boundary_execution = true
  []
[]
[LinearFVBCs]
  [no-slip-wall-u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left right top bottom'
    variable = vel_x
    functor = '0'
  []
  [no-slip-wall-v]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left right top bottom'
    variable = vel_y
    functor = '0'
  []
  [pressure-extrapolation]
    type = LinearFVExtrapolatedPressureBC
    boundary = 'left right top bottom'
    variable = pressure
    use_two_term_expansion = true
  []
[]
[Functions]
  [exact_u]
    type = ParsedFunction
    expression = 'x^2*(1-x)^2*(2*y-6*y^2+4*y^3)'
  []
  [exact_v]
    type = ParsedFunction
    expression = '-y^2*(1-y)^2*(2*x-6*x^2+4*x^3)'
  []
  [exact_p]
    type = ParsedFunction
    expression = 'x*(1-x)'
  []
  [mu]
    type = ParsedFunction
    expression = '1+(x-1)*x*(y-1)*y'
  []
  [forcing_u]
    type = ParsedFunction
    expression = '-(2*x-1)*y*(y-1)*(2*x-6*x^2+4*x^3)*(2*y-6*y^2+4*y^3)'
                 '-(1+x*(x-1)*y*(y-1))*(2*y-6*y^2+4*y^3)*(2-12*x+12*x^2)'
                 '-(2*y-1)*x*(x-1)*(x^2*(1-x)^2*(2-12*y+12*y^2))'
                 '-(1+x*(x-1)*y*(y-1))*(x^2*(1-x)^2*(-12+24*y))'
                 '+1-2*x+rho*4*x^3*y^2*(2*y^2-2*y+1)*(y-1)^2*(-1+2*x)*(x-1)^3'
    symbol_names = 'rho'
    symbol_values = '${rho}'
  []
  [forcing_v]
    type = ParsedFunction
    expression = '(2*y-1)*x*(x-1)*(2*y-6*y^2+4*y^3)*(2*x-6*x^2+4*x^3)'
                 '+(1+x*(x-1)*y*(y-1))*(2-12*y+12*y^2)*(2*x-6*x^2+4*x^3)'
                 '+(2*x-1)*y*(y-1)*(y^2*(1-y)^2*(2-12*x+12*x^2))'
                 '+(1+x*(x-1)*y*(y-1))*(y^2*(1-y)^2*(-12+24*x))'
                 '+rho*4*y^3*x^2*(2*x^2-2*x+1)*(x-1)^2*(-1+2*y)*(y-1)^3'
    symbol_names = 'rho'
    symbol_values = '${rho}'
  []
  [forcing_u_deviatoric]
    type = ParsedFunction
    expression = '-2*(2*x-1)*y*(y-1)*(2*x-6*x^2+4*x^3)*(2*y-6*y^2+4*y^3)'
                 '-2*(1+x*(x-1)*y*(y-1))*(2*y-6*y^2+4*y^3)*(2-12*x+12*x^2)'
                 '-(2*y-1)*x*(x-1)*(x^2*(1-x)^2*(2-12*y+12*y^2)-y^2*(1-y)^2*(2-12*x+12*x^2))'
                 '-(1+x*(x-1)*y*(y-1))*(x^2*(1-x)^2*(-12+24*y)-(2*y-6*y^2+4*y^3)*(2-12*x+12*x^2))'
                 '+1-2*x+rho*4*x^3*y^2*(2*y^2-2*y+1)*(y-1)^2*(-1+2*x)*(x-1)^3'
    symbol_names = 'rho'
    symbol_values = '${rho}'
  []
  [forcing_v_deviatoric]
    type = ParsedFunction
    expression = '2*(2*y-1)*x*(x-1)*(2*y-6*y^2+4*y^3)*(2*x-6*x^2+4*x^3)'
                 '+2*(1+x*(x-1)*y*(y-1))*(2-12*y+12*y^2)*(2*x-6*x^2+4*x^3)'
                 '-(2*x-1)*y*(y-1)*(x^2*(1-x)^2*(2-12*y+12*y^2)-y^2*(1-y)^2*(2-12*x+12*x^2))'
                 '-(1+x*(x-1)*y*(y-1))*(-y^2*(1-y)^2*(-12+24*x)+(2*x-6*x^2+4*x^3)*(2-12*y+12*y^2))'
                 '+rho*4*y^3*x^2*(2*x^2-2*x+1)*(x-1)^2*(-1+2*y)*(y-1)^3'
    symbol_names = 'rho'
    symbol_values = '${rho}'
  []
[]
[Executioner]
  type = SIMPLE
  momentum_l_abs_tol = 1e-14
  pressure_l_abs_tol = 1e-14
  momentum_l_max_its = 30
  pressure_l_max_its = 30
  momentum_l_tol = 0.0
  pressure_l_tol = 0.0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  momentum_equation_relaxation = 0.8
  pressure_variable_relaxation = 0.3
  num_iterations = 2000
  pressure_absolute_tolerance = 1e-8
  momentum_absolute_tolerance = 1e-8
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  print_fields = false
  pin_pressure = true
  pressure_pin_value = 0.25
  pressure_pin_point = '0.5 0.5 0.0'
[]
[Outputs]
  exodus = true
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    outputs = 'csv'
    execute_on = FINAL
  []
  [L2u]
    type = ElementL2FunctorError
    approximate = vel_x
    exact = exact_u
    outputs = 'csv'
    execute_on = FINAL
  []
  [L2v]
    type = ElementL2FunctorError
    approximate = vel_y
    exact = exact_v
    outputs = 'csv'
    execute_on = FINAL
  []
  [L2p]
    approximate = pressure
    exact = exact_p
    type = ElementL2FunctorError
    outputs = 'csv'
    execute_on = FINAL
  []
[]
(test/tests/linearfvbcs/robin/advection-2d-robin.i)
##################################################################
c     = 0.01 # advection velocity (+x direction)
amp   = 1.0 # sinusoid amplitude, for u_exact
u0    = 1.2 # any positive constant > 1.0
x_l   = ${fparse 0.0*pi}
x_r   = ${fparse pi}
y_l   = ${fparse 0.0*pi}
y_r   = ${fparse 1.0*pi}
alpha = 5.000 # robin BC coeff for gradient term
beta  = 2.000 # robin BC coeff for variable term
nx = 2
ny = 2
##################################################################
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim  = 2
    nx   = ${nx}
    ny   = ${ny}
    xmin = ${x_l}
    xmax = ${x_r}
    ymin = ${y_l}
    ymax = ${y_r}
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 0.01
  []
[]
[Functions]
  [u_exact]
    type = ParsedFunction
    expression = '${amp} * (${u0} - cos(x)) * sin(y)'
  []
  [source_fn]
    type = ParsedFunction
    expression = '${fparse c*amp} * sin(x) * sin(y)'
  []
  [gamma_fn]
    type = ParsedFunction
    expression = '(${fparse -amp*alpha}*sin(x)*sin(y)) + (${beta} * u_e)'
    symbol_names = 'u_e'
    symbol_values = 'u_exact'
  []
[]
[LinearFVKernels]
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "${c} 0 0"
    advected_interp_method = average
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_fn
  []
[]
[LinearFVBCs]
  [rob_l]
    type = LinearFVAdvectionDiffusionFunctorRobinBC
    variable = u
    boundary = "left"
    alpha = ${alpha}
    beta  = ${beta}
    gamma = gamma_fn
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "right"
    use_two_term_expansion = true
  []
  [top]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "top"
    functor = 0.0
  []
  [bottom]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "bottom"
    functor = 0.0
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = u_exact
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 4
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-7
  petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = 'lu       NONZERO               1e-10'
  linear_convergence = linear
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/multisystem/picard/linearfv_nonlinearfv/same_input.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 6
  []
[]
[Problem]
  nl_sys_names = 'v_sys'
  linear_sys_names = 'u_sys'
[]
[Variables]
  [v]
    type = MooseVariableFVReal
    initial_condition = 2.0
    solver_sys = v_sys
  []
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[FVKernels]
  [diffusion]
    type = FVDiffusion
    variable = v
    coeff = u
  []
  [source]
    type = FVBodyForce
    variable = v
    function = 3
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = v
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = 1
  []
[]
[FVBCs]
  [dir]
    type = FVFunctorDirichletBC
    variable = v
    boundary = "left right"
    functor = 2
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = 1
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 6
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = 'v_sys u_sys'
  l_abs_tol = 1e-12
  l_tol = 1e-10
  nl_abs_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  exodus = true
  execute_on = timestep_end
[]
(test/tests/linearfvbcs/robin/diffusion-2d-robin.i)
#################################################################
k     = 7.0 # diffusion coeff.
amp   = 3.6 # sinusoid amplitude, for u_exact
x1   = ${fparse 0.1*pi}
x2   = ${fparse 1.0*pi}
y1   = ${fparse 0.0*pi}
y2   = ${fparse 1.0*pi}
alpha = 2.000 # robin BC coeff for gradient term
beta  = 5.000 # robin BC coeff for variable term
nx = 2
ny = 2
##################################################################
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim  = 2
    nx   = ${nx}
    ny   = ${ny}
    xmin = ${x1}
    xmax = ${x2}
    ymin = ${y1}
    ymax = ${y2}
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 0.0
  []
[]
[Functions]
  [u_exact]
    type = ParsedFunction
    expression = '(${amp}*sin(x)*sin(y))'
  []
  [source_fn]
    type = ParsedFunction
    expression = '${fparse k*amp}*2.0*sin(x)*sin(y)'
  []
  [gamma_fn]
    type = ParsedFunction
    expression = '${fparse -amp*alpha}*cos(x)*sin(y) + ${beta} * u_e'
    symbol_names = 'u_e'
    symbol_values = 'u_exact'
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = ${k}
    use_nonorthogonal_correction = true
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_fn
  []
[]
[LinearFVBCs]
  [right]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "right"
    functor = u_exact
  []
  [top]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "top"
    functor = u_exact
  []
  [bottom]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "bottom"
    functor = u_exact
  []
  [robin_left]
    type = LinearFVAdvectionDiffusionFunctorRobinBC
    variable = u
    boundary = "left"
    alpha = ${alpha}
    beta  = ${beta}
    gamma = gamma_fn
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = u_exact
    execute_on = FINAL
  []
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 4
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-7
  multi_system_fixed_point = true
  multi_system_fixed_point_convergence = linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
  linear_convergence = linear
[]
(test/tests/multisystem/restore_multiapp/parent.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 6
  []
[]
[Problem]
  nl_sys_names = 'v_sys'
  linear_sys_names = 'u_sys'
[]
[Variables]
  [v]
    type = MooseVariableFVReal
    initial_condition = 2.0
    solver_sys = v_sys
  []
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[FVKernels]
  [diffusion]
    type = FVDiffusion
    variable = v
    coeff = u
  []
  [source]
    type = FVBodyForce
    variable = v
    function = 3
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = v
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = 1
  []
[]
[FVBCs]
  [dir]
    type = FVFunctorDirichletBC
    variable = v
    boundary = "left right"
    functor = 2
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = 1
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 6
    converge_at_max_iterations = true
  []
[]
[MultiApps]
  [sub]
    type = FullSolveMultiApp
    input_files = sub.i
    execute_on = TIMESTEP_END
    keep_solution_during_restore = true
  []
[]
[Executioner]
  type = Steady
  system_names = 'v_sys u_sys'
  l_abs_tol = 1e-12
  l_tol = 1e-10
  nl_abs_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
  fixed_point_max_its = 2
  fixed_point_min_its = 2
  #disable_fixed_point_residual_norm_check = true
  accept_on_max_fixed_point_iteration = true
[]
(test/tests/linearfvkernels/reaction/reaction-1d.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 10
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '1+sin(x)'
  []
  [source_func]
    type = ParsedFunction
    expression = '(1+sin(x))*(1+cos(x))'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1+cos(x)'
  []
[]
[Postprocessors]
  [l2error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = TIMESTEP_END
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
[]
[Outputs]
  [exodus]
    type = Exodus
    execute_on = TIMESTEP_END
  []
[]
(test/tests/linearfvkernels/advection/advection-2d-rz.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 2
    ny= 1
    ymax = 0.5
  []
  coord_type = RZ
  rz_coord_axis = Y
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [advection]
    type = LinearFVAdvection
    variable = u
    velocity = "0.0 0.5 0"
    advected_interp_method = average
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [inflow]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right bottom"
    functor = analytic_solution
  []
  [outflow]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = u
    boundary = "top"
    use_two_term_expansion = true
  []
[]
[Functions]
  [source_func]
    type = ParsedFunction
    expression = '1.0*pi*sin(x*pi)*cos(2*y*pi)'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = 'sin(x*pi)*sin(2*y*pi) + 1.5'
  []
[]
[Postprocessors]
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = FINAL
  []
  [h]
    type = AverageElementSize
    execute_on = FINAL
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 10
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = 'lu       NONZERO               1e-10'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = FINAL
  []
[]
(test/tests/outputs/debug/show_execution_linear_fv_flux.i)
[Mesh]
  [gen_mesh]
    type = GeneratedMeshGenerator
    dim = 1
    xmin = 0
    xmax = 10
    nx = 50
  []
  [left]
    type = ParsedSubdomainMeshGenerator
    input = 'gen_mesh'
    combinatorial_geometry = 'x < 0.5'
    block_id = '1'
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
  []
[]
[LinearFVKernels]
  [diffusion_1]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = 1.5
    block = 0
  []
  [diffusion_2]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = 2.5
    block = 1
  []
  [source_1]
    type = LinearFVSource
    variable = u
    source_density = 3.5
    block = 0
  []
  [source_2]
    type = LinearFVSource
    variable = u
    source_density = 4.5
    block = 1
  []
[]
[LinearFVBCs]
  [left_bc]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = 'left right'
    functor = 0
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
[]
[Debug]
  show_execution_order = 'NONLINEAR'
[]
(test/tests/linearfvkernels/block-restriction/block-restricted-diffusion-react.i)
source=1
diff_coeff=2
reac_coeff=3
[Mesh]
  [cmg]
    type = CartesianMeshGenerator
    dim = 1
    dx = '0.5 0.5'
    ix = '20 20'
    subdomain_id = '1 2'
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = ${diff_coeff}
    use_nonorthogonal_correction = false
    block = 1
  []
  [reaction]
    type = LinearFVReaction
    variable = u
    coeff = ${reac_coeff}
    block = 2
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = ${source}
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left"
    functor = 0
  []
[]
[Functions]
  [analytic_solution]
    type = ParsedFunction
    expression = 'if(x<0.5, -x*x*S/2/D+(S/C+0.5*0.5/2/D*S)/0.5*x, S/C)'
    symbol_names = 'S D C'
    symbol_values = '${source} ${diff_coeff} ${reac_coeff}'
  []
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    execute_on = TIMESTEP_END
    block = 2
  []
  [error]
    type = ElementL2FunctorError
    approximate = u
    exact = analytic_solution
    execute_on = TIMESTEP_END
    block = 2
  []
[]
[Convergence]
  [linear]
    type = IterationCountConvergence
    max_iterations = 1
    converge_at_max_iterations = true
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  multi_system_fixed_point=true
  multi_system_fixed_point_convergence=linear
  petsc_options_iname = '-pc_type -pc_hypre_type -ksp_rtol'
  petsc_options_value = 'hypre boomeramg 1e-10'
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = TIMESTEP_END
  []
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h_fp.i)
L = 30
nx = 600
bulk_u = 0.01
p_ref = 101325.0
T_in = 860.
q_source = 20000.
advected_interp_method = 'upwind'
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    xmin = 0
    xmax = ${L}
    nx = ${nx}
  []
  allow_renumbering = false
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
  advected_interp_method = ${advected_interp_method}
  u = vel_x
[]
[Problem]
  linear_sys_names = 'u_system pressure_system energy_system'
  previous_nl_solution_required = true
[]
[UserObjects]
  [rc]
    type = RhieChowMassFlux
    u = vel_x
    pressure = pressure
    rho = 'rho'
    p_diffusion_kernel = p_diffusion
  []
[]
[Variables]
  [vel_x]
    type = MooseLinearVariableFVReal
    solver_sys = u_system
    initial_condition = ${bulk_u}
  []
  [pressure]
    type = MooseLinearVariableFVReal
    solver_sys = pressure_system
    initial_condition = ${p_ref}
  []
  [h]
    type = MooseLinearVariableFVReal
    solver_sys = energy_system
    initial_condition = ${fparse 860.*240.}
  []
[]
[AuxVariables]
  [rho_var]
    type = MooseLinearVariableFVReal
  []
  [cp_var]
    type = MooseLinearVariableFVReal
  []
  [mu_var]
    type = MooseLinearVariableFVReal
  []
  [k_var]
    type = MooseLinearVariableFVReal
  []
  [T]
    type = MooseLinearVariableFVReal
    initial_condition = ${T_in}
  []
  [h_aux]
    type = MooseLinearVariableFVReal
  []
[]
[LinearFVKernels]
  [u_advection_stress]
    type = LinearWCNSFVMomentumFlux
    variable = vel_x
    mu = 'mu'
    momentum_component = 'x'
    use_nonorthogonal_correction = false
  []
  [u_pressure]
    type = LinearFVMomentumPressure
    variable = vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [p_diffusion]
    type = LinearFVAnisotropicDiffusion
    variable = pressure
    diffusion_tensor = Ainv
    use_nonorthogonal_correction = false
  []
  [HbyA_divergence]
    type = LinearFVDivergence
    variable = pressure
    face_flux = HbyA
    force_boundary_execution = true
  []
  [temp_advection]
    type = LinearFVEnergyAdvection
    variable = h
  []
  [source]
    type = LinearFVSource
    variable = h
    source_density = source_func
  []
[]
[LinearFVBCs]
  [inlet_u]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'left'
    variable = vel_x
    functor = ${bulk_u}
  []
  [inlet_h]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = h
    boundary = 'left'
    functor = 'h_from_p_T'
  []
  [inlet_T]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = T
    boundary = 'left'
    functor = ${T_in}
  []
  [outlet_p]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    boundary = 'right'
    variable = pressure
    functor = ${p_ref}
  []
  [outlet_h]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = h
    use_two_term_expansion = false
    boundary = 'right'
  []
  [outlet_u]
    type = LinearFVAdvectionDiffusionOutflowBC
    variable = vel_x
    use_two_term_expansion = false
    boundary = 'right'
  []
[]
[FluidProperties]
  [lead]
    type = LeadFluidProperties
  []
[]
[FunctorMaterials]
  [enthalpy_material]
    type = LinearFVEnthalpyFunctorMaterial
    pressure = ${p_ref}
    T_fluid = T
    h = h
    fp = lead
  []
  [fluid_props_to_mat_props]
    type = GeneralFunctorFluidProps
    fp = lead
    pressure = ${p_ref}
    T_fluid = 'T'
    speed = 1
    porosity = 1
    characteristic_length = 1
  []
  [source_func]
    type = ADParsedFunctorMaterial
    property_name = source_func
    functor_names = 'rho'
    expression = ${q_source}
  []
[]
[AuxKernels]
  [rho_out]
    type = FunctorAux
    functor = 'rho'
    variable = 'rho_var'
    execute_on = 'NONLINEAR'
  []
  [cp_out]
    type = FunctorAux
    functor = 'cp'
    variable = 'cp_var'
    execute_on = 'NONLINEAR'
  []
  [mu_out]
    type = FunctorAux
    functor = 'mu'
    variable = 'mu_var'
    execute_on = 'NONLINEAR'
  []
  [k_out]
    type = FunctorAux
    functor = 'k'
    variable = 'k_var'
    execute_on = 'NONLINEAR'
  []
  [T_from_h_functor_aux]
    type = FunctorAux
    functor = 'T_from_p_h'
    variable = 'T'
    execute_on = 'NONLINEAR'
  []
  [h_from_T_functor_aux]
    type = FunctorAux
    functor = 'h_from_p_T'
    variable = 'h_aux'
    execute_on = 'NONLINEAR'
  []
[]
[Postprocessors]
  [T_out_sim]
    type = ElementalVariableValue
    variable = T
    elementid = ${fparse nx-1}
  []
[]
[Executioner]
  type = SIMPLE
  momentum_l_abs_tol = 1e-12
  pressure_l_abs_tol = 1e-12
  energy_l_abs_tol = 1e-12
  momentum_l_tol = 0
  pressure_l_tol = 0
  energy_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system'
  pressure_system = 'pressure_system'
  energy_system = 'energy_system'
  momentum_equation_relaxation = 0.7
  pressure_variable_relaxation = 0.3
  energy_equation_relaxation = 0.95
  num_iterations = 100
  pressure_absolute_tolerance = 1e-8
  momentum_absolute_tolerance = 1e-8
  energy_absolute_tolerance = 1e-6
  print_fields = false
  momentum_l_max_its = 200
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  energy_petsc_options_iname = '-pc_type -pc_hypre_type'
  energy_petsc_options_value = 'hypre boomeramg'
  continue_on_max_its = true
[]
[Outputs]
  [out]
    type = CSV
  []
[]
(test/tests/transfers/multiapp_copy_transfer/linear_sys_to_aux/linear_sub.i)
[Mesh]
  [gmg]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 10
  []
[]
[Problem]
  linear_sys_names = 'u_sys'
[]
[Variables]
  [u]
    type = MooseLinearVariableFVReal
    solver_sys = 'u_sys'
    initial_condition = 1.0
  []
[]
[LinearFVKernels]
  [diffusion]
    type = LinearFVDiffusion
    variable = u
    diffusion_coeff = coeff_func
  []
  [source]
    type = LinearFVSource
    variable = u
    source_density = source_func
  []
[]
[LinearFVBCs]
  [dir]
    type = LinearFVAdvectionDiffusionFunctorDirichletBC
    variable = u
    boundary = "left right"
    functor = analytic_solution
  []
[]
[Functions]
  [coeff_func]
    type = ParsedFunction
    expression = '0.5*x'
  []
  [source_func]
    type = ParsedFunction
    expression = '2*x'
  []
  [analytic_solution]
    type = ParsedFunction
    expression = '1-x*x'
  []
[]
[Executioner]
  type = Steady
  system_names = u_sys
  l_tol = 1e-10
  petsc_options_iname = '-pc_type -pc_hypre_type'
  petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
  exodus = true
  execute_on = TIMESTEP_END
[]