cold_plasma_permittivity_LRP

plasmapy.physics.dielectric.cold_plasma_permittivity_LRP(B: Unit("T"), species, n, omega: Unit("rad / s"))

Magnetized Cold Plasma Dielectric Permittivity Tensor Elements.

Elements (L, R, P) are given in the “rotating” basis, ie. in the basis \((\mathbf{u}_{+}, \mathbf{u}_{-}, \mathbf{u}_z)\), where the tensor is diagonal and with B // z.

The \(\exp(-i \omega t)\) time-harmonic convention is assumed.

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

  • species (list of str) – The plasma particle species (e.g.: ['e', 'D+'] or ['e', 'D+', 'He+'].

  • n (list of ~astropy.units.Quantity) – list of species density in units convertible to per cubic meter. The order of the species densities should follow species.

  • omega (Quantity) – Electromagnetic wave frequency in rad/s.

Returns

  • left (~astropy.units.Quantity) – L (“Left”) Left-handed circularly polarization tensor element.

  • right (~astropy.units.Quantity) – R (“Right”) Right-handed circularly polarization tensor element.

  • plasma (~astropy.units.Quantity) – P (“Plasma”) dielectric tensor element.

Notes

In the rotating frame defined by \((\mathbf{u}_{+}, \mathbf{u}_{-}, \mathbf{u}_z)\) with \(\mathbf{u}_{\pm}=(\mathbf{u}_x \pm \mathbf{u}_y)/\sqrt{2}\), the dielectric tensor takes a diagonal form with elements L, R, P with:

\[ \begin{align}\begin{aligned}L = 1 - \sum_s \frac{\omega_{p,s}^2}{\omega\left(\omega - \Omega_{c,s}\right)}\\R = 1 - \sum_s \frac{\omega_{p,s}^2}{\omega\left(\omega + \Omega_{c,s}\right)}\\P = 1 - \sum_s \frac{\omega_{p,s}^2}{\omega^2}\end{aligned}\end{align} \]

where \(\omega_{p,s}\) is the plasma frequency and \(\Omega_{c,s}\) is the signed version of the cyclotron frequency for the species \(s\).

References

  • T.H. Stix, Waves in Plasma, 1992.

Examples

>>> from astropy import units as u
>>> from plasmapy.constants import pi
>>> B = 2*u.T
>>> species = ['e', 'D+']
>>> n = [1e18*u.m**-3, 1e18*u.m**-3]
>>> omega = 3.7e9*(2*pi)*(u.rad/u.s)
>>> L, R, P = permittivity = cold_plasma_permittivity_LRP(B, species, n, omega)
>>> L
<Quantity 0.63333549>
>>> permittivity.left    # namedtuple-style access
<Quantity 0.63333549>
>>> R
<Quantity 1.41512254>
>>> P
<Quantity -4.8903104>