- momentum_componentThe component of the momentum equation that this kernel applies to.C++ Type:MooseEnum Controllable:No Description:The component of the momentum equation that this kernel applies to. 
- rhie_chow_user_objectThe rhie-chow user-objectC++ Type:UserObjectName Controllable:No Description:The rhie-chow user-object 
- rhoThe density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:The density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- variableThe name of the variable that this residual object operates onC++ Type:NonlinearVariableName Unit:(no unit assumed) Controllable:No Description:The name of the variable that this residual object operates on 
PINSFVMomentumFriction
This kernel adds the friction term to the porous media Navier Stokes momentum equations. This kernel must be used with the canonical PINSFV variable set, e.g. pressure and superficial velocity, and supports Darcy and Forchheimer friction models in two flavors:
Standard friction formulation
Set parameter: "standard_friction_formulation" = 'true'
Darcy drag model (1) Forchheimer drag model (2)
Simplified friction formulation
Set parameter: "standard_friction_formulation" = 'false'
Darcy drag model (3) Forchheimer drag model (4)
where is the i-th component of the friction force (denoted by in Finite Volume Incompressible Porous media Navier Stokes, Eq. (1)), the friction factor, which may be anisotropic, the fluid dynamic viscosity, the fluid density, and the i-th component of the fluid superficial velocity. We have used a negative sign to match the notation used in Finite Volume Incompressible Porous media Navier Stokes, Eq. (1) where the friction force is on the right-hand-side of the equation. When moved to the left-hand side, which is done when setting up a Newton scheme, the term becomes positive which is what is shown in the source code itself. Darcy and Forchheimer terms represent fundamentally different friction effects. Darcy is meant to represent viscous effects and as shown in Eq. (1),Eq. (3), it has a linear dependence on the fluid velocity. Meanwhile, Forchheimer is meant to represent inertial effects and as shown in Eq. (2), Eq. (4) it has a quadratic dependence on velocity.
For the non-porous medium version of the above equations set parameter "is_porous_medium"  to false. (epsilon = 1)
Computation of friction factors and pre-factors
To outline how friction factors for Darcy and Forchheimer may be calculated, let's consider a specific example. We'll draw from the Ergun equation, which is outlined here. Let's consider the form:
where is the bed length, is the fluid dynamic viscosity and is representative of the diameter of the pebbles in the pebble bed. We can divide the equation through by , recognize that denotes such that , multiply the equation through by , move all terms to the left-hand-side, and do some term manipulation in order to yield:
If we define the hydraulic diameter as , then the above equation can be rewritten as:
Let's introduce the interstitial fluid velocity to rewrite the above equation as:
Then dividing through by :
(5)
We are now very close the forms for Darcy and Forchheimer espoused by Holzmann and SimScale which is:
(6)
Looking at Eq. (6) we can rearrange Eq. (5):
and arrive at the Ergun expression for the Darcy coefficient:
and the Ergun expression for the Forchheimer coefficient:
where we have made the multiplication explicit to make the 1.75 factor from the Ergun wikipedia page more recognizable. We perform a similar separation in the implementation of the Ergun Forchheimer coefficient outlined in FunctorErgunDragCoefficients.
Input Parameters
- Darcy_nameName of the Darcy coefficients property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:Name of the Darcy coefficients property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- Forchheimer_nameName of the Forchheimer coefficients property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:Name of the Forchheimer coefficients property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- blockThe list of blocks (ids or names) that this object will be appliedC++ Type:std::vector<SubdomainName> Controllable:No Description:The list of blocks (ids or names) that this object will be applied 
- is_porous_mediumTrueBoolean to choose the type of medium.Default:True C++ Type:bool Controllable:No Description:Boolean to choose the type of medium. 
- 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) 
- muThe dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:The dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- porosity1The porosity. 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 porosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- speedThe magnitude of the interstitial velocity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:The magnitude of the interstitial velocity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- standard_friction_formulationTrueBoolean to choose the type of friction formulation.Default:True C++ Type:bool Controllable:No Description:Boolean to choose the type of friction formulation. 
- uThe velocity in the x direction. Superficial in the case of porous treatment, interstitial otherwise. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:The velocity in the x direction. Superficial in the case of porous treatment, interstitial otherwise. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- vThe velocity in the y direction. Superficial in the case of porous treatment, interstitial otherwise. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:The velocity in the y direction. Superficial in the case of porous treatment, interstitial otherwise. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. 
- wThe velocity in the z direction. Superficial in the case of porous treatment, interstitial otherwise. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.C++ Type:MooseFunctorName Unit:(no unit assumed) Controllable:No Description:The velocity in the z direction. Superficial in the case of porous treatment, interstitial otherwise. 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 contributionC++ Type:std::vector<TagName> Controllable:No Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution 
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fillC++ Type:std::vector<TagName> Controllable:No Description:The extra tags for the matrices this Kernel should fill 
- extra_vector_tagsThe extra tags for the vectors this Kernel should fillC++ Type:std::vector<TagName> Controllable:No Description:The extra tags for the vectors this Kernel should fill 
- matrix_tagssystemThe tag for the matrices this Kernel should fillDefault:system C++ Type:MultiMooseEnum Options:nontime, system Controllable:No Description:The tag for the matrices this Kernel should fill 
- vector_tagsnontimeThe tag for the vectors this Kernel should fillDefault:nontime C++ Type:MultiMooseEnum Options:nontime, time Controllable:No Description:The tag for the vectors this Kernel should fill 
Contribution To Tagged Field Data Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.C++ Type:std::vector<std::string> Controllable:No Description:Adds user-defined labels for accessing object parameters via control logic. 
- enableTrueSet the enabled status of the MooseObject.Default:True C++ Type:bool Controllable:Yes Description:Set the enabled status of the MooseObject. 
- implicitTrueDetermines whether this object is calculated using an implicit or explicit formDefault:True C++ Type:bool Controllable:No Description:Determines whether this object is calculated using an implicit or explicit form 
- seed0The seed for the master random number generatorDefault:0 C++ Type:unsigned int Controllable:No Description:The seed for the master random number generator 
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.Default:False C++ Type:bool Controllable:No Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used. 
Advanced Parameters
- 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
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.C++ Type:MaterialPropertyName Unit:(no unit assumed) Controllable:No Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character. 
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.Default:False C++ Type:bool Controllable:No Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction. 
Material Property Retrieval Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/segregated/2d/2d-segregated-velocity-rz-slip.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-advection-slip.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/natural_convection/natural_circulation_pipe.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux-transient.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/solidification/pipe_solidification.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/hydraulic-separators/separator-no-jump.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/hydraulic-separators/separator-mixing.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/mms/porosity_change/pressure-interpolation-corrected.i)
- (modules/navier_stokes/examples/solidification/gallium_melting.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/discontinuous-body-forces.i)
- (modules/navier_stokes/test/tests/finite_volume/materials/ergun/ergun.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/friction/2d-rc-friction.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc-friction.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/segregated/2d-momentum.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc-rz-by-parts.i)
standard_friction_formulation
Default:True
C++ Type:bool
Controllable:No
Description:Boolean to choose the type of friction formulation.
standard_friction_formulation
Default:True
C++ Type:bool
Controllable:No
Description:Boolean to choose the type of friction formulation.
is_porous_medium
Default:True
C++ Type:bool
Controllable:No
Description:Boolean to choose the type of medium.
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/segregated/2d/2d-segregated-velocity-rz-slip.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
velocity_interp_method = 'rc'
pressure_tag = "pressure_grad"
[Mesh]
  coord_type = 'RZ'
  rz_coord_axis = X
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1.25'
    dy = '0.2'
    ix = '30'
    iy = '7'
  []
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
[]
[Problem]
  nl_sys_names = 'u_system v_system pressure_system'
  previous_nl_solution_required = true
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolatorSegregated
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 0.5
    solver_sys = u_system
    two_term_boundary_expansion = false
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 0.0
    solver_sys = v_system
    two_term_boundary_expansion = false
  []
  [pressure]
    type = INSFVPressureVariable
    solver_sys = pressure_system
    initial_condition = 0.2
    two_term_boundary_expansion = false
  []
[]
[FVKernels]
  inactive = 'u_friction v_friction'
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = ${mu}
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
    extra_vector_tags = ${pressure_tag}
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = vel_x
    u = vel_x
    v = vel_y
    momentum_component = 'x'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    standard_friction_formulation = false
    rho = ${rho}
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = ${mu}
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
    extra_vector_tags = ${pressure_tag}
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = vel_y
    u = vel_x
    v = vel_y
    momentum_component = 'y'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    standard_friction_formulation = false
    rho = ${rho}
  []
  [p_diffusion]
    type = FVAnisotropicDiffusion
    variable = pressure
    coeff = "Ainv"
    coeff_interp_method = 'average'
  []
  [p_source]
    type = FVDivergence
    variable = pressure
    vector_field = "HbyA"
    force_boundary_execution = true
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_x
    function = '1.1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_y
    function = '0.0'
  []
  [walls-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = vel_x
    momentum_component = 'x'
  []
  [walls-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = vel_y
    momentum_component = 'y'
  []
  [symmetry_u]
    type = INSFVSymmetryVelocityBC
    variable = vel_x
    boundary = 'bottom'
    momentum_component = 'x'
    mu = ${mu}
    u = vel_x
    v = vel_y
  []
  [symmetry_v]
    type = INSFVSymmetryVelocityBC
    variable = vel_y
    boundary = 'bottom'
    momentum_component = 'y'
    mu = ${mu}
    u = vel_x
    v = vel_y
  []
  [symmetry_pressure]
    type = INSFVSymmetryPressureBC
    boundary = 'bottom'
    variable = pressure
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 1.4
  []
[]
[FunctorMaterials]
  [darcy]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coefficient Forchheimer_coefficient'
    prop_values = '0.1 0.1 0.1 0.1 0.1 0.1'
  []
[]
[Executioner]
  type = SIMPLENonlinearAssembly
  momentum_l_abs_tol = 1e-14
  pressure_l_abs_tol = 1e-14
  momentum_l_tol = 0
  pressure_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  pressure_gradient_tag = ${pressure_tag}
  momentum_equation_relaxation = 0.5
  pressure_variable_relaxation = 0.3
  num_iterations = 150
  pressure_absolute_tolerance = 1e-13
  momentum_absolute_tolerance = 1e-13
  print_fields = false
  continue_on_max_its = true
[]
[Outputs]
  exodus = true
  csv = false
  perf_graph = false
  print_nonlinear_residuals = false
  print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux.i)
mu = 1.0
rho = 10.0
mu_d = 0.1
rho_d = 1.0
l = 2
U = 1
dp = 0.01
inlet_phase_2 = 0.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
  rhie_chow_user_object = 'rc'
  density_interp_method = 'average'
  mu_interp_method = 'average'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[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
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [phase_2]
    type = INSFVScalarFieldVariable
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
    momentum_component = 'x'
  []
  [u_drift]
    type = WCNSFV2PMomentumDriftFlux
    variable = vel_x
    rho_d = ${rho_d}
    fd = 'phase_2'
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = 'mu_mixture'
    limit_interpolation = true
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction]
    type = PINSFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    is_porous_medium = false
    momentum_component = x
    mu = mu_mixture
    rho = rho_mixture
    variable = vel_x
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
    momentum_component = 'y'
  []
  [v_drift]
    type = WCNSFV2PMomentumDriftFlux
    variable = vel_y
    rho_d = ${rho_d}
    fd = 'phase_2'
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    momentum_component = 'x'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = 'mu_mixture'
    limit_interpolation = true
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction]
    type = PINSFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    is_porous_medium = false
    momentum_component = y
    mu = mu_mixture
    rho = rho_mixture
    variable = vel_y
  []
  [phase_2_advection]
    type = INSFVScalarFieldAdvection
    variable = phase_2
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = 'upwind'
  []
  [phase_2_src]
    type = NSFVMixturePhaseInterface
    variable = phase_2
    phase_coupled = phase_1
    alpha = 0.1
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_x
    functor = '${U}'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_y
    functor = '0'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_x
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_y
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = '0'
  []
  [inlet_phase_2]
    type = FVDirichletBC
    boundary = 'left'
    variable = phase_2
    value = ${inlet_phase_2}
  []
[]
[AuxVariables]
  [drag_coefficient]
    type = MooseVariableFVReal
  []
  [rho_mixture_var]
    type = MooseVariableFVReal
  []
  [mu_mixture_var]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [populate_cd]
    type = FunctorAux
    variable = drag_coefficient
    functor = 'Darcy_coefficient'
  []
  [populate_rho_mixture_var]
    type = FunctorAux
    variable = rho_mixture_var
    functor = 'rho_mixture'
  []
  [populate_mu_mixture_var]
    type = FunctorAux
    variable = mu_mixture_var
    functor = 'mu_mixture'
  []
[]
[FunctorMaterials]
  [populate_u_slip]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    slip_velocity_name = 'vel_slip_x'
    momentum_component = 'x'
    u = 'vel_x'
    v = 'vel_y'
    rho = ${rho}
    mu = 'mu_mixture'
    rho_d = ${rho_d}
    particle_diameter = ${dp}
    linear_coef_name = 'Darcy_coefficient'
    outputs = 'out'
    output_properties = 'vel_slip_x'
    ghost_layers = 5
  []
  [populate_v_slip]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    slip_velocity_name = 'vel_slip_y'
    momentum_component = 'y'
    u = 'vel_x'
    v = 'vel_y'
    rho = ${rho}
    mu = 'mu_mixture'
    rho_d = ${rho_d}
    particle_diameter = ${dp}
    linear_coef_name = 'Darcy_coefficient'
    outputs = 'out'
    output_properties = 'vel_slip_y'
    ghost_layers = 5
  []
  [compute_phase_1]
    type = ADParsedFunctorMaterial
    property_name = phase_1
    functor_names = 'phase_2'
    expression = '1 - phase_2'
    outputs = 'out'
    output_properties = 'phase_1'
  []
  [CD]
    type = NSFVDispersePhaseDragFunctorMaterial
    rho = 'rho_mixture'
    mu = mu_mixture
    u = 'vel_x'
    v = 'vel_y'
    particle_diameter = ${dp}
  []
  [mixing_material]
    type = NSFVMixtureFunctorMaterial
    phase_2_names = '${rho} ${mu}'
    phase_1_names = '${rho_d} ${mu_d}'
    prop_names = 'rho_mixture mu_mixture'
    phase_1_fraction = 'phase_2'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  nl_rel_tol = 1e-10
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_shift_type'
    petsc_options_value = 'lu       NONZERO'
  []
[]
[Outputs]
  print_linear_residuals = true
  print_nonlinear_residuals = true
  [out]
    type = Exodus
    hide = 'Re lin cum_lin'
  []
  [perf]
    type = PerfGraphOutput
  []
[]
[Postprocessors]
  [Re]
    type = ParsedPostprocessor
    expression = '${rho} * ${l} * ${U}'
  []
  [lin]
    type = NumLinearIterations
  []
  [cum_lin]
    type = CumulativeValuePostprocessor
    postprocessor = lin
  []
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-advection-slip.i)
mu = 1.0
rho = 10.0
mu_d = 0.1
rho_d = 1.0
l = 2
U = 1
dp = 0.01
inlet_phase_2 = 0.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
  rhie_chow_user_object = 'rc'
  mu_interp_method = 'average'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = '${fparse l * 5}'
    ymin = '${fparse -l / 2}'
    ymax = '${fparse l / 2}'
    nx = 10
    ny = 6
  []
  uniform_refine = 0
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [phase_2]
    type = INSFVScalarFieldVariable
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
    momentum_component = 'x'
  []
  [u_advection_slip]
    type = WCNSFV2PMomentumAdvectionSlip
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    rho_d = ${rho_d}
    fd = phase_2
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = 'mu_mixture'
    limit_interpolation = true
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction]
    type = PINSFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    is_porous_medium = false
    momentum_component = x
    mu = mu_mixture
    rho = rho_mixture
    variable = vel_x
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
    momentum_component = 'y'
  []
  [v_advection_slip]
    type = WCNSFV2PMomentumAdvectionSlip
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    rho_d = ${rho_d}
    fd = phase_2
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = 'mu_mixture'
    limit_interpolation = true
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction]
    type = PINSFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    is_porous_medium = false
    momentum_component = y
    mu = mu_mixture
    rho = rho_mixture
    variable = vel_y
  []
  [phase_2_advection]
    type = INSFVScalarFieldAdvection
    variable = phase_2
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = 'upwind'
  []
  [phase_2_src]
    type = NSFVMixturePhaseInterface
    variable = phase_2
    phase_coupled = phase_1
    alpha = 0.1
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_x
    functor = '${U}'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_y
    functor = '0'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_x
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_y
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = '0'
  []
  [inlet_phase_2]
    type = FVDirichletBC
    boundary = 'left'
    variable = phase_2
    value = ${inlet_phase_2}
  []
[]
[AuxVariables]
  [drag_coefficient]
    type = MooseVariableFVReal
  []
  [rho_mixture_var]
    type = MooseVariableFVReal
  []
  [mu_mixture_var]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [populate_cd]
    type = FunctorAux
    variable = drag_coefficient
    functor = 'Darcy_coefficient'
  []
  [populate_rho_mixture_var]
    type = FunctorAux
    variable = rho_mixture_var
    functor = 'rho_mixture'
  []
  [populate_mu_mixture_var]
    type = FunctorAux
    variable = mu_mixture_var
    functor = 'mu_mixture'
  []
[]
[FunctorMaterials]
  [phase_1]
    property_name = 'phase_1'
    type = ADParsedFunctorMaterial
    functor_names = 'phase_2'
    expression = '1 - phase_2'
    outputs = 'out'
    output_properties = 'phase_1'
  []
  [populate_u_slip]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    slip_velocity_name = 'vel_slip_x'
    momentum_component = 'x'
    u = 'vel_x'
    v = 'vel_y'
    rho = ${rho}
    mu = 'mu_mixture'
    rho_d = ${rho_d}
    particle_diameter = ${dp}
    linear_coef_name = 'Darcy_coefficient'
    outputs = 'out'
    output_properties = 'vel_slip_x'
  []
  [populate_v_slip]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    slip_velocity_name = 'vel_slip_y'
    momentum_component = 'y'
    u = 'vel_x'
    v = 'vel_y'
    rho = ${rho}
    mu = 'mu_mixture'
    rho_d = ${rho_d}
    particle_diameter = ${dp}
    linear_coef_name = 'Darcy_coefficient'
    outputs = 'out'
    output_properties = 'vel_slip_y'
  []
  [CD]
    type = NSFVDispersePhaseDragFunctorMaterial
    rho = 'rho_mixture'
    mu = mu_mixture
    u = 'vel_x'
    v = 'vel_y'
    particle_diameter = ${dp}
  []
  [mixing_material]
    type = NSFVMixtureFunctorMaterial
    phase_2_names = '${rho} ${mu}'
    phase_1_names = '${rho_d} ${mu_d}'
    prop_names = 'rho_mixture mu_mixture'
    phase_1_fraction = 'phase_2'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  nl_rel_tol = 1e-10
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_shift_type'
    petsc_options_value = 'lu       NONZERO'
  []
[]
[Outputs]
  [out]
    type = Exodus
    hide = 'Re lin cum_lin'
  []
[]
[Postprocessors]
  [Re]
    type = ParsedPostprocessor
    expression = '${rho} * ${l} * ${U}'
  []
  [lin]
    type = NumLinearIterations
  []
  [cum_lin]
    type = CumulativeValuePostprocessor
    postprocessor = lin
  []
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/natural_convection/natural_circulation_pipe.i)
# natural convection through a pipe
# Reference solution in "reference_pipe_natural_convection.py"
# Reference mdot: 0.0792 kg/s
# this input
# iy   mdot
# 10   8.302364e-02
# 20   8.111192e-02
# 40   8.007924e-02
# 80   7.954403e-02
# 160  7.927201e-02
# Convergence to the analytical result is observed
height = 10.0
gravity = 9.81
p0 = 1e5
molar_mass = 29.0e-3
T0 = 328
Ru = 8.3145
Ri = '${fparse Ru / molar_mass}'
density = '${fparse p0 / (Ri * T0)}'
head = '${fparse height * density * gravity}'
k = 25.68e-3
gamma = 1.4
[Mesh]
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '0.1'
    ix = '2'
    dy = '${height}'
    iy = '5'
  []
[]
[GlobalParams]
  rhie_chow_user_object = pins_rhie_chow_interpolator
[]
[FluidProperties]
  [air]
    type = IdealGasFluidProperties
    molar_mass = ${molar_mass}
    k = ${k}
    gamma = ${gamma}
  []
[]
[Modules]
  [NavierStokesFV]
    compressibility = 'weakly-compressible'
    add_energy_equation = true
    gravity = '0 -${gravity} 0'
    density = rho
    dynamic_viscosity = mu
    specific_heat = cp
    thermal_conductivity = k
    initial_velocity = '0 1e-6 0'
    initial_pressure = ${p0}
    initial_temperature = ${T0}
    inlet_boundaries = 'bottom'
    momentum_inlet_types = 'fixed-pressure'
    momentum_inlet_function = '${fparse p0 + head}'
    energy_inlet_types = 'fixed-temperature'
    energy_inlet_function = '${T0}'
    energy_scaling = 1e-5
    wall_boundaries = 'left right'
    momentum_wall_types = 'slip slip'
    energy_wall_types = 'heatflux heatflux'
    energy_wall_function = '300 300'
    outlet_boundaries = 'top'
    momentum_outlet_types = 'fixed-pressure'
    pressure_function = '${fparse p0}'
    momentum_advection_interpolation = 'upwind'
    mass_advection_interpolation = 'upwind'
    porous_medium_treatment = true
    porosity = porosity
    energy_advection_interpolation = 'average'
  []
[]
[FVKernels]
  [u_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_x
    Darcy_name = linear_friction_coeff
    momentum_component = 'x'
    standard_friction_formulation = false
    rho = rho
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_y
    Darcy_name = linear_friction_coeff
    momentum_component = 'y'
    standard_friction_formulation = false
    rho = rho
  []
[]
[Executioner]
  type = Transient
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'lu        NONZERO'
  nl_rel_tol = 1e-8
  nl_abs_tol = 1e-6
  end_time = 1e4
  [TimeStepper]
    type = IterationAdaptiveDT
    dt = 0.1
    growth_factor = 2
    iteration_window = 2
    optimal_iterations = 6
  []
[]
[Functions]
  [mu_rampdown_fn]
    type = PiecewiseLinear
    x = '0    0.5  1   5  10 100 1000 2000'
    y = '1000 1000 100 10 1  1   1    0'
  []
[]
[FunctorMaterials]
  [fluid_props_to_mat_props]
    type = GeneralFunctorFluidProps
    fp = air
    pressure = pressure
    T_fluid = T_fluid
    speed = speed
    force_define_density = true
    neglect_derivatives_of_density_time_derivative = false
    mu_rampdown = 'mu_rampdown_fn'
    characteristic_length = 1
    porosity = porosity
  []
  [scalar_props]
    type = ADGenericFunctorMaterial
    prop_names = 'porosity loss_coeff'
    prop_values = '1       1.3'
  []
  [linear_friction]
    type = ADParsedFunctorMaterial
    property_name = 'linear_friction'
    expression = 'loss_coeff * rho'
    functor_names = 'loss_coeff rho'
  []
  [linear_friction_coeff]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'linear_friction_coeff'
    prop_values = 'linear_friction linear_friction linear_friction'
  []
[]
[AuxVariables]
  [rho_var]
    type = MooseVariableFVReal
  []
  [cp_var]
    type = MooseVariableFVReal
  []
  [rho_cp_T_fluid_var]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [rho_var_aux]
    type = FunctorAux
    variable = rho_var
    functor = rho
  []
  [cp_var_aux]
    type = FunctorAux
    variable = cp_var
    functor = cp
  []
  [rho_cp_T_fluid_var_aux]
    type = ParsedAux
    variable = rho_cp_T_fluid_var
    coupled_variables = 'rho_var cp_var T_fluid'
    expression = 'rho_var * cp_var * T_fluid'
  []
[]
[Postprocessors]
  [inlet_mfr]
    type = VolumetricFlowRate
    vel_x = superficial_vel_x
    vel_y = superficial_vel_y
    advected_quantity = rho
    boundary = bottom
    advected_interp_method = average
  []
  [outlet_mfr]
    type = VolumetricFlowRate
    vel_x = superficial_vel_x
    vel_y = superficial_vel_y
    advected_quantity = rho
    boundary = top
    advected_interp_method = average
  []
  [inlet_energy]
    type = VolumetricFlowRate
    vel_x = superficial_vel_x
    vel_y = superficial_vel_y
    advected_quantity = rho_cp_T_fluid_var
    boundary = bottom
    advected_interp_method = average
  []
  [outlet_energy]
    type = VolumetricFlowRate
    vel_x = superficial_vel_x
    vel_y = superficial_vel_y
    advected_quantity = rho_cp_T_fluid_var
    boundary = top
    advected_interp_method = average
  []
[]
[Debug]
  show_var_residual_norms = true
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux-transient.i)
mu = 1.0
rho = 10.0
mu_d = 0.1
rho_d = 1.0
l = 2
U = 1
dp = 0.01
inlet_phase_2 = 0.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
  rhie_chow_user_object = 'rc'
  density_interp_method = 'average'
  mu_interp_method = 'average'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[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
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 0
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [phase_2]
    type = INSFVScalarFieldVariable
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
  []
  [u_time]
    type = INSFVMomentumTimeDerivative
    variable = vel_x
    rho = 'rho_mixture'
    momentum_component = 'x'
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
    momentum_component = 'x'
  []
  [u_drift]
    type = WCNSFV2PMomentumDriftFlux
    variable = vel_x
    rho_d = ${rho_d}
    fd = 'phase_2'
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = 'mu_mixture'
    limit_interpolation = true
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction]
    type = PINSFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    is_porous_medium = false
    momentum_component = x
    mu = mu_mixture
    rho = rho_mixture
    variable = vel_x
  []
  [v_time]
    type = INSFVMomentumTimeDerivative
    variable = vel_y
    rho = 'rho_mixture'
    momentum_component = 'y'
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = 'rho_mixture'
    momentum_component = 'y'
  []
  [v_drift]
    type = WCNSFV2PMomentumDriftFlux
    variable = vel_y
    rho_d = ${rho_d}
    fd = 'phase_2'
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = 'mu_mixture'
    limit_interpolation = true
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction]
    type = PINSFVMomentumFriction
    Darcy_name = Darcy_coefficient_vec
    is_porous_medium = false
    momentum_component = y
    mu = mu_mixture
    rho = rho_mixture
    variable = vel_y
  []
  [phase_2_time]
    type = FVFunctorTimeKernel
    variable = phase_2
    functor = phase_2
  []
  [phase_2_advection]
    type = INSFVScalarFieldAdvection
    variable = phase_2
    u_slip = 'vel_slip_x'
    v_slip = 'vel_slip_y'
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = 'upwind'
  []
  [phase_2_src]
    type = NSFVMixturePhaseInterface
    variable = phase_2
    phase_coupled = phase_1
    alpha = 0.1
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_x
    functor = '${U}'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_y
    functor = '0'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_x
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_y
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = '0'
  []
  [inlet_phase_2]
    type = FVDirichletBC
    boundary = 'left'
    variable = phase_2
    value = ${inlet_phase_2}
  []
[]
[AuxVariables]
  [drag_coefficient]
    type = MooseVariableFVReal
  []
  [rho_mixture_var]
    type = MooseVariableFVReal
  []
  [mu_mixture_var]
    type = MooseVariableFVReal
  []
  [vel_slip_x_var]
    type = MooseVariableFVReal
  []
  [vel_slip_y_var]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [populate_cd]
    type = FunctorAux
    variable = drag_coefficient
    functor = 'Darcy_coefficient'
  []
  [populate_rho_mixture_var]
    type = FunctorAux
    variable = rho_mixture_var
    functor = 'rho_mixture'
  []
  [populate_mu_mixture_var]
    type = FunctorAux
    variable = mu_mixture_var
    functor = 'mu_mixture'
  []
  [populate_vx_slip_var]
    type = FunctorAux
    variable = vel_slip_x_var
    functor = 'vel_slip_x'
  []
  [populate_vy_slip_var]
    type = FunctorAux
    variable = vel_slip_y_var
    functor = 'vel_slip_y'
  []
[]
[FunctorMaterials]
  [populate_u_slip]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    slip_velocity_name = 'vel_slip_x'
    momentum_component = 'x'
    u = 'vel_x'
    v = 'vel_y'
    rho = ${rho}
    mu = 'mu_mixture'
    rho_d = ${rho_d}
    particle_diameter = ${dp}
    linear_coef_name = 'Darcy_coefficient'
  []
  [populate_v_slip]
    type = WCNSFV2PSlipVelocityFunctorMaterial
    slip_velocity_name = 'vel_slip_y'
    momentum_component = 'y'
    u = 'vel_x'
    v = 'vel_y'
    rho = ${rho}
    mu = 'mu_mixture'
    rho_d = ${rho_d}
    particle_diameter = ${dp}
    linear_coef_name = 'Darcy_coefficient'
  []
  [compute_phase_1]
    type = ADParsedFunctorMaterial
    property_name = phase_1
    functor_names = 'phase_2'
    expression = '1 - phase_2'
  []
  [CD]
    type = NSFVDispersePhaseDragFunctorMaterial
    rho = 'rho_mixture'
    mu = mu_mixture
    u = 'vel_x'
    v = 'vel_y'
    particle_diameter = ${dp}
  []
  [mixing_material]
    type = NSFVMixtureFunctorMaterial
    phase_2_names = '${rho} ${mu}'
    phase_1_names = '${rho_d} ${mu_d}'
    prop_names = 'rho_mixture mu_mixture'
    phase_1_fraction = 'phase_2'
  []
[]
[Executioner]
  type = Transient
  solve_type = 'NEWTON'
  nl_rel_tol = 1e-10
  dt = 0.1
  end_time = 1.0
[]
[Preconditioning]
  [SMP]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_shift_type'
    petsc_options_value = 'lu       NONZERO'
  []
[]
[Outputs]
  exodus = false
  [CSV]
    type = CSV
    execute_on = 'TIMESTEP_END'
  []
[]
[Postprocessors]
  [Re]
    type = ParsedPostprocessor
    expression = '${rho} * ${l} * ${U}'
  []
  [rho_outlet]
    type = SideAverageValue
    boundary = 'right'
    variable = 'rho_mixture_var'
  []
  [vslip_x]
    type = SideExtremeValue
    boundary = 'left'
    variable = 'vel_slip_x_var'
  []
  [vslip_y]
    type = SideExtremeValue
    boundary = 'left'
    variable = 'vel_slip_y_var'
  []
  [vslip_value]
    type = ParsedPostprocessor
    expression = 'sqrt(vslip_x*vslip_x + vslip_y*vslip_y)*vslip_x/abs(vslip_x)'
    pp_names = 'vslip_x vslip_y'
  []
[]
(modules/navier_stokes/test/tests/finite_volume/ins/solidification/pipe_solidification.i)
mu = 8.8871e-4
rho_solid = 997.561
rho_liquid = 997.561
k_solid = 0.6203
k_liquid = 0.6203
cp_solid = 4181.72
cp_liquid = 4181.72
L = 3e5
T_liquidus = 285
T_solidus = 280
advected_interp_method = 'average'
velocity_interp_method = 'rc'
U_inlet = '${fparse 0.5 * mu / rho_liquid / 0.5}'
T_inlet = 300.0
T_cold = 200.0
Nx = 30
Ny = 5
[GlobalParams]
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[Mesh]
  coord_type = 'RZ'
  rz_coord_axis = 'X'
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 10
    ymin = 0
    ymax = '${fparse 0.5 * 1.0}'
    nx = ${Nx}
    ny = ${Ny}
    bias_y = '${fparse 1 / 1.2}'
  []
  [rename1]
    type = RenameBoundaryGenerator
    input = gen
    old_boundary = 'left'
    new_boundary = 'inlet'
  []
  [rename2]
    type = RenameBoundaryGenerator
    input = rename1
    old_boundary = 'right'
    new_boundary = 'outlet'
  []
  [rename3]
    type = RenameBoundaryGenerator
    input = rename2
    old_boundary = 'bottom'
    new_boundary = 'symmetry'
  []
  [rename4]
    type = RenameBoundaryGenerator
    input = rename3
    old_boundary = 'top'
    new_boundary = 'wall'
  []
  [rename5]
    type = ParsedGenerateSideset
    input = rename4
    normal = '0 1 0'
    combinatorial_geometry = 'x>2.0 & x<8.0 & y>0.49999'
    new_sideset_name = 'cooled_wall'
  []
[]
[AuxVariables]
  [U]
    type = MooseVariableFVReal
  []
  [fl]
    type = MooseVariableFVReal
    initial_condition = 1.0
  []
  [density]
    type = MooseVariableFVReal
  []
  [th_cond]
    type = MooseVariableFVReal
  []
  [cp_var]
    type = MooseVariableFVReal
  []
  [darcy_coef]
    type = MooseVariableFVReal
  []
  [fch_coef]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [mag]
    type = VectorMagnitudeAux
    variable = U
    x = vel_x
    y = vel_y
  []
  [compute_fl]
    type = NSLiquidFractionAux
    variable = fl
    temperature = T
    T_liquidus = '${T_liquidus}'
    T_solidus = '${T_solidus}'
    execute_on = 'TIMESTEP_END'
  []
  [rho_out]
    type = FunctorAux
    functor = 'rho_mixture'
    variable = 'density'
  []
  [th_cond_out]
    type = FunctorAux
    functor = 'k_mixture'
    variable = 'th_cond'
  []
  [cp_out]
    type = FunctorAux
    functor = 'cp_mixture'
    variable = 'cp_var'
  []
  [darcy_out]
    type = FunctorAux
    functor = 'Darcy_coefficient'
    variable = 'darcy_coef'
  []
  [fch_out]
    type = FunctorAux
    functor = 'Forchheimer_coefficient'
    variable = 'fch_coef'
  []
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 0.0
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 0.0
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [T]
    type = INSFVEnergyVariable
    initial_condition = '${T_inlet}'
    scaling = 1.0
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = rho_mixture
  []
  [u_time]
    type = INSFVMomentumTimeDerivative
    variable = vel_x
    rho = rho_mixture
    momentum_component = 'x'
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = rho_mixture
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = ${mu}
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = vel_x
    momentum_component = 'x'
    u = vel_x
    v = vel_y
    Darcy_name = 'Darcy_coeff'
    Forchheimer_name = 'Forchheimer_coeff'
    rho = ${rho_liquid}
    mu = ${mu}
    standard_friction_formulation = false
  []
  [v_time]
    type = INSFVMomentumTimeDerivative
    variable = vel_y
    rho = rho_mixture
    momentum_component = 'y'
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = rho_mixture
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = ${mu}
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = vel_y
    momentum_component = 'y'
    u = vel_x
    v = vel_y
    Darcy_name = 'Darcy_coeff'
    Forchheimer_name = 'Forchheimer_coeff'
    rho = ${rho_liquid}
    mu = ${mu}
    standard_friction_formulation = false
  []
  [T_time]
    type = INSFVEnergyTimeDerivative
    variable = T
    rho = rho_mixture
    dh_dt = dh_dt
  []
  [energy_advection]
    type = INSFVEnergyAdvection
    variable = T
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
  []
  [energy_diffusion]
    type = FVDiffusion
    coeff = k_mixture
    variable = T
  []
  [energy_source]
    type = NSFVPhaseChangeSource
    variable = T
    L = ${L}
    liquid_fraction = fl
    T_liquidus = ${T_liquidus}
    T_solidus = ${T_solidus}
    rho = 'rho_mixture'
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'inlet'
    variable = vel_x
    function = '${U_inlet}'
  []
  [sym_u]
    type = INSFVSymmetryVelocityBC
    boundary = 'symmetry'
    variable = vel_x
    u = vel_x
    v = vel_y
    mu = ${mu}
    momentum_component = 'x'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'inlet'
    variable = vel_y
    function = 0
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'wall'
    variable = vel_x
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'wall'
    variable = vel_y
    function = 0
  []
  [sym_v]
    type = INSFVSymmetryVelocityBC
    boundary = 'symmetry'
    variable = vel_y
    u = vel_x
    v = vel_y
    mu = ${mu}
    momentum_component = y
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'outlet'
    variable = pressure
    function = 0
  []
  [sym_p]
    type = INSFVSymmetryPressureBC
    boundary = 'symmetry'
    variable = pressure
  []
  [sym_T]
    type = INSFVSymmetryScalarBC
    variable = T
    boundary = 'symmetry'
  []
  [cooled_wall]
    type = FVFunctorDirichletBC
    variable = T
    functor = '${T_cold}'
    boundary = 'cooled_wall'
  []
[]
[FunctorMaterials]
  [ins_fv]
    type = INSFVEnthalpyFunctorMaterial
    rho = rho_mixture
    cp = cp_mixture
    temperature = 'T'
  []
  [eff_cp]
    type = NSFVMixtureFunctorMaterial
    phase_2_names = '${cp_solid} ${k_solid} ${rho_solid}'
    phase_1_names = '${cp_liquid} ${k_liquid} ${rho_liquid}'
    prop_names = 'cp_mixture k_mixture rho_mixture'
    phase_1_fraction = fl
  []
  [mushy_zone_resistance]
    type = INSFVMushyPorousFrictionFunctorMaterial
    liquid_fraction = 'fl'
    mu = '${mu}'
    rho_l = '${rho_liquid}'
  []
  [friction]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coeff Forchheimer_coeff'
    prop_values = 'darcy_coef darcy_coef darcy_coef fch_coef fch_coef fch_coef'
  []
[]
[Executioner]
  type = Transient
  dt = 5e3
  end_time = 1e4
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  nl_abs_tol = 1e-8
  nl_max_its = 12
[]
[Postprocessors]
  [average_T]
    type = ElementAverageValue
    variable = T
    outputs = csv
    execute_on = FINAL
  []
[]
[VectorPostprocessors]
  [sat]
    type = LineValueSampler
    warn_discontinuous_face_values = false
    start_point = '0.0 0 0'
    end_point = '10.0 0 0'
    num_points = '${Nx}'
    sort_by = x
    variable = 'T'
    execute_on = FINAL
  []
[]
[Outputs]
  exodus = true
  [csv]
    type = CSV
    execute_on = 'FINAL'
  []
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/hydraulic-separators/separator-no-jump.i)
# This test describes a test where three parallel channels are
# separated using flow separators that act as slip boundary conditions.
# The different channels have different friction factors
# meaning that we expect different pressure drops.
# Channel 1 expected drop (analytic, Forchheimer only): 5.50E-03 Pa
# Channel 2 expected drop (analytic, Forchheimer only): 4.40E-02 Pa
# Channel 3 expected drop (analytic, Forchheimer only): 1.49E-01 Pa
rho=1.1
mu=1.1
advected_interp_method='average'
velocity_interp_method='rc'
[Mesh]
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '1'
    dy = '0.25 0.25 0.25'
    ix = '5'
    iy = '2 2 2'
    subdomain_id = '1 2 3'
  []
  [separator-1]
    type = SideSetsBetweenSubdomainsGenerator
    new_boundary = 'separator-1'
    primary_block = 1
    paired_block = 2
    input = mesh
  []
  [separator-2]
    type = SideSetsBetweenSubdomainsGenerator
    new_boundary = 'separator-2'
    primary_block = 2
    paired_block = 3
    input = separator-1
  []
  [inlet-1]
    type = ParsedGenerateSideset
    input = separator-2
    combinatorial_geometry = 'y < 0.25 & x < 0.00001'
    replace = true
    new_sideset_name = inlet-1
  []
  [inlet-2]
    type = ParsedGenerateSideset
    input = inlet-1
    combinatorial_geometry = 'y > 0.25 & y < 0.5 & x < 0.00001'
    replace = true
    new_sideset_name = inlet-2
  []
  [inlet-3]
    type = ParsedGenerateSideset
    input = inlet-2
    combinatorial_geometry = 'y > 0.5 & x < 0.00001'
    replace = true
    new_sideset_name = inlet-3
  []
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
  porosity = porosity
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolator
    u = superficial_vel_x
    v = superficial_vel_y
    pressure = pressure
  []
[]
[Variables]
  [superficial_vel_x]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 0.1
  []
  [superficial_vel_y]
    type = PINSFVSuperficialVelocityVariable
  []
  [pressure]
    type = BernoulliPressureVariable
    u = u
    v = v
    rho = ${rho}
  []
[]
[FVKernels]
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [u_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_x
    momentum_component = 'x'
    mu = ${mu}
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_x
    momentum_component = 'x'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
  []
  [v_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_y
    momentum_component = 'y'
    mu = ${mu}
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_y
    pressure = pressure
    momentum_component = 'y'
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_y
    momentum_component = 'y'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
  []
[]
[FVBCs]
  [inlet-u-1]
    type = INSFVInletVelocityBC
    boundary = 'inlet-1'
    variable = superficial_vel_x
    function = '0.1'
  []
  [inlet-u-2]
    type = INSFVInletVelocityBC
    boundary = 'inlet-2'
    variable = superficial_vel_x
    function = '0.2'
  []
  [inlet-u-3]
    type = INSFVInletVelocityBC
    boundary = 'inlet-3'
    variable = superficial_vel_x
    function = '0.3'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'inlet-1 inlet-2 inlet-3'
    variable = superficial_vel_y
    function = 0
  []
  [walls-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top bottom'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [walls-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top bottom'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  [separator-u]
    type = INSFVVelocityHydraulicSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [separator-v]
    type = INSFVVelocityHydraulicSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  [separator-p]
    type = INSFVScalarFieldSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = pressure
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0.4
  []
[]
[FunctorMaterials]
  [const]
    type = ADGenericFunctorMaterial
    prop_names = 'porosity'
    prop_values = '1.0'
  []
  [darcy-1]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '1.0 1.0 1.0'
    block = 1
  []
  [darcy-2]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '2.0 2.0 2.0'
    block = 2
  []
  [darcy-3]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '3.0 3.0 3.0'
    block = 3
  []
  [speed]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = superficial_vel_x
    superficial_vel_y = superficial_vel_y
    porosity = porosity
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = ' lu       NONZERO               1e-10'
  line_search = 'none'
  nl_rel_tol = 1e-10
  nl_max_its = 10
[]
[Postprocessors]
  [inlet_p1]
    type = SideAverageValue
    variable = 'pressure'
    boundary = 'inlet-1'
  []
  [inlet_p2]
    type = SideAverageValue
    variable = 'pressure'
    boundary = 'inlet-2'
  []
  [inlet_p3]
    type = SideAverageValue
    variable = 'pressure'
    boundary = 'inlet-3'
  []
  [drop-1]
    type = ParsedPostprocessor
    expression = 'inlet_p1 - outlet'
    pp_names = 'inlet_p1'
    constant_names = 'outlet'
    constant_expressions = '0.4'
  []
  [drop-2]
    type = ParsedPostprocessor
    expression = 'inlet_p2 - outlet'
    pp_names = 'inlet_p2'
    constant_names = 'outlet'
    constant_expressions = '0.4'
  []
  [drop-3]
    type = ParsedPostprocessor
    expression = 'inlet_p3 - outlet'
    pp_names = 'inlet_p3'
    constant_names = 'outlet'
    constant_expressions = '0.4'
  []
[]
[Outputs]
  csv = true
  execute_on = final
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/hydraulic-separators/separator-mixing.i)
# This test is designed to check for energy conservation
# in separated channels. The three inlet temperatures should be
# preserved at the outlets.
rho=1.1
mu=1e-4
k=2.1
cp=5.5
advected_interp_method='upwind'
velocity_interp_method='rc'
[Mesh]
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '0.25 1.0 0.25'
    dy = '0.25 0.25 0.25'
    ix = '4 20 4'
    iy = '5 5 5'
    subdomain_id = '1 2 5 1 3 5 1 4 5'
  []
  [separator-1]
    type = SideSetsBetweenSubdomainsGenerator
    input = mesh
    primary_block = '2'
    paired_block = '3'
    new_boundary = 'separator-1'
  []
  [separator-2]
    type = SideSetsBetweenSubdomainsGenerator
    input = separator-1
    primary_block = '3'
    paired_block = '4'
    new_boundary = 'separator-2'
  []
  [jump-1]
    type = SideSetsBetweenSubdomainsGenerator
    input = separator-2
    primary_block = '1'
    paired_block = '2'
    new_boundary = jump-1
  []
  [jump-2]
    type = SideSetsBetweenSubdomainsGenerator
    input = jump-1
    primary_block = '1'
    paired_block = '3'
    new_boundary = jump-2
  []
  [jump-3]
    type = SideSetsBetweenSubdomainsGenerator
    input = jump-2
    primary_block = '1'
    paired_block = '4'
    new_boundary = jump-3
  []
  [outlet-1]
    type = SideSetsBetweenSubdomainsGenerator
    input = jump-3
    primary_block = '2'
    paired_block = '5'
    new_boundary = outlet-1
  []
  [outlet-2]
    type = SideSetsBetweenSubdomainsGenerator
    input = outlet-1
    primary_block = '3'
    paired_block = '5'
    new_boundary = outlet-2
  []
  [outlet-3]
    type = SideSetsBetweenSubdomainsGenerator
    input = outlet-2
    primary_block = '4'
    paired_block = '5'
    new_boundary = outlet-3
  []
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
  porosity = porosity
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolator
    u = superficial_vel_x
    v = superficial_vel_y
    pressure = pressure
  []
[]
[Variables]
  [superficial_vel_x]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 0.1
  []
  [superficial_vel_y]
    type = PINSFVSuperficialVelocityVariable
  []
  [pressure]
    type = BernoulliPressureVariable
    u = superficial_vel_x
    v = superficial_vel_y
    rho = ${rho}
    pressure_drop_sidesets = 'jump-1 jump-2 jump-3 outlet-1 outlet-2 outlet-3'
    pressure_drop_form_factors = '0.1 0.2 0.3 0.1 0.2 0.3'
  []
  [T_fluid]
    type = INSFVEnergyVariable
    initial_condition = 300
  []
[]
[FVKernels]
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [u_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_x
    momentum_component = 'x'
    mu = ${mu}
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_x
    pressure = pressure
    momentum_component = 'x'
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_x
    momentum_component = 'x'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
  []
  [v_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_y
    momentum_component = 'y'
    mu = ${mu}
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_y
    pressure = pressure
    momentum_component = 'y'
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_y
    momentum_component = 'y'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
  []
  [temp_conduction]
    type = FVDiffusion
    coeff = ${k}
    variable = T_fluid
  []
  [temp_advection]
    type = INSFVEnergyAdvection
    variable = T_fluid
  []
  [temp_source]
    type = FVBodyForce
    variable = T_fluid
    function = heating
    block = '2 3 4'
  []
[]
[Functions]
  [heating]
    type = ParsedFunction
    expression = 'if(y<0.25, 10, if(y<0.5, 20, 30))'
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = superficial_vel_x
    function = '0.1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = superficial_vel_y
    function = 0
  []
  [inlet-T]
    type = FVDirichletBC
    variable = T_fluid
    boundary = 'left'
    value = 300
  []
  [walls-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top bottom'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [walls-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top bottom'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  [separator-u]
    type = INSFVVelocityHydraulicSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [separator-v]
    type = INSFVVelocityHydraulicSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  [separator-p]
    type = INSFVScalarFieldSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = pressure
  []
  [separator-T]
    type = INSFVScalarFieldSeparatorBC
    boundary = 'separator-1 separator-2'
    variable = T_fluid
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0.4
  []
[]
[FunctorMaterials]
  [porosity]
    type = ADPiecewiseByBlockFunctorMaterial
    prop_name = porosity
    subdomain_to_prop_value = '1 0.8
                               2 0.7
                               3 0.6
                               4 0.5
                               5 0.8'
  []
  [darcy-1]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '1.0 1.0 1.0'
    block = '1 5'
  []
  [darcy-2]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '3.0 3.0 3.0'
    block = 2
  []
  [darcy-3]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '1.5 1.5 1.5'
    block = 3
  []
  [darcy-4]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Forchheimer_coefficient'
    prop_values = '0.75 0.75 0.75'
    block = 4
  []
  [speed]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = superficial_vel_x
    superficial_vel_y = superficial_vel_y
    porosity = porosity
  []
  [ins_fv]
    type = INSFVEnthalpyFunctorMaterial
    temperature = 'T_fluid'
    rho = ${rho}
    cp = ${cp}
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = ' lu       NONZERO               1e-10'
  line_search = 'none'
  nl_rel_tol = 1e-10
[]
[Postprocessors]
  [outlet_T1]
    type = SideAverageValue
    variable = 'T_fluid'
    boundary = 'right'
  []
[]
[Outputs]
  csv = true
  execute_on = final
[]
(modules/navier_stokes/test/tests/finite_volume/pins/mms/porosity_change/pressure-interpolation-corrected.i)
mu = 1.1
rho = 1.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
darcy = 1.1
forch = 1.1
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 2
    ymin = -1
    ymax = 1
    nx = 2
    ny = 2
  []
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
  Darcy_name = 'Darcy_coefficient'
  Forchheimer_name = 'Forchheimer_coefficient'
  porosity = porosity
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolator
    u = u
    v = v
    porosity = porosity
    pressure = pressure
    smoothing_layers = 2
  []
[]
[Variables]
  [u]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [v]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
[]
[AuxVariables]
  [eps_out]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [eps_out]
    type = FunctorAux
    variable = eps_out
    functor = porosity
    execute_on = 'timestep_end'
  []
[]
[FVKernels]
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [mass_forcing]
    type = FVBodyForce
    variable = pressure
    function = forcing_p
  []
  [u_advection]
    type = PINSFVMomentumAdvection
    variable = u
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = u
    mu = ${mu}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = u
    pressure = pressure
    porosity = porosity
    momentum_component = 'x'
  []
  [u_drag]
    type = PINSFVMomentumFriction
    variable = u
    momentum_component = 'x'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [u_correction]
    type = PINSFVMomentumFrictionCorrection
    variable = u
    momentum_component = 'x'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [u_forcing]
    type = INSFVBodyForce
    variable = u
    functor = forcing_u
    momentum_component = 'x'
  []
  [v_advection]
    type = PINSFVMomentumAdvection
    variable = v
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = v
    mu = ${mu}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = v
    pressure = pressure
    porosity = porosity
    momentum_component = 'y'
  []
  [v_drag]
    type = PINSFVMomentumFriction
    variable = v
    momentum_component = 'y'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [v_correction]
    type = PINSFVMomentumFrictionCorrection
    variable = v
    momentum_component = 'y'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [v_forcing]
    type = INSFVBodyForce
    variable = v
    functor = forcing_v
    momentum_component = 'y'
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    functor = 'exact_u'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    functor = 'exact_v'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = u
    function = 'exact_u'
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = v
    function = 'exact_v'
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 'exact_p'
  []
[]
[FunctorMaterials]
  [darcy]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coefficient Forchheimer_coefficient'
    prop_values = '${darcy} ${darcy} ${darcy} ${forch} ${forch} ${forch}'
  []
  [speed]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = u
    superficial_vel_y = v
    porosity = porosity
  []
[]
[Functions]
  [porosity]
    type = ParsedFunction
    expression = '.5 + .1 * sin(pi * x / 4) * cos(pi * y / 4)'
  []
  [exact_u]
    type = ParsedFunction
    expression = 'sin((1/2)*y*pi)*cos((1/2)*x*pi)'
  []
  [forcing_u]
    type = ParsedFunction
    expression = 'darcy*mu*sin((1/2)*y*pi)*cos((1/2)*x*pi) + (1/2)*forch*rho*sqrt(sin((1/4)*x*pi)^2*cos((1/2)*y*pi)^2 + sin((1/2)*y*pi)^2*cos((1/2)*x*pi)^2)*sin((1/2)*y*pi)*cos((1/2)*x*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) - mu*(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)*(0.1*pi^2*sin((1/4)*x*pi)*sin((1/4)*y*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.025*pi^2*sin((1/4)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*x*pi)*cos((1/4)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.01*pi^2*sin((1/4)*x*pi)^2*sin((1/4)*y*pi)^2*sin((1/2)*y*pi)*cos((1/2)*x*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^3 - 1/4*pi^2*sin((1/2)*y*pi)*cos((1/2)*x*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)) - mu*(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)*(0.025*pi^2*sin((1/4)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*x*pi)*cos((1/4)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.1*pi^2*sin((1/2)*x*pi)*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/4)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.01*pi^2*sin((1/2)*y*pi)*cos((1/4)*x*pi)^2*cos((1/2)*x*pi)*cos((1/4)*y*pi)^2/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^3 - 1/4*pi^2*sin((1/2)*y*pi)*cos((1/2)*x*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)) + 0.025*pi*mu*(0.1*pi*sin((1/4)*x*pi)*sin((1/4)*y*pi)*sin((1/2)*y*pi)*cos((1/2)*x*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + (1/2)*pi*cos((1/2)*x*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5))*sin((1/4)*x*pi)*sin((1/4)*y*pi) - 0.025*pi*mu*(-0.1*pi*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/4)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - 1/2*pi*sin((1/2)*x*pi)*sin((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5))*cos((1/4)*x*pi)*cos((1/4)*y*pi) + 0.1*pi*rho*sin((1/4)*x*pi)^2*sin((1/4)*y*pi)*sin((1/2)*y*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - 0.1*pi*rho*sin((1/2)*y*pi)^2*cos((1/4)*x*pi)*cos((1/2)*x*pi)^2*cos((1/4)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - 1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi)^2*cos((1/2)*x*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) + (1/2)*pi*rho*sin((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi)^2/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) - pi*rho*sin((1/2)*x*pi)*sin((1/2)*y*pi)^2*cos((1/2)*x*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) - 1/4*pi*(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)*sin((1/4)*x*pi)*sin((3/2)*y*pi)'
    symbol_names = 'mu rho darcy forch'
    symbol_values = '${mu} ${rho} ${darcy} ${forch}'
  []
  [exact_v]
    type = ParsedFunction
    expression = 'sin((1/4)*x*pi)*cos((1/2)*y*pi)'
  []
  [forcing_v]
    type = ParsedFunction
    expression = 'darcy*mu*sin((1/4)*x*pi)*cos((1/2)*y*pi) + (1/2)*forch*rho*sqrt(sin((1/4)*x*pi)^2*cos((1/2)*y*pi)^2 + sin((1/2)*y*pi)^2*cos((1/2)*x*pi)^2)*sin((1/4)*x*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) - mu*(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)*(-0.1*pi^2*sin((1/4)*x*pi)^2*sin((1/4)*y*pi)*sin((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.025*pi^2*sin((1/4)*x*pi)^2*cos((1/4)*y*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.01*pi^2*sin((1/4)*x*pi)^3*sin((1/4)*y*pi)^2*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^3 - 1/4*pi^2*sin((1/4)*x*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)) - mu*(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)*(0.025*pi^2*sin((1/4)*x*pi)^2*cos((1/4)*y*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - 0.05*pi^2*cos((1/4)*x*pi)^2*cos((1/4)*y*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + 0.01*pi^2*sin((1/4)*x*pi)*cos((1/4)*x*pi)^2*cos((1/4)*y*pi)^2*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^3 - 1/16*pi^2*sin((1/4)*x*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)) + 0.025*pi*mu*(0.1*pi*sin((1/4)*x*pi)^2*sin((1/4)*y*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - 1/2*pi*sin((1/4)*x*pi)*sin((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5))*sin((1/4)*x*pi)*sin((1/4)*y*pi) - 0.025*pi*mu*(-0.1*pi*sin((1/4)*x*pi)*cos((1/4)*x*pi)*cos((1/4)*y*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 + (1/4)*pi*cos((1/4)*x*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5))*cos((1/4)*x*pi)*cos((1/4)*y*pi) + 0.1*pi*rho*sin((1/4)*x*pi)^3*sin((1/4)*y*pi)*cos((1/2)*y*pi)^2/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - 0.1*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/4)*y*pi)*cos((1/2)*y*pi)/(0.2*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 1)^2 - pi*rho*sin((1/4)*x*pi)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) - 1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) + (1/4)*pi*rho*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi)/(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5) + (3/2)*pi*(0.1*sin((1/4)*x*pi)*cos((1/4)*y*pi) + 0.5)*cos((1/4)*x*pi)*cos((3/2)*y*pi)'
    symbol_names = 'mu rho darcy forch'
    symbol_values = '${mu} ${rho} ${darcy} ${forch}'
  []
  [exact_p]
    type = ParsedFunction
    expression = 'sin((3/2)*y*pi)*cos((1/4)*x*pi)'
  []
  [forcing_p]
    type = ParsedFunction
    expression = '-1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi) - 1/2*pi*rho*sin((1/2)*x*pi)*sin((1/2)*y*pi)'
    symbol_names = 'rho'
    symbol_values = '${rho}'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
  petsc_options_value = 'asm      100                lu           NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = false
  csv = true
[]
[Postprocessors]
  [h]
    type = AverageElementSize
    outputs = 'console csv'
    execute_on = 'timestep_end'
  []
  [L2u]
    type = ElementL2FunctorError
    approximate = u
    exact = exact_u
    outputs = 'console csv'
    execute_on = 'timestep_end'
  []
  [L2v]
    type = ElementL2FunctorError
    approximate = v
    exact = exact_v
    outputs = 'console csv'
    execute_on = 'timestep_end'
  []
  [L2p]
    type = ElementL2FunctorError
    approximate = pressure
    exact = exact_p
    outputs = 'console csv'
    execute_on = 'timestep_end'
  []
[]
(modules/navier_stokes/examples/solidification/gallium_melting.i)
##########################################################
# Simulation of Gallium Melting Experiment
# Ref: Gau, C., & Viskanta, R. (1986). Melting and solidification of a pure metal on a vertical wall.
# Key physics: melting/solidification, convective heat transfer, natural convection
##########################################################
mu = 1.81e-3
rho_solid = 6093
rho_liquid = 6093
k_solid = 32
k_liquid = 32
cp_solid = 381.5
cp_liquid = 381.5
L = 80160
alpha_b = 1.2e-4
T_solidus = 302.93
T_liquidus = '${fparse T_solidus + 0.1}'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
T_cold = 301.15
T_hot = 311.15
Nx = 100
Ny = 50
[GlobalParams]
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 88.9e-3
    ymin = 0
    ymax = 63.5e-3
    nx = ${Nx}
    ny = ${Ny}
  []
[]
[AuxVariables]
  [U]
    type = MooseVariableFVReal
  []
  [fl]
    type = MooseVariableFVReal
    initial_condition = 0.0
  []
  [density]
    type = MooseVariableFVReal
  []
  [th_cond]
    type = MooseVariableFVReal
  []
  [cp_var]
    type = MooseVariableFVReal
  []
  [darcy_coef]
    type = MooseVariableFVReal
  []
  [fch_coef]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [mag]
    type = VectorMagnitudeAux
    variable = U
    x = vel_x
    y = vel_y
  []
  [compute_fl]
    type = NSLiquidFractionAux
    variable = fl
    temperature = T
    T_liquidus = '${T_liquidus}'
    T_solidus = '${T_solidus}'
    execute_on = 'TIMESTEP_END'
  []
  [rho_out]
    type = FunctorAux
    functor = 'rho_mixture'
    variable = 'density'
  []
  [th_cond_out]
    type = FunctorAux
    functor = 'k_mixture'
    variable = 'th_cond'
  []
  [cp_out]
    type = FunctorAux
    functor = 'cp_mixture'
    variable = 'cp_var'
  []
  [darcy_out]
    type = FunctorAux
    functor = 'Darcy_coefficient'
    variable = 'darcy_coef'
  []
  [fch_out]
    type = FunctorAux
    functor = 'Forchheimer_coefficient'
    variable = 'fch_coef'
  []
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 0.0
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 0.0
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [lambda]
    family = SCALAR
    order = FIRST
  []
  [T]
    type = INSFVEnergyVariable
    initial_condition = '${T_cold}'
    scaling = 1e-4
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = rho_mixture
  []
  [mean_zero_pressure]
    type = FVIntegralValueConstraint
    variable = pressure
    lambda = lambda
    phi0 = 0.0
  []
  [u_time]
    type = INSFVMomentumTimeDerivative
    variable = vel_x
    rho = rho_mixture
    momentum_component = 'x'
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = rho_mixture
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = ${mu}
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = vel_x
    momentum_component = 'x'
    u = vel_x
    v = vel_y
    Darcy_name = 'Darcy_coeff'
    Forchheimer_name = 'Forchheimer_coeff'
    rho = ${rho_liquid}
    mu = ${mu}
    standard_friction_formulation = false
  []
  [u_buoyancy]
    type = INSFVMomentumBoussinesq
    variable = vel_x
    T_fluid = T
    gravity = '0 -9.81 0'
    rho = '${rho_liquid}'
    ref_temperature = ${T_cold}
    momentum_component = 'x'
  []
  [u_gravity]
    type = INSFVMomentumGravity
    variable = vel_x
    gravity = '0 -9.81 0'
    rho = '${rho_liquid}'
    momentum_component = 'x'
  []
  [v_time]
    type = INSFVMomentumTimeDerivative
    variable = vel_y
    rho = rho_mixture
    momentum_component = 'y'
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = rho_mixture
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = ${mu}
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = vel_y
    momentum_component = 'y'
    u = vel_x
    v = vel_y
    Darcy_name = 'Darcy_coeff'
    Forchheimer_name = 'Forchheimer_coeff'
    rho = ${rho_liquid}
    mu = ${mu}
    standard_friction_formulation = false
  []
  [v_buoyancy]
    type = INSFVMomentumBoussinesq
    variable = vel_y
    T_fluid = T
    gravity = '0 -9.81 0'
    rho = '${rho_liquid}'
    ref_temperature = ${T_cold}
    momentum_component = 'y'
  []
  [v_gravity]
    type = INSFVMomentumGravity
    variable = vel_y
    gravity = '0 -9.81 0'
    rho = '${rho_liquid}'
    momentum_component = 'y'
  []
  [T_time]
    type = INSFVEnergyTimeDerivative
    variable = T
    rho = rho_mixture
    dh_dt = dh_dt
  []
  [energy_advection]
    type = INSFVEnergyAdvection
    variable = T
    velocity_interp_method = ${velocity_interp_method}
    advected_interp_method = ${advected_interp_method}
  []
  [energy_diffusion]
    type = FVDiffusion
    coeff = k_mixture
    variable = T
  []
  [energy_source]
    type = NSFVPhaseChangeSource
    variable = T
    L = ${L}
    liquid_fraction = fl
    T_liquidus = ${T_liquidus}
    T_solidus = ${T_solidus}
    rho = 'rho_mixture'
  []
[]
[FVBCs]
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'left right top bottom'
    variable = vel_x
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'left right top bottom'
    variable = vel_y
    function = 0
  []
  [hot_wall]
    type = FVDirichletBC
    variable = T
    value = '${T_hot}'
    boundary = 'left'
  []
  [cold_wall]
    type = FVDirichletBC
    variable = T
    value = '${T_cold}'
    boundary = 'right'
  []
[]
[FunctorMaterials]
  [ins_fv]
    type = INSFVEnthalpyFunctorMaterial
    rho = rho_mixture
    cp = cp_mixture
    temperature = 'T'
  []
  [eff_cp]
    type = NSFVMixtureFunctorMaterial
    phase_2_names = '${cp_solid} ${k_solid} ${rho_solid}'
    phase_1_names = '${cp_liquid} ${k_liquid} ${rho_liquid}'
    prop_names = 'cp_mixture k_mixture rho_mixture'
    phase_1_fraction = fl
  []
  [mushy_zone_resistance]
    type = INSFVMushyPorousFrictionFunctorMaterial
    liquid_fraction = 'fl'
    mu = '${mu}'
    rho_l = '${rho_liquid}'
    dendrite_spacing_scaling = 1e-1
  []
  [friction]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coeff Forchheimer_coeff'
    prop_values = 'darcy_coef darcy_coef darcy_coef fch_coef fch_coef fch_coef'
  []
  [const_functor]
    type = ADGenericFunctorMaterial
    prop_names = 'alpha_b'
    prop_values = '${alpha_b}'
  []
[]
[Executioner]
  type = Transient
  # Time-stepping parameters
  start_time = 0.0
  end_time = 200.0
  num_steps = 2
  [TimeStepper]
    type = IterationAdaptiveDT
    # Raise time step often but not by as much
    # There's a rough spot for convergence near 10% fluid fraction
    optimal_iterations = 15
    growth_factor = 1.5
    dt = 0.1
  []
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  nl_rel_tol = 1e-6
  nl_max_its = 30
  line_search = 'none'
[]
[Postprocessors]
  [ave_p]
    type = ElementAverageValue
    variable = 'pressure'
    execute_on = 'INITIAL TIMESTEP_END'
  []
  [ave_fl]
    type = ElementAverageValue
    variable = 'fl'
    execute_on = 'INITIAL TIMESTEP_END'
  []
  [ave_T]
    type = ElementAverageValue
    variable = 'T'
    execute_on = 'INITIAL TIMESTEP_END'
  []
[]
[VectorPostprocessors]
  [vel_x]
    type = ElementValueSampler
    variable = 'vel_x fl'
    sort_by = 'x'
  []
[]
[Outputs]
  exodus = true
  csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/discontinuous-body-forces.i)
mu = 1.1
rho = 1.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
  two_term_boundary_expansion = true
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = u
    v = v
    pressure = pressure
  []
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 10
    ymin = -1
    ymax = 1
    nx = 100
    ny = 9
  []
  [subdomain]
    type = SubdomainBoundingBoxGenerator
    bottom_left = '5 -1 0'
    top_right = '10 1 0'
    block_id = 1
    input = gen
  []
[]
[Variables]
  [u]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [v]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
[]
[FVKernels]
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = u
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = u
    mu = ${mu}
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = u
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction_linear]
    type = PINSFVMomentumFriction
    variable = u
    Darcy_name = friction_coefficient
    momentum_component = 'x'
    block = '1'
    standard_friction_formulation = false
    rho = ${rho}
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = v
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = v
    mu = ${mu}
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction_linear]
    type = PINSFVMomentumFriction
    variable = v
    Darcy_name = friction_coefficient
    momentum_component = 'y'
    block = '1'
    standard_friction_formulation = false
    rho = ${rho}
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = '0'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = u
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = v
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = '0'
  []
[]
[FunctorMaterials]
  [friction_coefficient]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'friction_coefficient'
    prop_values = '25 25 25'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/materials/ergun/ergun.i)
# This file simulates flow of fluid in a porous elbow for the purpose of verifying
# correct implementation of the various different solution variable sets. This input
# tests correct implementation of the primitive superficial variable set. Flow enters on the top
# and exits on the right. Because the purpose is only to test the equivalence of
# different equation sets, no solid energy equation is included.
porosity_left = 0.4
porosity_right = 0.6
pebble_diameter = 0.06
mu = 1.81e-5 # This has been increased to avoid refining the mesh
M = 28.97e-3
R = 8.3144598
# inlet mass flowrate, kg/s
mdot = -10.0
# inlet mass flux (superficial)
mflux_in_superficial = ${fparse mdot / (pi * 0.5 * 0.5)}
# inlet mass flux (interstitial)
mflux_in_interstitial = ${fparse mflux_in_superficial / porosity_left}
p_initial = 201325.0
T_initial = 300.0
rho_initial = ${fparse p_initial / T_initial * M / R}
vel_y_initial = ${fparse mflux_in_interstitial / rho_initial}
vel_x_initial = 0.0
superficial_vel_y_initial = ${fparse mflux_in_superficial / rho_initial}
superficial_vel_x_initial = 1e-12
# Computation parameters
velocity_interp_method = 'rc'
advected_interp_method = 'upwind'
# ==============================================================================
# GEOMETRY AND MESH
# ==============================================================================
[Mesh]
  [fmg]
    type = FileMeshGenerator
    file = 'ergun_in.e'
  []
  coord_type = RZ
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolator
    u = superficial_vel_x
    v = superficial_vel_y
    pressure = pressure
    porosity = porosity
  []
[]
[GlobalParams]
  porosity = porosity
  pebble_diameter = ${pebble_diameter}
  fp = fp
  # rho for the kernels. Must match fluid property!
  rho = ${rho_initial}
  fv = true
  velocity_interp_method = ${velocity_interp_method}
  advected_interp_method = ${advected_interp_method}
  # behavior at time of test creation
  two_term_boundary_expansion = false
  rhie_chow_user_object = 'rc'
[]
# ==============================================================================
# VARIABLES AND KERNELS
# ==============================================================================
[Variables]
  [pressure]
    type = INSFVPressureVariable
    initial_condition = ${p_initial}
  []
  [superficial_vel_x]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = ${superficial_vel_x_initial}
  []
  [superficial_vel_y]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = ${superficial_vel_y_initial}
  []
[]
[FVKernels]
  # Mass Equation.
  [mass]
    type = PINSFVMassAdvection
    variable = 'pressure'
  []
  # Momentum x component equation.
  [vel_x_time]
    type = PINSFVMomentumTimeDerivative
    variable = 'superficial_vel_x'
    momentum_component = 'x'
  []
  [vel_x_advection]
    type = PINSFVMomentumAdvection
    variable = 'superficial_vel_x'
    momentum_component = 'x'
  []
  [vel_x_viscosity]
    type = PINSFVMomentumDiffusion
    variable = 'superficial_vel_x'
    momentum_component = 'x'
    mu = 'mu'
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = 'superficial_vel_x'
    pressure = pressure
    momentum_component = 'x'
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = 'superficial_vel_x'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    momentum_component = 'x'
    speed = speed
    mu = 'mu'
  []
  # Momentum y component equation.
  [vel_y_time]
    type = PINSFVMomentumTimeDerivative
    variable = 'superficial_vel_y'
    momentum_component = 'y'
  []
  [vel_y_advection]
    type = PINSFVMomentumAdvection
    variable = 'superficial_vel_y'
    momentum_component = 'y'
  []
  [vel_y_viscosity]
    type = PINSFVMomentumDiffusion
    variable = 'superficial_vel_y'
    momentum_component = 'y'
    mu = 'mu'
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = 'superficial_vel_y'
    pressure = pressure
    momentum_component = 'y'
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = 'superficial_vel_y'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    momentum_component = 'y'
    mu = 'mu'
    speed = speed
  []
  [gravity]
    type = PINSFVMomentumGravity
    variable = 'superficial_vel_y'
    gravity = '0 -9.81 0'
    momentum_component = 'y'
  []
[]
# ==============================================================================
# AUXVARIABLES AND AUXKERNELS
# ==============================================================================
[AuxVariables]
  [T_fluid]
    initial_condition = ${T_initial}
    order = CONSTANT
    family = MONOMIAL
  []
  [vel_x]
    initial_condition = ${fparse vel_x_initial}
    order = CONSTANT
    family = MONOMIAL
  []
  [vel_y]
    initial_condition = ${fparse vel_y_initial}
    order = CONSTANT
    family = MONOMIAL
  []
  [porosity_out]
    type = MooseVariableFVReal
  []
[]
[AuxKernels]
  [vel_x]
    type = FunctorAux
    variable = vel_x
    functor = vel_x_mat
  []
  [vel_y]
    type = FunctorAux
    variable = vel_y
    functor = vel_y_mat
  []
  [porosity_out]
    type = FunctorAux
    variable = porosity_out
    functor = porosity
  []
[]
# ==============================================================================
# FLUID PROPERTIES, MATERIALS AND USER OBJECTS
# ==============================================================================
[FluidProperties]
  [fp]
    type = IdealGasFluidProperties
    k = 0.0
    mu = ${mu}
    gamma = 1.4
    molar_mass = ${M}
  []
[]
[FunctorMaterials]
  [enthalpy]
    type = INSFVEnthalpyMaterial
    temperature = 'T_fluid'
  []
  [speed]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = 'superficial_vel_x'
    superficial_vel_y = 'superficial_vel_y'
    porosity = porosity
    vel_x = vel_x_mat
    vel_y = vel_y_mat
  []
  [kappa]
    type = FunctorKappaFluid
  []
  [const_Fdrags_mat]
    type = FunctorErgunDragCoefficients
    porosity = porosity
  []
  [fluidprops]
    type = GeneralFunctorFluidProps
    mu_rampdown = mu_func
    porosity = porosity
    characteristic_length = ${pebble_diameter}
    T_fluid = 'T_fluid'
    pressure = 'pressure'
    speed = 'speed'
  []
[]
d = 0.05
[Functions]
  [mu_func]
    type = PiecewiseLinear
    x = '1 3 5 10 15 20'
    y = '1e5 1e4 1e3 1e2 1e1 1'
  []
  [real_porosity_function]
    type = ParsedFunction
    expression = 'if (x < 0.6 - ${d}, ${porosity_left}, if (x > 0.6 + ${d}, ${porosity_right},
        (x-(0.6-${d}))/(2*${d})*(${porosity_right}-${porosity_left}) + ${porosity_left}))'
  []
  [porosity]
    type = ParsedFunction
    expression = 'if (x < 0.6 - ${d}, ${porosity_left}, if (x > 0.6 + ${d}, ${porosity_right},
        (x-(0.6-${d}))/(2*${d})*(${porosity_right}-${porosity_left}) + ${porosity_left}))'
  []
[]
# ==============================================================================
# BOUNDARY CONDITIONS
# ==============================================================================
[FVBCs]
  [outlet_p]
    type = INSFVOutletPressureBC
    variable = 'pressure'
    function = ${p_initial}
    boundary = 'right'
  []
  ## No or Free slip BC
  [free-slip-wall-x]
    type = INSFVNaturalFreeSlipBC
    boundary = 'bottom wall_1 wall_2 left'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [free-slip-wall-y]
    type = INSFVNaturalFreeSlipBC
    boundary = 'bottom wall_1 wall_2 left'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  ## Symmetry
  [symmetry-x]
    type = PINSFVSymmetryVelocityBC
    boundary = 'left'
    variable = superficial_vel_x
    u = superficial_vel_x
    v = superficial_vel_y
    mu = 'mu'
    momentum_component = 'x'
  []
  [symmetry-y]
    type = PINSFVSymmetryVelocityBC
    boundary = 'left'
    variable = superficial_vel_y
    u = superficial_vel_x
    v = superficial_vel_y
    mu = 'mu'
    momentum_component = 'y'
  []
  [symmetry-p]
    type = INSFVSymmetryPressureBC
    boundary = 'left'
    variable = 'pressure'
  []
  ## inlet
  [inlet_vel_x]
    type = INSFVInletVelocityBC
    variable = 'superficial_vel_x'
    function = ${superficial_vel_x_initial}
    boundary = 'top'
  []
  [inlet_vel_y]
    type = INSFVInletVelocityBC
    variable = 'superficial_vel_y'
    function = ${superficial_vel_y_initial}
    boundary = 'top'
  []
[]
# ==============================================================================
# EXECUTION PARAMETERS
# ==============================================================================
[Executioner]
  type = Transient
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
  petsc_options_value = 'asm      lu           NONZERO                   200'
  line_search = 'none'
  # Problem time parameters
  dtmin = 0.01
  dtmax = 2000
  end_time = 3000
  # must be the same as the fluid
  # Iterations parameters
  l_max_its = 50
  l_tol     = 1e-8
  nl_max_its = 25
  # nl_rel_tol = 5e-7
  nl_abs_tol = 2e-7
  # Automatic scaling
  automatic_scaling = true
  verbose = true
  [TimeStepper]
    type = IterationAdaptiveDT
    dt                 = 0.025
    cutback_factor     = 0.5
    growth_factor      = 2.0
  []
  # Steady state detection.
  steady_state_detection = true
  steady_state_tolerance = 1e-7
  steady_state_start_time = 400
[]
# ==============================================================================
# POSTPROCESSORS DEBUG AND OUTPUTS
# ==============================================================================
[Postprocessors]
  [mass_flow_in]
    type = VolumetricFlowRate
    boundary = 'top'
    vel_x = 'superficial_vel_x'
    vel_y = 'superficial_vel_y'
    advected_quantity = ${rho_initial}
    execute_on = 'INITIAL TIMESTEP_END'
  []
  [mass_flow_out]
    type = VolumetricFlowRate
    boundary = 'right'
    vel_x = 'superficial_vel_x'
    vel_y = 'superficial_vel_y'
    advected_quantity = ${rho_initial}
    execute_on = 'INITIAL TIMESTEP_END'
  []
  [p_in]
    type = SideAverageValue
    variable = pressure
    boundary = 'top'
  []
  [dP]
    type = LinearCombinationPostprocessor
    pp_names = 'p_in'
    pp_coefs = '1.0'
    b = ${fparse -p_initial}
  []
[]
[Outputs]
  exodus = true
  print_linear_residuals = false
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/friction/2d-rc-friction.i)
mu = 1.1
rho = 1.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = INSFVRhieChowInterpolator
    u = vel_x
    v = vel_y
    pressure = pressure
  []
[]
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 5
    ymin = -1
    ymax = 1
    nx = 50
    ny = 10
  []
[]
[Variables]
  [vel_x]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [vel_y]
    type = INSFVVelocityVariable
    initial_condition = 1
  []
  [pressure]
    type = INSFVPressureVariable
  []
[]
[FVKernels]
  inactive = 'u_friction_quad v_friction_quad'
  [mass]
    type = INSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [u_advection]
    type = INSFVMomentumAdvection
    variable = vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'x'
  []
  [u_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_x
    mu = ${mu}
    momentum_component = 'x'
  []
  [u_pressure]
    type = INSFVMomentumPressure
    variable = vel_x
    momentum_component = 'x'
    pressure = pressure
  []
  [u_friction_linear]
    type = PINSFVMomentumFriction
    variable = vel_x
    Darcy_name = friction_coefficient
    momentum_component = 'x'
    rho = ${rho}
    standard_friction_formulation = false
  []
  [u_friction_quad]
    type = PINSFVMomentumFriction
    variable = vel_x
    speed = speed
    Forchheimer_name = friction_coefficient
    momentum_component = 'x'
    rho = ${rho}
    standard_friction_formulation = false
  []
  [v_advection]
    type = INSFVMomentumAdvection
    variable = vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    momentum_component = 'y'
  []
  [v_viscosity]
    type = INSFVMomentumDiffusion
    variable = vel_y
    mu = ${mu}
    momentum_component = 'y'
  []
  [v_pressure]
    type = INSFVMomentumPressure
    variable = vel_y
    momentum_component = 'y'
    pressure = pressure
  []
  [v_friction_linear]
    type = PINSFVMomentumFriction
    variable = vel_y
    Darcy_name = friction_coefficient
    momentum_component = 'y'
    rho = ${rho}
    standard_friction_formulation = false
  []
  [v_friction_quad]
    type = PINSFVMomentumFriction
    variable = vel_y
    speed = speed
    Forchheimer_name = friction_coefficient
    momentum_component = 'y'
    rho = ${rho}
    standard_friction_formulation = false
  []
[]
[FVBCs]
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_x
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = vel_y
    function = '0'
  []
  [walls-u]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_x
    function = 0
  []
  [walls-v]
    type = INSFVNoSlipWallBC
    boundary = 'top bottom'
    variable = vel_y
    function = 0
  []
  [outlet_p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = '0'
  []
[]
[FunctorMaterials]
  inactive = exponential_friction_coefficient
  [friction_coefficient]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'friction_coefficient'
    prop_values = '25 25 25'
  []
  [speed_material]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = vel_x
    superficial_vel_y = vel_y
    porosity = 1
    vel_x = vel_x_mat
    vel_y = vel_y_mat
  []
  [Re_material]
    type = ReynoldsNumberFunctorMaterial
    speed = speed
    characteristic_length = 2
    rho = ${rho}
    mu = ${mu}
  []
  [exponential_coeff]
    type = ExponentialFrictionMaterial
    friction_factor_name = 'exponential_coeff'
    Re = Re
    c1 = 0.25
    c2 = 0.55
  []
  [exponential_friction_coefficient]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'friction_coefficient'
    prop_values = 'exponential_coeff exponential_coeff exponential_coeff'
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  nl_rel_tol = 1e-12
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc-friction.i)
mu = 1.1
rho = 1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[Mesh]
  [mesh]
    type = CartesianMeshGenerator
    dim = 2
    dx = '2.5 2.5'
    dy = '1.0'
    ix = '20 20'
    iy = '20'
    subdomain_id = '1 2'
  []
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolator
    u = superficial_vel_x
    v = superficial_vel_y
    pressure = pressure
    porosity = porosity
  []
[]
[Variables]
  inactive = 'lambda'
  [superficial_vel_x]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [superficial_vel_y]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1e-6
  []
  [pressure]
    type = INSFVPressureVariable
  []
  [lambda]
    family = SCALAR
    order = FIRST
  []
[]
[AuxVariables]
  [porosity]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_condition = 0.5
  []
[]
[FVKernels]
  inactive = 'mean-pressure'
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [u_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_x
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_x
    mu = ${mu}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_x
    momentum_component = 'x'
    pressure = pressure
    porosity = porosity
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_x
    momentum_component = 'x'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    mu = ${mu}
    rho = ${rho}
    speed = speed
  []
  [v_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_y
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_y
    mu = ${mu}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_y
    momentum_component = 'y'
    pressure = pressure
    porosity = porosity
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_y
    momentum_component = 'y'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [mean-pressure]
    type = FVIntegralValueConstraint
    variable = pressure
    lambda = lambda
    phi0 = 0.01
  []
[]
[FVBCs]
  inactive = 'free-slip-u free-slip-v'
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = superficial_vel_x
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = superficial_vel_y
    function = 0
  []
  [no-slip-u]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = superficial_vel_x
    function = 0
  []
  [no-slip-v]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = superficial_vel_y
    function = 0
  []
  [free-slip-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [free-slip-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  [symmetry-u]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = superficial_vel_x
    u = superficial_vel_x
    v = superficial_vel_y
    mu = ${mu}
    momentum_component = 'x'
  []
  [symmetry-v]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = superficial_vel_y
    u = superficial_vel_x
    v = superficial_vel_y
    mu = ${mu}
    momentum_component = 'y'
  []
  [symmetry-p]
    type = INSFVSymmetryPressureBC
    boundary = 'bottom'
    variable = pressure
  []
  [outlet-p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0
  []
[]
[FunctorMaterials]
  [darcy]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coefficient Forchheimer_coefficient'
    prop_values = '0.1 0.1 0.1 0.1 0.1 0.1'
  []
  [speec]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = superficial_vel_x
    superficial_vel_y = superficial_vel_y
    porosity = porosity
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  nl_rel_tol = 1e-11
  nl_abs_tol = 1e-14
[]
# Some basic Postprocessors to visually examine the solution
[Postprocessors]
  [inlet-p]
    type = SideAverageValue
    variable = pressure
    boundary = 'left'
  []
  [outlet-u]
    type = SideIntegralVariablePostprocessor
    variable = superficial_vel_x
    boundary = 'right'
  []
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/segregated/2d-momentum.i)
mu = 1.1
rho = 1.1
pressure_tag = "pressure_grad"
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 5
    ymin = 0
    ymax = 1
    nx = 40
    ny = 6
  []
[]
[GlobalParams]
  advected_interp_method = 'average'
  velocity_interp_method = 'rc'
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolatorSegregated
    u = superficial_vel_x
    v = superficial_vel_y
    pressure = pressure
    porosity = porosity
  []
[]
[Problem]
  nl_sys_names = 'u_system v_system pressure_system'
  previous_nl_solution_required = true
[]
[Variables]
  [superficial_vel_x]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
    solver_sys = u_system
    two_term_boundary_expansion = false
  []
  [superficial_vel_y]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1e-6
    solver_sys = v_system
    two_term_boundary_expansion = false
  []
  [pressure]
    type = INSFVPressureVariable
    two_term_boundary_expansion = false
    solver_sys = pressure_system
  []
[]
[AuxVariables]
  [porosity]
    type = MooseVariableFVReal
    initial_condition = 0.5
  []
[]
[FVKernels]
  inactive = "u_friction v_friction"
  [u_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_x
    rho = ${rho}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_x
    mu = ${mu}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_x
    momentum_component = 'x'
    pressure = pressure
    porosity = porosity
    extra_vector_tags = ${pressure_tag}
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_x
    momentum_component = 'y'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [v_advection]
    type = PINSFVMomentumAdvection
    variable = superficial_vel_y
    rho = ${rho}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = superficial_vel_y
    mu = ${mu}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_pressure]
    type = PINSFVMomentumPressure
    variable = superficial_vel_y
    momentum_component = 'y'
    pressure = pressure
    porosity = porosity
    extra_vector_tags = ${pressure_tag}
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = superficial_vel_y
    momentum_component = 'y'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [p_diffusion]
    type = FVAnisotropicDiffusion
    variable = pressure
    coeff = "Ainv"
    coeff_interp_method = 'average'
  []
  [p_source]
    type = FVDivergence
    variable = pressure
    vector_field = "HbyA"
    force_boundary_execution = true
  []
[]
[FVBCs]
  inactive = 'slip-u slip-v'
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = superficial_vel_x
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = superficial_vel_y
    function = 0
  []
  [no-slip-u]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = superficial_vel_x
    function = 0
  []
  [no-slip-v]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = superficial_vel_y
    function = 0
  []
  [symmetry-u]
    type = INSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = superficial_vel_x
    u = superficial_vel_x
    v = superficial_vel_y
    mu = ${mu}
    momentum_component = 'x'
  []
  [symmetry-v]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = superficial_vel_y
    u = superficial_vel_x
    v = superficial_vel_y
    mu = ${mu}
    momentum_component = 'y'
  []
  [symmetry-p]
    type = INSFVSymmetryPressureBC
    boundary = 'bottom'
    variable = pressure
  []
  [outlet-p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0.4
  []
  ### Are disabled by default but we switch it on for certain tests ###
  [slip-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = superficial_vel_x
    momentum_component = 'x'
  []
  [slip-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = superficial_vel_y
    momentum_component = 'y'
  []
  #####################################################################
[]
[FunctorMaterials]
  [darcy]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coefficient Forchheimer_coefficient'
    prop_values = '0.01 0.02 0.03 0.01 0.02 0.03'
  []
  [speed]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = superficial_vel_x
    superficial_vel_y = superficial_vel_y
    porosity = porosity
  []
[]
[Executioner]
  type = SIMPLENonlinearAssembly
  momentum_l_abs_tol = 1e-14
  pressure_l_abs_tol = 1e-14
  momentum_l_tol = 0
  pressure_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  pressure_gradient_tag = ${pressure_tag}
  momentum_equation_relaxation = 0.85
  pressure_variable_relaxation = 0.45
  num_iterations = 150
  pressure_absolute_tolerance = 1e-13
  momentum_absolute_tolerance = 1e-13
  print_fields = false
  continue_on_max_its = true
[]
[Outputs]
  exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/pins/channel-flow/2d-rc-rz-by-parts.i)
mu = 1.1
rho = 1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    xmin = 0
    xmax = 5
    ymin = 0
    ymax = 1
    nx = 40
    ny = 10
  []
  coord_type = 'RZ'
  rz_coord_axis = 'X'
[]
[GlobalParams]
  rhie_chow_user_object = 'rc'
[]
[UserObjects]
  [rc]
    type = PINSFVRhieChowInterpolator
    u = u
    v = v
    pressure = pressure
    porosity = porosity
  []
[]
[Variables]
  [u]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1
  []
  [v]
    type = PINSFVSuperficialVelocityVariable
    initial_condition = 1e-6
  []
  [pressure]
    type = INSFVPressureVariable
  []
[]
[AuxVariables]
  [porosity]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_condition = 0.5
  []
[]
[FVKernels]
  inactive = 'v_pressure_volumetric'
  [mass]
    type = PINSFVMassAdvection
    variable = pressure
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
  []
  [u_advection]
    type = PINSFVMomentumAdvection
    variable = u
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_viscosity]
    type = PINSFVMomentumDiffusion
    variable = u
    mu = ${mu}
    porosity = porosity
    momentum_component = 'x'
  []
  [u_pressure]
    type = PINSFVMomentumPressureFlux
    variable = u
    momentum_component = 'x'
    pressure = pressure
    porosity = porosity
  []
  [u_friction]
    type = PINSFVMomentumFriction
    variable = u
    momentum_component = 'x'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
  [v_advection]
    type = PINSFVMomentumAdvection
    variable = v
    advected_interp_method = ${advected_interp_method}
    velocity_interp_method = ${velocity_interp_method}
    rho = ${rho}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_viscosity]
    type = PINSFVMomentumDiffusion
    variable = v
    mu = ${mu}
    porosity = porosity
    momentum_component = 'y'
  []
  [v_pressure_volumetric]
    type = PINSFVMomentumPressure
    variable = v
    momentum_component = 'y'
    pressure = pressure
    porosity = porosity
  []
  [v_pressure_by_parts_flux]
    type = PINSFVMomentumPressureFlux
    variable = v
    momentum_component = 'y'
    pressure = pressure
    porosity = porosity
  []
  [v_pressure_by_parts_volume_term]
    type = PNSFVMomentumPressureFluxRZ
    variable = v
    pressure = pressure
    porosity = porosity
    momentum_component = 'y'
  []
  [v_friction]
    type = PINSFVMomentumFriction
    variable = v
    momentum_component = 'y'
    Darcy_name = 'Darcy_coefficient'
    Forchheimer_name = 'Forchheimer_coefficient'
    rho = ${rho}
    speed = speed
    mu = ${mu}
  []
[]
[FVBCs]
  inactive = 'free-slip-u free-slip-v'
  [inlet-u]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = u
    function = '1'
  []
  [inlet-v]
    type = INSFVInletVelocityBC
    boundary = 'left'
    variable = v
    function = 0
  []
  [no-slip-u]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = u
    function = 0
  []
  [no-slip-v]
    type = INSFVNoSlipWallBC
    boundary = 'top'
    variable = v
    function = 0
  []
  [free-slip-u]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = u
    momentum_component = 'x'
  []
  [free-slip-v]
    type = INSFVNaturalFreeSlipBC
    boundary = 'top'
    variable = v
    momentum_component = 'y'
  []
  [symmetry-u]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = u
    u = u
    v = v
    mu = ${mu}
    momentum_component = 'x'
  []
  [symmetry-v]
    type = PINSFVSymmetryVelocityBC
    boundary = 'bottom'
    variable = v
    u = u
    v = v
    mu = ${mu}
    momentum_component = 'y'
  []
  [symmetry-p]
    type = INSFVSymmetryPressureBC
    boundary = 'bottom'
    variable = pressure
  []
  [outlet-p]
    type = INSFVOutletPressureBC
    boundary = 'right'
    variable = pressure
    function = 0
  []
[]
[FunctorMaterials]
  [darcy]
    type = ADGenericVectorFunctorMaterial
    prop_names = 'Darcy_coefficient Forchheimer_coefficient'
    prop_values = '0.1 0.1 0.1 0.1 0.1 0.1'
  []
  [speed]
    type = PINSFVSpeedFunctorMaterial
    superficial_vel_x = u
    superficial_vel_y = v
    porosity = porosity
  []
[]
[Executioner]
  type = Steady
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -pc_factor_shift_type'
  petsc_options_value = 'lu NONZERO'
  line_search = 'none'
  nl_rel_tol = 1e-11
  nl_abs_tol = 1e-14
[]
# Some basic Postprocessors to visually examine the solution
[Postprocessors]
  [inlet-p]
    type = SideAverageValue
    variable = pressure
    boundary = 'left'
  []
  [outlet-u]
    type = SideIntegralVariablePostprocessor
    variable = u
    boundary = 'right'
  []
[]
[Outputs]
  exodus = true
[]