{ "cells": [ { "cell_type": "markdown", "id": "dc53f4b3", "metadata": { "lines_to_next_cell": 0 }, "source": [ ".. index:: Solvers; Sicpy (Uzawa Saddle Point)\n", "\n", ".. index::\n", " triple: Equations; Stokes; Saddle-Point;\n", "\n", "# (Quasi) Stokes equations (Scipy based Uzawa scheme)" ] }, { "cell_type": "code", "execution_count": 1, "id": "44da68d0", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:01:25.926952Z", "iopub.status.busy": "2024-02-29T12:01:25.926480Z", "iopub.status.idle": "2024-02-29T12:12:02.256785Z", "shell.execute_reply": "2024-02-29T12:12:02.255469Z" } }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "0.05083519140846924\n", "0.004520464458157929\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "0.0001443814586563057\n", "6.536111368146219e-06\n", "3.445323821267375e-07\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "2.8767428662099426e-08\n", "1.6235031965105895e-09\n", "2.8276576181054936e-12\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "1.252191030596945e-13\n", "9.279809692348484e-15\n" ] } ], "source": [ "from matplotlib import pyplot\n", "import numpy\n", "from scipy.sparse import bmat, linalg\n", "from dune.grid import cartesianDomain\n", "from dune.alugrid import aluCubeGrid\n", "from dune.fem.space import lagrange, dgonb\n", "from dune.fem.operator import galerkin as galerkinOperator\n", "from ufl import SpatialCoordinate, CellVolume, TrialFunction, TestFunction,\\\n", " inner, dot, div, grad, dx, as_vector, transpose, Identity\n", "from dune.ufl import Constant, DirichletBC\n", "import dune.fem\n", "\n", "muValue = 1\n", "nuValue = 0.1\n", "order = 2\n", "\n", "grid = aluCubeGrid(constructor=cartesianDomain([0,0],[3,1],[30,10]))\n", "spcU = lagrange(grid, dimRange=grid.dimension, order=order, storage=\"numpy\")\n", "spcP = lagrange(grid, order=order-1, storage=\"numpy\")\n", "\n", "x = SpatialCoordinate(spcU)\n", "mu = Constant(muValue, \"mu\")\n", "nu = Constant(nuValue, \"nu\")\n", "u = TrialFunction(spcU)\n", "v = TestFunction(spcU)\n", "p = TrialFunction(spcP)\n", "q = TestFunction(spcP)\n", "\n", "exact_u = as_vector( [x[1] * (1.-x[1]), 0] )\n", "exact_p = (-2*x[0] + 2)*mu\n", "f = as_vector( [0,]*grid.dimension )\n", "f += nu*exact_u\n", "\n", "mainModel = (nu*dot(u,v) + mu*inner(grad(u)+grad(u).T, grad(v)) - dot(f,v)) * dx\n", "gradModel = -inner( p*Identity(grid.dimension), grad(v) ) * dx\n", "divModel = -div(u)*q * dx\n", "massModel = p*q * dx\n", "preconModel = inner(grad(p),grad(q)) * dx\n", "\n", "mainOp = galerkinOperator( (mainModel, DirichletBC(spcU,exact_u)), spcU)\n", "gradOp = galerkinOperator( gradModel, spcP,spcU)\n", "divOp = galerkinOperator( divModel, spcU,spcP)\n", "massOp = galerkinOperator( massModel, spcP)\n", "\n", "A = mainOp.linear()\n", "A = A.as_numpy\n", "def Ainv(rhs,target): target[:] = linalg.spsolve(A,rhs)\n", "G = gradOp.linear()\n", "G = G.as_numpy\n", "D = divOp.linear()\n", "D = D.as_numpy\n", "M = massOp.linear()\n", "M = M.as_numpy\n", "def Minv(rhs,target): target[:] = linalg.spsolve(M,rhs)\n", "\n", "if mainOp.model.nu > 0: # same as nu.value>0\n", " preconOp = galerkinOperator( (preconModel, DirichletBC(spcP,0)), spcP)\n", " P = preconOp.linear().as_numpy\n", " def Pinv(rhs,target): target[:] = linalg.spsolve(P,rhs)\n", "\n", "# discrete functions\n", "velocity = spcU.interpolate(spcU.dimRange*[0], name=\"velocity\")\n", "pressure = spcP.interpolate(0, name=\"pressure\")\n", "sol_u = velocity.as_numpy\n", "sol_p = pressure.as_numpy\n", "\n", "# aux. function\n", "rhsVelo = velocity.copy()\n", "rhsPress = pressure.copy()\n", "rhs_u = rhsVelo.as_numpy\n", "rhs_p = rhsPress.as_numpy\n", "r = numpy.copy(rhs_p)\n", "d = numpy.copy(rhs_p)\n", "precon = numpy.copy(rhs_p)\n", "xi = numpy.copy(rhs_u)\n", "\n", "# right hand side for Shur complement problem\n", "mainOp(velocity,rhsVelo)\n", "rhs_u *= -1\n", "xi[:] = G*sol_p\n", "rhs_u -= xi\n", "mainOp.setConstraints(rhsVelo)\n", "Ainv(rhs_u[:], sol_u[:])\n", "rhs_p[:] = D*sol_u\n", "Minv(rhs_p, r)\n", "if mainOp.model.nu > 0:\n", " precon.fill(0)\n", " Pinv(rhs_p, precon)\n", " r *= mainOp.model.mu\n", " r += mainOp.model.nu*precon\n", "d[:] = r[:]\n", "delta = numpy.dot(r,rhs_p)\n", "\n", "# cg type iteration\n", "for m in range(100):\n", " xi.fill(0)\n", " rhs_u[:] = G*d\n", " mainOp.setConstraints([0,]*grid.dimension, rhsVelo)\n", " Ainv(rhs_u[:], xi[:])\n", " rhs_p[:] = D*xi\n", " rho = delta / numpy.dot(d,rhs_p)\n", " sol_p += rho*d\n", " sol_u -= rho*xi\n", " rhs_p[:] = D*sol_u\n", " Minv(rhs_p[:],r[:])\n", " if mainOp.model.nu > 0:\n", " precon.fill(0)\n", " Pinv(rhs_p,precon)\n", " r *= mainOp.model.mu\n", " r += mainOp.model.nu*precon\n", " oldDelta = delta\n", " delta = numpy.dot(r,rhs_p)\n", " print(delta)\n", " if delta < 1e-14: break\n", " gamma = delta/oldDelta\n", " d *= gamma\n", " d += r" ] }, { "cell_type": "markdown", "id": "c9403f27", "metadata": { "lines_to_next_cell": 0 }, "source": [ "Now we can easily plot the velocity magnitude and pressure" ] }, { "cell_type": "code", "execution_count": 2, "id": "508259a4", "metadata": { "execution": { "iopub.execute_input": "2024-02-29T12:12:02.264268Z", "iopub.status.busy": "2024-02-29T12:12:02.263135Z", "iopub.status.idle": "2024-02-29T12:12:02.979170Z", "shell.execute_reply": "2024-02-29T12:12:02.978476Z" } }, "outputs": [ { "data": { "image/jpeg": "/9j/4AAQSkZJRgABAQEAMgAyAAD/2wBDAAgGBgcGBQgHBwcJCQgKDBQNDAsLDBkSEw8UHRofHh0aHBwgJC4nICIsIxwcKDcpLDAxNDQ0Hyc5PTgyPC4zNDL/2wBDAQkJCQwLDBgNDRgyIRwhMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjIyMjL/wAARCADbA+MDASIAAhEBAxEB/8QAHwAAAQUBAQEBAQEAAAAAAAAAAAECAwQFBgcICQoL/8QAtRAAAgEDAwIEAwUFBAQAAAF9AQIDAAQRBRIhMUEGE1FhByJxFDKBkaEII0KxwRVS0fAkM2JyggkKFhcYGRolJicoKSo0NTY3ODk6Q0RFRkdISUpTVFVWV1hZWmNkZWZnaGlqc3R1dnd4eXqDhIWGh4iJipKTlJWWl5iZmqKjpKWmp6ipqrKztLW2t7i5usLDxMXGx8jJytLT1NXW19jZ2uHi4+Tl5ufo6erx8vP09fb3+Pn6/8QAHwEAAwEBAQEBAQEBAQAAAAAAAAECAwQFBgcICQoL/8QAtREAAgECBAQDBAcFBAQAAQJ3AAECAxEEBSExBhJBUQdhcRMiMoEIFEKRobHBCSMzUvAVYnLRChYkNOEl8RcYGRomJygpKjU2Nzg5OkNERUZHSElKU1RVVldYWVpjZGVmZ2hpanN0dXZ3eHl6goOEhYaHiImKkpOUlZaXmJmaoqOkpaanqKmqsrO0tba3uLm6wsPExcbHyMnK0tPU1dbX2Nna4uPk5ebn6Onq8vP09fb3+Pn6/9oADAMBAAIRAxEAPwD3yV2jhd1ieVlGRGhGW9hkgfmRWTL4ltItOs73ybh47m2N3hVXMcIClnbLdt65AyeeAa1pY1mieNiwVwVJRip/AjkfUVmJ4a0pLOC0EMxggXZGr3MrYTABTJbJQhVyp+U45FAE+pa5pOjBDqmqWViH4U3VwkW76biM1n/8Jz4R/wChq0P/AMGEX/xVcX8eoGufBWnRoMk6oh/8gzV8+f2NP/cr6HK8iWPoe19py62tb0OariYUpcrPrf8A4Tnwj/0NWh/+DCL/AOKo/wCE58I/9DVof/gwi/8Aiq+Rzo8+5Rt60v8AY0/9yvRXCaf/AC+X3f8ABM/r1PufW/8AwnPhH/oatD/8GEX/AMVR/wAJz4R/6GrQ/wDwYRf/ABVfI/8AY8+/G3tS/wBjT/3KFwmn/wAvvw/4IfXqfc+t/wDhOfCP/Q1aH/4MIv8A4qj/AITnwj/0NWh/+DCL/wCKr5H/ALHn3429qX+xp/7lC4TT/wCX34f8EPr1PufW/wDwnPhH/oatD/8ABhF/8VR/wnPhH/oatD/8GEX/AMVXyONHn3MNvSl/saf+5QuE0/8Al8vu/wCCH16n3Prf/hOfCP8A0NWh/wDgwi/+Ko/4Tnwj/wBDVof/AIMIv/iq+Rxo8+5ht6Uv9jT/ANyhcJp/8vl93/BD69T7n1v/AMJz4R/6GrQ//BhF/wDFUf8ACc+Ef+hq0P8A8GEX/wAVXyONHn3MNvSl/saf+5QuE0/+Xy+7/gh9ep9z63/4Tnwj/wBDVof/AIMIv/iqP+E58I/9DVof/gwi/wDiq+R10ec5+XvS/wBjT/3KFwmmv434f8EPr1PufW//AAnPhH/oatD/APBhF/8AFUf8Jz4R/wChq0P/AMGEX/xVfI66POc/L3pf7Gn/ALlC4TTX8b8P+CH16n3Prf8A4Tnwj/0NWh/+DCL/AOKo/wCE58I/9DVof/gwi/8Aiq+R10ecqDtpf7Gn/uULhNNX9svu/wCCH16n3Prf/hOfCP8A0NWh/wDgwi/+Ko/4Tnwj/wBDVof/AIMIv/iq+R10ecqDtpTo84U/LQuE1a/tvw/4IfXqfc+t/wDhOfCP/Q1aH/4MIv8A4qj/AITnwj/0NWh/+DCL/wCKr5IGjzlR8tB0ecKflo/1TVr+2/D/AIIfXqfc+t/+E58I/wDQ1aH/AODCL/4qj/hOfCP/AENWh/8Agwi/+Kr5IGjzlR8tI2jzhSdtD4TVr+2/D/gh9ep9z64/4Tnwj/0NWh/+DCL/AOKo/wCE58I/9DVof/gwi/8Aiq+SP7Gn/uUjaPOFJ20PhNJX9svu/wCCH16n3Prj/hOfCP8A0NWh/wDgwi/+Ko/4Tnwj/wBDVof/AIMIv/iq+SP7Gn/uUjaPOFJ20PhNJX9svu/4IfXqfc+uP+E58I/9DVof/gwi/wDiqP8AhOfCP/Q1aH/4MIv/AIqvkj+xp/7lI2jzjHy96Hwmkr+2/D/gh9ep9z64/wCE58I/9DVof/gwi/8AiqP+E58I/wDQ1aH/AODCL/4qvkj+xp/7lI2jzjHy96Hwmkr+2/D/AIIfXqfc+uP+E58I/wDQ1aH/AODCL/4qj/hOfCP/AENWh/8Agwi/+Kr5I/saf+5SHR59yjb1ofCaX/L5fd/wQ+vU+59cf8Jz4R/6GrQ//BhF/wDFUf8ACc+Ef+hq0P8A8GEX/wAVXyR/Y0/9ykOjz7lG3rQ+E0v+Xy+7/gh9ep9z64/4Tnwj/wBDVof/AIMIv/iqP+E58I/9DVof/gwi/wDiq+SP7Gn/ALlJ/Y8+/G3tQ+E0v+X34f8ABD69T7n1x/wnPhH/AKGrQ/8AwYRf/FUf8Jz4R/6GrQ//AAYRf/FV8kf2NP8A3KT+x59+Nvah8Jpf8vvw/wCCH16n3Prj/hOfCP8A0NWh/wDgwi/+Ko/4Tnwj/wBDVof/AIMIv/iq+SP7Gn/uUn9jz78be1D4TS/5ffh/wQ+vU+59cf8ACc+Ef+hq0P8A8GEX/wAVR/wnPhH/AKGrQ/8AwYRf/FV8kf2NP/cpBo8+5ht6UPhNf8/l93/BD69T7n1x/wAJz4R/6GrQ/wDwYRf/ABVH/Cc+Ef8AoatD/wDBhF/8VXyR/Y0/9ykGjz7mG3pQ+E1/z+X3f8EPr1PufXH/AAnPhH/oatD/APBhF/8AFUf8Jz4R/wChq0P/AMGEX/xVfJH9jT/3KRdHnOfl70f6pq/8b8P+CH16n3Prj/hOfCP/AENWh/8Agwi/+Ko/4Tnwj/0NWh/+DCL/AOKr5I/saf8AuUi6POc/L3o/1TV/434f8EPr1PufXH/Cc+Ef+hq0P/wYRf8AxVH/AAnPhH/oatD/APBhF/8AFV8kf2NP/cpF0ecqDto/1TV7e2X3f8EPr1PufXH/AAnPhH/oatD/APBhF/8AFUf8Jz4R/wChq0P/AMGEX/xVfJH9jT/3KRdHnKg7aP8AVNXt7Zfd/wAEPr1PufXH/Cc+Ef8AoatD/wDBhF/8VR/wnPhH/oatD/8ABhF/8VXyR/Y0/wDcpF0ecqDto/1TV7e2X3f8EPr1PufXH/Cc+Ef+hq0P/wAGEX/xVbkUsc8KTQyLJFIoZHQ5DA8ggjqK+LDo8+PuV9ieHv8AkWtK/wCvOH/0AV42b5Ssu5Pf5ua/4W/zNqNeNW/L0L8rtHC7rE8rKMiNCMt7DJA/MismXxLaRadZ3vk3Dx3NsbvCquY4QFLO2W7b1yBk88A1rSxrNE8bFgrgqSjFT+BHI+orMTw1pSWcFoIZjBAuyNXuZWwmACmS2ShCrlT8pxyK8Y3J9R1vSdHCHU9UsrEP9z7TcJFu+m4jNUP+E38J/wDQ0aL/AODCL/4quJ+Pn/Ik6d/2FE/9EzV8/JWFWtyO1j1cFlqxMOdytr2Prf8A4Tfwn/0NGi/+B8X/AMVS/wDCbeFP+hn0X/wPi/8Aiq+Sx95fxqwvWsXi2uh6UOHoy/5efh/wT6s/4TXwp/0M+i/+B8X/AMVS/wDCa+FP+hm0b/wPi/8Aiq+Vl/1n/AasJ0qHjmvsnTT4WjL/AJe/h/wT6h/4TTwqf+Zm0b/wPi/+Kpf+Ez8Lf9DLo/8A4HRf/FV8yQf67/gP9a0I6ylmTj9n8TZcIwf/AC9/D/gn0X/wmXhf/oZNH/8AA6L/AOKpf+Ex8Mf9DHpH/gdF/wDFV89J/wAfH/Af61cj6VjLN2vsfj/wDJ8KxX/L38P+Ce8/8Jh4YP8AzMekf+B0X/xVL/wl/hn/AKGLSP8AwNj/APiq8Lg/10n4fyq6lZSzuS+x+P8AwDF8NRX/AC8/D/gns/8Awl3hk/8AMxaT/wCBsf8A8VS/8Jb4b/6GHSf/AANj/wAa8dg/18v4fyq7H0rKWfyX/Lv8f+AYvIIr/l5+H/BPVR4t8Nnp4g0r/wADY/8AGl/4Svw5/wBB/Sv/AAMj/wAa8utOsn/XQ1ejrKXEcov+H+P/AADF5LFfb/D/AIJ6GPFfhw9Nf0r/AMDI/wDGl/4Snw9/0HtL/wDAyP8Axrz+06S/9dDWhH1rKXFEov8Ahfj/AMAxeVJfa/A7AeKfDxGRr2ln/t8j/wAaX/hKPD//AEHdM/8AAuP/ABrkLD/j1T8f51ox1lPiyUXb2X4/8AweXpfaN4eJ/D5GRrumEf8AX3H/AI0f8JNoH/Qc03/wLj/xrD0//j0j/H+dXlrOXF807ex/H/gGE8MorcvDxPoBGRrmmf8AgXH/AI0f8JPoA665pn/gXH/jWVZ/8eifj/OluP8Aj3l/3T/Kto8VSbt7L8f+AcslY1P+En0D/oOaZ/4Fx/40f8JP4fHXXNM/8C4/8ayIf+PeP/dH8qbcf6iT/cP8q6I8Ryk7ez/H/gGEqtuhs/8ACT+H/wDoO6Z/4Fx/40h8UeHwMnXdMH/b3H/jWJD/AMe8f+6P5VDe/wDHq/4fzrphnkpO3J+P/AOeWMa6HRf8JR4f/wCg7pn/AIFx/wCNIfFHh4DJ13TB/wBvcf8AjWGaq3v/AB6v+H866IZq5fZ/EwlmTX2fxOm/4Sjw/wD9B3TP/AuP/Gg+KfDw669pf/gZH/jXOmqt1/yx/wCuq11Qxrl0OZ5y19j8f+AdZ/wlPh7/AKDul/8AgXH/AI0HxV4dHXXtL/8AAyP/ABrmDVW6/wCWX/XQV0wq83QwfEEl/wAu/wAf+Adj/wAJT4e/6D2l/wDgZH/jSHxV4dHXXtL/APAyP/GuTNVrn/ll/wBdBXRGPMc74mkv+Xf4/wDAO2/4Srw7/wBB7S//AAMj/wAaQ+K/DgODr+lf+Bkf+Nccap3P/HzD/wAC/lXTDCqXUw/1sle3svx/4B3v/CV+HP8AoYNK/wDAyP8AxpP+Et8Nj/mYdJ/8DY/8a88eqFx/x8Rfj/Ku2OVqS+L8AhxZKTt7L8f+Aepf8Jb4a/6GHSf/AANj/wDiqT/hL/DOcf8ACRaT/wCBsf8A8VXlMnSqL/8AH1/wD+tbLJk/t/h/wTaHE8pf8uvx/wCAex/8Jf4Z/wChi0j/AMDY/wD4qj/hMfDH/Qx6R/4HRf8AxVeLSVRl/wCPn/gH9a1WQxf2/wAP+Cbw4hlL/l3+P/APdv8AhMfC/wD0Mmj/APgdF/8AFUf8Jn4W/wChl0f/AMDov/iq8Ck6VSP+tk/Ctlw7F/8ALz8P+CbxzuUvsfj/AMA+if8AhM/C3/Qy6N/4Hxf/ABVH/CaeFf8AoZtG/wDA+L/4qvm6SqT/AOsf8K1XDEX/AMvfw/4JvHNXL7H4n07/AMJr4U/6GfRf/A+L/wCKo/4Tbwn/ANDPov8A4Hxf/FV8tydKpt99vwrVcJxdv3v4f8E2jmDf2T6v/wCE38Jf9DRon/gwi/8AiqP+E48Jf9DTon/gwi/+Kr5Hkqo3U/WtlwfBv+M/u/4JvHFN9D7D/wCE48Jf9DTon/gwi/8AiqT/AITjwj/0NOif+DCL/wCKr48pnr9aU+D4R/5ffh/wS1iPI+xv+E58I/8AQ06J/wCDCL/4qj/hOfCP/Q06J/4MIv8A4qvjg0ztXNPheMXb2v4f8Eft/I+yv+E58I/9DTon/gwi/wDiqP8AhOfCP/Q1aH/4MIv/AIqvjM0ztXNPh+MXb2n4f8Er2vkfZ/8AwnXhH/oatD/8GEX/AMVW5FLHPCk0MiyRSKGR0OQwPIII6ivhM19u+Hv+Ra0r/rzh/wDQBXk4/BLCuNpXvcuMuY0qKKK88sKKKKAOI+J1r9s0jSoSpbOog4Az/wAsJq8+/wCEeX/n3l/75r1Pxgm9NHX/AKfz/wCiJqyfs3tXVRzz6jH2V/M+Qz2hXnik6e1l+p56+gKJ4h9nk5z/AA+1S/8ACPL/AM+8v/fNdtLbf6VBx/e/lU/2b2rT/WrzPHeFxdlqee/2Av2rH2eT7mfu+9S/8I8v/PvL/wB81232b/Tun/LL+tT/AGb2o/1qt1B4XF9zz3+wF+1Y+zyfcz933qX/AIR5f+feX/vmu2+zf6d0/wCWX9an+ze1H+tVuoPC4vueepoCmeUfZ5OMfw+1S/8ACPL/AM+8v/fNdtFbf6VPx/d/lU/2b2o/1qt1B4XF33PPU0BTPKPs8nGP4fapf+EeX/n3l/75rtorb/Sp+P7v8qn+ze1H+tVuoPC4u+556mgKZ5R9nk4x/D7VL/wjy/8APvL/AN8120Vt/pU/H93+VT/Zvaj/AFqt1B4XF33PPYdAU+Z/o8nDkfdqX/hHl/595f8Avmu2trb/AF3H/LVqn+ze1H+tVuoSwuLvueew6Ap8z/R5OHI+7Uv/AAjy/wDPvL/3zXbW1t/ruP8Alq1T/Zvaj/Wq3UJYXF33PPbfQFaBT9nk7/w+9S/8I8v/AD7y/wDfNdtZ23+ipx6/zqf7N7Uf61W0uEsLi7vU89t9AVoFP2eTv/D70+Tw+oic/Z5fun+Gu4s7b/RU49f50+a2/cScfwn+VH+tXS4PC4vm3OEj8PqYkP2eX7o/hok8PqInP2eX7p/hru4bb9xHx/CP5UTW37iTj+E/yo/1q6XD6ri77nCR+H1MSH7PL90fw0y40BVgY/Z5O38PvXfQ237iPj+Efypl5bf6K/Hp/Oj/AFq6XBYXF825xP8Awjy/8+8v/fNRXGgKsDH7PJ2/h969C+ze1QXlt/or8en86P8AWq+lwjhcXdanE/8ACPL/AM+8v/fNRXGgKsDH7PJ2/h969C+ze1QXlt/or8en86P9ar6XCOFxd1qcT/wjy/8APvL/AN81FNoCjy/9Hk5cD7tehfZvaoLm2/1PH/LVaP8AWq/UI4XF33OJ/wCEeX/n3l/75qKbQFHl/wCjycuB92vQvs3tUFzbf6nj/lqtH+tV+oRwuLvucT/wjy/8+8v/AHzUT6AoniH2eTnP8PtXoX2b2qCW2/0qDj+9/Kj/AFqv1COFxd9zif8AhHl/595f++aifQFE8Q+zyc5/h9q9C+ze1QS23+lQcf3v5Uf61X6hHC4u+5xP/CPL/wA+8v8A3zUX9gL9qx9nk+5n7vvXoX2b2qD7N/p3T/ll/Wj/AFqv1BYXF9zif+EeX/n3l/75qL+wF+1Y+zyfcz933r0L7N7VB9m/07p/yy/rR/rVfqCwuL7nE/8ACPL/AM+8v/fNRf2Av2rH2eT7mfu+9ehfZvaoPs3+ndP+WX9aP9ar9QWFxfc4n/hHl/595f8Avmok0BTPKPs8nGP4favQvs3tUEVt/pU/H93+VH+tXmCwuLs9Tif+EeX/AJ95f++aiTQFM8o+zycY/h9q9C+ze1QRW3+lT8f3f5Uf61eYLC4uz1OJ/wCEeX/n3l/75qKHQFPmf6PJw5H3a9C+ze1QW1t/ruP+WrUf61eYLC4u25xP/CPL/wA+8v8A3zUUOgKfM/0eThyPu16F9m9qgtrb/Xcf8tWo/wBavMFhcXbc4n/hHl/595f++ait9AVoFP2eTv8Aw+9ehfZvaoLO2/0VOPX+dH+tXW4fVcXbc4n/AIR5f+feX/vmorfQFaBT9nk7/wAPvXoX2b2qCztv9FTj1/nR/rV1uH1XF23OJ/4R5f8An3l/75qK30BWgU/Z5O/8PvXoX2b2qCztv9FTj1/nR/rV1uH1XF23OIbw+uw/6PL0/u17P4e/5FrSv+vOH/0AVyj237tuOxrq/D3/ACLWlf8AXnD/AOgCsK+a/wBoW1+H9f8Ahj6Ph6lVp+09r5fqaVFFFYH0h5V8fTjwRp3/AGFE/wDRM1fPqN7H8q+g/j5/yJOnf9hRP/RM1fPyVxYn4j6fJr+x+f8AkODfMvynv2qwrf7LflUI+8v41YXrXHI+iop9xyv+8+633fSp1f8A2W/Kol/1n/AasJ0rGR6NBPuSQyfvvuN930960I5P9h/yqlB/rv8AgP8AWtCOuaodcU7PUekv7/7j/d9PerscvH+rf8qrJ/x8f8B/rVyPpXLUsc0k9dR0Ev76T93J2/h9qupN/wBM5P8Avmq0H+uk/D+VXUrlqNHJJO24sE37+X91J2/h9quxz/8ATKX/AL5qvB/r5fw/lV2PpXNUaOOSdhLWfmT91L/rD/DV+O4/6ZS/981WtOsn/XQ1ejrmqNX2OSSdhtpccS/uZf8AWH+Gr8dxz/qZv++arWnSX/roa0I+tc1Vq+xxyvYjsbn/AEZP3M3f+H3q+lz/ANMJv++KrWH/AB6p+P8AOtGOuaq1zPQ45bENhdf6In7ibv8Awe9Xluv+mE//AHxUOn/8ekf4/wA6vLXPNrmehw1dihZ3X+iJ+4n7/wAHvTri6/0eT9xP90/we1TWf/Hon4/zpbj/AI95f90/yrrg1zbHm1CrDdfuI/3E/wB0fwe1NuLr9xJ+4n+6f4ParMP/AB7x/wC6P5U24/1En+4f5V30mr7HDUKsN1+4j/cT/dH8HtUN7c/6K/7ibt/B71dh/wCPeP8A3R/Kob3/AI9X/D+dehSaucNTcYbn/phN/wB8VWvbn/RX/cTdv4PetE1Vvf8Aj1f8P516FG10cM9xhuf+mE3/AHxVa6uf9T+5m/1g/hrQNVbr/lj/ANdVr0qR587XGm5/6Yzf98VWubn/AFX7mb/WD+GtA1Vuv+WX/XQV6FI4ZNXGm5/6Yzf981Wubj/VfuZv9YP4avmq1z/yy/66CvRpnDJq+w03H/TGb/vmqdzcf6TD+5l/i/h9q0jVO5/4+Yf+Bfyr0KRxNq+xTe45/wBVL/3zVC4n/wBIi/dS9/4fatR6oXH/AB8Rfj/KvXp3sFJq+xVkn4P7qX/vmqLzf6V/q5Puf3fetOTpVF/+Pr/gH9a7I3Oyk12Kkk3/AEzk/wC+apSS/wCkf6t/uenvWlJVGX/j5/4B/WuqKZ2UmuxUkl4+4/5VSMn7yT5H7dq0ZOlUj/rZPwrqindanbTasVJJP9h/yqm7/vH+Vu3ar8lUn/1j/hXVFO61O2nYrSP/ALLflVNm+dvlbt2q9J0qm332/CuqKd1qdlMqyN/sn8qqM3J4PWrklVG6n611RTutTtpjc+xpmevB61JTPX60qqd1qaoaT7GmA8dDUhpnavOqp33LQwn2NNzx0NPNM7V59VO+5SGk+xr7e8Pf8i1pX/XnD/6AK+IjX274e/5FrSv+vOH/ANAFfK55vD5/ob0zSooorwTUqale/wBn6fLdCPzGTAVN20FiQBk9hkjJ7ViXnix7Oyt5jYB3eSdJUWViIxE+x2yEPy5/iYKo43FciujlijnieKWNZI3Uq6OMhgeoI7iqh0XSjGsZ0yy2I+9V8hcK2AMgY64AGfYUAZHiu4SO40aJklZjeFvljYjHkS/xYxn2zmqnnx/88Zv++Km8cT/Z7fSJM4xf4/8AIE1c7/a3+0a+F4jdVYxcm3KvzZ6GFy2OJhztGrLOn2u3/czfxfwe1T+fH/zxm/74rnZNV/0mA7j/ABfyqb+1v9o14TlXstTf+xIdjV89Pt/+pm/1X9z3qfz4/wDnjN/3xXO/2r/pudx/1f8AWpv7W/2jQ5V9NQ/sSHY1fPT7f/qZv9V/c96n8+P/AJ4zf98Vzv8Aav8Apudx/wBX/Wpv7W/2jQ5V9NQ/sSHY1Yp0+13H7mb+H+D2qfz4/wDnjN/3xXOx6r/pM53H+H+VTf2t/tGiUq99x/2JDsasU6fa7j9zN/D/AAe1T+fH/wA8Zv8Aviudj1X/AEmc7j/D/Kpv7W/2jRKVe+4f2JDsasU6fa7j9zN/D/B7VP58f/PGb/viudj1X/SZzuP8P8qm/tb/AGjRKVe+4f2JDsattOn779zN/rW/gqfz4/8AnjN/3xXO2+q4835j/rDU39rf7RolKvfcFkkH0NW2nT99+5m/1rfwVP58f/PGb/viudt9Vx5vzH/WGpv7W/2jRKVe+4LJIPoatnOn2RP3M3f+D3qfz4/+eM3/AHxXO2mq4tkG49/51N/a3+0aJSr8z1Eskg+hq2c6fZE/czd/4PepJ54/s8n7mb7p/g9qwLTVcWyDce/86fNquYJPmP3T/Km5V+bcP7Eha9jcgnj+zx/uZvuj+D2onnj+zyfuZvun+D2rDh1XEEfzH7o/lRNquYJPmP3T/KlzV+bcP7Eh2NyCeP7PH+5m+6P4Pao7ydPsj/uZu38HvWPDquII/mP3R/KmXeq5tnG49v501Kvzbh/YkLXsdF58f/PGb/vioLydPsj/ALmbt/B71lf2t/tGobvVc2zjce386UZV+ZajeSQXQ6Lz4/8AnjN/3xUF5On2R/3M3b+D3rK/tb/aNQ3eq5tnG49v50RlX5lqDySC6HRefH/zxm/74qC5nT9z+5m/1q/wVlf2t/tGobjVc+V8x/1gojKvfcP7Eguh0Xnx/wDPGb/vioLmdP3P7mb/AFq/wVlf2t/tGobjVc+V8x/1gojKvfcP7Eguh0Xnx/8APGb/AL4qCWdPtdv+5m/i/g9qyv7W/wBo1DJqv+kwHcf4v5URlXvuL+xIdjovPj/54zf98VBLOn2u3/czfxfwe1ZX9rf7RqGTVf8ASYDuP8X8qIyr33D+xIdjovPj/wCeM3/fFQeen2//AFM3+q/ue9ZX9rf7RqH+1f8ATc7j/q/60RlX7h/YkOx0Xnx/88Zv++Kg89Pt/wDqZv8AVf3Pesr+1v8AaNQ/2r/pudx/1f8AWiMq/cP7Eh2Oi8+P/njN/wB8VB56fb/9TN/qv7nvWV/a3+0ah/tX/Tc7j/q/60RlX7j/ALEh2Oi8+P8A54zf98VBFOn2u4/czfw/we1ZX9rf7RqGPVf9JnO4/wAP8qFKvZ6h/YkOx0Xnx/8APGb/AL4qCKdPtdx+5m/h/g9qyv7W/wBo1DHqv+kzncf4f5UKVez1D+xIdjovPj/54zf98VBbTp++/czf61v4Kyv7W/2jUNvquPN+Y/6w0KVez1D+xIdjovPj/wCeM3/fFQW06fvv3M3+tb+Csr+1v9o1Db6rjzfmP+sNClXs9Q/sSHY6Lz4/+eM3/fFQWc6fZE/czd/4Pesr+1v9o1Daari2Qbj3/nRzV+XcX9iQ7HRefH/zxm/74qCznT7In7mbv/B71lf2t/tGobTVcWyDce/86Oavy7h/YkOx0Xnx/wDPGb/vioLOdPsifuZu/wDB71lf2t/tGobTVcWyDce/86Oavy7h/YkOx0Dzx+W37mbof4K6Dw9/yLWlf9ecP/oArhH1X5G+Y9DXd+Hv+Ra0r/rzh/8AQBX1fC7qP2vP/d/U5MXgVhbWW/6E+pXv9n6fLdCPzGTAVN20FiQBk9hkjJ7ViXnix7Oyt5jYB3eSdJUWViIxE+x2yEPy5/iYKo43FciujlijnieKWNZI3Uq6OMhgeoI7iqh0XSjGsZ0yy2I+9V8hcK2AMgY64AGfYV9YcZ5r8fZlHg/TYcPuOpI2Qh248qUfexjPtnNeAI3sfyr6D+Pn/Ik6d/2FE/8ARM1fPyVxYn4j6fJr+x+f+Q4N8y/Ke/arCt/st+VQj7y/jVhetccj6Kin3HK/7z7rfd9KnV/9lvyqJf8AWf8AAasJ0rGR6NBPuSQyfvvuN930960I5P8AYf8AKqUH+u/4D/WtCOuaodcU7PUekv7/AO4/3fT3q7HLx/q3/Kqyf8fH/Af61cj6Vy1LHNJPXUdBL++k/dydv4farqTf9M5P++arQf66T8P5VdSuWo0ckk7biwTfv5f3Unb+H2q7HP8A9Mpf++arwf6+X8P5Vdj6VzVGjjknYS1n5k/dS/6w/wANX47j/plL/wB81WtOsn/XQ1ejrmqNX2OSSdhtpccS/uZf9Yf4avx3HP8AqZv++arWnSX/AK6GtCPrXNVavsccr2I7G5/0ZP3M3f8Ah96vpc/9MJv++KrWH/Hqn4/zrRjrmqtcz0OOWxDYXX+iJ+4m7/we9Xluv+mE/wD3xUOn/wDHpH+P86vLXPNrmehw1dihZ3X+iJ+4n7/we9OuLr/R5P3E/wB0/wAHtU1n/wAeifj/ADpbj/j3l/3T/KuuDXNsebUKsN1+4j/cT/dH8HtTbi6/cSfuJ/un+D2qzD/x7x/7o/lTbj/USf7h/lXfSavscNQqw3X7iP8AcT/dH8HtUN7c/wCiv+4m7fwe9XYf+PeP/dH8qhvf+PV/w/nXoUmrnDU3GG5/6YTf98VWvbn/AEV/3E3b+D3rRNVb3/j1f8P516FG10cM9xhuf+mE3/fFVrq5/wBT+5m/1g/hrQNVbr/lj/11WvSpHnztcabn/pjN/wB8VWubn/VfuZv9YP4a0DVW6/5Zf9dBXoUjhk1cabn/AKYzf981Wubj/VfuZv8AWD+Gr5qtc/8ALL/roK9GmcMmr7DTcf8ATGb/AL5qnc3H+kw/uZf4v4fatI1Tuf8Aj5h/4F/KvQpHE2r7FN7jn/VS/wDfNULif/SIv3Uvf+H2rUeqFx/x8Rfj/KvXp3sFJq+xVkn4P7qX/vmqLzf6V/q5Puf3fetOTpVF/wDj6/4B/WuyNzspNdipJN/0zk/75qlJL/pH+rf7np71pSVRl/4+f+Af1rqimdlJrsVJJePuP+VUjJ+8k+R+3atGTpVI/wCtk/CuqKd1qdtNqxUkk/2H/Kqbv+8f5W7dqvyVSf8A1j/hXVFO61O2nYrSP/st+VU2b52+Vu3ar0nSqbffb8K6op3Wp2UyrI3+yfyqozcng9auSVUbqfrXVFO61O2mNz7GmZ68HrUlM9frSqp3WpqhpPsaYDx0NSGmdq86qnfctDCfY03PHQ080ztXn1U77lIaT7Gvt7w9/wAi1pX/AF5w/wDoAr4iNfbvh7/kWtK/684f/QBXyuebw+f6G9M0qKKK8E1CiiigDg/irP8AZtD0uQAnGogcf9cZq80/tg/3H/KvQPjPL5PhfTH/AOokv/omavGv7Qr57NKCqV07dD7vhugqmDbf8z/JHQPq586I7H79vapf7YP9x/yrl21D94n40/8AtCvOeEXY99YSN3qdB/a5+052P9z096l/tg/3H/KuX/tD99n/AGaf/aFDwi7BHCR7nQf2uftOdj/c9Pepf7YP9x/yrl/7Q/fZ/wBmn/2hQ8IuwRwke50CaufOlOx+3b2qX+2D/cf8q5ddQ/eP+FP/ALQoeEXYI4SNtzoE1c+dKdj9u3tUv9sH+4/5Vy66h+8f8Kf/AGhQ8IuwRwkbbnQJq586U7H7dvapf7YP9x/yrl11D94/4U/+0KHhF2COEjbc6CHVyPM+R/vntUv9sH+4/wCVcvHqGN3+8af/AGhQ8Ir7BHCRtudBDq5HmfI/3z2qX+2D/cf8q5ePUMbv940/+0KHhFfYI4SNtzoINXIhUbH79vepf7YP9x/yrl4tQxGKf/aFDwivsEcJHlWp0EGrkQqNj9+3vT5NXPlv8j9D2rmYtQxGKV9Q+RvpTeEV9hLCR5dzpY9XPlJ8j9B2ok1c+W/yP0PauaTUPkX6UPqHyN9KX1RX2H9Ujy7nSx6ufKT5H6DtTJ9XJhYbH7dveudTUPkX6UkuoZjNNYRX2E8JHl3Oo/tg/wBx/wAqin1cmFhsft2965/+0KZLqGYzSWEV9hywkeV6nUf2wf7j/lUU+rkwsNj9u3vXP/2hTJdQzGaFhFfYJYSPK9TqP7YP9x/yqKbVyfL+R/vjtXP/ANoUyTUM7f8AeFCwivsEsJG251H9sH+4/wCVRTauT5fyP98dq5/+0KZJqGdv+8KFhFfYJYSNtzqP7YP9x/yqJ9XPnRHY/ft7Vz/9oUxtQ/eJ+NCwi7BLCRtudR/bB/uP+VRPq586I7H79vauf/tCmNqH7xPxoWEXYJYSNtzqP7YP9x/yqL+1z9pzsf7np71z/wDaFM/tD99n/ZoWEXYHhI6anUf2wf7j/lUX9rn7TnY/3PT3rn/7Qpn9ofvs/wCzQsIuwPCR01Oo/tg/3H/Kov7XP2nOx/uenvXP/wBoUz+0P32f9mhYRdgeEjpqdR/bB/uP+VRJq586U7H7dvauf/tCmLqH7x/woWEXYHhI3Wp1H9sH+4/5VEmrnzpTsft29q5/+0KYuofvH/ChYRdgeEjdanUf2wf7j/lUUOrkeZ8j/fPauf8A7QpkeoY3f7xoWEVtgeEjdanUf2wf7j/lUUOrkeZ8j/fPauf/ALQpkeoY3f7xoWEVtgeEjdanUf2wf7j/AJVFBq5EKjY/ft71z/8AaFMi1DEYo+qK2wfVI8251H9sH+4/5VFBq5EKjY/ft71z/wDaFMi1DEYo+qK2wfVI8251H9sH+4/5VFBq5EKjY/ft71z/APaFMi1DEYo+qK2wfVI82507awdp+R+npXu3h7/kWtK/684f/QBXzIdQ+U/Svpvw9/yLWlf9ecP/AKAK9nKKPs+fTt+p8lxTRVP2Vv736GlRRRXsnyR5V8fTjwRp3/YUT/0TNXz6jex/KvoP4+f8iTp3/YUT/wBEzV8/JXFifiPp8mv7H5/5Dg3zL8p79qsK3+y35VCPvL+NWF61xyPoqKfccr/vPut930qdX/2W/Kol/wBZ/wABqwnSsZHo0E+5JDJ+++433fT3rQjk/wBh/wAqpQf67/gP9a0I65qh1xTs9R6S/v8A7j/d9PerscvH+rf8qrJ/x8f8B/rVyPpXLUsc0k9dR0Ev76T93J2/h9qupN/0zk/75qtB/rpPw/lV1K5ajRySTtuLBN+/l/dSdv4farsc/wD0yl/75qvB/r5fw/lV2PpXNUaOOSdhLWfmT91L/rD/AA1fjuP+mUv/AHzVa06yf9dDV6Ouao1fY5JJ2G2lxxL+5l/1h/hq/Hcc/wCpm/75qtadJf8Aroa0I+tc1Vq+xxyvYjsbn/Rk/czd/wCH3q+lz/0wm/74qtYf8eqfj/OtGOuaq1zPQ45bENhdf6In7ibv/B71eW6/6YT/APfFQ6f/AMekf4/zq8tc82uZ6HDV2KFndf6In7ifv/B7064uv9Hk/cT/AHT/AAe1TWf/AB6J+P8AOluP+PeX/dP8q64Nc2x5tQqw3X7iP9xP90fwe1NuLr9xJ+4n+6f4ParMP/HvH/uj+VNuP9RJ/uH+Vd9Jq+xw1CrDdfuI/wBxP90fwe1Q3tz/AKK/7ibt/B71dh/494/90fyqG9/49X/D+dehSaucNTcYbn/phN/3xVa9uf8ARX/cTdv4PetE1Vvf+PV/w/nXoUbXRwz3GG5/6YTf98VWurn/AFP7mb/WD+GtA1Vuv+WP/XVa9KkefO1xpuf+mM3/AHxVa5uf9V+5m/1g/hrQNVbr/ll/10FehSOGTVxpuf8ApjN/3zVa5uP9V+5m/wBYP4avmq1z/wAsv+ugr0aZwyavsNNx/wBMZv8Avmqdzcf6TD+5l/i/h9q0jVO5/wCPmH/gX8q9CkcTavsU3uOf9VL/AN81QuJ/9Ii/dS9/4fatR6oXH/HxF+P8q9enewUmr7FWSfg/upf++aovN/pX+rk+5/d9605OlUX/AOPr/gH9a7I3Oyk12Kkk3/TOT/vmqUkv+kf6t/uenvWlJVGX/j5/4B/WuqKZ2UmuxUkl4+4/5VSMn7yT5H7dq0ZOlUj/AK2T8K6op3Wp202rFSST/Yf8qpu/7x/lbt2q/JVJ/wDWP+FdUU7rU7aditI/+y35VTZvnb5W7dqvSdKpt99vwrqindanZTKsjf7J/KqjNyeD1q5JVRup+tdUU7rU7aY3PsaZnrwetSUz1+tKqndamqGk+xpgPHQ1IaZ2rzqqd9y0MJ9jTc8dDTzTO1efVTvuUhpPsa+3vD3/ACLWlf8AXnD/AOgCviI19u+Hv+Ra0r/rzh/9AFfK55vD5/ob0zSooorwTUpaveSWGlXF1EFLxrkFxlV5wWb2HU9OBXN3Hi66jsbWWE2MhMs6zSgna4jcBVjXd/rJFIZVycj1BzXY0UAeWfHaZ08KaYghcr/aKnzMrtB8qXjrnP4Y4614P57/AN39a95+PJx4L04/9RRP/RM1fP8Avrz8VG8z6/IqrjhWk+r/AEJjO+9fl9e9L57/AN39arF/mWnb65+VHsKvK71JvPfzPu9vWl89/wC7+tVt/wA/4U7fRyoFXl3JvPfzPu9vWl89/wC7+tVt/wA/4U7fRyoFXl3JhO+9vl9O9L57/wB39arB/manb6OVBGvK25MJ33t8vp3pfPf+7+tVg/zNTt9HKgjXlbcmE772+X070vnv/d/Wqwf5mp2+jlQRrytuTLO/zfL39aXz3/u/rVZX6/Wnb6HFBGvK25Ms7/N8vf1pfPf+7+tVlfr9advocUEa8rbkyTvsHy/rS+e/939arI/yinb6HFBGvKy1JknfYPl/Whp32n5e3rVdH+UUpf5T9KOVXEq8uXcnWd9o+Xt60NO+0/L29agD/KPpQX+U/SjlQ/bytuTrO+0fL29aHnfYfl/WoA/yj6Ujv8po5VcTry5dyz57/wB39aR532H5f1qHfTXf5TQoocq8rPUs+e/939aR532H5f1qHfTXf5TQooJV5WepZ89/7v60jTv8vy9/Wod9NZ+n1oUUEq8rblnz3/u/rSNO/wAvy9/Wod9NZ+n1oUUEq8rblnz3/u/rSGd96/L696h300v8y0KKCVeVtyz57/3f1pDO+9fl9e9Q76aX+ZaFFBKvK25Z89/7v60nnv5n3e3rUO+m7/n/AAo5UDry7lnz3/u/rSee/mfd7etQ76bv+f8ACjlQOvLuWfPf+7+tJ57+Z93t61Dvpu/5/wAKOVA68u5Z89/7v60gnfe3y+neod9ND/M1HKgdeV1qWfPf+7+tIJ33t8vp3qHfTQ/zNRyoHXldalnz3/u/rSLO/wA3y9/Wod9NV+v1o5UDryutSz57/wB39aRZ3+b5e/rUO+mq/X60cqB15XWpZ89/7v60iTvsHy/rUO+mo/yijlQe3lfcs+e/939aRJ32D5f1qHfTUf5RRyoPbyvuWfPf+7+tIk77B8v61DvpqP8AKKOVB7eV9yyZ3wfl/Wvrnw9/yLWlf9ecP/oAr4/L8GvsDw9/yLWlf9ecP/oArswitc+c4hqOfs7vv+hLq95JYaVcXUQUvGuQXGVXnBZvYdT04Fc3ceLrqOxtZYTYyEyzrNKCdriNwFWNd3+skUhlXJyPUHNdjRXYfNnk/wAfZHHg/TYxC5U6kjeZkbQfKl465z+GOOteAIW/u/rX0H8fP+RJ07/sKJ/6Jmr5+SuLE/EfT5Mv3L9f8hwLbl+X171YUt/d/WoR95fxqwvWuOR9FRTvuOUt5n3f4fWp1Zv7n61Ev+s/4DVhOlYyPRoJ9ySFn877n8Pr71oRs/8Azz/8eqlB/rv+A/1rQjrmqM64p2eo9Gk8/wD1f8P973q7G8mP9V/49VZP+Pj/AID/AFq5H0rlqM5pJ66joHk86T916fxD0q6jy/8APH/x4VWg/wBdJ+H8qupXLUfkckk+4sDy+fL+59P4h6Vdjkl/54/+Piq8H+vl/D+VXY+lc1R+RxyWm4lrJNmT9x/y0P8AGKvxyTf88P8Ax8VWtOsn/XQ1ejrmqNX2OSS03G2kk2Jf9H/5aH+MVfjlnz/x7/8Aj4qtadJf+uhrQj61zVWr7HHJaEdjLP8AZk/0f1/jHrV9JZ/+fb/x8VWsP+PVPx/nWjHXNVa5nocctiGwluPsif6N6/8ALQetXlluP+fX/wAiCodP/wCPSP8AH+dXlrnm1zPQ4auxQs5bj7In+jev/LQetOuJbj7PJ/o38J/5aD0qaz/49E/H+dLcf8e8v+6f5V1wa5tjzahVhluPIj/0b+Ef8tB6U24luPIk/wBG/hP/AC0HpVmH/j3j/wB0fyptx/qJP9w/yrvpNX2OGoVYZbjyI/8ARv4R/wAtB6VDeyz/AGV/9G9P4x61dh/494/90fyqG9/49X/D+dehSfvbHDU3GGWf/n2/8fFVr2Wf7K/+jen8Y9a0TVW9/wCPV/w/nXoUXqjhm9Rhln/59v8Ax8VWupZ/3P8Ao/8Ay0H8YrQNVbr/AJY/9dVr0qR583qNMs//AD7/APj4qtcyz/uv9H/5aD+MVoGqt1/yy/66CvQpHDJ67DTLP/z7/wDj4qtcyT/uv9H/AOWg/jFXzVa5/wCWX/XQV6NI4ZNX2GmSf/n3/wDHxVO5km+0w/6P/e/jHpWkap3P/HzD/wAC/lXoUjibV9im8k2f9R/4+KoXEk32iL9x6/xj0rUeqFx/x8Rfj/KvXprQKTV9irJJLg/uf/HhVF3l+1f6n+D+8PWtOTpVF/8Aj6/4B/WuyK8zspNdipI8v/PL/wAeFUpHk+0f6r+D+971pSVRl/4+f+Af1rqivM7KTXYqSPJj/Vf+PVSLP5kn7v0/irRk6VSP+tk/CuqKd1qdtN6bFSRn/wCef61Tdn8x/k9O9X5KpP8A6x/wrqindanbTZWkZv7n61TYtvb5fTvV6TpVNvvt+FdUU7rU7KbKshb+7+tVGJyeO/rVySqjdT9a6op3Wp20xuT6UzJ5471JTPX60qqd1qaoaSfSmAnHSpDTO1edVTvuWhhJ9KbzjpTzTO1edVTvuUhpz6V9veHv+Ra0r/rzh/8AQBXxEa+3fD3/ACLWlf8AXnD/AOgCvls83h8/0N6ZpUUUV4JqFFFFAHlfx8JHgnTsf9BRP/RM1fPO5vT9a+iPjyM+C9O/7Cif+iZq+ftlcWIdpn0+Twbw7a7v9CuWbcPl/Wl3N6frUpT5lp2ysLo9RU5a6lfc2/7vb1pdzen61Ls+f8KdsougVOXcr7m3/d7etLub0/Wpdnz/AIU7ZRdAqcu5XDNuPy/rS7m9P1qUJ8zU7ZRdAqcu5XDNuPy/rS7m9P1qUJ8zU7ZRdAqcu5XDNuPy/rS7m9P1qUJ8zU7ZRdAqcu5XVm5+Xv60u5vT9alVOv1p2yi6BU5W3K6s3Py9/Wl3N6frUqp1+tO2UXQKnK25XVm2j5f1pdzen61KifKKdsougjTlbcrqzbR8v60pZsHj9alRPlFKU+U/Si6uCpy5dyEM2Bx+tBZsHj9amCfKPpQU+U/Si6D2crbkIZsDj9aRmbafl/Wpwnyj6Ujp8pourg6cuXci3N6frSMzbT8v61Y2U10+U0JoJU5W3Itzen60jM20/L+tWNlNdPlNCaCVOVtyLc3p+tIzNx8vf1qxsprJ0+tCaB05W3Itzen60jM3Hy9/WrGymsnT60JoHTlbci3N6frSFm3D5f1qxsppT5lougdOXci3N6frSFm3D5f1qxsppT5lougdOXci3N6frSbm3/d7etWNlN2fP+FF0Dpy7kW5vT9aTc2/7vb1qxspuz5/wougdOXci3N6frSbm3/d7etWNlN2fP8AhRdA6cu5Fub0/WkDNuPy/rVjZTQnzNRdA6ctNSLc3p+tIGbcfl/WrGymhPmai6B05aakW5vT9aRWbn5e/rVjZTVTr9aLoHTlfci3N6frSKzc/L39asbKaqdfrRdA6cr7kW5vT9aRWbaPl/WrGymonyii6D2cr7kW5vT9aRWbaPl/WrGymonyii6D2cr7kW5vT9aRWbaPl/WrGymonyii6D2cr7kRZsdP1r7I8Pf8i1pX/XnD/wCgCvj0pxX2F4e/5FrSv+vOH/0AV1YZ3ueFncXH2d/P9DSooorqPBPKvj7/AMiRp2B/zFE/9EzV8+oW/u/rX0H8fP8AkSdO/wCwon/omavn5K4sT8R9Pky/c/P/ACHAtuX5fXvVhS3939ahH3l/GrC9a45H0VFO+45S3mfd/h9anVm/ufrUS/6z/gNWE6VjI9Ggn3JIWfzvufw+vvWhGz/88/8Ax6qUH+u/4D/WtCOuaozrinZ6j0aTz/8AV/w/3versbyY/wBV/wCPVWT/AI+P+A/1q5H0rlqM5pJ66joHk86T916fxD0q6jy/88f/AB4VWg/10n4fyq6lctR+RyST7iwPL58v7n0/iHpV2OSX/nj/AOPiq8H+vl/D+VXY+lc1R+RxyWm4lrJNmT9x/wAtD/GKvxyTf88P/HxVa06yf9dDV6Ouao1fY5JLTcbaSTYl/wBH/wCWh/jFX45Z8/8AHv8A+Piq1p0l/wCuhrQj61zVWr7HHJaEdjLP9mT/AEf1/jHrV9JZ/wDn2/8AHxVaw/49U/H+daMdc1Vrmehxy2IbCW4+yJ/o3r/y0HrV5Zbj/n1/8iCodP8A+PSP8f51eWuebXM9Dhq7FCzluPsif6N6/wDLQetOuJbj7PJ/o38J/wCWg9Kms/8Aj0T8f50tx/x7y/7p/lXXBrm2PNqFWGW48iP/AEb+Ef8ALQelNuJbjyJP9G/hP/LQelWYf+PeP/dH8qbcf6iT/cP8q76TV9jhqFWGW48iP/Rv4R/y0HpUN7LP9lf/AEb0/jHrV2H/AI94/wDdH8qhvf8Aj1f8P516FJ+9scNTcYZZ/wDn2/8AHxVa9ln+yv8A6N6fxj1rRNVb3/j1f8P516FF6o4ZvUYZZ/8An2/8fFVrqWf9z/o//LQfxitA1Vuv+WP/AF1WvSpHnzeo0yz/APPv/wCPiq1zLP8Auv8AR/8AloP4xWgaq3X/ACy/66CvQpHDJ67DTLP/AM+//j4qtcyT/uv9H/5aD+MVfNVrn/ll/wBdBXo0jhk1fYaZJ/8An3/8fFU7mSb7TD/o/wDe/jHpWkap3P8Ax8w/8C/lXoUjibV9im8k2f8AUf8Aj4qhcSTfaIv3Hr/GPStR6oXH/HxF+P8AKvXprQKTV9irJJLg/uf/AB4VRd5ftX+p/g/vD1rTk6VRf/j6/wCAf1rsivM7KTXYqSPL/wA8v/HhVKR5PtH+q/g/ve9aUlUZf+Pn/gH9a6orzOyk12KkjyY/1X/j1Uiz+ZJ+79P4q0ZOlUj/AK2T8K6op3Wp203psVJGf/nn+tU3Z/Mf5PTvV+SqT/6x/wAK6op3Wp202VpGb+5+tU2Lb2+X071ek6VTb77fhXVFO61OymyrIW/u/rVRicnjv61ckqo3U/WuqKd1qdtMbk+lMyeeO9SUz1+tKqndamqGkn0pgJx0qQ0ztXnVU77loYSfSm846U80ztXnVU77lIac+lfb3h7/AJFrSv8Arzh/9AFfERr7d8Pf8i1pX/XnD/6AK+WzzeHz/Q3pmlRRRXgmpS1eS6i0q4ezDeeF+Uou5gM8lV5yQMkDByR0Nc/Pqevf2XZm1tb2aQX2yeU26pIYRMoXKvt+9GwJYD5cN908r1tFAHlfx385vCmmKqR+V/aKksXO7d5UvGMdMZ5z+FeDeXJ6LX0J8aojL4V0xB/0E1/9EzV4r9gavLxtVRqWfY+54cw3tcG5X+0/yRhmOTevC96d5cnota7WLeYg+tP+wNXJ7dHurA6vUw/Lk8zovSneXJ6LWv8AYW87H+zT/sDUe3QRwPmYflyeZ0XpTvLk9FrX+wt52P8AZp/2BqPboI4HzMMRyb24XtTvLk9FrXWxbzHH0p/2BqHXQRwOm5hiOTe3C9qd5cnota62LeY4+lP+wNQ66COB03MMRyb24XtTvLk9FrXWxbzHH0p/2BqHXQRwOm5hrHJzwvWneXJ6LWvHYsd3+8af9gah10EcDpuYaxyc8L1p3lyei1rx2LHd/vGn/YGoddBHA6bmGkcmwcLTvLk9FrXisWMYNP8AsDUOurhHA6LUw0jk2DhaVo5Np4XpWxFYsYwaV7BtjfSj26uJYH3dzGWOTaOF6UNHJtPC9K2UsG2L9KHsG2N9KPbq4/qPu7mMscm0cL0pHjk2Hha2ksG2L9KSWxYRk0e3VxPA+7uZHlyei0145Nh4Wtz7A1MlsWEZNCrq45YH3XqZHlyei0145Nh4Wtz7A1MlsWEZNCrq4SwPuvUyPLk9FprRyccL1rc+wNTJLFht/wB4UKuglgdNzI8uT0WmtHJxwvWtz7A1MksWG3/eFCroJYHTcyPLk9Fppjk3rwvetz7A1MaxbzEH1oVdBLA6bmR5cnotNMcm9eF71ufYGpjWLeYg+tCroJYHTcyPLk9FpvlyeZ0XpW59gamfYW87H+zQq6B4HbUyPLk9FpvlyeZ0XpW59gamfYW87H+zQq6B4HbUyPLk9FpvlyeZ0XpW59gamfYW87H+zQq6B4HbUyPLk9Fpojk3twvatz7A1MWxbzHH0o9ugeB1WpkeXJ6LTRHJvbhe1bn2BqYti3mOPpR7dA8DqtTI8uT0WmrHJzwvWtz7A1MjsWO7/eNHt0DwOq1Mjy5PRaascnPC9a3PsDUyOxY7v940e3QPA6rUyPLk9FpqRybBwtbn2BqZFYsYwaPbqwfUddzI8uT0WmpHJsHC1ufYGpkVixjBo9urB9R13Mjy5PRaakcmwcLW59gamRWLGMGj26sH1HXcyDHJg8LX134e/wCRa0r/AK84f/QBXy6bBsGvqLw9/wAi1pX/AF5w/wDoArvwNRT5reR8pxPQ9l7Lz5v0JdXkuotKuHsw3nhflKLuYDPJVeckDJAwckdDXPz6nr39l2ZtbW9mkF9snlNuqSGETKFyr7fvRsCWA+XDfdPK9bRXoHyh5P8AH0zf8IfpoCJ5X9pIS5c7t3lS8Yx0xnnP4V4Am72r6D+Pn/Ik6d/2FE/9EzV8/JXFifiPp8mX7n5/5Dhu3L071YXf/s1CPvL+NWF61xyPoqKHLv8AM/h+7U67/wDZqJf9Z/wGrCdKxkz0aESSHzPO/h+7/WtCPzf9j9apQf67/gP9a0I65qjOuMdHqPTzfP8A4Pu+/rV2PzsfwfrVZP8Aj4/4D/WrkfSuWozmlHcdB53nSf6vt6+lXU8//pn+tVoP9dJ+H8qupXLUZySWgsHn+fL/AKvt6+lXY/tH/TL9arwf6+X8P5Vdj6VzVGccloJa/aMyf6r/AFh9avx/af8Apl+tVrTrJ/10NXo65qj1OSS0G2n2nEv+q/1h9avx/as/8sf1qtadJf8Aroa0I+tc1WWpxyWhHY/avsyf6nv6+tX0+1/9Mf1qtYf8eqfj/OtGOuarL3noccloQ2H2v7ImPJ7+vrV5ftn/AEw/WodP/wCPSP8AH+dXlrnnL3nocNXYoWf2v7In+o7+vrTrj7X9nk/1H3T6+lTWf/Hon4/zpbj/AI95f90/yrrhL3jzahVh+1+RH/qPuj19Kbcfa/Ik/wBR90+vpVmH/j3j/wB0fyptx/qJP9w/yrvpPU4ahVh+1+RH/qPuj19Khvftf2V/9T29fWrsP/HvH/uj+VQ3v/Hq/wCH869Ck9ThqbjD9r/6Y/rVa9+1fZX/ANT29fWtE1Vvf+PV/wAP516FF6o4ZvUYftX/AEx/Wq119q/c/wCp/wBYPWtA1Vuv+WP/AF1WvSpHnzeo0/av+mP61WuftX7r/U/6wetaBqrdf8sv+ugr0KRwyeo0/av+mP61WuftP7r/AFP+sHrV81Wuf+WX/XQV6NM4ZPUaftP/AEx/Wqdz9p+0w/6r+L19K0jVO5/4+Yf+Bfyr0KRxOWuxTf7Tn/ll+tULj7R9oi/1Xf19K1Hqhcf8fEX4/wAq9emtApS12Ksn2jB/1X61Rfz/ALV/yz+57+tacnSqL/8AH1/wD+tdkUdlKXkVJPP/AOmf61Sk877R/B9z39a0pKoy/wDHz/wD+tdUUdlJlSTzcfwfrVI+b5kn3O3rWjJ0qkf9bJ+FdUVqjtpvQqSeZ/sfrVN/M8x/u9qvyVSf/WP+FdUVqjtpsrSeZ/s1Tbfvb7var0nSqbffb8K6ox1Wp2U2VZN3tVRt2T061ckqo3U/WuqMdVqdtMbz7Uznnp1qSmev1pVY6rU1Q059qYM47VIaZ2rzqsdS0MOfam847U80ztXn1VqUhpz7V9veHv8AkWtK/wCvOH/0AV8RGvt3w9/yLWlf9ecP/oAr5XPFrD5/ob0zSooorwTUKKKKAOB+LMTTaDpaIAT/AGiOv/XGavLf7PuP7qfrXsfj2H7RZ6THjOdQz/5Amrlf7J/2TXyOeYr2WJUfJfqfXZHi1Rwzjfq/0OCfT7jzo/lTv61J/Z9x/dT9a7OTSv8ASYRtPO7+VTf2T/smvIeYLQ9dZitdTgv7PuPtGNqfc9/WpP7PuP7qfrXZ/wBlf6Zjaf8AV/1qb+yf9k0PMECzFdzgv7PuPtGNqfc9/WpP7PuP7qfrXZ/2V/pmNp/1f9am/sn/AGTQ8wQLMV3OCTT7jzpPlTt61J/Z9x/dT9a7OPSv9JmG08bf5VN/ZP8Asmh5grgsxXc4JNPuPOk+VO3rUn9n3H91P1rs49K/0mYbTxt/lU39k/7JoeYK4LMV3OCTT7jzpPlTt61J/Z9x/dT9a7OPSv8ASZhtPG3+VTf2T/smh5grgsxXc4KLT7j5/lT759ak/s+4/up+tdnb6VnzflP+sNTf2T/smiWYK4RzFW3OCi0+4+f5U++fWpP7PuP7qfrXZ2+lZ835T/rDU39k/wCyaJZgrhHMVbc4KDT7gwr8qfrUn9n3H91P1rs7XSs2yHae/wDOpv7J/wBk0SzBXYRzFWWpwUGn3BhX5U/WnSafceW3yp0PrXa2ulZtkO09/wCdPl0rEL/Kfumm8euYFmK5dzh49PuPLX5U6D1ok0+48tvlTofWu4i0rMKfKfuiiXSsQv8AKfuml/aCuH9oq25w8en3Hlr8qdB602fT7gQt8qfrXdRaVmFPlP3RTLrSsWznae386ax65geYrl3OM/s+4/up+tRz6fcCFvlT9a73+yf9k1DdaVi2c7T2/nSjmCuglmKs9TjP7PuP7qfrUc+n3Ahb5U/Wu9/sn/ZNQ3WlYtnO09v50RzBXQSzFWepxn9n3H91P1qOXT7gbPlT749a73+yf9k1DcaVjyvlP+sFEcwVwlmKtucZ/Z9x/dT9ajl0+4Gz5U++PWu9/sn/AGTUNxpWPK+U/wCsFEcwVwlmKtucZ/Z9x/dT9ajfT7jzo/lTv613v9k/7JqGTSv9JhG087v5ULMFcHmK7nGf2fcf3U/Wo30+486P5U7+td7/AGT/ALJqGTSv9JhG087v5ULMFcHmK7nGf2fcf3U/Wo/7PuPtGNqfc9/Wu9/sn/ZNQ/2V/pmNp/1f9aFmCB5iu5xn9n3H91P1qP8As+4+0Y2p9z39a73+yf8AZNQ/2V/pmNp/1f8AWhZggeYrucZ/Z9x/dT9aj/s+4+0Y2p9z39a73+yf9k1D/ZX+mY2n/V/1oWYIHmK7nGf2fcf3U/Wo00+486T5U7etd7/ZP+yahj0r/SZhtPG3+VCzBWYPMVpqcZ/Z9x/dT9ajTT7jzpPlTt613v8AZP8AsmoY9K/0mYbTxt/lQswVmDzFaanGf2fcf3U/Wo4tPuPn+VPvn1rvf7J/2TUNvpWfN+U/6w0LMFYHmKutTjP7PuP7qfrUcWn3Hz/Kn3z613v9k/7JqG30rPm/Kf8AWGhZgrA8xV1qcZ/Z9x/dT9ajg0+4MK/Kn613v9k/7JqG10rNsh2nv/Oj+0FYP7RV9zjP7PuP7qfrUcGn3BhX5U/Wu9/sn/ZNQ2ulZtkO09/50f2grB/aKvucZ/Z9x/dT9ajg0+4MK/Kn613v9k/7JqG10rNsh2nv/Oj+0FYP7RV9zi20+42n5U6e9e/+Hv8AkWtK/wCvOH/0AV5y+lfI3ynpXo3h7/kWtK/684f/AEAV9Fw/iPbe08rfqfO8QYn23s/K/wChpUUUV9GfOHlXx9z/AMIRp2P+gon/AKJmr59Td7V9B/Hz/kSdO/7Cif8Aomavn5K4sT8R9Pky/c/P/IcN25enerC7/wDZqEfeX8asL1rjkfRUUOXf5n8P3anXf/s1Ev8ArP8AgNWE6VjJno0IkkPmed/D93+taEfm/wCx+tUoP9d/wH+taEdc1RnXGOj1Hp5vn/wfd9/WrsfnY/g/Wqyf8fH/AAH+tXI+lctRnNKO46DzvOk/1fb19Kup5/8A0z/Wq0H+uk/D+VXUrlqM5JLQWDz/AD5f9X29fSrsf2j/AKZfrVeD/Xy/h/KrsfSuaozjktBLX7RmT/Vf6w+tX4/tP/TL9arWnWT/AK6Gr0dc1R6nJJaDbT7TiX/Vf6w+tX4/tWf+WP61WtOkv/XQ1oR9a5qstTjktCOx+1fZk/1Pf19avp9r/wCmP61WsP8Aj1T8f51ox1zVZe89DjktCGw+1/ZEx5Pf19avL9s/6YfrUOn/APHpH+P86vLXPOXvPQ4auxQs/tf2RP8AUd/X1p1x9r+zyf6j7p9fSprP/j0T8f50tx/x7y/7p/lXXCXvHm1CrD9r8iP/AFH3R6+lNuPtfkSf6j7p9fSrMP8Ax7x/7o/lTbj/AFEn+4f5V30nqcNQqw/a/Ij/ANR90evpUN79r+yv/qe3r61dh/494/8AdH8qhvf+PV/w/nXoUnqcNTcYftf/AEx/Wq179q+yv/qe3r61omqt7/x6v+H869Ci9UcM3qMP2r/pj+tVrr7V+5/1P+sHrWgaq3X/ACx/66rXpUjz5vUaftX/AEx/Wq1z9q/df6n/AFg9a0DVW6/5Zf8AXQV6FI4ZPUaftX/TH9arXP2n91/qf9YPWr5qtc/8sv8AroK9GmcMnqNP2n/pj+tU7n7T9ph/1X8Xr6VpGqdz/wAfMP8AwL+VehSOJy12Kb/ac/8ALL9aoXH2j7RF/qu/r6VqPVC4/wCPiL8f5V69NaBSlrsVZPtGD/qv1qi/n/av+Wf3Pf1rTk6VRf8A4+v+Af1rsijspS8ipJ5//TP9apSed9o/g+57+taUlUZf+Pn/AIB/WuqKOykypJ5uP4P1qkfN8yT7nb1rRk6VSP8ArZPwrqitUdtN6FSTzP8AY/Wqb+Z5j/d7Vfkqk/8ArH/CuqK1R202VpPM/wBmqbb97fd7Vek6VTb77fhXVGOq1OymyrJu9qqNuyenWrklVG6n611RjqtTtpjefamc89OtSUz1+tKrHVamqGnPtTBnHapDTO1edVjqWhhz7U3nHanmmdq8+qtSkNOfavt7w9/yLWlf9ecP/oAr4iNfbvh7/kWtK/684f8A0AV8rni1h8/0N6ZpUUUV4JqUtXkuotKuHsw3nhflKLuYDPJVeckDJAwckdDWC2pakbCyVZNTW6a8ePf9gO1og+VaUeXwDHt+7tyxI4wdvV0UAc14sjmkn0bbMiQ/bDlTHlt3kS85zjHtj8ap/ZW/5+B/37FXvFz7F0c/9Px/9ETVl/afevlc6yqpisSqkeyX5nFiM3+qT9nfzGy2rfa7f/SB/F/yzHpVj7K3/PwP+/YqlLcf6Vb8/wB7+VT/AGn3ryXkFXQw/wBYvMb9lb7f/wAfA/1X/PMetWPsrf8APwP+/Yql9o/0/r/yy/rU/wBp96HkFUP9YvMb9lb7f/x8D/Vf88x61Y+yt/z8D/v2KpfaP9P6/wDLL+tT/afeh5BVD/WLzGxWrfa7j/SB/D/yzHpVj7K3/PwP+/YqlFcf6Vcc/wB3+VT/AGn3oeQVWH+sVuo2K1b7Xcf6QP4f+WY9KsfZW/5+B/37FUorj/Srjn+7/Kp/tPvQ8gqsP9YrdRsVq32u4/0gfw/8sx6VY+yt/wA/A/79iqUVx/pVxz/d/lU/2n3oeQVWH+sVuo22tW/ff6QP9a3/ACzFWPsrf8/A/wC/YqlbXH+u5/5atU/2n3oeQVWw/wBYrdRttat++/0gf61v+WYqx9lb/n4H/fsVStrj/Xc/8tWqf7T70PIKrYf6xW6jbO1b7In+kDv/AMsx61Y+yt/z8D/v2KpWdx/oqc+v86n+0+9Esgqt3D/WK2lxtnat9kT/AEgd/wDlmPWpJ7Vvs8n+kD7p/wCWY9KrWdx/oqc+v86fPcfuJOf4T/Kh5BV5rh/rF0uTQWreRH/pA+6P+WY9KJ7Vvs8n+kD7p/5Zj0qGC4/cR8/wj+VE9x+4k5/hP8qP7Aq3uH+sXS5NBat5Ef8ApA+6P+WY9KjvLVvsj/6QO3/LMetNguP3EfP8I/lTLy4/0V+fT+dCyCrzXD/WLpcu/ZW/5+B/37FV7y1b7I/+kDt/yzHrTvtPvUF5cf6K/Pp/OiOQVU0w/wBYr6XLv2Vv+fgf9+xVe8tW+yP/AKQO3/LMetO+0+9QXlx/or8+n86I5BVTTD/WK+ly79lb/n4H/fsVXubVv3P+kD/Wr/yzFO+0+9QXNx/qef8AlqtCyCqmH+sV+pd+yt/z8D/v2Kr3Nq37n/SB/rV/5ZinfafeoLm4/wBTz/y1WhZBVTD/AFiv1Lv2Vv8An4H/AH7FV5bVvtdv/pA/i/5Zj0p32n3qCW4/0q35/vfyoWQVUH+sV+pd+yt/z8D/AL9iq8tq32u3/wBIH8X/ACzHpTvtPvUEtx/pVvz/AHv5ULIKqD/WK/Uu/ZW/5+B/37FV/srfb/8Aj4H+q/55j1p32n3qD7R/p/X/AJZf1oWQVQ/1i8y79lb/AJ+B/wB+xVf7K32//j4H+q/55j1p32n3qD7R/p/X/ll/WhZBVD/WLzLv2Vv+fgf9+xVf7K32/wD4+B/qv+eY9ad9p96g+0f6f1/5Zf1oWQVQ/wBYvMu/ZW/5+B/37FV4rVvtdx/pA/h/5Zj0p32n3qCK4/0q45/u/wAqFkFXUP8AWLzLv2Vv+fgf9+xVeK1b7Xcf6QP4f+WY9Kd9p96giuP9KuOf7v8AKhZBV1D/AFi8y79lb/n4H/fsVXtrVv33+kD/AFrf8sxTvtPvUFtcf67n/lq1CyCrYP8AWLzLv2Vv+fgf9+xVe2tW/ff6QP8AWt/yzFO+0+9QW1x/ruf+WrULIKtg/wBYvMu/ZW/5+B/37FV7O1b7In+kDv8A8sx6077T71BZ3H+ipz6/zo/sCrawf6xdbl37K3/PwP8Av2Kr2dq32RP9IHf/AJZj1p32n3qCzuP9FTn1/nR/YFW1g/1i63Lv2Vv+fgf9+xVeztW+yJ/pA7/8sx6077T71BZ3H+ipz6/zo/sCrawf6xdblt7VvLb/AEgdD/yzFdD4e/5FrSv+vOH/ANAFcw9z+7bnsa6fw9/yLWlf9ecP/oAr6DIcvng/ac3W34XOrDZj9dvrt+v/AAxLq8l1FpVw9mG88L8pRdzAZ5KrzkgZIGDkjoawW1LUjYWSrJqa3TXjx7/sB2tEHyrSjy+AY9v3duWJHGDt6uivoDqPJ/j6sv8Awh+mkSIIf7SQFdh3bvKl5znGMZ4x+NeAIG/vfpX0H8fP+RJ07/sKJ/6Jmr5+SuLE/EfT5Mv3L9f8hwDbl+b17VYUN/e/SoR95fxqwvWuOTPoqMVccobzPvfw+lTqrf3/ANKiX/Wf8BqwnSsZM9GhFEkKv533/wCH0960I1f/AJ6f+O1Sg/13/Af61oR1zVGdcYqzHosnn/6z+H+771djSTH+t/8AHarJ/wAfH/Af61cj6Vy1Gc0orUdAknnSfvfT+EelXUSX/nt/46KrQf66T8P5VdSuWozklFCwJL58v770/hHpV2OOX/nt/wCOCq8H+vl/D+VXY+lc1RnHJaCWsc2ZP3//AC0P8Aq/HHN/z3/8cFVrTrJ/10NXo65qknc5JLQbaRzYl/0j/lof4BV+OKfP/Hx/44KrWnSX/roa0I+tc1WTuccloR2MU/2ZP9I9f4B61fSKf/n5/wDHBVaw/wCPVPx/nWjHXNVk+ZnHJaENhFcfZE/0n1/5Zj1q8sVx/wA/X/kMVDp//HpH+P8AOry1zzk+ZnDV2KFnFcfZE/0n1/5Zj1p1xFcfZ5P9J/hP/LMelTWf/Hon4/zpbj/j3l/3T/KuuEnzHm1CrDFceRH/AKT/AAj/AJZj0ptxFceRJ/pP8J/5Zj0qzD/x7x/7o/lTbj/USf7h/lXfSk7nDUKsMVx5Ef8ApP8ACP8AlmPSob2Kf7K/+k+n8A9auw/8e8f+6P5VDe/8er/h/OvQpP3jhqPUYYp/+fn/AMcFVr2Kf7K/+k+n8A9a0TVW9/49X/D+dehReqOGb1GGKf8A5+f/ABwVWuop/wBz/pH/AC0H8ArQNVbr/lj/ANdVr0qTPPm9Rpin/wCfj/xwVWuYp/3X+kf8tB/AK0DVW6/5Zf8AXQV6FI4ZPUaYp/8An4/8cFVrmOf91/pH/LQfwCr5qtc/8sv+ugr0aRwyk7jTHP8A8/H/AI4Kp3Mc32mH9/8A3v4B6VpGqdz/AMfMP/Av5V6FI4nJ3KbxzZ/1/wD44KoXEc32iL9/6/wD0rUeqFx/x8Rfj/KvXprQKUncqyRy4P77/wAdFUXSX7V/rv4P7o9a05OlUX/4+v8AgH9a7Io7KUmVJEl/56/+OiqUiSfaP9b/AAf3fetKSqMv/Hz/AMA/rXVFHZSkypIkmP8AW/8AjtUir+ZJ+89P4a0ZOlUj/rZPwrqjFXR203oVJFf/AJ6fpVN1fzH+f07Vfkqk/wDrH/CuqMVdHbTZWkVv7/6VTYNvb5vTtV6TpVNvvt+FdUYq6OymyrIG/vfpVRgcnnv6Vckqo3U/WuqMVdHbTY3B9aZg8896kpnr9aVWKujVDSD60wA461IaZ2rzqsVctDCD603nHWnmmdq86rFXKQ059a+3vD3/ACLWlf8AXnD/AOgCviI19u+Hv+Ra0r/rzh/9AFfLZ4tYfP8AQ3pmlRRRXgmoUUUUAcb8Rrn7Jp2kyhtuNQxnGf8AlhNXFf28f+fgf98itj46XTWXg7TZl6jU0H/kGavBv+Ejl9/zr6zJMphi8M6kl1a/I+XzjLZ4nEKcX0S/M9bk1wm4hP2gcbv4R6VN/bx/5+B/3yK8cPiKUup54z3p3/CRy+/51664epa6I8t5HV01PW/7cP2vP2gf6vH3R61N/bx/5+B/3yK8c/4SKXzM89Mdad/wkcvv+dC4epdkDyOr3PW/7cP2vP2gf6vH3R61N/bx/wCfgf8AfIrxz/hIpfMzz0x1p3/CRy+/50Lh6l2QPI6vc9bj1wi4mP2gc7f4R6VN/bx/5+B/3yK8cHiKUOx55x3p3/CRy+/50Lh6k+iB5HV7nrceuEXEx+0Dnb/CPSpv7eP/AD8D/vkV44PEUodjzzjvTv8AhI5ff86Fw9SfRA8jq9z1uPXCLiY/aBzt/hHpU39vH/n4H/fIrxweIpQ7HnnHenf8JHL7/nQuHqT6IHkdXuetwa4R5n+kDmQn7oqb+3j/AM/A/wC+RXjieIpRu68nPWnf8JHL7/nQuHqTWyB5HVb3PW4NcI8z/SBzIT90VN/bx/5+B/3yK8cTxFKN3Xk5607/AISOX3/OhcPUmtkDyOq3uet22uEW6j7QO/8ACPWpv7eP/PwP++RXjieIpVQDn86d/wAJHL7/AJ0R4epNJ2QPI6re563ba4RbqPtA7/wj1p8uukwuPtA+6f4RXj6eIpVQDn86VvEUpUjnp60lw9S5b2G8kq817nr8WukQoPtA+6P4RRLrpMLj7QPun+EV5AviKUKBz09aG8RSlSOenrT/ANXqVr2Qv7Dq33PX4tdIhQfaB90fwimXOuE27D7QO38I9a8jXxFKFA56etI/iKVkI5/Ok+HqXLew1klXmvc9j/t4/wDPwP8AvkVDc64TbsPtA7fwj1ryT/hI5ff86a/iKVkI5/OnLh6klshRyOqnuex/28f+fgf98iobnXCbdh9oHb+EeteSf8JHL7/nTX8RSshHP50S4epJbII5HVT3PY/7eP8Az8D/AL5FQz64T5f+kDiQH7oryT/hI5ff86a/iKU7evBz1ofD1JLZAsjqp7nsf9vH/n4H/fIqGfXCfL/0gcSA/dFeSf8ACRy+/wCdNfxFKdvXg560Ph6klsgWR1U9z2P+3j/z8D/vkVDJrhNxCftA43fwj0ryT/hI5ff86afEUpdTzxnvQ+HqS6IFkdVdT2P+3j/z8D/vkVDJrhNxCftA43fwj0ryT/hI5ff86afEUpdTzxnvQ+HqS6IFkdVdT2P+3j/z8D/vkVD/AG4fteftA/1ePuj1ryT/AISOX3/Om/8ACRS+ZnnpjrQ+HqXZAsjq9z2P+3j/AM/A/wC+RUP9uH7Xn7QP9Xj7o9a8k/4SOX3/ADpv/CRS+ZnnpjrQ+HqXZAsjq9z2P+3j/wA/A/75FQ/24fteftA/1ePuj1ryT/hI5ff86b/wkUvmZ56Y60Ph6l2QLI6vc9j/ALeP/PwP++RUMeuEXEx+0Dnb/CPSvJP+Ejl9/wA6aPEUodjzzjvQ+HqWmiBZHV11PY/7eP8Az8D/AL5FQx64RcTH7QOdv8I9K8k/4SOX3/OmjxFKHY88470Ph6lpogWR1ddT2P8At4/8/A/75FQwa4R5n+kDmQn7oryT/hI5ff8AOmp4ilG7ryc9aHw9SvsgWR1bbnsf9vH/AJ+B/wB8ioYNcI8z/SBzIT90V5J/wkcvv+dNTxFKN3Xk560Ph6lfZAsjq23PY/7eP/PwP++RUNtrhFuo+0Dv/CPWvJP+Ejl9/wA6aniKVUA5/Oj/AFepXtZB/YdW257H/bx/5+B/3yKhttcIt1H2gd/4R615J/wkcvv+dNTxFKqAc/nR/q9SvayD+w6ttz2P+3j/AM/A/wC+RUNtrhFuo+0Dv/CPWvJP+Ejl9/zpqeIpVQDn86P9XqV7WQf2HVtuewvrp2N/pA6f3RXq/h7/AJFrSv8Arzh/9AFfJJ8Ry4PX86+tvD3/ACLWlf8AXnD/AOgCvns/y6GD9ny9b/hb/M97I8DLC+05utvwuaVFFFfOnvHlXx9/5EjTsH/mKJ/6Jmr59QN/e/SvoP4+f8iTp3/YUT/0TNXz8lcWJ+I+nyZfufn/AJDgG3L83r2qwob+9+lQj7y/jVhetccmfRUYq45Q3mfe/h9KnVW/v/pUS/6z/gNWE6VjJno0IokhV/O+/wDw+nvWhGr/APPT/wAdqlB/rv8AgP8AWtCOuaozrjFWY9Fk8/8A1n8P933q7GkmP9b/AOO1WT/j4/4D/WrkfSuWozmlFajoEk86T976fwj0q6iS/wDPb/x0VWg/10n4fyq6lctRnJKKFgSXz5f33p/CPSrsccv/AD2/8cFV4P8AXy/h/KrsfSuaozjktBLWObMn7/8A5aH+AVfjjm/57/8AjgqtadZP+uhq9HXNUk7nJJaDbSObEv8ApH/LQ/wCr8cU+f8Aj4/8cFVrTpL/ANdDWhH1rmqydzjktCOxin+zJ/pHr/APWr6RT/8APz/44KrWH/Hqn4/zrRjrmqyfMzjktCGwiuPsif6T6/8ALMetXliuP+fr/wAhiodP/wCPSP8AH+dXlrnnJ8zOGrsULOK4+yJ/pPr/AMsx6064iuPs8n+k/wAJ/wCWY9Kms/8Aj0T8f50tx/x7y/7p/lXXCT5jzahVhiuPIj/0n+Ef8sx6U24iuPIk/wBJ/hP/ACzHpVmH/j3j/wB0fyptx/qJP9w/yrvpSdzhqFWGK48iP/Sf4R/yzHpUN7FP9lf/AEn0/gHrV2H/AI94/wDdH8qhvf8Aj1f8P516FJ+8cNR6jDFP/wA/P/jgqtexT/ZX/wBJ9P4B61omqt7/AMer/h/OvQovVHDN6jDFP/z8/wDjgqtdRT/uf9I/5aD+AVoGqt1/yx/66rXpUmefN6jTFP8A8/H/AI4KrXMU/wC6/wBI/wCWg/gFaBqrdf8ALL/roK9CkcMnqNMU/wDz8f8Ajgqtcxz/ALr/AEj/AJaD+AVfNVrn/ll/10FejSOGUncaY5/+fj/xwVTuY5vtMP7/APvfwD0rSNU7n/j5h/4F/KvQpHE5O5TeObP+v/8AHBVC4jm+0Rfv/X+Aelaj1QuP+PiL8f5V69NaBSk7lWSOXB/ff+OiqLpL9q/138H90etacnSqL/8AH1/wD+tdkUdlKTKkiS/89f8Ax0VSkST7R/rf4P7vvWlJVGX/AI+f+Af1rqijspSZUkSTH+t/8dqkVfzJP3np/DWjJ0qkf9bJ+FdUYq6O2m9CpIr/APPT9Kpur+Y/z+nar8lUn/1j/hXVGKujtpsrSK39/wDSqbBt7fN6dqvSdKpt99vwrqjFXR2U2VZA3979KqMDk89/SrklVG6n611Riro7abG4PrTMHnnvUlM9frSqxV0aoaQfWmAHHWpDTO1edVirloYQfWm846080ztXnVYq5SGnPrX294e/5FrSv+vOH/0AV8RGvt3w9/yLWlf9ecP/AKAK+Wzxaw+f6G9M0qKKK8E1KWrx3UulXEdmWE5X5QjbWIzyFbjBIyAcjBPUVzs1pr0mnRRxx3yOhuVjzdruVmYGB3YP86IpKsCWJI6Nwa6+igDj/iD4Lk8c6dY6a16LW0iuhPKyx7nyI3UYJOMZYDGPxrhf+GddO/6GK7/78L/jXeQ+Gr+KK9h22RScRbv3jf6QySlyZBs/5aKxVuuMD7wPHQaPZSafpkdtLsDKzsEjOUjDOWCLwOFBCjgcDoOldNLGV6MeWnNpCcU9zyL/AIZ107/oYrv/AL8L/jS/8M66d/0MV3/34X/GvaqK0/tLF/8APx/eLkj2PFP+GddOz/yMV3/34X/Gl/4Z107/AKGK7/78L/jXtVFH9pYv/n4/vDkj2PFP+GddOz/yMV3/AN+F/wAaX/hnXTv+hiu/+/C/417VRR/aWL/5+P7w5I9jxT/hnXTv+hiu/wDvwv8AjS/8M66d/wBDFd/9+F/xr2qqWr2cl/pVxaxFQ8i4Ac4VuclW9j0PsaP7Sxf/AD8f3hyR7Hjy/s7Wf2h1OvXQiCKVfyVyWycjGewC/n7VJ/wzrp3/AEMV3/34X/Gu9l8MXM+n2kL2lis9rdSTW8q3D4tlZ942rsw23O0KcDCjpnA62j+0sX/z8f3hyR7Hin/DOunf9DFd/wDfhf8AGl/4Z107/oYrv/vwv+Ne1UUf2li/+fj+8OSPY8U/4Z107/oYrv8A78L/AI0v/DOunf8AQxXf/fhf8a9qoo/tLF/8/H94ckex4p/wzrp3/QxXf/fhf8aX/hnXTv8AoYrv/vwv+Ne1UUf2li/+fj+8OSPY8U/4Z104f8zFd/8Afhf8aX/hnXTv+hiu/wDvwv8AjXtVZOv6fcaha2y2sFrLNDdwzhrhymwJIrEqQrHJAI7detH9pYv/AJ+P7w5I9jyeH9nazaMmXXrpG3sABCp+Xcdp69xg/jUn/DOunf8AQxXf/fhf8a9Li0a8PiFNTlW2UlxI7o5aQDydhgHyjMe7585HP8Pet+j+0sX/AM/H94ckex4r/wAM66d/0MV3/wB+F/xo/wCGddO/6GK7/wC/C/417VRR/aWL/wCfj+8OSPY8V/4Z107/AKGK7/78L/jSf8M66cf+Ziu/+/C/417XRR/aWL/5+P7w5I9jxX/hnXTv+hiu/wDvwv8AjSf8M66cf+Ziu/8Avwv+Ne10Uf2li/8An4/vDkj2PFf+GddO/wChiu/+/C/41Hcfs7Wa28rQ69dSShCURoVAZscAnPHNe3Vj3Nhd/wDCUWupW9rZmJbZ4JpHlKytueM9AhyFCNjJ6t260f2li/8An4/vDkj2PLf+GddO/wChiu/+/C/40n/DOunf9DFd/wDfhf8AGvUPDuitoovEWCC2gmkDpBDM0oBx8zFmAPPHHQY9626P7Sxf/Px/eHJHseK/8M66d/0MV3/34X/Gk/4Z107/AKGK7/78L/jXtdFH9pYv/n4/vDkj2PFf+GddO/6GK7/78L/jSf8ADOunf9DFd/8Afhf8a9roo/tLF/8APx/eHJHseK/8M66d/wBDFd/9+F/xpP8AhnXTv+hiu/8Avwv+Ne10Uf2li/8An4/vDkj2PFf+GddO/wChiu/+/C/41Gf2drMXKKNeujEUYs/krkNkYGM98t+XvXtcyloZFCI5KkBHOFb2PB4/A1ycfhy/OkR2s1jpeVuDcSWyyt9nkLI6ldvl/Kq5Rhwclcnn5qP7Sxf/AD8f3hyR7HCf8M66d/0MV3/34X/Gk/4Z107P/IxXf/fhf8a9ksoHtbC3t5JTNJFEqNI3VyAASfrU9H9pYv8A5+P7w5I9jxX/AIZ107/oYrv/AL8L/jSf8M66dn/kYrv/AL8L/jXtdFH9pYv/AJ+P7w5I9jxX/hnXTv8AoYrv/vwv+NJ/wzrp3/QxXf8A34X/ABr2uij+0sX/AM/H94ckex4r/wAM66d/0MV3/wB+F/xpP+GddO/6GK7/AO/C/wCNe10Uf2li/wDn4/vDkj2PFf8AhnXTv+hiu/8Avwv+NRx/s7WZeYSa9dKofEZEKncu0cnnjnI/CvW9fsp9Q0mS2traynkcgYvPuKM8sBtbLAdMjGfyrO/sCVtW0q+FrbRTWaCJp/tTyP5YDgIFKANkEZY4IJPXHJ/aWL/5+P7w5I9jzn/hnXTv+hiu/wDvwv8AjSf8M66d/wBDFd/9+F/xr2uij+0sX/z8f3hyR7Hiv/DOunf9DFd/9+F/xpP+GddOH/MxXf8A34X/ABr2uij+0sX/AM/H94ckex4r/wAM66d/0MV3/wB+F/xpP+GddOH/ADMV3/34X/Gva6KP7Sxf/Px/eHJHseK/8M66d/0MV3/34X/Gk/4Z104f8zFd/wDfhf8AGva6KP7Sxf8Az8f3hyR7Hiv/AAzrp3/QxXf/AH4X/GvXdHtpbLRLC1n2+dDbRxybem4KAcfiKoa9pV1qFzYz2kFo0lrKsoeZ9p4YHaD5bMAcdQw9wwypTTNGuLPW7i8k8gI/nZkQkyT75Ay7xgY8sAovJ4PbpWNbE1q9vaybt3GopbGhq8d1LpVxHZlhOV+UI21iM8hW4wSMgHIwT1Fc7Naa9Jp0Uccd8joblY83a7lZmBgd2D/OiKSrAliSOjcGuvorAZyPj7wdL42sbDTjeC2tYrrz5GWLc+RG6jBJxjLDjH41xY+ANkP+Y/cf+A6/412kPhq/iivYdtkUnEW7943+kMkpcmQbP+WisVbrjA+8Dx0Gj2Umn6ZHbS7Ays7BIzlIwzlgi8DhQQo4HA6DpUyhGWrRvSxVakrQlZHln/Cg7PIP9v3HH/Tuv/xVPHwItB/zH5//AAGH/wAVXr1FR7Gn2N1mWLW1RnkY+BVqGz/b0/TH/HsP/iqePgdbj/mPTf8AgMP/AIqvWaKX1el/KWs2xq2qM8oX4JQI24a9LnGP+PYf/FVMPg3GvTXX/wDAUf8AxVeo0UnhaL3iV/bOP/5+s8xHwfQPu/tx84x/x6j/AOKqQfCbb01w/wDgKP8A4uvSqpavZyX+lXFrEVDyLgBzhW5yVb2PQ+xqXg8O/sIn+18d/wA/WcBH8Kplupf+JwBEUUq/2bktzkY3cYAX8/arI+GDj/mNj/wE/wDs61pfDFzPp9pC9pYrPa3Uk1vKtw+LZWfeNq7MNtztCnAwo6ZwOtqXgMM94In+1MZ/z8Z58nwzkR2Ya0MtjP8Aon/2dSj4dTr01pP/AAE/+zrvKKl5dhXvBE/2jiv52cLH8PbiPO3WU5OebT/7OpR4Fux01mP/AMAz/wDF12tFJ5Zg3vTRP17E/wA7OLj8D3cW7brEXzHJzZn/AOLqYeD75emrw/8AgGf/AI5XXUVLyrBPemifrdf+ZnJReEb+GMIurwYHrZn/AOOVMPDWpDpq1v8A+AR/+OV09ZOv6fcaha2y2sFrLNDdwzhrhymwJIrEqQrHJAI7detS8nwL3pIj6xV/mMuz8O6vHBsfU7VCHYAfZCcruOD/AKzuMH2zVkaHqo/5i1r/AOAR/wDjlOi0a8PiFNTlW2UlxI7o5aQDydhgHyjMe7585HP8Pet+k8lwD19kiXVm92c3H4f1SKMIurW2B62R/wDjlK+gao6Mp1a1wRg/6Ef/AI5XR0VSyjArX2SM277nNr4e1RUVRq1tgDH/AB5H/wCOUj+HdTdGU6tbYIx/x5H/AOOV0tFUsswa2polwi90cyvhvU1RVGrW2AMf8eR/+OU2XwxqMsZRtXt8H0sj/wDHK6iitFgMMtoIh0ab3RzP/CNakf8AmLW3/gEf/jlV7zwvqr2cwi1S2eTYSiG0KhmHIGfM45rrqx7mwu/+EotdSt7WzMS2zwTSPKVlbc8Z6BDkKEbGT1bt1qlhKC2iiXhaL+yZ3/CL6j/0F7f/AMAj/wDHKZJ4Sv5NudXg+VgwxZn/AOOVo+HdFbRReIsEFtBNIHSCGZpQDj5mLMAeeOOgx71t1oqFNbIh4HDv7COU/wCEU1D/AKC8H/gEf/jlMfwffSbc6xD8rBhizP8A8crrqKtQitkR/Z2F/kRyf/CI3/8A0F4P/AM//HKjfwbeybc6xD8rBhizP/xyuwoq02tiP7Kwf/PtHI/8Ifff9BiH/wAAz/8AHKifwReSOrnWIsrnGLM9/wDtpXZ0VaqzWzJ/sfA/8+kcUfAl2f8AmMx/+Af/ANnVWf4fXjTQMmrRsNxDn7Ljau08/f55wPxrv65RNAv2g1BGs9Nh+1TpL5cErCN1SRWCMvlj767gzcnnoRgDRYuutpMFk+BW1JGafhzO3XWl/wDAT/7Ooj8MZDJvOtjOMf8AHp/9nXaaPZyafpcNtKU3JuO2M5VAWJCLwPlUEKOBwBwKvVX17EfzspZVgltTR52fhc7ddb/8lP8A7Oom+E259x1w5xj/AI9f/s69Joqv7QxX87KWW4RbU0eZt8Iw3XXG/wDAUf8AxVRf8KbjLE/26/P/AE6j/wCKr1Giq/tPF/8APxlrAYZbQR5YfgvE3XXZP/AYf/FVEfghblif7el5/wCnYf8AxVesUVX9q43/AJ+MpYOgtoo8kb4GWzddem/8Bh/8VVcfAeA3MgbXJhFsUq/kLktk5GM9gF/P2r1jU4ZbjTLiKCC1nldCEius+Ux7bsA8fhXPSeGbiWw0uE2dmtxYS7op/tTnyssrMyoIwp5BAU4AAGMdq/tjHf8AP1lLDUltE4k/AGybrr9x/wCA6/41Gf2e7A/8zBc/+A6/417NRVf21mH/AD9ZapQWyPGP+GetP/6GC6/8B1/xpP8AhnjT/wDoYLr/AMB1/wAa9oooedZg96rH7OHY8X/4Z40//oYLr/wHX/Gk/wCGdtO/6GC6/wDAdf8AGvaaKzebY171GHs49jxX/hnXTv8AoYbr/wAB1/xpP+GdNN/6GG7/APAdf8a9roqHmWLe9Rj5I9jxP/hnPTT/AMzDd/8Afhf8a9f0e2lstEsLWfb50NtHHJt6bgoBx+Iqhr2lXWoXNjPaQWjSWsqyh5n2nhgdoPlswBx1DD3DDKlNM0a4s9buLyTyAj+dmRCTJPvkDLvGBjywCi8ng9ulYVa9Stb2jvYaSWxu0UUViMKKKKACiiigDI8T+IrTwn4dutbv455La22b0gUFzudUGASB1Yd687/4aG8J/wDQN1z/AL8w/wDx2t740f8AJJdb+tv/AOlEdfJor2MtwFLEwbnfRmc5NbH0r/w0L4U/6Bmt/wDfmH/47S/8NCeFP+gZrf8A35h/+O182CnCvbp5BhZbt/f/AMAzdWR9Jf8ADQfhX/oGa3/35h/+O0f8NBeFf+gZrf8A35h/+O183iniuyHDOClu5fev8iXWkfR3/DQPhb/oGa3/AN+Yf/jtH/DQHhb/AKBmt/8AfmH/AOO185inCupcJ4C28vvX+RPt5H0aPj34Ybppmtf9+of/AI7Ug+O3htuml61/36h/+O187x9qtx1nLhXArrL71/kYTxdRdj6AHxx8PN00vWf+/cP/AMdp4+NegN00rWf+/cP/AMdrweLtVuKsJcNYJdZfev8AI5p5jWXY9xHxm0Jumlax/wB+4f8A47Tx8YNFPTSdY/74h/8AjteLRdKuR9a55cPYRdX9/wDwDmnm2IXb+vmewD4t6Q3TSdX/AO+If/jtSD4raU3TSNX/AO+Yf/jteSx9KuR1hLIsKur+/wD4BzTzzErZL7v+Ceoj4oaa3TSNW/75h/8AjtSD4l6e3TSNW/KD/wCO15rH1q3F2rCWT4ddX9//AADmnxFi1sl93/BPQh8RrJumj6r+UH/x2nj4gWh6aPqn5Qf/AB2uEj61cj6VzzyugtrnNPifGrpH7n/mdmvju3bpo2qf+QP/AI7T/wDhNof+gNqf/kD/AOO1ysParI6VwVcLCL0JXFGOfSP3P/M6L/hNIv8AoDan/wCQP/jtL/wmcf8A0BtT/OD/AOO1z4p4riqRUdjphxHjJbqP3f8ABN7/AITKP/oC6n+cH/x2l/4TBP8AoC6n+cH/AMdrDFPFcc6so7HVDPcU90vu/wCCbP8Awl6/9AXU/wA4P/jtL/wlqn/mC6n+cH/x2sgU9a46mMqR2sdUM2ry3t/XzNX/AISwf9AXUvzg/wDjtL/wlf8A1BdS/OD/AOO1mCniuKeaV47WOyGPqy3saP8AwlP/AFBdS/OD/wCO0v8AwlB/6Aupf99Qf/HazxTxXHPO8THZL7v+CdUMTN7l7/hJm/6Ampf99Qf/AB2l/wCElb/oCal/31B/8dqmKeK5J8RYyOyj93/BOqE29yz/AMJI/wD0BNS/76g/+O0p8SOP+YJqX/fUH/x2q4pW6Vxz4pxy2Ufuf+Z106aluSnxOR10XUv++oP/AI7TD4sC9dF1L84P/jtV5OlVZO9EeKcc+kfuf+Z2wwdN73Lx8YIOui6n+cH/AMdph8awjro2p/8AkD/47WVJVSSto8TY19I/c/8AM6Y5dRfc3T46tx10bVP/ACB/8dph8f2o66Pqn/kD/wCO1zcneqsnWt48Q4t9F93/AATojlNB9/6+R1Z+IlmOuj6r+UH/AMdph+JNgOuj6r+UH/x2uNkqpJW8c9xT6L7v+CdEckwz6v7/APgHdH4nacOukat/3zB/8dph+KeljrpGrf8AfMP/AMdrz6TvVWSto5ziX0X3f8E6Y8P4R9X9/wDwD0c/FjSR10nV/wDvmH/47TD8XdHXrpOr/wDfEP8A8drzGSqknet45rXfY6IcNYJ9Zfev8j1Y/GLRF66TrH/fEP8A8dqNvjToK9dK1n/v3D/8dryOTpVGfvXRDMaz3sb/AOq2BtvL71/kezH43eHh/wAwvWf+/cP/AMdpD8cfDo/5hes/9+of/jteHv0qBq3WMqPsYVOGsHHZy+9f5Hup+Ovhwf8AML1r/v1D/wDHaafjx4aH/ML1r/v1D/8AHa8GaoWrRYqbOKpkWFjs39//AAD30/HvwyP+YZrX/fqH/wCO16ZZ3SX1jb3cQYRzxrKoYcgMMjP518Xv1r7F8Pf8i1pX/XnD/wCgCumjUc73PDzLB08M48l9b/oaVFFFbnlhRRRQAUUUUAZHifxFaeE/Dt1rd9HPJbW2zekCgudzqgwCQOrDvXnf/DQ3hP8A6Buuf9+Yf/jtb3xo/wCSS639bf8A9KI6+TRXsZbgaWJg3O+jM5ya2PpX/hoXwp/0DNb/AO/MP/x2l/4aE8Kf9AzW/wDvzD/8dr5sFOFe1TyDCy3b+/8A4Bm6sj6R/wCGg/Cv/QM1v/vzD/8AHaX/AIaC8K/9AzW/+/MP/wAdr5vFPFdsOGcFLdy+9f5EutI+jv8AhoHwt/0DNb/78w//AB2j/hoDwt/0DNb/AO/MP/x2vnMU4V1LhPANby+9f5E+3kfRo+Pfhg9NM1r/AL9Q/wDx2pB8dvDbdNL1r/v1D/8AHa+d4+1W4+tZy4VwK6y+9f5GE8XUXY+gB8cfDrdNL1n/AL9w/wDx2nj416A3TStZ/wC/cP8A8drwePpVuKsJcNYJdZfev8jmnmNZdj3EfGbQm6aVrH/fuH/47Tx8YdFbppOsf98Q/wDx2vFoqtx9q55cPYRdX9//AADmnm2IXb+vmewj4t6Q3TSdX/74h/8AjtSD4raU3TSNX/75h/8AjteSx1bj6VhLIsKur+//AIBzTzzErovu/wCCepD4o6aemkat/wB8w/8Ax2pB8S9PbppGrflB/wDHa80j6ircfSsJZPh11f3/APAOafEOLXRfd/wT0MfEayPTR9V/KD/47Tx8QbQ9NH1T8oP/AI7XCR1bjrnnldBLS5zS4nxq6R+5/wCZ2i+O7dumjap/5A/+O0//AITaH/oDan/5A/8AjtcpF0FWhXBVwsI7ErijHPpH7n/mdF/wmkX/AEBtT/8AIH/x2l/4TOP/AKA2p/nB/wDHa58U8VxTio7HTDiPGPdR+7/gm9/wmUf/AEBdT/OD/wCO0v8AwmCf9AXU/wA4P/jtYYp4rjnVlHY6oZ7inul93/BNn/hL1/6Aup/nB/8AHaX/AIS5f+gLqf5wf/HayBTxXHUxlSO1jqhm2Ilvb+vmav8Awlg/6AupfnB/8dpf+Er/AOoLqX5wf/HazBThXFPNK8drHXDMKst7Gl/wlP8A1BdS/OD/AOO0v/CUH/oC6l/31B/8drPFPFck87xMdkvu/wCCdcMTOW5e/wCEmb/oCal/31B/8dpf+Elb/oCal/31B/8AHapiniuOfEWMjso/d/wTqhNvcs/8JI//AEBNS/76g/8AjtL/AMJI4/5gmpf99Qf/AB2q4pT0rknxTjlso/c/8zrpwUtyU+JyOui6l/31B/8AHaYfFYXroupfnB/8dqvJ0qrJSjxVjn0j9z/zO2GDpvuXz4wQddF1P84P/jtRnxrEOujan/5A/wDjtZUlVJOlbR4mxr6R+5/5nRHLqL7m6fHVuOujap/5A/8AjtMPj+1HXR9U/wDIH/x2ubk61VkrePEWLfRfd/wTpjlNB9zqz8RLNeuj6r+UH/x2mH4k2C9dH1X8oP8A47XGyVVk71vHPcU+i+7/AIJ0RyTDPq/v/wCAdyfidpw66Rq3/fMH/wAdph+KeljrpGrf98w//Ha8+kqrJW0c5xL6L7v+CdMeH8I+r+//AIB6OfixpI66Tq//AHzD/wDHaYfi7o466Tq//fEP/wAdrzGSqslbxzWu+x0R4awT6y+9f5Hqp+MWiL10nWP++If/AI7UbfGnQV66VrP/AH7h/wDjteRydKoz9K3hmNZ72N/9VsDbeX3r/I9nPxu8Pj/mF6z/AN+4f/jtNPxx8Oj/AJhes/8AfqH/AOO14e9Qt0rdY2o+xhU4awcdnL71/ke6H46eHB/zC9a/79Q//Haafjx4aH/ML1r/AL9Q/wDx2vBmqF61WKmziqZFhY7N/f8A8A99Px78Mj/mGa1/36h/+O16ZZ3SX1jb3cQYRzxrKoYcgMMjP518XvX2L4e/5FrSv+vOH/0AV00ajne54eZYOnhnHkvrf9DSooorc8sKKKKACiiigAooooA4L40f8kl1v62//pRHXyaK+0PGvh3/AISzwpd6IZWiW6eENIuMqqyoxIz7Ka8z/wCGc9N/6GG7/wC/C/417OWY6lh4NVO5nOLex8/CnCvoD/hnXTf+hguv/Adf8aX/AIZ207/oYLr/AMB1/wAa9ynnuEju39xk6UjwAU8V75/wzvp3/QwXX/gOv+NL/wAM86f/ANDBdf8AgOv+NdtPiPAx3b+4l0ZHggpwr3ab9ny0RAYdcuHbeoIMKj5Sw3Hr2GT+FS/8M96f/wBDBc/+A6/411rinL7bv7ifYTPDY+1W469qX4AWK9NfuP8AwHX/ABqRfgNZr016f/wHH/xVZy4ny99X9xhPCVHseORdqtxV66vwMtV6a7N/4DD/AOKqRfgnbr012X/wGH/xVc8uI8C+r+45p5fWex5TF0q5H1r0qT4NCOMGHWWdt6ghrcD5Sw3H73YZP4VYHwfRemtt/wCAv/2dc8s/wb6v7jlnlWIeyX3nm8fSrkdegL8JQvTWj/4C/wD2dSL8KyvTWv8AyV/+zrCWd4R9X9xzTyXFPZL7zho+tW4u1dmvwydemsr/AOAv/wBnUi/DiVemsJ/4CH/4uueWb4Z9X9xyzyDGPZL7zko+tXI+ldDJ8P7yNAYdUidt6ghrUj5Sw3H7/YZP4VYHgS5Xpq8X/gIf/jlc88yoPa/3HLPhvHvZL7zCh7VZHSthfBV2vTVoP/AM/wDxypP+EQvf+grb/wDgG3/xyvPq4mnJ6ErhrHrovvMYU8Vr/wDCJXo/5itv/wCAbf8Axyl/4RO+H/MVt/8AwDb/AOOVw1JKWx1Q4fxq3S+8yhTxV+fwtqaW8jQajayShCURrRgGbHAJ8zjmpf8AhFr/AP6Clt/4Bt/8crjqU5S2OqGSYtbpfeZwp61f/wCEY1D/AKClt/4Bt/8AHKcPDWoD/mKWv/gG3/xyuKphKstjshlWIW6X3lEU8Vc/4RvUf+gpa/8AgE3/AMdpf+Ed1Ef8xO1/8Am/+O1w1Msry2t951wwNZblQU8VLPoGrpbyNBqFnJKEJRGtGAZscAnzOOal/sDUv+gnaf8AgE3/AMdrjnk2Klsl951ww01uQCnipf7C1P8A6Cdp/wCATf8Ax2l/sPU/+gnZ/wDgE3/x2uOpw/jZbJfedcINbkYpW6VJ/Yuqf9BKz/8AAJv/AI7QdF1Q/wDMSs//AACb/wCO1xz4Zx72S+866c1HcqydKqyd6vXGh6wLeVodQsnlCEohtGAZscAnzOKG8Oai3XVLX/wCb/47Sjwxj10X3ndDF01uY8lVJK3z4Wv266pbf+Abf/HKjPhC9brqtv8A+Abf/HK3jw5jl0X3nTDMKK3Oak71Vk611h8F3bddWg/8Az/8cqM+BbluurQ/+Ah/+OVvHIMYui+86Y5ph11f3HGyVUkrtbjwDeC3laHVIXlCEohtSAzY4BO/jmg/DqVuusJ/4CH/AOLreOSYtdF950wznCrdv7jgJO9VZK9FPw0kbrrK/wDgJ/8AZ1Gfhcx66yP/AAF/+zreOT4ldF950wz7Brdv7jzSSqknevU2+FO7rrX/AJK//Z1G3wjVuutn/wABf/s63jleIXRfedMOI8Ct2/uPJpOlUZ+9ev3PwfIt5Wh1gvKFJRDb4DNjgE7uOaRvgrA/XXJP/AYf/FV0Qy+stzf/AFnwFt39x4s/SoGr2s/A+2P/ADHZf/AYf/FU0/Au1P8AzHpv/AYf/FVusHVRz1OIsDLZv7jxBqhavcz8CLQ/8x6f/wABx/8AFUw/ASzP/Men/wDAcf8AxVarDVEcVTO8JLZv7jwh+tfYvh7/AJFrSv8Arzh/9AFeXn4A2J/5j9x/4Dr/AI16rpFtLZ6LYWs2PNht44329NwUA4/KuqjTlC9zwczxdPEOPs+l/wBC7RRRW55QUUUUAFFFFAHBfGj/AJJLrf1t/wD0ojr5NFfaHjXw7/wlnhS70QytEt08IaRcZVVlRiRn2U15n/wznpv/AEMN3/34X/GvZyzG0sPBqp3M5xb2Pn4U4V9Af8M66b/0MF1/4Dr/AI0v/DO2nf8AQwXX/gOv+Ne5TzzCR3b+4ydKR4AKeK98/wCGd9O/6GC6/wDAdf8AGl/4Z50//oYLr/wHX/Gu2nxHgY7t/cS6MjwQU4V7tN+z5aKgMOuXDtvUEGFR8pYbj17DJ/Cpf+Ge9P8A+hguf/Adf8a61xTl9t39xPsJnhsfarcfWvah8ALFemv3H/gOv+NSL8BrNemvT/8AgOP/AIqs5cT5e+r+4wnhKj2PHI+lW4q9dX4GWq9Ndm/8Bh/8VUi/BK3Xprsv/gMP/iq55cR4F9X9xzTy+s9jymKrcfavS5Pg0I0Bh1lnbeoIa3A+XcNx+92GT+FWF+D8a9Nbb/wF/wDsq55Z/g31f3HLPKsQ9kvvPN46tx9K9BX4SqvTWj/4C/8A2dSL8KyvTWv/ACV/+zrCWd4R9X9xzTyXFPZL7zhY+oq3H0rsx8MXXprK/wDgL/8AZ1Ivw3lXprCf+Ah/+Lrnlm+GfV/ccs8gxj2S+85OOrcddFN8P7uO3kaDVIpJQhKI1sQGbHAJ38c1YXwJcr01eL/wEP8A8crnnmVBrS/3HNPhvHvZL7zBi6CrQrYXwVdr01aD/wAAz/8AHKk/4RC9/wCgrb/+AZ/+OV59XFU5PQhcM49dF95jCnitf/hEr0f8xW3/APANv/jlL/wid9/0Fbf/AMA2/wDjlcNSSlsdUOH8at0vvMoU8Vfn8LamlvI0Go2skoQlEa0YBmxwCfM45qX/AIRa/wD+gpbf+Abf/HK4505S2OqGSYtbpfeZwp4q/wD8IxqH/QUtv/ANv/jlL/wjWoD/AJilr/4Bt/8AHK4qmEqS2OuGVYhbpfeUhThV3/hG9R/6Clr/AOATf/HaX/hHdRH/ADE7X/wCb/47XFUy2vLa33nZDA1luVBTxUs+gaulvK0GoWckoQlEa0YBmxwCfM45qX+wNS/6Cdp/4BN/8drink2Klsl951ww1RbkAp4qX+wtT/6Cdp/4BN/8dpf7D1P/AKCdn/4BN/8AHa46nD+Nlsl951wg1uRilPSpP7F1T/oJWf8A4BN/8doOi6of+YlZ/wDgE3/x2uOfDOPeyX3nXTmo7lWTpVWSr1xoesC3laHULJ5QhKIbRgGbHAJ8zihvDmot11O1/wDAJv8A47Sjwxj10X3ndDF01uY8lVJOlb58LX5/5ilt/wCAbf8Axyoz4QvW66rb/wDgG3/xyt48OY5dF950wzCitzmpOtVZK6w+C7tuurQf+AZ/+OVGfAty3XVof/AQ/wDxyt45BjF0X3nTDNMOt2/uONkqrJ3rtLnwDeC3laHVIXlCEohtSAzY4BO/ihvh1M3XWE/8BD/8XW8ckxa6L7zphnOFW7f3HASVVkr0VvhpI3XWF/8AAT/7Ooz8LmPXWR/4C/8A2dbxyfErovvOmGfYNbt/ceaSVVkr1JvhTu661/5K/wD2dRn4Rq3XWz/4C/8A2dbxyvELovvOiHEeBXV/ceTSdKoz9K9jb4Oxt11tv/Ab/wCyqtF8F1lt42uNZeOUoC6LbghWxyAd3NdEMvrLc6P9Z8Bbd/ceNPULdK9qPwPtj/zHZf8AwGH/AMVTT8C7U/8AMem/8Bh/8VW6wdVHPU4iwMtm/uPEGqF69zPwItD/AMx6f/wHH/xVMPwEsz/zHp//AAHH/wAVWqw1RHFUzvCS2b+48IevsXw9/wAi1pX/AF5w/wDoAry8/AGxP/MfuP8AwHX/AOKr1XSLaWz0WwtZsebDbxxvt6bgoBx+VdVGnKF7ng5ni6eIcfZ9L/oXaKKK3PKCiiigAooooAKwLmHVF8Tw3EcdzJaeYoJScCMRFGUgoWGWDlWzjO0YB7HfooA4d9O8SGyK4vy21RIq3oDPP5cgaRW3/LGWMZCZHT7vY2ZbDXzPKW+2OpkJn8u62idPOUoIhvHl4jDA/dyT1PWuvooAp6THcxaRZx3pJuVhUSlm3HdjnJ7n371coooAKwLmHVF8Tw3EcdzJaeYoJScCMRFGUgoWGWDlWzjO0YB7HfooA4lrHXxp7RNBqUk0giRpI74ApKA/mzL+8AKH5NqHAyOVAzVzVtL1V5NTmsXvF8x4EijF25BUNvkdV8xdud23aGQ4Q44PPVUUAU9JS5i0axjvN32pLeNZt77m3hRnJ7nPerlFFABWBcw6ovieG4jjuZLTzFBKTgRiIoykFCwywcq2cZ2jAPY79FAHDzaX4gk0xl3aqt0nlunlXqgSXAWQOzEvkQMWT5Bgjbwoqe80vWjc6gLWbVhazTI7brpPMxvYsLf5sIuCPvYOBj2rsaKAKumrcrpdot5gXQhQTYORv2jdz35zVqiigArAuYdUXxPDcRx3Mlp5iglJwIxEUZSChYZYOVbOM7RgHsd+igDgb3TPFL6Z5cBvhN8nmEXeWeURyB3X94uIyxjwuRgrnZ2OtcWmqPf6i6wX32eSAKircn/WBh8y/vhx16eVxxhs8dRRQBT0iO5i0axjvSxukt0WYs+4lwozk9znvVyiigArA1yHVG1G2nso7mWKPYwWCcIARIpkDAsA25AVGcgH0zmt+igDltTttVkh8RxW1tqOZoF+xsl4FzLgj5D5gKDO0n7o4PXvcdbyTW3lNlqAtJbAq4FyoAkzkKAJPlfGRuXAz/F3rdooAztAW7Tw9p0d/HNHeR26JMJpA7lwoBJYEg5IznJ61o0UUAFYGuQ6o2o209lHcyxR7GCwThACJFMgYFgG3ICozkA+mc1v0UAc5PBqBudaRbXUDDKIzbsLkYLD72wCZWUdMjK5wayotM8StMgmkvkT7CIneO5DHf5OCRmQAtv/ANkHPPmAfLXcUUAZfh+G/t9ISPUWYz+ZIQHbcyoXJUE7m5Ax/Efqa1KKKACsDXIdUbUbaeyjuZYo9jBYJwgBEimQMCwDbkBUZyAfTOa36KAOW1G01iWfVorVNQihlkt2SRLlMsob96Isv8hK+oUVFYafr39oQvftdvCYFjdUucIqeUARw27zvMH3umCfmziuuooAztCtLiy0a2hu3le4wXl8yVpSrMSxUMxJIXO0c9AK0aKKACsDXIdUbUbaeyjuZYo9jBYJwgBEimQMCwDbkBUZyAfTOa36KAOYvbTU5Z9YSGPUY4naF4HW4BDkH5wg81SoPGRlO/41rex8XR3Vq7ywPs05oSzXTBVlIjO502nc2Q+DyMDtk7uwooAzPD0F9baJbw6iG+1IWDbpzMcbjty55PGK06KKACsLxNbanc20Y03zywWQAQzeUVkK/u3Y5GVB6jnORwcVu0UAczeLqMuqT3CaXqYhFmf3aXyoJ5CBhABLiMrjG4Dkk84HzQW9rqqz6aJbTVSsUMj3EgvRhid+Idvm84yMOcnhOc5K9bRQBleHoruLSx9tiuIZmkZvJnm80xgnhQ+5iwHqT68AYA1aKKACsLxNbanc20Y03zywWQAQzeUVkK/u3Y5GVB6jnORwcVu0UAcxe29/falcEWur2lu9mVLwXiBmdgOFUybUZf7wHJzzgfNTi0/xIDbhXuldYwsUktyGEeHlyZlDEOzIYugYBgcYFdnRQBj+HYL23s51vFuEzKDElzP5rqvloDltzdXDnr0PboNiiigAooooAKoXcWptOGt7m2W2AG+FrZmkb1CuJFAJHTjir9FAHLW9trKaJcW7RXnmPciRM3A3iEPHuQMZGKsy+YB856ZymQBTbTvExx5UlykvkuIne6BSNSkwRHG47nDNFlsHp1457WigDM0SG5htZhcJPGhmJhjuJvNkRMDhm3Nn5tx6nAIHatOiigArPvotTaQtbXFv9mC/PB9nYyv6hZPNUKSOASODzWhRQBzllFqK+HpYbiz1BbiXe4jN0HaHkbUEnmhmHOc7hna33cqKrW1lrKSaUZ4r5hFbFLzF1/rD5bAAfvcbs7e2c8+Zx83WUUAZPh2C+t9JEeoeb5wlkI82Te20sSuSWft23N9ew1qKKACqN1FqT3Ktb3Vqlrgb4ntmaRvXDiRQMjpxx71eooA5dIdZfwsbNYL+C7MrfO9wjyiLzgdofefnMRIUk8EckcZgsdO1/wDtaza7kvDaCLy3X7SMKm2QZfDZMufL+YZA5weBXX0UAZmg2lxZ6TGl00xnkZpXWadpWj3MSE3EnO0EL1xxmtOiigAqjdRak9yrW91apa4G+J7Zmkb1w4kUDI6cce9XqKAOZgi1Y+Ghata38dwJmY77lDIYhODt3iQnc0ZIBzxjkjgmg2neJjjypLlJfJcRO90CkalJgiONx3OGaLLYPTrxz2tFAGZokNzDazC4SeNDMTDHcTebIiYHDNubPzbj1OAQO1adFFABWZqkGtSpL/Zt9ZwAxkKstqztuwedwkAHb+E49+ladFAHO2cGoR6TeQS2uoYkJa3V7oPLGAsYwX80EkvvYYfoMZHAqrbWWspLpBmivmEMDrdkXXDjY4Uf6372SvYnOP3hxz1lFAGT4dgvrfTGjv8AzfM85yglk3sEJyM5d8fTc34dBrUUUAFZutw3M1lGtss7gTKZY7eXypHTuFbcuOcHqOAR3rSooA5gW+svHoReK7FzAEF4/ngI+OGJxJg5wSMo2cgfIeRTNhrn9nTpHFqCzG98y33XZyke1OGPnnjIbqXHU7OcV2dFABRRRQAVDdLvtJV2SyZQ/JC+x29lbIwffI+tTUUAcg1v4guNC0yMwX0V3DaGKTN0gb7TsTZK5V/nQHfkEknP3TTTYa79omaVL+S3M+ZY47wK0y7pseUd42AbocjK5Ckc9+xooAraelzHptql44e6WFBM46M+BuP55qzRRQAVDdLvtJV2SyZQ/JC+x29lbIwffI+tTUUAcbLbeI7jRdHQR30NzbQBZ0NxGGknCJtZ2DndHneGwdx9O9Kul69Jb6xFNJeNJO/+jyfa9gRjI4DJsYEIqFCVOCSD14rsaKAGRRiKFIwzMEUKC7FmOPUnkn3p9FFABVTU/M/sy48q2muZNh2wwzeU7nsA+Rt+uat0UAcpNa6u+iWMEUOoC9imUs73AHy5DEEibledoLbzhTlWzksu7HWJF1kRQ6gFm2fZALo5VsvuYfvxgcr0ZR0+Tg566igCK2EgtYRN/rQi7/mz82OecDP5VLRRQAVXvVZrGZUSV2KEBYm2sT7HcuP++h9asUUAccbPxC+kW0Oy9W4iE6gi6AIkLAwyMfMJZFG4FSzH2fG4hsNd+0TNKl/JbmfMscd4FaZd02PKO8bAN0ORlchSOe/Y0UAVtPS5j021S8cPdLCgmcdGfA3H881ZoooAKKKKAP/Z", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig = pyplot.figure(figsize=(20,10))\n", "velocity.plot(colorbar=\"horizontal\", figure=(fig, 121))\n", "pressure.plot(colorbar=\"horizontal\", figure=(fig, 122))" ] } ], "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 }