MaxwellianCollisionFrequencies

class plasmapy.formulary.collisions.frequencies.MaxwellianCollisionFrequencies(test_particle: str | int | ~numpy.integer | ~plasmapy.particles.particle_class.Particle | ~plasmapy.particles.particle_class.CustomParticle | ~astropy.units.quantity.Quantity, field_particle: str | int | ~numpy.integer | ~plasmapy.particles.particle_class.Particle | ~plasmapy.particles.particle_class.CustomParticle | ~astropy.units.quantity.Quantity, *, v_drift: ~astropy.units.quantity.Quantity = <Quantity 0. m / s>, T_a: ~astropy.units.quantity.Quantity, n_a: ~astropy.units.quantity.Quantity, T_b: ~astropy.units.quantity.Quantity, n_b: ~astropy.units.quantity.Quantity, Coulomb_log: ~astropy.units.quantity.Quantity)[source]

Bases: object

Compute collision frequencies between two slowly flowing Maxwellian populations.

The condition of “slowly flowing”, outlined by Eq. 2.133 in Callen [n.d.] requires that

\[v_{drift} ≪ \sqrt{ v_{T_a}^2 + v_{T_b}^2 }\]

where \(v_{drift}\) is the relative drift between the two species, \(v_{T_a}\) is the thermal velocity of species \(a\), and \(v_{T_b}\) is the thermal velocity of species \(b\).

Parameters:
  • test_particle (particle-like) – The test particle streaming through a background of field particles.

  • field_particle (particle-like) – The background particle being interacted with.

  • v_drift (Quantity, default: 0 m/s) – The relative drift between the test and field particles.

  • T_a (Quantity) – The temperature of the test particles in units convertible to kelvin or eV per particle.

  • n_a (Quantity) – The number density of the test particles in units convertible to m-3.

  • T_b (Quantity) – The temperature of the background field particles in units convertible to kelvin or eV per particle.

  • n_b (Quantity) – The number density of the background field particles in units convertible to m-3.

  • Coulomb_log (Quantity) – The value of the Coulomb logarithm for the interaction.

Raises:

ValueError – If the specified v_drift and T_a arrays do not have equal size.

Attributes Summary

Lorentz_collision_frequency

The Lorentz collision frequency.

Maxwellian_avg_ei_collision_freq

Average momentum relaxation rate for a slowly flowing Maxwellian distribution of electrons, relative to a population of stationary ions.

Maxwellian_avg_ii_collision_freq

Average momentum relaxation rate for a slowly flowing Maxwellian distribution of ions, relative to a population of stationary ions.

Attributes Documentation

Lorentz_collision_frequency

The Lorentz collision frequency.

The Lorentz collision frequency (see Ch. 5 of Chen [2016]) is given by

\[ν = n σ v \ln{Λ}\]

where \(n\) is the particle number density, \(σ\) is the collisional cross-section, \(v\) is the mean thermal velocity between particle species (see Equation 2.133 in Callen [n.d.]), and \(\ln{Λ}\) is the Coulomb logarithm accounting for small angle collisions.

See Equation (2.86) in Callen [n.d.].

This form of the Lorentz collision frequency differs from the form found in SingleParticleCollisionFrequencies in that \(v\) is the mean thermal velocity between particle species in this method (as opposed to the drift velocity between species).

Maxwellian_avg_ei_collision_freq

Average momentum relaxation rate for a slowly flowing Maxwellian distribution of electrons, relative to a population of stationary ions.

This function assumes that both populations are Maxwellian, and \(T_i ≲ T_e\).

Callen [n.d.] provides a derivation of this as an average collision frequency between electrons and ions for a Maxwellian distribution. It is thus a special case of the collision frequency with an averaging factor, and is on many occasions in transport theory the most relevant collision frequency that has to be considered. It commonly occurs in relation to diffusion and resistivity in plasmas.

Raises:
  • PhysicsError – If the test particles are not ‘slowly flowing’ relative to the field particles (see notes).

  • ValueError – If the specified interaction is not electron-ion.

Examples

>>> import astropy.units as u
>>> v_drift = 1 * u.m / u.s
>>> n_a = 1e26 * u.m**-3
>>> T_a = 1 * u.eV
>>> n_b = 1e26 * u.m**-3
>>> T_b = 1e3 * u.eV
>>> Coulomb_log = 10 * u.dimensionless_unscaled
>>> electron_ion_collisions = MaxwellianCollisionFrequencies(
...     "e-",
...     "Na+",
...     v_drift=v_drift,
...     n_a=n_a,
...     T_a=T_a,
...     n_b=n_b,
...     T_b=T_b,
...     Coulomb_log=Coulomb_log,
... )
>>> electron_ion_collisions.Maxwellian_avg_ei_collision_freq
<Quantity 2.8053078...e+15 Hz>
Maxwellian_avg_ii_collision_freq

Average momentum relaxation rate for a slowly flowing Maxwellian distribution of ions, relative to a population of stationary ions.

This function assumes that both populations are Maxwellian, and \(T_i ≲ T_e\).

Callen [n.d.] provides a derivation of this as an average collision frequency between ions and ions for a Maxwellian distribution. It is thus a special case of the collision frequency with an averaging factor.

Raises:
  • PhysicsError – The test particles are not ‘slowly flowing’ relative to the field particles (see notes).

  • ValueError – If the specified interaction is not ion-ion.

Examples

>>> import astropy.units as u
>>> v_drift = 1 * u.m / u.s
>>> n_a = 1e26 * u.m**-3
>>> T_a = 1e3 * u.eV
>>> n_b = 1e26 * u.m**-3
>>> T_b = 1e3 * u.eV
>>> Coulomb_log = 10 * u.dimensionless_unscaled
>>> ion_ion_collisions = MaxwellianCollisionFrequencies(
...     "Na+",
...     "Na+",
...     v_drift=v_drift,
...     n_a=n_a,
...     T_a=T_a,
...     n_b=n_b,
...     T_b=T_b,
...     Coulomb_log=Coulomb_log,
... )
>>> ion_ion_collisions.Maxwellian_avg_ii_collision_freq
<Quantity 1.1223822...e+08 Hz>