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]
  [annular]
    type = AnnularMeshGenerator
    nr = 10
    rmin = 1.0
    rmax = 10
    growth_r = 1.4
    nt = 4
    dmin = 0
    dmax = 90
  []
  [make3D]
    type = MeshExtruderGenerator
    extrusion_vector = '0 0 12'
    num_layers = 3
    bottom_sideset = 'bottom'
    top_sideset = 'top'
    input = annular
  []
  [shift_down]
    type = TransformGenerator
    transform = TRANSLATE
    vector_value = '0 0 -6'
    input = make3D
  []
  [aquifer]
    type = SubdomainBoundingBoxGenerator
    block_id = 1
    bottom_left = '0 0 -2'
    top_right = '10 10 2'
    input = shift_down
  []
  [injection_area]
    type = ParsedGenerateSideset
    combinatorial_geometry = 'x*x+y*y<1.01'
    included_subdomains = 1
    new_sideset_name = 'injection_area'
    input = 'aquifer'
  []
  [rename]
    type = RenameBlockGenerator
    old_block = '0 1'
    new_block = 'caps aquifer'
    input = 'injection_area'
  []
[]

[GlobalParams]
  PorousFlowDictator = dictator
[]

[Variables]
  [porepressure]
  []
[]

[PorousFlowBasicTHM]
  porepressure = porepressure
  coupling_type = Hydro
  gravity = '0 0 0'
  fp = the_simple_fluid
[]

[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1E6
    boundary = injection_area
  []
[]

[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    viscosity = 1.0E-3
    density0 = 1000.0
  []
[]

[Materials]
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [biot_modulus]
    type = PorousFlowConstantBiotModulus
    biot_coefficient = 0.8
    solid_bulk_compliance = 2E-7
    fluid_bulk_modulus = 1E7
  []
  [permeability_aquifer]
    type = PorousFlowPermeabilityConst
    block = aquifer
    permeability = '1E-14 0 0   0 1E-14 0   0 0 1E-14'
  []
  [permeability_caps]
    type = PorousFlowPermeabilityConst
    block = caps
    permeability = '1E-15 0 0   0 1E-15 0   0 0 1E-16'
  []
[]

[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]

[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 1E6
  dt = 1E5
  nl_abs_tol = 1E-13
[]

[Outputs]
  exodus = true
[]
(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]
  [annular]
    type = AnnularMeshGenerator
    nr = 10
    rmin = 1.0
    rmax = 10
    growth_r = 1.4
    nt = 4
    dmin = 0
    dmax = 90
  []
  [make3D]
    type = MeshExtruderGenerator
    extrusion_vector = '0 0 12'
    num_layers = 3
    bottom_sideset = 'bottom'
    top_sideset = 'top'
    input = annular
  []
  [shift_down]
    type = TransformGenerator
    transform = TRANSLATE
    vector_value = '0 0 -6'
    input = make3D
  []
  [aquifer]
    type = SubdomainBoundingBoxGenerator
    block_id = 1
    bottom_left = '0 0 -2'
    top_right = '10 10 2'
    input = shift_down
  []
  [injection_area]
    type = ParsedGenerateSideset
    combinatorial_geometry = 'x*x+y*y<1.01'
    included_subdomains = 1
    new_sideset_name = 'injection_area'
    input = 'aquifer'
  []
  [rename]
    type = RenameBlockGenerator
    old_block = '0 1'
    new_block = 'caps aquifer'
    input = 'injection_area'
  []
[]

[GlobalParams]
  PorousFlowDictator = dictator
[]

[Variables]
  [porepressure]
  []
  [temperature]
    initial_condition = 293
    scaling = 1E-8
  []
[]

[PorousFlowBasicTHM]
  porepressure = porepressure
  temperature = temperature
  coupling_type = ThermoHydro
  gravity = '0 0 0'
  fp = the_simple_fluid
[]

[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1E6
    boundary = injection_area
  []
  [constant_injection_temperature]
    type = DirichletBC
    variable = temperature
    value = 313
    boundary = injection_area
  []
[]

[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    viscosity = 1.0E-3
    density0 = 1000.0
    thermal_expansion = 0.0002
    cp = 4194
    cv = 4186
    porepressure_coefficient = 0
  []
[]

[Materials]
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [biot_modulus]
    type = PorousFlowConstantBiotModulus
    biot_coefficient = 0.8
    solid_bulk_compliance = 2E-7
    fluid_bulk_modulus = 1E7
  []
  [permeability_aquifer]
    type = PorousFlowPermeabilityConst
    block = aquifer
    permeability = '1E-14 0 0   0 1E-14 0   0 0 1E-14'
  []
  [permeability_caps]
    type = PorousFlowPermeabilityConst
    block = caps
    permeability = '1E-15 0 0   0 1E-15 0   0 0 1E-16'
  []

  [thermal_expansion]
    type = PorousFlowConstantThermalExpansionCoefficient
    biot_coefficient = 0.8
    drained_coefficient = 0.003
    fluid_coefficient = 0.0002
  []
  [rock_internal_energy]
    type = PorousFlowMatrixInternalEnergy
    density = 2500.0
    specific_heat_capacity = 1200.0
  []
  [thermal_conductivity]
    type = PorousFlowThermalConductivityIdeal
    dry_thermal_conductivity = '10 0 0  0 10 0  0 0 10'
    block = 'caps aquifer'
  []
[]

[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]

[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 1E6
  dt = 1E5
  nl_abs_tol = 1E-10
[]

[Outputs]
  exodus = true
[]
(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]
  [annular]
    type = AnnularMeshGenerator
    nr = 10
    rmin = 1.0
    rmax = 10
    growth_r = 1.4
    nt = 4
    dmin = 0
    dmax = 90
  []
  [make3D]
    type = MeshExtruderGenerator
    extrusion_vector = '0 0 12'
    num_layers = 3
    bottom_sideset = 'bottom'
    top_sideset = 'top'
    input = annular
  []
  [shift_down]
    type = TransformGenerator
    transform = TRANSLATE
    vector_value = '0 0 -6'
    input = make3D
  []
  [aquifer]
    type = SubdomainBoundingBoxGenerator
    block_id = 1
    bottom_left = '0 0 -2'
    top_right = '10 10 2'
    input = shift_down
  []
  [injection_area]
    type = ParsedGenerateSideset
    combinatorial_geometry = 'x*x+y*y<1.01'
    included_subdomains = 1
    new_sideset_name = 'injection_area'
    input = 'aquifer'
  []
  [rename]
    type = RenameBlockGenerator
    old_block = '0 1'
    new_block = 'caps aquifer'
    input = 'injection_area'
  []
[]

[GlobalParams]
  displacements = 'disp_x disp_y disp_z'
  PorousFlowDictator = dictator
  biot_coefficient = 1.0
[]

[Variables]
  [porepressure]
  []
  [temperature]
    initial_condition = 293
    scaling = 1E-8
  []
  [disp_x]
    scaling = 1E-10
  []
  [disp_y]
    scaling = 1E-10
  []
  [disp_z]
    scaling = 1E-10
  []
[]

[PorousFlowBasicTHM]
  porepressure = porepressure
  temperature = temperature
  coupling_type = ThermoHydroMechanical
  gravity = '0 0 0'
  fp = the_simple_fluid
  eigenstrain_names = thermal_contribution
  use_displaced_mesh = false
[]

[BCs]
  [constant_injection_porepressure]
    type = DirichletBC
    variable = porepressure
    value = 1E6
    boundary = injection_area
  []
  [constant_injection_temperature]
    type = DirichletBC
    variable = temperature
    value = 313
    boundary = injection_area
  []

  [roller_tmax]
    type = DirichletBC
    variable = disp_x
    value = 0
    boundary = dmax
  []
  [roller_tmin]
    type = DirichletBC
    variable = disp_y
    value = 0
    boundary = dmin
  []
  [roller_top_bottom]
    type = DirichletBC
    variable = disp_z
    value = 0
    boundary = 'top bottom'
  []
  [cavity_pressure_x]
    type = Pressure
    boundary = injection_area
    variable = disp_x
    component = 0
    factor = 1E6
    use_displaced_mesh = false
  []
  [cavity_pressure_y]
    type = Pressure
    boundary = injection_area
    variable = disp_y
    component = 1
    factor = 1E6
    use_displaced_mesh = false
  []
[]

[AuxVariables]
  [stress_rr]
    family = MONOMIAL
    order = CONSTANT
  []
  [stress_pp]
    family = MONOMIAL
    order = CONSTANT
  []
[]

[AuxKernels]
  [stress_rr]
    type = RankTwoScalarAux
    rank_two_tensor = stress
    variable = stress_rr
    scalar_type = RadialStress
    point1 = '0 0 0'
    point2 = '0 0 1'
  []
  [stress_pp]
    type = RankTwoScalarAux
    rank_two_tensor = stress
    variable = stress_pp
    scalar_type = HoopStress
    point1 = '0 0 0'
    point2 = '0 0 1'
  []
[]

[FluidProperties]
  [the_simple_fluid]
    type = SimpleFluidProperties
    bulk_modulus = 2E9
    viscosity = 1.0E-3
    density0 = 1000.0
    thermal_expansion = 0.0002
    cp = 4194
    cv = 4186
    porepressure_coefficient = 0
  []
[]

[Materials]
  [porosity]
    type = PorousFlowPorosity
    porosity_zero = 0.1
  []
  [biot_modulus]
    type = PorousFlowConstantBiotModulus
    solid_bulk_compliance = 2E-7
    fluid_bulk_modulus = 1E7
  []
  [permeability_aquifer]
    type = PorousFlowPermeabilityConst
    block = aquifer
    permeability = '1E-14 0 0   0 1E-14 0   0 0 1E-14'
  []
  [permeability_caps]
    type = PorousFlowPermeabilityConst
    block = caps
    permeability = '1E-15 0 0   0 1E-15 0   0 0 1E-16'
  []

  [thermal_expansion]
    type = PorousFlowConstantThermalExpansionCoefficient
    drained_coefficient = 0.003
    fluid_coefficient = 0.0002
  []
  [rock_internal_energy]
    type = PorousFlowMatrixInternalEnergy
    density = 2500.0
    specific_heat_capacity = 1200.0
  []
  [thermal_conductivity]
    type = PorousFlowThermalConductivityIdeal
    dry_thermal_conductivity = '10 0 0  0 10 0  0 0 10'
    block = 'caps aquifer'
  []

  [elasticity_tensor]
    type = ComputeIsotropicElasticityTensor
    youngs_modulus = 5E9
    poissons_ratio = 0.0
  []
  [strain]
    type = ComputeSmallStrain
    eigenstrain_names = thermal_contribution
  []
  [thermal_contribution]
    type = ComputeThermalExpansionEigenstrain
    temperature = temperature
    thermal_expansion_coeff = 0.001 # this is the linear thermal expansion coefficient
    eigenstrain_name = thermal_contribution
    stress_free_temperature = 293
  []
  [stress]
    type = ComputeLinearElasticStress
  []
[]

[Preconditioning]
  active = basic
  [basic]
    type = SMP
    full = true
    petsc_options = '-ksp_diagonal_scale -ksp_diagonal_scale_fix'
    petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
    petsc_options_value = ' asm      lu           NONZERO                   2'
  []
  [preferred_but_might_not_be_installed]
    type = SMP
    full = true
    petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
    petsc_options_value = ' lu       mumps'
  []
[]

[Executioner]
  type = Transient
  solve_type = Newton
  end_time = 1E6
  dt = 1E5
  nl_abs_tol = 1E-15
  nl_rel_tol = 1E-14
[]

[Outputs]
  exodus = true
[]
(modules/porous_flow/examples/tutorial/04.i)

Input Parameters

  • porepressureThe name of the porepressure variable

    C++ Type:VariableName

    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

    Controllable:No

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

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

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

    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>

    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>

    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>

    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

    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>

    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>

    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