def alfven_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_ci = pfp.gyrofrequency(B=B, particle=ion, signed=False, Z=z_mean)

    #Grid/vector creation for k?

    #Parameters kz

    kz = np.cos(theta.value) * k
    kx = np.sqrt(k**2 - kz**2)

    #Parameters sigma, D, and F to simplify equation 3
    A = (kz * v_A)**2
    F = ((kx * c_s) / omega_ci)**2

    omega = np.sqrt(A * (1 + F))
    print(omega_ci)
    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 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
예제 #4
0
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