{ "cells": [ { "cell_type": "markdown", "id": "1df232ee", "metadata": {}, "source": [ ".. index:: Equations; Linear Elasticity\n", "\n", "# Linear Elasticity: a deformed beam\n", "We first setup the domain and solution space, together with the\n", "vector valued discrete function $u$ which describes the displacement\n", "field:" ] }, { "cell_type": "code", "execution_count": 1, "id": "f845848b", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:01:25.979943Z", "iopub.status.busy": "2024-02-29T12:01:25.979507Z", "iopub.status.idle": "2024-02-29T12:07:45.435394Z", "shell.execute_reply": "2024-02-29T12:07:45.434051Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "from matplotlib import pyplot\n", "from dune.fem.plotting import plotPointData as plot\n", "from dune.grid import structuredGrid as leafGridView\n", "from dune.fem.space import lagrange as solutionSpace\n", "from dune.fem.scheme import galerkin as solutionScheme\n", "\n", "from ufl import *\n", "import dune.ufl\n", "\n", "gridView = leafGridView([0, 0], [1, 0.15], [100, 15])\n", "space = solutionSpace(gridView, dimRange=2, order=2, storage=\"istl\")\n", "displacement = space.interpolate([0,0], name=\"displacement\")" ] }, { "cell_type": "markdown", "id": "fc7c59a3", "metadata": { "lines_to_next_cell": 2 }, "source": [ "We want clamped boundary conditions on the left, i.e., zero displacement" ] }, { "cell_type": "code", "execution_count": 2, "id": "68035825", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:07:45.442870Z", "iopub.status.busy": "2024-02-29T12:07:45.441933Z", "iopub.status.idle": "2024-02-29T12:07:45.447693Z", "shell.execute_reply": "2024-02-29T12:07:45.446635Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "x = SpatialCoordinate(space)\n", "dbc = dune.ufl.DirichletBC(space, as_vector([0,0]), x[0]<1e-10)" ] }, { "cell_type": "markdown", "id": "09ec3520", "metadata": { "lines_to_next_cell": 2 }, "source": [ "Next we define the variational problem starting with a few constants\n", "describing material properties ($\\mu,\\lambda,\\rho$) and the gravitational\n", "force" ] }, { "cell_type": "code", "execution_count": 3, "id": "6858fe94", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:07:45.452802Z", "iopub.status.busy": "2024-02-29T12:07:45.452056Z", "iopub.status.idle": "2024-02-29T12:07:45.456938Z", "shell.execute_reply": "2024-02-29T12:07:45.455893Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "lamb = 0.1 # 1st Lamé coefficient (\\lambda)\n", "mu = 1 # 2nd Lamé coefficient (\\mu)\n", "rho = 1/1000. # material density\n", "g = 9.8 # gravitational acceleration" ] }, { "cell_type": "markdown", "id": "4fea81b3", "metadata": { "lines_to_next_cell": 2 }, "source": [ "Next we define the strain and stress\n", "\\begin{align*}\n", "\\epsilon(u) &= \\frac{1}{2}(\\nabla u + \\nabla u^T) \\\\\n", "\\sigma(u) &= 2\\mu\\epsilon(u) + \\lambda\\nabla\\cdot u I\n", "\\end{align*}\n", "where $I$ is the identity matrix." ] }, { "cell_type": "code", "execution_count": 4, "id": "b6b9803b", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:07:45.461861Z", "iopub.status.busy": "2024-02-29T12:07:45.461122Z", "iopub.status.idle": "2024-02-29T12:07:45.466391Z", "shell.execute_reply": "2024-02-29T12:07:45.465372Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "epsilon = lambda u: 0.5*(nabla_grad(u) + nabla_grad(u).T)\n", "sigma = lambda u: lamb*nabla_div(u)*Identity(2) + 2*mu*epsilon(u)" ] }, { "cell_type": "markdown", "id": "d826d017", "metadata": { "lines_to_next_cell": 2 }, "source": [ "Finally we define the variational problem\n", "\\begin{align*}\n", "\\int_\\Omega \\sigma(u)\\colon\\epsilon(v) = \\int_\\Omega (0,-\\rho g)\\cdot v\n", "\\end{align*}\n", "and solve the system" ] }, { "cell_type": "code", "execution_count": 5, "id": "bb769ccf", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:07:45.471279Z", "iopub.status.busy": "2024-02-29T12:07:45.470498Z", "iopub.status.idle": "2024-02-29T12:08:54.487048Z", "shell.execute_reply": "2024-02-29T12:08:54.485384Z" }, "lines_to_next_cell": 2 }, "outputs": [], "source": [ "u = TrialFunction(space)\n", "v = TestFunction(space)\n", "equation = inner(sigma(u), epsilon(v))*dx == dot(as_vector([0,-rho*g]),v)*dx\n", "\n", "scheme = solutionScheme([equation, dbc], solver='cg',\n", " parameters = {\"newton.linear.preconditioning.method\": \"ilu\"} )\n", "info = scheme.solve(target=displacement)" ] }, { "cell_type": "markdown", "id": "371d8410", "metadata": { "lines_to_next_cell": 2 }, "source": [ "We can directly plot the magnitude of the displacement field and the stress" ] }, { "cell_type": "code", "execution_count": 6, "id": "88449f12", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:08:54.497959Z", "iopub.status.busy": "2024-02-29T12:08:54.497658Z", "iopub.status.idle": "2024-02-29T12:09:35.440737Z", "shell.execute_reply": "2024-02-29T12:09:35.439809Z" }, "lines_to_next_cell": 2 }, "outputs": [ { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQEAMgAyAAD/2wBDAAgGBgcGBQgHBwcJCQgKDBQNDAsLDBkSEw8UHRofHh0aHBwgJC4nICIsIxwcKDcpLDAxNDQ0Hyc5PTgyPC4zNDL/2wBDAQkJCQwLDBgNDRgyIRwhMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjL/wAARCACFA+MDASIAAhEBAxEB/8QAHwAAAQUBAQEBAQEAAAAAAAAAAAECAwQFBgcICQoL/8QAtRAAAgEDAwIEAwUFBAQAAAF9AQIDAAQRBRIhMUEGE1FhByJxFDKBkaEII0KxwRVS0fAkM2JyggkKFhcYGRolJicoKSo0NTY3ODk6Q0RFRkdISUpTVFVWV1hZWmNkZWZnaGlqc3R1dnd4eXqDhIWGh4iJipKTlJWWl5iZmqKjpKWmp6ipqrKztLW2t7i5usLDxMXGx8jJytLT1NXW19jZ2uHi4+Tl5ufo6erx8vP09fb3+Pn6/8QAHwEAAwEBAQEBAQEBAQAAAAAAAAECAwQFBgcICQoL/8QAtREAAgECBAQDBAcFBAQAAQJ3AAECAxEEBSExBhJBUQdhcRMiMoEIFEKRobHBCSMzUvAVYnLRChYkNOEl8RcYGRomJygpKjU2Nzg5OkNERUZHSElKU1RVVldYWVpjZGVmZ2hpanN0dXZ3eHl6goOEhYaHiImKkpOUlZaXmJmaoqOkpaanqKmqsrO0tba3uLm6wsPExcbHyMnK0tPU1dbX2Nna4uPk5ebn6Onq8vP09fb3+Pn6/9oADAMBAAIRAxEAPwD3e+vI7CzkuZVZlQD5UHzMScADPckgfjVMa5E1taz/AGW52z3JtSPkzE4cod3zdNykZXP5VfubaG8tpLedN8UgwwyR+RHIPuKoHw/pxt44NlwI45/tCgXUoPmZznO7J55IPBJJ6kmgB+m6xHqV1e2wtbq3ltJCjC4QIXG5lDqMklDtOD3/ADqrceMfC9pO8Fz4k0eGZDho5L6JWU+4LZFaNtp1taXE9xEJDNPgSPJK8hIBJAG4nABZuBgcmvlT4y/8la136wf+k8ddWDwyxFTkbsTJ2Vz6W/4Tnwj/ANDVof8A4MIv/iqP+E58I/8AQ1aH/wCDCL/4qvi6lr1lkkf5/wAP+CR7TyPtD/hOfCP/AENWh/8Agwi/+Ko/4Tnwj/0NWif+DCL/AOKr4woq1kUX9v8AD/gi9r5H2f8A8Jz4R/6GrRP/AAYRf/FUf8Jz4R/6GnRP/BhF/wDFV8Y0tWuH4v8A5efh/wAEPa+R9m/8Jz4R/wChp0T/AMGEX/xVH/Cc+Ef+hp0T/wAGEX/xVfGdLWi4bi/+Xn4f8EXtvI+y/wDhOfCP/Q06J/4MIv8A4qj/AITnwj/0NOif+DCL/wCKr41pwrWPC8H/AMvfw/4Ivb+R9k/8Jx4R/wChp0T/AMGEX/xVH/CceEv+hp0T/wAGEX/xVfG4pwraPCUH/wAvX93/AARe3fY+xv8AhOPCX/Q06J/4MIv/AIqj/hOPCX/Q06J/4MIv/iq+OhS1suDYP/l8/u/4IvrD7H2J/wAJv4S/6GjRP/BhF/8AFUf8Jv4T/wCho0X/AMD4v/iq+P0PNWENJ8HQX/L5/d/wSZYprofXP/CbeFP+hn0X/wAD4v8A4ql/4TXwp/0M2jf+B8X/AMVXychqzGayfCUF/wAvfw/4JjLHtfZPqj/hNPCv/QzaN/4Hxf8AxVL/AMJn4W/6GXR//A6L/wCKr5fQ1ajPSsnwtFf8vfw/4JhLNWvs/ifS3/CZeF/+hk0f/wADov8A4ql/4THwx/0Mekf+B0X/AMVXzlGatRmspcNxX/Lz8P8AgmMs6kvsfj/wD6D/AOEw8M/9DHpH/gbH/wDFUv8Awl/hn/oYtJ/8DY//AIqvBozVuM1jLIIr/l5+H/BMJcQSX/Lv8f8AgHt3/CW+Gv8AoYdJ/wDA2P8A+Kpf+Es8N/8AQwaV/wCBsf8AjXjcZq3GaxlkiX2/w/4Jzy4mmv8Al1+P/APWv+Er8Of9B/Sv/AyP/Gj/AISrw7/0HtL/APAyP/GvMIz0q5EeawnlSj9r8DCXFs1/y6X3/wDAPRP+Ep8Pf9B7S/8AwMj/AMaP+Eo8P/8AQd0z/wAC4/8AGuIjPFTCuGeG5eo48WTf/Lr8f+Adj/wlHh//AKDumf8AgXH/AI0v/CT+H/8AoO6Z/wCBcf8AjXICniuWfunTDiWUv+Xf4/8AAOs/4Sfw/wD9BzTP/AuP/Gj/AISfQP8AoOaZ/wCBcf8AjXLiniuSeJceh0wz2Uvsfj/wDpf+En0D/oOaZ/4Fx/40f8JNoH/Qc03/AMC4/wDGudWpBXHPM3H7P4nVDNHL7P4m/wD8JNoH/Qb03/wLj/xo/wCEl0H/AKDem/8AgXH/AI1hiniuSeeSj9j8f+AdUMY5dDZ/4SXQf+g3pv8A4Fp/jR/wkug/9BvTf/ApP8ayhTxXJPiacf8Al3+P/AOqFTmNL/hJNC/6DWm/+BSf40f8JLoP/Qb03/wKT/GqApr1yy4vmn/CX3/8A6qdJT6mj/wk2g/9BvTf/AuP/Gk/4Sfw/wD9BzTP/AuP/GsWXvVOUVUeLZv/AJdL7/8AgHZHAKX2jpf+Eo8Pf9B3TP8AwLj/AMaP+Eq8O/8AQe0v/wADI/8AGuOlHWqUveto8USf/Lr8f+AdEcpUvt/gd7/wlfhz/oP6V/4GR/40n/CWeG/+hg0r/wADY/8AGvN5BVOUda3jxFJ/8u/x/wCAdEcijL7f4f8ABPVP+Et8Nf8AQw6T/wCBsf8A8VSf8Jf4Z/6GLSf/AANj/wDiq8fl71TkFbRz2T+x+P8AwDpjw3GX/Lz8P+Ce1/8ACYeGP+hj0j/wOj/+KpP+Ex8L/wDQyaP/AOB0X/xVeESjrVOWt45u39j8TojwpCX/AC9/D/gn0H/wmfhYf8zLo/8A4HRf/FUn/CaeFf8AoZtG/wDA+L/4qvm+ccVnvW8cxcvsly4Rgv8Al6/u/wCCfT//AAmvhT/oZtG/8D4v/iqP+E28Kf8AQz6L/wCB8X/xVfLD1XetVjG+hyT4ajH/AJefh/wT6v8A+E38J/8AQ0aL/wCB8X/xVaen6nYatbfadOvra8g3FfNt5VkXI6jKkjPNfGz19F/BD/kn3/b5L/Ja6KVbndrHk4/LVhYKSlfW2x2+razBo8SvNFNLlXkYRAEqiDLOckcAY6ZPPANJd6wtncXMJsruV4LVrpRCqsZgvVUG7JbOAAQM561NqGl2eqRrHeRGRVzjDsuQRgg7SMgjqDwe9Rtotm1ytwTdecsH2cOLuUHZ/wB9dePvdfetzyh1nqttd6UuosTbQHdv+0EIYyrFSG5wCCCCM8VRXxn4WYZXxLo5HtfRf/FVq2dnBYWy29upWMFm+ZixJYliSSSSSSTk+teS6N4RU6NYOsYANvG3T1UV5eZ5nDARi5K9ztwWHo1nJVZ8tvI9IHi7wyeniLSf/A2P/wCKpf8AhLPDf/QwaV/4Gx/41x0XhhV6qBVyPQYk64rw3xXHpTv8/wDgHTPBYWP/AC9/D/gnTf8ACV+Hf+g/pX/gZH/jS/8ACU+Hv+g9pf8A4GR/41gppVunpUy2tsn8Oan/AFqqP4aP4/8AAOOcMNH/AJefgbP/AAlHh/8A6Dumf+Bcf+NH/CTaB/0HNN/8C4/8ayh5K9EFHnKOgFaR4jxMtqH4/wDAOSdbDx+0a3/CS6D/ANBvTf8AwKT/ABo/4STQv+g1pv8A4FJ/jWObj3phuPeumGdYqX/Llff/AMA5J4+hHqbf/CS6D/0G9N/8Ck/xo/4SXQf+g3pv/gWn+NYBuPemG4967IZhiZf8uvx/4ByzzejE6L/hJtB/6Dem/wDgXH/jSf8ACTaB/wBBzTf/AALj/wAa5sz+9MM/vXXCviJfY/H/AIByzz+lHp+J0/8Awk2gf9BzTf8AwLj/AMaP+En0D/oOaZ/4Fx/41yhn96YZ/euuCry+yc0+JoR+z+P/AADrv+En8P8A/Qc0z/wLj/xo/wCEn8P/APQd0z/wLj/xrjWuVHVhUL3sY/irrhhK0uhhLitLan+P/AO4/wCEo8P/APQd0z/wLj/xpP8AhKPD/wD0HdM/8C4/8a4B9SQdKryap6GuqGWVZbkf62Te1H8f+Aej/wDCUeHv+g7pn/gXH/jSf8JV4d/6D2l/+Bkf+NeXSann+Kqsmpf7VdMclk95fgVHies/+XP4/wDAPW/+Er8Of9B/Sv8AwMj/AMaT/hLPDf8A0MGlf+Bkf+NeNyaj/tVWk1D3rdZC39v8P+CbR4irP/l1+P8AwD23/hLfDf8A0MOk/wDgbH/jSf8ACXeGv+hi0n/wNj/+Krwl7/3qu9971quHb/b/AA/4JtHPar/5d/j/AMA99/4S/wAM/wDQxaR/4Gx//FUf8Jh4Y/6GPSP/AAOi/wDiq+eXvfeqz3vvWi4Zv/y8/D/gm0c4qP8A5d/j/wAA+jv+Ex8L/wDQyaP/AOB0X/xVH/CZeFv+hk0f/wADov8A4qvml7z3qu9371quFk/+Xn4f8E2jmlR/Y/H/AIB9O/8ACZ+Fv+hl0b/wOi/+Ko/4TTwr/wBDNo3/AIHxf/FV8tPde9QPde9aLhJP/l7+H/BNVj5v7H4n1X/wmvhT/oZtG/8AA+L/AOKo/wCE28Kf9DPov/gfF/8AFV8mPce9QPce9aLg+L/5ev7v+CarGSf2T66/4Tbwn/0NGi/+B8X/AMVSf8Jv4S/6GjRP/BhF/wDFV8gNMTUZYnvVrg2P/P78P+CaLEvsfYf/AAnHhL/oadE/8GEX/wAVR/wnHhL/AKGnRP8AwYRf/FV8d5pM0nwdTX/L5/d/wR/WH2PsX/hOPCX/AENOif8Agwi/+Ko/4Tjwl/0NOif+DCL/AOKr45zTc1m+Eqa/5fP7v+CP277H2R/wnHhH/oadE/8ABhF/8VSf8Jx4R/6GnRP/AAYRf/FV8b5pM1k+FoL/AJe/h/wR+3fY+yf+E58I/wDQ06J/4MIv/iqP+E58I/8AQ06J/wCDCL/4qvjXNJms3wzBf8vfw/4I/beR9l/8Jz4R/wChp0T/AMGEX/xVH/Cc+Ef+hp0T/wAGEX/xVfGmaSofDkF/y8/D/gh7byPsz/hOfCP/AENOif8Agwi/+Ko/4Tnwj/0NWif+DCL/AOKr4ypKzfD8V/y8/D/gj9r5H2d/wnPhH/oatE/8GEX/AMVR/wAJz4R/6GrRP/BhF/8AFV8Y0lS8hj/z8/D/AIIe18j7P/4Tnwj/ANDVof8A4MIv/iqP+E58I/8AQ1aH/wCDCL/4qvjCkqHkcf5/w/4I/a+R9of8Jz4R/wChq0P/AMGEX/xVH/Cc+Ef+hq0P/wAGEX/xVfF9IaX9iR/n/D/gh7TyPtH/AITnwj/0NWh/+DCL/wCKo/4Trwj/ANDVof8A4MIv/iq+LqKX9ir+f8P+CHtD7R/4Trwj/wBDVof/AIMIv/iqP+E68I/9DVof/gwi/wDiq+LqKX9jR/n/AA/4Ie0PtH/hOvCP/Q1aH/4MIv8A4qtyKWOeFJoZFkikUMjochgeQQR1FfCVfb3h7/kWtK/684f/AEAVwY7BLDctne9y4yuaVFFFcBQUVU1K9/s/T5boR+YyYCpu2gsSAMnsMkZPaqK63M1paT/Y4x5t41pMvnH92wkaMlfl+YZXvt4/KgDZr5K+Mv8AyVrXfrB/6Tx19P6VrTaleXEDWwiWMFo2Em4sokdPmGBtOUzjnr7Gvl34wSrP8VdckQOAWhGHQoeIIx0IBr08p/3j5Miexw9LSUV9QjAWlpKBWkWA6gUlLW0WIUU4U2lrogxDqWminCumDJFFOFNFKK7IMljhTqaKcK7IPQliqeasIaripozTkRNFtDVmM1UQ1ZjNc80ck0XENWkNU0NWozXLNHHNFyM9KtRmqcZ4q1Ga5Zo4qiLsZq3GapR1bjrkmcNRF2M1bjNUoz0q3Ga5JnDURdiNXIz0qlFVuOuKqtDhqF+I1OKrwmrArxq61FDckFPFRinivKrI9CkyQU8UwU8V5VZHpUmSLTxUYqQV5NdHp0WPFSCoxUgryqyPSpMeKkFRinivJrI9OkSCkegUrV5dXc9Sg9SrJVKSrstU5R1qqZ69MpS96pyVdkqlLXZTO+kUpKpy96uyVTl6120z0KZRlqnJV2WqUldlM9CkUpe9U5e9XZe9U5K7KZ30yhPWe/etGbpWdJ1Nd9I6J7Fd6rvVh6rvXVE8usV3r6L+CH/JPv8At8l/ktfOr19FfBD/AJJ//wBvkv8AJa7sL8R8pnn8Fev6M9IorJ13Wm0eFXS2E58uWZw0mzCRrubBwct0wOM+oxTdR1uTTri6R7VXjhs2uUYS8uQcFSNvy9Rzk/Su4+XNivPdEu8aDpwz0to//QRXZWWpG40+S5nhMbxSSRSJFul5Rip24XLdPSvG9K8Rqmk2SbhxAg/8dFH9jLM/da+H9Twc9xdbDRg6S3uegfbPekN371xn/CSr6ikPiRf7wojwZFdD5t5xjH0Z2Juvemm69641vEY/vCom8R/7dddPhGK6GMswxkuh2hufemG5964hvEX+3+tQN4g/267afC8F0MnWxkjumu1HVh+dRteoOriuDbX/APbqFtd/2/1rshw5BdCPZ4uW7O9bUYh/FUTapGO9cA2uZ/iqFtb/ANquyGQ010F9SxEt2d82roOmKgfWR2Irg21r/aqFtZ/2q64ZNBfZLWWTe53T6z/tVWfWP9quIbV/9qom1b/arqhlaWyNo5V5HaPq3+1Vd9V/2q41tV/2qibVP9quiOX+R0Ryu3Q699U/2qrvqfvXJtqf+1UTal/tVvHA+R0Ry3yOpfUv9qoH1H/armG1H3qJtQ962jgjojl/kdI+oe9QPf8AvXOtf+9Rtek962jhDeOB8joHvveoHvvesI3ZPrTDcsa1WGSN44NI2XvfeoWvPeskzMe9NMhPetFRijaOGSNJrv3qFrr3qkW96TdV8sEaqjFFlrgmozKxqHdRuo54ItQSHlie9JmmZpM1LrJF2H5pM0zNGaylXHYdmkzTc0ZrCVYdhc0maTNJmsZVB2FzRmkzSVk5hYWkopKychi0lFFZtjEooorNsBDRRRUMYUlLSVDGFFFFQwEopaSpYBRRRSYwr7e8Pf8AItaV/wBecP8A6AK+Ia+3vD3/ACLWlf8AXnD/AOgCvBzr7Hz/AENafU0qKKK8M0GSxRzxPFLGskbqVdHGQwPUEdxVQ6JpTRJEdMsjGkgkVDbrhXAADAY4OABn2FXqKAIYbO1tpZZILaGKSY7pWRApc88sR1PJ6+tfKHxl/wCSta79YP8A0njr61r5K+Mv/JWtd+sH/pPHXp5T/vHyZE9jhKKKK+nRiLRQKKtCFpaSlraIhaWkpa6ICFpRSUorqgSOFKKQUorqgSOFKKQUortpksUVNHUNTR1pLYiWxZQ1ZjNVkqzHXPM5ZlqM1ajqrHVqPtXLM46hajq1HVWKrUdcsziqFuOrcfSqkdW4+lckziqFuPtVuOqkfarcfSuSZwVC5H1q3HVSPrVuOuOpscNQvQ1ZFVYqtCvFr7kQ3HinimCnivKrHoUh4qQVGKkFeXWPRpDxUgqMVIK8iuepRHinimCnivKrHp0iQU8UwU8V5NY9KkPFK3SkFK3SvLq7nqUCtLVOWrkvSqctED16ZSlqnLVyWqctdtM9CkUpepqnL1q7L3qlLXZTO+kU5KpS96uy1Sl7120z0KRTlqlLV2WqUtdlM9CkUp+hrOl+8a0Z+hrOl+8a7qR0T2Kz1XerD1XeuuJ5lYgevov4I/8AJPz/ANfkv8lr50evov4I/wDJPz/1+S/yWu7C/EfKZ7/BXr+jPQbmztbxEW6toZ1Rg6iVAwVvUZ6GoP7I0zzBJ/Z1p5gh8gN5C5EeMbM4+7jjHSrtFdx8sRwW8NrCsNvDHDEv3UjUKo78AV8b2+tulvGu48KB+lfZlfD8VtK0SEL1Ar7DhGEZVKvN2X6nFjacJpc5uDXm/vU7+3W/v1jCzmP8NO+xTelfb+wonmvD0DW/txv79Idab+/WV9il9KT7HL6U/Y0gWHoGmdZb+9TDrDf3qzjaSjtSG1kHamqVIpUKJfOrN/eNNOqOf4jVE28g7Unkyf3apUqZao0i4dTc9zTTqLnuaqGJ/wC7SeW3oar2cOxSpUyyb+Q9zTTeyHvVfY3oaTafQ0cq7FKnDsTm7k9aablz3qHaaMUW8ilCPYlM7n+KmmZ/71MxSbaNew+VD/Mb+8aTefU03FGDUuUuw7IXdRupuKMVLnMdkLuo3UmKTFQ5zCyF3UbqTFJioc5jshc0ZpMUYqHOQBmjNJiiocpDFzSZoxSYqHKQC5pM0UVDbGGaM0mKKhtgFFFFQ2xiZoooqW2AUmaWkqGMKKKKlgFJQaKhgFFFFSAlFFFSMKDRRUMBKKKKljCkNLSGpYBRRRUsAr7e8Pf8i1pX/XnD/wCgCviGvt7w9/yLWlf9ecP/AKAK8HOvsfP9Dan1NKiiivDNClq95JYaVcXUQUvGuQXGVXnBZvYdT04Fc7N4ruF06KaKSxeQG5DMASszRsAkaDdw8ikMvLcdA3WuvooAxNE1W7vdS1S0upbGVbaQGJ7MlgEZnCq5J/1gCDcBwMj14+YvjBI0nxV1xnieI7oRscgkfuI/Qkc9evevrivkr4y/8la136wf+k8denlP+8fJkT2OEooor6dGItFFLWkUIKWkFLW0UIWlpKUV0QQhaUUClFdUESKKUUgpwrrgiRRThTRTq7IIlgKnjFQip0FVIiRYSrMdV0FWYxXPM5JsspVmOq6CrUYrlmcdRliOrcdVoxVqMVyzOKoyzHVyMVVjFW465JnFULUfarcdVYxVuMdK5JnDULUdXI6qRCrkY6Vx1XocFQuQ1ZFQQipxXi13qKG48U8U0U8V5VZnoUh4p4pgp4ry6x6NIeKkFMFPFeTXPTojxUgqMVIK8msenSHinimCnivKrHp0h4pWpBQ9eXV3PUobleWqUverktUpadM9amVJKpy1blqnLXZTPQpFOXvVOXrVuSqcp6120z0KZTlqnJVyWqUldlM76RUl6mqUtXJe9Upa7KZ6FMpTdDWfJ1NaE9Z0nU130jonsV3qu9WHqu9dUTy6xA9fRfwR/wCSfn/r8l/ktfOb19F/BD/kn/8A2+S/yWu7C/EfKZ5/BXr+jOr8TazcaNbRyQeQCyyNmZSQ7KuViXBHzMeB16Hg1TvfEl0L24hsp9NC7G8gzk7ceR5i3DMG/wBTuwhwOv8AF2rqaK7j5cytH1Oe+0GO9kiSafLqVtSNshVyuULNgqcZBJ5BBr54sPDhk0+2fZ96JT+gr6cry/RrFTolgdg5to+3+yK9DAZrHL3JyfxfofO8QSrqNP2K7/oeejw0f+ef6Up8Nn/nn+len/YV/uD8qQ2K/wBwflXpLiuk+p8w545dDy8+HP8AYqNvDx/uV6ibBP7g/Kozp0R/gFbw4npvqR9Yxkd0eWt4fP8AcqJtAP8Acr1NtLhP8FRNpER7V1Q4hpv7QfX8RHdM8sbQj/cqJtDP9z9K9SbRoz2/SoX0Ne2K6YZ5Tf2ilmtRbpnl7aJ/sVEdGP8Acr059C9FFQPoR/uV0xziL+0aLOH1PNG0f/ZqJtIP92vSH0TH8H6VXfRf9muiOZp9TaOceZ522kn+7UR0rH8NehPo3+zUD6R/s1tHMfM3jmy7nANpn+zTG00/3a7t9J/2agbSv9mto4/zN45nfqcOdOI7VGbBh2Ndq+l/7NQNpmP4a1WORvHMU+pxxsmHrTDauK65tN9qhbTvatVi0zaOPTOVNu47U0wsO1dK2n/7NQtYe1WsTFmqxiZzxjPoaTb7VutY+1QtZe1WqsGarExZj7aTbWm1n7VE1p7VV4M0VaLKO2k21bNsRUZhYUezgy1NMgxSYqYoR2puKl0E9i+YixRipMUmKzeHC5HijFPxSYrKVAdxmKMU/FJisXSHcZijFOxSYrJ0x3G0U7FJWbgO43FFLiis3EBKSlorNoYlFFFZsBKKDRUsYUhpaSoYwoooqGAUlLSVLAKKKKTGFfb3h7/kWtK/684f/QBXxDX294e/5FrSv+vOH/0AV4OdfY+f6GtPqaVFFFeGaBRVLV7ySw0q4uogpeNcguMqvOCzew6npwK52bxXcLp0U0Uli8gNyGYAlZmjYBI0G7h5FIZeW46ButAHX18lfGX/AJK1rv1g/wDSeOvpjRNZuNS1C+gm8grCThYlIaI+ZIm18k5OEB7denQ14t49+GXirxd4+1nVtOsES2eSJFE86KxKwxgngnj8a9DLasKdbmm7KxM1dHitLXo//CjPHH/Ppaf+BS0v/CjPHH/Pnaf+BS19Csdhv50Y8rPN6WvR/wDhRvjf/nztP/ApaX/hRvjf/nztP/Apa0WYYX+dC5Zdjzilr0b/AIUb43/587X/AMClpf8AhR3jf/nztf8AwKWtY5jhP+fiFyS7HnNLXo3/AAo/xv8A8+dr/wCBS0f8KP8AG3/Pna/+BS1vHM8Gv+XiFyS7HnVKK9F/4Uh42/587X/wKWl/4Uj42/587X/wKWuiOa4Ff8vEL2cux52KUV6BP8F/GdtbyzyWdvsjQu224UnAGTgDrUv/AApLxt/z52v/AIErXTDOMAv+XqF7OfY88FLXof8AwpPxr/z523/gSlH/AApTxr/z523/AIEpXSs7y5L+KifZz7Hn6jmrCCu7X4LeNB1srb/wJT/GpV+DfjIf8uVv/wCBKf41Ms7y/wD5+ozlSqPocOgq1GK7Rfg/4xH/AC4wf+BKf41Mvwl8Xr1sIf8AwIT/ABrCWc4B/wDL1HNPD1X9k45BVqMV1i/Crxav/LhF/wCBCf41Mvwv8Vr10+P/AMCI/wDGueWb4J/8vEcs8JXf2WcvGKtRiulX4a+Kl66cn/gRH/8AFVMnw68UL105f+/8f/xVc8s0wb/5eI5Z4HEv7DOfjFW4x0rYHgbxGlwkB0073RnGJo8YBAPO7/aFW08D+JF66Yf+/wDF/wDFVzSzLCv7aOSeXYt7U2Y8Yq3GK008G+Il66W3/f8Ai/8AiqnXwn4gXrpb/wDf6L/4uuaePwz+2jknlWNe1NlCIVcjFWk8M66vXSpf+/0X/wAXVhNA1teukzf9/of/AIuuSrjKDWkjklk+Pf8Ay6ZHEOKmFTLo+sqP+QRP/wB/Yf8A4unjStY/6BFx/wB/Yf8A4uvKq1oPZjhk2OW9JkIp4qUaZrH/AEB7j/v7D/8AF04adq//AEB7n/v7D/8AF159XXY7aeVYxb02RininjT9W/6A9z/39h/+OU4WOrf9Ae5/7+w//HK86rSm9kd1PL8St4MRaeKjt4dSuLeOePR7vZIgddzwg4IyMgvxUwtdU/6A11/39h/+OV5tXC13tFnoUsHXW8Rwp4pgt9UH/MGuv+/sP/xynCHUx/zBrv8A7+Q//HK86rgMU9oM76dCot0SCpBUIj1Mf8wa7/7+Qf8AxynhNSH/ADBbz/v5B/8AHK86rlWNe1NnfTg1uTCkeowNRH/MFvP+/kH/AMcpGXUj/wAwW8/7+Qf/AByvPnkuYN/wmejRlFPVkUp61Slq49vqjdNGuv8Av7D/APHKrvYas3TR7n/v7D/8cq4ZNj1/y6Z6NPEUlvIz5e9U5T1rUfStZbppFx/39h/+LqpPpOsI8SHSbjMrlFxJEedpbn5+OFNdUMpxq3ps7aeMoLeSMmSqcvet5/Dmut00qX/v9F/8XVd/CviBumlSf9/ov/i66oZZi1/y7Z208wwq3mjnZT1qlJXTv4O8RN00tv8Av/F/8VUD+B/EjdNMP/f+L/4quqGX4lfYZ2U80wa3qI5OXvVOWuwfwB4mbppv/keP/wCKqu/w58Ut005f/AiP/wCKrqhgsQvsM7aecYFb1UcROeDWe/U130nwx8Vv009PxuI/8aqt8KPFx/5cIv8AwIT/ABrrp4aqt4ms85wDWlVHBvVd69Ab4SeMD/y4w/8AgSn+NRN8IPGJ/wCXGD/wJT/GuiNCp2PPq5rg3tUR569fRfwQ/wCSff8Ab5L/ACWvMW+DfjM/8uVv/wCBKf41698K9Ev/AA74Wn0vUofKuYbxywDBgQyIwII69fzBrsw8JRlqj5zN8VRrUkqcr6/ozuKKwvE2s3GjW0ckHkAssjZmUkOyrlYlwR8zHgdeh4NUNR8U3FvdXItDaTRIhaKMAszJ5HmifIbmPd8nA6/xdq6z586yvPdFv7NdC05S3Itowf8AvkV12l6lNeaS10yJcuskiKbbAWYK5UMu5scgd2/GvHdC8G+NY9LtmeyUhokZf36AgYGAQTxXh55gq2KhD2T2uehgMNha7l9Zklba56Mt1av0epAYH6OK5OLw34sTrpx/CeP/AOKq5HovihOumP8A9/4v/i6+XeUZlH4bm9XLcue00dD5KHoRSG2rJj07xKvXS5fwmi/+LqzHbeIV66VP/wB/of8A4umsHm8Notnn1cqwT2mi0bb2phtvalRNdH3tIuD/ANtYf/i6lA1f+LRbn8JYf/jlbw/teG9Jnn1cmwz2kisbb2qM2/tWhs1E9dEu/wAJYP8A45SGLUD/AMwW8/7+Qf8Axyu2liM0W9KR59XI6T2sZpt/amm39q0jBqB/5gt5/wB/IP8A45TDbaif+YLd/wDfyD/45XoUsZmC3ps4KmQxeyMwwe1RtbKeqj8q0ZYr+N4kbRrzMr7Fw8J52luf3nHCmnGy1M/8wa6/7+w//HK9Glj8Yt4M4KnD0uiMdrKM/wANQvp0Z7VuGw1M/wDMGuv+/sP/AMcpp0/VP+gPc/8Af2H/AOOV6NLM8St0zknw/XXwxZzz6Wp6Yqu+k/7NdOdN1X/oD3P/AH9h/wDjlNOl6sf+YPc/9/Yf/i676eb1VujB5Hjo/DFnISaVj+Gqsmmf7NdudJ1Y/wDMHuP+/sP/AMXUbaLqrddGn/7+w/8Axdd1POO5P9l5lHamzg30z/Zqs+m+1egN4e1Rv+YPcf8Af2H/AOLqu/hjV26aRN/39h/+Lrthm9LqzSOCzJb0mefPp3tVd9P9q9DfwlrTdNJl/wC/0X/xdQN4N1xumlSf9/ov/i66Y5vQ/nOiGFzDrSZ529h7VWex9q9DXwbrczyomlSZifY2ZYxztDcfNzwwpjeAvEDdNLP/AH/i/wDiq6I5zh/50dUMNjetNnnD2XtVd7P2r0lvh54iPTTP/I8f/wAVULfDjxIemmj/AL/x/wDxVbxzrC/8/EdEaGL6wZ5q9p7VA9p7V6Y3wz8TH/mGr/4ER/8AxVQt8LvFB6acn/gRH/8AFVtHPMJ/z8R0xpYnrBnmT2vtUD23tXp7fCnxUf8AmHx/+BEf+NRN8JfFh/5h8X/gQn+NbRz3Bf8APxG8aeI/lZ5e1vUTREV6g3wh8WnpYQ/+BCf41E3wd8XnpYwf+BKf41ss+wP/AD9RvGNbrE8xKkdRTcV6Y3wZ8Yn/AJcbf/wJT/Gq0/wb8XwoHeygwXVOLhTyzBR+pFX/AG5lz3qo1UJ9jzzFJivRf+FJ+Nf+fO2/8CUo/wCFJ+Nf+fK2/wDAlKl5xlr/AOXqL9nPsecYpMV6P/wpLxr/AM+dt/4EpSf8KS8a/wDPna/+BK1jLNsv/wCfqH7OfY85xSYr0f8A4Uj42/587X/wJWk/4Uj42/587X/wKWsJZpgP+fqHyS7HnGKSvR/+FIeNv+fO1/8AApaP+FH+Nv8Anztf/ApaxlmeC/5+IfJLseb0lekf8KP8b/8APna/+BS0n/CjvG//AD52v/gUtZPMsH/z8Q+SXY84pK9I/wCFHeN/+fO1/wDApajHwU8aNcPALO23oiuf9IXGCSBz/wABP+TWTzHCf8/EHJLsed0lekf8KN8b/wDPnaf+BS0f8KN8b/8APnaf+BS1m8wwv86HyS7Hm9JXpH/CjfG//Pnaf+BS0f8ACjPHH/Pnaf8AgUtQ8ww386Hyy7Hm9Fekf8KM8cf8+dp/4FLSf8KM8cf8+lp/4FLUfX8N/Og5Web0V6R/wozxx/z6Wn/gUtH/AAozxx/z6Wn/AIFLU/X8P/Og5Geb0V6R/wAKM8cf8+lp/wCBS0f8KM8cf8+lp/4FLS+vYf8AnQ+Vnm9fb3h7/kWtK/684f8A0AV80/8ACjPHH/Ppaf8AgUtfS/h8MvhrSg6MjiziDK3VTsHBryM1r06vJ7N3tf8AQuCa3NGiiivINAooooAwvE1zqdtbRnTfPDFZCDDD5paQL+7Rhg4UnqeMYHIzTYr3UZvEKBRcizdxiN7YrH5Jh3byxXIk8z5dpOcfw9636KACiisSfxl4WtbiW3uPEmjwzxOUkjkvolZGBwQQWyCD2ppN7AbdFYH/AAnPhH/oatE/8GEX/wAVR/wnPhH/AKGnRP8AwYRf/FU+SXYDforA/wCE58I/9DTon/gwi/8AiqP+E48I/wDQ06J/4MIv/iqfs5dgub9FYP8AwnHhL/oadE/8GEX/AMVR/wAJx4S/6GnRP/BhF/8AFUezn2YXN6isH/hOPCX/AENGif8Agwi/+Kpf+E28J/8AQ0aL/wCB8X/xVHs59mK6G3N7qMPieGLFz9jaRUKJbF4zGUb5y4UkN5mwYzgLyRjJE3h+6urtbx55LxohNtgF5beTJtAxnGxRgnJHBOOpycCP/hNfCn/Qz6L/AOB8X/xVL/wmnhX/AKGbRv8AwPi/+Ko9nPsF0blFYf8Awmfhb/oZdH/8Dov/AIql/wCEy8L/APQyaP8A+B0X/wAVS5Jdg5l3NuisX/hMfDH/AEMekf8AgdF/8VR/wmHhn/oY9I/8DY//AIqjkl2Fzx7m1RWN/wAJd4Z/6GLSf/A2P/4ql/4S3w1/0MOk/wDgbH/8VRyy7Bzx7mxRWP8A8JZ4b/6GDSv/AANj/wAaX/hK/Dn/AEH9K/8AAyP/ABpcr7C9pDua9UtXkuotKuHsw3nhflKLuYDPJVeckDJAwckdDVX/AISrw7/0HtL/APAyP/Gl/wCEp8Pf9B7S/wDwMj/xosw9pDujIn1PXv7Lsza2t7NIL7ZPKbdUkMImULlX2/ejYEsB8uG+6eV62sn/AISjw/8A9B3TP/AuP/Gj/hKPD/8A0HdM/wDAuP8AxpB7SHc1qKyv+Eo8P/8AQd0z/wAC4/8AGj/hJ/D/AP0HdM/8C4/8aB88e5q0Vlf8JPoH/Qc0z/wLj/xo/wCEn0D/AKDmmf8AgXH/AI0roOaPc1aKyv8AhJtA/wCg5pv/AIFx/wCNH/CTaB/0HNN/8C4/8aOZdx8y7mrRWX/wk2gf9BvTf/AuP/Gj/hJdB/6Dem/+Bcf+NLmj3C6NSsDXL3UbTUbY24ufs42MywWxlEg8xRIGwpIxHkjGCT69Kt/8JLoP/Qb03/wLT/Gj/hJdB/6Dem/+BSf40c8e4yppNzq8mvXsd6JPs2ZPLTyNiRBXATD4/eF1O44Py4IPWt+sv/hJdB/6Dem/+BSf40f8JJoX/Qa07/wKT/Gl7SHdAalFZf8Awkug/wDQb03/AMCk/wAaP+Em0H/oN6b/AOBcf+NHtId0OzNSisr/AISfQP8AoOaZ/wCBcf8AjR/wlHh//oO6Z/4Fx/40e0h3DlZq0Vk/8JT4e/6D2l/+Bkf+NJ/wlXh3/oPaX/4GR/40+ePcOV9jXorI/wCEr8Of9B/Sv/AyP/Gk/wCEs8N/9DBpX/gZH/jRzx7j5JdjYrjxqWuiz1FSt6bgRo0ebQgRv5rLIqEIQyhShBw5IyQHwRWx/wAJb4b/AOhh0n/wNj/xpP8AhLvDX/Qw6T/4Gx//ABVHPHuHJLsWNCmu59Gt5L4SfaTu3+Ym0nDEA4wvbH8K567V6DRrG/4S/wAM/wDQxaR/4Gx//FUn/CYeGP8AoY9I/wDA6L/4qjmXcfs59jaorE/4THwv/wBDJo//AIHRf/FUf8Jl4W/6GXR//A6L/wCKp8y7h7OfZm3RWH/wmnhX/oZtG/8AA+L/AOKo/wCE18K/9DNo3/gfF/8AFUXQezn2NyisL/hNfCn/AEM+i/8AgfF/8VR/wm3hP/oZ9F/8D4v/AIqi6FyS7G7RWF/wm/hP/oaNF/8AA+L/AOKrQ07WNL1hHfTNSs75Izh2tp1lCn0O0nFO4nFrdE90dtpKfMljwhO+FN7r7quDk+2D9K5ZtV1qbQtMmhW+S7NoTMDYsrPdBEKo6snyoxL5YAAYxuFdfRQIKKKKAMC5vdRh8TwxYufsbSKhRLYvGYyjfOXCkhvM2DGcBeSMZId4dudTuPO/tDzziOMt50Pl7Jju8xE4G5FwuG5zk/Ma3aKACio554bW3luLiVIYIkLySSMFVFAySSeAAO9Yn/Cc+Ef+hq0P/wAGEX/xVNJvYDforA/4Tnwj/wBDVon/AIMIv/iqP+E58I/9DTon/gwi/wDiqfJLsBv0Vgf8Jz4R/wChp0T/AMGEX/xVL/wnHhH/AKGnRP8AwYRf/FU/Zz7Bc3qKwf8AhOPCX/Q06J/4MIv/AIqj/hOPCX/Q06J/4MIv/iqPZz7MV0b1FYP/AAm/hL/oaNE/8GEX/wAVS/8ACbeE/wDoZ9F/8D4v/iqPZz7MLooW+qamkGom+/tTKQrInkWOXWTeweOP5MMMeXgnOQSc8HG/pL3D6TatduXuDGDIzRlCT7ggfyH0HSs//hNfCn/QzaN/4Hxf/FUv/CaeFf8AoZtG/wDA+L/4qj2c+wcy7m5RWJ/wmfhb/oZdH/8AA6L/AOKo/wCEy8L/APQyaP8A+B0X/wAVS5Jdhcy7m3RWL/wmHhj/AKGPSP8AwOi/+Ko/4S/wz/0MWkf+Bsf/AMVRyS7Bzx7m1RWN/wAJd4a/6GLSf/A2P/4ql/4S3w3/ANDDpP8A4Gx/40csuwc8e5sUVj/8JZ4c/wChg0r/AMDI/wDGl/4Svw5/0H9K/wDAyP8Axpcr7C9pDua9YXia51O2tozpvnhishBhh80tIF/dowwcKT1PGMDkZqf/AISrw7/0HtL/APAyP/Gj/hKPD3/Qd0z/AMC4/wDGizD2kO6K8d9qs3iuBBa3KaQ9q5DmNQrP+7IZs/Op5ZdpA6E8/wAO9WT/AMJR4f8A+g7pn/gXH/jS/wDCUeH/APoO6Z/4Fx/40h88e5q0Vlf8JR4f/wCg7pn/AIFx/wCNH/CT+H/+g5pn/gXH/jQHPHuatFZX/CT6B/0HNM/8C4/8aP8AhJ9A/wCg5pn/AIFx/wCNK6Dmj3NWisr/AISbQP8AoOab/wCBcf8AjR/wk2gf9BzTf/AuP/GjmXcfMjVorL/4SbQf+g3pv/gXH/jR/wAJLoP/AEG9N/8AAuP/ABpc8e4XRqVy8moa1Fqd+kcVxKfKuPIia3IhEihTDiQLzuG/OWIBwODgHU/4SXQf+g3pv/gWn+NH/CS6D/0G9N/8Ck/xo549xjPDs19Np7/b3mlkWUqk0sHktIuAc7CAVwSV5HO3Petesv8A4SXQv+g3pv8A4FJ/jR/wkmhf9BrTf/ApP8aXtId0BqUVl/8ACS6D/wBBvTf/AAKT/Gk/4SbQP+g5pv8A4Fx/40e0h3Q7M1aKyv8AhJ9A/wCg5pn/AIFx/wCNJ/wlHh//AKDumf8AgXH/AI0e0h3DlfY1qKyf+Ep8Pf8AQe0v/wADI/8AGk/4Srw7/wBB/S//AAMj/wAafPHuPlfY16KyP+Er8Of9B/Sv/AyP/Gk/4Szw3/0MGlf+Bsf+NHPHuHJLsW9XkuotKuHsw3nhflKLuYDPJVeckDJAwckdDXMXuqa8NOtDai+MqzzrIxsyGdQ37rI8s8lCDjCAnILoRitz/hLfDX/Qw6T/AOBsf/xVJ/wl3hr/AKGLSf8AwNj/APiqOaPcOSXY2aKxf+Ev8M/9DFpH/gbH/wDFUf8ACY+GP+hj0j/wOi/+Ko5l3H7OfY2qKxP+Ey8L/wDQyaP/AOB0X/xVH/CZ+Fv+hl0f/wADov8A4qnzLuHs59mbdFYf/CaeFf8AoZtG/wDA+L/4qk/4TXwp/wBDNo3/AIHxf/FUXQezn2N2isL/AITXwp/0M+i/+B8X/wAVR/wm3hP/AKGfRf8AwPi/+Kouhckuxu0Vg/8ACb+E/wDoaNF/8GEX/wAVWxa3dtfWsdzaXEVxbyDKSwuHRh6gjg0xNNbmfr9zcWlrbS2zXW/7XCrJb25m3RmRQ4YBWIGzccjHTrVaK91GbxCgUXIs3cYje2Kx+SYd28sVyJPM+XaTnH8Pet+igQUUUUAFFFFABRRRQAV8WeNf+R+8Sf8AYVuv/RrV9p18WeNf+R+8Sf8AYVuv/RrV7GS/xpen6mdTYxRThTRThX2NE52OFPFMFPFenSIY4U6minV6UNiGSx9atx9KqR9atx9KmZhULcXarkdU4u1XI65JnDULUfarkdU4+1XI65JnDULUVXIqpxVcirjmcNQtx1bi7VUjq3F2rlmcNQtx1bjqpHVuOuOZwVC3F1FXY/u1Si6irsf3a8rEGKJBTxTBTxXlVTrpDxTxTBTxXm1T0KQ4U8UwU8V5tU76RIKcKaKcK8uqelSHinimCnivNqnfSHinimCnivMqnoUhwp4pgp4rzqp6FIjaq0nerLVWk71jE9SlsVZKqSVbkqpJXTA76ZVk6GqslWpOhqrJXXA7aZVeqknSrb1Uk6V1QOymVZOtVZKtSdaqyV1QO2mVZKqyd6tSVVk711QO6mVJKqS1bkqpLXXA647FB+pqBqnfqagauyJjVIWqBqnaoGreJ5lUhavcfgB/yDNb/wCu0X/oLV4c1e4/AD/kGa3/ANdov/QWrrw/xnzmc/7u/VHsdFFFd58oFFFFABRRRQBgeOv+SfeJf+wVdf8Aopq+LxX2h46/5J94l/7BV1/6Kavi8V9Bke0/l+plVHCnimCnivrKJgxwpwpopwr06JDHilHWkFKOtehHYhliPpVtO1VI+lW07VnM5qhciq1H0qrFVqPpXHM4ahbj61bj6VUj61bj6VyTOGoW4u1XI+lU4u1XI+lckzhqFuPtVuOqkfarcdckzhqFuOrcdVI6tx1x1DhqFuKra9BVSKra9BXkYjcyjuSCnimCnivLqnbSHinCminCvNqnoUh4p4pgp4rzax6FIeKeKYKeK8uqejSHinimCnivNqnfSHinimCnivNqnoUxwpx6U0U49K86qehSIXqrJ0q09VZOlRA9WnsVZKqyVakqrJXVA7aZVk6VUk71bk6VUk711QO2mVpOlVJKtydKqSV1QO2mVZOhqpJVuToaqSV1wO2mVZKqyVakqrJXVA7aZVk71Smq7J3qlNXXTOtbFJqgbvU7VA3euuJzViFqhapmqFq2ieZVIHr6p+GP/JN9E/64n/0Nq+Vnr6p+GP8AyTfRP+uJ/wDQ2ruw27Plc7+CPqdbRRRXYfOhRRRQAUUUUAFFFZut38mnWUcsckEO+ZY2nuBmOIH+JhkfTqOSKANKvmjX/hB4u1rxRrWp2VrbtbXOo3MkZecKSplbsa9rXxVcyjTnjS2zcQQO1vyXld3KyJGc8GLGWyD7461Tj8Y37WQl3WDHYzh1RtsriONhbr83+sJdgDz937vUDow+JqYeXNT3E0nueIj4G+N/+fO0/wDApaUfA7xv/wA+dr/4FLX0Roms3GpahfQTeQVhJwsSkNEfMkTa+ScnCA9uvToa3a7453i47NfcR7OJ8uD4IeNv+fO1/wDApaUfBHxt/wA+dr/4FLX1FVLV7ySw0q4uogpeNcguMqvOCzew6npwK2jxFjo7NfcL2MT5sHwS8bf8+dr/AOBK0yH4NeMp0LpZ2+A7JzcKOVYqf1Br3SbxXcLp0U0Uli8gNyGYAlZmjYBI0G7h5FIZeW46ButJfeJdQim1WO2vNGAtHjKSXJ2IQfM3ISshy42DGdp6/LjBrVcUZguq+4XsIHiafBbxovWytv8AwJT/ABqdfg74yXrY2/8A4Ep/jXuul61d3mvzWczWbW7QCaFIc+dDgJlZwThWJc4A6gH056Ch8UZg+q+4l4am9z5xT4R+MFxmxh/8CU/xqdPhT4uXrYRf+BCf419D1S1e8ksNKuLqIKXjXILjKrzgs3sOp6cCs3xJjn1X3GbwFF7nhifC/wAWL10+P/wIj/xp9t8PPE8qF004YDsnMyDlWKnqfUGvTpvFdwunRTRSWLyA3IZgCVmaNgEjQbuHkUhl5bjoG61Pda5qNtd3itc6fHaJcpbrdSwMqwnY0jb/AN5hh/q0BBX5mPpiofEGNfVfcZPK8O90/vPNU+HfiheunD/v/H/8VU6eAvEy9dN/8jx//FV6Lo/iDUL7XEs7m2hjjktROUQgvCSsZw/zZBJdhgqvQYJ5rpqzeeYt9V9xm8mwr3T+88bTwR4kXrph/wC/8X/xVTp4O8RLjOlt/wB/4v8A4qvXapaveSWGlXF1EFLxrkFxlV5wWb2HU9OBUPN8S+q+4zeQ4N7p/eear4U8QL10uT/v9F/8XUkGga27youlTZifY2ZYhzgNx8/PBFdWPE5GnWNy97pwL6gbSU9A43lQyfP8pK7XGd3DD1zVKfxbqKNcJENPDLP5fmTFkjtfnlGJju6kRrjpy49s5vM8Q+q+4yfDmBe6f3mWmg62p50mb/v9D/8AF1YXR9ZUc6RP/wB/Yf8A4ut7S/EVxqOuxWbWksFvJZ+evnW7oxYbM8n5cfPjHJyp5ro6xli6styf9WcB2f3nAjStY/6BFx/39h/+Lpw0zWP+gRcf9/Yf/i67yqWr3klhpVxdRBS8a5BcZVecFm9h1PTgVhKo5bmi4dwS2T+85Eabq4/5g9z/AN/Yf/i6SK11SR5UXR7rMT7GzJCOdobj5+eGFag8TkadY3L3unAvqBtJT0DjeVDJ8/ykrtcZ3cMPXNV7nxbcwnVWcQ2cdo8Bj+2QOm5GdlOSxUEnbkEdAeQaxlBS3NFkeEWyf3kAsNWH/MHuf+/sP/xynCy1Yf8AMHuf+/sP/wAcrWs9au5/FdxpjLbvbLE0qPCwYoB5e3cQxOW3scFV6cbuTW/WUsLTluarKcMtk/vOMFpqo/5g91/39h/+OU4WuqD/AJg11/39h/8AjldjVLV7ySw0q4uogpeNcguMqvOCzew6npwKxll1CW9zZZfRWxzgt9U/6A11/wB/Yf8A45TYv7RkeVF0a8zE+xsvCOdobj95zwwq0PE5GnWNy97pwL6gbSU9A43lQyfP8pK7XGd3DD1zTrvVtWhi1qWKaydLR0its2z5klIz5Z/ec53xKCMcluKxlk+Flun95qsLTWxCItT/AOgNd/8AfyD/AOOU4JqX/QGu/wDv5B/8cq9bardHxP8A2VLNayqLUyv5KYZJAUGD87HB3E8qvbBbmtysZZBg5bp/eaqnFbHLBdS/6At5/wB/IP8A45Tsaj/0Bbz/AL+Qf/HK6eorqVoLSaZIzK8aMyxr1YgZwPrWUuGsBLdP7zWMnHY5opqR/wCYNef9/IP/AI5VYLqM7zImjXmYn2Nl4RztDcfvOeGFP/4S1hoV3cvqWlrcQypHG4QlJt8YZQFMgx8xZclsfIxOMHFm48Sra3eqw3Op6bALW0t7kNtLmPeWDBgHBfouMbfvr1yMyuGMvXR/ebRxNSOxQay1Y9NHuf8Av7D/APHKhbTNYbppFx/39h/+Lq9YeJ76417TtPmFntuIEeTyyCdxjdztIc8ZUDhSvX584B66rXDeBXR/eaLH1lseeto2stnGkT/9/Yf/AIuoH0DW26aTN/3+h/8Ai69JqK6laC0mmSMyvGjMsa9WIGcD61a4fwS6P7zRZpiFs19x5q3hvXT/AMwqX/v9F/8AF1VPhvXXuHgGlS70RXOZYsYYkDnd/sn/ACa6n/hLWGhXdy+paWtxDKkcbhCUm3xhlAUyDHzFlyWx8jE4wcJqPi6+tbq6iW3jjSKCNkeSF9rMZERmDkqhX5zjkD5clgDxayPCLo/vNFnGKWzX3HKN4P8AELdNLf8A7/Rf/FVC3grxG3TTD/3/AIv/AIqu60jxLNqOoafbG6sWkntvOnt0Xa8fBwQd5DZI6AHjJ3EYLdRVrJ8Kuj+81We4xbNfceLP4E8St003/wAjx/8AxVQv8PvE7dNNH/f+P/4qvb6iupWgtJpkjMrxozLGvViBnA+taLKsOuj+80XEWOXVfceHt8OfFLdNOX/wIj/+Kqm3w48USXDwLp670RXOZ0xgkgc5/wBk16j/AMJaw0K7uX1LS1uIZUjjcISk2+MMoCmQY+YsuS2PkYnGDhV8U3bazqNm8tjDFb25kRzg4YCLqTIuQS7AFgg+7gnnFrLqC7lribHrqvuPKG+FXi4/8uEX/gQn+NRn4TeLz/y4Q/8AgQn+NezaZrGpahcaR89p5dzbvPcoIWyoX5cq28gZZlxwwIDYY8GukrRYOkiZcSY6W7X3HzifhF4wP/LjD/4Ep/jUbfB/xif+XGD/AMCU/wAa+kqiupWgtJpkjMrxozLGvViBnA+tUsLTRjLPcZLdr7j5uPwb8Zn/AJcrf/wJT/GvTPhH4X1XwpHq9nq0CxSymKVNjhgV+cdR34/lW7H4omk0G+uRfaa1zbOgR1AEcoZFYYDSrjkuuS/VG78B58S3LS38Udxp/nRWNrdpE/yNGHLCTducZwACM7PvAEjOauFKMHdHJiMwrV48k9jqqK5XS/E09/qOnWrTWQkuEd5bbbiQIC4WRTvIIbaCFAPGTuIGT1VanEFFFFABRRWBrmvSaVqNtEXtooG2FzPkGQGRVbacgDYpLnIPHp1oAm8X20l74K161hGZZtOuI0HqxjYD+dfOI+Bnjj/nztP/AAKWverjxReRXVxEGs0VZWjYuh/0QCZYw0vzDIdWLr93gdxzVf8A4S/UPs/m+Va+Z5G/7PtbfjyPN8/73+q3fJjHX+LPFdWGxlXD39n1JcU9zw4fA3xv/wA+dp/4FLTh8D/G/wDz52v/AIFLX0not9JqOnefI0UhEskYlhGEkCuVDKMnggeprRrujnmLjs19xPsony4Pgh42/wCfO1/8Clpw+CPjb/nztf8AwJWvqGsDXNek0rUbaIvbRQNsLmfIMgMiq205AGxSXOQePTrW0eI8dHZr7hexifPU/wAGfGdtbyTyWdvsjUu224UnAGTgDrUo+CnjX/nztv8AwJSvcbnxVLbS3vn3djbRwmQfPEztCVlVEDjeMmUHcv3fx60l14ovrZ4lDWMrtp/2iRIvnET7QckhySnJOdoBGPmzwdVxTmC6r7hewgeKJ8GPGijmyt//AAJT/Gp1+D3jIf8ALjB/4Ep/jX0Hot9JqOmLcStE7ebKgkhGEkCuyhl5PBAB6nr1rQqXxPmD6r7iHhab3PnNPhL4vXrYQ/8AgSn+NTp8KvFq9bCL/wACE/xr6FrC8TazcaNbRyQeQCyyNmZSQ7KuViXBHzMeB16Hg1D4jxz6r7jN4Ci9zxqT4b+KbW3lnk09dkaF22zoTgDJwAeatr8N/FK9dOX/AMCI/wD4qvULnxPLFqkkcCRXUA/1cECl5pU8nzBMpBOULfu/u9ec9qoWvjDUpZrRXtrchgwkURukk7ebJHiJGOfl8sFs5wGzWbz/ABj6r7jN5Xh3un95wifD3xOvXTR/3/j/APiqsJ4D8Sr103/yPH/8VXqHhvVptXsppZnt5CkiqJLdSEbMaORyTyCxU89u3Stmoed4t9V9xk8lwr3T+88cTwT4jXrph/7/AMX/AMVVhPCHiFeulv8A9/ov/iq9brC8TazcaNbRyQeQCyyNmZSQ7KuViXBHzMeB16Hg1m83xL6r7jN5Bg3un95w6eFtfXrpUn/f6L/4un2uha1PbxTppU+yRA65kiBwRkZBfium1PxVc2txdG0W1uY442aKFcl5E8jzRNuBx5e75OB179qjt/E+pPNpqstjMLiXyWFuc+YRKUZ0y+dgUBshWzhs7cVDzOu97fcZvhvAvdP7zITQtbXrpM3/AH+h/wDi6nGkayB/yCLj/v7D/wDF11Xh3VTq+ltO1xazyR3E8DvbcIdkjKONzYyoU9e9a1YSxVSW5K4ZwHZ/ecCNK1gf8wi4/wC/sP8A8XThpmsf9Ai4/wC/sP8A8XXeVheJtZuNGto5IPIBZZGzMpIdlXKxLgj5mPA69DwawlNy3NFw9glsn95hDTtX/wCgPc/9/Yf/AIum29rqlxbxTx6PdbJEDrukhBwRkZBfitO98SXS308Vk9jIgQmFXPLJ5HmictuA8rd8nQDP8QqPT/FFxeXOmQNPYpNdSShoCuGeNWZRKrCQgA4GAN2TnBxkjKVOMtzVZJhFsn95VFjq3/QHuf8Av7D/APHKcLLVR/zB7r/v7D/8crf8O6r/AGxpRuGuLWd0uJoWe2+4dkjKDjc2MqFPU9a1qxlhKctzVZVh1sn95xgtNV/6A91/39h/+OU4W2qf9Aa6/wC/sP8A8crsawvE2s3GjW0ckHkAssjZmUkOyrlYlwR8zHgdeh4NYyy2hLe/3mqwFFbGaLfVP+gNd/8Af2H/AOOU23/tG4t4p49GvNkiB13PCDgjIyDJxWo+uIdbazj1HTxFJYG6iY8lcEfMfn+ZMHPGOB1rPsvE93djw/OlxpckGoW8bTRRMTN5jA5KLu+4COc5IwfSsZZNhZbp/eaxw1NbCiPUx/zBrv8A7+Qf/HKcE1If8wa8/wC/kH/xytTw7qv9saUbhri1ndLiaFntvuHZIyg43NjKhT1PWtasZZBgpbp/eaqCWxy4XUv+gLef9/IP/jlLjUcf8gW8/wC/kH/xyunrC8TazcaNbRyQeQCyyNmZSQ7KuViXBHzMeB16Hg1k+GcA90/vNYzcdii0epH/AJg13/38g/8AjlVY01K6t4549GvNkiB13PCDgjIyDJxWvNryLrL2iajp/kvYG6iYkEryMMTvAZSDnt061j2njC+mu9LhdbTFy5WQhcbh5zxgr+8PZQfk8wZPJVcMUuGMvXR/ebrFVFsI9hqzdNHuf+/sP/xyoG0rWG6aRcf9/Yf/AIuuo8O6qdX0tp2uLWeSO4ngd7bhDskZRxubGVCnr3rWq1w5gV0f3mix9ZbHnbaLrTdNIn/7+w//ABdQv4f1xs40mb/v9F/8XXpVZut38mnWUcsckEO+ZY2nuBmOIH+JhkfTqOSKtZBgl0f3mizTELZr7jgG8M683TSpf+/0X/xdVY/DWu3KF49LlwHZPmljHKkqf4vUGu2h8SmVtCYzWi/2jCxe3z+88wJkqp3dmBUgjOfSqNh4vuptNtLq8m02Bpb2KB1GdpV1VtqsWHzruweD06d60WR4RdH95os5xS2a+45RvB3iJhxpbf8Af+L/AOKqFvBHiRummH/v/F/8VXp2kaqb++1i1e4tZXsrvylEHBVDGjDcNx5BZhnjO3oK1qtZPhV0f3miz7GLZr7jxVvAfiVumm/+R4//AIqoH+Hvidummj/v/H/8VXuNZut38mnWUcsckEO+ZY2nuBmOIH+JhkfTqOSKtZXh10f3mi4ixy6r7jxpvhv4pbppy/8AgRH/APFVTX4beKbqMvHp64DsnzToOVYqe/qDXr8HiXzW0Jnns0GowsZIDzIJAmSFO7swKkbSc+9Y9n431CbSYLh7aJ5HvYYXkhgeSNI28rcDsZ/n/eELkjpkgEbTosvorYv/AFmx/dfceaH4U+Lj/wAuEX/gQn+NRn4S+Lz/AMuEP/gSn+Ne5aJq7apd6ii3VpdQW0gjWS3G0hssGUjc3TAG7jJ3ccZOzVrB0kTLiPHS3a+4+cD8IvGB/wCXGD/wJT/Goz8HvGR/5cYP/AlP8a+k6zdbv5NOso5Y5IId8yxtPcDMcQP8TDI+nUckVSw1NGMs8xct2vuPnpvg14zP/Llb/wDgSn+Ne5+AbG40vwRpthdx+XcWweKRc5wyyMDg9xUK+KXlOj+TNZNPeW/my2X/AC0GFO7DbuCHG3btJJB9CRHa+IdUvbGCS3n095pLyK3UrA7JIGRJHx+84KKZM9c7Ox4rSFOMNjixONq4hJVOh11FY2iau2qXeoot1aXUFtII1ktxtIbLBlI3N0wBu4yd3HGTs1ocgUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAf/2Q==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig = pyplot.figure(figsize=(20,10))\n", "displacement.plot(gridLines=None, figure=(fig, 121), colorbar=\"horizontal\")\n", "s = sigma(displacement) - (1./3)*tr(sigma(displacement))*Identity(2)\n", "von_Mises = sqrt(3./2*inner(s, s))\n", "plot(von_Mises, gridView=gridView, gridLines=None, figure=(fig, 122), colorbar=\"horizontal\")" ] }, { "cell_type": "markdown", "id": "2eaa81f7", "metadata": { "lines_to_next_cell": 2 }, "source": [ "Finally we can plot the actual displaced beam using a grid view that\n", "allows us to add a transformation of the geometry of each entity in the\n", "grid by prociding a grid function to the constructor. Note that this also\n", "allows for higher order transformation like in this case where the\n", "transformation is given by a second order Lagrange discrete function.\n", "We will highlight the flexibility of the # `GeometryGridView` in further\n", "examples:" ] }, { "cell_type": "code", "execution_count": 7, "id": "6f49d091", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:09:35.451032Z", "iopub.status.busy": "2024-02-29T12:09:35.450658Z", "iopub.status.idle": "2024-02-29T12:11:46.691143Z", "shell.execute_reply": "2024-02-29T12:11:46.690200Z" } }, "outputs": [ { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQEAMgAyAAD/2wBDAAgGBgcGBQgHBwcJCQgKDBQNDAsLDBkSEw8UHRofHh0aHBwgJC4nICIsIxwcKDcpLDAxNDQ0Hyc5PTgyPC4zNDL/2wBDAQkJCQwLDBgNDRgyIRwhMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjL/wAARCACFARcDASIAAhEBAxEB/8QAHwAAAQUBAQEBAQEAAAAAAAAAAAECAwQFBgcICQoL/8QAtRAAAgEDAwIEAwUFBAQAAAF9AQIDAAQRBRIhMUEGE1FhByJxFDKBkaEII0KxwRVS0fAkM2JyggkKFhcYGRolJicoKSo0NTY3ODk6Q0RFRkdISUpTVFVWV1hZWmNkZWZnaGlqc3R1dnd4eXqDhIWGh4iJipKTlJWWl5iZmqKjpKWmp6ipqrKztLW2t7i5usLDxMXGx8jJytLT1NXW19jZ2uHi4+Tl5ufo6erx8vP09fb3+Pn6/8QAHwEAAwEBAQEBAQEBAQAAAAAAAAECAwQFBgcICQoL/8QAtREAAgECBAQDBAcFBAQAAQJ3AAECAxEEBSExBhJBUQdhcRMiMoEIFEKRobHBCSMzUvAVYnLRChYkNOEl8RcYGRomJygpKjU2Nzg5OkNERUZHSElKU1RVVldYWVpjZGVmZ2hpanN0dXZ3eHl6goOEhYaHiImKkpOUlZaXmJmaoqOkpaanqKmqsrO0tba3uLm6wsPExcbHyMnK0tPU1dbX2Nna4uPk5ebn6Onq8vP09fb3+Pn6/9oADAMBAAIRAxEAPwD3+iiigArPuNUaG/ezh0+6uXSJJWaIxgAMWAHzuvPyGtCs2D/kZb7/AK87f/0OagA/tO7/AOgHqH/fcH/x2j+07v8A6Aeof99wf/HaTX5pbfR5JISwfzIlG1XJwZFB4T5uhPSvNB4o16PTFnSW4LiFfn2XH3vsryE/NEw+9tPpx6ZNAHpn9p3f/QD1D/vuD/47R/ad3/0A9Q/77g/+O1wU3i/XY7sZleOB52RTKjKABLDGfma2Hdn6+vPQArB8RrlJIxLLbPE4D7jJA5CHzmHSVOdsadQOSfUAgHef2nd/9APUP++4P/jtH9p3f/QD1D/vuD/47XLRePbuKN1vNLczpGWZEilHzBIm2gqrgktKFHPoTjIrRtvHul3E/wBnI/f+YI9kc0ZJYyNGAAWDZypONvA5OKANj+07v/oB6h/33B/8do/tO7/6Aeof99wf/HaZB4l0ieJZReCONgMPOjRLym8DLADO3nHpzWlFNFPGJIZEkQ9GRgR+YoAof2nd/wDQD1D/AL7g/wDjtH9p3f8A0A9Q/wC+4P8A47WlRQBm/wBp3f8A0A9Q/wC+4P8A47R/ad3/ANAPUP8AvuD/AOO1pUUAZv8Aad3/ANAPUP8AvuD/AOO0f2nd/wDQD1D/AL7g/wDjtaVFAGb/AGnd/wDQD1D/AL7g/wDjtH9p3f8A0A9Q/wC+4P8A47WlRQBm/wBp3f8A0A9Q/wC+4P8A47R/ad3/ANAPUP8AvuD/AOO1pUUAZv8Aad3/ANAPUP8AvuD/AOO0f2nd/wDQD1D/AL7g/wDjtaVFAGb/AGnd/wDQD1D/AL7g/wDjtH9p3f8A0A9Q/wC+4P8A47WlRQBm/wBp3f8A0A9Q/wC+4P8A47R/ad3/ANAPUP8AvuD/AOO1pUUAZv8Aad3/ANAPUP8AvuD/AOO0f2nd/wDQD1D/AL7g/wDjtaVFAGb/AGnd/wDQD1D/AL7g/wDjtH9p3f8A0A9Q/wC+4P8A47WlRQBm/wBp3f8A0A9Q/wC+4P8A47R/ad3/ANAPUP8AvuD/AOO1pUUAUbLUTd3U9tJZXFrLCiSETFDuViwGNjN/cPWio4P+Rlvv+vO3/wDQ5qKANKsDXNek0rUbaIvbRQNsLmfIMgMiq205AGxSXOQePTrW/RQBk2Gqm513U9Oa4tZDa+WypFw6hgchhuOSMdcDr0otJGfxPqQaF4wtrbgFiuHG+bkYJ4+uDWtWbB/yMt9/152//oc1AFfxWgfw7OpVWBlh4ZdwP71OxR//AEE14usVr/ZJD/ZFIgZ+fIU8WaAffgXu/r19/mPtHitwnh2diwUCWHkttA/ep33pj/vofWvG7S4aSGKFLokyRJFtW6/vLaL0+1D1x/8Aq+QAsxz2q3pSKS2jIuBtZHgjI/0qRzzHOg6RDt6YxwSkf2p9MEf+ls7W+Ewbg5/0aONVG0yoebg8YHXHGcU6d7tYpbhWvMukkikfaRjcl43UNKp4YHn1P+0aSSKAagXuI4I2FwMiaOFCwE6qOJYYwQFtmz8w/D5hQAH7Kt8TMLZJhcFlMiwIf+PhmYtuWBhhbcYIPfHGcUm+7FgG/wBI88QBjvEzqcQZyciZGLSXQ9Mf7O6o7YzXOnqsRlb9wI0ELSkFjCiBf3csgHz3RPKgduMlafJPbpfGecwrMZvOEmIlfb5zyFuRBJnbboBg9G7bgKAJIYYpLiS2sQpQu0Mi2xBO1pFh48h1Zh5cMhz5WRuOQc5qSLVZIZUumnKu5EwUuhlQNvuWGH8mQkjyl3BiCCcE5AqvcwzLbC3vN6yKnl7rgsFMgjWPAM6MoIluXJIkA4J44NThi7My/ubWYnld6wmFpDn7vnRBfKtcYwOvcCgDZi8Xa1pUawzzlnRSimclBNIiIhAEyqSTM+GxJwEbAHWuig8fLG4W+thGPM2l2zF1lMa4LZRs7WbIfGFPXFef200kMIubZG3IA7RwDCtIq+ac+QXTiWaIANGAdqjjBpbeAoJItPJPlHyVmtSch8C3Qlrc5xkzuN0WSB33UAeuWnijSruCKU3HkLIqspnG1fmTeo3fdyV+bGc4rXVldQyMGUjIIOQRXhsYhm3yxbVikOZXtsAMkhJIZrfBXFvCfvRE5dsgZq5Z63f2u65iklfbme5EWCwbmeQFoe3zxJmSM4LMPloA9norzuy8d3ttCftax3ix7lM0eGWR1CK2HjztzI+0bkUY7109j4t0q9Qt53lqG2+YSGjPzsgO9cqAWU4yRn0oA3aKZFLHNEskTrJGwyrKcgj2NPoAKKKKACiiigAooooAKKKKACiiigAooooAKKKKAM2D/kZb7/rzt/8A0OaiiD/kZb7/AK87f/0OaigDSoorA1zRru+1G2u7VbaRodhTz3KGMrIrnaQp++BsJ4wPXpQBv1mwf8jLff8AXnb/APoc1VNM0a4s9buLyTyAj+dmRCTJPvkDLvGBjywCi8ng9ulT2gmHifUvNdGU2tvsCoVKjfNweTk+/H0oAPEu7+w5Nu/d5sONnmZ/1qdPL+b8q8gtGuitooW6I/cdfPYf8uXZoGHp69uvAb2vULCHUrJ7S4GYnKsRtVvusGHDAg8gdRXMp8N9Ei8soCTHtwHghIO3y8Z+QZ/1S5+retAHl0dvA1vEk8cMfmRImZI4oz80dumR5kCg8Ssfve+e5lWZriKZ7NioaN5NlscY3JOw+WKf/pui8L3xjghvRY/h1DbmL7PqDxLGUwqeYgO3yeyyAZPkj8WJ7DFGT4dXg8mMXcc0MRj4lffnb5AOBIj4yITxnHzAdAcgHH3ckS3jT3yDZHKzDz1UkqskjYHnxKW+W2RQC+e3cgEMNx9mSAGZ4mUW4MTSMg4ityRsaZCOZwBtHcjsK2j4I1m1hWIW+2LasbtCCmcpHGxPlSDPD3BzsJ55B3YGTd6ffpKZ5rRo76XLIJAhcOwcjl1ifPmXKAEMeUY87RgAgtpl+a5slVtqfaJEtyM5G+cKfs7o2N0kAyU4YdtuKctskl19kQq0zE25GFaZF+W3zgCGUtsimbGG4IyMkmi4eKNle6VhbRP5i+cCVKhi42rOrKcx20YAEvII5OSKa0Mn2IW05YRlREZCzCMHAhBy/mREb5LhgQy5KnleBQA57nzHS6lCvOP36RyMDKn3rk48zy5ephUMHY9cE8CnSxP9ojtzIZLpB9ntmk/1m8YhDBZtki5kknf5XbDL0+Wk89RC1wYyto/787Rsj2Z88rj97DnYlsuPlwSox3pAksdmFQhFK+UJVOyHOPJDMf3kODLLOcnb9wkEZoAc/wBnnlEbhdjA+Wswz+6b+JVm2yIBbwDG2Rj84AIzTnkkeVTcoxuoiZTDJ800ZH+kSYEmyYZZoUGHY8PjPFNJiEOT/o+mz9eBHH5TDP8Atw7zBCP7jbpKGkeCDN0qrHgym3YhIyRidljD74SDI8C4RlOExjmgBxEnn4DefewcKeZH3IdoJH7u4UtcTM2eenQgUsYaeUGwZ2nTbHFKhLyxrgwxDfGEmxtWaQgoeAuQc0jRNBtss7pFYRxRzKoDyKTDGPLmLIR5skzZjkPQnAwKTalzIlqMhWX9wJFLMiMpRZFSUiQBYI5JBtkYZZcelAF7TtUubcC4srny4XI/eKwMYDY25kiGF8uCLdiSPgsxI5ro9N+IUohWW/iX7ORuMr4AQbWkbLplDtUxLyFJZwMCuPG2eZ5JwdxTdKCHeSGJgHdW+7cIBCkceRvXDnnmhpmWSS5kRpJVBa6aN9z5GJZEMsY3Y3NDHiVGAIIzxQB67Z+ItNuyq+cIJG6JMQMn5cgMCVbG5fuk4zitWvD8vbj7SJFJBLNMjBI53R8s3mIGiO+4bA3oo2xLyK1bDxBqWl7ILWcucCOFHxH5pz5UY2sTG2+UsxaOT7qHpQB63RXJaZ44t7q4S2uIHWWTBjCKVd0LMFYRt8xBClspvAAySK6Oy1Gz1GES2dzHMuA3ynkAjIyOoyOeaALVFFFABRRRQAUUUUAFFFFABRRRQBmwf8jLff8AXnb/APoc1FEH/Iy33/Xnb/8Aoc1FAGlRRVO51SztLuK1mlKyy42gIxAydq7iBhcngZIyeBQBcrNg/wCRlvv+vO3/APQ5qtLfW76hJYhn+0RxiQqY2A2k4BDEYP4HiqsH/Iy33/Xnb/8Aoc1AGlRRRQAUUUUAFIyhlKsAQeCD3paKAMqfw5pFwpBski3Ag+QTFnIAOduM8ADntx0rBvfh9bu8s9jdvb3L7mMu3Y5ciXndGVP3pSec9B24rs6KAPLL7wbqdlcvdLb+cgO5jAPnZQ24IGj2PysUKcq3LM3asJrK+hncEP8AbgQrXMSkvvP7tSWj2SjMkszEsjE7F+8K9wqC5srW9VVubeOULyu9QSp9Qex9xQB4jGYZC1xHjySMyvBzmNhuKM0IDgCCFRh4jy7ZxmnI8sLNdoy+eG3PJGRskmUiRsyQ5XHnzRjEkeMIAWGCK9L1HwLpd6VZMxlcABxvAGY8gHIcZWJV4YYGcda5fUPBerWpE8ZNzNGMpMMyNuG4g7gVlBMsjucFgAq8HGKAObVI4Itls/lwuPJjmVlRSOYFYsu6F8/6S/zKh4zupJWi+zgOogsLj5ihRY49jqGACNuh3i3jU/IyNmXFTzWjRSSK0TxgfujIc7kQgxqPNjAkG2CKVz5kZHPIOaiWdTFLekFEcM8zx4TIOJZQzxgxuAvkRYdFyZMEDFACzqxkCXkWZEYu9q68hhiWRUjlIyAfJhHlyZwhwDmnkyedjm6urbt88j70bk87LiNnuZB0z0PUColLw24SRYwD92P5UimZXDEKDuhbfcOudpX5YTwM09YHR0s4yXlGFt45EyxYboo2EUhJI3NPJuhdhkA4GKACB3LtLbShpI2VFuQc/OpMcRkkjAON/mzHzYySI1yTTY0ihiDW7+RBMAsUm5EBVgYkO8ZhciMTPhxG2W655pSwlkT7OxyqgROrO7wRspRCDxcR4hSSTB3qA445oylxG1ySIopMiSVWwEBQFlaSIceXbrGnzxn5pXz3oAJVSWAQbVggnO5YWQIg3JxiKQ7MpbqWyjg7plqzbXVzDMxaWWGWMF5VKu726kK7ko2J49sSQxhl3qMnnmqxeOCFjcZiilDNMq7EDqcSygDmCXI8mMbSjcY20TxyCRLeeFXkjZm+yhMgurBpRHE5BAMpSIeU/IjbANAHVWHju+tI2a8jE0UYYzYPmCIqhd13qMgJvhT51zuLZbiuwsvE+mXioTMIS5YIZGGx8EKdrglW+YgcHrXlILhy5zdTWrYZsu7bkc4BYYnhaS4dj8wxhRyQKjEksDtcW06sdy/6VvUJM6sUjZ5EHluGmMjgyqCQg55oA9yoryfTNd1DTTBbWbuQ4VbaGTCNIpxHDhWPlyBtrSs0bAnPUV1eleO7O9lEE8bJIwDKUU7ih3lWMZ+YZRA/G4AMOaAOtoqC1vLa+hE1rPHNGcfMjZxxnn0PNT0AFFFFABRRRQBmwf8AIy33/Xnb/wDoc1FEH/Iy33/Xnb/+hzUUAXpoIrmF4Z4klicbXR1DKw9CD1rKHhfS4rhJbSBbJNyPJDaxpGkpRt6FgFzkMM8EZ6HIrZooAz4tLMeqLftfXMsgtxAUcR7WAOdxwoOSeeCB7VDaQQw+J9SaKJEaS1t2cqoBY75uT6mtas2D/kZb7/rzt/8A0OagDSooooAKKKKACiiigAooooAKKKKACiiigCre6bZalGUvLaOYbSoLDkAjBweoyCRxXOX/AIHt5bg3dpKRPu3nzGKsSHMmPNX5sF9pO7d9xQMYrraKAPJL7w5rGkOzIhLMVBuMBN7jCq8hUGNx5kskn7xAcKPmzWNJBH5KQL/o1rKFWNJNsa7WXanDkws4gVpDhkO6bjmvdKyb/wAN6ZqCyh4fKaVWWRofl3hgA2R0OQMEkE4oA8lcLMfJuo9oALyQNGWMKlQ8h8mQiRAkKRRgxtj0PNNaZmmkluInlljB+0BXZ5ByHkj8xcSrl2ihxIrgYIzxXXal4EuoAXsxHdQBjJ9lKrsY72kK+W/ygM4iB2lflU1zdzp93aTBXikmngIMbMruxYNhWKkieMyTuXyjEYH3SBQBWYPEXujKkjbiz3Csqxzur5JMiAxnfcNgCVFBWMcinKvkhLSFQ5bYsMDoqrOVYpEDE5MThpmeQtG/QE4BpIJd0gls5RvUqkV0XztIykRaVFBKjEs582MHCDJPWokiiEKmJjbWlwAsbEpGCroVT5uYHcQB3Iby2zJ60AOWJ5biGG1kd3jVfspKu8ka4KRSBGxMox5s3yFwMrxUsbLLtubaRYgcJFcLKcQ7k2Rj7Qg3KI4FZysq9ScjmkkCTK1tMohRsu1u0RxCGTczeRIdy+XbrgGJjyygGo2Zpp5prpHZo1Y3Q3O8kSEB5UZ1AnQBBHDlgwG5uaAGFo/K+bNrZ3JJIwkS7XTJ45gZ1twBxsbMtWJUEhe3u4hGfmee3MRPlZUPKfIc7gEiVIlaJh7daHk8pZLq4cHzA5nnjdVE+CHnxKgMUoZzHEFkX0HFRYlhdLWWFJGDHbalAEmdXDMqROShD3BVT5TDKxnigCzBqt7bXT3DPKZISXuArvI6YKySKXGJV58qH94Hx8wzxXVWPjq4tLdTqQSaNVbM6ldkhQAOySL8hzIdiqQpPc1yO078R7ruW2OUJDu+5XwpK5E8LSXDFsqxGAOCBUcUzpIZ7WZSwZVS6Mg2sysUiaSVBhh5nmTHzUyQgyTQB7DZa5p98/lRzhJslfKk+ViQxU47NyrDIJHFaNeIRIIRHHZ7oRKEFsrBELhgY4PlbMEp275TjYctwRW3pXi++tJ0gRy0J2kRurNsjOCGMbnzAFhRnO1nBLD1xQB6pRXN6V4ysdQjXzlMD8B8HeI2KbyrcZUqvLbgAOma6SgDNg/5GW+/687f/wBDmoog/wCRlvv+vO3/APQ5qKAL00ywQvK4cqg3EIhdvwVQSfoBWI/idGvEit7OYxBolme4jkt3XzZBGm1HQFuep4AHqeK36hls7WeeKea2hkmhz5UjoCyZ67SeR+FAFUahN/b39n/Z42i8ky+akpLJyAN67cLk7sfMc7Tx1xDaTLL4n1JVDgx2tup3RsoJ3zdMjke44q3BpWnWtybi3sLWGcqEMscKq20AADIGcAADHoBUEH/Iy33/AF52/wD6HNQBpUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABUNxa292gW4hSQKdy7hkqfUHsfcVNRQByOpeArK5bzLR2ibbsALEFVKqmFkGGX5FKjO4DcflNcpqugX+mtPPdxsI5A5uJosIWVstLllGx/kRYlDopJkJwK9ZooA8N8mS3YQ3UCklmLWuxQsrKweRUjY+XIDL5UeInDYQjbUh37gnN1NasSF+d23q4yQpxPC0lwyr8uRhT1Ar1W+8MaZeq+IVhLlS4RRsfBLDchG1huOeR1rjr/wAC39mqCzlE0UQUQkjzPJIUIjbGO4bd80nyN97bheKAOWhkljke5hmWR0ZV+17/AJXdWKxmSVB0MxaT98nIjXJNSCMRIsNquxZgggSRUTzhkxwblOYJtxMkpI2ngnNTTW88MoYI9s8eI45C7k2wKlUUTLiWPZCssjLIGA7jmqcckZiaRj9mtbjO51KRLh0+b5lHks6242gMqHdN60APjJnlit4WZcKptgUZ2hQqVjdY3PmqEiEkv7t2ALLinrsmxcRnyhgIJVlY/ZwyYC+cgEsfl26kkSK65JyKJY12vb3UKxAlnmtjFgJlfMl/cOcEJEqoGhZT8w65phlaaaVruN5ZIlb7QoLySKOHlTdxOm4+VAN28AFuSKAG+fGY3eX/AEa2uNxlZdkY2umZOVBgdxBtjHyod0vrUssQDPBeQqpyzTW3l8LwJJgsLnBVVCRAxSKeOAaGYwb7uWQSNljLcI4CzlWDS4lQFGDTFIws0ZyABkU1S0ZW1ESSHICW7KqJcOsgwBGxMR8y4IyY3XKxHgUAbvhazfVfECi4w8luT5xZ97Ich5RubbKuX2R7XDjarDJr1Sub8F6b9h0jeXkfedqs+8FlUnLbXJKlmLMRnHNdJQBmwf8AIy33/Xnb/wDoc1FEH/Iy33/Xnb/+hzUUAaVFRzGVYXMCI8oHyK7lVJ9CQDj8jXPXF54gGp26S24gUvEFjtM3EUgMgEnmOY1KbU5HTJ7npQB0tZsH/Iy33/Xnb/8Aoc1VVvLyXxNbrCbs2EluWkWS3KKjYBXlowQf+Bk5ONoxkTWhmPifUvNRFUWtvsKuWLDfNyeBg+3P1oA1qKKKACiiigAooooAKKKKACiiigAooooAKKKKACiiigAooooAKKKKACiiigCpfaXZalGUvLZJQVZN3RgpGCAw5GRwcGuY1HwOslxJdWc7ebISZSTskcF97gsBtcHCrtkVgFGOK7KigDxy60fUdMxBdWe9cjEIRVjmfcr7RG2YzvnaMEIVO2I8Cqe0qyxxZuZbc5jDK7NuV9qt5ZImiLzvnMbkYX7vFe2uiSoySKrowwVYZBFc/qng7T7+HbAfszqP3Xy70jIRkUqp5XaHYqFKgHmgDzCGVmkM1tKpdSipcs+RuBKxF5UAON/mTnzU5Ea5Jq7o+mrd3sFvaq0cMwUx7lUbkIKRE4DQzHZ5kpBVG5zmtXUvB9/Zy+ZFGSgyiSI7sYkPyALKuJUCxAj5vMUlula/gPTmLXGpTQiJnOdu1B8zAddh2sVTYoO1Ty+RyaAO1hhjt4I4IUCRRqERR0AAwBT6KKAM2D/kZb7/AK87f/0OaiiD/kZb7/rzt/8A0OaigDSooooAKyJmurXXbi4j064uYZbaFA8LxjDK0pIIZ1PRhWvRQBkzazcQIHfQ9RwXVOGgPLMFH/LT1IqT+07v/oB6h/33B/8AHayZdQ1uPU7+OOO4kPk3HkRtbHylkUKYcOF53DfnLEA4HBwDqaHNdTW05uHuJEWbEEtzD5UjptUksu1cfMXHQcAfUgDv7Tu/+gHqH/fcH/x2j+07v/oB6h/33B/8drSooAzf7Tu/+gHqH/fcH/x2j+07v/oB6h/33B/8drSrN1uW6hso2tmnQGZRLJbxeZIidyq7WzzgdDwSe1AEces3Eryouh6jmJ9jZaAc7Q3H7znhhUn9p3f/AEA9Q/77g/8AjtZVvqOoz3Gho/2+KeSENeobMiAEAhgW2ZDFunzAYGeQRu6egDN/tO7/AOgHqH/fcH/x2j+07v8A6Aeof99wf/Ha0qKAM3+07v8A6Aeof99wf/Hajm1m4gQO+h6jguqcNAeWYKP+WnqRWtXNpqGrRa5dI0F1PHsm8uDyNse5dpiCybcfMu/JLEA4HHAIBpf2nd/9APUP++4P/jtH9p3f/QD1D/vuD/47UXhufUp9JLatHMl2txMp81FQlRI20gKSMbcAcnIGcnqdegDN/tO7/wCgHqH/AH3B/wDHaP7Tu/8AoB6h/wB9wf8Ax2tKigDN/tO7/wCgHqH/AH3B/wDHajh1m4nQumh6jgOyctAOVYqf+WnqDUHia51O2tozpvnhishBhh80tIF/dowwcKT1PGMDkZqKG91WXxQBiYaa4AjjNsVGwx7vMZyOG35TZwehxwaAND+07v8A6Aeof99wf/HaP7Tu/wDoB6h/33B/8drSooAzf7Tu/wDoB6h/33B/8do/tO7/AOgHqH/fcH/x2tKigDJk1m4ieJG0PUcyvsXDQHnaW5/eccKak/tO7/6Aeof99wf/AB2sEalros9RUrem4EaNHm0IEb+ayyKhCEMoUoQcOSMkB8EVv6FNdz6NbyXwk+0ndv8AMTaThiAcYXtj+Fc9dq9AAJ/ad3/0A9Q/77g/+O0f2nd/9APUP++4P/jtaVFAGb/ad3/0A9Q/77g/+O0f2ldD/mBX/wD33B/8drSrJ1+5uLS1tpbZrrf9rhVkt7czbozIocMArEDZuORjp1oAIdZuJ0Lpoeo4DsnLQDlWKn/lp6g1J/ad3/0A9Q/77g/+O1SivdRm8QoFFyLN3GI3tisfkmHdvLFciTzPl2k5x/D3rfoAy7D7TNrF5dzWM1rG9vDGnnNGSxVpSfuM395etFalFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFABRRRQAUUUUAFFFFAH/2Q==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "from dune.fem.view import geometryGridView\n", "position = space.interpolate( x+displacement, name=\"position\" )\n", "beam = geometryGridView( position )\n", "beam.plot()" ] } ], "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 }