lower_hybrid_frequency¶
-
plasmapy.formulary.
lower_hybrid_frequency
(B: Unit(‘T’), n_i: Unit(‘1 / m3’), ion: plasmapy.particles.particle_class.Particle, to_hz=False)¶ Return the lower hybrid frequency.
Aliases:
wlh_
- Parameters
B (Quantity) – The magnetic field magnitude in units convertible to tesla.
n_i (Quantity) – Ion number density.
ion (Particle) – Representation of the ion 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 ions are assumed to be singly charged.
- Returns
omega_lh – The lower hybrid frequency in radians per second.
- Return type
- Raises
TypeError – If either of
B
orn_i
is not aQuantity
, or ion is of an inappropriate type.UnitConversionError – If either of
B
orn_i
is in incorrect units.ValueError – If either of
B
orn_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.78372...e+08 rad / s> >>> lower_hybrid_frequency(0.2*u.T, n_i=5e19*u.m**-3, ion='D+', to_hz = True) <Quantity 92050879.3... Hz>
- Other Parameters
to_hz (bool) – Set
True
to to convert function output from angular frequency to Hz