stix
- plasmapy.dispersion.analytical.stix_.stix( )[source]
Calculate the cold plasma dispersion function presented by Stix [1992], and discussed by Bellan [2012]. This is an analytical solution of equation 8 in Bellan [2012]. See the Notes section below for additional details.
- Parameters:
B (
Quantity
) – The magnetic field magnitude in units convertible to T.w (
Quantity
, single valued or 1-D array) – Wavefrequency in units convertible to rad/s. Either singled valued or 1-D array of length \(N\).ions (a single or
list
of particle-like object(s)) – A list or single instance of particle-like objects representing the ion species (e.g.,"p"
for protons,"D+"
for deuterium,["H+", "He+"]
for hydrogen and helium, etc.). All ions must be positively charged.n_i (
Quantity
, single valued or 1-D array) – Ion number density in units convertible to m-3. Must be single valued or equal length toions
.theta (
Quantity
, single valued or 1-D array) – The angle of propagation of the wave with respect to the magnetic field, \(\cos^{-1}(k_z / k)\), in units convertible to radians. Either single valued or 1-D array of size \(M\).
- Returns:
k – An array of wavenubmers in units rad/m (shape \(N × M × 4\)). The first dimension maps to the
w
array, the second dimension maps to thetheta
array, and the third dimension maps to the four roots of the Stix polynomial.k[0]
is the square root of the positive quadratic solutionk[1] = -k[0]
k[2]
is the square root of the negative quadratic solutionk[3] = -k[2]
- Return type:
Quantity
of shape(N, M, 4)
- Raises:
TypeError – If applicable arguments are not instances of
Quantity
or cannot be converted into one.ValueError – Particles in
ions
are not positively charged ions.ValueError – The size of
n_i
is not the same as the length ofions
.ValueError – If of
B
orn_i
is negative.ValueError – If
w
is negative or zero.ValueError – If
w
ortheta
are not single valued or a 1-D array.UnitTypeError – If applicable arguments do not have units convertible to the expected units.
Notes
The cold plasma dispersion function is defined by Stix [1992] in section 1-3 (and present by Bellan [2012] in equation 8) to be
\[(S\sin^{2}(θ) + P\cos^{2}(θ))(ck/ω)^{4} - [ RL\sin^{2}(θ) + PS(1 + \cos^{2}(θ)) ](ck/ω)^{2} + PRL = 0\]where,
\[\begin{split}\mathbf{B_o} &= B_{o} \mathbf{\hat{z}} \\ \cos θ &= \frac{k_z}{k} \\ \mathbf{k} &= k_{\rm x} \hat{x} + k_{\rm z} \hat{z}\end{split}\]\[\begin{split}S &= 1 - \sum_{s} \frac{ω^{2}_{p,s}}{ω^{2} - ω^{2}_{c,s}}\\ P &= 1 - \sum_{s} \frac{ω^{2}_{p,s}}{ω^{2}}\\ D &= \sum_{s} \frac{ω_{c,s}}{ω} \frac{ω^{2}_{p,s}}{ω^{2} - ω_{c,s}^{2}}\end{split}\]\[R = S + D \hspace{1cm} L = S - D\]\(ω\) is the wave frequency, \(k\) is the wavenumber, \(θ\) is the wave propagation angle with respect to the background magntic field \(\mathbf{B_o}\), \(s\) corresponds to plasma species \(s\), \(ω_{p,s}\) is the plasma frequency of species \(s\), and \(ω_{c,s}\) is the gyrofrequency of species \(s\). The derivation of this dispersion relation assumed:
zero temperature for all plasma species (\(T_{s}=0\))
quasi-neutrality
a uniform background magntic field \(\mathbf{B_o} = B_{o} \mathbf{\hat{z}}\)
no D.C. electric field \(\mathbf{E_o}=0\)
zero-order quantities for all plasma parameters (densities, electric-field, magnetic field, particle speeds, etc.) are constant in time and space
first-order perturbations in plasma parameters vary like \(\sim e^{\left [ i (\textbf{k}\cdot\textbf{r} - ω t)\right ]}\)
Due to the cold plasma assumption, this equation is valid for all \(ω\) and \(k\) given \(\frac{ω}{k_{z}} \gg v_{Th}\) for all thermal speeds \(v_{Th}\) of all plasma species and \(k_{x} r_{L} \ll 1\) for all gyroradii \(r_{L}\) of all plasma species.
The relation predicts \(k \to 0\) when any one of P, R or L vanish (cutoffs) and \(k \to ∞\) for perpendicular propagation during wave resonance \(S \to 0\).
Examples
>>> import astropy.units as u >>> from plasmapy.particles import Particle >>> from plasmapy.dispersion.analytical.stix_ import stix >>> inputs = { ... "B": 8.3e-9 * u.T, ... "w": 0.001 * u.rad / u.s, ... "ions": [Particle("H+"), Particle("He+")], ... "n_i": [4.0e5,2.0e5] * u.m**-3, ... "theta": 30 * u.deg, ... } >>> stix(**inputs) <Quantity [ 6.03817661e-09-0.j, -6.03817661e-09+0.j, 6.97262784e-09-0.j, -6.97262784e-09+0.j] rad / m>