{ "cells": [ { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "# Cold Magnetized Plasma Dielectric Permittivity Tensor" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "This notebook shows how to calculate the values of the cold plasma tensor elements for various electromagnetic wave frequencies." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "import matplotlib.pyplot as plt\n", "import numpy as np\n", "\n", "from astropy import units as u\n", "from astropy.visualization import quantity_support\n", "\n", "from plasmapy.formulary import (\n", " cold_plasma_permittivity_LRP,\n", " cold_plasma_permittivity_SDP,\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "[astropy.units]: https://docs.astropy.org/en/stable/units/index.html\n", "[cold_plasma_permittivity_SDP]: https://docs.plasmapy.org/en/stable/api/plasmapy.formulary.dielectric.cold_plasma_permittivity_SDP.html\n", "[cold_plasma_permittivity_LRP]: https://docs.plasmapy.org/en/stable/api/plasmapy.formulary.dielectric.cold_plasma_permittivity_LRP.html\n", "\n", "For more information, check out the documentation on [cold_plasma_permittivity_SDP] and [cold_plasma_permittivity_LRP].\n", "\n", "Let's use [astropy.units] to define some parameters, such as the magnetic field magnitude, plasma species, and densities." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "B = 2 * u.T\n", "species = [\"e\", \"D+\"]\n", "n = [1e18 * u.m**-3, 1e18 * u.m**-3]" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "Let's pick some frequencies in the radio frequency range." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "f_min = 1 * u.MHz\n", "f_max = 200 * u.GHz\n", "\n", "f = np.geomspace(f_min, f_max, num=5000).to(u.Hz)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "[equivalency]: https://docs.astropy.org/en/stable/units/equivalencies.html\n", "[angular frequencies]: https://docs.plasmapy.org/en/stable/development/code_guide.html#angular-frequencies\n", "\n", "Next we convert these frequencies to angular frequencies. To do this we must specify an [equivalency] between a cycle per second and a hertz. The reasons why are described in the section of PlasmaPy's documentation on [angular frequencies]." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "ω_RF = f.to(u.rad / u.s, equivalencies=[(u.Hz, u.cycle / u.s)])" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "[Stix 1992]: https://link.springer.com/book/9780883188590\n", "\n", "In a Jupyter notebook, we can type `\\omega` and press tab to get \"ω\".\n", "\n", "Now we are ready to calculate the $S$ (sum), $D$ (difference), and $P$ (plasma) components of the dielectric tensor in the \"Stix\" frame with $B ∥ ẑ$. This notation is from [Stix 1992]." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "S, D, P = cold_plasma_permittivity_SDP(B=B, species=species, n=n, omega=ω_RF)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "Next we will filter the negative and positive values so that they can be plotted separately. We'll be doing this multiple times, so let's create a function using the `numpy.ma` module for working with masked arrays." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "def filter_negative_and_positive_values(arr):\n", " \"\"\"\n", " Return an array with only negative values of ``arr``\n", " and an array with only positive values of ``arr``.\n", " Each element of these arrays that does not meet the\n", " condition are replaced with `~numpy.nan`.\n", " \"\"\"\n", " arr_neg = np.ma.masked_greater_equal(arr, 0).filled(np.nan)\n", " arr_pos = np.ma.masked_less_equal(arr, 0).filled(np.nan)\n", " return arr_neg, arr_pos" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "Let's apply this function to `S`, `D`, and `P`." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "S_neg, S_pos = filter_negative_and_positive_values(S)\n", "D_neg, D_pos = filter_negative_and_positive_values(D)\n", "P_neg, P_pos = filter_negative_and_positive_values(P)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "[Quantity]: https://docs.astropy.org/en/stable/units/quantity.html\n", "[quantity_support]: https://docs.astropy.org/en/stable/api/astropy.visualization.quantity_support.html?highlight=quantity_support#quantity-support\n", "\n", "Because we are using [Quantity] objects, we will need to call [quantity_support] before plotting the results." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "quantity_support()" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "While we're at it, let's specify a few colorblind friendly colors to use in the plots." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "red, blue, orange = \"#920000\", \"#006ddb\", \"#db6d00\"" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "Let's plot the elements of the cold plasma dielectric tensor in the Stix frame. " ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "ylim = (1e-4, 1e8)\n", "\n", "plt.figure(figsize=(12, 6))\n", "\n", "plt.semilogx(f, abs(S_neg), red, lw=2, ls=\"--\", label=\"S < 0\")\n", "plt.semilogx(f, abs(D_neg), blue, lw=2, ls=\"--\", label=\"D < 0\")\n", "plt.semilogx(f, abs(P_neg), orange, lw=2, ls=\"--\", label=\"P < 0\")\n", "\n", "plt.semilogx(f, S_pos, red, lw=2, label=\"S > 0\")\n", "plt.semilogx(f, D_pos, blue, lw=2, label=\"D > 0\")\n", "plt.semilogx(f, P_pos, orange, lw=2, label=\"P > 0\")\n", "\n", "plt.ylim(ylim)\n", "plt.xlim(f_min, f_max)\n", "\n", "plt.title(\n", " \"Cold plasma dielectric permittivity tensor components in Stix frame\", size=16\n", ")\n", "\n", "plt.xlabel(\"Frequency [Hz]\", size=16)\n", "plt.ylabel(\"Absolute value\", size=16)\n", "\n", "plt.yscale(\"log\")\n", "\n", "plt.grid(True, which=\"major\")\n", "plt.grid(True, which=\"minor\")\n", "\n", "plt.tick_params(labelsize=14)\n", "\n", "plt.legend(fontsize=18, ncol=2, framealpha=1)\n", "\n", "plt.show()" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "[cold_plasma_permittivity_LRP]: https://docs.plasmapy.org/en/stable/api/plasmapy.formulary.dielectric.cold_plasma_permittivity_LRP.html\n", "\n", "Next let's get the dielectric permittivity tensor elements in the \"rotating\" basis where the tensor is diagonal and with $B ∥ ẑ$. We will use [cold_plasma_permittivity_LRP] to get the left-handed circular polarization tensor element $L$, the right-handed circular polarization tensor element $R$, and the plasma component $P$." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "L, R, P = cold_plasma_permittivity_LRP(B, species, n, ω_RF)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "Let's use the function from earlier to prepare for plotting." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "L_neg, L_pos = filter_negative_and_positive_values(L)\n", "R_neg, R_pos = filter_negative_and_positive_values(R)\n", "P_neg, P_pos = filter_negative_and_positive_values(P)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% raw\n" } }, "source": [ "Now we can plot the tensor components rotating frame too." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "plt.figure(figsize=(12, 6))\n", "\n", "plt.semilogx(f, abs(L_neg), red, lw=2, ls=\"--\", label=\"L < 0\")\n", "plt.semilogx(f, abs(R_neg), blue, lw=2, ls=\"--\", label=\"R < 0\")\n", "plt.semilogx(f, abs(P_neg), orange, lw=2, ls=\"--\", label=\"P < 0\")\n", "\n", "plt.semilogx(f, L_pos, red, lw=2, label=\"L > 0\")\n", "plt.semilogx(f, R_pos, blue, lw=2, label=\"R > 0\")\n", "plt.semilogx(f, P_pos, orange, lw=2, label=\"P > 0\")\n", "\n", "plt.ylim(ylim)\n", "plt.xlim(f_min, f_max)\n", "\n", "plt.title(\n", " \"Cold plasma dielectric permittivity tensor components in the rotating frame\",\n", " size=16,\n", ")\n", "\n", "plt.xlabel(\"Frequency [Hz]\", size=16)\n", "plt.ylabel(\"Absolute value\", size=16)\n", "\n", "plt.yscale(\"log\")\n", "\n", "plt.grid(True, which=\"major\")\n", "plt.grid(True, which=\"minor\")\n", "\n", "plt.tick_params(labelsize=14)\n", "\n", "plt.legend(fontsize=18, ncol=2, framealpha=1)\n", "\n", "plt.show()" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3" }, "widgets": { "application/vnd.jupyter.widget-state+json": { "state": {}, "version_major": 2, "version_minor": 0 } } }, "nbformat": 4, "nbformat_minor": 4 }