{ "cells": [ { "cell_type": "raw", "metadata": {}, "source": [ "Text provided under a Creative Commons Attribution license, CC-BY. All code is made available under the FSF-approved BSD-3 license. (c) Lorena A. Barba, Gilbert F. Forsyth 2017. Thanks to NSF for support via CAREER award #1149784." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "[@LorenaABarba](https://twitter.com/LorenaABarba)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "12 steps to Navier–Stokes\n", "======\n", "***" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hello! Welcome to the **12 steps to Navier–Stokes**. This is a practical module that is used in the beginning of an interactive Computational Fluid Dynamics (CFD) course taught by [Prof. Lorena Barba](http://lorenabarba.com) since Spring 2009 at Boston University. The course assumes only basic programming knowledge (in any language) and of course some foundation in partial differential equations and fluid mechanics. The practical module was inspired by the ideas of Dr. Rio Yokota, who was a post-doc in Barba's lab, and has been refined by Prof. Barba and her students over several semesters teaching the course. The course is taught entirely using Python and students who don't know Python just learn as we work through the module.\n", "\n", "This [Jupyter notebook](https://jupyter-notebook.readthedocs.io/en/stable/) will lead you through the first step of programming your own Navier–Stokes solver in Python from the ground up. We're going to dive right in. Don't worry if you don't understand everything that's happening at first, we'll cover it in detail as we move forward and you can support your learning with the videos of [Prof. Barba's lectures on YouTube](http://www.youtube.com/playlist?list=PL30F4C5ABCE62CB61).\n", "\n", "For best results, after you follow this notebook, prepare your own code for Step 1, either as a Python script or in a clean Jupyter notebook.\n", "\n", "To execute this Notebook, we assume you have invoked the notebook server using: `jupyter notebook`." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Step 1: 1-D Linear Convection\n", "-----\n", "***" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The 1-D Linear Convection equation is the simplest, most basic model that can be used to learn something about CFD. It is surprising that this little equation can teach us so much! Here it is:\n", "\n", "$$\\frac{\\partial u}{\\partial t} + c \\frac{\\partial u}{\\partial x} = 0$$\n", "\n", "With given initial conditions (understood as a *wave*), the equation represents the propagation of that initial *wave* with speed $c$, without change of shape. Let the initial condition be $u(x,0)=u_0(x)$. Then the exact solution of the equation is $u(x,t)=u_0(x-ct)$.\n", "\n", "We discretize this equation in both space and time, using the Forward Difference scheme for the time derivative and the Backward Difference scheme for the space derivative. Consider discretizing the spatial coordinate $x$ into points that we index from $i=0$ to $N$, and stepping in discrete time intervals of size $\\Delta t$.\n", "\n", "From the definition of a derivative (and simply removing the limit), we know that:\n", "\n", "$$\\frac{\\partial u}{\\partial x}\\approx \\frac{u(x+\\Delta x)-u(x)}{\\Delta x}$$\n", "\n", "Our discrete equation, then, is:\n", "\n", "$$\\frac{u_i^{n+1}-u_i^n}{\\Delta t} + c \\frac{u_i^n - u_{i-1}^n}{\\Delta x} = 0 $$\n", "\n", "Where $n$ and $n+1$ are two consecutive steps in time, while $i-1$ and $i$ are two neighboring points of the discretized $x$ coordinate. If there are given initial conditions, then the only unknown in this discretization is $u_i^{n+1}$. We can solve for our unknown to get an equation that allows us to advance in time, as follows:\n", "\n", "$$u_i^{n+1} = u_i^n - c \\frac{\\Delta t}{\\Delta x}(u_i^n-u_{i-1}^n)$$\n", "\n", "Now let's try implementing this in Python. \n", "\n", "We'll start by importing a few libraries to help us out.\n", "\n", "* `numpy` is a library that provides a bunch of useful matrix operations akin to MATLAB\n", "* `matplotlib` is a 2D plotting library that we will use to plot our results\n", "* `time` and `sys` provide basic timing functions that we'll use to slow down animations for viewing" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# Remember: comments in python are denoted by the pound sign\n", "import numpy #here we load numpy\n", "from matplotlib import pyplot #here we load matplotlib\n", "import time, sys #and load some utilities\n" ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "collapsed": false }, "outputs": [], "source": [ "#this makes matplotlib plots appear in the notebook (instead of a separate window)\n", "%matplotlib inline " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now let's define a few variables; we want to define an evenly spaced grid of points within a spatial domain that is 2 units of length wide, i.e., $x_i\\in(0,2)$. We'll define a variable `nx`, which will be the number of grid points we want and `dx` will be the distance between any pair of adjacent grid points. " ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [], "source": [ "nx = 41 # try changing this number from 41 to 81 and Run All ... what happens?\n", "dx = 2 / (nx-1)\n", "nt = 25 #nt is the number of timesteps we want to calculate\n", "dt = .025 #dt is the amount of time each timestep covers (delta t)\n", "c = 1 #assume wavespeed of c = 1" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We also need to set up our initial conditions. The initial velocity $u_0$ is given as \n", "$u = 2$ in the interval $0.5 \\leq x \\leq 1$ and $u = 1$ everywhere else in $(0,2)$ (i.e., a hat function).\n", "\n", "Here, we use the function `ones()` defining a `numpy` array which is `nx` elements long with every value equal to 1." ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "[ 1. 1. 1. 1. 1. 1. 1. 1. 1. 1. 2. 2. 2. 2. 2. 2. 2. 2.\n", " 2. 2. 2. 1. 1. 1. 1. 1. 1. 1. 1. 1. 1. 1. 1. 1. 1. 1.\n", " 1. 1. 1. 1. 1.]\n" ] } ], "source": [ "u = numpy.ones(nx) #numpy function ones()\n", "u[int(.5 / dx):int(1 / dx + 1)] = 2 #setting u = 2 between 0.5 and 1 as per our I.C.s\n", "print(u)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now let's take a look at those initial conditions using a Matplotlib plot. We've imported the `matplotlib` plotting library `pyplot` and the plotting function is called `plot`, so we'll call `pyplot.plot`. To learn about the myriad possibilities of Matplotlib, explore the [Gallery](http://matplotlib.org/gallery.html) of example plots.\n", "\n", "Here, we use the syntax for a simple 2D plot: `plot(x,y)`, where the `x` values are evenly distributed grid points:" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAXkAAAEACAYAAABWLgY0AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAEu1JREFUeJzt3W2MXOV5xvHr8svG9q6DhJFc8WbSlAilxXEqlRcHNUNV\nyUCrAFGrFEdUoWqaVKRB6odQoSJPFJS0UhIlCFrk1nFEFAckkGpDQInUZD6YAEpSHFMwFEIVUydx\n1AQ+7NrExr77YWbdZT27c3Z93p5n/z9p5RnP49lbw/HF5WfnnHFECACQp2VNDwAAqA4hDwAZI+QB\nIGOEPABkjJAHgIwR8gCQsZEhb/t829+x/ZztZ21/co51d9t+yfY+25vKHxUAsFArCqx5U9LfRsQ+\n2xOSfmj72xHxwvQC29dKemdEXGz7ckn3SbqimpEBAEWNbPIR8fOI2De4PSnpgKTzZi27XtL9gzVP\nSzrL9vqSZwUALNCC9uRtXyRpk6SnZz10nqRXZ9w/pNP/RwAAqFnhkB9s1Twk6bZBowcAtFyRPXnZ\nXqF+wH8tInYPWXJI0gUz7p8/+L3Zz8OFcgBgESLCi/lzhUJe0lckPR8RX57j8T2SbpX0oO0rJL0e\nEYeHLeSCaOW44QbpbW/r6sEHu02PkoWNG6XNm7u6775u06NkodvtqtvtNj1GNuxF5bukAiFv+32S\nPizpWdvPSApJd0jaICkiYntEPGb7OtsvS5qSdMuiJ0Ihk5PS2rVNT5GPiQnp2LGmpwDKNzLkI+IJ\nScsLrPtEKROhkKkpaWys6SnyMT4uHT/e9BRA+TjjNVGTk9LmzZ2mx8jGxIR08cWdpsfIRqfTaXoE\nDBDyiZqa4i9SmcbHpXe8o9P0GNng2GwPQj5Rk5P99olyTEz0X1MgN4R8oqam+u0T5Rgf77+mQG4I\n+QSdPCkdPSqtWdP0JPmgySNXhHyCjhzpB/wy/uuVZmKCJo88ERMJmpxkq6Zs4+M0eeSJkE/Q1BQ/\ndC0bTR65IuQTRJMvH00euSLkE0STLx9NHrki5BNEky8fTR65IuQTRJMvH00euSLkE0STLx9NHrki\n5BNEky8fTR65IuQTRJMvH00euSLkE0STL9/4eP9MYj64DLkh5BNEky/f8uX9D2E5erTpSYByEfIJ\noslXg3155IiQTxBNvhrsyyNHhHyC+MCQanC5YeSIkE8QHxhSDT44BDki5BNEk68GTR45IuQTRJOv\nBk0eOSLkE0STrwZNHjki5BNEk68GTR45IuQTRJOvBk0eOSLkExPRP/2eJl8+mjxyRMgn5ujR/un3\ny5c3PUl+aPLI0ciQt73D9mHb++d4/O2299jeZ/tZ2x8pfUqcwn58dWjyyFGRJr9T0pZ5Hr9V0nMR\nsUnS1ZK+YHtFGcPhdOzHV4cmjxyNDPmI2CvptfmWSFo7uL1W0i8j4s0SZsMQNPnq0OSRozIa9z2S\n9tj+qaQJSR8q4TkxB5p8dWjyyFEZP3jdIumZiDhX0nsl3WubGKoIlxmuDpcaRo7KaPK3SPqcJEXE\nj23/t6RLJP1g2OJut3vqdqfTUafTKWGEpYPLDFeHSw2jLXq9nnq9XinP5SjweWe2L5L0SERcOuSx\neyX9IiI+bXu9+uH+noj41ZC1UeT7YW67dkmPPtr/FeV68UXpAx/o/wq0iW1FhBfzZ0c2edu7JHUk\nrbN9UNI2SWOSIiK2S7pL0ldnvMXyU8MCHuWgyVeHJo8cjQz5iNg64vGfaf63WKJE7MlXhz155Igz\nXhNDk6/OdJNnRxE5IeQTQ5OvzsqV/ctFHDvW9CRAeQj5xNDkq8W+PHJDyCeGk6GqxQlRyA0hnxgu\na1AtLm2A3BDyiaHJV4smj9wQ8omhyVeLJo/cEPKJoclXiyaP3BDyiaHJV4smj9wQ8omhyVeLJo/c\nEPKJoclXiyaP3BDyCYngZKiq0eSRG0I+IceOScuWSWNjTU+SL5o8ckPIJ4T9+OrR5JEbQj4h7MdX\njyaP3BDyCaHJV48mj9wQ8gmhyVePJo/cEPIJoclXjyaP3BDyCaHJV48mj9wQ8gmhyVePJo/cEPIJ\n4aP/qseHeSM3hHxCONu1enz8H3JDyCeEJl89mjxyQ8gnhCZfvbEx6cSJ/iUkgBwQ8gnhB6/Vs2nz\nyAshnxDeQlkP3kaJnBDyCaHJ14O3USInhHxCaPL1oMkjJ4R8Qmjy9aDJIycjQ972DtuHbe+fZ03H\n9jO2/9P2d8sdEdNo8vWgySMnRZr8Tklb5nrQ9lmS7pX0xxHxO5L+tKTZMAtNvh40eeRkZMhHxF5J\nr82zZKukhyPi0GD9/5Y0G2ahydeDJo+clLEn/y5JZ9v+ru3v2765hOfEEDT5etDkkZMVJT3H70r6\nA0njkp60/WREvDxscbfbPXW70+mo0+mUMMLSQJOvB00eTev1eur1eqU8lyNi9CJ7g6RHImLjkMdu\nl7QqIj49uP+vkh6PiIeHrI0i3w+nO35cWr26/6vd9DR5+8xnpF//WrrrrqYnAfpsKyIW9Te/6HaN\nB1/D7JZ0le3lttdIulzSgcUMg7lNt3gCvno0eeRk5HaN7V2SOpLW2T4oaZukMUkREdsj4gXb35K0\nX9IJSdsj4vkKZ16S2I+vD3vyyMnIkI+IrQXWfF7S50uZCEOxH18fmjxywhmviaDJ14cmj5wQ8omg\nydeHJo+cEPKJoMnXhyaPnBDyiaDJ14cmj5wQ8omgydeHJo+cEPKJoMnXhyaPnBDyiaDJ14cmj5wQ\n8omYnKTJ12X16v5lDU6caHoS4MwR8omYmqLJ18Vmywb5IOQTwXZNvdiyQS4I+UTwg9d60eSRC0I+\nETT5etHkkQtCPhE0+XrR5JELQj4RNPl60eSRC0I+ETT5etHkkQtCPhE0+XrR5JELQj4RNPl60eSR\nC0I+ETT5etHkkQtCPgEnTvRPs1+9uulJlg6aPHJByCdgeqvGbnqSpYMmj1wQ8glgP75+NHnkgpBP\nAPvx9aPJIxeEfAJo8vWjySMXhHwCaPL1o8kjF4R8Amjy9aPJIxeEfAJo8vWjySMXhHwCaPL1o8kj\nF4R8Amjy9aPJIxcjQ972DtuHbe8fse73bB+3/cHyxoPEh3g3YXyckEceijT5nZK2zLfA9jJJ/yDp\nW2UMhbfiQ7zrt2aNdPSodPJk05MAZ2ZkyEfEXkmvjVj2N5IekvSLMobCW9Hk67d8ubRqlXTkSNOT\nAGfmjPfkbZ8r6YaI+GdJXF2lAjT5ZkxM8MNXpG9FCc/xJUm3z7g/b9B3u91TtzudjjqdTgkj5I0m\n34zpffn165ueBEtNr9dTr9cr5bkcEaMX2RskPRIRG4c89sr0TUnnSJqS9FcRsWfI2ijy/fBWN94o\n3Xyz9EF+pF2rSy+Vvv51aeNpRz1QL9uKiEXtlBRt8tYcDT0ifnPGIDvV/5/BaQGPxaPJN4N32CAH\nI0Pe9i5JHUnrbB+UtE3SmKSIiO2zllPTK8CefDPYk0cORoZ8RGwt+mQR8RdnNg6G4WSoZnBCFHLA\nGa8J4LIGzeDSBsgBIZ8AmnwzaPLIASGfAJp8M2jyyAEh33InT/bPulyzpulJlh6aPHJAyLfc0aP9\n0+uXL296kqWHJo8cEPItx358c2jyyAEh33LsxzeHJo8cEPItR5NvDk0eOSDkW44m3xyaPHJAyLcc\nTb45NHnkgJBvOZp8c2jyyAEh33I0+ebQ5JEDQr7laPLNockjB4R8y9Hkm0OTRw4I+ZbjA0OaM93k\n+TAzpIyQbzk+MKQ5K1ZIK1dKb7zR9CTA4hHyLUeTbxYfAYjUEfItR5NvFh8BiNQR8i1Hk28WTR6p\nI+RbjibfLJo8UkfItxxNvlk0eaSOkG85mnyzaPJIHSHfcjT5ZtHkkTpCvuVo8s2iySN1hHzLcVmD\nZnFpA6SOkG+xCC5Q1jQuUobUEfIt9sYb/dPqV6xoepKliyaP1BHyLUaLbx5NHqkbGfK2d9g+bHv/\nHI9vtf2jwdde25eWP+bSxH5882jySF2RJr9T0pZ5Hn9F0u9HxHsk3SXpX8oYDDT5NqDJI3Ujd3sj\nYq/tDfM8/tSMu09JOq+MwUCTbwOaPFJX9p78X0p6vOTnXLJo8s2jySN1pb1vw/bVkm6RdNV867rd\n7qnbnU5HnU6nrBGyQ5NvHk0eTej1eur1eqU8l6PAZ5sNtmseiYiNczy+UdLDkq6JiB/P8zxR5Puh\n7xvfkHbvlh54oOlJlq4DB6Qbb5ReeKHpSbCU2VZEeDF/tuh2jQdfw775heoH/M3zBTwWjibfPJo8\nUjdyu8b2LkkdSetsH5S0TdKYpIiI7ZLulHS2pH+ybUnHI+Ky6kZeOrg4WfO4QBlSV+TdNVtHPP5R\nSR8tbSKcwsXJmjd9gbIIyYv6xzLQLM54bTGafPPGxvq/HjvW7BzAYhHyLUaTbwcuN4yUEfItRpNv\nB/blkTJCvsVo8u1Ak0fKCPkWo8m3A00eKSPkW4wm3w40eaSMkG8xmnw70OSRMkK+xWjy7UCTR8oI\n+RajybcDTR4pI+RbjCbfDjR5pIyQbzGafDvQ5JEyQr6lpk+jnz6tHs2hySNlhHxLTbd4LorVPJo8\nUkbItxT78e1Bk0fKCPmW4gND2oMPDkHKCPmW4kO824MP80bKCPmWosm3B00eKSPkW4om3x40eaSM\nkG8pmnx70OSRMkK+pTgRqj14CyVSRsi3FG+hbA/eQomUEfItRZNvD5o8UkbItxRNvj1WrZLefFM6\nfrzpSYCFI+RbiibfHjbvsEG6CPmWosm3C/vySBUh31I0+XZhXx6pIuRbiibfLjR5pIqQbymafLvQ\n5JGqkSFve4ftw7b3z7Pmbtsv2d5ne1O5Iy5NNPl2ockjVUWa/E5JW+Z60Pa1kt4ZERdL+pik+0qa\nbUmjybcLTR6pGhnyEbFX0mvzLLle0v2DtU9LOsv2+nLGW7po8u1Ck0eqytiTP0/SqzPuHxr8Hs4A\nTb5daPJI1Yq6v+G559b9HdM0OSmtXdv0FJh2zjnSHXdIn/1s05MAC1NGyB+SdMGM++cPfm+om27q\nnrp95ZUdbd7cKWGE/Kxa1f9CO9x5p/Txjzc9BZaK732vpyef7J26/8UvLv65HBGjF9kXSXokIi4d\n8th1km6NiD+yfYWkL0XEFXM8TxT5fgCA/2dbEeHF/NmRTd72LkkdSetsH5S0TdKYpIiI7RHxmO3r\nbL8saUrSLYsZBABQvkJNvrRvRpMHgAU7kybPGa8AkDFCHgAyRsgDQMYIeQDIGCEPABkj5AEgY4Q8\nAGSMkAeAjBHyAJAxQh4AMkbIA0DGCHkAyBghDwAZI+QBIGOEPABkjJAHgIwR8gCQMUIeADJGyANA\nxgh5AMgYIQ8AGSPkASBjhDwAZIyQB4CMEfIAkDFCHgAyRsgDQMYIeQDIWKGQt32N7Rds/5ft24c8\n/nbbe2zvs/2s7Y+UPikAYMFGhrztZZLukbRF0m9Lusn2JbOW3SrpuYjYJOlqSV+wvaLsYfFWvV6v\n6RGywutZHl7L9ijS5C+T9FJE/CQijkt6QNL1s9aEpLWD22sl/TIi3ixvTAzDX6Ry8XqWh9eyPYqE\n/HmSXp1x/38GvzfTPZLebfunkn4k6bZyxgMAnImyfvC6RdIzEXGupPdKutf2REnPDQBYJEfE/Avs\nKyR1I+Kawf2/kxQR8Y8z1jwq6XMR8cTg/r9Luj0ifjDrueb/ZgCAoSLCi/lzRX44+n1Jv2V7g6Sf\nSfozSTfNWvMTSX8o6Qnb6yW9S9IrZQ0JAFickSEfESdsf0LSt9Xf3tkREQdsf6z/cGyXdJekr9re\nP/hjn4qIX1U2NQCgkJHbNQCAdFVyxuuok6cGa+62/dLgBKpNVcyRiwIno73f9uu2/2Pw9fdNzJkC\n2ztsH57xr85hazg2Cxj1WnJcLozt821/x/Zzg5NKPznHuoUdnxFR6pf6/+N4WdIGSSsl7ZN0yaw1\n10r65uD25ZKeKnuOXL4Kvp7vl7Sn6VlT+JJ0laRNkvbP8TjHZnmvJcflwl7P35C0aXB7QtKLZWRn\nFU2+yMlT10u6X5Ii4mlJZw1+YIvTFXk9JYkfahcQEXslvTbPEo7Nggq8lhLHZWER8fOI2De4PSnp\ngE4/J2nBx2cVIV/k5KnZaw4NWYO+Iq+nJF05+OfbN22/u57RssSxWS6Oy0WwfZH6/0p6etZDCz4+\nub5MHn4o6cKIOGL7Wkn/pv7bWIEmcVwuwuBE0ock3TZo9GekiiZ/SNKFM+6fP/i92WsuGLEGfSNf\nz4iYjIgjg9uPS1pp++z6RswKx2ZJOC4XbnBhx4ckfS0idg9ZsuDjs4qQP3XylO0x9U+e2jNrzR5J\nfy6dOqP29Yg4XMEsORj5es7ck7N9mfpvjeU8hblZc+8Vc2wuzJyvJcflonxF0vMR8eU5Hl/w8Vn6\ndk0UOHkqIh6zfZ3tlyVNSbql7DlyUeT1lPQntv9a0nFJRyV9qLmJ2832LkkdSetsH5S0TdKYODYX\nbNRrKY7LBbH9PkkflvSs7WfUv7rvHeq/s27RxycnQwFAxvj4PwDIGCEPABkj5AEgY4Q8AGSMkAeA\njBHyAJAxQh4AMkbIA0DG/g9WD4D/jFcM7gAAAABJRU5ErkJggg==\n", "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "pyplot.plot(numpy.linspace(0, 2, nx), u);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Why doesn't the hat function have perfectly straight sides? Think for a bit." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now it's time to implement the discretization of the convection equation using a finite-difference scheme. \n", "\n", "For every element of our array `u`, we need to perform the operation $u_i^{n+1} = u_i^n - c \\frac{\\Delta t}{\\Delta x}(u_i^n-u_{i-1}^n)$\n", "\n", "We'll store the result in a new (temporary) array `un`, which will be the solution $u$ for the next time-step. We will repeat this operation for as many time-steps as we specify and then we can see how far the wave has convected. \n", "\n", "We first initialize our placeholder array `un` to hold the values we calculate for the $n+1$ timestep, using once again the NumPy function `ones()`.\n", "\n", "Then, we may think we have two iterative operations: one in space and one in time (we'll learn differently later), so we'll start by nesting one loop inside the other. Note the use of the nifty `range()` function. When we write: `for i in range(1,nx)` we will iterate through the `u` array, but we'll be skipping the first element (the zero-th element). *Why?*" ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "collapsed": false }, "outputs": [], "source": [ "un = numpy.ones(nx) #initialize a temporary array\n", "\n", "for n in range(nt): #loop for values of n from 0 to nt, so it will run nt times\n", " un = u.copy() ##copy the existing values of u into un\n", " for i in range(1, nx): ## you can try commenting this line and...\n", " #for i in range(nx): ## ... uncommenting this line and see what happens!\n", " u[i] = un[i] - c * dt / dx * (un[i] - un[i-1])\n", " \n", " " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Note**—We will learn later that the code as written above is quite inefficient, and there are better ways to write this, Python-style. But let's carry on.\n", "\n", "Now let's try plotting our `u` array after advancing in time." ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAXkAAAEACAYAAABWLgY0AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAHLJJREFUeJzt3XuUVNWZxuHfh9wEFC9RFAXUKBIUB6IoCIMFMXIfLipI\nNBqMCWY0Q0YzGpPMyIxOIkk00aW5mFEczcJIuDZ3BC0RVIIXbFA0QUAIICpoECQDtHv+2IV2muqu\n6upTtatOvc9atazqOl31Wev0y67v7LOPOecQEZF4ahS6ABERyR+FvIhIjCnkRURiTCEvIhJjCnkR\nkRhTyIuIxFjGkDezk83sKTN7zcxWm9m/1LLdfWb2ZzNbZWZdoy9VRETqq3EW2xwAbnLOrTKzVsBL\nZrbIOffGwQ3MbCDweefcGWZ2AfBroEd+ShYRkWxlHMk7595xzq1K3d8NrAVOqrHZMODR1DYrgNZm\n1ibiWkVEpJ7q1ZM3s1OArsCKGk+dBGyu9ngLh/5DICIiBZZ1yKdaNVOB8akRvYiIFLlsevKYWWN8\nwD/mnJuVZpMtQLtqj09O/azm62ihHBGRHDjnLJffy3Yk/zDwunPu3lqerwCuBjCzHsCHzrnt6TZ0\nzukW0e32228PXkOcbvo89VkW660hMo7kzawXcCWw2sxeARzwfaCDz2z3oHNunpkNMrN1wB5gbIOq\nEhGRSGQMeefccuCwLLa7MZKKREQkMln15KU4JRKJ0CXESjl8ni+/DDNnQl0dgObN4etfhxNOyP19\nyuGzLBXW0H5Pvd7MzBXy/UTE27kTfvhDmD4drr0WDj+89m23boWpU+EHP4Abb4TGGgoGZ2a4HA+8\nKuRFYuyTT+Dhh33AX3op3HEHHHNM5t9buxa+/W149124/37o0yf/tUrtFPIicogXX4QbboDDDoMH\nHoBu3er3+875Ef3NN/uQ/+lP4cQT81Or1K0hIa9VKEViZscOuP56GDLE/3fZsvoHPIAZXH45vP46\nnHwydOkCP/857N8ffc2SPwp5kRiZMgU6d4YmTXzLZexYaNTAv/JWreCuu2D5cpg/3/+DUVkZTb2S\nf2rXiMTEihUwdKgP4nPPzc97OOd7/Hfe6WfqHH10ft5H/p568iJlbscOH+y/+AUMH57/9/vOd2DD\nBj8d03KKHqkPhbxIGfvkEz+C79QJ7r67MO+5b58/GHvZZfDd7xbmPcuZQl6kjE2cCLNmwTPP+F58\noWzaBN27+7n3vXoV7n3LkUJepEwtXQqjRsHKldCuXebtozZ3rp/B8/LLcNxxhX//cqEplCJlaPt2\n+MpX4JFHwgQ8wODBcNVV/lZVFaYGqZtCXqQEVVXBlVfC174GAwaEreWOO2DvXvjRj8LWIempXSNS\ngm6/3bdqnnyyONaW2boVzjsPHnsMvvSl0NXEj3ryImVk0SJ/ktNLLzVspcioLVkCX/2qX06hbdvQ\n1cSLevIiZWLLFrjmGvjd74or4MGP4K+/HsaMgQMHQlcjBynkRUpEVRVccYVf/rdv39DVpPeDH0Cz\nZr6dJMVB7RqREjFtml9DZsWKhq9Hk0/bt8MXvgCvvaZVK6Oido1IzDnnT3q67bbiDniANm38zJ97\n7w1diYBG8iIlIZmEceP8sr+HZbzicngbN/q1dNavh9atQ1dT+jSSF4m5iRP9GjGlEPAAp5wC/fvD\nb34TuhLRSF6kyFVW+hOe1q/3F9kuFatW+TNi16/3B2MldxrJi8TYT34C48eXVsADdO3qryb1u9+F\nrqS8aSQvUsRKvbf99NPwrW/5YwnFfsC4mGkkLxJT99wD111XmgEPkEjAkUf6pZAlDI3kRYrU++9D\nx46wZk1pLxMwdSr87Gfw/PO6ilSuNJIXiaEHHoCRI0s74AFGjICdO+HZZ0NXUp40khcpQnv2wKmn\n+mA888zQ1TTcgw/6ls3cuaErKU15Hcmb2UNmtt3MKmt5/kgzqzCzVWa22sy+lkshIvKZSZOgd+94\nBDzA1Vf7q0etXh26kvKTcSRvZr2B3cCjzrlz0jx/G3Ckc+42M/sc8CbQxjl3yDp0GsmLZHbgAJxx\nBjz+OPToEbqa6Pz4x7B2LTz6aOhKSk9DRvIZLzfgnFtmZh3q2gQ4InX/CGBHuoAXkexMmQLt28cr\n4MFPpfz85/0FwNu3D11N+YjiwOv9QGcz2wq8CoyP4DVFypJz/uSnW28NXUn0jjoKrr3WTwuVwoki\n5PsDrzjn2gLdgAfMrFUErytSdhYt8uvGDxwYupL8+M53fLtmx47QlZSPKK4OORb4MYBz7i0z2wB0\nAl5Mt/GECRM+vZ9IJEgkEhGUIBIPEyfCLbfEdz75SSfB8OHwy1/Cv/976GqKVzKZJJlMRvJaWU2h\nNLNTgNnOuS5pnnsAeNc5959m1gYf7v/gnNuZZlsdeBWpxcqVcNllsG4dNGkSupr8WbvWnwm7YQO0\naBG6mtKQ7ymUk4HngI5mtsnMxprZODP7ZmqTO4ELU1MsnwRuSRfwIlK3hx7yByfjHPDgrxrVtSvM\nnh26kvKgk6FEisCBA/7M1hUr/ElQcTdpkg/56dNDV1IatKyBSIl76ikf7uUQ8OD78kuWwK5doSuJ\nP4W8SBF44gkYPTp0FYVz9NHQpw9UVISuJP4U8iKB7dsHM2fC5ZeHrqSwRo3yJ35JfinkRQJbvBg6\ndYJ27UJXUljDhsEzz8CHH4auJN4U8iKBlVur5qAjj4R+/fy3GMkfhbxIQH/7m+9LX3ZZ6ErCGD3a\n/yMn+aOQFwlo0SI455zSvzBIroYMgeee0zIH+aSQFwmoXFs1B7VqBZdcAjNmhK4kvhTyIoHs3euv\nlHTppaErCUstm/xSyIsEMn8+nHcetGkTupKwBg3y6/a8917oSuJJIS8SyBNP+Lni5a5FC7+08rRp\noSuJJ4W8SAB79sCCBTByZOhKioNaNvmjkBcJYM4c6NkTPve50JUUhwEDYNUq2LYtdCXxo5AXCWDK\nlPKeVVNT8+YwdKhaNvmgkBcpsI8+8ksZDB8eupLiMmqUWjb5oJAXKbCKCvjHf/QrMcpnLrkEXn8d\n/vKX0JXEi0JepMDK/QSo2jRt6hct+8MfQlcSLwp5kQL68EO/8uKwYaErKU6jR2v54agp5EUKaNYs\n6NvXr8Aoh+rXz1/I/O23Q1cSHwp5kQJSq6ZuTZrAiBEazUdJIS9SIDt2wPLlfqqg1E4nRkVLIS9S\nIDNnwpe/7FdelNpddBFs3gxvvRW6knhQyIsUiE6Ayk7jxn5lTs2yiYZCXqQAPvrIXxxj4MDQlZSG\n4cNh9uzQVcSDQl6kAJYsgR491KrJVp8+sGaNrhgVBYW8SAHMmweDB4euonQ0bw6JhL88ojSMQl4k\nz5zzIT9oUOhKSsugQf7KWdIwCnmRPKus9CPTM84IXUlpGTjQr7lfVRW6ktKmkBfJs4OtGrPQlZSW\n9u2hbVt/aUDJXcaQN7OHzGy7mVXWsU3CzF4xszVm9nS0JYqUtrlz1arJlVo2DZfNSH4S0L+2J82s\nNfAAMMQ5dzZweUS1iZS8nTt9u+aii0JXUpoGDfLfhCR3GUPeObcM+KCOTb4CTHPObUlt/35EtYmU\nvIUL/SyR5s1DV1KaLrwQNmzQZQEbIoqefEfgGDN72sxWmtlXI3hNkVjQrJqGadzYLwUxf37oSkpX\n44he44tAP6Al8LyZPe+cW5du4wkTJnx6P5FIkEgkIihBpPhUVfnZIf/936ErKW2DBvmzX6+9NnQl\nhZNMJkkmk5G8ljnnMm9k1gGY7Zw7J81ztwLNnXP/mXr8P8B859whl+Q1M5fN+4nEwQsvwDe+AatX\nh66ktL37LnTs6P/btGnoasIwM5xzOc3PyrZdY6lbOrOA3mZ2mJm1AC4A1uZSjEicqFUTjeOP9yG/\nfHnoSkpTNlMoJwPPAR3NbJOZjTWzcWb2TQDn3BvAQqASeAF40Dn3ej6LFikFCvnoaJZN7rJq10T2\nZmrXSJnYtg06d/YthiZNQldT+lauhGuugdfLdPhYiHaNiNTDggV+VogCPhrnnutXpNywIXQlpUch\nL5IHWnUyWo0a+bVsNJWy/hTyIhHbvx8WL4YBA0JXEi9a4iA3CnmRiC1fDqefDm3ahK4kXi65BJ59\nFvbuDV1JaVHIi0RMrZr8OOoo6NYNIjpHqGwo5EUiplUn80ctm/pTyItEaONGeO89OO+80JXE08GQ\n10zs7CnkRSI0f76fBdJIf1l5cfbZfk2gN98MXUnp0K4oEiG1avLLTC2b+lLIi0Rk715YutTPApH8\n0RIH9aOQF4lIMgldu8LRR4euJN6+9CX44x9h167QlZQGhbxIRLQgWWG0bOmvGLV4cehKSoNCXiQC\nzinkC0ktm+wp5EUi8OabsG8fdOkSupLyMHiwD3lNpcxMIS8SgYOjeMtpMVipr9NPh1atYNWq0JUU\nP4W8SATUqik8tWyyo5AXaaCPPoIVK/ysDymcwYM1Xz4bCnmRBlq8GHr29O0DKZw+fWDNGnj//dCV\nFDeFvEgDqVUTRrNm0LcvLFoUupLippAXaYCDUye1tHAYatlkppAXaYBXX4UWLeCMM0JXUp4GDoSF\nC/2iZZKeQl6kAdSqCatdO2jb1i9zIOkp5EUaYO5ctWpCU8umbgp5kRzt2AGrV/tZHhKO5svXTSEv\nkqNFiyCRgObNQ1dS3nr29Ffk2ro1dCXFSSEvkiO1aopD48Z+Df/580NXUpwU8iI5qKqCBQv87A4J\nTy2b2inkRXKwciWceCK0bx+6EgEYMACWLPErgcrfyxjyZvaQmW03s8oM23U3s/1mNjK68kSKk1o1\nxeX44+HMM2HZstCVFJ9sRvKTgP51bWBmjYC7gIVRFCVS7DQ/vvioZZNexpB3zi0DPsiw2beBqcC7\nURQlUsy2bYMNG/ysDikeCvn0GtyTN7O2wHDn3K8AXTJBYm/+fPjyl6FJk9CVSHXnnuvPXdiwIXQl\nxaVxBK/xC+DWao/rDPoJEyZ8ej+RSJBIJCIoQaRw5s2DoUNDVyE1NWrkZzvNmwc33BC6moZJJpMk\nk8lIXstcFhdJNLMOwGzn3Dlpnlt/8C7wOWAP8E3nXEWabV027ydSrPbvh+OO89d0bdMmdDVS05Qp\n8L//G79lDswM51xOnZJs2zVGLSN059xpqdup+L78P6cLeJE4WLYMOnZUwBerSy6BZ5+Fjz8OXUnx\nyGYK5WTgOaCjmW0ys7FmNs7Mvplmcw3TJdY0q6a4HXUUdOsGEXU6YiGrdk1kb6Z2jZS4s86CSZPg\n/PNDVyK1mTgRNm+G++8PXUl0CtGuESl7Gzf664med17oSqQuB5ce1njSU8iLZGnePH/6fCP91RS1\ns87yawu98UboSoqDdleRLKkfXxrMdGJUdQp5kSzs3QtLl/rZG1L8dLWozyjkRbKQTELXrnD00aEr\nkWz06+dXCt21K3Ql4SnkRbIwe7ZWnSwlLVtC795+zf9yp5AXyeDAAZg2DS6/PHQlUh+XXebPgC13\nCnmRDJJJf3GQ004LXYnUx4gR8OST8NFHoSsJSyEvksGUKTBqVOgqpL6OOQZ69YI5c0JXEpZCXqQO\n+/fDjBkK+VI1ejQ88UToKsJSyIvUYckSOP106NAhdCWSi+HD4emn4a9/DV1JOAp5kTpMmeJHg1Ka\nWreGRAIqynhdXIW8SC327YNZs/wsDSldo0aVd8tGIS9Si0WLoHNnOPnk0JVIQ/zTP/k15j/IdKXq\nmFLIi9RCrZp4OOIIuPhimDkzdCVhKORF0vjb3/xZrpdeGroSiUI5t2wU8iJpLFjg16o58cTQlUgU\nhgyB55/31wMoNwp5kTSeeEKtmjhp2dJfC2D69NCVFJ5CXqSGjz+G+fNh5MjQlUiURo0qz7VsFPIi\nNcybB927w/HHh65EojRoELz4ImzfHrqSwlLIi9SgVk08HX64Xy562rTQlRSWQl6kmt27/fz4ESNC\nVyL5MHp0+bVsFPIi1cyZAxdeCMceG7oSyYf+/aGyErZuDV1J4SjkRapRqybemjWDoUNh6tTQlRSO\nQl4kZdcueOopv3KhxFe5tWwU8iIpFRXQpw8cdVToSiSfLr4Y1q6FzZtDV1IYCnmRFLVqykPTpv7b\n2h/+ELqSwlDIi+BXKFy61K9YKPFXTleMyhjyZvaQmW03s8panv+Kmb2aui0zsy7RlymSX7NmQb9+\ncOSRoSuRQujbF9avh40bQ1eSf9mM5CcB/et4fj3Qxzn3D8CdwG+jKEykkNSqKS9NmvhlK8rhAKw5\n5zJvZNYBmO2cOyfDdkcBq51z7Wp53mXzfiKF9NZbcMEF8PbbfiErKQ/Ll8M118Cbb8Jhh4Wupm5m\nhnPOcvndqHvy1wHzI35Nkby65x4YN04BX24uvNCvTzRjRuhK8qtxVC9kZn2BsUDvurabMGHCp/cT\niQSJRCKqEkTq7b334PHH/ZQ6KS9mcMst8KMf+YvDWE7j5PxIJpMkk8lIXiuSdo2ZnQNMAwY4596q\n43XUrpGi8h//Ae++C7/+dehKJIRPPvHX8f3Vr/zB2GJViHaNpW7p3rw9PuC/WlfAixSb3bv9H/fN\nN4euREJp1Aj+7d/gJz8JXUn+ZBzJm9lkIAEcC2wHbgeaAs4596CZ/RYYCbyN/4dgv3Pu/FpeSyN5\nKRr33efnxpfTOiZyqP/7Pzj1VH/Jx3PqnFoSTkNG8lm1a6KikJdisX8/nHGGn0J3ftohiZSTiRNh\nzRp47LHQlaTXkJCP7MCrSCmZMsWP3hTwAnD99XDaaX4abYcOoauJlpY1kLLjnO/B3nJL6EqkWLRu\nDV//Ovz856EriZ5CXsrOokV+VsWAAaErkWIyfjw8+ijs2BG6kmgp5KXsTJzoR/HFNC9awjvpJL86\n5S9/GbqSaOnAq5SVlSv9iS9vveXXLxGpbu1aSCT8wmWHHx66ms8U07IGIkXtpz+Fm25SwEt6X/gC\n9OgBjzwSupLoaCQvZWPdOujZEzZsgFatQlcjxWr5crj6avjTn4pn4TKN5EWycPfdfiEyBbzUpVcv\nOOEEmDYtdCXR0EheysL27dCpk19W9vjjQ1cjxa6iAv7rv/wxnGI4QK+RvEgG998PV1yhgJfsDBkC\ne/bA00+HrqThNJKX2Nu925/d+vzzcPrpoauRUvHww/6KYQsXhq5EI3mROt12GwwcqICX+rnySj/V\nttR781q7RmJt+nSYOxdefjl0JVJqmjXzF5QZPBi++EX/bbAUqV0jsbVxo1+AbM4cLUQmubvnHr+g\n3bPPhju/QksNi9Swfz/06ePPbv3ud0NXI6XMORg61F9BKtTFRRTyIjV873tQWelH8Y105Eka6P33\noVs3ePBBf3yn0BTyItUsWADf+Ibvwx93XOhqJC6WLoVRo+Cll/xiZoWkkBdJ2boVzj0Xfv97uOii\n0NVI3NxxBzz1FCxeXNglDzSFUgSoqoKrrvJX+VHASz58//u+/XfnnaEryZ5G8hIboUZZUl62bfNT\nKh9/3C9LXAhq10jZe+YZv2zBSy9B27ahq5G4W7jQXy7wlVcKc9xHIS9l7eDMh9/+Vpf0k8Ip5Awu\n9eSlbG3c6M9IHDNGAS+Fdccd8OGHcN11fjGzYqWQl5I1eTJ07+6ntd11V+hqpNw0aeKn6x444Gd0\nFevSGWrXSMnZtQtuuMGv9T15sj8IJhLS44/D+PH+AvE33RR9+0btGikbL7zg++8tWviDrAp4KQZj\nxsAf/wgzZkD//v58jWKhkJeSUFXle6DDhsHPfga/+Q20bBm6KpHPnHKKn+XVu7cffFRUhK7IU7tG\nit6mTf4kp8aN4bHHCn9KuUh9LV/u99kBA/y1hVu0aNjr5bVdY2YPmdl2M6usY5v7zOzPZrbKzLrm\nUohIdXv2+LXgr7rKt2cGD4Ynn1TAS2no1QtWrfLHj848E/71X/1SxVVVha8l40jezHoDu4FHnXPn\npHl+IHCjc26wmV0A3Ouc61HLa2kkL7X64AM/53j6dH/m6vnnw4gRMHy4TnCS0rV6te/VT5/uz5Yd\nNgxGjoR+/aBp0+xeI+8nQ5lZB2B2LSH/a+Bp59wTqcdrgYRzbnuabRXyAvg1unfu9K2YFSv8H8CK\nFdC3r/8DGDIEjjkmdJUi0Vq/3gf+jBnw2mswaJAP/bPPhvbtoVWr9L/XkJCP4vJ/JwGbqz3ekvrZ\nISEv8eMc7NsHe/emv338Mbzzjg/zmrdmzfyO3aULjBvnd3wdTJU4O+00uPlmf3vnHZg1yx9nWrfu\n7/8mat4aouDXeNXX7sLJ9KWp+vO13f/kk/S3qqrP/tu4MRx+eO23E0/0O2r37v5KTe3bQ7t2cMQR\n0f7/ipSSE07wg5tx4/zj6t9uq99efLFh7xNFyG8B2lV7fHLqZ2mNGTPh0/s9eya48MJEBCVIbSzD\nF7zqz6e7b+ZXdGzU6NDbwZ9rxUeRhjODY4/1t7/+NcmOHUlatvSXHWzQ62bZkz8F35Pvkua5QcAN\nqQOvPYBf6MCriEh08tqTN7PJQAI41sw2AbcDTQHnnHvQOTfPzAaZ2TpgDzA2l0JERCR6OhlKRKTI\nae0aERFJSyEvIhJjCnkRkRhTyIuIxJhCXkQkxhTyIiIxppAXEYkxhbyISIwp5EVEYkwhLyISYwp5\nEZEYU8iLiMSYQl5EJMYU8iIiMaaQFxGJMYW8iEiMKeRFRGJMIS8iEmMKeRGRGFPIi4jEmEJeRCTG\nFPIiIjGmkBcRiTGFvIhIjCnkRURiTCEvIhJjCnkRkRhTyIuIxFhWIW9mA8zsDTP7k5ndmub5I82s\nwsxWmdlqM/ta5JWKiEi9ZQx5M2sE3A/0B84CxphZpxqb3QC85pzrCvQF7jazxlEXK38vmUyGLiFW\n9HlGR59l8chmJH8+8Gfn3NvOuf3A74FhNbZxwBGp+0cAO5xzB6IrU9LRH1K09HlGR59l8cgm5E8C\nNld7/JfUz6q7H+hsZluBV4Hx0ZQnIiINEdWB1/7AK865tkA34AEzaxXRa4uISI7MOVf3BmY9gAnO\nuQGpx98DnHNuYrVt5gA/ds4tTz1eAtzqnHuxxmvV/WYiIpKWc85y+b1sDo6uBE43sw7ANuAKYEyN\nbd4GLgaWm1kboCOwPqoiRUQkNxlD3jlXZWY3Aovw7Z2HnHNrzWycf9o9CNwJPGJmlalfu8U5tzNv\nVYuISFYytmtERKR05eWM10wnT6W2uc/M/pw6gaprPuqIiyxORrvIzD40s5dTtx+GqLMUmNlDZra9\n2rfOdNto38xCps9S+2X9mNnJZvaUmb2WOqn0X2rZrn77p3Mu0hv+H451QAegCbAK6FRjm4HA3NT9\nC4AXoq4jLrcsP8+LgIrQtZbCDegNdAUqa3le+2Z0n6X2y/p9nicAXVP3WwFvRpGd+RjJZ3Py1DDg\nUQDn3AqgdeqArRwqm88TQAe1s+CcWwZ8UMcm2jezlMVnCdovs+ace8c5typ1fzewlkPPSar3/pmP\nkM/m5Kma22xJs4142XyeAD1TX9/mmlnnwpQWS9o3o6X9Mgdmdgr+W9KKGk/Ve//U+jLx8BLQ3jn3\nsZkNBGbip7GKhKT9MgepE0mnAuNTI/oGycdIfgvQvtrjk1M/q7lNuwzbiJfx83TO7XbOfZy6Px9o\nYmbHFK7EWNG+GRHtl/WXWthxKvCYc25Wmk3qvX/mI+Q/PXnKzJriT56qqLFNBXA1fHpG7YfOue15\nqCUOMn6e1XtyZnY+fmqszlOonVF7r1j7Zv3U+llqv8zJw8Drzrl7a3m+3vtn5O0al8XJU865eWY2\nyMzWAXuAsVHXERfZfJ7AZWb2LWA/sBcYHa7i4mZmk4EEcKyZbQJuB5qifbPeMn2WaL+sFzPrBVwJ\nrDazV/Cr+34fP7Mu5/1TJ0OJiMSYLv8nIhJjCnkRkRhTyIuIxJhCXkQkxhTyIiIxppAXEYkxhbyI\nSIwp5EVEYuz/AYvongtzVccXAAAAAElFTkSuQmCC\n", "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "pyplot.plot(numpy.linspace(0, 2, nx), u);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "OK! So our hat function has definitely moved to the right, but it's no longer a hat. **What's going on?**" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Learn More\n", "-----\n", "***" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "For a more thorough explanation of the finite-difference method, including topics like the truncation error, order of convergence and other details, watch **Video Lessons 2 and 3** by Prof. Barba on YouTube." ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQAAAQABAAD/2wCEABALDA4MChAODQ4SERATGCgaGBYWGDEjJR0oOjM9PDkz\nODdASFxOQERXRTc4UG1RV19iZ2hnPk1xeXBkeFxlZ2MBERISGBUYLxoaL2NCOEJjY2NjY2NjY2Nj\nY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY//AABEIAWgB4AMBIgACEQED\nEQH/xAAaAAEAAwEBAQAAAAAAAAAAAAAAAQMEAgUG/8QAPRAAAgICAAMDCAcHBQEBAQAAAAECAwQR\nEiExE0FRBRQiMjNhcbFyc4GRocHRI0JiY5Ky4RU0UvDxQ1MG/8QAFgEBAQEAAAAAAAAAAAAAAAAA\nAAEC/8QAHBEBAQEAAwEBAQAAAAAAAAAAAAEREiExYQJB/9oADAMBAAIRAxEAPwD78AAAAABBIAAA\nAAAAAAAAAQSAAAAAAACAJAAAAAAAAAAAAAACAJBVC+M4TnFPgjvT/wCWvAivIhZbKtcpRSk0/B/+\nAXAquuVNfaSW4L1mu5eJYntbQEgAAAAAAAAg4ttjVW5y6L8fACwEb5c+QAkghSi5OKfpLqjjIujj\n0TtkpS4Vvhitt+5AWgrpthdTC2uScJpST9wnbGE4Rf770n7wLACAJBxZZGuuU5PlFbZKfopvltdG\nB0Cu62NNfHL1V1fh7ywAAAAAAGfEzKMytzx5qajJweu5p6aNB5XkPyLX5J85mpuduTbKyb7lt8kg\nO8nMya7r3B09jQk5Rlvily29Gt5dMbIVynqc1tRZkv8AJ3bZFt/BDtYyjOqT9y6Mqn5PsnkXWTjc\n1ZOM+GNqS6Lk/hoD0PO6POVj8f7V9xEc2iyyyque7IJ7WivGhfRxVurcXbJqfEuje9/jozVY2U8q\nq25Tb4Zce7FwptdyA1YmVK+zhlFL9lCfLxe/0KpeVK5edRr9fHkova5d36neFj2U2cU0kuxhDr3r\ne/mc2Y90o5kOBNW2KcXtc+UV+TAvpzKLo2SrsTVb1LfLR3RkVZMOOqXEk9GTOxLMh38HSVcEtS1t\nqTevcd4OM6LLZtW7ml6Vlik2Ay/KFVE1VGSdrlFcOvGSRcsuh3SqVi44rbRjlRkql0KraeQp8fEv\nV4lLf5HGL5OnTdW59tLgtlJPtVw8989faBfgZrzK67FOvUpSWtNPl0LoZ2PbN112xc1HiXhoqqxr\nYvEckv2bnxc+m96Kq8K6ONg18KTpi1Nb6eg182BdHyhX5zVjSads6uPcd8LO5ZtNNVUr7I7musea\nf+CrHxrqpY/FHfDjKuWn0fI4xqMnF4ZKpWPsIwa4ktSW+X27/ADVfmY+Ot22JdH48i3jiq+0clwa\n3v3GLEw7aJx4tNRxo1796b/U6eNb/pNVGk7K4Q9HfJuOnr8NARk+U6YYFuTVNPh5LafX3o1dtGGP\n205ejw7b0ZLqL753T7PgcoRilxLnz2astXOlqj1m1vT563z179ATRkVZNfaUy4o9DNd5Rqjl041c\nlKydnDJeC03+Q8m41tHnPaKa47Nxc58Ta0lv8CqujJUcOp0pRps4pT4lzSTW/t2Bshl0TtnXGxOU\nOqOfOq5uLrsg4qWpN9/LfI8/F8mW0xqU1dOVanFOVqcejSeveaoYs40YEFFLsNKaT6eg18wNburS\ni+JeknJfAy2+VcaGHPJjJyhCSi+T3zZXXjZNMbnCMXKuLhRzXOO9/Z3L7CqODkJ5XKxqcYKDssUm\n9SbfwA3xyYOTk5x4GoOPj6W+p3ZkVVKTnNLgW2Z54srrclzXDG2qEU9801xfLaM1mDffj1zsUo39\nsrJquenrTjyfw5gejXfVYoOEk1P1S0x4VM8aPAoT1KcnKU5qT+JsAGXylbKnBslF6k9Qi/Byaivm\najD5X5YcX3K6pv8AriB3k40lgxox+UY8KaT6xTW19q5GezCnKVl/Dqc516inzUY93x6/eejYtwe2\n4657XUog4zsdatt4lFS09dH/AOARh1TjjShek5SlJvnve2c+TG1jSql1pslWt+CfL8NGqEeFa4pS\n+JkwOduZJdHe9fYkn8gNoBmyMaV01KNsoaWtKTQGkrttjVwcW/Sko/ayMep018Mpub31b2VZcHd2\nca9ScLYykt9EBfCyM48Uem9c0dHlxxMnlxQ5JXL1l3y3E0YdV1NtjubVfZV63LfpJPi/IC/Krnbi\n2wraU5Rai34mB49tcJzmnwxsraW9vhjr/LPTlqceTa33o47H+bZ96/QDLdjSvzYTlWnXFNet1bcd\nP8C7GsldVY23F8cor3dxZ2P82z71+g7H+bZ96/QDHTg3QzJ2Sy7nHUeuvS1v3FuXO6N1fBju6KTf\nKSWmX9j/ADbPvX6Dsf5tn3r9APMpoy/M8iiNHZxlKTpUp74E1+rZplRGFdFVcOFytjN899Obb+4s\ndlSy44rvt7WUXJLu0vfrrzLux/m2fev0A5yMedzThkWVa7oa5kwplCh1ytlbLn6U/wDBPY/zbPvX\n6E9l/Ns+9AeXHEvjjz40/Rqriue98L2/v5mnJoeTfW3WpVc9vi9ba0vm/wADV2T/AP1s+9foOx/m\n2fev0AyXzlbiZlaT3zqgvHlr5m9ckV9l/Ns+9foWJaSW2/iBIAAAAAAAAAAAAACu62NNM7JdIrb0\nV02T1u+ytS74J+qBeSVTvqrnGE5qMpJtJnbsgk25LSW+oEgzRvSzJxlNdm64yjt8ur/wWqUKKZOV\nm1Dbk2/tCasOa7IWR4oSUltra8URVfVbCM4TTUuhT5P0sXXepyT+O2BqKsjtOyap0ptpbfd7y0xZ\n+csNxcluHC5S8eTS/MQtae1irlV+848SLDy7MpR7HJuXA67HXNLn1X/hbDyiuK5WVTg6+aXDt8Pi\nXDW9AwV5XBhLOutbrlFS4UuSTLK81ydqnROuVceLUtekhhrUyITjZBSg9xfNM8yzOuq81lYlrJe2\ntcoR/Xn+BzK/Iwo1VtJqM+CEIrnNL/1DE5PQd/DkyrktRjXxuX2v9C2EozipRaaa2mjzrsvWRjXV\nVStd1bXAtb7n3/aU5DlX2VtMrKYWuSlDuhPrv8GXicnsg8qOTl35EuzTj2MVPg16+0uXzPSU918a\ni3y3rvJmLrsGDz3Jb1/p9y9/FH9RiZd7yFjZFMlJxclPl0Xj7xlNbwDz/KGWuxuqqU5TitS1y13/\nAC8CTstx6BTl0+cYtlW9OUdJ+D7mZrcyalZ2HA401xskn3p/4R1/qMFJ8UXGtScePa6630LlNiJz\nuyfJ3obhdyjZFdVz9JL363opsx7n5xdVXJTsgoVxctcKSfN/a/ka8LIjlQlbGpw29c+80E8WXWaH\nDh4spaly56b22+gxYeaYW7PXe7J6/wCTe38y/dc5a3Fyi9630M1cfOnernxV8bioeGgjRVarYKST\nW+aT6nZVGiuFkZRWuGLil7mRGyTzp1t+iq4yS+1/ogq59HrqYY0ZEMjtoxinOLjZ6XV9U/s6fBm5\nNNvTXLqRKcYKTlJJRW2BRXDI4Gp2alttct8tv/Bxdj5E/wD7cUeKL4eHW9S2/wADTGyEoxkprUvV\ne+pzK+ClKCfFYoOfCurQHOJVKnHUJvb4pP4JybS+5l55mB5TruoipOcrm3uPA99XyNtFlticrKuy\nXdFvbLiauABFCH0JAHj5SvspjOvDt85hPtIybWt+G99NcjTOmd2XGUoTjWk+XH1e09/gbgAJIJAA\nAAQSAAAAAAAAAAAAAACnKrlbjyjH1uTW/FczzM3DeRGFlmJwNTbsUWm3HXX/AL4HsEllxLNeXPFs\nyrW0kseyMG+LlJcO+X4lEPItkapzdjds6uCS316d/wAUe0SXlWeEecsCORZCd9EYxhF1qttSWu5/\nE4h5NtqhJxu45ycZSUuknFrh+HJaPTBOVXjHjryfl2eUpZMuzqXKS1z596/z8D0cPHljwmpz45Tm\n5t61zei8kW2k/MgcWVV2644KWvFHYI047OHpeivSe5e861y0SAMGfTCjyRbCCfDXHiS+D2aKMeuq\nDUW58XWUntyLZRU4uMluLWmhCKhCMI9IrSKmJSSSSWkuSIcU5KTS2uj8DoEVXGiqLTjBJp8XLxJs\nqhbHhsipR3vTR2AIJAAgqpp7Oy2cnuU5b37u5FwAHm5t1N1rolOmCg1uyUucX15L4d56R5V12DTn\n5Dyex5wi/Sim+/f5Fnxm/SUfJ8rG45VSjKKjOKkuaXRFSjh3ZuVxZFfZy09cS6uOmzbjRxMmMmsW\nEHHXKVaT5rZbHDoja7I1QW4paUVrlv8AUbYcZVGLaqo6tyq7EvRgq460jVVfC5yUVJcPXcWhN10V\nysaUYxW3pFUs+pVRsalwyrlZ05pLW/mPV8V4+D2OXK70VtybaXOW/H4HGNXkPKyIStjCHFxcEVz5\nrx+8vvynGUaqFGd81xRi3pa8WUTzo0dlK2EVZZZ2Vmn6uk3/AN+Je6z09BdDDm1cWZQ+KcYzi4+g\n9Ntc1z+xl+HlQzKFdWpKDbS4lrZdOEZ8PEt8L2viTxr14nYXqtzhTKCmpcSXWU2uXLwWtfaW+bWZ\nNk8xWdip65WLui1ra+89WyyNUHOb0jNmzjOicYtS4Jx7RLuW03v7DXJniyPHuflGqM0rlRFTilqO\nm9r8juGM8HOV74rI2cUfRW+Hb2l8+Z1j5dd+ZffD2NdfC7O5tNs1Tt7XtK6JanFetraT8BbSSK+z\nvqscaVF1zlxNt+p4/ebCrHt7alTa1LpJeDXUtMNQAAUAAAAAAAAAAAAAAAAAOZzjCLlNpJd7A6AA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAADLdiq69tpcE6nCXj/AN6moFlxLNUY2OseDXHKcn60\npdWXgEVXfUrqZ1vlxLW/AxvAsvrSybUmuWquS4e9faegCy4lmstGFXTbGzcpTjDs05d0d70cZfk2\nnKujZNyWt7SfJvWtm0DaZGNVZcaOzjKqL6ca3yXw8TYARVGZXK3GlGCTktSjvxT2ZsSq6jJnKdfE\nsh8c2n6r/wDNL7D0AXUsYMaEo4nmrq1OMXHbXov3nWJjX4rVcZ1yp6tvfG2bSJJuLUXp65MWmYzY\nnt8tr1e0WvjwrZqK6KlTWoRbfe2+rb6lhCAACgAAAAAAAAAAAAAAAB5vlq914jrhFynLT0u5bR6J\n5/lJQpnDKlCc1BPiSkkuXNbL+fWf144xcyVmfVXKT9W3fufEtfgmeoeV2u8mmOROuNkpqVari+fJ\n9W/ds9UVYAAigAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAB5ebjzzcmVFlbcOWm/VUe9/F9D1CCy4lmsEsGyzGq4mu0qr1Bdyly5/ga6Z2Ti3bX2b3\nyW98ibao3QcJ717nooWLZSo+b3SUY/uT9JS+3qh6ni+yzglXHW+OXD8OW/yOzDPIfbY8L4Ouamuf\n7ren0ZuI0kAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAABmzIRslRCaTi5tNP6LM9uTHydJwlLjr16Md+lF9y+BpyPa4/1n5MjMwKM1Vq+LfZy4\notPT2WZ/UvxpBxC2ucpRhOMpQ5SSe9HZFAQSAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAABnyPa4/1n5MtnZwOPoye3raW9FWR7XH+s/JmgDwa/Jqr\njOeBkJylPilPfPae9fbzREM/yjTl3K2nhq3w1dq9Lq+e/gvxR7VmPVbrjgnp7XcVzotSm67eJvpG\n1bSA5w8qeROxSr4Yw4VxJ7TbW2vs2azLCTonPix3FP0p2Q5pvXh1Lar67Ypwmnvp4gWggASAAAAA\nAAAAAAAAAEASAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAz5Htcf6z8mX\nlGR7XH+s/Jl4EgAAY82/EplF3R47f3YxjxT8NpDKypqxY2NFTyJLb30rXi/yXed4uFDGTk27Lpev\nbL1pfovcB5dOT5XlCTx8BJKxr9tenvT19hq/1S/GpdnlHBsoSfWt9qvw5mvD9nP62f8Acy6a3CW+\na0BRh5+JnVqeLkV2xa36LNB5NXkrHuwqLaUsfI7JJW1x0+a713neNO6N0caybpsjtqL9KNq8U3z8\nN/ED1AZu3urS7Whtt6/ZekWV31WSlGFico9V4AWgEASAAABAAwTsnnZHZUTcMeuX7S2P7zX7q/N/\nYLrJ51ssXHk41R5XWr+2Pv8AF9xtqqhTXGuuKjCK0kgOwAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAGfI9rj/AFn5MvKMj2uP9Z+TLwJMubkuiMYVR477Xw1x9/i/cizJyIY1\nErbHyXRLq34L3lOFRZxSyslft7V6u+Vce6K/MCzExljQe3x2TfFZN9ZM0PoQAKMP2c/rZ/3Mun6k\nvgU4fs5/Wz/uZdP1H8AKcD/YUfVr5HPlDHjk4lkZNxlFOUJx6wkujR1gf7DH+rXyLLvY2fRfyAw+\nTcm2Ma8XNlvIUNqzutXivf4o22U12xcZwTTKp41eVh112bXJOMovTi9dU/Eppy7ceao8oNJt6hcu\nUZ/HwYGmNLhYnG2aglrgfNFwJAAHNlkKoSnZJRhFbcm9JAdHn2XTzrJUYsuGmL1bcv7Y+/pz95y5\nXeU21DjpxO+fSVnw8F7+8311wqrjXXFRhFaUV0SAiimuiqNVUVGEVpJFgAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAABnyPa4/1n5MsutroqlZbJQhFbbfcUZ1sKOyttkowj\nNtt/RZVTVPOtjk5EXGmL3VS/7pe/wXvAnHpsyro5mTFwS9jTJer/ABP+L5G4kACCQBnw/Zz+tn/c\ny6fqS+BTh+zn9bP+5l0/Ul8GBTgf7DH+rXyLLvYWfRfyK8D/AGGP9WvkWX+ws+i/kAo9hX9FE2Vw\ntg4WQU4vqmtoij2Ff0V8iwDzuDJ8neyU8nGX/wA+s4L3f8jv/V8DmpZVcJR9aMnpx+KNxw6629uE\nW/egMb8qQtUfMq55PHyU4L0E/exDCsyJK3Pkptc1TH1Ivv8ApfabVFRWopJe4kAuRIAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAOZS4YuT7lvkdADyXCzMyMbIyYuFSs3V\nS/g2pS9/TS7j1jPke1x/rPyZoAAAAAAM+H7Oz62f9zLbfZT+iyrD9nZ9bP8AuZbb7Kf0WBT5O5+T\n8f6tfItv9hZ9F/IzYVtdXk/FVk1Higkt95pu9hZ9F/IBR7Cv6KLCuj2Ff0UWAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA5nJQi5Pel4I6AGLIyK3ZR63r/APF+\nDL/Oa/4v6Wc5Htcf6z8mXgVec1/xf0sjzmv+L+llwAq85r/i/pZXkZirx7J1pucYtpSi0jUZL08j\nIVH/AM4alZ7/AAX5liV53/8APeUb8im/zurg4bHwyjF+k223+JvzM2NdD7NNzl6MU4vWzXCMYLUU\nkuvIlpPW1vT2LZpJcfN13R8teQVVVOePdTKKjOMdtOOmj1rcmUfJ0uFSncquji+b0d+S4Rhg18MU\nuJbeu9+JrLc1J528byBmX2Ycnmce+LUPR7j0/Oa/4v6GVf7XI/k3P+mf6P5/E1kvdPz1FPnNf8X9\nDHnNfjL+llwI0p85r/i/pZZXZGxbjvXvWjsAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAABy4qTTa6PaDek2+46K77OxpnZwuTitpLvA4xsmOSpTr51p6UvF95eYMKFl\nWNVTZuucJek0uUuf5l1+ZXTxJbnNfuR6suJK0nMYRi5OK05Pb94rnGyuM484yW0dEUIJKsiuFlMo\n2NqHV6egFMa6oqmD9RdO/RYeRh2q3Dmp5Eq8mWnOzh00uSXX/vM9ddPEtSObK421yhNbjJaZjxst\nQkseyXHPi4Yy/wCS/VaezeZbcGmzKryUuG2D3xLv5a1+IhfibsyuDdcGpXbSUN97NB5uZOuHlKiy\nyD4KYvckt6cta+TLnkunAldKxT4JPbfet8/w+QxJW0FNV/aX2V8DSilJS7pJlxGgAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAABXdUrYpNtOMlJNeKOnH0WvE6IAqx\nK5VYlVc/WjFJlwAA5nCNkHCa3F8mvE6AFd9Mb6XXJtJ66deuya4dnBR4pS13yfM7AAAAQVwohDj0\nuU3tp9C0AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAENpNLfN9CHKMVttJHNlXaST4taTX3lSxEozjx7Ulrmt6A0EOSWtvryRRLE3Pbtlrny+O/1LFQl\nGKTfoy4gLDPHOolbXXFy3ZtRfC9NrfLf2M0JaM1OGqrIyc3JQcmlro5Pe/yA1AAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAACCQABBIAAAD//2Q==\n", "text/html": [ "\n", " \n", " " ], "text/plain": [ "" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from IPython.display import YouTubeVideo\n", "YouTubeVideo('iz22_37mMkk')" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQAAAQABAAD/2wCEABALDA4MChAODQ4SERATGCgaGBYWGDEjJR0oOjM9PDkz\nODdASFxOQERXRTc4UG1RV19iZ2hnPk1xeXBkeFxlZ2MBERISGBUYLxoaL2NCOEJjY2NjY2NjY2Nj\nY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY//AABEIAWgB4AMBIgACEQED\nEQH/xAAaAAEAAwEBAQAAAAAAAAAAAAAAAQMEAgUG/8QAQxAAAgIBAwEEBwUFBgQHAQAAAAECAxEE\nEiExBRNBURQiMmFxgZFCUqGxwSMz0eHwFTRTYnOScoKDoiQ1Q0RUk/EG/8QAFwEBAQEBAAAAAAAA\nAAAAAAAAAAECA//EAB8RAQEBAAICAwEBAAAAAAAAAAABEQIhMUEDElEiE//aAAwDAQACEQMRAD8A\n+/AAAAAAQ2l1aGQJBzvis+svV689CVyBIAAgkEASCCFJNtJrK6ryA6BAAkEEgAQmmsp5AEgjOOoT\nUkmnlPxAkAiUlH2ml4cgSCABIIDaXVgSAQBIBAEggkAAQBIAAAEASCCQAAAAgkACABIAAAAAAAAA\nAAAAAAAAAg8L/wDm+19R2jdrqL6XjTXyhG7wkk+F8T3jiqqFMdtUFGOW8JeL6gYlp6tTrNXG6Cmv\nVXPwKq5qu5adx3wq1KhW2/ZzXn8MtG2elUrLJxtsg7MbtrXgRDRVQhVFbv2ct6bfLeGsvz6geSlr\nH2c/Wo3S1fk8Nd54l9na7hq3WnBxhNVyioy3N9G14f8A4egtJWqlX623vO86+O7P5keipTm4WTgp\ny3SSfiBms7Qtqr1F84xVNMnBdcyllJP4FNXa1sqbcqFk4OGJRi4xe6W3HPiehLSVSpsqe7bZJyfP\nKb54IlpFZFxttnNPb7ujygM61mood0NTGuVkVFw7vOHueEufeRfqtTotM7NW6XKUlGGxPGX5mjU0\naeTm7p7e+iq87sdG2se/kmWkVkMWW2S5TTzjbjyAr7N1j1lc92HKEsboppPjPGTJprdVUr0+6dl2\npnCt4eFjL5+SPUqrcMuVkpt+ZVPRVzjjMk1Z3qknypAZ/Tb665wsrhK+E4R9XiLUn1/ryFt+u3um\nvuO8hVGc5Szh5b4X0NC0cPWcpSlKU1Nyb8V0+XBZ3Ee9nZzunBQfwWf4sDHRq9VKVFlsalTqI5go\n53Re3PJbXq5So0djis6iST92Yt/oWR0tcYURWcULEefdjkrr7PhCyqXe2ONLzXDPC4x8+oGRaqzT\n9nQlCyqEnOX7xN59Z9EiLO1LnTpJ193X3y5lYntz0xx0+ZuWhhGUJVznFwUkmvKTyyI6GMNPGiFt\nigouL5Tyn+vIFfa/ePsq3unDc4rl9MCd2slY6aFTvqinY5Zw2/BfQt1KohpfR5tqGxpJdcJZEtOr\npd9Gc65WQUZbfFfx5Aqo1l2qugqoRjDZCyTl1w88fHgjXK6XaGhVbr2bpOSkn5eBa4UaCPeLKjiF\nXwSeF+ZfKqM7a7HndXnHzWAMFeptUaadNVCLsc3mTeI4kdVay+ztC2jdTFQziEk9z8n5M1V6SuqU\nJRzmCklz5vLKrY0q+Lu1DzCW+MW1w3lfqBz2Q9Q9DF6mUJSy8bc9Mszam/U6iqNsO6jpe+gsPO94\nmln6npUUqivZGTcctrPhl5KHoIPEe8sVas7xQT4znP0z4AZP7Yb1yrik6++7lx2vdnO3OemMll+p\ntso1Up1QdFe6GG3mTT/I1LSKNjcLbIxlPe4p8Z/gdS0lcqLanu22tuXPmBTZrJwd8Nq7yM4xgvvK\nXR/n9Ci7tGyntCFLnVKE5qG2Ke5Z45fTqbZ6Wuerr1Lzvgmlh8P4/j9TNZptNXOEbNQ4Lvu+jByS\nzLOfmsgZ6Z3rTaiWq7uyPpDUUs/fNctZNb61Bd8rlXGPmnzn6Z+h1HT1ShbFSk4ys3teUs5FVVOo\ntq10d2514is8Yf68v6gYqe2O91kYR2uuU3WoqL3LnGc9McGjQay7Uam2NkqlGOV3aTU4v356/I7t\n0/cUWyrtsilusUc8Z6/TPgTpKO8jRqJ2znJQ9XOOMpZA2AHDshGexzipYzjPOAOyCQBX31fe93vW\n/wAjs8/TqqzR2Svaj+1nLOcNYk+TTHVQlNKv1opJymnws9C4mtBxKyMZRi+suEcR1VEre7jbFzzj\nbnkzW6muev0yi84ck/y/NMYWtxJjlrktbDTxjlNuLlnGGlnp4l1V8bZzjDL2PbJ44yQ1cAAoAAAA\nAAAAAAAAAAAAAABD6EkAed2bN3SnZON6mnJOU36r5fQt0l853OMpOUJxc4ZWMJPBdRR3VMq89ZSf\n1bf6nGk0zoy5yUpYUY48IroVGkw9q6qWmrgoS2OecSxnouhvMeq7Or1d6sunOUFHb3ecITN7LudL\ntJYrNPB95G2SWJSj0bLirT0V6amNVUdsI9EWkI8vtWmV2poXGO7tSb6KWFh/mVK2zfLVtyjsqhLb\nnjHOUelqNP38obpyUFndBdJfE4t0ztt9ZpVcZilzLHg/ca3pmxqIfCYJMtvHjdOFtN1cbbbLW5Wx\nhziHgseHh+JKWo1MdG56qdfeNuccJcrnB6ddVdSargo5eXgW012xUZxyk8rw5LrOJ3wVnd7luazg\n6K6tPVTKUq4YcnlvOWy0iqtTa6NPOxLLiuEZLrbdLQ5T1CsnN8PakoLGX8Te8NYZR3Wlq421rndz\n5lM3ww23vUwr3LbbsnCcH1Tcf5F/Zmrnq4SnsxVFKMX5vxNLsobzurynnOUZXpqIy/Y6l01+NcJL\nDGw+vJm7VunfVqKYycVH1YqKy5y65+CPVoe6iD3b8xXreZXH0WHsutcY6+BNPo1NeyqUIx8twtiz\njfa8w19nx9IlZa9/OVnxfm/0Rr72v/Ej9SHdV/iR+pNMZdZOyd8Ka593GK3TmnyuuPxMa1GssnTO\nalCvfCTljCaeFj6t/gatY6LLIN0K7wclNJL48nT1eGoSpThjrGaaRqM2U0V9s4zu1EoquyaVUfJd\nFn4mi+V8Uu4hCT8d0sHmu1JJVRnOmM9+OFtw8+fQ1+n7Iud9M6ocYl7Sf0Fns7nVaq3NwTsSUvFL\noedOyL7UlGMYzk2oyjJZcePaXu6GhdoQtg3pq53NPDSWPzJjqZ7pOeksi/BrDyIW+mW2cpqENW5w\n9TO2HDsn8vh+JmhXqY6LSaa1OquOYTm88P7PxPQ9Mvc8Ps+3C8cx/iX9/Fr1q7F7nAu4z5YNVKyX\neaaM7E4QUYQS5m8e035fzNnZrb0FCfWMdr+K4f5HNmqtU/U0Vsl97KRxpLZ1qxWaeyEXNyT4a5+B\nLOll7bzz6+yaY6yWolKU5N59Z58fy9xrptlbuzXKCTwt32veWE7jWaB8ppPD8yQRXk0aDNl0VdKN\n0LMubWcqS8vqaKuza6dPfTXJ7LV0f2XjBtwiS6mMFfZdcLa7VJ95BRW7HlnP1yyyGgqjY5vMvW3R\nT+zznj5msDaZFFmlpskpSgtye7cuHn4meGf7QUqqrIKWe9cliL8vnk3nFtkKa5WWSUYRWW34DTHR\nJCakk10ayiSKAAAAAAAAAAAAAAAAAAAAAAAAFF2qqpbTzKS6xissvM9Fco33zkvbkse9YCzPZHUu\nxZhTNr38ETt1CXqaVyfvmkaMJeBzZZCqO6ySjFeLYLZ6Z43ayXtaWMf+pk6nHVyj6s64P4Nmjqsg\nEuXWJVa/POor+Vf8y1VXS/eXS/5cIvjKM1mLTWcEkxq87WSWhUuuo1PymQuz0v8A3Gp/+w2kFYZV\noor/ANW9/wDUZ0tJBf8AqXf72aQBn9Ep3bmm375Mieg0s3mVEH8jQcW2xq2Z+3JRQWW8e4o/s3Sf\n4EPoP7N0n/x4fQvouhfUrIey84+uDqE42LMJKSzjKJka/wBef6zf2bpc/uYfQsjpKYLEa0vkXgZC\n/Jzvmq+4r+6jiWlpl1gn8i5yS6tL4kQnGazCSkvNDE+1/VHoOm/wo/Qh9n6WXWmP0OdTro6aN0pw\nco1bc46vJpnbCuCnOaUW0k2/MuLPk5fryu19BVDs+x6elKzosLnnj9T1ao7aoR8kkdNJ9QTF5fJy\n5cZKkAq1FqopnZL7K6eb8iua0GbS3ynp27V+0g3GaXmv5YLNPfDUQc687VJx580BY2optvCXVme/\nUul1NQUq5yUXLd0z0x5lt8I2UThP2ZRafwPPpljTXJQWIQ76ncs7U08fR5Klr1AYtHqXdq9RDdlR\nUGvpyWOTj2lGOeJ1Pj4P+Yw1pAIIqSG0k23hIKSfRoy9pOfou2Fc7FKSUlBZe3xCVorthbBTrkpR\nfijJ2pa1o591Ke9SS/ZrLTKYWWft416W6Nc03tksYl7vcaezYzr03d2QlGyLxKT+2/MuJu9LdNdG\n+pSi3xw1JYafwMmolLU3z0l8P/DWpwXg20st/wBeR6GEm2l16nlzSn2lXKcn30LXtjnhV4fOP1EK\n16bUxtlGuMWv2an8nwvyNEnti35LJg0tHouvaV3eRshhLj1EnlL8X9DfNZg15oVYr0tkrdNVZPG6\nUU3guM3Z3/l+n/01+RpJSIJKb7u52Lhuc1FJstTz0ChRXq67NXZp453VrLfgXnlae2vT6i+ye6U7\nLWntXEFuws/MRK9YAr76vvu53LvMbtvuCu5SUYuUnhJZbEZKUVJdGsoz9ownZoL4VvEnB4OaI3uu\nE1fGxNJ8wxwEawQSFCH0JIfRgYH2rCOleolFKKzxuWfoadPc7HOE47bIPDS6PyaM9XZ8H2dKlwVd\nlle2ckucmqihUp8uUpPMpPqzVz0zNRffKqyuKpnOMnhyXSPJccWKxuPdySxL1srOUdGWkggkAZNe\no7a5TTlGM/YX2200l9Wayu+pXUyreOej8n4MJXlW+lz0sdPVuWFOuck+i4w/oepppq3T1zTzuimZ\nK9Lq6a9lN9a85Sg5N/Hk16auVVEYTcXJZy4rC6+RqpFei47+P3bpfo/1NJl0uFdq0/C3P/ajSmmk\n08p+JmrEkEgKzajVumTiqZT48JJfmyzT29/RGxwcG/st5wJaemdveSri5/ea5LOhesTvUmfWUzuh\nDu3FThJSW7oaCq+6NFe+UZSXlGOWRaz6DT3aOuNDcZ1rlSzzzy/xycQvs9H0s6afVnJb1HpFP+Z3\nHtKmTxsvXxqkXLVU+LcfjForKKLLZanURnjZFx2YXuNBRLV6ePW+tfGRC1ulfTUVP/mRGlk6arJK\nU4KTXRs7SUVhJJe4p9LoxxbGX/DyVPtGpS2qu9v3VMZamyI1OjldqN2Y91Lbvi/HGRCm2EqarUra\not4klysY25/E79Nhj91d/sLqbldDdGMorOMSWGXLPKbL4WHDtrTw7Ip+9nZmloNJOTlLT1uT6txE\nW76aDL2jXZOqEq5VxUJqUt+cPHTp7zV0WARXkrUxqttlfqK4q2OH3Sbafn08jVoIwSc9LJejzedr\nTyn04+hF07o6yHr4g5KMIr7Sw22zcVmIMD099Nk1po6aMJcJScs4N55eqlK/VRSm4uu5QjGPV8Ze\nfdj8hFqzS6W7TTi3HTQhFP2E846s7Uo3doUW1tODok0/c2iqiMNRrdUpym4TSUFJ8NeOPoWd1Vpt\nWlCChWqJYS9zQRui1JZi017gV6Vxlpq3CKhFxTUV4FpGmaOg0sbu9jVied2cvqTqtPZfju9TZRj7\nmOfqaAB5no9qvjStdqZSay/Z9VfQ1afTTplmeott90mjLqNLqZam+Ud7hLa4d3NReUsYbN9KnGqK\nslumlyzVrMnaww6yjRyuVmoocpYw58+qvebjzNcrNRqnp0rFHakscLnq2/h4GYtbUqNPHjZBJFVW\noslrrqJxShGKlCS8clGs0j1N8Zwr9bTYlGT43SXOPh/E3w3OCco7ZPqvIqM2jsVfZVVjTajWm0uv\nCNaeVlHn6e6VWnpr4TdkoZfRJN/wPQTysrlMVYz63TK+l4jF2JqSz448C2iqFNe2uGxdcFgIqDz3\no5Q0ur2xzZbOU17+eD0QXUs1TdG6VtLrmowUs2L7yw+Prgqekl6cr1JbM7mvHOMfQ1gaYqjVJXzs\ndkpKSSUX0iUx091TVdNkY05zhr1o+5GsE0xxOUoyilFyy8NrwOwAoQSAIJAAAAAAABzZONcJTm0o\npZbfgdHFtcbqpVy6SWGB51dlkoauNPe1bcTinHMnw+En54NXZ+96SErZN2y5nnwfkWafTqiL9eU5\nS5lKXVluMFrMlZasR1mqz02xk/o/4GiqUJVQdaSi0mkvId3DdKW1ZksN+a/pnNFFVCaqhtT97ZFW\ngAKEEgAQSABDWepIArlTVL2q4P4xRg1lftUaeilWyw4NJcebfu8D0jB/Z9srZXPV2wsn7WzGPguC\nxK7heoaKN0I7tvt5WGsdfmX6a2V9XeOOIyeYLxx4Mzw7Ngk+9tst3PL3Pr9DTTTGmOyDezwWenuL\ncSaswCQZaAAAIJAGL+za3qJXyssc92Yvdjb7jVXBVwUE28eLeWdgamYg57uG/fsjv+9jk7AVyoxW\nOFx0ObKoWcyXOGvk+pnujXqr4whc42US3NL4fzNUU1FKTy/MqK9NTOmChK3fGKwvVwXAEUAIAkA5\nnNQhKUnhRWWB0QcVWK2qFiWFJJ4ZzqLoV1yUm8uLaUeuEBaSYOybK5abZCTbTy1z6ufDk6u1unal\nXJzW7MU4rq+jwMTWpwhOOHCLWc4aOklFYSwl0SKNDGVekqhNbWljD648M+8vbSTb6LkpEgiMlKKl\nF5TWUySKAENqKy3hICQVUairURbqlnHmsFoAAAAAAAI6ASCqnUU3pumyM1Hh4fQsAkEEgAAAIJK7\nozlTONcts2mk/JgWA82mx6OdlCbsSnWlulnG54/RmrR6pauqU1Fx2zlHD9zLiavBTqLpVWUJLKsn\ntfu4f8C4ipAAAHM5xhFym1GK8WVw1NM7XVC2MprnamBcAAMtmtVffZjxVOMX784/iajyddTO+6/T\nQTzc4y3eSS/il9TTRer9RXmDjJV7pZeMPOMY+RqxmXttBzLLg1F4fg8FVFV9cm7dR3qx02YMtLpP\nbFyfRLJRp9XXqHtjlTUIza8lLp+RfJboteaPO0mku09tdssSlju5KP3Ulh/h+JU1vrsVibimsPHK\nwTuju25W5rODNTO+2FtcozrmpNKbjxjPGPPgmrRKGoWosslZco7HJ8cfADUACKjct23PPXBJ505a\npa6tyjVDcpRi8t+/9DfDdtW5pvxwWzEl1TrLnTp2482Se2C85PocaBP0Z0TnKU6m4OTfL8n9MF2o\nrdteIzVck8qTingyVOGinOzU61TlNJPKUenwEKaGcpau+p2OUKPVjnq8+L+GMG8wVx097h6HqNso\nLrHnKb8TeKkeNqLJO2U6YNJbuI8uXrRy/wAGbuz53ONkdRnvFLPTonzj5dCqFepqvsdNFCcnmTdj\nz+Rofex2ztuhCKfMUuvuyavhJ5RrPb03+svyZqME9RHUS0+2MouOo2tP3RZuyZrSQDM6L3ZKXpUl\nFviKiuCK0mDtCyLtqoshOVUvWnti2njon/XgbIRcY4lNyfmyLqo3V7JOSX+WWGWJVOjbjGVa5qj7\nEvHHk/gU6NX+kytlU3C18SlLmEfBY/rqFotNK6Valc5RWX+1lwXU6Gmm1WQdjkvvWNotxMrSeZbL\ndq4bWo7L0owj1fHrN/JnpmPUKUNQnGdMJT4i5Q5fzJCs2k1He62mve5OMbJS9zz0+h6dqzVNecWY\n/Q7Yamq+NsMRk90FBRTyuXkvs1MI6qGmae+yDkuOMIENC86Kl/5F+ReZ+zv7jUvJY+nBfKO6LWWs\n+RKscV295bbDH7tpZ8+MlfaP9xtXmsfiV6XSKGpusk7M7/Vbk8NYR32lGUuz71B4lsbTL7PSNP62\ntvmvZSjWvist/mazyey/2dzeZKu6vvIqT9/8MGnV6lei1zjNwrsfNn3V5ls7SXptMdmphDtCuErN\ni2yWJcJvjB3XqWtLTO2L3zxHGOr/AK5OdZXHvqLbIKcItqWVnGV1JitZBlWtTujCNU9kuFLGM/D3\neZqIqSGsrDJIA+fzPTUXqvDVzlGCT5hFSwn9T1uz4qMLYw/dRsar5zxhfrkmOgojCcNram8vL565\nNEIRrgoQioxisJI1brMmOalCO6MHn1m38XyWEKKWcLGeWSZaAAAObG1XLb7WODo5sUnBqElGXg2s\n4A8vRaZaimcluhlQW7HO5ct/VmnR1LS32ULO14nH6JP8fzI9Erstk3q73Ne1FWYx8jTXCNFbW+Ti\nuW5yya1nFes9rT/6y/Jmky632tP/AK0f1NRlQABXFkd8JRzjKxnyK6NO6puTscsrGNqSReAmAACq\nbL6654kpZ90Wyt3USkrFXOUo9GoPKNQCMz1tcIuU42RS84M4j2ppJdLX84tGwp1N3cVOSWZPiMfN\n+CApfaeki8d42/dFsS1tmfU0tkodd7aiiIqzTXVynZKyM1tn7peePw+hm1Si+01vnP0dxStWfVz9\nn+fyNSJb03Kd8llQrWfORNXpPeN2913eONuc5LFOG7ZuW7GdvuOzKwAAVGE8ZXQkADmcI2R2zipL\nyZ5WrjbVbOUZQopgkoRWM2Sb8fd4HrmB9m126yy/UYtT9iLXs8f19TXFnk2RiklhJfA7KqaI052u\nTy/tPOPcWGVjzYWSnrqb4qMYWylHHjJYfL+aIpn33abdk3KnDdKl5p4eP68zbVpaqpucY8vzfReS\nONXGFdUJqEc1yShlezlpGtiSVRZo6NPbp501KL73lr3pmrTOuyLuri13jy8+OOC7wK6qK6ZNwWM+\n/p8CKtAOZSjCLlJ4S6siujmyarrlN9IptnRzOKnCUZLKksNAY9FqKk1VOX/iJ+vJNeL56/DBuM0N\nJGM63Kyc+69hN8L+ZpFSB5naObo6lOTjCmGUl4y6/wAD0zPZpK7L1bJyyusc+q/LKLLhZrFZiesr\n0tc5Kucc2xzlLxxnzf5G+dCndXNv930SIu0ldsUvWg1LcpQeHk7oqVNexTnPnOZvLFqSMdLktFKE\nLFXYrZRi357nwb1056mezSKUW65bZ953ib5xLGC+CkoJTe6Xi0sZFWOjmUVOLjJZTWGiSSK5jCMY\nqKikksJe4KEVFRUVtXCR0AY4nXCyO2cU0nk7AAjC8iQAAAAAAAAAAAAgkz2zktXTBPhqTa8+hflJ\ncsFeZFekdpV37motyjGKeMqPi/PkntGedRCKbdUMO/y2+H9eR3UtNo7P30rJTzGC67VnOFgsjfZK\n2UFpXHP2pNcl0nC11qpKfo7i013y5XwZ1fqFROpSwozbTk/Dgq2a3KS9GjBeGGdOGqfWVD+MWRO2\npAyuOt+zOj5xZG3Xff0/0YVrBnhHVfbnV8os6cNQ1xbFf8oF4Myr1S66iP8AsOtt/wDix/2hV4KG\ntR4WV/OJxJav7M6fnFhGoz6rSU6tRVqb2PKxJrDKsdoZ4npsfCROO0Pvab/uG4Yso0lVHsJt+cm2\n/wASx1QaknGLU/ayupmce0Pvab6SLIR1ePXnVn/LFl1JFldFdTzCCT8/EtM7pulHnUST/wAsUVeg\n3Z/v1/4fwI3OM/WxlGhk5aflttSkv+5lK0N2f77d+H8DTpqVp6tik5cttvxbeQtkkyVaV3XV0xUr\nJKK95YcyipLlZDMz2zz1kUswqtsXnCIr1cpv+63x/wCJJGlRSWEMBLnpS7ppc0T+WDJDWaivWQrs\nqnKu2WFJrDj7veeic93DvFZtW9LGfcF42TzHRj7VshXobN0sNr1V4trk2lV9MbqpQaTzFrkLxze2\nZdoJ12SVNiVXt7uMcZOpdo0Rrc+XibhhLlteBdHTwdbjOKbklufng5Wi06xitcNNfFdCdk+iXq6l\ns9biUXPPgo+ZlhrK+09JfGnKzFqO7x9/1NHoVDjBSrUlD2M+B1PSUTxurTwmljjhjtf5irT6uF6d\nkLU4qCbj5P8Ar8jmPaCVO6yuW/LWyPPTk0LS0qDgq0ouO1r3eRMKa64KEIJRXgC3gsjJTipLlNZR\nJzCEa4KEFiK4SOisAAAAAAAAIJAAEEgCCQAAAAAAAAVVX13KbrkmoScX7mgLQV2XQrjCTfE5KKa9\n/QsAAACm/Txv2tynGUc4cHhlcNFXHO+dtifhOeUaiAOKaa6Y7a4qMfJFgAW3QABAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAFd1ippnZLpCLkzytHbXoKLoam\nzlqM2vFuS5/FHoa+meooVUV6s5JT/wCHxKNVoG7Fdp0lPKck/tcmpjNYIapQ7F0kZqW9WRSXilGa\nXJ7x5a7GipY7x7JRxOPnLz+ryenWnGEVJ5aXL8xyz0cddAHM03CSTw2uH5GWnQPM7E1V19FlWphK\nNlM3HdJ+2s8P6HpgAAAAAAAACCQAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAEEgAAAAIJAHKhFNtJJvqdAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAgkAQS\nAAAAAAAAAAAAAAAAAAAAAAAAAAIJAAAAAAAAAAAACCQAAAAAAAAAAAAgABlZxnnqQ5RXVpHFtPeN\nvftzBx4OI6RLPrZWU+VnHOQLyHOMcZfV4RR6JnObZPjCLO5woJP1YvOMdQLCuV9cb40ttTksrjh/\nMsXBRbpXZq6b3ZhVZajjrlNfqBeSAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAQCQAIJAEEgAAAB//Z\n", "text/html": [ "\n", " \n", " " ], "text/plain": [ "" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "YouTubeVideo('xq9YTcv-fQg')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "For a careful walk-through of the discretization of the linear convection equation with finite differences (and also the following steps, up to Step 4), watch **Video Lesson 4** by Prof. Barba on YouTube." ] }, { "cell_type": "code", "execution_count": 10, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQAAAQABAAD/2wCEABALDA4MChAODQ4SERATGCgaGBYWGDEjJR0oOjM9PDkz\nODdASFxOQERXRTc4UG1RV19iZ2hnPk1xeXBkeFxlZ2MBERISGBUYLxoaL2NCOEJjY2NjY2NjY2Nj\nY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY2NjY//AABEIAWgB4AMBIgACEQED\nEQH/xAAaAAEAAwEBAQAAAAAAAAAAAAAAAQIDBAUG/8QAOBAAAgIBAwEFBgUDBAIDAAAAAAECAxEE\nEiExBRNBUXEUIjJhcoEjNJGhsRXB0UJS8PEz4SRik//EABYBAQEBAAAAAAAAAAAAAAAAAAABAv/E\nABcRAQEBAQAAAAAAAAAAAAAAAAABESH/2gAMAwEAAhEDEQA/APvwAAAK7knhtfqBYEEgAQM84AkE\nbl5r9SQAIbS8SHKKaTay+i8wLArGSlna08cPBIEggkAAAAAAAAAAABBJD4WQBJww16d1UZ7UrFJ8\nPmLXOP5OyM4SxtknlZWH1QFgAAIJMtRb3NTnjPK8fmBoDzb+0YR77daod20k0/jyuPTxOzT25pqV\nlkZWSgm2vHjloDcERalFSTynymiQAAAAEAAc1/4mrpqeduHN4fljH8l9PqIajvNieITcG34tDBuQ\nMpdWY13Z1Ftcmko7WufB/wDQG4II3LLWVleAFgZV31WqTrnGSi8PD6FXqqtrkpZSSlleOegGwOO2\n5WQpvjujstSlGXDWeOf1NoylTXZPUTWNzaflHwCNpPCb6/IrVZG2uM4PMZLKEZwnFSjJNNZTMKPw\ntVbSvhaVkflnOV+q/cDqAAUAAAAAAAAAAFLJONcpJZaTaR4uvTrhTOF7bltdiTWGm+rzzjw4PdOO\nfZtE5qTgnj4cpNx9PIC+leO9rzmNcsR+Sx0KPtLTJ4bn/wDnI6Kq41QUI9F4vq/my4GNGqq1Emq3\nLK65i0YXVQXa2nsS9+UJpvPodpSVMZXQted0E0vv/wBAcOl0dF2ihOccTWWpp4cXl8m1Wra7Khqr\nYvd3ak0vFk+wR7pVO23u087U0s/J/I3sqhbTKqS9ySw0vIDzpvVe1J6rupR7ixxjDP8A9eorlfLt\nKKr7pQWmTimnlHWtFF2b52WTfdutZfRP+/BpDTwharFncoKH2QHB2RddHTaSq7Y3Otzbj8sf5Lan\nX2wpUoKEPxZV7pptLDwunmdK0VcI1KuUouqLjF/Jk1aRU1qFdtiSbbbeW8vLyBjPUqqbtnGM5qjc\n5QfXnoiPadTp961Srk9kpxdecceDybLQ0qvu+XHu+75fgS9HCWe8nOeYOvLfRf5Ax0+p1Er61cq+\n7ui5Q2ZyunX9TuMvZ4KVUlnNScY+nH+DUCQAAAAAAADn19Vl+itqqltlOOM/ydBncrHBqqSjLwbW\nQOGWjlfb3+xVOupwqi/Bvxf6L9y2h01lN6bjtrjTGuKznGM5Jzq9232urOcf+Px8uprTDVKxOy+u\ncF1ShhgdQOTXQumodzGcms522bTPRV6iN2ba7Ixx1lapfsB3HJrKrJ3aeyKc41SbcF4vGE/sU1E7\nv6io1Qm9tPH+1tvx9MfuXlXfVp6oQ1CTisSlOO5yYE6XSKt3W2xi7Lp75LHC4wl+iMtbprZ3KdCW\ne7lFfKT6MtBayeduqqeOv4fT9yl8tXRFOzVV89Eqst/uB16WvutNXXjGyKil5JGxjpnKVEZTsjY5\nLKlFYTRsAAAAAAcl25a6CjjdKqai30zlHn6eWo0movhulYlJxjHHDk1uT/do9W+l2OEoS2zg8p4z\nx4omU4wU5SWFFZb8zUrNnXkNXOmNSslZbNb08NKU/wCyWMnaoR9srusqw7IKPK+GS/7/AGN9LK+U\nXK9RW73o48M+BuS0kDzY6WcoapZxfZPLcs4254Xpj+TuV9Tlt3JNycUn4tdS1ltdbSnOMc9MvqRb\n15Sq1N2pdNk3VCcGpKKwmuOhvDSXSjzsi47IpPo9rOydW7UV25+GLXrnH+DGWpjVqLY2ywlt2rxe\nfI1us5jP2Kz2a2Dscp2Pf5JS+XyLSoulVbUpf6ouDlzxw8G89TTW8Tsinz1flyytmsorvjTKWJSW\nenCXzZNq8c1fZrV0LJ2fBBJJeecv/B0xjJ66c8Yiq1HPm8tmEe0q+6tusi41QntUlzuecfyWhrbF\nqY031Rqc1mOJ5F0464yUs4aeOCxzaL4LF5Wz/lnSRoAAAAAAAAAAAAAAAAAAAAACCQBBIAAAAAAA\nAAAAAAABBE5bYSljOFnBYgDx27tRXoVFpTvsVs4xXwR6v98L7nR2NGfcWW2XOx3Tc+VjC6L9kjuh\nVCttwilnyJjGMViKSXyAsAAIOTW2S73T0QS/Fnlya+FR5/XodhSdcbFicU18wPHeonT7fqXaoxbf\ndvHM9sef34+x1XQ1k9Nit1W2t87vd2rHKTR3OuDSThFpdOC0YqKwlhAUoh3dFcGoxcYpYj0XHgaE\nEgAAAAAEHna+UsaxR5xSuG/U9EpKmubm5RT3rbL5osSs9NLUyb7+uuEfDbLJuEklhEkpJjyp0XPU\nNquWKZuyL/3ZeePtk7dRRK9QcbO7a5+FP+TcDTBJqKTeWvHzOW3Qwnf38X+LuTy+emVj9zqAXHBf\n2XDUWynZY8N52pefDX3KWaTOqVSllSpcZSfXblYR6RjKlvV13KWFGMoteecY/gu1nIotDQpZUOM5\n2593PoaWaauxPMcSbT3Lrx0NQNXI8zR6mdfa2p0VsJNNKyFmPdfmvXoemVdcHNTcVuXRliKkAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAgkgAASBAAAkAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAgkAQCQBBIAAAAAAAAAAFLZbK5S8lktF5i\nn5gSAQBIIJAAgkACABIIAEggkAAAAAAAAAAAAAAAAAAAAAAEEkACSCQAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAy1H5ez6WXh8EfQpqPy9n0svD4I+gFjl1NlkbIxqznDfCydRXbHduxzjGQMPbIJ\nRzF5fh4on2qDxhPnD/5+hdaepPOznzyV9kow13aw/m/LAGdesXdxdkXuf7v/AIy8tXCGdyaxw388\nZLvTUv8A0IPT1POYJ54Ayjq13soSi854X6f5Fts67nl4io5isfEzWOnqg04wSa+ZMqa5yzKOX6gY\nvWRTw4yclhPHnx/kvK1y07nDcnz0WXwyZ6aubzjHKbw+uC3dQ2KO33V0WQJqblXGUsZay8FyEklh\ncJACQAAAAAAAAAAAAAAAAAAAIAkAAAAAAAAAAQSAAAAAAAAAAAAAAAAAAAAAAAAAAAAGWo/L2fSy\n8Pgj6FNR+Xs+ll4fBH0AsQSeZ2jotRqNVXdTLb3cVj32udyfT0yB6QPL7M0Wr0+tvt1E1Kuz4I72\n9nPTn9TW/sx3doQ1S1VkFHH4a6PH3A7bLYVR3WTjCPm3gd7W7O73x34ztzzg4+1qNRqdL3WmhVJu\nS3d55eOPmZ1aG2Gvjc1BQS3ZT97O3G308QPSBWqU5VxdkVGb6pPoeQuy9Q/aYyn7t8Zp/iSfLeY4\n8gPZJOTQ02U6SNdsVGS44m5fuzPs/s16K2yb1Nl2/wAJeAHY7a1aq3OO9rKjnkV2QsjurkpLplPJ\n52p7Pst7VhqYQhHENrs3Pc1h8Y+5t2XpbdLXNWqCcmsKHThYz9wO4AgCQQSAAAAAAAAAAAAAAQCQ\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAABlqPy9n0svD4F6FNR+Xs+ll4fAvQCx5fatus\njOMNPuUXteYxy5PdyvlweoeZ2r2nLRThCuEW3tcpSkkknLH3A20k7vaL4XSlLE245jhKPhyTbZq1\nrIxrrzTxmXBloO1FrdTfQq1F09XuTzz4Fr9fbVro6eOllOLx+InwsgT2m7lVB0uaW579nXo8fvg0\nnK72FxWfae5zlf7sf5Kdp6v2PTKashCTkox3LO5+ReOonPU11RSwob7X5Z6Jfv8AoBxQs1EaapOV\n+O/4TXLhjnP3PUjYpTnFZzDGTy7O151UwtdKkpuxvEktsYvGfma6fteF+ruo7vmnOWpJ558F4gaa\nyesjqao0JOqfEnj4MPl/oYy1GpjptQ47nNW+63D/AEm9utUL9LHcoRuk47ZrEuhrpLnfCzckpQsl\nB48cMBTO6eihOcNtrhlxfmU0M9VOM/aq9jT93pydYA86p6r+qT3Ofdc9fhxhYx885L6qyz2jRyql\nZsc3vSXGMPr98HcAIJAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\ny1H5ez6WXh8EfQpqPy9n0svD4I+gFjOyiq7He1Qnjpuing0Mrb6qMd7ZGG54WX1YEwoqrea6oQfT\nMYpFyN8efeXHXklNNZTygKWVV24VtcJpcrdFPBMa4QbcYpN9cEwnGcd0GmvNDdHONyz5ZAzWloim\nu6g023ys9epaNFUJboVQjLnlRSfJM7IQcVKSTk8L5ss5JdWlgCk6a7JRlOuEpR+FuOWvQmqqFUds\nFhZb+7J3x/3L9SQJAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAABlqPy9n0svD4F6FNT+Xs+ll4fCvQCxxazQvU6rTXfh4pbeJxy+fI7SAPIt7HsnXdHvofieO18\n+9nL55fgenTUqaI1xUYqKxiKwjyr5ayen7QqjXqlKUvw5JLLXC4/c9HR1uvQ11p2ZUetnxfcCvZ2\nms0lDqsnCS3OS2Rwll5HsFXtntXPeHP2VC7T6ex2wt+JJRk8tvCTfpkjubv6730Xd3bhiSl8C46r\n5gdWsotvVfdWRg4TUnujnOCut0b1VFtWYLe4vLXk0+Rr46lxg9O8xUvxILhyXyZlo4dox1c5amyD\n07ztiuq8gOaXYlruomr4KNM21HZw05Zw/kl0PZPM7Rlr1rtOtKp9xx3uIp8Z8CmuWq7nX93G12Sc\nVUoeWF0++QPWByVSsnr9+2cYOlNxl4PL/ctro6l1RelksqXvR8ZR8Un4MDpB52nr7R9t32Tj7K+k\nH8SLa9WrWaKddds1Gx73Dok01z98Ad4AAkAAAAAAAEAkAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAGWp/LWfSy8fgXoU1P5az6WXj8C9ALEEkAc+tvemqjYkmt8Yy+Sbxkwp107oaOW\nFFX7nL7J8FF2zprYan2Zu2enTc49MpdcP7Guo1+kr1NFF2e8kt8Pdylw+c+mQModoWPRK3MJTV6r\nkl0acscfbk0v1llV9ilFRqhFvfhvHBbQ67R9qVuene9Vy8Y4w/M7AOXs6+eo0yss8W8cc48MnUCQ\nIBIAgEgCASAIJAAAAAAAAAAAAAAAAAAAAAAQBIAAAAAAAAAAAAAAAAAAAAAAAAAAAEASCCQAAAy1\nP5az6WXh8C9Cmp/LWfSy8fhXoBYgkAebPsfs+KtTr2K94niTWec4J9g0He6ex+/JLZXum3nr/wCz\nq1mneopUYy2yjJSi8Z5Tyc9eidMdIlYm6MpuSxuTQDQ6fQ6GrOl2whZLbnPDfTB1u6uMnFzjuSy1\nnk4v6fZ7N3Hexa75WN7cY97c8fc2loYvVPUKct7XR9AOiqyF1cbK5KUJcprxLnPoqHpdNGlz37ej\nxg6AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAgASAAIBJAEgAAAAAAAAAAAAAAAAAAAAAAAy1\nP5az6WXj8K9Cmp/LW/Sy8PhXoBYgk83tDX36XWU1V0boTx776Zz0/TkDzoaHtl23q6+cq7LE04WJ\nYju5x5cGktH2nO/SSsjKaqlFv8ZJLGctrxeDb+papblKuCTk8Sw8RSk1l/bk7ey7bbuz6rL3mx53\nPGPFgdgAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAIBIAAAAAAAAAAAAAAAAAAEASCCQAAA\nAAAAAMtV+Wt+ll4fAvQz1X5W36WaQ+BegFiGk+qJObVa2jSSqjdJqVstsElltgdGF5DhcIxes0yU\n831rZ8XvfD6msJwsgp1yUovo0+GBYHM9dplSrpWxVbnsUn55wbucIw3uSUcZz4AWBSu2u2O6ucZr\nzTyZ1amu2y2Ecp1PEtyx4AbgxlqaIzjCV0FKfMU319BRqadQpdzZGe1tPHg08AbAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAQSAIJAAAAAAAAAAy1X5a36WXh8C9Cmq/LW/Sy8\nPgXoBY5tRpXdfTarXDum3jannJ0kAeZZ2RXKFkFc0pLC91Pas5a+fPmehVDu6owznC64weVLsnUS\nm4ysh3Tkv9csuO/OP04PR01M69HCm2SlJRw2mBzx7NfskqZaiUn3neRm4r3XnPQ6p0Rt0zotbnGU\ndsn0yZaHR+yKa3uW71MnpLv6grlOPdqbk1l55jjGPVIDfR6KjRVuvTxcYt55eSK9LKGovtd8pK3H\nuuKxE6SQOG/s6N9kLJWPdGMY/CvBp/2NdLpPZp2ONjcJyclDC91t5fJ0gAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAADHVflbfoZpD4F6Geq/K2/SzSHwr0AsQ\nSQB4ne6yOo1tsbLntmoxi6vdjHzS8TSnVa+zWaaueYRlDNn4XV5458OD0dVqFpoRk1lSmo9cYy8G\nNeu732Vwikr8t5fRJAR2tdbTpVKne5b48QjubWef2N7pWy0kp6ZJ2uOYKXHJzvXyWlVzhHKuVcop\n543bTrus7qmdmM7VnAHP2c9a6Ze3xhGzdxt6YOfS2WS7Q1b/APkRiliKsWYt+a/wdeh1XtdLs2bM\nPGMnNZ2tGG78LMo543Lwko/3Apq9XrK76e5TlXsjKa7ptyy+fTgdn6ntC93q6Ci1/wCLdDapLL5f\nz+RvDXuyFDUMSstdUk30azn+Cs9fKGl1FzjFuiza4p5yuP7MB3s12rVCSt96p7sJ7E/Dn9TvIk9s\nG8N4XReJy9n6162qc3p7ads3HFixnDwB2EHEte5a2zSqn360225cYwsP7/2M32pinSzdWe+Scmnx\nHov5YHog856/UV2ahX6XZXCcYVT3fHl4+xXR9p26iVKlTCKsU23v6bZ7ePMD0wedd2jKtTkoKSc3\nCvLwntWW2/1X2KaTtmN99FEq9s7ao2ZzwsrOAPUJPJ1Has6dPK2NaluUp1pvHux8fub6ftOF+us0\nqhiUFlyzw+nT9QO8EZXHK5AEgAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAQBIIJAAAAAAAAAx1\nX5a36WaQ+Behnq/ytv0M0h8C9ALAEAROEZxcZxUk/BlXRU9i2L8N5jjjD/4y+V5lZThDG6SW54WX\n1YHDXqdDqNJPUVxTqpsefdx7yOy26FVErbXthFZl8kVv00LtNOjG2M14eBedSsolVZzGcdssfNAc\nk+1dDVjfbtzJxXuvqi992mqndvqTdcFOb259P4MV2Lp1CUXbqJblJNynzykn/CN7dFG2Vv4k4q2C\nhLa8PjxyBErdI46dT2R7xqdafHPXJSVum7i2dlO2FdiU014prkvV2dRXXRFqVncfBKby/wBSn9Lr\nVV9attavlunvln9AOx2RUtu5bv8Abnk4v6pRXRO26M6lC3u2ms8v09TssphZltYk1jcuqMNFoKtF\np+5hKycM5/EluArZ2hoouzfNJxT3Zj4Lr/JEdRo7NPVbXGM4Oe2GI45yVu7I0999lsrLouxNOMZ4\nXKSeF9kavQQUMQnPPeq1OTzyAlq6LKdS3FyhQ3GxNeSyyZWaOudVb7uMmnKEcGf9MgnqnG65vU8T\nUpZS9F4cG70tMra7JQTnWsRfkBxT7Q7Lr08YylHuuMLY315X9zoc9LvrrVcX3kHJYj/px1/cwh2J\npoWb1bflPKTnwuGsenLN1oIRdGycl3VbqTzy48ePnwApt0mp0lVm2CqfEFNeXH9g7qYWXruMOmKk\n2kuU/L9CdL2fTpaO6TnbHc5LvXuabKz7Pi777o22qd0NjTlmKXyQHVBwnCEo4ccZiUpthbKzasSh\nLbJfP/ovXBV1xhH4YpJFNPQqXY85lZNyb/58gNgAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAY6v8rb9DNIfAvQz1f5S36GaQ+BegFiCTye1P6p7dp3osdwubPnz0f2ApDsi6Dkpahz\ni5qWXnlKW7n+C39Ks30ylbGXdyjLLT4xnhfqYY7WjGxN3zVjfPu5gtzxt+2Op6PZMbodm0x1MZRt\nSe5TeX1fUDtIJAEAkAQCQBAJAEAkAQCQBAJAEAkAAAAAAAAAAAAAAAgkgASQSABAAkEEgAQAJAAA\ngkgASAAAAAAAACAJAAAAAAABjqvytv0s0j8C9CmpTemtSWW4svH4V6AWAAAAAAAAAAEAACQAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAQCQBBIAEAkAAQAJAAAAAACAJBBIEAkAAAAAAAAAQCQAAAAAAAA\nAAAEAkAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAQCQAAAAAACCSABJBIAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAACCSABJBIAAAAAAAAAAAAAAAAAAAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAIJIAZWcZ56kOSXVpFLKd7b3bcxcSi0qWfeys55XTnIG+SHOMcZf\nUw9lb62yfGEadzxFbvdTzjzA0M/aK+/7nd7/AKft6mmDlnoYS7QhrN0lKCxtT4b6Zf2YHUAAJAAA\nAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAQAAJIJAAgk\nAAAB/9k=\n", "text/html": [ "\n", " \n", " " ], "text/plain": [ "" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "YouTubeVideo('y2WaK7_iMRI')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Last but not least" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Remember** to rewrite Step 1 as a fresh Python script or in *your own* Jupyter notebook and then experiment by changing the discretization parameters. Once you have done this, you will be ready for [Step 2](./02_Step_2.ipynb).\n", "\n", "\n", "***" ] }, { "cell_type": "code", "execution_count": 11, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "\n", "\n", "\n", "\n", "\n" ], "text/plain": [ "" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from IPython.core.display import HTML\n", "def css_styling():\n", " styles = open(\"../styles/custom.css\", \"r\").read()\n", " return HTML(styles)\n", "css_styling()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "> (The cell above executes the style for this notebook. We modified a style we found on the GitHub of [CamDavidsonPilon](https://github.com/CamDavidsonPilon), [@Cmrn_DP](https://twitter.com/cmrn_dp).)" ] } ], "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.5.2" } }, "nbformat": 4, "nbformat_minor": 0 }