The Solution Class

pyEQL Solution Class

copyright:2013-2018 by Ryan S. Kingsbury
license:LGPL, see LICENSE for more details.
class pyEQL.solution.Solution(solutes=[], **kwargs)

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

Parameters:
solutes : list of lists, optional

See add_solute() documentation for formatting of this list. Defaults to empty (pure solvent) if omitted

volume : str, optional

Volume of the solvent, including the unit. Defaults to ‘1 L’ if omitted. Note that the total solution volume will be computed using partial molar volumes of the respective solutes as they are added to the solution.

temperature : str, optional

The solution temperature, including the unit. Defaults to ‘25 degC’ if omitted.

pressure : Quantity, optional

The ambient pressure of the solution, including the unit. Defaults to ‘1 atm’ if omitted.

pH : number, optional

Negative log of H+ activity. If omitted, the solution will be initialized to pH 7 (neutral) with appropriate quantities of H+ and OH- ions

Returns:
Solution

A Solution object.

See also

add_solute

Examples

>>> s1 = pyEQL.Solution([['Na+','1 mol/L'],['Cl-','1 mol/L']],temperature='20 degC',volume='500 mL')
>>> print(s1)
Components: 
['H2O', 'Cl-', 'H+', 'OH-', 'Na+']
Volume: 0.5 l
Density: 1.0383030844030992 kg/l

Methods

add_amount(solute, amount) Add the amount of ‘solute’ to the parent solution.
add_solute(formula, amount[, parameters]) Primary method for adding substances to a pyEQL solution
add_solvent(formula, amount) Same as add_solute but omits the need to pass solvent mass to pint
copy() Return a copy of the solution
get_activity(solute[, scale, verbose]) Return the thermodynamic activity of the solute in solution on the molal scale.
get_activity_coefficient(solute[, scale, …]) Return the activity coefficient of a solute in solution.
get_alkalinity() Return the alkalinity or acid neutralizing capacity of a solution
get_amount(solute, units) Return the amount of ‘solute’ in the parent solution.
get_bjerrum_length() Return the Bjerrum length of a solution
get_charge_balance() Return the charge balance of the solution.
get_chemical_potential_energy([…]) Return the total chemical potential energy of a solution (not including pressure or electric effects)
get_conductivity() Compute the electrical conductivity of the solution.
get_debye_length() Return the Debye length of a solution
get_density() Return the density of the solution.
get_dielectric_constant() Returns the dielectric constant of the solution.
get_hardness() Return the hardness of a solution.
get_ionic_strength() Return the ionic strength of the solution.
get_lattice_distance(solute) Calculate the average distance between molecules
get_mass() Return the total mass of the solution.
get_mobility(solute) Calculate the ionic mobility of the solute
get_molar_conductivity(solute) Calculate the molar (equivalent) conductivity for a solute
get_mole_fraction(solute) Return the mole fraction of ‘solute’ in the solution
get_moles_solvent() Return the moles of solvent present in the solution
get_osmolality([activity_correction]) Return the osmolality of the solution in Osm/kg
get_osmolarity([activity_correction]) Return the osmolarity of the solution in Osm/L
get_osmotic_coefficient([scale]) Return the osmotic coefficient of an aqueous solution.
get_osmotic_pressure() Return the osmotic pressure of the solution relative to pure water.
get_pressure() Return the hydrostatic pressure of the solution.
get_property(solute, name) Retrieve a thermodynamic property (such as diffusion coefficient) for solute, and adjust it from the reference conditions to the conditions of the solution
get_salt() Determine the predominant salt in a solution of ions.
get_salt_list() Determine the predominant salt in a solution of ions.
get_solute(i) Return the specified solute object.
get_solvent() Return the solvent object.
get_solvent_mass() Return the mass of the solvent.
get_temperature() Return the temperature of the solution.
get_total_amount(element, units) Return the total amount of ‘element’ (across all solutes) in the solution.
get_total_moles_solute() Return the total moles of all solute in the solution
get_transport_number(solute[, …]) Calculate the transport number of the solute in the solution
get_viscosity_dynamic() Return the dynamic (absolute) viscosity of the solution.
get_viscosity_kinematic() Return the kinematic viscosity of the solution.
get_viscosity_relative() Return the viscosity of the solution relative to that of water
get_volume() Return the volume of the solution.
get_water_activity() Return the water activity.
list_activities([decimals]) List the activity of each species in a solution.
list_concentrations([unit, decimals, type]) List the concentration of each species in a solution.
list_solutes() List all the solutes in the solution.
p(solute[, activity]) Return the negative log of the activity of solute.
set_amount(solute, amount) Set the amount of ‘solute’ in the parent solution.
set_pressure(pressure) Set the hydrostatic pressure of the solution.
set_temperature(temperature) Set the solution temperature.
set_volume(volume) Change the total solution volume to volume, while preserving all component concentrations
list_salts  
add_amount(solute, amount)

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.

See also

Solute.add_moles

add_solute(formula, amount, parameters={})

Primary method for adding substances to a pyEQL solution

Parameters:
formula : str

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’).

amount : str

The amount of substance in the specified unit system. The string should contain both a quantity and a pint-compatible representation of a unit. e.g. ‘5 mol/kg’ or ‘0.1 g/L’

parameters : dictionary, optional

Dictionary of custom parameters, such as diffusion coefficients, transport numbers, etc. Specify parameters as key:value pairs separated by commas within curly braces, e.g. {diffusion_coeff:5e-10,transport_number:0.8}. The ‘key’ is the name that will be used to access the parameter, the value is its value.

add_solvent(formula, amount)

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

copy()

Return a copy of the solution

TODO - clarify whether this is a deep or shallow copy

get_activity(solute, scale='molal', verbose=False)

Return the thermodynamic activity of the solute in solution on the molal scale.

Parameters:
solute : str

String representing the name of the solute of interest

scale : str, optional

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 : bool, optional

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)

Notes

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

References

[1]Robinson, R. A.; Stokes, R. H. Electrolyte Solutions: Second Revised Edition; Butterworths: London, 1968, p.32.
get_activity_coefficient(solute, scale='molal', verbose=False)

Return the activity coefficient of a solute in solution.

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

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 [3]

\[f_\pm = \gamma_\pm * (1 + M_w \sum_i \nu_i \m_i)\]
\[y_\pm = 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 : str

String representing the name of the solute of interest

scale : str, optional

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.

verbose : bool, optional

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 mean ion activity coefficient of the solute in question on the selected scale.

See also

get_ionic_strength, get_salt, activity_correction.get_activity_coefficient_debyehuckel, activity_correction.get_activity_coefficient_guntelberg, activity_correction.get_activity_coefficient_davies, activity_correction.get_activity_coefficient_pitzer

Notes

For multicomponent mixtures, pyEQL implements the “effective Pitzer model” presented by Mistry et al. [4]. 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 = 2 I \over (\nu_+ z_+^2 + \nu_- z_- ^2)\]

References

[2]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
[3]Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 165. Wiley Interscience, 1996.
[4]Robinson, R. A.; Stokes, R. H. Electrolyte Solutions: Second Revised Edition; Butterworths: London, 1968, p.32.
get_alkalinity()

Return the alkalinity or acid neutralizing capacity of a solution

Returns:
Quantity :

The alkalinity of the solution in mg/L as CaCO3

Notes

The alkalinity is calculated according to: [5]

\[Alk = F \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 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

[5]Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 165. Wiley Interscience, 1996.
get_amount(solute, units)

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’) 1. substance per mass of solvent (e.g., ‘mol/kg’) 1. mass of substance (e.g., ‘kg’) 1. moles of substance (‘mol’) 1. mole fraction (‘fraction’) 1. percent by weight (%)

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_bjerrum_length()

Return the Bjerrum length of a solution

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

\[\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.

Parameters:
None
Returns:
Quantity

The Bjerrum length, in nanometers.

References

[6]https://en.wikipedia.org/wiki/Bjerrum_length

Examples

>>> s1 = pyEQL.Solution()
>>> s1.get_bjerrum_length()
<Quantity(0.7152793009386953, 'nanometer')>
get_charge_balance()

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.

Returns:
float :

The charge balance of the solution, in equivalents.

Notes

The charge balance is calculated according to:

\[CB = F \sum_i n_i z_i\]

Where \(n_i\) is the number of moles, \(z_i\) is the charge on species i, and \(F\) is the Faraday constant.

get_chemical_potential_energy(activity_correction=True)

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: [7]

\[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

[7]Koga, Yoshikata, 2007. Solution Thermodynamics and its Application to Aqueous Solutions: A differential approach. Elsevier, 2007, pp. 23-37.
get_conductivity()

Compute the electrical conductivity of the solution.

Parameters:
None
Returns:
Quantity

The electrical conductivity of the solution in Siemens / meter.

Notes

Conductivity is calculated by summing the molar conductivities of the respective solutes, but they are activity-corrected and adjusted using an empricial exponent. This approach is used in PHREEQC and Aqion models [8] [9]

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

Where:

\[\begin{split}\alpha = \begin{cases} {0.6 \over \sqrt{|z_i|}} & {I < 0.36|z_i|} \\ {\sqrt{I} \over |z_i|} & otherwise \end{cases}\end{split}\]

Note: PHREEQC uses the molal rather than molar concentration according to http://wwwbrr.cr.usgs.gov/projects/GWC_coupled/phreeqc/phreeqc3-html/phreeqc3-43.htm

References

[8]http://www.aqion.de/site/77
[9]http://www.hydrochemistry.eu/exmpls/sc.html
get_debye_length()

Return the Debye length of a solution

Debye length is calculated as [10]

\[\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.

Parameters:
None
Returns:
Quantity

The Debye length, in nanometers.

References

[10]https://en.wikipedia.org/wiki/Debye_length#Debye_length_in_an_electrolyte
get_density()

Return the density of the solution.

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

Returns:
Quantity: The density of the solution.
get_dielectric_constant()

Returns the dielectric constant of the solution.

Parameters:
None
Returns:
Quantity: the dielectric constant of the solution, dimensionless.

Notes

Implements the following equation as given by [11]

\[\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

[11][1] 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.

get_hardness()

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.

Parameters:
None
Returns:
Quantity

The hardness of the solution in mg/L as CaCO3

get_ionic_strength()

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:
Quantity :

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

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.get_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.get_ionic_strength()
<Quantity(1.0000001004383303, 'mole / kilogram')>
get_lattice_distance(solute)

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:
Quantity : The average distance between solute molecules

Notes

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

\[d = ( C_i N_A ) ^ {-{1\over3}}\]

Examples

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

Return the total mass of the solution.

The mass is calculated each time this method is called. Parameters ———- None

Returns:
Quantity: the mass of the solution, in kg
get_mobility(solute)

Calculate the ionic mobility of the solute

Parameters:
solute : str

String identifying the solute for which the mobility is to be calculated.

Returns:
float : The ionic mobility. Zero if the solute is not charged.

Notes

This function uses the Einstein relation to convert a diffusion coefficient into an ionic mobility [12]

\[\mu_i = {F |z_i| D_i \over RT}\]

References

[12]Smedley, Stuart I. The Interpretation of Ionic Conductivity in Liquids. Plenum Press, 1980.
get_molar_conductivity(solute)

Calculate the molar (equivalent) conductivity for a solute

Parameters:
solute : str

String identifying the solute for which the molar conductivity is to be calculated.

Returns:
float

The molar or equivalent conductivity of the species in the solution. Zero if the solute is not charged.

Notes

Molar conductivity is calculated from the Nernst-Einstein relation [13]

\[\kappa_i = {z_i^2 D_i F^2 \over RT}\]

Note that the diffusion coefficient is strongly variable with temperature.

References

[13]Smedley, Stuart. The Interpretation of Ionic Conductivity in Liquids, pp 1-9. Plenum Press, 1980.

Examples

TODO

get_mole_fraction(solute)

Return the mole fraction of ‘solute’ in the solution

Notes

This function is DEPRECATED and will raise a warning when called. Use get_amount() instead and specify ‘fraction’ as the unit type.

get_moles_solvent()

Return the moles of solvent present in the solution

Parameters:
None
Returns:
Quantity

The moles of solvent in the solution.

get_osmolality(activity_correction=False)

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)

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='molal')

Return the osmotic coefficient of an aqueous solution.

Osmotic coefficient is calculated using the Pitzer model.[#]_ 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 [14]

\[g = - \phi * M_w \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 : str, optional

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
——-
Quantity :

The osmotic coefficient

Notes

For multicomponent mixtures, pyEQL adopts the “effective Pitzer model” presented by Mistry et al. [15]. 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

[14]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
[15]Robinson, R. A.; Stokes, R. H. Electrolyte Solutions: Second Revised Edition; Butterworths: London, 1968, p.32.
[16]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.9235996615888572, 'dimensionless')>
>>> s1 = pyEQL.Solution([['Mg+2','0.3 mol/kg'],['Cl-','0.6 mol/kg']],temperature='30 degC')
>>> s1.get_osmotic_coefficient()
<Quantity(0.891154788474231, 'dimensionless')>
get_osmotic_pressure()

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

Parameters:
None
Returns:
Quantity

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

Notes

Osmotic pressure is calculated based on the water activity [16] [17] :

\[\Pi = {RT \over 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

[17]Sata, Toshikatsu. Ion Exchange Membranes: Preparation, Characterization, and Modification. Royal Society of Chemistry, 2004, p. 10.
[18]http://en.wikipedia.org/wiki/Osmotic_pressure#Derivation_of_osmotic_pressure

Examples

>>> s1=pyEQL.Solution()
>>> s1.get_osmotic_pressure()
0.0
>>> s1 = pyEQL.Solution([['Na+','0.2 mol/kg'],['Cl-','0.2 mol/kg']])
>>> soln.get_osmotic_pressure()
<Quantity(906516.7318131207, 'pascal')>        
get_pressure()

Return the hydrostatic pressure of the solution.

Returns:
Quantity: The hydrostatic pressure of the solution, in atm.
get_property(solute, name)

Retrieve a thermodynamic property (such as diffusion coefficient) for solute, and adjust it from the reference conditions to the conditions of the solution

Parameters:
solute: str

String representing the chemical formula of the solute species

name: str

The name of the property needed, e.g. ‘diffusion coefficient’

Returns:
Quantity: The desired parameter
get_salt()

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 direclty 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., 1 M MgCl2 yields 1 M Mg+2 and 2 M Cl-).

Parameters:
None
Returns:
Salt

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_list()

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 direclty model these quantities.

The get_salt_list() method examines the ionic composition of a solution and simplifies it into a list of salts. The method retuns 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-).

Parameters:
None
Returns:
dict

A dictionary of Salt objects, keyed to the salt formula

get_solute(i)

Return the specified solute object.

get_solvent()

Return the solvent object.

get_solvent_mass()

Return the mass of the solvent.

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

Parameters:
None
Returns:
Quantity: the mass of the solvent, in kg

See also

get_amount

get_temperature()

Return the temperature of the solution.

Parameters:
None
Returns:
Quantity: The temperature of the solution, in Kelvin.
get_total_amount(element, units)

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

Parameters:
element : str

String representing the name of the element of interest

units : str

Units desired for the output. 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

See also

get_amount

Notes

There is currently no way to distinguish between different oxidation states of the same element (e.g. TOTFe(II) vs. TOTFe(III)). This is planned for a future release. (TODO)

get_total_moles_solute()

Return the total moles of all solute in the solution

get_transport_number(solute, activity_correction=False)

Calculate the transport number of the solute in the solution

Parameters:
solute : str

String identifying the solute for which the transport number is to be calculated.

activity_correction: bool

If True, the transport number will be corrected for activity following the same method used for solution conductivity. Defaults to False if omitted.

Returns:
float

The transport number of solute

Notes

Transport number is calculated according to [18] :

\[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.

If activity_correction is True, the contribution of each ion to the transport number is corrected with an activity factor. See the documentation for get_conductivity() for an explanation of this correction.

References

[19]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.

get_viscosity_dynamic()

Return the dynamic (absolute) viscosity of the solution.

Calculated from the kinematic viscosity

get_viscosity_kinematic()

Return the kinematic viscosity of the solution.

Notes

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

\[\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 :math: 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

[20]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.
get_viscosity_relative()

Return the viscosity of the solution relative to that of water

This is calculated using a simplified form of the Jones-Dole equation:

\[\eta_{rel} = 1 + \sum_i B_i m_i\]

Where \(m\) is the molal concentration and \(B\) is an empirical parameter.

See <http://downloads.olisystems.com/ResourceCD/TransportProperties/Viscosity-Aqueous.pdf> <http://www.nrcresearchpress.com/doi/pdf/10.1139/v77-148> <http://apple.csgi.unifi.it/~fratini/chen/pdf/14.pdf>

get_volume()

Return the volume of the solution.

Parameters:
None
Returns:
Quantity: the volume of the solution, in L
get_water_activity()

Return the water activity.

Returns:
Quantity :

The thermodynamic activity of water in the solution.

Notes

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

\[\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

[21]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')>
list_activities(decimals=4)

List the activity of each species in a solution.

Parameters:
decimals: int

The number of decimal places to display. Defaults to 4.

Returns:
dict

Dictionary containing a list of the species in solution paired with their activity

list_concentrations(unit='mol/kg', decimals=4, type='all')

List the concentration of each species in a solution.

Parameters:
unit: str

String representing the desired concentration unit.

decimals: int

The number of decimal places to display. Defaults to 4.

type : str

The type of component to be sorted. Defaults to ‘all’ for all solutes. Other valid arguments are ‘cations’ and ‘anions’ which return lists of cations and anions, respectively.

Returns:
dict

Dictionary containing a list of the species in solution paired with their amount in the specified units

list_solutes()

List all the solutes in the solution.

p(solute, activity=True)

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.

Examples

TODO

set_amount(solute, amount)

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.

See also

Solute.set_moles

set_pressure(pressure)

Set the hydrostatic pressure of the solution.

Parameters:
pressure : str

String representing the temperature, e.g. ‘25 degC’

set_temperature(temperature)

Set the solution temperature.

Parameters:
temperature : str

String representing the temperature, e.g. ‘25 degC’

set_volume(volume)

Change the total solution volume to volume, while preserving all component concentrations

Parameters:
volume : str quantity

Total volume of the solution, including the unit, e.g. ‘1 L’

Examples

>>> mysol = Solution([['Na+','2 mol/L'],['Cl-','0.01 mol/L']],volume='500 mL')
>>> print(mysol.get_volume())
0.5000883925072983 l
>>> mysol.list_concentrations()
{'H2O': '55.508435061791985 mol/kg', 'Cl-': '0.00992937605907076 mol/kg', 'Na+': '2.0059345573880325 mol/kg'}
>>> mysol.set_volume('200 mL')
>>> print(mysol.get_volume())
0.2 l
>>> mysol.list_concentrations()
{'H2O': '55.50843506179199 mol/kg', 'Cl-': '0.00992937605907076 mol/kg', 'Na+': '2.0059345573880325 mol/kg'}