- PorousFlowDictatorThe UserObject that holds the list of PorousFlow variable names
C++ Type:UserObjectName
Controllable:No
Description:The UserObject that holds the list of PorousFlow variable names
- m(1-porosity) exponent
C++ Type:double
Controllable:No
Description:(1-porosity) exponent
- nPorosity exponent
C++ Type:double
Controllable:No
Description:Porosity exponent
PorousFlowPermeabilityKozenyCarman
This Material calculates the permeability tensor from a form of the Kozeny-Carman equation, k = k_ijk * A * phin / (1 - phi)m, where k_ijk is a tensor providing the anisotropy, phi is porosity, n and m are positive scalar constants and A is given in one of the following forms: A = k0 * (1 - phi0)^m / phi0n (where k0 and phi0 are a reference permeability and porosity) or A = f * d2 (where f is a scalar constant and d is grain diameter.
Permeability is calculated from porosity using where and are user-defined constants.
Input can be entered in one of two forms depending on the value of poroperm_function
poroperm_function | input format |
---|---|
kozeny_carman_fd2 | |
kozeny_carman_fd2 |
The parameters and are then converted to the correct form internally.
Input Parameters
- at_nodesFalseEvaluate Material properties at nodes instead of quadpoints
Default:False
C++ Type:bool
Controllable:No
Description:Evaluate Material properties at nodes instead of quadpoints
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- computeTrueWhen false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
Default:True
C++ Type:bool
Controllable:No
Description:When false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
Default:NONE
C++ Type:MooseEnum
Options:NONE, ELEMENT, SUBDOMAIN
Controllable:No
Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
- dThe grain diameter, required for kozeny_carman_fd2
C++ Type:double
Controllable:No
Description:The grain diameter, required for kozeny_carman_fd2
- declare_suffixAn optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
- fThe multiplying factor, required for kozeny_carman_fd2
C++ Type:double
Controllable:No
Description:The multiplying factor, required for kozeny_carman_fd2
- k0The permeability scalar value (usually in m^2) at the reference porosity, required for kozeny_carman_phi0
C++ Type:double
Controllable:No
Description:The permeability scalar value (usually in m^2) at the reference porosity, required for kozeny_carman_phi0
- k_anisotropyA tensor to multiply the calculated scalar permeability, in order to obtain anisotropy if required. Defaults to isotropic permeability if not specified.
C++ Type:libMesh::TensorValue<double>
Controllable:No
Description:A tensor to multiply the calculated scalar permeability, in order to obtain anisotropy if required. Defaults to isotropic permeability if not specified.
- phi0The reference porosity, required for kozeny_carman_phi0
C++ Type:double
Controllable:No
Description:The reference porosity, required for kozeny_carman_phi0
- poroperm_functionkozeny_carman_fd2Function relating porosity and permeability. The options are: kozeny_carman_fd2 = f d^2 phi^n/(1-phi)^m (where phi is porosity, f is a scalar constant with typical values 0.01-0.001, and d is grain size). kozeny_carman_phi0 = k0 (1-phi0)^m/phi0^n * phi^n/(1-phi)^m (where phi is porosity, and k0 is the permeability at porosity phi0)
Default:kozeny_carman_fd2
C++ Type:MooseEnum
Options:kozeny_carman_fd2, kozeny_carman_phi0
Controllable:No
Description:Function relating porosity and permeability. The options are: kozeny_carman_fd2 = f d^2 phi^n/(1-phi)^m (where phi is porosity, f is a scalar constant with typical values 0.01-0.001, and d is grain size). kozeny_carman_phi0 = k0 (1-phi0)^m/phi0^n * phi^n/(1-phi)^m (where phi is porosity, and k0 is the permeability at porosity phi0)
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)
C++ Type:std::vector<std::string>
Controllable:No
Description:List of material properties, from this material, to output (outputs must also be defined to an output type)
- outputsnone Vector of output names where you would like to restrict the output of variables(s) associated with this object
Default:none
C++ Type:std::vector<OutputName>
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
Input Files
- (modules/porous_flow/examples/tutorial/07.i)
- (modules/porous_flow/test/tests/heterogeneous_materials/vol_expansion_poroperm.i)
- (modules/porous_flow/examples/coal_mining/fine_with_fluid.i)
- (modules/porous_flow/examples/tutorial/11_2D.i)
- (modules/porous_flow/test/tests/poroperm/PermFromPoro01_fv.i)
- (modules/porous_flow/examples/tutorial/11.i)
- (modules/porous_flow/test/tests/jacobian/fflux08.i)
- (modules/porous_flow/test/tests/poroperm/PermFromPoro02.i)
- (modules/porous_flow/examples/multiapp_fracture_flow/3dFracture/fracture_only_aperture_changing.i)
- (modules/porous_flow/test/tests/jacobian/basic_advection5.i)
- (modules/porous_flow/test/tests/poro_elasticity/vol_expansion_poroperm.i)
- (modules/porous_flow/test/tests/jacobian/basic_advection6.i)
- (modules/porous_flow/test/tests/poroperm/PermFromPoro01.i)
- (modules/porous_flow/examples/coal_mining/coarse_with_fluid.i)
(modules/porous_flow/examples/tutorial/07.i)
# Darcy flow with a tracer that precipitates causing mineralisation and porosity changes and permeability changes
[Mesh]
[annular]
type = AnnularMeshGenerator
nr = 10
rmin = 1.0
rmax = 10
growth_r = 1.4
nt = 4
dmin = 0
dmax = 90
[]
[make3D]
input = annular
type = MeshExtruderGenerator
extrusion_vector = '0 0 12'
num_layers = 3
bottom_sideset = 'bottom'
top_sideset = 'top'
[]
[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]
[]
[tracer_concentration]
[]
[]
[PorousFlowFullySaturated]
porepressure = porepressure
coupling_type = Hydro
gravity = '0 0 0'
fp = the_simple_fluid
mass_fraction_vars = tracer_concentration
number_aqueous_kinetic = 1
temperature = 283.0
stabilization = none # Note to reader: try this with other stabilization and compare the results
[]
[AuxVariables]
[eqm_k]
initial_condition = 0.1
[]
[mineral_conc]
family = MONOMIAL
order = CONSTANT
[]
[initial_and_reference_conc]
initial_condition = 0
[]
[porosity]
family = MONOMIAL
order = CONSTANT
[]
[permeability]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[mineral_conc]
type = PorousFlowPropertyAux
property = mineral_concentration
mineral_species = 0
variable = mineral_conc
[]
[porosity]
type = PorousFlowPropertyAux
property = porosity
variable = porosity
[]
[permeability]
type = PorousFlowPropertyAux
property = permeability
column = 0
row = 0
variable = permeability
[]
[]
[Kernels]
[precipitation_dissolution]
type = PorousFlowPreDis
mineral_density = 1000.0
stoichiometry = 1
variable = tracer_concentration
[]
[]
[BCs]
[constant_injection_of_tracer]
type = PorousFlowSink
variable = tracer_concentration
flux_function = -5E-3
boundary = injection_area
[]
[constant_outer_porepressure]
type = DirichletBC
variable = porepressure
value = 0
boundary = rmax
[]
[]
[FluidProperties]
[the_simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2E9
viscosity = 1.0E-3
density0 = 1000.0
[]
[]
[Materials]
[porosity_mat]
type = PorousFlowPorosity
porosity_zero = 0.1
chemical = true
initial_mineral_concentrations = initial_and_reference_conc
reference_chemistry = initial_and_reference_conc
[]
[permeability_aquifer]
type = PorousFlowPermeabilityKozenyCarman
block = aquifer
k0 = 1E-14
m = 2
n = 3
phi0 = 0.1
poroperm_function = kozeny_carman_phi0
[]
[permeability_caps]
type = PorousFlowPermeabilityKozenyCarman
block = caps
k0 = 1E-15
k_anisotropy = '1 0 0 0 1 0 0 0 0.1'
m = 2
n = 3
phi0 = 0.1
poroperm_function = kozeny_carman_phi0
[]
[precipitation_dissolution_mat]
type = PorousFlowAqueousPreDisChemistry
reference_temperature = 283.0
activation_energy = 1 # irrelevant because T=Tref
equilibrium_constants = eqm_k # equilibrium tracer concentration
kinetic_rate_constant = 1E-8
molar_volume = 1
num_reactions = 1
primary_activity_coefficients = 1
primary_concentrations = tracer_concentration
reactions = 1
specific_reactive_surface_area = 1
[]
[mineral_concentration]
type = PorousFlowAqueousPreDisMineral
[]
[]
[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/test/tests/heterogeneous_materials/vol_expansion_poroperm.i)
# Apply an increasing porepressure, with zero mechanical forces,
# and observe the corresponding volumetric expansion and porosity increase.
# Check that permeability is calculated correctly from porosity.
#
# P = t
# With the Biot coefficient being 1, the effective stresses should be
# stress_xx = stress_yy = stress_zz = t
# With bulk modulus = 1 then should have
# vol_strain = strain_xx + strain_yy + strain_zz = t.
#
# With the biot coefficient being 1, the porosity (phi) # at time t is:
# phi = 1 - (1 - phi0) / exp(vol_strain)
# where phi0 is the porosity at t = 0 and P = 0.
#
# The permeability (k) is
# k = k_anisotropic * f * d^2 * phi^n / (1-phi)^m
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 1
zmin = 0
zmax = 1
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
block = 0
PorousFlowDictator = dictator
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[p]
[]
[]
[BCs]
[p]
type = FunctionDirichletBC
boundary = 'bottom top'
variable = p
function = t
[]
[xmin]
type = DirichletBC
boundary = left
variable = disp_x
value = 0
[]
[ymin]
type = DirichletBC
boundary = bottom
variable = disp_y
value = 0
[]
[zmin]
type = DirichletBC
boundary = back
variable = disp_z
value = 0
[]
[]
[Kernels]
[p_does_not_really_diffuse]
type = Diffusion
variable = p
[]
[TensorMechanics]
displacements = 'disp_x disp_y disp_z'
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 1
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 1
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 1
variable = disp_z
component = 2
[]
[]
[AuxVariables]
[poro0]
order = CONSTANT
family = MONOMIAL
[]
[poro]
order = CONSTANT
family = MONOMIAL
[]
[perm_x]
order = CONSTANT
family = MONOMIAL
[]
[perm_y]
order = CONSTANT
family = MONOMIAL
[]
[perm_z]
order = CONSTANT
family = MONOMIAL
[]
[]
[ICs]
[poro0]
type = RandomIC
seed = 0
variable = poro0
max = 0.15
min = 0.05
[]
[]
[AuxKernels]
[poromat]
type = PorousFlowPropertyAux
property = porosity
variable = poro
[]
[perm_x]
type = PorousFlowPropertyAux
property = permeability
variable = perm_x
row = 0
column = 0
[]
[perm_y]
type = PorousFlowPropertyAux
property = permeability
variable = perm_y
row = 1
column = 1
[]
[perm_z]
type = PorousFlowPropertyAux
property = permeability
variable = perm_z
row = 2
column = 2
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'p'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[elasticity_tensor]
type = ComputeIsotropicElasticityTensor
bulk_modulus = 1
shear_modulus = 1
[]
[strain]
type = ComputeSmallStrain
[]
[stress]
type = ComputeLinearElasticStress
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = p
capillary_pressure = pc
[]
[p_eff]
type = PorousFlowEffectiveFluidPressure
[]
[porosity]
type = PorousFlowPorosity
fluid = true
mechanical = true
porosity_zero = poro0
solid_bulk = 1
biot_coefficient = 1
[]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
k_anisotropy = '1 0 0 0 2 0 0 0 0.1'
poroperm_function = kozeny_carman_fd2
f = 0.1
d = 5
m = 2
n = 7
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_atol -ksp_rtol'
petsc_options_value = 'gmres bjacobi 1E-10 1E-10 10 1E-15 1E-10'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
start_time = 0
dt = 0.1
end_time = 1
[]
[Outputs]
exodus = true
execute_on = 'timestep_end'
[]
(modules/porous_flow/examples/coal_mining/fine_with_fluid.i)
#################################################################
#
# NOTE:
# The mesh for this model is too large for the MOOSE repository
# so is kept in the the large_media submodule
#
#################################################################
#
# Strata deformation and fluid flow aaround a coal mine - 3D model
#
# A "half model" is used. The mine is 400m deep and
# just the roof is studied (-400<=z<=0). The mining panel
# sits between 0<=x<=150, and 0<=y<=1000, so this simulates
# a coal panel that is 300m wide and 1000m long. The outer boundaries
# are 1km from the excavation boundaries.
#
# The excavation takes 0.5 years.
#
# The boundary conditions for this simulation are:
# - disp_x = 0 at x=0 and x=1150
# - disp_y = 0 at y=-1000 and y=1000
# - disp_z = 0 at z=-400, but there is a time-dependent
# Young modulus that simulates excavation
# - wc_x = 0 at y=-1000 and y=1000
# - wc_y = 0 at x=0 and x=1150
# - no flow at x=0, z=-400 and z=0
# - fixed porepressure at y=-1000, y=1000 and x=1150
# That is, rollers on the sides, free at top,
# and prescribed at bottom in the unexcavated portion.
#
# A single-phase unsaturated fluid is used.
#
# The small strain formulation is used.
#
# All stresses are measured in MPa, and time units are measured in years.
#
# The initial porepressure is hydrostatic with P=0 at z=0, so
# Porepressure ~ - 0.01*z MPa, where the fluid has density 1E3 kg/m^3 and
# gravity = = 10 m.s^-2 = 1E-5 MPa m^2/kg.
# To be more accurate, i use
# Porepressure = -bulk * log(1 + g*rho0*z/bulk)
# where bulk=2E3 MPa and rho0=1Ee kg/m^3.
# The initial stress is consistent with the weight force from undrained
# density 2500 kg/m^3, and fluid porepressure, and a Biot coefficient of 0.7, ie,
# stress_zz^effective = 0.025*z + 0.7 * initial_porepressure
# The maximum and minimum principal horizontal effective stresses are
# assumed to be equal to 0.8*stress_zz.
#
# Material properties:
# Young's modulus = 8 GPa
# Poisson's ratio = 0.25
# Cosserat layer thickness = 1 m
# Cosserat-joint normal stiffness = large
# Cosserat-joint shear stiffness = 1 GPa
# MC cohesion = 2 MPa
# MC friction angle = 35 deg
# MC dilation angle = 8 deg
# MC tensile strength = 1 MPa
# MC compressive strength = 100 MPa
# WeakPlane cohesion = 0.1 MPa
# WeakPlane friction angle = 30 deg
# WeakPlane dilation angle = 10 deg
# WeakPlane tensile strength = 0.1 MPa
# WeakPlane compressive strength = 100 MPa softening to 1 MPa at strain = 1
# Fluid density at zero porepressure = 1E3 kg/m^3
# Fluid bulk modulus = 2E3 MPa
# Fluid viscosity = 1.1E-3 Pa.s = 1.1E-9 MPa.s = 3.5E-17 MPa.year
#
[GlobalParams]
perform_finite_strain_rotations = false
displacements = 'disp_x disp_y disp_z'
Cosserat_rotations = 'wc_x wc_y wc_z'
PorousFlowDictator = dictator
biot_coefficient = 0.7
[]
[Mesh]
[file]
type = FileMeshGenerator
file = fine.e
[]
[xmin]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
new_boundary = xmin
normal = '-1 0 0'
input = file
[]
[xmax]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
new_boundary = xmax
normal = '1 0 0'
input = xmin
[]
[ymin]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
new_boundary = ymin
normal = '0 -1 0'
input = xmax
[]
[ymax]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
new_boundary = ymax
normal = '0 1 0'
input = ymin
[]
[zmax]
type = SideSetsAroundSubdomainGenerator
block = 30
new_boundary = zmax
normal = '0 0 1'
input = ymax
[]
[zmin]
type = SideSetsAroundSubdomainGenerator
block = 2
new_boundary = zmin
normal = '0 0 -1'
input = zmax
[]
[excav]
type = SubdomainBoundingBoxGenerator
input = zmin
block_id = 1
bottom_left = '0 0 -400'
top_right = '150 1000 -397'
[]
[roof]
type = SideSetsBetweenSubdomainsGenerator
primary_block = 3
paired_block = 1
input = excav
new_boundary = roof
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[wc_x]
[]
[wc_y]
[]
[porepressure]
scaling = 1E-5
[]
[]
[ICs]
[porepressure]
type = FunctionIC
variable = porepressure
function = ini_pp
[]
[]
[Kernels]
[cx_elastic]
type = CosseratStressDivergenceTensors
use_displaced_mesh = false
variable = disp_x
component = 0
[]
[cy_elastic]
type = CosseratStressDivergenceTensors
use_displaced_mesh = false
variable = disp_y
component = 1
[]
[cz_elastic]
type = CosseratStressDivergenceTensors
use_displaced_mesh = false
variable = disp_z
component = 2
[]
[x_couple]
type = StressDivergenceTensors
use_displaced_mesh = false
variable = wc_x
displacements = 'wc_x wc_y wc_z'
component = 0
base_name = couple
[]
[y_couple]
type = StressDivergenceTensors
use_displaced_mesh = false
variable = wc_y
displacements = 'wc_x wc_y wc_z'
component = 1
base_name = couple
[]
[x_moment]
type = MomentBalancing
use_displaced_mesh = false
variable = wc_x
component = 0
[]
[y_moment]
type = MomentBalancing
use_displaced_mesh = false
variable = wc_y
component = 1
[]
[gravity]
type = Gravity
use_displaced_mesh = false
variable = disp_z
value = -10E-6 # remember this is in MPa
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
use_displaced_mesh = false
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
use_displaced_mesh = false
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
use_displaced_mesh = false
component = 2
variable = disp_z
[]
[poro_vol_exp]
type = PorousFlowMassVolumetricExpansion
use_displaced_mesh = false
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
variable = porepressure
fluid_component = 0
[]
[mass0]
type = PorousFlowMassTimeDerivative
use_displaced_mesh = false
fluid_component = 0
variable = porepressure
[]
[flux]
type = PorousFlowAdvectiveFlux
use_displaced_mesh = false
variable = porepressure
gravity = '0 0 -10E-6'
fluid_component = 0
[]
[]
[AuxVariables]
[saturation]
order = CONSTANT
family = MONOMIAL
[]
[darcy_x]
order = CONSTANT
family = MONOMIAL
[]
[darcy_y]
order = CONSTANT
family = MONOMIAL
[]
[darcy_z]
order = CONSTANT
family = MONOMIAL
[]
[porosity]
order = CONSTANT
family = MONOMIAL
[]
[wc_z]
[]
[stress_xx]
order = CONSTANT
family = MONOMIAL
[]
[stress_xy]
order = CONSTANT
family = MONOMIAL
[]
[stress_xz]
order = CONSTANT
family = MONOMIAL
[]
[stress_yx]
order = CONSTANT
family = MONOMIAL
[]
[stress_yy]
order = CONSTANT
family = MONOMIAL
[]
[stress_yz]
order = CONSTANT
family = MONOMIAL
[]
[stress_zx]
order = CONSTANT
family = MONOMIAL
[]
[stress_zy]
order = CONSTANT
family = MONOMIAL
[]
[stress_zz]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_xx]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_xy]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_xz]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_yx]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_yy]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_yz]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_zx]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_zy]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_zz]
order = CONSTANT
family = MONOMIAL
[]
[perm_xx]
order = CONSTANT
family = MONOMIAL
[]
[perm_yy]
order = CONSTANT
family = MONOMIAL
[]
[perm_zz]
order = CONSTANT
family = MONOMIAL
[]
[mc_shear]
order = CONSTANT
family = MONOMIAL
[]
[mc_tensile]
order = CONSTANT
family = MONOMIAL
[]
[wp_shear]
order = CONSTANT
family = MONOMIAL
[]
[wp_tensile]
order = CONSTANT
family = MONOMIAL
[]
[wp_shear_f]
order = CONSTANT
family = MONOMIAL
[]
[wp_tensile_f]
order = CONSTANT
family = MONOMIAL
[]
[mc_shear_f]
order = CONSTANT
family = MONOMIAL
[]
[mc_tensile_f]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[saturation_water]
type = PorousFlowPropertyAux
variable = saturation
property = saturation
phase = 0
execute_on = timestep_end
[]
[darcy_x]
type = PorousFlowDarcyVelocityComponent
variable = darcy_x
gravity = '0 0 -10E-6'
component = x
[]
[darcy_y]
type = PorousFlowDarcyVelocityComponent
variable = darcy_y
gravity = '0 0 -10E-6'
component = y
[]
[darcy_z]
type = PorousFlowDarcyVelocityComponent
variable = darcy_z
gravity = '0 0 -10E-6'
component = z
[]
[porosity]
type = PorousFlowPropertyAux
property = porosity
variable = porosity
execute_on = timestep_end
[]
[stress_xx]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_xx
index_i = 0
index_j = 0
execute_on = timestep_end
[]
[stress_xy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_xy
index_i = 0
index_j = 1
execute_on = timestep_end
[]
[stress_xz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_xz
index_i = 0
index_j = 2
execute_on = timestep_end
[]
[stress_yx]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yx
index_i = 1
index_j = 0
execute_on = timestep_end
[]
[stress_yy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yy
index_i = 1
index_j = 1
execute_on = timestep_end
[]
[stress_yz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yz
index_i = 1
index_j = 2
execute_on = timestep_end
[]
[stress_zx]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zx
index_i = 2
index_j = 0
execute_on = timestep_end
[]
[stress_zy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zy
index_i = 2
index_j = 1
execute_on = timestep_end
[]
[stress_zz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zz
index_i = 2
index_j = 2
execute_on = timestep_end
[]
[total_strain_xx]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_xx
index_i = 0
index_j = 0
execute_on = timestep_end
[]
[total_strain_xy]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_xy
index_i = 0
index_j = 1
execute_on = timestep_end
[]
[total_strain_xz]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_xz
index_i = 0
index_j = 2
execute_on = timestep_end
[]
[total_strain_yx]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_yx
index_i = 1
index_j = 0
execute_on = timestep_end
[]
[total_strain_yy]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_yy
index_i = 1
index_j = 1
execute_on = timestep_end
[]
[total_strain_yz]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_yz
index_i = 1
index_j = 2
execute_on = timestep_end
[]
[total_strain_zx]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_zx
index_i = 2
index_j = 0
execute_on = timestep_end
[]
[total_strain_zy]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_zy
index_i = 2
index_j = 1
execute_on = timestep_end
[]
[total_strain_zz]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_zz
index_i = 2
index_j = 2
execute_on = timestep_end
[]
[perm_xx]
type = PorousFlowPropertyAux
property = permeability
variable = perm_xx
row = 0
column = 0
execute_on = timestep_end
[]
[perm_yy]
type = PorousFlowPropertyAux
property = permeability
variable = perm_yy
row = 1
column = 1
execute_on = timestep_end
[]
[perm_zz]
type = PorousFlowPropertyAux
property = permeability
variable = perm_zz
row = 2
column = 2
execute_on = timestep_end
[]
[mc_shear]
type = MaterialStdVectorAux
index = 0
property = mc_plastic_internal_parameter
variable = mc_shear
execute_on = timestep_end
[]
[mc_tensile]
type = MaterialStdVectorAux
index = 1
property = mc_plastic_internal_parameter
variable = mc_tensile
execute_on = timestep_end
[]
[wp_shear]
type = MaterialStdVectorAux
index = 0
property = wp_plastic_internal_parameter
variable = wp_shear
execute_on = timestep_end
[]
[wp_tensile]
type = MaterialStdVectorAux
index = 1
property = wp_plastic_internal_parameter
variable = wp_tensile
execute_on = timestep_end
[]
[mc_shear_f]
type = MaterialStdVectorAux
index = 6
property = mc_plastic_yield_function
variable = mc_shear_f
execute_on = timestep_end
[]
[mc_tensile_f]
type = MaterialStdVectorAux
index = 0
property = mc_plastic_yield_function
variable = mc_tensile_f
execute_on = timestep_end
[]
[wp_shear_f]
type = MaterialStdVectorAux
index = 0
property = wp_plastic_yield_function
variable = wp_shear_f
execute_on = timestep_end
[]
[wp_tensile_f]
type = MaterialStdVectorAux
index = 1
property = wp_plastic_yield_function
variable = wp_tensile_f
execute_on = timestep_end
[]
[]
[BCs]
[no_x]
type = DirichletBC
variable = disp_x
boundary = 'xmin xmax'
value = 0.0
[]
[no_y]
type = DirichletBC
variable = disp_y
boundary = 'ymin ymax'
value = 0.0
[]
[no_z]
type = DirichletBC
variable = disp_z
boundary = zmin
value = 0.0
[]
[no_wc_x]
type = DirichletBC
variable = wc_x
boundary = 'ymin ymax'
value = 0.0
[]
[no_wc_y]
type = DirichletBC
variable = wc_y
boundary = 'xmin xmax'
value = 0.0
[]
[fix_porepressure]
type = FunctionDirichletBC
variable = porepressure
boundary = 'ymin ymax xmax'
function = ini_pp
[]
[roof_porepressure]
type = PorousFlowPiecewiseLinearSink
variable = porepressure
pt_vals = '-1E3 1E3'
multipliers = '-1 1'
fluid_phase = 0
flux_function = roof_conductance
boundary = roof
[]
[roof]
type = StickyBC
variable = disp_z
min_value = -3.0
boundary = roof
[]
[]
[Functions]
[ini_pp]
type = ParsedFunction
symbol_names = 'bulk p0 g rho0'
symbol_values = '2E3 0.0 1E-5 1E3'
expression = '-bulk*log(exp(-p0/bulk)+g*rho0*z/bulk)'
[]
[ini_xx]
type = ParsedFunction
symbol_names = 'bulk p0 g rho0 biot'
symbol_values = '2E3 0.0 1E-5 1E3 0.7'
expression = '0.8*(2500*10E-6*z+biot*(-bulk*log(exp(-p0/bulk)+g*rho0*z/bulk)))'
[]
[ini_zz]
type = ParsedFunction
symbol_names = 'bulk p0 g rho0 biot'
symbol_values = '2E3 0.0 1E-5 1E3 0.7'
expression = '2500*10E-6*z+biot*(-bulk*log(exp(-p0/bulk)+g*rho0*z/bulk))'
[]
[excav_sideways]
type = ParsedFunction
symbol_names = 'end_t ymin ymax minval maxval slope'
symbol_values = '0.5 0 1000.0 1E-9 1 10'
# excavation face at ymin+(ymax-ymin)*min(t/end_t,1)
# slope is the distance over which the modulus reduces from maxval to minval
expression = 'if(y<ymin+(ymax-ymin)*min(t/end_t,1),minval,if(y<ymin+(ymax-ymin)*min(t/end_t,1)+slope,minval+(maxval-minval)*(y-(ymin+(ymax-ymin)*min(t/end_t,1)))/slope,maxval))'
[]
[density_sideways]
type = ParsedFunction
symbol_names = 'end_t ymin ymax minval maxval'
symbol_values = '0.5 0 1000.0 0 2500'
expression = 'if(y<ymin+(ymax-ymin)*min(t/end_t,1),minval,maxval)'
[]
[roof_conductance]
type = ParsedFunction
symbol_names = 'end_t ymin ymax maxval minval'
symbol_values = '0.5 0 1000.0 1E7 0'
expression = 'if(y<ymin+(ymax-ymin)*min(t/end_t,1),maxval,minval)'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'porepressure disp_x disp_y disp_z'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1 # MPa^-1
[]
[mc_coh_strong_harden]
type = TensorMechanicsHardeningExponential
value_0 = 1.99 # MPa
value_residual = 2.01 # MPa
rate = 1.0
[]
[mc_fric]
type = TensorMechanicsHardeningConstant
value = 0.61 # 35deg
[]
[mc_dil]
type = TensorMechanicsHardeningConstant
value = 0.15 # 8deg
[]
[mc_tensile_str_strong_harden]
type = TensorMechanicsHardeningExponential
value_0 = 1.0 # MPa
value_residual = 1.0 # MPa
rate = 1.0
[]
[mc_compressive_str]
type = TensorMechanicsHardeningCubic
value_0 = 100 # Large!
value_residual = 100
internal_limit = 0.1
[]
[wp_coh_harden]
type = TensorMechanicsHardeningCubic
value_0 = 0.05
value_residual = 0.05
internal_limit = 10
[]
[wp_tan_fric]
type = TensorMechanicsHardeningConstant
value = 0.26 # 15deg
[]
[wp_tan_dil]
type = TensorMechanicsHardeningConstant
value = 0.18 # 10deg
[]
[wp_tensile_str_harden]
type = TensorMechanicsHardeningCubic
value_0 = 0.05
value_residual = 0.05
internal_limit = 10
[]
[wp_compressive_str_soften]
type = TensorMechanicsHardeningCubic
value_0 = 100
value_residual = 1
internal_limit = 1.0
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2E3
density0 = 1000
thermal_expansion = 0
viscosity = 3.5E-17
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[eff_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = porepressure
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity_for_aux]
type = PorousFlowPorosity
at_nodes = false
fluid = true
mechanical = true
ensure_positive = true
porosity_zero = 0.02
solid_bulk = 5.3333E3
[]
[porosity_bulk]
type = PorousFlowPorosity
fluid = true
mechanical = true
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
ensure_positive = true
porosity_zero = 0.02
solid_bulk = 5.3333E3
[]
[porosity_excav]
type = PorousFlowPorosityConst
block = 1
porosity = 1.0
[]
[permeability_bulk]
type = PorousFlowPermeabilityKozenyCarman
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
poroperm_function = kozeny_carman_phi0
k0 = 1E-15
phi0 = 0.02
n = 2
m = 2
[]
[permeability_excav]
type = PorousFlowPermeabilityConst
block = 1
permeability = '0 0 0 0 0 0 0 0 0'
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 4
s_res = 0.4
sum_s_res = 0.4
phase = 0
[]
[elasticity_tensor_0]
type = ComputeLayeredCosseratElasticityTensor
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
young = 8E3 # MPa
poisson = 0.25
layer_thickness = 1.0
joint_normal_stiffness = 1E9 # huge
joint_shear_stiffness = 1E3 # MPa
[]
[elasticity_tensor_1]
type = ComputeLayeredCosseratElasticityTensor
block = 1
young = 8E3 # MPa
poisson = 0.25
layer_thickness = 1.0
joint_normal_stiffness = 1E9 # huge
joint_shear_stiffness = 1E3 # MPa
elasticity_tensor_prefactor = excav_sideways
[]
[strain]
type = ComputeCosseratIncrementalSmallStrain
eigenstrain_names = ini_stress
[]
[ini_stress]
type = ComputeEigenstrainFromInitialStress
eigenstrain_name = ini_stress
initial_stress = 'ini_xx 0 0 0 ini_xx 0 0 0 ini_zz'
[]
[stress_0]
type = ComputeMultipleInelasticCosseratStress
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
inelastic_models = 'mc wp'
cycle_models = true
relative_tolerance = 2.0
absolute_tolerance = 1E6
max_iterations = 1
tangent_operator = nonlinear
perform_finite_strain_rotations = false
[]
[stress_1]
type = ComputeMultipleInelasticCosseratStress
block = 1
inelastic_models = ''
relative_tolerance = 2.0
absolute_tolerance = 1E6
max_iterations = 1
tangent_operator = nonlinear
perform_finite_strain_rotations = false
[]
[mc]
type = CappedMohrCoulombCosseratStressUpdate
warn_about_precision_loss = false
host_youngs_modulus = 8E3
host_poissons_ratio = 0.25
base_name = mc
tensile_strength = mc_tensile_str_strong_harden
compressive_strength = mc_compressive_str
cohesion = mc_coh_strong_harden
friction_angle = mc_fric
dilation_angle = mc_dil
max_NR_iterations = 100000
smoothing_tol = 0.1 # MPa # Must be linked to cohesion
yield_function_tol = 1E-9 # MPa. this is essentially the lowest possible without lots of precision loss
perfect_guess = true
min_step_size = 1.0
[]
[wp]
type = CappedWeakPlaneCosseratStressUpdate
warn_about_precision_loss = false
base_name = wp
cohesion = wp_coh_harden
tan_friction_angle = wp_tan_fric
tan_dilation_angle = wp_tan_dil
tensile_strength = wp_tensile_str_harden
compressive_strength = wp_compressive_str_soften
max_NR_iterations = 10000
tip_smoother = 0.05
smoothing_tol = 0.05 # MPa # Note, this must be tied to cohesion, otherwise get no possible return at cone apex
yield_function_tol = 1E-11 # MPa. this is essentially the lowest possible without lots of precision loss
perfect_guess = true
min_step_size = 1.0E-3
[]
[undrained_density_0]
type = GenericConstantMaterial
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
prop_names = density
prop_values = 2500
[]
[undrained_density_1]
type = GenericFunctionMaterial
block = 1
prop_names = density
prop_values = density_sideways
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[min_roof_disp]
type = NodalExtremeValue
boundary = roof
value_type = min
variable = disp_z
[]
[min_roof_pp]
type = NodalExtremeValue
boundary = roof
value_type = min
variable = porepressure
[]
[min_surface_disp]
type = NodalExtremeValue
boundary = zmax
value_type = min
variable = disp_z
[]
[min_surface_pp]
type = NodalExtremeValue
boundary = zmax
value_type = min
variable = porepressure
[]
[max_perm_zz]
type = ElementExtremeValue
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30'
variable = perm_zz
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason'
# best overall
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
# best if you don't have mumps:
#petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -ksp_type -ksp_gmres_restart'
#petsc_options_value = ' asm 2 lu gmres 200'
# very basic:
#petsc_options_iname = '-pc_type -ksp_type -ksp_gmres_restart'
#petsc_options_value = ' bjacobi gmres 200'
line_search = bt
nl_abs_tol = 1e-3
nl_rel_tol = 1e-5
l_max_its = 200
nl_max_its = 30
start_time = 0.0
dt = 0.0025
end_time = 0.5
[]
[Outputs]
time_step_interval = 1
print_linear_residuals = true
exodus = true
csv = true
console = true
[]
(modules/porous_flow/examples/tutorial/11_2D.i)
# Two-phase borehole injection problem in RZ coordinates
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
xmin = 1.0
xmax = 10
bias_x = 1.4
ny = 3
ymin = -6
ymax = 6
[]
[aquifer]
input = gen
type = SubdomainBoundingBoxGenerator
block_id = 1
bottom_left = '0 -2 0'
top_right = '10 2 0'
[]
[injection_area]
type = ParsedGenerateSideset
combinatorial_geometry = 'x<1.0001'
included_subdomains = 1
new_sideset_name = 'injection_area'
input = 'aquifer'
[]
[rename]
type = RenameBlockGenerator
old_block = '0 1'
new_block = 'caps aquifer'
input = 'injection_area'
[]
coord_type = RZ
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pwater pgas T disp_r'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 1E-6
m = 0.6
[]
[]
[GlobalParams]
displacements = 'disp_r disp_z'
gravity = '0 0 0'
biot_coefficient = 1.0
PorousFlowDictator = dictator
[]
[Variables]
[pwater]
initial_condition = 20E6
[]
[pgas]
initial_condition = 20.1E6
[]
[T]
initial_condition = 330
scaling = 1E-5
[]
[disp_r]
scaling = 1E-5
[]
[]
[Kernels]
[mass_water_dot]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pwater
[]
[flux_water]
type = PorousFlowAdvectiveFlux
fluid_component = 0
use_displaced_mesh = false
variable = pwater
[]
[vol_strain_rate_water]
type = PorousFlowMassVolumetricExpansion
fluid_component = 0
variable = pwater
[]
[mass_co2_dot]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = pgas
[]
[flux_co2]
type = PorousFlowAdvectiveFlux
fluid_component = 1
use_displaced_mesh = false
variable = pgas
[]
[vol_strain_rate_co2]
type = PorousFlowMassVolumetricExpansion
fluid_component = 1
variable = pgas
[]
[energy_dot]
type = PorousFlowEnergyTimeDerivative
variable = T
[]
[advection]
type = PorousFlowHeatAdvection
use_displaced_mesh = false
variable = T
[]
[conduction]
type = PorousFlowHeatConduction
use_displaced_mesh = false
variable = T
[]
[vol_strain_rate_heat]
type = PorousFlowHeatVolumetricExpansion
variable = T
[]
[grad_stress_r]
type = StressDivergenceRZTensors
temperature = T
variable = disp_r
eigenstrain_names = thermal_contribution
use_displaced_mesh = false
component = 0
[]
[poro_r]
type = PorousFlowEffectiveStressCoupling
variable = disp_r
use_displaced_mesh = false
component = 0
[]
[]
[AuxVariables]
[disp_z]
[]
[effective_fluid_pressure]
family = MONOMIAL
order = CONSTANT
[]
[mass_frac_phase0_species0]
initial_condition = 1 # all water in phase=0
[]
[mass_frac_phase1_species0]
initial_condition = 0 # no water in phase=1
[]
[sgas]
family = MONOMIAL
order = CONSTANT
[]
[swater]
family = MONOMIAL
order = CONSTANT
[]
[stress_rr]
family = MONOMIAL
order = CONSTANT
[]
[stress_tt]
family = MONOMIAL
order = CONSTANT
[]
[stress_zz]
family = MONOMIAL
order = CONSTANT
[]
[porosity]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[effective_fluid_pressure]
type = ParsedAux
coupled_variables = 'pwater pgas swater sgas'
expression = 'pwater * swater + pgas * sgas'
variable = effective_fluid_pressure
[]
[swater]
type = PorousFlowPropertyAux
variable = swater
property = saturation
phase = 0
execute_on = timestep_end
[]
[sgas]
type = PorousFlowPropertyAux
variable = sgas
property = saturation
phase = 1
execute_on = timestep_end
[]
[stress_rr_aux]
type = RankTwoAux
variable = stress_rr
rank_two_tensor = stress
index_i = 0
index_j = 0
[]
[stress_tt]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_tt
index_i = 2
index_j = 2
[]
[stress_zz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zz
index_i = 1
index_j = 1
[]
[porosity]
type = PorousFlowPropertyAux
variable = porosity
property = porosity
execute_on = timestep_end
[]
[]
[BCs]
[pinned_top_bottom_r]
type = DirichletBC
variable = disp_r
value = 0
boundary = 'top bottom'
[]
[cavity_pressure_r]
type = Pressure
boundary = injection_area
variable = disp_r
postprocessor = constrained_effective_fluid_pressure_at_wellbore
use_displaced_mesh = false
[]
[cold_co2]
type = DirichletBC
boundary = injection_area
variable = T
value = 290 # injection temperature
use_displaced_mesh = false
[]
[constant_co2_injection]
type = PorousFlowSink
boundary = injection_area
variable = pgas
fluid_phase = 1
flux_function = -1E-4
use_displaced_mesh = false
[]
[outer_water_removal]
type = PorousFlowPiecewiseLinearSink
boundary = right
variable = pwater
fluid_phase = 0
pt_vals = '0 1E9'
multipliers = '0 1E8'
PT_shift = 20E6
use_mobility = true
use_relperm = true
use_displaced_mesh = false
[]
[outer_co2_removal]
type = PorousFlowPiecewiseLinearSink
boundary = right
variable = pgas
fluid_phase = 1
pt_vals = '0 1E9'
multipliers = '0 1E8'
PT_shift = 20.1E6
use_mobility = true
use_relperm = true
use_displaced_mesh = false
[]
[]
[FluidProperties]
[true_water]
type = Water97FluidProperties
[]
[tabulated_water]
type = TabulatedBicubicFluidProperties
fp = true_water
temperature_min = 275
pressure_max = 1E8
fluid_property_file = water97_tabulated_11.csv
[]
[true_co2]
type = CO2FluidProperties
[]
[tabulated_co2]
type = TabulatedBicubicFluidProperties
fp = true_co2
temperature_min = 275
pressure_max = 1E8
fluid_property_file = co2_tabulated_11.csv
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = T
[]
[saturation_calculator]
type = PorousFlow2PhasePP
phase0_porepressure = pwater
phase1_porepressure = pgas
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'mass_frac_phase0_species0 mass_frac_phase1_species0'
[]
[water]
type = PorousFlowSingleComponentFluid
fp = tabulated_water
phase = 0
[]
[co2]
type = PorousFlowSingleComponentFluid
fp = tabulated_co2
phase = 1
[]
[relperm_water]
type = PorousFlowRelativePermeabilityCorey
n = 4
s_res = 0.1
sum_s_res = 0.2
phase = 0
[]
[relperm_co2]
type = PorousFlowRelativePermeabilityBC
nw_phase = true
lambda = 2
s_res = 0.1
sum_s_res = 0.2
phase = 1
[]
[porosity]
type = PorousFlowPorosity
fluid = true
mechanical = true
thermal = true
porosity_zero = 0.1
reference_temperature = 330
reference_porepressure = 20E6
thermal_expansion_coeff = 15E-6 # volumetric
solid_bulk = 8E9 # unimportant since biot = 1
[]
[permeability_aquifer]
type = PorousFlowPermeabilityKozenyCarman
block = aquifer
poroperm_function = kozeny_carman_phi0
phi0 = 0.1
n = 2
m = 2
k0 = 1E-12
[]
[permeability_caps]
type = PorousFlowPermeabilityKozenyCarman
block = caps
poroperm_function = kozeny_carman_phi0
phi0 = 0.1
n = 2
m = 2
k0 = 1E-15
k_anisotropy = '1 0 0 0 1 0 0 0 0.1'
[]
[rock_thermal_conductivity]
type = PorousFlowThermalConductivityIdeal
dry_thermal_conductivity = '2 0 0 0 2 0 0 0 2'
[]
[rock_internal_energy]
type = PorousFlowMatrixInternalEnergy
specific_heat_capacity = 1100
density = 2300
[]
[elasticity_tensor]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 5E9
poissons_ratio = 0.0
[]
[strain]
type = ComputeAxisymmetricRZSmallStrain
eigenstrain_names = 'thermal_contribution initial_stress'
[]
[thermal_contribution]
type = ComputeThermalExpansionEigenstrain
temperature = T
thermal_expansion_coeff = 5E-6 # this is the linear thermal expansion coefficient
eigenstrain_name = thermal_contribution
stress_free_temperature = 330
[]
[initial_strain]
type = ComputeEigenstrainFromInitialStress
initial_stress = '20E6 0 0 0 20E6 0 0 0 20E6'
eigenstrain_name = initial_stress
[]
[stress]
type = ComputeLinearElasticStress
[]
[effective_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[volumetric_strain]
type = PorousFlowVolumetricStrain
[]
[]
[Postprocessors]
[effective_fluid_pressure_at_wellbore]
type = PointValue
variable = effective_fluid_pressure
point = '1 0 0'
execute_on = timestep_begin
use_displaced_mesh = false
[]
[constrained_effective_fluid_pressure_at_wellbore]
type = FunctionValuePostprocessor
function = constrain_effective_fluid_pressure
execute_on = timestep_begin
[]
[]
[Functions]
[constrain_effective_fluid_pressure]
type = ParsedFunction
symbol_names = effective_fluid_pressure_at_wellbore
symbol_values = effective_fluid_pressure_at_wellbore
expression = 'max(effective_fluid_pressure_at_wellbore, 20E6)'
[]
[]
[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 = 1E3
[TimeStepper]
type = IterationAdaptiveDT
dt = 1E3
growth_factor = 1.2
optimal_iterations = 10
[]
nl_abs_tol = 1E-7
[]
[Outputs]
exodus = true
[]
(modules/porous_flow/test/tests/poroperm/PermFromPoro01_fv.i)
# Testing permeability from porosity
# Trivial test, checking calculated permeability is correct
# k = k_anisotropic * f * d^2 * phi^n / (1-phi)^m
[Mesh]
[mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 3
xmin = 0
xmax = 3
[]
[]
[GlobalParams]
block = 0
PorousFlowDictator = dictator
[]
[Variables]
[pp]
type = MooseVariableFVReal
[FVInitialCondition]
type = FVConstantIC
value = 0
[]
[]
[]
[FVKernels]
[flux]
type = FVPorousFlowAdvectiveFlux
gravity = '0 0 0'
variable = pp
[]
[]
[FVBCs]
[ptop]
type = FVDirichletBC
variable = pp
boundary = right
value = 0
[]
[pbase]
type = FVDirichletBC
variable = pp
boundary = left
value = 1
[]
[]
[AuxVariables]
[poro]
type = MooseVariableFVReal
[]
[perm_x]
type = MooseVariableFVReal
[]
[perm_y]
type = MooseVariableFVReal
[]
[perm_z]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[poro]
type = ADPorousFlowPropertyAux
property = porosity
variable = poro
[]
[perm_x]
type = ADPorousFlowPropertyAux
property = permeability
variable = perm_x
row = 0
column = 0
[]
[perm_y]
type = ADPorousFlowPropertyAux
property = permeability
variable = perm_y
row = 1
column = 1
[]
[perm_z]
type = ADPorousFlowPropertyAux
property = permeability
variable = perm_z
row = 2
column = 2
[]
[]
[Postprocessors]
[perm_x_bottom]
type = PointValue
variable = perm_x
point = '0 0 0'
[]
[perm_y_bottom]
type = PointValue
variable = perm_y
point = '0 0 0'
[]
[perm_z_bottom]
type = PointValue
variable = perm_z
point = '0 0 0'
[]
[perm_x_top]
type = PointValue
variable = perm_x
point = '3 0 0'
[]
[perm_y_top]
type = PointValue
variable = perm_y
point = '3 0 0'
[]
[perm_z_top]
type = PointValue
variable = perm_z
point = '3 0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
# unimportant in this fully-saturated test
m = 0.8
alpha = 1e-4
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2.2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[permeability]
type = ADPorousFlowPermeabilityKozenyCarman
k_anisotropy = '1 0 0 0 2 0 0 0 0.1'
poroperm_function = kozeny_carman_fd2
f = 0.1
d = 5
m = 2
n = 7
[]
[temperature]
type = ADPorousFlowTemperature
[]
[massfrac]
type = ADPorousFlowMassFraction
[]
[eff_fluid_pressure]
type = ADPorousFlowEffectiveFluidPressure
[]
[ppss]
type = ADPorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[simple_fluid]
type = ADPorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = ADPorousFlowPorosityConst
porosity = 0.1
[]
[relperm]
type = ADPorousFlowRelativePermeabilityCorey
n = 0 # unimportant in this fully-saturated situation
phase = 0
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = Newton
type = Steady
l_tol = 1E-5
nl_abs_tol = 1E-3
nl_rel_tol = 1E-8
l_max_its = 200
nl_max_its = 400
[]
[Outputs]
file_base = 'PermFromPoro01_out'
csv = true
execute_on = 'timestep_end'
[]
(modules/porous_flow/examples/tutorial/11.i)
# Two-phase borehole injection problem
[Mesh]
[annular]
type = AnnularMeshGenerator
nr = 10
rmin = 1.0
rmax = 10
growth_r = 1.4
nt = 4
dmin = 0
dmax = 90
[]
[make3D]
input = annular
type = MeshExtruderGenerator
extrusion_vector = '0 0 12'
num_layers = 3
bottom_sideset = 'bottom'
top_sideset = 'top'
[]
[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'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pwater pgas T disp_x disp_y'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 1E-6
m = 0.6
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
gravity = '0 0 0'
biot_coefficient = 1.0
PorousFlowDictator = dictator
[]
[Variables]
[pwater]
initial_condition = 20E6
[]
[pgas]
initial_condition = 20.1E6
[]
[T]
initial_condition = 330
scaling = 1E-5
[]
[disp_x]
scaling = 1E-5
[]
[disp_y]
scaling = 1E-5
[]
[]
[Kernels]
[mass_water_dot]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pwater
[]
[flux_water]
type = PorousFlowAdvectiveFlux
fluid_component = 0
use_displaced_mesh = false
variable = pwater
[]
[vol_strain_rate_water]
type = PorousFlowMassVolumetricExpansion
fluid_component = 0
variable = pwater
[]
[mass_co2_dot]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = pgas
[]
[flux_co2]
type = PorousFlowAdvectiveFlux
fluid_component = 1
use_displaced_mesh = false
variable = pgas
[]
[vol_strain_rate_co2]
type = PorousFlowMassVolumetricExpansion
fluid_component = 1
variable = pgas
[]
[energy_dot]
type = PorousFlowEnergyTimeDerivative
variable = T
[]
[advection]
type = PorousFlowHeatAdvection
use_displaced_mesh = false
variable = T
[]
[conduction]
type = PorousFlowHeatConduction
use_displaced_mesh = false
variable = T
[]
[vol_strain_rate_heat]
type = PorousFlowHeatVolumetricExpansion
variable = T
[]
[grad_stress_x]
type = StressDivergenceTensors
temperature = T
variable = disp_x
eigenstrain_names = thermal_contribution
use_displaced_mesh = false
component = 0
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
variable = disp_x
use_displaced_mesh = false
component = 0
[]
[grad_stress_y]
type = StressDivergenceTensors
temperature = T
variable = disp_y
eigenstrain_names = thermal_contribution
use_displaced_mesh = false
component = 1
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
variable = disp_y
use_displaced_mesh = false
component = 1
[]
[]
[AuxVariables]
[disp_z]
[]
[effective_fluid_pressure]
family = MONOMIAL
order = CONSTANT
[]
[mass_frac_phase0_species0]
initial_condition = 1 # all water in phase=0
[]
[mass_frac_phase1_species0]
initial_condition = 0 # no water in phase=1
[]
[sgas]
family = MONOMIAL
order = CONSTANT
[]
[swater]
family = MONOMIAL
order = CONSTANT
[]
[stress_rr]
family = MONOMIAL
order = CONSTANT
[]
[stress_tt]
family = MONOMIAL
order = CONSTANT
[]
[stress_zz]
family = MONOMIAL
order = CONSTANT
[]
[porosity]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[effective_fluid_pressure]
type = ParsedAux
coupled_variables = 'pwater pgas swater sgas'
expression = 'pwater * swater + pgas * sgas'
variable = effective_fluid_pressure
[]
[swater]
type = PorousFlowPropertyAux
variable = swater
property = saturation
phase = 0
execute_on = timestep_end
[]
[sgas]
type = PorousFlowPropertyAux
variable = sgas
property = saturation
phase = 1
execute_on = timestep_end
[]
[stress_rr]
type = RankTwoScalarAux
variable = stress_rr
rank_two_tensor = stress
scalar_type = RadialStress
point1 = '0 0 0'
point2 = '0 0 1'
execute_on = timestep_end
[]
[stress_tt]
type = RankTwoScalarAux
variable = stress_tt
rank_two_tensor = stress
scalar_type = HoopStress
point1 = '0 0 0'
point2 = '0 0 1'
execute_on = timestep_end
[]
[stress_zz]
type = RankTwoAux
variable = stress_zz
rank_two_tensor = stress
index_i = 2
index_j = 2
execute_on = timestep_end
[]
[porosity]
type = PorousFlowPropertyAux
variable = porosity
property = porosity
execute_on = timestep_end
[]
[]
[BCs]
[roller_tmax]
type = DirichletBC
variable = disp_x
value = 0
boundary = dmax
[]
[roller_tmin]
type = DirichletBC
variable = disp_y
value = 0
boundary = dmin
[]
[pinned_top_bottom_x]
type = DirichletBC
variable = disp_x
value = 0
boundary = 'top bottom'
[]
[pinned_top_bottom_y]
type = DirichletBC
variable = disp_y
value = 0
boundary = 'top bottom'
[]
[cavity_pressure_x]
type = Pressure
boundary = injection_area
variable = disp_x
component = 0
postprocessor = constrained_effective_fluid_pressure_at_wellbore
use_displaced_mesh = false
[]
[cavity_pressure_y]
type = Pressure
boundary = injection_area
variable = disp_y
component = 1
postprocessor = constrained_effective_fluid_pressure_at_wellbore
use_displaced_mesh = false
[]
[cold_co2]
type = DirichletBC
boundary = injection_area
variable = T
value = 290 # injection temperature
use_displaced_mesh = false
[]
[constant_co2_injection]
type = PorousFlowSink
boundary = injection_area
variable = pgas
fluid_phase = 1
flux_function = -1E-4
use_displaced_mesh = false
[]
[outer_water_removal]
type = PorousFlowPiecewiseLinearSink
boundary = rmax
variable = pwater
fluid_phase = 0
pt_vals = '0 1E9'
multipliers = '0 1E8'
PT_shift = 20E6
use_mobility = true
use_relperm = true
use_displaced_mesh = false
[]
[outer_co2_removal]
type = PorousFlowPiecewiseLinearSink
boundary = rmax
variable = pgas
fluid_phase = 1
pt_vals = '0 1E9'
multipliers = '0 1E8'
PT_shift = 20.1E6
use_mobility = true
use_relperm = true
use_displaced_mesh = false
[]
[]
[FluidProperties]
[true_water]
type = Water97FluidProperties
[]
[tabulated_water]
type = TabulatedFluidProperties
fp = true_water
temperature_min = 275
pressure_max = 1E8
interpolated_properties = 'density viscosity enthalpy internal_energy'
fluid_property_file = water97_tabulated_11.csv
[]
[true_co2]
type = CO2FluidProperties
[]
[tabulated_co2]
type = TabulatedFluidProperties
fp = true_co2
temperature_min = 275
pressure_max = 1E8
interpolated_properties = 'density viscosity enthalpy internal_energy'
fluid_property_file = co2_tabulated_11.csv
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = T
[]
[saturation_calculator]
type = PorousFlow2PhasePP
phase0_porepressure = pwater
phase1_porepressure = pgas
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'mass_frac_phase0_species0 mass_frac_phase1_species0'
[]
[water]
type = PorousFlowSingleComponentFluid
fp = tabulated_water
phase = 0
[]
[co2]
type = PorousFlowSingleComponentFluid
fp = tabulated_co2
phase = 1
[]
[relperm_water]
type = PorousFlowRelativePermeabilityCorey
n = 4
s_res = 0.1
sum_s_res = 0.2
phase = 0
[]
[relperm_co2]
type = PorousFlowRelativePermeabilityBC
nw_phase = true
lambda = 2
s_res = 0.1
sum_s_res = 0.2
phase = 1
[]
[porosity_mat]
type = PorousFlowPorosity
fluid = true
mechanical = true
thermal = true
porosity_zero = 0.1
reference_temperature = 330
reference_porepressure = 20E6
thermal_expansion_coeff = 15E-6 # volumetric
solid_bulk = 8E9 # unimportant since biot = 1
[]
[permeability_aquifer]
type = PorousFlowPermeabilityKozenyCarman
block = aquifer
poroperm_function = kozeny_carman_phi0
phi0 = 0.1
n = 2
m = 2
k0 = 1E-12
[]
[permeability_caps]
type = PorousFlowPermeabilityKozenyCarman
block = caps
poroperm_function = kozeny_carman_phi0
phi0 = 0.1
n = 2
m = 2
k0 = 1E-15
k_anisotropy = '1 0 0 0 1 0 0 0 0.1'
[]
[rock_thermal_conductivity]
type = PorousFlowThermalConductivityIdeal
dry_thermal_conductivity = '2 0 0 0 2 0 0 0 2'
[]
[rock_internal_energy]
type = PorousFlowMatrixInternalEnergy
specific_heat_capacity = 1100
density = 2300
[]
[elasticity_tensor]
type = ComputeIsotropicElasticityTensor
youngs_modulus = 5E9
poissons_ratio = 0.0
[]
[strain]
type = ComputeSmallStrain
eigenstrain_names = 'thermal_contribution initial_stress'
[]
[thermal_contribution]
type = ComputeThermalExpansionEigenstrain
temperature = T
thermal_expansion_coeff = 5E-6 # this is the linear thermal expansion coefficient
eigenstrain_name = thermal_contribution
stress_free_temperature = 330
[]
[initial_strain]
type = ComputeEigenstrainFromInitialStress
initial_stress = '20E6 0 0 0 20E6 0 0 0 20E6'
eigenstrain_name = initial_stress
[]
[stress]
type = ComputeLinearElasticStress
[]
[effective_fluid_pressure_mat]
type = PorousFlowEffectiveFluidPressure
[]
[volumetric_strain]
type = PorousFlowVolumetricStrain
[]
[]
[Postprocessors]
[effective_fluid_pressure_at_wellbore]
type = PointValue
variable = effective_fluid_pressure
point = '1 0 0'
execute_on = timestep_begin
use_displaced_mesh = false
[]
[constrained_effective_fluid_pressure_at_wellbore]
type = FunctionValuePostprocessor
function = constrain_effective_fluid_pressure
execute_on = timestep_begin
[]
[]
[Functions]
[constrain_effective_fluid_pressure]
type = ParsedFunction
symbol_names = effective_fluid_pressure_at_wellbore
symbol_values = effective_fluid_pressure_at_wellbore
expression = 'max(effective_fluid_pressure_at_wellbore, 20E6)'
[]
[]
[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 = 1E3
[TimeStepper]
type = IterationAdaptiveDT
dt = 1E3
growth_factor = 1.2
optimal_iterations = 10
[]
nl_abs_tol = 1E-7
[]
[Outputs]
exodus = true
[]
(modules/porous_flow/test/tests/jacobian/fflux08.i)
# 1phase, 1component, constant viscosity, Kozeny-Carman permeability
# density with constant bulk, Corey relative perm, nonzero gravity, unsaturated with vanGenuchten
[Mesh]
type = GeneratedMesh
dim = 3
[]
[GlobalParams]
PorousFlowDictator = dictator
displacements = 'disp_x disp_y disp_z'
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[pp]
[]
[]
[ICs]
[disp_x]
type = RandomIC
variable = disp_x
min = -0.1
max = 0.1
[]
[disp_y]
type = RandomIC
variable = disp_y
min = -0.1
max = 0.1
[]
[disp_z]
type = RandomIC
variable = disp_z
min = -0.1
max = 0.1
[]
[pp]
type = RandomIC
variable = pp
min = -1
max = 1
[]
[]
[Kernels]
[grad_stress_x]
type = StressDivergenceTensors
variable = disp_x
component = 0
[]
[grad_stress_y]
type = StressDivergenceTensors
variable = disp_y
component = 1
[]
[grad_stress_z]
type = StressDivergenceTensors
variable = disp_z
component = 2
[]
[flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
gravity = '-1 -0.1 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp disp_x disp_y disp_z'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1.5
density0 = 1
thermal_expansion = 0
viscosity = 1
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[elasticity_tensor]
type = ComputeElasticityTensor
C_ijkl = '0.5 0.75'
# bulk modulus is lambda + 2*mu/3 = 0.5 + 2*0.75/3 = 1
fill_method = symmetric_isotropic
[]
[strain]
type = ComputeSmallStrain
[]
[stress]
type = ComputeLinearElasticStress
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[porosity]
type = PorousFlowPorosity
fluid = true
mechanical = true
porosity_zero = 0.1
biot_coefficient = 0.5
solid_bulk = 1
[]
[p_eff]
type = PorousFlowEffectiveFluidPressure
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
poroperm_function = kozeny_carman_phi0
k_anisotropy = '1 0 0 0 2 0 0 0 3'
phi0 = 0.1
n = 1.0
m = 2.0
k0 = 2
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[]
[]
[Preconditioning]
active = check
[andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'bcgs bjacobi 1E-15 1E-10 10000'
[]
[check]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -snes_type'
petsc_options_value = 'bcgs bjacobi 1E-15 1E-10 10000 test'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]
(modules/porous_flow/test/tests/poroperm/PermFromPoro02.i)
# Testing permeability from porosity
# Trivial test, checking calculated permeability is correct
# k = k_anisotropic * k0 * (1-phi0)^m/phi0^n * phi^n/(1-phi)^m
[Mesh]
type = GeneratedMesh
dim = 1
nx = 3
xmin = 0
xmax = 3
[]
[GlobalParams]
block = 0
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[InitialCondition]
type = ConstantIC
value = 0
[]
[]
[]
[Kernels]
[flux]
type = PorousFlowAdvectiveFlux
gravity = '0 0 0'
variable = pp
[]
[]
[BCs]
[ptop]
type = DirichletBC
variable = pp
boundary = right
value = 0
[]
[pbase]
type = DirichletBC
variable = pp
boundary = left
value = 1
[]
[]
[AuxVariables]
[poro]
order = CONSTANT
family = MONOMIAL
[]
[perm_x]
order = CONSTANT
family = MONOMIAL
[]
[perm_y]
order = CONSTANT
family = MONOMIAL
[]
[perm_z]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[poro]
type = PorousFlowPropertyAux
property = porosity
variable = poro
[]
[perm_x]
type = PorousFlowPropertyAux
property = permeability
variable = perm_x
row = 0
column = 0
[]
[perm_y]
type = PorousFlowPropertyAux
property = permeability
variable = perm_y
row = 1
column = 1
[]
[perm_z]
type = PorousFlowPropertyAux
property = permeability
variable = perm_z
row = 2
column = 2
[]
[]
[Postprocessors]
[perm_x_bottom]
type = PointValue
variable = perm_x
point = '0 0 0'
[]
[perm_y_bottom]
type = PointValue
variable = perm_y
point = '0 0 0'
[]
[perm_z_bottom]
type = PointValue
variable = perm_z
point = '0 0 0'
[]
[perm_x_top]
type = PointValue
variable = perm_x
point = '3 0 0'
[]
[perm_y_top]
type = PointValue
variable = perm_y
point = '3 0 0'
[]
[perm_z_top]
type = PointValue
variable = perm_z
point = '3 0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
# unimportant in this fully-saturated test
m = 0.8
alpha = 1e-4
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2.2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
k_anisotropy = '1 0 0 0 2 0 0 0 0.1'
poroperm_function = kozeny_carman_phi0
k0 = 1e-10
phi0 = 0.05
m = 2
n = 7
[]
[temperature]
type = PorousFlowTemperature
[]
[massfrac]
type = PorousFlowMassFraction
[]
[eff_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosity
porosity_zero = 0.1
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 0 # unimportant in this fully-saturated situation
phase = 0
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = Newton
type = Steady
l_tol = 1E-5
nl_abs_tol = 1E-3
nl_rel_tol = 1E-8
l_max_its = 200
nl_max_its = 400
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -ksp_type -ksp_gmres_restart'
petsc_options_value = ' asm 2 lu gmres 200'
[]
[Outputs]
csv = true
execute_on = 'timestep_end'
[]
(modules/porous_flow/examples/multiapp_fracture_flow/3dFracture/fracture_only_aperture_changing.i)
# Cold water injection into one side of the fracture network, and production from the other side
injection_rate = 10 # kg/s
[Mesh]
uniform_refine = 0
[cluster34]
type = FileMeshGenerator
file = 'Cluster_34.exo'
[]
[injection_node]
type = BoundingBoxNodeSetGenerator
input = cluster34
bottom_left = '-1000 0 -1000'
top_right = '1000 0.504 1000'
new_boundary = injection_node
[]
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 -9.81E-6' # Note the value, because of pressure_unit
[]
[Variables]
[frac_P]
scaling = 1E6
[]
[frac_T]
initial_condition = 473
[]
[]
[ICs]
[frac_P]
type = FunctionIC
variable = frac_P
function = insitu_pp
[]
[]
[PorousFlowFullySaturated]
coupling_type = ThermoHydro
porepressure = frac_P
temperature = frac_T
fp = water
pressure_unit = MPa
[]
[Kernels]
[toMatrix]
type = PorousFlowHeatMassTransfer
variable = frac_T
v = transferred_matrix_T
transfer_coefficient = heat_transfer_coefficient
save_in = joules_per_s
[]
[]
[AuxVariables]
[heat_transfer_coefficient]
family = MONOMIAL
order = CONSTANT
initial_condition = 0.0
[]
[transferred_matrix_T]
initial_condition = 473
[]
[joules_per_s]
[]
[normal_dirn_x]
family = MONOMIAL
order = CONSTANT
[]
[normal_dirn_y]
family = MONOMIAL
order = CONSTANT
[]
[normal_dirn_z]
family = MONOMIAL
order = CONSTANT
[]
[enclosing_element_normal_length]
family = MONOMIAL
order = CONSTANT
[]
[enclosing_element_normal_thermal_cond]
family = MONOMIAL
order = CONSTANT
[]
[aperture]
family = MONOMIAL
order = CONSTANT
[]
[perm_times_app]
family = MONOMIAL
order = CONSTANT
[]
[density]
family = MONOMIAL
order = CONSTANT
[]
[viscosity]
family = MONOMIAL
order = CONSTANT
[]
[insitu_pp]
[]
[]
[AuxKernels]
[normal_dirn_x_auxk]
type = PorousFlowElementNormal
variable = normal_dirn_x
component = x
[]
[normal_dirn_y]
type = PorousFlowElementNormal
variable = normal_dirn_y
component = y
[]
[normal_dirn_z]
type = PorousFlowElementNormal
variable = normal_dirn_z
component = z
[]
[heat_transfer_coefficient_auxk]
type = ParsedAux
variable = heat_transfer_coefficient
coupled_variables = 'enclosing_element_normal_length enclosing_element_normal_thermal_cond'
constant_names = h_s
constant_expressions = 1E3 # should be much bigger than thermal_conductivity / L ~ 1
expression = 'if(enclosing_element_normal_length = 0, 0, h_s * enclosing_element_normal_thermal_cond * 2 * enclosing_element_normal_length / (h_s * enclosing_element_normal_length * enclosing_element_normal_length + enclosing_element_normal_thermal_cond * 2 * enclosing_element_normal_length))'
[]
[aperture]
type = PorousFlowPropertyAux
variable = aperture
property = porosity
[]
[perm_times_app]
type = PorousFlowPropertyAux
variable = perm_times_app
property = permeability
row = 0
column = 0
[]
[density]
type = PorousFlowPropertyAux
variable = density
property = density
phase = 0
[]
[viscosity]
type = PorousFlowPropertyAux
variable = viscosity
property = viscosity
phase = 0
[]
[insitu_pp]
type = FunctionAux
execute_on = initial
variable = insitu_pp
function = insitu_pp
[]
[]
[BCs]
[inject_heat]
type = DirichletBC
boundary = injection_node
variable = frac_T
value = 373
[]
[]
[DiracKernels]
[inject_fluid]
type = PorousFlowPointSourceFromPostprocessor
mass_flux = ${injection_rate}
point = '58.8124 0.50384 74.7838'
variable = frac_P
[]
[withdraw_fluid]
type = PorousFlowPeacemanBorehole
SumQuantityUO = kg_out_uo
bottom_p_or_t = 10.6 # 1MPa + approx insitu at production point, to prevent aperture closing due to low porepressures
character = 1
line_length = 1
point_file = production.xyz
unit_weight = '0 0 0'
fluid_phase = 0
use_mobility = true
variable = frac_P
[]
[withdraw_heat]
type = PorousFlowPeacemanBorehole
SumQuantityUO = J_out_uo
bottom_p_or_t = 10.6 # 1MPa + approx insitu at production point, to prevent aperture closing due to low porepressures
character = 1
line_length = 1
point_file = production.xyz
unit_weight = '0 0 0'
fluid_phase = 0
use_mobility = true
use_enthalpy = true
variable = frac_T
[]
[]
[UserObjects]
[kg_out_uo]
type = PorousFlowSumQuantity
[]
[J_out_uo]
type = PorousFlowSumQuantity
[]
[]
[FluidProperties]
[true_water]
type = Water97FluidProperties
[]
[water]
type = TabulatedBicubicFluidProperties
fp = true_water
temperature_min = 275 # K
temperature_max = 600
interpolated_properties = 'density viscosity enthalpy internal_energy'
fluid_property_file = water97_tabulated.csv
[]
[]
[Materials]
[porosity]
type = PorousFlowPorosityLinear
porosity_ref = 1E-4 # fracture porosity = 1.0, but must include fracture aperture of 1E-4 at P = insitu_pp
P_ref = insitu_pp
P_coeff = 1E-3 # this is in metres/MPa, ie for P_ref = 1/P_coeff, the aperture becomes 1 metre
porosity_min = 1E-5
[]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
k0 = 1E-15 # fracture perm = 1E-11 m^2, but must include fracture aperture of 1E-4
poroperm_function = kozeny_carman_phi0
m = 0
n = 3
phi0 = 1E-4
[]
[internal_energy]
type = PorousFlowMatrixInternalEnergy
density = 2700 # kg/m^3
specific_heat_capacity = 0 # basically no rock inside the fracture
[]
[aq_thermal_conductivity]
type = PorousFlowThermalConductivityIdeal
dry_thermal_conductivity = '0.6E-4 0 0 0 0.6E-4 0 0 0 0.6E-4' # thermal conductivity of water times fracture aperture. This should increase linearly with aperture, but is set constant in this model
[]
[]
[Functions]
[kg_rate]
type = ParsedFunction
symbol_values = 'dt kg_out'
symbol_names = 'dt kg_out'
expression = 'kg_out/dt'
[]
[insitu_pp]
type = ParsedFunction
expression = '10 - 0.847E-2 * z' # Approximate hydrostatic in MPa
[]
[]
[Postprocessors]
[dt]
type = TimestepSize
outputs = 'none'
[]
[kg_out]
type = PorousFlowPlotQuantity
uo = kg_out_uo
[]
[kg_per_s]
type = FunctionValuePostprocessor
function = kg_rate
[]
[J_out]
type = PorousFlowPlotQuantity
uo = J_out_uo
[]
[TK_out]
type = PointValue
variable = frac_T
point = '101.705 160.459 39.5722'
[]
[P_out]
type = PointValue
variable = frac_P
point = '101.705 160.459 39.5722'
[]
[P_in]
type = PointValue
variable = frac_P
point = '58.8124 0.50384 74.7838'
[]
[]
[VectorPostprocessors]
[heat_transfer_rate]
type = NodalValueSampler
outputs = none
sort_by = id
variable = joules_per_s
[]
[]
[Preconditioning]
[entire_jacobian]
type = SMP
full = true
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = ' asm lu NONZERO 2 '
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
[TimeStepper]
type = IterationAdaptiveDT
dt = 1
optimal_iterations = 10
growth_factor = 1.5
[]
dtmax = 1E8
end_time = 1E8
nl_abs_tol = 1E-3
nl_max_its = 20
[]
[Outputs]
print_linear_residuals = false
csv = true
[ex]
type = Exodus
sync_times = '1 10 100 200 300 400 500 600 700 800 900 1000 1100 1200 1300 1400 1500 1600 1700 1800 1900 2000 2100 2200 2300 2400 2500 2600 2700 2800 2900 3000 3100 3200 3300 3400 3500 3600 3700 3800 3900 4000 4100 4200 4300 4400 4500 4600 4700 4800 4900 5000 5100 5200 5300 5400 5500 5600 5700 5800 5900 6000 6100 6200 6300 6400 6500 6600 6700 6800 6900 7000 7100 7200 7300 7400 7500 7600 7700 7800 7900 8000 8100 8200 8300 8400 8500 8600 8700 8800 8900 9000 10000 11000 12000 13000 14000 15000 16000 17000 18000 19000 20000 30000 50000 70000 100000 200000 300000 400000 500000 600000 700000 800000 900000 1000000 1100000 1200000 1300000 1400000 1500000 1600000 1700000 1800000 1900000 2000000 2100000 2200000 2300000 2400000 2500000 2600000 2700000 2800000 2900000'
sync_only = true
[]
[]
(modules/porous_flow/test/tests/jacobian/basic_advection5.i)
# Basic advection with 1 porepressure as a PorousFlow variable
[Mesh]
type = GeneratedMesh
dim = 1
nx = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[u]
[]
[P]
[]
[]
[ICs]
[P]
type = RandomIC
variable = P
min = -1
max = 1
[]
[u]
type = RandomIC
variable = u
[]
[]
[Kernels]
[dummy_P]
type = NullKernel
variable = P
[]
[u_advection]
type = PorousFlowBasicAdvection
variable = u
phase = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = P
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 1
m = 0.6
sat_lr = 0.1
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 3
density0 = 4
thermal_expansion = 0
viscosity = 150.0
[]
[]
[Materials]
[temperature_qp]
type = PorousFlowTemperature
[]
[ppss_qp]
type = PorousFlow1PhaseP
porepressure = P
capillary_pressure = pc
[]
[simple_fluid_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[effective_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[porosity]
type = PorousFlowPorosity
porosity_zero = 0.1
fluid = true
biot_coefficient = 0.5
solid_bulk = 1
[]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
poroperm_function = kozeny_carman_phi0
k0 = 5
m = 2
n = 2
phi0 = 0.1
[]
[relperm_qp]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[]
[darcy_velocity_qp]
type = PorousFlowDarcyVelocityMaterial
gravity = '0.25 0 0'
[]
[]
[Preconditioning]
[check]
type = SMP
full = true
#petsc_options = '-snes_test_display'
petsc_options_iname = '-snes_type'
petsc_options_value = ' test'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1
[]
(modules/porous_flow/test/tests/poro_elasticity/vol_expansion_poroperm.i)
# Apply an increasing porepressure, with zero mechanical forces,
# and observe the corresponding volumetric expansion and porosity increase.
# Check that permeability is calculated correctly from porosity.
#
# P = t
# With the Biot coefficient being 1, the effective stresses should be
# stress_xx = stress_yy = stress_zz = t
# With bulk modulus = 1 then should have
# vol_strain = strain_xx + strain_yy + strain_zz = t.
#
# With the biot coefficient being 1, the porosity (phi) # at time t is:
# phi = 1 - (1 - phi0) / exp(vol_strain)
# where phi0 is the porosity at t = 0 and P = 0.
#
# The permeability (k) is
# k = k_anisotropic * f * d^2 * phi^n / (1-phi)^m
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 1
zmin = 0
zmax = 1
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
block = 0
PorousFlowDictator = dictator
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[p]
[]
[]
[BCs]
[p]
type = FunctionDirichletBC
boundary = 'bottom top'
variable = p
function = t
[]
[xmin]
type = DirichletBC
boundary = left
variable = disp_x
value = 0
[]
[ymin]
type = DirichletBC
boundary = bottom
variable = disp_y
value = 0
[]
[zmin]
type = DirichletBC
boundary = back
variable = disp_z
value = 0
[]
[]
[Kernels]
[p_does_not_really_diffuse]
type = Diffusion
variable = p
[]
[TensorMechanics]
displacements = 'disp_x disp_y disp_z'
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 1
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 1
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 1
variable = disp_z
component = 2
[]
[]
[AuxVariables]
[poro]
order = CONSTANT
family = MONOMIAL
[]
[perm_x]
order = CONSTANT
family = MONOMIAL
[]
[perm_y]
order = CONSTANT
family = MONOMIAL
[]
[perm_z]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[poro]
type = PorousFlowPropertyAux
property = porosity
variable = poro
[]
[perm_x]
type = PorousFlowPropertyAux
property = permeability
variable = perm_x
row = 0
column = 0
[]
[perm_y]
type = PorousFlowPropertyAux
property = permeability
variable = perm_y
row = 1
column = 1
[]
[perm_z]
type = PorousFlowPropertyAux
property = permeability
variable = perm_z
row = 2
column = 2
[]
[]
[Postprocessors]
[poro]
type = PointValue
variable = poro
point = '0 0 0'
[]
[perm_x]
type = PointValue
variable = perm_x
point = '0 0 0'
[]
[perm_y]
type = PointValue
variable = perm_y
point = '0 0 0'
[]
[perm_z]
type = PointValue
variable = perm_z
point = '0 0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'p'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[elasticity_tensor]
type = ComputeIsotropicElasticityTensor
bulk_modulus = 1
shear_modulus = 1
[]
[strain]
type = ComputeSmallStrain
[]
[stress]
type = ComputeLinearElasticStress
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = p
capillary_pressure = pc
[]
[p_eff]
type = PorousFlowEffectiveFluidPressure
[]
[porosity]
type = PorousFlowPorosity
fluid = true
mechanical = true
porosity_zero = 0.1
solid_bulk = 1
biot_coefficient = 1
[]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
k_anisotropy = '1 0 0 0 2 0 0 0 0.1'
poroperm_function = kozeny_carman_fd2
f = 0.1
d = 5
m = 2
n = 7
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_atol -ksp_rtol'
petsc_options_value = 'gmres bjacobi 1E-10 1E-10 10 1E-15 1E-10'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
start_time = 0
dt = 0.1
end_time = 1
[]
[Outputs]
file_base = vol_expansion_poroperm
csv = true
execute_on = 'timestep_end'
[]
(modules/porous_flow/test/tests/jacobian/basic_advection6.i)
# Basic advection with 2 porepressure as PorousFlow variables
[Mesh]
type = GeneratedMesh
dim = 1
nx = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[u]
[]
[P0]
[]
[P1]
[]
[]
[ICs]
[P0]
type = RandomIC
variable = P0
min = -1
max = 0
[]
[P1]
type = RandomIC
variable = P1
min = 0
max = 1
[]
[u]
type = RandomIC
variable = u
[]
[]
[Kernels]
[dummy_P0]
type = NullKernel
variable = P0
[]
[dummy_P1]
type = NullKernel
variable = P1
[]
[u_advection]
type = PorousFlowBasicAdvection
variable = u
phase = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'P0 P1'
number_fluid_phases = 2
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 1
m = 0.6
sat_lr = 0.1
[]
[]
[FluidProperties]
[simple_fluid0]
type = SimpleFluidProperties
bulk_modulus = 3
density0 = 4
thermal_expansion = 0
viscosity = 150.0
[]
[simple_fluid1]
type = SimpleFluidProperties
bulk_modulus = 4
density0 = 3
thermal_expansion = 0
viscosity = 130.0
[]
[]
[Materials]
[temperature_qp]
type = PorousFlowTemperature
[]
[ppss_qp]
type = PorousFlow2PhasePP
phase0_porepressure = P0
phase1_porepressure = P1
capillary_pressure = pc
[]
[simple_fluid0_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid0
phase = 0
[]
[simple_fluid1_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid1
phase = 1
[]
[effective_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[porosity]
type = PorousFlowPorosity
porosity_zero = 0.1
fluid = true
biot_coefficient = 0.5
solid_bulk = 1
[]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
poroperm_function = kozeny_carman_phi0
k0 = 5
m = 2
n = 2
phi0 = 0.1
[]
[relperm0_qp]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
s_res = 0.1
sum_s_res = 0.1
[]
[relperm1_qp]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
s_res = 0.0
sum_s_res = 0.1
[]
[darcy_velocity_qp]
type = PorousFlowDarcyVelocityMaterial
gravity = '0.25 0 0'
[]
[]
[Preconditioning]
[check]
type = SMP
full = true
petsc_options = '-snes_test_display'
petsc_options_iname = '-snes_type'
petsc_options_value = ' test'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1
[]
(modules/porous_flow/test/tests/poroperm/PermFromPoro01.i)
# Testing permeability from porosity
# Trivial test, checking calculated permeability is correct
# k = k_anisotropic * f * d^2 * phi^n / (1-phi)^m
[Mesh]
type = GeneratedMesh
dim = 1
nx = 3
xmin = 0
xmax = 3
[]
[GlobalParams]
block = 0
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[InitialCondition]
type = ConstantIC
value = 0
[]
[]
[]
[Kernels]
[flux]
type = PorousFlowAdvectiveFlux
gravity = '0 0 0'
variable = pp
[]
[]
[BCs]
[ptop]
type = DirichletBC
variable = pp
boundary = right
value = 0
[]
[pbase]
type = DirichletBC
variable = pp
boundary = left
value = 1
[]
[]
[AuxVariables]
[poro]
order = CONSTANT
family = MONOMIAL
[]
[perm_x]
order = CONSTANT
family = MONOMIAL
[]
[perm_y]
order = CONSTANT
family = MONOMIAL
[]
[perm_z]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[poro]
type = PorousFlowPropertyAux
property = porosity
variable = poro
[]
[perm_x]
type = PorousFlowPropertyAux
property = permeability
variable = perm_x
row = 0
column = 0
[]
[perm_y]
type = PorousFlowPropertyAux
property = permeability
variable = perm_y
row = 1
column = 1
[]
[perm_z]
type = PorousFlowPropertyAux
property = permeability
variable = perm_z
row = 2
column = 2
[]
[]
[Postprocessors]
[perm_x_bottom]
type = PointValue
variable = perm_x
point = '0 0 0'
[]
[perm_y_bottom]
type = PointValue
variable = perm_y
point = '0 0 0'
[]
[perm_z_bottom]
type = PointValue
variable = perm_z
point = '0 0 0'
[]
[perm_x_top]
type = PointValue
variable = perm_x
point = '3 0 0'
[]
[perm_y_top]
type = PointValue
variable = perm_y
point = '3 0 0'
[]
[perm_z_top]
type = PointValue
variable = perm_z
point = '3 0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
# unimportant in this fully-saturated test
m = 0.8
alpha = 1e-4
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2.2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[permeability]
type = PorousFlowPermeabilityKozenyCarman
k_anisotropy = '1 0 0 0 2 0 0 0 0.1'
poroperm_function = kozeny_carman_fd2
f = 0.1
d = 5
m = 2
n = 7
[]
[temperature]
type = PorousFlowTemperature
[]
[massfrac]
type = PorousFlowMassFraction
[]
[eff_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosity
porosity_zero = 0.1
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 0 # unimportant in this fully-saturated situation
phase = 0
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = Newton
type = Steady
l_tol = 1E-5
nl_abs_tol = 1E-3
nl_rel_tol = 1E-8
l_max_its = 200
nl_max_its = 400
petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -ksp_type -ksp_gmres_restart'
petsc_options_value = ' asm 2 lu gmres 200'
[]
[Outputs]
csv = true
execute_on = 'timestep_end'
[]
(modules/porous_flow/examples/coal_mining/coarse_with_fluid.i)
# Strata deformation and fluid flow aaround a coal mine - 3D model
#
# A "half model" is used. The mine is 400m deep and
# just the roof is studied (-400<=z<=0). The mining panel
# sits between 0<=x<=150, and 0<=y<=1000, so this simulates
# a coal panel that is 300m wide and 1000m long. The outer boundaries
# are 1km from the excavation boundaries.
#
# The excavation takes 0.5 years.
#
# The boundary conditions for this simulation are:
# - disp_x = 0 at x=0 and x=1150
# - disp_y = 0 at y=-1000 and y=1000
# - disp_z = 0 at z=-400, but there is a time-dependent
# Young modulus that simulates excavation
# - wc_x = 0 at y=-1000 and y=1000
# - wc_y = 0 at x=0 and x=1150
# - no flow at x=0, z=-400 and z=0
# - fixed porepressure at y=-1000, y=1000 and x=1150
# That is, rollers on the sides, free at top,
# and prescribed at bottom in the unexcavated portion.
#
# A single-phase unsaturated fluid is used.
#
# The small strain formulation is used.
#
# All stresses are measured in MPa, and time units are measured in years.
#
# The initial porepressure is hydrostatic with P=0 at z=0, so
# Porepressure ~ - 0.01*z MPa, where the fluid has density 1E3 kg/m^3 and
# gravity = = 10 m.s^-2 = 1E-5 MPa m^2/kg.
# To be more accurate, i use
# Porepressure = -bulk * log(1 + g*rho0*z/bulk)
# where bulk=2E3 MPa and rho0=1Ee kg/m^3.
# The initial stress is consistent with the weight force from undrained
# density 2500 kg/m^3, and fluid porepressure, and a Biot coefficient of 0.7, ie,
# stress_zz^effective = 0.025*z + 0.7 * initial_porepressure
# The maximum and minimum principal horizontal effective stresses are
# assumed to be equal to 0.8*stress_zz.
#
# Material properties:
# Young's modulus = 8 GPa
# Poisson's ratio = 0.25
# Cosserat layer thickness = 1 m
# Cosserat-joint normal stiffness = large
# Cosserat-joint shear stiffness = 1 GPa
# MC cohesion = 2 MPa
# MC friction angle = 35 deg
# MC dilation angle = 8 deg
# MC tensile strength = 1 MPa
# MC compressive strength = 100 MPa
# WeakPlane cohesion = 0.1 MPa
# WeakPlane friction angle = 30 deg
# WeakPlane dilation angle = 10 deg
# WeakPlane tensile strength = 0.1 MPa
# WeakPlane compressive strength = 100 MPa softening to 1 MPa at strain = 1
# Fluid density at zero porepressure = 1E3 kg/m^3
# Fluid bulk modulus = 2E3 MPa
# Fluid viscosity = 1.1E-3 Pa.s = 1.1E-9 MPa.s = 3.5E-17 MPa.year
#
[GlobalParams]
perform_finite_strain_rotations = false
displacements = 'disp_x disp_y disp_z'
Cosserat_rotations = 'wc_x wc_y wc_z'
PorousFlowDictator = dictator
biot_coefficient = 0.7
[]
[Mesh]
[file]
type = FileMeshGenerator
file = mesh/coarse.e
[]
[xmin]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
new_boundary = xmin
normal = '-1 0 0'
input = file
[]
[xmax]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
new_boundary = xmax
normal = '1 0 0'
input = xmin
[]
[ymin]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
new_boundary = ymin
normal = '0 -1 0'
input = xmax
[]
[ymax]
type = SideSetsAroundSubdomainGenerator
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
new_boundary = ymax
normal = '0 1 0'
input = ymin
[]
[zmax]
type = SideSetsAroundSubdomainGenerator
block = 16
new_boundary = zmax
normal = '0 0 1'
input = ymax
[]
[zmin]
type = SideSetsAroundSubdomainGenerator
block = 2
new_boundary = zmin
normal = '0 0 -1'
input = zmax
[]
[excav]
type = SubdomainBoundingBoxGenerator
input = zmin
block_id = 1
bottom_left = '0 0 -400'
top_right = '150 1000 -397'
[]
[roof]
type = SideSetsBetweenSubdomainsGenerator
primary_block = 3
paired_block = 1
input = excav
new_boundary = roof
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[wc_x]
[]
[wc_y]
[]
[porepressure]
scaling = 1E-5
[]
[]
[ICs]
[porepressure]
type = FunctionIC
variable = porepressure
function = ini_pp
[]
[]
[Kernels]
[cx_elastic]
type = CosseratStressDivergenceTensors
use_displaced_mesh = false
variable = disp_x
component = 0
[]
[cy_elastic]
type = CosseratStressDivergenceTensors
use_displaced_mesh = false
variable = disp_y
component = 1
[]
[cz_elastic]
type = CosseratStressDivergenceTensors
use_displaced_mesh = false
variable = disp_z
component = 2
[]
[x_couple]
type = StressDivergenceTensors
use_displaced_mesh = false
variable = wc_x
displacements = 'wc_x wc_y wc_z'
component = 0
base_name = couple
[]
[y_couple]
type = StressDivergenceTensors
use_displaced_mesh = false
variable = wc_y
displacements = 'wc_x wc_y wc_z'
component = 1
base_name = couple
[]
[x_moment]
type = MomentBalancing
use_displaced_mesh = false
variable = wc_x
component = 0
[]
[y_moment]
type = MomentBalancing
use_displaced_mesh = false
variable = wc_y
component = 1
[]
[gravity]
type = Gravity
use_displaced_mesh = false
variable = disp_z
value = -10E-6 # remember this is in MPa
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
use_displaced_mesh = false
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
use_displaced_mesh = false
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
use_displaced_mesh = false
component = 2
variable = disp_z
[]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = porepressure
[]
[flux]
type = PorousFlowAdvectiveFlux
use_displaced_mesh = false
variable = porepressure
gravity = '0 0 -10E-6'
fluid_component = 0
[]
[poro_vol_exp]
type = PorousFlowMassVolumetricExpansion
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
variable = porepressure
fluid_component = 0
[]
[]
[AuxVariables]
[saturation]
order = CONSTANT
family = MONOMIAL
[]
[darcy_x]
order = CONSTANT
family = MONOMIAL
[]
[darcy_y]
order = CONSTANT
family = MONOMIAL
[]
[darcy_z]
order = CONSTANT
family = MONOMIAL
[]
[porosity]
order = CONSTANT
family = MONOMIAL
[]
[wc_z]
[]
[stress_xx]
order = CONSTANT
family = MONOMIAL
[]
[stress_xy]
order = CONSTANT
family = MONOMIAL
[]
[stress_xz]
order = CONSTANT
family = MONOMIAL
[]
[stress_yx]
order = CONSTANT
family = MONOMIAL
[]
[stress_yy]
order = CONSTANT
family = MONOMIAL
[]
[stress_yz]
order = CONSTANT
family = MONOMIAL
[]
[stress_zx]
order = CONSTANT
family = MONOMIAL
[]
[stress_zy]
order = CONSTANT
family = MONOMIAL
[]
[stress_zz]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_xx]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_xy]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_xz]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_yx]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_yy]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_yz]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_zx]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_zy]
order = CONSTANT
family = MONOMIAL
[]
[total_strain_zz]
order = CONSTANT
family = MONOMIAL
[]
[perm_xx]
order = CONSTANT
family = MONOMIAL
[]
[perm_yy]
order = CONSTANT
family = MONOMIAL
[]
[perm_zz]
order = CONSTANT
family = MONOMIAL
[]
[mc_shear]
order = CONSTANT
family = MONOMIAL
[]
[mc_tensile]
order = CONSTANT
family = MONOMIAL
[]
[wp_shear]
order = CONSTANT
family = MONOMIAL
[]
[wp_tensile]
order = CONSTANT
family = MONOMIAL
[]
[wp_shear_f]
order = CONSTANT
family = MONOMIAL
[]
[wp_tensile_f]
order = CONSTANT
family = MONOMIAL
[]
[mc_shear_f]
order = CONSTANT
family = MONOMIAL
[]
[mc_tensile_f]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[saturation_water]
type = PorousFlowPropertyAux
variable = saturation
property = saturation
phase = 0
execute_on = timestep_end
[]
[darcy_x]
type = PorousFlowDarcyVelocityComponent
variable = darcy_x
gravity = '0 0 -10E-6'
component = x
[]
[darcy_y]
type = PorousFlowDarcyVelocityComponent
variable = darcy_y
gravity = '0 0 -10E-6'
component = y
[]
[darcy_z]
type = PorousFlowDarcyVelocityComponent
variable = darcy_z
gravity = '0 0 -10E-6'
component = z
[]
[porosity]
type = PorousFlowPropertyAux
property = porosity
variable = porosity
execute_on = timestep_end
[]
[stress_xx]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_xx
index_i = 0
index_j = 0
execute_on = timestep_end
[]
[stress_xy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_xy
index_i = 0
index_j = 1
execute_on = timestep_end
[]
[stress_xz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_xz
index_i = 0
index_j = 2
execute_on = timestep_end
[]
[stress_yx]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yx
index_i = 1
index_j = 0
execute_on = timestep_end
[]
[stress_yy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yy
index_i = 1
index_j = 1
execute_on = timestep_end
[]
[stress_yz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yz
index_i = 1
index_j = 2
execute_on = timestep_end
[]
[stress_zx]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zx
index_i = 2
index_j = 0
execute_on = timestep_end
[]
[stress_zy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zy
index_i = 2
index_j = 1
execute_on = timestep_end
[]
[stress_zz]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_zz
index_i = 2
index_j = 2
execute_on = timestep_end
[]
[total_strain_xx]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_xx
index_i = 0
index_j = 0
execute_on = timestep_end
[]
[total_strain_xy]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_xy
index_i = 0
index_j = 1
execute_on = timestep_end
[]
[total_strain_xz]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_xz
index_i = 0
index_j = 2
execute_on = timestep_end
[]
[total_strain_yx]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_yx
index_i = 1
index_j = 0
execute_on = timestep_end
[]
[total_strain_yy]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_yy
index_i = 1
index_j = 1
execute_on = timestep_end
[]
[total_strain_yz]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_yz
index_i = 1
index_j = 2
execute_on = timestep_end
[]
[total_strain_zx]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_zx
index_i = 2
index_j = 0
execute_on = timestep_end
[]
[total_strain_zy]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_zy
index_i = 2
index_j = 1
execute_on = timestep_end
[]
[total_strain_zz]
type = RankTwoAux
rank_two_tensor = total_strain
variable = total_strain_zz
index_i = 2
index_j = 2
execute_on = timestep_end
[]
[perm_xx]
type = PorousFlowPropertyAux
property = permeability
variable = perm_xx
row = 0
column = 0
execute_on = timestep_end
[]
[perm_yy]
type = PorousFlowPropertyAux
property = permeability
variable = perm_yy
row = 1
column = 1
execute_on = timestep_end
[]
[perm_zz]
type = PorousFlowPropertyAux
property = permeability
variable = perm_zz
row = 2
column = 2
execute_on = timestep_end
[]
[mc_shear]
type = MaterialStdVectorAux
index = 0
property = mc_plastic_internal_parameter
variable = mc_shear
execute_on = timestep_end
[]
[mc_tensile]
type = MaterialStdVectorAux
index = 1
property = mc_plastic_internal_parameter
variable = mc_tensile
execute_on = timestep_end
[]
[wp_shear]
type = MaterialStdVectorAux
index = 0
property = wp_plastic_internal_parameter
variable = wp_shear
execute_on = timestep_end
[]
[wp_tensile]
type = MaterialStdVectorAux
index = 1
property = wp_plastic_internal_parameter
variable = wp_tensile
execute_on = timestep_end
[]
[mc_shear_f]
type = MaterialStdVectorAux
index = 6
property = mc_plastic_yield_function
variable = mc_shear_f
execute_on = timestep_end
[]
[mc_tensile_f]
type = MaterialStdVectorAux
index = 0
property = mc_plastic_yield_function
variable = mc_tensile_f
execute_on = timestep_end
[]
[wp_shear_f]
type = MaterialStdVectorAux
index = 0
property = wp_plastic_yield_function
variable = wp_shear_f
execute_on = timestep_end
[]
[wp_tensile_f]
type = MaterialStdVectorAux
index = 1
property = wp_plastic_yield_function
variable = wp_tensile_f
execute_on = timestep_end
[]
[]
[BCs]
[no_x]
type = DirichletBC
variable = disp_x
boundary = 'xmin xmax'
value = 0.0
[]
[no_y]
type = DirichletBC
variable = disp_y
boundary = 'ymin ymax'
value = 0.0
[]
[no_z]
type = DirichletBC
variable = disp_z
boundary = zmin
value = 0.0
[]
[no_wc_x]
type = DirichletBC
variable = wc_x
boundary = 'ymin ymax'
value = 0.0
[]
[no_wc_y]
type = DirichletBC
variable = wc_y
boundary = 'xmin xmax'
value = 0.0
[]
[fix_porepressure]
type = FunctionDirichletBC
variable = porepressure
boundary = 'ymin ymax xmax'
function = ini_pp
[]
[roof_porepressure]
type = PorousFlowPiecewiseLinearSink
variable = porepressure
pt_vals = '-1E3 1E3'
multipliers = '-1 1'
fluid_phase = 0
flux_function = roof_conductance
boundary = roof
[]
[roof_bcs]
type = StickyBC
variable = disp_z
min_value = -3.0
boundary = roof
[]
[]
[Functions]
[ini_pp]
type = ParsedFunction
symbol_names = 'bulk p0 g rho0'
symbol_values = '2E3 0.0 1E-5 1E3'
expression = '-bulk*log(exp(-p0/bulk)+g*rho0*z/bulk)'
[]
[ini_xx]
type = ParsedFunction
symbol_names = 'bulk p0 g rho0 biot'
symbol_values = '2E3 0.0 1E-5 1E3 0.7'
expression = '0.8*(2500*10E-6*z+biot*(-bulk*log(exp(-p0/bulk)+g*rho0*z/bulk)))'
[]
[ini_zz]
type = ParsedFunction
symbol_names = 'bulk p0 g rho0 biot'
symbol_values = '2E3 0.0 1E-5 1E3 0.7'
expression = '2500*10E-6*z+biot*(-bulk*log(exp(-p0/bulk)+g*rho0*z/bulk))'
[]
[excav_sideways]
type = ParsedFunction
symbol_names = 'end_t ymin ymax minval maxval slope'
symbol_values = '0.5 0 1000.0 1E-9 1 60'
# excavation face at ymin+(ymax-ymin)*min(t/end_t,1)
# slope is the distance over which the modulus reduces from maxval to minval
expression = 'if(y<ymin+(ymax-ymin)*min(t/end_t,1),minval,if(y<ymin+(ymax-ymin)*min(t/end_t,1)+slope,minval+(maxval-minval)*(y-(ymin+(ymax-ymin)*min(t/end_t,1)))/slope,maxval))'
[]
[density_sideways]
type = ParsedFunction
symbol_names = 'end_t ymin ymax minval maxval'
symbol_values = '0.5 0 1000.0 0 2500'
expression = 'if(y<ymin+(ymax-ymin)*min(t/end_t,1),minval,maxval)'
[]
[roof_conductance]
type = ParsedFunction
symbol_names = 'end_t ymin ymax maxval minval'
symbol_values = '0.5 0 1000.0 1E7 0'
expression = 'if(y<ymin+(ymax-ymin)*min(t/end_t,1),maxval,minval)'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'porepressure disp_x disp_y disp_z'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1 # MPa^-1
[]
[mc_coh_strong_harden]
type = TensorMechanicsHardeningExponential
value_0 = 1.99 # MPa
value_residual = 2.01 # MPa
rate = 1.0
[]
[mc_fric]
type = TensorMechanicsHardeningConstant
value = 0.61 # 35deg
[]
[mc_dil]
type = TensorMechanicsHardeningConstant
value = 0.15 # 8deg
[]
[mc_tensile_str_strong_harden]
type = TensorMechanicsHardeningExponential
value_0 = 1.0 # MPa
value_residual = 1.0 # MPa
rate = 1.0
[]
[mc_compressive_str]
type = TensorMechanicsHardeningCubic
value_0 = 100 # Large!
value_residual = 100
internal_limit = 0.1
[]
[wp_coh_harden]
type = TensorMechanicsHardeningCubic
value_0 = 0.05
value_residual = 0.05
internal_limit = 10
[]
[wp_tan_fric]
type = TensorMechanicsHardeningConstant
value = 0.26 # 15deg
[]
[wp_tan_dil]
type = TensorMechanicsHardeningConstant
value = 0.18 # 10deg
[]
[wp_tensile_str_harden]
type = TensorMechanicsHardeningCubic
value_0 = 0.05
value_residual = 0.05
internal_limit = 10
[]
[wp_compressive_str_soften]
type = TensorMechanicsHardeningCubic
value_0 = 100
value_residual = 1
internal_limit = 1.0
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2E3
density0 = 1000
thermal_expansion = 0
viscosity = 3.5E-17
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[eff_fluid_pressure]
type = PorousFlowEffectiveFluidPressure
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = porepressure
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity_bulk]
type = PorousFlowPorosity
fluid = true
mechanical = true
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
ensure_positive = true
porosity_zero = 0.02
solid_bulk = 5.3333E3
[]
[porosity_excav]
type = PorousFlowPorosityConst
block = 1
porosity = 1.0
[]
[permeability_bulk]
type = PorousFlowPermeabilityKozenyCarman
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
poroperm_function = kozeny_carman_phi0
k0 = 1E-15
phi0 = 0.02
n = 2
m = 2
[]
[permeability_excav]
type = PorousFlowPermeabilityConst
block = 1
permeability = '0 0 0 0 0 0 0 0 0'
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 4
s_res = 0.4
sum_s_res = 0.4
phase = 0
[]
[elasticity_tensor_0]
type = ComputeLayeredCosseratElasticityTensor
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
young = 8E3 # MPa
poisson = 0.25
layer_thickness = 1.0
joint_normal_stiffness = 1E9 # huge
joint_shear_stiffness = 1E3 # MPa
[]
[elasticity_tensor_1]
type = ComputeLayeredCosseratElasticityTensor
block = 1
young = 8E3 # MPa
poisson = 0.25
layer_thickness = 1.0
joint_normal_stiffness = 1E9 # huge
joint_shear_stiffness = 1E3 # MPa
elasticity_tensor_prefactor = excav_sideways
[]
[strain]
type = ComputeCosseratIncrementalSmallStrain
eigenstrain_names = ini_stress
[]
[ini_stress]
type = ComputeEigenstrainFromInitialStress
eigenstrain_name = ini_stress
initial_stress = 'ini_xx 0 0 0 ini_xx 0 0 0 ini_zz'
[]
[stress_0]
type = ComputeMultipleInelasticCosseratStress
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
inelastic_models = 'mc wp'
cycle_models = true
relative_tolerance = 2.0
absolute_tolerance = 1E6
max_iterations = 1
tangent_operator = nonlinear
perform_finite_strain_rotations = false
[]
[stress_1]
type = ComputeMultipleInelasticCosseratStress
block = 1
inelastic_models = ''
relative_tolerance = 2.0
absolute_tolerance = 1E6
max_iterations = 1
tangent_operator = nonlinear
perform_finite_strain_rotations = false
[]
[mc]
type = CappedMohrCoulombCosseratStressUpdate
warn_about_precision_loss = false
host_youngs_modulus = 8E3
host_poissons_ratio = 0.25
base_name = mc
tensile_strength = mc_tensile_str_strong_harden
compressive_strength = mc_compressive_str
cohesion = mc_coh_strong_harden
friction_angle = mc_fric
dilation_angle = mc_dil
max_NR_iterations = 100000
smoothing_tol = 0.1 # MPa # Must be linked to cohesion
yield_function_tol = 1E-9 # MPa. this is essentially the lowest possible without lots of precision loss
perfect_guess = true
min_step_size = 1.0
[]
[wp]
type = CappedWeakPlaneCosseratStressUpdate
warn_about_precision_loss = false
base_name = wp
cohesion = wp_coh_harden
tan_friction_angle = wp_tan_fric
tan_dilation_angle = wp_tan_dil
tensile_strength = wp_tensile_str_harden
compressive_strength = wp_compressive_str_soften
max_NR_iterations = 10000
tip_smoother = 0.05
smoothing_tol = 0.05 # MPa # Note, this must be tied to cohesion, otherwise get no possible return at cone apex
yield_function_tol = 1E-11 # MPa. this is essentially the lowest possible without lots of precision loss
perfect_guess = true
min_step_size = 1.0E-3
[]
[undrained_density_0]
type = GenericConstantMaterial
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
prop_names = density
prop_values = 2500
[]
[undrained_density_1]
type = GenericFunctionMaterial
block = 1
prop_names = density
prop_values = density_sideways
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Postprocessors]
[min_roof_disp]
type = NodalExtremeValue
boundary = roof
value_type = min
variable = disp_z
[]
[min_roof_pp]
type = NodalExtremeValue
boundary = roof
value_type = min
variable = porepressure
[]
[min_surface_disp]
type = NodalExtremeValue
boundary = zmax
value_type = min
variable = disp_z
[]
[min_surface_pp]
type = NodalExtremeValue
boundary = zmax
value_type = min
variable = porepressure
[]
[max_perm_zz]
type = ElementExtremeValue
block = '2 3 4 5 6 7 8 9 10 11 12 13 14 15 16'
variable = perm_zz
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-snes_converged_reason'
# best overall
# petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
# petsc_options_value = ' lu mumps'
# best if you do not have mumps:
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu superlu_dist'
# best if you do not have mumps or superlu_dist:
#petsc_options_iname = '-pc_type -pc_asm_overlap -sub_pc_type -ksp_type -ksp_gmres_restart'
#petsc_options_value = ' asm 2 lu gmres 200'
# very basic:
#petsc_options_iname = '-pc_type -ksp_type -ksp_gmres_restart'
#petsc_options_value = ' bjacobi gmres 200'
line_search = bt
nl_abs_tol = 1e-3
nl_rel_tol = 1e-5
l_max_its = 200
nl_max_its = 30
start_time = 0.0
dt = 0.014706
end_time = 0.014706 #0.5
[]
[Outputs]
time_step_interval = 1
print_linear_residuals = true
exodus = true
csv = true
console = true
[]