def collision_rate_ion_ion(T_i, n_i, ion_particle, coulomb_log=None, V=None): r""" momentum relaxation ion-ion collision rate From [3]_, equations (2.36) and (2.122) Considering a Maxwellian distribution of "test" ions colliding with a Maxwellian distribution of "field" ions. Note, it is assumed that electrons are present in such numbers as to establish quasineutrality, but the effects of the test ions colliding with them are not considered here. This result is an ion momentum relaxation rate, and is used in many classical transport expressions. It is equivalent to: * 1/tau_i from ref [1]_ eqn (1) pp. #, * 1/tau_i from ref [2]_ eqn (1) pp. #, * nu_i\i_S from ref [2]_ eqn (1) pp. #, Parameters ---------- T_i : ~astropy.units.Quantity The electron temperature of the Maxwellian test ions n_i : ~astropy.units.Quantity The number density of the Maxwellian test ions ion_particle: string String signifying a particle type of the test and field ions, including charge state information. This function assumes the test and field ions are the same species. coulomb_log : float or dimensionless ~astropy.units.Quantity, optional option to specify a Coulomb logarithm of the electrons on the ions. If not specified, the Coulomb log will is calculated using the .transport/Coulomb_logarithm function. References ---------- .. [1] Braginskii .. [2] Formulary .. [3] Callen Chapter 2, http://homepages.cae.wisc.edu/~callen/chap2.pdf """ from plasmapy.physics.transport import Coulomb_logarithm T_i = T_i.to(units.K, equivalencies=units.temperature_energy()) if coulomb_log is not None: coulomb_log_val = coulomb_log else: particles = [ion_particle, ion_particle] coulomb_log_val = Coulomb_logarithm(T_i, n_i, particles, V) Z_i = atomic.integer_charge(ion_particle) m_i = atomic.ion_mass(ion_particle) nu_i = 4 / 3 * np.sqrt(np.pi / m_i) / (4 * np.pi * eps0)**2 * e**4 * \ n_i * Z_i**4 * coulomb_log_val / (k_B * T_i)**1.5 return nu_i.to(1 / units.s)
def grab_charge(ion, z_mean=None): """Utility function to merge two possible inputs for particle charge. Parameters ---------- ion : str or `plasmapy.atomic.Particle` a string representing a charged particle, or a Particle object. z_mean : float An optional float describing the average ionization of a particle species. Returns ------- float if `z_mean` was passed, `z_mean`, otherwise, the integer charge of the `ion`. """ if z_mean is None: # warnings.warn("No z_mean given, defaulting to atomic charge", # PhysicsWarning) Z = atomic.integer_charge(ion) else: # using average ionization provided by user Z = z_mean return Z
def collision_rate_electron_ion(T_e, n_e, ion_particle, coulomb_log=None, V=None): r""" Momentum relaxation electron-ion collision rate From [3]_, equations (2.17) and (2.120) Considering a Maxwellian distribution of "test" electrons colliding with a Maxwellian distribution of "field" ions. This result is an electron momentum relaxation rate, and is used in many classical transport expressions. It is equivalent to: * 1/tau_e from ref [1]_ eqn (1) pp. #, * 1/tau_e from ref [2]_ eqn (1) pp. #, * nu_e\i_S from ref [2]_ eqn (1) pp. #, Parameters ---------- T_e : ~astropy.units.Quantity The electron temperature of the Maxwellian test electrons n_e : ~astropy.units.Quantity The number density of the Maxwellian test electrons ion_particle: str String signifying a particle type of the field ions, including charge state information. coulomb_log : float or dimensionless ~astropy.units.Quantity, optional Option to specify a Coulomb logarithm of the electrons on the ions. If not specified, the Coulomb log will is calculated using the `~plasmapy.physics.transport.Coulomb_logarithm` function. References ---------- .. [1] Braginskii .. [2] Formulary .. [3] Callen Chapter 2, http://homepages.cae.wisc.edu/~callen/chap2.pdf """ from plasmapy.physics.transport import Coulomb_logarithm T_e = T_e.to(u.K, equivalencies=u.temperature_energy()) if coulomb_log is not None: coulomb_log_val = coulomb_log else: particles = ['e', ion_particle] coulomb_log_val = Coulomb_logarithm(T_e, n_e, particles, V) Z_i = atomic.integer_charge(ion_particle) nu_e = 4 / 3 * np.sqrt(2 * np.pi / m_e) / (4 * np.pi * eps0) ** 2 * \ e ** 4 * n_e * Z_i * coulomb_log_val / (k_B * T_e) ** 1.5 return nu_e.to(1 / u.s)
def _boilerPlate(T, particles, V): """ Some boiler plate code for checking if inputs to functions in collisions.py are good. Also obtains reduced in mass in a 2 particle collision system along with thermal velocity. """ # checking temperature is in correct units T = T.to(u.K, equivalencies=u.temperature_energy()) # extracting particle information if not isinstance(particles, (list, tuple)) or len(particles) != 2: raise ValueError("Particles input must be a " "list or tuple containing representations of two " f"charged particles. Got {particles} instead.") masses = np.zeros(2) * u.kg charges = np.zeros(2) * u.C for particle, i in zip(particles, range(2)): try: masses[i] = ion_mass(particles[i]) except Exception: raise ValueError("Unable to find mass of particle: " f"{particles[i]}.") try: charges[i] = np.abs(e * integer_charge(particles[i])) if charges[i] is None: raise ValueError("Unable to find charge of particle: " f"{particles[i]}.") except Exception: raise ValueError("Unable to find charge of particle: " f"{particles[i]}.") # obtaining reduced mass of 2 particle collision system reduced_mass = masses[0] * masses[1] / (masses[0] + masses[1]) # getting thermal velocity of system if no velocity is given if np.isnan(V): V = np.sqrt(2 * k_B * T / reduced_mass).to(u.m / u.s) _check_relativistic(V, 'V') return T, masses, charges, reduced_mass, V
def collision_rate_ion_ion(T_i, n_i, ion_particle, coulomb_log=None, V=None, coulomb_log_method="classical"): r""" Momentum relaxation ion-ion collision rate From [3]_, equations (2.36) and (2.122) Considering a Maxwellian distribution of "test" ions colliding with a Maxwellian distribution of "field" ions. Note, it is assumed that electrons are present in such numbers as to establish quasineutrality, but the effects of the test ions colliding with them are not considered here. This result is an ion momentum relaxation rate, and is used in many classical transport expressions. It is equivalent to: * 1/tau_i from ref [1]_ eqn (1) pp. #, * 1/tau_i from ref [2]_ eqn (1) pp. #, * nu_i\i_S from ref [2]_ eqn (1) pp. #, Parameters ---------- T_i : ~astropy.units.Quantity The electron temperature of the Maxwellian test ions n_i : ~astropy.units.Quantity The number density of the Maxwellian test ions ion_particle: str String signifying a particle type of the test and field ions, including charge state information. This function assumes the test and field ions are the same species. V : ~astropy.units.Quantity, optional The relative velocity between particles. If not provided, thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T` where `mu` is the reduced mass. coulomb_log : float or dimensionless ~astropy.units.Quantity, optional Option to specify a Coulomb logarithm of the electrons on the ions. If not specified, the Coulomb log will is calculated using the ~plasmapy.physics.transport.Coulomb_logarithm function. coulomb_log_method : string, optional Method used for Coulomb logarithm calculation (see that function for more documentation). Choose from "classical" or "GMS-1" to "GMS-6". References ---------- .. [1] Braginskii .. [2] Formulary .. [3] Callen Chapter 2, http://homepages.cae.wisc.edu/~callen/chap2.pdf Examples -------- >>> from astropy import units as u >>> collision_rate_ion_ion(0.1*u.eV, 1e6/u.m**3, 'p') <Quantity 2.97315582e-05 1 / s> >>> collision_rate_ion_ion(100*u.eV, 1e6/u.m**3, 'p') <Quantity 1.43713193e-09 1 / s> >>> collision_rate_ion_ion(100*u.eV, 1e20/u.m**3, 'p') <Quantity 66411.80316364 1 / s> >>> collision_rate_ion_ion(100*u.eV, 1e20/u.m**3, 'p', coulomb_log_method='GMS-1') <Quantity 66407.00859126 1 / s> >>> collision_rate_ion_ion(100*u.eV, 1e20/u.m**3, 'p', V = c/100) <Quantity 6.53577473 1 / s> >>> collision_rate_ion_ion(100*u.eV, 1e20/u.m**3, 'p', coulomb_log=20) <Quantity 95918.76240877 1 / s> """ from plasmapy.physics.transport.collisions import Coulomb_logarithm T_i = T_i.to(u.K, equivalencies=u.temperature_energy()) m_i = atomic.ion_mass(ion_particle) if V is not None: V = V else: # ion thermal velocity (most probable) V = np.sqrt(2 * k_B * T_i / m_i) if coulomb_log is not None: coulomb_log_val = coulomb_log else: particles = [ion_particle, ion_particle] coulomb_log_val = Coulomb_logarithm(T_i, n_i, particles, V, method=coulomb_log_method) Z_i = atomic.integer_charge(ion_particle) # this is the same as b_perp in collisions.py, using most probable thermal velocity for V # and using ion mass instead of reduced mass bperp = (Z_i * e)**2 / (4 * np.pi * eps0 * m_i * V**2) # collisional cross-section sigma = np.pi * (2 * bperp)**2 # collisional frequency with Coulomb logarithm to correct for small angle collisions nu = n_i * sigma * V * coulomb_log_val # this coefficient is the constant that pops out when comparing this definition of # collisional frequency to the one in collisions.py coeff = np.sqrt(8 / np.pi) / 3 # collisional frequency modified by the constant difference nu_i = coeff * nu return nu_i.to(1 / u.s)
def lower_hybrid_frequency(B, n_i, ion='p+'): r""" Return the lower hybrid frequency. Parameters ---------- B : ~astropy.units.Quantity The magnetic field magnitude in units convertible to tesla. n_i : ~astropy.units.Quantity Ion number density. ion : str, optional Representation of the ion species (e.g., 'p' for protons, 'D+' for deuterium, or 'He-4 +1' for singly ionized helium-4), which defaults to protons. If no charge state information is provided, then the ions are assumed to be singly charged. Returns ------- omega_lh : ~astropy.units.Quantity The lower hybrid frequency in radians per second. Raises ------ TypeError If either of `B` or `n_i` is not a `~astropy.units.Quantity`, or ion is of an inappropriate type. ~astropy.units.UnitConversionError If either of `B` or `n_i` is in incorrect units. ValueError If either of `B` or `n_i` contains invalid values or are of incompatible dimensions, or ion cannot be used to identify an ion or isotope. Warns ----- ~astropy.units.UnitsWarning If units are not provided, SI units are assumed Notes ----- The lower hybrid frequency is given through the relation .. math:: \frac{1}{\omega_{lh}^2} = \frac{1}{\omega_{ci}^2 + \omega_{pi}^2} + \frac{1}{\omega_{ci}\omega_{ce}} where :math:`\omega_{ci}` is the ion gyrofrequency, :math:`\omega_{ce}` is the electron gyrofrequency, and :math:`\omega_{pi}` is the ion plasma frequency. Example ------- >>> from astropy import units as u >>> lower_hybrid_frequency(0.2*u.T, n_i=5e19*u.m**-3, ion='D+') <Quantity 5.78372733e+08 rad / s> """ # We do not need a charge state here, so the sole intent is to # catch invalid ions. try: atomic.integer_charge(ion) except Exception: raise ValueError("Invalid ion in lower_hybrid_frequency.") omega_ci = gyrofrequency(B, particle=ion) omega_pi = plasma_frequency(n_i, particle=ion) omega_ce = gyrofrequency(B) omega_lh = ((omega_ci * omega_ce)**-1 + omega_pi**-2)**-0.5 # TODO possibly optimize the above line via np.sqrt omega_lh = omega_lh return omega_lh.to(u.rad / u.s)
def plasma_frequency(n, particle='e-', z_mean=None): r"""Calculate the particle plasma frequency. Parameters ---------- n : ~astropy.units.Quantity Particle number density in units convertible to per cubic meter particle : str, optional Representation of the particle species (e.g., 'p' for protons, 'D+' for deuterium, or 'He-4 +1' for singly ionized helium-4), which defaults to electrons. If no charge state information is provided, then the particles are assumed to be singly charged. z_mean : ~astropy.units.Quantity, optional The average ionization (arithmetic mean) for a plasma where the a macroscopic description is valid. If this quantity is not given then the atomic charge state (`int`) of the ion is used. This is effectively an average plasma frequency for the plasma where multiple charge states are present. Returns ------- omega_p : ~astropy.units.Quantity The particle plasma frequency in radians per second. Raises ------ TypeError If n_i is not a `~astropy.units.Quantity` or particle is not of an appropriate type. UnitConversionError If `n_i` is not in correct units ValueError If `n_i` contains invalid values or particle cannot be used to identify an particle or isotope. Warns ----- ~astropy.units.UnitsWarning If units are not provided, SI units are assumed Notes ----- The particle plasma frequency is .. math:: \omega_{pi} = Z e \sqrt{\frac{n_i}{\epsilon_0 m_i}} At present, astropy.units does not allow direct conversions from radians/second for angular frequency to 1/second or Hz for frequency. The dimensionless_angles equivalency allows that conversion, but does not account for the factor of 2*pi. The alternatives are to convert to cycle/second or to do the conversion manually, as shown in the examples. Example ------- >>> from astropy import units as u >>> plasma_frequency(1e19*u.m**-3, particle='p') <Quantity 4.16329453e+09 rad / s> >>> plasma_frequency(1e19*u.m**-3, particle='D+') <Quantity 2.94462452e+09 rad / s> >>> plasma_frequency(1e19*u.m**-3) <Quantity 1.78398636e+11 rad / s> """ try: m = atomic.ion_mass(particle) if z_mean is None: # warnings.warn("No z_mean given, defaulting to atomic charge", # PhysicsWarning) try: Z = atomic.integer_charge(particle) except Exception: Z = 1 else: # using user provided average ionization Z = z_mean Z = np.abs(Z) # TODO REPLACE WITH Z = np.abs(grab_charge(particle, z_mean)), some bugs atm except Exception: raise ValueError(f"Invalid particle, {particle}, in " "plasma_frequency.") omega_p = u.rad * Z * e * np.sqrt(n / (eps0 * m)) return omega_p.si
def coupling_parameter(T, n_e, particles, z_mean=np.nan * u.dimensionless_unscaled, V=np.nan * u.m / u.s, method="classical"): r"""Coupling parameter. Coupling parameter compares Coulomb energy to kinetic energy (typically) thermal. Classical plasmas are weakly coupled Gamma << 1, whereas dense plasmas tend to have significant to strong coupling Gamma >= 1. Parameters ---------- T : ~astropy.units.Quantity Temperature in units of temperature or energy per particle, which is assumed to be equal for both the test particle and the target particle n_e : ~astropy.units.Quantity The electron density in units convertible to per cubic meter. particles : tuple A tuple containing string representations of the test particle (listed first) and the target particle (listed second) z_mean : ~astropy.units.Quantity, optional The average ionization (arithmetic mean) for a plasma where the a macroscopic description is valid. This is used to recover the average ion density (given the average ionization and electron density) for calculating the ion sphere radius for non-classical impact parameters. V : ~astropy.units.Quantity, optional The relative velocity between particles. If not provided, thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T` where `mu` is the reduced mass. method: str, optional Selects which theory to use when calculating the Coulomb logarithm. Defaults to classical method. Returns ------- coupling : float or numpy.ndarray The coupling parameter for a plasma. Raises ------ ValueError If the mass or charge of either particle cannot be found, or any of the inputs contain incorrect values. UnitConversionError If the units on any of the inputs are incorrect TypeError If the n_e, T, or V are not Quantities. Warns ----- ~astropy.units.UnitsWarning If units are not provided, SI units are assumed ~plasmapy.utils.RelativityWarning If the input velocity is greater than 5% of the speed of light. Notes ----- The coupling parameter is given by .. math:: \Gamma = \frac{E_{Coulomb}}{E_{Kinetic}} The Coulomb energy is given by .. math:: E_{Coulomb} = \frac{Z_1 Z_2 q_e^2}{4 \pi \epsilon_0 r} where :math:`r` is the Wigner-Seitz radius, and 1 and 2 refer to particle species 1 and 2 between which we want to determine the coupling. In the classical case the kinetic energy is simply the thermal energy .. math:: E_{kinetic} = k_B T_e The quantum case is more complex. The kinetic energy is dominated by the Fermi energy, modulated by a correction factor based on the ideal chemical potential. This is obtained more precisely by taking the the thermal kinetic energy and dividing by the degeneracy parameter, modulated by the Fermi integral [1]_ .. math:: E_{kinetic} = 2 k_B T_e / \chi f_{3/2} (\mu_{ideal} / k_B T_e) where :math:`\chi` is the degeneracy parameter, :math:`f_{3/2}` is the Fermi integral, and :math:`\mu_{ideal}` is the ideal chemical potential. The degeneracy parameter is given by .. math:: \chi = n_e \Lambda_{deBroglie} ^ 3 where :math:`n_e` is the electron density and :math:`\Lambda_{deBroglie}` is the thermal deBroglie wavelength. See equations 1.2, 1.3 and footnote 5 in [2]_ for details on the ideal chemical potential. Examples -------- >>> from astropy import units as u >>> n = 1e19*u.m**-3 >>> T = 1e6*u.K >>> particles = ('e', 'p') >>> coupling_parameter(T, n, particles) <Quantity 5.80330315e-05> >>> coupling_parameter(T, n, particles, V=1e6*u.m/u.s) <Quantity 5.80330315e-05> References ---------- .. [1] Dense plasma temperature equilibration in the binary collision approximation. D. O. Gericke et. al. PRE, 65, 036418 (2002). DOI: 10.1103/PhysRevE.65.036418 .. [2] Bonitz, Michael. Quantum kinetic theory. Stuttgart: Teubner, 1998. """ # boiler plate checks T, masses, charges, reduced_mass, V = _boilerPlate(T=T, particles=particles, V=V) if np.isnan(z_mean): # using mean charge to get average ion density. # If you are running this, you should strongly consider giving # a value of z_mean as an argument instead. Z1 = np.abs(atomic.integer_charge(particles[0])) Z2 = np.abs(atomic.integer_charge(particles[1])) Z = (Z1 + Z2) / 2 # getting ion density from electron density n_i = n_e / Z # getting Wigner-Seitz radius based on ion density radius = Wigner_Seitz_radius(n_i) else: # getting ion density from electron density n_i = n_e / z_mean # getting Wigner-Seitz radius based on ion density radius = Wigner_Seitz_radius(n_i) # Coulomb potential energy between particles if np.isnan(z_mean): coulombEnergy = charges[0] * charges[1] / (4 * np.pi * eps0 * radius) else: coulombEnergy = (z_mean * e) ** 2 / (4 * np.pi * eps0 * radius) if method == "classical": # classical thermal kinetic energy kineticEnergy = k_B * T elif method == "quantum": # quantum kinetic energy for dense plasmas lambda_deBroglie = thermal_deBroglie_wavelength(T) chemicalPotential = chemical_potential(n_e, T) fermiIntegral = Fermi_integral(chemicalPotential.si.value, 1.5) denom = (n_e * lambda_deBroglie ** 3) * fermiIntegral kineticEnergy = 2 * k_B * T / denom if np.imag(kineticEnergy) == 0: kineticEnergy = np.real(kineticEnergy) else: raise ValueError("Kinetic energy should not be imaginary." "Something went horribly wrong.") coupling = coulombEnergy / kineticEnergy return coupling.to(u.dimensionless_unscaled)
def Alfven_speed(B, density, ion="p+", z_mean=None): r""" Returns the Alfven speed. Parameters ---------- B : ~astropy.units.Quantity The magnetic field magnitude in units convertible to tesla. density : ~astropy.units.Quantity Either the ion number density in units convertible to 1 / m**3, or the mass density in units convertible to kg / m**3. ion : str, optional Representation of the ion species (e.g., `'p'` for protons, `'D+'` for deuterium, or `'He-4 +1'` for singly ionized helium-4), which defaults to protons. If no charge state information is provided, then the ions are assumed to be singly charged. z_mean : ~astropy.units.Quantity, optional The average ionization (arithmetic mean) for a plasma where the a macroscopic description is valid. If this quantity is not given then the atomic charge state (integer) of the ion is used. This is effectively an average Alfven speed for the plasma where multiple charge states are present. Returns ------- V_A : ~astropy.units.Quantity with units of velocity The Alfven velocity of the plasma in units of meters per second. Raises ------ TypeError The magnetic field and density arguments are not Quantities and cannot be converted into Quantities. ~astropy.units.UnitConversionError If the magnetic field or density is not in appropriate units. ~plasmapy.utils.RelativityError If the Alfven velocity is greater than or equal to the speed of light ValueError If the density is negative, or the ion mass or charge state cannot be found. UserWarning if units are not provided and SI units are assumed. Warns ----- ~plasmapy.utils.RelativityWarning If the Alfven velocity exceeds 10% of the speed of light Notes ----- The Alfven velocity :math:`V_A` is the typical propagation speed of magnetic disturbances in a plasma, and is given by: .. math:: V_A = \frac{B}{\sqrt{\mu_0\rho}} where the mass density is :math:`\rho = n_i m_i + n_e m_e`. This expression does not account for relativistic effects, and loses validity when the resulting speed is a significant fraction of the speed of light. This function switches B and density when B has units of number density or mass density and density has units of magnetic field strength. Examples -------- >>> from astropy import units as u >>> from plasmapy.constants import m_p, m_e >>> B = 0.014*u.T >>> n = 5e19*u.m**-3 >>> rho = n*(m_p+m_e) >>> ion = 'p' >>> Alfven_speed(B, n, ion) <Quantity 43173.87029559 m / s> >>> Alfven_speed(B, rho, ion) <Quantity 43173.87029559 m / s> >>> Alfven_speed(B, rho, ion).to(u.cm/u.us) <Quantity 4.31738703 cm / us> """ utils._check_quantity(B, 'B', 'Alfven_speed', u.T) utils._check_quantity(density, 'density', 'Alfven_speed', [u.m**-3, u.kg / u.m**3], can_be_negative=False) B = B.to(u.T) density = density.si if density.unit == u.m**-3: try: m_i = atomic.ion_mass(ion) if z_mean is None: # warnings.warn("No z_mean given, defaulting to atomic charge", # PhysicsWarning) try: Z = atomic.integer_charge(ion) except AtomicError: Z = 1 else: # using average ionization provided by user Z = z_mean # ensuring positive value of Z Z = np.abs(Z) except AtomicError: raise ValueError("Invalid ion in Alfven_speed.") rho = density * m_i + Z * density * m_e elif density.unit == u.kg / u.m**3: rho = density try: V_A = (np.abs(B) / np.sqrt(mu0 * rho)).to(u.m / u.s) except Exception: raise ValueError("Unable to find Alfven speed") return V_A
def gyrofrequency(B, particle='e-', signed=False, z_mean=None): r"""Calculate the particle gyrofrequency in units of radians per second. Parameters ---------- B : ~astropy.units.Quantity The magnetic field magnitude in units convertible to tesla. particle : str, optional Representation of the particle species (e.g., 'p' for protons, 'D+' for deuterium, or 'He-4 +1' for singly ionized helium-4), which defaults to electrons. If no charge state information is provided, then the particles are assumed to be singly charged. z_mean : float or ~astropy.units.Quantity, optional The average ionization (arithmetic mean) for a plasma where the a macroscopic description is valid. If this quantity is not given then the atomic charge state (integer) of the ion is used. This is effectively an average gyrofrequency for the plasma where multiple charge states are present, and should not be interpreted as the gyrofrequency for any single particle. signed : bool, optional The gyrofrequency can be defined as signed (negative for electron, positive for ion). Default is `False` (unsigned, i.e. always positive). Returns ------- omega_c : ~astropy.units.Quantity The particle gyrofrequency in units of radians per second Raises ------ TypeError If the magnetic field is not a Quantity or particle is not of an appropriate type ValueError If the magnetic field contains invalid values or particle cannot be used to identify an particle or isotope UserWarning If units are not provided and SI units are assumed Notes ----- The particle gyrofrequency is the angular frequency of particle gyration around magnetic field lines and is given by: .. math:: \omega_{ci} = \frac{Z e B}{m_i} The particle gyrofrequency is also known as the particle cyclotron frequency or the particle Larmor frequency. The recommended way to convert from angular frequency to frequency is to use an equivalency between cycles per second and Hertz, as Astropy's dimensionles_angles() equivalency does not account for the factor of 2*pi needed during this conversion. The dimensionless_angles() equivalency is appropriate when dividing a velocity by an angular frequency to get a length scale. Examples -------- >>> from astropy import units as u >>> gyrofrequency(0.1*u.T) <Quantity 1.75882002e+10 rad / s> >>> gyrofrequency(0.1*u.T, signed=True) <Quantity -1.75882002e+10 rad / s> >>> gyrofrequency(0.01*u.T, 'p') <Quantity 957883.32241481 rad / s> >>> gyrofrequency(0.01*u.T, 'p') <Quantity 957883.32241481 rad / s> >>> gyrofrequency(0.01*u.T, particle='T+') <Quantity 319964.54975911 rad / s> >>> omega_ce = gyrofrequency(0.1*u.T) >>> print(omega_ce) 17588200236.02124 rad / s >>> f_ce = omega_ce.to(u.Hz, equivalencies=[(u.cy/u.s, u.Hz)]) >>> print(f_ce) 2799249007.6528206 Hz """ try: m_i = atomic.ion_mass(particle) if z_mean is None: # warnings.warn("No z_mean given, defaulting to atomic charge", # PhysicsWarning) try: Z = atomic.integer_charge(particle) except AtomicError: Z = 1 else: # using user provided average ionization Z = z_mean except Exception: raise ValueError( "Invalid particle {} in gyrofrequency".format(particle)) if not signed: Z = abs(Z) omega_ci = u.rad * (Z * e * np.abs(B) / m_i).to(1 / u.s) return omega_ci
def ion_sound_speed(*ignore, T_e=0 * u.K, T_i=0 * u.K, gamma_e=1, gamma_i=3, ion='p+', z_mean=None): r""" Returns the ion sound speed for an electron-ion plasma. Parameters ---------- T_e : ~astropy.units.Quantity, optional, keyword-only Electron temperature in units of temperature or energy per particle. If this is not given, then the electron temperature is assumed to be zero. T_i : ~astropy.units.Quantity, optional, keyword-only Ion temperature in units of temperature or energy per particle. If this is not given, then the ion temperature is assumed to be zero. gamma_e : float or int, keyword-only The adiabatic index for electrons, which defaults to 1. This value assumes that the electrons are able to equalize their temperature rapidly enough that the electrons are effectively isothermal. gamma_i : float or int, keyword-only The adiabatic index for ions, which defaults to 3. This value assumes that ion motion has only one degree of freedom, namely along magnetic field lines. ion : str, optional, keyword-only Representation of the ion species (e.g., `'p'` for protons, `'D+'` for deuterium, or 'He-4 +1' for singly ionized helium-4), which defaults to protons. If no charge state information is provided, then the ions are assumed to be singly charged. z_mean : ~astropy.units.Quantity, optional, keyword-only The average ionization (arithmetic mean) for a plasma where the a macroscopic description is valid. If this quantity is not given then the atomic charge state (integer) of the ion is used. This is effectively an average ion sound speed for the plasma where multiple charge states are present. Returns ------- V_S : ~astropy.units.Quantity The ion sound speed in units of meters per second. Raises ------ TypeError If any of the arguments are not entered as keyword arguments or are of an incorrect type. ValueError If the ion mass, adiabatic index, or temperature are invalid. ~plasmapy.utils.PhysicsError If an adiabatic index is less than one. ~astropy.units.UnitConversionError If the temperature is in incorrect units. UserWarning If the ion sound speed exceeds 10% of the speed of light, or if units are not provided and SI units are assumed. Notes ----- The ion sound speed :math:`V_S` is approximately given by .. math:: V_S = \sqrt{\frac{\gamma_e Z k_B T_e + \gamma_i k_B T_i}{m_i}} where :math:`\gamma_e` and :math:`\gamma_i` are the electron and ion adiabatic indices, :math:`k_B` is the Boltzmann constant, :math:`T_e` and :math:`T_i` are the electron and ion temperatures, :math:`Z` is the charge state of the ion, and :math:`m_i` is the ion mass. This function assumes that the product of the wavenumber and the Debye length is small. In this limit, the ion sound speed is not dispersive (e.g., frequency independent). When the electron temperature is much greater than the ion temperature, the ion sound velocity reduces to :math:`\sqrt{\gamma_e k_B T_e / m_i}`. Ion acoustic waves can therefore occur even when the ion temperature is zero. Example ------- >>> from astropy import units as u >>> ion_sound_speed(T_e=5e6*u.K, T_i=0*u.K, ion='p', gamma_e=1, gamma_i=3) <Quantity 203155.0764042 m / s> >>> ion_sound_speed(T_e=5e6*u.K) <Quantity 203155.0764042 m / s> >>> ion_sound_speed(T_e=500*u.eV, T_i=200*u.eV, ion='D+') <Quantity 229586.01860212 m / s> """ if ignore: raise TypeError("All arguments are required to be keyword arguments " "in ion_sound_speed to prevent mixing up the electron " "and ion temperatures. An example call that uses the " "units subpackage from astropy is: " "ion_sound_speed(T_e=5*u.K, T_i=0*u.K, " "ion='D+')") try: m_i = atomic.ion_mass(ion) if z_mean is None: # warnings.warn("No z_mean given, defaulting to atomic charge", # PhysicsWarning) try: Z = atomic.integer_charge(ion) except AtomicError: Z = 1 else: # using average ionization provided by user Z = z_mean except AtomicError: raise ValueError("Invalid ion in ion_sound_speed.") if not isinstance(gamma_e, (float, int)): raise TypeError("The adiabatic index for electrons (gamma_e) must be " "a float or int in ion_sound_speed") if not isinstance(gamma_i, (float, int)): raise TypeError("The adiabatic index for ions (gamma_i) must be " "a float or int in ion_sound_speed") if not 1 <= gamma_e <= np.inf: raise utils.PhysicsError( "The adiabatic index for electrons must be between " "one and infinity") if not 1 <= gamma_i <= np.inf: raise utils.PhysicsError( "The adiabatic index for ions must be between " "one and infinity") T_i = T_i.to(u.K, equivalencies=u.temperature_energy()) T_e = T_e.to(u.K, equivalencies=u.temperature_energy()) try: V_S_squared = (gamma_e * Z * k_B * T_e + gamma_i * k_B * T_i) / m_i V_S = np.sqrt(V_S_squared).to(u.m / u.s) except Exception: raise ValueError("Unable to find ion sound speed.") return V_S
def inertial_length(n, particle='e-'): r"""Calculate the particle inertial length, Parameters ---------- n_i : ~astropy.units.Quantity Particle number density in units convertible to m**-3. particle : str, optional Representation of the particle species (e.g., 'p' for protons, 'D+' for deuterium, or 'He-4 +1' for singly ionized helium-4), which defaults to electrons. If no charge state information is provided, then the particles are assumed to be singly charged. Returns ------- d_i : ~astropy.units.Quantity Particles inertial length in meters. Raises ------ TypeError If n_i not a Quantity or particle is not a string. ~astropy.units.UnitConversionError If n_i is not in units of a number density. ValueError The particle density does not have an appropriate value. UserWarning If units are not provided and SI units are assumed. Notes ----- The particle inertial length is also known as an particle skin depth and is given by: .. math:: d_i = \frac{c}{\omega_{pi}} Example ------- >>> from astropy import units as u >>> inertial_length(5*u.m**-3, particle='He+') <Quantity 2.02985802e+08 m> >>> inertial_length(5*u.m**-3) <Quantity 2376534.75601976 m> """ try: Z = atomic.integer_charge(particle) except AtomicError: raise ValueError(f"Invalid particle {particle} in inertial_length.") if Z: Z = abs(Z) omega_p = plasma_frequency(n, particle=particle) d = (c / omega_p).to(u.m, equivalencies=u.dimensionless_angles()) return d