- PorousFlowDictatorThe UserObject that holds the list of PorousFlow variable names
C++ Type:UserObjectName
Controllable:No
Description:The UserObject that holds the list of PorousFlow variable names
 - advective_flux_calculatorPorousFlowAdvectiveFluxCalculator UserObject. This determines whether the advection describes a movement of a fluid component in a fluid phase, or movement of heat energy in a fluid phase
C++ Type:UserObjectName
Controllable:No
Description:PorousFlowAdvectiveFluxCalculator UserObject. This determines whether the advection describes a movement of a fluid component in a fluid phase, or movement of heat energy in a fluid phase
 - variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
 
PorousFlowFluxLimitedTVDAdvection
This Kernel implements the advective flux  where  and  are defined via one of the PorousFlowAdvectiveFluxCalculator UserObjects.  See the PorousFlow Kuzmin-Turek page for details of to use this Kernel in PorousFlow simulations.
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
 - displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
 - 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)
 
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
 - extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
 - extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
 - matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
 - vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default: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.
 - diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
 - enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
 - implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
 - save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
 - search_methodnearest_node_connected_sidesChoice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes).
Default:nearest_node_connected_sides
C++ Type:MooseEnum
Options:nearest_node_connected_sides, all_proximate_sides
Controllable:No
Description:Choice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes).
 - seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
 - use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
 
Advanced Parameters
- 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/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_02.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_2D_angle.i)
 - (modules/porous_flow/test/tests/numerical_diffusion/pffltvd.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_01.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_03.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_2D.i)
 - (modules/porous_flow/test/tests/recover/pffltvd.i)
 - (modules/porous_flow/test/tests/pressure_pulse/pressure_pulse_1d_2phasePS_KT.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_05.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_3D.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_1D_adaptivity.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_2D_trimesh.i)
 - (modules/porous_flow/test/tests/heat_advection/heat_advection_1d_KT.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_1D.i)
 - (modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_04.i)
 
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_02.i)
# Checking the Jacobian of Flux-Limited TVD Advection, 1 phase, 3 components, unsaturated, using flux_limiter_type = none
[Mesh]
  type = GeneratedMesh
  dim = 3
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = -1
  ymax = 2
[]
[GlobalParams]
  gravity = '1 2 -0.5'
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
  [tracer0]
  []
  [tracer1]
  []
[]
[ICs]
  [pp]
    variable = pp
    type = RandomIC
    min = -1
    max = 0
  []
  [tracer0]
    variable = tracer0
    type = RandomIC
    min = 0
    max = 1
  []
  [tracer1]
    variable = tracer1
    type = RandomIC
    min = 0
    max = 1
  []
[]
[Kernels]
  [fluxpp]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = pp
    advective_flux_calculator = advective_flux_calculator_0
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer0
    advective_flux_calculator = advective_flux_calculator_1
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer1
    advective_flux_calculator = advective_flux_calculator_2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1
    density0 = 0.4
    viscosity = 1.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp tracer0 tracer1'
    number_fluid_phases = 1
    number_fluid_components = 3
  []
  [pc]
    type = PorousFlowCapillaryPressureVG
    alpha = 1
    m = 0.5
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    fluid_component = 1
  []
  [advective_flux_calculator_2]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    fluid_component = 2
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = pp
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'tracer0 tracer1'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityCorey
    phase = 0
    n = 2
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1.21 0 0  0 1.5 0  0 0 0.8'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-snes_type'
    petsc_options_value = 'test'
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 1
  num_steps = 1
  dt = 1
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_2D_angle.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, mploying PorousFlow Kernels and UserObjects, with superbee flux-limiter
# 2D version with velocity = (0.1, 0.2, 0)
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 10
  xmin = 0
  xmax = 1
  ny = 10
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x - 2 * y'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1 | x > 0.3 | y < 0.1 | y > 0.3, 0, 1)'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_boundary_porepressure]
    type = FunctionDirichletBC
    variable = porepressure
    function = '1 - x - 2 * y'
    boundary = 'left right top bottom'
  []
  [no_tracer_at_boundary]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = 'left right top bottom'
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 0.3
  dt = 0.1
[]
[Outputs]
  [out]
    type = Exodus
    execute_on = 'initial final'
  []
[]
(modules/porous_flow/test/tests/numerical_diffusion/pffltvd.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, employing PorousFlow Kernels and UserObjects, with superbee flux-limiter
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 100
  xmin = 0
  xmax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.3,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0 0'
    num_points = 101
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 6
  dt = 6E-2
  nl_abs_tol = 1E-8
  timestep_tolerance = 1E-3
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_01.i)
# Checking the Jacobian of Flux-Limited TVD Advection, 1 phase, 1 component, full saturation, using flux_limiter_type = none
# This is quite a heavy test, but we need a fairly big mesh to check the upwinding is happening correctly
[Mesh]
  type = GeneratedMesh
  dim = 3
  nx = 3
  xmin = 0
  xmax = 1
  ny = 4
  ymin = -1
  ymax = 2
  bias_y = 1.5
  nz = 4
  zmin = 1
  zmax = 2
  bias_z = 0.8
[]
[GlobalParams]
  gravity = '1 2 -0.5'
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
[]
[ICs]
  [pp]
    variable = pp
    type = RandomIC
    min = 1
    max = 2
  []
[]
[Kernels]
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = pp
    advective_flux_calculator = advective_flux_calculator
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1
    density0 = 0.4
    viscosity = 1.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp'
    number_fluid_phases = 1
    number_fluid_components = 1
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator]
    type = PorousFlowAdvectiveFluxCalculatorSaturated
    flux_limiter_type = None
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = pp
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1.21 0 0  0 1.5 0  0 0 0.8'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-snes_type'
    petsc_options_value = 'test'
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 1
  num_steps = 1
  dt = 1
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_03.i)
# Checking the Jacobian of Flux-Limited TVD Advection, 2 phases, 2 components, using flux_limiter_type = None
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 2
  xmin = 0
  xmax = 1
  ny = 2
  ymin = -1
  ymax = 2
  bias_y = 1.5
[]
[GlobalParams]
  gravity = '1 2 -0.5'
  PorousFlowDictator = dictator
[]
[Variables]
  [ppwater]
  []
  [ppgas]
  []
  [massfrac_ph0_sp0]
  []
  [massfrac_ph1_sp0]
  []
[]
[ICs]
  [ppwater]
    type = RandomIC
    variable = ppwater
    min = -1
    max = 0
  []
  [ppgas]
    type = RandomIC
    variable = ppgas
    min = 0
    max = 1
  []
  [massfrac_ph0_sp0]
    type = RandomIC
    variable = massfrac_ph0_sp0
    min = 0
    max = 1
  []
  [massfrac_ph1_sp0]
    type = RandomIC
    variable = massfrac_ph1_sp0
    min = 0
    max = 1
  []
[]
[Kernels]
  [flux_ph0_sp0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = ppwater
    advective_flux_calculator = advective_flux_calculator_ph0_sp0
  []
  [flux_ph0_sp1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = ppgas
    advective_flux_calculator = advective_flux_calculator_ph0_sp1
  []
  [flux_ph1_sp0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = massfrac_ph0_sp0
    advective_flux_calculator = advective_flux_calculator_ph1_sp0
  []
  [flux_ph1_sp1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = massfrac_ph1_sp0
    advective_flux_calculator = advective_flux_calculator_ph1_sp1
  []
[]
[FluidProperties]
  [simple_fluid0]
    type = SimpleFluidProperties
    bulk_modulus = 1.5
    density0 = 1
    thermal_expansion = 0
    viscosity = 1
  []
  [simple_fluid1]
    type = SimpleFluidProperties
    bulk_modulus = 0.5
    density0 = 0.5
    thermal_expansion = 0
    viscosity = 1.4
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'ppwater ppgas massfrac_ph0_sp0 massfrac_ph1_sp0'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureVG
    alpha = 1
    m = 0.5
  []
  [advective_flux_calculator_ph0_sp0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    phase = 0
    fluid_component = 0
  []
  [advective_flux_calculator_ph0_sp1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    phase = 0
    fluid_component = 1
  []
  [advective_flux_calculator_ph1_sp0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    phase = 1
    fluid_component = 0
  []
  [advective_flux_calculator_ph1_sp1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = None
    phase = 1
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow2PhasePP
    phase0_porepressure = ppwater
    phase1_porepressure = ppgas
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
  []
  [simple_fluid0]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid0
    phase = 0
  []
  [simple_fluid1]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid1
    phase = 1
  []
  [relperm0]
    type = PorousFlowRelativePermeabilityCorey
    n = 2
    phase = 0
  []
  [relperm1]
    type = PorousFlowRelativePermeabilityCorey
    n = 3
    phase = 1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1.21 0 0  0 1.5 0  0 0 0.8'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-snes_type'
    petsc_options_value = 'test'
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 1
  num_steps = 1
  dt = 1
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_2D.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, employing PorousFlow Kernels and UserObjects, with superbee flux-limiter
# 3D version
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 10
  xmin = 0
  xmax = 1
  ny = 4
  ymin = 0
  ymax = 0.5
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.3,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0.5 0'
    num_points = 11
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 6
  dt = 6E-2
  nl_abs_tol = 1E-8
  timestep_tolerance = 1E-3
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/recover/pffltvd.i)
# Tests that PorousFlow can successfully recover using a checkpoint file.
# This test contains stateful material properties, adaptivity, integrated
# boundary conditions with nodal-sized materials, and TVD flux limiting.
#
# This test file is run three times:
# 1) The full input file is run to completion
# 2) The input file is run for half the time and checkpointing is included
# 3) The input file is run in recovery using the checkpoint data
#
# The final output of test 3 is compared to the final output of test 1 to verify
# that recovery was successful.
[Mesh]
  [mesh]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 10
    xmin = 0
    xmax = 1
  []
  # To get consistent ordering of results with distributed meshes
  allow_renumbering = false
[]
[Adaptivity]
  initial_steps = 1
  initial_marker = tracer_marker
  marker = tracer_marker
  max_h_level = 1
  [Markers]
    [tracer_marker]
      type = ValueRangeMarker
      variable = tracer
      lower_bound = 0.02
      upper_bound = 0.98
    []
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '2 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.3,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 2
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = NodalValueSampler
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 0.2
  dt = 0.05
[]
[Outputs]
  csv = true
[]
(modules/porous_flow/test/tests/pressure_pulse/pressure_pulse_1d_2phasePS_KT.i)
# Pressure pulse in 1D with 2 phases, 2components - transient
# Using KT stabilization
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 10
  xmin = 0
  xmax = 100
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [ppwater]
    initial_condition = 2e6
  []
  [sgas]
    initial_condition = 0.3
  []
[]
[AuxVariables]
  [massfrac_ph0_sp0]
    initial_condition = 1
  []
  [massfrac_ph1_sp0]
    initial_condition = 0
  []
  [ppgas]
    family = MONOMIAL
    order = FIRST
  []
[]
[Kernels]
  [mass_component0]
    type = PorousFlowMassTimeDerivative
    variable = ppwater
    fluid_component = 0
  []
  [flux_component0_phase0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = ppwater
    advective_flux_calculator = afc_component0_phase0
  []
  [flux_component0_phase1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = ppwater
    advective_flux_calculator = afc_component0_phase1
  []
  [mass_component1]
    type = PorousFlowMassTimeDerivative
    variable = sgas
    fluid_component = 1
  []
  [flux_component1_phase0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = sgas
    advective_flux_calculator = afc_component1_phase0
  []
  [flux_component1_phase1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = sgas
    advective_flux_calculator = afc_component1_phase1
  []
[]
[AuxKernels]
  [ppgas]
    type = PorousFlowPropertyAux
    property = pressure
    phase = 1
    variable = ppgas
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'ppwater sgas'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
    pc = 1e5
  []
  [afc_component0_phase0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    fluid_component = 0
    phase = 0
    flux_limiter_type = superbee
  []
  [afc_component0_phase1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    fluid_component = 0
    phase = 1
    flux_limiter_type = superbee
  []
  [afc_component1_phase0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    fluid_component = 1
    phase = 0
    flux_limiter_type = superbee
  []
  [afc_component1_phase1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    fluid_component = 1
    phase = 1
    flux_limiter_type = superbee
  []
[]
[FluidProperties]
  [simple_fluid0]
    type = SimpleFluidProperties
    bulk_modulus = 2e9
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
  [simple_fluid1]
    type = SimpleFluidProperties
    bulk_modulus = 2e7
    density0 = 1
    thermal_expansion = 0
    viscosity = 1e-5
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow2PhasePS
    phase0_porepressure = ppwater
    phase1_saturation = sgas
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
  []
  [simple_fluid0]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid0
    phase = 0
  []
  [simple_fluid1]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid1
    phase = 1
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-15 0 0 0 1e-15 0 0 0 1e-15'
  []
  [relperm_water]
    type = PorousFlowRelativePermeabilityCorey
    n = 1
    phase = 0
  []
  [relperm_gas]
    type = PorousFlowRelativePermeabilityCorey
    n = 1
    phase = 1
  []
[]
[BCs]
  [leftwater]
    type = DirichletBC
    boundary = left
    value = 3e6
    variable = ppwater
  []
  [rightwater]
    type = DirichletBC
    boundary = right
    value = 2e6
    variable = ppwater
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
    petsc_options_value = 'bcgs bjacobi 1E-15 1E-20 10000'
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1e3
  end_time = 1e4
[]
[VectorPostprocessors]
  [pp]
    type = LineValueSampler
    warn_discontinuous_face_values = false
    sort_by = x
    variable = 'ppwater ppgas'
    start_point = '0 0 0'
    end_point = '100 0 0'
    num_points = 11
  []
[]
[Outputs]
  file_base = pressure_pulse_1d_2phasePS_KT
  print_linear_residuals = false
  [csv]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_05.i)
# Checking the Jacobian of Flux-Limited TVD Advection, 2 phases, 2 components, using flux_limiter_type != None
#
# Here we use snes_check_jacobian instead of snes_type=test.  The former just checks the Jacobian for the
# random initial conditions, while the latter checks for u=1 and u=-1
#
# The Jacobian is correct for u=1 and u=-1, but the finite-difference scheme used by snes_type=test gives the
# wrong answer.
# For u=constant, the Kuzmin-Turek scheme adds as much antidiffusion as possible, resulting in a central-difference
# version of advection (flux_limiter = 1).  This is correct, and the Jacobian is calculated correctly.
# However, when computing the Jacobian using finite differences, u is increased or decreased at a node.
# This results in that node being at a maximum or minimum, which means no antidiffusion should be added
# (flux_limiter = 0).  This corresponds to a full-upwind scheme.  So the finite-difference computes the
# Jacobian in the full-upwind scenario, which is incorrect (the original residual = 0, after finite-differencing
# the residual comes from the full-upwind scenario).
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 10
  xmin = 0
  xmax = 5
[]
[GlobalParams]
  gravity = '1.1 2 -0.5'
  PorousFlowDictator = dictator
[]
[Variables]
  [ppwater]
  []
  [ppgas]
  []
  [massfrac_ph0_sp0]
  []
  [massfrac_ph1_sp0]
  []
[]
[ICs]
  [ppwater]
    type = FunctionIC
    variable = ppwater
    function = 'if(x<1,0,if(x<4,sin(x-1),1))'
  []
  [ppgas]
    type = FunctionIC
    variable = ppgas
    function = 'x*(6-x)/6'
  []
  [massfrac_ph0_sp0]
    type = FunctionIC
    variable = massfrac_ph0_sp0
    function = 'x/6'
  []
  [massfrac_ph1_sp0]
    type = FunctionIC
    variable = massfrac_ph1_sp0
    function = '1-x/7'
  []
[]
[Kernels]
  [flux_ph0_sp0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = ppwater
    advective_flux_calculator = advective_flux_calculator_ph0_sp0
  []
  [flux_ph0_sp1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = ppgas
    advective_flux_calculator = advective_flux_calculator_ph0_sp1
  []
  [flux_ph1_sp0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = massfrac_ph0_sp0
    advective_flux_calculator = advective_flux_calculator_ph1_sp0
  []
  [flux_ph1_sp1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = massfrac_ph1_sp0
    advective_flux_calculator = advective_flux_calculator_ph1_sp1
  []
[]
[FluidProperties]
  [simple_fluid0]
    type = SimpleFluidProperties
    bulk_modulus = 1.5
    density0 = 1
    thermal_expansion = 0
    viscosity = 1
  []
  [simple_fluid1]
    type = SimpleFluidProperties
    bulk_modulus = 0.5
    density0 = 0.5
    thermal_expansion = 0
    viscosity = 1.4
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'ppwater ppgas massfrac_ph0_sp0 massfrac_ph1_sp0'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureVG
    alpha = 1
    m = 0.5
  []
  [advective_flux_calculator_ph0_sp0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = minmod
    phase = 0
    fluid_component = 0
  []
  [advective_flux_calculator_ph0_sp1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = vanleer
    phase = 0
    fluid_component = 1
  []
  [advective_flux_calculator_ph1_sp0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = mc
    phase = 1
    fluid_component = 0
  []
  [advective_flux_calculator_ph1_sp1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = superbee
    phase = 1
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow2PhasePP
    phase0_porepressure = ppwater
    phase1_porepressure = ppgas
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
  []
  [simple_fluid0]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid0
    phase = 0
  []
  [simple_fluid1]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid1
    phase = 1
  []
  [relperm0]
    type = PorousFlowRelativePermeabilityCorey
    n = 2
    phase = 0
  []
  [relperm1]
    type = PorousFlowRelativePermeabilityCorey
    n = 3
    phase = 1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1.21 0 0  0 1.5 0  0 0 0.8'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options = '-snes_check_jacobian'
  []
[]
[Executioner]
  type = Transient
  solve_type = Linear # this is to force convergence even though the nonlinear residual is high: we just care about the Jacobian in this test
  end_time = 1
  num_steps = 1
  dt = 1
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_3D.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, employing PorousFlow Kernels and UserObjects, with superbee flux-limiter
# 3D version
[Mesh]
  type = GeneratedMesh
  dim = 3
  nx = 10
  xmin = 0
  xmax = 1
  ny = 4
  ymin = 0
  ymax = 0.5
  nz = 3
  zmin = 0
  zmax = 2
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.3,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0.5 2'
    num_points = 11
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 0.3
  dt = 6E-2
  nl_abs_tol = 1E-8
  timestep_tolerance = 1E-3
[]
[Outputs]
  print_linear_residuals = false
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_1D_adaptivity.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, mploying PorousFlow Kernels and UserObjects, with superbee flux-limiter
# 1D version with adaptivity
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 10
  xmin = 0
  xmax = 1
[]
[Adaptivity]
  initial_steps = 1
  initial_marker = tracer_marker
  marker = tracer_marker
  max_h_level = 1
  [Markers]
    [tracer_marker]
      type = ValueRangeMarker
      variable = tracer
      lower_bound = 0.02
      upper_bound = 0.98
    []
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.3,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0 0'
    num_points = 11
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 6
  dt = 6E-2
  nl_abs_tol = 1E-8
  timestep_tolerance = 1E-3
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_2D_trimesh.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, mploying PorousFlow Kernels and UserObjects, with superbee flux-limiter
# 2D version
[Mesh]
  type = FileMesh
  file = trimesh.msh
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
  block = '50'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.305,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0.04 0'
    num_points = 101
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 6
  dt = 6E-2
  nl_abs_tol = 1E-8
  timestep_tolerance = 1E-3
[]
[Outputs]
  print_linear_residuals = false
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/heat_advection/heat_advection_1d_KT.i)
# 1phase, heat advecting with a moving fluid
# Using the Kuzmin-Turek stabilization scheme
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 50
  xmin = 0
  xmax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [temp]
    initial_condition = 200
  []
  [pp]
  []
[]
[ICs]
  [pp]
    type = FunctionIC
    variable = pp
    function = '1-x'
  []
[]
[BCs]
  [pp0]
    type = DirichletBC
    variable = pp
    boundary = left
    value = 1
  []
  [pp1]
    type = DirichletBC
    variable = pp
    boundary = right
    value = 0
  []
  [spit_heat]
    type = DirichletBC
    variable = temp
    boundary = left
    value = 300
  []
  [suck_heat]
    type = DirichletBC
    variable = temp
    boundary = right
    value = 200
  []
[]
[Kernels]
  [mass_dot]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [fluid_advection]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = pp
    advective_flux_calculator = fluid_advective_flux
  []
  [energy_dot]
    type = PorousFlowEnergyTimeDerivative
    variable = temp
  []
  [heat_advection]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = temp
    advective_flux_calculator = heat_advective_flux
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'temp pp'
    number_fluid_phases = 1
    number_fluid_components = 1
  []
  [pc]
    type = PorousFlowCapillaryPressureVG
    m = 0.6
    alpha = 1.3
  []
  [fluid_advective_flux]
    type = PorousFlowAdvectiveFluxCalculatorSaturated
    flux_limiter_type = superbee
  []
  [heat_advective_flux]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedHeat
    flux_limiter_type = superbee
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 100
    density0 = 1000
    viscosity = 4.4
    thermal_expansion = 0
    cv = 2
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
    temperature = temp
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.2
  []
  [rock_heat]
    type = PorousFlowMatrixInternalEnergy
    specific_heat_capacity = 1.0
    density = 125
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1.1 0 0 0 2 0 0 0 3'
  []
  [relperm]
    type = PorousFlowRelativePermeabilityCorey
    n = 2
    phase = 0
  []
  [massfrac]
    type = PorousFlowMassFraction
  []
  [PS]
    type = PorousFlow1PhaseP
    porepressure = pp
    capillary_pressure = pc
  []
[]
[Preconditioning]
  [andy]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
    petsc_options_value = 'gmres bjacobi 1E-15 1E-10 10000'
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 0.01
  end_time = 0.6
[]
[VectorPostprocessors]
  [T]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0 0'
    num_points = 51
    sort_by = x
    variable = temp
  []
[]
[Outputs]
  file_base = heat_advection_1d_KT
  [csv]
    type = CSV
    sync_times = '0.1 0.6'
    sync_only = true
  []
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/pffltvd_1D.i)
# Using flux-limited TVD advection ala Kuzmin and Turek, mploying PorousFlow Kernels and UserObjects, with superbee flux-limiter
# 1D version
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 10
  xmin = 0
  xmax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
  []
  [tracer]
  []
[]
[ICs]
  [porepressure]
    type = FunctionIC
    variable = porepressure
    function = '1 - x'
  []
  [tracer]
    type = FunctionIC
    variable = tracer
    function = 'if(x<0.1,0,if(x>0.3,0,1))'
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = tracer
  []
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = tracer
    advective_flux_calculator = advective_flux_calculator_0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [flux1]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = porepressure
    advective_flux_calculator = advective_flux_calculator_1
  []
[]
[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1
    boundary = left
  []
  [no_tracer_on_left]
    type = DirichletBC
    variable = tracer
    value = 0
    boundary = left
  []
  [remove_component_1]
    type = PorousFlowPiecewiseLinearSink
    variable = porepressure
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 1
    use_mobility = true
    flux_function = 1E3
  []
  [remove_component_0]
    type = PorousFlowPiecewiseLinearSink
    variable = tracer
    boundary = right
    fluid_phase = 0
    pt_vals = '0 1E3'
    multipliers = '0 1E3'
    mass_fraction_component = 0
    use_mobility = true
    flux_function = 1E3
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    thermal_expansion = 0
    viscosity = 1.0
    density0 = 1000.0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
  [advective_flux_calculator_0]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 0
  []
  [advective_flux_calculator_1]
    type = PorousFlowAdvectiveFluxCalculatorSaturatedMultiComponent
    flux_limiter_type = superbee
    fluid_component = 1
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = porepressure
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = tracer
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = the_simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-2 0 0   0 1E-2 0   0 0 1E-2'
  []
[]
[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]
[VectorPostprocessors]
  [tracer]
    type = LineValueSampler
    start_point = '0 0 0'
    end_point = '1 0 0'
    num_points = 11
    sort_by = x
    variable = tracer
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 6
  dt = 6E-2
  nl_abs_tol = 1E-8
  timestep_tolerance = 1E-3
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/flux_limited_TVD_pflow/jacobian_04.i)
# Checking the Jacobian of Flux-Limited TVD Advection, 1 phase, 1 component, unsaturated, using flux_limiter_type != none
# This is quite a heavy test, but we need a fairly big mesh to check the flux-limiting+TVD is happening correctly
#
# Here we use snes_check_jacobian instead of snes_type=test.  The former just checks the Jacobian for the
# random initial conditions, while the latter checks for u=1 and u=-1
#
# The Jacobian is correct for u=1 and u=-1, but the finite-difference scheme used by snes_type=test gives the
# wrong answer.
# For u=constant, the Kuzmin-Turek scheme adds as much antidiffusion as possible, resulting in a central-difference
# version of advection (flux_limiter = 1).  This is correct, and the Jacobian is calculated correctly.
# However, when computing the Jacobian using finite differences, u is increased or decreased at a node.
# This results in that node being at a maximum or minimum, which means no antidiffusion should be added
# (flux_limiter = 0).  This corresponds to a full-upwind scheme.  So the finite-difference computes the
# Jacobian in the full-upwind scenario, which is incorrect (the original residual = 0, after finite-differencing
# the residual comes from the full-upwind scenario).
[Mesh]
  type = GeneratedMesh
  dim = 3
  nx = 3
  xmin = 0
  xmax = 1
  ny = 4
  ymin = -1
  ymax = 2
  bias_y = 1.5
  nz = 4
  zmin = 1
  zmax = 2
  bias_z = 0.8
[]
[GlobalParams]
  gravity = '1 2 -0.5'
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
[]
[ICs]
  [pp]
    variable = pp
    type = RandomIC
    min = -1
    max = 0
  []
[]
[Kernels]
  [flux0]
    type = PorousFlowFluxLimitedTVDAdvection
    variable = pp
    advective_flux_calculator = advective_flux_calculator
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1
    density0 = 0.4
    viscosity = 1.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp'
    number_fluid_phases = 1
    number_fluid_components = 1
  []
  [pc]
    type = PorousFlowCapillaryPressureVG
    alpha = 1
    m = 0.5
  []
  [advective_flux_calculator]
    type = PorousFlowAdvectiveFluxCalculatorUnsaturated
    flux_limiter_type = minmod
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseP
    porepressure = pp
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [relperm]
    type = PorousFlowRelativePermeabilityCorey
    phase = 0
    n = 2
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1.21 0 0  0 1.5 0  0 0 0.8'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options = '-snes_check_jacobian'
  []
[]
[Executioner]
  type = Transient
  solve_type = Linear # this is to force convergence even though the nonlinear residual is high: we just care about the Jacobian in this test
  end_time = 1
  num_steps = 1
  dt = 1
[]