{ "cells": [ { "cell_type": "code", "execution_count": 1, "source": [ "%pylab inline" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "Populating the interactive namespace from numpy and matplotlib\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2021-01-13T06:16:40.844266Z", "start_time": "2021-01-13T06:16:39.992092Z" } } }, { "cell_type": "markdown", "source": [ "## Notebook magic" ], "metadata": {} }, { "cell_type": "code", "execution_count": 1, "source": [ "from IPython.core.magic import Magics, magics_class, line_cell_magic\n", "from IPython.core.magic import cell_magic, register_cell_magic, register_line_magic\n", "from IPython.core.magic_arguments import argument, magic_arguments, parse_argstring\n", "import subprocess\n", "import os" ], "outputs": [], "metadata": { "ExecuteTime": { "end_time": "2021-01-13T06:16:40.857076Z", "start_time": "2021-01-13T06:16:40.852721Z" } } }, { "cell_type": "code", "execution_count": 2, "source": [ "@magics_class\n", "class PyboardMagic(Magics):\n", " @cell_magic\n", " @magic_arguments()\n", " @argument('-skip')\n", " @argument('-unix')\n", " @argument('-pyboard')\n", " @argument('-file')\n", " @argument('-data')\n", " @argument('-time')\n", " @argument('-memory')\n", " def micropython(self, line='', cell=None):\n", " args = parse_argstring(self.micropython, line)\n", " if args.skip: # doesn't care about the cell's content\n", " print('skipped execution')\n", " return None # do not parse the rest\n", " if args.unix: # tests the code on the unix port. Note that this works on unix only\n", " with open('/dev/shm/micropython.py', 'w') as fout:\n", " fout.write(cell)\n", " proc = subprocess.Popen([\"../../micropython/ports/unix/micropython\", \"/dev/shm/micropython.py\"], \n", " stdout=subprocess.PIPE, stderr=subprocess.PIPE)\n", " print(proc.stdout.read().decode(\"utf-8\"))\n", " print(proc.stderr.read().decode(\"utf-8\"))\n", " return None\n", " if args.file: # can be used to copy the cell content onto the pyboard's flash\n", " spaces = \" \"\n", " try:\n", " with open(args.file, 'w') as fout:\n", " fout.write(cell.replace('\\t', spaces))\n", " printf('written cell to {}'.format(args.file))\n", " except:\n", " print('Failed to write to disc!')\n", " return None # do not parse the rest\n", " if args.data: # can be used to load data from the pyboard directly into kernel space\n", " message = pyb.exec(cell)\n", " if len(message) == 0:\n", " print('pyboard >>>')\n", " else:\n", " print(message.decode('utf-8'))\n", " # register new variable in user namespace\n", " self.shell.user_ns[args.data] = string_to_matrix(message.decode(\"utf-8\"))\n", " \n", " if args.time: # measures the time of executions\n", " pyb.exec('import utime')\n", " message = pyb.exec('t = utime.ticks_us()\\n' + cell + '\\ndelta = utime.ticks_diff(utime.ticks_us(), t)' + \n", " \"\\nprint('execution time: {:d} us'.format(delta))\")\n", " print(message.decode('utf-8'))\n", " \n", " if args.memory: # prints out memory information \n", " message = pyb.exec('from micropython import mem_info\\nprint(mem_info())\\n')\n", " print(\"memory before execution:\\n========================\\n\", message.decode('utf-8'))\n", " message = pyb.exec(cell)\n", " print(\">>> \", message.decode('utf-8'))\n", " message = pyb.exec('print(mem_info())')\n", " print(\"memory after execution:\\n========================\\n\", message.decode('utf-8'))\n", "\n", " if args.pyboard:\n", " message = pyb.exec(cell)\n", " print(message.decode('utf-8'))\n", "\n", "ip = get_ipython()\n", "ip.register_magics(PyboardMagic)" ], "outputs": [], "metadata": { "ExecuteTime": { "end_time": "2021-01-13T06:16:40.947944Z", "start_time": "2021-01-13T06:16:40.865720Z" } } }, { "cell_type": "markdown", "source": [ "## pyboard" ], "metadata": {} }, { "cell_type": "code", "execution_count": 57, "source": [ "import pyboard\n", "pyb = pyboard.Pyboard('/dev/ttyACM0')\n", "pyb.enter_raw_repl()" ], "outputs": [], "metadata": { "ExecuteTime": { "end_time": "2020-05-07T07:35:35.126401Z", "start_time": "2020-05-07T07:35:35.105824Z" } } }, { "cell_type": "code", "execution_count": 9, "source": [ "pyb.exit_raw_repl()\n", "pyb.close()" ], "outputs": [], "metadata": { "ExecuteTime": { "end_time": "2020-05-19T19:11:18.145548Z", "start_time": "2020-05-19T19:11:18.137468Z" } } }, { "cell_type": "code", "execution_count": 58, "source": [ "%%micropython -pyboard 1\n", "\n", "import utime\n", "import ulab as np\n", "\n", "def timeit(n=1000):\n", " def wrapper(f, *args, **kwargs):\n", " func_name = str(f).split(' ')[1]\n", " def new_func(*args, **kwargs):\n", " run_times = np.zeros(n, dtype=np.uint16)\n", " for i in range(n):\n", " t = utime.ticks_us()\n", " result = f(*args, **kwargs)\n", " run_times[i] = utime.ticks_diff(utime.ticks_us(), t)\n", " print('{}() execution times based on {} cycles'.format(func_name, n, (delta2-delta1)/n))\n", " print('\\tbest: %d us'%np.min(run_times))\n", " print('\\tworst: %d us'%np.max(run_times))\n", " print('\\taverage: %d us'%np.mean(run_times))\n", " print('\\tdeviation: +/-%.3f us'%np.std(run_times)) \n", " return result\n", " return new_func\n", " return wrapper\n", "\n", "def timeit(f, *args, **kwargs):\n", " func_name = str(f).split(' ')[1]\n", " def new_func(*args, **kwargs):\n", " t = utime.ticks_us()\n", " result = f(*args, **kwargs)\n", " print('execution time: ', utime.ticks_diff(utime.ticks_us(), t), ' us')\n", " return result\n", " return new_func" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2020-05-07T07:35:38.725924Z", "start_time": "2020-05-07T07:35:38.645488Z" } } }, { "cell_type": "markdown", "source": [ "__END_OF_DEFS__" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "# numpy.linalg\n", "\n", "Functions in the `linalg` module can be called by prepending them by `numpy.linalg.`. The module defines the following seven functions:\n", "\n", "1. [numpy.linalg.cholesky](#cholesky)\n", "1. [numpy.linalg.det](#det)\n", "1. [numpy.linalg.eig](#eig)\n", "1. [numpy.linalg.inv](#inv)\n", "1. [numpy.linalg.norm](#norm)\n", "1. [numpy.linalg.qr](#qr)" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## cholesky\n", "\n", "`numpy`: https://docs.scipy.org/doc/numpy-1.17.0/reference/generated/numpy.linalg.cholesky.html\n", "\n", "The function of the Cholesky decomposition takes a positive definite, symmetric square matrix as its single argument, and returns the *square root matrix* in the lower triangular form. If the input argument does not fulfill the positivity or symmetry condition, a `ValueError` is raised." ], "metadata": {} }, { "cell_type": "code", "execution_count": 18, "source": [ "%%micropython -unix 1\n", "\n", "from ulab import numpy as np\n", "\n", "a = np.array([[25, 15, -5], [15, 18, 0], [-5, 0, 11]])\n", "print('a: ', a)\n", "print('\\n' + '='*20 + '\\nCholesky decomposition\\n', np.linalg.cholesky(a))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "a: array([[25.0, 15.0, -5.0],\n", "\t [15.0, 18.0, 0.0],\n", "\t [-5.0, 0.0, 11.0]], dtype=float)\n", "\n", "====================\n", "Cholesky decomposition\n", " array([[5.0, 0.0, 0.0],\n", "\t [3.0, 3.0, 0.0],\n", "\t [-1.0, 1.0, 3.0]], dtype=float)\n", "\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2020-03-10T19:25:21.754166Z", "start_time": "2020-03-10T19:25:21.740726Z" }, "scrolled": true } }, { "cell_type": "markdown", "source": [ "## det\n", "\n", "`numpy`: https://docs.scipy.org/doc/numpy/reference/generated/numpy.linalg.det.html\n", "\n", "The `det` function takes a square matrix as its single argument, and calculates the determinant. The calculation is based on successive elimination of the matrix elements, and the return value is a float, even if the input array was of integer type." ], "metadata": {} }, { "cell_type": "code", "execution_count": 495, "source": [ "%%micropython -unix 1\n", "\n", "from ulab import numpy as np\n", "\n", "a = np.array([[1, 2], [3, 4]], dtype=np.uint8)\n", "print(np.linalg.det(a))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "-2.0\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2019-10-19T13:27:24.246995Z", "start_time": "2019-10-19T13:27:24.228698Z" }, "scrolled": true } }, { "cell_type": "markdown", "source": [ "### Benchmark\n", "\n", "Since the routine for calculating the determinant is pretty much the same as for finding the [inverse of a matrix](#inv), the execution times are similar:" ], "metadata": {} }, { "cell_type": "code", "execution_count": 557, "source": [ "%%micropython -pyboard 1\n", "\n", "from ulab import numpy as np\n", "\n", "@timeit\n", "def matrix_det(m):\n", " return np.linalg.inv(m)\n", "\n", "m = np.array([[1, 2, 3, 4, 5, 6, 7, 8], [0, 5, 6, 4, 5, 6, 4, 5], \n", " [0, 0, 9, 7, 8, 9, 7, 8], [0, 0, 0, 10, 11, 12, 11, 12], \n", " [0, 0, 0, 0, 4, 6, 7, 8], [0, 0, 0, 0, 0, 5, 6, 7], \n", " [0, 0, 0, 0, 0, 0, 7, 6], [0, 0, 0, 0, 0, 0, 0, 2]])\n", "\n", "matrix_det(m)" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "execution time: 294 us\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2019-10-20T07:14:59.778987Z", "start_time": "2019-10-20T07:14:59.740021Z" } } }, { "cell_type": "markdown", "source": [ "## eig\n", "\n", "`numpy`: https://docs.scipy.org/doc/numpy/reference/generated/numpy.linalg.eig.html\n", "\n", "The `eig` function calculates the eigenvalues and the eigenvectors of a real, symmetric square matrix. If the matrix is not symmetric, a `ValueError` will be raised. The function takes a single argument, and returns a tuple with the eigenvalues, and eigenvectors. With the help of the eigenvectors, amongst other things, you can implement sophisticated stabilisation routines for robots." ], "metadata": {} }, { "cell_type": "code", "execution_count": 18, "source": [ "%%micropython -unix 1\n", "\n", "from ulab import numpy as np\n", "\n", "a = np.array([[1, 2, 1, 4], [2, 5, 3, 5], [1, 3, 6, 1], [4, 5, 1, 7]], dtype=np.uint8)\n", "x, y = np.linalg.eig(a)\n", "print('eigenvectors of a:\\n', y)\n", "print('\\neigenvalues of a:\\n', x)" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "eigenvectors of a:\n", " array([[0.8151560042509081, -0.4499411232970823, -0.1644660242574522, 0.3256141906686505],\n", " [0.2211334179893007, 0.7846992598235538, 0.08372081379922657, 0.5730077734355189],\n", " [-0.1340114162071679, -0.3100776411558949, 0.8742786816656, 0.3486109343758527],\n", " [-0.5183258053659028, -0.292663481927148, -0.4489749870391468, 0.6664142156731531]], dtype=float)\n", "\n", "eigenvalues of a:\n", " array([-1.165288365404889, 0.8029365530314914, 5.585625756072663, 13.77672605630074], dtype=float)\n", "\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2020-11-03T20:25:26.952290Z", "start_time": "2020-11-03T20:25:26.930184Z" } } }, { "cell_type": "markdown", "source": [ "The same matrix diagonalised with `numpy` yields:" ], "metadata": {} }, { "cell_type": "code", "execution_count": 6, "source": [ "a = array([[1, 2, 1, 4], [2, 5, 3, 5], [1, 3, 6, 1], [4, 5, 1, 7]], dtype=np.uint8)\n", "x, y = eig(a)\n", "print('eigenvectors of a:\\n', y)\n", "print('\\neigenvalues of a:\\n', x)" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "eigenvectors of a:\n", " [[ 0.32561419 0.815156 0.44994112 -0.16446602]\n", " [ 0.57300777 0.22113342 -0.78469926 0.08372081]\n", " [ 0.34861093 -0.13401142 0.31007764 0.87427868]\n", " [ 0.66641421 -0.51832581 0.29266348 -0.44897499]]\n", "\n", "eigenvalues of a:\n", " [13.77672606 -1.16528837 0.80293655 5.58562576]\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2020-11-03T20:13:27.236159Z", "start_time": "2020-11-03T20:13:27.069967Z" } } }, { "cell_type": "markdown", "source": [ "When comparing results, we should keep two things in mind: \n", "\n", "1. the eigenvalues and eigenvectors are not necessarily sorted in the same way\n", "2. an eigenvector can be multiplied by an arbitrary non-zero scalar, and it is still an eigenvector with the same eigenvalue. This is why all signs of the eigenvector belonging to 5.58, and 0.80 are flipped in `ulab` with respect to `numpy`. This difference, however, is of absolutely no consequence. " ], "metadata": {} }, { "cell_type": "markdown", "source": [ "### Computation expenses\n", "\n", "Since the function is based on [Givens rotations](https://en.wikipedia.org/wiki/Givens_rotation) and runs till convergence is achieved, or till the maximum number of allowed rotations is exhausted, there is no universal estimate for the time required to find the eigenvalues. However, an order of magnitude can, at least, be guessed based on the measurement below:" ], "metadata": {} }, { "cell_type": "code", "execution_count": 559, "source": [ "%%micropython -pyboard 1\n", "\n", "from ulab import numpy as np\n", "\n", "@timeit\n", "def matrix_eig(a):\n", " return np.linalg.eig(a)\n", "\n", "a = np.array([[1, 2, 1, 4], [2, 5, 3, 5], [1, 3, 6, 1], [4, 5, 1, 7]], dtype=np.uint8)\n", "\n", "matrix_eig(a)" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "execution time: 111 us\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2019-10-20T07:18:52.520515Z", "start_time": "2019-10-20T07:18:52.499653Z" } } }, { "cell_type": "markdown", "source": [ "## inv\n", "\n", "`numpy`: https://docs.scipy.org/doc/numpy-1.17.0/reference/generated/numpy.linalg.inv.html\n", "\n", "A square matrix, provided that it is not singular, can be inverted by calling the `inv` function that takes a single argument. The inversion is based on successive elimination of elements in the lower left triangle, and raises a `ValueError` exception, if the matrix turns out to be singular (i.e., one of the diagonal entries is zero)." ], "metadata": {} }, { "cell_type": "code", "execution_count": 5, "source": [ "%%micropython -unix 1\n", "\n", "from ulab import numpy as np\n", "\n", "m = np.array([[1, 2, 3, 4], [4, 5, 6, 4], [7, 8.6, 9, 4], [3, 4, 5, 6]])\n", "\n", "print(np.linalg.inv(m))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "array([[-2.166666666666667, 1.500000000000001, -0.8333333333333337, 1.0],\n", " [1.666666666666667, -3.333333333333335, 1.666666666666668, -0.0],\n", " [0.1666666666666666, 2.166666666666668, -0.8333333333333337, -1.0],\n", " [-0.1666666666666667, -0.3333333333333333, 0.0, 0.5]], dtype=float64)\n", "\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2021-01-13T06:17:13.053816Z", "start_time": "2021-01-13T06:17:13.038403Z" } } }, { "cell_type": "markdown", "source": [ "### Computation expenses\n", "\n", "Note that the cost of inverting a matrix is approximately twice as many floats (RAM), as the number of entries in the original matrix, and approximately as many operations, as the number of entries. Here are a couple of numbers: " ], "metadata": {} }, { "cell_type": "code", "execution_count": 552, "source": [ "%%micropython -pyboard 1\n", "\n", "from ulab import numpy as np\n", "\n", "@timeit\n", "def invert_matrix(m):\n", " return np.linalg.inv(m)\n", "\n", "m = np.array([[1, 2,], [4, 5]])\n", "print('2 by 2 matrix:')\n", "invert_matrix(m)\n", "\n", "m = np.array([[1, 2, 3, 4], [4, 5, 6, 4], [7, 8.6, 9, 4], [3, 4, 5, 6]])\n", "print('\\n4 by 4 matrix:')\n", "invert_matrix(m)\n", "\n", "m = np.array([[1, 2, 3, 4, 5, 6, 7, 8], [0, 5, 6, 4, 5, 6, 4, 5], \n", " [0, 0, 9, 7, 8, 9, 7, 8], [0, 0, 0, 10, 11, 12, 11, 12], \n", " [0, 0, 0, 0, 4, 6, 7, 8], [0, 0, 0, 0, 0, 5, 6, 7], \n", " [0, 0, 0, 0, 0, 0, 7, 6], [0, 0, 0, 0, 0, 0, 0, 2]])\n", "print('\\n8 by 8 matrix:')\n", "invert_matrix(m)" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "2 by 2 matrix:\n", "execution time: 65 us\n", "\n", "4 by 4 matrix:\n", "execution time: 105 us\n", "\n", "8 by 8 matrix:\n", "execution time: 299 us\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2019-10-20T07:10:39.190734Z", "start_time": "2019-10-20T07:10:39.138872Z" } } }, { "cell_type": "markdown", "source": [ "The above-mentioned scaling is not obeyed strictly. The reason for the discrepancy is that the function call is still the same for all three cases: the input must be inspected, the output array must be created, and so on. " ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## norm\n", "\n", "`numpy`: https://numpy.org/doc/stable/reference/generated/numpy.linalg.norm.html\n", "\n", "The function takes a vector or matrix without options, and returns its 2-norm, i.e., the square root of the sum of the square of the elements." ], "metadata": {} }, { "cell_type": "code", "execution_count": 11, "source": [ "%%micropython -unix 1\n", "\n", "from ulab import numpy as np\n", "\n", "a = np.array([1, 2, 3, 4, 5])\n", "b = np.array([[1, 2, 3], [4, 5, 6], [7, 8, 9]])\n", "\n", "print('norm of a:', np.linalg.norm(a))\n", "print('norm of b:', np.linalg.norm(b))" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "norm of a: 7.416198487095663\n", "norm of b: 16.88194301613414\n", "\n", "\n" ] } ], "metadata": { "ExecuteTime": { "end_time": "2020-07-23T20:41:10.341349Z", "start_time": "2020-07-23T20:41:10.327624Z" } } }, { "cell_type": "markdown", "source": [ "## qr\n", "\n", "`numpy`: https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html\n", "\n", "\n", "The function computes the QR decomposition of a matrix `m` of dimensions `(M, N)`, i.e., it returns two such matrices, `q`', and `r`, that `m = qr`, where `q` is orthonormal, and `r` is upper triangular. In addition to the input matrix, which is the first positional argument, the function accepts the `mode` keyword argument with a default value of `reduced`. If `mode` is `reduced`, `q`, and `r` are returned in the reduced representation. Otherwise, the outputs will have dimensions `(M, M)`, and `(M, N)`, respectively." ], "metadata": {} }, { "cell_type": "code", "execution_count": 3, "source": [ "%%micropython -unix 1\n", "\n", "from ulab import numpy as np\n", "\n", "A = np.arange(6).reshape((3, 2))\n", "print('A: \\n', A)\n", "\n", "print('complete decomposition')\n", "q, r = np.linalg.qr(A, mode='complete')\n", "print('q: \\n', q)\n", "print()\n", "print('r: \\n', r)\n", "\n", "print('\\n\\nreduced decomposition')\n", "q, r = np.linalg.qr(A, mode='reduced')\n", "print('q: \\n', q)\n", "print()\n", "print('r: \\n', r)\n" ], "outputs": [ { "output_type": "stream", "name": "stdout", "text": [ "A: \n", " array([[0, 1],\n", " [2, 3],\n", " [4, 5]], dtype=int16)\n", "complete decomposition\n", "q: \n", " array([[0.0, -0.9128709291752768, 0.408248290463863],\n", " [-0.447213595499958, -0.3651483716701107, -0.8164965809277261],\n", " [-0.8944271909999159, 0.1825741858350553, 0.408248290463863]], dtype=float64)\n", "\n", "r: \n", " array([[-4.47213595499958, -5.813776741499454],\n", " [0.0, -1.095445115010332],\n", " [0.0, 0.0]], dtype=float64)\n", "\n", "\n", "reduced decomposition\n", "q: \n", " array([[0.0, -0.9128709291752768],\n", " [-0.447213595499958, -0.3651483716701107],\n", " [-0.8944271909999159, 0.1825741858350553]], dtype=float64)\n", "\n", "r: \n", " array([[-4.47213595499958, -5.813776741499454],\n", " [0.0, -1.095445115010332]], dtype=float64)\n", "\n", "\n" ] } ], "metadata": {} } ], "metadata": { "kernelspec": { "name": "python3", "display_name": "Python 3.8.5 64-bit ('base': conda)" }, "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.5" }, "toc": { "base_numbering": 1, "nav_menu": {}, "number_sections": true, "sideBar": true, "skip_h1_title": false, "title_cell": "Table of Contents", "title_sidebar": "Contents", "toc_cell": false, "toc_position": { "height": "calc(100% - 180px)", "left": "10px", "top": "150px", "width": "382.797px" }, "toc_section_display": true, "toc_window_display": true }, "varInspector": { "cols": { "lenName": 16, "lenType": 16, "lenVar": 40 }, "kernels_config": { "python": { "delete_cmd_postfix": "", "delete_cmd_prefix": "del ", "library": "var_list.py", "varRefreshCmd": "print(var_dic_list())" }, "r": { "delete_cmd_postfix": ") ", "delete_cmd_prefix": "rm(", "library": "var_list.r", "varRefreshCmd": "cat(var_dic_list()) " } }, "types_to_exclude": [ "module", "function", "builtin_function_or_method", "instance", "_Feature" ], "window_display": false }, "interpreter": { "hash": "ce9a02f9f7db620716422019cafa4bc1786ca85daa298b819f6da075e7993842" } }, "nbformat": 4, "nbformat_minor": 4 }