gyrofrequency

plasmapy.formulary.frequencies.gyrofrequency(B: Quantity, particle: str | Integral | Particle | CustomParticle | Quantity, signed: bool = False, Z: Real | None = None, mass_numb: Integral | None = None, *, to_hz=False) Quantity[source]

Calculate the particle gyrofrequency in units of radians per second.

Aliases: oc_, wc_

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

  • particle (particle-like) – Representation of the particle species (e.g., 'p+' for protons, 'D+' for deuterium, or 'He-4 1+' for singly ionized helium-4).

  • signed (bool, default: False) – The gyrofrequency can be defined as signed, where its sign is the sign of the charge number. Defaults to unsigned (i.e., always positive).

  • Z (real number, optional) – The charge number of an ion or neutral atom, if not provided in particle.

  • mass_numb (integer, optional) – The mass number of an isotope, if not provided in particle.

Returns:

The particle gyrofrequency in units of radians per second.

Return type:

Quantity

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

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

Warns:

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

Notes

The particle gyrofrequency is the angular frequency of particle gyration around magnetic field lines and is given by:

\[ω_{c} = \frac{|Z| e B}{m}\]

If signed is True, then \(|Z|\) is replaced with \(Z\). A particle’s gyrofrequency is also known as its cyclotron frequency or Larmor frequency.

The recommended way to convert from angular frequency to frequency is to use an equivalency between cycles per second and hertz, as Astropy’s dimensionless_angles equivalency does not account for the factor of \(2π\) needed during this conversion. The dimensionless_angles equivalency is appropriate when dividing a velocity by an angular frequency to get a length scale.

Examples

>>> import astropy.units as u
>>> gyrofrequency(0.1 * u.T, "e-")
<Quantity 1.7588...e+10 rad / s>
>>> gyrofrequency(0.1 * u.T, "e-", to_hz=True)
<Quantity 2.79924...e+09 Hz>
>>> gyrofrequency(0.1 * u.T, "e-", signed=True)
<Quantity -1.75882...e+10 rad / s>
>>> gyrofrequency(0.01 * u.T, "p+")
<Quantity 957883.32... rad / s>
>>> gyrofrequency(0.01 * u.T, "p+", signed=True)
<Quantity 957883.32... rad / s>
>>> gyrofrequency(0.01 * u.T, particle="T+")
<Quantity 319964.5... rad / s>
>>> gyrofrequency(0.01 * u.T, particle="T+", to_hz=True)
<Quantity 50923.9... Hz>
>>> gyrofrequency(250 * u.G, particle="Fe", mass_numb=56, Z=13)
<Quantity 560682.34875287 rad / s>
>>> omega_ce = gyrofrequency(0.1 * u.T, "e-")
>>> print(omega_ce)
1758820... rad / s
>>> f_ce = omega_ce.to(u.Hz, equivalencies=[(u.cy / u.s, u.Hz)])
>>> print(f_ce)
279924... Hz
Parameters:

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