AbstractMHDWave

class plasmapy.dispersion.analytical.mhd_waves_.AbstractMHDWave(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)[source]

Bases: ABC

Abstract base class for magnetohydrodynamic waves.

Attributes Summary

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.

sound_speed

The sound speed of the plasma.

Methods Summary

angular_frequency(k, theta)

Calculate the angular frequency of magnetohydrodynamic waves.

group_velocity(k, theta)

Calculate the group velocities of magnetohydrodynamic 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

abstract angular_frequency(k: Quantity, theta: Quantity) Quantity[source]

Calculate the angular frequency of magnetohydrodynamic waves.

Parameters:
  • k (Quantity) – 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:

omega – An \(N × M\) array of computed wave frequencies in units rad/s.

Return type:

Quantity

Raises:
  • UnitTypeError – If applicable arguments do not have units convertible to the expected units.

  • ValueError – If k is negative or zero.

  • ValueError – If k or theta 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.

abstract group_velocity(k: Quantity, theta: Quantity) Quantity[source]

Calculate the group velocities of magnetohydrodynamic waves.

Parameters:
  • k (Quantity) – 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 increasing theta, the second dimension maps to the k array, and the third dimension maps to the theta 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 or theta 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{domega}{dmathbf{k}}
= hat{mathbf{k}} frac{partialomega}{partial k}
  • hat{mathbf{theta}} frac{partial v_{ph}}{partialtheta}

where \(ω\) is the angular frequency, \(\mathbf{k}\) is the wavevector, \(θ\) is the angle between \(\mathbf{k}\) and the unperturbed magnetic field, and \(v_{ph}\) is the phase velocity.

phase_velocity(k: Quantity, theta: Quantity) Quantity[source]

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:

Quantity

Raises:
  • UnitTypeError – If applicable arguments do not have units convertible to the expected units.

  • ValueError – If k is negative or zero.

  • ValueError – If k or theta 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.