plasma_frequency¶
-
plasmapy.formulary.parameters.
plasma_frequency
(n: Unit("1 / m3"), particle: plasmapy.particles.particle_class.Particle, z_mean=None, to_hz=False) -> Unit("rad / s")¶ Calculate the particle plasma frequency.
Aliases:
wp_
Parameters: - n (Quantity) – Particle number density in units convertible to per cubic meter
- 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.
- z_mean (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 (
int
) of the ion is used. This is effectively an average plasma frequency for the plasma where multiple charge states are present.
Returns: omega_p – The particle plasma frequency in radians per second.
Return type: Raises: TypeError
– If n_i is not aQuantity
or particle is not of an appropriate type.UnitConversionError
– Ifn_i
is not in correct unitsValueError
– Ifn_i
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 plasma frequency is
\[\omega_{pi} = Z e \sqrt{\frac{n_i}{\epsilon_0 m_i}}\]At present, astropy.units does not allow direct conversions from radians/second for angular frequency to 1/second or Hz for frequency. The dimensionless_angles equivalency allows that conversion, but does not account for the factor of 2*pi. The alternatives are to convert to cycle/second or to do the conversion manually, as shown in the examples.
Example
>>> from astropy import units as u >>> plasma_frequency(1e19*u.m**-3, particle='p') <Quantity 4.16329...e+09 rad / s> >>> plasma_frequency(1e19*u.m**-3, particle='p', to_hz=True) <Quantity 6.62608...e+08 Hz> >>> plasma_frequency(1e19*u.m**-3, particle='D+') <Quantity 2.94462...e+09 rad / s> >>> plasma_frequency(1e19*u.m**-3, 'e-') <Quantity 1.78398...e+11 rad / s> >>> plasma_frequency(1e19*u.m**-3, 'e-', to_hz=True) <Quantity 2.83930...e+10 Hz>
Other Parameters: to_hz (bool) – Set True
to to convert function output from angular frequency to Hz