- nodal_void_volume_uoThe name of the NodalVoidVolume UserObject.
C++ Type:UserObjectName
Controllable:No
Description:The name of the NodalVoidVolume UserObject.
- variableThe name of the variable that this object applies to
C++ Type:AuxVariableName
Controllable:No
Description:The name of the variable that this object applies to
NodalVoidVolumeAux
This AuxKernel extracts the void volume, , associated with each node from a NodalVoidVolume UserObject. Its variable
must therefore be a nodal variable (not a monomial, for instance).
is allows translation between moles of reactant at a node, , and concentration (mol/volume-of-solute) at the node:
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- 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
- check_boundary_restrictedTrueWhether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
Default:True
C++ Type:bool
Controllable:No
Description:Whether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
- execute_onLINEAR TIMESTEP_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, PRE_DISPLACE, ALWAYS.
Default:LINEAR TIMESTEP_END
C++ Type:ExecFlagEnum
Options:NONE, INITIAL, LINEAR, NONLINEAR, TIMESTEP_END, TIMESTEP_BEGIN, FINAL, CUSTOM, PRE_DISPLACE, ALWAYS
Controllable:No
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, PRE_DISPLACE, ALWAYS.
- 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.
- 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
Input Files
- (modules/geochemistry/test/tests/nodal_void_volume/nodal_void_volume.i)
- (modules/combined/examples/geochem-porous_flow/geotes_2D/aquifer_geochemistry.i)
- (modules/geochemistry/test/tests/nodal_void_volume/except.i)
- (modules/geochemistry/test/tests/nodal_void_volume/nodal_void_volume_adaptive.i)
- (modules/combined/examples/geochem-porous_flow/forge/aquifer_geochemistry.i)
- (modules/combined/examples/geochem-porous_flow/geotes_weber_tensleep/aquifer_geochemistry.i)
(modules/geochemistry/test/tests/nodal_void_volume/nodal_void_volume.i)
# Computes nodal void volume and compares with the Postprocessor hand-calculated values
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1 1 2 2'
dy = '1 4'
[]
[]
[Variables]
[u]
[]
[]
[Kernels]
[u]
type = Diffusion
variable = u
[]
[]
[Executioner]
type = Transient
end_time = 1
[]
[Outputs]
csv = true
[]
[UserObjects]
[nodal_void_volume]
type = NodalVoidVolume
porosity = porosity
concentration = u
[]
[]
[AuxVariables]
[porosity]
family = MONOMIAL
order = CONSTANT
[]
[vol]
[]
[]
[AuxKernels]
[porosity]
type = FunctionAux
variable = porosity
function = 'if(x<4, 1, 2)'
[]
[vol]
type = NodalVoidVolumeAux
variable = vol
nodal_void_volume_uo = nodal_void_volume
[]
[]
[Postprocessors]
[quarter]
type = PointValue
point = '0 0 0'
variable = vol
[]
[half]
type = PointValue
point = '1 0 0'
variable = vol
[]
[three_quarters]
type = PointValue
point = '2 0 0'
variable = vol
[]
[one_and_half]
type = PointValue
point = '4 0 0'
variable = vol
[]
[one]
type = PointValue
point = '6 0 0'
variable = vol
[]
[one_and_quarter]
type = PointValue
point = '0 1 0'
variable = vol
[]
[two_and_half]
type = PointValue
point = '1 1 0'
variable = vol
[]
[three_and_three_quarters]
type = PointValue
point = '2 1 0'
variable = vol
[]
[seven_and_half]
type = PointValue
point = '4 1 0'
variable = vol
[]
[five]
type = PointValue
point = '6 1 0'
variable = vol
[]
[]
(modules/combined/examples/geochem-porous_flow/geotes_2D/aquifer_geochemistry.i)
# Simulates geochemistry in the aquifer. This input file may be run in standalone fashion but it does not do anything of interest. To simulate something interesting, run the porous_flow.i simulation which couples to this input file using MultiApps.
# This file receives pf_rate_H2O, pf_rate_Na, pf_rate_Cl, pf_rate_SiO2 and temperature as AuxVariables from porous_flow.i.
# The pf_rate quantities are kg/s changes of fluid-component mass at each node, but the geochemistry module expects rates-of-changes of moles at every node. Secondly, since this input file considers just 1 litre of aqueous solution at every node, the nodal_void_volume is used to convert pf_rate_* into rate_*_per_1l, which is measured in mol/s/1_litre_of_aqueous_solution.
# This file sends massfrac_Na, massfrac_Cl and massfrac_SiO2 to porous_flow.i. These are computed from the corresponding transported_* quantities.
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 14 # for better resolution, use 56 or 112
ny = 8 # for better resolution, use 32 or 64
xmin = -70
xmax = 70
ymin = -40
ymax = 40
[]
[]
[GlobalParams]
point = '0 0 0'
reactor = reactor
[]
[SpatialReactionSolver]
model_definition = definition
geochemistry_reactor_name = reactor
charge_balance_species = "Cl-"
constraint_species = "H2O Na+ Cl- SiO2(aq)"
# ASSUME that 1 litre of solution contains:
constraint_value = " 1.0 0.1 0.1 0.000555052386"
constraint_meaning = "kg_solvent_water bulk_composition bulk_composition free_concentration"
constraint_unit = " kg moles moles molal"
initial_temperature = 50.0
kinetic_species_name = QuartzLike
# Per 1 litre (1000cm^3) of aqueous solution (1kg of solvent water), there is 9000cm^3 of QuartzLike, which means the initial porosity is 0.1.
kinetic_species_initial_value = 9000
kinetic_species_unit = cm3
temperature = temperature
source_species_names = 'H2O Na+ Cl- SiO2(aq)'
source_species_rates = 'rate_H2O_per_1l rate_Na_per_1l rate_Cl_per_1l rate_SiO2_per_1l'
ramp_max_ionic_strength_initial = 0 # max_ionic_strength in such a simple problem does not need ramping
add_aux_pH = false # there is no H+ in this system
evaluate_kinetic_rates_always = true # implicit time-marching used for stability
execute_console_output_on = ''
[]
[UserObjects]
[rate_quartz]
type = GeochemistryKineticRate
kinetic_species_name = QuartzLike
intrinsic_rate_constant = 1.0E-2
multiply_by_mass = true
area_quantity = 1
activation_energy = 72800.0
[]
[definition]
type = GeochemicalModelDefinition
database_file = "small_database.json"
basis_species = "H2O SiO2(aq) Na+ Cl-"
kinetic_minerals = "QuartzLike"
kinetic_rate_descriptions = "rate_quartz"
[]
[nodal_void_volume_uo]
type = NodalVoidVolume
porosity = porosity
execute_on = 'initial timestep_end' # "initial" means this is evaluated properly for the first timestep
[]
[]
[Executioner]
type = Transient
dt = 1E5
end_time = 7.76E6 # 90 days
[]
[AuxVariables]
[temperature]
initial_condition = 50.0
[]
[porosity]
initial_condition = 0.1
[]
[nodal_void_volume]
[]
[pf_rate_H2O] # change in H2O mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Na] # change in H2O mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Cl] # change in H2O mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_SiO2] # change in H2O mass (kg/s) at each node provided by the porous-flow simulation
[]
[rate_H2O_per_1l] # rate per 1 litre of aqueous solution that we consider at each node
[]
[rate_Na_per_1l]
[]
[rate_Cl_per_1l]
[]
[rate_SiO2_per_1l]
[]
[transported_H2O]
[]
[transported_Na]
[]
[transported_Cl]
[]
[transported_SiO2]
[]
[transported_mass]
[]
[massfrac_Na]
[]
[massfrac_Cl]
[]
[massfrac_SiO2]
[]
[massfrac_H2O]
[]
[]
[AuxKernels]
[porosity]
type = ParsedAux
args = free_cm3_QuartzLike
function = '1000.0 / (1000.0 + free_cm3_QuartzLike)'
variable = porosity
execute_on = 'timestep_end'
[]
[nodal_void_volume_auxk]
type = NodalVoidVolumeAux
variable = nodal_void_volume
nodal_void_volume_uo = nodal_void_volume_uo
execute_on = 'initial timestep_end' # "initial" to ensure it is properly evaluated for the first timestep
[]
[rate_H2O_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_H2O nodal_void_volume'
variable = rate_H2O_per_1l
# pf_rate = change in kg at every node
# pf_rate * 1000 / molar_mass_in_g_per_mole = change in moles at every node
# pf_rate * 1000 / molar_mass / (nodal_void_volume_in_m^3 * 1000) = change in moles per litre of aqueous solution
function = 'pf_rate_H2O / 18.0152 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Na_per_1l]
type = ParsedAux
args = 'pf_rate_Na nodal_void_volume'
variable = rate_Na_per_1l
function = 'pf_rate_Na / 22.9898 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Cl_per_1l]
type = ParsedAux
args = 'pf_rate_Cl nodal_void_volume'
variable = rate_Cl_per_1l
function = 'pf_rate_Cl / 35.453 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_SiO2_per_1l]
type = ParsedAux
args = 'pf_rate_SiO2 nodal_void_volume'
variable = rate_SiO2_per_1l
function = 'pf_rate_SiO2 / 60.0843 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[transported_H2O_auxk]
type = GeochemistryQuantityAux
variable = transported_H2O
species = H2O
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Na]
type = GeochemistryQuantityAux
variable = transported_Na
species = Na+
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Cl]
type = GeochemistryQuantityAux
variable = transported_Cl
species = Cl-
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_SiO2]
type = GeochemistryQuantityAux
variable = transported_SiO2
species = 'SiO2(aq)'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_mass_auxk]
type = ParsedAux
args = 'transported_H2O transported_Na transported_Cl transported_SiO2'
variable = transported_mass
function = 'transported_H2O * 18.0152 + transported_Na * 22.9898 + transported_Cl * 35.453 + transported_SiO2 * 60.0843'
execute_on = 'timestep_end'
[]
[massfrac_H2O]
type = ParsedAux
args = 'transported_H2O transported_mass'
variable = massfrac_H2O
function = 'transported_H2O * 18.0152 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Na]
type = ParsedAux
args = 'transported_Na transported_mass'
variable = massfrac_Na
function = 'transported_Na * 22.9898 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Cl]
type = ParsedAux
args = 'transported_Cl transported_mass'
variable = massfrac_Cl
function = 'transported_Cl * 35.453 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_SiO2]
type = ParsedAux
args = 'transported_SiO2 transported_mass'
variable = massfrac_SiO2
function = 'transported_SiO2 * 60.0843 / transported_mass'
execute_on = 'timestep_end'
[]
[]
[Postprocessors]
[cm3_quartz]
type = PointValue
variable = free_cm3_QuartzLike
[]
[porosity]
type = PointValue
variable = porosity
[]
[solution_temperature]
type = PointValue
variable = solution_temperature
[]
[massfrac_H2O]
type = PointValue
variable = massfrac_H2O
[]
[massfrac_Na]
type = PointValue
variable = massfrac_Na
[]
[massfrac_Cl]
type = PointValue
variable = massfrac_Cl
[]
[massfrac_SiO2]
type = PointValue
variable = massfrac_SiO2
[]
[]
[Outputs]
exodus = true
csv = true
[]
(modules/geochemistry/test/tests/nodal_void_volume/except.i)
# Exception test: the nodal void volume AuxVariable is a constant monomial, ooops!
[Mesh]
type = GeneratedMesh
dim = 1
[]
[Variables]
[u]
[]
[]
[Kernels]
[u]
type = Diffusion
variable = u
[]
[]
[Executioner]
type = Transient
end_time = 1
[]
[UserObjects]
[nodal_void_volume]
type = NodalVoidVolume
porosity = 1
[]
[]
[AuxVariables]
[vol]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[vol]
type = NodalVoidVolumeAux
variable = vol
nodal_void_volume_uo = nodal_void_volume
[]
[]
(modules/geochemistry/test/tests/nodal_void_volume/nodal_void_volume_adaptive.i)
# Computes nodal void volume, when using adaptivity, and compares with the Postprocessor hand-calculated values
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1 1 2 2'
dy = '1 4'
[]
[]
[Adaptivity]
initial_marker = u_marker
marker = u_marker
max_h_level = 1
[Markers]
[u_marker]
type = ValueRangeMarker
variable = u
invert = true
lower_bound = 0.02
upper_bound = 0.98
[]
[]
[]
[Variables]
[u]
[]
[]
[ICs]
[u]
type = FunctionIC
variable = u
function = 'if(x<2,0,1)'
[]
[]
[Kernels]
[dot]
type = TimeDerivative
variable = u
[]
[u]
type = Diffusion
variable = u
[]
[]
[Executioner]
type = Transient
dt = 1
end_time = 2
[]
[Outputs]
csv = true
[]
[UserObjects]
[nodal_void_volume]
type = NodalVoidVolume
porosity = porosity
[]
[]
[AuxVariables]
[porosity]
family = MONOMIAL
order = CONSTANT
[]
[vol]
[]
[]
[AuxKernels]
[porosity]
type = FunctionAux
variable = porosity
function = 'if(x<4, 1, 2)'
[]
[vol]
type = NodalVoidVolumeAux
variable = vol
nodal_void_volume_uo = nodal_void_volume
[]
[]
[Postprocessors]
[quarter]
type = PointValue
point = '0 0 0'
variable = vol
[]
[half]
type = PointValue
point = '1 0 0'
variable = vol
[]
[three_quarters]
type = PointValue
point = '2 0 0'
variable = vol
[]
[one_and_half_to_34s]
type = PointValue
point = '4 0 0'
variable = vol
[]
[one_to_14]
type = PointValue
point = '6 0 0'
variable = vol
[]
[one_and_quarter]
type = PointValue
point = '0 1 0'
variable = vol
[]
[two_and_half]
type = PointValue
point = '1 1 0'
variable = vol
[]
[three_and_three_quarters]
type = PointValue
point = '2 1 0'
variable = vol
[]
[seven_and_half_to_334]
type = PointValue
point = '4 1 0'
variable = vol
[]
[five_to_54]
type = PointValue
point = '6 1 0'
variable = vol
[]
[]
(modules/combined/examples/geochem-porous_flow/forge/aquifer_geochemistry.i)
# Simulates geochemistry in the aquifer. This input file may be run in standalone fashion, which will study the natural kinetically-controlled mineral changes in the same way as natural_reservoir.i. To simulate the FORGE injection scenario, run the porous_flow.i simulation which couples to this input file using MultiApps.
# This file receives pf_rate_H pf_rate_Na pf_rate_K pf_rate_Ca pf_rate_Mg pf_rate_SiO2 pf_rate_Al pf_rate_Cl pf_rate_SO4 pf_rate_HCO3 pf_rate_H2O and temperature as AuxVariables from porous_flow.i
# The pf_rate quantities are kg/s changes of fluid-component mass at each node, but the geochemistry module expects rates-of-changes of moles at every node. Secondly, since this input file considers just 1 litre of aqueous solution at every node, the nodal_void_volume is used to convert pf_rate_* into rate_*_per_1l, which is measured in mol/s/1_litre_of_aqueous_solution.
# This file sends massfrac_H massfrac_Na massfrac_K massfrac_Ca massfrac_Mg massfrac_SiO2 massfrac_Al massfrac_Cl massfrac_SO4 massfrac_HCO3 to porous_flow.i. These are computed from the corresponding transported_* quantities.
# The results depend on the kinetic rates used and these are recognised to be poorly constrained by experiment
[UserObjects]
[rate_Albite]
type = GeochemistryKineticRate
kinetic_species_name = Albite
intrinsic_rate_constant = 1E-17
multiply_by_mass = true
area_quantity = 10
activation_energy = 69.8E3
one_over_T0 = 0.003354
[]
[rate_Anhydrite]
type = GeochemistryKineticRate
kinetic_species_name = Anhydrite
intrinsic_rate_constant = 1.0E-7
multiply_by_mass = true
area_quantity = 10
activation_energy = 14.3E3
one_over_T0 = 0.003354
[]
[rate_Anorthite]
type = GeochemistryKineticRate
kinetic_species_name = Anorthite
intrinsic_rate_constant = 1.0E-13
multiply_by_mass = true
area_quantity = 10
activation_energy = 17.8E3
one_over_T0 = 0.003354
[]
[rate_Calcite]
type = GeochemistryKineticRate
kinetic_species_name = Calcite
intrinsic_rate_constant = 1.0E-10
multiply_by_mass = true
area_quantity = 10
activation_energy = 23.5E3
one_over_T0 = 0.003354
[]
[rate_Chalcedony]
type = GeochemistryKineticRate
kinetic_species_name = Chalcedony
intrinsic_rate_constant = 1.0E-18
multiply_by_mass = true
area_quantity = 10
activation_energy = 90.1E3
one_over_T0 = 0.003354
[]
[rate_Clinochl-7A]
type = GeochemistryKineticRate
kinetic_species_name = Clinochl-7A
intrinsic_rate_constant = 1.0E-17
multiply_by_mass = true
area_quantity = 10
activation_energy = 88.0E3
one_over_T0 = 0.003354
[]
[rate_Illite]
type = GeochemistryKineticRate
kinetic_species_name = Illite
intrinsic_rate_constant = 1E-17
multiply_by_mass = true
area_quantity = 10
activation_energy = 29E3
one_over_T0 = 0.003354
[]
[rate_K-feldspar]
type = GeochemistryKineticRate
kinetic_species_name = K-feldspar
intrinsic_rate_constant = 1E-17
multiply_by_mass = true
area_quantity = 10
activation_energy = 38E3
one_over_T0 = 0.003354
[]
[rate_Kaolinite]
type = GeochemistryKineticRate
kinetic_species_name = Kaolinite
intrinsic_rate_constant = 1E-18
multiply_by_mass = true
area_quantity = 10
activation_energy = 22.2E3
one_over_T0 = 0.003354
[]
[rate_Quartz]
type = GeochemistryKineticRate
kinetic_species_name = Quartz
intrinsic_rate_constant = 1E-18
multiply_by_mass = true
area_quantity = 10
activation_energy = 90.1E3
one_over_T0 = 0.003354
[]
[rate_Paragonite]
type = GeochemistryKineticRate
kinetic_species_name = Paragonite
intrinsic_rate_constant = 1E-17
multiply_by_mass = true
area_quantity = 10
activation_energy = 22E3
one_over_T0 = 0.003354
[]
[rate_Phlogopite]
type = GeochemistryKineticRate
kinetic_species_name = Phlogopite
intrinsic_rate_constant = 1E-17
multiply_by_mass = true
area_quantity = 10
activation_energy = 22E3
one_over_T0 = 0.003354
[]
[rate_Laumontite]
type = GeochemistryKineticRate
kinetic_species_name = Laumontite
intrinsic_rate_constant = 1.0E-15
multiply_by_mass = true
area_quantity = 10
activation_energy = 17.8E3
one_over_T0 = 0.003354
[]
[rate_Zoisite]
type = GeochemistryKineticRate
kinetic_species_name = Zoisite
intrinsic_rate_constant = 1E-16
multiply_by_mass = true
area_quantity = 10
activation_energy = 66.1E3
one_over_T0 = 0.003354
[]
[definition]
type = GeochemicalModelDefinition
database_file = '../../../../geochemistry/database/moose_geochemdb.json'
basis_species = 'H2O H+ Na+ K+ Ca++ Mg++ SiO2(aq) Al+++ Cl- SO4-- HCO3-'
remove_all_extrapolated_secondary_species = true
kinetic_minerals = 'Albite Anhydrite Anorthite Calcite Chalcedony Clinochl-7A Illite K-feldspar Kaolinite Quartz Paragonite Phlogopite Zoisite Laumontite'
kinetic_rate_descriptions = 'rate_Albite rate_Anhydrite rate_Anorthite rate_Calcite rate_Chalcedony rate_Clinochl-7A rate_Illite rate_K-feldspar rate_Kaolinite rate_Quartz rate_Paragonite rate_Phlogopite rate_Zoisite rate_Laumontite'
[]
[nodal_void_volume_uo]
type = NodalVoidVolume
porosity = porosity
execute_on = 'initial timestep_end' # "initial" means this is evaluated properly for the first timestep
[]
[]
[SpatialReactionSolver]
model_definition = definition
geochemistry_reactor_name = reactor
charge_balance_species = 'Cl-'
constraint_species = 'H2O H+ Na+ K+ Ca++ Mg++ SiO2(aq) Al+++ Cl- SO4-- HCO3-'
# Following numbers are from water_60_to_220degC_out.csv
constraint_value = ' 1.0006383866109 9.5165072498215e-07 0.100020379171 0.0059389061065 0.011570884507621 4.6626763057447e-06 0.0045110404925255 5.8096968688789e-17 0.13500708594394 6.6523540147676e-05 7.7361407898089e-05'
constraint_meaning = 'kg_solvent_water free_concentration free_concentration free_concentration free_concentration free_concentration free_concentration free_concentration bulk_composition free_concentration free_concentration'
constraint_unit = ' kg molal molal molal molal molal molal molal moles molal molal'
initial_temperature = 220
temperature = temperature
kinetic_species_name = ' Albite Anorthite K-feldspar Quartz Phlogopite Paragonite Calcite Anhydrite Chalcedony Illite Kaolinite Clinochl-7A Zoisite Laumontite'
kinetic_species_initial_value = '4.324073236492E+02 4.631370307325E+01 2.685015418378E+02 7.720095013956E+02 1.235192062541E+01 7.545461404965E-01 4.234651808835E-04 4.000485907930E-04 4.407616361072E+00 1.342524904876E+01 1.004823151125E+00 4.728132387707E-01 7.326007326007E-01 4.818116116598E-01'
kinetic_species_unit = ' moles moles moles moles moles moles moles moles moles moles moles moles moles moles'
evaluate_kinetic_rates_always = true # otherwise will easily "run out" of dissolving species
source_species_names = 'H2O H+ Na+ K+ Ca++ Mg++ SiO2(aq) Al+++ Cl- SO4-- HCO3-'
source_species_rates = 'rate_H2O_per_1l rate_H_per_1l rate_Na_per_1l rate_K_per_1l rate_Ca_per_1l rate_Mg_per_1l rate_SiO2_per_1l rate_Al_per_1l rate_Cl_per_1l rate_SO4_per_1l rate_HCO3_per_1l'
ramp_max_ionic_strength_initial = 0 # max_ionic_strength in such a simple problem does not need ramping
execute_console_output_on = ''
add_aux_molal = false # save some memory and reduce variables in output exodus
add_aux_mg_per_kg = false # save some memory and reduce variables in output exodus
add_aux_free_mg = false # save some memory and reduce variables in output exodus
add_aux_activity = false # save some memory and reduce variables in output exodus
add_aux_bulk_moles = false # save some memory and reduce variables in output exodus
adaptive_timestepping = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 15
ny = 10
xmin = -100
xmax = 200
ymin = -100
ymax = 100
[]
[injection_node]
input = gen
type = ExtraNodesetGenerator
new_boundary = injection_node
coord = '0 0 0'
[]
[]
[Executioner]
type = Transient
[TimeStepper]
type = FunctionDT
function = 'max(1E6, 0.3 * t)'
[]
end_time = 4E12
[]
[AuxVariables]
[temperature]
initial_condition = 220.0
[]
[porosity]
initial_condition = 0.01
[]
[nodal_void_volume]
[]
[free_cm3_Kfeldspar] # necessary because of the minus sign in K-feldspar which does not parse correctly in the porosity AuxKernel
[]
[free_cm3_Clinochl7A] # necessary because of the minus sign in Clinochl-7A which does not parse correctly in the porosity AuxKernel
[]
[pf_rate_H] # change in H mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Na]
[]
[pf_rate_K]
[]
[pf_rate_Ca]
[]
[pf_rate_Mg]
[]
[pf_rate_SiO2]
[]
[pf_rate_Al]
[]
[pf_rate_Cl]
[]
[pf_rate_SO4]
[]
[pf_rate_HCO3]
[]
[pf_rate_H2O] # change in H2O mass (kg/s) at each node provided by the porous-flow simulation
[]
[rate_H_per_1l]
[]
[rate_Na_per_1l]
[]
[rate_K_per_1l]
[]
[rate_Ca_per_1l]
[]
[rate_Mg_per_1l]
[]
[rate_SiO2_per_1l]
[]
[rate_Al_per_1l]
[]
[rate_Cl_per_1l]
[]
[rate_SO4_per_1l]
[]
[rate_HCO3_per_1l]
[]
[rate_H2O_per_1l]
[]
[transported_H]
[]
[transported_Na]
[]
[transported_K]
[]
[transported_Ca]
[]
[transported_Mg]
[]
[transported_SiO2]
[]
[transported_Al]
[]
[transported_Cl]
[]
[transported_SO4]
[]
[transported_HCO3]
[]
[transported_H2O]
[]
[transported_mass]
[]
[massfrac_H]
[]
[massfrac_Na]
[]
[massfrac_K]
[]
[massfrac_Ca]
[]
[massfrac_Mg]
[]
[massfrac_SiO2]
[]
[massfrac_Al]
[]
[massfrac_Cl]
[]
[massfrac_SO4]
[]
[massfrac_HCO3]
[]
[massfrac_H2O]
[]
[]
[AuxKernels]
[free_cm3_Kfeldspar]
type = GeochemistryQuantityAux
variable = free_cm3_Kfeldspar
species = 'K-feldspar'
quantity = free_cm3
execute_on = 'timestep_begin timestep_end'
[]
[free_cm3_Clinochl7A]
type = GeochemistryQuantityAux
variable = free_cm3_Clinochl7A
species = 'Clinochl-7A'
quantity = free_cm3
execute_on = 'timestep_begin timestep_end'
[]
[porosity_auxk]
type = ParsedAux
args = 'free_cm3_Albite free_cm3_Anhydrite free_cm3_Anorthite free_cm3_Calcite free_cm3_Chalcedony free_cm3_Clinochl7A free_cm3_Illite free_cm3_Kfeldspar free_cm3_Kaolinite free_cm3_Quartz free_cm3_Paragonite free_cm3_Phlogopite free_cm3_Zoisite free_cm3_Laumontite'
function = '1000.0 / (1000.0 + free_cm3_Albite + free_cm3_Anhydrite + free_cm3_Anorthite + free_cm3_Calcite + free_cm3_Chalcedony + free_cm3_Clinochl7A + free_cm3_Illite + free_cm3_Kfeldspar + free_cm3_Kaolinite + free_cm3_Quartz + free_cm3_Paragonite + free_cm3_Phlogopite + free_cm3_Zoisite + free_cm3_Laumontite)'
variable = porosity
execute_on = 'timestep_end'
[]
[nodal_void_volume_auxk]
type = NodalVoidVolumeAux
variable = nodal_void_volume
nodal_void_volume_uo = nodal_void_volume_uo
execute_on = 'initial timestep_end' # "initial" to ensure it is properly evaluated for the first timestep
[]
[rate_H_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_H nodal_void_volume'
variable = rate_H_per_1l
function = 'pf_rate_H / 1.0079 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Na_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Na nodal_void_volume'
variable = rate_Na_per_1l
function = 'pf_rate_Na / 22.9898 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_K_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_K nodal_void_volume'
variable = rate_K_per_1l
function = 'pf_rate_K / 39.0983 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Ca_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Ca nodal_void_volume'
variable = rate_Ca_per_1l
function = 'pf_rate_Ca / 40.08 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Mg_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Mg nodal_void_volume'
variable = rate_Mg_per_1l
function = 'pf_rate_Mg / 24.305 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_SiO2_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_SiO2 nodal_void_volume'
variable = rate_SiO2_per_1l
function = 'pf_rate_SiO2 / 60.0843 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Al_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Al nodal_void_volume'
variable = rate_Al_per_1l
function = 'pf_rate_Al / 26.9815 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_Cl_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Cl nodal_void_volume'
variable = rate_Cl_per_1l
function = 'pf_rate_Cl / 35.453 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_SO4_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_SO4 nodal_void_volume'
variable = rate_SO4_per_1l
function = 'pf_rate_SO4 / 96.0576 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_HCO3_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_HCO3 nodal_void_volume'
variable = rate_HCO3_per_1l
function = 'pf_rate_HCO3 / 61.0171 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[rate_H2O_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_H2O nodal_void_volume'
variable = rate_H2O_per_1l
function = 'pf_rate_H2O / 18.01801802 / nodal_void_volume'
execute_on = 'timestep_begin'
[]
[transported_H_auxk]
type = GeochemistryQuantityAux
variable = transported_H
species = 'H+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_Na_auxk]
type = GeochemistryQuantityAux
variable = transported_Na
species = 'Na+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_K_auxk]
type = GeochemistryQuantityAux
variable = transported_K
species = 'K+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_Ca_auxk]
type = GeochemistryQuantityAux
variable = transported_Ca
species = 'Ca++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_Mg_auxk]
type = GeochemistryQuantityAux
variable = transported_Mg
species = 'Mg++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_SiO2_auxk]
type = GeochemistryQuantityAux
variable = transported_SiO2
species = 'SiO2(aq)'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_Al_auxk]
type = GeochemistryQuantityAux
variable = transported_Al
species = 'Al+++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_Cl_auxk]
type = GeochemistryQuantityAux
variable = transported_Cl
species = 'Cl-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_SO4_auxk]
type = GeochemistryQuantityAux
variable = transported_SO4
species = 'SO4--'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_HCO3_auxk]
type = GeochemistryQuantityAux
variable = transported_HCO3
species = 'HCO3-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_H2O_auxk]
type = GeochemistryQuantityAux
variable = transported_H2O
species = 'H2O'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_begin'
[]
[transported_mass_auxk]
type = ParsedAux
args = ' transported_H transported_Na transported_K transported_Ca transported_Mg transported_SiO2 transported_Al transported_Cl transported_SO4 transported_HCO3 transported_H2O'
variable = transported_mass
function = 'transported_H * 1.0079 + transported_Cl * 35.453 + transported_SO4 * 96.0576 + transported_HCO3 * 61.0171 + transported_SiO2 * 60.0843 + transported_Al * 26.9815 + transported_Ca * 40.08 + transported_Mg * 24.305 + transported_K * 39.0983 + transported_Na * 22.9898 + transported_H2O * 18.01801802'
execute_on = 'timestep_end'
[]
[massfrac_H_auxk]
type = ParsedAux
args = 'transported_H transported_mass'
variable = massfrac_H
function = 'transported_H * 1.0079 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Na_auxk]
type = ParsedAux
args = 'transported_Na transported_mass'
variable = massfrac_Na
function = 'transported_Na * 22.9898 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_K_auxk]
type = ParsedAux
args = 'transported_K transported_mass'
variable = massfrac_K
function = 'transported_K * 39.0983 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Ca_auxk]
type = ParsedAux
args = 'transported_Ca transported_mass'
variable = massfrac_Ca
function = 'transported_Ca * 40.08 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Mg_auxk]
type = ParsedAux
args = 'transported_Mg transported_mass'
variable = massfrac_Mg
function = 'transported_Mg * 24.305 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_SiO2_auxk]
type = ParsedAux
args = 'transported_SiO2 transported_mass'
variable = massfrac_SiO2
function = 'transported_SiO2 * 60.0843 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Al_auxk]
type = ParsedAux
args = 'transported_Al transported_mass'
variable = massfrac_Al
function = 'transported_Al * 26.9815 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Cl_auxk]
type = ParsedAux
args = 'transported_Cl transported_mass'
variable = massfrac_Cl
function = 'transported_Cl * 35.453 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_SO4_auxk]
type = ParsedAux
args = 'transported_SO4 transported_mass'
variable = massfrac_SO4
function = 'transported_SO4 * 96.0576 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_HCO3_auxk]
type = ParsedAux
args = 'transported_HCO3 transported_mass'
variable = massfrac_HCO3
function = 'transported_HCO3 * 61.0171 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_H2O_auxk]
type = ParsedAux
args = 'transported_H2O transported_mass'
variable = massfrac_H2O
function = 'transported_H2O * 18.01801802 / transported_mass'
execute_on = 'timestep_end'
[]
[]
[GlobalParams]
point = '0 0 0'
reactor = reactor
[]
[Postprocessors]
[temperature]
type = PointValue
variable = 'solution_temperature'
[]
[porosity]
type = PointValue
variable = porosity
[]
[solution_temperature]
type = PointValue
variable = solution_temperature
[]
[massfrac_H]
type = PointValue
variable = massfrac_H
[]
[massfrac_Na]
type = PointValue
variable = massfrac_Na
[]
[massfrac_K]
type = PointValue
variable = massfrac_K
[]
[massfrac_Ca]
type = PointValue
variable = massfrac_Ca
[]
[massfrac_Mg]
type = PointValue
variable = massfrac_Mg
[]
[massfrac_SiO2]
type = PointValue
variable = massfrac_SiO2
[]
[massfrac_Al]
type = PointValue
variable = massfrac_Al
[]
[massfrac_Cl]
type = PointValue
variable = massfrac_Cl
[]
[massfrac_SO4]
type = PointValue
variable = massfrac_SO4
[]
[massfrac_HCO3]
type = PointValue
variable = massfrac_HCO3
[]
[massfrac_H2O]
type = PointValue
variable = massfrac_H2O
[]
[cm3_Albite]
type = PointValue
variable = 'free_cm3_Albite'
[]
[cm3_Anhydrite]
type = PointValue
variable = 'free_cm3_Anhydrite'
[]
[cm3_Anorthite]
type = PointValue
variable = 'free_cm3_Anorthite'
[]
[cm3_Calcite]
type = PointValue
variable = 'free_cm3_Calcite'
[]
[cm3_Chalcedony]
type = PointValue
variable = 'free_cm3_Chalcedony'
[]
[cm3_Clinochl-7A]
type = PointValue
variable = 'free_cm3_Clinochl-7A'
[]
[cm3_Illite]
type = PointValue
variable = 'free_cm3_Illite'
[]
[cm3_K-feldspar]
type = PointValue
variable = 'free_cm3_K-feldspar'
[]
[cm3_Kaolinite]
type = PointValue
variable = 'free_cm3_Kaolinite'
[]
[cm3_Quartz]
type = PointValue
variable = 'free_cm3_Quartz'
[]
[cm3_Paragonite]
type = PointValue
variable = 'free_cm3_Paragonite'
[]
[cm3_Phlogopite]
type = PointValue
variable = 'free_cm3_Phlogopite'
[]
[cm3_Zoisite]
type = PointValue
variable = 'free_cm3_Zoisite'
[]
[cm3_Laumontite]
type = PointValue
variable = 'free_cm3_Laumontite'
[]
[cm3_mineral]
type = LinearCombinationPostprocessor
pp_names = 'cm3_Albite cm3_Anhydrite cm3_Anorthite cm3_Calcite cm3_Chalcedony cm3_Clinochl-7A cm3_Illite cm3_K-feldspar cm3_Kaolinite cm3_Quartz cm3_Paragonite cm3_Phlogopite cm3_Zoisite cm3_Laumontite'
pp_coefs = '1 1 1 1 1 1 1 1 1 1 1 1 1 1'
[]
[pH]
type = PointValue
variable = 'pH'
[]
[]
[Outputs]
[exo]
type = Exodus
execute_on = final
[]
csv = true
[]
(modules/combined/examples/geochem-porous_flow/geotes_weber_tensleep/aquifer_geochemistry.i)
#########################################
# #
# File written by create_input_files.py #
# #
#########################################
# Simulates geochemistry in the aquifer. This input file may be run in standalone fashion but it does not do anything of interest. To simulate something interesting, run the porous_flow.i simulation which couples to this input file using MultiApps.
# This file receives pf_rate_H pf_rate_Cl pf_rate_SO4 pf_rate_HCO3 pf_rate_SiO2aq pf_rate_Al pf_rate_Ca pf_rate_Mg pf_rate_Fe pf_rate_K pf_rate_Na pf_rate_Sr pf_rate_F pf_rate_BOH pf_rate_Br pf_rate_Ba pf_rate_Li pf_rate_NO3 pf_rate_O2aq pf_rate_H2O and temperature as AuxVariables from porous_flow.i
# The pf_rate quantities are kg/s changes of fluid-component mass at each node, but the geochemistry module expects rates-of-changes of moles at every node. Secondly, since this input file considers just 1 litre of aqueous solution at every node, the nodal_void_volume is used to convert pf_rate_* into rate_*_per_1l, which is measured in mol/s/1_litre_of_aqueous_solution.
# This file sends massfrac_H massfrac_Cl massfrac_SO4 massfrac_HCO3 massfrac_SiO2aq massfrac_Al massfrac_Ca massfrac_Mg massfrac_Fe massfrac_K massfrac_Na massfrac_Sr massfrac_F massfrac_BOH massfrac_Br massfrac_Ba massfrac_Li massfrac_NO3 massfrac_O2aq to porous_flow.i. These are computed from the corresponding transported_* quantities.
[UserObjects]
[definition]
type = GeochemicalModelDefinition
database_file = '../../../../geochemistry/database/moose_geochemdb.json'
basis_species = 'H2O H+ Cl- SO4-- HCO3- SiO2(aq) Al+++ Ca++ Mg++ Fe++ K+ Na+ Sr++ F- B(OH)3 Br- Ba++ Li+ NO3- O2(aq)'
equilibrium_minerals = 'Siderite Pyrrhotite Dolomite Illite Anhydrite Calcite Quartz K-feldspar Kaolinite Barite Celestite Fluorite Albite Chalcedony Goethite'
[]
[nodal_void_volume_uo]
type = NodalVoidVolume
porosity = porosity
execute_on = 'initial timestep_end' # initial means this is evaluated properly for the first timestep
[]
[]
[SpatialReactionSolver]
model_definition = definition
geochemistry_reactor_name = reactor
charge_balance_species = 'Cl-'
swap_out_of_basis = 'NO3- H+ Fe++ Ba++ SiO2(aq) Mg++ O2(aq) Al+++ K+ Ca++ HCO3-'
swap_into_basis = ' NH3 Pyrrhotite K-feldspar Barite Quartz Dolomite Siderite Calcite Illite Anhydrite Kaolinite'
# ASSUME that 1 litre of solution contains:
constraint_species = 'H2O Quartz Calcite K-feldspar Siderite Dolomite Anhydrite Pyrrhotite Illite Kaolinite Barite Na+ Cl- SO4-- Li+ B(OH)3 Br- F- Sr++ NH3'
constraint_value = ' 0.99778351 322.177447 12.111108 6.8269499 6.2844304 2.8670301 1.1912027 0.51474767 0.3732507 0.20903322 0.0001865889 1.5876606 1.5059455 0.046792579 0.013110503 0.006663119 0.001238987 0.00032108 0.000159781 0.001937302'
constraint_meaning = 'kg_solvent_water bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition bulk_composition'
constraint_unit = "kg moles moles moles moles moles moles moles moles moles moles moles moles moles moles moles moles moles moles moles"
prevent_precipitation = 'Fluorite Albite Goethite'
initial_temperature = 92
temperature = temperature
source_species_names = 'H+ Cl- SO4-- HCO3- SiO2(aq) Al+++ Ca++ Mg++ Fe++ K+ Na+ Sr++ F- B(OH)3 Br- Ba++ Li+ NO3- O2(aq) H2O'
source_species_rates = ' rate_H_per_1l rate_Cl_per_1l rate_SO4_per_1l rate_HCO3_per_1l rate_SiO2aq_per_1l rate_Al_per_1l rate_Ca_per_1l rate_Mg_per_1l rate_Fe_per_1l rate_K_per_1l rate_Na_per_1l rate_Sr_per_1l rate_F_per_1l rate_BOH_per_1l rate_Br_per_1l rate_Ba_per_1l rate_Li_per_1l rate_NO3_per_1l rate_O2aq_per_1l rate_H2O_per_1l'
ramp_max_ionic_strength_initial = 0 # max_ionic_strength in such a simple problem does not need ramping
execute_console_output_on = '' # only CSV and exodus output for this simulation
add_aux_molal = false # save some memory and reduce variables in output exodus
add_aux_mg_per_kg = false # save some memory and reduce variables in output exodus
add_aux_free_mg = false # save some memory and reduce variables in output exodus
add_aux_activity = false # save some memory and reduce variables in output exodus
add_aux_bulk_moles = false # save some memory and reduce variables in output exodus
adaptive_timestepping = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 3
xmin = -75
xmax = 75
ymin = 0
ymax = 40
zmin = -25
zmax = 25
nx = 15
ny = 4
nz = 5
[]
[aquifer]
type = ParsedSubdomainMeshGenerator
input = gen
block_id = 1
block_name = aquifer
combinatorial_geometry = 'z >= -5 & z <= 5'
[]
[injection_nodes]
input = aquifer
type = ExtraNodesetGenerator
new_boundary = injection_nodes
coord = '-25 0 -5; -25 0 5'
[]
[production_nodes]
input = injection_nodes
type = ExtraNodesetGenerator
new_boundary = production_nodes
coord = '25 0 -5; 25 0 5'
[]
[]
[GlobalParams]
point = '-25 0 0'
reactor = reactor
[]
[Executioner]
type = Transient
solve_type = Newton
end_time = 7.76E6 # 90 days
[TimeStepper]
type = FunctionDT
function = 'min(3E4, max(1E4, 0.2 * t))'
[]
[]
[AuxVariables]
[temperature]
initial_condition = 92.0
[]
[porosity]
initial_condition = 0.1
[]
[nodal_void_volume]
[]
[free_cm3_Kfeldspar] # necessary because of the minus sign in K-feldspar which does not parse correctly in the porosity AuxKernel
[]
[pf_rate_H] # change in H mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Cl] # change in Cl mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_SO4] # change in SO4 mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_HCO3] # change in HCO3 mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_SiO2aq] # change in SiO2aq mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Al] # change in Al mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Ca] # change in Ca mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Mg] # change in Mg mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Fe] # change in Fe mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_K] # change in K mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Na] # change in Na mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Sr] # change in Sr mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_F] # change in F mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_BOH] # change in BOH mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Br] # change in Br mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Ba] # change in Ba mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_Li] # change in Li mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_NO3] # change in NO3 mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_O2aq] # change in O2aq mass (kg/s) at each node provided by the porous-flow simulation
[]
[pf_rate_H2O] # change in H2O mass (kg/s) at each node provided by the porous-flow simulation
[]
[rate_H_per_1l]
[]
[rate_Cl_per_1l]
[]
[rate_SO4_per_1l]
[]
[rate_HCO3_per_1l]
[]
[rate_SiO2aq_per_1l]
[]
[rate_Al_per_1l]
[]
[rate_Ca_per_1l]
[]
[rate_Mg_per_1l]
[]
[rate_Fe_per_1l]
[]
[rate_K_per_1l]
[]
[rate_Na_per_1l]
[]
[rate_Sr_per_1l]
[]
[rate_F_per_1l]
[]
[rate_BOH_per_1l]
[]
[rate_Br_per_1l]
[]
[rate_Ba_per_1l]
[]
[rate_Li_per_1l]
[]
[rate_NO3_per_1l]
[]
[rate_O2aq_per_1l]
[]
[rate_H2O_per_1l]
[]
[transported_H]
[]
[transported_Cl]
[]
[transported_SO4]
[]
[transported_HCO3]
[]
[transported_SiO2aq]
[]
[transported_Al]
[]
[transported_Ca]
[]
[transported_Mg]
[]
[transported_Fe]
[]
[transported_K]
[]
[transported_Na]
[]
[transported_Sr]
[]
[transported_F]
[]
[transported_BOH]
[]
[transported_Br]
[]
[transported_Ba]
[]
[transported_Li]
[]
[transported_NO3]
[]
[transported_O2aq]
[]
[transported_H2O]
[]
[transported_mass]
[]
[massfrac_H]
[]
[massfrac_Cl]
[]
[massfrac_SO4]
[]
[massfrac_HCO3]
[]
[massfrac_SiO2aq]
[]
[massfrac_Al]
[]
[massfrac_Ca]
[]
[massfrac_Mg]
[]
[massfrac_Fe]
[]
[massfrac_K]
[]
[massfrac_Na]
[]
[massfrac_Sr]
[]
[massfrac_F]
[]
[massfrac_BOH]
[]
[massfrac_Br]
[]
[massfrac_Ba]
[]
[massfrac_Li]
[]
[massfrac_NO3]
[]
[massfrac_O2aq]
[]
[massfrac_H2O]
[]
[]
[AuxKernels]
[free_cm3_Kfeldspar]
type = GeochemistryQuantityAux
variable = free_cm3_Kfeldspar
species = 'K-feldspar'
quantity = free_cm3
execute_on = 'timestep_end'
[]
[porosity_auxk]
type = ParsedAux
args = 'free_cm3_Siderite free_cm3_Pyrrhotite free_cm3_Dolomite free_cm3_Illite free_cm3_Anhydrite free_cm3_Calcite free_cm3_Quartz free_cm3_Kfeldspar free_cm3_Kaolinite free_cm3_Barite free_cm3_Celestite free_cm3_Fluorite free_cm3_Albite free_cm3_Chalcedony free_cm3_Goethite'
function = '1000.0 / (1000.0 + free_cm3_Siderite + free_cm3_Pyrrhotite + free_cm3_Dolomite + free_cm3_Illite + free_cm3_Anhydrite + free_cm3_Calcite + free_cm3_Quartz + free_cm3_Kfeldspar + free_cm3_Kaolinite + free_cm3_Barite + free_cm3_Celestite + free_cm3_Fluorite + free_cm3_Albite + free_cm3_Chalcedony + free_cm3_Goethite)'
variable = porosity
execute_on = 'timestep_end'
[]
[nodal_void_volume_auxk]
type = NodalVoidVolumeAux
variable = nodal_void_volume
nodal_void_volume_uo = nodal_void_volume_uo
execute_on = 'initial timestep_end' # initial to ensure it is properly evaluated for the first timestep
[]
[rate_H_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_H nodal_void_volume'
variable = rate_H_per_1l
function = 'pf_rate_H / 1.0079 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Cl_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Cl nodal_void_volume'
variable = rate_Cl_per_1l
function = 'pf_rate_Cl / 35.453 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_SO4_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_SO4 nodal_void_volume'
variable = rate_SO4_per_1l
function = 'pf_rate_SO4 / 96.0576 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_HCO3_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_HCO3 nodal_void_volume'
variable = rate_HCO3_per_1l
function = 'pf_rate_HCO3 / 61.0171 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_SiO2aq_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_SiO2aq nodal_void_volume'
variable = rate_SiO2aq_per_1l
function = 'pf_rate_SiO2aq / 60.0843 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Al_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Al nodal_void_volume'
variable = rate_Al_per_1l
function = 'pf_rate_Al / 26.9815 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Ca_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Ca nodal_void_volume'
variable = rate_Ca_per_1l
function = 'pf_rate_Ca / 40.08 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Mg_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Mg nodal_void_volume'
variable = rate_Mg_per_1l
function = 'pf_rate_Mg / 24.305 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Fe_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Fe nodal_void_volume'
variable = rate_Fe_per_1l
function = 'pf_rate_Fe / 55.847 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_K_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_K nodal_void_volume'
variable = rate_K_per_1l
function = 'pf_rate_K / 39.0983 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Na_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Na nodal_void_volume'
variable = rate_Na_per_1l
function = 'pf_rate_Na / 22.9898 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Sr_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Sr nodal_void_volume'
variable = rate_Sr_per_1l
function = 'pf_rate_Sr / 87.62 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_F_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_F nodal_void_volume'
variable = rate_F_per_1l
function = 'pf_rate_F / 18.9984 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_BOH_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_BOH nodal_void_volume'
variable = rate_BOH_per_1l
function = 'pf_rate_BOH / 61.8329 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Br_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Br nodal_void_volume'
variable = rate_Br_per_1l
function = 'pf_rate_Br / 79.904 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Ba_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Ba nodal_void_volume'
variable = rate_Ba_per_1l
function = 'pf_rate_Ba / 137.33 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_Li_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_Li nodal_void_volume'
variable = rate_Li_per_1l
function = 'pf_rate_Li / 6.941 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_NO3_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_NO3 nodal_void_volume'
variable = rate_NO3_per_1l
function = 'pf_rate_NO3 / 62.0049 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_O2aq_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_O2aq nodal_void_volume'
variable = rate_O2aq_per_1l
function = 'pf_rate_O2aq / 31.9988 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[rate_H2O_per_1l_auxk]
type = ParsedAux
args = 'pf_rate_H2O nodal_void_volume'
variable = rate_H2O_per_1l
function = 'pf_rate_H2O / 18.01801802 / nodal_void_volume'
execute_on = 'timestep_end'
[]
[transported_H_auxk]
type = GeochemistryQuantityAux
variable = transported_H
species = 'H+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Cl_auxk]
type = GeochemistryQuantityAux
variable = transported_Cl
species = 'Cl-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_SO4_auxk]
type = GeochemistryQuantityAux
variable = transported_SO4
species = 'SO4--'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_HCO3_auxk]
type = GeochemistryQuantityAux
variable = transported_HCO3
species = 'HCO3-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_SiO2aq_auxk]
type = GeochemistryQuantityAux
variable = transported_SiO2aq
species = 'SiO2(aq)'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Al_auxk]
type = GeochemistryQuantityAux
variable = transported_Al
species = 'Al+++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Ca_auxk]
type = GeochemistryQuantityAux
variable = transported_Ca
species = 'Ca++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Mg_auxk]
type = GeochemistryQuantityAux
variable = transported_Mg
species = 'Mg++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Fe_auxk]
type = GeochemistryQuantityAux
variable = transported_Fe
species = 'Fe++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_K_auxk]
type = GeochemistryQuantityAux
variable = transported_K
species = 'K+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Na_auxk]
type = GeochemistryQuantityAux
variable = transported_Na
species = 'Na+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Sr_auxk]
type = GeochemistryQuantityAux
variable = transported_Sr
species = 'Sr++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_F_auxk]
type = GeochemistryQuantityAux
variable = transported_F
species = 'F-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_BOH_auxk]
type = GeochemistryQuantityAux
variable = transported_BOH
species = 'B(OH)3'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Br_auxk]
type = GeochemistryQuantityAux
variable = transported_Br
species = 'Br-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Ba_auxk]
type = GeochemistryQuantityAux
variable = transported_Ba
species = 'Ba++'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_Li_auxk]
type = GeochemistryQuantityAux
variable = transported_Li
species = 'Li+'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_NO3_auxk]
type = GeochemistryQuantityAux
variable = transported_NO3
species = 'NO3-'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_O2aq_auxk]
type = GeochemistryQuantityAux
variable = transported_O2aq
species = 'O2(aq)'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_H2O_auxk]
type = GeochemistryQuantityAux
variable = transported_H2O
species = 'H2O'
quantity = transported_moles_in_original_basis
execute_on = 'timestep_end'
[]
[transported_mass_auxk]
type = ParsedAux
args = ' transported_H transported_Cl transported_SO4 transported_HCO3 transported_SiO2aq transported_Al transported_Ca transported_Mg transported_Fe transported_K transported_Na transported_Sr transported_F transported_BOH transported_Br transported_Ba transported_Li transported_NO3 transported_O2aq transported_H2O'
variable = transported_mass
function = 'transported_H * 1.0079 + transported_Cl * 35.453 + transported_SO4 * 96.0576 + transported_HCO3 * 61.0171 + transported_SiO2aq * 60.0843 + transported_Al * 26.9815 + transported_Ca * 40.08 + transported_Mg * 24.305 + transported_Fe * 55.847 + transported_K * 39.0983 + transported_Na * 22.9898 + transported_Sr * 87.62 + transported_F * 18.9984 + transported_BOH * 61.8329 + transported_Br * 79.904 + transported_Ba * 137.33 + transported_Li * 6.941 + transported_NO3 * 62.0049 + transported_O2aq * 31.9988 + transported_H2O * 18.01801802'
execute_on = 'timestep_end'
[]
[massfrac_H_auxk]
type = ParsedAux
args = 'transported_H transported_mass'
variable = massfrac_H
function = 'transported_H * 1.0079 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Cl_auxk]
type = ParsedAux
args = 'transported_Cl transported_mass'
variable = massfrac_Cl
function = 'transported_Cl * 35.453 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_SO4_auxk]
type = ParsedAux
args = 'transported_SO4 transported_mass'
variable = massfrac_SO4
function = 'transported_SO4 * 96.0576 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_HCO3_auxk]
type = ParsedAux
args = 'transported_HCO3 transported_mass'
variable = massfrac_HCO3
function = 'transported_HCO3 * 61.0171 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_SiO2aq_auxk]
type = ParsedAux
args = 'transported_SiO2aq transported_mass'
variable = massfrac_SiO2aq
function = 'transported_SiO2aq * 60.0843 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Al_auxk]
type = ParsedAux
args = 'transported_Al transported_mass'
variable = massfrac_Al
function = 'transported_Al * 26.9815 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Ca_auxk]
type = ParsedAux
args = 'transported_Ca transported_mass'
variable = massfrac_Ca
function = 'transported_Ca * 40.08 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Mg_auxk]
type = ParsedAux
args = 'transported_Mg transported_mass'
variable = massfrac_Mg
function = 'transported_Mg * 24.305 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Fe_auxk]
type = ParsedAux
args = 'transported_Fe transported_mass'
variable = massfrac_Fe
function = 'transported_Fe * 55.847 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_K_auxk]
type = ParsedAux
args = 'transported_K transported_mass'
variable = massfrac_K
function = 'transported_K * 39.0983 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Na_auxk]
type = ParsedAux
args = 'transported_Na transported_mass'
variable = massfrac_Na
function = 'transported_Na * 22.9898 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Sr_auxk]
type = ParsedAux
args = 'transported_Sr transported_mass'
variable = massfrac_Sr
function = 'transported_Sr * 87.62 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_F_auxk]
type = ParsedAux
args = 'transported_F transported_mass'
variable = massfrac_F
function = 'transported_F * 18.9984 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_BOH_auxk]
type = ParsedAux
args = 'transported_BOH transported_mass'
variable = massfrac_BOH
function = 'transported_BOH * 61.8329 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Br_auxk]
type = ParsedAux
args = 'transported_Br transported_mass'
variable = massfrac_Br
function = 'transported_Br * 79.904 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Ba_auxk]
type = ParsedAux
args = 'transported_Ba transported_mass'
variable = massfrac_Ba
function = 'transported_Ba * 137.33 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_Li_auxk]
type = ParsedAux
args = 'transported_Li transported_mass'
variable = massfrac_Li
function = 'transported_Li * 6.941 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_NO3_auxk]
type = ParsedAux
args = 'transported_NO3 transported_mass'
variable = massfrac_NO3
function = 'transported_NO3 * 62.0049 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_O2aq_auxk]
type = ParsedAux
args = 'transported_O2aq transported_mass'
variable = massfrac_O2aq
function = 'transported_O2aq * 31.9988 / transported_mass'
execute_on = 'timestep_end'
[]
[massfrac_H2O_auxk]
type = ParsedAux
args = 'transported_H2O transported_mass'
variable = massfrac_H2O
function = 'transported_H2O * 18.01801802 / transported_mass'
execute_on = 'timestep_end'
[]
[]
[Postprocessors]
[memory]
type = MemoryUsage
outputs = 'console'
[]
[porosity]
type = PointValue
variable = porosity
[]
[solution_temperature]
type = PointValue
variable = solution_temperature
[]
[massfrac_H]
type = PointValue
variable = massfrac_H
[]
[massfrac_Cl]
type = PointValue
variable = massfrac_Cl
[]
[massfrac_SO4]
type = PointValue
variable = massfrac_SO4
[]
[massfrac_HCO3]
type = PointValue
variable = massfrac_HCO3
[]
[massfrac_SiO2aq]
type = PointValue
variable = massfrac_SiO2aq
[]
[massfrac_Al]
type = PointValue
variable = massfrac_Al
[]
[massfrac_Ca]
type = PointValue
variable = massfrac_Ca
[]
[massfrac_Mg]
type = PointValue
variable = massfrac_Mg
[]
[massfrac_Fe]
type = PointValue
variable = massfrac_Fe
[]
[massfrac_K]
type = PointValue
variable = massfrac_K
[]
[massfrac_Na]
type = PointValue
variable = massfrac_Na
[]
[massfrac_Sr]
type = PointValue
variable = massfrac_Sr
[]
[massfrac_F]
type = PointValue
variable = massfrac_F
[]
[massfrac_BOH]
type = PointValue
variable = massfrac_BOH
[]
[massfrac_Br]
type = PointValue
variable = massfrac_Br
[]
[massfrac_Ba]
type = PointValue
variable = massfrac_Ba
[]
[massfrac_Li]
type = PointValue
variable = massfrac_Li
[]
[massfrac_NO3]
type = PointValue
variable = massfrac_NO3
[]
[massfrac_O2aq]
type = PointValue
variable = massfrac_O2aq
[]
[massfrac_H2O]
type = PointValue
variable = massfrac_H2O
[]
[free_cm3_Siderite]
type = PointValue
variable = free_cm3_Siderite
[]
[free_cm3_Pyrrhotite]
type = PointValue
variable = free_cm3_Pyrrhotite
[]
[free_cm3_Dolomite]
type = PointValue
variable = free_cm3_Dolomite
[]
[free_cm3_Illite]
type = PointValue
variable = free_cm3_Illite
[]
[free_cm3_Anhydrite]
type = PointValue
variable = free_cm3_Anhydrite
[]
[free_cm3_Calcite]
type = PointValue
variable = free_cm3_Calcite
[]
[free_cm3_Quartz]
type = PointValue
variable = free_cm3_Quartz
[]
[free_cm3_K-feldspar]
type = PointValue
variable = free_cm3_K-feldspar
[]
[free_cm3_Kaolinite]
type = PointValue
variable = free_cm3_Kaolinite
[]
[free_cm3_Barite]
type = PointValue
variable = free_cm3_Barite
[]
[free_cm3_Celestite]
type = PointValue
variable = free_cm3_Celestite
[]
[free_cm3_Fluorite]
type = PointValue
variable = free_cm3_Fluorite
[]
[free_cm3_Albite]
type = PointValue
variable = free_cm3_Albite
[]
[free_cm3_Chalcedony]
type = PointValue
variable = free_cm3_Chalcedony
[]
[free_cm3_Goethite]
type = PointValue
variable = free_cm3_Goethite
[]
[]
[Outputs]
exodus = true
csv = true
[]