lower_hybrid_frequencyο
- plasmapy.formulary.frequencies.lower_hybrid_frequency(
- B: Annotated[Quantity, Unit('T')],
- n_i: Annotated[Quantity, Unit('1 / m3')],
- ion: str | int | integer | Particle | CustomParticle | Quantity,
- *,
- to_hz: bool = 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
Born_iis not aQuantity, or ion is of an inappropriate type.UnitConversionError β If either of
Born_iis in incorrect units.ValueError β If either of
Born_icontains invalid values or are of incompatible dimensions, or ion cannot be used to identify an ion or isotope.
- Warns:
UnitsWarningβ If units are not provided, SI units are assumed.
Notes
The lower hybrid frequency is given through the relation
\[\frac{1}{Ο_{lh}^2} = \frac{1}{Ο_{ci}^2 + Ο_{pi}^2} + \frac{1}{Ο_{ci}Ο_{ce}}\]where \(Ο_{ci}\) is the ion gyrofrequency, \(Ο_{ce}\) is the electron gyrofrequency, and \(Ο_{pi}\) is the ion plasma frequency.
The lower hybrid frequency constitutes a resonance for electromagnetic waves in magnetized plasmas, namely for the X-mode. These are waves with their wave electric field being perpendicular to the background magnetic field. For the lower hybrid frequency, ion and electron dynamics both play a role. As the name suggests, it has a lower frequency compared to the upper hybrid frequency. It can play an important role for heating and current drive in fusion plasmas.
Examples
>>> import astropy.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>