- brine_fpThe name of the user object for brine
C++ Type:UserObjectName
Description:The name of the user object for brine
- capillary_pressureName of the UserObject defining the capillary pressure
C++ Type:UserObjectName
Description:Name of the UserObject defining the capillary pressure
- co2_fpThe name of the user object for CO2
C++ Type:UserObjectName
Description:The name of the user object for CO2
PorousFlowBrineCO2
Fluid state class for brine and CO2
A high precision equation of state for brine and CO, including the mutual solubility of CO into the liquid brine and water vapor into the CO-rich gas phase using the accurate fugacity-based formulation of Spycher et al. (2003) and Spycher et al. (2005).
This model is suitable for simulations of geological storage of CO in saline aquifers.
For more details, see the documentation of the brine and CO equation of state.
Input Parameters
- execute_onTIMESTEP_ENDThe list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM.
Default:TIMESTEP_END
C++ Type:ExecFlagEnum
Options:NONE INITIAL LINEAR NONLINEAR TIMESTEP_END TIMESTEP_BEGIN FINAL CUSTOM
Description:The list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM.
- liquid_fluid_component0The fluid component number of the primary liquid component
Default:0
C++ Type:unsigned int
Options:
Description:The fluid component number of the primary liquid component
- liquid_phase_number0The phase number of the liquid phase
Default:0
C++ Type:unsigned int
Options:
Description:The phase number of the liquid phase
- salt_component2The component number of salt
Default:2
C++ Type:unsigned int
Options:
Description:The component number of salt
Optional Parameters
- allow_duplicate_execution_on_initialFalseIn the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
Default:False
C++ Type:bool
Options:
Description:In the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
- force_preauxFalseForces the GeneralUserObject to be executed in PREAUX
Default:False
C++ Type:bool
Options:
Description:Forces the GeneralUserObject to be executed in PREAUX
- 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
Options:
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
Input Files
- modules/porous_flow/test/tests/jacobian/brineco2_twophase.i
- modules/porous_flow/examples/co2_intercomparison/1Dradial/1Dradial.i
- modules/porous_flow/examples/lava_lamp/2phase_convection.i
- modules/porous_flow/test/tests/fluidstate/brineco2_2.i
- modules/porous_flow/test/tests/fluidstate/brineco2_ic.i
- modules/porous_flow/test/tests/jacobian/brineco2_liquid_2.i
- modules/porous_flow/test/tests/jacobian/brineco2_twophase_nonisothermal.i
- modules/porous_flow/test/tests/fluidstate/brineco2_nonisothermal.i
- modules/porous_flow/test/tests/jacobian/brineco2_gas.i
- modules/porous_flow/test/tests/fluidstate/theis_brineco2.i
- modules/porous_flow/examples/co2_intercomparison/1Dradial/properties.i
- modules/porous_flow/examples/lava_lamp/1phase_convection.i
- modules/porous_flow/test/tests/fluidstate/brineco2.i
- modules/porous_flow/test/tests/fluidstate/brineco2_hightemp.i
- modules/porous_flow/test/tests/jacobian/brineco2_liquid.i
- modules/porous_flow/test/tests/fluidstate/theis_brineco2_nonisothermal.i
modules/porous_flow/test/tests/jacobian/brineco2_twophase.i
# Tests correct calculation of properties derivatives in PorousFlowFluidState
# for conditions that are appropriate for two phases
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.05
[../]
[]
[Variables]
[./pgas]
[../]
[./zi]
[../]
[]
[ICs]
[./pgas]
type = RandomIC
min = 1e6
max = 4e6
variable = pgas
seed = 1
[../]
[./z]
type = RandomIC
min = 0.2
max = 0.8
variable = zi
seed = 2
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = zi
fluid_component = 1
[../]
[./adv0]
type = PorousFlowAdvectiveFlux
variable = pgas
fluid_component = 0
[../]
[./adv1]
type = PorousFlowAdvectiveFlux
variable = zi
fluid_component = 1
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1e1
pc_max = 1e4
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[./water]
type = Water97FluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = 50
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[AuxVariables]
[./sgas]
family = MONOMIAL
order = CONSTANT
[../]
[]
[AuxKernels]
[./sgas]
type = PorousFlowPropertyAux
property = saturation
phase = 1
variable = sgas
[../]
[]
[Postprocessors]
[./sgas_min]
type = ElementExtremeValue
variable = sgas
value_type = min
[../]
[./sgas_max]
type = ElementExtremeValue
variable = sgas
value_type = max
[../]
[]
modules/porous_flow/examples/co2_intercomparison/1Dradial/1Dradial.i
# Intercomparison problem 3: Radial flow from an injection well
#
# From Pruess et al, Code intercomparison builds confidence in
# numerical simulation models for geologic disposal of CO2, Energy 29 (2004)
#
# A variation with zero salinity can be run by changing the initial condition
# of the AuxVariable xnacl
[Mesh]
type = GeneratedMesh
dim = 1
nx = 500
xmax = 10000
bias_x = 1.01
[]
[Problem]
type = FEProblem
coord_type = 'RZ'
rz_coord_axis = Y
[]
[GlobalParams]
PorousFlowDictator = 'dictator'
gravity = '0 0 0'
[]
[AuxVariables]
[pressure_liquid]
order = CONSTANT
family = MONOMIAL
[]
[saturation_gas]
order = CONSTANT
family = MONOMIAL
[]
[x1]
order = CONSTANT
family = MONOMIAL
[]
[y0]
order = CONSTANT
family = MONOMIAL
[]
[xnacl]
initial_condition = 0.15
[]
[]
[AuxKernels]
[pressure_liquid]
type = PorousFlowPropertyAux
variable = pressure_liquid
property = pressure
phase = 0
execute_on = 'timestep_end'
[]
[saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = 'timestep_end'
[]
[x1]
type = PorousFlowPropertyAux
variable = x1
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = 'timestep_end'
[]
[y0]
type = PorousFlowPropertyAux
variable = y0
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = 'timestep_end'
[]
[]
[Variables]
[pgas]
initial_condition = 12e6
[]
[zi]
initial_condition = 0
scaling = 1e4
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pgas
[]
[flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pgas
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = zi
[]
[flux1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = zi
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 3
[]
[pc]
type = PorousFlowCapillaryPressureVG
alpha = 5.099e-5
m = 0.457
sat_lr = 0.0
pc_max = 1e7
[]
[fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[]
[]
[Modules]
[FluidProperties]
[co2sw]
type = CO2FluidProperties
[]
[co2]
type = TabulatedFluidProperties
fp = co2sw
[]
[water]
type = Water97FluidProperties
[]
[watertab]
type = TabulatedFluidProperties
fp = water
temperature_min = 273.15
temperature_max = 573.15
fluid_property_file = water_fluid_properties.csv
save_file = false
[]
[brine]
type = BrineFluidProperties
water_fp = watertab
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = '45'
[]
[brineco2]
type = PorousFlowFluidState
gas_porepressure = 'pgas'
z = 'zi'
temperature_unit = Celsius
xnacl = 'xnacl'
capillary_pressure = pc
fluid_state = fs
[]
[porosity]
type = PorousFlowPorosityConst
porosity = '0.12'
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-13 0 0 0 1e-13 0 0 0 1e-13'
[]
[relperm_water]
type = PorousFlowRelativePermeabilityVG
m = 0.457
phase = 0
s_res = 0.3
sum_s_res = 0.35
[]
[relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 1
s_res = 0.05
sum_s_res = 0.35
[]
[]
[BCs]
[rightwater]
type = PorousFlowPiecewiseLinearSink
boundary = 'right'
variable = pgas
use_mobility = true
PorousFlowDictator = dictator
fluid_phase = 0
multipliers = '0 1e9'
PT_shift = '12e6'
pt_vals = '0 1e9'
mass_fraction_component = 0
use_relperm = true
[]
[rightco2]
type = PorousFlowPiecewiseLinearSink
variable = zi
boundary = 'right'
use_mobility = true
PorousFlowDictator = dictator
fluid_phase = 1
multipliers = '0 1e9'
PT_shift = '12e6'
pt_vals = '0 1e9'
mass_fraction_component = 1
use_relperm = true
[]
[]
[DiracKernels]
[source]
type = PorousFlowSquarePulsePointSource
point = '0 0 0'
mass_flux = 1
variable = zi
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'gmres bjacobi lu NONZERO'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 8.64e8
nl_max_its = 25
l_max_its = 100
dtmax = 5e6
[TimeStepper]
type = IterationAdaptiveDT
dt = 100
[]
[]
[VectorPostprocessors]
[vars]
type = NodalValueSampler
sort_by = x
variable = 'pgas zi xnacl'
execute_on = 'timestep_end'
outputs = spatial
[]
[auxvars]
type = ElementValueSampler
sort_by = x
variable = 'saturation_gas x1 y0'
execute_on = 'timestep_end'
outputs = spatial
[]
[]
[Postprocessors]
[pgas]
type = PointValue
point = '25.25 0 0'
variable = pgas
outputs = time
[]
[sgas]
type = PointValue
point = '25.25 0 0'
variable = saturation_gas
outputs = time
[]
[zi]
type = PointValue
point = '25.25 0 0'
variable = zi
outputs = time
[]
[massgas]
type = PorousFlowFluidMass
fluid_component = 1
outputs = time
[]
[x1]
type = PointValue
point = '25.25 0 0'
variable = x1
outputs = time
[]
[y0]
type = PointValue
point = '25.25 0 0'
variable = y0
outputs = time
[]
[xnacl]
type = PointValue
point = '25.25 0 0'
variable = xnacl
outputs = time
[]
[]
[Outputs]
print_linear_residuals = false
perf_graph = true
sync_times = '2.592e6 8.64e6 8.64e7 8.64e8'
[time]
type = CSV
[]
[spatial]
type = CSV
sync_only = true
[]
[]
modules/porous_flow/examples/lava_lamp/2phase_convection.i
# Two phase density-driven convection of dissolved CO2 in brine
#
# Initially, the model has a gas phase at the top with a saturation of 0.29
# (which corresponds to an initial value of zi = 0.2).
# Diffusion of the dissolved CO2
# component from the saturated liquid to the unsaturated liquid below reduces the
# amount of CO2 in the gas phase. As the density of the CO2-saturated brine is greater
# than the unsaturated brine, a gravitational instability arises and density-driven
# convection of CO2-rich fingers descend into the unsaturated brine.
#
# The instability is seeded by a random perturbation to the porosity field.
# Mesh adaptivity is used to refine the mesh as the fingers form.
[GlobalParams]
PorousFlowDictator = 'dictator'
gravity = '0 -9.81 0'
[]
[Adaptivity]
max_h_level = 2
marker = marker
initial_marker = initial
initial_steps = 2
[./Indicators]
[./indicator]
type = GradientJumpIndicator
variable = zi
[../]
[../]
[./Markers]
[./marker]
type = ErrorFractionMarker
indicator = indicator
refine = 0.8
[../]
[./initial]
type = BoxMarker
bottom_left = '0 1.95 0'
top_right = '2 2 0'
inside = REFINE
outside = DO_NOTHING
[../]
[../]
[]
[Mesh]
type = GeneratedMesh
dim = 2
ymax = 2
xmax = 2
ny = 40
nx = 40
bias_y = 0.95
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.01
[../]
[./saturation_gas]
order = FIRST
family = MONOMIAL
[../]
[./xco2l]
order = FIRST
family = MONOMIAL
[../]
[./density_liquid]
order = FIRST
family = MONOMIAL
[../]
[./porosity]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = 'timestep_end'
[../]
[./xco2l]
type = PorousFlowPropertyAux
variable = xco2l
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = 'timestep_end'
[../]
[./density_liquid]
type = PorousFlowPropertyAux
variable = density_liquid
property = density
phase = 0
execute_on = 'timestep_end'
[../]
[]
[Variables]
[./pgas]
[../]
[./zi]
scaling = 1e4
[../]
[]
[ICs]
[./pressure]
type = FunctionIC
function = 10e6-9.81*1000*y
variable = pgas
[../]
[./zi]
type = BoundingBoxIC
variable = zi
x1 = 0
x2 = 2
y1 = 1.95
y2 = 2
inside = 0.2
outside = 0
[../]
[./porosity]
type = RandomIC
variable = porosity
min = 0.25
max = 0.275
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pgas
[../]
[./flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pgas
[../]
[./diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pgas
disp_long = '0 0'
disp_trans = '0 0'
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = zi
[../]
[./flux1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = zi
[../]
[./diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = zi
disp_long = '0 0'
disp_trans = '0 0'
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2sw]
type = CO2FluidProperties
[../]
[./co2]
type = TabulatedFluidProperties
fp = co2sw
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = '45'
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = 'pgas'
z = 'zi'
temperature_unit = Celsius
xnacl = 'xnacl'
capillary_pressure = pc
fluid_state = fs
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = porosity
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-11 0 0 0 1e-11 0 0 0 1e-11'
[../]
[./relperm_water]
type = PorousFlowRelativePermeabilityCorey
phase = 0
n = 2
s_res = 0.1
sum_s_res = 0.2
[../]
[./relperm_gas]
type = PorousFlowRelativePermeabilityCorey
phase = 1
n = 2
s_res = 0.1
sum_s_res = 0.2
[../]
[./diffusivity]
type = PorousFlowDiffusivityConst
diffusion_coeff = '2e-9 2e-9 2e-9 2e-9'
tortuosity = '1 1'
[../]
[]
[Preconditioning]
active = basic
[./mumps_is_best_for_parallel_jobs]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[./basic]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = 'gmres asm lu NONZERO 2 '
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 1e6
nl_max_its = 25
l_max_its = 100
dtmax = 1e4
nl_abs_tol = 1e-6
[./TimeStepper]
type = IterationAdaptiveDT
dt = 10
growth_factor = 2
cutback_factor = 0.5
[../]
[]
[Functions]
[./flux]
type = ParsedFunction
vals = 'delta_xco2 dt'
vars = 'dx dt'
value = 'dx/dt'
[../]
[]
[Postprocessors]
[./total_co2_in_gas]
type = PorousFlowFluidMass
phase = 1
fluid_component = 1
[../]
[./total_co2_in_liquid]
type = PorousFlowFluidMass
phase = 0
fluid_component = 1
[../]
[./numdofs]
type = NumDOFs
[../]
[./delta_xco2]
type = ChangeOverTimePostprocessor
postprocessor = total_co2_in_liquid
[../]
[./dt]
type = TimestepSize
[../]
[./flux]
type = FunctionValuePostprocessor
function = flux
[../]
[]
[Outputs]
print_linear_residuals = false
perf_graph = true
exodus = true
csv = true
[]
modules/porous_flow/test/tests/fluidstate/brineco2_2.i
# Injection of supercritical CO2 into a single brine saturated cell. The CO2 initially fully
# dissolves into the brine, increasing its density slightly. After a few time steps,
# the brine is saturated with CO2, and subsequently a supercritical gas phase of CO2 saturated
# with a small amount of H2O is formed. Salt is included as a nonlinear variable.
[Mesh]
type = GeneratedMesh
dim = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
temperature = 30
[]
[Variables]
[./pgas]
initial_condition = 20e6
[../]
[./z]
[../]
[./xnacl]
initial_condition = 0.1
[../]
[]
[DiracKernels]
[./source]
type = PorousFlowSquarePulsePointSource
variable = z
point = '0.5 0.5 0'
mass_flux = 2
[../]
[]
[BCs]
[./left]
type = DirichletBC
value = 20e6
variable = pgas
boundary = left
[../]
[./right]
type = DirichletBC
value = 20e6
variable = pgas
boundary = right
[../]
[]
[AuxVariables]
[./pressure_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./pressure_water]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_water]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_gas]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./pressure_water]
type = PorousFlowPropertyAux
variable = pressure_water
property = pressure
phase = 0
execute_on = 'initial timestep_end'
[../]
[./pressure_gas]
type = PorousFlowPropertyAux
variable = pressure_gas
property = pressure
phase = 1
execute_on = 'initial timestep_end'
[../]
[./saturation_water]
type = PorousFlowPropertyAux
variable = saturation_water
property = saturation
phase = 0
execute_on = 'initial timestep_end'
[../]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = 'initial timestep_end'
[../]
[./density_water]
type = PorousFlowPropertyAux
variable = density_water
property = density
phase = 0
execute_on = 'initial timestep_end'
[../]
[./density_gas]
type = PorousFlowPropertyAux
variable = density_gas
property = density
phase = 1
execute_on = 'initial timestep_end'
[../]
[./viscosity_water]
type = PorousFlowPropertyAux
variable = viscosity_water
property = viscosity
phase = 0
execute_on = 'initial timestep_end'
[../]
[./viscosity_gas]
type = PorousFlowPropertyAux
variable = viscosity_gas
property = viscosity
phase = 1
execute_on = 'initial timestep_end'
[../]
[./x1_water]
type = PorousFlowPropertyAux
variable = x1_water
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = 'initial timestep_end'
[../]
[./x1_gas]
type = PorousFlowPropertyAux
variable = x1_gas
property = mass_fraction
phase = 1
fluid_component = 1
execute_on = 'initial timestep_end'
[../]
[./x0_water]
type = PorousFlowPropertyAux
variable = x0_water
property = mass_fraction
phase = 0
fluid_component = 0
execute_on = 'initial timestep_end'
[../]
[./x0_gas]
type = PorousFlowPropertyAux
variable = x0_gas
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = 'initial timestep_end'
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = z
fluid_component = 1
[../]
[./mass2]
type = PorousFlowMassTimeDerivative
variable = xnacl
fluid_component = 2
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas z xnacl'
number_fluid_phases = 2
number_fluid_components = 3
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = z
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 10
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
[./density_water]
type = ElementIntegralVariablePostprocessor
variable = density_water
execute_on = 'initial timestep_end'
[../]
[./density_gas]
type = ElementIntegralVariablePostprocessor
variable = density_gas
execute_on = 'initial timestep_end'
[../]
[./viscosity_water]
type = ElementIntegralVariablePostprocessor
variable = viscosity_water
execute_on = 'initial timestep_end'
[../]
[./viscosity_gas]
type = ElementIntegralVariablePostprocessor
variable = viscosity_gas
execute_on = 'initial timestep_end'
[../]
[./x1_water]
type = ElementIntegralVariablePostprocessor
variable = x1_water
execute_on = 'initial timestep_end'
[../]
[./x0_water]
type = ElementIntegralVariablePostprocessor
variable = x0_water
execute_on = 'initial timestep_end'
[../]
[./x1_gas]
type = ElementIntegralVariablePostprocessor
variable = x1_gas
execute_on = 'initial timestep_end'
[../]
[./x0_gas]
type = ElementIntegralVariablePostprocessor
variable = x0_gas
execute_on = 'initial timestep_end'
[../]
[./sg]
type = ElementIntegralVariablePostprocessor
variable = saturation_gas
execute_on = 'initial timestep_end'
[../]
[./sw]
type = ElementIntegralVariablePostprocessor
variable = saturation_water
execute_on = 'initial timestep_end'
[../]
[./pwater]
type = ElementIntegralVariablePostprocessor
variable = pressure_water
execute_on = 'initial timestep_end'
[../]
[./pgas]
type = ElementIntegralVariablePostprocessor
variable = pressure_gas
execute_on = 'initial timestep_end'
[../]
[./xnacl]
type = ElementIntegralVariablePostprocessor
variable = xnacl
execute_on = 'initial timestep_end'
[../]
[./x0mass]
type = PorousFlowFluidMass
fluid_component = 0
phase = '0 1'
execute_on = 'initial timestep_end'
[../]
[./x1mass]
type = PorousFlowFluidMass
fluid_component = 1
phase = '0 1'
execute_on = 'initial timestep_end'
[../]
[./x2mass]
type = PorousFlowFluidMass
fluid_component = 2
phase = '0 1'
execute_on = 'initial timestep_end'
[../]
[]
[Outputs]
csv = true
file_base = brineco2_2
execute_on = 'initial timestep_end'
perf_graph = true
[]
modules/porous_flow/test/tests/fluidstate/brineco2_ic.i
# Tests correct calculation of z (total mass fraction of NCG summed over all
# phases) using the PorousFlowFluidStateIC initial condition. Once z is
# calculated by the initial condition, the thermophysical properties are calculated
# and the resulting gas saturation should be equal to that given in the intial condition
[Mesh]
type = GeneratedMesh
dim = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
temperature_unit = Celsius
[]
[Variables]
[./pgas]
initial_condition = 1e6
[../]
[./z]
[../]
[]
[ICs]
[./z]
type = PorousFlowFluidStateIC
saturation = 0.5
gas_porepressure = pgas
temperature = 50
variable = z
xnacl = 0.1
fluid_state = fs
[../]
[]
[AuxVariables]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_water]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./saturation_water]
type = PorousFlowPropertyAux
variable = saturation_water
property = saturation
phase = 0
execute_on = timestep_end
[../]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = timestep_end
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = z
fluid_component = 1
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas z'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = 50
[../]
[./waterncg]
type = PorousFlowFluidState
gas_porepressure = pgas
z = z
fluid_state = fs
capillary_pressure = pc
xnacl = 0.1
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
[./sg]
type = ElementIntegralVariablePostprocessor
variable = saturation_gas
execute_on = 'initial timestep_end'
[../]
[./sw]
type = ElementIntegralVariablePostprocessor
variable = saturation_water
execute_on = 'initial timestep_end'
[../]
[./z]
type = ElementIntegralVariablePostprocessor
variable = z
execute_on = 'initial timestep_end'
[../]
[]
[Outputs]
csv = true
[]
modules/porous_flow/test/tests/jacobian/brineco2_liquid_2.i
# Tests correct calculation of properties derivatives in PorousFlowFluidState
# for conditions that give a single liquid phase, including salt as a nonlinear variable
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[./pgas]
[../]
[./zi]
[../]
[./xnacl]
[../]
[]
[ICs]
[./pgas]
type = RandomIC
min = 5e6
max = 8e6
variable = pgas
[../]
[./z]
type = RandomIC
min = 0.01
max = 0.03
variable = zi
[../]
[./xnacl]
type = RandomIC
min = 0.01
max = 0.15
variable = xnacl
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = zi
fluid_component = 1
[../]
[./mass2]
type = PorousFlowMassTimeDerivative
variable = xnacl
fluid_component = 2
[../]
[./adv0]
type = PorousFlowAdvectiveFlux
variable = pgas
fluid_component = 0
[../]
[./adv1]
type = PorousFlowAdvectiveFlux
variable = zi
fluid_component = 1
[../]
[./adv2]
type = PorousFlowAdvectiveFlux
variable = xnacl
fluid_component = 2
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi xnacl'
number_fluid_phases = 2
number_fluid_components = 3
[../]
[./pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
pc_max = 1e3
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = 50
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
modules/porous_flow/test/tests/jacobian/brineco2_twophase_nonisothermal.i
# Tests correct calculation of properties derivatives in PorousFlowFluidState
# for nonisothermal two phase conditions, including salt as a nonlinear variable
[Mesh]
[./mesh]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 2
[../]
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[./pgas]
[../]
[./zi]
[../]
[./xnacl]
[../]
[./temperature]
[../]
[]
[ICs]
[./pgas]
type = RandomIC
min = 1e6
max = 4e6
variable = pgas
seed = 1
[../]
[./z]
type = RandomIC
min = 0.2
max = 0.8
variable = zi
seed = 1
[../]
[./xnacl]
type = RandomIC
min = 0.01
max = 0.15
variable = xnacl
seed = 1
[../]
[./temperature]
type = RandomIC
min = 20
max = 80
variable = temperature
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = zi
fluid_component = 1
[../]
[./mass2]
type = PorousFlowMassTimeDerivative
variable = xnacl
fluid_component = 2
[../]
[./adv0]
type = PorousFlowAdvectiveFlux
variable = pgas
fluid_component = 0
[../]
[./adv1]
type = PorousFlowAdvectiveFlux
variable = zi
fluid_component = 1
[../]
[./adv2]
type = PorousFlowAdvectiveFlux
variable = xnacl
fluid_component = 2
[../]
[./energy]
type = PorousFlowEnergyTimeDerivative
variable = temperature
[../]
[./heat]
type = PorousFlowHeatAdvection
variable = temperature
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi xnacl temperature'
number_fluid_phases = 2
number_fluid_components = 3
[../]
[./pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
pc_max = 1e3
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = temperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature = temperature
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[./rock_heat]
type = PorousFlowMatrixInternalEnergy
specific_heat_capacity = 1000
density = 2500
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
automatic_scaling = true
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
modules/porous_flow/test/tests/fluidstate/brineco2_nonisothermal.i
[Mesh]
[./mesh]
type = GeneratedMeshGenerator
dim = 2
[../]
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[./pgas]
initial_condition = 20e6
[../]
[./z]
initial_condition = 0.2
[../]
[./temperature]
initial_condition = 70
[../]
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.1
[../]
[./pressure_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./pressure_water]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_water]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./enthalpy_water]
order = CONSTANT
family = MONOMIAL
[../]
[./enthalpy_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./internal_energy_water]
order = CONSTANT
family = MONOMIAL
[../]
[./internal_energy_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_gas]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./pressure_water]
type = PorousFlowPropertyAux
variable = pressure_water
property = pressure
phase = 0
execute_on = timestep_end
[../]
[./pressure_gas]
type = PorousFlowPropertyAux
variable = pressure_gas
property = pressure
phase = 1
execute_on = timestep_end
[../]
[./saturation_water]
type = PorousFlowPropertyAux
variable = saturation_water
property = saturation
phase = 0
execute_on = timestep_end
[../]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = timestep_end
[../]
[./density_water]
type = PorousFlowPropertyAux
variable = density_water
property = density
phase = 0
execute_on = timestep_end
[../]
[./density_gas]
type = PorousFlowPropertyAux
variable = density_gas
property = density
phase = 1
execute_on = timestep_end
[../]
[./viscosity_water]
type = PorousFlowPropertyAux
variable = viscosity_water
property = viscosity
phase = 0
execute_on = timestep_end
[../]
[./viscosity_gas]
type = PorousFlowPropertyAux
variable = viscosity_gas
property = viscosity
phase = 1
execute_on = timestep_end
[../]
[./enthalpy_water]
type = PorousFlowPropertyAux
variable = enthalpy_water
property = enthalpy
phase = 0
execute_on = timestep_end
[../]
[./enthalpy_gas]
type = PorousFlowPropertyAux
variable = enthalpy_gas
property = enthalpy
phase = 1
execute_on = timestep_end
[../]
[./internal_energy_water]
type = PorousFlowPropertyAux
variable = internal_energy_water
property = internal_energy
phase = 0
execute_on = timestep_end
[../]
[./internal_energy_gas]
type = PorousFlowPropertyAux
variable = internal_energy_gas
property = internal_energy
phase = 1
execute_on = timestep_end
[../]
[./x1_water]
type = PorousFlowPropertyAux
variable = x1_water
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = timestep_end
[../]
[./x1_gas]
type = PorousFlowPropertyAux
variable = x1_gas
property = mass_fraction
phase = 1
fluid_component = 1
execute_on = timestep_end
[../]
[./x0_water]
type = PorousFlowPropertyAux
variable = x0_water
property = mass_fraction
phase = 0
fluid_component = 0
execute_on = timestep_end
[../]
[./x0_gas]
type = PorousFlowPropertyAux
variable = x0_gas
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = timestep_end
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = z
fluid_component = 1
[../]
[./heat]
type = TimeDerivative
variable = temperature
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas z temperature'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = z
temperature = temperature
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
[./density_water]
type = ElementIntegralVariablePostprocessor
variable = density_water
[../]
[./density_gas]
type = ElementIntegralVariablePostprocessor
variable = density_gas
[../]
[./viscosity_water]
type = ElementIntegralVariablePostprocessor
variable = viscosity_water
[../]
[./viscosity_gas]
type = ElementIntegralVariablePostprocessor
variable = viscosity_gas
[../]
[./enthalpy_water]
type = ElementIntegralVariablePostprocessor
variable = enthalpy_water
[../]
[./enthalpy_gas]
type = ElementIntegralVariablePostprocessor
variable = enthalpy_gas
[../]
[./internal_energy_water]
type = ElementIntegralVariablePostprocessor
variable = internal_energy_water
[../]
[./internal_energy_gas]
type = ElementIntegralVariablePostprocessor
variable = internal_energy_gas
[../]
[./x1_water]
type = ElementIntegralVariablePostprocessor
variable = x1_water
[../]
[./x0_water]
type = ElementIntegralVariablePostprocessor
variable = x0_water
[../]
[./x1_gas]
type = ElementIntegralVariablePostprocessor
variable = x1_gas
[../]
[./x0_gas]
type = ElementIntegralVariablePostprocessor
variable = x0_gas
[../]
[./sg]
type = ElementIntegralVariablePostprocessor
variable = saturation_gas
[../]
[./sw]
type = ElementIntegralVariablePostprocessor
variable = saturation_water
[../]
[./pwater]
type = ElementIntegralVariablePostprocessor
variable = pressure_water
[../]
[./pgas]
type = ElementIntegralVariablePostprocessor
variable = pressure_gas
[../]
[./x0mass]
type = PorousFlowFluidMass
fluid_component = 0
phase = '0 1'
[../]
[./x1mass]
type = PorousFlowFluidMass
fluid_component = 1
phase = '0 1'
[../]
[]
[Outputs]
csv = true
execute_on = timestep_end
[]
modules/porous_flow/test/tests/jacobian/brineco2_gas.i
# Tests correct calculation of properties derivatives in PorousFlowFluidState
# for conditions that give a single gas phase
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.05
[../]
[]
[Variables]
[./pgas]
[../]
[./zi]
[../]
[]
[ICs]
[./pgas]
type = RandomIC
min = 5e4
max = 1e5
variable = pgas
[../]
[./z]
type = RandomIC
min = 0.9
max = 0.99
variable = zi
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = zi
fluid_component = 1
[../]
[./adv0]
type = PorousFlowAdvectiveFlux
variable = pgas
fluid_component = 0
[../]
[./adv1]
type = PorousFlowAdvectiveFlux
variable = zi
fluid_component = 1
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
pc_max = 1e3
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[./water]
type = Water97FluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = 50
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[AuxVariables]
[./sgas]
family = MONOMIAL
order = CONSTANT
[../]
[]
[AuxKernels]
[./sgas]
type = PorousFlowPropertyAux
property = saturation
phase = 1
variable = sgas
[../]
[]
[Postprocessors]
[./sgas_min]
type = ElementExtremeValue
variable = sgas
value_type = min
[../]
[./sgas_max]
type = ElementExtremeValue
variable = sgas
value_type = max
[../]
[]
modules/porous_flow/test/tests/fluidstate/theis_brineco2.i
# Two phase Theis problem: Flow from single source.
# Constant rate injection 2 kg/s
# 1D cylindrical mesh
# Initially, system has only a liquid phase, until enough gas is injected
# to form a gas phase, in which case the system becomes two phase.
#
# This test takes a few minutes to run, so is marked heavy
[Mesh]
type = GeneratedMesh
dim = 1
nx = 2000
xmax = 2000
[]
[Problem]
type = FEProblem
coord_type = RZ
rz_coord_axis = Y
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[AuxVariables]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x1]
order = CONSTANT
family = MONOMIAL
[../]
[./y0]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = timestep_end
[../]
[./x1]
type = PorousFlowPropertyAux
variable = x1
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = timestep_end
[../]
[./y0]
type = PorousFlowPropertyAux
variable = y0
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = timestep_end
[../]
[]
[Variables]
[./pgas]
initial_condition = 20e6
[../]
[./zi]
initial_condition = 0
[../]
[./xnacl]
initial_condition = 0.1
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pgas
[../]
[./flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pgas
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = zi
[../]
[./flux1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = zi
[../]
[./mass2]
type = PorousFlowMassTimeDerivative
fluid_component = 2
variable = xnacl
[../]
[./flux2]
type = PorousFlowAdvectiveFlux
fluid_component = 2
variable = xnacl
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi xnacl'
number_fluid_phases = 2
number_fluid_components = 3
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2sw]
type = CO2FluidProperties
[../]
[./co2]
type = TabulatedFluidProperties
fp = co2sw
[../]
[./water]
type = Water97FluidProperties
[../]
[./watertab]
type = TabulatedFluidProperties
fp = water
temperature_min = 273.15
temperature_max = 573.15
fluid_property_file = water_fluid_properties.csv
save_file = false
[../]
[./brine]
type = BrineFluidProperties
water_fp = watertab
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./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 = 2
phase = 0
s_res = 0.1
sum_s_res = 0.1
[../]
[./relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 1
[../]
[]
[BCs]
[./rightwater]
type = DirichletBC
boundary = right
value = 20e6
variable = pgas
[../]
[]
[DiracKernels]
[./source]
type = PorousFlowSquarePulsePointSource
point = '0 0 0'
mass_flux = 2
variable = zi
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 1e5
[./TimeStepper]
type = IterationAdaptiveDT
dt = 1
growth_factor = 1.5
[../]
[]
[VectorPostprocessors]
[./line]
type = LineValueSampler
sort_by = x
start_point = '0 0 0'
end_point = '2000 0 0'
num_points = 10000
variable = 'pgas zi xnacl x1 saturation_gas'
execute_on = 'timestep_end'
[../]
[]
[Postprocessors]
[./pgas]
type = PointValue
point = '4 0 0'
variable = pgas
[../]
[./sgas]
type = PointValue
point = '4 0 0'
variable = saturation_gas
[../]
[./zi]
type = PointValue
point = '4 0 0'
variable = zi
[../]
[./massgas]
type = PorousFlowFluidMass
fluid_component = 1
[../]
[./x1]
type = PointValue
point = '4 0 0'
variable = x1
[../]
[./y0]
type = PointValue
point = '4 0 0'
variable = y0
[../]
[./xnacl]
type = PointValue
point = '4 0 0'
variable = xnacl
[../]
[]
[Outputs]
print_linear_residuals = false
perf_graph = true
[./csvout]
type = CSV
execute_on = timestep_end
execute_vector_postprocessors_on = final
[../]
[]
modules/porous_flow/examples/co2_intercomparison/1Dradial/properties.i
# Liquid and gas properties for code intercomparison problem 3
#
# From Pruess et al, Code intercomparison builds confidence in
# numerical simulation models for geologic disposal of CO2, Energy 29 (2004)
#
# This test simply calculates density and viscosity of each phase for
# various pressures and salinities, as well as mass fractions of CO2 in the
# liquid phase and H2O in the gas phase.
#
# Four versions of this are run:
# 1) No CO2, 0 salt mass fraction (pure water)
# 2) Enough CO2 to form gas phase, 0 salt mass fraction (pure water)
# 3) No CO2, 0.15 salt mass fraction
# 4) Enough CO2 to form gas phase, 0.15 salt mass fraction
#
# These results compare well with detailed results presented in Pruess et al,
# Intercomparison of numerical simulation codes for geologic disposal of CO2,
# LBNL-51813 (2002)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
xmax = 4
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[AuxVariables]
[density_liquid]
order = CONSTANT
family = MONOMIAL
[]
[density_gas]
order = CONSTANT
family = MONOMIAL
[]
[viscosity_liquid]
order = CONSTANT
family = MONOMIAL
[]
[viscosity_gas]
order = CONSTANT
family = MONOMIAL
[]
[x1]
order = CONSTANT
family = MONOMIAL
[]
[y0]
order = CONSTANT
family = MONOMIAL
[]
[xnacl]
initial_condition = 0.0
[]
[]
[AuxKernels]
[density_liquid]
type = PorousFlowPropertyAux
variable = density_liquid
property = density
phase = 0
execute_on = timestep_end
[]
[density_gas]
type = PorousFlowPropertyAux
variable = density_gas
property = density
phase = 1
execute_on = timestep_end
[]
[viscosity_liquid]
type = PorousFlowPropertyAux
variable = viscosity_liquid
property = viscosity
phase = 0
execute_on = timestep_end
[]
[viscosity_gas]
type = PorousFlowPropertyAux
variable = viscosity_gas
property = viscosity
phase = 1
execute_on = timestep_end
[]
[x1]
type = PorousFlowPropertyAux
variable = x1
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = timestep_end
[]
[y0]
type = PorousFlowPropertyAux
variable = y0
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = timestep_end
[]
[]
[Variables]
[pgas]
order = CONSTANT
family = MONOMIAL
[]
[zi]
initial_condition = 0.0
[]
[]
[Functions]
[pic]
type = ParsedFunction
value = 'if(x<1,12e6,if(x<2,16e6,if(x<3,20e6,24e6)))'
[]
[]
[ICs]
[pic]
type = FunctionIC
function = pic
variable = pgas
[]
[]
[Kernels]
[diffusionp]
type = NullKernel
variable = pgas
[]
[diffusionz]
type = NullKernel
variable = zi
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[]
[fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[]
[]
[Modules]
[FluidProperties]
[co2]
type = CO2FluidProperties
[]
[brine]
type = BrineFluidProperties
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = 45
[]
[brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = 'gmres asm lu NONZERO 2 '
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
perf_graph = true
csv = true
execute_on = timestep_end
file_base = properties_water
[]
[VectorPostprocessors]
[vpp]
type = ElementValueSampler
variable = 'pgas density_liquid density_gas viscosity_liquid viscosity_gas x1 y0'
sort_by = x
[]
[]
modules/porous_flow/examples/lava_lamp/1phase_convection.i
# Two phase density-driven convection of dissolved CO2 in brine
#
# The model starts with CO2 in the liquid phase only. The CO2 diffuses into the brine.
# As the density of the CO2-saturated brine is greater
# than the unsaturated brine, a gravitational instability arises and density-driven
# convection of CO2-rich fingers descend into the unsaturated brine.
#
# The instability is seeded by a random perturbation to the porosity field.
# Mesh adaptivity is used to refine the mesh as the fingers form.
[GlobalParams]
PorousFlowDictator = 'dictator'
gravity = '0 -9.81 0'
[]
[Adaptivity]
max_h_level = 2
marker = marker
initial_marker = initial
initial_steps = 2
[./Indicators]
[./indicator]
type = GradientJumpIndicator
variable = zi
[../]
[../]
[./Markers]
[./marker]
type = ErrorFractionMarker
indicator = indicator
refine = 0.8
[../]
[./initial]
type = BoxMarker
bottom_left = '0 1.95 0'
top_right = '2 2 0'
inside = REFINE
outside = DO_NOTHING
[../]
[../]
[]
[Mesh]
type = GeneratedMesh
dim = 2
ymin = 1.5
ymax = 2
xmax = 2
ny = 20
nx = 40
bias_y = 0.95
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.01
[../]
[./saturation_gas]
order = FIRST
family = MONOMIAL
[../]
[./xco2l]
order = FIRST
family = MONOMIAL
[../]
[./density_liquid]
order = FIRST
family = MONOMIAL
[../]
[./porosity]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = 'timestep_end'
[../]
[./xco2l]
type = PorousFlowPropertyAux
variable = xco2l
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = 'timestep_end'
[../]
[./density_liquid]
type = PorousFlowPropertyAux
variable = density_liquid
property = density
phase = 0
execute_on = 'timestep_end'
[../]
[]
[Variables]
[./pgas]
[../]
[./zi]
scaling = 1e4
[../]
[]
[ICs]
[./pressure]
type = FunctionIC
function = 10e6-9.81*1000*y
variable = pgas
[../]
[./zi]
type = ConstantIC
value = 0
variable = zi
[../]
# [./zi]
# type = BoundingBoxIC
# variable = zi
# x1 = 0
# x2 = 2
# y1 = 1.95
# y2 = 2
# inside = 0.1
# outside = 0
# [../]
[./porosity]
type = RandomIC
variable = porosity
min = 0.25
max = 0.275
[../]
[]
[BCs]
[./top]
type = DirichletBC
value = 0.04
variable = zi
boundary = top
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pgas
[../]
[./flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pgas
[../]
[./diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pgas
disp_long = '0 0'
disp_trans = '0 0'
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = zi
[../]
[./flux1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = zi
[../]
[./diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = zi
disp_long = '0 0'
disp_trans = '0 0'
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2sw]
type = CO2FluidProperties
[../]
[./co2]
type = TabulatedFluidProperties
fp = co2sw
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = '45'
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = 'pgas'
z = 'zi'
temperature_unit = Celsius
xnacl = 'xnacl'
capillary_pressure = pc
fluid_state = fs
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = porosity
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-11 0 0 0 1e-11 0 0 0 1e-11'
[../]
[./relperm_water]
type = PorousFlowRelativePermeabilityCorey
phase = 0
n = 2
s_res = 0.1
sum_s_res = 0.2
[../]
[./relperm_gas]
type = PorousFlowRelativePermeabilityCorey
phase = 1
n = 2
s_res = 0.1
sum_s_res = 0.2
[../]
[./diffusivity]
type = PorousFlowDiffusivityConst
diffusion_coeff = '2e-9 2e-9 2e-9 2e-9'
tortuosity = '1 1'
[../]
[]
[Preconditioning]
active = basic
[./mumps_is_best_for_parallel_jobs]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = ' lu mumps'
[../]
[./basic]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = 'gmres asm lu NONZERO 2 '
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 1e6
nl_max_its = 25
l_max_its = 100
dtmax = 1e4
nl_abs_tol = 1e-6
[./TimeStepper]
type = IterationAdaptiveDT
dt = 100
growth_factor = 2
cutback_factor = 0.5
[../]
[]
[Functions]
[./flux]
type = ParsedFunction
vals = 'delta_xco2 dt'
vars = 'dx dt'
value = 'dx/dt'
[../]
[]
[Postprocessors]
[./total_co2_in_gas]
type = PorousFlowFluidMass
phase = 1
fluid_component = 1
[../]
[./total_co2_in_liquid]
type = PorousFlowFluidMass
phase = 0
fluid_component = 1
[../]
[./numdofs]
type = NumDOFs
[../]
[./delta_xco2]
type = ChangeOverTimePostprocessor
postprocessor = total_co2_in_liquid
[../]
[./dt]
type = TimestepSize
[../]
[./flux]
type = FunctionValuePostprocessor
function = flux
[../]
[]
[Outputs]
print_linear_residuals = false
perf_graph = true
exodus = true
csv = true
[]
modules/porous_flow/test/tests/fluidstate/brineco2.i
# Tests correct calculation of properties in PorousFlowBrineCO2
[Mesh]
[./mesh]
type = GeneratedMeshGenerator
dim = 2
[../]
[]
[GlobalParams]
PorousFlowDictator = dictator
temperature = 30
[]
[Variables]
[./pgas]
initial_condition = 20e6
[../]
[./z]
initial_condition = 0.2
[../]
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.1
[../]
[./pressure_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./pressure_water]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_water]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./enthalpy_water]
order = CONSTANT
family = MONOMIAL
[../]
[./enthalpy_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./internal_energy_water]
order = CONSTANT
family = MONOMIAL
[../]
[./internal_energy_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_gas]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./pressure_water]
type = PorousFlowPropertyAux
variable = pressure_water
property = pressure
phase = 0
execute_on = timestep_end
[../]
[./pressure_gas]
type = PorousFlowPropertyAux
variable = pressure_gas
property = pressure
phase = 1
execute_on = timestep_end
[../]
[./saturation_water]
type = PorousFlowPropertyAux
variable = saturation_water
property = saturation
phase = 0
execute_on = timestep_end
[../]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = timestep_end
[../]
[./density_water]
type = PorousFlowPropertyAux
variable = density_water
property = density
phase = 0
execute_on = timestep_end
[../]
[./density_gas]
type = PorousFlowPropertyAux
variable = density_gas
property = density
phase = 1
execute_on = timestep_end
[../]
[./viscosity_water]
type = PorousFlowPropertyAux
variable = viscosity_water
property = viscosity
phase = 0
execute_on = timestep_end
[../]
[./viscosity_gas]
type = PorousFlowPropertyAux
variable = viscosity_gas
property = viscosity
phase = 1
execute_on = timestep_end
[../]
[./enthalpy_water]
type = PorousFlowPropertyAux
variable = enthalpy_water
property = enthalpy
phase = 0
execute_on = timestep_end
[../]
[./enthalpy_gas]
type = PorousFlowPropertyAux
variable = enthalpy_gas
property = enthalpy
phase = 1
execute_on = timestep_end
[../]
[./internal_energy_water]
type = PorousFlowPropertyAux
variable = internal_energy_water
property = internal_energy
phase = 0
execute_on = timestep_end
[../]
[./internal_energy_gas]
type = PorousFlowPropertyAux
variable = internal_energy_gas
property = internal_energy
phase = 1
execute_on = timestep_end
[../]
[./x1_water]
type = PorousFlowPropertyAux
variable = x1_water
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = timestep_end
[../]
[./x1_gas]
type = PorousFlowPropertyAux
variable = x1_gas
property = mass_fraction
phase = 1
fluid_component = 1
execute_on = timestep_end
[../]
[./x0_water]
type = PorousFlowPropertyAux
variable = x0_water
property = mass_fraction
phase = 0
fluid_component = 0
execute_on = timestep_end
[../]
[./x0_gas]
type = PorousFlowPropertyAux
variable = x0_gas
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = timestep_end
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = z
fluid_component = 1
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas z'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = z
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
[./density_water]
type = ElementIntegralVariablePostprocessor
variable = density_water
[../]
[./density_gas]
type = ElementIntegralVariablePostprocessor
variable = density_gas
[../]
[./viscosity_water]
type = ElementIntegralVariablePostprocessor
variable = viscosity_water
[../]
[./viscosity_gas]
type = ElementIntegralVariablePostprocessor
variable = viscosity_gas
[../]
[./enthalpy_water]
type = ElementIntegralVariablePostprocessor
variable = enthalpy_water
[../]
[./enthalpy_gas]
type = ElementIntegralVariablePostprocessor
variable = enthalpy_gas
[../]
[./internal_energy_water]
type = ElementIntegralVariablePostprocessor
variable = internal_energy_water
[../]
[./internal_energy_gas]
type = ElementIntegralVariablePostprocessor
variable = internal_energy_gas
[../]
[./x1_water]
type = ElementIntegralVariablePostprocessor
variable = x1_water
[../]
[./x0_water]
type = ElementIntegralVariablePostprocessor
variable = x0_water
[../]
[./x1_gas]
type = ElementIntegralVariablePostprocessor
variable = x1_gas
[../]
[./x0_gas]
type = ElementIntegralVariablePostprocessor
variable = x0_gas
[../]
[./sg]
type = ElementIntegralVariablePostprocessor
variable = saturation_gas
[../]
[./sw]
type = ElementIntegralVariablePostprocessor
variable = saturation_water
[../]
[./pwater]
type = ElementIntegralVariablePostprocessor
variable = pressure_water
[../]
[./pgas]
type = ElementIntegralVariablePostprocessor
variable = pressure_gas
[../]
[./x0mass]
type = PorousFlowFluidMass
fluid_component = 0
phase = '0 1'
[../]
[./x1mass]
type = PorousFlowFluidMass
fluid_component = 1
phase = '0 1'
[../]
[]
[Outputs]
csv = true
file_base = brineco2
execute_on = 'TIMESTEP_END'
perf_graph = false
[]
modules/porous_flow/test/tests/fluidstate/brineco2_hightemp.i
# Tests correct calculation of properties in PorousFlowBrineCO2 in the elevated
# temperature regime (T > 110C)
[Mesh]
type = GeneratedMesh
dim = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
temperature = 250
[]
[Variables]
[./pgas]
initial_condition = 20e6
[../]
[./z]
initial_condition = 0.2
[../]
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.1
[../]
[./pressure_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./pressure_water]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./saturation_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_water]
order = CONSTANT
family = MONOMIAL
[../]
[./density_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_water]
order = CONSTANT
family = MONOMIAL
[../]
[./viscosity_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x0_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_water]
order = CONSTANT
family = MONOMIAL
[../]
[./x1_gas]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./pressure_water]
type = PorousFlowPropertyAux
variable = pressure_water
property = pressure
phase = 0
execute_on = timestep_end
[../]
[./pressure_gas]
type = PorousFlowPropertyAux
variable = pressure_gas
property = pressure
phase = 1
execute_on = timestep_end
[../]
[./saturation_water]
type = PorousFlowPropertyAux
variable = saturation_water
property = saturation
phase = 0
execute_on = timestep_end
[../]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = timestep_end
[../]
[./density_water]
type = PorousFlowPropertyAux
variable = density_water
property = density
phase = 0
execute_on = timestep_end
[../]
[./density_gas]
type = PorousFlowPropertyAux
variable = density_gas
property = density
phase = 1
execute_on = timestep_end
[../]
[./viscosity_water]
type = PorousFlowPropertyAux
variable = viscosity_water
property = viscosity
phase = 0
execute_on = timestep_end
[../]
[./viscosity_gas]
type = PorousFlowPropertyAux
variable = viscosity_gas
property = viscosity
phase = 1
execute_on = timestep_end
[../]
[./x1_water]
type = PorousFlowPropertyAux
variable = x1_water
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = timestep_end
[../]
[./x1_gas]
type = PorousFlowPropertyAux
variable = x1_gas
property = mass_fraction
phase = 1
fluid_component = 1
execute_on = timestep_end
[../]
[./x0_water]
type = PorousFlowPropertyAux
variable = x0_water
property = mass_fraction
phase = 0
fluid_component = 0
execute_on = timestep_end
[../]
[./x0_gas]
type = PorousFlowPropertyAux
variable = x0_gas
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = timestep_end
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = z
fluid_component = 1
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas z'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = z
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
[./density_water]
type = ElementIntegralVariablePostprocessor
variable = density_water
[../]
[./density_gas]
type = ElementIntegralVariablePostprocessor
variable = density_gas
[../]
[./viscosity_water]
type = ElementIntegralVariablePostprocessor
variable = viscosity_water
[../]
[./viscosity_gas]
type = ElementIntegralVariablePostprocessor
variable = viscosity_gas
[../]
[./x1_water]
type = ElementIntegralVariablePostprocessor
variable = x1_water
[../]
[./x0_water]
type = ElementIntegralVariablePostprocessor
variable = x0_water
[../]
[./x1_gas]
type = ElementIntegralVariablePostprocessor
variable = x1_gas
[../]
[./x0_gas]
type = ElementIntegralVariablePostprocessor
variable = x0_gas
[../]
[./sg]
type = ElementIntegralVariablePostprocessor
variable = saturation_gas
[../]
[./sw]
type = ElementIntegralVariablePostprocessor
variable = saturation_water
[../]
[./pwater]
type = ElementIntegralVariablePostprocessor
variable = pressure_water
[../]
[./pgas]
type = ElementIntegralVariablePostprocessor
variable = pressure_gas
[../]
[./x0mass]
type = PorousFlowFluidMass
fluid_component = 0
phase = '0 1'
[../]
[./x1mass]
type = PorousFlowFluidMass
fluid_component = 1
phase = '0 1'
[../]
[]
[Outputs]
csv = true
execute_on = 'TIMESTEP_END'
perf_graph = false
[]
modules/porous_flow/test/tests/jacobian/brineco2_liquid.i
# Tests correct calculation of properties derivatives in PorousFlowFluidState
# for conditions that give a single liquid phase
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[AuxVariables]
[./xnacl]
initial_condition = 0.05
[../]
[]
[Variables]
[./pgas]
[../]
[./zi]
[../]
[]
[ICs]
[./pgas]
type = RandomIC
min = 5e6
max = 8e6
variable = pgas
[../]
[./z_liquid]
type = RandomIC
min = 0.01
max = 0.03
variable = zi
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
variable = pgas
fluid_component = 0
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
variable = zi
fluid_component = 1
[../]
[./adv0]
type = PorousFlowAdvectiveFlux
variable = pgas
fluid_component = 0
[../]
[./adv1]
type = PorousFlowAdvectiveFlux
variable = zi
fluid_component = 1
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi'
number_fluid_phases = 2
number_fluid_components = 2
[../]
[./pc]
type = PorousFlowCapillaryPressureVG
m = 0.5
alpha = 1
pc_max = 1e4
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[./water]
type = Water97FluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = 50
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-12 0 0 0 1e-12 0 0 0 1e-12'
[../]
[./relperm0]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 0
[../]
[./relperm1]
type = PorousFlowRelativePermeabilityCorey
n = 3
phase = 1
[../]
[./porosity]
type = PorousFlowPorosityConst
porosity = 0.1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 1
end_time = 1
nl_abs_tol = 1e-12
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[AuxVariables]
[./sgas]
family = MONOMIAL
order = CONSTANT
[../]
[]
[AuxKernels]
[./sgas]
type = PorousFlowPropertyAux
property = saturation
phase = 1
variable = sgas
[../]
[]
[Postprocessors]
[./sgas_min]
type = ElementExtremeValue
variable = sgas
value_type = min
[../]
[./sgas_max]
type = ElementExtremeValue
variable = sgas
value_type = max
[../]
[]
modules/porous_flow/test/tests/fluidstate/theis_brineco2_nonisothermal.i
# Two phase nonisothermal Theis problem: Flow from single source.
# Constant rate injection 2 kg/s of cold CO2 into warm reservoir
# 1D cylindrical mesh
# Initially, system has only a liquid phase, until enough gas is injected
# to form a gas phase, in which case the system becomes two phase.
[Mesh]
[./mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 40
xmin = 0.1
xmax = 200
bias_x = 1.05
[../]
[]
[Problem]
type = FEProblem
coord_type = RZ
rz_coord_axis = Y
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[AuxVariables]
[./saturation_gas]
order = CONSTANT
family = MONOMIAL
[../]
[./x1]
order = CONSTANT
family = MONOMIAL
[../]
[./y0]
order = CONSTANT
family = MONOMIAL
[../]
[]
[AuxKernels]
[./saturation_gas]
type = PorousFlowPropertyAux
variable = saturation_gas
property = saturation
phase = 1
execute_on = timestep_end
[../]
[./x1]
type = PorousFlowPropertyAux
variable = x1
property = mass_fraction
phase = 0
fluid_component = 1
execute_on = timestep_end
[../]
[./y0]
type = PorousFlowPropertyAux
variable = y0
property = mass_fraction
phase = 1
fluid_component = 0
execute_on = timestep_end
[../]
[]
[Variables]
[./pgas]
initial_condition = 20e6
[../]
[./zi]
initial_condition = 0
[../]
[./xnacl]
initial_condition = 0.1
[../]
[./temperature]
initial_condition = 70
[../]
[]
[Kernels]
[./mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pgas
[../]
[./flux0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pgas
[../]
[./mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = zi
[../]
[./flux1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = zi
[../]
[./mass2]
type = PorousFlowMassTimeDerivative
fluid_component = 2
variable = xnacl
[../]
[./flux2]
type = PorousFlowAdvectiveFlux
fluid_component = 2
variable = xnacl
[../]
[./energy]
type = PorousFlowEnergyTimeDerivative
variable = temperature
[../]
[./heatadv]
type = PorousFlowHeatAdvection
variable = temperature
[../]
[./conduction]
type = PorousFlowHeatConduction
variable = temperature
[../]
[]
[UserObjects]
[./dictator]
type = PorousFlowDictator
porous_flow_vars = 'pgas zi xnacl temperature'
number_fluid_phases = 2
number_fluid_components = 3
[../]
[./pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[../]
[./fs]
type = PorousFlowBrineCO2
brine_fp = brine
co2_fp = co2
capillary_pressure = pc
[../]
[]
[Modules]
[./FluidProperties]
[./co2]
type = CO2FluidProperties
[../]
[./brine]
type = BrineFluidProperties
[../]
[../]
[]
[Materials]
[./temperature]
type = PorousFlowTemperature
temperature = temperature
[../]
[./brineco2]
type = PorousFlowFluidState
gas_porepressure = pgas
z = zi
temperature = temperature
temperature_unit = Celsius
xnacl = xnacl
capillary_pressure = pc
fluid_state = fs
[../]
[./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 = 2
phase = 0
s_res = 0.1
sum_s_res = 0.1
[../]
[./relperm_gas]
type = PorousFlowRelativePermeabilityCorey
n = 2
phase = 1
[../]
[./rockheat]
type = PorousFlowMatrixInternalEnergy
specific_heat_capacity = 1000
density = 2500
[../]
[./rock_thermal_conductivity]
type = PorousFlowThermalConductivityIdeal
dry_thermal_conductivity = '50 0 0 0 50 0 0 0 50'
[../]
[]
[BCs]
[./cold_gas]
type = DirichletBC
boundary = left
variable = temperature
value = 20
[../]
[./gas_injecton]
type = PorousFlowSink
boundary = left
variable = zi
flux_function = -0.159155
[../]
[./rightwater]
type = DirichletBC
boundary = right
value = 20e6
variable = pgas
[../]
[./righttemp]
type = DirichletBC
boundary = right
value = 70
variable = temperature
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
petsc_options_iname = '-ksp_type -pc_type -sub_pc_type -sub_pc_factor_shift_type -pc_asm_overlap'
petsc_options_value = 'gmres asm lu NONZERO 2'
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 1e4
automatic_scaling = true
nl_abs_tol = 1e-7
nl_rel_tol = 1e-5
[./TimeStepper]
type = IterationAdaptiveDT
dt = 1
growth_factor = 1.5
[../]
[]
[Postprocessors]
[./pgas]
type = PointValue
point = '2 0 0'
variable = pgas
[../]
[./sgas]
type = PointValue
point = '2 0 0'
variable = saturation_gas
[../]
[./zi]
type = PointValue
point = '2 0 0'
variable = zi
[../]
[./temperature]
type = PointValue
point = '2 0 0'
variable = temperature
[../]
[./massgas]
type = PorousFlowFluidMass
fluid_component = 1
[../]
[./x1]
type = PointValue
point = '2 0 0'
variable = x1
[../]
[./y0]
type = PointValue
point = '2 0 0'
variable = y0
[../]
[]
[Outputs]
print_linear_residuals = false
perf_graph = true
csv = true
[]
References
- N. Spycher, K. Pruess, and J. Ennis-King.
CO$_2$-H$_2$O mixtures in the geological sequestration of CO$_2$. I. Assessment and calculation of mutual solubilities from 12 to 100C and up to 600 bar.
Geochimica et Cosmochimica Acta, 67:3015–3031, 2003.[BibTeX]
@article{spycher2003, author = "Spycher, N. and Pruess, K. and Ennis-King, J.", title = "{CO$\_2$-H$\_2$O mixtures in the geological sequestration of CO$\_2$. I. Assessment and calculation of mutual solubilities from 12 to 100C and up to 600 bar}", journal = "Geochimica et Cosmochimica Acta", volume = "67", pages = "3015--3031", year = "2003" }
- N. Spycher, K. Pruess, and J. Ennis-King.
CO$_2$-H$_2$O mixtures in the geological sequestration of CO$_2$. II. Partitioning in chloride brine at 12-100C and up to 600 bar.
Geochimica et Cosmochimica Acta, 69:3309–3320, 2005.[BibTeX]
@article{spycher2005, author = "Spycher, N. and Pruess, K. and Ennis-King, J.", title = "{CO$\_2$-H$\_2$O mixtures in the geological sequestration of CO$\_2$. II. Partitioning in chloride brine at 12-100C and up to 600 bar}", journal = "Geochimica et Cosmochimica Acta", volume = "69", pages = "3309--3320", year = "2005" }