Coulomb_logarithm

plasmapy.physics.transport.collisions.Coulomb_logarithm(T, n_e, particles, z_mean=<Quantity nan>, V=<Quantity nan m / s>, method='classical')

Estimates the Coulomb logarithm.

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 density in units convertible to per cubic meter.
  • particles (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) for a plasma where the a macroscopic description is valid. This is used to recover the average ion density (given the average ionization and electron density) for calculating the ion sphere radius for non-classical impact parameters.
  • V (Quantity, optional) – The relative velocity between particles. If not provided, thermal velocity is assumed: \(\mu V^2 \sim 2 k_B T\) where mu is the reduced mass.
  • method (str, optional) – Selects which theory to use when calculating the Coulomb logarithm. Defaults to classical method.
Returns:

lnLambda – An estimate of the Coulomb logarithm that is accurate to roughly its reciprocal.

Return type:

float or numpy.ndarray

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.

    If the n_e, T, or V are not Quantities.

  • PhysicsError – If the result is smaller than 1.

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

Warns:
  • ~astropy.units.UnitsWarning – If units are not provided, SI units are assumed
  • ~plasmapy.utils.RelativityWarning – If the input velocity is greater than 5% of the speed of light.

Notes

The classical Coulomb logarithm is given by

\[\ln{\Lambda} \equiv \ln\left( \frac{b_{max}}{b_{min}} \right)\]

where \(b_{min}\) and \(b_{max}\) are the inner and outer impact parameters for Coulomb collisions [1].

The outer impact parameter is given by the Debye length: \(b_{min} = \lambda_D\) which is a function of electron temperature and electron density. At distances greater than the Debye length, electric fields from other particles will be screened out due to electrons rearranging themselves.

The choice of inner impact parameter is either the distance of closest approach for a 90 degree Coulomb collision or the thermal deBroglie wavelength, whichever is larger. This is because Coulomb-style collisions cannot occur for impact parameters shorter than the deBroglie wavelength because quantum effects will change the fundamental nature of the collision [2], [3].

Errors associated with the classical Coulomb logarithm are of order its inverse. If the Coulomb logarithm is of order unity, then the assumptions made in the standard analysis of Coulomb collisions are invalid.

For dense plasmas where the classical Coulomb logarithm breaks down there are various extended methods. These can be found in D.O. Gericke et al’s paper, which has a table summarizing the methods [4]. The GMS-1 through GMS-6 methods correspond to the methods found it that table.

It should be noted that GMS-4 thru GMS-6 modify the Coulomb logarithm to the form:

\[\ln{\Lambda} \equiv 0.5 \ln\left(1 + \frac{b_{max}^2}{b_{min}^2} \right)\]

This means the Coulomb logarithm will not break down for Lambda < 0, which occurs for dense, cold plasmas.

plasmapy.physics.transport.collisions.Classical()

classical Landau-Spitzer approach. Fails for large coupling parameter where Lambda can become less than zero.

GMS-1

1st method listed in Table 1 of reference [3] Landau-Spitzer, but with interpolated bmin instead of bmin selected between deBroglie wavelength and distance of closest approach. Fails for large coupling parameter where Lambda can become less than zero.

GMS-2

2nd method listed in Table 1 of reference [3] Another Landau-Spitzer like approach, but now bmax is also being interpolated. The interpolation is between the Debye length and the ion sphere radius, allowing for descriptions of dilute plasmas. Fails for large coupling parameter where Lambda can become less than zero. 3rd method listed in Table 1 of reference [3] classical Landau-Spitzer fails for argument of Coulomb logarithm Lambda < 0, therefore a clamp is placed at Lambda_min = 2

GMS-4

4th method listed in Table 1 of reference [3] Spitzer-like extension to Coulomb logarithm by noting that Coulomb collisions take hyperbolic trajectories. Removes divergence for small bmin issue in classical Landau-Spitzer approach, so bmin can be zero. Also doesn’t break down as Lambda < 0 is now impossible, even when coupling parameter is large.

GMS-5

5th method listed in Table 1 of reference [3] Similar to GMS-4, but setting bmin as distance of closest approach and bmax interpolated between Debye length and ion sphere radius. Lambda < 0 impossible.

GMS-6

6th method listed in Table 1 of reference [3] Similar to GMS-4 and GMS-5, but using interpolation methods for both bmin and bmax.

Examples

>>> from astropy import units as u
>>> n = 1e19*u.m**-3
>>> T = 1e6*u.K
>>> particles = ('e', 'p')
>>> Coulomb_logarithm(T, n, particles)
14.545527226436974
>>> Coulomb_logarithm(T, n, particles, V=1e6 * u.m / u.s)
11.363478214139432

References

[1]Physics of Fully Ionized Gases, L. Spitzer (1962)
[2]Francis, F. Chen. Introduction to plasma physics and controlled fusion 3rd edition. Ch 5 (Springer 2015).
[3]Comparison of Coulomb Collision Rates in the Plasma Physics and Magnetically Confined Fusion Literature, W. Fundamenski and O.E. Garcia, EFDA–JET–R(07)01 (http://www.euro-fusionscipub.org/wp-content/uploads/2014/11/EFDR07001.pdf)
[4]Dense plasma temperature equilibration in the binary collision approximation. D. O. Gericke et. al. PRE, 65, 036418 (2002). DOI: 10.1103/PhysRevE.65.036418