{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Simple FEM-BEM coupling for the Helmholtz equation with FEniCS" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Background" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "For this problem, you will need FEniCS installed alongside Bempp.\n", "\n", "In this tutorial, we will solve the problem of a wave travelling through a unit cube, $\\Omega = [0,1]^3$ with different material parameters inside and outside the domain. The incident wave is given by\n", "\n", "$$\n", "u^\\text{inc}(\\mathbf{x})=\\mathrm{e}^{\\mathrm{i} k \\mathbf{x}\\cdot\\mathbf{d}},\n", "$$\n", "\n", "where $\\mathbf{x}=(x,y,z)$ and $\\mathbf{d}$ is the direction of the incident wave. In the implementation we use, $\\mathbf{d} = \\frac{1}{\\sqrt{3}}(1,1,1)$.\n", "\n", "The PDE is\n", "\n", "$$\n", "\\Delta u + n(\\mathbf{x})^2 k^2 u = 0, \\quad \\text{ in } \\Omega\\\\\n", "\\Delta u + k^2 u = 0, \\quad \\text{ in } \\mathbb{R}^3 \\backslash \\Omega\n", "$$\n", "\n", "In this example, we use \n", "\n", "$$\n", "n(\\mathbf{x}) = 0.5\n", "$$\n", "Since the interior wavenumber is constant one could have also used a BEM/BEM coupling approach. However, here we demonstrate the use of FEM for the interior problem using the FEniCS finite element package." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### FEM Part" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In $\\Omega$, the FEM part is formulated as\n", "\n", "$$\n", "\\int_\\Omega \\nabla u\\cdot\\nabla v -k^2\\int_\\Omega n^2uv - \\int_{d\\Omega} v\\frac{\\partial u}{\\partial \\nu} = 0,\n", "$$\n", "\n", "or\n", "\n", "$$\n", "\\langle\\nabla u,\\nabla v\\rangle_\\Omega - k^2\\langle n^2u,v\\rangle_\\Omega - \\langle \\lambda,v\\rangle_\\Gamma=0,\n", "$$\n", "\n", "where $\\lambda=\\frac{\\partial u}{\\partial \\nu}$.\n", "\n", "Later, we will write this as the following operator equation\n", "\n", "$$\n", "\\mathsf{A}u-k^2 \\mathsf{M}u-\\mathsf{M}_\\Gamma \\lambda = 0\n", "$$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### BEM Part" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In $\\mathbb{R}^3 \\backslash \\Omega$, we let $u = u^\\text{inc}+u^\\text{s}$, where $u^\\text{inc}$ is the incident wave and $u^\\text{s}$ is the scattered wave. As given in Integral equation methods in scattering theory by Colton & Kress,\n", "\n", "$$\n", "0 = \\mathcal{K}u^\\text{inc}-\\mathcal{V}\\frac{\\partial u^{inc}}{\\partial \\nu},\\\\[2mm]\n", "u^\\text{s} = \\mathcal{K}u^\\text{s}-\\mathcal{V}\\frac{\\partial u^{s}}{\\partial \\nu},\n", "$$\n", "where $\\mathcal{K}$ and $\\mathcal{V}$ are the double single layer potential operators. Adding these, we get\n", "\n", "$$\n", "u^\\text{s} = \\mathcal{K}u-\\mathcal{V}\\lambda.\n", "$$\n", "\n", "This representation formula will be used to find $u^\\text{s}$ for plotting later.\n", "\n", "Taking the trace on the boundary gives\n", "\n", "$$\n", "u-u^\\text{inc} = \\left(\\tfrac{1}{2}\\mathsf{Id}+\\mathsf{K}\\right)u -\\mathsf{V}\\lambda.\n", "$$\n", "\n", "This rearranges to\n", "\n", "$$\n", "u^\\text{inc} = \\left(\\tfrac{1}{2}\\mathsf{Id}-\\mathsf{K}\\right)u+\\mathsf{V}\\lambda.\n", "$$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Full Formulation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The full blocked formulation is\n", "\n", "$$\n", "\\begin{bmatrix}\n", " \\mathsf{A}-k^2 \\mathsf{M} & -\\mathsf{M}_\\Gamma\\\\\n", " \\tfrac{1}{2}\\mathsf{Id}-\\mathsf{K} & \\mathsf{V}\n", "\\end{bmatrix}\n", "\\begin{bmatrix}\n", " u\\\\\n", " \\lambda\n", "\\end{bmatrix}=\\begin{bmatrix}\n", " 0\\\\\n", " u^\\text{inc}\n", "\\end{bmatrix}.\n", "$$\n", "\n", "This formulation is not stable for all frequencies due to the possibility of interior resonances. But it is sufficient for this example and serves as a blueprint for more complex formulations." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Implementation" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We begin by importing DOLFIN, the FEniCS python library, Bempp and NumPy." ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [], "source": [ "import dolfin\n", "import bempp.api\n", "import numpy as np" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we set the wavenumber ``k`` and the direction ``d`` of the incoming wave." ] }, { "cell_type": "code", "execution_count": 21, "metadata": {}, "outputs": [], "source": [ "k = 6.\n", "d = np.array([1., 1., 1])\n", "d /= np.linalg.norm(d)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We create a DOLFIN mesh. Later, the boundary mesh will be extracted from this.\n", "\n", "A mesh could be created from a file changing this line to ``mesh = dolfin.Mesh('/path/to/file.xml')``." ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [], "source": [ "mesh = dolfin.UnitCubeMesh(10, 10, 10)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we make the DOLFIN and Bempp function spaces.\n", "\n", "The function ``fenics_to_bempp_trace_data`` will extract the trace space from the DOLFIN space and create the matrix ``trace_matrix``, which maps between the dofs (degrees of freedom) in DOLFIN and Bempp." ] }, { "cell_type": "code", "execution_count": 23, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "FEM dofs: 1331\n", "BEM dofs: 1200\n" ] } ], "source": [ "from bempp.api.external import fenics\n", "\n", "fenics_space = dolfin.FunctionSpace(mesh, \"CG\", 1)\n", "trace_space, trace_matrix = \\\n", " fenics.fenics_to_bempp_trace_data(fenics_space)\n", "bempp_space = bempp.api.function_space(trace_space.grid, \"DP\", 0)\n", "\n", "print(\"FEM dofs: {0}\".format(mesh.num_vertices()))\n", "print(\"BEM dofs: {0}\".format(bempp_space.global_dof_count))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We create the boundary operators that we need." ] }, { "cell_type": "code", "execution_count": 24, "metadata": {}, "outputs": [], "source": [ "id_op = bempp.api.operators.boundary.sparse.identity(\n", " trace_space, bempp_space, bempp_space)\n", "mass = bempp.api.operators.boundary.sparse.identity(\n", " bempp_space, bempp_space, trace_space)\n", "dlp = bempp.api.operators.boundary.helmholtz.double_layer(\n", " trace_space, bempp_space, bempp_space, k)\n", "slp = bempp.api.operators.boundary.helmholtz.single_layer(\n", " bempp_space, bempp_space, bempp_space, k)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We create the DOLFIN function spaces and the function (or in this case constant) ``n``." ] }, { "cell_type": "code", "execution_count": 25, "metadata": {}, "outputs": [], "source": [ "u = dolfin.TrialFunction(fenics_space)\n", "v = dolfin.TestFunction(fenics_space)\n", "n = 0.5" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We make the vectors on the right hand side of the formulation." ] }, { "cell_type": "code", "execution_count": 26, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ ":3: NumbaPerformanceWarning: np.dot() is faster on contiguous arrays, called on (array(float64, 1d, A), readonly array(float64, 1d, C))\n", " result[0] = np.exp(1j * k * np.dot(x, d))\n" ] } ], "source": [ "@bempp.api.complex_callable\n", "def u_inc(x, n, domain_index, result):\n", " result[0] = np.exp(1j * k * np.dot(x, d))\n", "u_inc = bempp.api.GridFunction(bempp_space, fun=u_inc)\n", "\n", "# The rhs from the FEM\n", "rhs_fem = np.zeros(mesh.num_vertices())\n", "# The rhs from the BEM\n", "rhs_bem = u_inc.projections(bempp_space)\n", "# The combined rhs\n", "rhs = np.concatenate([rhs_fem, rhs_bem])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We are now ready to create a ``BlockedLinearOperator`` containing all four parts of the discretisation of\n", "$$\n", "\\begin{bmatrix}\n", " \\mathsf{A}-k^2 \\mathsf{M} & -\\mathsf{M}_\\Gamma\\\\\n", " \\tfrac{1}{2}\\mathsf{Id}-\\mathsf{K} & \\mathsf{V}\n", "\\end{bmatrix}.\n", "$$" ] }, { "cell_type": "code", "execution_count": 27, "metadata": {}, "outputs": [], "source": [ "from bempp.api.assembly.blocked_operator import BlockedDiscreteOperator\n", "from bempp.api.external.fenics import FenicsOperator\n", "from scipy.sparse.linalg.interface import LinearOperator\n", "blocks = [[None,None],[None,None]]\n", "\n", "trace_op = LinearOperator(trace_matrix.shape, lambda x:trace_matrix @ x)\n", "\n", "A = FenicsOperator((dolfin.inner(dolfin.nabla_grad(u),\n", " dolfin.nabla_grad(v)) \\\n", " - k**2 * n**2 * u * v) * dolfin.dx)\n", "\n", "blocks[0][0] = A.weak_form()\n", "blocks[0][1] = -trace_matrix.T * mass.weak_form().to_sparse()\n", "blocks[1][0] = (.5 * id_op - dlp).weak_form() * trace_op\n", "blocks[1][1] = slp.weak_form()\n", "\n", "blocked = BlockedDiscreteOperator(np.array(blocks))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we solve the system, then split the solution into the parts assosiated with u and λ. For an efficient solve, preconditioning is required." ] }, { "cell_type": "code", "execution_count": 28, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "/home/betcke/.conda/envs/bempp_default/lib/python3.7/site-packages/scipy/sparse/linalg/dsolve/linsolve.py:296: SparseEfficiencyWarning: splu requires CSC matrix format\n", " warn('splu requires CSC matrix format', SparseEfficiencyWarning)\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "Number of iterations: 272\n" ] } ], "source": [ "from bempp.api.assembly.discrete_boundary_operator import InverseSparseDiscreteBoundaryOperator\n", "from scipy.sparse.linalg import LinearOperator\n", "\n", "# Compute the sparse inverse of the Helmholtz operator\n", "# Although it is not a boundary operator we can use\n", "# the SparseInverseDiscreteBoundaryOperator function from\n", "# BEM++ to turn its LU decomposition into a linear operator.\n", "P1 = InverseSparseDiscreteBoundaryOperator(\n", " blocked[0,0].to_sparse().tocsc())\n", "\n", "# For the Laplace slp we use a simple mass matrix preconditioner. \n", "# This is sufficient for smaller low-frequency problems.\n", "P2 = InverseSparseDiscreteBoundaryOperator(\n", " bempp.api.operators.boundary.sparse.identity(\n", " bempp_space, bempp_space, bempp_space).weak_form())\n", "\n", "# Create a block diagonal preconditioner object using the Scipy LinearOperator class\n", "def apply_prec(x):\n", " \"\"\"Apply the block diagonal preconditioner\"\"\"\n", " m1 = P1.shape[0]\n", " m2 = P2.shape[0]\n", " n1 = P1.shape[1]\n", " n2 = P2.shape[1]\n", " \n", " res1 = P1.dot(x[:n1])\n", " res2 = P2.dot(x[n1:])\n", " return np.concatenate([res1, res2])\n", "\n", "p_shape = (P1.shape[0] + P2.shape[0], P1.shape[1] + P2.shape[1])\n", "P = LinearOperator(p_shape, apply_prec, dtype=np.dtype('complex128'))\n", "\n", "# Create a callback function to count the number of iterations\n", "it_count = 0\n", "def count_iterations(x):\n", " global it_count\n", " it_count += 1\n", "\n", "from scipy.sparse.linalg import gmres\n", "soln, info = gmres(blocked, rhs, M=P, callback=count_iterations)\n", "\n", "soln_fem = soln[:mesh.num_vertices()]\n", "soln_bem = soln[mesh.num_vertices():]\n", "\n", "print(\"Number of iterations: {0}\".format(it_count))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we make DOLFIN and Bempp functions from the solution." ] }, { "cell_type": "code", "execution_count": 29, "metadata": {}, "outputs": [], "source": [ "# Store the real part of the FEM solution\n", "u = dolfin.Function(fenics_space)\n", "u.vector()[:] = np.ascontiguousarray(np.real(soln_fem))\n", "\n", "# Solution function with dirichlet data on the boundary\n", "dirichlet_data = trace_matrix * soln_fem\n", "dirichlet_fun = bempp.api.GridFunction(trace_space, coefficients=dirichlet_data)\n", "\n", "# Solution function with Neumann data on the boundary\n", "neumann_fun = bempp.api.GridFunction(bempp_space, coefficients=soln_bem)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now evaluate the solution on the slice $z=0.5$ and plot it. For the exterior domain, we use the respresentation formula\n", "\n", "$$\n", "u^\\text{s} = \\mathcal{K}u-\\mathcal{V}\\frac{\\partial u}{\\partial \\nu}\n", "$$\n", "\n", "to evaluate the solution." ] }, { "cell_type": "code", "execution_count": 31, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAj8AAAHwCAYAAABJ+g7LAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjEsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8QZhcZAAAgAElEQVR4nOy9a7Bty1Ue9vWce+2z70PmGmRHT5CdgGNDQQBFGBOX5TKOQRYIEmJjHF55qCAogQopOyEOdhzbwVUpV1BkUBQDijDmkSCwICJEJCZI5QgwsiRbyA8F5NJFCliSdaX7OHuvtWbnxxxj9Biju+daa5+9z7537fFVnTPX7Nnds+dzz/H1N8ZIOWcEAoFAIBAI3BYMNz2AQCAQCAQCgfuJ+PgJBAKBQCBwqxAfP4FAIBAIBG4V4uMnEAgEAoHArUJ8/AQCgUAgELhViI+fQCAQCAQCtwrx8RMIBCqklF6aUnpUrb8npfTSa9rXX0opfTil9P9dR/+HIqX0+pTSX7pk27+QUvqb93u/gUDgMMTHT+AZh5TS+1NKT6WUHlf/npdSelFKKbvyx1NKf5LavZ62f7nr77+n8m9Y2GdOKT1B/X04pfTDKaVH1PafTynddfv9Kdr2Umr/Rtfn51D5zy/s95T+oP5T2v/7U0rfn1J60WXO3WWRc/7MnPPPX3W/KaUXAvh2AL8v5/ycK+ozp5T+FVd26Y+Sm4D/+AwEAleL+PgJPFPxZTnnh9W/D6ptj7htP6q2/RMAX88rKaUTAP8OgP93j31+Ts75YQC/G8BvB/AX3PZXuf1+mdr2zwH8gZTSp6iyr6fxLOF/BfDlAL4GwCcB+BwAvwLgj+wx3mcCPg3AR3LOv3VoQ7p2gUAgcDDi4ydw2/BTAL4opfTbaf1LALwbwN5TLjnnjwN4E4Dfd8B+LwD8JICvBoCU0gjgTwD4oV6DlNIXA/ijAF6Rc/7lnPMm5/xYzvmv55y/j+o8L6X0ppTSR1NK70sp/YeqvZlGaUxlvT+l9F+klH41pfQvUko/kFI664zl/TQeZlF+LKX0hpTSJ2hK7MWq7uellP4+bftfUko/2prOof7eAuB5xJS9nsq/nPr8GDFqv9eN48+mlN4N4InLfgCllP7VlNJb6Lz945TSn+jUe2lK6dGU0p9JKf1WSulDKaWvSCm9LKX0T6j9d7hmpwvn5vfSMX2Mtn25a4uU0kMAfkadF2Y2P6bWnyCG60WXOf5A4LYjPn4Ctw13MX+4fDWtfx2ANxzSAX04fQWAtx+47zfQ/gDgjwF4D4AP9qvjiwH8Us75Awt1fhjAowCeB+CrAPyVlNIhrNCfprH8ywA+A8Cf27PdlwP4EQCPYD6frwHmaToAPwHg9QA+mcb3la0Ocs4/B+BLAXyQmLJvSCl9BrX5NgC/A8CbAfwU9cv4UwD+OGaGb7P3kRLo4+ItAP4WgN9J/X1PSukzO02eA+AMwPMBfCeA/wnAvwvg8wH8QQDfmVL63ap+79ysMH98/x+03/8YwA+llH6POy9PuPPycM75gzlnYTQBfDeAtwL4jUOPPxAIxMdP4JmLnyRL+GMppZ902z6stn1MMweENwD4upTSJwH4Q5gZmX3wjpTSxwB8GMCnAvgf3fZXu/3+N3pjzvnvAvhk+mO3z0fXpwD4UG8j6WX+DQB/Nud8N+f8TgB/A8DX7nk8APCanPMHcs4fBfCXMX8I7IO35ZzfnHPeAvhBzNNxAPD7AZwAeHXOeZ1zfiOAXzpgPH8SwP+Wc35LznkN4L8D8ACAP6DqvJrG/NRCP+/Q1wLAf662vRzA+3POP0Bs2jsA/Djmj8cW1gD+Mo3nRwA8G8B355w/kXN+D+aP2M9W9ZfOzcMAvivnfJFz/r8A/DT2P+cAgDRr2L4GwL9NYwoEAgciPn4Cz1R8BVnCj+Scv8Jte7ba9kjO+b16Y875bZhZhT8H4Kf9H1GajuDphT+oNn1ezvkRzCzA9wJ4q5sm+k/cfv+rxrh/EMCrAPxhzAzJEj4C4LkL258H4KM550+osn+GmaHYF5pV+mfU5z7Q04RPAjijKajnAfiNbDMmLzFXHs+jcQAAcs4TtdfHtE9/n6evBYDvUts+DcAXuI+jP42Z4WnhI/QhAwB8r/ym2v4U5o8axtK5+QAdE+Og65VS+lzMTNJX5pz/+b7tAoGARXz8BG4r/iZmL6OKfSHPJp5ueGtj+xozw/K7AHzWgfv9QQD/EYA355yf3FH35wC8JKX0gs72D2Jmkp6lyj4VZSrkCQAPqm2tP+4vdG2XpuH2wYcAPD+llDr72IUPYv44AQBQPy+End7JvtGB+ACA/9t9qD6cc/7me+x3Fz4I4IUpJf3e1ddLozrGlNLvwPzB/Kqc89+/niEGArcD8fETuK14NWYx8S8c2pDEyt+I2eL/tUPa5px/HfNU23+5R92fw6xN+YmU0uenlE5SSs9KKX1TSunfIy3Q3wXw36aUzlJKnw3g30cRUb8TwMtSSp+cUnoOZh2Nx7eklF6QUvpkAN8B4EcbdQ7B/wNgC+BVNN5XAHjJAe1/DMAfTyn9EdLIfDuAc8zHeVX4aQCfkVL62pTSiv79643p0avGL2L+IP0ztM+XAvgyzFNpHr8J4FNoapY9234cwA8578VAIHAJxMdP4BihvWIeTyn9p75CzvmjOef/003P7MK7UkqPA/gXmN3Uv5K0MozXuP3+SquTnPPbnGv+Er4Ks+j3RwE8BuAfAngxZlYImPUiL8LMKvwEgD+fc34LbftBAO8C8H7MItvWH82/Rdt+jf7dU5C9nPMFgH8L80fYxzALg38a8wfMPu3/MbX5HzBrq74Mc1iDi3sZl9vHJwD8m5hF7x/EPE31VwHcuap9dPZ7gVkM/aWYj+17AHxdzvkfNer+I8zC71+jabmXYBZXf5u7xz71OsccCBwr0mHv/kAgcCxIKb0fwH9ADNN17ucXAbw25/wD17mfQCAQ2BfB/AQCgStFSukPpZSeQ9NeX4/ZE+p/v+lxBQKBAOPGPn5Io/BLKaV3kXfNf92ok1JKr05z8LZ3p5Q+7ybGGggEDsLvwTzd9hhmzc5X5Zy7LvuBQCBwv3Fj017kxfFQzvlxEja+DcC35pzfruq8DHMgsJcB+ALMsTW+4EYGHAgEAoFA4ChwY8xPnvE4ra7on/8SewWAN1DdtwN4JKW0FPckEAgEAoFAYBE3qvlJKY0ppXcC+C0Ab8k5/6Kr8nzYgGaP4rAAboFAIBAIBAIGN5oVmaKm/msppUcwxzL5rJzzP1RVUqtZq6+U0isBvBIAhvH08x/4bb8ToCm9pFvQbynjaT/fa+6uHIbkDiHZ8tzcRnvltrJe15Wywa0f0KbsN9t1VZaGeTnwkBKvlyMYYMuGNNl11G0S7LZyyNksNfxN0bs6WdXk31O2532ikzCpE1XqzMttdm3demubzCZzncmtq0GUe5HWJ7u+T53qfm7Wcfe6VFVtGv001w+BfwZUWTkd9l7XF1meg8Gtu/s5a1Nuz+fDtJO67jng+3dQ97q7/3nbSCddlurc+m0DeD2b9eYzJc+DLddIrfOsQdew8TrEBP9c2Ht8Unby1j0zWzqBXL5VF2Lr62TbdprscwMAeeKLx4Ox64vPx2TXu88JgDTl9rbW3wT/t6R6PqofC3D3utlUnou75x/DxebJHRf16vDH/vBD+SMf3e6ueCB+5d3nP5tz/pIr7/gecaMfP4yc88dSSj+POcO2/vh5FDY67AvQiUCbc34dgNcBwMOf/ML82V/8rRjW9OJYq5fJhsvmpyLRetrS+sTr6iae6AnSQelbMC9fuonpayGPAy1p/cSuA8C0ohfDCb0YVna5XakhcRmle5xO/TrXK222p9luo/UsSzoHp+UBGFdz2ep0zh95ZzUvz07nlEJnJyWv5AMnc9mDJxdm+cC4NssHhxKy5c5A/Q3ztlXauuW8Xf9BGN2F2DoCk1+s61xu73UeAQB36YSc0zZef0qd3KfoJHLZkxtezuV3eX1d2txdz/1dXMzLDa1vL2hs5/P+07pc7+F8/j1czMuRlgNFxBlVZBs+ZeNFbq7LUt/r5/ZeH/mepyWXY1POZ9rO5z1x2ZbvfVryH4Gp/5LP5QuZBqKuDz8Hcv/P5wW0zs9AXqk/oPxcnNIf0Dtc197zmzvqWaKyLZVt7/jyMn7Zdmc+xkzLRMuRnoc7ZyWN1gN0/z90Ol+Ih0/ni/bbVncBAM+i5W87uSttPulkzozxSeO8fBYvh7nOQ3ThHxxKSKSH0tz/mTwPE63zR1HB6Y6Pny1dO/3njW+Xu/RRcs7PCT0fT+T5hD0xlRBIn5gemJfbObPLY9s5iPhjm7n847Scf8/tPn5Bbdbz+hMX1O/5vDw/L8/qhn7nu/TMnM9jG2g58vNxtxwvPyvy7NC9z3VPzu1zousMF3SdL9zzcVHOlPx9WNPzwX8fNlSHP4amhT8SKbWXJ6NUkb8PqxFv/9XX9fu6Bnzko1v80s9efdio8bn/9NlX3ukV4MY+fihU+5o+fB7AnMH6r7pqb8IcKfZHMAueH9vLayTRB4T/+gfAWXXkBc33nZiP3IU2Efgh46//zotfvVW4u+QsBs8M5AYTwI1Tl6LRle3YStvWAF1dtu5k+2C22kPa/1bxVjGDLVz/8WIbcx9sHbPJrtpQdgDup+qPtq/QSPhN3Y2Z+6/H5JmrXcelyxyxBz5vfB4zyouuOgueBVGoGJLq+i9hPmgeNt9NE+z63C19CA+d3unlrlkW/yFU1mumoVjQ/MfClg+01F2OzpovlrqnfkqbrZR1ngvz3PGYBtOWn81NnwiQcU6ODSzMSTlRNVPSWVdnfaJ7eUt/1c/kTqKPoQaPw3fY6D6GZF2/v6TKDstuD4GEZ3p12S4k1YY/F/npFfJUxsqDab2jpcDVWXpO+GYnI4RK9ciHKsQmfQRxbf4I0h/65SVvlx4b9ZEle9923uHXhwxg2mnhHw9ukvl5LoD/mVIFDAB+LOf80ymlbwKAnPNrMUe2fRmA92FOEPiNNzXYQCAQCAQCx4Eb+/jJOb8bwOc2yl+rfmcA33Jw3wCmEyCRqTCduI2AmBMyvc/TT2TK6ekonhoTg0PM+P6neXIWZyaLWaxVmW8ulsIgH91E69MasxLJWKue4eH1fSydZQYoKyuKmQo+5PlbtVhqaYEF8RgXTBluI4xM9myOqgs7Lz0665LPl54OE9bG9T84HUarH9FkpNPu+FvnwVaYz9tkLiFf5x7Nkpo/bcHhDJBfmxpl0juZ3/5+tli2bHnKAFCMK7MtssH2Mai+Mk3FeC0GMzVNPYewvva52DZYovIsMY1jWRt+VLfq2M8rjZdjfhZ0YaKZ4edutJoZzRZ5GmJKNPVGVbbm6jFFxsdFz1SPAVJ1+L4828fy38UC7fFXpcUSeQgD5CVArYEk+xzsw5Q2b2UAfNKHxjTisPFlHQaohV0MkGqfMC7XuxZkbP1794gREZ4DgUAgEAjcKjwtBM+BQCAQCARuDrPm536zTTeH4/z4EcFzrVQU4bG4WNI6c6BjQ2BLVD3PjHgx8z70JNfNLD7lDZuaZpSph2zp9/3AtO/lp780IcgCOJ7+qljfBjyVvSgUbri/6/KR5iX1QznJdFRn/57+Ry2OLks3pbkwJn5axsbx+GP0oQB4vfgLqdlTEXS662ymyDoHmw653lxlMPsd1SYvgu5OGjSmsopbsTtv2v2evWR8qIfBCaEbz9TghdVSl56pSYtN/bItgNZ1kpsSK1NydvoLAPgwznmmzLlsy24b017ejdwLnvVUmbiPD7zkm2VenKVa1M/PBznDiZdXTwgNzBFmAchBn15i+qvpyND5C7M0Be6nkMkpUp6devoLC4Jnvz136/g/AUsTI/UWN/0F9KfADpj+ut+4TYLnmPYKBAKBQCBwq3CUzE9OFN8mO4sUtbusD6Q1eEtxLqW21NHQsB4B+yUvJrSzVtlV2DNAQIMFovgWbafkJvpCP6AOLMdthMoy9ebfZPXSSWSRrtibOmAhd+Pcvqsgh2pMhWWxAkiOwzM6YfLcZjBt2eL0QeNGJYxmsaoXOA8isC39n3O/ru7JZN3jzXH4oHcd5kdbtWv6vSExdLkVLTNj4a4hkxUty9YFAyzXtX8/eTffnQyQQmZhvmeAlhwDiEHKrm7z0Pnxc923BM/CBrltzO7ocBVJGBjfn1PaKga2CJzn9XX25Q3Bc8X4eGE1rZ/0n3MWx2+dC7ytw/vk2ECuwsL14D2XNvcmgF4SNF8WngECGqwp11lihFrvuy7a/Vel6j1V2P1LMkD3eQYqIwtLeBsQzE8gEAgEAoFbhaNkfpCA7apofrSLeC9NwOTYndTQBLD7O2/LzrBufsFXZA6NSSLmamFKzwKhsTVYFrGCvZXPAeaMVdO22+sUA/W8ODM+nhjbpjLHzcGLi96Fxlilt2gEFNxapoTZnDX1r5kZ1unw/DTrInxANa1B8Kd2K+wNa3/KozBSzAGOEF0CLtLYTogJGopF12O3fFqC1jWEMEDzaondVzNw0tqzORXjB3j339KamaUFBsBZwz0GyO2RtjkGyGxsW7uip+NrpqQsPiCi6OeKsKZ1CFSXT6Z1i7eCKsfk+u4a5SW4IYUwkKjirP2pu6pZoYZru0OdOoJYo6HRhorOXHBPdoeXgIhGC2fh16+KAfJohZbYBX+faf1cySLTG0Tr3boPn+mxJwPU2vMuBsj/vgHxcQieA4FAIBAI3BpkFKPwNuAoP35Y81Ppe3ijgoT8Z0tuVX/D9zQ3iT2h6IPeJBYcrLUq6GiB5m3ELPGgnAaoFaLOB0KsbPQGi8NlNVvE1UpBFhaHrFXRDbEWqDRj5sLrW7zn0zg0mB+yBFcUkVJ0PI5Rmdu7NtjtGVGlqhDGh1JlqPO0osh+a7qwzDox47RKKzM2O17LCvF6MxFsdb7m9Q2tT4NKhSFMjw3Fn51uwSaytWV+WdgwHSyOi9z9mZJeuHvdlRHxcBAD5Fja1iu4Z6cPDeansEWe8WnY6DtSYRQtUJ9BZvaGu+C4qHd1ws5eIMQqiWhpU+WtG11gTL1NPMPOqS3l3+PnI7k8VCj3Ffc2uiXD6oaWWZux9W7zddwVXtIGtZK4esgzQ3/S/PmptZAaTxMGCNjPEyxwJTjKj59AIBAIBAKHIaa9nulIc/bm5C061N4f3RjnRvPjvK54G8c34Wrqu180PQcwQF1PsMkPGhg4zoowMHNtyWDf0I3UHkId1qDBHnjPiKnRP3tcbIbRNuUxNzyevFfU6JieE9HfFFuUPcBWxMxsJdEpLy3rMpdx/z5tBvWhdBdbTnDI8VF4fy7jvO6Ly/g4VuIpZo/DsEVUNrpt56wBUmbkljUNct4tA8T3kdUJ+WX7WtrrbfVA/p7YR9Yh3W35vlIeMFt3/jsxgUxKDL/Nt2WPPdUNPzI7k6ICkibDvxsqxkcNtfYEY/bGDm3TYH56yybz09nmk6KaMmaS+B5xHmE6JQbH86k8whz0ub0uHVC3i07cMA1+ZlgHVDFAzePz9//VMUBLWw5igALXhuP8+AkEAoFAILA3MnCrXN2P8+MnUTLThmdHnaDRag3a2gCrpxm8wUPaHG0nCHvjGaB9hu+8WSQqbiMFZdGLsLaExkgmkNWA8NIzPY410IFyuWhwbWV9qOpOdKzsczLQCbvYzOyN1vycDOxRRQwJa2Yktg6xLMZzy7JBRZNjGZoWBmetSrwfddAc7bZofOYjYabJ63v0Prk/XmftDx+Ptlo9G+R1Qefqnrmgsq2LCcQRfznK8pLmZ3TXv7ruzbpeO0EMo5b88Aqb3aIPshogPe59GSBTtxcVGjV80taWkq/syzG6/r3RiApdMceiWbJeWTrCM+uApKnohXhpWR4NHwuIGZ9prI9MvMc6SVFbUaG5kmh/9mD6eNvpPp5bPKRLeHnto/nxqBigThygGf58ewaoVcfjihmgG/oIuT3xnSPOTyAQCAQCgVuG42R+AoFAIBAI7I2MHK7uz3RkEjwzbNj7jqgtu3UFSTAqJTS15HkzlflQ3OAleGLbtbd9AMyD81isuJk6nBcbnjphcWM2TQfNcItYlpeW25aprMa0l3etHpzweW7PUzAcgn8uX3sBtBb9eoGzmy7yru8AsKJtPO3lhc885aBFnJ5ul9QVIoSuqWeuw1MMpyJ4pmmwqZxcL4aW5dSfKitTfHZqTM6Jmh6UxKgihqYpKzq3QusvTHMyJN4fC7v1A1KlxOA23H/jvuWpMN+HG8Y8cOpPBlMLnHV5q0xSYrg+m9NfMoXsnLeNeNmFsPC7bqTE2PopMv9eaSQl5iPktixD7gmf9W+/ZBf4qXHU3j2+SoraAPezkmiT1nFC91iuMx/AvNhr+ovhfemvCV0BdBO96S9gfxH0AdNfIuBvJEMdhj32FbgXHOXHTyAQCAQCgQOQlYfjLcBxfvwMwHSasfzl3HNrbH3hc/oELzwe1P8oika9FxaVkvBZCAa2jhdD81NbGktSzFLlQr+xolkhdRrMDBvxmcfbCZg3l3Fb6o9zWAyNutRxYRbm5ZbYCWmqGI27G+fqPnCQQMeCKFOdE4wyu1IYoJkCYBZkMm7AvB/L+JRkqLXVupK2NrCctFGC5FNin+4SKzQ44XMRQGuRNNUZrJt8OQeFvuTreFfOOzFAcn1IEK3oSLkeLGIWhs+zO/oaWhF0LYCmeoaF5EFydwtMEIuht/QgOBbH3OM97GCA9JA8ZEzN546dFGxpcX1vMAEdkXTtWAF5MCbpb379XnRc33W3PXd4zRKVVBg9logZWfVckBu8T4mBDgOk8XRmgLw7vDBAKh3P1N156+5Z+vvQHcWepYV1FhZosw3i55pxnB8/gUAgEAgE9kbG7fL2OsqPn5yA7SlQWWdNdL7om2kheN3pH9ZcWruii0aGtQ4ugtrSyESTU3JYFEjgQ+dqO9n9TcrqHjbWisve5Vk0O+XYBylzgxPXal2X9snb2OIUVoJ0K2uVRFTcuplNmbcxQ+Jd34Hi/s5pJgoDRG7zOZtyoJ8Cw2t/5kOzrwDuh5kk1vxw+gsAuHD6I27zRLpj1jWD1XOHbwVETI4V8u7w64Fdesu4902JYXJ8Dr17wt4bRifEjBKJWAbX79BgfqqUGC4gotEA7XL7FY2carOgA1K7NyhsjRO/Cf1iaruxOZ1TKzCi2/tkXwVNOeCuVBhL8G7wwlyZ94gZUqkr7zy+99XzwT/cEIQdzLpsXrmzKwXNfdYAAeXQuwxQap3j62SAAGaB0v06IQapXPdbgHB1DwQCgUAgcKtwlMzPnN4iqxD2+wSrWvqitzqIuiurrQBU4DryxmIdCq+zxZunBtHoixqfqCXImtP+sBnJhpv2QOP8hqzfcRoQYaf0/joBEDOzXY0ElyUAIq87BkixBhekBxqcxqfyhBoUi0MeYOfJMkCis5HEpMV6Eg8ttGFTYdRs0Hw85FUmnjFqTLQvZoN4DLzfJ4kBMoERp43pR5aOCdK/TzopMQY6txeKtduKJ9i87r3xJCWGEoYJWVAFs+yzRSUlBrVZT7aOukfkp7un66SoCo7VrJKhsjzFxCWdTH8yggaLNFC/LYK17thvsxoffueU1BilavJsjSRDZVZnxlq1yTsYHqNrcx5mnh3yXmBA0Yj5pKhnhiMBVjolhpzweTG6IZp1pxkq/bQCLXIH/U1XCdYg7mSAAPQ9wa6RAWoyT9eHDOOwfPQI5icQCAQCgcCtwpEyPxn5dFKWTsNLo/tV3WCLOok/vVGWtAfMwN5FtO49qyb2WDHDdj8uAedlkpTvovcAE8YnWetMMwGljNadJ5eJCcSEAscAcsySxP/RqT4SMz/2fHmGQ3tJFXaoHVvnLjNByspnLc6pnPC+BmHseIZJX9JDOU+n1N8Fj4WYoFO3X8MW9WIDtfRBg2XAxk5MIB1DqRsTSK6lu6/nDkxZzQTVLGg3JYb0pdQufA+ymxpvW7PWgaDFS+zVxcfmmSDerVqX+1UyjDotRcMjabiwTWrUtiIPadNjApp92W2S3oLGphmbtXvJ9Ly/NIT5kUSnHP+qTonRjXvjixuv0F0MkMZK3rtVbqB+o/seC4iXN8MAAYWFBHZIVa8Jt0nzc5wfP4FAIBAIBPZGRnz8PPORAJxOJfmj+tIuVmj7S91HNJ7Lst3mbxC2ZvUHvSQWtTGC2Dpl7Y+J2lxMwHnTAgFUWQVsWPW0QCiaA87jyN5fzMwwKTFstPVB7JDThxS9UOlfEly6usyIcd1JteHDX6+ZbWGvJors3GB+ih5ojoPjtTItzQ//Zq8s9szjZSsZqkR/hh0bQ+uHJu6HGKA1WY/MBA2kQRlMnJ95nKfJan9GF/8HAE42D9hjd+elnLcyxqd6MYFY6yBaHXUgnVhAhQnics2yWH2Q5NUUlrARE2jke4PZR+tdJINVdZgBYj2Pl+I0YwR5GofjbS0QDoxlTQBp+Fz3JfKzZWCbaHiEAVbHw92s3fpiNOgdHmGtqNAyfmGLqA+OmZWKBogTo3ajQiswG8QE9HgAAyR/iPdggAZ5VtvO2mkPNv1qGaBWHY8FBugC913zc9twnB8/gUAgEAgEDsI+IRSOBSF4DgQCgUAgcKtwnMzPkDGcboW+1B+zkxNBdwO+qTY+GKB34fXpHOYyovPXduqhSkSqU2Jw/zxtJE0vIYDmg9faZQ6smGwdmV1pCJ5LKgxeb09/mX05cSyXS/BDRedO1MHGCZ9HF/RQTzmx0JnLRAzs0l7oqawLCoDohc98P+ipBjSmwIBa+KwFyRxMkqfCJHUF74/GeKoCI/bE0LLclrq8r547vCRF1cLwbkBEcu0mAXQezNzlXOaE7jLdJfd6uYajT3XSuw/UNu5/cCr8kvxWC55dGYuht3a6pemR7p8dn0AVu6fA2gERedlQ/s9beABVmYtzWk+NqSZ8X3LRxrVpucL3xNDbBZG0THf5k7hoHrdTYpjnxJ3+wW9opsLoXJB7EED7Kesl+Omv+Xdn59W5bCjDD5j+4nf9tE+zK0ZofgKBQCAQCNwqZKRmLKhjxVF+/KQEnJxuRTOpvWb5W3wqJiYtaXtLZOYtWVm3LIh1GW4vxR1+YO5jjkUAACAASURBVBZGuaKL+7sVRzeJH1+2zwe7E0Wz8LkKfqiID7baRczsAyE2jllc3h1LILpXk3xzXk4iyiWhsDs/Op3CSELnwvzMfAszQqOImVUb+s3iYklHQeva4pEghjxcZzWKuDLV/TNWrtwzQfPvjV0muzxVlETPHf7EM0FaJJ3OzDiZCbrrAiJyMERg/4CIJnOCu/97wmcAchNI/y4gohBBOtwCP8hrToZqxdKeAap+t2DE0fuJoJdD1Lk8EU1WZ5cI2jJDuq6w2Oxav6DN6AVGLEJo5eouQRmpzP1F4Lo6Ye52sP1vfdgIw+bQ+8ON5aYYoENwUCqMe0iGag99LP+H4PlacZQfP4FAIBAIBA7DbRI8H+XHT0oZp6cbYQ1aX/CSdLPD0BjNz4VnfGjpknsa1qjDComRTYNKyoqS1Besf3BMkLZme1PYVeBFYw3TQrzsmRrj/VC56mNwrs3e9d2EBHCMT639IcvTBEbk82MDIK7pRKWGCzczMd7dm3UvXhcDACvW/EwUAJHuhFNhi0r/3J4ZoJFO6uhYKJMMFW2Mvq8GGyXj5dQY8AEZNeND7NDAS6t/0mzX4Bgf7w4vKTHUzdQNiFixnuq6009heFwS1GZAROnXBkQs95t6Lrw7PIeJqLRAyoTupMJoMkKij7NinOTdprUW7sJ30mGADJgac0yAsEQ8jJZuhLVp3MVcR3MjWXRAMHV9ctR9wIEQF6dBaJNPhaGxFcaV2dhmFzBU9r6BEO9j7s/93eAPZ4Ca7/I83EiQw9uEo/z4CQQCgUAgsD9C8HwEGFLGA6frZmArlg9shclgBsh5uUBbq21vFqnb8Grx2gipu7FtdEDBRJQJMz2caLQYDpqOYrOuOsQZDa8c2Y9Yj9wHM00wSz0Wsb6d95c+x4MLdlfpNxwbNtex552Zhi0HP2RGQ+tr2BNsszLbmPXgtBS6jfeW4tQXpxOt6+CDmefdbQqJXlLUuY5d9xYuh9xcKX3EWnRAc+1TZqM4hUUuY2Kmise0KzUGUPQ/vaSovHxSHfs56aeYAZrSCS1toEo02LtLBUTkuiNvrwMv8jPELG2VHJVv2K2+r7hsso2WtECSNoOCKHIgQxl/7QXIsf/K6XAMUNMFzTIB5fQssAaSOJUZGR5jYSD4VbLLr6kZGJG1P8wwLQZCtMe0FbFgXXflUmDsZoBU5WNJhVFd313aLwAY7ru3F5Bqb78jxu050kAgEAgEAgEcMfPz4GpdUkqobWI1klfR5KxVCBPUaOR0L967BQ3Nj/ci854wJj0Ex9JhFoS9pLz2B2jG8Zl329AfdVBSYXCfbPmqOp4VYlbKxX0BChskxE61rBmyXgwgsSZd/B9AJS51rBB7Pp0Qe7FSByKeYVur8blLGppBsSzMuHhNDlPCg9MAafDp8PYgr2tHmdGxQcwEMSulx8Rj8DGBRqdzGhrapZUwP9ZD7GSYvcGMxyEtz1kvJ9fHxQTSHnvuObhMTCBhSn0sLUBOGt9zrFXzyVGTut6ia/N6oB5TChSRDN96/CxxSgx9szsGo1b6+C12K3Xi1he8guQFZTVAW8VkeX7qnMtdTKDLiFr1dIhvz7G62nInpojtIL1uzjgPyq+nMQPkiPdLpcLw2i9Y/U8eDr9O94KMhUS314iU0vcDeDmA38o5f1ZjewLw3QBeBuBJAN+Qc37Hve43mJ9AIBAIBAI3hdcD+JKF7V8K4NPp3ysBfO9V7PSImZ8LsWiNLoXZAtI5sKZkI5Yt61J0LBoXm4StV6d1aOkgCsOTm+XGq2ywbVkPVOLwKAuh8gSDXTZQsUGTs7LJotMBjkvsH2+x1wzT4LVRjt0SfYSJ+2LPoddeMQO3Nk34uq7M8Ivmx3plzWVWMzOKLsgmOJ3rECPDnlR5Zcq5X20NC1vjTvHoznnLPmSPmJJIlTy51JhYD8Q6HR8TaOViHM3jJ1ZLvOMsE1TOkzr2PWMCmdhZ8hwIZUJLd6/PO53LKraItkt5acRjkufCxeiqvMGAIu5zeiBhh/aJmp7t85fVufUxgXYzQO2tQDldm1ZUaH4mvW6kkTiVtTjCANG2iz2YHs8K+eSoRgvS+6vRILs4PpB4hMkh2vhO+pi9RX5IMlTV6L5ANEBy3pb4hA4D1EmGehPeXjcheM45/0JK6UULVV4B4A15Fqm+PaX0SErpuTnnD93LfoP5CQQCgUAgcF14dkrp76l/rzyw/fMBfECtP0pl94SjZH4CgUAgEAjsj5yvzdvrwznnF99D+xYdtQd1u4yj/PgZ0oSHV+eVIFb/HihhZhFAE62fGrS+EyDXYk1aXtSUpZ/KKtNfvF81deJc2yVoIE97KZbXBySsUlZc5tZgdlkHUxRhNU+JtV3f9biraQ8vfF7r88SidBazcl07/bVVB7RO9told539OlBPifE0np8OA9TUDwufaXmXHpci0sxVG4+awtfTFDwGEvv6sTbGxNNrHBCxTMXVx7FyARz91JgPFtkq6wVE3OiwDvTMbOW+tdNd41A/FxIaoHou6Fw07qu85j6s2724wpsQCvRDnhMWL9O6dh7wU2C9df1S6ARELAlZWy7dvcmxuc0JXaetme/oiaD7UyackoJntXkWm+My6vQXIobG4ZhGe3Qybj3l1Iv52J3+KqPxZ6s//QVcRyDEqXGeeigu8OVPqqQRkfcqb3F99a73DUx7LYU4uEE8CuCFav0FAD54r53GtFcgEAgEAoGnK94E4OvSjN8P4LF71fsAR8r8jCnjWSfnJbibshT4tzM8SyqMgS3b0t9WwvSTRTVYy7O4viv7yYs9ky33rvBzXWsFF1de21SPW2KpuVhu3n29CafOLcEPVdnWbhTGx5UD5TzUjIxjytRuhQlxbJGkIZFzWy7IRs6lZ3r4sGrmh+HLCmNSs4MlQSqNhpkG1P1z2Trb61rYnBqrZEsnJ5pupcJgMbSs872S6+PoucOP4vpO5SoNCLNOvM07DYjzgE6J4ZhK/3wsBkT07FBDSF8HTZT8LNQ/s5BKkMz3Qkf4bNgDZm92iKD1ds5Py4wPB0RkJnFqqH/7bvD+JdEQwDbcoudyxXZVdThg4Yzt8uE1cYhbfCsVBgue5Xy4Ksw0aHF/uTbtAd9vBugQ2FRKFCC0wzGk5rkt98D9jjeYsSOdyTUhpfTDAF6KWRv0KIA/D4opm3N+LYA3Y3Zzfx9mV/dvvIr9HuXHTyAQCAQCgac/cs5/asf2DOBbrnq/R/nxM6YJz1rd7Wh+rEXrQ/4/xe7T+qPcB0R0Oh4RvmgmhS3yC/t1X2uA9DZr9Zb0GXY385howUawT4XRCIK4tw5IMz9sTYqrsO0sqZQC5Txb696zOhrMOQwVZcV9WOtV12B91jksWmlNemxQYXd2nxzPEukghGXnVmezFTZntwXNdYaFJKhFD2SDGhYmqJESw7nzD54JMuwEPxdbt26fE30ezyUZ7byUgIji+l6FrVMsIN/cdn1spIopTOtg6uQLfm6kSe0OzwESXVLUuQrrgTgQItU9xB3eBURk/Z6x+iUNy2SPo8sAAbXmx62bIbo6TngoDFA9fJUMNdllQwOylTqDqXsQ5FBZyNWo02BlNUz5DSdDbZ0DrwPqu8Fr9q6U3X9X90hvcV+QUnphSunvpJTem1J6T0rpWxt1XppSeiyl9E769503MdZAIBAIBI4ZGfMH2lX/e7riJpmfDYBvzzm/I6X0LAC/klJ6S875V129t+acX35Ix0PKeHg8L0kgtabB6R28RctWoLGgyTOsCojoUmFoXYrXAYn1Si4XrSCHVQBBzxKpJKjZaW+KxoH7ctYfoIQ0u7QN5bdUlWl2a3maAJJbO/6a6al1C7V10+ujnFs+397Wu7gmS2n0e2Lr0WQuOKdNzMSIKMQ0scHcmOHpDFwV3+HUGrDWsLBRkhKjoRPKlvkZvaapkQZEvOBcgETvNad/86ViZz6+HSbDblqW1AcMlfuh4SHWC54onmP6+ZM0GtSdPBciDCuD6gRClFO5mBLDjs0HRNQBQwtYH5RNFzUDpA5gFwNk0KtjvcD0vrVW5b7CM0BAzQLtYIDMtmsKhHgvWpj9GSDAXM+n73fDUeDGPn5Irf0h+v2JlNJ7MQcu8h8/gUAgEAgErhnb+z/XdmN4Wmh+KLT15wL4xcbmL0wpvQuzX/9/lnN+T6ePV2LO+4FHnnuGZ413i1fLVKwBTur4hAvx75c6bktylm0SjzDSAokWoYxnGpyF2fX6Km2q9BbOem0nTnUMjMTjgVlSz2age2mAOoyPxOJQmp9J9Cg9/Y7Z/Vy396z1NEC0h3l/MMsFG+9KIUzJ7swFgNPmnKprOFWeYfZkNBkh0YlYrU/R/pSz0PME8zGBWl5ro38+HOOjGT+JQSQpMeblBennOCUGoPVazPTYlBg+AfDcP/1wbNHY8RwDgJHNbY4JNPIzRYyMIhpKegxifHwsoJYbpIffNrnnci404/XeXzUDpGvtYoBa3kPLDJDeJz/e64VD9Jj20IeIhuSE99duM2nPR2JRn04M0L3Ev+EYQWtJMEt95sYfDpW25BZ9h9wIbvzjJ6X0MIAfB/BtOeePu83vAPBpOefHU0ovA/CTmJObVcg5vw7A6wDgBZ/5SQc8woFAIBAI3G5kpBtxdb8p3OiRpjkz5Y8D+KGc8xv99pzzx3POj9PvNwNYpZSefZ+HGQgEAoHA0WPKw5X/e7rixpifNHPN3wfgvTnnv9ap8xwAv5lzzimll2D+WPvIrr5nwfNdceVdqbmfOrT/vK0ERLSB4HRdnka4K6kxaPpLBJ5akMxiaJumIaMthDZ1nBtwJfhsbJNurIevde1mQSdx3JWYeYkvkzrZ9KUb8XTXzukvQ9XDbfP75ZOhB+cEnLTmp7/28VS+DLYiHNXzLG4QnekvG0dgXqzc6eHprqFpm9ipE9/WzhVYtW1vGkwLnr0IukyntaeJW2WDuwcvVN2N3Pd2yhg+9ENjPrQInJ1jALvHN5wHJPAoT2klFhuXysOaAhNyt95poBUE0fuIO4jru3ap5x8bOyXTm/6ax1nN/fAAfK+dMrWeW1vbzxDPCrbSX/gM8Ix9/tD5lBhuGGYQk8xNuimrp4EA+l4g55aW7bshd1+JgavBTU57fRGArwXwD1JK76Sy7wDwqQA4suNXAfjmlNIGwFMAvpoCHgUCgUAgELgiZNxMhOebwk16e70NO75tc86vAfCaQ/seMOGh4bzpyits0GRZIRZtslt8K9mjD/DGbSTIWyqnc+sSmJaw/izwZGpGMwHudCS/VNZqlR7DWsOFhSkYnGWYXOoL2aqHZDXShXjgcW/UmEe/710MkCoCW9ddBXQ90uIXPe+PhYRca0Ex6LdpK7YXtG0v0WOHAdpK22J5Ts4Ur1kcHXxwMEveJiwIu7ObIfL5t6JrGergreKWKLojgNbM6L4pMQC5oSQlBi09Q2puEe8cIAJnZnz4PCoWteMGX56bmoGr3OG39niycpyoAiHuYIB0FbnHJzuWYduiDXtUYs2i7nKDX2RZ8zKLelnwM3PQH9Vu1TYDBOwRCPEGU2H03ieeAQLsE/80njE6Cty44DkQCAQCgcDNIiOFq/szHUPKeHA4r1geoDA9XtMgGqBUW8lFH+SYn46lCxQWqFi4LiAiB3s715ofZ8E2rOCCluWn2jKT0kir4AMi9hggQHnUu92URKp9y5MZn+wTLhrthB+fY6dahrrXAbkw/qJTUH1LfEi3TdiiBeaH1731ukXdhutwQkeuMyWKbqm6OHM6oC0tVy61B+2BmlsGaOJAnlxXpdw4Ew2RVxfY5Kh6TEOlQPDrteZHXN07bvG6Lt8vfNuzBmhKu1Ni1ElQpVOzfa6Tm9tYcpJ1TD3HOnGICZacMMuZFFEgbJCcWscEtTDZH6mT4NSccxcKoJRXP9DX1i39QfPb7DN0v8NINIbSQGM0+6bCeBokQ/UjXDe2TRmh+blmHOXHTyAQCAQCgcPwdE5HcdU4yo+fARlnad1M9jiKpwtbzMQOba0GSFur3iNMNA6cJJXqJWPh8g+y6n1yVJYVaF+PTmC/JhPpdUAdrxDNYPH8+0A0yHTC69Riqq2mlg5IrxveZ8p2G3ttnLix5AVrVcx8z+pUQ5NtbG0n1hd4qkB1t3HMT5Fs1CyO92rh8ycsj+q/58UiLxNv5AOY6Npt6Z47E6vUMzXq2klqirZZaBKoVoZ/3W/dgV0dnaVc0l40PMQa3mN6zIBiSWn9nNY5JcbklnNllzi4k+ZiMVWMY1PNLSjSMepobVlT9gYzmjzWA7lAiLmZzqIDiXbH46gZJs+8LUcn6TE+nWes2XaG9o8FLBt6CAt0L8EBBffAAPnH7lLJUOeGB8O/L3rnQpOFm6zKhuaFujbkjEhsGggEAoFAIHCsOGrm55Qs6otcPttX2TI8rAeqvFwaHmJF22C9yFJD21BincxlnHSTNUDexgcKE+A1P0NDt1OjowFS1gZb5oXloPXR1jUMkO/OaRqs5McyVczIDNqacX0mz96IdgJmLElRAS50S80OiY6n7GcSPdC8LhYWb580S5TsNjovm2k+UVthhEqbdeZtA61TMtzRlrMWaP49h/GvPcEmKi8HcEpVejGBGGbd3RJ8Cs8OYYCKybw3auan9gzzeqAqKaqyeotOju+vTkoMdejFy4u3ZdPHUkwgXvr0MoOKz8NjECZra1ngpveXP91yz9sKWdukkriYmYxuICl9QLyDHetLbc3uzf44Qee6SUlfI9whm/useu3ZN2z7bPV0QDUDJIzRNccCuqtHcN81P+lqmLpnCIL5CQQCgUAgcKtwlMxPShlnw1os85WaRF8T8+O1P8z0cF0dFbqyVoUl6nu1+ASQfskKf61tYIuvb4/v8VWevZXX2MaMjHhjwbZRFp3omDqaH9u/tZTFnGcLt8H8ZEeBsRFW2J2GF43ESXFtZdnQNlA/02T1PMziZHUhSuJBqjMNZt0zQQCwzu06zOoI86Ou4ZrMyIfIE4yZoIksz5XR/JB2zHmEFS1Q/97gbZ41Wh3AAI3Z1tGJf0d4pid36/pt5Vaxz8eFjplFS+8x6ZnGKk4W0I0K3UoS3I0n5GMQ6ePw8baYqWQtkL5vd2g45NRs1flklmtjGZ+rYYB0Wa8tDcmsDaaXEg16YShXiYZ+blcy1BYDVPkT3lA06FY8srtAiaV2n5BxuzQ/R/nxEwgEAoFA4DDcpgjPt+dIA4FAIBAIBHCkzE9Cxim2AAueFUdZXCBpmmukKYZpa8rtFBbV3bI4up38sU5CuBs6wFUhaHdNfwFdujq77Wa6yLe1U1qShmIsbXx6i+oQG+7xagJpXmehKnUy6JD/LoabZBDoBTJU+yxTZMmsl9kcNX3nhNM8Jcr705kL1jIW66rK01+83ExaBOr6pTYshGZXeG1ZsSi6FxDxTE3Xend4nv5iITSf054LvAaP4LQKrTAfyRK24na+0H8zkFwbyU+DNQKGXsi+56VPidEMBtoLFCru8eoe9IlR3X5Eo66DTnKS2LWbfnKi6aSmsKrgnB4uZAOA+vly01/JiMl7F2VJ8LxLBF2XexE098B361300Uv1cCmYy8EvA3qjyjvNTUvq4LVVN73pL+AqAiEecuw5p+rZuG5kpKu9Pk9zBPMTCAQCgUDgVuFImR8vWNbCUU5wSEJnFj4PLGKe108nHearHTq9JXTubWu5w3uICJqWiwyQE/eKgei71+Jlsix9LsYS8M+73harOHmhaMN6bQmaTV1mX1T2TUn8yoYbWVLC6rjl/Nu5w3s3eaq7bbWRQIxMH5D1qlgiEUXLcjDlfK418yNs0CktxeV9PqA1iaPXKuxCcYNnlsgyQRdD4QUfIv5jIrubU2OwEHrFY9eiXOwHbbxy+1O665jROmPLtyk27XTceLv0UmD47Rp8y11IUlRKhSFjrgXJhb3xomgq1/egbLP3eJVktdm/TZgqyUn59aHZIhFDuwNcTInBbJDDls+fqsppM3hMa99oSfC8zACZhMNUdRcDBBQWaCnJ8JVCWLo1rfJLjl8kuqq79+TX9TJAS9DMy3CfBc/A7dL8HOXHTyAQCAQCgf2RYWOXHTuO8uMnIRuWRv/eumBhntURzc9Qf/X7UP9+/arQZYAa1lflnu7ZF2U8VMyMaGWIEWJtjr7/x2zbOphyx8RU5byqGRnWSIzJbPMaHcv8+CXVYfdiyZHQauP6a2h+JscGbYgFlHQXwgQp5oe2MePDbvDM+GxWxPIYzQ/VGZkVukt15r4ewnkZE1nza7Jo16QHukM0winr0tQFZzZo3CtI5gyvg1j17vE65+pBDFB3/00W1RUIAzSvlsdCXQ/3OFSanwaLc0hKDN5WubxzMlRe3+rjmUzbEhixf31KqhjLSkkgUpMKw+5HWO66V/XbMz67NED1NkkO7BggN7z7C3/QjZePhCrpNtUvzz3d4BcYIB/2QspNah2d7uj+Mz+3CUf58RMIBAKBQOAQJBOH7NhxtB8/va9mCVCY2cODvWR4Dp2pANNoXrhAb4d4tfTQUtdnYRG4DperduztIZqfhlcULDMzSZ4Arst12uUAMGRnKbvElHZfvt/2NUhqB9lZssI+OfZG6yR6eqC0tX0p2Ra2nvEpghFbrvqZuB824reWnZqUwGbjPMFEA5Ttcq0CI56fcJoAzwDxenk818P8+8FEKTFcaow1nagzpcRgjzCf5HEfeYJ4PkngSvYyI+tV7eeMxSWOAdq2dB6XeONUgRHdOrO5UypHJikwxAuL72OnXQPqm7lTt80W2aVogNbsEWqOZK7K9zqPzWviFgx+0e3VEhZ4WdZuBkj3cO8MUNlrzZ3sYoAuw3LYtCl7MJRAk/Ri1rTPAGmvMevhuwSdyuYyuK6ZhcCMo/34CQQCgUAgsB9C83NE8N5ZZpvM9yZTV7RASvMjHmH56pieSTyslNeJK/PEick3ysYi607Ee8l6M5kwIY7pYR3B5BgUk+phYE2Oq9PwgCnkkIgAzGCbxp2zrjmkexYdDx3HSdkR63YGEn1MJG5hpmc72Xrz70zbuF/XZlvXrTzChBEitkV7iG2tR9hW9ELE+GxJ+3NaXi6VPujEMj7s/TWXzdsuBqrjUmOckRZoUozGSnRAVg/Ewz7kNcfPx8p5gQEQzxrWJU0kfGlZvqJ52vHmWWICPPMjSYNVmRCHclNaBkjfjNlposoj2dEANcq4DXNPU6GLq3FXHmGOxcmN2FnlwLI5HKOfkzq0a9iC+80A6X32GKDLxLLx6VQAdO+nwZ0M8zfBHaJngHTyW5ZuiaaOtt2Bd90zO78nBPNzvTjqj59AIBAIBAL7ITQ/R4ApJ8Xu9K2LQaxh5x6ydA9c4ou+FzlzMkp/Z4HS+kUz+SazQ1aXUmmANPtBV1s8m1zsm8mtz/ukpfMEK0bJguVWJTStNUXFyGZhDa3T/tLIY9JslGWDigeXXZo4P8zsuFgrqekZZtkh8QgT5ofGoTU/tDPvEbbZ2lg+G0UrX0zzBSnaHxcbSMcEot8XjhXyWiAdG+iMGJhTuoieCfJaoEOgvcD42Jh9kvLKTcvqQYAS+XoRp/Oixwa12IONu6wlcLgEglEd2Ie+ignkNUDzTpt1GUJ6qVMixFHi6MzsseUsfMUwJacH8tC6OvnJ9zaPRQZFz7keZ90j99ZZ12W9tjV6DNBSNGhGHY+nPhc9lmR0Gs7RePFuqA6xwp3I93Md1zFHVN+Huerc4vK+b/wVnvKwX99XiJzTrZr2uj1HGggEAoFAIIAjZn4CgUAgEAjsDx+D6JhxlB8/meIVsJi5TZNm3wiAovP3YXkX7hOmD/lmmsZkyksizHoqqyTStMJnLY7eylQY75AFt9wZTwGpMfFvuurCCE9+XYt/qTt2zx15nbbrNAEi3PTnlqaapmUKX0OSobL7ukmJYae3eDoqkeg4ueOb67hj3PI55mmwVNX1wuntlo+DzrUKYMfncOv739prytNgc10WPDu3eBJA83QYUE+FnQ+zytsLoXVgxDW5fvN01MqlxhARc2PKwDsL8JQAPzf6vl05N/gVTSec8YOjnpOr1BQsiWV5mwQMlSksHodyx15yg4ee2tLlubltKSVGJyJA+bV1zwlUpAHsAfHGznq1DGVTT3P241MuCZ53aQOuZ/rrnlCU6AXdPLC73dlHd3pO9xAnt6aBK7A0IacQPF8zjvLjJxAIBAKBwP7IcF6cR46j/fjZIkkgQ+3+2/uaFybIM0CqrHdfbNnFVyUxZDdftnT5pmKLk9d10LvNydpsKzphK6LVZYUdsqJlbwXSAOcycZPnOtw/W56qiatTBNB2CRRrWwIXNoSDev80cLvRCbVTK8DcyMHiLCvkx2pc3VeeLbJ1t4rFKakvLHuWhC1y9aDYID5PVHfjkqNutzoZKt8DNgDixZaE0FN5PDfM+FDZw+PM8PSE0ADwEAVCXBMVVtzhiQkiRbdmY5jx2S5YvT1wW0m2yv2rm+TMucXLs3OPAeGAZfd4FkCzhl3XlPeDTwPi0lvo5z97l/mi3DeLfVJiJGYauK42+KskqNy0PlYJgAi+by0DNDiX97nM7pLXK3f/SyRDXdq25AJ/vuff332CzEr6Ij7nh6RlMX8DWODs+r8EA+T310omOp0MbZf+a0W6VdNet+dIA4FAIBAIBHCkzE8GWcQiOWhYhFTktT+eAQIaOiCnGxBrVlmvvE2Yno4GaNtwdWfmR1ImrDZ1XcfEFPdrXlomgnYwD98xGj3Xd13Xs0KyH2WZ8vCqPJq8zuyKtlpdQMRi0Xq6SLE4TL3Q+c6js3RZA7QZVBuyuokC2J4ye8Tby65KIERiMkS/Y4eWTGBEW7a9w4wca334PNbX0KfG8Bog/VsCIeb2krVAQNED8f35IDFBF6QFOqWD1i7qzNaUdC/ZrDNa2h1hjRIHBU1mP0Cx+CcSwHirt5X8UZ6L8g6zQwAAIABJREFUE7vPXvgIjV4NTahI4uDefbsHs1GGYtnPfQaTkuVfkmGqJ9WrZnPqzroJhfmx27LmpzAKPvAhb6lToCxpfg5ngHwyWj1yGzDhctgnXQYnr5Ygtmjf6wAw8j0sfzeoD67QIMh2skA7qIf7rfnJ2O+ZOhYE8xMIBAKBQOBW4SiZnykn3M3FArYhzd3XdIcB0hi8FUwW1inZjz6sP1C0DFuy4qaBNT+W3dGB7NjLR7x+TsgqEn2PsoY55YIwQH5JFqKethYWZ14v2ha7Pqk2PnkosypimKvw/czAiNdXaphDHtI/76Bd19gjbKKLeMkxNXRuhlEzWOyhldx+uW3N4niPMK8BaiZO3dq6nPYi36mZnwtXZ5JUGDYlBlDYoIsta3/oXvFM0FgHRhTGh9YfGiglBi0v1D3ISUr53pZ0LwtWsZwDZ917JggonmBc94xZD3Z4aukfrlmEyYlRM/EdMlq3W83m9L26vBZIt3EMCbMH4gVG9632IuRkqC4gYhU4FG02SNdJzFKp/r0OaDcDpAZ+BQxQeYSUnpGWmx3sjWZ3vP5n13q7Q1+wUZvond+ZBagYIOhzdzkG6DLJXu8VrefvWHGUHz+BQCAQCAT2R0a6VdNeR/nxkzEzP20rlb7mOwzQErw3C1u03poFgLOK6eHUBeSVk+Y2D5L1DQCbkS1PauMYH2aC5jLLJFSxgRzLM69YbYywIOIZRkvFgmRiT1ijlE6s9kfrMDhyfEmGyvtj84hZF8XIlCBGtK1TruH1FMIEkdcMM2RKg8UxeTgWEFvXfG4HkziVWSLT7WIyVLhzx222jk0zHntUd33HXUPHBAFF77UWdpCWk/UC07GBxCNsZHbIeoSdETvKXmG6fx8baHQM0F6WdAMlQSp5giXPBNG5UNfYe4J5HV2LGdr1Ere30An1OyMLk+us4NT6ndp1WikxHEPidUHCuqS6TEgIfgds7XOo+63eZdn+SGrMXgd0tQyQLuu0pYV1auPzP18XnwxVelhgRrpMvvqLV/Rs9l4eh7rt4LSgzAD589LiTVZyCg5jn+6/t9ftwlF+/AQCgUAgEDgM1Qf/EeP2HGkgEAgEAoEAjpT5mTDg7nQqn3YPNWttqC4HZHMCOUXdelrdC5+nhhh0FIHcTOiuaXmH1I2tzN0PjBdmfzLFQcEPtevz1k2JTdV0F7taq2kWl73du75zfLyW4LkSCI9uO6BSXzCXzfQ+d0Ztodq4smq6a0dWaz0mmeZiQedYBsdTYXlrlzIFqKYUq2N1Wd63bvpQ15EpMecO71NlzHWzOTSZsmy4xXt3eO8Wf+HE8vM2e4/ddSkx1tkGSgTKdBdni5fpL5qeKtNW9WTEPu7wVRu6gSrnAZWdvoSNoHPJx3jVby+6qX0gRB28VFW2tVrBOE29UrcKkOjc5YfGPE9a2CZ15Hbfb/pr3rMVQRdhNe2P6um3Y98NviV4XhZBt0ICeBH0RNfFH/pl1Cnafbwn3m+5mA8ld49pO7j8I607hcsOnf4a7vO0V842nMqxI5ifQCAQCAQCtwrHyfzkhCem01KgPvHOnL2wcsGrmNVph/xvfxWz+FN/qXsR9EoCy81LzwDp37x8YLSMz4VKhbF1ARAlEFzH9X3+zczPvJ6kDa07lgfQAt7crGtSPJzwfqhfEjyPjmnKKsGfuIZ7t3jP+LQCIzJ4DMkFiVOCYRGInnAIALuESW9hxdCSMNUd81YrU31qEMcWLQWQhAtbwEL0tWGJLMO3YXd4lxpjo+4RFj+v6f45lxQYlgl6cCyC5wed+zsnTGUGiJlM7R7PZdsDRNBsXUvwT7qIJeVGuXY+JYanHiafc0Bvc5ZsPsCyrRggkwyVty27eZvddRgfzx61IL4CvC59qkpMjYgDg3tOsluqFRFBb10bxwDpMexmgFplOwTQCoUBsiEI+DDPDwhPIuVLImlwmhEOoqra5c69LX8/LAOkm/vztC8DdBOu7uHtFQgEAoFA4NZgdnW/PZNBR/nxs8WAT0wPqPWW72ipCwAT2RnFtbev+elBuyayldjT/qycBgiomZ/NwCkLyPoelQ7Cp0KgZJ/Tihkf1oZo5ocHxyyCDBaAYiKUwZHctuR0Qza5Jy1Z+3Ni98Pu5UkxPxInnt3kvS6o5fLOQds6OiBhllSbJL75g11S0MBRB5abbP9bcStmBshqgfRvrwfauvOkr0dhhXi/tNxaJggANiQm8kluexogQAXLZObw5ILq0H3mAiQCSg80OPf4jgYIANZk264o18k+Yfl9AkVhgpgBUgoPDhR6hzVdcikXXN1FH3TvlmzjsUBxx2Y413afFNVgX72Qquo0JcIAaSEMvc2rOHxe66MPxLvB9y6d2k/PDb7sZ0nzs4sBqrfVIQhmLKXB6LEmrVRHhSXq37fy92Cw62Oh2WiQpf/BHWKPAeppe3pJuANXg6P8+AkEAoFAIHAY9nFSOBYc5cfPhIQnleanFbuA0094PQFffM38DBLQqu3NskQV8lf9IJqG2YRi7c86FXvA64A2ZNk+MLJVr9IdcB1JmMrMA62fWKZg3sYeYG32ZnLeX3MZTJ2e91erjiRhPLFslJHsMHtDx5o4U6oELmywOC0dkFpnLZP2EPNsEHuVseYoK8ZkoHPr9UDluFgLVHbd0wP1GCDVrWKN7Pqkxi9D4fQik133ATH1b06TwaxQYXwsaziXMfNDjM947sqJ5VEHLx5iLuhny5vGB9PrwXhOkhV/yh6awgRZ76+tDmrJ98AVvOGYYdD2uU/M6RmgfRiO3E3/svsPUAmIqJ6/zWS2Se98/7Y0OZUOyDFA/JwoXVUvEGL7ivbOx+EMkD/X+nqsXdUu89Mor4IJ0j2z5Bk2Dm22ZjRJfP1+2kt9mGeNOylwPTjKj59AIBAIBAL7IyMEz894bPOAxzYPYmIdjImPw3FE7HJNVuVZg/lZuZhA1f5YX7BgzbIVwbGBVk77o3/fIXao6C7IcleaH9FzeO1Hx/sLADjXq+QQ9SkXGh5JPpln0QXVGoGK+REGhtY58alOCzHaeDsV47MPerGBFGTcNDhheoRxUjPy3J4Yk4EZMedFkxTLUnRA1jNMND7eGwyozm11/o03mfMIq7z66PpvleaHPcJOrfbHs4bGQ2wiFqdKmMoM0PzKMFq1ZD3CmFXle15rGsbOdfUaIA22zIvnpPUI42SsWvsjTBiYPbUsrX7J+xd+7w+A1pgIMcIFrEeRpkueT37V6YUalXp/k6wuhJmRfRkgNb4eA3SJZKjXzgA57Q9QzjtfI3+biXSqcY49GzR2mH7ThvfuNED2ejDbz3XaWKnfRSc03YDm53YJnm/PkQYCgUAgEAjgSJmfKSc8vr1TrD+lafAsDccrYeuRrdZTxch4HVDPK8AwTB1LpjBAFNlWe82QBb0iF4473vtr1DGB2syPrJ80mB9mD1bMfjBrwLoXPo5ax1PFAso1S8RMjzAkvGQZD58fZdBwhFSJRO3Ym+TLUfQ7ycckYTS9wFhbxGPkdd6/Pg46z8JgcR3uyuqG5rr22Mq56DBBUHogz5DJuS3d86H6BKl8Db032FzHejx5z7AL1gKttOZnMEvWA/E9eE7PyYO51qrxs1N0bfOY9PPCzOc+HmGMormzzMbKRYXWbFQPV2XZeh2Q5kV0uaEg6IHIO6NBq3u9UBY7xyRRoPdkgOYyt88dDJBuv4sB0mVqlG4HezBk3XKl8XJLrwFaSoJ6SCwgPma5f2nRSobq2zCrI5qf1jVV8YJuYgJqX8/mY8CNMT8ppRemlP5OSum9KaX3pJS+tVEnpZRenVJ6X0rp3Smlz7uJsQYCgUAgEDge3CTzswHw7Tnnd6SUngXgV1JKb8k5/6qq86UAPp3+fQGA76VlIBAIBAKBK8Jty+11Yx8/OecPAfgQ/f5ESum9AJ4PQH/8vALAG/I8J/H2lNIjKaXnUtsutnnAJ9ZnmE5qITLTekx/s6iYqVum0NdGiMypKYjGz4767IxB9+shKTFU2HQvgi6BEEkAbaYaeCps3tZzfZ+UKFemSniMVTJUrlhPzVQcvXPTBtQ0UIn931nqNs79XcS+MMs0avdc3rkVX1dolft0Go0prCxibj6e0dTlJKlTLo/PKFNW3AlPg03Ul50SnIfixeM8TcV9qHF7cbRc32yqtoJaTpO9F7YyVcb3ihJJSzJdmy5DkqCOdioWKFNga5cKowT0LFO77LbuE6MekhRVph5c+AjdpkyF8Ta6Zo1UGD2B8z6eL3tPfwHYW+RrpkOW3eD1EPkc8jtu1/QXsOAG35n+0nV3TX/NZajK7Cha58KVdaYH3WTVvD+qu4XFktC5nvaaqvJuIER6DMZG+gtJgO1T9sih8xRwOcBRvZduZNorBM/3FymlFwH4XAC/6DY9H8AH1PqjVBYIBAKBQCBwKdy44Dml9DCAHwfwbTnnj/vNjSZNMz+l9EoArwSAB/6lh/HY+qyINpWVytYcW7BstfIX74WE8b+QNiUgobVkR7FA+2I3EVjv8UVdLFnaD1mtq8zJUIu1vBHmx7kxe9f3E+3mTxYTGyLMCHA5W09G/NtmJ5gUMUaTr9sTQBs3YysUHqQt92mZDWDBZpxq+7tCLymqDi7mkqoKE+QZINOt3VaEz1bw3EwtIOfJMmfNgIgugGNyrM6k2nh3eLjrL0xQw+3bp8bwz9JGBRRkxufB0TNAdP+qQJ7ipk4WLjNBE/osag/MvDKrs1L2vgQt5WCm7A7fYHNqV3dmI/u2t7/DdjFAetvl3LyXGSANn1ahxwCZ3vdmgMpKjwEy543vzwYr5EawUNZhu0wxr1jGjaNqXKCGZ4OECeLUFa2AiMmyjmOyx66Fz1Uy1FSnTprXa4wp7SVyv0rMub1uz7TXjTI/KaUV5g+fH8o5v7FR5VEAL1TrLwDwwVZfOefX5ZxfnHN+8Z1HHmhVCQQCgUAgELg55iellAB8H4D35pz/WqfamwC8KqX0I5iFzo/t0vsAsyX7+Lq4utsgh8T8kAW7da68Xr8A1MkcmQliK7YVzM2DLZ86oWOxAtb0s2KAGgERSyoMq13i5Kfl2JUGRCx/u9w4N2rjwS3MjNfXOF2P+e3qVhqg/VQButTYI95NnetIObMjjTa9de02z6ER6ISwZKXajwbHPPMBEeUcDHZp6vAQOudN7bwKNVCd49K/hDCgKht3TwhbZFKgtEMnbE5sQEQOgggAD4yzXc0MjARGpBN3plzQy728/Ay1rO5d0GkKhGEa+Fjn5Z1hU7WbxquzAfsMEHCZVBjl0hzCALX78wxQGVGDAVo6/Y6VTe7x1nv1vNduBqjVU+/8LOmEHANEywt1YD3mZwlVKgxXbnSg8mo7pzqsAeJngMam3ie9IKD3C7fJ1f0mp72+CMDXAvgHKaV3Utl3APhUAMg5vxbAmwG8DMD7ADwJ4BtvYJyBQCAQCASOCDfp7fU27IjjRF5e33Jo31NOePzijrAsWvPjQ/p7DZBngIDdXiytYG6MQ5Kheg2DpMJw2h8AmAbLYJVUHrZ/zTSVbfN6dhoQWVdjkniBubIRaVEuoTecakOqdbl9aLRlBsh2R+zHej4/7HwlzkXGtarDynnPMV13FwPU7nFx3PacuK1OA2SSuVbb2uu6f3GcYw2LJFC1WiAdCHPyrKDzCNvkOpVElTaDU8a49CxAYYEkKanc65btPEQDNAqrUPbDz6IPhMin3KS3cOkx/PIyaCZDLaNz2y6hAVpIhVHvsc0AzWPoBELkdXlZLAxJWMm+Pq/LALUVnJ3xL50fVyfZPYsXmNKfnbseKq2U9vZyeqDiKVZ7hgloV+IdXJ1E/ruhihoeYPcLGft5OB4LblzwHAgEAoFA4OZxm1zdj/LjZ8oJT16sSqybhjeLMD47rFZgtxfLUjwTtgD9XPFWNEBLXhucCoMYIKX5YebnTuZEpjaOyZL1WqbsLfPDLIJWRXCdyTMz2VtjQH9Ofg+dQtrBAClLaGiUAUDabs3eUovtcZ5ci+A6vEPWOGzsfuahLFtMS6+U4v1W8USqjhuDu29annTeU6cY8YNZ36gD8Szg1qXImBwjBBQ90Jr1QPS8eC0QUO5HZpKE1XTJSnXclBIvpc0GcRttfbOEi61t6T/b+D+tfnpYYoJ6nmGLyVB3MkC6zKP3bLXaWB2PPo5uLCAelFTVNwltcto3jkxjEiPzT0u8gY86L/6xvQIGSBLDDqYUUDogrsmkWspmOY/bMz32Xix6Hn3fWsZnHOy69Knv68T9Ba4bR/nxEwgEAoFA4ADk2+XqfpQfP9M04Knz08pzRf/uebGcj/MpeeikzAiLDmhse7F4DxZdVnQ87ZuqFwEaKFaE6BeUteq9V4TR8rqFBdZrnxt9I2wBR8L2zEzf48J7qpT4HMr+EnOLC8hCkySQjd24poKBy3l/Ks4rR0gWNsgmJtwLVYwg1XjTGeawuOpK2d2u4Qsju+5Yv47lmX8n21auu/UCM4a6aIisHsh7COr4UYUNsjq6lrelZ1Y5eS8/L/xMactadDteu+J0Fkv3c/UMaaccZlFd1GlJCMvH0YoO3XjOetg/GnTjmeoyi63nr7WtQGtPyi3NsY1oFPyXge9r85ywtofvr2yWhiTyj/dG700zTvfCAEGSxfo2Odk2mt3jLRIL6IDkp94bcSna/+gZnkWCsdwdueVVGrgyHOXHTyAQCAQCgf2REa7ugUAgEAgEbhli2usZjpyB9cWJovXLBe0ldbyYOEGoDecPAA+MHFCQpr+GC7PuKXtdxiLoXaLNJZSEhYUGLWkCiDJ3rsNC2Teo+stgI7MrvekvKGq+Q9UnXw8yJ5NlqmoulmGnoWrDdYeBl1s7SDo3RoTMImV2j2eum31ts5n7wT7Qwk6h12laTfbt6XDl10qRE4qeOjnhrjl9lxF/ttsmNx3VmvbauGcn+ymgxnRqSebLfbSeJRsIUVJTJDu1ZcT9cl68YPQenqXGtNp2wR3+KnEtyVCbdX05rTVubxbf8jHLDHtDeVsCbe5+TngqjGfTqukvefz0UffmhZbOhSur3jX1O4lTrPA7zUvgU8vV3bm4V9v11eT7tTMlxqkw9JRaEU7frimvlNKXAPhuzHfc38g5f5fb/lIAfxvAr1PRG3POf/Fe9nmUHz+BQCAQCAT2x03F+UkpjQD+OoA/ijml1S+nlN6Uc/5VV/WtOeeXX9V+j/PjJydsz0exNkwQN5cIki3ZIny2Qui5bLRL58oromNlifTC+A8uMd4+IdXFvVL3T0tvre5y113ez+6xsHXELMWknxXP0nhWhzYbMkoYH2eZMdPEbdRhidsql/FumQmS/SkJKe00bSa7jQMYarf44u+NveHbsBiaGCaxIpUVmNx5Kud/KeDjDgZoIehkFTyxlRKjEkPPddadgJhA3x2+JawvIRkGs37HJR7Vruh8j6/kOWMre657GStZW+5DtuyTjNVFvVuKgXLIHw1JN8JtubyMqPTrGYzkr3sL7jr78gXsFXK0Smthhc+tMBL8ahAGqErhUvY09PbeFf23ytoCaPtQ0N8AdrunQW7YcaWRCmNnMlTF240dplLK5Q9UNSQa1a2ZgnoJgPflnH8NACid1SsA+I+fK8XtiWgUCAQCgUCgi4nc3a/y3x54PoAPqPVHqczjC1NK70op/UxK6TPv9ViPk/mZAJwPJVGnoicyJ3P0zM/WJmzUaSE4pL/XMPDSa4KA2g2eLVsfsFDrFnbpgYyF23Hd1bqjq0CxeLhgXucw8Zr5mYRIYDaHGZq+u25hdpKpy30I46M1P9yfMD9WA8RtB00xsR6IGZ8NMz6kBRqUHeAyvKZOQMSsmYGOK3LRPFgmaB42hfrnocmxMpWlrzft0zEB9bulZQ1Xo3J1FTMj3bDrMw1bGJ8TWqrr4XRBnvHRARGrbcwAOVdkzWAKi+o0dhUTtOBu7KGfpanTbqJzLMtWUlRhsi7/3PUZIKDrBu8TDDexnwaohYoB0l0wK87P0D4pMDjhL1fyuSQ2Nf3RZYAEC5qfA5Khehabb+2NqnrupERlpF4L1NAJuftzlHJaVydXp8K436qfjGuL8/PslNLfU+uvyzm/Tq23duoP/x0APi3n/HhK6WUAfhLAp9/LoI7z4ycQCAQCgcDTAR/OOb94YfujAF6o1l8A4IO6Qs754+r3m1NK35NSenbO+cOXHdRxfvzkhOHuIIxPVtb2xiVzrDRATr8wt/EBEW1gxKIFKtYfB3Nj7cLW6Qq2ZPmsoLxaJEx/2xI15UyZsE7BmWrD0PcuY8tjrDwYautFZuTZ8uG0HdT/RhllW2Zt2COMx8uNed5fWX1etyOs0YDmEgA48wizNWw5cRA6ZpG0uIj1QLIcnXm3VccsLA3rg/bQAJXY+LZcdBC8f+Uh5pmljdUVJKMPonujO1HdYoL6jJtu04pKV8gumxZAGBu9l+yXqbnUvyWdjCTgZa0PJ+pVz5/zaJQlndPCopY2h+iA/PMmHmfDtlX92uAZIGDBE6yVXqZrtfcYoNKm1ofRXlu3Pgc+pGu0yCsx4+MqZWGCGmBGyR9PNw2M7sl7fdl102Wlo7LnWA+5x+sNTgtkmZ+2rtO/k403GP0ch3wjmp8b0hn9MoBPTyn9LgC/AeCrAXyNrpBSeg6A38w555TSSzBfrI/cy06P8+MnEAgEAoHA0x45501K6VUAfhazq/v355zfk1L6Jtr+WgBfBeCbU0obAE8B+Op8jyGwj/PjZwKG81S8vU61N4u1GtfMDjkmqJUSQ5Yn1lOllUS0p2FgTY54oGmpCSdhdIfTYoKkbAcD1DKSxNOsSrBnl0CxSk4Gjl9jLRxDKjgdEDMz02C9s0zuTafbqZgepwUCtD6I6pIFOqytFmhSk/YDaUuGNR0HWZeD8/oCoPRBc5vkPbgYLS1QocjsktGIJ5S2rD2w60YWwcfPJmiyS68FMliyfk1ngLecaw1QtuW6yQ7tjylDu840bqo2d6iMnyFOQ8FeYOINZmIDLSeg7KWbARrPA52/k0H3ZbmAfbQS++opNKvWjQVUeT7pWnY//dQorW1VlQrlFeO0Od6BC0rj4yD7ZcZVP0ss9JPH2ad9aTFDIphzB+DWmyxou86U1N8NIYrnOuf+PcijaHjMCtPuvL1ku36aiNUe8nT/WZh8c0EOc85vBvBmV/Za9fs1AF5zlfsMb69AIBAIBAK3CkfJ/KQMjOepGBMmzg9734CWpEGYnAU61daqt2hZryBLzRZRnQfGtoV7p+E54uOMeHZniQHiPUuEUNb8qMipo3iIzfs5max13NL8nAxbU8bMzzhMZgkAwzDfTmsq246s3xjpMAaznH/TuR2pX9YROAZIe25xXS5jJzvRAgnLo/Q14kVGdU+o7Zo0D1sVZ8SzQWKd0hg4VsmwYB7zftnSbHiteVRxUtTlZjZINEusu7AklRFXZBczqXi3LFl3bcvZ33lZx6JxWh8PfZZ6rFDNADVYVF+HtT/JrgNFS8cxfPjeF0J0D8+wpWSVJa4WeXNek5W+KxbQpaJBm+u0qNgppVqWwowMv0uzvceTZjf5+fbaHzdjkfR+WWtHHSZJhmrpbH3sY0V123NR7nm9X3usPsq83jyJRxgnt6YqsmxpfvqM+ry91gTpKOb3m4XJuDnm5yZwlB8/gUAgEAgEDsNt+viJaa9AIBAIBAK3CkfJ/KQJGM+BxO7Fp9qNmae1aArGTYOxq3s+1VNldtprVzh/W9YO59+E28TC5zHVLre9gIijE3oOatprled+OBijD8legseV/a2GeRQs9jylbadU/pSa9uKptvNxLrtw7vDiCq+mi/i8ZDfdJZpHmRZTB+nK+tNg5XqMVCfLdBfR8CyW3pQx8Vh4KoxTXySe1iwR/8qYdjkeCKV+gGWlp7142oCXrG8XcTnvRom8eylD/BCaQ9pv+muuYd3hJQVKw9Wd4ePi1dNfalqtmvZyDgfkRKADI3oRNPfhpx5aGCtnAp4yaxyHn9a8JuyVCkN+7Tn9pbYlPw3lROxax+sFz7znwbWp9wITUsJA7b/si6+DyInn/yc/xVWOvTf9dULnYrNPMtSF8zVJYucZZfrLygKA5RAiel2HZRBRNLI6v/cH1xjk8GmJYH4CgUAgEAjcKhwl84MMDBconotKvJxWLLxjBgi0dK682qh3cc89E9SycP0XdEvI2YX7JK3c2qEZnrYlNZAQUzNEbA0zG8T9esbnRLl0nyRmeja0zq7v1L9hfmjbeEJjnMd2QUzQmgTRnBgWKKxNFgbIMkHDCVlH50m1oWOku5cZHmF3uK1ii7womgXPA5luWTE/g4igmWXxDBDV1VbsZZKh9lJULKUuYOaHKR9aF7JON3EJZYWk6DBCS2M8iAFyp8cKnmm5aj8fzWcJ9tkpzE8y20+Vk7hPTcH3NjM/msXZhw2q6rlnR1W6VlwqFcYBDFDtr95vWwIg8nVvxNfgtBXSDbuG8zu1/9yk6iAtE2RPtT32HgPUev22BM7zuipwrGmWMAvz+kWT+YEp82kuSroL5TSiyrY3wML0nBaOEcf58RMIBAKBQOAg3KJM8sf58TO7uhfLwcSmExaH1235JIkoSxNxh3cuvUvB3DwWtT4OwurQIFmjM6Zi2bK+wTM8lXuuGpJPscH9TtTvipJWamuWf59PvG3u/ynS/JwqfZDXA62ICXpqPa+LFmgst92Gfm+ZtZFUFW0tEFAYn9FpfApbxNtrHc/gND5MQmUVeYDZJz4NEhCRzDzR3Si2qEqC2kmG2oRLyNp0h2eyRu5Ba6HLWNV5kuwiTME4JqhigIAqeGKJ5eeYANWm5wbvHzHABvAzbRb0QR6950xbyt4Nnp+XiqkBrB83Gi7JhUJWA6b++P73J2Hhcb8KXcVBqTAWGKCi7eH1Tl2jyXHvUBlAzcj4usV9fXeoAc8KJYleWj9bJdKD14wtRHvdwYy1U2EJhRUkAAAgAElEQVTwOru+z6vbNL90zltNKu0Pv6sb7vGiw7z/mp/bhqP8+AkEAoFAILA/cr5dru7H+fEzAeN5Lt5exmvGL50177zBgGJlFy0DWaluty1b399MPrjaYiI89nwiGmqNImK540OkC1vky0v/klSVxnTK89ZsxRL9caYYpjvkSXM+rWh9Q0taH09LXUpDcEplq3Frlk+urecYUNggrwfqaYHm38nUYQaI8suCiaWsJvhFH0QaItH6sAZIp8LgOhuuS2NwTFBSmqLiheWsVbGs92CCvAbImNCuqpdKCDmhLPTOPc4MFlutOlgje4v57BxFL8QF+njaOiDPAOk6lUfYqv/SrVJi+MChtN4KNFhSFBCF1TCmK0+wHUwQUHuEXYYB6un19gFfBp1kY3cgRM90oE4e6uvKfazrtZkSYV10IQdonTglhdXlCRGnn48uKUR98MFr9kiau5Qb++T78UlQU18fVAuDPANUXgrMAtXMj1u2vL3SJO/lwPXgOD9+AoFAIBAIHIQQPD/DUTQ/zOaobQ02SJfDaSkANcUsXg3zqjj7iF6hP6auV9YC88NWwJqsCc0ayZy20zD4JHqtMP7eyOa+LlgDpOzJO6QDOk/M+Mzbzpj52RaW6HQ4pSUxQIP1FGPr+O5YUrc+ScyP6IFOKEXGCY3lwjJBADBRYtl84dgbx/hob6/C/JBV5+L+DMqE5swjzA55Jiht5o6T8vbyzE9hXdhyblz/HSyBjgnkQ+9331G5/l2xncxgCfNTmhSGhwvssuTj7TMBS6kwKh1QLyVGM2ZWm/GRempMEhNosM8Hp6OY1LllRrTLxCykufAMEPfB6TUMLmHIH8IO8S28PwOEZc9CoOH1VerKPe30QjrrS5HteDqTBZnsJaX2zx3w/SovXDskf/Vb41ZJaxr1/AVZ1gDZ7h37n/w5Lp5g52l+3/lYQD5xrv49pHwDzE/E+QkEAoFAIBA4Whwl8xMIBAKBQOAwxLTXMx0ZGC96gmcrbOZtW9btZkfpAhKQqyTbtkLOQnDXp7N3K7VCnnsXSO+SzukpAGCb7BSA35HODszoiaJZFHpGxPmFyiVxlqksXQAA7uT5RPlpMKBMgd2h6a8ijr4z90WC6MfHO9KGxdBP0VTYU7R+TusXLkM8gJIl3rvD85Knvy7KMVbbeGqMdbBqWo3FzxIQkafBXODFtE2qzfK0F3xWa6CtkNfQsxMya7D/y8m7xcs9P7h1NUPD9yNPC/GUmBc+m/2ISJrv6fb0F9Bygy9h6Uy9S2iBTTiJzptNpsYanDdPVbFjgZ/SWszyzn0k28dVTX8dgv2nv9CY1kpmtdyv5cKnhjTA1G3d45ObEuO2HFbCxBZZvvgimt7WZTyGKvWG1NQnnx9Sd0FaGeArV3e/5GlpNcVLg9mQbOEuT3/xSBZTYkzY3qKYOzeB4/z4CQQCgUAgsDcywtX9GY8kzA+tNwXPznrxlo6+CTqGCFuaNQME8YP2lrL/2m+J3TiQ4IkTYprgg2T2bL27euKgh9ymCJJHxwZVDBBZRWepsDnMAp3xNmKC2IphJmhuN+/rwTw7eT5JjA8zQE8RvcYu8QDwxMhl7B4/L58SITQJrk/KrbomlmbLwmdigractPSCyhVLIq7z7OLOiVOZCTpRdde2zcRMEDNAvFQXnMmywqZYdij5TJ66jNcXDN5cWZi83m/j+2WDOjmD1xAaXqhN4/cCaCMb7WxjUfSkxyiP234MkIZPgeHLW6ljeikxTB2XCkOHYpjHbBkgoGaBeL0KJKqOp2KBbpgB0tsqse/SH0F5Zzo3+AbzI4FmYRmeKvihunZ1oEVe9c+LVlazAJmYno3db80A6TXnBi9DUedAHk7veMDr1k0eKO9mrsvRNO6KU0H2ezHu7+Hqfr04yo+fQCAQCAQCByBfbrr5mYqj/PhJU6Ygh7yutnlmx89ft+atS+vmNtbMZFM2L9n68kyPJAEdyiXgpKFPueShYkUq3+1KB7QHA+DZIGGAFrQMZ2R98PzzmqgSZofWSlhzl9zfn5hI40P7OSNhzZPM8myK5qfnFs9MELc5GcsY77o0GRtmaIgBmsSNXVm4F07Hw67vrtxsc27wnPzUM0GmX74cvOwEGpzrWH1Z9eKpvYur9dzTIuyDStcBSevi3eIlDhwvtVEqrvNW+1PaLLgM78EA7XohL23eJ1fRlqmrXfHwGuktGMKqdhggQLnXHxAI8SrQZ4DKzstt6TWPvK6aVDqghboEvpqiy6I6HPzQME10PaR3vr/2SB7MpztTv4M8W7Xeqe8GX1+Q7HRAubqna7YquxQYhQGa6zzlXOCBm3Z13+95ORYErxYIBAKBQOBW4SiZH2RgWE8yd6xDs3cT+fl0BI0v4Gqbq6IT0XG3k1P8CwNE9ILuYnSMz1NbShDK5SpDq6SqoB62ne9YrU1gjc9BDBAfBzNAwvjM4zeeYZhZIGaFmAm6Q0zQnWkuv2PSZ8y/HxjnbU+QJ9gpeX2duhQZAHBK+h/2EGM90JqZoJM6MKJPk1F0PbRUT4IPnjgxE8RBD2n4WSU2rXRBHQZI64QkEr73MHTanCWIceiZoNJtLd+Qe5+t4sbzwcyPuM8wc0XlWr6SbFnx9uL7qlSdOs9OnwECmDdY43B09UGqvEqL0SECdBsO2MmB6rbOI6zpGbZnKgyd7sAHOWx5ifptqSMe8wwQoHiKihnZQwNU6SRR191Rp7ybh0Ybfmg8O984Pg44yp6MEoSQtEYNzU/N8+xmgDZ0XkanAfJeX3QApmwSJojuZ9YAmfQy5fpOU+OcXyMy6uflmBHMTyAQCAQCgVuFo2R+Zm+vCRPHYVFJSksivVJ3Xu+Yx3Mtu6XlDeC76Cj9IcwPa3/Kfu5ubOJPXp5P86CZOQGK/oe1P7zksP6tqVuJ80NL7xlWGKCGZUV12frdEnuzVszPBVnowvzklVk/S6dmHVB6oMnHBpqXHxdNUGGLHicvshWdH06RUWuByu29lUSppNsRTzHQUjEmzjMseZaImZ91aZNG25bvPc8AmdgkXl/jRBlNA96XeQal4VlVxSjxXWjNT8WM0qobo2awJEEq9+8ZIK1pcG2qZ6iZDNXqgPju2cdKreoc8MZbTIpKQ/IM0BJ6yVBZl8LLsaKECpb0efuMATgsGWpqaFlq7WOPPVfkR+VVywwcpQUxcX6ItXHskLzakt2fQcmbQiNraIoIfDnrbvjdrfU783Ksninn5WWev7Yej5kffgYuNFkkzM9NuJ3frvQWR/nxEwgEAoFA4DCEt9czHTkjrScMuY4eOjkFfR2VtPXl25nfdU2y2sBkE8eayQNpZpixGWqWhT3AVtv5srD3F7Me51NJCMr6H9bPTLyE9c5qYay0P+31JqhbSRypxB+8z7scG0hiAnEyVNIEKQaLf4suSBKnrs36x4cHpA17hnGk6JXogihKNDNmykNsTRorjg3UixKtf0/swcWsDscIEpZHnRbRASW7zvF+ZF2xRT6mjo9I3rCgD3o79ZKh7mPcVQyQ1f5oA9GzQXLWvaU+d0SbPKNgn612MlTLAClnO1MPUNHYTxoRltGOCVRBgsXwc6ez38qgACgdzx7RoD0DVIm7zCvqMkqnNlpHXHmCyUlse4HNlajE3YrN+DyVV62rKxe39D8K8yNiK2qzW/vjEwgL0zrY57DVfrgwq44s6ilE7DG3NXft55C9wdQrQSVBvV0BB28Cx/nxEwgEAoFA4CCE4DkQCAQCgUDgSHGczE8Ghs2ETHRpXtVulIV73EVnuo6BajoBMrWlaFR2l3bC54nqbGguZVRTMxebueycpmiYSr9LLu93lMqUxc9c94zmYC44CCFNOZlAWZ2P+t7011zWbsOc91YdMh/Jin6taUrsjBTCdzMLn9W0V1q7JU1lJRvgUaccEFE0pcn4hAigqY0TQgNKDH3RDozI05JASZMx+MSpJyyEBi1bgmdQW1rnaTBeVydM3OAlBQbcUtWtXISxN/x0V4+Ob8GnxmhOeThNqeiEJfih7tCVMeef7LSBbrIrEKKf/rpyVNNfgBdBs1jZp8KwoSb8PBE/b53tKFPh0oc4LdSvbp0a4VD0AiFOLfvYT125cn1DlYSj9r6tXd+VUwp7rXNoj5PBlHMAQxs80x2z7NdOlWXjzO9u1InP3xJ8ElTufv4xGrW3Fzz7dRsEEQC29KI4B+6/q3u+XczPcX78BAKBQCAQOAi3SWd0ox8/KaXvB/ByAL+Vc/6sxvaXAvjbAH6dit6Yc/6LOzvOGWm9Bcb6G74YUjsCWi24LFbCUff1DxQWiAmFyQmfmQHSKSvY/X3FST1J+MwC36eGYnkK+zFtTD9rYYBsQlJAsUA7GaBSxq1PeykKtNVCVtYdsk7XmRmgudKK1lfK35uZHWF+SOC82m5MOQufdR1Ok1GE4ZYl0lYzi6JFDL2mwIjMYKl7ZbogV1Q+Xcz8SFJUFj4rB+GOGJoN9MEJooHC/AwsinYBBJO2/Byx4IWdexn7/j494D3XS44690dlwua4chMQkVg1GrAPqlfKk9uiGCBmUW+KAQJqEXQvFYZhXjnQXoehaYikJ/fcDQuC9/EQOrADfsqYoBQGaDHVghUvW1d3L3h2QufG/VSlwOBgnCOzR/5FjBIIsa1vL0yQTkbMDFLnD/4iA5Tsu9QnHqaOTVlVxzFAQGGBtglYpGQD94ybZn5eD+A1AN6wUOetOeeX35/hBAKBQCBwOxGu7vcJOedfSCm96Kr7TRlIm6mE6Ne6FLYmpKTNADWjufvAVsLmoG7jkiVKHWYRiGnYjoX5YR0K63iYwWAty+lYLM+nRAdE7uSTXW+mn2DXygNmtU/FUp8xsoXeYII4ncGWlqdU5YKeqBWYzSnm111igZ4kk+00s9aHghturfZH//Z6oJUwP7TUbYY5iCJbxxxckkPLr9fqOjBbx8EOB6f9cSky5rJ5WQIikpW6cdt1AlVmSFgXJMyPZYLMNhFlOO1Px/24hUqqoS/4DmNTUmJoZsaNocUOCSRRKo2fKQanAdIpGnhfwgC51AWtZKj3iwUaJAcKl7tqqUdFFHjXd92GGTBmhSUgYiPNhWd+PLPkkyvvAxcTlsZkmZlmWgu4bS5halcDpNow88ku7/KHmV/io74HLbvSTYKq1nspMFrovzF3M2JVgEQeiuxfl6kAi7foQ+Qm8Ezw9vrClNK7Uko/8/+z93axtmzZedCYVWvtve+93bbbaYiNbQUL9VMeQKhlg3gxAiOnFakRCsjJA06E1HKUFm+gSCAi8YAs3oJspWkZi/iBGIRI0g8NBvISXiJ1ZAkUO6C0TMCdNjbuGLvvz9l7rarJQ41vzDG+OavWWueec/a5a81P2qpVVbOqZlXNql3jG98YI6X0x5+7Mx0dHR0dHdeInNMr/3tb8dxur1P4NRH5Yznn91NKXxCRvyEin2s1TCl9SUS+JCLysP8ekSlLUhbBf0CvK33gwyVL1P02WQi+2AcwA5EB8r9h6YMpMdbAtECld4gAO0Dzc9SorzFGfYkUfQsSHz7OMRHiAxUgFXE6ICtSev7ABOOzV8tkaH03V1qMBSX6K2s/nLWqNjqYGOiBkO6+6JA883MMyyomyOZLP7hkyEhTsG3Lb2WDVLAz4T4/8b3zGq9lWu6zWtkgI1Aqw48RaHswRo6YBytS2pZEiBKmm1FY51qOXlq0EhlWoXUcZnxYGydSHjzonXTlTM/d4LRRYBSbugppMUAib0wHBNCLBePq4NgoA2l/gFZixK1IsHDc5ipigrZKY5xggzx/hX9oq3trFMqtygitaYDcb9PiIPnnDm1j9JdbZO/8mt7Uxi6C0l75iCo74x81ht5IC3D8wO6s6YFWNKNLn5aZKQ0xDOwNIMvb/bHyqvFWMz855z/MOb+vv78uIvuU0mdX2n415/z5nPPn73bvvdF+dnR0dHR0dHxy8FYzPymlHxCR38k555TSj8nysfadkxvmvHzN21dssVvw3b/KANnXuuOESNNTUgS1GaDlt06RCwZtVBMyk35EpBTfhP4E1uOLMUZ9iYi8UKYCEWC19kcZICdMYR0QM0CTi6thsMZnNC1Q6/tZo6VkDm3AAO2dzXiXIuODiLBRLTSc897lpYeFbPodHMeWx8Kty7ql7c7ypET9g48MM22E3rMnjf5CaZL5ACbCnbtF8ek+0F0aI57RQBvbDbYlRmhZFnVBpmez1P9xeQDfTsy3jLwVw+8sg5Cs/K3IMCF2KFG+n9ClqX3w+rRc1MybjgQDKO/P3r17ELnV4IKWbbbKyhDsGXCle0rZmjaz9LHgBkDRAenzQKxOjEQ8MaCaxVApqmuOTA8e1UDiGJOk7KxFcm1cg9Jx3UnU6bWAEhj2f8Peiy0mnKO82tFfsSSGe7c8g+bnlmRGzx3q/tdE5CdE5LMppW+JyF8Skb2ISM75KyLyp0Tkz6eUjiLykYj8dM63pEfv6Ojo6OjoeNV47mivP31i/c/LEgp/Oaa5fMmHr/KoAzIGyLQGeuzgu4XWg6K7VqbNZZYbRpdTVNDyW5kYzWh60IzPTxoR9jiW27VX/Q+YkftxiWa6J+3Pi1R0QqYDUovtSafGuqg9Mzes71m/OccNndBAVt6g+we7YpErQWkFOzJaq6MyWWbReaMY1jVbusYE5bAvv25g5kciIyRSbFXokMAEoRjtUbVAIQ2P3W8dc2McK0MjC/hgTI+uM82PrveaH9Mn6AJjfihSzJvDaxFg55gPlV4B8+v3v2J8GgE3FbnBuodGcFRhaUkzoWidzlok2OtigGzcYlw1bDSL1CJx1BoTtKzDflULp9sgF04Y40M85uqYbzzfl0SA4d8GihvbtaZIrmUlLVvRqE2NbVgPlFby/yybDGHbKhdQy2bmZXP8sckAHdvPQWCD9brbMooWbubbcm0uIAJfDbJ0zU9HR0dHR0dHx7Xirdb8dHR0dHR0dLwh3JCo5Do/frIslCbRmCLOBQZKEu4CUJwQ8joKGcUoBwhdVxMXum1smYQpkuCZG8S7vXTZrO6vo/pFDrtl4xdH5/ZC+DuEz+YGW9xf+1ndYa4sxAczFQ1FYkH1rxyoLIWIS8Wv3fSlHdcwpkgogtYHtR0odhLfDnavsCKGwovUYfAvgyohnHd7kXsgVVPtmVMkWzg8EiJS+nsTQgfXqF4PJDWEy8wSAZa27PZK0VtYxrHzxa0XJaXr1rqM7Fpid9cFAmi//6qgJYfFN0T3ufIOnCeAXpa9HvcXu4lW3UZeE78iV2Rh8n4jMeKaG0xESpi9psp4MZ/ztGofLnB78WlwIsSYRLUd2m7aX3Pb+ncCtcF04vXlOFUCRLvd7Io7w/011W3s/wO1LYEyLZfssrYuOwLBcxRAh+1Teh7Bc3d7dXR0dHR0dHRcJ66T+QFaYtmisFtmybzICB0+1kLCEq6suzfmB/OexdFtKgYoskQoiClSGJ+sDBNKXxx0/smFYyPE3UpfTGBxYgLAe6faw7IXWbdVha2VlNDpUxAka9FQ7SZYIYiWhwsMBWOEHLO0B3tCDJBUDFCxhi0kf3g8/+AncEkCuGKcleVJz2OCuJsEwiaEdhcMPxECnSCwpjB5ERf+jumqANr1syqJgSmJTy/BBfe7YnmkZqGYASoEkD/3tY6us1G1xvs1lcJA4doVJsizOhDm4/FCuDqE3MaQuhMx8f7K+PTWK9hRFDfm47QEz0UM3d5/uoARMmLGbTJT7Z+0ImaWhuB5mmObUtQ3xX2LY46QAJHmSxJE3znqcMUAhWqry/a4S3S9NlmERK04eMCXcnFh8K+gTu3FuKVY6s78dHR0dHR0dNwUrpv5acFCIcH06Bc8ksgNbIKWUGRoM8SKVcYQ5dkXrUSbQ1w3MFvki/NZ+Ltaqcr4oOzFkwt1h5X30biwOFa+wRL8LdP7wTE/c2SFRmN12iUlRIq/+oXFL+t5wM5zlkJhcZBssP1tHZZbW5hDsKyWCSz12R3ogeKhZ71wE9FQXnswX+DL9hZlCy0LujCI1H+c6xFaIN9vbWqJEaGJ0vVO52RtSBfETJBnOTMxP1vJB8tGjWWnwGHxWyAdBzM/dni3LyuBsda5Si9UHc7Nv5lSGNsJC2PKCStaeoLlEdnWubEOCM/zoKwFjjNKvf9al3I5WhUk7Ehc/XltHEhhdrgIqhExuWaQ06RMNMbRLg72Evrunu0ZDw/rkjauxXHZZm2oJ8dY4kicLsJSp3h9j7V1OqA3HOqe5bY0P7f38dPR0dHR0dERkeVMC+Y6cJ0fP0lULb/xBW9mStRbwGwJidmgD1IdULG+dR5FGp3pCEYH2yAplrnjD60IMT0e2ipNxGUvROrSF2B+oPn5QKO+di5THpd9gAWIeUsa6LRFJY/7UdusiTTErDdmgMrx9XxybdJAQzQRqwPrycetTHqgifr/oBb1pMdBYVKRwgpNdhywCdAPlBOZdq/+BVAM3qFaxon+jBHyuhfLkwZ2AGNEl5MR67epIsKiHKJdDJUeHSYGmu/ItcuWG234OCx/8AwWmjQ0EiJirG37QPWaZf7jlcIokYDLNngXFPamfvcURkYZHyQ9tASiMYHhsqwd+bWpUcO2ljBUp8ayNpgfXIMdjrt+HgxucXC/TZJmjVSTo3NV8kPxDM+JSLGgE4oJELHNsBL9tWw/xANWUZANsZqFeko4D2Bw98sOyYlukYyVn39xutGUn0Xzc0u4zo+fjo6Ojo6OjotwS4Ln6/34GS70mZK/N3mlP5geKP7hQ2ftjzPSsq1bmSedh4jIsFL6Ih/AABXm5wnW3ABrLlqPXMhTxOX3MeYnbtNEVQH2qLO1iY5lXNAUep5548maV27WqGbRIeiQsG5ZhrIcB90/NEwPqEIohelBcdfZ2KHlfCZfYNbS9WMao0uyzbtIvRN0Mdgcf5YW2bESBZJ8biOME9IDlZxA9QGM6SFtEQfSbZWfaJWoOImNS5GiEV8xPfZctIZDpe3JK8vduhMM0LLs8lxAZ0c5ujdsnVuKmB+LpHT3Hc/zyu69fsfKyOjFHSmSyzRALeaVmJ6toqin2CCvG0HQrAVfQb/DzF/j+Um0Du/dihmS8hzMYPLtWeV9OaadI8C42Opcn6eRctCKrjBAIu7VCZ1klSsrPvehC8+U5+eWcL0fPx0dHR0dHR3n44Y+uPrHT0dHR0dHx80j9Wivq8AwvBxnD1dBiMaGiyFO4X9h4bP/DbdWQpp13YZdXCJODE2lL1AqYXZiPZRWQPg7RI0W8k6h7yJ1deeLQC6TSU9wcvt/SBAiqzuKaPgtKh1g99d0xjacCM6E2+6+m4ssIdx/uX4PWv7DEsNJSQ+AFANHpB5gkbTjq9ntxfMHqTGb22tBSXpGU3EiaIwnoabmBnOuAFxKLoXBLqdzBMmyMu/3Y+fBDeptqv1Qn3Kjbe0Kawmd19ZFYe3HdX+d+hdh7iPnPjcxtE6f4KqG25ZC30WKKwwuZRZAt55luLXhvkVaChwnpLKgMPtNF/gJbLl+j1UbuJR11p/GimuM3V3+1DkMvpTP0HeQ1pfxJUaKu2tA53RFDCbw+8H0ZdxfjBJ97wXuOOCK67fjleF6P346Ojo6Ojo6zkd3e33CkdLCmpj47SXuaGsbEsLZ1JIguqYTTVnwjLB2nxgR4e+wxnQeRTEhfBYphTSfnjTUdqUIZ0xlD+uOhJdnMTIqulYr5T150uXlpA9qWT6AZdGTh+C5WKCnqVUwPkxeiNR1ByFm5gRwg2ORwPjc6QU/kPjbF5O00iDGCiFEfzkOCs4encU274jFye35yCKMeo5LPzMlhQxns8amgOkpscOlDYmhOelhVUHE76cSpMZpqE37CphyZnzCLueNdWGJ7xR1bmOI1+TWyxdDNV5jg12tSmFQeLxnJ0wEbawmnqUaxuykOcwDxgT5dwLtf5yj6LoZsn/Bf0h7DnJ7fuYQdalD2VPF+CjDHELdaVuLMNCgC4tm92MkegYs9H0r2aH9D4gHtpZukBQdvt6zitmN63lZD3V/vbjOj5+Ojo6Ojo6O85Gla34+6chJJI+jJP08Dx/w7Ec9o7oZh0Ja+CS0PiNbA1KXH0D5DCRDo6SHIsWPbHogK4YZ55eugIlRJmMg67HF/KxYbBbKrWaHL6w4aaegiXnKOJ5qZlJRszwou1JYFU1CaEnc6n640w8gmUpge8qyGHI+ncEo2XEpDHjvy4DoOT5qkkScxzFFPdWdo/qsDwP0WZMuX9YXi9dbqzo+x2hxtvUoBE74xwyQOG2PFfFli7Oxc2aDVhifJpn6cSxV7L+1K3T7bAZoqzPrY+RcBkikGPiXFP4sOqDcXt7Q5I1UogIaHU4rIeL0ItgvveyM3XEvvaIDQoLT2BeD+0/B65gdOucfKELgsWUoKbOS3LBof2i9FKYH2h9jfqzQqVTbVAkQse1YD/KSTLT9ACA1SliKPqzp9Boh8NlKX+SbckE9B67y46ejo6Ojo6PjQtzQB9eVfvwkkd1QNA7e0rJqc1HZX63n5X7vVCzPCkf6aC9oenTdwNofYnn8MkuDjqgPmx+qtrPm8TuSNWH9aFimXNJh2tO8s8IOmlixYn40WeB7w2NpqywQpvuKAYoJGUW2CzWKFDbH+/fRf9Pi2LxqaJCcMJQwOIPiI6xZ5DtjgMoN3+kxLekktAYUhecL2bIOYoYwQWKCNo/aqNaxAa2Xj/biH7B0NxIjlrY44ImpOP3PWvDVx3yhVokRTzJAbumq3q/Wi5R5sHQ57KJVDLUVxSfSfu7WWNkyvnZh+XIUXafjaC8oM5N0fh1rxUpHp6AbSCc0lJNF47De9+kSVMwnNHDEpou4oCvS4BQ2B+sdQz3HQVK2ncO8139aFK/RT5Rpc/CdooyhK6UwfHLctccimQayAR/t9Szplm/H7XX5f4SOjo6Ojo6Ojk8wrpP5SSLzbnCp/xsK/7UvXETGBJFeoQoAACAASURBVO3Eia9h9k1Lg/GBLsgiw9Sq9OUtTNuD+cgAOeLH5YhR9sPOcf2Wsi+etT6zMSmO+TGNz7LfF+Nia77ITzottif0Pygr8UAM0J3UkVVVhBZZq63cIcz0cOmKJ2KpRAqbZYVNrXTF+ve/Wd2VhV5rM5gVmpTpmUwLpPfdnQ+qlZTIsJhMKW/0LeM2Y3dbxT1naouWFhnYYEHOZYBay0iOdBFI+9Pq0ioD1OhTsjFOGg33e6YLk6kPuC9RGhWpEWaAhlSr2dZKR3AU5thifnTZYY55vSYXIWi5eoiZKfuoGRvofxBdycWNEQ22xdBuRbadKvtiUZC5XjabJifqeaq8P64NR2OBLTLt5Vyul+UAAsFuobhghLwOSd/XVNS6cUL1b+2oXYkhTgf3PyD8/3mOPD835PbqzE9HR0dHR0fHTeFqmZ+8G4pWxhWog092zVgpCUj91zits8Z0WH8cWI+mB1LLH+yOWiKDM3lgVawxQMFAhVFvbnCwIMBya4MholPT9syxcOdR5x9nV0B1t+znUTNJPyrl8Djsw1Sk5MV5Ly86IOTWAQOEgqM+CsWWIRKlUXRxOa/6Ox2sDRge6HpYlyRSGCpjsog1ajFArWN6tCLpOM+SWfOI0pn9/QajB22JZvPVviWX0dt+UeQL9lFYigZnYiwOC8Iw9doGrIoMqLEsNK5D34jx4eV8qNhoHVUOoDV2qkEyoCk0Mqbj8fm10JZYIlxLZoD8fqtcQBtERysCczl+bk5FRHb6fOx0bGDcjo0ipWO5eZsIGZ6hyxNkbE+6/xyOEwuons4FZLiLs6VI8DKfeayIi/Q0RlH1OjqbGmxSyQ2k+7Acb9hHY9u1HECTDg4fvcsRYKz9aehzmhFgQkzryrpxeCbNzw0xP9f58dPR0dHR0dFxPrK0IiquFt3t1dHR0dHR0XFTuErmJ6ck891QFyIVKa6wFt/qEVxMKUxZAN0qxmiaXgp555IYvqKBRU5TOHyVlM73SVG0gOrOaYSQcrhpJnfXQd1dR+/2mqPb66Nx4bFf7BbX1odD4bXfVaEzlhUBNIe+l4SCXHD0ksKKE4W0c0LGlnD7xQz31y4s94VNS8h8CtNzsB7ODBdXOT9Q8Fg2z1SqxJsmGGNjRaLr6uj+Wpaxu5T9RY1Zc/FCOEr7aLmYMLZX3F0h3kBo2ZorufFcrjVtSr1Xhdqt1hGWgFQfnpncX8vWUQTdSoQoEh/TtYSIpRRNQ0g/xfHDJSxQQkakiJ/xHJQwdvQsh21bGIXaWqi9c3uRO+acZ3benyd8FhF5Mrd8dJHZC9EuetnGxM9wd+lGeA+ii5MTPFelWvS8BguBb4iX4aLWxWmic29sU7m/uAiqGySDviPnJFtD9LXhWaLrnwmd+eno6Ojo6Oi4KVwl8yNJZLobzGpKTlRcmJd1oVq9P/0y32JiRNrWKjM+aplYCLz7/EQ/LQlWdbytAOAFxeCFmLmss5TysKh0epy0UCcYoKlYrw+Tho/vlmUfaaj7R9MyfWcsQb4fjgvzc6+lIsD43FfMT11E9C4heVudCFEkCi5ZiGxlLogBaobs0/RxjkLoVptjjoyYJVFsMEKXsER8W1n4nIIqN45XG7ZEs2Q3CAsTw0wPWCIEBLh1SOOP3WMnLIR2t8fasrgUQtI6dqAWgW7RONJexk0DwzTEZfXT4jpFlr/RH2DmIIBuFNI8lwFqYa0IcUyhoAwPQtw5/YLUbZkdimWBIzjJYVk+hOVBWA0xtJXGuJwuYFbVX9u6+GlcPtP4EiksqgmeqUQGJ0gM28wogRGfrfCvAfvBPNh5Y3c2rgEVxM4YI2D/Iz0oIlry4jlYmBtifq7z46ejo6Ojo6PjMtyQ4PkqP37yIDLfpeJbdWHrzPxYePo5X7ymh4hhwM3xYv5k+I+j7xnHHVw5gsIsYeO44xyaVnZvmAWrE9LGY9kdim6qb71ifhxjgvB3ZYDe2S0szgvTALlQd9UDgQ26UwYITBBYHjBCfplZuDoPa7KVmG0NM5W78OHrlghxbmt9PPPzqG3A9EALdSDGx7M8nCjv46BIzMrNM50Zj42Bxq8vfks6oIoBksY4Y00PSxpYC9TaHVUJaGp+aEFOdduTYM1GQ4d0Mllj2A/uKzFAzVIZ9GxewACdYnxayTNHYnpQLqXFWOIZMgY0kQbI9WVMbfYGIfATmCZfQBVtSTRxybO6BWZ+DsYA0XPtx6DRRLqMkhxOxryXTSZKgFhritx9WEuASMWIgwZorQSGMaJzfR6hTtEN0TDPgK756ejo6Ojo6JCUX/3fWcdN6adSSv9HSumbKaW/2FifUkr/ma7/31JK//zHPderZH4kJZnuUkka6DU/SHOONOWWQh1f9Kd3XyVCNKt8ayNtMtPUF0Pl3dHUW2wlAoJ82bR/z34g2SPc92B8skZiTEfVtuydNan6n/24dPRpH7U/92OJ3LrTNg/G/CzzYH7QdtfQ/MCSNeaHEql5TQKXBwCYffEFWq0Y6hwjwYztchY0ND5o+0hRcE+NqDjsx6zteV0XtIbUpDCwbpkaQ2L6INKjeZEXa9UoeZ81c9fT1mEckTFcNDQ1C7KWjNBrilb4yoYIyP1usTUnUF1KnjaYn0JogHlYFlj3c6tT5zJAYhfzcGjrgcprpWZ+rJgunoupTj6IZ+dg7CkiwsC8Ng8bUBc6HcK+/TpOgNhKUHqKDWpqflbWHVgLFDQ/UbcDFq+OAnMHZ+YQ+4BXwL87TQhE0XwsEGolMFxJhNhKgmhbDxeyoJ9gpJRGEfkFEflJEfmWiHwjpfS1nPNvuGZ/QkQ+p38/LiJ/Racvjc78dHR0dHR03Drya/o7jR8TkW/mnH8z5/wkIr8iIl+kNl8UkV/OC/6OiHxfSukHX/ZURa6U+clJlPnReR9RhWgrRFaxZWBWwBl3rVEEte6MhP1mij7wuo7SN9oYGha324GbsG6BigKKiJXUmPdqnaqFM6ufed5rvhmnQ5qOaq3ulnVPx2XIgAm62xWLEMzP/W5hhfbK/IABujPmp5zJPS3DNkXbEK1Yj820+hJZl0kiEzNbTqC6wOkR68Dw5MgSMQMkUpieIzE+E+mDuLjsFoLmJ8fxYrtZ0wKJFKmKlcCItEetASosEDNApQZE2LUuyrROz5XYo7As7m4zcmsNW6UyOJ2PXRaetrbHVBPkWLTX6JvivtbHXmZbxVAXHGNTp/0RnTrmB8+BaePQNjKlIuXZGZUSQUmVgz74JUqr7B9sEDM+JSdQHSnGUV4DFUFtYuU/zZZWzp4ZnccrmaNWRUQm1vxU70W0cweodJiYDmEbv+Fg71UJ0xLJVbCaA2hFC7QcRyM9Uzrvf9B14IdE5Lfc/LekZnVabX5IRH77ZQ96lR8/HR0dHR0dHZcgnbDkXxqfTSn9XTf/1ZzzV+OBK2w5wtfaXITr/PjZYH44v06ijKCmCWppAjaOJyIx8/NaYE0RAyzzrh3yEllGWbPKan1BlVOF/OCWQdppQGaLblCrdafWCqZYPjqtjK6b9suFOY6qh1HG52nnokDUAgQbtAMTpFOwRV6/A7ZoR9of0zboie2HYoGy5ucUAyTiGB+6MQfL3VOWow0zPcwAeU2RLbN1Kc5b9F2tbSgG4Uu8eJjxSY11VqRUxzhFip3DAFX6Hc+mCvqPJZFpCpofPJNnMkDhNE7BtVvLG2SaqVS3zcRulSmxCSJ2sia14mtoZFirGOoCMEAsE0kNzU8dGQbmp1zc3Rx1cqYBQhSYan8md0CMYSvAa3mE1lkcnJExPyusUViHs6b/OFu5snBNOSM97sexpRMCa6PzFu1lDco2lgW6Yv1J6OZ+cwFVYS3QsdpkOwu0xHe05dXibd4UXg/Z9Hs5589vrP+WiPyIm/9hEfn2S7S5CF3z09HR0dHR0fFc+IaIfC6l9KMppTsR+WkR+Rq1+ZqI/Nsa9fUviMgf5Jxf2uUlcq3MT0dHR0dHR8dleAaZUc75mFL6soj8qixBzb+Uc/71lNLP6vqviMjXReQLIvJNEflQRP7cxz3uVX78LEkOHXvpE79ZsVMVM8LNBRqeCpEuGzWWueUGx+BWRUl5G0677tqym8uiKB1dWlK8Q+CsK3YSt3Hi5YR16spKqnbMEDwfoxvM/54P0VU27dQF5NoOCIfXBIg7dX+9gCtrxNSLNCHO1HVDpNAx711bLP6041/g/uL5uREWX4uXaX5uuL1W3F0zrffLOKnbS7nBjOf3vh9y17AA2oqWuvs9xGOa+wvjaaskBiU3ZPeX365yf7FbKvi9JIJvM7sitto2BM/VMYkPT3Rflt8YN+r6yby81RUOg19Q3F/Lmke3Dbu9+Hnwbi8ub7FPGnhg7i8NPHDFUGdTo7dLYGyVsJhSdJFZqPs5/oSN/zwTPZNrpTD883Hk96AlQlwmqRHqXpXAqATQ7t2JMHh6SAeaD0OVS2CsebLcwArh77cjeJac89dl+cDxy77ifmcR+Quv8phX+fHT0dHR0dHRcSFu53vrSj9+ksh0n4zFyaOzbI8Q/eo8MUEc9hiWVSLj0yMlryRA5H2JSFVKgAsqemvSLE0KaS/npYtdeK5Z5rgGsEzABO3UohqdpQ4x9EjzyvxM7toaG6TLDloMFaGwoy4fx1okDTYI4bmjMUDR4hVpWMH0xKazFbLu2robtGZproWvixSmJxvzExkgFmQuv1OYsqDzY7+HKjaIB6FOvaWOkF1mgIzVWWdzTjFAfjtjfFJsa8+dZ2b4NM7BStvWPjiBJDO9FuLeEMAONn5yXN5IT1E2jWMC5GxKdfLDdcGzPgOB+YlCZwQJ7LOmp1DB85OL2d/riwIi6FKsZsFWmRkTPiOEviF4tuc5xefZcMF/oJr5KevsuWMGDgldjRFyA4sYnqoEhmubmFGi93eL7GKWCAVNy742BvRzsD5Z5JZqe3XBc0dHR0dHR8dN4SqZn5xEpnuRAZ92XvOjloiFvIPx0fnBCp7WX/1clJQ1GS+Vjtz7e1cSdBX/cm1tmwU9RYYGiQq9MamGn2jdzsISqZ5n2MX1fn8Z60awQy19kK7T6z5pGDzYITBCw662VocR7JCyOUO0cMfBs0XRGh4buiCRWiqyhdatW2N+WkxcYXzaTE9Lx8OMT15hgNa2PxvEaJTUDDWLw2UyTjFAfvvTDFDphG1D+6X6vMsiJq7oZjVZ1Hi4zWeTw99Zf2T7cOdho5EYIC6J4Z/ZWjI4hONw8sOlT9vMj2c9OS0EQt+NCcLUaX5KAkRlhyS2bcESIOp0sr6ksNzDCqOumNshvUal+dGprD9//JzNpLUr79b6fcth8IWZEdcWy9olMDIL9txvGwFg2s9ggFJLNPYGcCslNUQ689PR0dHR0dFxY7hO5mcQmXy0l9e9gLnQr3B8waOkRB7o6188G6TbEEPDpSv8slXhBrYNy+LCRL7t7JNhWbRMZFtg1CUwM760gF6HAdFdY1yOazK461WYHl1nGqAU1ou09EHEEul0clYeWCH0FwU7wfykITJBIoUdslx92sYkLGQte5yjB2J2pWJ+VtqJNNicKpLLs3fUxlY0NA0VK/gyDBBRGrbLmsU5xQDFcU1szgoDFLpNbEtho2JX/brNTIgMZnw22lZDgvvYMOpNczfGdYW1w3x9IqsMkM5PvkNK3ZYix9DQRC2NSCkFsx8W5c6OC5ymmtVBhBbYoIMe7860QNDVFRQdUGSfoBc6BDYKpS90bSNJ5imA8cHzx1GYIjXz8zTHsWjsjnshTsTwpKptg/23cUUlMFrEFt1o1ohWDNDbgM78vBmklH4ppfS7KaW/t7L+lZex7+jo6Ojo6LhtPDfz81+KyM+LyC+vrH+5MvZJZL7PpvUJ5S2I7QALUubVinUub7AcVvqCtD+tvBBVPh/+ot6wSG0vZk7qcpeW3qxrFLqkUh6Zziu2gc5mmTf2Bqkx3DbQAWVqa9fLs0TUNlf6IAnzYd0Q52cwP8YMFfPoaDmUwBJF02kc64uKS7mDDmnStP66X2+UwYpEnqLDYaS2KezDrztqWzBYRYOwtPW6JGxv0UbG8Ol484PGjy3fX8vVgxUXMEKNpplYHI4Ys3IXns1ZyeuzWQqjUbIlzHnWi1bys9Rki2jjrSKoFfiSknUvUvQ/RsZSRBjn/YkH3WaAwq3W3wdmfkgT539ziRgwPR/Od2He/36hzIxFhOkJ7VNLiRTBGh9fDPVgYXwH7bee9QYDNNFCzFuUJTFBYZ1FeSlLRPdhduy5sTir2h/3rK6VwGB9kHuRDKzbwf+JuLSUtHDbP5fm55ZwkvlJKX05pfSZ13HwnPPfFpF/vNHklZex7+jo6Ojo6KiR8qv/e1txDvPzAyLyjZTSr4nIL4nIr+b8xpIQnF3GPqX0JRH5kojI7vs+o5ofMByOMYGuhfRAlebHMT+m+YERYwxQtIL91/qaPqiKAGgxP+wrRlM3khL9KFoKsF1qqftCe5TxGvPQ07TYnMLIqOVUsTulLbNCmaLHmBFafmdal+PUju++00kfhCy1mJ+g5xrqizsdoaFAW9zUVtuY9WSeoq0QIremqNso9BTa6j5DBBdrfTCtLduyEY8nmoa27f1ehEpvg/tS6yGqDNJb2aBtt2BGU+yiZ4s498+WPkjWl4lIM/HOWr4u5Oji2rEi5X6WfFs6BUPa2N9gHcdG2wyQX8dZoDn6S6SwKh8e7/R4y7q16C//+04ZnoM+tHeWFXqZD+wO5fPBcaBD8mP8HtnYKULMVYSN87T9Ms+ZnmvND+vybJ1OD4j+8rd6ZvYmvs9DPygLdGEUMUZqHoHHIMZZxXk2ipg2dtfxinHyEuec/0NZ3E7/hYj8WRH5Byml/ySl9M+85r6JNIn5NmGdc/5qzvnzOefPD++995q71dHR0dHRcWXI6dX/vaU46/tSmZ7/R/+OIvIZEflvU0r/6Wvsm8hrKGPf0dHR0dHRcds46fZKKf27IvIzIvJ7IvKLIvLv5ZwPKaVBRP6BiPz7r7F/XxORL6eUfkUWofN5ZewHEjw7AayVvNB1xd0V1w/eXYQiqBBDV+4v3cYVES1cuc6TG6wIS90mNf+tu6rbMqyMBsXQhw9vc0ucmjp3jl0vuBB1uYWxl92XdRLWDewyc4kRZ0qMWImksf+WSNpcJO158Yn4+OLR/dkKgQelzm2aCQfXdtMIX19r09zXif1uur9Wj3dBW4AE0MtuSLR8TikMu+4rbb0GlNxda2Hym1hxQfh1lZtCaFz548Dlzc8O3YcgpLdhWfkSl/WtbutJwrV71OrEOPdH7wqvSl8sL6j31Q22a4S63w+LEHk/xyKoED5D8Dy5k8dvLoUB3DVitzkxokGv25hbZ78ACRBZ6BxC3UkMzeVmSjFiPxZ1HYWxW5tQ4gjbR5dlLWtw7mAdFGAYKJNJW/Jv5YqeIf49y8u9Ez6hOEfz81kR+Tdyzv+XX5hznlNKf/LjHDyl9NdE5CdE5LMppW+JyF8SfaZeVxn7jo6Ojo6Ojgb6x09Bzvk/2lj39z/OwXPOf/rE+iwvU8Y+ZZnvZ2MRhoP7GgerYSLfyABZZKcXXJIYeiBBLUpl+BBVK5SKw3HxzUa6dQMxPpyEyy8DBrWcckVx+IPiPIgdIgF0SIxI4e+Z2DQfFl+VwCDB80DrRZzwHDnQxnhfWDwd1lVMD9bXprz9qlgEsu43sNnm1EvjIkamZj9O7vcS3/rWfu3QJPDk3ft9NETQIo4B8s+F3bOVruS43u2+ce/afffHXBOdNksX8PFot/E89JxNsBvvWQl1dzuqlp3BACUql6KbHFIU7vvfnAARoe8fpJbg+T4s289gh5aXHJId+m3A/GDKYuaWlmJovsTapTAm0Gl1nddlvb6UpiAMX/py1BfS8Y6YH71wj+4lbSUw7MUHxkfHvmtbJT7MsU0Jm2+cvW0brw/ffRGRFCpqv716mWvAc+f56ejo6Ojo6HgL8DaHpr9qXOfHTxKRu1my0gr+Y7zS/JC2BPqUfCjbgJ0w4wfzSCSIkPFghekysEJkJRnCZ79OSDdg2qKGToH1CmnDd86dK1qGNhMk4pgXXWbMDOYbjMxMSQ1nm+o1cNtMd7ot/O5gligh4uz0VMZMsd4CliKxewFm5OGcV5iNFi55M2wxJafwuiIkuA+tw/DwhBaH++RnjZZYYYDcTk0fZPc77qN1vExtVhkgz8zINkIFCdZLrdyr1JgpCSr13HEeVF7B75d1QVUFVXckVwI4NJ2100dXufhJrx0SIFqyQ0p6uBvubRus209U/HRW5oc0QCIid1r8dF7L2ucwVs8BWKkYAh9AYfAT3uNIdmgMWa354USImbQ+s2NzDpDX6MvGtD/6D2Nq3Dsucl0nP/TbUAmMldMMj9KxXtbxenCdHz8dHR0dHR0dl6EzP59wDFmGh0nmAyzRQv1YqQowDMreDKZ7Ub+s172gDdggsvoGs8aclkU/9zmHV8UAheiA9siz5SH6ABayzsMpvaZx2EAm9iPk9gITpieJXIMl2stZX2BpKBIM8/MO68v+oceaTfOjxwVjRqVFlv3ivkro2xojFM4NrB0lg2yaWhWjAMbsjDfEqzbd1pikS47zMdrmjbfi2cVQxUVQEQNUJT/0RYI5AmyNAdrov2GL3WGp2Dn/CHAcYrIwrnxEUqUDqvbfVIEs+7FzjQzQ0Z8nJUAspTCY+fGaHyp6Oh/D/B1pgEJbY4DaGqDld8RQ3fYNpnolEYslPfQFki3xoUbFQfszD2E6OTcAWKGDsXRgh+K7VUTqEhis/Wno9Aobb/8gdB5MVn2aNgJmeR7Nzw19/PQ8kh0dHR0dHR03hatkflIS2d8d5aiJZ2afI+aglhNMEOSQgOVp7IFjZCyaKOqBBmIYBqcTErKC7CufrL9QIJIYpeq8WlYFSm/McVpStudqmyrCjA4Y5B0WWRUZoAFsjo/csiivyADNpP2Z9k7PQYVkE1gjm5cw749ZND5YnsN80C6xRmNozwesMQvGlG2YSdU2600vaVOxTufoeC457kswTFUx1JdhgHAc1ve0+sBtGttUEWJbWGF8TkWBiTTGFZnzg3vWTA90MpovqEDCAUolBjyXjv1Q/c+TMdKR+SmFT2vm535Anp8Y/YWyFw9OBPlCadq7YQr7WMv7I1IHbrE8r53lCI3QImp+fLSX6YC4+CktD+VlMJ0jA2QlMLxYdKUExlr+HxE/ntqhjbgGjeovIsf0xpMjv+21uF41OvPT0dHR0dHRcVO4Suano6Ojo6Oj40K8xbW4XjWu8uMnpSx3d0erVn506uVJXWGoEs6i2WQC37K/gcSxmJcDjtfqRPzB7i+rYu1LMEB4B8Ei79O7sNjddYzzlqirsU1rf7HvtbvIXH40HXa1mLyEtEPMHJd78fJ0TzSy9hGh7WjrRdJwhVkGeLge6N75siZcLbw6rw3XD6cEMBfHGQklN91FtN9NmHtlzSeKnfr9vwSHferd9zLuL+/64VIYJPxnt5HfvhKV8jV2u17Lzbh5fhzOzG6wxvNilx0V4Mn95UPdyzJykWwKrPnG4sTUXe/OB78hgn6y95W6vyj5oYgPh5/CPELfH9LykvsglfD4/bi4wl5oZtK9+qSthIV/wAkjjd/iDvMnv+ICIz/F7F7S5u7axQSIBxI+T26bkgAxVomfSAC9rFumVQkMEkD7MTJV7s22o8W7BOPwf4YPke726ujo6Ojo6Oi4Tlwl8zMMWd69O8iTWjGHsXxbH9QaAgM0oczFQEJoxxaZGBplMigpoIVr+06smJwQQM4koBPxouW1fdSohc5q7R11J65AXgmH55hbRcPqLgQWWBWkno8MkIhIVqtr2MdQVCSDTLrcG4a7F0sfDu+B6cE0XqcU0tLrMkvAqNOZ51sMVpy3pHo+tQFddwtfNqs7TMLMWuLIJqwtH7BuusoSrRADy8/tAbTJDF0ipCbmJVdlWdzGeN5WlNXNPttw3WZMwmW0+4s+bPR/7TJY/+vQ5yqlBBe9bY3B+IgWMfRZrga60Q3BM0VU27tugOBZs4uC3RER+XBUxsdC3mPo+8McBdAiIi/Skpl0r2zRQSMPsK0vgjpb6YsItrrv3M0rJS/iTTOh8/BYHWcSYny0T8ddZHV86gEOf+cSGHNgfqC6xkVuj8WQDJTubxX6Dn7HDZ2Q3/EZiJ8ueO7o6Ojo6OjouFJcJ/OTsrx39yQ7tXQejzmsExE5jsoKGQPU1gItv/XLHSxHilNY0ENLj8H+fPINe/KFWQ4zKoxp8IIOPo5ap1NkfNLRWU+mE5rDNsWK2fjsh6YBicVw7l7zgzQCR7VowARNamHZ8cs2x3eW3+OjXkM1ETlt/Dz7C6XHBvOD6dyeF/FMj3afw+EbepFSsDOccttQp6HBLBG3C79XWKJmQVAuTWKMzwUm2yUW5VrbzbhvWux1bdTxigGyMbjRSWZ8iPETkZJ0sGLkoM1x/WfijTQ+raSHpxIhlj45vRMxYfa803E8e1Bf5sg+e90ZkqxmQfqO5UE4kFZtdH0qCRA1bN20P1TuIhRDPYYpip+iuLIvVopEiJZZhM5m1PMJy7XtQ8UAHX0rmX2Y/wBNj17bMWp8CiPk3j36opiMAYqMz1OLtaMSGHbvLOTddbBKhIgVVL+jkWdjEHkW5qdrfjo6Ojo6Ojo6rhRXyfyMaZb39k8unbuLblD9z+NhOXVEhB0sMkwTIzqTK1tKeYqwMEOXzXC3EnNErpgUx5esmKCjiYxDWxuCHcVZ0ydA3+MivJJGcAizQ8z8tBggZrvAhh3K93PaKXumTE/e6zxqVhjzU3Y7KCs33eu1RaY0iq7wuhH018pnUNJDni79bU9bbI6xQhzBY9vkapvVyKNKD1NvwsVVS+K8unGVCJOZpi3T7eMwPmtao1abTI2CDIK1PdsaoOWnrrMEpLCkMRbrTeweEgN00TVYYZiWY2NdQw/kuuzH4Z6BvAAAIABJREFUbRl7UQtSl7tYf4+UxY37jUSINl61tAPEcWCoXWXhuvTFMr0bFpblflj0PQ9zSXJojI8lQozFTw/upXZINRvUgmfN99BE2aLIAM3KJkH7I+KLnQ5hCu1P0QCVl4KVwLjjaK+oARIRgfOAS2Bw0kOPMkYii11EbK17O9iqzWjS14F8W5qfq/z46ejo6Ojo6LgQ/ePnk40hZfnU7rH4scfirx6PC7WAnBdPU0z9brmBUrEQEBEGqYr51HV9+UCvmRlYp3NlRarV57QsxmRYzpu4PHuX95pRQNRSctFeYHzScQptjAlqMT9oM0QPabKoL7cc+1XrCuxTBuOjTJC3LjLYLjNk6TjkN/eNzarXa1gILLrm4jU+7WlwALNehCKGmBFaZkLXSh2TNSaIf4vThDSYvooNWmWC6m1keIVvtHMYFLuXzASVlVxQtvAYqd6GmRE8D9j/WG9jw4WZPruHrkcnDOyixXHL5nhuq0F4bkXRAYEt0HHLkZ/h3LdYAn8kr2taJjMVQcV77MlrfrR/e+T3ORH95X8/aD2fku/nWLfVCw/tz6Qn19RHKoZE1wH91wUPUucRmtPT0kbP2TRA0P5QmQuRov85QhdkuYDqtriIVQkMbEssT1xG860xXmF4Hs3PDeEqP346Ojo6Ojo6LkRnfj7ZGNMsn96/kJ0W63uai6VgFo1GM7wAE6SWz6P6w7119GQ+82XCxQUt22roBfy8ZLmxfqChZUkczYTu+8gRs+plGy0WB1MwNaoFyrbem9Bc3VEPaBFwrgOqpwLblI19QnXSXHX5uF+u/3DgyI54Yrmhd2KGB5dnxnF8EUNieNaYoGUlLSMLdJPFwZAAO2XL88ltNvcPUEZq0x81xgETMPV67OuCN96W1me1A0H0o9O4zrJDVxom3zjF3Y5xdcjVNMdldp+ZKWv1dw2BfdRNKGKryv/TmKnGHEWApibLE69BO+Iw0mglwlGfE9UzHt11etSHZqfRr6wButdszr4YKhif/fTOMrWIMGV+XBXiO32eWfuDW8fRX8uyig5cjltpgMpxJu0Ta39Y83MYy7880wGp9gdMDzNBIkX/w1mgLfNzcQsUVPmoiMU+kd+pWWy545XhKj9+Ojo6Ojo6Oi7DLQme+7dlR0dHR0dHx03hKpmfQbJ8anw0sd5jw+3FhfwghB4sIZgTKupvC4eHAA+CQndkwHTHLG5jutx5lVjoPFSuGefGQXj3MboCys4w9YJLKO8gdFZ3lwmgp2ob8x+Ybym6+nyfkKEwTZju4rYN7D9a6Orj9z3oEhDhbfeXHnw5NtHH7N7JbgFfwzU3WFjHrpMtVyO7qipXWS2Srkps0D6a4ffcJ7iCEP7f6luVEXGlr34Zg913l5hNDXdRqcUQx1Upjuo3iR3N7E7QArY5lNGgQ289S+YvjV3b1HSb0Fn7i1fMGZazuSrJbZfZLeI6U4VJt/x2Nm7Y5Yqkh8vs5NznB60S/ILdXlzw1L2o7ocobL63Ehi63Lmj7gYkQIS7C20gfF7aDY2rbQkQqVrsXeX+Epn0mLMWYjXBM0piUPkLkfJuhvsL83CDebeXFTvVi/pIbi8Li/eieBsjfA/XXhKENx3qfmO4yo+fjo6Ojo6OjgtxQ26vq/z4GdMsn9o9yv28WBkfzne2DkLnvYqhYeGA8YGFMzrmZzTLnEfGcvlaRp8VwmPGByLdBvNjQmcU9bQin8v6YfRMRgpT4dIbTQYgr0yJCZoadBTTKnYcl+RQw2Rz3um50jYQSTuLZvrUYjUOTxy+GhOztcpPsKWcybJqh7pHhsQYtHDzaBsioVqFZ6uElOeIpLE/S9aXw/kFlkjiMk7ex+HTvk8l92DiBXWfGGsFWv31OhVKH2gcYjewbo6NsztAMfR1HZ4DYoCyez6qdAjE9IXyNXy9W/dqDRWj2xY+B+B4xPQiLqCZTHH1Bnh6UNdU54HxxUx1YYEQ6IH3IBggvC/vnOCZmZ+9JTuMyQ9FSukLiKCf9AbsTQC9YHYP4LjCeOwr+rOsQwHVSY/9niyh7zOFvh9c5lMUvDbmx8LZVQjtBokVRqUEiE+WsgTvHPc+PLMIanxPOkayEz+vFVf58dPR0dHR0dFxAfJtCZ6v8uNnSFk+Nb4wf/beWS0fTgsLhASIpgHS6YdpYYQ8yzPQiCCDSpgBEpEqDLukQYeTGPNuvxNNNakhin2ijINIScbIIbxchmILmRkg0wS5M0HY+kq8tGdxUMbC9CdWykDbmq6jYPwDTVn/rjJAVfFY3bcT5Zi1SKUwMmtbNkLdbSc4L6cBWWV8eLdRmBK6XTFB3A+/P2KArI+NN5G1jYet2R0pFmVVK3SL8VllrIip8cVKuc0mSKcF67himLwFzGwKsUfGCJWjWELBkZo2klqW4raReTvnNNZC3VND51aPn8gMWNULvykTrqu6Eb8s7r+wnfpOCnonFHZe2jwpe/aRMT+NEkEnip82i6BqAkSEut/pwzWQ9scvg+ZnsKSyYOmT7rNsYwkQU7sEBrQ/h1z+5aHoqTFAYHz2NfNzXEmACAaIy1+IiExcAmMl7UkccO4BfA7m54Y+fnq0V0dHR0dHR8dN4SqZn1Fm+fTwQl6kmHZdpFglYIAwP6SZpjXzAzaImaBiqXmrAiwHiubpPsgPjGKmyzoJy4Zd1P4MrrzFvFNrQqO98k79+UdNLGiaHK8JWDEliKkRVxLD9D+c7BDrw0w0Xc02tUinWvNTroNeO0SvQcOCc3cFWo0sQL5FS2yG+wSL0elfKCLMZDvWZd8WKxO1wXLBgd2564RYImEGyF9Gvh3VfMsipONgnrUtvr+mMdDZTY3Pxjq3Ppb2yLTutPlYSnkQe9M6PoYgv61Y8+M3MTYlan/ABIVClFwCw9hG9CmO5yboGtsl8NGWRG7Zfoc4NsOjtpeIKnKvHiOZWNPCohI9KSLzgEgwZT/wTCnT8yFFgYkU1vyemB8kP0T017IuRoBBP7k35keXu7tnv4kIHexXrfnZW1SXMkBWTkMZH40Ce88VQ11LgGjzoQgqmJ+YANFKYFD5C5Fyh2am53mMBJT7+yxJDjvz09HR0dHR0dFxnbhK5meQLO8OjxZhgMJ7Io7p0U9c6IHA5uwb0V7DCuOzBZbRWLoctRgswqeW14gSVjJZtJcyQO5uGSOCaCWK/kIUky8/kRvRVqdPBKUqWmapBKohqbVVZAnETqH8hYsmy+/caxs9zgiGDIyQTo/Ogh7AhOl1McYBLNEyH3saLWWrLEDsSPhtIgwsJ1PUb8O6IG7T0AAZ0zDRPB1uWdbSeDTYisAsRa1B5vNolFPIp0y/1jNgTByYt+1dLJtwnyIrEcqZTHGHXtsTlvsLZhVV9DgWlaPLG3mdVoveMoMisq5vogc/Rm6RLgi7QwFgenZDW9MW2d5iA7+MGJ5CPMT5sM4Kgar2R19C0BZ65sciwMYpzHMUmEgdAcbRX6bvcRcKzA80PTPYeH2IKgZIREqJIWh+IgPE+X9ERJ50ICH3D/RAzAiJFP0PcgAZE0TRX74YatEH6QLdh+P6ZR1vXvOT5EzZ3pWgMz8dHR0dHR0dN4WrZH5SyvIwHJyVUcQyL/Ki9SnWhea3sKgDTItVwevWEMgDK5K3zM+WIVTZCsvA7DaiCDBofMD4QOcT1u3BgsBqNKe9zjs7A1YPMUAhS/Pqya3k+3HUFc7Z9kb5RRKiyI5OvPT//uNl+kc+Q21I++PPAxolY5bC4cx8aTJclXBnHRWLw1Fkl2h+eD4cSKdzbNrM87PGAGF9q09r2p/WcF67PDSfxrlaV11/zLai1jBWUpxv5nVqdDMsJ5ZVxGl9EO2FawuWzbFHZR3Ga+wDyTo2O8OanxQ0P+hoexdgzkKmantvNApn+k42+slRiSVDeb0N3hc4dIn+0gzQjvkBK76jIqi7ZrQXR4ItD+uQqdCpG4wcAVa0e5EBaoG1P5jeQ5ckB2sL/Y8VQeW8P+7agvE5UCQYR/V69vEJRVB5fJKYJxSydePnWfL83BDzc5UfPx0dHR0dHR0XIHe3V0dHR0dHR0fH1eIqmZ9Bsjykg6U2H5wKFPQqhHh7zZRnbrCE0MsVbnoDXuzGodUsuJz3cAW5sNN9DsuQxA/zkwv3LgkRta26xIYd5vXcj4XfR/kJK0mhNC/cUIncVMtCVvJyGYoaRRwN0SfKZ+h058yLz3zvMtU2EDijxAcLn0WK+JnLEZgr0ETNZZu58uOQC6sheK5cWGtCaJHiSqxcZbGPTZjrwTpb7aNyVVWiZVrv+0dJ1i4y7ugaV24qKaLYlOK0uB9buyW3nT0ng++y/oabs14n4t1f7vmz9A36PGBM7ur7ncgVlkmAzkkQ/aFYcM6dCqULpnivWOTdquVbEh9S1VVZc4OJrAmf4RH3Y9GWccCEhr4ftQMvhhI0gjD4/YiEscs8BM/3Q51axELeKfQdZTP27h19qMLfOY9D/W4e7Zz1uOT+ghD6waWindNSAgMFqp/0nJ8o5F1E5LD7SPcDOUMMdT+a8NmVxOB3Pksf+L3iT4B+vjF05qejo6Ojo6Oj4zpx3cyPFdErlsgBBfb0q36wlOm10Nn2dyLUfTaWp8X84OsfoZJx+exVmroOodogSrjshYjIfATDo/MqfJ6PyihBab1zoe5gejR5VzrqDrEcIlZf3gJWLxYV018qVMUqiY5A2LwLdU/f/WD58e472kQZIGOL9Hyc4BLi58L0YIVOLIOhYyckMgFFe7ougF5jcdrLSYQdF1dh7JvrGgbuai/5djSYn6o6xDnWXZWtLy5PvpgpMT6pYoJc0ypBaDwjhDXn8FjwBakLdPJGLEQd4hAMBGbF+CB9xJjjcn+/iRE7B8z4pLbpX1THIqXYMVNNei2CWDbxQG0zQJFViyJvK4GhlNAMEbArGvs4LizQR2MsfXE3HnW+IXhGGPwUWR2EuHt2ngXOuD6c/HDcuPhYd2ek2rJtZH6WlU+anPE9WQTQFvI+FubHSmFQAkQOeZ/cxeX/Cwd+92xmMUw9yeFrxlV+/HR0dHR0dHRchlsSPF/lx0+SrEm1YopzEa/tiYX2oAWC1iekW2+wQR5gXcIysjxnsuAmJMUKhU1hRWC+rf0REZmR4A8J/YwJUubBNEXFfID+J03EACl9lNUqS6OPA45xmsViboS+n8pu1yj2aMVPq6yQLRGLHoaSwhmzQUxZfJJh9UatT7m1dds1wzJvrYY+RdfyJQlHYY1PPHp0SrPkA0bwHM+r9fJqskJrnVpDVcLCrcKUGB/WAvF2rYNDbua1E8mKVSIcGzcazxTCjt2+6JxNc6d6Mx9Obs+gvg056WSV9FCcNsZOPod527t/dTDjM7fvmddDWaoK3TEI0GTvl8YgkUglcpmOMVYR1SkYIKGpvpOG8k540v18OMTSF/thSVi6c+9LJI3dV0VQowYovG+NVdR1Kb6T99WDUrNAKH466DWHYmlyx0HfEP5+SMs5IgQ+aH4oDB7TSvvj7kdhHyPjY0VQtV24h37cPofm54ZwlR8/HR0dHR0dHReiMz+fbKSUg87HJ8WC9QA2CFbGaEm3UN7C+aBPMD/A3PD3WopzSoPOEQAiIjOiQRDddRfno+ZnmU5U9HS2ttAPOUtE9UAZzI8WE00WaaUMULCgc2hji5kBEjEqI60lTWwsT6r1OYnGQ2nRS+gCSnpYOQcvbtBtTA9E+w2nF5mlaltYxXWev1rsU6/gLtUMEPfNbc5D0ZgGYoCq7X1jY9VaBzoPfogww4P5MnV6rVP7pX2KlGdkArs1gQGKJQ2i3ClGQ9nzBsbPaVhM40MFVO1ZMu2P0/RViRCh9TrjWlrpi3V2EwDDg0jMIhPS5zu0RjQc6aqMnSINkOs3yIexZBbUaUx+KCJyRAkMjYr6yJifqAESEbmHDui4POec9BAJE/37Fklp6wSI8ZrvG8/HGkZEvHnmR8fNnd74hwQGaDmvd10R1CfSARnzM5MGyLFF9j8AEYyUffJIjJBIYTdds47XhKv8+Ono6Ojo6Oi4AFk683MNGCUbY7P3Cn/9sjbmB6nUEfU11NEHp2D6Hmd+FOZnoHkwQbA2nTUJjY/luMG86Hw5Jhc9nYjxScryzL4gqKZmR7SVaWeM+al1PFZKAm345FvmCagMy4/CwhdX2PQ7v7800fIWFVgTJJ6ZivvlXDeBJSGtTzGhydIV8c742IZPtSaWahYH8gveV+NwG2lMSlsql8FMVjPPDxv8rIdZP1wFsDujK2+RKqYnMj5j0PzUOqAWvOZnAquSMI8V2IeWI2iMEZaSgbCZQ0SVbgN2SI33hCgwTJsFQanNRBfbjytOkUXlLlJL50aFXgc7r6FqWmNpY5VhiKUScfoa0jCVc9blrrxMflS2w3IAQUuEaRkbd0pJcy6gtbIXIuV9jf3hnEfVAA0Vq+rbxOcL87O95wv2pv+E9mfZ/0EZoKdUWBzWAR0GZYJ2YII0mtcNkqr4Kd1e03+6PplOiM7vTeGWBM89z09HR0dHR0fHTeFqmZ+Ojo6Ojo6OC3BDzM+zfvyklH5KRP6yLLz1L+acf47W/4SI/E0R+T910X+Xc/6PT+5XFpHy2LiTo9KsA7u9EE6Z4/QczI2w02kfl1UUKERwjt4/qtsJ2mILU6fQ9+W3rqLSF0XorP1wgucSYqtCZ3Ov6TCAO8y7DyDOxDxKYeAATghtYeug1a2KPHw/VFZDRNJnPy0B51SYtwPqPshNZbR+ywXE2skNd1RJNhgF0Cx8XtrE3Veevsr9VR9zoyerAu0q9L11zpj9GC82LiXixf1waRTBc3R3ecGzuTJOub2c4HmA+5eTKep4OkqNkh4vuoe4zIxISTSKZweC50zPmBdJI9EmPCOcXuESsWpJ8wC/m1+rfTN3J/xs5M+TmsYvQmcOfXdtKMQd5S7K8jj2RURmiKAHcn8h6eFYSmHs0sMyNVF0LB/UKifECRDrsh+44+Xc7xI/IG24JB4m4obwGVKIh4bbC+dsFd+RBBLlLsy15cpbUKg7Jz20qesTPLlZxtvyQT0Dnu3jJ6U0isgviMhPisi3ROQbKaWv5Zx/g5r+LznnP/nGO9jR0dHR0XFDuKXvredkfn5MRL6Zc/5NEZGU0q+IyBdFhD9+Xgqj5CbzA4tzhNgM1qklN2wwPivKKHzlTyNE1O6r3yyCmAyLkxx6ZmY2gTNE0bDuIrsj4sLhqfTFzMudsBOsE9iHWUPeB7IivdHK9pStswRzXhxN106tpTRGszJ5WuSjF8v0nYfQB7sNNnXHMfYG4s8oYrZQ+ybLEs+kCIj9eZDZbrMstEyNNivz2IMbS2usUzMCnUP0uY1NvVlPiuY1NqL1wuPjUKmKEL5OunZsYgyQY3FGZoVOlIwREZkGPDO6jY4zJKUDE3Twp8TaYTCtmuQwJBe1hKPxeUtIeqhEQ3Zvy7r4Kaa4QFGovPRzbVCsCJ/DiYCxapS1IFSvq7VwdrdMiOGpC55KtU3Ge0+Z3EcVAft7ivu816LKYAl3VvZCmXgXnWBJDbEMBLKlI8H+fSCLHsf2kfymTWCdMU1U3PphKCMKTD3C3xH6/i6Vu/D/A47E+jPTwwzQMuNE0F2R+1rxnJf3h0Tkt9z8t3QZ419MKf2vKaX/PqX0x99M1zo6Ojo6Om4M+TX8vaV4TuanZbrwpfo1EfljOef3U0pfEJG/ISKfa+4spS+JyJdERH7gh0aZJDk9j7NEzJ8Mxkctnc2aBYoTn4rNUHdof8g3fECh08n5ldXiPEB7AJc3hbEvy7S7qP9Joe+wZpMvE2A6IL0eemxmUgILgqlZj7oEBUddkVJhnRSHvFsB1VEYbDGXJGsSl78iJLKo/cjLEtcx8WNEimeWjG3SNqbxWe9DFbZOep7YYVq3yvy4/a+ta7S9FK1ipQOVOwDjA6tfpGZ8fFi0R9AUISmcbgMmCJf86Yz+Fr0F3URxLBCeGbA6mFLSQ5Gi/xkoDL4KEfcsy4lXTDPUnZaVsTjUbZl4teOChaoHJZYVxicyQNke2ZotAvODpIpHZX6enDbqIzA+VPqiLnfhxggxMIWVz/HEAvBC1DkLfT8N7A4h79D+3Lnqtwel/6AH4tB3Y4BcGZCjvudY/4n/BawFiti9eR/UW/6x8qrxnMzPt0TkR9z8D4vIt32DnPMf5pzf199fF5F9SumzrZ3lnL+ac/58zvnz3/f95wz5jo6Ojo6OjlvEczI/3xCRz6WUflRE/pGI/LSI/BnfIKX0AyLyOznnnFL6MVk+1r7zKg6ORFmYzuYjjv7flY0DkMjQF1Dl5Iam/dGCo7aNY2bMOgUDRIkQZ6cP4uKnxuoY4wMGqPST2SCrVQmDirQHy8rIfphuQaO+skt2V7J4oXOR+Ukc/SUi+V3V+uDY1KZYl7UFnSuRycdgh9ztXmd6cJxqk1oORAxQMxrLLp1a2xKvtTcILRFfas9XjJD4qLs4lkua/ZoFOckKIZrJJwncRVaCk935JIcjlUDgosFghLw1jP6OQ3xmBor+GtKG0YOhbUSKG4NgWqGbw/nQs+XLy2SUk1FGxBggY37AoLg+2P0APbHWV8/mNNggtzzl0/arPRY2Ldvg3rBmaS354bIu7hiax6yMz2Es/1ZeEAtY7n/UAHl2fjDmJ8bxrUd/FUyIImskQmSskbLMPC37W37fQQ9EDBBYnYNjtUsRVGX5VRfEzE9uMD85SxUx+rqRZNvxcW14to+fnPMxpfRlEflVWdjJX8o5/3pK6Wd1/VdE5E+JyJ9PKR1F5CMR+emcP07AbkdHR0dHR8et41nz/Kgr6+u07Cvu98+LyM9fvF9ZvsRLJFdtVQBmlfInVbBIW5lERCYUVoRF6kwRWAD4yq+sADBBPo0/FT9FhMFkGiCn+YHWh/UKrAEKzA+mRF3QuQ+rM1IitjA9enMY5nU76sumY9lp+sP3l00/8z3LdDeEaSmR4fZHjM/HKgDYMBBLtBUYjRil0c7zE7USzADR4ZY2xApx9Jq/d2uMz1rJCv+7qgu5MuQ3Qbdhty/3fSCr3kpgNHQ9d8jzQgwQP5ceeIb2xqIuF27UZ+is28+5VXxZGUQ/WrQVMUEN5iftosYOhYYxxK1YqtfKYBmezTXT3mvJWAdktVuIvZP16zBQmZnwvGA3YFqrPD+6D0eqlXPT8fqoDBDeg+65eEp77QMxP5T/x0d7mR6MWEGLxEVf/Go95oPEm4WSJxbZdcG7wrP/d7pfaH/2yjJaTqAhRn+JFEaM88Cdo/nJUp6jN4obohZ6hueOjo6Ojo6Onufnk44sSZ5kLFEC7suatT2AWanQAgXBBX4QA8TaH5/nh772wQQdMxifmP9HpDA/nAPIor4882PRVlHjMzEDFJgAmHmYb4lX6BQpiUs6wELU4zv9jpnO7JmkjM/ZMT/5e97VbdGlIWyTOYmMSNHcvIRT3JgTHIb0Njqj6zDbZoCaOqGZNh7omns2h7etOuDWsKaHC5qWztYH2GKHTmHlEvv8VChymojp2dHUrwMDtLM8LzlMPdhCPqrVzdsE3Qjth3OqTIH5wYHAvC6zSlqUzM+e5OQs0Dvo8rCNTo/uOFzotzW2GRztNceb6KPJ8HpK/PwdcBjoedy6KsqL+trI8Cwr7FC2zM+lKdiPJ9UBvSAd2Nhg/hB1ZdFeDV3QsrH7bWzsY2iC7M2jaclOD37TgwYdErREUftzMA3QEnPoo72M7UfUlzATVDM/2S07lQX9FpBS+n4R+a9F5J8WkX8oIv9Wzvn3G+3+oYh8VxYl3THn/PlT++5plDo6Ojo6OjrkLczz8xdF5G/lnD8nIn9L59fwL+ec/7lzPnxE+sdPR0dHR0dHx9uJL4rIX9Xff1VE/vVXtePrdHvlhXKESG0K9OWK24s/UYMCdqZl0f3FwmcRkSlD6LxcYtC/B0qH7t1eVvzUqPkogI4iTYiiEfK64u5ylKq5xtDWTnn9GxibG2PP7i7PoU8rbi/rQKq2SY8aMv+wD11hDWDMAP8S7i6IiSVeJ7v83stJguA191dMCaBT8iRWJSz8pcZxaHhV4et+5VqSQ2nMnxA42/Jwsdv3LlFY+ehSHJRQc9FpdGl44Sq7u+4qt9d6VkgrEgx3mvqYnrQwbwhsIJfBWlFJEVdMEi4HSnaY9nruPskhucQQ7m/CZwt9924vne6i8FngGoNr2buYLEiB1es4vpuH6B7rJAJDb07ltT9Q8tKxCtWPbjF/TnyOJfmhP2d9z6H46QABNFyidzot931HJS+s3AWNzTBW7ORwjouv7wECZSteWvYx0hjBe3Fu+HpLCQy45Ej4PMT3vYgLhOG0J2M75J3xLG6vt8/T9kdzzr8tIpJz/u2U0j+50i6LyP+YlpfRf55z/uqpHV/lx09HR0dHR0fHBcivTfD82ZTS33XzX/UfJyml/1lEfqCx3X9wwTH+pZzzt/Xj6H9KKf3vOee/vbXBVX78ZElyyDt5IbU1CfEzhM9ryQzDcrOyo4k+SSyANzul39OA0PZHbQvLYAxT/9VfmJ+YzA3T7PZ/VJbFiAUOXwdL5E/DmCNYkSu2YSMU1qxStfgHWHfHWvCcpraVaiJmZxmirYW4jzxNYX3sE7odrdeLQt9hdLu+VskG6XAVA+TbrCVCzLTeA225jX8T4d7xLay2cfslRsPYLxI+hxeeFV/cfgv6MNxEy7iERcuqZwYIIlfMtzCTBX1Qy3o316JZDhM25kfnYwWJ5TU4V2VllOUEE+QIX2OF9A06oAQGhM9ggHZuXGGMY8wjLQKeB5z6+iVYZ1UdkH7CxiQJq31RWiHmB6LlVQG0FMZnJNaUhc8i5bU0432o8y8sxB37cqHu0h5HDF+EGmMDImi8bzF9EDA1fpvL/9NzEVQIn2d9iH0xVAic+Z3/3k4TI6KPjRfWnNPzhLq/Hvzelg4n5/yvrq1LKf1OSukHlfX5QRH53ZU06UXFAAAgAElEQVR9fFunv5tS+uuyFE7f/Pjpmp+Ojo6Ojo6Ot1Hw/DUR+Rn9/TMi8je5QUrpvZTSp/FbRP41Efl7p3Z8tczPi7y38MQh78tK9QVbOCNVG2wxQbYMUg+1GO50/7OaavAzi5Tii7OVvogh7y3mxxIg5pgAkbU/IsWaN5IFGiBiAlJD21BMtRV1QKpZFvPzP0URUNo57cdx+Q3LlkNuTavjZULvLUMQBVIrBgj6iFBJM/aJQ265XdierNSNePMacRgEdmQ1EaLtN1r7YX9rMpfQJzSO7JZdFmKPWvvP3P9mfDwWxTHCKfh9qPvdrh2uDq3GLnnmJzI892B+rNxB3JcHp4/Y6bQUKfaan3XtEJ+PaeruwPzoVLU+Vlh47+63LdN9HEn70yiGOmhSQJTEQKLQjOSQlgw00LVrJ7C+DCytsk81A1SoJa64UZe50PmQ5FDCRlz8NLvCpiUhor7LdF2d/PDOtimMIQqbRr3NVkJMYLJrqfNgZqS8o5nFaZVWWYOF32f0KWqBRETuVA90r2wQWCgrco0+uv/Cnt3soe4iIvJzIvLfpJT+HRH5v0Xk3xQRSSn9UyLyiznnL4jIHxWRv67jeyci/1XO+X84teOr/Pjp6Ojo6OjouAxv2/dWzvk7IvKvNJZ/W0S+oL9/U0T+2Uv3fZUfP7MkeTHvzbc95mI9wXqwr34KheEIAI9xxVSH9scXNn2AP5dKX5j1StFfIo0EiMT8TA3mhxMgmu97jpb78pt1QG1VS4ysYg1A1O0MLokbLNlSdDGH+RbGj9Qqem+x/MD04Po0S1iwxmeD8eFlZ9SBrDel7rfIorVEiIk3CteCdDxoisWN8hYW/cOMD5b78+OEjnO7bYstYrD+YL+fqnU8ZQtepLBA98b8HEObfYP5GS1SJz5D9zp9bDBMa1Zzq5hkeZb0OWMGCCUsPIujp2RJDhGsiIKnpP1Z2igjg+SlR32mTN+m/ffPFGlyKiKo9WxZ5JMySmCAGkkVk+mAdNkhPt/o2+gYy8IO6X2maK+xkRARCU7xWpp020cVS4X7TYkPtxJgMiZSckCL825SnY17QErElpZc2dAATSvUcIlEw/+aso89FUGdLSp4SYg4jfXLaNqV42z9L+r4+LjKj5+Ojo6Ojo6OC3FD31tX+fEz50E+mO9LnojA/OTmFFEA+MJfiwLzGGhbn//HfMz61X/QsJCDRYFpOvSxMD8WxUI5gGbT9zj9DumADmbR6r5M8+OsC0uVz1YMLMKoGfCrzKqzPCB6fGepDAe1/MH8nIj6EhE5fuoubGMWpzEbNRtVIlRiHxkX5QPyMgvoc1Z23FLKFIaHhDXch8DMkA6INTr+NnAEGDF/ZtC2iqHabVjR+JzB/JSOLA2OxzJu92Nb89Mqb7FnjQ8xQLCW92kr5GkBWNMd5f1Z+rCt+ZnDs6QsDRgfPG/K9Niz5HROVsYCWh/W/hhb5I6p+4P2pyp+armzXN8tXIoovy3ND567kbadls6kFjODfpvWh9lVp6caMdW20Pg0mFgQ27mKKtOSD8YAeb0WWFSc8sv/R8Z7GLnYHpR1EXHsvLZhr0B7fysMkP0vOFZtJz135ATC9F3tSxiLxm4O9h56k7glsqlHe3V0dHR0dHTcFK6T+ZEkH873zcygPlOpiLMQ+Ys31WxRNa8f7Ij+ml27O4sAixaHWQMoeOo0PxwJBuaH8/+IOOaHGJ8jkS3Rfj4V4gSWxUVrcNQHZXwF2yNSojxg9doUTRpWRcX4MPND0U2h25VWhhY0jDQmBJoaILu9bTPI8gA1qB8KHixaoLnVGTpgFXXXODgfhxif5vkYCxinL5HmxOBzxbBlbgyQxKlI0eUgqguMD/KjtJifNRbnoCFVeF78NmbFr2SODpofWlYyqi/LwdjMrWgvG+too+tJ+7P81ilpfwpLpH0NlUF1/1bwl0VfGzcRGddHHEcjKkO0FwYS2LvI0Jh+x7M5KFBM7wTMDz4yjJgjPN+zFUHVd55nOU2a1Nb8XBIFZaz5AJbHsSz68NzZmFtuWon+Op1tfAusK93bmI/RX+82HvR5fIZoryw35fbqzE9HR0dHR0fHTeEqmZ+Ojo6Ojo6OC3FDzM9VfvxMeZDvTg8lEZU06EvLJwYaOSY/DP6DlcJ6ZVcsfK7D301cp1TrYQBl/+i2iS4xhLxPSpN7qvVYub1A1cNFFpcvy1jwTAJnc6l4DjqHZcX9pefuKG6EvQ/HKHxmt4uHJUZEckNyaZkY2DPDK16ilwlj30TlCtWJZZL0TVNcxH1MjU5XImVqs+H2MlcTXDTWN9eUEi2am4vdYF4kzS5KFkujVMaGF28gt4Gn77GM3VvsEogurLYbYlY3BdxeL+aSzBQu72GKN7EkkXPPEpWTsWcLIel3Okad63Imd5e5uUgIPbvjz+QKQ2mMhHmEuLvyECZ+RukNuLC2fJaU7NDVwAl9FikusBJST+Vr4KZy4fd4vkcLi9d9qdt7HOsxXhU/tSncX2WTA1zq+o5cc39twUpHkKh5ds6OEnoe3VA23nygzIoLbKKXjm83W/h7TIBYgmFQA6Xe3yTDG3d7JemC546Ojo6Ojo6Oq8VVMj+L4PlOxikmoPIwyxBhoPQZ+OBDbTlLXLUvFSw6lqWEzoMBiqHvEPxNjTBHY4U05N3Soc8+IaJapcz86Jf7k7Vzp0HHqUWTDXExh69WyQ7d/vU3hJ0DrGJYni2m4WHZCJebw9fBOG1Zuix03rSKq40bfVoxCTj6OpCDEs+xrIjXL20kOTQyhxk5f8xEq8BAmJXv9r6S5DCZoJeZQLcNs4R0Xr68hezjOrbU9y7JYSlgGhkgMD4Plnju6LbRJHScXJSCBjxbtNdipeeUQsDzZeJYSiMxU+i7iEiGwBlMD5IdVqHvrr8mdI7bDFocuJTTcEwZKJFJnxOuX+PD1vn1xM+BMYE+lB7jBqyQdmol9F2kFjFzMIQPi8cyY4hZQN1gfiYNBcc7jEtwnAMLcQfTN9aM35O+Z590/EB0DwG09xjcUehIfbzI2ntYyhXt052mREFx7ZAcF0WyJbU9Fq8bnfnp6Ojo6Ojo6LhOXCfzk5O8P93bfCtk0VKHk+5iHOq2D0RLDBTCaMdxn834mq+Ln6q1KjHEV0TkoNbqYdRU7GyR7p3mgNaVqXYVTJDrH2yXigEi6iEU34xSD2fJIdTdaQGQ/Azp9GHhzm0GSERkfFRrFwnmSPMDcGHNFkqlknPanm/i2PDJK8tlgy3i4rFBu0RjECxXJRwqG1YyIdL65AbzU5W1oJD3luYnGZMYGSDch3F8Oau0WMFr2p/6uTDmhy1hPXdY7odcrHPoNdZ0ep5x5WeIy8lMpv3xmh9oYZSRgZ5H63Pimnq2yELaoQdCEVQwQlYctQymjHD4HWhB7cvc0v7goHTOPNbdvIW/o7iq0Sw6j7IarhTDcIhsOXRIuGWD0/yUoqc4oC7HtpbKwp2zDkxcuqeV8imXwO6tOw8rNKp9gB7oQGPRr/s4bEwp05J0/9ErsPRl+X2fjmfVWn7VuIg1/4TjKj9+Ojo6Ojo6Oi7AjeX5ucqPnykP8t3Dg817rUCJHCHWBoyPTt4bHoVxZ9bqSmJEB06ECIthD8aHor9EnA5Ion8ayQ+9pgj6n3kfrdZW4UYAPvQ1BqgUPnQnRkkG1xKbiZRoD0tkBk0DdClWbLVsM91HHQTXjmgaeRzNYuyELt4yzs54uKvtrS+ZFnur/tRe6wgxY3oswJD0O6FTkUGyJsbu0DUO6+LU2KEGo1W1YdJAp8djufF3u7aNyskO/bJS5iKWteBkhyJFJ2c6PbK+oe+ZxLFFMxLWnbbU8ZzVSUVTmJ99eYu9Miam/UGbZT2XvRBxjM8R2h7SDek1DRF7YJ0scnKIy8dybRM/F+sn7A6A/caoMtEozDSClSqdskKmxv5C+1O/EyqN4Io+yLPNeA/inXNUDdAloUgVI76DvqZ0DloxMO1g3u9N++MS0K4kQgQwNluanzVgm73TE1nh1eGyZI4dl+MqP346Ojo6Ojo6LsMtfW9d5cfPlJN893hvX87e+luL/rA2lajCsUDESpT9nx4xaINIHli2d87CRb4JaBhmSsnurQrkB2LtT0vrw1hjgOaSV97aVmUmLEqDpuJyecAiNO0PND8SpiJFB5FpJHIJhvhQRgYpIUrNbgydWLX9CawwPZXm5yL/eEPzM2RaBl0PxkrLgo/sUF5hd5rrWOvTuLZrp7RRm/YkvAWL58x0c/Q8GAOUvOZnGUitKBy/3EfNcGmBLSuaWQIwPQdlVw97zbflmR+wQ8i/YxogMJl6n5zmx/L86HSyEhkS2nr2bmDmB1FfO3pARMp4AlDe4oycQFZmBuwjmCDk7nJaSDznwxCfd3vufQSoagLHMxmgZQfQ4jADtH4ajJmY7/KedPcQDLWuO1geHsy74r32vo6RhazjafeFcwHp862n7P8/eZ3RcxQ2vSVc5cdPR0dHR0dHx4W4oe+tq/z4mWWQD493LiOoY36qIowrmgD/sa5NLAJFrdI7M63rzdeKoU6IdslQ/JcvfcsGbbkeSAPkIhWKtRr1CWzxeLAOqGKAkPHU6ZBc2t5lYoyPHs9neNbLUxgg3YVlfNZ5p48ZKLt0pfVpsDBV9JVNiRHy2oZTD3VgiVYYny0d0ikYheWOA8Mc1wD6jRwZoGVzWKkS1yFSxUg7r6+hdbgfYBosksdfW21s1xjzYJrq8XUO21j6BEYG2Z9j9FeJAivPBVggY35WnllvYZumKJ/W/pjlb1Fj0AAp44PMz475mUzjQ8zPXWRv3GmYtmeioqhYjn0NjuKYEfmlbUcwpRbtVR5ARObZ+B3oXm4B27D2Bw/8XHYCvdGASDHtb36C9sftlhieku8nrpeNCFMxBgjvUM0DdAYTVEXyuXcg7jPWHYaYV+2+EXEIdvxJ03Pf6Q226MILXg7GYLrzwLt+L8dnYX5uye3V8/x0dHR0dHR03BSukvnp6Ojo6OjouBA3xPxc5cfPPCd5/8klOdy4o5wyv0mPswgaTLA6jBCqGFxdK5TsmvBZpAiarQSG0r2TFeBztPsYXWJraItmI+D+mpludn0odS2jQHFo0NVMaSNcdjBhZ7lOqEU5OPfAsrG2xdQni4MbR6n6ShTNYfOuTyVLQQr7aKISOGMfZ7whqE0yYWctJge9XULedRehcqpeb9O5kp+wcZ9nSwpHwlq6FmHIs3DaCmpGl6JP3seu1nMSUgJFmByTHvpAgDsKMy5iU3WDQMzqngULSc6nfT5VQWF1JSGdRBFAl/1z+PvhiHkNEUfY+p07DhdBXRU+u2eaUhnMM8YBqok6l6UusxFxtBVV2wprwmds40PdVfw8k9t0gKjZPcsc/g7hs1jZDPTRdYVc4ZO9e/Qaazt+ZcTTie4unl9OCSlElv2+My5vQisz5FyK9yiiq64xlGGBG8zcuY1iqGuJEefGu9v+/6TulHnduMqPn46Ojo6Ojo4LkG9L83OVHz9zTvLhYV8VWPS4RJhmWGGAzBr2yRRzDOWtdkXCZ5Fi4cIiQOj7JDXLU6zd7ZBODrP04DTxkPdNbjmXtTABJDFAvk1hfHQTmGg6TY79GFQkCQaIL5eriek6rk0z2Bt0HEwH2jmWhZIMcmmHLQYokfXNrE5oyyJs60C9rfUO1w3DB9fN67V32A3uO7E6Jngu2wxmMes6Zn4gYvcFO8cc2nKyw6bgWXBqNBY3GCCwsQNZxy1r2UpeUII5PNcPFkxQi1mNWeIw8EY/EVBgwmdlBp7mOC8icphi+DsLoOc7XEc31sH4rAqfU1i/bKPL9kiACMEzXjou4oDGKcYMSlh40XIFzl3AAuhAB8/h3FDiBs/9cHClKsAQg6hCMMShHfoeDmVkIwYwAjK0G667eMVUtVztsavfl8d8CMuOzPw5ph2Mz72lIVFhMoW+ezF+8QiQd2EzC2vHm8JVfvx0dHR0dHR0XIjO/HyyMeckLw47F+pe39E65P38u46vfiQ6m1Wjs/de6NP1DkQkMkOTpTuPybYeTPvjdBYIY0ZJDDA947q1fQpggg6pDIujsQfR6sL84Ikl0vpY0dMUl/tLg7Tzxg7l2AYMhz+MMRmkWbFqF2Q5Lj+jRuYcrQ8zPqZVamyauMQGl8Koipa6ddR/swtD+L2u01szcNJEo7L8GNG20GKAaQBZYJSNO46dM3WGBGGz079Ywj88D8xGnqEBsmSEekCv+TH2Rqf30PzQM7sPv5f9HMAkbWh/mFkF83PcxbIXT3N5LpD4EOHvRy1weuDQd6dvmziZIYW4T5YYsfStSnyI646crFss5GF5qCpV2CXJOXnsixRKErogaH4Qhn8sbQdmfji5IYopp/Cw6jS2KWM7an/4t4gIihNl0vz4Mz+Z3HIsrBr0QCUdgmp/qByL13BCL8npHDDOttIvjGepNV8tktyW26urqjo6Ojo6OjpuClfJ/OSc5PFxv9lmKwJsDZX2ZuD1qWp7J9uJ2TxgCUDzA5/xBM3DsG4LrKVXn86wuouhRZoZETMFkFiMI8JmH72Ey4IdUGRHMs1P2b3pIEiXwr2b3f3i6hUWzGLsStT3LDtYud8Njc4a42Oanxa7Y9qYGLlVdtqii7S/ukMrG4DDuKcTluuACDeSYiQ+d7+MGDdMoafyeouKTaNAodwICbR1XAi0MfZOsUCc/FDEW9VghaALwjb1fiZiiwxaqsY/L5ZMFMnuxmjdP6og7XFXbsiTam0Od8oW6PxsTFAsfCoiVvpiotIXiPpiJkikVKhIxPiUMVrYiaI3ww3RdUhQKbRepIpczBvvmNKI9EBIuHiMUV9LF3QZEiESO8wFk5fGtIwLLhMDJOIISutiHKc8RkVc4kMksdzFcRsiw6AHGhem5xQD5H/zdLD3vD7TLU3qWZkpXwMuKtnzyUZnfjo6Ojo6OjpuClfK/IgcD2NhMty6rQiwlwXy8DSZGV2EqK5ziqAWH/GyDdgjnxdittw/MSLM1kttvaweb6NPuF5IJQ8NUCEAvK8+6oLYZ2+BYo0KhWbJgjQii9BHJLEOyIwVYikC2Wb9pHNtLebILNp/Uqs+tbZhVgjsUWNscD4Ti7DSqCMfWWUaiR3u67Icdn/piz9OtJQzMXF2bX2kjWlMov4kkZbFS2gs4gkMAxggRMa0chDRsm39A+53vO97PFuNbco6MEAxUmxOT6Wx7hC6uUN+sUyVAXqxU+bHaX44Auxprzo9jvry4xalL4jxgRYID47X10DrgwgwK5uhfW1p1yqLFsVWMT+Va53X3oNr0V/+t5WrgQApjh1/LoNFdwlNiSV2x+ZXV66eYd8ArE1sAfkRa3+Wtsr86L0rJTCg9SrlLcD8gB28HxARFiPD7l3CMvzG/4dJz3VP73Mf2ejz/DwHB3NLmp+r/Pjp6Ojo6OjouABZbiraq7u9Ojo6Ojo6Om4K18n8zEnmx1EOp1u+GnDyQ7cMFCpCdycKe9zCWPwsy7ae5qVjTlSRfR6U1j/jDrMY27sE2VsDN5i5v5zPhBMi5oHC4UGLp5rirsTQJNINlxaugEYYvEhxF3lKP9EPE0c3BMqmqyQXlrm7pto8qtxcllgunk/Yhl0LoPuRpC4kH4yuMJQMsSSWLT21qaHhcojuL5QmcVn87T5wUjobR+b+8qHuJBwlwbNPDmhBA5R8k922recDoe3o0h1KC9g2Tihs7i49H7tA+hy6it041GFYTpLFrI/j4vaCIFZE5NESIOo2E4e+6/1y16mEv8NNuCxfC31f2kBsj50kmvc+S0zjc1E5i/y4m+g6m8+Y/M8NcHqHbH11748jxh6Ez7HblahZXPg7jl1p7EkA7ZdRGHw2VxaJvqWM2yJ8jm1jKYw4to/jUacxvUNMRItEiMuNfYAUIcXlg3+v2816HsHzLeVf7MxPR0dHR0dHx03hOpmfLCJPRUbm9bXFoPn4zs1SQkK/+gdveS6/UaQUMOu0YVCxtQtGZlILdC8++5kex8TWZ/BcZ95tLnvhl3E4/CEVqx6FNLMuy9r/mYS93tLJJqTECWFfaKtr3aUpTAUxQHNYHK5xpstubI4pJH3YOomXjxTyPq0LnsEOlTBgiSeyBRucKmZ1Mdx5F9kgMEHF2o+Wr0cmASmHFw+O+bF6mVR6gctdBFaKmB5Yzq2QdysmSWVYtmBFI4npK8Lneh8DPWDlbWBSWFuH8fkuhcE/6cX4lDI/vtAlwt/B/Dxp0kMLgUfou0sGOSHcm8TQVgrDGDTXcTB9FOI+NZiftVeatUjxGVsWYTA0ngO/TSs7J2DsKcaoY1zB+GCZMT9RdO87zwzlmgA69gn9j4Mb7xecViBtOQEiTSd/73YxkaeVwFgRQosUMfS9vpvxHOw1xwSKpe6rqs6LGDq3/km8btyQ5uc6P346Ojo6Ojo6LkKP9npDSCn9lIj8ZVlc87+Yc/45Wp90/RdE5EMR+bM55187ueOcZHgcjD2YzviAvqi8BSyD8QxznprcScziF0LfV/pZtD/FJITFjKKPfDzWUATQXbfSApZgrhE+SykCPkKXXduD7hjXGyUwrDRGoytgJRAKi5BwEFkmF3Ebc43KRIwPMxxLP3HA9rSl+RnA8KxMY1kIXVcxPzG8fKsoqmmilPnJPnufMglgg8AEFUu9caF4DLAmSxPQ+eOAyLOpFcdMYerD8E07QdqfEjLsk9GlsG6y+fOt3NGmKUyHrTFvu2cGSGRCSDLOQ8frYdSQd8026csdcPg7GJ/HPTRAen5eG3UHHVBkfqzAKbE7/ndJdoibqPNh3K6zfyLuGfCrcQAOX18pz9IE6938GOdEoXi+TQvEnWs8v3aSW30gBstYYawf4mpxOiBcgpVx7H9zyZM1LdDyu60HupcYAv/gzsvrgZ6F+bkhPNvHT0ppFJFfEJGfFJFvicg3Ukpfyzn/hmv2J0Tkc/r34yLyV3Ta0dHR0dHR8aqQpXZ7XjGek/n5MRH5Zs75N0VEUkq/IiJfFBH/8fNFEfnlvDhu/05K6ftSSj+Yc/7tzT1nkeEpFebBrTqv3OjLoVXeAoms6qSE6v91PSqW7OkBWApAtvvwYPqRjZ3g7h9pvnU80vzwco+ktIHVKjVjlZIguu6V5GNqEZKF6yOSqiSGlRagcQIcwUXlG0JiOVip0PocifkBuxOS0emyI0RL2paZnq2XC7QNajEmF0mXRi2XgIgjZRrqwpY+5V9JmCZSWDQrawENyuj0FkiiiLx1VKCTkx2KlHvFyQ7ZKhZxRSMpoqqUmAAj5BgTFmwRwPiMjj0Y7PmDWQ8qS+JyKZqxgwqbDqrTezctGqAXw8LyvMilZM6nlRVCyYuPJmiAtACmRX+5SDeNADsQe7YW/SXi2SCwEVHbEh8/DjuND//QiuAyHRCebx1nRhO2wgihwTlDD8SlYih5ZqvECoY9xqeryBqP2wRRvMTueI0Z9EAzMz9WvNcxMhYBtsxPpAtiLdDyu60HQvHpuRFJBzZIhpv6DnkWPGe01w+JyG+5+W/pskvbdHR0dHR0dHxMpPzq/95WPCfzs2GfX9RmaZjSl0TkSyIiu+/9jCxpbqKv1W/8OhmgFiaKSLGU506zAz3QpObSuJF0wdZlzKewj63cQxV0FAyN/DWnsBUZBkw6X9QWXgPC20L/EH33Pu0FjDdUG0BKI5tvBL7BIBs+UOsRWpaP6tw90MSMHyw7nt7Zh3mwO/m+PD7pO3+4LPueT+l+H7VT2vF3Hpb13/n9ss0f+cyy7MOPwrnLp99bpq6t/BPfv/T//YVxyKotmb73HRER2f1/qk/57Lu2ye4FogWX+elOowb1Gsy7qK8ScQTJnU4PNG+lGNw2lL8GbMeWNVwiwZbpE3Q1Og1RYGdKH7zmZzAdkDJkVqG1tc/lZMD8IAfQk277blYGyDE/j8oGvauD773d0gZlL6D9eXLMD3IAIQJsIhataIC8ZsaoqqXbVt5iaRO0UqZjW2GAQKS4+12K34IaiUyQtXM0RD7FKoeIKmKqSPtjeaX8/TDNz5rW57QGKJEGqDBZpQ2YHS56Cpbn6EthkB7I7iEXQ/X6HdrfPB5pvo6G9GxQqyxMx6vDczI/3xKRH3HzPywi336JNiIiknP+as758znnz4/vvfdKO9rR0dHR0XH1yK/h7y3FczI/3xCRz6WUflRE/pGI/LSI/Blq8zUR+bLqgX5cRP7gpN5HxDQ/ltnW+YiR0wPsw2Sij1Z5xPPQygjKeUx8JmSR7bw8puMxVufybNAGf9hLcs7w/inai6ceFRtE19afjemAqGu4lukQLUcPRH2hYCQixsAAJbcN9jPdK6vzBAZIdWFz6dVw0PxEypQMT8s1RVSUWd8fueKY9ws1kj58ETsJy/ODD5f5u8Ie5Pc/jE31HqbvfrAs2Je28r7G1yHCUKthju8vjEPWPDO7D8q4Qv+Tanp2ynJN96qReayfD+icRq1ka7mAkJ9lBwveaX6mOP7B+Bx06tmPx0kZnl1b+8M5dvyyMj3/jTqCcVXqcI+x6BkAnbnXcTspe4ocXe8NNfMDHRByAFnUl2qAXqgG6LCvmR8ufsrMT3Jak/X8Pg1NTvWMRAYIGsjBa6MOsSipHWCIGp3NyKMNDY69CizfTmR2mxFuIK85RJPyhYVoTtL0gaAc7T3ceI8Y44N3gPbRrrUb47j1pgeKeh7ODi1S57kqjE9kidYYnjcd7ZXk7XZTvWo828dPzvmYUvqyiPyqLCrNX8o5/3pK6Wd1/VdE5OuyhLl/U5ZQ9z/3XP3t6Ojo6OjouA48a56fnPPXZfnA8cu+4n5nEfkLb7pfHR0dHR0dN4WcbyrE7CozPKcsMjrBc3Q8xhDYLKDbaR+viv9bucJW4NFRxg8DaNCYwBCus6GSBxcM5DiCaNqHws8mun7yTc/CcIGr4aKEkTo19xdtCuFq9oJLCBUpitmWI/FIFGcAACAASURBVHmfv6kpTjerKvALwEKQwalTskO3TI7qIoMbTd1Ttr5VsJASIWYNiRafRBOh8yh9oTy8ST8bLy1LI0DnbIJuK15atoFrb4YQFZH7RxKsttxex+jWMeGzu9gmeNaDvphj6QgTPLvBaS6xrM+DuafWx9lIAQbs/vJh8Xs91h1C3vVaPlDoO9xfIi78HW4v7TfKXiAEPrj81AVm7pBpJfTdpQYoqRmiC6S4k7wbZ0UIbO4u8g25pnZfMQ8XkG3i/YT8kPLxZBVWAoME3Eh6GLa3KTrF51e2KV2KB7flVCZk+a1ucrvuaIuEiG7/J8Lh5z22KceHGJrLZZgUYuNCTTIEGUXHq8dVfvx0dHR0dHR0XIau+fmkIy/hz+XDOcWVInKKAXJS1gprX+RbX+oohcFp/VuC5Nklulra1izOFgskUofCiywp0yMuZ4BeFyoGSGyBNqhFoJa4DKM45vQL7A6zHkUMj52u9y2VWNg4dSJp/DbGRxkgMD/ZSmL4GPGVNw0VOBURSTp+sjIKVWHWVtkMlA5B+YoDmDGdqmU6uMq/EDbDEgfjg+WzDqGgqyfhLkK4j8ZsFPbjSQXPp4TPXvCMkPODJUJEwVztk954b0lDeLrGAM3u2oPV3AtE0cr8CJif5WR9keL3NMT9RV6mCH3/cKGc5R1lhN7dlTcJEjwaI3anTJYlyosM0LIOP3RqLALmy7nV/7hWBMKeLAIrxOUncACwnf7dxvVlmJBplKqoembPdUt0H3czVAdoMFyZfpTMqgtUsBweORaNc0JEL3gm5sfu1T7FXbj9T2uCZxI6+/8bGLfvjOl5Qt1v6OPnLfi319HR0dHR0dHx5nCVzI9pflofzpaa/TwGSKSwQPkE49P6Up8b6frXgDD82ap6LpOasSmaBYS4gwli7U9wcguFwRtOM0DnaH4u0foYYG1ltqCi1eRPw6x2LvKIVPk7MBHOZ69sx6BWqzFAmPe3rlUGoNl3r4OgQqawaM2CBgPk7uVakVMkgPMWtLIGRYuhlvrKtiJOH4QQfVj1xiLp2NmVfhTmp834YPnsToNLXsxU1uLRMz9IAjjvaKoaGkyHElZeWCFlkvS6TEj0h7QL7n7sTty6wV05PEtr2p+9Dj7P/DwguaEmQHy00Hed38XzW85d2S5ofxD6vodObGk3Oz2V6X9Ms6L3vfkuamkcxb3zqJmIDLg1SGZpujNdj+fDM0xTg+ZYw1oYPDE/yT9ruszGHhIxTnjXtc6TjpPjdcKrafb/8UhPZRogZo2kFQ4fNT9Fh+gSuO6nsBtOoshTxnNofm7J7dWZn46Ojo6Ojo6bwv/f3tfG3PJddf3WnOc89/5bMEgKfQFUTAhRSVSCBSQxoLw0/VIg1hATIYaEaEKi32wk8S0mVj+YIMFIU0lqIiJRKw2Ul6ISPhixhZSXUtBKCNR/QwMitv333uc5Z7YfZq2111p7z5w59z73Pvecs37JvXNmz549e8/eM8/s3/6ttc6S+cEIDDfluAB4wrpoqnWy9jRVWf9VHx0hVs1P+42qgRlVfNH/jnVBUpU+eXIGSKvWmSIcYxE2Bw2GqkwQ6znMjEqnsOooTTQsnKzWWfUU+S3BPcXX3SAz3U3tD9XESOBPmXkujic5OQqS4rY0eUtkgKTedjYcRAU6ekPdyFiIFZ46K/OzEX0H31Nld+r1hcVR3YVs9/39KU3uF7N1vN3t2uCeonu50S07BRx98FCx+vJpHDS0iBaHtx0LrlG1H3JsflDXsBh73p/Acg5s1UFibbSwQLJ9wGztq5j5eWWcnF5K2AvbZtFA3V57R3nKIuxNH6p1EVv7qXit1eLURzIyPUFX4zQ/cMdkLAhrWNS6z4wRYSa1DvD7HSy/i017TKWaIKiqlYkMkKlEYG1KvE9j55RgOVcdL1prL+mjcC7Eomvat1FBKuMjW+lLf7yHsZAP8fI8UDDPRp8hzvPjJ5FIJBKJxHG4nG+f8/z4ocLBGmcD47nc3VQ/Qej7AnoS6DpvR5QQWSJhfFRrYCYC28DexGCom84oriEwVjJA4ZoAMMRQG08wgizrVRrNj+g55Ph07y07MupM08/GtGrK/JhZcfBPo/5rrtq8VQ/EW2GChEHRGbQVT4jGqj9bqzPQev+0TcoWlXCOqb+kaYIwNDLLF31PHaUyewczMMShDOiK9WESlNEMA5n5D3qfRH/htT49Vk3ue7T6utkZXzc7ZnquJjbn0xwG4vHGa38+NT7Qc2T8P6JpXG6Z+VFGhu/Kxjy1Az+tNZzF4RAxwgBJOfK8XIvmx7wBHnGjG+3P6AOeSrsA4IYtv8Ti7SZaf2194FPAPAf8DhrHwGQs6V5iuvr5MWNdiJJgITko0yd6McP8BEaGArOxigEKVXVEcmBk6ljzjM9o2q4sUPB1VFmdeaYsvj8qu9PqkOI7Z4w+x4yFWLTem3vXubBIJi39/DxbnOXHTyKRSCQSieNwSYLn8/z4KcDmtqBv/XBgdtRB1AHNMUBLX+prvuKjtqfxBWQnr4FgiMFQcScMENDogOSUXhzYMJp6fizWonBhYwkzOECnp9EHCgXrr8EaVomWhafSkcGwvm7GrbAeXoNB0fOy1WhFCzHyeYQtcr0xw/jU48VkFa2PVJzLGwLj46zJuJN2UxqxlRHtxEN1q+cYwn3Smb/cp979Un1Q0P4Io2FMw8TjsWzF38+n95NG5pXNtH3VWNlH9abM2p8tmyI9ZAboEd+TDVkmQEz/ps32iMDFg1pS8jgowizV8q+FWQpBUKP2x1l7CbslDNAo98f7+enp2/SYMhmDO+4a22h+ZCvpZlzpOArFqrUXn2qGlQa5Ve0P5wmBQS3iH9XFV0Jgaar2xjNYVvOjOqDAPkXLLctKjeG9EfU8Xi4Z8jRskbcGs82QOoi0rucTqId78fNzQTjPj59EIpFIJBLH4YJie6WpeyKRSCQSiYvCWTI/KngOZuwTZqjhevZsubPLXwscbk/U1tuf0sKyVzR17ITC0MCofGxu+QvoL4EBJhTGM1r+WsI4I/6T/dtA+0/HmNIOwQrHYPo+Vj95LU2tyzto8or4ebyS5S6m86+8gBhX9d6SBB4VZ4a8NKYOAHuC6OrjXxqGWRR/w9XJYTChLyZ8BoWQGrrctS1ufzCm1WMQOA9heTBubR4xna+BTpnu39XyRfwsy12PdNnLOzkUU3EAeMDiZ1lakoC/j3htRlwsbNx6i3Q47/Kt3axxVxAgXhCsE0Uxe78mEVZLCAwWOhMLnzd1+U4dOvLy182VX/a61eUvY+oeNPFtNBP7UgiNDeLiKm624n4eNyH8y9zyl83TCJ/j8rNrAI4GhVN1yUqWXDdmmVPeGxLUOISu0OfeiZglS99M3jlWjctd+r7y5zpLcXGMqEtwcr3QIfaUjvj5eSI1P4lEIpFIJC4HBWnqfvJYIXguHbNPh4UZ4rwA+uluZyNwWyhOgp8+kClZCIpaHbIZNieYqQsTJMLqp2KApgJ9HTeHZy5zs5vWQZg9x6eNQfxZg322szxlhUTQq4EKa9bq7l7YGylvaqAGgTQsS1HmJ3ha3HhWx56j5+oAWsEANSfPUAOmDjEkBrT+wgCZGfQuHBOTdwl0Ghgz+7s6SPTC570xdb8NgudHezF5n8bcp3bM8gw1lMSDURwKTmyQsCwyxjfK/FjBs2fVRBC7XbHSP87cfyuolmtJXa7V9J3N8gsLoI2zxpc2U9rjUbbTsd3VtB8F0FNdRGg7Yybds4Kg+N6bEUDbvNK2OQaoI3hWpidE0nFjIwzphllYMdYpDHE1u9/bhgRRdBEBtOz69KluXsSsY74T2FTbJq8CZYI4TyfsSPRrqiEyhMDqtNUyP3PhlBJ3g/P8+EkkEolEIrEaBCz6aDo3nOXHD5WC4baYJVXzNd18TJcmT6fEbuqyCfx0aw+aM674uu8xKHOhMJyDwgbTXEMYHmV8Zhmgek6LTigMmenEui0gtn9O+2NnQUJUNPoHYWp6wVCD5mcfZ2wmrxyrWgbPAJWtmLVWRkP1NVcShmDLecREndtnNTnw0Nl1zVwPLoRnsO3o2hkr4yNMkGeEyLFF00a1Pg0D5I8DVdsjjJiGBVHHiKbvmAV6dDs9H9eb6T49ZFZEtD+v7I3mR/U0U57roLfZDDLwTJu5ouoIj6uwL8IW1TqJabswPvsjgvjOMkBB+wMAj5kVemkjYTo4DAhvhfnpzf6jTkTDwBi6tbJAMlaCaXgnUO+c1meIWiDj2kB0Qjomgjm8k7VFNihqgELg0VWQMt1JnvUt4ttS71/Q9cA8m1EPtBDeQm+x1RPCsEidMB1VD8R1k8PBBN6mTfmReIY4y4+fRCKRSCQSR+KwI/SzwXl+/JRJx9B1a9Zhg/Skbnr3ZIfIAAE1WGGnatM5HUYD4dgSZh1grWCAZIYos9U5BghYowOaD4UxV1e7H+/TOMP8WMRAgbodvVWFj1HiZ5jVMVtbj8oOyb4wQTIN3rgyXKWkDNmG2Z5TKew48Kh6DvRsVNeKUMJYxBAbnVm90l1ROFXFUq6OQMsGqdNDncH7Wf/0m7cx6KlYfd2a54KZsaj9eWU3MT3X7Kzxwb6OtwfM4mzHh9NWLawkEKnQVbbtvOVz98rQTMnXriP6U+ybFUsAEu5FGCCpizBAD6x2ifU/D9hpo7RLtEA7CTdiWc5R0vpVdayBhMBQ6iJogHrhLaJzw8AECQNk720JLKAwQcIsuSCo8uzM6IL0+XDPkm9jfGZ7kT2UpdFj/GPjT7IWezVGrNcDNQ4NgUYPFNlltRbtvEe0XH23SZ+O3IxeEO0rFyrjeeGSlr3Sz08ikUgkEomLwnkyP5j0Bzr5W/zEu3sGaCqFLV2WLr0Swkoc5e58lQZo6v45BsilPYkl2B3Crn83+odZ6696jrI0wZeHWmt0rECqS/xQXsdyRWeUcbYnZUqyYWhIft8GFkfLWuCgh8FvKWztsRlo/e1MXZkd3g8aINVO2WCockytvSSP1/4AwJ4twG5Z8/OImZ4t63YebaaTPzlUzc/14K271MorBtm1CNaID/gcYUV8UMw+pPTbFZNhYYCE8Xmk+/VGPeTfov253XBgU2YSbyXQqdGSSbDTqBPpYaePBTMLYb+GsjBtl8ZH6y4lOX06YLRvMgR3Po+L+hK0Y8IkyuDTYKjWQizogvSNPGcxZvPqo1RcXn3szXKA6pyif7CODjCGyYgsUa9uNRRJ1ANJWdxP5i+H048+b2uvgoZ1O2ck85NIJBKJROKicLbMTyKRSCQSibUoh82Tzwjn+fHDgmeltW6XMguYAlWmsTWjnDunXteYKULEbH7563mNrSYiPIBxdsmqv/wF3JEjxCfAXNgLiyh83us5UkYbHXuerrYibM+V1yWzdjktYog/eKlBl7+s80yNsi3UP3PyEiJjYbBQFDrHrSm/cdjZCDHNKSoEL2HL6SJ2tabuGund71MIdwEAIy977VgEfXMly1/TBT7FUd2vzAUeDGKezuJiERmv4Oj3vAbzah6TtyL2Ny4BxHnhXHky9m/NeDpklKDLX2ZJQ54vEUHf8hrMjusoy2Cje4/Mj/85VDN43oqOviOOj0thag4vWcIyGNAuc+m5Idr7VF4QyMsmuJ5wg1BF9r4/4vJXD/VZlfZxOrVjvVniC0YQ3slheG+EkDr1JpvKVEsYn6cxu7fhTIyBxD18h1xSeItc9kokEolEInFROEvmh0rBsB8x8rfdcV94vRnWIRF0x/YSXsxWdH/CKjKKscr0XfIs9ejBG+EZIGCNI8S7YYBi/WOQ1x4aU3epc2CEpvK8uSka8aEteSaP4nDdhjlT9E09l27ZbJ1n/GL6jr04a+ypJ/UCvnwxgXdq05inZQMbRFF3mKGrMNWaAQdTd2UAhCG4rfdPAr6OYvLODNDjYer4KxY+bw3zI2kiKhYHg8OKqbG8AyTsSw2OWsvfYh2TZBnRvTIy/mYObupfGSAAuFUHiDve5+0w3QNxfriG+Vl8J/D92RMzSUpYDm4fMPp4FTyLmNmzO+SEwnzKzLnDzuYln5cCgyKBdL0TCMnMeVcwQDMsacNkWBZU2By5P+IZQtpqWaLA4sQgrqSMsr12YI7rkVBn29+C4fkLnoGLWvZK5ieRSCQSicQLByJ6KxF9iIhGIvqyhXxvIqJfI6KPENHb1pR9lswPAGBfdBY2mm+82a+95iP7CRig0pu9rGeA5sxYV633S0/uwn6vvIMMQD35kCPEJw6FMXdJIT+u4kzX6C1En7V8NedmQFmh6Po/igMcIqO3wACR30pxVYrDs0DDzJRB2JqdPyb25Pvagjn9T+PscGOm6MyylJ4eyMKa+Qcz4zmtw2DM48foCHHnZ/VkQyNIyAtmfHbshO6GNT9Xm2kgvGJM3SPTIxodSe+5gNjzeBFdzauGx9N1mMKwgVOvOeTFcIyWiK95iAEazPOxIZ+25X5/UMS8WZzftWP9GM2PDAUh3OT2j2Je7sYgubTK5shx3jXs3ZweaAhsEVAD4tbYw6IZ47Gh49cUL0yiDFsNpBr6pdNNDcvSMKY2s5RDLm8vQOt45RmeMbrI6Gl+4rFgDl91hq2+dIzJzwPFM1cvCH4ZwDcD+L65DES0AfC9AL4OwEcBvJ+I3lNK+ZWlgs/34yeRSCQSicTJopTyYSAYirR4I4CPlFJ+nfP+IIC3ALjAj58yzRIq91I/Z4UFWiN/qIg3/pAGyB5bxwABLW/SC+oZ0cwI1zBAwQHc8k047AhRsTYUxhoG6AkQ75N176jWdrqvU9FwBGj7VWanMmUjfypQhQOB+YE4iWMWZjBBakl+iw5oN9WSduI10AROnYtvoBlaTVHrCHFGA9Qrbq32B9bJIe+L5ofJlWLaHHVAI9f3ltjqS6yDzBR+o5ofz/gsoTIzg9vKdW5NSIEtfKDUygDNT4MjS6MM0MLglvdQdNYoLJSwPLuhZX7mwuW48snfu7jdDaIBqvdPGB9hZsYY4FRflOac4PhwaBggV6mpDuKJVBlSshvY528I705hXSIDZG/JHOPTpFuHnnIVuU/q9FArotiMvi5tmIugKTTXpqg3bDCjd7oPFubZaH5eQ0QfMPvvKKW84w7L/zwAv2X2Pwrgyw+ddJ4fP4lEIpFIJI7Ds1lq+51SypJe56cAvK5z6LtKKT+8ovwljcoszvfjp5gZgkmWWZAs+a8iP+rZ7UW66b1zAgNU2rrJ5H7Jfmo1IgNk0wQyU1sTwJEzL4UU2ESa4CmCod71yIw6oIYBos7sq9Ej+L4s9hxhfHjau5FZJM+SNx3dDW28tRfdcu2ueLs395p1QLMMkPprMXWSKbhqf3pt7eyb8qvWJ+ghrA8UsdgRVof7TrU/JvrmEBifwh2zZ2bh5uaKq2TOCdP3aOWl4V/MzFq0PuI7R/bF2uvVJs5BZWKm8XodAqYujfm9hig4zACpP6HiGSCxRNtqMFRrbTkde2mz9I6ZoMyP7k/bygBN+/apFBZILcKU6fGWYT3fPU0Q1KAFsnmq9aNngKh5tmoeWepQ3dDG53XWkDE0TERgYSz0qZaxLUVZCzfRuhW/3zCjzgLNa32idVefwKzvmEvxuVNK+dqnLOKjAL7A7H8+gJcPnXS+Hz+JRCKRSCRW40Sjur8fwBcR0RcC+N8AvgXAXz500tl+/FApKD0/ETv+RL9iDYZ4peXDOnvpTYbDerU5Mn9SrVE3r50hypE1HnRWw/bwIUswmaUtLjbrwjuAvvZnjgFSa5diGQ1/7QZPMUKpM3WaZYBcnwpD0tf4yEx3M7QzzzgbFuOreo6ZQYt+gH3eDDwmy46ZIGPtpcEkxRIsWrFooXYGLXUSBsizUbpdGrYHtD9T3cJW/P3Ivej4faFwX8Ybtv4KbAVgZSK+zdESylpJCdOjVl/sR+h289gdB6rmZpYBKusZoLl9C/VQXYL1V4f5OQYNQxY0P8I83Qy1HeJnSSRrY9AAaX85RjGkyXMgbKcxYa3+fXx51RoyaIGmCk/HRPMT3qvK/PWsazvFAWisGLvHIiFqnjGtd2B+KGiBXPnBuqv5G7D44NHzt/Z6AUFE3wTgewB8DoAfJaIPllK+gYjeAOCdpZQ3l1J2RPSdAH4CwAbA95dSPnSo7LP9+EkkEolEInEEXjDmp5TybgDv7qS/DODNZv+9AN57TNn58ZNIJBKJxKWj4H4szO4J5/vxM5bqttzwmUpB7gNFLKeF/WXMLX/1jvl06ro2H1wpd7r8Bax3hOgafzOTyS9/Ae0SWLP8xbvXnRWmY5a/hkDj13ReRljBFwszL/EWC9UL65IY+eWvEgXQZqmsLonFrQigRbRplr1Y/DuI+Jcd/REvzQ47syQqgucau8Nvl+6fLGXItWW/Y/JeZszfKVD41tQ9mr+L6buKRG14C1kiuZHlD6kLN0NEue3qXS0jCJxluzNLTfJ7x+4CHrPzxMesxn7VUMf1I16niaEvmuWvhbEe0XN+GM3VZRyLqb0In/dmLEp/7mcejE1HGF6dQXoXAY/FZYBZ9hLx8u1t4WvzcqGKmqV/jPm9PBZhaawGRzUVvI3LXJwen6Guutcvf2m/S92dc85QntRhP99PFJ4din1q26GvMrkvfuwN4VVn6xRD6rTODZefucSzwfl+/CQSiUQikVgFQjlVwfMT4SI+fmyHqlBtJm8kImzawhU6aYdE0PHrv6bNMkB3NS4PMEBuVquN9wzQJrA5ABoR9BDoiA3PoPc9blU0udxX0cGcnXXPMTvDgiA1imZlqwJoM81SFkJd/XsxM4JoE6i+BasDuDBDZH2tiz4hAmcO8VDEmZvMlg3zQ3tv6j7LAPUQzODF4WKNvTF/qqIx6TWzbqbPJGSBODIUUbNlboZwfzSih9RNTMbdtLf/mpJZtzACOxOxU47djmLqzsJnMYEfakc8GLYAgMcsfH7QCKDFCWLlYKsZ/NM/lJvg/HBvxcV8v204DneuZX70ueNtYH42LPq25wgjKXl2zD7K/dkvOkYU4bN38OmZRCk/qoljSzri5eaYtM+LjwFreh6EzZGsdcGC/Tmto0RbBW/8ID4yh/D4WTsEvWVzgueeyT7Mu+VyvkPuBRfx8ZNIJBKJROIAkvk5A8wEcZSv7LUMEHCMDuhpGKA27xwDBDwjR4hLoyEwQGL+e92p/sYKQtAznTd3kmfTOpOWQsLNdubxsbSgdYjpQGV6aogESZ+2N4YJEL1J1aF4XUo1263Xqqbsfl9YnSGYl/tjvGVqb1B3/qb+zAIp88O6IG3inOm7RZiZV60GzeZptGnB+eH0m7fB5F0ttu19YgJjiNoMudeSz4SfqKPpylZhkfkRU/Obq+mCogG6GacybjZ1sL/E+p/HzIwIE9SawNeaRPP0yFQeYxbflGlD8gYqWjQ6G31eWkY0hgG5Yo2Pbk3droQN4mOPb6/cdW55uzPPxxj1QDOhMaY6SD8H5qTR/tRzKnMUxnTQznizcq+9iWbm+vpwQz0wPkvm8Hqe5OHng/wYhHkuCt+I+gypCNVf11SqGIbpXjQ/F/Txs07Xm0gkEolEInEmOFvmpxAtireUAeqo9CPu1hJsjgGyx5YZIFuXu2SAhgXLCMXCzRAWaB8ORuZncLPVCXPaCZmRbgzzswks0dwsu2sJE7ZVA7TVvBoGgI+NPOstIVCotYwq0dGbOPGT8Ba6X+s0yjG5/2L9deWZoOmY1/pQDLCogRU7mgYtxO/2mJ/o+LAEBqgGcqzlRMeHOo4k3IVjAuSksA3Pha26sEB7ZXo4z+hn3Xsz694HNki2Nxtv/WV/v7Txmh9hgrZsvvbAaH62TGtV5sfrbKoDw+Nn0k5zpzQa74fnbuiO8eA8kdkuYXyuDIOlz1fQA93sQzBZo/nZBYuwudAYgCFC1FFoCftN8w0C0xPS7a0t8X0eGJ/loOBSiJzb6bOYpK9ozwQNZrAXvrgwlPLemg9wCri/Ac/b7Lzgokzdk/lJJBKJRCJxUThb5mfyDt5T0wcU/+MUGCB75GkYoCeZlQoqI7NQhup5wrm9m8zFCHtUQwDwEDWnWD8lQOvfR2ZYPQuxOeanN4O+Ud3DlC4+gWQ222VMxGot7AsTNBjNz0asvDhNZBWildkY/zhjYIOUFYrMj2HvqpYBHmHoOR1ST4PhMrdl1muL9Q+niyWRc9oTt/2xbqGxXCUosBga8mApTD3st0bzI4zPnq27ZLudtjdjfWJuxBfQ6DU/kQl6NFR28GHp64G2wSeQfT6iNWLUCXXzyQAShkk0P8UzTEB9ZmJdHivjI9ofo6Ka0QM1WqCh/qmoz4UwQFN6DI0BVDZIWNMh9HdjBdZD9IvTaGZqHjVCVctM3hfLQ2eOxdu5d/7YeZZigrwDtE61MNUBafBbf9keKil0P5qfNHV/xiCizwbwbwH8EQC/AeAvlVJ+r5PvNwB8ApPmcVdK+bLnV8tEIpFIJBLniPtift4G4D+VUt5ORG/j/b81k/drSim/c1TpNM3K1Run8/Dsv2xV+xPYlheZAZqu/Yy9QQc0vkNkhtNxiFT9/Ih/H9mfv6nK9AhrozNeX7a7tu4L4+M1QJbNkVmwzHqlPXWma2bQ/PvTvBV26Ja3e/HIbC1gNn7WK9tB/ft4LdD0WzQ+nFfYnds2r7BBco4yQaIBEtbFzmzVJ49ccGZW1/FU3bXCsac45idsRaaiuqH2UtKd6vl3KeijkFyi+ZHZfdD6FHOO6IGEAZI8ov25vao395Z1QDdX001+yIyPMEEvbYQRqk/XjpkMGTdzlmHOKzTffzkmz/CSf6ohUhlqERaOAxhGr89prL1GrttYX/vbyAopw9TXAk1pch3xCTTtqxbI+cwKHrwRnusZzZdD8XnknT26MeI1aWo0Wtxhb9kYdDuRhPTMUvy7wT/keesNdunf0f+lWCRXbCXug4S5IObnvjQ/bwHwLv79LgDfeE/1Bz5tsAAAIABJREFUSCQSiUQigTJ9/Nz1vxcU9/Xx89pSyscAgLefO5OvAPhJIvo5IvqOpQKJ6DuI6ANE9IGb20/dcXUTiUQikUicC57ZshcR/RSA13UOfdcRxXxVKeVlIvpcAO8jol8tpfxML2Mp5R0A3gEAf+AzP69gQ5HVnPIFhrma8PaXv4DOEphQ9cyC383yV+9YTG/zzi5/rTHtXIGDQUKP8AapYubO8lcTKDIsg9nj6hBRxZ99obNbcqCHnCZLY7Ik0FmekGtrWIBp++lBhNDe1BdozeGLBiuVJSwuy4qLgxizhKUxaxYvomEZczr2QmgJ73xQlsTkgisMADRra6Y+ncvHXVBJv+Q2cJ2kKl5YL+N2giw+6WpdWLawv6WLRl66KrKszdfb72vbxy1fZ+/F0DsWOt/ua8Me8xLYQxY+izn8w82UVwTRshwGAA8GvxR2aBnM/pagp3os3OTeMlgV9XP7QqgHoI7TLS9+3/LyloT0qAJoI3gOy8FXLOq+3kz7j3ZbzmcdI/olscdqIMDLYFT/rKg5fFgVGheFzgdeXtHZoUnTMT/6LLpUummfj1mnit1ryxos76seXdfVTGa5KNmsZnm2fVGqK4D7EDwXvNBMzV3jmX38lFK+du4YEf02Eb2+lPIxIno9gI/PlPEybz9ORO8G8EYA3Y+fRCKRSCQSiTW4L8HzewB8G4C38/aHYwYiejWAoZTyCf799QD+wZrCCxHGzaBCVezqrIVEmCazxsaJWyuSbkTQIRxEZICANSxQb3ZzSATd1ukgA7RUgxnHfxbR+eCiebwGJw2sDd+4fRQzm2O1/DDDVaFnL6ikN/fV64X96bcXm8aQGLadMrOtASFHt/20BsKs50gM0r2Yr8uMdyGI6BjE0KJD3QQmCABGNX+XcoXx4bx83DkflHATG39s1gS+B52x+zHphkF0JCcCUhFhm1OHMG6biWanbsJmjTJTF4GzMkF8inUwx9feCQMkAmhJN4Ln3X7PWzaHv/LbB8II7evr8uGVOEKcjr3E/fyAO0LEuI75Cc4SJWxKI4427RioHcuAMWs3aXsuT56ZLbM3GsyVB1iPjapbNm1ncXTPPH6zuwZQ3x8bNQzY8rbW6YaMah+GARK3BZiHMIpjHBNdQTKnjeSPSbo8U86LgFLR7lwhZqj3gMSVAwnmq402f2tUzR/f5975q4jCp6oYSuk+SJgLcnJ4Xx8/bwfwQ0T07QB+E8BbAYCI3gDgnaWUNwN4LYB3c+yUKwA/UEr58XuqbyKRSCQSZ4308/OMUUr5XQB/oZP+MoA38+9fB/Ann+gCBJQrMgyNmW6LNqJedNoEM0o7CGbN4GcYIOBZB0O1x9YxQEDLAtHTLCoLS0Etu6KzlxmfAI5NKpLmg6EKIhPUu46a1gtjI6EFbBiN6Po/MEHO8VvjFM5vIxMEAI83YvY73ZidaoCm43sNe2Fm9SH46ebG62yMXzmwhKSax7OTQ2F8IhMEGJNzYYCCSfrSDLrB0lCMwU6V5GrHcX2sit8PTuosgzVqmtd1CJsjDNNoHDyWLZcvrBHnKcoE1QvsRRekjhF9SAwxhX9wZUzduZMebiTvlGfHuqAdn3NtXgoP5Ank/nUBTFFZnLlnAVhmXgcpT0zpA/t0S54JsseElboag1n8INvKMUUT+k+zLkh6eWMYUXnH8NBWLWIdggsMkBYTdGA6blv9jo6VwAQJq2pdQcj9UXP4wf8NcC4aoiRRxUTh74U5qfjHwa1ATBBXB6ZOcoTKRX2I3AfO18NzIpFIJBKJ9bigD67z/PghYLyiyhYsCF/qDLPPADmEWavOQAIDBBgdEInS/wiBRXPBp2eAprQJzzsUxhIjE4OhaqDTFQ4R1YmiWr54LdC1mUFfjxKAkrUSMuNlMY1lsJrZsGzpgcvbc/z2iLcS2mG3EY3ONNseDYWlLE5kfESfYMJbNJZhwREiidbIdO4QAotGa6xqEWNm6ku6CqCJNABYFpW3oevs2B8RZtt6PWGCAjVg6qlkjTI+kh6YJwCjaKGuPONTLcNM8M2rkc+JlmHiEFGsp+o5wgbtroQBYqsvvtm3vH3J6rZUwzJt1WniAj1c2aCg+QljvYc9FywM07bD/AgbtC0+UOuW9TtiGWb7MDoIVQeoIYSMRQ0gPG3Fbi4yQDatDWsxrzvTvm8C8IZ3tRkjOgZFt7Mhn6fnEFFOmiPn3LMkf1MGtx8x9Pao3I/m54Jwnh8/iUQikUgk1qPAxTM7d5zlx0+h6ucDCF/WM7THHAPUv4CcI5Y9vG/LFh3QLc8qVeF/PwyQPSITG6nuk2h/YrgLwPjZaUJhiJUO5xva6w1BvyPN6fn5qddjvynBr0+19jKz1cHnUcanU360CLsKljBRAzS1yeuAYnDUHW/3RvNTmD3Ya2BTbxlm/fxoWgifoQFTw/HpN8+CVfsTGKaOb6ASxnZz25eGStT+BBYJqGNB9BU6GRYWSlgcYySk2p6gB1J9h2h+jOmTskLCAMn+1vsImo6xBZiUEyzDZF+swQBgv2WmJwRVFeansjyG0dgMTZqDdr/RqglTpu+PBT2QskFez7YX7Q/nuzbWRTel7wsoMkDWQuwqMKJVAzTt2/dJDSBs7dIqIgMEGB1QYHHGoKvphVhp/f34sWMie1Tm80rGBFtf9UJVSHndVmBxuYjUHFFEQJy30QABg3TNLe6B+SkXtex1Xx6eE4lEIpFIJO4FZ8n8gAjjVf8bXQMrrmSA4k+HGQYIaHVAd8sA9Y7F9JYtmrMEq35rF2qy4AtoNdQPUIfFiYFNl4qJvoGkDPH0LOyO9WfCM1thhUQPFGfJ9rdsPxG8Q1dGqNYjWoYJA/SIfa083kwz3xsTDFW93w5eDyQzz8FoflQKEDU/MZCqvXFivcKWYVqGFEt+ljz99seqKCPMfG0XerLO0Ee8sRYwqvUpPqv68uHjTvPD90m0St6oCeXK709pnHcbz12wDAtMT91ndueqXkCYnb1oijhdmKCo73Fpm1n+oIU8M3yFUdkdOWzYR/G7Ezydb8PlLPMUdUDCBM0xQICxAIvPEFuEDV3mZ/n9Yd9BkXQU7VJ1D90ylo0l2BjTAxtp8uj7W/pS85ryN6Lb4UM8UJsxb6GmW+FYtDi2DJA8k7S49vDskMxPIpFIJBKJxHniPJmfRCKRSCQSx+GCmJ/z/PgJgmeLGb97bRFuZWVG/BkO2+WDRgR9p8tftjLLy1/F0cxMu0vQyuBgbGn5a46u7qXPhctonCDCCIWZbm8DnMq5rdAz5t2oWS5xmcUc88LmGy5PhNBbo1avy2Y7bkcwjw+m8EArAtUydBmM74nxXFjF0NO+LIOJu/vRiqNVDCvrHbKU1V/+cudIHrGsJr9MO5gouFp8uLfq+h+Bsp8S/VbOKfEHzLMS6y1LGbxvg5QGc/5G6CzLYOZtFvPIqs3YEXnL+fLOKKEOew2nYc7RZbqwNMZC6NJZ9loLa74+8DtFx4+Gteg9D0H4r6FbgrDWrA/uJRQGj/Vt8ctgccxPeaYbHU3bo8m7xTHL5Y0ZvNhLjH75y5mON04yg+l7x3lmDfcSlnj5XBuGh3hQqGZZ1h1lSZn89RzC8pcIoNUE3ubdj7W4C/oQuQ+c58dPIpFIJBKJ9ShIU/dTR4GYNM7PutYyQEBHBH2AAZrO8bPrQwzQVIdDA+/JGSCfVxgfPwNRJ4idoo6Zv87N8pYcsum5YioeTND3zvw3lBsZID3HCC+FDQpm+GoWbxS2c0FQWxGoDYnRN4f/pAZ9lBl8K1CVAKm3TF0oA3RT2ySBEwsrmos6RpT2yJTU3Ajy2yh0lnbZXhEyS50RzjBAvVlp0+2dYaATZBHCBitgNYG3IRICE6NCZ2aCdN+wRWrqrmENJF3EzKZO23AOOztsQmRYlwB6H7xItid0joim7k3oFTtu+T5ssXPnShDTngF5ZHyi808LIciFARJDgBthfvhGXRtmdMs+BeJ1IvPUa2NE7z7JvRRSZQxbeamSHexB6HzI9L2XV5hDCZBbTAgUNX9ns3ji97c6RpTxtMTWHGCAAHhnkM/9O6R4S4MzRwqeE4lEIpFIXBTOkvnBAOyvCYdZkWfEAJljaxmgqQ7P2xFinwGad6PWomfWGrHkij8yMMPMzMM5U1ST9n3I469v76PMkEeE2XDPyaGY8pboCHFeSxHNf9UBXAgFcDVcN/VVtotnhGIOvzfsxxgCpY4U6JzA8thjMrkemjwtuxOfg4YB6gyzOcanBntEAwrDtDJbcrI5JtqeyPSEWXcx5uvynA1svjw27I65uLABwg4Ie6POFIXdMXUaIzsRzKQ7iCxHZHp6zE98dkT7I2NzNA4L98JYip4tnhvCwFjIuJV3kGiAbotnguz5c8GC17wTBL37VYq/3yXcW3FQ6VZootZHC/Ppzq1DcJYYQ7tYU3cK5u8jM0DiMkVN3+3AlrF8gExxQbSlnGNewneJC9IZJfOTSCQSiUTionCWzM8U3gKoFk+Hz7lvBggwOqAis+01eHYMENCyQHHGvhQaoxcCA4Bvc+OwcLRVqrNW5/Uu5A3X64XCiOc+xN7ldaEwoh4ohM9Qyxijg5ib/Yqup4YCsJqf0eURZkbu6Y2516LDqpEjBrcPTTe6FC03dhowlyD1j/7ZSpgVW8eF1fOb7IZ9200zXRPHkSEaqoWNMDwhtIc6frSaH9FviHWXWIRtg0bD1K/RA13741arI+yDWiJxJcVSaYH0WrCc9OMBMIF4B2+FdSvOCI14ScNZcH9qUNTA+FjGsn1WfBlV/7ZrzlEGdoERPYQlzU90NigsW+n0hwa7DSFW1OqvY+3V6oQopBvrPg5fosSkWH8J+yh1tvoduXhkgOLDZSHtoPL8WZiCFDwnEolEIpG4MFzQstd5fvxQ1PyEgwuIDJBNa0q6SwYIqDqgna/DfTFA0xGxLpE8dztkor5pE9gcnZl2boLmFTf7ShGEc811KrPBM1vaN3lFM/EI3veQ6C1iaAygtQjbRD2Hzorbc+KWwtb+FgZI+kN9AvHWkzyxf/t6CN8FzDqG8BN1Fs7pY3NKM9KaEAPwmov22qYsG6pi8MeKMj2cvPHp9nwJc6PMj+h6LEsUWAE18hH24Lq4fZdHLIM4i8qOjvDvE/vfjRHRjI3COkYGqDZadEAzLs4a7c90Ta9nq9cdXN696V15DmoA4T5DugY9f2zVYm7aj2Nwpxos2x9er1VPDroeFxKD2xYtGLt5eSvjSa/tX0zUeXc2DNACKDx3iWeH8/z4SSQSiUQicRwu6KsrBc+JRCKRSCQuCmfJ/BQCxmugpf0t1i1/AceEwjDXmfuAnln+Ajpm8He6/NU7JsnhPpm1k7pywebXeqQdOrNCzo4J7Nw5bciKIBaEod2DWHkflph6qJGv+Rymr7eGqh/0Pu1cXhU1j2KaboWjfXP4KsJulwTm7ksvAnbbc+wQkffsQomgXsn3r5pyN2XW5YFRlwt4mVBWH0XYaZ0QSpos18mBKIB2abJf/L4kWz21rCWJOXwQPneXvfi3Om1UZ4eyDGbqHwSvNeyBb58zk4ZfmhETccnRe06iV4K5iOdXQztGrtSVwtalW0ebYpZ+zdtqBh9dQvTGrR/jsvQr7XLLXkFAra4ansJBnhU+R0eRjak7b/fmnCo8V18QAFoBNKqnCTQOEOdM31HHhDhA1OdBnqLGn4R96iQifLOOh1nci6/BexBZ3yPO8uMnkUgkEonEESgwYrbzx3l+/BCwvw4JihmF5QLuPRhqYICW6mBqM3/RVaEwPGYZINr0sjssOTrbLBxz55oGb8IMU2eynFXZI+tvbKYjIhMEVHYjHlNTYQ3GuuC0MYg+owA6/rb7a4JAVnJujgFCI8YtylbwdYIzOaC++9Q6V/blFvcCqApTMldZO4OO415n3T5vV7fNB1V7vREmSGbdhplhYaqaw+/jtid49ixBI9juCWADhHloGSA0z8os89NhB6+CE80ofLa/Y2iKPXnTd4sm3EvDvPpgwd1zOiFi1kJZHlN+w/hwuoYQ6TidVIeU2lfCtsz0pUmjwHI2/R/OmzC49O77eMcOC6UIGXu9vInnjvP8+EkkEolEInEcctnrtFGdHPYwpwO6OwYIOMIM3s1EwrR3hgFaWwdfWO+id8gA3TFWmc1y4welJ2Qz7V+bmaiGgxB2qMP4KDqm8lO5C+eIU8YV4TnuEvVdxQyQ00GwXkcTZEthv5ZXgzvyVmbZolFTBsiKcvy4lVAVXQeYVZwx5dHZt39OaOElLOWLXbkGQ90YNicGP73ys3nbHa3mR+ro71OfkAv3EnLq0CRLcfIYPw7Mz6bH/GhYlOmFJtq0ygDVF11lg/x2K8FQNXRF+9RG/VyjxTPjqmrr+mbySy8lDcwq2pkNuXQA2G992k7OGeXc1tT9Npi/awiMBfau0fY0Ju/te3A+BAY7Qew95jLWxYlsuN7SWE88O5zlx08ikUgkEokjcUEfYuf58UPimGyJzXnBGCB3zM+KdapoJAPH1KGWP2cJtuY+eUQGCDjMAq3RshyTdxOsTZqbQbc1Mxe3ndE2+HLDMWGhghPF2jEdPGMnEmPQQwicDkLSxHllw/x4KxegN1OeNiqrClvA6oBaq6g51Fm3jHVhW9qGxKFAQfsjwSQt86NsEDM+wR9mn/mZvT8hHeiwQXFsyGVNSAwhrJTxYY2OWC8OLesi+jKxANsqE8SszlCfuq2GvJjGZdX+TOninNOFhTgQQbOGrDBaNXEQ2jwPHax8DvYm46jatCHsh+1ow41MaeoAUR1WyuDsWOzNMKANuzOldtKApTfxID93/qRqvcjj2NGDsfzEs8J5fvwkEolEIpE4AiVje506yiB+fo5hNk6AAXraOsRZanPBJ2eAplL6OqCl4KdPg1mtASdv7AK8GATp/RfTuiOYIB0qPfOmGRboCRigXrDHOd8ngmgJMyXyucFHSdXv8DlXli7yvk7EmEjDRATtD4DG745OqOUWL5ELDcMkFy5NnkYbEfxT2XAEyvxIe0Rbou2w1l6+zfK87SMzsOAzq3bHXDoA0QGJ9keOLYQ1iczPVWR+zM2NlmDboP1RJsjE1BH/PcrmzOjarLWkhogJOrql18c+hGHZD95/kLVaEx2Q3H/R/OxE86Nbq2/jc9UHkKTLc9OOkajtmWP+7LHo+wlhP9ih8v/+7az+sPZsDWZPkttwH8tPBShP4avp1JAenhOJRCKRSFwUzpL5aTU/8zO2Nr33xb3MAt1VMNQ4sy2hTnfFAJkrhv0nZ4DstSMDdLOqLneIToe8enjssjykQFc4r7f9WdcsA2TTnoABmrOAcXnmdA8y0+0Ee7xRy5dlSxir+dFj7AdHrb9EQ8F1sx6SNS16LaewXULweut0Q40vILGekfK5XfY6G/bzolZfEmSSXJ1thXu+YNqCETLNNE6JAaOnkjSpr7AfAwcpXaH5ueIGifXXA2MC+go7N2sZIA54ytofYYKADtPDN0EYme2CCKV6mfZ+tnq3RJif+ByI1mfctJofYYN247Qdt/2xD1jm09OP++ATyOqEKuMTWKHI7pg8dSzExjLL41huecACA6S3X+6JYZicN/M1D84d44KWvZL5SSQSiUQicVE4T+YnkUgkEonEcUhT9xMHFYwPRhibZ3Pw0NLO0y9/AU8WDFUocqVheRuXv/z5/Tp0rz87rueWv3rHJMuTL39Fse4zQ6dDDi5/AQdF0L3wGSKorMsGfZWv0v+YF39qXmv+q0tjfeFzz/GbiJ8bsbR0EC8BuHA+wbW/HpPAoLqMVE9pBM8iNg4CaMBolOcGowqg7VgPIuhm7PPz4rpQluDCUoa078ouacglfciC1uR5xfKXLvXJUpxd1/bLXSEGJ3ay7GWWQR4P041Xp5zkl8E+ZV4Essz1ythf/roOwufp2LQkJmOvCqB5+TPs27rEW9AsfwF6Px4O7HZCxhXfA6nLaJxmynLXfiMuGkRMzAJodXZoHCPqMlcI6RLCXhQjeK4+ZaOwHe3+3HJXCfsOsqTuC6aont7XPlQHoaU8/1WvUnBJsb1y2SuRSCQSicRF4UyZHwDb0UwQh3AQWM8A2by9Y308TTDUwwxQzbuWAQKcrm6uNny9tk7r7pPH8wqFEbHv3e1wQ3QmyrjuCDv3wcHbkjl8DJyq5ZGv00PjeHFUxsdvhYHomf/eFp9nJ6LQ4mfJQJ0FCwOkIQD2XvhcTJDP6hxOzHE9W8SaWWdSLyJNSSMNO8Hn9FiiQJREuDEoz0G1X+6fZLtQKCaZSTdMkHFGF4XOIWTBTA3dnp5Ccd/SXr4Oo7JDvM/C551pfHSAuN0wu7Ob2B0bCuPBhtPGGAJj2j5ilmcwQn0NShqEzz3nhnPQPJEBMmkjP/3Cdo6YngNhPa0zyNsyccTC7NxuZIyzMPzK7wMtE1odIUaxfyt0r84/p03j/BDouAmJ78XOezIM7iLC51sJ3cMHBjPO9nrD7gcXtOyVzE8ikUgkEomLwnkyP0PB8GCvH889x1PrGaBe3pi+UJVQh2McIc4xQPbYWgZosQ7Nx/4TmPb2julaPQf06zFAM/qfnqO/OazJq6yKzLKOCIWxpulx9ismrxvus+sw852qwLNf0f4oA+D1F1PeqN8R5sczPntjyrvf7vw5o5/9xqClUx5uquiFJBwEb8eg/XFpcrt4Qq4Mkx1wg94gPukJZpqz2h+XaSperi2VELP81pLeBiaQI4hHFJHZUR2Hp7Rs2wfVQvl7MAYN0G6ojAapDoi1PpvpZovm59FQX+GfVDYomLrzVoKkDubGSciLG8zkWWBzIisUn4GpPGasuPyoLXo4CKNpwuTMOECcc34ItOO/an88FbfrOkb0Wq+WAexh5u+HuSeF/PlV9uRN3MmEv1A90P4eND8AygVpfs7z4yeRSCQSicQRKBe17HWWHz9EwGZbp6bOmEV/nQADJFofawEgeYqfZa2KMRgmzE/lCDG6hO+eIwdlXd8zQIAJIXCHqM4CF1ooBJDoYpwVFjMykJnzVG+x6Ir6nh5EMzHyds9lbNGGFngYmB6ZGd6W+njesl7jllmBlzaTLiK6/t9d1XZomljNbIXp8eEBnEf7PedVJ4fF74smyLw5VBe04bzCrmxEq1PzNpZgEgw1Mig9BMaHgmPE/jnSRg4lUPz1gcXQlP6IGSOBUKiz/MAA9Yag6o70Hvj0YjQgwgKJBdgjtfry1l9AawE2hDyRCQKAa6b2hFmqGqDRXWfovVlmGCCL+BzItVXnxo23Gry9jn/W+GwecxmeCRLnh0BlfnTMy3MdhojTxOkxcnmj88MJhzSPvb8ffizIpbV7b93h6Zi8HigcSNw5zvLjJ5FIJBKJxBEouCgPz2f58UNU8ODBDo87x1od0D0zQAtjTUtbMSDVIizG3LSzlyuxNmgPHYZvGy36t/Dn1LyeAQKAopZOfSwG+Qw6mP328P2vgQ95K9Owzgy96hPYTwrCjLoX7DHcj6j9Gc050bJG9Be3vH3gZsO+zTr73Ux5bsd2NixWMdEHyv6K9UfB+gswFmAaAoOPifVX8AMEtLogEs2Pbo12STVDxR0rGiLDz5IBM8Jo5rlbZH643hqEk7cbq5+bcJAB6jEcFLde11OM/xpleG78OTHchWWlRi5HGKAbZmiGwWt/7G8ZczH4adT+AJXhkbTNbGBLE7YlCgxXvEcOaeBE/wZUDZywnI0VJKePxleTWEEK86nWXXxcGd5OeAv1kaXXWWBxZpnqztjshV2xWYUs3NnBLs9ZWZYcJZ4aZ/nxk0gkEolE4khcUFT3s/z4Gajg4bbOVJ6OAbJpEXMM0NI5/mpu0hpjYsqsVbQANo7mzGy3YYBsNXase9h4bckq7Y8sX+tlj7hPjXdUa6VR3JHIAPWaudZDtGWNRmVxDp+75/v9UM+Rmf80O+2xOMOM7iEyV3bWLRY1ogPSYJJMnTw0nXfLacIGvYp9nwgD9BIzQNYCRn7vVQ+0d3WK1l+AsdBirY/65VGWiNtlOkqDn+6lXNH6eO3PlMZ5N35MV72LMDO1HRQdscihOFgWGVJvYUVu1O9tjsoAzWg0AFSfLcJuRQ1QtOxyaXIhr5UqN3z/LFukXsCnPLdBA/SKKX8TmR8Kuh3V73TYx8jmDP64B7+owjnHaOH2warMauHq+PfWkNJ2GfvWn5f4/BkbP1jeGtIHQ/Xb6O7J+h6i5h028x4h23ei2fTH4lgZbFHCBu3nL5G4G5zlx08ikUgkEon1KIBze3HuyI+fRCKRSCQuHaXkstezBhG9FcDfA/DHALyxlPKBmXxvAvDdADYA3llKefu68gte2t52j8UlsMPLX8BhEfT88le7TDQPJVmFARazZVmuMkVEB4jtdYvLB9Tlg0rvswl3XP6ymr3mWQhtLU93n4ROFgpbl78WrJhjUqS0qxDaOPwLprViBq9iSti8XCem24V+H8P+3gpHSyuCPoRNuP9R+GzN7+WauvzFpu+vGtjknZcEHhsHedf8+4GEB7jyJu77KxGB1nbs1RxelrCmdFnCGq/8kpbPIwk8NvdeJD39jlsek1f+uvacuhQmY94vH1BPCH0wBEZtM4V1NB3FYRiTSZAQBbImVkXL6G6Xjg03YYnsxtxbKf+Wl4u4T2/F+aFZbm2XvbzAWfY3pu0bNYv3y1CbsMT4arRYYwgwBx3rInw2z18U/osZvDo7HLwpPADsJG0jy1/TPYzLXnbJPNZS/mJE8/iprTOyiDBI/Cq3DNDiduNYsa9LHRO7BW114k5wX8zPLwP4ZgDfN5eBiDYAvhfA1wH4KID3E9F7Sim/8nyqmEgkEonE5SCXvZ4xSikfBvxMqoM3AvhIKeXXOe8PAngLgIMfPwMVvGp7s5hnPQMErDeDX2JB1piGhyuLrlDEofZKyuIEBigMXrIzT8m7lgFCvR8yx1oleNZwHOSydFmniDxBAAAH1ElEQVSwYFZaBYv++pbHi4Ln3qzOpsffFj2HiJo2hK2EoxBhslHcygx2H6gyFXh22KgI6QfdmrK2zPQoA8QMwJYdIYoQ+iXTDjF7j8FPG0dwlmXRsBbCxIiDROkQP4Z8GtxWxdFWqN8wP3xt2ZexbpglET9r9AS1PA9RUq1yNDgKrQ0MA8sWLO2Qykn3ikh6Z94O/GCVnWdgWnbHMlje/D3O/DXdNIO0HBmL0/5O7oG7TV7gvAnODrdhHwCGvdT/wHKHGbYP4ztMiLkjDAGa4q0hQBGnjN4QQBggYWFeZUNiRCZ3BfMT04qew8dN/apBhmf6Ft/rkelR4XNgjWx/iNk7Aenk8NniRdb8fB6A3zL7HwXw5fdUl0QikUgkzhup+Xl6ENFPAXhd59B3lVJ+eE0RnbTZaQQRfQeA7+Ddx+/7mu/+5RXXODW8BsDv3HclngHOtV3A+bYt23V6ONe2nWu7vvh5XuwT+L2f+Kny717zDIp+IfvmmX38lFK+9imL+CiALzD7nw/g5YXrvQPAOwCAiD5QSvmyp7z+C4ds1+nhXNuW7To9nGvbzrldz/N6pZQ3Pc/r3TdW+be7J7wfwBcR0RcS0TWAbwHwnnuuUyKRSCQSiRPHvXz8ENE3EdFHAXwlgB8lop/g9DcQ0XsBoJSyA/CdAH4CwIcB/FAp5UP3Ud9EIpFIJBLng/uy9no3gHd30l8G8Gaz/14A732CS7zjyWv3QiPbdXo417Zlu04P59q2bFfiaFA5LrR3IpFIJBKJxEnjRdb8JBKJRCKRSNw5Tv7jh4jeSkQfIqKRiGYV/0T0JiL6NSL6CBG97XnW8UlBRJ9NRO8jov/J2z84k+83iOiXiOiDz9tC4Bgc6gOa8M/4+C8S0ZfeRz2PxYp2fTUR/T73zweJ6O/cRz2PBRF9PxF9nIi6biNOtb+AVW07uT4joi8gov9CRB/md+Lf6OQ5yT5b2bZT7LOHRPTfiegXuF1/v5PnJPvshUcp5aT/YYoP9sUAfhrAl83k2QD4XwD+KIBrAL8A4I/fd91XtO2fAHgb/34bgH88k+83ALzmvut7oC0H+wCT3uvHMPl4+goAP3vf9b6jdn01gB+577o+Qdv+HIAvBfDLM8dPrr+OaNvJ9RmA1wP4Uv79mQD+xzk8Y0e07RT7jAB8Bv/eAvhZAF9xDn32ov87eeanlPLhUsqvHcimoTJKKTcAJFTGi463AHgX/34XgG+8x7o8Ldb0wVsA/Ksy4b8B+Cwiev3zruiRONWxdRCllJ8B8H8WspxifwFY1baTQynlY6WUn+ffn8BkJft5IdtJ9tnKtp0cuB8+ybtb/heFuCfZZy86Tv7jZyV6oTJO4cF5bSnlY8D08AP43Jl8BcBPEtHPsafrFxFr+uAU+2ltnb+Sqe0fI6I/8Xyq9sxxiv11DE62z4jojwD405iYBIuT77OFtgEn2GdEtCGiDwL4OID3lVLOrs9eRLzIsb0U9JxDZTxPLLXtiGK+qpTyMhF9LoD3EdGv8sz2RcKaPnhh+2kBa+r88wD+cCnlk0T0ZgD/EcAXPfOaPXucYn+txcn2GRF9BoB/D+BvllL+XzzcOeVk+uxA206yz0opewB/iog+C8C7iehLSilWi3bSffai4iQ+fspzDpXxPLHUNiL6bSJ6fSnlY0xzfnymjJd5+3EiejempZgX7eNnTR+8sP20gIN1ti/pUsp7ieifE9FrSikvZMybI3CK/bUKp9pnRLTF9HHwr0sp/6GT5WT77FDbTrXPBKWU/0tEPw3gTQDsx8/J9tmLjEtZ9jrVUBnvAfBt/PvbADQsFxG9mog+U34D+Hr4B+dFwZo+eA+Ab2Xrhq8A8Puy7PcC42C7iOh1RET8+42Ynrvffe41vXucYn+twin2Gdf3XwL4cCnln85kO8k+W9O2E+2zz2HGB0T0EoCvBfCrIdtJ9tmLjpNgfpZARN8E4HsAfA6mUBkfLKV8AxG9AcA7SylvLqXsiEhCZWwAfH85jVAZbwfwQ0T07QB+E8BbgSkMCLhtAF6LiSoFpv78gVLKj99TfWcx1wdE9Nf4+L/A5M37zQA+AuAVAH/1vuq7Fivb9RcB/HUi2gH4NIBvKaW88LQ1Ef0bTBY0r6EpHM3fxSTIPNn+Eqxo2yn22VcB+CsAfok1JADwtwH8IeDk+2xN206xz14P4F1EtMH0sfZDpZQfOfX34ikgPTwnEolEIpG4KFzKslcikUgkEokEgPz4SSQSiUQicWHIj59EIpFIJBIXhfz4SSQSiUQicVHIj59EIpFIJBIXhfz4SSQSiUQicVHIj59EIpFIJBIXhfz4SSQSq0BEf4aIfpGIHrJn8Q8R0Zfcd70SiUTiWKSTw0QisRpE9A8BPATwEoCPllL+0T1XKZFIJI5GfvwkEonV4Phl7wfwCMCf5YjUiUQicVLIZa9EInEMPhvAZwD4TEwMUCKRSJwckvlJJBKrQUTvAfCDAL4QwOtLKd95z1VKJBKJo3HyUd0TicTzARF9K4BdKeUHOAr1fyWiP19K+c/3XbdEIpE4Bsn8JBKJRCKRuCik5ieRSCQSicRFIT9+EolEIpFIXBTy4yeRSCQSicRFIT9+EolEIpFIXBTy4yeRSCQSicRFIT9+EolEIpFIXBTy4yeRSCQSicRFIT9+EolEIpFIXBT+P4K1iKpFSZXGAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "# The next command ensures that plots are shown within the IPython notebook\n", "%matplotlib inline\n", "\n", "Nx=200\n", "Ny=200\n", "xmin, xmax, ymin, ymax=[-1,3,-1,3]\n", "plot_grid = np.mgrid[xmin:xmax:Nx*1j,ymin:ymax:Ny*1j]\n", "points = np.vstack((plot_grid[0].ravel(),\n", " plot_grid[1].ravel(),\n", " np.array([0.5]*plot_grid[0].size)))\n", "plot_me = np.zeros(points.shape[1], dtype=np.complex128)\n", "\n", "x,y,z = points\n", "bem_x = np.logical_not((x>0) * (x<1) * (y>0) * (y<1) * (z>0) * (z<1))\n", "\n", "slp_pot= bempp.api.operators.potential.helmholtz.single_layer(\n", " bempp_space, points[:, bem_x], k)\n", "dlp_pot= bempp.api.operators.potential.helmholtz.double_layer(\n", " trace_space, points[:, bem_x], k)\n", "\n", "plot_me[bem_x] += np.exp(1j * k * (points[0, bem_x] * d[0] \\\n", " + points[1, bem_x] * d[1] \\\n", " + points[2, bem_x] * d[2]))\n", "plot_me[bem_x] += dlp_pot.evaluate(dirichlet_fun).flat\n", "plot_me[bem_x] -= slp_pot.evaluate(neumann_fun).flat\n", "\n", "fem_points = points[:, np.logical_not(bem_x)].transpose()\n", "fem_val = np.zeros(len(fem_points))\n", "for p,point in enumerate(fem_points):\n", " result = np.zeros(1)\n", " u.eval(result, point)\n", " fem_val[p] = result[0]\n", "\n", "plot_me[np.logical_not(bem_x)] += fem_val\n", "\n", "plot_me = plot_me.reshape((Nx, Ny))\n", "\n", "plot_me = plot_me.transpose()[::-1]\n", "\n", "# Plot the image\n", "from matplotlib import pyplot as plt\n", "fig=plt.figure(figsize=(10, 8))\n", "plt.imshow(np.real(plot_me), extent=[xmin, xmax, ymin, ymax])\n", "plt.xlabel('x')\n", "plt.ylabel('y')\n", "plt.colorbar()\n", "plt.title(\"FEM-BEM Coupling for Helmholtz\")\n", "plt.show()" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.8.5" } }, "nbformat": 4, "nbformat_minor": 4 }