{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "Langmuir probe data analysis\n", "============================\n", "\n", "Let's analyze a few Langmuir probe characteristics using the\n", "`diagnostics.langmuir` subpackage. First we need to import the module and some\n", "basics.\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "%matplotlib inline\n", "\n", "import astropy.units as u\n", "import matplotlib.pyplot as plt\n", "import numpy as np\n", "import os\n", "\n", "from pprint import pprint\n", "\n", "from plasmapy.diagnostics.langmuir import Characteristic, swept_probe_analysis" ] }, { "cell_type": "raw", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The first characteristic we analyze is a simple single-probe measurement in\n", "a low (ion) temperature, low density plasma with a cylindrical probe. This\n", "allows us to utilize OML theory implemented in :func:`~plasmapy.diagnostics.langmuir.swept_probe_analysis`.\n", "The data has been preprocessed with some smoothing, which allows us to obtain\n", "a Electron Energy Distribution Function (EEDF) as well." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "# Load the bias and current values stored in the .p pickle file.\n", "path = os.path.join(\"langmuir_samples\", \"Beckers2017.npy\")\n", "bias, current = np.load(path)\n", "\n", "# Create the Characteristic object, taking into account the correct units\n", "characteristic = Characteristic(u.Quantity(bias, u.V),\n", " u.Quantity(current, u.A))\n", "\n", "# Calculate the cylindrical probe surface area\n", "probe_length = 1.145 * u.mm\n", "probe_diameter = 1.57 * u.mm\n", "probe_area = (probe_length * np.pi * probe_diameter +\n", " np.pi * 0.25 * probe_diameter**2)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we can actually perform the analysis. Since the plasma is in Helium an\n", "ion mass number of 4 is entered. The results are visualized and the obtained\n", "EEDF is also shown.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false }, "nbsphinx-thumbnail": { "output-index": 2 } }, "outputs": [], "source": [ "pprint(swept_probe_analysis(characteristic,\n", " probe_area, 'He-4+',\n", " visualize=True,\n", " plot_EEDF=True))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The cyan and yellow lines indicate the fitted electron and ion currents,\n", "respectively. The green line is the sum of these and agrees nicely with the\n", "data. This indicates a successful analysis.\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The next sample probe data is provided by David Pace. It is also obtained\n", "from a low relatively ion temperature and density plasma, in Argon.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "# Load the data from a file and create the Characteristic object\n", "path = os.path.join(\"langmuir_samples\", \"Pace2015.npy\")\n", "bias, current = np.load(path)\n", "characteristic = Characteristic(u.Quantity(bias, u.V),\n", " u.Quantity(current, u.A))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Initially the electrons are assumed to be Maxwellian. To check this the fit\n", "of the electron growth region will be plotted.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "swept_probe_analysis(characteristic,\n", " 0.738 * u.cm**2,\n", " 'Ar-40 1+',\n", " bimaxwellian=False,\n", " plot_electron_fit=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "It can be seen that this plasma is slightly bi-Maxwellian, as there are two\n", "distinct slopes in the exponential section. The analysis is now performed\n", "with bimaxwellian set to True, which yields improved results.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "pprint(swept_probe_analysis(characteristic,\n", " 0.738 * u.cm**2,\n", " 'Ar-40 1+',\n", " bimaxwellian=True,\n", " visualize=True,\n", " plot_electron_fit=True))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The probe current resolution of the raw data is relatively poor, but the\n", "analysis still performs well in the ion current region. The bi-Maxwellian\n", "properties are not significant but do make a difference. Check this analysis\n", "without setting `bimaxwellian` to True!\n", "This is reflected in the results, which indicate that the temperatures of\n", "the cold and hot electron population are indeed different, but relatively\n", "close.\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This Helium plasma is fully bi-Maxwellian.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "# Import probe data and calculate probe surface area.\n", "path = os.path.join(\"langmuir_samples\", \"Beckers2017b.npy\")\n", "bias, current = np.load(path)\n", "characteristic = Characteristic(u.Quantity(bias, u.V),\n", " u.Quantity(current, u.A))\n", "probe_length = 1.145 * u.mm\n", "probe_diameter = 1.57 * u.mm\n", "probe_area = (probe_length * np.pi * probe_diameter +\n", " np.pi * 0.25 * probe_diameter**2)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "`plot_electron_fit` is set to True to check the bi-Maxwellian properties.\n", "The fit converges nicely to the two slopes of the electron growth region.\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "pprint(swept_probe_analysis(characteristic,\n", " probe_area,\n", " 'He-4+',\n", " bimaxwellian=True,\n", " plot_electron_fit=True,\n", " visualize=True))" ] } ], "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 }