ion_sound_speed

plasmapy.formulary.speeds.ion_sound_speed(T_e: Quantity, T_i: Quantity, ion: str | Integral | Particle | CustomParticle | Quantity, n_e: Quantity = None, k: Quantity = None, gamma_e=1, gamma_i=3, Z=None) Quantity[source]

Return the ion sound speed for an electron-ion plasma.

Aliases: cs_

Parameters:
  • T_e (Quantity) – Electron temperature in units of temperature or energy per particle. If this is not given, then the electron temperature is assumed to be zero.

  • T_i (Quantity) – Ion temperature in units of temperature or energy per particle. If this is not given, then the ion temperature is assumed to be zero.

  • ion (particle-like) – Representation of the ion species (e.g., 'p+' for protons, 'D+' for deuterium, or 'He-4 +1' for singly ionized helium-4).

  • n_e (Quantity) – Electron number density. If this is not given, then ion_sound_speed will be approximated in the non-dispersive limit (\(k^2 λ_{D}^2\) will be assumed zero). If n_e is given, a value for k must also be given.

  • k (Quantity) – Wavenumber (in units of inverse length, e.g. m-1). If this is not given, then ion_sound_speed will be approximated in the non-dispersive limit (\(k^2 λ_{D}^2\) will be assumed zero). If k is given, a value for n_e must also be given.

  • gamma_e (float or int) – The adiabatic index for electrons, which defaults to 1. This value assumes that the electrons are able to equalize their temperature rapidly enough that the electrons are effectively isothermal.

  • gamma_i (float or int) – The adiabatic index for ions, which defaults to 3. This value assumes that ion motion has only one degree of freedom, namely along magnetic field lines.

  • Z (Quantity, optional) – The average ionization (arithmetic mean) for a plasma where a macroscopic description is valid. If this quantity is not given then the charge number of the ion is used. This is effectively an average ion sound speed for the plasma where multiple charge states are present.

  • mass_numb (positive integer, optional) – The mass number of an isotope corresponding to ion.

Returns:

V_S – The ion sound speed in units of meters per second.

Return type:

Quantity

Raises:
  • TypeError – If any of the arguments are not entered as keyword arguments or are of an incorrect type.

  • ValueError – If the ion mass, adiabatic index, or temperature are invalid.

  • PhysicsError – If an adiabatic index is less than one.

  • UnitConversionError – If the temperature, electron number density, or wavenumber is in incorrect units.

Warns:
  • RelativityWarning – If the ion sound speed exceeds 5% of the speed of light.

  • UnitsWarning – If units are not provided, SI units are assumed.

  • PhysicsWarning – If only one of k or n_e is given, the non-dispersive limit is assumed.

Notes

The ion sound speed \(V_S\) is given by

\[V_S = \sqrt{\frac{γ_e Z k_B T_e + γ_i k_B T_i}{m_i (1 + k^2 λ_{D}^2)}}\]

where \(γ_e\) and \(γ_i\) are the electron and ion adiabatic indices, \(k_B\) is the Boltzmann constant, \(T_e\) and \(T_i\) are the electron and ion temperatures, \(Z\) is the charge state of the ion, \(m_i\) is the ion mass, \(λ_D\) is the Debye length, and \(k\) is the wavenumber.

In the non-dispersive limit (\(k^2 λ_D^2\) is small) the equation for \(V_S\) is approximated (the denominator reduces to \(m_i\)).

When the electron temperature is much greater than the ion temperature, the ion sound velocity reduces to \(\sqrt{γ_e k_B T_e / m_i}\). Ion acoustic waves can therefore occur even when the ion temperature is zero.

Examples

>>> import astropy.units as u
>>> n = 5e19 * u.m**-3
>>> k_1 = 3e1 * u.m**-1
>>> k_2 = 3e7 * u.m**-1
>>> ion_sound_speed(
...     T_e=5e6 * u.K,
...     T_i=0 * u.K,
...     ion="p+",
...     gamma_e=1,
...     gamma_i=3,
... )
<Quantity 203155... m / s>
>>> ion_sound_speed(
...     T_e=5e6 * u.K,
...     T_i=0 * u.K,
...     n_e=n,
...     k=k_1,
...     ion="p+",
...     gamma_e=1,
...     gamma_i=3,
... )
<Quantity 203155... m / s>
>>> ion_sound_speed(
...     T_e=5e6 * u.K,
...     T_i=0 * u.K,
...     n_e=n,
...     k=k_2,
...     ion="p+",
...     gamma_e=1,
...     gamma_i=3,
... )
<Quantity 310.31... m / s>
>>> ion_sound_speed(T_e=5e6 * u.K, T_i=0 * u.K, n_e=n, k=k_1, ion="p+")
<Quantity 203155... m / s>
>>> ion_sound_speed(T_e=500 * u.eV, T_i=200 * u.eV, n_e=n, k=k_1, ion="D+")
<Quantity 229585... m / s>