{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Shan-Chen Two-Component Lattice Boltzmann" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from lbmpy.session import *\n", "from lbmpy.updatekernels import create_stream_pull_with_output_kernel\n", "from lbmpy.macroscopic_value_kernels import macroscopic_values_getter, macroscopic_values_setter\n", "from lbmpy.maxwellian_equilibrium import get_weights" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This is based on section 9.3.3 of Krüger et al.'s \"The Lattice Boltzmann Method\", Springer 2017 (http://www.lbmbook.com).\n", "Sample code is available at [https://github.com/lbm-principles-practice/code/](https://github.com/lbm-principles-practice/code/blob/master/chapter9/shanchen.cpp)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Parameters" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "N = 64 # domain size\n", "omega_a = 1. # relaxation rate of first component\n", "omega_b = 1. # relaxation rate of second component\n", "\n", "# interaction strength\n", "g_aa = 0.\n", "g_ab = g_ba = 6.\n", "g_bb = 0.\n", "\n", "rho0 = 1.\n", "\n", "stencil = LBStencil(Stencil.D2Q9)\n", "weights = get_weights(stencil, c_s_sq=sp.Rational(1, 3))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Data structures\n", "\n", "We allocate two sets of PDF's, one for each phase. Additionally, for each phase there is one field to store its density and velocity.\n", "\n", "To run the simulation on GPU, change the `default_target` to gpu" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "dim = stencil.D\n", "dh = ps.create_data_handling((N, ) * dim, periodicity=True, default_target=ps.Target.CPU)\n", "\n", "src_a = dh.add_array('src_a', values_per_cell=stencil.Q)\n", "dst_a = dh.add_array_like('dst_a', 'src_a')\n", "\n", "src_b = dh.add_array('src_b', values_per_cell=stencil.Q)\n", "dst_b = dh.add_array_like('dst_b', 'src_b')\n", "\n", "ρ_a = dh.add_array('rho_a')\n", "ρ_b = dh.add_array('rho_b')\n", "u_a = dh.add_array('u_a', values_per_cell=stencil.D)\n", "u_b = dh.add_array('u_b', values_per_cell=stencil.D)\n", "u_bary = dh.add_array_like('u_bary', u_a.name)\n", "\n", "f_a = dh.add_array('f_a', values_per_cell=stencil.D)\n", "f_b = dh.add_array_like('f_b', f_a.name)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Force & combined velocity\n", "\n", "The two LB methods are coupled using a force term. Its symbolic representation is created in the next cells.\n", "The force value is not written to a field, but directly evaluated inside the collision kernel." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The force between the two components is\n", "$\\mathbf{F}_k(\\mathbf{x})=-\\psi(\\rho_k(\\mathbf{x}))\\sum\\limits_{k^\\prime\\in\\{A,B\\}}g_{kk^\\prime}\\sum\\limits_{i=1}^{q}w_i\\psi(\\rho_{k^\\prime}(\\mathbf{x}+\\mathbf{c}_i))\\mathbf{c}_i$\n", "for $k\\in\\{A,B\\}$\n", "and with \n", "$\\psi(\\rho)=\\rho_0\\left[1-\\exp(-\\rho/\\rho_0)\\right]$." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "def psi(dens):\n", " return rho0 * (1. - sp.exp(-dens / rho0));" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "zero_vec = sp.Matrix([0] * dh.dim) \n", "\n", "force_a = zero_vec\n", "for factor, ρ in zip([g_aa, g_ab], [ρ_a, ρ_b]):\n", " force_a += sum((psi(ρ[d]) * w_d * sp.Matrix(d)\n", " for d, w_d in zip(stencil, weights)), \n", " zero_vec) * psi(ρ_a.center) * -1 * factor\n", "\n", "force_b = zero_vec\n", "for factor, ρ in zip([g_ba, g_bb], [ρ_a, ρ_b]):\n", " force_b += sum((psi(ρ[d]) * w_d * sp.Matrix(d)\n", " for d, w_d in zip(stencil, weights)), \n", " zero_vec) * psi(ρ_b.center) * -1 * factor\n", " \n", "f_expressions = ps.AssignmentCollection([\n", " ps.Assignment(f_a.center_vector, force_a),\n", " ps.Assignment(f_b.center_vector, force_b)\n", "])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The barycentric velocity, which is used in place of the individual components' velocities in the equilibrium distribution and Guo force term, is\n", "$\\vec{u}=\\frac{1}{\\rho_a+\\rho_b}\\left(\\rho_a\\vec{u}_a+\\frac{1}{2}\\vec{F}_a+\\rho_b\\vec{u}_b+\\frac{1}{2}\\vec{F}_b\\right)$." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "u_full = list(ps.Assignment(u_bary.center_vector,\n", " (ρ_a.center * u_a.center_vector + ρ_b.center * u_b.center_vector) / (ρ_a.center + ρ_b.center)))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Kernels" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "lbm_config_a = LBMConfig(stencil=stencil, relaxation_rate=omega_a, compressible=True,\n", " velocity_input=u_bary, density_input=ρ_a, force_model=ForceModel.GUO, \n", " force=f_a, kernel_type='collide_only')\n", "\n", "lbm_config_b = LBMConfig(stencil=stencil, relaxation_rate=omega_b, compressible=True,\n", " velocity_input=u_bary, density_input=ρ_b, force_model=ForceModel.GUO, \n", " force=f_b, kernel_type='collide_only')\n", "\n", "\n", "\n", "collision_a = create_lb_update_rule(lbm_config=lbm_config_a,\n", " optimization={'symbolic_field': src_a})\n", "\n", "collision_b = create_lb_update_rule(lbm_config=lbm_config_b,\n", " optimization={'symbolic_field': src_b})" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [], "source": [ "stream_a = create_stream_pull_with_output_kernel(collision_a.method, src_a, dst_a, \n", " {'density': ρ_a, 'velocity': u_a})\n", "stream_b = create_stream_pull_with_output_kernel(collision_b.method, src_b, dst_b, \n", " {'density': ρ_b, 'velocity': u_b})\n", "\n", "config = ps.CreateKernelConfig(target=dh.default_target)\n", "\n", "stream_a_kernel = ps.create_kernel(stream_a, config=config).compile()\n", "stream_b_kernel = ps.create_kernel(stream_b, config=config).compile()\n", "collision_a_kernel = ps.create_kernel(collision_a, config=config).compile()\n", "collision_b_kernel = ps.create_kernel(collision_b, config=config).compile()\n", "\n", "force_kernel = ps.create_kernel(f_expressions, config=config).compile()\n", "u_kernel = ps.create_kernel(u_full, config=config).compile()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Initialization" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [], "source": [ "init_a = macroscopic_values_setter(collision_a.method, velocity=(0, 0), \n", " pdfs=src_a.center_vector, density=ρ_a.center)\n", "init_b = macroscopic_values_setter(collision_b.method, velocity=(0, 0), \n", " pdfs=src_b.center_vector, density=ρ_b.center)\n", "init_a_kernel = ps.create_kernel(init_a, ghost_layers=0).compile()\n", "init_b_kernel = ps.create_kernel(init_b, ghost_layers=0).compile()\n", "\n", "\n", "dh.run_kernel(init_a_kernel)\n", "dh.run_kernel(init_b_kernel)" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [], "source": [ "def init():\n", " dh.fill(ρ_a.name, 0.1, slice_obj=ps.make_slice[:, :0.5])\n", " dh.fill(ρ_a.name, 0.9, slice_obj=ps.make_slice[:, 0.5:])\n", "\n", " dh.fill(ρ_b.name, 0.9, slice_obj=ps.make_slice[:, :0.5])\n", " dh.fill(ρ_b.name, 0.1, slice_obj=ps.make_slice[:, 0.5:])\n", " \n", " dh.fill(f_a.name, 0.0)\n", " dh.fill(f_b.name, 0.0)\n", "\n", " dh.run_kernel(init_a_kernel)\n", " dh.run_kernel(init_b_kernel)\n", " \n", " dh.fill(u_a.name, 0.0)\n", " dh.fill(u_b.name, 0.0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Timeloop" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [], "source": [ "sync_pdfs = dh.synchronization_function([src_a.name, src_b.name])\n", "sync_ρs = dh.synchronization_function([ρ_a.name, ρ_b.name])\n", "\n", "def time_loop(steps):\n", " dh.all_to_gpu()\n", " for i in range(steps):\n", " sync_ρs() # force values depend on neighboring ρ's\n", " dh.run_kernel(force_kernel)\n", " dh.run_kernel(u_kernel)\n", " dh.run_kernel(collision_a_kernel)\n", " dh.run_kernel(collision_b_kernel)\n", " \n", " sync_pdfs()\n", " dh.run_kernel(stream_a_kernel)\n", " dh.run_kernel(stream_b_kernel)\n", " \n", " dh.swap(src_a.name, dst_a.name)\n", " dh.swap(src_b.name, dst_b.name)\n", " dh.all_to_cpu()" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [], "source": [ "def plot_ρs():\n", " plt.figure(dpi=200)\n", " plt.subplot(1,2,1)\n", " plt.title(\"$\\\\rho_A$\")\n", " plt.scalar_field(dh.gather_array(ρ_a.name), vmin=0, vmax=2)\n", " plt.colorbar()\n", " plt.subplot(1,2,2)\n", " plt.title(\"$\\\\rho_B$\")\n", " plt.scalar_field(dh.gather_array(ρ_b.name), vmin=0, vmax=2)\n", " plt.colorbar()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Run the simulation\n", "### Initial state" ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAACeQAAAQKCAYAAAArG+jpAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjQuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/MnkTPAAAACXBIWXMAAB7CAAAewgFu0HU+AACX2ElEQVR4nOzdfZhtZ1kf/u99MiSE12h4CSTRAFETFAFNKAgSFYWgUQIVUFSMEBFrqeIbtGCJoK0itmjAipUS5IcgIBBMQcQXIkTBxOILJqAi0XMICQFMCElIzOzn98esc85mmJk9M3vt2XvNfD7Xta79rL2e9ax7OLmutt/eez3VWgsAAAAAAAAAAAAwnX3zLgAAAAAAAAAAAAB2Aw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAQJKkqu5UVc+rqj+vqk9U1S3dcWVVvbKqvnLeNQIAAAAAzIJ8FACAvlRrbd41AABzVlWPTfK/kxy/wbTPJfm+1tqbdqYqAAAAAIDZk48CANAnDXkAsMdV1XlJXjD21T8kuTTJLUm+OsnXjl27KcnXtNY+vGMFAgAAAADMiHwUAIC+acgDgD2sqv5rkp/tTg8k+cHW2u+vmnNWkjclOar76oLW2g/sXJUAAAAAAP2TjwIAMAsa8gBgj6qqb0ryh0kqyb8k+brW2sfWmfszSV7Ynd6Q5JjW2mhHCgUAAAAA6Jl8FACAWdGQBwB7UFUdkeTyJF+epGUlbHrfBvNPSLJ/7KtTbMsAAAAAAAyRfBQAgFnaN+8CAIC5eFJWwqYk+Z2NwqYkaa0dSPKZsa/uOqvCAAAAAABmTD4KAMDMaMgDgL3pB8bGL9/kPf82Np74it2qentVte54zJaqAwAAAACYnV7z0ao6dSwLXX3cVFVXVNX/qKp7TF05AAALb2neBQAAO6uq7pzkm7rTTya5ZBP3LOXzf/W5f7253fzHJnns2FcPSvLOLRUKAAAAANCzGeWjX9N93pTkA2PfH5XkS5Kc0h2Pr6oHtdau30bpAAAMhDfkAcDe83VJjujG722tTXzbXZKTc7iR/4Yk16w3sapul+R/dqfXdZ8P3HqZAAAAAAC9m0U+erAh749ba48YO05PclySH+mun5TkKduuHACAQdCQBwB7z78bG1+xyXu+fmz8pxNCqmcl+YokH0ny37rvHrTp6gAAAAAAZmcW+ejBhrwPrPo+bcWv5XAT30mbfCYAAANly1oA2HsePDb+xCbvefzYeN2tZ6vqbkn+a3f6U0kObr3w5VV1dGvt5k1XCQAAAADQv1nkow/qPr+gIW/MqPv8+CafCQDAQHlDHgDsPeOB0+0mTa6qk5M8pjv9XJLXbjD955PcNcnFrbW3JPnb7vsjknzV1ksFAAAAAOhVr/loVd0vyTHd6ZoNeVX1mCT3SnJLkjdvrVwAAIZGQx4A7CFV9cVJvnTsq/tv4rZfzuH/O8OrWmufXmftByY5Nyu/9Hx2krTWrs3hX5k+cDs1AwAAAAD0YUb56MHtaq9rrV059qyjqurLq+q/JnlTkluTnNta+5dtFQ8AwGBoyAOAveXBq84fX1VftN7kqnpuku/oTj+d5Gc2WPtXsvJ/t3h1a238l6AH35L3oK2VCgAAAADQq1nkowcb8o6pqnbwyMrb9D6c5Lwkb0zy71pr/980xQMAMAwa8gBgbxkPnG7Nyvayr6+qY8YnVdVdq+rlSf5799Vyku9rrX1qrUWr6juTnJHkxiTPW3X5g92nN+QBAAAAAPM0i3z04Jr/nOSSseP/ZaWJr5L8+yQP6+lvAABgwS3NuwAAYEeNB07PS/ILSR6d5J+r6o+SfDLJCVlprrtDN285yTNba29fa8Gqun2SX+pOf6G19vFVUw6+Ie+rq6paa236PwMAAAAAYMt6z0fH1nz+6jfgVdURSZ6dlfz016rq71prf9rLXwIAwMIq/3/iALB3VNXlSU7tTk/MynYL52f9t+YeSPKDrbXf32DN5yd5UZL9Sb6itXbzqusPSfL+7vR+rbV/2v5fAAAAAACwPX3no1V1QlZy0ST56tba364z74NJvjLJr7XWfmSb5QMAMBC2rAWAPaKq7pDky7vTT7fWDrTWfi3JI5O8MStbKtya5NokFyf5T1lpsNuoGe/4JM/tTk9MclNVtfEjh5vxkuRBff5NAAAAAACbMYt8NMnXdJ+3JLlig3kf7T7vuc3yAQAYEFvWAsDe8cAkR3Tjvz74ZWvtkiSXbHPNX0hyxySfTXL9BvPunuTIroY3b/NZAAAAAADbNYt89GBD3t+11m7bYN4J3efV23wOAAADoiEPAPaOB4+N/2raxarq3yX5niQtyWNaa3+2wdw3JHlivCEPAAAAAJiPXvPRzsGGvL9eb0KXoz6wO317T88FAGCB2bIWAPaO8cBp3YBoM6qqkvxqkkry2o2a8Tof7D4fuOEsAAAAAIDZ6C0fHXOwIe+vVl+oqn1V9dQkF2YlR31na01DHgDAHuANeQCwd/QZOD01yUOS3JjkOZuYf7Ah70ur6pjW2nVTPh8AAAAAYCt6bcirqrsnOb47/eGqetLY5S9Kcp8kR3fnb0nyfdM+EwCAYdCQBwB7QFUtJfmq7vTfklw+xVp3SvLfu9P/1lq7ahO3fXBs/KAk797u8wEAAAAAtqLPfHTM14yNT1l17bNJrkzyviS/1Vp7dw/PAwBgIDTkAcDecP8kR3XjD7XWbt3uQq21zya59xbv+fusbMsAAAAAALDTestHD2qtvTMyTwAA1rBv3gUAADui1+0YAAAAAAAGRD4KAMCOqdbavGsAAAAAWDhV9TVJzkzy9VnZ3uoeWdne6qokf5bkla219/T8zO9K8gNJvjrJFyW5Osl7kry8tfa+Ta5xbJL/lOTsJCdl5a0dH03y1iS/2lr7VJ81AwAAAAC7j3x0+zTkAQAAAKxSVRcneeQmpr4mybnTbnlVVbdP8sYkZ60zZZTkvNbaiyasc3qSC5Pca50pVyV5XGvtsu3WCgAAAADsbvLR6WjIAwAAAFilqv4xyf2yEtC8MSu/wvyXJEckeViSn0hyfDf9da21p0z5vNcmObjGnyT5le7ZD0jyX7pakuQHW2u/uc4axyf5yyT3THJbkv+R5KLu8llJfjzJUpJrknxta+1j09QMAAAAAOxO8tHpaMgDAAAAWKWqLkryW0l+t7W2vMb1uyW5JMmXd189crvbM1TVGUne3Z3+XpLHjz+ze9ZfJvmSJP+a5L6ttevWWOeCJN/fnT6ptfbGVdefmOQN3emrWmtP2069AAAAAMDuJh+djoY8AAAAgG2oqrOyEhAlya+21n50m+v83yTfmmQ5yUmttQNrzPmuJK/rTn+ytfbLq67fM8nHsvIL1Xe21s5c51m/n+Qx3bOOb61ds52aAQAAAIC9TT66vn2zWBQAAABgD3j32Ph+603aSFXdKcmjutN3rRU2dd6c5DPd+AlrXP+OrIRNSfKqDR55Qfd5RHcPAAAAAMB2vHtsLB8doyEPAAAAYHuOHBuPtrnGQ5Ic1Y0vXm9Sa+3WJO87eE9V3W7VlK8fG6+7zqprj9hskQAAAAAAq8hH16EhDwAAAGB7zhgbf2iba5y6hTUOXl9K8mXrrHN9a+3q9RZorX08h39Jeup68wAAAAAAJpCPrmNpVgvvJlV1VJIHdKfXZmUfYQAAAPaOI5LcvRv/bWvtlnkW07eqWkpy3Lzr2Kbjson/t/oGWx1sS1XtS/Lcsa/esM2lThwbT6px/6r7Ll9jnc38nfuTfOWqZwOsSz4KAACw5+3qfDSRkW6VfHRjGvI25wFJLp13EQAAACyE05NcNu8ienZcPj/M2I2q5/WenZXtFJLkLa217f43ceex8WcnzL1xbHynddaZtMb4OqvXAFiPfBQAAICDdmM+mshIt0o+ugFb1gIAAABsQVWdkeQXutNPJPnhKZa7/dj41glzx395fPQ660xaY3yd1WsAAAAAAGxIPjqZN+RtzrUHByc880ezdOe7zLMWAAAAdthtN3wmB379Vw6eXrvR3KF73ztOzL3uccS8y5jo459YzkMfe+gHq6cnuXonnltVX5nkLVnJVG5J8qTW2jVTLPm5sfGRE+YeNTa+eY117rCJNcbXWb0GwHoO/Z99D73L43LUvjvMsxYAAAB22C2jm/K+z1x48HRX56OJjHQj8tHN0ZC3OYf2WF66812ydNdj5lgKAAAAc7Y8ecpw3eseR+SEe99u3mVs1dWttQOzfkhV3SfJHyT5oqz8d/DdrbWLp1z2hrHxpC0S7jg2Xr31wg1ZCZw2s83CwXU2s30DQDL2f/Ydte8Ouf0+O14DAADsYbs6H01kpOuRj26ehjwAAADgkFFaRhnNu4yJRmk7+ryquneSP0xy7yQtydNaa2/pYenxkOyEJJdtMPfEsfH+VdcOJLlnt8YkB9dZvQYAAAAA7Hky0i8kH92afbNaGAAAAGA3qKq7JXlXkvt2Xz2rtfZbPS1/+dj4lAlzD16/Lck/rrPOXavquPUWqKp7JblLd3rFZosEAAAAAPYm+ejWacgDAAAAWEdV3TXJO5Pcv/vqua21l/f4iEuT3NqNz9igjiOTPPTgPa21W1dNee/YeN11Vl27ZLNFAgAAAAB7j3x0ezTkAQAAAKyhqu6Q5P8m+Zruq59vrf1in89ord2Q5I+602+uqvW2VHhCDv9yc62tIN6WHNpH4wc2eOQ53eeouwcAAAAA4AvIR7dPQx4AAABwyHIbDeaYpe4Xl29J8vDuq19prT1/G+ucU1WtO85bZ9pLus+lJC+vqiNWrXG3JAeDruuS/ObqBVprVyd5bXf6mKr6zjVqeWKSx3Snr+nuAQAAAADGzDv3XISMVD46naVZLQwAAAAwYK9L8uhu/MdJXllVX7XB/Ftba3+/nQe11v64ql6f5LuSfEeSd1XVS5NcleQBSZ6X5Eu66c9trf3rOks9L8mZSe6e5HVVdVqSi7prZyX5iW58bZIth2cAAAAAwJ4hH52ChjwAAACAL/SEsfE3JfmbCfP/OclJUzzvaVnZcuFbk3xjd4wbJXlRa+0V6y3QWttfVd+e5K1JjkvynO4Yd3WSs1trB6aoFQAAAADY3eSjU9CQBwAAABwySssobd5lTDSEGreitXZzkm+rqqckOSfJA5Mck+SaJO9J8rLW2p9vYp33V9UDkvxokrNzOAT7aJILk7y0tfapnssHAAAAgF1DRrrzdls+qiEPAAAAYJXWWvW0zgVJLtjC/N9O8ttTPvOTSX6mOwAAAAAAtkQ+Op1983owAAAAAAAAAAAA7CYa8gAAAAAAAAAAAKAHtqwFAAAADmkZZZTRvMuYqA2gRgAAAABgeGSkTMsb8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB7asBQAAAA5Zbi3Lrc27jImGUCMAAAAAMDwyUqblDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANCDpXkXAAAAACyOUVpGafMuY6Ih1AgAAAAADI+MlGl5Qx4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9MCWtQAAAMAhoyTLA9jqYDTvAgAAAACAXUlGyrS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgMUxSssobd5lTDSEGgEAAACA4ZGRMi1vyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHtiyFgAAADhkubUst8Xf6mAINQIAAAAAwyMjZVrekAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD1YmncBAAAAwOIYdceiG0KNAAAAAMDwyEiZljfkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPbFkLAAAAHDJKy3LavMuYaDSAGgEAAACA4ZGRMq0de0NeVd2tqn66qi6pqqur6paquqqq3l9Vv1RVD9vEGmdW1Zur6kB3/4Hu/Myd+BsAAAAAALZLRgoAAACw++3IG/Kq6olJ/leSY1dduld3PCTJlyU5e537K8mvJ3nGqkvHJ3l8ksdX1W8keWZrTfsnAAAAALBQZKQAAAAAe8PMG/Kq6qlJXpWVt/F9Iiuh03uTfDrJcUnul+Tbk/zbBsv8XA4HTR9I8uIkH+nu/ekkD+6uX5vk+b3/EQAAAAAA2yQjBQAAANg7ZtqQV1WnJvmNrARN70ny7a2169eYen5VHbnOGidnJVBKksuSPLK1dnN3fmlVvS3JxUlOS/KcqnpVa+0jff4dAAAAsFcsJ1kewHuVluddAMAmyUgBAABgWGSkTGvfjNc/P8lRST6Z5AnrBE1JktbaretcenYONw4+ayxoOnjfTUme1Z0uJfmxaQoGAAAAAOiRjBQAAABgD5lZQ15VnZLkUd3py1prn9zGGpXkcd3ph1pr71trXvf9h7vTs7v7AAAAAADmRkYKAAAAsPfM8g15Txwbv/HgoKq+qKq+rKqO3cQa90lyfDe+eMLcg9dPSHLSZosEAAAAAJgRGSkAAADAHjPLhryHdp/XJ7miqr6nqv46yaeT/H2ST1bVP1XVC6rqTuuscerY+EMTnjd+/dR1ZwEAAADrGg3oABgAGSkAAAAMzLxzTxnp8C3NcO37d59XJjk/yY+sMec+Sc5L8p1V9ZjW2lWrrp84Nj4w4Xn717lvoqo6YcKU47ayHgAAAABABpKRykcBAAAA+jPLhrwv7j5PSfLAJNcleW6SNyf5TJIHJHlhkscm+aokb6yqr2+tjTdw3nls/NkJz7txbLzer0nXs3/yFAAAAACALRlKRiofBQAAAOjJLBvy7th9HpVkOcljW2vvG7t+WVWdleSirAROX5fkCUneNDbn9mPjWyc875ax8dHbqhgAAAD2uFEqy6l5lzHRaAA1AkRGCgAAAIMjI2Vas2zI+1wOB05vXBU0JUlaa6Oq+qmshE1J8t35/LDpc2PjIyc876ix8c1brHXS9g3HJbl0i2sCAAAAAHvbUDJS+SgAAABAT2bZkHdDDodN71hvUmvt76rqY0mOT3L6GmscNGmLhTuOjSdt3bC6hgMbXa/SUQoAAAAAbNkgMlL5KAAAAEB/9s1w7f1j4w0DnbG591j1/fh9J0xYY/xXnPvXnQUAAAAAsDNkpAAAAAB7zCzfkPd3OfxrziMmzD14/bZV318+Nj5lwhrj16+YMBcAAABYw6itHItuCDUCREYKAAAAgyMjZVqzfEPen46N7zdh7n27z4+t+v6jSa7qxmdMWOORY2tcOak4AAAAAIAZk5ECAAAA7DGzbMh7W5J/68ZPWG9SVZ2R5Nju9D3j11prLcmF3ekpVfXQddZ4aA7/+vPC7j4AAAAAgHmSkQIAAADsMTNryGutfSrJb3an31JV37V6TlXdOclLx756xRpLvTSHt2k4v6qOXrXG0UnO705vW7UeAAAAsAXLqcEcAItORgoAAADDM+/cU0Y6fLN8Q16SvCDJv3Tj11TV+VX1jVX1tVV1TpK/SPKg7vr/aq1dunqB1trfJ3lJd3pakkuq6slVdVpVPTnJJd33SfJLrbV/mNHfAgAAAACwVTJSAAAAgD1kaZaLt9auraozs7I1w8lJ/mN3rPZ/kvzoBks9L8k9kjwtyYOTvH6NOa9M8vypCgYAAAAA6JGMFAAAAGBvmfUb8tJauyIrv/D8qSTvT/LpJLcmOZDkd5J8U2vt6a21f9tgjVFr7elJvi3JhUmu6ta4qjv/1tbaua210Sz/FgAAAACArZKRAgAAAOwdM31D3kGttRuzsqXCSybNnbDO25O8vZeiAAAAgC+wnMpyat5lTDSEGgHGyUgBAABgGGSkTGvmb8gDAAAAAAAAAACAvUBDHgAAAAAAAAAAAPRgR7asBQAAAIahtWTUFn+rg9bmXQEAAAAAsBvJSJmWN+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4llNZTs27jImGUCMAAAAAMDwyUqblDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0ANb1gIAAACHLGdflgfw+70h1AgAAAAADI+MlGn5lwEAAAAAAAAAAIAeaMgDAAAAAAAAAACAHmjIAwAAAAAAAAAAgB4szbsAAAAAYHG0Vhm1mncZE7UB1AgAAAAADI+MlGl5Qx4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9MCWtQAAAMAhy6ksZ/G3OhhCjQAAAADA8MhImZY35AEAAAAAAAAAAEAPNOQBAAAAAAAAAABADzTkAQAAAAAAAAAAQA+W5l0AAAAAsDiW274st8X//d4QagQAAAAAhkdGyrT8ywAAAAAAAAAAAEAPNOQBAAAAAAAAAABAD2xZCwAAABwySmU0gN/vjVLzLgEAAAAA2IVkpExr8f/rAQAAAAAAAAAAgAHQkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPViadwEAAADA4hilspyadxkTjQZQIwAAAAAwPDJSpuUNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0IOleRcAAAAALI7lti/LbfF/vzeEGgEAAACA4ZGRMi3/MgAAAAAAAAAAANADDXkAAAAAAAAAAADQA1vWAgAAAIeMUhml5l3GREOoEQAAAAAYHhkp0/KGPAAAAAAAAAAAAOiBhjwAAAAAAAAAAADogYY8AAAAAAAAAAAA6MHSvAsAAAAAFsco+7I8gN/vjQZQIwAAAAAwPDJSpuVfBgAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6YMtaAAAA4JDlti/LbfF/vzeEGgEAAACA4ZGRMi3/MgAAAABrqKp7VNVZVfXCqnpHVX2yqlp3XNDTM75hbM3NHu9eZ60rN3n/lX3UDgAAAADsXvLR7fOGPAAAAIC1XTPvAtbx4XkXAAAAAADsevLRbdKQBwAAADDZ/iRXJHl0z+temuQBm5j3siRndONXT5h7YZLnb3D91k08DwAAAADgIPnoFmjIAwAAAA4ZpTLKvnmXMdEotROPeWFWAqFLW2vXVNVJST7a5wNaazcm+eBGc6rqmCQP7U7/sbX2ZxOWva61tuGaAAAAAMDaZKSHyEe3SUMeAAAAwBpaay+Ydw2dJyc5qhu/Zp6FAAAAAAB7g3x0+xa/nRMAAABgb3tq99kykMAJAAAAAKAng8tHvSEPAAAAOGTUKsttR7aDncpoADX2oarul+TrutP3tNZ63RICAAAAAPh8MtLFMdR81BvyAAAAABbXU8fGr97kPY+sqr+pqhur6qaq+mhV/U5VnV1Vuz+lAwAAAAB2i0Hmo96QBwAAAAzdcZNylNbagR2qpW/f233enORNm7znPqvOT+qOJyW5pKqe3Fr7WC/VAQAAAACLYLdmpIPMRzXkAQAAAEN36SbmDO7NcFX19Unu252+pbX2mQm33JrkbUn+IMkHk1yf5JgkD0vyw0lOTPLwJO+qqoe11q6fRd0AAAAAwI7bdRnpkPNRDXkAAADAIcvZl+Xsm3cZEw2hxh5839j4tzYx/yGttevW+P7dVfWyrPyC9NFJTk3ygiQ/PnWFAAAAALDLyEgXxmDzUQ15AAAAwNCdnuTqeRfRp6o6KskTu9OrkvzhpHvWCZsOXruhqp6U5CNJjk3yjKp6bmvt1h7KBQAAAADma1dlpEPPRzXkAQAAAEN3dWvtwLyL6NnjsrKdQpK8trW2PO2CrbXrq+r1SX4kyR2TnJbkz6ZdFwAAAACYu92WkQ46H9WQBwAAABwyavsyaou/1cEQapzSU8fGm9mOYbMuHxsf3+O6AAAAALAryEgXwqDz0V39LwMAAAAwNFV1jySP6U7/X2vtg30u3+NaAAAAAAC92g35qIY8AAAAgMXylBze1aDPX38myf3Hxlf1vDYAAAAAwLQGn49qyAMAAABYLAe3Y7gtyW/3tWhV3TXJk7vTm5Jc1tfaAAAAAAA9GXw+qiEPAAAAOGQ5+wZzDEFVnVNVrTvO28T8r0zy4O70Ha21azf5nDOr6ugNrt85yRuSHNt99crW2i2bWRsAAAAA9pJ55567KSPdq/no0uQpAAAAAHtPVT0iycljX91tbHxyVZ0zPr+1dkEPj/3+sfGrt3Dfc5O8tqrenOS9ST6S5LNJjknysCQ/nOTEbu6Hk5w3baEAAAAAwO4lH90+DXkAAAAAazs3nx8AjXt4d4y7YJqHVdW+JE/pTv81yUVbXOKLs1LzuRvM+dMkT2mtfXrrFQIAAAAAe4h8dJs05AEAAACHjJIst5p3GRON5l3AbDwqyfHd+He2uGXCT3b3PyzJV2Tl16rHJLkpyVVJ3p/kdUn+oLXW+ioYAAAAAHYbGenc7Jp8VEMeAAAAwBpaa+ckOWfKNS7IJn8Z2lp7V5JtJX2ttcuSXLadewEAAAAAVpOPbt++eRcAAAAAAAAAAAAAu4GGPAAAAAAAAAAAAOiBLWsBAACAQ0bZl9EAfr83hBoBAAAAgOGRkTIt/zIAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANCDpXkXAAAAACyO5bYvy23xf783hBoBAAAAgOGRkTIt/zIAAAAAAAAAAADQAw15AAAAAAAAAAAA0ANb1gIAAACHjFIZpeZdxkRDqBEAAAAAGB4ZKdPyhjwAAAAAAAAAAADogYY8AAAAAAAAAAAA6IGGPAAAAAAAAAAAAOjB0rwLAAAAABbHqO3Lclv83++NBlAjAAAAADA8MlKm5V8GAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHpgy1oAAADgkOVUlgfw+73l1LxLAAAAAAB2IRkp01r8/3oAAAAAAAAAAABgADTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4Rq0yajXvMiYaQo0AAAAAwPDISJmWN+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA9sWQsAAAAcMsq+LA/g93ujAdQIAAAAAAyPjJRp+ZcBAAAAAAAAAACAHmjIAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeLM27AAAAAGBxjNq+jNri/35vCDUCAAAAAMMjI2Va/mUAAAAAAAAAAACgBxryAAAAAAAAAAAAoAe2rAUAAAAOWU5lOTXvMiYaQo0AAAAAwPDISJmWN+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4Rm1fRm3xf783hBoBAAAAgOGRkTIt/zIAAAAAAAAAAADQAw15AAAAAAAAAAAA0ANb1gIAAACHjJIsp+ZdxkSjeRcAAAAAAOxKMlKm5Q15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQg6V5FwAAAAAsjlHbl1Fb/N/vDaFGAAAAAGB4ZKRMy78MAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPTAlrUAAADAIcttX5YHsNXBEGoEAAAAAIZHRsq0/MsAAAAAAAAAAABADzTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4Wiqj1LzLmKgNoEYAAAAAYHhkpEzLG/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHS/MuAAAAAFgcy21fltvi/35vCDUCAAAAAMMjI2Va/mUAAAAAAAAAAACgBxryAAAAAAAAAAAAoAe2rAUAAAAOGaUyajXvMiYaZfFrBAAAAACGR0bKtGb6hryqaps83r2Jtc6sqjdX1YGquqX7fHNVnTnLvwEAAAAAYLtkpAAAAAB7y8JvWVsrXpHkHUken+T4JEd2n49P8o6qekVVafsEAAAAAHYdGSkAAADAcOzUlrX/K8mvbXD9xg2u/VySZ3TjDyR5cZKPJLlfkp9O8uDu+rVJnj91pQAAAAAA/ZORAgAAAOwBO9WQ94nW2ge3elNVnZyVQClJLkvyyNbazd35pVX1tiQXJzktyXOq6lWttY/0UjEAAADsQcupLC/+C/WzHC+BAgZHRgoAAAADICNlWov+X8+zc7hp8FljQVOSpLV2U5JndadLSX5s50oDAAAAAJg5GSkAAADAgCxsQ15VVZLHdacfaq29b6153fcf7k7P7u4DAAAAABg0GSkAAADA8OzUlrXbcZ8kx3fjiyfMvTjJVyQ5IclJST46u7IAAABg92qtMmqL38fRBlAjQA9kpAAAALDDZKRMa6fekPfEqvpwVd1cVTdU1T9U1aur6hs3uOfUsfGHJqw/fv3UdWcBAAAAAMyHjBQAAABgD9ipN+Tdf9X5yd3x1Kp6a5JzWmvXr5pz4tj4wIT1969z36ZU1QkTphy31TUBAAAAAMYsbEYqHwUAAADoz6wb8m5K8rYkf5SVX2h+Nsndk5yR5JlJjk1ydpILq+pbWmv/NnbvncfGn53wnBvHxnfaRp37J08BAAAAANiyIWSk8lEAAACAnsy6Ie/41tp1a3z/rqo6P8k7kjw4K+HTDyf51bE5tx8b3zrhObeMjY/eRp0AAABAklH2ZZR98y5joiHUCNCRkQIAAMCAyEiZ1kwb8tYJmg5eu6aqvjPJFUmOTPKsfH7Y9Lmx8ZETHnXU2PjmLZaZTN7C4bgkl25jXQAAAABgDxtIRiofBQAAAOjJrN+Qt6HW2j9V1buSfFuSk6vq3q21q7rLN4xNnbTFwh3HxpO2blirjgMbXa+qrS4JAAAAADDRImSk8lEAAACA/sy1Ia9zeVbCpiQ5PsnBsGk8BDphwhrjv+Dc31NdAAAAsOcst8pyW/zGiyHUCLAFMlIAAABYEDJSprUImwmv91/H5WPjUyasMX79iunKAQAAAADYUTJSAAAAgF1iERry7j82vmps/NGx8zMmrPHI7vNjSa7spywAAAAAgB0hIwUAAADYJebakFdV903yLd3pP7XWPnbwWmutJbmwOz2lqh66zhoPzeFff17Y3QcAAAAAsPBkpAAAAAC7y8wa8qrq26tqaYPr90zypiS36756+RrTXprktm58flUdvWqNo5Oc353e1s0HAAAAtmnUajAHwKKTkQIAAMDwzDv3lJEO37phUA/OT3K7qvrdJH+elW0Sbk5ytyTfkOSZSY7t5r43a4RNrbW/r6qXJHluktOSXFJVv5jkI0nul+Q5SR7cTf+l1to/zOqPAQAAAADYIhkpAAAAwB4zy4a8JLl3kmd1x3p+N8m5rbVb1rn+vCT3SPK0rARLr19jziuTPH+KOgEAAAAAZkFGCgAAALCHzLIh7/uTnJHkYUnum5Vffd4lyWeT7E/yZ0le3Vr7840Waa2Nkjy9+xXpM5Kc3q31ySSXJnlFa+0ds/ojAAAAYC9pbV9Gbd+8y5ioDaBGgMhIAQAAYHBkpExrZg15rbWLk1zc43pvT/L2vtYDAAAAAJglGSkAAADA3qNVEgAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6MLMtawEAAIDhWU5lOTXvMiYaQo0AAAAAwPDISJmWN+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA9sWQsAAAAcMmrJqC3+VgejNu8KAAAAAIDdSEbKtLwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAerA07wIAAACAxTFq+zJqi//7vSHUCAAAAAAMj4yUafmXAQAAAAAAAAAAgB5oyAMAAAAAAAAAAIAe2LIWAAAAOKSlMkrNu4yJ2gBqBAAAAACGR0bKtLwhDwAAAGANVXWPqjqrql5YVe+oqk9WVeuOC3p8znlj6046vmET6x1bVT9bVX9dVddX1We68c9W1bF91Q0AAAAA7F7y0e3zhjwAAACAtV0z7wK2qqpOT3JhknutuvTV3XFuVT2utXbZjhcHAAAAAAyJfHSbNOQBAAAATLY/yRVJHj3j5zxgwvWPrnehqo5P8ntJ7pnktiT/I8lF3eWzkvx4knsnuaiqvra19rHpywUAAAAA9gD56BZoyAMAAAAOWW6V5VbzLmOiHarxhUkuTXJpa+2aqjopGwQ+fWitfXCK238+K2FTkjyltfbGsWvvqarLkryhm/OiJE+b4lkAAAAAsCvJSA+Rj27TvlktDAAAADBkrbUXtNYuaq0t/NYMVXXPJN/bnb5zVdiUJOm+e2d3+tTuHgAAAACALyAf3T4NeQAAAADD9x1JjujGr9pg3gXd5xHdPQAAAAAAQ7dQ+aiGPAAAAIDh+/qx8cUbzBu/9ogZ1QIAAAAAsJMWKh9dmtXCAAAAwPCM2r6M2uL/fm8INW5HVb0rydckuXOS65JcnuT3k7yitfavG9x6avd5fWvt6vUmtdY+XlWfSXKXsXsAAAAAgI6MdH52Sz6qIQ8AAAAYuuOqasMJrbUDO1TLtL55bHz3JGd0x3Oq6pzW2oXr3Hdi97mZv3N/kq8cuwcAAAAAGLbdkpHuinxUQx4AAAAwdJduYs7GadT8/W2Styb5iyRXJbldkq9I8j1JHp3kmCS/W1Xf3lp7xxr337n7/OwmnnVj93mnKeoFAAAAABbH0DPSXZWPasgDAAAADhmlMmqLnMusGC10drRlL22tnbfG9+9P8ltV9UNJfj3JEUl+s6pObq3dvGru7bvPWzfxvFu6z6O3UywAAAAA7GYy0h236/JRDXkAAADA0J2e5Op5F7FdrbXrJlx/RVWdluTcJPdO8oQkr1017XNJ7pDkyE088qjuc3VoBQAAAAAM02Az0t2Yj2rIAwAAAIbu6tbagXkXMWOvyErglCRn5AsDpxuyEjhtZpuFO3afm9m+AQAAAABYfLs9Ix1UPrpvVgsDAAAA0JvLx8bHr3H9YNh2wibWOrH73D9VRQAAAAAAO2NQ+ag35AEAAACHtFRGqXmXMVEbQI09m/QHX57ka5PctaqOa62tuT1FVd0ryV260yt6rA8AAAAAdgUZ6UIaVD7qDXkAAAAAi+/+Y+Or1rj+3rHxGRusM37tkqkqAgAAAADYGYPKRzXkAQAAACy+HxobX7zG9bclGXXjH9hgnXO6z1F3DwAAAADAohtUPqohDwAAADhk1JJRqwEc8/5fanOq6pyqat1x3hrXH1BVJ09Y44eSPL07vTrJW1bP6bZgeG13+piq+s411nliksd0p69Zb9sGAAAAANjLZKT92av56NKsFgYAAAAYsqp6RJLxMOhuY+OTq+qc8fmttQu28ZivTfKbVfUnSd6R5G+TfCormc0pSb43ybd0c5eT/FBr7cZ11npekjOT3D3J66rqtCQXddfOSvIT3fjaJM/fRq0AAAAAwB4hH90+DXkAAAAAazs3yfevc+3h3THugm0+54gk39wd6/lUkqe31tbdRqG1tr+qvj3JW5Mcl+Q53THu6iRnt9YObLNWAAAAAGBvkI9uk4Y8AAAAgPl5e1a2W3hYkgcnuWeSY5NUkk8n+eskv5/kgtbaZyYt1lp7f1U9IMmPJjk7yUndpY8muTDJS1trn+r3TwAAAAAA2JZdmY9qyAMAAAAOGbV9GbV98y5jop2osbV2TpJzplzjgmzwy9DW2ieS/J/u6EVr7ZNJfqY7AAAAAIAtkJGukI9u3+L/1wMAAAAAAAAAAAADoCEPAAAAAAAAAAAAemDLWgAAAOCQUauMWs27jImGUCMAAAAAMDwyUqblDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANCDpXkXAAAAACyOUSqj1LzLmGgINQIAAAAAwyMjZVrekAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPbBlLQAAAHBIa5VRW/ytDtoAagQAAAAAhkdGyrS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgMUxapVRq3mXMdEQagQAAAAAhkdGyrS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAemDLWgAAAOAQ2zEAAAAAAHuZjJRpeUMeAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPRAQx4AAAAAAAAAAAD0YGneBQAAAACLY9Qqo1bzLmOiIdQIAAAAAAyPjJRpeUMeAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPTAlrUAAADAIS3JKIu/1UGbdwEAAAAAwK4kI2Va3pAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9WJp3AQAAAMDiGLXKqNW8y5hoCDUCAAAAAMMjI2Va3pAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9WJp3AQAAAMDiGKUyajXvMiYaZfFrBAAAAACGR0bKtLwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6YMtaAAAA4JBRG8h2DAOoEQAAAAAYHhkp0/KGPAAAAAAAAAAAAOiBhjwAAAAAAAAAAADogYY8AAAAAAAAAAAA6MHSvAsAAAAAFseoVUat5l3GREOoEQAAAAAYHhkp0/KGPAAAAAAAAAAAAOiBhjwAAAAAAAAAAADogS1rAQAAgMNapQ1hq4Mh1AgAAAAADI+MlCl5Qx4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPRgad4FAAAAAItjlMooNe8yJhpCjQAAAADA8MhImZY35AEAAAAAAAAAAEAPNOQBAAAAAAAAAABAD2xZCwAAABwyapVRW/ytDoZQIwAAAAAwPDJSpuUNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0IOleRcAAAAALI7WKq3VvMuYaAg1AgAAAADDIyNlWt6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9sGUtAAAAcMioJaMBbHUwavOuAAAAAADYjWSkTMsb8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgBxryAAAAAAAAAAAAoAdL8y4AAAAAWBytVVqreZcx0RBqBAAAAACGR0bKtLwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6YMtaAAAA4JDWKqMBbHVgOwYAAAAAYBZkpEzLG/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHS/MuAAAAAFgcLUlr865isgGUCAAAAAAMkIyUaXlDHgAAAAAAAAAAAPRAQx4AAAAAAAAAAAD0wJa1AAAAwCGjVEapeZcx0RBqBAAAAACGR0bKtLwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAerA07wIAAACAxdFapbWadxkTDaFGAAAAAGB4ZKRMyxvyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB0vzLgAAAABYHKNWGbWadxkTDaFGAAAAAGB4ZKRMyxvyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHtqwFAAAADmlt5Vh0Q6gRAAAAABgeGSnT8oY8AAAAAAAAAAAA6IGGPAAAAAAAAAAAAOiBhjwAAAAAAAAAAADowVwa8qrqxVXVxo5v2MQ9Z1bVm6vqQFXd0n2+uarOnH3FAAAAsEe0ShvAkVbz/l8KYCoyUgAAAFhQC5B/ykiHbccb8qrqgUmevYX5VVWvSPKOJI9PcnySI7vPxyd5R1W9oqr8VwYAAAAALDwZKQAAAMDutaMNeVW1L8n/TrKU5BObvO3nkjyjG38gyXcneUj3+YHu+2ckeVF/lQIAAAAA9E9GCgAAALC7Le3w8/5TktOTfCjJW5L8540mV9XJSX66O70sySNbazd355dW1duSXJzktCTPqapXtdY+MpPKAQAAYA84tN3BghtCjQDrkJECAADAApORMq0de0NeVZ2Yw7/Q/OEkt27itmfncNPgs8aCpiRJa+2mJM/qTpeS/Nj0lQIAAAAA9E9GCgAAALD77eSWtb+W5E5JXt1ae/ekyVVVSR7XnX6otfa+teZ133+4Oz27uw8AAAAAYNHISAEAAAB2uR1pyKuqJyU5K8mnk/zUJm+7T5Lju/HFE+YevH5CkpO2Wh8AAAAAwCzJSAEAAAD2hqXJU6ZTVcck+ZXu9DmttWs3eeupY+MPTZg7fv3UJB/d5DMAAACAMaNWGbXFf7HSTtRYVfdI8pDuOL07ju0uv7q1dk5Pz7lLkm9N8qgkX5vkvknukOT6JH+X5KIkv9lau27COlcm+dJNPPKfW2snbb9iYKtkpAAAADAcMtIV8tHtm3lDXpIXJzkuyZ8leeUW7jtxbHxgwtz969y3KVV1woQpx211TQAAAGDwrpn1A6rqsUnekuSoNS7fLckZ3fGTVfXdrbU/mXVNwEwsdEYqHwUAAADWIB/dppk25FXVI5Kcm+S2JM9srbUt3H7nsfFnJ8y9cWx8py0846D9k6cAAAAAe9j+JFckeXTP6x6blbBplORdSX4/yV8nuS4r205+T5InJ7lnkouq6uGttb+asOaFSZ6/wfVbpysZ2IqBZKTyUQAAAGAj8tEtmFlDXlUdmeQ3klSS/9la+9stLnH7sfGk/yFuGRsfvcXnAAAAAJ3WVo5Ft0M1vjDJpUkuba1dU1Unpf8tIP8tySuS/LfW2r+suvaBJL9XVZck+dWsbNPwy1nZumEj17XWPthzncA2yEgBAABgeGSkh8hHt2mWb8j7L0lOTfIvSX52G/d/bmx85IS5468tvHkbz5q0hcNxWfkPDAAAANgjWmsv2IFn/E6S35kw5/yqemqS05J8Q1Ud21r71KxrA3oxlIxUPgoAAAB8Hvno9s2kIa+qTknyn7vTZ7XWbtxo/jpuGBtP2mLhjmPjSVs3fIHW2oGNrlfVVpcEAAAA6NO7sxI47UtynyQLHTgBw8pI5aMAAADAgnt3BpSPzuoNec/Oyi82/ynJHarqu9aY81Vj42+qquO68e914dR4CHTChOeN/4Jz/1aLBQAAAFhw42++Gs2tCmArZKQAAAAA/RhUPjqrhryD/yPcN8nrNjH/Z8bG90lyY5LLx747ZcL949ev2MTzAAAAgDW0lrS2+G9Cam3eFey4M7rP25L844S5j6yqv0lyvySV5Jokf5GVjObC1vbg/3owHzJSAAAAGCAZ6UIaVD46q4a8Pnw0yVVJ7p3D/6Ou55Hd58eSXDnDmgAAAIDFc9yk7RQnbce4yKrq25J8dXf6ztbaZybccp9V5yd1x5OSXFJVT26tfazXIoFZkZECAAAAm7FrM9Ih5qP7ZrFoa+2c1lptdCT52bFbvnHs2pXdGi3Jhd31U6rqoWs9q/v+4K8//cobAAAA9p5Ls7I940bHIFXVFyd5eXe6nM9/g9ZqtyZ5W5L/mOQbkjw4yTcm+S85/L/Bw5O8q6ruOot6gcNkpAAAAMAO2pUZ6VDz0UV+Q16SvDTJD2alzvOr6pGttZsPXqyqo5Oc353e1s0HAAAAtqmlhrEdQxa/xmlV1RFJXpvkS7uvfq619oENbnlIa+26Nb5/d1W9LMmbkjw6yalJXpDkx3ssF5idl0ZGCgAAADtGRroYhpyPzuQNeX1prf19kpd0p6ele21gVZ1WVU9Ockn3fZL8UmvtH+ZRJwAAADBXpyc5ccIxRL+W5Mxu/H+TvGijyeuETQev3ZCVLRk+1X31jKo6socagRmTkQIAAACbsBsz0sHmo4v+hrwkeV6SeyR5WlZeJfj6Nea8Msnzd7IoAAAAYGFc3Vo7MO8i+lRV/z3JM7rT9yZ5YmtteZo1W2vXV9Xrk/xIkjtmpYHnz6YqFNgpMlIAAABgI7sqIx16PrrQb8hLktbaqLX29CTfluTCJFdlZc/fq7rzb22tndtaG82xTAAAAIBeVNVzkjy3O/1/Sc4a355ySpePjY/vaU1gxmSkAAAAwF6xG/LRub0hr7V2XpLztjD/7UnePqt6AAAAgKR1x6IbQo3bUVX/IckvdKdXJHlMa+36Ph/R41rAlGSkAAAAsHhkpPOzW/LRhX9DHgAAAMBeUFXfl+Rl3ek/Jfnm1tone37M/cfGV/W8NgAAAADAtuymfFRDHgAAAMCcVdUTkrwqK7/QPJDkUa21XgOhqrprkid3pzcluazP9QEAAAAAtmO35aMa8gAAAIBDWqvBHENQVedUVeuO89aZ8+gkr0tyRJJPZOWXn1du8TlnVtXRG1y/c5I3JDm2++qVrbVbtvIMAAAAANgL5p177qaMdK/mo0uzWhgAAABgyKrqEUlOHvvqbmPjk6vqnPH5rbULtvGMhyZ5S5Ijk/xbkmcnuV1VfdUGtx1orV236rvnJnltVb05yXuTfCTJZ5Mck+RhSX44yYnd3A8nOW+rtQIAAAAAe4d8dPs05AEAAACs7dwk37/OtYd3x7gLtvGMM5PcoRvfLslrN3HPD6zzrC/OSs3nbnDvnyZ5Smvt01uoEQAAAADYe+Sj26QhDwAAAGD4fjLJo7Lya8+vyMqvVY9JclOSq5K8PyvbPvxBa63NqUYAAAAAgFlYqHxUQx4AAABwWOuORbcDNbbWzklyzpRrXJANfhnaWjsvPWyP0Fq7LMll064DAAAAAHuejHRlefnotu2bdwEAAAAAAAAAAACwG2jIAwAAAAAAAAAAgB7YshYAAAA4rFVaq3lXMdkQagQAAAAAhkdGypS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgMXR2sqx6IZQIwAAAAAwPDJSpuUNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0IOleRcAAAAALI7WKq3VvMuYaAg1AgAAAADDIyNlWt6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9sGUtAAAAcFirlWPRDaFGAAAAAGB4ZKRMyRvyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB0vzLgAAAABYHK2tHItuCDUCAAAAAMMjI2Va3pAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD2wZS0AAABwWOuORTeEGgEAAACA4ZGRMiVvyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHmjIAwAAAAAAAAAAgB4szbsAAAAAYHG0Vmmt5l3GREOoEQAAAAAYHhkp0/KGPAAAAAAAAAAAAOiBhjwAAAAAAAAAAADogS1rAQAAgM/X5l0AAAAAAMAcyUiZgjfkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPNOQBAAAAAAAAAABAD5bmXQAAAACwOFqrtFbzLmOiIdQIAAAAAAyPjJRpeUMeAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPTAlrUAAADAYa07Ft0QagQAAAAAhkdGypS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgEVS3bHohlAjAAAAADA8MlKm4w15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADW9YCAAAAh7XuWHRDqBEAAAAAGB4ZKVPyhjwAAAAAAAAAAADogYY8AAAAAAAAAAAA6IGGPAAAAAAAAAAAAOjB0rwLAAAAABZI645FN4QaAQAAAIDhkZEyJW/IAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAe2LIWAAAAOKwlaTXvKiazHQMAAAAAMAsyUqbkDXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANCDpXkXAAAAACyO1laORTeEGgEAAACA4ZGRMi1vyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHmjIAwAAAAAAAAAAgB4szbsAAAAAYIG07lh0Q6gRAAAAABgeGSlT8oY8AAAAAAAAAAAA6IGGPAAAAAAAAAAAAOiBLWsBAACAw1qtHItuCDUCAAAAAMMjI2VK3pAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9WJp3AQAAAMDiqCTV5l3FZDXvAgAAAACAXUlGyrS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAemDLWgAAAOCw1h2Lbgg1AgAAAADDIyNlSt6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPViadwEAAADAAmm1ciy6IdQIAAAAAAyPjJQpeUMeAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPTAlrUAAADAYa07Ft0QagQAAAAAhkdGypS8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgAXSumPRDaFGAAAAAGB4ZKRMyRvyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHtqwFAAAADrMdAwAAAACwl8lImZI35AEAAAAAAAAAAEAPNOQBAAAAAAAAAABADzTkAQAAAAAAAAAAQA+W5l0AAAAAsEBarRyLbgg1AgAAAADDIyNlSt6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9sGUtAAAAcEi1lWPRDaFGAAAAAGB4ZKRMyxvyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB0vzLgAAAABYIK07Ft0QagQAAAAAhkdGypS8IQ8AAABgDVV1j6o6q6peWFXvqKpPVlXrjgtm9Mzvqqp3VtXHq+pzVXVlVb2mqh66hTWOraqfraq/rqrrq+oz3fhnq+rYWdQNAAAAAOwu8tHt84Y8AAAAgLVds1MPqqrbJ3ljkrNWXfrS7nhKVZ3XWnvRhHVOT3JhknutuvTV3XFuVT2utXZZP5UDAAAAALuUfHSbvCEPAAAAYLL9Sf5ghuu/MofDpj9JcnaShyR5epKPZCXDeWFVnbveAlV1fJLfy0rYdFuSFyd5ZHe8uPvu3kku6uYCAAAAAGyGfHQLvCEPAAAAYG0vTHJpkktba9dU1UlJPtr3Q6rqjCRP6U5/L8njW2vL3fmlVfW2JH+Z5EuSvLiq3tRau26NpX4+yT278VNaa28cu/aeqrosyRu6OS9K8rR+/xIAAAAAYBeRj26TN+QBAAAArKG19oLW2kWttVlvzfDT3edykv8wFjYdrOOTSZ7TnX5RVn4V+nmq6p5Jvrc7feeqsOngOm9M8s7u9KndPQAAAAAAX0A+un0a8gAAAADmpKrulORR3em7WmsH1pn65iSf6cZPWOP6dyQ5ohu/aoNHXtB9HtHdAwAAAAAwF7s1H9WQBwAAABxSSaoN4Jj3/1D9eUiSo7rxxetNaq3dmuR9B++pqtutmvL1Y+N111l17RGbLRIAAAAA9goZ6Y7alfmohjwAAACA+Tl1bPyhCXMPXl9K8mXrrHN9a+3q9RZorX08h39Jeup68wAAAAAAdsCuzEeXZrUwAAAAwA45rmrj34NusNXBvJ04Np5U4/5V912+xjqb+Tv3J/nKVc8GAAAAAIZrqBnprsxHNeQBAAAAQ3fpJuYs6g4Odx4bf3bC3BvHxndaZ51Ja4yvs3oNAAAAAGCYhpqR7sp8VEPeFn3p//7H3H6fvBoAAGAv+dzos7ly3kXslFYrx6IbQo2bc/ux8a0T5t4yNj56nXUmrTG+zuo1ACb65x88OUt3PWbeZQAAALCDbrv+uuSX5l3FDpKR7qRdmY9qyAMAAACG7vQkV8+7iG363Nj4yAlzjxob37zGOnfYxBrj66xeAwAAAAAYpqFmpLsyH9WQBwAAAAzd1a21A/MuYptuGBtPeiX/HcfGq7deuCErgdNmXut/cJ3NbN8AAAAAACy+oWakuzIf3TerhQEAAIABagM6dofxkOyECXNPHBvvX2edSWuMr7N6DQAAAABg3rnn3spId2U+qiEPAAAAYH4uHxufMmHuweu3JfnHdda5a1Udt94CVXWvJHfpTq/YbJEAAAAAADOwK/NRDXkAAAAA83Npklu78RnrTaqqI5M89OA9rbVbV01579h43XVWXbtks0UCAAAAAMzArsxHNeQBAAAAzElr7YYkf9SdfnNVrbelwhNy+Jebb1nj+tuSjLrxD2zwyHO6z1F3DwAAAADAXOzWfFRDHgAAAHBYG9AxAFV1TlW17jhvnWkv6T6Xkry8qo5Ytcbdkvxid3pdkt9cvUBr7eokr+1OH1NV37lGLU9M8pju9DXdPQAAAADAuHnnnrsoI92r+ejSrBYGAAAAGLKqekSSk8e+utvY+OSqOmd8fmvtgu08p7X2x1X1+iTfleQ7kryrql6a5KokD0jyvCRf0k1/bmvtX9dZ6nlJzkxy9ySvq6rTklzUXTsryU9042uTPH87tQIAAAAAe4N8dPs05AEAAACs7dwk37/OtYd3x7gLpnjW07Ky5cK3JvnG7hg3SvKi1tor1lugtba/qr49yVuTHJfkOd0x7uokZ7fWDkxRKwAAAACw+8lHt0lDHgAAAHBItZVj0Q2hxq1ord2c5Nuq6ilJzknywCTHJLkmyXuSvKy19uebWOf9VfWAJD+a5OwkJ3WXPprkwiQvba19qufyAQAAAGDXkJHuvN2Wj86sIa+qDnYtnp7ktCTHZ+WVgEdnZT/fy5O8PckrN/OHVtWZSZ6R5CHdOtcm+Yskv9Fa+/0Z/AkAAADAHtZaOycr4c80a1yQLfwytLX220l+e8pnfjLJz3QHMEcyUgAAAGCo5KPbN8s35D0kyevWuXb3JGd0x09V1fe21t651sSqqiS/npWgadzxSR6f5PFV9RtJntla20W9nwAAAADAwMlIAQAAAPaYWW9Zuz/JnyT5y2788ST7kpyQ5DuTPCHJ3ZK8rapOb639zRpr/FwOB00fSPLiJB9Jcr8kP53kwd31a5M8f2Z/CQAAAADA1slIAQAAAPaQWTbk/Ulr7Us2uP6Gqjo7yVuSHJnkBUn+/fiEqjo5K4FSklyW5JHdnsFJcmlVvS3JxVnZ7uE5VfWq1tpHevwbAAAAYG9p3bHohlAjgIwUAAAAhkdGypT2zWrh1tryJua8NcmHutNHrjHl2TncNPissaDp4P03JXlWd7qU5Me2UysAAAAAQN9kpAAAAAB7z8wa8rbgxu7z9uNfVlUleVx3+qHW2vvWurn7/sPd6dndfQAAAAAAQyEjBQAAANgl5tqQV1WnJnlQd/qhVZfvk+T4bnzxhKUOXj8hyUl91AYAAAB7UhvQAbALyEgBAABgwcw795SRDt6ON+RV1R2q6suq6seT/EmSI7pLv7Jq6qlj49VB1Grj109ddxYAAAAAwJzJSAEAAAB2r6WdeEhVnZPkVRtMeUmS16767sSx8YEJj9i/zn2bUlUnTJhy3FbXBAAAAAA4aJEzUvkoAAAAQH92pCFvA3+V5Jmttfevce3OY+PPTljnxrHxnbZRx/7JUwAAAAAAevdXmX9GKh8FAAAA6MlONeS9Ncll3fjoJPdL8qQkj0/y2qr6sdbaRavuuf3Y+NYJ698yNj56ijoBAABgT6u2ciy6IdQIsMpbIyMFAACAhScjZVo70pDXWrsuyXVjX12a5PVV9X1JXp3kwqp6emvtgrE5nxsbHznhEUeNjW/eRomTtnA4Lis1AwAAAABs2YJnpPJRAAAAgJ7Mdcva1tprquqsrPwS9GVVdWFr7V+7yzeMTZ20xcIdx8aTtm5Yq44DG12vqq0uCQAAAAAw0SJkpPJRAAAAgP7sm3cBSS7sPu+Y5LFj34+HQCdMWGP8F5z7+ygKAAAA9qRWwzkAdg8ZKQAAACyKeeeeMtLBW4SGvGvHxl86Nr58bHzKhDXGr18xdUUAAAAAADtHRgoAAACwSyxCQ97xY+PxrRQ+muSqbnzGhDUe2X1+LMmV/ZQFAAAAALAjZKQAAAAAu8QiNOQ9cWz8twcHrbWWw1s1nFJVD13r5u77g7/+vLC7DwAAAABgKGSkAAAAALvEzBryquqcqrr9hDnPTvKt3emVSd67aspLk9zWjc+vqqNX3X90kvO709u6+QAAAMA02gAOgAGQkQIAAMBAzTv/lJEO2tIM1z4vyS9X1e9mJUT6SFa2W7hzkgck+Z4kD+/m3prkB1trt40v0Fr7+6p6SZLnJjktySVV9YvdWvdL8pwkD+6m/1Jr7R9m+PcAAAAAAGzFeZGRAgAAAOwps2zIS5IvTvKD3bGeA0me1lr7w3WuPy/JPZI8LSvB0uvXmPPKJM+fok4AAAAAgFmQkQIAAADsIbNsyHtUkm9O8o1JTk1yzyTHJvlckmuS/FWSi5K8obV203qLtNZGSZ7e/Yr0GUlOT3K3JJ9McmmSV7TW3jG7PwMAAAD2jmorx6IbQo0AkZECAADA4MhImdbMGvJaax/JyrYJr+hpvbcneXsfawEAAAAAzJqMFAAAAGDv2TfvAgAAAAAAAAAAAGA30JAHAAAAAAAAAAAAPZjZlrUAAADAALXuWHRDqBEAAAAAGB4ZKVPyhjwAAAAAAAAAAADogYY8AAAAAAAAAAAA6IEtawEAAIDDWlJD2OpgCDUCAAAAAMMjI2VK3pAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9WJp3AQAAAMACad2x6IZQIwAAwP/f3r1HW3fW9aH//sJLyIWLlCipCQhCD4lWPQhhkAoEq0VKsEYOyqWWBohU2sER0RaqdghYUSmchkJbaeUQUaoFhYSLFOgBUkCUxEOtlmAMF08CBUUI5EISwvs7f6z5Zq/s7Nvaa+53r/nuz2eMOdYz13zmM38xj8mbL8+aDwAwPTJSluQNeQAAAAAAAAAAADACC/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAgO7XcBAAAAwArp4Vh1U6gRAAAAAJgeGSlL8oY8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADACW9YCAAAAt6meHatuCjUCAAAAANMjI2VZ3pAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADACC/IAAAAAAAAAAABgBIf2uwAAAABghfRwrLop1AgAAAAATI+MlCV5Qx4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADCCQ/tdAAAAALA6qmfHqptCjQAAAADA9MhIWZY35AEAAAAAAAAAAMAILMgDAAAAAAAAAACAEdiyFgAAALg9Wx0AAAAAAAeZjJQleEMeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwQno4Vt0UagQAAAAApkdGypK8IQ8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjMCWtQAAAMBtqmfHqptCjQAAAADA9MhIWZY35AEAAAAAAAAAAMAILMgDAAAAAAAAAACAEViQBwAAAAAAAAAAACM4tN8FAAAAACukh2PVTaFGAAAAAGB6ZKQsyRvyAAAAAAAAAAAAYAQW5AEAAAAAAAAAAMAIbFkLAAAA3KZ6dqy6KdQIAAAAAEyPjJRleUMeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwQno4Vt0UagQAAAAApkdGypK8IQ8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjMCWtQAAAMAa2zEAAAAAAAeZjJQleUMeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwOqpnx6qbQo0AAAAAwPTISFmWN+QBAAAAAAAAAADACCzIAwAAAAAAAAAAgBFYkAcAAAAAAAAAAAAjOLTfBQAAAAArpIdj1U2hRgAAAABgemSkLMkb8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACGxZCwAAANyerQ4AAAAAgINMRsoSvCEPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzAgjwAAAAAAAAAAAAYwaH9LgAAAABYHdWzY9VNoUYAAAAAYHpkpCzLG/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAhsWQsAAACs6eFYdVOoEQAAAACYHhkpS/KGPAAAAIAtVNV9q+plVXVFVd1QVV+oqg9X1U9V1UlLjPvoquoFj/dtMtandnj/p3ZbLwAAAABwMMlIF+MNeQAAAACbqKpzk7w+yT3mvj4pyVnDcUFVPa67P3GUSvrTo/QcAAAAAAAZ6S5YkAcAAACwgar6jiRvyCxcuj7JLyZ5b5ITkzw5yY8meVCSt1fVWd19/YKPuCzJt+2g36uSnDO0f22bvpck+dktrt+yg+cBAAAAAMhId8mCPAAAAOA21bNj1R2lGi/MLGi6NcljuvtDc9feU1V/luSlSc5I8rwkL15k8O6+IcmfbNWnqr4uycOH06u6+/e2Gfba7t5yTAAAAABgczLS27kwMtKFHbefDwcAAABYRVV1VpJHD6evWRc0HfHyJFcM7edW1Z33oJQnJbnL0P71PRgfAAAAAOAOZKS7Z0EeAAAAwB2dN9d+7UYduvtwktcNp/fMWjg1pqcdeVwmEjYBAAAAAMeE8+baMtIFWJAHAAAArOkJHXvrkcPnDUn+cIt+l861HzFmAVX1gCR/azh9f3d/cszxAQAAAIAN7HfuKSO9zVQzUgvyAAAAAO7ozOHzqu6+dYt+H9vgnrE8ba79azu851FV9T+q6oaqurGqPllV/7mqzquqGrk+AAAAAODYJSPdpUNH4yEAAAAAe+jU7XKU7r5mp4NV1QlJThlOt7yvu79YVTckOTnJfXb6jB36keHzK0l+e4f33H/d+f2G44eTfLCqntTdnx6lOgAAAABgVchI1+x7RmpBHgAAADB1l+2gzyK/fLzbXPv6HfQ/EjbddYFnbKmqHpnkm4fTN3f3l7e55ZYkb0nyriR/kuRLSb4uydlJnp1ZEPZdSd5dVWd395fGqhUAAAAA2Hcy0hXKSC3IAwAAANb0cKy6va3xhLn2LTvof/PweeKINfyDufbrdtD/Yd197Qbfv6+qXpXZr0cfk9mWET+X5HlLVwgAAAAAxyIZaSIjXYoFeQAAAMDUnZXksyOOd9Nc+/gd9L/L8PmVMR5eVXdJ8kPD6WeS/Nft7tkkaDpy7bqq+uEkH09yryTPqqoXdPdOgjQAAAAAYPXJSFcoI7UgDwAAAJi6z3b3NSOOd91ceydbLJw8fO5k64ad+IHMtlJIktd399eWHbC7v1RVv5Xkn2RW70OT/N6y4wIAAAAAK0FGuo2jmZFakAcAAADcpoZj1e1ljd19U1V9PskpSU7fso6qe2YtbLp6pBKeNtfeyVYMO/XRufZpI44LAAAAAMcMGamMdFnH7dXAAAAAABN2xfD5wKra6geNZ2xwz65V1Tck+b7h9P/t7j9Zdsz54UccCwAAAAA4tslId8mCPAAAAIA7+sDweXKSh2zR75y59gdHeO5Ts7ajwZi//EySb5lrf2bksQEAAACAY4uMdJcsyAMAAAC4o4vn2k/fqENVHZe1rROuTfLeEZ57ZLxbk/ynEcZLklTVPZI8aTi9McnlY40NAAAAAByTLp5ry0gXYEEeAAAAsKYndOyh7v5wkvcPp8+sqrM36PaTSc4c2q/o7q/OX6yq86uqh+OF2z2zqr41yYOH03d091/upNaqemxVnbjF9bsleUOSew1fvaa7b97J2AAAAABw4Ox37ikjnXxGutX+vgAAAAAH2Y9ntsXCiUneVVUvyewXnicmeXKSZw39rkzy8hGe9w/n2r+2wH0vSPL6qnpTZttIfDzJ9Um+LsnZSZ6d5D5D3z9N8sJlCwUAAAAADgQZ6S5YkAcAAACwge7+SFU9KclvJLl7kpds0O3KJOd293XLPGvY2uGpw+kXk7xtwSH+WpILhmMz/y3JU7v7C4tXCAAAAAAcNDLS3bEgDwAAAFjTSe3xVgejOEo1dvdbq+rbM/sl6LlJTk9yS5Krkrwxyau6+8YRHvU9SU4b2v95we0Sfmq4/+wkD0pySma//LwxyWeS/EGS30zyru6ewt9dAAAAANg/MtLbP0ZGujAL8gAAAAC20N1/nuR5w7HIfRcluWiHfd+dpBatbbj38iSX7+ZeAAAAAIDtyEgXc9x+FwAAAAAAAAAAAADHAgvyAAAAAAAAAAAAYAS2rAUAAADW9HCsuinUCAAAAABMj4yUJXlDHgAAAAAAAAAAAIzAgjwAAAAAAAAAAAAYgQV5AAAAAAAAAAAAMIJD+10AAAAAsGJ6vwsAAAAAANhHMlKW4A15AAAAAAAAAAAAMAIL8gAAAAAAAAAAAGAEtqwFAAAAblM9O1bdFGoEAAAAAKZHRsqyvCEPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzAgjwAAAAAAAAAAAAYwaH9LgAAAABYIT0cq24KNQIAAAAA0yMjZUnekAcAAAAAAAAAAAAjsCAPAAAAAAAAAAAARmDLWgAAAOA21bNj1U2hRgAAAABgemSkLMsb8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACCzIAwAAAAAAAAAAgBEc2u8CAAAAgBXSw7HqplAjAAAAADA9MlKW5A15AAAAAAAAAAAAMAIL8gAAAAAAAAAAAGAEtqwFAAAAblM9O1bdFGoEAAAAAKZHRsqyvCEPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzAgjwAAAAAAAAAAAAYwZ4uyKuq76yqn66qd1TV1VV1c1VdX1VXVtVFVfXIBcd7bFW9qaquGca6Zjh/7F79NQAAAMCB0hM6ACZARgoAAAATs9+5p4x08g7t1cBVdWmSR21w6fgkf2M4/mFV/XqSC7r7li3GqiS/kuRZ6y6dluQHk/xgVf2HJD/W3aYbAAAAALDvZKQAAAAAB89eviHvtOHzM0lekeSJSR6W5Owkz0vy6eH6P0hy0TZj/cusBU0fSfKUYaynDOcZrv/8CHUDAAAAAIxBRgoAAABwwOzZG/KSfCzJTyf5ne7+2rprvz/86vODSf63JE+pqn/f3e9fP0hVPTDJPxtOL0/yqO7+ynB+WVW9JcmlSR6a5PlV9dru/vge/PUAAADAsW8qWx1MoUYAGSkAAABMj4yUJe3ZG/K6+/Hd/YYNgqYj1z+f5CfnvnriJkP9RNYWDj5nLmg6Ms6NSZ4znB5K8txdFw0AAAAAMBIZKQAAAMDBs5db1u7E++baD1h/saoqyQ8Mpx/r7t/faJDh+z8dTs8b7gMAAAAAWHXvm2vLSAEAAAAmbr8X5B0/1z68wfX7JzltaF+6zVhHrp+e5H7LlQUAAAAAcFTISAEAAACOIYe277Knzplrf2yD62ducz2bXD8zySd3WxQAAAAcVNWzY9VNoUaAHZKRAgAAwAqRkbKsfVuQV1XHJXnB3Fdv2KDbfeba12wz5NWb3LeTWk7fpsupi4wHAAAAALCdVclI5aMAAAAA49nPN+T9RJKHDe03d/flG/S521z7+m3Gu2GufdcFa7l6+y4AAAAAAKNalYxUPgoAAAAwkn1ZkFdV5yT5peH0L5I8e5OuJ8y1b9lm2Jvn2ifusjQAAADAVgcAe05GCgAAACtMRsoSjvqCvKr61iRvHp59c5If7u7PbdL9prn28dsMfZe59lcWLGu77RtOTXLZgmMCAAAAANzBCmak8lEAAACAkRzVBXlVdf8k70pyzyRfS/KU7r50i1uum2tvt8XCyXPt7bZuuJ3uvmar61W1yHAAAAAAABtaxYxUPgoAAAAwnuOO1oOq6huT/Nck35jZix2f0d1v3ua2+SDo9G36zv+K8+rFKwQAAAAA2DsyUgAAAIBj31F5Q15VnZLk3Um+efjqOd39uh3c+tG59hnb9J2/fsUC5QEAAACD6k5173cZ25pCjQDzZKQAAAAwDTJSlrXnb8irqnskeWeSbxm+ekF3/9sd3v7JJJ8Z2uds0/dRw+enk3xqkRoBAAAAAPaKjBQAAADg4NjTBXlVdVKStyf5zuGrX+juX97p/d3dSS4ZTs+oqodv8pyHZ+3Xn5cM9wEAAAAA7CsZKQAAAMDBsmcL8qrq+CRvTvJdw1ev6O6f3cVQFya5dWi/sqpOXPecE5O8cji9degPAAAA7EZP6ABYcTJSAAAAmKD9zj1lpJN3aA/H/s0kjxna70nymqr6m1v0v6W7r1z/ZXdfWVUvS/KCJA9N8sGq+uUkH0/ygCTPT/Lgofu/6u4/G+svAAAAAABgCTJSAAAAgANmLxfkPWGu/beT/I9t+v95kvttcu1nknxDkmdkFiz91gZ9XpNkN78uBQAAAADYCzJSAAAAgANmz7asHVN3H+7uZyY5N8klST6T5Jbh85Ikj+vuC7r78D6WCQAAAACwJ2SkAAAAANOwZ2/I6+7agzF/N8nvjj0uAAAAMFM9O1bdFGoEkJECAADA9MhIWdYk3pAHAAAAAAAAAAAAq86CPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwQno4Vt0UagQAAAAApkdGypK8IQ8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjMCWtQAAAMBtqmfHqptCjQAAAADA9MhIWZY35AEAAAAAAAAAAMAILMgDAAAAAAAAAACAEViQBwAAAAAAAAAAACM4tN8FAAAAACukh2PVTaFGAAAAAGB6ZKQsyRvyAAAAAAAAAAAAYAQW5AEAAAAAAAAAAMAIbFkLAAAA3KZ6dqy6KdQIAAAAAEyPjJRleUMeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwQno4Vt0UagQAAAAApkdGypK8IQ8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjMCWtQAAAMDtlK0OAAAAAIADTEbKMrwhDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGMGh/S4AAAAAWCHds2PVTaFGAAAAAGB6ZKQsyRvyAAAAAAAAAAAAYAQW5AEAAAAAAAAAAMAIbFkLAAAA3KZ6dqy6KdQIAAAAAEyPjJRleUMeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwQno4Vt0UagQAAAAApkdGypK8IQ8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjMCWtQAAAMBt6vDsWHVTqBEAAAAAmB4ZKcvyhjwAAAAAAAAAAAAYgQV5AAAAAAAAAAAAMAIL8gAAAAAAAAAAAGAEh/a7AAAAAGCF9HCsuinUCAAAAABMj4yUJXlDHgAAAAAAAAAAAIzAgjwAAAAAAAAAAAAYgS1rAQAAgNtUz45VN4UaAQAAAIDpkZGyLG/IAwAAANhCVd23ql5WVVdU1Q1V9YWq+nBV/VRVnbTk2C+sqt7h8egdjHevqnpRVf1RVX2pqr48tF9UVfdaplYAAAAA4GCSkS7GG/IAAAAANlFV5yZ5fZJ7zH19UpKzhuOCqnpcd39iP+qbV1VnJbkkyV9fd+nbh+OCqvqB7r78qBcHAAAAAEySjHRxFuQBAAAAbKCqviPJGzILl65P8otJ3pvkxCRPTvKjSR6U5O1VdVZ3X7/kI79tm+uf3KLW05K8Ncm9k9ya5P9K8rbh8uOTPC/JNyZ5W1U9pLs/vWStAAAAAMAxTka6OxbkAQAAAGu6Z8eqOzo1XphZ0HRrksd094fmrr2nqv4syUuTnJFZmPPiZR7W3X+yxO2/kFnQlCRP7e43zl17f1Vdnllwdu8kP5/kGUs8CwAAAACOXTLSeRdGRrqw4/ZqYAAAAICpGrY2ePRw+pp1QdMRL09yxdB+blXd+WjUtl5V3TvJjwyn71wXNCVJhu/eOZw+bbgHAAAAAGBDMtLdsyAPAAAA4I7Om2u/dqMO3X04yeuG03tmLZw62v5ekjsN7Q1rHVw0fN5puAcAAAAAYDPnzbVlpAuwIA8AAADgjh45fN6Q5A+36HfpXPsRe1fOlh451750016rUSsAAAAAMA0y0l06tFcDAwAAANNTPTtW3VGo8czh86ruvnWLfh/b4J5dqap3J/nOJHdLcm2Sjyb5L0le3d1f3OLWI8/9Und/drNO3f2/qurLSe6+bK0AAAAAcKySkd5GRrpLFuQBAAAAU3dqVW3Zobuv2elgVXVCklOG0y3v6+4vVtUNSU5Ocp+dPmMT3zvX/vok5wzH86vq/O6+ZJP7jjx3J3+NVyf51ixfKwAAAACwOmSkO6h1sOcZqQV5AAAAwNRdtoM+W6dRt3e3ufb1O+h/JGy66wLPmPfHSS5O8uEkn0ly5yQPSvL3kzwmydcl+Z2q+v7ufscW9e601ixRKwAAAACwemSki9WaJWrdlgV5AAAAwO1NYDuGPXbCXPuWHfS/efg8cRfPurC7X7jB93+Q5HVV9Y+S/EqSOyX51ap6YHd/ZV3fI/Xuda0AAAAAcDDISGWkS7AgDwAAAJi6s5J8dsTxbpprH7+D/ncZPteHQNvq7mu3uf7qqnpokguSfGOSJyR5/bpuNyU5KXtcKwAAAACwsmSkK5SRWpAHAAAATN1nu/uaEce7bq69k20LTh4+d7Idwm68OrOwKUnOyR3DpusyC5tWoVYAAAAA4OiTka5QRnrcXg0MAAAAMEXdfVOSzw+np2/Vt6rumbUA5+o9Kumjc+3TNrh+JGjbstbBfYbPvaoVAAAAAJg4GelyLMgDAAAAblM9nWOPXTF8PrCqttph4IwN7hlbbXP9SBh1j6o6ddNBqv56krsPp3tVKwAAAABM2n7nnjLSDU0qI7UgDwAAAOCOPjB8npzkIVv0O2eu/cE9quVb5tqf2eD6B+ba52xwfaNre1UrAAAAAHBskJHukgV5AAAAAHd08Vz76Rt1qKrjkjxtOL02yXv3qJZ/NNe+dIPrb0lyeGhvWOvg/OHz8HAPAAAAAMBmLp5ry0gXYEEeAAAAsKZ7Osee/p+hP5zk/cPpM6vq7A26/WSSM4f2K7r7q/MXq+r8qurheOH6m6vq26rqgVvVUVX/KMkzh9PPJnnzBrV+Nsnrh9Pvq6onbjDODyX5vuH014d7AAAAAID19jv3lJGu7zO5jHSr/X0BAAAADrIfz2zbghOTvKuqXpLZLzxPTPLkJM8a+l2Z5OW7GP8hSX61qt6b5B1J/jjJX2WW15yR5EeS/J2h79eS/KPuvmGTsX4myWOTfH2S36yqhyZ523Dt8ZkFY0nyl0l+dhe1AgAAAAAHj4x0FyzIAwAAANhAd3+kqp6U5DeS3D3JSzbodmWSc7v7ul0+5k5Jvnc4NvNXSZ7Z3ZtuodDdV1fV92e2jcSpSZ4/HPM+m+S87r5ml7UCAAAAAAeIjHR3LMgDAAAA2ER3v7Wqvj2zX4Kem+T0JLckuSrJG5O8qrtv3OXwv5vZVgtnJ3lwknsnuVeSSvKFJH+U5L8kuai7v7yDWv+gqr5tqPW8JPcbLn0yySVJLuzuv9plrQAAAADAASQjXZwFeQAAAMBtqmfHqjuaNXb3nyd53nAsct9FSS7a4vpfJPm/h2MU3f35JP9iOAAAAACABclI70hGupjj9uvBAAAAAAAAAAAAcCyxIA8AAAAAAAAAAABGYMtaAAAAYE0Px6qbQo0AAAAAwPTISFmSN+QBAAAAAAAAAADACCzIAwAAAAAAAAAAgBFYkAcAAAAAAAAAAAAjOLTfBQAAAACro3p2rLop1AgAAAAATI+MlGV5Qx4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEtawEAAIA1h3t2rLop1AgAAAAATI+MlCV5Qx4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADCCQ/tdAAAAALBCejhW3RRqBAAAAACmR0bKkrwhDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwJa1AAAAwG0qSU1gq4Pa7wIAAAAAgGOSjJRleUMeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwgkP7XQAAAACwQjpJ935Xsb0JlAgAAAAATJCMlCV5Qx4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEtawEAAIA1ndQUtjqYQo0AAAAAwPTISFmSN+QBAAAAAAAAAADACCzIAwAAAAAAAAAAgBFYkAcAAAAAAAAAAAAjOLTfBQAAAAArpIdj1U2hRgAAAABgemSkLMkb8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACCzIAwAAAAAAAAAAgBEc2u8CAAAAgNVR3anu/S5jW1OoEQAAAACYHhkpy/KGPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwAlvWAgAAAGsOD8eqm0KNAAAAAMD0yEhZkjfkAQAAAAAAAAAAwAgsyAMAAAAAAAAAAIARWJAHAAAAAAAAAAAAIzi03wUAAAAAq6O6U937Xca2plAjAAAAADA9MlKW5Q15AAAAAAAAAAAAMAIL8gAAAAAAAAAAAGAEtqwFAAAA1vRwrLop1AgAAAAATI+MlCV5Qx4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADCCQ/tdAAAAALBKOune7yJ2YAo1AgAAAADTIyNlOd6QBwAAAAAAAAAAACOwIA8AAAAAAAAAAABGYMtaAAAA4DbVs2PVTaFGAAAAAGB6ZKQsyxvyAAAAAAAAAAAAYAQW5AEAAAAAAAAAAMAILMgDAAAAAAAAAACAERza7wIAAACAFdI9O1bdFGoEAAAAAKZHRsqSvCEPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzAlrUAAADAberw7Fh1U6gRAAAAAJgeGSnL2tM35FXVN1TV46vqxVX1jqr6fFX1cFy0i/EeW1Vvqqprqurm4fNNVfXYPSgfAAAAAGApMlIAAACAg2Wv35D3uTEGqapK8itJnrXu0mlJfjDJD1bVf0jyY93dYzwTAAAAAGAEMlIAAACAA2RP35C3ztVJ3rXLe/9l1oKmjyR5SpKHDZ8fGb5/VpKfX6ZAAAAAAIA9JCMFAAAAOMbt9RvyXpzksiSXdffnqup+ST65yABV9cAk/2w4vTzJo7r7K8P5ZVX1liSXJnlokudX1Wu7++OjVA8AAAAHTffsWHVTqBFgRkYKAAAAUyIjZUl7+oa87v657n5bdy+zLcNPZG3h4HPmgqYjz7gxyXOG00NJnrvEswAAAAAARiMjBQAAADhYjuaWtQurqkryA8Ppx7r79zfqN3z/p8PpecN9AAAAAACTJiMFAAAAmJaVXpCX5P5JThval27T98j105Pcb68KAgAAgGNaT+gAOBhkpAAAAHA07XfuKSOdvFVfkHfmXPtj2/Sdv37mpr0AAAAAAKZDRgoAAAAwIYf2u4Bt3Geufc02fa/e5L5tVdXp23Q5dZHxAAAAAABGsucZqXwUAAAAYDyrviDvbnPt67fpe8Nc+64LPufq7bsAAAAAABx1RyMjlY8CAAAAjGTVF+SdMNe+ZZu+N8+1T9yDWgAAAOCYV92p7v0uY1tTqBFgJDJSAAAAOIpkpCxr1Rfk3TTXPn6bvneZa39lwedst33DqUkuW3BMAAAAAIBlHY2MVD4KAAAAMJJVX5B33Vx7uy0WTp5rb7d1w+109zVbXa+qRYYDAAAAABjLnmek8lEAAACA8az6grz5IOj0bfrO/4rz6j2oBQAAAI593bNj1U2hRoBxyEgBAADgaJKRsqTj9ruAbXx0rn3GNn3nr1+xB7UAAAAAABxtMlIAAACACVn1BXmfTPKZoX3ONn0fNXx+Osmn9qogAAAAAICjSEYKAAAAMCErvSCvuzvJJcPpGVX18I36Dd8f+fXnJcN9AAAAAACTJiMFAAAAmJaVXpA3uDDJrUP7lVV14vzF4fyVw+mtQ38AAABgNzrJ4QkclpkAB8uFkZECAADA0SEjZUmH9nLwqnpEkgfOfXXKXPuBVXX+fP/uvmj9GN19ZVW9LMkLkjw0yQer6peTfDzJA5I8P8mDh+7/qrv/bLS/AAAAAACAJchIAQAAAA6WPV2Ql+SCJP9wk2vfNRzzLtqk788k+YYkz8gsWPqtDfq8JsnPLl4iAAAAAMCekZECAAAAHCBT2LI23X24u5+Z5NwklyT5TJJbhs9Lkjyuuy/o7sP7WCYAAAAAwJ6QkQIAAABMw56+Ia+7z09y/ojj/W6S3x1rPAAAAOD2qjvVvd9lbGsKNQIkMlIAAACYGhkpy5rEG/IAAAAAAAAAAABg1VmQBwAAAAAAAAAAACPY0y1rAQAAgInpJFPY6mACJQIAAAAAEyQjZUnekAcAAAAAAAAAAAAjsCAPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzg0H4XAAAAAKyQ7tmx6qZQIwAAAAAwPTJSluQNeQAAAAAAAAAAADACC/IAAAAAAAAAAABgBLasBQAAANYcHo5VN4UaAQAAAIDpkZGyJG/IAwAAAAAAAAAAgBFYkAcAAAAAAAAAAAAjsCAPAAAAAAAAAAAARnBovwsAAAAAVkd1p7r3u4xtHc0aq+q+Sf7PJOcmuW+Sm5NcleQNSf5dd9+4xNh3T/K4JN+T5CFJvjnJSUm+lOR/Jnlbkl/t7mu3GedTSb5pB4/88+6+327rBQAAAIBjnYx0g2fJSBdiQR4AAADAJqrq3CSvT3KPua9PSnLWcFxQVY/r7k/sYuy/m+TNSe6yweVTkpwzHD9VVU/p7vcu+gwAAAAAgGXISBdnQR4AAADABqrqOzL7hedJSa5P8otJ3pvkxCRPTvKjSR6U5O1VdVZ3X7/gI+6VWdB0OMm7k/yXJH+U5Nokpyf5+0melOTeSd5WVd/V3f99mzEvSfKzW1y/ZcEaAQAAAIADSka6OxbkAQAAAGu6Z8eqOzo1XphZ0HRrksd094fmrr2nqv4syUuTnJHkeUlevOD4X03y6iQv6e7/b921jyR5a1V9MMm/Gep4eWbbNmzl2u7+kwXrAAAAAACOkJHOuzAy0oUdt58PBwAAAFhFVXVWkkcPp69ZFzQd8fIkVwzt51bVnRd5Rnf/5+7+sQ2Cpvk+r0xy+XD66Kq61yLPAAAAAADYDRnp7lmQBwAAAHBH5821X7tRh+4+nOR1w+k9sxZOje19w+dxSe6/R88AAAAAAJh33lxbRroAC/IAAAAA7uiRw+cNSf5wi36XzrUfsUe13GWufXiPngEAAAAAME9GukuH9rsAAAAAYJV00r3fRezAntd45vB5VXffukW/j21wz9jOGT5vTXLVNn0fVVX/I8kDklSSzyX5cJLfTHJJ9yT+5gIAAADAPpKRDmSku2RBHgAAADB1p1bVlh26+5qdDlZVJyQ5ZTjd8r7u/mJV3ZDk5CT32ekzFqjl3CTfPpy+s7u/vM0t67druN9w/HCSD1bVk7r706MWCQAAAADsNxnpmn3PSC3IAwAAAKbush302TqNur27zbWv30H/I2HTXRd4xraq6q8l+bfD6deS/Istut+S5C1J3pXkT5J8KcnXJTk7ybMzC8K+K8m7q+rs7v7SmLUCAAAAAPtKRrpCGakFeQAAAMCansh2DHtb4wlz7Vt20P/m4fPEsQqoqjsleX2Sbxq++pfd/ZEtbnlYd1+7wffvq6pXJfntJI/JbMuIn0vyvLFqBQAAAIBjiow0kZEuxYI8AAAAYOrOSvLZEce7aa59/A7632X4/MqINfy7JI8d2m9P8vNbdd4kaDpy7bqq+uEkH09yryTPqqoXdPdOgjQAAAAAYPXJSFcoI7UgDwAAAJi6z3b3NSOOd91ceydbLJw8fO5k64ZtVdUvJnnWcPqBJD/U3V9bZszu/lJV/VaSf5JZvQ9N8ntLFQoAAAAArAoZ6TaOZkZ63F4MCgAAADBV3X1Tks8Pp6dv1beq7pm1sOnqZZ9dVc9P8oLh9P9N8vjuHutXpR+da5820pgAAAAAwDFGRrocb8gDAAAA1hwejlW39zVekeSRSR5YVYe6+9ZN+p2x7p5dq6p/nOSX5sb6vu7+0jJjrn/EiGMBAAAAwLFJRnqEjHSXvCEPAAAA4I4+MHyenOQhW/Q7Z679wd0+rKr+QZJXDaefSPK93f35LW7ZjW+Za39m5LEBAAAAgGOLjHSXLMgDAAAAuKOL59pP36hDVR2X5GnD6bVJ3rubB1XVE5K8NrNfZ16T5Hu6e9QwqKrukeRJw+mNSS4fc3wAAAAA4Jhz8VxbRroAC/IAAACA21T3ZI691N0fTvL+4fSZVXX2Bt1+MsmZQ/sV3f3V+YtVdX5V9XC8cKPnVNVjkvxmkjsl+YvMfvX5qUVqrarHVtWJW1y/W5I3JLnX8NVruvvmRZ4BAAAAAAfFfueeMtLpZ6SH9mpgAAAAgIn78cy2WDgxybuq6iWZ/cLzxCRPTvKsod+VSV6+6OBV9fAkb05yfJKvJvmJJHeuqr+5xW3XdPe16757QZLXV9WbMttG4uNJrk/ydUnOTvLsJPcZ+v5pkhcuWisAAAAAcCDJSHfBgjwAAACADXT3R6rqSUl+I8ndk7xkg25XJjm3u6/bxSMem+SkoX3nJK/fwT1PT3LRBt//tSQXDMdm/luSp3b3FxaoEQAAAAA4oGSku2NBHgAAAMAmuvutVfXtmf0S9Nwkpye5JclVSd6Y5FXdfeM+lpgkP5XkezL7peeDkpyS2S8/b0zymSR/kNmWD+/q3uN9LAAAAACAY4qMdHEW5AEAAABrumfHqjuKNXb3nyd53nAsct9F2fiXmkeuvzAjbI3Q3ZcnuXzZcQAAAACAyEg3fJSMdBHH7XcBAAAAAAAAAAAAcCywIA8AAAAAAAAAAABGYMtaAAAAYM3hnh2rbgo1AgAAAADTIyNlSd6QBwAAAAAAAAAAACOwIA8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjODQfhcAAAAArJBO0r3fVWxvAiUCAAAAABMkI2VJ3pAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACM4NB+FwAAAACskk6697uIHZhCjQAAAADA9MhIWY435AEAAAAAAAAAAMAILMgDAAAAAAAAAACAEdiyFgAAAFjTE9mOYQo1AgAAAADTIyNlSd6QBwAAAAAAAAAAACOwIA8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjODQfhcAAAAArJDDPTtW3RRqBAAAAACmR0bKkrwhDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwJa1AAAAwJo+PDtW3RRqBAAAAACmR0bKkrwhDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGMGh/S4AAAAAWCHds2PVTaFGAAAAAGB6ZKQsyRvyAAAAAAAAAAAAYAQW5AEAAAAAAAAAAMAIbFkLAAAArDncs2PVTaFGAAAAAGB6ZKQsyRvyAAAAAAAAAAAAYAQW5AEAAAAAAAAAAMAILMgDAAAAAAAAAACAERza7wIAAACAFdI9O1bdFGoEAAAAAKZHRsqSvCEPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzAlrUAAADA7dnqAAAAAAA4yGSkLMEb8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACCzIAwAAAAAAAAAAgBEc2u8CAAAAgBXSPTtW3RRqBAAAAACmR0bKkrwhDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwJa1AAAAwJrDh5M6vN9VbO/wBGoEAAAAAKZHRsqSvCEPAAAAAAAAAAAARmBBHgAAAAAAAAAAAIzAgjwAAAAAAAAAAAAYwaH9LgAAAABYId2zY9VNoUYAAAAAYHpkpCzJG/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAhsWQsAAACssR0DAAAAAHCQyUhZkjfkAQAAAAAAAAAAwAgsyAMAAAAAAAAAAIARWJAHAAAAAAAAAAAAIzi03wUAAAAAK+RwJ9X7XcX2Dk+gRgAAAABgemSkLMkb8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACCzIAwAAAAAAAAAAgBEc2u8CAAAAgNXRfTjdh/e7jG1NoUYAAAAAYHpkpCzLG/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAhsWQsAAACs6SSHe7+r2N4ESgQAAAAAJkhGypK8IQ8AAAAAAAAAAABGYEEeAAAAAAAAAAAAjMCCPAAAAAAAAAAAABjBof0uAAAAAFgh3bNj1U2hRgAAAABgemSkLMkb8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACGxZCwAAAKw5fDjJ4f2uYnuHJ1AjAAAAADA9MlKW5A15AAAAAAAAAAAAMAIL8gAAAAAAAAAAAGAEFuQBAAAAAAAAAADACA7tdwEAAADACumeHatuCjUCAAAAANMjI2VJ3pAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgy1oAAADgNn34cDqH97uMbfXh1a8RAAAAAJgeGSnL8oY8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADACC/IAAAAAAAAAAABgBIf2uwAAAABghXTPjlU3hRoBAAAAgOmRkbIkb8gDAAAAAAAAAACAEViQBwAAAAAAAAAAACOwZS0AAACwpjs5PIGtDmzHAAAAAADsBRkpS/KGPAAAAAAAAAAAABiBBXkAAAAAAAAAAAAwAgvyAAAAAAAAAAAAYASH9rsAAAAAYIV0Jzm831Vsr3u/KwAAAAAAjkUyUpY0uTfkVdV9q+plVXVFVd1QVV+oqg9X1U9V1Un7XR8AAAAAwF6RjwIAAACstkm9Ia+qzk3y+iT3mPv6pCRnDccFVfW47v7EftQHAAAAALBX5KMAAAAAq28yC/Kq6juSvCGzgOn6JL+Y5L1JTkzy5CQ/muRBSd5eVWd19/X7VSsAAABMVR/udK3+VgdtOwbggJGPAgAAwNEhI2VZk1mQl+TCzMKmW5M8prs/NHftPVX1Z0lemuSMJM9L8uKjXiEAAAAAwN64MPJRAAAAgJV33H4XsBNVdVaSRw+nr1kXNh3x8iRXDO3nVtWdj0ZtAAAAAAB7ST4KAAAAMB2TWJCX5Ly59ms36tDdh5O8bji9Z9YCKgAAAACAKTtvri0fBQAAAFhhU9my9pHD5w1J/nCLfpfOtR+R5N17VhEAAAAci/pwksP7XcX2egI1AoxHPgoAAABHi4yUJU3lDXlnDp9XdfetW/T72Ab3bKuqTt/qSHLqbooGAAAApq+q7ltVL6uqK6rqhqr6QlV9uKp+qqpOGvE5T66qd1bV/6qqm6rqU1X161X18AXGuFdVvaiq/qiqvlRVXx7aL6qqe41VK3DUyUcBAACAfSMjXczKvyGvqk5Icspwes1Wfbv7i1V1Q5KTk9xngcdcvcvyAAAAgGNYVZ2b5PVJ7jH39UlJzhqOC6rqcd39iSWecUKSNyZ5/LpL3zQcT62qF3b3z28zzllJLkny19dd+vbhuKCqfqC7L99trcDRJx8FAAAA9pOMdHFTeEPe3eba1++g/w3D5133oBYAAAA4pvXhnsyx16rqO5K8IbOg6fokP5PkbyX5niT/cej2oCRvr6plcojXZC1oem+S85I8LMkzk3w8s/zmxVV1wRa1npbkrZkFTbcmeWmSRw3HS4fvvjHJ24a+wHTIRwEAAOAo2u/cU0Y6/Yx05d+Ql+SEufYtO+h/8/B54gLP2O7XoqcmuWyB8QAAAIDpuzCzX3remuQx3f2huWvvqao/yyzIOSPJ85K8eNEHVNU5SZ46nL41yQ9299eG88uq6i1J/jDJfZO8tKp+u7uv3WCoX0hy76H91O5+49y191fV5ZkFZ/dO8vNJnrForcC+kY8CAAAA++XCyEgXNoU35N001z5+B/3vMnx+ZacP6O5rtjqSfHaRggEAAIBpG7Y2ePRw+pp1QdMRL09yxdB+blXdeReP+mfD59eS/OO5oClJ0t2fT/L84fSemf0idH2t907yI8PpO9cFTUfGeWOSdw6nTxvuAaZBPgoAAAAcdTLS3ZvCgrzr5to7ebXhycPnTrZvAAAAANjIeXPt127UobsPJ3ndcHrPrIVTOzJs4fA9w+m7h0UvG3lTki8P7SdscP3vJbnTVrUOLho+7zTcA0yDfBQAAADYD+fNtWWkC1j5BXndfVOSzw+np2/Vt6rumbXA6eq9rAsAAACOSX14OsfeeuTweUNm2yFs5tK59iMWfMbDsvYmq0s369TdtyT5/SP3bPAr00fOtTcdJ8vVCuwT+SgAAAAcZfude8pI72BqGenKL8gbHHm14QOr6tAW/c7Y4B4AAACARZ05fF7V3bdu0e9jG9yz6DPWj7PVcw4l+RubjPOl7t50W8nu/l9Z+xXporUC+0s+CgAAABxtMtJd2iq8WSUfyGwl48lJHpLkDzbpd85c+4MjPv/IKw1z8+EbRxwWAACAKVj334J32qzfseDm3JT0flexvZtz0/zpqVW1Zf8ttjq4g6o6Ickpw+mW93X3F6vqhswyi/vs9BmD+f7b1Tf/pqv7JPnoBuPs5K/x6iTfmsVrBfbXyuSjt1735a36AQAAcAxa99+Cx3Q+mshIExnpsqayIO/iJP98aD89GwROVXVckqcNp9cmee+Iz//6I43f//IlIw4LAADABH19kj/f7yL2ymV5z36XsBuX7aDP1mnU7d1trn39DvofCZvuusAzFn3ODXPt9c85Ms5Oa91oDGC1XZwVyUev+ZVXjDgsAAAAE3RM56OJjHQgI13CJLas7e4PJ3n/cPrMqjp7g24/mbVXCb6iu796VIoDAAAAjjUnzLVv2UH/m4fPE/fwOTfPtdc/58g4e1krsI/kowAAAMBRJiNdwlTekJckP57ZNgsnJnlXVb0ks195npjkyUmeNfS7MsnLR372Hyc5a2j/ZZKvzV07NWurTM9Ksuk+xDAwZ1iUOcMizBcWZc6wKHOGRR0rc+ZOWXs70B/vZyF75LOZ7hamp+aO/62+rPm9Ho7fQf+7DJ9f2cPn3GWuvf45NyU5aQdjzI+zaK3A/pOPcqwwZ1iUOcOizBkWZc6wKHOGRR0Lc+ZYz0cTGel6MtIlTGZBXnd/pKqelOQ3ktw9yUs26HZlknO7+7qRn31zkss3urZu/+XPLrLfMgeTOcOizBkWYb6wKHOGRZkzLOoYmzPH7DYM3X1rkqn+vdmLuudzhZ1sW3Dy8LmT7RB2+5yT59rrn3NdZmHTXtYK7DP5KMcKc4ZFmTMsypxhUeYMizJnWNQxNGeO2Xw0kZFuQEa6hElsWXtEd781ybcn+deZhUs3Jrk2szDo+Uke3N1X7VuBAAAAwOR1901JPj+cnr5V36q6Z9YCnKsXfNR8ULblc3L7X+euf86RcbYbY36cRWsFVoB8FAAAADgaZKTLmcwb8o7o7j9P8rzhAAAAANgLVyR5ZJIHVtWh4ReyGzlj3T2L+Ogm42z1nFuTrF9s89EkD0lyj6o6tbs33Pakqv56Zm/V2k2twIqQjwIAAABHiYx0lyb1hjwAAACAo+QDw+fJmQU5mzlnrv3BBZ9xWZJbNhjndqrq+CQPP3JPd9+yrssH5tqbjpPlagUAAAAADhYZ6S5ZkAcAAABwRxfPtZ++UYeqOi7J04bTa5O8d5EHdPd1Sf6f4fR7q2qz7RSekLVfbb55g+tvSXJ4q1oH5w+fh4d7AAAAAAA2c/FcW0a6AAvyAAAAANbp7g8nef9w+syqOnuDbj+Z5Myh/Yru/ur8xao6v6p6OF64yaNeNnweSvJvq+pO68Y4JckvD6fXJvnVDWr9bJLXD6ffV1VPXN+nqn4oyfcNp7++2ZYNAAAAAACJjHQZFuQBAAAAbOzHk3wlsyDoXVX1z6vq4VX13VX16iQvHfpdmeTlu3lAd78nyW8Np38vybur6u9V1UOr6ulJfj/JfYfrL+juL24y1M8k+cuh/ZtV9UtV9Yjh+KUk/2m49pdJfnY3tQIAAAAAB46MdBcO7eXgAAAAAFPV3R+pqicl+Y3MtkN4yQbdrkxy7rC1wm49Yxj/cUm+ezjmHU7y89396i1qvbqqvj+zbSROTfL84Zj32STndfc1S9QKAAAAABwQMtLdqe7ey/EBAAAAJq2qvimzX4Kem+T0JLckuSrJG5O8qrtv3OS+85O8djh9UXe/cJvnPDXJ+Um+I8nXJflcZltCvKq7P7TDWk8Zaj0vyf2Grz+Z5JIkF3b3X+1kHAAAAACAI2Ski7EgDwAAAAAAAAAAAEZw3H4XAAAAAAAAAAAAAMcCC/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAgsyAMAAAAAAAAAAIARWJAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeUuoqvtW1cuq6oqquqGqvlBVH66qn6qqk/a7PvZeVX1DVT2+ql5cVe+oqs9XVQ/HRbsY77FV9aaquqaqbh4+31RVj92D8tkHVfWdVfXTw3y5evj7fH1VXVlVF1XVIxccz5w5hlXV3avqyVX18qq6tKquqqovVdUtVfUXVfW+qvpnVXWvHY5nvhxgVfXSuX9HdVU9egf3mDMHwLp5sdXxvh2MZc4cMFV1yvDvog9W1WeHv++fqao/qKp/VVVn72AM8wYAJkxGioyURclIWYSMlDHJSNmMjJTdko8Cm6nu3u8aJqmqzk3y+iT32KTLnyZ5XHd/4uhVxdFWVVv9P9Cvdff5OxynkvxKkmdt0e0/JPmx9v+0k1VVlyZ51A66/nqSC7r7li3GMmcOgKr63iTv3kHXzyf5ke5+5ybjmC8HXFV9R5LLkxya+/q7u/t9m/Q3Zw6Qbf48M+/S7n70JmOYMwdQVf1Qkn+fZKv/0eOS7j5vk/vNGwCYOBkpiYyUxchIWZSMlLHISNmKjJTdkI8CW/GGvF0Y/sD2hsyCpuuT/EySv5Xke5L8x6Hbg5K8varuui9Fsh+uTvKuXd77L7P2L9qPJHlKkocNnx8Zvn9Wkp9fpkD23WnD52eSvCLJEzP7+3x2kucl+fRw/R8kuWibscyZg+PqJK9L8uNJnpDZfPmuJE9K8sYkX0tySpK3VNW3bzKG+XKAVdVxmf355FCSv9jhbebMwfTvk3zbFsfTt7jXnDlgquppSX4rs7DpL5K8KMnfSfKQJOcm+T8z+x9MvrrFMOYNAEyYjJRNyEjZjoyU3ZCRshQZKQuQkbIj8lFgO96QtwtV9d4kj05ya5JHdfeH1l3/p0leOpz+XHe/+OhWyNFSVS9KclmSy7r7c1V1vySfHC7v6NefVfXAJFdk9h8Bl2c2p74yd/2kJJcmeWhmc+6M7v74mH8dHB1V9bbMQoPf6e6vbXD9lCQfTPK/DV89qrvfv0E/c+aAqKo7bTRX1vU5L8mbh9M3dff/se66+XLAVdVzk/zrJB/LbK788+HShr/+NGcOnrlff76ou1+4i/vNmQOmqs7MLBC6S5L3J/n+7v7SJn2P3+iNFuYNAEyfjJQjZKQsQkbKomSkjEFGynZkpCxCPgrshDfkLaiqzsosaEqS16wPmgYvz+wfnkny3Kq689GojaOvu3+uu9/W3Z9bYpifyNrrsZ8z/y/a4Rk3JnnOcHooyXOXeBb7qLsf391v2Cw86O7PJ/nJua+euMlQ5swBsV3QNPS5OLMQIdl4uw/z5QCrqvtk7ZdTz06y6TYvc8wZFmXOHDyvzCxs+nySJ2wWNiXJFttLmTcAMGEyUubJSFmEjJRFyUhZloyUo8ScOVjko8C2LMhb3Hlz7ddu1KG7D2f2C68kuWfWwim4nWFf+B8YTj/W3b+/Ub/h+z8dTs8b7uPY9L659gPWXzRn2MQNw+cJ81+aLyT5d0numtkbCd63XWdzhkWZMwdPVZ2R2TZ0SfKq4X8sW3QM8wYApu+8ubaMlKX48yEbeN9cW0bKTslI2YyMlD1lzhws8lFgpyzIW9wjh88bkvzhFv0unWs/Yu/KYeLun+S0oX3pVh3nrp+e5H57VRD77vi59uENrpsz3M7wWuz/fTj92LrL5ssBVlU/nOTxSb6Q5J/u8DZzhkWZMwfPD82133ikUVX3rKq/UVX32sEY5g0ATJ+MlDH58yHryUhZiIyUzchIOUrMmYNFPgrsiAV5iztz+Lyqu2/dot/8H/jP3LQXB9383Fj/H4nrmVMHwzlz7Y3mhDlDquqk4Q/1z0vy3iR3Gi69Yl1X8+WAqqqvy9p8eH53/+UObzVnDrYfqqo/raqvVNV1VfVnVfVrVfXdW9xjzhw8Dx8+v5Tkiqr6+1X1R5kF21cm+XxVfaKqfq6q7rrJGOYNAEyfjJQx+fMh68lI2ZaMlO3ISNklGSnbkY8CO3Jo+y4cUVUnJDllOL1mq77d/cWquiHJyUnus9e1MVnzc2PLOZXk6k3u4xhRVcclecHcV2/YoJs5c0BV1fnZZBugwcuSvH7dd+bLwfXSJKcm+b0kr1ngPnPmYPuWdecPHI6nVdXFSc7v7i+t62POHDxH5smnkrwyyT/ZoM/9k7wwyROr6vu6+zPrrps3ADBhMlL2gD8fchsZKVuRkbIgGSm7ISNlO/JRYEe8IW8xd5trX7+D/jcMn5utfIZF5tQNc21z6tj0E0keNrTf3N2Xb9DHnGG9/57k4d39T7u7110zXw6gqnpEkguS3JrkxzaYF1sxZw6mG5P8VpIfzWzrsQcneUySX0jyV0Of85JcUlV3XnevOXPw/LXh84zMwqZrk/xYkm9IckKSs5K8Y+jzN5O8cfgf1OaZNwAwbTJSxubPh8yTkbIb/z0yUubISNkFGSk7JR8FdsQb8hZzwlz7lh30v3n4PHEPauHYsMicunmubU4dY6rqnCS/NJz+RZJnb9LVnDm4Lk5yJIA8MckDkvxwkh9M8vqqem53v23dPebLAVNVxyf5D0kqyb/u7j9ecAhz5mA6rbuv3eD7d1fVKzMLDx6c2ZZBz07yb+b6mDMHz8nD512SfC3J3+3u35+7fnlVPT7J25L83SR/K8kTkvz2XB/zBgCmTUbK2Pz5kCQyUnbk4shI2YaMlF2SkbJT8lFgR7whbzE3zbWP30H/uwyfX9mDWjg2LDKn7jLXNqeOIVX1rUnenNki6ZuT/HB3f26T7ubMAdXd13b3nwzHZd39W939hCRPS/LNmf0q6/x1t5kvB89PJzkzyf+X5EW7uN+cOYA2CZqOXPtckidmLRR4zrou5szBM//3/I3rwqYkSXcfTvJP5756yhZjmDcAMD0yUsbmz4fISNkRGSk7JCNlYTJSFiAfBXbEgrzFXDfX3snrQI+sjt7J1g0cTIvMqZPn2ubUMaKq7p/kXUnumdmvKJ7S3ZducYs5w+10968neWNm/05/VVXdc+6y+XKAVNUZSf75cPqc7r5hq/6bMGe4g+7+RJJ3D6cPrKpvnLtszhw883/P37FZp+7+n0k+PZyetcUY5g0ATI+MlLH58+EBJyNlWTJSjpCRsldkpMyRjwI7YsvaBXT3TVX1+SSnJDl9q77DH/aP/MPx6r2ujcm6Zq695ZxKcp+5tjl1DBj+sP5fk3xjkk7yjO5+8za3mTNs5JLMtmY4ObPXX/+n4Xvz5WD5icx+SfWJJCdV1ZM36PM359p/u6pOHdpvHcIpc4bNfDTJuUP7tCSfGdrmzMFzdZIj/+y4ZquOQ9/TknzDuu/NGwCYMBkpe8CfDw8wGSkjkpGSyEjZWzJSEvkosEMW5C3uiiSPzGzl+6HuvnWTfmesuwc28tG59hmb9rrjdXNq4qrqlMx+SfPNw1fP6e7X7eBWc4aN/OVc+5vm2ubLwXLkteXfnOQ3d9D/X8y175/khpgzbK42+d6cOXj+Z9Z+0Xmnbfoeub7+v5nMGwCYPhkpY/LnwwNKRsrIZKQkMlL2loyURD4K7JAtaxf3geHz5CQP2aLfOXPtD+5dOUzcJ7P264lztuqY5FHD56eTfGqvCmLvVdU9krwzybcMX72gu//tDm83Z9jIaXPt+ddVmy8sypxhM98y1/7MXNucOXj+21z7Adv0PfI/qn163ffmDQBMn4yUMfnz4QEkI2UPyEgZiznDZmSkJPJRYIcsyFvcxXPtp2/UoaqOS/K04fTaJO/d25KYqu7uzF6jniRnVNXDN+o3fH9k9fslw31MUFWdlOTtSb5z+OoXuvuXd3q/OcMmfmiu/cdHGubLwdLd53d3bXUkedHcLd89d+1TwxjmDHdQVd+c5O8Mp5/o7tvCA3PmQHpLkq8O7Sds1qmqzklyr+H0/fPXzBsAOCZcPNeWkbIUfz48eGSk7BEZKTJS9oyMlDnyUWBHLMhbUHd/OGv/wHxmVZ29QbefTHLm0H5Fd391gz5wxIVZe03tK6vqxPmLw/krh9Nbh/5MUFUdn+TNSb5r+OoV3f2zuxjqwpgzB0JVnV9VJ2zT5yeSPG44/VTW3lJwxIUxX1jMhTFnDoyq+v6qOrTF9Xsn+e0kdx6+2uhtBRfGnDkwuvuvkvzqcPp3qurJ6/tU1d1y+7/Pr95gqAtj3gDAZMlI2QMXxp8PDwQZKYuSkbJPLow5c2DISFmEfBTYqbKIdnFV9eDMtlg4MbPXXr8ks194npjkyUmeNXS9MslDu/u6/aiTvVdVj0jywLmvTknyr4b2B7P2L+MkSXdftMk4v5jkBcPpR5L8cpKPZ/aa2+cnefBw7Re7+6fHqJ2jr6p+J2u/lHhPkucm2eofwrd095WbjGXOHABV9akkd0vyO5mFSB/P7N87d0vybUn+ftbCy1uSnNvd/3WDccwXkiRV9cIkPzecfnd3v2+TfubMATH8c+bOmf1z5kOZhdZfyezPNI9O8mNZ+xXfB5J8b3ffvME45swBUlVfn+TyJPfNLAz6lSRvSvLlzP799Pys/XLz33f3P95kHPMGACZMRsoRMlIWISNlUTJSxiYjZT0ZKYuSjwI7YUHeLlXV9yf5jSR336TLlZn9of+qo1cVR1tVXZTkH+60//Aq7I3GOS7Jf0zyjC1uf02SZ3X34UVqZHVU1aL/wP3z7r7fJmOZMwfA8B+B37SDrtckeUZ3v3uTccwXkiwUNpkzB8QC/5z5nSQXdPe1m4xjzhwwVXVmZtszPHCLbv93kh/b7G045g0ATJ+MlERGymJkpCxKRsrYZKSsJyNlN+SjwHYsyFtCVX1Tkh9Pcm6S0zP75c1VSd6Y5FXdfeM+lsdRMFbYNDfe4zL79fBZmf3q4vNJLkvy6u5+x+4rZRWMGTbNjWnOHMOq6gFJvjfJd2e2zc+9M/sV1k1JPpfkvyd5W5I37OTfOeYLOw2b5vqbM8e4qjonyTlJzk7yzZn9fb57Zr80vzrJ7yX5te7+0A7HM2cOkKo6Ocmzkzwxyd9Ictckf5HZW1Be3d3v3eE45g0ATJiMFBkpi5CRsigZKWOTkbKejJTdko8CW7EgDwAAAAAAAAAAAEZw3H4XAAAAAAAAAAAAAMcCC/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAgsyAMAAAAAAAAAAIARWJAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADACC/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAgsyAMAAAAAAAAAAIARWJAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADACC/IAAAAAAAAAAABgBBbkAQAAAAAAAAAAwAgsyAMAAAAAAAAAAIARWJAHAAAAAAAAAAAAI7AgDwAAAAAAAAAAAEZgQR4AAAAAAAAAAACMwII8AAAAAAAAAAAAGIEFeQAAAAAAAAAAADCC/x8Tfd0LLPOa2gAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "init()\n", "plot_ρs()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Run the simulation until converged" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAACeQAAAQKCAYAAAArG+jpAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjQuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/MnkTPAAAACXBIWXMAAB7CAAAewgFu0HU+AACZZElEQVR4nOzde5RsZ1kn/u9z0kkI12i4BJJogKgJyk0JA4KJikLQKIERVFSMEFHHYRRvMAMOEXRGEWdQwBkcGYL8EAQEghkQ8UIEFCZx8IIJqDHoOcQEAiaEEHJI1/v7o/c5p2i6u7q7dnXV7v581tqr3l37rXc/ffqsFc6Xp/ZbrbUAAAAAAAAAAAAA09k37wIAAAAAAAAAAABgN9CQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHACRJquqOVfWcqvrzqvpYVd3aHR+pqldU1VfOu0YAAAAAgFmQjwIA0Jdqrc27BgBgzqrqsUn+V5KTNpj22STf11p7485UBQAAAAAwe/JRAAD6pCEPAPa4qrowyfPG3vr7JJcluTXJA5J8zdi1zyT56tbah3esQAAAAACAGZGPAgDQNw15ALCHVdV/TvJz3emBJD/YWvv9VXPOTfLGJMd2b13UWvuBnasSAAAAAKB/8lEAAGZBQx4A7FFV9Y1J/jBJJfnnJF/bWvvoOnN/Nsnzu9ObkhzfWhvtSKEAAAAAAD2TjwIAMCsa8gBgD6qqo5JckeTLk7SshE3v22D+yUn2j711um0ZAAAAAIAhko8CADBL++ZdAAAwF0/KStiUJL+zUdiUJK21A0k+NfbWXWZVGAAAAADAjMlHAQCYGQ15ALA3/cDY+GWb/MznxsYTH7FbVW+rqtYdj9lSdQAAAAAAs9NrPlpVZ4xloauPz1TVlVX136rq7lNXDgDAwluadwEAwM6qqjsl+cbu9Pok793EZ5by+d/63L/e3G7+Y5M8duytByV5x5YKBQAAAADo2Yzy0a/uXj+T5ANj7x+b5EuSnN4dj6+qB7XWbtxG6QAADIQn5AHA3vO1SY7qxu9prU182l2S03Kkkf+mJNetN7Gqjk7y37vTG7rXB269TAAAAACA3s0iHz3UkPfHrbVHjh1nJjkxyY92109N8uRtVw4AwCBoyAOAveffjI2v3ORnvm5s/KcTQqpnJPmKJFcl+S/dew/adHUAAAAAALMzi3z0UEPeB1a9n7bi13Okie/UTd4TAICBsmUtAOw9Dx4bf2yTn3n82HjdrWer6q5J/nN3+tNJDm298OVVdVxr7ZZNVwkAAAAA0L9Z5KMP6l6/oCFvzKh7/ZdN3hMAgIHyhDwA2HvGA6ejJ02uqtOSPKY7/WyS12ww/ReS3CXJpa21Nyf5m+79o5J81dZLBQAAAADoVa/5aFXdN8nx3emaDXlV9Zgk90xya5I3ba1cAACGRkMeAOwhVfXFSb507K37beJjv5Ij/5vhla21T66z9gOTXJCVb3o+M0laax/PkW+ZPnA7NQMAAAAA9GFG+eih7WpvaK19ZOxex1bVl1fVf07yxiQHk1zQWvvnbRUPAMBgaMgDgL3lwavOH19VX7Te5Kp6dpJv704/meRnN1j7V7Pyvy1e1Vob/ybooafkPWhrpQIAAAAA9GoW+eihhrzjq6odOrLyNL0PJ7kwyRuS/JvW2v83TfEAAAyDhjwA2FvGA6eDWdle9nVVdfz4pKq6S1W9LMl/7d5aTvJ9rbVPrLVoVX1HkrOT3JzkOasuf7B79YQ8AAAAAGCeZpGPHlrzn5K8d+z4f1lp4qsk/zbJw3v6GQAAWHBL8y4AANhR44HTc5L8YpJHJ/mnqvqjJNcnOTkrzXW37+YtJ/nh1trb1lqwqm6X5Je7019srf3LqimHnpD3gKqq1lqb/scAAAAAANiy3vPRsTWfu/oJeFV1VJJnZiU//fWq+tvW2p/28pMAALCwyv8nDgB7R1VdkeSM7vSUrGy38JKs/9TcA0l+sLX2+xus+dwkL0iyP8lXtNZuWXX9oUne353et7X2j9v/CQAAAAAAtqfvfLSqTs5KLpokD2it/c068z6Y5CuT/Hpr7Ue3WT4AAANhy1oA2COq6vZJvrw7/WRr7UBr7deTnJXkDVnZUuFgko8nuTTJf8hKg91GzXgnJXl2d3pKks9UVRs/cqQZL0ke1OfPBAAAAACwGbPIR5N8dfd6a5IrN5h3dfd6j22WDwDAgNiyFgD2jgcmOaob/9WhN1tr703y3m2u+YtJ7pDk00lu3GDe3ZIc09Xwpm3eCwAAAABgu2aRjx5qyPvb1tptG8w7uXu9dpv3AQBgQDTkAcDe8eCx8V9Ou1hV/Zsk35OkJXlMa+3PNpj7+iRPjCfkAQAAAADz0Ws+2jnUkPdX603octQHdqdv6+m+AAAsMFvWAsDeMR44rRsQbUZVVZJfS1JJXrNRM17ng93rAzecBQAAAAAwG73lo2MONeT95eoLVbWvqp6S5OKs5KjvaK1pyAMA2AM8IQ8A9o4+A6enJHlokpuTPGsT8w815H1pVR3fWrthyvsDAAAAAGxFrw15VXW3JCd1pz9SVU8au/xFSe6d5Lju/M1Jvm/aewIAMAwa8gBgD6iqpSRf1Z1+LskVU6x1xyT/tTv9L621azbxsQ+OjR+U5F3bvT8AAAAAwFb0mY+O+eqx8emrrn06yUeSvC/Jb7XW3tXD/QAAGAgNeQCwN9wvybHd+EOttYPbXai19ukk99riZ/4uK9syAAAAAADstN7y0UNaa++IzBMAgDXsm3cBAMCO6HU7BgAAAACAAZGPAgCwY6q1Nu8aAAAAABZOVX11knOSfF1Wtre6e1a2t7omyZ8leUVr7d093/O7kvxAkgck+aIk1yZ5d5KXtdbet8k1TkjyH5Kcl+TUrDy14+okb0nya621T/RZMwAAAACw+8hHt09DHgAAAMAqVXVpkrM2MfXVSS6YdsurqrpdkjckOXedKaMkF7bWXjBhnTOTXJzknutMuSbJ41prl2+3VgAAAABgd5OPTkdDHgAAAMAqVfUPSe6blYDmDVn5FuY/JzkqycOT/GSSk7rpr22tPXnK+70myaE1/iTJr3b3vn+S/9TVkiQ/2Fr7zXXWOCnJXyS5R5Lbkvy3JJd0l89N8hNJlpJcl+RrWmsfnaZmAAAAAGB3ko9OR0MeAAAAwCpVdUmS30ryu6215TWu3zXJe5N8effWWdvdnqGqzk7yru7095I8fvye3b3+IsmXJPnXJPdprd2wxjoXJfn+7vRJrbU3rLr+xCSv705f2Vp76nbqBQAAAAB2N/nodDTkAQAAAGxDVZ2blYAoSX6ttfZj21zn/yT5liTLSU5trR1YY853JXltd/pTrbVfWXX9Hkk+mpVvqL6jtXbOOvf6/SSP6e51Umvtuu3UDAAAAADsbfLR9e2bxaIAAAAAe8C7xsb3XW/SRqrqjkke1Z2+c62wqfOmJJ/qxk9Y4/q3ZyVsSpJXbnDLi7rXo7rPAAAAAABsx7vGxvLRMRryAAAAALbnmLHxaJtrPDTJsd340vUmtdYOJnnfoc9U1dGrpnzd2HjddVZde+RmiwQAAAAAWEU+ug4NeQAAAADbc/bY+EPbXOOMLaxx6PpSki9bZ50bW2vXrrdAa+1fcuSbpGesNw8AAAAAYAL56DqWZrXwblJVxya5f3f68azsIwwAAMDecVSSu3Xjv2mt3TrPYvpWVUtJTpx3Hdt0Yjbxb/UNtjrYlqral+TZY2+9fptLnTI2nlTj/lWfu2KNdTbzc+5P8pWr7g2wLvkoAADAnrer89FERrpV8tGNacjbnPsnuWzeRQAAALAQzkxy+byL6NmJ+fwwYzeqntd7Zla2U0iSN7fWtvt34k5j409PmHvz2PiO66wzaY3xdVavAbAe+SgAAACH7MZ8NJGRbpV8dAO2rAUAAADYgqo6O8kvdqcfS/IjUyx3u7HxwQlzx795fNw660xaY3yd1WsAAAAAAGxIPjqZJ+RtzscPDR7w9f8+x9zuzvOsBQAAgB128LOfyl+/66WHTj++0dyhe9/bT8k9737UvMuY6F8+tpyHPfbwF1bPTHLtTty3qr4yyZuzkqncmuRJrbXrpljys2PjYybMPXZsfMsa69x+E2uMr7N6DYD1HP5v38NO+b4cu+QBmwAAAHvJrbd9Ou/b/+pDp7s6H01kpBuRj26OhrzNObzH8jG3u3OOPe74OZYCAADAnC1PnjJc97z7UTn5XkfPu4ytura1dmDWN6mqeyf5gyRflJW/B9/dWrt0ymVvGhtP6nC5w9h49dYLN2UlcNpMl8yhdTazfQNAMvbfvmOX7pjbHX2njeYCAACwu+3qfDSRka5HPrp5GvIAAACAw0ZpGWU07zImGqXt6P2q6l5J/jDJvZK0JE9trb25h6XHQ7KTk1y+wdxTxsb7V107kOQe3RqTHFpn9RoAAAAAsOfJSL+QfHRr9s1qYQAAAIDdoKrumuSdSe7TvfWM1tpv9bT8FWPj0yfMPXT9tiT/sM46d6mqE9dboKrumeTO3emVmy0SAAAAANib5KNbpyEPAAAAYB1VdZck70hyv+6tZ7fWXtbjLS5LcrAbn71BHcckedihz7TWDq6a8p6x8brrrLr23s0WCQAAAADsPfLR7dGQBwAAALCGqrp9kv+T5Ku7t36htfZLfd6jtXZTkj/qTr+pqtbbUuEJOfLNzbW2gnhrcngfjR/Y4Jbnd6+j7jMAAAAAAF9APrp9GvIAAACAw5bbaDDHLHXfuHxzkkd0b/1qa+2521jn/Kpq3XHhOtNe1L0uJXlZVR21ao27JjkUdN2Q5DdXL9BauzbJa7rTx1TVd6xRyxOTPKY7fXX3GQAAAABgzLxzz0XISOWj01ma1cIAAAAAA/baJI/uxn+c5BVV9VUbzD/YWvu77dyotfbHVfW6JN+V5NuTvLOqXpzkmiT3T/KcJF/STX92a+1f11nqOUnOSXK3JK+tqockuaS7dm6Sn+zGH0+y5fAMAAAAANgz5KNT0JAHAAAA8IWeMDb+xiR/PWH+PyU5dYr7PTUrWy58S5Jv6I5xoyQvaK29fL0FWmv7q+rbkrwlyYlJntUd465Ncl5r7cAUtQIAAAAAu5t8dAoa8gAAAIDDRmkZpc27jImGUONWtNZuSfKtVfXkJOcneWCS45Ncl+TdSV7aWvvzTazz/qq6f5IfS3JejoRgVye5OMmLW2uf6Ll8AAAAANg1ZKQ7b7floxryAAAAAFZprVVP61yU5KItzP/tJL895T2vT/Kz3QEAAAAAsCXy0ensm9eNAQAAAAAAAAAAYDfRkAcAAAAAAAAAAAA9sGXtFt1476UcfWd/bAAAAHvJ5z61d/4d2DLKKKN5lzFRG0CNALvS9Z9M9t067yoAAADYSaOb513BjpKRMi1PyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHuydPXcAAACAiZZby3Jr8y5joiHUCAAAAAAMj4yUaXlCHgAAAAAAAAAAAPRAQx4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9GBp3gUAAAAAi2OUllHavMuYaAg1AgAAAADDIyNlWp6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9sGUtAAAAcNgoyfIAtjoYzbsAAAAAAGBXkpEyLU/IAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHizNu4ChOXh8y+j4xd8nGgAAgP7ctm/v/DtwlJZRFv/nHUKNALtR+9zBtDp63mUAAACwg1o7OO8SdpSMlGl5Qh4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9MCWtQAAAMBhy61luS3+VgdDqBEAAAAAGB4ZKdPyhDwAAAAAAAAAAADogYY8AAAAAAAAAAAA6IGGPAAAAAAAAAAAAOjB0rwLAAAAABbHqDsW3RBqBAAAAACGR0bKtDwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6YMtaAAAA4LBRWpbT5l3GRKMB1AgAAAAADI+MlGnt2BPyququVfUzVfXeqrq2qm6tqmuq6v1V9ctV9fBNrHFOVb2pqg50nz/QnZ+zEz8DAAAAAMB2yUgBAAAAdr8deUJeVT0xyf9IcsKqS/fsjocm+bIk563z+UryP5M8fdWlk5I8Psnjq+o3kvxwa037JwAAAACwUGSkAAAAAHvDzBvyquopSV6ZlafxfSwrodN7knwyyYlJ7pvk25J8boNlfj5HgqYPJHlhkqu6z/5Mkgd31z+e5Lm9/xAAAAAAANskIwUAAADYO2bakFdVZyT5jawETe9O8m2ttRvXmPqSqjpmnTVOy0qglCSXJzmrtXZLd35ZVb01yaVJHpLkWVX1ytbaVX3+HONGSy2jJV8wBQAA2Ev20r8Dl5MsD+DHXZ53AQCbtNsy0taSlgH8hwIAAIDe7LXnsMtImda+Ga//kiTHJrk+yRPWCZqSJK21g+tcemaONA4+YyxoOvS5zyR5Rne6lOTHpykYAAAAAKBHMlIAAACAPWRmDXlVdXqSR3WnL22tXb+NNSrJ47rTD7XW3rfWvO79D3en53WfAwAAAACYGxkpAAAAwN4zyyfkPXFs/IZDg6r6oqr6sqo6YRNr3DvJSd340glzD10/Ocmpmy0SAAAAAGBGZKQAAAAAe8wsG/Ie1r3emOTKqvqeqvqrJJ9M8ndJrq+qf6yq51XVHddZ44yx8Ycm3G/8+hnrzgIAAADWNRrQATAAMlIAAAAYmHnnnjLS4Vua4dr3614/kuQlSX50jTn3TnJhku+oqse01q5Zdf2UsfGBCffbv87nJqqqkydMOXEr6wEAAAAAZCAZqXwUAAAAoD+zbMj74u719CQPTHJDkmcneVOSTyW5f5LnJ3lskq9K8oaq+rrW2ngD553Gxp+ecL+bx8brfZt0PfsnTwEAAAAA2JKhZKTyUQAAAICezLIh7w7d67FJlpM8trX2vrHrl1fVuUkuyUrg9LVJnpDkjWNzbjc2PjjhfreOjY/bVsUAAACwx41SWU7Nu4yJRgOoESAyUgAAABgcGSnTmmVD3mdzJHB6w6qgKUnSWhtV1U9nJWxKku/O54dNnx0bHzPhfseOjW/ZYq2Ttm84McllW1wTAAAAANjbhpKRykcBAAAAejLLhrybciRsevt6k1prf1tVH01yUpIz11jjkElbLNxhbDxp64bVNRzY6HqVjlIAAAAAYMsGkZHKRwEAAAD6s2+Ga+8fG28Y6IzNvfuq98c/d/KENca/xbl/3VkAAAAAADtDRgoAAACwx8zyCXl/myPf5jxqwtxD129b9f4VY+PTJ6wxfv3KCXMBAACANYzayrHohlAjQGSkAAAAMDgyUqY1yyfk/enY+L4T5t6ne/3oqvevTnJNNz57whpnja3xkUnFAQAAAADMmIwUAAAAYI+ZZUPeW5N8rhs/Yb1JVXV2khO603ePX2uttSQXd6enV9XD1lnjYTny7c+Lu88BAAAAAMyTjBQAAABgj5lZQ15r7RNJfrM7/eaq+q7Vc6rqTklePPbWy9dY6sU5sk3DS6rquFVrHJfkJd3pbavWAwAAALZgOTWYA2DRyUgBAABgeOade8pIh29pxus/L8m3JvmSJK+uqkckeVOSTyW5f5Jn5ci3Nv9Ha+2y1Qu01v6uql6U5NlJHpLkvVX1S0muyso2D89K8uBu+i+31v5+hj9PUt0BAADA3uHfgQBs3+7LSAEAAABY10wb8lprH6+qc7KyNcNpSf59d6z2v5P82AZLPSfJ3ZM8NSvB0uvWmPOKJM+dqmAAAAAAgB7JSAEAAAD2lpltWXtIa+3KJA9K8tNJ3p/kk0kOJjmQ5HeSfGNr7Wmttc9tsMaotfa0rHyT9OIk13RrXNOdf0tr7YLW2miWPwsAAAAAwFbJSAEAAAD2jllvWZskaa3dnORF3THNOm9L8rZeigIAAAC+wHIqywPYo3cINQKMk5ECAADAMMhImdbMn5AHAAAAAAAAAAAAe4GGPAAAAAAAAAAAAOjBjmxZCwAAAAxDa8moLf5WB63NuwIAAAAAYDeSkTItT8gDAAAAAAAAAACAHmjIAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeLM27AAAAAGBxLKeynJp3GRMNoUYAAAAAYHhkpEzLE/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAe2rAUAAAAOW86+LA/g+3tDqBEAAAAAGB4ZKdPymwEAAAAAAAAAAIAeaMgDAAAAAAAAAACAHmjIAwAAAAAAAAAAgB4szbsAAAAAYHG0Vhm1mncZE7UB1AgAAAAADI+MlGl5Qh4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9MCWtQAAAMBhy6ksZ/G3OhhCjQAAAADA8MhImZYn5AEAAAAAAAAAAEAPNOQBAAAAAAAAAABADzTkAQAAAAAAAAAAQA+W5l0AAAAAsDiW274st8X//t4QagQAAAAAhkdGyrT8ZgAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB7asBQAAAA4bpTIawPf3Rql5lwAAAAAA7EIyUqa1+H97AAAAAAAAAAAAYAA05AEAAAAAAAAAAEAPNOQBAAAAAAAAAABAD5bmXQAAAACwOEapLKfmXcZEowHUCAAAAAAMj4yUaXlCHgAAAAAAAAAAAPRAQx4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9GBp3gUAAAAAi2O57ctyW/zv7w2hRgAAAABgeGSkTMtvBgAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6YMtaAAAA4LBRKqPUvMuYaAg1AgAAAADDIyNlWp6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPViadwEAAADA4hhlX5YH8P290QBqBAAAAACGR0bKtPxmAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHtqwFAAAADltu+7LcFv/7e0OoEQAAAAAYHhkp0/KbAQAAAFhDVd29qs6tqudX1dur6vqqat1xUU/3+PqxNTd7vGudtT6yyc9/pI/aAQAAAIDdSz66fZ6QBwAAALC26+ZdwDo+PO8CAAAAAIBdTz66TRryAAAAACbbn+TKJI/ued3Lktx/E/NemuTsbvyqCXMvTvLcDa4f3MT9AAAAAAAOkY9ugYY8AAAA4LBRKqPsm3cZE41SO3Gb52clELqstXZdVZ2a5Oo+b9BauznJBzeaU1XHJ3lYd/oPrbU/m7DsDa21DdcEAAAAANYmIz1MPrpNGvIAAAAA1tBae968a+h8Z5Jju/Gr51kIAAAAALA3yEe3b/HbOQEAAAD2tqd0ry0DCZwAAAAAAHoyuHzUE/IAAACAw0atstx2ZDvYqYwGUGMfquq+Sb62O313a63XLSEAAAAAgM8nI10cQ81HPSEPAAAAYHE9ZWz8qk1+5qyq+uuqurmqPlNVV1fV71TVeVW1+1M6AAAAAGC3GGQ+6gl5AAAAwNCdOClHaa0d2KFa+va93estSd64yc/ce9X5qd3xpCTvrarvbK19tJfqAAAAAIBFsFsz0kHmoxryAAAAgKG7bBNzBvdkuKr6uiT36U7f3Fr71ISPHEzy1iR/kOSDSW5McnyShyf5kSSnJHlEkndW1cNbazfOom4AAAAAYMftuox0yPmohjwAAADgsOXsy3L2zbuMiYZQYw++b2z8W5uY/9DW2g1rvP+uqnppVr5B+ugkZyR5XpKfmLpCAAAAANhlZKQLY7D5qIY8AAAAYOjOTHLtvIvoU1Udm+SJ3ek1Sf5w0mfWCZsOXbupqp6U5KokJyR5elU9u7V2sIdyAQAAAID52lUZ6dDzUQ15AAAAwNBd21o7MO8ieva4rGynkCSvaa0tT7tga+3Gqnpdkh9NcockD0nyZ9OuCwAAAADM3W7LSAedj2rIAwAAAA4btX0ZtcXf6mAINU7pKWPjzWzHsFlXjI1P6nFdAAAAANgVZKQLYdD56K7+zQAAAAAMTVXdPcljutP/11r7YJ/L97gWAAAAAECvdkM+qiEPAAAAYLE8OUd2Nejz259Jcr+x8TU9rw0AAAAAMK3B56Ma8gAAAAAWy6HtGG5L8tt9LVpVd0nynd3pZ5Jc3tfaAAAAAAA9GXw+qiEPAAAAOGw5+wZzDEFVnV9VrTsu3MT8r0zy4O707a21j2/yPudU1XEbXL9TktcnOaF76xWttVs3szYAAAAA7CXzzj13U0a6V/PRpclTAAAAAPaeqnpkktPG3rrr2Pi0qjp/fH5r7aIebvv9Y+NXbeFzz07ymqp6U5L3JLkqyaeTHJ/k4Ul+JMkp3dwPJ7lw2kIBAAAAgN1LPrp9GvIAAAAA1nZBPj8AGveI7hh30TQ3q6p9SZ7cnf5rkku2uMQXZ6XmCzaY86dJntxa++TWKwQAAAAA9hD56DZpyAMAAAAOGyVZbjXvMiYazbuA2XhUkpO68e9sccuEn+o+//AkX5GVb6sen+QzSa5J8v4kr03yB6211lfBAAAAALDbyEjnZtfkoxryAAAAANbQWjs/yflTrnFRNvnN0NbaO5NsK+lrrV2e5PLtfBYAAAAAYDX56Pbtm3cBAAAAAAAAAAAAsBtoyAMAAAAAAAAAAIAe2LIWAAAAOGyUfRkN4Pt7Q6gRAAAAABgeGSnT8psBAAAAAAAAAACAHmjIAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeLM27AAAAAGBxLLd9WW6L//29IdQIAAAAAAyPjJRp+c0AAAAAAAAAAABADzTkAQAAAAAAAAAAQA9sWQsAAAAcNkpllJp3GRMNoUYAAAAAYHhkpEzLE/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHS/MuAAAAAFgco7Yvy23xv783GkCNAAAAAMDwyEiZlt8MAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPTAlrUAAADAYcupLA/g+3vLqXmXAAAAAADsQjJSprX4f3sAAAAAAAAAAABgADTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4Rq0yajXvMiYaQo0AAAAAwPDISJmWJ+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA9sWQsAAAAcNsq+LA/g+3ujAdQIAAAAAAyPjJRp+c0AAAAAAAAAAABADzTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4Rm1fRm3xv783hBoBAAAAgOGRkTItvxkAAAAAAAAAAADogYY8AAAAAAAAAAAA6IEtawEAAIDDllNZTs27jImGUCMAAAAAMDwyUqblCXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANCDpXkXAAAAACyOUduXUVv87+8NoUYAAAAAYHhkpEzLbwYAAAAAAAAAAAB6oCEPAAAAAAAAAAAAemDLWgAAAOCwUZLl1LzLmGg07wIAAAAAgF1JRsq0PCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6sDTvAgAAAIDFMWr7MmqL//29IdQIAAAAAAyPjJRp+c0AAAAAAAAAAABADzTkAQAAAAAAAAAAQA9sWQsAAAActtz2ZXkAWx0MoUYAAAAAYHhkpEzLbwYAAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgMXRUhml5l3GRG0ANQIAAAAAwyMjZVqekAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD1YmncBAAAAwOJYbvuy3Bb/+3tDqBEAAAAAGB4ZKdPymwEAAAAAAAAAAIAeaMgDAAAAAAAAAACAHtiyFgAAADhslMqo1bzLmGiUxa8RAAAAABgeGSnTmukT8qqqbfJ41ybWOqeq3lRVB6rq1u71TVV1zix/BgAAAACA7ZKRAgAAAOwtC79lba14eZK3J3l8kpOSHNO9Pj7J26vq5VWl7RMAAAAA2HVkpAAAAADDsVNb1v6PJL++wfWbN7j280me3o0/kOSFSa5Kct8kP5Pkwd31jyd57tSVAgAAAAD0T0YKAAAAsAfsVEPex1prH9zqh6rqtKwESklyeZKzWmu3dOeXVdVbk1ya5CFJnlVVr2ytXdVLxQAAALAHLaeyvPgP1M9yPAQKGBwZKQAAAAyAjJRpLfrfnmfmSNPgM8aCpiRJa+0zSZ7RnS4l+fGdKw0AAAAAYOZkpAAAAAADsrANeVVVSR7XnX6otfa+teZ173+4Oz2v+xwAAAAAwKDJSAEAAACGZ6e2rN2Oeyc5qRtfOmHupUm+IsnJSU5NcvXsygIAAIDdq7XKqC1+H0cbQI0APZCRAgAAwA6TkTKtnXpC3hOr6sNVdUtV3VRVf19Vr6qqb9jgM2eMjT80Yf3x62esOwsAAAAAYD5kpAAAAAB7wE49Ie9+q85P646nVNVbkpzfWrtx1ZxTxsYHJqy/f53PbUpVnTxhyolbXRMAAAAAYMzCZqTyUQAAAID+zLoh7zNJ3prkj7LyDc1PJ7lbkrOT/HCSE5Kcl+Tiqvrm1trnxj57p7Hxpyfc5+ax8R23Uef+yVMAAAAAALZsCBmpfBQAAACgJ7NuyDuptXbDGu+/s6pekuTtSR6clfDpR5L82tic242ND064z61j4+O2UScAAACQZJR9GWXfvMuYaAg1AnRkpAAAADAgMlKmNdOGvHWCpkPXrquq70hyZZJjkjwjnx82fXZsfMyEWx07Nr5li2Umk7dwODHJZdtYFwAAAADYwwaSkcpHAQAAAHoy6yfkbai19o9V9c4k35rktKq6V2vtmu7yTWNTJ22xcIex8aStG9aq48BG16tqq0sCAAAAAEy0CBmpfBQAAACgP3NtyOtckZWwKUlOSnIobBoPgU6esMb4Nzj391QXAAAA7DnLrbLcFr/xYgg1AmyBjBQAAAAWhIyUaS3CZsLr/e24Ymx8+oQ1xq9fOV05AAAAAAA7SkYKAAAAsEssQkPe/cbG14yNrx47P3vCGmd1rx9N8pF+ygIAAAAA2BEyUgAAAIBdYq4NeVV1nyTf3J3+Y2vto4eutdZakou709Or6mHrrPGwHPn258Xd5wAAAAAAFp6MFAAAAGB3mVlDXlV9W1UtbXD9HknemOTo7q2XrTHtxUlu68YvqarjVq1xXJKXdKe3dfMBAACAbRq1GswBsOhkpAAAADA88849ZaTDt24Y1IOXJDm6qn43yZ9nZZuEW5LcNcnXJ/nhJCd0c9+TNcKm1trfVdWLkjw7yUOSvLeqfinJVUnum+RZSR7cTf/l1trfz+qHAQAAAADYIhkpAAAAwB4zy4a8JLlXkmd0x3p+N8kFrbVb17n+nCR3T/LUrARLr1tjziuSPHeKOgEAAAAAZkFGCgAAALCHzLIh7/uTnJ3k4Unuk5Vvfd45yaeT7E/yZ0le1Vr7840Waa2Nkjyt+xbp05Oc2a11fZLLkry8tfb2Wf0QAAAAsJe0ti+jtm/eZUzUBlAjQGSkAAAAMDgyUqY1s4a81tqlSS7tcb23JXlbX+sBAAAAAMySjBQAAABg79EqCQAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9mNmWtQAAAMDwLKeynJp3GRMNoUYAAAAAYHhkpEzLE/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAe2rAUAAAAOG7Vk1BZ/q4NRm3cFAAAAAMBuJCNlWp6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPViadwEAAADA4hi1fRm1xf/+3hBqBAAAAACGR0bKtPxmAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHtqwFAAAADmupjFLzLmOiNoAaAQAAAIDhkZEyLU/IAwAAAFhDVd29qs6tqudX1dur6vqqat1xUY/3uXBs3UnH129ivROq6ueq6q+q6saq+lQ3/rmqOqGvugEAAACA3Us+un2ekAcAAACwtuvmXcBWVdWZSS5Ocs9Vlx7QHRdU1eNaa5fveHEAAAAAwJDIR7dJQx4AAADAZPuTXJnk0TO+z/0nXL96vQtVdVKS30tyjyS3JflvSS7pLp+b5CeS3CvJJVX1Na21j05fLgAAAACwB8hHt0BDHgAAAHDYcqsst5p3GRPtUI3PT3JZkstaa9dV1anZIPDpQ2vtg1N8/BeyEjYlyZNba28Yu/buqro8yeu7OS9I8tQp7gUAAAAAu5KM9DD56Dbtm9XCAAAAAEPWWntea+2S1trCb81QVfdI8r3d6TtWhU1Jku69d3SnT+k+AwAAAADwBeSj26chDwAAAGD4vj3JUd34lRvMu6h7Par7DAAAAADA0C1UPqohDwAAAGD4vm5sfOkG88avPXJGtQAAAAAA7KSFykeXZrUwAAAAMDyjti+jtvjf3xtCjdtRVe9M8tVJ7pTkhiRXJPn9JC9vrf3rBh89o3u9sbV27XqTWmv/UlWfSnLnsc8AAAAAAB0Z6fzslnxUQx4AAAAwdCdW1YYTWmsHdqiWaX3T2PhuSc7ujmdV1fmttYvX+dwp3etmfs79Sb5y7DMAAAAAwLDtlox0V+SjGvIAAACAobtsE3M2TqPm72+SvCXJ/01yTZKjk3xFku9J8ugkxyf53ar6ttba29f4/J26109v4l43d693nKJeAAAAAGBxDD0j3VX5qIY8AAAA4LBRKqO2yLnMitFCZ0db9uLW2oVrvP/+JL9VVT+U5H8mOSrJb1bVaa21W1bNvV33enAT97u1ez1uO8UCAAAAwG4mI91xuy4f1ZAHAAAADN2ZSa6ddxHb1Vq7YcL1l1fVQ5JckOReSZ6Q5DWrpn02ye2THLOJWx7bva4OrQAAAACAYRpsRrob81ENeQAAAMDQXdtaOzDvImbs5VkJnJLk7Hxh4HRTVgKnzWyzcIfudTPbNwAAAAAAi2+3Z6SDykf3zWphAAAAAHpzxdj4pDWuHwrbTt7EWqd0r/unqggAAAAAYGcMKh/1hDwAAADgsJbKKDXvMiZqA6ixZ5N+4CuSfE2Su1TVia21NbenqKp7Jrlzd3plj/UBAAAAwK4gI11Ig8pHPSEPAAAAYPHdb2x8zRrX3zM2PnuDdcavvXeqigAAAAAAdsag8lENeQAAAACL74fGxpeucf2tSUbd+Ac2WOf87nXUfQYAAAAAYNENKh/VkAcAAAAcNmrJqNUAjnn/SW1OVZ1fVa07Llzj+v2r6rQJa/xQkqd1p9cmefPqOd0WDK/pTh9TVd+xxjpPTPKY7vTV623bAAAAAAB7mYy0P3s1H12a1cIAAAAAQ1ZVj0wyHgbddWx8WlWdPz6/tXbRNm7zNUl+s6r+JMnbk/xNkk9kJbM5Pcn3Jvnmbu5ykh9qrd28zlrPSXJOkrsleW1VPSTJJd21c5P8ZDf+eJLnbqNWAAAAAGCPkI9un4Y8AAAAgLVdkOT717n2iO4Yd9E273NUkm/qjvV8IsnTWmvrbqPQWttfVd+W5C1JTkzyrO4Yd22S81prB7ZZKwAAAACwN8hHt0lDHgAAAMD8vC0r2y08PMmDk9wjyQlJKsknk/xVkt9PclFr7VOTFmutvb+q7p/kx5Kcl+TU7tLVSS5O8uLW2if6/REAAAAAALZlV+ajGvIAAACAw0ZtX0Zt37zLmGgnamytnZ/k/CnXuCgbfDO0tfaxJP+7O3rRWrs+yc92BwAAAACwBTLSFfLR7Vv8vz0AAAAAAAAAAAAwABryAAAAAAAAAAAAoAe2rAUAAAAOG7XKqNW8y5hoCDUCAAAAAMMjI2VanpAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9WJp3AQAAAMDiGKUySs27jImGUCMAAAAAMDwyUqblCXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0ANb1gIAAACHtVYZtcXf6qANoEYAAAAAYHhkpEzLE/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHS/MuAAAAAFgco1YZtZp3GRMNoUYAAAAAYHhkpEzLE/IAAAAAAAAAAACgBxryAAAAAAAAAAAAoAe2rAUAAAAOsx0DAAAAALCXyUiZlifkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPNOQBAAAAAAAAAABAD5bmXQAAAACwOEatMmo17zImGkKNAAAAAMDwyEiZlifkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPbFkLAAAAHNaSjLL4Wx20eRcAAAAAAOxKMlKm5Ql5AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQg6V5FwAAAAAsjlGrjFrNu4yJhlAjAAAAADA8MlKm5Ql5AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADDXkAAAAAAAAAAADQg6V5FwAAAAAsjlEqo1bzLmOiURa/RgAAAABgeGSkTMsT8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB7asBQAAAA4btYFsxzCAGgEAAACA4ZGRMi1PyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHmjIAwAAAAAAAAAAgB4szbsAAAAAYHGMWmXUat5lTDSEGgEAAACA4ZGRMi1PyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHtiyFgAAADiiVdoQtjoYQo0AAAAAwPDISJmSJ+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA805AEAAAAAAAAAAEAPluZdAAAAALA4RqmMUvMuY6Ih1AgAAAAADI+MlGl5Qh4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9MCWtQAAAMBho1YZtcXf6mAINQIAAAAAwyMjZVqekAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD1YmncBAAAAwOJordJazbuMiYZQIwAAAAAwPDJSpuUJeQAAAAAAAAAAANADDXkAAAAAAAAAAADQA1vWAgAAAIeNWjIawFYHozbvCgAAAACA3UhGyrQ8IQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqwNO8CAAAAgMXRWqW1mncZEw2hRgAAAABgeGSkTMsT8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB7asBQAAAA5rrTIawFYHtmMAAAAAAGZBRsq0PCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6sDTvAgAAAIDF0ZK0Nu8qJhtAiQAAAADAAMlImZYn5AEAAAAAAAAAAEAPNOQBAAAAAAAAAABAD2xZCwAAABw2SmWUmncZEw2hRgAAAABgeGSkTMsT8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgBxryAAAAAAAAAAAAoAdL8y4AAAAAWBytVVqreZcx0RBqBAAAAACGR0bKtDwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAerA07wIAAACAxTFqlVGreZcx0RBqBAAAAACGR0bKtDwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6YMtaAAAA4LDWVo5FN4QaAQAAAIDhkZEyLU/IAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHsylIa+qXlhVbez4+k185pyqelNVHaiqW7vXN1XVObOvGAAAAPaIVmkDONJq3n9SAFORkQIAAMCCWoD8U0Y6bDvekFdVD0zyzC3Mr6p6eZK3J3l8kpOSHNO9Pj7J26vq5VXlbxkAAAAAsPBkpAAAAAC714425FXVviT/K8lSko9t8mM/n+Tp3fgDSb47yUO71w907z89yQv6qxQAAAAAoH8yUgAAAIDdbWmH7/cfkpyZ5ENJ3pzkP240uapOS/Iz3enlSc5qrd3SnV9WVW9NcmmShyR5VlW9srV21UwqBwAAgD3g8HYHC24INQKsQ0YKAAAAC0xGyrR27Al5VXVKjnxD80eSHNzEx56ZI02DzxgLmpIkrbXPJHlGd7qU5MenrxQAAAAAoH8yUgAAAIDdbye3rP31JHdM8qrW2rsmTa6qSvK47vRDrbX3rTWve//D3el53ecAAAAAABaNjBQAAABgl9uRhryqelKSc5N8MslPb/Jj905yUje+dMLcQ9dPTnLqVusDAAAAAJglGSkAAADA3rA0ecp0qur4JL/anT6rtfbxTX70jLHxhybMHb9+RpKrN3kPAAAAYMyoVUZt8R+stBM1VtXdkzy0O87sjhO6y69qrZ3f033unORbkjwqydckuU+S2ye5McnfJrkkyW+21m6YsM5HknzpJm75T621U7dfMbBVMlIAAAAYDhnpCvno9s28IS/JC5OcmOTPkrxiC587ZWx8YMLc/et8blOq6uQJU07c6poAAADA4F036xtU1WOTvDnJsWtcvmuSs7vjp6rqu1trfzLrmoCZWOiMVD4KAAAArEE+uk0zbcirqkcmuSDJbUl+uLXWtvDxO42NPz1h7s1j4ztu4R6H7J88BQAAANjD9ie5Msmje173hKyETaMk70zy+0n+KskNWdl28nuSfGeSeyS5pKoe0Vr7ywlrXpzkuRtcPzhdycBWDCQjlY8CAAAAG5GPbsHMGvKq6pgkv5Gkkvz31trfbHGJ242NJ/1B3Do2Pm6L9wEAAAA6ra0ci26Hanx+ksuSXNZau66qTk3/W0B+LsnLk/yX1to/r7r2gSS/V1XvTfJrWdmm4VeysnXDRm5orX2w5zqBbZCRAgAAwPDISA+Tj27TLJ+Q95+SnJHkn5P83DY+/9mx8TET5o4/tvCWbdxr0hYOJ2blLxgAAACwR7TWnrcD9/idJL8zYc5LquopSR6S5Our6oTW2idmXRvQi6FkpPJRAAAA4PPIR7dvJg15VXV6kv/YnT6jtXbzRvPXcdPYeNIWC3cYG0/auuELtNYObHS9qra6JAAAAECf3pWVwGlfknsnWejACRhWRiofBQAAABbcuzKgfHRWT8h7Zla+sfmPSW5fVd+1xpyvGht/Y1Wd2I1/rwunxkOgkyfcb/wbnPu3WiwAAADAght/8tVoblUAWyEjBQAAAOjHoPLRWTXkHfpDuE+S125i/s+Oje+d5OYkV4y9d/qEz49fv3IT9wMAAADW0FrS2uI/Cam1eVew487uXm9L8g8T5p5VVX+d5L5JKsl1Sf5vVjKai1vbg396MB8yUgAAABggGelCGlQ+OquGvD5cneSaJPfKkT/U9ZzVvX40yUdmWBMAAACweE6ctJ3ipO0YF1lVfWuSB3Sn72itfWrCR+696vzU7nhSkvdW1Xe21j7aa5HArMhIAQAAgM3YtRnpEPPRfbNYtLV2fmutNjqS/NzYR75h7NpHujVakou766dX1cPWulf3/qFvf/qWNwAAAOw9l2Vle8aNjkGqqi9O8rLudDmf/wSt1Q4meWuSf5/k65M8OMk3JPlPOfJn8Igk76yqu8yiXuAIGSkAAACwg3ZlRjrUfHSRn5CXJC9O8oNZqfMlVXVWa+2WQxer6rgkL+lOb+vmAwAAANvUUsPYjiGLX+O0quqoJK9J8qXdWz/fWvvABh95aGvthjXef1dVvTTJG5M8OskZSZ6X5Cd6LBeYnRdHRgoAAAA7Rka6GIacj87kCXl9aa39XZIXdacPSffYwKp6SFV9Z5L3du8nyS+31v5+HnUCAAAAc3VmklMmHEP060nO6cb/J8kLNpq8Tth06NpNWdmS4RPdW0+vqmN6qBGYMRkpAAAAsAm7MSMdbD666E/IS5LnJLl7kqdm5VGCr1tjziuSPHcniwIAAAAWxrWttQPzLqJPVfVfkzy9O31Pkie21panWbO1dmNVvS7Jjya5Q1YaeP5sqkKBnSIjBQAAADayqzLSoeejC/2EvCRprY1aa09L8q1JLk5yTVb2/L2mO/+W1toFrbXRHMsEAAAA6EVVPSvJs7vT/5fk3PHtKad0xdj4pJ7WBGZMRgoAAADsFbshH53bE/JaaxcmuXAL89+W5G2zqgcAAABIWncsuiHUuB1V9e+S/GJ3emWSx7TWbuzzFj2uBUxJRgoAAACLR0Y6P7slH134J+QBAAAA7AVV9X1JXtqd/mOSb2qtXd/zbe43Nr6m57UBAAAAALZlN+WjGvIAAAAA5qyqnpDklVn5huaBJI9qrfUaCFXVXZJ8Z3f6mSSX97k+AAAAAMB27LZ8VEMeAAAAcFhrNZhjCKrq/Kpq3XHhOnMeneS1SY5K8rGsfPPzI1u8zzlVddwG1++U5PVJTujeekVr7dat3AMAAAAA9oJ55567KSPdq/no0qwWBgAAABiyqnpkktPG3rrr2Pi0qjp/fH5r7aJt3ONhSd6c5Jgkn0vyzCRHV9VXbfCxA621G1a99+wkr6mqNyV5T5Krknw6yfFJHp7kR5Kc0s39cJILt1orAAAAALB3yEe3T0MeAAAAwNouSPL961x7RHeMu2gb9zgnye278dFJXrOJz/zAOvf64qzUfMEGn/3TJE9urX1yCzUCAAAAAHuPfHSbNOQBAAAADN9PJXlUVr7t+RVZ+bbq8Uk+k+SaJO/PyrYPf9Baa3OqEQAAAABgFhYqH9WQBwAAABzRumPR7UCNrbXzk5w/5RoXZYNvhrbWLkwP2yO01i5Pcvm06wAAAADAnicjXVlePrpt++ZdAAAAAAAAAAAAAOwGGvIAAAAAAAAAAACgB7asBQAAAI5oldZq3lVMNoQaAQAAAIDhkZEyJU/IAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHizNuwAAAABgcbS2ciy6IdQIAAAAAAyPjJRpeUIeAAAAAAAAAAAA9EBDHgAAAAAAAAAAAPRAQx4AAAAAAAAAAAD0YGneBQAAAACLo7VKazXvMiYaQo0AAAAAwPDISJmWJ+QBAAAAAAAAAABADzTkAQAAAAAAAAAAQA9sWQsAAAAc0WrlWHRDqBEAAAAAGB4ZKVPyhDwAAAAAAAAAAADogYY8AAAAAAAAAAAA6IGGPAAAAAAAAAAAAOjB0rwLAAAAABZHayvHohtCjQAAAADA8MhImZYn5AEAAAAAAAAAAEAPNOQBAAAAAAAAAABAD2xZCwAAABzRumPRDaFGAAAAAGB4ZKRMyRPyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB0vzLgAAAABYHK1VWqt5lzHREGoEAAAAAIZHRsq0PCEPAAAAAAAAAAAAeqAhDwAAAAAAAAAAAHpgy1oAAADg87V5FwAAAAAAMEcyUqbgCXkAAAAAAAAAAADQAw15AAAAAAAAAAAA0AMNeQAAAAAAAAAAANCDpXkXAAAAACyO1iqt1bzLmGgINQIAAAAAwyMjZVqekAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPbBlLQAAAHBE645FN4QaAQAAAIDhkZEyJU/IAwAAAAAAAAAAgB5oyAMAAAAAAAAAAIAeaMgDAAAAAAAAAACAHizNuwAAAABgkVR3LLoh1AgAAAAADI+MlOl4Qh4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9MCWtQAAAMARrTsW3RBqBAAAAACGR0bKlDwhDwAAAAAAAAAAAHqgIQ8AAAAAAAAAAAB6oCEPAAAAAAAAAAAAerA07wIGZyj7RAMAANCfvfTvwKH8u3cINQIAAAAAwyMjZUqekAcAAAAAAAAAAAA90JAHAAAAAAAAAAAAPbBlLQAAAHBES9Jq3lVMZjsGAAAAAGAWZKRMyRPyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHGvIAAAAAAAAAAACgB0vzLgAAAABYHK2tHItuCDUCAAAAAMMjI2VanpAHAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9WJp3AQAAAMACad2x6IZQIwAAAAAwPDJSpuQJeQAAAAAAAAAAANADDXkAAAAAAAAAAADQA1vWAgAAAEe0WjkW3RBqBAAAAACGR0bKlDTkbdG+2yr7bvMXGgAAYC/x70AAWFGVVPnvIgAAwF7iX4GwNbasBQAAAAAAAAAAgB5oyAMAAAAAAAAAAIAe2LIWAAAAOKySVJt3FZPZJgMAAAAAmAUZKdPyhDwAAAAAAAAAAADogYY8AAAAAAAAAAAA6IEtawEAAIAjWncsuiHUCAAAAAAMj4yUKXlCHgAAAAAAAAAAAPRAQx4AAAAAAAAAAAD0QEMeAAAAAAAAAAAA9GBp3gUAAAAAC6TVyrHohlAjAAAAADA8MlKm5Al5AAAAAAAAAAAA0AMNeQAAAAAAAAAAANADW9YCAAAAR7TuWHRDqBEAAAAAGB4ZKVPSkLdFx9xQOXpkD2YAAIC9ZN+n/DsQAJKkjj4mte+YeZcBAADADqrR55LPzrsKGA5b1gIAAAAAAAAAAEAPNOQBAAAAAAAAAABAD2xZCwAAABzRumPRDaFGAAAAAGB4ZKRMyRPyAAAAAAAAAAAAoAca8gAAAAAAAAAAAKAHtqwFAAAAjrAdAwAAAACwl8lImZIn5AEAAAAAAAAAAEAPNOQBAAAAAAAAAABADzTkAQAAAAAAAAAAQA+W5l0AAAAAsEBarRyLbgg1AgAAAADDIyNlSp6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9sGUtAAAAcFi1lWPRDaFGAAAAAGB4ZKRMS0PeFt3l6tty7HG3zbsMAAAAdtCtt/h3IAAkSe76xcnRd5p3FQAAAOykzx2b3DTvImA4bFkLAAAAAAAAAAAAPdCQBwAAAAAAAAAAAD2wZS0AAABwROuORTeEGgEAAACA4ZGRMiVPyAMAAABYQ1XdvarOrarnV9Xbq+r6qmrdcdGM7vldVfWOqvqXqvpsVX2kql5dVQ/bwhonVNXPVdVfVdWNVfWpbvxzVXXCLOoGAAAAAHYX+ej2eUIeAAAAwNqu26kbVdXtkrwhybmrLn1pdzy5qi5srb1gwjpnJrk4yT1XXXpAd1xQVY9rrV3eT+UAAAAAwC4lH90mT8gDAAAAmGx/kj+Y4fqvyJGw6U+SnJfkoUmeluSqrGQ4z6+qC9ZboKpOSvJ7WQmbbkvywiRndccLu/fuleSSbi4AAAAAwGbIR7fAE/IAAAAA1vb8JJcluay1dl1VnZrk6r5vUlVnJ3lyd/p7SR7fWlvuzi+rqrcm+YskX5LkhVX1xtbaDWss9QtJ7tGNn9xae8PYtXdX1eVJXt/NeUGSp/b7kwAAAAAAu4h8dJs8IQ8AAABgDa2157XWLmmtzXprhp/pXpeT/LuxsOlQHdcneVZ3+kVZ+Vbo56mqeyT53u70HavCpkPrvCHJO7rTp3SfAQAAAAD4AvLR7dOQBwAAADAnVXXHJI/qTt/ZWjuwztQ3JflUN37CGte/PclR3fiVG9zyou71qO4zAAAAAABzsVvzUQ15AAAAwGGVpNoAjnn/QfXnoUmO7caXrjeptXYwyfsOfaaqjl415evGxuuus+raIzdbJAAAAADsFTLSHbUr81ENeQAAAADzc8bY+EMT5h66vpTky9ZZ58bW2rXrLdBa+5cc+SbpGevNAwAAAADYAbsyH12a1cIAAAAAO+TEqo2/D7rBVgfzdsrYeFKN+1d97oo11tnMz7k/yVeuujcAAAAAMFxDzUh3ZT6qIQ8AAAAYuss2MWdRd3C409j40xPm3jw2vuM660xaY3yd1WsAAAAAAMM01Ix0V+ajGvK26A4f/tfc7ujb5l0GAAAAO2jpczfNu4Sd02rlWHRDqHFzbjc2Pjhh7q1j4+PWWWfSGuPrrF4DYKKbv+KLcttxx8+7DAAAAHbQrbcsJVfPu4odJCPdSbsyH9WQBwAAAAzdmUmunXcR2/TZsfExE+YeOza+ZY11br+JNcbXWb0GAAAAADBMQ81Id2U+qiEPAAAAGLprW2sH5l3ENo0/fnHSFgl3GBuv3nrhpqwETpvZZuHQOpvZvgEAAAAAWHxDzUh3ZT66b1YLAwAAAAPUBnTsDuMh2ckT5p4yNt6/zjqT1hhfZ/UaAAAAAMC8c8+9lZHuynxUQx4AAADA/FwxNj59wtxD129L8g/rrHOXqjpxvQWq6p5J7tydXrnZIgEAAAAAZmBX5qMa8gAAAADm57IkB7vx2etNqqpjkjzs0GdaawdXTXnP2HjddVZde+9miwQAAAAAmIFdmY9qyAMAAACYk9baTUn+qDv9pqpab0uFJ+TINzffvMb1tyYZdeMf2OCW53evo+4zAAAAAABzsVvzUQ15AAAAwBFtQMcAVNX5VdW648J1pr2oe11K8rKqOmrVGndN8kvd6Q1JfnP1Aq21a5O8pjt9TFV9xxq1PDHJY7rTV3efAQAAAADGzTv33EUZ6V7NR5dmtTAAAADAkFXVI5OcNvbWXcfGp1XV+ePzW2sXbec+rbU/rqrXJfmuJN+e5J1V9eIk1yS5f5LnJPmSbvqzW2v/us5Sz0lyTpK7JXltVT0kySXdtXOT/GQ3/niS526nVgAAAABgb5CPbp+GPAAAAIC1XZDk+9e59ojuGHfRFPd6ala2XPiWJN/QHeNGSV7QWnv5egu01vZX1bcleUuSE5M8qzvGXZvkvNbagSlqBQAAAAB2P/noNmnIAwAAAA6rtnIsuiHUuBWttVuSfGtVPTnJ+UkemOT4JNcleXeSl7bW/nwT67y/qu6f5MeSnJfk1O7S1UkuTvLi1tonei4fAAAAAHYNGenO22356Mwa8qrqUNfimUkekuSkrDwS8Lis7Od7RZK3JXnFZn7QqjonydOTPLRb5+NJ/m+S32it/f4MfgQAAABgD2utnZ+V8GeaNS7KFr4Z2lr77SS/PeU9r0/ys90BzJGMFAAAABgq+ej2zfIJeQ9N8tp1rt0tydnd8dNV9b2ttXesNbGqKsn/zErQNO6kJI9P8viq+o0kP9xa20W9nwAAAADAwMlIAQAAAPaYWW9Zuz/JnyT5i278L0n2JTk5yXckeUKSuyZ5a1Wd2Vr76zXW+PkcCZo+kOSFSa5Kct8kP5Pkwd31jyd57sx+EgAAAACArZORAgAAAOwhs2zI+5PW2pdscP31VXVekjcnOSbJ85L82/EJVXVaVgKlJLk8yVndnsFJcllVvTXJpVnZ7uFZVfXK1tpVPf4MX+j6Tyb7bp3pLQAAAFgwo5vnXcHOad2x6IZQI8AuzEhvvPdSjr7zrL/nDQAAwCL53Kf22L8DZaRMad+sFm6tLW9izluSfKg7PWuNKc/MkabBZ4wFTYc+/5kkz+hOl5L8+HZqBQAAAADom4wUAAAAYO+ZWUPeFhx6zMDtxt+sqkryuO70Q62196314e79D3en53WfAwAAAAAYChkpAAAAwC4x14a8qjojyYO60w+tunzvJCd140snLHXo+slJTu2jNgAAANiT2oAOgF1ARgoAAAALZt65p4x08Ha8Ia+qbl9VX1ZVP5HkT5Ic1V361VVTzxgbrw6iVhu/fsa6swAAAAAA5kxGCgAAALB7Le3ETarq/CSv3GDKi5K8ZtV7p4yND0y4xf51PrcpVXXyhCknbnVNAAAAAIBDFjkjlY8CAAAA9GdHGvI28JdJfri19v41rt1pbPzpCevcPDa+4zbq2D95CgAAAABA7/4y889I5aMAAAAAPdmphry3JLm8Gx+X5L5JnpTk8UleU1U/3lq7ZNVnbjc2Pjhh/VvHxsdNUScAAADsadVWjkU3hBoBVnlLZKQAAACw8GSkTGtHGvJaazckuWHsrcuSvK6qvi/Jq5JcXFVPa61dNDbns2PjYybc4tix8S3bKHHSFg4nZqVmAAAAAIAtW/CMVD4KAAAA0JO5blnbWnt1VZ2blW+CvrSqLm6t/Wt3+aaxqZO2WLjD2HjS1g1r1XFgo+tVtdUlAQAAAAAmWoSMVD4KAAAA0J998y4gycXd6x2SPHbs/fEQ6OQJa4x/g3N/H0UBAADAntRqOAfA7iEjBQAAgEUx79xTRjp4i9CQ9/Gx8ZeOja8YG58+YY3x61dOXREAAAAAwM6RkQIAAADsEovQkHfS2Hh8K4Wrk1zTjc+esMZZ3etHk3ykn7IAAAAAAHaEjBQAAABgl1iEhrwnjo3/5tCgtdZyZKuG06vqYWt9uHv/0Lc/L+4+BwAAAAAwFDJSAAAAgF1iaVYLV9X5SV7XWvvsBnOemeRbutOPJHnPqikvTvKDWanzJVV1VmvtlrHPH5fkJd3pbd38mWqfO5hWR8/6NgAAACyQ1g7Ou4SdpY0DoBe7MSM9eHzL6Hj/oQAAANhLbtu3B/8duAd/ZPozs4a8JBcm+ZWq+t2shEhXZWW7hTsluX+S70nyiG7uwSQ/2Fq7bXyB1trfVdWLkjw7yUOSvLeqfqlb675JnpXkwd30X26t/f0Mfx4AAAAAgK24MDJSAAAAgD1llg15SfLFWfn25g9uMOdAkqe21v5wnevPSXL3JE/NSrD0ujXmvCLJc6eoEwAAAABgFmSkAAAAAHvILBvyHpXkm5J8Q5IzktwjyQlJPpvkuiR/meSSJK9vrX1mvUVaa6MkT+u+Rfr0JGcmuWuS65NcluTlrbW3z+7HAAAAgL2j2sqx6IZQI0BkpAAAADA4MlKmNbOGvNbaVVnZNuHlPa33tiRv62MtAAAAAIBZk5ECAAAA7D375l0AAAAAAAAAAAAA7AYa8gAAAAAAAAAAAKAHM9uyFgAAABig1h2Lbgg1AgAAAADDIyNlSp6QBwAAAAAAAAAAAD3QkAcAAAAAAAAAAAA9sGUtAAAAcERLaghbHQyhRgAAAP7/9u482rarrhP99xdOegJEokQTkCYWiQ0+hDCIAoFSkSKokUevYoCIUDUo6RRsakijqAivQoGlWFJEMCUFCgmNFOAjpCCCJBZ2RSCE7iWkAkZISB/Cne+PvW7Ozslp9j57nbP3OvvzGWONPddec8/1u9zJved+M/eaADA8MlJm5Al5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHFuQBAAAAAAAAAABAD1bmXcDQtJY0mzADAAAslbZM/wxs3bHohlAjwB60b6Vl34o/hAEAAJbJ0v07UEbKjDwhDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6MHKvAsAAAAAFkjrjkU3hBoBAAAAgOGRkTIjT8gDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wZS0AAABwq2qjY9ENoUYAAAAAYHhkpMzKE/IAAAAAAAAAAACgBxbkAQAAAAAAAAAAQA8syAMAAAAAAAAAAIAeWJAHAAAAAAAAAAAAPbAgDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0YGXeBQAAAAALpHXHohtCjQAAAADA8MhImZEFeQAAAAAATKa6AwAAgOXh34EwFVvWAgAAAAAAAAAAQA8syAMAAAAAAAAAAIAe2LIWAAAAuFW10bHohlAjAAAAADA8MlJm5Ql5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHtqwFAAAAbstWBwAAAADAMpORMgNPyAMAAAAAAAAAAIAeWJAHAAAAAAAAAAAAPbAgDwAAAAAAAAAAAHqwMu8CAAAAgAXSumPRDaFGAAAAAGB4ZKTMyBPyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPbFkLAAAA3Kra6Fh0Q6gRAAAAABgeGSmz8oQ8AAAAAAAAAAAA6IEFeQAAAAAAAAAAANADC/IAAAAAAAAAAACgByvzLgAAAABYIK07Ft0QagQAAAAAhkdGyow8IQ8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCWtQAAAMCtqo2ORTeEGgEAAACA4ZGRMitPyAMAAAAAAAAAAIAeWJAHAAAAAAAAAAAAPbAgDwAAAAAAAAAAAHqwMu8CAAAAgAXSumPRDaFGAAAAAGB4ZKTMyBPyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPbFkLAAAArLIdAwAAAACwzGSkzMgT8gAAAAAAAAAAAKAHFuQBAAAAAAAAAABADyzIAwAAAAAAAAAAgB6szLsAAAAAYHFUGx2Lbgg1AgAAAADDIyNlVp6QBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9GBl3gUAAAAAC6R1x6IbQo0AAAAAwPDISJmRJ+QBAAAAAAAAAABADyzIAwAAAAAAAAAAgB7YshYAAAC4LVsdAAAAAADLTEbKDDwhDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6MHKvAsAAAAAFke10bHohlAjAAAAADA8MlJm5Ql5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHtqwFAAAAVrXuWHRDqBEAAAAAGB4ZKTPyhDwAAACATVTVParqVVV1UVVdV1VfraqPV9ULq+qwGcZ9eFW1KY8PbTDWFyb8/Be2Wy8AAAAAsJxkpNPxhDwAAACADVTVKUnOSnLnsbcPS3Jid5xeVY9urX1ul0r69C7dBwAAAABARroNFuQBAAAArKOqvj/JWzMKl65N8ttJzk1yaJInJfn5JPdN8p6qOrG1du2Ut7ggyfdN0O91SU7u2n+yRd9zkvz6JtdvnuB+AAAAAAAy0m2yIA8AAAC4VbXRseh2qcYzMgqabknyyNbaR8eufbCqPpPklUmOT/L8JC+bZvDW2nVJ/mmzPlV1lyQP7k4vaa399RbDXtVa23RMAAAAAGBjMtLbOCMy0qkdMM+bAwAAACyiqjoxycO70zesCZr2e3WSi7r2c6vqwB0o5YlJDu7ab96B8QEAAAAAbkdGun0W5AEAAADc3qlj7Teu16G1ti/Jm7rTI7MaTvXpqftvl4GETQAAAADAnnDqWFtGOgUL8gAAAIBVbUDHznpo93pdkr/dpN95Y+2H9FlAVd0nyQ92px9urX2+z/EBAAAAgHXMO/eUkd5qqBmpBXkAAAAAt3dC93pJa+2WTfp9ap3P9OWpY+0/mfAzD6uqf6iq66rq+qr6fFX996o6taqq5/oAAAAAgL1LRrpNK7txEwAAAIAddPRWOUpr7bJJB6uqQ5Ic1Z1u+rnW2teq6rokhye5+6T3mNDPdK83JPnzCT9zrzXn9+yOJyQ5v6qe2Fr7Ui/VAQAAAACLQka6au4ZqQV5AAAAwNBdMEGfab75eMRY+9oJ+u8Pm+44xT02VVUPTXLv7vQdrbWvb/GRm5O8M8n7k/xTkquT3CXJSUmenVEQ9kNJPlBVJ7XWru6rVgAAAABg7mSkC5SRWpAHAAAArGrdseh2tsZDxto3T9D/pu710B5r+Nmx9psm6P+g1tpV67z/oap6XUbfHn1kRltG/EaS589cIQAAAADsRTLSREY6EwvyAAAAgKE7MckVPY5341j7oAn6H9y93tDHzavq4CSP704vT/JXW31mg6Bp/7VrquoJST6b5K5JnllVL26tTRKkAQAAAACLT0a6QBmpBXkAAADA0F3RWrusx/GuGWtPssXC4d3rJFs3TOInM9pKIUnOaq19c9YBW2tXV9Vbkvy7jOp9YJK/nnVcAAAAAGAhyEi3sJsZqQV5AAAAwK2qOxbdTtbYWruxqq5MclSSYzeto+rIrIZNl/ZUwlPH2pNsxTCpT461j+lxXAAAAADYM2SkMtJZHbBTAwMAAAAM2EXd63FVtdkXGo9f5zPbVlXfluTHutP/1Vr7p1nHHB++x7EAAAAAgL1NRrpNFuQBAAAA3N5HutfDkzxgk34nj7XP7+G+T8nqjgZ9fvMzSb57rH15z2MDAAAAAHuLjHSbLMgDAAAAuL2zx9pPW69DVR2Q1a0Trkpybg/33T/eLUn+Ww/jJUmq6s5JntidXp/kwr7GBgAAAAD2pLPH2jLSKViQBwAAAKxqAzp2UGvt40k+3J0+o6pOWqfbC5Kc0LVf01r7xvjFqjqtqlp3vGSre1bV9yS5f3f63tbaP09Sa1U9qqoO3eT6EUnemuSu3VtvaK3dNMnYAAAAALB05p17ykgHn5Futr8vAAAAwDL7xYy2WDg0yfur6hUZfcPz0CRPSvLMrt/FSV7dw/1+bqz9J1N87sVJzqqqt2e0jcRnk1yb5C5JTkry7CR37/p+OslLZi0UAAAAAFgKMtJtsCAPAAAAYB2ttU9U1ROT/GmSOyV5xTrdLk5ySmvtmlnu1W3t8JTu9GtJ3j3lEN+S5PTu2Mj/TPKU1tpXp68QAAAAAFg2MtLtsSAPAAAAWNWS2uGtDnqxSzW21t5VVffL6JugpyQ5NsnNSS5J8rYkr2utXd/DrX44yTFd+79PuV3CC7vPn5TkvkmOyuibn9cnuTzJ3yT5syTvb60N4XcXAAAAAOZHRnrb28hIp2ZBHgAAAMAmWmtfTPL87pjmc2cmOXPCvh9IUtPW1n32wiQXbuezAAAAAABbkZFO54B5FwAAAAAAAAAAAAB7gQV5AAAAAAAAAAAA0ANb1gIAAACrWncsuiHUCAAAAAAMj4yUGXlCHgAAAAAAAAAAAPTAgjwAAAAAAAAAAADogQV5AAAAAAAAAAAA0IOVeRcAAAAALJg27wIAAAAAAOZIRsoMPCEPAAAAAAAAAAAAemBBHgAAAAAAAAAAAPTAlrUAAADAraqNjkU3hBoBAAAAgOGRkTIrT8gDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6sDLvAgAAAIAF0rpj0Q2hRgAAAABgeGSkzMgT8gAAAAAAAAAAAKAHFuQBAAAAAAAAAABAD2xZCwAAANyq2uhYdEOoEQAAAAAYHhkps/KEPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcr8y4AAAAAWCCtOxbdEGoEAAAAAIZHRsqMPCEPAAAAAAAAAAAAemBBHgAAAAAAAAAAAPTAlrUAAADAraqNjkU3hBoBAAAAgOGRkTIrT8gDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6sKML8qrqB6rqV6vqvVV1aVXdVFXXVtXFVXVmVT10yvEeVVVvr6rLurEu684ftVO/BgAAAFgqbUAHwADISAEAAGBg5p17ykgHb2WnBq6q85I8bJ1LByX5ru74uap6c5LTW2s3bzJWJfnDJM9cc+mYJD+V5Keq6o+SPKu1ZroBAAAAAHMnIwUAAABYPjv5hLxjutfLk7wmyeOSPCjJSUmen+RL3fWfTXLmFmP9ZlaDpk8keXI31pO783TXX95D3QAAAAAAfZCRAgAAACyZHXtCXpJPJfnVJH/RWvvmmmsf6771eX6Sf5XkyVX1B621D68dpKqOS/LL3emFSR7WWruhO7+gqt6Z5LwkD0zyoqp6Y2vtszvw6wEAAIC9byhbHQyhRgAZKQAAAAyPjJQZ7dgT8lprj2mtvXWdoGn/9SuTvGDsrcdtMNTzsrpw8DljQdP+ca5P8pzudCXJc7ddNAAAAABAT2SkAAAAAMtnJ7esncSHxtr3WXuxqirJT3ann2qtfWy9Qbr3P92dntp9DgAAAABg0X1orC0jBQAAABi4eS/IO2isvW+d6/dKckzXPm+LsfZfPzbJPWcrCwAAAABgV8hIAQAAAPaQla277KiTx9qfWuf6CVtczwbXT0jy+e0WBQAAAMuq2uhYdEOoEWBCMlIAAABYIDJSZjW3BXlVdUCSF4+99dZ1ut19rH3ZFkNeusHnJqnl2C26HD3NeAAAAAAAW1mUjFQ+CgAAANCfeT4h73lJHtS139Fau3CdPkeMta/dYrzrxtp3nLKWS7fuAgAAAADQq0XJSOWjAAAAAD2Zy4K8qjo5ye90p19J8uwNuh4y1r55i2FvGmsfus3SAAAAAFsdAOw4GSkAAAAsMBkpM9j1BXlV9T1J3tHd+6YkT2itfXmD7jeOtQ/aYuiDx9o3TFnWVts3HJ3kginHBAAAAAC4nQXMSOWjAAAAAD3Z1QV5VXWvJO9PcmSSbyZ5cmvtvE0+cs1Ye6stFg4fa2+1dcNttNYu2+x6VU0zHAAAAADAuhYxI5WPAgAAAPTngN26UVV9R5K/SvIdGT3Y8emttXds8bHxIOjYLfqOf4vz0ukrBAAAAADYOTJSAAAAgL1vV56QV1VHJflAknt3bz2ntfamCT76ybH28Vv0Hb9+0RTlAQAAAJ1qLdXavMvY0hBqBBgnIwUAAIBhkJEyqx1/Ql5V3TnJ+5J8d/fWi1trvz/hxz+f5PKuffIWfR/WvX4pyRemqREAAAAAYKfISAEAAACWx44uyKuqw5K8J8kPdG/9Vmvtdyf9fGutJTmnOz2+qh68wX0enNVvf57TfQ4AAAAAYK5kpAAAAADLZccW5FXVQUnekeSHurde01r79W0MdUaSW7r2a6vq0DX3OTTJa7vTW7r+AAAAwHa0AR0AC05GCgAAAAM079xTRjp4Kzs49p8leWTX/mCSN1TV927S/+bW2sVr32ytXVxVr0ry4iQPTHJ+Vf1uks8muU+SFyW5f9f991prn+nrFwAAAAAAMAMZKQAAAMCS2ckFeY8da//rJP+wRf8vJrnnBtd+Lcm3JXl6RsHSW9bp84Yk2/l2KQAAAADATpCRAgAAACyZHduytk+ttX2ttWckOSXJOUkuT3Jz93pOkke31k5vre2bY5kAAAAAADtCRgoAAAAwDDv2hLzWWu3AmH+Z5C/7HhcAAAAYqTY6Ft0QagSQkQIAAMDwyEiZ1SCekAcAAAAAAAAAAACLzoI8AAAAAAAAAAAA6IEFeQAAAAAAAAAAANCDlXkXAAAAACyQ1h2Lbgg1AgAAAADDIyNlRp6QBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YMtaAAAA4FbVRseiG0KNAAAAAMDwyEiZlSfkAQAAAAAAAAAAQA8syAMAAAAAAAAAAIAeWJAHAAAAAAAAAAAAPViZdwEAAADAAmndseiGUCMAAAAAMDwyUmbkCXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAe2rAUAAABuVW10LLoh1AgAAAAADI+MlFl5Qh4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6IEFeQAAAAAAAAAAANCDlXkXAAAAACyQ1h2Lbgg1AgAAAADDIyNlRp6QBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YMtaAAAA4DbKVgcAAAAAwBKTkTILT8gDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6sDLvAgAAAIAF0troWHRDqBEAAAAAGB4ZKTPyhDwAAAAAAAAAAADogQV5AAAAAAAAAAAA0ANb1gIAAAC3qjY6Ft0QagQAAAAAhkdGyqw8IQ8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOjByrwLAAAAABZI645FN4QaAQAAAIDhkZEyI0/IAwAAAAAAAAAAgB5YkAcAAAAAAAAAAAA9sGUtAAAAcKvaNzoW3RBqBAAAAACGR0bKrDwhDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6MHKvAsAAAAAFkjrjkU3hBoBAAAAgOGRkTIjT8gDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wZS0AAABwq2qjY9ENoUYAAAAAYHhkpMzKE/IAAAAANlFV96iqV1XVRVV1XVV9tao+XlUvrKrDZhz7JVXVJjwePsF4d62ql1bV31fV1VX19a790qq66yy1AgAAAADLSUY6HU/IAwAAANhAVZ2S5Kwkdx57+7AkJ3bH6VX16Nba5+ZR37iqOjHJOUm+fc2l+3XH6VX1k621C3e9OAAAAABgkGSk07MgDwAAAGAdVfX9Sd6aUbh0bZLfTnJukkOTPCnJzye5b5L3VNWJrbVrZ7zl921x/fOb1HpMkncluVuSW5L8P0ne3V1+TJLnJ/mOJO+uqge01r40Y60AAAAAwB4nI90eC/IAAACAVa2NjkW3OzWekVHQdEuSR7bWPjp27YNV9Zkkr0xyfEZhzstmuVlr7Z9m+PhvZRQ0JclTWmtvG7v24aq6MKPg7G5JXp7k6TPcCwAAAAD2LhnpuDMiI53aATs1MAAAAMBQdVsbPLw7fcOaoGm/Vye5qGs/t6oO3I3a1qqquyX5me70fWuCpiRJ9977utOndp8BAAAAAFiXjHT7LMgDAAAAuL1Tx9pvXK9Da21fkjd1p0dmNZzabT+R5A5de91aO2d2r3foPgMAAAAAsJFTx9oy0ilYkAcAAABwew/tXq9L8reb9DtvrP2QnStnUw8da5+3Ya/FqBUAAAAAGAYZ6Tat7NTAAAAAwPBUGx2LbhdqPKF7vaS1dssm/T61zme2pao+kOQHkhyR5Kokn0zyP5K8vrX2tU0+uv++V7fWrtioU2vt/1TV15PcadZaAQAAAGCvkpHeSka6TRbkAQAAAEN3dFVt2qG1dtmkg1XVIUmO6k43/Vxr7WtVdV2Sw5PcfdJ7bOBHxtrfmuTk7nhRVZ3WWjtng8/tv+8kv8ZLk3xPZq8VAAAAAFgcMtIJau3seEZqQR4AAAAwdBdM0GfzNOq2jhhrXztB//1h0x2nuMe4f0xydpKPJ7k8yYFJ7pvkp5M8MsldkvxFVf14a+29m9Q7aa2ZoVYAAAAAYPHISKerNTPUuiUL8gAAAIDbGsB2DDvskLH2zRP0v6l7PXQb9zqjtfaSdd7/myRvqqpfSPKHSe6Q5I+r6rjW2g1r+u6vd6drBQAAAIDlICOVkc7AgjwAAABg6E5MckWP49041j5ogv4Hd69rQ6Attdau2uL666vqgUlOT/IdSR6b5Kw13W5Mclh2uFYAAAAAYGHJSBcoI7UgDwAAABi6K1prl/U43jVj7Um2LTi8e51kO4TteH1GYVOSnJzbh03XZBQ2LUKtAAAAAMDuk5EuUEZ6wE4NDAAAADBErbUbk1zZnR67Wd+qOjKrAc6lO1TSJ8fax6xzfX/Qtmmtnbt3rztVKwAAAAAwcDLS2ViQBwAAANyq2nCOHXZR93pcVW22w8Dx63ymb7XF9f1h1J2r6ugNB6n69iR36k53qlYAAAAAGLR5554y0nUNKiO1IA8AAADg9j7SvR6e5AGb9Dt5rH3+DtXy3WPty9e5/pGx9snrXF/v2k7VCgAAAADsDTLSbbIgDwAAAOD2zh5rP229DlV1QJKndqdXJTl3h2r5hbH2eetcf2eSfV173Vo7p3Wv+7rPAAAAAABs5Oyxtox0ChbkAQAAAKtaG86xo/8ztI8n+XB3+oyqOmmdbi9IckLXfk1r7RvjF6vqtKpq3fGStR+uqu+rquM2q6OqfiHJM7rTK5K8Y51ar0hyVnf6Y1X1uHXGeXySH+tO39x9BgAAAABYa965p4x0bZ/BZaSb7e8LAAAAsMx+MaNtCw5N8v6qekVG3/A8NMmTkjyz63dxkldvY/wHJPnjqjo3yXuT/GOSf8korzk+yc8k+dGu7zeT/EJr7boNxvq1JI9K8q1J/qyqHpjk3d21x2QUjCXJPyf59W3UCgAAAAAsHxnpNliQBwAAALCO1tonquqJSf40yZ2SvGKdbhcnOaW1ds02b3OHJD/SHRv5lyTPaK1tuIVCa+3SqvrxjLaRODrJi7pj3BVJTm2tXbbNWgEAAACAJSIj3R4L8gAAAAA20Fp7V1XdL6Nvgp6S5NgkNye5JMnbkryutXb9Nof/y4y2Wjgpyf2T3C3JXZNUkq8m+fsk/yPJma21r09Q699U1fd1tZ6a5J7dpc8nOSfJGa21f9lmrQAAAADAEpKRTs+CPAAAAOBW1UbHotvNGltrX0zy/O6Y5nNnJjlzk+tfSfJfu6MXrbUrk/yH7gAAAAAApiQjvT0Z6XQOmNeNAQAAAAAAAAAAYC+xIA8AAAAAAAAAAAB6YMtaAAAAYFXrjkU3hBoBAAAAgOGRkTIjT8gDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6sDLvAgAAAIDFUW10LLoh1AgAAAAADI+MlFl5Qh4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6IEtawEAAIBV+9roWHRDqBEAAAAAGB4ZKTPyhDwAAAAAAAAAAADogQV5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHK/MuAAAAAFggrTsW3RBqBAAAAACGR0bKjDwhDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0wJa1AAAAwK0qSQ1gq4OadwEAAAAAwJ4kI2VWnpAHAAAAAAAAAAAAPbAgDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0YGXeBQAAAAALpCVpbd5VbG0AJQIAAAAAAyQjZUaekAcAAAAAAAAAAAA9sCAPAAAAAAAAAAAAemDLWgAAAGBVS2oIWx0MoUYAAAAAYHhkpMzIE/IAAAAAAAAAAACgBxbkAQAAAAAAAAAAQA8syAMAAAAAAAAAAIAerMy7AAAAAGCBtO5YdEOoEQAAAAAYHhkpM/KEPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcr8y4AAAAAWBzVWqq1eZexpSHUCAAAAAAMj4yUWXlCHgAAAAAAAAAAAPTAgjwAAAAAAAAAAADogS1rAQAAgFX7umPRDaFGAAAAAGB4ZKTMyBPyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHqzMuwAAAABgcVRrqdbmXcaWhlAjAAAAADA8MlJm5Ql5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHtqwFAAAAVrXuWHRDqBEAAAAAGB4ZKTPyhDwAAAAAAAAAAADogQV5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHK/MuAAAAAFgkLWlt3kVMYAg1AgAAAADDIyNlNp6QBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YMtaAAAA4FbVRseiG0KNAAAAAMDwyEiZlSfkAQAAAAAAAAAAQA8syAMAAAAAAAAAAIAeWJAHAAAAAAAAAAAAPViZdwEAAADAAmltdCy6IdQIAAAAAAyPjJQZeUIeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOiBLWsBAACAW9W+0bHohlAjAAAAADA8MlJmtaNPyKuqb6uqx1TVy6rqvVV1ZVW17jhzG+M9qqreXlWXVdVN3evbq+pRO1A+AAAAAMBMZKQAAAAAy2Wnn5D35T4GqapK8odJnrnm0jFJfirJT1XVHyV5Vmut9XFPAAAAAIAeyEgBAAAAlsiOPiFvjUuTvH+bn/3NrAZNn0jy5CQP6l4/0b3/zCQvn6VAAAAAAIAdJCMFAAAA2ON2+gl5L0tyQZILWmtfrqp7Jvn8NANU1XFJfrk7vTDJw1prN3TnF1TVO5Ocl+SBSV5UVW9srX22l+oBAABg2bQ2OhbdEGoEGJGRAgAAwJDISJnRjj4hr7X2G621d7fWZtmW4XlZXTj4nLGgaf89rk/ynO50JclzZ7gXAAAAAEBvZKQAAAAAy2U3t6ydWlVVkp/sTj/VWvvYev269z/dnZ7afQ4AAAAAYNBkpAAAAADDstAL8pLcK8kxXfu8Lfruv35sknvuVEEAAACwp7UBHQDLQUYKAAAAu2neuaeMdPAWfUHeCWPtT23Rd/z6CRv2AgAAAAAYDhkpAAAAwICszLuALdx9rH3ZFn0v3eBzW6qqY7focvQ04wEAAAAA9GTHM1L5KAAAAEB/Fn1B3hFj7Wu36HvdWPuOU97n0q27AAAAAADsut3ISOWjAAAAAD1Z9AV5h4y1b96i701j7UN3oBYAAADY86q1VGvzLmNLQ6gRoCcyUgAAANhFMlJmtegL8m4cax+0Rd+Dx9o3THmfrbZvODrJBVOOCQAAAAAwq93ISOWjAAAAAD1Z9AV514y1t9pi4fCx9lZbN9xGa+2yza5X1TTDAQAAAAD0ZcczUvkoAAAAQH8WfUHeeBB07BZ9x7/FeekO1AIAAAB7X2ujY9ENoUaAfshIAQAAYDfJSJnRAfMuYAufHGsfv0Xf8esX7UAtAAAAAAC7TUYKAAAAMCCLviDv80ku79onb9H3Yd3rl5J8YacKAgAAAADYRTJSAAAAgAFZ6AV5rbWW5Jzu9PiqevB6/br393/785zucwAAAAAAgyYjBQAAABiWhV6Q1zkjyS1d+7VVdej4xe78td3pLV1/AAAAYDtakn0DOCwzAZbLGZGRAgAAwO6QkTKjlZ0cvKoekuS4sbeOGmsfV1WnjfdvrZ25dozW2sVV9aokL07ywCTnV9XvJvlskvskeVGS+3fdf6+19pnefgEAAAAAADOQkQIAAAAslx1dkJfk9CQ/t8G1H+qOcWdu0PfXknxbkqdnFCy9ZZ0+b0jy69OXCAAAAACwY2SkAAAAAEtkCFvWprW2r7X2jCSnJDknyeVJbu5ez0ny6Nba6a21fXMsEwAAAABgR8hIAQAAAIZhR5+Q11o7LclpPY73l0n+sq/xAAAAgNuq1lKtzbuMLQ2hRoBERgoAAABDIyNlVoN4Qh4AAAAAAAAAAAAsOgvyAAAAAAAAAAAAoAc7umUtAAAAMDAtyRC2OhhAiQAAAADAAMlImZEn5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHliQBwAAAAAAAAAAAD1YmXcBAAAAwAJpbXQsuiHUCAAAAAAMj4yUGXlCHgAAAAAAAAAAAPTAgjwAAAAAAAAAAADogS1rAQAAgFX7umPRDaFGAAAAAGB4ZKTMyBPyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHqzMuwAAAABgcVRrqdbmXcaWdrPGqrpHkn+f5JQk90hyU5JLkrw1yX9urV0/w9h3SvLoJD+c5AFJ7p3ksCRXJ/nfSd6d5I9ba1dtMc4XknznBLf8YmvtntutFwAAAAD2OhnpOveSkU7FgjwAAACADVTVKUnOSnLnsbcPS3Jid5xeVY9urX1uG2P/myTvSHLwOpePSnJyd7ywqp7cWjt32nsAAAAAAMxCRjo9C/IAAAAA1lFV35/RNzwPS3Jtkt9Ocm6SQ5M8KcnPJ7lvkvdU1YmttWunvMVdMwqa9iX5QJL/keTvk1yV5NgkP53kiUnuluTdVfVDrbW/22LMc5L8+ibXb56yRgAAAABgSclIt8eCPAAAAGBVa6Nj0e1OjWdkFDTdkuSRrbWPjl37YFV9Jskrkxyf5PlJXjbl+N9I8vokr2it/X9rrn0iybuq6vwk/6mr49UZbduwmataa/80ZR0AAAAAwH4y0nFnREY6tQPmeXMAAACARVRVJyZ5eHf6hjVB036vTnJR135uVR04zT1aa/+9tfasdYKm8T6vTXJhd/rwqrrrNPcAAAAAANgOGen2WZAHAAAAcHunjrXfuF6H1tq+JG/qTo/MajjVtw91rwckudcO3QMAAAAAYNypY20Z6RQsyAMAAAC4vYd2r9cl+dtN+p031n7IDtVy8Fh73w7dAwAAAABgnIx0m1bmXQAAAACwSFrS2ryLmMCO13hC93pJa+2WTfp9ap3P9O3k7vWWJJds0fdhVfUPSe6TpJJ8OcnHk/xZknNaG8RvLgAAAADMkYy0IyPdJgvyAAAAgKE7uqo27dBau2zSwarqkCRHdaebfq619rWqui7J4UnuPuk9pqjllCT3607f11r7+hYfWbtdwz274wlJzq+qJ7bWvtRrkQAAAADAvMlIV809I7UgDwAAABi6Cybos3kadVtHjLWvnaD//rDpjlPcY0tV9S1Jfr87/WaS/7BJ95uTvDPJ+5P8U5Krk9wlyUlJnp1REPZDST5QVSe11q7us1YAAAAAYK5kpAuUkVqQBwAAAKxqA9mOYWdrPGSsffME/W/qXg/tq4CqukOSs5J8Z/fWb7bWPrHJRx7UWrtqnfc/VFWvS/LnSR6Z0ZYRv5Hk+X3VCgAAAAB7iow0kZHOxII8AAAAYOhOTHJFj+PdONY+aIL+B3evN/RYw39O8qiu/Z4kL9+s8wZB0/5r11TVE5J8Nsldkzyzql7cWpskSAMAAAAAFp+MdIEyUgvyAAAAgKG7orV2WY/jXTPWnmSLhcO710m2bthSVf12kmd2px9J8vjW2jdnGbO1dnVVvSXJv8uo3gcm+euZCgUAAAAAFoWMdAu7mZEesBODAgAAAAxVa+3GJFd2p8du1reqjsxq2HTprPeuqhcleXF3+r+SPKa11te3Sj851j6mpzEBAAAAgD1GRjobT8gDAAAAVu3rjkW38zVelOShSY6rqpXW2i0b9Dt+zWe2rar+bZLfGRvrx1prV88y5tpb9DgWAAAAAOxNMtL9ZKTb5Al5AAAAALf3ke718CQP2KTfyWPt87d7s6r62SSv604/l+RHWmtXbvKR7fjusfblPY8NAAAAAOwtMtJtsiAPAAAA4PbOHms/bb0OVXVAkqd2p1clOXc7N6qqxyZ5Y0bfzrwsyQ+31noNg6rqzkme2J1en+TCPscHAAAAAPacs8faMtIpWJAHAAAA3KpaG8yxk1prH0/y4e70GVV10jrdXpDkhK79mtbaN8YvVtVpVdW64yXr3aeqHpnkz5LcIclXMvrW5xemqbWqHlVVh25y/Ygkb01y1+6tN7TWbprmHgAAAACwLOade8pIh5+RruzUwAAAAAAD94sZbbFwaJL3V9UrMvqG56FJnpTkmV2/i5O8etrBq+rBSd6R5KAk30jyvCQHVtX3bvKxy1prV61578VJzqqqt2e0jcRnk1yb5C5JTkry7CR37/p+OslLpq0VAAAAAFhKMtJtsCAPAAAAYB2ttU9U1ROT/GmSOyV5xTrdLk5ySmvtmm3c4lFJDuvaByY5a4LPPC3Jmeu8/y1JTu+OjfzPJE9prX11ihoBAAAAgCUlI90eC/IAAAAANtBae1dV3S+jb4KekuTYJDcnuSTJ25K8rrV2/RxLTJIXJvnhjL7ped8kR2X0zc/rk1ye5G8y2vLh/a3t8D4WAAAAAMCeIiOdngV5AAAAwKrWRsei28UaW2tfTPL87pjmc2dm/W9q7r/+kvSwNUJr7cIkF846DgAAAAAQGem6t5KRTuOAeRcAAAAAAAAAAAAAe4EFeQAAAAAAAAAAANADW9YCAAAAq/a10bHohlAjAAAAADA8MlJm5Al5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHFuQBAAAAAAAAAABAD1bmXQAAAACwQFqS1uZdxdYGUCIAAAAAMEAyUmbkCXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPVuZdAAAAALBIWtLavIuYwBBqBAAAAACGR0bKbDwhDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0wJa1AAAAwKo2kO0YhlAjAAAAADA8MlJm5Al5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHFuQBAAAAAAAAAABAD1bmXQAAAACwQPa10bHohlAjAAAAADA8MlJm5Al5AAAAAAAAAAAA0AML8gAAAAAAAAAAAKAHtqwFAAAAVrV9o2PRDaFGAAAAAGB4ZKTMyBPyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHqzMuwAAAABggbQ2OhbdEGoEAAAAAIZHRsqMPCEPAAAAAAAAAAAAemBBHgAAAAAAAAAAAPTAlrUAAADAqn1tdCy6IdQIAAAAAAyPjJQZeUIeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQg5V5FwAAAAAskNZGx6IbQo0AAAAAwPDISJmRJ+QBAAAAAAAAAABADyzIAwAAAAAAAAAAgB7YshYAAAC4LVsdAAAAAADLTEbKDDwhDwAAAAAAAAAAAHpgQR4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6MHKvAsAAAAAFkhro2PRDaFGAAAAAGB4ZKTMyBPyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPbFkLAAAArNq3L6l9865ia/sGUCMAAAAAMDwyUmbkCXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPVuZdAAAAALBAWhsdi24INQLsRa07AAAAWB7L9u9AGSkz8oQ8AAAAAAAAAAAA6IEFeQAAAAAAAAAAANADW9YCAAAAq2zHAAAAAAAsMxkpM/KEPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcr8y4AAAAAWCD7WlJt3lVsbd8AagQAAAAAhkdGyow8IQ8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOjByrwLAAAAABZHa/vS2r55l7GlIdQIAAAAAAyPjJRZeUIeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOiBLWsBAACAVS3JvjbvKrY2gBIBAAAAgAGSkTIjC/KmVJVU1bzLAAAAYBf5VyAAjBxwS+WAW/zNCAAAsEz8OxCmY8taAAAAAAAAAAAA6IEFeQAAAAAAAAAAANADW9YCAAAAq1obHYtuCDUCAAAAAMMjI2VGnpAHAAAAAAAAAAAAPbAgDwAAAAAAAAAAAHpgy1oAAABg1b59SfbNu4qt7RtAjQAAAADA8MhImZEn5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHliQBwAAAAAAAAAAAD1YmXcBAAAAwAJpbXQsuiHUCAAAAAAMj4yUGXlCHgAAAAAAAAAAAPTAgjwAAAAAAAAAAADogS1rAQAAgFu1ffvSsm/eZWyp7Vv8GgEAAACA4ZGRMisL8qZUBx6UOuCgeZcBAADALqp930hunHcVADB/B11VOXBfzbsMAAAAdtEBX/fvQJiGLWsBAAAAAAAAAACgBxbkAQAAAAAAAAAAQA9sWQsAAACsam10LLoh1AgAAAAADI+MlBl5Qh4AAAAAAAAAAAD0wII8AAAAAAAAAAAA6IEtawEAAIBVrSX7BrDVge0YAAAAAICdICNlRp6QBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9GBl3gUAAAAAC6S1JPvmXcXWWpt3BQAAAADAXiQjZUaDe0JeVd2jql5VVRdV1XVV9dWq+nhVvbCqDpt3fQAAAAAAO0U+CgAAALDYBvWEvKo6JclZSe489vZhSU7sjtOr6tGttc/Noz4AAAAAgJ0iHwUAAABYfINZkFdV35/krRkFTNcm+e0k5yY5NMmTkvx8kvsmeU9Vndhau3ZetQIAAMBQtX0trRZ/q4NmOwZgychHAQAAYHfISJnVYBbkJTkjo7DpliSPbK19dOzaB6vqM0lemeT4JM9P8rIdqeKob0kOPGJHhgYAAGBBfePg5Jp5FwHAkjsjC5CP3vnzt+TgQ2/ZiaEBAABYUDfd4N+BMI0D5l3AJKrqxCQP707fsCZs2u/VSS7q2s+tqgN3ozYAAAAAgJ0kHwUAAAAYjkEsyEty6lj7jet1aK3tS/Km7vTIrAZUAAAAAABDdupYWz4KAAAAsMCGsmXtQ7vX65L87Sb9zhtrPyTJB3asIgAAANiL2r4k++ZdxdbaAGoE6I98FAAAAHaLjJQZDeUJeSd0r5e01jbbmPpT63xmS1V17GZHkqO3UzQAAAAwfFV1j6p6VVVdVFXXVdVXq+rjVfXCqjqsx/s8qareV1X/p6purKovVNWbq+rBU4xx16p6aVX9fVVdXVVf79ovraq79lUrsOvkowAAAMDcyEins/BPyKuqQ5Ic1Z1etlnf1trXquq6JIcnufsUt7l0m+UBAAAAe1hVnZLkrCR3Hnv7sCQndsfpVfXo1trnZrjHIUneluQxay59Z3c8pape0lp7+RbjnJjknCTfvubS/brj9Kr6ydbahdutFdh98lEAAABgnmSk0xvCE/KOGGtfO0H/67rXO+5ALQAAALCntX1tMMdOq6rvT/LWjIKma5P8WpIfTPLDSf5L1+2+Sd5TVbPkEG/IatB0bpJTkzwoyTOSfDaj/OZlVXX6JrUek+RdGQVNtyR5ZZKHdccru/e+I8m7u77AcMhHAQAAYBfNO/eUkQ4/I134J+QlOWSsffME/W/qXg+d4h5bfVv06CQXTDEeAAAAMHxnZPRNz1uSPLK19tGxax+sqs9kFOQcn+T5SV427Q2q6uQkT+lO35Xkp1pr3+zOL6iqdyb52yT3SPLKqvrz1tpV6wz1W0nu1rWf0lp729i1D1fVhRkFZ3dL8vIkT5+2VmBu5KMAAADAvJwRGenUhvCEvBvH2gdN0P/g7vWGSW/QWrtssyPJFdMUDAAAAAxbt7XBw7vTN6wJmvZ7dZKLuvZzq+rAbdzql7vXbyb5t2NBU5KktXZlkhd1p0dm9I3QtbXeLcnPdKfvWxM07R/nbUne150+tfsMMAzyUQAAAGDXyUi3bwgL8q4Za0/yaMPDu9dJtm8AAAAAWM+pY+03rtehtbYvyZu60yOzGk5NpNvC4Ye70w90i17W8/YkX+/aj13n+k8kucNmtXbO7F7v0H0GGAb5KAAAADAPp461ZaRTWPgFea21G5Nc2Z0eu1nfqjoyq4HTpTtZFwAAAOxJbd9wjp310O71uoy2Q9jIeWPth0x5jwdl9UlW523UqbV2c5KP7f/MOt8yfehYe8NxMlutwJzIRwEAAGCXzTv3lJHeztAy0oVfkNfZ/2jD46pqZZN+x6/zGQAAAIBpndC9XtJau2WTfp9a5zPT3mPtOJvdZyXJd20wztWttQ23lWyt/Z+sfot02lqB+ZKPAgAAALtNRrpNm4U3i+QjGa1kPDzJA5L8zQb9Th5rn9/j/fc/0jA33WKnBwAAgGWz5t+Cd9io315wU25M2ryr2NpNuXH89Oiq2rT/Jlsd3E5VHZLkqO5008+11r5WVddllFncfdJ7dMb7b1Xf+JOu7p7kk+uMM8mv8dIk35PpawXma2Hy0Ztv/Ppm/QAAANiD1vxbcE/no4mMNJGRzmooC/LOTvIrXftpWSdwqqoDkjy1O70qybk93v9b9zc+dumbexwWAACAAfrWJF+cdxE75YJ8cN4lbMcFE/TZPI26rSPG2pN8M29/2HTHKe4x7X2uG2uvvc/+cSatdb0xgMV2dhYkH/2HD72ux2EBAAAYoD2djyYy0o6MdAaD2LK2tfbxJB/uTp9RVSet0+0FWX2U4Gtaa9/YleIAAACAveaQsfbNE/S/qXs9dAfvc9NYe+199o+zk7UCcyQfBQAAAHaZjHQGQ3lCXpL8YkbbLBya5P1V9YqMvuV5aJInJXlm1+/iJK/u+d7/mOTErv3PSb45du3orK4yPTHJhvsQQ8ecYVrmDNMwX5iWOcO0zBmmtVfmzB2y+nSgf5xnITvkigx3C9Ojc/t/q89qfK+Hgybof3D3esMO3ufgsfba+9yY5LAJxhgfZ9pagfmTj7JXmDNMy5xhWuYM0zJnmJY5w7T2wpzZ6/loIiNdS0Y6g8EsyGutfaKqnpjkT5PcKckr1ul2cZJTWmvX9Hzvm5JcuN61NfsvXzHNfsssJ3OGaZkzTMN8YVrmDNMyZ5jWHpsze3YbhtbaLUmG+nuzE3WP5wqTbFtwePc6yXYI273P4WPttfe5JqOwaSdrBeZMPspeYc4wLXOGaZkzTMucYVrmDNPaQ3Nmz+ajiYx0HTLSGQxiy9r9WmvvSnK/JP8xo3Dp+iRXZRQGvSjJ/Vtrl8ytQAAAAGDwWms3JrmyOz12s75VdWRWA5xLp7zVeFC26X1y22/nrr3P/nG2GmN8nGlrBRaAfBQAAADYDTLS2QzmCXn7tda+mOT53QEAAACwEy5K8tAkx1XVSvcN2fUcv+Yz0/jkBuNsdp9bkqxdbPPJJA9IcueqOrq1tu62J1X17Rk9VWs7tQILQj4KAAAA7BIZ6TYN6gl5AAAAALvkI93r4RkFORs5eax9/pT3uCDJzeuMcxtVdVCSB+//TGvt5jVdPjLW3nCczFYrAAAAALBcZKTbZEEeAAAAwO2dPdZ+2nodquqAJE/tTq9Kcu40N2itXZPk/+1Of6SqNtpO4bFZ/dbmO9a5/s4k+zartXNa97qv+wwAAAAAwEbOHmvLSKdgQR4AAADAGq21jyf5cHf6jKo6aZ1uL0hyQtd+TWvtG+MXq+q0qmrd8ZINbvWq7nUlye9X1R3WjHFUkt/tTq9K8sfr1HpFkrO60x+rqset7VNVj0/yY93pmzfasgEAAAAAIJGRzsKCPAAAAID1/WKSGzIKgt5fVb9SVQ+uqkdU1euTvLLrd3GSV2/nBq21DyZ5S3f6E0k+UFU/UVUPrKqnJflYknt011/cWvvaBkP9WpJ/7tp/VlW/U1UP6Y7fSfLfumv/nOTXt1MrAAAAALB0ZKTbsLKTgwMAAAAMVWvtE1X1xCR/mtF2CK9Yp9vFSU7ptlbYrqd34z86ySO6Y9y+JC9vrb1+k1ovraofz2gbiaOTvKg7xl2R5NTW2mUz1AoAAAAALAkZ6fZUa20nxwcAAAAYtKr6zoy+CXpKkmOT3JzkkiRvS/K61tr1G3zutCRv7E5f2lp7yRb3eUqS05J8f5K7JPlyRltCvK619tEJaz2qq/XUJPfs3v58knOSnNFa+5dJxgEAAAAA2E9GOh0L8gAAAAAAAAAAAKAHB8y7AAAAAAAAAAAAANgLLMgDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcW5M2gqu5RVa+qqouq6rqq+mpVfbyqXlhVh827PnZeVX1bVT2mql5WVe+tqiurqnXHmdsY71FV9faquqyqbupe315Vj9qB8pmDqvqBqvrVbr5c2v0+X1tVF1fVmVX10CnHM2f2sKq6U1U9qapeXVXnVdUlVXV1Vd1cVV+pqg9V1S9X1V0nHM98WWJV9cqxv6NaVT18gs+YM0tgzbzY7PjQBGOZM0umqo7q/i46v6qu6H7fL6+qv6mq36uqkyYYw7wBgAGTkSIjZVoyUqYhI6VPMlI2IiNlu+SjwEaqtTbvGgapqk5JclaSO2/Q5dNJHt1a+9zuVcVuq6rN/g/0J6210yYcp5L8YZJnbtLtj5I8q/k/7WBV1XlJHjZB1zcnOb21dvMmY5kzS6CqfiTJByboemWSn2mtvW+DccyXJVdV35/kwiQrY28/orX2oQ36mzNLZIufZ8ad11p7+AZjmDNLqKoen+QPkmz2Hz3Oaa2dusHnzRsAGDgZKYmMlOnISJmWjJS+yEjZjIyU7ZCPApvxhLxt6H5ge2tGQdO1SX4tyQ8m+eEk/6Xrdt8k76mqO86lSObh0iTv3+ZnfzOrf9F+IsmTkzyoe/1E9/4zk7x8lgKZu2O618uTvCbJ4zL6fT4pyfOTfKm7/rNJztxiLHNmeVya5E1JfjHJYzOaLz+U5IlJ3pbkm0mOSvLOqrrfBmOYL0usqg7I6OeTlSRfmfBj5sxy+oMk37fJ8bRNPmvOLJmqemqSt2QUNn0lyUuT/GiSByQ5Jcm/z+g/mHxjk2HMGwAYMBkpG5CRshUZKdshI2UmMlKmICNlIvJRYCuekLcNVXVukocnuSXJw1prH11z/ZeSvLI7/Y3W2st2t0J2S1W9NMkFSS5orX25qu6Z5PPd5Ym+/VlVxyW5KKN/BFyY0Zy6Yez6YUnOS/LAjObc8a21z/b562B3VNW7MwoN/qK19s11rh+V5Pwk/6p762GttQ+v08+cWRJVdYf15sqaPqcmeUd3+vbW2v+95rr5suSq6rlJ/mOST2U0V36lu7Tutz/NmeUz9u3Pl7bWXrKNz5szS6aqTsgoEDo4yYeT/Hhr7eoN+h603hMtzBsAGD4ZKfvJSJmGjJRpyUjpg4yUrchImYZ8FJiEJ+RNqapOzChoSpI3rA2aOq/O6A/PJHluVR24G7Wx+1prv9Fae3dr7cszDPO8rD4e+znjf9F297g+yXO605Ukz53hXsxRa+0xrbW3bhQetNauTPKCsbcet8FQ5syS2Cpo6vqcnVGIkKy/3Yf5ssSq6u5Z/ebUs5NsuM3LGHOGaZkzy+e1GYVNVyZ57EZhU5Jssr2UeQMAAyYjZZyMlGnISJmWjJRZyUjZJebMcpGPAluyIG96p46137heh9bavoy+4ZUkR2Y1nILb6PaF/8nu9FOttY+t1697/9Pd6and59ibPjTWvs/ai+YMG7iuez1k/E3zhST/OckdM3oiwYe26mzOMC1zZvlU1fEZbUOXJK/r/mPZtGOYNwAwfKeOtWWkzMTPh6zjQ2NtGSmTkpGyERkpO8qcWS7yUWBSFuRN76Hd63VJ/naTfueNtR+yc+UwcPdKckzXPm+zjmPXj01yz50qiLk7aKy9b53r5gy30T0W+//qTj+15rL5ssSq6glJHpPkq0l+acKPmTNMy5xZPo8fa79tf6Oqjqyq76qqu04whnkDAMMnI6VPfj5kLRkpU5GRshEZKbvEnFku8lFgIhbkTe+E7vWS1totm/Qb/4H/hA17sezG58bafySuZU4th5PH2uvNCXOGVNVh3Q/1z09ybpI7dJdes6ar+bKkquouWZ0PL2qt/fOEHzVnltvjq+rTVXVDVV1TVZ+pqj+pqkds8hlzZvk8uHu9OslFVfXTVfX3GQXbFye5sqo+V1W/UVV33GAM8wYAhk9GSp/8fMhaMlK2JCNlKzJStklGylbko8BEVrbuwn5VdUiSo7rTyzbr21r7WlVdl+TwJHff6doYrPG5semcSnLpBp9jj6iqA5K8eOytt67TzZxZUlV1WjbYBqjzqiRnrXnPfFler0xydJK/TvKGKT5nziy3715zflx3PLWqzk5yWmvt6jV9zJnls3+efCHJa5P8u3X63CvJS5I8rqp+rLV2+Zrr5g0ADJiMlB3g50NuJSNlMzJSpiQjZTtkpGxFPgpMxBPypnPEWPvaCfpf171utPIZpplT1421zam96XlJHtS139Fau3CdPuYMa/1dkge31n6ptdbWXDNfllBVPSTJ6UluSfKsdebFZsyZ5XR9krck+fmMth67f5JHJvmtJP/S9Tk1yTlVdeCaz5ozy+dbutfjMwqbrkryrCTfluSQJCcmeW/X53uTvK37D2rjzBsAGDYZKX3z8yHjZKRsx99FRsoYGSnbICNlUvJRYCKekDedQ8baN0/Q/6bu9dAdqIW9YZo5ddNY25zaY6rq5CS/051+JcmzN+hqziyvs5PsDyAPTXKfJE9I8lNJzqqq57bW3r3mM+bLkqmqg5L8UZJK8h9ba/845RDmzHI6prV21Trvf6CqXptReHD/jLYMenaS/zTWx5xZPod3rwcn+WaSf9Na+9jY9Qur6jFJ3p3k3yT5wSSPTfLnY33MGwAYNhkpffPzIUlkpEzk7MhI2YKMlG2SkTIp+SgwEU/Im86NY+2DJuh/cPd6ww7Uwt4wzZw6eKxtTu0hVfU9Sd6R0SLpm5I8obX25Q26mzNLqrV2VWvtn7rjgtbaW1prj03y1CT3zuhbWaet+Zj5snx+NckJSf6/JC/dxufNmSW0QdC0/9qXkzwuq6HAc9Z0MWeWz/jv+dvWhE1JktbaviS/NPbWkzcZw7wBgOGRkdI3Px8iI2UiMlImJCNlajJSpiAfBSZiQd50rhlrT/I40P2royfZuoHlNM2cOnysbU7tEVV1ryTvT3JkRt+ieHJr7bxNPmLOcButtTcneVtGf6e/rqqOHLtsviyRqjo+ya90p89prV23Wf8NmDPcTmvtc0k+0J0eV1XfMXbZnFk+47/n792oU2vtfyf5Und64iZjmDcAMDwyUvrm58MlJyNlVjJS9pORslNkpIyRjwITsWXtFFprN1bVlUmOSnLsZn27H/b3/+F46U7XxmBdNtbedE4luftY25zaA7of1v8qyXckaUme3lp7xxYfM2dYzzkZbc1weEaPv/5v3fvmy3J5XkbfpPpcksOq6knr9Pnesfa/rqqju/a7unDKnGEjn0xyStc+JsnlXducWT6XJtn/Z8dlm3Xs+h6T5NvWvG/eAMCAyUjZAX4+XGIyUnokIyWRkbKzZKQk8lFgQhbkTe+iJA/NaOX7Smvtlg36Hb/mM7CeT461j9+w1+2vm1MDV1VHZfRNmnt3bz2ntfamCT5qzrCefx5rf+dY23xZLvsfW37vJH82Qf//MNa+V5LrYs6wsdrgfXNm+fzvrH6j8w5b9N1/fe2/mcwbABg+GSl98vPhkpKR0jMZKYmMlJ0lIyWRjwITsmXt9D7SvR6e5AGb9Dt5rH3+zpXDwH0+q9+eOHmzjkke1r1+KckXdqogdl5V3TnJ+5J8d/fWi1trvz/hx80Z1nPMWHv8cdXmC9MyZ9jId4+1Lx9rmzPL53+Ote+zRd/9/1HtS2veN28AYPhkpPTJz4dLSEbKDpCR0hdzho3ISEnko8CELMib3tlj7aet16GqDkjy1O70qiTn7mxJDFVrrWX0GPUkOb6qHrxev+79/avfz+k+xwBV1WFJ3pPkB7q3fqu19ruTft6cYQOPH2v/4/6G+bJcWmuntdZqsyPJS8c+8oixa1/oxjBnuJ2quneSH+1OP9dauzU8MGeW0juTfKNrP3ajTlV1cpK7dqcfHr9m3gDAnnD2WFtGykz8fLh8ZKTsEBkpMlJ2jIyUMfJRYCIW5E2ptfbxrP6B+YyqOmmdbi9IckLXfk1r7Rvr9IH9zsjqY2pfW1WHjl/szl/bnd7S9WeAquqgJO9I8kPdW69prf36NoY6I+bMUqiq06rqkC36PC/Jo7vTL2T1KQX7nRHzhemcEXNmaVTVj1fVyibX75bkz5Mc2L213tMKzog5szRaa/+S5I+70x+tqiet7VNVR+S2v8+vX2eoM2LeAMBgyUjZAWfEz4dLQUbKtGSkzMkZMWeWhoyUachHgUmVRbTTq6r7Z7TFwqEZPfb6FRl9w/PQJE9K8syu68VJHthau2YedbLzquohSY4be+uoJL/Xtc/P6l/GSZLW2pkbjPPbSV7cnX4iye8m+WxGj7l9UZL7d9d+u7X2q33Uzu6rqr/I6jclPpjkuUk2+0P45tbaxRuMZc4sgar6QpIjkvxFRiHSZzP6e+eIJN+X5KezGl7enOSU1tpfrTOO+UKSpKpekuQ3utNHtNY+tEE/c2ZJdH/OHJjRnzMfzSi0viGjn2kenuRZWf0W30eS/Ehr7aZ1xjFnlkhVfWuSC5PcI6Mw6A+TvD3J1zP6++lFWf3m5h+01v7tBuOYNwAwYDJS9pORMg0ZKdOSkdI3GSlryUiZlnwUmIQFedtUVT+e5E+T3GmDLhdn9EP/JbtXFbutqs5M8nOT9u8ehb3eOAck+S9Jnr7Jx9+Q5JmttX3T1MjiqKpp/8D9YmvtnhuMZc4sge4fgd85QdfLkjy9tfaBDcYxX0gyVdhkziyJKf6c+Yskp7fWrtpgHHNmyVTVCRltz3DcJt3+a5JnbfQ0HPMGAIZPRkoiI2U6MlKmJSOlbzJS1pKRsh3yUWArFuTNoKq+M8kvJjklybEZffPmkiRvS/K61tr1cyyPXdBX2DQ23qMz+vbwiRl96+LKJBckeX1r7b3br5RF0GfYNDamObOHVdV9kvxIkkdktM3P3TL6FtaNSb6c5O+SvDvJWyf5O8d8YdKwaay/ObPHVdXJSU5OclKSe2f0+3ynjL5pfmmSv07yJ621j044njmzRKrq8CTPTvK4JN+V5I5JvpLRU1Be31o7d8JxzBsAGDAZKTJSpiEjZVoyUvomI2UtGSnbJR8FNmNBHgAAAAAAAAAAAPTggHkXAAAAAAAAAAAAAHuBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOiBBXkAAAAAAAAAAADQAwvyAAAAAAAAAAAAoAcW5AEAAAAAAAAAAEAPLMgDAAAAAAAAAACAHliQBwAAAAAAAAAAAD2wIA8AAAAAAAAAAAB6YEEeAAAAAAAAAAAA9MCCPAAAAAAAAAAAAOjB/w+8ewg8LMcbVQAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "init()\n", "time_loop(10000)\n", "plot_ρs()" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [], "source": [ "assert np.isfinite(dh.gather_array(ρ_a.name)).all()\n", "assert np.isfinite(dh.gather_array(ρ_b.name)).all()" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.9.9" } }, "nbformat": 4, "nbformat_minor": 2 }