{ "cells": [ { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "Analysing ITER parameters\n", "=========================\n", "\n", "Let's try to look at ITER plasma conditions using the `physics` subpackage.\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "import matplotlib.pyplot as plt\n", "import numpy as np\n", "\n", "from astropy import units as u\n", "from mpl_toolkits.mplot3d import Axes3D\n", "\n", "from plasmapy import formulary" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The radius of electric field shielding clouds, also known as the :func:`~plasmapy.formulary.parameters.Debye_length`,\n", "would be\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "electron_temperature = 8.8 * u.keV\n", "electron_concentration = 10.1e19 / u.m ** 3\n", "print(formulary.Debye_length(electron_temperature, electron_concentration))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that we can also neglect the unit for the concentration, as\n", "1/m^3 is the a standard unit for this kind of Quantity:\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "print(formulary.Debye_length(electron_temperature, 10.1e19))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Assuming the magnetic field as 5.3 Teslas (which is the value at the major\n", "radius):\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "B = 5.3 * u.T\n", "\n", "print(formulary.gyrofrequency(B, particle=\"e\"))\n", "\n", "print(formulary.gyroradius(B, T_i=electron_temperature, particle=\"e\"))" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The electron :func:`~plasmapy.formulary.parameters.inertial_length` would be\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "print(formulary.inertial_length(electron_concentration, particle=\"e\"))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In these conditions, they should reach thermal velocities of about\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "print(formulary.thermal_speed(T=electron_temperature, particle=\"e\"))" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "And the Langmuir wave :func:`~plasmapy.formulary.parameters.plasma_frequency` should be on the order of\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "print(formulary.plasma_frequency(electron_concentration, particle=\"e-\"))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let's try to recreate some plots and get a feel for some of these quantities.\n", "We will also compare our data to real-world plasma situations.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false }, "tags": [ "nbsphinx-thumbnail" ] }, "outputs": [], "source": [ "n_e = np.logspace(4, 30, 100) / u.m ** 3\n", "plt.plot(n_e, formulary.plasma_frequency(n_e, particle=\"e-\"))\n", "plt.scatter(\n", " electron_concentration,\n", " formulary.plasma_frequency(electron_concentration, particle=\"e-\"),\n", " label=\"Our Data\",\n", ")\n", "\n", "# IRT1 Tokamak Data\n", "# http://article.sapub.org/pdf/10.5923.j.jnpp.20110101.03.pdf\n", "n_e = 1.2e19 / u.m ** 3\n", "T_e = 136.8323 * u.eV\n", "B = 0.82 * u.T\n", "plt.scatter(n_e, formulary.plasma_frequency(n_e, particle=\"e-\"), label=\"IRT1 Tokamak\")\n", "\n", "# Wendelstein 7-X Stellerator Data\n", "# https://nucleus.iaea.org/sites/fusionportal/Shared%20Documents/FEC%202016/fec2016-preprints/preprint0541.pdf\n", "n_e = 3e19 / u.m ** 3\n", "T_e = 6 * u.keV\n", "B = 3 * u.T\n", "plt.scatter(\n", " n_e, formulary.plasma_frequency(n_e, particle=\"e-\"), label=\"W7-X Stellerator\"\n", ")\n", "\n", "# Solar Corona\n", "# Estimated by Nick Murphy\n", "n_e = 1e15 / u.m ** 3\n", "T_e = 1e6 * u.K\n", "B = 0.005 * u.T\n", "T_e.to(u.eV, equivalencies=u.temperature_energy())\n", "plt.scatter(n_e, formulary.plasma_frequency(n_e, particle=\"e-\"), label=\"Solar Corona\")\n", "\n", "# Interstellar (warm neutral) Medium\n", "# Estimated by Nick Murphy\n", "n_e = 1e6 / u.m ** 3\n", "T_e = 5e3 * u.K\n", "B = 0.005 * u.T\n", "T_e.to(u.eV, equivalencies=u.temperature_energy())\n", "plt.scatter(\n", " n_e, formulary.plasma_frequency(n_e, particle=\"e-\"), label=\"Interstellar Medium\"\n", ")\n", "\n", "# Solar Wind at 1 AU\n", "# Estimated by Nick Murphy\n", "n_e = 7e6 / u.m ** 3\n", "T_e = 1e5 * u.K\n", "B = 5e-9 * u.T\n", "T_e.to(u.eV, equivalencies=u.temperature_energy())\n", "plt.scatter(\n", " n_e, formulary.plasma_frequency(n_e, particle=\"e-\"), label=\"Solar Wind (1AU)\"\n", ")\n", "\n", "\n", "plt.xlabel(\"Electron Concentration (m^-3)\")\n", "plt.ylabel(\"Langmuir Wave Plasma Frequency (rad/s)\")\n", "plt.grid()\n", "plt.xscale(\"log\")\n", "plt.yscale(\"log\")\n", "plt.legend()\n", "plt.title(\"Log-scale plot of plasma frequencies\")\n", "plt.show()" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "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", "version": "3.8.2" }, "widgets": { "application/vnd.jupyter.widget-state+json": { "state": {}, "version_major": 2, "version_minor": 0 } } }, "nbformat": 4, "nbformat_minor": 4 }