thermal_speed¶
-
plasmapy.formulary.
thermal_speed
(T: Unit("K"), particle: plasmapy.particles.particle_class.Particle, method='most_probable', mass: Unit("kg") = <Quantity nan kg>, ndim=3) -> Unit("m / s")¶ Return the most probable speed for a particle within a Maxwellian distribution.
Aliases:
vth_
- Parameters
T (Quantity) – The particle temperature in either kelvin or energy per particle
particle (Particle) – Representation of the particle species (e.g.,
'p'
for protons,'D+'
for deuterium, or'He-4 +1'
for singly ionized helium-4). If no charge state information is provided, then the particles are assumed to be singly charged.method (str, optional) – Method to be used for calculating the thermal speed. Options are
'most_probable'
(default),'rms'
, and'mean_magnitude'
.mass (Quantity) – The particle’s mass override. Defaults to NaN and if so, doesn’t do anything, but if set, overrides mass acquired from
particle
. Useful with relative velocities of particles.ndim (int) – Dimensionality of space in which to calculate thermal velocity. Valid values are 1,2,3.
- Returns
V – particle thermal speed
- Return type
- Raises
TypeError – The particle temperature is not a ~astropy.units.Quantity
UnitConversionError – If the particle temperature is not in units of temperature or energy per particle
ValueError – The particle temperature is invalid or particle cannot be used to identify an isotope or particle
- Warns
RelativityWarning – If the ion sound speed exceeds 5% of the speed of light, or
~astropy.units.UnitsWarning – If units are not provided, SI units are assumed.
Notes
The particle thermal speed is given by:
\[V_{th,i} = \sqrt{\frac{N k_B T_i}{m_i}}\]where the value of N depends on the dimensionality and the definition of \(v_{th}\): most probable, root-mean-square (RMS), or mean magnitude. The value of N in each case is
Values of constant N¶ Dim.
Most-Probable
RMS
Mean-Magnitude
1D
0
1
\(2/\pi\)
2D
1
2
\(\pi/2\)
3D
2
3
\(8/\pi\)
The definition of thermal velocity varies by the square root of two depending on whether or not this velocity absorbs that factor in the expression for a Maxwellian distribution. In particular, the expression given in the NRL Plasma Formulary [1] is a square root of two smaller than the result from this function.
Examples
>>> from astropy import units as u >>> thermal_speed(5*u.eV, 'p') <Quantity 30949.6... m / s> >>> thermal_speed(1e6*u.K, particle='p') <Quantity 128486... m / s> >>> thermal_speed(5*u.eV, particle='e-') <Quantity 132620... m / s> >>> thermal_speed(1e6*u.K, particle='e-') <Quantity 550569... m / s> >>> thermal_speed(1e6*u.K, "e-", method="rms") <Quantity 674307... m / s> >>> thermal_speed(1e6*u.K, "e-", method="mean_magnitude") <Quantity 621251... m / s>