two_fluid

plasmapy.dispersion.analytical.two_fluid_.two_fluid(
B: Quantity,
ion: str | int | integer | Particle | CustomParticle | Quantity,
k: Quantity,
n_i: Quantity,
theta: Quantity,
*,
T_e: Quantity,
T_i: Quantity,
gamma_e: float = 1,
gamma_i: float = 3,
mass_numb: int | None = None,
Z: float | None = None,
)[source]

Using the solution provided by Bellan [2012], calculate the analytical solution to the two fluid, low-frequency (\(ω/kc ≪ 1\)) dispersion relation presented by Stringer [1963]. This dispersion relation also assumes a uniform magnetic field \(\mathbf{B}_0\), no D.C. electric field \(\mathbf{E}_0=0\), and quasi-neutrality. For more information see the Notes section below.

Parameters:
  • B (Quantity) – The magnetic field magnitude in units convertible to T.

  • ion (particle-like) – Representation of the ion species (e.g., 'p' for protons, 'D+' for deuterium, 'He-4 1+' for singly ionized helium-4, etc.).

  • k (Quantity) – Wavenumber in units convertible to rad/m. May be either single valued or a 1D array of length \(N\).

  • n_i (Quantity) – Ion number density in units convertible to m-3.

  • theta (Quantity) – The angle of propagation of the wave with respect to the magnetic field, \(\cos^{-1}(k_z / k)\), in units must be convertible to radians. May be either single valued or a 1D array of size \(M\).

  • T_e (Quantity, keyword-only) – The electron temperature in units of K or eV.

  • T_i (Quantity, keyword-only) – The ion temperature in units of K or eV.

  • gamma_e (float or int, keyword-only, default: 1) – The adiabatic index for electrons. The default value assumes that the electrons are able to equalize their temperature rapidly enough that the electrons are effectively isothermal.

  • gamma_i (float or int, keyword-only, default: 3) – The adiabatic index for ions. The default value assumes that ion motion has only one degree of freedom, namely along magnetic field lines.

  • mass_numb (int, keyword-only, optional) – The mass number of an isotope corresponding to ion.

  • Z (real number, keyword-only, optional) – The charge number corresponding to ion.

Returns:

omega – A dictionary of computed wave frequencies in units rad/s. The dictionary contains three keys: 'fast_mode' for the fast mode, 'alfven_mode' for the Alfvén mode, and 'acoustic_mode' for the ion-acoustic mode. The value for each key will be a \(N × M\) array.

Return type:

Dict[str, Quantity]

Raises:
Warns:

PhysicsWarning – When the computed wave frequencies violate the low-frequency (\(ω/kc ≪ 1\)) assumption of the dispersion relation.

Notes

The complete dispersion equation presented by Stringer [1963] (equation 1 of Bellan [2012]) is:

\[\begin{split}\left( \cos^2 θ - Q \frac{ω^2}{k^2 {v_A}^2} \right) & \left[ \left( \cos^2 θ - \frac{ω^2}{k^2 {c_s}^2} \right) - Q \frac{ω^2}{k^2 {v_A}^2} \left( 1 - \frac{ω^2}{k^2 {c_s}^2} \right) \right] \\ &= \left(1 - \frac{ω^2}{k^2 {c_s}^2} \right) \frac{ω^2}{{ω_{ci}}^2} \cos^2 θ\end{split}\]

where

\[\begin{split}Q &= 1 + k^2 c^2/{ω_{pe}}^2 \\ \cos θ &= \frac{k_z}{k} \\ \mathbf{B}_0 &= B_0 \mathbf{\hat{z}}\end{split}\]

\(ω\) is the wave frequency, \(k\) is the wavenumber, \(v_A\) is the Alfvén velocity, \(c_s\) is the sound speed, \(ω_{ci}\) is the ion gyrofrequency, and \(ω_{pe}\) is the electron plasma frequency. This relation does additionally assume low-frequency waves \(ω/kc ≪ 1\), no D.C. electric field \(\mathbf{E}_0=0\) and quasi-neutrality.

Following section 5 of Bellan [2012], the exact roots of the above dispersion equation can be derived and expressed as one analytical solution (equation 38 of Bellan [2012]):

\[\frac{ω}{ω_{ci}} = \sqrt{ 2 Λ \sqrt{-\frac{P}{3}} \cos\left( \frac{1}{3} \cos^{-1}\left( \frac{3q}{2p} \sqrt{-\frac{3}{p}} \right) - \frac{2π}{3}j \right) + \frac{Λ A}{3} }\]

where \(j = 0\) represents the fast mode, \(j = 1\) represents the Alfvén mode, and \(j = 2\) represents the acoustic mode. Additionally,

\[\begin{split}p &= \frac{3B-A^2}{3} \; , \; q = \frac{9AB-2A^3-27C}{27} \\ A &= \frac{Q + Q^2 β + Q α + α Λ}{Q^2} \; , \; B = α \frac{1 + 2 Q β + Λ β}{Q^2} \; , \; C = \frac{α^2 β}{Q^2} \\ α &= \cos^2 θ \; , \; β = \left( \frac{c_s}{v_A}\right)^2 \; , \; Λ = \left( \frac{k v_{A}}{ω_{ci}}\right)^2\end{split}\]

Examples

>>> import astropy.units as u
>>> from plasmapy.dispersion.analytical import two_fluid
>>> inputs = {
...     "k": 0.01 * u.rad / u.m,
...     "theta": 30 * u.deg,
...     "B": 8.3e-9 * u.T,
...     "n_i": 5e6 * u.m ** -3,
...     "T_e": 1.6e6 * u.K,
...     "T_i": 4.0e5 * u.K,
...     "ion": "p+",
... }
>>> omegas = two_fluid(**inputs)
>>> omegas
{'fast_mode': <Quantity 1520.57... rad / s>,
 'alfven_mode': <Quantity 1261.75... rad / s>,
 'acoustic_mode': <Quantity 0.688152... rad / s>}
>>> inputs = {
...     "k": [1e-7, 2e-7] * u.rad / u.m,
...     "theta": [10, 20] * u.deg,
...     "B": 8.3e-9 * u.T,
...     "n_i": 5e6 * u.m ** -3,
...     "T_e": 1.6e6 * u.K,
...     "T_i": 4.0e5 * u.K,
...     "ion": "He+",
... }
>>> omegas = two_fluid(**inputs)
>>> omegas['fast_mode']
<Quantity [[0.00767..., 0.00779... ],
           [0.01534..., 0.01558...]] rad / s>