pygcc.pygcc_utils

Created on Wed Mar 17 16:02:22 2021

@author: Adedapo Awolayo and Ben Tutolo, University of Calgary

Copyright (c) 2020 - 2021, Adedapo Awolayo and Ben Tutolo, University of Calgary

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

Module Contents

Classes

calcRxnlogK

This class implemetation calculates logK values for any reaction with the option for extrapolation where rho < 350kg/m3

write_database

Class to write the new database for either GWB, EQ3/6, ToughReact, Pflotran into

Functions

calc_elem_count_molewt(formula, **kwargs)

This function calculates the molecular mass and the elemental composition of a substance given

importconfile(filename, *Rows)

This function imports numeric data from a text file as column vectors.

var_name(var)

feval(funcName, *args)

roundup_tenth(x)

roundup_hundredth(x)

read_specific_lines(file, lines_to_read)

derivative(f, a[, method, h])

Compute the derivative of f, f'(a) with step size h.

info(name, dic)

This function checks for naming convention of the species in the direct-access or source database

drummondgamma(TK, I)

This function models solubility of CO2 gas in brine using Drummond equation

Henry_duan_sun(TK, P, I)

This function evaluates the solubility of CO2 gase in brine using Duan_Sun Formulation

co2fugacity(TK, P[, poy])

This function computes the fugacity and density of CO2 by Duan and Sun 2003

gamma_correlation(TC, P[, method])

This function calculates the CO2 activity correlation coefficients at

Helgeson_activity(TC, P, I[, Dielec_method])

This function calculates the solute activity coefficient, solvent osmotic coefficient,

aw_correlation(TC, P[, Dielec_method])

Calculates the water activity correlation coefficients at given temperature and pressure

outputfmt(fid, logK, Rxn, *T[, dataset, logK_form])

This function writes logK and Rxn data to any file using GWB, EQ36, Pflotran and ToughReact format

main_function_name(module)

print main_function name

Attributes

eps

J_to_cal

IAPWS95_COEFFS

pygcc.pygcc_utils.eps = 2.220446049250313e-16
pygcc.pygcc_utils.J_to_cal = 4.184
pygcc.pygcc_utils.IAPWS95_COEFFS
pygcc.pygcc_utils.calc_elem_count_molewt(formula, **kwargs)[source]

This function calculates the molecular mass and the elemental composition of a substance given by its chemical formula. This is modified from https://github.com/cgohlke/molmass/blob/master/molmass/molmass.py

Parameters:
  • formula (string) – chemical formula

  • Elementdic (dict) – dictionary, containing the atomic mass of element database

Returns:

  • elements (dict) – dictionary of elemental composition and their respective number of atoms

  • molewt (float) – Calculated Molecular Weights [g/mol]

Usage:

[elements, molewt] = calc_elem_count_molewt(formula)

Examples of valid formulas are “H2O”, “[2H]2O”, “CH3COOH”, “EtOH”, “CuSO4:5H2O”, “(COOH)2”, “AgCuRu4(H)2[CO]12{PPh3}2”, “CGCGAATTCGCG”, and, “MDRGEQGLLK” .

Examples

>>> elements, molewt = calc_elem_count_molewt("CuSO4:5H2O")
>>> elements, molewt
    {'O': 9.0, 'H': 10.0, 'S': 1, 'Cu': 1}, 249.6773
pygcc.pygcc_utils.importconfile(filename, *Rows)[source]
This function imports numeric data from a text file as column vectors.

[Var1, Var2] = importconfile(filename) Reads data from text file filename for the default selection.

[Var1, Var2] = importconfile(filename, StartRow, EndRow) Reads data from rows StartRow through EndRow of text file filename.

Example:

[Var1, Var2] = importconfile(‘100bar.con’, 6, 13);

pygcc.pygcc_utils.var_name(var)[source]
pygcc.pygcc_utils.feval(funcName, *args)[source]
pygcc.pygcc_utils.roundup_tenth(x)[source]
pygcc.pygcc_utils.roundup_hundredth(x)[source]
pygcc.pygcc_utils.read_specific_lines(file, lines_to_read)[source]
pygcc.pygcc_utils.derivative(f, a, method='central', h=0.001)[source]

Compute the derivative of f, f’(a) with step size h. :param f: Vectorized function of one variable :type f: function :param a: Compute derivative at x = a :type a: number :param method: Difference formula: ‘forward’, ‘backward’ or ‘central’ :type method: string :param h: Step size in difference formula :type h: number

Returns:

Difference formula:

central: f(a + h) - f(a - h))/2h forward: f(a + h) - f(a))/h backward: f(a) - f(a-h))/h

Return type:

float

pygcc.pygcc_utils.info(name, dic)[source]

This function checks for naming convention of the species in the direct-access or source database

Parameters:
  • name (string) – species name

  • dic (dict) – dictionary of species from direct-access or source database

Returns:

lst – resulting search of all species with the input name

Return type:

list

Examples

>>> from pygcc.pygcc_utils import db_reader
>>> ps = db_reader() # utilizes the default direct-access database, speq21
>>> info('ss_', ps.dbaccessdic)
    ['ss_Anorthite',  'ss_Albite_high', 'ss_K-feldspar', 'ss_Ferrosilite',
     'ss_Enstatite', 'ss_Clinoenstatite', 'ss_Hedenbergite', 'ss_Diopside',
     'ss_Forsterite', 'ss_Fayalite', 'ss_Annite',  'ss_Phlogopite',
     'ss_Anorthite',  'ss_Albite_high', 'ss_K-feldspar', 'ss_Ferrosilite',
     'ss_Enstatite', 'ss_Clinoenstatite', 'ss_Hedenbergite', 'ss_Diopside',
     'ss_Forsterite',  'ss_Fayalite', 'ss_Annite', 'ss_Phlogopite']
pygcc.pygcc_utils.drummondgamma(TK, I)[source]

This function models solubility of CO2 gas in brine using Drummond equation

Parameters:
  • TK (float, vector) – Temperature [K]

  • I – Ionic strength

Examples

>>> log10_gamma = drummondgamma( 500, 0.5)
>>> log10_gamma
    0.0781512920184902
pygcc.pygcc_utils.Henry_duan_sun(TK, P, I)[source]

This function evaluates the solubility of CO2 gase in brine using Duan_Sun Formulation

Parameters:
  • TK (float, vector) – Temperature [K]

  • P (float, vector) – Pressure [bar]

  • I (float) – Ionic strength

Returns:

  • log10_co2_gamma (float, vector) – co2 aqueous activity coefficients in log10

  • mco2 (float, vector) – co2 aqueous molalities

Usage

log10_co2_gamma, mco2 = Henry_duan_sun( TK, P, I)

Examples

>>> log10_co2_gamma, mco2 = Henry_duan_sun( 500, 250, 0.5)
>>> log10_co2_gamma, mco2
    0.06202532, 1.6062557

References

  1. Duan, Z. and Sun, R., 2003. An improved model calculating CO2 solubility in pure water and aqueous NaCl solutions from 273 to 533 K and from 0 to 2000 bar. Chemical geology, 193(3-4), pp.257-271.

pygcc.pygcc_utils.co2fugacity(TK, P, poy=True)[source]

This function computes the fugacity and density of CO2 by Duan and Sun 2003 Also Calculate the Saturation Index SI and partial pressure of CO2(g) at any given T, P using the Duan equation of state. A Poynting correction factor is also applied. :param TK: Temperature [K] :type TK: float, vector :param P: pressure [bar] :type P: float, vector :param poy: Option to apply a Poynting correction factor [True is the default] :type poy: bool

Returns:

  • fCO2 (float, vector) – co2 fugacity [bar]

  • denCO2 (float, vector) – co2 density [g/cm3]

  • SI (float, vector) – co2 saturation index

  • pCO2 (float, vector) – co2 partial pressure [bar]

Usage

[fCO2, denCO2, SI, pCO2] = co2fugacity( TK, P)

Examples

>>> fCO2, denCO2, SI, pCO2 = co2fugacity( np.array([400,420]), np.array([150, 200]))
>>> fCO2
    array([114.50387423, 150.87672517])
>>> denCO2
    array([2.6741242 , 3.33134397])
>>> SI
    array([1.73291553, 1.84597591])
>>> pCO2
    array([54.78127571, 71.0710152])
pygcc.pygcc_utils.gamma_correlation(TC, P, method=None)[source]

This function calculates the CO2 activity correlation coefficients at given temperature T and pressure P

Parameters:
  • TC (float, vector) – Temperature [°C]

  • P (float, vector) – Pressure [bar]

  • method – specify the activity model ['Duan_Sun' or 'Drummond']

Usage

Examples

>>> log10_co2_gamma, mco2 = Henry_duan_sun( 500, 250, 0.5)
>>> log10_co2_gamma, mco2
    0.06202532, 1.6062557

References

  1. Segal Edward Drummond, 1981, Boiling and Mixing of Hydrothermal Fluids: Chemical Effects on Mineral Precipitation, page 19

  2. Wolery, T. J., Lawrence Livermore National Laboratory, United States Dept. of Energy, 1992. EQ3/6: A software package for geochemical modeling of aqueous systems: package overview and installation guide (version 7.0)

pygcc.pygcc_utils.Helgeson_activity(TC, P, I, Dielec_method=None, **rhoEDB)[source]

This function calculates the solute activity coefficient, solvent osmotic coefficient, and solvent activity at given temperature and pressure using equations 298, 190 and 106 in Helgeson, Kirkham and Flowers, 1981, A.J.S. p.1249-1516 Parameters: ———-

TC : Temperature [°C] P : pressure [bar] I : ionic strength Dielec_method : specify either ‘FGL97’ or ‘JN91’ or ‘DEW’ as the method to calculate dielectric

constant (optional), if not specified default - ‘JN91’

rhoEDBdictionary of water properties like density (rho), dielectric factor (E) and

Debye–Hückel coefficients (optional)

Returns:

aw : solvent activity phi : solvent osmotic coefficient mean_act : solute activity coefficient

Usage

[aw, phi, mean_act] = Helgeson_activity( TC, P, I)

pygcc.pygcc_utils.aw_correlation(TC, P, Dielec_method=None, **rhoEDB)[source]

Calculates the water activity correlation coefficients at given temperature and pressure Parameters: ———-

TC : Temperature [°C] P : pressure [bar] Dielec_method : specify either ‘FGL97’ or ‘JN91’ or ‘DEW’ as the method to calculate dielectric

constant (optional), if not specified default - ‘JN91’

rhoEDBdictionary of water properties like density (rho), dielectric factor (E) and

Debye–Hückel coefficients (optional)

Returns:

ch20 : water activity correlation coefficients

Usage

[ch20] = aw_correlation( TC, P)

class pygcc.pygcc_utils.calcRxnlogK(**kwargs)[source]

This class implemetation calculates logK values for any reaction with the option for extrapolation where rho < 350kg/m3

Parameters:
  • T (float, vector) – Temperature [°C]

  • P (float, vector) – Pressure [bar]

  • Specie (string) – specify the species for logK calculation, either the Product species of any reaction or solid solutions or clay like ‘AnAb’ or ‘AbOr’ or ‘FoFa’ or ‘EnFe’ or ‘DiHedEnFe’ or ‘clay’

  • Specie_class (string, optional) – specify the class of species like ‘aqueous’, ‘minerals’, ‘liquids’, or ‘gases’

  • elem (list) – list containing nine parameters with clay names and elements compositions with the following format [‘Montmorillonite_Lc_MgK’, ‘Si’, ‘Al’, ‘FeIII’, ‘FeII’, ‘Mg’, ‘K’, ‘Na’, ‘Ca’, ‘Li’]

  • dbaccessdic (dict) – direct-acess database dictionary

  • rhoEGextrap (dict) – dictionary of water properties like density (rho), dielectric factor (E) and Gibbs Energy for density region 350-550kg/m3

  • group (string) – specify the structural layering of the phyllosilicate, for layers composed of 1 tetrahedral + 1 octahedral sheet (1:1 layer) - specify '7A', 2 tetrahedral + 1 octahedral sheet (2:1 layer) - specify '10A', or the latter with a brucitic sheet in the interlayer (2:1:1 layer) - specify '14A' (optional), if not specified, default is ‘10A’ for smectites, micas, et cetera

  • ClayMintype (string) – specify either ‘Smectite’ or ‘Chlorite’ or ‘Mica’ as the clay type, if not specified default - ‘Smectites’

  • X (float) – volume fractions of any (Anorthite, Albite, Forsterite, Enstatite) or mole fraction of Mg

  • cpx_Ca (float) – number of moles of Ca in formula unit (=1 for Di, Hed), must be greater than zero

  • sourcedic (dict) – source database reactions dictionary

  • specielist (list of list, optional) – source database species grouped into categories [element, basis, redox, aqueous, minerals, gases, oxides]

  • Dielec_method (string) – specify either ‘FGL97’ or ‘JN91’ or ‘DEW’ as the method to calculate dielectric constant, default is ‘JN91’

  • heatcap_method (string) – specify either ‘SUPCRT’ or ‘Berman88’ or ‘HP11’ or ‘HF76’ as the method to calculate thermodynamic properties of any mineral or gas, default is ‘SUPCRT’

  • ThermoInUnit (string) – specify either ‘cal’ or ‘KJ’ as the input units for species properties (optional), particularly used to covert KJ data to cal by supcrtaq function if not specified default - ‘cal’

  • Al_Si (string) – specify either ‘pygcc’ or ‘Arnórsson_Stefánsson’ as the input to express Al and Si species in solid solution (optional), ‘Arnórsson_Stefánsson’ expresses them as ‘Al(OH)4-’ and ‘H4SiO4(aq)’, respectively while pygcc uses ‘Al3+’ and ‘SiO2(aq)’, respectively if not specified default - ‘pygcc’

  • sourceformat (string) – specify the source database format, either ‘GWB’ or ‘EQ36’

  • densityextrap (float, vector) – specify the extrapolation option for density-logK, ‘Yes’/True or ‘No’/False

Returns:

  • log_K (float, vector) – logarithmic K value(s)

  • dGrxn (float, vector) – Total reaction Gibbs energy [cal/mol]

Usage

The general usage of water_dielec is as follows:

  1. For water dielectric properties at any Temperature and Pressure:

    calclogK = calcRxnlogK(T = T, P = P, Dielec_method = ‘JN91’, **kwargs),

    where T is temperature in celsius and P is pressure in bar

  2. For water dielectric properties at any Temperature and density :

    calclogK = calcRxnlogK(T = T, rho = rho, Dielec_method = ‘JN91’, **kwargs),

    where T is temperature in celsius and rho is density in kg/m³

  3. For water dielectric properties at any Temperature and Pressure on steam saturation curve:

    calclogK = calcRxnlogK(T = T, P = ‘T’, Dielec_method = ‘JN91’, **kwargs),

    where T is temperature in celsius, followed with a quoted character ‘T’ to reflect steam saturation pressure

    calclogK = calcRxnlogK(P = P, T = ‘P’, Dielec_method = ‘JN91’, **kwargs),

    where P is pressure in bar, followed with a quoted character ‘P’ to reflect steam saturation temperature

Examples

>>> ps = db_reader(sourcedb = './default_db/thermo.com.dat',
                   sourceformat = 'gwb', dbaccess = './default_db/speq21.dat')
>>> calclogK = calcRxnlogK(T = 100, P = 50, Specie = 'H2S(aq)',
                           Specie_class = 'aqueous',
                           dbaccessdic = ps.dbaccessdic,
                           sourcedic = ps.sourcedic, specielist = ps.specielist)
>>> calclogK.logK, calclogK.dGrxn
     -6.4713,       11049.3168
kwargs
__calc__(**kwargs)[source]
AllRxnslogK(TC, P, rhoEG)[source]

This function calculates logK values of all reactions including solid-solution and clay minerals

Parameters:
  • TC (temperature [°C]) –

  • P (pressure [bar]) –

  • rhoEG (dictionary of water properties like density (rho),) – dielectric factor (E) and Gibbs Energy (optional)

Returns:

  • logK (logarithmic K value(s))

  • dGrxn (Total reaction Gibbs energy [cal/mol])

  • dGP (Product specie Gibbs energy [cal/mol])

  • dGRs (Reactant species Gibbs energy [cal/mol])

  • Rxn (dict)

  • The calculated dictionary of reaction thermodynamic properties has the following properties

    • type: solid-solution mineral type, [-]

    • name: solid-solution name, [K]

    • formula: solid-solution mineral formula, [-]

    • MW: Molecular weight, [g/mol]

    • min: solid-solution mineral properties, [‘formula’, ‘source date’, dG[cal/ml], dH[cal/mol], S[cal/mol-K], V[cm3/mol], a[cal/mol-K], b[10^3 cal/mol/K^2], c[10^-5 cal/mol/K]]

    • spec: list of species, [-]

    • coeff: list of corresponding coefficients of species above, [-]

    • nSpec: Total number of species, [-]

    • V: Molar volume, [cm3/mol]

    • source: Source of thermo data, [kJ/kg·K]

    • elements: list of elements and their total numbers, [-]

Usage

The general usage of calcRxnlogK without the optional argument is as follows:

  1. Not on steam saturation curve:

    [logK, dGrxn, dGP, dGRs] = AllRxnslogK(TC, P, Prod, dbaccessdic, sourcedic, specielist),

    where TC is temperature in celsius and P is pressure in bar;

  2. On steam saturation curve:

    [logK, dGrxn, dGP, dGRs] = AllRxnslogK(TC, ‘T’, Prod, dbaccessdic, sourcedic, specielist),

    where TC is temperature in celsius, followed with a quoted char ‘T’

    [logK, dGrxn, dGP, dGRs] = AllRxnslogK(P, ‘P’, Prod, dbaccessdic, sourcedic, specielist),

    where P is pressure in bar, followed with a quoted char ‘P’.

  3. Meanwhile, usage with any specific dielectric constant method (‘FGL97’) for condition not on steam saturation curve is as follows. Default method is ‘JN91’

    [logK, dGrxn, dGP, dGRs] = AllRxnslogK( TC, P, Prod, dbaccessdic, sourcedic, specielist, Dielec_method = ‘FGL97’)

RxnlogK(TC, P, Prod, rhoEG)[source]

This function calculates logK values of any reaction

Parameters:
  • TC (temperature [°C]) –

  • P (pressure [bar]) –

  • Prod (Product species of the reaction) –

  • dbaccessdic (direct-acess database dictionary) –

  • sourcedic (source database reactions dictionary) –

  • specielist (source database species grouped into) – [element, basis, redox, aqueous, minerals, gases, oxides]

  • Dielec_method (specify either 'FGL97' or 'JN91' or 'DEW' as the method) – to calculate dielectric constant, default is ‘JN91’

  • sourceformat (source database format, either 'GWB' or 'EQ36', default is 'GWB') –

  • heatcap_method (specify either 'SUPCRT' or 'Berman88' or 'HP11' or 'HF76' as the method to calculate thermodynamic properties) – of any mineral, default is ‘SUPCRT’

  • rhoEG (dictionary of water properties like density (rho),) – dielectric factor (E) and Gibbs Energy (optional)

Returns:

  • logK (logarithmic K value(s))

  • dGrxn (Total reaction Gibbs energy [cal/mol])

  • dGP (Product specie Gibbs energy [cal/mol])

  • dGRs (Reactant species Gibbs energy [cal/mol])

Usage

The general usage of calcRxnlogK without the optional argument is as follows:

  1. Not on steam saturation curve:

    [logK, dGrxn, dGP, dGRs] = RxnlogK(TC, P, Prod, dbaccessdic, sourcedic, specielist),

    where TC is temperature in celsius and P is pressure in bar;

  2. On steam saturation curve:

    [logK, dGrxn, dGP, dGRs] = RxnlogK(TC, ‘T’, Prod, dbaccessdic, sourcedic, specielist),

    where TC is temperature in celsius, followed with a quoted char ‘T’

    [logK, dGrxn, dGP, dGRs] = RxnlogK(P, ‘P’, Prod, dbaccessdic, sourcedic, specielist),

    where P is pressure in bar, followed with a quoted char ‘P’.

  3. Meanwhile, usage with any specific dielectric constant method (‘FGL97’) for condition not on steam saturation curve is as follows. Default method is ‘JN91’

    [logK, dGrxn, dGP, dGRs] = RxnlogK( TC, P, Prod, dbaccessdic, sourcedic, specielist, Dielec_method = ‘FGL97’)

densitylogKextrap(TC, P, rhoEGextrap)[source]

This function calculates logK values extrapolation for conditions where rho < 350kg/m3

Parameters:
  • TC (temperature [°C]) –

  • P (pressure [bar]) –

  • rhoEGextrap (dictionary of water properties like density (rho), dielectric factor (E) and) – Gibbs Energy for density region 350-550kg/m3

Returns:

logK

Return type:

extrapolated logarithmic K value(s)

Usage

[logK] = densitylogKextrap(TC, P, ‘H2S(aq)’, dbaccessdic = dbaccessdic, sourcedic = sourcedic, specielist = specielist),

pygcc.pygcc_utils.outputfmt(fid, logK, Rxn, *T, dataset=None, logK_form=None)[source]

This function writes logK and Rxn data to any file using GWB, EQ36, Pflotran and ToughReact format

Parameters:
  • fid (string) – file ID

  • logK (float, vector) – logarithmic K value(s)

  • Rxn (dict) – dictionary of reaction thermodynamic properties

  • T (float, vector) – Temperature value(s), optional, required when ‘polycoeffs’ is specified for logK_form

  • dataset (string) – specify the dataset format, either ‘GWB’, ‘EQ36’, ‘Pflotran’ or ‘ToughReact’

  • logK_form (string) – specify the format of logK either as a set of eight values one for each of the dataset’s principal temperatures, or blocks of polynomial coefficients, [values, polycoeffs], default is ‘a set of eight values’ (optional)

Return type:

Output data to the file with filename described in fid with any format mentioned above.

Examples

>>> fid = open('./logK_details.txt', 'w')
>>> logKRxn = calcRxnlogK(T = 100, P = 'T', X = 0.634,
                          Specie = 'Plagioclase', densityextrap = True)
>>> # output in EQ36 format
>>> outputfmt(fid, logKRxn.logK, logKRxn.Rxn, dataset = 'EQ36')
>>> fid.close()
class pygcc.pygcc_utils.write_database(**kwargs)[source]

Class to write the new database for either GWB, EQ3/6, ToughReact, Pflotran into a new folder called “output”

Parameters:
  • T (float, vector) – Temperature [°C]

  • P (float, vector) – Pressure [bar]

  • cpx_Ca (string) – number of moles of Ca in solid solution of clinopyroxene (optional) if it is ommitted solid solution of clinopyroxene will not be included, 0 < nCa >=1

  • solid_solution (string, bool) – specify the inclusion of solid-solution [Yes/True or No/False], default is ‘No’

  • clay_thermo (float, vector) – specify the inclusion of clay thermodynamic properties [Yes/True or No/False], default is ‘No’

  • logK_form (float, vector) – specify the format of logK either as a set of eight values one for each of the dataset’s principal temperatures, or blocks of polynomial coefficients, [values, polycoeffs] default is ‘a set of eight values’

  • densityextrap (float, vector) – specify the utilization of density extrapolation [Yes/True or No/False], default is ‘Yes’

  • dbaccess (string) – direct-access database filename and location (optional)

  • dbBerman_dir (string) – filename and location of the Berman mineral database (optional)

  • dbHP_dir (string) – filename and location of the supcrtbl mineral and gas database, optional

  • dbaccessformat (string, optional) – specify the direct-access/sequential-access database format, either ‘speq’ or ‘supcrtbl’, default is ‘speq’

  • sourcedb (string) – source database filename and location (optional)

  • sourceformat (string) – source database format, either ‘GWB’ or ‘EQ36’, default is ‘GWB’

  • sourcedb_codecs (string) – specify the name of the encoding used to decode or encode the sourcedb file, optional

  • objdb (string) – new database filename and location (optional)

  • co2actmodel (string) – co2 activity model equation [Duan_Sun or Drummond] (optional), if not specified, default is ‘Drummond’

  • Dielec_method (string) – specify either ‘FGL97’ or ‘JN91’ or ‘DEW’ as the method to calculate dielectric constant, default is ‘JN91’ (optional)

  • heatcap_method (string) – specify either ‘SUPCRT’ or ‘Berman88’ or ‘HP11’ or ‘HF76’ the method to calculate thermodynamic properties of any mineral, default is ‘SUPCRT’

  • ThermoInUnit (string) – specify either ‘cal’ or ‘KJ’ as the input units for species properties (optional), particularly used to covert KJ data to cal by supcrtaq function if not specified default - ‘cal’

  • dataset (string) – specify the dataset format, either ‘GWB’, ‘EQ36’, ‘Pflotran’ or ‘ToughReact’, default is old GWB database [‘GWB’] (optional)

  • print_msg (string, bool) – print debug message [True or False], default is False

Return type:

Output the new database to an ASCII file with filename described in ‘objdb’ if specified.

Usage

With any Temperature and Pressure:
  1. General format with default dielectric constant and CO2 activity model and exclusions of solid solutions for GWB

    write_database(T = T, P = P, cpx_Ca = nCa, dataset = ‘GWB’, sourceformat = ‘GWB’)

  2. Inclusion of solid solutions and exclusion of solid solution of clinopyroxene and clay thermo

    write_database(T = T, P = P, solid_solution = ‘Yes’, clay_thermo = ‘Yes’, dataset = ‘GWB’, sourceformat = ‘GWB’)

  3. Inclusion of all solid solutions and clay thermo with emph{‘Duan_Sun’} CO2 activity model and ‘FGL97’ dielectric constant calculation

    write_database(T = T, P = P, cpx_Ca = nCa, solid_solution = ‘Yes’, clay_thermo = ‘Yes’, co2actmodel = ‘Duan_Sun’, Dielec_method = ‘FGL97’, dataset = ‘GWB’, sourceformat = ‘GWB’)

With any Temperature or Pressure on the steam saturation curve:
  1. General format with default dielectric constant and CO2 activity model

    write_database(T = T, P = ‘T’, cpx_Ca = nCa, dataset = ‘GWB’, sourceformat = ‘GWB’),

    where T is temperature in celsius and P is assigned a quoted character ‘T’ to reflect steam saturation pressure

    write_database(P = P, T = ‘P’, cpx_Ca = nCa, dataset = ‘GWB’, sourceformat = ‘GWB’),

    where P is pressure in bar and T is assigned a quoted character ‘P’ to reflect steam saturation temperature

Examples

>>> write_database(T=np.array([ 0.010, 25., 60., 100., 150.,  200.,
                               250.,  300.]), P=200, cpx_Ca = 0.5,
                   sourceformat = 'GWB', solid_solution = True,
                   clay_thermo = True, dataset = 'GWB')
>>> write_database(T=np.array([ 0.010, 25., 60., 100., 150.,  200.,  250.,
                               300.]), P=200, cpx_Ca = 0.5,
                   sourceformat = 'GWB', solid_solution = True,
                   clay_thermo = True, dataset = 'Pflotran',
                   sourcedb = './pygcc/default_db/thermo.2021.dat')
kwargs
__calc__(**kwargs)[source]
write_GWBdb(T, P)[source]

This function writes the new GWB database into a new folder called “output”

Parameters:
  • T (temperature [°C]) –

  • P (pressure [bar]) –

Return type:

Outputs the new database to an ASCII file with filename described in ‘objdb’.

Usage
Example:
  1. General format with default dielectric constant and CO2 activity model and exclusions of solid solutions

    write_GWBdb(T, P)

  2. Inclusion of solid solutions and clay thermo and exclusion of solid solution of clinopyroxene

    write_GWBdb(T, P)

  3. Inclusion of all solid solutions and clay thermo with emph{‘Duan_Sun’} CO2 activity model and ‘JN91’ dielectric constant calculation

    write_GWBdb(T, P)

write_EQ36db(T, P)[source]

This function writes the new EQ3/6 database into a new folder called “output”

Parameters:
  • T (temperature [°C]) –

  • P (pressure [bar]) –

Return type:

Outputs the new database to an ASCII file with filename described in ‘objdb’.

Usage
Example:
  1. General format with default dielectric constant and CO2 activity model and exclusions of solid solutions

    write_EQ36db(T, P )

  2. Inclusion of solid solutions and clay thermo and exclusion of solid solution of clinopyroxene

    write_EQ36db(T, P )

  3. Inclusion of all solid solutions and clay thermo with emph{‘Duan_Sun’} CO2 activity model and ‘FGL97’ dielectric constant calculation

    write_EQ36db(T, P )

write_pflotrandb(T, P)[source]

This function writes the new pflotran database into a new folder called “output”

Parameters:
  • T (temperature [°C]) –

  • P (pressure [bar]) –

Return type:

Outputs the new database to an ASCII file with filename described in ‘objdb’.

Usage
Example:
  1. General format with default dielectric constant and CO2 activity model and exclusions of solid solutions

    write_pflotrandb(T, P )

  2. Inclusion of solid solutions and clay thermo and exclusion of solid solution of clinopyroxene

    write_pflotrandb(T, P )

  3. Inclusion of all solid solutions and clay thermo with emph{‘Duan_Sun’} CO2 activity model and ‘FGL97’ dielectric constant calculation

    write_pflotrandb(T, P )

write_ToughReactdb(T, P)[source]

This function writes the new ToughReact database into a new folder called “output”

Parameters:
  • T (temperature [°C]) –

  • P (pressure [bar]) –

Return type:

Outputs the new database to an ASCII file with filename described in ‘objdb’.

Usage
Example:
  1. General format with default dielectric constant and CO2 activity model and exclusions of solid solutions

    write_ToughReactdb(T, P, dbaccess = ‘location’, sourcedb = ‘location’,

    objdb = ‘location’, sourceformat = ‘GWB’)

  2. Inclusion of solid solutions and clay thermo and exclusion of solid solution of clinopyroxene

    write_ToughReactdb(T, P, solid_solution = ‘Yes’, clay_thermo = ‘Yes’, dbaccess = ‘location’,

    sourcedb = ‘location’, objdb = ‘location’, sourceformat = ‘GWB’)

  3. Inclusion of all solid solutions and clay thermo with emph{‘Duan_Sun’} CO2 activity model and ‘FGL97’ dielectric constant calculation

    write_ToughReactdb(T, P, nCa, solid_solution = ‘Yes’, clay_thermo = ‘Yes’, dbaccess = ‘location’,

    sourcedb = ‘location’, objdb = ‘location’, co2actmodel = ‘Duan_Sun’, Dielec_method = ‘FGL97’, sourceformat = ‘GWB’)

pygcc.pygcc_utils.main_function_name(module)[source]

print main_function name