def test_upper_hybrid_frequency(): r"""Test the upper_hybrid_frequency function in parameters.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_lower_hybrid_frequency(): r"""Test the lower_hybrid_frequency function in parameters.py.""" ion = "He-4 1+" omega_ci = gyrofrequency(B, particle=ion) omega_pi = plasma_frequency(n=n_i, particle=ion) omega_ce = gyrofrequency(B) 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) == lower_hybrid_frequency( 1.3 * u.T, 1e19 * u.m ** -3 ) assert_can_handle_nparray(lower_hybrid_frequency)
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,\omega) = \sum_e \frac{2\pi}{k} \bigg |1 - \frac{\chi_e}{\epsilon} \bigg |^2 f_{e0,e} \bigg (\frac{\omega}{k} \bigg ) + \sum_i \frac{2\pi Z_i}{k} \bigg |\frac{\chi_e}{\epsilon} \bigg |^2 f_{i0,i} \bigg ( \frac{\omega}{k} \bigg ) where :math:`\chi_e` is the electron component susceptibility of the plasma and :math:`\epsilon = 1 + \sum_e \chi_e + \sum_i \chi_i` is the total plasma dielectric function (with :math:`\chi_i` being the ion component of the susceptibility), :math:`Z_i` is the charge of each ion, :math:`k` is the scattering wavenumber, :math:`\omega` 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 `Sheffield`_. 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^-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`, shape (Ni, ), optional A list or single instance of `~plasmapy.particles.Particle`, or strings convertible to `~plasmapy.particles.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 degree 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 "Plasma Scattering of Electromagnetic Radiation" by Sheffield et al. `ISBN 978\\-0123748775`_. This code is a modified version of the program described therein. For a concise summary of the relevant physics, see Chapter 5 of Derek Schaeffer's thesis, DOI: `10.5281/zenodo.3766933`_. .. _`ISBN 978\\-0123748775`: https://www.sciencedirect.com/book/9780123748775/plasma-scattering-of-electromagnetic-radiation .. _`10.5281/zenodo.3766933`: https://doi.org/10.5281/zenodo.3766933 .. _`Sheffield`: https://doi.org/10.1016/B978-0-12-374877-5.00003-8 """ 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 electon 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 " f"{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 " f"{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.integer_charge * 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)) # Normal vector along k k_vec = (scatter_vec - probe_vec) * u.dimensionless_unscaled # 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, fract in enumerate(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_plasma_frequency(): r"""Test the plasma_frequency function in parameters.py.""" assert plasma_frequency(n_e, "e-").unit.is_equivalent(u.rad / u.s) assert plasma_frequency(n_e, "e-", to_hz=True).unit.is_equivalent(u.Hz) assert np.isclose(plasma_frequency(1 * u.cm**-3, "e-").value, 5.64e4, rtol=1e-2) assert np.isclose(plasma_frequency(1 * u.cm**-3, particle="N").value, 3.53e2, rtol=1e-1) assert np.isclose( plasma_frequency(1 * u.cm**-3, particle="N", to_hz=True).value, 56.19000195094519, ) with pytest.raises(TypeError): plasma_frequency(u.m**-3, "e-") with pytest.raises(u.UnitTypeError): plasma_frequency(5 * u.m**-2, "e-") assert np.isnan(plasma_frequency(np.nan * u.m**-3, "e-")) with pytest.warns(u.UnitsWarning): assert plasma_frequency(1e19, "e-") == plasma_frequency( 1e19 * u.m**-3, "e-") assert plasma_frequency(n_i, particle="p").unit.is_equivalent(u.rad / u.s) # Case where Z=1 is assumed assert plasma_frequency(n_i, particle="H-1+") == plasma_frequency(n_i, particle="p") assert np.isclose(plasma_frequency(mu * u.cm**-3, particle="p").value, 1.32e3, rtol=1e-2) with pytest.raises(ValueError): plasma_frequency(n=5 * u.m**-3, particle="sdfas") with pytest.warns(u.UnitsWarning): plasma_freq_no_units = plasma_frequency(1e19, particle="p") assert plasma_freq_no_units == plasma_frequency(1e19 * u.m**-3, particle="p") plasma_frequency(1e17 * u.cm**-3, particle="p") # testing for user input z_mean testMeth1 = plasma_frequency(1e17 * u.cm**-3, particle="p", z_mean=0.8).si.value testTrue1 = 333063562455.4028 errStr = f"plasma_frequency() gave {testMeth1}, should be {testTrue1}." assert np.isclose(testMeth1, testTrue1, atol=0.0, rtol=1e-6), errStr assert_can_handle_nparray(plasma_frequency)
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, ie. with B // z. The :math:`\exp(-i \omega 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 \omega t)` time-harmonic convention as :math:`\varepsilon = \varepsilon_0 A`, with :math:`A` being .. math:: \varepsilon = \varepsilon_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{\omega_{p,s}^2}{\omega^2 - \Omega_{c,s}^2} D = \sum_s \frac{\Omega_{c,s}}{\omega} \frac{\omega_{p,s}^2}{\omega^2 - \Omega_{c,s}^2} P = 1 - \sum_s \frac{\omega_{p,s}^2}{\omega^2} where :math:`\omega_{p,s}` is the plasma frequency and :math:`\Omega_{c,s}` is the signed version of the cyclotron frequency for the species :math:`s`. References ---------- - T.H. Stix, Waves in Plasma, 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) >>> 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 = parameters.gyrofrequency(B=B, particle=s, signed=True) omega_p = parameters.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""" 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. This is often modulated by the dispersion of the plasma or by relativistic effects. See em_wave.py for ways to calculate this. 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 : str 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 [1]_ .. math:: \chi_e(k, \omega) = - \frac{\alpha_e^2}{2} Z'(x_e) \chi_i(k, \omega) = - \frac{\alpha_i^2}{2}\frac{Z}{} Z'(x_i) \alpha = \frac{\omega_p}{k v_{Th}} x = \frac{\omega}{k v_{Th}} :math:`chi_e` and :math:`chi_i` are the electron and ion permittivities respectively. :math:`Z'` is the derivative of the plasma dispersion function. :math:`\alpha` 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 EM wave propagating through the plasma. References ---------- .. [1] J. Sheffield, D. Froula, S. H. Glenzer, and N. C. Luhmann Jr, Plasma scattering of electromagnetic radiation: theory and measurement techniques. Chapter 5 Pg 106 (Academic press, 2010). Example ------- >>> from astropy import units as u >>> from numpy import pi >>> from astropy.constants import c >>> T = 30 * 11600 * u.K >>> n = 1e18 * u.cm**-3 >>> particle = 'Ne' >>> z_mean = 8 * u.dimensionless_unscaled >>> vTh = parameters.thermal_speed(T, particle, method="most_probable") >>> omega = 5.635e14 * 2 * pi * u.rad / u.s >>> kWave = omega / vTh >>> permittivity_1D_Maxwellian(omega, kWave, T, n, particle, z_mean) <Quantity -6.72809...e-08+5.76037...e-07j> """ # thermal velocity vTh = parameters.thermal_speed(T=T, particle=particle, method="most_probable") # plasma frequency wp = parameters.plasma_frequency(n=n, particle=particle, z_mean=z_mean) # scattering parameter alpha. # explicitly removing factor of sqrt(2) to be consistent with Froula alpha = np.sqrt(2) * (wp / (kWave * vTh)).to(u.dimensionless_unscaled) # The dimensionless phase velocity of the propagating EM wave. zeta = (omega / (kWave * vTh)).to(u.dimensionless_unscaled) chi = alpha ** 2 * (-1 / 2) * plasma_dispersion_func_deriv(zeta.value) return chi
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, ie. in the basis :math:`(\mathbf{u}_{+}, \mathbf{u}_{-}, \mathbf{u}_z)`, where the tensor is diagonal and with B // z. The :math:`\exp(-i \omega 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{\omega_{p,s}^2}{\omega\left(\omega - \Omega_{c,s}\right)} R = 1 - \sum_s \frac{\omega_{p,s}^2}{\omega\left(\omega + \Omega_{c,s}\right)} P = 1 - \sum_s \frac{\omega_{p,s}^2}{\omega^2} where :math:`\omega_{p,s}` is the plasma frequency and :math:`\Omega_{c,s}` is the signed version of the cyclotron frequency for the species :math:`s`. References ---------- - T.H. Stix, Waves in Plasma, 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 = parameters.gyrofrequency(B=B, particle=s, signed=True) omega_p = parameters.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 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 : `str` or `~plasmapy.particles.Particle`, optional An instance of `~plasmapy.particles.Particle`, or a string convertible to `~plasmapy.particles.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 "Radiation Processes in Plasmas" by Bekefi. `ISBN 978\\-0471063506`_. .. _`ISBN 978\\-0471063506`: https://ui.adsabs.harvard.edu/abs/1966rpp..book.....B/abstract """ # Default n_i is n_e/Z: if n_i is None: n_i = n_e / ion_species.integer_charge # 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: " "hbar*ω/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.integer_charge 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 time_plasma_frequency(self): plasma_frequency(1e19 * u.m**-3, particle='p', to_hz=True)
def hollweg_dispersion_solution( *, 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.deg, gamma_e: Union[float, int] = 1, gamma_i: Union[float, int] = 3, z_mean: Union[float, int] = None, ): # 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.integer_charge) 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 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 not (k.ndim == 0 or k.ndim == 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() theta = theta.to(u.radian) if not (theta.ndim == 0 or theta.ndim == 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 c_s = pfp.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 = pfp.Alfven_speed(B, n_i, ion=ion, z_mean=z_mean) omega_ci = pfp.gyrofrequency(B=B, particle=ion, signed=False, Z=z_mean) omega_pe = pfp.plasma_frequency(n=n_e, particle="e-") # Parameters kx and kz kz = np.cos(theta.value) * k kx = np.sqrt(k**2 - kz**2) # Bellan2012JGR beta param equation 3 beta = (c_s / v_A)**2 # Parameters D, F, sigma, and alpha to simplify equation 3 D = (c_s / omega_ci)**2 F = (c / omega_pe)**2 sigma = (kz * v_A)**2 alpha = (k * v_A)**2 # Polynomial coefficients: c3*x^3 + c2*x^2 + c1*x + c0 = 0 c3 = (F * kx**2 + 1) / sigma c2 = -((alpha / sigma) * (1 + beta + F * kx**2) + D * kx**2 + 1) c1 = alpha * (1 + 2 * beta + D * kx**2) c0 = -beta * alpha * sigma omega = {} fast_mode = [] alfven_mode = [] acoustic_mode = [] # If a single k value is given if np.isscalar(k.value) == True: w = np.emath.sqrt(np.roots([c3.value, c2.value, c1.value, c0.value])) fast_mode = np.max(w) alfven_mode = np.median(w) acoustic_mode = np.min(w) # If mutliple k values are given else: # a0*x^3 + a1*x^2 + a2*x^3 + a3 = 0 for (a0, a1, a2, a3) in zip(c3, c2, c1, c0): w = np.emath.sqrt( np.roots([a0.value, a1.value, a2.value, a3.value])) fast_mode.append(np.max(w)) alfven_mode.append(np.median(w)) acoustic_mode.append(np.min(w)) omega['fast_mode'] = fast_mode * u.rad / u.s omega['alfven_mode'] = alfven_mode * u.rad / u.s omega['acoustic_mode'] = acoustic_mode * u.rad / u.s return omega
def hirose_dispersion_solution( *, 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.deg, gamma_e: Union[float, int] = 1, gamma_i: Union[float, int] = 3, z_mean: Union[float, int] = None, ): # 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.integer_charge) 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 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 not (k.ndim == 0 or k.ndim == 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() theta = theta.to(u.radian) if not (theta.ndim == 0 or theta.ndim == 1): raise ValueError( f"Argument 'theta' needs to be a single valued or 1D array astropy " f"Quantity, got array of shape {k.shape}." ) n_e = z_mean * n_i c_s = pfp.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 = pfp.Alfven_speed(B, n_i, ion=ion, z_mean=z_mean) omega_pi = pfp.plasma_frequency(n=n_i, particle=ion) #Grid/vector creation for k? #Parameters kz kz = np.cos(theta.value) * k #Parameters sigma, D, and F to simplify equation 3 A = (kz * v_A) ** 2 B = (k * c_s) ** 2 C = (k * v_A) ** 2 D = ((k * c) / omega_pi ) ** 2 #Polynomial coefficients where x in 'cx' represents the order of the term c3 = 1 c2 = A * (1 + D) + B + C c1 = A * (2 * B + C + B * D) c0 = -B * A ** 2 [L1, L2, L3] = np.roots([c3, c2.value, c1.value, c0.value]) [omega1, omega2, omega3] = [np.emath.sqrt(L1), np.emath.sqrt(L2), np.emath.sqrt(L3)] return omega1, omega2, omega3
def two_fluid_dispersion_solution( *, 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.deg, gamma_e: Union[float, int] = 1, gamma_i: Union[float, int] = 3, z_mean: Union[float, int] = None, ): r""" Using the solution provided by Bellan 2012, calculate the analytical solution to the two fluid, low-frequency (:math:`\omega/kc \ll 1`) dispersion relation presented by Stringer 1963. This dispersion relation also assummes 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 :math:`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 :math:`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 :math:`m^{-3}`. T_e : `~astropy.units.Quantity` The electron temperature in units of :math:`K` or :math:`eV`. T_i : `~astropy.units.Quantity` The ion temperature in units of :math:`K` or :math:`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 :math:`deg`. 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 :math:`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`. 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 Springer 1963 [2]_ (equation 1 of Bellan 2012 [1]_) 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 Bellan 2012 [1]_ the exact roots of the above dispersion equation can be derived and expressed as one analytical solution (equation 38 of Bellan 2012 [1]_): .. 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 References ---------- .. [1] PM Bellan, Improved basis set for low frequency plasma waves, 2012, JGR, 117, A12219, doi: `10.1029/2012JA017856 <https://agupubs.onlinelibrary.wiley.com/doi/10.1029/2012JA017856>`_. .. [2] TE Stringer, Low-frequency waves in an unbounded plasma, 1963, JNE, Part C, doi: `10.1088/0368-3281/5/2/304 <https://doi.org/10.1088/0368-3281/5/2/304>`_ Examples -------- >>> from astropy import units as u >>> from plasmapy.dispersion import two_fluid_dispersion >>> 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_dispersion_solution(**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_dispersion_solution(**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 or ion.is_category("element")): raise ValueError( f"The particle passed for 'ion' must be an ion or element.") # validate z_mean if z_mean is None: try: z_mean = abs(ion.integer_charge) 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 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 not (k.ndim == 0 or k.ndim == 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(f"Argument 'k' can not be a or have negative values.") # validate argument theta theta = theta.squeeze() theta = theta.to(u.radian) if not (theta.ndim == 0 or theta.ndim == 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 = pfp.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 = pfp.Alfven_speed(B, n_i, ion=ion, z_mean=z_mean) omega_ci = pfp.gyrofrequency(B=B, particle=ion, signed=False, Z=z_mean) omega_pe = pfp.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