- PorousFlowDictatorThe UserObject that holds the list of PorousFlow variable names
C++ Type:UserObjectName
Controllable:No
Description:The UserObject that holds the list of PorousFlow variable names
- diffusion_coeffList of diffusion coefficients. Order is i) component 0 in phase 0; ii) component 1 in phase 0 ...; component 0 in phase 1; ... component k in phase n (m^2/s
C++ Type:std::vector<double>
Controllable:No
Description:List of diffusion coefficients. Order is i) component 0 in phase 0; ii) component 1 in phase 0 ...; component 0 in phase 1; ... component k in phase n (m^2/s
- tortuosityList of tortuosities. Order is i) phase 0; ii) phase 1; etc
C++ Type:std::vector<double>
Controllable:No
Description:List of tortuosities. Order is i) phase 0; ii) phase 1; etc
PorousFlowDiffusivityConst
This Material provides constant tortuosity and diffusion coefficients
Tortuosity is defined as the ratio of the shortest path to the effective path, such that .
Input Parameters
- at_nodesFalseEvaluate Material properties at nodes instead of quadpoints
Default:False
C++ Type:bool
Controllable:No
Description:Evaluate Material properties at nodes instead of quadpoints
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this boundary condition applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this boundary condition applies
- computeTrueWhen false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
Default:True
C++ Type:bool
Controllable:No
Description:When false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
Default:NONE
C++ Type:MooseEnum
Options:NONE, ELEMENT, SUBDOMAIN
Controllable:No
Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
- declare_suffixAn optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)
C++ Type:std::vector<std::string>
Controllable:No
Description:List of material properties, from this material, to output (outputs must also be defined to an output type)
- outputsnone Vector of output names were you would like to restrict the output of variables(s) associated with this object
Default:none
C++ Type:std::vector<OutputName>
Controllable:No
Description:Vector of output names were you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
Input Files
- (modules/porous_flow/test/tests/dispersion/diff01.i)
- (modules/porous_flow/test/tests/jacobian/disp01.i)
- (modules/porous_flow/test/tests/jacobian/disp02.i)
- (modules/porous_flow/examples/flow_through_fractured_media/coarse_3D.i)
- (modules/porous_flow/test/tests/chemistry/2species_equilibrium_2phase.i)
- (modules/porous_flow/examples/lava_lamp/1phase_convection.i)
- (modules/porous_flow/examples/flow_through_fractured_media/coarse.i)
- (modules/porous_flow/test/tests/jacobian/diff02.i)
- (modules/porous_flow/examples/lava_lamp/2phase_convection.i)
- (modules/porous_flow/examples/flow_through_fractured_media/fine_thick_fracture_transient.i)
- (modules/porous_flow/test/tests/dispersion/disp01.i)
- (modules/porous_flow/examples/flow_through_fractured_media/fine_transient.i)
- (modules/porous_flow/test/tests/chemistry/2species_equilibrium.i)
- (modules/porous_flow/test/tests/jacobian/disp03.i)
- (modules/porous_flow/test/tests/dispersion/diff01_action.i)
- (modules/porous_flow/test/tests/dispersion/disp01_heavy.i)
- (modules/porous_flow/test/tests/jacobian/disp04.i)
- (modules/porous_flow/test/tests/chemistry/2species_predis.i)
- (modules/porous_flow/test/tests/jacobian/diff01.i)
(modules/porous_flow/test/tests/dispersion/diff01.i)
# Test diffusive part of PorousFlowDispersiveFlux kernel by setting dispersion
# coefficients to zero. Pressure is held constant over the mesh, and gravity is
# set to zero so that no advective transport of mass takes place.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
xmax = 10
bias_x = 1.1
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[AuxVariables]
[velocity]
family = MONOMIAL
order = FIRST
[]
[]
[AuxKernels]
[velocity]
type = PorousFlowDarcyVelocityComponent
variable = velocity
component = x
[]
[]
[ICs]
[pp]
type = ConstantIC
variable = pp
value = 1e5
[]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[]
[BCs]
[left]
type = DirichletBC
value = 1
variable = massfrac0
boundary = left
[]
[right]
type = DirichletBC
value = 0
variable = massfrac0
boundary = right
[]
[pright]
type = DirichletBC
variable = pp
boundary = right
value = 1e5
[]
[pleft]
type = DirichletBC
variable = pp
boundary = left
value = 1e5
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
disp_trans = 0
disp_long = 0
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = massfrac0
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = massfrac0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
disp_trans = 0
disp_long = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 1000
viscosity = 0.001
thermal_expansion = 0
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = massfrac0
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.3
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1 1'
tortuosity = 0.1
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
[]
[]
[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
dt = 1
end_time = 20
[]
[VectorPostprocessors]
[xmass]
type = NodalValueSampler
sort_by = id
variable = massfrac0
[]
[]
[Outputs]
[out]
type = CSV
execute_on = final
[]
[]
(modules/porous_flow/test/tests/jacobian/disp01.i)
# Test the Jacobian of the dispersive contribution to the diffusive component of
# the PorousFlowDisperiveFlux kernel. By setting disp_long and disp_trans to the same
# non-zero value, and diffusion to zero (by setting tortuosity to zero), the purely
# dispersive component of the flux is zero, and the only flux is due to the contribution
# from disp_trans on the diffusive flux.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 3
xmin = 0
xmax = 1
ny = 1
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[ICs]
[pp]
type = RandomIC
variable = pp
max = 2e1
min = 1e1
[]
[massfrac0]
type = RandomIC
variable = massfrac0
min = 0
max = 1
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
gravity = '1 0 0'
disp_long = 0.1
disp_trans = 0.1
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
gravity = '1 0 0'
disp_long = 0.1
disp_trans = 0.1
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 10
thermal_expansion = 0
viscosity = 1
[]
[]
[]
[Materials]
[temp]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac0'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-2 1e-1'
tortuosity = 1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[relperm]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[]
[Preconditioning]
active = smp
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]
(modules/porous_flow/test/tests/jacobian/disp02.i)
# Test the Jacobian of the dispersive contribution to the diffusive component of
# the PorousFlowDisperiveFlux kernel along with a non-zero diffusion.
# By setting disp_long and disp_trans to the same non-zero value, the purely
# dispersive component of the flux is zero, and the only flux is due to diffusion
# and its contribution from disp_trans.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 3
xmin = 0
xmax = 1
ny = 1
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[ICs]
[pp]
type = RandomIC
variable = pp
max = 2e1
min = 1e1
[]
[massfrac0]
type = RandomIC
variable = massfrac0
min = 0
max = 1
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
gravity = '1 0 0'
disp_long = 0.1
disp_trans = 0.1
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
gravity = '1 0 0'
disp_long = 0.1
disp_trans = 0.1
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 10
thermal_expansion = 0
viscosity = 1
[]
[]
[]
[Materials]
[temp]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac0'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-2 1e-1'
tortuosity = '0.1'
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[relperm]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[]
[Preconditioning]
active = smp
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]
(modules/porous_flow/examples/flow_through_fractured_media/coarse_3D.i)
# Flow and solute transport along 2 2D eliptical fractures embedded in a 3D porous matrix
# the model domain has dimensions 1 x 1 x 0.3m and the two fracture have r1 = 0.45 and r2 = 0.2
# The fractures intersect each other and the domain boundaries on two opposite sides
# fracture aperture = 6e-4m
# fracture porosity = 6e-4m
# fracture permeability = 1.8e-11 which is based in k=3e-8 from a**2/12, and k*a = 3e-8*6e-4;
# matrix porosity = 0.1;
# matrix permeanility = 1e-20;
[Mesh]
type = FileMesh
file = coarse_3D.e
block_id = '1 2 3'
block_name = 'matrix f1 f2'
boundary_id = '1 2 3 4'
boundary_name = 'rf2 lf1 right_matrix left_matrix'
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[pp]
[]
[tracer]
[]
[]
[AuxVariables]
[velocity_x]
family = MONOMIAL
order = CONSTANT
block = 'f1 f2'
[]
[velocity_y]
family = MONOMIAL
order = CONSTANT
block = 'f1 f2'
[]
[velocity_z]
family = MONOMIAL
order = CONSTANT
block = 'f1 f2'
[]
[]
[AuxKernels]
[velocity_x]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_x
component = x
aperture = 6E-4
[]
[velocity_y]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_y
component = y
aperture = 6E-4
[]
[velocity_z]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_z
component = z
aperture = 6E-4
[]
[]
[ICs]
[pp]
type = ConstantIC
variable = pp
value = 1e6
[]
[tracer]
type = ConstantIC
variable = tracer
value = 0
[]
[]
[BCs]
[top]
type = DirichletBC
value = 0
variable = tracer
boundary = rf2
[]
[bottom]
type = DirichletBC
value = 1
variable = tracer
boundary = lf1
[]
[ptop]
type = DirichletBC
variable = pp
boundary = rf2
value = 1e6
[]
[pbottom]
type = DirichletBC
variable = pp
boundary = lf1
value = 1.02e6
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
disp_trans = 0
disp_long = 0
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = tracer
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = tracer
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = tracer
disp_trans = 0
disp_long = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp tracer'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'tracer'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro1]
type = PorousFlowPorosityConst
porosity = 6e-4 # = a * phif
block = 'f1 f2'
[]
[diff1]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1.e-9 1.e-9'
tortuosity = 1.0
block = 'f1 f2'
[]
[poro2]
type = PorousFlowPorosityConst
porosity = 0.1
block = 'matrix'
[]
[diff2]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1.e-9 1.e-9'
tortuosity = 0.1
block = 'matrix'
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[permeability1]
type = PorousFlowPermeabilityConst
permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11' # 1.8e-11 = a * kf
block = 'f1 f2'
[]
[permeability2]
type = PorousFlowPermeabilityConst
permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
block = 'matrix'
[]
[]
[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 = 20
dt = 1
[]
[VectorPostprocessors]
[xmass]
type = LineValueSampler
start_point = '-0.5 0 0'
end_point = '0.5 0 0'
sort_by = x
num_points = 41
variable = tracer
outputs = csv
[]
[]
[Outputs]
[csv]
type = CSV
execute_on = 'final'
[]
[]
(modules/porous_flow/test/tests/chemistry/2species_equilibrium_2phase.i)
# Using a two-phase system (see 2species_equilibrium for the single-phase)
# The saturations, porosity, mass fractions, tortuosity and diffusion coefficients are chosen so that the results are identical to 2species_equilibrium
#
# PorousFlow analogy of chemical_reactions/test/tests/aqueous_equilibrium/2species.i
#
# Simple equilibrium reaction example to illustrate the use of PorousFlowMassFractionAqueousEquilibriumChemistry
#
# In this example, two primary species a and b are transported by diffusion and convection
# from the left of the porous medium, reacting to form two equilibrium species pa2 and pab
# according to the equilibrium reaction:
#
# reactions = '2a = pa2 rate = 10^2
# a + b = pab rate = 10^-2'
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
[]
[Variables]
[a]
order = FIRST
family = LAGRANGE
[InitialCondition]
type = BoundingBoxIC
x1 = 0.0
y1 = 0.0
x2 = 1.0e-10
y2 = 1
inside = 1.0e-2
outside = 1.0e-10
[]
[]
[b]
order = FIRST
family = LAGRANGE
[InitialCondition]
type = BoundingBoxIC
x1 = 0.0
y1 = 0.0
x2 = 1.0e-10
y2 = 1
inside = 1.0e-2
outside = 1.0e-10
[]
[]
[]
[AuxVariables]
[eqm_k0]
initial_condition = 1E2
[]
[eqm_k1]
initial_condition = 1E-2
[]
[pressure0]
[]
[saturation1]
initial_condition = 0.25
[]
[a_in_phase0]
initial_condition = 0.0
[]
[b_in_phase0]
initial_condition = 0.0
[]
[pa2]
family = MONOMIAL
order = CONSTANT
[]
[pab]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[pa2]
type = PorousFlowPropertyAux
property = secondary_concentration
secondary_species = 0
variable = pa2
[]
[pab]
type = PorousFlowPropertyAux
property = secondary_concentration
secondary_species = 1
variable = pab
[]
[]
[ICs]
[pressure0]
type = FunctionIC
variable = pressure0
function = 2-x
[]
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Kernels]
[mass_a]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = a
[]
[flux_a]
type = PorousFlowAdvectiveFlux
variable = a
fluid_component = 0
[]
[diff_a]
type = PorousFlowDispersiveFlux
variable = a
fluid_component = 0
disp_trans = '0 0'
disp_long = '0 0'
[]
[mass_b]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = b
[]
[flux_b]
type = PorousFlowAdvectiveFlux
variable = b
fluid_component = 1
[]
[diff_b]
type = PorousFlowDispersiveFlux
variable = b
fluid_component = 1
disp_trans = '0 0'
disp_long = '0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'a b'
number_fluid_phases = 2
number_fluid_components = 3
number_aqueous_equilibrium = 2
aqueous_phase_number = 1
[]
[pc]
type = PorousFlowCapillaryPressureConst
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9 # huge, so mimic chemical_reactions
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow2PhasePS
capillary_pressure = pc
phase0_porepressure = pressure0
phase1_saturation = saturation1
[]
[massfrac]
type = PorousFlowMassFractionAqueousEquilibriumChemistry
mass_fraction_vars = 'a_in_phase0 b_in_phase0 a b'
num_reactions = 2
equilibrium_constants = 'eqm_k0 eqm_k1'
primary_activity_coefficients = '1 1'
secondary_activity_coefficients = '1 1'
reactions = '2 0
1 1'
[]
[simple_fluid0]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[simple_fluid1]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 1
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.8
[]
[permeability]
type = PorousFlowPermeabilityConst
# porous_flow permeability / porous_flow viscosity = chemical_reactions conductivity = 1E-4
permeability = '1E-7 0 0 0 1E-7 0 0 0 1E-7'
[]
[relp0]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[relp1]
type = PorousFlowRelativePermeabilityConst
phase = 1
[]
[diff]
type = PorousFlowDiffusivityConst
# porous_flow diffusion_coeff * tortuousity * porosity = chemical_reactions diffusivity = 1E-4
diffusion_coeff = '5E-4 5E-4 5E-4
5E-4 5E-4 5E-4'
tortuosity = '0.25 0.25'
[]
[]
[BCs]
[a_left]
type = DirichletBC
variable = a
boundary = left
value = 1.0e-2
[]
[b_left]
type = DirichletBC
variable = b
boundary = left
value = 1.0e-2
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 10
end_time = 100
[]
[Outputs]
print_linear_residuals = true
exodus = true
perf_graph = true
[]
(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.
#
# Note: this model is computationally expensive, so should be run with multiple cores.
[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
[]
[porosity]
type = RandomIC
variable = porosity
min = 0.25
max = 0.275
seed = 0
[]
[]
[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/examples/flow_through_fractured_media/coarse.i)
# Flow and solute transport along a fracture embedded in a porous matrix
# The fracture is represented by lower dimensional elements
# fracture aperture = 6e-4m
# fracture porosity = 6e-4m = phi * a
# fracture permeability = 1.8e-11 which is based on k=3e-8 from a**2/12, and k*a = 3e-8*6e-4
# matrix porosity = 0.1
# matrix permeanility = 1e-20
[Mesh]
type = FileMesh
file = 'coarse.e'
block_id = '1 2 3'
block_name = 'fracture matrix1 matrix2'
boundary_id = '1 2'
boundary_name = 'bottom top'
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[AuxVariables]
[velocity_x]
family = MONOMIAL
order = CONSTANT
block = 'fracture'
[]
[velocity_y]
family = MONOMIAL
order = CONSTANT
block = 'fracture'
[]
[]
[AuxKernels]
[velocity_x]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_x
component = x
aperture = 6E-4
[]
[velocity_y]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_y
component = y
aperture = 6E-4
[]
[]
[ICs]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[pp_matrix]
type = ConstantIC
variable = pp
value = 1E6
[]
[]
[BCs]
[top]
type = DirichletBC
value = 0
variable = massfrac0
boundary = top
[]
[bottom]
type = DirichletBC
value = 1
variable = massfrac0
boundary = bottom
[]
[ptop]
type = DirichletBC
variable = pp
boundary = top
value = 1e6
[]
[pbottom]
type = DirichletBC
variable = pp
boundary = bottom
value = 1.002e6
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = pp
disp_trans = 0
disp_long = 0
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = massfrac0
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = massfrac0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = massfrac0
disp_trans = 0
disp_long = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = massfrac0
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro_fracture]
type = PorousFlowPorosityConst
porosity = 6e-4 # = a * phif
block = 'fracture'
[]
[poro_matrix]
type = PorousFlowPorosityConst
porosity = 0.1
block = 'matrix1 matrix2'
[]
[diff1]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-9 1e-9'
tortuosity = 1.0
block = 'fracture'
[]
[diff2]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-9 1e-9'
tortuosity = 0.1
block = 'matrix1 matrix2'
[]
[permeability_fracture]
type = PorousFlowPermeabilityConst
permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11' # 1.8e-11 = a * kf
block = 'fracture'
[]
[permeability_matrix]
type = PorousFlowPermeabilityConst
permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
block = 'matrix1 matrix2'
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[]
[Preconditioning]
[basic]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
end_time = 10
dt = 1
# controls for nonlinear iterations
nl_max_its = 15
nl_rel_tol = 1e-14
nl_abs_tol = 1e-12
[]
[VectorPostprocessors]
[xmass]
type = LineValueSampler
start_point = '-0.5 0 0'
end_point = '0.5 0 0'
sort_by = x
num_points = 41
variable = massfrac0
outputs = csv
[]
[]
[Outputs]
[csv]
type = CSV
execute_on = 'final'
[]
[]
(modules/porous_flow/test/tests/jacobian/diff02.i)
# Test the Jacobian of the diffusive component of the PorousFlowDisperiveFlux kernel for two phases.
# By setting disp_long and disp_trans to zero, the purely diffusive component of the flux
# can be isolated. Uses constant tortuosity and diffusion coefficients
[Mesh]
type = GeneratedMesh
dim = 2
nx = 3
xmin = 0
xmax = 1
ny = 1
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[sgas]
[]
[massfrac0]
[]
[]
[AuxVariables]
[massfrac1]
[]
[]
[ICs]
[sgas]
type = RandomIC
variable = sgas
max = 1
min = 0
[]
[massfrac0]
type = RandomIC
variable = massfrac0
min = 0
max = 1
[]
[massfrac1]
type = RandomIC
variable = massfrac1
min = 0
max = 1
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = sgas
gravity = '1 0 0'
disp_long = '0 0'
disp_trans = '0 0'
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
gravity = '1 0 0'
disp_long = '0 0'
disp_trans = '0 0'
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'sgas massfrac0'
number_fluid_phases = 2
number_fluid_components = 2
[]
[pc]
type = PorousFlowCapillaryPressureConst
pc = 0
[]
[]
[Modules]
[FluidProperties]
[simple_fluid0]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 10
thermal_expansion = 0
viscosity = 1
[]
[simple_fluid1]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 1
thermal_expansion = 0
viscosity = 0.1
[]
[]
[]
[Materials]
[temp]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow2PhasePS
phase0_porepressure = 1
phase1_saturation = sgas
capillary_pressure = pc
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac0 massfrac1'
[]
[simple_fluid0]
type = PorousFlowSingleComponentFluid
fp = simple_fluid0
phase = 0
[]
[simple_fluid1]
type = PorousFlowSingleComponentFluid
fp = simple_fluid1
phase = 1
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-2 1e-1 1e-2 1e-1'
tortuosity = '0.1 0.2'
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[relperm0]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[relperm1]
type = PorousFlowRelativePermeabilityConst
phase = 1
[]
[]
[Preconditioning]
active = smp
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]
(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.
#
# Note: this model is computationally expensive, so should be run with multiple cores,
# preferably on a cluster.
[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
[]
[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'
[]
[]
[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
seed = 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/examples/flow_through_fractured_media/fine_thick_fracture_transient.i)
# Using a single-dimensional mesh
# Transient flow and solute transport along a fracture in a porous matrix
# advective dominated flow in the fracture and diffusion into the porous matrix
#
# Note that fine_thick_fracture_steady.i must be run to initialise the porepressure properly
[Mesh]
file = 'gold/fine_thick_fracture_steady_out.e'
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[pp]
initial_from_file_var = pp
initial_from_file_timestep = 1
[]
[massfrac0]
[]
[]
[AuxVariables]
[velocity_x]
family = MONOMIAL
order = CONSTANT
block = fracture
[]
[velocity_y]
family = MONOMIAL
order = CONSTANT
block = fracture
[]
[]
[AuxKernels]
[velocity_x]
type = PorousFlowDarcyVelocityComponent
variable = velocity_x
component = x
[]
[velocity_y]
type = PorousFlowDarcyVelocityComponent
variable = velocity_y
component = y
[]
[]
[ICs]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[]
[BCs]
[top]
type = DirichletBC
value = 0
variable = massfrac0
boundary = top
[]
[bottom]
type = DirichletBC
value = 1
variable = massfrac0
boundary = bottom
[]
[ptop]
type = DirichletBC
variable = pp
boundary = top
value = 1e6
[]
[pbottom]
type = DirichletBC
variable = pp
boundary = bottom
value = 1.002e6
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
disp_trans = 0
disp_long = 0
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = massfrac0
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = massfrac0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
disp_trans = 0
disp_long = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = massfrac0
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro_fracture]
type = PorousFlowPorosityConst
porosity = 1.0 # this is the true porosity of the fracture
block = 'fracture'
[]
[poro_matrix]
type = PorousFlowPorosityConst
porosity = 0.1
block = 'matrix1 matrix2'
[]
[diff1]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-9 1e-9'
tortuosity = 1.0
block = 'fracture'
[]
[diff2]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-9 1e-9'
tortuosity = 0.1
block = 'matrix1 matrix2'
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[permeability1]
type = PorousFlowPermeabilityConst
permeability = '3e-8 0 0 0 3e-8 0 0 0 3e-8' # this is the true permeability of the fracture
block = 'fracture'
[]
[permeability2]
type = PorousFlowPermeabilityConst
permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
block = 'matrix1 matrix2'
[]
[]
[Functions]
[dt_controller]
type = PiecewiseConstant
x = '0 30 40 100 200 83200'
y = '0.01 0.1 1 10 100 32'
[]
[]
[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 = 86400
#dt = 0.01
[TimeStepper]
type = FunctionDT
function = dt_controller
[]
# controls for nonlinear iterations
nl_max_its = 15
nl_rel_tol = 1e-14
nl_abs_tol = 1e-9
[]
[VectorPostprocessors]
[xmass]
type = LineValueSampler
start_point = '0.4 0 0'
end_point = '0.5 0 0'
sort_by = x
num_points = 167
variable = massfrac0
[]
[]
[Outputs]
perf_graph = true
console = true
csv = true
exodus = true
[]
(modules/porous_flow/test/tests/dispersion/disp01.i)
# Test dispersive part of PorousFlowDispersiveFlux kernel by setting diffusion
# coefficients to zero. A pressure gradient is applied over the mesh to give a
# uniform velocity. Gravity is set to zero.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 100
xmax = 10
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[AuxVariables]
[velocity]
family = MONOMIAL
order = FIRST
[]
[]
[AuxKernels]
[velocity]
type = PorousFlowDarcyVelocityComponent
variable = velocity
component = x
[]
[]
[ICs]
[pp]
type = FunctionIC
variable = pp
function = pic
[]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[]
[Functions]
[pic]
type = ParsedFunction
value = 1.1e5-x*1e3
[]
[]
[BCs]
[xleft]
type = DirichletBC
value = 1
variable = massfrac0
boundary = left
[]
[xright]
type = DirichletBC
value = 0
variable = massfrac0
boundary = right
[]
[pright]
type = DirichletBC
variable = pp
boundary = right
value = 1e5
[]
[pleft]
type = DirichletBC
variable = pp
boundary = left
value = 1.1e5
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
variable = pp
disp_trans = 0
disp_long = 0.2
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = massfrac0
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = massfrac0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
disp_trans = 0
disp_long = 0.2
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e9
density0 = 1000
viscosity = 0.001
thermal_expansion = 0
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = massfrac0
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.3
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '0 0'
tortuosity = 0.1
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
[]
[]
[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 = 1e3
dtmax = 50
[TimeStepper]
type = IterationAdaptiveDT
growth_factor = 1.5
cutback_factor = 0.5
dt = 1
[]
[]
[VectorPostprocessors]
[xmass]
type = NodalValueSampler
sort_by = id
variable = massfrac0
[]
[]
[Outputs]
[out]
type = CSV
execute_on = final
[]
[]
(modules/porous_flow/examples/flow_through_fractured_media/fine_transient.i)
# Using a mixed-dimensional mesh
# Transient flow and solute transport along a fracture in a porous matrix
# advective dominated flow in the fracture and diffusion into the porous matrix
#
# Note that fine_steady.i must be run to initialise the porepressure properly
[Mesh]
file = 'gold/fine_steady_out.e'
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Variables]
[pp]
initial_from_file_var = pp
initial_from_file_timestep = 1
[]
[massfrac0]
[]
[]
[AuxVariables]
[velocity_x]
family = MONOMIAL
order = CONSTANT
block = fracture
[]
[velocity_y]
family = MONOMIAL
order = CONSTANT
block = fracture
[]
[]
[AuxKernels]
[velocity_x]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_x
component = x
aperture = 6E-4
[]
[velocity_y]
type = PorousFlowDarcyVelocityComponentLowerDimensional
variable = velocity_y
component = y
aperture = 6E-4
[]
[]
[ICs]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[]
[BCs]
[top]
type = DirichletBC
value = 0
variable = massfrac0
boundary = top
[]
[bottom]
type = DirichletBC
value = 1
variable = massfrac0
boundary = bottom
[]
[ptop]
type = DirichletBC
variable = pp
boundary = top
value = 1e6
[]
[pbottom]
type = DirichletBC
variable = pp
boundary = bottom
value = 1.002e6
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
disp_trans = 0
disp_long = 0
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = massfrac0
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = massfrac0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
disp_trans = 0
disp_long = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = massfrac0
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro_fracture]
type = PorousFlowPorosityConst
porosity = 6e-4 # = a * phif
block = 'fracture'
[]
[poro_matrix]
type = PorousFlowPorosityConst
porosity = 0.1
block = 'matrix1 matrix2'
[]
[diff1]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-9 1e-9'
tortuosity = 1.0
block = 'fracture'
[]
[diff2]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-9 1e-9'
tortuosity = 0.1
block = 'matrix1 matrix2'
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[permeability_fracture]
type = PorousFlowPermeabilityConst
permeability = '1.8e-11 0 0 0 1.8e-11 0 0 0 1.8e-11' # kf=3e-8, a=6e-4m. 1.8e-11 = kf * a
block = 'fracture'
[]
[permeability_matrix]
type = PorousFlowPermeabilityConst
permeability = '1e-20 0 0 0 1e-20 0 0 0 1e-20'
block = 'matrix1 matrix2'
[]
[]
[Functions]
[dt_controller]
type = PiecewiseConstant
x = '0 30 40 100 200 83200'
y = '0.01 0.1 1 10 100 32'
[]
[]
[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 = 86400
[TimeStepper]
type = FunctionDT
function = dt_controller
[]
# controls for nonlinear iterations
nl_max_its = 15
nl_rel_tol = 1e-14
nl_abs_tol = 1e-9
[]
[VectorPostprocessors]
[xmass]
type = LineValueSampler
start_point = '0.4 0 0'
end_point = '0.5 0 0'
sort_by = x
num_points = 167
variable = massfrac0
[]
[]
[Outputs]
perf_graph = true
console = true
csv = true
exodus = true
[]
(modules/porous_flow/test/tests/chemistry/2species_equilibrium.i)
# PorousFlow analogy of chemical_reactions/test/tests/aqueous_equilibrium/2species.i
#
# Simple equilibrium reaction example to illustrate the use of PorousFlowMassFractionAqueousEquilibriumChemistry
#
# In this example, two primary species a and b are transported by diffusion and convection
# from the left of the porous medium, reacting to form two equilibrium species pa2 and pab
# according to the equilibrium reaction:
#
# reactions = '2a = pa2 rate = 10^2
# a + b = pab rate = 10^-2'
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
[]
[Variables]
[a]
order = FIRST
family = LAGRANGE
[InitialCondition]
type = BoundingBoxIC
x1 = 0.0
y1 = 0.0
x2 = 1.0e-10
y2 = 1
inside = 1.0e-2
outside = 1.0e-10
[]
[]
[b]
order = FIRST
family = LAGRANGE
[InitialCondition]
type = BoundingBoxIC
x1 = 0.0
y1 = 0.0
x2 = 1.0e-10
y2 = 1
inside = 1.0e-2
outside = 1.0e-10
[]
[]
[]
[AuxVariables]
[eqm_k0]
initial_condition = 1E2
[]
[eqm_k1]
initial_condition = 1E-2
[]
[pressure]
[]
[pa2]
family = MONOMIAL
order = CONSTANT
[]
[pab]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[pa2]
type = PorousFlowPropertyAux
property = secondary_concentration
secondary_species = 0
variable = pa2
[]
[pab]
type = PorousFlowPropertyAux
property = secondary_concentration
secondary_species = 1
variable = pab
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 2-x
[]
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Kernels]
[mass_a]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = a
[]
[flux_a]
type = PorousFlowFullySaturatedDarcyFlow
variable = a
fluid_component = 0
[]
[diff_a]
type = PorousFlowDispersiveFlux
variable = a
fluid_component = 0
disp_trans = 0
disp_long = 0
[]
[mass_b]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = b
[]
[flux_b]
type = PorousFlowFullySaturatedDarcyFlow
variable = b
fluid_component = 1
[]
[diff_b]
type = PorousFlowDispersiveFlux
variable = b
fluid_component = 1
disp_trans = 0
disp_long = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'a b'
number_fluid_phases = 1
number_fluid_components = 3
number_aqueous_equilibrium = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9 # huge, so mimic chemical_reactions
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pressure
[]
[massfrac]
type = PorousFlowMassFractionAqueousEquilibriumChemistry
mass_fraction_vars = 'a b'
num_reactions = 2
equilibrium_constants = 'eqm_k0 eqm_k1'
primary_activity_coefficients = '1 1'
secondary_activity_coefficients = '1 1'
reactions = '2 0
1 1'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.2
[]
[permeability]
type = PorousFlowPermeabilityConst
# porous_flow permeability / porous_flow viscosity = chemical_reactions conductivity = 1E-4
permeability = '1E-7 0 0 0 1E-7 0 0 0 1E-7'
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[diff]
type = PorousFlowDiffusivityConst
# porous_flow diffusion_coeff * tortuousity * porosity = chemical_reactions diffusivity = 1E-4
diffusion_coeff = '5E-4 5E-4 5E-4'
tortuosity = 1.0
[]
[]
[BCs]
[a_left]
type = DirichletBC
variable = a
boundary = left
value = 1.0e-2
[]
[b_left]
type = DirichletBC
variable = b
boundary = left
value = 1.0e-2
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 10
end_time = 100
[]
[Outputs]
print_linear_residuals = true
exodus = true
perf_graph = true
[]
(modules/porous_flow/test/tests/jacobian/disp03.i)
# Test the Jacobian of the dispersive contribution to the PorousFlowDisperiveFlux
# kernel by setting the diffusive component to zero (tortuosity = 0).
[Mesh]
type = GeneratedMesh
dim = 2
nx = 3
xmin = 0
xmax = 1
ny = 1
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[ICs]
[pp]
type = RandomIC
variable = pp
max = 2e1
min = 1e1
[]
[massfrac0]
type = RandomIC
variable = massfrac0
min = 0
max = 1
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
gravity = '1 0 0'
disp_long = 0.2
disp_trans = 0.1
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
gravity = '1 0 0'
disp_long = 0.2
disp_trans = 0.1
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 10
thermal_expansion = 0
viscosity = 1
[]
[]
[]
[Materials]
[temp]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac0'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-2 1e-1'
tortuosity = 1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[relperm]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[]
[Preconditioning]
active = smp
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]
(modules/porous_flow/test/tests/dispersion/diff01_action.i)
# Test diffusive part of PorousFlowDispersiveFlux kernel by setting dispersion
# coefficients to zero. Pressure is held constant over the mesh, and gravity is
# set to zero so that no advective transport of mass takes place.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
xmax = 10
bias_x = 1.1
[]
[GlobalParams]
PorousFlowDictator = andy_heheheh
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[ICs]
[pp]
type = ConstantIC
variable = pp
value = 1e5
[]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[]
[BCs]
[left]
type = DirichletBC
value = 1
variable = massfrac0
boundary = left
[]
[right]
type = DirichletBC
value = 0
variable = massfrac0
boundary = right
[]
[pright]
type = DirichletBC
variable = pp
boundary = right
value = 1e5
[]
[pleft]
type = DirichletBC
variable = pp
boundary = left
value = 1e5
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = massfrac0
disp_trans = 0
disp_long = 0
gravity = '0 0 0'
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = pp
disp_trans = 0
disp_long = 0
gravity = '0 0 0'
[]
[]
[Modules]
[FluidProperties]
[the_simple_fluid]
type = SimpleFluidProperties
thermal_expansion = 0.0
bulk_modulus = 1E7
viscosity = 0.001
density0 = 1000.0
[]
[]
[]
[PorousFlowUnsaturated]
porepressure = pp
gravity = '0 0 0'
fp = the_simple_fluid
dictator_name = andy_heheheh
relative_permeability_type = Corey
relative_permeability_exponent = 0.0
mass_fraction_vars = massfrac0
[]
[Materials]
[poro]
type = PorousFlowPorosityConst
porosity = 0.3
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1 1'
tortuosity = 0.1
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
[]
[]
[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
dt = 1
end_time = 20
[]
[VectorPostprocessors]
[xmass]
type = NodalValueSampler
sort_by = id
variable = massfrac0
[]
[]
[Outputs]
[out]
type = CSV
execute_on = final
[]
[]
(modules/porous_flow/test/tests/dispersion/disp01_heavy.i)
# Test dispersive part of PorousFlowDispersiveFlux kernel by setting diffusion
# coefficients to zero. A pressure gradient is applied over the mesh to give a
# uniform velocity. Gravity is set to zero.
# Mass fraction is set to 1 on the left hand side and 0 on the right hand side.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 200
xmax = 10
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
compute_enthalpy = false
compute_internal_energy = false
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[AuxVariables]
[velocity]
family = MONOMIAL
order = FIRST
[]
[]
[AuxKernels]
[velocity]
type = PorousFlowDarcyVelocityComponent
variable = velocity
component = x
[]
[]
[ICs]
[pp]
type = FunctionIC
variable = pp
function = pic
[]
[massfrac0]
type = ConstantIC
variable = massfrac0
value = 0
[]
[]
[Functions]
[pic]
type = ParsedFunction
value = 1.1e5-x*1e3
[]
[]
[BCs]
[xleft]
type = DirichletBC
value = 1
variable = massfrac0
boundary = left
[]
[xright]
type = DirichletBC
value = 0
variable = massfrac0
boundary = right
[]
[pright]
type = DirichletBC
variable = pp
boundary = right
value = 1e5
[]
[pleft]
type = DirichletBC
variable = pp
boundary = left
value = 1.1e5
[]
[]
[Kernels]
[mass0]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = pp
[]
[adv0]
type = PorousFlowAdvectiveFlux
fluid_component = 0
variable = pp
[]
[diff0]
type = PorousFlowDispersiveFlux
variable = pp
disp_trans = 0
disp_long = 0.2
[]
[mass1]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = massfrac0
[]
[adv1]
type = PorousFlowAdvectiveFlux
fluid_component = 1
variable = massfrac0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
disp_trans = 0
disp_long = 0.2
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e9
density0 = 1000
viscosity = 0.001
thermal_expansion = 0
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = massfrac0
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.3
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '0 0'
tortuosity = 0.1
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1e-9 0 0 0 1e-9 0 0 0 1e-9'
[]
[]
[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 = 1e3
dtmax = 10
[TimeStepper]
type = IterationAdaptiveDT
growth_factor = 1.5
cutback_factor = 0.5
dt = 1
[]
[]
[VectorPostprocessors]
[xmass]
type = NodalValueSampler
sort_by = id
variable = massfrac0
[]
[]
[Outputs]
[out]
type = CSV
execute_on = final
[]
[]
(modules/porous_flow/test/tests/jacobian/disp04.i)
# Test the Jacobian of the PorousFlowDisperiveFlux kernel
[Mesh]
type = GeneratedMesh
dim = 2
nx = 3
xmin = 0
xmax = 1
ny = 1
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[ICs]
[pp]
type = RandomIC
variable = pp
max = 2e1
min = 1e1
[]
[massfrac0]
type = RandomIC
variable = massfrac0
min = 0
max = 1
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
gravity = '1 0 0'
disp_long = 0.2
disp_trans = 0.1
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
gravity = '1 0 0'
disp_long = 0.2
disp_trans = 0.1
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 10
thermal_expansion = 0
viscosity = 1
[]
[]
[]
[Materials]
[temp]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac0'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-2 1e-1'
tortuosity = '0.1'
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[relperm]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[]
[Preconditioning]
active = smp
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]
(modules/porous_flow/test/tests/chemistry/2species_predis.i)
# PorousFlow analogy of chemical_reactions/test/tests/solid_kinetics/2species_without_action.i
#
# Simple equilibrium reaction example to illustrate the use of PorousFlowAqueousPreDisChemistry
#
# In this example, two primary species a and b diffuse towards each other from
# opposite ends of a porous medium, reacting when they meet to form a mineral
# precipitate. The kinetic reaction is
#
# a + b = mineral
#
# where a and b are the primary species (reactants), and mineral is the precipitate.
# At the time of writing, the results of this test differ from chemical_reactions because
# in PorousFlow the mineral_concentration is measured in m^3 (precipitate) / m^3 (porous_material)
# in chemical_reactions the mineral_concentration is measured in m^3 (precipitate) / m^3 (fluid)
# ie, PorousFlow_mineral_concentration = porosity * chemical_reactions_mineral_concentration
[Mesh]
type = GeneratedMesh
dim = 2
xmax = 1
ymax = 1
nx = 40
[]
[Variables]
[a]
order = FIRST
family = LAGRANGE
initial_condition = 0
[]
[b]
order = FIRST
family = LAGRANGE
initial_condition = 0
[]
[]
[AuxVariables]
[eqm_k]
initial_condition = 1E-6
[]
[pressure]
[]
[mineral]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[mineral]
type = PorousFlowPropertyAux
property = mineral_concentration
mineral_species = 0
variable = mineral
[]
[]
[GlobalParams]
PorousFlowDictator = dictator
gravity = '0 0 0'
[]
[Kernels]
[mass_a]
type = PorousFlowMassTimeDerivative
fluid_component = 0
variable = a
[]
[diff_a]
type = PorousFlowDispersiveFlux
variable = a
fluid_component = 0
disp_trans = 0
disp_long = 0
[]
[predis_a]
type = PorousFlowPreDis
variable = a
mineral_density = 1000
stoichiometry = 1
[]
[mass_b]
type = PorousFlowMassTimeDerivative
fluid_component = 1
variable = b
[]
[diff_b]
type = PorousFlowDispersiveFlux
variable = b
fluid_component = 1
disp_trans = 0
disp_long = 0
[]
[predis_b]
type = PorousFlowPreDis
variable = b
mineral_density = 1000
stoichiometry = 1
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'a b'
number_fluid_phases = 1
number_fluid_components = 3
number_aqueous_kinetic = 1
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 2e9 # huge, so mimic chemical_reactions
density0 = 1000
thermal_expansion = 0
viscosity = 1e-3
[]
[]
[]
[Materials]
[temperature]
type = PorousFlowTemperature
temperature = 298.15
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pressure
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'a b'
[]
[chem]
type = PorousFlowAqueousPreDisChemistry
primary_concentrations = 'a b'
num_reactions = 1
equilibrium_constants = eqm_k
primary_activity_coefficients = '1 1'
reactions = '1 1'
specific_reactive_surface_area = '1.0'
kinetic_rate_constant = '1.0e-8'
activation_energy = '1.5e4'
molar_volume = 1
gas_constant = 8.314
reference_temperature = 298.15
[]
[mineral_conc]
type = PorousFlowAqueousPreDisMineral
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[porosity]
type = PorousFlowPorosityConst
porosity = 0.4
[]
[permeability]
type = PorousFlowPermeabilityConst
# porous_flow permeability / porous_flow viscosity = chemical_reactions conductivity = 4E-3
permeability = '4E-6 0 0 0 4E-6 0 0 0 4E-6'
[]
[relp]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[diff]
type = PorousFlowDiffusivityConst
# porous_flow diffusion_coeff * tortuousity * porosity = chemical_reactions diffusivity = 5E-4
diffusion_coeff = '12.5E-4 12.5E-4 12.5E-4'
tortuosity = 1.0
[]
[]
[BCs]
[a_left]
type = DirichletBC
variable = a
boundary = left
value = 1.0e-2
[]
[a_right]
type = DirichletBC
variable = a
boundary = right
value = 0
[]
[b_left]
type = DirichletBC
variable = b
boundary = left
value = 0
[]
[b_right]
type = DirichletBC
variable = b
boundary = right
value = 1.0e-2
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 5
end_time = 50
[]
[Outputs]
print_linear_residuals = true
exodus = true
perf_graph = true
[]
(modules/porous_flow/test/tests/jacobian/diff01.i)
# Test the Jacobian of the diffusive component of the PorousFlowDisperiveFlux kernel.
# By setting disp_long and disp_trans to zero, the purely diffusive component of the flux
# can be isolated.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 3
xmin = 0
xmax = 1
ny = 1
ymin = 0
ymax = 1
[]
[GlobalParams]
PorousFlowDictator = dictator
[]
[Variables]
[pp]
[]
[massfrac0]
[]
[]
[ICs]
[pp]
type = RandomIC
variable = pp
max = 2e1
min = 1e1
[]
[massfrac0]
type = RandomIC
variable = massfrac0
min = 0
max = 1
[]
[]
[Kernels]
[diff0]
type = PorousFlowDispersiveFlux
fluid_component = 0
variable = pp
gravity = '1 0 0'
disp_long = 0
disp_trans = 0
[]
[diff1]
type = PorousFlowDispersiveFlux
fluid_component = 1
variable = massfrac0
gravity = '1 0 0'
disp_long = 0
disp_trans = 0
[]
[]
[UserObjects]
[dictator]
type = PorousFlowDictator
porous_flow_vars = 'pp massfrac0'
number_fluid_phases = 1
number_fluid_components = 2
[]
[]
[Modules]
[FluidProperties]
[simple_fluid]
type = SimpleFluidProperties
bulk_modulus = 1e7
density0 = 10
thermal_expansion = 0
viscosity = 1
[]
[]
[]
[Materials]
[temp]
type = PorousFlowTemperature
[]
[ppss]
type = PorousFlow1PhaseFullySaturated
porepressure = pp
[]
[massfrac]
type = PorousFlowMassFraction
mass_fraction_vars = 'massfrac0'
[]
[simple_fluid]
type = PorousFlowSingleComponentFluid
fp = simple_fluid
phase = 0
[]
[poro]
type = PorousFlowPorosityConst
porosity = 0.1
[]
[diff]
type = PorousFlowDiffusivityConst
diffusion_coeff = '1e-2 1e-1'
tortuosity = '0.1'
[]
[permeability]
type = PorousFlowPermeabilityConst
permeability = '1 0 0 0 2 0 0 0 3'
[]
[relperm]
type = PorousFlowRelativePermeabilityConst
phase = 0
[]
[]
[Preconditioning]
active = smp
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = Newton
dt = 1
end_time = 1
[]
[Outputs]
exodus = false
[]