Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 2 additions & 0 deletions changelog/1664.feature.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@
Add the `~plasmapy.formulary.densities` submodule
Add the `~plasmapy.formulary.densities.critical_density` function to the `~plasmapy.formulary.densities` module. This function calculates the critical density of a plasma for a given frequency of radiation.
1 change: 1 addition & 0 deletions docs/about/credits.rst
Original file line number Diff line number Diff line change
Expand Up @@ -68,6 +68,7 @@ in parentheses are `ORCID author identifiers <https://orcid.org>`__.
* Jakub Polak
* :user:`James Kent <jdkent>`
* :user:`Jasper Beckers <jasperbeckers>`
* :user:`Jayden Roberts <JaydenR2305>`
* :user:`Joao Victor Martinelli <JvPy>`
* :user:`Joshua Munn <jams2>`
* :user:`Isaias McHardy <jota33>` (:orcid:`0000-0001-5394-9445`)
Expand Down
9 changes: 9 additions & 0 deletions docs/formulary/densities.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,9 @@
.. _densities:

******************************************************************
Density Plasma Parameters (`plasmapy.formulary.densities`)
******************************************************************

.. currentmodule:: plasmapy.formulary.densities

.. automodapi:: plasmapy.formulary.densities
3 changes: 3 additions & 0 deletions docs/formulary/index.rst
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,9 @@ physical quantities helpful for plasma physics.
| .. toctree:: Collisions <collisions> | `plasmapy.formulary.collisions` |
| :maxdepth: 1 | |
+--------------------------------------------------------+-----------------------------------------+
| .. toctree:: Densities <densities> | `plasmapy.formulary.densities` |
| :maxdepth: 1 | |
+--------------------------------------------------------+-----------------------------------------+
| .. toctree:: Dielectrics <dielectric> | `plasmapy.formulary.dielectric` |
| :maxdepth: 1 | |
+--------------------------------------------------------+-----------------------------------------+
Expand Down
2 changes: 2 additions & 0 deletions plasmapy/formulary/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@

from plasmapy.formulary.braginskii import *
from plasmapy.formulary.collisions import *
from plasmapy.formulary.densities import *
from plasmapy.formulary.dielectric import *
from plasmapy.formulary.dimensionless import *
from plasmapy.formulary.distribution import *
Expand Down Expand Up @@ -37,6 +38,7 @@
for modname in (
"braginskii",
"collisions",
"densities",
"dielectric",
"dimensionless",
"distribution",
Expand Down
57 changes: 57 additions & 0 deletions plasmapy/formulary/densities.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,57 @@
"""Functions to calculate plasma density parameters."""
__all__ = [
"critical_density",
]

import astropy.units as u

from astropy.constants.si import e, eps0, m_e

from plasmapy.utils.decorators import validate_quantities


@validate_quantities(
omega={"can_be_negative": False},
validations_on_return={
"units": [u.m**-3],
},
)
def critical_density(omega: u.rad / u.s) -> u.m**-3:
r"""Calculate the plasma critical density for a radiation of a given frequency.

Parameters
----------
omega: `~astropy.units.Quantity`
The radiation frequency in units of angular frequency.

Returns
-------
n_c : `~astropy.units.Quantity`
The plasma critical density.

Notes
-----
The critical density for a given frequency of radiation is
defined as the value at which the electron plasma frequency equals
the frequency of the radiation.

The critical density is given by the formula

.. math::
n_{c}=\frac{m_{e}\varepsilon_0\omega^{2}}{e^{2}}

where :math:`m_{e}` is the mass of an electron,
:math:`\varepsilon_0` is the permittivity of free space, :math:`\omega`
is the radiation frequency, and :math:`e` is the elementary charge.

Examples
--------
>>> from astropy import units as u
>>> critical_density(5e15 * u.rad/u.s)
<Quantity 7.85519457e+27 1 / m3>

"""

n_c = m_e * eps0 * omega**2 / (e**2)

return n_c.to(u.m**-3, equivalencies=u.dimensionless_angles())
35 changes: 35 additions & 0 deletions plasmapy/formulary/tests/test_densities.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,35 @@
"""Tests for functionality contained in `plasmapy.formulary.frequencies`."""
import astropy.units as u
import numpy as np
import pytest

from plasmapy.formulary.densities import critical_density
from plasmapy.formulary.frequencies import plasma_frequency


class TestCriticalDensity:
"""Test the plasma_critical_density function in densities.py."""
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

missed a spot


n_i = 5e19 * u.m**-3
omega = plasma_frequency(n=n_i, particle="e-")

@pytest.fixture()
def n_c(self):
"""Get the critical density for the example frequency"""
return critical_density(self.omega)

def test_units(self, n_c):
"""Test the return units"""

assert n_c.unit.is_equivalent(u.m**-3)

def test_value(self, n_c):
"""
Compare the calculated critical density with the expected value.

The plasma critical density for a given frequency of radiation is
defined as the value at which the electron plasma frequency equals
the frequency of the radiation.
"""

assert np.isclose(n_c.value, self.n_i.value)