cold_plasma_permittivity_LRP

plasmapy.formulary.dielectric.cold_plasma_permittivity_LRP(
B: Quantity,
species: ParticleList | Sequence[str | int | integer | Particle | CustomParticle | Quantity],
n: list[Quantity] | Quantity,
omega: Quantity,
) RotatingTensorElements[source]

Magnetized cold plasma dielectric permittivity tensor elements.

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

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

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

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

  • n ((k,) list of 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 (Quantity) – L (“Left”) Left-handed circularly polarization tensor element.

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

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

Notes

In the rotating frame defined by \((\mathbf{u}_{+}, \mathbf{u}_{-}, \mathbf{u}_z)\) with \(\mathbf{u}_{±}=(\mathbf{u}_x ± \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{ω_{p,s}^2}{ω\left(ω - Ω_{c,s}\right)}\\R = 1 - \sum_s \frac{ω_{p,s}^2}{ω\left(ω + Ω_{c,s}\right)}\\P = 1 - \sum_s \frac{ω_{p,s}^2}{ω^2}\end{aligned}\end{align} \]

where \(ω_{p,s}\) is the plasma frequency and \(Ω_{c,s}\) is the signed version of the cyclotron frequency for the species \(s\) [Stix, 1992].

Examples

>>> import astropy.units as u
>>> from numpy 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.63333...>
>>> permittivity.left  # namedtuple-style access
<Quantity 0.63333...>
>>> R
<Quantity 1.41512...>
>>> P
<Quantity -4.8903...>