def test_raises(self, args, kwargs, _error): """ Test scenarios that cause plasma_frequency to raise an Exception. """ with pytest.raises(_error): plasma_frequency(*args, **kwargs)
def test_to_hz(self, args, kwargs): """Test behavior of the ``to_hz`` keyword.""" wp = plasma_frequency(*args, **kwargs) fp = plasma_frequency(*args, to_hz=True, **kwargs) assert isinstance(fp, u.Quantity) assert fp.unit == u.Hz assert fp.value == wp.value / (2.0 * np.pi)
def test_normal_vs_lite_values(self, inputs): """ Test that plasma_frequency and plasma_frequency_lite calculate the same values. """ particle = Particle(inputs["particle"]) inputs_unitless = { "n": inputs["n"].to(u.m**-3).value, "mass": particle.mass.value, } if "z_mean" in inputs: inputs_unitless["z_mean"] = inputs["z_mean"] else: try: inputs_unitless["z_mean"] = np.abs(particle.charge_number) except Exception: inputs_unitless["z_mean"] = 1 if "to_hz" in inputs: inputs_unitless["to_hz"] = inputs["to_hz"] lite = plasma_frequency_lite(**inputs_unitless) pylite = plasma_frequency_lite.py_func(**inputs_unitless) assert pylite == lite normal = plasma_frequency(**inputs) assert np.allclose(normal.value, lite)
def test_values(self, args, kwargs, expected, rtol): """Test various expected values.""" wp = plasma_frequency(*args, **kwargs) assert isinstance(wp, u.Quantity) assert wp.unit == u.rad / u.s assert np.allclose(wp.value, expected, rtol=rtol)
def test_upper_hybrid_frequency(): r"""Test the upper_hybrid_frequency function in frequencies.py.""" omega_uh = upper_hybrid_frequency(B, n_e=n_e) omega_uh_hz = upper_hybrid_frequency(B, n_e=n_e, to_hz=True) omega_ce = gyrofrequency(B, "e-") omega_pe = plasma_frequency(n=n_e, particle="e-") assert omega_ce.unit.is_equivalent(u.rad / u.s) assert omega_pe.unit.is_equivalent(u.rad / u.s) assert omega_uh.unit.is_equivalent(u.rad / u.s) assert omega_uh_hz.unit.is_equivalent(u.Hz) left_hand_side = omega_uh**2 right_hand_side = omega_ce**2 + omega_pe**2 assert np.isclose(left_hand_side.value, right_hand_side.value) assert np.isclose(omega_uh_hz.value, 69385868857.90918) with pytest.raises(ValueError): upper_hybrid_frequency(5 * u.T, n_e=-1 * u.m**-3) with pytest.warns(u.UnitsWarning): assert upper_hybrid_frequency(1.2, 1.3) == upper_hybrid_frequency( 1.2 * u.T, 1.3 * u.m**-3) with pytest.warns(u.UnitsWarning): assert upper_hybrid_frequency(1.4 * u.T, 1.3) == upper_hybrid_frequency( 1.4, 1.3 * u.m**-3) assert_can_handle_nparray(upper_hybrid_frequency)
def test_normal_vs_lite_values(self, kwargs, expected): """ Test that `permittivity_1D_Maxwellian_lite` and `permittivity_1D_Maxwellian` calculate the same values. """ wp = plasma_frequency(kwargs["n"], kwargs["particle"], kwargs["z_mean"]) vth = thermal_speed(kwargs["T"], kwargs["particle"], method="most_probable") kwargs["kWave"] = kwargs["omega"] / vth val = permittivity_1D_Maxwellian(**kwargs) val_lite = permittivity_1D_Maxwellian_lite( kwargs["omega"].value, kwargs["kWave"].to(u.rad / u.m).value, vth.value, wp.value, ) assert ( np.isclose(val, val_lite, rtol=1e-6, atol=0.0), "'permittivity_1D_Maxwellian' and 'permittivity_1D_Maxwellian_lite' " "do not agree.", )
def test_warns(self, args, kwargs, _warning, expected): """ Test scenarios the cause plasma_frequency to issue a warning. """ with pytest.warns(_warning): wp = plasma_frequency(*args, **kwargs) assert isinstance(wp, u.Quantity) assert wp.unit == u.rad / u.s if expected is not None: assert np.allclose(wp, expected)
def test_proton_electron_plasma(self): """ Test proton-electron plasma against the (approximate) analytical formulas """ B = 1 * u.T n = [1, 1] * 1 / u.m ** 3 omega = 1 * u.rad / u.s omega_ce = gyrofrequency(B, particle="e", signed=True) omega_pe = plasma_frequency(n[0], particle="e") omega_cp = abs(omega_ce) / 1860 omega_pp = omega_pe / 43 S_analytical = ( 1 - omega_pe ** 2 / (omega ** 2 - omega_ce ** 2) - omega_pp ** 2 / (omega ** 2 - omega_cp ** 2) ) D_analytical = +omega_ce / omega * omega_pe ** 2 / ( omega ** 2 - omega_ce ** 2 ) + omega_cp / omega * omega_pp ** 2 / (omega ** 2 - omega_cp ** 2) P_analytical = 1 - (omega_pe ** 2 + omega_pp ** 2) / omega ** 2 species = ["e", "p"] S, D, P = tuple_result = cold_plasma_permittivity_SDP(B, species, n, omega) assert tuple_result.sum is S assert tuple_result.difference is D assert tuple_result.plasma is P assert isinstance(tuple_result, StixTensorElements) assert np.isclose(S, S_analytical) assert np.isclose(D, D_analytical) assert np.isclose(P, P_analytical) L, R, P = rotating_tuple_result = cold_plasma_permittivity_LRP( B, species, n, omega ) assert rotating_tuple_result.left is L assert rotating_tuple_result.right is R assert rotating_tuple_result.plasma is P assert isinstance(rotating_tuple_result, RotatingTensorElements)
def test_lower_hybrid_frequency(): r"""Test the lower_hybrid_frequency function in frequencies.py.""" ion = "He-4 1+" omega_ci = gyrofrequency(B, particle=ion) omega_pi = plasma_frequency(n=n_i, particle=ion) omega_ce = gyrofrequency(B, "e-") omega_lh = lower_hybrid_frequency(B, n_i=n_i, ion=ion) omega_lh_hz = lower_hybrid_frequency(B, n_i=n_i, ion=ion, to_hz=True) assert omega_ci.unit.is_equivalent(u.rad / u.s) assert omega_pi.unit.is_equivalent(u.rad / u.s) assert omega_ce.unit.is_equivalent(u.rad / u.s) assert omega_lh.unit.is_equivalent(u.rad / u.s) left_hand_side = omega_lh**-2 right_hand_side = (1 / (omega_ci**2 + omega_pi**2) + omega_ci**-1 * omega_ce**-1) assert np.isclose(left_hand_side.value, right_hand_side.value) assert np.isclose(omega_lh_hz.value, 299878691.3223296) with pytest.raises(ValueError): lower_hybrid_frequency(0.2 * u.T, n_i=5e19 * u.m**-3, ion="asdfasd") with pytest.raises(ValueError): lower_hybrid_frequency(0.2 * u.T, n_i=-5e19 * u.m**-3, ion="asdfasd") with pytest.raises(ValueError): lower_hybrid_frequency(np.nan * u.T, n_i=-5e19 * u.m**-3, ion="asdfasd") with pytest.warns(u.UnitsWarning): assert lower_hybrid_frequency(1.3, 1e19, "p+") == lower_hybrid_frequency( 1.3 * u.T, 1e19 * u.m**-3, "p+") assert_can_handle_nparray(lower_hybrid_frequency)
def inertial_length(n: u.m**-3, particle: Particle) -> u.m: r""" Calculate a charged particle's inertial length. **Aliases:** `cwp_` Parameters ---------- n : `~astropy.units.Quantity` Particle number density in units convertible to m\ :sup:`-3`\ . particle : `~plasmapy.particles.particle_class.Particle` Representation of the particle species (e.g., 'p+' for protons, 'D+' for deuterium, or 'He-4 +1' for singly ionized helium-4). Returns ------- d : `~astropy.units.Quantity` The particle's inertial length in meters. Raises ------ `TypeError` If ``n`` is not a `~astropy.units.Quantity` or ``particle`` is not a string. `~astropy.units.UnitConversionError` If ``n`` is not in units of a number density. `ValueError` The particle density does not have an appropriate value. Warns ----- : `~astropy.units.UnitsWarning` If units are not provided and SI units are assumed. Notes ----- The inertial length of a particle of species :math:`s` is given by .. math:: d = \frac{c}{ω_{ps}} The inertial length is the characteristic length scale for a particle to be accelerated in a plasma. The Hall effect becomes important on length scales shorter than the ion inertial length. The inertial length is also known as the skin depth. Examples -------- >>> from astropy import units as u >>> inertial_length(5 * u.m ** -3, 'He+') <Quantity 2.02985...e+08 m> >>> inertial_length(5 * u.m ** -3, 'e-') <Quantity 2376534.75... m> """ omega_p = frequencies.plasma_frequency(n, particle=particle) return c / omega_p
def thermal_bremsstrahlung( frequencies: u.Hz, n_e: u.m**-3, T_e: u.K, n_i: u.m**-3 = None, ion_species: Particle = "H+", kmax: u.m = None, ) -> np.ndarray: r""" Calculate the bremsstrahlung emission spectrum for a Maxwellian plasma in the Rayleigh-Jeans limit :math:`ℏ ω ≪ k_B T_e` .. math:: \frac{dP}{dω} = \frac{8 \sqrt{2}}{3\sqrt{π}} \bigg ( \frac{e^2}{4 π ε_0} \bigg )^3 \bigg ( m_e c^2 \bigg )^{-\frac{3}{2}} \bigg ( 1 - \frac{ω_{pe}^2}{ω^2} \bigg )^\frac{1}{2} \frac{Z_i^2 n_i n_e}{\sqrt(k_B T_e)} E_1(y) where :math:`E_1` is the exponential integral .. math:: E_1 (y) = - \int_{-y}^∞ \frac{e^{-t}}{t}dt and :math:`y` is the dimensionless argument .. math:: y = \frac{1}{2} \frac{ω^2 m_e}{k_{max}^2 k_B T_e} where :math:`k_{max}` is a maximum wavenumber approximated here as :math:`k_{max} = 1/λ_B` where :math:`λ_B` is the electron de Broglie wavelength. Parameters ---------- frequencies : `~astropy.units.Quantity` Array of frequencies over which the bremsstrahlung spectrum will be calculated (convertible to Hz). n_e : `~astropy.units.Quantity` Electron number density in the plasma (convertible to m\ :sup:`-3`\ ). T_e : `~astropy.units.Quantity` Temperature of the electrons (in K or convertible to eV). n_i : `~astropy.units.Quantity`, optional Ion number density in the plasma (convertible to m\ :sup:`-3`\ ). Defaults to the quasi-neutral condition :math:`n_i = n_e / Z`\ . ion_species : `str` or `~plasmapy.particles.particle_class.Particle`, optional An instance of `~plasmapy.particles.particle_class.Particle`, or a string convertible to `~plasmapy.particles.particle_class.Particle`. kmax : `~astropy.units.Quantity` Cutoff wavenumber (convertible to radians per meter). Defaults to the inverse of the electron de Broglie wavelength. Returns ------- spectrum : `~astropy.units.Quantity` Computed bremsstrahlung spectrum over the frequencies provided. Notes ----- For details, see :cite:t:`bekefi:1966`\ . """ # Default n_i is n_e/Z: if n_i is None: n_i = n_e / ion_species.charge_number # Default value of kmax is the electrom thermal de Broglie wavelength if kmax is None: kmax = (np.sqrt(const.m_e.si * const.k_B.si * T_e) / const.hbar.si).to( 1 / u.m) # Convert frequencies to angular frequencies ω = (frequencies * 2 * np.pi * u.rad).to(u.rad / u.s) # Calculate the electron plasma frequency ω_pe = plasma_frequency(n=n_e, particle="e-") # Check that all ω < wpe (this formula is only valid in this limit) if np.min(ω) < ω_pe: raise PhysicsError( "Lowest frequency must be larger than the electron " f"plasma frequency {ω_pe:.1e}, but min(ω) = {np.min(ω):.1e}") # Check that the parameters given fall within the Rayleigh-Jeans limit # hω << kT_e rj_const = (np.max(ω) * const.hbar.si / (2 * np.pi * u.rad * const.k_B.si * T_e)).to( u.dimensionless_unscaled) if rj_const.value > 0.1: raise PhysicsError("Rayleigh-Jeans limit not satisfied: " f"ℏω/kT_e = {rj_const.value:.2e} > 0.1. " "Try lower ω or higher T_e.") # Calculate the bremsstrahlung power spectral density in several steps c1 = ((8 / 3) * np.sqrt(2 / np.pi) * (const.e.si**2 / (4 * np.pi * const.eps0.si))**3 * 1 / (const.m_e.si * const.c.si**2)**1.5) Zi = ion_species.charge_number c2 = (np.sqrt(1 - ω_pe**2 / ω**2) * Zi**2 * n_i * n_e / np.sqrt(const.k_B.si * T_e)) # Dimensionless argument for exponential integral arg = 0.5 * ω**2 * const.m_e.si / (kmax**2 * const.k_B.si * T_e) / u.rad**2 # Remove units, get ndarray of values arg = (arg.to(u.dimensionless_unscaled)).value return c1 * c2 * exp1(arg)
def two_fluid( *, B: u.T, ion: Union[str, Particle], k: u.rad / u.m, n_i: u.m**-3, T_e: u.K, T_i: u.K, theta: u.rad, gamma_e: Union[float, int] = 1, gamma_i: Union[float, int] = 3, z_mean: Union[float, int] = None, ): r""" Using the solution provided by :cite:t:`bellan:2012`, calculate the analytical solution to the two fluid, low-frequency (:math:`\omega/kc \ll 1`) dispersion relation presented by :cite:t:`stringer:1963`. This dispersion relation also assumes a uniform magnetic field :math:`\mathbf{B_o}`, no D.C. electric field :math:`\mathbf{E_o}=0`, and quasi-neutrality. For more information see the **Notes** section below. Parameters ---------- B : `~astropy.units.Quantity` The magnetic field magnitude in units convertible to T. ion : `str` or `~plasmapy.particles.particle_class.Particle` Representation of the ion species (e.g., ``'p'`` for protons, ``'D+'`` for deuterium, ``'He-4 +1'`` for singly ionized helium-4, etc.). If no charge state information is provided, then the ions are assumed to be singly ionized. k : `~astropy.units.Quantity`, single valued or 1-D array Wavenumber in units convertible to rad/m`. Either single valued or 1-D array of length :math:`N`. n_i : `~astropy.units.Quantity` Ion number density in units convertible to m\ :sup:`-3`. T_e : `~astropy.units.Quantity` The electron temperature in units of K or eV. T_i : `~astropy.units.Quantity` The ion temperature in units of K or eV. theta : `~astropy.units.Quantity`, single valued or 1-D array The angle of propagation of the wave with respect to the magnetic field, :math:`\cos^{-1}(k_z / k)`, in units must be convertible to radians. Either single valued or 1-D array of size :math:`M`. gamma_e : `float` or `int`, optional 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`, optional 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. z_mean : `float` or int, optional The average ionization state (arithmetic mean) of the ``ion`` composing the plasma. Will override any charge state defined by argument ``ion``. Returns ------- omega : Dict[str, `~astropy.units.Quantity`] A dictionary of computed wave frequencies in units rad/s. The dictionary contains three keys: ``'fast_mode'`` for the fast mode, ``'alfven_mode'`` for the Alfvén mode, and ``'acoustic_mode'`` for the ion-acoustic mode. The value for each key will be a :math:`N × M` array. Raises ------ TypeError If applicable arguments are not instances of `~astropy.units.Quantity` or cannot be converted into one. TypeError If ``ion`` is not of type or convertible to `~plasmapy.particles.particle_class.Particle`. TypeError If ``gamma_e``, ``gamma_i``, or ``z_mean`` are not of type `int` or `float`. ~astropy.units.UnitTypeError If applicable arguments do not have units convertible to the expected units. ValueError If any of ``B``, ``k``, ``n_i``, ``T_e``, or ``T_i`` is negative. ValueError If ``k`` is negative or zero. ValueError If ``ion`` is not of category ion or element. ValueError If ``B``, ``n_i``, ``T_e``, or ``T_I`` are not single valued `astropy.units.Quantity` (i.e. an array). ValueError If ``k`` or ``theta`` are not single valued or a 1-D array. Warns ----- : `~plasmapy.utils.exceptions.PhysicsWarning` When the computed wave frequencies violate the low-frequency (:math:`\omega/kc \ll 1`) assumption of the dispersion relation. Notes ----- The complete dispersion equation presented by :cite:t:`stringer:1963` (equation 1 of :cite:t:`bellan:2012`) is: .. math:: \left( \cos^2 \theta - Q \frac{\omega^2}{k^2 {v_A}^2} \right) & \left[ \left( \cos^2 \theta - \frac{\omega^2}{k^2 {c_s}^2} \right) - Q \frac{\omega^2}{k^2 {v_A}^2} \left( 1 - \frac{\omega^2}{k^2 {c_s}^2} \right) \right] \\ &= \left(1 - \frac{\omega^2}{k^2 {c_s}^2} \right) \frac{\omega^2}{{\omega_{ci}}^2} \cos^2 \theta where .. math:: Q &= 1 + k^2 c^2/{\omega_{pe}}^2 \\ \cos \theta &= \frac{k_z}{k} \\ \mathbf{B_o} &= B_{o} \mathbf{\hat{z}} :math:`\omega` is the wave frequency, :math:`k` is the wavenumber, :math:`v_A` is the Alfvén velocity, :math:`c_s` is the sound speed, :math:`\omega_{ci}` is the ion gyrofrequency, and :math:`\omega_{pe}` is the electron plasma frequency. This relation does additionally assume low-frequency waves :math:`\omega/kc \ll 1`, no D.C. electric field :math:`\mathbf{E_o}=0` and quasi-neutrality. Following section 5 of :cite:t:`bellan:2012` the exact roots of the above dispersion equation can be derived and expressed as one analytical solution (equation 38 of :cite:t:`bellan:2012`): .. math:: \frac{\omega}{\omega_{ci}} = \sqrt{ 2 \Lambda \sqrt{-\frac{P}{3}} \cos\left( \frac{1}{3} \cos^{-1}\left( \frac{3q}{2p} \sqrt{-\frac{3}{p}} \right) - \frac{2 \pi}{3}j \right) + \frac{\Lambda A}{3} } where :math:`j = 0` represents the fast mode, :math:`j = 1` represents the Alfvén mode, and :math:`j = 2` represents the acoustic mode. Additionally, .. math:: p &= \frac{3B-A^2}{3} \; , \; q = \frac{9AB-2A^3-27C}{27} \\ A &= \frac{Q + Q^2 \beta + Q \alpha + \alpha \Lambda}{Q^2} \; , \; B = \alpha \frac{1 + 2 Q \beta + \Lambda \beta}{Q^2} \; , \; C = \frac{\alpha^2 \beta}{Q^2} \\ \alpha &= \cos^2 \theta \; , \; \beta = \left( \frac{c_s}{v_A}\right)^2 \; , \; \Lambda = \left( \frac{k v_{A}}{\omega_{ci}}\right)^2 Examples -------- >>> from astropy import units as u >>> from plasmapy.dispersion.analytical import two_fluid >>> inputs = { ... "k": 0.01 * u.rad / u.m, ... "theta": 30 * u.deg, ... "B": 8.3e-9 * u.T, ... "n_i": 5e6 * u.m ** -3, ... "T_e": 1.6e6 * u.K, ... "T_i": 4.0e5 * u.K, ... "ion": "p+", ... } >>> omegas = two_fluid(**inputs) >>> omegas {'fast_mode': <Quantity 1520.57... rad / s>, 'alfven_mode': <Quantity 1261.75... rad / s>, 'acoustic_mode': <Quantity 0.688152... rad / s>} >>> inputs = { ... "k": [1e-7, 2e-7] * u.rad / u.m, ... "theta": [10, 20] * u.deg, ... "B": 8.3e-9 * u.T, ... "n_i": 5e6 * u.m ** -3, ... "T_e": 1.6e6 * u.K, ... "T_i": 4.0e5 * u.K, ... "ion": "He+", ... } >>> omegas = two_fluid(**inputs) >>> omegas['fast_mode'] <Quantity [[0.00767..., 0.00779... ], [0.01534..., 0.01558...]] rad / s> """ # validate argument ion if not isinstance(ion, Particle): try: ion = Particle(ion) except TypeError: raise TypeError( f"For argument 'ion' expected type {Particle} but got {type(ion)}." ) if not ion.is_ion and not ion.is_category("element"): raise ValueError( "The particle passed for 'ion' must be an ion or element.") # validate z_mean if z_mean is None: try: z_mean = abs(ion.charge_number) except ChargeError: z_mean = 1 elif isinstance(z_mean, (int, np.integer, float, np.floating)): z_mean = abs(z_mean) else: raise TypeError( f"Expected int or float for argument 'z_mean', but got {type(z_mean)}." ) # validate arguments for arg_name in ("B", "n_i", "T_e", "T_i"): val = locals()[arg_name].squeeze() if val.shape != (): raise ValueError( f"Argument '{arg_name}' must a single value and not an array of " f"shape {val.shape}.") locals()[arg_name] = val # validate arguments for arg_name in ("gamma_e", "gamma_i"): if not isinstance(locals()[arg_name], (int, np.integer, float, np.floating)): raise TypeError( f"Expected int or float for argument '{arg_name}', but got " f"{type(locals()[arg_name])}.") # validate argument k k = k.squeeze() if k.ndim not in [0, 1]: raise ValueError( f"Argument 'k' needs to be a single valued or 1D array astropy Quantity," f" got array of shape {k.shape}.") if np.any(k <= 0): raise ValueError("Argument 'k' can not be a or have negative values.") # validate argument theta theta = theta.squeeze() if theta.ndim not in [0, 1]: raise ValueError( f"Argument 'theta' needs to be a single valued or 1D array astropy " f"Quantity, got array of shape {k.shape}.") # Calc needed plasma parameters n_e = z_mean * n_i with warnings.catch_warnings(): warnings.simplefilter("ignore", category=PhysicsWarning) c_s = ion_sound_speed( T_e=T_e, T_i=T_i, ion=ion, n_e=n_e, gamma_e=gamma_e, gamma_i=gamma_i, z_mean=z_mean, ) v_A = Alfven_speed(B, n_i, ion=ion, z_mean=z_mean) omega_ci = gyrofrequency(B=B, particle=ion, signed=False, Z=z_mean) omega_pe = plasma_frequency(n=n_e, particle="e-") # Bellan2012JGR params equation 32 alpha = np.cos(theta.value)**2 beta = (c_s / v_A).to(u.dimensionless_unscaled).value**2 alphav, kv = np.meshgrid(alpha, k.value) # create grid Lambda = (kv * v_A.value / omega_ci.value)**2 # Bellan2012JGR params equation 2 Q = 1 + (kv * c.value / omega_pe.value)**2 # Bellan2012JGR params equation 35 A = ((1 + alphav) / Q) + beta + (alphav * Lambda / Q**2) B = alphav * (1 + 2 * Q * beta + Lambda * beta) / Q**2 C = beta * (alphav / Q)**2 # Bellan2012JGR params equation 36 p = (3 * B - A**2) / 3 q = (9 * A * B - 2 * A**3 - 27 * C) / 27 # Bellan2012JGR params equation 38 R = 2 * Lambda * np.emath.sqrt(-p / 3) S = 3 * q / (2 * p) * np.emath.sqrt(-3 / p) T = Lambda * A / 3 omega = {} for ind, key in enumerate(("fast_mode", "alfven_mode", "acoustic_mode")): # The solution corresponding to equation 38 w = omega_ci * np.emath.sqrt( R * np.cos(1 / 3 * np.emath.arccos(S) - 2 * np.pi / 3 * ind) + T) omega[key] = w.squeeze() # check for violation of dispersion relation assumptions # (i.e. low-frequency, w/kc << 0.1) wkc_max = np.max(w.value / (kv * c.value)) if wkc_max > 0.1: warnings.warn( f"The {key} calculation produced a high-frequency wave (w/kc == " f"{wkc_max:.3f}), which violates the low-frequency (w/kc << 1) " f"assumption of the dispersion relation.", PhysicsWarning, ) return omega
def hollweg( *, B: u.T, ion: Union[str, Particle], k: u.rad / u.m, n_i: u.m**-3, T_e: u.K, T_i: u.K, theta: u.rad, gamma_e: Union[float, int] = 1, gamma_i: Union[float, int] = 3, z_mean: Union[float, int] = None, ): r""" Calculate the two fluid dispersion relation presented by :cite:t:`hollweg:1999`, and discussed by :cite:t:`bellan:2012`. This is a numberical solver of equation 3 in :cite:t:`bellan:2012`. See the **Notes** section below for additional details. Parameters ---------- B : `~astropy.units.Quantity` The magnetic field magnitude in units convertible to T. ion : `str` or `~plasmapy.particles.particle_class.Particle` Representation of the ion species (e.g., ``'p'`` for protons, ``'D+'`` for deuterium, ``'He-4 +1'`` for singly ionized helium-4, etc.). If no charge state information is provided, then the ions are assumed to be singly ionized. k : `~astropy.units.Quantity`, single valued or 1-D array Wavenumber in units convertible to rad/m. Either single valued or 1-D array of length :math:`N`. n_i : `~astropy.units.Quantity` Ion number density in units convertible to m\ :sup:`-3`. T_e : `~astropy.units.Quantity` The electron temperature in units of K or eV. T_i : `~astropy.units.Quantity` The ion temperature in units of K or eV. theta : `~astropy.units.Quantity`, single valued or 1-D array The angle of propagation of the wave with respect to the magnetic field, :math:`\cos^{-1}(k_z / k)`, in units convertible to radians. Either single valued or 1-D array of size :math:`M`. gamma_e : `float` or `int`, optional 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`, optional 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. z_mean : `float` or int, optional The average ionization state (arithmetic mean) of the ``ion`` composing the plasma. Will override any charge state defined by argument ``ion``. Returns ------- omega : Dict[str, `~astropy.units.Quantity`] A dictionary of computed wave frequencies in units rad/s. The dictionary contains three keys: ``'fast_mode'`` for the fast mode, ``'alfven_mode'`` for the Alfvén mode, and ``'acoustic_mode'`` for the ion-acoustic mode. The value for each key will be a :math:`N x M` array. Raises ------ TypeError If applicable arguments are not instances of `~astropy.units.Quantity` or cannot be converted into one. TypeError If ``ion`` is not of type or convertible to `~plasmapy.particles.particle_class.Particle`. TypeError If ``gamma_e``, ``gamma_i``, or ``z_mean`` are not of type `int` or `float`. ~astropy.units.UnitTypeError If applicable arguments do not have units convertible to the expected units. ValueError If any of ``B``, ``k``, ``n_i``, ``T_e``, or ``T_i`` is negative. ValueError If ``k`` is negative or zero. ValueError If ``ion`` is not of category ion or element. ValueError If ``B``, ``n_i``, ``T_e``, or ``T_I`` are not single valued `astropy.units.Quantity` (i.e. an array). ValueError If ``k`` or ``theta`` are not single valued or a 1-D array. Warns ----- : `~plasmapy.utils.exceptions.PhysicsWarning` When :math:`\omega / \omega_{\rm ci} > 0.1`, violation of the low-frequency assumption. : `~plasmapy.utils.exceptions.PhysicsWarning` When :math:`c_{\rm s} / v_{\rm A} > 0.1`, violation of low-β. : `~plasmapy.utils.exceptions.PhysicsWarning` When :math:`|θ - π/2| > 0.1`, violation of quasi-perpendicular propagation. Notes ----- The dispersion relation presented in :cite:t:`hollweg:1999` (equation 3 in :cite:t:`bellan:2012`) is: .. math:: \left( \frac{\omega^2}{k_{\rm z}^2 v_{\rm A}^2} - 1 \right) & \left[ \omega^2 \left( \omega^2 - k^2 v_{\rm A}^2 \right) - \beta k^2 v_{\rm A}^2 \left( \omega^2 - k_{\rm z}^2 v_{\rm A}^2 \right) \right] \\ &= \omega^2 \left(\omega^2 - k^2 v_{\rm A}^2 \right) k_{\rm x}^2 \left( \frac{c_{\rm s}^2}{\omega_{\rm ci}^2} - \frac{c^2}{\omega_{\rm pe}^2} \frac{\omega^2}{k_{\rm z}^2v_{\rm A}^2} \right) where .. math:: \mathbf{B_o} &= B_{o} \mathbf{\hat{z}} \\ \cos \theta &= \frac{k_z}{k} \\ \mathbf{k} &= k_{\rm x} \hat{x} + k_{\rm z} \hat{z} :math:`\omega` is the wave frequency, :math:`k` is the wavenumber, :math:`v_{\rm A}` is the Alfvén velocity, :math:`c_{\rm s}` is the sound speed, :math:`\omega_{\rm ci}` is the ion gyrofrequency, and :math:`\omega_{\rm pe}` is the electron plasma frequency. In the derivation of this relation Hollweg assumed low-frequency waves :math:`\omega / \omega_{\rm ci} \ll 1`, no D.C. electric field :math:`\mathbf{E_o}=0`, and quasi-neutrality. :cite:t:`hollweg:1999` asserts this expression is valid for arbitrary :math:`c_{\rm s} / v_{\rm A}` (β) and :math:`k_{\rm z} / k` (θ). Contrarily, :cite:t:`bellan:2012` states in §1.7 that due to the inconsistent retention of the :math:`\omega / \omega_{\rm ci} \ll 1` terms the expression can only be valid if both :math:`c_{\rm s} \ll v_{\rm A}` (low-β) and the wave propgation is nearly perpendicular to the magnetic field. This routine solves for ω for given :math:`k` values by numerically solving for the roots of the above expression. Examples -------- >>> from astropy import units as u >>> from plasmapy.dispersion.numerical import hollweg_ >>> inputs = { ... "k": np.logspace(-7, -2, 2) * u.rad / u.m, ... "theta": 88 * u.deg, ... "n_i": 5 * u.cm ** -3, ... "B": 2.2e-8 * u.T, ... "T_e": 1.6e6 * u.K, ... "T_i": 4.0e5 * u.K, ... "ion": Particle("p+"), ... } >>> omegas = hollweg(**inputs) >>> omegas {'fast_mode': <Quantity [2.62911663e-02+0.j, 2.27876968e+03+0.j] rad / s>, 'alfven_mode': <Quantity [7.48765909e-04+0.j, 2.13800404e+03+0.j] rad / s>, 'acoustic_mode': <Quantity [0.00043295+0.j, 0.07358991+0.j] rad / s>} """ # validate argument ion if not isinstance(ion, Particle): try: ion = Particle(ion) except TypeError: raise TypeError( f"For argument 'ion' expected type {Particle} but got {type(ion)}." ) if not (ion.is_ion or ion.is_category("element")): raise ValueError( "The particle passed for 'ion' must be an ion or element.") # validate z_mean if z_mean is None: try: z_mean = abs(ion.charge_number) except ChargeError: z_mean = 1 else: if not isinstance(z_mean, (int, np.integer, float, np.floating)): raise TypeError( f"Expected int or float for argument 'z_mean', but got {type(z_mean)}." ) z_mean = abs(z_mean) # validate arguments for arg_name in ("B", "n_i", "T_e", "T_i"): val = locals()[arg_name].squeeze() if val.shape != (): raise ValueError( f"Argument '{arg_name}' must be single valued and not an array of " f"shape {val.shape}.") locals()[arg_name] = val # validate arguments for arg_name in ("gamma_e", "gamma_i"): if not isinstance(locals()[arg_name], (int, np.integer, float, np.floating)): raise TypeError( f"Expected int or float for argument '{arg_name}', but got " f"{type(locals()[arg_name])}.") # validate argument k k = k.squeeze() if not (k.ndim == 0 or k.ndim == 1): raise ValueError( f"Argument 'k' needs to be single valued or a 1D array astropy Quantity," f" got array of shape {k.shape}.") if np.any(k <= 0): raise ValueError("Argument 'k' can not be a or have negative values.") # validate argument theta theta = theta.squeeze() if theta.ndim not in (0, 1): raise ValueError( f"Argument 'theta' needs to be a single valued or 1D array astropy " f"Quantity, got array of shape {theta.shape}.") # Single k value case if np.isscalar(k.value): k = np.array([k.value]) * u.rad / u.m # Calc needed plasma parameters with warnings.catch_warnings(): warnings.simplefilter("ignore", category=PhysicsWarning) n_e = z_mean * n_i c_s = ion_sound_speed( T_e=T_e, T_i=T_i, ion=ion, n_e=n_e, gamma_e=gamma_e, gamma_i=gamma_i, z_mean=z_mean, ).value v_A = Alfven_speed(B, n_i, ion=ion, z_mean=z_mean).value omega_ci = gyrofrequency(B=B, particle=ion, signed=False, Z=z_mean).value omega_pe = plasma_frequency(n=n_e, particle="e-").value cs_vA = c_s / v_A thetav, kv = np.meshgrid(theta.value, k.value) # Parameters kx and kz kz = np.cos(thetav) * kv kx = np.sin(thetav) * kv # Define helpful parameters beta = (c_s / v_A)**2 alpha_A = (kv * v_A)**2 alpha_s = (kv * c_s)**2 # == alpha_A * beta sigma = (kz * v_A)**2 D = (c_s / omega_ci)**2 F = (c_si_unitless / omega_pe)**2 # Polynomial coefficients: c3*x^3 + c2*x^2 + c1*x + c0 = 0 c3 = F * kx**2 + 1 c2 = -alpha_A * (1 + beta + F * kx**2) - sigma * (1 + D * kx**2) c1 = sigma * alpha_A * (1 + 2 * beta + D * kx**2) c0 = -alpha_s * sigma**2 # Find roots to polynomial coefficients = np.array([c3, c2, c1, c0], ndmin=3) nroots = coefficients.shape[0] - 1 # 3 nks = coefficients.shape[1] nthetas = coefficients.shape[2] roots = np.empty((nroots, nks, nthetas), dtype=np.complex128) for ii in range(nks): for jj in range(nthetas): roots[:, ii, jj] = np.roots(coefficients[:, ii, jj]) roots = np.sqrt(roots) roots = np.sort(roots, axis=0) # Warn about NOT low-beta if c_s / v_A > 0.1: warnings.warn( f"This solver is valid in the low-beta regime, " f"c_s/v_A << 1 according to Bellan, 2012, Sec. 1.7 " f"(see documentation for DOI). A c_s/v_A value of {cs_vA:.2f} " f"was calculated which may affect the validity of the solution.", PhysicsWarning, ) # Warn about theta not nearly perpendicular theta_diff_max = np.amax(np.abs(thetav - np.pi / 2)) if theta_diff_max > 0.1: warnings.warn( f"This solver is valid in the regime where propagation is " f"nearly perpendicular to B according to Bellan, 2012, Sec. 1.7 " f"(see documentation for DOI). A |theta - pi/2| value of " f"{theta_diff_max:.2f} was calculated which may affect the " f"validity of the solution.", PhysicsWarning, ) # dispersion relation is only valid in the regime w << w_ci w_max = np.max(roots) w_wci_max = w_max / omega_ci if w_wci_max > 0.1: warnings.warn( f"This solver is valid in the regime w/w_ci << 1. A w " f"value of {w_max:.2f} and a w/w_ci value of " f"{w_wci_max:.2f} were calculated which may affect the " f"validity of the solution.", PhysicsWarning, ) omegas = { "fast_mode": roots[2, :].squeeze() * u.rad / u.s, "alfven_mode": roots[1, :].squeeze() * u.rad / u.s, "acoustic_mode": roots[0, :].squeeze() * u.rad / u.s, } return omegas
theta = theta.value.squeeze() if theta.ndim not in (0, 1): raise TypeError( "Argument 'theta' needs to be a single value or 1D array " f" astropy Quantity, got array of shape {theta.shape}.") elif np.isscalar(theta): theta = np.array([theta]) # Generate mesh grid of w x theta w, theta = np.meshgrid(w, theta, indexing="ij") # Generate the plasma parameters needed wps = [] wcs = [] for par, dens in zip(species, densities): wps.append(plasma_frequency(n=dens * u.m**-3, particle=par).value) wcs.append(gyrofrequency(B=B, particle=par, signed=True).value) # Stix method implemented S = np.ones_like(w, dtype=np.float64) P = np.ones_like(S) D = np.zeros_like(S) for wc, wp in zip(wcs, wps): S -= (wp**2) / (w**2 - wc**2) P -= (wp / w)**2 D += ((wp**2) / (w**2 - wc**2)) * (wc / w) R = S + D L = S - D # Generate coefficients to solve, a * k**4 + b * k**2 + c = 0
def cold_plasma_permittivity_SDP(B: u.T, species, n, omega: u.rad / u.s): r""" Magnetized cold plasma dielectric permittivity tensor elements. Elements (S, D, P) are given in the "Stix" frame, i.e. with :math:`B ∥ \hat{z}` :cite:p:`stix:1992`. The :math:`\exp(-i ω t)` time-harmonic convention is assumed. Parameters ---------- B : `~astropy.units.Quantity` Magnetic field magnitude in units convertible to tesla. species : `list` of `str` List of the plasma particle species, e.g.: ``['e', 'D+']`` or ``['e', 'D+', 'He+']``. n : `list` of `~astropy.units.Quantity` `list` of species density in units convertible to per cubic meter The order of the species densities should follow species. omega : `~astropy.units.Quantity` Electromagnetic wave frequency in rad/s. Returns ------- sum : `~astropy.units.Quantity` S ("Sum") dielectric tensor element. difference : `~astropy.units.Quantity` D ("Difference") dielectric tensor element. plasma : `~astropy.units.Quantity` P ("Plasma") dielectric tensor element. Notes ----- The dielectric permittivity tensor is expressed in the Stix frame with the :math:`\exp(-i ω t)` time-harmonic convention as :math:`ε = ε_0 A`, with :math:`A` being .. math:: ε = ε_0 \left(\begin{matrix} S & -i D & 0 \\ +i D & S & 0 \\ 0 & 0 & P \end{matrix}\right) where: .. math:: S = 1 - \sum_s \frac{ω_{p,s}^2}{ω^2 - Ω_{c,s}^2} D = \sum_s \frac{Ω_{c,s}}{ω} \frac{ω_{p,s}^2}{ω^2 - Ω_{c,s}^2} P = 1 - \sum_s \frac{ω_{p,s}^2}{ω^2} where :math:`ω_{p,s}` is the plasma frequency and :math:`Ω_{c,s}` is the signed version of the cyclotron frequency for the species :math:`s`. Examples -------- >>> from astropy import units as u >>> from numpy import pi >>> B = 2*u.T >>> species = ['e', 'D+'] >>> n = [1e18*u.m**-3, 1e18*u.m**-3] >>> omega = 3.7e9*(2*pi)*(u.rad/u.s) >>> permittivity = S, D, P = cold_plasma_permittivity_SDP(B, species, n, omega) >>> S <Quantity 1.02422...> >>> permittivity.sum # namedtuple-style access <Quantity 1.02422...> >>> D <Quantity 0.39089...> >>> P <Quantity -4.8903...> """ S, D, P = 1, 0, 1 for s, n_s in zip(species, n): omega_c = gyrofrequency(B=B, particle=s, signed=True) omega_p = plasma_frequency(n=n_s, particle=s) S += -(omega_p**2) / (omega**2 - omega_c**2) D += omega_c / omega * omega_p**2 / (omega**2 - omega_c**2) P += -(omega_p**2) / omega**2 return StixTensorElements(S, D, P)
def permittivity_1D_Maxwellian( omega: u.rad / u.s, kWave: u.rad / u.m, T: u.K, n: u.m**-3, particle, z_mean: u.dimensionless_unscaled = None, ) -> u.dimensionless_unscaled: r""" Compute the classical dielectric permittivity for a 1D Maxwellian plasma. This function can calculate both the ion and electron permittivities. No additional effects are considered (e.g. magnetic fields, relativistic effects, strongly coupled regime, etc.). Parameters ---------- omega : `~astropy.units.Quantity` The frequency, in rad/s, of the electromagnetic wave propagating through the plasma. kWave : `~astropy.units.Quantity` The corresponding wavenumber, in rad/m, of the electromagnetic wave propagating through the plasma. T : `~astropy.units.Quantity` The plasma temperature — this can be either the electron or the ion temperature, but should be consistent with density and particle. n : `~astropy.units.Quantity` The plasma density — this can be either the electron or the ion density, but should be consistent with temperature and particle. particle : `str` The plasma particle species. z_mean : `~numbers.Real` The average ionization of the plasma. This is only required for calculating the ion permittivity. Returns ------- chi : `~astropy.units.Quantity` The ion or the electron dielectric permittivity of the plasma. This is a dimensionless quantity. Notes ----- The dielectric permittivities for a Maxwellian plasma are described by the following equations (see p. 106 of :cite:t:`froula:2011`): .. math:: χ_e(k, ω) = - \frac{α_e^2}{2} Z'(x_e) χ_i(k, ω) = - \frac{α_i^2}{2}\frac{Z}{} Z'(x_i) α = \frac{ω_p}{k v_{Th}} x = \frac{ω}{k v_{Th}} :math:`χ_e` and :math:`χ_i` are the electron and ion permittivities, respectively. :math:`Z'` is the derivative of the plasma dispersion function. :math:`α` is the scattering parameter which delineates the difference between the collective and non-collective Thomson scattering regimes. :math:`x` is the dimensionless phase velocity of the electromagnetic wave propagating through the plasma. Examples -------- >>> from astropy import units as u >>> from numpy import pi >>> from plasmapy.formulary import thermal_speed >>> T = 30 * 11600 * u.K >>> n = 1e18 * u.cm**-3 >>> particle = 'Ne' >>> z_mean = 8 * u.dimensionless_unscaled >>> vth = thermal_speed(T, particle, method="most_probable") >>> omega = 5.635e14 * 2 * pi * u.rad / u.s >>> k_wave = omega / vth >>> permittivity_1D_Maxwellian(omega, k_wave, T, n, particle, z_mean) <Quantity -6.72809...e-08+5.76037...e-07j> For user convenience `~plasmapy.formulary.dielectric.permittivity_1D_Maxwellian_lite` is bound to this function and can be used as follows: >>> from plasmapy.formulary import plasma_frequency >>> wp = plasma_frequency(n, particle, z_mean=z_mean) >>> permittivity_1D_Maxwellian.lite( ... omega.value, k_wave.value, vth=vth.value, wp=wp.value ... ) (-6.72809...e-08+5.76037...e-07j) """ vth = thermal_speed(T=T, particle=particle, method="most_probable").value wp = plasma_frequency(n=n, particle=particle, z_mean=z_mean).value chi = permittivity_1D_Maxwellian_lite( omega.value, kWave.value, vth, wp, ) return chi * u.dimensionless_unscaled
def cold_plasma_permittivity_LRP(B: u.T, species, n, omega: u.rad / u.s): r""" Magnetized cold plasma dielectric permittivity tensor elements. Elements (L, R, P) are given in the "rotating" basis, i.e. in the basis :math:`(\mathbf{u}_{+}, \mathbf{u}_{-}, \mathbf{u}_z)`, where the tensor is diagonal and with :math:`B ∥ z`\ . The :math:`\exp(-i ω t)` time-harmonic convention is assumed. Parameters ---------- B : `~astropy.units.Quantity` Magnetic field magnitude in units convertible to tesla. species : `list` of `str` The plasma particle species (e.g.: ``['e', 'D+']`` or ``['e', 'D+', 'He+']``. n : `list` of `~astropy.units.Quantity` `list` of species density in units convertible to per cubic meter. The order of the species densities should follow species. omega : `~astropy.units.Quantity` Electromagnetic wave frequency in rad/s. Returns ------- left : `~astropy.units.Quantity` L ("Left") Left-handed circularly polarization tensor element. right : `~astropy.units.Quantity` R ("Right") Right-handed circularly polarization tensor element. plasma : `~astropy.units.Quantity` P ("Plasma") dielectric tensor element. Notes ----- In the rotating frame defined by :math:`(\mathbf{u}_{+}, \mathbf{u}_{-}, \mathbf{u}_z)` with :math:`\mathbf{u}_{\pm}=(\mathbf{u}_x \pm \mathbf{u}_y)/\sqrt{2}`, the dielectric tensor takes a diagonal form with elements L, R, P with: .. math:: L = 1 - \sum_s \frac{ω_{p,s}^2}{ω\left(ω - Ω_{c,s}\right)} R = 1 - \sum_s \frac{ω_{p,s}^2}{ω\left(ω + Ω_{c,s}\right)} P = 1 - \sum_s \frac{ω_{p,s}^2}{ω^2} where :math:`ω_{p,s}` is the plasma frequency and :math:`Ω_{c,s}` is the signed version of the cyclotron frequency for the species :math:`s` :cite:p:`stix:1992`. Examples -------- >>> from astropy import units as u >>> from numpy import pi >>> B = 2*u.T >>> species = ['e', 'D+'] >>> n = [1e18*u.m**-3, 1e18*u.m**-3] >>> omega = 3.7e9*(2*pi)*(u.rad/u.s) >>> L, R, P = permittivity = cold_plasma_permittivity_LRP(B, species, n, omega) >>> L <Quantity 0.63333...> >>> permittivity.left # namedtuple-style access <Quantity 0.63333...> >>> R <Quantity 1.41512...> >>> P <Quantity -4.8903...> """ L, R, P = 1, 1, 1 for s, n_s in zip(species, n): omega_c = gyrofrequency(B=B, particle=s, signed=True) omega_p = plasma_frequency(n=n_s, particle=s) L += -(omega_p**2) / (omega * (omega - omega_c)) R += -(omega_p**2) / (omega * (omega + omega_c)) P += -(omega_p**2) / omega**2 return RotatingTensorElements(L, R, P)
def spectral_density( wavelengths: u.nm, probe_wavelength: u.nm, n: u.m**-3, Te: u.K, Ti: u.K, efract: np.ndarray = None, ifract: np.ndarray = None, ion_species: Union[str, List[str], Particle, List[Particle]] = "H+", electron_vel: u.m / u.s = None, ion_vel: u.m / u.s = None, probe_vec=np.array([1, 0, 0]), scatter_vec=np.array([0, 1, 0]), ) -> Tuple[Union[np.floating, np.ndarray], np.ndarray]: r""" Calculate the spectral density function for Thomson scattering of a probe laser beam by a multi-species Maxwellian plasma. This function calculates the spectral density function for Thomson scattering of a probe laser beam by a plasma consisting of one or more ion species and a one or more thermal electron populations (the entire plasma is assumed to be quasi-neutral) .. math:: S(k,ω) = \sum_e \frac{2π}{k} \bigg |1 - \frac{χ_e}{ε} \bigg |^2 f_{e0,e} \bigg (\frac{ω}{k} \bigg ) + \sum_i \frac{2π Z_i}{k} \bigg |\frac{χ_e}{ε} \bigg |^2 f_{i0,i} \bigg ( \frac{ω}{k} \bigg ) where :math:`χ_e` is the electron component susceptibility of the plasma and :math:`ε = 1 + \sum_e χ_e + \sum_i χ_i` is the total plasma dielectric function (with :math:`χ_i` being the ion component of the susceptibility), :math:`Z_i` is the charge of each ion, :math:`k` is the scattering wavenumber, :math:`ω` is the scattering frequency, and :math:`f_{e0,e}` and :math:`f_{i0,i}` are the electron and ion velocity distribution functions respectively. In this function the electron and ion velocity distribution functions are assumed to be Maxwellian, making this function equivalent to Eq. 3.4.6 in :cite:p:`sheffield:2011`\ . Parameters ---------- wavelengths : `~astropy.units.Quantity` Array of wavelengths over which the spectral density function will be calculated. (convertible to nm) probe_wavelength : `~astropy.units.Quantity` Wavelength of the probe laser. (convertible to nm) n : `~astropy.units.Quantity` Mean (0th order) density of all plasma components combined. (convertible to cm\ :sup:`-3`\ .) Te : `~astropy.units.Quantity`, shape (Ne, ) Temperature of each electron component. Shape (Ne, ) must be equal to the number of electron components Ne. (in K or convertible to eV) Ti : `~astropy.units.Quantity`, shape (Ni, ) Temperature of each ion component. Shape (Ni, ) must be equal to the number of ion components Ni. (in K or convertible to eV) efract : array_like, shape (Ne, ), optional An array-like object where each element represents the fraction (or ratio) of the electron component number density to the total electron number density. Must sum to 1.0. Default is a single electron component. ifract : array_like, shape (Ni, ), optional An array-like object where each element represents the fraction (or ratio) of the ion component number density to the total ion number density. Must sum to 1.0. Default is a single ion species. ion_species : `str` or `~plasmapy.particles.particle_class.Particle`, shape (Ni, ), optional A list or single instance of `~plasmapy.particles.Particle`, or strings convertible to `~plasmapy.particles.particle_class.Particle`. Default is ``'H+'`` corresponding to a single species of hydrogen ions. electron_vel : `~astropy.units.Quantity`, shape (Ne, 3), optional Velocity of each electron component in the rest frame. (convertible to m/s) Defaults to a stationary plasma [0, 0, 0] m/s. ion_vel : `~astropy.units.Quantity`, shape (Ni, 3), optional Velocity vectors for each electron population in the rest frame (convertible to m/s) Defaults zero drift for all specified ion species. probe_vec : float `~numpy.ndarray`, shape (3, ) Unit vector in the direction of the probe laser. Defaults to [1, 0, 0]. scatter_vec : float `~numpy.ndarray`, shape (3, ) Unit vector pointing from the scattering volume to the detector. Defaults to [0, 1, 0] which, along with the default `probe_vec`, corresponds to a 90° scattering angle geometry. Returns ------- alpha : `float` Mean scattering parameter, where ``alpha`` > 1 corresponds to collective scattering and ``alpha`` < 1 indicates non-collective scattering. The scattering parameter is calculated based on the total plasma density ``n``. Skw : `~astropy.units.Quantity` Computed spectral density function over the input ``wavelengths`` array with units of s/rad. Notes ----- For details, see :cite:t:`sheffield:2011`\ . This code is a modified version of the program described therein. For a concise summary of the relevant physics, see Chapter 5 of :cite:t:`schaeffer:2014`\ . """ if efract is None: efract = np.ones(1) else: efract = np.asarray(efract, dtype=np.float64) if ifract is None: ifract = np.ones(1) else: ifract = np.asarray(ifract, dtype=np.float64) # If electron velocity is not specified, create an array corresponding # to zero drift if electron_vel is None: electron_vel = np.zeros([efract.size, 3]) * u.m / u.s # If ion drift velocity is not specified, create an array corresponding # to zero drift if ion_vel is None: ion_vel = np.zeros([ifract.size, 3]) * u.m / u.s # Condition ion_species if isinstance(ion_species, (str, Particle)): ion_species = [ion_species] if len(ion_species) == 0: raise ValueError("At least one ion species needs to be defined.") for ii, ion in enumerate(ion_species): if isinstance(ion, Particle): continue ion_species[ii] = Particle(ion) # Condition Te if Te.size == 1: # If a single quantity is given, put it in an array so it's iterable # If Te.size != len(efract), assume same temp. for all species Te = np.repeat(Te, len(efract)) elif Te.size != len(efract): raise ValueError( f"Got {Te.size} electron temperatures and expected {len(efract)}.") # Condition Ti if Ti.size == 1: # If a single quantity is given, put it in an array so it's iterable # If Ti.size != len(ion_species), assume same temp. for all species Ti = [Ti.value] * len(ion_species) * Ti.unit elif Ti.size != len(ion_species): raise ValueError( f"Got {Ti.size} ion temperatures and expected {len(ion_species)}.") # Make sure the sizes of ion_species, ifract, ion_vel, and Ti all match if ((len(ion_species) != ifract.size) or (ion_vel.shape[0] != ifract.size) or (Ti.size != ifract.size)): raise ValueError( f"Inconsistent number of species in ifract ({ifract}), " f"ion_species ({len(ion_species)}), Ti ({Ti.size}), " f"and/or ion_vel ({ion_vel.shape[0]}).") # Make sure the sizes of efract, electron_vel, and Te all match if (electron_vel.shape[0] != efract.size) or (Te.size != efract.size): raise ValueError( f"Inconsistent number of electron populations in efract ({efract.size}), " f"Te ({Te.size}), or electron velocity ({electron_vel.shape[0]}).") # Ensure unit vectors are normalized probe_vec = probe_vec / np.linalg.norm(probe_vec) scatter_vec = scatter_vec / np.linalg.norm(scatter_vec) # Define some constants C = const.c.si # speed of light # Calculate plasma parameters vTe = thermal_speed(Te, particle="e-") vTi, ion_z = [], [] for T, ion in zip(Ti, ion_species): vTi.append(thermal_speed(T, particle=ion).value) ion_z.append(ion.charge_number * u.dimensionless_unscaled) vTi = vTi * vTe.unit zbar = np.sum(ifract * ion_z) ne = efract * n ni = ifract * n / zbar # ne/zbar = sum(ni) # wpe is calculated for the entire plasma (all electron populations combined) wpe = plasma_frequency(n=n, particle="e-") # Convert wavelengths to angular frequencies (electromagnetic waves, so # phase speed is c) ws = (2 * np.pi * u.rad * C / wavelengths).to(u.rad / u.s) wl = (2 * np.pi * u.rad * C / probe_wavelength).to(u.rad / u.s) # Compute the frequency shift (required by energy conservation) w = ws - wl # Compute the wavenumbers in the plasma # See Sheffield Sec. 1.8.1 and Eqs. 5.4.1 and 5.4.2 ks = np.sqrt(ws**2 - wpe**2) / C kl = np.sqrt(wl**2 - wpe**2) / C # Compute the wavenumber shift (required by momentum conservation) scattering_angle = np.arccos(np.dot(probe_vec, scatter_vec)) # Eq. 1.7.10 in Sheffield k = np.sqrt(ks**2 + kl**2 - 2 * ks * kl * np.cos(scattering_angle)) # Normalized vector along k k_vec = (scatter_vec - probe_vec) * u.dimensionless_unscaled k_vec = k_vec / np.linalg.norm(k_vec) # Compute Doppler-shifted frequencies for both the ions and electrons # Matmul is simultaneously conducting dot product over all wavelengths # and ion components w_e = w - np.matmul(electron_vel, np.outer(k, k_vec).T) w_i = w - np.matmul(ion_vel, np.outer(k, k_vec).T) # Compute the scattering parameter alpha # expressed here using the fact that v_th/w_p = root(2) * Debye length alpha = np.sqrt(2) * wpe / np.outer(k, vTe) # Calculate the normalized phase velocities (Sec. 3.4.2 in Sheffield) xe = (np.outer(1 / vTe, 1 / k) * w_e).to(u.dimensionless_unscaled) xi = (np.outer(1 / vTi, 1 / k) * w_i).to(u.dimensionless_unscaled) # Calculate the susceptibilities chiE = np.zeros([efract.size, w.size], dtype=np.complex128) for i in range(len(efract)): chiE[i, :] = permittivity_1D_Maxwellian(w_e[i, :], k, Te[i], ne[i], "e-") # Treatment of multiple species is an extension of the discussion in # Sheffield Sec. 5.1 chiI = np.zeros([ifract.size, w.size], dtype=np.complex128) for i, ion in enumerate(ion_species): chiI[i, :] = permittivity_1D_Maxwellian(w_i[i, :], k, Ti[i], ni[i], ion, z_mean=ion_z[i]) # Calculate the longitudinal dielectric function epsilon = 1 + np.sum(chiE, axis=0) + np.sum(chiI, axis=0) econtr = np.zeros([efract.size, w.size], dtype=np.complex128) * u.s / u.rad for m in range(efract.size): econtr[m, :] = efract[m] * (2 * np.sqrt(np.pi) / k / vTe[m] * np.power( np.abs(1 - np.sum(chiE, axis=0) / epsilon), 2) * np.exp(-xe[m, :]**2)) icontr = np.zeros([ifract.size, w.size], dtype=np.complex128) * u.s / u.rad for m in range(ifract.size): icontr[m, :] = ifract[m] * ( 2 * np.sqrt(np.pi) * ion_z[m] / k / vTi[m] * np.power(np.abs(np.sum(chiE, axis=0) / epsilon), 2) * np.exp(-xi[m, :]**2)) # Recast as real: imaginary part is already zero Skw = np.real(np.sum(econtr, axis=0) + np.sum(icontr, axis=0)) return np.mean(alpha), Skw
def test_nans(self): assert np.isnan(plasma_frequency(np.nan * u.m**-3, "e-"))
class TestPlasmaFrequency: """ Test class for `plasmapy.formulary.frequencies.plasma_frequency`. Note: Testing of `plasma_frequency_lite` is done in a separate test class. """ @pytest.mark.parametrize( "bound_name, bound_attr", [("lite", plasma_frequency_lite)], ) def test_lite_function_binding(self, bound_name, bound_attr): """Test expected attributes are bound correctly.""" assert hasattr(plasma_frequency, bound_name) assert getattr(plasma_frequency, bound_name) is bound_attr def test_lite_function_marking(self): """ Test plasma_frequency is marked as having a Lite-Function. """ assert hasattr(plasma_frequency, "__bound_lite_func__") assert isinstance(plasma_frequency.__bound_lite_func__, dict) for bound_name, bound_origin in plasma_frequency.__bound_lite_func__.items(): assert hasattr(plasma_frequency, bound_name) attr = getattr(plasma_frequency, bound_name) origin = f"{attr.__module__}.{attr.__name__}" assert origin == bound_origin @pytest.mark.parametrize( "args, kwargs, _error", [ ((u.m**-3, "e-"), {}, TypeError), (("not a density", "e-"), {}, TypeError), ((5 * u.s, "e-"), {}, u.UnitTypeError), ((5 * u.m**-2, "e-"), {}, u.UnitTypeError), ((), {"n": 5 * u.m**-3, "particle": "not a particle"}, ValueError), ], ) def test_raises(self, args, kwargs, _error): """ Test scenarios that cause plasma_frequency to raise an Exception. """ with pytest.raises(_error): plasma_frequency(*args, **kwargs) @pytest.mark.parametrize( "args, kwargs, _warning, expected", [ ( (1e19, "e-"), {}, u.UnitsWarning, plasma_frequency(1e19 * u.m**-3, "e-"), ), ((1e19, "p"), {}, u.UnitsWarning, plasma_frequency(1e19 * u.m**-3, "p")), ], ) def test_warns(self, args, kwargs, _warning, expected): """ Test scenarios the cause plasma_frequency to issue a warning. """ with pytest.warns(_warning): wp = plasma_frequency(*args, **kwargs) assert isinstance(wp, u.Quantity) assert wp.unit == u.rad / u.s if expected is not None: assert np.allclose(wp, expected) @pytest.mark.parametrize( "args, kwargs, expected, rtol", [ ((1 * u.cm**-3, "e-"), {}, 5.64e4, 1e-2), ((1 * u.cm**-3, "N"), {}, 3.53e2, 1e-1), ((1e17 * u.cm**-3, "p"), {"z_mean": 0.8}, 333063562455.4028, 1e-6), ( (5e19 * u.m**-3, "p"), {}, plasma_frequency(5e19 * u.m**-3, particle="H-1+").value, 1e-5, ), ((m_p.to(u.u).value * u.cm**-3,), {"particle": "p"}, 1.32e3, 1e-2), ], ) def test_values(self, args, kwargs, expected, rtol): """Test various expected values.""" wp = plasma_frequency(*args, **kwargs) assert isinstance(wp, u.Quantity) assert wp.unit == u.rad / u.s assert np.allclose(wp.value, expected, rtol=rtol) @pytest.mark.parametrize( "args, kwargs", [((1 * u.cm**-3, "N"), {}), ((1e12 * u.cm**-3,), {"particle": "p"})], ) def test_to_hz(self, args, kwargs): """Test behavior of the ``to_hz`` keyword.""" wp = plasma_frequency(*args, **kwargs) fp = plasma_frequency(*args, to_hz=True, **kwargs) assert isinstance(fp, u.Quantity) assert fp.unit == u.Hz assert fp.value == wp.value / (2.0 * np.pi) def test_nans(self): assert np.isnan(plasma_frequency(np.nan * u.m**-3, "e-")) def test_can_handle_numpy_arrays(self): assert_can_handle_nparray(plasma_frequency)