Hall_parameter
- plasmapy.formulary.dimensionless.Hall_parameter(
- n: Quantity,
- T: Quantity,
- B: Quantity,
- ion: str | int | integer | Particle | CustomParticle | Quantity,
- particle: str | int | integer | Particle | CustomParticle | Quantity,
- coulomb_log: float | None = None,
- V: Quantity | None = None,
- coulomb_log_method: str = 'classical',
Calculate the
particle
Hall parameter for a plasma.The Hall parameter for plasma species \(s\) (
particle
) is given by:\[β_{s} = \frac{Ω_{c s}}{ν_{s s^{\prime}}}\]where \(Ω_{c s}\) is the gyrofrequncy for plasma species \(s\) (
particle
) and \(ν_{s s^{\prime}}\) is the collision frequency between plasma species \(s\) (particle
) and species \(s^{\prime}\) (ion
).Aliases:
betaH_
- Parameters:
n (
Quantity
) – The number density associated withparticle
.T (
Quantity
) – The temperature of associated withparticle
.B (
Quantity
) – The magnetic field.ion (
Particle
) – The type of ionparticle
is colliding with.particle (
Particle
) – The particle species for which the Hall parameter is calculated for. 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.coulomb_log (
float
, optional) – Preset value for the Coulomb logarithm. Used mostly for testing purposes.V (
Quantity
) – The relative velocity betweenparticle
andion
. If not provided, then theparticle
thermal velocity is assumed (thermal_speed
).coulomb_log_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 ofCoulomb_logarithm
for more information about these methods.
- Returns:
Hall parameter for
particle
.- Return type:
See also
gyrofrequency
,fundamental_electron_collision_freq
,fundamental_ion_collision_freq
,Coulomb_logarithm
Notes
For calculating the collision frequency
fundamental_electron_collision_freq
is used whenparticle
is an electron andfundamental_ion_collision_freq
whenparticle
is an ion.The collision frequencies are calculated assuming a slowly moving Maxwellian distribution.
Examples
>>> import astropy.units as u >>> Hall_parameter(1e10 * u.m**-3, 2.8e2 * u.eV, 2.3 * u.T, "He-4 +1", "e-") <Quantity 2.500...e+15>