{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Variable Coefficient Poisson" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We want to solve an equation of the form $\\nabla \\cdot (\\alpha \\nabla \\phi) = f$\n", "\n", "We'll do this with periodic boundary conditions. Consider the coefficient $\\alpha$ of the form:\n", "\n", "$$\\alpha = 2 + \\cos(2 \\pi x) \\cos(2\\pi x)$$\n", "\n", "and the source, $f$:\n", "\n", "$$f = -16 \\pi^2 \\left ( \\cos(2 \\pi x) \\cos(2\\pi y) + 1 \\right ) \\sin(2\\pi x) \\sin(2\\pi y)$$\n", "\n", "The solution to this (with periodic BCs) is:\n", "\n", "$$\\phi = \\sin(2\\pi x) \\sin(2\\pi y)$$" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false }, "tags": [] }, "outputs": [], "source": [ "import numpy as np\n", "import matplotlib.pyplot as plt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Setting up the solver" ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "tags": [] }, "outputs": [], "source": [ "import pyro.multigrid.variable_coeff_MG as MG" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "tags": [] }, "outputs": [], "source": [ "def true(x, y):\n", " return np.sin(2.0*np.pi*x)*np.sin(2.0*np.pi*y)" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "tags": [] }, "outputs": [], "source": [ "def alpha(x, y):\n", " return 2.0 + np.cos(2.0*np.pi*x)*np.cos(2.0*np.pi*y)" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "tags": [] }, "outputs": [], "source": [ "def f(x, y):\n", " return -16.0*np.pi**2*(np.cos(2*np.pi*x)*np.cos(2*np.pi*y) + 1) * \\\n", " np.sin(2*np.pi*x)*np.sin(2*np.pi*y)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let's create a patch to store the coefficient $\\alpha$" ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "tags": [] }, "outputs": [], "source": [ "from pyro.mesh import patch\n", "import pyro.mesh.boundary as bnd" ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "tags": [] }, "outputs": [], "source": [ "N = 128\n", "\n", "g = patch.Grid2d(N, N, ng=1)\n", "d = patch.CellCenterData2d(g)\n", "bc_alpha = bnd.BC(xlb=\"periodic\", xrb=\"periodic\",\n", " ylb=\"periodic\", yrb=\"periodic\")\n", "d.register_var(\"alpha\", bc_alpha)\n", "d.create()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we can fill the coefficient" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "tags": [] }, "outputs": [], "source": [ "a = d.get_var(\"alpha\")\n", "a[:, :] = alpha(g.x2d, g.y2d)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "With periodic BCs, solvability requires that $f$ sum to zero over the domain. Let's check that." ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "tags": [ "nbval-ignore-output" ] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "rhs sum: 2.07633187e-12\n" ] } ], "source": [ "rhs = f(g.x2d, g.y2d)\n", "print(f\"rhs sum: {np.sum(rhs[g.ilo:g.ihi+1, g.jlo:g.jhi+1]):20.10g}\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we can create the multigrid object" ] }, { "cell_type": "code", "execution_count": 10, "metadata": { "tags": [] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "cc data: nx = 2, ny = 2, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n", "cc data: nx = 4, ny = 4, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n", "cc data: nx = 8, ny = 8, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n", "cc data: nx = 16, ny = 16, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n", "cc data: nx = 32, ny = 32, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n", "cc data: nx = 64, ny = 64, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n", "cc data: nx = 128, ny = 128, ng = 1\n", " nvars = 4\n", " variables:\n", " v: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " f: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " r: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", " coeffs: min: 0.0000000000 max: 0.0000000000\n", " BCs: -x: periodic +x: periodic -y: periodic +y: periodic \n", "\n" ] } ], "source": [ "mg = MG.VarCoeffCCMG2d(N, N,\n", " xl_BC_type=\"periodic\", yl_BC_type=\"periodic\",\n", " xr_BC_type=\"periodic\", yr_BC_type=\"periodic\",\n", " coeffs=a, coeffs_bc=bc_alpha,\n", " verbose=1, vis=0, true_function=true)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Initialize the solution to 0" ] }, { "cell_type": "code", "execution_count": 11, "metadata": { "tags": [] }, "outputs": [], "source": [ "mg.init_zeros()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now initialize the RHS" ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "tags": [ "nbval-ignore-output" ] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Source norm = 81.3868428575047\n" ] } ], "source": [ "rhs = f(mg.x2d, mg.y2d)\n", "mg.init_RHS(rhs)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Solving the system" ] }, { "cell_type": "code", "execution_count": 13, "metadata": { "tags": [ "nbval-ignore-output" ] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "source norm = 81.3868428575047\n", "<<< beginning V-cycle (cycle 1) >>>\n", "\n", " level = 6, nx = 128, residual change: 81.3868 → 112.091\n", " level = 5, nx = 64, residual change: 79.2048 → 101.02\n", " level = 4, nx = 32, residual change: 71.2411 → 68.0705\n", " level = 3, nx = 16, residual change: 47.6577 → 14.4578\n", " level = 2, nx = 8, residual change: 9.8448 → 0.0171409\n", " level = 1, nx = 4, residual change: 0.0106141 → 1.69329e-16\n", " bottom solve\n", " level = 1, nx = 4, residual change: 3.07629e-16 → 3.07629e-16\n", " level = 2, nx = 8, residual change: 0.0168243 → 0.0168243\n", " level = 3, nx = 16, residual change: 14.4901 → 14.4901\n", " level = 4, nx = 32, residual change: 69.371 → 69.371\n", " level = 5, nx = 64, residual change: 103.884 → 103.884\n", " level = 6, nx = 128, residual change: 115.511 → 115.511\n", "cycle 1: relative err = 1.0000000000000007, residual err = 0.025573219961900512\n", "\n", "<<< beginning V-cycle (cycle 2) >>>\n", "\n", " level = 6, nx = 128, residual change: 2.08132 → 2.02861\n", " level = 5, nx = 64, residual change: 1.43441 → 1.83684\n", " level = 4, nx = 32, residual change: 1.29871 → 1.24144\n", " level = 3, nx = 16, residual change: 0.877452 → 0.256041\n", " level = 2, nx = 8, residual change: 0.180588 → 0.000272285\n", " level = 1, nx = 4, residual change: 0.000187447 → 6.55867e-17\n", " bottom solve\n", " level = 1, nx = 4, residual change: 1.19217e-16 → 1.19217e-16\n", " level = 2, nx = 8, residual change: 0.000298554 → 0.000298554\n", " level = 3, nx = 16, residual change: 0.33718 → 0.33718\n", " level = 4, nx = 32, residual change: 1.72809 → 1.72809\n", " level = 5, nx = 64, residual change: 2.60971 → 2.60971\n", " level = 6, nx = 128, residual change: 2.89721 → 2.89721\n", "cycle 2: relative err = 2.1803634390217064, residual err = 0.0006486396426301177\n", "\n", "<<< beginning V-cycle (cycle 3) >>>\n", "\n", " level = 6, nx = 128, residual change: 0.0527907 → 0.0515129\n", " level = 5, nx = 64, residual change: 0.0364241 → 0.0468113\n", " level = 4, nx = 32, residual change: 0.033097 → 0.0318323\n", " level = 3, nx = 16, residual change: 0.0224975 → 0.00656631\n", " level = 2, nx = 8, residual change: 0.00463131 → 6.95479e-06\n", " level = 1, nx = 4, residual change: 4.7921e-06 → 3.47155e-16\n", " bottom solve\n", " level = 1, nx = 4, residual change: 6.30821e-16 → 6.30821e-16\n", " level = 2, nx = 8, residual change: 7.63309e-06 → 7.63309e-06\n", " level = 3, nx = 16, residual change: 0.00864876 → 0.00864876\n", " level = 4, nx = 32, residual change: 0.0442789 → 0.0442789\n", " level = 5, nx = 64, residual change: 0.0665472 → 0.0665472\n", " level = 6, nx = 128, residual change: 0.0736819 → 0.0736819\n", "cycle 3: relative err = 0.04844393523115633, residual err = 1.659245815001406e-05\n", "\n", "<<< beginning V-cycle (cycle 4) >>>\n", "\n", " level = 6, nx = 128, residual change: 0.00135041 → 0.00131762\n", " level = 5, nx = 64, residual change: 0.000931668 → 0.00119765\n", " level = 4, nx = 32, residual change: 0.000846751 → 0.000816239\n", " level = 3, nx = 16, residual change: 0.000576837 → 0.000168502\n", " level = 2, nx = 8, residual change: 0.00011884 → 1.78399e-07\n", " level = 1, nx = 4, residual change: 1.22925e-07 → 2.1827e-16\n", " bottom solve\n", " level = 1, nx = 4, residual change: 3.96622e-16 → 3.96622e-16\n", " level = 2, nx = 8, residual change: 1.95801e-07 → 1.95801e-07\n", " level = 3, nx = 16, residual change: 0.000221902 → 0.000221902\n", " level = 4, nx = 32, residual change: 0.0011347 → 0.0011347\n", " level = 5, nx = 64, residual change: 0.00170278 → 0.00170278\n", " level = 6, nx = 128, residual change: 0.00188597 → 0.00188597\n", "cycle 4: relative err = 0.0012759605329324085, residual err = 4.2728946362388976e-07\n", "\n", "<<< beginning V-cycle (cycle 5) >>>\n", "\n", " level = 6, nx = 128, residual change: 3.47757e-05 → 3.38943e-05\n", " level = 5, nx = 64, residual change: 2.39659e-05 → 3.073e-05\n", " level = 4, nx = 32, residual change: 2.1726e-05 → 2.09217e-05\n", " level = 3, nx = 16, residual change: 1.47845e-05 → 4.32098e-06\n", " level = 2, nx = 8, residual change: 3.04737e-06 → 4.57173e-09\n", " level = 1, nx = 4, residual change: 3.15043e-09 → 4.90471e-17\n", " bottom solve\n", " level = 1, nx = 4, residual change: 8.91242e-17 → 8.91242e-17\n", " level = 2, nx = 8, residual change: 5.01821e-09 → 5.01821e-09\n", " level = 3, nx = 16, residual change: 5.68972e-06 → 5.68972e-06\n", " level = 4, nx = 32, residual change: 2.90707e-05 → 2.90707e-05\n", " level = 5, nx = 64, residual change: 4.36992e-05 → 4.36992e-05\n", " level = 6, nx = 128, residual change: 4.85557e-05 → 4.85557e-05\n", "cycle 5: relative err = 3.301203447716335e-05, residual err = 1.1068868945958364e-08\n", "\n", "<<< beginning V-cycle (cycle 6) >>>\n", "\n", " level = 6, nx = 128, residual change: 9.0086e-07 → 8.76274e-07\n", " level = 5, nx = 64, residual change: 6.19593e-07 → 7.90474e-07\n", " level = 4, nx = 32, residual change: 5.58852e-07 → 5.36009e-07\n", " level = 3, nx = 16, residual change: 3.78756e-07 → 1.10732e-07\n", " level = 2, nx = 8, residual change: 7.80911e-08 → 1.17095e-10\n", " level = 1, nx = 4, residual change: 8.06977e-11 → 5.69861e-16\n", " bottom solve\n", " level = 1, nx = 4, residual change: 1.0355e-15 → 1.0355e-15\n", " level = 2, nx = 8, residual change: 1.28541e-10 → 1.28541e-10\n", " level = 3, nx = 16, residual change: 1.45795e-07 → 1.45795e-07\n", " level = 4, nx = 32, residual change: 7.4452e-07 → 7.4452e-07\n", " level = 5, nx = 64, residual change: 1.12439e-06 → 1.12439e-06\n", " level = 6, nx = 128, residual change: 1.25658e-06 → 1.25658e-06\n", "cycle 6: relative err = 8.544249588823554e-07, residual err = 2.88200772432267e-10\n", "\n", "<<< beginning V-cycle (cycle 7) >>>\n", "\n", " level = 6, nx = 128, residual change: 2.34558e-08 → 2.27531e-08\n", " level = 5, nx = 64, residual change: 1.6088e-08 → 2.03781e-08\n", " level = 4, nx = 32, residual change: 1.44068e-08 → 1.37252e-08\n", " level = 3, nx = 16, residual change: 9.69812e-09 → 2.83563e-09\n", " level = 2, nx = 8, residual change: 1.99971e-09 → 2.99732e-12\n", " level = 1, nx = 4, residual change: 2.06576e-12 → 1.74907e-17\n", " bottom solve\n", " level = 1, nx = 4, residual change: 3.17826e-17 → 3.17826e-17\n", " level = 2, nx = 8, residual change: 3.29051e-12 → 3.29051e-12\n", " level = 3, nx = 16, residual change: 3.73325e-09 → 3.73325e-09\n", " level = 4, nx = 32, residual change: 1.90594e-08 → 1.90594e-08\n", " level = 5, nx = 64, residual change: 2.89959e-08 → 2.89959e-08\n", " level = 6, nx = 128, residual change: 3.26638e-08 → 3.26638e-08\n", "cycle 7: relative err = 2.210681933627904e-08, residual err = 7.534885150074738e-12\n", "\n" ] } ], "source": [ "mg.solve(rtol=1.e-11)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Visualizing the solution" ] }, { "cell_type": "code", "execution_count": 14, "metadata": { "tags": [] }, "outputs": [], "source": [ "v = mg.get_solution()" ] }, { "cell_type": "code", "execution_count": 15, "metadata": { "tags": [ "nbval-ignore-output" ] }, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAgUAAAGiCAYAAAB3bbXGAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjYuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/P9b71AAAACXBIWXMAAA9hAAAPYQGoP6dpAACNWklEQVR4nO2de5xcVZXvf+dUv1DpVsikEySE6EVeuaOQDKSTG8YHBiKCOirx4g3ihGg+UUPIoBKByUNnMugI4ZUI3mBkiJCPQgQ/xkDrSAATUGKCCgxyeZjAdBMTJQEl/ahz7h9dp2qtPmv12aequrse68unPpzs2mefffbZVb1r/dZa2wvDMIRhGIZhGHWPP9odMAzDMAyjMrBFgWEYhmEYAGxRYBiGYRhGDlsUGIZhGIYBwBYFhmEYhmHksEWBYRiGYRgAbFFgGIZhGEYOWxQYhmEYhgHAFgWGYRiGYeSwRYFhGIZhGACKWBQ8+OCDOPfcc3HUUUfB8zz88Ic/TDxn69atmDJlClpaWvC2t70N3/rWt4rpq2EYhmFUHMP1d/Guu+7CSSedhObmZpx00knYtGnTMPSek3pR8Je//AXvfOc7ceONNzrVf/755/GBD3wAM2fOxM6dO/GVr3wFixYtwl133ZW6s4ZhGIZRaQzH38Xt27djzpw5mDt3Lh5//HHMnTsX559/Ph599NHhug0AgFfKhkie52HTpk348Ic/rNb58pe/jHvvvRdPPfVUvmzBggV4/PHHsX379mIvbRiGYRgVR7n+Ls6ZMwcHDx7ET37yk3yds88+G295y1twxx13DFv/G4at5Rzbt2/HrFmzWNlZZ52FdevWoa+vD42NjbFzenp60NPTk/93EAT405/+hCOPPBKe5w13lw3DMIwyEoYhXn31VRx11FHw/eFzZTt06BB6e3vL0lYYhrG/N83NzWhubi65bZe/i9u3b8ell14aq7N69eqSrz8Uw74o6O7uRnt7Oytrb29Hf38/9u3bh/Hjx8fOWbVqFVasWDHcXTMMwzBGkD179uDoo48elrYPHTqESRPfhO692bK096Y3vQmvvfYaK1u2bBmWL19ectsufxe1Ot3d3SVffyiGfVEAILbaihQL7Vf/0qVLsWTJkvy/Dxw4gGOOOQZHL7sSfktLefpUtGiiUO72SiGsIGtK2Qe6BMo8LGUf5gpvr4KeZGXNqwr6vJW9J2Ua5uDQIby4/Gs4/PDDy9OgQG9vL7r3ZvH8joloPbw0a8TBVwNMmvIH7NmzB62trfnyclgJIlz+Lkp1httaPuyLgnHjxsVWNnv37kVDQwOOPPJI8RzNROO3tMQWBep3Q9JkVj7IxbeX8D6K+MAW+4GsoO/Lor+lUp7ndMtJbSrvq9/5SX+YUrdXpnby5yn9S3v96G2HZxKm/WM9QvNjWBmhz6nnsvgo9ntKLR/6mmm/L0dC/m093C95UZBvq7WVLQrKhcvfRa3OYOtBuRn2PAUdHR3o7OxkZffffz+mTp0q+hMYhmEYRrFkw6Asr+HE5e+iVmf69OnD2rfUi4LXXnsNu3btwq5duwAMhFbs2rULu3fvBjBg+r/wwgvz9RcsWIA//OEPWLJkCZ566inceuutWLduHS677LLy3IFhGIZh5AgQluWVhuH4u3jJJZfg/vvvx9VXX43/+q//wtVXX42f/vSnWLx4ccljNBSp5YPHHnsM73nPe/L/jrT/T33qU1i/fj26urryAwEAkyZNwubNm3HppZfipptuwlFHHYXrr78eH/3oR4vutJPJKmf2cqtL2tYuqtRP1SeHcr0dd7Nb1cqtSsddzO2eUp7Yjqcca7D68gmF68gTi95mWimB1Y8acmmDHKsm/rRjMbgfQ7XnVJ5i4larfKA9cIfviRAuEyehbfV7x5PraBHrYfQ/on/Tql687kgRIECpv/PTtjAcfxenT5+OO++8E1deeSWuuuoqvP3tb8fGjRtx+umnl3h3Q1NSnoKR4uDBg2hra8Mxq74Gv6XFFgUJ1OOiwKW82EUB65VL/ajcS55YodaeLQoSruleddgZoUWBUztJbaReFCS0GSqLgkF1g0OHsPvyK3HgwIFh0eiBwt+J/3766LI4Gh51/IvD2t9KZUSiDwzDMAxjJMiGIbIl/tYt9fxqpvoXBYJ1ACCrVheLgFZHWzFLdRxXzOIxQb1mwnlqG6OM0w8aL3Yw6Jdt0nlw+mXt5QbG5dc5P09pWzGXRoehZkNNG1khWQdoOy4WARcrQFJ9B8sHtL7S4qQ6DuPjVdAkD4uNCiDnqX+DtF/wXFeI11e/uxTJK5TbUy0IQqEqJYwwxfgESG3UK7ZLomEYhmEYAGrBUmAYhmEYOQKEyJqloGiqc1GQJBkAiJxHNZnAC+RyTT7wgqEdcVJJDYPrKCTVLzrR0lBI55bi1FWk018omONd6zNTPqvjxcpYG76DuZ0eK3a2yBLrkffTSgmJkgG5vpNM4FJHvf+BctXsL0g0sTboqUo74pRQnRiLn+SSHFSShKya5OPlXAFwkQ+U9pTvwLwMoNV1+N4T2xtUJy8rMCf9ypASTD4oDZMPDMMwDMMAUK2WAoL2azq/Tg2Uukq5ZhHQLAv58pSOg8VaFsoV+khJs5JP5zjoVl4I4VOuk9IZUOtjVCekS2H2A96L1R2oH0pVtFPzFoJQmo9D9E8lyZFQW9prv/yVcs+XJ1++XHO+VH7tc6sBqaNYFqJyzTqgl4vFTudKuDgOquH7CVYD9sNbsxSwcnKs5Cmgv+aj2wwD+QPBuse+D2n9UC6nVYLIeqR9qOXikcCiD0qj6hcFhmEYhhERAGVIXlS/mHxgGIZhGAaAWrMUJJnbqWRAzGIeK5fbU8vDhPcTztP6rdVxOY/iuZjB0ljKXEKyNZNikgxQgnzAHQYTykP5fb2c/IM5D1ITOykvt7Om1o7UJpUA2DE5LUkmAL+fqI6vyAE+PU+VBtzrMKVDaYPiIg34KcY/cPnIODgXBsz0H39fkw9c6gSKPBBdk81HekOKzsZuR/uZLEkJ2udnFMmWIfqg1POrmdpaFBiGYRh1TTYceJXaRr1iiwLDMAyjZjCfgtKoqkWBF+aslkoMLrOsRk81rWSglPMIhXh9ra4e5UDqaPIBMftJaZtHKhIhfcRBcqy6ZPpnj9UXTJWDjqlMIHn/x67jC+9TM3kol+u3L3uG50szpD3Fu5zKO6Fgsh/cATEnQVrJgMkE8rEvyApcJtDkA7kNTQbICOe6SAa+JiWIpcnnBi4RB0q5dq4kJdCyLPk+0mSHQJESAuXc6DuORh8EnvxZCgOxWM+pQYmuw54Vvc4gyaJc8pkx7FTVosAwDMMwhiKAh2yJq5CgjlcxtigwDMMwaoYgdHMWTWqjXqn+RUGCN75msmfHWbkNtb507FSX2qeVOknSg0MEw7CkP04ihUyglXNTfyjWZWZ9RW5gZnhaP0rsos0ZLUDXZdyk3C/lGm81EiOMlymSgZcJSDmpw2SCQh0qCfiCfOCTuhlPLqcKUIaUM/lAkBg0+UCTIyj+ME5y7dejk7k/d0zfzyp1s4FPziPXIeVZMri0PDLhM8mASg1kgrDPlSd/EJiUwKJzormnyHxG1VL9iwLDMAzDyJEtg3xQ6vnVjC0KDMMwjJrBFgWlUZWLAnUjNkk+UCIOtN3CNFkhqVyTCcoW5ZCXQxQJYjSkhDJIBgPlXrxMiSDQpAG2PwErF2QF2gaJENDGhF1T8dimZtS8KVbJGQP5EXKSIg7osRJBACUZETPJZ2TJIMPKB441mYBGEGQUWaFBkQ9oeWT6d5EM2LEyij77ABVHEMqaEpUSNPlAkhLoef3E7E/lA1qeZbICuWcWXVAoj2QFJjVkiRxBxiQg5VTrYXspFGrIETTKPFX3RDAqnqpcFBiGYRiGRBB6TiGmSW3UK7YoMAzDMGoGkw9Ko6YWBWJqdjUxkHJMpAGfRiUEQ5f7JC+mk0zgEKEgSgmqZCCf5yQlUKQ6Dp8PN8lAybWe65guGciRBWqUgSIrBBmPXm6gLN6NWLmGtg2tNPeGxSdeyj+v7WWgSAYZTT5g8sDAcUMmXgboMgGt0+DJ5VQGiOpQs39a+UCLSqBIsoImE/A6xGSfUj6I2u8n12n0qTRQKM/4vljeL8gEA3XIc85JAvQriqKZ9amUoO7jYZEGdUFNLQoMwzCM+iYLH9kSNwDWFlX1gC0KDMMwjJohLINPgbYLZj1QnYuChIRFQCHqQDXTO0gJmmTgUakgK7ShRS2w9hS5gfUrXkeVBpSIBwYzcZdu0A4dzImquV+QGLSoAXrTtE4kBwx5HRJd4OfuOaRl1AxMu6cmDCLHgVYnlxfeV/SDlMmQNHkgii7Q9i/gx0MnJgIGmfszhUncmJMN6PuNDjJBI2mjgUkCpL5wLpcJSHuKrECh9TNOetnQZAUJYDD0D1AfqROw44E6VD6gUQYBkQxonb5sYbJqskKfIKv0qdv80YgH8rkiGyF4fH9wchzGD7UPzSj+TTWfgtIozcZiGIZhGEbNUJ2WAsMwDMMQyIY+sg6Oo0O3UabOVCHVvyhIkhLKFH0gSQb02FcTHclmfV5fS3YUNz+za1MJQNs/gaJFKBQJT6ojm9v0/QyI6T+ykGrJizRZQdl2OKSyAh0imqioUKPQPUVKULduVqIPkuZeapKkDBZxUHibmeGVxERJkgFQMPEzOYBJCXJ5k98vljOJgZwbRRE0kDIW2UAmvyYfZJAsMaSB7VWgGFZpnX6iTWVZQqKBchq10EcmpCYr0HvIMMmgcC6NLOnTbiSBMJRDDgIa+RMIn3dtbo4iATwEJRrBg2HdKKayMfnAMAzDMAwAtWApMAzDMIwc5mhYGjW1KJA88522JXbZn0CQDICCDMDLFKmB5i7X5AOlPOqLltwIWrkWcVAO65i0f0GsnNSnpkgaUZDNlWuSQUYuZ8mBlNztNEIhOlUNzqC3oEQZsDrKnghRH70yD/dAo/Q4Hn2gJS9Stz2mMoAgGQyUD0zKJj8bKxsol2UCVt+XIxF4+cBxoyITMFO6IhNkaIRCGUadmvuzQjQBwGWFIOwX6/T5A5OYygt0HFjEAZEG6L31e7Ks4EPUxURoqB0dHV/5/DDJkX4mo++bCtzjoDw+BSYfGIZhGIZR51SXpSAkrxz8V7EXL0+Zj8All4AvOhrKv/adrAb9isMgsyCEwj1Q64Dchp7+uAwrYTVtMfkHS11M8wrQ4yh/gHweveeAzFif/urJEOsAqHUgXs5+4NNLkmcV0L7SMXexIOSaCdl8FGK8AW7tcHHaYjvRxetqzoXM0ZA67ylOh9wSkBXKCr+Im5RymlegOUOsCeSDQK0CeYdG5X3NCpBRciBkymApoGZkmneAWQ0Ua0IftQrkJlF/QC0j1NGQWAdY2ubC/fR6xX1d01FoYE645N6oEzCdq+Q7Jksmf3TInI3pRBxkLRzJ390DjoYlbohk8oFhGIZhVD9BGdIcW/SBYRiGYRh1T01ZCjQHO6lMT3Os5AxQ8hBEJmcXyYDKBDzvgXIs9Iu/rzgEMQdFh9wEaaQEp3wEskzgafJBzqRJzfQ0RTCVFVjmYOJflSQZ0HIqDdDnw6yfbB4ocog6nrnzRsrRUNslkUkJgXispS6WHQ2TJYNmoq01s/JCBD1rW5AKWBlpg8oB3BmRSgayG2manAVa3nzuUCjLBFmWh6AhVoc6ETaQCU9lBSYfBMV9RdN5Hyo7N9KdFqkUF9D5LslVILKBORrmWbNmDb7xjW+gq6sLJ598MlavXo2ZM2eKdS+66CJ897vfjZWfdNJJeOKJJwAA69evx6c//elYnddffx0tLS2p++eKWQoMwzCMmiGAX5ZXGjZu3IjFixfjiiuuwM6dOzFz5kzMnj0bu3fvFutfd9116Orqyr/27NmDI444Ah//+MdZvdbWVlavq6trWBcEQI1ZCgzDMIz6Jht6LJtksW2k4ZprrsG8efNw8cUXAwBWr16N++67D2vXrsWqVati9dva2tDW1pb/9w9/+EP8+c9/jlkGPM/DuHHjiriD4qn+RUEKE65LxIFTHUEeUCWDPlka8F3kA0kqYBEJQfx9YFBuX9JvJf64aFhuAiVyQJEMkClUiszzkqQwuKvs1sixz255aCnBYyZRJcpAOVbrJA1tKcOtmb6FNMcuEQcZtVzOHxDJA5pkcFimIA1wyYBGIsiygiQVaNEHLBKBTGwelTCMuyRCjj6gsgKVEvpIf/PygSApAEAP+SrWdon0U0QfBIpkQFMl03kQkA9ToMyhQEpzTNEiEaqYgwcPsn83NzejubmZlfX29mLHjh24/PLLWfmsWbOwbds2p+usW7cOZ555JiZOnMjKX3vtNUycOBHZbBbvete78NWvfhWnnHJKEXfijskHhmEYRs2QzUUflPoCgAkTJuR/1be1tYm/+vft24dsNov29nZW3t7eju7u7sT+dnV14Sc/+UneyhBxwgknYP369bj33ntxxx13oKWlBTNmzMAzzzxTwugkU/2WAsMwDMPIEYQ+cwQtro0Bq8eePXvQ2tqaLx9sJaB4g6wnYRjGyiTWr1+PN7/5zfjwhz/MyqdNm4Zp06bl/z1jxgyceuqpuOGGG3D99de73EZRVOWiQB1mKbpAiUjgsoLm8e9yLEQfKFEGTDLop+WkQVVKCGJ9ZXUDauOuHPmAl1PJgEQXRLJBhmU6KtRlplCS6pVchqo+VKWgUkIkG9BoAjfpiEZwyLKGNLe0IdY2V6SESZIByP3Qt4X3AT3iQNu9sEGo4xZlIB+3kOgDLiXEy+n7TZ4sNfA0x3LCogx7oMXBZQIlkZEiH/SGUvRB4R56gkaxr4dIuYb2Ry+SCmhkQ5Z8N2TJc80qkQgekRiU7OX5uaVJ73T+elUsJbS2trJFgcSYMWOQyWRiVoG9e/fGrAeDCcMQt956K+bOnYumpqYh6/q+j7/7u78bdkuByQeGYRhGzVBO+cCFpqYmTJkyBZ2dnay8s7MT06dPH/LcrVu34v/9v/+HefPmJV4nDEPs2rUL48ePd+5bMVSlpcAwDMMwJAKkjx6Q2kjDkiVLMHfuXEydOhUdHR245ZZbsHv3bixYsAAAsHTpUrz00ku47bbb2Hnr1q3D6aefjsmTJ8faXLFiBaZNm4bjjjsOBw8exPXXX49du3bhpptuKva2nKj+RYGWsEiKPtAiEeh5Dvsg+CxRUVS3eMnA64/LBIPbQX8kHyhSA5MM6GYBimRQZvmAaWe+LBmwBxDE63APfmqmp5JB4d6oN7i2nwHd2c2PdmNUkv2w6AdlfmgRB2L0gTY306L0Nyr3tIRF7LhwmhpxoO1tEO19QPcySCkZUJmghUkC8fo84kBLXkTKFfnAT/3VHifwZPmA73FQ+Bql/WoM49EHVBqg3v+HwnRfxYGyY2OQe1asjETy0N0YaSQC+8gq88ljcy+afAkRCXXCnDlzsH//fqxcuRJdXV2YPHkyNm/enI8m6OrqiuUsOHDgAO666y5cd911YpuvvPIKPvOZz6C7uxttbW045ZRT8OCDD+K0004b1nup/kWBYRiGYeQoJvmQ1EZaFi5ciIULF4rvrV+/PlbW1taGv/71r2p71157La699trU/SgVWxQYhmEYNUN50hzXr7td7S8KEuQFoPiIA3rssq+Bk2TQLx9HkgCTFJiUIEcfeEo5owx7H7BtlGkUARk4j0UXkONo22Mmb2gfSkVK8GSPaRaJEFk82bNMjkTwtCEslzxQKg7ygVPEgRJ9ECUyotsfS1EDgJtkwMt788dN0t4H5Lwmut8BefZNyj4IFG1PBAnNyYya7HtRiDKg1+8l0QcsEiJ3fRoFckjb7lyZ+lwmKIwhjbCJjqlM0BDIz5hGIvTTvTHYZ0nRyIyapfYXBYZhGEbdEMBjC6Vi26hXbFFgGIZh1AwmH5RG1S8KEq1bLvKBQ1IjVW6ILIfq3ghKsiEXyYCWRyZAKgcQs6CnlKubBaSREhwkA7DkJzRso/DhCmnbxPPay30A6efQZZ0eKmZOtsUr7VcUKUI/71oSK5Y8Sa6fJB+o0QnxqkOj5IjKO4Cr7yd7kWeY2XhoWaGBmMmpvKDtZeAiGVC5oSAfkPaQjb0/uA6LOFAiEShSeVaZcSxhEbHr0+vT6ANf2ZMhk6vj07AWTTKgVchkpeNMJY4GEuXQkNuamT4/NiY+/U5Jnh/63IqXlRgFWDbS5hnQ2qhX6vfODcMwDMNgVL2lwDAMwzAigtBjTpnFtlGvVOeiQE0YXziMLGBO0oDDdrhSxAE9VvdPcNkWOUkyAID+nInQRTLQIg4UWSFMEX3gqfsd+HJ5htwnvWaG2tNzY0hmY1iwlHLDLkuYJD9cj0kZ8eesRhlQySDlPghSV9RhdRlu7TtJuJAuDZBjUt7AJAO5nG6d3OhLUQHyMY9KIImJEiSDgTp9sfPYPgiQr8nlg1AuTyHaUCczLh8EYrkqX4QJUoYqGRTabiZzmerc9P77icTQn3tWLPrAL0RE+IH8vLWttfXkRcIXLGUU/6YGZZAPSs1zUM3U750bhmEYhsGoTkuBYRiGYQiUZ+vk+v29XP2LAnXvzkH/H3SsRi04JDWS9kpgcoAmNWRlmYBJApJkAMAT5AMtEgHZwnmhGn3gkNRIIGTyAd1YgNwDc0nOkDqylBCV014wKYGOfT+VD0i5IhnQLa3zkQgpn7HLXElMZMTmacr4A8VEG7XoEnHA9kEg19eSF/E6A4PI5AUqH2jRB4pMIEkG9FiXD+j15b5mVPnAHToNsqqUoMgHNPoACfIBvSb5LNHERFmyPwEbZ/K56iPyQBSJ4JOvdu0Z9ymyi1MkgnQTbJ6Onn6QhadGk6Rpo16p3+WQYRiGYRiM6rcUUFLEiusWBCU3AT0mvz4Ljob0ffnXObcgKLshar/+o2NiPWAWAW1nxCArl7Of35DrRPjKqplZCkiaY/LLhaVcps6I1IIQzUK6vWGWWgToTxTmkVVoTtmlUnI69GjbpBv6PPCUOvKxVLdsCI+C/Zoj5b5iKeAWBM3pMO7U5+Jo6HYs75IYHfNdFIl1gOUACEk5tRTkD9VfPBmhLCuUAUBA2s6S4z4y0j6zFMiWikIbsuUhy3ZjlB0KXY4jp0P6/PyA5FGgacfVHTXl+eRJk7wCf1CbfFAatbUoMAzDMOqaLEo3/2uLxHqgfpdDhmEYhmEwatZSkGTOdXI0dJAY8uG6Ls6FSspjUSaInRulOVYkg2w2XnegEmmb5AnQ8hdI0GUzkQCYOZGa2KkEQXZGVN2QcvIAlQmYQ6MiKzDHReIEpzod5vqoSQMu80CVoKDUGUYkc66Lw5hmKmYx9kIdF+fCpgRpAOCpiyUzuCYZNNFyRSag0kDG88Q6Eo3kmDkahnQc6HWolCBLjnSSZ3P3TKUBmncg69EUxoUTm4i+1UeONafDnpyHbsZBIko7VyiilFAhmHxQGjW7KDAMwzDqD9sQqTSKuvM1a9Zg0qRJaGlpwZQpU/DQQw8NWX/Dhg145zvfiTe84Q0YP348Pv3pT2P//v1FddgwDMMwNMLc1smlvMJK9KAcIVJbCjZu3IjFixdjzZo1mDFjBm6++WbMnj0bTz75JI455phY/YcffhgXXnghrr32Wpx77rl46aWXsGDBAlx88cXYtGlTWW5iSBzMvbpkIEsCzL4YebQzk7RDxAGtn1UkA0FWcJIM1DwFDumPk6Cpivl2hPlDF9Mis7JG1ycyAY0QoNdhEQxsfEh9Ju8IkSCBLC/weSBHHKg3MYrWVDUVLUGXDIgXP4tnz8bqZJQdCDVPeFpHy1NAUxdHeQhcJANq7tdkggyNENB2+hQIVMmgAJUV+N8P+hmjpV7u/4X7ZbkJSM/prosuY5sR5IEGdl6yfOAiE1SyZGCUj9SWgmuuuQbz5s3DxRdfjBNPPBGrV6/GhAkTsHbtWrH+I488gmOPPRaLFi3CpEmT8L/+1//CZz/7WTz22GMld94wDMMwKJF8UOqrXkl15729vdixYwdmzZrFymfNmoVt27aJ50yfPh0vvvgiNm/ejDAM8fLLL+MHP/gBzjnnHPU6PT09OHjwIHsZhmEYRhLRLomlvuqVVPLBvn37kM1m0d7ezsrb29vR3d0tnjN9+nRs2LABc+bMwaFDh9Df34/zzjsPN9xwg3qdVatWYcWKFfE3Qi/3KhSl9QyX3nfyOtd2T4yOlZ30oMgKTqZ8SRKgbbtIBlRiIG2r6Y+ToKZ5Jh8QWSFTMLR6WSXilyY+iupoOzCGimQQyten4x8Kz83lubrMDxXpOlob9B5YZyAeaztDSriYil2Oo3S9viIZZJTUvjTiwFfKGwUzN0tM5CAZNCoyAZMPUvz+CZQUyoEmGSjlNNlRJIlkiXxAowx6SQQBHZ8+Us7Hlsg7EJ5VmZ69irBDZ6jMWXjhyIXkGCVTlI3EG6TPhWEYK4t48sknsWjRIvzzP/8zduzYgS1btuD555/HggUL1PaXLl2KAwcO5F979uwpppuGYRhGnZHNbZ1c6qteSWUpGDNmDDKZTMwqsHfv3pj1IGLVqlWYMWMGvvjFLwIA/vZv/xZvfOMbMXPmTHzta1/D+PHjY+c0Nzejubk5TdcMwzAMoyzmf5MPHGlqasKUKVPQ2dmJj3zkI/nyzs5OfOhDHxLP+etf/4qGBn6ZTM68HKbYna9omLU1+XrFmoqhmeNZpALEOiziQJEb8mNF9zIIFTlAkwy0pEZpoJcnEoD2EWJSAkukROQOX0oqJI9PyKIPoBxr8o3SSeF9NwsqbXv0vkS0XRI1tFz9an0h+oAnOkoub9I854W+NLKyAi6SQSOJEcikiDigZEgbPnnGfVTqYOEzhUMalUDvI5IStHun49NHr59ynJNM/2mfvbZLolG7pA5JXLJkCebOnYupU6eio6MDt9xyC3bv3p2XA5YuXYqXXnoJt912GwDg3HPPxfz587F27VqcddZZ6OrqwuLFi3HaaafhqKOOKu/dGIZhGHVNAJ+FeBbbRr2SelEwZ84c7N+/HytXrkRXVxcmT56MzZs3Y+LEiQCArq4u7N69O1//oosuwquvvoobb7wR//RP/4Q3v/nNeO9734urr766fHdhGIZhGACyoYdsiZa7Us+vZopKc7xw4UIsXLhQfG/9+vWxsi984Qv4whe+UMylhg8HmUD3GCcmtdyx5hTuaYl01KgEJRIhn3jHpW6yZBBqEQeSrODJq2a630BIt5NmlZSIAk+4D3YPdFtmRVZQx1O+TL6+y94HlLSSUoXAtlFWOq7ufYC4eVp732fRB4o0oJi+WZRDXqYo9I8mDKKz0EUy0CIOfEHsCtQvAfr5ofsQEO//kN4zaVPYK4Gb75UxDJUxVKI8JClBe1+TF2i/6vfPoQHY3geGYRhGDWGOhqVhiwLDMAyjZgjLsEtiWMcZDW1RMBhVSqBma6G+FtngIhk41Q/4/zEo4kCKVBhcP5DL1evnkeuGxMTPpARmkpf7xSIR8vdGPojDMp7RtQtFPCLFYb+DGoCakynSdslAwbScUQYlwyQA2STO6isSQ8HETtuWkxGxY0UykGSCgfrCF74ajUPqknujvyQDcp1AiUTIJEQf8CvKkQUZxfQvRXBIskzsOso8qHay8FhiqGLbqFfqdzlkGIZhGAbDLAWGYRhGzRCEpfsEpMn8XmvU/KIgdTKicpA2KZNSX5MH8igRB6lndFJ/WW53hw8bvT7L50/dseNSAZMX0vRvMCORFAtInexouPCY2X9krumUH5+QUaQE7o0vvU+vWdxeBoAiGSjvZx0SezGZgo6/Mvf8/P81CSCdKT/t+BeLz4KHquMvZVAGn4JSz69m6vfODcMwDMNg1LylwDAMw6gfAnjM+bPYNuqV2rUUhBhSLoh28xzSIhYisZ2yE4aFl8SAYFac6BUG5JVwHZf+0fbSknQfxfavFByed6p5M8x43tBqju+F+ZdeJ8i/RooMwvxLfr/wSosPL//KeH7+lap/5DzaXlqS7iNpHIYDl+ftMm+S5t5oEmU0LPWVljVr1mDSpEloaWnBlClT8NBDD6l1H3jgAXieF3v913/9F6t311134aSTTkJzczNOOukkbNq0KXW/0lK7iwLDMAzDGAE2btyIxYsX44orrsDOnTsxc+ZMzJ49m6X8l3j66afR1dWVfx133HH597Zv3445c+Zg7ty5ePzxxzF37lycf/75ePTRR4f1XmxRYBiGYdQMkaNhqS8AOHjwIHv19PSI17zmmmswb948XHzxxTjxxBOxevVqTJgwAWvXrh2yr2PHjsW4cePyrwzZVXb16tV4//vfj6VLl+KEE07A0qVL8b73vQ+rV68u21hJ2KLAMAzDqBkCePlUx0W/cnLRhAkT0NbWln+tWrUqdr3e3l7s2LEDs2bNYuWzZs3Ctm3bhuzrKaecgvHjx+N973sffv7zn7P3tm/fHmvzrLPOSmyzVMzR0DAMwzAE9uzZg9bW1vy/m5ubY3X27duHbDaL9vZ2Vt7e3o7u7m6x3fHjx+OWW27BlClT0NPTg//4j//A+973PjzwwAM444wzAADd3d2p2iwXtbsoSPATcfIjGQ1HmiTvHRo4nNWryW1TwxBNM5zC0YntdFiCoSkpoH40vJgcLllJ8ybpsbkkcBmNeOykFLJ0WjembJvtdhjSFMHu90nzFKi7J7q0k/j+yM9xl+ftMm9G0v83LWEZog/C3Pmtra1sUTAU3qDvrDAMY2URxx9/PI4//vj8vzs6OrBnzx78+7//e35RkLbNcmHygWEYhlEzlCwdpNxlccyYMchkMrFf8Hv37o390h+KadOm4Zlnnsn/e9y4cSW3WQy2KDAMwzBqhnI6GrrQ1NSEKVOmoLOzk5V3dnZi+vTpzu3s3LkT48ePz/+7o6Mj1ub999+fqs1iqF35IAdb8I2UtS6teUepT81EoVTHV3YVpMcuEgNtW7ILpr0fKg3QPvrKBy3XvmoWK9N4lh06bKMYs013pRypnO1pc8tnld8fbIdB8X16zcLNBXT3QIeMBtkEKcEltTHvF90xkfRLrR/9Xx43bXzU64/QhGMbq47mJK9wlixZgrlz52Lq1Kno6OjALbfcgt27d2PBggUAgKVLl+Kll17CbbfdBmAgsuDYY4/FySefjN7eXtx+++246667cNddd+XbvOSSS3DGGWfg6quvxoc+9CHcc889+OlPf4qHH354WO+l5hcFhmEYRv2Q1vyvtZGGOXPmYP/+/Vi5ciW6urowefJkbN68GRMnTgQAdHV1sZwFvb29uOyyy/DSSy/hsMMOw8knn4wf//jH+MAHPpCvM336dNx555248sorcdVVV+Htb387Nm7ciNNPP72ke0vCFgWGYRhGzTBaaY4XLlyIhQsXiu+tX7+e/ftLX/oSvvSlLyW2+bGPfQwf+9jHUvelFGxRMBhlLoTMlE+3x4v+72D61o6d6vv8/+C7ltH+cdmB1PcLxs0w0MyVggFU8dz2fKF/g64PpV+szfy9lTI+LvUH/kd/BIiyDKlbi2h6aZbJENSsP3CsecvT87Kk7UCZN7Qd6Zh6/GeJaZ6pUnR3RbaLJ52/yhwXpAItyoDKBLQvWSjHipSQFcZQG8+A9JuOp5Z6N4v4swrYM1GuU8c7ARo6tigwDMMwaobRkA9qCVsUGIZhGDWDLQpKo34XBQ5RCWrkghAVwE3Scl3dlE7toprnfhB/n11UjkSgXQlJJAKVEjgp9qajkkFGiTJgx9p9ekPXVcYqVI+hHHtCe3JdxmhEsJQBahDXNFImE1BTtWDC1t5n5m7FPE7PDdi5Qaw+NcfTaUBnJvX47/PIxA5JLWUnwDRxBlQO6COhPFrEAQ32ydIgoDBu1mdjEjqMIZTxD+PPSntf+2NH50cF5yUyRoD6XRQYhmEYNYdZCkrDFgWGYRhGzWCLgtKo/UWB6mkuG8lSJzuSog+4m7RyTMx1xPTuBcQYKcgNoU/MoySzCItEyMgSALt7eqNpMt748T4NlJObI9dXIw58oY4mE2Q0CQLKsSY9YGhSJiMK07Q9jIQpE8wEgrf6kPVzbWY1c7dDeS8x/jd6/YU6Ql/6SFlG8eyn3W4kxVRKoF/sforkQFzSIDIBkxKSIw76hHvTpIFeInuUMs5Jf8zSPns6nyp5vwOjfNT+osAwDMOoG0IUl2dgcBv1ii0KDMMwjJrB5IPSqK5FgRfmXsme5kWbil2c+6U69H3NxK6Z26mXNDPDU7twzrxIbXgZ2STqZQsmVCYlUDlCTfiSAOu3fD/snun1NRkgqqNFMKjjpsgNisQQPSuX5+oyP1Sk66iRLKFciWkCtNhd4+Ce7sUfR57s9NdXHzV3K17xfWHh66WJmPVpOZUHomOamIia6QfpX2K5z0z5RD7w3H/7BUqSIhZlwKSEAn1CxMFAuZ/7f4aUFcYhUMZNH9tCOwGEZ1WmZ6+Sq8PnI5Rjb0Q3B7FFQWlYSivDMAzDMABUm6XAMAzDMIbALAWlUfuLAhfTr2rmpSZpUu4T21jOVK3uPaAlIyLmTJr4J6SJh6iJNHdNj7ThYhD1QjmaIWTm1xQfAE0mUPZBYJIAjUqQEhzRe8vIbbAkRbS+EnEQSn1UpB5NlnKZN6MbfeCJxxQXU3EfmeT91DwdDkQLZJUEO9wkLkcZUO96nxjifSJdReU+ZOmEGfDJbVJTPktwpMgKSWgyQVapQyWDXjJGfWSi9eaOaRQGHRMWfeAwtiyKQYg+6GfnydEJLpKBy9yqNGxRUBomHxiGYRiGAaAeLAWGYRhG3RCGXslWjWqxigwHNbsoEJ+pFlmQ0utc2kaZmqGZSZolJiJ2RmrbpNoEkSbCjCAlELMgc8Ymx8x8H9AND2h7xABKTfVSIiNf+YCwZETUVpsRyz0lwVF0/dBl/4SMJx7Tc7VnUYg+oJIBYu8PVe4iK4zU94n0xRUo5t7QwWwcKNv0RnWYKTsofHXQhEG9xHO+kWy2waIMiGSQIfMwkzPKZ9SIg8IhjQpopBECiqxAkVJ7ZYUyYLCUUDjuY+Z7Yu5nEQLxaAE2hpAjEehxLytXxp9JPV6sT9reB2nnCqWS/2gG8ErOU1Dq+dWMyQeGYRiGYQCoYUuBYRiGUX+Yo2Fp1NaiQHiObglkaB1qWiamfFqH2B/DbC76gJn9icmc2DPVSAQmE8hbIOdlBfrEaMQD3TOBJi+i16TSAJUPWEIkDI26l4EcfcFkAl+OKJCiD1jEAZUJtIgDtpeEHDUS1WfPz2keaHWUdqTzyoVgEWemX1LuYirWTMv9ATFV5+SgBjJnNK947TgDKhkMnbyIkqXbFZMBbVT2J9DaoyZRmmxIgm+FLMsEPJFTPMpgoDwuD2gyQdrxTDqmz0+TiFwkJZaDSJrkFZgP2HwKSsPkA8MwDMMwANSCpYClLxV+3ZXgPMZ/Zcp1omUV+3VKf4RrTod8CS5exyNPJ1+DWB7YL3LqYaWlMw5Yx8jl3Zf7bmmb5XtW8xBE5Q2FXzeh5mjYILehOhdSx8SoespnnNbaJFoIUqTZjaE5e0X/V9Igq78ESQf7WW4CX6kT/fokqYq9ZKfDQ0Fj/tj3ST6CcOi02nwnQdqnLCmXnRh91VKQJk8Bvb48hlougV5ov+AHxuhQWBgTLU8BHTcX58L+IJ5fIu0zdnE6VDJvk8LK+HVt8kFpVP+iwDAMwzBymHxQGrYoMAzDMGqGsAyWAlsUVBva8xLMvJrplzupyeZ7rb5ocmZlstMhFLMcbZCFZReyxBakhCxpg8oEpJyZ75l8IDsxsuGUpAQtDbKa5liWFdQ8BFGeAk0mYJKB4lCY4FxIj5PyGAxuQ0uLrM4tSbqiuHzfaNZu4UKauTdLJCXqJEfNzY2+XN5HHNUacnkFGohkoDkd9hDTd4Z8Dg6xeUjuR/BqUlP4EsdWWofLBIFYTpEdGuWHwuUDLc0zzSswdB4CTSY4FDaJ7fVoUoJ2nHtu9PnR5xoo8yCrzBs1zbH0BUupQAdEw43qXBQYhmEYhkAI+bdN2jbqFVsUGIZhGDVDAA+eZTQsmqpfFGgmXC92oB87SQyeIjHkrHQedeynpmyasZWZ4pKjQZmUkGtfizhgcfVUjsiSjvnK+jfNstpBSmAygVMkghcvU+UDOeWxJiVAeFZqlIGDNJAm4iAxj4ErLLogfiyVDRwnm4Gz1LRMxl/yWO9nHu8F7/8e8jVCUxgfImZwVTIg5ZE8kKUyAdlpMaC7DpK+NJGIBx9yZENGKZfIKpHaVD6g19f6JeUhYPIBkQx6mJQgRxn0kOgDKeJg4Nhn/x/cb/q8XeaHPrfiZXX987qGqPpFgWEYhmFEWPRBadiiwDAMw6gZgtBjVtli26hXan9R4CIfKAls9OO4pzux7DErPXO6Jo34SS7YADe99+eOqU5BIg7YDoQkSRHz6NdsgcWiJjJSog+YlEDHMEpepEUZUCmBeINrqZDJcdBAryNFHyDVcaqERSOFQ8IiajampueMJhnQSIRc1IGvJAPyibRGPeoZqmRA+usPJCCm5nvNy7+RyArUyz8TypEIxRIofaFREbRfWn8jGYAnJiq8f0iJONASGdE6PVlynSBKNDW0FDS4TlaJUBAjDoyapvYXBYZhGEbdEIZliD6oY/8IWxQYhmEYNYP5FJRG9S8KEsy5emQBqUrPU8zGAUtCFE8C7lOTG9mFj25MyC2ospQQEs9rj5hl832kexxQnYLta0AjFAqH3jDKB0wa8JU6LGGRYNZX5ABaN1ClBC2REa3D/x87T0lYpEpNSREK5ZIXmIe3Fyt32u+AtJFVZIW+bGFgqCQQHftEuvKJXkbLXaBJiJr9eF8aSTRBE0mYRMupKZ1en+2SmLJfElklYRFPZESlBLqfQVw+0CSAPiWpEY040I7ZLo25yd1L3yfPNaskLAqYzOmSvCj6v1BmVDXVvygwDMMwjBxmKSgNWxQYhmEYNYNFH5RGVS4KlM2SxYREnmLC1fc1kO39SZEIYaZwHs2GxSz8pCtMSlAkA7bdbnQdYufzfObGTTqlXHTU5QMlwZEUFeCQmEiTDFjEAZEKCtEHtMzlWE5qlJgMS/lecRl5lugKcbkKKPyaoe9ryYsC5mleqNRP5hiNRMhQ7/Xc/Oz15K8LXymn0C/ZZr8QOUDN8425ciYTKPIBTUZEZQK290Ep21Xn+0fGkEYiOMgH0rG2/XGPQzk/JjKFkNRIjT7QIg7UpEYoHAPxOsoQl/pHuRRGy9FwzZo1+MY3voGuri6cfPLJWL16NWbOnCnWvfvuu7F27Vrs2rULPT09OPnkk7F8+XKcddZZ+Trr16/Hpz/96di5r7/+OlpaWtJ30JHktHqGYRiGYahs3LgRixcvxhVXXIGdO3di5syZmD17Nnbv3i3Wf/DBB/H+978fmzdvxo4dO/Ce97wH5557Lnbu3Mnqtba2oquri72Gc0EAVKmlwDAMwzAkBiwFpfoUpKt/zTXXYN68ebj44osBAKtXr8Z9992HtWvXYtWqVbH6q1evZv/+13/9V9xzzz340Y9+hFNOOSVf7nkexo0bl7r/pVD9iwLt2UvRBzTiwMFsrNZhexsMzB4mGUCREmhQgJL3hx2TE7ycydfLxiUF2g9ad+BCiv2PUoa9D9hnUNmaWJcSEpIKuUgJimTAJQahDXWLZPlYrZMUaVBS9AE9WXieTCaQj7NkbH2lPCskLAJ4JEIS1CufHecSEwHc3E4lgcbcA2JlRGqgkQW0Do8+kCMO0t2D/LBov4NQlgxY8iJBEtDkBWrW71GjD2TJoDcbL6dbJ2sRB9rW2tocSow0qBAdvpyOhgcPHmTlzc3NaG5uZmW9vb3YsWMHLr/8clY+a9YsbNu2zel6QRDg1VdfxRFHHMHKX3vtNUycOBHZbBbvete78NWvfpUtGoYDkw8MwzAMQ2DChAloa2vLv6Rf/fv27UM2m0V7ezsrb29vR3d3t9N1vvnNb+Ivf/kLzj///HzZCSecgPXr1+Pee+/FHXfcgZaWFsyYMQPPPPNMaTeVQPVbCgzDMAwjR4jSUyZE5+/Zswetra358sFWAoo3yIoahmGsTOKOO+7A8uXLcc8992Ds2LH58mnTpmHatGn5f8+YMQOnnnoqbrjhBlx//fWOd5KemloUsK2TE0y43NxNyrVIBOrFTq3zOVO0Khko5VSaoBb2UJAMACCylkqSAsATE4VKxAGzoJY9+kAuZ2OryQe5OlrSIdqGtJfBwLmkjiAZ0HJp6+tYuSZlJMkEpHxYrKlCIqOQuYjLyWaod3lAJgKNROhjybLoHgJkkFLAkieRwWoIC6b/fvJB6PMF+YA8IB5lICcp0hIZFYuasIhKI+qeCHGpgO0pwbY8JhEHWS0x0dCSAVBIVMSu4xBxoEUfQIlEQEL0wWhSTvmgtbWVLQokxowZg0wmE7MK7N27N2Y9GMzGjRsxb948fP/738eZZ545ZF3f9/F3f/d3w24pMPnAMAzDMIqkqakJU6ZMQWdnJyvv7OzE9OnT1fPuuOMOXHTRRfje976Hc845J/E6YRhi165dGD9+fMl9HoqashQYhmEYdU459QNHlixZgrlz52Lq1Kno6OjALbfcgt27d2PBggUAgKVLl+Kll17CbbfdBmBgQXDhhRfiuuuuw7Rp0/JWhsMOOwxtbW0AgBUrVmDatGk47rjjcPDgQVx//fXYtWsXbrrpphJvbmiKWhSkSdIAAD09PVi5ciVuv/12dHd34+ijj8YVV1yBf/zHf0x3YY+8cnDJgJjQc6ZqasnWvMipWV3b78AntjOaOz+q7iQZUJMs2faYRSIULKcs8VBkzuaSAQ1n0GQC0raWiaRY2HPwlHJSX0kClDfVa+Z7IQERwJ9DGikhUCQDl30QXLZULiQvCuNlg+rSYzYP1UiReMQJPY+a7L2Azj054sDLkj04SEOFWIF0MMnAL0xmJh+QyIZ+UieSFaik0OPRPRaIZIBALlckhmLRpAG2LTXkcnqclw+opBBQ+cAXy3sDOUJB29sgqq/td9CfpRJIcsQBi8SQ5qdmph+Ug21E4xLKIB+k1f7mzJmD/fv3Y+XKlejq6sLkyZOxefNmTJw4EQDQ1dXFchbcfPPN6O/vx+c+9zl87nOfy5d/6lOfwvr16wEAr7zyCj7zmc+gu7sbbW1tOOWUU/Dggw/itNNOK+3eEki9KIiSNKxZswYzZszAzTffjNmzZ+PJJ5/EMcccI55z/vnn4+WXX8a6devwP/7H/8DevXvR398v1jUMwzCMYhmtjIYLFy7EwoULxfeiP/QRDzzwQGJ71157La699tr0HSmR1IuCtEkatmzZgq1bt+K5557Lx2Aee+yxpfVaQXI0ZL/m6OqV/uKjE0CzICi+VtHvEp62WLYOeMwkQA7pToa+XB713SMWBvbLn/4qdMhHUIYMsLpzoVZHczqMxlZxuNQcALVf/Foegqg+szzQ55o2T4FmhRppR0P67Gn/6C8+j/4SpCl6C/Rl3SdFAPnXZAPxcqXl9JdwA5nYtNzPfR01EOsBTVXcwHITyH3VLAjFolkEtDrUEkB/iffnJl/A8hjIlgLNIsCdB2ULguhoSK0DqqMh+cwouQnE8gp0NDRKI5WjYZSkYdasWax8qCQN9957L6ZOnYqvf/3reOtb34p3vOMduOyyy/D666+r1+np6cHBgwfZyzAMwzCSiKIPSn3VK6ksBcUkaXjuuefw8MMPo6WlBZs2bcK+ffuwcOFC/OlPf8Ktt94qnrNq1SqsWLEiTdcMwzAMY8CKMcI+BbVEUY6GaZI0BEEAz/OwYcOGvFflNddcg4997GO46aabcNhhh8XOWbp0KZYsWZL/98GDBzFhwgSlMwnHyvuarEBDsvmuhrShuM1Mcxak48LkAE1WoBJHED/2qLmbOZWFsboxaP0y5CnQnAtZHcX0LznmOdV1cAxMykPA5QNFjlDzFxSOi517qVEcR6U0x3QXTSoZ0AmazaYyEMqXVhzqsgGRD8jYNgRyLoEGvzDQkazgBw1i3UYvWRpgeQrKvUtimCwf9DFnxLgDoiYT0LrM6TArSwzZhDwEVDLoV5wO6TzQnAt5zhPhe0qbm0bVkmpRUEyShvHjx+Otb31rfkEAACeeeCLCMMSLL76I4447LnaOlF/aMAzDMJIYLUfDWiHVT4ZikjTMmDED//3f/43XXnstX/b73/8evu/j6KOPLqLLhmEYhqEQlulVp6SWD9Imabjgggvw1a9+FZ/+9KexYsUK7Nu3D1/84hfxj//4j6J04IRDzHc+RbG2s5cScTA4vjZCkxKi5n3quU0vKUgAA/2TZQUtD0G+jtK/1OWqvV8oczB9a3H4WoSCVF9LJ6xHIsjHWhRB3pwtRCQMdZ7eL+04ChUhbaeVElTJgEoC8SJo+SJoLD2ZcFz2k38jRFJBA5VayEUz1CTty5EFNOIgQ+r4ZL7n5QNFJlCPIZdr+IK+pkkDvA65fyX6Qj8eaL9fkBSAQbkEFJmgX93tMJ6HIKumNpbzEYQ0zbESicCiD6Jjl+9Xo6pIvShIm6ThTW96Ezo7O/GFL3wBU6dOxZFHHonzzz8fX/va18p3F4ZhGIYBlCV6wKIPUpImSQMwsAXkYMnBMAzDMIYFs1QUTU3tfSCZpD3FbOtgLWSSgRe32rJylrBIiwSg8oEiK4iSAZCf5FqqYk8x4+mJjJTyNGgKhJOUQMu9eJlLkiAXc7+0G2NCdMKQ5Q7RB6FgWR0WRA9w2XOcqQ5K9EGoeFflL0NN1mRcM2Q30ayDrEA/Hw2kvC93Jc9FMlDkA4okE6RFjThIKR9EY0fP61d2JnSRCbJMBojXD1yiDLLy9bXkRSYP1Ac1tSgwDMMw6huTD0rDFgWGYRhG7VCO6IE6toRU5aJAW8RJUkHoU3O7N/jtXDk51GQFzXs88gBnckByZIGapIhORqmcJSAidUdKMtAog5Sg75NAjrWdCR0iFCDJByl3QNQlCyH6RLO8pow+8JgZWphD1ENcDj5ASCeFr5mQA3JM9iTI1Q/IPfpshz0iH5C2af6pfiofsF0N45IAlQZcpARKOXZG1AiUSZ4kGdByFnGg1GXRAjQCSpEVpD0MnKIMtCRFSREHgChdaRLmyOOh9H0Z69dSUFxqM8MwDMMwao6qtBQYhmEYhojJByVR/YsCzcqTYJJW0RyWE8zgWg6Y0CGyINRkAJbISHg/rUzgMNGlNp3M3QnPYcg2pTwoyt4HLhEKidJDysREaZMajdzeB9HAyVm2QkX3oBIE26qb3WgYq0/PY+ZpIit4xDzNzf2FciYf+HGpQJMPKJp84DLM0rmBwyTXPj7auZJ8wGWCoaWGwcesjnZurjxU9zJIlgZ4Ob0heuwJZagMbFFQEiYfGIZhGIYBoBYsBYZhGIYRYVsnl0RVLQqiZ81MijSfvrCfAfPAZjYh+aFzj23yDyURTN4qqsgBmkygRQ4kefO6nEdx2iI5janM4bMSKttoJ0kwLomOnJIhqREKpZ0XP6Ye/cJ1SHP8Hohpnk9QiG9Qc78UiRBoD4Vch3SAXVKZQ4Efvz7dlpneg88kAFJFkQGS6jDRw0FK0MopvsO8jQhcPjIpJAOAOOsrcgALKnGoEygJhvLX1JJYKeepkkFC9IEacTC4fATN8bZLYmmYfGAYhmEYBoAqsxQYhmEYxpCYo2FJ1NaiQLDoiaZ+DJISmGlXMYdpZl5p+1p6mkuCoRR10kccJNtNHayvaZorPhKhXPJBUkSKg0wgbsM9uO000QflQporTObSBpkcsnktyxRUKog+N1m6vwe1KtO6rGl5MuvygXCe1h4rF4udzpVwSXGrmZeTZAVdJqCVFfkAWh1anvufi0zAvsccytn1pU5VCOZTUBImHxiGYRiGAaDWLAWGYRhGXeOF6ayfWhv1StUvCtRc75G5Xc7Honpg0xzxnmpmJYe+UFcx1XpKuYYqMUjvU0qZ0NK5pVjSkiQDpW5amSKV3KBGHygyQdpohaT+pUWNLMmZpANlYlPTL5t88odGkxWi9tlHQIkA0iKDWL+VSAxxSmiTvIRvbalbJXmbO0Qi5MuU89TrK3KDul12GC9zkgmGihwQjvNRBw4RUCOO+RSURNUvCgzDMAwjj/kUlER1LgrY81JisaNYcW3FlxGbUI/5boukK5Fjj7Z7IW3OYVWdaE1Ia2EYZYp2TEzhOKi2Mfjc3MCU4sSo7Xwo1ZGuPVRfWXNsuilerjknV4957pGa2q927fo0Pl1wvlXbo73T7jnenF6nzI6Dw42LY6L4uXWyFGgTTmk76Re81tcki4BWR2m7gh6PkZLqXBQYhmEYhoTJByVhiwLDMAyjdrBFQUlU/6IgQUrwlLqauU5LU6DKCtF5NMZbO8+hXHcedNe4Ksl0l0qaUzruJB84lCc5N6rnpa2fd2hMJxloOEkJQpGTrMAulCQxKDfh0h5Bf54pcglUkuSb5vPmYr53Kk8xAE6ypYMDotSmSQY1R/UvCgzDMAwjwiwFJWGLAsMwDKN2sOiDkqjKRQELs2ax1axWrK7aSPw09/LobZeogOQqQ0gPKZaulbTKLfazlfI8p1tOajON7AAk20vTtscqye2IUoLWDdX7X5MBEvqkvk+86NPakDV9L9V5o0yxn7e0Q1VslIPL+0XKFMOSK8UYVapyUWAYhmEYEpbRsDRsUWAYhmHUDuZTUBJVvyhI7Zkeva3Y+4ueCw4njtg8qyQ9rJKW3GUelhEb5mIz/bpsH5imGy5JhUbKrl9J86qCPm9lH5YKGuZKZ82aNfjGN76Brq4unHzyyVi9ejVmzpyp1t+6dSuWLFmCJ554AkcddRS+9KUvYcGCBazOXXfdhauuugrPPvss3v72t+Nf/uVf8JGPfGRY78N2STQMwzCMEti4cSMWL16MK664Ajt37sTMmTMxe/Zs7N69W6z//PPP4wMf+ABmzpyJnTt34itf+QoWLVqEu+66K19n+/btmDNnDubOnYvHH38cc+fOxfnnn49HH310WO/FC8OStgIZEQ4ePIi2tjYcs+pr8FtaytJmTa+oK+iXS0X9oqt0S0GFt1dJ06qi5lUFDUylfq8Fhw5h9+VX4sCBA2htbS1Po4OI/k5MvLr0vxPBoUP4w5fd+3v66afj1FNPxdq1a/NlJ554Ij784Q9j1apVsfpf/vKXce+99+Kpp57Kly1YsACPP/44tm/fDgCYM2cODh48iJ/85Cf5OmeffTbe8pa34I477ijl9oakuuQDD2X7ogvLbCOpoO8FwAuS64wUFTQwlfqFOVztlXsaVNLf4bJ/gEuhkgamUheCIzlEZQxJPHjwICtubm5Gc3MzK+vt7cWOHTtw+eWXs/JZs2Zh27ZtYvPbt2/HrFmzWNlZZ52FdevWoa+vD42Njdi+fTsuvfTSWJ3Vq1cXc0fOVNAnyzAMwzAqhwkTJqCtrS3/kn7179u3D9lsFu3t7ay8vb0d3d3dYrvd3d1i/f7+fuzbt2/IOlqb5aK6LAWGYRiGMRRljD7Ys2cPkw8GWwkog/N/hGGo5wRR6g8uT9tmOaj+RUGxCWLUfO3pruP8PkpN7DIC5w0HxbqsDENil8SuVGhil6LbUfpXfHsJ76OIqTdCiX+GlWK/pKsgQVfiFtWlJOgaLsq4KGhtbU30KRgzZgwymUzsF/zevXtjv/Qjxo0bJ9ZvaGjAkUceOWQdrc1yYfKBYRiGYRRJU1MTpkyZgs7OTlbe2dmJ6dOni+d0dHTE6t9///2YOnUqGhsbh6yjtVkuqt9SYBiGYRg5RiOj4ZIlSzB37lxMnToVHR0duOWWW7B79+583oGlS5fipZdewm233QZgINLgxhtvxJIlSzB//nxs374d69atY1EFl1xyCc444wxcffXV+NCHPoR77rkHP/3pT/Hwww+XdnMJVOeiwMVkFT1VFznAaVvZ5OvL55VJpkgzSytKPkhTt/j9KEL2j5HaVlY5IYz+p2wrm7D7cay+0PbAsZeirnx5l/qp+uRQbtuDp2hD2826zNuD6/WTtr+WJ5ZX5MexLIxCRsM5c+Zg//79WLlyJbq6ujB58mRs3rwZEydOBAB0dXWxnAWTJk3C5s2bcemll+Kmm27CUUcdheuvvx4f/ehH83WmT5+OO++8E1deeSWuuuoqvP3tb8fGjRtx+umnl3hzQ1NdeQr+LRd/aouChGu6Vx12RmhR4NROUhupFwUJbbrsNa/94bZFwZDU46LApbzYRQHrlUv9qNxLnlihl8tTsHRk8hQc+7V/KUuegheuvGJY+1upVKelwDAMwzAkbO+Dkqj6RYFoHQAKq1YXi4CLFSCpvsOKGVpfaXFSHZf88xX0MyosdrtXcp5qy9J+wXNdIV5f+0Ws9JUZ09ivee368UJVSnBB6WO+HReLgItFIqlOCZYPinrNhPPUNkYZJ0uBFzsY9H2UdB6G/CWer5IbmFA7T7mmFkDBPntC10PWuIs5avixXRJLw6IPDMMwDMMAUAOWAsMwDMPIU8Y0x/VIVS4KEiUDIG8DcZIJXOr4Wp2BctXsL5j2Ym3QU5V2RP8h1YmxeNuXZEYsyRVVNcnHy7kC4CIfKO0p5va8DKDVpXsFOEgMibIC23sgpZSQJBmQ9jWZwHO4HzYlg6EdKlNJDYPrKCTVL1cCqMRzS/kbUKTTHzPfO7StSQIeq+PFylgb6veYcqzYkqOp75H3K0ZKMJ+CkqjKRYFhGIZhSJhPQWlU/6IgyZFQ85rQVsxKuefLP4fy5ZrTjvJrn1sNSB3FshCVa9YBvVwsdjpXwsVxUA3fT7AasB/emqWAlZNjJU8B/TUf3WYYyD+/QuWHDlj9UC6nVYLIepQyxFJB+zWdbz1Q6irlmkVAsyzky1M6DhZrWShX6CMlzZd8OsdBt/JCCJ9ynZTOgHpa4tz79HuPzWsvVnegfihV0U7NWwhCaT4O0T+j8qn+RYFhGIZhRJh8UBK2KDAMwzBqhzLIB7YoqBWSYnqpBMCOyWlJMgG4c01Ux1fkAJ+ep0oD7nWY0qG0QXGRBvwUpr7A4cOSJBMMtENN//H3NfnApU6gyAPRNZlzFL0hxT7Lboc5D9JyweaqyRGlkGRup5IB6ZPHyuX21PIw4f2E87R+a3VczqN4Lp6wacbf4fMQatJQkgxQgnzAHQYTykP5fb2c/IM5D8rfe+KQm2RQE9TWosAwDMOob0w+KAlbFBiGYRi1gy0KSqK6FgVe9JJdXsWcBGklAyYTyMe+ICtwmUCTD+Q2NBkgI5zrIhn4mpQgliafGzi4EmufIe1cSUqgZVli+tZkh0CREgLl3MicTqMPAqbXkP4FYrEei02JrsOeFb2OLFlwM7wSIUCq5OunlQyUch6hEK+v1dWjHEgdTT4g8o2UtnmkIhHSRxwoD0Vr04uXhb489zSZQPL+j13HF96nUkMol+u3L0f45EszpD0lSsgLw7oO8as2qmtRYBiGYRhDYHkKSsP2PjAMwzAMA0AtWApUD94wXqZIBl4mIOWkDpMJCnWoJOAL8oFP6mY8uZxaDjOknMkHgsSgyQeaHEHxh1EoCxQDpJO5P3dM388qdbOBT84j1yHlWTK4tDwy4TPJgEoNZIIwMy8bT1KH2XlJlfzcU8zDaUnwxtdM9uw4K7eh1peOnepS+7RSJ0l6cIhgGJb0x0mkkAm0cm7qD8W6zKyvyA1MEqD1c02G2pzRfga6jBvti/CsjNqg+hcFhmEYhhFhjoYlYYsCwzAMo2Ywn4LSqM5FQVLEAT1WIgigJCNiJvmMLBlkWPnAsSYT0AiCjCIrNCjyAS2PTP8ukgE7Vpa8PnMNL44glG2RVErQ5ANJSqDn9ROzP5UPaHmWyQrknll0QaE8khWY1JAlcgQZk4CUU62H7aVQqCF7XivzVN0TgaBuqCnJB0rEgbbroyYrJJVrMkHZohzycogiQYyGlFAGyWCgXJgTSgSBJg2w/QlYuSAr0DZIhIA2JuyaSuQN1TLy81kOpOFRDvIlh5c6/qNeKuZoaBiGYRgGgGq1FBiGYRiGhPkUlETtLgqk/PPaXgaKZJDR5AMmDwwcN2TiZYAuE9A6DZ5cTmWAqA41+6eVD7SoBIokK2gyAa9DTPYp5YOo/X5ynUafSgOF8ozvi+X9gkwwUIc855wkQB3xKZpZn0oJav73ckUaJCBusaEmBlKOyQD4dDCCocv9rCwZuEUl0GMhYRGto0oG8nlOUgJFquPwzNwkA08pH7ioLhnIkQVqlIEiKwQZj15uoCzejVi5hraduDT3KuXvqPkUlIbJB4ZhGIZhAKhlS4FhGIZRf5h8UBLVuShISlgE5KMLtP0L+PHQiYmAQeb+TMG22piTDej7jQ4yQSNpo4FJAqS+cC6XCUh7iqxAofUzZbCPZQUJYDBUJugjdQJ2PFCHygc0yiAgkgGt05ctuFVrskKfIKv0ZbV7pxEPxNObbITg8X1lyXEYP2R225TygpZ8hpnbc6bilF7+mumfSgYelQqyQhta1AJrT5EbWL/idVRpQIl4YDATd+lzPHR4bqq5X5AYtKgBetO0TiQHDHkdEl3g5+45pGVUzqPdU79HyXGg1cnt70EjurQ5O8KYfFAaJh8YhmEYhgGgyEXBmjVrMGnSJLS0tGDKlCl46KGHnM77xS9+gYaGBrzrXe8q5rKGYRiGMTRhmV51Smr5YOPGjVi8eDHWrFmDGTNm4Oabb8bs2bPx5JNP4phjjlHPO3DgAC688EK8733vw8svv1xSpxlJJjAWcVB4m5nhlcRESZIBUDDxMzmASQlyeZPfL5YziYGcG0URNJAyFtlA7LmafJBBssSQBrZXgbK+pHX6iU0zyxISDZTTqIW+IEPOk2UFeg8ZJhkUzqWRJX3ajSQQhnLIQUA9xmnOfyFRzXDufVCu6ANJMqDHvproSDbr8/pyHWmvBHZtKgFo+ydQtAiFIuHJ0eSHqO9nQOZHNCW15EWarEBlKVqHygp0iGiiokKNQvcUKUHdulmJPkiae6OK+RSURGpLwTXXXIN58+bh4osvxoknnojVq1djwoQJWLt27ZDnffazn8UFF1yAjo6OxGv09PTg4MGD7GUYhmEYxvCSalHQ29uLHTt2YNasWax81qxZ2LZtm3red77zHTz77LNYtmyZ03VWrVqFtra2/GvChAlpumkYhmHUKZGjYamveiWVfLBv3z5ks1m0t7ez8vb2dnR3d4vnPPPMM7j88svx0EMPoaHB7XJLly7FkiVL8v8+ePCg28JAiERg+x0oyYvUbY+pDCBIBgPlA7bOJj8bKxsol2UCVt+XIxF4+cBxoyITMFO6IhNkaIRCGexj1NyfFaIJAC4rBGG/WKfPH7B5UnmBjgOLOCDSAL23fk+WFXyI9lQRuscCHR1f2fuAmaqpaTcyiTvscZAWyTPfaVtil/0JBMkAKMgAvEyRGugeFJp8oJRHfdGSG0Er1yIOyvHFLu1fECsn9elcoREF2Vy5Jhlk5HKWHEiZhzRCITpVDc6gt6BEGbA6yp4IUR+9Mg93WTD5oCSKCkkcnPktDEMxG1w2m8UFF1yAFStW4B3veIdz+83NzWhubi6ma4ZhGEY9Y4uCkki1KBgzZgwymUzMKrB3796Y9QAAXn31VTz22GPYuXMnPv/5zwMAgiBAGIZoaGjA/fffj/e+973O1w+RW8VqP8DYTnTRQeFtzbmQORpS5z3F6ZBbArJCWeEXcZNSTvMKNGeINYH8jKJWgbxDo/K+ZgXIKDkQMmWY9VmWzpjuXkisA4o1oY9aBXJ5APoDahmhjobEOsDSNhfup9crLuUGHYUG5rxF7o06j1HHL/LLNUsWxdEhc1JD8q9Mp9S95AQp1WzaHQu1XAK+6Ggo/9p3shr0Kw6DzIIQCvdArQNyG3r64zJ8s6tpi+kXC7kkrc9SFw/0hToIsrTF5H4CMpV9ar3KEOsAqHUgXs5+4NNLkmcV0L7SMXexIOSaCdl8pJ8T0kY5/kgbI0aqb9KmpiZMmTIFnZ2d+MhHPpIv7+zsxIc+9KFY/dbWVvz2t79lZWvWrMF//ud/4gc/+AEmTZpUZLcNwzAMI44lLyqN1D+vlixZgrlz52Lq1Kno6OjALbfcgt27d2PBggUABvwBXnrpJdx2223wfR+TJ09m548dOxYtLS2xcsMwDMMoGZMPSiL1omDOnDnYv38/Vq5cia6uLkyePBmbN2/GxIkTAQBdXV3YvXt32TvqhORoqO2SyKSEQDzWUhfLjobJkkEzsck2s/JCBD1rW5AKWBlpg8oB3BmRSgay+1GanAVBKGs33KFQlgmyLA9BQ6wOdSJsIPZMKisw+SAoTjKgptdQ2bmR7rRITbhBINv7PSYJjKyjIUORHVSTvEPq4sjk7CIZUJmA5z1QjoV+8fcVx07moKjIB5Q0UoJTPgJZJvA0+SAnG1AzPU0RTGUFljmYOCAmSQa0nEoD9PkwFYvNA0UOUcczd14lOhpWOH/+85+xaNEi3HvvvQCA8847DzfccAPe/OY3i/X7+vpw5ZVXYvPmzXjuuefQ1taGM888E//2b/+Go446Kl/v3e9+N7Zu3crOnTNnDu68885U/Ssqo+HChQvxwgsvoKenBzt27MAZZ5yRf2/9+vV44IEH1HOXL1+OXbt2FXNZwzAMwxiSSg9JvOCCC7Br1y5s2bIFW7Zswa5duzB37ly1/l//+lf8+te/xlVXXYVf//rXuPvuu/H73/8e5513Xqzu/Pnz0dXVlX/dfPPNqftXnRsiGYZhGIZEBcsHTz31FLZs2YJHHnkEp59+OgDg29/+Njo6OvD000/j+OOPj53T1taGzs5OVnbDDTfgtNNOw+7du1km4Te84Q0YN25cSX2s/kWBtqQT0hy7RBxk1HI5f0AkD2iSwWGZgjTAJQMaiSDLCpJUoEUfsEgEIhPwqIRh3CURcvQBlRWolNBH+puXDwRJAQB6yDTVdon0U0QfBIpkQFMl03kQEDNroMyhQEpzTNEiEVxIYcJ1iThwqiPIA6pk0CdLA76LfCBJBSwiIYi/DwzK7Uv6reSRKBqWm0CJHFAkA2QKlSLzvCQpDO4quzVy7LNbHlpK8Ji0pUQZKMdqnaShrUH9YHA23VLD5bdv3462trb8ggAApk2bhra2Nmzbtk1cFEgcOHAAnufFJIcNGzbg9ttvR3t7O2bPno1ly5bh8MMPT9XH6l8UGIZhGEZEGS0Fg5PmLVu2DMuXLy+62e7ubowdOzZWPnbsWDUB4GAOHTqEyy+/HBdccAFaW1vz5Z/85CcxadIkjBs3Dr/73e+wdOlSPP744zErQxK2KDAMwzBqhqFS2aRpAwD27NnD/vBqVoLly5djxYoVQ7b5q1/9aqBtwZKoJQAcTF9fHz7xiU8gCAKsWbOGvTd//vz88eTJk3Hcccdh6tSp+PWvf41TTz01se2IqlwUhEmSAQrmMy3fiBZ9kFGiDGja3QahjluUgXzcQqIPuJQQL6fvN3my1MDTHMsJizLMhlwcXCZQEhkp8kFvKEUfFO6hJ2gU+3qIlGsEoew/G0kFNLIhG5Bj8lyzSiSCRyQGJettfm4pgRps/nrK15f69SBFFygRCVxW0Dz+XY6F6AMlyoBJBv20nDSoSglBrK+sbkBt3JUjH/ByKhmQ6IJoDmVYpqNC3VD+/ChZhmnwA5MS8t97AZ2zEI/VJFKhLGtIc0sbYm1zxWqjtbWVLQo0Pv/5z+MTn/jEkHWOPfZY/OY3vxF3Cf7jH/8oJgCk9PX14fzzz8fzzz+P//zP/0zs16mnnorGxkY888wztb8oMAzDMAyRUXA0HDNmDMaMGZNYr6OjAwcOHMAvf/lLnHbaaQCARx99FAcOHMD06dPV86IFwTPPPIOf//znOPLIIxOv9cQTT6Cvrw/jx493vxEUGZJoGIZhGJVIJYcknnjiiTj77LMxf/58PPLII3jkkUcwf/58fPCDH2ROhieccAI2bdoEAOjv78fHPvYxPPbYY9iwYQOy2Sy6u7vR3d2N3t5eAMCzzz6LlStX4rHHHsMLL7yAzZs34+Mf/zhOOeUUzJgxI1Ufq99SoCQnisq5Fy71XKfHhdPUiANtb4No7wO6l0FKyYDKBC1MEojX5xEHWvIiUq7IB766j5o7gSfLB3yPg8IUo/1qDOPRB1QaoN7/h8J00zRQdmwMcs+KlREPcLobI41EoPPDV+aTx+ZeNPkSIhJcUUy4YvSBFolAz3PYB8FniYqiusVLBh6TG0g5ORf9kXygSA1MMqCbBSiSQZnlA6b5+rJkwB5AEK/DPfipmZ5KBoV7o1E92n4GdCdDP9qNUfleZNEPyvzQIg7E6ANtbo4mo2ApSMOGDRuwaNEizJo1C8BA8qIbb7yR1Xn66adx4MABAMCLL76YT3T0rne9i9X7+c9/jne/+91oamrCz372M1x33XV47bXXMGHCBJxzzjlYtmwZMhn3nWKBWlgUGIZhGEaVcMQRR+D2228fsg7deO3YY4/lW7YLTJgwIZbNsFhsUWAYhmHUFpVitahCan9R4CAfOEUcKNEHUSIjuv2xFDUAuEkGvLw3f9wk7X1Azmui+x0Q82OTsg8CRdsTQSKruKFQk30vCuYqev1eEn3AIiFy16dRIIe0bXIVLxguExTGkHpmR8dUJmgI5GdMIxH66d4YHo0+UGyrI02CvAAUH3FAj132NXCSDPrl40gSYJICkxLk6ANPKWeUYe8Dto0yjSIgA+ex6AJyHG17zOQNzaVLkRI8OfKFRSJE6gF7lsmRCJ42hJUoDyRguySWhjkaGoZhGIYBoB4sBYZhGEb9UOGOhpVO9S8KlNwieQdw9f1kL/IMMxsPLSs0EDM5lRe0vQxcJAMqNxTkA9IesrH3B9dhEQdKJAJFKs8qqXRYwiJidKLXp9EHvrInQyZXx6fu0JpkQKsQcyodZypxNJAoh4bc1sz0+bExIc+Nmc/VuYIhj/l2yiiaRJXCRT5QIxRCub60pbK6NwI5UUhGBAwhGdDySCrIxiWFwe1BTWTkEH0glTtIBmBJrGjYRmG+MYcwEkHj5aQCqhi4TIlQkavY3KL9iiJF6GdGed4eS54k10+SD9TohHjVEcPkg9Iw+cAwDMMwDAC1YCkwDMMwjAiTD0qiOhcFmt1NsGXp0gA5JuUNTDKQy+nWyY2+FBUgH/OoBJKYKEEyGKjTFzuP7YMA+ZpcPgjl8hSfAOrNz+WDQCxX5YswQcpQJYNC283EbEoTJtH77yd21P7cs2LRB34hIsIP5Oetba2tJy8SsgpRnOzGyeXSZVRpwGE7XCnigB6r+ye4bIucJBkAQH/uublIBlrEgSIrJMV5Uzx1vwNfLs+Q+6TXJOVRX+gO32Hho8ynBEuYJD9cj0kZ8eesRhlQySDlPghSV9RhHcU/qiYflIbJB4ZhGIZhAKhWS4FhGIZhSJh8UBLVvyhQ7DyR0csl4oDtg0Bmg5a8iNcZMHkyeYHKB1r0gSITSJIBPdblA3p9ua8ZVT5wh1ocs6qUoMgHNPoACfIBvSaRAGhioizZn4CNM0mS1EfkgSgSwSfTXnvGfYrs4hSJIN0Em6cpQxHUPZgH/X/QsRq1oEUZJOyVwOQATWrIyjIBkwQkyQCAJ8gHWiQCsoXzQjX6wCGpkUDIHibdWIDcA3P/J7nlWcRBPBKBzQIqJdCx76fyASlXJAO6pXU+EiHlM3aZK4mJjNg8HcW/qrYoKInqXxQYhmEYRg7zKSiN2loUCD+o2K85Uu4rlgJuQdCcDuNOfS6Ohm7H8i6J0THfRZFYB1gOgJCUU0tB/lC1Dkj7aWWFMgAISNtZctxHRtpnlgLZUlFoQ7Y8ZNlujLJDoctx5HRIn58fkDwKNF2tuqOmPJ886VukhNwEKtJlUv7KY79gtV+L5KEXHA3p+/Kvc25BUHZD1H79R8fEesAsAtrOiEFWLmc/vyHXifCVh8UsBSTNMbFAsZTL1BmRWhCiaUa3N8xSiwC1TjDP2kJzyi6VktOhR9umH2p1HnhKHflYqmvUBrW1KDAMwzDqG5MPSsIWBYZhGEbN4IUht4IU2Ua9UrOLAsmc6+IwppmKWYy9UMfFubApQRoAeOpiyQyuSQZNtFyRCagVMeN5Yh2JRnLMHA1DOg70OlRKkE3V1LSezd0zlQZo3oGsR1MYF05sInbRPnKsOR325Dy7Mg4SUdq5QhGlhGEgyZzr5GjoIDHk8yG4OBcqKY9FmSB2bpTmWJEMstl43YFKpO3CNUMtf4EE1ch8ObUwM7FTCYLsjKi6k+Y+b1QmYA6NiqzAHBeJo6PqdJjroyYNuMwDVYKCUseoKWp2UWAYhmHUISYflIQtCgzDMIyawaIPSqPmFwVqKlqCLhkQL34Wz56N1ckoOxBqnvC0jpangKYujvIQuEgG1NyvyQQZGiGg7RAnEKiSQQEqK3AbKjX50lIv9//C/bLcBKTndNdFl7HNCPJAAzsvWT5wkQlGSjIQcTD36pKBLAkwnSjyaGcmaYeIA1o/q0gGgqzgJBmoeQoc0h8nQVMV8+0I84cuz5tJCdH1iUxAIwTodVgEAxsfUp/JO0IkSCDLC3weyBEH6k3U8R/KeqLmFwWGYRhGHWHyQUnYosAwDMOoGUw+KI3qWhRET5tljSGHiZk2CriYil2Oo3S9viIZZJTUvjTiwFfKGwUzN0tM5CAZNCoyAZMPUiQ6DpQUyoEmGSjlNNlRJIlkiXxAowx6SQQBHZ8+Us7Hlsg7EJ5VmZ69irBDZ6jMWTZPQ5r4Rq7i4hme1IbatrZ7YnSs7KQHRVZwMuVLkgBt20UyoBIDaVtNf5wENc0z+YA8z0xh7nn0+hSa+Ciqo+3AGCqSQShfn45/KDw3l+fqMj9UpOtobYSenqrbqDiqa1FgGIZhGENh8kFJ2KLAMAzDqBlMPiiNml8UaLskami5+tX6QvQBT3SUXN6kec4LfWlkZQVcJINGEiOQSRFxQMmQNnxizuyjUgdzuy4c0qgEeh+RlKDdOx2fPnr9lOOcZPpP++y1XRJHHGbCdeh3kaZiaOZ4FqkAsQ6LOFDkhry5n+5lECpygCYZaEmN0kAvTyQA7REzKYElUiJyhy8lFZLHJ2TRB1CONflG6aTwvpsSRtuuEgnALAUlkWbXXMMwDMMwapiatxQYhmEY9UU9m/9LpW4XBWwbZcVWpO59gLh5WnvfZ9EHijSgmL5ZlENepij0jyYMoiYfF8lAizjwBSNpoNnSqFtzSPchIN7/Ib1n0qawVwI33ytjGCpjqER5SFKC9r4mL9B+VYkBdQAHmUD3GCf3nDvWPM09LZGOGpWgRCLkE++41E2WDEIt4kCSFTz580D3GwjpdtKskhJR4An3we6BbsusyArqeMqXydd32fuAklZSqmTCkN9/sW3UKSYfGIZhGIYBoI4tBYZhGEbtYdEHpWGLgkFQczJF2i4ZKJiWM4rNLcMkANkkzuorEkPBxE7blpMRsWNFMpBkgoH6gvFI9eImdcm9BcRcGZDrBEokQiYh+oBfUY4syCifYimCQ5JlYtdR5kFNoEoJ1Gwt1NdMqi6SgVP9gP8fgyIOpEiFwfUDuVy9fh65bkhM/ExKYCZ5uV8sEiF/bzRB1XCMZ3TtQhGPSHHY76DaseiDkjD5wDAMwzAMAGYpMAzDMGoIL+A+0MW2Ua/U/KLAY2b/kbmmU358QkYxXXJvfOl9es3i9jIAFMlAeT/rkBCGyRR0/BXzp5//vyYBpPuEph3/YvGZ0/no2RtTJyMqB2m9s5X6mjyQR4k4SLWXwRDXF993yURFr882F6BhNXGpgMkLafo3mJHyjk+Z7KgiMPmgJEw+MAzDMAwDQB1YCgzDMIz6waIPSqNmLQWeN7QV0PfC/EuvE+RfI0UGYf4lv194pcWHl39lPD//StU/ch5tLy1J95E0DsOBy/N2mTdJc69sJJhJ8zuNDzWEIcpjbk1DlFxGM4EHYeGVuu2AvBKu49I/2l5aku6j2P6VgsPzTjVvKhH2/Ep4DRN//vOfMXfuXLS1taGtrQ1z587FK6+8MuQ5F110ETzPY69p06axOj09PfjCF76AMWPG4I1vfCPOO+88vPjii6n7V7OLAsMwDKP+oIuaUl7DxQUXXIBdu3Zhy5Yt2LJlC3bt2oW5c+cmnnf22Wejq6sr/9q8eTN7f/Hixdi0aRPuvPNOPPzww3jttdfwwQ9+EFm6eZgDJh8YhmEYhsDBgwfZv5ubm9Hc3Fx0e0899RS2bNmCRx55BKeffjoA4Nvf/jY6Ojrw9NNP4/jjj1fPbW5uxrhx48T3Dhw4gHXr1uE//uM/cOaZZwIAbr/9dkyYMAE//elPcdZZZzn30SwFhmEYRu0QlukFYMKECXkzf1tbG1atWlVS17Zv3462trb8ggAApk2bhra2Nmzbtm3Icx944AGMHTsW73jHOzB//nzs3bs3/96OHTvQ19eHWbNm5cuOOuooTJ48ObHdwZilwDAMw6gZyulouGfPHrS2tubLS7ESAEB3dzfGjh0bKx87diy6u7vV82bPno2Pf/zjmDhxIp5//nlcddVVeO9734sdO3agubkZ3d3daGpqwlve8hZ2Xnt7+5DtStTsoiDJT4Sm4tXrjLwhJZvgtEfVocaUbbPdDkOaItj9PmmeAnX3RJd2Et8f+f0IXZ63y7wZMb+xhK44dHV0tn1M8sKkCSDSyaGDdjukaYZTPBS202EJ3wFJiVFGxBt18DWTq1TsvBkFWltb2aJAY/ny5VixYsWQdX71q18BADzhuYdhKJZHzJkzJ388efJkTJ06FRMnTsSPf/xj/MM//IN6XlK7EjW7KDAMwzDqkFHYOvnzn/88PvGJTwxZ59hjj8VvfvMbvPzyy7H3/vjHP6K9vd35euPHj8fEiRPxzDPPAADGjRuH3t5e/PnPf2bWgr1792L69OnO7QK2KDAMwzBqiNHIUzBmzBiMGTMmsV5HRwcOHDiAX/7ylzjttNMAAI8++igOHDiQ6o/3/v37sWfPHowfPx4AMGXKFDQ2NqKzsxPnn38+AKCrqwu/+93v8PWvfz3VvdS8o2EYevlXKeHPaQhCL/9yIQs//2LtwCOv+D5uAX2FYeFF/nO6fhjkX8W8Pxh6fd4vJNxH4X7Z9ZXxUa+fcvyLhc4nOs9GmtArvEBfw0mUjMHVNKnUp3HXYnu+r7y8wqsc/U17P/T6Wh+F9tn9lnL9tPWLhcwnNs+MojjxxBNx9tlnY/78+XjkkUfwyCOPYP78+fjgBz/IIg9OOOEEbNq0CQDw2muv4bLLLsP27dvxwgsv4IEHHsC5556LMWPG4CMf+QgAoK2tDfPmzcM//dM/4Wc/+xl27tyJ//N//g/+5//8n/loBFfMUmAYhmHUDhW+98GGDRuwaNGifKTAeeedhxtvvJHVefrpp3HgwAEAQCaTwW9/+1vcdttteOWVVzB+/Hi85z3vwcaNG3H44Yfnz7n22mvR0NCA888/H6+//jre9773Yf369chk0qW6s0WBYRiGUTNUeprjI444ArfffvuQdejGWYcddhjuu+++xHZbWlpwww034IYbbiipf7YoGITmgZ4lNjNqlo5M3Zq3PD0vS9oOFK9m2o50TD3+s2TiUCsq212R7f5GDfeKKV6QCLQoAypP0L5koRyTOvQqWWEMtfEMSL/peGYVmyZtJ3pWAXsmynVGIfJkxFDMvyExRYf0WzEqHsr0nnTsVN/n/wfffZL2z2N9JfX9wswKA+0ZCjKY8nn0fKF/g64PpV+szfy9lTI+LvUH/kendaieJxcb9Y0tCgzDMIzaoRyOY8PteFbB2KLAMAzDqB0q3Keg0qnbRQF95oO93fPlTDIgpmrBhK29z8zdinmcnhuwc4NYfWqOp5ZN6koSEDN9n0eyv4SklrITYJq94Kgc0EeyzASqTEDPpcdxsz4bk9BhDKGMfxh/Vtr7WrQCnR9V9V3hKceEUKvDzPNerG6o1NVN6VTfosfExB6Z/n3lQlTSIXOMdiUkk4xKCZwUjldUMsjQvmrH2n16Q9dVxipUj6Ece0J7cl2Gw1ypFjyUwaegLD2pTmpYODUMwzAMIw11aykwDMMwapBRyGhYS9T8ooA+W5fEMoHgrT5k/VybWc3c7VDeS8yZjV5/oY7Qlz5SllE8+2m3G0kxlRKoqdxPYTDikgaRCZiUkBxx0CfcmyYN9BLZo5RxTkpmlPbZ0/k0qt8hqqe53ClVMkhqX5UAoByT8SGmdy8gM0GQG0KfmPeJwxeLRFBir9nd0xtN4zjmx/s0UE5ujlxfjTjwhTqaTJDRJAgox5r0gKFxkRJA66Rou0Ko9JDESqco+WDNmjWYNGkSWlpaMGXKFDz00ENq3bvvvhvvf//78Td/8zdobW1FR0eHU8ylYRiGYRgjS+pFwcaNG7F48WJcccUV2LlzJ2bOnInZs2dj9+7dYv0HH3wQ73//+7F582bs2LED73nPe3Duuedi586dJXfeMAzDMBhhmV51Smr54JprrsG8efNw8cUXAwBWr16N++67D2vXrsWqVati9VevXs3+/a//+q+455578KMf/QinnHKKeI2enh709PTk/33w4MGBg3zybaoJkEPVPTpOoCUjSnkcebJTc3MfNXcrXvF9YWHom4hZn5ZTeSA6pomJqJl+kN1ULPeZKZ/IBylsZYGSpIhFGTApoUCfEHEwUO7n/p8hZYVxCJRx08e20E4A4VmV6dmr5Orw+QjlmJpn6XNL9jQv2lTs4twv1aHvayZ2zdxOo12YGZ5+gHPPjc7fjPy7xcsWPjNMSqByhJq4KwHWb/l+2D3T62syQFRHi2BQx02RGxSJIXpWLs/VZX6oSNfR2iiHPT8FXhjyZ19kG/VKKktBb28vduzYkc/ZHDFr1ixs27bNqY0gCPDqq6/iiCOOUOusWrUKbW1t+deECRPSdNMwDMMwjCJItSjYt28fstlsbN/n9vZ2dHd3O7XxzW9+E3/5y1/y2ztKLF26FAcOHMi/9uzZk6abhmEYRr0SlOlVpxQVfTB4288wDONbgQrccccdWL58Oe655x6MHTtWrdfc3Izm5uZiuhaDe4jLfXQxFfcRG1w/NU+HA9ECWSXBDjeJy1EG1LveJ7PRJybPqNyHLJ2wWUxuk5ryWYIjRVZIQpMJskodKhn0kjHqI+vR3twxjcKgY8KiDxzGlkUxCNEH/ew8OTrBRTJwmVsjgovpVzXzUpM0KffJg8uZqtW9B7RkRMRkTBP/0M1emJk2d02PtOEyM71QjmYIQ1mOSUSTCZR9EJgkQKMSpARH9N4ychssSRGtr0QchFIfFalHk6Vc5k31RB+YfFAKqRYFY8aMQSaTiVkF9u7dG7MeDGbjxo2YN28evv/976fe39kwDMMwjOEnlXzQ1NSEKVOmoLOzk5V3dnZi+vTp6nl33HEHLrroInzve9/DOeecU1xPDcMwDCMJiz4oidTywZIlSzB37lxMnToVHR0duOWWW7B7924sWLAAwIA/wEsvvYTbbrsNwMCC4MILL8R1112HadOm5a0Mhx12GNra2sp4KxzJnBso5t7QwWwcKNv0RnWYKTsoDCtNGNRLPOcbSZJ2FmVAJIMMMYtmckb5jBpxUDikUQGNNEJAkRUoUkqYrFAGDJYSCsd9zHxPzP0sQiAeLcDGEHIkAj3uZeXK+DOpx4v1Sdv7IO1coYyUlCBeRossSOl1Lm2jTM3QzCTNEhORiUAnHNUmiDQRkuiDvMmWPDNlinPzfUA3PKDtkRlKTfVSIiNfeWYsGRE102fEck9JcBRdP3TZPyHjicf0XO1ZFKIPqGSA2PtDlbvICqOpliViGQ1LIvWiYM6cOdi/fz9WrlyJrq4uTJ48GZs3b8bEiRMBAF1dXSxnwc0334z+/n587nOfw+c+97l8+ac+9SmsX7++9DswDMMwjByW0bA0inI0XLhwIRYuXCi+N/gP/QMPPFDMJQzDMAzDGGFqa+8DYXXHTL+k3MVUrJmW+wNiqs6ZERuIqVLziteOM6CSwdDJiyhZul0xsfM1KvsTaO1RxxKabEiCb4UsywQ8kVM8ymCgPC4PaDJB2vFMOqbPT5OIXCQlloNIsqcOx68N6TJOCWRoHWpaJqZ8WodYvsNsLvqAmf2JyZzIBGokApMJ5C2Q87IC/VaiEQ90zwSavIhek0oDVD5gCZEwNOpeBnL0BZMJfDmiQIo+YBEHVCbQIg7YXhJy1EhUnz0/p3mg1VHakc6rFEw+KInaWhQYhmEYdY0XpEtgqbVRr1T/okBz9or+r6RBVn8JkqVvP8tN4Ct1ol+fJFWxl+x0eChozB/7PslHEA49G/lOgrRPWVIuOzH6qqUgTZ4Cen15DLVcAr3QfsEPjNGhsDAmWp4COm4uzoX9QTy/RNpn7OJ0qGTeJoUl/KRiAqfw664E5zH+K1OuE00z9uuU/gjXnA65KUW8jke+gfI1iOWB/SKnjotaOuOAdYxc3n2Ou6Vtlu9ZzUMQlTcU5maoORo2yG2ozoXUMTGqnvIZp7U2iRaCehbia4jqXxQYhmEYRoTJByVhiwLDMAyjdihHnoH6XRNU6aJAe2CCDUwz92aJKTLLnAgL9rVGXy7vI45qDTnxqYFIBprTYQ8xfWeI09YhZv4k9yOkllJT+BKHKFqHywSBWE6RHRpl0zeXD7Q0zzSvwNB5CDSZ4FDYJLbXo0kJ2nHuudHnR59roMyDrDJv1DTH0lZ1FKfcvcnl4o54zGxMzcqy+V6rL5qcWZnsdAh1TIisAFJcyPZdkBKyNO6eyASknJnv2edHdmJkwyn9CtTSIKtpjmVZQc1DEOUp0GQCJhkoDoUJzoX0OCmPweA2tLTI6tySpCtKJTogGk5U56LAMAzDMARs74PSsEWBYRiGUTuYT0FJVP+igEUXxI+lsoHjZDNwlpqWidlP8ljvZx7vBe//HjLENIXxIWIGVyUDUh7JA1kqE5CdFgO66yDpSxOJePCV/UAzKfYJzSrbZVD5gF5f65eUh4DJB0Qy6GFSghxl0EOiD6SIg4Fjn/1/cL/p83aZH/rcipeVolFqJlwvdqAfO0kMniIx5IaQhmkxUzbNvM3GJ3lrFSYl5NrXIg5YXD2VI7KkY74y0Gm+5B2kBCYTOEUiePEyVT6QUx5rUgKEZ6VGGThIA2kiDhLzGBhVR/UvCgzDMAwjIgRS/M7R26hTbFFgGIZh1AzmU1Aatb8ocEhYRM3G1PSc0SQDGomQizrwlWRAPjHJUo96hioZkP76AwmIqfle8/JvJLIC9fLPhLSPpafsCpS+0KgI2i+tv5EMwBMTFd4/pEQcaImMaJ2eLLlOECWaGloKGlwnq0QoiBEHo4GLfKAksNGP457uRKFhVnrm/E8a8ZNCaQBueu/PHVOdgkQcsB0ISZIi5tGvaTrFoiYyUqIPmJRAxzBKXqRFGVApgXzutVTI5DhooNeRog+Q6jhVwqJKJEQZfArK0pOqJFn0MwzDMAyjLqh9S4FhGIZRP1j0QUlU/6KAeXh7sXKn/Q5IG1lFVujLFszTVBKIjn1i8vSJndVPubMGTULU7Mf70kiiCZpIwiRaTk3p9Ppsl8Qy7PiRVRIW8URG8s6IvYJ8oEkAfUpSIxpxoB2zXRpz8kEvfZ8816ySsChg5nGX5EXR/4WyYkgw5+qRBaQqPU8xGwcsCVF8Mwef3i/ZhY+OD1fCZCkhJHPcI5+lfB/pHgdUp2D7GtAIhcKhN4zyAZMGfKUOS1gkmPUVOYDWDVQpQUtkROvw/8fOUxIWqVJTUoRCJcoLAUrvSx1viGTygWEYhmEYAGrBUmAYhmEYOSz6oDSqclHAEqQgbuYECuZc+r6WvChgnuaFSv3EtEkjETLUez1nCu315KH0lXIKNUk3+4XIAWqeb8yVM5lAkQ9oMiIqE7C9D8qwzSkzsdNIBAf5QDrWtj/ucSjnx0SmEJIaqdEHWsSBmtQIhWMgXkcZYs8hUkHZLFlMSOQpJlx9XwPZ3p8UiRBmCufRraWZhZ90hUkJimTAttuNrkM+g57PwnFIp5SLjrp8oCQ4kqICHBITaZIBizggUkEh+oCWuRzLSY0Sk2EpU3lU/6RWuE/Bn//8ZyxatAj33nsvAOC8887DDTfcgDe/+c3qOZ6SUOvrX/86vvjFLwIA3v3ud2Pr1q3s/Tlz5uDOO+9M1b+qXBQYhmEYRjVywQUX4MUXX8SWLVsAAJ/5zGcwd+5c/OhHP1LP6erqYv/+yU9+gnnz5uGjH/0oK58/fz5WrlyZ//dhhx2Wun+2KDAMwzBqhzJaCg4ePMiKm5ub0dzcXHSzTz31FLZs2YJHHnkEp59+OgDg29/+Njo6OvD000/j+OOPF88bN24c+/c999yD97znPXjb297Gyt/whjfE6qal+hcFzNYl2HOZTCAfZ4npzFfKs0LCIoBHIiRBvfLZcS4xEcDN7VQSaMzZCFkZkRpoZAGtw6MPZJfadPcgm7Fov4NQlgxY8iJBEtDkBWrW71GjD2TJoDcbL6dbJ2sRB9rW2tocSow0KCW5kXaqFH1AIw4czMZqHba3wcANMckAipRAgwKUvD/smJzg5WQDLxuXFGg/aN2BCyk6DqUMex+wR6hsTaxLCQlJhVykBEUy4BKD0IYScaAmL9LqJEUaVEr0QRkXBRMmTGDFy5Ytw/Lly4tudvv27Whra8svCABg2rRpaGtrw7Zt29RFAeXll1/Gj3/8Y3z3u9+NvbdhwwbcfvvtaG9vx+zZs7Fs2TIcfvjhqfpY/YsCwzAMwxgG9uzZg9bW1vy/S7ESAEB3dzfGjh0bKx87diy6u7ud2vjud7+Lww8/HP/wD//Ayj/5yU9i0qRJGDduHH73u99h6dKlePzxx9HZ2Zmqj7YoMAzDMGqHMuYpaG1tZYsCjeXLl2PFihVD1vnVr34FQHYaDMNQdSYczK233opPfvKTaGlpYeXz58/PH0+ePBnHHXccpk6dil//+tc49dRTndoGam1RICQyCpk3spxshnqXB8SUTiMR+liSFbqHALHjpYAlTyIzuCEsmP77iW23zxfkA2JD5FEGcpIiLZFRsagJi6g0ou6JEJcK2J4SbMtjEnGQ1RITDS0ZAIVERew6DhEHWvQBlEgEJEQflEIoW6pJoVyXm4SVSATqxU6t8zlTtCoZKOVUmqAW9lCQDAAgmtqSpADwMLFQiThgSljZow/kcja2mnyQq6MlHaJtSHsZDJxL6giSAS2Xtr6OlWtShktCIkG6qhRGIyTx85//PD7xiU8MWefYY4/Fb37zG7z88sux9/74xz+ivb098ToPPfQQnn76aWzcuDGx7qmnnorGxkY888wzdbwoMAzDMOqbUQhJHDNmDMaMGZNYr6OjAwcOHMAvf/lLnHbaaQCARx99FAcOHMD06dMTz1+3bh2mTJmCd77znYl1n3jiCfT19WH8+PHJN0CwjIaGYRiGMQKceOKJOPvsszF//nw88sgjeOSRRzB//nx88IMfZE6GJ5xwAjZt2sTOPXjwIL7//e/j4osvjrX77LPPYuXKlXjsscfwwgsvYPPmzfj4xz+OU045BTNmzEjVx6qyFHgYMA+qizhmzs15Mit56z3iGk21HBpx4GVJ7nbSUCFWIB1MMvALMgCTD0hkQz+pE8kKVFLo8egeC0QyQCCXKxJDsWjSANuWGnI5Pc7LB1RSCKh84IvlvYEcoaDtbRDV1/Y76M9SCSQ54iBQIg7y81OzrSq5exITxYBLV5GpmlqyNS9yeh1tvwOffLBo7vyoupNkQKU1su0xi0QoTGuWeCgyZ3PJgIYzaDIBaVvLKFUs7Dl4SjmpryQBypvqNfO9kIAI4M8hjZQQKJKByz4ILlsqF5IXhfGyQXUHvrgxcgThoAlSZBvDxIYNG7Bo0SLMmjULwEDyohtvvJHVefrpp3HgwAFWdueddyIMQ/zv//2/Y202NTXhZz/7Ga677jq89tprmDBhAs455xwsW7YMmUw6ibuqFgWGYRiGMSQVntHwiCOOwO23355w+fj1P/OZz+Azn/mMWH/ChAmxbIbFUluLAsnRkFoE6AqY/uIjvwBYymPSXF/WfZIEzAGPOhEGYjn9JdxAfs3Tcj/3qBqI9YCmKm5guQnkvmoWhGLRLAJaHWoJoL/E+3M/awKWx0C2FGgWAe48KFsQREdDah1QHQ3JryslN4FYPkKOhuzXHP3VTOc77YtmQVB+UESzhqctlq0DHjMJkEO6k6Evl0d994iFgf3yp59lh3wEZZjiunOhVkdzOozGVnG41BwAtV/8Wh6CqD6zPNDnmjZPgWaFqmBHQ6M0amtRYBiGYdQ5ZbAUjO7uDaOKLQoMwzCM2qHC5YNKp/oXBYrDkZTmmO6+Rs2f1CyZzcpm8FTdUBzqsgGRD4j5ryGQcwk0+AW7XyQr+EGDWLfRS5YGWJ6Ccu+SGCbLB33MGTHugKjJBLQuczrMyhJDNiEPAZUM+hWnQzoPNOdCHisvOB1qczMtQzlwDfG+JivQfAR8V0PaULzDmrMgddRlcoAmK1CJI4gfe9TczZyDw1jdGLR+Gb7YNedCVkcx/UuOeU51HRwDk/IQcPlAkSPU/AWF42LnnlG9VP+iwDAMwzAighAlm/+HMfqg0rFFgWEYhlE7hMEgM16RbdQp1bkoUCUDKgnEi6DFGdNYemKX5LmoZfN4JBU0UBMduWiGmqR9ObKARhxkSB2frFbz8oEiE6jHkMs1fMEuq0kDvA65fyX6Qj8eaL9fkBSAQbkEFJmgX93tMJ6HIKumNpbzEYQ0zbESicCiD6JjbedElx8hSWZboJCi2OU6WoQCde5n1ePx9j6NwKGXFCSAgf7JsoKWhyBfR8vjkLZctfcLZQ6mby0OX/uOkepr6YT1SAT5WIsiyMsNQkTCUOfp/dKOo1AR0rZJCTVBdS4KDMMwDEPCHA1LwhYFhmEYRu1gPgUlUbuLAtEDXPYcZ6qDEn0gZZhil6Ema2IKzJBd6LIOsgI1rTaQ8r7clTwXyUCRDyiSTJAWNeIgpXwQjR09r1/ZmdBFJsgyGSBeP3CJMsjK19eSF5Ut0iABySTtKWZbB9WHSQYsIEcoZwmLtEgA+rlSZAVRMgAKyYuUVMWeMt56IiOlPA2aAuEkJdByL17mkiTIxdwv7cbosktiCdEH+WiKSpQJzFJQErYhkmEYhmEYAGrZUmAYhmHUHyHKYCkoS0+qkupcFCgmxVAyO1IPcSUPSUjtj75mQg7IMdmTIFc/IDKBz3bYI/IBaZvmLemn8gHb1TAuCVBpwEVKoJRjZ0SNQLGzJkkGtJxFHCh1WbQA9ZxXZAVpDwOnKAMtSVFSxAEgSlea6VtDM8tKUkHo03nvDX47dk1VVtC8x6NIHiYHJEcWqEmK6P1L5SwBEak7UpKBRhmkBH2fBHKs7UzoEKEAST5IuQOiLlkI0SeagjaasoLJByVh8oFhGIZhGACq1VJgGIZhGBJBAO4qW2wb9Un1LwqSEhmxrVdpVdleRiUItsUrs6mFsfr0PGaeJiY3j5inubm/UM7kA3pugnxA0eQDF4uedG7gYAvUjG3auZJ8wGWCoaWGwcesjnZurjxU9zJIlgZ4Ob0heuwJZSgebfgTTNIq2vddghlcy+UVOkQWhJoMwBIZCe+nlQkcxllq08ncnfAchmxTeFba3gcuEQqJ0kPKxERpkxpV9N4HJh+UhMkHhmEYhmEAqAVLgWEYhmFEmKWgJKprURAlqmI2MDmPet6cH2g2LSIBELsXi0pQTJ6BH78+3ZaZ2id9JgGQKooMkFSHWRAdpAStnOKnMPu5JPoKU0gGAHHWV+QAFlTiUCdQEgzlr6klsVLOUyWDhOgDNeJAKQ+VZ8/y6dN2cnY+NmfZheTnwCNvyD+UhF55dUuRAzSZQIscSIrKcDmP4rRFcprveIfPA9tSWTtXmk4J0QmudfQIhdLOix/TyCzhOqQ5fg/k+9Ub4WgEy2hYEiYfGIZhGIYBoNosBYZhGIYxBGEYsLwyxbZRr9TWokAyNTLzqOYmTA6Z9VWWKahUEJlWszQvPLUq07qsadnOqssHwnlae6xcLHY6V0KTBniddOdG5bpMQCsr8gG0OrQ89z8XmYAl2HEoZ9eXOlUmhEuKpn4MkhL45KOV5DpCgiE1oZJLgqEUddJHHCTPyRRTfHgjEcolHyRFpDjIBOI23IPbThN9UCmEYenmf/MpMAzDMIwaICyDT0EdLwrMp8AwDMMwDAC1YClQPZJzJulAMaFS0y+1LSobJGiyQtQ+s8gqnuOaRznrt+LBK9VWzf5pbKVDXT9HSYtmh0iEfJlynnp9RW5Qt8sO42VOMoFD5IDoUV+m5EXqnh2RuV3Oq6VG0tD56SnzWtorQZMDtIgDl3mjSgzS+5SS5qTUkRLaS5IMlLppZYpUcoMafaDIBGmjFZL6N5oEAQ+JKQbzKTAMwzCMGsDkg5KoykUB+zWieUflFnqeJ//M0n7Nq78Y6MJRcNpS26O9U6wAFM1SkNg/1kblTGgXx0Tx8+tkKVB+tWttJ/2C1/qaIsaeHSttOz0e1hXy647m1IhixbX2MmIT6jHfbZF0JXLQdLgHLbcHJdGakNbCMMoU7ZiYwnFQbWPwubmBKcWJUdv5UKqjfqdVitXASE1VLgoMwzAMQyIMAoQlygcWkmgYhmEYtYDJByVR9YsCJylBKHKSFdiFkiQGlyBupT2CbkZMkUugkkx3aT5bLuZ7p/IUA+Bk7nZwQJTaTCsZaCRICZ5SV/te09IUqLJCdB7N1aGd51CuOw+6P7eqkw8ilI47yQcO5UnOjep5aevnHRpNMqg1qn5RYBiGYRh5grD0VaNZCgzDMAyjBghDcM/wYtuoT6pzUcBiq8mhJCVoz1b1/tdkgIQ+qe8TL/q0q1fNLpzqvFGm2M9W2qEqNsrB5f0izd2lxNiHyryWpAS1G+WSY6K3XaICkqsMIT2keOiV9J1d7Oct5XlOt1zk91TRsmXa9oyKpzoXBYZhGIYhEAZh+h9gg9swS4FhGIZh1ABhgNLlAwtJTMWaNWvwjW98A11dXTj55JOxevVqzJw5U62/detWLFmyBE888QSOOuoofOlLX8KCBQuK7jQjrYdzvkJ57Vsu5jJvpOz6VeuaPbyUfVhGaJhTe6ZHbyu/dorutovsUWzbaamgeVVRn7cyD0slDbMrZikojdQbIm3cuBGLFy/GFVdcgZ07d2LmzJmYPXs2du/eLdZ//vnn8YEPfAAzZ87Ezp078ZWvfAWLFi3CXXfdVXLnDcMwDMMoH6ktBddccw3mzZuHiy++GACwevVq3HfffVi7di1WrVoVq/+tb30LxxxzDFavXg0AOPHEE/HYY4/h3//93/HRj35UvEZPTw96enry/z5w4AAAIDh0KG13dcq9oi5vc6VRSb9cKuinRtl7UuGWh2q1jDhRQfOqoj5vFWopiL67R+IXeH/YU7L5vx99ZepNFRKmoKenJ8xkMuHdd9/NyhctWhSeccYZ4jkzZ84MFy1axMruvvvusKGhIezt7RXPWbZsWZSSyl72spe97FUjr2effTbNn5xUvP766+G4cePK1tdx48aFr7/++rD1t1JJZSnYt28fstks2tvbWXl7ezu6u7vFc7q7u8X6/f392LdvH8aPHx87Z+nSpViyZEn+36+88gomTpyI3bt3o62tLU2X64qDBw9iwoQJ2LNnD1pbW0e7OxWLjZMbNk5u2Dglc+DAARxzzDE44ogjhu0aLS0teP7559Hb21uW9pqamtDS0lKWtqqJohwNB8fyh2Gox/cr9aXyiObmZjQ3N8fK29ra7EPnQGtrq42TAzZObtg4uWHjlIzvp3ZjS0VLS0td/iEvJ6me0JgxY5DJZGJWgb1798asARHjxo0T6zc0NODII49M2V3DMAzDMIaLVIuCpqYmTJkyBZ2dnay8s7MT06dPF8/p6OiI1b///vsxdepUNDY2puyuYRiGYRjDRWpbzpIlS/B//+//xa233oqnnnoKl156KXbv3p3PO7B06VJceOGF+foLFizAH/7wByxZsgRPPfUUbr31Vqxbtw6XXXaZ8zWbm5uxbNkyUVIwCtg4uWHj5IaNkxs2TsnYGFUPXhimjxFZs2YNvv71r6OrqwuTJ0/GtddeizPOOAMAcNFFF+GFF17AAw88kK+/detWXHrppfnkRV/+8pfLl7zIMAzDMIyyUNSiwDAMwzCM2mN4XUENwzAMw6gabFFgGIZhGAYAWxQYhmEYhpHDFgWGYRiGYQCooEXBmjVrMGnSJLS0tGDKlCl46KGHhqy/detWTJkyBS0tLXjb296Gb33rWyPU09ElzTjdfffdeP/734+/+Zu/QWtrKzo6OnDfffeNYG9Hh7RzKeIXv/gFGhoa8K53vWt4O1ghpB2nnp4eXHHFFZg4cSKam5vx9re/HbfeeusI9Xb0SDtOGzZswDvf+U684Q1vwPjx4/HpT38a+/fvH6Hejg4PPvggzj33XBx11FHwPA8//OEPE8+p1+/wimc0N16IuPPOO8PGxsbw29/+dvjkk0+Gl1xySfjGN74x/MMf/iDWf+6558I3vOEN4SWXXBI++eST4be//e2wsbEx/MEPfjDCPR9Z0o7TJZdcEl599dXhL3/5y/D3v/99uHTp0rCxsTH89a9/PcI9HznSjlHEK6+8Er7tbW8LZ82aFb7zne8cmc6OIsWM03nnnReefvrpYWdnZ/j888+Hjz76aPiLX/xiBHs98qQdp4ceeij0fT+87rrrwueeey586KGHwpNPPjn88Ic/PMI9H1k2b94cXnHFFeFdd90VAgg3bdo0ZP16/Q6vBipiUXDaaaeFCxYsYGUnnHBCePnll4v1v/SlL4UnnHACK/vsZz8bTps2bdj6WAmkHSeJk046KVyxYkW5u1YxFDtGc+bMCa+88spw2bJldbEoSDtOP/nJT8K2trZw//79I9G9iiHtOH3jG98I3/a2t7Gy66+/Pjz66KOHrY+VhsuioF6/w6uBUZcPent7sWPHDsyaNYuVz5o1C9u2bRPP2b59e6z+WWedhcceewx9fbW5D3Yx4zSYIAjw6quvDutOZaNJsWP0ne98B88++yyWLVs23F2sCIoZp3vvvRdTp07F17/+dbz1rW/FO97xDlx22WV4/fXXR6LLo0Ix4zR9+nS8+OKL2Lx5M8IwxMsvv4wf/OAHOOecc0aiy1VDPX6HVwtF7ZJYTkZqO+Zqp5hxGsw3v/lN/OUvf8H5558/HF0cdYoZo2eeeQaXX345HnroITQ0jPrHYUQoZpyee+45PPzww2hpacGmTZuwb98+LFy4EH/6059q1q+gmHGaPn06NmzYgDlz5uDQoUPo7+/HeeedhxtuuGEkulw11ON3eLUw6paCiOHejrlWSDtOEXfccQeWL1+OjRs3YuzYscPVvYrAdYyy2SwuuOACrFixAu94xztGqnsVQ5q5FAQBPM/Dhg0bcNppp+EDH/gArrnmGqxfv76mrQVAunF68sknsWjRIvzzP/8zduzYgS1btuD555+3tO4C9fodXumM+k8j247ZjWLGKWLjxo2YN28evv/97+PMM88czm6OKmnH6NVXX8Vjjz2GnTt34vOf/zyAgT9+YRiioaEB999/P9773veOSN9HkmLm0vjx4/HWt74VbW1t+bITTzwRYRjixRdfxHHHHTesfR4NihmnVatWYcaMGfjiF78IAPjbv/1bvPGNb8TMmTPxta99zX4B56jH7/BqYdQtBbYdsxvFjBMwYCG46KKL8L3vfa/mdc20Y9Ta2orf/va32LVrV/61YMECHH/88di1axdOP/30ker6iFLMXJoxYwb++7//G6+99lq+7Pe//z1838fRRx89rP0dLYoZp7/+9a/wff61mslkABR+CRv1+R1eNYySgyMjCvtZt25d+OSTT4aLFy8O3/jGN4YvvPBCGIZhePnll4dz587N14/CWS699NLwySefDNetW1cX4Sxpx+l73/te2NDQEN50001hV1dX/vXKK6+M1i0MO2nHaDD1En2QdpxeffXV8Oijjw4/9rGPhU888US4devW8Ljjjgsvvvji0bqFESHtOH3nO98JGxoawjVr1oTPPvts+PDDD4dTp04NTzvttNG6hRHh1VdfDXfu3Bnu3LkzBBBec8014c6dO/Ohm/YdXj1UxKIgDMPwpptuCidOnBg2NTWFp556arh169b8e5/61KfCv//7v2f1H3jggfCUU04Jm5qawmOPPTZcu3btCPd4dEgzTn//938fAoi9PvWpT418x0eQtHOJUi+LgjBMP05PPfVUeOaZZ4aHHXZYePTRR4dLliwJ//rXv45wr0eetON0/fXXhyeddFJ42GGHhePHjw8/+clPhi+++OII93pk+fnPfz7kd419h1cPtnWyYRiGYRgAKsCnwDAMwzCMysAWBYZhGIZhALBFgWEYhmEYOWxRYBiGYRgGAFsUGIZhGIaRwxYFhmEYhmEAsEWBYRiGYRg5bFFgGIZhGAYAWxQYhmEYhpHDFgWGYRiGYQCwRYFhGIZhGDn+P7BHUtYMBpqWAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots()\n", "\n", "im = ax.imshow(np.transpose(v.v()),\n", " interpolation=\"nearest\", origin=\"lower\",\n", " extent=[mg.xmin, mg.xmax, mg.ymin, mg.ymax])\n", "fig.colorbar(im, ax=ax)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Comparing to the exact solution" ] }, { "cell_type": "code", "execution_count": 16, "metadata": { "tags": [] }, "outputs": [], "source": [ "phi = true(mg.x2d, mg.y2d)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "With periodic BCs all around, there is nothing to normalize the solution, so we subtract off the average of $\\phi$ from the MG solution to ensure it is normalized (we'll do the same with the true solution, just to be sure)" ] }, { "cell_type": "code", "execution_count": 17, "metadata": { "tags": [] }, "outputs": [], "source": [ "e = v - np.sum(v.v()) / N**2 - (phi - np.sum(phi[mg.ilo:mg.ihi+1, mg.jlo:mg.jhi+1]) / N**2)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we can look at the norm of the error:" ] }, { "cell_type": "code", "execution_count": 18, "metadata": { "tags": [] }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "error = 9.754984685e-05\n" ] } ], "source": [ "error_norm = e.norm()\n", "print(f\"error = {error_norm:20.10g}\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and we can plot the error" ] }, { "cell_type": "code", "execution_count": 19, "metadata": { "tags": [ "nbval-ignore-output" ] }, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 19, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAh8AAAGiCAYAAABH4aTnAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjYuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/P9b71AAAACXBIWXMAAA9hAAAPYQGoP6dpAACjzElEQVR4nO29e5QV1Zn3/606fbobjPQoRBoiKubFiCGTaLNCgBd1okI0RnOZJfkxQ2KijCwmIqBjJJgBHUdGk0FiFC8ZDBpRWRPDxKyXGNp53zASiFECZrwsk0lQvHRLMKQBhT6Xqt8fferUd/fZT1fVuXT3aZ5PVoXtrqd27VOX0/s83/082/F934eiKIqiKEo/4Q50BxRFURRFObrQwYeiKIqiKP2KDj4URVEURelXdPChKIqiKEq/ooMPRVEURVH6FR18KIqiKIrSr+jgQ1EURVGUfkUHH4qiKIqi9Cs6+FAURVEUpV/RwYeiKIqiKP1K4sHHf/3Xf+Ezn/kMxo4dC8dx8B//8R+Rx2zZsgVtbW1obm7GqaeeinvvvbecviqKoiiKMgRIPPh499138dGPfhR33XVXLPvdu3fjoosuwowZM7Bz50584xvfwMKFC/H4448n7qyiKIqiKPWPU8nCco7jYOPGjfjsZz8r2nz961/HE088gZdffrlYN3/+fDz//PPYvn17uadWFEVRFKVOaaj1CbZv346ZM2cadbNmzcLatWuRzWaRTqdLjunu7kZ3d3fxvz3Pw5/+9CeMHDkSjuPUusuKoihKFfF9HwcPHsTYsWPhurWbanjkyBFkMpmqtNXY2Ijm5uaqtKWUUvPBR2dnJ0aPHm3UjR49GrlcDvv27cOYMWNKjlm5ciVuuummWndNURRF6Udef/11nHjiiTVp+8iRIxh/8vvQuTdflfZaW1uxe/duHYDUiJoPPgCUeCsCpUfyYixduhRLliwp/ndXVxdOOukknLj8RrgxHwSnXDGpbBGK2xhi3pmyL2bvdipvoqaXdqBvW5Uuc9kM8mtbtctTred5sFCFl6LsFmJeSu/IEbyx4hYce+yx5Z4pkkwmg869eezecTJGHFuZd+XAQQ/j215DJpPRwUeNqPngo7W1FZ2dnUbd3r170dDQgJEjR1qPaWpqQlNTU0m929wMt7lZ/u7gessLaRwXow3xhYx64YT9g/U7L/K7S9ov1IsfM0E7Yp+kiyj1Jcm3qg4+Kratyn0T2pAUVz/hMzHg97nalP195ETaIOo7k9qI873cH7L5iGPdigcfSu2p+eBj6tSp+MlPfmLUbd68GZMnT7bO91AURVGUcsn7HvIVDuTzvledzigiiYeHhw4dwq5du7Br1y4APaG0u3btwp49ewD0SCZf+tKXivbz58/Ha6+9hiVLluDll1/GAw88gLVr1+K6666rzidQFEVRlAIe/KpsSm1J7Pl47rnn8Fd/9VfF/w7mZnz5y1/GunXr0NHRURyIAMD48eOxadMmLF68GHfffTfGjh2LO++8E1/4wheq0H0AnuD283v9i17eVsGdKEozok3l8k4coiSbxLJvuRKI0BG2MTyrjt0GUjmqT9S40R53i/poa3qoTcmpOWU+E+I9Np4Jy7Fxng1HeLGT2ic5LgZRqkLixAbSw5pADpFsxb6I33uWNqkRbs/xBu4l8+ChUr9F5S0oUSQefJx77rnoKzXIunXrSurOOecc/PrXv056KkVRFEVRhiD9Eu2iKIqiKP1B3veRLz93ZrENpbbU5+AjhkxSrPcEW6meZ2+z5y3qnHGkG4FY0TtJqHLkhy/IHqKNa6936IPaJBjfZVt7twzXruTKNy5Aqbudr3dFEox0bD9/b0nPT7U/myi1JJVXRHs/3v7e9S7JbDHsHdt5+LCEEky5l1l8TIQb5wvfK75VDiF5UpKJPeFFML4nqR3eUfhulOQV473v5/ehGnM2dM5H7dF4JEVRFEVR+pX69HwoiqIoigUPPvLq+Rj01NXgw/ELLjzBRWiTVVg6MctOiW1vG6NtYfJzUJ88Skaoh90mEXFklygbvsTsHxOkFkPoEOQT3yVJyybN+HY3eSw5xt5FuwQTRy6phmRR4++vqrizY3zO8JmIlk58N9omUj6RJBWhbSeG7GJUu16JrSNJKrFsyrsRviiv2CUTqd6UWJy+2+DvPamey9xF/g4svMssxRjfqcY9dip7nxKiskt9oLKLoiiKoij9Sl15PqxIk0gDj0TePoHUKPM6RIK3QyoXf0zHmJwq2RiIqZDjj8R9yT2QIP9GnAmkhkeEbTy7jSPUB8dyHgEzh4dfYlvSNptDqi91SRj7q/3rTOpUJU3W8geZ+Ez0PRE0lreDbST7Qr0jeDscwwPm28tk77p2m6BebIO7JNgwSTwfcbwdjEf1nsXD0VMme8+11JGtF15Eow2afOrnyUa6n8F3pssPgvUj9Dsa7VIf1P/gQ1EURVEKeEAVkowptUZlF0VRFEVR+pW69HyI894sEockrxhlX6hPINOYtn7J/hKbGBNUbcTL6RBjYmBEPU8OlSZ/irJLjLLtuvgpe3vcAV+STBJIMLHklWpNPh2MJP081knIFUgtfK9YGkkVJoJK8oobvkCuIK+Y5dA+ZZFPUpIsI0gtriDNMK7wgnoRD52hHgvySp4kE/6K8ag+WBDNOI5kFI8OZDnGcwR9NM/vHr/MwX6qMq7bwL00+SpEu1R6vBJNXQ4+FEVRFMVG3kcVVrWtTl8UGR18KIqiKEMGnfNRH9Tn4EOSKSx5PETZRah3WVLx4tT7Je1JqdsTpWtHsggXJk60ixi1UpRdfOt+s0yNSDYpod5mwx+XJRjjmpCruFwJRlgOV4x8kcNn+p2yI3Kk46IiXMgmsdSSsssuptRSWu/GkFdSKbKhvqZYanGl+p5ySpBU2JbrXeGGVyPPh0cXUYpwMWUXqrfY5Ehq4c9gtOGxNBP2y6OX2fgjbFvBVsg3MuSkSqXq1OfgQ1EURVEseHCQr3D04+noqebo4ENRFEUZMni+OSG33DaU2lL/gw9JvigmGbPvl6QWMSKGZiDZ7OWEZL61vuy065IU49j1gKioFkCSXRzrflMusUszXsqe2M1jKcUmq3AdLPvNIlwkk2BsfTKko0okmMFIJVKLLaV+JVKLRV4BeskqBSlFklds0gkANAg2acm+KLvY5ZUGll3ohruCPZMk2sWUV1jqsMsuOSOqxbXXF8pZl/bnU1QO+5dzLC9HL4wIF7rOfnB+KRW7okRQ/4MPRVEURSmQr4LsUunxSjQ6+FAURVGGDDr4qA/qa/Dh01bAlCws7v6kUkuObXy7jaUdlmWSSjCxol0iNUh7JIIR+RJDdglkCDnahdeDsLfBn5OlFldaabNgkzS8LakEE3iQDeVkqEkw1ZZaQM9CBVKLy7ILR5CwZFIos9TSIMgrDanwJZTklUZ6yW2ySwO9yMZ5BHmFZZdUjCgYGyyp5CXZhaUWQV7JURgZ12cLL5xLUouZNC2sd2IkszDWiKFzBtKZ7wgvQa/nTadQKL2pr8GHoiiKovSB5zuRGWXjtKHUFh18KIqiKEMGlV3qg7ocfIjLmEQsdS8mAjNkFHJRGhIMlXOlcowh43il+0vqY6wFExn5IhCVLKrHhiQQSyKwyIRgADyWLFJ2OcaIVCEfrvHLolDtGrbCOWHHcP7yZ7PJJ5YlKgDzUvlCNrFYEozNuBZEtZ8ggViJvS2hmLhWS7TUYiQOs0gtXJakljRJLY0pu2SSpheRZZfGVPgyB7JKA+1POyzLkFwjyCtxIl9sxIlwyZJWyfZZss944Vc3yy5BPUfyOPnQVnpkfEMSDcsceeTz91TwjhlSnT3KztE/5IqFuhx8KIqiKIqNPFzkK1ywPR9tolSIDj4URVGUIYNfhTkfvs75qDn1P/iQZArP/Ld32RUSiEkRMSy1uCTHBMdK7fFxEKQWM/lZ3xEu0noviddzYZmEls0OXOhGVIsgtTi8JDd9Ht9IICbJMaXrRMRJacy/Z+JFx5RKMLYImN79cwwJIlqCYYofOekU/yqv21J2VEvvNoN6IcKlGlILEMoqaapjeSVtRLjYJZhGejlZdmmyyC7cBke+pIUomBSipZaUY38q85Zsd6akQlKLa6/nCJfufHiebif8Gg/6FScCx/jqFNaT4b54dC2cICrOsUgxA4zO+agPKvNNKYqiKIqiJKT+PR+KoiiKUiDvu1ZPU7I2qtQZRaQ+Bx9CFIjNE2rILjyLSJRgqF6QWqzRLkadvWzYCBEupuxil2aiMN3tLKmQTY7Xbildx8VvCI19vj4swdDTY6zhwus9GPKK3ZUZmthtJXlFlGBEj2nBVRyRhMzsSbxEZFGnr5qEbGlHbLsaUguXhfVZzLVaoqUWKVlYY0NBdqE2mhrCF49lkuZUNjxOkF2aWHah+iCahevSguxirP9C99tcFyZ+ejyOcOE/kIbsIkS+dNMLZ0belEbhxInAMd4xuoZ5llNJnuX3Oiga8qS09k8/48GBV6FT30usmSpJUdlFURRFUZR+pT49H4qiKIpiQSec1gf1P/iwJKsC7Gu7cFlat0Vc80Wyz/ol+90suWTzgtQi1Bs2QkRMFJIb3k+Fji4xmqPgZjUSiLEcQXKMsYZLg3AfjFn09nUggmZMGYWTO0VLMMkm2kckIevVGUmC4ZOa1zzIcEemrKYl/V6LimaJkTSsbKkFKMotjpRMjGWXMqUWAGgqRrvYo1QMqYVe1GFUz/amBFMqsTS72ZK63mWWVNKGbhuSSrAqEeef8ATZxSiTzmlE53j8HcMSUOlXuiGdGGvL2NeTydO94gRmxjJRwTMn6t4D98e7OnM+VHapNSq7KIqiKIrSrwytwYcfbk6weQm3vF/c3DzsW84vbkXbrBdueb+4OVnPurmZvHVzsuFmrc/kIjepbTeTo004f8aDkzH7l8p6xc3NhJuT84ubmw23VCbc+Fq5OYQb2Qd1Dm1sa9wToz7cjHrjftJmrXOKGzwUN6MN375ZnzcfhWU8HWN/1R5xJ9zCk8K6+a5v3Qw717dvqXBz3MKWCjfXDbdUQz7cUl5xazC2fHFrbAi3phRvOTSlchjWkC1uzSnecsXtmIZMcRuWyoabmylu70t1F7djU0dKtvfF2I51eTts3yxt996Gu90Y7nYbxwV1w91u8fzD3Yx1a3azxY0/f5ObK2z54tZobLnilnbzxS3levbN8Yub63rFzXH8Hu8HP0uDhJ4Jp5VvSVmzZg3Gjx+P5uZmtLW14emnn+7TfsuWLWhra0NzczNOPfVU3HvvvSU2jz/+OM444ww0NTXhjDPOwMaNGxOf90c/+hFmzZqFUaNGwXEc7Nq1q6SNc889F47jGNsXv/jFZBcgIUNr8KEoiqIc1XiF9OqVbEmjZTZs2IBFixZh2bJl2LlzJ2bMmIELL7wQe/bssdrv3r0bF110EWbMmIGdO3fiG9/4BhYuXIjHH3+8aLN9+3bMnj0bc+fOxfPPP4+5c+fisssuwzPPPJPovO+++y6mT5+Of/mXf+nzM8ybNw8dHR3F7b777kt0DZKigw9FURRFqYBVq1bhiiuuwJVXXomJEydi9erVGDduHO655x6r/b333ouTTjoJq1evxsSJE3HllVfiq1/9Kr797W8XbVavXo0LLrgAS5cuxemnn46lS5fivPPOw+rVqxOdd+7cufjHf/xHnH/++X1+huHDh6O1tbW4tbS0VHZRIqivCafszi4gznfye/2L3pM5uV4qx8jdUZhcakwspQmnZp4PapzKjjG7034e2CZASZOizJlhdM7SiaU9JjRxsqGnTU65zpNMXb4+Rp4AKqfJ3uivUA5bDNug2qQTUQ37RE949ERUnscmZfwIs6tHp2WXc3QINrbJpXEmkyaYWArYJ5fyxNIUTRQtd2JpT5kmhRZyevDEUi43GZNMM1Qme5pE2kTLUhv1hTJPIG12wvZ4Jds0tZGS8nwIE05tv57N3B6UWwP2yac8mfWIn7aex3puYTJpjuqN1YM5jTznDeGU6pzkp2DjkK2Y56OoR/YP1ZxweuDAAaO+qakJTU1NRl0mk8GOHTtwww03GPUzZ87Etm3brO1v374dM2fONOpmzZqFtWvXIpvNIp1OY/v27Vi8eHGJTTD4KOe8fbF+/Xo8/PDDGD16NC688EIsX74cxx57bOJ24lJfgw9FURRF6QOvDNmktI2ewce4ceOM+uXLl2PFihVG3b59+5DP5zF69GijfvTo0ejs7LS239nZabXP5XLYt28fxowZI9oEbZZzXom/+Zu/wfjx49Ha2ooXXngBS5cuxfPPP4/29vZE7SRBBx+KoijKkCHvO8hXmFI4OP7111/HiBEjivW9vR6M0yve3/f9kroo+971cdpMel4b8+bNK5YnTZqECRMmYPLkyfj1r3+Ns846K1Fbcan/wYfkyreEoYtyjZFng+qFVOuGlFIoi1JLNmzEYdnFs9sbCTbyFjkmTvy58OA5Lv0aMFaypfrg/KnQxcrXJ5BlANO16/qmOBJ2W0iYIYoWpftrKsFIypWR3r08CcasE5J+xCEij0diqYXkFQi5O8yU6T3lWkotQCixsNQyTChz3o7hbiiZsLzC9U1U31iQUpodKc8HSS2c5wNCno8Y6dUDGcBzwocmY+T2CB9OY/VcTlefIPdEjt5pI107PbQZklFSVDbTtYdtssRS/IoZJCnVa8WIESOMwYeNUaNGIZVKlXgb9u7dW+KVCGhtbbXaNzQ0YOTIkX3aBG2Wc964nHXWWUin0/jd735Xs8GHTjhVFEVRhgyVRroEW1waGxvR1tZWIlG0t7dj2rRp1mOmTp1aYr9582ZMnjwZ6XS6T5ugzXLOG5cXX3wR2WwWY8aMqaidvqh/z4eiKIqiFPB81/DKltdGMu/kkiVLMHfuXEyePBlTp07F/fffjz179mD+/PkAgKVLl+LNN9/EQw89BACYP38+7rrrLixZsgTz5s3D9u3bsXbtWjz66KPFNq+55hqcffbZuO2223DppZfixz/+MZ566ils3bo19nkB4E9/+hP27NmDt956CwDwyiuvAEAxquX3v/891q9fj4suugijRo3CSy+9hGuvvRZnnnkmpk+fXt4FjEFdDj7ieK2L2X9jRbVw2beWjVVwLavTsi1HtThxIlxy5M5lOcYT5JiAONEuQj2nWmfZxSmkU2YpiKNd+JwuyzGGpFWaOr1whL1fFB/S9/7qSDA+r8wrnJF7YrYXX4KJpTJFpU7vjS2ypUpSi2ukTC9dnbaWUguXJamFZZThKbukYkgwjj0KJqhneaVRkF2MFW6FpyUlPLfG+iCFIssrfB6OfOGoFtdnqZaK9C6zlBP0l2Upjp5J05edGO0ilF1DdgnmJ0ABMHv2bLzzzju4+eab0dHRgUmTJmHTpk04+eSTAQAdHR1G7o3x48dj06ZNWLx4Me6++26MHTsWd955J77whS8UbaZNm4bHHnsMN954I775zW/igx/8IDZs2IApU6bEPi8APPHEE/jKV75S/O8geVgwebaxsRH/+Z//ie985zs4dOgQxo0bh09/+tNYvnw5UiS/VxvH9wd/EvsDBw6gpaUFJ9/yz3Cbm0HfO3AzDpXD+sAm1c22fsn+njLP0fDt9cax9EIWBhQ854PXdnFozkc9DD4Q1FMdDz6MMrXhp1NkQ1+MZM8huLwWjJcurCdDQ2Fjv1EObXzDnupTgn3K/BcAaMVycw0bfud4SozxR9x+bHEuhiPt94X6aHvznDr4AOpk8FGABx95Yy5GWM+DD54XcsQL69/zmoRyY8+/+cZi3aF8uP/dXFg+SOX3cqH9u9mwfDgbnvNwJix3Z3r6m8vS5+mmlyZDP2gyDrzDR/D69d9EV1dX5ByKcgn+Tnzv120YfmxlfzTfO5jHvLN21LS/Rzt16flQFEVRFBseUHG0S/ylApVyqf/Bh+S3sSUZi4iMAcqQY4JoF0FeMT0Zdm+HI3g+yo52YYSEY0Zf2MvhlXo+HHYD8PlFzUJIMmZg2xMVAWPaVCLB2I+zU4kEU9zP+eKSfi+KyZssdVX2dgChx6OW3g4g9HLE8XZIkSzHuKGrM4nngz0cjRTVwqvaGknGEi7YE3g5+DxmMrHo9lhe8ejz2BKU8Qq4aY8/T3ieBk6UZkS4lMorpfWWDvZjIjGl/qn/wYeiKIqiFKhOkjENBK01OvhQFEVRhgzVSa+ug49aM6QGH+I6LwFSkrE4ycdEmaY02sWREojFkVq4niWOfL5QJbg2WQpy7X59R1rzxUgo5pbU8TkdIVTEIRtJmRksEowkrySVYIw1cYyHIthPxxnPD8lffN2kSabGOi+WemmtlipILUAot9RSauFyHKllOMkrXLbJK73LweRSOcmYXWpxheeQ10KR5hoEib44yVcqoUyRJ9mFpRaz7/mSvrJcwnKMIbUI9q7Qx6Ico1KLUiZDavChKIqiHN14cIysruW2odQWHXwoiqIoQwaVXeqDoT/4iBHhEqfM8olRDuQO9qsbuoM98kWUWvJh2TfsC23SzHURQTPwDdmFo1m47ynzfCApBr0kGE4yRo+SEXwhdMvl/yg2E5WEzDAWbaIkmKgIGKCP5EniWhala7dUGO0Xfc7A5W1ILbSbZRcus1s9QmrpKfc8c7WUWgDgfYWkPPJaLdFSixjtYpFYGunus9Qiyy7RpAUZIlu4cSxvZI0Wc7DBuUDSDq8Fw5Ey3Pe88W/vsimvePayEMXl1InEkjQ9utSGUlv0CiuKoiiK0q8Mfc+HoiiKctTg+Q68SpOMVc1tqUgMrcGHJJ9YcIwIDnu9LMewrFKIdhFkGU4UZiYnsyQQQy+pJZ8vtfctUkxcXLvsYix7XziP08CJxSgahpozIj/4PL49tbEhgdiVAsGaqVyCMZKDCfpKrDUrBDmkeDv5uRK6Ld7BqAgXLhu3VYhwcaOlllREErFaSi1AKLcMT9mlE5ZajjGkGSGlulHmhGJeoY5Tp1NECF1PfgqTJuxmgTRoPyvcfM+QVzgRWfgVzSngs0I6eNfpSYGeEhOI0f2OkVhMori2S6Rl/+NVQXbRPB+1R6+woiiKoij9ytDyfCiKoihHNZ7vwqswWqXS45Vo6n7wEeUhdJKuhRKVqAy9ZRpbGzEiXywJxHpsJGmmp+xznW2l276gF8pwrbKUU0gu5vPKvMLncSzyExBrxXi4dGxwpmRJyHq3Hl+CMVvjftAKydS0cZX5lJwPzsjfFmQZExKIVYKxwm1ptItD0omTUGpJRyQRq6XUAoRyi7kybViWpJZjHEl2KZVaeuqDRFyg/WGZ5ZUU6WVJbyG3kyk873yePEs99JRxojKWVHiFW5ZVzHP21LvCfimxWByi5BgjqV6ilqtLHo51ReGkbSi1RYd3iqIoiqL0K3Xv+UiCOHCP4R2J9KBI+2PUGynTOY8HeSSKHg/JexLHC8I5RFxOCMGr3Rba5PTqnIOAPR/CaQx/RAz7cCJoaV3pgRV4QQqf0xOyzDu2yZy96yVvh6XzZkr16O6ZfY1fNvN5UDc4dTrZGBNLXfZ8UNkt9YhU4u0wU6OXeju4XvJ2xMnhwd6OJl6p1sjjUfgXoLrwwqYMD1iyX8CecNMbC+1kfPZ2UFr2GCvmGt4On+v7fqAkLwkjpVGX6gczKrvUB0fV4ENRFEUZ2uRRuWwSI5WjUiE6vFMURVEUpV9Rz0ctSCrBSHgWWUWQWnxh8qcBz6IkCcZYnTVox5hYSmNUlouo6VgSTIQ9u7hZGnHtmacRS4KxpSOXeksuZo+lKEmO8ez1gVrm2C9bciTZpSCliKnTY0gtnEa9kfN4kMQSTDRtdO2ySyypJWVPmW5KJj021ZJami1SCxBKIGlBXkkJz4crPG+esbwz7TDuuV9oOySLaCQJxhUXBlBUdqkPdPChKIqiDBl0Ybn6oKwrvGbNGowfPx7Nzc1oa2vD008/3af9+vXr8dGPfhTDhw/HmDFj8JWvfAXvvPNOWR1WFEVRFAkfDrwKN19DbWtOYs/Hhg0bsGjRIqxZswbTp0/HfffdhwsvvBAvvfQSTjrppBL7rVu34ktf+hLuuOMOfOYzn8Gbb76J+fPn48orr8TGjRsr/gBJUvCLtlLua+PYGj6MSdKkx5FafMElK8yi8i36QBK5JK5NIgmGymIUjHRExMqzvpFKmuQVIwCIU0+THMQ2HHxgqffpgXOM6CbqEt82sd9cz5EtgexCu6UIF0Fq4dweRtkl+aIgpTRzng/aH0dqMaNd7KvTBnJLLaUWIJRb0nThWFJJJX3XfeFpZZUzQnczIlaqHGAS51e8tJZJkjVOKpIWlaOOxJ6PVatW4YorrsCVV16JiRMnYvXq1Rg3bhzuueceq/0vf/lLnHLKKVi4cCHGjx+P//2//zeuuuoqPPfccxV3XlEURVGYQHapdFNqS6IrnMlksGPHDsycOdOonzlzJrZt22Y9Ztq0aXjjjTewadMm+L6Pt99+Gz/84Q/x6U9/WjxPd3c3Dhw4YGyKoiiKEkWwqm2lm1JbEsku+/btQz6fx+jRo4360aNHo7Oz03rMtGnTsH79esyePRtHjhxBLpfDJZdcgu9+97vieVauXImbbropSdd6sLitfUFSkerN1UTtp/EtNo6UrUrMYkWIubwjEKQWP+Fqt+yqD3QFnzQah6Uez77CLZNYgilcl6gkZIAswZjpzd0+bYzEUUaufHuUTJ5lDZKunAgJxvCkVxTt0veqtnEiXEzZhZJvRUgtANBYqDejWkIZpYmiYMqVWoAwZXq1pJZmllroZgXRLJLUIkW1iAj6Wz5Btoi88GXD+Sr41zivusrHBqu5StEa/EfVE944X+gL1wdlVVqUcinLt9R7GXLf98WlyV966SUsXLgQ//iP/4gdO3bgySefxO7duzF//nyx/aVLl6Krq6u4vf766+V0U1EURTnKyMOtyqbUlkSej1GjRiGVSpV4Ofbu3VviDQlYuXIlpk+fjn/4h38AAPzlX/4ljjnmGMyYMQO33HILxowZU3JMU1MTmpqaknRNURRFUaoim6jsUnsSDT4aGxvR1taG9vZ2fO5znyvWt7e349JLL7Ue895776GhwTxNKlg5tb+nR8eSUSSXPLcT8WDGmS1fy+gZRop8YRP2Dhc+s8MuVlrnxeH1YfgwQepJIsEYLXAUir3a+A9zDYrwM3sswQRee77fJHk5grzCUobH9hzYwI9KEO3Cl154lMTlhqRnz1jHpXCvBNmFr0mDsG4LJxxrFOSYQG5heUWKcGmiBGJcZqnFSDJmWZ02qdTSTNcnHSG19NT3PM8sr7gJQyu9GgoOcVKDR02INCQallp8u1wj/bFVWUWpFYlDbZcsWYK5c+di8uTJmDp1Ku6//37s2bOnKKMsXboUb775Jh566CEAwGc+8xnMmzcP99xzD2bNmoWOjg4sWrQIH//4xzF27NjqfhpFURTlqMaDa8yJKbcNpbYkHnzMnj0b77zzDm6++WZ0dHRg0qRJ2LRpE04++WQAQEdHB/bs2VO0v/zyy3Hw4EHcdddduPbaa/EXf/EX+OQnP4nbbrutep9CURRFUdDj0ZEm8CZpQ6ktZaVXX7BgARYsWGDdt27dupK6q6++GldffXU5p0qELQolydLkfZWNZ5FzcqUK0SHk1kaOXfn2pesNWYHXVmEb19AHev7N5+37ky7D6AtRHgWtwKdsWg41zgnJxKRhbrQL12rBkoZxXPiYShObzcsW9tHPswxRqKPz8O3xDUnDLq8gTrmgUhnKSSU5pCKCp3htFymxmFTmCBdeu4XLwdouTYZcYi8Pl9ZlIZuoJGJJo1okqSVNK6nYolkkqSXlsDSRbA0VTwhXC2rzlrqeMkeslEavAECWPg/bZP0GKqcs7XGUjF2CMcqGTWmEi1HWP9JKmejaLoqiKMqQQSec1gc6+FAURVGGDH4VVrX1NcNpzan/wUdUCIUgncQrc8QHlVlWcXoXAJAt8twGSRYeh0qwL921loOIBp91gnxCrUWKLrLWc2Ix6ge57DlKRpRgpK5YbAz5ifezjMPfCVR26Z54QhRMUO/mKNkct0Hn8e23AU6KpBmv78gXMZgrju4iPtelWo4R4WLIK2E9R7g0GlKLvdxkkWM4wqWJpBGWWuQEYXHKBXmnggRiUVJLT9kp7I/+A1ORBEMPQLC2C7eQF6QW/sMZyCg9bXB9+NVtyCoFG7bNefwZ7JJKHKnFrEcpg8RbkIcTK2Ioqg2ltujwTlEURVGUfqX+PR+KoiiKUsDzK5+zkXCFCqUM6tLz4TvhlsiWNt+NsaVoc53YGxzaXDd6S6WKm5Nyixt4c52ejdp2aIPjFjfHdYqbAfcr8sL5tHnFzff84mbU+35xQz5f3Py8V9y43rrlws2hzaz3ord8uLm0hfv94ubSxvWOxxuiNz/cYN0SPLS94cbpGS4WHb+4pXhzvcitMZUrbk28uaVbs5u1bk20STbNjrTlilva8ZB2PDTSlnZQ3Bodp7ilHbe4peCEmxNurvE/p7ilHDeW5AL0SC3Bxnj0v7zvhxukrSfSJevz5lq3DFLFLes3WLcgnDTvO8j6qdhbzqPNd8PNC7c8bbzQGj/Otq+JwYJXmPNR6ZaUNWvWYPz48WhubkZbWxuefvrpPu23bNmCtrY2NDc349RTT8W9995bYvP444/jjDPOQFNTE8444wxs3Lgx8Xl/9KMfYdasWRg1ahQcx8GuXbtK2uju7sbVV1+NUaNG4ZhjjsEll1yCN954I9kFSEhdDj4URVEUZbCwYcMGLFq0CMuWLcPOnTsxY8YMXHjhhUbOK2b37t246KKLMGPGDOzcuRPf+MY3sHDhQjz++ONFm+3bt2P27NmYO3cunn/+ecydOxeXXXYZnnnmmUTnfffddzF9+nT8y7/8i9j/RYsWYePGjXjsscewdetWHDp0CBdffDHySecVJsDx+z3HeXIOHDiAlpYWnHzLP8NtbgalDECqmyaVZag+Y6vzBVs/suxyfbdH9Xnj395lJ0s3L5uj+py1Hrmw7OdK63225VTneV5ZkyeFCg9PktsuTIg1PCucgp3t0+mwnibcglLuO0GZ0/A3hO35aaqnsp9OWcteI9fTZLtC2WukCXqNNOmPZjTmqd4z6mGt96g+sPGo27zfa6RJsGkuhzZ+I93PdFh2qD7V0FNON1JejLS9PCxNK9I2hOX3pWmF2QbK0UGr0x7T0GNzbOpIsY4nmb6P6o8xVqyllWwdez3n9GguTDRtpgmnTfQoNdJz1RSROh2Q06fH9XoA8iTTHGXsyNO7lKX6LB17pGDTTa/dEZpMyuV3/fBhec8L17h61yiHNge9YcXyoXxzT13hXwA4kAvL7+bCNg5S+VA2LL+XpfNnw4fySDZ8oDOZnnI2E9Z5mfAzIEPfExkH3uEjeP36b6KrqwsjRoxALQj+Tsz9f/8fGt/XGH1AH2QOZfCDv3o0dn+nTJmCs846C/fcc0+xbuLEifjsZz+LlStXlth//etfxxNPPIGXX365WDd//nw8//zz2L59O4CehJ4HDhzAT3/606LNpz71KRx33HF49NFHE5/31Vdfxfjx47Fz50587GMfK9Z3dXXh/e9/P37wgx9g9uzZAIC33noL48aNw6ZNmzBr1qzIz18O9T/nQ/Bi25KMGR5voewLEQ+8poZRX4hs4UgWv4G+YSiqxaFZ50YECZVBX0LgAUXhjzuvoWKMG7nM7RnZrZLN1re2LUXBcCKyOGvBUBSDj1yhzp5Ny4iwoT9MRnQMRQE5ObLnZGFB8i+OTOGIHbr0ho1Qhv2yFG2kJGOm35o/MzcCa5nXcQmuC0e78P6UsIaLUXaEMkW7pAsXKc2DAvoF0GgMIOyRLI1OaXs9ZV6jJTgP6DgaWPAAooYDDglez4UHHJxYjCNcska5cJyRHMwe1cLljFB/hAYotmO5Lmfst0fBeGJUS3Q5rCytGgiqmeH0wIEDRr1t0dNMJoMdO3bghhtuMOpnzpyJbdu2Wdvfvn07Zs6cadTNmjULa9euRTabRTqdxvbt27F48eISm9WrV5d9Xhs7duxANps1+jN27FhMmjQJ27Ztq9ngQ2UXRVEURbEwbtw4tLS0FDebF2Pfvn3I5/MlK7uPHj26ZAX4gM7OTqt9LpfDvn37+rQJ2iznvFJfGhsbcdxxx1XUTlLq3/OhKIqiKAXKnTDauw0AeP311w3ZpbfXg+m99IPv++JyEJJ97/o4bSY9b1yq1Y7E0Bp8WFzVRtIwQ1LxY9QLxzaQu7LgcXUp+ZRPicWMuRC0zoixqAjZcIyXw/Mesl6pLbl1zSRfJE2wq5TlGNYJksz/iCPBkOvf6Au3Y1nDxs+RrZRMjeeZsKTFNpTkja+5W5CxfJo24xiJ34SkYXxPDBuqNx60Xv9WE4syZUotVOZ6Q14Jr5UhrxhJxkhKKUgsLLWkjURgXJ+zlg0JBqVSCxC6YWnqiyG1cDIxU16prtTC8zxYajHkFSrz3I4s2dPUtKI1z+3IUL+P+OGnNtdqaRBswna6abLQkUK5myYcZT2SYDy7HMNJyViyyNPzbsouQQGDDg9VSK9eeI5GjBgROedj1KhRSKVSJV6CvXv3lnglAlpbW632DQ0NGDlyZJ82QZvlnFfqSyaTwf79+w3vx969ezFt2rTY7SRFZRdFURRFKZPGxka0tbWhvb3dqG9vbxf/eE+dOrXEfvPmzZg8eTLShUn6kk3QZjnntdHW1oZ0Om2009HRgRdeeKGmg4/68nwUcx3QCFyYUFocyPPwyrY/bpkdDnmLTcruGTEmiHI9NcgTR6XJp07BJpic2VNHE1jDo4zU7cbkT54UWg0vCMOTWflGGJNf7VE4xagZ9oawJ0OafEplw1OS43q6LwWPFE8shejhIBubVwPyJFLHdgkr+YUo/IgLPB7GJFQqu0KZJ5ymBY9I2lI2vR0Zq608sdQ+yTRFFyaYaMqr1LpCinQzdXrlruE43o4svfh5w8MRljOWSaZA6PGQJplmjMmkaWs5Q14Q9nbYJpyaHg7XXvbs5TxPjOeJqOwFKdj4wrtRMqm6H1Ov+3CMdPHltpGEJUuWYO7cuZg8eTKmTp2K+++/H3v27MH8+fMBAEuXLsWbb76Jhx56CEBPZMtdd92FJUuWYN68edi+fTvWrl1bjGIBgGuuuQZnn302brvtNlx66aX48Y9/jKeeegpbt26NfV4A+NOf/oQ9e/bgrbfeAgC88sorAHo8Hq2trWhpacEVV1yBa6+9FiNHjsTxxx+P6667Dh/5yEdw/vnnl3cBY1Bfgw9FURRF6YOBWNV29uzZeOedd3DzzTejo6MDkyZNwqZNm3DyyScD6PEkcO6N8ePHY9OmTVi8eDHuvvtujB07FnfeeSe+8IUvFG2mTZuGxx57DDfeeCO++c1v4oMf/CA2bNiAKVOmxD4vADzxxBP4yle+UvzvL37xiwCA5cuXY8WKFQCAO+64Aw0NDbjssstw+PBhnHfeeVi3bh1SLPNXmfrK8/HPtxTyfHBuD3uej0CyTnEdl+lnSUqspzLVu0b+D6+wn35tU1nK+REr/wfNgSjm+bDk/gB65fnIhmqzcXs57JXzB/tlej6EyUgOP7CcF4TzfLBHJl2a58Mxcn5IeT6i839wzo+gzHX5JiHnR4yymdujtN7ICcK2TfbcHl4TeXUoFwgot4fbSB6EQn6PRsrz0Uy5PYZTbo/hacrLQbk9jm3gPB9heURDmLsjyO/B+TyOdQ+Hx9GLZeT24JwfDi8+l7eWg5weSReNa4D9C7LcfB5JPR/ddCx7Pmw5PczcHuEzyzk83vOl3B6Ul4PyeLzHOT8K9ZzP4116EN/NhWXO7XE4Fz6Ihym3Rzfl9uBykN8jnyMvCeX2QJbzfLg9eT6+fmO/5Pn4XPtXkD6msjwf2Xcz2HjB92va36Od+vd8GMkUeDZezz9WKabPMrkWeR4ou+dZgilIKbzf5b/lnPNDyMshSjBAqb3UBpuSy9XIs2EkCLNPFi3qEHEGIb792vPAxjgP3QzHyEsSJEHgfB40mONkHDyA4YvOZbqhppTii3UA4NBE4aR5PqxSi0AS25JjLeM9rnLj5PmIyOchlbmOJ7NKk0wbYZdazHLY91TxX3s+D2mVWqbaAw7O52FILcKAg6WWjCGx9JSPCBNIjbJlAilgn1haatNg/AsA3XlKDkblvCC15IVVcI0JpyUF9Ku0otQ/9T/4UBRFUZQCAyG7KMnRwYeiKIoyZPCqMOG00uOVaOpz8CE9F5bIl3jyCkWkCJ58zu/Aa3YEbnte54NzPnCyG5fnWbD7Xoi4MCJYggQ0VCd57yUbo56DSawSjJCKXZBaEsMROcFnM+ae2GUUIw8KyzR0D8H5VAzZq3Acy1Vk65CtY0hEiCyLM/37qouLoNMEl9+NE+0Ce31KiILhvCBuoZwCyyU5qy1LM64Q1cJlc6UCp3Bc+c8VSykswdjWaIkjtUSt1dJTD6oPz9ltmd8hSi0xyjy3o5vkm8N5iogpyi4USROjnDVkF45w4ZxGti9Y4V4N+pmEykBTn4MPRVEURbGgskt9oIMPRVEUZcigg4/6oO4HH1KSsWLRWJnWXpYTi7H70Z6AymsouIpZUaHjXE4EJrgiWQbw6ZYYj79XenAsCUaScfhYlhgKbmYjAoZd1hVJLdQDukZFCYbb5mtiRAYJ0UOcat0ILy6VaYxU8HGkk4TloM1KolpELBJLPKmlVEbpbcOyCteH0S4ktfB+kia4DVNeEaQWKtviVKSoljjYpBYglFtiJRDj1OllSi095bTxL2CGznLEihF2S1JLnCiYILIlQ9pwRki1HifaxUwsVppe3ZfeB0WJoO4HH4qiKIoSoJ6P+kAHH4qiKMqQQQcf9UFdDj4M+cDIsEQ2hbIUNGEEUPAqqCxT8LFGwrHS9QwMKSZtXzfFETKeGQEkdCSf1LftJ2JJMBzhw/WciKyQeZRd+axn+Rb5p19Jmvwsar8ondgjX8x2orvS30jPhxThYsoxbCNEO8XcD5iySy3xYpzHFs1i1AlSS7lRLT1lklL8HvnkiCCjvEdSS7cgqbAEwxEunFAsqM/kKaqFyryGC0e45PIktUjrufgWCcbyXQj08c4oSoG6HHwoiqIoig0flefpGIS/K4YcOvhQFEVRhgwqu9QH9T/4kKJdopKM8dpnLJnw9Htp3QKLS54XCGPRxffZnWm3YQe2NLc/OLu49ovFFoiZiIyWow8kGJ8iGMDRIW60u908EX0it8wXWpJ6pMgXcQ0d89/eiJEvho09ZsiMwgkePnsblcBX0BbtAsv+3kgJxwwbI0FYIcmYEUljT0gmSS2pMkN/WBoxXuAYsg/DicOCNo1IFuq3tFZLUqnFXK+lRzJ5l6ST94yoFnuEy3u0KJwpr9glmCBx2BGq4wgXQ4LJC5EveSHaxbdILNJ35ACig4/6oPw4NkVRFEVRlDKof8+HoiiKohRQz0d9UF+DD4e2AqbUQpEqBRe/EyvaJSy7xhLSQuSLsXaLJQ5FdN/bHU1JJJg467bEsY+yMdqji+hLL6UkjZDU4nBokpvA6VauXBMDJ3H0jH3Nl+jzxLctOTbBx5dllOgOxIlgKRdemj0d0RdPDP8SXmbDwv4ZOJolaL/WUgsnCwvklqRSixzhQn3JU18KNiy1ZPP2Mke4yInFKCrPWHop/touji9/P9UCHXzUByq7KIqiKIrSr9SX50NRFEVR+sD3HdlDm6ANpbbU/+BDkmCCJGNG6AkVhQUmjOVChAgXawekxFVy2id7LUfecHSM5fxVk2AclkZK11nxKQmZIVPQeirGWi0SLLXYJJik8oqkR1Sy/owSC094fvMJHewUU4XSlHoQH+y8caQdlm/MJGKB7AKqC8uZKkstPTaNpcdxArEYUguXjxgRLmG/grVdunO0tosgteQSRrjw2i7F5GKDcG0XD07FeT4qPV6JRmUXRVEURVH6lfr3fCiKoihKAZ1wWh/U5+DDiS4Xnx0hL5EvJBMT64VyURoRH9Y4ggiH4Vhri3P4o5KQ9XTFsdZL9obEgsKy6RylwnIIr30jSTBx4MRmKTc4Ee0Pb4QRJZNUaomQYPw4Eo1gMxi/n6QvzThfpnkhgiTcb2/DOM6QRigqjMpZKqcsESeNxvUm6SRGZBI/hXlDdgkJJBa2PUIvfpY+T9RaLYC5XotNagFCKeW9vBDhklBq4TJLLEE0C6/bkk0otXhRES6AkGRMKPczOuejPlDZRVEURVGUfqU+PR8C5uTSQuppnjAlpU7nFWup2uUf9sJ0Tbe4n6mSF8RSa/RP+EUuTiyVcm44pfV+Lm/djzyv0uvZbeLA3pTA45Fy7ftTnAtfyBUSwztSzHoupSaQPGoxiPKgVPJDKnqR3mSN87Mse0roF7LlmWRvh0ceqyx5CtJOnupDG07vzrk20oWucM4NfmXjwB4Ofld4QmngkUnq7bClSwfklOlsH3g82MNxiLwglXg7uvNc7vlMmZw9twd7OwzPR569HTzZnBMi2SacDj4Pgcou9cGQGnwoiqIoRzcqu9QHOvhQFEVRhgx+FTwfOvioPXU5+JCeC8fiNjdshdwe4nnExVR5Ip1f0nQtJRjjPE5flkHT9smn4tkL9g5LGpTnwzdm7Ro556UeRBNMLjUmuQqSCkkzvivUU84RLhfbNyQaWMsso/hCvYjT698qYszpi/iCNJcBCMt5o94VbNwSG5ZUsn74TGRIvkg7ObIJ63lV26wwwTrIxZE2ZJlkGBNOeZKr75aUM+A6zucRncPDyNEhpEw3bXrkllpKLUAot+QSSi1Gbo+8kNvDNqHUSKMu2CqKhbocfCiKoiiKDR+V/RYK2lBqiw4+FEVRlCGDBweOZjgd9NT/4EOKXAg8isKKpHFi0n1BmjFXnnUKdaURML1tqyLBsIoQerjF1XBNqSW83Zy7w8gLEuR3z9ujXRyKYODwIl/6qWHkiBeic4L2paiWBl522I0uC5EvgUwjyjIx8sfEKvczRgpsqjejWkpllN42WSP6o7TMUS95Q7KwyzEp8fdj+ODy+dOFpzgvHMfSjZRzhNuTZZeez8OSSiaG7BJHaolanbaWUgsQyi25OFEtQup0M406qGzoj4V/ab+6C5QE1P/gQ1EURVEKaLRLfaCDD0VRFGXI4PmOOfm1zDaU2lJXgw/f6dkch6dYk7vQVm1MrA8NHMFVzUhShs3GZRd3TSWYvpOQ9T6PqUDYw2McI7Kj4DY3Ik949VpOOsTXU8g/H4fg/EJUC0sqPkswVDakFLZneakY7UJd5bIhRQk2MaSWqGRmFcHPreUEnhDhIiUWy3rhNfRcu30gU0iyTIqeOJYp4sCJyAIpRZRrhOq8ILWwvJRBad+zUjIxLzrapdsSyQJEr05bS6mlp9zzmeWoFq63R7hYk4n1LhfuhRjh0rteJRmlF3U1+FAURVGUvvD9KkS76GCp5ujgQ1EURRky6JyP+qDuBx/8jNjyPwkLbvaSRiqXYGwRMKXnsdvXUoIR11xx7TZFCYZWqfVZdqF6Jx8jyZhUb+uXIJcYkSwstUjRLg0s01B9QZqxSjEwnxVjnSAxEVm0TZ91cRG+CINL6wlSjCTB5ITEYoasQnJM1ukps9TAUotLz4frm0+2DTMKhdeCyZXYShJMXrigeSFxGJ8nU6bswjKKIbtESC1crqXUAoRyS0VSi1QfkWRMpRUlCXU/+FAURVGUAPV81AdiVm5FURRFqTeCVW0r3ZKyZs0ajB8/Hs3NzWhra8PTTz/dp/2WLVvQ1taG5uZmnHrqqbj33ntLbB5//HGcccYZaGpqwhlnnIGNGzcmPq/v+1ixYgXGjh2LYcOG4dxzz8WLL75o2Jx77rlwHMfYvvjFLya+BkmoT8+HmTnLWu/bTGsgwQTqgSSpVF+CkfaT1MCKBUUGGcEuOSkRV087DiUZM9ZZYTkm5VvrzW7H8MUmiXYx1mqxyyteg12+Ccpegz3qJ6kcI5Yjvrcq+VFlBhUVPg+vyUIudi7nhDIn38qRTNHtkRzi5gu29kRcLneKFTqH+kLltBO20+iUrhGTcmJIN779d1OWolryRmIxToRmSzJGUocQyWJbq6XHPjw2KolYLaUWAMgX7D0haVhVpBYATmAjqaoDKMEMxITTDRs2YNGiRVizZg2mT5+O++67DxdeeCFeeuklnHTSSSX2u3fvxkUXXYR58+bh4Ycfxi9+8QssWLAA73//+/GFL3wBALB9+3bMnj0b//RP/4TPfe5z2LhxIy677DJs3boVU6ZMiX3e22+/HatWrcK6detw2mmn4ZZbbsEFF1yAV155Bccee2yxT/PmzcPNN99c/O9hw4YlvWyJUM+HoiiKolg4cOCAsXV3d1vtVq1ahSuuuAJXXnklJk6ciNWrV2PcuHG45557rPb33nsvTjrpJKxevRoTJ07ElVdeia9+9av49re/XbRZvXo1LrjgAixduhSnn346li5divPOOw+rV6+OfV7f97F69WosW7YMn//85zFp0iQ8+OCDeO+99/DII48YfRo+fDhaW1uLW0tLS4VXr2908KEoiqIMGXo8H06FW09b48aNQ0tLS3FbuXJlyfkymQx27NiBmTNnGvUzZ87Etm3brH3cvn17if2sWbPw3HPPIZvN9mkTtBnnvLt370ZnZ6dh09TUhHPOOaekb+vXr8eoUaPw4Q9/GNdddx0OHjxo7Xu1qE/ZhSH/nk+ShFOsQ0kdEE+CcaRIEUs1eV5rLMEIn8hQUewtulTvGeu5UDlXsGd5I8dRLdR2PkaSsSTRLiz/GGuu2BOIQZJaOOFYmmyCaBfazxIMr+XjS3JMksiXGAnJYiG4voNLa0a4gMokxwjySsZjdz89K3SiBklSs8DPledmw/NLUSgkuxTlFl4OCPZze7zODH02KZKG14I54jeW9IOlFpZXujkixrNLKofzLMGQNGWRXWoptQCh3MLyih9HXkkitXC9lFhsAKnmhNPXX38dI0aMKNY3NTWV2O7btw/5fB6jR4826kePHo3Ozk5r+52dnVb7XC6Hffv2YcyYMaJN0Gac8wb/2mxee+214n//zd/8DcaPH4/W1la88MILWLp0KZ5//nm0t7db+18N6n/woSiKoig1YMSIEcbgoy96/1j1fV/+ASvY966P02Y1bObNm1csT5o0CRMmTMDkyZPx61//GmeddZb4GSpBZRdFURRlyOBXaYvLqFGjkEqlSrwce/fuLfE4BLS2tlrtGxoaMHLkyD5tgjbjnLe1tRUAEvUNAM466yyk02n87ne/E20qpS49H0ZiMXMPlZyS/UklGF88om/6TYIR3PpGGyS1+CSTcBSMa5EYHJZUWJZhv34DlfPC65ok2oUPE2UXKktRLRapBQC8Qr0hy6SEtqWEY4IcY7sXsdaE4bAA4wCu5+rSRo1gE4524cRiMSJfup3wwTWej3wS2YUlEI5wIXmFnmxrYrEYoRJ5wa1uSC1iwrHSpGlGJI8gr7C9FOFyhKQUQ9IqyC21lFqAUG6pqdTCZboPjtBGf9PfeT4aGxvR1taG9vZ2fO5znyvWt7e349JLL7UeM3XqVPzkJz8x6jZv3ozJkycjnU4Xbdrb27F48WLDZtq0abHPG0gp7e3tOPPMMwH0zBXZsmULbrvtNvEzvfjii8hmsxgzZkzs65CUuhx8KIqiKMpgYcmSJZg7dy4mT56MqVOn4v7778eePXswf/58AMDSpUvx5ptv4qGHHgIAzJ8/H3fddReWLFmCefPmYfv27Vi7di0effTRYpvXXHMNzj77bNx222249NJL8eMf/xhPPfUUtm7dGvu8juNg0aJFuPXWWzFhwgRMmDABt956K4YPH445c+YAAH7/+99j/fr1uOiiizBq1Ci89NJLuPbaa3HmmWdi+vTpNbtmOvhQFEVRhg7VWEU34fGzZ8/GO++8g5tvvhkdHR2YNGkSNm3ahJNPPhkA0NHRgT179hTtx48fj02bNmHx4sW4++67MXbsWNx5553FHB8AMG3aNDz22GO48cYb8c1vfhMf/OAHsWHDhmKOjzjnBYDrr78ehw8fxoIFC7B//35MmTIFmzdvLub4aGxsxH/+53/iO9/5Dg4dOoRx48bh05/+NJYvX45UimbhVxnH95OnY1mzZg2+9a1voaOjAx/+8IexevVqzJgxQ7Tv7u7GzTffjIcffhidnZ048cQTsWzZMnz1q1+Ndb4DBw6gpaUFJ628BW5zc68PIByUxC3IbkbD5Ug2HtdTm/lSW/IwG/Vuzm5j1vv2+kLZNY4TbI027GWjjxTNEvTXcLWzBOPZ5QBHcM0nSTZkeDrFhF+ciMxeb5NagFBi8dIc7QIq83GCDa0Yn290rPVB2VYHAF6jT/W+1cZvpOuZpvtD5VSh3JAOb2Y6HT4IzVxuCMtNVB7WEEakNKfs5WGFchM9ZMPYlqJaWF6JU3bpBUkJkS028sJ0NZZPPL9v2eWIILuwpGJGA9mllowQ4cL22YKUUkupBSC5JaHUEimvlJQL8mzU929w2iNHsOeGG9HV1RV7AmdSgr8Tp65bBnd4c/QBfeC9dwR/uPyfa9rfo53Eno+kmdwA4LLLLsPbb7+NtWvX4n/9r/+FvXv3Ipcr1XoVRVEUpRIGIsOpkpzEgw/OqAb0ZGH72c9+hnvuuceagOXJJ5/Eli1b8Ic//AHHH388AOCUU06prNeEMfnUOj/Ut9oaHhFOY8GNV3kiqi94sJJMRPWM03H+CZ5Map8Iyd4BwwtCM06LqRbYk0L9Zs8H5/bwPZ6JyTaIjTkpU8r5Qfbs7Wiw29tyehj5PAzbsN4zbLi+70mmPf2y7a/g20yYuepbJpzyRDn+BZ2le+zSr2xOZc65PVxLfz1jNVyatEpts3ckTW46bjtty+0B0wsSBZ+f83x4dNHzEZ4P9nDkPLvng/N2sBckk7d7O3g1YM7pkS28H0m9HebqtBHeDiD0eFTi7eDbIEy8tD7O+gdbSUCiUNtyMrk98cQTmDx5Mm6//XZ84AMfwGmnnYbrrrsOhw8fFs/T3d1dktZWURRFUaKoPLtp5dEySjSJPB/lZHL7wx/+gK1bt6K5uRkbN27Evn37sGDBAvzpT3/CAw88YD1m5cqVuOmmm5J0TVEURVF6vDWVDh508FFzyop2SZLJzfM8OI6D9evXFxeqWbVqFf76r/8ad999t3XlvKVLl2LJkiXF/z5w4ADGjRsX2S/b8yI/QiQZGHksbBaIIcEIZ+IcDPbmDMwVc0uzlIgZ3w2pxf7ZzHwe1BeefBpMoOVVZY3JtjxpVZBXpIzqgpDqW1Ot036WRgTZhS+oJMEEHnFTipGOo3qWy/gaCunYrX2V8nxIRKRU5zK749llz9IJu/iNlY7z0V8BgdRkyi4sD4ZlljfSHk8stcsuNnlHWtVWWsnWSCNPN8hcvZckjkKZJZIctd1t5OqIll24nSxPKGWbwn3h/flYUgtNmBZWqrVOLq2S1BJ3QqmiJCXR4KOcTG5jxozBBz7wAWOFvIkTJ8L3fbzxxhuYMGFCyTFNTU3WHPqKoiiK0hc64bQ+SDTngzOqMe3t7cWsa72ZPn063nrrLRw6dKhY99vf/hau6+LEE08so8uKoiiKItDf+dWVskgsuyTN5DZnzhz80z/9E77yla/gpptuwr59+/AP//AP+OpXv2qVXPrEQfzgksLDI0l3vYSjsOTaXY6+oJkE1kkjYOJIMKZ9Ia5eiKAwFrI1VqwNq32XU6pzjgw6tOAVNnKZGFIL1XN69Rh5U6KQU5CTjZgC3S6ZmPk/SvcbuT2ECBebdFN6fks5qdQiYcmvAISud4408uhe5TmNOktkFL4Up1vBs8fyRs6IJOE06hy9QlKLEPniCvaRfTJW8pWiXai/LI0U7HOC7JIRZBdDpvHY3i67cARLYJ9UavE4qoXllag8HpLUYsgrXBaklqhbEveW6R9zpReJBx9JM7m9733vQ3t7O66++mpMnjwZI0eOxGWXXYZbbrmlep9CURRFUYCqRKtotEvtKWvC6YIFC7BgwQLrvnXr1pXUnX766SVSjaIoiqLUBPW0DHqG7toutoErqwRCrihDgqE9LHfYomCkJGQOu8mFrpoRLvZyUT0wPKzcP/azUpEjWSKkFiCM9jHSxUtl335RHGFF1kiM6BAhsZgtmRf6kmOoHERtiPJKWPaE5GOIklqov7Y6IN4lMZ4bIyqh9NryrzSOfOFyjvQ3R1iB2PC2cwRJIFNQuFSDGz4IGXqAuL5BkFSMlOpcn+BhkeQVKSKHZZVAAsn5dhmF5SqWXbg+a0hafUstQCixSFKLcQ8lqUVKHGaJYImXLj2h1KJ/0JUqMnQHH4qiKMpRh8ou9YEOPhRFUZShQzWiVdTLU3OOrsFHac6unmIcCSYqEZmQhMxcnyVagpEITinLMuSG5gRiFOFiSCbcX6/0s5mr19o/u+SqdaRfDVLwvCXJmBT5ItYnkGDkNVyEiJkYa77Y13YR1hUSykZ0lej6tkW7hMaew65+oQ1uzrJWDAB4qfDgQL6QZJeUILu4gqQSFeEiRb14wnMVJ/KFJZDAPm+RYgBZXjGjh+z1eYvUAoQSmLhWi7BuC6KSiQH2yJakCcSGlNSSJCyyrzaUWpIoz4eiKIqiKEqlHF2eD0VRFGVoo7JLXXD0Dj5qIMGE++1txJJgJG9fUJ+31AHgFe05wsXj9U9YMkkJckxRdpHculROvJ5LfFdmnGgXUXZhG4vswpEsxpowQlSLGB3j2OuDYxOv5yIhXedAdskbGlq4O4Zj05BaOMKF14tJ9bSZcoUoFSqnJNlFkGCcJFnohH4zLLWYn8cmu9j3m+W+2wDs8kppucdGkld8KXolhtRiLR91Uguhg4+6QGUXRVEURVH6laPX86EoiqIMPXwH4roaSdpQasrQGnwkeV4kNSChBBNEaojNxZBgHEFiMMp5S5+k5F92L7wpu1giXHpOENRxorTS/X3VSzqS5GGPfM9jRL6ICb1skS+CXCInJxPOKR5b+KBCn2I9pxFSC7djXD9y05u3lRNahfUcDeUZ6+CURooYUoshr4TtseziCLJL1MdPGu3iCzaSpBTUs6RiJlUrlUt6ypKkwvKJPXFYsDYUS2S+tOZKRAKx3mWrRFolqaXstZkGEF3Vtj5Q2UVRFEVRlH5laHk+FEVRlKMbnXBaF+jgozeCBCMRuBqlw+JIML50hEViMCQSQXYQI1LENVpKbeR1W2Ctjxft0je+JdlYTyPCAXGiXWw20nERa7UAMWUaS7SLsKSGjJC0zXhWipnn+DDuCNmmKAqGPoSRoIwkE5YSgtuSNySVFJWFqBbBZy9FuCSJfJGiXUSpxWIjrYkj1gs2NnmlxD4oSzJKNaSW3uUohuofWJ3zUReo7KIoiqIoSr+ing9FURRlyOD4ySbKSm0otWVoDT7Ky2cVr+mIyBefozpidEly6xvrgRhBNQUj9qoLMookwYAkAzGCxbPUxZghL9uUeSOEw5Ku+WKz8QV/nxwxQ/W2qJZeNsWydI+TIq3TEYSZ5O22vlEOO+O4xoNF9SQ38PNeuLmOEUbF0Vqltn3VM9VOMmbKIX0fa1wfjlIxjO3t+YLsISYL8/uo63WeiqSWQjuxIlyqwKBUJ3TOR10wtAYfiqIoytGNzvmoC4bu4CMYuVbyDMWZfOqUGsTxghhNSL/gDbdJwcPCv4SMuYX2uH5fmFjqS7+u/NK6OMTxjiRrUGguzv2M8ojE8JiYXhC6t2wfkdJd9sYI7Ym/hMmbwJMYgwM40YbgBeFf3D7fLOMz2L0Wxf5KngyjHnYkz4dgHoX4WIkekb4b8UXPg+BJKdeDEcfDEcNGOmekI6lKv+j1b7NSKUN38KEoiqIcfajsUhfo4ENRFEUZOujgoy4Y+oOPGjxExuRTq7xjl2BgOw7yxDjHJgMIx3E+jUhJpTdRbtsqSSpRLuHErtw4+T9sNpLsUclkVtv547QhEOv5CFz/bGscKM2SFh5EeuDM/pZ2INbnSTiZVJRsgnMmfZelByrqGY9zXCX2BRLLK0bbZb6zCVF5RakVQ3/woSiKohw9qOejLtDBh6IoijJ00GiXukAHHxVie0ZNDze7wYXjynUPc7XkmpcOlf4jYsQ/WJPviN8VUd8hcS59FeSdSqKujAiXCDlEXCHZMAqL8a5b1ENeHQb00UooLTpJ3lmpPoakYp4zQduKMsjRwYeiKIoyZNAMp/WBDj4URVGUoYPO+agLdGG5GhBIjiWba9/g+tbNT1m2hmSbJ21paUOfW75xcG5yn/v6rPL1iXV96b5Y7x/fZ+GZkB8i++Z4vDmlW562nH0DbU7WtW+ZqM0ZYpvwOYXrg17XsXg9pS1v2Wz3z3N63eNwk54JZXCwZs0ajB8/Hs3NzWhra8PTTz/dp/2WLVvQ1taG5uZmnHrqqbj33ntLbB5//HGcccYZaGpqwhlnnIGNGzcmPq/v+1ixYgXGjh2LYcOG4dxzz8WLL75o2HR3d+Pqq6/GqFGjcMwxx+CSSy7BG2+8UcZViI8OPhRFURSlAjZs2IBFixZh2bJl2LlzJ2bMmIELL7wQe/bssdrv3r0bF110EWbMmIGdO3fiG9/4BhYuXIjHH3+8aLN9+3bMnj0bc+fOxfPPP4+5c+fisssuwzPPPJPovLfffjtWrVqFu+66C88++yxaW1txwQUX4ODBg0WbRYsWYePGjXjsscewdetWHDp0CBdffDHyeU6ZXF0c308cPd/vHDhwAC0tLTjpX26B29w8MJ2o4VWqir5Yrf4NtVne1bi4ZV6SquUtqZZ9garc4qEmilfpuR9U73J/t90H3pEj2LP0RnR1dWHEiBE1OUfwd+Lk2yr/O+EdOYLXvn4jXn/9daO/TU1NaGpqKrGfMmUKzjrrLNxzzz3FuokTJ+Kzn/0sVq5cWWL/9a9/HU888QRefvnlYt38+fPx/PPPY/v27QCA2bNn48CBA/jpT39atPnUpz6F4447Do8++mis8/q+j7Fjx2LRokX4+te/DqDHyzF69GjcdtttuOqqq9DV1YX3v//9+MEPfoDZs2cDAN566y2MGzcOmzZtwqxZs8q6hlHU15yPgXQz1vBvsrTKaqI2qtU/YxncIUAVLkzZf1ASJ8VKZl7urarOH8gh5jSt1mCqCu9hTcf/A/Xboj/PW8VQ23HjxhnVy5cvx4oVK4y6TCaDHTt24IYbbjDqZ86ciW3btlmb3759O2bOnGnUzZo1C2vXrkU2m0U6ncb27duxePHiEpvVq1fHPu/u3bvR2dlpnKupqQnnnHMOtm3bhquuugo7duxANps1bMaOHYtJkyZh27ZtOvhQFEVRlP7E5vnozb59+5DP5zF69GijfvTo0ejs7LS229nZabXP5XLYt28fxowZI9oEbcY5b/Cvzea1114r2jQ2NuK4446L3f9qoIMPRVEUZehQxWiXESNGxJaJnF7rA/i+X1IXZd+7Pk6b1bLpTRybStDBRy2Q8g8lXQPDmrgq6TkTtN1Xfb0SNaUpRhIpsYmoNTjiLHUunj/OOUtvlhPVp1714u2O+vJO2u8BJtILL33JJswlluS9kvrkJHxnh9o0rYrp51DbUaNGIZVKlXgJ9u7dW+JxCGhtbbXaNzQ0YOTIkX3aBG3GOW9rayuAHu/GmDFjRJtMJoP9+/cb3o+9e/di2rRp8S5CGQwx4VZRFEVR+o/Gxka0tbWhvb3dqG9vbxf/eE+dOrXEfvPmzZg8eTLS6XSfNkGbcc47fvx4tLa2GjaZTAZbtmwp2rS1tSGdThs2HR0deOGFF2o6+FDPh6IoijJkGIgMp0uWLMHcuXMxefJkTJ06Fffffz/27NmD+fPnAwCWLl2KN998Ew899BCAnsiWu+66C0uWLMG8efOwfft2rF27thjFAgDXXHMNzj77bNx222249NJL8eMf/xhPPfUUtm7dGvu8juNg0aJFuPXWWzFhwgRMmDABt956K4YPH445c+YAAFpaWnDFFVfg2muvxciRI3H88cfjuuuuw0c+8hGcf/75lVzGPtHBR6VEuVPjuFDFpdxt63D3fe6Sc1Zin+S4GAzIkukRcohkm1RqsS6PTo1we7xWSyw8Qb7xe/2LXrdK7GtSm8rlnThEPU5VC10u851lG0d8Z4XzJFljiBo314Cy99HW9FEtxQxAhtPZs2fjnXfewc0334yOjg5MmjQJmzZtwsknnwygx5PAuTfGjx+PTZs2YfHixbj77rsxduxY3HnnnfjCF75QtJk2bRoee+wx3HjjjfjmN7+JD37wg9iwYQOmTJkS+7wAcP311+Pw4cNYsGAB9u/fjylTpmDz5s049thjizZ33HEHGhoacNlll+Hw4cM477zzsG7dOqRSqaRXLjb1ledj5QDm+ajyF5kOPkoZTIOPpG3bBx+CrRd9fuPS6uBD6kZEg/Hrkw4+5Hc22qZoK/VJas+w7/udLXuhxRrhHTmCPTf0T56PU27556rk+Xj1xmU17e/Rjno+FEVRlKHDAHg+lOQM+cGH9GuqIrdkkl9OST0cor0fb3/vepfcszHsHdt5+LCEXpByL7P47gs3zhd+iftWj4QQycK2nm+vFzwSPu8oJP+S5BXDfR8nIiXKO+EJtlI9R+FworKoc8bxnghUfTn4akR+iB4LJ9rGtddzpIrNC8K52aR3yZDopO8M4wKUvsuO8PgmRjp2EP9x1lVt6wONdlEURVEUpV8Z8p4PRVEU5SiiiunVldoxZAcfVXGbxXj+is9oDOnEd6NtIuUTSVIR2nZiyC5GteuV2EpJj+LZlHcjfFFesUsmUr0psTh9tyFMBDXqucxdZPnCLZyHGjEkGMfuVzclEHtfbLIKH2eWnRLb0vMI9SitTz5RVaiH3SYRcWSXKBvj3bQf5wu3W5JPfJfup02a8e3vaSw5xt5FuwQTRy6phhwzGOUJnfNRFwzZwYeiKIpy9KFzPuqDITX4qOkDI05O63siaCxvB9tI9oV6R/B2OMYvJ99eJnvXtdsE9WIb3CXBhkni+Yjj7WA8qvcsHo6eMtl7rqWObL3wIhpt0ORTP0820v3MF/51+UGwfoR4SJNIA49E3u49Mcr5sCx5O6Ry8cd0jMmpko2BUO8kiLX246RDjwiBjTOB1PCIsI1nt3GE+uBYI22/4bDwS2xL2mZzSPWlLgljf7XVBKlTihLBkBp8KIqiKEc5KrvUBTr4UBRFUYYOVZBddPBRe3Tw0Zukbknr5LUKpBZ2s7I0kipMBJXkFTf097qCvGKWQ/uURT5JSbKMILW4gjTDuMI3ghfhCzZUB0FeyZNkwik6PKrPF3zlxnEko3h0IMsxniP41Unu8A1febCfqozrFv2QifOXLRKHJK8YZV+oTyDTmLZ+yf4SmxgTVG3Ey80TY4J3RD1PDpUmf4qyS4yy7br4nK3aSHRgzxkjTnK1nyYUXeJ8j1Vr8qmilIEOPhRFUZShg8oudYEOPhRFUZShgw4+6oIhNfgoeya3dFxUhAvZJJZaUnbZxZRaSuvdGPJKKkU21NcUSy2uVN9TTgmSCttyvSu8rdXI8+HRRZQiXEzZheotNjmSWvgzGG14LM2E/fLIV26kxbCtYGtIDULCCAlJprDk8RBlF6HeZUnFi1Pvl7QnpW5PlK4dySJcmDjRLmLUiuWdlWUUakSySQn1Nhv+uCzBGNeEn/0yJRhhOVwx8kUOn1GUmjCkBh+KoijK0Y3m+agPdG0XRVEURVH6laPX81GJ1GJLxVyJ1GKRV4BeskpBSpHkFZt0AgANgk1asi/KLnZ5pYFlF/LJuoI9kyTaxZRXWOqwyy45I6rFtdcXylmX9udTVA77l3Oix+VGhAtdZz84v5SKPSmSfFFMMmbfL0ktYkQMfX6bvZyQzLfWl512XZJiHLseEBXVAkiyi2Pdb8oldmnGS9kTu3kspdhkFa6DZb9ZhItkEoytT4Z0VIkEoyhV5OgdfCiKoihDD51wWhfo4ENRFEUZMuicj/rg6Bp8VFtqAc2Yr0BqcVl24QgSlkwKZZZaGgR5pSEV+swleaWR/Oo22aWB/PHGeQR5hWWXVIwoGBssqeQl2YWlFkFeyVH4AddnCz5xl6QWM2laWM8ShISxRgydM5DOfEfwXwvPWyxpwre4+5NKLTm28e02lnb4miSVYGJFu0RecntEmRH5EkN2CWQIOdqF1/Wxt8Gfk6UWV1oxuWAjLBwsklSCCb6njCfvaJRgdPAw6NEJp4qiKIqi9CtHl+dDURRFGdronI+6YOgOPqJchAkSiJXY2xKKiWu1REstRuIwi9TCZUlqSZPU0piySyZp8sOz7NKYCv3wgazSQPvTDssyJNcI8kqcyBcbcSJcsuTjZvss2We88LFm2SWo50geJx/aSo+Mb7jSwzJHHvm81kkgAxhSnT06wxHOKi5jErHUvZgIzJBR6Jk0JBgq50rlGEPG8Ur3l9THWAsmMvJFIM47y3KMLRFYZEIwAB5LFim7HGNEqpDWYkRxFapdw1Y4J+wYCgh/Npt8YllqCDAvlS9kE4slwdiMBwk656M+UNlFURRFUZR+Zeh6PhRFUZSjD5Vd6oL6Gnw4SO7mE+zLjmrp3aZr8ady9EoVpBYglFXSVMfyStqIcLFLMI1ujsphfZNFduE2OPIlLUTBpBAttaQcu0M5b8mSZEoqJLW49nqOcOnOh+fpdsJHPOhXnAgcQxkQ1pPhvnh0LZwgmsKQAyrwT0syhWf+27vsCgnEpIgYllroUSkeK7XHx0GQWszkZ31HuEjrvSRez4Vlkjyv0VKIRnLtsoxnfAa6x/R5fCOBmCTH+CVFL8YXGL8N8aJjSiUYWwRM7/45xvdbtATDFD9y3D/S/fjHXGWX+kBlF0VRFEVR+pWyBh9r1qzB+PHj0dzcjLa2Njz99NOxjvvFL36BhoYGfOxjHyvntIqiKIrSN36VNqWmJJZdNmzYgEWLFmHNmjWYPn067rvvPlx44YV46aWXcNJJJ4nHdXV14Utf+hLOO+88vP322xV1WsTiI7QsG1KwrYLUwmVhfRZzrZZoqUVKFtbYUJBdqI2mhtA3zjJJcyobHifILk0su1B9EM3CdWlBdjHWf6G31VwXJn5aJY5wYSnGkF2EyJduinAxI29Ko3DiROAYUgtdwzy74VPskqdyoegIckAs6VCIArF13ZBdSBqBKMFQvSC1WKNdjDp72bARIlxM2cUuzURhXk9h2fscr91Suo6L3xAa+3x9WIKhb0hjDRdet8eQV+w3NzSx20pviSjBiM9QQfKLSEJm9iReIrKo04vftf2NzvmoCxJ7PlatWoUrrrgCV155JSZOnIjVq1dj3LhxuOeee/o87qqrrsKcOXMwderUyHN0d3fjwIEDxqYoiqIoytAg0eAjk8lgx44dmDlzplE/c+ZMbNu2TTzu+9//Pn7/+99j+fLlsc6zcuVKtLS0FLdx48Yl6aaiKIpylBJMOK10U2pLItll3759yOfzGD16tFE/evRodHZ2Wo/53e9+hxtuuAFPP/00GhrinW7p0qVYsmRJ8b8PHDggD0CiolliJA0rW2oBinKLIyUTY9mlTKkFAJqK0S72KBVDaiEf+zCqZ3tTgimVWJrdbEld7zJLKmnD3x+SSrCaRZ7Gwp4guxhl8o8b0TkUruDmWQIqff4M6cRYW8a+nkye7hUnMDOWFwmeOVEvSeiftiSrAmBd24XL0rot4povkn3WL9nvZvkaC1KLUG/YCBExUUgufj9lvyeGlFC450YCMZYjSI4x1nBpEO6DEQ1lv89BM6aMQtcT0RJMsoCpiCRkvTojSTB8UvOaB884mbKaNpASjMoudUFZobZOr7fA9/2SOgDI5/OYM2cObrrpJpx22mmx229qakJTU1M5XVMURVGOZnTwURckGnyMGjUKqVSqxMuxd+/eEm8IABw8eBDPPfccdu7cia997WsAAM/z4Ps+GhoasHnzZnzyk5+soPsm1twdST0cCbwdQOjxkLwdKfJelOvt6Cn3eCd4kil7L5rJq2F4OwwPR9ZqY/NyNMXwfLBXI815uolUjJ+zecvPpKwfPprsEWHPR7eXttZHYa6SK6yY67rWMq8AzJ8tz3k+nIJ9LX79WZwp8gqzQpk9GDFyfhTzfEjeDqo3PBw5e70xydS2erCQ54Ox/dgBAD/FH1SYcFrwjnDeHRiTh8NqL03p9G1eJwD0GML8kil1BbBTIakXxLBP9M2dzAvCSP66MDW78D2qKBEkeoQbGxvR1taG9vZ2fO5znyvWt7e349JLLy2xHzFiBP77v//bqFuzZg3+7//9v/jhD3+I8ePHl9ltRVEURSlFk4zVB4lllyVLlmDu3LmYPHkypk6divvvvx979uzB/PnzAfTM13jzzTfx0EMPwXVdTJo0yTj+hBNOQHNzc0m9oiiKolSMyi51QeLBx+zZs/HOO+/g5ptvRkdHByZNmoRNmzbh5JNPBgB0dHRgz549Ve+ogZg/weICrKHUAoRySy2lFiCUW3hiabMhr/Ak0wyV7fJKE8kkRn2hzPJKsxO2x1IDSy0pKc+HMH3Os/h5zdwelFtDkF1Y9jniG77vPhGlFqo3Vg8WpBbXkFpI9inYOGQr5vkwZunZZwCK81Z9i60xmZPrpbJdArFNLhWlFiPPBzVuyC48u1OQXWxyiyTBmDN86ZylE0t7TEh6aOhpk1Ou8yRTl6+Pke+FymmyN/orlMMWwzaodjBJMPwYShk/QtElOi37oMn/oQwqyspwumDBArz66qvo7u7Gjh07cPbZZxf3rVu3Dj//+c/FY1esWIFdu3aVc1pFURRF6ZPBHmq7f/9+zJ07t5hKYu7cufjzn//c5zG+72PFihUYO3Yshg0bhnPPPRcvvviiYdPd3Y2rr74ao0aNwjHHHINLLrkEb7zxRuJz79mzB5/5zGdwzDHHYNSoUVi4cCEymfAH6KuvvgrHcUq2J598MtF10LVdFEVRlKGDX6WtRsyZMwe7du3Ck08+iSeffBK7du3C3Llz+zzm9ttvx6pVq3DXXXfh2WefRWtrKy644AIcPHiwaLNo0SJs3LgRjz32GLZu3YpDhw7h4osvRj4fepGjzp3P5/HpT38a7777LrZu3YrHHnsMjz/+OK699tqSPj311FPo6OgobkmDR+prVVsbEXk8EkstPANeyN1hpkzvKddSagFCiYWllmFCmSNchrvhiJXlFa7nyJbGgpTS7EjRLiS1cJ4PCHk+YqRXD+QWzwnHwhkjt0f4mBqr53K6+hgREgEcvWKkayd/c4ZklBSVzXTtYZsssTgWya+iSADJlR9Euwj7jVfDyLNB9RERLkAot4hSSzZshCNcOGzEkFdYgslb5Jg491KIdnHo3sJYyZbqg/OnwvvK1yeQZQBTonN9UxwJuy0kzBBFi9L9NZVgJOXKSO9engRjj4CBztgUePnll/Hkk0/il7/8JaZMmQIA+N73voepU6filVdewYc+9KGSY3zfx+rVq7Fs2TJ8/vOfBwA8+OCDGD16NB555BFcddVV6Orqwtq1a/GDH/wA559/PgDg4Ycfxrhx4/DUU09h1qxZsc69efNmvPTSS3j99dcxduxYAMC//uu/4vLLL8c///M/Y8SIEcV+jRw5Eq2trWVfC/V8KIqiKEOHKno+ei/z0d3dXVHXtm/fjpaWluIffwD4xCc+gZaWFjFL+O7du9HZ2WlkFm9qasI555xTPGbHjh3IZrOGzdixYzFp0qSiTZxzb9++HZMmTSoOPABg1qxZxSkWzCWXXIITTjgB06dPxw9/+MPE10IHH4qiKMqQwanSBgDjxo0zlvpYuXJlRX3r7OzECSecUFJ/wgkniFnCg/q+Mot3dnaisbERxx13XJ82Uefu7OwsOc9xxx2HxsbGos373vc+rFq1Cj/84Q+xadMmnHfeeZg9ezYefvjhyM/P1KXsIq9USzaBNFIlqcU1koiVrk5bS6mFy5LUwjLK8JRdUjEkGMceBRPUs7zSKMguUsIxJiX4fPOGi7bnH5ZX+Dwc+cJRLa7PLn4qkkuepZwwgRq3TSnaSY8Qo12EsmvILkG0C8omjte67CRjYqrzsGhbnZZtOarFiRPhkiN9h+UY4/yWZyhOtItQz6nWWXZxCmnxWQriaBc+p8tyjCFplaZOLxxh75ctKZd1f3UkGM65Jwmf3BOzvfgSTCyVqY6jXV5//XVDapAyb69YsQI33XRTn209++yzAGBNkCdlCWfiZhbvyybOuaNsRo0ahcWLFxf3TZ48Gfv378ftt9+Ov/3bv+2zP0xdDj4URVEUxUo1JowWjh8xYoQx+JD42te+hi9+8Yt92pxyyin4zW9+g7fffrtk3x//+EdrlnAAxXkVnZ2dGDNmTLGeM4u3trYik8lg//79hvdj7969mDZtWtEm6tytra145plnjP379+9HNpsV+wf0yDf/9m//Ju63obKLoiiKMmQYiFDbUaNG4fTTT+9za25uxtSpU9HV1YVf/epXxWOfeeYZdHV1FQcJvRk/fjxaW1vR3t5erMtkMtiyZUvxmLa2NqTTacOmo6MDL7zwQtEmzrmnTp2KF154AR0dHUWbzZs3o6mpCW1tbeLn37lzpzEwikP9ez7E5E2WuipLLUAot9RSagFCiSWO1CJFshzjhpOlksguLK80UlQLr2prJBlL+LMjkGD4PGYysej2WF7x6PPYEpTxCrhpjz9PeJ4Gxy6p2OSV0npLByuZ/S8d6pfulyJfrKvhlpQFOSaIdhHkFVNGsUstjiC7lB3twggJx4y+kOziB6sRsxTDOgWfX9QshCRjBrY9UREwpk0lEoz9ODuVSDDF/dSIrmprZ+LEifjUpz6FefPm4b777gMA/N3f/R0uvvhiI9Ll9NNPx8qVK/G5z30OjuNg0aJFuPXWWzFhwgRMmDABt956K4YPH445c+YAAFpaWnDFFVfg2muvxciRI3H88cfjuuuuw0c+8pFi9Eucc8+cORNnnHEG5s6di29961v405/+hOuuuw7z5s0reoAefPBBpNNpnHnmmXBdFz/5yU9w55134rbbbkt0Lep/8KEoiqIodcL69euxcOHCYmTKJZdcgrvuusuweeWVV9DV1VX87+uvvx6HDx/GggULsH//fkyZMgWbN2/GscceW7S544470NDQgMsuuwyHDx/Geeedh3Xr1iFFIeVR506lUvg//+f/YMGCBZg+fTqGDRuGOXPm4Nvf/rbRv1tuuQWvvfYaUqkUTjvtNDzwwAOJ5nsAgOP7SX9m9D8HDhxAS0sLTvqXW+A2N5ujbiMdeli0TjhVz0fY9iD0fPAk1KyQ54MnnBplz17/rhdOEDuUbwYAvOc1hnW5cP+B3LDwuJzd5t1sWP9els6ZKS1nM7Qyb4Yezm7KL5GhVVO7afXc8BbC7Sab8FYh1W3+23McPbMZe71UdjN0P6kc5PdwM/mSOgBw8kKejwH2fBiJWNjzkSr1fKAhvCc8+dSnL26fUqobNlTvUb1n1BdWtaVVcjknR7C/tBza+IY91adK7Q1HDk8+5XrXbsPfo15K+n4t/Gt8v/J+37D1jhzBnqU3oqurK9YcinII/k58+KpbkWpsrqitfOYIXrzvGzXt79FOfXk+inFQvllXwLfVS2u1VGHAAYSDjloOOLgcZ8AxnAYZXLYNMnqXg8gWOcmYfcDhCgMOY9l5wRcbJPriJF+phDJFnmQXllrMvudL+spyCcsxhtQi2LtCH4tyTA0SLYnrvAQI+8UySy2iTFMa7eJICcTiDDi4ngcahUGM+HuIpSDX/iw50kDESCjmltTxOR0hVMQhG0mZGSwSjCSvJJVgjDVxjIci2E/HGc8PyV++3685x3RV2/pAJ5wqiqIoitKv1JfnQ1EURVH6YhBPOFVC6n/wIUa7FJ4eYU4Ir89ilNmtHiG19JR7XMW1lFoA4H0FYV9eqyVaahHnfFgklkZyyrLUIssu0aQFX2a2cONY3sgaLeZgg+eIpB2eI8KRMtz3vPFv77Ipr3j2sjD73xlIP22MCJc4ZZZPjHIgd7Bf3dAd7JEvotRCc0R8w77QJkUgiQiagW/ILhzNwn1PmecDSTHoJcFwkjH6uhSmOpjSCP9HsZmoJGSGsWgTJcFERcAAcp42eU2i0rVbBjSqRUBll/pAZRdFURRFUfqV+vd8KIqiKEqAyi51QX0OPqIiXLhseGGFCBeql6SWVEQobS2lFiCUW4an7NIJSy3HGNKMsJ6LUeawWq9Qx+u2UEQIRy+GRXCkXhzYsR60n+U33ogmYHmFE5GFjy+vP5MV1qJxnZ4Q2JSYQIzud4zEYhLFtV0iLctAkk9s/TAiOOz1shzDN6DweQRZhsNlzeRkljBa9JJa8vlSe98ixcTFtcsuxrL3hfM4DfZ41FhLlPj2J96QQOyKr2DNVC7BGMnBBH0l1tpDggRTvJ38XAnd7u+/4yq71AcquyiKoiiK0q/Up+dDURRFUWyo7FIXDK3BB/vKAinFkFrIPZxQaklHJBGrpdQChHKLEe3ihGVJajnGkWSXUqmlpz5IxAXaH5bZ2ZwiP2tSFxq3kyn4cPk8eZZ6yLHMicpYUslwhItjn98fRL64wn4psVgcouQYIxlTopajXcBO0oygUYnK0FumsbURI/LFkkCsx0aSZnrKPtd5UqyGAMkrxj1hKaeQXMznLK3C53Es8hMQb8V4l44NzpQsCVnv1uNLMGZr3A96Z6lp4yrzKTkfnJG/LcgyRtdksPjRdfBRFwytwYeiKIpyVKNzPuqD+hx8RK1kS2Uzn0e420idTjbGxFKXPR9Udks9IpV4O8zU6KXeDq6XvB1xcniwt6OJ12sx8ngU/gWoLrywKeOXU7IplZ6Zf5nO39NOhvanaT97QaR1Ywxvh8/1fX+LSF4SRkqjLtUPJGKXYnhHIj0o0v4Y9UbKdM7jQR6JosdD8p7E8YJwDhFeu8VY7bbQJqdX51wy7PkQTmP4I2LYhxNBS+tKD6zAC1L4nJ6QZd4Rvi+dON4OS+eFV1o9B0ok9Tn4UBRFURQbKrvUBTr4UBRFUYYMPQvZVTZ6qPR4JZr6H3xIsktBShFTp8eQWjiNeiPn8SCJJZho2ujaZZdYUkvKnjLdlEx6bKoltTRbpBYglEDSgrySEhzLruAq9oxlQWmH8W77hbZDsohGkmBcMaG0UjWSSjASnkVWEaQWX5j8acCzKEmCMVZnDdoxJpbyevDJJpbGygtS7B5NDuWcOfYVBBBLgrEtKyH1lr4DPZaiJDnGs9cHapljv2yKEkn9Dz4URVEUJUBll7pABx+KoijKkEGjXeqDuhp8+CjMvhb9nxzZEsgutFuKcBGkFs7tYZRdki8KUkoz5/mg/XGkFjPaxb46bSC31FJqAUK5JU0XjiWVVMIIFzP43y7B5CP8tUbESpW/FPIxkhN4wtKdUr2NSlzSSVYOFW2l3NfGsTVcojRJmvQ4UosvSGvCgri+RR9IIpfEtUkkwVBZjIKRjohYedY3lgQgecUIACIZmuUgtuHAI0u9Tw+cY0Q3UZeq4YlQhhx1NfhQFEVRlD5R2aUu0MGHoiiKMmRQ2aU+qP/BR8SqtnEiXEzZhZJvRUgtANBYqDejWkIZpYmmsZcrtQBhyvRqSS3NLLWQLzaIZpGkFimqRUTw2+Yl/7iFvKAl5OmGs3ziUR/52Hyh3hOkFpZRPMGx7gt94fqgXJPvL8sz7guSilRvrgptP41vsXGkbFViFitCzOUdgSC1+AlXu2XJNdAVfHoGHZZ6PPsKt0xiCaZwXaKSkAGyBGOmN3f7tDESABrfkfYomTxLKvRqOhESjKGI6h9sJQH1P/hQFEVRlACVXeoCHXwoiqIoQwaVXeqDuhx8+JLUYqzjUpjRLsguvC5Hg7BuCyccaxTkmEBuYXlFinBpogRiXGapxUgyZlmdNqnU0kzXJx0htfTU97icWV5xRWexHa+GPxvyMfoSFcFiSDQstfh2uUaKahmU30+xZBTJJc/tRFznOJExtYyeYaTIFzZhla/wmR2WymidF4fXh+HDBKkniQRjtMBRKPZq4z/MtYTCz+yxBBN8JfD9JsnLEeQV/m702J4D1PhRCaJd+NILj1K/vyfq+agLBssiyIqiKIqiHCXUpedDURRFUSRUNhn81P/gI2LSPa/tIiUWk8oc4cJrt3A5WNulyZBL7OXh0rosZBOVRCxpVIsktaRpJRVbNIsktaQcliaSraHiCWEOQW3eUtdT5oiV0ugVAMjS52GbrN9A5ZSlPY6SsUswRhl2OcYW7ZIoO1hMbFEoYsKphGWju5yTK1WIDiF5Ejl25duXrjdkBV5bhW04Cia4zvm8fX/8AKlCe0KUR+G59SmblkONc0IyMWmYa7+3kRIMSxrGceFz6ghylXnZwj76eZaTC3V0Hr49viFN2+UVxCkXXjdDiqldLsBk+H7loTcaulNzVHZRFEVRFKVfqX/Ph6IoiqIU0GiX+qA+Bx/ilPJSH6AR4WLIK2E9R7g0GlKLvdxkkWM4wqWJpBGWWuQEYXHKBXmnggRiUVJLT9kp7I92ilUkwZBbM1jbhVvIC1ILSyCBjNLTBteHj7UhqxRs2DbnRUe4xJFazHqUUokEExVCIUgn8coc8UFlllWc3gUAZIs8t0GShcehEizHuNZyEJnms06QT6i1SO5yaz0nFqN+0PcBR8mIEozUFYuNIT/xfpZx+NWjskv3xBOiYIJ6N0fJ5rgNOo9vvw1wUiTNeH1HvojqxED+8dZol7pAZRdFURRFUfqV+vR8KIqiKIoFx+s1gbbMNpTaUv+DDyFJUujmpAgXLseIdmlMlUa1AGZCsaAsRbjEiYKJs15LuvA2NApSS2OMBGJRUkuPTXxnmCS1cFRLnuUVSOUesnQrs4aMEpYzJB2ZkSwktZCuwNJMULbVAUCO1vTISdKMxzKO3YVerKvAdZtEpTFsLcEjfZZTVGb5xC0tG8mqjNALN7rMSbyMpde5M17QuPU8RhSKII2YUS0RN8DYH0OC4eXjOREZfU7jtlnD74RoF+Ezi9FDUjRLLthP7zddH0NqoXsvyitCkrHiA2/UDWhqMfPUKrsMelR2URRFURSlX6l/z4eiKIqiFNBol/qgvgYfwVMluJlZYgncpRzt4giyS4NUdoQy+THTheiTtMPRMKF00kjSiSSvNDql7fWUeY2W4Dyg40hqYenEslYLUB2pRYLXc2GphSUYjnDJGuXCcUZyMHtUC5czQv0Rv7HPYw2pxdhvl1o8MaoluhxWllbFRpBgbEnGJAlGWuvDWAaH3hWjvhDZwpEsfgN9IIpqcei6GfKFJ2g9efLrF6QZXkPFZ2mEy9yekd2qTLE+jgTDicjirAVD76+PXKFOklHIlr4DDNmJ5Z0c2bMcFiT/YrmEI3YMqSW6DPtlKdpIScZK5JgaJNwT0SRjdYHKLoqiKMqQofgbtcKtVuzfvx9z585FS0sLWlpaMHfuXPz5z3/u8xjf97FixQqMHTsWw4YNw7nnnosXX3zRsOnu7sbVV1+NUaNG4ZhjjsEll1yCN954I/G5r7nmGrS1taGpqQkf+9jHrP357//+b5xzzjkYNmwYPvCBD+Dmm282fyjEQAcfiqIoitJPzJkzB7t27cKTTz6JJ598Ert27cLcuXP7POb222/HqlWrcNddd+HZZ59Fa2srLrjgAhw8eLBos2jRImzcuBGPPfYYtm7dikOHDuHiiy9Gnjxzcc7t+z6++tWvYvbs2da+HDhwABdccAHGjh2LZ599Ft/97nfx7W9/G6tWrUp0HepLdonCOrmcpRYh8sWQV8IbZcgrRpIxklIKEgtLLWkjERjX56xlQ4JBqdQChKPENELShrzCkopUrlxq4QgXlloMeYXKWbLPkn34iUNv7hGWVKjfR/zwU0sRLqZN2E63F9YfKZS7PWrDk6Jd7AnMOJIm70myS1BA9bFIKb4oqfgx6oVjG0h2Ktwgl5JPcWSMEclC64xI0S4gWcVpCOv9rFdqS7+mzFgKkib42rMcwzpBkl9lcSQY+v4w+sLtWNaw8XNkKyVT44RjLGmxDSV542vuFmQsP/x6gWMkfpOiWqTEYlRvPGi9/h1MDOJol5dffhlPPvkkfvnLX2LKlCkAgO9973uYOnUqXnnlFXzoQx8q7YrvY/Xq1Vi2bBk+//nPAwAefPBBjB49Go888giuuuoqdHV1Ye3atfjBD36A888/HwDw8MMPY9y4cXjqqacwa9as2Oe+8847AQB//OMf8Zvf/KakP+vXr8eRI0ewbt06NDU1YdKkSfjtb3+LVatWYcmSJeLaRL1Rz4eiKIoyZKim7HLgwAFj6+7u7vvkEWzfvh0tLS3FP/4A8IlPfAItLS3Ytm2b9Zjdu3ejs7MTM2fOLNY1NTXhnHPOKR6zY8cOZLNZw2bs2LGYNGlS0aacc0uf4ZxzzkFTU1OxbtasWXjrrbfw6quvxm6n/j0fwiAr+GViTEKlsiuUecJpWvCIpC1l09uRsdrKE0vtk0xTNPwOJpqmjbTK9rwdZj6Pyid6xfF2ZCnZQt7wcITljGWSKRB6PKRJphljMmnaWs6QF4S9HbYJp2I+D2HCqZTnw0iv7pX++vZtvxRLyuyu49l7jtXE9HIUCvwTImmeDzHnh8UmZfeMGBNEuZ4aZD1Ymnwa5P8IJmf21NH1Do8yUrcbkz95Umg1vCAMT2blG2FMfuX8I+zlKM3zYXg7pMmnVDY8JTmup/tS8EjxxFKIHg6yEZ5VaRKpdU7EYPSCVMi4ceOM/16+fDlWrFhRdnudnZ044YQTSupPOOEEdHZ2iscAwOjRo4360aNH47XXXivaNDY24rjjjiuxCY4v59xSf0455ZSS8wT7xo8fH6ud+h98KIqiKEpAFaNdXn/9dYwYMaJYzb/2mRUrVuCmm27qs8lnn30WAKyyhO/7kXJF7/1xjultU+654/RFal9CBx+KoijKkKGaeT5GjBhhDD4kvva1r+GLX/xinzannHIKfvOb3+Dtt98u2ffHP/6xxLMR0NraCqDHqzBmzJhi/d69e4vHtLa2IpPJYP/+/Yb3Y+/evZg2bVrRJum5pf709pTs3bsXQKl3pi+G1ODDNujiKteYfBojt4cln4dU5jqezCpNMm2EXWpJC+nTU8V/7fk8pNTpTLmp0yWpRUqdzpNMJaklY0mffkSYQGqULRNIAfvE0lKbBuNfAOjOh+UMlfOC1JIXVsH1banWJXklKYIcY59wGqdMfeV5oOyeZwmmIKXwfpcVCM75IeTlECUYoNReaoNNPSHPhpSC3SbBxPl17Nuvvc/yhZSC3chLUihTXhOf3nWHk3HQBFHjonPZt8teQdlWBwAOTRROmucjyR/0oy1J16hRozBq1KhIu6lTp6Krqwu/+tWv8PGPfxwA8Mwzz6Crq6s4SOjN+PHj0draivb2dpx55pkAgEwmgy1btuC2224DALS1tSGdTqO9vR2XXXYZAKCjowMvvPACbr/99rLPLX2Gb3zjG8hkMmhs7MmptHnzZowdO7ZEjukLnXCqKIqiDB38Km01YOLEifjUpz6FefPm4Ze//CV++ctfYt68ebj44ouNSJfTTz8dGzduBNAjZSxatAi33norNm7ciBdeeAGXX345hg8fjjlz5gAAWlpacMUVV+Daa6/Ff/7nf2Lnzp3427/9W3zkIx8pRr/EPff//M//YNeuXejs7MThw4exa9cu7Nq1C5lMz1zGOXPmoKmpCZdffjleeOEFbNy4EbfeemuiSBdgiHk+FEVRlKObwZ5eff369Vi4cGExMuWSSy7BXXfdZdi88sor6OrqKv739ddfj8OHD2PBggXYv38/pkyZgs2bN+PYY48t2txxxx1oaGjAZZddhsOHD+O8887DunXrkKIJ2XHOfeWVV2LLli3F/w68Lbt378Ypp5yClpYWtLe34+///u8xefJkHHfccViyZAmWLFmS6Do4ftK0ZAPAgQMH0NLSgnG3/xPcYc3wG6nLTSR3NJJk0VhYbbaRUppT+X2NYUTKMQ1h+X3pMJRqRPpIsXxsA5VTYfl9hfKxLu8/HLbthu2xDef2OMbI+UH5RCyp1JuEFWvjpFFnWIKxrU4bR2rJknTEUssRMaolPGe3EcHSMwZ+l9Kis9TyrhdO8nqPylx/KN9MNmE7h3KhzeFC/aFcuP89oXw4R9JNjmSaLMk0VJ/NUl6QXE/Zy1KkBpWRJdc31VOqGLg5ktEyVM/H5izHCeUU3QjThqRIqT7jFc5XWtfzGcKym8lb650s11MSirylnnJh+DmyzdmP84X2fE7d7rME45fUGQhSi4Eh7wg5T/jYdM/zxDk30BA+Pw6VA9seG5KrGsN6Px3a+02hjdfYULCl9ztNOVuaqJ7LjY69TPMrbfW8nxROeE0kyTYC3pEjeO3GZejq6oo1h6Icgr8T02bdjIZ0c/QBfZDLHsG2n/1jTft7tKOeD0VRFGXo4PlGIruy21Bqig4+FEVRlKHDIM5wqoQMrcGHJaGYmFgMXO/Zy2STgr0+jHbJkS3tJ2mC22Aboy/0cThPkG1msBTVEgeb1AKEckusBGKcOr1MqaWnnDb+BUxJhSNWTAmm0WojRcEEkS0ZinbJCKnW40S7mInFStOrG4JmBV9mUpKxYtFYmdZelhOLcaI0ewIqr6HHhrKyw6PjXE4EJnxOjvzw6WvHEDUsvzbNlOpC25zkTDqWo1AK76QRAcPvQ8KcBwb8GegaFSNVuG2+JkZkkBA9xKnWuZ7SqztBenVOBS8muCu/HLQ5GKNaHFRhzkdVeqL0hUa7KIqiKIrSrwwtz4eiKIpydFPFDKdK7TiqBh+SK81IPmbINJ5gI8ySj7kfMGWXWuLFOI8tmsWoE6SWcqNaesokpRSiXI4IMgpHuHQLkgpLMIfzpYnFuD6Tp3VjqMxruGR5bZc8SS3Sei6+RYLx2MUeFp0YCccM+cDIlEc2bnC+0rreZZZJDJmCjzUSjpX23ZBi0uw0ZTnAnvHMCCChI/mkvm0/EUuC4QRqXM9RMIWIFJZkYdy/Af7DkzT5WdR+UTphqUe46nX4N3iwh9oqPZQlu6xZswbjx49Hc3Mz2tra8PTTT4u2P/rRj3DBBRfg/e9/P0aMGIGpU6fiZz/7WdkdVhRFURSlvkk8+NiwYQMWLVqEZcuWYefOnZgxYwYuvPBC7Nmzx2r/X//1X7jggguwadMm7NixA3/1V3+Fz3zmM9i5c2fFnVcURVEUg0Gc4VQJSSy7rFq1CldccQWuvPJKAMDq1avxs5/9DPfccw9WrlxZYr969Wrjv2+99Vb8+Mc/xk9+8pNi5rTedHd3o7s7TM514MCBWH1jx6Et2gWW/b0xo2AEG2PZ+56yIddwVItjj3BhUmX6+DwI/vYYsg/DicOCNo1IFuq3tFZLUqnFXK+lRzJ5l6ST94TEYiyvvJcPy6a8YpdgMoVoliNUxxEuhgSTFyJf8kK0i21JcsPdXcnaLvZmAjlGjmQhW5ZMOIxK6qPFJc9JpPhXi++zLGW34SdS+sUTnF1c+8ViC/QhwbANLUcfSDA+RaKBo0PcZO8PJxyDW+Z9lqQeKfJFXEPH/Lc3YuSLYWOPGTKjcIKHz97GQOL4fq/PUF4bSm1J5PnIZDLYsWNHMTVrwMyZM7Ft27ZYbXieh4MHD+L4448XbVauXImWlpbiNm7cuCTdVBRFURRlEJNo8LFv3z7k8/mSZXNHjx5dssSuxL/+67/i3XffLa68Z2Pp0qXo6uoqbq+//nqSbiqKoihHK16VNqWmlBXt0nvlOt/3Y61m9+ijj2LFihX48Y9/jBNOOEG0a2pqQlNTk7jf7EssMwB9ySjRLrY4ESzlwkuzpyP64olhA4IEY1jYPwNHswTt11pqeddYo6VHJkkqtcgRLtSXPPWlYMNSSzZvL3OEi5xYjKI5jGVELGvdM5ZkTT3/YS+bUgtFqhRc/I5w681ol7DscmQHJ7ridjiqp2hjj94xXe/2Zy+JBCNJKlWTYGzt0UX0pfsmSSMktRjfg26C33blyjUxiCUhiFJLkvPEt60lKrvUB4kGH6NGjUIqlSrxcuzdu7fEG9KbDRs24IorrsC///u/F5f4VRRFURTl6COR7NLY2Ii2tja0t7cb9e3t7Zg2bZp43KOPPorLL78cjzzyCD796U+X11NFURRFiUKjXeqCxLLLkiVLMHfuXEyePBlTp07F/fffjz179mD+/PkAeuZrvPnmm3jooYcA9Aw8vvSlL+E73/kOPvGJTxS9JsOGDUNLS0sVP8rgwRPGdPmEKwbQXHyUpmKC6JPOG0faYfnGTCIWyC6gurCcqbLU0mPTWHocJxCLIbVw+YgR4RL2K1jbpTtHa7sIUksuYYQLr+1STC4WI7IgFpIEEyQZM0JPqCgsFGQsFyJEuFg7ICWuktP32Ws58oajYyznr5oE47A0UrrOik9JyAyXO62nYqzVIsFSi02CSSqvSLpyJevPDHU0w2ldkHjwMXv2bLzzzju4+eab0dHRgUmTJmHTpk04+eSTAQAdHR1Gzo/77rsPuVwOf//3f4+///u/L9Z/+ctfxrp16yr/BIqiKIpSQDOc1gdlTThdsGABFixYYN3Xe0Dx85//vJxTKIqiKIoyRDmq1nbxhFnsUj2TFyJIwv32NozjDGmEogmonKUyJyULIk4aHbtjOR/DTchRBnlDdgkJJBa2PULySpY+T9RaLYC5XotNagFCKeW9vBDhklBq4TJLLEE0C6/bkk0otXhRES6AkGRMKEsIkS9WCUbIL+cLycTEeqFclEbE9ySOIMJhONba4jMXlYSspyuOtV6yNyQW5AonojZYDuG1byQJJg6c2CzlBiei/eGNMKJkkkotERKMH0eiEWwqyY03YKjsUhccVYMPRVEUZWjjeImTPFvbUGpL3Q8+ohd3TDZ0Zy+E7CmhX8iW32ns7fDol06WPAVpJ0/1oQ2nd+dcG+lCVzjnBv+AjQN7OPjd4gmlgUcmqbfDli4dkFOms33g8WAPxyHyglTi7ejOc7nnM2Vy9twe7O0wPB959nbwJEVOpGGbcFr9n43m5NKeG8crkkpeDX5YjDwb/MNemK7pFvczVfKCWGqN/gm/yMWJpVLODae03s/lrfuR51V6PbtNHNibEnwPpFz7/hTnwhdyhcTwjvgRj55Rn/DjRHlQ6tJLogwYdT/4UBRFUZQiKrvUBTr4UBRFUYYO1cjToWOPmjOkBh/GnL4IH6CZPjos5416V7BxS2xYUsn6oTs3Q37wtJMjm7CeV7XNChPzglwcaUOWSYYx4ZQnufpuSTkDruN8HtE5PIwcHULKdNOmR26ppdQChHJLLqHUYuT2yAu5PWwTSo006oKtgPT4Oha3uWEr5PYQzyP0xXg/Ch2WUqRXW4IxzuP0ZRk0bZ98Kp69YO+wpEF5Pnxj1q6Rc17qQTTB5FJjkqsgqZA047tCPeUc4XKxfUOigbXMMoov1Is4vf5VlIQMqcGHoiiKcnSja7vUBzr4UBRFUYYOOuejLhjygw8jBTbVm1EtpTJKb5usEf1RWuaol7whWdjlmJToew+lGT5/uuDozgvHsXQj5Rzh9mTZpefzsKSSiSG7xJFaolanraXUAoRySy5OVIuQOt1Mow4qWxJw1CC9OlO8bcKKpHHO7wvSjLnyrFOoK42A6W1bFQmGVYTwdRBXwzWllvA54NwdRl6QIL973h7t4lAkGocX+dIfJCNHvBCdE7QvRbU08LLDbnRZiHwJZBpRlomRPyZWWVEqZMgPPhRFUZSjCB+9R8TltaHUFB18KIqiKEMGnfNRHwytwQe7yi3SgydEuEiJxbJe6Ar1XLt9IFNIskyKhuAsU8SBE5EFUooo1wjVeUFqYXkpg9K+Z6VkYl50tEu3JZIFiF6dtpZSS0+55zPLUS1cb49wsSYT610u3AsxwkWoN1as5ZWt2K1uqzYCpEIDR5AcGUnKsNm4/PzUVILpOwlZ7/OYCoQ9PMYxIjsK75URecKr13LyOL6eQv75OATnF6JaWFLxWYKhsiGluPaImGLZogKW2Eo2MaSWqGRmA4qPKsz5qEpPlD7oe8ESRVEURVGUKjO0PB+KoijK0Y1Gu9QF9T/4EPx+wbPjCVKMJMHkhMRihqxCckzW6Smz1MBSi0vuXNePngVlRqHwWjC5EltJgskLvtK8kDiMz5MpU3ZhGcWQXSKkFi7XUmoBQrmlIqlFqo9IMlaJG9eUYyz77Qsn95JGKpdgbBEwpeex29dSghHXXHHtNkUJhlap9Vl2oXonHyPJmFRv65cglxiRLCy1SNEuDSzTUH1BmrFKMTCfFWOdIDERWbRNn3UDgYfK+6ILy9UclV0URVEURelX6t/zoSiKoigFNNqlPhhSgw9zMnohMRKvyUIudi7nhDIn38qRTNHtkRzi5gu29kRcLneKPbsO9YXKaSdsp9EpXSMm5cSQbny7QytLUS15I7EYJ0KzJRkjqUOIZLGt1dJjHx4blUSsllILAOQL9p6QNKwqUgsAJ7CRvPFxvtfMzFnWet9mWgMJJlAPJEml+hJMRBIymHKASxfUCHbJSYm4etpxKMmYsc4KyzEp31pvdjvGDU0S7WKs1WKXV7wGu3wTlL0Ge9RPUjlGLEfIGgMaBaNzPuoClV0URVEURelXhpTnQ1EURTnKUc9HXVD/gw/B9R08O2aEC6hMcowgr2Q8dvdT1AqdqEFyxVrwyOXqudnw/FIUCskuRbmFl5EQpmR7vM4MfTYpkobXgjniN5b0g6UWlle6OSLGs0sqh/MswZA0ZZFdaim1AKHcwvKKH0deSSK1cL2UWCwpJCv4dA+dYh1K6oB4EowjRYpYqul1qLEEI3wiQ0Wxt+hSvWes50LlXMGe5Y0cR7VQ2/kYScaSRLuw/GOsuWJPIAZJauGEY2myCaJdaD9LMLyWjy/JMUkiX2IkJOt3BvngY//+/Vi4cCGeeOIJAMAll1yC7373u/iLv/iLPrrj46abbsL999+P/fv3Y8qUKbj77rvx4Q9/uGjT3d2N6667Do8++igOHz6M8847D2vWrMGJJ56Y6NzXXHMNtm7dihdeeAETJ07Erl27jL68+uqrGD9+fEkff/rTn+JTn/pU7OugsouiKIqi9BNz5szBrl278OSTT+LJJ5/Erl27MHfu3D6Puf3227Fq1SrcddddePbZZ9Ha2ooLLrgABw8eLNosWrQIGzduxGOPPYatW7fi0KFDuPjii5GnuU1xzu37Pr761a9i9uzZffbpqaeeQkdHR3H75Cc/meg61L/nQ1EURVECqpjn48CBA0Z1U1MTmpqaLAfE4+WXX8aTTz6JX/7yl5gyZQoA4Hvf+x6mTp2KV155BR/60IdKjvF9H6tXr8ayZcvw+c9/HgDw4IMPYvTo0XjkkUdw1VVXoaurC2vXrsUPfvADnH/++QCAhx9+GOPGjcNTTz2FWbNmxT73nXfeCQD44x//iN/85jfiZxk5ciRaW1vLvhb1NfjwncLGLk8qWhYoMIJNONqFE4vFiHzpdsJLxbPr3XwS2YUlEI5wIXmFXMjWxGIxQiXygovfkFrEhGOlSdOMSB5BXmF7KcLlCEkphqRVkFtqKbUAodxSU6mFy3QfHKENCSOxmLmHSk7J/qQSjC8e0Tf9JsEIbn2jDXqvfJJJjPfUIjE4LKmwLMP6bAOV88KNSxLtwoeJsguVpagWi9QCAF6h3pBlUkLbUsIxQY6x3YtYa8I4fswQr+pQzVDbcePGGfXLly/HihUrym53+/btaGlpKf7xB4BPfOITaGlpwbZt26yDj927d6OzsxMzZ84s1jU1NeGcc87Btm3bcNVVV2HHjh3IZrOGzdixYzFp0iRs27YNs2bNKuvcfXHJJZfgyJEjmDBhAhYvXoy//uu/TnR8fQ0+FEVRFKUvqjjn4/XXX8eIESOK1ZV4PQCgs7MTJ5xwQkn9CSecgM7OTvEYABg9erRRP3r0aLz22mtFm8bGRhx33HElNsHx5Zzbxvve9z6sWrUK06dPh+u6eOKJJzB79mw8+OCD+Nu//dvY7ejgQ1EURVEsjBgxwhh8SKxYsQI33XRTnzbPPvssAPskb9/35cnfBXrvj3NMb5tyz82MGjUKixcvLv735MmTsX//ftx+++1H8eDDlmTMkFpCJ6+UcIzXbWFpwHDhOjy7Pv4lZHmFpQyu57Lr2NeLiSIvzCPmc3p+37LLEUF2YUnFjAaySy0ZIcKF7bMFKaWWUgtAcktCqSVSXikpF9z6tfA0G98RfuH/SyNgenfJYLBLMPxF6Nj7ZCQZI3XSNW4n2ZNNkMePJVN+TgzZhZ8DQWJNcp8NaUJM+EVGgmRik1qAUGLhCBd6BY2yKccgWdkW7TJY8Kog83jJjv/a176GL37xi33anHLKKfjNb36Dt99+u2TfH//4xxLPRkAwr6KzsxNjxowp1u/du7d4TGtrKzKZDPbv3294P/bu3Ytp06YVbZKeOy6f+MQn8G//9m+JjtFoF0VRFGXoEMgulW4JGDVqFE4//fQ+t+bmZkydOhVdXV341a9+VTz2mWeeQVdXV3GQ0Jvx48ejtbUV7e3txbpMJoMtW7YUj2lra0M6nTZsOjo68MILLxRtyjl3XHbu3GkMjOJQ/54PYcaTb5lwyhNS+Rd0lmajufQrm1OZc24P1zKq9ozVcGnSKrXdRD/RgrTsvdtO23J7wPS2RMHn5zwfHv1MyUd4PtjDkfPsno9uwUuUydu9HexV4pwe2YLnKam3w1ydNsLbAYQej0q8HXwbhIm91h9d1VrV1uoI8K22xoRXTmPBjVfZC8J5JJgkXhDPOB17O/gdtE+EZO+Am+MJpeEHLabMsXhDevbbc3v4Hv/0ZxvExpyUaZ8UKn0eM1+HvRzY+JKHQ/SCcH3fk0x7+mXbr4m5opg4cSI+9alPYd68ebjvvvsAAH/3d3+Hiy++2Jjwefrpp2PlypX43Oc+B8dxsGjRItx6662YMGECJkyYgFtvvRXDhw/HnDlzAAAtLS244oorcO2112LkyJE4/vjjcd111+EjH/lIMfol7rn/53/+B4cOHUJnZycOHz5czPNxxhlnoLGxEQ8++CDS6TTOPPNMuK6Ln/zkJ7jzzjtx2223JboW9T/4UBRFUZQiVZhwWsmvhQjWr1+PhQsXFiNTLrnkEtx1112GzSuvvIKurq7if19//fU4fPgwFixYUEwytnnzZhx77LFFmzvuuAMNDQ247LLLiknG1q1bh1QqHFnGOfeVV16JLVu2FP/7zDPPBNATdXPKKacAAG655Ra89tprSKVSOO200/DAAw8kmu8BAI7vD/48sgcOHEBLSwvG3XYL3GHN8BvpdxSFwjlp8hQUyg1p8jakaW5FKiw3kU0T1Tc3hFlIm1M5a31jwZvRRJ4M9nA0USZT9XzUzvPB9XXr+YgTsZmkbcHzYXxO4zOTjcf1QWiq3Zbr3by93rBhjwOHxuaS2/ack74DjGO5nsqFvjuU1dT4PILnw7z39eP5MOZ/0M/NfJpDeqPraRoYgsTF5n7fauulfXhHjuC1ZTeiq6sr1gTOcgj+Tpw//mo0uJVFpeS8bjy1+7s17e/RTn16PoSJfraMx+yO5z9QLJ3wHzpjhcx89OUJXJSm7MJu5bBsTDL1eGKpffBhk3ekVW2llWyNNPLkYzdX76U/9IUyDxRy1Ha3kasjevDB7WR5QinbFO4L78/HklpowCGsVGudXFqlAYf4R6eGw3nbeEcWRUgycO0DkWQSjHAm4e+zNKHMXDG3dIqsmPFdmIhqSC3CBEnboMjhVWWlwUdeGGQI91jKL+FbU63Tfh5MCIMPvqDSQCR43cwBiXQc1bNcxtdQSMdu7asg0SiKjfocfCiKoiiKDc9Hxb8AEka7KMnRwYeiKIoydPC9nq3SNpSaUv+DD0t+BSB0vfMMdY9Gs5zbI8euVZr2HsdzGMgqLG/kjPkUnEad53CQ1CLM/3AF+8g+GSv5SnM+qL8sjRTsc4LskhFkF0Om8djeLrvwPI7APqnU4vHcDpZXovJ4SFKLMC+i7DTpSX88JXVVF9qXFsw1qwUJhqVKQTMJrJNGwMSRYEz7wtwSIYLCWMjWWLE2rPZden+EVOLBK8733pRaqJ7Tqwv3vvw5H/aynALdLpmY+T9K90s5P6SIGE+QWqxllVqUMqn/wYeiKIqiBFQxvbpSO3TwoSiKogwddM5HXVCXgw/HSCZGOyxpkTmxGEe+cDlHfltHWLnS8LhyBEkgU9A0+wY3dDhnSMbh+gZBUmGphVewdRO8TJK8IkXksKwSSCA53y6jsFzFsouZop4lrb6lFiCUWCSpxbiHktQihc9aIljipUtPKLUMxHeVzc3Nr4CQK8qQYDg1u3BZAs3EHgHT630UumpGuNjLRfXAuH3cP9Y6qMiRLBFSC0Chtl5pXUnZt18UI6olyb03okOE8FpDguGyJMdQOYi+E+WVsOwJyccQJbVQf211wMC8DuHJ1fNRD2h6dUVRFEVR+pW69HwoiqIoihUfVfB8VKUnSh/U1eDDQY8LXHyurNEulHXPYVd/9Pl8Q97h2eDhwYF8IckuKUF2cQVJJSrCRYp68YSQhziRLyyBBPZ5ixQDyPKKGT1kr89bpBYglMDEtVqEdVsQlUwMsEe2JE0gNpiklihKc3b1FONIMFGJyIQkZOb6LNESjERwSlmWITmRE4hRhIshmXB/vdLPZq5ea//s0nPgSCFG0peTbSlzKdpFjILhct8SjLyGixAxE2PNF/vaLvT89BXJ05+RMCq71AUquyiKoiiK0q/UledDURRFUfrE82C6NMttQ6kl9T/4kNzjgeySN3yv4e4YTh9DauEIF14vJtXTZsoVolSonJJkF0GCcZJkLxL6zbDUYn4em+xi32+W+24DsMsrpeUeG0le8aXolRhSi7U81KQWiRpIMOF+exuxJBjJBR/U5y11AHhFe45w8Xj9E5ZMUoIcU5RdJHmOyonXc4mvL8SJdhFlF7axyC4cyWKsCSNEtYjRMY69Pjh2UK7norJLXaCyi6IoiqIo/Ur9ez4URVEUJUA9H3VBfQ4+IqQWAGESHHYFkpveXMaDE1qF9TyL3jPWTyiNFDGkFkNeCdtj2cURZJcoz2XSaBdfsJEkpaCeJRUzqVqpXNJTliQVlk/sicOCNUVYIvOlNVciEoj1Lltd61WSWspe06MSkrQjqQEJJZggUkNsLoYE4wgSg1HOW/okJf+yq6mm7GKJcOk5QVDHidJK9/dVL+lI0jMRef9jRL6ICb1skS+CXCInJxPOKR5b+KBCnwZUgtEMp3WByi6KoiiKovQr9en5UBRFURQLvu/B9yuLVqn0eCWa+hx8CMl+jCW/ixmL+DD2G5JtiqJgyM9oJCgjyYSlhMCbnDcklRSVhagWwT8rRbgkiXyRol1EqcViI62JI9YLNjZ5pcQ+KEsySjWklt7lKIay11WQYCSCWygdFkeC8aUjLBKDIZEIsoMYkSKu0VJqI6/bAmt9vGiXvvEtycZ6GhEOiBPtYrORjotYqwWIKdNYol2MSzWQsovvVy6b6JyPmlOfgw9FURRFseFXYc6HDj5qjs75UBRFURSlX6l/z4e0TkcwAzxvtzVWxCZ/ouOy/5XkFY52Ya9xwUfrGNPveZZ/qW1f9Uy1k4yZckjfxxrXh6NUDGN7e74ge4jJwvw+6nqdpyKppdBOrAiXKlAT13N5+aziNR0R+eJzVEeMLklufWM9EON1Kxixe1+QUSQJBiQZiBEsnqUuRqSTbFPmjRAOS7rmi83GF35WyhEzVG+LaullUyxL93gg8Tx7Zrwk6JyPmlP/gw9FURRFCVDZpS6or8FH8EwZv4Q5hbJl8ikn2hC8IPyL2+efN8YEL7vXovhzSPJkGPWwI3k+BPMoxNdG9Ij03YgfY1aZL1zPRB6MOB6OGDbSOSMdSVX6vum3yXZBfys5X5zJp06pQRwviNGE9AvecJsUPCx8v4054vZ76QsTS33xGYrorEAc70iyBoXm4tzPKI9IDI+J6QWhe8v2ESndZW+M2d6ATkBVBiX1NfhQFEVRlD7wPQ9+hbKLhtrWHh18KIqiKEMHlV3qgrocfIiuVbYJXP9saxwoza7jevvkU19wL9r32/uXKDe3eXorid8VyQ8aIcHEOq4S+wKJ5RWj7Qippd7klShq8D3p214D4/PaJRjp9RHfU5sMIBzH+TQiJZXelPtMVEuaKe1GzAaj661tSjJKJZNZbeeP04aiWKjLwYeiKIqiWPH8xD/uSlDPR83RwYeiKIoydPB9mHkXym1DqSV1P/gwIlwi5BBxZU3DKCyKLlJJpolquwIG9FVI6JJ2kkg6Un3CnApVj0RQitguufkKsE4iHFeuzMfV0vMmHSr9R8QzUemP5loR7/so/n6rjNLHsVHyjkotShLqfvChKIqiKAG+55vz+8ppQz0fNUcHH4qiKMrQwfdQueyioba1pqzBx5o1a/Ctb30LHR0d+PCHP4zVq1djxowZov2WLVuwZMkSvPjiixg7diyuv/56zJ8/v+xOy+5+/o/yfIBmEEy5aZOH2Ki5SmEdVbksQ+zS1itJJQBp5ddEt1OSY5K0AQyiMKUqUY0Xq8xLMhgvpXo+6oPEC8tt2LABixYtwrJly7Bz507MmDEDF154Ifbs2WO13717Ny666CLMmDEDO3fuxDe+8Q0sXLgQjz/+eMWdVxRFURSl/kjs+Vi1ahWuuOIKXHnllQCA1atX42c/+xnuuecerFy5ssT+3nvvxUknnYTVq1cDACZOnIjnnnsO3/72t/GFL3zBeo7u7m50d3cX/7urqwsA4B05Er+j5Y7kyzus17mH2Ki5Wp6PajRSy0s70LdtoH9F1vDzDyqv12D8uV4JdeD5CL67+8OjkPO7K5ZNcshWqTeKiJ+A7u5uP5VK+T/60Y+M+oULF/pnn3229ZgZM2b4CxcuNOp+9KMf+Q0NDX4mk7Ees3z58iBFnW666aabbkNk+/3vf5/kT04iDh8+7Le2tlatr62trf7hw4dr1t+jnUSej3379iGfz2P06NFG/ejRo9HZ2Wk9prOz02qfy+Wwb98+jBkzpuSYpUuXYsmSJcX//vOf/4yTTz4Ze/bsQUtLS5IuH1UcOHAA48aNw+uvv44RI0YMdHcGLXqd4qHXKR56naLp6urCSSedhOOPP75m52hubsbu3buRyWSq0l5jYyOam5ur0pZSSlkTTp1ekzl93y+pi7K31Qc0NTWhqamppL6lpUVf7hiMGDFCr1MM9DrFQ69TPPQ6ReO6iacZJqK5uVkHDHVCoidh1KhRSKVSJV6OvXv3lng3AlpbW632DQ0NGDlyZMLuKoqiKIpS7yQafDQ2NqKtrQ3t7e1GfXt7O6ZNm2Y9ZurUqSX2mzdvxuTJk5FOpxN2V1EURVGUeiexD2zJkiX4t3/7NzzwwAN4+eWXsXjxYuzZs6eYt2Pp0qX40pe+VLSfP38+XnvtNSxZsgQvv/wyHnjgAaxduxbXXXdd7HM2NTVh+fLlVilGCdHrFA+9TvHQ6xQPvU7R6DVSeuP4fvLYpzVr1uD2229HR0cHJk2ahDvuuANnn302AODyyy/Hq6++ip///OdF+y1btmDx4sXFJGNf//rXK0sypiiKoihK3VLW4ENRFEVRFKVcajv1WFEURVEUpRc6+FAURVEUpV/RwYeiKIqiKP2KDj4URVEURelXBs3gY82aNRg/fjyam5vR1taGp59+uk/7LVu2oK2tDc3NzTj11FNx77339lNPB5Yk1+lHP/oRLrjgArz//e/HiBEjMHXqVPzsZz/rx94ODEmfpYBf/OIXaGhowMc+9rHadnCQkPQ6dXd3Y9myZTj55JPR1NSED37wg3jggQf6qbcDR9LrtH79enz0ox/F8OHDMWbMGHzlK1/BO++800+9HRj+67/+C5/5zGcwduxYOI6D//iP/4g85mj9DlcKDOTCMgGPPfaYn06n/e9973v+Sy+95F9zzTX+Mccc47/22mtW+z/84Q/+8OHD/WuuucZ/6aWX/O9973t+Op32f/jDH/Zzz/uXpNfpmmuu8W+77Tb/V7/6lf/b3/7WX7p0qZ9Op/1f//rX/dzz/iPpNQr485//7J966qn+zJkz/Y9+9KP909kBpJzrdMkll/hTpkzx29vb/d27d/vPPPOM/4tf/KIfe93/JL1OTz/9tO+6rv+d73zH/8Mf/uA//fTT/oc//GH/s5/9bD/3vH/ZtGmTv2zZMv/xxx/3AfgbN27s0/5o/Q5XQgbF4OPjH/+4P3/+fKPu9NNP92+44Qar/fXXX++ffvrpRt1VV13lf+ITn6hZHwcDSa+TjTPOOMO/6aabqt21QUO512j27Nn+jTfe6C9fvvyoGHwkvU4//elP/ZaWFv+dd97pj+4NGpJep29961v+qaeeatTdeeed/oknnlizPg424gw+jtbvcCVkwGWXTCaDHTt2YObMmUb9zJkzsW3bNusx27dvL7GfNWsWnnvuOWSz2Zr1dSAp5zr1xvM8HDx4sKYrSw4k5V6j73//+/j973+P5cuX17qLg4JyrtMTTzyByZMn4/bbb8cHPvABnHbaabjuuutw+PDh/ujygFDOdZo2bRreeOMNbNq0Cb7v4+2338YPf/hDfPrTn+6PLtcNR+N3uGJS1qq21WTfvn3I5/MlC9ONHj26ZEG6gM7OTqt9LpfDvn37MGbMmJr1d6Ao5zr15l//9V/x7rvv4rLLLqtFFweccq7R7373O9xwww14+umn0dAw4K9Dv1DOdfrDH/6ArVu3orm5GRs3bsS+ffuwYMEC/OlPfxqy8z7KuU7Tpk3D+vXrMXv2bBw5cgS5XA6XXHIJvvvd7/ZHl+uGo/E7XDEZcM9HgOM4xn/7vl9SF2Vvqx9qJL1OAY8++ihWrFiBDRs24IQTTqhV9wYFca9RPp/HnDlzcNNNN+G0007rr+4NGpI8S57nwXEcrF+/Hh//+Mdx0UUXYdWqVVi3bt2Q9n4Aya7TSy+9hIULF+If//EfsWPHDjz55JPYvXu3Lidh4Wj9Dld6GPCfeqNGjUIqlSr5JbF3796SkXFAa2ur1b6hoQEjR46sWV8HknKuU8CGDRtwxRVX4N///d9x/vnn17KbA0rSa3Tw4EE899xz2LlzJ772ta8B6Pkj6/s+GhoasHnzZnzyk5/sl773J+U8S2PGjMEHPvABtLS0FOsmTpwI3/fxxhtvYMKECTXt80BQznVauXIlpk+fjn/4h38AAPzlX/4ljjnmGMyYMQO33HKL/qIvcDR+hysmA+75aGxsRFtbG9rb24369vZ2TJs2zXrM1KlTS+w3b96MyZMnI51O16yvA0k51wno8XhcfvnleOSRR4a87pz0Go0YMQL//d//jV27dhW3+fPn40Mf+hB27dqFKVOm9FfX+5VynqXp06fjrbfewqFDh4p1v/3tb+G6Lk488cSa9negKOc6vffee3Bd82s1lUoBCH/ZK0fnd7jSiwGa6GoQhLOtXbvWf+mll/xFixb5xxxzjP/qq6/6vu/7N9xwgz937tyifRCmtXjxYv+ll17y165de1SEaSW9To888ojf0NDg33333X5HR0dx+/Of/zxQH6HmJL1GvTlaol2SXqeDBw/6J554ov/Xf/3X/osvvuhv2bLFnzBhgn/llVcO1EfoF5Jep+9///t+Q0ODv2bNGv/3v/+9v3XrVn/y5Mn+xz/+8YH6CP3CwYMH/Z07d/o7d+70AfirVq3yd+7cWQxJ1u9wpTeDYvDh+75/9913+yeffLLf2Njon3XWWf6WLVuK+7785S/755xzjmH/85//3D/zzDP9xsZG/5RTTvHvueeefu7xwJDkOp1zzjk+gJLty1/+cv93vB9J+iwxR8vgw/eTX6eXX37ZP//88/1hw4b5J554or9kyRL/vffe6+de9z9Jr9Odd97pn3HGGf6wYcP8MWPG+H/zN3/jv/HGG/3c6/7l//2//9fnd41+hyu9cXxffYGKoiiKovQfAz7nQ1EURVGUowsdfCiKoiiK0q/o4ENRFEVRlH5FBx+KoiiKovQrOvhQFEVRFKVf0cGHoiiKoij9ig4+FEVRFEXpV3TwoSiKoihKv6KDD0VRFEVR+hUdfCiKoiiK0q/o4ENRFEVRlH7l/wfKKjDVtQhRiQAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots()\n", "\n", "im = ax.imshow(np.transpose(e.v()),\n", " interpolation=\"nearest\", origin=\"lower\",\n", " extent=[mg.xmin, mg.xmax, mg.ymin, mg.ymax])\n", "fig.colorbar(im, ax=ax)" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.11.2" } }, "nbformat": 4, "nbformat_minor": 4 }