{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Accuracy of Newton-Cotes\n", "\n", "Copyright (C) 2020 Andreas Kloeckner\n", "\n", "
\n", "MIT License\n", "Permission is hereby granted, free of charge, to any person obtaining a copy\n", "of this software and associated documentation files (the \"Software\"), to deal\n", "in the Software without restriction, including without limitation the rights\n", "to use, copy, modify, merge, publish, distribute, sublicense, and/or sell\n", "copies of the Software, and to permit persons to whom the Software is\n", "furnished to do so, subject to the following conditions:\n", "\n", "The above copyright notice and this permission notice shall be included in\n", "all copies or substantial portions of the Software.\n", "\n", "THE SOFTWARE IS PROVIDED \"AS IS\", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR\n", "IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,\n", "FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE\n", "AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER\n", "LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,\n", "OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN\n", "THE SOFTWARE.\n", "
" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "import numpy as np\n", "import numpy.linalg as la\n", "import matplotlib.pyplot as pt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A function to make Vandermonde matrices:\n", " \n", "(Note that the ordering of this matrix matches the convention in our class but *disagrees* with `np.vander`.)" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "def vander(nodes, ncolumns=None):\n", " if ncolumns is None:\n", " ncolumns = len(nodes)\n", " result = np.empty((len(nodes), ncolumns))\n", " for i in range(ncolumns):\n", " result[:, i] = nodes**i\n", " return result" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Fix a set of nodes:" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "# nodes = [0.5] # Midpoint\n", "# nodes = [0]\n", "#nodes = [0, 1] # Trapezoidal\n", "nodes = [0, 0.5, 1] # Simpson's\n", "#nodes = [0, 1/3, 1]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Find the weights for the Newton-Cotes rule for the given nodes on $[0,1]$:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "(a, b) = (0, 1)\n", "nodes = np.array(nodes)\n", "n = len(nodes)\n", "\n", "degs = np.arange(n)\n", "rhs = 1/(degs+1)*(b**(degs+1 - a**(degs+1)))\n", "weights = la.solve(vander(nodes).T, rhs)\n", "print(weights)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here is a function and its definite integral from $0$ to $x$:\n", "\n", "$$\\text{int_f}(x)=\\int_0^x f(\\xi)d\\xi$$" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "fdeg = 9\n", "def f(x):\n", " return sum((-x)**i for i in range(fdeg + 1))\n", "def int_f(x):\n", " return sum(\n", " (-1)**i*1/(i+1)*(\n", " (x)**(i+1)-0**(i+1)\n", " )\n", " for i in range(fdeg + 1))\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Plotted:" ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "plot_x = np.linspace(0, 1, 200)\n", "\n", "pt.plot(plot_x, f(plot_x), label=\"f\")\n", "pt.fill_between(plot_x, 0*plot_x, f(plot_x),alpha=0.3)\n", "pt.plot(plot_x, int_f(plot_x), label=\"$\\int f$\")\n", "pt.grid()\n", "pt.legend(loc=\"best\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This here plots the function, the interpolant, and the area under the interpolant:" ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "# fix nodes\n", "h = 1\n", "x = nodes * h\n", "\n", "# find interpolant\n", "coeffs = la.solve(vander(x), f(x))\n", "\n", "# evaluate interpolant\n", "plot_x = np.linspace(0, h, 200)\n", "interpolant = vander(plot_x, len(coeffs)) @ coeffs\n", "\n", "# plot\n", "pt.plot(plot_x, f(plot_x), label=\"f\")\n", "pt.plot(plot_x, interpolant, label=\"Interpolant\")\n", "pt.fill_between(plot_x, 0*plot_x, interpolant, alpha=0.3, color=\"green\")\n", "pt.plot(x, f(x), \"og\")\n", "pt.grid()\n", "pt.legend(loc=\"best\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Compute the following:\n", "\n", "* The true integral as `true_val` (from `int_f`)\n", "* The quadrature result as `quad` (using `x` and `weights` and `h`)\n", "* The error as `err` (the difference of the two)\n", "\n", " (Do not be tempted to compute a relative error--that has one order lower.)\n", "\n", "Compare the error for $h=1,0.5,0.25$. What order of accuracy do you observe?" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Estimate the order of accuracy:\n", "\n", "We assume that the error depends on the mesh spacings $h$ as\n", "$E(h)\\approx C h^p$ for some unknown power $p$. Taking the $\\log$\n", "of this approximate equality reveals a linear function in $p$:\n", "$$\n", "E(h) \\approx C h^p \\quad \\iff \\quad \\log E(h) \\approx \\log(C) +\n", "p\\log(h).\n", "$$\n", "You can now either do a least-squares fit for $\\log C$ and $p$ from\n", "a few data points $(h,E(h))$ (more accurate, more robust), or you\n", "can use just two grid sizes $h_1$ and $h_2$, and estimate the slope:\n", "(less accurate, less robust)\n", "$$\n", " p \\approx \\frac{ \\log(\\frac{E(h_2)}{E(h_1)}) } {\\log(\\frac{h_2}{h_1})}.\n", "$$\n", "This is called the *empirical order of convergence* or EOC.\n", "\n" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "for i in range(len(errors)-1):\n", " print(np.log(errors[i+1]/errors[i])/np.log(1/2))" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "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.3" } }, "nbformat": 4, "nbformat_minor": 4 }