Capillary pressure test descriptions

The details of the capillary pressure curve implementations can be found in here.

Brooks-Corey

Figure 1: Brooks-Corey Capillary Pressure test cases.

There are two test cases for Brooks-Corey including the log-extension turning on, and log-extension turning off at low saturation.

The input file for log-extension off test:

# Test Brooks-Corey capillary pressure curve by varying saturation over the mesh
# lambda = 2, sat_lr = 0.1, log_extension = false

[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 1
  nx = 500
[]

[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
  PorousFlowDictator = dictator
[]

[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
  [p0]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e6
  []
  [s1]
  []
[]

[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
  [s0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [s1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
[]

[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
  [s0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s0aux
  []
  [s1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s1aux
  []
  [p0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p0aux
  []
  [p1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p1aux
  []
[]

[Functions<<<{"href": "../../../../syntax/Functions/index.html"}>>>]
  [s1]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = x
  []
[]

[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
  [s1]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = s1
    function<<<{"description": "The initial condition function."}>>> = s1
  []
[]

[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
  [p0]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p0
  []
  [s1]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
  []
[]

[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
  [dictator]
    type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
    porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation.  Jacobian entries involving derivatives wrt these variables will be computed.  In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'p0 s1'
    number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
    number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureBC<<<{"description": "Brooks-Corey capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureBC.html"}>>>
    lambda<<<{"description": "Brooks-Corey exponent lambda"}>>> = 2
    log_extension<<<{"description": "Use a logarithmic extension for low saturation to avoid capillary pressure going to infinity. Default is true.  Set to false if your capillary pressure depends on spatially-dependent variables other than saturation, as the log-extension C++ code for this case has yet to be implemented"}>>> = false
    pe<<<{"description": "Brooks-Corey entry pressure. Must be positive"}>>> = 1e5
    sat_lr<<<{"description": "Liquid residual saturation.  Must be between 0 and 1. Default is 0"}>>> = 0.1
  []
[]

[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
  [temperature]
    type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
  []
  [ppss]
    type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
    phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = p0
    phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = s1
    capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
  []
  [kr0]
    type = PorousFlowRelativePermeabilityVG<<<{"description": "This Material calculates relative permeability of a phase using the van Genuchten model", "href": "../../../../source/materials/PorousFlowRelativePermeabilityVG.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 0
    m<<<{"description": "The van Genuchten exponent of the phase"}>>> = 0.5
  []
  [kr1]
    type = PorousFlowRelativePermeabilityCorey<<<{"description": "This Material calculates relative permeability of the fluid phase, using the simple Corey model ((S-S_res)/(1-sum(S_res)))^n", "href": "../../../../source/materials/PorousFlowRelativePermeabilityCorey.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 1
    n<<<{"description": "The Corey exponent of the phase."}>>> = 2
  []
[]

[VectorPostprocessors<<<{"href": "../../../../syntax/VectorPostprocessors/index.html"}>>>]
  [vpp]
    type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = 's0aux s1aux p0aux p1aux'
    start_point<<<{"description": "The beginning of the line"}>>> = '0 0 0'
    end_point<<<{"description": "The ending of the line"}>>> = '1 0 0'
    num_points<<<{"description": "The number of points to sample along the line"}>>> = 500
    sort_by<<<{"description": "What to sort the samples by"}>>> = id
  []
[]

[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = NEWTON
  nl_abs_tol = 1e-6
[]

[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
  [sleft]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = left
  []
  [sright]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 1
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = right
  []
[]

[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
[]
(moose/modules/porous_flow/test/tests/capillary_pressure/brooks_corey1.i)

The input file for log-extension on test:

# Test Brooks-Corey capillary pressure curve by varying saturation over the mesh
# lambda = 2, sat_lr = 0.1, log_extension = true

[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 1
  nx = 500
[]

[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
  PorousFlowDictator = dictator
[]

[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
  [p0]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e6
  []
  [s1]
  []
[]

[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
  [s0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [s1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
[]

[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
  [s0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s0aux
  []
  [s1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s1aux
  []
  [p0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p0aux
  []
  [p1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p1aux
  []
[]

[Functions<<<{"href": "../../../../syntax/Functions/index.html"}>>>]
  [s1]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = x
  []
[]

[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
  [s1]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = s1
    function<<<{"description": "The initial condition function."}>>> = s1
  []
[]

[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
  [p0]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p0
  []
  [s1]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
  []
[]

[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
  [dictator]
    type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
    porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation.  Jacobian entries involving derivatives wrt these variables will be computed.  In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'p0 s1'
    number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
    number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureBC<<<{"description": "Brooks-Corey capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureBC.html"}>>>
    lambda<<<{"description": "Brooks-Corey exponent lambda"}>>> = 2
    log_extension<<<{"description": "Use a logarithmic extension for low saturation to avoid capillary pressure going to infinity. Default is true.  Set to false if your capillary pressure depends on spatially-dependent variables other than saturation, as the log-extension C++ code for this case has yet to be implemented"}>>> = true
    pe<<<{"description": "Brooks-Corey entry pressure. Must be positive"}>>> = 1e5
    sat_lr<<<{"description": "Liquid residual saturation.  Must be between 0 and 1. Default is 0"}>>> = 0.1
  []
[]

[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
  [temperature]
    type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
  []
  [ppss]
    type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
    phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = p0
    phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = s1
    capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
  []
  [kr0]
    type = PorousFlowRelativePermeabilityVG<<<{"description": "This Material calculates relative permeability of a phase using the van Genuchten model", "href": "../../../../source/materials/PorousFlowRelativePermeabilityVG.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 0
    m<<<{"description": "The van Genuchten exponent of the phase"}>>> = 0.5
  []
  [kr1]
    type = PorousFlowRelativePermeabilityCorey<<<{"description": "This Material calculates relative permeability of the fluid phase, using the simple Corey model ((S-S_res)/(1-sum(S_res)))^n", "href": "../../../../source/materials/PorousFlowRelativePermeabilityCorey.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 1
    n<<<{"description": "The Corey exponent of the phase."}>>> = 2
  []
[]

[VectorPostprocessors<<<{"href": "../../../../syntax/VectorPostprocessors/index.html"}>>>]
  [vpp]
    type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = 's0aux s1aux p0aux p1aux'
    start_point<<<{"description": "The beginning of the line"}>>> = '0 0 0'
    end_point<<<{"description": "The ending of the line"}>>> = '1 0 0'
    num_points<<<{"description": "The number of points to sample along the line"}>>> = 500
    sort_by<<<{"description": "What to sort the samples by"}>>> = id
  []
[]

[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = NEWTON
  nl_abs_tol = 1e-6
[]

[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
  [sleft]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = left
  []
  [sright]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 1
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = right
  []
[]

[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
[]
(moose/modules/porous_flow/test/tests/capillary_pressure/brooks_corey2.i)

van Genuchten

Figure 2: van Genuchten Capillary Pressure test cases.

There are three test cases for van Genuchten including the log-extension turning on, log-extension turning off at low saturation, and changing the scale factor.

The input file for log-extension off test:

# Test van Genuchten relative permeability curve by varying saturation over the mesh
# van Genuchten exponent m = 0.5 for both phases
# No residual saturation in either phase

[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 1
  nx = 500
[]

[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
  PorousFlowDictator = dictator
[]

[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
  [p0]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e6
  []
  [s1]
  []
[]

[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
  [s0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [s1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
[]

[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
  [s0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s0aux
  []
  [s1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s1aux
  []
  [p0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p0aux
  []
  [p1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p1aux
  []
[]

[Functions<<<{"href": "../../../../syntax/Functions/index.html"}>>>]
  [s1]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = x
  []
[]

[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
  [s1]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = s1
    function<<<{"description": "The initial condition function."}>>> = s1
  []
[]

[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
  [p0]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p0
  []
  [s1]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
  []
[]

[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
  [dictator]
    type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
    porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation.  Jacobian entries involving derivatives wrt these variables will be computed.  In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'p0 s1'
    number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
    number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
    alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1e-5
    m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
    sat_lr<<<{"description": "Liquid residual saturation.  Must be between 0 and 1. Default is 0"}>>> = 0.1
    log_extension<<<{"description": "Use a logarithmic extension for low saturation to avoid capillary pressure going to infinity. Default is true.  Set to false if your capillary pressure depends on spatially-dependent variables other than saturation, as the log-extension C++ code for this case has yet to be implemented"}>>> = false
  []
[]

[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
  [temperature]
    type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
  []
  [ppss]
    type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
    phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = p0
    phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = s1
    capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
  []
  [kr0]
    type = PorousFlowRelativePermeabilityVG<<<{"description": "This Material calculates relative permeability of a phase using the van Genuchten model", "href": "../../../../source/materials/PorousFlowRelativePermeabilityVG.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 0
    m<<<{"description": "The van Genuchten exponent of the phase"}>>> = 0.5
  []
  [kr1]
    type = PorousFlowRelativePermeabilityCorey<<<{"description": "This Material calculates relative permeability of the fluid phase, using the simple Corey model ((S-S_res)/(1-sum(S_res)))^n", "href": "../../../../source/materials/PorousFlowRelativePermeabilityCorey.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 1
    n<<<{"description": "The Corey exponent of the phase."}>>> = 2
  []
[]

[VectorPostprocessors<<<{"href": "../../../../syntax/VectorPostprocessors/index.html"}>>>]
  [vpp]
    type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = 's0aux s1aux p0aux p1aux'
    start_point<<<{"description": "The beginning of the line"}>>> = '0 0 0'
    end_point<<<{"description": "The ending of the line"}>>> = '1 0 0'
    num_points<<<{"description": "The number of points to sample along the line"}>>> = 500
    sort_by<<<{"description": "What to sort the samples by"}>>> = id
  []
[]

[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = NEWTON
  nl_abs_tol = 1e-6
[]

[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
  [sleft]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = left
  []
  [sright]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 1
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = right
  []
[]

[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
[]
(moose/modules/porous_flow/test/tests/capillary_pressure/vangenuchten1.i)

The input file for log-extension on test:

# Test van Genuchten relative permeability curve by varying saturation over the mesh
# van Genuchten exponent m = 0.5 for both phases
# No residual saturation in either phase

[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 1
  nx = 500
[]

[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
  PorousFlowDictator = dictator
[]

[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
  [p0]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e6
  []
  [s1]
  []
[]

[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
  [s0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [s1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
[]

[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
  [s0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s0aux
  []
  [s1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s1aux
  []
  [p0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p0aux
  []
  [p1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p1aux
  []
[]

[Functions<<<{"href": "../../../../syntax/Functions/index.html"}>>>]
  [s1]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = x
  []
[]

[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
  [s1]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = s1
    function<<<{"description": "The initial condition function."}>>> = s1
  []
[]

[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
  [p0]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p0
  []
  [s1]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
  []
[]

[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
  [dictator]
    type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
    porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation.  Jacobian entries involving derivatives wrt these variables will be computed.  In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'p0 s1'
    number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
    number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
    alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1e-5
    m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
    sat_lr<<<{"description": "Liquid residual saturation.  Must be between 0 and 1. Default is 0"}>>> = 0.1
    log_extension<<<{"description": "Use a logarithmic extension for low saturation to avoid capillary pressure going to infinity. Default is true.  Set to false if your capillary pressure depends on spatially-dependent variables other than saturation, as the log-extension C++ code for this case has yet to be implemented"}>>> = true
  []
[]

[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
  [temperature]
    type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
  []
  [ppss]
    type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
    phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = p0
    phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = s1
    capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
  []
  [kr0]
    type = PorousFlowRelativePermeabilityVG<<<{"description": "This Material calculates relative permeability of a phase using the van Genuchten model", "href": "../../../../source/materials/PorousFlowRelativePermeabilityVG.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 0
    m<<<{"description": "The van Genuchten exponent of the phase"}>>> = 0.5
  []
  [kr1]
    type = PorousFlowRelativePermeabilityCorey<<<{"description": "This Material calculates relative permeability of the fluid phase, using the simple Corey model ((S-S_res)/(1-sum(S_res)))^n", "href": "../../../../source/materials/PorousFlowRelativePermeabilityCorey.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 1
    n<<<{"description": "The Corey exponent of the phase."}>>> = 2
  []
[]

[VectorPostprocessors<<<{"href": "../../../../syntax/VectorPostprocessors/index.html"}>>>]
  [vpp]
    type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = 's0aux s1aux p0aux p1aux'
    start_point<<<{"description": "The beginning of the line"}>>> = '0 0 0'
    end_point<<<{"description": "The ending of the line"}>>> = '1 0 0'
    num_points<<<{"description": "The number of points to sample along the line"}>>> = 500
    sort_by<<<{"description": "What to sort the samples by"}>>> = id
  []
[]

[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = NEWTON
  nl_abs_tol = 1e-6
[]

[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
  [sleft]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = left
  []
  [sright]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 1
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = right
  []
[]

[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
[]
(moose/modules/porous_flow/test/tests/capillary_pressure/vangenuchten2.i)

The input file for applied-scaling test:

# Test van Genuchten relative permeability curve by varying saturation over the mesh
# van Genuchten exponent m = 0.5 for both phases
# No residual saturation in either phase

[Mesh<<<{"href": "../../../../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 1
  nx = 500
[]

[GlobalParams<<<{"href": "../../../../syntax/GlobalParams/index.html"}>>>]
  PorousFlowDictator = dictator
[]

[Variables<<<{"href": "../../../../syntax/Variables/index.html"}>>>]
  [p0]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e6
  []
  [s1]
  []
[]

[AuxVariables<<<{"href": "../../../../syntax/AuxVariables/index.html"}>>>]
  [s0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [s1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p0aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
  [p1aux]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
[]

[AuxKernels<<<{"href": "../../../../syntax/AuxKernels/index.html"}>>>]
  [s0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s0aux
  []
  [s1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = saturation
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = s1aux
  []
  [p0]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 0
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p0aux
  []
  [p1]
    type = PorousFlowPropertyAux<<<{"description": "AuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.", "href": "../../../../source/auxkernels/PorousFlowPropertyAux.html"}>>>
    property<<<{"description": "The fluid property that this auxillary kernel is to calculate"}>>> = pressure
    phase<<<{"description": "The index of the phase this auxillary kernel acts on"}>>> = 1
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = p1aux
  []
[]

[Functions<<<{"href": "../../../../syntax/Functions/index.html"}>>>]
  [s1]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = x
  []
[]

[ICs<<<{"href": "../../../../syntax/ICs/index.html"}>>>]
  [s1]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = s1
    function<<<{"description": "The initial condition function."}>>> = s1
  []
[]

[Kernels<<<{"href": "../../../../syntax/Kernels/index.html"}>>>]
  [p0]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p0
  []
  [s1]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
  []
[]

[UserObjects<<<{"href": "../../../../syntax/UserObjects/index.html"}>>>]
  [dictator]
    type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../../../../source/userobjects/PorousFlowDictator.html"}>>>
    porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation.  Jacobian entries involving derivatives wrt these variables will be computed.  In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'p0 s1'
    number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 2
    number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 2
  []
  [pc]
    type = PorousFlowCapillaryPressureVG<<<{"description": "van Genuchten capillary pressure", "href": "../../../../source/userobjects/PorousFlowCapillaryPressureVG.html"}>>>
    alpha<<<{"description": "van Genuchten parameter alpha. Must be positive"}>>> = 1e-5
    m<<<{"description": "van Genuchten exponent m. Must be between 0 and 1, and optimally should be set to >0.5"}>>> = 0.5
    sat_lr<<<{"description": "Liquid residual saturation.  Must be between 0 and 1. Default is 0"}>>> = 0.1
    s_scale<<<{"description": "CapillaryPressure = f(Seff * s_scale) - f(s_scale), where f is the van Genuchten expression.  Setting s_scale<1 is unusual but sometimes helps fully saturated, 2-phase PP simulations converge as the zero derivative (1/f'(S=1)=0) is removed"}>>> = 0.8
    log_extension<<<{"description": "Use a logarithmic extension for low saturation to avoid capillary pressure going to infinity. Default is true.  Set to false if your capillary pressure depends on spatially-dependent variables other than saturation, as the log-extension C++ code for this case has yet to be implemented"}>>> = false
  []
[]

[Materials<<<{"href": "../../../../syntax/Materials/index.html"}>>>]
  [temperature]
    type = PorousFlowTemperature<<<{"description": "Material to provide temperature at the quadpoints or nodes and derivatives of it with respect to the PorousFlow variables", "href": "../../../../source/materials/PorousFlowTemperature.html"}>>>
  []
  [ppss]
    type = PorousFlow2PhasePS<<<{"description": "This Material calculates the 2 porepressures and the 2 saturations in a 2-phase situation, and derivatives of these with respect to the PorousFlowVariables.", "href": "../../../../source/materials/PorousFlow2PhasePS.html"}>>>
    phase0_porepressure<<<{"description": "Variable that is the porepressure of phase 0 (the liquid phase)"}>>> = p0
    phase1_saturation<<<{"description": "Variable that is the saturation of phase 1 (the gas phase)"}>>> = s1
    capillary_pressure<<<{"description": "Name of the UserObject defining the capillary pressure"}>>> = pc
  []
  [kr0]
    type = PorousFlowRelativePermeabilityVG<<<{"description": "This Material calculates relative permeability of a phase using the van Genuchten model", "href": "../../../../source/materials/PorousFlowRelativePermeabilityVG.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 0
    m<<<{"description": "The van Genuchten exponent of the phase"}>>> = 0.5
  []
  [kr1]
    type = PorousFlowRelativePermeabilityCorey<<<{"description": "This Material calculates relative permeability of the fluid phase, using the simple Corey model ((S-S_res)/(1-sum(S_res)))^n", "href": "../../../../source/materials/PorousFlowRelativePermeabilityCorey.html"}>>>
    phase<<<{"description": "The phase number"}>>> = 1
    n<<<{"description": "The Corey exponent of the phase."}>>> = 2
  []
[]

[VectorPostprocessors<<<{"href": "../../../../syntax/VectorPostprocessors/index.html"}>>>]
  [vpp]
    type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = 's0aux s1aux p0aux p1aux'
    start_point<<<{"description": "The beginning of the line"}>>> = '0 0 0'
    end_point<<<{"description": "The ending of the line"}>>> = '1 0 0'
    num_points<<<{"description": "The number of points to sample along the line"}>>> = 500
    sort_by<<<{"description": "What to sort the samples by"}>>> = id
  []
[]

[Executioner<<<{"href": "../../../../syntax/Executioner/index.html"}>>>]
  type = Steady
  solve_type = NEWTON
  nl_abs_tol = 1e-6
[]

[BCs<<<{"href": "../../../../syntax/BCs/index.html"}>>>]
  [sleft]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = left
  []
  [sright]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = s1
    value<<<{"description": "Value of the BC"}>>> = 1
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = right
  []
[]

[Outputs<<<{"href": "../../../../syntax/Outputs/index.html"}>>>]
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
[]
(moose/modules/porous_flow/test/tests/capillary_pressure/vangenuchten3.i)