The Solution Class

pyEQL Solution Class

copyright:2013-2015 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) Return the thermodynamic activity of the solute in solution
get_activity_coefficient(solute) Routine to determine the activity coefficient of a solute in solution.
get_amount(solute, units) Return the amount of ‘solute’ in the parent solution
get_chemical_potential_energy([...]) Return the total chemical potential energy of a solution (not including
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_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_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_osmotic_coefficient() Calculate the osmotic coefficient
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)
get_salt() Match ions in the solution to a parent salt.
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_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() List the activity of each species in a solution.
list_concentrations([unit]) 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
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)

Return the thermodynamic activity of the solute in solution

Parameters:

solute : str

String representing the name of the solute of interest

temperature : Quantity, optional

The temperature of the solution. Defaults to 25 degrees C if omitted

Returns:

The thermodynamic activity of the solute in question (dimensionless)

Notes

The thermodynamic activity is independent of the concentration scale used. However, the concentration and the activity coefficient must use corresponding scales. [1] [2] In this module, ionic strength, activity coefficients, and activities are all calculated based on the molal (mol/kg) concentration scale.

References

[1]http://adsorption.org/awm/utils/Activity.htm
[2]http://en.wikipedia.org/wiki/Thermodynamic_activity#Activity_coefficient
get_activity_coefficient(solute)

Routine to determine the activity coefficient of a solute in solution. The correct function is chosen based on the ionic strength of the parent solution.

Parameters:

solute : str

String representing the name of the solute of interest

Returns:

The molal (mol/kg) scale mean ion activity coefficient of the solute in question

See also

get_activity_coefficient_debyehuckel, get_activity_coefficient_guntelberg, get_activity_coefficient_davies, get_activity_coefficient_pitzer

get_amount(solute, units)

Return the amount of ‘solute’ in the parent solution

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.

Returns:

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

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

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

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

See also

get_ionic_strength, get_molar_conductivity, get_activity_coefficient

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 [4] [5]

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

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

Return the Debye length of a solution

Debye length is calculated as [6]

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

NOTE: The influence of ionic strength on the dielectric constant is not currently accounted for. The dielectric constant of pure water is used in the calculation.

Parameters:

None

Returns:

Quantity

The Debye length.

See also

get_ionic_strength, h2o.water_dielectric_constant

References

[6]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_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, but the returned value does not carry units.

Returns:

float :

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

TODO

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_mole_fraction(solute)

Return the mole fraction of ‘solute’ in the solution

Parameters:

solute : str

String representing the name of the solute of interest

Returns:

float

The mole fraction of ‘solute’ in the parent Solution object

See also

get_solvent_mass

Notes

This function assumes water is the solvent with MW = 18

Examples

TODO

get_moles_solvent()

Return the moles of solvent present in the solution

Parameters:

None

Returns:

Quantity

The moles of solvent in the solution.

get_osmotic_coefficient()

Calculate the osmotic coefficient

Returns:

Quantity :

The osmotic coefficient

Notes

For ionic strengths below 0.5 mol/kg, the osmotic coefficient is assumed to equal 1.0. 1.0 will also be returned at higher ionic strengths if appropriate Pitzer parameters are not supplied.

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

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

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

Examples

If ‘soln’ is pure water, return 0 >>> soln.get_osmotic_pressure() 0.0

If ‘soln’ is 0.5 mol/kg NaCl >>> soln.get_osmotic_pressure() 2262808... 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()

Match ions in the solution to a parent salt.

Parameters:

None

Returns:

Salt

Salt object containing information about the parent salt.

See also

salt_ion_match.py

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_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

See also

get_conductivity

Notes

Transport number is calculated according to [9] :

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

Where C is the concentration in mol/L.

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

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

See also

get_density_dynamic, get_viscosity_relative

Notes

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

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

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

float :

The thermodynamic activity of water in the solution.

Notes

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

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

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

If ‘soln’ is a 0.5 mol/kg NaCl solution at 25 degC: >>> soln.get_water_activity() 0.9835...

If ‘soln’ is a 5.11 mol/kg NaHCO2 (sodium formate) solution at 25 degC: (literature value from Cabot specialty fluids is 0.82) >>> soln.get_water_activity() 0.8631...

list_activities()

List the activity of each species in a solution.

Returns:

dict

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

list_concentrations(unit='mol/kg')

List the concentration of each species in a solution.

Parameters:

unit: str

String representing the desired concentration unit.

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'}