- PorousFlowDictatorThe UserObject that holds the list of PorousFlow variable namesC++ Type:UserObjectName Controllable:No Description:The UserObject that holds the list of PorousFlow variable names 
- diffusion_coeffList of diffusion coefficients. Order is i) component 0 in phase 0; ii) component 1 in phase 0 ...; component 0 in phase 1; ... component k in phase n (m^2/sC++ Type:std::vector<double> Unit:(no unit assumed) Controllable:No Description:List of diffusion coefficients. Order is i) component 0 in phase 0; ii) component 1 in phase 0 ...; component 0 in phase 1; ... component k in phase n (m^2/s 
- tortuosityList of tortuosities. Order is i) phase 0; ii) phase 1; etcC++ Type:std::vector<double> Unit:(no unit assumed) Controllable:No Description:List of tortuosities. Order is i) phase 0; ii) phase 1; etc 
PorousFlowDiffusivityConst
This Material provides constant tortuosity and diffusion coefficients
Tortuosity is defined as the ratio of the shortest path to the effective path, such that .
Input Parameters
- at_nodesFalseEvaluate Material properties at nodes instead of quadpointsDefault:False C++ Type:bool Controllable:No Description:Evaluate Material properties at nodes instead of quadpoints 
- blockThe list of blocks (ids or names) that this object will be appliedC++ Type:std::vector<SubdomainName> Controllable:No Description:The list of blocks (ids or names) that this object will be applied 
- boundaryThe list of boundaries (ids or names) from the mesh where this object appliesC++ Type:std::vector<BoundaryName> Controllable:No Description:The list of boundaries (ids or names) from the mesh where this object applies 
- computeTrueWhen false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.Default:True C++ Type:bool Controllable:No Description:When false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies. 
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skippedDefault:NONE C++ Type:MooseEnum Options:NONE, ELEMENT, SUBDOMAIN Controllable:No Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped 
- declare_suffixAn optional suffix parameter that can be appended to any declared 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 declared properties. The suffix will be prepended with a '_' character. 
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.C++ Type:std::vector<std::string> Controllable:No Description:Adds user-defined labels for accessing object parameters via control logic. 
- enableTrueSet the enabled status of the MooseObject.Default:True C++ Type:bool Controllable:Yes Description:Set the enabled status of the MooseObject. 
- 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 
- search_methodnearest_node_connected_sidesChoice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes).Default:nearest_node_connected_sides C++ Type:MooseEnum Options:nearest_node_connected_sides, all_proximate_sides Controllable:No Description:Choice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes). 
- seed0The seed for the master random number generatorDefault:0 C++ Type:unsigned int Controllable:No Description:The seed for the master random number generator 
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.Default:False C++ Type:bool Controllable:No Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used. 
Advanced Parameters
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)C++ Type:std::vector<std::string> Controllable:No Description:List of material properties, from this material, to output (outputs must also be defined to an output type) 
- outputsnone Vector of output names where you would like to restrict the output of variables(s) associated with this objectDefault:none C++ Type:std::vector<OutputName> Controllable:No Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object 
Outputs 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/chemistry/2species_predis.i)
- (modules/porous_flow/test/tests/chemistry/2species_equilibrium.i)
- (modules/porous_flow/test/tests/jacobian/diff02.i)
- (modules/porous_flow/test/tests/dispersion/disp01_fv.i)
- (modules/porous_flow/examples/lava_lamp/2phase_convection.i)
- (modules/porous_flow/examples/flow_through_fractured_media/coarse_3D.i)
- (modules/porous_flow/examples/fluidflower/fluidflower.i)
- (modules/porous_flow/examples/solute_tracer_transport/solute_tracer_transport.i)
- (modules/porous_flow/test/tests/jacobian/disp04.i)
- (modules/porous_flow/examples/flow_through_fractured_media/fine_thick_fracture_transient.i)
- (modules/porous_flow/test/tests/dispersion/diff01_fv.i)
- (modules/porous_flow/examples/flow_through_fractured_media/coarse.i)
- (modules/porous_flow/test/tests/jacobian/disp02.i)
- (modules/porous_flow/test/tests/dispersion/diff01.i)
- (modules/porous_flow/test/tests/jacobian/diff01.i)
- (modules/porous_flow/test/tests/chemistry/2species_equilibrium_2phase.i)
- (modules/porous_flow/examples/lava_lamp/1phase_convection.i)
- (modules/porous_flow/test/tests/dispersion/disp01.i)
- (modules/porous_flow/test/tests/jacobian/disp01.i)
- (modules/porous_flow/test/tests/jacobian/disp03.i)
- (modules/porous_flow/test/tests/dispersion/diff01_action.i)
- (modules/porous_flow/examples/flow_through_fractured_media/fine_transient.i)
- (modules/porous_flow/test/tests/dispersion/disp01_heavy.i)
- (modules/porous_flow/examples/solute_tracer_transport/solute_tracer_transport_2D.i)
(modules/porous_flow/test/tests/chemistry/2species_predis.i)
# PorousFlow analogy of chemical_reactions/test/tests/solid_kinetics/2species_without_action.i
#
# Simple equilibrium reaction example to illustrate the use of PorousFlowAqueousPreDisChemistry
#
# In this example, two primary species a and b diffuse towards each other from
# opposite ends of a porous medium, reacting when they meet to form a mineral
# precipitate. The kinetic reaction is
#
# a + b = mineral
#
# where a and b are the primary species (reactants), and mineral is the precipitate.
# At the time of writing, the results of this test differ from chemical_reactions because
# in PorousFlow the mineral_concentration is measured in m^3 (precipitate) / m^3 (porous_material)
# in chemical_reactions the mineral_concentration is measured in m^3 (precipitate) / m^3 (fluid)
# ie, PorousFlow_mineral_concentration = porosity * chemical_reactions_mineral_concentration
[Mesh]
  type = GeneratedMesh
  dim = 2
  xmax = 1
  ymax = 1
  nx = 40
[]
[Variables]
  [a]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0
  []
  [b]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0
  []
[]
[AuxVariables]
  [eqm_k]
    initial_condition = 1E-6
  []
  [pressure]
  []
  [mineral]
    family = MONOMIAL
    order = CONSTANT
  []
[]
[AuxKernels]
  [mineral]
    type = PorousFlowPropertyAux
    property = mineral_concentration
    mineral_species = 0
    variable = mineral
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Kernels]
  [mass_a]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = a
  []
  [diff_a]
    type = PorousFlowDispersiveFlux
    variable = a
    fluid_component = 0
    disp_trans = 0
    disp_long = 0
  []
  [predis_a]
    type = PorousFlowPreDis
    variable = a
    mineral_density = 1000
    stoichiometry = 1
  []
  [mass_b]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = b
  []
  [diff_b]
    type = PorousFlowDispersiveFlux
    variable = b
    fluid_component = 1
    disp_trans = 0
    disp_long = 0
  []
  [predis_b]
    type = PorousFlowPreDis
    variable = b
    mineral_density = 1000
    stoichiometry = 1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'a b'
    number_fluid_phases = 1
    number_fluid_components = 3
    number_aqueous_kinetic = 1
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9 # huge, so mimic chemical_reactions
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
    temperature = 298.15
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pressure
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'a b'
  []
  [chem]
    type = PorousFlowAqueousPreDisChemistry
    primary_concentrations = 'a b'
    num_reactions = 1
    equilibrium_constants = eqm_k
    primary_activity_coefficients = '1 1'
    reactions = '1 1'
    specific_reactive_surface_area = '1.0'
    kinetic_rate_constant = '1.0e-8'
    activation_energy = '1.5e4'
    molar_volume = 1
    gas_constant = 8.314
    reference_temperature = 298.15
  []
  [mineral_conc]
    type = PorousFlowAqueousPreDisMineral
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.4
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    # porous_flow permeability / porous_flow viscosity = chemical_reactions conductivity = 4E-3
    permeability = '4E-6 0 0 0 4E-6 0 0 0 4E-6'
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [diff]
    type = PorousFlowDiffusivityConst
    # porous_flow diffusion_coeff * tortuousity * porosity = chemical_reactions diffusivity = 5E-4
    diffusion_coeff = '12.5E-4 12.5E-4 12.5E-4'
    tortuosity = 1.0
  []
[]
[BCs]
  [a_left]
    type = DirichletBC
    variable = a
    boundary = left
    value = 1.0e-2
  []
  [a_right]
    type = DirichletBC
    variable = a
    boundary = right
    value = 0
  []
  [b_left]
    type = DirichletBC
    variable = b
    boundary = left
    value = 0
  []
  [b_right]
    type = DirichletBC
    variable = b
    boundary = right
    value = 1.0e-2
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 5
  end_time = 50
[]
[Outputs]
  print_linear_residuals = true
  exodus = true
  perf_graph = true
  hide = eqm_k
[]
(modules/porous_flow/test/tests/chemistry/2species_equilibrium.i)
# PorousFlow analogy of chemical_reactions/test/tests/aqueous_equilibrium/2species.i
#
# Simple equilibrium reaction example to illustrate the use of PorousFlowMassFractionAqueousEquilibriumChemistry
#
# In this example, two primary species a and b are transported by diffusion and convection
# from the left of the porous medium, reacting to form two equilibrium species pa2 and pab
# according to the equilibrium reaction:
#
#      reactions = '2a = pa2     rate = 10^2
#                   a + b = pab  rate = 10^-2'
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 10
[]
[Variables]
  [a]
    order = FIRST
    family = LAGRANGE
    [InitialCondition]
      type = BoundingBoxIC
      x1 = 0.0
      y1 = 0.0
      x2 = 1.0e-10
      y2 = 1
      inside = 1.0e-2
      outside = 1.0e-10
    []
  []
  [b]
    order = FIRST
    family = LAGRANGE
    [InitialCondition]
      type = BoundingBoxIC
      x1 = 0.0
      y1 = 0.0
      x2 = 1.0e-10
      y2 = 1
      inside = 1.0e-2
      outside = 1.0e-10
    []
  []
[]
[AuxVariables]
  [eqm_k0]
    initial_condition = 1E2
  []
  [eqm_k1]
    initial_condition = 1E-2
  []
  [pressure]
  []
  [pa2]
    family = MONOMIAL
    order = CONSTANT
  []
  [pab]
    family = MONOMIAL
    order = CONSTANT
  []
[]
[AuxKernels]
  [pa2]
    type = PorousFlowPropertyAux
    property = secondary_concentration
    secondary_species = 0
    variable = pa2
  []
  [pab]
    type = PorousFlowPropertyAux
    property = secondary_concentration
    secondary_species = 1
    variable = pab
  []
[]
[ICs]
  [pressure]
    type = FunctionIC
    variable = pressure
    function = 2-x
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Kernels]
  [mass_a]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = a
  []
  [flux_a]
    type = PorousFlowFullySaturatedDarcyFlow
    variable = a
    fluid_component = 0
  []
  [diff_a]
    type = PorousFlowDispersiveFlux
    variable = a
    fluid_component = 0
    disp_trans = 0
    disp_long = 0
  []
  [mass_b]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = b
  []
  [flux_b]
    type = PorousFlowFullySaturatedDarcyFlow
    variable = b
    fluid_component = 1
  []
  [diff_b]
    type = PorousFlowDispersiveFlux
    variable = b
    fluid_component = 1
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'a b'
    number_fluid_phases = 1
    number_fluid_components = 3
    number_aqueous_equilibrium = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9 # huge, so mimic chemical_reactions
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pressure
  []
  [massfrac]
    type = PorousFlowMassFractionAqueousEquilibriumChemistry
    mass_fraction_vars = 'a b'
    num_reactions = 2
    equilibrium_constants = 'eqm_k0 eqm_k1'
    primary_activity_coefficients = '1 1'
    secondary_activity_coefficients = '1 1'
    reactions = '2 0
                 1 1'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.2
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    # porous_flow permeability / porous_flow viscosity = chemical_reactions conductivity = 1E-4
    permeability = '1E-7 0 0 0 1E-7 0 0 0 1E-7'
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [diff]
    type = PorousFlowDiffusivityConst
    # porous_flow diffusion_coeff * tortuousity * porosity = chemical_reactions diffusivity = 1E-4
    diffusion_coeff = '5E-4 5E-4 5E-4'
    tortuosity = 1.0
  []
[]
[BCs]
  [a_left]
    type = DirichletBC
    variable = a
    boundary = left
    value = 1.0e-2
  []
  [b_left]
    type = DirichletBC
    variable = b
    boundary = left
    value = 1.0e-2
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 10
  end_time = 100
[]
[Outputs]
  print_linear_residuals = true
  exodus = true
  perf_graph = true
  hide = eqm_k0
[]
(modules/porous_flow/test/tests/jacobian/diff02.i)
# Test the Jacobian of the diffusive component of the PorousFlowDisperiveFlux kernel for two phases.
# By setting disp_long and disp_trans to zero, the purely diffusive component of the flux
# can be isolated. Uses constant tortuosity and diffusion coefficients
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
[]
[Variables]
  [sgas]
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [massfrac1]
  []
[]
[ICs]
  [sgas]
    type = RandomIC
    variable = sgas
    max = 1
    min = 0
  []
  [massfrac0]
    type = RandomIC
    variable = massfrac0
    min = 0
    max = 1
  []
  [massfrac1]
    type = RandomIC
    variable = massfrac1
    min = 0
    max = 1
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = sgas
    gravity = '1 0 0'
    disp_long = '0 0'
    disp_trans = '0 0'
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    gravity = '1 0 0'
    disp_long = '0 0'
    disp_trans = '0 0'
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'sgas massfrac0'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
    pc = 0
  []
[]
[FluidProperties]
  [simple_fluid0]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 10
    thermal_expansion = 0
    viscosity = 1
  []
  [simple_fluid1]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 1
    thermal_expansion = 0
    viscosity = 0.1
  []
[]
[Materials]
  [temp]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow2PhasePS
    phase0_porepressure = 1
    phase1_saturation = sgas
    capillary_pressure = pc
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac0 massfrac1'
  []
  [simple_fluid0]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid0
    phase = 0
  []
  [simple_fluid1]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid1
    phase = 1
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [diff]
    type = PorousFlowDiffusivityConst
     diffusion_coeff = '1e-2 1e-1 1e-2 1e-1'
     tortuosity = '0.1 0.2'
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1 0 0 0 2 0 0 0 3'
  []
  [relperm0]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [relperm1]
    type = PorousFlowRelativePermeabilityConst
    phase = 1
  []
[]
[Preconditioning]
  active = smp
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1
  end_time = 1
[]
[Outputs]
  exodus = false
[]
(modules/porous_flow/test/tests/dispersion/disp01_fv.i)
# Test dispersive part of FVPorousFlowDispersiveFlux kernel by setting diffusion
# coefficients to zero. A pressure gradient is applied over the mesh to give a
# uniform velocity. Gravity is set to zero.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 20
  xmax = 10
  bias_x = 1.1
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
    type = MooseVariableFVReal
  []
  [massfrac0]
    type = MooseVariableFVReal
  []
[]
[AuxVariables]
  [velocity]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [velocity]
    type = ADPorousFlowDarcyVelocityComponent
    variable = velocity
    component = x
  []
[]
[ICs]
  [pp]
    type = FunctionIC
    variable = pp
    function = pic
  []
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[Functions]
  [pic]
    type = ParsedFunction
    expression = '1.1e5-x*1e3'
  []
[]
[FVBCs]
  [xleft]
    type = FVDirichletBC
    value = 1
    variable = massfrac0
    boundary = left
  []
  [xright]
    type = FVDirichletBC
    value = 0
    variable = massfrac0
    boundary = right
  []
  [pright]
    type = FVDirichletBC
    variable = pp
    boundary = right
    value = 1e5
  []
  [pleft]
    type = FVDirichletBC
    variable = pp
    boundary = left
    value = 1.1e5
  []
[]
[FVKernels]
  [mass0]
    type = FVPorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = FVPorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = FVPorousFlowDispersiveFlux
    variable = pp
    disp_trans = 0
    disp_long = 0.2
  []
  [mass1]
    type = FVPorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [adv1]
    type = FVPorousFlowAdvectiveFlux
    fluid_component = 1
    variable = massfrac0
  []
  [diff1]
    type = FVPorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0.2
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e9
    density0 = 1000
    viscosity = 0.001
    thermal_expansion = 0
  []
[]
[Materials]
  [temperature]
    type = ADPorousFlowTemperature
  []
  [ppss]
    type = ADPorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = ADPorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = ADPorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = ADPorousFlowPorosityConst
    porosity = 0.3
  []
  [diff]
    type = ADPorousFlowDiffusivityConst
    diffusion_coeff = '0 0'
    tortuosity = 0.1
  []
  [relp]
    type = ADPorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability]
    type = ADPorousFlowPermeabilityConst
    permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type'
    petsc_options_value = 'gmres      asm      lu           NONZERO'
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 3e2
  dtmax = 100
  nl_abs_tol = 1e-12
  [TimeStepper]
    type = IterationAdaptiveDT
    growth_factor = 2
    cutback_factor = 0.5
    dt = 10
  []
[]
[VectorPostprocessors]
  [xmass]
    type = ElementValueSampler
    sort_by = id
    variable = 'massfrac0 velocity'
  []
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/examples/lava_lamp/2phase_convection.i)
# Two phase density-driven convection of dissolved CO2 in brine
#
# Initially, the model has a gas phase at the top with a saturation of 0.29
# (which corresponds to an initial value of zi = 0.2).
# Diffusion of the dissolved CO2
# component from the saturated liquid to the unsaturated liquid below reduces the
# amount of CO2 in the gas phase. As the density of the CO2-saturated brine is greater
# than the unsaturated brine, a gravitational instability arises and density-driven
# convection of CO2-rich fingers descend into the unsaturated brine.
#
# The instability is seeded by a random perturbation to the porosity field.
# Mesh adaptivity is used to refine the mesh as the fingers form.
#
# Note: this model is computationally expensive, so should be run with multiple cores,
# preferably on a cluster.
[GlobalParams]
  PorousFlowDictator = 'dictator'
  gravity = '0 -9.81 0'
[]
[Adaptivity]
  max_h_level = 2
  marker = marker
  initial_marker = initial
  initial_steps = 2
  [Indicators]
    [indicator]
      type = GradientJumpIndicator
      variable = zi
    []
  []
  [Markers]
    [marker]
      type = ErrorFractionMarker
      indicator = indicator
      refine = 0.8
    []
    [initial]
      type = BoxMarker
      bottom_left = '0 1.95 0'
      top_right = '2 2 0'
      inside = REFINE
      outside = DO_NOTHING
    []
  []
[]
[Mesh]
  type = GeneratedMesh
  dim = 2
  ymax = 2
  xmax = 2
  ny = 40
  nx = 40
  bias_y = 0.95
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pgas
  []
  [flux0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pgas
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pgas
    disp_long = '0 0'
    disp_trans = '0 0'
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = zi
  []
  [flux1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = zi
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = zi
    disp_long = '0 0'
    disp_trans = '0 0'
  []
[]
[AuxVariables]
  [xnacl]
    initial_condition = 0.01
  []
  [saturation_gas]
    order = FIRST
    family = MONOMIAL
  []
  [xco2l]
    order = FIRST
    family = MONOMIAL
  []
  [density_liquid]
    order = FIRST
    family = MONOMIAL
  []
  [porosity]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  [saturation_gas]
    type = PorousFlowPropertyAux
    variable = saturation_gas
    property = saturation
    phase = 1
    execute_on = 'timestep_end'
  []
  [xco2l]
    type = PorousFlowPropertyAux
    variable = xco2l
    property = mass_fraction
    phase = 0
    fluid_component = 1
    execute_on = 'timestep_end'
  []
  [density_liquid]
    type = PorousFlowPropertyAux
    variable = density_liquid
    property = density
    phase = 0
    execute_on = 'timestep_end'
  []
[]
[Variables]
  [pgas]
  []
  [zi]
    scaling = 1e4
  []
[]
[ICs]
  [pressure]
    type = FunctionIC
    function = 10e6-9.81*1000*y
    variable = pgas
  []
  [zi]
    type = BoundingBoxIC
    variable = zi
    x1 = 0
    x2 = 2
    y1 = 1.95
    y2 = 2
    inside = 0.2
    outside = 0
  []
  [porosity]
    type = RandomIC
    variable = porosity
    min = 0.25
    max = 0.275
    seed = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pgas zi'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
    pc = 0
  []
  [fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2
    capillary_pressure = pc
  []
[]
[FluidProperties]
  [co2sw]
    type = CO2FluidProperties
  []
  [co2]
    type = TabulatedBicubicFluidProperties
    fp = co2sw
  []
  [brine]
    type = BrineFluidProperties
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
    temperature = '45'
  []
  [brineco2]
    type = PorousFlowFluidState
    gas_porepressure = 'pgas'
    z = 'zi'
    temperature_unit = Celsius
    xnacl = 'xnacl'
    capillary_pressure = pc
    fluid_state = fs
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = porosity
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-11 0 0 0 1e-11 0 0 0 1e-11'
  []
  [relperm_water]
    type = PorousFlowRelativePermeabilityCorey
    phase = 0
    n = 2
    s_res = 0.1
    sum_s_res = 0.2
  []
  [relperm_gas]
    type = PorousFlowRelativePermeabilityCorey
    phase = 1
    n = 2
    s_res = 0.1
    sum_s_res = 0.2
  []
  [diffusivity]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '2e-9 2e-9 2e-9 2e-9'
    tortuosity = '1 1'
  []
[]
[Preconditioning]
  active = basic
  [mumps_is_best_for_parallel_jobs]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
  [basic]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 1e6
  nl_max_its = 25
  l_max_its = 100
  dtmax = 1e4
  nl_abs_tol = 1e-6
  [TimeStepper]
    type = IterationAdaptiveDT
    dt = 10
    growth_factor = 2
    cutback_factor = 0.5
  []
[]
[Functions]
  [flux]
    type = ParsedFunction
    symbol_values = 'delta_xco2 dt'
    symbol_names = 'dx dt'
    expression = 'dx/dt'
  []
[]
[Postprocessors]
  [total_co2_in_gas]
    type = PorousFlowFluidMass
    phase = 1
    fluid_component = 1
  []
  [total_co2_in_liquid]
    type = PorousFlowFluidMass
    phase = 0
    fluid_component = 1
  []
  [numdofs]
    type = NumDOFs
  []
  [delta_xco2]
    type = ChangeOverTimePostprocessor
    postprocessor = total_co2_in_liquid
  []
  [dt]
    type = TimestepSize
  []
  [flux]
    type = FunctionValuePostprocessor
    function = flux
  []
[]
[Outputs]
  print_linear_residuals = false
  perf_graph = true
  exodus = true
  csv = true
[]
(modules/porous_flow/examples/flow_through_fractured_media/coarse_3D.i)
# Flow and solute transport along 2 2D eliptical fractures embedded in a 3D porous matrix
# the model domain has dimensions 1 x 1 x 0.3m and the two fracture have r1 = 0.45 and r2 = 0.2
# The fractures intersect each other and the domain boundaries on two opposite sides
# fracture aperture = 6e-4m
# fracture porosity = 6e-4m
# fracture permeability = 1.8e-11 which is based in k=3e-8 from a**2/12, and k*a = 3e-8*6e-4;
# matrix porosity = 0.1;
# matrix permeanility = 1e-20;
[Mesh]
  type = FileMesh
  file = coarse_3D.e
  block_id = '1 2 3'
  block_name = 'matrix f1 f2'
  boundary_id = '1 2 3 4'
  boundary_name = 'rf2 lf1 right_matrix left_matrix'
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
  []
  [tracer]
  []
[]
[AuxVariables]
  [velocity_x]
    family = MONOMIAL
    order = CONSTANT
    block = 'f1 f2'
  []
  [velocity_y]
    family = MONOMIAL
    order = CONSTANT
    block = 'f1 f2'
  []
  [velocity_z]
    family = MONOMIAL
    order = CONSTANT
    block = 'f1 f2'
  []
[]
[AuxKernels]
  [velocity_x]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_x
    component = x
    aperture = 6E-4
  []
  [velocity_y]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_y
    component = y
    aperture = 6E-4
  []
  [velocity_z]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_z
    component = z
    aperture = 6E-4
  []
[]
[ICs]
  [pp]
    type = ConstantIC
    variable = pp
    value = 1e6
  []
  [tracer]
    type = ConstantIC
    variable = tracer
    value = 0
  []
[]
[BCs]
  [top]
    type = DirichletBC
    value = 0
    variable = tracer
    boundary = rf2
  []
  [bottom]
    type = DirichletBC
    value = 1
    variable = tracer
    boundary = lf1
  []
  [ptop]
    type = DirichletBC
    variable = pp
    boundary =  rf2
    value = 1e6
  []
  [pbottom]
    type = DirichletBC
    variable = pp
    boundary = lf1
    value = 1.02e6
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    disp_trans = 0
    disp_long = 0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = tracer
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = tracer
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = tracer
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp tracer'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'tracer'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro1]
    type = PorousFlowPorosityConst
    porosity = 6e-4   # = a * phif
    block = 'f1 f2'
  []
  [diff1]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1.e-9 1.e-9'
    tortuosity = 1.0
    block = 'f1 f2'
  []
  [poro2]
    type = PorousFlowPorosityConst
    porosity = 0.1
    block = 'matrix'
  []
  [diff2]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1.e-9 1.e-9'
    tortuosity = 0.1
    block = 'matrix'
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability1]
    type = PorousFlowPermeabilityConst
    permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11'   # 1.8e-11 = a * kf
    block = 'f1 f2'
  []
  [permeability2]
    type = PorousFlowPermeabilityConst
    permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
    block = 'matrix'
  []
[]
[Preconditioning]
  active = basic
  [mumps_is_best_for_parallel_jobs]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
  [basic]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 20
  dt = 1
[]
[VectorPostprocessors]
  [xmass]
    type = LineValueSampler
    start_point = '-0.5 0 0'
    end_point = '0.5 0 0'
    sort_by = x
    num_points = 41
    variable = tracer
    outputs = csv
  []
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = 'final'
  []
[]
(modules/porous_flow/examples/fluidflower/fluidflower.i)
# FluidFlower International Benchmark study model
# CSIRO 2023
#
# This example can be used to reproduce the results presented by the
# CSIRO team as part of this benchmark study. See
# Green, C., Jackson, S.J., Gunning, J., Wilkins, A. and Ennis-King, J.,
# 2023. Modelling the FluidFlower: Insights from Characterisation and
# Numerical Predictions. Transport in Porous Media.
#
# This example takes a long time to run! The large density contrast
# between the gas phase CO2 and the water makes convergence very hard,
# so small timesteps must be taken during injection.
#
# This example uses a simplified mesh in order to be run during the
# automated testing. To reproduce the results of the benchmark study,
# replace the simple layered input mesh with the one located in the
# large_media submodule.
#
# The mesh file contains:
# - porosity as given by FluidFlower description
# - permeability as given by FluidFlower description
# - subdomain ids for each sand type
#
# The nominal thickness of the FluidFlower tank is 19mm. To keep masses consistent
# with the experiment, porosity and permeability are multiplied by the thickness
thickness = 0.019
#
# Properties associated with each sand type associated with mesh block ids
#
# block 0 - ESF (very fine sand)
sandESF = '0 10 20'
sandESF_pe = 1471.5
sandESF_krg = 0.09
sandESF_swi = 0.32
sandESF_krw = 0.71
sandESF_sgi = 0.14
# block 1 - C - Coarse lower
sandC = '1 21'
sandC_pe = 294.3
sandC_krg = 0.05
sandC_swi = 0.14
sandC_krw = 0.93
sandC_sgi = 0.1
# block 2 - D - Coarse upper
sandD = '2 22'
sandD_pe = 98.1
sandD_krg = 0.02
sandD_swi = 0.12
sandD_krw = 0.95
sandD_sgi = 0.08
# block 3 - E - Very Coarse lower
sandE = '3 13 23'
sandE_pe = 10
sandE_krg = 0.1
sandE_swi = 0.12
sandE_krw = 0.93
sandE_sgi = 0.06
# block 4 - F - Very Coarse upper
sandF = '4 14 24 34'
sandF_pe = 10
sandF_krg = 0.11
sandF_swi = 0.12
sandF_krw = 0.72
sandF_sgi = 0.13
# block 5 - G - Flush Zone
sandG = '5 15 35'
sandG_pe = 10
sandG_krg = 0.16
sandG_swi = 0.1
sandG_krw = 0.75
sandG_sgi = 0.06
# block 6 - Fault 1 - Heterogeneous
fault1 = '6 26'
fault1_pe = 10
fault1_krg = 0.16
fault1_swi = 0.1
fault1_krw = 0.75
fault1_sgi = 0.06
# block 7 - Fault 2 - Impermeable
# Note: this fault has been removed from the mesh (no elements in this region)
# block 8 - Fault 3 - Homogeneous
fault3 = '8'
fault3_pe = 10
fault3_krg = 0.16
fault3_swi = 0.1
fault3_krw = 0.75
fault3_sgi = 0.06
# Top layer
top_layer = '9'
# Boxes A, B an C used to report values (sg, sgr, xco2, etc)
boxA = '10 13 14 15 34 35'
boxB = '20 21 22 23 24 26'
boxC = '34 35'
# Furthermore, the seal sand unit in boxes A and B
seal_boxA = '10'
seal_boxB = '20'
# CO2 injection details:
# CO2 density ~1.8389 kg/m3 at 293.15 K, 1.01325e5 Pa
# Injection in Port (9, 3) for 5 hours.
# Injection in Port (17, 7) for 2:45 hours.
# Injection of 10 ml/min = 0.1666 ml/s = 1.666e-7 m3/s = ~3.06e-7 kg/s.
# Total mass of CO2 injected ~ 8.5g.
inj_rate = 3.06e-7
[Mesh]
  [mesh]
    type = FileMeshGenerator
    file = 'fluidflower_test.e'
    # file = '../../../../large_media/porous_flow/examples/fluidflower/fluidflower.e'
    use_for_exodus_restart = true
  []
[]
[Debug]
  show_var_residual_norms = true
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 -9.81 0'
  temperature = temperature
  log_extension = false
[]
[Variables]
  [pgas]
    family = MONOMIAL
    order = CONSTANT
    fv = true
  []
  [z]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    scaling = 1e4
  []
[]
[AuxVariables]
  [xnacl]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_condition = 0.0055
  []
  [temperature]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_condition = 20
  []
  [porosity]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_from_file_var = porosity
  []
  [porosity_times_thickness]
    family = MONOMIAL
    order = CONSTANT
    fv = true
  []
  [permeability]
    family = MONOMIAL
    order = CONSTANT
    fv = true
    initial_from_file_var = permeability
  []
  [permeability_times_thickness]
    family = MONOMIAL
    order = CONSTANT
    fv = true
  []
  [saturation_water]
    family = MONOMIAL
    order = CONSTANT
  []
  [saturation_gas]
    family = MONOMIAL
    order = CONSTANT
  []
  [pressure_water]
    family = MONOMIAL
    order = CONSTANT
  []
  [pc]
    family = MONOMIAL
    order = CONSTANT
  []
  [x0_water]
    order = CONSTANT
    family = MONOMIAL
  []
  [x0_gas]
    order = CONSTANT
    family = MONOMIAL
  []
  [x1_water]
    order = CONSTANT
    family = MONOMIAL
  []
  [x1_gas]
    order = CONSTANT
    family = MONOMIAL
  []
  [density_water]
    order = CONSTANT
    family = MONOMIAL
  []
  [density_gas]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  [porosity_times_thickness]
    type = ParsedAux
    variable = porosity_times_thickness
    coupled_variables = porosity
    expression = 'porosity * ${thickness}'
    execute_on = 'initial'
  []
  [permeability_times_thickness]
    type = ParsedAux
    variable = permeability_times_thickness
    coupled_variables = permeability
    expression = 'permeability * ${thickness}'
    execute_on = 'initial'
  []
  [pressure_water]
    type = ADPorousFlowPropertyAux
    variable = pressure_water
    property = pressure
    phase = 0
    execute_on = 'initial timestep_end'
  []
  [saturation_water]
    type = ADPorousFlowPropertyAux
    variable = saturation_water
    property = saturation
    phase = 0
    execute_on = 'initial timestep_end'
  []
  [saturation_gas]
    type = ADPorousFlowPropertyAux
    variable = saturation_gas
    property = saturation
    phase = 1
    execute_on = 'initial timestep_end'
  []
  [density_water]
    type = ADPorousFlowPropertyAux
    variable = density_water
    property = density
    phase = 0
    execute_on = 'initial timestep_end'
  []
  [density_gas]
    type = ADPorousFlowPropertyAux
    variable = density_gas
    property = density
    phase = 1
    execute_on = 'initial timestep_end'
  []
  [x1_water]
    type = ADPorousFlowPropertyAux
    variable = x1_water
    property = mass_fraction
    phase = 0
    fluid_component = 1
    execute_on = 'initial timestep_end'
  []
  [x1_gas]
    type = ADPorousFlowPropertyAux
    variable = x1_gas
    property = mass_fraction
    phase = 1
    fluid_component = 1
    execute_on = 'initial timestep_end'
  []
  [x0_water]
    type = ADPorousFlowPropertyAux
    variable = x0_water
    property = mass_fraction
    phase = 0
    fluid_component = 0
    execute_on = 'initial timestep_end'
  []
  [x0_gas]
    type = ADPorousFlowPropertyAux
    variable = x0_gas
    property = mass_fraction
    phase = 1
    fluid_component = 0
    execute_on = 'initial timestep_end'
  []
  [pc]
    type = ADPorousFlowPropertyAux
    variable = pc
    property = capillary_pressure
    execute_on = 'initial timestep_end'
  []
[]
[FVKernels]
  [mass0]
    type = FVPorousFlowMassTimeDerivative
    variable = pgas
    fluid_component = 0
  []
  [flux0]
    type = FVPorousFlowAdvectiveFlux
    variable = pgas
    fluid_component = 0
  []
  [diff0]
    type = FVPorousFlowDispersiveFlux
    variable = pgas
    fluid_component = 0
    disp_long = '0 0'
    disp_trans = '0 0'
  []
  [mass1]
    type = FVPorousFlowMassTimeDerivative
    variable = z
    fluid_component = 1
  []
  [flux1]
    type = FVPorousFlowAdvectiveFlux
    variable = z
    fluid_component = 1
  []
  [diff1]
    type = FVPorousFlowDispersiveFlux
    variable = z
    fluid_component = 1
    disp_long = '0 0'
    disp_trans = '0 0'
  []
[]
[DiracKernels]
  [injector1]
    type = ConstantPointSource
    point = '0.9 0.3 0'
    value = ${inj_rate}
    variable = z
  []
  [injector2]
    type = ConstantPointSource
    point = '1.7 0.7 0'
    value = ${inj_rate}
    variable = z
  []
[]
[Controls]
  [injection1]
    type = ConditionalFunctionEnableControl
    enable_objects = 'DiracKernels::injector1'
    conditional_function = injection_schedule1
  []
  [injection2]
    type = ConditionalFunctionEnableControl
    enable_objects = 'DiracKernels::injector2'
    conditional_function = injection_schedule2
  []
[]
[Functions]
  [initial_p]
    type = ParsedFunction
    symbol_names = 'p0 g H rho0'
    symbol_values = '101.325e3 9.81 1.5 1002'
    expression = 'p0 + rho0 * g * (H - y)'
  []
  [injection_schedule1]
    type = ParsedFunction
    expression = 'if(t >= 0 & t <= 1.8e4, 1, 0)'
  []
  [injection_schedule2]
    type = ParsedFunction
    expression = 'if(t >= 8.1e3 & t <= 1.8e4, 1, 0)'
  []
[]
[ICs]
  [p]
    type = FunctionIC
    variable = pgas
    function = initial_p
  []
[]
[FVBCs]
  [pressure_top]
    type = FVPorousFlowAdvectiveFluxBC
    boundary = top
    porepressure_value = 1.01325e5
    variable = pgas
  []
[]
[FluidProperties]
  [water]
    type = Water97FluidProperties
  []
  [watertab]
    type = TabulatedBicubicFluidProperties
    fp = water
    save_file = false
    pressure_min = 1e5
    pressure_max = 1e6
    temperature_min = 290
    temperature_max = 300
    num_p = 20
    num_T = 10
  []
  [co2]
    type = CO2FluidProperties
  []
  [co2tab]
    type = TabulatedBicubicFluidProperties
    fp = co2
    save_file = false
    pressure_min = 1e5
    pressure_max = 1e6
    temperature_min = 290
    temperature_max = 300
    num_p = 20
    num_T = 10
  []
  [brine]
    type = BrineFluidProperties
    water_fp = watertab
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pgas z'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [sandESF_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${sandESF_pe}
    lambda = 2
    block = ${sandESF}
    pc_max = 1e4
    sat_lr = ${sandESF_swi}
  []
  [sandC_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${sandC_pe}
    lambda = 2
    block = ${sandC}
    pc_max = 1e4
    sat_lr = ${sandC_swi}
  []
  [sandD_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${sandD_pe}
    lambda = 2
    block = ${sandD}
    pc_max = 1e4
    sat_lr = ${sandD_swi}
  []
  [sandE_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${sandE_pe}
    lambda = 2
    block = ${sandE}
    pc_max = 1e4
    sat_lr = ${sandE_swi}
  []
  [sandF_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${sandF_pe}
    lambda = 2
    block = ${sandF}
    pc_max = 1e4
    sat_lr = ${sandF_swi}
  []
  [sandG_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${sandG_pe}
    lambda = 2
    block = ${sandG}
    pc_max = 1e4
    sat_lr = ${sandG_swi}
  []
  [fault1_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${fault1_pe}
    lambda = 2
    block = ${fault1}
    pc_max = 1e4
    sat_lr = ${fault1_swi}
  []
  [fault3_pc]
    type = PorousFlowCapillaryPressureBC
    pe = ${fault3_pe}
    lambda = 2
    block = ${fault3}
    pc_max = 1e4
    sat_lr = ${fault3_swi}
  []
  [top_layer_pc]
    type = PorousFlowCapillaryPressureConst
    pc = 0
    block =  ${top_layer}
  []
  [sandESF_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = sandESF_pc
  []
  [sandC_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = sandC_pc
  []
  [sandD_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = sandD_pc
  []
  [sandE_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = sandE_pc
  []
  [sandF_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = sandF_pc
  []
  [sandG_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = sandG_pc
  []
  [fault1_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = fault1_pc
  []
  [fault3_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = fault3_pc
  []
  [top_layer_fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2tab
    capillary_pressure = top_layer_pc
  []
[]
[Materials]
  [temperature]
    type = ADPorousFlowTemperature
    temperature = temperature
  []
  [sandESF_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = sandESF_fs
    capillary_pressure = sandESF_pc
    block = ${sandESF}
  []
  [sandC_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = sandC_fs
    capillary_pressure = sandC_pc
    block = ${sandC}
  []
  [sandD_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = sandD_fs
    capillary_pressure = sandD_pc
    block = ${sandD}
  []
  [sandE_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = sandE_fs
    capillary_pressure = sandE_pc
    block = ${sandE}
  []
  [sandF_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = sandF_fs
    capillary_pressure = sandF_pc
    block = ${sandF}
  []
  [sandG_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = sandG_fs
    capillary_pressure = sandG_pc
    block = ${sandG}
  []
  [fault1_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = fault1_fs
    capillary_pressure = fault1_pc
    block = ${fault1}
  []
  [fault3_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = fault3_fs
    capillary_pressure = fault3_pc
    block = ${fault3}
  []
  [top_layer_brineco2]
    type = ADPorousFlowFluidState
    gas_porepressure = pgas
    z = z
    temperature_unit = Celsius
    xnacl = xnacl
    fluid_state = top_layer_fs
    capillary_pressure = top_layer_pc
    block = ${top_layer}
  []
  [porosity]
    type = ADPorousFlowPorosityConst
    porosity = porosity_times_thickness
  []
  [permeability]
    type = ADPorousFlowPermeabilityConstFromVar
    perm_xx = permeability_times_thickness
    perm_yy = permeability_times_thickness
    perm_zz = permeability_times_thickness
  []
  [diffcoeff]
    type = ADPorousFlowDiffusivityConst
    tortuosity = '1 1'
    diffusion_coeff = '2e-9 2e-9 0 0'
  []
  [sandESF_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${sandESF_swi}
    sum_s_res = ${fparse sandESF_sgi + sandESF_swi}
    scaling = ${sandESF_krw}
    block = ${sandESF}
  []
  [sandESF_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${sandESF_sgi}
    sum_s_res = ${fparse sandESF_sgi + sandESF_swi}
    scaling = ${sandESF_krg}
    block = ${sandESF}
  []
  [sandC_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${sandC_swi}
    sum_s_res = ${fparse sandC_sgi + sandC_swi}
    scaling = ${sandC_krw}
    block = ${sandC}
  []
  [sandC_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${sandC_sgi}
    sum_s_res = ${fparse sandC_sgi + sandC_swi}
    scaling = ${sandC_krg}
    block = ${sandC}
  []
  [sandD_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${sandD_swi}
    sum_s_res = ${fparse sandD_sgi + sandD_swi}
    scaling = ${sandD_krw}
    block = ${sandD}
  []
  [sandD_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${sandD_sgi}
    sum_s_res = ${fparse sandD_sgi + sandD_swi}
    scaling = ${sandD_krg}
    block = ${sandD}
  []
  [sandE_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${sandE_swi}
    sum_s_res = ${fparse sandE_sgi + sandE_swi}
    scaling = ${sandE_krw}
    block = ${sandE}
  []
  [sandE_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${sandE_sgi}
    sum_s_res = ${fparse sandE_sgi + sandE_swi}
    scaling = ${sandE_krg}
    block = ${sandE}
  []
  [sandF_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${sandF_swi}
    sum_s_res = ${fparse sandF_sgi + sandF_swi}
    scaling = ${sandF_krw}
    block = ${sandF}
  []
  [sandF_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${sandF_sgi}
    sum_s_res = ${fparse sandF_sgi + sandF_swi}
    scaling = ${sandF_krg}
    block = ${sandF}
  []
  [sandG_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${sandG_swi}
    sum_s_res = ${fparse sandG_sgi + sandG_swi}
    scaling = ${sandG_krw}
    block = ${sandG}
  []
  [sandG_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${sandG_sgi}
    sum_s_res = ${fparse sandG_sgi + sandG_swi}
    scaling = ${sandG_krg}
    block = ${sandG}
  []
  [fault1_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${fault1_swi}
    sum_s_res = ${fparse fault1_sgi + fault1_swi}
    scaling = ${fault1_krw}
    block = ${fault1}
  []
  [fault1_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${fault1_sgi}
    sum_s_res = ${fparse fault1_sgi + fault1_swi}
    scaling = ${fault1_krg}
    block = ${fault1}
  []
  [fault3_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    s_res = ${fault3_swi}
    sum_s_res = ${fparse fault3_sgi + fault3_swi}
    scaling = ${fault3_krw}
    block = ${fault3}
  []
  [fault3_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    s_res = ${fault3_sgi}
    sum_s_res = ${fparse fault3_sgi + fault3_swi}
    scaling = ${fault3_krg}
    block = ${fault3}
  []
  [top_layer_relperm0]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 0
    lambda = 2
    block = ${top_layer}
  []
  [top_layer_relperm1]
    type = ADPorousFlowRelativePermeabilityBC
    phase = 1
    nw_phase = true
    lambda = 2
    block = ${top_layer}
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options = '-ksp_snes_ew'
    petsc_options_iname = '-ksp_type -pc_type -pc_factor_mat_solver_package -sub_pc_factor_shift_type'
    petsc_options_value = 'gmres lu mumps NONZERO'
    # petsc_options_iname = '-ksp_type -pc_type -pc_hypre_type -sub_pc_type -sub_pc_factor_shift_type -sub_pc_factor_levels -ksp_gmres_restart'
    # petsc_options_value = 'gmres hypre boomeramg lu NONZERO 4 301'
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  dtmax = 60
  start_time = 0
  end_time = 4.32e5
  nl_rel_tol = 1e-6
  nl_abs_tol = 1e-8
  nl_max_its = 15
  l_tol = 1e-5
  l_abs_tol = 1e-8
  # line_search = none # Can be a useful option for this problem
  [TimeSteppers]
    [time]
      type = FunctionDT
      growth_factor = 2
      cutback_factor_at_failure = 0.5
      function = 'if(t<1.8e4, 2, if(t<3.6e4, 20, 60))'
    []
  []
[]
[Postprocessors]
  [p_5_3]
    type = PointValue
    variable = pgas
    point = '0.5 0.3 0'
    execute_on = 'initial timestep_end'
  []
  [p_5_3_w]
    type = PointValue
    variable = pressure_water
    point = '0.5 0.3 0'
    execute_on = 'initial timestep_end'
  []
  [p_5_7]
    type = PointValue
    variable = pgas
    point = '0.5 0.7 0'
    execute_on = 'initial timestep_end'
  []
  [p_5_7_w]
    type = PointValue
    variable = pressure_water
    point = '0.5 0.7 0'
    execute_on = 'initial timestep_end'
  []
  [p_9_3]
    type = PointValue
    variable = pgas
    point = '0.9 0.3 0'
    execute_on = 'initial timestep_end'
  []
  [p_9_3_w]
    type = PointValue
    variable = pressure_water
    point = '0.9 0.3 0'
    execute_on = 'initial timestep_end'
  []
  [p_15_5]
    type = PointValue
    variable = pgas
    point = '1.5 0.5 0'
    execute_on = 'initial timestep_end'
  []
  [p_15_5_w]
    type = PointValue
    variable = pressure_water
    point = '1.5 0.5 0'
    execute_on = 'initial timestep_end'
  []
  [p_17_7]
    type = PointValue
    variable = pgas
    point = '1.7 0.7 0'
    execute_on = 'initial timestep_end'
  []
  [p_17_7_w]
    type = PointValue
    variable = pressure_water
    point = '1.7 0.7 0'
    execute_on = 'initial timestep_end'
  []
  [p_17_11]
    type = PointValue
    variable = pgas
    point = '1.7 1.1 0'
    execute_on = 'initial timestep_end'
  []
  [p_17_11_w]
    type = PointValue
    variable = pressure_water
    point = '1.7 1.1 0'
    execute_on = 'initial timestep_end'
  []
  [x0mass]
    type = FVPorousFlowFluidMass
    fluid_component = 0
    phase = '0 1'
    execute_on = 'initial timestep_end'
  []
  [x1mass]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '0 1'
    execute_on = 'initial timestep_end'
  []
  [x1gas]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '1'
    execute_on = 'initial timestep_end'
  []
  [boxA]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '0 1'
    block = ${boxA}
    execute_on = 'initial timestep_end'
  []
  [imm_A_sandESF]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandESF_sgi}
    block = 10
    execute_on = 'initial timestep_end'
  []
  [imm_A_sandE]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandE_sgi}
    block = 13
    execute_on = 'initial timestep_end'
  []
  [imm_A_sandF]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandF_sgi}
    block = '14 34'
    execute_on = 'initial timestep_end'
  []
  [imm_A_sandG]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandG_sgi}
    block = '15 35'
    execute_on = 'initial timestep_end'
  []
  [imm_A]
    type = LinearCombinationPostprocessor
    pp_names = 'imm_A_sandESF imm_A_sandE imm_A_sandF imm_A_sandG'
    pp_coefs = '1 1 1 1'
    execute_on = 'initial timestep_end'
  []
  [diss_A]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 0
    block = ${boxA}
    execute_on = 'initial timestep_end'
  []
  [seal_A]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '0 1'
    block = ${seal_boxA}
    execute_on = 'initial timestep_end'
  []
  [boxB]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '0 1'
    block = ${boxB}
    execute_on = 'initial timestep_end'
  []
  [imm_B_sandESF]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandESF_sgi}
    block = 20
    execute_on = 'initial timestep_end'
  []
  [imm_B_sandC]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandC_sgi}
    block = 21
    execute_on = 'initial timestep_end'
  []
  [imm_B_sandD]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandD_sgi}
    block = 22
    execute_on = 'initial timestep_end'
  []
  [imm_B_sandE]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandE_sgi}
    block = 23
    execute_on = 'initial timestep_end'
  []
  [imm_B_sandF]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${sandF_sgi}
    block = 24
    execute_on = 'initial timestep_end'
  []
  [imm_B_fault1]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 1
    saturation_threshold = ${fault1_sgi}
    block = 26
    execute_on = 'initial timestep_end'
  []
  [imm_B]
    type = LinearCombinationPostprocessor
    pp_names = 'imm_B_sandESF imm_B_sandC imm_B_sandD imm_B_sandE imm_B_sandF imm_B_fault1'
    pp_coefs = '1 1 1 1 1 1'
    execute_on = 'initial timestep_end'
  []
  [diss_B]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = 0
    block = ${boxB}
    execute_on = 'initial timestep_end'
  []
  [seal_B]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '0 1'
    block = ${seal_boxB}
    execute_on = 'initial timestep_end'
  []
  [boxC]
    type = FVPorousFlowFluidMass
    fluid_component = 1
    phase = '0'
    block = ${boxC}
    execute_on = 'initial timestep_end'
  []
[]
[Outputs]
  print_linear_residuals = false
  perf_graph = true
  # exodus = true
  [csv]
     type = CSV
  []
[]
(modules/porous_flow/examples/solute_tracer_transport/solute_tracer_transport.i)
# Longitudinal dispersivity
disp = 0.7
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 100
    xmin = 0
    xmax = 100
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
    initial_condition = 1e5
  []
  [C]
    initial_condition = 0
  []
[]
[AuxVariables]
  [Darcy_vel_x]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  [Darcy_vel_x]
    type = PorousFlowDarcyVelocityComponent
    variable = Darcy_vel_x
    component = x
    fluid_phase = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure C'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[Kernels]
  [mass_der_water]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [adv_pp]
    type = PorousFlowFullySaturatedDarcyFlow
    variable = porepressure
    fluid_component = 1
  []
  [diff_pp]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = porepressure
    disp_trans = 0
    disp_long = ${disp}
  []
  [mass_der_C]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = C
  []
  [adv_C]
    type = PorousFlowFullySaturatedDarcyFlow
    fluid_component = 0
    variable = C
  []
  [diff_C]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = C
    disp_trans = 0
    disp_long = ${disp}
  []
[]
[FluidProperties]
  [water]
    type = Water97FluidProperties
  []
[]
[Materials]
  [ps]
    type = PorousFlow1PhaseFullySaturated
    porepressure = porepressure
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.25
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-11 0 0   0 1E-11 0   0 0 1E-11'
  []
  [water]
    type = PorousFlowSingleComponentFluid
    fp = water
    phase = 0
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = C
  []
  [temperature]
    type = PorousFlowTemperature
    temperature = 293
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '0 0'
    tortuosity = 0.1
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    kr = 1
    phase = 0
  []
[]
[BCs]
  [constant_inlet_pressure]
    type = DirichletBC
    variable = porepressure
    value = 1.2e5
    boundary = left
  []
  [constant_outlet_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1e5
    boundary = right
  []
  [inlet_tracer]
    type = DirichletBC
    variable = C
    value = 0.001
    boundary = left
  []
  [outlet_tracer]
    type = PorousFlowOutflowBC
    variable = C
    boundary = right
    mass_fraction_component = 0
  []
[]
[Preconditioning]
  [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'
  []
[]
[Executioner]
  type = Transient
  end_time = 17280000
  dtmax = 86400
  nl_rel_tol = 1e-6
  nl_abs_tol = 1e-12
  [TimeStepper]
    type = IterationAdaptiveDT
    dt = 1000
  []
[]
[Postprocessors]
  [C]
    type = PointValue
    variable = C
    point = '50 0 0'
  []
  [Darcy_x]
    type = PointValue
    variable = Darcy_vel_x
    point = '50 0 0'
  []
[]
[Outputs]
  file_base = solute_tracer_transport_${disp}
  csv = true
[]
(modules/porous_flow/test/tests/jacobian/disp04.i)
# Test the Jacobian of the PorousFlowDisperiveFlux kernel
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[ICs]
  [pp]
    type = RandomIC
    variable = pp
    max = 2e1
    min = 1e1
  []
  [massfrac0]
    type = RandomIC
    variable = massfrac0
    min = 0
    max = 1
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    gravity = '1 0 0'
    disp_long = 0.2
    disp_trans = 0.1
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    gravity = '1 0 0'
    disp_long = 0.2
    disp_trans = 0.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 10
    thermal_expansion = 0
    viscosity = 1
  []
[]
[Materials]
  [temp]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac0'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-2 1e-1'
    tortuosity = '0.1'
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1 0 0 0 2 0 0 0 3'
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
[]
[Preconditioning]
  active = smp
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1
  end_time = 1
[]
[Outputs]
  exodus = false
[]
(modules/porous_flow/examples/flow_through_fractured_media/fine_thick_fracture_transient.i)
# Using a single-dimensional mesh
# Transient flow and solute transport along a fracture in a porous matrix
# advective dominated flow in the fracture and diffusion into the porous matrix
#
# Note that fine_thick_fracture_steady.i must be run to initialise the porepressure properly
[Mesh]
  file = 'gold/fine_thick_fracture_steady_out.e'
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
    initial_from_file_var = pp
    initial_from_file_timestep = 1
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [velocity_x]
    family = MONOMIAL
    order = CONSTANT
    block = fracture
  []
  [velocity_y]
    family = MONOMIAL
    order = CONSTANT
    block = fracture
  []
[]
[AuxKernels]
  [velocity_x]
    type = PorousFlowDarcyVelocityComponent
    variable = velocity_x
    component = x
  []
  [velocity_y]
    type = PorousFlowDarcyVelocityComponent
    variable = velocity_y
    component = y
  []
[]
[Problem]
  # massfrac0 has an initial condition despite the restart
  allow_initial_conditions_with_restart = true
[]
[ICs]
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[BCs]
  [top]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = top
  []
  [bottom]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = bottom
  []
  [ptop]
    type = DirichletBC
    variable = pp
    boundary = top
    value = 1e6
  []
  [pbottom]
    type = DirichletBC
    variable = pp
    boundary = bottom
    value = 1.002e6
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    disp_trans = 0
    disp_long = 0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = massfrac0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro_fracture]
    type = PorousFlowPorosityConst
    porosity = 1.0 # this is the true porosity of the fracture
    block = 'fracture'
  []
  [poro_matrix]
    type = PorousFlowPorosityConst
    porosity = 0.1
    block = 'matrix1 matrix2'
  []
  [diff1]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-9 1e-9'
    tortuosity = 1.0
    block = 'fracture'
  []
  [diff2]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-9 1e-9'
    tortuosity = 0.1
    block = 'matrix1 matrix2'
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability1]
    type = PorousFlowPermeabilityConst
    permeability = '3e-8 0 0 0 3e-8 0 0 0 3e-8' # this is the true permeability of the fracture
    block = 'fracture'
  []
  [permeability2]
    type = PorousFlowPermeabilityConst
    permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
    block = 'matrix1 matrix2'
  []
[]
[Functions]
  [dt_controller]
    type = PiecewiseConstant
    x = '0    30   40 100 200 83200'
    y = '0.01 0.1  1  10  100 32'
  []
[]
[Preconditioning]
  active = basic
  [mumps_is_best_for_parallel_jobs]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
  [basic]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 86400
  #dt = 0.01
  [TimeStepper]
    type = FunctionDT
    function = dt_controller
  []
  # controls for nonlinear iterations
  nl_max_its = 15
  nl_rel_tol = 1e-14
  nl_abs_tol = 1e-9
[]
[VectorPostprocessors]
  [xmass]
    type = LineValueSampler
    start_point = '0.4 0 0'
    end_point = '0.5 0 0'
    sort_by = x
    num_points = 167
    variable = massfrac0
  []
[]
[Outputs]
  perf_graph = true
  console = true
  csv = true
  exodus = true
[]
(modules/porous_flow/test/tests/dispersion/diff01_fv.i)
# Test diffusive part of FVPorousFlowDispersiveFlux kernel by setting dispersion
# coefficients to zero. Pressure is held constant over the mesh, and gravity is
# set to zero so that no advective transport of mass takes place.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
  [mesh]
    type = GeneratedMeshGenerator
    dim = 1
    nx = 20
    xmax = 10
    bias_x = 1.2
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
    type = MooseVariableFVReal
  []
  [massfrac0]
    type = MooseVariableFVReal
  []
[]
[AuxVariables]
  [velocity]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [velocity]
    type = ADPorousFlowDarcyVelocityComponent
    variable = velocity
    component = x
  []
[]
[ICs]
  [pp]
    type = ConstantIC
    variable = pp
    value = 1e5
  []
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[FVBCs]
  [left]
    type = FVDirichletBC
    value = 1
    variable = massfrac0
    boundary = left
  []
  [right]
    type = FVDirichletBC
    value = 0
    variable = massfrac0
    boundary = right
  []
  [pright]
    type = FVDirichletBC
    variable = pp
    boundary = right
    value = 1e5
  []
  [pleft]
    type = FVDirichletBC
    variable = pp
    boundary = left
    value = 1e5
  []
[]
[FVKernels]
  [mass0]
    type = FVPorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [diff0_pp]
    type = FVPorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    disp_trans = 0
    disp_long = 0
  []
  [mass1]
    type = FVPorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [diff1_x]
    type = FVPorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 1000
    viscosity = 0.001
    thermal_expansion = 0
  []
[]
[Materials]
  [temperature]
    type = ADPorousFlowTemperature
  []
  [ppss]
    type = ADPorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = ADPorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = ADPorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = ADPorousFlowPorosityConst
    porosity = 0.3
  []
  [diff]
    type = ADPorousFlowDiffusivityConst
    diffusion_coeff = '1 1'
    tortuosity = 0.1
  []
  [relp]
    type = ADPorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability]
    type = ADPorousFlowPermeabilityConst
    permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  dt = 1
  end_time = 20
[]
[VectorPostprocessors]
  [xmass]
    type = ElementValueSampler
    sort_by = id
    variable = massfrac0
  []
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/examples/flow_through_fractured_media/coarse.i)
# Flow and solute transport along a fracture embedded in a porous matrix
# The fracture is represented by lower dimensional elements
# fracture aperture = 6e-4m
# fracture porosity = 6e-4m = phi * a
# fracture permeability = 1.8e-11 which is based on k=3e-8 from a**2/12, and k*a = 3e-8*6e-4
# matrix porosity = 0.1
# matrix permeanility = 1e-20
[Mesh]
  type = FileMesh
  file = 'coarse.e'
  block_id = '1 2 3'
  block_name = 'fracture matrix1 matrix2'
  boundary_id = '1 2'
  boundary_name = 'bottom top'
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [velocity_x]
    family = MONOMIAL
    order = CONSTANT
    block = 'fracture'
  []
  [velocity_y]
    family = MONOMIAL
    order = CONSTANT
    block = 'fracture'
  []
[]
[AuxKernels]
  [velocity_x]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_x
    component = x
    aperture = 6E-4
  []
  [velocity_y]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_y
    component = y
    aperture = 6E-4
  []
[]
[ICs]
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
  [pp_matrix]
    type = ConstantIC
    variable = pp
    value = 1E6
  []
[]
[BCs]
  [top]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = top
  []
  [bottom]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = bottom
  []
  [ptop]
    type = DirichletBC
    variable = pp
    boundary =  top
    value = 1e6
  []
  [pbottom]
    type = DirichletBC
    variable = pp
    boundary = bottom
    value = 1.002e6
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = pp
    disp_trans = 0
    disp_long = 0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = massfrac0
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = massfrac0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = massfrac0
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro_fracture]
    type = PorousFlowPorosityConst
    porosity = 6e-4   # = a * phif
    block = 'fracture'
  []
  [poro_matrix]
    type = PorousFlowPorosityConst
    porosity = 0.1
    block = 'matrix1 matrix2'
  []
  [diff1]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-9 1e-9'
    tortuosity = 1.0
    block = 'fracture'
  []
  [diff2]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-9 1e-9'
    tortuosity = 0.1
    block = 'matrix1 matrix2'
  []
  [permeability_fracture]
    type = PorousFlowPermeabilityConst
    permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11'   # 1.8e-11 = a * kf
    block = 'fracture'
  []
  [permeability_matrix]
    type = PorousFlowPermeabilityConst
    permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
    block = 'matrix1 matrix2'
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
[]
[Preconditioning]
  [basic]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 10
  dt = 1
# controls for nonlinear iterations
  nl_max_its = 15
  nl_rel_tol = 1e-14
  nl_abs_tol = 1e-12
[]
[VectorPostprocessors]
  [xmass]
    type = LineValueSampler
    start_point = '-0.5 0 0'
    end_point = '0.5 0 0'
    sort_by = x
    num_points = 41
    variable = massfrac0
    outputs = csv
  []
[]
[Outputs]
  [csv]
    type = CSV
    execute_on = 'final'
  []
[]
(modules/porous_flow/test/tests/jacobian/disp02.i)
# Test the Jacobian of the dispersive contribution to the diffusive component of
# the PorousFlowDisperiveFlux kernel along with a non-zero diffusion.
# By setting disp_long and disp_trans to the same non-zero value, the purely
# dispersive component of the flux is zero, and the only flux is due to diffusion
# and its contribution from disp_trans.
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[ICs]
  [pp]
    type = RandomIC
    variable = pp
    max = 2e1
    min = 1e1
  []
  [massfrac0]
    type = RandomIC
    variable = massfrac0
    min = 0
    max = 1
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    gravity = '1 0 0'
    disp_long = 0.1
    disp_trans = 0.1
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    gravity = '1 0 0'
    disp_long = 0.1
    disp_trans = 0.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 10
    thermal_expansion = 0
    viscosity = 1
  []
[]
[Materials]
  [temp]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac0'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-2 1e-1'
    tortuosity = '0.1'
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1 0 0 0 2 0 0 0 3'
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
[]
[Preconditioning]
  active = smp
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1
  end_time = 1
[]
[Outputs]
  exodus = false
[]
(modules/porous_flow/test/tests/dispersion/diff01.i)
# Test diffusive part of PorousFlowDispersiveFlux kernel by setting dispersion
# coefficients to zero. Pressure is held constant over the mesh, and gravity is
# set to zero so that no advective transport of mass takes place.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 20
  xmax = 10
  bias_x = 1.1
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [velocity]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [velocity]
    type = PorousFlowDarcyVelocityComponent
    variable = velocity
    component = x
  []
[]
[ICs]
  [pp]
    type = ConstantIC
    variable = pp
    value = 1e5
  []
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[BCs]
  [left]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = left
  []
  [right]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = right
  []
  [pright]
    type = DirichletBC
    variable = pp
    boundary = right
    value = 1e5
  []
  [pleft]
    type = DirichletBC
    variable = pp
    boundary = left
    value = 1e5
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    disp_trans = 0
    disp_long = 0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = massfrac0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 1000
    viscosity = 0.001
    thermal_expansion = 0
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.3
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1 1'
    tortuosity = 0.1
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  dt = 1
  end_time = 20
[]
[VectorPostprocessors]
  [xmass]
    type = NodalValueSampler
    sort_by = id
    variable = massfrac0
  []
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/jacobian/diff01.i)
# Test the Jacobian of the diffusive component of the PorousFlowDisperiveFlux kernel.
# By setting disp_long and disp_trans to zero, the purely diffusive component of the flux
# can be isolated.
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[ICs]
  [pp]
    type = RandomIC
    variable = pp
    max = 2e1
    min = 1e1
  []
  [massfrac0]
    type = RandomIC
    variable = massfrac0
    min = 0
    max = 1
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    gravity = '1 0 0'
    disp_long = 0
    disp_trans = 0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    gravity = '1 0 0'
    disp_long = 0
    disp_trans = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 10
    thermal_expansion = 0
    viscosity = 1
  []
[]
[Materials]
  [temp]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac0'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-2 1e-1'
    tortuosity = '0.1'
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1 0 0 0 2 0 0 0 3'
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
[]
[Preconditioning]
  active = smp
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1
  end_time = 1
[]
[Outputs]
  exodus = false
[]
(modules/porous_flow/test/tests/chemistry/2species_equilibrium_2phase.i)
# Using a two-phase system (see 2species_equilibrium for the single-phase)
# The saturations, porosity, mass fractions, tortuosity and diffusion coefficients are chosen so that the results are identical to 2species_equilibrium
#
# PorousFlow analogy of chemical_reactions/test/tests/aqueous_equilibrium/2species.i
#
# Simple equilibrium reaction example to illustrate the use of PorousFlowMassFractionAqueousEquilibriumChemistry
#
# In this example, two primary species a and b are transported by diffusion and convection
# from the left of the porous medium, reacting to form two equilibrium species pa2 and pab
# according to the equilibrium reaction:
#
#      reactions = '2a = pa2     rate = 10^2
#                   a + b = pab  rate = 10^-2'
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 10
[]
[Variables]
  [a]
    order = FIRST
    family = LAGRANGE
    [InitialCondition]
      type = BoundingBoxIC
      x1 = 0.0
      y1 = 0.0
      x2 = 1.0e-10
      y2 = 1
      inside = 1.0e-2
      outside = 1.0e-10
    []
  []
  [b]
    order = FIRST
    family = LAGRANGE
    [InitialCondition]
      type = BoundingBoxIC
      x1 = 0.0
      y1 = 0.0
      x2 = 1.0e-10
      y2 = 1
      inside = 1.0e-2
      outside = 1.0e-10
    []
  []
[]
[AuxVariables]
  [eqm_k0]
    initial_condition = 1E2
  []
  [eqm_k1]
    initial_condition = 1E-2
  []
  [pressure0]
  []
  [saturation1]
    initial_condition = 0.25
  []
  [a_in_phase0]
    initial_condition = 0.0
  []
  [b_in_phase0]
    initial_condition = 0.0
  []
  [pa2]
    family = MONOMIAL
    order = CONSTANT
  []
  [pab]
    family = MONOMIAL
    order = CONSTANT
  []
[]
[AuxKernels]
  [pa2]
    type = PorousFlowPropertyAux
    property = secondary_concentration
    secondary_species = 0
    variable = pa2
  []
  [pab]
    type = PorousFlowPropertyAux
    property = secondary_concentration
    secondary_species = 1
    variable = pab
  []
[]
[ICs]
  [pressure0]
    type = FunctionIC
    variable = pressure0
    function = 2-x
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Kernels]
  [mass_a]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = a
  []
  [flux_a]
    type = PorousFlowAdvectiveFlux
    variable = a
    fluid_component = 0
  []
  [diff_a]
    type = PorousFlowDispersiveFlux
    variable = a
    fluid_component = 0
    disp_trans = '0 0'
    disp_long = '0 0'
  []
  [mass_b]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = b
  []
  [flux_b]
    type = PorousFlowAdvectiveFlux
    variable = b
    fluid_component = 1
  []
  [diff_b]
    type = PorousFlowDispersiveFlux
    variable = b
    fluid_component = 1
    disp_trans = '0 0'
    disp_long = '0 0'
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'a b'
    number_fluid_phases = 2
    number_fluid_components = 3
    number_aqueous_equilibrium = 2
    aqueous_phase_number = 1
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9 # huge, so mimic chemical_reactions
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow2PhasePS
    capillary_pressure = pc
    phase0_porepressure = pressure0
    phase1_saturation = saturation1
  []
  [massfrac]
    type = PorousFlowMassFractionAqueousEquilibriumChemistry
    mass_fraction_vars = 'a_in_phase0 b_in_phase0 a b'
    num_reactions = 2
    equilibrium_constants = 'eqm_k0 eqm_k1'
    primary_activity_coefficients = '1 1'
    secondary_activity_coefficients = '1 1'
    reactions = '2 0
                 1 1'
  []
  [simple_fluid0]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [simple_fluid1]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 1
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.8
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    # porous_flow permeability / porous_flow viscosity = chemical_reactions conductivity = 1E-4
    permeability = '1E-7 0 0 0 1E-7 0 0 0 1E-7'
  []
  [relp0]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [relp1]
    type = PorousFlowRelativePermeabilityConst
    phase = 1
  []
  [diff]
    type = PorousFlowDiffusivityConst
    # porous_flow diffusion_coeff * tortuousity * porosity = chemical_reactions diffusivity = 1E-4
    diffusion_coeff = '5E-4 5E-4 5E-4
                       5E-4 5E-4 5E-4'
    tortuosity = '0.25 0.25'
  []
[]
[BCs]
  [a_left]
    type = DirichletBC
    variable = a
    boundary = left
    value = 1.0e-2
  []
  [b_left]
    type = DirichletBC
    variable = b
    boundary = left
    value = 1.0e-2
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 10
  end_time = 100
[]
[Outputs]
  print_linear_residuals = true
  exodus = true
  perf_graph = true
[]
(modules/porous_flow/examples/lava_lamp/1phase_convection.i)
# Two phase density-driven convection of dissolved CO2 in brine
#
# The model starts with CO2 in the liquid phase only.  The CO2 diffuses into the brine.
# As the density of the CO2-saturated brine is greater
# than the unsaturated brine, a gravitational instability arises and density-driven
# convection of CO2-rich fingers descend into the unsaturated brine.
#
# The instability is seeded by a random perturbation to the porosity field.
# Mesh adaptivity is used to refine the mesh as the fingers form.
#
# Note: this model is computationally expensive, so should be run with multiple cores.
[GlobalParams]
  PorousFlowDictator = 'dictator'
  gravity = '0 -9.81 0'
[]
[Adaptivity]
  max_h_level = 2
  marker = marker
  initial_marker = initial
  initial_steps = 2
  [Indicators]
    [indicator]
      type = GradientJumpIndicator
      variable = zi
    []
  []
  [Markers]
    [marker]
      type = ErrorFractionMarker
      indicator = indicator
      refine = 0.8
    []
    [initial]
      type = BoxMarker
      bottom_left = '0 1.95 0'
      top_right = '2 2 0'
      inside = REFINE
      outside = DO_NOTHING
    []
  []
[]
[Mesh]
  type = GeneratedMesh
  dim = 2
  ymin = 1.5
  ymax = 2
  xmax = 2
  ny = 20
  nx = 40
  bias_y = 0.95
[]
[AuxVariables]
  [xnacl]
    initial_condition = 0.01
  []
  [saturation_gas]
    order = FIRST
    family = MONOMIAL
  []
  [xco2l]
    order = FIRST
    family = MONOMIAL
  []
  [density_liquid]
    order = FIRST
    family = MONOMIAL
  []
  [porosity]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  [saturation_gas]
    type = PorousFlowPropertyAux
    variable = saturation_gas
    property = saturation
    phase = 1
    execute_on = 'timestep_end'
  []
  [xco2l]
    type = PorousFlowPropertyAux
    variable = xco2l
    property = mass_fraction
    phase = 0
    fluid_component = 1
    execute_on = 'timestep_end'
  []
  [density_liquid]
    type = PorousFlowPropertyAux
    variable = density_liquid
    property = density
    phase = 0
    execute_on = 'timestep_end'
  []
[]
[Variables]
  [pgas]
  []
  [zi]
    scaling = 1e4
  []
[]
[ICs]
  [pressure]
    type = FunctionIC
    function = 10e6-9.81*1000*y
    variable = pgas
  []
  [zi]
    type = ConstantIC
    value = 0
    variable = zi
  []
  [porosity]
    type = RandomIC
    variable = porosity
    min = 0.25
    max = 0.275
    seed = 0
  []
[]
[BCs]
  [top]
    type = DirichletBC
    value = 0.04
    variable = zi
    boundary = top
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pgas
  []
  [flux0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pgas
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pgas
    disp_long = '0 0'
    disp_trans = '0 0'
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = zi
  []
  [flux1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = zi
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = zi
    disp_long = '0 0'
    disp_trans = '0 0'
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pgas zi'
    number_fluid_phases = 2
    number_fluid_components = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureConst
    pc = 0
  []
  [fs]
    type = PorousFlowBrineCO2
    brine_fp = brine
    co2_fp = co2
    capillary_pressure = pc
  []
[]
[FluidProperties]
  [co2sw]
    type = CO2FluidProperties
  []
  [co2]
    type = TabulatedBicubicFluidProperties
    fp = co2sw
  []
  [brine]
    type = BrineFluidProperties
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
    temperature = '45'
  []
  [brineco2]
    type = PorousFlowFluidState
    gas_porepressure = 'pgas'
    z = 'zi'
    temperature_unit = Celsius
    xnacl = 'xnacl'
    capillary_pressure = pc
    fluid_state = fs
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = porosity
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-11 0 0 0 1e-11 0 0 0 1e-11'
  []
  [relperm_water]
    type = PorousFlowRelativePermeabilityCorey
    phase = 0
    n = 2
    s_res = 0.1
    sum_s_res = 0.2
  []
  [relperm_gas]
    type = PorousFlowRelativePermeabilityCorey
    phase = 1
    n = 2
    s_res = 0.1
    sum_s_res = 0.2
  []
  [diffusivity]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '2e-9 2e-9 2e-9 2e-9'
    tortuosity = '1 1'
  []
[]
[Preconditioning]
  active = basic
  [mumps_is_best_for_parallel_jobs]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
  [basic]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 1e6
  nl_max_its = 25
  l_max_its = 100
  dtmax = 1e4
  nl_abs_tol = 1e-6
  [TimeStepper]
    type = IterationAdaptiveDT
    dt = 100
    growth_factor = 2
    cutback_factor = 0.5
  []
[]
[Functions]
  [flux]
    type = ParsedFunction
    symbol_values = 'delta_xco2 dt'
    symbol_names = 'dx dt'
    expression = 'dx/dt'
  []
[]
[Postprocessors]
  [total_co2_in_gas]
    type = PorousFlowFluidMass
    phase = 1
    fluid_component = 1
  []
  [total_co2_in_liquid]
    type = PorousFlowFluidMass
    phase = 0
    fluid_component = 1
  []
  [numdofs]
    type = NumDOFs
  []
  [delta_xco2]
    type = ChangeOverTimePostprocessor
    postprocessor = total_co2_in_liquid
  []
  [dt]
    type = TimestepSize
  []
  [flux]
    type = FunctionValuePostprocessor
    function = flux
  []
[]
[Outputs]
  print_linear_residuals = false
  perf_graph = true
  exodus = true
  csv = true
[]
(modules/porous_flow/test/tests/dispersion/disp01.i)
# Test dispersive part of PorousFlowDispersiveFlux kernel by setting diffusion
# coefficients to zero. A pressure gradient is applied over the mesh to give a
# uniform velocity. Gravity is set to zero.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 100
  xmax = 10
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [velocity]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [velocity]
    type = PorousFlowDarcyVelocityComponent
    variable = velocity
    component = x
  []
[]
[ICs]
  [pp]
    type = FunctionIC
    variable = pp
    function = pic
  []
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[Functions]
  [pic]
    type = ParsedFunction
    expression = 1.1e5-x*1e3
  []
[]
[BCs]
  [xleft]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = left
  []
  [xright]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = right
  []
  [pright]
    type = DirichletBC
    variable = pp
    boundary = right
    value = 1e5
  []
  [pleft]
    type = DirichletBC
    variable = pp
    boundary = left
    value = 1.1e5
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    variable = pp
    disp_trans = 0
    disp_long = 0.2
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = massfrac0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0.2
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e9
    density0 = 1000
    viscosity = 0.001
    thermal_expansion = 0
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.3
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '0 0'
    tortuosity = 0.1
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 1e3
  dtmax = 50
  [TimeStepper]
    type = IterationAdaptiveDT
    growth_factor = 1.5
    cutback_factor = 0.5
    dt = 1
  []
[]
[VectorPostprocessors]
  [xmass]
    type = NodalValueSampler
    sort_by = id
    variable = massfrac0
  []
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/test/tests/jacobian/disp01.i)
# Test the Jacobian of the dispersive contribution to the diffusive component of
# the PorousFlowDisperiveFlux kernel. By setting disp_long and disp_trans to the same
# non-zero value, and diffusion to zero (by setting tortuosity to zero), the purely
# dispersive component of the flux is zero, and the only flux is due to the contribution
# from disp_trans on the diffusive flux.
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[ICs]
  [pp]
    type = RandomIC
    variable = pp
    max = 2e1
    min = 1e1
  []
  [massfrac0]
    type = RandomIC
    variable = massfrac0
    min = 0
    max = 1
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    gravity = '1 0 0'
    disp_long = 0.1
    disp_trans = 0.1
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    gravity = '1 0 0'
    disp_long = 0.1
    disp_trans = 0.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 10
    thermal_expansion = 0
    viscosity = 1
  []
[]
[Materials]
  [temp]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac0'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-2 1e-1'
    tortuosity = 1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1 0 0 0 2 0 0 0 3'
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
[]
[Preconditioning]
  active = smp
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1
  end_time = 1
[]
[Outputs]
  exodus = false
[]
(modules/porous_flow/test/tests/jacobian/disp03.i)
# Test the Jacobian of the dispersive contribution to the PorousFlowDisperiveFlux
# kernel by setting the diffusive component to zero (tortuosity = 0).
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 3
  xmin = 0
  xmax = 1
  ny = 1
  ymin = 0
  ymax = 1
[]
[GlobalParams]
  PorousFlowDictator = dictator
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[ICs]
  [pp]
    type = RandomIC
    variable = pp
    max = 2e1
    min = 1e1
  []
  [massfrac0]
    type = RandomIC
    variable = massfrac0
    min = 0
    max = 1
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    gravity = '1 0 0'
    disp_long = 0.2
    disp_trans = 0.1
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    gravity = '1 0 0'
    disp_long = 0.2
    disp_trans = 0.1
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e7
    density0 = 10
    thermal_expansion = 0
    viscosity = 1
  []
[]
[Materials]
  [temp]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = 'massfrac0'
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.1
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-2 1e-1'
    tortuosity = 1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1 0 0 0 2 0 0 0 3'
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
[]
[Preconditioning]
  active = smp
  [smp]
    type = SMP
    full = true
  []
[]
[Executioner]
  type = Transient
  solve_type = Newton
  dt = 1
  end_time = 1
[]
[Outputs]
  exodus = false
[]
(modules/porous_flow/test/tests/dispersion/diff01_action.i)
# Test diffusive part of PorousFlowDispersiveFlux kernel by setting dispersion
# coefficients to zero. Pressure is held constant over the mesh, and gravity is
# set to zero so that no advective transport of mass takes place.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 20
  xmax = 10
  bias_x = 1.1
[]
[GlobalParams]
  PorousFlowDictator = andy_heheheh
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[ICs]
  [pp]
    type = ConstantIC
    variable = pp
    value = 1e5
  []
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[BCs]
  [left]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = left
  []
  [right]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = right
  []
  [pright]
    type = DirichletBC
    variable = pp
    boundary = right
    value = 1e5
  []
  [pleft]
    type = DirichletBC
    variable = pp
    boundary = left
    value = 1e5
  []
[]
[Kernels]
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = massfrac0
    disp_trans = 0
    disp_long = 0
    gravity = '0 0 0'
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = pp
    disp_trans = 0
    disp_long = 0
    gravity = '0 0 0'
  []
[]
[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    thermal_expansion = 0.0
    bulk_modulus = 1E7
    viscosity = 0.001
    density0 = 1000.0
  []
[]
[PorousFlowUnsaturated]
  porepressure = pp
  gravity = '0 0 0'
  fp = the_simple_fluid
  dictator_name = andy_heheheh
  relative_permeability_type = Corey
  relative_permeability_exponent = 0.0
  mass_fraction_vars = massfrac0
[]
[Materials]
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.3
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1 1'
    tortuosity = 0.1
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  dt = 1
  end_time = 20
[]
[VectorPostprocessors]
  [xmass]
    type = NodalValueSampler
    sort_by = id
    variable = massfrac0
  []
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/examples/flow_through_fractured_media/fine_transient.i)
# Using a mixed-dimensional mesh
# Transient flow and solute transport along a fracture in a porous matrix
# advective dominated flow in the fracture and diffusion into the porous matrix
#
# Note that fine_steady.i must be run to initialise the porepressure properly
[Mesh]
  file = 'gold/fine_steady_out.e'
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [pp]
    initial_from_file_var = pp
    initial_from_file_timestep = 1
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [velocity_x]
    family = MONOMIAL
    order = CONSTANT
    block = fracture
  []
  [velocity_y]
    family = MONOMIAL
    order = CONSTANT
    block = fracture
  []
[]
[AuxKernels]
  [velocity_x]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_x
    component = x
    aperture = 6E-4
  []
  [velocity_y]
    type = PorousFlowDarcyVelocityComponentLowerDimensional
    variable = velocity_y
    component = y
    aperture = 6E-4
  []
[]
[Problem]
  # massfrac0 has an initial condition despite the restart
  allow_initial_conditions_with_restart = true
[]
[ICs]
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[BCs]
  [top]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = top
  []
  [bottom]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = bottom
  []
  [ptop]
    type = DirichletBC
    variable = pp
    boundary = top
    value = 1e6
  []
  [pbottom]
    type = DirichletBC
    variable = pp
    boundary = bottom
    value = 1.002e6
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = pp
    disp_trans = 0
    disp_long = 0
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = massfrac0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2e9
    density0 = 1000
    thermal_expansion = 0
    viscosity = 1e-3
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro_fracture]
    type = PorousFlowPorosityConst
    porosity = 6e-4 # = a * phif
    block = 'fracture'
  []
  [poro_matrix]
    type = PorousFlowPorosityConst
    porosity = 0.1
    block = 'matrix1 matrix2'
  []
  [diff1]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-9 1e-9'
    tortuosity = 1.0
    block = 'fracture'
  []
  [diff2]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '1e-9 1e-9'
    tortuosity = 0.1
    block = 'matrix1 matrix2'
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability_fracture]
    type = PorousFlowPermeabilityConst
    permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11' # kf=3e-8, a=6e-4m.  1.8e-11 = kf * a
    block = 'fracture'
  []
  [permeability_matrix]
    type = PorousFlowPermeabilityConst
    permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
    block = 'matrix1 matrix2'
  []
[]
[Functions]
  [dt_controller]
    type = PiecewiseConstant
    x = '0    30   40 100 200 83200'
    y = '0.01 0.1  1  10  100 32'
  []
[]
[Preconditioning]
  active = basic
  [mumps_is_best_for_parallel_jobs]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
  [basic]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 86400
  [TimeStepper]
    type = FunctionDT
    function = dt_controller
  []
  # controls for nonlinear iterations
  nl_max_its = 15
  nl_rel_tol = 1e-14
  nl_abs_tol = 1e-9
[]
[VectorPostprocessors]
  [xmass]
    type = LineValueSampler
    start_point = '0.4 0 0'
    end_point = '0.5 0 0'
    sort_by = x
    num_points = 167
    variable = massfrac0
  []
[]
[Outputs]
  perf_graph = true
  console = true
  csv = true
  exodus = true
[]
(modules/porous_flow/test/tests/dispersion/disp01_heavy.i)
# Test dispersive part of PorousFlowDispersiveFlux kernel by setting diffusion
# coefficients to zero. A pressure gradient is applied over the mesh to give a
# uniform velocity. Gravity is set to zero.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
  type = GeneratedMesh
  dim = 1
  nx = 200
  xmax = 10
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
  compute_enthalpy = false
  compute_internal_energy = false
[]
[Variables]
  [pp]
  []
  [massfrac0]
  []
[]
[AuxVariables]
  [velocity]
    family = MONOMIAL
    order = FIRST
  []
[]
[AuxKernels]
  [velocity]
    type = PorousFlowDarcyVelocityComponent
    variable = velocity
    component = x
  []
[]
[ICs]
  [pp]
    type = FunctionIC
    variable = pp
    function = pic
  []
  [massfrac0]
    type = ConstantIC
    variable = massfrac0
    value = 0
  []
[]
[Functions]
  [pic]
    type = ParsedFunction
    expression = 1.1e5-x*1e3
  []
[]
[BCs]
  [xleft]
    type = DirichletBC
    value = 1
    variable = massfrac0
    boundary = left
  []
  [xright]
    type = DirichletBC
    value = 0
    variable = massfrac0
    boundary = right
  []
  [pright]
    type = DirichletBC
    variable = pp
    boundary = right
    value = 1e5
  []
  [pleft]
    type = DirichletBC
    variable = pp
    boundary = left
    value = 1.1e5
  []
[]
[Kernels]
  [mass0]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = pp
  []
  [adv0]
    type = PorousFlowAdvectiveFlux
    fluid_component = 0
    variable = pp
  []
  [diff0]
    type = PorousFlowDispersiveFlux
    variable = pp
    disp_trans = 0
    disp_long = 0.2
  []
  [mass1]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = massfrac0
  []
  [adv1]
    type = PorousFlowAdvectiveFlux
    fluid_component = 1
    variable = massfrac0
  []
  [diff1]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = massfrac0
    disp_trans = 0
    disp_long = 0.2
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'pp massfrac0'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[FluidProperties]
  [simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 1e9
    density0 = 1000
    viscosity = 0.001
    thermal_expansion = 0
  []
[]
[Materials]
  [temperature]
    type = PorousFlowTemperature
  []
  [ppss]
    type = PorousFlow1PhaseFullySaturated
    porepressure = pp
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = massfrac0
  []
  [simple_fluid]
    type = PorousFlowSingleComponentFluid
    fp = simple_fluid
    phase = 0
  []
  [poro]
    type = PorousFlowPorosityConst
    porosity = 0.3
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '0 0'
    tortuosity = 0.1
  []
  [relp]
    type = PorousFlowRelativePermeabilityConst
    phase = 0
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
  []
[]
[Preconditioning]
  [smp]
    type = SMP
    full = true
    petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = 'gmres      asm      lu           NONZERO                   2             '
  []
[]
[Executioner]
  type = Transient
  solve_type = NEWTON
  end_time = 1e3
  dtmax = 10
  [TimeStepper]
    type = IterationAdaptiveDT
    growth_factor = 1.5
    cutback_factor = 0.5
    dt = 1
  []
[]
[VectorPostprocessors]
  [xmass]
    type = NodalValueSampler
    sort_by = id
    variable = massfrac0
  []
[]
[Outputs]
  [out]
    type = CSV
    execute_on = final
  []
[]
(modules/porous_flow/examples/solute_tracer_transport/solute_tracer_transport_2D.i)
# Longitudinal dispersivity
disp = 5
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 100
    xmin = -50
    xmax = 50
    ny = 60
    ymin = 0
    ymax = 50
  []
[]
[GlobalParams]
  PorousFlowDictator = dictator
  gravity = '0 0 0'
[]
[Variables]
  [porepressure]
    initial_condition = 1e5
  []
  [C]
    initial_condition = 0
  []
[]
[AuxVariables]
  [Darcy_vel_x]
    order = CONSTANT
    family = MONOMIAL
  []
  [Darcy_vel_y]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[AuxKernels]
  [Darcy_vel_x]
    type = PorousFlowDarcyVelocityComponent
    variable = Darcy_vel_x
    component = x
    fluid_phase = 0
  []
  [Darcy_vel_y]
    type = PorousFlowDarcyVelocityComponent
    variable = Darcy_vel_y
    component = y
    fluid_phase = 0
  []
[]
[UserObjects]
  [dictator]
    type = PorousFlowDictator
    porous_flow_vars = 'porepressure C'
    number_fluid_phases = 1
    number_fluid_components = 2
  []
[]
[Kernels]
  [mass_der_water]
    type = PorousFlowMassTimeDerivative
    fluid_component = 1
    variable = porepressure
  []
  [adv_pp]
    type = PorousFlowFullySaturatedDarcyFlow
    variable = porepressure
    fluid_component = 1
  []
  [diff_pp]
    type = PorousFlowDispersiveFlux
    fluid_component = 1
    variable = porepressure
    disp_trans = 0
    disp_long = ${disp}
  []
  [mass_der_C]
    type = PorousFlowMassTimeDerivative
    fluid_component = 0
    variable = C
  []
  [adv_C]
    type = PorousFlowFullySaturatedDarcyFlow
    fluid_component = 0
    variable = C
  []
  [diff_C]
    type = PorousFlowDispersiveFlux
    fluid_component = 0
    variable = C
    disp_trans = 0
    disp_long = ${disp}
  []
[]
[FluidProperties]
  [water]
    type = Water97FluidProperties
  []
[]
[Materials]
  [ps]
    type = PorousFlow1PhaseFullySaturated
    porepressure = porepressure
  []
  [porosity]
    type = PorousFlowPorosityConst
    porosity = 0.25
  []
  [permeability]
    type = PorousFlowPermeabilityConst
    permeability = '1E-11 0 0   0 1E-11 0   0 0 1E-11'
  []
  [water]
    type = PorousFlowSingleComponentFluid
    fp = water
    phase = 0
  []
  [massfrac]
    type = PorousFlowMassFraction
    mass_fraction_vars = C
  []
  [temperature]
    type = PorousFlowTemperature
    temperature = 293
  []
  [diff]
    type = PorousFlowDiffusivityConst
    diffusion_coeff = '0 0'
    tortuosity = 0.1
  []
  [relperm]
    type = PorousFlowRelativePermeabilityConst
    kr = 1
    phase = 0
  []
[]
[DiracKernels]
  [source_P]
    type = PorousFlowSquarePulsePointSource
    point = '0 0 0'
    mass_flux = 1e-1
    variable = porepressure
  []
  [source_C]
    type = PorousFlowSquarePulsePointSource
    point = '0 0 0'
    mass_flux = 1e-7
    variable = C
  []
[]
[BCs]
  [constant_outlet_porepressure_]
    type = DirichletBC
    variable = porepressure
    value = 1e5
    boundary = 'top left right'
  []
  [outlet_tracer_top]
    type = PorousFlowOutflowBC
    variable = C
    boundary = top
    mass_fraction_component = 0
  []
  [outlet_tracer_right]
    type = PorousFlowOutflowBC
    variable = C
    boundary = right
    mass_fraction_component = 0
  []
  [outlet_tracer_left]
    type = PorousFlowOutflowBC
    variable = C
    boundary = left
    mass_fraction_component = 0
  []
[]
[Preconditioning]
  [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'
  []
[]
[Executioner]
  type = Transient
  end_time = 17280000
  dtmax = 100000
  nl_rel_tol = 1e-6
  nl_abs_tol = 1e-12
  [TimeStepper]
    type = IterationAdaptiveDT
    dt = 1000
  []
[]
[Postprocessors]
  [C]
    type = PointValue
    variable = C
    point = '0 25 0'
  []
  [Darcy_x]
    type = PointValue
    variable = Darcy_vel_x
    point = '0 25 0'
  []
  [Darcy_y]
    type = PointValue
    variable = Darcy_vel_y
    point = '0 25 0'
  []
[]
[Outputs]
  file_base = solute_tracer_transport_2D_${disp}
  csv = true
  exodus = true
[]