impact_parameter

plasmapy.formulary.collisions.lengths.impact_parameter(
T: ~astropy.units.quantity.Quantity,
n_e: ~astropy.units.quantity.Quantity,
species,
z_mean: float = nan,
V: ~astropy.units.quantity.Quantity = <Quantity nan m / s>,
method: str = 'classical',
)[source]

Impact parameters for classical and quantum Coulomb collision.

Parameters:
  • T (Quantity) – Temperature in units of temperature or energy per particle, which is assumed to be equal for both the test particle and the target particle.

  • n_e (Quantity) – The electron number density in units convertible to per cubic meter.

  • species (tuple) – A tuple containing string representations of the test particle (listed first) and the target particle (listed second).

  • z_mean (Quantity, optional) – The average ionization (arithmetic mean) of a plasma for which a macroscopic description is valid. This parameter is used to compute the average ion density (given the average ionization and electron density) for calculating the ion sphere radius for non-classical impact parameters. z_mean is a required parameter if method is "ls_full_interp", "hls_max_interp", or "hls_full_interp".

  • V (Quantity, optional) – The relative velocity between particles. If not provided, thermal velocity is assumed: \(μ V^2 \sim 2 k_B T\) where \(μ\) is the reduced mass.

  • method (str, optional) – The method by which to compute the Coulomb logarithm. The default method is the classical straight-line Landau-Spitzer method ("classical" or "ls"). The other 6 supported methods are "ls_min_interp", "ls_full_interp", "ls_clamp_mininterp", "hls_min_interp", "hls_max_interp", and "hls_full_interp". Please refer to the docstring of Coulomb_logarithm for more information about these methods.

Returns:

bmin, bmax – The minimum and maximum impact parameters (distances) for a Coulomb collision.

Return type:

tuple of floats

Raises:
  • ValueError – If the mass or charge of either particle cannot be found, or any of the inputs contain incorrect values.

  • UnitConversionError – If the units on any of the inputs are incorrect.

  • TypeError – If any of n_e, T, or V is not a Quantity.

  • RelativityError – If the input velocity is same or greater than the speed of light.

Warns:
  • UnitsWarning – If units are not provided, SI units are assumed.

  • RelativityWarning – If the input velocity is greater than 5% of the speed of light.

Notes

The minimum and maximum impact parameters may be calculated in a variety of ways. The maximum impact parameter is typically the Debye length.

For quantum plasmas the maximum impact parameter can be the quadratic sum of the debye length and ion radius (Wigner_Seitz_radius) [Gericke et al., 2002]:

\[b_{max} = \left(λ_{De}^2 + a_i^2\right)^{1/2}\]

The minimum impact parameter is typically some combination of the thermal de Broglie wavelength and the distance of closest approach for a 90° Coulomb collision. A quadratic sum is used for all GMS methods, except for GMS-5, where b_min is simply set to the distance of closest approach [Gericke et al., 2002].

\[b_{min} = \left(Λ_{de Broglie}^2 + ρ_⟂^2\right)^{1/2}\]

Examples

>>> import astropy.units as u
>>> n = 1e19 * u.m**-3
>>> T = 1e6 * u.K
>>> species = ("e", "p")
>>> impact_parameter(T, n, species)
(<Quantity 1.051...e-11 m>, <Quantity 2.182...e-05 m>)
>>> impact_parameter(T, n, species, V=1e6 * u.m / u.s)
(<Quantity 2.534...e-10 m>, <Quantity 2.182...e-05 m>)