Source code for torx.normalization.unit_conversions_m
"""
Small functions converting from a physical params dict to quantities.
These are mostly designed to give the value of a particular quantity when
calculated from the reference parameters.
"""
import numpy as np
from torx import Quantity as Q_, unit_registry as ureg
elementary_charge = Q_(1, "elementary_charge")
electron_mass = Q_(1, "electron_mass")
proton_mass = Q_(1, "proton_mass")
atomic_mass_units = Q_(1, "amu")
speed_of_light = Q_(1, "speed_of_light")
epsilon0 = Q_(1, "vacuum_permittivity")
mu0 = Q_(1, "vacuum_permeability")
@property
def rho_s0(self):
"""
Reference Ion Larmor Radius [m].
This is also the scale length for perpendicular quantities.
"""
return (
np.sqrt(self.Ti0 * self.Mi) / (self.elementary_charge * self.B0)
).to_base_units()
@property
def c_s0(self):
"""Reference Sound speed [m/s]."""
return (np.sqrt(self.Te0 / self.Mi)).to_base_units()
@property
def v_A0(self):
"""Reference Alfven speed [m/s]."""
return (self.B0 / np.sqrt(self.mu0 * self.n0 * self.Mi)).to_base_units()
@property
def tau_0(self):
"""Time normalization [s]."""
return (self.R0 / self.c_s0).to("s")
@ureg.wraps(ret="", args=["cm**-3", "eV"], strict=True)
def Coulomb_logarithm_ee(n, Te):
"""
Coulomb logarithm (unitless) for thermal electron-electron collisions.
You can pass n and Te in any units, and these will be converted to
cgs+eV such that the following formula applies.
"""
if Te > 10:
return 24.0 - np.log(np.sqrt(n) / Te)
else:
return 23.0 - np.log(np.sqrt(n) / Te ** (1.5))
@property
def lnLambda0(self):
"""
Thermal electron-electron Coulomb logarithm for reference parameters.
Note that we usually the Coulomb logarithm to be constant to ensure that the
collision operator is bilinear.
"""
return Coulomb_logarithm_ee(self.n0, self.Te0)
def tau_ee(n, Te, Z):
"""
Electron-electron collision time.
This is the mean time required for the direction of motion of an individual electron
to be deflected by approximately 90 degrees due to collisions with electrons.
Using 3.182 from Fitzpatrick (in SI units), but taking the parametric dependency
from equation 2.5e from Braginskii (n = n_e = Z n_i, so n_i = n/Z.
Therefore Z**2 * n_i = Z n) -- since Braginskii is c.g.s. and conversion to s.i. is a pain.
"""
lnLambda = Coulomb_logarithm_ee(n, Te)
return (
(
6.0
* np.sqrt(2.0)
* np.pi**1.5
* epsilon0**2
* np.sqrt(electron_mass)
* Te**1.5
)
/ (lnLambda * elementary_charge**4 * Z * n)
).to("s")
def tau_ii(n, Te, Ti, Mi, Z):
"""
Ion-ion collision time.
This is the mean time required for the direction of motion of an individual ion
to be deflected by approximately 90 degrees due to collisions with ions.
Using 3.184 from Fitzpatrick (in SI units), but taking the parametric dependency from
equation 2.5i from Braginskii (n = n_e = Z n_i, so n_i = n/Z.
Therefore Z**4 * n_i = Z**3 n).
"""
lnLambda = Coulomb_logarithm_ee(n, Te)
return (
(12.0 * np.pi**1.5 * epsilon0**2 * np.sqrt(Mi) * Ti**1.5)
/ (lnLambda * elementary_charge**4 * Z**3 * n)
).to("s")
@property
def tau_e0(self):
"""Reference electron-electron collision time using Z_eff."""
return tau_ee(self.n0, self.Te0, self.Z_eff)
@property
def tau_i0(self):
"""Reference ion-ion collision time."""
return tau_ii(self.n0, self.Te0, self.Ti0, self.Mi, self.Z)
@property
def delta(self):
"""
Ratio of major radius to reference ion Larmor radius.
Conversion from perpendicular length scale to parallel length scale (>>1).
"""
return self.R0 / self.rho_s0
@property
def rhostar(self):
"""
Ratio of reference ion Larmor radius to major radius (inverse delta).
Conversion from parallel length scale to perpendicular length scale (<<1).
"""
return self.rho_s0 / self.R0
@property
def beta_0(self):
"""Reference electron dynamical beta."""
return self.c_s0**2 / (self.v_A0**2)
@property
def mu(self):
"""Electron to ion mass ratio."""
return (self.electron_mass / self.Mi).to("")
@property
def zeta(self):
"""Ti0/Te0: Ion to electron reference temperature ratio."""
return self.Ti0 / self.Te0
@property
def tau_e_norm(self):
"""Normalized reference electron collision time."""
return self.tau_e0 * self.c_s0 / self.R0
@property
def tau_i_norm(self):
"""Normalized reference ion collision time."""
return self.tau_i0 * self.c_s0 / self.R0
@property
def nu_e0(self):
"""Normalized reference electron collision frequency."""
return self.tau_e_norm**-1
@property
def nu_i0(self):
"""Normalized reference ion collision frequency."""
return self.tau_i_norm**-1
@property
def chipar_e(self):
"""Parallel electron heat conductivity coefficient."""
return 3.16 * self.tau_e_norm * self.mu**-1
@property
def chipar_i(self):
"""Parallel ion heat conductivity coefficient."""
return 3.90 * self.zeta * self.tau_i_norm
@property
def etapar_e(self):
"""Parallel electron resistivity coefficient."""
return 0.51 * self.nu_e0 * self.mu
@property
def etapar_i(self):
"""Parallel ion viscosity coefficient."""
return 0.96 * self.tau_i_norm