FastMagnetosonicWave

class plasmapy.dispersion.analytical.mhd_waves_.FastMagnetosonicWave(B: ~astropy.units.quantity.Quantity, density: (Unit("1 / m3"), Unit("kg / m3")), ion: str | ~numbers.Integral | ~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: ~numbers.Integral | None = None, Z: ~numbers.Real | None = None)[source]

Bases: AbstractMHDWave

A class to represent fast magnetosonic 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 or int, keyword-only, default: 5/3) – The adiabatic index for the plasma.

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

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

Raises:

Notes

Fast magnetosonic waves are also referred to as fast magnetoacoustic waves.

Examples

>>> import astropy.units as u
>>> from plasmapy.dispersion.analytical import FastMagnetosonicWave
>>> fast = FastMagnetosonicWave(1e-3 * u.T, 1e16 * u.m**-3, "p+", T=2.5e6 * u.K)
>>> fast.angular_frequency(1e-5 * u.rad / u.m, 0 * u.deg)
<Quantity 2.18060973 rad / s>
>>> fast.phase_velocity(1e-5 * u.rad / u.m, 0 * u.deg)
<Quantity 218060.97295233 m / s>
>>> fast.alfven_speed
<Quantity 218060.97295233 m / s>

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 a fast magnetosonic waves.

group_velocity(k, theta)

Calculate the group velocities of fast magnetosonic 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(k: Quantity, theta: Quantity)[source]

Calculate the angular frequency of a fast magnetosonic 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.

Notes

The angular frequency \(ω\) of a fast magnetosonic wave is given by the equation

\[ω^2 = \frac{k^2}{2} \left( c_{ms}^2 + \sqrt{c_{ms}^4 - 4 v_A^2 c_s^2 \cos^2 θ} \right)\]

where \(k\) is the wavenumber, \(v_A\) is the Alfvén speed, \(c_s\) is the ideal sound speed, \(c_{ms} = \sqrt{v_A^2 + c_s^2}\) is the magnetosonic speed, and \(θ\) is the angle between the wavevector and the equilibrium magnetic field.

Examples

>>> import astropy.units as u
>>> from plasmapy.dispersion.analytical import FastMagnetosonicWave
>>> fast = FastMagnetosonicWave(1e-3 * u.T, 1e16 * u.m**-3, "p+", T=2.5e6 * u.K)
>>> fast.angular_frequency(1e-5 * u.rad / u.m, 0 * u.deg)
<Quantity 2.18060973 rad / s>
>>> fast.angular_frequency([1e-5, 2e-4] * (u.rad / u.m), 0 * u.deg)
<Quantity [ 2.18060973, 43.61219459] rad / s>
>>> fast.angular_frequency(1e-5 * u.rad / u.m, [0, 45, 90] * u.deg)
<Quantity [2.18060973, 2.65168984, 2.86258485] rad / s>
>>> fast.angular_frequency([1e-5, 2e-4] * (u.rad / u.m), [0, 45, 90] * u.deg)
<Quantity [[ 2.18060973,  2.65168984,  2.86258485],
           [43.61219459, 53.03379678, 57.251697  ]] rad / s>
group_velocity(k: Quantity, theta: Quantity)[source]

Calculate the group velocities of fast magnetosonic 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 \times N \times 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}} v_{ph}
  • 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

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.