permittivity_1D_Maxwellian

plasmapy.formulary.dielectric.permittivity_1D_Maxwellian(omega: Quantity, kWave: Quantity, T: Quantity, n: Quantity, particle, z_mean=None) Quantity[source]

Compute the classical dielectric permittivity for a 1D Maxwellian plasma.

This function can calculate both the ion and electron permittivities. No additional effects are considered (e.g. magnetic fields, relativistic effects, strongly coupled regime, etc.).

Parameters:
  • omega (Quantity) – The frequency, in rad/s, of the electromagnetic wave propagating through the plasma.

  • kWave (Quantity) – The corresponding wavenumber, in rad/m, of the electromagnetic wave propagating through the plasma.

  • T (Quantity) – The plasma temperature — this can be either the electron or the ion temperature, but should be consistent with density and particle.

  • n (Quantity) – The plasma density — this can be either the electron or the ion density, but should be consistent with temperature and particle.

  • particle (str) – The plasma particle species.

  • z_mean (Real) – The average ionization of the plasma. This is only required for calculating the ion permittivity.

Returns:

chi – The ion or the electron dielectric permittivity of the plasma. This is a dimensionless quantity.

Return type:

Quantity

Notes

The dielectric permittivities for a Maxwellian plasma are described by the following equations (see p. 106 of Froula et al. [2011]):

\[ \begin{align}\begin{aligned}χ_e(k, ω) = - \frac{α_e^2}{2} Z'(x_e)\\χ_i(k, ω) = - \frac{α_i^2}{2}\frac{Z}{} Z'(x_i)\\α = \frac{ω_p}{k v_{Th}}\\x = \frac{ω}{k v_{Th}}\end{aligned}\end{align} \]

\(χ_e\) and \(χ_i\) are the electron and ion permittivities, respectively. \(Z'\) is the derivative of the plasma dispersion function. \(α\) is the scattering parameter which delineates the difference between the collective and non-collective Thomson scattering regimes. \(x\) is the dimensionless phase velocity of the electromagnetic wave propagating through the plasma.

Examples

>>> import astropy.units as u
>>> from numpy import pi
>>> from plasmapy.formulary import thermal_speed
>>> T = 30 * 11600 * u.K
>>> n = 1e18 * u.cm**-3
>>> particle = "Ne"
>>> Z = 8
>>> vth = thermal_speed(T, particle, method="most_probable")
>>> omega = 5.635e14 * 2 * pi * u.rad / u.s
>>> k_wave = omega / vth
>>> permittivity_1D_Maxwellian(omega, k_wave, T, n, particle, Z)
<Quantity -6.72955...e-08+5.76163...e-07j>

For user convenience permittivity_1D_Maxwellian_lite is bound to this function and can be used as follows:

>>> from plasmapy.formulary import plasma_frequency
>>> wp = plasma_frequency(n, particle, Z=Z)
>>> permittivity_1D_Maxwellian.lite(
...     omega.value, k_wave.value, vth=vth.value, wp=wp.value
... )
(-6.72955...e-08+5.76163...e-07j)