The Solution Class

The Solution class defines a pythonic interface for creating, modifying, and estimating properties of electrolyte solutions. It is the core feature of pyEQL and the primary user-facing class.

Creating a solution

A Solution created with no arguments will default to pure water at pH=7, T=25 degrees Celsius, and 1 atm pressure.

>>> from pyEQL import Solution
>>> s1 = Solution()
>>> s1.pH
6.998877352386266

Alternatively, you can use the autogenerate() function to easily create common solutions like seawater:

>>> from pyEQL.functions import autogenerate
>>> s2 = autogenerate('seawater')
<pyEQL.solution.Solution object at 0x7f057de6b0a0>

You can inspect the solutes present in the solution via the components attribute. This comprises a dictionary of solute formula: moles, where ‘moles’ is the number of moles of that solute in the Solution. Note that the solvent (water) is present in components, too.

>>> s2.components
{'H2O': 55.34455401423017,
 'H+': 7.943282347242822e-09,
 'OH-': 8.207436858780226e-06,
 'Na+': 0.46758273714962967,
 'Mg+2': 0.052661180523467986,
 'Ca+2': 0.010251594148212318,
 'K+': 0.010177468379526856,
 'Sr+2': 9.046483353663286e-05,
 'Cl-': 0.54425785619973,
 'SO4-2': 0.028151873448454025,
 'HCO3-': 0.001712651176926199,
 'Br-': 0.0008395244921424563,
 'CO3-2': 0.00023825904349479546,
 'B(OH)4': 0.0001005389715937341,
 'F-': 6.822478260456777e-05,
 'B(OH)3': 0.0003134669156396757,
 'CO2': 9.515218476861175e-06
 }

To get the amount of a specific solute, use get_amount() and specify the units you want:

>>> s2.get_amount('Na+', 'g/L')
<Quantity(10.6636539, 'gram / liter')>

Finally, you can manually create a solution with any list of solutes, temperature, pressure, etc. that you need:

>>> from pyEQL import Solution
>>> s1 = Solution(solutes={'Na+':'0.5 mol/kg', 'Cl-': '0.5 mol/kg'},
                  pH=8,
                  temperature = '20 degC',
                  volume='8 L'
                  )

Class reference

The remainder of this page contains detailed information on each of the methods, attributes, and properties in Solution. Use the sidebar on the right for easier navigation.

class pyEQL.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, solvent: str | list = 'H2O', engine: Literal['native', 'ideal'] = 'native', database: str | Path | Store | None = None)[source]

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

__init__(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, solvent: str | list = 'H2O', engine: Literal['native', 'ideal'] = 'native', database: str | Path | Store | None = None)[source]
Parameters:
  • solutes

    dict, optional. Keys must be the chemical formula, while values must be str Quantity representing the amount. For example:

    {“Na+”: “0.1 mol/L”, “Cl-”: “0.1 mol/L”}

    Note that an older “list of lists” syntax is also supported; however this will be deprecated in the future and is no longer recommended. The equivalent list syntax for the above example is

    [[“Na+”, “0.1 mol/L”], [“Cl-”, “0.1 mol/L”]]

    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

  • pe – the pE value (redox potential) of the solution. Lower values = more reducing, higher values = more oxidizing. At pH 7, water is stable between approximately -7 to +14. The default value corresponds to a pE value typical of natural waters in equilibrium with the atmosphere.

  • solvent – Formula of the solvent. Solvents other than water are not supported at this time.

  • engine

  • database – path to a .json file (str or Path) or maggma Store instance that contains serialized SoluteDocs. None (default) will use the built-in pyEQL database.

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

Return the total mass of the solution.

The mass is calculated each time this method is called. :param None:

Returns:

Quantity

Return type:

the mass of the solution, in kg

property solvent_mass

Return the mass of the solvent.

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

Returns:

Quantity

Return type:

the mass of the solvent, in kg

See also

get_amount()

property volume: Quantity

Return the volume of the solution.

Returns:

Quantity: the volume of the solution, in L

property temperature: Quantity

Return the temperature of the solution in Kelvin.

property pressure: Quantity

Return the hydrostatic pressure of the solution in atm.

property pH: Quantity

Return the pH of the solution.

p(solute, activity=True)[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:

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

Return type:

Quantity

property density: Quantity

Return the density of the solution.

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

Returns:

Quantity

Return type:

The density of the solution.

property dielectric_constant: Quantity

Returns the dielectric constant of the solution.

Parameters:

None

Returns:

Quantity

Return type:

the dielectric constant of the solution, dimensionless.

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

Return the dynamic (absolute) viscosity of the solution.

Calculated from the kinematic viscosity

See Also:

viscosity_kinematic

property 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

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

See Also:

viscosity_dynamic()

property conductivity

Compute the electrical conductivity of the solution.

Parameters:

None

Returns:

The electrical conductivity of the solution in Siemens / meter.

Return type:

Quantity

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 [aq] [hc]

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

Where:

\[\alpha = \begin{cases} {\frac{0.6}{\sqrt{| z_{i} | }}} & {I < 0.36 | z_{i} | } {\frac{\sqrt{I}}{| z_i |}} & otherwise \end{cases}\]

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

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:

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')>
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 n_i z_i\]

where \(n_i\) is the number of moles, and \(z_i\) is the charge on species i.

Returns:

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

Return type:

float

property alkalinity

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]

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

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

The hardness of the solution in mg/L as CaCO3

Return type:

Quantity

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 .. [wk3] https://en.wikipedia.org/wiki/Debye_length#Debye_length_in_an_electrolyte

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.

Parameters:

None

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 osmotic_pressure

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

Returns

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

See also

get_water_activity get_osmotic_coefficient get_salt

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
0.0
>>> s1 = pyEQL.Solution([['Na+','0.2 mol/kg'],['Cl-','0.2 mol/kg']])
>>> soln.osmotic_pressure
<Quantity(906516.7318131207, 'pascal')>
get_amount(solute, units)[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’) 2. substance per mass of solvent (e.g., ‘mol/kg’) 3. mass of substance (e.g., ‘kg’) 4. moles of substance (‘mol’) 5. mole fraction (‘fraction’) 6. percent by weight (%) 7. number of molecules (‘count’)

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

Return type:

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

get_total_amount(element, units)[source]

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’

Return type:

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

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.

See also

get_amount

add_solute(formula, amount)[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 (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’

add_solvent(formula, amount)[source]

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

add_amount(solute, amount)[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

Return type:

Nothing. The concentration of solute is modified.

set_amount(solute, amount)[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.

Return type:

Nothing. The concentration of solute is modified.

get_total_moles_solute() Quantity[source]

Return the total moles of all solute in the solution.

get_moles_solvent() Quantity[source]

Return the moles of solvent present in the solution.

Returns

The moles of solvent in the solution.

get_osmolarity(activity_correction=False)[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_osmolality(activity_correction=False)[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_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., 1 M MgCl2 yields 1 M Mg+2 and 2 M Cl-).

Parameters:

None

Returns:

Salt object containing information about the parent salt.

Return type:

Salt

See also

get_activity(), get_activity_coefficient(), get_water_activity(), get_osmotic_coefficient(), get_osmotic_pressure(), get_viscosity_kinematic()

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() dict[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.

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

Parameters:

None

Returns:

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

Return the activity coefficient of a solute in solution.

The model used to calculte 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

Quantity: the activity coefficient as a dimensionless pint Quantity

get_activity(solute: str, scale: Literal['molal', 'molar', 'rational'] = 'molal', verbose: bool = False)[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)

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_osmotic_coefficient(scale: Literal['molal', 'molar', 'rational'] = 'molal')[source]

Return the osmotic coefficient of an aqueous solution.

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

get_water_activity()[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')>
get_chemical_potential_energy(activity_correction=True)[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:

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

Return type:

Quantity

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_property(solute: str, name: str) Quantity | None[source]

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

Return type:

The desired parameter or None if not found

get_transport_number(solute, activity_correction=False)[source]

Calculate the transport number of the solute in the solution.

Parameters:
  • solute – String identifying the solute for which the transport number is to be calculated.

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

  • Returns – The transport number of solute

  • 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.

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

  • 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.

get_molar_conductivity(solute)[source]

Calculate the molar (equivalent) conductivity for a solute.

Parameters:

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

Returns

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

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

Note that the diffusion coefficient is strongly variable with temperature.

References

[smed]

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

get_mobility(solute)[source]

Calculate the ionic mobility of the solute.

Parameters:

solute (str) – String identifying the solute for which the mobility is to be calculated.

Returns:

float

Return type:

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

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

References

[smed]

Smedley, Stuart I. The Interpretation of Ionic Conductivity in Liquids. Plenum Press, 1980.

get_lattice_distance(solute)[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:

Quantity

Return type:

The average distance between solute molecules

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}}\]
copy()[source]

Return a copy of the solution.

as_dict() dict[source]

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

classmethod from_dict(d: dict) Solution[source]

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

list_solutes()[source]

List all the solutes in the solution.

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

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:

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

Return type:

dict

list_activities(decimals=4)[source]

List the activity of each species in a solution.

Parameters:

decimals (int) – The number of decimal places to display. Defaults to 4.

Returns:

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

Return type:

dict

to_json() str

Returns a json string representation of the MSONable object.

unsafe_hash()

Returns an hash of the current object. This uses a generic but low performance method of converting the object to a dictionary, flattening any nested keys, and then performing a hash on the resulting object

classmethod validate_monty(v)

pydantic Validator for MSONable pattern