{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "```\n", "This notebook sets up and runs a set of benchmarks to compare\n", "different numerical discretizations of the SWEs\n", "\n", "Copyright (C) 2016 SINTEF ICT\n", "\n", "This program is free software: you can redistribute it and/or modify\n", "it under the terms of the GNU General Public License as published by\n", "the Free Software Foundation, either version 3 of the License, or\n", "(at your option) any later version.\n", "\n", "This program is distributed in the hope that it will be useful,\n", "but WITHOUT ANY WARRANTY; without even the implied warranty of\n", "MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the\n", "GNU General Public License for more details.\n", "\n", "You should have received a copy of the GNU General Public License\n", "along with this program. If not, see .\n", "```" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Import modules and set up environment" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "#Lets have matplotlib \"inline\"\n", "%matplotlib inline\n", "%config InlineBackend.figure_format = 'retina'\n", "\n", "#Import packages we need\n", "import numpy as np\n", "from matplotlib import animation, rc\n", "from matplotlib import pyplot as plt\n", "\n", "import os\n", "import pyopencl\n", "import datetime\n", "import sys\n", "\n", "#Set large figure sizes\n", "rc('figure', figsize=(6.0, 4.0))\n", "rc('animation', html='html5')\n", "\n", "#Import our simulator\n", "from SWESimulators import FBL, CTCS,KP07, CDKLM16, PlotHelper, Common\n", "#Import initial condition and bathymetry generating functions:\n", "from SWESimulators.BathymetryAndICs import *" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "#Make sure we get compiler output from OpenCL\n", "os.environ[\"PYOPENCL_COMPILER_OUTPUT\"] = \"1\"\n", "\n", "#Set which CL device to use, and disable kernel caching\n", "if (str.lower(sys.platform).startswith(\"linux\")):\n", " os.environ[\"PYOPENCL_CTX\"] = \"0\"\n", "else:\n", " os.environ[\"PYOPENCL_CTX\"] = \"1\"\n", "os.environ[\"CUDA_CACHE_DISABLE\"] = \"1\"\n", "os.environ[\"PYOPENCL_COMPILER_OUTPUT\"] = \"1\"\n", "os.environ[\"PYOPENCL_NO_CACHE\"] = \"1\"\n", "\n", "#Create OpenCL context\n", "cl_ctx = pyopencl.create_some_context()\n", "print \"Using \", cl_ctx.devices[0].name" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "#Create output directory for images\n", "imgdir='images_convergence_' + datetime.datetime.now().strftime(\"%Y_%m_%d-%H_%M_%S\")\n", "os.makedirs(imgdir)\n", "print \"Saving images to \" + imgdir" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "def setBwStyles(ax):\n", " from cycler import cycler\n", "\n", " ax.set_prop_cycle( cycler('marker', ['.', 'x', 4, '+', '*', '1']) +\n", " cycler('linestyle', ['-.', '--', ':', '-.', '--', ':']) +\n", " #cycler('markersize', [15, 15, 15, 15, 15, 15]) +\n", " cycler('color', ['k', 'k', 'k', 'k', 'k', 'k']) )" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "def rebin(a, *args):\n", " '''rebin ndarray data into a smaller ndarray of the same rank whose dimensions\n", " are factors of the original dimensions. eg. An array with 6 columns and 4 rows\n", " can be reduced to have 6,3,2 or 1 columns and 4,2 or 1 rows.\n", " example usages:\n", " >>> a=rand(6,4); b=rebin(a,3,2)\n", " >>> a=rand(6); b=rebin(a,2)\n", " '''\n", " shape = a.shape\n", " lenShape = len(shape)\n", " factor = np.asarray(shape)/np.asarray(args)\n", " evList = ['a.reshape('] + \\\n", " ['args[%d],factor[%d],'%(i,i) for i in range(lenShape)] + \\\n", " [')'] + ['.sum(%d)'%(i+1) for i in range(lenShape)] + \\\n", " ['/factor[%d]'%i for i in range(lenShape)]\n", " #print ''.join(evList)\n", " return eval(''.join(evList))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Global parameters" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "width = 512000\n", "height = 512000\n", "\n", "domain_sizes = [16, 32, 64, 128, 256]#, 512, 1024, 2048, 4096]\n", "reference_domain_size = 4 * max(domain_sizes)\n", "\n", "\n", "#schemes = [\"FBL\"] \n", "schemes = [\"FBL\", \"CTCS\", \"KP\", \"CDKLM\"]\n", "\n", "#Timestep size \n", "dt = 8000/reference_domain_size\n", " \n", "g = 9.81\n", "r = 0.0\n", "\n", "# Coriolis parameters: f + beta * y\n", "f = 8.0e-5\n", "\n", "timesteps = 5\n", "\n", "end_time = (timesteps - 0.01)*dt\n", "make_netCDF = False\n", "\n", "print(\"Timesteps = \" + str(end_time / dt))" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "def initDataBump(h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f):\n", " \n", " waterHeight = 50\n", " \n", " def my_exp(i, j):\n", " size = 0.3\n", " x = (i + 0.5 - reference_domain_size/2.0) / float(reference_domain_size)\n", " y = (j + 0.5 - reference_domain_size/2.0) / float(reference_domain_size)\n", " return np.exp(-10*(x*x/(size*size)+y*y/(size*size))) * (np.sqrt(x**2 + y**2) < size)\n", " \n", " def my_cos(i, j):\n", " size = 0.6\n", " x = 2*(i + 0.5 - reference_domain_size/2.0) / float(reference_domain_size)\n", " y = 2*(j + 0.5 - reference_domain_size/2.0) / float(reference_domain_size)\n", " r = np.sqrt(x**2 + y**2)\n", " return 0.5*(1.0 + np.cos(np.pi*r/size)) * (r < size)\n", " \n", " #Generate disturbance at reference scale and downsample \n", " disturbance = np.fromfunction(lambda i, j: my_cos(i,j), (reference_domain_size, reference_domain_size)) \n", " disturbance = rebin(disturbance, nx, ny)\n", " \n", " validCells = [ghosts[2], eta0.shape[0] - ghosts[0], ghosts[3], eta0.shape[1] - ghosts[1]]\n", " \n", " eta0.fill(0.0)\n", " eta0[validCells[0]:validCells[1], validCells[2]:validCells[3]] += (0.01*disturbance)\n", " h0.fill(waterHeight)\n", " u0.fill(0.0)\n", " v0.fill(0.0)\n", "\n", "def initDataBalancedBump(h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f):\n", " bump_posx = 0.5\n", " bump_posy = 0.5\n", " bump_height = 0.25\n", " bump_width_factor = 20*nx\n", " waterHeight = 50 \n", " initializeBalancedBumpOverPoint(eta0, u0, v0, # allocated buffers to be filled with data (output)\n", " nx, ny, dx, dy, ghosts, # grid data\n", " bump_posx, bump_posy, # relative placement of bump center\n", " bump_height, bump_width_factor, # bump information\n", " f, waterHeight, # parameters defined at the bump centre (coriolis force, water depth)\n", " g)\n", " \n", " # Scale eta to be out of geostrophic balance\n", " eta0 *= 1.1\n", " h0.fill(waterHeight);\n", " \n", "def initData(h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f):\n", " initDataBump(h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f)\n", " \n", "def testInitData(domain_size):\n", " \n", " nx = domain_size\n", " ny = domain_size\n", " \n", " dx = float(width/nx)\n", " dy = float(height/ny)\n", " \n", " ghosts = [1, 1, 1, 1] # north, east, south, west\n", " dataShape = (ny + ghosts[0]+ghosts[2], \n", " nx + ghosts[1]+ghosts[3])\n", "\n", " h0 = np.zeros(dataShape, dtype=np.float32);\n", " eta0 = np.zeros(dataShape, dtype=np.float32);\n", " u0 = np.zeros((dataShape[0], dataShape[1]+1), dtype=np.float32);\n", " v0 = np.zeros((dataShape[0]+1, dataShape[1]), dtype=np.float32);\n", " \n", " initData(h0, eta0, u0, v0, nx, ny, dx, dy, ghosts, g, f)\n", " \n", " return eta0\n", " \n", "plt.figure()\n", "for i, domain_size in enumerate(domain_sizes):\n", " eta0 = testInitData(domain_size)\n", " plt.subplot(1, len(domain_sizes)+1, i+1)\n", " plt.imshow(eta0, interpolation='nearest')\n", " print(\"Max={:.05f}, min={:.05f}, sum={:.010f}\".format(np.max(eta0), np.min(eta0), np.sum(eta0/(domain_size*domain_size))))" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "def plotData(eta0, u0, v0, eta1, u1, v1):\n", " fig, axarr = plt.subplots(2, 3)\n", " axarr[0, 0].imshow(eta0, interpolation=\"nearest\")\n", " axarr[0, 1].imshow(u0, interpolation=\"nearest\")\n", " axarr[0, 2].imshow(v0, interpolation=\"nearest\")\n", " axarr[1, 0].imshow(eta1, interpolation=\"nearest\")\n", " axarr[1, 1].imshow(u1, interpolation=\"nearest\")\n", " axarr[1, 2].imshow(v1, interpolation=\"nearest\")\n", " print(\"Eta0: Maximum = {:.05f}, minimum = {:.05f}\".format(np.max(eta0), np.min(eta0)))\n", " print(\"Eta1: Maximum = {:.05f}, minimum = {:.05f}\".format(np.max(eta1), np.min(eta1)))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Forward Backward Linear" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "def runFBL(domain_size):\n", " #Clean up old simulator if any:\n", " if 'fbl_sim' in globals():\n", " fbl_sim.cleanUp()\n", " \n", " nx = domain_size\n", " ny = domain_size\n", " \n", " dx = float(width/nx)\n", " dy = float(height/ny)\n", " \n", " ghosts = [0, 0, 0, 0] # north, east, south, west\n", " dataShape = (ny + ghosts[0]+ghosts[2], \n", " nx + ghosts[1]+ghosts[3])\n", "\n", " h0 = np.zeros(dataShape, dtype=np.float32);\n", " eta0 = np.zeros(dataShape, dtype=np.float32);\n", " u0 = np.zeros((dataShape[0], dataShape[1]+1), dtype=np.float32);\n", " v0 = np.zeros((dataShape[0]+1, dataShape[1]), dtype=np.float32);\n", "\n", " # Generate bump in geostrophic balance\n", " initData(h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f)\n", "\n", " #Initialize simulator\n", " reload(FBL)\n", " fbl_sim = FBL.FBL(cl_ctx, \\\n", " h0, eta0, u0, v0, \\\n", " nx, ny, \\\n", " dx, dy, dt, \\\n", " g, f, r, \\\n", " write_netcdf=make_netCDF)\n", "\n", " t = fbl_sim.step(end_time)\n", " eta1, u1, v1 = fbl_sim.download()\n", " print \"\\t\\tt=\" + str(t) + \"\\tMax eta: \" + str(np.max(eta1))\n", " \n", " return [eta0, u0, v0, eta1, u1, v1]\n", "\n", "[eta0, u0, v0, eta1, u1, v1] = runFBL(16)\n", "plotData(eta0, u0, v0, eta1, u1, v1)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "if make_netCDF:\n", " fbl_sim.cleanUp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Centered in time, centered in space" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "scrolled": false }, "outputs": [], "source": [ "#Centered in time, centered in space\n", "\n", "def runCTCS(domain_size):\n", " #Clean up old simulator if any:\n", " if 'ctcs_sim' in globals():\n", " ctcs_sim.cleanUp()\n", " \n", " nx = domain_size\n", " ny = domain_size\n", " \n", " dx = float(width/nx)\n", " dy = float(height/ny)\n", " \n", " ghosts = [1,1,1,1] # north, east, south, west\n", " validDomain = np.array([1,1,1,1])\n", " dataShape = (ny + ghosts[0]+ghosts[2], \n", " nx + ghosts[1]+ghosts[3])\n", "\n", " h0 = np.zeros(dataShape, dtype=np.float32);\n", " eta0 = np.zeros(dataShape, dtype=np.float32);\n", " u0 = np.zeros((dataShape[0], dataShape[1]+1), dtype=np.float32);\n", " v0 = np.zeros((dataShape[0]+1, dataShape[1]), dtype=np.float32); \n", "\n", " initData(h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f)\n", " \n", " # Eddy viscocity parameter\n", " A = 0.5*dx\n", " \n", " reload(CTCS)\n", " ctcs_sim = CTCS.CTCS(cl_ctx, \\\n", " h0, eta0, u0, v0, \\\n", " nx, ny, dx, dy, dt, \\\n", " g, f, r, A, \\\n", " write_netcdf=make_netCDF)\n", "\n", " t = ctcs_sim.step(end_time)\n", " eta1, u1, v1 = ctcs_sim.download()\n", " \n", " # Remove ghost cells\n", " eta1 = eta1[validDomain[3]:-validDomain[1], validDomain[2]:-validDomain[0]]\n", " \n", " print \"\\t\\tt=\" + str(t) + \"\\tMax eta: \" + str(np.max(eta1))\n", " \n", " return [eta0, u0, v0, eta1, u1, v1]\n", "\n", "[eta0, u0, v0, eta1, u1, v1] = runCTCS(16)\n", "plotData(eta0, u0, v0, eta1, u1, v1)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "if make_netCDF:\n", " ctcs_sim.cleanUp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## CDKLM 16" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "scrolled": false }, "outputs": [], "source": [ "def runCDKLM(domain_size):\n", " #Clean up old simulator if any:\n", " if 'cdklm_sim' in globals():\n", " cdklm_sim.cleanUp()\n", "\n", " #Coriolis well balanced reconstruction scheme\n", " \n", " nx = domain_size\n", " ny = domain_size\n", " \n", " dx = float(width/nx)\n", " dy = float(height/ny)\n", "\n", " ghosts = np.array([2,2,2,2]) # north, east, south, west\n", " validDomain = np.array([2,2,2,2])\n", " dataShape = (ny + ghosts[0]+ghosts[2], \n", " nx + ghosts[1]+ghosts[3])\n", "\n", " Hi = np.zeros((dataShape[0]+1, dataShape[1]+1), dtype=np.float32)\n", " eta0 = np.zeros(dataShape, dtype=np.float32)\n", " u0 = np.zeros(dataShape, dtype=np.float32)\n", " v0 = np.zeros(dataShape, dtype=np.float32)\n", "\n", " initData(Hi, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f)\n", "\n", " #Initialize simulator\n", " reload(CDKLM16)\n", " cdklm_sim = CDKLM16.CDKLM16(cl_ctx, \\\n", " eta0, u0, v0, Hi, \\\n", " nx, ny, dx, dy, dt, \\\n", " g, f, r, \\\n", " rk_order=2, \n", " write_netcdf=make_netCDF)\n", "\n", "\n", " t = cdklm_sim.step(end_time)\n", " eta1, u1, v1 = cdklm_sim.download()\n", " \n", " # Remove ghost cells\n", " eta1 = eta1[validDomain[3]:-validDomain[1], validDomain[2]:-validDomain[0]]\n", " \n", " print \"\\t\\tt=\" + str(t) + \"\\tMax eta: \" + str(np.max(eta1))\n", " \n", " return [eta0, u0, v0, eta1, u1, v1]\n", "\n", "[eta0, u0, v0, eta1, u1, v1] = runCDKLM(16)\n", "plotData(eta0, u0, v0, eta1, u1, v1)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "if make_netCDF:\n", " cdklm_sim.cleanUp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Kurganov-Petrova 2007" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "scrolled": false }, "outputs": [], "source": [ "def runKP(domain_size):\n", " #Clean up old simulator if any:\n", " if 'kp07_sim' in globals():\n", " kp07_sim.cleanUp()\n", " \n", " # Kurganov-Petrova 2007\n", " \n", " nx = domain_size\n", " ny = domain_size\n", " \n", " dx = float(width/nx)\n", " dy = float(height/ny)\n", " \n", " ghosts = np.array([2,2,2,2]) # north, east, south, west\n", " validDomain = np.array([2,2,2,2])\n", " dataShape = (ny + ghosts[0]+ghosts[2], \n", " nx + ghosts[1]+ghosts[3])\n", "\n", " Hi = np.zeros((dataShape[0]+1, dataShape[1]+1), dtype=np.float32)\n", " eta0 = np.zeros(dataShape, dtype=np.float32)\n", " u0 = np.zeros(dataShape, dtype=np.float32)\n", " v0 = np.zeros(dataShape, dtype=np.float32)\n", "\n", " initData(Hi, eta0, u0, v0, \\\n", " nx, ny, dx, dy, ghosts, \\\n", " g, f)\n", "\n", " #Initialize simulator\n", " reload(KP07)\n", " kp07_sim = KP07.KP07(cl_ctx, \\\n", " eta0, Hi, u0, v0, \\\n", " nx, ny, dx, dy, dt, \\\n", " g, f, r, \\\n", " write_netcdf=make_netCDF,\\\n", " use_rk2=True)\n", "\n", " t = kp07_sim.step(end_time)\n", " eta1, u1, v1 = kp07_sim.download()\n", " \n", " # Remove ghost cells\n", " eta1 = eta1[validDomain[3]:-validDomain[1], validDomain[2]:-validDomain[0]]\n", " \n", " print \"\\t\\tt=\" + str(t) + \"\\tMax eta: \" + str(np.max(eta1))\n", " \n", " return [eta0, u0, v0, eta1, u1, v1]\n", "\n", "[eta0, u0, v0, eta1, u1, v1] = runKP(16)\n", "plotData(eta0, u0, v0, eta1, u1, v1)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "if make_netCDF:\n", " kp07_sim.cleanUp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Control " ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "scrolled": false }, "outputs": [], "source": [ "for scheme in schemes:\n", " print \"Scheme: \" + scheme\n", " \n", " data = {};\n", " \n", " # Make reference solution\n", " print \"\\tDomain size (reference solution): \" + str(reference_domain_size)\n", " [_, _, _, eta1_ref, _, _] = eval(\"run\" + scheme + \"(\" + str(reference_domain_size) + \")\")\n", " \n", " data[str(reference_domain_size)] = eta1_ref\n", "\n", " # Run all domain sizes\n", " for domain_size in domain_sizes:\n", " print \"\\tDomain size: \" + str(domain_size)\n", " [_, _, _, eta1, _, _] = eval(\"run\" + scheme + \"(\" + str(domain_size) + \")\")\n", " \n", " data[str(domain_size)] = eta1\n", " \n", " \n", " out_filename = imgdir + \"/\" + scheme + \"_data.npz\"\n", " np.savez(out_filename, **data)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "error = np.zeros([len(schemes), len(domain_sizes)])\n", "\n", "for k, scheme in enumerate(schemes):\n", " print \"Scheme: \" + scheme\n", " \n", " in_filename = imgdir + \"/\" + scheme + \"_data.npz\"\n", " npzfile = np.load(in_filename)\n", " \n", " #Get reference\n", " eta1_ref = npzfile[str(reference_domain_size)].astype(np.float64)\n", " \n", " # Run all domain sizes\n", " for l, domain_size in enumerate(domain_sizes):\n", " eta1 = npzfile[str(domain_size)].astype(np.float64)\n", " \n", " print(\"Max={:.05f}, min={:.05f}, sum={:.010f}\".format(np.max(eta1), np.min(eta1), np.sum(eta1/(domain_size*domain_size))))\n", "\n", " #ver 1 : downsample til minste opplk\u00f8sning\n", " \"\"\"\n", " eta1_ref_downsampled = rebin(eta1_ref, min(domain_sizes), min(domain_sizes))\n", " eta1_downsampled = rebin(eta1, min(domain_sizes), min(domain_sizes))\n", " tmp =eta1_ref_downsampled - eta1_downsampled\n", " error[k, l] = np.linalg.norm(tmp.flatten(), ord=2)\n", " \"\"\"\n", " \n", " #\"\"\"\n", " #ver 2: downsample til current oppl\u00f8sning\n", " eta1_ref_downsampled = rebin(eta1_ref, domain_size, domain_size)\n", " eta1_downsampled = eta1\n", " tmp =eta1_ref_downsampled - eta1_downsampled\n", " error[k, l] = np.linalg.norm(tmp, ord='fro') / (domain_size*domain_size)\n", " #\"\"\"\n", " \n", " \"\"\"\n", " #ver 3: upsample til refereanseoppl\u00f8sning\n", " eta1_ref_downsampled = eta1_ref\n", " upsampling = np.ones(np.divide(eta1_ref.shape, eta1.shape))\n", " eta1_downsampled = np.kron(eta1, upsampling)\n", " tmp =eta1_ref_downsampled - eta1_downsampled\n", " error[k, l] = np.linalg.norm(tmp.flatten(), ord=2)\n", " \"\"\"\n", " \n", " \n", "fig = plt.figure()\n", "setBwStyles(fig.gca())\n", "\n", "x = np.linspace(domain_sizes[0], domain_sizes[-1], 100);\n", "\n", "#scaling = np.min(error[:,0]) * domain_sizes[0]**0.5 * 0.5\n", "#plt.loglog(x, scaling/(np.sqrt(x)), '-', color='gray', label='Order 0.5')\n", "\n", "scaling = np.max(error[:,0]) * domain_sizes[0] * 2\n", "plt.loglog(x, scaling/x, '-', color='gray', label='Order 1')\n", "\n", "scaling = np.min(error[:,0]) * domain_sizes[0]**2 * 0.5\n", "plt.loglog(x, scaling/(x*x), '-', color='gray', label='Order 2')\n", "\n", "for k in range(len(schemes)):\n", " print \"Scheme \" + str(schemes[k])\n", " for l in range(len(domain_sizes)):\n", " print \"\\tDomain size: \" + str(domain_sizes[l]) + \": \" + str(error[k,l])\n", " plt.loglog(domain_sizes, error[k,:], label=schemes[k], markersize=15)\n", "#plt.loglog(domain_sizes, np.abs(error[0,:]-error[1,:]), label=\"Diff\", markersize=15)\n", " \n", "plt.xlabel('Number of cells')\n", "plt.ylabel('Error')\n", "plt.legend(markerscale=0.5)\n", "\n", "plt.savefig(imgdir + \"/\" + \"convergence.pdf\")" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [] } ], "metadata": { "git": { "suppress_outputs": true }, "kernelspec": { "display_name": "Python 2", "language": "python", "name": "python2" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 2 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython2", "version": "2.7.6" } }, "nbformat": 4, "nbformat_minor": 1 }