AlfvenWave
- class plasmapy.dispersion.analytical.mhd_waves_.AlfvenWave(
- B: ~astropy.units.quantity.Quantity,
- density: (Unit("1 / m3"),
- Unit("kg / m3")),
- ion: str | int | ~numpy.integer | ~plasmapy.particles.particle_class.Particle | ~plasmapy.particles.particle_class.CustomParticle | ~astropy.units.quantity.Quantity,
- *,
- T: ~astropy.units.quantity.Quantity = <Quantity 0. K>,
- gamma: float = 1.6666666666666667,
- mass_numb: int | None = None,
- Z: float | None = None,
Bases:
AbstractMHDWave
A class to represent magnetohydrodynamic Alfvén waves.
- Parameters:
B (
Quantity
) – The magnetic field magnitude in units convertible to T.density (
Quantity
) – Either the ion number density \(n_i\) in units convertible to m-3 or the total mass density \(ρ\) in units convertible to kg m-3.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.T (
Quantity
, keyword-only, optional) – The plasma temperature in units of K or eV, which defaults to zero.gamma (
float
orint
, keyword-only, default: 5/3) – The adiabatic index for the plasma.mass_numb (
int
, keyword-only, optional) – The mass number corresponding toion
.Z (
float
orint
, keyword-only, optional) – The charge number corresponding toion
.
- Raises:
TypeError – If applicable arguments are not instances of
Quantity
or cannot be converted into one.TypeError – If
ion
is not particle-like.UnitTypeError – If applicable arguments do not have units convertible to the expected units.
ValueError – If any of
B
,density
, orT
is negative.ValueError – If
ion
is not of category ion or element.ValueError – If
B
,density
, orT
are not single valuedastropy.units.Quantity
(i.e. an array).
See also
Examples
>>> import astropy.units as u >>> from plasmapy.dispersion.analytical import AlfvenWave >>> alfven = AlfvenWave(1e-3 * u.T, 1e16 * u.m**-3, "p+") >>> alfven.angular_frequency(1e-5 * u.rad / u.m, 0 * u.deg) <Quantity 2.18060973 rad / s> >>> alfven.phase_velocity(1e-5 * u.rad / u.m, 0 * u.deg) <Quantity 218060.97295233 m / s> >>> alfven.alfven_speed <Quantity 218060.97295233 m / s>
Attributes Summary
The Alfvén speed of the plasma.
The ratio of thermal pressure to magnetic pressure.
The magnetosonic speed of the plasma.
The sound speed of the plasma.
Methods Summary
angular_frequency
(k, theta)Calculate the angular frequency of magnetohydrodynamic Alfvén waves.
group_velocity
(k, theta)Calculate the group velocities of magnetohydrodynamic Alfvén waves.
phase_velocity
(k, theta)Calculate the phase velocities of magnetohydrodynamic waves.
Attributes Documentation
- alfven_speed
The Alfvén speed of the plasma.
- beta
The ratio of thermal pressure to magnetic pressure.
- magnetosonic_speed
The magnetosonic speed of the plasma.
Defined as \(c_{ms} = \sqrt{v_A^2 + c_s^2}\) where \(v_A\) is the Alfvén speed and \(c_s\) is the sound speed.
- sound_speed
The sound speed of the plasma.
Defined as \(c_s = \sqrt{γ k_B T / m_i}\) where \(gamma\) is the adiabatic index of the fluid, \(k_B\) is the Boltzmann constant, \(T\) is the temperature of the fluid, and \(m_i\) is the mass of the ion species in the fluid.
Methods Documentation
- angular_frequency( )[source]
Calculate the angular frequency of magnetohydrodynamic Alfvén waves.
- Parameters:
k (
Quantity
, single valued or 1-D array) – Wavenumber in units convertible to rad/m`. Either single valued or 1-D array of length \(N\).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 must be convertible to radians. Either single valued or 1-D array of size \(M\).
- Returns:
omega – An \(N × M\) array of computed wave frequencies in units rad/s.
- Return type:
- Raises:
UnitTypeError – If applicable arguments do not have units convertible to the expected units.
ValueError – If
k
is negative or zero.ValueError – If
k
ortheta
are not single valued or a 1-D array.
- Warns:
PhysicsWarning
– When the computed wave frequencies violate the low-frequency (\(ω ≪ ω_c,ω_p\)) assumption of the dispersion relation.
Notes
The angular frequency \(ω\) of a magnetohydrodynamic Alfvén wave is given by
\[ω = k v_A \cosθ\]where \(k\) is the wavenumber, \(v_A\) is the Alfvén speed, and \(θ\) is the angle between the wavevector and the equilibrium magnetic field.
Examples
>>> import astropy.units as u >>> from plasmapy.dispersion.analytical import AlfvenWave >>> alfven = AlfvenWave(1e-3 * u.T, 1e16 * u.m**-3, "p+") >>> alfven.angular_frequency(1e-5 * u.rad / u.m, 0 * u.deg) <Quantity 2.18060973 rad / s> >>> alfven.angular_frequency([1e-5, 2e-4] * (u.rad / u.m), 0 * u.deg) <Quantity [ 2.18060973, 43.61219459] rad / s> >>> alfven.angular_frequency(1e-5 * u.rad / u.m, [0, 45, 90] * u.deg) <Quantity [2.18060973e+00, 1.54192393e+00, 1.33523836e-16] rad / s> >>> alfven.angular_frequency([1e-5, 2e-4] * (u.rad / u.m), [0, 45, 90] * u.deg) <Quantity [[2.18060973e+00, 1.54192393e+00, 1.33523836e-16], [4.36121946e+01, 3.08384785e+01, 2.67047673e-15]] rad / s>
- group_velocity( )[source]
Calculate the group velocities of magnetohydrodynamic Alfvén waves.
- Parameters:
k (
Quantity
, single valued or 1-D array) – Wavenumber in units convertible to rad/m`. Either single valued or 1-D array of length \(N\).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\).
- Returns:
group_velocity – An array of group_velocities in units m/s with shape \(2 × N × M\). The first dimension maps to the two coordinate arrays in the direction of
k
and in the direction of increasingtheta
, the second dimension maps to thek
array, and the third dimension maps to thetheta
array.- Return type:
Quantity
of shape(2, N, M)
- Raises:
UnitTypeError – If applicable arguments do not have units convertible to the expected units.
ValueError – If
k
is negative or zero.ValueError – If
k
ortheta
are not single valued or a 1-D array.
- Warns:
PhysicsWarning
– When the computed wave frequencies violate the low-frequency (\(ω ≪ ω_c,ω_p\)) assumption of the dispersion relation.
Notes
The group velocity \(\mathbf{v}_g\) is given by
\[\]- mathbf{v}_g = frac{dω}{dmathbf{k}}
= ± hat{mathbf{B}} v_{A}
where \(\hat{\mathbf{B}}\) is the unit vector in the direction of the unperturbed magnetic field and \(v_A\) is the Alfvén speed.
- phase_velocity( ) Quantity
Calculate the phase velocities of magnetohydrodynamic waves.
- Parameters:
k (
Quantity
) – Wavenumber in units convertible to rad/m`. Either single valued or 1-D array of size \(N\).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\).
- Returns:
phase_velocity – An \(N × M\) array of computed phase velocities in units of m/s.
- Return type:
- Raises:
UnitTypeError – If applicable arguments do not have units convertible to the expected units.
ValueError – If
k
is negative or zero.ValueError – If
k
ortheta
are not single valued or a 1-D array.
- Warns:
PhysicsWarning
– When the computed wave frequencies violate the low-frequency (\(ω ≪ ω_c,ω_p\)) assumption of the dispersion relation.