lower_hybrid_frequency

plasmapy.physics.parameters.lower_hybrid_frequency(B, n_i, ion='p+')

Return the lower hybrid frequency.

Parameters
  • B (Quantity) – The magnetic field magnitude in units convertible to tesla.

  • n_i (Quantity) – Ion number density.

  • ion (str, optional) – Representation of the ion species (e.g., ‘p’ for protons, ‘D+’ for deuterium, or ‘He-4 +1’ for singly ionized helium-4), which defaults to protons. If no charge state information is provided, then the ions are assumed to be singly charged.

Returns

omega_lh – The lower hybrid frequency in radians per second.

Return type

Quantity

Raises
  • TypeError – If either of B or n_i is not a Quantity, or ion is of an inappropriate type.

  • UnitConversionError – If either of B or n_i is in incorrect units.

  • ValueError – If either of B or n_i contains invalid values or are of incompatible dimensions, or ion cannot be used to identify an ion or isotope.

Warns

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

Notes

The lower hybrid frequency is given through the relation

\[\frac{1}{\omega_{lh}^2} = \frac{1}{\omega_{ci}^2 + \omega_{pi}^2} + \frac{1}{\omega_{ci}\omega_{ce}}\]

where \(\omega_{ci}\) is the ion gyrofrequency, \(\omega_{ce}\) is the electron gyrofrequency, and \(\omega_{pi}\) is the ion plasma frequency.

Example

>>> from astropy import units as u
>>> lower_hybrid_frequency(0.2*u.T, n_i=5e19*u.m**-3, ion='D+')
<Quantity 5.78372733e+08 rad / s>