{ "cells": [ { "cell_type": "markdown", "id": "c42971d8", "metadata": {}, "source": [ "# Time Dependent Reaction Diffuison: spiral wave\n", "\n", "This demonstrates the simulation of spiral waves in an excitable media. It consists of system of reaction diffusion equations with two components. Both the model parameters and the approach for discretizing the system are taken from http://www.scholarpedia.org/article/Barkley_model.\n", "\n", "We use the _Barkley model_ in its simplest form:\n", "\\begin{align*}\n", "\\frac{\\partial u}{\\partial_t} &= \\frac{1}{\\varepsilon}f(u,v) + \\Delta u \\\\\n", "\\frac{\\partial v}{\\partial_t} &= h(u,v)\n", "\\end{align*}\n", "where\n", "\\begin{equation*}\n", "f(u,v)=u\\Big(1-u\\Big)\\Big(u-\\frac{v+b}{a}\\Big)\n", "\\end{equation*}\n", "The function $h$ can take different forms, e.g., in its simplest form\n", "\\begin{equation*}\n", "h(u,v) = u - v~.\n", "\\end{equation*}\n", "Finally, $\\varepsilon,a,b$ for more details on how to chose these parameters check the web page provided above.\n", "\n", "We employ a carefully constructed linear time stepping scheme for this model: let $u^n,v^n$ be given functions approximating the solution at a time $t^n$. To compute approximations $u^{m+1},v^{m+1}$ at a later time\n", "$t^{n+1}=t^n+\\tau$ we first split up the non linear function $f$ as follows:\n", "\\begin{align*}\n", "f(u,v) = f_I(u,u,v) + f_E(u,v)\n", "\\end{align*}\n", "where using $u^*(V):=\\frac{V+b}{a}$:\n", "\\begin{align*}\n", "f_I(u,U,V) &= \\begin{cases}\n", "u\\;(1-U)\\;(\\;U-U^*(V)\\;) & U < U^*(V) \\\\\n", "-u\\;U\\;(\\;U-U^*(V)\\;) & U \\geq U^*(V)\n", "\\end{cases} \\\\\n", "\\text{and} \\\\\n", "f_E(U,V) &= \\begin{cases}\n", "0 & U < U^*(V) \\\\\n", "U\\;(\\;U-U^*(V)\\;) & U \\geq U^*(V)\n", "\\end{cases} \\\\\n", "\\end{align*}\n", "Thus $f_I(u,U,V) = -m(U,V)u$ with\n", "\\begin{align*}\n", "m(U,V) &= \\begin{cases}\n", "(U-1)\\;(\\;U-U^*(V)\\;) & U < U^*(V) \\\\\n", "U\\;(\\;U-U^*(V)\\;) & U \\geq U^*(V)\n", "\\end{cases}\n", "\\end{align*}\n", "Note that $u,v$ are assumed to take values only between zero and one so that therefore $m(u^n,v^n) > 0$. Therefore, the following time discrete version of the Barkley model has a linear, positive definite elliptic operator on its left hand side:\n", "\\begin{align*}\n", "-\\tau\\Delta u^{n+1} +\n", "(1+\\frac{\\tau}{\\varepsilon} m(u^n,v^n))\\; u^{n+1}\n", "&= u^n + \\frac{\\tau}{\\varepsilon} f_E(u^n,v^n) \\\\\n", "v^{n+1} &= v^n + \\tau h(u^n,v^n)\n", "\\end{align*}\n", "Which can now be solved using a finite element discretization for $u^n,v^n$.\n", "\n", "Note that by taking the slow reaction $h(u,v)$ explicitly, the equation for $v^{n+1}$ is purely algebraic. We will therefore construct a scalar model for computing $u^{n+1}$ only and compute $v^{{n+1}}$ be using the interpolation method on the space applied to\n", "$v^n + \\tau h(u^n,v^n)$.\n", "\n", "Let's get started by importing some standard python packages, ufl, and some part of the dune-fempy package:" ] }, { "cell_type": "code", "execution_count": 1, "id": "dea44946", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:01:25.976164Z", "iopub.status.busy": "2024-02-29T12:01:25.975909Z", "iopub.status.idle": "2024-02-29T12:01:26.868492Z", "shell.execute_reply": "2024-02-29T12:01:26.867093Z" } }, "outputs": [], "source": [ "import ufl\n", "import dune.ufl\n", "import dune.grid\n", "import dune.fem" ] }, { "cell_type": "markdown", "id": "2662bb57", "metadata": {}, "source": [ "In our attempt we will discretize the model as a 2x2 system. Here are some possible model parameters and initial conditions (we even have two sets of model parameters to choose from):" ] }, { "cell_type": "code", "execution_count": 2, "id": "cafd16c4", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:01:26.873839Z", "iopub.status.busy": "2024-02-29T12:01:26.873458Z", "iopub.status.idle": "2024-02-29T12:01:26.880191Z", "shell.execute_reply": "2024-02-29T12:01:26.879315Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "dimRange = 1\n", "dt = 0.25\n", "linearSpiral = True\n", "\n", "if linearSpiral:\n", " spiral_a = 0.75\n", " spiral_b = 0.02\n", " spiral_eps = 0.02\n", " spiral_D = 1./100\n", " def spiral_h(u,v): return u - v\n", "else:\n", " spiral_a = 0.75\n", " spiral_b = 0.0006\n", " spiral_eps = 0.08\n", " def spiral_h(u,v): return u**3 - v\n", "\n", "x = ufl.SpatialCoordinate(ufl.triangle)\n", "initial_u = ufl.conditional(x[1]>1.25,1,0)\n", "initial_v = ufl.conditional(x[0]<1.25,0.5,0)" ] }, { "cell_type": "markdown", "id": "fa3c4ce1", "metadata": {}, "source": [ "Now we set up the reference domain, the Lagrange finite element space (second order), and discrete functions for $(u^n,v^n($, $(u^{n+1},v^{n+1})$:" ] }, { "cell_type": "code", "execution_count": 3, "id": "29cfa0cb", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:01:26.883507Z", "iopub.status.busy": "2024-02-29T12:01:26.883297Z", "iopub.status.idle": "2024-02-29T12:07:45.676066Z", "shell.execute_reply": "2024-02-29T12:07:45.675108Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "gridView = dune.grid.structuredGrid([0,0],[2.5,2.5],[30,30])\n", "space = dune.fem.space.lagrange( gridView, dimRange=dimRange, order=1 )\n", "\n", "uh = space.interpolate( initial_u, name=\"u\" )\n", "uh_n = uh.copy()\n", "vh = space.interpolate( initial_v, name=\"v\" )\n", "vh_n = vh.copy()" ] }, { "cell_type": "markdown", "id": "01ded985", "metadata": { "lines_to_next_cell": 0 }, "source": [ "We define the model in two steps:\n", "- first we define the standard parts, not involving $f_E,f_I$:\n", "- then we add the missing parts with the required _if_ statement directly using C++ code" ] }, { "cell_type": "code", "execution_count": 4, "id": "e9bdc53e", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:07:45.680052Z", "iopub.status.busy": "2024-02-29T12:07:45.679728Z", "iopub.status.idle": "2024-02-29T12:07:45.687917Z", "shell.execute_reply": "2024-02-29T12:07:45.687076Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "u = ufl.TrialFunction(space)\n", "phi = ufl.TestFunction(space)\n", "\n", "# right hand side (time derivative part\n", "a_ex = ufl.inner(uh_n, phi) * ufl.dx\n", "# left hand side (heat equation in first variable + backward Euler in time)\n", "a_im = (dt * spiral_D * ufl.inner(ufl.grad(u), ufl.grad(phi)) +\n", " ufl.inner(u,phi)) * ufl.dx\n", "\n", "ustar = (vh_n[0]+spiral_b)/spiral_a\n", "a_ex += ufl.conditional(uh_n[0]\n", " \n", " Your browser does not support the video tag.\n", "" ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQEAMgAyAAD/2wBDAAgGBgcGBQgHBwcJCQgKDBQNDAsLDBkSEw8UHRofHh0aHBwgJC4nICIsIxwcKDcpLDAxNDQ0Hyc5PTgyPC4zNDL/2wBDAQkJCQwLDBgNDRgyIRwhMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjL/wAARCADQANgDASIAAhEBAxEB/8QAHwAAAQUBAQEBAQEAAAAAAAAAAAECAwQFBgcICQoL/8QAtRAAAgEDAwIEAwUFBAQAAAF9AQIDAAQRBRIhMUEGE1FhByJxFDKBkaEII0KxwRVS0fAkM2JyggkKFhcYGRolJicoKSo0NTY3ODk6Q0RFRkdISUpTVFVWV1hZWmNkZWZnaGlqc3R1dnd4eXqDhIWGh4iJipKTlJWWl5iZmqKjpKWmp6ipqrKztLW2t7i5usLDxMXGx8jJytLT1NXW19jZ2uHi4+Tl5ufo6erx8vP09fb3+Pn6/8QAHwEAAwEBAQEBAQEBAQAAAAAAAAECAwQFBgcICQoL/8QAtREAAgECBAQDBAcFBAQAAQJ3AAECAxEEBSExBhJBUQdhcRMiMoEIFEKRobHBCSMzUvAVYnLRChYkNOEl8RcYGRomJygpKjU2Nzg5OkNERUZHSElKU1RVVldYWVpjZGVmZ2hpanN0dXZ3eHl6goOEhYaHiImKkpOUlZaXmJmaoqOkpaanqKmqsrO0tba3uLm6wsPExcbHyMnK0tPU1dbX2Nna4uPk5ebn6Onq8vP09fb3+Pn6/9oADAMBAAIRAxEAPwD3uaRooXdIXmZRkRxlQzew3ED8zWQ/ia3S30+4+xXphvHKeYsalIfmC5kbdgAk8EEggZGeK1L2zh1CymtLgM0MylHCOyEg9RlSCPwNVDoVi0EcLC5eKNgyo93KynGMAgt8w+UcHI68cnIBJqGs6XpOz+0dSs7Pf9z7TOse76biM1S/4THwv/0Mmj/+B0X/AMVXM/F5S3h3TAP+gkv/AKJlryqOBielezgMqjiqPtHO2ttjnrYiNLc98/4TDwx/0Mekf+B0X/xVKPF3ho9PEOkn/t9j/wDiq8StbKSQgKhJrqtK0OTKs6UsVl1HDq8qn4f8E8atn3LLlpw5n6/8A9GHijw+3TXdMP0u4/8AGnf8JLoP/Qb03/wLj/xrAhhjto+cCql7qKRqQGFfL4jMIU5Wgrnu5bHE43Vw5UdQfFHh8ddd0wf9vcf+NJ/wlfhwdfEGlf8AgZH/AI15reaiXPWsuW5LHrXOs0l/L+J9RT4fclrP8P8Agnrv/CW+G/8AoYNK/wDA2P8Axo/4S3w3/wBDDpP/AIGx/wCNeNGak82n/aj/AJfxNv8AVtf8/Pw/4J7N/wAJb4b/AOhg0r/wNj/xpR4r8OHpr+lf+Bkf+NeNLNg1ahuip60v7Uf8v4ky4cttP8P+CeujxP4fPTXdMP8A29x/40v/AAk2g/8AQb03/wAC4/8AGvOrHU9uMmujtbyOZByKqGaXdpRsePj8trYaPND3jfPirw6v3tf0sfW8j/xpp8X+GR18RaR/4Gx//FVzt9p63CEqAa4vUtHmhdj5Z219DgaNDFL47HxtTO6lGpyVKdvn/wAA9V/4THwx/wBDHpH/AIHRf/FUn/CY+F/+hk0f/wADov8A4qvCp7dlzxVJ0IPSvchw/CW1T8P+CerQxsKyuj6C/wCEx8L/APQyaP8A+B0X/wAVR/wmPhf/AKGTR/8AwOi/+Kr55IpKv/VuP/Pz8P8AgnVzH0P/AMJj4X/6GTR//A6L/wCKo/4THwv/ANDJo/8A4HRf/FV88UUf6tx/5+fh/wAEfMfQ/wDwmPhf/oZNH/8AA6L/AOKqSDxV4duriOC31/S5ppGCpHHeRszE9AADkmvnStfwt/yNukf9fkX/AKEKyr5BGlSlU9peyb27L1DmPouiiivmigooooA4n4m27XWkaVEgyTqIP/kGWuNsvDk0j/MMCvRvFuzbpHmY2/bj1/64S1kS6ha2ynDKKt508FT9jHfc8fGZRi8fiF7L4bWILHR4bNAXxmpLrUobVMKVrC1LXlkBWNyPpXOXF20rZ3k/jXzmKx9bEyvJn12TcH0cOlOrqzdvtcmdiFfA9qxJr6WQnL1TLknqabya47H29HCU6StFErTMerUwuTTkgL1bi0x5OmaDVzhHco7iaMmtmPw/Iw6mpf8AhG5fU0XRk8XRXUwdxFKJCO9bEmgSJ61Tl0x4+uaLoqNelLZkCXDr0NXrfU54yMPWc8JSmZxRYqVOE1sdtp+uMzBZHFbZ8i9iI+U15jHKynO4/nWxp2qm3fJc4z3NaUq06MuaDPms34bw+Mg9NTS1fw8cF4v5Vx95YtE5BGCK9MsdSgvIQrMCT71n6zoaTKZIxzX3OTcQc9qdV6n5Tj8txWT1bpe6eYyRlTURFbN7ZGKRlI5FZTpg19zSqqauj1cFjI14JoiopTSVsegFa/hb/kbdI/6/Iv8A0IVkVr+Fv+Rt0j/r8i/9CFcuO/3Wp/hf5DW59F0UUV+cmgUVU1OY2+mXEwN0CiFh9li8yX6KuDk/hXOS6hrZ0+yMTXhn8uTeRZEF5soY43ygwpVmDOAoyvVelAEXxIuDbaXpcg6/2gB/5BlrzubUWkHNdx8VWlGlaUAieT9vGX3ndu8mXA24xj3z+FebZJrwsxX775H2nD1KMsK5NdX+hI8hao+TTkRmPAzWtYaTNcFcR8VwHuzqQpK8jOht2lYDaT+FblhoJmwWjOPcV0dlpMVqgeUKCPWotQ1+3sVKxAFh6V24TL6+Llamj4rOOMKGFfJT1fkPt9AtIVBdB+NTO2m2g5VOPauLvvFl1JkA7R7Vhz6zPK2WJP419ZhOEW1esz4yvnOcY13p+6vU9GfxDpsRwAOPQVH/AMJTYf3f/Ha80a/kam/bJK9mPCuES1RyfV8zlrKq/vPUo9f02YgFRz6irSrpt4OETn2ryVdQkWrlvrlxCflYj8a5cRwlQkv3ehUKmc4d3hUv8z0S78PW0iExov4Cuav9HaBuIzj2Wnaf4tuEwJMMPrXW2l/a6jCMldx7V8pj8gxOE95ao+gy3jOpSmqWMi0zziSJkbBUj8KarEGu21XQy5Z4lGMVyNzaSQs25cYNeFqtGfpGDx1LFQUoO5YstQa3K4zxXbabqaXcao3UjvXnIJFaenX8kNwhB4pqTi+aJy5pllLGUnGSOj13RlkVpY159q4K7tihb5eleqWV0l7AVbBOK5TxDpRi3SKvB61+gcPZv7Rezm9T8Vx2Eq5Ri+V/CzhHXFMq5cxFTVQivuYSuj6DD1VUimhK1/C3/I26R/1+Rf8AoQrIrX8Lf8jbpH/X5F/6EK58d/utT/C/yOhbn0XRRRX5yaBRRRQBwfxV/wCQFpf/AGER/wCiZa82iiMjYFenfE2LzdI0pPXUR/6Jlrl9J0YSTAspxXhZk/33yPsMlxMKGBcpPq/yRFpGivP8x6V1JNtpVtl2AIHem3FxbaPaE8AgV5/reuNeyE7jjsK9XJcini5KdRe6fD57xFXxtZ4bCbdWaWteJJJpCsUm1PbvXJz3ssjMS5OagkkLtnNR1+l4XBUsPBRgjzsHlkKXvT96XmOLsepptFFdp6iSWwUUUUDClzSUUCHrK69DWlY6nPDKpEhGD2rKpQcGs6lKM1Zo5cRg6VeNpI9U0TXo7tRDMw3Y70/V9IWaF5Iz6nivNrC7MEwbcRjHevSPD+rpqFuYZCC3Svz/AIhyFQvXoo48szGvk2JVOTvTfU4u6tjAcHNQoxVgQa63xDpiou9V4JrkpV2Pivh/Jn7Lg8TDE0lOL3Oi0G/Mc5DOcHHeuov4UvrBiOeK88tJfLkzXeaFci4s9p9xXTg67oVlNHynF2UxxGGlNLVHn+q2hglKkdDWJIuGrvfFNkEcOB1NcRcptc1+v5ZiVWoqR+cZPXavSluipWv4W/5G3SP+vyL/ANCFZJrW8Lf8jbpH/X5F/wChCurHf7rU/wAL/I+kR9F0UUV+cmgUUyWaOCF5ppEjiQbmd2ACj1JPSoH1OwjW2Z762VbrHkFpVAlzjG3n5uo6etAGL4vhWdNIRhkfbs/+QJay7ieDTrdn4BAq549u/sNjpc/pfY/8gS15nqmtTXOQX+X0rfCZP9crKq9keJmuYYlP6nR2av8AeR63rcl5I+fu54Fc47lzzTpZGdiSair9Bw+HhQgoxRvgcFHDw82FFFFdB6AUUUUAFFFFABRRRQAUUUUAKDg1uaJqT2kgZSeKwqlhkKdKxr0o1YOLOHH4ZYii49T2AFdT0tXxkkZrhtWtzDdlcV03hK782zETHOKo+JrbbdBgOor8dzjC/VsVKK2PoOB8xlOm8PN6xOXVsNXU+HLooNpPeuVbhyK1tJlKEc15bPvcbSVSk4s63X7cT2YbFea30e2UivVp187TgfavM9Wj23Lj3r9I4XxDnR5X0PwepD6tmkorzMNuprV8Lf8AI26R/wBfkX/oQrLccmtTwt/yNukf9fkX/oQr6nG/7rU/wv8AI+ng7pH0XRRRX50bDJRIYnETKkhB2s67gD2JGRn8xXPL4ZuG0XT9Onu7ORrJURLgWREiqoUZQmQ7H+U/Nz1HHHPSUUAeefFtFXRtLn3OHF+Exvbbgwy/w5xn3xmvJpJSx616z8YP+Rb0z/sJL/6JlryEmvs8givqrfm/0OedKLqcwlFFFe6WFFFFAwpQM0gqzbwl3AxSlKyuZVaqpxcmRpEW7VYWzLfwmuk0fQnumBKYX1rrrfQLSFQXRc14GOz2hhXaT1PAeY4jET5cPG55a1mR/CarvEVPSvWLrw/ayRsUQZx2rgtWshbTsnpWuX5xSxjtAdHMa0KqpV42bOfopzjDGm17Z9EndXClBpKKAO48JzFZkXPXFbviRMlW/wBmue8KRlriMgdK6LxIwG0Z/hr8v4sSVdWOXg66zKqltqcJJxK31q7YnBFUpOZW+tXbEZNfJs/YKvwHocPzaaP92vOtcXF7J9a9Fg+XTV/3a86105vZPrX3PCd7M/Cc5t/azt3Odf7xrT8Lf8jbpH/X5F/6EKzH+8a0/C3/ACNukf8AX5F/6EK+4xv+61P8L/I96lsj6Looor86NwooooA88+MP/ItaZ/2El/8ARMteQV6/8Yf+Ra0z/sJL/wCiZa8gr7Th/wD3R+r/AEIluFFFFe4IKKKKAHoMsK6LQ7A3F2i446msCAZkX616D4Tth5nmEdq8rNsT7ChKS7Hz2c1ZNxor7R0WI9NssgAECuWu9duXPySAVo+IrvCugPSuRgBc9a/Ha9aVao5yZ+ncO5PRwuFTa1PQtKlnktVkmIxiuI8TzxNqMmzsMV2V1L9j0RMHb8teW6nM0l2xLE5r7fhHBOzrN7n57ms1jc2dOKsoXKTHLGkoor9APUSsrBT40LHimgZrW0iwa7kCqCSTWdWoqcHJnLjMTHD0nNnaeELMrAJTTPE1yDchQegrcgjj0rS1U8HFcRrF1594WzX47nWL+s4qTWyPa4Hy6UU8TNayMw8uT71q6XEXIwO9ZajJrp/Dtr5nJHevJZ+g4yoqdJtnUSnydOGfSvNNXk3XUh969A164FvZBc9q80v5N0pNfo/C+HcKPM+p+D1p/Wc0lJbamc/U1qeFv+Rt0j/r8i/9CFZTda1fC3/I26R/1+Rf+hCvqsb/ALrU/wAL/I+mgrJH0XRRRX50bDJTIIXMKI8oHyq7FQT7kA4/I1gzeJ2j0mwvksg32ixa/kQzY8uNVQsAdvzN84wOM4PIrelijnieKWNZI3BVkcZDA9QR3qpHomlRQRwR6ZZJFG4kSNbdQqsAAGAxwcADPtQBw3ximUaDpcOH3HUA2dh248mX+LGM+2c15JXr/wAYf+Ra0z/sJL/6JlryCvtOH/8AdH6v9CJbhRRRXuCCiiigCxbf6xfrXpHhQ/J+FeawnEi/Wu68K3QW5CE9RXg59Sc8NJI+ZzdctenU6Jk/iJT5shwa5mB9ld1r1oZIHdR2rg3Vk61+Q2s2mftOUV418LFx8j0aJU1DS1Q8nbXm2taZLBeuvltge1dHo+sSQOiEjbnmunlt7TU4ssFLEda+lyHO/qMuSpsfnfEGSYjCYp43DK990eONFICQUNCxSE/cNeiXHhB2djGy4PQU+18IlGBlYYHWvuXxBhFDm5keB/a9b4VSdzi9O0ya5mCiImvQ9F0eLTLbzZVAfrV2K3s9NjyoXIrC1bXXZGSMgD2r4/OeJJYlOlR0R6uWZBi8zrKtiVaK6CeINXEgEadM1ykr+Y+TT553lOWOabGhdgAK+R82frWEw0MNSUIKyRLaxGSTGK7vQbUW9nuIx3rD0HTWebcyfLXR6lcpY2LAcHFdWDw8sRWUEfI8W5xHDYd04vVnOeKb4O4QHpXEXL7nNaep3bTSFietYsjZav17LcMqFJRPzrKKDd6st2MNa3hb/kbdI/6/Iv8A0IVkVr+Fv+Rt0j/r8i/9CFdWO/3Wp/hf5H0iPouiiivzk0CiiigDzz4w/wDItaZ/2El/9Ey15BXr/wAYf+Ra0z/sJL/6JlryCvtOH/8AdH6v9CJbhRRRXuCCiiigByHDCt7R7wwXSPnpXPirNvLsfOawxFJVIOLPNzHC+3pNdT1+3mTUbEqCCSK5DV9Je2BPak8P60LZ9rt8pNdfLDb6tbA8HPNflWd5VPDVXOK0OvhbiB4SX1XEOx5vl42GCQRV221OeFhiZh+NaeqaIIH+VeKwZbcoxGK+f3P1anVpYiF1qjoYvEUqgZmP50k3iKVgQJjXMFPajZRYj+z8PfmsvuNG51KaYnMrH8aoM7N1JNCxk1bgszKcYNGx0pQprQrJGXrc0vR3mlRj0q/pvh9ZArODiugP2bTYewIFa0aM60uWCPm864ioYKm9dRUSLTrckkZxXE67qxuGZQ3yg1Y1rXfPZkRvlrkri43k81+jZFkqw8VOa1PySpVr5tifa1PhIp5Sx61WJpWOabX18Y2R9HRpKnGyCtfwt/yNukf9fkX/AKEKyK1/C3/I26R/1+Rf+hCufHf7rU/wv8jdbn0XRRRX5yaBRVLV47qXSriOzLCcr8oRtrEZ5CtxgkZAORgnqK5+Wy1iXT7QKmqx3VvdSNGBdpsaMvuQTHfudQpCnBJyrdeCQDK+MAz4b0z/ALCS/wDomWvICK9g+LEU0ukaWA6eT9vGU2HcW8mXnOcY9sfjXlE0JU19jkE0sM15v9DkqV4xq+zb1KtFB4NFe+bBRRRQMKVTg0lFAmrly3uShrqNI8QSWwVTkrXGA4qeOdlxg1xYrB068eWSuePjstjV9+OjPXLXU7a+jG4DJ9adPpVtPllVefSvMrbUpYiCr4rbtfEt1HgFwRXxeO4WUm5UtDDC5tmWXaP3kb8vh4c7UqMeHjn/AFdQReK2/ixUp8VfSvFfDWKTsexHjuslZwLsGgRqQWQVfjs7S1XJVRXMz+KpSDtIrIuvEFzLnMmPpXZhuFasneozhxPFePxS5aUbHY32uwWiEIOfauQ1PXZLljzgelY09/I5JLZqlJMzHrX1+AyWjhkrLU8yngK+Jnz4mVyae5Lk1UZs0hYmm17sYKK0PoKGHjSjZIKKKKs6ArX8Lf8AI26R/wBfkX/oQrIrX8Lf8jbpH/X5F/6EK5cd/utT/C/yGtz6Looor85NAooooA5bxzaLe2mlQsMg32f/ACBLXmmr6MbUk7flr1HxbKIho7Hp9uI/8gS1lXtrFfWrDAJIrbCZs8JXVJ7M8HN8DXusZS6afceNSptdhjvUVburae9vK4KkYJrDIIr9Eo1VUgpI7MBi44ilzISiiitjuCiiigAooooEODEd6eJWHc1FRSaTIlTjLcsC4YfxH86X7S394/nVaip5EZPDU30JzOx/iNMMhPeo6KaikXGjCOyHFvem5ooqjVJIKKKKBhRRRQAVr+Fv+Rt0j/r8i/8AQhWRWv4W/wCRt0j/AK/Iv/QhXLjv91qf4X+QLc+i6KKK/OTQKKiubmGztpLid9kUYyxwT+QHJPsOtUn1/TY7aKd52EcpcD90+V2HDlhjKhTwS2AO+KAOe+JEph0zSXBwRqI/9Ey1j6LqRlfY7Z4rQ+KZxomln/qIj/0TLXFaROUuRyfzrwsx0rJrsfT5dhIYjLpKS6v9DovE2mrNaNKi8gZrza6gMZr1/i405gecivNddt/KmZcd6+74Xx8qtLkl0PzKUHgMydBfC2c/RStwaSvsj6FO6CiiigYUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFa/hb/kbdI/6/Iv/AEIVkVr+Fv8AkbdI/wCvyL/0IVy47/dan+F/kC3Pouiiivzk0K99Zx39nJbSsyq4HzIfmUg5BGe4IB/CstvC9tJb+S95dsCZg7Epl0mIaVDhejMM8YIzwQOK3KKAOXu7zQ9YxdXmltO9m0bwCWKNjMJWeJNuSeGbOA2OxOKvWOlaDf2i3Eei2SAsyMj2se5WVirKcAjIII4JHFTL4d0xTNtgkHnDDDz5MABtw2jd8uGyVxjaScYyavWtrDZW6wQIVjXJALFiSTkkk8kkkkk8kmgCn/wj2if9AfT/APwGT/Cj/hHtE/6A+n/+Ayf4VpUUAYd/puhWEKOdCtJnkkEccUVtFudj2G7A6AnkjpWbDH4aF5ZhPD6tFqiJPHctZp5KbkG1ST90kKPlHc+9dNeWUF/B5NwrFQwYFHZGUjoQykEH6Gq0eh2EQtVWOXZa48mM3EhRcEkfKWwcZ4yOMDGMDAAn/CPaJ/0B9P8A/AZP8KP+Ee0T/oD6f/4DJ/hWlRQBm/8ACPaJ/wBAfT//AAGT/CqFzpGgzXx0iTQ7dDcW0kizJbxgYBVWAI+YMN69vxroapT6TaXOpQ6hJ5/2mFCiMlxIihSwYjaGCnJVc5HOMHigDJ0u38O6tu8nQrWPEaTJ5trF+8ifOxxjPB2twcHjkCtH/hHtE/6A+n/+Ayf4VNY6XZ6bv+yxFN+AcuzYAzhRknCjJwowBk4FXKAM3/hHtE/6A+n/APgMn+FB8P6KAT/Y1gfYWyf4VpUHkYoA5CePw/Hp0es/8Iw0gR5bZreGzid0IcqxYA7cAxnnPQ+9bcehaFLEkiaTpzI6hlIt0IIP4UkXh3T4bSS1jN4sMjl2H26fJJLMed+RksxPr36CtOONIYkijUJGihVVRgADoBQBn/8ACPaJ/wBAfT//AAGT/Cj/AIR7RP8AoD6f/wCAyf4VpUUAYGp2Gh6XbC4k8P200QIDmK2h+TJwOGIJ5OMLknsKrwWPhyLU30xdEgJ85k857aMp5pXzSg/i+62RxjAxnjFbd/pdrqTW7XKyk28nmR+XO8eGxjPykZ4J6+tNh0eygvFu0jczqoG55nbJ2hdxBJBbaANx+YjjNAEf/CPaJ/0B9P8A/AZP8KP+Ee0T/oD6f/4DJ/hWlRQBm/8ACPaJ/wBAfT//AAGT/Csue00RLyaD/hHhHPaILyNoLeLdJ5bBhs2nIJIwA23POOhx01Zy6JZJd3F0n2lZrlt0pW7lAY7dvTdgYHTHTgjkCgCXS9Sj1WxW6jimhyzK0U6hZEIJBDAE4PHSipLOygsIDDbqwUsXJd2dmYnJJZiST9TRQB//2Q==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "import matplotlib.pyplot as plt\n", "from matplotlib import animation\n", "\n", "from numpy import linspace\n", "fig, ax = plt.subplots()\n", "ax.set_xlim(( 0, 2.5))\n", "ax.set_ylim(( 0, 2.5))\n", "triangulation = gridView.triangulation(1)\n", "levels = linspace(-0.1, 1.1, 256)\n", "ax.set_aspect('equal');\n", "t = 0.\n", "stepsize = 0.5\n", "nextstep = 0.\n", "anim = animation.FuncAnimation(fig, animate, init_func=init, frames=20, interval=100, blit=True)\n", "\n", "try:\n", " movie = anim.to_html5_video()\n", " from IPython.display import HTML, display\n", " display( HTML(movie) )\n", "except: # ffmpeg probably missing\n", " try:\n", " anim.save(\"spiral.html\")\n", " from IPython.display import IFrame\n", " IFrame(src='./nice.html', width=700, height=600)\n", " except:\n", " pass" ] }, { "cell_type": "markdown", "id": "035fa12b", "metadata": {}, "source": [ "... if the movie is not showing you might have to rerun the notebook ..." ] } ], "metadata": { "jupytext": { "cell_metadata_filter": "-all" }, "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.10.12" } }, "nbformat": 4, "nbformat_minor": 5 }