PorousFlowBasicTHM

This Action allows simple simulation of fully-saturated, single-phase, single-component simulations. It is mainly used in poro-elastic and thermo-poro-elastic simulations. The fluid flow may be optionally coupled to mechanics and/or heat flow using the coupling_type flag.

This Action offers simplified physics compared with the PorousFlowFullySaturated and PorousFlowUnsaturated Actions.

Fluid flow

The fluid equation is a simplified form of the full PorousFlow fluid equation (see PorousFlowFullySaturatedMassTimeDerivative and PorousFlowFullySaturatedDarcyBase for the derivation): (1) In this equation, the fluid density, , appears in parentheses, because the user has the option of including it or not using the multiply_by_density flag. Note that the fluid-mass time derivative is close to linear, and is perfectly linear if multiply_by_density=false, and this also almost linearises the flow term. Extremely good nonlinear convergence should therefore be expected, but there are some knock-on effects that are documented in PorousFlowFullySaturatedMassTimeDerivative.

The term only appears if the coupling_type includes "Mechanical". The term only appears if the coupling_type includes "Thermo".

To simulate Eq. (1) the PorousFlowBasicTHM Action employs the following Kernels:

For isothermal simulations, the fluid properties may still depend on temperature, so the temperature input parameter may be set to any real number, or to an AuxVariable if desired.

Upwinding and mass lumping are not used by these Kernels. These numerical schemes are typically unnecessary in fully-saturated, single-phase simulations, but the lack of stabilization means the results from PorousFlowBasicTHM will typically be slightly different to the remainder of PorousFlow.

commentnote

Even though there is only one fluid phase in this fully saturated action, some objects may require a relative permeability material to work. Examples include PorousFlowDarcyVelocityComponent AuxKernels which requires relative permeability, or PorousFlowPeacemanBorehole wells. By default, this action adds a constant relative permeability of one, so that the fluid is perfectly mobile.

Heat flow

For anisothermal simulations, the energy equation reads where the final term is only used if coupling with mechanics is also desired. To simulate this DE, PorousFlowBasicTHM uses the following kernels:

Solid Mechanics

For simulations that couple fluid flow to mechanics, the equations are already written in governing equations, and PorousFlowBasicTHM implements these by using the following kernels:

Required Materials

PorousFlowBasicTHM adds many Materials automatically, however, to run a simulation you will need to provide more Materials for each mesh block, depending on your simulation type, viz:

commentnote

Before choosing non-standard pressure_unit or time_unit, you should read the essay on these non-standard choices here

Examples

A simple example of using PorousFlowBasicTHM is documented in the PorousFlow tutorial with input file:

# Darcy flow
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [annular]
    type = AnnularMeshGenerator<<<{"description": "For rmin>0: creates an annular mesh of QUAD4 elements. For rmin=0: creates a disc mesh of QUAD4 and TRI3 elements. Boundary sidesets are created at rmax and rmin, and given these names. If dmin!=0 and dmax!=360, a sector of an annulus or disc is created. In this case boundary sidesets are also created at dmin and dmax, and given these names", "href": "../meshgenerators/AnnularMeshGenerator.html"}>>>
    nr<<<{"description": "Number of elements in the radial direction"}>>> = 10
    rmin<<<{"description": "Inner radius.  If rmin=0 then a disc mesh (with no central hole) will be created."}>>> = 1.0
    rmax<<<{"description": "Outer radius"}>>> = 10
    growth_r<<<{"description": "The ratio of radial sizes of successive rings of elements"}>>> = 1.4
    nt<<<{"description": "Number of elements in the angular direction"}>>> = 4
    dmin<<<{"description": "Minimum degree, measured in degrees anticlockwise from x axis"}>>> = 0
    dmax<<<{"description": "Maximum angle, measured in degrees anticlockwise from x axis. If dmin=0 and dmax=360 an annular mesh is created. Otherwise, only a sector of an annulus is created"}>>> = 90
  []
  [make3D]
    type = MeshExtruderGenerator<<<{"description": "Takes a 1D or 2D mesh and extrudes the entire structure along the specified axis increasing the dimensionality of the mesh.", "href": "../meshgenerators/MeshExtruderGenerator.html"}>>>
    extrusion_vector<<<{"description": "The direction and length of the extrusion"}>>> = '0 0 12'
    num_layers<<<{"description": "The number of layers in the extruded mesh"}>>> = 3
    bottom_sideset<<<{"description": "The boundary that will be applied to the bottom of the extruded mesh"}>>> = 'bottom'
    top_sideset<<<{"description": "The boundary that will be to the top of the extruded mesh"}>>> = 'top'
    input<<<{"description": "the mesh we want to extrude"}>>> = annular
  []
  [shift_down]
    type = TransformGenerator<<<{"description": "Applies a linear transform to the entire mesh.", "href": "../meshgenerators/TransformGenerator.html"}>>>
    transform<<<{"description": "The type of transformation to perform (TRANSLATE, TRANSLATE_CENTER_ORIGIN, TRANSLATE_MIN_ORIGIN, ROTATE, SCALE)"}>>> = TRANSLATE
    vector_value<<<{"description": "The value to use for the transformation. When using TRANSLATE or SCALE, the xyz coordinates are applied in each direction respectively. When using ROTATE, the values are interpreted as the Euler angles phi, theta and psi given in degrees."}>>> = '0 0 -6'
    input<<<{"description": "The mesh we want to modify"}>>> = make3D
  []
  [aquifer]
    type = SubdomainBoundingBoxGenerator<<<{"description": "Changes the subdomain ID of elements either (XOR) inside or outside the specified box to the specified ID.", "href": "../meshgenerators/SubdomainBoundingBoxGenerator.html"}>>>
    block_id<<<{"description": "Subdomain id to set for inside/outside the bounding box"}>>> = 1
    bottom_left<<<{"description": "The bottom left point (in x,y,z with spaces in-between)."}>>> = '0 0 -2'
    top_right<<<{"description": "The bottom left point (in x,y,z with spaces in-between)."}>>> = '10 10 2'
    input<<<{"description": "The mesh we want to modify"}>>> = shift_down
  []
  [injection_area]
    type = ParsedGenerateSideset<<<{"description": "A MeshGenerator that adds element sides to a sideset if the centroid of the side satisfies the `combinatorial_geometry` expression.", "href": "../meshgenerators/ParsedGenerateSideset.html"}>>>
    combinatorial_geometry<<<{"description": "Function expression encoding a combinatorial geometry"}>>> = 'x*x+y*y<1.01'
    included_subdomains<<<{"description": "A set of subdomain names or ids whose sides will be included in the new sidesets. A side is only added if the subdomain id of the corresponding element is in this set."}>>> = 1
    new_sideset_name<<<{"description": "The name of the new sideset"}>>> = 'injection_area'
    input<<<{"description": "The mesh we want to modify"}>>> = 'aquifer'
  []
  [rename]
    type = RenameBlockGenerator<<<{"description": "Changes the block IDs and/or block names for a given set of blocks defined by either block ID or block name. The changes are independent of ordering. The merging of blocks is supported.", "href": "../meshgenerators/RenameBlockGenerator.html"}>>>
    old_block<<<{"description": "Elements with these block ID(s)/name(s) will be given the new block information specified in 'new_block'"}>>> = '0 1'
    new_block<<<{"description": "The new block ID(s)/name(s) to be given by the elements defined in 'old_block'."}>>> = 'caps aquifer'
    input<<<{"description": "The mesh we want to modify"}>>> = 'injection_area'
  []
[]

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

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [porepressure]
  []
[]

[PorousFlowBasicTHM<<<{"href": "../../syntax/PorousFlowBasicTHM/index.html"}>>>]
  porepressure<<<{"description": "The name of the porepressure variable"}>>> = porepressure
  coupling_type<<<{"description": "The type of simulation.  For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient.  For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material"}>>> = Hydro
  gravity<<<{"description": "Gravitational acceleration vector downwards (m/s^2)"}>>> = '0 0 0'
  fp<<<{"description": "The name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid"}>>> = the_simple_fluid
[]

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

[FluidProperties<<<{"href": "../../syntax/FluidProperties/index.html"}>>>]
  [the_simple_fluid]
    type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../fluidproperties/SimpleFluidProperties.html"}>>>
    bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 2E9
    viscosity<<<{"description": "Constant dynamic viscosity (Pa.s)"}>>> = 1.0E-3
    density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1000.0
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [porosity]
    type = PorousFlowPorosity<<<{"description": "This Material calculates the porosity PorousFlow simulations", "href": "../materials/PorousFlowPorosity.html"}>>>
    porosity_zero<<<{"description": "The porosity at zero volumetric strain and reference temperature and reference effective porepressure and reference chemistry.  This must be a real number or a constant monomial variable (not a linear lagrange or other type of variable)"}>>> = 0.1
  []
  [biot_modulus]
    type = PorousFlowConstantBiotModulus<<<{"description": "Computes the Biot Modulus, which is assumed to be constant for all time.  Sometimes 1 / BiotModulus is called storativity", "href": "../materials/PorousFlowConstantBiotModulus.html"}>>>
    biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.8
    solid_bulk_compliance<<<{"description": "Reciprocal of the drained bulk modulus of the porous skeleton.  If strain = C * stress, then solid_bulk_compliance = de_ij de_kl C_ijkl.  If the grain bulk modulus is Kg then 1/Kg = (1 - biot_coefficient) * solid_bulk_compliance."}>>> = 2E-7
    fluid_bulk_modulus<<<{"description": "Fluid bulk modulus"}>>> = 1E7
  []
  [permeability_aquifer]
    type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../materials/PorousFlowPermeabilityConst.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = aquifer
    permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '1E-14 0 0   0 1E-14 0   0 0 1E-14'
  []
  [permeability_caps]
    type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../materials/PorousFlowPermeabilityConst.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = caps
    permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '1E-15 0 0   0 1E-15 0   0 0 1E-16'
  []
[]

[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
  active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = basic
  [basic]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
  []
  [preferred_but_might_not_be_installed]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = ' lu       mumps'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Transient
  solve_type = Newton
  end_time = 1E6
  dt = 1E5
  nl_abs_tol = 1E-13
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(moose/modules/porous_flow/examples/tutorial/01.i)

A TH example that uses PorousFlowBasicTHM is also documented in the PorousFlow tutorial with input file:

# Darcy flow with heat advection and conduction
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [annular]
    type = AnnularMeshGenerator<<<{"description": "For rmin>0: creates an annular mesh of QUAD4 elements. For rmin=0: creates a disc mesh of QUAD4 and TRI3 elements. Boundary sidesets are created at rmax and rmin, and given these names. If dmin!=0 and dmax!=360, a sector of an annulus or disc is created. In this case boundary sidesets are also created at dmin and dmax, and given these names", "href": "../meshgenerators/AnnularMeshGenerator.html"}>>>
    nr<<<{"description": "Number of elements in the radial direction"}>>> = 10
    rmin<<<{"description": "Inner radius.  If rmin=0 then a disc mesh (with no central hole) will be created."}>>> = 1.0
    rmax<<<{"description": "Outer radius"}>>> = 10
    growth_r<<<{"description": "The ratio of radial sizes of successive rings of elements"}>>> = 1.4
    nt<<<{"description": "Number of elements in the angular direction"}>>> = 4
    dmin<<<{"description": "Minimum degree, measured in degrees anticlockwise from x axis"}>>> = 0
    dmax<<<{"description": "Maximum angle, measured in degrees anticlockwise from x axis. If dmin=0 and dmax=360 an annular mesh is created. Otherwise, only a sector of an annulus is created"}>>> = 90
  []
  [make3D]
    type = MeshExtruderGenerator<<<{"description": "Takes a 1D or 2D mesh and extrudes the entire structure along the specified axis increasing the dimensionality of the mesh.", "href": "../meshgenerators/MeshExtruderGenerator.html"}>>>
    extrusion_vector<<<{"description": "The direction and length of the extrusion"}>>> = '0 0 12'
    num_layers<<<{"description": "The number of layers in the extruded mesh"}>>> = 3
    bottom_sideset<<<{"description": "The boundary that will be applied to the bottom of the extruded mesh"}>>> = 'bottom'
    top_sideset<<<{"description": "The boundary that will be to the top of the extruded mesh"}>>> = 'top'
    input<<<{"description": "the mesh we want to extrude"}>>> = annular
  []
  [shift_down]
    type = TransformGenerator<<<{"description": "Applies a linear transform to the entire mesh.", "href": "../meshgenerators/TransformGenerator.html"}>>>
    transform<<<{"description": "The type of transformation to perform (TRANSLATE, TRANSLATE_CENTER_ORIGIN, TRANSLATE_MIN_ORIGIN, ROTATE, SCALE)"}>>> = TRANSLATE
    vector_value<<<{"description": "The value to use for the transformation. When using TRANSLATE or SCALE, the xyz coordinates are applied in each direction respectively. When using ROTATE, the values are interpreted as the Euler angles phi, theta and psi given in degrees."}>>> = '0 0 -6'
    input<<<{"description": "The mesh we want to modify"}>>> = make3D
  []
  [aquifer]
    type = SubdomainBoundingBoxGenerator<<<{"description": "Changes the subdomain ID of elements either (XOR) inside or outside the specified box to the specified ID.", "href": "../meshgenerators/SubdomainBoundingBoxGenerator.html"}>>>
    block_id<<<{"description": "Subdomain id to set for inside/outside the bounding box"}>>> = 1
    bottom_left<<<{"description": "The bottom left point (in x,y,z with spaces in-between)."}>>> = '0 0 -2'
    top_right<<<{"description": "The bottom left point (in x,y,z with spaces in-between)."}>>> = '10 10 2'
    input<<<{"description": "The mesh we want to modify"}>>> = shift_down
  []
  [injection_area]
    type = ParsedGenerateSideset<<<{"description": "A MeshGenerator that adds element sides to a sideset if the centroid of the side satisfies the `combinatorial_geometry` expression.", "href": "../meshgenerators/ParsedGenerateSideset.html"}>>>
    combinatorial_geometry<<<{"description": "Function expression encoding a combinatorial geometry"}>>> = 'x*x+y*y<1.01'
    included_subdomains<<<{"description": "A set of subdomain names or ids whose sides will be included in the new sidesets. A side is only added if the subdomain id of the corresponding element is in this set."}>>> = 1
    new_sideset_name<<<{"description": "The name of the new sideset"}>>> = 'injection_area'
    input<<<{"description": "The mesh we want to modify"}>>> = 'aquifer'
  []
  [rename]
    type = RenameBlockGenerator<<<{"description": "Changes the block IDs and/or block names for a given set of blocks defined by either block ID or block name. The changes are independent of ordering. The merging of blocks is supported.", "href": "../meshgenerators/RenameBlockGenerator.html"}>>>
    old_block<<<{"description": "Elements with these block ID(s)/name(s) will be given the new block information specified in 'new_block'"}>>> = '0 1'
    new_block<<<{"description": "The new block ID(s)/name(s) to be given by the elements defined in 'old_block'."}>>> = 'caps aquifer'
    input<<<{"description": "The mesh we want to modify"}>>> = 'injection_area'
  []
[]

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

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [porepressure]
  []
  [temperature]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 293
    scaling<<<{"description": "Specifies a scaling factor to apply to this variable"}>>> = 1E-8
  []
[]

[PorousFlowBasicTHM<<<{"href": "../../syntax/PorousFlowBasicTHM/index.html"}>>>]
  porepressure<<<{"description": "The name of the porepressure variable"}>>> = porepressure
  temperature<<<{"description": "For isothermal simulations, this is the temperature at which fluid properties (and stress-free strains) are evaluated at.  Otherwise, this is the name of the temperature variable.  Units = Kelvin"}>>> = temperature
  coupling_type<<<{"description": "The type of simulation.  For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient.  For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material"}>>> = ThermoHydro
  gravity<<<{"description": "Gravitational acceleration vector downwards (m/s^2)"}>>> = '0 0 0'
  fp<<<{"description": "The name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid"}>>> = the_simple_fluid
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [constant_injection_porepressure]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
    value<<<{"description": "Value of the BC"}>>> = 1E6
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = injection_area
  []
  [constant_injection_temperature]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
    value<<<{"description": "Value of the BC"}>>> = 313
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = injection_area
  []
[]

[FluidProperties<<<{"href": "../../syntax/FluidProperties/index.html"}>>>]
  [the_simple_fluid]
    type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../fluidproperties/SimpleFluidProperties.html"}>>>
    bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 2E9
    viscosity<<<{"description": "Constant dynamic viscosity (Pa.s)"}>>> = 1.0E-3
    density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1000.0
    thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0.0002
    cp<<<{"description": "Constant specific heat capacity at constant pressure (J/kg/K)"}>>> = 4194
    cv<<<{"description": "Constant specific heat capacity at constant volume (J/kg/K)"}>>> = 4186
    porepressure_coefficient<<<{"description": "The enthalpy is internal_energy + P / density * porepressure_coefficient.  Physically this should be 1.0, but analytic solutions are simplified when it is zero"}>>> = 0
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [porosity]
    type = PorousFlowPorosity<<<{"description": "This Material calculates the porosity PorousFlow simulations", "href": "../materials/PorousFlowPorosity.html"}>>>
    porosity_zero<<<{"description": "The porosity at zero volumetric strain and reference temperature and reference effective porepressure and reference chemistry.  This must be a real number or a constant monomial variable (not a linear lagrange or other type of variable)"}>>> = 0.1
  []
  [biot_modulus]
    type = PorousFlowConstantBiotModulus<<<{"description": "Computes the Biot Modulus, which is assumed to be constant for all time.  Sometimes 1 / BiotModulus is called storativity", "href": "../materials/PorousFlowConstantBiotModulus.html"}>>>
    biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.8
    solid_bulk_compliance<<<{"description": "Reciprocal of the drained bulk modulus of the porous skeleton.  If strain = C * stress, then solid_bulk_compliance = de_ij de_kl C_ijkl.  If the grain bulk modulus is Kg then 1/Kg = (1 - biot_coefficient) * solid_bulk_compliance."}>>> = 2E-7
    fluid_bulk_modulus<<<{"description": "Fluid bulk modulus"}>>> = 1E7
  []
  [permeability_aquifer]
    type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../materials/PorousFlowPermeabilityConst.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = aquifer
    permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '1E-14 0 0   0 1E-14 0   0 0 1E-14'
  []
  [permeability_caps]
    type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../materials/PorousFlowPermeabilityConst.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = caps
    permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '1E-15 0 0   0 1E-15 0   0 0 1E-16'
  []

  [thermal_expansion]
    type = PorousFlowConstantThermalExpansionCoefficient<<<{"description": "Computes the effective thermal expansion coefficient, (biot_coeff - porosity) * drained_coefficient + porosity * fluid_coefficient.", "href": "../materials/PorousFlowConstantThermalExpansionCoefficient.html"}>>>
    biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.8
    drained_coefficient<<<{"description": "Volumetric coefficient of thermal expansion of the drained porous skeleton (ie the porous rock without fluid, or with a fluid that is free to move in and out of the rock)"}>>> = 0.003
    fluid_coefficient<<<{"description": "Volumetric coefficient of thermal expansion for the fluid"}>>> = 0.0002
  []
  [rock_internal_energy]
    type = PorousFlowMatrixInternalEnergy<<<{"description": "This Material calculates the internal energy of solid rock grains, which is specific_heat_capacity * density * temperature.  Kernels multiply this by (1 - porosity) to find the energy density of the porous rock in a rock-fluid system", "href": "../materials/PorousFlowMatrixInternalEnergy.html"}>>>
    density<<<{"description": "Density of the rock grains"}>>> = 2500.0
    specific_heat_capacity<<<{"description": "Specific heat capacity of the rock grains (J/kg/K)."}>>> = 1200.0
  []
  [thermal_conductivity]
    type = PorousFlowThermalConductivityIdeal<<<{"description": "This Material calculates rock-fluid combined thermal conductivity by using a weighted sum.  Thermal conductivity = dry_thermal_conductivity + S^exponent * (wet_thermal_conductivity - dry_thermal_conductivity), where S is the aqueous saturation", "href": "../materials/PorousFlowThermalConductivityIdeal.html"}>>>
    dry_thermal_conductivity<<<{"description": "The thermal conductivity of the rock matrix when the aqueous saturation is zero"}>>> = '10 0 0  0 10 0  0 0 10'
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'caps aquifer'
  []
[]

[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
  active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = basic
  [basic]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    petsc_options<<<{"description": "Singleton PETSc options"}>>> = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = ' lu       mumps'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Transient
  solve_type = Newton
  end_time = 1E6
  dt = 1E5
  nl_abs_tol = 1E-10
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(moose/modules/porous_flow/examples/tutorial/03.i)

A THM example that uses PorousFlowBasicTHM is also documented in the PorousFlow tutorial with input file:

# Darcy flow with heat advection and conduction, and elasticity
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [annular]
    type = AnnularMeshGenerator<<<{"description": "For rmin>0: creates an annular mesh of QUAD4 elements. For rmin=0: creates a disc mesh of QUAD4 and TRI3 elements. Boundary sidesets are created at rmax and rmin, and given these names. If dmin!=0 and dmax!=360, a sector of an annulus or disc is created. In this case boundary sidesets are also created at dmin and dmax, and given these names", "href": "../meshgenerators/AnnularMeshGenerator.html"}>>>
    nr<<<{"description": "Number of elements in the radial direction"}>>> = 10
    rmin<<<{"description": "Inner radius.  If rmin=0 then a disc mesh (with no central hole) will be created."}>>> = 1.0
    rmax<<<{"description": "Outer radius"}>>> = 10
    growth_r<<<{"description": "The ratio of radial sizes of successive rings of elements"}>>> = 1.4
    nt<<<{"description": "Number of elements in the angular direction"}>>> = 4
    dmin<<<{"description": "Minimum degree, measured in degrees anticlockwise from x axis"}>>> = 0
    dmax<<<{"description": "Maximum angle, measured in degrees anticlockwise from x axis. If dmin=0 and dmax=360 an annular mesh is created. Otherwise, only a sector of an annulus is created"}>>> = 90
  []
  [make3D]
    type = MeshExtruderGenerator<<<{"description": "Takes a 1D or 2D mesh and extrudes the entire structure along the specified axis increasing the dimensionality of the mesh.", "href": "../meshgenerators/MeshExtruderGenerator.html"}>>>
    extrusion_vector<<<{"description": "The direction and length of the extrusion"}>>> = '0 0 12'
    num_layers<<<{"description": "The number of layers in the extruded mesh"}>>> = 3
    bottom_sideset<<<{"description": "The boundary that will be applied to the bottom of the extruded mesh"}>>> = 'bottom'
    top_sideset<<<{"description": "The boundary that will be to the top of the extruded mesh"}>>> = 'top'
    input<<<{"description": "the mesh we want to extrude"}>>> = annular
  []
  [shift_down]
    type = TransformGenerator<<<{"description": "Applies a linear transform to the entire mesh.", "href": "../meshgenerators/TransformGenerator.html"}>>>
    transform<<<{"description": "The type of transformation to perform (TRANSLATE, TRANSLATE_CENTER_ORIGIN, TRANSLATE_MIN_ORIGIN, ROTATE, SCALE)"}>>> = TRANSLATE
    vector_value<<<{"description": "The value to use for the transformation. When using TRANSLATE or SCALE, the xyz coordinates are applied in each direction respectively. When using ROTATE, the values are interpreted as the Euler angles phi, theta and psi given in degrees."}>>> = '0 0 -6'
    input<<<{"description": "The mesh we want to modify"}>>> = make3D
  []
  [aquifer]
    type = SubdomainBoundingBoxGenerator<<<{"description": "Changes the subdomain ID of elements either (XOR) inside or outside the specified box to the specified ID.", "href": "../meshgenerators/SubdomainBoundingBoxGenerator.html"}>>>
    block_id<<<{"description": "Subdomain id to set for inside/outside the bounding box"}>>> = 1
    bottom_left<<<{"description": "The bottom left point (in x,y,z with spaces in-between)."}>>> = '0 0 -2'
    top_right<<<{"description": "The bottom left point (in x,y,z with spaces in-between)."}>>> = '10 10 2'
    input<<<{"description": "The mesh we want to modify"}>>> = shift_down
  []
  [injection_area]
    type = ParsedGenerateSideset<<<{"description": "A MeshGenerator that adds element sides to a sideset if the centroid of the side satisfies the `combinatorial_geometry` expression.", "href": "../meshgenerators/ParsedGenerateSideset.html"}>>>
    combinatorial_geometry<<<{"description": "Function expression encoding a combinatorial geometry"}>>> = 'x*x+y*y<1.01'
    included_subdomains<<<{"description": "A set of subdomain names or ids whose sides will be included in the new sidesets. A side is only added if the subdomain id of the corresponding element is in this set."}>>> = 1
    new_sideset_name<<<{"description": "The name of the new sideset"}>>> = 'injection_area'
    input<<<{"description": "The mesh we want to modify"}>>> = 'aquifer'
  []
  [rename]
    type = RenameBlockGenerator<<<{"description": "Changes the block IDs and/or block names for a given set of blocks defined by either block ID or block name. The changes are independent of ordering. The merging of blocks is supported.", "href": "../meshgenerators/RenameBlockGenerator.html"}>>>
    old_block<<<{"description": "Elements with these block ID(s)/name(s) will be given the new block information specified in 'new_block'"}>>> = '0 1'
    new_block<<<{"description": "The new block ID(s)/name(s) to be given by the elements defined in 'old_block'."}>>> = 'caps aquifer'
    input<<<{"description": "The mesh we want to modify"}>>> = 'injection_area'
  []
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  displacements = 'disp_x disp_y disp_z'
  PorousFlowDictator = dictator
  biot_coefficient = 1.0
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [porepressure]
  []
  [temperature]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 293
    scaling<<<{"description": "Specifies a scaling factor to apply to this variable"}>>> = 1E-8
  []
  [disp_x]
    scaling<<<{"description": "Specifies a scaling factor to apply to this variable"}>>> = 1E-10
  []
  [disp_y]
    scaling<<<{"description": "Specifies a scaling factor to apply to this variable"}>>> = 1E-10
  []
  [disp_z]
    scaling<<<{"description": "Specifies a scaling factor to apply to this variable"}>>> = 1E-10
  []
[]

[PorousFlowBasicTHM<<<{"href": "../../syntax/PorousFlowBasicTHM/index.html"}>>>]
  porepressure<<<{"description": "The name of the porepressure variable"}>>> = porepressure
  temperature<<<{"description": "For isothermal simulations, this is the temperature at which fluid properties (and stress-free strains) are evaluated at.  Otherwise, this is the name of the temperature variable.  Units = Kelvin"}>>> = temperature
  coupling_type<<<{"description": "The type of simulation.  For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient.  For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material"}>>> = ThermoHydroMechanical
  gravity<<<{"description": "Gravitational acceleration vector downwards (m/s^2)"}>>> = '0 0 0'
  fp<<<{"description": "The name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid"}>>> = the_simple_fluid
  eigenstrain_names<<<{"description": "List of all eigenstrain models used in mechanics calculations. Typically the eigenstrain_name used in ComputeThermalExpansionEigenstrain.  Only needed for thermally-coupled simulations with thermal expansion."}>>> = thermal_contribution
  use_displaced_mesh<<<{"description": "Use displaced mesh computations in mechanical kernels"}>>> = false
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [constant_injection_porepressure]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = porepressure
    value<<<{"description": "Value of the BC"}>>> = 1E6
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = injection_area
  []
  [constant_injection_temperature]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
    value<<<{"description": "Value of the BC"}>>> = 313
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = injection_area
  []

  [roller_tmax]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = dmax
  []
  [roller_tmin]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = dmin
  []
  [roller_top_bottom]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
  []
  [cavity_pressure_x]
    type = Pressure<<<{"description": "Applies a pressure on a given boundary in a given direction", "href": "../bcs/Pressure.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = injection_area
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
    component<<<{"description": "The component for the pressure"}>>> = 0
    factor<<<{"description": "The magnitude to use in computing the pressure"}>>> = 1E6
    use_displaced_mesh<<<{"description": "Whether to use the displaced mesh."}>>> = false
  []
  [cavity_pressure_y]
    type = Pressure<<<{"description": "Applies a pressure on a given boundary in a given direction", "href": "../bcs/Pressure.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = injection_area
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
    component<<<{"description": "The component for the pressure"}>>> = 1
    factor<<<{"description": "The magnitude to use in computing the pressure"}>>> = 1E6
    use_displaced_mesh<<<{"description": "Whether to use the displaced mesh."}>>> = false
  []
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  [stress_rr]
    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
  []
  [stress_pp]
    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"}>>>]
  [stress_rr]
    type = RankTwoScalarAux<<<{"description": "Compute a scalar property of a RankTwoTensor", "href": "../auxscalarkernels/RankTwoScalarAux.html"}>>>
    rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_rr
    scalar_type<<<{"description": "Type of scalar output"}>>> = RadialStress
    point1<<<{"description": "Start point for axis used to calculate some cylindrical material tensor quantities"}>>> = '0 0 0'
    point2<<<{"description": "End point for axis used to calculate some material tensor quantities"}>>> = '0 0 1'
  []
  [stress_pp]
    type = RankTwoScalarAux<<<{"description": "Compute a scalar property of a RankTwoTensor", "href": "../auxscalarkernels/RankTwoScalarAux.html"}>>>
    rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_pp
    scalar_type<<<{"description": "Type of scalar output"}>>> = HoopStress
    point1<<<{"description": "Start point for axis used to calculate some cylindrical material tensor quantities"}>>> = '0 0 0'
    point2<<<{"description": "End point for axis used to calculate some material tensor quantities"}>>> = '0 0 1'
  []
[]

[FluidProperties<<<{"href": "../../syntax/FluidProperties/index.html"}>>>]
  [the_simple_fluid]
    type = SimpleFluidProperties<<<{"description": "Fluid properties for a simple fluid with a constant bulk density", "href": "../fluidproperties/SimpleFluidProperties.html"}>>>
    bulk_modulus<<<{"description": "Constant bulk modulus (Pa)"}>>> = 2E9
    viscosity<<<{"description": "Constant dynamic viscosity (Pa.s)"}>>> = 1.0E-3
    density0<<<{"description": "Density at zero pressure and zero temperature"}>>> = 1000.0
    thermal_expansion<<<{"description": "Constant coefficient of thermal expansion (1/K)"}>>> = 0.0002
    cp<<<{"description": "Constant specific heat capacity at constant pressure (J/kg/K)"}>>> = 4194
    cv<<<{"description": "Constant specific heat capacity at constant volume (J/kg/K)"}>>> = 4186
    porepressure_coefficient<<<{"description": "The enthalpy is internal_energy + P / density * porepressure_coefficient.  Physically this should be 1.0, but analytic solutions are simplified when it is zero"}>>> = 0
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [porosity]
    type = PorousFlowPorosity<<<{"description": "This Material calculates the porosity PorousFlow simulations", "href": "../materials/PorousFlowPorosity.html"}>>>
    porosity_zero<<<{"description": "The porosity at zero volumetric strain and reference temperature and reference effective porepressure and reference chemistry.  This must be a real number or a constant monomial variable (not a linear lagrange or other type of variable)"}>>> = 0.1
  []
  [biot_modulus]
    type = PorousFlowConstantBiotModulus<<<{"description": "Computes the Biot Modulus, which is assumed to be constant for all time.  Sometimes 1 / BiotModulus is called storativity", "href": "../materials/PorousFlowConstantBiotModulus.html"}>>>
    solid_bulk_compliance<<<{"description": "Reciprocal of the drained bulk modulus of the porous skeleton.  If strain = C * stress, then solid_bulk_compliance = de_ij de_kl C_ijkl.  If the grain bulk modulus is Kg then 1/Kg = (1 - biot_coefficient) * solid_bulk_compliance."}>>> = 2E-7
    fluid_bulk_modulus<<<{"description": "Fluid bulk modulus"}>>> = 1E7
  []
  [permeability_aquifer]
    type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../materials/PorousFlowPermeabilityConst.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = aquifer
    permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '1E-14 0 0   0 1E-14 0   0 0 1E-14'
  []
  [permeability_caps]
    type = PorousFlowPermeabilityConst<<<{"description": "This Material calculates the permeability tensor assuming it is constant", "href": "../materials/PorousFlowPermeabilityConst.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = caps
    permeability<<<{"description": "The permeability tensor (usually in m^2), which is assumed constant for this material"}>>> = '1E-15 0 0   0 1E-15 0   0 0 1E-16'
  []

  [thermal_expansion]
    type = PorousFlowConstantThermalExpansionCoefficient<<<{"description": "Computes the effective thermal expansion coefficient, (biot_coeff - porosity) * drained_coefficient + porosity * fluid_coefficient.", "href": "../materials/PorousFlowConstantThermalExpansionCoefficient.html"}>>>
    drained_coefficient<<<{"description": "Volumetric coefficient of thermal expansion of the drained porous skeleton (ie the porous rock without fluid, or with a fluid that is free to move in and out of the rock)"}>>> = 0.003
    fluid_coefficient<<<{"description": "Volumetric coefficient of thermal expansion for the fluid"}>>> = 0.0002
  []
  [rock_internal_energy]
    type = PorousFlowMatrixInternalEnergy<<<{"description": "This Material calculates the internal energy of solid rock grains, which is specific_heat_capacity * density * temperature.  Kernels multiply this by (1 - porosity) to find the energy density of the porous rock in a rock-fluid system", "href": "../materials/PorousFlowMatrixInternalEnergy.html"}>>>
    density<<<{"description": "Density of the rock grains"}>>> = 2500.0
    specific_heat_capacity<<<{"description": "Specific heat capacity of the rock grains (J/kg/K)."}>>> = 1200.0
  []
  [thermal_conductivity]
    type = PorousFlowThermalConductivityIdeal<<<{"description": "This Material calculates rock-fluid combined thermal conductivity by using a weighted sum.  Thermal conductivity = dry_thermal_conductivity + S^exponent * (wet_thermal_conductivity - dry_thermal_conductivity), where S is the aqueous saturation", "href": "../materials/PorousFlowThermalConductivityIdeal.html"}>>>
    dry_thermal_conductivity<<<{"description": "The thermal conductivity of the rock matrix when the aqueous saturation is zero"}>>> = '10 0 0  0 10 0  0 0 10'
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'caps aquifer'
  []

  [elasticity_tensor]
    type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../materials/ComputeIsotropicElasticityTensor.html"}>>>
    youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 5E9
    poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.0
  []
  [strain]
    type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../materials/ComputeSmallStrain.html"}>>>
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = thermal_contribution
  []
  [thermal_contribution]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../materials/ComputeThermalExpansionEigenstrain.html"}>>>
    temperature<<<{"description": "Coupled temperature"}>>> = temperature
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 0.001 # this is the linear thermal expansion coefficient
    eigenstrain_name<<<{"description": "Material property name for the eigenstrain tensor computed by this model. IMPORTANT: The name of this property must also be provided to the strain calculator."}>>> = thermal_contribution
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 293
  []
  [stress]
    type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../materials/ComputeLinearElasticStress.html"}>>>
  []
[]

[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
  active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = basic
  [basic]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    petsc_options<<<{"description": "Singleton PETSc options"}>>> = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = ' lu       mumps'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Transient
  solve_type = Newton
  end_time = 1E6
  dt = 1E5
  nl_abs_tol = 1E-15
  nl_rel_tol = 1E-14
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(moose/modules/porous_flow/examples/tutorial/04.i)

Input Parameters

  • porepressureThe name of the porepressure variable

    C++ Type:VariableName

    Unit:(no unit assumed)

    Controllable:No

    Description:The name of the porepressure variable

Required Parameters

  • active__all__ If specified only the blocks named will be visited and made active

    Default:__all__

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:If specified only the blocks named will be visited and made active

  • add_darcy_auxTrueAdd AuxVariables that record Darcy velocity

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Add AuxVariables that record Darcy velocity

  • add_stress_auxTrueAdd AuxVariables that record effective stress

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Add AuxVariables that record effective stress

  • base_nameThe base_name used in the TensorMechanics strain calculator. This is only relevant for mechanically-coupled models.

    C++ Type:std::string

    Controllable:No

    Description:The base_name used in the TensorMechanics strain calculator. This is only relevant for mechanically-coupled models.

  • biot_coefficient1The Biot coefficient (relevant only for mechanically-coupled simulations)

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The Biot coefficient (relevant only for mechanically-coupled simulations)

  • coupling_typeHydroThe type of simulation. For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient. For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material

    Default:Hydro

    C++ Type:MooseEnum

    Options:Hydro, ThermoHydro, HydroMechanical, ThermoHydroMechanical

    Controllable:No

    Description:The type of simulation. For simulations involving Mechanical deformations, you will need to supply the correct Biot coefficient. For simulations involving Thermal flows, you will need an associated ConstantThermalExpansionCoefficient Material

  • dictator_namedictatorThe name of the dictator user object that is created by this Action

    Default:dictator

    C++ Type:std::string

    Controllable:No

    Description:The name of the dictator user object that is created by this Action

  • displacementsThe name of the displacement variables (relevant only for mechanically-coupled simulations)

    C++ Type:std::vector<VariableName>

    Unit:(no unit assumed)

    Controllable:No

    Description:The name of the displacement variables (relevant only for mechanically-coupled simulations)

  • eigenstrain_namesList of all eigenstrain models used in mechanics calculations. Typically the eigenstrain_name used in ComputeThermalExpansionEigenstrain. Only needed for thermally-coupled simulations with thermal expansion.

    C++ Type:std::vector<MaterialPropertyName>

    Unit:(no unit assumed)

    Controllable:No

    Description:List of all eigenstrain models used in mechanics calculations. Typically the eigenstrain_name used in ComputeThermalExpansionEigenstrain. Only needed for thermally-coupled simulations with thermal expansion.

  • fluid_properties_typePorousFlowSingleComponentFluidType of fluid properties to use. For 'PorousFlowSingleComponentFluid' you must provide a fp UserObject. For 'PorousFlowBrine' you must supply a nacl_name. For 'Custom' your input file must include a Material that provides fluid properties such as density, viscosity, enthalpy and internal energy

    Default:PorousFlowSingleComponentFluid

    C++ Type:MooseEnum

    Options:PorousFlowSingleComponentFluid, PorousFlowBrine, Custom

    Controllable:No

    Description:Type of fluid properties to use. For 'PorousFlowSingleComponentFluid' you must provide a fp UserObject. For 'PorousFlowBrine' you must supply a nacl_name. For 'Custom' your input file must include a Material that provides fluid properties such as density, viscosity, enthalpy and internal energy

  • flux_limiter_typeVanLeerType of flux limiter to use if stabilization=KT. 'None' means that no antidiffusion will be added in the Kuzmin-Turek scheme

    Default:VanLeer

    C++ Type:MooseEnum

    Options:MinMod, VanLeer, MC, superbee, None

    Controllable:No

    Description:Type of flux limiter to use if stabilization=KT. 'None' means that no antidiffusion will be added in the Kuzmin-Turek scheme

  • fpThe name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid

    C++ Type:UserObjectName

    Controllable:No

    Description:The name of the user object for fluid properties. Only needed if fluid_properties_type = PorousFlowSingleComponentFluid

  • gravity0 0 -10Gravitational acceleration vector downwards (m/s^2)

    Default:0 0 -10

    C++ Type:libMesh::VectorValue<double>

    Unit:(no unit assumed)

    Controllable:No

    Description:Gravitational acceleration vector downwards (m/s^2)

  • inactiveIf specified blocks matching these identifiers will be skipped.

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:If specified blocks matching these identifiers will be skipped.

  • mass_fraction_varsList of variables that represent the mass fractions. With only one fluid component, this may be left empty. With N fluid components, the format is 'f_0 f_1 f_2 ... f_(N-1)'. That is, the N^th component need not be specified because f_N = 1 - (f_0 + f_1 + ... + f_(N-1)). It is best numerically to choose the N-1 mass fraction variables so that they represent the fluid components with small concentrations. This Action will associated the i^th mass fraction variable to the equation for the i^th fluid component, and the pressure variable to the N^th fluid component.

    C++ Type:std::vector<VariableName>

    Unit:(no unit assumed)

    Controllable:No

    Description:List of variables that represent the mass fractions. With only one fluid component, this may be left empty. With N fluid components, the format is 'f_0 f_1 f_2 ... f_(N-1)'. That is, the N^th component need not be specified because f_N = 1 - (f_0 + f_1 + ... + f_(N-1)). It is best numerically to choose the N-1 mass fraction variables so that they represent the fluid components with small concentrations. This Action will associated the i^th mass fraction variable to the equation for the i^th fluid component, and the pressure variable to the N^th fluid component.

  • multiply_by_densityFalseIf true, then the Kernels for fluid flow are multiplied by the fluid density. If false, this multiplication is not performed, which means the problem linearises, but that care must be taken when using other PorousFlow objects.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:If true, then the Kernels for fluid flow are multiplied by the fluid density. If false, this multiplication is not performed, which means the problem linearises, but that care must be taken when using other PorousFlow objects.

  • nacl_nameName of the NaCl variable. Only required if fluid_properties_type = PorousFlowBrine

    C++ Type:VariableName

    Unit:(no unit assumed)

    Controllable:No

    Description:Name of the NaCl variable. Only required if fluid_properties_type = PorousFlowBrine

  • number_aqueous_equilibrium0The number of secondary species in the aqueous-equilibrium reaction system. (Leave as zero if the simulation does not involve chemistry)

    Default:0

    C++ Type:unsigned int

    Controllable:No

    Description:The number of secondary species in the aqueous-equilibrium reaction system. (Leave as zero if the simulation does not involve chemistry)

  • number_aqueous_kinetic0The number of secondary species in the aqueous-kinetic reaction system involved in precipitation and dissolution. (Leave as zero if the simulation does not involve chemistry)

    Default:0

    C++ Type:unsigned int

    Controllable:No

    Description:The number of secondary species in the aqueous-kinetic reaction system involved in precipitation and dissolution. (Leave as zero if the simulation does not involve chemistry)

  • pressure_unitPaThe unit of the pressure variable used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid

    Default:Pa

    C++ Type:MooseEnum

    Options:Pa, MPa

    Controllable:No

    Description:The unit of the pressure variable used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid

  • save_component_rate_inList of AuxVariables into which the rate-of-change of each fluid component at each node will be saved. There must be exactly N of these to match the N fluid components. The result will be measured in kg/s, where the kg is the mass of the fluid component at the node (or m^3/s if multiply_by_density=false). Note that this saves the result from the MassTimeDerivative Kernels, but NOT from the MassVolumetricExpansion Kernels.

    C++ Type:std::vector<AuxVariableName>

    Unit:(no unit assumed)

    Controllable:No

    Description:List of AuxVariables into which the rate-of-change of each fluid component at each node will be saved. There must be exactly N of these to match the N fluid components. The result will be measured in kg/s, where the kg is the mass of the fluid component at the node (or m^3/s if multiply_by_density=false). Note that this saves the result from the MassTimeDerivative Kernels, but NOT from the MassVolumetricExpansion Kernels.

  • stabilizationFullNumerical stabilization used. 'Full' means full upwinding. 'KT' means FEM-TVD stabilization of Kuzmin-Turek

    Default:Full

    C++ Type:MooseEnum

    Options:None, Full, KT

    Controllable:No

    Description:Numerical stabilization used. 'Full' means full upwinding. 'KT' means FEM-TVD stabilization of Kuzmin-Turek

  • strain_at_nearest_qpFalseOnly relevant for models in which porosity depends on strain. If true, then when calculating nodal porosity that depends on strain, the strain at the nearest quadpoint will be used. This adds a small extra computational burden, and is only necessary for simulations involving: (1) elements that are not linear lagrange or (2) certain PorousFlow Dirac Kernels (as specified in their documentation). If you set this to true, you will also want to set the same parameter to true for related Kernels and Materials (which is probably easiest to do in the GlobalParams block)

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Only relevant for models in which porosity depends on strain. If true, then when calculating nodal porosity that depends on strain, the strain at the nearest quadpoint will be used. This adds a small extra computational burden, and is only necessary for simulations involving: (1) elements that are not linear lagrange or (2) certain PorousFlow Dirac Kernels (as specified in their documentation). If you set this to true, you will also want to set the same parameter to true for related Kernels and Materials (which is probably easiest to do in the GlobalParams block)

  • temperature293.0For isothermal simulations, this is the temperature at which fluid properties (and stress-free strains) are evaluated at. Otherwise, this is the name of the temperature variable. Units = Kelvin

    Default:293.0

    C++ Type:std::vector<VariableName>

    Unit:(no unit assumed)

    Controllable:No

    Description:For isothermal simulations, this is the temperature at which fluid properties (and stress-free strains) are evaluated at. Otherwise, this is the name of the temperature variable. Units = Kelvin

  • temperature_unitKelvinThe unit of the temperature variable

    Default:Kelvin

    C++ Type:MooseEnum

    Options:Kelvin, Celsius

    Controllable:No

    Description:The unit of the temperature variable

  • time_unitsecondsThe unit of time used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid

    Default:seconds

    C++ Type:MooseEnum

    Options:seconds, hours, days, years

    Controllable:No

    Description:The unit of time used everywhere in the input file except for in the FluidProperties-module objects. This can be set to the non-default value only for fluid_properties_type = PorousFlowSingleComponentFluid

  • use_displaced_meshFalseUse displaced mesh computations in mechanical kernels

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Use displaced mesh computations in mechanical kernels

Optional Parameters

  • blockThe list of block ids (SubdomainID) on which the porous flow is defined on

    C++ Type:std::vector<SubdomainName>

    Controllable:No

    Description:The list of block ids (SubdomainID) on which the porous flow is defined on

  • 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.

Advanced Parameters