matscipy.electrochemistry.poisson_boltzmann_distribution

Calculate ionic densities consistent with the Poisson-Boltzmann equation.

Copyright 2019, 2020 IMTEK Simulation University of Freiburg

Authors:

Functions

charge_density(x, c, z, u[, T, ...])

Charge density due to Poisson-Boltzmann distribution :param x: distance from the surface [m] :type x: (N,) ndarray :param c: bulk concentrations (i.e. far from the surface) [mol/m^-3] :type c: (M,) ndarray :param z: number charge of each ionic species [1] :type z: (M,) ndarray :param u: electrostatic potential at the metal/liquid interface against bulk [V] :type u: float :param T: temperature of the solution [K] (default: 298.15) :type T: float :param relative_permittivity: relative permittivity of the ionic solution, 80 for water [1] :type relative_permittivity: float.

concentration(x, c, z, u[, T, ...])

The concentration of ions near a charged surface.

debye(c, z[, T, relative_permittivity, ...])

Calculate the Debye length (in SI units per default).

gamma(u[, T])

Calculate term from Gouy-Chapmann theory.

ionic_strength(c, z)

Compute ionic strength from charges and concentrations

main()

Do stuff.

potential(x, c, z, u[, T, relative_permittivity])

The potential near a charged surface in an ionic solution.

test()

Run docstring unittests

matscipy.electrochemistry.poisson_boltzmann_distribution.ionic_strength(c, z)

Compute ionic strength from charges and concentrations

Parameters:
  • c ((M,) ndarray) – M bulk concentrations [concentration unit, i.e. mol m^-3]

  • z ((M,) ndarray) – M number charges [number charge unit, i.e. 1]

Returns:

I – ionic strength ( 1/2 * sum(z_i^2*c_i) ) [concentration unit, i.e. mol m^-3]

Return type:

float

matscipy.electrochemistry.poisson_boltzmann_distribution.debye(c, z, T=298.15, relative_permittivity=79, vacuum_permittivity=8.8541878128e-12, R=8.314462618, F=96485.33212)

Calculate the Debye length (in SI units per default).

The Debye length indicates at which distance a charge will be screened off.

Parameters:
  • c ((M,) ndarray) – bulk concentrations of each ionic species [mol/m^3]

  • z ((M,) ndarray) – charge of each ionic species [1]

  • T (float) – temperature of the solution [K] (default: 298.15)

  • relative_permittivity (float) – relative permittivity of the ionic solution [1] (default: 79)

  • vacuum_permittivity (float) – vacuum permittivity [F m^-1] (default: 8.854187817620389e-12 )

  • R (float) – molar gas constant [J mol^-1 K^-1] (default: 8.3144598)

  • F (float) – Faraday constant [C mol^-1] (default: 96485.33289)

Returns:

lambda_D – Debye length, sqrt( epsR*eps*R*T/(2*F^2*I) ) [m]

Return type:

float

matscipy.electrochemistry.poisson_boltzmann_distribution.gamma(u, T=298.15)

Calculate term from Gouy-Chapmann theory.

Parameters:
  • u (float) – electrostatic potential at the metal/solution boundary in Volts, e.g. 0.05 [V]

  • T (float) – temperature of the solution in Kelvin [K] (default: 298.15)

Return type:

float

matscipy.electrochemistry.poisson_boltzmann_distribution.potential(x, c, z, u, T=298.15, relative_permittivity=79)

The potential near a charged surface in an ionic solution.

A single value is returned, which specifies the value of the potential at this distance.

The decimal package is used for increased precision. If only normal float precision is used, the potential is a step function. Steps in the potential result in unphysical particle concentrations.

Parameters:
  • x ((N,) ndarray) – z-distance from the surface [m]

  • c ((M,) ndarray) – bulk concentrations of each ionic species [mol/m^3]

  • z ((M,) ndarray) – charge number of each ionic species [1]

  • u (float) – electrostatic potential at the metal/solution boundary in Volts, e.g. 0.05 [V]

  • T (float) – temperature of the solution [Kelvin] (default: 298.15)

  • relative_permittivity (float) – relative permittivity of the ionic solution [] (default: 79)

Returns:

phi – Electrostatic potential [V]

Return type:

(N,) ndarray

matscipy.electrochemistry.poisson_boltzmann_distribution.concentration(x, c, z, u, T=298.15, relative_permittivity=79)

The concentration of ions near a charged surface.

Calculates the molar concentration of ions of a species, at a distance x away from a substrate/solution interface. Potential difference u between substrate and bulk solution leads to non-neutrality close to the interface, with concentrations converging to their bulk values at high distances.

Parameters:
  • x ((N,) ndarray) – distance from the surface [m]

  • c ((M,) ndarray) – bulk concentrations (i.e. far from the surface) [mol/m^-3]

  • z ((M,) ndarray) – number charge of each ionic species [1]

  • u (float) – electrostatic potential at the metal/liquid interface against bulk [V]

  • T (float) – temperature of the solution [K] (default: 298.15)

  • relative_permittivity (float) – relative permittivity of the ionic solution, 80 for water [1]

Returns:

c – molar concentrations of ion species [mol/m^3]

Return type:

(M,N) ndarray

matscipy.electrochemistry.poisson_boltzmann_distribution.charge_density(x, c, z, u, T=298.15, relative_permittivity=79)

Charge density due to Poisson-Boltzmann distribution :param x: distance from the surface [m] :type x: (N,) ndarray :param c: bulk concentrations (i.e. far from the surface) [mol/m^-3] :type c: (M,) ndarray :param z: number charge of each ionic species [1] :type z: (M,) ndarray :param u: electrostatic potential at the metal/liquid interface against bulk [V] :type u: float :param T: temperature of the solution [K] (default: 298.15) :type T: float :param relative_permittivity: relative permittivity of the ionic solution, 80 for water [1] :type relative_permittivity: float

Returns:

c – charge density [C/m^3]

Return type:

(N,) ndarray

matscipy.electrochemistry.poisson_boltzmann_distribution.test()

Run docstring unittests

matscipy.electrochemistry.poisson_boltzmann_distribution.main()

Do stuff.