{ "cells": [ { "cell_type": "markdown", "id": "d25d6269", "metadata": {}, "source": [ "# The OSRC preconditioner for high-frequency scattering\n" ] }, { "cell_type": "markdown", "id": "bcc4387b", "metadata": {}, "source": [ "### Background\n", "\n", "Solving acoustic scattering problems with the BEM incurs increasing computational costs for higher frequencies. The mesh needs to be fine enough to cover several elements per wavelength, thus increasing the number of degrees of freedom and the size of the discretization matrix. Furthermore, the number of iterations it takes for GMRES to converge to an accurate solution often increases with frequency as well. This slower convergence strongly depends on the specific boundary integral formulation and can be alleviated with well-designed preconditioners. This tutorial presents a preconditioning strategy that uses On-Surface Radiation Condition (OSRC) operators and which is especially effective at high frequencies.\n", "\n", "Let us consider acoustic scattering of an incoming wave $u^\\text{inc}$ at a sound-hard obstacle $\\Omega$ with surface $\\Gamma$. The corresponding Helmholtz problem reads\n", "$$\n", "\\begin{cases}\n", "\\Delta u^+ + k^2 u^+ = 0, &\\text{in }\\Omega^+; \\\\\n", "\\frac{\\partial u^+}{\\partial\\nu} = -\\frac{\\partial u^\\text{inc}}{\\partial\\nu}, &\\text{on }\\Gamma; \\\\\n", "\\lim_{|\\mathbf{x}|\\to\\infty}|\\mathbf{x}|\\left(\\nabla u^+ \\cdot \\frac{\\mathbf{x}}{|\\mathbf{x}|} - \\imath k u^\\text{+}\\right) = 0\n", "\\end{cases}\n", "$$\n", "where $u^+$ is the scattered field and $k$ the wavenumber in the exterior domain $\\Omega^+ = \\mathbb{R}^3\\setminus\\Omega$. Hence, $u = u^+ + u^\\text{inc}$ is the total field.\n", "\n", "The scattered field $u^\\text{+}$ can be represented as\n", "$$\n", "u^+ = \\mathcal{K}\\phi \\quad\\text{in } \\Omega^+,\n", "$$\n", "where $\\mathcal{K}$ is the double layer potential operator defined by\n", "$$\n", "\\left[\\mathcal{K}\\phi\\right](\\mathbf{x}) = \\int_{\\Gamma} \\frac{\\partial}{\\partial\\nu_\\mathbf{y}} G(\\mathbf{x},\\mathbf{y})\\phi(\\mathbf{y}) \\,\\mathrm{d}\\mathbf{y} \\quad\\text{for } \\mathbf{x}\\in\\Omega^+,\n", "$$\n", "with\n", "$$\n", "G(\\mathbf{x},\\mathbf{y}) = \\frac{\\mathrm{e}^{\\imath k|\\mathbf{x}-\\mathbf{y}|}}{4\\pi|\\mathbf{x}-\\mathbf{y}|}\n", "$$\n", "Green's function of the three-dimensional Helmholtz equation. Here,\n", "$$\\phi = \\gamma_\\text{D} u \\quad\\text{on } \\Gamma\n", "$$\n", "is the unknown Dirichlet trace of the total acoustic field.\n", "\n", "The Burton-Miller formulation to compute $\\phi$ is given by\n", "$$\n", "\\left(\\tfrac{1}{2}\\mathsf{Id}-\\mathsf{K}+\\eta\\mathsf{W}\\right)\\phi = \\gamma_0 u^\\text{inc}+\\eta \\gamma_\\text{N} u^\\text{inc}\n", "$$\n", "for some $\\eta \\neq 0$. Here, $\\mathsf{Id}$, $\\mathsf{K}$ and $\\mathsf{W}$ are the identity, double-layer and hypersingular boundary operators, respectively, and $\\gamma_\\text{N}$ the Neumann trace operator.\n", "\n", "A typical choice for the Burton-Miller parameter is $\\eta = \\imath/k$ but one could also use a boundary integral operator for $\\eta$. A perfect choice is the exterior Neumann-to-Dirichlet map $\\mathsf{NtD}$ since\n", "$$\n", "\\tfrac{1}{2}\\mathsf{Id}-\\mathsf{K}-\\mathsf{NtD}\\circ\\mathsf{W} = \\mathsf{Id}\n", "$$\n", "and the linear system becomes trivially to solve. However, the NtD-map is a complicated non-local operator that cannot be calculated explicitly for general scattering problems. The idea of the OSRC preconditioner is to use the pseudodifferential operator approximation\n", "$$\n", "\\mathsf{NtD}\\approx \\frac{1}{\\mathrm{i}k}\\left(1+\\frac{\\Delta_{\\Gamma}}{(k+\\mathrm{i}\\epsilon)^2}\\right)^{-\\frac{1}{2}}\n", "$$\n", "for a regularization parameter $\\epsilon>0$. We localize the square-root operator by a Padé approximation. Details of this OSRC operator are given in Antoine & Darbas (2007).\n", "\n", "The OSRC-approximated NtD operator and its inverse the Dirichlet-to-Neumann operator are available in Bempp. This notebook demonstrates how to use these for high-frequency scattering computations.\n" ] }, { "cell_type": "markdown", "id": "0403b4cf", "metadata": {}, "source": [ "### Implementation\n", "\n", "Let us start with importing the relevant libraries." ] }, { "cell_type": "code", "execution_count": 1, "id": "9f8df394", "metadata": {}, "outputs": [], "source": [ "import bempp.api\n", "import numpy as np" ] }, { "cell_type": "markdown", "id": "6b344dda", "metadata": {}, "source": [ "The following defines the wavenumber $k$ and the Dirichlet and Neumann data of the incident plane wave travelling in the positive $x$ direction. Notice that we need both traces in the Burton-Miller formulation.\n", "\n", "The wavenumber should be adapted to the computing resources available. A higher wavenumber requires more memory and compute time." ] }, { "cell_type": "code", "execution_count": 2, "id": "efa54b9f", "metadata": {}, "outputs": [], "source": [ "k = 7" ] }, { "cell_type": "code", "execution_count": 3, "id": "5d4d326f", "metadata": {}, "outputs": [], "source": [ "@bempp.api.complex_callable\n", "def dirichlet_fun(x, n, domain_index, result):\n", " result[0] = np.exp(1j * k * x[0])\n", "\n", "@bempp.api.complex_callable\n", "def neumann_fun(x, n, domain_index, result):\n", " result[0] = 1j * k * n[0] * np.exp(1j * k * x[0])" ] }, { "cell_type": "markdown", "id": "6eb09184", "metadata": {}, "source": [ "For this example we will use an elongated ellipsoid. The element size is chosen to roughly correspond to six elements per wavelength. The function space consists of continuous, piecewise-linear basis functions." ] }, { "cell_type": "code", "execution_count": 4, "id": "418c7518", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "The discrete function space has 1807 degrees of freedom.\n" ] } ], "source": [ "wavelength = 2*np.pi / k\n", "h = wavelength / 6\n", "grid = bempp.api.shapes.ellipsoid(3, 1, 1, h=h)\n", "\n", "space = bempp.api.function_space(grid, \"P\", 1)\n", "print(\"The discrete function space has {0} degrees of freedom.\".format(space.global_dof_count))" ] }, { "cell_type": "markdown", "id": "e5faf8b2", "metadata": {}, "source": [ "We can now form the Burton-Miller operator and its OSRC-preconditioned variant. The NtD operator only works for continuous function spaces as the OSRC approximation assembles a Laplace-Beltrami operator. Optionally, the damped wavenumber $k + \\imath\\epsilon$ in the OSRC approximation, and the number of expansion terms and branch cut angle for the Padé series can be specified for the NtD operator to fine-tune its accuracy." ] }, { "cell_type": "code", "execution_count": 5, "id": "3d9270b0", "metadata": {}, "outputs": [], "source": [ "identity = bempp.api.operators.boundary.sparse.identity(space, space, space)\n", "dlp = bempp.api.operators.boundary.helmholtz.double_layer(space, space, space, k)\n", "hyp = bempp.api.operators.boundary.helmholtz.hypersingular(space, space, space, k)\n", "ntd = bempp.api.operators.boundary.helmholtz.osrc_ntd(space, k)\n", "\n", "burton_miller = .5 * identity - dlp + (1j/k) * hyp\n", "osrc_bm = .5 * identity - dlp - ntd * hyp" ] }, { "cell_type": "markdown", "id": "a871f392", "metadata": {}, "source": [ "We next assemble the right-hand side, which also includes the OSRC operator." ] }, { "cell_type": "code", "execution_count": 6, "id": "6c0d2537", "metadata": {}, "outputs": [], "source": [ "dirichlet_grid_fun = bempp.api.GridFunction(space, fun=dirichlet_fun)\n", "neumann_grid_fun = bempp.api.GridFunction(space, fun=neumann_fun)\n", "\n", "rhs_fun_bm = dirichlet_grid_fun + (1j/k) * neumann_grid_fun\n", "rhs_fun_osrc = dirichlet_grid_fun - ntd * neumann_grid_fun" ] }, { "cell_type": "markdown", "id": "83216d5c", "metadata": {}, "source": [ "We can now solve the Burton-Miller formulation using GMRES. We use a strong form discretisation: this automatically performs mass matrix preconditioning of the entire linear system." ] }, { "cell_type": "code", "execution_count": 7, "id": "9d5fea88", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "The linear system was solved in\n", " 16 iterations for the Burton-Miller formulation, and\n", " 6 iterations with the OSRC preconditioner.\n" ] } ], "source": [ "from bempp.api.linalg import gmres\n", "\n", "sol_bm, info, it_count_bm = bempp.api.linalg.gmres(\n", " burton_miller, rhs_fun_bm, use_strong_form=True,\n", " return_iteration_count=True)\n", "\n", "sol_osrc, info, it_count_osrc = bempp.api.linalg.gmres(\n", " osrc_bm, rhs_fun_osrc, use_strong_form=True,\n", " return_iteration_count=True)\n", "\n", "print(\"The linear system was solved in\")\n", "print(\" {0} iterations for the Burton-Miller formulation, and\".format(it_count_bm))\n", "print(\" {0} iterations with the OSRC preconditioner.\".format(it_count_osrc))" ] }, { "cell_type": "markdown", "id": "b6f5e6c6", "metadata": {}, "source": [ "The number of iterations is indeed smaller when using the OSRC preconditioner, without incurring much computational overhead.\n", "\n", "We now want to plot the radar cross section in the $z=0$ plane. To compute it we use the far-field operators implemented in Bempp." ] }, { "cell_type": "code", "execution_count": 8, "id": "5d4e2f8d", "metadata": {}, "outputs": [], "source": [ "theta = np.linspace(0, 2 * np.pi, 400)\n", "points = np.array([np.cos(theta), np.sin(theta), np.zeros(len(theta))])\n", "\n", "dlp_far_field = bempp.api.operators.far_field.helmholtz.double_layer(space, points, k)\n", "far_field = dlp_far_field * sol_osrc\n", "\n", "max_incident = np.max(np.abs(dirichlet_grid_fun.coefficients))\n", "radiation_pattern = (np.abs(far_field/max_incident)**2).ravel()\n", "db_pattern = 10 * np.log10(4 * np.pi * radiation_pattern)" ] }, { "cell_type": "code", "execution_count": 9, "id": "58842986", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAe8AAAH3CAYAAACSDNc7AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjQuMiwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8rg+JYAAAACXBIWXMAAAsTAAALEwEAmpwYAAEAAElEQVR4nOydd1hb59nG71cDCYEGCCSxMXsZbMCA8czee+80ezZJm522SZu0zejI1yZtVtvs2TRx9nKGbTwxBozZe4qtCdrv94cQxjZb88D5XeGKLem85xGWzn2e8T4PoZSChYWFhYWFhTlw/G0ACwsLCwsLy8JgxZuFhYWFhYVhsOLNwsLCwsLCMFjxZmFhYWFhYRiseLOwsLCwsDAMVrxZWFhYWFgYBiveLCzLEEJIFiGkfJbnXyWEPDHPtQSEkHpCiMJzFrKwsMwGK94sLAEAIaSdEDJOCDEQQtQT4hk65XkJIeRZQkjnxGuaJ/4eMfH8ekLITkKIlhAyQggpI4SsmeWUjwP40zxt20wIcUyc10AI6SGE/Nb1PKXUDODfAB5Y3LtnYWFZKKx4s7AEDmdRSkMBrAKwGsBDAEAICQKwFUA2gFMBSACUAhgGUEQIkQD4DMDfAYQDiAHwWwDm6U5CCIkCcByAjxdgWy+lNHTCvvUArieEnDvl+bcBXEMIESxgTRYWlkXCijcLS4BBKVUD+BpOEQeAqwHEAziPUlpLKXVQSgcopY9TSr8AkDZx3DuUUjuldJxS+g2ltHqGU5wEoIJSanI9QAhZTQipIIToCSHvARDOYl8bgJ0AsqY81g1gFEDJYt83CwvL/GHFm4UlwCCExAI4DUDzxEMnAviKUmqY4ZBGAHZCyGuEkNMIIWFznGIlgIYp5wuC0wt/A07P/QMAF8xiXyqAdQB2H/VUHYC8Oc7NwsLiAVjxZmEJHD4mhOgBdAEYAPDoxONyAH0zHUQp1cEZyqYAXgYwSAj5hBCinOEQGQD9lL+XAOADeJZSaqWU/hfAvqOOiSaEaAghOjhvFvYA2HHUa/QTa7OwsHgZVrxZWAKHcymlYgCbAWQAiJh4fBhA1GwHUkrrKKXXUkpjAeQAiAbw7AwvHwUgnvL3aAA99MgpRR1HHdNLKZVRSiVwCvQ4gNeOeo0YgGY2O1lYWDwDK94sLAEGpfQnAK/icDX4dwBOIYSEzPP4+onjc2Z4STUm8uQT9AGIIYSQKY/Fz7K+Fs4CtbOOeioTQNV8bGRhYXEPVrxZWAKTZwGcRAhZBWcuugvAh4SQDEIIhxAiJ4Q8TAg5feKxX07kykEIiQNwGY7NSbv4FkA+IcRVlLYLgA3AzwkhPELI+QCKZjJsYgvbpQAOTXksBs58+UznZGFh8SCseLOwBCCU0kEArwP49cQ+6hMB1MMpvDoAe+EMq++BM9dcDGAPIcQIp4DWAPjlDGv3A/gewDkTf7cAOB/AtXCG1C8B8L+jDot27fOGM6QeDuCKKc9fDuC1CVtZWFi8DDkyzcXCwrIcIIRkwZmzLqJuXgQm9nZXAdhIKR3whH0sLCyzw4o3CwsLCwsLw2DD5iwsLCwsLAyDFW8WFhYWFhaGwYo3CwsLCwsLw2DFm4WFhYWFhWGw4s3CskQghNxFCKkhhBwihNw98Vg4IeRbQkjTxP/Dprz+GUJIOSFkk9+MZmFhWRSseLOwLAEIITkAboSzuUoegDMnBog8CGArpTQVzrGiD068PmPi0I0Abve9xSwsLO7AijcLy9IgE8BuSukYpdQG4CcA58HZiMXVg/w1AOdO/JkLwAHnMBMCFhYWRsGKNwvL0qAGwMaJtqkiAKcDiAOgpJT2AcDE/xUTfz4EQATnZLB/+sdkFhaWxcLztwEsLCzuQymtI4Q8BWf7VAOcHc9scxxzpy9sY2Fh8Tys583CskSglP6LUppPKd0IYARAE4B+QkgUAEz8n21fysKyBGDFm4VliUAIUUz8Px7OQSPvAPgEwDUTL7kGwBb/WMfCwuJJ2N7mLCxLBELIdgByAFYAv6CUbiWEyAG8D+d87k4AF1FKR/xoJgsLiwdgxZuFhYWFhYVhsGFzFhYWFhYWhsGKNwsLCwsLC8NgxZuFhYWFhYVhsOLNwsLCwsLCMFjxZmFhYWFhYRiseLOwsLCwsDAMVrxZWFhYWFgYBiveLCwsLCwsDIMVbxYWFhYWFobBThVjYWEQhBAOgBAAYgChE/8/4s/BwcESHo/H53A4XA6HwyWEcAFwCCHU4XA4KKV2h8NhdzgcdrPZPGa1WnUA9BM/hqP/TCk1+/6dsrCwzAbbHpWFxc8QQngAlACiAERxOJwYuVyeIhQKEwHEWq1WJYfDEfB4PC6PxyMhISE0NDSUSqVSIpVKOVKplCuTyXgymYwvk8kEISEhhMfjgcPhTP4QQtDQ0IDU1FQ4HI7JH7PZDJ1OZ9NoNBaNRmPVaDQ2jUbj0Ol0Dp1ORwwGAzGZTLDb7Q6bzWbjcDg6Ho/Xa7PZOvR6fYter+8A0Augb+JHR9mLCguL12HFm4XFBxBCQgCkAEiTy+WrQkJCVlut1lRCSKhQKOSqVCoaFxfHSUxMFCYkJIiio6M5UVFRiI6OhlKpBJ/Pd9uGH3/8EZs3b1708ZRS6HQ69PX1obe3F319feju7ja3t7ePtbe3W3t6eqDVamGxWKw8Hq8PwKGBgYF9FoulHkAjgB5KqcPtN8LCwsKKNwuLJyGEqACsDg0NXRkWFlbocDgyKaVyiUTCy8zMJHl5eSHZ2dnBaWlpSElJgUgk8plt7or3fKGUoq+vD42NjWhoaHBUV1frqqurLZ2dncRisYzz+fx2q9Va3d/fX04prQFQw4bmWVgWBiveLCyLhBASDaAgIiJig1Ao3Gi32xNiYmI469evF65atUqcnp5OUlNTIZfL/W0qAN+J92zYbDa0t7ejsbERtbW1lt27d+sOHDjgGB8fN3C53MrBwcHvzGbzPgAHWUFnYZkZVrxZWOYBISQSwNrIyMiNAoFgo91uj4+LiyMbNmwILi0tFRcUFCA+Ph6EEH+bOiOBIN4zMT4+jurqauzbt8+2bds2TUVFBR0bGzNwudyqoaGhrSaTaTeASkqpzd+2srAEAqx4s7BMw4RYb4yJiTnPbrevVyqVwSeddNKkUMfFxQW0UE9HIIv3dJhMJlRXV6O8vNy2detWzb59+xwOh6PTYDBs0Wq1XwM4wIo5y3KFFW8WFhwj1huUSmXw6aefLjr55JNDSkpKIBQK/W2i2zBNvKejo6MDP/74o+Ozzz4b3b17t4NS2m4wGD7RarVfgfXMWZYRrHizLEsIIXwAG6Oioq6klB6vUCiEZ5xxRsjJJ58cUlxcjODgYH+b6HGWgngfjUvMP/3005E9e/ZQSmmbRqN512g0/o9S2uFv+1hYvAUr3izLBkJIGJfLPT0qKupaSmneySefzL/kkktkGzduXJJifTRLUbyPpq2tDZ9++qn1rbfe0nR2dhptNttHQ0ND7wDYz25TY1lKsOLNsqQhhCRLJJILQ0JCrpRIJMqLL7445IILLhDl5uYyLmftLstBvKei0+nw9ddf03feeWd4165ddi6Xu6Onp+dVAFsppeP+to+FxR1Y8WZZchBCchUKxY0cDufslJSU4CuuuCLs7LPP5kVHR/vbNL+y3MR7KjabDbt27cIHH3yg++STT8w2m615cHDwJYvF8iGlVO9v+1hYFgor3ixLAkJIXFhY2HVBQUHX5uTkhN56663y0047jfiyCUqgs5zF+2gaGxvx1ltvjb322mtjNpvtQE9Pz/8B+IZSavW3bSws84EVbxbGQgiRiUSiS6VS6e0qlUp5yy23hF100UW8sLAwf5sWkLDifSyUUuzfvx+vvPKKdsuWLWZCyLd9fX3PA9jN9mhnCWRY8WZhFIQQAYfDOUOlUt0lFAqzrrvuupCrrroqOD4+3t+mBTyseM+OzWbD1q1b8cILLwzv2rXLZLfb3x0aGnqJUtrob9tYWI6GFW8WRkAIyVAqlfdxudyzLrzwQsH1118vyc3N9bdZjIIV7/ljNBqxZcsWx/PPPz/c2to6pNFonjGZTO+yhW4sgQIr3iwBCyFEIBAILggPD39wxYoVUffff7/8jDPOIDweO4Z+KpRS2O122Gy2yf87HA64vtuu/x84cAD5+fkAMFlpz+VywefzweVy4RojynIk3d3d+Oc//2l89dVXx+x2+6f9/f3PUErr/W0Xy/KGFW+WgIMQEhcZGXkvj8e79PLLLxfdcccdoYmJif42yydQSmGxWGAymWA2mzE+Pg6z2QyTyQSTyQSr1QqbzYajv7cu8eXxeOByueByuQAOizQhBN3d3YiJiTlC1F1i7/qZaV2BQAChUDj5/6k/rnMtdWw2Gz7//HP69NNPD7W1tfUODg4+brPZtrBd3Vj8ASveLAEBcarMpujo6MfCw8OzHnjggfCLL76YGxQU5G/TPI7FYoHRaITBYIDBYIDRaMTY2BgcDmcPkaCgoElhFAgECA4OnhTNoKAgcLncRXnICw2bU0rhcDhgtVphNpuPuIlw/ZjNZtjt9km7Q0JCEBoaitDQUISEhEAkEi1Jb76lpQV//vOfdf/73/+MVqv1lZGRkb9RSof8bRfL8oEVbxa/QgjhCYXCy6RS6WPr16+XPfTQQ+EFBQX+NsttKKUYGxuDVqudFGiDwQC73Q4+nz+tyHnbg/VmztsVMTj6psRoNAIAhELh5PsVi8WQSCTg8/lescWXGI1GvP7665Y///nPOqPR+INarX6YUtrsb7tYlj6seLP4BUKIUCwW3ygSiR645JJLJA899JBYpVL526xFQSmF0WiEVquFRqOBRqOBxWKBSCSCVCqFWCyeFGt/5uv9VbBGKYXJZJoUdb1eD61WC7vdDrFYDKlUCplMBqlUylhBp5Tiyy+/pI888shwf3//gb6+vvsppZX+totl6cKKN4tPIYRIwsPD7xEIBLfcdNNNkrvvvlskk8n8bdaCMJlMGB4exujoKLRaLSwWC0JCQiCTySZFSCAQ+NvMYwi0anOHwwGDwTB5w+MS9NDQUMhkMoSHh0MmkzEu7L5jxw489NBDwy0tLc19fX33UUq3+9smlqUHK94sPoEQEhkZGflIUFDQ5b/85S/FN998s5Ap3c9cYj00NITR0VHw+XxEREQgLCwsYIV6OgJNvKeDUgq9Xg+NRoORkRFoNBoEBQUhIiICERERjBLzqqoqPPLIIyP79+/vGxgYuN/hcHzJNn5h8RSseLN4FUJIrEqleiI4OPiMX/3qV7Irr7ySF+hFaDOJtUs8mFpdzQTxno7x8fHJfw+NRgOBQAC5XM4YMW9ubsZvf/tbzbfffjui0+l+PT4+/i474YzFXVjxZvEKhJBwhULxB7FYfMGTTz4Zfv7553MC9SJLKYVGo0F/fz/6+/vB4/GWhFgfDVPF+2jGx8cxNDQ0mboIDQ2FSqWCUqlEIN8Y9vb24oknntD973//GxoeHv65zWb7gvXEWRYLK94sHoUQIpLL5Q8FBwff/Pjjj8uuuuoqfiCKn91ux+DgINRqNUZHRyGRSBAVFYXIyEjGFk3NxVIR76m4wuxqtRr9/f3gcDhQKpVQqVQIDQ31t3nT0tnZiXvvvXd0+/btXWq1+hZK6S5/28TCPFjxZvEIhBC+VCq9VSQSPfSLX/xCdueddwoDLRdsNpuhVquhVqsxPj6OyMhIqFQqhIWFBXzo1RMsRfE+GpPJhP7+/mP+jcPDwwNufvuhQ4dw5513DtfV1dVOiHitv21iYQ6seLO4BSGEBAcHXyqRSJ665pprwh555JFQiUTib7Mmsdls6OvrQ3d3N2w2G6KioqBUKiEWi/1tms9ZDuI9laOjK5GRkYiNjYVUKg0oIS8rK8Mdd9wxrFart6nV6rsopV3+tokl8GHFm2XRcLncEyIjI/9x5plnKn//+99LlUqlv00C4NyCNDg4iK6uLhgMBkRFRSE2NhYhISH+Ns2vLDfxnorD4cDAwAC6u7thMBgQHR2N2NhYBMqOB0opPv/8c8c999wzotPp/jcwMPAQpXTE33axBC6seLMsGEJIvFKpfDU/P3/Vc889F5aUlORvk0ApxejoKLq7uzE8PIyIiAjExsZCJpMFlJflT5azeE/FarVORmMcDgdiY2MRHR0dEMVudrsdb775pvXhhx8eNRgMv9XpdC+wleks08GKN8u8IYQIIiIifi2RSG5+8cUX5SeeeKLfVdFkMqGzsxO9vb0Qi8WIi4tDRETEsshhLxRWvI9lfHwcPT096O3thVAoRGJiIiIjI/1+w6fX6/Hwww/rPvjgg97+/v6rKKXlfjWIJeBgxZtlXgQFBZ0aHh7+4s9//vOIe++9V+RPL4VSiqGhIbS3t2N8fBzx8fGIiYlZslXinoIV79nRarVob2/HyMgIYmJiEB8fD6FQ6Febamtrcc0114x0dXV909/ffweldNivBrEEDKx4s8wKISRBqVS+WlhYmPfCCy+ExcbG+s0Wq9WKrq4udHZ2QiqVIjExkQ2LLwBWvOeH1WpFT08POjo6EBISgqSkJISFhfntc0YpxTvvvGO7//77R/R6/e90Ot0/2VA6CyveLNNCCBHI5fLfSKXSm/wdItfr9Whra8Pw8DDi4uIQHx8fEPlJpsGK98Jw1VG0trZibGwMCQkJiI2N9VvTnolQuv7999/vHRgYuIpSus8vhrAEBKx4sxwDj8fbEBER8ebPf/7zyHvvvTfYX0I5PDyMpqYmOBwOJCUlQalUsl62G7DivXhMJhPa29vR19eH6OhoJCUl+S1NU1tbi2uvvXaks7Pzy/7+/lsopQa/GMLiV1jxZpmEECJSKBTPJSUlnfPuu++GJyQk+NwGSikGBgbQ3NyMoKAgpKamgmlTxwIVVrzdx263o7OzE+3t7VAoFEhOTvZLXpxSildeecXy61//enBwcPBKu93+o8+NYPErrHizAJj0tt/6zW9+o7j11lsFvvZwKaXo7e1FS0sLJBIJUlJSAra9JVNhxdtzOBwO9PT0oLW1FTKZDCkpKX7pI9DV1YXLLrtspLm5+fP+/v7bWC98+cCK9zJnwtv++4oVK8599913wxMTE316fofDga6uLrS1tSEiIgLJyckIDg72qQ3LBVa8PQ+lFGq1Gs3NzRCJREhNTYWvOwxSSvHSSy9ZHn300YGhoaErbTbbTz41gMUvsOK9jOHxeOsiIiLe/vWvf6287bbbfOptOxwOdHR0oL29HVFRUUhKSmKL0LwMK97eg1I6WaPB5XKRkZHhcxHv7OzEpZdeOtrS0vLpwMDAbZRSo08NYPEprHgvQya87f9LTEw8/7333vOpt+0Kjzc1NUGlUiE5OZndn+0jWPH2DSMjI6irq4NIJEJGRoZPI0mUUrz44ovmxx57bHBoaOgKm822zWcnZ/EprHgvMwghq5RK5ccPP/yw8o477hD6shPZwMAA6uvrIZPJkJaW5vcGGMsNVrx9B6UU/f39aGhogFwuR1pamk8jS52dnbjkkktGWltbPxgYGPg5pdTis5Oz+ARWvJcJhBASHh7+i8jIyIe2bNkiT09P99m5R0dHUVdXB4FAgIyMjGU/IMRfsOLteyil6O7uRnNzM6Kjo5GcnAwej+eTczscDjz11FNjzz77bPvAwMCZlNI2n5yYxSew4r0MIITIFArFf88444zCf/zjH1JfebwGgwF1dXWw2+3IzMyEVCr1yXlZpocVb/9ht9vR3t6Ozs5OJCYmIiEhwWf99/fs2YOLLrpocGRk5E6DwfCeT07K4nVY8V7iEEKKlUrlh3/729+UF198sU9u+W02GxoaGjA8PIysrCxERET44rQsc8CKt/+xWq1oaWlBf38/srKyEBkZ6ZPzjo6O4tJLL9VUVVV92t/ffxOl1OSTE7N4DVa8lyiEEI5cLn84Kirqnk8++SR8xYoVXj+nqxitsbERK1asQEJCAtsRLYBgxTtwGBsbQ01NDTgcDrKzs31S1EYpxd/+9jfTH/7wh66JMHqj10/K4jVY8V6CEEIiFArFxxdddNHKv/zlLxJfFMro9XocPHgQIpEIWVlZ7LavAIQV78Cjv78fdXV1iI2NRVJSkk9C6QcOHMAFF1wwNDQ0dJ9Op3vV6ydk8QqseC8xeDzeusjIyA9eeuklxVlnneX1CQpWqxWNjY0YGRlBTk4OwsLCvH1KlkXCindgYrfb0dzcjL6+PmRnZ/sklK7X63HVVVdpdu3atXVgYOBKNozOPFjxXkKEhYXdGhMT88SXX34ZHhcX59VzsSFy5sGKd2Dj61A6pRTPP/+86Xe/+13b4ODgiZTSXq+ekMWjsOK9BCCE8BQKxUslJSXnvfvuuzJvf+nHx8dRVVUFoVDIhsgZBCvezGBgYAC1tbWTVenevin+6aef6GWXXdbf19d3DqV0r1dPxuIxfNehg8UrEELCFQpF2R133HHxxx9/7FXhppSio6MDe/bsQXJyMlatWsUKNwuLh1EoFFi/fj30ej12796NsbExr55v06ZNZOfOnaq0tLQvZDLZz7x6MhaPwYr3AiGE/JsQMkAIqZny2OOEkGpCSCUh5BtCSPTE44mEkPGJxysJIS9MOWYzIaScEPK0G7ZkKpXKAy+//HL+r3/96xBv3qGPj49jz5490Gg0WL9+vc+2uLCwLEd4PB5WrlyJ1NRU7N27F+3t7fBmlDQxMRH79++XFxcX/0WpVP6DEOKRehlCCJcQcoAQ8tnE3x8jhPRMuSaePuW1z0xcEzd54txLHVa8F86rAE496rFnKKW5lNJVAD4D8Jspz7VQSldN/Nwy5fFbAWwAwCWEZCzUCJFIdFZCQsK2rVu3xp999tle279NKUVnZ+ekt52Xl+ezDlEsLMudiIgIrF+/HjqdzuteeGhoKL788kvZddddd1VkZOSPhBBPdFW6C0DdUY/9dco18QsAmHIN3Ajgdg+cd8nDivcCoZRuAzBy1GO6KX8NATCfW2TOxOscAObtMhNCSGRk5GO5ubmv7t+/PyI7O3u+hy4Yl7c9OjqKdevWsd42C4sf4PF4yM3NnfTCOzo6vOaFczgc/PGPfwx97rnnihUKRSUhJHWxaxFCYgGcAeCVebycC+e1kGIB18PlDCveHoIQ8ntCSBeAK3Ck571iImz0EyFkw5THXwGwEwCHUnr0nelM5+ArFIqPzjnnnHu2b98eLpfLPfcGjqK7u/sIb5ud/MXC4l9cXrhWq8Xu3bthMnlvd9fFF1/M//rrrxPj4uJ2BAUFbV7kMs8CuB9OUZ7KHRNpxn8TQsIAgFJ6CIAIwA4A/1zk+ZYVbLX5IiCEJAL4jFKaM81zDwEQUkofJYQIAIRSSocJIQUAPgaQfZSnPt9zhkZGRn5733335d53330iN9/CjNhsNhw8eBAOhwO5ubmsaC8h2GrzpcPg4CBqamqQnZ0NhULhtfOo1Wocd9xxw93d3bfp9fr353scIeRMAKdTSm8jhGwGcC+l9ExCiBLAEJwe9uMAoiil13nD9qUO63l7nrcBXAAAlFIzpXR44s/7AbQASFvogoSQiMjIyD1/+ctf8r0p3DqdDmVlZQgPD0d+fj4r3CwsAUpkZCTWrl2LlpYW1NbWwuE42rn1DCqVCrt375ZnZmb+Mzw8/K4FHLoOwNmEkHYA7wI4nhDyJqW0n1Jqp5Q6ALwMoMgLZi8LWPH2AEflhc4GUD/xeKSrapMQkgQgFUDrAtdOUCgU+954442MK6+80iv7slxbwA4cOIDVq1ezDVdYWBiAUChESUkJeDwedu7cifHxca+cRyqVYtu2beHFxcWPKRSKp8k8Lg6U0ocopbGU0kQAlwL4nlJ6JSEkasrLzgNQM+0CLHPClg0vEELIOwA2A4gghHQDeBTA6YSQdDhzOx0AXFXlGwH8jhBiA2AHcAuldOTYVWc818qoqKivt2zZErVmzRpPvo1JrFYrqqqqwOPxsG7dOraSnIWFQRBCkJaWBrlcjt27dyMzMxMqlcrj5xEKhfjss89kN9xww81ffPFFDCHkakqpfRFLPU0IWQVn2LwdwM2etHM5wea8AxQej7chJibmw2+//TYyLW3BkfZ5odFoUFlZiZSUFMTGxnrlHCyBA5vzXtpYLBYcOHAAIpEI2dnZXhlyQinFb37zG8OLL764Z3Bw8Ey2J7r/YMU7ABGLxRdGR0e/8MMPP8ijo6O9co729nZ0dnYiPz8foaGhXjkHi2ehlMJqtcJkMsFsNsNiscBms834Qymd/AGcM53Dw8MBOD02LpcLHo83449AIIBQKIRQKPTJtCsW96GUorW1Fb29vSgsLPRaf/R//OMfpscee6x2cHDweEqp1isnYZkVVrwDjPDw8FuTk5Of+Pbbb8NlMpnH13c4HDh48CDsdjvy8vLA5Xp98BjLPHE4HBgfH4fBYIDBYIDRaJwUarvdGaHk8/mTghoUFAQejzejCHM4nMnaBUIIduzYgXXr1k0Kut1unxT6qX+22WywWq0wm82T53cVRAUFBUEoFCI4OBihoaEICQlBaGgo2yY3wBgZGUFVVRXy8vImb9g8zYcffmi79dZb2wYHBzdQSvu9chKWGWHFO4CQy+X35OTk/Prrr78OEwqFHl/fbDajvLwcSqUSycnJbFGan3A4HNDpdNBqtZNC7So2EolECA0NnRRGoVAIgUDgkVoEd8PmlFJYLBaYTKbJmwyj0QiDwQCr1Qoulzsp5hKJBFKp1OuTsVhmZnx8HPv27UNCQgISEhK8co6tW7c6Lr/88o6BgYF1lNI+r5yEZVpY8Q4Q5HL5fXl5eQ9/+eWXMoFA4PH1dTodKioqkJmZCaVS6fH1WabHbrdPCrVGo4FOpwOlFGKxGFKpFGKxGKGhoQgODvb6zZS3c942m21SzLVaLbRaLUwmE4KDgyGVSiGTySCTySAUCtkbRx9hs9lw4MABBAcHIysryyvpjx9//NFxySWXdA4MDKynlPZ4/AQs08KKdwAQERHxcH5+/n2fffaZzBvhx76+PjQ0NKCgoABisdjj67Mcxmq1YmRkBENDQxgZGYHD4Zj0QmUyGSQSid8q+v1RsEYphclkgkajmfwxmUwIDQ1FREQE5HI5xGIxK+ZehFKKxsZGjIyMoKCgwCspju3bt9MLL7ywa8ID7/b4CViOgRVvP6NSqX6TnZ39wFdffSXydFMUX3xplztHizWlFHK5fPInkBrdBEq1OaUUBoMBQ0NDGBoagl6vh1gsRkREBCIiIhAaGsqKuRfw9k38jz/+SK+44oqh3t7eAkppl8dPwHIErHj7kcjIyAcLCgoeeOyxx2Th4eHw5JYwu92OiooKr4bLlitGoxFqtRpqtRo2mw1yuXzSiwwksT6aQBHvo6GUQq/XY3h4GENDQzAYDAgPD4dKpUJERARbVOlBvJU+M5lMrn7r9JprrukaGBgoZUPo3oUVbz8RERHxy9WrV//q888/l/F4PFRUVEAikXhEwC0WC/bt24fY2FivFaosJyilGB0dhVqtxsDAAAQCAVQqFVQqFaMKsgJVvI/G4XBgZGQEarUaQ0NDEIlEUCqVUKlU8EY9yHLDbDZj3759iI+PR3x8vNvruYR75cqVkMvlrhB6x4SAs0VsXoIVbz8gl8vvzs3NffSrr76aLE5zOBweEfDx8XHs3bsXaWlpiIqKmvsAlmmhlGJkZATd3d0YGRmBVCqFSqWCQqFgbBc6poj3VFwhdrVajf5+526k6OhoxMTEsELuBjabDeXl5ZDL5UhNXfTUz2OE28VPP/3kuPjiizsGBgbWstvIvAMr3j4mLCzs+uzs7Ge+++67Y7aDuSvger0e+/fvP+aLxDJ/9Ho9urq60N/fD5lMhtjYWERERCyJHCwTxftoTCYTenp60NvbCz6fj9jYWKhUKsbeUPkTh8OByspKBAUFITs7e8Gf8ZmE28X333/vuPTSS1sHBwfXUEo1HjKbZQJWvH2ISCQ6IyUl5fVdu3aFh4SETPuaxQq4qylDQUEBJBKJp0xeFphMJnR3d6O3txcCgWBSEJZarnUpiPdUDAYDuru7oVarIRaLERcXh8jIyCVxo+UrKKWora2FyWTC6tWr510bM5dwu/joo49sN998c83g4GAppdQ7k1OWKax4+whCSHFycvLne/bskc/lFS9UwPv7+1FXV4eioiKIRF6bGLqkoJRicHAQbW1tsFgsiI2NRUxMzJKuyF9q4u3CVZPQ3d2N4eFhREdHIz4+nlH1CP6mpaUFAwMDWLNmzZxRjPkKt4sXX3zR/Jvf/GbnwMDASYscZsIyDax4+wBCSEZsbOxPZWVlivkWiMxXwLu6utDe3o6ioiI2BzgPLBYLOjo60N3djfDwcCQmJkIqlfrbLJ+wVMV7KjabDb29vejo6IBAIMCKFSuWTNrD23R3d6OtrW3Wa8lChdvFo48+anzhhRc+HRgYuJyyouMRWPH2MoSQGJVKtefbb7+NycnJWdCxcwl4S0sLBgcHUVhYyOb85kCr1aK1tRVarRYJCQmIi4tbdr+z5SDeU3H9m+t0OsTHxy/Lf/OFMjAwgNraWhQXFx8TuViscAPO6MiNN96o++STT14eGBi415M2L1dY8fYihJCwyMjIfR9++GHShg0bFnXrP5OANzc3Y2RkBIWFhewe7hmglGJoaAjNzc3gcDhISkpa1l7YchNvF1OjLVFRUUhKSlrS6RF3GR4exsGDB48QcHeE24Xdbsc555yj2bVr1+PDw8N/8aTNyxFWvL0EISQ4MjJy18svv5x9zjnnuHW7f7SANzU1QaPRoKCggBXuaaCUQq1Wo7m5GSKRCKmpqWwRH5aveLtwOBzo6upCW1sb5HI5UlJS2Lz4DIyMjKC6uhpFRUXgcDhuC7cLs9mMTZs2jdbW1t6m0+ne9ZC5yxJWvL0AIYSnUCi+ffzxx9fedNNNHklEuwTcYrEgKCgI+fn5rHAfhcPhQE9PD1pbWxEWFoaUlBS2gG8Ky128XVBK0dfXh+bmZojFYqSmprIz7adhdHQUBw4cAADk5eV5bPupTqdDSUnJSHNz80UWi+V7jyy6DGETQF5AoVC8fNNNN63xlHADAIfDgVgsRkdHBxISEljhngKlFF1dXWhpaYFSqURxcTG8MVKVZWlACEF0dDSioqIwODiIqqoqCAQCZGRksCI+heDg4MnZ756MUEgkEvzwww/hRUVF7xJCNlJK6z22+DKC9bw9THh4+O0bNmx44uOPP5Z5Mrfa1NQErVaLVatWobKy0mOtVJmMKzze2NiIiIgIpKamsrnMWWA975kZHBxEfX09JBIJ0tPTl/3N39QcN4fDQVVV1bRFbO5QU1ODE044oXNgYGAVpXTUYwsvE1jx9iA8Hm99RkbGln379oV78kPe0tIyORmMw+F4rJUqkxkeHkZdXR1CQ0ORnp7O5i7nASves+MKpzc2NkKpVCIlJSWgB814i+mK01w5cE8L+Mcff2y76aabKgYHB9dRSm0eW3gZwMZePQQhJF6pVP73q6++8qhwt7W1YXh4+IjiNA6Hg/z8fOh0OjQ2NnrsXExAp9Nhz549aGlpQV5eHlatWsUKN4tHcIXTN27ciODgYOzYsQPNzc2w25dPX5GZqsrDw8ORm5uLPXv2wGw2e+x85557Lu/WW2/NjoyM/KfHFl0msOLtAQghIQqF4rsPP/xQGRsb67F1e3t70dfXN21V+XITcKvViurqalRXVyM1NRVFRUVemUnMwsLhcJCYmIgNGzaAUort27dDrVb72yyvM9d2sPDwcGRlZWHv3r2wWq0eO+9jjz0WUlxcfGFYWNjNHlt0GcCKt5sQQohCofjoj3/8Y0JJSYnH1h0cHERzczPWrFkzY4/t5SDglFJ0dHRgx44dCAsLw7p16xAeHu5vs1iWATweD6mpqSguLkZPTw/27NkDo9Hob7O8wnz3cSsUCiQlJaG8vBwOh8Mj5yaE4L333pOpVKo/8Hi8Uo8sugxgxdtNIiMj/3jRRReVXHfddR6rlNJqtTh06BCKiormzLktZQHXaDQoKyuDTqfD+vXrERcXt2wbrLD4j+DgYBQUFCA5ORnl5eVoaGhYUqH0hTZgiYmJgVKpREVFBTxVMyUSifDNN9+EK5XKDwkhcR5ZdInDircbiMXii9LT02/6v//7P4/Fb41GIyoqKlBYWDjvitelJuAWiwXV1dU4dOgQcnNzsXLlymVZOMQSWERERGDDhg3g8XhLJpS+2M5pSUlJEIlEqKmp8ZiAx8XF4f3331cqFIpvCSFsg4Y5YMV7kRBCMiMjI//x6aefhnlqdKTZbMa+ffuwevXqBe83XSoC3tfXh7KyMoSFhaG0tJTtjMYSUHA4HCQnJ6OkpATd3d0oLy/3aAGXL3G35WlmZiZsNhuam5s9ZtO6devIE088kahQKN4nbJhtVljxXgSEEJFCofhsy5YtETKZzCNrWq1W7N27F9nZ2VjsmkwWcIvFgv3796O7uxulpaVsiJwloBEKhSgsLERMTAx27tyJ3t5ef5u0IDzRq5wQgry8PIyMjKCjo8Njtt14442C448/fr1MJrvFY4suQVjxXgQKheI/jz32WMzKlSs9sp7D4UB5eTmSkpIQGRnp1lpMFHC1Wo2ysjKoVCqsWbOGHW3KwhiioqKwbt069Pb2ory8HBaLxd8mzYknhNsFh8NBQUEBurq6PJpGeOWVV6RyufxxQki2xxZdYrDivUCkUulVRUVFJ99yyy0eU5iDBw8iMjISMTExHlmPKQJusVhQUVGBrq4ulJaWeuz9s7D4kqCgIBQWFiI6OhplZWUB7YV7Urhd8Hg8FBUVob6+HjqdziNrhoSEYMuWLXKFQvEpm/+eHla8FwAhJCUsLOwvb731lsdan7a1tcFutyM5Odkj67kIdAEfHBxEWVkZlEol622zLAmio6NRWlqK3t5eVFRUwGYLrIZh3hBuF0FBQSgoKJgcnuQJsrOz8fjjj0crlcrXPLLgEoMV73lCCBEoFIrPP/zwwwhPFVENDQ2hp6cHeXl5XsnvBqKAU0pRX1+PxsZGlJSUsN42y5JCIBCgsLAQERER2LFjB7Rarb9NAuBd4XYhFouRmZnp0T3gN954o6CkpOREiURyjUcWXEKw4j1PFArFCw888EBcQUGBR9YzGo04ePAgCgsLZ2zC4gkCScDHx8exc+dOAEBpaSnb1pRlyRIfH4+CggJUVlaira3NY9upFoMvhNuFUqlEZGQkampqPLIeIQRvvPGGTC6X/4kQkuqRRZcIrHjPg5CQkAtyc3PPueeeezyiNlarFeXl5Vi1apVPphcFgoD39/dj9+7dSE9PR0ZGBltJzrLkEYvFWL9+PXQ6HcrLyz3aUnS++FK4XaSkpMBqtaK9vd0j64nFYnz44YcRCoXic0LI8h73NgVWvOeAEJIQHh7+z/feey/ME4JDKUVFRQWSk5MRFhbmAQvnh78E3OFw4NChQ2htbUVpaSkiIiJ8dm4WFn/D5XKRl5c3Wcw2Ouq7yZf+EG7A6S2vWrUKXV1dGBoa8sia+fn5eOihh+IUCsVLHllwCcCK9ywQQjgKheLjd999N9JT/bTr6+shFovhyQEm88XXAm42m7Fr1y7w+XyUlJSwRWksy5aYmBisWbMGBw8e9Oie6Jnwl3C74HK5KCwsxMGDBzE2NuaRNe+66y5hdnb2mUKh8HSPLMhwWPGehbCwsF9cdNFFSevWrfPIer29vdDpdMjMzPTIeovBVwKu0+mwa9cupKSkIC0tjQ2Tsyx7QkJCUFpaioGBAdTU1HisqOto/C3cLoKDg7Fq1SqUl5d7pBc8IQRvvfVWmEwme4UQIvWAiYyGFe8ZIIQky2SyB5555hmPlJYbjUY0NjYiPz/f70LmbQHv6+tDRUUFCgoKoFQqPb4+CwtT4fF4KCwsBJ/Px549ezze1CVQhNtFWFgY4uLiPFbAFhUVhWeeeSZCqVT+yyMLMhhWvKdhIlz+4VtvvRXhiYpoh8OBiooK5ObmBsyADW8IOKUUDQ0NaG9vR2lpKTtvm4VlGgghSE9PR2JiInbu3Am9Xu+RdQNNuF0kJibCYrGgp6fHI+tdeeWV/JycnOOFQuEZHlmQobDiPQ1hYWG/vPjii5PWrl3rkfXq6uoQFRUVcHOoPSngNpttsqK2uLgYQUEem5DKwrIkiYqKQn5+Pvbv34/+/n631gpU4QYOF7A1NTV5ZB46IQRvvvlmWFhY2MuEEJn7FjITVryPYiJcfv/TTz/tEbdRrVZDr9d7vIOap/CEgFssFuzevRsKhQI5OTngcNiPFQvLfJBIJFi7di2am5sXXcgWyMLtgs/nIy8vDxUVFR7J9atUKlf4/N8eMI+RsFfZKbjC5W+//bZHwuXj4+Ooq6vD6tWr/Z7nng13BHx8fHyyMC0hIcFLFrKwLF0EAgFKSkqgVqvR2Ni4oIYuTBBuF2FhYYiOjkZtba1H1rviiiv4ubm5m4VC4ZkeWZBhsOI9hbCwsF9eeumlSSUlJW6v5cpzr1y5khFbpBYj4DqdbvLCoVKpvGwhC8vShcvlYs2aNRgbG8PBgwfnJeBMEm4XSUlJMBqNHplANtF9LSwsLOyl5Rg+J/5s2xdIEEKSk5OTd9fU1ER4outZXV0dOBwO0tPTPWCd73DddEgkEqSlpc34upGREVRVVaGwsJAtTAsQKKWwWCywWq2glE7+OByOyeZAhYWFIISAEAIOhwMOhwOBQODVFr0s88fV+99oNGL16tUz/rswUbhdWCwWlJWVoaSkxCMtkt9++23rL37xi8/VavV5HjCPMbDiDYAQQpRK5a7//e9/xaWlpW6vNzQ0hMbGRqxduzagw+UzMZeAq9VqNDQ0oKioiO1P7gMopTCbzTCZTDCZTEf82fV3Vx4xKCgIfD5/UpxdQk0IQW9vL6KioibFnFIKu91+xPF8Ph9CoXDyRyAQHPF3Ho/nz1/FsqG1tRVqtRpr1qw5ZocKk4XbxfDwMOrr61FaWur2NZJSihNOOGFk27Zt59tstp88ZGLAw4o3gODg4AvOPvvsV9577z2Zu2vZbDbs2LEDRUVFEImYO4Z2JgHv7u5Ge3s7ioqK2IpyL+BwOKDX66HRaKDRaKDT6WC3248Q0aMFdb6e848//ojNmzfP+DylFFardcabBJPJBJvNBpFIBJlMBplMBqlU6pP+/MuRnp4etLa2HrF7YykIt4tDhw5BKBR6pJi3ra0NJSUlbQMDA+mUUt83kfcDy168CSEipVLZVF1dHa1QKNxer7q6GhKJBImJie4b52eOFvCuri50dnaiuLiY9cA8gMPhgE6ng1arhUajgVarhcPhgFgsPkIcPfW7nku85wOlFGNjY5M2azQamM1miEQiSKXSSbtZQfcMriK2kpISOByOJSPcAGC327Fjxw4UFBQgNDTU7fUeeeQR4wsvvPD74eHhP3rAvIBn2Yu3QqF49uGHH7757rvvdvtqMzQ0hKamJpSUlDAyXD4dLgF3OByw2WwoKipihXuRUEqh1WqhVqsxODgIu90OiUQyKXqeFOrp8IR4TwelFOPj45Ni7hJ0sVgMlUoFhULBRmncoL+/H3V1dXA4HMjLy1sSwu1idHQUtbW1Hgmfj4+PIyMjY6izszOXUtrnIRMDlmV9FSaEpKalpV15xx13uC3cNpsNNTU1KC4uXjLCDTir0CMjI1FXV4fExERWuBeI3W7H0NAQ1Go1RkZGIJFIoFKpkJycHDDd9tyFEAKRSASRSITo6GgATkHX6XRQq9XYvXs3uFwuVCoVlEqlR7ys5YRUKoXVagWXy11yxaFhYWGQyWRobW11O3weHByM5557LvzGG298GcCS3z62bD1vTxepVVdXQyqVLrm9zt3d3ejo6MCaNWsmUwKzVaGzOKep9ff3Q61Ww2g0IjIyEiqVCuHh4X5tYOMtz3s+jI+PT/5OTCbTEb+TpXSz62mm5rhtNhsaGhpQUlKypCIZng6fH3fccSPbt29f8sVry1a8g4ODzz/rrLP+9f7778vcXWtwcBDNzc1LKlwOOAtm2tvbJ3Pc891GthwxmUzo7u6e3L+qVCqhVCohFosD5jPhT/Geis1mw+DgINRqNUZHRyGTyRAbG4vIyMiA+V0FAtMVp6nV6snU3FKJ3ACeDZ+3traitLS0rb+/f0kXry1L8fZkkZqrury4uHhJbZsaGBiYvMufepFgBfwwlFIMDAygo6MDJpMJcXFxiI6ODtimPIEi3lOhlGJ0dBRdXV0YGRlBdHQ04uPjl9R3aTHMVlXe29s7eVO9lPbn19bWQiAQeKT6/JFHHjH+85//fGJkZORJD5gWkCxL8VYoFH99+OGHb/FEkdpSDJdrNBpUVlZi7dq10wrRchdwk8mEjo4O9Pb2Ijw8HAkJCZDJZP42a04CUbynYrPZ0Nvbi87OTvB4PCQmJkKpVC47b3w+28Ha2towNDQ02XRnKeAKnxcWFiIkJMSttSaK1wY7Ozvzlmrx2rITb0JIXHJyckV9fX2Eu8VXGo0Ghw4d8kioJ1AwGAzYt28fioqKZv0CLUcBHx0dRWtrK4xGIxISEhATE8OoAr5AF++p6HQ6tLe3Y2RkBHFxcYiPj19SYeKZWMg+7vr6elgsFqxcuXLJXH+Gh4fR1NTkkcLfjz/+2HHrrbd+2NfXd7GHzAsoll1vc5VK9exf/vKXcHcvupRS1NTULKkvjslkQnl5OfLz8+e88/XGPPBAxOFwoKenBzt27EBTUxMSExOxYcMGJCQkMEq4mYZEIkFubi7WrVsHSil27NiBgwcPemSkZKCy0AYs6enpoJQuqe+fXC5HUFCQR3qfn3POORyFQnE8IYRZParnybISb0JIllKp3HTWWWe5/b47Ozshk8kgkUg8YZrfsVqt2Lt3L3JyciCVSud1zFIWcEopent7sW3bNmg0GuTn56OoqAhyuXzJ3KwxAT6fj5SUFGzevBkRERE4cOAAKisrMT4+7m/TPMpiOqcRQpCbmwutVov29nbvGuhDsrKy0NDQALvd7tY6hBA899xzcpVK9U8PmRZQLCvxVqlU/3j++efl7l58LRYLWltbGTd0ZCYcDgf27duHlJQURERELOjYpSjgg4OD2LFjBwYHB1FSUoLs7GxGt7pdChBCEBUVhXXr1kGlUmHv3r04dOgQLBaLv01zG3danhJCUFBQgJ6eHo94q4GAUChEXFwcmpqa3F5rw4YNSEtLyyWEFHnAtIBi2Yg3IaQkIyMjZ926dW6vVV9fj5SUlCWTgzt48CAUCsVkg42FslQEXKPRYNeuXejo6MDq1auRl5fHtvkMMAghUKlU2LhxI8RiMcrKytDY2AibzeZv0xaFJ3qVu8aJ1tfXQ6fTedhC/7BixQr09/d7JE3y3HPPyZVK5YtkiYXMloV4TzRkeeHvf/+7230FNRoN9Ho9YmNjPWGa32lra4Pdbnd7ewaTBdxgMKC8vBx1dXXIzMxEYWEh2wUswCGEID4+Hhs3bgSXy8X27dvR3t4+OR2NCXhyyEhQUBAKCgpQUVGxJKIRHA4HOTk5OHjwoNtrrVy5EqWlpQkcDudkD5gWMCwL8ebxeKdt2LAhPicnx611llqR2tDQEHp6epCXl+eR98M0ATeZTKiqqsKBAweQkJCAtWvXMmLLF8thuFwukpOTsX79eoyPj2Pbtm3o6elBoO+i8cZ0MLFYjIyMDJSXlzPqJmYm5HI5+Hy+R9IBf/7zn8MiIyOfI4QsGc1bMm9kJgghHLlc/rc//elPYe6u1dXVBalUuiSK1IxGIw4ePIjCwkKPNnpggoA7HA40NjZi9+7diIyMxPr16xEZGelvs1jcgM/nIzMzEyUlJRgeHsaOHTug1Wr9bda0eHOsp0qlQmRkJGpqajy6rr/Izs5GfX2928VrK1aswNlnn60QiUSXecg0v7PkxVskEl1x3nnnRbrbRMVms6GlpQUZGRkessx/WK1WlJeXY9WqVV7J6QaygOt0OpSVlYFSio0bNyI6OnpJRFFYnAiFQuTm5iIvLw/V1dVoaGgIKC/UF/O4U1JSYLVal0QFulAoRGxsLFpbW91e64knnpCIxeInCSFLojH8khZvQghPLBb/4fHHH3fbVW5tbV0SjSIopThw4ACSk5MRFuZ2MGJGAk3AXd52ZWUlcnNzkZ6e7tchISzeRSKRYN26dSCEeMQLN9vsqOnR4oPyLryxqx0Ox8LD8r4QbsBZD7Bq1Sp0dXVheHjYa+fxFStWrEB3d7fbuXyFQoHrr78+TCwW3+Qh0/zKkr56CYXCyy+99FLpQrc/HY3FYkFPTw8SExM9Y5gfaW1tRXBwsE8K7gJFwKd62+vXr5/3PnYWZsPhcJCWloZVq1Yt2gvXmax45KODyP7N1zjz7ztw33+r8esth/Bj48CC1vGVcLvgcrkoLCxEdXU14wvYXHUNntg6dt9994WIRKIHCSGM77C0ZMWbEMKRSqWPPfjgg24PwG1sbERycjLjhwCMjo6ip6cHWVlZPjunPwWc9baXPgN605zFaUd74fPdTvVVTR9O/PNPeGdvJy5eE4fnLl+Nb+/ZiCipEC9va5u3jb4WbhfBwcHIyMjAgQMHAr6Aby7i4uIwNDTkdnMemUyGSy65RLIUct9L9krG4/HOPPXUU8NUKpVb64yNjWF4eBhxcXEessw/WK1WVFVVIT8/3+c3If4QcJe37XA4WG97CWKy2vHQ/6pR9PutuPjFXTjQOTrr66d64ZWVlbN64WqtCTe9Xo5b3qxARKgAH9++Dn84byXOzI1GqlKM4hXhqFPP7wbAX8LtIioqCiKRyCM5Y39CCEFaWhoaGhrcXuuhhx4Si8Xi3zG98pzRxs8EIYREREQ8+eijj8rcXauhoQHp6emMLmqilKKqqgopKSl+27/sKwGnlKK5uXnS287IyGC97SVG5/AYLvjnTryztwvnrY5BvVqPJ7+sn9exEokE69evBwDs2LEDer3+iOe3VPbgpL/8hJ8aB/HgaRnYcsc65MbKJp93OCh2NA9hfcrcqTh/C7eL7Oxs9Pb2YnR09hucQEelUsFgMBzzb7aYdU455ZQwHo93podM8wuMj/vPwMbi4mLlihUr3FpEp9PBaDRCqVR6yCz/0NHRAR6P5/fGMi4Br6ioQGNjo8enkdlsNhw4cADBwcFYv349K9pLkO1Ng7j9rQoAwL+uKcQJmUqc8OcfER4y/wJiDoeD9PR0qFQq7N+/HxkZGVAqlfjrd03429YmFCWG45mLcpEgP3Y4z4EuDYYMFpyUNfs1wRvCTSnFwR4t9raNYNBgxrDBgiGDGWarA9etXzGjTa7v3b59+7Bu3TrGFt0SQpCZmYm6ujoUFbnX7fTRRx+VfvPNN08SQj6lDM0pLMmrW1RU1DNPPPFEuLvr1NXVISsri9Fet06nQ0dHB1auXOlvUwB4zwMfGxvDzp07oVQqkZOTwwr3EuRQrxY3v7Ef0bJgfP7zDTgh0ylW/TozlJKFb3mUSqUoLS1FfVMLrn9lO/62tQkXF8bizRuKpxVuAPi2th88DsHmdMWM63pauNuGjHj2u0ac8OefcPZzZXji8zr8e0cbypqHMGywoFc7jhtfL8f1r+6DzmSddo2QkBCkpqaisrKS0flvuVwOh8PhdhQhKSkJRUVFSgAbPWOZ71lynjchJH/Tpk0r3O2mNjIyAgAID3f7HsBvOBwOVFZWYvXq1QFVbOdpD3x4eBjV1dXIy8tj9L8Xy8wM6E248bVySIR8vH5dERQTYm0w22Aw2xYl3gBgsAJ/ryHY36HHFTmh+O252eDxpr/xs9kd+OJgH4qTwiENnt579ZRwD+rN+KSqF59U9qCqWwtCgJIVcty0MQknZikhDwmadCqsdgde2taKZ75uwGdVfbi8OH7aNWNiYjAwMICenh6/R+HcISsrCzU1NSgtLXVrnSeeeCJ83759zwBg5NCSJSfeUVFRT//xj390b28YnBXmmZmZnjDJbzQ1NUGlUgVkRzhPCXhHRwc6OztRUlKC4OBgD1vJ4knahozgEoJ4+cImtJmsdtz4+n6MjlnxwS1rJ4UbAGp7nYVjCQtcEwBaBg247tV96NOa8Nzlq5EjsWDXrl0oLCyc9rP0fnk3OkfG8PDp018XPCHcg3ozXvipBW/u7oDZ5kB2tASPnJ6JM/OiECWd/vPN53JwWVE8nvm6AWbb7J3IcnJysGPHDkRERDB26I5EIgGfz8fw8LBbN0grV65EampqEiEkn1Ja4UETfcKSEm9CSFp+fn7e2rVr3VrH1dCByRXKWq0WAwMD8MQUNW/hjoA7HA7U1NTAarWitLQ0oCILLEfSPmTEH7+sw9eH+pGuFOPre+YfqaSU4t4PqlDdrcELVxYgJ+bI7+SnVb0Q8jnYlLaw9ra7WoZxy5v7weMQvHtTCfLjnQ2LxGIxdu/ejVWrVh3RxMhotuGv3zWiMCEMp2Qfm1t2V7hHjBa8uK0Fr+/sgNlmx/n5sbh5YxJSlfPb6cqd8MLtczSP4fP5yM7ORlVVFYqKihibEkxNTUV9fb3baYknn3xSfsEFFzwD4ATPWOY7lpR4q1Sq3/7ud79zO8nU1NSE1NRUT5jkFxwOB6qqqrBq1aqAz/0uRsAtFgvKy8sRGRm5ZIbELFVMVjuuf20fBnRmrIyR4lCvFiarHUL+/G62nv2uCZ9V9+HB0zJwSvaR2z5dYewTMpUIEcz/Uralsgf3flCFBHkI/nPtGsSFH/ba5XI5iouLUV5ejhUrVkxuEX1le5vTK76y4JjPmzvCrRmz4JXtbfhPWRvGrHackxeNn5+QiqTIhe0KKe9wpvkU80gfKBQK9PT0oLu7m7FbYGUyGSil0Gq1bjlZa9euhUKhyCWErKCUzn/zfgAQ2Ff2BUAIkQUFBZ142mmnuXUlNxgMMJlMft3a4S6NjY2IiooKyHD5dCykiE2n02Hnzp1ISkpCamoqK9wBzl+/bUTLoBH/uDIfN21MgoMCrYPzm9H8WXUv/m9rEy4qcHqhR1PWMoxhowVn5c5/Dv1HB7pxz3uVyI8Pw4e3lh4h3C5EIhFKS0vR19eHQ4cOYUBnwovbWnBajgoFCUe2FF6scOtMVjz7XSM2PPUDnvuhGZszFPjm7o149tLVCxZuSin++l0j4sNFOC1nfn0tcnJy0NzcDJPJtKBzBRKpqake6br28MMPhysUigc8YJJPWTLiLZPJbr7rrrvE7nqazc3NSElJ8ZBVvker1WJwcNDt+dy+Zj4CPjo6iv379yM/Px/uNt9h8T77O0bx8vZWXFYUjw2pkZBPbOcaHZu7XWf36Bge+vAg8uNl+P1500dXPq3qhVjAw+b0+YXMP9zfjV+8X4WSJDle/VnRjEVnAMDj8bBmzRpwOBz86t2dsNgcuP/UI4cSLUa4HQ6KN3d3YMNTP+DZ75pQmiLHV3dvwPOX5887RH4039b2o6ZHhzuPTwGfO7/r39TwOVOrz+VyOUwmE4zGY28Gr7vuOigUCkwtXH7ssccQExODVatWYdWqVfjiiy8AAOeddx6Hx+OdTwhZeOGEH1kS4k0I4QgEgjuuv/56gTvrjI+PQ6vVMnZft6sZS15eXsCHy6djNgEfHh5GVVUViouLGRNRWM6YrHbc90EVoqTBePh0p+i1DDkvsisipt+G5cLuoPjF+1WgAP7v0tUImqb622S14+saNU7JUc0rBP/f/d24979VKE2W41/XrEFw0NzHEELAl8fhu7ZxnJgoQEL44YKxxQj3oV4tzv/nTvzq4xpkRUnw2Z3r8eJVhchQLf7z7HBQ/PW7JiTKRThvdcyCjlUoFBAIBOjt7V30+f0JIQQpKSlobm4+5rlrr70WX3311TGP33PPPaisrERlZSVOP/10AM4bmRtuuCEkJCTkaq8b7UGYd4WfBg6Hc8qpp54a6m6BWWtrK5KTkxkbim1vb0dERASjxW06AR8aGsLBgwdRXFwMkYhRN8fLlr9+14jWISOeuiAXYqHTw20ZMCAkiIso6ex52Ze2tWJv2wgeOzt72rA2APzYMAi92Yaz8uYOmX9Q3oX7/luFdckR8xZuwHkz/NtPD0EUxMONpbHYt28f7Hb7goXbYLbh8c9qcdbfd6B7dAzPXrIKb99YfEzx3WJ4r7wLdX06/PyEVPDm6XVPJSsrC42NjbBap98fHugolUpoNJpjwv8bN25c0LbR22+/XRQaGno/YdDFf0mIt0ql+s19990nc2cNi8WCwcFBREfPP38WSJjNZrS3t3u8a5k/mCrgFRUVOHToELsVjEGotSb8p6wdFxbEYn3q4V2bzQMGpChCZ705runR4i/fNuD0lSpckD+zJ/nO3k5EigVYlzy7eL5f3oX7P6zG+pQIvHJN4bwL5QDg7b2d2N40hAdOy0BBTjqUSiV27dqFXbt2zUu4KaWTw03+XdaGy4risfUXm3Hu6hiPOAiVXRo8uuUQ1qdE4JxVC/O6XQQFBWHFihUBMbZ3MRBCkJycjJaWlnm9/rnnnkNubi6uu+66Ixq9KBQKlJSUSAC4t1XJhzBevAkh8RERESnZ2dlurdPW1obExERGhpsBoLa2Funp6eDxlsYGAg6Hg7i4OKjVakRGRjJ2T+pSoWtkDLta5jcb+oWfWuBwUNx1wpE7NpoG9EhWzFyMZbLacfd7lQgTBeH35868i6CxX4+fGgdxdUnCrN7m+/u68MCEcL989cKEu2tkDL//vA7rUyJw5UTTE5VKhbGxMXA4HMhksjmPv/4153CTsJAgfHhrKX5/3kpIRZ5pTTqoN+OWN/ZDIRHg75etBpez+JuBhIQEjIyMzHvaWqARHR2NwcHBOUef3nrrrWhpaUFlZSWioqLwy1/+8ojn7733Xnl0dPRD3rTVkzBTqaYQGRl513333Rc29ytnxuFwoLe3l7HbJkZGRmAymRAVFeVvUzzG0NAQ6uvrcdxxx2FsbIyxnsFS4D9lbTjprz/hild2QzNHsZlaa8LbeztxYUHsESFvncmKfp0ZKbOI95Nf1qN5wIA/XZSHsFl6lf97RxsEPA6uKEmY8TXv7evE/R9WY0Nq5IKF2+GguO+/VeAQgqcuzAUhZDJUXlBQgISEBJSXl087lcxic+AfPzbjpL/+hD2tw/j1mVn49I51k/vIPYHV7sDtb1VAM27Bi1cVzPq7mg+EEOTk5KCmpoaRxWscDgfx8fHo7Oyc9XVKpRJcLhccDgc33ngj9u7de8Tz69atg1AoLCaEMKJNI6PFmxDC53K5V1x44YVudejo7e2FQqFgZKMPSilqamqW1H7nkZER1NTUoLi4GMHBwX6bB84CNKj1+O2ntUgID4GDOkO1s+Hyum8/7sgdGw1q5ySoVMX0FdXbGgfx6s52XFuaiI2zNFwZNpjxvwM9uKAgdsZhJO/u7cQDHx7EprRIvHRVwYKEGwBe39WO3a0j+PWZmYiRBR+T405MTERkZOQxAr6/YxRn/G07nv6qAcelK/DdLzfh+vUrFpWLno3ff16Hve0jeOqCXGRHe6aRVFhYGEJCQhhbvBYXF4eurq5Zbz76+vom//zRRx/h6BbahBDceeedEplMdqPXDPUgjBZvHo931gUXXBDsbki1vb0d7k4g8xeuIjV/jfr0NBqNBtXV1SguLp4MlftjHjiLk/+UtUHI5+CVawpByOziPZPXDQB7Wp0h96P3SQPAqNGCez+oQqoiFA+elnHM81N5c3cnLDYHrls3/ff14wM9ePB/B7E5PRIvLkK424eMePKremxOj8TFhXEzFqclJSUhLCwMFRUVMFvtePqrelz0wk6MWez497WF+OeVBTO2M3WHD/d349Wd7bh+/YpF57lnIjMzk7HFa3w+H3K5HP39/QCAyy67DGvXrkVDQwNiY2Pxr3/9C/fffz9WrlyJ3Nxc/PDDD/jrX/96zDrXXnutQCAQ3M6EwjVGJ0gVCsV9d9xxh1ul1RqNBkFBQYysYrZarWhvb8eGDRv8bYpHGBsbw4EDB1BUVHRMcZq3x4myHMuI0YKPDvTg/HynGMeFidAyS4OVmbxuAChrHkZWlGRab/lXH9dgdMyC//xszaxia7La8cbudhyXHjlt+P272n788oMqrE2S44UrFy7cdoezFSufy8GT5+fCbDbPWlWempqKr3ZV4bS/bEXrqBWXFMbh12dlIXQB3d4Wws7mITz80UGUJIXjoTluchZDUFAQEhMT0dLSgowMz6/vbVasWIGamhqoVCq88847xzx//fXXz7mGTCbD2rVrRR9//PFaADu9YKbHYKznTQiJEIvFye5+yNra2hjrdTc3N2PFihVLokjNZrOhvLwceXl5CAmZfh8w64H7lnf2dsJsc+Bn6xIBOMVNMMPErdm8bpPVjv2doyidpjL8qxo1Pj/Yh7tPTJszBPxJVS+GDBbcsOHYbmu7WoZx29sVyImW4OUFVpW7+E9ZG8o7RvHYWdmQCTCrcDscFK9sb8XPP+/F0JgNfzg1Hk9dmOtV4b7utX1IlIfgH1cUeDwU7yIhIQFqtRpms9kr63sTsVgMSikMBoNb69x6663yqKio2z1kltdgrHiHhoZeccMNN7jldVutVmi1WkREuD2EzOeYTCb09/cjPn768X9MglKKiooKrFixYs69mayA+war3YHXd7VjQ2oE0iY6f41b7QieQRRn87orOkZhsTlQmnKkCGrHrfjNlhpkRklw0zTtT6dCKcW/d7QhQyU+5iaguluDG17bh4RwEV79WdGiBLR5wIBnvm7AiZlKnJ4ln1W4u0fHcPkru/HE53XYlBaJb+/ZhFjOKIaH51eNv1DKJoQ7ITwEb99YPGOu3xNwOBykpqaioaHBa+fwJomJiejo6HBrjeOPPx4ATiKEeGZrgJdgrHiLxeJbrrzySrc6qnV1dSE2NpaRhV6u0DFTt7ZNpaGhASEhIfOu9mcF3Pvsbh1Gv86MK4oPV3SPWWzTNjjRjFnw7r5OnLs6ZtqmKmUtQ+ByCIpWHCmET31VjyGDGU9dsHLOtp5lzcOoV+tx/foVR3xfmwf0uObfexEWEoQ3ri9eVOW1yWrHne8cgCiIi0dPT8WePXumFW5KKf67vxunPbsdB7u1ePrCXLx0VQGUshCsWbMG1dXVGBsbW/D5Z6OseQjXTxFueahbl7x5ER0dDY1GM23b0UBHpVJhYGAAdvvso1Fng8fj4ZxzzgkCcIrnLPM8jIy3EkKSioqK5O70t6aUoqurCyUlJR60zDcYjUZotVqsXLnS36a4TW9vLzQaDYqLixd0HJsDd+JwOGA2m2E2m2EymSZ/zGYzLBYLKKWglMLhcMBoNGLPnj0ghIDL5UIgEEAoFEIoFB7xZx6Ph611AxDwDo/adDgoTFbHtJ73++VdMFlnLiLb2TKMvFjpER7xntZhvL2nEzduWIHcWNmc7/PFbS2IFAtw9qrDTZS6RsZw5St7weVw8Ob1xVDN0bltJp78sh51fTr847JctBw6MK1wDxvMePijg/j6UD+KEsPx54vzjrhRCQ4ORl5eHsrLy1FaWuqRVFZZ8xCue3UfVkSE4K0bfCPcgLPqOiMjA/X19SgoKHB7PUopho0W9OtMsNrp5GN08nmAQwCVVAilWAiOG3vWORwOVCoV+vr6EBsbu+h1brjhBumnn356F4DPFr2Il2GkeMvl8htuu+02t/bijY6OIjQ0FAKBb74QnqSurg4ZGRmMjBhMRavVorGxEevWrVvUe1luAm6xWKDRaKDRaKDVamEwGMDhcI4RYalUCqFQCD6fDw6HA0IICCHYuXMncnNz0dSvh1TIhYjrmBR7V4tJk8kEq9WKzw6MI1chwKC6BzKZDJwgpzAe7XnbHRSv7exA8YpwZEUfm8XSm6yo7tbi1k2HB+WYrHY89L+DiAsPxj0nzf1vVtOjdXY6OzUDAp7z/AN6E6761x6MWWx47+a1SJyjX/pMfFvbj1d3tuOq4liIRlumFe6tdf144MOD0I1b8dBpGbhhQ9K0TVHCw8ORmJiIyspKFBQcOzp0IexocnrcvhZuFwqFAs3NzfMauUkpRffoONqGjOjVjKNXa3L+XzOOvok/m23H7omfDgGPg7hwERLCRYiXi5CuFOP4TAUU4vnfmCUkJODAgQNuiXd+fj74fH4eIURCKQ3I7jWME29CCFEqlVeef/75bm3K7uzsZGS+WKvVwmq1IjJyfpOUAhWz2YyKigqsWbMGfP7iU0tLWcANBgPUajVGRkZgNBoRFBQEqVQKmUyGqKgohIbO3mr0aDgcDkbNwCX/qkCaSoyPbyud9vgGtR5DX2/DTeujYLFY0NTUhJ5h5/VrdFCNvj4RIiMjwePx8F1dP3o04/j1mZnTnnNv2wjsDnpEvvvv3zehdciIN64vgiho7kvQCz+1QCzg4YoS5/dVO2bF1f/ai36dGW/eUIzMqMWVvvRpx3Hff6uQpRJjo3QEK1fmHiHcYxZnT/J39nYhQyXGG9cXzXmu+Pj4yXROenr6ouyaKtxv31ji1Rz3bGRmZqKuru6Y6OSI0YKqbg2quiZ+urUYMR5u3sMhgFIiRJRUiOxoCU7OUiJKKoRKKpwcMkNAMPEfAMBBKXo1JnSOjKFj2IiO4THsbBnGuNUOQoDChDCcmhOFU3NUiJHNvgVPJBKBEAKj0Thj8etcEEJw9dVXhzz99NMXAvj3ohbxMowTbwAFBQUFwWLx4sbnAYDdbsfo6Cjy8vI8aJZvaGhoWPRFIVBwOBwoLy9Hdna2R/anLxUBp5RiZGQEarUag4ODCA4OhlKpRFZWFkJCQjwSafn9F3UYt9pR1aXB9/UDOCHz2Al639U598qesyYJSonT45H064HvtyFOKcfIyAgaGhoQHByMF/aMI1oqxInTrAMA25uGIOBxJjuM1fbq8OJPrbggPxYbUue+Ae0cHsMXB/tw48YkSIR8jFlsuO61fWgZNODf166Zdt/4fLA7KO5+txIWqwNXp9qQn5d3hHDX9Gjx83cOoG3YiJs3JeEXJ6VNev1zkZWVhT179qC3t3fBsxK2Nw3ihtfK/S7cgLNxi9kObK3uQJvWjsouDaq6NegaGQcAEAKkKcQ4MVOBvDgZUhViRMuEUEqE8x5NOhuUUjT2G/BlTR++qlHj8c9qnQNe8qLx8OkZs+6jj42NRXd3t1vXymuuuUb08ssv3wlWvD2DSqW67bbbbnOrPLy/vx8KhYJxYWe9Xg+r1bqgaTmBSHV1NVQqFRQKhcfWZKqAU0rR39+P3t5eaLVahIWFQaVSeaVPvcbkwOfVfbhlUzI+PtCDDyu6pxXv7+sHsDJGOincANCrdU5tSo9TIDvR+fmraFXjQO9+XJopxK6dZVAoFIiLizuiZ8K2pkEUJ8kh5HNhd1A8+L9qSIP5+NUZ03vqR/Py9lbwOBxct24Fxi123Ph6OQ50juL5y/PnJf4z8dz3zdjTNoIbVgbhpJLDwu1wUPxrRxue/roe8hAB3rqhGKXJC7vccDgcFBYWoqysDKGhofOe8ve/im488GE1kiND/SbclFI09OvxQ/0gfqgfwP6OEdips4o+RhaMVXEyXFmcgLw4Gf7x+4fw1b8+Ro9CgadragA4uyOefuo5aG9vR2JiIt5//32EhS3uBosQgnSVGOkqMe4+MQ1tQ0b8d38XXt7ehu9q+3HH8Sm4fv2KabcFRkVFoaysDGlpaYu+ziclJUEul8cQQqIppQHXeo5R4k0I4alUqjNOPvlkt9Zx947MXzQ1NSE1NXXuFwYwPT09sNlsSEqafWvQYmCSgI+NjaGzsxN9fX2IiIjAihUrIJPJvHpD2aZz5h1PzFSgtk836UFNZdhgRkXn6DFDRfo0ztdOHef5QeUAhHwOHrxoA0Q8QK1W48CBA+BwOEhISIBNIEHroHGyYv0/ZW2o7tbi75etnldV+JDBjPfLu3De6hhIg/m44bVy7GwZxp8vysNpKxffx39v2wj+b2sj1kZxccspqyeFe0Bvwi/fr8L2piGcnKXEUxfkLrpvOJ/PR0FBASoqKrB+/fpZWy9TSvHsd034v61NKE2W459XFkAa7LtdSmMWG3Y2D+P7hgH8WD8weaOWFSXBzZuSIRpT44y1K7Ei6shaAMs1V+Lu22/B1VcfHoP95JNP4oQTTsCDDz6IJ598Ek8++SSeeuopj9i5IiIE952SgUvXxOOJz2vxzNcN+KpGjTevLz5m4Aufz0doaCi0Wu2cQ2Rm4+abb5b96le/ugbAH92z3vMwSrwBrD/xxBP57uRIrVYrxsfH5yzCCDTGxsZgNBoZnes2m81uFajNh0AWcEophoeH0draCovFgoSEBGzcuNFnPfXbtQ5wCJAVLUGMLBg1PdpjXvNDwyAoxTFh8F6tCWQilwk4t4d9dKAH562OgUzkFLj4+HjEx8dDr9ejo6MD71c4vbG1iVJ0DBvx528acUKGAmfmzk94X9vZDovdgatLE3DDa+UoaxnCny7Mw/n5iy9E0oxZcNc7FYgIJnjq4vzJHg/f1/fjvg+qYbTY8PvzcnB5Ubzbn1GxWIzY2Fg0NDQgKytr2teYbXY8+OFBfHSgBxcVxOL3562czAt7k45hI76vH8APDYPY3ToMi82BkCAu1qdG4OcnpGJzumKyen9gIBzd3R3HiPfGjRvR3t5+xGNbtmzBjz/+CAC45pprsHnzZo+Jt4u4cBFevKoQX9Wo8fN3DuCqf+/BG9cXH3PD4wqduyPel156Kf8Pf/gDK97uolKprrr88svdGs+zmDxUINDc3Izk5GTGhfqnUl1djczMTAQFeTccGGgCTilFX18fmpqaIBaLkZaW5tYFZbG06xxIUYRCFMRDbFgwRowWjFlsRxSN/dAwAKVEgOyjKsfV2nFEhgomc5nv7XNuD7umNPGY84jFYuTk5OBvFeNQho5C3VSFP1XYwOUAj5+bM6/PsMlqx+u7OrAxNRJPflmPspYhPHNhHi4oWLxwA8AD/63CgN6MV6/MQXyUAiarHU9+WY9Xd7YjQyXGu5eVIFW5+Hqao0lKSkJZWRlGR0ePCR+PGi24+Y392Ns+gvtOScdtm733/bY7KPa0DWNr3QB+qB9A65BzD3dSZAiuKknA8RkKFCaGTZvXj4yMRENDA8bGxuZsI93f3z853TAqKgoDAwOefzMTnJqjwj+vzMctb+7H1f/ei/dvLjnCfoVCgbq6OlBKF/17jYiIQERERDghREkp7feU7Z6AMeI9UWV+ynHHHefWOt3d3cjPz/eQVb7BbDZjZGSE0fu6e3p6Jvdg+oJAEfDBwUHU19dDIpFM27PdV1BK0aZ14KQcGQBMVuz2asaRMjHpi1KKPa3D2JAaeczFrk9rQtTEMTa7A6/v6sDaJDkyVNPnc612B3a2jOCsvCi08kJRO1iHG1YKMNrTighR6pw3cAM6M7TjVlR0jsJgtuHpC3JxoZvC/XV1N76uHcBt6+OwITsBjf16/PydA6hX6/GzdYl44NSMRbVVnQ1CCFatWoX9+/cfET5vHzLiZ6/uQ49mHH+7bDXOzvOOQ1Hbq8PHlT3YUtmDfp0ZQTwOSpLkuHptAo7LUCBBPnc1NiEEycnJaGlpCbhr0AmZSjx7yWrc/nYF3tjVcUTrXA6HA7lcjsHBQbfqay677DJxY2Pj2QBe9oDJHoMx4g0gIysrK8idCWLj4868nb8uoIultbUVSUlJjPW6p4bLfYk/BVyj0aCurg58Ph+rV6/2+9S3fp0ZOgvFyhin2MaEOb8D3aOHxbtl0IghgwXFK44tiOzVjE+2Sf2ubgA9mnH85qzpQ8EAcKBTA4PZhtgwEZ75ugEnZCjw8GUF6O7uRllZGWJjY2fty8/jOj/rBrMNT12Qi4sK59d9bya0eiN+/XE1EsKFuOvUbLyxuwNPfFaLUAEP/7l2DY7L8Fzx5NGEhoYeET4vbx/Bja+XAwDevqEYhYmeLUDt045jS2UvPj7Qg3q1HjwOweZ0BX59ZjSOz1DMa3ve0URFRaGxsRFms3nW3hhKpRJ9fX2IiopCX1+fR4tSZ+KM3Ci8Vx6Jv3/fjAsLYifTOIAzdN7e3u6WHeeff77w+eefvwYBJt6M6a0ZFhZ28VVXXeXWp7y7u9utjfv+wG63Q61WM87uqfgqXD4dvm6lajQaUV5ejrq6OmRmZqKwsNDvwg0ATQPOedrpE56yfKIQSzt+ePzjnjZnVXFJ0rFtQfu0pskc6Ks72xAjC55xexjgnM8NOCuohXwu/nj+SnA4HMTHx2Pjxo0ghGD79u1ob2+fdgZzeEgQ8uJkePqCXFzspnCbTCY89l4ZBsYoTsmJxoX/3IVff1yD4iQ5vrx7g1eF20VSUhJGR0fxdlkjLn95D8JEQfjotnUeE26dyYr393Xhspd2o/TJ7/Hkl/UQBXHx+Lk52PvIiXjlmkKcmRu9KOEGnN73ihUrjslxH83ZZ5+N1157DQDw2muv4ZxzzlnU+RbKQ6dlQGey4uXtrUc8HhYWBp1OB5vNtui109PTweVy0wghAeX1MUa8hULhpWeeeaZbMS3XHSGT6O7uRnR0NGN7mPf09IDL5fosXD4dvhBwSilaW1uxb98+JCQkYO3atX7Ja89E30QFsStcLpgID5ush3tA724dgVIiQIL8yLymzmTDmMWOaGkw6vp02N06gqvXJkzbZcyFWuc8X8ugEb87JxuKKdvOuFwuUlJSsH79ehiNRpSVlR0zCUrI52LL7evc9rhNJhM+/X4nvmh3Xrxf2taKEaMFz1yYi1evXbOgzl3uslsvw8OfNmFVnBT/u6100V3hXFhsDnxX24/b367Amie+w/0fVkOtM+HuE9Lw032b8b/b1uGqkgSPbTmLiYlBb28vHA7nroXpZmY/+OCD+Pbbb5Gamopvv/0WDz74oEfOPReZURKUJsuxte7IHDshBEql0u3c+5lnnhkE4Hi3FvEwjAibE0Iis7Oz5e5UWo+Pj4PL5frF+1sslFJ0dHSgqKjI36YsCn+Fy6fDmyF0o9GIyspKyGQybNiwwWfV41Ox2R247a0KZEVLcNvmlGMqlvs0TjFVSp0hT+HE8yar80LsyneXJMmnyXdPbBOTCfHaznYI+RxcsmZ2UV0xIUzrUuQz5nP5fD6ys7MxOjqK/fv3IzY21qPpIZPJhN27d+PHoRCYbeNQiAW48/gUXLIm3icV3S4sNgce/ugg/ru/GyemSHBHkeSI0O5CoJSiolODjw/04LPqXoyOWREeEoTLiuJx7uoY5MVKvZZe4/F4UCgUUKvViI6OnnZmNgBs3brVK+efi7VJcvzpm0aMGC1H3LCoVCq0tbW5Vah8ySWXSLds2XINgM89YKpHYIR4BwUFnXnZZZe5FXtUq9V+9f4Ww+joKEJCQuBOnt9fUEpRVVXlt3D5dHhawCmlaGtrQ2dnJ3Jzc/3aPKd1yIhvavvxTa2zIPbuE498b2rdOCRBZLIa11WYZbY5Pe/24TEM6M0oTjr2PbiEP5jPxUcHenB+fuyc4nN2XjQa1Ho8elbWnGISFhaG9evXo6GhAWVlZVi1apXbqQaXcK9cuRJamQkFSQpcVhQ/7VQ0bzKgM+G2typQ3jGKu09Mxc+PT8GuXbswMjKyoM/LqNGCD/Z34e09nWgfHoOAx8HJ2SqctzoaG1IjPdLRbD64ercH4o6dtROjYve2jeDUnMPXeplMBq1W61bVeWlpKex2+wZCCIdSOr9G7V6GEeIdGRn5s/PPP9+tfEN/f3/AVUrOhatLERPp6+sDj8cLuBsmTwm4y9uWSqV+87an4tqzHREqwOfVfceId6/GhHDh4QuX4CjPe3fr9Plu4HDIfVvjIMw2B66dZnvY0cSFi/C3y1bP234ul4usrCyMjo6ivLwccXFxi/bCpwq3XC7HOce+JZ+wr30Et71VAaPZhr9fthpnTUQg8vLysH//fmzYsGHOdFhVlwZv7O7Ap1W9MNscKEoMx+3HpeDUHBXEQt+Pmw4JCQGXy4Ver4c7Laq9gWvK26DedMTjhBCEhYVhZGRk2vns84HH46GoqIj3ySef5AMod9dWTxDwiVRCiJDL5aZnZGQseg2bzQaTybToJvX+wGKxQK/XM7IVqsPhQGNj44yNKfyNuznwwcFB7N27F5mZmcjJyfG7cANATY8OQj4Ht2xKQtOAAW1DR85iVmuPFG8elwMeh0zmvPe0DiMiVICkafKwrrD5V4fUKE2WI13lvYt2WFgYNmzYgPHxcezbt2/BhUZHC7c/oJTiP2VtuOyl3QgV8PDx7esmhRtwVp9HRkaiq6tr2uNNVjs+KO/C2c/twDnPl+GLg324qDAWX929Ae/fshYXFcb5RbhdJCYmzlm45g9cI2eNlmNneatUKqjVarfWv+KKK+QRERGXurWIBwl48Qaw4ZRTTuG7k8dxd5+fP+jq6kJcXBwjt4d1dnZCqVQGdLh/MQLuKkpraGjA2rVrA+rGSjtuhTxEMBku/Lb2yAtVr3YcYcIjP0tCPhdmm8OZ724bQXFS+LSft96JsHm/zjwvr9tduFwucnJyEBUVhZ07d2JsbGxexwWCcI9b7Lj7vUr89tNabE5XYMsd6ya32E0lNTUVra2tsNsPC03HsBF/+KIOJX/civv+W40xix2/Oycbex4+AU+cu3LGPfW+RqlUYmho6AjbA4FgPheEAEbzsTd8ERERGBwcdGv9U089lfB4vHPdWsSDBHzYXKFQnHnmmWe61VVNrVYjISHBUyZ5HUopuru7UVpa6m9TFozNZkNbWxvWr1/vb1PmZCEhdIfDgerqalBKsXbt2oDwtqcjNkyEzCgJvqsbwE0bnTO0DWYb9CYb5Ed5awIeByarHZ0jY+jTmqYNmQPA3nZnSD0pImTaQSbeIi4uDiEhIdizZw9yc3NnFeRAEO6OYSNufmM/Gvr1uPfkNNy2OQWcGSry+Xw+4uLi0NTcgh6HFK/v6sBPjYPgcghOyVbiqpJElMxwM+VvCCGIiopCb28v4uLc2w3gSSx2BygFgqbJ//N4PAQHB8NgMCy6nkIikUAmk0kIIVJK6bG9hX1MwHveXC73lA0bNiz6eEopNBrNoifb+AOdTgeRSOTWnGt/0drairi4OMbYPh8P3Gw2Y9euXRCLxVi1alXACreL0mQ5qrs1sNmd+Wz1RNg7THjk113I58JkdWBP6wgAoGSa5iwAcGZuNE5fqcJ7N6+ddXuYNwgPD0dJSQlqa2tnDNUGgnB/X9+Ps/6+A2qdCa/+rAh3HJ86o3ADzgEwX3dSXPFOM65/rRx1fTrcfWIqdj54PP5xRQHWJh9b9R9IuPqGBxKDejMAQCGZvomMJ0Lnp5xyigDA4gXJgwS0500ICUlJSQl3R3hHR0e9Pq3J0zCxmQzgzNP39PRg48aN/jZlQczmgev1euzfvx9ZWVkBnXqhONzoJCdGApPVgdYhI9KUYgzqLQAAqeDI74CAx4HZZsfutmHIQ4KQopjeI3ng1MXXm3iC4OBglJaWorKyEnq9Hjk5h/uj+1u4HQ6Kv33fhGe/a0JWlAQvXlUwWTh1NK5tXm/u7sDn1X2w2B1YHR2Cm1ND8LOTC3xWMe4JQkNDYbPZ5uy45ktc4h0pnt4epVKJ/fv3IyUlZdHnOO200yTvv//+OQA+W/QiHiKgxRtA6QknnOCWC8e0LWKUUgwMDCAzc37zjgOJpqYmJCcnB7xnOh3TCbhOp8P+/ftRUFAw75nM/sR1f7oyxjkx72C3FmlKMTRjTvEOPeqbJJjwvA90alC0IjBDtC64XC7y8/NRW1uLqqoq5OXlwWw2+1W4tWNW3PN+Jb6vH8AF+bH4/Xk50/ZGt9gc+LSqF/8ua8OhXh1CBTxcVhSHK0sSkBwZgu3bt8NmMYM/TdvmxMREiMVicLlc8Hg8lJcHRKEzAGfTlp6eHq+M910MXaPOCJNKMv3GJKFQCIfDAYvFsujtq+vWrQOlNCCatQS0eCuVyrNOP/10mTtrDA8PM2oG9tDQEMLDwxnXUW18fBxDQ0MBW2E+H6YKeFVVFUZGRrBmzZqAaG+6EFZEhCKYz0VNrxYXFMRidMzZAjU06OiCNQ5aBw3o0Yzjpo2BcQGeDUIIsrKy0NDQgH379sFoNM6ZC/cWdX063PzGfvRpx/H4uTm4svjYEaKjRgve2tOB13Z1YFBvRqoiFE+cm4NzV8dMVkYDQFpaGhobG5GXlzftuX744YfJ0aWBRExMDPbu3Rsw4n2oVws+l8wYQQKchWvDw8OL7rQZGhqKsLAwMSFERinVLNJUjxDQCsHhcE52JwRrtVpBKWVM/hVgbsi8oaEBaWlpAe29zQcOh4Pk5GT09PQgMjIyoISbUopPqnqxrXHwmH7gXEJgszsf43IIsqIlk3u/Ryc972PD5q7RkNM1ZwlECCFITEzEyMgIBAKBXyr+t1T24Lx/lMFss+Pdm9biqpKEIz73LYMGPPLRQax9civ+9E0jMqMkeO26Inxzz0ZcWZJwhHADzlysTqc7pkVsoCMQCMDn8wPG7tpeHdKU4lm758nlcgwNDbl1nom8t98rcgNWvAkhIaGhoeHu9Id2Z1O+P7Db7dBoNAG1BWk+6PV6GAwGRqUnZkKn06GyshIbNmyAyWTyyTCT+dI0YMDP3zmAq/+9F18fOrLwRhLMh950eMjIyhgpDvXq4HBQaMYsEPI5COIeu1UMAGQiPtIUgdVwYyZcOe41a9YgPDwclZWV0w428QYOB8VTX9XjrncrkRsjw6d3rkdBgrMeh1KKnc1DuO7VfTjhzz/hg/3dOCcvBl/fvRGvX1eETWnHjll1QQhBZmYm6urqpn3u5JNPRkFBAV566SWvvr/FECiFa5RSHOrVHTOH/mjCw8MxMjLi1rlOPfVUSVRUlG8mrsxCIIfNS0888US3XOahoSFGiXd/fz+USiXjvNf6+npkZGQwzu6jMRqN2L9/PwoLCyEWiwNiHvhU7I7DIvXXb5twcpZqsqJZGsyH0WKH1e4An8tBdrQEYxY7WoeMGB2zImyadqaSia1jxSvCZ62MDhSOLk4LDw9HbW0tampqvN490WC24e53D+C7ugFcXhyPx87KRhCPM5nPfmVHG+r6dJCHBOHuE1NxZUkCIkLnX8gVERGB5uZmaDSaIwbalJWVITo6GgMDAzjppJOQkZERUAWhKpUKO3bsQHp6ul+//2qdCSNGC3Im6j1mgsfjgcvlLom8d8B63kql8qzTTjtN5s4aTPO8mTj1bGxsDGazOSBzcgvBarWivLwcq1evnmz76OtxonPBmxDYk7KUaOjX44uavsnnJELnfbje5GxQ4Wro0divh2bMMm0v8kfOyMTvz8vBg6cFfnHkdFXlrhy41WpFW1ub187dOTyG8/9Rhh8aBvG7c7Lx+3NzYDTb8Nz3TVj31Pf45QdVsDscePqCXJQ9eDzuPjFtQcLtIiUlBa2tR460dPUQVygUOO+887B3716PvCdPwePxIBKJ/B46P9SjA4A5PW/AGTofHh5e9Lmm5r0XvYgHCFjx5nA4J7uzv5tp+W6HwwGdThdQYyTnQ1tbG1asWOFvM9yCUoqKigqkpKQc8/sPJAF3hbmPS1cgVRGK575vngwZS0XOz7luYj53iiIUhABN/YYJz/vY74FSIsQVxQmTE8ACldm2gxFCkJeXh56eHrc7aE3HrpZhnPP8DvTrzHj9uiKsS4nArz6uOSKf/fp1Rfj67o24eE3ctNXm80Uul8NgMMBkcna0MxqN0Ov1k3/+5ptvkJOT45H35Uk8sX/aXQ716kAI5tWFLiIiYknkvQNSvAkh/KCgoGWV7x4dHUVYWBijQs82mw0DAwOMixYcTV1dHSQSCWJiYqZ9PlAEPFoWjCAeB21DBlxdmoh6tR6Hep0ehysErp0Q7+AgLuLCRGgc0GN0zDJt2JwJzGcfN5fLRWFhIWpqamA0Gqd9zWJ4c3cHrvrXHphtDpydF41//Ng8mc8+d1UMvrnHmc/eOEs+eyG4ivFczWj6+/uxfv165OXloaioCGeccQZOPfVUt8/jaZRKJfr7+/1qQ2O/HvHhIoQI5s4Eu4aUuMOmTZskCoXiOLcWcZOAFG8Ame7msIaGhhgVymXafnTA2X89JiaGcdvaptLd3Q2DwYC5Bt8EgoBzOQSpilA09htwVm4U+FyCjw70AHAWrAGAbkrRWpoyFE39eowYLQgLYUYEaioLacAiFAqxevVqlJeXw2q1zvraubDaHfj1xzX41cc1sDkoxix2vLG7A/06M+45MQ07HzweT16QO23PcneJiYlBX18f7HY7kpKSUFVVhaqqKhw6dAiPPPKIx8/nCQQCASilMJvNfrOhaUCP1HkWXU7Ney+WgoICBAUF+bXTWkBedblcbuGmTZtmrzyYg+HhYUZVbQ8ODiIyMtLfZswbSik6OjoYO7IUcEY7WlpakJ+fPy/PKRAEPF0pRmO/HjJREE7IUGJLZQ9sdgekLvEePzyUIUUhRmO/AZoxK1SSwB0SMx2L6Zwmk8mQnJyMiooKtyrQtzUO4o3dHRAFcXFiphJPnJuD7fcfh+9+sQl3nZi6qHz2fOFyuVCpVOjp6fHaObyBP71vm92BtiHjrPu7j8bdvHdsbCwcDodfG7sHpHirVKoTioqKFh3ns9lsjMp3GwwGBAcHM6oz2fDwMCQSyaIrNv2N1WpFZWUlCgsLwePNf9OFvwU8TSVGn9YE7bgV5+fHYMhgwfamoWPC5oDT83ahYJB4u9PyNDY2FmKxGE1NTYs+/4bUSHx6x3oc+M1JeOWaQlxZkjBjy1NvkJiYiM7OTp+dzxP4M+89bLTAaqeIDZu+s9p0uCvehBDExcVxCCHRc7/aOwSkeDscjjWrV69e9PE6nQ5SqVuOu09hYsic6V53TU0NUlNTFzXj3Z8Cnj4Rqm3q12NzugJhIj7+d6BnMiw+YjwcupxavMMUz9sTvcozMjKgVquh0+kWdXwQj4OVsVIIeP65mQ4ODgaPx1u0/f5ALBZjbGzML2NCDRMjQMXC+d+Ey2QyaLXuDQbbsGFDMIACtxZxg4ATb0IIj8/ny1zbdRbD0XslA52BgYGAHnpxNBaLBQaDgVGT2qbS398Pq9U6Y4HafPCFgD/7XSOe/qr+iOYrGVHO70VNjxZBPA7OzI3GN4fUsDsoJEIehgyH83iZUYe/Q0oGiLenhoxwOBzk5eWhsrISDofDgxb6joSEBHR0dPjbjHlDCEF4eDhGR0d9fm7X/O6QoPmLd1BQECwWi1vpldLSUnFkZKTfNt0HnHjDA8VqWq2WMeLtapQfPM1QgkClq6sLsbGxjKqMd2G1WlFbW4vc3Fy37fe2gD/7XRP+8WMLbn5jPxwTDVqipMGIkgpR3uG8SJ6+MgpmmwM/NgwiIlSAQcNhz5sQggvyna12o2SBLd6eng4mlUqhVCrdCp/7E6VSiaGhIb94sovFE1uwFgOB83tscyxMiENCQtzanVBQUACBQMCKtwsul1u4ceNGt2LeOp0O7njuvoRpUQJKKbq7uxEX59dajUVz8OBBpKamQij0jJh5U8BzY51fg50tw3hrz2EvrCAhDPsnxLtoRTjkIUH4qkaNiFABhvRHVvw+fWEutt9/3GROPBDx1ljP1NRU9Pf3Myr87ILD4UClUqGvr2/uFwcInugbvhhcue7u0bEFHSeTyaDRaBZ93ri4ONjt9vhFL+AmASfe7har2e12EEIYU/zFtC1tOp0OIpGIkYVqarUaNpvNrXD5dHhLwNOUYoiFPKxNkuOPX9ajc9h5cVqTGI4+rQk9mnFwOQQnZSnxff0AJME8DBmOFG8uh/i02GqheHMeN4fDwapVqxgbPo+LiwuIvuHzRSAQwGaz+TxaIBPxIRbw0DWycPF2J+9NCEFsbCyHEOKXRhcBJ97uFqtptVpGzF52wbT+62q1mpFNWex2O+rq6jwSLp8Obwj4Kdkq6E02nLs6GlxC8JtPagBgchhGebuz0cSpOSoYzDbU9ekxbFz83lVf403hdiGRSKBUKo9pO+riq6++Qnp6OlJSUvDkk096xYbFEhoaCpPJBJvNNveLAwRPDP5YKIQQZEZLsLNleEE5bKlU6pbnDQDr16/3W9FaQIk3IYRwuVyZO5XibL7buzCtuM5Fa2sr4uLiPBYunw5PC/iG1AiIgrio7tbijuNT8GPDIHa2DCFDJUZIEHcydF6aHAGxkIcezTg0Y1ZY7YHvZfpCuF2kpKSgq6vrmOYtdrsdt99+O7788kvU1tbinXfeQW1trVdtWSgKhcIrbV+9hWtetq85d1UMmgYMqOmZf4pEIBC4XbRWUlIiDg8PL170Am4QUOINQOVuSFOj0TBmmxjT8t0mkwkcDodxIXOLxYLu7m6f9GD3pIAL+Vwcl67A14f6cfXaRERLhXjyy3pwCMHq+DDsa3eKdxCPgxMzlZPHjQa49+1L4QacjU+SkpKOKV7bu3cvUlJSkJSUhKCgIFx66aXYsmWL1+1ZCIHQN3wh+CvvfcbKKARxOXivfGH7490tWktLS0NoaOjiQ8VuEGjinZaTk+NWZY1Op2NM2Jxp+W61Wg2lUjn3CwOMpqYmJCUl+awOwpMCfmqOCkMGM/Z3jOIXJ6ejuluLzw/2oSAhDA1q3WQ71HNWHe4VYQlgz9vXwu0iLi4OAwMDGB8fn3ysp6fniMLL2NjYgOtsFhYWBo1G47OZ5e7ir7y3VMTHBQUxeHdvF+rV8/e+pVKpW3nvlJQUWK1Wv8wLDijx5vF46Xl5eYsuE6eUglLKmGI1pnne/f39jGsmMz4+jsHBQZ9Xx3tKwE/KUkIazMc7+zpx3uoYZKjEeObrBuQnhMFBgX1tzvzihtTDrXUD9TrvL+EGnP8eaWlpR/xbTCeIgbb9kRACmUzml/3Ti0Uqlfqlwv/+UzIgFvLwq49qJrdWzoVEIpmc3LYYQkNDweFw/OItBpR4KxSKwoyMjPnvtD+K8fFxiESBW1l7NEajcVEdvvyBzWbD+Pg4QkPn3z84EGhoaEB6erpfhqcsRsAruzSo6jrsaQn5XFyQH4tvDqkxOmbBA6dloHNkDE39egTxONjV4swvcjkEaxKdhWyOAFRvfwq3i6ioKOh0usnZ07Gxsejq6pp8vru7e3J+diDBtNC5J7qXLYawkCA8fHomyjtG8fL26QsUjyY0NNTtWeQRERGEEOLzQRoBJd4cDmdlWtriIxAGg4ExYmixWBAUFBRwd/ozwbTBKYDz86DX6/0aLViogF/60i6c83wZHvm4ZlLALyuKg9VO8eH+bmxOi0RBQhj+taMNOdES7Go9XBz06s+K8NdL8hA/sTXMYrHAaDRCp9NBo9HAbrdDr9djbGzMp1unAkG4gYmq5MxM1NfXAwDWrFmDpqYmtLW1wWKx4N1338XZZ5/tN/tmIjIyEgMDA/42Y954oop7sVxYEIvTV6rw5Ff12Fo396AUd3PeAJCTk8MFkOrWIosgoMTbarXGuRPeNBgMjPEM2ZC592lpaUFqaqrfb5AWIuCuMZNv7+nEY58cAqUUqUox1iSG4Z29zmKcu09MRZ/WhM6RcdT26aAZs8BsNsOgGUaWyIh9+/bhhx9+wN69e1FbW4umpia0t7fDbDajoaEBBw8exLZt27Bt2zZUVVWho6MDWq3WK4IeKMLtIiIiAiaTCWNjY+DxeHjuuedwyimnIDMzExdffDGys7P9beIx8Hg8CIVCj84q9ybu5pHdgRCCP1+0CjnRUvz8nQOo7Z09fM/hcCbTrYslLy9PTAjxed570SFqT0MI4SUmJga5E940GAweb8DhLZhUFU8pxcjICPLy8vxtyryxWCwYHR1Fbm6uv00BcFjAKyoq0NjYiJkiTBcXxqG6W4sMlRiv7eqAVBSEX5yUhsuK4vGL96uwq2UY61MiUJgQNtki9V+fbkNJbDBkMhlkMhliYmIQEhJyzE2LRqNBYWHh5N+tVuukV97c3Ay9Xo/g4GDEx8dDqVS6nWoINOF2sWLFCrS1tSE7Oxunn346Tj/9dH+bNCeu0HlycrK/TZkTV82R3W73S/1RcBAXr1xTiHOeK8NV/9qDN28oRmbUzGlpgUAAs9m86G2kmZmZQUqlsgDAG4s0eVEEkuedkJSU5Fayjkk5ZCZ53qOjo5DJZH73YBdCZ2cn4uPjA8rm+XjgFxXGIloqhJDPxYUFsfjb1iZ8fUiN01dGQRrMx+u7OzA4OIiTow837tAKlFi3bh2ys7MRExOD0NDQeb1vPp8PuVyO5ORkFBQUYPPmzcjMzMTw8DB++ukn1NfXw2w2z7nOdASqcAPO3PfAwACjmp8olUpG5b3FYrFbhWDuopQI8c5NJQjicXDZy7tR0zNzJMDdvHdaWhr4fL7PPZtAEu+03NxctzpomM1mCAQCT9njVcbGxhhzozE0NMSofDelFF1dXQHZf30uARfwuLjzhFRUdmlwQoYCebFS/PL9KnSPjuGcbDm+rlGjoqEDl2zOmyxQ29nquY5WEokEOTk52LRpE0QiEXbv3o2qqipYLPPfOx7Iwg04/w1iYmICblvYbAQHB8NqtTKmzau7fcM9wYqIELx/81qEBPFw+cu7Udk1vT3uiveKFStgtVq930TiKAJGvEUiUWZOTs6it4m5epoHkqc1E1arFTwejxG2AszqWgc4bzbCwsLA5wfmMI65BPzCglgkyEV45psG/P2yfAh4HFz9chnyZGZwCEG5NgRSqRQPnZ4JAGgeMMx7a8xCbIyPj8fGjRsRERGBsrKyeQ3JCHThdhEfH4/OzoU19PA3EomEMUNW/LVd7GjiwkV4/5a1kImCcOUreyZbCk/F3aI1Pp8PPp8fTHx8QQ8Y8ZbJZBkJCQmLfvNGo5ExxWpMCu8DzCoEBICOjg4kJCT424xZmU3A+VwOHjs7G62DRryzsxE3ZRP0G+3Y2i/AaStVeG9fFwxmG/Ljw8CZ+Mb0aManOYv7EEIQExOD0tJSdHd3Y//+/TN64UwRbgAQCoUQCAR+K6xaDP7agrUYPLEFy1PEyILx/s1roRALcPW/92J365HtWz1h68TnPcytRRZIwIg3l8tNcGfgBZMEkUliaLFYwOfzGRMlsFqtMBqNjIgUzCbg65PCUBwjxL929+CEtQW475QMfHFQDVEQF3qzDf8td+5P3v+rk/DgaRmIlnm3P75AIMCaNWugUqlQVlZ2TP6VScLtIiEhgVHedyCEoudLUFDQglIt3kYlFeLdm0sQIwvGtf/Zix1Nh1u4BgcHH9F5bzHExsYSAD6d2BQw4m2326PdEe/x8XHGDPhgkngzqbAOOLyljSk3Gy4B12q1+O+2KowaLTCZTNi5cyfuOyERXA4XT33bgps2JGFDagS2VPYiJIiL/+xsh8NBERYShFs2JYPL8c37dXnhbW1taG5uBsBM4Qac+6eHhoYY03pUIpEwxvMmhIDD4fi8TepsKMTOIrZEeQiue20ffmhw7p33xLUiISGBj+Uq3g6HQxYWtviog8lk8urEKE/CpCiBVqtlzJY2gJn91zkcDmRx6bj3i26sfvxb3P/GNqRnZqEwKxl3nZiKb2v78U1tP/5y8SqIhXwYLXZ0DI/h+3r/NO4QCAQoLi6GVqvFwYMHsWvXLsYJN+D8vYeGhvq1Knoh8Hg8OBwOxhSteaIBiqeJCBXgnRtLkKYMxc2v78e3tc5GLlwu163dB4mJiSIsV/Hm8/l8d+6A3Nmn52uY1AmOSZ63w+GAXq9n1M2GC7n48Gf3kxYrbvtvI7pGxnD9+hXIipLgVx8fBI9D8NdLDu9I+XdZmz9MBeAUvqysLHR3d0MikTBOuF0wrfWou724fUkg5b2nEhYShLduKEFmtAS3vbUf39f3QygUwmQyLXrN2NjYIJlMluRBM+ckIMSbEBIsEoncssVkMjFimxjThqcwKcQ/PDyM8PBwxoTMpxLCo7giwzlqNScyCPV9Opz7fBkO9mjxp4vyoBmz4tFPDmFDaiRu3exs1LGzZRh1ff6p6DWZTNizZw/WrFkDAJMh9IXy2GOPISYmBqtWrcKqVavwxRdfeNLMOVEqlejvn7uNZqDgz9ajCyVQxRsApMF8vHF9ETJUEtzyZgUaNHBLvKOiohASErL8xBtAVHR0tFuJJ6bs8TabzYyZh820YjW1Wg136ib8hcPhwL59+3DnabnYkBqBhhErrs4NQRCH4rKXdqN92Ig7j0/FJ1W9+KpGjV+clIakCGfk5nef1vrc3qk57oiICKxevRrDw8Po7e1d1Hr33HMPKisrUVlZ6fNuZ675Au5cuH0Jk4rWAjFsPhWJkI/XryvCCnkInvhpEHvbhuc+aAaioqLA4XDiPWjenASMeCcmJrrdqpUJIsOkyWdMyndTSjE0NMTI8G1zczMiIyOhUirx3GX5iAsX4b06E+4sFCM5PAi3vVUBIZ+D7Ghn+NxgsuG164oAAGabbwuCpitO43A4WL16NRoaGhgjglNRqVSM8b792Td8oYhEIoyNjfnbjFkJCwnCmzcUIzKUjwc+a0Xr4OIiBVFRUbDb7csy5x2dkJCw6FJxh8PBCOEGmFVYx6R8t9FohEgkYkw6woVOp4NarUZqqnMokVTEx7+vWQMK4KVqE+4sCMWGxFD88ct6xIWJoB13hs/jwkXY+eDx+M/Pinxm62xV5UFBQcjMzER1dfWCq7efe+455Obm4rrrrvPL3Gomhc6ZVLTG5/NhtVr9bcacRIoF+Nv56SAAbnitHNqxhdssFotBKfVpfjEgxFsoFMbGx8cvWryZVKzGJPHW6/WQSPwyZ37BaDQauLNbwR84HA5UVlZi1apVRwwBSYwIwb+uWYN+nRnPHjDj2pXBOCVVgq8OqWG108nwebQsGNJg33SRm892MJVKBR6Pd0zb0RNPPBE5OTnH/GzZsgW33norWlpaUFlZiaioKPzyl7/0xds5gtDQ0IAO7x4NU+xlikMFAMkqGe4vlaF7dBy3vrUfVvvCb454PB7fl13WAkK8pVJpokKhWPTxTClWA5hnK1NuNJgUJXDR1NQEpVI57Q1SQUIYXr66EG2DRvy90oIrsgQ4N0s2+fzDHx3EiNE3TTAWso975cqVaGpqOiJ8/t1336GmpuaYn3POOQdKpRJcLhccDgc33ngj9u7d6+23cwyEkIBrKjIb7lZG+5JA2+s9E0KhECtCHfjj+Suxs2UYf9vatOA1xGIxAeCznGhAiDefzw93x8NjkufNJFstFgtjiuuY1n/dYDCgv79/MlwOAM0DepQ1D8Fodu43XZ8agecuX42DPTr8vcqGC1L5uDwvHAAwYrTgVx8f9LqdC23AwufzkZWVhZqamnmtP7Vf+kcffYScnJxF2+oOTKriZpJ4M8VWV4j/goJYXFQQi+d/aJ62D/psiMViCmDR8zkWSkCIN4fDkYrFi3/PrkEfTIBJ3izAjNAXpZRRNxoA0NLSgvT09MlwOaUU5zxXhite2YP8x7/Fne8cwMFuLU7OVuH/Ll2Fik4N/lJhxSkJXNxQ6BTRLw6q8dGBbq/ZuNjOaUqlEmazeV6h3fvvvx8rV65Ebm4ufvjhB/z1r391x+RFw6Qqbtf8aSYgFAoZYevU69yjZ2cjJiwY97xfCb1p/vnvieJen+W9A0K8KaUSd8TbZrMxRryZsqXNbrcfkYcNZJjUsQ5wfgY0Gg2mpooIIciOdlb2R0mF+LF+AGc9twNXvLIbifIQPH95Pmp6tXim3Ix1URzcXuIc0XrPe1Ww2DxfvORuy9OkpCS0tLTM+bo33ngDBw8eRHV1NT755BO/bfVjkngzxZsFnDcaTLHVRaiAh79evArdo+P4v+/mHz6XSqUEy83zppSK3WkEwiTxppQyQhSZcpMBMC/f3d7ejsTExGOiGv+8Mh+ZURJ0jY7jls3JeOi0DNT36XHWczvwY8MA/nDeSjQOGPDUPhPyI4B7Nyhx6Zo4UHi2N7cnepWrVCqMjIwwotoYCPw9yVNhkngzyVYAkzslChPDcVFBLF7f1YGukfltd5NKpTwsN/F2OByh7njedrudMeLNFJgU3tfpdIzaj97X14eYmJhjnpOHCvDBLWuxKS0Sz3zdgMouDT6+fR2uX7cC/93fjd99Wou1SXI0Dxjwxz1jSJfYcd1KIQQ8z22P89SQEUIIoqOjF924xdcQQhiztYlJgsiU3yng7G8+dQveL05KB4cDPP11w7yOl8lkXCy3sLndbhe6MxGMKZ43U6YXAcwSbyZNlBsZGYFEIpnx8xoq4OGVqwvx0GkZ+La2H+f/cycKE8Pw9T0bUZgYhp8aB2FzUDQPGvHE7jG0qkePGSe6WDw9HSwuLg5dXV0esMw3eGI0pC/g8XhuDdHwJTwejxHV5sCxw0lUUiFu3JCET6t60dg/dz95mUwWhOXmeXO5XK47hVFMEW+m2AkwS7yZZGtfXx9iY2NnfQ2HQ3DzpmR8fPs6RIYKcMubFXjqy3r88fxc/OfaNZOtUTuGx3DLF0No7B1xW8C9MdYzODgYhBBGFCwBzCyuCnSYdqNxtK0/W7cCQVwO3tk799z3ZSneHDeTwDabjRGdtZgk3kzKeVutVvD5vmlW4i4Lyc/nxEix5Y51ePC0DPzUOIiT/vITejTj+OKuDfjVGZmTr7vjq2G3BNyb87jZQjDv4O4IS1/BdPEODwnCydlKfHSgBybr7BEEiUTCCQkJ8VmnqIAQb66byssUUWSKnQCzvFmAGd4IpRRWq3VBW9r4XA5u2ZSMr+/eiJwYKX71cQ0uemEXVsXJUP6rE5EfLwMAjAQpodPpFizg3hRugFnizaTKaKbcaDBdvAHgsqJ4aMask7O/ZyI0NBTBwcHh3rLvaPyuJIQQbnJystueNxNEkUmFdUzxvJm0pW0h41XNNjt6RsfB5RAEB3Ghkgrx9o3F2FLZiye/rMeFL+zCWXnR+Pvl+eBxCOQhQeCQBFRUVKCxsRFpaWlznsPbwg049752d3tvL7onEQqFfumtvhhce70DfVwvUyIEwMz5+bVJckiEPOxsGcJZedEzHi8WixEUFLR8xBuAICgoyK1KLqaIIlNuMgDm/E6ZFCGYTxe4/R0j+M2WQ6jr08Fx1LdCJRFiRUQISpPl2N48hE+revHlwT7cdlwKbtmUBFEQD/n5+fMScF8IN+DcghXok6VcMCXnDThFkQmFYEzyvGe60eBwCAoTw7G3bfaOawKBAIQQn1XOzumyEEL+TQgZIITUTHlsFSFkNyGkkhBSTggpmvLcQ4SQZkJIAyHklCmPb5547dNH28Dlct0Sb0opI8KmTMnNA8z5nTJJvDUazZxb2u55rwrDBgvuOD4Vf74oD3++KA9PnJuDX56UhtIUOcw2O7bWD2BQ7xQZm4Pib1ubkPWbr7GndRgcDgf5+fmzhtB9JdzA4b7hTBBFJoXNCSGM2L1y9ParQGa232lhYhhaBo0YNsz8OeZyuSCEzOsCTwg5dUIjmwkhD048Fk0I+Z4QsoUQMmdIZT6u1asAngPw+pTHngbwW0rpl4SQ0yf+vpkQkgXgUgDZAKIBfEcISaOU2gHcCmADgCcIIRmU0vqJtdytV2MMTBFEgDljVu12O2NuiAwGA5KTk2d9jSiIixABD3edkAouZ/rfP6UUnSNjqOzSoLJLg/+UtQMADnRpUJwknxTw6TxwXwq3C7FYDIPBEPBpGCZta+JwOIwQbyYx2/UuK8o5e6NtyAh56PSfYw6HMy/xnnjN8wBOAtANYB8h5BMAVwO4E0ASgCsBvDDbOnOqJqV0G4Cj4wUUgGuSiBSAqxPDOQDepZSaKaVtAJoBuLxyzsRxDgBTf0vLRrwBZhRWAczpBMekG6L5pCJu3ZyMuj4d3tjVPuNrCCFIkIfgnFUxePSsbLQ/eQba/ng6btl0+MZgOg+cUupz4QaYE+JlijcLMMtWpjDb71Qe4hTs2Sb5TVwv53PRLALQTCltpZRaALwLp3Zy4dTHozVy+vPN40TTcTeAZwghXQD+BOChicdjAEztytA98RgAvAJgJwAOpbRuqg1MEAlPwKQvG1NsZcpNBjC/4rqz86InO6xpx+bfmWq6G5ipAl5dXQ2j0ehz4XbZxpS8J1NCvAAYcUMEMMdOSumM//4ykXMr6q7W4VmPt9vt88nhzaSTzwF4EcAtAN6ca5HFViTdCuAeSumHhJCLAfwLwImY/m6BAgCl9GsAX0/zPHE4HCE//vjjIk1xhiPdOd5XWK1W2O32I8YgBip6vR4//fSTv82YE9fvVK1W+9uUOTEYDNi2bduckYIEnhU/Wez4bOt2xIjdvzFxOBzQ650dog4e9P4Y0aPR6/Xg8Xge6wTnTfR6PSOuJXq9HlwuF83Nzf42ZU6Ycn02GAzgcDhobW095rlRk1PU/7uvHZvEg9Mev3XrVvT29p41j1NNq5OU0g4AG+dr72LF+xoAd038+QM4vWrAeQcRN+V1sTgcUp8JByHEuHnz5kWPhfrxxx+xefPmxR7uM3p7e6HX65Genu5vU+bkxx9/xKZNmwI+JK1WqzE6OorMzMy5X+xntm/fjtLS0llz9AM6E+7dsQNpylBc9v/tnXd4W+XZxu+jYctbnpL3TmLHjp3YSewsRqCsAmVDGWWV3ZZRaIH2K5s2EKCMQtk7zDADDQGynMTxivfeU8PakrX1fn/IMk7iraNxnPO7Ll2XLUvveSRL5z7vM8/ZBM40ce+54opxr1u3DpWVlUhISJhTGRmdtLS0IDIyEmKx2KvHXQhMOZe0t7cjNDQUCQnTly75C0x5Tzs7OyEQCKbsgNg7agD27MFfz8nDyWtTp3x+SEgI3nnnnc/mcKiF6ORxLPSyfhjASeM/nwrANTftawCXUxQVSFFUOoBsABWzrOVgkqvKHfxdCCfDlJgakxJ3OBzOjC5Eh4Pgrk9qoTdb8eJvV9Em3Hl5eYiNjUVQUNCCGrm4i8PhYERSIVM+RwCzcj2Ywkz/f8V4rDtBOH0l2LiOzUXMKgFkUxSVTlFUAJxJ3l/Pw1QAc9h5UxS1DcDJAGIoihoE8A8Avwfwb4qieABMAG4CAEJIE0VRnwBoBmADcPt4pvlMnFDizZQTBFNEkUnvqWvC0nQd1sp7FDjQqcCjv8nDEtH0LZINZhtaJTq0jGixq1kKcbgAT1yYf1R2+mThjomJAeB8r+ZaB04nFouFET0DmCSITKkGYRrTvadynbOEMCp4+u6IDocDc9A7EEJsFEXdAWcYmQvgTUJI03xtnfUbRQi5Ypo/FU3z+McBPD4PG2gRbyZ88SiKYkxCjMtWf98x8Xg8xowcDA8Ph0ajQUjI1BEis9X52fihSQKz1Q4BnwtCCPRmO4bUYxhUGdE7akCfcgyTr1cCeRz837m5CAl0fp2nEm4XM5WReQqdTgd3Rv56C5vNxpge+UxJ1GTKhTUws4bUDmjA51IzXlTPY+cNQsh3AL5bgJkT+MPlsMNut7uluq7dl7+LN5O6DTFlR8uUHs+As8+3UqmcNk558tJY3H36ErxX3of9HaNHPzeYj6TIIOQmhCNMwEfDkAYAcH5hAh48J2dOwu3CmwLucDjgcDgYsfNmUsMfJpzvAGb1YZjpPa3sVSIvMQJBAdO/FrvdPqedN134wzfKZDKZ3PoUukRxPgMffAGTmkAw5T119XhmAkKhED09PdP+naIo/HFzNv5wahY0RissdgcoOHubhwby0DtqwEPfNKFhSIOlojA8cv5yrM34pexrLsLtwlsCrtVqER4ePvsD/QCTyeT3jWRcMKVbI5NaQk/3no5ZbKgfVOP69ekzPn9sbAx2u332wd804fN3lRBCUlNT3VI0pogik3berh1tcHCwr02ZEabE5oFf3tPZdk0URUE4KbZmtNjx9M42vLqvGwE8Dv7+61xcU5oKPvcXt+l8hNuFNwR8Li1h/QUm7byZYivTxHsqW79rkMBqJzh1WdyMz9fr9bBYLDM3QKcRvwiauBv0Zkrck0kTdpjkjgaYE1sLCgqC0Wic8+N/bJbitGf24sXdnThnRTx+vuck3LAhHRU9SmyvcU7rWohwu5hLL3R3mMswFn+BKYIIMMdWpgw4AqYX720V/ciIDcGa9JkHhul0Ouj1+um7uNCMX4i33c1tM1N2tEyxE2CWeDPJdR4ZGQmFYvbvt2bMirs/rsWN71YhNJCHj28qwbOXFUJttOKSVw7iytcP4+5P6mA0Ghcs3C48JeCEEKhUKnbn7QGYkrDG9J13u1SH6j4VLl+dPGuOgUajsVssFo0nbZyMX/z3CSFu6TdT3OZM6fEMMEsQmXShkZSUhIGBgRkfU9WrxJn/3oev6obxx83Z+OYPG5ARG4r7tzfgV8/uQ2Wvc+b0sxfn4fDhw24JtwtPCLharUZoaChjTt5ms5kR4s0ULxPAfPH+948dCOJzcdGq4xu3HItGozED8FrM2y/Em8PhjBkMhgU/nyk7WiZkh7pgkiAGBwfDnc+PNwkJCQEhZErXOSEEr+7rwmWvliOAx8EXt63DrSdl4r97u3DyU7uxraJ/4rHPXpyHSH0PLcLtgm4BHxgYQHJy8uwP9BPGxsYYId4Wi4VRiXVMFe+afhV2NIzg95sypp0kNhm1Wm0FoPegiUfhL+Kt0+sX/pqZIt5MgkniLRQKoVarfW3GnElOTkZvb+9R99nsDjzwRSOe+K4Vv8oV4bNb1qGqV4VNT+3G1l3tWJMehbToYFAUsOWCXNqF2wVdAm61WqFQKBAXN3OSj7/g8ogxIYObSVnxFouFkbXzhBA8saMFMaGBuHlTxpyer1arbTjRdt4UReldgxMWApPEe7YWmf4Ck9zmERER0Gi8Fmpym8TEREgkkqM+s4/taMG2in5cWpyEpeIwnP38fjzybTMyY0Pw4e/Xgs/loFcxhkd+vQxxxj4YwpJRI/XMZ54OAe/v70dy8uxxQn+BaSVtTPAQAMyydTLvH+5HVZ8K9/xqyUQPhdlQq9UOnGg7bwAad3bebBY3/TDlpAs431Oz2cyYWCCXy0VSUhL6+51ucInGhLcP9gIAPqsexHM/diAnPhwf/n4t3r9hLT443I8fmqX4+9lLkGTpx4+jYbj901bc8n61x2x0R8AdDgf6+/uRmjr1AAd/RKPRsIl1HoApeQST6ZLr8fiOZmxaEovLiuce9tFqtQRe3Hn7RTDC4XCo3dl5BwQEQKvV0miR53CJ93QtMv0JV4tUJmS1BgcHw2g0+n1duou0tDSUlZUhMTERwmA+zs4XQ64zozQjGr9ZmYiM2FA4HAR//rQOO+pHcN+vspBuH8TXI6H4ot45/vS1a4o9auNC68A7OzuRkJDAGHcp4EyuS0lJ8bUZc8JkMjGi3SzAHBe/q/eC1e7AXR/XIojPxVMXr5jXcKATUrwtFsuoOzHLwMBARuxmAebsvIFfXOdBQdNP0vEXXHFvpog3n89HTk4O6uvrsXr1avznyqNHBRBC8OCXDdh+ZAh3npqBbAzjgx4BfmqXAQDOK0jA5hyRx+2cr4BrtVpIJBJs2LDB47bRCZN23mazGbGxsb42Y06YzWZGiLfZbEZAQAD+76sm1A9q8MpVqyAKn5/HQKPRcHCixbwVCkX3yMjIghu1ME0QmWLrfBuK+JKIiAhGJa0BgFgsBpfLxdDQ0FH3E0Lw8DfN2FYxgFs3pWEpR4LXWnn4qd3ZvCkmNAAPnbfca3bO1YXucDhQW1uLgoICRnhrXDApWQ0AjEYjIy6oAebUo5tMJvzQZ8O2in7ccUoWzsyLn/caRqPRQQjxWqKQX7yrDodjqK+vb2yhzw8ICIDFYqHTJI/his8yASYJYmRkJJRKr3UmpI28vDy0t7cf9ZnYsrMNbx/sxXWlKcjlyfBCA4XqQR0yYp2hlkfPz0NUiHd7zs9FwDs6OiASiRizg3XBpBauAHNK2pgyPAUAdjVL8HatFmfni3H36fNvE0wIgc1m86oI+YV4Axjp6+tb8HaUKR8QgFleAiZlcQcGBoIQwpiLOBcBAQHIyclBXV0dCCH4vHoQL+/pwmVFicjmyrGlxo4uhRH3/GoJBpVGnLMiHmflx8Nic2DM4t0kzZkEXKvVQiqVIjs726s20YFEIoFYLPa1GXPCFUNmwjmPKfHuQ10KPLijG0tiBdh6SeG84twulEoluFyuygPmTYvfiHd/f79b/c2ZMiubSQ1FwsPDGZMICAAikQhSqdTXZswbsVgMgUCAr8vqcP8XDShJj0QSRvFkpQVjVoJ3r1+DnY0ShAp4eOS85bDaHcj7x05c+J+DXrd1KgE3Go2oqalBYWEhI1ykxyKXy2mvl/cUTOoVPzY25vc5KDX9KtzwTiXEYXw8+5vsGUd+zsTw8DA4HM4wzebNiL9806QSicStS0mm7Gj5fD5jytpcMUAm1KUDThGUSCS+NmNBZC3NxeN7JBAGchBL6fFcjRnxwiB8eft61PSrUTeowSPnL0d0aCC2/K8VFrsDrRKv5cYcxWQBb25uRkVFBfLz8xlTJz0ZvV6PoKAgxnQBY5KLX6/XIzQ01NdmTEvjkAa/e7MCsWGBePjkaIgjF27ryMgIbDZbH43mzYpfiDchxGY2m91SCKaIN+AUcCZMQQOY5ToPCwuDXq9nzMXGZF7a0wmZwQ7lmA3fdFtxyjIRPr91Hcw2O57d1Y6z8sQ4Jz8ePzRJ8Np+50zwJy7I95m9HA4Hy5cvR29vL8LCwhAdHT37k/wQJrnMAad4M2Xn7c/iXd2nwhWvlSNcwMcHN65FMGV1K49gZGQEGo2mg0YTZ8UvxBsAbDab2Z2TLpPEOzQ0FO40pfEmTBJviqIQExMzp6ld/sYPTU6PAY/LwU0rw/DntaEI4nPx50/rERLIxSPn52FQZcSfP60DAATxubhgZaLP7B0bG8Phw4dRXFwMh8PhkXGi3kAqlUIk8nzJHV0wwRXtwmAw+GU/i4Ndo7j6jcOIDgnAxzeXICky2O1mMn19fcaxsbGZJw7RjN+IN5/PH5XJZAt+PpNKsEJDQxkT92Za33Amus5NJhPOSbLhD5uSsfvek3H/pRthMpnwr88OoHZAjX+cuxzhQTzc8WENtCZnyOWadakLjs+5i1QqxeHDh7FixQrExcV5dB64J7FYLHA4HIzI3AZ+qUVmQrIa4J8XGrvbZLjurUokCoPwyc2lSIp02uduM6q+vj4jgBGazJwTfiPeFEUNjYws/LW7OmwxgZCQEMbsvJmWtBYdHQ2FQsGI5EXAKdzl5eW46tRC3HP2CsRHBIGiKCRn5eCjZgNyojhYK+bgye9aUTeoQUZMCLgcCteUpnndVqvViiNHjqCvrw+lpaWIiooC4Ll54J5meHiYdZl7CFerYn+60Phf4whuercKWXGh+PjmUsSNN2Gx2+1u29nf32/DiSreVqu1Z3h44cl6THJFM8lWpiWtcTgcxMbGMmL37RLuqaaDPf1DG4xWB56+Yg0+Le/E2wd7cXlxIhQGC85YLkKi0Nmkw1v93GUyGcrKyhAbG4vVq1cft1tlooC7hqcwBSaJt781kvnyyBBu//AI8hMj8OHvS47qkzA2Nua2e39oaIjCiSrecrm8qaenZ8Fp2EzqBhYSEsIYtznAvN13amoq+vq8mvg5b2YS7roBNT6qHMC169LA4fHxco0B+fEh4Blk0BituKbEOfBDZbAg/f7v8FXt0FSHoIWxsTEcOXIEPT09KC0tRVJS0rS7FCYJuEajgUAgYIzLHGBWC1edTuc3yWqfVg3grk9qsSYtCu/dsBYRQUf33KcjsU6j0TgIIV49qfuNeNvt9ra6uroFK4TrhMKEyVKu2ApTdrNCoRAqlVf7D7hFWFgYHA6H314gzSTcAPDEdy2IDgnE79al4ab3qhAexMMb15XgsCIAmVEBGOurR1dXF/64rQYA0DBIb0IhIQQqlQpVVVWoqamBWCzGmjVr5iR0TBHw3t5eRk09I4T4bQLYVPhLPfpn1YO47/N6bMiKwVvXrZ5yvKe74q3X60EI8XpWr9+IN4D2hoYGtwqgg4ODMTa24C6rXoVJu9m4uDjGNT9JS0tDT0+Pr804jtmE+3C3Aod7lLjlpAz8dXs9pFoz/nt1MTpkenTIDLh1cw42btyIdrkJ+zudWfUlSYFuXwgSQjA2NobOzk7s27cPXV1dyMjIwIYNGxAfHz+vmKC/C7jFYoFKpUJcXJyvTZkzWq0WYWFhfhVDngl/qEffXjOIez+rw/rMGLx2TTEE/KkTPN29KOrs7ASPx/P6B92fOhNI3Yl5A7+4o5lwdeoqwYqMjPS1KbMSEhICi8UCm83GmGYW8fHxaG9v9yubZxNuAHjh505EhwSgpl+FA50KPH1JAQqThbjxnUpEhwTg1yviwedz8VHzL16FELseBw4cACEEYWFhEAqFEAqFCA8Pn/K1E0JgNBqhVquh0WigVqthMpkQFBQEsViMdevWuT3Oc6HjRL1Bf38/UlJSGCOEAPPq0X19Ht7fIce9n9WjNCN6RuEG3N95t7e3Q6/X1yx4gQXiH2c1AIQQEh8fb3BnhJwrEYwJV9RCoRD9/f2+NmPOxMbGQiaTISEhwdemzAkOh4OkpCT09/cjIyPD1+bMSbh7Rg0o6xwFAHzXIMG9ZyzFxUVJGFCO4adWGe44JQsCPhfVfUrsbZcDAJaJw1BaVADAWe6i1WqhVqsxMDAAjUZzVNa9TqfDnj17ADj7IgiFQkRGRiItLc0jyUX+KOAOhwMDAwOMG1kqlUpRUlLiazPmhMVi8WlJW6dMj9s+qEF2XChevaZ41pJKm83m1sVqU1PTmEqlql/wAgvEb8QbAPh8fld3d3d6Tk7Ogp4fGhp63HhFf4VJzU8AZ/10X18fY8QbcLrOy8rKkJqa6tNxj3MRbsDZrtHFnadl4/ZTsgA4E24A4PI1KQCAZ3a1IyKID43RipKMXzqbcTiciV33VOzZswcnn3yym69mfvibgPf390MsFrvtWfAmRqMRXC4XAQHenSS3UHzpMtearLjhnUoE8jh4/XfFCJ0ixj0Zq9XqtmeutrbWAMDrbnN/innDZDLVuRMjY0uwPEdkZCTUajUjEgJd8Pl8JCcno7u722c2zFW4AWB5QjguWJmI569YiT9tdk7nsjsIPq0exKbsWCQKg1DercCBTgU2ZDnXWpseNfH88m4F/rOn03MvZoH4SwzcZrOhp6cHWVlZPrNhIUilUka5zH1Z0vbU/9owoBzDK1cVTTRgmQk63Putra0A4PUvnl+Jt1wur2lubl5wmzQmzfUGnFnROp1vhkvMF4qiGJd1DgDp6ekYHBz0yediPsINABmxoXj2skKcV5Aw4XLc1yHHiMaEy1YngxCCrT+0QRQeiMw4Z4xuzSTxvvzVcmz5X5tnXoyb+IOAd3d3IyUlhVG7boB58W5fiXftgBrvH+7DNaVpKE6Lmv0JcD/eTQiBTqezEUK8XqfsV+INoL2urm7BW2eKosDj8Rgj4JGRkVAqlb42Y84wsfUol8tFRkYGOjq8OjNg3sI9HZ9UDiAqJACn5Yiwq1mKyl4V/nBqNmoH1FgiCkV0qDM/ZLLL3V/xpYBbLBYMDQ0hLS3Nq8d1F5vNBpPJxIgkXMC3JW0Pfd2EuLBA3POruYdmNBqNW9PwFAoFOByOfMELuIG/iXdHU1OTW30thUIhY2LJTBui4UpaYxopKSkYHR31WhMfuoR7VG/GrmYpLlyZCIoC/vm/VmTGhuDioiRU9SqxNv2XePdHlc7kR46fJ1D7SsA7OjqQmZnp09yHhSCXyxEbG+trM+aMwWBAcHCw15PVWiVa1A6occtJmQgTzN2z4q6XoL29HRRFNS94ATfwK/EmhGgUCgVxJ64aERHBmEEarhg9U+LIPB4PAoHAb5ufTAdFUViyZAna2jzvUqZLuAHg+0YJbA6CS4qT8VHlALrlBvz1rBx0yvQYs9ixetxlbrM7sKPe2ZkxJnRhlRrexNsCbjQaMTo6yqhWqC4kEgni4+N9bcacUSgUbn/uF8IXNUPgcSicVzD3hFpCCCwWy4KrmwCgpaXFoVKpqha8gBv4lXgDAIfDGXKn3ptJU7AoimJU3BtgpusccNptMBg8GrOnU7gB4MdmKdJjQpAYGYR//9iONelROC0nDtV9ztdQlOrsEXC4RwnVmHM+PBPEG/CugDc2NmLZsmWMqusGful0x4ReEC5GR0d9It7fNY7gpCWxE2GkuUDH1LOysjK1wWCocGuRBeJ34m21Wsuqq6sX/HwmjdsEnFOwRkdHfW3GnBGJRHBn+puvoCgKBQUFqK+v90iGP93CrTfbcKhLgdNy4vD8Tx0Y1Vvw4Nk5oCgK1X0qiMMFSIhwtiv9vnEEQXwulonDEB3KjHIiwDsCPjQ0BC6Xy6iZ3S6USiUiIyMZc9FBCHE7hrwQzDY7BlVG5CfNrzyNjsS68vJyO4BatxZZIH4n3jKZbO+hQ4cWrL6upDWr1UqnWR6DaXHvoKAgcLlcRnkLXISGhiIpKclV2kEbdAs3AOxvl8NidyAyJACv7e/GFWuSUZAsBADU9KuwKlUIiqLgcBDsbJLi5KWxMFhsjNl5u/CkgJvNZrS3tyMvL4/Wdb1FX18fo1z9rkQ1b19sDKmMIARIiZrfLtpd8bbZbFCr1WZCiE9Ohn4n3gCq9+3b59bWmY17e5bU1FRGdYebTEZGBtRqNW1Z/p4QbsBZs83lUNhW0Y9EYRAePCcXACDVmjCoMmJVitOVWt2vglxnxpl5YozqLIgOOXrn3S7V4Zkf2vz68+UJASeEoK6uDjk5OYxpbjIZq9UKrVaL6Ojo2R/sJ/jKZS7ROKuLxRHzmxDn7pS2lpYW8Pn8lgUv4Cb+KN797g6UYOPenkUsFkMmkzGmwcxk6HSfe0q4AWBIbYLdQTCgNGLLxSsmOkXVHBPv/r5BggAuB6UZ0TBa7cfF/H717D48/3MnFAb/Lp+kW8CHh4fB4/EYVR89mcHBwRnHr/ojvhJvPs8pYzb73C9QCSEwm81ujYStqqoiKpXq5wUv4CZ+J96EEMLhcIbciasyrfVoXFwco5LAOBwORCIRo2yeTGhoKJKTk91yn3tSuIFfXID/vDAf6zJ/Wb+mX4UAHgfLEyJgdxB83ziCjdkxMNucFZaTY96qSYLt2p34M3QJuMlkQnt7O/Lz82m0znsQQjAwMMAol7mrr763490AEBLgvLA1mOc+lJKmZDWVXq8/5NYibuB34g0ANpvtgDtJa0zbyYpEIsaN3ExNTUVfX5+vzVgwGRkZ0Gg0C3rfPS3cAPDgOTnofPysiX7mLppHtMgRhyGAx8GhLgVGNCZcsCoRo3ozACBmkni7stIBZog34L6AOxwO1NTUIDc3d8Gd1DRGK17d14XybgWsdrfaTiwItVqNoKAgt0qYvI1SqUR0dLRPPAVhAqd4a4xzz3OiKVnNBuCIW4u4gV+Kt1Qq3eNu0ppAIPBaUw53CQwMBEVRMJvNvjZlzrgSU5h0kTQZiqJQVFSElpaWeb0Gbwg3AHA5FHjc47+eHVI9skVhAJzzisMEPJyWI4J6vFRMGPyLeLdJf3ldEi0zxBtwT8AbGxsRExPjVnZ5n8KAJ79vxeWvlmPVI7tw2wfV+KRyADIvvYc9PT1IT0/3yrHowpctXOMjBAjkcdApm3tzztHRUbfyCex2O5RKpYUQol3wIm7il+INGpLWYmJiGFeCxTQ3dEZGhk+HfrhLYGAgVq1aherq6jm11PWWcE+HxmiFTGdGVlwoDGYbvm+U4NcrEiDgc6E1OcU7fFJ3qQ6pDqLwQHA5FGN23i4WIuC9vb2wWq3Izs5269grkoR47epihARwoTPb8F2DBPd9Xo81T/yES145iO8bRmB3eCYB0GQyQa/XMypRjRACuVzuk+8EAPC4HCwVh6FFMncddbd+vrW1FXw+n96ylXnir+Ld564oMK1+monNT+Li4qBSqRjTS34qwsPDsXTpUlRXVx81+/pYfC3cACZ2Ftlxofi+UQKj1Y6LViUCALQmZ7wvPOiX8YbtUj2WicMRExoAKYN23i7mI+AKhQIDAwMoLCykxXV7Wq4IX9y+HilRwQjgcnDzpgzce8ZSjGhMuPWDGpz89G68WdYD/TzirHOht7cXaWlpjEpU0+v1CAkJ8Wnr2RxxOFpGdHOqqjCbzeDz+W7ZW1VVRdRqtc+S1QA/FW9CCKEoqre3t3fBazCpxzngjNMbjUZGZXBTFIWUlBRGx74BID4+HtHR0Whqapry7/4g3ADQKXO6wbPjwrC9ZhApUcETWeda49E7b7uDoFOux1JxGCKDAyY6sDGNuQj42NgY6uvrUVxcTKuALBGF4cvb12NVqhD/3dcNncmGn+85GS9fuQqiMAEe+bYZpU/8hMd3NGNI7X6Izm63Y2RkBImJiTRY7z0kEonPm+Asiw+D0mCBTDd76JGOrPhdu3YpdTrdPrcWcRO/FG8A0Ol0X+3Zs2fB2SIcDgeBgYGMiXsDzHP1A86hH4ODgzPuWplAdnY2zGYzjr1g9BfhBoBuuQEBPA7URgsOdilwafEvpUQ6kw0BXA4EfKd49SkMsNgcyI4LRWRwANRjx3tHLDYHdrfKMKga8+rrmC8zCbjVakVVVRUKCgoQFBRE+7GjQgLw3g1rceXaFLyytwu3vl+NDdkx+OzWdfjy9vU4eVkc3jzQi01bduOOD2twpH/h7XeHhoYgFosZNzzFH+aN58Q7s9xbRmZ3nbsb7waA/fv3OwBUurWIm/iteGu12p3ffPONW500mCaGCQkJGBoa8rUZ84LH40EkEmFgYMDXprgFRVFYuXIlhoaGMDg4CMC/hBsATFY7ggO4eP6nDoQLeLhmXdrE37Qm63Euc8C5e4wKCYByCvH+5/etuO7tSvzfV1N7HPyJqQTcZrPh8OHDyMrKQlTU3OY3LwQ+l4PHL8jHo+cvx552OS56+SD6FWMoTBbihStWYt99p+CGDenY2ybHBf85iItenn9c3OFwoLu7m3GJakajERRF+TwzPkfsEu/Zk0/djXcPDw/DbrcPE0J8Gi/0W/EGUFtRUeHWdo5p4h0ZGQmNRsMo1zng3LV2d3czzu5j4XK5WLt2LXp7e9HX1+dXwu1CPWbFjy0y3LQp46jkNK3RetTvTcMacCineAuD+UfVfLuQ6pxx8PJuBWw+KImaL5MFvLW1FRUVFUhLS0NCwtwnSbnD1aVpePf6NZBqzTjvpTLsbXeOcU4UBuGBs3Nw6IHN+Me5uZDpnHHxk57ajTfKeqAzzR6y6O/vh0gkcqtpiC8YGhryCzd/RDAficKgWXfeJpPJ7Xj3nj17iNFo/GbBC9CE34o3IcQOwK24N5PapALO3R8TE9f4fD6SkpIYnXnugsfjoaCgAA0NDYiPj/cr4eaPl45FhQTg2vVH79B0JttEvSsANAxpkB0XhqAALqJCAqAxWuE4JpnHbHVebI1Z7Gieg7vRH+BwOFixYgV6e3vB4/GQlJTk1eOvz4rB13eshzhcgGvfqsCLP3fAMb7DDg3k4br16djz51PwylVFiI8Q4NFvm7HuyZ/x2LfN04YnbDYbenp6kJWV5c2XQgtDQ0Neu3iajWXisFnFm46RpTt27FAplcrv3VqEBvxWvAFAq9V+tXv3brfi3gKBAGNj/h3Tm0xSUtKE25ZJZGRkYHBwkDEDYabDZDKhuroaxcXFkMvlfvW/uH5DOv51UT6+un39RLtUFwazDSHj9xFC0DikQV6is29zZHAAHAQ4NmfNbHMgKdIZJy7vZsZwHKvVisrKSuTk5IDD4XhlHvixpEaHYPtt63BeQQKe/qEdN79fPVGqBzhr9M/ME+PTW9bhq9vX45RlcXjrYC9OemoPbn6vCgc6R4/Kiu7p6UFycvKCm8r4Cq1WC4FA4De943Piw9E9aoDJOr0HkI5ktbKyMjsAn8zwnoy/i7fbce+4uDjIZDK6TPI4YWFhMJvNjGrYAjhdzhkZGejo6PC1KQtmcoxbLBZj7dq16Ovr85ts+gRhEC5bnYLkKaYnEQCu6iKJ1oRRvQX5ic44YGSIUxT01qN33iarHcmRwciIDUF5Nz2DWjyJxWLB4cOHkZaWhtTUVK/NA5+K4AAenrusEP84Nxe7W2U4/8UDaJceH28tSBbi+StWYv99p+D3GzNQ0aPEla8fxuZn9uKtAz0Y1RgwODjIuFg34Oy/7k8tXJeKw2B3EPSMTt0ihBAyMWZ1ofhLvBvwc/GGM+7tVjcEJrqhExMTMTw87Gsz5k1ycjJkMhmjMvxdTJWcxufzsXbtWkilUjQ1Nfn1ZC4KgMu8hkFniaRrvnHkeNc1neVY8XZAwOdgbXo0KnuUHms8QgdarRYHDx5EVlbWRIzVG/PAZ4KiKFy3Ph0f/r4EOpMNv3npAL6tn/p7myAMwl/PWoZD92/GM5cWIFzAx8PfNGP9lr3YPhCIVuncu4P5A4QQSKVSn5eITSZ+fKrYdD0NdDodQkNDweEsXPb8Jd4N+Ll4j8e9u92Je4eEhMBkMsFmo7eZgidJTExkXNY54DyZLlmyBG1tbb42ZV7MlFXO4/GwevVqcLlcHD582G/DAhT1i3g3DjmT1XLjneIdHeLMBNYeI95mmx2BPC5KMqKgM9vQPDx1vPC1fd245b1qr7UHPRaJRIKamhqsWrXquJIkXws4AKxJj8KOP25ATnw47vjwCB7f0TxtAqCAz8WFq5Lw5e3r8dnvi7E2notdHRqc83wZLnr5IL48MgSzzf8TP5VKJYRCoV+VtcWGOT/n8mlqvelo4frtt98q/SHeDfi5eAOATqf7evfu3W59mmNjYyGXy+kyyeMIBAJwuVzo9cy6GgecDU90Oh1jep7PpRyMoigsW7YMKSkpOHDggF/+XyhQIHCKc93gL8lqACCKcJ7UVKapd94lGc6a1+ni3j+2SPG/JgnOe/HAjPFEuiGEoKOjA93d3Vi3bt20E6v8QcBF4QJs+30JrilNxWv7e3DVG4cnhsVMB083gn9dtAKH7z8NfzsnBwq9GXd+XIt1T/6MLf9r9ev6+/7+fr9ymQOTxHua950OT0FZWZkDfhDvBhgg3lqt9vvt27e7FZBjouucqVO7KIpCTk6OW+M2vcV867gTEhKwcuVKVFVV+V8exfjO22C24VC3AuuyfmlCERMSCD6XgvI48bZDwOdCFC5ARkzItOK9Ytz9LtGa8HWtd8I5drsdNTU1MBqNKCkpmTUpyh8EPIDHwSPn5+GZSwtwpF+NXz9fhpppmrbodDro9XqIxWJEBPNx48YM/HzPyXj3+jVYlRqJV/Z2YdOW3bjxnUrsbZdPZLT7A1arFRqNxu/6r7veIh7n+NayJpMJHA7HreS6/v5+2O32AX+IdwMMEG8AtVVVVXZ33N5RUVFQq9V+HbM8FrFYDJlMxsja6ZiYGNhsNqhUC+825WkW2oAlIiICpaWl6OjoQGtrq990luNQzpao+zvksNgc+FXuL+5BDoeCOEIApeloW802x0RHtrUZUajonTruHR/hzEgXhwvw5oEej3+PtFotDhw4gOjoaKxYsWLOMUp/EHAAuHBVEj6/dR14XAqX/fcQPjjcd9x71tLSgpycnKN6mHM4FDYticVr1xRj332n4NaTM3GkX43fvVmBU7buwWv7uqfslOdtBgcHkZSU5Hf918csTo0IDuAd9zc6dt1ffvmlVaPRvOfWIjTi9+JNCCFcLnffgQMHFrwGRVGIiIjwazE5Fg6HA5FIhJGREV+bsiCWL1+OhoYGvxG3ybjbOS0wMBClpaXgcDgoKyvzix76KVHB6B41YGeTFMJgPlanHZ1RGx8RNOXOO5DnPAWUZERDZ7JNWSfrSgQ6Y7kIrRIdDvd4JjPd4XCgvb0dtbW1KCgoQFpa2rzX8BcBz0uMwLd/2IB1mTF48ItG/OXz+olYtmuG/Ew716TIYNx7xjIcvP9U/PvyQsSGBuLx71qw9omf8KePjmB3m8wnjXUIIX7pMgcAg9n5/gYHHB+HpyPe/cEHH6gNBsMXbi1CI34v3gAwNDT09scff+zWGZKJrvO0tLTjem0zhfDwcIhEIr8rHaOr5akrOa+wsBD19fVoa2vz6YVKfpIQSoMFXxwZwuZlouNmgSdECKAw/iLehBCYbQ4Eunbe6dPHvcXj4r0mPRrCYD7ePtA7qz0jGiMufvnglOVTU+HabRNCsGHDBkRERMzpeVNxrIA/s6sdN75TOZGF7y2EwQF489rV+MOpWfikahDXvFEBuWYMzc3NWLFixZzWCORxcX5hIj67dR2+++NGXFyUhD1tclz3ViVKnvwJD33dhLoB73kVlUolQkNDfd4OdSqahp3/3/SYkKPut9vtMBqNCA0NXfDaWq0W/f39BkKI38QyGSHeAHZ/++23Fnc+oHFxcYxKWgOA4OBg8Pl8v9jZLYTs7GxIpVJotf7RvcsTvcrDw8Oxfv16UBTl0114QdIvYndVScpxf48XBkFtJhNucbPNeaEh4DtPAeIIAdKig6es93a5zVVjFlyxJgU/NEtmTab6qUWGqj4VLv3vIVhn2CFO3m2vWLECS5cudauUx8VkAVcpFfixRYZzXyzDLe9Vo2OOFxR0wOVQuOdXS/HcZYU40q/G+S/uQ3BcyoLaoOYmhOPxC/JR8eBm/PfqIqxOi8KHh/tx/ksHsHnrXjz/Uwf6FZ5Ncuvu7kZGRoZHj7FQDnUpEBrIQ37i0Rd+crkcsbGxbq29c+dOYrfbt7u1CM0wQrwJISaHw9HmTgkSj8dDQECAX2YKz0RGRgZj245yOBwUFhaitrbW5+5zTw4ZOXYX3tjY6PUZ5/mJEXjuskIc/OupWJlyfBOKhAgB7AQTGdBmq/P/Ecj7xcVYkhGNih7FcXHv2LBAcDkUJBoTripJBUVReK985g2Ia8etHrPipd2dUz5GJpOhrKxs2t22xmh1a0fpEvAzU3kQhTrjoHvaZfjVc/tw9ye1s2aD08lvVibi+QuzoDU7cOv2blT3LTz0EMjj4ozlYrx8VREqHzwN/7wwH7FhgXhmVzs2PbUbF798EO+X903Zz94djEYjTCaTW01OPIXDQXCgcxSr0yKP8zqNjIzQ4TJXyOXyj9xahGYYId4AMDo6+s727dvdKjRNSkpiXP10TEwMNBoN4zquufAH97m3poOFh4dPiNCBAwfQ3t7utf4CFEXhNysTkSCceiyma/c8PD532jQef3XtvAGneGuniHtzORREYYEY1hiRKAzCGctF+KhiAEbL9MmUmbFOF+XKFCFe/LkTjUO/eCTUajUOHjyIgYEBFBUVTbnbNlrsWPvEj7j7kzq3msdwOByUrinC7/KdrtQbN2Tgpo0Z+LZ+BGc+t39iuIinsVgsCNQO4rObSxAu4OGK1w7jmzr3M/cjgvm4fE0KPr65FGV/OQX3nbkUGqMVf/uyEWue+BG/f7cK3zWM0FLi19vbu6A8BG/wdd0wehVjOL/w6CEpdrsdarXaralz49Pr7ACq3TSTVhgj3maz+ett27a55e+Kj4/HyMgIo7LOKYpCeno6Y3ffgNN9LpPJfOJS9vZYT4qikJycjE2bNoHL5WL//v3o6enxuechLcbZUrVL7mwdqTU6m81M7pG+NsN5gpsqIU0cIYBE47x2vnZdOjRGK744Mv2F8KYlTjfl5mVxiAoJwD2f1EGp1qKyshItLS1Yvnw5ioqKEBISMuXzDRYbTFYHvjgyhHs/dV/Ab/51KQpFAXizrBs3bEjHN3dsQFQIH797swKPfdvs8cYojY2NWLJkCZYlRuGL29ajICkCf9h2BC/t7qTtfJQUGYzbTs7CD3dtwo4/bsC169JQN6DGbR/UYPXjP+Ivn9VjT5tsQUJutVohkUj8YoLYsZisdjy1sw15ieE4r+DoISmuLHN3MuMPHToELpe7nxDiV9m3jBFvQohMoVAo3Ylb83g8hIWFMWrSGOBsOyqRSPy2u9dscDgcFBQUoK6uzqsi5st53FwuF5mZmdiwYQPMZjP27duHgYEBn4l4ekwoBFygYVANABgc34G7BpMAzt15anTwlElr8RFBE+K9Oi0SyxPC8fbB6cvG0qKDkRIVjNoBNR4+ZwnapDr849PDSE9PR2lp6awJaa5l02NCsP3IEO79zD0B53K52HplCSx2B/784WEsEYXi6zs24JrSVLxe1oMLXjqITplnQmoSiQR2u31i+lZkSADev3Etzi9MwFM72/CXz+tnzAuYLxRFYXlCBB48JxeH7t+M925Yg9NzRfimfhjXvlWJwkd+wA1vV+K98r45N4Lp6+tDSkoKLfkIdEIIwcPfNGNIbcT9Z+WAc0yNt6uszR0++eQT7dDQ0DtuLeIB/Os/MQsmk+njHTt2uHWJzMSpXRwOBykpKYxs2uIiPDwcYrHYa+5zXwr3ZPh8PpYtW4aSkhJotVrs3bsX7e3tXg+DcDkUUsM5qB93Xw+qnOKdKDx6yElpRjTKuxXHiUl8hAAjGhMIIaAoCteuS0O7VI9DXdNPI1uTHIqyDjmCdf34dW40dnRbMDA2t3aaMaEBWCIKRSCPg7tPX4LtNUO477N6twQ8My4Md562BPt69XhhRxUEfC4eOT8Pr19TjBGNEee+UIaPKvpp9cxZLBa0tLRgxYoVR+3+AnlcPHdZIf64ORufVA3i2rcqoDHSf3HO5VDYmB2LZy4tRM3fT8db163GZcXJaJfp8PcvG7HhX7vxq2f34snvW6b8vwNO1/PAwABSU1Npt88dCCH4x9dN2FbRj1tPzsT6rKO/5xaLBUajcdrOfHPlm2++MQP4ya1FPACjxFulUn38zjvvuFWs7WqV6ms35nxJTU3FwMAAI5u2uMjKyvKK+9xfhHsyAoEAy5cvx8aNGxEQEIDy8nJUV1djdHTUa2GctAgOmoe1sNodGFIZwedSiAs7uuTn5KVx0JlsqOw92nUujhDAaLVPCMy5BQmICgnAWwd7j3qcxWJBV1cX9u7di6xQK0w2goCEHDx5aRHE4QL85fN6WGyzf/coisLNmzLRKtEhPzECd522BJ/XDOIvn7sn4LeenIXSjGi8dEiGHysaAQCn5Yrwvzs3oSg1En/d3oDbP6yB5tj5qQvE5S6fqrSKoijcffoSbL2kABU9Slz08kEMKD2XLS7gc3HK0jg8fH4e9t17Cn68+yT87ZwcxIQG4o39Pbj81XKsenQXbv+gBp9WDUz0CB8YGEBCQgJ4vOObn/gKzZgV931Wj3cP9eGmTRm474ylxz1meHjYbTd/a2srrFZrByHE76YtMUq8CSHNra2tenearXA4HMTExDCubIzH4yEhIQEDAwO+NmXBuNzntbW1Hkvk8kfhngyPx0NaWho2bdqE9PR09Pf3Y+/evWhra/N4F8D0cC7MNgc6pHoMqY1IEAYd52bcmB2DAC4HP7Uc3f7VlfA2Mu46F/C5+O2aFPzYIkW3TIOhoSFUVVXh4MGDAIB169bhytNXg8ehsK9DjjABH49dkId2qR6v7O2ak73nFiQgPkKAl/d24U+nZeOu05bgs+pB/PXz+gW3C+VyKDx3eSFCAvl49KchNDY72/iKwgV49/o1uP+sZfihSYqz/r0Ph92ccT40NAS73T6rgFxUlIR3r18LmdaEC/5zAHUDareOOxcoikJWXChu3JiBD39fgiP/dzpeuWoVzs6LR2WvEvd+Vo/Vj/+IX7+wH49914Z6rQBdcr3P27Q6HAQfVfTjlK178HnNIG4/JRP3n7Vsypj20NCQ2+L91ltv6eVy+YtuLeIhGCXeAGAymd785JNP3DrzM9F1DgDp6eno7e1lnNdgMuHh4UhPT8eRI0doFyp/F+7JUBSFqKgorFq1CuvXr0dISAg6OzuxZ88e1NXVQSqV0u5lSY9wft2r+pToHTUgcYrM9JBAHkozo/Fji/So/0+80FmX7Ip7j42NYWOCcxTpli8OQ6fTISsrCyeddBIyMzMREBCAMAEfq1IjsbfNeaF86jIRzitIwIs/d6JTNnvuaQCPgxs2pKOiR4mafhX+dFo27jwtG59WD+KRb5sX/PkRhQuw9dJC9GlseOnA8EQnNg6Hws0nZWL7besQyOfiitfKsfWHtgV1MtNoNOjo6EBhYeGcHl+aGY3tt61HUAAXV75+eNqe6J4iTMDHmXnx+NfFK3D4gc349g8b8OdfLQGP2FE2bMe9nzdi89a9KHjkB1z1+mE8tbMVO5sk047fpJthtRFvlvXg1y+U4a/bG5ARE4Jv/rAB954xtXCPjTk9GEFBU1dfzAVCCD788EOj1Wr9asGLeBCKSZnXAEBRVFJhYWHNkSNHFlx1TwjBnj17sHHjRr9yBc2FpqYmREREuJ2E4WsaGxvB5/OxdOnx7q6FwCThngmHwwGFQgGJRAKFQoHg4GBER0dDKBQiIiLCrc/r7t27sbWRh37FGLQmG/561jLcclLmcY9771Av/v5VE368exOy4sJACEGPVIVTnzuEO9ZGY3WUGQEBARCLxdhSNoqyLiXK79+MkMDjbXtpdyee2tmGygdPQ2xYIEb1Zpz2zF5kxYbik5tLj9v5H4vBbMO6f/6Mkowo/PfqYhBC8PiOFrxe1oM/bc7GXacvWfD78cg3zXjzQA/+si4Cm5fFYcmSX9YymG146OsmfFo9iIJkIZ69tAAZsXPr0GU2m3Ho0CEUFxfPu6vXiMaIK14tx6jegneuX4OiVN/VVBNCsG/fPhQVr8aQzo66ATVqB9WoH1SjdUQH2/guXBwuwIqkCCwVhyE+IggJQgEShEFIEAYdVc0wn+PKdGb0jhpwZECN7xslE96IpaIw3HJyBn5TmDhjBnl7ezsCAwPditMfPHgQl1xyybdDQ0PnLngRD8Is5QJACBlMSEiQ9Pf3x6akHN9Jai5QFDVRNuaPPXpnIisrC4cOHUJCQoLfZX7Oh9zcXBw+fBjDw8MTWbgLZbEIN+AMLcTGxiI2NhaEEOj1eiiVSgwODqKpqQkOhwNhYWEQCoUQCoUIDQ1FQEDAnEphKIrCI+fn4cL/HERMaACuKT3+xGa321GS4hScD/Y04PQkpxjxAp07GAdPgHXrisDn8wEAN3KE+L7pILYfGcLVJcevtz4rBk/tbEN5twLnFiQgJjQQfz8nF/d8WocPDvfh6tK0GW0OCeThmtJUvLi7E11yPTJjQ/HgOTnQmqz4908diAji4/oN6bO+9qn4y1lLcbhHgf8eGUOGUAWgfULAQwJ5eOqSApy8NA4PftmAs5/fjwfPzploUjMdDocDVVVVyM3NXVA7zviIIHx0Uykuf/UQfvdmBd65fjWKUhdeo+wOQ0NDiIqKQmhIMJaGAEvFYbh0tfN8abLa0TyiRd2AGnUDatQPavBTq+y4fIQwAQ+JwiDERwggjghCANf53rneQ9dbSQgwpDaiT2FAv3IMJusv3o4VSRG494ylOCtPPKcLKEIIhoeHsX79erde/6uvvqoaHh5+3q1FPAjjdt4AEBwcfNMDDzzw3N/+9rcF+0TGxsZw5MgRt//BvqC1tRWBgYFIT1/YSctfsFqtOHDgAFauXLngXtaLSbjngsPhgE6ng1qthlqthsFgmOjmRlEUBAIBBAIBAgMDJ0Sdw+GAoii0trZiyZIl+L5VhdhgLnKiuTCZTDCZTBNliBwOB0FBQfjLbjXCBAH4+Ka1EAgEIIQg/f7v8MfN2bh70m6XEILzXzoAg9mGH+8+6Thhs9kdWPnoLvx6RTyevHDFxHOuebMCR/rV+OGuTdM2lnExqjdj/T9/xm8KE/Gvi51r2B0Ed3xYg+8bJXjq4hW4pHhhF+Hdcj1+/UIZ8hMjcPcqHiIjIo7agQOAVGvCfZ/VY2+7HBuzY/DUxQUT/d4nQwhBXV0dwsLCkJl5vEdjPkg0JlzxWjlkWhPevWGN1wXc4XBg3759KC0tnXMfc5vdAZnOjBGNEUNqE4bVRoyonT+PaIyQak2w2p16Q4hr+jwAAoByVjSkRIUgLToYqdHBSIkOwRJR6ES+xVxRKpXo6elBUVHRvJ43GavViuTkZKlUKk0khPhlljAjxZuiKGFGRkZ7V1eXWw1ry8vLsXz5coSFhdFlmlewWq0oKytjpNv/WHQ6HaqqqrBu3bp5Dzs40YR7NhwOx4QYm81mWCwW50mSEDgcDnR2dmLp0qUTgh4YGDgh9nw+/yjh3fpDG17a3Ynqv52OyBDnDOScv/8PV5em4oGzc4467vaaQdz9SR3eu2ENNmYf/5W88Z0qdMh02HvvKRP3DSjH8Ktn92FdZjRe/13xrJ6Dv3/ZiI8q+3HgL6ciLtwpnGabHTe+U4UDnaP4z5VFODNvYS0wP68exD2f1uHOzVnYEKlDeHj4cQJOCMEHh/vx+I4WBPA4eOw3eTj3mIYg3d3d0Gg0KCwspGVc5mQBf+f6NShO856Ad3d3w2q10hbW8iY1NTVITU11a974t99+ixtvvPENiURyI42m0Qoj/a6EELXJZGqur693a53U1FRGTu3i8/lISUlBV9fcsnb9mbCwMOTm5qKqqmpeiXiscB8Ph8NBcHAwoqKiEB8fj9TUVKSlpSE9PR2ZmZkIDAxEWloaUlNTkZycjLi4OISHh0/pdt+cI4KDALuapRP3BQVwp2yJes6KeMSEBuCNsp4p7VqXGY0+xdhRDUGSo4Jxz6+W4KdWGb6tn33s7Y0b02FzELx76JdeB4E8Ll65qggFyUL8cdsRlHWMzrrOVFxUlIQLVibi+Z87wYnNmHKcKEVRuKokFd/9aSPSY0Lwh21H8MdtRyZKyuRyOYaHh4+r53YHcYQA235fgrhwAX73ZgWqehfeD30+2Gw29PX1ue098AUWiwU6nc6tdqgA8J///GdUKpW+TJNZHoGR4g0AEonkuddff92tcVUikQijo6OMrJ1OS0vDyMiI1wdgeAKRSASRSISGhoY5PZ4Vbs9TkBSB7LhQvHOodyKrO4jPxdgU4h3I4+K69enY0yafsszJ1Tzj4DENXa5bn46CpAg89HXTrEM0UqND8KtcEd4/3HfUBURIIA9vX7sGGbEhuOm9qgVnaT98/nKIwwX482f1yM0vmHYeeHpMCD67pRT3nL4E3zWM4Izn9uGH+n40NjaiuLgYXO7cmtDMFXGEAB/dVALRuIAfW3/vCbq6upCamspIr97AwACSk5PduoDSarWoqakZA1BDn2X0w1jxdjgc33/++edmd8qmOBwOEhISGFk25mq/OdUJholkZmbCbrejp2fq3ZsLVri9A0VRuG59OpqGtagY73UeFMCdti/2NaWpEAbz8e+fju+gt0QUiuiQgOO6sXE5FJ68cAU0Rise3dE8q003bMiAesyK7UeO/r5GBPPx7vVrEBsWiOveqkSrZP7X9OECPrZcXIBuuQHP/thx1DzwY+FxOfjD5mx8cdt6hARycdOHDdgpjwDh8Od93LkgChdg2yQBrx9vcesJzGYzRkZG/HYAyUwQQibE2x0+//xzu9VqfZ/4eUyZseJNCDETQnbv3r3brXVSU1PR19fHqGElLpKSkqBQKGAwGHxtittQFIWCggIMDQ1hZGRqNyor3N7lgpWJEAbz8daBXgDOk6NtmovlMAEfv9+YgZ9bZag9ZvdNURRKM6NxsOv4bnK5CeG49eRMbK8Zwv8aZ3afr06LxIqkCLxR1nNcs5C4cAHev2EtBHwOrn6jAr2j8/9ObMiOwZVrU/B6WQ+q+9UzCjgA5MaH4u9r+LisMBYf1YzgnOf3H/fa6cIl4JHBAbj5vWqPjTNta2tDVlYWIytZ5HI5hELhRCXEQnnppZeUSqXyDZrM8hjM+w9NYmRk5Kmnn356YYGucQQCAYKDg+FO1zZfQVEUli9fjsbGRl+bQgtcLhdr165FR0cHpFLpUX9jhdv7BAVwceXaFPzQLEGrRIueUQOWiqfvE/27dWmIDObj3z8eL3brs2Ig1ZrRPYWo/nFzNlYkReCv2xtmbPpBURRu2JCObrkBe9plx/09OSoY79+wFja7A1e8Vr4gAX/g7BwkCoNw76d1MNkc0wq4w+FAZWUl0pIT8a/L1+DDG9fCZLXjwv8cwCPfNENvpr+DoChcgP9eXQSlwYLbP6ihdZgJ4Gwso9Vq/XJy2Fzo7u5GRkaGW2t0dnZiZGRkkBDi92McGS3ehJCq2tpa9bEn+vmSkZHB2JGbMTEx4PF4x4kdU+Hz+SgpKUFra+tEC1tWuH3H1SVp4FAUfvdmBRzEOZ97OkIDefj9pgzsbpMftwNdPZ4pXdN3/EUyn8vBs5cVwmS1497P6mf0gp2dHw9xuACv7586vJItCsMHN5bAZLXj8lfnL+AhgTw8dXEBehVj2PK/NnA4nOME3FXLHRcXN+FeXpcVg//dtQlXrk3FWwd7sHnrHnzXQP/44bzECDx5YT4O9yjxxHcttK1LCEFDQwPy8/NpS7jzJgaDAXa73e0hJM8//7xWJpM9QZNZHoXR4g0AOp3umVdffdWtbv6RkZEYGxuDyeSdVn90k5ubi5aWFkYm3k1FQEAA1q5di6amJoyMjLDC7UPEEQK8cpWzXjaAy0FhknDGx19T6tx9P3fM7jstOhgBXA465VOP3XQ2X8nFvnb5URnlx8LncnDt+jQc7FKgaXjqATe5CeH48PfO8Z+Xv1qOnnkKeGlmNK5dl4a3D/biYNfoUQLe1taGmpoaREVFHbfLCxfw8ehv8vDFbesRHRKI2z6owXVvV6JfQe+wkQtXJeG69Wl460AvttfQk68zODiI8PDwBfdb8DU9PT1u970wm8349NNPjTabzS/boR4L48XbYDC899prrxncSVyjKArp6emMLb0KCgpCUlISOjs7fW0KbQgEAhQWFqKqqgopKSmscPuQ03JF+Omek/H9nRsnar6nIzSQh5s2ZWJPmxxHJmV+87gcpMeEoGuGmdlXrU3BKUtj8cR3LeiQTt/7/IrVKQjic/FmWe+0j8mJD8eHv187LuCH0D3NRcN03HfmUqRFB+O+z+qhN9vA4XBQWFiIvr4+mM1mZGVlTfvcwmQhvr5jPf7v17mo7FHi9Gf34sWfO2C20Xdx/cDZOVibHoX7tzegcci9KX1WqxWdnZ3IycmZ/cF+iNVqhVwuh1i8sDp/F59//rndbrd/RAihfzarB2C8eBNC9Far9cedO3e6tU5iYiJkMtlEpymmkZGRAYlEMtGQn+mYTCbU1tZi1apV6O/vZ9wUuMVGaCAPmXPs7X1NaSoig/nY+kP7UW7jrLhQdMwg3hRFYcvFBQgN5OFPH9VOOzo0IpiPS4uT8HXdEJQzlJgtEzsF3GonuPzV8nkJeHAAD09fUoAhtRFPftcCh8OBI0eOIDU1FYGBgbNWefC4HFy/IR0/3XMyNufE4ekf2nH2v/fPOP98PvC5HLx05SpEhTgT2LSmhZ+3WltbkZmZ6Xail6/o7u5GWlqa20l2W7ZsUcrl8ufoscrzMF68AUAikTzxxBNPuJW4xuFwkJaWxtjYN4fDwfLly+dcK+3PTI5xJyYmoqSkBM3NzZDJjk9SYvE/QgJ5uPO0JSjrHMWn1b+4dbPiQjGgHJu23AwAYsMC8c+LVqB5RItnp0h8c3FJcTKsdoLvZ8lQXyYOx7bfl8DucAp41zwEvDgtCtetS8eHFf349MdDiIqKwtKlS2fNQp+MOEKA/1xZhLeuWw3LeCLd3Z/U0pItHhMaiBd/uwpDaiPemsELMRMajQYajYZxMx5c2Gw2DA8PY6FzLlzU19dDJpN1EkJ66bHM8ywK8SaENHZ1dQ256/ZOSUnB8PCwx2ZNexpX8trw8LCvTVkwUyWnCQSCiSS26crIWPyLq0tSsTY9Co9+04xhtRGAU7wdBOiWzxyDPj1XhCvWJOOVvV3TztRenhCOjNgQfF07+2d9qTgM224qgYM4Bbxzht3/sdxxSgaCeRQ+a7dOdBybKoltNk5ZGocf7jwJd5yShW/qhrF56178d2/XlB3r5kNRaiROzxXh9bJuaIzz230TQlBfX8/YJDUA6O/vR1JSktvNcf71r3+pRkZG/kGTWV5hUYg3AMjl8ke3bt3qVsc1LpeLpKQk9Pf302WW18nPz0dbWxvMZs/UgXqSmbLKAwMDUVJSgq6uLvT1TZ/QxOIfcDgUnrq4AHZC8JfPnRnkKVHBAIAB1eyhnb+dk4vUqGDc9fHUu1SKonBeQQIqepUTM8ZnYokoDNt+XwJCgCtem5uAm81mtNRV49rVIlQOGrC/45fQzUIEPCiAiz+fsRTf/2kjCpKFePL7VmzcshtvlPXM6I2YjTtPy4bOZJu2Pe10dHZ2IiYmhrFJag6HA319fW43lFGpVPjpp590AH6ixTAvsWjE22azfbV9+3aDuw1L0tLS0NfXN68+2/5EQEAAli1bxjj3+VzKwQICAlBaWgq5XI76+nrG/o9OFFKig3H/2TnY3zGKDyv6J3bgCXOYEhUSyMMLV6yCcsyCm9+rnlLczitIACHAt/Vz8zRli8Lw0U1rQQhw2X8PzdhKVaPR4NChQ8jOzsYfzy5EUmQQnviu9ajmMAsRcADIigvDu9evwWe3lGKpOBSPftuMTVt2452DvQsS8eUJETgrT4w3y3qgHptbu2SdTofh4WFGDh5xMTg4CLFY7Has/tVXXzUZjcZ/E0IYdUJZNOJNCLHZbLZX33zzTbeaffP5fIjFYka2THURHx8PiqIY4z6fTx03l8tFUVERBAIBysvLF0Vv98XMlWtSsD4rGn/7shGPj9clZ8aFzOm5+UkReObSQlT3qfDXz4+v/86IDUVeYji+qZv75zwrLgyf3FyCkEAerni1HN81HB+GGR4expEjR1BcXIy4uDgE8ri494ylaBnR4svaoaMeu1ABB5wx9Q9uLMFHN5UgLSYE//i6Cac8vQfvl/dNm6w3HbefkgW92YbvGiSzPtbhcKC2thYFBQWM7KQGOF3+dDRlsdlseOmll3RarfZ1mkzzGsz8z02DQqF47qmnntK4G7N2NW1hYstUF0xxny+kAQtFUViyZAkyMjJw8OBBaLVuRUtYPAiHQ+E/vy3CrSdlQmO0Ypk4DMEBcx94cXZ+PO49Yym+rB3Giz8fXwp5XkEC6gY182rGkhEbii9uW4e8xAjc9kEN/rOnc2J0amtrK/r7+7F+/XqEhv6SXX/uigSsSIrA0zvbjtsduyPgAFCSEY2PbyrBBzeuRYIwCH/7shGnPL0Hb5T1QKabW+8J10z0uXRd6+rqQmxsLIRC4bxt9RdGRkYQHR097zHCx7Jt2za7xWL5lBDCuJPIohJvQojaYrF8um3bNreyQAIDAxEVFcXo5CiX+9zdsamexN3OaWKxGKtWrUJNTQ0kktl3HCy+ISKYj/vOXIbDD2zGJ7eUzvv5t52ciQtXJmLrrvbjXORn58cDAH5unV8lQnRoID64cS3OK0jAlv+14b7P6nDocAXsdjvWrl17nCuWw6Hw17OWYVhjwntTNJFxV8ApisL6rBh8dksp3rl+DeLCA/Hot80oeeInXP3GYXxWPQjdDOVgjvGNBmeWvDOtVouRkZHj5pUzCUIIurq63B5Z6nA48Mgjj6ikUumjNJnmVRaVeAOAVCp99JFHHlG5Gw/Nzs5GR0cHo3ff8fHx4HA4GBoamv3BXoaulqfh4eFYt24duru7Gf//WuwEB/AQLph/fJKiKDx5UT6KUyNxzyd1R40dTRQGQcDnYGg8nj4fBHwu/n15IW7dmIpPq4fwdKUJSRlLps28XpcZg6LUSHxxZOrvk7sCDjhf60lLYvHFbeux665NuO3kLPQqDPjzp3UofuxH3PZBNb44MohDXQp0SHVQGizolOnwx21HADhL7abD4XCgrq6O0e5ywLnrDg8PR3BwsFvrfPvttw69Xr+LEMLIK3/m/gengRAi0ev1P3z77bduqXdQUBBiYmIwMDBAl2k+YcWKFejo6PCryWN09yoPCAhASUkJjEYjampqFk2bWJZfCORx8d+rixAXHojr3q7EgU5nWweKoiAOF0Ayw0CTmVAqlVgTLMcjZ2ehQWLARS8fxIBy+mz4M5aL0DyixeA0GfN0CLiLbFEY/nzGUuy79xRsv20dLl+djMPdStz1cR2ueK0cpz+7D6se3YXTntmHhiENHv1NHn6VO32XsZaWFojFYsZmlwPOC5D29na3PQeEEDz44INKiUTyIE2meR1qMe5UKIpKz8vLq6ivr49xp37RYrHgwIED2LRpk9t1hL5EpVKhsbER69ev9/kVt6eHjPT29qK/vx8rV65EWFgY7eszmT179uDkk0/2tRlu0TNqwO/frUK3XI+7TluC20/JwhWvlcPuIPjs1nVzXseV8DQ8PIzi4mIEBQWhvFuBm9+rBgA8eHYOLilOOm4X3i3X49Ste/HQubm4dv30vbQdDgdqamoQHh5Oq4vaZnegS26AQm/GqMEChd4Mk9WBi4uSZtx1SyQS9Pb2Yu3atYyt6QaAvr4+GAwG5ObmurXOnj178Nvf/va74eHhc2gyzessSvEGgISEhB3btm07+6STTnJrnY6ODlAUNWMvYybQ2dkJk8mEvLw8n9ngrelgarUadXV1SExMRGZmJqNPVnSyGMQbAAxmGx78ogFf1g6jIFmIDqkOKVHB+N+dm+b0fL1ej7q6OgiFQixbtuyoC/PeUQPu+6weFb1KlGZE44kL85Eec3R2/OateyCOEOCDG0tmPI6nBHy+GI1GlJeXY926dW4nePkSu92Offv2YcOGDW6Xh5WWlirKy8tPJoQwdp7yonObuxgZGfnL/fff71bLVABIT0/HwMAAY3ueu8jMzIRer/dZYpc3x3oKhUJs2LABVqsVBw4cgF4/v6EULP5NSCAPz15WiMd+k4cOqQ6nLIvDS1eumvV5rkSnqqoq5OTkYPny5cd51NJiQvDRTSV44oJ8NA5rcMZz+/DS7s6jsriXicPRNDx7cjKdLvSF4rqAyM/PZ7RwA84e5snJyW4Ld3V1Nfr6+lqZLNzAIt55A0BCQkLZt99+u37Vqtm/2DPR29sLo9HI2Kk7LsxmMw4ePIiSkhIEBc3eKIMufDmPW6VSob6+nt2FY/HsvCdDCJnT/9RgMKC2tnbK3fZ0yLQm/OPrJnzfKEFqdDCKU6OwVByKp39ox0WrkvDkhflzstGXO/CWlhZwOBxGN2MBnJPDysrKaAlhnnnmmcqdO3eeQwgpp8k8n7Bod94AMDIycu8DDzzg9hiflJQUSKVSxs77dhEYGIgVK1agpqbGa93JfCncgHNW+4YNGybyF9hd+OJiNuF2xbYrKyun3W1PR1y4AC9fVYRXry5CanQI9nfI8cR3rbA7CG7YkDZnG321A5fL5VCpVIwuC3PR0dGBjIwMt4W7ra0N9fX1A0wXbmCR77wBID4+vnrnzp2rVqxY4dY6w8PDkMvlKCgooMky39HR0QGj0Qh335PZ8LVwH4trF56UlISMjIwTbhe+GHfeM7GQ3fZsKPRm6M02pEbPrUvcZLy5AzcYDKioqEBpaSkEAoFHj+VpjEYjDh8+jE2bNrmdcHveeeepduzYcYndbmdUH/OpWNQ7bwCQSCS33HHHHW7vvuPj46HVaqHT6egwy6dkZWXBarWit7fXY8fwN+EGftmFm81mHDhwYFH8L1mOx+FwLHi3PRvRoYELEm7Aeztwq9WKqqoqFBYWMl64AeduecmSJW4Ld21tLSorK7sXg3ADJ4B4E0IqOzo66vbt2+fWOhRFTczLZrq3gqIoFBYWYmBgAKOjbuf0HYc/CrcLLpeL3NxcLF++HHV1dThy5AiMxvk3+GDxPwghGB4exr59+2A2m7Fx40ZERUX52qyj8LSAE0JQU1ODzMxMREZG0r6+t1GpVBgbG0N8fLzba91+++0KiURyMw1m+QWLXrwBQCKR3HrHHXeMuiu6UVFRCAoKYnTbVBdcLhfFxcVoaGjA2NjsIxrnij8L92QiIyOxfv16iMViVFRUoKmpiR1ywmDkcjnKysogl8tRUlKCnJwcv+3N4EkBb21tRVhYGJKSkmhd1xcQQtDY2EjLvPE9e/agu7v7CCGkmibzfM4JId6EkHa5XP7zV1995XaWVm5uLtra2uDu8BN/ICgoCIWFhaiqqqLl9TBFuF1QFIX4+Hhs2rQJYWFhOHDgANrb2xfF//ZEQa1W49ChQ+jr68PKlStRUFDACFexJwR8aGgIWq2W8VUxLvr6+hAVFeV2syVCCP7whz+MSiSS22gyzS84IcQbACQSyd1//vOfle62zgwMDERaWprP6jbpJjIyEunp6aipqXErHMA04Z4MRVFISUmZKEPZv38/ent72Xnhfoxer0dVVRVaWlqQk5OD4uLio6aAMQE6BVytVqOzsxOrVq1aFImYZrMZPT09tJS4bd++3TE6OrqLENJBg2l+wwkj3oSQIZ1O98nbb7/tdreVtLQ0jI6OLpqEp+TkZISFhaGxsXFBAs5k4Z4Ml8tFZmYmNmzYAKPRiH379mFoaIjxOQ6LCZPJNJGrkJqaitLSUkaPtqRDwA0Gw8T8cXcbmPgLLS0tWLJkCXi8uY+PnQqbzYZ7771XKZFI7qHJNL/hhBFvAJDJZH976KGHVO7Wa1MUhby8vAWLnT+ybNky2Gw2dHTM7+J0sQj3ZPh8PnJyclBSUoLR0VHs27cPfX197MATH6LValFbW4vy8nLExsZiw4YNiI2N9bVZtOCOgJtMJlRWVmLlypUICVlYFry/4UpSS0hIcHutN99806LX6z8khDA/UekYTijxJoSoxsbG/vPvf//b7fTiqKgoCASCRZG8BjgvSAoKCqBSqdDXd/y84qlYjMI9GYFAgIKCApSUlMBkMmHfvn1obm5ms9O9BCEEEokEBw8eRFNTE+Lj43HSSSchISFhUbiGJ7MQAbdaraioqMDy5csZ7X2YDJ1JakajEY888ohaLpf/gybz/IoTSrwBQKlUbnn22WfVCoXbpd+LKnkNcJ5AiouLMTAwMGsP9MUu3JMJDAzE0qVLcdJJJyE8PBxVVVWorKyEXC5fNJ4Xf8JsNqOzsxN79+6FVCpFfn4+SktLIRKJFp1oT2Y+Am6321FZWYnMzMxF44EA6EtSA4Cnn356zGg0vkgIUbtvmf+x6DusTUVISMiVF1100Uvvvvuu24Nt+/r6oNFoPN6tzJtYLBYcOnQI+fn5U9bJnkjCPR0uD4VarUZCQgJSUlL8PsvZnzusEUIwOjo6MUcgOTkZSUlJiyaGOx9m68RGCEF1dTWioqKQkZHhAws9g6uT2oYNG9yOdQ8MDGD16tX9Uqk0mxCyKGtAT0jxpiiKiouLq96xY8fK4uJit9YihKC8vBzZ2dmLSshcYwSLiooQHh4+cT8r3Edjs9kwNDSE/v5+BAQEICEhASKRCAEBAb427Tj8TbwJIdBoNJBIJJBIJBAKhUhNTYVQKFzUO+y5MJ2AE0LQ0NAwkZexWCCE4PDhw7R5Es4++2zVjz/+eIXFYtlJg3l+yQkp3gBAUVTu8uXL99XX10e723ZvbGwMFRUVtFwx+hN6vR6VlZUTAs4K98xotVqMjIxAJpOBw+FAJBJBLBb7TQmTP4i33W7H6OgoJBIJlEolwsPDIRaLIRKJFtV3hw6OFXBCCJqamuBwOGiJCfsTLi8WHbMjdu3a5bj66qv3SCSSzTSY5recsOINACKR6JVHHnnkuptvvtntbdJidJ8DgE6nQ1VVFfLz89HY2MgK9xwxmUyQSqWQSCQwGo2IjY2FWCxGVFSUz066vhJvs9k88V6MjY0hJiZm4r1w98J5seMS8LCwMNhsNthsNqxYsWJRCTed7nKLxYJly5aN9vT0FBFC+mky0S85ocWboqgwsVjc3tTUJHa3B/JidZ8DwOjoKA4dOoSCggKkpKT42hzGYbPZjtptCoVCiESiiYoFb52IvSXedrsdWq124jVTFHWUF2IxCY83sNvt2LNnD3g8HjZt2rSo3j+6z5uPP/742LPPPvvM6Ojo32kwz685ocUbAEJCQq665JJLXnj77beF7q5F5xWkv+Bylaenp6O7u/u4GDjL/CCEQKVSQSaTQa1Ww2g0QiAQQCgUIiIiAkKhEEFBQR45QXtCvF1CrVaroVarodVqAQDh4eGIjo6GSCRCYGAgrcc8kSCEoLm5GVarFVarFREREYtiPrcLOj2Wg4ODKC4u7pdKpUsIIWYazPNrTnjxpiiKEolE1Tt27FhZVFTk9nqLyX1+bIzb5UJnBZw+CCEwmUxQq9XQaDRHCbpLzOkSdHfF22azQavVTtg5WahddoaHh/vtQBCmcWyM2zUxzBvzwL0B3Zudc845R7Vr167fWiyW/9Fgnt9zwos3AFAUtTwvL29vXV2d28lrdGdN+orpktNcAl5QUOB34xYXE0ajcUIkNRoNTCYTCCGgKAqBgYEQCAQQCARH/ez6fbrP8FTiTQiBzWaDyWSauJnN5qN+d/Ux4HK5CAsLm7ioYIXaczgcDtTX14PL5SIvL2/iwm22MjKmQLe7fNeuXWQ8Se1UGsxjBKx4jyMSif778MMP/+6WW25x28fnKrNav369X5YMzcZsWeVGoxEVFRVYunQpxGKxDyw8cXE4HBPieqzIun6f7jut0+mmbH7B4/GOuwCY/DuPx1tUcVZ/x263o6qqCpGRkcjOzj7uvV8MAt7d3Q2DwYD8/Hy31zKZTMjNzT0hktQmw4r3OBRFhYhEotbKysqk5ORkt9cbGRlBf38/1qxZw6gT31zLwSwWCyoqKpCSksImsTEEfygVY5kZi8WCyspKJCUlITU1ddrHMVnA1Wo16uvrsX79elo8N3feeafuvffee1yhUPyLBvMYA1unMQ4hxCCXy6++/PLLlXRc0MTHxyM4OBjd3d00WOcd5lPHHRAQgJKSEgwPD6Ojo4NtE8rC4iZGoxGHDh1CRkbGjMINeGYeuDewWq2ora3FqlWraBHuw4cPY9u2bd1KpfIpGsxjFKx4T8Jut+/p7Oz89tVXX6UlU3H58uUYHh6GSqWiYzmPspAGLDweD2vWrIFOp0NTUxMr4CwsC0Sn0+Hw4cPIy8tDfHz8nJ7DNAEnhKCurg5ZWVm0NC4ymUz47W9/OyqTyS4mhDhoMJFRsOJ9DDKZ7LZ//OMf8v5+90Mnri9XXV0drFa3x4h7DHc6p3E4HKxcuRIcDgdVVVWLZkgLC4u3GB0dRVVVFVatWoXo6Oh5PZdJAt7f3w8ej4ekpCRa1vvLX/6iU6vVTxFCOmlZkGGw4n0MhBDD6OjoVXS5z0NCQpCdnY3a2lq/3JnS0fKUoijk5uYiLi4OBw8eZEdmsrDMkZ6eHrS0tKC0tHTB5ZdMEHCtVove3l7k5eXRsl55eTk++uijbqVS+TQtCzIQVrynwGaz7e3q6vrmv//9Ly3u88TERAQEBMx5Tra3oLtXeWpqKvLy8lBeXg6lUkmDhSwsixOHw4G6ujoolUqsW7fO7Yl0/izgNpsNR44cwcqVK2mp5zaZTLjyyitHZTLZRSeiu9wFK97TIJPJbn/ooYdkdLjPASAvL2+igYs/4KkhI1FRUSgpKUFjYyPoeu9YWBYTZrMZhw4dQkhICG2JW4D/CnhDQwPS0tJoa+x033336VQq1RZCSBctCzIUVryngW73OZfLxapVq3DkyBFYLL4dL+vp6WBBQUFYt24dpFIpGhsb4XCcsBfHLCxHodVqcejQIWRlZSErK4v2MlJ/E/Cenh4QQmgrJy0vL8fHH3/cpVKpttKyIINhxXsGbDbbvu7u7q+ef/55Ex3rhYWFYdmyZaiqqvKZoHlrrCePx0NxcTH4fD7Ky8thMtHyFrKwMJbBwUHU1NSgqKgIIpHIY8fxFwEfHR3F4OAgCgoKaLlI0ev1J3R2+bGw4j0LUqn0tieeeGKwtraWlvXEYjFiY2PR2NhIy3rzwdvzuCmKwtKlS5GVlYVDhw5BJpN5/JgsLP6GK+YrlUqxfv36Kbvc0Y2vBdxgMKChoQGrV6+mLSxw7bXXqlUq1f0nurvcBSves0AIMclksl9feOGFozqdjpY1s7KyYLVa0dvbS8t6c8Hbwj2ZuLg4lJaWorOzE83NzawbneWEQavV4sCBA4iKisKqVavA5/O9dmxfCbjVakVVVRUKCwvdTsRz8eabb1r279+/V6VSvU7LgosAVrznACGkTalU3nfNNdeo6Yh/UxSFwsJCDAwMYHR0lAYLZ8aXwu1CIBCgtLQUPB4Phw4dYsvJWBY1hBD09fVNZFmnpqb6pE2ytwXcNfksMzMTkZGRtKzZ2tqK+++/f0Qmk11Fy4KLBFa854harX7r4MGDu9944w1ass24XC6Ki4vR0NCAsbExOpacEn8QbhcURWHJkiVYtmwZysvLMTIy4lN7WFg8gdVqRU1NDZRKJdavX+/z8bneFPDW1laEhYXR1ojFaDTi/PPPV8hksvMIIXpaFl0ksOI9D2Qy2dUPPvjgcHNzMy3rBQUFoaCgwGOdyfxJuCcTHR2N9evXY2BgADU1NX7dfY6FZT7IZDKUlZVBJBLRVtdMB94Q8KGhIWi1WuTk5NC25s0336yRy+WPEELqaVt0kcCK9zwghBhkMtl5v/nNbxR07ZajoqKQnp6OmpoaWjuw+atwuwgICMDq1asRFxeHsrIySKVSX5vEwrJgXAM3enp6UFpaStvOk048KeAqlQqdnZ1YtWoVbeGBjz76yLpz587DKpXqBVoWXGSw4j1PCCENo6OjD9100020dVtJTk5GeHg46uvraRFwfxduFxRFISkpCaWlpRPxQXYXzsI0XLvtqKgorFmzhrYkLU/gCQHX6/Wora2dKA2lg66uLtx5551SmUx2KfHHvtJ+ACveC0ClUr20a9euQx988AFtSrN06VIAcPsLxRThnoxAIMDq1asRGxvL7sJZGIPVakVdXR26u7tRUlKClJQUnySlzRc6BdxkMk0MVQkJCaHFPrPZjPPOO08hlUovJIT4R0tKP4QV7wVACCEymeyyu+++e6iuro6WNSmKwooVK6DRaBZcQsZE4XZx7C68urqabezC4pcQQjA8PIyysjJERkZi7dq1CAoK8rVZ84IOAbdaraioqEBeXh4iIiJosYsQguuuu04jlUq3EEIqaVl0kcKK9wIhhGhlMtkZ5557rlwul9OyJkVRKCoqwtDQ0LwzsZks3JMRCARYs2YNEhMTUV5ejq6uLrYunMVv0Ov1KC8vh1Qqxbp16xiz254KdwTcbrejoqICWVlZtJ5vnnvuOeNPP/20S6FQPEXboosUVrzdgBDSLpPJfnfWWWcp6epXzuVysWbNGnR0dMy5BnyxCPdkxGIxNm7cCJvNhv3793ulHp6FZTpsNhuam5tRU1ODpUuXYuXKlQgMDPS1WW6zEAEnhKC6uhoJCQlISEigzZZdu3Y5/vnPf3bIZLIr2Tj37LDi7SYmk+n73t7erddff72Grs8bn8/H6tWr0djYCK1WO9vxF51wu+ByuVi6dCmKi4vR1dXFutJZvI7LRb5//34EBQVh48aNiIqK8rVZtDIfASeEoL6+HuHh4UhPT6fNho6ODlx99dUSmUx2OiHEt5ObGAIr3jSgUCie3LVr1890DTABnDXgxcXFqK6unraJy2IW7smEhIRg7dq1E6701tZWNiudxeMolUocPHhwwkWenp7OWBf5bMxVwNva2gD8kmBLBxqNBmeddZZCKpWeQwhhByDMEYr1TtADRVGBsbGxVdu2bcvdvHkzbRdFKpUKdXV1xyXFnCjCfSwOhwO9vb3o6+tDamoq0tLSwOGw16BzYc+ePTj55JN9bYbfo9Pp0NLSAkIIcnJyfN4hzZs4HA7U1NQgPDwcS5YsOepvnZ2dUKlUKC4upu0ixm63Y/PmzaqamprbtVrtNloWPUFgz3o0QQgxy+Xy06+66ippVxd9Q28iIyOxYsUKHD58eKIf+Ikq3IBzh5CRkYGNGzfCarVi3759GBgYoLXBDcuJidFoxJEjR1BXV4fMzEysXbv2hBJuYPodeFdXF5RKJYqKimj1Ptxzzz265ubmN1nhnj/szptmKIpalZWV9UN1dXU0nV98pVKJ+vp6rFy5EkeOHDkhhXsqzGYzOjo6oFQqsWTJEohEokXr2nQXduc9NWazGZ2dnRgdHcXSpUvZzxCO3oHzeDzI5XKsXr2aVi/Xe++9Z/nzn/9cNh7nZktK5gkr3h4gIiLiyhUrVrzw888/R9I5AnBkZATV1dVYtWoVrVmei4GxsTG0t7dDq9UiMzMTCQkJJ/wJ+FhY8T4ao9GIzs5OKBQKZGRkIDk5mf3MTMLhcGDfvn2w2+045ZRTaBXuvXv3Oi655JIuuVy+ih04sjBYt7kH0Gg0H7S0tDx9+eWXa+iqUTaZTGhra0NOTg7a29vZkZrHEBwcjMLCQqxevRoqlQp79+5FX18f7Ha7r01j8TP0ej2OHDmCyspKREdH46STTmJ0vban6OnpgUAgQHh4ODo7O2lbt6GhAZdddtmwXC4/mRXuhcPuvD2ISCT6z6WXXnrVCy+8EObOOsfGuBUKBRoaGhjZ2clbWCwWdHd3Y2RkBMnJyUhLS/ObCU++4kTfeavVanR0dMBisSA7OxuxsbGsYE9DV1cXFAoFiouLAWDaJLb50tvbiw0bNkiHhoY2EUI8P2B8EcOKtwehKIqKi4vbfuedd/7q/vvvD17IGtMlp7li4GvWrEFw8IKWPiGw2Wzo7e3FwMAA4uLikJaWRlsPZqZxIoq3w+GARCJBT08PeDwesrOzF12dNt10dHRArVajqKhowlU+Uxb6XBkdHcXatWtHu7u7z2Zbn7oPK94ehqIofmxs7O4tW7YUX3vttfNqyTRbVrlKpUJtbS2KiopOuKzY+eJwODAyMoLe3l7weDykp6efcDuvE0m8zWYz+vr6MDQ0hNjYWKSlpSE0NNTXZvk1hBA0NzfDZDJh5cqVx8W43RFwvV6P0tJSZUdHx5Umk+l/dNp9osKKtxegKCokNjb28Ntvv73s7LPP5s7lOXMtB9PpdKiurkZ+fj6io6Nps3kx4xr+olQqkZSUhJSUlEXR6nI2Frt4E0IwOjqKvr4+jI2NISUlBUlJSSd8uGQuOBwO1NbWIiAgAMuXL5/2onYhAm61WnHqqaeqGhoa7lKr1e/QafeJDCveXoKiqOi4uLiqr7/+Om3t2rUzPna+ddxGoxEVFRVYsmQJ4uPj6TJ50WO1WjE0NIT+/n4EBgYiKSkJYrEYXO6crq8Yx2IVb71ej8HBQUgkEkRERCAtLQ1CofCE8qq4g81mQ1VVFWJiYpCVlTXr4+cj4A6HAxdffLF6//79W+Ry+ZN02czCirdXoSgqRSwWl+/evTt+2bJlUz5moQ1YLBYLKisrkZSUhNTUVLpMPmHQarUYHByEVCqFUChEUlISYmJiFpUALCbxNpvNGBoawtDQEPh8/sSFF7vLnh8WiwUVFRVISUlBSkrKnJ83VwG/4447tJ988sl7MpnsDjrsZfkFVry9DEVReUlJST/t3bs3LiMj46i/uds5zW63o6qqCpGRkcjOzl5UwuMtCCFQKpUYHByEUqlEXFwcEhISFsVOjunibbFYIJVKMTg4CJvNhsTERCQmJp4QIQ9P4PLYLVu2DCKRaN7Pn03A//a3v+lfe+21H2Qy2cXslDD6YcXbB1AUtTI5OXnnvn37YtPS0gDQ1/LU4XCgvr4eHA4H+fn5jBccX+JwOCCVSjEyMgKNRoPIyEiIxWLExsYy0rXORPE2GAyQSCSQSCSw2+0QiURITExkk8/cRKvVorq6GgUFBW5l308n4A899JDh5Zdf/lkmk11ACGGbLXgAVrx9BEVRxSkpKd/v378/Ji4ujtZe5YQQtLa2Qq/XY+XKlawrkQYcDgdUKhUkEgnkcjmCgoIgFoshEokgEAh8bd6cYIJ4E0Im3meZTAaBQDDxPrM9DehhdHQUDQ0NKC4uRliYWy0oABwv4I8//vjY888/v1cmk51HCLHRYDLLFLDi7UMoilqbmpq646mnnoo+5ZRTaO9V3tfXh76+PqxevZo98dGMXq+HRCKBVCqF1WpFVFQUYmJiEB0d7bduXH8Ub0IINBoNRkdHoVAoMDY2BqFQOOHhYC886aWnpweDg4NYvXo1rRedLgH/9NNPLW+//fZemUx2NivcnoUVbx/D5XJLU1JSduzduzdyPgkjc0WhUKC+vt5t9xjL9NjtdqhUqgkBmizmUVFRfrMz9wfxdjgc0Gg0UCgUE2IdERGB6OhoxMbGsg2HPITD4UBDQwNsNhsKCws9EvZ58sknx1544YWqkZGR0wghVtoPwHIUrHj7ARRFrU1JSfl27969Ma4YOJ2MjY2hqqoKaWlp88ooZVkYk8VcpVLBZDIhKCgIQqEQEREREAqFEAgEXs9H8LZ42+126HQ6qNVqqNVqaLVaOBwOhIeHIzo6GjExMQgODmbzMjyM2WxGVVUVRCIRMjMzPfJ+P/roo4YXX3yxTCaTncsKt3dgxdtPoCiqODk5+bs9e/bEHpuFTgc2mw1HjhxBUFAQcnNzaZ0QxDIzhBAYjUao1WpoNBqo1eoJQQ8NDZ24hYSEeFTUPSXedrsdBoMBer0eer0eBoMBOp0OhBCEh4dPXLC4xkuyeA+tVouamhrk5OQsKKN8Lvzf//2f4ZVXXtkrl8vPZ13l3oMVbz+CoqhVSUlJ3//0009x7g4AmApCCNrb26FUKlFUVISAgADaj8EyNwghMJlME4LnEj2TyQSKohAcHDwh5pNvgYGBC3Z5LkS8CSGwWq0wmUwTN7PZjLGxMRgMBpjNZnC5XISEhBx1ERIWFsbIjPzFxMjICNra2lBUVERLYtqxEEJw//33G956662fZDLZhWxWuXdhxdvPoCgqPz4+fudXX30Vv3r1ao8cw9Nfahb3cDgcMBgME+J4rHC6xszyeLwZbxRFTdwAoKWlBa7mQIQQ2O122Gy2aW+u4/D5/OMuIFxeg4CAANbt7Wd44yLdbrfjxhtv1H733XffymSya1jh9j6seM8TiqKSAbwLQAzAAeBVQsi/KYp6CsC5ACwAugBcRwhRUxSVBqAFQNv4EuWEkFvG1zoZwNMAfiaE3DfpGGlxcXF73nnnneQzzzzTI/5tlzstIyODjYMzEELIlII7WZAdDgdc329CCLq6upCVlTUhtlwuF3w+HzweD1wu9zjxZ3fOzMNsNk+UbeXk5HgkPGYymXDBBReoq6urX5PL5X9xNWChKEoAYB+AQAA8AJ8RQv5BUdQlAB4CkANgDSGkavzxaZjnuZFlEoQQ9jaPG4B4AKvGfw4D0A4gF8CvAPDG7/8XgH+N/5wGoHGatT4GEARgK4Blx/wtNjY2tvndd981Ew9htVpJdXU1qaqqIlar1VOHYfETdu/e7WsTWDyIXC4nP//8M5FIJB47hlqtJqtXr1ZGRkb+iRx/PqMAhI7/zAdwGEAJnKK9FMAeAMWTHr+gcyN7c97YrKV5QggZIYTUjP+sg/PKMZEQ8gP5JVmjHEDSHJbjACBw7uCP8j0SQuRyuXzNPffcU7Nly5Yx+l7BL/B4PKxatQqxsbEoKyuDRqPxxGFYWFg8CCHOpkxtbW0oKSnxWGLayMgI1q5dq2xpablFqVT+ewo7CCFEP/4rf/xGCCEthJC2Yx8/C9OeG1mcsOLtBuNun5VwXmFO5noA30/6PZ2iqCMURe2lKGrjpPtfB3AQAIcQ0nLs+oQQvVwu37R169af7rzzTt34FSntpKSkoKioCLW1tejp6YGnjsPCwkIvRqMRBw8eBACsW7fOY82YOjo6UFJSMtrd3X2RTqf7ZLrHURTFpSiqFoAMwC5CyLHnxmNZ0LmRhY15LxiKokIB7AXwOCFk+6T7HwRQDOBCQgihKCoQTleSgqKoIgBfAlhOCNHO41icuLi4l0866aTLPvjggwg+n0/vixnHbrejsbERFosFhYWF8NRxWHyDPzRpYaEPqVSK5uZm5Ofn096dcTJVVVU477zzJCMjI2cSQurm8hyKooQAvgDwB0JI4/h9ewD8mfwS83b73Hgiw+68FwBFUXwAnwP44Bjh/h2AXwO4koxfFRFCzIQQxfjP1XAms82rDowQ4pBKpTfv3r37mdNPP12t1+tnf9IC4HK5KCgoQEJCAg4cOAClUumR47CwsCwcu92OpqYmdHd3Y926dR4V7p07dzrOOeecvpGRkXVzFW4AIISo4YxxnznDY9w+N57IsOI9Tyhnqu4bAFoIIc9Muv9MAH8BcB4hZGzS/bEURXHHf84AkA2geyHHlsvlj9TW1t5ZXFys6O3tdeNVzExiYiJWr16NlpYWNDU1wW5nq0BYWPwBlUqFsrIyCAQClJSUeKyPPiEEzz77rPHqq69ukslkqwkhPbM9Z/xcJxz/OQjAaQBaZ3k8LefGExFWvOfPegBXAziVoqja8dvZAF6EM/t81/h9r4w/fhOAeoqi6gB8BuAWQsiCt7Rqtfqdtra2s9etWyfZt2+fx2IeISEhEzG0/fv3s7twFhYfYrfb0dzcjKamJhQVFXmszSngnJt+1VVXaf75z39+J5fLVxNC5HN8ajyA3RRF1QOohDPm/S1FURdQFDUIoBTADoqido4/ntZz44kGG/NmKBRFJcbFxe16+OGHM2655RaPjrEyGAyora2FUCjEsmXL2PpfhsLGvJmJWq1GXV0dkpKSkJGR4dGmODKZDGeddZayr69vi0Kh2EJYgfBbWPFmMBRFBcXFxX1y7rnnbnr55ZfDPZlgRghBd3c3BgcHsWLFCkRGRnrsWCyegRVvZmG329He3g6FQoGCggKPd0Osra3FueeeK5fL5VebTKadsz+DxZewbnMGQwgxymSy87766qtnN27cqFIoFB47FkVRyMzMRFFR0YT7jo2Fs7B4BrVajbKyMvD5fKxfv97jwv3xxx9bzzjjjJ7BwcH1rHAzA3bnvUgIDg4+Ny4u7s0dO3bELF++3KPHcu3CBwYGkJubi7i4OI8ej4Ue2J23/2O1WtHa2gqNRuOV3bbD4cADDzygf+utt+pkMtk5hBC2UxNDYMV7EUFRVE5cXNzOl156Kf7iiy/2+OxFo9GIpqYmOBwO5OXlITg42NOHZHEDVrz9F0IIBgcH0dnZiczMTCQnJ3t84Itarcall16qrqur+1Amk/2RsMNFGAUr3osMiqIi4+LiPj/33HOLXnzxxXCBQODxY8rlcjQ1NSEhIQGZmZlsQpufwoq3f6LRaNDY2IiwsDDk5OR4pTlSZWUlLr74YrlCofiTXq/f5vEDstAOK96LEIqiqKioqD/HxcX95euvv47Ozs72+DEdDge6urowNDSEnJwcj/VXZlk4rHj7F5Nd5Pn5+YiIiPD4MQkheOqpp8a2bt3aP+4mZ+uqGQor3osYiqKKRCLRF1u3bhVfeeWVXul16nKl2+125Ofns650P4IVb/+AEIKBgQF0dXV5zUUOAAqFAhdffLGqubn5c5lMdjshxOLxg7J4DFa8FzkURYXHxcV9fNppp5W+9tprEd4SU7lcjubmZsTExCA7OxsBAQFeOS7L9LDi7XtGR0fR0tKCiIgIr7nIAWD//v3kt7/9rXx0dPRmo9H4pVcOyuJRWPE+AaAoioqMjLwtJibm4S+//DI6NzfXK8d1OBwYGBhAd3f3RIMJNh7uO1jx9h0ajQbNzc3g8XjIyclBaGioV47rcDjw2GOPGV566aVumUz2a0JIv1cOzOJxWPE+gaAoKi8uLu6bxx57LOHGG28M8IarDnA2m+jp6cHAwAAyMjKQnJwMDodtMeBtWPH2PgaDAa2trTCbzcjJyfFqcyOpVIoLL7xQ1dnZ+Z5MJruHEqpaDwAAESRJREFUEGLz2sFZPA4r3icYFEWFiESiN5cvX376+++/HxkfH++1Y1utVnR0dEAmk2HJkiWIj4/3SqyPxQkr3t7DZDKhvb0darUay5Yt82ovBEIItm3bZrvnnntGVSrVdSaT6X9eOziL12DF+wRFIBCcFRkZ+fqWLVtir7rqKr43RdRkMqGtrQ1arRbLli1DTEwMK+JegBVvz2O1WtHV1QWJRILs7GwkJCR49bMtk8lwzTXXqGpra/dJpdLrCCEqrx2cxauw4n0CQ1FUhEgkej0vL2/z+++/HykWi716fL1ej7a2NoyNjSErKwtisZgVcQ/CirfnMJlM6OrqgkwmQ3p6OlJSUrweGvroo49sd91116hKpbrJZDJ949WDs3gdVrxZIBAIzo6MjHz9qaeeirnyyiu9ugsHgLGxMXR0dECtViMjIwOJiYlsTNwDsOJNPwaDAZ2dnT797MpkMvzud79THTlyhN1tn0Cw4s0CYGIX/kZ+fv6p77//fqQvmqyYTCZ0d3dDKpUiNTUVqampbHY6jbDiTR9arRYdHR0+9xp9/PHHtrvuumtUrVbfPDY29rXXDWDxGax4sxyFQCA4JzIy8jVf7cIBZ9ywp6cHQ0NDSExMRHp6utfqYRczrHi7j1KpREdHBxwOB7KzsxEdHe0T0ZbJZLj22mtVNTU1ZVKp9HfsbvvEgxVvluMY34W/kpqaevrbb78dnZOT4xM77HY7+vr60N/fj6ioKKSnp3t8ytJihhXvhWG32zE8PIze3l4IBAJkZ2dDKBT6xBabzYYXXnjB9K9//Uup0WhuMxqNX/nEEBafw4o3y7RQFLVaJBK9e+mllyY+/vjjYb4STkIIpFIpenp6QAhBeno6m9y2AFjxnh9GoxG9vb2QSCQQi8VIS0tDUFCQz+w5cOAArr/++lGVSvW+XC5/kBAy5jNjWHwOK94sM0JRFCc8PPyWsLCwf2zZsiXqiiuu4PlSNPV6PXp6ejA6OoqEhASkpKT49ITKJFjxnh1CCORyOXp6emCxWJCWloaEhASf5l5IpVLccccd6v3797dIpdJrCCGdPjOGxW9gxZtlTlAUFS0SiV5ISUn51dtvv+21FqvTYbPZMDw8jL6+PgQEBCAtLQ2xsbFslvoMsOI9PUajEQMDAxgaGkJUVBTS0tK8MuVrJlwu8i1btijUavUfTCbTl4Q9YbOMw4o3y7ygKKpYJBK9d8kllyQ+8cQTPnOlT0aj0aCvrw8KhQIxMTFISkqCUChk3erHwIr30VitVoyMjGBwcBAOhwNJSUlITEz0i+TIsrIy3HDDDQqVSvWOXC7/O+siZzkWVrxZ5o3LlR4aGvqPv//978Ibb7wxgMfj+dosOBwOyGQyDA4OQq/XIz4+HklJSQgJCfG1aX4BK97+/xnp6urCXXfdpaqoqGhlXeQsM8GKN8uCoSgqIjY29uGQkJArn3rqqciLLrqI6y+7XavVColEgoGBAdjtdiQlJSEhIQGBgYG+Ns1nnKjiTQiBSqXCwMAAlEolYmNjkZSUhIiICL/xzkilUjzwwAOaHTt2SOVy+e12u/1HX9vE4t+w4s3iNhRFxYvF4qejo6N/9fzzz0efeuqp/nFGHMdkMmFoaAhDQ0PgcrkQi8UQi8V+s9vyFieSeNvtdigUCkgkEigUCkRERCApKQmxsbF+I9iAs9nL448/rn/nnXdUWq32PqPR+AkhxOFru1j8H1a8WWiDoqgssVj8Unp6etFLL70UvXLlSl+bdBxGoxESiQRSqRQmkwlxcXEQi8WIjIz0q5O6J1js4m2xWCCVSiGRSKDX6xETEwOxWIzo6Gi/S2Q0m8144YUXTFu3blUbjcYnNBrNK4QQq6/tYmEOrHiz0A5FUUVisfiV4uLijOeeey4qMzPT1yZNic1mg1wuh0QigVqthlAohFgsRmxsLPwhhk83i028CSHQ6/UTF2OEEIhEIojFYoSFhfnlxZjdbsd7771n/fvf/642Go3/VSgUT7LJaCwLgRVvFo/B5XJPjY2N/c/mzZvjHn744cisrCxfmzQtrrioRCLB6OgoKIpCdHQ0YmJiEBUVtSjEnOniTQiBwWCAQqHA6OgotFotQkJCIBKJIBKJIBAIfG3itNhsNnz00Uf2hx56SKXX6z+TSqUPEkKUvraLhbmw4s3iUSiKojgczplxcXFbVq1alfD4449HFRYW+tqsWbFarRMioVAowOFwGC/mTBNvQgjGxsYwOjo6IdbBwcGIiYlBTEwMwsPD/XJ3PRmTyYQ33njD8q9//UtjNps/k8lkjxJCRnxtFwvzYcWbxWtQFLUhPj5+S0ZGRvaTTz4Zs3HjRl+bNGcsFgsUCsXEjcPhQCgUTtxCQ0P9Xkj8XbytVivUajU0Gg3UajV0Oh3jxNqFVqvFCy+8YHzppZd0FovldYVC8TQ7PISFTljxZvE6FEWtiI+P3xIXF1f02GOPRZ9zzjkUU07KLlxC47oZDAbweDxERERAKBQiIiLC7+Ku/iTek98/jUYDvV4PLpc78d4JhUK/e//mglwux5YtW/Tvv/++1mQyPatWq//DxrRZPAEr3iw+g6KozPj4+MdDQ0M3P/TQQ5GXXHIJ1x+6Wy2UY3eOLkEKCQlBaGjoxC0kJMQnvbK9Ld6EEJjNZuj1euj1ehgMBuj1ehiNxqMudJjiuZiJvr4+PPbYY5pvvvlGrdfrHzUYDO+y2eMsnoQVbxafQ1FUfFxc3N/5fP5FN9xwQ+jtt98eHBcX52uzaMFms02I1mQRs9vtCAwMnBD24OBgCAQCCAQCBAYGeqS0iW7xJoTAZrPBZDJN3Fyv1WAwgBCCwMDAoy5aQkNDERQUxGihdkEIwc8//4zHH398tLW1VaZQKP5hsVi2s3XaLN6AFW8Wv4GiqODQ0NCrQ0JC7ispKRE+8MADUatXr14UJ/pjIYTAYrFMiN3Y2NiEAJrNZri+l4GBgUeJekBAAHg83nE3LpcLHo8343s1k3g7HA7YbLZpb2az+Sj7bDYbAIDH403YJxAIJgQ6JCTE72qr6UKv1+Ptt9+2PPPMM1qz2bxveHj4MULIEV/bxXJiwYo3i98xHgBfn5CQ8FBoaOiKu+66S3jVVVfxQ0NDfW2aV3G5nV2CaTKZYLFYJgTVbrcfJ7IzodfrMd17yOFwjroIOPY2+SJCIBAwMtveXerq6vDMM8+od+7cOWa1Wl9VKpUvEkIUvraL5cSEFW8Wv4aiKFFkZORtAQEBN55xxhnBd999t7CgoMDXZjESf0pYYwpjY2P46KOPbFu3blWp1erW4eHhxwHsYl3jLL6GFW8WRkBRFAfAaQkJCQ8EBQXlXn/99aFXX311UHJysq9NYwyseM8Nm82Gn3/+Ga+88ori4MGDJrvd/t7o6OgLhJBhX9vGwuKCFW8WxkFRVGRwcPBlERERt4vFYtGtt94aeckll/CEQqGvTfNrWPGeHkIIqqur8frrr2u++uorM0VRu0ZGRl4CUE7YkySLH8KKNwujoSgqJTIy8rqAgIBr8/LyQm+77bboc845hzqRR39OByvex9Pd3Y133nln7N133zVYrda6oaGh5wD8wJZ5sfg7J17WCcuighDSD+BhAA9TFLWisbHxFgAXnHLKKQFXX3111KmnnurXPa9ZvE9PTw++/PJLy5tvvqlRKpWDSqXyBZPJ9BkhROdr21hY5gq782ZZdIzHxzeKxeKrCSFn5ubmBlx11VVRv/71r7mLpX58IZyoO2+Hw4GKigp89tln+u3btxstFkufSqV6d2xs7DO2zzgLU2HFm2VRM152tiwyMvISgUBweVRUVMwVV1wResEFFwTl5OQsyhry6TiRxHtsbAw//vgjPvzwQ8XevXvtXC738MjIyNsOh+MHQoje1/axsLgLK94sJxQURcUGBAT8Oi4u7lqKonLOOecc/nnnnSfcsGEDwsLCfG2eR1nM4k0IQXt7O3766Sfbtm3bVJ2dnUZCyLdSqfQDAIcJIXZf28jCQieseLOcsFAUFQjgZLFYfBGAU4VCYfgZZ5wReOaZZ4avX79+0Yn5YhJvl1j//PPPtm+++UZVW1tLuFxuq0ql+spgMHxNCOn0tY0sLJ6EFW8WlnEoiooAsEEsFp+PRSjmTBZvQgg6OjqOEmsOh9OqUqm+NBgMuwA0s41TWE4kWPFmYZmGY8U8NDQ0vLi4mLNp0ybh6tWrufn5+WBSSRqTxHtkZATV1dU4dOiQYf/+/YbOzk6M76xZsWZhASveLCxzhqKoYAAFQUFBa6OjozdbrdYVoaGhQccKur+WpvmreA8PD08I9b59+wzd3d3gcDjDNpvtgFQq3QOgGkAv2yyFheUXWPFmYXEDiqKCAKwQCARrYmJiNttstoLg4OCQpUuXYsWKFUH5+fmhS5YsQXZ2NnzdAc6X4m2329Hf34/29na0trba6urqtA0NDbahoSFwOJxBq9VaJpPJ9sIp1P2sULOwzAwr3iwsNDOeCJcBYEl4eHh+REREkd1uXwogUigUcpcvX04VFBSE5ubmCtLT05GQkICYmBiPj9D0tHgbjUaMjIxgcHAQ7e3tpL6+XltfX2/u7e2lzGazic/nD9hstgaZTFZlt9vbALQDkLFCzcIyf1jxZmHxIhRFhQHIBrAkNjZ2pUAgyHM4HAk2my2az+fzBQIBNz4+niQnJ3PT0tICU1NTQxISEqj4+HgIhUKEhYUhLCwMAoFg3jXq8xVvu90OvV4PnU4HnU4HmUzmEmdLX1/fWG9vr3VwcJAolUrKarVaCSFjfD5fQgjpUygUVUajsRlOgR5gS7VYWOiFFW8WFj+CoqgAAGIA8QDi+Xx+YlRUVHZAQEAagEhCSKjdbg91OBxBXC6Xw+FwuBwOh8vn86nQ0FBHeHg4QkNDwePxKA6HAy6XCw6HQ3G5XOh0urCgoCCtw+GA3W6Hw+GAyWQiWq0WOp2OGhsboxwOh3385iCE2DgcjoHL5eoA6BwOx4jBYOhWq9XdAEYm3VTs7pmFxbuw4s3CsgigKIoHIBRAGIAQAJxjblwAjiluVgC68ZuRFWEWFmbAijcLCwsLCwvD8GyGDAsLCwsLCwvtsOLNwsLCwsLCMFjxZmFhYWFhYRiseLOwsLCwsDAMVrxZWFhYWFgYBiveLCwsLCwsDIMVbxYWFhYWFobBijcLC0OhKCqZoqjdFEW1UBTVRFHUn8bv/5iiqNrxWy9FUbWTnnM/RVGdFEW1URR1xqT7T6YoqoqiqC0+eCksLCzzhOdrA1hYWBaMDcA9hJCa8Z7p1RRF7SKEXOZ6AEVRWwFoxn/OBXA5gOUAEgD8SFHUkvG+47cC2AjgMYqilhFCWr39YlhYWOYOu/NmYWEohJARQkjN+M86AC0AEl1/p5yTSy4FsG38rvMBfEQIMRNCegB0Algz/jcOAAJny9T5TTxhYWHxOqx4s7AsAiiKSgOwEsDhSXdvBCAlhHSM/54IYGDS3wfxi9i/DuAgAA4hpMWz1rKwsLgL6zZnYWE4FEWFAvgcwJ2EEO2kP12BX3bdwNQ7agIAhJCdAHZ6zEgWFhZaYcWbhYXBUBTFh1O4PyCEbJ90Pw/AhQCKJj18EEDypN+TAAx7w04WFhZ6Yd3mLCwMZTym/QaAFkLIM8f8+TQArYSQwUn3fQ3gcoqiAimKSgeQDaDCO9aysLDQCbvzZmFhLusBXA2gYVI52AOEkO/gzCqf7DIHIaSJoqhPADTDmal++3imOQsLC8Ng53mzsLCwsLAwDNZtzsLCwsLCwjBY8WZhYWFhYWEYrHizsLCwsLAwDFa8WVhYWFhYGAYr3iwsLCwsLAyDFW8WFhYWFhaGwYo3CwsLCwsLw2DFm4WFhYWFhWH8PzeyBtn4C7DVAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "%matplotlib inline\n", "from matplotlib import pyplot as plt\n", "\n", "fig = plt.figure(figsize=(8, 8))\n", "\n", "plt.polar(theta, db_pattern)\n", "\n", "plt.ylim(db_pattern.min()-1, db_pattern.max()+1)\n", "plt.title('RCS (dB)')\n", "\n", "plt.show()" ] }, { "cell_type": "markdown", "id": "66e8b790", "metadata": {}, "source": [ "### Further information\n", "\n", "More information about the OSRC-preconditioned Burton-Miller formulation, including convergence plots and combination with fast matrix algebra, can be found in Betcke, van 't Wout & Gélat (2017).\n", "\n", "The OSRC preconditioner can also be used to speed up acoustic BEM simulations for multiple and penetrable domains. Benchmark results for various boundary integral formulations are reported in van 't Wout et al. (2021).\n", "More information about the OSRC-preconditioned PMCHWT formulation, with application to focused ultrasound propagation in the human body, can be found in Haqshenas et al. (2021)." ] } ], "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.10" } }, "nbformat": 4, "nbformat_minor": 5 }