plasma_frequency

plasmapy.formulary.frequencies.plasma_frequency(n: Unit('1 / m3'), particle: str | Integral | Particle | CustomParticle | Quantity, z_mean=None, to_hz=False)

Calculate the particle plasma frequency.

Aliases: wp_

Lite Version: plasma_frequency_lite

Parameters:
  • n (Quantity) – Particle number density in units convertible to m-3.

  • particle (Particle) – Representation of the particle species (e.g., "p" for protons, "D+" for deuterium, or "He-4 +1" for singly ionized helium-4). If no charge state information is provided, then the particles are assumed to be singly charged.

  • z_mean (Real, optional) – The average ionization (arithmetic mean) for the particle species in the plasma. Typically, the charge state will be derived from the particle argument, but this keyword will override that behavior.

Returns:

omega_p – The particle plasma frequency in radians per second. Setting keyword to_hz=True will apply the factor of \(1/2π\) and yield a value in Hz.

Return type:

Quantity

Raises:
  • TypeError – If n is not a Quantity or particle is not of an appropriate type.

  • UnitConversionError – If n is not in correct units.

  • ValueError – If n contains invalid values or particle cannot be used to identify a particle or isotope.

Warns:

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

Notes

The particle plasma frequency is

\[ω_{p} = Z |e| \sqrt{\frac{n}{\epsilon_0 m}}\]

where \(m\) is the mass of the particle, \(e\) is the fundamental unit of charge, \(Z\) is the average charge state z_mean of the particle species, \(n\) is the particle number density. This form of the plasma frequency has units of radians / s, but using the to_hz will apply the factor of \(1/2π\) to give a value in Hz.

Examples

>>> from astropy import units as u
>>> plasma_frequency(1e19 * u.m**-3, particle='p')
<Quantity 4.16329...e+09 rad / s>
>>> plasma_frequency(1e19 * u.m**-3, particle='p', to_hz=True)
<Quantity 6.62608...e+08 Hz>
>>> plasma_frequency(1e19 * u.m**-3, particle='D+')
<Quantity 2.94462...e+09 rad / s>
>>> plasma_frequency(1e19 * u.m**-3, 'e-')
<Quantity 1.78398...e+11 rad / s>
>>> plasma_frequency(1e19 * u.m**-3, 'e-', to_hz=True)
<Quantity 2.83930...e+10 Hz>

For user convenience plasma_frequency_lite is bound to this function and can be used as follows.

>>> from plasmapy.particles import Particle
>>> mass = Particle("p").mass.value
>>> plasma_frequency.lite(n=1e19, mass=mass, z_mean=1)
416329...
>>> plasma_frequency.lite(n=1e19, mass=mass, z_mean=1, to_hz=True)
662608...
Parameters:

to_hz (bool) – Set True to to convert function output from angular frequency to Hz