kinetic_alfven
- plasmapy.dispersion.numerical.kinetic_alfven_.kinetic_alfven(
- 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,
Using the equation provided in Bellan [2012], this function calculates the numerical solution to the kinetic Alfvén dispersion relation presented by Hirose et al. [2004].
- 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. Either single valued or 1-D array of length \(N\).n_i (
Quantity
) – Ion number density in units convertible to m-3.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.theta (
Quantity
) – The angle of propagation of the wave with respect to the magnetic field, \(\cos^{-1}(k_z / k)\), in units convertible to rad. Either single valued or 1-D array of size \(M\).gamma_e (real number, 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 (real number, 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 (integer, keyword-only, optional) – The mass number 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 a key for each:
theta
value provided. The value for each key will be an \(N × M\) array.- Return type:
Dict[str,
Quantity
]- Raises:
TypeError – If applicable arguments are not instances of
Quantity
or cannot be converted into one.TypeError – If
ion
is not particle-like.TypeError – If
gamma_e
,gamma_i
, orZ
are not real numbers.UnitTypeError – If applicable arguments do not have units convertible to the expected units.
ValueError – If any of
B
,k
,n_i
,T_e
, orT_i
is negative.ValueError – If
k
is negative or zero.ValueError – If
ion
is not of category ion or element.ValueError – If
B
,n_i
,T_e
, orT_i
are not single valuedastropy.units.Quantity
(i.e. an array).ValueError – If
k
ortheta
are not single valued or a 1-D array.
Notes
Using the 2 × 2 matrix approach method from Bellan [2012], this function computes the corresponding wave frequencies in units of rad / s. This approach comes from Hasegawa and Uberoi [1982], Morales and Maggs [1997] and Lysak and Lotko [1996]; who argued that a 3 × 3 matrix that describes warm plasma waves can be represented as a 2 × 2 matrix because the compressional (i.e., fast) mode can be factored out. The result is that the determinant, when in the limit of \(ω ≫ k_{z}^{2} c^{2}_{\rm s}\), reduces to the kinetic Alfvén dispersion relation.
\[ω^2 = k_{\rm z}^2 v_{\rm A}^2 \left(1 + \frac{k_{\rm x}^2 c_{\rm s}^2}{ω_{\rm ci}^2} \right)\]With \(c_{\rm s}\) being the wave speed and \(ω_{\rm ci}\) as the gyrofrequency of the respective ion. The regions in which this is valid are \(ω ≪ ω_{\rm ci}\) and \(ν_{\rm Te} ≫ \frac{ω}{k_{z}} ≫ ν_{\rm Ti}\), with \(ν_{\rm Ti}\) standing for the thermal speed of the respective ion. There is no restriction on propagation angle.
Examples
>>> import numpy as np >>> import astropy.units as u >>> from plasmapy.particles import Particle >>> from plasmapy.dispersion.numerical import kinetic_alfven_ >>> inputs = { ... "B": 8.3e-9 * u.T, ... "ion": Particle("p+"), ... "k": np.logspace(-7, -2, 2) * u.rad / u.m, ... "n_i": 5 * u.m**-3, ... "T_e": 1.6e6 * u.K, ... "T_i": 4.0e5 * u.K, ... "theta": 30 * u.deg, ... "gamma_e": 3, ... "gamma_i": 3, ... "Z": 1, ... } >>> kinetic_alfven(**inputs) {30.0: <Quantity [1.24901116e+00, 3.45301796e+08] rad / s>}