API Reference
Dimensionless numbers
PlasmaFormulary.plasma_beta — Function
plasma_beta(T, n, B)Compute the plamsa beta (β), the ratio of thermal pressure to magnetic pressure.
Arguments
T: The temperature of the plasma.n: The particle density of the plasma.B: The magnetic field in the plasma.
Example
julia> plasma_beta(1e6u"K", 1e19u"m^-3", 0.2u"T")
0.008674873511172188References
Frequencies
PlasmaFormulary.plasma_frequency — Function
plasma_frequency(n::NumberDensity, [q::Charge, mass::Mass]; to_hz = false)
plasma_frequency(n::NumberDensity, p::ParticleLike; kw...)Calculate the plasma frequency of a species.
\[ω_p = \sqrt{\frac{n q^2}{m ε₀}}\]
The plasma frequency is a characteristic frequency of the plasma. More often, it refers to the frequency at which electrons oscillate in the plasma.
Set to_hz = true to convert output from angular frequency to Hz.
Examples
julia> plasma_frequency(1e19u"m^-3") # plasma frequency
1.7839863654934653e11 rad s⁻¹
julia> plasma_frequency(1e19u"m^-3", :p) # proton plasma frequency
4.1632945624883513e9 rad s⁻¹References
PlasmaFormulary.gyrofrequency — Function
gyrofrequency(B::BField, p::ParticleLike; kw...)
gyrofrequency(B::BField, mass::Mass, q::Charge; to_hz = false)Calculate the gyrofrequency (or cyclotron frequency) of a charged particle's circular motion in a magnetic field. The gyrofrequency is the angular frequency of a charged particle's gyromotion around magnetic field lines.
Set to_hz = true to convert output from angular frequency to Hz.
References:
Examples
julia> gyrofrequency(0.01u"T", :p) # proton gyrofrequency
957883.3292211705 rad s⁻¹
julia> gyrofrequency(0.1u"T", :e; to_hz = true) # electron gyrofrequency as frequency
2.799248987233304e9 Hz
julia> gyrofrequency(250u"Gauss", "Fe"; z=13) # Fe2+ ion gyrofrequency
560682.3520611608 rad s⁻¹PlasmaFormulary.lower_hybrid_frequency — Function
lower_hybrid_frequency(B::BField, n_i::NumberDensity, ion::ParticleLike; kw...)Calculate the lower hybrid frequency.
\[\frac{1}{ω_{lh}^2} = \frac{1}{ω_{ci}^2 + ω_{pi}^2} + \frac{1}{ω_{ci} ω_{ce}}\]
where ωci is the ion gyrofrequency, ωce is the electron gyrofrequency, and ω_pi is the ion plasma frequency.
Set to_hz = true to convert output from angular frequency to Hz.
Examples
julia> lower_hybrid_frequency(0.2u"T", 5e19u"m^-3", "D+")
5.783727350182709e8 rad s⁻¹References
PlasmaFormulary.upper_hybrid_frequency — Function
upper_hybrid_frequency(B::BField, n_e::NumberDensity; to_hz = false)Calculate the upper hybrid frequency.
\[ω_{uh}^2 = ω_{ce}^2 + ω_{pe}^2\]
where ωce is the electron gyrofrequency and ωpe is the electron plasma frequency.
Set to_hz = true to convert output from angular frequency to Hz.
Examples
julia> upper_hybrid_frequency(0.2u"T", 5e19u"m^-3")
4.004594195988481e11 rad s⁻¹
julia> upper_hybrid_frequency(0.2u"T", 5e19u"m^-3"; to_hz = true)
6.37350961368681e10 HzReferences
Lengths
PlasmaFormulary.gyroradius — Function
Calculate the radius of circular motion for a charged particle in a uniform magnetic field
References: PlasmaPy API Documentation
Examples
julia> gyroradius(0.2u"T", Unitful.me, Unitful.q, 1e6u"K")
0.00015651672339994665 mPlasmaFormulary.inertial_length — Function
inertial_length(n::NumberDensity, q::Charge, mass::Mass)
inertial_length(n::NumberDensity, p::ParticleLike; kw...)Calculate the inertial length, the characteristic length scale for a particle to be accelerated in a plasma:
\[dᵢ = \frac{c}{ω_p}\]
where $ω_p$ is the plasma frequency.
The Hall effect becomes important on length scales shorter than the ion inertial length.
See also: plasma_frequency.
References: PlasmaPy API Documentation
PlasmaFormulary.Debye_length — Function
Debye_length(n::NumberDensity, T::EnergyOrTemp)Calculate the Debye length, exponential scale length for charge screening in an electron plasma with stationary ions:
\[λ_D = \sqrt{\frac{ε₀ k_B T}{n q^2}}\]
where $T$ is the electron temperature, $n$ is the electron density, and $q$ is the elementary charge.
Velocites
PlasmaFormulary.Alfven_velocity — Function
Alfven_velocity(B::BField, ρ)
Alfven_velocity(𝐁::Vector{BField}, ρ)Calculate the Alfven velocity for magnetic field vector. See also Alfven_speed.
PlasmaFormulary.diamagnetic_drift — Function
diamagnetic_drift(∇p, B::BFields, n, q)Calculate the diamagnetic drift velocity given by:
\[𝐯 = - (∇p × 𝐁) / (q n |𝐁|^2)\]
where $∇p$ is the pressure gradient.
PlasmaFormulary.ExB_drift — Function
ExB_drift(𝐄::EFields, 𝐁::BFields)Calculate the $E × 𝐁$ drift velocity given by:
\[𝐯 = (𝐄 × 𝐁) / |𝐁|^2\]
PlasmaFormulary.force_drift — Function
force_drift(𝐅, 𝐁::BFields, q)Calculate the general force drift for a particle in a magnetic field given by:
\[𝐯 = (𝐅 × 𝐁) / (q |𝐁|^2)\]
Speeds
PlasmaFormulary.Alfven_speed — Function
Alfven_speed(B::BFieldOrBFields, ρ)
Alfven_speed(B::BFieldOrBFields, n::NumberDensity, mass_number = 1)Alfvén speed $V_A$, the typical propagation speed of magnetic disturbances in a quasineutral plasma.
Note that this is different from the Alfven velocity, see also Alfven_velocity. References: PlasmaPy API Documentation
PlasmaFormulary.ion_sound_speed — Function
ion_sound_speed(T_e, T_i, m_i, Z; γ_e=1, γ_i=3, n_e=nothing, k=nothing)Calculate the ion sound speed for an electron-ion plasma given by:
\[V_S = \sqrt{\frac{γ_e Z k_B T_e + γ_i k_B T_i}{m_i (1 + k^2 λ_{D}^2)}}\]
If both n_e and k are given, includes dispersive correction.
Arguments
T_e: Electron temperature (Kor energy per particle)T_i: Ion temperature (Kor energy per particle)m_i: Ion massZ: Ion charge state (default: 1)γ_e: Electron adiabatic index (default: 1)γ_i: Ion adiabatic index (default: 3)n_e: Electron number density (optional)k: Wavenumber (optional)
PlasmaFormulary.thermal_speed — Function
thermal_speed(T::EnergyOrTemp, mass::Mass, method::ThermalVelocityMethod = MostProbable(), ndim = 3)Calculate the speed of thermal motion for particles with a Maxwellian distribution.
\[v_{th} = C_o \sqrt{\frac{k_B T}{m}}\]
where $T$ is the temperature associated with the distribution, $m$ is the particle's mass, and $C_o$ is a constant of proportionality. ```
PlasmaFormulary.thermal_temperature — Function
thermal_temperature(V::Velocity, mass::Mass, method::ThermalVelocityMethod = MostProbable(), ndim = 3)Current density
PlasmaFormulary.Alfven_current_density — Function
Alfven_current_density(Va::Velocity, n::NumberDensity)
Alfven_current_density(B::BFieldOrBFields, n::NumberDensity)
Alfven_current_density(B::BFieldOrBFields, dᵢ::Length)Calculate the Alfvén current density $J_A$, a natural scaling for current density:
\[J_A = e n V_A = \frac{B}{μ₀ dᵢ}\]
where $V_A$ is the Alfven speed, $dᵢ$ is ion inertial length.
See also: Alfven_speed, inertial_length.
Aliases
PlasmaFormulary.Aliases — Module
This module provides aliases of the most common plasma functionality for user convenience.
Aliases are denoted with a trailing underscore (e.g., alias_).
Constants
PlasmaFormulary.thermal_velocity_coefficients — Function
thermal_velocity_coefficients(method::ThermalVelocityMethod, ndim::Int)Get the thermal speed coefficient corresponding to the desired thermal speed definition.
Arguments
method::ThermalVelocityMethod: Method to be used for calculating the thermal speed. Valid values areMostProbable(),RMS(),MeanMagnitude(), andNRL().ndim::Val{Int}: Dimensionality (1D, 2D, 3D) of space in which to calculate thermal speed. Valid values areVal(1),Val(2), orVal{3}.
Otherwise undocumented functions
This section will be removed once the documentation is complete.
PlasmaFormulary.magnetic_pressure — Method
magnetic_pressure(B)Calculate the magnetic pressure.
PlasmaFormulary.thermal_pressure — Method
thermal_pressure(T, n)Calculate the thermal pressure for a Maxwellian distribution.
Arguments
T: The particle temperature or energy.n: The particle number density.