cold_plasma_permittivity_SDP

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

Magnetized cold plasma dielectric permittivity tensor elements.

Elements (S, D, P) are given in the “Stix” frame, i.e. with \(B ∥ \hat{z}\) [Stix, 1992].

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

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

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

  • n (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:

  • sum (Quantity) – The “sum” dielectric tensor element, \(S\).

  • difference (Quantity) – The “difference” dielectric tensor element, \(D\).

  • plasma (Quantity) – The “plasma” dielectric tensor element, \(P\).

Notes

The dielectric permittivity tensor is expressed in the Stix frame with the \(\exp(-i ω t)\) time-harmonic convention as \(ε = ε_0 A\), with \(A\) being

\[\begin{split}ε = ε_0 \left(\begin{matrix} S & -i D & 0 \\ +i D & S & 0 \\ 0 & 0 & P \end{matrix}\right)\end{split}\]

where:

\[ \begin{align}\begin{aligned}S = 1 - \sum_s \frac{ω_{p,s}^2}{ω^2 - Ω_{c,s}^2}\\D = \sum_s \frac{Ω_{c,s}}{ω} \frac{ω_{p,s}^2}{ω^2 - Ω_{c,s}^2}\\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\).

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)
>>> permittivity = S, D, P = cold_plasma_permittivity_SDP(B, species, n, omega)
>>> S
<Quantity 1.02422...>
>>> permittivity.sum   # namedtuple-style access
<Quantity 1.02422...>
>>> D
<Quantity 0.39089...>
>>> P
<Quantity -4.8903...>