upper_hybrid_frequency

plasmapy.physics.parameters.upper_hybrid_frequency(B: Unit("T"), n_e: Unit("1 / m3"), to_hz=False)

Return the upper hybrid frequency.

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

  • n_e (Quantity) – The electron number density.

Returns

omega_uh – The upper hybrid frequency in radians per second.

Return type

Quantity

Raises
  • TypeError – If either of B or n_e is not a Quantity.

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

  • ValueError – If either of B or n_e contains invalid values or are of incompatible dimensions.

Warns

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

Notes

The upper hybrid frequency is given through the relation

\[\omega_{uh}^2 = \omega_{ce}^2 + \omega_{pe}^2\]

where \(\omega_{ce}\) is the electron gyrofrequency and \(\omega_{pe}\) is the electron plasma frequency.

Example

>>> from astropy import units as u
>>> upper_hybrid_frequency(0.2*u.T, n_e=5e19*u.m**-3)
<Quantity 4.00459419e+11 rad / s>
>>> upper_hybrid_frequency(0.2*u.T, n_e=5e19*u.m**-3, to_hz = True)
<Quantity 63735096112.815445 Hz>
Other Parameters

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