Point and line sources and sinks
This page may be read in conjunction with the theory behind Dirac Kernels in PorousFlow.
Geometric tests
The test suite contains many core tests that demonstrate:
when a point sink is placed at a node, it withdraws fluid (or heat) from only that node, for instance
# fully-saturated situation with a poly-line sink at one
# of the nodes. Because there is no fluid flow, the
# other nodes should not experience any change in
# porepressure.
# The poly-line sink has length=2 and weight=0.1, and
# extracts fluid at a constant rate of 1 kg.m^-1.s^-1.
# Therefore, in 1 second it will have extracted a total
# of 0.2 kg.
# The porosity is 0.1, and the elemental volume is 2,
# so the fluid mass at the node in question = 0.2 * density / 4,
# where the 4 is the number of nodes in the element.
# In this simulation density = dens0 * exp(P / bulk), with
# dens0 = 100, and bulk = 20 MPa.
# The initial porepressure P0 = 10 MPa, so the final (after
# 1 second of simulation) is
# P(t=1) = 0.950879 MPa
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 1
xmin = 0
xmax = 2
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
initial_condition = 1E7
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pls_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1e-7
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e7
density0 = 100
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[]
[DiracKernels]
[pls]
type = PorousFlowPolyLineSink
fluid_phase = 0
point_file = pls01_21.bh
line_length = 2
SumQuantityUO = pls_total_outflow_mass
variable = pp
p_or_t_vals = '0 1E7'
fluxes = '1 1'
[]
[]
[Postprocessors]
[pls_report]
type = PorousFlowPlotQuantity
uo = pls_total_outflow_mass
[]
[fluid_mass0]
type = PorousFlowFluidMass
execute_on = timestep_begin
[]
[fluid_mass1]
type = PorousFlowFluidMass
execute_on = timestep_end
[]
[zmass_error]
type = FunctionValuePostprocessor
function = mass_bal_fcn
execute_on = timestep_end
indirect_dependencies = 'fluid_mass1 fluid_mass0 pls_report'
[]
[p00]
type = PointValue
variable = pp
point = '0 0 0'
execute_on = timestep_end
[]
[p01]
type = PointValue
variable = pp
point = '0 1 0'
execute_on = timestep_end
[]
[p20]
type = PointValue
variable = pp
point = '2 0 0'
execute_on = timestep_end
[]
[p21]
type = PointValue
variable = pp
point = '2 1 0'
execute_on = timestep_end
[]
[]
[Functions]
[mass_bal_fcn]
type = ParsedFunction
expression = abs((a-c+d)/2/(a+c))
symbol_names = 'a c d'
symbol_values = 'fluid_mass1 fluid_mass0 pls_report'
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 1
dt = 1
solve_type = NEWTON
[]
[Outputs]
file_base = pls01
exodus = false
csv = true
execute_on = timestep_end
[]
(modules/porous_flow/test/tests/dirackernels/pls01.i)when a point sink that is proportional to mobility (or relative permeability, etc) is placed in an element where some nodes have zero mobility (or relative permeability, etc), then fluid (or heat) is not extracted from those nodes, for instance
# Test that the upwinding works correctly.
#
# A poly-line sink sits at the centre of the element.
# It has length=4 and weight=0.5, and extracts fluid
# at a constant rate of
# (1 * relative_permeability) kg.m^-1.s^-1
# Since it sits at the centre of the element, it extracts
# equally from each node, so the rate of extraction from
# each node is
# (0.5 * relative_permeability) kg.s^-1
# including the length and weight effects.
#
# There is no fluid flow.
#
# The initial conditions are such that all nodes have
# relative_permeability=0, except for one which has
# relative_permeaility = 1. Therefore, all nodes should
# remain at their initial porepressure, except the one.
#
# The porosity is 0.1, and the elemental volume is 2,
# so the fluid mass at the node in question = 0.2 * density / 4,
# where the 4 is the number of nodes in the element.
# In this simulation density = dens0 * exp(P / bulk), with
# dens0 = 100, and bulk = 20 MPa.
# The initial porepressure P0 = 10 MPa, so the final (after
# 1 second of simulation) is
# P(t=1) = 8.748592 MPa
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 1
xmin = 0
xmax = 2
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[]
[ICs]
[pp]
type = FunctionIC
variable = pp
#function = if((x<1)&(y<0.5),1E7,-1E7)
function = if((x<1)&(y>0.5),1E7,-1E7)
#function = if((x>1)&(y<0.5),1E7,-1E7)
#function = if((x>1)&(y>0.5),1E7,-1E7)
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pls_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1e-7
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e7
density0 = 100
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[relperm]
type = PorousFlowRelativePermeabilityFLAC
phase = 0
m = 2
s_res = 0.99
sum_s_res = 0.99
[]
[]
[DiracKernels]
[pls]
type = PorousFlowPolyLineSink
fluid_phase = 0
point_file = pls03.bh
use_relative_permeability = true
line_length = 4
SumQuantityUO = pls_total_outflow_mass
variable = pp
p_or_t_vals = '0 1E7'
fluxes = '1 1'
[]
[]
[Postprocessors]
[pls_report]
type = PorousFlowPlotQuantity
uo = pls_total_outflow_mass
[]
[fluid_mass0]
type = PorousFlowFluidMass
execute_on = timestep_begin
[]
[fluid_mass1]
type = PorousFlowFluidMass
execute_on = timestep_end
[]
[zmass_error]
type = FunctionValuePostprocessor
function = mass_bal_fcn
execute_on = timestep_end
indirect_dependencies = 'fluid_mass1 fluid_mass0 pls_report'
[]
[p00]
type = PointValue
variable = pp
point = '0 0 0'
execute_on = timestep_end
[]
[p01]
type = PointValue
variable = pp
point = '0 1 0'
execute_on = timestep_end
[]
[p20]
type = PointValue
variable = pp
point = '2 0 0'
execute_on = timestep_end
[]
[p21]
type = PointValue
variable = pp
point = '2 1 0'
execute_on = timestep_end
[]
[]
[Functions]
[mass_bal_fcn]
type = ParsedFunction
expression = abs((a-c+d)/2/(a+c))
symbol_names = 'a c d'
symbol_values = 'fluid_mass1 fluid_mass0 pls_report'
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 1
dt = 1
solve_type = NEWTON
[]
[Outputs]
file_base = pls03
exodus = false
csv = true
execute_on = timestep_end
[]
(modules/porous_flow/test/tests/dirackernels/pls03.i)Basic point sources/sinks
The following input file extracts fluid at the rate (kg.s)
# Test PorousFlowSquarePulsePointSource DiracKernel
[Mesh]
type = GeneratedMesh
dim = 2
bias_x = 1.1
bias_y = 1.1
ymax = 1
xmax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = pp
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1e-7
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.2
[]
[]
[Postprocessors]
[total_mass]
type = PorousFlowFluidMass
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
nl_abs_tol = 1e-14
dt = 200
end_time = 2000
[]
[Outputs]
perf_graph = true
file_base = squarepulse1
csv = true
execute_on = 'initial timestep_end'
[con]
output_linear = true
type = Console
[]
[]
[ICs]
[PressureIC]
variable = pp
type = ConstantIC
value = 20e6
[]
[]
[DiracKernels]
[sink1]
type = PorousFlowSquarePulsePointSource
start_time = 100
end_time = 300
point = '0.5 0.5 0'
mass_flux = -0.1
variable = pp
[]
[sink]
type = PorousFlowSquarePulsePointSource
start_time = 600
end_time = 1400
point = '0.5 0.5 0'
mass_flux = -0.1
variable = pp
[]
[source]
point = '0.5 0.5 0'
start_time = 1500
mass_flux = 0.2
end_time = 2000
variable = pp
type = PorousFlowSquarePulsePointSource
[]
[]
(modules/porous_flow/test/tests/dirackernels/squarepulse1.i)MOOSE produces the expected result:
Theis tests
In the fully-saturated, isothermal case with no mechanical coupling and constant, large fluid bulk modulus, the fluid flow equations reduce to the form conventionally used in groundwater flow: where
is the "head", . Here is the fluid density.
is the "specific storage". If the rock compressibility is ignored then , where is the fluid bulk modulus
is the hydraulic conductivity: , where is the permeability and is the fluid viscosity
This equation is identical to the physics described by the PorousFlowFullySaturatedMassTimeDerivative and the PorousFlowFullySaturatedDarcyBase Kernels, which may be used explicitly, or through the PorousFlowBasicTHM Action. When using these, the specific storage in the above equation is replaced by (reciprocal of the Biot modulus), and by , and by . Finally, remember that this formulation is volume based, not mass based like the remainder of PorousFlow.
Place a constant volumetric sink of strength (m.m.s) acting along an infinite line in an isotropic 3D medium. The situation is therefore two dimensional. Theis provided the solution for the head (1) which is frequently used in the groundwater literature. "Ei" is the exponential integral function, with values for small (where is the Euler number ), and for large .
This is checked using a number of PorousFlow tests (the extraction specified using mass rate, with volume rate, and from single and 2-phase systems). For instance:
# Theis problem: Flow to single sink using BasicTHM
# SinglePhase
# RZ mesh
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
xmax = 100
bias_x = 1.05
coord_type = RZ
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
initial_condition = 20E6
[]
[]
[PorousFlowBasicTHM]
dictator_name = dictator
add_darcy_aux = false
fp = simple_fluid
gravity = '0 0 0'
multiply_by_density = false
porepressure = pp
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
viscosity = 0.001
[]
[]
[Materials]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.05
[]
[biot_mod]
type = PorousFlowConstantBiotModulus
fluid_bulk_modulus = 2E9
biot_coefficient = 1.0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-14 0 0 0 1E-14 0 0 0 1E-14'
[]
[]
[DiracKernels]
[sink]
type = PorousFlowSquarePulsePointSource
point = '0 0 0'
mass_flux = -0.16E-3 # recall this is a volumetric flux because multiply_by_density = false in the Action, so this corresponds to a mass_flux of 0.16 kg/s/m because density=1000
variable = pp
[]
[]
[VectorPostprocessors]
[pp]
type = LineValueSampler
num_points = 25
start_point = '0 0 0'
end_point = '100 0 0'
sort_by = x
variable = pp
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 200
end_time = 1E3
nl_abs_tol = 1e-10
[]
[Outputs]
perf_graph = true
[csv]
type = CSV
execute_on = final
[]
[]
(modules/porous_flow/test/tests/dirackernels/theis_rz.i)MOOSE agrees with the Theis solution:
Flow from a source in a radial 1D problem admits a similarity solution, where the fluid pressure (and saturation, mass fraction, etc) is a function of the similarity variable . This may be observed in Eq. (1). Using this fact, two-phase immiscible flow may be tested. The following input file describes injection of a gas phase into a fully liquid-saturated model at a constant rate:
# Two phase Theis problem: Flow from single source
# Constant rate injection 0.5 kg/s
# 1D cylindrical mesh
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmax = 2000
bias_x = 1.05
coord_type = RZ
rz_coord_axis = Y
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[ppwater]
initial_condition = 20e6
[]
[sgas]
initial_condition = 0
[]
[]
[AuxVariables]
[massfrac_ph0_sp0]
initial_condition = 1
[]
[massfrac_ph1_sp0]
initial_condition = 0
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = ppwater
[]
[flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = ppwater
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = sgas
[]
[flux1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = sgas
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'ppwater sgas'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureConst
pc = 1e5
[]
[]
[FluidProperties]
[simple_fluid0]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
viscosity = 1e-3
thermal_expansion = 0
[]
[simple_fluid1]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 10
viscosity = 1e-4
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow2PhasePS
phase0_porepressure = ppwater
phase1_saturation = sgas
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac_ph0_sp0 massfrac_ph1_sp0'
[]
[simple_fluid0]
type = PorousFlowSingleComponentFluid
fp = simple_fluid0
phase = 0
compute_enthalpy = false
compute_internal_energy = false
[]
[simple_fluid1]
type = PorousFlowSingleComponentFluid
fp = simple_fluid1
phase = 1
compute_enthalpy = false
compute_internal_energy = false
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.2
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[]
[relperm_water]
type = PorousFlowRelativePermeabilityCorey
n = 1
phase = 0
[]
[relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 1
phase = 1
[]
[]
[BCs]
[rightwater]
type = DirichletBC
boundary = right
value = 20e6
variable = ppwater
[]
[]
[DiracKernels]
[source]
type = PorousFlowSquarePulsePointSource
point = '0 0 0'
mass_flux = 0.5
variable = sgas
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options = '-snes_converged_reason -ksp_diagonal_scale -ksp_diagonal_scale_fix -ksp_gmres_modifiedgramschmidt -snes_linesearch_monitor'
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap -snes_atol -snes_rtol -snes_max_it'
petsc_options_value = 'gmres asm lu NONZERO 2 1E-8 1E-10 20'
[]
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 1e4
[TimeStepper]
type = IterationAdaptiveDT
dt = 10
growth_factor = 2
[]
[]
[VectorPostprocessors]
[line]
type = NodalValueSampler
sort_by = x
variable = 'ppwater sgas'
execute_on = 'timestep_end'
[]
[]
[Postprocessors]
[ppwater]
type = PointValue
point = '4 0 0'
variable = ppwater
[]
[sgas]
type = PointValue
point = '4 0 0'
variable = sgas
[]
[massgas]
type = PorousFlowFluidMass
fluid_component = 1
[]
[]
[Outputs]
file_base = theis3
print_linear_residuals = false
perf_graph = true
[csv]
type = CSV
execute_on = timestep_end
execute_vector_postprocessors_on = final
[]
[]
(modules/porous_flow/test/tests/dirackernels/theis3.i)Figure 3 shows the comparison of similarity solutions calculated with either fixed radial distance or fixed time. In this case, good agreement is observed between the two results for both liquid pressure and gas saturation.
Further tests of this kind, involving multi-component, multi-phase flow, including mutual dissolution of gas and liquid phases may be found here.
Peaceman borehole fluxes
The test suite that checks that the Peaceman flux (2) is correctly implemented. A vertical borehole is placed through the centre of a single element, and fluid flow to the borehole as a function of porepressure is measured. The tests are
A production borehole with , with a fully-saturated medium.
# fully-saturated
# production
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = -1
xmax = 1
ymin = -1
ymax = 1
zmin = -1
zmax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
initial_condition = 1E7
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[borehole_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1e-7
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-12 0 0 0 1E-12 0 0 0 1E-12'
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[]
[]
[DiracKernels]
[bh]
type = PorousFlowPeacemanBorehole
# Because the Variable for this Sink is pp, and pp is associated
# with the fluid-mass conservation equation, this sink is extracting
# fluid mass (and not heat energy or something else)
variable = pp
# The following specfies that the total fluid mass coming out of
# the porespace via this sink in this timestep should be recorded
# in the pls_total_outflow_mass UserObject
SumQuantityUO = borehole_total_outflow_mass
# The following file defines the polyline geometry
# which is just two points in this particular example
point_file = bh02.bh
# First, we want Peacemans f to be a function of porepressure (and not
# temperature or something else). So bottom_p_or_t is actually porepressure
function_of = pressure
fluid_phase = 0
# The bottomhole pressure
bottom_p_or_t = 0
# In this example there is no increase of the wellbore pressure
# due to gravity:
unit_weight = '0 0 0'
# PeacemanBoreholes should almost always have use_mobility = true
use_mobility = true
# This is a production wellbore (a sink of fluid that removes fluid from porespace)
character = 1
[]
[]
[Postprocessors]
[bh_report]
type = PorousFlowPlotQuantity
uo = borehole_total_outflow_mass
[]
[fluid_mass0]
type = PorousFlowFluidMass
execute_on = timestep_begin
[]
[fluid_mass1]
type = PorousFlowFluidMass
execute_on = timestep_end
[]
[zmass_error]
type = FunctionValuePostprocessor
function = mass_bal_fcn
execute_on = timestep_end
indirect_dependencies = 'fluid_mass1 fluid_mass0 bh_report'
[]
[p0]
type = PointValue
variable = pp
point = '0 0 0'
execute_on = timestep_end
[]
[]
[Functions]
[mass_bal_fcn]
type = ParsedFunction
expression = abs((a-c+d)/2/(a+c))
symbol_names = 'a c d'
symbol_values = 'fluid_mass1 fluid_mass0 bh_report'
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 0.5
dt = 1E-2
solve_type = NEWTON
[]
[Outputs]
file_base = bh02
exodus = false
csv = true
execute_on = timestep_end
[]
(modules/porous_flow/test/tests/dirackernels/bh02.i)An injection borehole with MPa, with a fully-saturated medium.
# fully-saturated
# injection
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = 1
xmax = 3
ymin = -1
ymax = 1
zmin = -1
zmax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
initial_condition = 0
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[borehole_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1e-7
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-12 0 0 0 1E-12 0 0 0 1E-12'
[]
[relperm]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[]
[]
[DiracKernels]
[bh]
type = PorousFlowPeacemanBorehole
variable = pp
SumQuantityUO = borehole_total_outflow_mass
point_file = bh03.bh
function_of = pressure
fluid_phase = 0
bottom_p_or_t = 'insitu_pp'
unit_weight = '0 0 0'
use_mobility = true
character = -1
[]
[]
[Postprocessors]
[bh_report]
type = PorousFlowPlotQuantity
uo = borehole_total_outflow_mass
[]
[fluid_mass0]
type = PorousFlowFluidMass
execute_on = timestep_begin
[]
[fluid_mass1]
type = PorousFlowFluidMass
execute_on = timestep_end
[]
[zmass_error]
type = FunctionValuePostprocessor
function = mass_bal_fcn
execute_on = timestep_end
indirect_dependencies = 'fluid_mass1 fluid_mass0 bh_report'
[]
[p0]
type = PointValue
variable = pp
point = '2 0 0'
execute_on = timestep_end
[]
[]
[Functions]
[mass_bal_fcn]
type = ParsedFunction
expression = abs((a-c+d)/2/(a+c))
symbol_names = 'a c d'
symbol_values = 'fluid_mass1 fluid_mass0 bh_report'
[]
[insitu_pp]
type = ParsedFunction
expression = '0.5e7*x' #bh is located at x=2
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 0.5
dt = 1E-2
solve_type = NEWTON
[]
[Outputs]
file_base = bh03
exodus = false
csv = true
execute_on = timestep_end
[]
(modules/porous_flow/test/tests/dirackernels/bh03.i)A production borehole with MPa with a fully-saturated medium, with
use_mobility=true
# fully-saturated
# production
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = -1
xmax = 1
ymin = -1
ymax = 1
zmin = -1
zmax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Functions]
[dts]
type = PiecewiseLinear
y = '1E-2 1E-1 1 1E1 1E2 1E3'
x = '0 1E-1 1 1E1 1E2 1E3'
[]
[]
[Variables]
[pp]
initial_condition = 0
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[borehole_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.8
alpha = 1e-5
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-12 0 0 0 1E-12 0 0 0 1E-12'
[]
[relperm]
type = PorousFlowRelativePermeabilityFLAC
m = 2
phase = 0
[]
[]
[DiracKernels]
[bh]
type = PorousFlowPeacemanBorehole
variable = pp
SumQuantityUO = borehole_total_outflow_mass
point_file = bh02.bh
fluid_phase = 0
bottom_p_or_t = -1E6
unit_weight = '0 0 0'
use_mobility = true
character = 1
[]
[]
[Postprocessors]
[bh_report]
type = PorousFlowPlotQuantity
uo = borehole_total_outflow_mass
[]
[fluid_mass0]
type = PorousFlowFluidMass
execute_on = timestep_begin
[]
[fluid_mass1]
type = PorousFlowFluidMass
execute_on = timestep_end
[]
[zmass_error]
type = FunctionValuePostprocessor
function = mass_bal_fcn
execute_on = timestep_end
indirect_dependencies = 'fluid_mass1 fluid_mass0 bh_report'
[]
[p0]
type = PointValue
variable = pp
point = '0 0 0'
execute_on = timestep_end
[]
[]
[Functions]
[mass_bal_fcn]
type = ParsedFunction
expression = abs((a-c+d)/2/(a+c))
symbol_names = 'a c d'
symbol_values = 'fluid_mass1 fluid_mass0 bh_report'
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 1E3
solve_type = NEWTON
[TimeStepper]
type = FunctionDT
function = dts
[]
[]
[Outputs]
file_base = bh04
exodus = false
csv = true
execute_on = timestep_end
[]
(modules/porous_flow/test/tests/dirackernels/bh04.i)An injection borehole with with an unsaturated medium, with
use_mobility=true
# unsaturated
# injection
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = 1
xmax = 3
ymin = -1
ymax = 1
zmin = -1
zmax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Functions]
[dts]
type = PiecewiseLinear
y = '500 500 1E1'
x = '4000 5000 6500'
[]
[]
[Variables]
[pp]
initial_condition = -2E5
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[borehole_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.8
alpha = 1e-5
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-12 0 0 0 1E-12 0 0 0 1E-12'
[]
[relperm]
type = PorousFlowRelativePermeabilityFLAC
m = 2
phase = 0
[]
[]
[DiracKernels]
[bh]
type = PorousFlowPeacemanBorehole
variable = pp
SumQuantityUO = borehole_total_outflow_mass
point_file = bh03.bh
fluid_phase = 0
bottom_p_or_t = 0
unit_weight = '0 0 0'
use_mobility = true
character = -1
[]
[]
[Postprocessors]
[bh_report]
type = PorousFlowPlotQuantity
uo = borehole_total_outflow_mass
[]
[fluid_mass0]
type = PorousFlowFluidMass
execute_on = timestep_begin
[]
[fluid_mass1]
type = PorousFlowFluidMass
execute_on = timestep_end
[]
[zmass_error]
type = FunctionValuePostprocessor
function = mass_bal_fcn
execute_on = timestep_end
indirect_dependencies = 'fluid_mass1 fluid_mass0 bh_report'
[]
[p0]
type = PointValue
variable = pp
point = '2 0 0'
execute_on = timestep_end
[]
[]
[Functions]
[mass_bal_fcn]
type = ParsedFunction
expression = abs((a-c+d)/2/(a+c))
symbol_names = 'a c d'
symbol_values = 'fluid_mass1 fluid_mass0 bh_report'
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 6500
solve_type = NEWTON
[TimeStepper]
type = FunctionDT
function = dts
[]
[]
[Outputs]
file_base = bh05
exodus = false
csv = true
execute_on = timestep_end
[]
(modules/porous_flow/test/tests/dirackernels/bh05.i)Further commentary, and the results demonstrating that MOOSE is correct may be found in the page discussing the theory behind Dirac Kernels in PorousFlow.
Comparison with a steady-state 2D analytic solution
The test
# Comparison with analytical solution for cylindrically-symmetric situation
[Mesh]
type = FileMesh
file = bh07_input.e
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Functions]
[dts]
type = PiecewiseLinear
y = '1000 10000'
x = '100 1000'
[]
[]
[Variables]
[pp]
initial_condition = 1E7
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[fflux]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
gravity = '0 0 0'
[]
[]
[BCs]
[fix_outer]
type = DirichletBC
boundary = perimeter
variable = pp
value = 1E7
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp'
number_fluid_phases = 1
number_fluid_components = 1
[]
[borehole_total_outflow_mass]
type = PorousFlowSumQuantity
[]
[pc]
type = PorousFlowCapillaryPressureVG
m = 0.8
alpha = 1e-5
[]
[]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
viscosity = 1e-3
density0 = 1000
thermal_expansion = 0
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseP
porepressure = pp
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1E-11 0 0 0 1E-11 0 0 0 1E-11'
[]
[relperm]
type = PorousFlowRelativePermeabilityFLAC
m = 2
phase = 0
[]
[]
[DiracKernels]
[bh]
type = PorousFlowPeacemanBorehole
variable = pp
SumQuantityUO = borehole_total_outflow_mass
point_file = bh07.bh
fluid_phase = 0
bottom_p_or_t = 0
unit_weight = '0 0 0'
use_mobility = true
re_constant = 0.1594 # use Chen and Zhang version
character = 2 # double the strength because bh07.bh only fills half the mesh
[]
[]
[Postprocessors]
[bh_report]
type = PorousFlowPlotQuantity
uo = borehole_total_outflow_mass
execute_on = 'initial timestep_end'
[]
[fluid_mass]
type = PorousFlowFluidMass
execute_on = 'initial timestep_end'
[]
[]
[VectorPostprocessors]
[pp]
type = LineValueSampler
variable = pp
start_point = '0 0 0'
end_point = '300 0 0'
sort_by = x
num_points = 300
execute_on = timestep_end
[]
[]
[Preconditioning]
[usual]
type = SMP
full = true
petsc_options = '-snes_converged_reason'
petsc_options_iname = '-ksp_type -pc_type -snes_atol -snes_rtol -snes_max_it -ksp_max_it'
petsc_options_value = 'bcgs bjacobi 1E-10 1E-10 10000 30'
[]
[]
[Executioner]
type = Transient
end_time = 1E3
solve_type = NEWTON
[TimeStepper]
# get only marginally better results for smaller time steps
type = FunctionDT
function = dts
[]
[]
[Outputs]
file_base = bh07
[along_line]
type = CSV
execute_on = final
[]
[exodus]
type = Exodus
execute_on = 'initial final'
[]
[]
(modules/porous_flow/test/tests/dirackernels/bh07.i)checks whether the Peaceman borehole extracts fluid at the correct rate and the corresponding reduction in porepressure agrees with an analytic solution of Laplace's equation. Details and results demonstrating the correctness of PorousFlow can be found in the page discussing the theory behind Dirac Kernels in PorousFlow.