hollweg
- plasmapy.dispersion.numerical.hollweg_.hollweg(
- 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,
Calculate the two-fluid dispersion relation presented by Hollweg [1999], and discussed by Bellan [2012].
This is a numerical solver of equation 3 in Bellan [2012]. See the Notes section below for additional details.
- Parameters:
B (
Quantity
) – The magnetic field magnitude in units convertible to tesla.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.). If no charge state information is provided, then the ions are assumed to be singly ionized.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.theta (
Quantity
) – 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\).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 (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 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 an \(N x 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 of type or convertible toParticle
.TypeError – If
gamma_e
,gamma_i
, orz_mean
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.
- Warns:
PhysicsWarning
– When \(ω / ω_{\rm ci} > 0.1\), violation of the low-frequency assumption.PhysicsWarning
– When \(c_{\rm s} / v_{\rm A} > 0.1\), violation of low-β.PhysicsWarning
– When \(|θ - π/2| > 0.1\), violation of quasi-perpendicular propagation.
Notes
The dispersion relation presented in Hollweg [1999] (equation 3 in Bellan [2012]) is:
\[\begin{split}\left( \frac{ω^2}{k_{\rm z}^2 v_{\rm A}^2} - 1 \right) & \left[ ω^2 \left( ω^2 - k^2 v_{\rm A}^2 \right) - β k^2 v_{\rm A}^2 \left( ω^2 - k_{\rm z}^2 v_{\rm A}^2 \right) \right] \\ &= ω^2 \left(ω^2 - k^2 v_{\rm A}^2 \right) k_{\rm x}^2 \left( \frac{c_{\rm s}^2}{ω_{\rm ci}^2} - \frac{c^2}{ω_{\rm pe}^2} \frac{ω^2}{k_{\rm z}^2v_{\rm A}^2} \right)\end{split}\]where
\[\begin{split}\mathbf{B}_0 &= B_0 \mathbf{\hat{z}} \\ \cos θ &= \frac{k_z}{k} \\ \mathbf{k} &= k_{\rm x} \hat{x} + k_{\rm z} \hat{z}\end{split}\]\(ω\) is the wave frequency, \(k\) is the wavenumber, \(v_{\rm A}\) is the Alfvén velocity, \(c_{\rm s}\) is the sound speed, \(ω_{\rm ci}\) is the ion gyrofrequency, and \(ω_{\rm pe}\) is the electron plasma frequency. In the derivation of this relation Hollweg assumed low-frequency waves \(ω / ω_{\rm ci} ≪ 1\), no D.C. electric field \(\mathbf{E}_0=0\), and quasineutrality.
Hollweg [1999] asserts this expression is valid for arbitrary \(c_{\rm s} / v_{\rm A}\) (β) and \(k_{\rm z} / k\) (\(θ\)). Contrarily, Bellan [2012] states in §1.7 that due to the inconsistent retention of the \(ω / ω_{\rm ci} ≪ 1\) terms the expression can only be valid if both \(c_{\rm s} ≪ v_{\rm A}\) (low-β) and the wave propagation is nearly perpendicular to the magnetic field.
This routine solves for \(ω\) for given \(k\) values by numerically solving for the roots of the above expression.
Examples
>>> import astropy.units as u >>> from plasmapy.dispersion.numerical import hollweg_ >>> inputs = { ... "k": np.logspace(-7, -2, 2) * u.rad / u.m, ... "theta": 88 * u.deg, ... "n_i": 5 * u.cm ** -3, ... "B": 2.2e-8 * u.T, ... "T_e": 1.6e6 * u.K, ... "T_i": 4.0e5 * u.K, ... "ion": "p+", ... } >>> omegas = hollweg(**inputs) >>> omegas {'fast_mode': <Quantity [2.62911663e-02+0.j, 2.27876968e+03+0.j] rad / s>, 'alfven_mode': <Quantity [7.48765909e-04+0.j, 2.13800404e+03+0.j] rad / s>, 'acoustic_mode': <Quantity [0.00043295+0.j, 0.07358991+0.j] rad / s>}