{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Scattering from a sphere using a combined direct formulation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Background" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In this tutorial, we will solve the problem of scattering from the unit sphere $\\Omega$ using a combined integral formulation and an incident wave defined by\n", "\n", "$$\n", "u^{\\text{inc}}(\\mathbf x) = \\mathrm{e}^{\\mathrm{i} k x}.\n", "$$\n", "\n", "where $\\mathbf x = (x, y, z)$.\n", "\n", "The PDE is given by the Helmholtz equation:\n", "\n", "$$\n", "\\Delta u + k^2 u = 0, \\quad \\text{ in } \\mathbb{R}^3 \\backslash \\Omega,\n", "$$\n", "\n", "where $u=u^\\text{s}+u^\\text{inc}$ is the total acoustic field and $u^\\text{s}$ satisfies the Sommerfeld radiation condition\n", "\n", "$$\n", "\\frac{\\partial u^\\text{s}}{\\partial r}-\\mathrm{i}ku^\\text{s}=o(r^{-1})\n", "$$\n", "\n", "for $r:=|\\mathbf{x}|\\rightarrow\\infty$.\n", "\n", "From Green's representation formula, one can derive that\n", "\n", "$$\n", "u(\\mathbf x) = u^\\text{inc}-\\int_{\\Gamma}g(\\mathbf x,\\mathbf y)\\frac{\\partial u}{\\partial\\nu}(\\mathbf y)\\mathrm{d}\\mathbf{y}.\n", "$$\n", "\n", "Here, $g(\\mathbf x, \\mathbf y)$ is the acoustic Green's function given by\n", "\n", "$$\n", "g(\\mathbf x, \\mathbf y):=\\frac{\\mathrm{e}^{\\mathrm{i} k |\\mathbf{x}-\\mathbf{y}|}}{4 \\pi |\\mathbf{x}-\\mathbf{y}|}.\n", "$$\n", "\n", "The problem has therefore been reduced to computing the normal derivative $u_\\nu:=\\frac{\\partial u}{\\partial\\nu}$ on the boundary $\\Gamma$. This is achieved using the following boundary integral equation formulation.\n", "\n", "$$\n", "(\\tfrac12\\mathsf{Id} + \\mathsf{K}' - \\mathrm{i} \\eta \\mathsf{V}) u_\\nu(\\mathbf{x}) = \\frac{\\partial u^{\\text{inc}}}{\\partial \\nu}(\\mathbf{x}) - \\mathrm{i} \\eta u^{\\text{inc}}(\\mathbf{x}), \\quad \\mathbf{x} \\in \\Gamma.\n", "$$\n", "\n", "where $\\mathsf{Id}$, $\\mathsf{K}'$ and $\\mathsf{V}$ are identity, adjoint double layer and single layer boundary operators. More details of the derivation of this formulation and its properties can be found in the article Chandler-Wilde et al (2012).\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Implementation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First we import the Bempp module and NumPy." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import bempp.api\n", "import numpy as np" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We define the wavenumber" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "k = 15." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The following command creates a sphere mesh." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "grid = bempp.api.shapes.regular_sphere(3)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "As basis functions, we use piecewise constant functions over the elements of the mesh. The corresponding space is initialised as follows." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "piecewise_const_space = bempp.api.function_space(grid, \"DP\", 0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now initialise the boundary operators.\n", "A boundary operator always takes at least three space arguments: a domain space, a range space and the test space (dual to the range). In this example we only work on the space $\\mathcal{L}^2(\\Gamma)$ and we can choose all spaces to be identical." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "identity = bempp.api.operators.boundary.sparse.identity(\n", " piecewise_const_space, piecewise_const_space, piecewise_const_space)\n", "adlp = bempp.api.operators.boundary.helmholtz.adjoint_double_layer(\n", " piecewise_const_space, piecewise_const_space, piecewise_const_space, k)\n", "slp = bempp.api.operators.boundary.helmholtz.single_layer(\n", " piecewise_const_space, piecewise_const_space, piecewise_const_space, k)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Standard arithmetic operators can be used to create linear combinations of boundary operators." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "lhs = 0.5 * identity + adlp - 1j * k * slp" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now form the right-hand side by defining a GridFunction using Python callable." ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "@bempp.api.complex_callable\n", "def combined_data(x, n, domain_index, result):\n", " result[0] = 1j * k * np.exp(1j * k * x[0]) * (n[0]-1)\n", "\n", "grid_fun = bempp.api.GridFunction(piecewise_const_space, fun=combined_data)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can now use GMRES to solve the problem." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [], "source": [ "from bempp.api.linalg import gmres\n", "neumann_fun, info = gmres(lhs, grid_fun, tol=1E-5)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "`gmres` returns a grid function `neumann_fun` and an integer `info`. When everything works fine info is equal to 0.\n", "\n", "At this stage, we have the surface solution of the integral equation. Now we will evaluate the solution in the domain of interest. We define the evaluation points as follows." ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [], "source": [ "Nx = 200\n", "Ny = 200\n", "xmin, xmax, ymin, ymax = [-3, 3, -3, 3]\n", "plot_grid = np.mgrid[xmin:xmax:Nx * 1j, ymin:ymax:Ny * 1j]\n", "points = np.vstack((plot_grid[0].ravel(),\n", " plot_grid[1].ravel(),\n", " np.zeros(plot_grid[0].size)))\n", "u_evaluated = np.zeros(points.shape[1], dtype=np.complex128)\n", "u_evaluated[:] = np.nan" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This will generate a grid of points in the $x$-$y$ plane.\n", "\n", "Then we create a single layer potential operator and use it to evaluate the solution at the evaluation points. The variable ``idx`` allows to compute the solution only at points located outside the unit circle of the plane. We use a single layer potential operator to evaluate the solution at the observation points." ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [], "source": [ "x, y, z = points\n", "idx = np.sqrt(x**2 + y**2) > 1.0\n", "\n", "from bempp.api.operators.potential import helmholtz as helmholtz_potential\n", "slp_pot = helmholtz_potential.single_layer(\n", " piecewise_const_space, points[:, idx], k)\n", "res = np.real(np.exp(1j *k * points[0, idx]) - slp_pot.evaluate(neumann_fun))\n", "u_evaluated[idx] = res.flat" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now plot the slice of the domain solution." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "Text(0.5, 1.0, 'Scattering from the unit sphere, solution in plane z=0')" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjYAAAHwCAYAAAC17yUBAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8li6FKAAAgAElEQVR4nOy9eZBv21Ue9u1zflMP9943aHqSnp4AIUPAIDPGZVeCE4opThmqsA1xSMDMthNcgdhhFMYYcCqOw+TIAhOhGEEYDMGYKkxisBiMCAghJGSDAEnv6T296U49/sadP/b61l57nfPr2/3u0Le796rqOt2/3zlnz6fP+va3vhVijKhWrVq1atWqVTsP1px2BapVq1atWrVq1e6U1RebatWqVatWrdq5sfpiU61atWrVqlU7N1ZfbKpVq1atWrVq58bqi021atWqVatW7dxYfbGpVq1atWrVqp0bqy821V6QhRC+IYTwg3fp3n8uhPCHIYTdEMLn3o0yTlifV4cQYghhcB/UZTeE8OH3oJwYQnjN3S7nblkI4YtDCL96G9e/IYTwzXeyTnLfu7luzvSYVat2p6y+2NxHFkL48yGEXw8h3AghXA0h/FoI4ZNv856dB3wI4U0hhG+/nfvGGL8jxvhlt3OPI+zbAHxfjHE7xvgzd6mMtRZCeF8I4dPvdbnHMemTPwbuzDhW618jMcavijH+/Ttd1l1eN2fSxHH4pRDCfgjh39+va6/a2bH6YnOfWAjhMoCfA/C9AB4C8AoAfw/A9DTr1Wf3ALl4DMC715QdQgh13p4Rux9Qrmr3vf0ogN8B8DCAbwTwkyGEF59ulaqdaYsx1p/74AfAJwG4fotzvhzAewDsAPh9AJ8gn/9PAP7IfP558vlHAzgEsASwC+A6gK8AMAcwk8/+pZz7cgA/BeBZAH8C4L835X4rgJ8E8M8B3ATwZfLZP5fvXw0gAvhvAXwAwHMAvtFcvwHghwFck/r/HQBPrGnjHwFYATiQ+o0B/DKAfwDg1+Tz10h9fxbAVQDvBfDlrr4/IfXdAfB7AF4L4OsBPAPgcQCfsab8/9OV/3eO0b7GjMHzAH4cwENr7v/FAH7VfRYBvEZ+fxOA7wfwr6TubwPwEf7cdePo7hsA/GNp8w0A7wTwsaacNwD4RSnn3wJ4zJXzVQD+UMbt+wEE8/1fl7G8BuAXeq79m3Ltn8hnHyVlXQXwHwD8lWOui6PacAXAm5Hm7PsBfBOAxvezGb+Bue8vI83jzhox/fPtbu29V+r/swBefty+cu35Vhxz3fRce5wx4zz6L5BeFm4izfdvNecdWS5ONp//pfQbf1YAvvgEz73XIjlvl8xnvwLgq+7Ec7X+XMyfU69A/ZGBAC7LQ+SHAXw2gAfd938ZwAcBfLI87F/Dh5p893J5IP1VAHsAHpHv9AFv7uUf2g2A3wbwLQBGAD4cwB8D+Ez5/luR/ol+rpy7seYB/QPy3cfLw+qj5fvvkofwgwBeifTPqffFRs5/H4BPN3//sjyAPwbAAMBQ7vdPAEwAvA7pn9t/bup7COAz5fw3I72sfaNc++WQf7jHLP9W7fvbAH5D2jYG8E8B/Oiae/eNh3+xuQrgU6TuPwLgx44499uPaMdnyrg+IHPmo828eBPSP8f/ROr83bZeUs7PybWvkv79LPnuc5H+yX+01PGbAPy6u/YXkZDHDQBbSP9cv0TO/wSkf6Yfc4x1cVQb3gzg/wZwScboDwB8qe9nHPFic5w1AuA/k/p+gvTV9wJ463H6qqc934pjrpuea48zZpwbnwbgTyOt148D8DSAz73T89nV77MAPAngUfn7nUjOVN/PP5FzPg/Ae9x9vg/A996tZ239Of8/FdK/TyzGeBPAn0d+4DwbQvjZEMJL5ZQvA/A/xxj/v5jsvTHG98u1PxFjfDLGuIox/l9InuOnnKD4Twbw4hjjt8UYZzFxOH4AwBeYc/5djPFnpIyDNff5ezHGgxjj7wL4XaQHJgD8FQDfEWO8FmN8AsD3nKButDfFGN8dY1wAeBlSX/3dGONhjPEdAH4QwBeZ838lxvgLcv5PAHgxgO+KMc4B/BiAV4cQHjhhHda17yuRPN4nYoxTpH9en38b2zD/Isb4m1L3H0F6cXshNkf6p/9RSAjCe2KMT5nv/1WM8a1S528E8GdDCI+a778rxng9xvgBAL9k6vGVAL5T7rcA8B0AXhdCeMxc+50xxqsyV/4igPfFGP+PGOMixvh2JHTw819oG0IILdJL/NfHGHdijO8D8I9QzoE7ZX8NwA/FGN8uffX1SH31anPOur46jq2bV312qzEDAMQYfznG+HuyXt+JtN3znx6z3BPP5xDCa5FeNP9qjPFxqcPHxRgfWPPzN+TSbSQkztoNpDGvVu0FWX2xuY9MHtpfHGN8JYCPRUJh/jf5+lEkaLhjIYT/JoTwjhDC9RDCdbn2RSco+jEAL+f1co9vAPBSc87jx7jPh8zv+0gPLUg77PXHuZc3e83LAVyNMe6Yz96PxEuiPW1+PwDwXIxxaf6Gqd9xbV37HgPw06bv3oO0tfFSvDBbV86JLMb4b5C83+8H8HQI4Y3C5aI9bs7dRUKKXn6MejwG4LtNe68ioSm2/+14PQbgU938+mtIL6gvtA0vQkIX329O93PgTtnLbTnSV8+7sm5nzE5y7a3GDAAQQvhUIeQ+G0K4gbRV5p8Jd2Q+hxCuICFn3xxj/JUj6t5nu0hotbXLSMhUtWovyOqLzX1qMcZ/jwQ9f6x89DiAj/DniZf8AwD+FoCHY4wPAHgX0j8aICFAndu7vx9H2pqxHtWlGOPnHHHNSewpJFib1vEwj2G2/CcBPBRCsF7dq5C26u6EnbStjwP4bNd/kxhjX332AGzyjxDCLf+53049Y4zfE2P8RKRtvNcC+B/N1zoOIYRtpK2jJ49R7uMAvtK1dyPG+Otr6vY4gH/rzt+OMX71Mcpa14bnkNAcixKtmwN7ctw0n9l+v1U/PmnLCSFsIRFd79R8O4kdd8zegsQFejTGeAWJmxN6zuuzY89nIfK/BcAvxRj/qfvu3SJP0PfzBjnt3QA+3K3lj8ea4IFq1Y5j9cXmPrEQwkeFEL42hPBK+ftRAF+ItNcNpK2WrwshfKJEBr1GXmq2kB7Mz8p1X4L8MgQk5OKVIYSR+8xqofwmgJshhL8bQtgIIbQhhI+93VBzYz8O4OtDCA+GEF6B9BL2gk2g7l8H8J0hhEkI4eMAfCnSts2dMN8/t7I3APgH3IoJIbw4hPCX1pz7uwA+JoTwuhDCBAnmvyv1DCF8snjuQ6R/7iTJ0j5HJAZGAP4+gLdxG+EW9gak8fwYKedKCOEvH3H+zwF4bQjhi0IIQ/n55BDCR8v1XxxCeN9J2iDo248j9fsl6fv/AYkwXliM8Vmkl5D/Wub2X0fpJPStEWtvAfAlMmZjpK23t8n21722447ZJSRU8zCE8CkA/qsTlHGS+fwPkJ5BX+O/iDF+jLzA9v18lZzzBwDeAeD1spY/D4kT9FMnqG+1aoXVF5v7x3YAfCqAt4UQ9pBeaN4F4GuBxKNBeoi8Rc79GaRIhd9H4hb8O6QH9J9Gih6i/Rsk7+dDIYTn5LN/BuA/Eqj5Z+SfxH+JxAv4EyRv+AeRok7uhH0bgCfk3v8PUoTV7YaxfyESCfJJAD8N4PUxxl+8zXvSvhPAN0n/fN0xzv9uJO/4X4cQdpDG7lP7TpQH+bch9cMfAnjBInJw49jz/WUkNO8a0lbK8wD+F/P9WwC8Hmk74xORtoduaTHGnwbwDwH8WAjhJtI8/ewjzt8B8BlInK0nkbZA/iESMRVIKMSv9V99ZBv+O6SXnT9G6se3APihNff5ciSk53kk5MeiS31rxNb//wXwzUj/bJ9Cein6An/ePbLjjtnfAPBtMh+/Bekl8Lh27PmMtA7/YwDXDBpzrHlk7AuQokKvIQUafL68jFar9oIsxHg7OwzVqp3cQghfDeALYoyezFjtHlkI4U1IkWnfdB/U5V8D+JoY43tOuy73s91PY1at2v1sFbGpdtcthPBISGkSmhDCn0JCoX76tOtV7f6wGONn1JeaatWq3Sk7NVVQ4Re8FQmOHgD4yRjj60+rPtXuqo2QtDA+DEnD4seQNGiqVatWrVq1O2qnthUVQggAtmKMu0IM/FUkOPo3bnFptWrVqlWrVq1ar50aYhPTG9Wu/DmUn0r4qVatWrVq1aq9YDtVjo2EXr4DKQ/ML8YY33aa9alWrVq1atWqnW071cy7Emb8upCk7X86hPCxMcZ32XNCCF+BlPAPYTL8xMkrH8bmYA4AGDWLfC/RnpquUpMO5kMAwGrWAgCaWTrPXIJARQ+RrVq1chzKPYcJQGpHWfpj1KbfB2EFAFhGKXeZyl3O003CLGthtVq2AFIr+UJeK5fDfO5KlDTiOJ07Hqa2ss1DqTTLBYCD5ahoc5w1vW0OBg/j5XHg2jxKlRsNc5vHcoNGbrCI6f6Hi3TRUsvLddKy51HKTsfYBCmPx1ynlbS5lbLZ1y3LXaVypss2t8ONbzuXtrKvbZs5vgMpW/u6bPPYTBK2eS6T43CROqwzr+a5HB1nLVfKY18btRTf37ZsAJhJubNFXqocX51XPC7Z5liUC+Q5pm0epXOGw1TesLHyNijaPF9wTrt5ZdocWDbXkutjHgej3L4NmdMTeyMAU5mU+/N00WKWxzvMQ1F249ewcdXWreMx13DDtZQv4jpmmVzHfg2HlV6ia4lzeun6eCRt3hzM9Bo+u6JcfCCVZZtXps2NCCO0XMfLcn6t2rKvbdsbafNokI58fqykw2ZmLc3nqe2+j0Pf86Nh2WV5OGI+B7eWplxLc7eWcjeZOS0HaetSymObm3GevxvDcl41crHOq0XuqMWUZZfjHPxzy66lketveX5sDNPFm/IQaswkmcvi35PJsfsHTz8XY7xnGcw/8y9sxeevdtf47dhvv3P6CzHGz7qjN71DdqovNrQY4/UQwi8jJVF7l/vujQDeCACbH/lIfO0//lJ8/EuS0OarNq7qeXw4vXcvzZV3P5OERfceT4KWW4+nCbz5dF6do9008bg4p1fSPfYfkYfNK9KifPDlOZXJqx9IZT44Sqr8O4skxfFH1x4GAFx9Mkm/bDye/2NvP5HK3HxWXhCm8lK0kQree0l+uOy+Sh44H3YIAPjIVzwDAHjdQ08AAF46vAkAuLHc0GveffMRAMC7nkrHxQe2UpufSO3ZeDaVNzjMbefiPHw4nbP38vTd8tFU7qtelvv2wy89n+43SE/YZ6apT//Dcy8BAFz/YFJE335fnk7bj6cyN5+RRT5Ni2qxJQ+Xl6b+2X1lfmDsP5bOfUj6+7Er19I95B/g1WkSjn3f8w/pNYcfTCrw2x9I7dj+YCp3cjX1dTAvGfNL0t8vTcfdV6XPF4+VbWZ7AWBDHlIfOizbvPeEzKsPtEW5ADC+ltrKh/KM5b5Mjo/mOq1elebRYy9JZb9qO7WZD+MP7qf59MfPPKzXRBnf7Q+kvtt+Uv6JXE115cN4djnPwd1HyjbPH01j+fKXXk/HbZ+uB/jQXhrXDz6TUmq1j09SeSIHt/V0bvPwpvyjHshL70NS3itlfj2W6vjwY9f0mte9OAnZvmbzmaLcPz5Iyv+/82wSq3728Qf1u43H0/zhOh7tlC9U08vmhf+l/ev4wx9M4/vwOAkS78wnuewbqZ+ffjK1efKB9I/o0gdkDT+T7jXYz/8klmN5brwk9feOrOGDD0v/6B57VZLG+cSHP6DXPDpJ4z2VN4J376a1+zsfSm3ef1/OMrD9/nT/rQ/JS8PNVDadhOkDXMP5BW3v0XTO5ivSbv+rH0rlvXiS/j6QN4PHd3K6tKeeTr8PnkzPtM2n0v0nV1O5rXnhWIz5/EjH/UfkBeCVaS09KvP51Zfyc4Qvc88epjX7h8+nZ/XeB8u1tPVknlcTWUt0BrmGdx9Jbd19TNb7YzkDw8e/LP1/+FPbKasKX3Def5jG9u3PZhH0p9+fniXbf5Lm1aUPyEvKs/LcWqS/Z5fMWnpFqsPOh6W/w4elefTxr0jz+eMup+N2e6jXfHCa5vBvPZ8W4Fs//R/ZVCB33Z6/usRv/sKr7ug920f+8CRpe+6pnWZU1IsBzOWlZgPApyOJdh1xTfI8uEDGIXsDc6TJRiSF3kGGJuQeq/xPhZ5PI+cQwdHbymJaGDdwIW7/XF6kVnJtw+cpi7Oeo3o1pVvZh1rEgdRpkAoftqWnxWPTQ0ei9xdWPKbP6dVaJIH1I6JCx2o+L5EJID8EaURq2C9hyXblc9RZuQVryvYTEay2kX9a0na+XEzkOBgYz0P6i32oyIQcQ5P/0Wl/0+OiJy99vSUe16VhfiDR49wcpIf9RBCO3aErd2TKkYd+XITiO57Lcm1btOxBKptI0eYgvcAODYJ2qGUTCSqPIZZ/AxmtYn8H6eMmrB8gzu3o5hPXCecOkOdWXK3kuxLdIQpwOM/zak88551lerFopQDOt+XqiAwAfp3x7+wjKEKGYbrveFgiJ5xXi1W+iHMO7r5arPZBt+35KBctSoRxahb6XJ8j6biShvC5ZdeFIo3uOaLoJ8fftH1d8gSiJVy7c4t+LqUOigyVbUYRaBKKgz5PGulrgZfYxwAwCETKZI74uafPaFOK1kEQlGWJmPEYY7fBnE98ZvryUyXK/u7chuVE+3+jPK6k37RP3ZhaO2q93U2LAFb2AX3O7TQRm0cA/HBIWXobAD8eY/y5U6xPtWrVqlWrVu2M22lGRb0TwJ85yTUBEYNmpV700Ozfrpx3l9GL9Hfj3rIB42XS26Cnol6ZeEIGsVm5V3q+gbf0AlruyZq6kFtBj3pZ/m33b1dt6UkTgWpD+ba9Mi4ZPc6OZ70o22k9bHW4eM6cXqZshS1yAw6Fd0BOAjk22ud9Toj3qNuyrfQubT+x7+g1bwn2zS0wevJETQBgR7zx1SDdiPwGbg9Yz3odcjJxHIgrgwO9hqjgniA23L8nCpARm9wMokahdeUq38Rsh0rZ29LGy4LYcLy3h2n7bWzafCDXa1t53zGJD2xfF61SRFD6mojgSMbWzm/1Lj1q0Yd+CmQfZHwbnXvyvcwvIoJA5jrsLlPftlLxQxlnoglhkSugy0A9afmTHrdB6KK0MciR3JptGeftdtpp83ggW5iC4qkn7z1668FLP/hnDJ8fXEsW+dyXNmfvXsaHBRjPXstWrlYoPicyVaI8Uqem30vnM2OxNBdJPzdLjiFRErYzn5oRjLK+waGtlgfJ5zbnmj4zPUpiH1NEavgMW7o6LV2/IT8biWqT66L/N1rzT4DPa/e86iA3sed3hxbNZCBIi1j2wGanhdgAEct4cRCbqjxcrVq1atWqVTs3dl+Qh6tVq1atWrVqd8cSx+biyMSdqRebgET80u2fnoEipLtaR3hcWPhcfudWlMKcJMMS5uzWZSg3XnCrSGBVhbAHJvLFkzu5VdC3HaMk5LJQwptKOjRMQYYELwVWliAitCQG8zizWLKEDLuQTra53H5riqPewhGmLf6nkO6A18q9HGHaEmkZmkoiLbeguGVwMJCtqIEJieaW0KifRNyYOi1d2QwDHrRCHpYtik0TbzqWiIoN2RbzWxWsvyXqKnFZtxzLo50bY5KHpexLLcnD6f6XpQ82RpmEea2z/SYkZW6/rcqtKttm1pfhuAy5Zl9bIi3bym2r7rZMD0NVSZ0l2ZPbNEuzxXkgJPS9xRjWdueyTcPwY7MVpVtcbnsk9GHP3CVpSyL+RlP2td2KYn80bh3rGj0C4/ZbUtza4ZbanhmQfdmXXOk2hmwn656I6Vu/TdLZmpK/7dNcxoxE2cGareyCdMtnpiMPN24bKLXV70/KuQ230Lvh3oxO4metEnf7tzxTnUgW5pYUirros9TUR/8HyA1b3ZIq+yR9mA6ZXO/oARzvYnuMRylH/l+oHIUMiH1eapDJKb5cXCTycN2KqlatWrVq1aqdGztTiA1CQkbanjdPIhokthK9IPqSvY98jZKH5W09OEKiFms8CXo+JNKOxSsgWa5h+HFP2Ck/azqeb9lG+yG9DiI19PT2jRoXRcVWDLFd0PNyiI0N926dJ6QqY3LoCZ8kiuAJ07yXbfPShSLH4EPcy5BrAGgHJYpAIi0964N2WHwPAI2iF9IuRcccuRHdsFnWe+BEF4eGYe7D7EnqDd7L7EGr1Dlz480QU1s2vVgiRPQy+ffIEh6JIhCFkbYraXlV/g30iMYJGXl7WKJiyzY3hP3McdFwYxLBC4RgTXxxLI/WsyZK4eUEiOQsVLjNIDYU5nNzmuiRJTR757jVtUsyqYhftnniekmBpUPbdD4Pcj91kCsXtsxnEWUSbJvp1ZMwrc8t+4jzUdEM89Z57Oc8dI4QjRy1pRAjkYNinbt66zOT6EifuOkaWQc+I+xaIhGfoeCc+30SGVoO67R0z6s1IddARk6W7ob8v2EDMYjmRYeeK3LTh0ryWiWyyzOagpYulB8AlvIQWK2Lw7/LFhGxPKW8kKdhFbGpVq1atWrVqp0bO1uIjTP7RszfZ0QvFLFJ33sRLSC//ecQTldA6IZMjpzwFOX+6YU04oVY7ym6UGcvKlaYQ0woO04Pj16H5SWomN6iFNjSUFvlDpm2O086f5EOrWmzolNypFdLmXYvkmd/X47L8F/Pc7GIzUR4H5uO68Ij+9yGaypy4lGRHi+wB4Q6sfUJgaUv1vx+QmtdiGqfmFlw3IR13mYp/MijKKsKZ+eSCzO3HiVD0CkOyDBzChwWqUCIYETWhXO9G75MW4pnTeRG+SbCSVnKfLaUqkb5Y/L3rLxv0xMaHpUrR85FicyNTUqHjNgIh45zm1yxHomGHCpczg2P2MwMpMk0KESYNSXLonxu2fvkD8rydNzN0zwoUlOG8480PUoP38I9E5TH5GQj0ncuxN1xbvoQG3JsVGCV69iFXPc+F2PJqfEh6CUSmDqC/xMYds069a8lFMf8POH8zVXxz07OrzyfU/mHZgEu3Bw/Davk4WrVqlWrVq3aubAIYFlfbO5v69vDJBO9EyHkRbOsqBj3a91beRaVSwfLoqdHRxSB3t9EEJtWvA+re7V2W9Xtadvfl4vSkyNCQ2/D8hI0QSL5RGsiU3rDu3wQBqMaWtvm1DaNnJEOYtRMFqszkWCKyJR9m9MASLmGl9O2JX+JAoxDcRnb0EUvOua8qb6+zQJzEs2gkveSTmJVcj6AzGlSdEHRMXK4rDeL4jPP84LxLln23EW9sQbkCxzl6ZEHECmSxs/7eD8uWobRXuQxWVExRmoxEeu+IjbSDMvh0Ygs+U65KKxjOgbDL/LjyDbS840qFJfPaXyk35wDXKKV6TupAyXvVTjNcy9MlBrnnqAJ+y51Rk502d3FXyvmJ7Y0kYZTTW1QJoTkc6tAnjqpDdaUa8oj149rlGiniuQxwWyPgJ8XYNR5vDDcFEXESZ6Sg0NuLJ9lrVCeEzU9EmXtIDccY7Om1qQ2yOlvLGIjt3V14PxtPHJj6uSRGy131f3/5FPwVLu7diZfbKpVq1atWrVqx7e6FXW/WkxeDxNRWs+ae9ZEOKJ6yS5aomdsO3La7i1+ZJIujpVvwr1qOYccmx6ZcP+SnvUu6JEaL1Y8zugS6BGhIZrBNAeAkUV3ieE6wIZh+OfIChfhIp6e1YrxejL0MplegJFJ0XJsfHQSo7y8Hojx4L0vs3Ietk/gB+T9bR8l0Zf4k14wKRWU+Wf/eYl/IHtdRMyUz8QINHI+TObjdlZGcPBvjVab55ZOF+X9mRCS3IQDgUeKZIWLcm57TQ9aL7AlRROx8TwTq71BNIcpLMjP8fo5wHoNHc+lskjg0EWjrXwEimpRWfSiRGzaWdnodpbrr/o3i5Kv1qczQtPIP02RAmm7VMlF3wH5edFBIY8IzeAaIgJIXhEjG22e18ZxXHIko7a0Ux6TaTIKimNJ5GbktGR666tcm1gcgdz/Hk3ieuxDJnyEIevmI5NWlr/E3xn1tiZyKxqIvJPoU7k268OudAx9VJRwq3q1xrSS6bAkL5Io3Mo8o7mTsDpiUlS7Y3a2XmyqVatWrVq1aieyCFyocO8z9WITkd7CicpM21x9vh3Ts/UJIfuRmv5oJX07d0kZgez5bAp6wWR2RG76ZA+yd1F6OZqkze6nOy+TCfRyIkrR+DARFvQUVIvGI1CM6DDObUZU5G/xqIOgL1YrxuvJ0AMaK6+oTBgImP1ydsg6Xoz5mF7ewu1Rk9/CMT4weiCq3SMoSDt3KIlJ/En0iLooPM5n6b5UvL252NBriJDtCKKyP011CR3EJpfTTkk4kL9H/BxFuQAw65Q9kXLTxbuCIh3MzVJdlGX3Ra0AGdEBupobqvXhvFgbLZO1e+S+1CvSOZOvU74NvVflosjHws9hhBWQ5xhRBNra5Ijo0TVhtCO1lcytMjIn6KeLMLRaUNoOh+KQt7Ry3AuLVtGWHf6Y1K3hs8HMRbi5TsSBCHPxTHDPDWq3HPFsC00ZFUWEedKUfT6w+khEL/ySdYkoU1347OLfUn2ZV4pQ9KBinFec40H5LSXPJf1eolFrk1OaZ9tcOXNpnInsqyJxH89lXaTZERybHE1bcsMWjj8FZI7N8lSjoi6OVVysWrVq1apVq3Zu7EwhNtWqVatWrVq1k1lErOHe96vFGDBbtgolH5ikclOFPhlX7Bm7cg8LKbqEZ357BoMyeR5gQhWdxD63lbiN2RdmrITWuTu3SPMgvwgkTVh14UMXDU6fw3tLKHnl2mVxXIqsrRzJsx10t6I2mQBSYGwvtMUQXt/lha0h/VnIfeG23Uik5ZYUt4iUwAsAQhZl3koeuR1UbkXJdpVuU6XP5zPZ8pqPpJxMHmayU24VzZiYccZ7oTgCQCNla+gzRQrdFhgALOYM52fZqY0kVt6cpb+n87zvw7I9cVm3w4jem22Z4ObVwomJsY9t6DO3qTQ8ltV2REvArBnC/Uw0ym0Z2eJk6DiQBQBJTs9EfEcqbdY/kDMRX/rcEly5TTInebgkavdtRfE5sk6I0acxAMzWhCMPM4Em0zNY2QifDFHL6wtxl999UtvgUmhY4w6wT679ZqEAACAASURBVEbpxS7HdiuqLcd5Xah1KtNtj2nCT9lC7wl59tt8rFuzRmgSMM9oleQoK5W3y/Lnc5eM0suDHLUVled2P03BnrNOxoPPZkse5jqr5OF7Y2fqxaZatWrVqlWrdkKLwPLiADZn68UmIoXmzl5AOJ16VdbL9AkZncdFUhu9Z8AgNc7jWpd8E+iGa1JUjF6a9c5y2HIozvFehvX4lGxJz8cn7hs5dwTdJJUr510yfB3IyevY9q63KXXuRalIMnRkQw29tuiFeNSCjuwIkXbaCILDkOuZQS/WkHhVwM2EA7dDllkSjhlaTySIYd9AlqInYXm5KK/VpIyGWNlSNI5tn5V1LNssoeaCyLCNDIEmkjOb5YmbkZ/092AqxNlpOS6DiamTXEO050D6cFdD3Mswc8CgOVxTTpDMWvQebwf9FBmBoZURSA0gOZ1lEy1UGYGCTCq/eG/ZzTfAzrV0nM85vmx7RuZoTEZJRCt0tBqkHuYxQxTHpwthm4ea1sCsKVlfRCejZ6SataTPDyI1Lj1MJyEl0NHi9GH9mnTVrPOuUN4adASGPOxEKLmW+kTqlm7QNOCiKcfZIoH6LOtDTgBDHs73XipiU/6fUNTkCAXAjvRHD0KXQ8LlA5KuNfihFIQEcjj/fFkRm3thZ+rFplq1atWqVat2Mou4WFFRZ+rFJsaAxarRN2/7RrwuHfxRb96eW5PPoSdRhpJao/dBj0Sl713yTaDrYdHzUjF4m5ySnocrsnGCWxZR6Saj5BFF+6xlobGyzZrWoG/jntcyVJV71tJ2q84enJdJJIXj4XkuADAX9GKPXJd5QhHIByCSwhD4VE6JenlULEvuA8ElGNQxWjJcswzBBXK/L11qgxw+LfcyqTqggnmxONcfAaiQJLlUKtTWkBfC5Ii5Tq0TGtQ+nppJh1KsTkPNp8KtmZdh5jfI7TETd1f72yVZ9RL/yN689oJ6vtIHPdIJ5NZQRoB8te1h+pzoodGiNKHmjnvR5whraLiMrxz7hBhp/I58L8oudBIfmiXVTUYpJwlfjetzYtesTFgNeXbPmN5UIEs3tzkPegQao0Ou2bccX/LkLH+wGZTIiX9m2javE8rz3ECLqivHxT2rSZvRR59FRzQ1hwjlKecGa82vY+Xa9Al8+kd7KI/KU+zl/ZScJD+GNrTbPyvvvYUOYnaereJi1apVq1atWrVzY2cKsQESh0XFxWyUDz1pvjUr077cN7boBaMKFLlpyzfxPrE9lfd3aR24N6/JN4tkfDyWERyxKRMHFua8AHp4XswMMFEkTEbpvNrl0LvTPWx/11brVa0cOsU2E01gQkgTSKWITOu4Lm3jopgM3yQ6ITXyDzi2/HxhEZs1qRRyCg2DhsXyeBLjOGhtnWcX7WRRUUJ+d/xycsI+Qd3k4likFUjHxqUXCPNyIjVGNFDnnvJNSj7T9fkmgJJPRjSHfBxGF2V0zJQlv9OzLRAsQDvdRgb5SB2KA5J7Mxql7/dGJlJLkk8qh2foBPWans6Wy1cuSaFNJEsjl4p8L8/T0fnVN4ccj48ICNOTsF1AXscHKpS3HiH1ySg98ht6xmNBvgfTCci8YtJeIjeW99OoAKMcPcrd90D0a0oQG6aCmRnEhoJ57H//zO4IpMJGJ8l3RJiPiExSrksso6OmPTzMuAYhX5dmp6hfT1oY2y4bBaa8m1NCbCJUZ/FCWEVsqlWrVq1atWp33EIIPxRCeCaE8K41339aCOFGCOEd8vMtd6LcM4XYxBiwWjW9OgQqqU7Peg1Ss2qNF+v3T9ehF3avVBn3glrQG2B0A6NmjD5L41IpKC+jZ7/ecxO4B879+Y3GCKaIMYqEGjQdPY2eiIJ10uT0amcm6eL+stT7oA4IvduoWjKmzUQIGLEzI0pVev1Ww4ccDialpLfHsaUHVuyLr/O06F326BZ1ohoY/dZShj5XaqJ6H5I+QvgSC3q1jp8FZBRBkTmXaDSaOch5ysg2pnBgm3OqjnzNOgQo80BKjg9gdYOkjo7PxOirgUFsdoTjxIgsH43V9iBCGg2n/I/1sFVGDxbF56qxIn29O8zlkDfGlAY++WbJhWDjy3KXTsPHru/MaXLRby51QMkvkmcNU7Twcxlbzp0NI3bEhLIHbWqQcl2clow1RRwX1Ozhs4afG5R1UUboeF6LIjaWY9O654fO1xLZTm11dXPpHRbLLipGpDfryXA990cmAXkM8/O7fJb1cW060UmrEj2aW9TEIzY6eDw6Lpctm+hkU65hrtlVz47C6kixr7trp8SxeROA7wPw5iPO+ZUY41+8k4WeqRebatWqVatWrdrJLOJ0XmxijG8NIbz6Xpd75l5srLfeGi+We8X0fBqvydDDsVkt1ntHABDFC5kbN5BIDVEMeiT09CKjKGxUlHozJe/D7+MCxuNkFIl4rfRiqdZqI6nGsodPj2vJtpO9f5T3p1EN4gFT08V4WoweIb9oj22nGu6i9BwBm4xSvEtGciiaQQ//1ouNHhAjtRqzp71cM748NsFEzpFTtS7xp6jiXpI+BrK3vS+ROkziOBPeB8VrbVJEogiMeFqOyXWSuhnOCLVaNqXsTeFhcG6rpotBedSL9fwDGe9ezpaOd8nZIU+NfKbGzHWOLyOyWiYa7VFbptIzvXqfoDOqMq1NDGjgFWSODRVpmaBRo4zQRWyo0cR2aTJOdBEyj355jx4wUWnkjWl0n7RTuUqmTm1/v7Mm5BVZfpxXAB45xKZEVxmiI/UnjdBHWxrEpqMnQ36JT3pqIE+uq4WbXxmVNPdX3o1rMyOSBE2aLi3HZlAcZ6ryzIvZsD50JLhj+X0f6YmIEBNQcmwXZg5yjQaHouv3Pc/onOBXrnHzq8/IcTq9qKi7Yi8KIfyW+fuNMcY3voD7/NkQwu8CeBLA18UY3327FTtzLzbVqlWrVq1atZNZbyqJ27PnYoyfdJv3eDuAx2KMuyGEzwHwMwA+8nYrdq5eH6tVq1atWrVqZ8NijDdjjLvy+88DGIYQXnS79z1ziE0IGdq1InKNYKKUwG8E2l3pllSXQNsbvogMSxKut/A5t55IsjxwEuxY9sDSHqHswKoW4pVfKEGvoaJpK2SzhzzMbTiK6y08PKzbED1kUibdk7+5FUViMJDJwoR2dyQMmEkFGQbcl8wzi9SRnXzrmMNGxpft4pFbbtwOAoCZbCctXTLPlW4H5fvqltC43BIaSFjx1jD17QPDfb2GIckHw3TjjVHaOtgdybbfSODtce70wUTKli0CLU/04JZmK2rUKfug6IutoaQ6GOXOnXIbzKcEGZZE2r555aUMogtNtcPDOc1tjeAFAW26iimlBrhdVZJvMS/TVgA5ie1hTH3r05Q0TrrBtsNLGvDSYkvQzQmGVA88ydPsp1DcDS4IoHEJKItQeg1xlyPrKkc+pyZGjXJTyMObss3LLcfgAgAA88xqdNCkPKlLD6GZDwFukTPh59yFf9vtk0ZTG8j8cWkx7DZ+n2ifbTS3OGd2a3PJOohcB+eXpuzoPo+9GN7KbT13Qq/RFcrj3NbQd7v9zef1spxH+Wash60Ty5ZxYKBHU86vPjst8vBpcWxuZSGElwF4OsYYQwifgjSSz9/ufc/ci021atWqVatW7fgWETr8qnthIYQfBfBpSHycJwC8HsAQAGKMbwDw+QC+OoSwAHAA4Ati7OhBn9jO1otNiAghqgc0bhadU+jVt/I2veyQzHruy2508uD0VGdHIDYk2Wpys16mLuvAkEXxVHyiQGTvkoRWkleJ1PBoiZfeQ/BCdDlhnWmyS3lAIIjCbUoMRiYL07P1YoShF6Vyc5Nhk4qclSTA1JB+CfrLglYxXNMmUtylKOGwRE40DNiQrPU7evJjIeiOUnmXh4dS3qFeQ8+aHu/WKHXUs4LYKHnYqPMrsXVQ/s1zokFsJoIAsewHDVoEADujSXEeABxKvYkAsV0DF/ps0YtMnKYcQuldcj737cPnJKdE5kr0AsgCjDxHBRqFhx0EwTmYdefVvnQMZf8XjlRs16wPtSUhm3N9ZTT3NG2IzBEifRQhVDK6cdM19Ng9C9al7ACAZtgW3/Ea77fbchhuvSkMbMoKULJhZQjTnaSUoUQXOmkNAA1jnjmhvGkn5Dp3bnDoREcuwiI2jVvHbtrkkGsbeFEK5WXy8BGkd/fs7JCJieYZVK+L2LAuPSJ5q7JsnwbF1yOVxTqkI8O8256dBG+3/y/7bFmM8Qtv8f33IYWD31E7Wy821apVq1atWrUT210gD9+3dqZebAJSSOKoR/CKg+b5Jj5RWWG6V02eifxNkoocbdJFhsWOW5H35z7xMUSmvOeTuRHmXPHUVI5d2kp0qunxMr34kyIo4mUyytTyAmgquqZciJ42C9+GZaro1xGhix2PzoWMMprctj240GfyTS4P0pF78wzLBgAoctIffm292BwiLHUUfs5Y7kfJeyZlBIDNJkEOFFQjgtYoUlTyEexnOQRZPncIQipbUClBbK60qa0c52uDlOpg0yA21xy/RxGhiUNsDO9npXwTQVIEvWBfcp7Zh19O0Cgf+HBWK1K3JDHNc1LSx5xnTL4JZK7WjWWZgJOoqCbf7Em26bkPrPbShntzbg0onSDo56AMtbbJdNumv42K3DiuDZC5Ro3wiPj8YCjxokdFjm0dN2X/UwCyQHG94KYiN2xoF2UI7tnFPiXyqOHHvUqALEfK75HK6CT+dIi4T2tgf1ehPIbWSz8pSNWDavhnqOfaWMTGjyHbmJ+TFrGRavvEuMepixwZJk+khojg8gK9SNxvdqZebKpVq1atWrVqJ7P7lTx8t+zMvdi0zUojnyzHJiM2EhVF3ska1j5gRZlKpEZvS4/LvOFnWe6SB+CTb1qelno8ygcgakEuhvGsKfbkuA+t27e1HBtFUJTzIlVw3qWNYmEdWvUy17d5sca7y8nrGIWTv9PoFYnUUY9aeCCKIBi+SStIBiOEiNRcUcQm3XRzmJG61nFdcuJP6YM+YS0X1UAPj/NpbKJX6Fl7CfrgxrmMmijjYrreZT534CLAiBAxWoVIkUWpNNkpxeoYbTUuUSobqcVzlFc0TuVdEqSI6NjUQAUsUyN1jko94iLvuuJxUqd5niRMsrmzSDwiCvMRyVHUcGkKctL9WnwPZ4voFKS/yM0iYsNIw7kZPIoCep5l5rF0OTZE6xSdYtoKpghZk9YAyJFgRAIH0teHRTtKhDdHJLn1aB8Rq3IdK2JDvhyjlqxY3Tpxurb7nPJoTkbOGHUXi3KA/MxU3o8ciWwperKeopKRcLfughGw7KQnkT6ODtlOZZVrJrg69Gnu5TrI80N2BygOm7lbuRxF2E8NxQmK0l0EuzgtrVatWrVq1aqdeztTiE0IEYN2qZ611YXwugz6Xuz2Sm0qAnqTncSMPs3AMr//rSNgNe4tvi/SySM2fQkU+arpJXbYPqIW1rPWKCXxEIeqN5KOmQNgvMxQRg2p5yJ185yhVDXZS6ZHNCgl722bcySQ9K1EgunnynfJdaKmC9NGEK3YFs4L+QETk6RStTc854JH6/3xnDVOU96Lz+O9DKW3vfCaG3178e6z0HeOs9bF0PDvnNzVoG2arNDrjUgB9LBN7lLOwSj9TT2eKyMiNikay86rrcE2AGAwLHkfmfNhPHimFfDBcM4DXhr0hdytm5Jsk9wEIjnk44R5Ho/GcSE0kazj3gCGA0G+2qCMsmP0m43CGvtklA6dypGGZjxU88Zz9Mq0AoeGUMbniOfM6Xy2aSRcmhBFThyPzRrryWcX0ZGpJv4s+S6A4Z6s4zEV6Ce5c/xbvtDnV3eyE/llfzPtQtaSkWuP0ADrrG8p16YcuZWeTPFsW5bzNtdhjaBNz0fk2KjelhwtEhh8qOo9tghgdYFwjIvT0mrVqlWrVq3aubczhtgAw3al+9GWCzEXz4MekL4XR+cNmCgZRWz45k30Qvd8j0AvpBxyLrg3H9SLNjoUzqPWRIQ+osDUhUYdCHpa1PzYN+44VYJXwl9ghBMT9lF3JFjtDZcw0TsSVpGUnoiPIhkLYrMzLCOSgKwVQ2RGo6GI1CiKZTwtl/CTiNxEJG/Jc2mK8I9+D8gn+UyFuc/k6BP22b4l4Y7cBJ6zEhRhoLyKLhKofasaKLxpHmQiQFlfpFySqx4EbV0Ennru1HQxiIo65tLfRC8uCWpxZZAQm7nx4C+PEu9mPE6N3BN+znIk6Jvh8ERyqYigeH0Tfm74DXNpOyMNqWuyPxeFWs7nmeEqMAGn47NEhW5M3xBlaUseC6MpGf1mUaqJVwB2CQ9zg/OvjaA3WeNGvhDkhmrLByaxLNWWicBmZI58Nfv8SEdFbIalP+ojw2z9iMIQmTlwHBurqr5co77ruTaA4Q26pJTKO5F29CWG5JxfOrX2xnFtgD7kROqiaswl0gX0cF0cacdyiTxS07ijn7+FKXDmkGxebIo9So34XtlFIg9XxKZatWrVqlWrdm7sTCE21apVq1atWrWTWYwXKyrqTL3YBEQMmlUW4TP44Nyd62W6szR6vka3aEJ5zlGhfoQdhzyJhDGBP1W4rUdgi7Bto9AxCzZtdJC9wvOyPcLymZgSMIkFXTJKH3LbC+t6gqBA4Br2ikzWZXgsibQaDiwpCVYjQ0R0JGGmcFDiKU81a611yQkJITPHydKRD4EMtTflcHTaDgDt3IXlisz/VEiqO9KnNxYbes24SVsG1+dJKG9fUk1EJnmUbZJ2mssZTMv+ppAey2XSUCCnGLi5KMXquM26K1uP3AIDcnhsh2x5BC9Rk18yYV9bbvtdatK2zDyYrSgmXhWi8c5Y5va4JIIDwGIibVzG4jslyPcIZWapexFqQxmCvJql48AsbmpyMmUDpQy4lo4innI7gGkMmJ6EJHUgh4L71AY6b3uIutxy8NsYXZG8vMXJrUcVyuPzSmUjTDNc4k8rBmm/LxP8llXMW55C3OUW97JLHl6XELIQHe2Qh0lcL7eG+raiVChP01eU89kGeKyd27o9Jte0ZiuqRw4E6Cc0++1pTalAWkIov7cW3XabT/VjAwL43Zq8y/fEesUYz6ldnFe4atWqVatWrdq5t7OF2ITkafaSSMXoiWiCRoruOW8KsORheth0fXy5xhuXMunpDsRNIumwEU94Zcl/mrANxVGth+BK74kelUVo0t/Z+/NCZs06T964Cz48WvUGpd4bJtHkpZ7kkEAWyjsKpdK2R4bWl+UX4blOhj2Htpdy8NOF8TI9iVeTL8pxZseB30k/yXE2FdRklsKOr4839RqKxt2QkOS9aapDmDZyLxTlAUB7yISQcg/pn1ZU15rDPAGms4GUnZCaa6NUNr0+itftm6SkIFpEkjjbrOgFivLT7xwHFEav0gsRAhmhI9GYmvcEHqwA4FB+J5rkSeL07G1Ybus86ExoFSIy165BuHR8dZzparNPjLdPAuiq9FQZWk1UrEihIROJIe4LhzBmkmyP9+tDhhWx4RrOa3ZnmcaVYfwMDVZieE/Y+sqlYlHAo0c24lZJKWeK2Bj0U8LT23XIdU+d9Cj9Q9L1kElWi8ycpWkKGC+SZ5/RmuKgrIQ+v1zINdBD4uW9ehAbH2TQqUPfWlrz/4Hl6RoqhDgpHHt64d6nkd37tOzitLRatWrVqlWrdu7tbCE2wrHpMxVQExdi6dILBA3FtB4d94VLkSZadJLZQEZmGDI6EDdJk9iJR7q03pOXoFcvQepkwn/pZMwFiegLFQVyiCyQQzbp+azdjzZeJkNG1duT25MvYyX8idTQs6VnndssXBiLUql4mKuKD8M2TpWXf98Vr3YeSvRib5493yg8DKIvRE4GhxQltIgN5BxBAoSvMZd77MwSKnZ9ljk29LSI5kwF3cnlyXm9iI3UZVjWzeRuxXyaxvEG0SLh8rDca1Ph9kxzm7tokRxnJWLT2jBpDUEWNEzmDFExenOWF6BoDl148k4Yqm/C+xfCu/EJOFdM5SCpMyjqCABjQYKY1Hbl5rjnXgCGQ9VBbNhOI+Y3L1GkhUvaOpJ5ZVEqIrFDQSxnTK7K9UK0xKYXWJPaQFNbCKJq12x3PZNz013DHXR1UKIWPkkmAEVeNbWBez4y/HthUyosOIZl/ftEKHMIeFk3htZrqhDD1Rt0wq7L8O5OuTC8wD75BkAfKBYJ8UJ5FCPsJHUtKuPrUHLCLD/Ri5kGx7vkfLKimiy7ObWw70oerlatWrVq1aqdE4u4WMrDZ/rFxiaCpFz5jAJqgth0IhXsm/eSvA8XZaDiT+kwbLseHSMrmNht5NCL3pfjNUhNKCK1xBsjYiMRO/siwkdekZVC7wgIeg+Pgl7W4xqQC8FzpE5MRDnIsIJPbcB+Z5vpLS17PKFOUjmNTpOj4U/MhG9CSf1r85JvQiRlz6IXwjchCkLkhKhJYxJ/MiotR9YIB2JW9vWu4TPR+6P8/1K870b5H10ksJ1SoEu+G7OOXcE5Rv4cCIeGXB56eJoQcprHuxHEpnH8nvagRC8Wk1yOojfS30QC9xh1tSrnF9Ddk2dEVdS0IRYBJE+tnFc6vySiimJ/QE52aucakJOS9vIbXHQjxzejVHkc2M9cS+SrKUolC2QYMp+Mc41rXtOFMOGoJq7tLvAOOsm+kIVxaFAaRruRZ6SIzbK7VvX+bl1nBKc8AplPxrVJXlGORJMIQ5MuRjl6XJsuSumohJCxKcvzInlAjlbqcF2Ua1OWCxjEPZZ18Lwfi4QQJcqITZC6iIhqb0PkdvpsLv839Al9+jHKEbPCszSkJBULPCWOzUWzM/1iU61atWrVqlW7tS29A3yO7Uy+2Gj0hHFRKIt+qIiNeOee5d63xek9Ie4bD+h95IuI2IxDiVZknYKjvIF0VBl48ULo/QOGC7FkhAijGCg7L95A3yTllrXqXoTiaLniyoHQ6JWSAzEx4iFMbaDRaM5V0YAFE33S8awF0WiZmHDWjXhZCHqxS67LnJoukhxRkBTqzgBAcHwZRVDEkw8Lk0ZiXkrfa7JCqTf5BpaL0UC4KPKZJtU8ghcA1cKQOmii0Z5zl4ycIQ8jee6DVTqZyMpintFJDo3yiYRnMjgUlEFIBAMTfcXAH0ZkkbPDBJRXF1sAysSy5IFwDvroojJKpuRUkYuieR8FCZyYaDsmO73kou2Ul9GjfROcZ83x5bwqEBs2ZQ1fzfOLgIxYDV2iUZ/IlmklbJ189CMtytyZr/IYsg5EbBjpl7mBZi05dCJHMpZ8H5vGhRwbakN5jSiP3KTCyWkqkeQO1+aIOhGxzukr8nzShJ/uGal/en4LSh0qwHJuHK/J/O71ZNhWziuLmnTAbo1sk18Uwcljt05bR+eO4/bYuvRFZlW783YmX2yqVatWrVq1aseziHChwr3ri021atWqVat2zm1Vo6LuT4sIWMWQhdtidyuKcG90kO66LLGACZv00uQC59qQxQ3Z81DInttLHl41W17Nmm0ZbmsY9L8j7kWhPh+qasvj7yQMciuNXEVuO9lw0KXbKqB0/FhCekdGjtwLIWroqNRJidpGwVzTFnALalr2D/9m2DSQw64pRkcS8Uz6n1nM7bZM0xFg5LZPuR0EGMj7BGgw265icm5XrxPKny6SL0P59zFMtwiYRoJbYFYSwGULZ+hzQ9Jy31YUQ8RJmBYyMgnZ18YlURvIW4EqDjgnIb8M07Xmw4CVVCxbUZtGRuAByR7uhR81VcewzLANmP7248t50DMHmcJCyeFC3N2XWHQ7v/nw96kNVrqmyi1cIK91H3bttfbsGubWMsebxGaGhpchz2x8eV8vrlmEe8t6ZmqUwZoUB0XgAbeCOuRh1qMn3pvGLeamfGb6tAbA+m17H2gAlNtS6W+W74ovnodl2DUHIBN4bQGucm6rk2sp9KV5UNHL8iZM9dOaQaR8Q3sfZPm+CHamXmyqVatWrVq1aieziIulPHy2Xmxi8mCVKLzK4ZNMMNchOh7hnedEbo37W04Q9IJCYkD2fOgNzIN1k4xMeE/Iok/MqIiNIcgpkVi+IzHXJ1qziMqIomdrBNRIeCwSENLzdN5eH/l55VIcsN81RFWEvYZFskKp7yHDryXUMhI5YBJJUxDDcuckgpciZjq2RTY+d+SfPKXNi1nHV488R7xM6cdNE348YkI9mQOtnKPJEYfdhJD8LDQleZvNiXbVadJRjmspA8/xOJao2KL0BoNJPZETf6bjzBG1r86I2OTymeaBiTrDvByz9gixy/yFHKSdlpRO0vBDgz25NNWfKTyY1oCoD2Dnq0PDSPKc27XENgsqIsjTjhCmmdbACvR1UxuUKKiGVg+7c1DJ+lpHqRPvYeathlvLySQ2a1oDE1DQSXbKYp0sheEma2oDzmk+L7oieeYaF3ZN5MaT33vr4sKuGVo9tnC0GGUi/JzuiHcCHcSmkxrHE9qxXihPk1OaZ9w6AVGOabMq06OkOpXnesRGUV7zUGKd2hrufU/sbL3YVKtWrVq1atVOZBGhhnvfrxaRwm4ZKmn3b6erMixX94t9dKB5Q6fXl8OiURzpZQ6NF+sTcFLkq5vKwYZwC2rhQpEZptn07qeXk5BIkYYwGs+YaAI5CSuVgWe7XGwmjHfnQ1NdsjwgIzRss3IUJLXBalqmNQAMh0baqtL3DcO8GUZpPVOpbyyP+r2iF9aDp0ddjmWkKKFJ85BDdeVv6SemkdgaJSjiyjBzPhi2ShmB8Tj9vTdO10pXKI8JABYb0h/ibS42muKc5digCmNBiQTuIlpEL//GIKEK9MBTGyFtlqN67iV6EQpvXI4UgBRkYB2fCchpHqaHMv6HLrTeoG1E6HTtcFxLqkKREoU8te225NhsD1MY+HAkaQ0McEcUUkOcHX8pGAhCkQZp81y4WfuSjJLz2CI2fLYsV96Vl3KdlIJt41JDwstzKWxYzFuZ20xGmRN/lqlgUv3l2MNBAboyFQA0PeWEFgAAIABJREFU2ehYxpPpT/gs601N4+7vuYElYsPPHFpxBMeGc3qkCSFLcUV/b8CMoee+uOd7n/wFw66JnHAtD4x8R5+kQLovB7VE2VMdWLf+slme5W6RYzNoS0S22t2xM/ViU61atWrVqlU7udWUCvepxRgwX7YaUTA1An1EGNTTkrdofZmmM2s8raDCVvxOLuXeuEvoBuT9Ws87mQn/I0vu53rn9AFEKeQeYb18OlwCTu6RU/DKRgEwiqR1USQZgToiMsztaxNxolcLZM+WbSVHYW8m5wj3ou3h2LROKC+4RHt9badn6xPp0dNrDXqxYOSXSvgL4jRhYtBcwGLSj5wMBIW5PBbOx2hPr6HHyWi7pyYJTdgdb8k9xOPeyPOK0UgrKZupDYTSUSA2o07Z+0Vf3BwJD2SUJ9RCOVRlhM5qyIgapvXIdeo4tJT7lyP5TDZVB9M8UDhxMGO0FduZ28HfOQ4616Xayx6AgEjJJCToZyX1vSzCfZNRmlAHpr8ojLeUNBWKzJGTZogbnrPBtcnoupuSVNWiCuwHTTWwJl2J5bPoZ45DFWU9cg33ibNpFBwF4DRVR18UXImcdNFogwQSMXFCeYrYaGRSp0qdRJCKVpjqN47345vmRfKALtdFn6vrBPuQUU99TnU4Np2qaT8TOeE8U56RRU1cdKP2pataH4rENZQR5vLFwSaUzXyf0+HYxIgLlQTz4rS0WrVq1apVq3bu7WwhNkieFNGRRZvdpsU6ervbG4+GmxJdyoF8jlza83a9dEhNTuWQ3DVN5WCjGpw0uf7dw5Bfl4CT3Ist8WYHK6OtI5EmjNiZC6dGZe5ZFcu5cB4QeQjUiKFXC3Q9WyZqnFJzQ7zLYFEq5Tew7c5ld1ocQEbKhsoLSO3akqyVM0Ek6MkDwKEkV6QHT3RksEFkyCI2DjmZpO+2J+l+Dwha8uAgoyabkgCU47wtPJwwWco9RArfcmwmjIwLxd/k46wmJkXHON3v4XFCiV482oG1PYFANidX9LPrcj1RKSJRg3GpM2Mjd6Lj5cDNvT5ek+d9aBoMRrxNDZI5ZfSIzP9pec1qTi5JXrPzQnglc242ZLyZfiGODL+IiM2IfSpzkJGGNhGkQyMYhUjNmD1J0TFr8kXKNeLcdlGKei+LhtHb5/OEUVxEfGVdWq6eR2+036X/Sj2eEv3KWk1lA217yXUh/459qjwXTWDbJ0YEub8cj0ggvD4hpNOSQea8dLgua8Cxoi4OPfKRWjYyiX2Z67Ao6lJERZGjt87Fd5weW6bqF626ayeVb6OijuA23RMLRVLO824VsalWrVq1atWqnRs7W4hNDFisGuUBzM1rdocVr0iNRGv0RDNoFI7TN/F7voUX67Rc9iVMg14gdSgKj8t7PMK4Jx+gRC3KelOHYlsQm0sSQTIMGVEhssG9/LlTw+3Th8iaOqWHuJCL9+aZY0NOzUw8LfJv5qI3Q6+2SO7oPLgoejKdBH6FWmqq4IZ46oxOemBY8k6eGedwnOsSVbQSz52oyGJMXlEXSSEPJAraQ2SA3I4rBrHZaoRTI9FJWzIO7RrNoHT/sj9UT4jnGARiU9CnK8OkwktNF3qbRMu2TZuvTYhSEbERZGhL+F3zMhornUs+jni6Q6dzIiiZVceNDkVQfgM1cYyadHsgOkWC8gym5FjJPUS3yCKBXEOe1Ejvdkjv1iR3XCmnimMp661hVJZBUtx6ZnsYOclIt5lNTin1YwSVzu0jOGE6px1XD0OiJmVkEtD13Fk3VUq3NJB1HBuNfmM7zTWMTiJK0pSohfLVitA5OXjEhvy44vkhxzUaO7zv0DwIPdeFKrx81vUBCoxOWp9QVuZ1j56N59hoH1iOjVOYXofcFGhVB7Ehii9oojTEcmyU93NKiE3ExeLYnKkXm2rVqlWrVq3aye0iKQ+fWktDCI+GEH4phPCeEMK7Qwhfc1p1qVatWrVq1aqdDztNxGYB4GtjjG8PIVwC8NshhF+MMf7+URctl01OvtgDrRHyCy4hpBL7jGBb8KRhF33No4XnSSLdD2lfgYJelERHn8CWEtAIuZZbUHY7RlFx2SqYCPlvS0isFDNrDCbeCd1cQ/6zIlPEfTWBomwZRAntPZjnqcGtJ5LPGPqs6SN4W4sGcytQtwZkzJT8mT4txOqEkLs9Sm1l2LUl89rvAaCVrah8v3K7qTGS7BrWz92QjphY2mOZhExO1m0RJf+JuJhs+y17IGw/LXMKh5JUCmRy52bDrca0JdUIjM0EkZvDvBUV2GYhLstulZKIG9kCsSHoCxdqPpQwc5KhSU6fGXKvkju5poja65aUgdrnshUlc7yVJJuyk4cg82vPqO3tLJjaIKVuIMmTnqXC9oXIorTDjTfHtkhtQcI0t1ikAVkcr/v44zpezktivCbUVDn97lry67kZluk4bDoJhl/vNuOyAr3bxjy6IITl+m0yL5mwKeRhzmOWX4rVyW3WbEk1JgCgQ+JdsyVk0wpwXY2VRFwmqeyEXAMdcnIur6zbyjyj5z7s2hGZi21ATx528e+6FWb7di15eD1OoGX7Pbt7ZEwgfVHs1BCbGONTMca3y+87AN4D4BWnVZ9q1apVq1at2tm3+4JjE0J4NYA/A+BtPd99BYCvAIDBi65gtWo6SceAjNTwbTy0TqzOhXIDxkPoV0/XN/G5JRcKQkNiGP/WMFYnJV/8zvBrHiX5ZkF4FI+U5E6GctPjotczNeThnGRPrlXSXzoylQNF8+w5rYbwyhdClpyZBIokWZLQyjd/eoV9yNNS0xcIgiWe0dKJ5JHMCmQRusujhFIQqSGZl8Rtm6SShOklQ9wdKdnOlY5DtcaBsXvRzK+iqTNWZZLVPpn7zmfHcNI4hkRqRg4pssTTRojL0Qn1aTi5eLMMa7e/rzYECRwLaXkkpGUJdT8wiUevivBj0L5Nn3uPHsgEUyhiw1Qagr6IaOGhQWyuS5LNG+N0JGpFVHTVU5AGAzgScVyG4m97ThyWKG7jkJsCkRWkkuHpTOyqhGmX3qAwhzxoWgNFXfO8ZcgzibQ+7LovESRJ4Tn0uSQa92m/DRxaMXYh9VasLvr0AkQrSFK2aQUccpKT9nK9lCHXtg4M61cS7y1Qk1SOIw8raiXFG6Ro5Ui8NEVubA4bhxJ5uQ1z01wXt675LKDkyFHIzWmRh4GLxbE59RebEMI2gJ8C8LdjjDf99zHGNwJ4IwBMPuIVPUu3WrVq1apVq7bOIo5+4TpvdqovNiGEIdJLzY/EGP/Fca6xb+ZFkjF5y8+ITRnmnZEbc70fZ7428U1cvMCpQS8UsRH3qCvB3lNp9QYo4CVv9kxSaZP8iXfJcGKGco8FfumTZe8k4HQh3I1LvgnkMOhmXoan81rbz1G9L0HFQikVPxuUXjRgwnEnJeeIYdgroguGY7MhiM0lCbtmaPulJqEK+20iVoyMl8mw1qV38jxRCj1IilxEXglTR+waqGPepO9uCA+E4cCLNR59+r0MSW31nMCb6rmcW4r8xXJJ9olqqUOtkgZyLhNQ6jzL1yqXSULctyU1xIOC2FAYcHeROR/PjrZT/SVJKOemT9VR1G3lUAXXP7OpEX4UGYFr800AwHyQGkLhPEVBexBaXUvk0XAemLB7DbOXeg9lvnLtNi4MHDDojUttwPq3TiwPAAx4Ko2Vc6Q8ojMU1wS6SFxHar9n3iqvZVHWIaMm6zkURGqIinV4LqbeOsECOXSOIwiLGpV15DOo75+o57oocuIQm6OoID5Bpy8X6CInREHJS7QcG4/s5zp4rk23DkHTkhDNLcVbrbHtfc/vanfeTu3FJoQQAPwzAO+JMf6vp1WPatWqVatW7Xxb6GzNnWc7TcTmzwH4IgC/F0J4h3z2DTHGn7/VheR2DA1i03I/20mFLzU6qoxcSDeS45ooIkrJ20gRIjQsm2/pfbwfGlEi9XBj6VEvjWe9aktvj95F6zb1rdjSugScmoiSYmlGAl+9GY1u6OEGibG/KepFzgDF3faHXW+ZfA+mOqAn76OWyH9I9y0jNhgl471Zy+7XyCwnbNY4fhEArIRPRKoDI3UOhPdxY55QmauLLb2G3uXVWfrs5iw1LB6mvm4PBSk8MFwh/k5HV8a91QihPHb7UvZ1QS2uLlM5E2n7jqBHVtiObdaAFu/sqwii+UxWOiN1tiXKipFnLxrsAsgePQA8O06IzVhSThyOUv8oP2qY2+ETcHZF3mTcDVq1K33JtmuSVUFsFCm1cJwTo/OpUlaWY8M5JjyuoYtSGsjYTs1jUNfxskQ9FanRRLa5Svq7l9QXBISo64bJEpvTR5TRSbEnzUruQ5a9KuvioqWAMkoI6IrVaZSURWw8yq1kQym/NyoKxTGuSS8AZOTERyetQ03konQO55Pn/RAhskKcq5JjQwFIPkMHpqOY1sbzLW3KjPRB/jWXLf9TVixXnsPK8ckNUbTolKKiLpqd2otNjPFXsZa+Wa1atWrVqlW7E1Y5NmfAcgTUsvMdPZCBcFQW3guwr1LqBchxVUYbcJ99vugiNkQvFuJmrNZ4kvb3zHmQN/phqcFhz6W3x7b65Jv7RsOfdVoKwjRwyQqpN2J1R1Q+3TsQoSwfyDwA63ECJhmlS2sAZM+ZlA16vsqJcKkjAHSSjtLjyQlHh0V7gYxSsWpEp7Kkf25gDIzQEcTjQCJ1DlJfXjtMyMGzw229hsjc89P03c5BQlDCIe+VzrOIzWCfEvRlNNFgX+q6nx8wewdjuX9Cap4ZXQaQPWoiRbsmxcVqenTETl/UEi345KqC0FxpE3IzNmQhRkptSjqHfYmo0vQVE8PDGpXcFO/xarLCRW47NWOI0FB/5OZMUCpGUM3zvXwUjpbv0hoAJhmlRBiOpcMYaUiu1qFBIjS1AZFA5amlI5GbIuGrS6GhzxWvJWPQMI/YjAYsiA1FxxQp4HpelHWxa9nryZCPOJKTydkbm2g7Rr9ljg0vlj97UItGuTayVsk36YkQ0jQLLrVCR3OsiFztR0466Lppr6bccTpFTQ/PpVO2Rq7K3z0RWr4OPipq7v8n9NSh2t21M/liU61atWrVqlU7vlWOzX1sTRMz18OEJjTKAxFeBgkIHW2GfC99cXe8DB9lYHUuFs4joak+QVN6AIBJlDgsJ5ZHL+z1XiV1KrAOkRpG8ADAvnjzUbgb2YP3+hfGW/CwpGrRSD8Oct9uDxM5hMkoB6KWStVaqv+WUVFlG+HVWXvGYyn9rerO0lZ6fTcWieOxM8ttXwp6MSIKcyhesqAmFrGhzQ94rpR3mMq7MU1IwfOjzLGht31DNFemhxIVJ2jPIAUVYXhgIvT2xBNdMYqoRIiI8gDATMq+Pk33f054LZzHzx6mv28e5EitQH7PQdlmolWqo2IAtuzZln3htUVGMaOgV6RxWzLOz0pE1WJSJt8EgMWGaB1Jf3Ou69zuGe+58COobD0TvSRyb2ZULzacJJ3bVJTmoQcpJVdENaEk2em2RCcRLTlsjCbUGj0ZXUuzMuor1ank4xDtIfpDXsXEDAjRG0ZKjdZputi6OH6JIrGqTJzPZZSQf07liKRSRwe4BdclNSifqzwTFMeVRhqKgnPPP1OiFlxbTYfnYpBAqrQrCllGHCpqYnhYPjrJJ3+0UVEaGcckqq1+IQV0qp//b3DuUe9MuI5a7n2kGxNjuFBbURenpdWqVatWrVq1c29nDrGpVq1atWrVqp3MPHJ1nu1MvdiEENE0KyUN925FCaSrBDEHJRYEOMpzE3X0EuUMUTXiT0sXRskwvpZbXgyZHOSCuBWhn/HQF4LOa1YlpHu4KlXA9oyQ2roEnBka72ORsg5SHiFYqSNhewC4PEx7J0zISCh3S1IbDCRMe27avK5tnTBgAyHPHJn02jxtCVHIS4m0U5M4kCReCaVmqPXgQObB1EjGh3T/lls51EuTrQ6GXu/M8rYPSbZ7s7RdwnBlEncZUt9MzdbEtNyKYnkkNDfT3CkMG795KCTicWojtwi4RUWCMwA0+34bTEiZsv3G+baY5r5t5feZCAtOZV55MTFLbiTBlGTbQKG+sYgUFltRMg5yO36n27Btdz1yq4YQ/kKYzQdzuf+M42TaIfxbvyUV2+6eAZtGscsNGbRLMp+5FbVndAo4t6NL/OkTyVoi/opbwItymyS6ZTc0Ycabkh1UBfNI4mW6jD5ZCg37LtMLZBKxCb932+jchmAqGMoJFGKX+uySe+ga7tmPWZZt1a3OpSPSGgIv0yy0LvS50S0wtx0E0w/rwq+VTpC/ZwqcY20JaWbXsjxt8zHIw/BtVvKwGQ8mHa4CfffEztSLTbVq1apVq1btZBbRr2J+Xu3Mvdg0TdTw402TVI6mSeXWSJTbaLtmWZ7TdEjEJQmwOFfluSW8XDyf4BIUAl0JevXo+uYZRd0EJWICSiI0C8kDsbfIXiYT93UQJ98FxvvwooH0rMNI0Jhh7tsHXDLKVkjIE/HkSdSe9zhGWQJdvDGVqpfPjTc+nZZE2ueHgthoyHX6e9/I8jcu7FrFCA8FsZkZxEaSjmZZfIb0lsjcbFWiGEBG6nzyS/WiiyR59KDFI9UEfuwLe2NBUgSt2hMi+ExSORApWkxznUaCNA0cUjPYW0o7pT1jg3QIURqCFu3KfW8uEjq1LzkumGTVmpI8BZFbiuCdvT8RmtgKadUnO2UiyoER1SRhk0Jq0sczCeGPTPdhqkSUTZEy+a5nyDqpDYg8bbcledgmp6TopE+OqI49CbuGlK5CeVIXhl+vC7kGMpF4U+rCxK6a4LQXteAHfn7J14Y8rMgJRePQXxdLHvbISU56uh690Lnu0AtPpAWyUB6NSKymRekL95Y1uy7FQV6H68nDHVTSMNibDmlYyu0gN1hvJ0itcHoWLtRW1MVpabVq1apVq1bt3NuZQmxCSN6jJnAzKln0+tYlk8t7wTZM0yE2ukdOAoocLC/HiQMSGWIyuUZCS/vk7NclebMoElEiogeHIqVPHgBFzPYNYqPh6Gyj7wKKTg3Mnq8IqmlYtnjhTHi4NcwJ++jhXmoSLMK9Y+6RqzNVCHilo0/ESSoEeRPNzHhagtgwnPvabLMoZ2fOMGCD2FCMkPcXTzrMmXXTJLxzIcLrzIaDMoyfRy/BDu/Z2d89R+EIry/nGRSPVyYQ0QsLhzVOgLGVUPP2gJ0tYcYb+ZqB8HIoDrgrnJ6rwmO6ukhh5eR+AFliYOVgCxXD60l6SieVeUR5ZLLTdpj7lpIClLpfkROhiSibor1Fm8lXoighQ7uL9AtyTVNKQTDce1tgvt1h5mwp52SNTEQWyTPcFOECZZkFjmF/yDWQ04WQ60IU2vNcAMN18aiFCoo61ARZCJHooxerazW9gU1LU4qaamqFQRlyDfRx5eQLIpCrLlriEQMiJwyxJ0WoH61yfeDCvi1io5wtxy/qC3fW28kxJ00OrORa82XPfZi5WfCnjZZE9AsGnleriE21atWqVatW7Y5bCOGHQgjPhBDeteb7EEL4nhDCe0MI7wwhfMKdKPdsITaIGA6WGkEwNPvDOUW9c8c1SWL60+5D26Rx9m99E6dXYl506c1zf57ler7J4iQRQqYePgEnI57IqZkL98LyQDRqK5aeFU9RYcCQK8DklPS0yYEYSKLAzYERExPEZmh15I2xfHJVAIDq8eRCDA5LD4t/K/cDwFx4JPtTESEcJ0+afc6+WBmvvHWeY65A1+NSz7PjkYpHzyScpu1MYMgklEQciFZoUs+JQcPGwnkSbs2S6NiQvCZTT2nAwCU9VW+2J4oi85bkSI7HYSnLPzjIBbWuvxllRd7S0+OUyuGSUQ+8JskpGUGl/CL53s7nLDZZ8nuoI7maCC9rbOaVRCmRX7I3N9FuyOhLwbFx84oIDrk9xfp2HBcKezISiQgkhQhTnaSAoYtO8pFJNiGkT3Eggn2rRcm5KNILyPiOXWqFdkCelH1+yPFW6QXs8pS16KPfGCHUl5ZGuXKO66Ll9kRHreO6EDWZmsneh2QARqBPo6JM1OCtUE92k1n/7G+mYPFCfUVKBYdCKrrtUyz0gB1eqE85NrGLkt0PfJtTEgx8E4DvA/DmNd9/NoCPlJ9PBfC/y/G27Ey92FSrVq1atWrVTmYR4VS2omKMbw0hvPqIU/4SgDfHGCOA3wghPBBCeCTG+NTtlHumXmxCiBi2S007b3Uh6NSRa7PUZHZyrUakWI4NNSvkHO7t8xQ6KoZzwX16q6ED5D1y1XSxSIHn1ChSw+gGg3SIt8eIEHIslFMj3q3X00kXl9EMOeJJIgtMwklFahzHZjjstq91KBg9LkZsMflm6VmLRzgl/0OixlbUnSF/Il8TXJvpbQKpLr0Lk3vk5H2oZpCgSMbLXEryxgUjdQgQSEqIS5OETD083tNr6EnTC3xuIyEce5sSpbYh99wwyNO2pAKQ8V1s8Rypx4bxGCep7G0pmxortPFQSCp9HnxnXsk5Kv9v0jxotJj0saRyuHqQUBmmcjgYZOIMdYN2RDdoJVySgUc2gcw1kj8zoiWfSDqGTYPYsK1M2UFPd11aA8BwqaZlG4nGEcEB8pr30UlEScgn2jYolUYnEZlz/Di1ghcnc5zIGVFXzufV+vQCI5eMksl7TRBcRhYVPZJnm0+x0BNt55ETjxxYjo2mF9CoqH70AkAnJYrXpZr1IjbkbJXISeNQEwssrBRlXcMz0vLNOleui9OVYfSdWTiaeDe4Njfu2KOTpOYQmz6NqHPKb3lRCOG3zN9vjDG+8YT3eAWAx83fT8hnF+fFplq1atWqVat2cvPh9nfAnosxftJt3qPvje+2VQwrebhatWrVqlWrdhr2BIBHzd+vBPDk7d70zCE2bYi6TWK3SAhwcwuKYl+eNEyCX/qO8aCEdl1huhWVv2A4KLco8ueL8tym+9LpQ84bSoybsHNuRZH8R+E2QsqrHkjUk+805Fa2mUjkLOBaJQ2juLZ18DBgtp4iM4uLqJtsj1H6fmyl70lWZSiyCObxtllgzbRnUY4d4VvWRetmttRW2laSVYWou0kCbz53sSnkPtk2WsiW0HAjjeVDkyRA+JLxjl6z2ZQikM9sXEp9sLkp9xSi4JZpu0jsc7xnW1KunLPYNGkLNlPZD04SgfXFk92iPBJqn5vkrcGlpDTQbcRhCZdngcbutowSTF0aCZKID5ZmK2qa2sjQcOj2oYyL3Y2l6GSHmJ2OzAC/Pc7h5A+OJHu4hF9TTI5pLGLPGlLCtGxBtVNmlpe1NMvwP7dGSW5XSQh5bmw5EjGQt8U8iXflwn+Lvl2U5GHdknUE3qmJj+eaaqTjdHtd2n5otx65HUO5Br8d4zJsp9/7t4RmPoTfGEVGO1vZLN+Sh71YH6tLYcAlCbzHIA+7tAZ9KRVuSR42U2XJdDRu+62PwKvh3k6YrxP2fWRqhVCUyy2wvraf1pZUjPl/431mPwvgb4UQfgyJNHzjdvk1wBl8salWrVq1atWq3f8WQvhRAJ+GxMd5AsDrAQwBIMb4BgA/D+BzALwXwD6AL7kT5Z6pF5sQkkeTSWcmEZ3KWYtXwBDkJb1YEnXN/ejNe+dDhaLEozeJ4ojU0JOnFzhygn1FvR3BTj1dJtoz5GGfgJPER5/no+0RkdMEnOrJpyPRDEsuzOG58oHblLSiYlOXgPOGsGB3RUgvTl1SSRiy6pQS9FJ468bDdhfHVdpDQiVDrreGqZyRIaAeTFIdlkLinW/Qc3dsVgDzTSIm6e+lICdXNtJYEi15ZHRdr/Hj/NBGQnU+tJHCoxcbw+Leqc1EbGJRHo+rzTwQlzYTWvDSSUKJHp1cS3WLTJ2R2vfU5iW9Zm9jImULOiXH1STVheJxRdisIg3yAdEkJh6l+KFxl29MUzmHh+m+zQGJ32xn7lwCmNoyesItUTYKP5pUHcPUlxTMI3l4wgSsXiQPyKJ783Jesc9bQx5uhLw71/QC5SSnMODEsN63hTw8FNmDOZN4amLIHpG3VflsyUfx3BfdRLZeLI6BEAz7j+bJvPIyBQNXhx7Jg6Dk4f7QZ1pTkIfldj6lAgm1ptzgMnwqGi1HEmlnBrXwiXw9eqEInQ280LI9U96h7H0pFWQud8o1EyqTh1GUTYRIkTpb/BrUyLe5jzzcJxJ4r+yUoqK+8BbfRwB/806Xe6ZebKpVq1atWrVqJ7MU7n1xKLVn6sUmIKIJUeXA7UDx7ZhIw1KRmvQ90Yq+ZIXRCdt5xEYT48EgNuK2UoRpWMRaokAKMs+He/Cll1BwbDSRZYnY0CjgNjLljenhjigeR+Sm5NFYeKQTMsxIYfE6LNeCCM2wSVDQjXn6e0/4GWHWDd3OHBomgBRXjg5i4L69aZzUif1NvgO5GETDro439ZJ9EX5bTmT85av5vOS5AJlbQ5n/KCHIW+NU8YdGKcz7oTbzXMjD2B+kthNdGI1Sn08lTH5lEkJSlI5eJBNELifi2U/y2F2ZCGIzvgkAeGSY0CKKae3KzT6w8aBes7OZQrPnwt0hh2ewL2HmMh6FaCDFAelEcq7LuJPDZWUEKJS4FKG/oYj7SS5UmChpg5SU606/l7XEcGoAeFAQmyttGl962NtEdUZOJM/UWyX1BZ1q+hAbuQ3D1A+XJWpBfotN/LkhF3EOTonUkIum6QW63m9OcSAfyFqmsKTlXBxKllzl2Mi6Vp6L4dgoWkShzcYNnp5o6iL9z3HleiZy04TuMzQ4rovn2hQhz7qOyyp4xCZLNuQUHRSfJHLiURPLscmh7g7d8eWarli5VBbHEcfzHBuPZNvx7gAfUjYRkb4kmIrYXKAM26dpZ+rFplq1atWqVat2cuvTUTqvduZebFrjgtuB4tvxTDwEvrW3LmKgI70PGJG3ci+bnJWR5diIR8fkdZ6bon6rkf3vRGbNuSHdjS5RTzcWp2i7KRDYGI9uIskEGyf3v3L8gL62e2GthYjj7Rp5+5tlcDZ4AAAgAElEQVSC2NDLuzET7oUiNtIOEyim/KVV6VXS61Mv1MzAOKKIW+pjIjUvHiX+CaPhrm1kxObaZormWUxEBE3QkcWU9cj3z6KEUifpr41hyZvaMtDTRMc7nTPuRL9J3a2oGBMyKmeh/LwxiSCZvoGy/g8PyqioG6PU9w9Msuz/B7dSHRYS+cWIrPmBIBGzMvoLABYuGSUF/xoXeWR5KIxsodok0zEw4m14YHheyrehWGOJPNIs0khhvCttQsrmwqG6NEqft6NyPgOlN58+EK7NnFwbw/vhMM5L9ICcC+p6WKFPTbArKOgNInIickn0Apbn4hIyatSjrOv5gqhJTlzLOgwdtNU6kTwg820UOWEdfEqYIgpOUEIXnaTlOoFRwArllRyblUNubNtcflRFKVmu5WwRkWM5nSidHoG+DoJSAoPmxPwrUe6ZS4ZJ83/bBnSiozy3p7fsEl0nz9Pyi5YOzal2d+3MvdhUq1atWrVq1Y5vEedW/bjXzuSLDfkHdg+Tngg9BCagG5JrQ+TA70uji9So9z3oRkXRo2NCSL9/6/VzgIxkaLJC5djQ08v38EFVjBBiygZGbczNNUzcxwSWTMCZZcg7Te6kmmDivrlkwGPEEwBcH5WIDbVVeK5GnvWhYexb8cZXTAhJ3ZmJOVW0Wq6Mk8dOpOYlo8Q/IY/pyihHCA2F6zITJEJRGaJUpu1+fD1NgvvfdkzbyAgH4UkwjQTRDNePgIk+c5EimR5iPF9xNYeaJkTml7jj1Fi5bFItsM3k8iwkIotcm1ZQBnJ7gKyhQ54P+5oRSGNB/RjBA5i5vCgRuYFUZbBvEZuSD9OKphGBzeWyOwnZ5suC3LCPidSxnYqwwXjQ5HtoAkgiN12ODef2voQ4keux7CFTcn0TpY2Kggpio4lM87XNgiI+0rcusnGmiE2Gnsix8UaujdXBWjnEhmU38BPL3IhcF5mnTH8yVfiHh/URQlnLhcd8ro78Gq7LymnJABk1YhvXkVntx50kmHyMH4HYcH0xyq7LKwqdc31qFp9Soaxf/wtCLreMyrJlnt7LxcUiD1+cllarVq1atWrVzr2dOcQmxmB0CrregLLwV6VH15d9wic4W6lmQvo+UMfG8AKGTq+GPJ+518+xiS1dNJQm3xw4NKmoWzq2omsxES7GlkRjrYzXsDlInBNNwOmRJ7FS50LqNi+Pc9GkoXcLZL7NQNrOPo706OlF9XhaRGhUnVOUgcn50CSJAEaCIjw0TtEyLxomvomNUkrtzRwYtnnalt4mjxYB83yiKGN0IG1l9Nf1ZebwDFfp/s8tUiTSzXmq+HyW+oDIRBEhRH7PqtQVauTchcmQerBIZe/LSUQtWgfd2TnIRIlzQTKoLM0km/TsF0ZbR7V7tmQ+idryAxsJHWFCygYZqePco7EKjDxqDw1ic5jGjg7p4IDaRrI+yHcwk6R1UUmXJDrqsvCNJqP0+U6mpmTeGPVkeDvyXBa5TrrehAN2OE99uyvaQERNhkbcilw2RWmpDaW8tW6EUPReuEPqyFvbMxybfYHbGC3W8eR7eCYeueEXfY64IifKdRlIuVIHeUbMO6Qlc7EiNz1RUQ458dGVUfkmfXpYkliW9deo1FjeC+hwXTSSVaOVyioDXY4NlYBbTYJ5hD+/Bq0qEKGmPNejVAsXlQVkjk0fSniv7CJFZFXEplq1atWqVat2buzMITbVqlWrVq1atePbfZwr6q7YmXuxWSEouXNqdMdVxlrgx+gIrQpVWqEltwXlQ5CDQPE23Ltxe1qENUlaXgpp2er1aWS4bEmpWF3TBcwUqXTpHLYkdnWbW1EGViSxWLcOSDx0UutFhKdA9wyJ1S0DCRUmbA/k7ZLxoAsvpzoLPGxCI3NSyhLqzmkA5PuNvHVwScK8HximrYgrbdqSIrl0f5Xg+6HZU1PCY3TjraH1uexGt1DkXNl225Ekj89N03bTU4MH9BoSXD80TSkUnj9IezrLvdQ/o710r+FenhfDvVVRF255yg4Lpge5T64dpv2jZ2aJEP3MMJWz1aRx3lmlrS9LPPWijR6u53hYYvZiU8ZoM3XMZUnl8PAkhVpfkT4fmFQdzw1TW/facs7rfJqZdXGQ7stZMzjkVp2I43GLc5G3Y7jtxq0gtvmysJMZ9n9jbCQeRiUJXdeQEnet2KUchTw8lS1Hih7uyXza6sGtdeuPSTB9yPXQpKtwzxbdyZFbKIF3kceQz66hS0qZRfK6Yporl16A5faSWXV7RJ5PLhkm5/WxiLR94d5HlW3ua7e6rEAhkJ8jK1eu3THx4dd5/6e32HQut71juQ3HZ7fdIvKxJB2R1h4xxM72m96rLPf+Ig+fbjqHe20Xp6XVqlWrVq1atXNvZwqxSfkugobx2eRmGobLN2Lv1bo38fR7SYrLYY2CrDCBn0FpWhcyyDpQrp3e2WjeJQ+T2KgEOFePog7iKTIcd0uk/C9paGz2Bigap3dZg9TYUFh6k0RqcmisEN8W+f5Eo4ZREvWJV98IQkSeseFGYiGE1najfHdmkkqGKkfjjVMUbUNgJCYnZFhrDvM3InILitJJO4S4OxARuYERbAsMxxRSbbuX7rO/nyrz9EFCTcYG2iJh+qmDKwCAG3sS+r4rfSK85tGuIUHfkPFQTz7NzdlOumawk/v25l6CVZ6U+zPNAMXrnpklBOf6dEOvIXFZky26MHtN4GeF7YRoPBSCNgX/XjZJofQPSp4Ei4Y9M079cZUhz45nGiw6Mp0X3w0lFLwlOiUpD/bmhkArhGl6uAxxZ7qSTRFOjDbcW+YNCdOUEWBqhb5UJsyYMJPUCjtCAN8VSKuF60Dk4AA+AwgOL4U8vBzlOajr2Yc+E7ET9PjQpBcggZnISUf236IWipjwWAZGHEUejnIO0QOP2Ni15JETj9jYKma0yNVXH2pddJdISRMdUuOsjzysxGUV0nPITV+4t0uGyXBvWyePUnWTYcJ90dNmvReK+896iNnLUxLo4//Oi2IVsalWrVq1atWqnRs7U4gNYnrjJTozanL16XnoG7GPoFYvxLx5U/BKRbfkc0rfi3fQ2DQOUg5DRTXMXITN4oLh3rlo7vt3POtAjk/+TH8fMMxbhPnoxQokQb4JkBElBQjEQyTwwCSbFFEDcj9oskom35SjCtAZY3+Q9zMUQcAZpe+NF0t+x2JWegkMfVY0wfA3hkSCnOQ6+5oe9s1ZRi9m0/Qdw62H++nakfBcKP8PAGHJVAOC5uzL/ffSPa7up/sSAQOAgQwauTWHeyMpp+TWDHeMJMDNhDg5+g+Gl1L5g93cJ4e7aRw/tJeQmUuS+JOcqquzVO7VgxyCvpBklyOmOBCUiuOc51MXCaSI4wPjhNg8MroBAHhwkLg2VuL/yUmq0xOSgFPHjsvODq3wxpqZ8L0kFJzIWXMoiM0sIzY7Mp4cX4Z9M10J000wGSaQkb6cHkO8cbrLpk4ZsZTw30XJ82H5Y0PE0gSZitisijb7ZJjpOwlbbsu5TomDyLQGBrE5kM4cL0UYUQqIPfHL67guHKpeBMGFXftkmETmrIBeJ2zdITervnBvn0xX75WOFqGw7QeygF5H7KIvpYJDTny4d1E2yrJnitjI/4jYg9is4dpoehT7vSt7HdemL33CaaZUqOHe1apVq1atWrVqZ9DOFGITAcxXjXo3Y+Nt0PNYt4epb/imxYwu8ZLhFIpihIJlk8+dTDd5Aipatyi9KWA9UpMTQloUSa5hgkbxWjcd78RKsvNNvCOBz4inWQ9iI2U3c5cCgo6vcVHIqaEwHhGiDRFQ2xOeDOX6gSznn8Xq0pEchT5eABEaepU7q4SgEJ16bp6ilizfZCnoxVgRFBFYEwSlPcwDEcRjX2wIynMgYyVowr5ER10fWo5Nut+ufBenpfAcEaGBKafZn0ubBbWQCJqh1HVgoqLCfvr9uqBFTw4T14bjzrbe2M1tDhKRNXARWYMDQbqGHFtk00A8GbtWkm9K5NnLBtcBACMDNb5MkJo/2EjjPpskhINoSWx7BnEhApZTJqVMHzMSbW+aEZvr84RCcZwZ1UeEVBO+DnPfrpTPRa4W03qI6KWJemRSSCKYK5nre4LYMDqKnB4g8zE6fBMK6ZFrYzk2RE48ekFUQ8pnZBKQ5/iGoM4LV26RCsQnpZTnBZFg5QoWOULkOyI2il4MpPwu38TrhHq0pOTYcB1LmX4qELExcIbnnLA/Vj7Ksodj00ksu4bnYm/gE08ShSs5Nr2XHh0Vta5sV579v8Fn9GmFXEfUXFHVqlWrVq1atXNkFync+0y92MQYMF+2+uZv92zzW7K7xnkddp/4/2fvXWNt2bLzoDGr1mu/zuve27df7r5tu7EgQRGR5RDxJwgSmciSISLI5AdWQtRyHGPLbtttCTASQopj4kR27OA0tiEIQbCQDJawMI8/iZAQthDBaRvL7ftw39u3773nuZ/rVTX5Mcc3H2POWWvt89pnnT0+6ajWqlVzjjmratep+Y1vjGFkfojKm3jsmwxRBW7FhZXXAqsxaGzicUgnMqJWfFG7yNaYNQljFL90K+upiBCKgVXmGjl0oK0RTA1W0UTk83DU9D9NlEcDmpMjjszCSvrR1K3gH3CkTTeJWCShpfFEgIxiiIojLlmndML1Fu6uQrFLIqK7S8fYPLoICVrAeHBQT4jGOefV/jxiX3jbLpBjhefKOqA116KIV9bSX49oOx9x0yPiLdLy9Milz+wFl9BofYr/MCeUWVjMuawDsyLIHYS5Lk8D0zFhjQ7XBqXJCVgjZs6YTYhLKmCOUjsFzdYtZm4mEdX4Ohs42uMIrT13PcDGoTwGEdF4hPvftUfEVLgHWZt2ETM2aQkL5LPBAxhsWTuK89gwKzVNx9DwtYu1L4C/t/nvAyU0zpixOY0S/uDv2zMnFbYk+ZtFXhk8W+T/H0LnQhSeGxftJLFbXNGLyKDAHqWMc1FYAPZC5LMZIVIvmkgv9YkiAsk2OWuRaV08wZHaJQqsVNDQsTapoodMbfs9RftJG6Hv8XlsYDdmqUT0rNRierasEBUV5pxOAN/ifDm4nzIdk+KZYKdebBQKhUKhUFwS9nqFe+/Ui40lonXXhPwEpRwMAL5muRjiaAbxW5aTIVe3L8TbP1Zea6zGfB6Hwvix6Ge2JGzDMfgMvcGE86i0osM47wXyYyCny0gUtmzA2EQFAnvvrE7njNXgeBRW7je4QOItzrGCc/Bwwitt1kAs4nwjo3Sl5Rc3vYjCivL9nHOE073FAREFHQiAPDNn0aofOhlE37Q+CocZmyi/imSp5DWCrqiNsu+2cl+Tai38Cj7ORMtRMv428oxgPYQDKzlE7mBVOedzYubhenttDefOmRxzfqELt13POMfOXrA34vO0XEBj4frt+Q8EEUmTKGX2q+MTIiK6zTlvvrHP/e9zzqa9KLPqhP8uoGtBJmCwVMyOdcvQ5iFHtz1YHyRj8JFJfHKaqBhnJ3Im9czc9JwxO2ZkrbzHmR1EBOMpa232mhBhiGglubL22gsfHRU9R5gU9Pd8pr1wm3XEli1FlJLPXN7nLIbUfWCOVcaZKGNGOxSG5Bt3xKxizE7W9B/lqKg0Ai/zchSioqTGZu3nKq5TDKnv8fpEfC8O2Y2torFZxZqemm2p5YlQ1RVV7BKVGSzFs8NOvdgoFAqFQqG4HCxdr3BvfbFRKBQKheIlh7qiXlBYa6jrG0+n9ja4GSA49Gw/XAYQ2vmCl7HIk/utpAXfppAbxgKa2UgxWtQ/XBSW85yHxIARfc7jG/nQV+mCcvbOu0CfI+FYv0zLC7RLuKAgcI0UwqPUR+Qp5XEayk1EdGviXBFIu48CfndHTsw75vDoRVvgknkXPByWhbMQ0CJxGxHRYu7mgWR0k6RqJ9Hxwok814twDcbS7bZC2YpcZO1F276YoNsNwfZ0wskQJ0HdC+E0aPrTPQ553uPrzlHY6/1wj7QHUx5Llx4Lt0m4dN716IWyuG9FgrXYZYdzN75AckBOinfK4fhsr4vKWVyccRkSDjlHGPlJt8f2OCQ9Eg9DUIzr33ABTYTLr/ZD/1N2T+G8ww0ky3rYRWiD6/mAw773eWJIoYAVZixkX6MoJQukIZRuphD55iJPDx6LTNR3MQquTeka8qgkyXOfhVtGZtKHKypyQ6AEy6Trku9FV4V3E8sxiP0FD4sfI9z37HpaGoiV68nq8tDnqH8plK648eNgDsxtLUS8NQFvzXay36TbePydKIZpbUE8XMuyN+iKEk2EYFpuiYKQuCv8/6B4+tipFxuFQqFQKBSXg+axecHR9U1e7j6CTyyHVZ5IzJckmeog9hP9VNKRE4UVHQDxXSa8K4QsBnEhRMMpc+DGx4JJHj+SlXXcIRLzYVVLFEKD7SoN98bWMzWxItGXc4D4lfdzmPlBxFrcYqbmzuiUbTvKAeLe1uRMjb8MIvTcJ+pjMStKIRARLZlNOGHGY2/MBSf5glys0rIVRERN55dLbuPZGDBS4UJ0M2ZOOFyZo8rJ7rlB3dh3g3l978S3wRzBHh0fukZnh+78r454BXwjGhMzATjv60N3zOqQ7YbqCGRZkLvHSfAOp05li3vujJm7dXQ/BRYEgmkWSp+58TdcJHOyF26syQnfG8zcPJy7c3t37Vi3497NCywNUQi/PmLx+GTK7NU+izCjcPL1HopR8jWS5QV8Qcpwns6WXMpi5cTDh5xOAIVlS+nnPVshRMQNh32nRW7TtmBTwa6CsYGImKiQ6BOrb8FslsK9g2CeD5bsRTQfsMB4fvjnCFb08fMkEw+n8ykKeKV4GIJWixBovidj8XAnqQhhPy79gq1ISlhjimLbXOvUi3jB2JjC8zxjjfx3IZgu6fIrpQ1i1sQzopVEfVlxzJJNOWdxrp1NIZRWPFPs3IuNQqFQKBSKy0EZmxcUltzbPd56Y5U3mA1fvE4k1vIrvWgl6X33YnXjixciPDRa1cwFY+MLuVXCQxPbSHWP1d8o3RLlq0wktUIJB2ghztZRiConHCNfyJLngeR7PqQ7Wqn4MGVKts3ULadQjJGI6ObIaSyOYnqFiEZRaDARJcn2oHlpF2nyOJx/JFpDIUoiohXrLxbMOJzzvFAcU6Zed/3ZZPzQXHSsA4lXYusDXqUe8IrqkDUqh44teX3fMVKf2bvv26D46MHIhZo/4hDlt25w4jxmahan+QVvmFmCFmVxi0Nubwb9z+SG6/9jR872KzNXjNKHAzMDspiGkgpeqwWTrGsxC2bZVu4GGJ2Fe2R86s4HCnA+5BIN31i4QpcfccHLVmZqJKIps1UT1lKdcemMdaSxASvVgD0UK2p0G2uF5szYPFqy1mbkqKzeswr548mzFr60AfRk6e9EERMq2QtmJsB0nkeMDf7eViLsGkxwSJIX7HitnHjWxAntiFItCZiSecN6JSTZrBXxpYitwPOkT78n2p4NGptS6HOmdYHdks4Fz8iMORFzjh6EXlPDzyHPYgzMuaZ1yZISDjA2UrcUP0dqjE1mN+6/ovsZsovzcHUlFa5XHhsNqlcoFAqFQvHSYKcYG7KG+t4U4/GxAkF0iWnLqbeHIgfw1u6L2vGKex0xNktOhjcyeRE5139ql6jAzAh9S8k3jpUi/LTnLCbAPI8jxmbBq3qshkGkeOYJfY6iFcQYUSU8Fo4qaTnZ3mHE2BwyU4MEamfG2cbK2q8Co3x6KFeAxHmjc8HYTLHSDhfBIO2/SFI3aTkBHetNUCDUzcNtO9bNrA44Ok0wa0REixvQwzBjc+T6feXIMVKfOXBMzbfMPvRtjlr3Gwol3t13mpRvHDoGZ3404T7DufV6Emaw1swQLW8wu3QjRHu9etMxNW8cOtuvT10ZA2i5wBA+PDjwbdaIepqCPRDrE2ZsmvOgkxpzkc3xKTN+p+4avnd+i4iI3t27436PimBCS4V7bsJJG0/4XlkHEolWHIHVLGPaoK61ISJaMTN3snKMDUosAEikl7ChssSB0NrEejmp//DRjkhW16XRUURRGn6+B/2KHl0Ixiiek2ROcr1JxNiwHZSFwVy9zqXEnojILCsZm/jZJpr6MgOI9OxyjY0fn0jaORghJBPlVfQmRGHO0ORB61SLxoptVqOiCrIgT2RWWJLkeyWhahbxVNJu+f8w0raIBEvYqoLO6nnjOuWxUcZGoVAoFArFS4PdYmzIvd3jTbiJXrPB1GDb+JWdyG8Rv7RKbU2fbqEZweqNKEQTTMQSwoi3d5tETWCbjmEoLTje7LFiP+3S1XO8ysSqy3TllQ+Ymj6OEPJ5QNger3ynE7f8PBgFxgZMDfQXSHkPVmHFBQjbRei/FUzN6IwLHDJTBH1GRBCE8aMPvpY+lwy3HU9DoxXnT1mxfmZ5CJ95m5wDIqLlERgU9x0RSXf2XCTQx7no46fGD3ybo8YxNmCnXp06huVgxjlXZhy1Nov0JhwthPuIg35ofcDMx2FgUj5x4Gx+675jiT49cczNGbMl0Fa9f+OGb3NyNE3m2h24Y9opR2NduGsHjRURUcsap9EF3yscHfXhmWOgvnbgGJtpk+eGgtYF18Ny3pYolRKtuXxDy1FPPvoNkX9eaxPdg7xiP126cR+POAqOG0MDU9J+yHwynqCLF+P4e5OrfaE3iYvp+lwrYDJqupM4X47PlUXJNrMbsxfcv2cvfITQFhob6FpkMczYXsU25oeRxGUeMq1Lls8mPhcVRlwwzraosUnzyYB5khrHeKCyHEleBqck0HGQBShLGptGNhdMTWy/yiJhW2CKqjl7nhfs9RIPK2OjUCgUCoXipcEOMjbhc5yVF+zNmPUYKJy3Fm/2pQyVKMzotRFIGtKnq5z4s3z79blc/CoqjM1HIPmznetwPLDSFQX7EAWFOcfRWSH/RNoVVpVY0ZlolQl9Brb9hItujhwbMo6iY8DULHmZhBw6pys3phVnAp7NQ/+jCzAEfD3mzAQwldXwKjkJwkHuHr52+5xW+NbUsSYzZm7O9sPc7x+4saxZW7PyjA1lAJPS7TGbNOM8LROnIbrJ+XrA0hARHTBbtc9VHBEhhPuMBDPoPqcraJ9DiSPBZtPAirwydVFQYGq+ZeyYmyWl0XBv3njFt3l40+l7ljc5S+5N1n0dO60KbquYVUDhT0SptefuqOMz1+brF44R2msjXQ7/QSDPi7/ncY8H0tAzNg1yKbE9r32R2XgpaNigETtZpfc4NDB9l/+hhIKylI0FkJl5ZYbYFT8c5m14DPocJNB/iBW2zCVDRGS8HTAp6bG+bZKFl8fg2SNmULv02ZNAsAj9gD25OK9FCMXfcUyWPR1fk/4rtkFKC+YmtmUyxgYHDzAKkhHyrElZ50IUHocyh0zCmmxgMbxmq8T0+2GnD5vSub7qPDaaoE+hUCgUCsVLhev0YqOuKIVCoVAoFC8Ndo+xscZT1WOZII6C26Rld8YKroISswsGEaLhTmzXOXW5FvwyxtKwsJJGEFiGY7JidTIHWlzpgKnoHgnaWEB5xpw+wsyLYZroTtDlKBgYw7ugmMqHIDQU3wyDXPpqkW5z2jn3BVwHdg7xcOh/xC6PZuFOpuGCkEaEJqeJzlybA3bV3J4619DHpq7EwWLsxoFzQkR0su8Ep+t9LpfA7qYGyQqjKnxwV+C84JrVCo0SEc354iD02Ydhs/sESRGbKPGcT5Ao7i9fEDSyB9cW3F+vtcENRkT0cPIREYVwcCKit245txRcUPNbbkLjE3dd/NlpIjcDix9xb6OQ5vLCHf3huXNvwd1HFETbx5yU0AtNDZIihnl41yaXqehxHw+5ohB2vU4TMuL8XCwL4mGfGC7tF0nyYvh7WyhD8fe39qHdsXiY5+HDrtFIuCSSdA5lt0z2zImLIvpEee47/p6L4mEp4hXnspigT7hofLg3yjD4ooyxkDa15yHcP4O2Ky6w2DYOwnX1At5C3dp6SQVhr+AqksUo/dzj+2nDnEvJ+GpFMH2XvgBoPverEg9rgj6FQqFQKBSKHcXuMTYUhMKjiFUAw4Bkbo0oglmECC/0TI2v8Mabwlu2TwjIB4+RPE4kBiSqC+yAWECLMXTMBEBAOWeWAmUMJHPkBpXawyq2mxYYG6SBx6q2TZcs62j5B7HwgrmAByuX+v5k4VgMs2TBaMTYQKTagKnpPF3h7EJQHYk+7dSdiKOZ6+jjMw6/nj4kolAcEQndiIi+seeYhuMZJw3kkHCwZE2XXztcb6yOETr/gOPAP+TyCUREE74gX1/ddvbmTmR7eu7G0HLY9CjUjqTxGZbYlMy1mXMI/zIwThCFgyXCleJ8g/Sx1rFVn+RzQER069AZ++Cmuw6L28zYnHFCO+PGFt9XSGCIew+MFti2RxeuzQfR3MHYnK84CeE6fVzEtyBCv9dgyvg+9uLeUtI6hMOvIOLlMhhI948keev4jwnzSEOdfeK/eDG+QTyMv+s45BksDpJEIuVDKEuSzocokCPeHk6TCEWO2QspUgZDhDIumYA3Qi1ZHRVOk/8uBK14fnVxuHctCGGIsREi3prd+LNPnCeTEQ6VVKgwKNsUwexEEEhRPFwpiTOY2LXyf0sQLceMTVmM/jxxnQpw7uSLjUKhUCgUiu1xnTIP796LjbE+UVissQGLM2HNgte8SHYk1rMgzBtkAq+WGpESv/Sm69PM86rWM0VeY1MP/w32YSfqlz+vV2ko7JzDo8FVxP5SX6CvwSrW7e98uvkh/QE6oWSucdIyaGo6Pvj+0jEb0EA0C+hMormtfYYut/Gh50gMCDYparPnJv/qnguBBlODUOgTHsdHk8Aq7E2c0WMwTumCPoyDiBou2QBmaXnu5njvwjEfX5s7VmZscu3WO3Ona/naiStBsHzoBr730PU5fRDpZh6lYgFjkZTQHXtxEib9PodZgxF6Y3yXiIjusAgGTORhVID05tR9/sahu/eWXCpifhsaGL4u0Ug5kbkAACAASURBVNzXzNh0Qoti+D67uHCszP1m3/825hIKYBWWKJcAbUSij3JbXE/c257FAEkV612EnsWXBkHSOmZybFxcNUsal9pJWIua/gN6E7YbyaP8vr4SglxarWelFCpFMONvni0iwdiU2AvBFmV/syV7lVBzqf+I2YugdZEPTWG3tM9vZehz+IyEdV77AjuyEGWJuakxNDLsm3LyRjJFNg55xzNYnGPJFBU1NtKgSADYR//ZlHQ3imeH3XuxUSgUCoVCsTXsNcs8vHMvNsbE+pawssYqGyUVkDBPJk+KVwOe8JEaG7HfFlYQXt/D9woStrUiMSBRYZXhEwOKcVCUHHCNxF1crA6lHNqcTfDlI0aIVnH7e1EuIZ57bRWL1TmisIiCpqbj5dExa1wWC1l8s3CiODIHK/ZuxvPhhG7r/dBmtu9YitdmrmzB6+NH7nvrtDYtL6/iJHK+lAWSi3FQj09EF4iOkDSOI3jWJ278948dA/X2hKON+vzP4t0zx9TcfeDYotEDd8yUqy/s3QvXZXofyQh5DFwyYD1zbZa3Qv/v3rlJRES/f/hxN9eRm+snuawDorFWET3ii4JO3XZ1yAntbgodTURFeEaF60z6e5LZkG7OiSBNYJNaZh89u7BEdJygxeL+vUaE7aDJuJCUEpcOxUJXabhPB21NrJOqsBc+CWJ86Egcw+PFIX5ekYYHyQC9rqcSNbOV3iT7u4+uh0/YZlO7vDUDEUI1NiFm0KoRO4KNSZIfVkobFDU2NQ2Kt5s/C3qhcfERYF7TWPiPtzbXDZFJROG53Ql2KmHgRQFOudsUNDaBncptluwSFRg5xTPFzr3YKBQKhUKhuBxUPPycYIz5FSL6LiL60Fr7x7drY0MkUsTY+CKY/JsprBiIRASSKKWAFT0iEkqrJm+baRZoexBBIvPnEEWrCsHUYEdcCBKrbOh78KYvU6GPoqQPI9ZCkGdsmJmAjmWM+YX22ZhEwU+USyAietg6xgbiszOOkkHkSFta1WIxNsEy0m3XXKyyY8am2wvzuMmFJe9MnMbmTssFJ1kUc9I7uiFmLzBew/qZETM0iEwaXUSsAiKxuDkKcS4O3Fzfn3Kpgj5lDoiI7p46Vmd97Oa+/8j1NXvgxj+7G1ik0V03bizdGo62Wu+587i8Ga7lyS332+8dvu7mysVH701dcUrkG/lwGYpgIlLONCjAyffcEa8KOfKsie4rnDJoYPpJqi+zHNm2CllwaA3mwfodrg0YlAI70k/TZb7XfZXynGAlvS4vfftVGpnk2qT9+AhA/J7kl0m1ZzUtREzJeialS/8OM7vxkKBxkVFYmX4tGluW0yVlp+KoqIx0rmpsooNkVUcx51KRysHopHg+obtg2zMqKSuWRkWxGZ+/Ju2s9Mj2zMkmhqgAmceml0U3iepzrWmI4jFUbJeKbfo5F6I0nw80j83zxH9JRN95xWNQKBQKhULxkuBKGRtr7T8yxrxxlWNQKBQKheJlh7qiXlQYS01jfYj1NOLa4RrCVrqivNY3dkUh9T1cUhWhoInuB4SaI6wcIuVJk4qHi1yYcEV58XDkMvDetUI5B9et6yQOdfeVpsccau7LJPAWdHGBdvVjYXsrFoieLoMratY6MSyozDknassqH0cenI7LOEAsjHDv9R5Cn91x/SxckMOpc8McseJ3IsKu5zwxJLUjIjqfu88jrlY9PuXrccIuw5PQBxIW2obLLxzwmA7dGM/2QuI/ADO8OHd22rOG7bn9k1MWqz8KKmXzyCXVs2u+R/jEzw7c+Ge3gzB7ecvt+9qhEyfjut7dd64oCOXvLUIY9n0OT0fZDSRXRDkDn2wxFosLgSmE5p7qh9C4y91w/iqjPwiDE7cPpf2hgj2aSrdQDJRfkCWcfbK/yC3Tp4dkyfe2CPf2bhiIWKM/8B7ifVFl2x9REtJW7BXnimOFkLX39ijdEuVusEqCvuT/LSHiDW4Z7r4XdonCXKULvuR+Ey6aLBih8H+oFA+HcG90KrYRtgnzztoI0XBmd8imONexdzrMXaishV2buDgLthXPDC/8i40x5gtE9AUiotGrN694NAqFQqFQ7BYsabj3CwVr7ZeJ6MtERHvf+knbNDYT8BIFIXEjXr2NEKbFIclI8++TuRVYHaKU/YFoF2wR7IJF8gUOS5GLPrw7XVHHYbnNOj1WvuGPCoU/Z2O2zYwNRJNZcb44rFyEmiN5HUJuz5aBVZiNpsncfAFOnBafJC0Kb0RI9R6zCnynIczbJ+abhJMNATbOaccnEaLh+2vHYnw0P/RtFufMgpy4YyfHfH045Lo9DaJeO+VCg1x2YXHm2rQXLJzmkOfFJAhoEUrfcyK7FkUvfWoACNDDybUrZ9su3ba5cEzU6Mx9n5yE/sfHnCDvoaNb/mjkEvU9nLs54947X4Trcc7sUccJBn3hVBYE/+GPfpGuCp/72Z9xHyor4AS49IId8ZAsRqE/KQy2hZV1xq6IBHHxUK1I8y+fBTIMON7Xy3IqgjWJDfn/ZDAWyRCViAjBhgQWju0mcxcsAuyKUOtEzCoLfsrkhAUBLeYsi27642z+2YryAmZItCznXLnupRMW2KnUXvxMNZKRQ1vkusQ5SdiqynWu2HU20/tK8Wzxwr/YKBQKhUKheALYshThZcVVh3v/t0T0Z4joVWPMu0T0H1lrf3moTdP0nrWIw73B1DS1hHwFNsYzJ15jI2sewGboE9qajLFB+HchMWDGCPlimzYdB4XVt1wpBm1Pao8oMCoth313WWKyfO4+xH2VbtcLTs+/CrcGCnBCywMfshGhvHF5BKTwxzn1x7AOpEO48SguZMqMEC89Tzib3JK/v7t0bMZHFwfB0CmPjSOsp8ecvO7YaV7MWdC+WA5hbw+55ACzUzKkN77e0Ex5NmwiQulnPL9ZYFTaSVTZk4ioTe+rWDrUcnmH9sSN5WLEyQ8v3Bi9v34eluOGi2m+8wM/Si8a3vqhMlv0xt/9mWyfkYyJWMyaoYKQgqHxfSQaG8Em+B9463Uu0W+1MG+vc8FqPYwp0y9VCvCmodW8T5YVEN+T9pm+JN0m88woCLAVIlFcib2osVTxTqmxqTEnSVJC3mWF7TrJvXHO2xTBhIFe2iXK9D0ZWVVg/TbpiaSeiSjSSl4hY3OdakVdabi3tfbfttZ+wlo7ttZ+etNLjUKhUCgUit2AMeY7jTG/b4z5qjHmJwq//xljzCNjzP/D/37yadjdKVeUMS4KaVxI0NeUsukRhVW40EQknyvKeLy1p4yN62if0/rLKCVoIkwSQoAxMEOzwjH8PY6KkhECvnwEyglElSYZ0KaMmP1Yj9IlkGSMiIjQTccSFGhsUBFwFaW3X/EqD6wRtg2iccBiTMOc1zPY5NUL8o9NoPvJeVGUjXi0dkzNN1onFl9wvv73Llzk0IPTECHUnqTRUONTdy7MuaNCzDxobKhlhkZqgzihoZnxtZ0tfJM91i+d8JyPj7gkxA3Hysw52d7kToioMt0tts3RZHtcUoGjopKipPDl4/68YE0SMzRvf/+Lx8o8Dt7+9+q6nzd+/m+5DxnLILbxT0LrUmZsaLBfK1kZokzXU2Mv0oRtKUMjSzhIu0QFxkTaHdCb1CKSSiUPvD1Bi1lpN7JdjRCK+8vYoprdiBGSY5HnYMBVUmVQijRPaqdWdDPteNhuqaSCrUZF5X0Xz/dzhCV67uHexpiWiH6BiP4sEb1LRL9ljPl1a+3vikP/sbX2u56m7atO0KdQKBQKheLlw3cQ0VettW9aa5dE9A+J6Lufh+GdYmyIbFVjA4C5yfy4pTw2vYygEp15xiaP3JnGdRCoUMohdneL4prNOl1SJCySWCGCLYK+BYxNzFDtjxwrgdIKMiIi5MvJl0StYGwMR/90URVP+KZ98VE+H+OJOwfrKVL7hxUBPmeFMaXuIFrBnK8co3F3meZwQd6aDy5cyYP5edCwTDmiaXzB53/J23V+b9jxiMfGuXRQiPPAtTk4cnqcT9049m1ujN2+Y85x8xafizNmVFpmWNpViHTqx26c7YI1TxyNNb/ttovbYc4oYIlSBG//9ZeDobkM3hZaoZIex0OuoIU2ItaZVFf1vognvse/iWPkKldG5VCBManZLbC4MhpqkL0Q/fkorBKrICFZjEIkWi1CqKQlkfq9LEqpgD5jTrjJUASduAyyXMWQXSs0k1neHCow5JvsEuWsXbUYZsyGbWalni2eSUmFV40xvx19/zJHMQOfIqKvRd/fJaI/VejnTxtj/gkRfZ2IftRa+5UnHdiOvdgoFAqFQqG4LJ5BVNRda+23D/w+4Cj0+L+J6LPW2lNjzJ8nov+BiD7/pAPbqRcbQ0RtEzQ2MWuBYoHYomhkluU3YUfKV1qugOKCk4iKGjdr2cy1sfVVARgTs2ZWoaB9CY24Da9IkP13D9qeaJkzg8amkYIA7orJi3aVz7djpgbTwbZ0ZnzeICYnZhM3poup299Pw7IGEVKeCZLsGAiVZVQQkrMII08NNDfnHJV1/9xpb+xFsNMs0/6gW7IjPmYUMU+HjulZsk7G15W85Tr51M1HRET0x2983bf52MSxNw+4kOWIDf0ua5Aulgc8vzAPRIS1nMUZ52Jxy+1f3Ann4qtf+hFSpJB6nM/9XGBwfH4RkdlY/u4OqjzJpXbHFliLDWzCkOaixiIk/7FsYpQHdEVZ9l3PYsQ0sdR9CMagxCAM2E7sU64rqjE1SVdZNJQYQ2nOsuAnzrE8F6U8NsJumaUqN881PTkTWNPYeHsFjU0xwu/lxbtE9E3R90+TY2U8rLXH0effMMb8PWPMq9bau09ieKdebBQKhUKhUFweV1Ar6reI6PPGmM8R0XtE9D1E9JfiA4wxHyeiD6y11hjzHeReJe89qWF9sVEoFAqFQvFUYa1dG2N+gIh+k4haIvoVa+1XjDHfx7//IhH9m0T014wxayK6IKLvsfbJnWa79WJjXKhxSTwsg6B9tLRMihdrSuGWassCQVC7oza4eOAGm7F4GGn/pWg5phy9aFiKh/H7AO2MBHFwN8kwc6KQLDAr/Am2GG6myBWF89JOEHJeoYkpFPqEOwzbCy49cAJX1CT0L8s5yOvQwkW1CJw+Ck3em7pwbpRuWLJLar509sw6tJEpz/uxO7Y5cGJfOwrHLm67/ucs3l3edo1v3T4jIqI/dvN9IiL6jsM3fZvXWseU3utY0Mzn/WTl+v/DuRvTPCrMiXBuXG8kJfz9n/xhUlweb/1gHiruSzcIDIU8VzMBlsTDm1xQJSHtBjFpGmac2hssK+DHIMTC0gVVEN1m9oQ7KHlO1dJelNxN0h3jjxEusMKc7YaQ+gTivG8bZu4OgliZv5fcQbXzLh+DxcSPfKhwh5US9Mnr/bxh7ZUwNmSt/Q0i+g2x7xejzz9PRD//tO3u1ouNQqFQKBSKS0OLYL6gMGSpMZbagddeXLxeFB0zvmxCxFpgBYFU33J1wG/mbRzuzbUHxgj3tukp7AviYV/YEgUT11iq5KsFuRqDIBhh3odtKBEAjHw5B8yL5ywYomaZzz2MLd0fr7ymI3fQ4dglrgODc752YtwHUze286iUgGdKMCYIqJk1ahcssD2PRLdc1PF4ymwLriVv1whBj88XN19zcsDVkevDMgvXT0L/8ztpuHV3x7Ffn7n5kIiI/vmDd4mI6J+bfODbvMYn9aP2gohCQc63D14hIqIPbjgm5yQSNFuT3hNvflEFwk8bKN0QC4s3osJImCjlgBcPS9bCHyC2tIWQVqT2J6I8vFuwF0kXNdYi22+zNnIMWchzwlZV/uOTId3FsVROWPy3WisfITW48bmtzd0XHC2Ihn0jW7SbBHbUhMuZeDj6TTJym0Tqse0rC/e+XtipFxuFQqFQKBSXhxbBfIHRGutZjJi5WVindUChL1xEGe6drFCsWCJIPzpvxxFjA9uNSMQXws3BlsQaG+hYUtaomFBL6HvGIyTmc+zCzORh5mBQPBnlw8t5C5ZkGU3eMzZpWLwfRjTnvZFjZG6MHGMDhuiMGZoJj/G8MJ/A1PAO3vZ8540uwnnqOOndYuGuJbRNsrCpHYXv0K+sDlg3w0U8uymH/UdlHha33b7lTdd+/4Zjvz5zcJ+IiN4YuwjDT7aBfTlsZjxsp8M5ahxzAwbNh9hHJSJ6TlgoE88pnj6gvxlM6uclKDL+Wmzl5xiyaZN/DqHINfYi13aYjL0YcBdUWITh8gL8PPQPKjGWgsZGsrYWz9R4zlm4dXkMscamxhbJEOj4FBh5bmvamuLchd1CuHctSaBMPJiEuotndPW2Ktm5XuHeV4ade7FRKBQKhUJxOVyFePiqsFMvNoajooY0NkjqBo2NTIJnBlZnearsXGODSCw5hqDtSfUtRHEBTowFg0qZIqLchzzhUgqIfILGJxaCIVLHCrYI5A4S8zXL/Lz5sWDjmaK8ZMONkWMrkKjvZDTjYyHmiR3qvIttjxbSjtt2e2EeqwtmvZbMunDpg3aEIp/OznIaTu56H0Up05UvNDx9kP3Q8ojbHLn2r+w5Bur2+JyIiGYmLzC6sjxXPk8frV1Wvw85ieDZhYuGMvPA8ihT8/yBpH4JcyPZl9r3QqkDWarByEPje30Te1BigWSUTCXXW2Izi9Sx6TaJihId1TQ9W7FVBYZ5E1tU+k9UsFJZJNjQ3GXhT5kcr5Sgr1KmwhQ1NuX/9EtRUTWmJmtUonmuyB1kyVyrFxstgqlQKBQKheKlwU4xNkRObyE1F0Qhn4xkTvJU5YVXZr8KEL5kfjOfNIEh8IyNEKVA2yOZIqKIqUE0FFgSBPkUGBtoNlD8ct9rbNx2ToGKgG1fRsKXR0ijobx9oiAr6tLviMbyLAwRHY3AbJzxsF1/91rHWnhGqzBnlDwYzUVUGjvP1/OIeeJxd2Cc+CROx2u2476DpSGKyiNwiQOsArsFDym6w7s9tj1KmTiwfA97lz/ngy4kvuyISygsP05ERP/kzGUI/8P7r7q2H3GU1PcrS/MiIC7HENibWtRKOQLK7fQHxYfmOrzCbzU7SfdZhNDAUCu2B7U13rbYSk1PXzhWDrvAWljJFkmzJSYEf9cij00ekRRFX2Vz3hCBVpqAZ2rEuY5/q0VFlaRPImePzB/mC42W2LArFPBeI+2wMjYKhUKhUCheHuwUY4M8NiX4qCQo+CsrobIPW0RFiRVKrLEBU9PwcgNMEVb9nm2INDYyMksEVKWrP5ARzNhMOePwlGmYCXc8j+bRWcnYuO+IRPLMzSoMKs7IG48BEUd7k6A3ORq76KGbnMsFgN4nRIjl7Au0NaNzLvzpI8IwxjARyR6BrTqcpNmW4xXSfc5ts/bFPFmXI7NJU3RuufkFZzJ+f3GTiIj+YOJYmeN+5tuc905D8ztnnyYiov/rw88SEdHD953W5h1lal5YFHU3Eao5TJKD0q9FdmaT5gIoZDi+zBiC7VTzUizGKHPnSD1Lrchn1NRHJJUkI5muqEL3lPLkCNuDp03Yxt+wZIxM3IlkTLJorHBwLQVNHiEbs0iVgUtmrqDdGsyy/Cxhr5d4WBkbhUKhUCgULw12irFRKBQKhULxGLhGIpude7FpKlenY85w3aeJ8kaScu0L7SuiPMP0Y+z+kvbhAlvz1rL9ppCgzwhVWRDGhWOleBjh3mP27TQFLnPt3WCwzWP1Cfq4TReNXSS88hTv2O2fjUIiwJsc5n2Dk9Mt+WCcF+kCIyJqWbwL0fDoAspg970fp66jGIYT8+2zO+yVmRMtT3hCcFERES1WXIaBSxp0HDION1yczxDjQ2j28akT/r41e4WH5n6/Nf6Yb/Nw5Y75g4evERHRB+/eJiKid77wY9m4FS8mNrmkSqiy9s2G35NO0m2SnC0L9047TJ4yVdHwFr4cmYhPinoH/rOTLvlyuLcYg+yvOOdy6LkMsU9sXqb4pZyAvA5biHqrIfaxzapYPO/TvABuIHVFKRQKhUKhUOwgdo6xMcZ6AW8XvZdBNAz2QiaBMiWmhnCoDPPGltkGk4uHez6o83aZKepTpiixnYUUCrsUmJOGO/AshU8MyIxH9Pa94g46Zj8aiIc5vNv4Ug5hUJbHD7YIJQ4w8P3x0h+Lwpv7DdMwLKj19pkpMsuYseGQ6jmLrZmxwWppvR5hIGFMfB5GEzfXWzPHEH1q5opUoozBNAq/P126sZzOuHDmCP2mAmo3f9hhwXfr2r5Ht4iI6GThQuhn40DznM7dMcf3DoiI6J2/+uOk2E1UmZttFrIVAW+xvfg7NwV2pFZss9intCmfUyVI1qAW8hw9R0z2QRySFMHcwBaVWAthu5g0VfY5dN7jY4sqaGxFiHuJsan1V2JnKnP2hLxkpApjuQpcp1pRytgoFAqFQqF4abBTjI2hvCAi0GXh3vzDNmGUlQR9QWNTj9HzGhtminyYeRzu7X27IpFTW2dsDGtsUGRxHHdIKVu1ZLoFyQFHPkEf99WBqog1NiaxB8bGjJmxGQXG5qhxjA1KDsyNC5NeceNVx5qbhLFx29GFG3czR+x5msAwnnvvuqW9mTv29dkJERF9auoYG88YRXhvfJO7Zc0QSirw8MdnUTg56lUiNJzHveBKmvfn7iQ043Cuey6q+c5f+VJmW7GbkMxN8kiprcZrK/kCMoZGaktKvw3Yz5LUSRbBiAPjjiVbsYktifqT6ShKGpt6GYmCOCVLllqwXRkLZYn5thAHVVirwSSIsgsZ9p20kcKciv3o89CwnyUsXS+NzU692CgUCoVCobgkLF2pG+x5YydfbDw7Y0saG0QnVd7aYyCrE9gL8XZeioqSYwBrgaioLEW6+Bzb9QxRlEwOxS9hE/qeViYEtFFZAWYeekT9+KKb0Nbk44euqAdrxHdCU2BsZixUkWUkwBStODIp1rNAY9Ms3GDMwtFHdpzecv0ozL3n4paHM8fMvD49JiKiT09ciQMkJzyJEuih3AVWI2BqRi6QiqaPoog2TgbYnUF/g2vG15J1P/0snNt3/pom4HtZUdLc1PLMbQX5Z7ZNNE6tq9iuZCsqZQyKHfikdDyGbQpPerupPTybiKicFHDArvucak+qyQlLxFaNISokB5XlHOSctwmkyhmp0oniZ6eMaCtpbEr6HsUzw06+2CgUCoVCodge10k8vNMvNl307o2SBr1wjm7lihXaGqxG/Mt7sehmythAZ5IxRfJzZLdYUI8/N5zLZdyk2hogYWyYObGITgJj41dNBcamTaOhsG3ZLiKQiIK2BhFZK6Z3LjonillzAcrxItLYLDkaauErW7oxjfg6wX6o5Ul26trcRjTU9AEREX189CgZ+35z23/2BUCXbu6Tc/d9+sj1NbsXIpzahdu33kPkHObODOA41zwprheqDM2Qtqb2jBnQXNSeS6Voo5q2xg4+3GBH2K7kzUmGV9OXJCyS2CeZE//8jaKuhnQ9BbslO9VyBjEka7SFviUrV+HnXrC3SVtT0BBdlbbmumKnX2wUCoVCoVBsgWv0cqUvNgqFQqFQvNQwGhX1wsKkbqFYPAzXEES9QUDGTUsORoh4JeUK8fAAfwiXF1xgQTwMewPz8EnxIOANP8HD1EK4LF6z4QZa9eHSwQ1GUjxcES07A6lwGVW9RyPXeNoEF44MNZ9zXDZcUR2LcKexeHjJ7jyuKI6Q8yBWZhdSlOvP7LljX52dElFwQb3CpRzmFu7GcN0v1hwjzmUSxq4pzR7wPO6GiuTN3M2pPXJG+4kTIS9upr6nt37oi6S4PoCImKiUvK/idhh0hWC7QSRbslNyx8iQ5yy0utBdFt6dBjUMjknYyxKXlo6tlWgoueQ3nY/YDSfPRyWs3CaZPitjGpxzZa41u1E/QX4gxNGJeHjAtuKpY7debBQKhUKhUFwe1+ilaidfbMCWxOJh7JPi4exiRgUna2/l1fTdFId5j3gMzBT1abj3kHAsiJURbh79xrvAFknhMgo1Lmy4dEtQPlnK8rStbaLw+JYFs0I8PGLxcMzYhBIWrv+5dSzJ+dopf+3S2W9DhDg1YGzWgjbiMXTTXDw84cR8dybnRET0SuvolyMWUM+ZmZpHcz9ZOPalPXW/Ibx78sANpr13EgysUEiUw71fYbqIT0sXosgVimEM/SdxKaam2DR9JtSYmiHmxouFpaB1IOxYPjozQe026mdpv9BkILw7sRs32iQaHgqtloxKwXYtCWJROC2RhbiL76Ti4eeNnXyxUSgUCoVCsSVsnm/nZcZOv9gkCfoyxqbcJr62VibKE2/nZuA+yIpfSv9qMZ8T/+YTA3KTSGODfY1IwhXCy6GxCY28xkamS5camzgRIPQ9qBmJEg4tl0AoTAC2F9DYQN+CopuRxqZZiTIOONdjN9ZuwqzbLNjZZ5HOq2MwNa6Uw9jbd20frA58m5MLx7qMTznM+9iNf3TPZeizD9NQcSIic+OA5846Ke7uzS/+SHas4nohS9rnKdRtYoXFoVswNtn/NTK0m6jAKEPckTaNQ6trpRO2Yg5qBSAvoSsqFn3cYLvEjmxkTErVbqSuaIghkmQOzn+bam5StgqDE9uuwtxEn6/Pq8XVYqdfbBQKhUKhUGyBa+QO28kXGxkBRRRW81asGAZXKLUVycBrNRLC9aKkQieXXkM+bKGxsZHux3LRSxmR5e0RF26MoqLAGpku1dh4wE7cn0/Ql0ZFjZEYMIqEwnmeM3dyzsKYi5X7blZsdx0syOSAFpqeCUeRTREVFTM2Thdzs3WRTDMeA0ZyZp3dexFjMz93+w44Gmry0OlozLHb0Z2ehTFNWVPTclLFAzem5e1r9BeveDxsQeNvpamR3dX0fU2BIajazTUd1YigrHHhc01fMhAZlLFHA8y1tDdYcLIaDSU7iX+TzAlVxyLnmBX8lOUsksZldso/h+NCyC9EVNT14Ys0x6pCoVAoFIqXBrvF2NhIy0KBPSEKjEYWFSUR9lTTzgAAIABJREFUCWe8xgYZ9vG2zr8PFb/sBHPjo6IGVkheY+OLysFQdAxy6Ei7fknBGpuopMJQRFZsN9HYtGLu7FNum9xp7QtvsrbmvHMsCYpvNszYNMkKBSsd6Ik4B80EzA0POWJsDseOsdlvXBHMMS9zuHYlPewcU3NvERib/syNaXzKjNOJ68Oeucgquw7RXc2ha9cdOuOLm25sX/1x1dYoUpQKZMYoMhGXWY3XmJpS9I/U1tTsFoKWjGATBgtObjOWqm3BkpSeg/KZLIdUKPJZ1dhkOqZwgPFaF8EabRMJNvBsziCirgJTk26L470KXCNiWhkbhUKhUCgULw12i7Fh+Dw2cVSUZ1CGWYtStEH4bXMERO/ZCzAnkikaGLjU2AyMKbMLO/w91th0tayiQJN3an00FI6pm18Jlggam8XK7Tci47EbAxgbzgTdQlPDmiTW2NhxGOzeyEVFzTi8CrqlFQ/uo/URERHdmwfGpj1zv43POO/POTM2yyipDsa0t0dERMtbTmszv3N9fM6Kp4PHykdSus0ehx3JxvIYEUhSChh/l/lyBpmasqYw6H0GnofyJ8mWxMtt7MuSgQn7MdGcFQFOj016ktoaPA8RFdUOiDWzqCg+tJD5/YXIY/MijOE5YSdfbBQKhUKhUGwJS1sJ4F8WqCtKoVAoFArFS4OdZGz6wpsn3BYZ27ZNcqxLhHvD/dWbVDyMMHNTGFsuykOYNyXb2HZWSkEMat0XxMMVgZq3G4mH+4p4uCyYdsfK4pfrjsPMCwkBgzAbLicWD4/ZJYVSCpPQaNJykUrmleGCOme799eHRET06CLUPmgvXH+jC+5n6dxYll1hPsSbiOyRc2EtbrmxLO+UsnspFDkGXQnbUvwFoe5lkuFlzxbx956McdOzIHSa/WbFcreYEBD9wDbcPcINUzxvchpCsOvdP0R52RsJ2ImEuj7UXLqEhtxiwvWUuaDicyLDu2vi4dgV9QI8akp1oF9WKGOjUCgUCoXipcFOMjZAX3iN3yjiTcRyBfFuyU6RIUpFw31lhZTAJ8oTIrnSSk7YXvXppVoVyknkxTb5A1iTOBGg3yeHaJM+nS1nG8Uvwdh0YGwgGi7NHeJhUUqBSRgyo4ixiePFiWjOS7iH/T4RhcR85/NQOROMTbtghmbt+kDZCgiGiYi6W66f+W3325s/rGHeimEg7PtzP1cO+x7EkFD3SVALdR5kR2QGvcKY/POoImQuMEIyOaCRouFi2ovUjky3sVUBUMEUJYwNZ3jwzMkAa5WFd/tgCsHUJJ2IwVREw3EwhT8PVylzuUaMzU6/2CgUCoVCodgC10g8vFMvNpYck1AqqdBL6kGi4LvemNa8AM/QQGPDjQYrp4p+DZgTH14eH5u+Vvswdt4Gpiiaey200tstMDa11OGMLlo2galB8csls0c9MzZNyZ8O2yMeL7assUG4eRMxNuMGJRSaxO7DzjEt95mxWS3CbbvPUd2y6KaZOW2NOQyh4SsO817cuj5/4IrngA23U/HRsCm8eygBoNR4bLMSz3Qt+X6bsRQYamGQPjkdH4NnkNTbldqCqZHPIKlvKYzb75b24lQTkqnBWEr/B4jnn7fdinMRQzDUnjWS2p4uH7xtrhFtcoXYqRcbhUKhUCgUl8cLkUvnOWGnX2yS8gqIiroM3WbSY2tNYztgExpmBjppdzBBn9DYbMEQebs+QZ87eB0xNrbSkWRqbFxSQa6WKtoeopCYD1uUUsAxjT82jrAAMwONDZgbSralKKyQCNAxLGBsHi1dNFS/CBFhnMuPTIflE4+Fo6Hs0b4/dnmD53HjGv2FK64MQxFI9UZb/CaZmm30LH4MvL+kHdkUgRT3L5gaz5gMFL/MGPIaU1MqAFrRE8nyCW5M6bZqP/qc2UbCUjA5xax+/FXYKyXoy+ajeKbY6RcbhUKhUCgUG2BJxcMvOkpRShkqx9iYpTF26NBB26GEw+aIeUmoGO+8pnQrP1NghEIph5QtIdqCpfKrtJyxqSFlw1jzwhqbNTMqthJ1EPcPm71nbMAesRam4HNesgDnjBmbR8zYnK04J806DL5BBARSbYzc2MyeY3fWhyGPzfLQtfvqlzQaSnE5vPWDm6Ojqn+GpT+QS2lrtmdoqnb8M0Boe+IIpAqjJCOfiAb0JHJsJXZEjKXGlrBxtl3uvxSBJPcVmRp0IyOx8Fxq0Qn6jJ6dfmhlLU/GYkX92EJ5G8XTx06+2CgUCoVCodgWRqOidgVx5E5WBBPYRsdSik7a0nbIY7N9Y8ngXEZjI+0SBf9vtuAasGNF9mPZOD63UmOz7tNsy+X8NWBmeLys7+lFzoq2wNjAzhkX2zzpHPtysXaMkVlF/vROtGfGxnL+nO5gHPrdvz5/2Ipni+HnSYVKKbWRh5YinYYyDG8an2RqGvF9gC0mmbE3nrTM+FvJNFwcqmRJ0O1QJFJ2njAORCSFn4LGRTzXwRDF/RhxUhswNqnBYpvKGEoamxfineIauaI087BCoVAoFIqXBjvN2CgUCoVCodgC14ix2bkXm0uFc5cwRL36/SJJXqF8QSlMmYjKN09NvLbFVPqsdMP2ri9/CARrSaryTXbzcG/sqxXdTFOVixBzNPEpy93GbBHujRIO87W7XZPEV8I2SjfAvbjei4qFhuoKCsVj4VIuqKFjN7mgiuULKm1L9qT7ueKSStuLPybp1Y+jsOFmkaHOsov8TzUTD9dCrImiMGvfTzqxUmh3lrxPPq5K/wfAdi2owRY/umOFkDlz4UVmrlMumavEzr3YKBQKhUKhuCSu0UvVTr7YbBNifRk8Dgkkw74vdc9swdxsYmT6UoK+Sn9BPFzoU+wKjFDoH8JlLx72sdz1MYZVmRHbelOEtq9sm+xHKYe1LLpJ0QrIM0QsHubw8m4azSPUzlQonh4uw9RIVJiaJLlbLTFfursSVi4Ymk2FLmODMjggYWzSsRjBUpiUWClDlDMwksGJ+rci3HvwfAkWp5A/NJpIaggMsmTku9KzMxuDYJFKLM9VvVxYekEUzM8HKh5WKBQKhULx1GGM+U5jzO8bY75qjPmJwu/GGPNz/Pv/a4z5k0/D7k4yNkBS6uApvI1uFXaNkgrUbTiygAqTMgSva5Hh5dssC/1qZICpqel/4jFUQumtXH2UfNdYaLXl/bHGphcFPrlaAi24/gLCzOP06Zm+x5dwYF3QOAqL3+m7XfHC40keQUPJ94b0NzW78u9PsjmVNA9ElP59lcYYfc6YGpsfmo+tMpYmZU1cP8Nj8b9GzyY5puEyN+nWDJ2XGmp24u8vAFnyvPU9xpiWiH6BiP4sEb1LRL9ljPl1a+3vRof9a0T0ef73p4joP+PtE0EZG4VCoVAoFE8b30FEX7XWvmmtXRLRPySi7xbHfDcR/VfW4f8kolvGmE88qeGXbg27KWoqlueYxyBdXkTYbVYmRMmqYRNbVGKEZITWYCfiEJ8Q0G/rg+2EhgqaHlvRGLidaRSW5dINfXSH9y/d3a7YdWTynELoTrWUgkl/H2ZJRGcljc0GJmgwaaBgeQwK/V5iLKYwJt9P7VEjNTelfZdibjAW8XOJjZZf5aktna+rxNMfw6vGmN+Ovn/ZWvvl6PuniOhr0fd3KWdjSsd8iojef5KBXemj3hjznUT0s0TUEtEvWWt/6irHo1AoFAqFYivctdZ++8DvpdfRbRypT/wKdmUvNlv6357xIJ6bpaeG7jG8h08qP5IMyqWAHDpbaHkAzxYJu9vkMMqZoehHTSKheJ4YisaRGDpmk4Zj6M/iElq6y/yN1g6ReWwu9V9UKffNtu0f80/bFumVp4RSl9frEfQuEX1T9P3TRPT1xzjm0tj4P5Yx5geMMbef1FAB2/jfFAqFQqFQPCGMfbr/tsBvEdHnjTGfM8ZMiOh7iOjXxTG/TkT/DkdH/YtE9Mha+0RuKKLtxMMfJ8em/CqHbj0tnqPmW0tgjPmCMea3jTG/vXp08ZRMKxQKhUKheFaw1q6J6AeI6DeJ6PeI6FettV8xxnyfMeb7+LDfIKI3ieirRPSfE9H3Pw3bG11R1tr/wBjzHxLRnyOiv0xEP2+M+VUi+mVr7R8+ge2tfGssRvoyEdGNb3v9ehF5CoVCoVA8DVxBgj5r7W+Qe3mJ9/1i9NkS0V9/2na3Ek+w8W/wvzUR3Sai/94Y89NPYPuZ+NYUCoVCoVBcX2xkbIwxP0hE30tEd4nol4jox6y1K2NMQ0R/QEQ//pi2vf+NiN4j53/7S4/Z1+NhB/mflvrNBwlsFQ46ZNNc3qZHL4xvMYiGD5J2SwUzJQyrDcM2+vEapRRXvAB4DOHsYNp/2e82t/NlQp4v8TdaG0L2J3ap8hJbm8/xmH/ahgf8xMWVi51vue95wNJO/n/3uNgmKupVIvoL1tp34p3W2t4Y812Pa9hauzbGwP/WEtGvWGu/8rj9KRQKhUKhqEBfbAKstT858NvvPYnxkv/tSbFpVf8k5MOLiq0KzhGlRew23ORN4a+g4ZPXbBMimSWrGmBSBFpvx21HvMW1Lc4T/Xa8Xbtts47Gv85aKRRXCtzLMky6VFagmqPuUsyNKGw59L+d+Cm240NIfGaG9O96cEzZWFLWxESGbTZegcKzz4gmVox1mzHJMPPku/ytNpZLJERVPF1oLlaFQqFQKF5yXKc0Xjv9YhOXlm+fwlXbpovH0bh4CMf0NvYwR9jteTlSYlSq9qwwnPwmtqUx8I+NGHAW+F9a1YCx6cr7Y992I5iaMTeaMtUyavjcNzH1JMosrPt4NzWrcKxRxkbxLHEZ7YvEoDaFmcoNzE1xLJ4xZSYCO/AYK2Xv2MCOxJ+RQzPT8Q2yI0JQhLFwWYYkL+cmpsZ/j1geTgrq2fmtmBreoE1ziYtYY4SUpbky7PSLjUKhUCgUii2gjM2LjeYpC2Ueh+wBewE241Iv59Jewb5kR/LfwznIdEXiq9S5DB3r5xX13woGBZqXoRPnbXJUVNjWm4J1g52ZWRER0bRx21HrGts2Gj5OPPrvOrbjjm0X0TyWUUOF4mkhE8rw/m0eCtCNCTamxFrg73frMgPx2BCdCAONpDyIqgKZkmYErCmKX0rmZhtdi2Rq8GyIhy8Z5YqOJSlu3KRNBktFZNoat+3lORiKVhM6IzBGCbl/iXIVzwzX6MXmCYoAKRQKhUKhULxY2DnGZptcJoMYULeH/elrdcxebGJShlcFA2OpQEYi1fQuxaHgkF6srmhzdFjcPxgUr/dpxNKqoBnyDE2X+tG91oa/l/JHwN5+syAior3WMTazkRPJ2LausTErNsD2RxedP7S9GJemqlBsjcGonxpz4xsXPnsmYrOO5lISHhGVCObDgB3xB9QGSNn4kzmDoZHsyDB5nIxJMjX4m7WlGWYsUsqOxIyNj9DifipN0wHiedShv7TRUJ4bH7iFMWBrC8dcEWNzifpOLwWUsVEoFAqFQvHSYOcYG4VCoVAoFJfENUqms9MvNnHoddVFIyjZImyFOt7CdhDbbt9Y5rfbKsycfUcQtcX2qgn6Buz4RHl9uXF8buEa8uLhJk2YV+TGvSuKzxO7pBr2DMEV1vd1V9RBsyQioqN2TkREeyPnkrLjKLSzFe3XLB5eObdVe7YK/Z5PCwNVKC6P0t9siGKu+D4Soa7Yhl7yLpC4TriVtnpeVVxQCP+OOXtbexD5sUa/wwXUCrdPReyboBdbHNsNCZorY5KCXYqDC9JzOeg+9CJrngfG4oMTCtn2KmPw9kvX+yqhriiFQqFQKBSK3cNOMjZbsSOVY+KQ58swJtJ2KSx601C8HckQDQiaEQINRmrcpEJeoi0E1SL02n0ebhL3jznOEHbNjIppC6tL0T9sNitmbtZgctw7dYmxmXAmvQMWD99sz933sftOozD4nu9gv4hixsZeOJanPZ2Ffk/3iYjoW//m3yYioq9+6UfygSsUBXzu535m4zEZWSwZnPiAGqNRWtlDKCtCk418fpQeA+I3n6jPC1yjZ0IjmAiJmN3JJpAGKAw9kqR4mDo5/Dw4oMbY+HnEKSBqtgeeU17ADPEwq6JNI859BDBc1mcDxZj4YMkmE12pglfFwwqFQqFQKBQ7iJ1kbICUVWAm4DKvpeI1vNY0Kd3AjmFoXlppdzB9N/zbVNwOofVLC/cuOiol6DPpKiroaHh/iHyOGBWMLbVXCvfGdtIK1qgwAW8bJQ48Y0PJdr1FuPctZmxuThwL00zDRHqO4PZaG1zThWtrTs79sZPjG67/Y03Up3j2yELDk0qNG/Qsg0jDootNan+aGfMRaVPwI5a7Q2ySD/MWzA0M+udKKXSb0mM6kZywZrOw32t8IgbaSuZJstMFhtx0gskaCj2XY/NMDcZUsINdmqDvuUAZG4VCoVAoFC8NdoqxMeRYAsmaEG2hdSn4oWs+6iEGJUuUhyiDoUayXzAopWgsmRxQ2Gl8srxo7vDpVqOixCqKCoyN0LrEUVEobbBq3FJkwnRLgxIHWKkkkRzc/xpMDW9XYHLY7Dq8W6/6NrENu2Bs7ozPiIhoPA3VLLsJ9zNGxjBeec1Zj2NC/+OHbt/0oSbqUzxFbNB0FJP6yUKQW2lthD2UWNhC15KVWQEBHJNIiGwSzI1tCh1jbvh7A3Pt+4CdWFeURiAZUeahxED5ZJwydyCS8EEzFJOweN6BCRpYvgeNDR8LvY/U9pROrtf58HXwjI1gj5N+6mN5prDXS2OzUy82CoVCoVAoHgP6YrMbaApXqtmkdYmbINfKhgteisJqZamDal6Kgr0BfYs0hf7HzTrZPy6UecjSZ3gfNlZT0W9dOpYwxDxPzpjpFTAoKHHQMmOzwmqpmLMCGhu3FGqXiJJyP9uIsVn2qfZlxsunW41jbF5hxmZ/tvTHnO25/ropr9xGrg+L1ex50Ni0D93n2YMDIiL65r/jIl3e/OEvFgauUBC98XfdPfJYC23JksQ/PcnKPWN38MyLc63IraQ80KQQi2QKjIO06w8RTI3XERYirCRzLVmrwlA9IyRZacGWUPLoEP3K/Djx89YzNeIZ2og+msK5Bdg2WCMfYRWN6SnXbVZswE6/2CgUCoVCodgCyti82CgzKDkrUd4R/VTR2AzdAJ6p8ZqXNAtvKXtnHgVVYW4i271YYbViUKMmKu5YKUop7dM69IFMwPApIypA2iUK5xZ5bMDYjJixWYqoAGcT/cMma4NYY9OCdFlGjE034iPdvjFrbY44A/Gd0SkREd3cm/s2x8zYrPd43BOnn8FKr5+HY82JY3ymD++4Q++r1kaxHUqpaDzKxEOhk1KTLRp7aYoVu434PfpNPgsq+sE46kdq8fD3HHLfbNaZVDMQR0PJ9T7yHESRWpKBl4+nLDor6sZrXWCnYB+/yUzDUjfYFybiC3HCHrfZJiJM8Uyxky82CoVCoVAotsd1Eg9ruLdCoVAoFIqXBjvJ2MAVlYR7yyKYmTumIBSWrqEtimGG0gZOUDsWIuLhBH18SJ/aK40ps4uEgHzANBITtzLcUI6hzztFuDVcUZ5+LtiGeBhh2PvsGpqO3f6zEdtvYwqZP0M8zP23CxYRL9ytZ1ahzcXauYbmnHXPu8DY7mujEyIiemV25tu8fcAC5gPHA/f7Lv67mXAceOSKshcXREQ04bDv2f2dvP0VV4hBl1S1Uf2nqktqC3GxD7E2hUGZSgdDguZe9DM0x+wEiLBvn+iu7o7xM/bf83n4MGzpcsqjLMJnuIRwKFxo6KswlpAskHdzf15HHYuTa+4x6YoqSAyuE2tyldAnu0KhUCgULzuu0UvVbr3YGFFGIbpSTcacVK5iUgRTCGiFBtcWhbQpc+KZlA0C3tge3uQlWxKPIY8oZKbIgCkqiIcrifqC3WjuFfFw1+feSdiesHh4n5W/KK3Qj93vfczYYAUFm8zcNEtseciL0OZ05ViW835KREQrpqD2WSh9q3VMzSvTwNg0B5w88NDdyqsj18fswBW8NFG4t11xYsFTZpwecVHMn+aimD+uRTEVDgjzrqFY+NX/uIWBCnPimZtC1UUjWJGscRyRLPaFsisDoc9ibP75V60qWWjk/+7zMcmSLxlz49tEgmaUNigTRDntE7qNkhCm7EtRPIzvMvy7UGQzGwMYc9+H2FIkkL6qBH3XDLv1YqNQKBQKheJysNfLDbaTLzalkgpgMGQxysFEWLUw7y00NmBqYLetraIK9rzGRmzdZ4RLph1Ie7HGZsSMDdKPW0m6FOwg3LvhEHBoblada7yKnMo4z0jQB43N3th9t8zY2FHM2PiBu/5RBHPpxj9asE4qYmzOF45tedTtERHRnMdwhx3fB8bZRaI+IqLZvtu3OnRtlrfcLT29cejssa4mQcdFPc/cmCYPdvLPQPE8MVgyxcSbrZibTMbi9SD8Ndam1EqmePPQtcSUjWBoas/B0hhFeLT/W46PrTHTmFcj7JeGIMl1r8uJDqk9k0XIdSkRoGePhPYwIZH68m+y1EKiU/RFQiXlxBs8OqNraLu0yZXgGr3YaFSUQqFQKBSKlwY7vVRtCmUFNkUnJYEDPlEetvyD9zHXbUN3ArbE24W/dUhjI7Q2Ua49v6+vFKX0Gpuo0Zi1LiQS5WXMTRcmZNZgbPg7a2zWHUcXFSYA21Mk6htxXQSOiuqjfHeyKCXmLEsrtPMwyIuF6+DuyrEtJ/2MiIhWrOkBW3U7YmyO9lyE091D19/ihutv7xVXNmG8COUXiDU2fdvy3Fm3xN198884rc2bX1StzXVFpq3ZZoktjikVc+QfNjX1OpO4cKNkLbzWppF2o7/vppIob5uoLqmt2YLJzhmUS0Rq4edSYU7sq9kusUmyeKcR16OksRFTxfPQoAhnNCgbwrjSLYlj4+uDMXXD5+CZQhkbhUKhUCgUit3DTjI2Po9NEhUlGZsKcxP7SmWKbfn2PhAVJXO7DEUm5cUpheYlVs8L21lpBT5gaoLGZgL2RuazadK2po8Yro6ZJu4G2zVrbBZ9uDU67hDn22ttRo4NMRNnv5sEyqaf8Ipx5NqGqAM+X6yxaSJCZT537e8vXbTSvc4xN6+1gaFx9sPcj6aOsfng0I1hcdONe3mb89ksjsKc58zYHLrfemiCkGMnpLxRKIYxuPBGlMwAU1DZl0UIUSEaSkb3lPQuImLHP4TwzJHRUdFnWWoGz6KYfAkMSmUC3n7UxkdopU2yLhI7MoqrQjsUyz0I5sTbjR/OuU3XhvdvYozixo2wFz1/S/qe5wlD10s8rIyNQqFQKBSKlwa7xdjYlMGINTYyUkoWgbMl/67I/Ou782/p9fdrMEPQuozQuJb9N0afRiTFGhvklYHGphfv+EFjE1gLr7HxGYB5DPK1NV52+Dw2iIpyB6/XrnHM2KxsmsRBFsNsx27uscamA2MzbnlMjbCHqKhoeBfu2Ltzx9R8Y32TiIheGx0n9uPr7nU+M3cOkM9mfrtlO3v+2Hbh2q333G/LQzemfpJ0T5/72aCzeOuHvkiKlxuDOWtkquFtltyCovHMTSlyR9rhH5JDcazMDFyLTCIqFKUULEIetJRraZDTBblk4ky6rTjWMyaSni58Hhh30idFU63pHn325ah5Ve8zYNdW5jqkL0K3nhXjjdT4UBSxdpV5bK4RY7NbLzYKhUKhUCguB6uuKIVCoVAoFIqdxE4zNrF4GC4hKR4ecgkFFxTEvPjO9C0fVwp9DuHXnOyN3UGmlExLsrNSrFwoqdAjUV5fyuWdllSYsFvKtEjUx92LJHnJGHyCPkq2S7Z70QW/0ty6zx1PBO4wuKJGYxYPT8P1gCuqn7pbrOEQa7gEfYLASDxsFs72g7lzH723uE1EwRU14TmfR74jJExsJ27u3T6Hfd9k11fkUmtWlsfpxrC45Y5Z3qCkrZ1co6WNIkEWfm2yD3xg3GhTr7nfp1ZM0z974pBnYdNHGyMxYDGOHHbgqkmT1YUgCiObRP3j2HR/vI+y84P9hedvk86tFg6fiId7+ZtwrZUgxbz+XGxu4t2GPhWIOOlJR9IlWLa7yfZzwzV6rCljo1AoFAqF4qXBTjI2rS94GYUvi4R5GXNSFIwh3LosIkY67SJjI0ocBPFwGnItP8d2PVMUJ85jxgQ219y4ozTkOmFsmC1qpHiYi1LKsG83VwiXIeZ1+/uVs3O+DqzInFXBnZgImKLx2G0vYvHwFIwNMzXM3FCT9gEBNRFRw4zN6dwVwfxg4aiUd8evEBHRfuOUxo/W+77NkhktlNLoeNhrl5+PFlGiQ6y+1vtuu7zJx9xmpuuIi2SOw7n97K/8TSIieuevfIkULxcgGh7SHviV/KVEw2UkLMyGY1PWQoqGRVjxUAdCLAymw5RKLchoclGU0sSJRLMkoIJuEXaT8foACyGqLgh1/eNbskZ+i3PRR23S8HgvpDYpEx/377/jPGXi4dKJQr/iOjSp3Xjfxgv/LHGNGJudfLFRKBQKhUKxPa6TeHinXmwsldkToojFQbikL1S2oUOKtTWpXzUwNnWPHVijUZNqbOIIaZ9jSryt+9DnWGPDqyKk3l73eVFKoqDxIQrMScOrFshKsEWodbJa8GwRz4MZG1tgbFDaAFobnA+pLzqLtClgbBBa3c5cW7kCjufOUeS04ER9H8xdcr1b41tEFDQ9YHKIiM5XnGwPoakNmBs2cBDpGvgUrrj8woqZmvaOY4JuHJ0TEdFsHELpwR599r9g5uYvK3Oz6xhkaioP/0AqSHqBsmeMTKRnpCbDHSQ7rtr33eNvBYyHbNOUxpQyyF6jUig4KZrkduPilLVw6CzUOppIk9reirzwOh8+uBcPkKE6FbUw7yJ7X+mioIP0h2a6Itk4/+2F0NpcA+zUi41CoVAoFIrHgDI2Ly6sNVGK/5LGJq0PbwuaF4mggOcdfpvqXIiCzqRBQkBZDLOosamEAUi7FBiUnhPmLZl2AWODyKQmWm6OuYN2BMbG8palXRcuAAAgAElEQVRXSNg2+UkwPlkg72C7YEKIiE47x9ic9469WAr2aIxorISxYW3QjDVDe3yrMUuFMcXnBOdhtXT9P0R0FDM2YKbuLw58mwfn7ph+7tqM1mm/SdLAGfd/g+fMTM2nXn1IRESfObpPRES3xhe+zcOV6/8Ppq8REdFnf+mniYjonb/646TYLVQT8dnK58IxYIJtfONWiyLiWKFroUh/AaZxaEx+n2At8NwoBU7KsYDoEOn+Y31IiAgSXXndSaxX43216KTCuZC2fSI7+SgoXI8QSSonJPdTzpiI/wsSqih7Jov+ZCFQ+TmC77aksZFsjuKZYudebBQKhUKhUFwClpSxeZEhSwwArde6IGIqZU78C30hQijX2vDuQlSUtO81Nijl0CCXTLQSQpSS9/mmKyMTFeb07A1rbJada7xi5qak9wFL1TJzwlURgtZmzG3aaOxYsWHOILpWbv98HW6NR1yW4LgP5QncWNyxKABqRzFj47brGbQ2o2SuYHT6wh1oOTLsfOnolntzx9CAlUOeGyKiszNHw5hz1vJwiYaoTmboF1FjXH7hxqFjZj534x4REf3Joz8iIqLXxw99G7BUr01OiYjo/+A5f/bL/ykREb3zhR/LDSleKAyWTKigtrCGtsTEjGy1k3Qb55vxuptKXpnkKeN1JqI/ySoUx5BGIuUlIuqiDxkhlJQtkDm/aqxVzFqA8RHRUf6ZLJ6/8T6bPZsxkOrwq1qbROdSYbayXGMlxqZmu8RWGZvvfM64TmSR5rFRKBQKhULx0mDnGBuFQqFQKBSXxDVibHbqxcaSqYZ7N1m4txDxDoTbGSEQk6nEuz4XD/cyYZ53B7GrJQ73FsmssijBgngYrqhF5y7RwrqtFO4SEbU8YOkS6sfsSmOhbj+O2lbGYFh8e7EMqtuTlXP3POqcCwhzXvQI/xb0NgUX0xph3/soceCOgYsKYySK3He8WfFJPF06IfOavx+z+4mIqDtxYxifcgj6KYur2SVlozu8w7llk3sTFz7+iekjIiL6/PQbbju+F9rwwbdaFwoOEfc/XrmO3/h7f4uIiN7+/h8lxYsF74KqcPB2yO1T+e49FvF+mUCvhpgfH3DzFAYqbMONnAqQEy+1/4OG7fRYKjyLMn2rHFvUf54wr+LLiV3/UjxskhZlJ40Mu4bLXKTkiEXMcBNi3DZzfUUygdqJl2HesXDal3UQbbJcFrl4eCiIRfH0sFMvNgqFQqFQKC6P66Sx2bkXm96WWRuwCBANI+y6k6m/S2/oWYI+7GcBb1SI0oddi1dvMEVIkreOzHjx8IjbCAFz8mIvxMNgLc65VsDcckK66Bz4QpCeseHhg6mZIOw7Xy54AkisiFbrMOeTtRPQPlg5ES/YKRTK9IxWYc6oVwmGBquZ1Z7bdoF88SyPYeYGK9MFsyPYLs8Dm9SeOEOTY3fs+JjPxSK1TxQl7eOLg3FjPrcax8q83oY/i6lxtlpybM69g0MiIvr6HVeP4XfmzgCYGyJlb64SiVC4lCgtAd9nhUKQMuw3W4wXEvRVNbyF/1F8kjr/fEL4tbAvP1OelmJIO5zNXTDZJn4WSRZXiGRtzCz7hHnlQfp5xOHkeA7iOVtNgxGzI/xBJE8NifqoDiEeDqk/8jlXQ85L4mEvXAZzA+bPxuZSsXVBUKx4dti5FxuFQqFQKBSXxAvG2Bhj7hDRf0dEbxDR20T0b1lrHxSOe5uIToioI6K1tfbbN/W9Uy821rpVdjcQzIXVt2ROivqW2gpFMDaxxiYkyhOMjWCKYikMCB+/YhE+5mTx58Ov03DvBdMZizjjnB8uVnuwzVvkxIPWZpLHqMrkgb5Q5Doci/IKxxz2jfByMDme3enzlRA0NHyoBwpRrqMI8n6PdVIT1iuBgUJBUNhZhJM7OgdT477PHoKxQVh5GBPGsmSW5+QCTJQrqomSEUQr32Zs3LFHfJ5eGzlDH5u58O+DPUcNPZgFauiNn2fdzQ8oc/O8UAzplqvkjBngbR//Acq/h/RYz9wkmSVlSPUGu/G+igYweUzJ8GIcK1mTmLqRD5dGMBzyGVQbZ3RQnJYiC7+u6U1iNOlvWcK8ov4xtRfY9ZQtKRWpNDJRXkHnEpiUstIn+z+BKDznNjF0CRW/iT18xrD0wr3YENFPENH/bq39KWPMT/D3Ws2af9lae3fbjlXKpFAoFAqF4nnju4noH/Dnf0BE//rT6ninGBsios4aXxgyZk2yqCjJ1HhFfNSZ1NtkSni3WUWMDWxLnY+3L1gT9xnJ6ODrLThhxRiwKgAbcuE1Nvkl63gsjVgFhgR9zDxFjA0WEGFMYhjRnC/Wjsk4ZtoFUVjHSyeQWXrGpjAdjIXJELBWHGBF673oPHHivOnUMSaIWpLn2kQCpnbutuMz18/k2PUxOnfbmKXqWTuzvMEJAI/d+P/o7A4REb19+CoREX1y9JZv85p1upuHPLcTTlIIfRHuB2iiiIiahdv3zT/zt4mI6M0v/ggpng0+93Nc0JK/F4Mm/aK+HOJkSqFBEjWdS9xGrOR9Qr6C3sR36xkUPsTbH3o2pMdWWZPIpnweBpYpjxCy4ryUWIuMRZJbHNcU5iyexaH0gLAX2UQ5CiPLSVTsRt1FhYnzOUvmKrtUKHkRR115FqnM3Pg+4nObzfX5wtALKe953Vr7PhGRtfZ9Y8zHKsdZIvpfjHNJ/H1r7Zc3dbxzLzYKhUKhUCiuHK8aY347+v5l+dJhjPnfiOjjhbb//iXs/EvW2q/zi8//aoz5/6y1/2iowU692CCPzZDGJmhd+JXe+1kb3ob31lBmwW1kiXroXGKNDTQuK2ZOOvEe3Eh/OwXmxHJJgyw6qeCuR3QSGAEwBChI2UTLhDWLeHyggGeNYB8amyjaoJdjQ1sSH4gWXF7hdJUKZU6XXBRz4cbWLCNGCCUNvJaH7YE94q66/bAMHO27Rjf2HQ1zY+L0K4hKOxs51moVny9uPmJNzfjE9TE6cW3jc923jm1Z3nAdLO67cf/RK67I5u8cfpqIQs4aIqLXWqepude5aKivzl8nIqJ3z1yb02PWHT0Kf0rjR8yg8TX8tv/47xAR0e//5A+T4ungcz9bKZMwtCytaF7iUgeSMc2alBbcPqdKJVKnpK+oaV6adHfSXLIUTbrfRBq3rDglmBsRhZU8SgVzEjrL+5cRpD63iz8g70pqXqTGpqg5lCyVeDZ7rU00XCNYqqwoZVFjQykqduN9eScCpet9lbTJ0yeL7m4S8lpr/9Xab8aYD4wxn2C25hNE9GGlj6/z9kNjzK8R0XcQ0eCLjWpsFAqFQqF4yWHs0/33FPDrRPS9/Pl7ieh/zMZszIEx5gifiejPEdE/3dTxTjE2ZB17AoZiFYUeyeKQ/uUfKxQRmUREaVFIYcd1ytE4Xeh7xR1KrQvs+8ik2L8K5kQwKDI3RjJ+/q1j23POQIx8Nk10ZyFiyoqO8FVqbWL4iClMx68sQv/Q8MyZNVrxRFCkcs1RSqNlpDNhxsZnUsaYoLGZ8HmahuXP3r5jWV7Zc4zJrYkrUnkOjc/E0Tzno9DGNriwbHflDJozx/rEM54yezO76drMH7jvDx+4/DxfOfyEO64JFTQ/NnGMzYO1O+Z3Hn2SiIi+ds8xNuauux6zj4Kl6QM3mHbJc2V26o99yTE3izvh3H71S6q/2QToaIgo18ENIPu7qmldkuy44tiK3WLGXlmgUSLWV4i0u3hemKEUxFLz0vEHPMeSqCjRVlJOksWgoLGpZSVPzqdkifryw8wU5uzzVEn9YwmSRZcMylAxTM9SgZ3K50wyMkvaLemLoKerZF9GIeTkWgrvgIKIiH6KiH7VGPPvEtEfEdFfJCIyxnySiH7JWvvnieh1Ivo1PpcjIvpvrLX/86aOd+vFRqFQKBQKxeVxNbrlKqy194joXyns/zoR/Xn+/CYR/YnL9q3vjwqFQqFQKF4a7BRjY4mo6413B0n3E1EIu/bJ3QTd2UfuJ9uUaV8p1ltH4uEl+2xWffnUwRWVCNR8+HUqHobAr0jFIqqxT91AF94VFbhRuKnicbrBcFdw/xRcUSg54N1VvM1LbQYX1MXKjWWOQpnsikLhSaJQ0qBdpcsEeKa8S3AS5nE0c41e4+R3H5ueEBHRGYeZoxjng73DaPzjpD/Q6Wbta0P4Y5tTN9nJietvcsxuvodu/3s3XJmESRv8ZzfGLhT8mG2/dd99X37okvrtf+D6OHg/cqnd5QSDC7ftpm5w49OW94fr9M/8Jy4kfHXkztNbP/RFuu6QyfaKJUf61C1TKnUQ/r5r8bjoKw/LrYUBl1xhXnKLMXkXRWo38XbAZV0RDw8lq8sKaOJ2jRNkSpeQEDDLpHVuDHXbiX0quaAqbaLPPtxeurtrAl7K3WBGhl8LN1DRNlxS3gWWX++sMLE8x0XhNF+zSmMTnbAsWeBV4AVjbJ4ldurFRqFQKBQKxSXx9AS/O4Ede7Ex1FfEw2AwpIiXRIhfkk4bK7qaiNezJlG4N7MjCxYPI+x6bcHC5AI8KVzOxMMJi5QOAYwNimEi7DsGSh6g5EC8uoj7hN14fCgMCebGjpnxGgXWAkJlJMoDM7RaunNgmIFo56H/ds4C2oWYlxBWmiiR4f7YJeR7deIYm09MHhIR0WLk5vxw5UKr39u/6dus9risAxfVtJyQz474XESMDT63czeI0YWzPTpzbc9OHCvz3uiGb3N/7JiZk7mzc3rPfZ/ddf3vfej6OHg/lGGYfHTm7M3dPrvHLBuXbuhHofxCx8VBmYijN36Bi2nyqbwOBTVRgiITWGaxw2FfCDMWncUFFLFKlsdItiKlFcrHUPqsSP/GBEtbE7QmlEfave3E82ooeScO7cAAg0HID83sYVtkL3hfJdy7WPpFnI+hgpyNCKyosSWlAqBBMM3brLRCbg8i3vBdsGLx55pYHNsoCCIvvZPbJkoZmxci3PsaYcdebBQKhUKhUFwayti8mLDWhR5DY5MwNrWYQbiCPWsSsSP43Jffpn3ZgWh1hmRxCLse86s89CdgbGxhdQZ9D8KuZZI8osIKAiyJTRP1xQCLhMKVvuSAYKLiU+Q1NdDYcPg1jd12PA5LlHGTlotAwsKeV4wNh3m3kcZmNMeW2R4+/wiBlmHgzo7beXPkwrw/PnpERERLvs7391zI9VuHd3yb94/cvtUhs3iHbmKjfY6xjpZydsxh8V5b4TYNZ/yzc76285CIcAE90QWzLieuD19085HrZHx/7ts0HzmmyS7dZJtDN8bR3pjPQXQNseLE/bnHRVynbvvGf/03krEREZm5O/+7VGTT62aSVazbgHHI/v5EWG3ym1g1lxbCuMeyKiSStUgYAv7S5av72G4abs7Hbpu0jiI2CTt8un8eRmlCwrbX9nSFMUlZkWdsEALNh7X59djEXhBRnriuwlIlrEUlcd4gS5U2zYpvDpZSEHPupb244wqVUrzetUKcaOOZwILGpi0NWPG0cSVRUcaYv2iM+YoxpjfGbCxBrlAoFAqF4vHxAiboe2a4KsbmnxLRXyCiv3/ZhjWNDZgTWTBR+u0T1qJFQiWszsq+5ZSxYY0Nb3t+PV/6sgbQz0TMkGBM/Oqc7SWRWm16LACWBHZwDoiI5lzyoIPGBsnxhB888S1D7zNOt4aZgtk4aFNmI6cVafyKi1cffF4QRJQyNqxfOefoNL7T1jOODGKWZ7mOSh7wicK1PGodc9PyhXg0eUBERK/tnfk27x+6RHlgbBY3OALphtPLtE3o33J0UsdjCIU5U8agT/QTzE6tmmTcYJ6g12nmy2CHmRq7dOfNdCnlECeJ9KUljtyc925xOQkuK4Fzfb4Iupzzc9fos7/80+4YZnDAcLz1g1cXWeVLHVQYg3juXt8hzr//HR/i5GhS6wK2wuvVYpuS8RF2kBQvtolIxRb3uhTDpGxr0q0vcSCYYD+efAWP7rtOaF6S6JlUkyeTE8pIsXhM4SAwNZw8Dg/EksamwlIlfQqmJhxTob2JqEH/vZhjekmLzf2zzIjvA8yNLK3gGeckQR8OFvak3SRBn7BZYaliAgpz768yKuoa4UpebKy1v0dENJhlU6FQKBQKxdPBNXqn2i2NDbnVNPQmq4i1wJt3L5Z9viCk0LkQ5RqbWorv2IeKiCAwNmvjxrBknYtnjAo3USmXDlFgDoiiopRY/Yk3fDA1y5ixWbFtZhXGrLGRZQ1i+jCMhbessRmxtuZgEhiIw/EimdsYFI3UqqyDgZaLUo4u3LE41+MZMx8LbhwVzoRWCEwcmJqjxjE3d0YuWgp5boiIpvuOFVkfcUFOFLi8407qeBwxQqxtWh65fesDjkjaY9Zl5k7YdBoinDDFFfQ5I2iGoI3BNtJ7jUF/8W8cubU+GLP9cP1XN1jTxEzNZ+44VupT+05fhFXmvcW+b/Pu2LFU963L59OxDqi9cPP6/N9wuXGocL09I4h5yFIaW6S396v0LswDui5E6hjRBPaSkiawBcYGW7nq7woLIL+iFgxKLJ8AKyI0KYG9yFfWvUX+q2GNTcoi8XTQRuZ2KURK+ggh9CGilIY0NhlbhXPehcnL4pBhjOkE+jZiND17IRm0AkuF/nyhYDFGtDT5Z88a+cLEZcYu6S9jUvCcrdwzkT0URO4bMPQFlkq2lfOJ7yt/z9nst5JdZxM7r+7t4kV3Hz1NPLMXm6Fy5dbarNjVQD9fIKIvEBGNXr254WiFQqFQKBTXGc/sxWaoXPkl+/kyEX2ZiGj2LZ+yfW88O7KIwok6SlmcrCBkKY8NmssICOnrj950g9YlzWMDBgVFK4u5DaTWxueXiWxhRcsr3Ea84YORitkq5LghZmwMEw5GMjaFMfmxsL3xxDU6nATBzM3xPLH9aOnyyUg2KY50apecD2eehj+N9pzB0bnbNvNwQU4XnGG4c/qYZbK8J5rxxA5GYWz7nK34wb5jNKC1ASsT5+7B58VNZm44XQ30LUdHbp6vHJyTBNiqU2Zd1vt8Hxy67fTmzB9rVkfuA9843R0XFTV/hSOs7kQM4CuOGfvmV10k1Z+4/R4REX1meo/tuv7/aBwiwZC36OGJuw49R+i1iEQ7Sxk7onBfQdOz3sd+Zm44Gi7OBJ3lgUHU3TK9z2JbDfYhIgjaAqz643sdS10eQzNKb9Ce81/bLl5hC71JJ7YFHZnURATmgO3GS1nOCt6BfZS5XSo6F9dIjEWwJjE8W0Q4x8x8+MzH8cFi/H26v5jTR7BFmd4Ec4+idLqKxmYoKiqwd9IuZQCDYZFrTNozYhvBs24iV5DM2ZW0MekzFNuEsanZFOe6iR5jvdQ0VbI8lzQ28pn53GDpWrmitFaUQqFQKBSKlwZXFe79bxhj3iWiP01E/5Mx5jevYhwKhUKhUFwL2Kf87wXGVUVF/RoR/drjtY3Ew7GrAgUrZWkD2b5QUkEKi2si4rj/RLhMQfgKV1RSNE0o+IKQk+3F7hIRgty2aVFPP46ozANKKRAEnCuEJHMfKEQZC+CEoBH2JlxKAYJhIqJbY+eagbvvYOxcKyi7UHThgcJdMofL/O3ojIXBF9iGRmdz51K5v3T93++cOBYuKFzvceTzwngti5/X7BFaHeSh9EhGuGRX1OqGG9PkppvrJ45c0c1PHjwiiRGf/7c5Ud/yphvL/DYLts9DOPa0SV1Ry9vut/NXWXj+SlT481UXuv7P3vyAiIj+hf133BjGTkR83jvf0TzKMjdpX3fd96k7b3zi5jU5Zuo9chXhvlod8nnhe269n7qgxnuhUcuuIV/WY8EC6i6n4OH2bBap0NTf23BJxQ9EuApGKONRzk3fr6MbS7pl+FaAq8BGXTRIiCmFoN5FwaLxgnsAolqZCqDklvGJNvFdpP2XLrDYpkyd0JVKHYjo0Szcuynsly4aMWegiwT/G9P+x4JpPD5EYIJ0DSUh7sIN5t2GIsCjMOzMNTSYoE/YNmLupuTmkxD24rFVy3kIu/G5hrDcXFGCPkPXSzysriiFQqFQKBQvDXYq3JvILYJ9McaIsYHIEkxGV0iU53aEjyBdZHIv3y1WJYW3ejA3GAsEvJ6xiYumZStG/ipCromCiJd49YrQapQbaApLFKyofcgtL7objtiOV+5+TFLsx2MaMUN00IZw79tgbHjAN1hMPJ26ZfoFC0/j0NFocK77lTMIMfEYBSjPQ5s5sx4fzR1T8wEXuwRDA1HxRRfYkVB0lFd9fEd3U7BxkTiSC06uDnm8R278d2441uSNIyfY/ea9u9k09lp3Ek+XzvY3ThwFsuDtKCp5kBQbJaLFTWZqWP/b3Qmq3k/fdOzQt+1/g4iIPj9xzM0dvngfMVsVs1RLvtfWC7ednTJT88jNa/aQV4zRanw9Sxksc8g/gBnkUPeD/cDUjZkNW/M9fcrdrTDX6G8rsITuu2cxsOKF3ZiywaXjMSDVQFjlcxqDNm+TreAhlI8Z2SxRZXpdwByMIqYIZFTXIewbD4mkaTFhW0geVwkHTpLVpWyR/xtuQXWmjHDSXLBGfn6FQo2yOKRkp5o2prjQmNumf1opIySFy1l4dP7Q9IwJnyh8D4Jp8TCK+stKaGSJATNz/nq0ovhmqcxDrZyOt1dgo+vFL1O7RES2BVt0hbSJMjYKhUKhUCgUu4edY2zIGh9yHetcGn6FRxh2Dw2KzZqHz9DY8HefrM4zN4U3fBzLHSH0HKvonkNi20gXUAu3lknTYttm7A6ecmmDWeu2IzA38WrAp3lPw3yZZPBamxhYYcuVFuY6bQOrcMhxxEiYd2fiGI69iTNwNuVVZ5Ro0GsrePmC5GHNuk/GFJdhsKy7uXvhNDbvLW4lY56zgfvLkKzuAgUlve6DN/ga65dQ8JNDntt9N8dX9hwj9U0zp2v53PRD32bCFw+MyV0uunn/ltsuj7no5jzYsU1670HTs7zJjMpRmPQn9lw1TWhqwNTMoPvim+S0C+Hkjxbuszl1ticoyPmA75kHbl4Jq8ClJqA9whLV8n22t+fs3olC3XEPnK9YV8SlO1ZY+UYMAVhBXE9/z6NkB540sQajSXVkuNe97gTMbBseUzUWwYe2lxjZCjsJuz7hZGR7veK/Z/G3X0zQJxkbEe7tj4s+e7YIyeN8eZeC3qTGUmErSzkQFco5sF3oPzDPiEEI4deSWk7tJf1K9iLT2ITPYDC8FkmUj6iGXkeN/E8Zg/P/s/dusZZ12XnQN9dae+9zqfov7b652wYixYpAlgmSZR5AQBQnblkhxgFLQUgE8mDlISI8IDmhJSKIIoEsISTDg1uJJZAcApJp2cIm2C1AgQeDTdQQO21HHcvB7b78/V+q6lz3Za3JwxzfmGOONfepqt9/1Tn7nDmk0jr7tuZlXWqNb37j+16g7NvZSqQ3Uba9BxUrjve+vsy4W4ZjIw10t0h0CR6+u8dxeA82LVq0aNGiRYsXj4gHtRR1kA82FPvalSVOAAxyokZ0N+zIrSXPDCjpE1ermohlH7ayUMyKkeoa/B4mf82cshOuzVKySRpRDrKToZY+ubX2bktUZP7dfltWUAUnWz+YdPxEUISlEBnIuXm8Sun5t6WiZlpl6Ikcl2khxpwbSamnkv/RZSoPwjpNxNl1glS+ef1G0ScKMr53farvXUmVEiuDlF+k/Ao72eTdyHckozoZUifeHhIS9cn+TH/BiqxrsUn49NF3AAD+0Unavn+a0JPdSb6UyHVicrY7lqEfpfaOlpn0RLHBhcwtj9Sl/PadMVVYfd2gV0/OE2K1eCrigB8IEvR+2sfiyTUHqL+ZloL4yH4nOb/CUTrObx6n33zqOI+diM0HXWrvmRwXWpDaU5BIzXBVojmU2ie/qUhYpXsL4dZwXtSskPYlhgOjnCm5aFQY0IgR6u7lkMwF+phZC2pieCb8bEt+kaBFPqO3gm16vcnriX3R6z0U+7ZjZNvKCaRg4g1cDM+pCbzs7D3H/Zxt98pzKdGytJ/Jdne2r5sQG+UV0UTUtVu06dt+AUsFPewerbqhOiqPuYRdivu5d718juFo8dV9HBvZ9rayLd4Bjs0DioN8sGnRokWLFi1avHg8pHLvw3uwCTFXIk3zygFW7ng9mdpBVaTGGwQ6o0D7hO+REvZFOT2UnTcZ5MyM0j/pW96P6tekLx8LUsMqJWZ6lhfgESUvA96JvYGdg35dIhxEGVgBMxkYiUjNSZfScnJuHovWDbP+cZVJNjsiNoLidBuRxx9KvnrJ00i/uV6n/Xwgxo8c8/Uuvf/uZebYbC/Te8vL9NtBKCJEDmzQTiBrnpSVbYyFOXhHbuysjjoSPgitL9SSAJnXo0vartrEtkcU6mxKsM63x7TdiJ3AP958AgDwu+fZUuH6SUJfHj1N+zl6UiI13VnaxqN8PIIQi3huk2+0FN2aT54kpOazx0/MPMhxlf4OvSBosUT7gIzY0NahG8mbICKIWTB7XQ4yx4s6YnPZrfQ3BBY9akF00kY3lNyz3G7akt/C9oF8bNZiLJuPnSNvWZrGruQE7dU5MW8Mcv1yS2BRK4UqHJu9NhJum3bo+DHkFZHXJO3sRoPYeL4Jo4Ze7OnL3BAyv9G7qiRvM1DTDZshJ+74Vy0cHDKXkRsU7aYX5X79PmbcHlSOrxuz5zGltuf/l7R4dXF4DzYtWrRo0aJFi5eLB/RMdXAPNpZlbzNfNYckx8bpyVSzKEVsJMNWJeCSrW9Vf/m3Z7czq2Wq1dmqKPJJuL3BlJIZRC/IzUqySXIx2O5RbxRi+0paAZhKpFi8BoB+HYttt5HKHVa+mPSJKBgrhMi5oTrxQvRsxqO8/92xIEDHgthsE3owLQQlcSrP6U1BeeTYXcpv+nAsrxPMcH5pDCdFyZjGjwsRW1meu0oL5Mocqh1fXwtaskn7eyrOkERPUsfThgrARFjUeJSIoEmxZ+ccEbtNiUgBwHvrxBf62iYjMgBwIe399uWnAQDfeC5v6akAACAASURBVJb5Rv3T1DZ1a5ZPUwNEasKVwCeLfHlzvsl9Gk9SJ98+Tb/5jFRnfXb1Qf6NHJwnuxI503Pc8KOIkC2upmLsWu1TQ2yGshqKCCArHNdyLlqtFe5GUUmn2WSvgU6meR/Hhjo9RwaxIUo07NEdqZkidkof43U259mldvPfzOAXNIaUeeJ1v61l9g69UFSsxut7DseGaJW9t2k1mm+7wikhd202VgJbDimybc34PloVdcP/vO4W57k9VY6Na5eIWpVjsyeUM1SU05Zbr4/kx1f77DbiIS1FNR2bFi1atGjRosW9iYNDbFq0aNGiRYsWLxkPCLE5uAebEKJC4pPBnWmlwC3LvbvRQeEVktnMnJI4ljOGBIClYNCeRKxQJeHJCnl4thQ17odRSSLkktOjvlyKYokykInEl86wT4l3IopHcTwA6IXMyy0h/Y28vthlwuZ1TJg+l6RIKiWRluW616s8kPFIjseJkJRHyuOXSyLG2zELInIfLPeV317J0hTNGAGgF0uGhdQgcwlqcTZf+6A8/uYNOScu0+v3r9JSyzc3abnn9xdv628e92mZ6pu7ZO/w7jp5EdCwk2XmJM0CwHBJmJ7vyLLciSz3nWdLiG9cpDa/uvwkAGsbkcb6O+eprPzZs7w8tjyTpcFzWU64SMchrOWcIHxuDED1eMhu4mnq3CdPzwEA3330PoAsFAjk5bdll5bJeDxINLfiioNaZJTXRb8sCbx2mYNLqI+Wqd9vLK+K315sU/s1wqWShuWaotGrNYzkEvBMNK4rhfmsGOUYuBTFBur/G9hlXV/urctUs6Wo/Bu1SlHyMJdJ6svKgCXSlueXLjUXAn31tv1SVFHq7tv29xEjnTATCXRL/SpWWBBo5XydLUWlz/UWasdeaRuoLYHtF+jTJdSKOWVVELFo172BPM+zpUZuXVl7arRSat7ilcXBPdi0aNGiRYsWLV4i4sPi2BzWg02I6LpYN4KkYN5UF8qbCTphfqC9UB/JjYMprV46a4PBEcRqT/NeRCwjN5VMy2U8RIZOpNybaMkbQ06XSb6MakZZjscjNwDQbSV7EmuDTsq/NyIlf77NiA3l/K/7VEut2WUoM99LY6nA0mqaL3ZjSRrm+8bPEnEpmY4gZJ70x1L0uM3UMJXyvxak5kL6dLYtxp6+nDabx1LSK8jH07MEY/x/jxMy8chAEfz7nW0SyvvGZUJYrs5X0k7a1/JZbmh55sjiIgFAlGr3KF927z5KCNDvSqnzmSBlFJr81nn6fLrIkzsIsEHidxhd6iik4ekkT+5W0CIagK4epXF99iSVd3/XMiE239Gfw4fah+wE3VuHoh9AJg0Pl+nYkZAfTksan7UPWSxJGk5w11uLK2lPSMtDabGQPiwRWD3+BKtMRjzyHuCQWWbNFL+06OdOLh6itHvJwxVElqCCksd9GbDZF0nDFi0CKqJ1NrwwH7dqqWD66X7qjRkVMZrMteSF8jxJtiY66pFwP+YwHzOPpxK0nRFl4VFZuW/X3q/9p93VkBP/2qFT+vZNY/co0WzMZfvFZxVCcYuPPg7rwaZFixYtWrRo8fLREJu7GQEpq1CBJ3OkmOUxu4wus9MMz66Ne5Ent8ZLefOVQWzIeSG/5Fx27Mv4iqTAiYjlMmwU26JPEkSGVpKastT61CA2x5LtUyyOvBUtX9cMzIydaALLZWU7rVPGemYQm6dCzHhT+CbkXuy8QGJvjge9KVcl2sKSay0HP8m/iSvJYpdOsA3ztXEGOVRqHyFihP2VIDbmeA8L8nGkRPw8vb4Szsu3LhM68miRS695nN8TscB3ztJ38EyEAZ+l9ldPc9+OPpBzgoafNFRcpu3OoBjXj9JcvrNK+yXiyHLyi6v0Oe0m0n5ly6GJdUJcCUIjiM3uNKM820dEi9KPPvEooSOfXqUy70/0afu4y2Sh98bUJ5a4b0S0rhfDT3KJAGA4F5TtcluMNchvlb9mhAxpofDmMrVJqw6anS4rPgkeMVHEkRwbAxCRC7TP/JLilxaxoYluRhfkg1i/n9i2Z99x3J4CsaFVipxfvG8NL8Kx0XuZNNuV76c3y992jle06Lg1RsL7+EQ1tFvfc6izI6t0FY4NUehr5djsH/MMMXFzAH8Pt/3eU3Zd49jsE+jzwoC2zapcx7525br26NHrioCHtRTVyr1btGjRokWLFvcmDgqxAdJTMLOa2hqmVic5E0zNcq15HbMN3U0oNxTRMogNM/iVpIxcp589iRfrxI5zIfyWcTu3XyDbn1VdRKJ62SFNGU+MOhozz7BI+yUqwq3aGFi1rFjyfIjYQPp0scn8DAq0kW+yFkiIlTtjpSJBD4P0YZTqGPmJVueMWWsP/YlwLo5SO2/OqmRSn4IxRZx8JRh5S4KWhG0+dv218KOuF7KV31ylMT8T4b9vLR/rb9QI8lpEAoVbQ7RnIXya5bPcp8VTsb+QtjshEm2PpeLm0lTuSGUWzTzPh1JiX4UmzfnEDH0Sy4DxiDBY6v90LAjLm/nyJq9oOk3jeesoze3Hh8SpeUOQmt6cuFs5zmfbtN+NVKOtpP8LUwE1XMlY18L36uvcmrjIvzkVxOZjUtL2Zp/6RB6Z568B+7k1NHq1lWA+oyb/htwOIjWPDGJzNbpsO5T3CLUpKTg2dR6OJ7lZpICICe8fu0heWVnBk/ZT7jdXVXK/lSodj1KR1yJfYnXn1tzbstBn/V5WoFREixyfaM4rMueI4+QtXKVQTaDPWyqo8OMMKZr9NHN52A54bA1yLX87cG1ul1Dh/TyPV1RUzuot+BZhk3iLbb/maIhNixYtWrRo0eLexEEhNiFE9N2ka++deSLm3yp34HRsFJmwiI23NvDVE6wgMD8iUkPkhn1hFjiTIzf7VxRBUZISybF9GseSM8SgQeOJqdxh5tktJduX6iKiJMzsCwPKUGY6mgWKseHVJvMznmwSWnFKxEZgmHOp4NnSMHCcIzdq/OjMF3V7lI/h8VEax8eOE9fik0dlhc7lLv3onZUxLVysZMvXwjcR3kmAOeBEzhTVkbHLmLcbGZdBq9advHct7Vyl18srVgZJdnZpMt8LMQclx0bmfRArhc7ovyi3SdbgRzlv1TCQ2azhpnDutsJT2j4WHsvAirM02evH+XjvTqU9QcXeWiV05LGgJLx+tqZs6cmYkDoe/+lykDGnzy1i0wm3hgiZStErupS2YZV/88aq5Na8KVV3uk85OSeTTpM3020cp0oQm3FpuEge8JGPWPHE6+bYekNILO2NAphVJNlrVtGbUH4ndyRtLALluXq8zmsy/PsMJ9mucmxuqPj0HBtWd16P+b+AfdVJVY4NuWsOGc/3ULnWzC3Bc2z6fdouN1RF7dXNsQAXUSpa4oQS4bdzrJpo+1CqyuHYZ4LpTTctOjOAY7491OQhcWwO6sGmRYsWLVq0aPGSEfGgqqLaUlSLFi1atGjR4t7EQSE2AWnJh5DuIswhRQ1dkpLfOsIdYKTI9zjUwpHPgLwURfIuoWQV8dsj+GT3T/Jw3pqSS1maYHn0RtZwRvglKVuCLsthJJ7qsoxsV7KsdZ2fY7ls4cW4OBfbXe7TMyGPng5paWItfTrbcClKhNs2uY+EyecETi6PyWuzxHKySoP/jlUik35q+az47dVxGtA3jrPT9dPj1KedWAaMIvw3Hc1P7biQsuvejZ2fC4Y9WtEyOQf0PS5xerHFrSkhpRDiTpZlVIixxnCUjbTDZRKSSjcrsYE4ystjLBdnCfdmzXJyih6WnwPA7lj2vyqXUju5QGibwVJ+AHhXRAlJnO6E6Mwy7/4qn4PdRiZExAJpbTAtZCu77Zf5N28JOfztIR1vLouxD7RMIYEamJOGO3Wnr5CHvcU05RtE+O+RkofN2iDP2z3yDXn50txH+Lcu4XgmatpY+wIeXx6HrVxTvS8zt20rebhcUlUCb0XoE25ZhkvnS0deBvaXIvt209/ldrZko+2aMXflMtjgBPuq5d7+/u1cxW9CITpPHlZqQWXN7nnl3ib0OvaM4z3tAs9ZanxNsa88/T5GQ2xatGjRokWLFvcmDguxCRGLftQnfytHPox8Kne/cWSzQqDPZT4zxEbCCgESsTnqSvKfivgNzFhtJ9xAfPm3LR0VIuskZdfXUh+9FvhlqjyLktw8SNvXvsRaEJt+lVEYZrZEL6LbrSVs0oSScv/Xu9SX83V6PUqDS8PBJB+TGfY+grYV9VtKJv1Y6rA/vjhLX5G07Vx8Gt48ziJyT06kXFbMNreyXZzIJHQWtRD0a8V5kbEK2Xq1mGexzLo00yKZl+XLHYnZuR0lLrPsXgjNI9ELYyMxiXHoyZGI1QmhliX8zP6u11ZsT8q5xcyz3xAFK9uxpfQcY9+XaRvRERKF3xdRPgD4lpiCnl2lHfUkTMv0k7ALQNEpnQOaneo5mNo/Os4nCS0U3hLS8GlHo9e0X4oVjrt8cq5krGybW5bWc66LvvAtuT4oaEky/CPjYErkhPPuhT6r1+yuPCe8ASVukI0gcZkk9SpqoqgF0YqyXUUgaxm5I8zy3Fa5CjOQmQEnBQcnd3/E/L66j0hrEYrBGQiT1Bv2oCVl22UfZoRm+xuHFuVikwp5eA9xedaeFTflOaFjLgfAV7bce5JO3W659+01/brjoB5sWrRo0aJFixYvH60q6o5GQFofrWUb/HvfGmYuuc7vZY6NfMaMyD2J27Jyclu8xQHRozC47A1zNIRP/9lawWQDUs5Ko8crQUcuReSNcvOT2am3HIjSh8y5KJGb1CdyIFD2l0vXJgth5sxya/ZJS8JF7p9S+wAwiCkljRoZvSI48/VpZqvMZh93KaMnOvZ0kVAFlioDwNekfHknCM32lCXQwjsxGfxOBPL4nVFEAqOUnJ+IMCDF64CcdW1lDs5W6UfjKrVHVIQl1gDQn6R54XmkfRM0ieKEABAFcXrrJLX5mdOnAIBTIZGwJPlymxGbb5+m40DLhO1VOZdaWm9MSQn0EYmjuOJTQWqYzZJXAwDfvEp/X4t44OKaCBFlCyrX2iDzviI6lt4ej9L+T1cWsUlIDY/zAhQnJL9MzDwrHDRaKfTr0r7CZtbk+fDc7kQc8JTcGkFsrNjleXdUzIfnnvFa5RykPonUhIzdm27SmmV5g9An0QxFigrUgqhReb9Q0dEF71vz/7lUZNShF2z/wrjQPg/tLsRNPco9U7hLG4tAeaSd8+GtUmommJ7Lsw9dtzvo93BsLPfFo0XzEnf3f8JNbc/aM+eIbG/LUuGhxUE92LRo0aJFixYtXjIiHpTy8GE92ISEyPCJf2UQG1ZI9X6dmOHWTgEYI8hy3VxF8gS1sOgI14VpbcA+sNKCcv/RzOw+HktQW4P8nv4tvILrnQjESepLToQVUmP1iC6NS4aodgaLGmIjn6mIn7w/1LkYQGb206CR1VBhk3Zmdc5IWyByw2B7auFgRP2IimgFgWQ3nGtyId4wVgtLqaTaHKVMe3silUJSOTQMuf2dVEwRsaEB53Aq0v6Cmnzq6Ex/w0yac/zBqVSGnYhNAiuUHltUQawfpOntG4LYOCNKAFg9SqjBp09TBdg/cfx+MdbjPvFc3hcTTgD44CSp7dFKgZwhW5Xmg+f0To7Zk3WCjb6xeRNARgK/vcmIzfvX0s5V+g3BFr3s7KEVC4VpSSRL0CmtUksdODVErMcyRh5fxjUtOwQZxMZcf1LApOaXmxKxsX3K4oByTjthvlr7RGTVTFcqsrQKjohjIdBHBMvZSFBkUa4lywkkt+ZEkDle38rBMIipR0wy0kzUpFIhRLTI8XsyUlQKjBZt6wDK9suqqLI6yf/U81wAg9QQpaKBMFzcxOVxVgo1AETpe4pSEQ0ruT2pn/WqKG+hUevfPhsJb+Vgo/Zei48+DuvBpkWLFi1atGjx0vGQVsEO6sEmIKILUZ/0j7qcaTEb0DVMx3b3a6fpb5f5KNeG2Zpkm+ZxfpQ0kOu17AMl0ol0TFYC35lRRreYXaBI7INwba7FruBCKpLOhNRh+8SqJd0LEzgnZ09NESDzMHJlkPxGeAjUxAFsxpNaIOeG+iJeWwQABuHWDFdlhsI+9GvJuNY5y11vaa5Z8olYLcPqKKtfxH5eS9XPJOgFkYJYVEWV2/E07efxSdr/p44TUvNdxx/obxZOH//bpwnF+OajtJPtI+GDGM2YIMeD3dw8Iqoj7T/O+3z7NKEGnzlO3JrPLj8oxsx4a/V27tMynetrnlceCWTV3dboCgmas5UKtieiTfPNVUKErhZprj/YZGTobE3vBjnOLnu1mjHTSrhfxyWfiDwmiJXC6SKfJOS2MINWLR2ik8IrClYfiYiNnF9qchrnmXc23kzbpWjoPF6kOScqZu8jOp5Y3gPIZ+m3c44N0ZsZ4EGrA0FQjwxiQ+0cIif+PCuqfJyFglZmeZ6L/bkiNmnLyklFmGVb1QLzPBNn5QAYfs++ascKaqEmmM7iwEeNz+I5PXrP1IGa37i2OXZyEUuOTTnmmf6Zjm9eFeU5NkSpBocUAXdDx+YhxUE92LRo0aJFixYtPkQ0xObuRoeYK5PCvCqKa8nMkpg16ZK1ffJ2OjY5I5KvUv13nPNZtD+hfDrvhWOzs1VRRExoRtnXkRvAZAhUAJa2Wb3AtXhrjkkjOxpnaoXCHuQGMLwbpzOCpXCIFoYPIJoqqkPhqximMqsFKlUrkmEpYnNVbgHgShCCJ9uU5rNih8afzOi3BqLQ6q3gjjfn3FzMWiUmlUxxJRm8VEN9fJVMN79z8UR/w3ONnKbfP34LAPDt04R07E7LaiwA6IQfxTFvxYCSSFE4NoiNVGB9cpnQok8OiWtD3gcROnIjgLn2hnIglDMm4zPn4HApSJYYWT67TPt9b5U6x3P8YpuVh9dETJy5qZ5Py9wAx8oKMFWCpuLwquS3APm4buVgXQh/7Nku9e1SzofOoHqZWyODnphal7pMgEUq029WC3K10vEmKrawJ64E9WxA002HlgSrNM3qMBp/8rqTXRDFPRoMwixtHoVSY0XPZ2v8qSq/rIZy2xsUrbs9+jULh8ICc/RiH7cnvXdD22a/NQNholS+WmmGmtg+7ENqKhGc+aVWzBLJKTg27rd7qqKKPj2HW8N27dhHrbC9naeLgIe1FNWUh1u0aNGiRYsW9yYODrFp0aJFixYtWrxExNjKve9yhBANKSxDyApzkqBHuJ6YVAW59OThDPWWxEGFpc3fSgZDCbkSdt6aMmNK3ashIEXjdKlo/5IUl5doY3Ah2L7tE4XzdrL+lQ3jyn1akml0Im7c6pLBIsPmjxal1L2KEcpx0JUvS8zmnDoxt14MG4drMf1b57FvhNjK0uZ3t4+Kdj/YnRbjLcfsCH0VxFqXqWSZkmKKlNh/cygl/oFMLGWZ/RtCPCWBdyNLeDTfBICRCv0Tl/nkuMt3KRQH5GU+Emm5LNI5fXx7vHc7lteTiJ3el+4rwdVCz+xff5F+u75KB/zJkVELBLDe5VvC6MqXlXCuS3rmc/GP2J2mLylBm2PuSzl9ABjl4uQS45kwjZ/Jj69FANKS0rnM1nH5h7vjMu+Q+6TXm1teJXHXl5kDuTiAS8BcEtZ2t+U9w/+dvpQ2PN8WQnA/NQM56ctlsNkSRaXkeU4elmvLnBsMnve8H3E5huXlXApb1OqZY3kt5dJuS6Atz7GZKGE3J+r6JSieC1zCoXloVaCPy2B+qahynXOplstvmaw8F86bhSvz9pxu27a/x3RuzPZc70AJi0Yefh3RlqJatGjRokWLex4hfrT//sD9CeHHQgi/GUKYQgjff8P3PhdC+O0QwldDCH/5RfZ9cIgNYMt+8+P0kQpOCSlOSLyTIxHPHdcscsPX8sGuJPACwFoYiSQ8jp5MTFsDk0RlVCQUW/alQFJc96aJZpgi1Ed0xiI2IghH8bVuRwl2aaZyEmb0QrZCsFyKEeSjZS43plkhM53zRUIvFvLda53jCizGuZWyXCI4Si42NgxrQWyerhN51IrFATmTf7rJKMN2k37TCfLTOxG5UgZetiXXU4PloPa8YkbNrcoKUMDQzSOQkQ2eajrX2hGD5sn5Q3Iyz6utXJpnU5qLZ9vsaMkxD0RqLuXYXUhGuuH5bErdj2h7Ia9lriknsF6IfYE5AWcWHQvuS7YnhsTds8yepOFQ/Hbo5ychx/xMCNIki38gx5fjXNhybwrxET1QsniJigJWwiB9lygk0Qpmz6O5AHl9k0wdWOruEBvtB2BKzUlgZp/kuAylOB6Q0SJf5k0pBUseniE1fg4cWpK+lDZKXKaoqSI1RAbzcYkzpEZ26yxnbF98pY1HL5YVVJ3o1D7kpIb8ehLv7J5WK/d29jfa9xchTDvzzbJ/tQm37d498jCAu1gV9RsA/gyAn973hRBCD+C/AvAnAHwNwK+FEH4hxvgPbtrxQT7YtGjRokWLFi0ON2KMXwGAUAEbTPwAgK/GGH9Hvvu3AfwIgPv1YGOzq948gjLzWUpm4pET3dYW3/RpvURuuL5uy72vJA2k7Pt2KqdQq8pNhqp8Fi0/LTk2dhdzM0rpijMGJOcm/U2zQLER0OxStl5My4TnnVDwjrwaAHhDyBvMQM6WKcM+WqYGriiOZ0wXVbzNixHuKa0HoNYMLPP9YFPyP54KavHkKr8/SvnyEQ0axcKBVg52zMzg1XqARqNS1nwpX7iYMoeHQR4IxRBpJjkT8oIRsvPlsrLdbvNJSNuAp4JGvbdLvCLyT97ZpLJyCuoBwHQhiI2UcC8EqVmcl6iY5V7Q+oHIFtHIyZkXFoJqRKWEE0QjS/J1tseGzyLnj/KJ1Fy1TBM3pk/ngtTwOian6myT3h/lfFiZhHtGCemI1NB003CdeBgXJTeM6C55TDy2QEZsKBbJcyUjNTIei9g4yxSe+7RyWChikzk2vg/aZ/JMDNoWHGKjdg4OcozmP4gox46IDdESbbdyMyBal003pd1aubfj2CjQx1Jr5fYY9NMJ5c0MISe3xbzkPA/Qvb6BY+PlKaq2Bjr2sl2O2c7tPvPLzLGR/4vMzY1ze5uWCgda7v1ZAL9nXn8NwD//vB8d3INNixYtWrRo0eLW4+MhhF83r78QY/yC/UII4UsAPl357edjjD//Am3U4JznPqId1INNRMAUg2azNusgYsP1zV4yRQrl6bq3QRD071LjzWQJgtgYtb2rsZR9V46Er0ywD/iK2HiOTdqMlhdAmfy+zAIYXIO3mS+zS5oFUkqfXAvNNmt8E7dczIzC8gFYLcS16qeLhB6cCmLzvlSdWME2NdcUdIpcBaWrxzmiAukT5/vciMUBOZO/XOcMO1yX3JFBCpoWlyTS5N+zOo28Hto5XAhC9P7mVMbzSH9zHVOWzYqsc7G22GxcRm8cEFilxKyWYoFsd3Od5+npdRrTt9YJmWHVChG6b1wnk8onV5lj012m3y8u0mtya4YLqVYTPpM91/M5UPIomIhS2NIaNa5XItp3nMZM082dGI0O2YtUdzQqUoMiRmnXVrR9sE2cGnJtnuzS63PhjPF8Ls4RPV8FJRGkZpK+8rwDDEInVgq0PVGjS1eVBWRElufgPvPLYNCSCCI2RGrkA0Fsjoa5aS85NhsZO6/rHbk9uwrHRo4h9oji2VIQzv8+80v9iZncLA4oY91jugmY8XuODbk9XYmSAaWgaupuLNt13B77t+e6zOwrCsSmrEqa2VVUYiYE6FFXy8vZgxaRP6NVu6Zd/j/R3RbRJWL/efPh490Y417SLwDEGH/wD9jG1wB8t3n9XQC+/rwftaqoFi1atGjR4r5H/Ij/vZ74NQDfE0L4QyGEJYA/C+AXnvejw0JsYqpQ4nrlZAgzzDwGryfDih3aGBR8Fsm0HHnJIza2KkoRG8fD2Dh9m6LfrjpGUZmOlSO2T7IdSsTGP+mPUzf/m5obrAxSc8o5OuIrp9xSf6HBwCyPWitEcFg5FUT7ZlzlzHe3ooaLZN2SARM1yXNvGuVxlfGQz8IM6EIy+Y1BbDoiNoIeLK4ka7qkfo7lOgnSIKhOfymIzVXa77evE1LzjeVb+huOmRVa718lVGG8Em6PWEIsLg3fi20TNAplu0RcAOBMrA2+dfpYxupMN69Sny4vM3pFG4pBxtqL0Wh3vS3GrFovRV/kDTllloImPBZOlbU84Ll8dpT6eHGc5onH1PJZmA0qUkNEiEayatmQr5snHrERThUrtVAzGySoQCRwIrdG+mQRG6mG6sgfc/wGomLXhhzGqsMtKwwdT83OKcNfx3krqMVQ8lsAg5RoNZwgNyOvl7z/jJjwNfmDJfJb6lQJv0c5Ng6tqtxDJ3fMfDVUsLwid89UfpFDbCzPhONnJamvEAqOt1a2zc6xQXJf5LVBtr35JdvlHN9UCZY5c7yGpB3LFdyD2PgxW8RmC2m76dhohBB+FMBPAfgEgF8MIXw5xvhDIYTPAPgbMcYfjjHuQgh/EcD/DKAH8DMxxt983r4P6sGmRYsWLVq0aPHycdfIwzHGLwL4YuX9rwP4YfP6lwD80svs++AebMYYtHLBa8gAmWnvq6K0ImmwlQMlajADW7ikbSpHWCFF5Ca/76tlKsiNr5qguKnl2HC3qoFRmnvWdBAmrWIoMy0WYfTruSZDv2S1R/nbUTO53CdWrRyFtMNHQmhhdj8sidjk/dNocpTKmTCSC1Fm1raSCkN57DhW9kX7ZipGeh2roBdr4Vhd7mZjHiTLH664Tfu5FhXe968TgvDOKuvnkJvwzav03lPhugRBXZTTc2HaoZ6MU6Tdiuovq5kAYH2Z2v5Aqp5shgtkTZ/xOl+qi3WdQ6XolGaz9lyXLXcjvKgT4Ul9xyoRdk6HTBbi/FOd+GKVthOr4My11LPQz/EkiDzsRHH6fJ2Rpw9kfzQ1JX+JCGmoXEO8ZqalU0Ve8rwy78m51QlqwUyd9w1Wv1ErCMi8GnYBiAAAIABJREFUrlEQxp4mmE7puJhbaugM7pyW85k8E8tv4TVFns+OyE2FY8M5VMRE+T1y/3KVn+nvsjKI/J6l7GxTufUrv8dxa7qK2easupTdJTrmqrEAg2A4AtbkuD0FF9AhKDOtrOC2AHqnoXMTx0YliHw1lEOrqp6bipTFol0qHdux69ju2tPFPY2De7Bp0aJFixYtWrxkPCCvqEYebtGiRYsWLVrcmzgoxCbGgO3YZ9jWQJokwXnxJ28dMPV22ad8rvPl33kpKn+PAnlcDmO7tDxQU0ZLgNtTmqhLYeYokPAYRFSMcuxc9qHgVV8TevLCWoSSt1PxGgD6FY0U5btStrzdluMD5kt+ujRFewEhSW7NshKXmnYkD4/l+1ItD8vBjm7MhPAHtzRVhIOSvew8S58BY8ApWxWrk9dn16lT7x2d6m9YIvyBGHNeC9GYxGOWPJPICwCLc1kGk+ULHm8ugZH8C+TS73Npm7L/PI+1lN8uTeyxytDlES51muUaJXOLyN5wlPr4seO0lvbx1TkA4DHr5k08lqW5b8uS49SnA10sRzhTQiWuyxxv1xRfzCfJmbGJALIBp73efDuUSlDLhlASmSdDHuZyTO8kE9ayVnQZ0pxTKBDI5egUu9wrblkQdaUPtEyhpcaiJJP2lVKSjVz8a71/SLvVkmc2WC41apm5XeFR8nC5FMVlGba7LcjD9SXtbEVi+h/Ke6UX+vRLYKnt9DfvJ5Nb9tbxFUUObpnVLQnlwgxzb9Ny67LsmmRxex+JjofgxzoTba2FIw/zvlWQh8PcqPR1x0NaBTuoB5sWLVq0aNGixUvG6y3RvvU4qAebiIDtrtfsxtoZUABp8o/0XVmGWrUvcK+zYaa8NicEhbRIHlbERqTxRy3XrGTYLuvLyM28D0FQEIp7nSpik15b8To1ZJRHck/gVHG2be6AkmxV1E/GJ2Wul8aygWKEj7t5Nm/b3wx5opg5E6EhEVSRnKMSQQBy2fjJSkjKUoLM0vOLIWXTnbWrqJS6pg4QbstvBV82qwTXtBOKsl2ZsevxFjhqlEx+oYiEoALXuaFuXSI2vYjH9SrLb85RyY5Z5sv2PHHahsoG8JwmUiACerSx2BnLg1HnWwQYj9MckzT8qcUzAMCJURokIsoMNDjko5QPkO1YXm9EbFiWb0v1L2VOaUa5cweR2f9UIIHclmiFktHt9e36S0kGSjWw9JZ2FqlP0sDWiQNyV4pMGKRjUZKHqfdHI96MOOYJI2pB9Ij3E09aBkzJsxdY82iJvZs74jLRAxoI87yyZrp679qH1Nj2ea/01x8JtFXEJr23jeV/O774wZLu/dijvx6cfAGQ59uXuBM1sddUdGjR3BBZXlcMfnU3DrGh+ebCjt3/Z9PilcZBPdi0aNGiRYsWLV4uAlCoZd/3OKgHmxhTRqtGdXHe/XHPE70vtQaMtQFfM+NyiE3hf+YyHWaBa5Yz04jSliy6dfoZN6IirEUU5EgcLU+ldptZyJkRUlPhMY7RTwGzDiOwRaPE3nNsNhRSy4jQuSA2J71YGkjGywzbVZ8WY8rGn2VmnTk2pkx6JWJxIvz39jLxP7hWTv7D0rgiXoshJ7PmUXgllNq3c+0FzfbFznA8tERY3otTmeF5QS/79yzTvuG+wnOMSA3Fxcg3opEjYMqtiV4cE3VLkx2HOWIjbgUYT9J+Hh2lOf6Y+DJ8bEgcG0r9A8D7YsjZOadD5TMZdX4tt9+W/A9SdshnGo0B6EY4NeOilPlX/thAq44KErgqr1XyW4rEmFXEcszIj2NZOTkXluujhrdeHFDvI0RlDL9vSSNOee1kC4aKKBtRC5rpkqNHKYOFFejjeeRNL9XCwd23kBFfcvIWTkaAtjSWY+MRm4xolkgwYKbHcRjJ7dF2Dc+kd/OgKPtUomOWz5KvMzcHwY3ZoHMLZ35JpKwmjqe78zIFRI1u0NOboVTOSsFyqvzYW7zaOKgHmxYtWrRo0aLFh4gH9Gx1cA820xQ0u7k0alzMKjdarXQzcgOYCilmdo6HQ/Qk1ETxUCI3FLeKO5o+5u8q/8Ax/ENN9YmZOxGbnhyblGFTHvzUIDYrSe9owxDd+rfu2mR8zAJz5ZRkQMK5KaTvJd3nejkRG88rql04mQ9SVoBpdmtE/Y4FsXlrlUqNvkPQBKIX5CH8/upN/c2lcEbGo7I6ZjzqZ3OgsvvkpKj0fYmSLU1JCrM/ogjkmeRsjeM06/b8m8JttPNwhqxpx5Lhkku1KLlU18s05ier3CeiXbtjQdlOykybGfz2JPdJC39WgtgsUztvDCK2aF08XRDBmoT3sRB0j/wiABiuZQ43ZSVYv3Ycrp21AhFOBSsZKUkvKBWrigqOjYw9WyeU53wRsv9Jzk9yp4hAErW4MPcRNbx16KpqyLlKyvReKPrJc9wb2Fr7Ato4kL+m1Vi78lgClWohohV6fsnH5m4e1EZijh7YvmwsT1HmqdeqKIeSVJYyFMHi0Fw1lkVsOneDoHHw5Ko5i7GzupGVYGx3D7cHmJtQeoE+y8P0bWcEltvno63oyvO3ZoJ5F+IhLUU1HZsWLVq0aNGixb2Jg0JsYgyYpk5REqu1oojN6JjvezRkAMwtFRSp4esykwfm+jF8+mf2qUaUpqoha6ug2GYOxhy5Yba37EsDPdoZWOl7Vk6Rk+ArwGYy5EDWf3Hr6dRLudrmuX0mHIRjEb1hxnshWSa1e7qXqQRzuheAMWSUsb0tiM2RkJQul6m9R6s89neJ2EgmvxUUY3Es1SWGU0Uuito9kEbk7AUeL3P1F7OvK0FOeociaDWOMYScVoOMXfg5giZRS8aiVOEoTdQbwnlhlRIRG57r7x2f6G8uTlLHdzJGjjlIxq28k0yTUruLblFqBM2yWYMqkMNGpCNuqP8jc3NtsmQacRKxkXNOTSTH+fVIJFQ5EXKyHAvnpmrV4aw4WElTBT+pnbQrK9vIqeHcXleq4Iim+spF5bMEgzxpRRbRCyK9ZX+s7tblJEiNoJ+s9AQRG3vd+IosZySr7ZqqRN4/PM+EKBW5PRtj8Ov5Pd2+aixAtZI8x0aRImrIGG4PtcVYEUaOnkeKiqoox7HJH7ituY8sHL+nv5EoQyiubC8jNZWxe7TIjznQ/sZwGmmhUfVmeA3xwMq9G2LTokWLFi1atLg3cVCITYsWLVq0aNHiZSPW0ad7Ggf3YDNNIZdYT/PukwyXCa31sk2gBjOXr+HgXPu3L+FkObC3NQAMHK9CeaVQW7EaIL9nWTmX2FayHHPS0WE7L8ecCOE0LMtlkuw2zFJYWw+KMrTUUkpjd3luKQa46tJyCEm855sEo49iSbAsyn9lzJ44zfmpHA9K35OkzLGeytrHm+JfwHJwAOiFVEuhPy4z7YRQaxUBtrIURdKttxd4c5WWoLgcBORSc0L2q6M0yMtjvxyU57Y/LZei2JfRtQsAi6OSMP3xZSq7PpEJ5NLQt4+zzcP50Ynsj0trLKVOn2eiqyE0kwjvMFouj5DMeh0zkZbCdRdynIMQgYVv7JaiRMRvx1J0R4KtHO+FLPHSLoT2Fby+l0Imv1iaUneK4clSSO9EKYvrTv7ebbkUVZZ7cxl5Y9jc++wc9Fqi4KRZdp08GV2XouT4U8SzsLJPcSF90SUhEndrRHxvoaDCgCi2QC4+8PcpHm+Km1rycPQCfb4P1tGcpGHnLN6ppcJ+Gwme00pKl3vCUJNOiOXxja7Mm9tglqI8edeXWls5By4NDSSya5n3ftKwF+YLTpTwhZbAbiEekqVCW4pq0aJFixYtWtybODjEBqjLgTP72krms9vtyT7sUysRGye+RYEtypIvjNoes0qiCoMSxcon/KJkcY+1Qb9l9jf/LrMYjrV3xDSiGbZPvRM0o4CbZnamVJV/z6wImLAYQjPL6xWpkSzzeiMZqDeVBNCJYNuwLgmIo0OpiBABc7l0ZjwcM8UJOd40ZjFmJJlXCLrSRcTOiNSprYB89zjt/1TsBd4+SoKAnxDUxLZNdPDxcYIrLoiaSFk5USAA2K2ZTqaNIkVOJC/tryQNf3KZrA1ImGa596NFPt7dimNeFGP3AGaNpxilaV4nJK8+HU+L1wDw/iZ1+HItFgRi4knycG8QmyCkYdpIBNcZzayHPPYjRxYnIZ4owrFAgBeL3E5GU8txZRPU/B4lDGiPsBbyMMnQvHZHk8EHNc+VrZMnoFSAOa3mSC/vJ3qfSl+2c8v7Fa8pRUhr9ym+peKAzsJBkaL8XW/8qeXdNN105psAEIlCaU21dja9tGi3Gm+W90wvSljYSMgOiBrxOMc9pptAhbg8Kz6Qt814Sbz3ooScg8JSwbet5pfy+gXIw2x7Th7ej1bdSjygpaiG2LRo0aJFixYt7k3cCmITQvhJAP8qgA2AfwTg340xPnnR30d98jfZuFobiEAfyzZnT+R2R+yQvPRZAGXJ+wpiI2vIXM+tifgx5mJ45ZYy9Om7aTs5U8TRPenbdWP2j+gFk7DRyc93Rs6esvua9bls0wbnW7lNkl1SzCxIZtxnUAGkANFsk/PfLyjkForfpv2KCaVksZSd9zYZNhNiWevYU6it5DvYAakFgWyjIFtEBt5cJJ7Lx4bMsSFKxGz7DeH3fOuI3J7UkEVstuuSJ+HLy6PhjLDE/K1FQos+1osooSA2T/vEczkeMkrVeZFAd8wUvajYekxicnouAozvbxJSQzkBonEA8P46fba+llJ3Ob4U5rOmqt22LBuflfVTOsFYQzySueTYHwlicyUHaEUZA8Of8Ip5PtM2QKbahATahWzL8/d4mKv66XVMA05eSyrCR2JF/o0X+tT3BfW8HkujS2BupksT1Jr0Q55D6QOIwPpr2F4XjgNIRFsudIqbbk25t0pVcE59mbk1/nT3DY+cVNEKtXEQlJM2JTdxexxa5HlGHHNXEegjv4fCgLyH7kzZvbdImZV7MypoVR4z2y05PZZjM7kS99cesc7buq9xW4jNrwD43hjj9wH4hwD+yi31o0WLFi1atGhxj+JWEJsY4y+bl78K4N940d/aJ3O7XrlzTHuunfZu3bYUf3KQDZ/KHdNeM0cAx4LYUKzuqE+Zz0ARv5r+Enk3U5npEkGxZoLkBWyFI8Rsj2viNYEnru1SSHA3MNsss0xWqpSf8bV0dZhnQMxildskc62meTsiQrlPaop4TcRG1qEFRSJPw/JyNpvUiTNBDc4E6uA6PeXnC1ExX/WmnXZbzCsp1M5AThJWmlEEEQCOQjrOb/YJzSFyMixKbo8teOF8K2Kj2b500qAWPLd4XhGpWXKdvkK20ETa8bm88GNxPDjPwv85v05z+cEmIULetgIAzjdiPSDCfIOcm1lwstI3zWbd+SVba2D6xiLN89uC2JA39nSQPvUeOsih1xL7QqSzQGxkK/3m+UoEkOezRVu1ApLX0IJb2YdWmhnej08PnZUDUeQLg4YR8SV/bZ95r90/r1kWJ+2r4qzF2lVkqZGwRWz8tUQbA2kwWNsQh4yyj7wP1gwniZjkyqx6JdhN1TvZoqV83RkRVc/vUYRLbXAMcr0HLZoZjlbvI45XdIOVwljh97z2eEAcm7tAHv7zAP67fR+GEH4cwI8DwPDxN/d9rUWLFi1atGixLx7Oc82re7AJIXwJwKcrH30+xvjz8p3PA9gB+Nl9+4kxfgHAFwDg+A9/JnbdVDeldBUI/km8qzLuy622yc/Vjj7/SBEbt9Xskuvc9sGcfyvjXrZi8GYzX3JOIrU3djTLE+PJQYwnTZqYM0/pv67Jp61Wcixzp5j9kW+iBRuSqS4NSkUjTiIbXD8n90bn0cwteUM0RVR9E0Fw1BzRZNhbQWyebRJS88Eu2wgAwLNdep+GgcAcNfKVFYUmhucv7LnQ7do4dVK4bs51dJ4bk0eBUMngX2DBl8eQ2eWGFVXUmTHVKzRKZHEYK9CU+zKWiB0AEIRiZdP1Op1HTwWxYbZp0bBLqSKKfm5Lmkv6W2X+BRWh4aiakqYfHxmxo7eWCQUjGka0ihWHXaUhcuaUn7bxY86TTaSG1xQtAyanEWU5dDSU7RSRq1eedfOk3CCzKNojp8eabWYNnZITqLuy55OiXyVyMuO5VMCAnUNJeGqTx7QzxzvMNL8cn8WcT2pH46rUgl4vRHnzQLxm0nafdk/NbFONZOvcHlsF1inHpdyPIkbWhdahRDPD0bzT3BcOaU9VVA2t0j5ML3AzaPEHjlf2YBNj/MGbPg8h/DkAfwrAH4/xAWFkLVq0aNGixWuOh+TufVtVUZ8D8BMA/uUY4+XL/Lbvp5mxmw1mYzN9ArcWD2SkRNeu99BkrPKwV8UlYrPPiBLIGY8+9VNNs4J0aP8EsbnWbC+tz58L78SqmDIr0/N2T5WXzeA1A6WBn1QIBTEePF5kxIbKsGqcKJ0kD0dFWO11wwyIisNEpwTByRU2hmNzLRU7wu1gxQ6P6ZPtSfE5kFWPV65ih5m8PUVYGaTcE8lQiVKsVX03z20vwi+5kkNMVrWKBfN23HFVs1Nq9hjtHmbzrLo6GxOCwnNbUaptzvajjLm/TvuhAvBwTe4Wq6ZMpZZ8l/O9XYsKr6gKkztk+Qes1CEUoIrGPJ8KXSSZl2WJ1ExaCZb6dGoQm8ciYfxYEBsGETPNbs188djNKrT0WjaZOzlHBH48UuORIeTruBfEZnQcm4yOzOERf9wnV4VlzTYZO5/Bu2s3jYlbIjYo+1Kpjptp6AhCQ/6Hr8ZKDchuHFCmvBYDPWa0iFs59xzPZzSd4nU14/c4tKQA5N3YZogNOTbmHr1w/y94bk8x589Bi7QK6wauXq+VYPuRGm37NnVsHlDcFsfmvwSwAvArIZ04vxpj/Au31JcWLVq0aNHifkdDbF5txBj/8G2026JFixYtWjy4iJjzhu5x3IWqqBeOECKGbsqS2Qb6mwlCaSksRafmJaq5VJRLBbH4bW3fhHIpaMYlKS7TdLQ1MDObSxNZVi7wZ6wslzirAYrgnQth9mxISxPWAJQl4dn4s5wDDymn/pXlmtxSQO3USPg/lrJcLsOdDWl9YSHLVWuB/+3y20w0TqX20+vOL88AgJBiaVb4dHtU7IOlyRfrvCzj7Ry4NOGtHNIYhXRL4TYhMF9Le8+kHPepIS1v+17ek7blOIzS14USeHOXdBlMjuUo7el3NvlAXIotxdNt2v/7Qza7BID3t+n1uRlzuJZxSMm8LkVdyBLhyONhyMNr9qEkp/P8sqRhbYfnfVceX13GXOX9j/J7LtnsSBpmt0WU0AoN8hqigJoXcNtyycAYTualqDo5vTPf9dczx7OQZYsjt7QKAOfSv4UsRVE6YZ8YImDO6V25xZbk4TkBfHDsY7Vy6OfXkl+K0uPbuWUSc9+iDALJybxfbAPtDNK2KDPfc78gMRzFtSRbvxwmofYJhTMnij6NTkS1bjgpn3nysC7HseR6/mMvSqhzYM91XV9Nm32rSbFCHvbLYAtHWh7t0p2008jDrycO6sGmRYsWLVq0aPFyERAbefiuRkASgGJ5ps16hlnNdpkFeCInYLO86L6TtjUxJW9GSRSDVgsdMy5DYlSTOmY3JN518/37ckM1KxSkgCJfFrFhSfgomSGz1hkZzyZnSsZjnySblUzVZtZvDCVic7FIfViJhcPZgmPO+6fMvI515zIVX16JnLlxzESiSN6+2DrzTWTUxaMX/RURIiO+Jn1i6TPJt9diGfBETB8/WGXEhqTeTFxOr+N1SeAlUmT7wLYVKSJqcp3n4kqQmPfEvuDU+lKY9y+vM2LD38thUdJwtq8QIufKlPKqeJ+84UrfiUoW6KQQQsPMXJWoTB6HolNq44HiNyyfXpjSai9kRkI8S5GJpAWDcBH1UsRmzYuV4zRjdvdxkm6JrtJ0c2XMEi8Wqe2lnNuXvI49x9fKCChRXr5KmxA5n3c7J0iH+b2lc+hYgfiqdUJ5IXuUxI6XiA2Jy5xTolMbZz1T7gjF/qNHnGHQwD1I1uTLzG2fpO2dQ5irwnyONDwjLct8WQmQzpWab+W/OfZlrMgH7OX9dm6cqCE282sntW/I+zXicotXFgf1YNOiRYsWLVq0+BDREJu7GSFELPpROTY207qiAVlXPv7rQzQTogrHZnJr43yKp8CWfcpm6WDnkJt9RpSAWY9m5uXWi6s2DJLdkDfDdemriqHelQipjcKbGJzEviJStazEZWe0ZWCJN5CtBk5key6Izekyfedd4eWYLuXMfUkBPaI6LtusjN2XqpL/obwAg/7Myn83LH0WIT2TWbPt3qE8Gyl9frZOnJ73V5nnwuz+vbUgNmJF4JGixVVuh2iRisZ5pMiWuPu2h9QOs78zKW2neCEwF+bT7VahxtnYZ9mwZJksVSXiaMXqiF4MIgFA9GV3FIptaqvMqD3KECqJqgq2xdTOpdSGk+tEZM7OV0ZqBE3alie1hds9V4TXJtFIIpEWOeI8ZIuUOnevoMjwPTkuio5tS+uRGgLM4zxw3p2VA2D5cEQXYjE+7aJB4SZKGIzlfUOtZ7xEhB2OIjW+3fwdz63xQ+P+LbLMsfI9va+6H9uXvHa0BJ08Ly8MaH6vfJYM9wDI99BdYSOBYmwz40+2W0WrKDtSL/e2HJsaWvTa4wE92DRcrEWLFi1atGhxb+KwEBukigZvZwDYNWSWKKTNPKsx+3OVUr6qgSJ/dm2cT95cv1WpfW9E2ZuMS4TMaG3QC4rhDd3sez5TZHbDbMdWWHAdnZUurHyhUJtHbgCgYyWCQ3G4Vm35S74CTBEcqZzqVpLRr/L+M8eCWSszUZfR27GzwqHC9wAykhNra+SKUglCsC3NNwErwy+vWSEkJo8Xwp+x1VjLUXhE8t5aEJZchQXZGo7N1Vi0PcyQIrP2LlVdtC+wdhFAFmi0HKXO8WUyL8vNV01UjAiK8GZOpKzr0YJ8k3zcr5epT+8epe9cHKU5mLwNB4Bx57gPPK5MmqVrFv0kp4aik09FnPB8W6JUhUmsHDsVflT0IhTboi8yVvLHiEbyPLaIDe8pfecuDMfVKwxA5c+eJqHb8rtxmp+3WSSQ3KOSz1QIfDqjyRm3piL0OdLOQflq5a2+yvVw+1chxgpioxVZbjfTTPRy/l8MkVfye2agWAUdif4zd3+fVcUiV2RN0vEZtyd9iKIBtsPK1Yq4qece+bbHisoihT13lerD1xIRD6rcuyE2LVq0aNGiRYt7EweF2CCkTEorkUwqdyxpOKsuaExW5a9wd9RS2U3utfxIsmQiIkDOdPKaafls2DnND/u3XyuvWh2wKTWVk/fJO4kl5wbIWiTsb+e0VTJyM0cvOod4xCoPoOQTcd7JR+hZNVPh2CgPg/Pl9E1s1g/h6tCI0PM+9NiaDEm767qtlSpjTlOC0yvyFgecR4uaTNIHngOsPOt3JSpG7kd6z1Un8TeOg2Hb3imXSirbnNFodXncV68M5UkTqxYaaUvezKNlQi0+tkzOJpa3tpa+HIsNwvmC1VFyLBeG+7LPiJEcBRmntRU4E7sInle0jyA6thNLgtU273SmOUWUSlDRaTG/lsgNWgk6RV2mx0J6stwIapHoORbL60PRz63hL8nPx5W/pkqkxp63g+Nn0MqBOlhFVaXTbvHtMhO3vJ/dzukUyWQMDomyzhCsMPImupkLaPq0pyKLKMxOtWPmNhKsyPLndE1vy2oxAQZFUo5NBalxFVlex6aoitpTTMtra+pKdD31qfyuN2Ummm93zT6N3oj3NcZDKvduiE2LFi1atGjR4t7EQSE2ARFdiFlDxqS+zDSp0xB6rlWXWUh1v67SQZ/imcmbddErSXnJD8icG8fst8u4HrFxhno2qfG6MlS19JmWbW+KZYbI7JIIgTcKTG2ykkbGLgiEZnoGEZo8KgXycCTrJK+o0N7gVhAazWpDsbW8nE5QhH28j8tBqnQWOTWl6jGrJ1Q/h2vxfe77PnTnpuDYNctzXK1snlf8qL59gdAqGXK2ZI6D5Ww5tWg1MhXdGvIDdkemmk+RMlEAXqXr5a1lMqB8e0iIjUVBrxbpIPJ4EFGb2JeaDpPnogjaQv4XuUQA8ESUpAf5MhEcVvnxN5Y7ovNNjgcrDJ2SNpDPwSjn+pEggcoR69LWVq8QvVHk0qEhPZFOg9Dp/cNzn9xht9cwKzt7+RKr7wap3NpWeXflfoPrWyELxOt5LPVT/H2qQBscUuPRkWDhHX/opS83cWw4t7kiq37PvAmx8dwtbb64H5aIjUe7C+0ef2kqakTO1rwf2r+Khg6QOTZThWNzq8rDDwixOagHmxYtWrRo0aLFy0Z8UA82bSmqRYsWLVq0aHFv4iARG8L1i26+FEVju86ZySmEWXuUc9CuQskC525NeSCFrtaC17IUfB+8mnYkzbCM0pHyavA5pP8kzBK65rZW3ujFwzpnGdFtMhSeTRHLJSnKv1uS52zZDR7Olj8MkW+fuBcJp1oybITIPKH1rUVaJlErByH1klwMANeylJX3x/0Lkdr0k2X3St5WonYpkmaXDFj2ThkBFZpzsHlhAEoYG6UQY/W7cpwp4a8CcTRqFFKpXX7jOcLlJVobdEdcM0gbK6AnFdWYjqTMe5WY5W8IkfZjw4V0PY+dpc8zQj7DvAwqCihTwCVOEn9lWYk2CQBwLkJ8XIbh8dWl30opbvRLBQOlFGS7nBOmuYTGudRxKWM+71+XMUYKS5Yl3JksPq+dDfyNl1CgSaJh92pZOZeipG+Ui9iYpUdvcqlzDBY7oOgrAF1G59LH6JaidMnTLnHrPVO2WuQg7VSOg+/TJHPghQHt+GfL9nr/kNdm2UfvH5GUAq7HFbsoVpS8OCDb43k1GQLvrNTcLT2xvbokR/lbPXfYrvmC78trj4iG2LRo0aJFixYtWhxiHCRiwyik0DuWBhOxIdERxbYqje1LhfcYUQI5A7mUemWiGCxhHPdkbbZPyRySAAAgAElEQVQdnwXYslwVUGOWKejEicjAk0i7NGiVJ6/NMjspP+6M/HwQkrCK1UlmvZXM+sKUPJ+L2h7nmGPPRnrzzFq7sgfZyDYT+UfHSyF3Cmn47YWUIEtKerEUK4dVtnt4ImTYUUqQdyTJHjnbCgCTKzUfpQw4LMu5fiTmiEBGbI4GlrbTMkOQokqZ8SgkXjXBpEihlrjPCdPHQtB9UxAUtkvkbLnKx/tyJef2ypXQH8uYWcW+Mn3i4SSRVhCC054EbRmfKVLtnaIXheZon1DKB6AYcy/ih3p+iTjhep0z+Itt6hRRKpqeju58qpXaKhFf5nZaVRAbGTOvJSJDvG9QXJNWDkDOrHnNe8RGycO7PDcq5uZLhxUJjMU4gUxgZqh0glo5YBaeNJxtYsotkAsJeG0SPfBIjVpHAIrYTDPyMIsS5hd47lPZHudvV1RRyMajFq7MvBTtlC1lD3hqOMTIEoI9cZlICkvQb0bVyyKEitYe9ht/ytgryHYmLt9eufdDEug76AebFi1atGjRosXz4yHp2BzUg01E2Gv7riXIztpgdKiIzax1HdVnARLBGVECOYO+UsRGuAOK2JQl12k/OgB5Q15W1m/5d7eHF3AqaIJFVMiBiPtECTWrmovIZQsCmQuxF7gyZbkUTiPXhQjOtUepxtyw8pX2CGBpBmTGziySY30kAmpHgtg8FtPC00VGbAINGlfSF0m+WercG14I39NSc0FOetnHY3J7pAQayKXIRBNWKxGrE26PtmfQkeGoRO3Idckl17lPC0Fi3lymsb21KIXyWHZ/YlCqC0FsRkVs0vvjpjyfRyN+qHwT2gtQ9JDoRSWdYxm0IpZ7BCABYKApJX04KVegXC4Z+y4fcIoeKvK3L5O314crcVcOyYrImR2znNt9hU8CYEOeXJwjslvhAs3sSbbktczLpBkenRyUP5UhlcyxkevccW7sNTwr76ahbODr+T0HThzQ3x+J+A7G9FTtHDwXLLit7ZOTyIjOONhKZdhrHTBI8w18NUXYY4nq5C+USBGQrxlv58A+vZDY5QyxmZd7RyeiuU8Y0LY93Wa59wOKg3qwadGiRYsWLVp8iGiIzd2MGFP26I0obZCbQMQmzgTc8nfVyoAS9HvWTu2aLNdKNbMjYiNIzsS1eYteuKxmxr+prN9SmI+8ACI1XJu/GDIvgBU1ukaufCI3HnNi7xP5YqZ3ZapXKHFPFIFCapfCkZiEl2OKlYxYmWS4fg7YFVtZ48TpVIjRmXCyggTIyFa2DCjREXuKjPqZfFdQlyNBYd4gaiJidXbMROjI7zmboSa2AomITdq/IkWuMgkAVsIreoNCeYrYCDIknVWRPABYcsyxHJdTry+zTIqJld9RVEazzHzcL6RqiUacFH0jb2Ywxp8dRSBHomAlN4WVO6O5LkaH0Gi1oxpCksdkxAk51mWZYo/OXBUw1TyumovcB8r9n43Z9PRSrmPaOeh57MT3CtFFV6njRTaX5DMNGeKinQMtWeZWDnn3Hqmxdg72/eKe4/4Pq5nbAsDKIDadM+DMyAT29km3NA7elVVlO4NaLFG2rfYzDpmLFQQ7OG6K3r+UY2Pu0bR1mEqIqFaR5BHkWR8qHBu/G/7/wPZ4LdnVhZ0af94ix+YBxUE92LRo0aJFixYtXjIiXkoB/dDjoB5syLHhGupYKx2QYIYyr8axFUjkJITiO7NMxcTk1lN1PZcaCTegF2rgR80PRS1Mn7huLlnMUjIsohXcnprKHfJwWPGikujOdLOrSODPNXwk+zBr4+TSMIMnv0el74WXQz5C+pvVMczk0/sj+QCcH4Neea0NInLkf5CPMBgyAZEtrxGkx9lcy6PaPMh2WZpuPpY5fXPIHBtFThYlcsJKKvJ0xgygKd9GOVqqsSNdKirBBC1i231qexFK+f+V4WcEr9E0Q+YwC/YluqoVIo9nU0IttoaX8EyO9+UmfYeVTeJEUBh/DtelJky/ke8STahYT3QOReBrmp6yYqw0V+W2rERSbSiLfrrTnRk0q/ouu7Q9NwePBqijIA8Dq6J25bVbhKvmmWaIjVS+9Rmx4XnF7N6bntrrwhtwFvwe83mBBLtukl/EqkoiRPa84rU0UkPHwXsWBfKGwd5Qdtpj4QDksaoWl0NspsK+wN0r3fiU02V4kNQWI8dmZpNRdIbtcMv/Cxz/skKNCXqfknum08+xiFHWOZvv5/VEfFBLUY3J1KJFixYtWrS4N3FYiE1Ma7dUtVybVI5r1czy9dk8uOy2wrHx2gVZpZi8hPmTLtsh658IR1REwqx3qwKw7NatidtMa26cR56J6NnQwM9kmceisdItmOEKikCEgONc5OfYPHb5zh5NBiCjUqrhI1ntmhUuW1a+5N8zOVXERsbcL8i9KLcAsN2VKAIVj1edNxy9YcHbvV3Tn8gmdmlDLQ9mryfdPLM+kQHxOzNuj0ECiR4wu88ImnzBqMrmTJpquKXWCrNaew768zHzpUp+hqVTsHKGx4q8GfKlnu5OAJSVJM/EpJLaM6xw4jHtDccmSJWd74t2lTQUw3ehmjOrhahMeyJVbws5n9emimxyfCJeU1VFWk6HXGdZDTedv2ddGvuFuZbIlSNXZJ+hpW1n6ktF61yBJhwuV9kIAEfBcKaQrzfyjgq+jLtvqB5VV94/ij4q8ltWQRE14jW0MBybXs7L0V1SuQLKcKqoizOW9zKqRY8OzbDBc54Gr5kbyBPXIjZEgtw5r1pjJUIEzDV0upnxZ/47Kw3L7tw2/2g2DD2/FLFx3J6tufnsq+Z9rdEQmxYtWrRo0aJFi8OL9mDTokWLFi1a3PeI8aP99weMEMKPhRB+M4QwhRC+/4bv/W4I4e+HEL4cQvj1F9n3gS1FBezGXuHya7MUNTnoU6c9lNupYrCmnzkCIstFeyM7XjWfhCGcq9x8/swb6OmWAmK2CpJCV66ZflYCnRtQca/BLY/QLJGGkKZMNBvc+SWpG5bf4CDXnZOdN2P2yxVBJOj7ZVd8bgnH2006ADRHJKlTpe9lCYFkZiBD0KpN6MtQ7TBi5T0TFDGzRpAkLqvAWed27M6v2nux9h0Xo8sx+DovUeTPozNo9EudXDIYLZmbZddC6r0SQvDTTVqOOe3nS1EkibP0edAlkbKd4m9fCuuW/zqzDEebChpxkhR+JkTtpZC6rxf5eHDZR5eW/dxagiuJxTtKMqSxkRTN43xuxC5ZBKDrMY70nJcuzPFwJH2tmB/KJU5akqS2ZfmFJFuUSyp2KdsvYatAX8dxzk9o9pP3LhKyed/geWWtWTo9t8sx12UqPHlY3tb7V2UJypH/1c5hHxne9oUv/fXNogSzfrbVJaFOxliWmRf3NrbtlqS4BLbvXpG+lDYs4fbCgHYZTufjtkww72b8BoA/A+CnX+C7fyzG+O6L7vigHmxatGjRokWLFi8ZEXeu3DvG+BUACF5c6yOIg3uwGacOG1fOB+Sn440KITFVdDuwxDFf4reH5Gml2JnhLDQLELKhlmnOyzX5VY/YZEO9OPvu6Ihoo6s3tAag7F/nbCQmJ1pG5AaYkzDjDWNmptW5yeSYWUpc2EjomGndIGiSSO93LAc2qMKGBpwi/Pdsd1y098wJAwLAKPsZ1MyznNsSVXBkS3lN08X1NCelk+CYzU6FJO6Oc0EAJ3HWiR/6rBYw0vPSNg0ZFU0Q1MpaXMRtaW1AP0VFx6R9K6DX09pAtuu1kIdFfPGJoCc2y9Ryfkek9ShG+lDeIhF/8Nv0ZRKCgWx2SquMhSI2qU9eeDINoGxPm+ecm8xdiekyX9ciOknyO0uer42yococzAwtZUuhT0vE1+tMXsu1RME7vT4rdfiUrCCxeUfbCkv89nYOz0HH0pfL65jEZZLgeT9ZWkuF55CG7X2Kx1nPbQUy5V5KKwcDefCeuezLtoOShivjcKG7c9ddtMKP7t7JOSBS03V2HHW0SPtQMY7MfZD7x1haj5C0bAUo1RDzRgjoVUYE4sG6YEYAvxzSAfzpGOMXnveDg3uwadGiRYsWLVrcenzccV6+4B86QghfAvDpym8/H2P8+Rds51+IMX49hPBJAL8SQvitGOPfvekHB/VgE5HWMzeubNOGiveRh8DMuvag7KwUtNxbsks4gS0gl6aqiZzLBnTXlbX+LGrluDam8pNr1soLGP26bbkmDxjBK8k+dlw3Hlxp7NJkNcwyVTxOOryojFmyeXIFNAPt9mcAGbUox8w56JV3lH9DJOJSkIIn2xKxeSqvLXoBh/x4YUCL2HheD40/17I/ckqsxD4zaiInzKxViHFb8lzsmDh2LW0nN2WbU9K1oAjPBDl5uizHrPYVmzxm8mR8SX1/PRXt2hL0ThEtOSe2pdnp9XLOC5iZBnptfVt1T+6LM5vNYohyXhnlyseLUpSQXCqiC0sl9Zh29vE/5HVveGQ8vrSC2OxKxIbci43hbI1e8t4LfC7m5dhEQlWAUe4bar4581DJCCDRQTWU3ZVWDvbvLMQn1zcnw1vCmD5kccDSloT8xBqK5OdWEecCsfF9cr8N87F7A06WmneO51Ljovl7aOfNhs1x88cwc+d4n7RwunR7j5XCjX0hh0vaI1KzndWK35H46Mu9340x7iX9pibjD/5BG4kxfl2274QQvgjgBwDc+GDTqqJatGjRokWLFncuQginIYTH/BvAn0QiHd8YB4XYIIZUGeWkq4Gcaa53RGzKJ/qapLsGn8rd03sYyooCIGc+3F7JlsZ91TVUzYCYbQjfpFJNxGyf6AURgkuBVrjdGqVBL/7kCfjKCzBcBbUeUOSGlVqS4RnTxVOBBjjmM1Zh0WiUGVdRzeAJACWCUuNEQFCqjWStzGKfBZH2l0z72hh0Bsc3Ia+kFy6P5f0MC6I56TX5JhuxhCBqQi4PkM8xokUXG4G4FCmCtJfb8aKERMPIcwmGV3QtSAzb/mB7AhusWrq2iI3yZfaMmYiN4VQpP4OCdjQp3GMYaIOHUvkI7tyxbREl8oajUSwoaCEB5POKopM9K89CWcVSNV8k/YYIDUEkI/joEbKtoFTMqHcV9FOR15kVC7dyDAuxy/IzNd8sR1G0w+s3V/ql4+utHNJYy2tGM289LmVfU/8FsZF7F1ESCkGSY2M5MLpbN8ceaQYAqDmlfDYTwSPHJx9LVmR1ityMxXer3Bp378x8vvK1RWx8RVa27ig5T+nDsm3lW1Kctfb/hpsf8jmV28O5rfxn8wp4si8WEXeOPBxC+FEAPwXgEwB+MYTw5RjjD4UQPgPgb8QYfxjApwB8UQjGA4C/FWP8O8/b92E92LRo0aJFixYtXj7umPJwjPGLAL5Yef/rAH5Y/v4dAP/sy+774B5sYjRZpnnEH31VlDDUwz7khjtL75ZtcP1YntoXFcSGa9VXfcq0BicPfpMWA1SGfH/lDnkB1N6gAeXlIMiN4RfNeEXu/M32AiarcfYRTNiphUObBiAbbnLsjxep7dXAKpYykwesJlAotqgZcbqYXKXQMKU+EY1j5g1kHoWiMOTusPrKzC0RBlYPEW2hfg6rsZ4Y1IS8orOtVCcJcuJ5LkUF0nXZdr8sEZa+ot2T207IUC8H8ULa5fdSv0M5VuVqlYhNeV7JlmaYjgvBLNpm8Mofk3Miqj7SvMquW7G6pERstMBMtGjUsBVZ3p/2Al7LZ/KVhmYcqu2yLTNra8zZu3maHErF/dvMumfFDDWhBqJuJXfIXmMzexKUwfuU1d26FCiLCOyVt3IwKG62D9CDJu0RiUXRfup/+u5qKG1CWHlGPpPlVKmGDu+ZTrfIcmzYpX1WEzx/rYYMj3cn1zUrDgONbL3/BixSI31wFYa1SsM5YlPq5/SmKipr6JRzm3fGfpghuqpKcmyep58DlBVZLV5dHNyDTYsWLVq0aNHiJeOOITavMhp5uEWLFi1atGhxb+LgEBuLFFoBJC5bUCRpEphwcJBlTQxv2jMLQR77FgUBrhS6IpRMyDcMhDbzfhQiZudDCWdbmFObGlmiWjpe027gajIy8LJEQ6h9JkRX4zPPiIeybDLQcTkvRT2io7iM+aJPfVCCscL2ef+jW65QQqsSmeclqvtC3XNJ9qwsTfjjqw7qO2OH4Qi0JJqSqM3S57NddnumQ+85l4RkGSyXl6fvFUsgugzmbCS4ZGRc0HeuxJ0EaS4JcYliMiTrXmF5LjnJHHhyqQ1HAOYyK5cMjrSU3ywdLGjV4VzjSQw2agvdqrQ6mJzwI68Lu6y7cCRhShnsZCcUyyvsBShKqKX7JbHWkoeDlknLOefE41S2wFwgJNt27pxWe5JFeW3ZMftzmXvl0sSVEQLk2Pnexlk5FEsfPJ5c1uUK3UwE0fxG+s97F7e9U5yzhGZeV8GdX52TagDMfdgve1MYsC/tJIDSBgaw5GH+mAMz4/BSGbyF+nmqWBVwCYpCjDy2gxUldNQBv5vaPdo7vkd3ftWCS383SWS82ogPCrE5uAebFi1atGjRosVLRAQwHazy8EvHQT3YhBDRdVPViDKbBUpW40jDmskXsv98z2UKbvc1QiWRG2YhzEyCIx0ChtxHFMdlB1WjRhLSxrL0mUiNzf6Y1edSUY6vHKc1y/PlmUzcmFHYEnciNY/6JH1/OqQSZBKMg5TysmQcMMJsKzkuJGwuy+zSEo6J/DDbUxFEJf9RDHHe7zwubitjnspzQLOwkXMtpOVxflmsnYCaN/604MM+Eq9HEIB8ns5l4AUlcSiDHaOPLHDmCNswqALl/oXMe7pIx/ZUPCksOvlYzCmXy/Te5ZJl5Onz0Qg+kkQ9+fJoJWXOO02y8FY6dR0pvJlOnjXn2pZwO8SG6Jh+bq/vPdczSatLhyAAGXnNhrIyZi3lll0WooEOgZX3ibKtK4KiisjJWLfefPMmlHUgOuavJYNGD2WJsy899hYhqb9Spszz1JmdFohNX154GZVOG6Iixwb5pQmokrb1ei77VhM3VaQ9lO9XtA9NqTnLu6XkPdIOZ17uPRPkYx/0PmL7VP7fQnuV0ZPSzUC8iGqLVxsH9WDTokWLFi1atPgQ0Zai7m70/aRP4n3lidibUeZy73K9GMicC6IrylmgsVrlPOAaNdfImYUcedE6W/rsymSnvsy4bhKmilrGXmazlP8HjDHjtsy4siGkvK6gVfssJ2yGx7/9mJndKh/BOB0oD0Ml6Lvi/ZmVA4AgKAL3S94Hs77zIfFcBmP3sHNoWC59lXYNvOOFC18meBy0IFUzOr6O9svF9mV873zZNbfBGvfpWMsS5LgoT6TJICp67jkzylOxNXhrcZneN9DTuZT1U1TvQtELFFv791xaIBYfWDHJtRp/luKTPLc3UuJujVK9+CHRMaKgoSZC5tBIolIWTWAQhRxkfq4dsqjnV+0c4rlA7h5tUVSywSI2pQHnbtxPNlPJBEVo3HGvHI9saVCWd2eUTBAbK3JKXiLRSI/U3MTd4nwIQrdywoBARrnZtr9ne2HAsg8cF4rv1pAtb+dAVJ0crt4gNsEjNnzfcWtqHBtF+L2Fg5NQAPaIA7Z4ZXFwDzYtWrRo0aJFi5eMhtjczQghYuimGfcC2GPmhnkWUIqWedY/35dMaNrPdu8diqFr2YLYjLYqylcvaKZVVggB87VqBvtAOXjLA2HlFEX91HRx67JbU7nTEUlxvA+K/O1uKFfqnPQ9M+HShA/FWKVAyFTUcC5MNZHjfbxBjocc7xNBbKxB59XgUQTud66WFj1PgseIFWGs6DDQFs+x3lXQ6FidFHv6conuZLn2cgvkqgxWC2nlRijRq96MeZ9o3LgqbREsB4bWBtNK+CViRvmIRpRDMqK0iM2TIQkVqhklhfqqPBMOiG/ISwoCssrPXBgUnTzrE2dLK/+kAm0n57UprKkI88nxqLoVSldoRinbk4GcsfXsJ/xMUUFWObrxlTwQ3j94TcmX9liEAMCqT39no9E9nBXztyI2i9K+IlcYmnubQwZooUC0hGKBW8uxoTjgTARvjjzOqyrlA3c+20oonlv+Xq1ikYp4mfuUq24MHasssTf8dcw+sNrOoiYziwOHxBIxioYX5+0dPLqu/bDcTHevbPFq46AebFq0aNGiRYsWLxvxznlFvco4qAebgFQxs+zm2cC6k+xuzxNxzgbye8z6mNH5Chc+pVsbemY8o8sQNbMXxGZXGE56xEaytMraeEYTJCtwXAvdp2mfbHzlFXkpdPIStnlueqI6u3LsrPqxiNBaOsqxT27szHqmSpbJbDJ4o0BFdPK4lo738XhIiA05PReC2BwZI8UnC/J75LgoEsQ+GtsNGjUq50nmmNo9lJ8fstAMNTDIvWC1zIZ8k0V5bAFglOMbyIVx9gJ2zKxOOhHNmMdiX0EeCNELViYB8+qknVSe9SvCCeX7aexpG+W35M08GsqKt95k00TKtIrEozIm9PR0VSydohdpaw1MyTl5pohNaS8wCWfMgFX5+vX3aHe+AfPKu8WCx5djnyM25MoNM5n/st1Ce0reU52ibTlmoiKWz8JKKfI+Qi7nK8Zjx8FzrdP3SwTS/sYjEbx213ISkt+0thWAM/sZh6yYfUavR+WMg1U/x9xwe3fQtOJPbhwESgsuIO9hRE787V0vc4P8auVq+f8E56AvCDNln2bcmgrS3+3R0JlVY5l2iMSqWfLrjgjE2eTd32jKwy1atGjRokWLexOHhdiEiGU/ajUDTdWAnIHouqbXC9AncbMO7XRG5lmmrEubdWi2o8gN1XDBp3ZpropeyHZRZjk180i/Vq1rtKyaMdmPVuoQsXHrw7mCxGQdqtwqrwUE2ciYmTUD2XCTRnpcn2cGGl271dizJm+1N8hrYDXUY0ERWE3xhvBArJEiVEMHxXakEq45D/bxewZBQ07I7TGZPNfpL2QOyE1ZC/IxqqZL7tK04gHnZ07DZ2mqJWR/1JHhGJn1vbEUZetFRqkuiNhQDbeiK2M/B8w51pdZpaIyTgUYyBWAyjHzqIxFP3lIYolaKIdL0Jf1NneKKtq8nr0Kb6SSdlF14lAENYIsKw0Bg9jIcaZpK01dT6wEtMRK1XDrJ3NXyeAVyNo5xENekytkK8LmqCdhn5IzZv8eF+Vv9ikeA5m3slOkhvfHUgdrbcmATkNnbnBpuYDlPUwrS6X/ywpi0zliDPtGxeNe59Y0qXyWemVWrCA2C0VsynNbOT5VKXbZzlS9HZ/GfGcPrVPvzS+zovBa4gEtRTXEpkWLFi1atGhxb+KgEJsWLVq0aNGixYeIVu59NyMgQei0L1gV4k8COwqeSeGlfeXTQBZVC06YT6HQSokqCY6XUj+rhFqVvp+348s1o8LNJJ6a73JpRrZLtzxD2P6yy2sfXH5ToqMnwEmppCUPdzsv5id9286Xop7tErmTc8wyXbUZkCUDy4vLRnFln2Zhjkun5LsSOuaSAdsnkRfI5Fs1aFyUy02W8DgTB5QlnaNlWfrM5aBizEJcPlmlvjzbQ1oG8jIY4XNdKmLJtbGeWM7avi7Gfkb7ioVZfnOCiNERsxlVI1YthRWyuzdFNOsaeblV9juScC67Mvp2vvyaRNpexPWCbGkiCuRz7Grg0qa0xyWb2rWkJdxOMoHH3Sy/6fw44UdeQzyvtmai/HIJrUeU2LrjuMxSsHRXycO6JCV9Y8m7ORlVMqFzpFKaMpol2ixTIPcNCojecG8jIZf3Li5FeSuHXaFLUS4p52Wecq5TX8p7F6ewczYGfaUum/Otpe5jOV81SQ6epvPlMfYxv5VLzUU8VbYUgLxpKSrTE/i6fN/2QXcjf/Ru7FY2gkUI/Y3r9a8wYsRD8opqS1EtWrRo0aJFi3sTB4XYIKSnYhJpj0zKSCKxmlEy83G29DeVqnr5bGYQKoAHYD2VGRC3LI+eKtnmrGJUCY/ptc20ScJj2eTKmclxe2RUy2gamcXiZDsTn5oTp5Uk58tyd7YstyR5kkh7tXVluTaDd2Xk3T7DUUPKm1XwhhJVWGgmNJdE96J7WqEf7dzKVjN5yZpl/mrCbTzHSDhV4vKiRE0seZjoDbO+0Rt/mnJvT2glYZqIzbGgCgtjShqG8tz2Am2zjNsGSaVTKdjGMmArY0B0kiXawaF8FrWgGSWvnV7QKV6iRARHYx1AUj7FIGdCmIq6GvTCEfBzCX9p2QEYwvRA8nCJ9Co6YxAbIlZ6HSsRGMXYg7FmUdDAG6LekCDbbB7I5cCKQBogxYtPqhyFr+k27UUn6GkNOIFM/B9rpFgGL1USss13ZwUQMl00eayhIjzXdDs5kriWe1cQm0hUumxPy8wNXOzvE7yWPDppd6hmm67UXYUBB4tsuX1In4jG+DJzALgignWb5OEHtBTVEJsWLVq0aNGixb2Jw0JsJHz5HpCfjllmqGaULpO30vdassi3HHLDbM2WadK0juu1zHyulW/ixJvMfmdxQ2l4p4gNy4HLElVr4KeS9w6d8uZ5c/3w3LcssCYlqqMdM038VrJNY18T1aH55jbvX+0cvK3DjIeQu6J2DhMzOimtF3W/cV9mb8OJyFWl6fUzIjbMsOZIINfpOe88v3h8oqInJptVcUWX1aqQmilxp3BXKLNLogmDs3KoDtmbbfrz2P4t5yfPaWbyZ2Pi8tgyZD3OgsyRJ5MlAuaIjUrRU56fKAapEjbrd8exc1yFUCt9djYSai9ATlXxXZmXvhRM85m7RamIyBJNyvYk5Zj7Td4HhxHkt4r4ukNms/WVclBKNIn3rclybLxZrjv9a+2Ncn3xGiWnhucVuTfjNM9tvZGptm8upnxOy5eJ1DiZjdH8hmiYSmaQ3zMTBjRj21NynhEbHmNzTXX1a4nnlz3vqFm3zyz5Jo6Nt9kYnMmq/f9pEZ5/Hb/qiA+IY3OQDzYtWrRo0aJFixeN+KCWog76waZ8IqaJoKyfOzNKb2sAGPSG2xknhVUNOesgl0bFxJw0+SQNKmcFmLfVliIAAA1BSURBVJvJTf6R3wQRG5fJaSWHcIlOjYjckTcpdFYNuRrLjN1lG/mDtLGZ3E55RTTzk7HPzDfNMOTvfl1WnvWLWHxuUR5WzFwKUkBUjJUVnHMrmEhxQK2acOJZN2Z/HyL2CbcV83gDoPS8GB3URGSlzDKJCpZZppfCt3OryAkFGDdpLs8EhWPl22hWpymgx+PSOcSmt4KPO1dNNDvX52PtVSSwzjfp3fkMGE6TbDuKLqq1Rf4us3rysHjsOEbPLwLyub1TxImVYC6TNzwQ3j9mFgSslurn4yQnkCgGkVmOuShW8qDKrIJH2jeFc+StbNTOIW0XU2nlUIS7J5Dno/cTc2JPaqnAPpZjV9New1/ywp5rj9g4xLHsm0dqULy2SJEXnaSVA/tUv5ZQbOGRGtOnmSFqVyKNbN8iv+RU3SZi85DioB9sWrRo0aJFixbPiYgHpTx8WA82MSEJY4XzzKx+nxml111IfzMLcNmAy7gtgreLJQ+ESA0zo0kyvEVFbr5za8jVtXinjUC9C6/pYuXgme1pRQV5H848r0Cr/NidHkQts5gcx4WoDjNGm3j3M24NrQ+I7pTIDZCl56lvQhSBcS5CMFZjhxVZC1et4quxav3k652rzil0TWQR3ltoKFKk3KTcz/3ZH1/n47BzOiNED8gL4GurpUSELBD1IsVqj5YMMNeT2UilE002zxdpuzPlOJxnGqNmk0JWjOQugfPseVwuuy0z67Lij0GdokEuoo3R/ZlxbOSnNQ6dh+hGh4Jxbi9NxVBGbEQC3xvKErGxCNXC3Y8cv4v3oqPBVnGmv3kt0XxTuYF2lzNOTVlhmHkh5ou0c1DkJh1vzjnRl+I631NxRHSmsHHxxpu+QKuC2NBMl/dM1dCpoKvaDvvUu/uV4611pipKNZrcDvn/huX9+Iosj3qGG5ZvvDaU59aUBqBzK5wWry4O68GmRYsWLVq0aPHy8YDcvQ/qwSYiZVuKGFTIDMxAOqc8nM3a9qMWqgfiqldqxURbZn+yJXoR3do8MF+fV7TC6VLoIIEZbNS7qpmFgUdmlWDKQ3DZ7dJUKDilVuUxaHWDUSlWRn95YWgycxNqQdVjVT9mBZXs2/BApk2a+IttyqBZhcVjytdWY0crssj/2JZoheVCKIpEg8ZtmdVeOG4PkLM7Vg+xyoTHWatmjK5JVuGV18rT4JyYKhxBBtg2q5N4vJ/JmMmJSf2u68loZRK1ZIx0iZ/v3bZExzi3ln/AsU6CDPQOeapms3sq8lgpRqNTICMz5Isxq326SHOwoOGoqRBSTRfHsbhJs4dZOblZV1M510QCgTwfkxzfQav3WNZVaYfNOX0q6udQg+jUmPZ61eOlN9+0SIpDNPJlWCI1JYImx3kq71MMIgcFZ4z3TKcm7BWt01j9G/JdVlU6800ba6ehEzz6WXTUoel6/Pk++S3m+nNoN/+f4Lltq1wVTfVIjT/eZgpy1am044xlZzpJpk/e/LTFq4mDerBp0aJFixYtWrxcRACxcWxatGjRokWLFvciYmxLUa86Qgh/DcCPIAF97wD4d2KMX3/e7yICdlOXBZ4MzDk6bJQQK9+eHNkMmJOHtTzaiT+VyzL1skYtj3ZiUwCMHDuJh/tLR7nCRONBXeraAyUDuZSQZYeTh22dUSCQl6kywVj258w3gWzfsHAEbSWC3rAMkIXhyrF2jqyXviPEVlkyoHUD5/xSlgko8Q/UlmVkS4l/s/9pT6n5djPI/lN7T3fH+hvCyrokJGJ1cEtqZpUhL4OxT0J+zeXSecI20va5K7vmmEnuXdsxy+91rBTM45hl7sdtPmc61ycKSW5dObANXZZyiZ4XcAOAOHTFe9kGQD6XuV8YVr23sODyG4m0XMIpDCGVsFn2JduizMniNFlkefG5HOelM3UF5kuNe4n+Zn3aL0FlYmt5LVmS9ImMOZPF3dKQWVvLFi/lfYPt6hKnXQpWA866uLwaNVrnWpnnyZn0qi2NtWzw6/M0AZ7Ka3hj7tE8n55rdlos+8gLvwTl5tguqfHeOOkS1B5hQGBWar5PXLEq50A5Afkyzyeex31lbe1GcdEWH1nclqXCT8YYvy/G+EcB/I8A/qNb6keLFi1atGhx7yNO8SP9d5fjVhCbGOMz8/IUddpY5XdJFE4FnuIcsSHCEWepnOzDZpleZMoJUvGxz5ouDlrK52W6uVNpzmQ3Wiar5FUKj5Hwmr9LYinJi8wuvIFcjYSWUaoyo6ubbTpisWRrnWTUKvqHnGmS5EkkhfPCrMk+Jt9kOlr02aKj8vfosr6rIEiNK8VNL8o57J2VQyhIvURMSuSGYmwkLZ9tc5n5TsZONEcNIZ3UvjWE7NZEpYTYuCm/awnTHIu2TcQGRKnS+xTJK8ah5d4OAYzlawAzuXpm1jdlkIoi7CPiL+bnE5GbTFiX3whic7TI59WpIDY0/mTwfOu7G24L/jrj5V6RWSC6RtSLc8p2KPwImLJ6J/SozeocVK4l3coHzsh2ZS50FgFMQZA/GQjvWyX6IuNx9xEVp+PxH82x3DN1LDrgfbIwV2XxgSMP52u5cq6440CEaO2ETAFg2ZXk3dm558rNiz6oPYkjEzvUxAb/T+A9k4TmkjxczvdsN4oM2v83yi1LzfdZOdhoiM3riVvj2IQQ/jqAfxvAUwB/7Lb60aJFixYtWtz7eEAcmxBfkX9ECOFLAD5d+ejzMcafN9/7KwCOYox/dc9+fhzAj8vL7wXwGx91X+9RfBzAu7fdiTsebY5ujjY/z482RzdHm5/nxx+JMT5+XY2FEP4O0nH5KOPdGOPnPuJ9fiTxyh5sXrgDIfyTAH4xxvi9L/DdX48xfv9r6NZBRpuf50ebo5ujzc/zo83RzdHm5/nR5ujVxq2Qh0MI32Ne/mkAv3Ub/WjRokWLFi1a3K+4LY7NfxpC+CNIFL1/DOAv3FI/WrRo0aJFixb3KG6rKupf/5A//cJH2pH7F21+nh9tjm6ONj/PjzZHN0ebn+dHm6NXGLfOsWnRokWLFi1atPio4rYE+lq0aNGiRYsWLT7yOLgHmxDCXwsh/L8hhC+HEH45hPCZ2+7TXYoQwk+GEH5L5uiLIYS3brtPdy1CCD8WQvjNEMIUQmiVCRIhhM+FEH47hPDVEMJfvu3+3LUIIfxMCOGdEEKTnKhECOG7Qwj/awjhK3J9/aXb7tNdihDCUQjh/woh/D8yP//xbffpvsbBLUWFEN6gcnEI4d8D8M/EGBv5WCKE8CcB/C8xxl0I4T8DgBjjT9xyt+5UhBD+aSTi+k8D+A9ijL9+y1269Qgh9AD+IYA/AeBrAH4NwL8ZY/wHt9qxOxQhhH8JwDmA/+ZF5CkeWoQQvhPAd8YY/14I4TGA/xvAv9bOoRQhhADgNMZ4HkJYAPg/APylGOOv3nLX7l0cHGLzYe0YHkrEGH85xkgx+V8F8F232Z+7GDHGr8QYf/u2+3HH4gcAfDXG+Dsxxg2Av41kVNtCIsb4dwG8f9v9uKsRY/xGjPHvyd9nAL4C4LO326u7EzHFubxcyL/2/9criIN7sAGSHUMI4fcA/FtoBpo3xZ8H8D/ddidaHER8FsDvmddfQ/tPqcWHjBDCPwXgnwPwf95uT+5WhBD6EMKXAbwD4FdijG1+XkHcyQebEMKXQgi/Ufn3IwAQY/x8jPG7AfwsgL94u719/fG8+ZHvfB7ADmmOHly8yBy1KKLmzteyyRYvHSGERwB+DsC/7xD2Bx8xxjHG+EeRkPQfCCG0Jc1XELdmgnlTxBh/8AW/+rcA/CKAqs/UfY3nzU8I4c8B+FMA/ng8NBLVRxQvcQ61SPE1AN9tXn8XgK/fUl9aHGgId+TnAPxsjPF/uO3+3NWIMT4JIfxvAD6H5n/4kcedRGxuimbHcHOEED4H4CcA/OkY4+Vt96fFwcSvAfieEMIfCiEsAfxZAL9wy31qcUAh5Ni/CeArMcb//Lb7c9cihPAJVqmGEI4B/CDa/1+vJA6xKurnABR2DDHG37/dXt2dCCF8FcAKwHvy1q+2qrEyQgg/CuCnAHwCwBMAX44x/tDt9ur2I4TwwwD+CwA9gJ+JMf71W+7SnYoQwn8L4F9Bckn+FoC/GmP8m7faqTsUIYR/EcD/DuDvI92fAeA/jDH+0u316u5ECOH7APzXSNdXB+C/jzH+J7fbq/sZB/dg0+L/b+8ObSoKgjCM/hMQjwAGh6QIDJJCqIAySOiEoCkAicGR0AFFICCDAOQTV93cyTkVjPwyu9kFAPbZ3FEUAMA+wgYAGEPYAABjCBsAYAxhAwCMIWwAgDGEDQAwhrABUlWXVfVWVbuqOq6qd//YAFvkgT4gSVJVd0l2SY6SfHT3/cojASwmbIAkyd8fUa9JPpNcdff3yiMBLOYoCvh3luQkyWl+NzcAm2NjAyRJquopyWOSiyTn3X278kgAix2uPQCwvqq6SfLV3Q9VdZDkpaquu/t57dkAlrCxAQDGcMcGABhD2AAAYwgbAGAMYQMAjCFsAIAxhA0AMIawAQDGEDYAwBg/dtpOYeTXS1cAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "%matplotlib inline\n", "\n", "u_evaluated = u_evaluated.reshape((Nx, Ny))\n", "\n", "from matplotlib import pyplot as plt\n", "fig = plt.figure(figsize=(10, 8))\n", "plt.imshow(np.real(u_evaluated.T), extent=[-3, 3, -3, 3])\n", "plt.xlabel('x')\n", "plt.ylabel('y')\n", "plt.colorbar()\n", "plt.title(\"Scattering from the unit sphere, solution in plane z=0\")" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "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.8.2" } }, "nbformat": 4, "nbformat_minor": 4 }