An acidic solution
A particularly simple geochemistry
example involves computing the molality, etc, of species in an acidic solution with fixed pH.
The MOOSE input file: step 1
The database and basis species must be specified for all geochemistry
simulations. In this case, the full moose_geochemdb.json
database is used, and the basis species are simply HO, H and Cl. A piecewise linear interpolation for the equilibrium constants and the Debye-Huckel activity coefficients is used so that the result agrees exactly with the analogous GWB simulation:
[UserObjects<<<{"href": "../../../syntax/UserObjects/index.html"}>>>]
[definition]
type = GeochemicalModelDefinition<<<{"description": "User object that parses a geochemical database file, and only retains information relevant to the current geochemical model", "href": "../../../source/userobjects/GeochemicalModelDefinition.html"}>>>
database_file<<<{"description": "The name of the geochemical database file"}>>> = "../../../database/moose_geochemdb.json"
basis_species<<<{"description": "A list of basis components relevant to the aqueous-equilibrium problem. H2O must appear first in this list. These components must be chosen from the 'basis species' in the database, the sorbing sites (if any) and the decoupled redox states that are in disequilibrium (if any)."}>>> = "H2O H+ Cl-"
piecewise_linear_interpolation<<<{"description": "If true then use a piecewise-linear interpolation of logK and Debye-Huckel parameters, regardless of the interpolation type specified in the database file. This can be useful for comparing with results using other geochemistry software"}>>> = true # to reproduce the GWB result
[]
[]
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl.i)The MOOSE input file: step 2
The next piece of the input file involves specifying the initial conditions and the simulation type. This is a time-independent solve (just the equilibrium configuration is sought). The charge-balance species is chosen to be Cl, there is 1kg of solvent water and the pH is fixed to 2 (via log10activity = -2
):
[TimeIndependentReactionSolver<<<{"href": "../../../syntax/TimeIndependentReactionSolver/index.html"}>>>]
model_definition<<<{"description": "The name of the GeochemicalModelDefinition user object (you must create this UserObject yourself)"}>>> = definition
geochemistry_reactor_name<<<{"description": "The name that will be given to the GeochemistryReactor UserObject built by this action"}>>> = reactor
charge_balance_species<<<{"description": "Charge balance will be enforced on this basis species. This means that its bulk mole number may be changed from the initial value you provide in order to ensure charge neutrality. After the initial swaps have been performed, this must be in the basis, and it must be provided with a bulk_composition constraint_meaning."}>>> = "Cl-"
constraint_species<<<{"description": "Names of the species that have their values fixed to constraint_value with meaning constraint_meaning. All basis species (after swap_into_basis and swap_out_of_basis) must be provided with exactly one constraint. These constraints are used to compute the configuration during the initial problem setup, and in time-dependent simulations they may be modified as time progresses."}>>> = "H2O H+ Cl-"
constraint_value<<<{"description": "Numerical value of the containts on constraint_species"}>>> = " 1.0 -2 1E-2"
constraint_meaning<<<{"description": "Meanings of the numerical values given in constraint_value. kg_solvent_water: can only be applied to H2O and units must be kg. bulk_composition: can be applied to all non-gas species, and represents the total amount of the basis species contained as free species as well as the amount found in secondary species but not in kinetic species, and units must be moles or mass (kg, g, etc). bulk_composition_with_kinetic: can be applied to all non-gas species, and represents the total amount of the basis species contained as free species as well as the amount found in secondary species and in kinetic species, and units must be moles or mass (kg, g, etc). free_concentration: can be applied to all basis species that are not gas and not H2O and not mineral, and represents the total amount of the basis species existing freely (not as secondary species) within the solution, and units must be molal or mass_per_kg_solvent. free_mineral: can be applied to all mineral basis species, and represents the total amount of the mineral existing freely (precipitated) within the solution, and units must be moles, mass or cm3. activity and log10activity: can be applied to basis species that are not gas and not mineral and not sorbing sites, and represents the activity of the basis species (recall pH = -log10activity), and units must be dimensionless. fugacity and log10fugacity: can be applied to gases, and units must be dimensionless"}>>> = "kg_solvent_water log10activity bulk_composition"
constraint_unit<<<{"description": "Units of the numerical values given in constraint_value. Dimensionless: should only be used for activity or fugacity constraints. Moles: mole number. Molal: moles per kg solvent water. kg: kilograms. g: grams. mg: milligrams. ug: micrograms. kg_per_kg_solvent: kilograms per kg solvent water. g_per_kg_solvent: grams per kg solvent water. mg_per_kg_solvent: milligrams per kg solvent water. ug_per_kg_solvent: micrograms per kg solvent water. cm3: cubic centimeters"}>>> = " kg dimensionless moles"
ramp_max_ionic_strength_initial<<<{"description": "The number of iterations over which to progressively increase the maximum ionic strength (from zero to max_ionic_strength) during the initial equilibration. Increasing this can help in convergence of the Newton process, at the cost of spending more time finding the aqueous configuration."}>>> = 0 # max_ionic_strength in such a simple problem does not need ramping
execute_console_output_on<<<{"description": "When to execute the geochemistry console output"}>>> = initial
abs_tol<<<{"description": "If the residual of the algebraic system (measured in mol) is lower than this value, the Newton process (that finds the aqueous configuration) is deemed to have converged"}>>> = 1E-15
[]
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl.i)The MOOSE input file: optional step 3
Postprocessors may be used to produce some CSV output. In this example, all the information is already obtained in the console output, but by way of example, the following postprocessors output the activity of water, the pH, the solvent-water mass, the Cl molality, etc:
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[pH]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'pH'
[]
[solvent_mass]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'kg_solvent_H2O'
[]
[molal_Cl-]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'molal_Cl-'
[]
[mg_per_kg_HCl]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'mg_per_kg_HCl'
[]
[activity_OH-]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'activity_OH-'
[]
[bulk_H+]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'bulk_moles_H+'
[]
[temperature]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'solution_temperature'
[]
[]
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl.i)In addition:
[GlobalParams<<<{"href": "../../../syntax/GlobalParams/index.html"}>>>]
point = '0 0 0'
[]
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl.i)[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[]
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl.i)The equivalent GWB input file
The results may be compared with those generated by Geochemists Workbench. The GWB input file is
# React script that is analogous to HCl.i
data = thermo.tdat verify
conductivity = conductivity-USGS.dat
temperature = 25 C
H2O = 1 free kg
Cl- = .01 molal
balance on Cl-
pH = 2
printout species = long
epsilon = 1e-15
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl.rea)Results
The geochemistry
module produces the same result as GWB (up to precision):
Summary:
Total number of iterations required = 9
Error in calculation = 6.888e-17mol
Charge of solution = 0mol (charge-balance species = Cl-)
Mass of solvent water = 1kg
Mass of aqueous solution = 1kg
pH = 2
Ionic strength = 0.01097mol/kg(solvent water)
Stoichiometric ionic strength = 0.01097mol/kg(solvent water)
Activity of water = 0.9996
Temperature = 25
Basis Species:
H+; bulk_moles = 0.01097mol; bulk_conc = 11.05mg/kg(solution); molality = 0.01097mol/kg(solvent water); free_conc = 11.06mg/kg(solvent water); act_coeff = 0.9114; log10(a) = -2
Cl-; bulk_moles = 0.01097mol; bulk_conc = 388.8mg/kg(solution); molality = 0.01097mol/kg(solvent water); free_conc = 389mg/kg(solvent water); act_coeff = 0.8956; log10(a) = -2.008
Equilibrium Species:
HCl; molality = 7.805e-11mol/kg(solvent water); free_conc = 2.846e-06mg/kg(solvent water); act_coeff = 1; log10(a) = -10.11; HCl = 1*H+ + 1*Cl-; log10K = 6.1
OH-; molality = 1.149e-12mol/kg(solvent water); free_conc = 1.954e-08mg/kg(solvent water); act_coeff = 0.8971; log10(a) = -11.99; OH- = 1*H2O - 1*H+; log10K = 13.99
Without an Action
Although most users will want to use the Action system, it is possible to build all geochemistry
models without using Actions for slightly more fine-grained control. The equivalent input file is:
# This is an example of an input file that does not utilize an action. Its functionality is the same as HCl.i
# This solves for molalities in a system just containing HCl
[GlobalParams<<<{"href": "../../../syntax/GlobalParams/index.html"}>>>]
point = '0 0 0'
[]
[Mesh<<<{"href": "../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 1
nx= 1
[]
[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
[u]
[]
[]
[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
[u]
type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../../source/kernels/Diffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
[]
[]
[AuxVariables<<<{"href": "../../../syntax/AuxVariables/index.html"}>>>]
[solution_temperature]
[]
[kg_solvent_H2O]
[]
[activity_H2O]
[]
[bulk_moles_H2O]
[]
[pH]
[]
[molal_H+]
[]
[molal_Cl-]
[]
[molal_HCl]
[]
[molal_OH-]
[]
[mg_per_kg_H+]
[]
[mg_per_kg_Cl-]
[]
[mg_per_kg_HCl]
[]
[mg_per_kg_OH-]
[]
[activity_H+]
[]
[activity_Cl-]
[]
[activity_HCl]
[]
[activity_OH-]
[]
[bulk_moles_H+]
[]
[bulk_moles_Cl-]
[]
[bulk_moles_HCl]
[]
[bulk_moles_OH-]
[]
[]
[AuxKernels<<<{"href": "../../../syntax/AuxKernels/index.html"}>>>]
[solution_temperature]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H+'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = solution_temperature
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = temperature
[]
[kg_solvent_H2O]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H2O'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = kg_solvent_H2O
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = molal
[]
[activity_H2O]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H2O'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = activity_H2O
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = activity
[]
[bulk_moles_H2O]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H2O'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = bulk_moles_H2O
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = bulk_moles
[]
[pH]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H+'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = pH
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = neglog10a
[]
[molal_H+]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H+'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'molal_H+'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = molal
[]
[molal_Cl-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'Cl-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'molal_Cl-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = molal
[]
[molal_HCl]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'HCl'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'molal_HCl'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = molal
[]
[molal_OH-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'OH-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'molal_OH-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = molal
[]
[mg_per_kg_H+]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H+'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'mg_per_kg_H+'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = mg_per_kg
[]
[mg_per_kg_Cl-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'Cl-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'mg_per_kg_Cl-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = mg_per_kg
[]
[mg_per_kg_HCl]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'HCl'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'mg_per_kg_HCl'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = mg_per_kg
[]
[mg_per_kg_OH-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'OH-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'mg_per_kg_OH-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = mg_per_kg
[]
[activity_H+]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H+'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'activity_H+'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = activity
[]
[activity_Cl-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'Cl-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'activity_Cl-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = activity
[]
[activity_HCl]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'HCl'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'activity_HCl'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = activity
[]
[activity_OH-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'OH-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'activity_OH-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = activity
[]
[bulk_moles_H+]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'H+'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'bulk_moles_H+'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = bulk_moles
[]
[bulk_moles_Cl-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'Cl-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'bulk_moles_Cl-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = bulk_moles
[]
[bulk_moles_HCl]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'HCl'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'bulk_moles_HCl'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = bulk_moles
[]
[bulk_moles_OH-]
type = GeochemistryQuantityAux<<<{"description": "Extracts information from the Reactor and records it in the AuxVariable", "href": "../../../source/auxkernels/GeochemistryQuantityAux.html"}>>>
species<<<{"description": "Species name"}>>> = 'OH-'
reactor<<<{"description": "The name of the Geochemistry*Reactor user object."}>>> = reactor
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'bulk_moles_OH-'
quantity<<<{"description": "Type of quantity to output. These are available for non-kinetic species: activity (which equals fugacity for gases); bulk moles (this will be zero if the species is not in the basis); neglog10a = -log10(activity); transported_moles_in_original_basis (will throw an error if species is not in original basis). These are available only for non-kinetic non-minerals: molal = mol(species)/kg(solvent_water); mg_per_kg = mg(species)/kg(solvent_water). These are available only for minerals: free_mg = free mg; free_cm3 = free cubic-centimeters; moles_dumped = moles dumped when special dump and flow-through modes are active. These are available for minerals that host sorbing sites: surface_charge (C/m^2); surface_potential (V). These are available for kinetic species: kinetic_moles; kinetic_additions (-dt * rate = mole increment in kinetic species for this timestep). If quantity=temperature, then the 'species' is ignored and the AuxVariable will record the aqueous solution temperature in degC"}>>> = bulk_moles
[]
[]
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[pH]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'pH'
[]
[solvent_mass]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'kg_solvent_H2O'
[]
[molal_Cl-]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'molal_Cl-'
[]
[mg_per_kg_HCl]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'mg_per_kg_HCl'
[]
[activity_OH-]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'activity_OH-'
[]
[bulk_H+]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'bulk_moles_H+'
[]
[temperature]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = 'solution_temperature'
[]
[]
[Executioner<<<{"href": "../../../syntax/Executioner/index.html"}>>>]
type = Steady
[]
[UserObjects<<<{"href": "../../../syntax/UserObjects/index.html"}>>>]
[definition]
type = GeochemicalModelDefinition<<<{"description": "User object that parses a geochemical database file, and only retains information relevant to the current geochemical model", "href": "../../../source/userobjects/GeochemicalModelDefinition.html"}>>>
database_file<<<{"description": "The name of the geochemical database file"}>>> = "../../../database/moose_geochemdb.json"
basis_species<<<{"description": "A list of basis components relevant to the aqueous-equilibrium problem. H2O must appear first in this list. These components must be chosen from the 'basis species' in the database, the sorbing sites (if any) and the decoupled redox states that are in disequilibrium (if any)."}>>> = "H2O H+ Cl-"
piecewise_linear_interpolation<<<{"description": "If true then use a piecewise-linear interpolation of logK and Debye-Huckel parameters, regardless of the interpolation type specified in the database file. This can be useful for comparing with results using other geochemistry software"}>>> = true # to reproduce the GWB result
[]
[reactor]
type = GeochemistryTimeDependentReactor<<<{"description": "UserObject that controls the time-dependent geochemistry reaction processes. Spatial dependence is not possible using this class", "href": "../../../source/userobjects/GeochemistryTimeDependentReactor.html"}>>>
model_definition<<<{"description": "The name of the GeochemicalModelDefinition user object."}>>> = definition
charge_balance_species<<<{"description": "Charge balance will be enforced on this basis species. This means that its bulk mole number may be changed from the initial value you provide in order to ensure charge neutrality. After the initial swaps have been performed, this must be in the basis, and it must be provided with a bulk_composition constraint_meaning."}>>> = "Cl-"
constraint_species<<<{"description": "Names of the species that have their values fixed to constraint_value with meaning constraint_meaning. All basis species (after swap_into_basis and swap_out_of_basis) must be provided with exactly one constraint. These constraints are used to compute the configuration during the initial problem setup, and in time-dependent simulations they may be modified as time progresses."}>>> = "H2O H+ Cl-"
constraint_value<<<{"description": "Numerical value of the containts on constraint_species"}>>> = " 1.0 -2 1E-2"
constraint_meaning<<<{"description": "Meanings of the numerical values given in constraint_value. kg_solvent_water: can only be applied to H2O and units must be kg. bulk_composition: can be applied to all non-gas species, and represents the total amount of the basis species contained as free species as well as the amount found in secondary species but not in kinetic species, and units must be moles or mass (kg, g, etc). bulk_composition_with_kinetic: can be applied to all non-gas species, and represents the total amount of the basis species contained as free species as well as the amount found in secondary species and in kinetic species, and units must be moles or mass (kg, g, etc). free_concentration: can be applied to all basis species that are not gas and not H2O and not mineral, and represents the total amount of the basis species existing freely (not as secondary species) within the solution, and units must be molal or mass_per_kg_solvent. free_mineral: can be applied to all mineral basis species, and represents the total amount of the mineral existing freely (precipitated) within the solution, and units must be moles, mass or cm3. activity and log10activity: can be applied to basis species that are not gas and not mineral and not sorbing sites, and represents the activity of the basis species (recall pH = -log10activity), and units must be dimensionless. fugacity and log10fugacity: can be applied to gases, and units must be dimensionless"}>>> = "kg_solvent_water log10activity bulk_composition"
constraint_unit<<<{"description": "Units of the numerical values given in constraint_value. Dimensionless: should only be used for activity or fugacity constraints. Moles: mole number. Molal: moles per kg solvent water. kg: kilograms. g: grams. mg: milligrams. ug: micrograms. kg_per_kg_solvent: kilograms per kg solvent water. g_per_kg_solvent: grams per kg solvent water. mg_per_kg_solvent: milligrams per kg solvent water. ug_per_kg_solvent: micrograms per kg solvent water. cm3: cubic centimeters"}>>> = " kg dimensionless moles"
ramp_max_ionic_strength_initial<<<{"description": "The number of iterations over which to progressively increase the maximum ionic strength (from zero to max_ionic_strength) during the initial equilibration. Increasing this can help in convergence of the Newton process, at the cost of spending more time finding the aqueous configuration."}>>> = 0 # max_ionic_strength in such a simple problem does not need ramping
abs_tol<<<{"description": "If the residual of the algebraic system (measured in mol) is lower than this value, the Newton process (that finds the aqueous configuration) is deemed to have converged"}>>> = 1E-15
[]
[nnn]
type = NearestNodeNumberUO<<<{"description": "Finds and outputs the nearest node number to a point", "href": "../../../source/userobjects/NearestNodeNumberUO.html"}>>>
[]
[]
[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[console_output]
type = GeochemistryConsoleOutput<<<{"description": "Outputs results from a GeochemistryReactor at a particular point", "href": "../../../source/outputs/GeochemistryConsoleOutput.html"}>>>
geochemistry_reactor<<<{"description": "The name of the GeochemistryReactor UserObject"}>>> = reactor
nearest_node_number_UO<<<{"description": "The NearestNodeNumber UserObject that defines the physical point at which to query the GeochemistryReactor"}>>> = nnn
solver_info<<<{"description": "Print information (to the console) from the solver including residuals, swaps, etc"}>>> = true
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = initial
[]
[]
(moose/modules/geochemistry/test/tests/equilibrium_models/HCl_no_action.i)