pyEQL package

Submodules

pyEQL.activity_correction module

pyEQL activity correction library.

This file contains functions for computing molal-scale activity coefficients of ions and salts in aqueous solution.

Individual functions for activity coefficients are defined here so that they can be used independently of a pyEQL solution object. Normally, these functions are called from within the get_activity_coefficient method of the Solution class.

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

pyEQL.activity_correction.get_activity_coefficient_davies(ionic_strength, z=1, temperature='25 degC')[source]

Return the activity coefficient of solute in the parent solution according to the Davies equation.

Parameters:
  • ionic_strength (Quantity) – The ionic strength of the parent solution, mol/kg.

  • z (int, optional) – The charge on the solute, including sign. Defaults to +1 if not specified.

  • temperature (str, Quantity, optional) – String representing the temperature of the solution. Defaults to ‘25 degC’ if not specified.

Returns:

The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless.

Return type:

Quantity

Notes

Activity coefficient is calculated according to:

\[\ln \gamma = A^{\gamma} z_i^2 ({\sqrt I \over (1 + \sqrt I)} + 0.2 I)\]

Valid for 0.1 < I < 0.5

References

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed,

pp 103. Wiley Interscience, 1996.

pyEQL.activity_correction.get_activity_coefficient_debyehuckel(ionic_strength, z=1, temperature='25 degC')[source]

Return the activity coefficient of solute in the parent solution according to the Debye-Huckel limiting law.

Parameters:
  • z (int, optional) – The charge on the solute, including sign. Defaults to +1 if not specified.

  • ionic_strength (Quantity) – The ionic strength of the parent solution, mol/kg.

  • temperature (str, Quantity, optional) – String representing the temperature of the solution. Defaults to ‘25 degC’ if not specified.

Returns:

The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless.

Return type:

Quantity

Notes

Activity coefficient is calculated according to:

\[\ln \gamma = A^{\gamma} z_i^2 \sqrt I\]

Valid only for I < 0.005

References

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed,

pp 103. Wiley Interscience, 1996.

pyEQL.activity_correction.get_activity_coefficient_guntelberg(ionic_strength, z=1, temperature='25 degC')[source]

Return the activity coefficient of solute in the parent solution according to the Guntelberg approximation.

Parameters:
  • z (int, optional) – The charge on the solute, including sign. Defaults to +1 if not specified.

  • ionic_strength (Quantity) – The ionic strength of the parent solution, mol/kg.

  • temperature (str, Quantity, optional) – String representing the temperature of the solution. Defaults to ‘25 degC’ if not specified.

Returns:

The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless.

Return type:

Quantity

Notes

Activity coefficient is calculated according to:

\[\ln \gamma = A^{\gamma} z_i^2 {\sqrt I \over (1 + \sqrt I)}\]

Valid for I < 0.1

References

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed,

pp 103. Wiley Interscience, 1996.

pyEQL.activity_correction.get_activity_coefficient_pitzer(ionic_strength, molality, alpha1, alpha2, beta0, beta1, beta2, C_phi, z_cation, z_anion, nu_cation, nu_anion, temperature='25 degC', b=1.2)[source]

Return the activity coefficient of solute in the parent solution according to the Pitzer model.

Parameters:
  • ionic_strength – The ionic strength of the parent solution, mol/kg

  • molality – The molal concentration of the parent salt, mol/kg

  • alpha1 – Coefficients for the Pitzer model. This function assigns the coefficients proper units of kg ** 0.5 / mol ** 0.5 after they are entered.

  • alpha2 – Coefficients for the Pitzer model. This function assigns the coefficients proper units of kg ** 0.5 / mol ** 0.5 after they are entered.

  • beta0 – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • beta1 – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • beta2 – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • C_phi – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • z_cation – The charge on the cation and anion, respectively

  • z_anion – The charge on the cation and anion, respectively

  • nu_cation – The stoichiometric coefficient of the cation and anion in the salt

  • nu_anion – The stoichiometric coefficient of the cation and anion in the salt

  • temperature – String representing the temperature of the solution. Defaults to ‘25 degC’ if not specified.

  • b – Coefficient. Usually set equal to 1.2 and considered independent of temperature and pressure. If provided, this coefficient is assigned proper units of kg ** 0.5 / mol ** 0.5 after entry.

Returns:

Quantity

The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless

Examples

>>> get_activity_coefficient_pitzer(0.5*ureg.Quantity('mol/kg'),0.5*ureg.Quantity('mol/kg'),1,0.5,-.0181191983,-.4625822071,.4682,.000246063,1,-1,1,1,b=1.2)
0.61915...
>>> get_activity_coefficient_pitzer(5.6153*ureg.Quantity('mol/kg'),5.6153*ureg.Quantity('mol/kg'),3,0.5,0.0369993,0.354664,0.0997513,-0.00171868,1,-1,1,1,b=1.2)
0.76331...

Notes: the examples below are for comparison with experimental and modeling data presented in the May et al reference below.

10 mol/kg ammonium nitrate. Estimated result (from graph) = 0.2725

>>> get_activity_coefficient_pitzer(10*ureg.Quantity('mol/kg'),10*ureg.Quantity('mol/kg'),2,0,-0.01709,0.09198,0,0.000419,1,-1,1,1,b=1.2)
0.22595 ...

5 mol/kg ammonium nitrate. Estimated result (from graph) = 0.3011

>>> get_activity_coefficient_pitzer(5*ureg.Quantity('mol/kg'),5*ureg.Quantity('mol/kg'),2,0,-0.01709,0.09198,0,0.000419,1,-1,1,1,b=1.2)
0.30249 ...

18 mol/kg ammonium nitrate. Estimated result (from graph) = 0.1653

>>> get_activity_coefficient_pitzer(18*ureg.Quantity('mol/kg'),18*ureg.Quantity('mol/kg'),2,0,-0.01709,0.09198,0,0.000419,1,-1,1,1,b=1.2)
0.16241 ...

References

Scharge, T., Munoz, A.G., and Moog, H.C. (2012). Activity Coefficients of Fission Products in Highly Salinary Solutions of Na+, K+, Mg2+, Ca2+, Cl-, and SO42- : Cs+. /Journal of Chemical& Engineering Data (57), p. 1637-1647.

Kim, H., & Jr, W. F. (1988). Evaluation of Pitzer ion interaction parameters of aqueous electrolytes at 25 degree C. 1. Single salt parameters. Journal of Chemical and Engineering Data, (2), 177-184.

May, P. M., Rowland, D., Hefter, G., & Königsberger, E. (2011). A Generic and Updatable Pitzer Characterization of Aqueous Binary Electrolyte Solutions at 1 bar and 25 °C. Journal of Chemical & Engineering Data, 56(12), 5066-5077. doi:10.1021/je2009329

Beyer, R., & Steiger, M. (2010). Vapor Pressure Measurements of NaHCOO + H 2 O and KHCOO + H 2 O from 278 to 308 K and Representation with an Ion Interaction (Pitzer) Model. Journal of Chemical & Engineering Data, 55(2), 830-838. doi:10.1021/je900487a

pyEQL.activity_correction.get_apparent_volume_pitzer(ionic_strength, molality, alpha1, alpha2, beta0, beta1, beta2, C_phi, V_o, z_cation, z_anion, nu_cation, nu_anion, temperature='25 degC', b=1.2)[source]

Return the apparent molar volume of solute in the parent solution according to the Pitzer model.

Parameters:
  • ionic_strength (Quantity) – The ionic strength of the parent solution, mol/kg.

  • molality (Quantity) – The molal concentration of the parent salt, mol/kg.

  • alpha1 (number) – Coefficients for the Pitzer model. This function assigns the coefficients proper units of kg ** 0.5 / mol ** 0.5 after they are entered.

  • alpha2 (number) – Coefficients for the Pitzer model. This function assigns the coefficients proper units of kg ** 0.5 / mol ** 0.5 after they are entered.

  • beta0 (number) – Pitzer coefficients for the apparent molar volume. These ion-interaction parameters are specific to each salt system.

  • beta1 (number) – Pitzer coefficients for the apparent molar volume. These ion-interaction parameters are specific to each salt system.

  • beta2 (number) – Pitzer coefficients for the apparent molar volume. These ion-interaction parameters are specific to each salt system.

  • C_phi (number) – Pitzer coefficients for the apparent molar volume. These ion-interaction parameters are specific to each salt system.

  • V_o (number) – The V^o Pitzer coefficient for the apparent molar volume.

  • z_cation (int) – The formal charge on the cation and anion, respectively.

  • z_anion (int) – The formal charge on the cation and anion, respectively.

  • nu_cation (int) – The stoichiometric coefficient of the cation and anion in the salt.

  • nu_anion (int) – The stoichiometric coefficient of the cation and anion in the salt.

  • temperature (str, Quantity) – String representing the temperature of the solution. Defaults to ‘25 degC’ if not specified.

  • b (number, optional) – Coefficient. Usually set equal to 1.2 and considered independent of temperature and pressure. If provided, this coefficient is assigned proper units of kg ** 0.5 / mol ** 0.5 after entry.

Returns:

The apparent molar volume of the solute, cm ** 3 / mol.

Return type:

Quantity

Examples

Notes: the example below is for comparison with experimental and modeling data presented in the Krumgalz et al reference below.

0.25 mol/kg CuSO4. Expected result (from graph) = 0.5 cm ** 3 / mol

>>> get_apparent_volume_pitzer(1.0*ureg.Quantity('mol/kg'),0.25*ureg.Quantity('mol/kg'),1.4,12,0.001499,-0.008124,0.2203,-0.0002589,-6,2,-2,1,1,b=1.2)
0.404...

1.0 mol/kg CuSO4. Expected result (from graph) = 4 cm ** 3 / mol

>>> get_apparent_volume_pitzer(4.0*ureg.Quantity('mol/kg'),1.0*ureg.Quantity('mol/kg'),1.4,12,0.001499,-0.008124,0.2203,-0.0002589,-6,2,-2,1,1,b=1.2)
4.424...

10.0 mol/kg ammonium nitrate. Expected result (from graph) = 50.3 cm ** 3 / mol

>>> get_apparent_volume_pitzer(10.0*ureg.Quantity('mol/kg'),10.0*ureg.Quantity('mol/kg'),2,0,0.000001742,0.0002926,0,0.000000424,46.9,1,-1,1,1,b=1.2)
50.286...

20.0 mol/kg ammonium nitrate. Expected result (from graph) = 51.2 cm ** 3 / mol

>>> get_apparent_volume_pitzer(20.0*ureg.Quantity('mol/kg'),20.0*ureg.Quantity('mol/kg'),2,0,0.000001742,0.0002926,0,0.000000424,46.9,1,-1,1,1,b=1.2)
51.145...

Notes: the examples below are for comparison with experimental and modeling data presented in the Krumgalz et al reference below.

0.8 mol/kg NaF. Expected result = 0.03

>>> get_apparent_volume_pitzer(0.8*ureg.Quantity('mol/kg'),0.8*ureg.Quantity('mol/kg'),2,0,0.000024693,0.00003169,0,-0.000004068,-2.426,1,-1,1,1,b=1.2)
0.22595 ...

References

May, P. M., Rowland, D., Hefter, G., & Königsberger, E. (2011). A Generic and Updatable Pitzer Characterization of Aqueous Binary Electrolyte Solutions at 1 bar and 25 °C. Journal of Chemical & Engineering Data, 56(12), 5066-5077. doi:10.1021/je2009329

Krumgalz, Boris S., Pogorelsky, Rita (1996). Volumetric Properties of Single Aqueous Electrolytes from Zero to Saturation Concentration at 298.15 K Represented by Pitzer’s Ion-Interaction Equations. Journal of Physical Chemical Reference Data, 25(2), 663-689.

pyEQL.activity_correction.get_osmotic_coefficient_pitzer(ionic_strength, molality, alpha1, alpha2, beta0, beta1, beta2, C_phi, z_cation, z_anion, nu_cation, nu_anion, temperature='25 degC', b=1.2)[source]

Return the osmotic coefficient of water in an electrolyte solution according to the Pitzer model.

Parameters:
  • ionic_strength (Quantity) – The ionic strength of the parent solution, mol/kg.

  • molality (Quantity) – The molal concentration of the parent salt, mol/kg.

  • alpha1 (number) – Coefficients for the Pitzer model. This function assigns the coefficients proper units of kg ** 0.5 / mol ** 0.5 after they are entered.

  • alpha2 (number) – Coefficients for the Pitzer model. This function assigns the coefficients proper units of kg ** 0.5 / mol ** 0.5 after they are entered.

  • beta0 – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • beta1 – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • beta2 – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • C_phi – Coefficients for the Pitzer model. These ion-interaction parameters are specific to each salt system.

  • z_cation (int) – The formal charge on the cation and anion, respectively.

  • z_anion (int) – The formal charge on the cation and anion, respectively.

  • nu_cation (int) – The stoichiometric coefficient of the cation and anion in the salt.

  • nu_anion (int) – The stoichiometric coefficient of the cation and anion in the salt.

  • temperature (str, Quantity) – String representing the temperature of the solution. Defaults to ‘25 degC’ if not specified.

  • b (number, optional) – Coefficient. Usually set equal to 1.2 and considered independent of temperature and pressure. If provided, this coefficient is assigned proper units of kg ** 0.5 / mol ** 0.5 after entry.

Returns:

The osmotic coefficient of water, dimensionless.

Return type:

Quantity

Examples

Experimental value according to Beyer and Stieger reference is 1.3550

>>> get_osmotic_coefficient_pitzer(10.175*ureg.Quantity('mol/kg'),10.175*ureg.Quantity('mol/kg'),1,0.5,-.0181191983,-.4625822071,.4682,.000246063,1,-1,1,1,b=1.2)
1.3552 ...

Experimental value according to Beyer and Stieger reference is 1.084

>>> get_osmotic_coefficient_pitzer(5.6153*ureg.Quantity('mol/kg'),5.6153*ureg.Quantity('mol/kg'),3,0.5,0.0369993,0.354664,0.0997513,-0.00171868,1,-1,1,1,b=1.2)
1.0850 ...

Notes: the examples below are for comparison with experimental and modeling data presented in the May et al reference below.

10 mol/kg ammonium nitrate. Estimated result (from graph) = 0.62

>>> get_osmotic_coefficient_pitzer(10*ureg.Quantity('mol/kg'),10*ureg.Quantity('mol/kg'),2,0,-0.01709,0.09198,0,0.000419,1,-1,1,1,b=1.2)
0.6143 ...

5 mol/kg ammonium nitrate. Estimated result (from graph) = 0.7

>>> get_osmotic_coefficient_pitzer(5*ureg.Quantity('mol/kg'),5*ureg.Quantity('mol/kg'),2,0,-0.01709,0.09198,0,0.000419,1,-1,1,1,b=1.2)
0.6925 ...

18 mol/kg ammonium nitrate. Estimated result (from graph) = 0.555

>>> get_osmotic_coefficient_pitzer(18*ureg.Quantity('mol/kg'),18*ureg.Quantity('mol/kg'),2,0,-0.01709,0.09198,0,0.000419,1,-1,1,1,b=1.2)
0.5556 ...

References

Scharge, T., Munoz, A.G., and Moog, H.C. (2012). Activity Coefficients of Fission Products in Highly Salinary Solutions of Na+, K+, Mg2+, Ca2+, Cl-, and SO42- : Cs+. /Journal of Chemical& Engineering Data (57), p. 1637-1647.

Kim, H., & Jr, W. F. (1988). Evaluation of Pitzer ion interaction parameters of aqueous electrolytes at 25 degree C. 1. Single salt parameters. Journal of Chemical and Engineering Data, (2), 177-184.

May, P. M., Rowland, D., Hefter, G., & Königsberger, E. (2011). A Generic and Updatable Pitzer Characterization of Aqueous Binary Electrolyte Solutions at 1 bar and 25 °C. Journal of Chemical & Engineering Data, 56(12), 5066-5077. doi:10.1021/je2009329

Beyer, R., & Steiger, M. (2010). Vapor Pressure Measurements of NaHCOO + H 2 O and KHCOO + H 2 O from 278 to 308 K and Representation with an Ion Interaction (Pitzer) Model. Journal of Chemical & Engineering Data, 55(2), 830-838. doi:10.1021/je900487a

pyEQL.engines module

pyEQL engines for computing aqueous equilibria (e.g., speciation, redox, etc.).

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

class pyEQL.engines.EOS[source]

Bases: ABC

Abstract base class for pyEQL equation of state classes.

The intent is that concrete implementations of this class make use of the standalone functions available in pyEQL.activity_correction and pyEQL.equilibrium as much as possible. This facilitates robust unit testing while allowing users to “mix and match” or customize the various models as needed.

abstract equilibrate(solution)[source]

Adjust the speciation and pH of a Solution object to achieve chemical equilibrium.

The Solution should be modified in-place, likely using add_moles / set_moles, etc.

Parameters:

solution – pyEQL Solution object

Returns

Nothing. The speciation of the Solution is modified in-place.

Raises:

ValueError if the calculation cannot be completed, e.g. due to insufficient number of parameters or lack of convergence.

abstract get_activity_coefficient(solution, solute)[source]

Return the molal scale activity coefficient of solute, given a Solution object.

Parameters:
  • solution – pyEQL Solution object

  • solute – str identifying the solute of interest

Returns

Quantity: dimensionless quantity object

Raises:

ValueError if the calculation cannot be completed, e.g. due to insufficient number of parameters.

abstract get_osmotic_coefficient(solution)[source]

Return the molal scale osmotic coefficient of a Solution.

Parameters:

solution – pyEQL Solution object

Returns

Quantity: dimensionless molal scale osmotic coefficient

Raises:

ValueError if the calculation cannot be completed, e.g. due to insufficient number of parameters.

abstract get_solute_volume()[source]

Return the volume of only the solutes.

Parameters:

solution – pyEQL Solution object

Returns

Quantity: solute volume in L

Raises:

ValueError if the calculation cannot be completed, e.g. due to insufficient number of parameters.

class pyEQL.engines.IdealEOS[source]

Bases: EOS

Ideal solution equation of state engine.

equilibrate(solution)[source]

Adjust the speciation of a Solution object to achieve chemical equilibrium.

get_activity_coefficient(solution, solute)[source]

Return the molal scale activity coefficient of solute, given a Solution object.

get_osmotic_coefficient(solution)[source]

Return the molal scale osmotic coefficient of solute, given a Solution object.

get_solute_volume(solution)[source]

Return the volume of the solutes.

class pyEQL.engines.NativeEOS(phreeqc_db: Literal['vitens.dat', 'wateq4f_PWN.dat', 'pitzer.dat', 'llnl.dat', 'geothermal.dat'] = 'llnl.dat')[source]

Bases: EOS

pyEQL’s native EOS. Uses the Pitzer model when possible, falls back to other models (e.g. Debye-Huckel) based on ionic strength if sufficient parameters are not available.

equilibrate(solution)[source]

Adjust the speciation of a Solution object to achieve chemical equilibrium.

get_activity_coefficient(solution, solute)[source]

Whenever the appropriate parameters are available, the Pitzer model [may] is used. If no Pitzer parameters are available, then the appropriate equations are selected according to the following logic: [stumm].

I <= 0.0005: Debye-Huckel equation 0.005 < I <= 0.1: Guntelberg approximation 0.1 < I <= 0.5: Davies equation I > 0.5: Raises a warning and returns activity coefficient = 1

The ionic strength, activity coefficients, and activities are all calculated based on the molal (mol/kg) concentration scale. If a different scale is given as input, then the molal-scale activity coefficient \(\gamma_\pm\) is converted according to [rbs]

\[f_\pm = \gamma_\pm * (1 + M_w \sum_i \nu_i m_i)\]
\[y_\pm = \frac{m \rho_w}{C \gamma_\pm}\]

where \(f_\pm\) is the rational activity coefficient, \(M_w\) is the molecular weight of water, the summation represents the total molality of all solute species, \(y_\pm\) is the molar activity coefficient, \(\rho_w\) is the density of pure water, \(m\) and \(C\) are the molal and molar concentrations of the chosen salt (not individual solute), respectively.

Parameters:
  • solute – String representing the name of the solute of interest

  • scale – The concentration scale for the returned activity coefficient. Valid options are “molal”, “molar”, and “rational” (i.e., mole fraction). By default, the molal scale activity coefficient is returned.

Returns:

The mean ion activity coefficient of the solute in question on the selected scale.

Notes

For multicomponent mixtures, pyEQL implements the “effective Pitzer model” presented by Mistry et al. [mistry]. In this model, the activity coefficient of a salt in a multicomponent mixture is calculated using an “effective molality,” which is the molality that would result in a single-salt mixture with the same total ionic strength as the multicomponent solution.

\[m_{effective} = \frac{2 I}{(\nu_{+} z_{+}^2 + \nu_{-}- z_{-}^2)}\]

References

[may]
May, P. M., Rowland, D., Hefter, G., & Königsberger, E. (2011).

A Generic and Updatable Pitzer Characterization of Aqueous Binary Electrolyte Solutions at 1 bar and 25 °C.

Journal of Chemical & Engineering Data, 56(12), 5066-5077. doi:10.1021/je2009329

[stumm]

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 165. Wiley Interscience, 1996.

[rbs]

Robinson, R. A.; Stokes, R. H. Electrolyte Solutions: Second Revised Edition; Butterworths: London, 1968, p.32.

[mistry]

Mistry, K. H.; Hunter, H. a.; Lienhard V, J. H. Effect of composition and nonideal solution behavior on desalination calculations for mixed electrolyte solutions with comparison to seawater. Desalination 2013, 318, 34-47.

get_osmotic_coefficient(solution)[source]

Return the molal scale osmotic coefficient of solute, given a Solution object.

Osmotic coefficient is calculated using the Pitzer model. [may] If appropriate parameters for the model are not available, then pyEQL raises a WARNING and returns an osmotic coefficient of 1.

If the ‘rational’ scale is given as input, then the molal-scale osmotic coefficient \(\phi\) is converted according to [rbs]

\[g = - \phi M_{w} \frac{\sum_{i} \nu_{i} m_{i}}{\ln x_{w}}\]

where \(g\) is the rational osmotic coefficient, \(M_{w}\) is the molecular weight of water, the summation represents the total molality of all solute species, and \(x_{w}\) is the mole fraction of water.

Parameters:

scale – The concentration scale for the returned osmotic coefficient. Valid options are “molal”, “rational” (i.e., mole fraction), and “fugacity”. By default, the molal scale osmotic coefficient is returned.

Returns:

The osmotic coefficient

Return type:

Quantity

Notes

For multicomponent mixtures, pyEQL adopts the “effective Pitzer model” presented by Mistry et al. [mstry]. In this approach, the osmotic coefficient of each individual salt is calculated using the normal Pitzer model based on its respective concentration. Then, an effective osmotic coefficient is calculated as the concentration-weighted average of the individual osmotic coefficients.

For example, in a mixture of 0.5 M NaCl and 0.5 M KBr, one would calculate the osmotic coefficient for each salt using a concentration of 0.5 M and an ionic strength of 1 M. Then, one would average the two resulting osmotic coefficients to obtain an effective osmotic coefficient for the mixture.

(Note: in the paper referenced below, the effective osmotic coefficient is determined by weighting using the “effective molality” rather than the true molality. Subsequent checking and correspondence with the author confirmed that the weight factor should be the true molality, and that is what is implemented in pyEQL.)

References

[may]

May, P. M., Rowland, D., Hefter, G., & Königsberger, E. (2011). A Generic and Updatable Pitzer Characterization of Aqueous Binary Electrolyte Solutions at 1 bar and 25 °C. Journal of Chemical & Engineering Data, 56(12), 5066-5077. doi:10.1021/je2009329

[rbs]

Robinson, R. A.; Stokes, R. H. Electrolyte Solutions: Second Revised Edition; Butterworths: London, 1968, p.32.

[mstry]

Mistry, K. H.; Hunter, H. a.; Lienhard V, J. H. Effect of composition and nonideal solution behavior on desalination calculations for mixed electrolyte solutions with comparison to seawater. Desalination 2013, 318, 34-47.

Examples

>>> s1 = pyEQL.Solution({'Na+': '0.2 mol/kg', 'Cl-': '0.2 mol/kg'})
>>> s1.get_osmotic_coefficient()
<Quantity(0.923715281, 'dimensionless')>
>>> s1 = pyEQL.Solution({'Mg+2': '0.3 mol/kg', 'Cl-': '0.6 mol/kg'},temperature='30 degC')
>>> s1.get_osmotic_coefficient()
<Quantity(0.891409618, 'dimensionless')>
get_solute_volume(solution)[source]

Return the volume of the solutes.

class pyEQL.engines.PhreeqcEOS(phreeqc_db: Literal['vitens.dat', 'wateq4f_PWN.dat', 'pitzer.dat', 'llnl.dat', 'geothermal.dat'] = 'phreeqc.dat')[source]

Bases: NativeEOS

Engine based on the PhreeqC model, as implemented via the phreeqpython package.

get_activity_coefficient(solution, solute)[source]

Return the molal scale activity coefficient of solute, given a Solution object.

get_osmotic_coefficient(solution)[source]

Return the molal scale osmotic coefficient of solute, given a Solution object.

PHREEQC appears to assume a unit osmotic coefficient unless the pitzer database is used. Unfortunately, there is no easy way to access the osmotic coefficient via phreeqcpython

get_solute_volume(solution)[source]

Return the volume of the solutes.

pyEQL.equilibrium module

pyEQL methods for chemical equilibrium calculations (e.g. acid/base, reactions, redox, complexation, etc.).

NOTE: these methods are not currently used but are here for the future.

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

pyEQL.equilibrium.adjust_temp_arrhenius(rate_constant, activation_energy, temperature, reference_temperature=<Quantity(25, 'degree_Celsius')>)[source]

Adjust a reaction equilibrium constant from one temperature to another.

Parameters:
  • rate_constant (Quantity) – The parameter value (usually a rate constant) being adjusted.

  • activation_energy (Quantity) – The activation energy of the process, in kJ/mol.

  • temperature (Quantity) – The desired reaction temperature.

  • reference_temperature (Quantity, optional) – The temperature at which equilibrium_constant is valid. Defaults to 25 degrees C if omitted.

Returns:

The adjusted reaction equilibrium constant.

Return type:

Quantity

Notes

This function implements the Arrhenius equation to adjust measured rate constants to other temperatures. TODO - add better reference

\[ln(\frac{K2}{K1}) = \frac{E_a}{R} ( \frac{1}{T_{1}} - {\frac{1}{T_2}} )\]

References

http://chemwiki.ucdavis.edu/Physical_Chemistry/Kinetics/Reaction_Rates/Temperature_Dependence_of_Reaction_Rates/Arrhenius_Equation

Examples

>>> adjust_temp_arrhenius(7,900*ureg.Quantity('kJ/mol'),37*ureg.Quantity('degC'),97*ureg.Quantity('degC')) 
1.8867225...e-24
pyEQL.equilibrium.adjust_temp_pitzer(c1, c2, c3, c4, c5, temp, temp_ref=<Quantity(298.15, 'kelvin')>)[source]

Calculate a parameter for the Pitzer model based on temperature-dependent coefficients c1,c2,c3,c4,and c5.

Parameters:
  • c1 (float) – Temperature-dependent coefficients for the pitzer parameter of interest.

  • c2 (float) – Temperature-dependent coefficients for the pitzer parameter of interest.

  • c3 (float) – Temperature-dependent coefficients for the pitzer parameter of interest.

  • c4 (float) – Temperature-dependent coefficients for the pitzer parameter of interest.

  • c5 (float) – Temperature-dependent coefficients for the pitzer parameter of interest.

  • temp (Quantity) – The temperature at which the Pitzer parameter is to be calculated.

  • temp_ref (Quantity, optional) – The reference temperature on which the parameters are based. Defaults to 298.15 K if omitted.

Note

As described in the PHREEQC documentation.

pyEQL.equilibrium.adjust_temp_vanthoff(equilibrium_constant, enthalpy, temperature, reference_temperature=<Quantity(25, 'degree_Celsius')>)[source]

Adjust a reaction equilibrium constant from one temperature to another.

Parameters:
  • equilibrium_constant (float) – The reaction equilibrium constant for the reaction.

  • enthalpy (Quantity) – The enthalpy change (delta H) for the reaction in kJ/mol. Assumed independent of temperature (see Notes).

  • temperature (Quantity) – The desired reaction temperature in degrees Celsius.

  • reference_temperature (Quantity, optional) – The temperature at which equilibrium_constant is valid. Defaults to 25 degrees C if omitted.

Returns:

The adjusted reaction equilibrium constant.

Return type:

float

Note

This function implements the Van’t Hoff equation to adjust measured equilibrium constants to other temperatures.

\[ln(K2 / K1) = {\delta H \over R} ( {1 \over T_1} - {1 \over T_2} )\]

This implementation assumes that the enthalpy is independent of temperature over the range of interest.

References

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 53. Wiley Interscience, 1996.

Examples

>>> adjust_temp_vanthoff(0.15,ureg.Quantity('-197.6 kJ/mol'),ureg.Quantity('42 degC'),ureg.Quantity(' 25degC')) 
0.00203566...

If the ‘ref_temperature’ parameter is omitted, a default of 25 C is used.

>>> adjust_temp_vanthoff(0.15,ureg.Quantity('-197.6 kJ/mol'),ureg.Quantity('42 degC')) 
0.00203566...
pyEQL.equilibrium.alpha(n, pH, pKa_list)[source]

Returns the acid-base distribution coefficient (alpha) of an acid in the n-deprotonated form at a given pH.

Parameters:
  • n (int) – The number of protons that have been lost by the desired form of the acid. Also the subscript on the alpha value. E.g. for bicarbonate (HCO3-), n=1 because 1 proton has been lost from the fully-protonated carbonic acid (H2CO3) form.

  • pH (float or int) – The pH of the solution.

  • pKa_list (list of floats or ints) – The pKa values (negative log of equilibrium constants) for the acid of interest. There must be a minimum of n pKa values in the list.

Returns:

The fraction of total acid present in the specified form.

Return type:

float

Notes

The acid-base coefficient is calculated as follows: [stm]

\[\alpha_n = \frac{term_n}{[H+]^n + k_{a1}[H+]^{n-1} + k_{a1}k_{a2}[H+]^{n-2} ... k_{a1}k_{a2}...k_{an} }\]

Where :math: ‘term_n’ refers to the nth term in the denominator, starting from 0

References

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 127-130. Wiley Interscience, 1996.

Examples

>>> alpha(1,8,[4.7]) 
0.999...

The sum of all alpha values should equal 1

>>> alpha(0,8,[6.35,10.33]) 
0.021...
>>> alpha(1,8,[6.35,10.33]) 
0.979...
>>> alpha(2,8,[6.35,10.33]) 
2.043...e-09

If pH is equal to one of the pKa values the function should return 0.5.

>>> alpha(1,6.35,[6.35,10.33])
0.5

pyEQL.functions module

pyEQL functions that take Solution objects as inputs or return Solution objects.

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

pyEQL.functions.donnan_eql(solution: Solution, fixed_charge: str)[source]

Return a solution object in equilibrium with fixed_charge.

Parameters:
  • solution – Solution object The external solution to be brought into equilibrium with the fixed charges

  • fixed_charge – str quantity String representing the concentration of fixed charges, including sign. May be specified in mol/L or mol/kg units. e.g. ‘1 mol/kg’

Returns:

A Solution that has established Donnan equilibrium with the external (input) Solution

Notes

The general equation representing the equilibrium between an external electrolyte solution and an ion-exchange medium containing fixed charges is

\[\big(\frac{a_{-}}{\bar a_{-}} \big)^(\frac{1}{z_{-}) \big(\frac{\bar a_{+}}{a_{+}}\big)^(\frac{1}{z_{+}) \exp(\frac{\Delta \pi \bar V}{RT z_{+} \nu_{+}})\]

Where subscripts \(+\) and \(-\) indicate the cation and anion, respectively, the overbar indicates the membrane phase, \(a\) represents activity, \(z\) represents charge, \(\nu\) represents the stoichiometric coefficient, \(V\) represents the partial molar volume of the salt, and \(\Delta \pi\) is the difference in osmotic pressure between the membrane and the solution phase.

In addition, electroneutrality must prevail within the membrane phase:

\[\bar C_{+} z_{+} + \bar X + \bar C_{-} z_{-} = 0\]

Where \(C\) represents concentration and \(X\) is the fixed charge concentration in the membrane or ion exchange phase.

This function solves these two equations simultaneously to arrive at the concentrations of the cation and anion in the membrane phase. It returns a solution equal to the input solution except that the concentrations of the predominant cation and anion have been adjusted according to this equilibrium.

NOTE that this treatment is only capable of equilibrating a single salt. This salt is identified by the get_salt() method.

References

Strathmann, Heiner, ed. Membrane Science and Technology vol. 9, 2004. Chapter 2, p. 51.

http://dx.doi.org/10.1016/S0927-5193(04)80033-0

See also

get_salt()

pyEQL.functions.entropy_mix(solution1: Solution, solution2: Solution)[source]

Return the ideal mixing entropy associated with mixing two solutions.

Parameters:
  • solution1 – The two solutions to be mixed.

  • solution2 – The two solutions to be mixed.

Returns:

The ideal mixing entropy associated with complete mixing of the Solutions, in Joules.

Notes

The ideal entropy of mixing is calculated as follows

\[\Delta_{mix} S = \sum_i {(n_c + n_d) R T \ln x_b} - \sum_i {n_c R T \ln x_c} - \sum_i {n_d R T \ln x_d}\]

Where \(n\) is the number of moles of substance, \(T\) is the temperature in kelvin, and subscripts \(b\), \(c\), and \(d\) refer to the concentrated, dilute, and blended Solutions, respectively.

Note that dissociated ions must be counted as separate components, so a simple salt dissolved in water is a three component solution (cation, anion, and water).

References

Koga, Yoshikata, 2007. *Solution Thermodynamics and its Application to Aqueous Solutions:

A differential approach.* Elsevier, 2007, pp. 23-37.

pyEQL.functions.gibbs_mix(solution1: Solution, solution2: Solution)[source]

Return the Gibbs energy change associated with mixing two solutions.

Parameters:
  • solution1 – a solution to be mixed.

  • solution2 – a solution to be mixed.

Returns:

The change in Gibbs energy associated with complete mixing of the Solutions, in Joules.

Notes

The Gibbs energy of mixing is calculated as follows

\[\Delta_{mix} G = \sum_i {(n_c + n_d) R T \ln a_b} - \sum_i {n_c R T \ln a_c} - \sum_i {n_d R T \ln a_d}\]

Where \(n\) is the number of moles of substance, \(T\) is the temperature in kelvin, and subscripts \(b\), \(c\), and \(d\) refer to the concentrated, dilute, and blended Solutions, respectively.

Note that dissociated ions must be counted as separate components, so a simple salt dissolved in water is a three component solution (cation, anion, and water).

References

Koga, Yoshikata, 2007. Solution Thermodynamics and its Application to Aqueous Solutions:

A differential approach. Elsevier, 2007, pp. 23-37.

pyEQL.salt_ion_match module

pyEQL salt matching library.

This file contains functions that allow a pyEQL Solution object composed of individual species (usually ions) to be mapped to a solution of one or more salts. This mapping is necessary because some parameters (such as activity coefficient data) can only be determined for salts (e.g. NaCl) and not individual species (e.g. Na+)

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

class pyEQL.salt_ion_match.Salt(cation, anion)[source]

Bases: MSONable

Class to represent a salt.

get_effective_molality(ionic_strength)[source]

Calculate the effective molality according to [mistry].

\[2 I \over (\nu_+ z_+^2 + \nu_- z_- ^2)\]
Parameters:

ionic_strength – Quantity The ionic strength of the parent solution, mol/kg

Returns:

the effective molality of the salt in the parent solution

Return type:

Quantity

References

[mistry]

Mistry, K. H.; Hunter, H. a.; Lienhard V, J. H. Effect of composition and nonideal solution behavior on desalination calculations for mixed electrolyte solutions with comparison to seawater. Desalination 2013, 318, 34-47.

pyEQL.solute module

pyEQL Solute class.

This file contains functions and methods for managing properties of individual solutes. The Solute class contains methods for accessing ONLY those properties that DO NOT depend on solution composition. Solute properties such as activity coefficient or concentration that do depend on compsition are accessed via Solution class methods.

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

class pyEQL.solute.Datum(value: str, reference: str | None = None, data_type: Literal['computed', 'experimental', 'fitted', 'unknown'] = 'unknown')[source]

Bases: object

Document containing data for a single computed or experimental property.

as_dict()[source]

Return a dictionary representation of the Datum.

data_type: Literal['computed', 'experimental', 'fitted', 'unknown'] = 'unknown'
property magnitude

Return the numerical value of a Datum.

reference: str | None = None
property uncertainty

Return the uncertainty of a Datum.

property unit

Return the unit of a Datum.

value: str
class pyEQL.solute.Solute(formula: str, charge: int, molecular_weight: str, elements: list, chemsys: str, pmg_ion: ~pymatgen.core.ion.Ion, formula_html: str, formula_latex: str, formula_hill: str, formula_pretty: str, oxi_state_guesses: dict[str, float], n_atoms: int, n_elements: int, size: dict = <factory>, thermo: dict = <factory>, transport: dict = <factory>, model_parameters: dict = <factory>)[source]

Bases: object

represent each chemical species as an object containing its formal charge, transport numbers, concentration, activity, etc.

Parameters:

formula – Chemical formula for the solute. Charged species must contain a + or - and (for polyvalent solutes) a number representing the net charge (e.g. ‘SO4-2’).

as_dict()[source]

Return a dictionary representation of the Solute.

charge: int
chemsys: str
elements: list
formula: str
formula_hill: str
formula_html: str
formula_latex: str
formula_pretty: str
classmethod from_formula(formula: str)[source]

Create an Ion document from a chemical formula. The formula is passed to pymatgen.core.Ion.from_formula() and used to populate the basic chemical informatics fields (e.g., formula, charge, molecular weight, elements, etc.) of the IonDoc.

model_parameters: dict
molecular_weight: str
n_atoms: int
n_elements: int
oxi_state_guesses: dict[str, float]
pmg_ion: Ion
size: dict
thermo: dict
transport: dict

pyEQL.solution module

pyEQL Solution Class.

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

class pyEQL.solution.Solution(solutes: list[list[str]] | dict[str, str] | None = None, volume: str | None = None, temperature: str = '298.15 K', pressure: str = '1 atm', pH: float = 7, pE: float = 8.5, balance_charge: str | None = None, solvent: str | list = 'H2O', engine: EOS | Literal['native', 'ideal', 'phreeqc'] = 'native', database: str | Path | Store | None = None, default_diffusion_coeff: float = 1.6106e-09, log_level: Literal['DEBUG', 'INFO', 'WARNING', 'ERROR', 'CRITICAL'] | None = 'ERROR')[source]

Bases: MSONable

Class representing the properties of a solution. Instances of this class contain information about the solutes, solvent, and bulk properties.

property TDS: Quantity

Alias of total_dissolved_solids().

add_amount(solute: str, amount: str)[source]

Add the amount of ‘solute’ to the parent solution.

Parameters:
  • solute – str String representing the name of the solute of interest

  • amount – str quantity String representing the concentration desired, e.g. ‘1 mol/kg’ If the units are given on a per-volume basis, the solution volume is not recalculated If the units are given on a mass, substance, per-mass, or per-substance basis, then the solution volume is recalculated based on the new composition

Returns:

Nothing. The concentration of solute is modified.

add_solute(formula: str, amount: str)[source]

Primary method for adding substances to a pyEQL solution.

Parameters:
  • formula (str) – Chemical formula for the solute. Charged species must contain a + or - and

  • charge ((for polyvalent solutes) a number representing the net) –

  • amount (str) – The amount of substance in the specified unit system. The string should contain

  • g/L'. (both a quantity and a pint-compatible representation of a ureg. e.g. '5 mol/kg' or '0.1) –

add_solvent(formula: str, amount: str)[source]

Same as add_solute but omits the need to pass solvent mass to pint.

property alkalinity: Quantity

Return the alkalinity or acid neutralizing capacity of a solution.

Returns:

The alkalinity of the solution in mg/L as CaCO3

Return type:

Quantity

Notes

The alkalinity is calculated according to [stm]

\[Alk = \sum_{i} z_{i} C_{B} + \sum_{i} z_{i} C_{A}\]

Where \(C_{B}\) and \(C_{A}\) are conservative cations and anions, respectively (i.e. ions that do not participate in acid-base reactions), and \(z_{i}\) is their signed charge. In this method, the set of conservative cations is all Group I and Group II cations, and the conservative anions are all the anions of strong acids.

References

[stm] (1,2)

Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 165. Wiley Interscience, 1996.

property anions: dict[str, float]

Returns the subset of components that are anions.

The returned dict is sorted by amount in descending order.

as_dict() dict[source]

Convert the Solution into a dict representation that can be serialized to .json or other format.

balance_charge

Standardized formula of the species used for charge balancing.

property bjerrum_length: Quantity

Return the Bjerrum length of a solution.

Bjerrum length represents the distance at which electrostatic interactions between particles become comparable in magnitude to the thermal energy.:math:lambda_B is calculated as

\[\lambda_B = {e^2 \over (4 \pi \epsilon_r \epsilon_o k_B T)}\]

where \(e\) is the fundamental charge, \(\epsilon_r\) and \(\epsilon_r\) are the relative permittivity and vacuum permittivity, \(k_B\) is the Boltzmann constant, and \(T\) is the temperature.

Returns:

The Bjerrum length, in nanometers.

Return type:

Quantity

References

https://en.wikipedia.org/wiki/Bjerrum_length

Examples

>>> s1 = pyEQL.Solution()
>>> s1.bjerrum_length
<Quantity(0.7152793009386953, 'nanometer')>
property cations: dict[str, float]

Returns the subset of components that are cations.

The returned dict is sorted by amount in descending order.

property charge_balance: float

Return the charge balance of the solution.

Return the charge balance of the solution. The charge balance represents the net electric charge on the solution and SHOULD equal zero at all times, but due to numerical errors will usually have a small nonzero value. It is calculated according to:

\[CB = \sum_i C_i z_i\]

where \(C_i\) is the molar concentration, and \(z_i\) is the charge on species i.

Returns:

The charge balance of the solution, in equivalents (mol of charge) per L.

Return type:

float

property chemical_system: str

Return the chemical system of the Solution as a “-” separated list of elements, sorted alphabetically. For example, a solution containing CaCO3 would have a chemical system of “C-Ca-H-O”.

components

Special dictionary where keys are standardized formula and values are the moles present in Solution.

property conductivity: Quantity

Compute the electrical conductivity of the solution.

Returns:

The electrical conductivity of the solution in Siemens / meter.

Notes

Conductivity is calculated by summing the molar conductivities of the respective solutes.

\[EC = {F^2 \over R T} \sum_i D_i z_i ^ 2 m_i = \sum_i \lambda_i m_i\]

Where \(D_i\) is the diffusion coefficient, \(m_i\) is the molal concentration, \(z_i\) is the charge, and the summation extends over all species in the solution. Alternatively, \(\lambda_i\) is the molar conductivity of solute i.

Diffusion coefficients \(D_i\) (and molar conductivities \(\lambda_i\)) are adjusted for the effects of temperature and ionic strength using the method implemented in PHREEQC >= 3.4. [aq] [hc] See get_diffusion_coefficient for further details.

References

See also

ionic_strength get_diffusion_coefficient() get_molar_conductivity()

database

Store instance containing the solute property database.

property debye_length: Quantity

Return the Debye length of a solution.

Debye length is calculated as [wk3]

\[\kappa^{-1} = \sqrt({\epsilon_r \epsilon_o k_B T \over (2 N_A e^2 I)})\]

where \(I\) is the ionic strength, \(\epsilon_r\) and \(\epsilon_r\) are the relative permittivity and vacuum permittivity, \(k_B\) is the Boltzmann constant, and \(T\) is the temperature, \(e\) is the elementary charge, and \(N_A\) is Avogadro’s number.

Returns The Debye length, in nanometers.

References

property density: Quantity

Return the density of the solution.

Density is calculated from the mass and volume each time this method is called.

Returns:

The density of the solution.

Return type:

Quantity

property dielectric_constant: Quantity

Returns the dielectric constant of the solution.

Parameters:

None

Returns:

the dielectric constant of the solution, dimensionless.

Return type:

Quantity

Notes

Implements the following equation as given by Zuber et al.

\[\epsilon = \epsilon_{solvent} \over 1 + \sum_i \alpha_i x_i\]

where \(\alpha_i\) is a coefficient specific to the solvent and ion, and \(x_i\) is the mole fraction of the ion in solution.

References

A. Zuber, L. Cardozo-Filho, V.F. Cabral, R.F. Checoni, M. Castier, An empirical equation for the dielectric constant in aqueous and nonaqueous electrolyte mixtures, Fluid Phase Equilib. 376 (2014) 116-123. doi:10.1016/j.fluid.2014.05.037.

property elements: list

Return a list of elements that are present in the solution.

For example, a solution containing CaCO3 would return [“C”, “Ca”, “H”, “O”]

equilibrate(**kwargs) None[source]

Update the composition of the Solution using the thermodynamic engine.

Any kwargs specified are passed through to self.engine.equilibrate()

Returns:

Nothing. The .components attribute of the Solution is updated.

classmethod from_dict(d: dict) Solution[source]

Instantiate a Solution from a dictionary generated by as_dict().

classmethod from_file(filename: str | Path) Solution[source]

Loading from a .yaml or .json file.

Parameters:

filename (str | Path) – Path to the .json or .yaml file (including extension) to load the Solution from. Valid extensions are .json or .yaml.

Returns:

A pyEQL Solution object.

Raises:

FileNotFoundError – If the given filename doesn’t exist on the file system.

classmethod from_preset(preset: Literal['seawater', 'rainwater', 'wastewater', 'urine', 'normal saline', 'Ringers lactate']) Solution[source]

Instantiate a solution from a preset composition.

Parameters:

preset (str) – String representing the desired solution. Valid entries are ‘seawater’, ‘rainwater’, ‘wastewater’, ‘urine’, ‘normal saline’ and ‘Ringers lactate’.

Returns:

A pyEQL Solution object.

Raises:

FileNotFoundError – If the given preset file doesn’t exist on the file system.

Notes

The following sections explain the different solution options:

  • ‘rainwater’ - pure water in equilibrium with atmospheric CO2 at pH 6

  • ‘seawater’ or ‘SW’- Standard Seawater. See Table 4 of the Reference for Composition [1]

  • ‘wastewater’ or ‘WW’ - medium strength domestic wastewater. See Table 3-18 of [2]

  • ‘urine’ - typical human urine. See Table 3-15 of [2]

  • ‘normal saline’ or ‘NS’ - normal saline solution used in medicine [3]

  • ‘Ringers lacatate’ or ‘RL’ - Ringer’s lactate solution used in medicine [4]

References

get_activity(solute: str, scale: Literal['molal', 'molar', 'rational'] = 'molal') Quantity[source]

Return the thermodynamic activity of the solute in solution on the chosen concentration scale.

Parameters:
  • solute – String representing the name of the solute of interest

  • scale – The concentration scale for the returned activity. Valid options are “molal”, “molar”, and “rational” (i.e., mole fraction). By default, the molal scale activity is returned.

  • verbose – If True, pyEQL will print a message indicating the parent salt that is being used for activity calculations. This option is useful when modeling multicomponent solutions. False by default.

Returns:

The thermodynamic activity of the solute in question (dimensionless Quantity)

Notes

The thermodynamic activity depends on the concentration scale used [rs] . By default, the ionic strength, activity coefficients, and activities are all calculated based on the molal (mol/kg) concentration scale.

References

[rs]

Robinson, R. A.; Stokes, R. H. Electrolyte Solutions: Second Revised Edition; Butterworths: London, 1968, p.32.

get_activity_coefficient(solute: str, scale: Literal['molal', 'molar', 'fugacity', 'rational'] = 'molal') Quantity[source]

Return the activity coefficient of a solute in solution.

The model used to calculate the activity coefficient is determined by the Solution’s equation of state engine.

Parameters:
  • solute – The solute for which to retrieve the activity coefficient

  • scale – The activity coefficient concentration scale

  • verbose – If True, pyEQL will print a message indicating the parent salt that is being used for activity calculations. This option is useful when modeling multicomponent solutions. False by default.

Returns:

the activity coefficient as a dimensionless pint Quantity

Return type:

Quantity

get_amount(solute: str, units: str = 'mol/L') Quantity[source]

Return the amount of ‘solute’ in the parent solution.

The amount of a solute can be given in a variety of unit types. 1. substance per volume (e.g., ‘mol/L’, ‘M’) 2. equivalents (i.e., moles of charge) per volume (e.g., ‘eq/L’, ‘meq/L’) 3. substance per mass of solvent (e.g., ‘mol/kg’, ‘m’) 4. mass of substance (e.g., ‘kg’) 5. moles of substance (‘mol’) 6. mole fraction (‘fraction’) 7. percent by weight (%) 8. number of molecules (‘count’) 9. “parts-per-x” units, where ppm = mg/L, ppb = ug/L ppt = ng/L

Parameters:
  • solute – str String representing the name of the solute of interest

  • units – str Units desired for the output. Examples of valid units are ‘mol/L’,’mol/kg’,’mol’, ‘kg’, and ‘g/L’ Use ‘fraction’ to return the mole fraction. Use ‘%’ to return the mass percent

Returns:

The amount of the solute in question, in the specified units

get_chemical_potential_energy(activity_correction: bool = True) Quantity[source]

Return the total chemical potential energy of a solution (not including pressure or electric effects).

Parameters:

activity_correction – bool, optional If True, activities will be used to calculate the true chemical potential. If False, mole fraction will be used, resulting in a calculation of the ideal chemical potential.

Returns:

Quantity

The actual or ideal chemical potential energy of the solution, in Joules.

Notes

The chemical potential energy (related to the Gibbs mixing energy) is calculated as follows: [koga]

\[E = R T \sum_i n_i \ln a_i\]

or

\[E = R T \sum_i n_i \ln x_i\]

Where \(n\) is the number of moles of substance, \(T\) is the temperature in kelvin, \(R\) the ideal gas constant, \(x\) the mole fraction, and \(a\) the activity of each component.

Note that dissociated ions must be counted as separate components, so a simple salt dissolved in water is a three component solution (cation, anion, and water).

References

[koga]

Koga, Yoshikata, 2007. *Solution Thermodynamics and its Application to Aqueous Solutions:

A differential approach.* Elsevier, 2007, pp. 23-37.

get_components_by_element() dict[str, list][source]

Return a list of all species associated with a given element.

Elements (keys) are suffixed with their oxidation state in parentheses, e.g.,

{“Na(1.0)”:[“Na[+1]”, “NaOH(aq)”]}

Species associated with each element are sorted in descending order of the amount present (i.e., the first species listed is the most abundant).

get_el_amt_dict()[source]

Return a dict of Element: amount in mol.

Elements (keys) are suffixed with their oxidation state in parentheses, e.g. “Fe(2.0)”, “Cl(-1.0)”.

get_lattice_distance(solute: str) Quantity[source]

Calculate the average distance between molecules.

Calculate the average distance between molecules of the given solute, assuming that the molecules are uniformly distributed throughout the solution.

Parameters:

solute – str String representing the name of the solute of interest

Returns:

The average distance between solute molecules

Return type:

Quantity

Examples

>>> soln = Solution([['Na+','0.5 mol/kg'],['Cl-','0.5 mol/kg']])
>>> soln.get_lattice_distance('Na+')
1.492964.... nanometer

Notes

The lattice distance is related to the molar concentration as follows:

\[d = ( C_i N_A ) ^ {-{1 \over 3}}\]
get_moles_solvent() Quantity[source]

Return the moles of solvent present in the solution.

Returns:

The moles of solvent in the solution.

get_osmolality(activity_correction=False) Quantity[source]

Return the osmolality of the solution in Osm/kg.

Parameters:

activity_correction – bool If TRUE, the osmotic coefficient is used to calculate the osmolarity. This correction is appropriate when trying to predict the osmolarity that would be measured from e.g. freezing point depression. Defaults to FALSE if omitted.

get_osmolarity(activity_correction=False) Quantity[source]

Return the osmolarity of the solution in Osm/L.

Parameters:

activity_correction – bool If TRUE, the osmotic coefficient is used to calculate the osmolarity. This correction is appropriate when trying to predict the osmolarity that would be measured from e.g. freezing point depression. Defaults to FALSE if omitted.

get_osmotic_coefficient(scale: Literal['molal', 'molar', 'rational'] = 'molal') Quantity[source]

Return the osmotic coefficient of an aqueous solution.

The method used depends on the Solution object’s equation of state engine.

get_salt() Salt[source]

Determine the predominant salt in a solution of ions.

Many empirical equations for solution properties such as activity coefficient, partial molar volume, or viscosity are based on the concentration of single salts (e.g., NaCl). When multiple ions are present (e.g., a solution containing Na+, Cl-, and Mg+2), it is generally not possible to directly model these quantities. pyEQL works around this problem by treating such solutions as single salt solutions.

The get_salt() method examines the ionic composition of a solution and returns an object that identifies the single most predominant salt in the solution, defined by the cation and anion with the highest mole fraction. The Salt object contains information about the stoichiometry of the salt to enable its effective concentration to be calculated (e.g., if a solution contains 0.5 mol/kg of Na+ and Cl-, plus traces of H+ and OH-, the matched salt is 0.5 mol/kg NaCl).

Returns:

Salt object containing information about the parent salt.

Examples

>>> s1 = Solution([['Na+','0.5 mol/kg'],['Cl-','0.5 mol/kg']])
>>> s1.get_salt()
<pyEQL.salt_ion_match.Salt object at 0x7fe6d3542048>
>>> s1.get_salt().formula
'NaCl'
>>> s1.get_salt().nu_cation
1
>>> s1.get_salt().z_anion
-1
>>> s2 = pyEQL.Solution([['Na+','0.1 mol/kg'],['Mg+2','0.2 mol/kg'],['Cl-','0.5 mol/kg']])
>>> s2.get_salt().formula
'MgCl2'
>>> s2.get_salt().nu_anion
2
>>> s2.get_salt().z_cation
2
get_salt_dict(cutoff: float = 0.01, use_totals: bool = True) dict[str, dict][source]

Returns a dict of salts that approximates the composition of the Solution. Like components, the dict is keyed by formula and the values are the total moles present in the solution, e.g., {“NaCl(aq)”: 1}. If the Solution is pure water, the returned dict contains only ‘HOH’.

Parameters:
  • cutoff – Lowest salt concentration to consider. Analysis will stop once the concentrations of Salts being analyzed goes below this value. Useful for excluding analysis of trace anions.

  • use_totals – Whether to base the analysis on total element concentrations or individual species concentrations.

Notes

Salts are identified by pairing the predominant cations and anions in the solution, in descending order of their respective equivalent amounts.

Many empirical equations for solution properties such as activity coefficient, partial molar volume, or viscosity are based on the concentration of single salts (e.g., NaCl). When multiple ions are present (e.g., a solution containing Na+, Cl-, and Mg+2), it is generally not possible to directly model these quantities.

The get_salt_dict() method examines the ionic composition of a solution and simplifies it into a list of salts. The method returns a dictionary of Salt objects where the keys are the salt formulas (e.g., ‘NaCl’). The Salt object contains information about the stoichiometry of the salt to enable its effective concentration to be calculated (e.g., 1 M MgCl2 yields 1 M Mg+2 and 2 M Cl-).

Returns:

dict

A dictionary of Salt objects, keyed to the salt formula

get_total_amount(element: str, units: str) Quantity[source]

Return the total amount of ‘element’ (across all solutes) in the solution.

Parameters:
  • element – The symbol of the element of interest. The symbol can optionally be followed by the oxidation state in parentheses, e.g., “Na(1.0)”, “Fe(2.0)”, or “O(0.0)”. If no oxidation state is given, the total concentration of the element (over all oxidation states) is returned.

  • units – str Units desired for the output. Any unit understood by get_amount can be used. Examples of valid units are ‘mol/L’,’mol/kg’,’mol’, ‘kg’, and ‘g/L’.

Returns:

The total amount of the element in the solution, in the specified units

get_total_moles_solute() Quantity[source]

Return the total moles of all solute in the solution.

get_transport_number(solute: str) Quantity[source]

Calculate the transport number of the solute in the solution.

Parameters:

solute – Formula of the solute for which the transport number is to be calculated.

Returns:

The transport number of solute, as a dimensionless Quantity.

Notes

Transport number is calculated according to :

\[t_i = {D_i z_i^2 C_i \over \sum D_i z_i^2 C_i}\]

Where \(C_i\) is the concentration in mol/L, \(D_i\) is the diffusion coefficient, and \(z_i\) is the charge, and the summation extends over all species in the solution.

Diffusion coefficients \(D_i\) are adjusted for the effects of temperature and ionic strength using the method implemented in PHREEQC >= 3.4. See get_diffusion_coefficient for further details.

References

Geise, G. M.; Cassady, H. J.; Paul, D. R.; Logan, E.; Hickner, M. A. “Specific ion effects on membrane potential and the permselectivity of ion exchange membranes.”” Phys. Chem. Chem. Phys. 2014, 16, 21673-21681.

See also

get_diffusion_coefficient() get_molar_conductivity()

get_water_activity() Quantity[source]

Return the water activity.

Returns:

The thermodynamic activity of water in the solution.

Return type:

Quantity

Notes

Water activity is related to the osmotic coefficient in a solution containing i solutes by:

\[\ln a_{w} = - \Phi M_{w} \sum_{i} m_{i}\]

Where \(M_{w}\) is the molar mass of water (0.018015 kg/mol) and \(m_{i}\) is the molal concentration of each species.

If appropriate Pitzer model parameters are not available, the water activity is assumed equal to the mole fraction of water.

References

Blandamer, Mike J., Engberts, Jan B. F. N., Gleeson, Peter T., Reis, Joao Carlos R., 2005. “Activity of water in aqueous systems: A frequently neglected property.” Chemical Society Review 34, 440-458.

Examples

>>> s1 = pyEQL.Solution([['Na+','0.3 mol/kg'],['Cl-','0.3 mol/kg']])
>>> s1.get_water_activity()
<Quantity(0.9900944932888518, 'dimensionless')>
property hardness: Quantity

Return the hardness of a solution.

Hardness is defined as the sum of the equivalent concentrations of multivalent cations as calcium carbonate.

NOTE: at present pyEQL cannot distinguish between mg/L as CaCO3 and mg/L units. Use with caution.

Returns:

The hardness of the solution in mg/L as CaCO3

Return type:

Quantity

property ionic_strength: Quantity

Return the ionic strength of the solution.

Return the ionic strength of the solution, calculated as 1/2 * sum ( molality * charge ^2) over all the ions.

Molal (mol/kg) scale concentrations are used for compatibility with the activity correction formulas.

Returns:

The ionic strength of the parent solution, mol/kg.

Return type:

Quantity

Notes

The ionic strength is calculated according to:

\[I = \sum_i m_i z_i^2\]

Where \(m_i\) is the molal concentration and \(z_i\) is the charge on species i.

Examples

>>> s1 = pyEQL.Solution([['Na+','0.2 mol/kg'],['Cl-','0.2 mol/kg']])
>>> s1.ionic_strength
<Quantity(0.20000010029672785, 'mole / kilogram')>
>>> s1 = pyEQL.Solution([['Mg+2','0.3 mol/kg'],['Na+','0.1 mol/kg'],['Cl-','0.7 mol/kg']],temperature='30 degC')
>>> s1.ionic_strength
<Quantity(1.0000001004383303, 'mole / kilogram')>
list_activities(**kwargs)
list_concentrations(**kwargs)
list_salts(**kwargs)
list_solutes(**kwargs)
property mass: Quantity

Return the total mass of the solution.

The mass is calculated each time this method is called.

Returns: The mass of the solution, in kg

property neutrals: dict[str, float]

Returns the subset of components that are neutral (not charged).

The returned dict is sorted by amount in descending order.

property osmotic_pressure: Quantity

Return the osmotic pressure of the solution relative to pure water.

Returns:

The osmotic pressure of the solution relative to pure water in Pa

Notes

Osmotic pressure is calculated based on the water activity [sata] [wk]

\[\Pi = -\frac{RT}{V_{w}} \ln a_{w}\]

Where \(\Pi\) is the osmotic pressure, \(V_{w}\) is the partial molar volume of water (18.2 cm**3/mol), and \(a_{w}\) is the water activity.

References

[sata]

Sata, Toshikatsu. Ion Exchange Membranes: Preparation, Characterization, and Modification. Royal Society of Chemistry, 2004, p. 10.

Examples

>>> s1=pyEQL.Solution()
>>> s1.osmotic_pressure
<Quantity(0.495791416, 'pascal')>
>>> s1 = pyEQL.Solution([['Na+','0.2 mol/kg'],['Cl-','0.2 mol/kg']])
>>> soln.osmotic_pressure
<Quantity(906516.7318131207, 'pascal')>
p(solute: str, activity=True) float | None[source]

Return the negative log of the activity of solute.

Generally used for expressing concentration of hydrogen ions (pH)

Parameters:
  • solute – str String representing the formula of the solute

  • activity – bool, optional If False, the function will use the molar concentration rather than the activity to calculate p. Defaults to True.

Returns:

Quantity

The negative log10 of the activity (or molar concentration if activity = False) of the solute.

property pH: float | None

Return the pH of the solution.

property pressure: Quantity

Return the hydrostatic pressure of the solution in atm.

print(mode: Literal['all', 'ions', 'cations', 'anions', 'neutrals'] = 'all', units: Literal['ppm', 'mol', 'mol/kg', 'mol/L', '%', 'activity'] = 'mol', places=4)[source]

Print details about the Solution.

Parameters:
  • mode – Whether to list the amounts of all solutes, or only anions, cations, any ion, or any neutral solute.

  • units – The units to list solute amounts in. “activity” will list dimensionless activities instead of concentrations.

  • places – The number of decimal places to round the solute amounts.

set_amount(solute: str, amount: str)[source]

Set the amount of ‘solute’ in the parent solution.

Parameters:
  • solute – str String representing the name of the solute of interest

  • amount

    str Quantity String representing the concentration desired, e.g. ‘1 mol/kg’ If the units are given on a per-volume basis, the solution volume is not recalculated and the molar concentrations of other components in the solution are not altered, while the molal concentrations are modified.

    If the units are given on a mass, substance, per-mass, or per-substance basis, then the solution volume is recalculated based on the new composition and the molal concentrations of other components are not altered, while the molar concentrations are modified.

Returns:

Nothing. The concentration of solute is modified.

solvent

Formula of the component that is set as the solvent (currently only H2O(aq) is supported).

property solvent_mass: Quantity

Return the mass of the solvent.

This property is used whenever mol/kg (or similar) concentrations are requested by get_amount()

Returns:

The mass of the solvent, in kg

See also

get_amount()

property temperature: Quantity

Return the temperature of the solution in Kelvin.

to_file(filename: str | Path) None[source]

Saving to a .yaml or .json file.

Parameters:

filename (str | Path) – The path to the file to save Solution. Valid extensions are .json or .yaml.

property total_dissolved_solids: Quantity

Total dissolved solids in mg/L (equivalent to ppm) including both charged and uncharged species.

The TDS is defined as the sum of the concentrations of all aqueous solutes (not including the solvent), except for H[+1] and OH[-1]].

property viscosity_dynamic: Quantity

Return the dynamic (absolute) viscosity of the solution.

Calculated from the kinematic viscosity

property viscosity_kinematic: Quantity

Return the kinematic viscosity of the solution.

Notes

The calculation is based on a model derived from the Eyring equation and presented in

\[\ln \nu = \ln {\nu_w MW_w \over \sum_i x_i MW_i } + 15 x_+^2 + x_+^3 \delta G^*_{123} + 3 x_+ \delta G^*_{23} (1-0.05x_+)\]

Where:

\[\delta G^*_{123} = a_o + a_1 (T)^{0.75}\]
\[\delta G^*_{23} = b_o + b_1 (T)^{0.5}\]

In which \(\nu\) is the kinematic viscosity, MW is the molecular weight, \(x_{+}\) is the mole fraction of cations, and \(T\) is the temperature in degrees C.

The a and b fitting parameters for a variety of common salts are included in the database.

References

Vásquez-Castillo, G.; Iglesias-Silva, G. a.; Hall, K. R. An extension of the McAllister model to correlate kinematic viscosity of electrolyte solutions. Fluid Phase Equilib. 2013, 358, 44-49.

property volume: Quantity

Return the volume of the solution.

Returns:

the volume of the solution, in L

Return type:

Quantity

water_substance

IAPWS instance describing water properties.

pyEQL.utils module

pyEQL utilities

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.

class pyEQL.utils.FormulaDict(dict=None, /, **kwargs)[source]

Bases: UserDict

Automatically converts keys on get/set using pymatgen.core.Ion.from_formula(key).reduced_formula.

This allows getting/setting/updating of Solution.components using flexible formula notation (e.g., “Na+”, “Na+1”, “Na[+]” all have the same effect)

pyEQL.utils.create_water_substance(temperature: float, pressure: float)[source]

Instantiate a water substance model from IAPWS.

Parameters:
  • temperature – the desired temperature in K

  • pressure – the desired pressure in MPa

Notes

The IAPWS97 model is much faster than IAPWS95, but the latter can do temp below zero. See https://github.com/jjgomera/iapws/issues/14. Hence, IAPWS97 will be used except when temperature is less than 0 degC.

Returns:

A IAPWS97 or IAPWS95 instance

pyEQL.utils.format_solutes_dict(solute_dict: dict, units: str)[source]

Formats a dictionary of solutes by converting the amount to a string with the provided units suitable for passing to use with the Solution class. Note that all solutes must be given in the same units.

Parameters:
  • solute_dict – The dictionary to format. This must be of the form dict{str: Number} e.g. {“Na+”: 0.5, “Cl-”: 0.9}

  • units – The units to use for the solute. e.g. “mol/kg”

Returns:

A formatted solute dictionary.

Raises:

TypeError if solute_dict is not a dictionary.

pyEQL.utils.interpret_units(unit: str) str[source]

Translate commonly used environmental units such as ‘ppm’ into strings that pint can understand.

Parameters:

unit – string representing the unit to translate

Returns: a unit that pint can understand

pyEQL.utils.standardize_formula(formula: str)[source]

Convert a chemical formula into standard form.

Parameters:

formula – the chemical formula to standardize.

Returns:

A standardized chemical formula

Raises:

ValueError if formula cannot be processed or is invalid.

Notes

Currently this method standardizes formulae by passing them through pymatgen.core.ion.Ion.reduced_formula(). For ions, this means that 1) the charge number will always be listed explicitly and 2) the charge number will be enclosed in square brackets to remove any ambiguity in the meaning of the formula. For example, ‘Na+’, ‘Na+1’, and ‘Na[+]’ will all standardize to “Na[+1]”

Module contents

pyEQL is a python package for calculating the properties of aqueous solutions and performing chemical thermodynamics computations.

copyright:

2013-2024 by Ryan S. Kingsbury

license:

LGPL, see LICENSE for more details.