gyrofrequency¶
-
plasmapy.formulary.
gyrofrequency
(B: Unit("T"), particle: plasmapy.particles.particle_class.Particle, signed=False, Z=None, to_hz=False) -> Unit("rad / s")¶ Calculate the particle gyrofrequency in units of radians per second.
Parameters: - B (Quantity) – The magnetic field magnitude in units convertible to tesla.
- particle (Particle) – Representation of the particle 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 particles are assumed to be singly charged.
- signed (bool, optional) – The gyrofrequency can be defined as signed (negative for electron,
positive for ion). Default is
False
(unsigned, i.e. always positive). - Z (float or Quantity, optional) – The average ionization (arithmetic mean) for a plasma where the
a macroscopic description is valid. If this quantity is not
given then the atomic charge state (integer) of the ion
is used. This is effectively an average gyrofrequency for the
plasma where multiple charge states are present, and should
not be interpreted as the gyrofrequency for any single particle.
If not provided, it defaults to the integer charge of the
particle
.
Returns: omega_c – The particle gyrofrequency in units of radians per second
Return type: Raises: TypeError
– If the magnetic field is not aQuantity
or particle is not of an appropriate typeValueError
– If the magnetic field contains invalid values or particle cannot be used to identify an particle or isotope
Warns: ~astropy.units.UnitsWarning – If units are not provided, SI units are assumed
Notes
The particle gyrofrequency is the angular frequency of particle gyration around magnetic field lines and is given by:
\[\omega_{ci} = \frac{Z e B}{m_i}\]The particle gyrofrequency is also known as the particle cyclotron frequency or the particle 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
dimensionles_angles
equivalency does not account for the factor of 2*pi needed during this conversion. Thedimensionless_angles
equivalency is appropriate when dividing a velocity by an angular frequency to get a length scale.Examples
>>> from astropy import 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> >>> 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
Other Parameters: to_hz (bool) – Set True
to to convert function output from angular frequency to Hz