{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# p08: Eigenvalues of harmonic oscillator \n", "\n", "$$\n", "-u'' + x^2 u = \\lambda u, \\qquad x \\in \\mathbb{R}\n", "$$" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "%config InlineBackend.figure_format='svg'\n", "from numpy import pi,arange,linspace,sin,zeros,diag,sort\n", "from scipy.linalg import toeplitz\n", "from numpy.linalg import eig" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "N = 6\n", " 4.614729169954720e-01\n", " 7.494134621050521e+00\n", " 7.720916053006558e+00\n", " 2.883248377834012e+01\n", "N = 12\n", " 9.781372812986103e-01\n", " 3.171605320647185e+00\n", " 4.455935291166790e+00\n", " 8.924529058119932e+00\n", "N = 18\n", " 9.999700014993121e-01\n", " 3.000644066795828e+00\n", " 4.992595324407724e+00\n", " 7.039571897981499e+00\n", "N = 24\n", " 9.999999976290406e-01\n", " 3.000000098410861e+00\n", " 4.999997965273280e+00\n", " 7.000024998156542e+00\n", "N = 30\n", " 9.999999999999772e-01\n", " 3.000000000000737e+00\n", " 4.999999999975596e+00\n", " 7.000000000508624e+00\n", "N = 36\n", " 9.999999999999927e-01\n", " 2.999999999999998e+00\n", " 5.000000000000000e+00\n", " 6.999999999999988e+00\n" ] } ], "source": [ "L = 8.0\n", "for N in range(6,37,6):\n", " h = 2.0*pi/N; x = h*linspace(1,N,N); x = L*(x-pi)/pi\n", " col = zeros(N)\n", " col[0] = -pi**2/(3.0*h**2) - 1.0/6.0\n", " col[1:] = -0.5*(-1.0)**arange(1,N)/sin(0.5*h*arange(1,N))**2\n", " D2 = (pi/L)**2 * toeplitz(col)\n", " evals,evecs = eig(-D2 + diag(x**2))\n", " eigenvalues = sort(evals)\n", " print(\"N = %d\" % N)\n", " for e in eigenvalues[0:4]:\n", " print(\"%24.15e\" % e)" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "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.12.8" } }, "nbformat": 4, "nbformat_minor": 4 }