- 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
- gravityGravitational acceleration vector downwards (m/s^2)
C++ Type:libMesh::VectorValue<double>
Controllable:No
Description:Gravitational acceleration vector downwards (m/s^2)
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Controllable:No
Description:The name of the variable that this residual object operates on
PorousFlowFullySaturatedDarcyBase
Darcy flux suitable for models involving a fully-saturated, single phase, single component fluid. No upwinding is used
Describes the differential term The nomenclature is described here. This is fully-saturated, single-component, single-phase Darcy flow.
The multiplication by is optional: this is indicated by the parenthases. The reason for this is that the time-derivative part described by PorousFlowFullySaturatedMassTimeDerivative linearises if the total differential equation is not multiplied by density. Please see that page for a full description of the effects of not multiplying by density.
No upwinding is performed, which means many nodal Material properties are not needed and numerical diffusion is reduced. However, the numerics are less well controlled: the whole point of full upwinding is to prevent over-shoots and under-shoots. Other Kernels implement Kuzmin-Turek TVD stabilization.
Input Parameters
- 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
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Controllable:No
Description:The displacements
- multiply_by_densityTrueIf true, then this Kernel is the fluid mass flux. If false, then this Kernel is the fluid volume flux (which is common in poro-mechanics)
Default:True
C++ Type:bool
Controllable:No
Description:If true, then this Kernel is the fluid mass flux. If false, then this Kernel is the fluid volume flux (which is common in poro-mechanics)
- 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.
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.
- diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Controllable:No
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- 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
- save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Controllable:No
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- 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
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- (modules/porous_flow/test/tests/poro_elasticity/mandel_fully_saturated_volume.i)
- (modules/porous_flow/test/tests/heat_advection/heat_advection_1d_fully_saturated.i)
- (modules/porous_flow/test/tests/gravity/fully_saturated_grav01b.i)
- (modules/porous_flow/test/tests/gravity/fully_saturated_grav01a.i)
- (modules/porous_flow/test/tests/poro_elasticity/mandel_fully_saturated.i)
- (modules/porous_flow/test/tests/poro_elasticity/terzaghi_fully_saturated_volume.i)
- (modules/porous_flow/test/tests/pressure_pulse/pressure_pulse_1d_fully_saturated_2.i)
- (modules/porous_flow/test/tests/pressure_pulse/pressure_pulse_1d_fully_saturated.i)
Child Objects
References
No citations exist within this document.(modules/porous_flow/test/tests/poro_elasticity/mandel_fully_saturated_volume.i)
# Mandel's problem of consolodation of a drained medium
# Using the FullySaturatedDarcyBase and FullySaturatedFullySaturatedMassTimeDerivative kernels
# with multiply_by_density = false, so that this problem becomes linear
#
# A sample is in plane strain.
# -a <= x <= a
# -b <= y <= b
# It is squashed with constant force by impermeable, frictionless plattens on its top and bottom surfaces (at y=+/-b)
# Fluid is allowed to leak out from its sides (at x=+/-a)
# The porepressure within the sample is monitored.
#
# As is common in the literature, this is simulated by
# considering the quarter-sample, 0<=x<=a and 0<=y<=b, with
# impermeable, roller BCs at x=0 and y=0 and y=b.
# Porepressure is fixed at zero on x=a.
# Porepressure and displacement are initialised to zero.
# Then the top (y=b) is moved downwards with prescribed velocity,
# so that the total force that is inducing this downwards velocity
# is fixed. The velocity is worked out by solving Mandel's problem
# analytically, and the total force is monitored in the simulation
# to check that it indeed remains constant.
#
# Here are the problem's parameters, and their values:
# Soil width. a = 1
# Soil height. b = 0.1
# Soil's Lame lambda. la = 0.5
# Soil's Lame mu, which is also the Soil's shear modulus. mu = G = 0.75
# Soil bulk modulus. K = la + 2*mu/3 = 1
# Drained Poisson ratio. nu = (3K - 2G)/(6K + 2G) = 0.2
# Soil bulk compliance. 1/K = 1
# Fluid bulk modulus. Kf = 8
# Fluid bulk compliance. 1/Kf = 0.125
# Soil initial porosity. phi0 = 0.1
# Biot coefficient. alpha = 0.6
# Biot modulus. M = 1/(phi0/Kf + (alpha - phi0)(1 - alpha)/K) = 4.705882
# Undrained bulk modulus. Ku = K + alpha^2*M = 2.694118
# Undrained Poisson ratio. nuu = (3Ku - 2G)/(6Ku + 2G) = 0.372627
# Skempton coefficient. B = alpha*M/Ku = 1.048035
# Fluid mobility (soil permeability/fluid viscosity). k = 1.5
# Consolidation coefficient. c = 2*k*B^2*G*(1-nu)*(1+nuu)^2/9/(1-nuu)/(nuu-nu) = 3.821656
# Normal stress on top. F = 1
#
# The solution for porepressure and displacements is given in
# AHD Cheng and E Detournay "A direct boundary element method for plane strain poroelasticity" International Journal of Numerical and Analytical Methods in Geomechanics 12 (1988) 551-572.
# The solution involves complicated infinite series, so I shall not write it here
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
ny = 1
nz = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 0.1
zmin = 0
zmax = 1
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
PorousFlowDictator = dictator
block = 0
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'porepressure disp_x disp_y disp_z'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[porepressure]
[]
[]
[BCs]
[roller_xmin]
type = DirichletBC
variable = disp_x
value = 0
boundary = 'left'
[]
[roller_ymin]
type = DirichletBC
variable = disp_y
value = 0
boundary = 'bottom'
[]
[plane_strain]
type = DirichletBC
variable = disp_z
value = 0
boundary = 'back front'
[]
[xmax_drained]
type = DirichletBC
variable = porepressure
value = 0
boundary = right
[]
[top_velocity]
type = FunctionDirichletBC
variable = disp_y
function = top_velocity
boundary = top
[]
[]
[Functions]
[top_velocity]
type = PiecewiseLinear
x = '0 0.002 0.006 0.014 0.03 0.046 0.062 0.078 0.094 0.11 0.126 0.142 0.158 0.174 0.19 0.206 0.222 0.238 0.254 0.27 0.286 0.302 0.318 0.334 0.35 0.366 0.382 0.398 0.414 0.43 0.446 0.462 0.478 0.494 0.51 0.526 0.542 0.558 0.574 0.59 0.606 0.622 0.638 0.654 0.67 0.686 0.702'
y = '-0.041824842 -0.042730269 -0.043412712 -0.04428867 -0.045509181 -0.04645965 -0.047268246 -0.047974749 -0.048597109 -0.0491467 -0.049632388 -0.050061697 -0.050441198 -0.050776675 -0.051073238 -0.0513354 -0.051567152 -0.051772022 -0.051953128 -0.052113227 -0.052254754 -0.052379865 -0.052490464 -0.052588233 -0.052674662 -0.052751065 -0.052818606 -0.052878312 -0.052931093 -0.052977751 -0.053018997 -0.053055459 -0.053087691 -0.053116185 -0.053141373 -0.05316364 -0.053183324 -0.053200724 -0.053216106 -0.053229704 -0.053241725 -0.053252351 -0.053261745 -0.053270049 -0.053277389 -0.053283879 -0.053289615'
[]
[]
[AuxVariables]
[stress_yy]
order = CONSTANT
family = MONOMIAL
[]
[tot_force]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[stress_yy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yy
index_i = 1
index_j = 1
[]
[tot_force]
type = ParsedAux
args = 'stress_yy porepressure'
execute_on = timestep_end
variable = tot_force
function = '-stress_yy+0.6*porepressure'
[]
[]
[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
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
component = 2
variable = disp_z
[]
[mass0]
type = PorousFlowFullySaturatedMassTimeDerivative
biot_coefficient = 0.6
multiply_by_density = false
coupling_type = HydroMechanical
variable = porepressure
[]
[flux]
type = PorousFlowFullySaturatedDarcyBase
multiply_by_density = false
variable = porepressure
gravity = '0 0 0'
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 8
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
[]
[eff_fluid_pressure_qp]
type = PorousFlowEffectiveFluidPressure
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = porepressure
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst # only the initial value of this is ever used
porosity = 0.1
[]
[biot_modulus]
type = PorousFlowConstantBiotModulus
biot_coefficient = 0.6
solid_bulk_compliance = 1
fluid_bulk_modulus = 8
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1.5 0 0 0 1.5 0 0 0 1.5'
[]
[]
[Postprocessors]
[p0]
type = PointValue
outputs = csv
point = '0.0 0 0'
variable = porepressure
[]
[p1]
type = PointValue
outputs = csv
point = '0.1 0 0'
variable = porepressure
[]
[p2]
type = PointValue
outputs = csv
point = '0.2 0 0'
variable = porepressure
[]
[p3]
type = PointValue
outputs = csv
point = '0.3 0 0'
variable = porepressure
[]
[p4]
type = PointValue
outputs = csv
point = '0.4 0 0'
variable = porepressure
[]
[p5]
type = PointValue
outputs = csv
point = '0.5 0 0'
variable = porepressure
[]
[p6]
type = PointValue
outputs = csv
point = '0.6 0 0'
variable = porepressure
[]
[p7]
type = PointValue
outputs = csv
point = '0.7 0 0'
variable = porepressure
[]
[p8]
type = PointValue
outputs = csv
point = '0.8 0 0'
variable = porepressure
[]
[p9]
type = PointValue
outputs = csv
point = '0.9 0 0'
variable = porepressure
[]
[p99]
type = PointValue
outputs = csv
point = '1 0 0'
variable = porepressure
[]
[xdisp]
type = PointValue
outputs = csv
point = '1 0.1 0'
variable = disp_x
[]
[ydisp]
type = PointValue
outputs = csv
point = '1 0.1 0'
variable = disp_y
[]
[total_downwards_force]
type = ElementAverageValue
outputs = csv
variable = tot_force
[]
[dt]
type = FunctionValuePostprocessor
outputs = console
function = if(0.15*t<0.01,0.15*t,0.01)
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'gmres asm lu 1E-14 1E-10 10000'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
start_time = 0
end_time = 0.7
[TimeStepper]
type = PostprocessorDT
postprocessor = dt
dt = 0.001
[]
[]
[Outputs]
execute_on = 'timestep_end'
file_base = mandel_fully_saturated_volume
[csv]
interval = 3
type = CSV
[]
[]
(modules/porous_flow/test/tests/heat_advection/heat_advection_1d_fully_saturated.i)
# 1phase, heat advecting with a moving fluid
# Using the FullySaturated Kernel
[Mesh]
type = GeneratedMesh
dim = 1
nx = 50
xmin = 0
xmax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[temp]
initial_condition = 200
[]
[pp]
[]
[]
[ICs]
[pp]
type = FunctionIC
variable = pp
function = '1-x'
[]
[]
[BCs]
[pp0]
type = DirichletBC
variable = pp
boundary = left
value = 1
[]
[pp1]
type = DirichletBC
variable = pp
boundary = right
value = 0
[]
[spit_heat]
type = DirichletBC
variable = temp
boundary = left
value = 300
[]
[suck_heat]
type = DirichletBC
variable = temp
boundary = right
value = 200
[]
[]
[Kernels]
[mass_dot]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[advection]
type = PorousFlowFullySaturatedDarcyBase
variable = pp
[]
[energy_dot]
type = PorousFlowEnergyTimeDerivative
variable = temp
[]
[convection]
type = PorousFlowFullySaturatedHeatAdvection
variable = temp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'temp pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 100
density0 = 1000
viscosity = 4.4
thermal_expansion = 0
cv = 2
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = temp
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.2
[]
[rock_heat]
type = PorousFlowMatrixInternalEnergy
specific_heat_capacity = 1.0
density = 125
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1.1 0 0 0 2 0 0 0 3'
[]
[massfrac]
type = PorousFlowMassFraction
[]
[PS]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'gmres bjacobi 1E-15 1E-10 10000'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 0.01
end_time = 0.6
[]
[VectorPostprocessors]
[T]
type = LineValueSampler
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 51
sort_by = x
variable = temp
[]
[]
[Outputs]
file_base = heat_advection_1d_fully_saturated
[csv]
type = CSV
sync_times = '0.1 0.6'
sync_only = true
[]
[]
(modules/porous_flow/test/tests/gravity/fully_saturated_grav01b.i)
# Checking that gravity head is established
# 1phase, constant and large fluid-bulk, constant viscosity, constant permeability
# fully saturated with fully-saturated Kernel
# For better agreement with the analytical solution (ana_pp), just increase nx
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmin = -1
xmax = 0
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[InitialCondition]
type = RandomIC
min = 0
max = 1
[]
[]
[]
[Kernels]
[flux0]
type = PorousFlowFullySaturatedDarcyBase
variable = pp
gravity = '-1 0 0'
[]
[]
[Functions]
[ana_pp]
type = ParsedFunction
vars = 'g B p0 rho0'
vals = '1 1E3 0 1'
value = '-B*log(exp(-p0/B)+g*rho0*x/B)' # expected pp at base
[]
[]
[BCs]
[z]
type = DirichletBC
variable = pp
boundary = right
value = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e3
density0 = 1
viscosity = 1
thermal_expansion = 0
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[]
[Postprocessors]
[pp_base]
type = PointValue
variable = pp
point = '-1 0 0'
[]
[pp_analytical]
type = FunctionValuePostprocessor
function = ana_pp
point = '-1 0 0'
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = Newton
[]
[Outputs]
execute_on = 'timestep_end'
file_base = fully_saturated_grav01b
[csv]
type = CSV
[]
[]
(modules/porous_flow/test/tests/gravity/fully_saturated_grav01a.i)
# Checking that gravity head is established
# 1phase, constant fluid-bulk, constant viscosity, constant permeability
# fully saturated
# For better agreement with the analytical solution (ana_pp), just increase nx
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmin = -1
xmax = 0
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[InitialCondition]
type = RandomIC
min = 0
max = 1
[]
[]
[]
[Kernels]
[flux0]
type = PorousFlowFullySaturatedDarcyBase
variable = pp
gravity = '-1 0 0'
[]
[]
[Functions]
[ana_pp]
type = ParsedFunction
vars = 'g B p0 rho0'
vals = '1 1.2 0 1'
value = '-B*log(exp(-p0/B)+g*rho0*x/B)' # expected pp at base
[]
[]
[BCs]
[z]
type = DirichletBC
variable = pp
boundary = right
value = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1.2
density0 = 1
viscosity = 1
thermal_expansion = 0
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[]
[Postprocessors]
[pp_base]
type = PointValue
variable = pp
point = '-1 0 0'
[]
[pp_analytical]
type = FunctionValuePostprocessor
function = ana_pp
point = '-1 0 0'
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = Newton
[]
[Outputs]
execute_on = 'timestep_end'
file_base = fully_saturated_grav01a
[csv]
type = CSV
[]
[]
(modules/porous_flow/test/tests/poro_elasticity/mandel_fully_saturated.i)
# Mandel's problem of consolodation of a drained medium
# Using the FullySaturatedDarcyBase and FullySaturatedMassTimeDerivative kernels
#
# A sample is in plane strain.
# -a <= x <= a
# -b <= y <= b
# It is squashed with constant force by impermeable, frictionless plattens on its top and bottom surfaces (at y=+/-b)
# Fluid is allowed to leak out from its sides (at x=+/-a)
# The porepressure within the sample is monitored.
#
# As is common in the literature, this is simulated by
# considering the quarter-sample, 0<=x<=a and 0<=y<=b, with
# impermeable, roller BCs at x=0 and y=0 and y=b.
# Porepressure is fixed at zero on x=a.
# Porepressure and displacement are initialised to zero.
# Then the top (y=b) is moved downwards with prescribed velocity,
# so that the total force that is inducing this downwards velocity
# is fixed. The velocity is worked out by solving Mandel's problem
# analytically, and the total force is monitored in the simulation
# to check that it indeed remains constant.
#
# Here are the problem's parameters, and their values:
# Soil width. a = 1
# Soil height. b = 0.1
# Soil's Lame lambda. la = 0.5
# Soil's Lame mu, which is also the Soil's shear modulus. mu = G = 0.75
# Soil bulk modulus. K = la + 2*mu/3 = 1
# Drained Poisson ratio. nu = (3K - 2G)/(6K + 2G) = 0.2
# Soil bulk compliance. 1/K = 1
# Fluid bulk modulus. Kf = 8
# Fluid bulk compliance. 1/Kf = 0.125
# Soil initial porosity. phi0 = 0.1
# Biot coefficient. alpha = 0.6
# Biot modulus. M = 1/(phi0/Kf + (alpha - phi0)(1 - alpha)/K) = 4.705882
# Undrained bulk modulus. Ku = K + alpha^2*M = 2.694118
# Undrained Poisson ratio. nuu = (3Ku - 2G)/(6Ku + 2G) = 0.372627
# Skempton coefficient. B = alpha*M/Ku = 1.048035
# Fluid mobility (soil permeability/fluid viscosity). k = 1.5
# Consolidation coefficient. c = 2*k*B^2*G*(1-nu)*(1+nuu)^2/9/(1-nuu)/(nuu-nu) = 3.821656
# Normal stress on top. F = 1
#
# The solution for porepressure and displacements is given in
# AHD Cheng and E Detournay "A direct boundary element method for plane strain poroelasticity" International Journal of Numerical and Analytical Methods in Geomechanics 12 (1988) 551-572.
# The solution involves complicated infinite series, so I shall not write it here
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
ny = 1
nz = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 0.1
zmin = 0
zmax = 1
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
PorousFlowDictator = dictator
block = 0
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'porepressure disp_x disp_y disp_z'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[porepressure]
[]
[]
[BCs]
[roller_xmin]
type = DirichletBC
variable = disp_x
value = 0
boundary = 'left'
[]
[roller_ymin]
type = DirichletBC
variable = disp_y
value = 0
boundary = 'bottom'
[]
[plane_strain]
type = DirichletBC
variable = disp_z
value = 0
boundary = 'back front'
[]
[xmax_drained]
type = DirichletBC
variable = porepressure
value = 0
boundary = right
[]
[top_velocity]
type = FunctionDirichletBC
variable = disp_y
function = top_velocity
boundary = top
[]
[]
[Functions]
[top_velocity]
type = PiecewiseLinear
x = '0 0.002 0.006 0.014 0.03 0.046 0.062 0.078 0.094 0.11 0.126 0.142 0.158 0.174 0.19 0.206 0.222 0.238 0.254 0.27 0.286 0.302 0.318 0.334 0.35 0.366 0.382 0.398 0.414 0.43 0.446 0.462 0.478 0.494 0.51 0.526 0.542 0.558 0.574 0.59 0.606 0.622 0.638 0.654 0.67 0.686 0.702'
y = '-0.041824842 -0.042730269 -0.043412712 -0.04428867 -0.045509181 -0.04645965 -0.047268246 -0.047974749 -0.048597109 -0.0491467 -0.049632388 -0.050061697 -0.050441198 -0.050776675 -0.051073238 -0.0513354 -0.051567152 -0.051772022 -0.051953128 -0.052113227 -0.052254754 -0.052379865 -0.052490464 -0.052588233 -0.052674662 -0.052751065 -0.052818606 -0.052878312 -0.052931093 -0.052977751 -0.053018997 -0.053055459 -0.053087691 -0.053116185 -0.053141373 -0.05316364 -0.053183324 -0.053200724 -0.053216106 -0.053229704 -0.053241725 -0.053252351 -0.053261745 -0.053270049 -0.053277389 -0.053283879 -0.053289615'
[]
[]
[AuxVariables]
[stress_yy]
order = CONSTANT
family = MONOMIAL
[]
[tot_force]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[stress_yy]
type = RankTwoAux
rank_two_tensor = stress
variable = stress_yy
index_i = 1
index_j = 1
[]
[tot_force]
type = ParsedAux
args = 'stress_yy porepressure'
execute_on = timestep_end
variable = tot_force
function = '-stress_yy+0.6*porepressure'
[]
[]
[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
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
component = 2
variable = disp_z
[]
[mass0]
type = PorousFlowFullySaturatedMassTimeDerivative
biot_coefficient = 0.6
coupling_type = HydroMechanical
variable = porepressure
[]
[flux]
type = PorousFlowFullySaturatedDarcyBase
variable = porepressure
gravity = '0 0 0'
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 8
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
[]
[eff_fluid_pressure_qp]
type = PorousFlowEffectiveFluidPressure
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = porepressure
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst # only the initial value of this is ever used
porosity = 0.1
[]
[biot_modulus]
type = PorousFlowConstantBiotModulus
biot_coefficient = 0.6
solid_bulk_compliance = 1
fluid_bulk_modulus = 8
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1.5 0 0 0 1.5 0 0 0 1.5'
[]
[]
[Postprocessors]
[p0]
type = PointValue
outputs = csv
point = '0.0 0 0'
variable = porepressure
[]
[p1]
type = PointValue
outputs = csv
point = '0.1 0 0'
variable = porepressure
[]
[p2]
type = PointValue
outputs = csv
point = '0.2 0 0'
variable = porepressure
[]
[p3]
type = PointValue
outputs = csv
point = '0.3 0 0'
variable = porepressure
[]
[p4]
type = PointValue
outputs = csv
point = '0.4 0 0'
variable = porepressure
[]
[p5]
type = PointValue
outputs = csv
point = '0.5 0 0'
variable = porepressure
[]
[p6]
type = PointValue
outputs = csv
point = '0.6 0 0'
variable = porepressure
[]
[p7]
type = PointValue
outputs = csv
point = '0.7 0 0'
variable = porepressure
[]
[p8]
type = PointValue
outputs = csv
point = '0.8 0 0'
variable = porepressure
[]
[p9]
type = PointValue
outputs = csv
point = '0.9 0 0'
variable = porepressure
[]
[p99]
type = PointValue
outputs = csv
point = '1 0 0'
variable = porepressure
[]
[xdisp]
type = PointValue
outputs = csv
point = '1 0.1 0'
variable = disp_x
[]
[ydisp]
type = PointValue
outputs = csv
point = '1 0.1 0'
variable = disp_y
[]
[total_downwards_force]
type = ElementAverageValue
outputs = csv
variable = tot_force
[]
[dt]
type = FunctionValuePostprocessor
outputs = console
function = if(0.15*t<0.01,0.15*t,0.01)
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'gmres asm lu 1E-14 1E-10 10000'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
start_time = 0
end_time = 0.7
[TimeStepper]
type = PostprocessorDT
postprocessor = dt
dt = 0.001
[]
[]
[Outputs]
execute_on = 'timestep_end'
file_base = mandel_fully_saturated
[csv]
interval = 3
type = CSV
[]
[]
(modules/porous_flow/test/tests/poro_elasticity/terzaghi_fully_saturated_volume.i)
# Terzaghi's problem of consolodation of a drained medium
# The FullySaturated Kernels are used, with multiply_by_density = false
# so that this becomes a linear problem with constant Biot Modulus
#
# A saturated soil sample sits in a bath of water.
# It is constrained on its sides, and bottom.
# Its sides and bottom are also impermeable.
# Initially it is unstressed.
# A normal stress, q, is applied to the soil's top.
# The soil then slowly compresses as water is squeezed
# out from the sample from its top (the top BC for
# the porepressure is porepressure = 0).
#
# See, for example. Section 2.2 of the online manuscript
# Arnold Verruijt "Theory and Problems of Poroelasticity" Delft University of Technology 2013
# but note that the "sigma" in that paper is the negative
# of the stress in TensorMechanics
#
# Here are the problem's parameters, and their values:
# Soil height. h = 10
# Soil's Lame lambda. la = 2
# Soil's Lame mu, which is also the Soil's shear modulus. mu = 3
# Soil bulk modulus. K = la + 2*mu/3 = 4
# Soil confined compressibility. m = 1/(K + 4mu/3) = 0.125
# Soil bulk compliance. 1/K = 0.25
# Fluid bulk modulus. Kf = 8
# Fluid bulk compliance. 1/Kf = 0.125
# Fluid mobility (soil permeability/fluid viscosity). k = 1.5
# Soil initial porosity. phi0 = 0.1
# Biot coefficient. alpha = 0.6
# Soil initial storativity, which is the reciprocal of the initial Biot modulus. S = phi0/Kf + (alpha - phi0)(1 - alpha)/K = 0.0625
# Consolidation coefficient. c = k/(S + alpha^2 m) = 13.95348837
# Normal stress on top. q = 1
# Initial porepressure, resulting from instantaneous application of q, assuming corresponding instantaneous increase of porepressure (Note that this is calculated by MOOSE: we only need it for the analytical solution). p0 = alpha*m*q/(S + alpha^2 m) = 0.69767442
# Initial vertical displacement (down is positive), resulting from instantaneous application of q (Note this is calculated by MOOSE: we only need it for the analytical solution). uz0 = q*m*h*S/(S + alpha^2 m)
# Final vertical displacement (down in positive) (Note this is calculated by MOOSE: we only need it for the analytical solution). uzinf = q*m*h
#
# The solution for porepressure is
# P = 4*p0/\pi \sum_{k=1}^{\infty} \frac{(-1)^{k-1}}{2k-1} \cos ((2k-1)\pi z/(2h)) \exp(-(2k-1)^2 \pi^2 ct/(4 h^2))
# This series converges very slowly for ct/h^2 small, so in that domain
# P = p0 erf( (1-(z/h))/(2 \sqrt(ct/h^2)) )
#
# The degree of consolidation is defined as
# U = (uz - uz0)/(uzinf - uz0)
# where uz0 and uzinf are defined above, and
# uz = the vertical displacement of the top (down is positive)
# U = 1 - (8/\pi^2)\sum_{k=1}^{\infty} \frac{1}{(2k-1)^2} \exp(-(2k-1)^2 \pi^2 ct/(4 h^2))
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 10
xmin = -1
xmax = 1
ymin = -1
ymax = 1
zmin = 0
zmax = 10
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
PorousFlowDictator = dictator
block = 0
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'porepressure disp_x disp_y disp_z'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[porepressure]
[]
[]
[BCs]
[confinex]
type = DirichletBC
variable = disp_x
value = 0
boundary = 'left right'
[]
[confiney]
type = DirichletBC
variable = disp_y
value = 0
boundary = 'bottom top'
[]
[basefixed]
type = DirichletBC
variable = disp_z
value = 0
boundary = back
[]
[topdrained]
type = DirichletBC
variable = porepressure
value = 0
boundary = front
[]
[topload]
type = NeumannBC
variable = disp_z
value = -1
boundary = front
[]
[]
[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
[]
[poro_x]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
variable = disp_x
component = 0
[]
[poro_y]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
variable = disp_y
component = 1
[]
[poro_z]
type = PorousFlowEffectiveStressCoupling
biot_coefficient = 0.6
component = 2
variable = disp_z
[]
[mass0]
type = PorousFlowFullySaturatedMassTimeDerivative
coupling_type = HydroMechanical
biot_coefficient = 0.6
multiply_by_density = false
variable = porepressure
[]
[flux]
type = PorousFlowFullySaturatedDarcyBase
multiply_by_density = false
variable = porepressure
gravity = '0 0 0'
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 8
density0 = 1
thermal_expansion = 0
viscosity = 0.96
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[elasticity_tensor]
type = ComputeElasticityTensor
C_ijkl = '2 3'
# bulk modulus is lambda + 2*mu/3 = 2 + 2*3/3 = 4
fill_method = symmetric_isotropic
[]
[strain]
type = ComputeSmallStrain
[]
[stress]
type = ComputeLinearElasticStress
[]
[eff_fluid_pressure_qp]
type = PorousFlowEffectiveFluidPressure
[]
[vol_strain]
type = PorousFlowVolumetricStrain
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = porepressure
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst # only the initial value of this is used
porosity = 0.1
[]
[biot_modulus]
type = PorousFlowConstantBiotModulus
biot_coefficient = 0.6
fluid_bulk_modulus = 8
solid_bulk_compliance = 0.25
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1.5 0 0 0 1.5 0 0 0 1.5'
[]
[]
[Postprocessors]
[p0]
type = PointValue
outputs = csv
point = '0 0 0'
variable = porepressure
use_displaced_mesh = false
[]
[p1]
type = PointValue
outputs = csv
point = '0 0 1'
variable = porepressure
use_displaced_mesh = false
[]
[p2]
type = PointValue
outputs = csv
point = '0 0 2'
variable = porepressure
use_displaced_mesh = false
[]
[p3]
type = PointValue
outputs = csv
point = '0 0 3'
variable = porepressure
use_displaced_mesh = false
[]
[p4]
type = PointValue
outputs = csv
point = '0 0 4'
variable = porepressure
use_displaced_mesh = false
[]
[p5]
type = PointValue
outputs = csv
point = '0 0 5'
variable = porepressure
use_displaced_mesh = false
[]
[p6]
type = PointValue
outputs = csv
point = '0 0 6'
variable = porepressure
use_displaced_mesh = false
[]
[p7]
type = PointValue
outputs = csv
point = '0 0 7'
variable = porepressure
use_displaced_mesh = false
[]
[p8]
type = PointValue
outputs = csv
point = '0 0 8'
variable = porepressure
use_displaced_mesh = false
[]
[p9]
type = PointValue
outputs = csv
point = '0 0 9'
variable = porepressure
use_displaced_mesh = false
[]
[p99]
type = PointValue
outputs = csv
point = '0 0 10'
variable = porepressure
use_displaced_mesh = false
[]
[zdisp]
type = PointValue
outputs = csv
point = '0 0 10'
variable = disp_z
use_displaced_mesh = false
[]
[dt]
type = FunctionValuePostprocessor
outputs = console
function = if(0.5*t<0.1,0.5*t,0.1)
[]
[]
[Preconditioning]
[andy]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
start_time = 0
end_time = 10
[TimeStepper]
type = PostprocessorDT
postprocessor = dt
dt = 0.0001
[]
[]
[Outputs]
execute_on = 'timestep_end'
file_base = terzaghi_fully_saturated_volume
[csv]
type = CSV
[]
[]
(modules/porous_flow/test/tests/pressure_pulse/pressure_pulse_1d_fully_saturated_2.i)
# Pressure pulse in 1D with 1 phase - transient
# using the PorousFlowFullySaturatedDarcyBase Kernel
# and the PorousFlowFullySaturatedMassTimeDerivative Kernel
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
xmin = 0
xmax = 100
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
initial_condition = 2E6
[]
[]
[Kernels]
[mass0]
type = PorousFlowFullySaturatedMassTimeDerivative
variable = pp
[]
[flux]
type = PorousFlowFullySaturatedDarcyBase
variable = pp
gravity = '0 0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature_qp]
type = PorousFlowTemperature
[]
[ppss_qp]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid_qp]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[biot_modulus]
type = PorousFlowConstantBiotModulus
fluid_bulk_modulus = 2E9
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-15 0 0 0 1E-15 0 0 0 1E-15'
[]
[]
[BCs]
[left]
type = DirichletBC
boundary = left
value = 3E6
variable = pp
[]
[]
[Preconditioning]
[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-20 10000'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1E3
end_time = 1E4
[]
[Postprocessors]
[p005]
type = PointValue
variable = pp
point = '5 0 0'
execute_on = 'initial timestep_end'
[]
[p015]
type = PointValue
variable = pp
point = '15 0 0'
execute_on = 'initial timestep_end'
[]
[p025]
type = PointValue
variable = pp
point = '25 0 0'
execute_on = 'initial timestep_end'
[]
[p035]
type = PointValue
variable = pp
point = '35 0 0'
execute_on = 'initial timestep_end'
[]
[p045]
type = PointValue
variable = pp
point = '45 0 0'
execute_on = 'initial timestep_end'
[]
[p055]
type = PointValue
variable = pp
point = '55 0 0'
execute_on = 'initial timestep_end'
[]
[p065]
type = PointValue
variable = pp
point = '65 0 0'
execute_on = 'initial timestep_end'
[]
[p075]
type = PointValue
variable = pp
point = '75 0 0'
execute_on = 'initial timestep_end'
[]
[p085]
type = PointValue
variable = pp
point = '85 0 0'
execute_on = 'initial timestep_end'
[]
[p095]
type = PointValue
variable = pp
point = '95 0 0'
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
file_base = pressure_pulse_1d_fully_saturated_2
print_linear_residuals = false
csv = true
[]
(modules/porous_flow/test/tests/pressure_pulse/pressure_pulse_1d_fully_saturated.i)
# Pressure pulse in 1D with 1 phase - transient
# using the PorousFlowFullySaturatedDarcyBase Kernel
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
xmin = 0
xmax = 100
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
initial_condition = 2E6
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[flux]
type = PorousFlowFullySaturatedDarcyBase
variable = pp
gravity = '0 0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-15 0 0 0 1E-15 0 0 0 1E-15'
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 0
phase = 0
[]
[]
[BCs]
[left]
type = DirichletBC
boundary = left
value = 3E6
variable = pp
[]
[]
[Preconditioning]
[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-20 10000'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1E3
end_time = 1E4
[]
[Postprocessors]
[p005]
type = PointValue
variable = pp
point = '5 0 0'
execute_on = 'initial timestep_end'
[]
[p015]
type = PointValue
variable = pp
point = '15 0 0'
execute_on = 'initial timestep_end'
[]
[p025]
type = PointValue
variable = pp
point = '25 0 0'
execute_on = 'initial timestep_end'
[]
[p035]
type = PointValue
variable = pp
point = '35 0 0'
execute_on = 'initial timestep_end'
[]
[p045]
type = PointValue
variable = pp
point = '45 0 0'
execute_on = 'initial timestep_end'
[]
[p055]
type = PointValue
variable = pp
point = '55 0 0'
execute_on = 'initial timestep_end'
[]
[p065]
type = PointValue
variable = pp
point = '65 0 0'
execute_on = 'initial timestep_end'
[]
[p075]
type = PointValue
variable = pp
point = '75 0 0'
execute_on = 'initial timestep_end'
[]
[p085]
type = PointValue
variable = pp
point = '85 0 0'
execute_on = 'initial timestep_end'
[]
[p095]
type = PointValue
variable = pp
point = '95 0 0'
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
file_base = pressure_pulse_1d_fully_saturated
print_linear_residuals = false
csv = true
[]
(modules/porous_flow/include/kernels/PorousFlowFullySaturatedHeatAdvection.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "PorousFlowFullySaturatedDarcyBase.h"
/**
* Advection of heat via flux via Darcy flow of a single phase
* fully-saturated fluid. No upwinding is used.
*/
class PorousFlowFullySaturatedHeatAdvection : public PorousFlowFullySaturatedDarcyBase
{
public:
static InputParameters validParams();
PorousFlowFullySaturatedHeatAdvection(const InputParameters & parameters);
protected:
virtual Real mobility() const override;
virtual Real dmobility(unsigned pvar) const override;
/// Enthalpy of each phase
const MaterialProperty<std::vector<Real>> & _enthalpy;
/// Derivative of the enthalpy wrt PorousFlow variables
const MaterialProperty<std::vector<std::vector<Real>>> & _denthalpy_dvar;
};
(modules/porous_flow/include/kernels/PorousFlowFullySaturatedDarcyFlow.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "PorousFlowFullySaturatedDarcyBase.h"
/**
* Darcy advective flux for a fully-saturated,
* single-phase, multi-component fluid.
* No upwinding or relative-permeability is used.
*/
class PorousFlowFullySaturatedDarcyFlow : public PorousFlowFullySaturatedDarcyBase
{
public:
static InputParameters validParams();
PorousFlowFullySaturatedDarcyFlow(const InputParameters & parameters);
protected:
/**
* The mobility of the fluid = mass_fraction * density / viscosity
*/
virtual Real mobility() const override;
/**
* The derivative of the mobility with respect to the PorousFlow variable pvar
* @param pvar Take the derivative with respect to this PorousFlow variable
*/
virtual Real dmobility(unsigned pvar) const override;
/// mass fraction of the components in the phase
const MaterialProperty<std::vector<std::vector<Real>>> & _mfrac;
/// Derivative of mass fraction wrt wrt PorousFlow variables
const MaterialProperty<std::vector<std::vector<std::vector<Real>>>> & _dmfrac_dvar;
/// The fluid component for this Kernel
const unsigned int _fluid_component;
};