{ "cells": [ { "cell_type": "markdown", "id": "cbbf7b8b", "metadata": {}, "source": [ "# Scaling laws of translational static friction\n", " \n", " How does the static friction $F_s$ of a crystalline interface scale with\n", " contact size $N$?\n", " \n", " ## Physics\n", " \n", " **Commensurate contact ($\\theta = 0$):** Every particle sits at an equivalent\n", " substrate site, so all $N$ particles contribute the same maximum force when\n", " the cluster is pushed by one lattice period. The forces add **coherently**:\n", " \n", " $$F_s = N \\cdot F_{1s} \\qquad \\text{(linear scaling)}$$\n", " \n", " where $F_{1s}$ is the maximum substrate force on a single particle.\n", " \n", " **Incommensurate contact ($\\theta \\neq 0$, small rotation):** Particles sit at\n", " different substrate phases. As the cluster slides, these forces partially\n", " **cancel** because they point in different directions. For a large, truly\n", " incommensurate contact the cancellation gives:\n", " \n", " $$F_s \\sim \\sqrt{N} \\qquad \\text{(sublinear — structural superlubricity)}$$\n", " \n", " A $N^{1/4}$ regime is sometimes observed when the cluster\n", " contains only a fraction of one Moiré tile, like in a cicular cluster.\n", " The crossover is controlled by the Moiré period:\n", " the length scale over which the local commensurability repeats.\n", " \n", " Note that the exact depenend on $N$ is not monotonic but it depends on the fraction of Moiré tiles included in the cluster, hence we expect oscillations.\n", " \n", " ## Method\n", " \n", " For each cluster size $N_l$ (controlling $N$ via `cluster_from_params`), we compute\n", " the translational energy map $E(x_\\text{cm}, y_\\text{cm})$ over one substrate period\n", " and extract $F_s = \\max |F_x|$ over the grid.\n", " We do this for $\\theta = 0$ (commensurate) and $\\theta = \\theta_\\text{inc}$\n", " (incommensurate). Both circle and hexagon shapes are tested to verify shape independence.\n", " \n" ] }, { "cell_type": "code", "execution_count": 1, "id": "93e5c60b", "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "from numpy import sqrt, pi\n", "import matplotlib.pyplot as plt\n", "from time import time\n", "\n", "from flake.substrate import substrate_from_params, get_ks\n", "from flake.cluster import cluster_from_params, rotate\n", "from flake.maps import translational_map" ] }, { "cell_type": "markdown", "id": "0a89d22b", "metadata": {}, "source": [ " ## System setup\n", " \n", " Sinusoidal triangular substrate ($\\epsilon = 1$, lattice spacing $a = 1$).\n", " The cluster lattice spacing equals the substrate spacing at $\\rho = 1$\n", " (commensurate contact).\n", "\n" ] }, { "cell_type": "code", "execution_count": 2, "id": "b04d06f5", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Single-particle E at origin: -1 (expected -1)\n" ] } ], "source": [ "# Sinusoidal triangular substrate, spacing=1.\n", "ks = get_ks(1, 3, 4./3., 0.)\n", "\n", "params = {\n", " 'sub_basis': [[0, 0]],\n", " 'epsilon': 1.0,\n", " 'well_shape': 'sin',\n", " 'ks': ks,\n", " 'a1': np.array([1.0, 0.0]),\n", " 'a2': np.array([0.5, -sqrt(3.)/2.]),\n", " 'cl_basis': [[0, 0]],\n", " 'cluster_shape': 'circle',\n", " 'N1': 25, 'N2': 25, # placeholder, overwritten in loop\n", "}\n", "\n", "pen_func, en_func, en_inputs = substrate_from_params(params)\n", "\n", "# Sanity check: single-particle energy at origin should be -epsilon.\n", "pos_test = np.zeros((1, 2))\n", "e_test = pen_func(pos_test, np.zeros(2), *en_inputs)[0][0]\n", "print('Single-particle E at origin: %.4g (expected %.4g)' % (e_test, -params['epsilon']))\n", "assert abs(e_test + params['epsilon']) < 1e-6, \\\n", " 'Substrate sanity check failed -- check ks values'" ] }, { "cell_type": "markdown", "id": "696d8981", "metadata": {}, "source": [ " ## Grid for translational map\n", " \n", " For the sinusoidal substrate we use a Cartesian grid (no unit cell).\n", " One substrate period $\\lambda = 2\\pi / |\\mathbf{k}|$ is sufficient to locate\n", " $\\max |F_x|$. \n", " \n", " The incommensurate angle $\\theta_\\text{inc} = 1.5°$ is small enough to stay\n", " near the commensurate contact but large enough to activate the superlubricity\n", " regime clearly at these cluster sizes.\n", " A larger angle (e.g. $30°$ for triangular symmetry) would produce deeper\n", " superlubricity but is harder to distinguish from a random phase distribution.\n" ] }, { "cell_type": "code", "execution_count": 3, "id": "ebbe0f85", "metadata": {}, "outputs": [], "source": [ "# Substrate period along x.\n", "period = 2. * pi / np.linalg.norm(ks[0])\n", "n_grid = 50 # 50x50 is coarse but sufficient to locate max |Fx|\n", "\n", "xx = np.linspace(0., period, n_grid)\n", "yy = np.linspace(0., period, n_grid)\n", "pos_cm_grid = np.array([[x, y] for x in xx for y in yy])\n", "\n", "# Incommensurate orientation: small rotation breaks commensurability.\n", "# 1.5 deg is small enough to stay near commensurate contact but large\n", "# enough to show the superlubricity regime clearly at these cluster sizes.\n", "theta_inc = 1.5 # degrees" ] }, { "cell_type": "markdown", "id": "b434f241", "metadata": {}, "source": [ "## Sweep over cluster sizes: circle and hexagon\n", " \n", " For each shape and each size $N_l$ (half-side of bounding box), compute:\n", " - $F_s^\\text{comm}$: $\\max |F_x|$ at $\\theta = 0$ (commensurate)\n", " - $F_s^\\text{inc}$: $\\max |F_x|$ at $\\theta = \\theta_\\text{inc}$ (incommensurate)\n", " " ] }, { "cell_type": "code", "execution_count": 4, "id": "26acf7cf", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Nl range: 4 to 100 (64 sizes)\n", "circle Nl= 4 N= 19 Fs_comm=53 Fs_inc=51.8\n", "circle Nl= 14 N= 199 Fs_comm=555.1 Fs_inc=428.9\n", "circle Nl= 24 N= 583 Fs_comm=1626 Fs_inc=700.6\n", "circle Nl= 34 N= 1159 Fs_comm=3233 Fs_inc=345.6\n", "circle Nl= 47 N= 2221 Fs_comm=6196 Fs_inc=737.5\n", "circle Nl= 65 N= 4243 Fs_comm=1.184e+04 Fs_inc=556.5\n", "circle Nl= 91 N= 8281 Fs_comm=2.31e+04 Fs_inc=1216\n", "circle done: 68s\n", "hexagon Nl= 4 N= 37 Fs_comm=103.2 Fs_inc=97.24\n", "hexagon Nl= 14 N= 547 Fs_comm=1526 Fs_inc=652.9\n", "hexagon Nl= 24 N= 1657 Fs_comm=4622 Fs_inc=1249\n", "hexagon Nl= 34 N= 3367 Fs_comm=9393 Fs_inc=789\n", "hexagon Nl= 47 N= 6487 Fs_comm=1.81e+04 Fs_inc=1154\n", "hexagon Nl= 65 N=12481 Fs_comm=3.482e+04 Fs_inc=324.7\n", "hexagon Nl= 91 N=24571 Fs_comm=6.854e+04 Fs_inc=2626\n", "hexagon done: 172s\n" ] } ], "source": [ "shapes = ['circle', 'hexagon']\n", "colors = {'circle': ('tab:blue', 'tab:orange'),\n", " 'hexagon': ('tab:cyan', 'tab:red')}\n", "markers = {'circle': ('o', 's'),\n", " 'hexagon': ('^', 'D')}\n", "\n", "# Log-spaced sizes: Nl ~ 4 to ~125.\n", "Nls = np.unique(np.round(10**np.linspace(0.6, 2, 100)).astype(int))\n", "print('Nl range: %i to %i (%i sizes)' % (Nls[0], Nls[-1], len(Nls)))\n", "\n", "results = {} # results[shape] = {'Ns_comm', 'Fs_comm', 'Ns_inc', 'Fs_inc'}\n", "\n", "for shape in shapes:\n", " params['cluster_shape'] = shape\n", " Ns_comm = np.zeros(len(Nls), dtype=int)\n", " Fs_comm = np.zeros(len(Nls))\n", " Ns_inc = np.zeros(len(Nls), dtype=int)\n", " Fs_inc = np.zeros(len(Nls))\n", "\n", " t0 = time()\n", " for i, Nl in enumerate(Nls):\n", " params['N1'] = params['N2'] = int(Nl)\n", " pos = cluster_from_params(params)\n", " N = pos.shape[0]\n", "\n", " # Commensurate: no rotation, direct reference-frame positions.\n", " r = translational_map(pos, en_func, en_inputs, None,\n", " n_grid, n_grid,\n", " pos_cm_grid=pos_cm_grid, n_jobs=1)\n", " Fs_comm[i] = r['force'][:, 0].max()\n", " Ns_comm[i] = N\n", "\n", " # Incommensurate: pre-rotate pos before calling translational_map.\n", " pos_rot = rotate(pos, theta_inc)\n", " r = translational_map(pos_rot, en_func, en_inputs, None,\n", " n_grid, n_grid,\n", " pos_cm_grid=pos_cm_grid, n_jobs=1)\n", " Fs_inc[i] = r['force'][:, 0].max()\n", " Ns_inc[i] = N\n", "\n", " if i % 10 == 0:\n", " print('%s Nl=%3i N=%5i Fs_comm=%.4g Fs_inc=%.4g'\n", " % (shape, Nl, N, Fs_comm[i], Fs_inc[i]))\n", "\n", " results[shape] = {'Ns_comm': Ns_comm, 'Fs_comm': Fs_comm,\n", " 'Ns_inc': Ns_inc, 'Fs_inc': Fs_inc}\n", " print('%s done: %is' % (shape, int(time() - t0)))" ] }, { "cell_type": "markdown", "id": "0ecef451", "metadata": {}, "source": [ "## Single-particle reference\n", " \n", " $F_{1s}$ is the maximum $F_x$ on a single particle: the atomic-scale friction unit.\n", " For a commensurate cluster, $F_s = N \\cdot F_{1s}$ should hold exactly.\n" ] }, { "cell_type": "code", "execution_count": 5, "id": "0f35c1b8", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "F1s (single particle) = 2.79\n" ] } ], "source": [ "pos_single = np.zeros((1, 2))\n", "r_single = translational_map(pos_single, en_func, en_inputs, None,\n", " n_grid, n_grid,\n", " pos_cm_grid=pos_cm_grid, n_jobs=1)\n", "F1s = r_single['force'][:, 0].max()\n", "print('F1s (single particle) = %.4g' % F1s)" ] }, { "cell_type": "markdown", "id": "a2ba5a51", "metadata": {}, "source": [ "## Scaling plot: all shapes together" ] }, { "cell_type": "code", "execution_count": 6, "id": "9e90c76d", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAA3UAAAJJCAYAAAANwR5CAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjksIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvJkbTWQAAAAlwSFlzAAAXEgAAFxIBZ5/SUgAA9TxJREFUeJzsnQd4U2Ubhp803S1Q9mjL3siUoSAgCoiLpQxBkSFOfsW9ByJDUUTBgTJUBEQRQRARkaGCiuytzELZowW61/mv5ysnJmmSpm3SJM17X9dpk5Oz5/d87zJomqZBEARBEARBEARB8EkCPL0BgiAIgiAIgiAIQuERUScIgiAIgiAIguDDiKgTBEEQBEEQBEHwYUTUCYIgCIIgCIIg+DAi6gRBEARBEARBEHwYEXWCIAiCIAiCIAg+jIg6QRAEQRAEQRAEH0ZEnSAIgiAIgiAIgg8jok4QBEEQBEEQBMGHEVEnCIIgCIIgCILgw4ioEwRBEARBEARB8GFE1AmCIAiCIAiCIPgwIuoEQRAEQRAEQRB8GBF1giAIgiAIgiAIPoyIOkHwEa6//noYDAa89tprnt6UEsORI0fUMeXAz67ks88+U8utWbMm3I2+D2vXrnXJ8rgcfZlC8V537rwm3U1xXvO+zowZM9Sxeuihhwo1vy9fJ+5CnlslA2+9to8ePYqQkBDUrVsXGRkZ8EZE1JVgNE3DN998gz59+qBGjRoICwtDZGQk6tSpg+uuuw5PPPEEvvvuO1y6dCnPvFOmTFHiYdu2bfBVEhMT1T5w4GdBEIoOnwm8p/iMEASh4CQlJeHll19WDcQXX3zRJ9+/vrCNJRW2Z5YsWYJXXnkFt912G6pWrWoSQexYKSpDhw41Lc/RkJWVBX+ievXqGDZsGA4ePIgPPvgA3kigpzdAcN9N37t3b6xbt840LjAwEOHh4aq34dChQ1i/fj3effddzJ49W93E1g/suLg41ePaokULnz0GY8aMUZ+5f1FRUZ7eJEFwCw0aNFD/eX+7Ai5HX6Y1bMTxvmJH0ejRo12yPuE/goKCTMeen4WSx6RJk3Dq1Ck88sgjiImJyfO7L7x/vXEbHT23ShKLFy9W4sLdhIaGokyZMnZ/90dPjhdeeAGzZs3CG2+8oc6Bt7UrRdSVUIYMGaIEndFoVA2vBx54QFnoAgICVO/Knj17sGLFCsybN8/TmyoIQhHZt2+fS49h27ZtXb5MwTmio6Pl2JdgUlNT8f7776vPhXW9FGzjT8+tKlWqoGXLlmjVqhWuvvpq9O3b1+XrGDBggEssfyXNWnfLLbcoSyldqJ966il4EyLqSiD79+/H0qVL1Wf2Jjz33HMWv9Ni16xZMzU888wz6iUjCIIgCIJ7WbBggfIioXWrSZMmcriFAnPPPffk8a4Sio+7775bibpPPvkETz75pFdZLCWmrgRi7uPeq1evfKdnrJ0OfeR5gdKtgtC8bO1HbS+YlX7G999/P2rVqqViBcyD5fk5P39v3Y/b0cNq7969ymWlcePGKFWqlIoRpLvFwIED8e233yInJ8eUVITbocPP5vvA36332XxcQQKwrefndnTv3h2VKlVSllHzxCa7du1S32+44QZlOeWxL126tOpxe+mll3Du3Dm4A1pn+QDiNlaoUEG5dZUvX14dO/bGzZw5M888RT1n5vNfvnwZzz//vFof95nbQPfgv/76K99tp5swH6J099PdQdgj++abb6rYlPy2i7Gl7FFjHCn32dm4A2eSPhQkoPvnn3/GzTffjIoVK6pjwAYdO13S0tJsTl+Q68qZRCkrV65U94keX1uuXDnVsfO///0Pf/zxh1PXO7/rbj98Rlg/G7hN2dnZyqWM39966y2Hx4TXHafjvcxrxBWwp57Pofr16yt3LF4zsbGxuOaaa5TrjL2efD47vv76a3Vd0lrGZxjPFXvBn332WXXvmpOZmYnvv/9erat169YqriU4OFidn5tuugnz589X115BcXRNWZ+XAwcOYPjw4Wr/uL087iNHjsTx48cdrmPnzp3qvmdvP49P7dq11XVw5swZtyWbKOzx4nXPbeG5sPU759W31/ockQkTJqjfOnbsaFNg8Z6sXLmyeibSlapevXro2bOnipmxd28WhU8//VT9HzRoUJ7fCvL+teb06dN47LHH1LuO55T7xPs9P8sVr/u5c+cq6wPn4TnhseazxtY5Keg2Fuc7z9G1a/0837x5M/r376+uQ947vAeYZyAhIcHhOpKTkzF58mR07txZvcd4vHjf8fs777yjzoO7oQeWr8NnFD3IzJ9dvJb4THMGJir58MMP0aVLF9N5qFKlimrz/vjjj/mew1dffRWNGjVS1yOfQbz+f/nlF6faPrfffrt6Z9GA4qrkZC5DE0ocX3/9NZ/Cali5cmWB5p00aZJWuXJlLSAgQM1funRp9d180Dl8+LBpPXPnztUiIyPV5/DwcC0iIkKrUaOGaVp+5m+zZ8+2u+57771XTcP/tpg4caJpuziEhoZq5cqVsxiXkJCgpu3Tp49WoUIF03h+Nt8H/q7z6quvqmk6d+5sd9vWrFljWpY15vM/8cQT6rPBYNDKli2rGY1G9bv1cTDffk6rj4uOjtb27dtncxu4fE5jvjxnyMrK0rp162ZaB4cyZcpoISEhFuOsKeo50+efPHmy1qBBA/U5ODhYXVP6OnnuZs6caXPZ2dnZ2qOPPmqxjbzGeEz171zukSNH7G7XkCFDtDvuuMO0Lp4T/tf3yfwa5mdzOA3Hm1/H1jg7/wcffGA6z1FRUVpgYKBpvpYtW2oXLlwo0nWlL4vXqTXJyclav379LI5jqVKl1DWgf2/evLlT1zvvHf388ThaPxv4/DDf9nr16mk5OTl2j1+7du3UdCNHjsxzPAt6nRM+78yv66CgIHW8zffd1nLPnj2rderUyWI6zqc/0zj06tXL7jHSn5U8rubjeNx5HRfkunH0m/k6V69ebdo+rtf8mqpWrZoWHx9v8xgtWrRIHRfze4rPIn6uWrWq6botTPPA0T1T2OPFc6PfO9u3b7f4LSMjQ71v9Pnfe++9POvt2rWr+u2VV16xGD9s2LA8zxbzZdk6/kUlMTHR9L76448/XPb+XbZsmVapUiXT+9f8HuAytm3bZnN7zp8/n+e6N38ucOjZs6eWnp5eqG0s6juvoDh6T5tfm2yv6PcA99e8DdGkSRPt8uXLNpe/efNmLTY21uL9xf0xP97vvvuu5gn09Tt6XztLfm2xosLjyPeYvs1hYWGmZxmvpwULFji8B/nO53nSp+H1ZH3dPvjggzbXffr0aa1x48Y23xFczkcffeRU2+fGG29U0zzzzDOaNyGirgTCm0B/aDZt2lT7559/CrwMZy5q85cKb0g20P7++2/T7+brLapA+PDDDy1eMlu3brVotLIxN2DAAO3ixYs2t8/Ry9lVok5/KD377LPamTNn1G9paWkWooMi47PPPtPi4uJM4/jCXLVqlda2bVs1f6tWrVwq6ubMmWN6oc6YMcP0wmJjmw84NvLuvPNOt4k6Pmz5AGdnQ2Zmpvptz549pv1hY5QPeWteeukl9TsbKxRFbIDoDTmeD4oh/XhZNwT17eI54fLffvtt07XB/T9x4kSxiTo2svjiYIP16NGj6reUlBT18tAbA+adDIW5rhyJuv79+5saIFzGsWPHLBrMbOBYvwCdbRzZg4JCF98UH7bYsWOHaR2bNm1yiairU6eOmrd79+7azp07TeNTU1O1Xbt2aWPGjMlzPfOa7NChg5qP5+PNN980HWdy/Phxbfr06drzzz9vMd9ff/2lPfDAA9rPP/9s8dzhdUpxoYtfW0LDFaKO9xSfhXv37jU9R9gY0oXSPffck2e9Bw8eNAkX3jf6ceezgPvBc2re2Coojq6NohyvZs2a2Www//bbb6aGoC3hzWPCBqP1vaHPx3uC51t/tpBz585pP/30k3qG8Ny7EoovvSHJe9hV71+eM17D+vuX1zSPM0U6f+/YsaPNzj79GdyiRQtt6dKl6l1KkpKStM8//9wkFEePHl2obSzqO6+gOPPc0kXvfffdZ3oec7+nTZtmEnovv/xynvk5rd5RTGH31VdfmY4X75/du3drr732mvbll19qJUXUxcTEqI45Hi8+V6666irtscce0/79999CL/vSpUta9erV1fL5n203veNvw4YNSqyZd8RZPwN5bTZs2FD9dv3112tr16413UuJiYmqE1l/Z06ZMiXP+nv06GESkuxQ1ufl+WUbkh3P+jPS0bF84YUX1DRs93oTIupKKOz5Nu/FYAP44YcfVhcxGzuOes8L81Lh9PZ6t4oqEGjF0BsqAwcOzHfbPSXqONCiUlh4/NjLyeWw0eEqUffQQw+p+e6///4CzecqUceBL3BrKGz4wuDvt9xyi8VvPF8UBXzw2utl5suBLx3O/91339ncLg7vv/++3e0vDlGnX1u2LDYU2fo0GzduLPR1ZU/U8bjrv7FjxFmKKupI7969TfesLUaNGpWnQVcUUccOCn1eXbQ7g34O+Jz84YcfNFfxzTffqOVSaLpD1HXp0sXmNcXrXW+06J0oOiNGjDB1lJgLGR1aTMytDgXF2WujoMeLooK/3X777RbjKdI5noJb73E3Pybr1q0zdWiZiygKOV38FycUC1wvG8eufP+ykcvnqTXff/+9aRrzzhzyxRdfmOZlY9gWFP28L9jQ5f1V0G0s6juvoDjz3LL3viK6R0TdunXz/Hb33Xer38qXL28Sg85gvt7CDLY66opL1HHge1j3DtHH8XooyPvEHP3e4zLYuWvNyZMnLTqWrJ+Br7/+uumdyg5eWyxatEhNQxFu/gzUO3M4sLPbGj47+Fx15ljqzyt2Gptbsj2NxNSVUOhrzDo4ERERyid+69atatyIESPQtGlT5XtM/3FX+X+PGjVKxbe5g4ULF6p4G8Y80Jfdm4JSzWGcE2NvCguPH/3yye+//+6y7dJT7jKFtifo0KEDbrzxxjzj6cv+9NNPq8/MxHrx4kXTb/RlZ2xWjx490Lx5c5vLpU8745/ITz/9ZHOasmXLKr99T8PYEV4f1jCGQE9p/tVXX7n8umLqZXLVVVcVe6Y9fX2shWkdN8PkTF9++aX6bH5+GMtwpbPRImbQGXg96Mf45MmTBT5GjKng4CpuvfVW9Z+xxu649xgfaOua0uOoeYwZ86HDY8r4NP3cMKbSGsa8Ms7IEzg6XoybIb/++qt6LuisWbNG/Wetrnbt2qkEJFu2bMnz+7XXXqvidqyfiWfPnrVYnrs5ceKE+s+YNVfCZA3msfE6jBdkrJEeR2mOHkfNa8Fe2nrGkzL2l/FL+rF0Je565znzPLaFfu8wrislJcUiBovxl4SJ5xgH5iw8L4xVLOygn7/ihBk1p02bpmJ609PTceHCBVXPmM8PxkXyenj44YdNz5OCoL/n+vXrp2LarGHb9MEHH7Q7v37dsv1qr+RL7969Vdwm3zuMndRh3Wb9PTN48OA88/F5au/asIZxfHq+Aj5HvAURdSUUZrh8/fXXVTDqnDlzcN9996nGsf6AYEA8a9Sxsbdx40aXNNzdxYYNG0wvGAY1eyt169ZVAbf5sWzZMpWkgIHZFN3mAeZM1EDi4+Ndtl1sqHLZTFLAlzyD3/XGRXHAAPn8fmOwvnljjMlR9OQefMjbG1hjkehB+9a0adPGIy9F63vRVpIG/SWiJ0LZtGlTka4rR/cOG73FTbdu3VQDgI2CL774Ik9HDRvgbNTZShhRGNh40jsP2BnAwrxMxMMGiD34Qv77779Nwe8FhZ1NrDnGhinPEa81/V42rxnoyvtZhyLGFtWqVTN9ZmNMh7VJecyJ3pC2haOEUUWlsMeL0zM5BDt+9EYak5gwwQ+vISZO0oXf6tWrTfPpn/XfdHidMJkIOzt5b7KhePjwYbgbvfFnS1C741rgs0cXkObXAoXsn3/+qT6z88TRM/aff/5x+Ix1huJ+5zmCx57P1PzuHfOEKXw2M9FPYZ4T3G92UhR2aN++PYqbRx99VCWkY1ItPSkL70+WTeAzVU9Cx86EgiSD4rNY71xwpl1gDduz+nVIA4W9a7Zq1aqmJGrm163exujUqZNd4wDbsrxv8sP8HhZRJxQb7IFj5kBm3GJWTL4UmYVPfzCxJ+OOO+4ocpavwjY6nUHvteUDxpvJ7xhQuLABy2PPFxkbEXzI0Zqk98qxoaH3DLoKZn1kpkg2oGgR4zYwux97G2kpckcPrDlclzO/saNBRxedPA60Jtsb9ONk3qtaXNels7BHz9xKYO8YmO+/q/bBk/cOX5rMdGie8U+HmVgJr0VXWviZ5ZSdV3zJjh07VmW8pAWP9wDFhHnDlpw/f97UWCvoMfr3339VFl6WhaEFietkzzEb0fr9rOPK+1mH+2UL8waJvm/WDQ/zxmtB7teiUJTjxfcYsyWaCzV2WLDDgKKM+6w3BPXfaanUhYu1qGNnA68VXnsUhuz0pODgvcZGONOVFyZzaX7o71lHzwNXXgvm14P5tcD7gMdOFy+OnrH6fPaesd74zivqsbI+XuaWY29vh7gbZpCml4AumNgx4iy87tiRlt9zRvdesca8M5ptV0fXbc6VTOjm163+DHT0/OO9qVvhHGFuGXdHltzCIpY6P4MP0K5duyqrzb333mvqIWNj31tT7Hqru2VBjwF7g2kl43S0ItA1Sndt0Hvl7rzzTjWtqxsUdHPkC5XWWbomsPHC8043RzaG6Aph/hLzNLpLFN0OdXc8R4O9tMIlIfVzUfbB0/cO0+3zJcm06mzIE37WXa100efKwrDsjeXzjL3NtO7z5U7LL8UEe+jNLTlFOT7sEOE9RFceuvVQILJhSnHOe9m8rIA7BEJR8MR1UdTjZS3a9P/6eLpY8v3Ga4vPMp5zCghaGGxZsuh+xUbpxx9/rIQcO7nY6NPLWtA6SJczVzeISX5p892Nucsp078784wtqDu0p995rsTTz1Fvg/eauQeAJ65blrdy5rodaqPckivOp3kHoX5fewMi6vwY8waV7mLhLvQeMEc9GuYxVebQnF5U9w93bVthfMnZKzxmzBjVwLSOiXFn3Bt7p0aPHq1inNiTtWPHDrUtujvcRx995Jbj4qhmlvlv5hYpd59zZ3DV/rNH0ZELoH4M3GFV9PRxZI8nPQHMrXX6fwouDq6G9xRrl7333nvKbYovX9bhouBjY5qWA/180IVGj8soyDE6duyYybWVjVY2TK1d6jwVw2oP8zguR+7X+dW4KwyuOF66tU0Xa9aijp0HdFWjUKSLmP47rbT2Ym+4DYzp5LP56NGjKpaKMVNs9P3222+FEjKOsOUK6QnYCNWfb+58Nnj6nefq52hhjhdj8Ry5t+Y36PdNSYD3m95J6Wy7wFXnwfz+c/T8Y6eDM7UTze9hV8fIFgURdX6MuduTtTuI/uB1Ve8Z3S30l7st2JtuL6ZI9ynn7wVJgGD+8nC0H/ltG3GmSHZ+6MvX3YisoQ+4K9bjLEyYwwa2Hg9Jt1xXnTNzHLl36r/xXJkfF32bVq1a5THXBn3/aUnQXZWsceZ80d2EDURb8Lpct26d+syCzK5Gv3eWLl3qsmUW9NmgJ0xhxwEbcHp8nautdI7crSjk9AB7dmjocR1s2DIeq6DHyPyesHc/89r1JuheqCcIcVQw1x3FdF1xvHRxRncqTstYSDYSW7RoYdOapz9brF0vHUG3TBYr1+M8rZ+JRYXupyS/+D1Xv3+t4XEszHVf0G30tndeYeGzWY/NLujxohuwIzfB/AZHHYKeQndrJnp8nTPwGDZr1izfdoG5N4U5tPLrbpuFuW5btWql/uvvXFuw00h3EXWEfg8zfk9/rnoDIupKILzYGL+QH59//nmei12HmYOIHlhfVPQMhrQS2XoJcFvsBUrTNZDbwxvt8ccfd/pFp+9Dfvuhbxt7b2y9YNiot44JKgx6hrHt27fb/J0xQEwk4GrsCRJr33DrHtSinDNz6A5lq6FIsfbOO++oz7SsmD8Y6bbHBjd7zF599VWHy+dLTw+KdiX6/nPfeQxsvazpzuoM48aNM/n4Wx9DveFDFzBXw2Bysnv37jyW2MJS0GcDG+NMyMTzzX3kOXVlghSd/Bo/5jEQ5te6foyWL1+uBmcwzxZo637mffzGG2/Am6D1iYkOCF0ObbkA0j1OT1zhSlxxvHjNMPERYRIwvg/oIml+LnUBx/ACvcPJlqgr7DOxqDBBg/5OYWZBe7j6/WsLvVPFmevelmXRmW301DvP1dCFd+DAgerzxIkTHXYAW0P3P2fcBO0N7kxcZIv82le8FsaPH68+02XZnmC3h/6eowu2LQ8x3ht8Ptlj5MiR6j876fKL57tgdd3qrr689+bNm2dz3/V9yw+9rajf016Dp2sqCK6HRURZVJW1v1hA1LzOB+t6bNmyRRs6dKipFgcLgFrXOxo8eLD6rX379qpOXFHqwFnXy2LRTxZ4JSxCy2KRrFlSrlw5uzVkPv74Y9P8LDBrXXycRV1ZiNe8qC2Jjo5W8/zvf//LU7NJh/uu19xp0KCBKuDKWngczxoxjRo1Mm2brVvGmTp35sW0WdeExYz12iasy6LXYWINHHvHoLB16lhsc9iwYdry5cu1hIQE03jWqRo7dqypUD23yZXnzLz4OKdjXRf9HLBg8g033GCqg2NetN66BpVeSNm8mDSXw2uA07AQrHWNI0f18wpyDV933XXqNxbxZTFfFuzV6zfx3jC/LvIrPs7CpnqtKBbD5vFm/Sz9mi7sdUUc1TRinTj+xmfCc889l6f4+KeffqoNHz7c6XpP+/fvN/3GYtfOMHXqVIvaS/ZqJhalTh23uWnTpuraZP0j/ZnGe3n9+vXqNy6XtQ3186hfS/p55vl466231HHRYQFqLvOZZ54xjeOy9QK6LJarF/HWC+iy9p5+L9s6L66oU+cIe+vludOLcbdu3Vq9C/Rj9Msvv2i1atVyS/Hxoh4vnRdffNHiOuJ1ZQ7PpV54mAPrm9p67vN51q9fP23hwoUW9ddYN+2jjz5SzzbOb11w3hXoBb1ZvNoernz/2qsnx3uga9euppphfBeYF1tnkefVq1er+rZ8hhdmG4v6zisoRa2v6eiY8rlpXnyczz69NiDvH76fnnrqKVX/rzjgM8p8ML8nzMfrBdLt1aGzhtvfp0+fPPcG95X1YOvXr2+a19E1bA+2H/T6sjVr1lTtDL328J9//qme046Kj/Me1Z/lvC65v3rbhLCNw7YO2wyNGzfWrOnWrZvpvcxrQq9fyfN71113OV18nG1FTlPYen3uQkRdCWTFihUWLz79oc0GqN6A1we+TM0f5OZFW/Vp2ehmo5YPQ/MHYkFEHeFNZr5u3rhsaOqiK7+G+Pjx403Tc2DjhPtkPs5ctBC+qPTfWFSXD2PuAxvY1seMDW99Wt7UeoObBbLnz59fZFHHbWOhV3053G4eA/04P/DAAw6PQWFFnT6fPpQuXVoN5uPuvPNOm4WMi3LO9MYEG8T6A5DngA9ifXnc908++cTmdvNBz2K95tcszzkbAeaFUDn8/vvvbhF1FI7mx4rXREREhPrMorksVm1vfvNGxLRp00z7wUaz+bXWvHlzi5eSq0UdX+p9+/bNcw2YnwduQ0HEw4033mjRcNafDe+++67dF7l+3DiYN+pdKerM95HHmNcKG5Tm+/3rr7/mmZeNn44dO1pcl7zWzQWCtfBm55n5svnM0BsD3FfzThFvEXWEnSvm283zp283O8H065b3akFx1HAuyvHSofA0P8e7d++22Yml/86OTVuYN2o58DybNyQ5UOhT2Liaxx57TC1/0KBBdqdx5fvXUZFw3pe33XZbnmeD+XuJA89bYbaxqO88bxJ1ZPPmzaaOYn2/+YzR2woc7D0DXY35OXM02HqOOhJ11sXSeW9av3P5bPjggw8Kve3sxDW/3/gc0J+1fB5RMDs6D2yzXnPNNXme16Wt2jW2isizQ8H8muR7Qt8WXp9sj+gdUGz32eKff/4xtavPnDmjeRPiflkCoSsb3WiYKEAv8MiYObpJ0I2gXr16qsAsg5gZl2ArvStNyj/88IPKlEm3OPp2MzC1KEHVzLTIbWIMBN1b6I7G2Cm6+7z//vv5zv/8888rNw6a3/U6M3S54v7cddddWLRokYXLJWHqXa6TPvGMI6C7IPfBOjibx4xxT6znxVgqZlmiawGD5lkXyTxAt7DwODLomclK6BvOgGG6GNK9gskDHLkcFIWpU6eqkgasV8djxfcBXQd53nv27KkKiNIVwparUVHPGeHxZC1EHksmq6DrE2NhmOaa/uu6O4UtdzG6WTGhCwud8jrWa1VxmYwXY1ZPHlN31UnkftPNgq43TGTC/WfyD9bwYYkQPUYmPzg9C6SzfhqPM4eGDRuq/WNKdXdmz+I9z3PMWlF9+vRR552ukLz2GN/ALJF6iQFnYXwcXaHr16+vMg3qzwZ7rli8L7t37+7WBCl0zeN1yRg+Lp/nidkLmRGR55HZL5kxzVbNQE5LF2EWRGctRwa+M+EGjx2XxWvX2i2Hzwpm9GTRbN7bdAfkcpjlkc8MvWaet0EXJLom8t3A/eT9yNTyjz32mHJn0l3mXB0n4orjxXtej//mM9nW/Wde48pePN3LL7+snl+8H3gf8l6gCzfvcdZXZEF6Xg+sqeZqmJiFsGyCvTT+7nj/2rsvGZtE90u6xenPZ8YtMnaJ9yxjDG25yTmzjZ5657kLhqnwGUIXTL1cCt1HeR9xnyZPnuxyt/LihvcMwwV4vzLGlO0mvnN5rfAZy4zUPAZ8JxcWtsf0RG28zvgs4HOHGdmZvViP97QH32EM6+A1xDYM49p4zWZkZKjrjG2LKVOmmDIum8PnBtu9fAY0aNBAvYt5TbJ9xFg+tkf0BGj2noFMvEX4/PCmJCnEQGXn6Y0QBKHkwYcrX/AsEG4rrbDgP7ChyJc309hPnz692JKkCAXnxRdfVAKW4uiXX36RQ+gGeGyZKIIxtUOGDJFjLAhewv79+1VnJWFGXHbum0PJxM7xgwcPqoQr3hZTJ5Y6QRAEwa2wR5WCjr29vt6TXZJhnTYW5Sa0KgvugZYQ8tZbb3l1fTZB8DcmTJig/tMLwFrQEXqDUNDRu8vbBB0RUScIgiC4Db4A6epCHnzwQYtSKkLxQ7dDuo+xJpueupuWVLrgsZHC7HN0KWIGWsF9xZvpBsustHR9FwSheNi3b59y+6RrpnnmVY6nKzg9iwhd7q1h+AVDJuiyOWnSJK88ZbnVJwVBEATBhbCUAcurMH6VL8OYmBgVFyt4lkOHDqk4WZ4LxjgxloWxh7rA43f2RrszzlOAahQ2adJExaQKglA8pKWlqXIIet1SPu94DzImT4dx5vfcc0+eeVn2irHIrM3HOr/eiMTUCYLgFiSmzr/Rzz/FAS1AdDXTExwJnoPJUJgQhj3Vx48fV26xTD7ChgpdipgwRS/wKwjFzYIFC9Q1WBCY5IUdFYKQH7TOMTHYqlWrVAIgeiawQ4tJkmhBZ7y3tya5cgYRdYIgCIIgCILHYcZlusEVBGZN5HyC4O+IqBMEQRAEQRAEQfBhJFGKIAiCIAiCIAiCDyOiThAEQRAEQRAEwYcRUScIgiAIgiAIguDDiKgTBEEQBEEQBEHwYUTUCYIgCIIgCIIg+DBSfNxLqFKlCpKTk1G9enVPb4ogCIIgCIIgCMXI0aNHERERgVOnThVqfrHUFQDWQTEYDHmGtWvXoqhQ0LGqvSAIgiAIgiAI/kVmZqbSA4VFLHWF4Pfff4fRaDR9b9y4MYqKbqHbvXt3kZclCIIgCIIgCILv0KRJkyLNL6KuELRr1w6BgXLoBEEQBEEQBEHwPOJ+KQiCIAiCIAiC4MOUGFG3efNmTJw4EX379kVMTIwp3i0/UlNT8corr6B+/foIDQ1FtWrVMHz4cBw/ftzuPNHR0cpS16xZMyxcuNDFeyIIgiAIgiAIguA8Bk3TNJQAevfujSVLluQZ72j30tLS0KVLF/z555+oWrUqOnbsiCNHjmDjxo2oWLGiGl+7dm3T9D/99JMSj3S/pBicOXMmFi9erIZevXq5xI82v5g67k8JOWWCIAiC4Nc42wEtCELJp4mTWqDEi7o333xTZYxp06aNGmrWrIn09HSHAuill17CuHHjcO2112LlypWIjIxU4ydPnownn3wSnTt3zjezJYVgRkYG/vrrL7eeSArQEydOqHWVkFMmCIIgCH4NBV1wcLDyEqK3kCAI/ksTEXW24cPRkaijOKpUqRIuXryILVu2oGXLlha/N2/eHDt27MCmTZtw9dVX2z0Bb7/9Nl588UW1LnedSAo61q7Izs4u0joEQRAEQfA+mFGbWbBF2AmC/9KkiKLOb1M4rl+/Xgm6OnXq5BF05M4771SibunSpQ5FXXFACx0FHR/2ejyfIAiCIAi+TVZWlorh171xzEM+BEEQCoLfqoPt27er/61atbL5uz6ews4etAJ+9913NkWhq+A6aFUkFHR00xAEQRAEwffhO53v9oMHD5rCKyTGThCEwuC3oo7ujISZMm2hj4+Li7Ow3rVt21ZlvaS75YwZM/DHH3/g+++/L3JhQT7QaTV0lBhFLHSCIAiCULLQ3+36+15EnSAIhcFvRV1SUpL6Hx4ebvP3iIgI9f/y5cumcSx7QCEXHx+vvtNCt2zZMtxyyy3Fss2CIAiCIAiCIAjW+K2oKwzjx49XQ1GwF/xoz4InCIIgCIIgCILgF8XHC4peviAlJcXm7yyPQEqVKlWs2yUIgiAIgiAIglAQ/FbUMXUw0V0prdHH16hRo1i3SxAEQcif4cOHmwo3mw+sPVoc6x09enSe38qUKYNJkyY5zHRYt25d/Pjjj/BGFi9ejFmzZjk17QsvvIAqVaqoYzFlyhS3b5uv0r9/f7z88sue3gxBEPwAvxV1rENHWKPOFvp4JkVxF4mJiThy5IgaMjMzkZOTA39m6NChuO6667xmOSUZR43LlStXok2bNggLC1P3yV9//eWWbdi2bRs6duyo1lOrVi1MmzbN4ndpDAmO2LlzJ3r27KmSVZkPDzzwgNvXGxISkufeYVKtS5cu4aqrrrI779y5c1Vpmptvvhm+LOo2bNiACRMmqBqtPOYDBw4slu3zRZ5++mm8//77SEhI8PSmCIJQwvFbUdehQwfVq8qsk2xcWrNw4UL1//bbb3fbNrB3k41ZDvv378f58+fhz7A385NPPvH0ZvgF9hqXzOTat29fJah++OEHVKhQAf369VMi0JWcPXsW3bp1Q+nSpVWyoYcfflhZPubMmWOaRhpDgj2YIXDPnj3qOX7NNddYDLxm3b3eu+++G//++y8OHTpk+m3Xrl3qvyNR9+GHH2LIkCHwdbjv5JFHHlHHnBa7wsDOTNZgLcmwg6xatWrqmSsIgneRlp2DaXGn1ZBeAgwrAf5cG2bUqFGmF5MeQ0cmT56s6tN17tzZrYXH2Yg9fPiwGurVq4fy5cvD3Rw+l4yJP+7DqHlb1H9+9xZY0qFx48Z+3QAoLmw1Lmk5pmsZr38KqhtuuEFNd+zYMZdb6z7++GPltvXNN9/gxhtvVOu7//77MXbsWNM00hgS7EExxXjoRo0aeWS9gwYNUrXFzK11tOBFRUUhNjbW5rzsuNu4caPqNLHVichsyuxoYTkd3gd6KRt29lE0mbt7sp4ZvUjMOx1///13lYm5UqVKqsOyU6dO2LRpU551zZ8/Hy1atFDr4nLvuusuVaKHHg6ff/451q9fb3Jl/eyzz/LMz+mGDRumPhuNRjUdvU3IO++8o4pn05LJZ/m8efNselFwPLNJcxv4fLGHvW3VcXZ9PL70TGAsPZ97PH7r1q1Ty2bcfO/evdXzr6jz2YPn/Msvv8x3OkEQipd5J8/jjUMn1TDv5AXfP/xaCWHZsmVau3btTIPBYOAb0WIcpzEnNTVVjed0VatW1fr372/6XrFiRe3gwYPFtv2NGzdWgzXZ2dnanj171MDPRWHB30e12s//oNV4dplp4HeOLy7mzZunNW/eXAsJCdEqV66sDRw4UEtLS1O/3XvvvVqHDh1M0+rf586dq9WrV08LCAjQDh8+XODlkG+++UZr0aKFmj46Olp7/fXXtZycnEJvK3n77be1WrVqacHBwVqjRo3UdpqjbwfXXadOHS0iIkK75557tPT0dG3t2rVq2ZGRkVqvXr20hIQEl82bH//++6+6xvfv328xfvz48Wp/srKyTOOSkpLUtNb7VlS4b9xHc1avXq3WZX7fvfDCC+qeFARzvvvuO3Wt7Nu3T8vMzDQN+d3Trlrv6dOnteHDh2u33nqr6bfBgwfnee6Y8/HHH2vly5e3+Zzhs2306NHaTz/9pL377rtaeHi4er7ofPvtt2qa3377TX1/6aWXtKioKO348eOmaXiPvvPOO2oZHIYOHaqWc/Tof8/3mTNnqu3nti9fvlwtl8+VS5cuaQcOHNBuueUWrVmzZtoff/yhhjNnzuTZXk73/PPPq+Xo0/G5OHXqVPXeffHFF7UVK1Zo999/v5qG69HhPV+hQgX1rluwYIH2ww8/aImJiTaPl6NtJc6ur1KlSuoZsmTJEm369Onqef7AAw+o9wG3Yf78+Vq5cuW0//3vf0Wezx5Lly7VjEajadvd+Z4XBME5UrOytRbrd2mVV29VAz+nefj+s6cFnKXEiLrZs2erB7qjgdNYk5KSor388suq4czGeZUqVdTL8NixY8W6/e4WdYfOJuURdObCjr+7m/xe0rZEna0GQEGX40yjyZsaFEWdNz/sNS6bNm2qPfvssxaN5BMnTqj9+vrrr/NMz+vRfFpbg71rlp0mEyZMsBinr8v8GObXGBL8k7Fjx9p8xm/fvr1Ay3n66ae1smXLqmt/69at2q5du9T91aZNG7vr5bVL2OHC54je0cOOFt6f9hg5cqTWqVMni3EUobGxsdojjzxiMX7SpEnq2ZeRkWEax05Hdm6tX79eCwwMtPk+s74369atq5alj+P77e6777Y7n60OMVt8+umn6njrsCOIy7bejx49emht27a1WD6fxRTjjshvWwuyvqCgIC0+Pt7iOHLbN23aZBr3zDPPqPNQ1PnsceTIETXvunXr7O6viDpBcK+Am3rklBp04Tbz2BmToNOHWfFnPXoaRNSVENwt6iYs32tT0OnDxB/3au6kMA0KWw2Agi6nII0mb2lQFHXe/LDVuDx37pzDDpENGzbkWc6rr76ab0eKtTVOh43Sjz76KI/l3NoqmF9jSHDPc8h60K26/G9vGh12ctj6/fbbbzdNM23atCK9vAYMGKC1bNlS+/vvvy0Gcytzfvz888/akCFDlChjhxEFEwUdOxIcrbdLly7qMzuYeB1zfymg2OnC/bIH9793794W4/hs4/W9Zs0ai84Q3uMcT6uYDq1mFJRcD5811vD3ESNGaNWqVTN5qnDQhSbfIfy+cuVKl4s6/T6ltd26c4zPMf28cPkUmvmR37YWZH1XXXWVxTS0MtLLwRx2mvFc6pbews5nj8uXL9vtHCMi6gTBvcw0E3AUbhR25lY6fWjpYWtdUUWdFB/3IPTF1/3xGS/G+AR3EZ+Qks/vqXAn//zzD06dOlXgJAGMl2jQoEGhl8OAfsZs3HnnnRbJPrp06aLiuI4ePapi+QqyrSx3wd/vuOMOi/FMKPLggw+quD/9XHLbGXujw3UxPsM8VpPjTp48qWJoGJ+iU5R5HcFtL1eunMW47du3q/+rV6+2qM3IGBVmbrOV/IExcLfddpvDdRU1aYUeZ3r69OkiLUcoWTB+jUlSWrduXaRlDBgwQMVjMRbtp59+Ugl8HF3TnIcxoIRxa0wSwrg6lshhrJejJClpaWl57rtz586Znke24LNLfz5VrFhRTff111/bzPB57733Yvfu3RgzZoyKVwsPD1exb1wv0RNxVa1aFa6GzyDCeD5zKleurN5t3E9+tjWNLfLb1oKsj+fJOp7e1ji+H/jsDgwMLNJ8tuA1RvRzIQhC8SZDmXr0jOn7+3GnkZWTg5PpmXmmPZGeifknL2BotPsSbrkTEXUehNkv+QLW4UvbXcSUDc/n9zC4k8I2KKxf2gVdTkEaTd7UoCjqvI6w1bjkcaAQtT5OTFvOhrO50NNh4oL8GmgBAbZzMZUtWxYXL160GKd3cPA3HWkMFS8UBY5g0oj8prnpppvynYbJqTgUBt5jTDriqHQBOyOY8IeJLdjhwUQ/FDnmMNEIE4Fcf/31qlNj7dq16j5icg6m6LfuJNHX+9hjj5nGMXssk4tQYBJHoo73nHVSDf0+5DJsJYky79D6+eefVWIhbvfzzz+PHj16qAQiJDU1VYnS2bNnW3RGmd9jegcJn1+OtrMw6M/KM2fOoEmTJqbx7IwJCgqy6NxxpvMpv20tyPq8Af28Wz93BUFwvYCbEX9WfR4ZWxEhAQEqGYq5gOPnSUfsdxRT9N1VtZya19fwvS0uQRRn9ssBbWJhDLD9MuX4/q1tZ2xzFeYv6YJg3QAo6HLMG01///13nsFWdtP81mHeoDDHWxsUzjQuKQytjzUbhCtWrMDgwYNtLuf1119X++toYDZNW9CSsG/fPotx+nfzhqw0hgRreJ1QYNkTJsxO+dJLL+G3335TWYx//fXXPIKO0OLGcjL0BqDV+dNPP8V3332nOttsWQBtrZfCit4AS5YsUc8FR89wPuNZy84cXutMd89OFa7TetA7U5KSkjBy5EhleaN447PmtddeMy2HVkLWOWUnj86ff/5pkVmS6+I2OsrCyPnNs0s6CzN2spNn0aJFFuOZPbJVq1YF9kLJb1tdvT53Q48Q/RoQBKH4slmm51ha6XQuZdnPpK5b63wRsdR5EKa/5kDYAHYntSpEYELfpnh+0U5k5+SmytYF3cS+TdXv7sT8Jc36ZMW1HPNGk7Mum/mtw7xBYW7Z8tYGhTVsWOh1GHWYhp3CjnUbdcvlxIkTldVsxIgRNpdTFPdLWnNYbJwWBhYfJ9wmbhsb2TrSGBJsuUASe6KO1mF2ULDT7J577lGp6e3BjgkO5tgr32FrvSxDQKv8ggUL7HoD6Fx77bUYP368Kp8TERFh2tZJkyYpscbi1F27dlXjKBRpmaNYJM8884y6P2l9pLX+vffeU/PQ5ZsdU3yPcFvo+cFls0PmlVdesfA24HK5fs5H8cY0+xRwXAddrCkg+ezjc4/j6PpN0etMZyOfeSxETismnxk85hTIdE21LtLuDM5sqyvX5yq++OIL1ZHF52iNGjVM47ds2aKOIzuzBEHwvJsleb1uNQyuZvv5ZoRz4Szehog6P4LWuDY1y+HrTcdUDB1dLjnO3YLO2Ze0O5bjbKOpJDQoHDUq8mtcduzYUblSPvTQQ8q165dfflEWCxYgt2XlIBTLHAoDYw95LFnknI3vrVu3Yvr06Zg1a5bFdNIYEqxhkW+6qttz/aVLIqdZvnw5XnjhBbRv3151UBQVLpOdH6VLlzaNo3iktY6eAPm5NLLuI58dq1atQq9evUzjWfOOy+Q9ybqQfOawc4P3BqFbKOs6fv/99yb3a4pVCkk+ozZv3qw6BVncmtY8xgnWrFlT3b9vvfWWxTaw/hqXT7dqPiv4/KLFUrfwsQNnw4YNKj6PwpDunJzHGVj3la7dH3zwgTre7ByiQOTxKQz5baur1+cKaC2lC69eY1CHHg+saScIgvuYV0A3y4+PncW90RV80s3SLi5N2yJ4dZ06b4CZDZlVjOUjWPtt0KBB+dapK+pyCDPaXXvttVpYWJhWpkwZrXXr1tpbb71V6G3VM2jWrFlTZVtr2LCh9uWXX1rMb2s7mDGSdfJsleNg1ruizqt/1+v52YJZJnkMFi9ebDGeGSZZb4/Z9a655hrt119/1dwJU8hzH7m+6tWra++//36eafr27asy+gmCs/zzzz+mzyw98uijj3rNwWPpE9azE/wHZr5k6QtHz9OS9p4XhOIuUWAvm2XlfIbZHi5h4Orslwb+8bSwFGAK9rZOMMCeP2ZjJHSNsZd4QhAKApNM0FLnKL7G0zCOiK5t7OWmJVEQnIHWK1qbmCmW1uoZM2Z4TYKKAwcOoEWLFsqSridTEko2dJWlazljPO0h73lBKBiz4s/ihf3H1ecJ9WPU/+f/jS+Um2WoMcDrtYCziPuln5Q0EARzWM6BjUsmXPDWxuXMmTNVjKIIOqEg0GXQW2H2UMaSsiyKt953gmth5wKFnSAI7omde+/IKTrDw6/cLO0gos5PShoIgq81LqUxJJREnI1RE0oG9hJNCYJQ8PIEtmLnTmX8V4PYFr5ee64giPull1jqunfvrix1e/futZhG3DIEQRAEoeQi73lByN/Fclh0BVWi4Jo/9+bJaFklOBBr2jZEsL3SXV7mZuku90vv38MSDNNQM0sZB2Yvk3g5QRAEQRAEwZ+xVZ6Ago6152yVKKC1bsmZREQYjTYHXxB0rsA/9lIQBEEQBEEQBK8TcNPiTquBws1eeYIvjp9X4s4e75vN769ITJ0gCIIgCIIgCMUOBdwbh06qzxGBRgyqWs7CSqfz9pGTuJhlX7Sd8KPYOXuIqBMEQRAEQRAEweNullk5OTZdLCnoHJUn0GPn/BkRdR5EShoIgiAIgiAI/ogtN8tJR+y7WPpTeYLCIEfFwyUNatWqpYb9+/fj/PnzntwcQRAEQRAEQXB73Bz/23KzvJSVna+LpWAbsdR5kNGjR5tqFuklDQRBEARBEAShpMbNsTyBvUyWxJGbpb+7WDpCRJ2HSxpwICxpIAiCIAiCIAglpWi4piFP3NwdlaMcZrIUN8vCIaJOEARBEARBEASXW+U0TcsTN/fUP/F2rXREMlkWDhF1giAIguAG7rrrLvzyyy8oV64c9u3bZ/Hbs88+i/bt26NRo0a47777cO7cOQQEBOD+++/Ho48+KudDEASfz2b53pFTgA13yb8Tk7D3uqsQHGDflVLcLAuOJEoRvAbGF1533XVesRxXbUtJIisrC3Xr1sWPP/6Y57eVK1eiTZs2CAsLQ/PmzfHXX3+5fP2bNm3CkCFD1DYYDAa89NJLTs332Wefqemth7Vr15qmycnJUY3pyMhI3HjjjThzxjJ4u3///nj55Zddvk9CyeaBBx7AihUrbP62atUqdO3aFSEhIfjwww+xZ88e/PHHH5g2bRp27txZ7NsqCILg6myWpzKycCojr0XuZEYWlpxJRITRaHcINYpEKShyxASvgY3mTz75BN6AN22LtzB37lyEhobi5ptvthj//fffo2/fvkr4/PDDD6hQoQL69eunRKArWb9+Pf78808ltsuUKVPg+X///XfVaNaHVq1amX778ssvVaN66dKlSjS++OKLFvM+/fTTeP/995GQkOCSfRGKzvDhw22K9XHjxhXLepnoyhpel5MmTTJ9v/7665WVzpqjR4+ifPnyiIiIQI0aNXDVVVep8aVKlULDhg0RHx/v1n0QBEEormyW9njfbF7BNYj7pT/AKNXMFMfTBIUDBs9mFKpTp47D3zMzM5V7UnFkCc1vW/wRWhNoKbOutchG7uTJk5Wli0RHR6uGKa11HTp0cNn6//e//+Gxxx5Tn2vWrFng+du1a4fAQNuPPG7rc889hy5duiiL4zXXXGPxO8dVq1ZNCdtRo0YVcg8EV0JrVs+ePfH8889bjKcod/d6aV2jxZplaXTi4uJw6dIlk0BzxPLly/N0jpBDhw5h8+bNuPbaa12+3YIgCJ7KZmkLiZtzPWKp8yBsEB85ckQNFCx0AXMLFHTjqzke8hN9LmL+/Plo0aKFsvhUqVJFxZykp6fbdHnUv8+bNw/169dX8xw7dizf5dhi4cKFaNmypZo+JiYGY8eOVcG79rC3LUuWLEGDBg1Uj/ott9yCkydPFmgfyTvvvIPatWurhmHjxo3V/tlaN7eZDVS6BFJMZWRkYN26dWrZXH/v3r3VNVTU+ZyBdRQ3btyoLHLmfPTRRyhdujRGjBhhGsfjqzdyXQkFvbugSPzqq6+UJW7WrFmqdqQ13Hda9ATPw3uXllV2GlCAmw+0FLt7vXfffTf+/fdfJcJ0du3apf47K+r4/DCHgpDXGIWinhVZEATB2+PmaHG7lJXlMJtlleBAFUN3sFNTi2FglbyeDELhEVHnQfyt+Dgby4MGDcLVV1+N7777Tll+WMqBosMe//zzj3KneuONN5RrXNmyZQu8HIqsAQMGKFcougo+9dRTmDhxorIuFYQDBw4oMThhwgTMmDFD9ahbW23y2zbGy9CVj9NwWzp27IjBgwfniVPj9fD222+rbeTw9ddfq+QJdPl64YUX8Omnn+K3337DK6+84pL58mP16tXKXczaCsJjS7dLNnbpbsmBjVN7ZTrYcaFPZ29wV+cGLYi01DVr1kwJX+vYJzbW6So3fvx4NVhD6wnj+i5fvuyW7ROch2IqJSVFJRnxxHp5//J6Mr9vacGjGIuNjXW4DHbwHD58WHUOmY/r1auX6pih67IgCIK3ulpax805k83ylJ0YOombcy3ifulB/Kn4OBvqjFNiD/fMmTNN460tP9ZcuHBBxULpDaCCLodig1nmHnroIbz77rumY03xQGFHweNsjUBuy99//21qtDHuhS573CZakfLbtuzsbCVQH374YSVSyU033aTia1577TULdyxajLZs2aIajoQZ9KZPn65EBQUj2bp1qxJVjPUq6nz5QQHbpEkTi3HshGBDlsObb76ZZx7dYmfO66+/jjFjxjhc17333quSm7iKqlWrquNO98vU1FR1bthwXrx4sWpIE1obeW4PHjyI6tWrKyurNU2bNlXnkMevU6dOLts+oeDoiURowTeP3eQzlPFu7l4vrXG8dynqHnnkEZOlzvoesQUT9HTu3Nn0nc8NikRen7bi9ARBELzF1TIkIAAfHssbNyfZLL0DEXUexJ+Kj9PidurUqTwxWflBN0XzHu2CLocuUnTZvPPOOy0af4ydosWMgsrZ+Dk2IM174Wkl4DKZKZFulvltG0Ugf7/jjjssxlNgPPjgg0ow6MKe+6wLM8JtpDulLsz0cXT/pHDVG7KFnS8/uN3WCR+2b99usuLRrVOH7qQUjLbc0Bh3d9tttzlcl6vd59j45qDD9dNCSmucLuoIjz3PsT1oqSSnT9t3MRGKB93VkbGb1tckLbHO8swzzyirOzsgvvjiC/UcZqcM/9Pd2NZ6K1asiEqVKqlOGHZA0MpGV2r+Zh2LyZg/dhawXAHX8cQTTyi3ZHPXSwpDWvW53XqmzFdffRV9+vQp8HERBEFwp6vlpCMncSkrx242y6HR7nN/F/JHRJ1QLOiupbSaFAQ2noqyHDamdBFnCwo+Z0WddZxLcHCw+p+WlubUtunxd9b7VLlyZRVTyW3lZ2Kd3ZHrsjWOopJiUE8AUtj58oP7aC3qeOwohKyPLd1TW7dubSH0dCh+rfe/OGPndCjmrDNc5gcb7ubnW/AcFFCMkbXOUOuMpcy8pAAFOu9LWrTpRsxslHSxttfxwPXqnRXdunVTbtW0vLEMBuvQjRw50mJ6ulhbQ/FGLwGdW2+91X3x1IIgCEXA2tXSlqDTYUzdXVXLKWue4BnkyAvFgm7lsJVYxBHWlqSCLkcXIp9//rnqMbcezC1YRSW/bdPFnnUNNDYsaRlwZ4KHosLjaJ1chcLQ+vxcvHhRWRsYJ2gLul9yXx0NzKbpjej7bytFfUmBoujxxx83faebMMfR9ViPK+X3Dz74wG5SoZ9++klNw/86/F13NSecn9NweYV1g2THgfVQEBd2LoOxthTrtJzR8kbruyNLMufRRR07S2iZo6WNsay02DmTJGXHjh2qnqMgCIK3licoTIkCPZul4DnEUicUC3QLpKhh9kD2cBfXcjg9U9HTqlRQ109Xbxvdr2ipWrRokYV1i0k7WDPNm2Mq69Wrlye5CF1RKewYh6ZbO2mBYDIb82yYnna/tIZup3R3o6WnINBVVz8WguegVZsiisltigItZozdZAIlum3S4kbrNeNNBw4cmKfDQl+vXlaDUAiyw0gv3eGMqBMEQfD28gS54+2XKHi9bjUMrpbbkW2OEZ4tjeXviKgTigW61DGGadiwYcr9j8lD2LPNEgGMv7LlqueK5XB6FgPm9Ewi0rVrVzWOsXY///yzmq8495Euf2wUUvjQekFxwZ5+6+yXxQljiWgdozhjIWRbMPMj9y05OVm5qBHGpdGVkkloWCuMLmzM6MoC5OHh4TaXQ4HNoTCcPXtWlWYgzEBIdzcKTW6PnmTG1r4wnrJt27aqEc/zwRgqFh+35RrnCCagoTXWUdydr7N7926L7y+99JIadJj91Hoa66Q2jF+0nobJjsxhchE9wUhB4XmnwHIkoBjXyeyvdI+kiGcdQutrki6TvJ4Yt8s4OSbQ4X3JmDpmkGVyoPzW26NHD3VP8x5nh45urRcEQfDV8gSDquZ6ozgqUfDxsbO4N7qCuFp6GSLq/AEWFn/hRP7TuBm6X1HsMOaKjW82oNiw0mPT3LUcZpZjdkOKEpYY4HS0tjCGxtXkt20sgcCYLLqf0apFCxcte2wcegrG89BC4ahu3w033KBEKeOQ9OQidCH75ptvVJIXiipavlauXKnEnjugUDBP9/7tt9+qgeKNtR7t7QtFGIUcE9UQbueyZcvy1AnLD7qVssaf4FnMM1DagoKfQpTXC69RdubY62SgOzAHcygAnV0vryXGwS5YsMBu3K4gCII3ibgZ8WfV55GxFW2WJ6CFTv9sDykc7p0YNEctOaHY0AP8rXu42UhlVkXdva84kkgIgi3o7kZLnT8W4E5KSlKNdwo7d4lWwTlYb5Ei3To2VYedJnye0ip/zz33WMT7FXW9vPZ1N1zzjhy6YDIWsaC1LwWByHteKC5mxZ/FC/uPq89j60ar8gTW4q1aSBDWtG2AwHyyY9PVUurMFY8WcBYRdV6CiDrB22FSixYtWijXRj1Lp7/w3nvvKVdPFm4XvB/WI1y+fLk6b+3bt7fINikI3oaIOqG4rHTX/rXXJOJKBwbYzWY5sX6MlCfwQVEnZh8PZ9Oj2xgHxmpIWmvBm2E81bRp00xujP4Ea/1RIAjeD+Nl6XbJepB0sabAEwRB8PeMlgUtT2CeCVPwDSSmzoMwqcSYMWNM3xmsLwjejHlaen/CXjZPwftgPOuGDRuUEGe8JV01BUEQ/DmjJZOfFKY8gRQT9y1E1HmQ0aNHmxrJ3bt39+qU9oIgCL7A7NmzPb0JgiAIXpXRMisnR8oT+AEi6jxIVFSUGgiLLguCIAiCIAhCYbGV0XLSESlP4A9ITJ0gCIIgCIIg+Gi8nA6/23KzvJSVna+rpeD7iKVOEARBEARBEHwwXm5YdAWz3y4U2M1SL08g+D5iqRMEQRAEQRAEH4yX0611/M/v9vj42FlVey7CaMwzSL25koGIOkEQBEEQBEHwsbIE/EzrXH5WOiJuliUfcb8UBEEQBEEQBB8sSzA17rT67a4q5dC/SlmHyxQ3y5KNiDpBEARBEARB8MGyBFJTTtAR90tBEARBEARB8GI3S0dlCcxj6wT/RUSdUGyw0Pp1110nR9wPyMrKQt26dfHjjz/m+W3lypVo06YNwsLC0Lx5c/z1119u2YZt27ahY8eOaj21atXCtGnTLH7v378/Xn75ZbesWxAEQRAK4mbJgXFxUpZAKCzifikIgsuZO3cuQkNDcfPNN1uM//777zFo0CC8+uqrePPNNzFu3Dj069cPhw4dQmCg6x5HZ8+eRbdu3dC2bVssW7YMW7ZswejRo1GmTBncc889apqnn34aXbt2xRNPPIGyZR3HIQiCIAiCJ90siZQlEBwhljpBEFzOhx9+iCFDhliMS0xMxPDhwzF58mQlqG644QY13bFjx1xurfv4449hMBjwzTff4MYbb1Tru//++zF27FjTNLQWVqtWTQlQQXAHd911FypVqoSGDRvm+e3ZZ5/FkiVLLMbddtttNqcVBKFkUlA3SylLIDhCRJ2fYct3u7hhQ6ZBgwYoVaoUbrnlFpw8mZvdSWfhwoVo2bKlsvTExMSohrimaeq38+fPo0qVKsrqopORkYFmzZrh9ttvV99///13tVw2pmiZ6dSpEzZt2mSxjuzsbDz++OPKQlO+fHm8+OKLeOGFF1CzZs082/vOO++gdu3aCAkJQePGjTFv3jybbqX57Zct5s+fjxYtWqh95X6xEZienl7gdfOY0d0xMjJSiSkek3Xr1qllc3t69+6tRFVR53OG/fv3Y+PGjejbt6/F+I8++gilS5fGiBEjTON4fklcXBxcyU8//aTOQXh4uGkcLYLcNloFdbiNX375pUvXLQg6DzzwAFasWGHzgKxatUpZinUWLFiAqKgoOXiC4CftL3GzFFyOJniMhIQE7fDhw2qoV6+e1rBhwzzTZGdna3v27FEDPxeVmcfOaJVXb1XDrPizWnFy7733apUrV9auvvpq7dtvv9W++uorrVKlSlrfvn1N08ybN08LCAjQRo8erf3000/au+++q4WHh2tvv/22aRrOy2l+++039f2ll17SoqKitOPHj6vvc+fO1d555x01P4ehQ4eqZRw9etS0jLFjx2pBQUHam2++qS1fvly7/fbbtejoaK1GjRoW2zx16lTNYDBoL774orZixQrt/vvvp7pU8xRkv2wxc+ZMtazhw4er5XHee+65R7t06VKB1s11tWvXTluyZIk2ffp0LSQkRHvggQe0Fi1aaAsWLNDmz5+vlStXTvvf//5X5Pmc4eOPP9bKly+fZ3zTpk21Z599VsvMzDQNJ06cUPv09ddf55me17v5tLYGe/dExYoVtQkTJliM09dlfvyWLl2qGY1G0zEXfIdhw4ap82k9vPHGG8Wy3sceeyzPb6VLl9beeusti3F8vjdo0MBiXFxcnNatWzeLd0H79u3Vc956WqHk4+r3vOA9HDqbpE1Yvld7ZO5mbeBPOy3aXxz079bD9KOntaSsLJtDapZcIyWVxo0bq6GwiKjzIK+++qpFY4QNUXc+7PkgaLF+l+mhwc9pxfgCoZCgkDIXVxRrgYGBat9ycnK02NhY7ZFHHrGYb9KkSVqFChW0jIwM07j+/fsrIbx+/Xo1/+zZs22uUxcGdevWVcsh/M5jTYGhw2VXrVrVQtRlZWVpVapUybM9PXr00Nq2bev0ftnbLi777rvvtvl7QdcdHx9vcWx4PW3atMk07plnnlHHtqjzOcPIkSO1Tp06WYw7d+6czQa4PmzYsCHf+8PWwP2wBY/9Rx99ZDEuNTVVzUPRr3PkyBE1bt26dQXaR8HztG7dWuvZs6f2xx9/WAxnz551+3rZAVK/fn2L8fq1ZN5pYE/U8dqcPHmyxT3D+WxNK5R8RNSVTBb8fVSr/fwPWo1nl2nVn1+mVV6+ydT+av77Tq357/+1x6yHlsXcPhNKhqiTRCkehC6EdIMj3bt3h9FoLHbfbWZaGhZdAcVF/fr1ERsba/reqFEjlSnxzJkzuHjxooqvuvPOO9U4nS5duqiYqKNHj6JOnTpqHDMZNmnSRMVl0YVJP456koznn39eZV6kC6TuunngwAH1n8vhNLfeeqtpnqCgINx0001Ys2aNaVx8fDxOnTqFO+64w2If6Mb34IMPKhdO/Zw52i+6VVrzzz//qGVbx50VZt10+YyOjjZNw2NEd8qrr77aYpx+LBhrVpT58oPbXa5cOYtx27dvV/9Xr16t3Dp16E76/vvv46qrrsqzHMbAMcbIERUqFO3apestOX3afgyD4H3wetyzZ4+6H6655ppiX+/dd9+NmTNnKldeukeTXbt2qf+2rmVrli9fjkmTJpncxenizKRCR44ccfMeCILgbhfLGfFnkZCSgdnf7UROTm77Izs6Agj9r413KuO/No4tpPacUBhE1HkQxk/oMRQUFe7Enu/21LjTGFS1HEICiie80jpmJDg4WP1PS0vDuXPnTCLOFhR8uqirWLGimu7rr79WcSvm3Hvvvdi9ezfGjBmjxBbjqoYNG6bWQSi0zBv09gSCHhPH2DxzKleujMzMTLW9/JzfftmCsYGkatWqNn8vyLoZN2i9blvjKDIpBvUsk4WdLz+4z9aijueOItT63E6YMAGtW7e2EHo6FMPW+29NgJ3rlrGS7CQwR48NNM90yVhFfZsF34FiKiUlRXWeeGK9zODKWDl2HD3yyCPqt507d6rngHnnji0YM3v48GHVqaKLuvXr16t4Xt5r7GBgEp+///67WPZJEISic/hcMhb8fQxr01OwrUzueymwWjgCDp5H8r51COw0BNbdolWCA7GmbUMEB9juMDXmmUMQHCOJUvwEWuRspcjVe4O8AV0IfP7556pBYz2YW5B+/vlnldmQCVJoldMb5ampqSpJBpOr3HfffSpJCkWDeQNfFwq6iNSx/q4LLl0E6rDRRRFeFCuRLijtJVNx57qL4zxaJ1dhY9Xa0sdzwobx4MGDbS7n9ddfV/vqaGA2TVtQzO/bt89inP5db0wTfTutRajg3VBA6eeZ15Y+6FZ5d6+X1jha9s3rMNJSR++B/Fi7di06d+5s+v7cc8/h+PHjykpHgceOKxF0guA7zP37KDp9uxlTj57GtpD/nkFZtUvh7JIJuPDjVGRs/yPPfLTWLTmTiAij0eYQapQmulAw5IrxA2ilY+0Te7zvwUyY5rCxzRTztOpQiFkPujUnKSkJI0eOVNY3CjgKnddee83UC56Tk2OylJE///xTLVOnevXqShTRBUqH1i8uyxxmZqS1aNGiRRbjmTGyVatWRXKX5b5SuNnLvOjOdbubevXq5clmSesFG90HDx40jZs4caKymplnw7R2v7Ql7s0H/bxbwwY3zy9Fvvmx47bp7nK6K66+zYLvoLs6Mv2/ucjXRZezPPPMM0rQs3OIxepp4WfmXdY3tLdeegmwY4juknTX1rPV8jdr18uePXvi2muvVdc972mW8+B1ycysgiCUDAvd838fQmb90siuXwZaYA7S//g198dQIyJuH4xSjzyH4Kuv8er2l1AyEPdLP7bSeZvvNl3pGGdCsZaQkKBi5Tju33//VZY5vaYTG2IUCGwg0V3wvffeU/MwvobWPDbK6HoZERGhrEGvvPKKhZsj3Qgfe+wxZQnSG3Ssl0ZLkrk7H4UTSx1wWooPlgD47rvvVO+8eQ99Yfd1/PjxarspQJlan41D7iNjzChg3bXuovDFF18o6xgbqTVq1LA5DRux3Lfk5GR1DkjHjh1VQ/ihhx5SltVffvkFU6ZMwQ8//GBRdsAcCnwOhYFxhzyO/fv3V7GrW7duxfTp0zFr1iyL6ViUnFZTWnwE34ECivf5J598YjHeGUuZeUkBdgjRWs7rkdcKr1da+e3FcpoLNxa3ZwkQWt5YC5GWYHY2mfP999/nWQafN+zQsAVdMK0tzIIgeJ+bZXxCCmLKhuN8ajoyapWClp2NtNUrkPz5dGSfOIayU2YiuFkrBN7cCYHBRq9vfwklAxF1fsBdVcqhf5X/4oi82XebsSqsZUZRQKFFwUMrChtchA0oFpZmY0mP/7rnnntUjScKpM2bN6ti0mxcDRgwQDWSKB7eeusti/XQ5YlxbVwPBRaFCi14dH8yZ9SoUcq184MPPlANMbpG0brWo0ePIu8rk7tw/xhXRrFE8cbGoW5ldOe6CwutoIyvc+TmxuQ1FKVsNPfq1UuNCwsLU+6yFFu0cLBBvnLlSiX23AGtKewI4DFkQhzGH7ITgNeKOXT/ZC0+wbegRa5Dhw7Kgl+UZfAZwbhKWs5oqWcCJUfJeTgP71HC5w+TtLCThc8Odso4kyRlx44dhd5mQRA8x9ebjuH5RTuRfSX5CcmKCUfaha1Imv0hso8chCEiEhHDHkZg3Stu/sFGGPcmIvB4ivpqNBjweq8m6NMq2uvaX4LvY2AKTE9vhPBfDzPdf6wb0cyUqLvs2UsMIRQN3gZ0a2SjbM6cOXI4iwiT19BS582FvenGS7FHYecucelrOGvpYrbed999V31+/PHHlUDXn10UR0888YRTy6HQpqtsQaCrNC1qb7/9Nh599FEUFlrnPvvsM2XBZXZWPassLeQDBw7MEwOqr5eZd+kaTNgpxBjgN954Q3U8MS7XOgGTIOSHvOe9n31nLqPbd1tVW8EYlwSDBmgG4NzBb5DyzedAaCjC+w5CRP8hCChtmXAsUgO6n9ZQo2wY+reORa0KuR4sguCsFnAWsdQJfgnd8ZYuXYp27dqpxhotZWzYscEmFB2WoGjRooVyb9OzdHobTElPIS+CzregeyLvWUdWMZbKoGCkeyQbYX/99VceN19a3NatW6diLGnZ5fVAaznLFVAw0uqf33ppNacIpNs0XbxF0AlCyXSz/DMnXcXNkYwDuxEWEIPsGqURVu92IDsT4YOGw1jOdodOkgFo2ylWXCwFtyOWOi9BLHXFC+P0mKCDrlB0cWTChVdffVXFtgmugVaQpk2bWmQt9SbYiKcbKIWd4DtQsDFjKjsMbJW8YMkBCi/2dNLtl/G55mUsirpexunSRZxQMFLM0ZWb5TposRSEgiKWOu92s9QCgPSOVZAZtw9Jsz5Ext8bUO62JxD08CCL2nMWpGVjQvkK6H91jMnFUrJZCvkhljpBKARMjPHbb7/JsXMj5gXhvRF7WTcF78Y8A6Ut6KJO10kmyGEMJZMMuWq9zOKqCzrC9dBaRxdMZ+LpBEHwraLh6jfDGVwcPxnpv6/mTY/QG2+GsUs7+4KOhBphiI5QpQkEobgQS52XIJY6QRAE18BSFiwdwMy47du3t5ttUhC8AbHUeZ5Z8Wfxwv7j6nPgnkQEHktWn5P2rMH5ZZNplkfIdTcgcthDCKxVF0jNQvCGMyq2jslPlj16HWqUt3TxFuucUFDEUicIgiAIZq7VtMTfcccdqlzB/v375dgIgmA3do7Fw2cGpVGdqXHppZJhyEyFMSgMQe2vQ8iFzYgYfB+CGpolkgoLRE7VcIQcT8HEPlehcaXcOrqC4EkkUYogCIJQYmCJkA0bNiAyMlLVUpwxY4anN0kQBC+OnUuPDkdW4yhkXziPlHmzkLL0G5S66R6UvaovtOaxiGqXm+nXmtBGUfih39VoIIJO8BJE1AmCIAglhtmzZ3t6EwRB8OK4OdItLFwJuixNQ0ZFDckzpiJl0TwgLQ2B9RvBeE0rZJeNcBg3x6yWf2SmowHESid4ByLqBEEQBEEQhBJdnmBtegq2lcmt9bs+M1Blt0zNOYZzw5+GlpwEY806iBz2MEKu66KSIGl7ExGy6oTdouFqvBQOF7wIEXWCIAiCIAhCibPMPfrnAfy46xQMRy8j47oqaryWnoZfk1JhDAAM7Zsj8OdGCLu5F0Jv6AGDWbbKsIZR6F62jBQNF3wGEXWCIAiCIAhCiWLqgZP4PiMFqF8aAaUCoRlzkPr9IiR/OQPBbduj7OCnYSgVhnKTP7E5vxQNF3wNEXWCIAiCIAhCibLSTT9+FvSO1LKzkbzjFyR/MR3ZJ48joGw5BNVtiJyoYAT/ckKVJTDHvESBuFcKvoSIOg+SmJioBpKZmQmjFKkUBEEQBEEoVNxcfEIKYsqGw1i7lLK0MaNlwpMjkR13GIbIUoi8bxTC+tyFgLArNeWqhcNwNLcmHTEGGKREgeCziKjzIFOmTMGYMWNM3ytWrOjJzREEQRAEQfAp9NIETHxCcqAhPbkUDOXLKKucsUIlhF53I8IHDEFApGWmyqhmFTCgdmWcSkhDTNkw9G8di1oVIjy0J4JQNETUeZDRo0dj6NCh6nP37t3FUicIgiAIguBk0fCVKck4dDYZRk1TeSjTju5Awl9zgZ8jUPbt6SqLZdRbH6n/tjiTmYXoplXwYnQFOeaCzyOizoNERUWpgQQFBXlyUwRBEARBEHyuaDjKlUHGwT24/M0nSIvbBhiNCGvSG1pmJgxBQZaCLjULYX+cxdie/5UnkLg5oaSQW7BDEIoBWiWvu+46OdZyPOySlZWFunXr4scff7QYv3LlSrRp0wZhYWFo3rw5/vrrL5dfR5s2bcKQIUPU+tkIeOmll5ya77PPPlPTWw9r1661mC4nJwf3338/IiMjceONN+LMmTOm3/r374+XX37Z5fskCIJQkpKfjN17DE9uPoxMg4as2qWgZWch8eUncPbth5EWtx1hbbuh/OeLUfrxl5Sgy0NYIEbd3Qx3t6mOCKNRDaGsbSAIJQCx1PkhmSdPqv9BVat6elP8Gjbi09PTPb0ZXsXcuXMRGhqKm2++2TTu+++/x6BBg/Dqq6/izTffxLhx49CvXz8cOnQIgYGue4StX78ef/75p+p4OHfuXIHn//333y1cqBs3bmzx+5dffok9e/Zg6dKl+Oqrr/Diiy/i008/Vb89/fTT6Nq1K5544gmULVvWBXsjCIJQsph38jw+OHVelSgwhOYAoVdsbMFBCOl0I8p0vBuhORWBPQZgT27R8JEda+HxbvUtliOWOaGkIqLODwVd3JB71ecaX3wuws6D1KlTx5Or90o+/PBDZS3TYXbY4cOHY/LkycrKRaKjo9GwYUNlrevQoYPL1v2///0Pjz32mPpcs2bNAs/frl07hyKT2/vcc8+hS5cuyup4zTXXmH7j92rVqilRO2rUqELugSAIQsnMZtn76mhMPXoG2adOIGnOJ8jcuwvlP/0KBmMgyrwwTv2nayV+Ow3DlYQpzGQ5uHWuRU4Q/AGxOfuhoMs8dkwN6vMVq11xsmTJEjRo0AClSpXCLbfcgpNW27Bw4UK0bNlSWWxiYmIwduxYaFruQ/r8+fOoUqWKSjKjk5GRgWbNmuH22283WUy43EqVKqFMmTLo1KmTcq0zJzs7G48//riyipQvX15ZTV544QWbjfl33nkHtWvXRkhIiLK+zJs3z6ZbaX77lZ87qrPLmT9/Plq0aKGOD4/FXXfdZWHxc3Z7eZzpakh3QAopHsd169apZXPdvXv3NpXcKMp8zrJ//35s3LgRffv2NY376KOPULp0aYwYMcI0jtcEiYuLgysJCHDv45DXFi10CQkJmDVrFmrVqmXxO/eb1jxBEAR/j5nrOnkdPl53EMt2nFT/b5z5M/6ZNBbn7u2NtB+XIKBMFHIu5r5nlKAjYYHIjon4rzRB36aSyVLwK8RS54eCzjTuirArTovdgQMHlEibMGGCqs336KOPKsvEt99+axIsd999txpPVzu6q1FwhYeH48knn1QCjNYcut/deeedSmRweceOHcOKFSvUMo4ePapc2XThx2V27twZ+/btQ2xsrBrH9X/wwQd444030LRpUyUetmzZksfSMm3aNOUaR8HXsWNHLFq0CIMHD1Zi0NxFML/9ctXxoRigwKH1itOkpqZi8eLFSlhRxDm7vRRQb7/9trKAnTp1Sq2Hx5jWJM7L+K9HHnkEr7zyCt5///0iz+cMq1evVueXglGH547xZhT1jLcjly5dsptciOvnkJ94c4eAowWRnQ4U0tx/Xp/mPPDAA7j++utRrlw5VK5cGT///LPF79dee6265i9fvqzEsSAIgj/Fy82IP4uElAzM/m4ncq5Y20ji+vm4+PdCID0dgQ2aIHLEKARf3c5mRsvA+mUwsk5lDGpdXQSd4H9oglfQuHFjNViTnZ2t7dmzRw38XBgyTpzQ9nftpu1p0NDmwN84jbu59957taCgIO3o0aOmcW+//bYWGBio9i0nJ0eLjY3VHnnkEYv5Jk2apFWoUEHLyMgwjevfv79Wr149bf369Wr+2bNn21wnl5uZmanVrVtXLYfwe8WKFbVnn33WNB2XXbVqVa1GjRqmcVlZWVqVKlXybE+PHj20tm3bOr1fjo5Hhw4dnF4OB27P3XffbXN5Bd3e+Ph4i+PJx8GmTZtM45555hl1Poo6n7OMHDlS69Spk+n7uXPn1LLtDRs2bMizjFdffdXhPBy4H/nB6+DFF190artXrFihjRs3Tlu1apW2dOlSrXfv3mo9ixcvtnmO/vnnHy01NTXPb0eOHFHzrVu3zqn1CoJQMnDFe96XOXQ2SRv4006t8uqtaoietkar/sz3Wo1nl6khsutALbBWXa3M2MlapV+2mKazN8yOP+vpXRIEl2oBZxFLnR9a6PJMU4wWu/r165usZaRRo0bKAsNMgBcvXlQWN1o4dKsMYQwSrU+0wOlxaLRINWnSBDfccIOyyun1/sjZs2fx/PPPqwyKdF3UXTdpBSNcDqe59dZbTfPQ6nPTTTdhzZo1pnHx8fHKGnXHHXdY7AOthA8++KBy4dQTYzjaL7pIuuL40G2P22Mec2ZOQbaX7p20LOnwuNKd8uqrr7YYpx8/vUe0sPM5A7edViyd7du3myx45pYrupPSCnjVVVflWQbj7m677TaH66lQwbX1iHjdcNDh+mklHT9+PHr16mUxLY8/z7EtaKUkp0+fdun2CYIgeKtlbkf8Rfy84gBS21cGgozQ0lJxcetipH23EFWHTIYhNALhj4xCeKkwGGx4WIRmaxiaFYL+rWJRo3y4GieJUAR/RUSdUKzodfl0goOD1f+0tDRTxkGKOFtQ8OmirmLFimq6r7/+Wrm1mXPvvfdi9+7dGDNmjGpA0z1w2LBhah1ETyWvN6LtNfb1WDbG5plD1zm6RnJ7+Tm//SoIjpZD1z5S1Y7wLsj2MtbQej22xlFQUgzqbqmFnc8ZuI/moo7nmyLI+nqg22nr1q1tuihSQFvvf3HHzhGKOboNFwS6zxbmmhEEQfDFTJZvHMp9ZwU0ioJmzEHq4oVInjsDOefPwVi2EjITTsLYpjkMZXLj5GyRZjSgZqNKaFxJXNYFQURdCYeWN1rgHFnrgmJjvSITpt6g//zzz/Okg9etRDqMR/rmm29UghRa5Xr06KEShzDG7KeffsLs2bMtLFq0AurojX7rtPXW33XxRBFIq6AOLSm07Lna4pMfugileLNlpfK27S3M+TdPsEJhaG3p43lk7OTEiRNtLuP1119XYt4RFP2sLedt6PtuLmwFQRBKopWOmSx1Uv7dgMtPvoucUycQULY8So16BmFdeyH4zwTgeDKMJ1LUdEaDAcsevc5kkdMRy5wg5CKizs+FnbcIOl20Ma07LTT2XAxJUlISRo4cqaxvrFlGAfjaa6+phj6zQDJRhm7hIqw9xmXqVK9eXQmc5cuXq8yYhJYsikHz+ZhlkZYfJhsxtxYx+2OrVq0sapIV1/GhcGOGxG7duuX53du2t6DUq1dPbasO3VAp7A4ePGiy0PIcM+mLeTZMT7tfWkO30++++05lcC0IdAvWj4MgCEJJEnGT/j2O7ccuouK5DCRWCMbJoP9CLHISL0BLTkLkff9DeJ+7YAgLU+OzYzIQeCxZhUOrbJZ9rhKLnCA4QESdHws7bxJ0ulvcpEmTlFhj/Bhj5Tju33//VZY5pvonzzzzjGrsMwMjXf/ee+89NQ9jxxjbxcY0rTURERHKssNMhOYui3QJZD0yWnVoFaG1jxk1aRUyd82jCKILHaelkGCmTTbWGavHwRPHh3Fa3FeKT6bAp4jlcWGMGd0RvWl7db744guVrZPirEaNGnanY/ZH7l9ycrI6d4xLo1X1oYceUtbYX375BVOmTMEPP/ygXGptwU4BDoWBcZYszUBSUlJUtlSKTG6LeeZQ6/1hDGjbtm3VdcTzMWPGDPzxxx+qaHpBYPZVWmPtxdwJgiD4Ik/9dRAL05OBICAg4SJS/92GpK9mIuqVN2GsGo2wW/og9MabERBp6UIZUK80bilXBjXKhqF/61jJZikI+SCizk+FHfEmQaczaNAgVZeMjXsKLYoXWi6Y1p6sXbsWH3/8sWow67Fc99xzDxYsWKDEzubNm1UBZ1ryBgwYoGqDUQi89dZbFuthEWjGqHE9FEtspNOCxxp35rCcAGOcWP6AViJajGgpo7unJ2BCGB4TxpVRXFC83XjjjSYLo7dtL6HllPF1esIaezDpDYXpqlWrVExaWFiYcrFlkheKKor1lStXKrHnDhiHyY4BHZaR4EDhduTIEbv7QxFGIcdENYTbuWzZMlVjsCDQrZQ1/gRBEEpK0fA2dcrh24uXgFAjMrZsRNLMacjcuxMIDkHmP7uVqDMEBanBmoygALTtFIOh0d4dOiAI3oKBKTA9vRECTDFQbFiawwbkP//8Y3K/c0WSB73guLcJOk/C24AuioxVmzNnjqc3x29h0hta6vytCDddipnEhsLOXaJVEATvxB3veU8VDX9+0U5km9WYy46NQErgSVz++F1kbvubrjIIu7UvIgaPgLGCVVKr1CwEbziDQBjweq8m6NMqWsXLhRp983gIgqu0gLOIpc4PETEHbN26FUuXLkW7du1UPB2tXkyhz1IJgudg6YoWLVqo5C56pk5/YObMmapTQQSdIAi+aqGzFnSaAcisXQo523Yic8dmhN7UE5FD7lfWOZuEBeK6TtUxrlUtcbUUhEIgok7wSxgnxTi9d955R7krNmzYUMVPdejQwdOb5tfUrVtXCWu6MvqTqGOtP8aGCoIg+JKQm/v3UazLSMXF1CxksTYpvYHOxyPx97kIv7kfjKHRCG53Hcp/sQSB1WJM8xr3JiLweApGdqyFx7v9F0csljlBKDwi6gS/hHFQv/32m6c3Q7CBeSF5f8FeNk9BEARvdrVMjw5HVuMooHQwEHYZSd/MRPLu1YCWA0OL+ohER5WEzFzQkeyakQg5noLBrasjwsszMwuCryCiThAEQRAEQci3NMGM+LNISMnA7O92IhsasmqXQvaFc0ie8ylSf1jEAqMIqd4Mpe4cCeNt19lfWFggbru1nrhZCoILEVEnCIIgCIIgOGTeyfN441BuorXAalfKyoQakXMiAanff4OghlehdNd7ERneEBpzm6w6YSoaPv2eVvjr8AUcT0xDdFQo+raKQe3yEXLEBcGFiKgTBEEQBEEQHFrpph49oz7nJF1G4oZ5CKxZB6GNb0JQ7Xoo99GXCKzXCIa0bGi/nYYhBxZFw7s1qKwGQRDch4g6QRAEQRAEwcS+M5fx0rY4XEzNQOfgMITULo0TiZeQsvgrJC/4HNqliyoBSmjXm9T0QfUb584YFohqTSugtSEYMVI0XBCKFRF1hWTnzp2qyHCVKlVMRYfdAQOMObCOWlZWlqnItCAIgiAIvg/f7ebve29IgvLUxoPIaBQFBAVg945jSPr+ZyQvmI2chPMIqFwVpR54HKHdb7U5f3L1CLzdvjFCfLTeniD4KiLqCsno0aNRvnx5uBs+4Cnk0tPTcfz4cURHRyMwUE6bIAiCIJQEQcd3O+G73tOijmUKnlu8Exkd/nOVzKxiwOXnpiEgIhKlHnseYbf0gSEoyO4yzmRmYf7JCxgaXaGYtloQBCLqoBAsXrwYhw4dwvDhwzFnzhy3X0nVqlXD0aNHVT21gwcPun19giAIgiAUH0ajUb3rPSHiFvx9DPEJKYgpG67cLdMrhyJ1/UogPU0JOGOliig78QME1W8EQ2iYxfwRGnDjGQ3VryQ/qVE+3FRvThCE4kVEXQHJyMjAU089hYkTJ2Lv3r0oDkJDQ1G9enWcOHFCrZ+umIIgCIIg+Da6Nw4FHd/1nqg1l52T26Zg2yL1wF9I2PIVso4cUG6WoTfdDoMxEMHNWtlcRrIBaN8xVqxyguAFlBhRt3nzZvz888/YuHGjGnR3hvwEUGpqKiZMmICvvvpKWcPKlSuHHj16YOzYscrV0ZopU6agYsWKGDBgAF577TUUF3zY165dW+2PiDpBEARB8H2KM45OrzNHuoWFmwQd2xRpR7Yh8bc5yDj5LxASivCBQxEx4F4l6HSMexPxUO3KeLxbfYvlilVOELyDEiPqKMKWLFlSoHnoznjDDTfgzz//RNWqVdGrVy8cOXIEs2fPxrJly9R4Cimd06dPY9y4cVixYgU8hbcEUguCIAiC4P3oLpZr01OwrUxu8pL1mYEmCx1ysnD+x/eQnZyIsNv7I2LIfTCWr5hnOTm1SqFfqxhEGI3FvQuCIPiTqLv22mvRrFkztGnTRg01a9ZUyUUc8cYbbyjhxnlXrlyJyMhINX7y5Ml48sknVczc2rVrTdO/8MILyorH6QVBEARBELwZ3cUyCxrSO1YxjV+9fzsyt65H6Za3wmAMQoXbngQa1gQ6NrS7LC3UiD8y09EApYpp6wVB8EtR9+yzzxZoesamTZs2TX3+4IMPTIKOPPHEE/j888+xbt065dZ59dVXY9euXfjyyy+VCExMTDRZ+ui2wO/h4eFSbkAQBEEQBI9DV8tJ/x7HJ5sPI0DTkB0TAYQakXXkIJI++wjpv/4CGAMRUb8DjBFRCK3eFFoagFUn1Pz9ro5GVHgwjiemIdosCYq4WgqC91JiRF1BWb9+PS5evIg6deqoenPW3HnnndixYweWLl2qRN2BAweUEGzVKm+wcNmyZfHRRx/hwQcfLKatFwRBEARBsM28k+fxwanzQP3SMGbnID30IpInvI20X34EcnIQ0r4zIgc9gID40sAVL0xDDv9qMAYY8EinuqhVIUIOryD4EH4r6rZv367+2xJp5uMp7Mh1112HNWvWWEzz2Wef4YcffsA333yD+vUtA4cFQRAEQRCKO3YuLiEFP1cJgF5VILtuaVx6YjQyd29H8NXtEDn8EQQ1aqp+0wyJMBxNNi2Dgm5i36Yi6ATBB/FbUcdMlyQmJsbm7/r4uLg49b9ChQq4/vrrLaZhvF1ISEie8Y5o0qSJzfGsP0eroSAIgiAIQkEyWu6Iv4hVyw8gJ1tDVmwEMsJykPbrLwjv2Q8ICkDkQ08A6ekIbtnGYv6oZhUwoHZlnEpIQ0zZMPRvHSuCThB8FL8VdUlJSeo/Y+FsERGR63Zw+fLlYt0uQRAEQRAEZ90s3zh0Un0OrBYOw4HTSPz9O6QsmQ8tNQWBdRsguHEzNdjiTGYWoptWwYvRFeSAC4KP47eizhWwTl1Ba9Xt3r27QBY8QRAEQRAER26WOakpuPj3QiR//Tm01CQYa9RC5NCHENTwKos6c4HHU2A0GPB6rybo0yq3Fq8kPxGEkoHfijo922VKSorN35OTc33MS5WS1L2CIAiCIHje1fLRPw/gx12nEHDksspomVU1Sv12edpbSPtxCYxVohHx2DMIvfEWGKzqyQXXL4P7alfGoNbVxcVSEEogfivqqlevrv7Hx8fb/F0fX6NGjWLdLkEQBEEQBGvL3E/JSdhbPlBltAzIyEBy3O8IbXCzEm8R/YcgqEEThN3cG4agIJsHL9VoQHSjiiLoBKGE4reirnnz5ur/li1bbP6uj2dBc3fB+nZ6zbvMzEwYrXrVBEEQBEHwX6yLh2vZ2UhbvQLJn32M7JPxQKABYd1uRWCN2mowd7NUnw0GLHv0OlVjTn3XU2IKglDi8FtR16FDB5QpU0Zlndy2bRtatGhh8fvChQvV/9tvv91t2zBlyhSMGTPG9L1ixYpuW5cgCIIgCL6T0TIhJQOzv9uJnBwNWTHhSNu4VhUOzz5yEIaISEQMewghHfJm386uGYnAY8kINBgwsc9VaFxJwkgEwR8waJp2pexkySI0NBTp6elwtHsvvfQSxo0bh/bt22PlypWmjJeTJ0/Gk08+ic6dO6uyBcVhqevevbuy1O3du9dt6xMEQRAEwXtdLOMTUnC+QjDWBGWp8YF7EmGMT8aFc6uQNGMKGzcI7ztIuVsGlC5jd3k3ZAZiXKta4mopCD6EnjTRXlJFv7HUsQj42LFjTd8zMjLU/2uuucY07uWXX8att95qIepWrVqFDRs2oF69eujYsaOqS/fXX38pq9msWbPcus1RUVFqIEF2fOAFQRAEQSj5LpbZORq0ACg3SwQZkbFjC3Ji6iDYEIGwdj2Rc/EcwgcMhbFceZsZLa3dLEONAR7cK0EQipsSI+rOnj2rxJg15uM4jbU1b82aNZgwYQLmzZuHxYsXo1y5chg6dKgSiPYKkwuCIAiCILjSxZJkR0cgM24fkmZ9iIy/NyBi6IMIuGskAoICUOqhJ226WoYcT8GbvcXNUhD8nRLrfulvJldBEARBEHyDWfFn8cL+4yYXS8bApZ+LQ8LehUjfsAYwGBB6Qw9E3PsAAmMcZ+F+ulpFPNkgt+acIAi+i7hf+jCS/VIQBEEQ/Ct2bu7fRzEzKI2pKdW4rNqlkPHbbzj79auApiHkuhuUhS6odj3TfNYZLaV4uCAIJdb90heR7JeCIAiC4F+xc+nR4chqHIXsM6cAoxHG8hUReG07hJ67BeF33IWghrmeO+aENYxC97JlUKNsGPq3jpUEKIIg5EHcLz2IZL8UBEEQBD+JnfsyNxlKStMgpHz7OVKWfoPQrregzNOvARnZQLDjWrUT68dgaHSFYtt2QRCKF3G/9GEk+6UgCIIglFzmnTyPNw6dVJ8DorJxedVXSJqyCEhLQ2D9Rgi9vnvuhMFGhy6WapwUDhcEwQHifikIgiAIguDimnNxCSn4uUoAqMWyjh/DhTGDoCUnwVijNiKHPYyQjjfAYMiNqyPiYikIQlEQUScIgiAIguACN8sd8RexavkB5GRryKwciIxAwFixMozVYhDS7joEX9MRoV1ugsGY19UyyQC07RQrLpaCIBQKEXWCIAiCIAhFsMqtTU/BtjK5xb6NlYKQ+uNiXNz0DQKr10LZd6Yri1yZlyYAqVkIXnsaBrNiUuaFw8XFUhCEwiKizoNISQNBEARB8O1sllnQkN6xCrTsbKT9shzJn01H9qnjMESVRfC1naDl5MAQkCv4EBaInKrhqi4dMQYYMLGPFA4XBKHoiKjzIFLSQBAEQRB800JHQcdsltmxEdAMmTg/8m5kHzkIQ2QpRA59BGH9BiEgLDzPvKGNotC9nJQnEATBtYio8yCjR4/G0KFD1efu3bvDaMPHXhAEQRAE74mbI4l7LiArOwfZqYnIql0FhhAjghs3Q0CHLggfMAQBkaXsLkdi5wRBcAci6jyIlDQQBEEQBN+Lm6u+eSNOz30P2ZmXUa7nIlVsoPRTr1jMZ12iQI+bU9+lPIEgCC5GRJ0gCIIgCIITcXOZe3chadY0nN78FxBgRNhNvaClpcEQEZln3uyakSp2LtAgcXOCILgfEXWCIAiCIAgA9p25jJe2xeFiagZaZBnx7aZ45GhAVkw4Lr79EtJWrwAMBoR2uxWlbhmCgOb17R+3sEB06FID41rVQq0KEXJ8BUFwKyLqBEEQBEHwe2iVe2rjQWQ0igKCArB3TwIC0lNhCAlDdp3SCKhcFSGduiJy6IMIrFkHIVkaAlafRA5V3xXoZvl6rybo0yo69zsMCDVeyXwpCILgRkTUeRApaSAIgiAI3hE399zincjoUFl9zz51Ahe/n46M339F5TFzgFAjIu/7n6o3p5MeaMDTdzdD1uHLiE9IRUzZMPRvHStWOUEQPIKIOg8iJQ0EQRAEwfMwEUpG1XBkpyQg+ZMZSP1hEZCZieB6zZFRLls1lswFnc7c84n4o3sjhOh16ARBEDyEiDoPIiUNBEEQBMEzMXOdg8MwuE11ZVmLS0jBxW1LkDTnYyA9DYENmiByxCgEN20NQ4j9ptKJ9EzMP3kBQ6MryGkUBMGjiKjzIFLSQBAEQRA8FTN3GrPW7sfEfi2RWCEYhshwBEbHImLYQwhpf73JMmdemsBe3JwgCIKnEVEnCIIgCIJfxcxp6Wm4uHUxUr76Ak9enoRSA9sg7NY+CLvtDhisXClzakViYJVySEnPlrg5QRC8FhF1giAIgiCU2KLh8QkpKpEJY+Y0Yw5Sl3yL5C8/Rc75czCWrYS0oMvIzM6BwWi7SaSFBqJRs2riYikIglcjok4QBEEQhBJDWnYOHv3zAH7cdQoBRy7DoAGaAUgOPYBLb7yFnFMnEFC2PEqNegZhXXshZGMCbj2t4e1+zRB3PgWLtsTjeGIaoqNC0bdVDGqUDxcXS0EQvB4RdYIgCIIglBjL3E/JSdhbPhCoXxrGzGwExacgOyYCuGyElpyEyJGPIrz3QBjCwtR8OdUyUDMqDBFGIxpXKoXGPRp5elcEQRAKjIg6QRAEQRB8PgnK84t2Igsa0jtWgaZpyPjjVyTN+hDl2o+AofONCA7pgArzlyMgPMJi3uzapdD7StITQRAEX0VEnQeR4uOCIAiCUHQLHQVddo6G7NgIZOzZjKRZ05C5ZycQHIL0kPMIDc3NUWmwEnRECzXij8x0NEApORWCIPgsIuo8iBQfFwRBEISCx8zNiD+rPo+MrahcLino0s/F4cKKz5Cx/W8gMBBhPfsh4u77YIwoh+BfTqjYOtI8pgyiy4ZLzJwgCCUKEXUeRIqPC4IgCELBMlmerxCMNUFZanxEoFGNIzlVwpCxaytCb7odkUMegLHqfy6VOVXDEXgsGcYAA97r10IVHBcEQShJiKjzIFJ8XBAEQRCci5ejNU4LgIqZQ5ARWUeP4Ilxz6Nrqy7QjE1g6NAMFRb8CGO5CnmWkVUrEsEnUjCxT1MRdIIglEhE1AmCIAiC4PXxciQ7OgLZiaeR9Pl0pP28DMjJwYnKlaF1aAuEGmEMzSvoFGGBeGRIc/RrIAlRBEEomYioEwRBEATBK9h35jJe2haHi6kZ6BwchuTUTJOgy0q9iAs/fo7Un74DsrIQ3LItIoY/jNRWrfF6YCRe/34PsrUrgXMAjAYDXu/VBH2uZLbMTZUiCIJQMhFRJwiCIAiCV7hZPrXxIDIaRQFBAdi75zSCjiWbfs+uGoa0T5YhqH4jRA4fheBWbdX4E+mZQI0I/DK6k1pGfEIqYsqGoX/rWHG1FATBbxBRJwiCIAiCx90sn1u8ExkdKpvGZVQ24PL8uTCGRKBU297QmldHuenzYYypDoPB0ur2ftxp3HVNIzzbo6EHtl4QBMHzBHh6AwRBEARB8G+Y2TKjariKi9NSU5E8fzbODe+JS+vnI3nvr8iKzv0tMLZGHkGnW+vmn7zgkW0XBEHwBsRSJwiCIAiCR8sU7D51GZl1Q5G6aD6S585ETsJ5BFSuilIjHkOZ4NYIOJGCwJOpduPl1DiJmRMEwY8RUScIgiAIQrEWDt8RfxGrlh9ATvaVJCixEUCIAckLPgcMBpR69DmE3dIHhuBg5OxJxOCq5RAVHizxcoIgCHYQUedBEhMT1UAyMzNhNBo9uTmCIAiC4FbmnTyPNw6dVJ+NVUKR8fMKZCacQOgLj8MQbETU+PcRGB0LQ2iYaZ7s2qUw9OpaaFCplJwdQRAEO4io8yBTpkzBmDFjTN8rVqzoyc0RBEEQBLda6aYePQNN05C+YR2SZ32IrMP7YQgvjeCchxCACATVqZ9nPi3UiD8y09EAIuoEQRDsYdD4dBU8bqnr3r27stTt3btXzoYgCIJQ4ph57AyeWrAISTM/QNa+XUBIKCI79kaZJj1hjChjmq5OxQg0qFIa0VGh6NsqBjXKh6t4uVCj5HYTBKHk0qRJE/V/9+7dhZpfLHUeJCoqSg0kKCjIk5siCIIgCC5PgBJTNhwD2sSiWrkwZaVTgu7gPwjrMxARg0fAGBYF42+nYbgSW0dualhZShMIgiAUEBF1giAIgiC4xL3y0T8P4MddpxBw5DIMGpBx+iDeHPcH7n7rdZzKyELpJ19GQEQkjFWqmebLjolA4JUi48YAgyoaLgiCIBQMEXWCIAiCIBTZMvdTchL2lg8E6peGFn8YSd98ipR/1qtpvtrZF8b69WzGzGXVioQxPhmBBgMm9m2KWhUi5GwIgiAUEBF1giAIgiAUiq83HcPzi3YiCxrSO1ZB1ol4JH8xHWmrlgM5OQir2w6l+o2EoX49+wsJC0SHLjUwrlUtEXSCIAiFRESdIAiCIAiFstBR0GXnaMhmnblQI5JmTkP6mp8QfHU7lOkyBOHGGlDRcqtO4JamVfF2v2Y2lyWJUARBEIqGiDpBEARBEAqcAOViagYykhJxefsKhFz3CAwAIoc9hPDb70Rwi9ZAahY0JkHJ4dwaakaFIULqsQqCILgFEXWCIAiCIBQoAUpOWhIubfwOlzYtgZaZhjI3NENop64IjKkBcCBhgaYkKJIARRAEwb2IqBMEQRAEwS5TD5zE9xkpKgFKQHIyUpd8hUt/fYuc9GQElo9BxAOjEHLdDTbnZRKU4BMpmNhHEqAIglAwMk+eVP+DqlaVQ+cEIuoEQRAEQbBrpZt+/CyUbyWAS39/i+Rfv0Bgmcoo2/V+hNx0K7KvKm//6IUF4pEhzdGvQbQcYUEQnBZvHBc35F71ucYXn4uwcwIRdYIgCIIg2IybQ2wozv74HUI6d1P15cIHDEaQVhqlYzvBYAyCdjINQadO5CZDuYLRYMDrvZqgT6toUxIUQRAEW9gSb/q4zGPH1Hh+FmGXPyLqBEEQBEEwlSdgNkstJxvJ+35F4qYFyD4Zj8iEC4gYPAIBZcoiZMgg4LfTzH2ikqDc1TYGUeHBiE9IRUzZMFU8XGrNCYKQH7bEW7W3J+HEU0+bxqnpjh0TYecEIuoEQRAEwQ/Zd+YyXtoWp7JYtsgy4ttN8UrQpe7/A4m/fYnMc0dhiIhExLCHENZn4H8zWiVAub9THRFxgiAUSdCZxNugwUB2dt7pRdjli4g6D5KYmKgGkpmZCaOkehYEQRCKySr31MaDyGgUBQQFYO+eBARqQMreX3Fu6SQYgkIQ3n8YIgbdi4DSZfLMLwlQBEFwpaAzYUPQmeYTYecQEXUeZMqUKRgzZozpe8WKFT25OYIgCIKfxM09t3gnMjpUVt8zdmxBdlRlGA1AeIP2KHN+IMJu6wftmtr2FyIJUARBKEQyFGJX0BVieZIZ8z9E1HmQ0aNHY+jQoepz9+7dxVInCIIguDXxyYA2sep7RtVwZMbtQ9KsD5Hx9waE9uiFkNtGK5fKqI53Q7sMYNUJ03JGdqyFx7vVt1i2JEARBKGgyVAYM1dYgmJjVcIUIpkx8yKizoNERUWpgQQFBXlyUwRBEIQSnPhE59PfDiHGcA7n1n6G9A1rAIMBoTfejIhBw5BVPhLG+GRVXJwJUFQmFIq3AAMGt66OCAkREAShiK6WTIJiKxmKBXzWWLlhWgs6yYyZlwAb4wRBEARB8OHacmP3HsOTmw8jSzMvNgCknjyIXycMVYIupEMXlPt0Acq8OB6BMTVMCVDMoaCb2FcKhwuCUDAhpw+2kqHowo5CzRol3ubNtfjNnqDTl6fGnfzPvdNfMWia1RO/gHzxxRcu25ghQ4bAX2nSpIn6v3v3bk9viiAIguDD7pYvbDmMNUFZ6nvgnkRg9xFomWkIKh+LHGg4v38BQrv2QFDD3PeOOYa0LAxJDUZKeraUJxAEocDoQk7LygJLVGadsC22KNSsLXa6eDOvVUfsCTrr5dXw8SLlRdUCRRZ1AQEBMBhcU1g020HGm5KOiDpBEAShqO6WTICSwgQooUZkXziPlDkzkbJsIUKiG6HKXeORFRuBrMa5bv/2mFg/BkOjK8jJEATBdVkt8xF2xFqUFTS5SpCPC7uiagGXxNQ1b94cvXr1KvT8ixcvxo4dO1yxKYIgCILglxY6xs9lRIcjJzMZyV9+jpRF84C0NATF1keZtn3BPlzj8WQYT6SoxCd9W8Vg0ZZ4HE9MQ3RUqPpeo3y4JEARBMHtgk4nqFIlkyXOWozp38W1shhFXYsWLfDqq68Wev4jR46IqBMEQRCEImS0ZPxcRkUN5+7pCe3SRRhr1EbksIcR0roTQn8/oxKgIIdxclCJT2pViEDjHo3kmAuCUOyCriBWNU7DaUu6+6XHRV3p0qURHh5epGWEhYWp5QiCIAiCULCMltNX70WtkBRkx0QjoGIUwm7pg8BadRF6Qw8YrmSsZAIUliuQxCeCIPiSoHNG2Imgc1FMneAaJKZOEARBsJfNckb8WSSkZGD2lzuRk5372tayM5G042dc3LAAhqAQlPtiEQwRITaXEZiRg4e0MAy6YqETBEHwVBxdUSxq1usrSYKuiTfE1AmCIAiC4B7mnTyPNw7lJgwIrBYOY9wlJO9Zh4u/z0XWxdMICC+DyBsGAsH2qxRlBQcgun5FEXSCUExcvHgR6enpqFSpUok+5s64RgZScBkAgzGwyALMfH2kpAg6VyCiThAEQRC82Eo39egZ0/fMmhE4O+5+ZJzcj4CQCER1GoJSV98OQ2gYmh7IwJ4T55Bt5oBjNBjweq8m6NMqWhKgCIIbSU1NRVBQEAIDA3HgwAE0aNAAjz/+ON5+++0Sedz15CUUVM64Rpq+u0CA6etz1fJKCm4RdSxNsGrVKhw7dgzlypVD/fr10bhxY1X+QBAEQRAE55KgGGuXwom0DOScPgljlWowhAch9OpOCD3fCqXb9oExNDJ3xhygY+3ymNq/hYq7i09IlTpzguBGGL2kl/SaM2cORo4ciZUrV6JTp06oXbs2Bg8ejHbt2pXIc2BdQ86esHOna6SIuWIQdQkJCejSpQt27txpMT4kJARXXXUVWrZsaRqaNWumkqQIgiAIgr9jnQRFMwBJO4/j8pyPkHXkICrM+wEBkaUQ9tCDCPntdG42yyswAUr/1rHKvfLZHg09txOC4Ac8/fTT+Oabb3Dw4EEYjUY0bNgQ119/vcl4wf9ffPEFSiLWMW38bEvYEXGN9HFRN2HCBFWe4LrrrkO/fv0wa9YsbN++HZGRkdi0aRM2b95smpY3QkZGhqs3QRAEQRB8ss6cLujST/yDhL/nI33fJsAYiLBbegHZWbkThwWaslkSyWgpCO4jMzNTWd2qVq2K9957T42LiIhA9erVce7cOVSuXBlt2rTBihUr/DIpCj/bEnaesqZlmrmF+hsuz37ZqFEjdZHT9TI0NBTDhg1TvRVZWVnKJfPll1/Gxo0bVY9GfHw8/v33X1eu3meR7JeCIAj+l83ScCQJpxJSlbvktmOJ6vfzKz9E0tblgMGA0K63IGLIAwiMjrVYRqQGdD+toUbZMJOFThCEorN//36MGzcOAwYMwM0336zGNW/eHLGxsVi2bJnfHuL8slx6QxbKTBtuob6E12W/PHr0KLp166YEnTn0O+Z4umYOGTIEhw8fxp49e1y9ekEQBEHwnWyWhxKVxS079RKMYbn1WoMr1UZoi06IePRRBNasY3MZSQagbadYDI2uUKzbLggljfXr12Pu3Ll48803UapUKeVF9vnnnytLnC7q/v77bwQHB8NfcaZsgbXFzpvcQv0Fl2cuoZuluaDTbwLdzZJZgaZPn64sdFOnTnX16gVBEATBq6107x4+bfqeHpmMcz++j+MfDkVmQq7Qi2x+Eyrd/AwiDoQhZNUJNYxKD8HBTk0thoFVynlwTwTBN/n5559VO1SHoUEfffQR/vzzT/W9Vq1aOHHiBF5//XXTNP4s6HzdLTTzijumP+ByUUfzNF0vdcqXL6/+8wbRYU9I586dVU+IIAiCIPiTle5sVhayL5zDpfcn4tzIPkjesRLBVetDy8owebYYcgBDtqaGQA0Y3Lo6IoxGiyHUKBmlBSG/eLi1a9eaBBuZOHGiKjXAGnJk0KBByuWya9eupvuP8XPCf+hxcnSxtIen3C8dWREz/UzYudz9kglSPvvsM2WZY88GM1yS33//HTVr1rSYlnU8/JnExEQ16A8emvwFQRCEklmeoM/V0armXMrShbj84TtAehoCGzRB5D0PodTlGghgdV4rJAmKIBSspNaWLVtUSQEaFVJSUnDjjTeiV69eWLRokSmhHzOv69a3ihUrqkFwjDO16LxJ0HmLW6hPi7q+ffuqG2f58uXo3bs3br31VpUlaPz48ejRowcqVKigekSYNIWf/ZkpU6ZgzJgxpu/yUBEEQSiZ5QlyMtPw4eFTyGgUperNBVaLQcTwhxHS/nplGcjZk4iAK9ksW1aPUiIwRpKgCIJDmOuPBoJ69eqp7wsXLsTAgQNV5nUm6itTpgw+/fRTlehEp23btnJUi5BRsjhr0Qkezn5p3luiW57ef/99jB49GuHh4cpXmTcgLXlPPfWUCkz1V8wtdd27d1fHa+/evZ7eLEEQBKGIFrquk9cpQUcxd3nLclz6ayHKjHkHwa1bq4YoWLj4Sk0rRWqWqj0XaDBg1ROdJZulINjh8uXLKoyH3H///Uq0nTlzRnWM8z/bnBR2rI0suCejpDdlmfSFrJw+m/1Sx9yV8NFHH1UvMaaI5YYyWcrQoUMtglD9kaioKDWQoKAgT2+OIAiCUMjkJ5P+PY7txy6i4rkMnExIRVZmBpK2r8TFPxYgO+kCjGUrQctKVdPTMsdyBRaEBUKLjcDEtnVE0AmCGWw/8p7h/1atWikDATNWEnqDMUEfQ1hIpUqV8MYbb8jxK4aMkp6sReftbqElzlJni5ycHJw6dUq5XUomIUukTp0gCIJvMmrDfixMz3WdDNyTiKzf/8C5H99D9sXTCIiIQplr+iOy1c0wBOXG8NzStCre7tcMcedTsGhLPI4npiE6KhR3topBw0q5FghB8Hfo1cVYOJbBevbZZ9W4UaNGKYFHa5zqHBGK3QLmrUIp02p7vXU7vdpSx2DUnj174rbbblMFxWmFs0dAQACqVatW1FUKgiAIglew78xlfHvxErTg3AZmVu1SCNgWBS0jFVGdh6JUq9sQEHylzE92bh9qzagwlb2ycaVSaNyjkSc3XxC8hqVLlyoPLmZGb9y4MWJiYlSGSrYddaZNm+bRbfQXnMko6W2CKcjMYke8bfuKgyLnQ6aIY4/JTTfdpCxw9GOeN28eEhISXLOFgiAIguClvLb1CNK2rseFB+5C2qrlQKgRAS0aIebhz1Hmmjv/E3Rm2Sz7t7afFlwQ/AXWhrvrrrssxh06dAhHjhxRn1nzmNa6p59+2kNb6J8UJKOkt5UKCLoi7PxR0LnEUsci4v/88w+WLFmielm+/fZbfPPNNyqmrn379sqKd/vtt5syEwmCIAiCL5cmGNAmVsW9rVi1CosffwLpe3cCwSHISbygps+qFQljfDJgFdwg5QkEf4XJ8WgAYDbKkSNHqnEsPcB247lz55RR4Oabb8bZs2ctLHOCUFCC/FDMuS2m7vz58/jhhx/w/fff4+eff1ZZiuj3TFFH32i6abKWnfhCWyIxdYIgCN5fmoDkJJ5Emc2zsWvTBrqrIOyWPoi4+z4YK1QyTcPYusBjyVKeQPBLTp8+rTr52d5jvWI2NVnQOzo6Gps3b1bTXLhwQWWxlERx3kdJyijpT1rArYlSmI1o9erVqidm2bJlOHr0qBJz5cqVwy233KKseHTbZOYif0dEnSAIgveWJjDPwpedchEnPrkfpW/oisC7R8JYNTrvzKlZCF9/Br88LuUJhJLP8ePHVf3hwYMHq7CcjRs3ol27dnjxxRdN2Sh37tyJOnXqqOyVgvfja4lSSgJeLeqs2bFjh3LTpMDbtGmTekEyC+bUqVNN5nh/RUSdIAiC50sTzIg/qz6PjK2Id3/6Fx+vO4jM88eQ+Ps8BFesiTLtB6jfM8pryGkd43B5/UMi8H57CT0QSh50kzx48CCuueYa9f2ll15SZaso5tq0aaNqFbOt17lzZ1PpJsH3KAkZJX0Jj2e/LAg0wXN4+eWXVWkDWvA4XLx4sTg3QxAEQRDyMO/kebxxKDfwPyLQiD3/HsC5H95D8u7VgJbDFM6maYMSDLjltKZKExDr8gR9W8WgdvkIOcpCieDSpUtITk5WLpSkW7duqh138uRJZb0eNGiQapDq+ROYV4EhN4JvU5IySmZeSeriq9vvEUvdypUrlUjr27evRQFywTFiqRMEQfCsle7av/biZHomcpIuQ/vsQyR8vxDZWVkIqd4MUR3vQWiMZfmBh66vg2d7NPTYNguCu0hLS1MJS+hNRfHG8gLDhw/HJ598on7/8ssvVbkB1pCTmLiSj68LoswrFkdvF6ZeZ6mbPHky/vrrL/Tr18/uNKmpqSpVrSRLEQRBELwho6Wxdikl6IghOATnf1+LmCZNkdn4TgRXb55nfilNIJQkcnJyTFknFy1apCxvCxcuVMntqlSpghEjRqBTp06m6e+++24Pbq334euiJz98eb8yrVxIvbHGnqtwuajbs2cPWrZs6XCaffv2qYxIkyZNwsMPP+zqTRAEQRAEpzNaZqcnI+HMKmjpKSj1yFMwBAej3NTPUS66Gp4KLo1Xvttlkf1SShMIJYlXX30VH3/8MeLi4lSHe6NGjZSACwsLM00zffp0+Dv2hJu1FUiH061fvx5//vknRo0ahZCQkGLeYiHTRrIXby2e7pWijvVGWJ9OZ8yYMdi/f78y1etQ9DGtLcseiKgTBEEQissyt+/UJaz756wqIZeTmYbLW5bh4sZF0FIuwVijNiIzMpSoM1ashJMZWcisGY5VT3RWQjA+IRUxZcNU8XDWqRMEX4MRN7SyUbzNnDlTjStdurSKhaObZc2aNZWoYyiNYF+45ZQrp0Twv5s2IWTaB4i5kFujsn2jRqgYEoL3mzVX082ZM0cJYtZrrl+/vhxSLynLkFlChZ3LY+rKly+vAmi/+uor9X3YsGH44osvVCYkc1jSgOltjzmoWO9PSEydIAhC8VjmtJxsXN66HJf++BrZyQkIqFQVkfc+gNDut8JgtOzrrBYShD+uaYQQKYgs+CAsJUVL3K233oo777xTjaOnFGPlWHJKsO+OeuLECRw+fBj/btmC7e+9j6PnzyE+MxPxOTk4k5GhBDJ5tEIFPFi+gvo84thRVAgMxJtVq6lMkUkvvoCkwEBcffXVFpZPwbN19nS8LZun18XUsZeHpQvMLXd65iT2BukwxS1T4gqCIAiCOy101sXDYQhA8u616mOZfo8iZMRgZZ2zxYn0TMw/eQFDo3MbbYL7KenxSe7k77//xueff66yjFeuXFkJic8++wxly5Y1ibo1a9aU/OQmFFyZKXlGZ546pX4LLFcaSenZKHWlXfrb7+vx5ezP8dSD96Nem9Y4cDgeDRrmTYIUZTQiJigIrSpUQExYGKqlp6NV2H9192bGVv9vXceOIXLceDShaBBBJxQDLhd1TJDCYpPMkNSlSxf18GBCFNYu6dq1q0WhShaoFARBEAR31Zqjy2VWdjZS9v6GtKM7Ub7HKPVOqtDzaRgjomAICQV+ze18VPN0rIXHu1m6SRlhkBNUTKLN2Sx1IvxyYRuLHemPPfaY+r5r1y588MEHqv11xx13oGLFijh9+jQqVapkOnYlXtCRzBSkjqmK/SeBIxdzcDQ1GwfPaNh9IBDHMjJxIjsDgUYDzj9TSk2+989MfPJTKq7+Yw1q9klFjVcO4uFhw1Dm9/WompyshByHSOus7maCzuZmlFA3P18qxZBpx1rnbVY6r3S/TExMRIsWLSzcKhkgumHDBmXqL1WqlHK7ZFwd/bb5WRD3S0EQBFcwK/4sXth/XH0eXy8aCyfPwaovpyLz7BEYgsNQbcSHCCxd0ea8TIDC+DmJlys8+YktR6LN2ULHzgi/kij6srKyVOINdkromSj79Omj8hNcuHABZcqUUW0wekHVrVvXLzKMswm7bt06ZSSgWyl5+umn8eWXc3Dq1Ok801OSVQkKQmxIIBrWzManfUKgpRqx6+dySEoKQHmjEcGR2ag2fQ5OPP9Kvu57zqJfy6bvJei69EU3zCAvFXRFdb90uagj9EN+5ZVXVCZMxs698MILaNeunQoqZfFxugckJSXhueeew7hx41y9ep9EYuoEQRBcV2suffOfyJj9IVL27IQhMASlrr4NpdvdAWPYf2EA5ugZLfu1jpXTUEjyE1uORJuzjS9nhJ+v1KTKD+Yi2LJlC2JjY1VZAdaFoxvltddei19++UVNs23bNlWK4KqrrjKVJCgpsHlKgcq4NvPh0KFDyljQu3dvNQ3DeRizpscI0vV05Q9LUfHIQVQzhCA2KAjRwYGICQxWgi7oitgNishCtWsv4MQf5ZCZbOk5Fli5MhBgQNbJUy7ZF16n1d6ehBNPPe3z16WvkelkZ5E34JWizha03A0YMECldiWdO3fG0qVLERkZWRyr93pE1AmCIBS+ztyANrFYk5ZistIljn0W6b+tRocBd+NYuRthCC+bZxls2nVpWAkNqpSSjJZubjg5Em16Yzc/NymSn/DzpQacNWyO/fvvv2jQoIH6vnz5cpXgZMqUKSb3yrlz56osim3atEFJgB38zMRJS1t8fLwqdWUu4JKTk/PMw3bjxIkT8cgjj5iOSdWqVXHDDTeo7+oauOceZMbnPgscYtAAzbZFMzAsSz0kslJshwoF8pri7ydyrcL2sHWN+9J1WRLI9JGOHp8RdeaxdFxlTExMca7W6xFRJwiCULhsliTz7EEkn/kVIY89pzJYZp86ocbH1qiBJ4NK2a01J5Y518XB2RNbxGEmOsYpWWXILmgD2pE4tCUw7e1LcXPx4kXlNkmefPJJTJ48WVmiatWqhcuXLyvxwlwFDGvxRTIzM1WnPveJIi0lJcUkUN98803lsbV161a1f/y9du3aSuDVqFFDfeZxsB4qVKhgci21PpfOZj10FnvCztlruyDXpeBeMr3ovi8xok6wjYg6QRAE5yx0XSevMwm0jHNHcfG3L5Hy7wb1vew7nyC4paUVY2L9GHQOCZNac4WkIHFwhbVm5IsTws/RNNaNcOKJBjVT5etukowBYwwck5uQVatWYdmyZXj88ceVqPEltm/frvbD2k2Sgo77rMO8ChSyFGU//fSTqmHMEB3mWKC7KTv+WcfYaJ2QxAaZJ04gbsgQ9bnGrNzi6HHDH3DOQucktlw07Vmhtawsi2u9IFZobxYagh+JOgbn0nQeHu44A5Aj2HOTlpaGcuXKwZthmuCpU6fiwIEDqveJLhLPPPMMBg4cWORli6gTBEHIn4k/7sPH6w4iM/EULq6fl1uaQMtByDWdEDHiEQTVyVvgV2rNFZ7CxMG5WpC5ShhaL6c4G9QUKzfddBN69uyJ8ePHq3G0UlHUMVulMyLGU1CEUaRVq1ZNZdFk+4fxbLSu6XkRmLBl8eLFpnlYwkq3tJlb3PiZbaeiJnGxdrF0FB9XWPdLLrPGDecQFJGDzOQAxK2uAJStgRpz5tiMF9XROw3yE3TeJOx8wYrlDzTxdJ06pssdOnQoZs6cWehl0C96zpw5KrOTN5OQkGB6kFHI8gF21113qc8cLwiCILgXxtCRpG0rkLxrNUJrtECpfiMRcEt7u/NIrTnXuVWap2h3lkCm0y+i66QrXDizzpyxmMZeunlXNHBpdaMFih3BTBTHRCeE7QUdulZ6A0zAwkR21lY2/TM778mnn36K++67T5VEYIIW3W2UPProoxg8eLBJvDGhi7uyb5quSzOLHIUcBV1BhJ0jIWgu6HK/56jv+N8Km9eF+biC3BvegK/EmwnFYKmjGwFF3axZswq9jGHDhuGLL75Qpndfg+4TDND95ptvirQcsdQJgiDksu/MZby0LQ4XUzPQOTgMg9tUR6SWouKNQq7ui882nUF26iVknjmC0BrNoNGb7UoD0ladOb3WXKixZGUH9KRbZUFiilyV5KQoyVaccc003wZbx8QRbAMtWbIE3333nWoX/frrr6qzlzV79aLfnoLNPF1gsbzUzz//jNGjRytRxuR17du3V9OYw+npBmluabv99tvRqlWrPMssTlK2b8fxxx+330HgyGJndg0ERVdTLpvqnJ86ZeG2aS3oLHjhBBAc4fT2OnsfeUpI+XJioZJIE0+7X/LhxXooen2QwvDbb7+pXiFfFHW9evVSDzZzt4PCIKJOEAQhNwnKUxsPIqNRlDocAZvikbp4LtK2LUVaSjJem/A2vrjUyCLpiY7UmStet8qCCLb81mm9Xa4si+CMFa8gGQoZI/buu+8qixUtVOSpp57CRx99pBpjNWvWNLVnCu1WyaZZZq5V2v5Gh6vODDbj6Elky9LG/xERESoZCaHr54svvqhKS7Vu3RqnTp3C//73Pwvxxv/Vq1dHSEgIvImUbdsQN/ju/M+lnTg4RyUFdHdOJMTZF3SFEHXeXCvNW7fLn2niDaLOFVAYFUXUbd68WfU+bdy4UQ30Xyf57V5qaiomTJiAr776CkePHlVxfT169MDYsWNVL5Ut6CbKOMAff/xRWRm//fZb3HzzzSgKIuoEQfB3mATlxinrkNKhMnK0dKR+Nx/JX30OLekygsrH4t23xuOhoYOwcMvxPNkvJZul8xQ5W2URSgm4onB4QSyMzsbj2ZuO+xE55V189+uvqhZa27Zt1XiWFSAsQUAYG8fcAsHBwXCJcMtIAd6ua/rKa/3f8zmgsbl++Vyh+OzFIVi55lccOnwEly5dyrMIunpSpDVu3BgLFy5U41g24Ny5c2jYsKGFK6i3w/N65K67kGWjoLg1uqWNxG29GjAEWFynahpb193Rg8DU1vYFXSFFnTdaxLzdguivNPG0qFu3bh1cBWvXFRa6OdD1wRpHu8fkLF26dFHuB3Sh7NixI44cOaJEIWMFOZ49VuawV4vT6j1wH374Ie6//34UFRF1giD4O0yCMu3QKWQ1jkLKsm9xefIbMFaNRumu96BUlWvx8A318WyPhiYBSKtefEIqYsqGSZ05K+w1Xl3tVlmYot+uiFlzJhZQ34789sVa0J3NysL65GTcVKoUwgICcJSdvX9sUFa59957T02zd+9eZc2iFSy/bbK7rRnJwPhqJtEWf0nD4cQcHE7IwaGEHBxO1PC/tkFoFxOIlEwNEeMv466rAjHvjtzEdIO+TcEf8dmoFRWQO/R+DrXqNjBZ2ypXruwRF0lXU5AyBdauk5nDtyiLplPXWgGso4Uhv84ItfhiEFEFtcSLsPOjRClFEWKu5Nprr0WzZs1UQU4OdH9g8K8j3njjDSXcOO/KlStNhdAZt8F6McOHD8fatWst5mF9FLotsH7MihUrMGrUKJQvXx533HGHW/dPEAShJMOMeisWzUV6+zYwIgphN/VU9eZCu94CQ5YB+O20EnA6tSpEmASe4Fzj0ZnGnHnyEA4Fcc/id11EOWoIuqKRmN/yrbfD3r5Q0F3IzMC/R4+h9ZUs3t9fuoh3zp5FpcBAtI+IQOz58/i4eQv0efBBU+ObKfgLcsy1nGwVwxVctSri4o5i/tff4MDOHTjyZwqOXMpG3EUNWTYMRF1qGtEuBggPMmBslxA0q/yfd9TcvmGWou35ZwtlRSqpgk6NY4IaZ48Jj6Ubj5+9+0MSlQiuosTWqaNbAUWdvd3LyMhQqXmZqnfLli1o2bKlxe/NmzfHjh07sGnTJuVyYY+RI0cqa6XuglFYxFInCEJJh9a1BX/TupaCmLLhGNAmFtXLhir399dee02Viwnr2Q+lR7+QZ97APYn4X50qIuSKMVauIG6VvoC+L4lxR3ApOwfVa9dWcVadr78eey5ewh/16iHIYMCxjAxsT0vFdRGRiLoSE+coJktf7sW4OBzPzMTJMqVR7tFH0fvmm9X497duxcwLF7DqqupofVMiNidk4NqZuRah8AADapczoHb5K9a2KANqlw1ArbIBqBkVgMjgAliFnjoABIe71LrkK4IuMDwLNW+0EQtXSHfJ4rLgeer+EvdL78TjljpfZf369UrQ1alTJ4+gI8xWRVG3dOlSh6KO5Q1mz57t5q0VBEHwbegqaR4Hxw63KTO+RNC2bxB/6F+Via/SiFHQ+tiu+5lduxR6t7Id5yw4V4JAtxQUxK3SfB7ii4KOoRYktGpVlP1gGho3bIjulSrhuyv78tSkSTgweTJysrJVwzw2OFgNOuaC7sihQ4jLyMCpG2/Epa5dceToUexbswbxycm4YJYXoMaoUWhxzTXKpZMWv2ZhobiclFvrrG7r81jYoByqZIWirNGI4MjsPBYm1kVDJoBgB/Fd1pjF4Llc2HgzBg3R7S/kHj9rYUuxVVQo6K64ydrFRQlU7JXZcDWOngW+3HHj7/itqNu+fbv6r6fntUYfT2HnCKYHpqunIAiCYN9CZ53YJP3YLpxeNA6GoBA8+OiTaDzyAYw7m8zQJptooUb8kZmOBiglh7mAPe/OCruiulV6C8xOqSdx++mnn1TM/YwZM1QdtfL16+PBESPQrHFj074MuO8+ZN56K47cMwQnjxxBZECAiqW7lJ2Nt1OS0fmaduhxJSvm0ydOYFvaFTfgvXvVv4rGQNQICkb78CDEBAchJigI1YOCTTF6d0RFqYFkJgNnf62ExrzSr7TAmKGRYi9PoWuKQ0eZGEs4znRCBIZlIbrDBYRXuFLnmILOBwRsQe7X4jy+Iuh8G78Vdcx0SWJiYmz+ro9nQU4dJlVh7ByzRrHnj4lZ5s2bp+rQFNS0as3BgweV1VAQBKGkQZdLCrq0Y7sQEByG4Mp1EBJ7FcrecB8iGnVGzZtbY2SjOhjeKFf0xZ1PwaIt8TiemIboqFD0bRWDGuXDVa05IS8FjZUrTGPOF8Qceeutt1RRb6byj4qKUu9cllxiMWzCLJEjRo1S6f4ZP29eAuDI4cNITUvDB9HR6BJZCpHVq+PbNauROWMmbrxSaPuuslG4Nbu0Em4cooOCEFrQLOBa3utYF3bWqfjNxZ4/4rATwlE9uRJ2v7qLkmCJF/7Db0VdUlKS+s8UxLbQM1oxIYp5nN3UqVNx7Ngx9TvTBNM987bbbiumrRYEQfBe0rJzMCP+rPo8MrYiQq40drds3oTTX7+LtMNbVLHwygPHqwQPpdv0Vr8zAYp5YfDGlUqhcQ/LRBT+jquz45WUxtx9992nhNrXX3+tvlO8MYnJL7/8grNnz6q4d5Y7WrVqlUpqduHChTzLYBkCetx06twZNSpXRrV9/yCodGnlcrnpsdEIvXLsye2lc8WdO1DC7peKFqLP2opXKFgewRY+Emtnulbvuce5AuFCoY6v/lnwXfxW1BWGKVOmqKEo2At+tGfBEwRB8BXmnTyPNw7lNoAjAo1om3gaL7/8Mr777jsGviC8cWdEdRiUZz6WJBDsk192vMLEypnPp3/2dlhS6PnnnzdluqaFjYnKWCOOmbirVauG+fPnKyHH8gPshL311lsRGxur0vvzPct0//qgF9vmfOY1d3UBTUoHBSGzOIWPPSvellao0WqLaVyBxEwJiLXLvVa/QNyQIeo7M4mqzJZ5JnRBDJ2bKez96u5tEnwfvxV1evkCFhG3RXJysvpfqpTEbwiCIDhjpZt69Izp+9ubduCfO7ojOzsb3W+5DXuq3gxjhRp55mPR8P6tY+UAO+mqZc8lq7CJD7yxMZeQkKCuGxbKZiwcBdf+/ftVmAL/22Lnzp0q2Y65BY+CjtY50rRpU/z6669Ord/8mOQb01WlsrJ2ZZ08ZTk+LCu39l2KVTPLoNkUbvmRefwEjmQ3V+syGI3/iRqrIuUlmaBq1VBjzhyvvW4LgiQqEdyB34o6Fg0lfGnYQh9fo0beRoggCIK/Y12ewFi7FOKPHUP2uTMIbtwM50uXwx3Pv4ynet6irCrW2S91QTexb1NVc04oenY8X0l8wJj0I0eOWMSzPfTQQzhx4gSWL1+uYuJGjx6Njh074oMPPkD9+vVVHDxdJG+88UYV125tbTMXdITWPA5uj+lqk5t0jS6Sehyc7hpoa7x1zJwJlk4wy55pi6yz50zTxI14KPe8li8Nf8KbruOi4iv3q+A7+K2oY3wcYY06W+jjXfFSsAddRjjohXeNV+rhCIIgeDPWAi0rJQEJR5cjZflCGCtVQfnPFqnC4ftvvhPNrs6NjaM1rk3NcmpextDR5ZLjRNC5NjueN8XKMf6csWzm4o3/T5q5N+q0b99eJR/74osvMGDAAPW9e/fuqgYsa8rSa8bcRdLjjW+rmC5+tpWx0mL8gu8QVL0Wapw6hbjhD/wXHxYTjWpvv40TTz9jP3mGlegzXQMzP0IQSgaujhv1BbzpfhV8Hyk+fvEitm7dqurNFab4eFFgsd0xY8aYvlesWBFnzvznviQIguCNFrquk9cpQZedloRLf32Ly1uWQstIQ2D9RogcMQrBra9ViVDIxPoxGBqd26gVnKMwRcJtLUNN44YGImsMnjt3Tgk0fm7Xrp0az/cZM0Lv2rULQUFBWLx4Mfr06WOaj0lM6P3SoEEDZXWbNm2aimf78ccflfcMBRw7ONmZ6ikB51Rso5aDGi0354lrU7XlbMS7mcaPizfFsNmKk7R73h1Y8dQ1QGFnHV/mjFumF8XU5Rc36lHcWHzcnwWtkBcpPl5ImO1q1KhRGDduHB555BGsXLnSlPGSaY4p6Bh47S5BR+heMnToUPWZPZJiqRMEwdvdLC+mZihBp2Vn4eSsUci+fA7G6rUROfxhhHS8wSTmdN6PO427qpYzZcIUiqcRV9TlMq6cLpK6dc3c0sZBzyDNcgG//fab+pyVlaXercwwSetahQoVVJIcluthkhO+5/jee+WVV9T0rVu3RnR0tKmcj7cnDDMllslMQdCsvDVu7SUvsTXeVpIam66e+bhlcroj945E9JR3EX7FA6mkxo16DD7T3Cx+vWZfBZ+mxFjqfvjhB4wdO9b0fePGjRY9iIRZ2Bg4be7Xf/311+Ovv/5C1apVlf8+69LxO61mf/75p/LX9wV1LgiC4G43Sy0rA5lnDiO4WgP1/dLmpUC18gi6p69K3mAPsdYV3CqRn7WuqLE3FGAU4OxMpIfIe++9hw4dOuCWW25Rv1NoUcSZw2lpUdPj2Pj/qquuwu23365+ZzgBa8OZe6LonjC0wPE7C4FTzPk0GcnA+GoFn89Jy5h+7rWsK8lWrhQxt4vRiMDKlVBz7tz/rgdnttELLHX2rnOJLRP8kSZF1AIlJqaO9WgoxqwxH8dpzAkNDcWaNWswYcIE5TJCV5Fy5cop6xkFor3C5IIgCP5godMFnZadiaSdq3Bx/VfQMtNQ7cGZMIZGovTVt0OjAW7NadN8IzvWwuPd6lssS4qGF182Sx12alKs2bO0Md7tjz/+UElsKLjGjx+vvFZ0UTdixAiVHVoXcBz4TgwM/K/ZwAyVuofJzTffrBoi7BilWKRVjiJPzzRNd8w33ngDfsNTB4Dg8EKl2ze34hGHrrhXLHkUfl5n4XJT3KggCCXcUufriKVOEARvYuKP+/DRmn+RvGctLv4+D1kXTyMgvAzKXNMfpVreDENgcJ55mM1y1ROdJfmJC60S1vNYT8s4tAMHDpjE2LJly/Dss88q10lbJXtKly5tEmmvvvqqsqLl5ORg3759Ks4tPDzcqZIDnTp1UiEKjIsjFGwUdO+//z7Cwkp43cFitoIVJNbOdH0wxs7NcWDFEjcaE21Zk85HCqYLQom21H311Veq1+7BBx9U31lrZtCgQdi7d69KUfzZZ5+pIGp/QrJfCoLgDfXlZsTnejGMjK1oin1jDN2ZhWOQdngLDCERiOp4D0q17omA4NwGO5tV5j2CUp7A9VaJjIwMxCUlYf/wYdj8+lgcS03FhbBQNPrwQ1O4AcMKvv76a6SmpirvE1rEaHljrJt5yn/9M9+z1nGPTErSuHFju9v++++/48knn1TCrVu3bsrFkjHo5nVcX3rpJZQoHCXHYBKSYqQgsXZOWbj0fXMk+rxFPCXEAVNbA3pcohe4jAoli8wSlKSm2EQdXRxZDFSHSUr4EmI9mqlTp+K5557D9OnT4U9MmTIlT/ZLQRCE4mTeyfN441DuSy3cGIC2iadVnBSTokQ0vh7BVeqhdNs+yt3SnIFtYxEVHizlCYog6HI0DVmahuCAAPXbC5064XSrVog7dUq5SLJeqrUzjeFoHG41G8f3Kl0f9eluuukmZb0rCnPnzsWXX36pSgww8QldKLlM1pFT22AwqJjzEg0FT37WOFsulgV0tSyosDsyaDCymCU7n+QpDoWdM/vmZvHkyL3YNI1V2QhB8Kusq94s6viC0s2KzIy1evVqldyEfvfMkvXEE0/A35Dsl4IgeNpKN/VobhmVjG2b8Ojo4UjbvUO54g1oE4tPm91oUSzc3Cp3f6c64mbpBHRV3L9jBzaePo24C+dxR5kolDEasTctDXcdjcOo8hVwX/nyatpvT5zA4QMHlDWNlrW2bdvmsbSxJEBISIhp+bScFRUmSWEnKztXycGDB7Fu3Tol5CjwWeKHJQwkQ7MVFHTFaDVig5NZLo8//nj+yVN8gHwLu4ugE/w566o3izq+DJhti1DQ8aXUpUsX9b1y5cpK6PkbdGHRM4XRZUYQBKG4ShNQtK1JS8HRbVuRNGsaMjb/BRgD0XngYOVaV7VCBCb0bWqR/dLn3CytXejU91TLaYLCLN3MdLczZ2pTBYYBWak4ffo0tmzdjkNMQnIkDoePHLnyP0652ZvTMiwMLcPCUSUwEO3Cw1H1yrOfsVDLZ89ClUaNTO8Fd8DtYThEo0aNVEwcoVWOVkHG4tEKxw5HftbFo54lU/A8LFvALJfOZEbVG67e3Eh1prC7P1OSXAO93XsiswQk5yk2Uccexw8++EDVpGGvIF1EdCHDjFwsKSAIgiC4vzQB+eT3Q7h4ZBEuzJupRExot1sRMeQBXKpdG+UqV1bT9G8dizY1y6n54xNSEVM2TI3zCUHnrJuZPbezK/PymJ24rCErB6hVNjfe8I1f0/HX8Wws3XxCFXhevCkDD/6QZloEJWJsGQOaRwWg1q0DUatOHdSqWRPVS5dG1NQPgVOnUTYwENNjYt2evp2C86efflJlB2gBZPKUhx56CEOGDDGJOmZ/ZhFwPdaOyVQE78WZzKjE293KzAWL2p977lExdCLo/js+7jiH/i4UM0tw1tViE3XvvPMOevbsiaZNmyI2NhYzZsww/bZgwQK0b9++uDZFEATBL0sTkKxL5xBYugIyqoUjp8zVCDl2GJFDH0ZgrdwC0CfSMzH/5AUMja6gvlPAPdujIbwCZ6xnBUzwwFi086kaDidoOJyYgwOvjsHh+JOIO3kCh7cmIS4xB5k5wB2NArGwf26c1L5zOfgzPlvVOg1lAe7KQRhbvRRigoPQrvNl1KkKBBv1bVie+y8u919mmwDEbWmFzOMn/svuN/MjBJUvnZtVsRD7YM758+exY8cOkycMS/UwQdnChQtxxx13KPH2448/qo5WnXr16hVqXYKXWbisBJ03u5XZEizMcsmkKD5joXPD88jdroElLYbMHVlXM31Y2BVZ1DGVct26dfOdjn75tMjxhcNacObZtyZNmqRcMAVBEATXQpdLCrqsS2dUnTnWm6t01zgYOndDSOi1CGlzbZ553o87jbuqljNlwvQaipjg4eTlHCzck4X2sUZcXc2oBF21yUk4lWQeNzhJ/Q0PC0OtSOCmuoGoFRWAa2P+cz/8rHcoAgMMLHaKzOQAhGyrgDvCcl+nQVvDYSh77r9sfVawwVqj1RbEJeaK5hotNyNoVivnapzZaBwmJSWpeDeWIyD33nsvVq5cqWL56EZ72223Yc6cOSarHOnRo4fjY+hv2GqcF3OGy6IKO2JL0HljI9WuYGHZAl8RdG5MOOMu18CSGEMmuFjUNWjQQGXfcjZzZfkrAeHm0Hrnj0hJA0EQ3M0/h4/iwqrpuLztRyA7CyGxV0GLLQ9DqP0YKWtrndf0etspq5qVoyH+koZDCTk4PPtzHNj7Lw7HxeVa23ZdxsaREaheJkC5UT66Ig1ju4QoUcfOxVvrBSp3yeoRRkTGl0HVnFDEBAWhco3qqNlqi02rgRJ0bCSdOo241RWQmfzfq5SfOc6RC5kSdjecM33Ow9t17TYO09PTVdFv1pNjyYMqVaqgY8eOyvpGHn74YfTq1cuUDZMhD3fffbfj4+rvFMZN10twtlC5hSigVdhDOBQstFjDv3GXa2BJjSFzS9bVWPe5w3u9qOOLgzVxbLFixQrl4kHLnJAXKWkgCII7mT9/PuY+MQKZ6akIrlpP1ZoLrdkSoChZlev+N7JjLTzerX6eeY1K6nhPw/pscg5O3bkMehfg4n2ZmLoxA4cTcnD0ooZsXe998YhpnsoVK6B22QCkZOb+2KhiAH4ZEo4mFf+zQM7oGaasbUqchfz3Ssw6fkJZ0+yJMzXPA49ZCLqCCjtnoJWVyWnI+vXrVbZLhjMwLo7lBh599FGVEVNHL0Iu+A9sfBbIrcxD4ilfwTLiIdSYueW/QuPFUCrCH1wDS3IMmbtiUoN89Fi4NaaOL5dhw4Zh5syZ7lyNzyIlDQRBcDWXL19WVhxmK2zRooWKlzrfoCdC6rT7z+1d6YlcsTC4dXVEuDOzoZNxJ9zuw6ezlUhjbBtj3JpWDsB9rYLVJHd+k4qdM3vjwqO5s1xI1fD38WyVvOS2+kYl3qqHGxF+vj6qXryE6KAglK4em+veeEVAhQcZcEMty9eeSdAVQJz9N0+uMLaFM8IuP6b+lYGX16Rh18ORiNE0NKlXE+2vaYfK5cuY4u/Gv/Zi7sT87i0Fo0si9lxifRFmbaXl1xEu3jenBQuFnQ83qr2Nkh5D5o6Y1CAfPgZuT5RiXThV+A8paSAIgqtKFPS6qgJ+WPAZJk6ciMmTJ6vshkxbv2vHNnyzOd5zpQmsLHAsuD1jS6YSb4eUeMvB4YyyOHfufJ5Ze9YNwn1Xws2GtwjCqVajkJ31ntr2e5sHYViLIJNQNQmttHTgSir+zPjjiEvIx9pmR9C5Upw5y/+Wp+JAQg5+HJx7TipFGNC4ohFnkzXEZKYi6sPGWNUJwJ6HcodiLhjt1xRzTbqCUtLdyvwBOYeeO96kJNwbxZb9UhAEQSh6sfAZ8WexI/4iVi0/gJxsDVp2JpK2r8RLfyxAdtIFlV2Yljodip7iKE2Qk5ODU6dO4fDhwyop1sCBA1XZmt/Xb8DgKZdVHNuQ5sHKqfOZn9NwMT13vmqlDGjQtD5uvqkmah37VlneqocGImh3eVQKCkRm8nklpu5tEQw88Sjw9vtqPt0lsbDWtqKgx8SZZ7HMO43jWlu0NPI4tKhixKi2udbIsykaTidpyMjWVPbMAVcFqUEQfN2tTASLZ86hHHfnY1J9XdAREXWCIAg+wryT5/HGodwaQ4HVwpG9cSvOLnkT2RdPIyAiCuW7PoCVn41Hw+i8ccyuKE1w+dIlHNi3K7e49uHcItuHDrPQ9hEciTuqUvzrsExNnTp1UK5cWZQPM5hS/FNkLhsUjgrhBtSMCkBooAF44efcmcb/8J9AywxETibyFWROW9u2Xo0aMz9GEF0Xr0C5VGPURcSNeFBZ9WyhSg7Mmp4nzkfNezYRcffem2deW4Ju95lsTN+cifuvDsJVlYwoHQJ8vTsTyZmaSdR92TfMlIRFEEqaW5k3i86SfA7luDumJF1zLhF1LFPAlMoVKhRjpjRBEAQ/s9JNPXoGWk4OkJ2NrNqlEPhvZRWzFtX5XpRqdTsCgkOxZOcZm6LOqXWkpamBruF6uZkTJ07g3XffVd/ffectvPr6OIt5QoxQ4uz6GANqRQWh9h2voVa9BqhYsaL6vXGjRtjyQKTFPNdVD3RaoFlY2lAEDAG57nORudulExRZETXmzClUQzMoOiLPvGqemR9h05E4rFz1C1558XklZE/+sgZTP7odde98FVclv6nE2+HHIlE+/L+kLUUWdI7S8Eu8nd/gzW5lNgWLrVqN6gf/jRF19Tn0ZrEvuA6DVsSgt4CAAFNMA7NcMoajcePGatATgcyaNctV21tiSxp0795dJTbYu3evpzdLEAQvZOaxM3hqzldImvUBQjregMghDyBwTyKMcZdgCPgv0cntzath6l0tbS6DqfAp0ugiqbtJmn/mb4888gimTZ2qYuG69rgNf2/egsQzx9VzfsOva7Fq3B2oFWVQbpJMTlIl0oAAs4ZX5vAtqjFmaiiwkeZELafMU6cQ16eHXYubsn4t+tFmA0TNO/wB+9Y2Jxov1gkFCtLgyThxAot69cLlrCwMXb5czcOslFOnTlXP9IYNGyqxHBcXh/r16sGQlWq1gBTbZQysE3XkN40jJN7OwQl07hr15pg6e9c0Mb+GbY0rjmLaebbt6EF1vxO7lnhvP+bFcLwcnq9CLs8bxb6QS5MmTdT/3bt3wyOijsJtx44d2L59uyp2alrwlYuY/2vXro2WLVuiVatW6j8HvRfXn3nttdcwZswY03cekzNnznh0mwRB8J7kJwPa5Ma9rVi1Cnf8v737gI6yzP44flMmgYSEQCihBggggnQQEEFERBQBRVfEglhQQVBsq7KIBZSyKqDCigKK+EdQFwEbClJcUaRLVUMndEJPgbT/uU+YMQlJSJnJzDvz/ZwzJ5nJG/KmDDO/uc9z79BnJHHL7yJBwRLa934pc98jIkmpEvy/w+KX5X/xgZ1ipH/zcvLTTz9Jq1atpHbt2ma2mc4D1VChs81ymx+q/0/rsTfeeKP0v+sf5knuscR0CQ/+e+nkpZhK2/qWpirmeMJQgCfMGgTzC2V25hX9mTNzD3Z5dHkrTDgr6BMeDcf6mFe2bFmzxFTVrllTAm02id2xw1yPjY2VhIQEadKkiXnx0ykz+kZXkyLz9CfI7lSCYcadLvn3XZLh9nyCpPyrunk3z72u/M26hLODIrwo1GW1b98+2bBhg3mws7/VV3/tX8LRTls3x1etakLe/PnzxVdRqQNgp01McnaozEg4LuEr35PNa34RCQiU0t1vlZB/3COSkiJphw5I2oE4Sd+yU9Lj4iT11BGpdt+bsuSf18u2Vcule/fuMmXKFHn44YfNv3X11VebZZX28GZ/q5ewsLDCP7m7xNJJR5jSJVUFCXUPDsy37Xa2fzOPJyPFqbZl/TfM52b5HH0M++uvv8zPSmfDbdq0yYS1p59+Wt544w1zzIIFC0w41r2EWR/rnKYIv5N8W/J7QUhBwRXovlHCoc4bq6OA14S63Jw9e9ZU8rKGvc2bN0tSUpJ54NNXPFH8XyQAa1foury13AS6jPQ0SUs8JYFlyktGaorEvXe/BISHSkZoGUmPPybpJy5u/a+VscDwivL2J/Nl4M1tTMV/8eLFJmDUqlXL5QEir2Yl9v1ltukX5hLkt/wy/rRT2rE7a3mRrjwpV66ceV9HRGiAW758uXTs2NF0+nz11VelW7du0rZtWykRxQ11OfGE2WcUuIpNqAPcyuNDXW70AVFf9dSQ16dPn5L+8h6JUAf4Bv0v9+jRo9n2tMWmVpClidUkJT5ODn78lGjf/5pDPzPHH9/0hZz59iPxL1deAqKqSUCVzIt/VFXzNrBKdbmuXDUZ3bqe80YUFCJAXKr7pCPY5egemf2gzKqRM5ZQFnV5kb7AqPua1T/+8Q/58ccfze9Jb9PHKt0b/uijj5p94x63RLAge/JyItT5hEsNn85233JhqLvoPkmlDnB6FnDLSAPdX6Abx/UCAN5I97F9//33FzUj2b17t9lrlVW9djfIsTNpkrBliUhGugSWq2oqdtoApWyDm+SOfwyQife2yfNrBYiflArwL/G9QgUaJ7Bvn1laWZBQ5qwObYU5VleNaLWtWbNmMmNG5ryi9u3bm/1yZ86cMUtWmzZtKhMnThS30t8RS9GKxkf2zBU20Dnun/3u+3uptAvPQ9GcA3Ad5tQBQCGlpKTI3r17zT7iTp06mdu0sjNs2DAZPXq0dO7c2VTkevXq5fgcHcQdHR1tlkTa97RVq1bN7Cue++U8SU9LleCaTaRcx3sluNrf1aAAW4jUiyovoReqSIU70cSCdZ6MP+0xG+dd3Y597dq1MnjwYHn22Weld+/eUrp0adOkyr7U0t4ADF6kgPcDQrPrg6WrAyTgywh1AJCDBrJDhw7l2vZf38bFxZll5EqrOWXKZM5h05B3/Phx836pUqVk1qxZpimUhjgNcPblffrv655ifTtp0iRp0qy5HK5/iwTVbHrR7yLA30/uaFXDZb8jMw7gwYGFD1DPbBdbUIhEO2GcQF7Bzv5+cXzxxRemYczcuXNNQxgdvfPHH3/I4cOHHcfo/kPA2+Q3dDrP5ZcurhQ6KoO6JNupXw0AoQ6AT8q6h2ry5Mnm+pAhQ8z1J554wswYy0mX4mlAs48K0IvdddddZ+a8ZdW3b99s10+dOmUGeWvzKA0ZGuy++eYbUyX6fG3cRd0vNdCN6d3YeXvlcls+mSWQFepVdO2kGBQqtpoxRR7enZ+ihrn//Oc/JpDbx8Xs379fVqxYIVu3bpU2bdqY39mxY8ccv3uvoEsHtdJUmD12DCoXXw92rhw+nd/ST8eS7GnrLr3XFoBnN0pBJkYaAK7d06Zz2fKqtukSuxEjRphjdYZbYmKi7LgwY0wDl1Zvcrb+z7pELzfJaekyNe6oeX9AjYoSfGE+me6he/fdd2XcuHGmkle/fn355ZdfTAv8nF0wdbRB3IkkqV6utKnQZQt0hd0blE8zgnw7Vs74SGwVyhb86zhpnEBhaZX0008/Nctab7jhBke4XrNmjcTHx0tgYKD52evb4OBg8RnF7ZTpbUsRacpRsPunk/YeFmQvX57nAPiwRlZslIJMEyZMuGj4OIDCtZ3XmWFXXHGFWVanFRjdJ6WhTSs0OV+z0sqYLoPU43VZpN2cOXNMYww7/Tf0UlizDsbLqJ2ZXd5CAwPk/moVzPLKkSNHmuV+Gj50rtm9995rgkZOGuCe69bA5XuD8mtwYl5Fv69/kZdMurIhgnajXLhwoXTp0kWqVKlilsAOHDhQevbs6Qh1U6dONf+X2n++oaFeFE6AIirQ/ZNmPIClUalzIyp1QN40kGloy1lp03CklTT14YcfygMPPGCajegTe21goqGtZs2a2Sps9oqb3u7sio1W1+as3id7TiTKoih/OXvhRewqwTZZ2fZyue+uu8x8s+HDh8tDDz1UvK9f2IpDLscXpGNlcV5FL8o4gbxoVVOrbl27djXXZ8+ebZa06niB+++/39GgRjtX5qx6+rTijj+gUue1nHn/vNTXccbcScCXNKJSZ126P0cv9s54gC/va/v222/NE/SsIe706cyujFlpt0JddqdNL7T1vA6GbtiwoeN+pIO3S4oulbTvg0utESoplcPl3NIfJGnxt5Ix8i2ZdfC4vP3226aRSkiIc/aHaChTttDMRi2epjhP0s6ePWv2w9Wtmxk6tEulLq/UIF+pUiVTofv4448dVTn7ckvkQMUFLrh/WmEvH+DLWH4JwKWhTTtF2oOaVlRatGjheEVKn6gvXbrUXP/hhx/MLDCdY6kVtZYtW160p00vlStXNssole5N04s7aIXOHujSJUPOHlknCW+9J6k7Y8UvJFRSd+2Qd0JKy11tL3fsrctVIfaxmE6VSyqYm6I7HytSsNPP0c+95MDwEnjSpfsez58/bwK6LqXU32+9evXMfkOllc1rrrnGUd2sUKGCWboKF8utiYoXznGDazlr7iSAgiHUASj2E/MNGzY4lkdmrbTpLLfU1FTHsS+++KIj1GkHSXulWj3zzDPy+OOPS40aNSxRudYllxrokvdslOMrZ0rK7m0iwaUkpM99Enpnf/EvGyEHzqXIpwePS/9qmUGsuLPkTKfKC0FMQ9klg13Oroj2mzUU5jGKwJVPujTka2jXUK5/MzqzT5vVPP/88+Z27TqadW+jzvvTC0pYbsszrbwkM4/7wUXHwPk/+hLYawsgE6EOwCWXxGmrfntFTCtq2slxzJgxZtmjLo1r27Ztts/R5YZadbnpppuyVduaNv17DtuMGZlzyOyqV6/ukb8J+565uBOJUr1ciPRpndmRUq+rpN3rJCUuVkr36iOh9zwoAZHZGx69veew9K1SPv9qXQFnyWUNYBruLhns8lmGl9soAlcGumnTppngruMF9O+mQYMG0qFDB9M8xk73HQJOx3JUt3Lm3EkAeSPUAT5Om4toRS23SptetOOg7gfTcKcVFt2z9t1338mjjz5qnpxrY5LXX3/dEdz0rTatsC+RtLKse+bs3v18kdSMWyTXD8gMIKG9+0nQ4/dLQFTu1bYCVesKMUsu+8cKEOzc9Cr6P//5T7OM8ueffzbXtVvlZZddZmb12Yezf//99075WgA8G2EOcD1CHeBD3SS16YQGtFtuucXcpsHsgw8+MPuZctIn4bq/SZtSaFjT8BcUFCR33HGH3HnnnY6W8dro5IUXXhBvk3XPnEo5tk9O/vyJJP65QjRe3dC9pwT4VxXbCX+xnRSRzX8v7wrw85OvH79aoiMzl3QFSNEC7t+dKvfnc4xzgp39/aLQWXA69y8mJsYspVQa3rShjTa7CQ8PN1VbvcADlyIWpCOmVThp1hoAWA2hDvAC+gQ6a6Uta7VNKyb9+/c3lTN94q1BzR7qLr/8cunVq1e2ZiT6vi6J0y6TudFg5wvse+ZSTh6SUytmScKWZSIZ6VK6bhuJ6HCPlLrsahndPjQz+KX9XckL8PeTMbdeIQ0rhZXcyZaLFhmyUCQqqkh7gwob5v766y8zf+/uu++WK6+80rxQ8PXXX5t9knaTJ092dDaFB/H2pYhOmuUIAFZDqAMsQDsE7t6923SF1GVrulxS9yfZg5vO88otfNWqVcs0p7DTSp12nLTTxhR6wcXse+aSd62ThM1LpFR0M4noeK8EV73swseTzKDw1rXKm2Waer16udJyR6vMPXfO4OhUub5lrssvS6qb3OrVq2XBggXy0ksvmQqtztjUUQ065FtDnb5gsHXrVilXrpzjcwh0AACUHEKdhwwf16VtPAnyXbr88eDBg9n2siUmJsrYsWPNx3UWmy5xXLlypbRp08Z0h/ziiy/MfjYdDZDboG1dPqkdBbNiptel6R5C/bmHNMychVamyfViq1BTStW4IstRGVI7PHO4t759rnON7P+IDv120hIvE+ymT8lslFJC855+/fVX8/fYu3dvc10D3ahRo6RHjx4mxOm4CQ1x2uzELmugAwAAJYtQ50YTJkyQV155xXFdX/WG91uyZImsWbMm21LJPXv2mNEAOStto0ePNsGsY8eO8uyzz5oGJCoqKkqSkpIc87tQfPoCy5tvvmnul9oU5rEn/SWgVCdJE1uOQCdSxv+8PLW6k8jqfIaDvxbntCVetqgol8170hcUNm7caCrA9pCmey21q+mtt95qqnAPPvigWbKrcwaVvgClS3cBAIBn8MvQ7glwe6Wua9eu5onStm3b+G1YlIYsXSKpIU27QtqXPmplQ2ezaUMS1bdvX5k9e7Z5X/et6XG5DdnWS9aZXXANDXDvvPOOjBs3ztwfNdiMHDnSVKm+WLf/ou6Xumfu3z1jpPfCK/NsbKKiv1xoxgZcUmGGjx88WOxOlfpffmxsrBkhoXvh9AUF/RvUQd/2v1HtShkaGmrmyOWs9sILeVNzEa2Ss6cOgAXpyiu1ZcuWIn0+lTo30sHL9uHLVhi27Os0oMXFxeXZ+l+Xq9lpSBg8eLAJ6lpxy/r71e6AOmRbQ1vlypW9ovW/J0tOS5epcUfN+wNqVLxoXpzOTRs2bJj5fUycONE0/wjQY1IS5Y4m5eXKaq1l7ro42X8ySapFlJbeLapLLe2BsjCvTpUXhoM/8IiZA3fJ4FWIxhVF7VR57NgxqVAhM2xOnTpVHn74YbOkUpdTalOc1157zcyMs9OOp/Ah3t48BQB8AJU6L0nnKD6tYGjVJiwss2vhokWLzL41rdpoc5FVq1aZ/Ww56bLZnJU2XS6pM7ngftPjjsqw2MwmI6PrV5d7KpWVjz76SBo3bmyGpuvv/b///a/07Nnz786eBXm1P59AV5JNTHKTmprqGDlx//33mwY5J06cMJVhrdJpkxOtzGUdBg94RRWPSh0Ai6JSBxRynlZerf/1ouFtx44d5ljdZ/T+++/Lvffea27XUQDaKTJrMxJdtlamTBl+Bx5cpXtn7xHzfkZamrz8/jQZ9X8fyM4dO8woh3nz5plK6e23317kr5FXoDMf27fPLJcsqWCn1WQNqtpAR7831b59e7OEUmfGaajTv2OtJAMFwogAALAEll/C62iVYvny5SZs2Str2vjhyy+/lCNHMp/gZ6UVDR0VoMfWr1/fcfsDDzwg/fr1cyxb0+5+2kQDnj0wXOfL6TiC6uVCJKBOmBxIPi/nfl4qZz+cLGm7d0jpsDB59dVXzcy+4sov0JVEsNPK/oABA8xySp1FqMt9da9c1rEVWpHTC+Czw9VzOwawMN1frUp6FQg8G6EOlqJL5bQrX84Km150X1SXLl3MMdp4Risxc+fONZ+nnf1iYmLMx3M2JNEnwfalalnRot1adFZc1qYmGX4i5yVK0o7GyamXnxEJDpaQvvdLTL+H5J9d2160t84KdMC3vrDwySefmA6oGt60uZLumbPTFy8An8X+QHg5ZzTMgnci1MHjnD592lTZdMnY9u3bTfMKe4jT7pLaZTKn8PBwR6MSbUqizSB0mZkdFTbvr9BlDXTJ+zZLarkACSxVTQJrREv4c69KUKu2ElC+gmit9tODx6V/tcwKrFOGg+dTrSvOvrpp06aZ0PbGG2+Y6xretLGLVug01Ol+Tr2NGZfwWVbb8wc4IdDZR9uU5PJ+eD4apXgIX2qUcv78edNG3V5h08qaLo9UOmB7zJgx5nbdr7Z+/XozDkAbWGiXvtza/uulfPnydJH0YWO++0PeW75Dzh2MlZP/mynJu9ZJQPVoiZzxZa5/F1WDbfJr28vzrtYVpNnCM9tFgjKXcaUcOmS6XabEZTZkKUqg0xcrtAKnQU1nwimdE/fdd99JfHy8GTGgx+iLHcwnRInx9MYjnn5+gIsCnbsbcsH5aJQCj6UDtrXKkHOppI4FyDoesUaNGo5Q17p1a7OPLesf+L59+6Rq1arMy0Keft+0SY58OV6S/vpVX6uS0q2uk5DHh+QZ9A+cSyl+tU4D3YUnijqPTscXFGY4uAa1hQsXmkYm+gKGhjVtxNOuXTtHqNMq9ccff2wCndJGJwAA35JXoHNHQy54Lip1HsJqlTptj26vpmmV7NSpU3LnnXeaVv5abVPXXnutLFu2LNsetZyt//V9vWRdKgkUtqNpZKUoOZd4VkrXbycRV98ttsq1si23GtChtjx5/d9NcFSA+EmpAH+nvvqf314Hvc/88ssv0r17d3NdK3A33XSTjB8/3tG0RRv86P8F9uY8gNt5eiXM088PcGGgy4qKnfVRqYNLJCcnm/1rOZuR2KtuGuLU559/btrB6x44fcKqIS/rkO0hQ4Y4QlzZsmX5bcEptHq7adMmE4pCQ0Jk9Nix8u6asxIUFZN5QHqy49gAfz+5q7pNgo4cKfirmEXsoJd1OPj58HDZ+eefjnmF+mLHlClTzMgMvU/osO8ZM2bI9ddf7/j8a665pmDnBwAAkAWVOg/hSZW677//Xrp163bR7bo8TDtFZq229e7d23HuuqQyr+VugDNo59PRo0fLe++9Zzqa6t7MsqUD832l3owdWN9SxM/fpctTdK+o7nmzv3ih95GQkBDZvHmzuf7rr7/Khg0bpE+fPqa6DViCpzcioVIHH3Cpah1VOu9ApQ5OV7duXbnjjjsuWiape9+0YUleCHRwFV26+O9//9vsMUtMTJSWLVvKqFGjTNfT/J5w/j1HLrOBiTP3Heigb/2bt3dpbdq0qTz++OMmdKrBgwdne7FD98rpBbAURgQAbn8hxL4KhEYpyA8jDdzo5MmT5qJSUlI8pi25znObM2eOu08DcISiTp06ycaNG6Vhw4YycuRI0xXyUi8i5DYY3FkbymfPni2DBg2Sb775xgQ1feGjc+fOjqWW6qmnnuI3CAC+RgOdC/Z55hbsqNAhK+tN3/UiOjvNXgmLjY013fAAZLb3/+GHH8yPQsPbiBEjZObMmSbY6ZLfogS6nMFOl7MU1EsvvSTNmjWT9PR0c12XIWtzH917qvQFma+++kr69+/Prw8A4BL2YKdhjkCHnNhT5yGVuq5du5onhjoCAPCVgeFzVu+TuBOJUr1ciPRpXUOqhdtk+vTpZmnloUOHzP3hkp1Rc+ypyS/QZZXXA6JWzR977DGpVKmSOQ/17LPPmvEDixcvlsqVKxfn2wa8g6fstfOU8wBKcJ+n/UVJRhh4l0bF7K/B8ks3ioiIMBdls9nceSpAifpszT55Ye4mSUvPnFeYkZ4mb03+QPzWfS5HDuwzA7h1D53u43Q17eiqe/V69eplxnDofVFHC1TJEvbGjh1rzgeAa5eYFRp7/uCDCHPIDaEOQIlX6LIGOhX/zXhJ2LpM/IND5ekXRsjLw541YzKKwhaaLtGdj+VZrduanCw/+vnJ6Pf+Yx4Yz//5pwl1wcHBJtSpVatWZRvBoc1QAAAAPBWhDoBLJaely9S4o+b9ATUqmiWXqWnpcv7AnxJcrYG5PbTJ9RJQtpKEX9lbKnZsWuRAl1uwW3ssRfalpEh37ZSpgc1mk/9sj5Vbd+6U62NipH79+rJ161Zp0CDzXBQzFQEAVsYSTd9DqAPgUrMOxsuonZnr/0MDA+S3X/4nh/9vgpzbv1Uq931dStVsIqWjm5qLijuRVLgvkGVQuDYy2bR5i2mk0qTxFRJ96JD0bdpKtp46LV3DwiSkZk0ZNG6s/CMtTZo3b24+R4+9/PLLnf1tA4D3YQ+jpebaKVfOZ4VnIdQBcGmV7p29R8z7Kds2yxPPDZIzq38V8Q+QiGbXS5lykWKTzA6SdtUjShVq3IHOiNOmJlpdOxEfL82vvEpuv/12+eyzz8RWM0ZGT50qR0aPkeCICMeDW7TTv1MA8AGespcSBR5U7sz5rPBshDoALutmuTQ5UQ6eS5GEzz6Ws++NN00N2rWtIzNbH5GY8r+JiF6y291se75f5+jRo6aRitLgduedd8onn3wid999t0RGRsobb7whLVq0cBzf5bbbJOWqq8z7PKgBAFwqy+qRfI8pgUDnzPms8HyEOgBO72apJi34WUrf3kK7jEhw2w6SsnWTlOn/qCRHhUn1VXeJZKTk+m/Vigy5aMSAvTvs448/LpMmTZJjx45JuXLlpEOHDmb8wKWGfvNABgAoEW7qyJpboHN8jGDnEwh1AJzazTL11BE5ueJTSdjyo4SXfUlK39BTqpQJF3n8eTlaPlL09ctPo26S/gfnX3Jp5TXXXCOBgYGyZMkSc9tVV10lZ8+eNRcNdVWrVpV3332X3yAAOBmNNqwjv0DnOIZg5/UIdQCKRZdcaqBLO3tCTq38TM5s+E4kLVVsV7SQwJq1peLxeBk//lVz7JNPjjDB7u2a90jfQ99KcJZq3c4T6dLvyyTpU/Y/MuTJZ0wDk1q1ajkCnl7XpZZ6AeC7S8zgejTaAKyH4UsAikX30CXtWCP7pzwkZ9Z+JUGVakvkY+Ok3MSpUrVyVRPoqh07Yi76voa8A6UqyfDjzaTrzATZfjzd/DuVQ/3kz/h0OXnqlOPf/vjjj81FAx0AD1xilt+F+63lqz726o69agfPZBqA6Z65GjXyPqZGDfbVeTkqdQCK5PTp0+Lv5ye1w/2kbNVoOVO5lpS/speUqddaJMBfIr9aL6N+miZRCcfN8QtOnZINhw/J3NEDpdJ1Z+S7Qwnyzp402XIkTeqW95fQID85/EwZ8R/2PL8RwGpode8VLtloIzJz3ic8N9jltgyTQOcb/DJ0XRPcrlGjRubtli1b3H0qQL6SkpJMs5IxY8bIkEGPyEv+F+9pSzjtJxP/GyRBKTbH0O9hBw/I16dPy/KYulKprEhUx6MSGJoupQJzVOFohw1Yz/kEWt17+b4sEwym/Uds0//uLpwr/g/3qN8jgc53sgCVOgAFcv78eZk6daqMGjVKDh48KDVq1JCYOnVEdoucSMqQhdtTpUUVf6lT2ib7l0XKv/fvlmqBf4e6JypUlGGVKkuZgABJSRA59FNFie58TCQwc/klAMDDG208OFCip60TW1TUpfdSUr11i6wVO8UoA99BqANwSX/88YfceOONsnv3bjPoe/To0RITEyO39egmMuYp2Xo0Te6amySvdAqWF1rZzLLMSdWqSbQtyPFvVL4wlsChXLTIkIUiOZ8c0FwBADyX/h9dkJb9DCp3e7Czvw/fQKMUALlKT083owNUVFSUlCpVyoS5nTt3yokTJ+SOO+6Q3zduMh+/slqAfNirlDzY3Ca20HRTgbuyQpCUD8z9dSOzHGTmTLHVjKG5AgC4GY02vPN3SqDzLYQ6ANnoNtt58+ZJ06ZNZciQIeY2nRenAe+5556T0NBQ6du3r0yePFmqVc18BdAW4Cf9mwVJtfDM/1Lswc4WmnrRT5f1/QBgrWDH/9uA52P5JQBTldMw99NPP8mzzz4ra9eulYCAAOnevbu5feDAgSbUpaamis1mk2bNmpmLaY6QB3uw27O+paTE7c+8jZbKAGCpDor8vw1YA6HOjU6ePGkuKiUlxTyJBkraggULpH///hIdHS0bNmyQwMBAs19uwIABpjKnHn300SL92ybYTZ9iNtcrNmwDgGej0QZgTYQ6N5owYYK88sorjusVK1Z05+nAR4wdO1amT58umzZl7oerVauW1K1bVxITE024GzFihNSuXdtpX0+7pLFhGwCsg0YbgPUQ6txo6NCh5km06tq1K5U6OJ0unRw8eLAEBQXJ+PHjHaMJdLnlrbfeKsnJybJ48WJZtWqVqRbr0spCd0HTmUSXOMZWpQCd0gBYVwH/L4B10GQDsBZCnRtFRESYiyr0k2kgF/v375c333xTrrvuOrMfzs/PT1avXm1CnNKRBDt27DAdLLdv326O0QqdNj8p0t+gn1/BWlsD8G78XwAAbkWoAyxMl1DOnj1bnnrqKYmMjBR/f39TkdMKnAY2pZW4sLAw08lyypQppiJ37bXXmiHiV111lbu/BQCAN6J6C5QoRhoAFqJdKWfMyBwoqn799Vd5/fXXZdmyZfLHkTMyZFO8tJs0V6rcPFh2HUswFbrw8HBTsYuPj5cWLVqYkLdkyRICHQDA9dXb/C56DACn8MvQTTdwu0aNGpm3W7ZscfepwENoINu8ebPZA9eqVStzW48ePWThwoVm+HeZMmXk6NGjmUsq0yvKP9fulvOXZy7n9V8bJ0nzZkn545vlz43rTEfLhIQECQkJMQEPAAAA3pMFqNQBHkJfX9F9blpRU7qEsnXr1jJ8+HDHMS+88IKptAUHBzs6plao3VCGfbVVztcOk4zkJEmY/ZEcevUeObniU9l3+Lis2PCHOVb3zRHoAAAAvA976gA3Onz4sFSuXNm8/91335l9cJMmTZJBgwaZqpruj6tfv77j+Iv2wGVkyNyVf0lGFZsk/jhfkqZNktSTx8W/chUp3/MhqVD1Sll5NECuKelvDAAAACWGSh1QgnQppd2wYcMkKipK9u3b5whsAwcOlObNmzuO0XDXpUuXvP/BlEQZvPZ6iYjJkOCdsRKecFYG16krDSZMk5p3d5XfQwfJkeMnXPtNAYAr6O6Q8wn5X9hBAgAGlTqghNx4441y5MgR0+xEtWvXTu677z5H0NPxFpMnT77o87ThyZzV+yTuRKJULxcifVrXkOjypeWzzz6TCePfknv7Xi+pCQEy7cBeqV6rlpT295de746RJ58cIZ9G3STVIkrzOwZgPSmJIq9Xzf8YnY3HWBUAINQBrnDgwAHp06eP3HDDDY49cXXq1JGyZctKWlqaGTSvTU/0kp/P1uyTF+ZukrT0DMe+uwnTZ0nw71/I3thtZlTBtGPdZPwXr0q1k8dF/DOL79WOHZHx41+VMYMGyeRm1fglAwAAeDGWXwJO8NNPP5kAt379enO9UqVKsmvXLjl79qzjGN0rpzPlNNBdUkaG7D54VF6du1qC0pOktCRLxr51cmTmU3L4i5Gyb9cOeXjwUBn33/kybtlyE+Jy0tuenzxZ1hzYz+8YAOASKQcPmgsA92L5JVAEc+bMMY1NPvzwQ9NRUgd6L1++XLZt22b2xOkIAd0rV+RukymJUmtKXdkc9PdN758+L4MPJ8vg1jYZ1iFIvmzeXa4b+5qk5hLosga76BeelZSZM8VWpUrRzgUAgNweqg4elD397jPvR388g8cZwI2o1AEFmBc3bdo0mT59uuO2X375RWbOnGlmxKmOHTua2XF33XWX45jijg9YdzBNen6aKDuOp5vr9zezSeyQMvLOTaWlSpi/HDqdLAX5CkylAwC4KtCl7NtnLuZ9KnaA2xDqgBxOnz4tn376qfz++++ZdxJ/fxk1apSMGTMmW+fK48ePS+3atc11m80mpUs7pyHJ1q1b5fY775aW7yfIV3+lyvc7UjO/RoCfREf8fZeNiK4m0dOniC008+O50Y+ZY6jSAQBcEOgctxHsALci1MHnnTp1ShYsWCDnzp0zP4u9e/eaitvHH3/s+Nlop8mff/7ZcV1ny2nTE2fSqp92w2zcuLH898v50qN+oGx4JFQGtc6yBjOL3i2qiy0qSqI7H8s12JlApx+LivL53zEAwHWBzvExgh3gNoQ6+JzExERHFc7ewKRXr16ycuVKc71Ro0YyY8YMefLJJx3HtG7d2jQ/caU///zTBMlrr71Wfv1piSzoGyJNo/JuqlIrMkTEFiK21+Ik+suFYqv+d5dLfd/c9lqcOQYAAFcGOscxBDvALWiUAq+nTUy0GlehQgVzXbtUbt68WY4dO2Y6Ufbu3dtU3erXr+/YC9evXz+Xn9fRo0fNks5bb71Vrr76aunatausWrXKBEgzVPfHAvwjum8vKFRsNWMkeuZMNqwD8B76gpTOobvUMQAAQh28s7GJXrQDpe6Pq169uvTs2VM++eQT8/EBAwbI4cOHzdBv3QfXoEEDcykpJ0+elDfffFMmTJhgRh6cOXPGhDoNkybQFZHum9PuY/b3AcDSLrxoBc9hf5zJr1pnq1GDTpiAG7D8El7lxx9/NMsk586da66Hh4ebStyVV17pOEarcM8++6zTGpsUVEJCgowePdoMIdfGKxo2da/ee++959QHXAIdAMDVwU7D20UfI9ABbsPyS1ja22+/LePHj5eNGzdKWFiYCUy1atUyVTq7jz76SDzBnj175F//+pdER0ebc7777ruznWc2LDsCAFioYkegA9yLUAdLeeqpp8ySyqlTpzpuCwkJMYO+GzZsaEYMrFmzRjxlL582XImMjDT75vT8fvjhBzPTLigo946WDiw7AgBYJNgpho8D7uWXkZGR4eZzwIWOi2rLli38PC44cuSIjB07Vtq0aSN33HGHua1Lly6yf/9+83PS+XH651vcId/Opvv5Zs+eLS+99JJs377dLP387bff3H1aAAA4nX3gOEv/AfdmAfbUFYLuf+revbtUqVLFdEvUikvW2WUonm3btsmLL74ocXFx5nqpUqVk4sSJ8tVXXzmO+fLLL81wbg10ypMCnQbMefPmSdOmTc3SSg2lr776qixevNjdpwYAgEuwlxvwDIS6QtBuhdoWX+eaff7551KtWjW57rrrss08Q8Ft2LBBpkyZ4riu++K0gciiRYscTU7ss9vsdN+cJwW5rI4fPy733nuv7Ny5U55//nnZtWuXCal6zgAAAICrsPyyEOLj483+qKzL7Bo3bizt27eX999/v1i/CG9ffqlVLP3etH1/u3btzG1azZo1a5YZL6AdK3WWnIa4Fi1a5N1AxMP873//k71795rvRX377bfm/KOiotx9agAAALCI4mYBQl0x9enTx1Ro7NWlovK2UKchTitVoaGhUrlyZUlLSzOB+PLLL5dff/3VHLN27VrT9ERD8SUbh3gYbcYyfPhw+f777833pUtGdbkoAAAAUFjsqbtAA8KYMWPMTDKd/6VL9AqyTC8pKUlGjBgh9evXN0/Kq1atKg888IBpxnEpGlRWr14tdevWLfQvzhsdPHjQhDm1YsUKiYmJkenTp5vrAQEB8tZbb5nmIXYtW7aUa6+91lKBTkO3/o3pkHDtZHnXXXeZkEqgAwAAgLtYY41bAYwcOVLmz59fqM9JTk6Wzp07y8qVK03zk169esnu3bvlww8/lK+//trcrnPP8vLuu++apXeDBg0SX3T+/HlHIHv99dfNDLbNmzebVxpatWplfi5ahbPTsGxlGuJvvvlm8zdyyy23mCYouvwWAAAAcCevCXW6T6tJkyamgqIXHUB97ty5fD9Hm3JocNPP1apLmTJlzO1aUXr66adNCFm2bFmun6st6rUZhi7B88Un9rfffrtpEBMbG2uu68+wX79+jq6UWrnShjJW8seRMzJu6Vo5k5wiTStVl7tb1xRb8glZunSpaYCi1cb33ntPypcvb/7GAAAAAE/gtXvqNFRoqMvr29Mqk705x7p166R58+bZPq5t6bUbo+6d0mWCWWmlpm3btmakwZw5c5zSjdGT99SdOHHCVDH1ex43bpy5TQPtjh07TGfK0qVLi2Xp30dKovx3XZyMW75RXvvqXXPzM53ukX0/LJAzv/8g6WlpZt6cvlAAAAAAOBt76opI93xpoNN9XzkDnb0SpbLOSFMnT540s+r0Cf6MGTM8tr1+caxatUq6desmP/30k7keEREhhw4dMstV7XT/oo51sHSgUymJIq9XlRvmdpSxCyZKmcMH5fNtW+TP0Y/K8dVfS5PIVJn28f9JdHS0u88UAAAA8O7ll4Vlny2n7edzY79dq3VZq3vaJCMxMVGWLFli/UCTBw2quuTwtttuM9VIva6jBrwxwKqUBH+J/amKVD4VLwMOHJAViQkSExQkj9asIo/3Oinvlq3vtd87AAAArM9nQ502OFHaKTM39tv37NnjuE0bfyxfvlw++OAD065fLyo4ODjXal9+pdWcdCmjVg09gS431Ypk1tDqraHmzO49Mv2/paRrUOb+y0cjI6VHeLh0Dw+XAD8/2bMsSE7W0U6oBfv9AgAAACXNZ0Pd2bNnzduQkJBcP67z1ZQOy7ZbvHixGTj+4IMPZjtWl+bpPjtvoc1OvLUKaZeSkiIfjB8vr44YIYfPnZOPatjkypBQaRkSIll3UKaeDZBe//dvSbm1pdiqVHHjGQMAAAC589lQVxTOCG55NULJq4IH548lmDVrlrz88suyc+dOKW+zyfMVK0nTUnmH2JAg7iYAAADwXJn9532QfXyB7o/LTUJCgnkbFhZWoucF1/rnP/9pRi8cP35cXnvtNRm3aKlcV7+BBF8YxZDT/gqVZP24N6jSAQAAwGP5bAmiZs2a5m1cXFyuH7ffTtdDa9ORFj/++KN06tRJAgMD5eGHHzbjLp555hkpV66cJCefkXPt/5IjP4aZpZZZBZZJk9bt/5LgxnXddv4AAADApfhspU7n0CmdUZcb++060NxVtBmJLunUi+7x0v16cB5tatOhQwe5/vrrZebMmea2yy67zFToNNCpUv7+UrZ0otS69qjYQlMdn6vv6236MT0GAAAA8FQ+W6lr3769lC1b1nSd3LBhgzRr1izbx7/44gvztkePHi47hwkTJsgrr7ziuF6xYkWXfS1fonP2hg8fLosWLXJU5zTY5coWIjLsgNi0KnvokOx54BFzc/T0KWKLivr7GAAAAMBD+WwJIigoSAYPHmzef+yxxxx76NRbb71l5tNdc801pr2/qwwdOtQxGqFevXoSGRnpsq/lK3QgfJs2bUyn0nvvvVf++OMPmTJlSp6jK0RHNQSFmoutZoxEz5xpLvq+/XZzDAAAAOCh/DJ005EX+Oabb2TkyJHZqjX6rekTfLsXX3xRunfv7rienJxs9lr99ttvUqVKFbNUT+fS6XWtmq1cuVLq1KlTIudv736ZV3dM5C02NlYqVapkKq/x8fHyxBNPyAsvvEBHUQAAAFhCcbOA11Tqjh49asKY/WLPqllv02Oy0oYZS5cuNWFP59XNmzfPhLr+/fubPXUlFehQ9AHyDz30kFx++eUyfvx4c5tWOz/55BMCHQAAAHyG11TqrI5KXcEdOnRIXn/9dbOs8vz589KxY0dzXfdJAgAAAL6WBXy2UYon0O6XelHa/TIgIHtLfVxMK66dO3c28wVbtWplOllqExQ/9r0BAADAR3nN8ksr0u6XtWvXNhfdF6b7wXCxM2fOyJEjR8z7zZs3N6Huyy+/NPsmu3btSqADAACAT2P5pYdU6jScaKVu27Zt7jwlj5KUlCSTJk2SMWPGyM033ywfffSRu08JAAAAcDqWX1pYRESEuSibTSelQek+uWnTpplupgcPHpQaNWqY8RIAAAAALsaeOngUDXFXXXWV7N69WypXriwTJ06URx55RIKDg919agAAAIBHItTB7dLT0+X48eNSoUIFiYqKksaNG5sgN2TIEAkNDXX36QEAAAAejVAHt9FpGjo0fvjw4WZw+LJly0zTkwULFvBbAQAAAAqI7pdupE1SdJmhXnSkgVasfMWSJUvMXLkePXrIn3/+Ka1bt5bU1FR3nxYAAABgOVTq3DzS4JVXXnFcr1ixong7Da833XSTLF68WAIDA2XgwIHyr3/9S6pVq+buUwMAAAAsiZEGbuRLIw0SEhIc++P69esn/v7+8tJLL5kZfQAAAIAva9SokXm7ZcuWIn0+lTo38oWRBrq0csSIEbJ+/XrzR6rfp86b01AHAAAAoPh4Zg2X0H2C999/vzRs2FA+++wzqVevnpw4cSLzj45ABwAAADgNoQ5O9/zzz0v9+vVNRa5jx46yYsUK0+WyUqVK/LQBAAAAJyPUwSmydq7UvYHNmzeXRYsWmS6XOkwcAAAAgGsQ6tzIG0YanD59Wl5++WWJiYlxLK/UBigrV66ULl26mLlzAAAAAFyHUOfmkQba/VEvsbGxEh8fL1aRmJgo48aNM+euYxlKlSole/fuNR8LCgoizAEAAAAlhFDnRkOHDpVdu3aZizYSiYyMFCv45JNPTGXuueeek7CwMJk+fbrpbNm0aVN3nxoAwJNlZIicT8j/oscAAAqFkQZuZKWRBhkZGY7qm1bp1LvvvisPPfSQBAcHu/nsAACWkJIo8nrV/I8ZdkAkKHOuKQCgYAh1ENG9fImZSz9/3XVM/rNspxw7nSwVwkvJIx1ryb71K2TkvyfK/PnzpW7dumZUwT333CMhISH89AAAAAA3I9QhM9C9Udf8JNpduGhl7qtfUmXQa+dk4+F0CQsrI5s2bTKhTquKnl5ZBAAAAHwFoQ4XWRmXKk8sTJZV+9OldKDIs1cFSecxP0i3Dhr3AAAAAHgSGqXgImfOiaw/mC6PtbbJjsfLyLjrS8mnvx/nJwUAAAB4ICp1uEiXOgGye2gZqRr2d+Y/evocPykAAADAAxHq3Dx8XC9Kh48HBASIJ9Aul1XDsg8NrxhOh0sAQAk4n9lhORtbiD448eMHgDwQ6tw8fFwHd9tVrFhRPNWj19Rx9ykAAHzBhcZd2TDmAADyRahz8/Dx/v37m/e7du3qMZW63LSrXcHdpwAAsDqtuGlAy1mZyy3IAQAKjFDnRlYaPg4AQLHpEkoGi1tWysGD5q2tShV3nwqAHAh1EAmJFHlm+6WPAQAAPhvo9vS7z7wf/fEMgh3gYQh1EPH3Fynjufv5AACA+wNdyr595rq+T7ADPAtz6gAAAFCgQGdu27cv87YLyzEBuB+hDgAAAAUKdI6PEewAj0KoAwAAQIEDneMYgh3gMdhTBwAAPGvMQW7HAADyRKUOAAC4f8xBfhc9BiVKxxaYZig1auR9TI0aNEwBPAShzo1Onjwpu3fvNpeUlBRJT0935+kAAAAUKNgR6ADPQqhzowkTJkjt2rXNJTY2VuLj4915OgAAAJcMdgQ6wPP4ZWRkZLj7JHy5UqcX1bVrVwkICJBt27a5+7QAAACyYfg44FqNGjUyb7ds2VKkz6dRihtFRESYi7LZbO48FQAAgEtW7OzvA/AshDoAAABcEmEO8FzsqQMAAAAACyPUAQAAAICFEeoAAAAAwMIIdQAAAABgYYQ6AAAAALAwQh0AAAAAWBihDgAAAAAsjFAHAAAAABZGqAMAAAAACyPUAQAAAICFBbr7BHzZyZMnzUWlpKRIQECAu08JAAAAgMVQqXOjCRMmSO3atc0lNjZW4uPj3Xk6AAAAACzILyMjI8PdJ+Grslbqunbtaip127Ztc/dpAQAAAChBjRo1Mm+3bNlSpM9n+aUbRUREmIuy2WzuPBUAAAAAFsXySwAAAACwMEIdAAAAAFgYoQ4AAAAALIxQBwAAAAAWRqgDAAAAAAsj1AEAAACAhRHqAAAAAMDCCHUAAAAAYGGEOgAAAACwMEIdAAAAAFgYoQ4AAAAALIxQBwAAAAAWRqgDAAAAAAsLdPcJAAAAwKIyMkRSEvM/xhYi4udXUmcE+CRCHQAAAIpGA93rVfM/ZtgBkaBQfsKAC7H8EgAAAAAsjEqdG508edJcVEpKigQEBLjzdAAAAABYEJU6N5owYYLUrl3bXGJjYyU+Pt6dpwMAAADAggh1bjR06FDZtWuXudSrV08iIyPdeToAAAAALIjll24UERFhLspms7nzVAAAAABYFJU6AAAAALAwQh0AAAAAWBjLLwEAAFA0Olhc59Bd6hgALkWoAwAAQNH4+TFYHPAALL8EAAAAAAsj1AEAAACAhRHqAAAAAMDCCHUAAAAAYGGEOgAAAACwMEIdAAAAAFgYoQ4AAAAALIxQBwAAAAAWRqgDAAAAAAsj1AEAAACAhRHqAAAAAMDCCHUAAAAAYGGEOgAAAACwMEIdAAAAAFiYX0ZGRoa7TwIiYWFhkpKSIjExMfw4AAAAAB+yY8cOsdlscubMmSJ9PpU6DxEaGmp+kZ4qPT1djh49at768nm54us5498szr9RlM8tzOcU9Fj9z0wvvs5T72slfW6u+lpWu7+56njub0X7+ZYUHtuK/3Pgsc2zeOp9zZMe2zQHaB4oMq3UAZeya9cureiat758Xq74es74N4vzbxTlcwvzOQU9tmHDhubi6zz1vlbS5+aqr2W1+5urjuf+VrSfb0nhsa34Pwce2zyLp97XvOWxTVGpAwAAAAALI9QBAAAAgIUR6lAgERER8tJLL5m3vnxervh6zvg3i/NvFOVzC/M5nvq346k8+edVkufmqq9ltfubq4/3dZ768+Kxrfg/Bx7bPIun3te85bFN0f0SgEdo1KiRebtlyxZ3nwrg9bi/AdzX4F2o1AEAAACAhVGpAwAAAAALo1IHAAAAABZGqAMAAAAACyPUAQAAAICFEeoAAAAAwMIIdQAAAABgYYQ6AAAAALAwQh0Ay/nss8+ke/fuUqVKFSlbtqx07NhRfv75Z3efFuB1ZsyYIa1atZKIiAgJDQ2VFi1ayOzZs919WoDX27RpkwQGBkr16tXdfSqwiEB3nwAAFNaECROkXr16MmnSJClTpox8+OGHct1118mqVaukadOm/EABJzlx4oTccsst0qxZMylVqpTMmzdP+vbta97X2wG4xtChQyUyMpIfLwqM4eMALCc+Pj7bg116ero0btxY2rdvL++//75bzw3wdldffbWpkn/++efuPhXAK+mLJ08++aTceeedMnPmTImLi3P3KcECWH4JwHJyvnrp7+8vV1xxhezatctt5wT40v0vJSXF3acBeKXz58/LM888I2PGjJHg4GB3nw4shFAHwKnWrl1rHox69+5t9gL4+fmZy6UkJSXJiBEjpH79+mZpV9WqVeWBBx6Q/fv3X/Jz09LSZPXq1VK3bl0nfReA5yvJ+1pqaqqcPn1a5syZI4sWLZJHHnnEyd8N4NlK6v6m2wsqVqwoffr0ccF3AW/G8ksATqX7bObPn3/R7RkZGXl+TnJyslx77bWycuVKs6yrQ4cOsnv3brNHTh/c9PY6derk+fkTJ06Up59+WtavX2+WYQK+oKTua4cOHTLHqoCAAJk8ebI8/PDDLviOAN++vx0+fNiEv4ULF0q7du3k5ZdflqlTp7L8EgVCoxQATqUPRE2aNJHWrVubS61ateTcuXP5fs6oUaPMg5t+7g8//GCan6i33nrLhDV9VXPZsmW5fu5vv/0mzz//vAwfPpxAB59SUve1ChUqmEr4mTNnzJPNwYMHmyWYt912m0u/P8DX7m/Dhg2Tbt26meOBwqJSB8CldLmJPvDl9Wqm7h+oVKmSnDp1StatWyfNmzfP9nHtZrlx40ZZs2aNtGzZMtvH9BXPtm3bmpEGuiysIEthAG/lyvtaVgMGDJDly5fLX3/95fTvAfDV+9vmzZvNWw2BtWvXNsfock9tlLJlyxYJCQmRoKCgEvneYE3sqQPgVitWrDAPejExMRc96Knbb7/dvP3qq6+y3X7y5Ekzq05fLdVZWgQ6wDX3tZx0vMHOnTv5cQNOvL9t377dBEGdBVmuXDlzGTt2rBw4cMC8P336dH7eyBfLLwG41e+//27e6gNZbuy36yuadvrAp5vVExMTZcmSJVK6dOkSOlvAt+5rufnll1/MiykAnHd/01EhS5cuzXbMRx99JN98840ZH6J77YD8EOoAuNXevXvNW+0mlhv77Xv27HHcNmjQILP864MPPjBjDOyjDLT9c26viAIo2n1Nmzzo3rkGDRqYpg/aKGLWrFnMgwSc/Nime1c7deqU7Rjdb6ePazlvB3JDqAPgVmfPnjVvdb9AbkJDQ81bbdJgt3jxYjNw/MEHH8x2bHR0tNlnB8A59zXd9/POO+/Ivn37zMcbNmxolovdfPPN/IgBJz+2AcVBqANgOQQ3oGTozCy9ACh5OtJAL0BB0CgFgFvZWzzr/rjcJCQkmLdhYWElel6At+G+BnB/g/ci1AFwq5o1a5q3cXFxuX7cfrsurQTAfQ2wAh7bUNIIdQDcSvfsKJ3jkxv77Tr0FQD3NcAKeGxDSSPUAXCr9u3bS9myZWXHjh2yYcOGiz7+xRdfmLc9evRww9kB3oP7GsD9Dd6LUAfArYKCgmTw4MHm/ccee8yxh0699dZbZobPNddcIy1btnTjWQLWx30N4P4G7+WXkZGR4e6TAOA9dFDqyJEjHddXrVol+t9MmzZtHLe9+OKL0r17d8d1nX+lc3h+++03qVKlinTo0MHM7tHrFStWlJUrV0qdOnVK/HsBPBn3NYD7G2DHSAMATnX06FETxnLKepsek1WpUqVk6dKlMnr0aDPYeN68eVK+fHnp37+/CYh5DW8FfBn3NYD7G2BHpQ4AAAAALIw9dQAAAABgYYQ6AAAAALAwQh0AAAAAWBihDgAAAAAsjFAHAAAAABZGqAMAAAAACyPUAQAAAICFEeoAAAAAwMIIdQAAAABgYYQ6AAAAALAwQh0AAAAAWBihDgAAAAAsjFAHAAAAABZGqAMAwA2GDh0qfn5+5vLmm2/mesyZM2fE399fQkNDJT09vcTPEQBgDYQ6AADc4Pfff3e8v2DBgjyPycjIkMaNG5twBwBAbniEAADAjaHusssukxUrVkh8fHyexzRt2rTEzw8AYB2EOgAAStjevXvlxIkTUqNGDbnnnnskLS1Nvvnmm4uO27Bhg3nbrFkzfkcAgDwR6gAAKGFZw1rPnj3zXIJJpQ4AUBCEOgAASljWsNakSROpVauWfP/993Lu3DnHMVq927x5s2mkoscAAJAXQh0AACUs57LKHj16yNmzZ2XJkiWOY/766y9JSkqSmJgYKVOmDL8jAECeCHUAALipUpc11Kn58+c7jmE/HQCgoAh1AACUIJ09t3PnTgkLC5M6deqY2zp16iTh4eHy9ddfmxEGiv10AICCItQBAFCCNm7caIKb7pPT/XLKZrNJt27dZP/+/bJ27VpzG5U6AEBBEeoAAChB9rCWc/ZcziWYVOoAAAVFqAMAwI376exuuukmCQwMNKMNjhw5IocOHZLy5cubWXYAAOSHUAcAgAdU6jTAtW/f3izPnDdvXq7HAACQG0IdAAAlxD57LiAgQBo3bnzRx+1LMEePHp1rNQ8AgNwQ6gAAKCH22XP16tWT0qVLX/Txnj17mre7d+82b6nUAQAKglAHAICb99PZadhr0KCB4zqVOgBAQfhl2AfiAAAAAAAsh0odAAAAAFgYoQ4AAAAALIxQBwAAAAAWRqgDAAAAAAsj1AEAAACAhRHqAAAAAMDCCHUAAAAAYGGEOgAAAACwMEIdAAAAAFgYoQ4AAAAALIxQBwAAAAAWRqgDAAAAAAsj1AEAAACAhRHqAAAAAMDCCHUAAAAAYGGEOgAAAACwMEIdAAAAAFgYoQ4AAAAAxLr+H6q6GdngNUMtAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots(dpi=150, figsize=(6, 4))\n", "\n", "for shape in shapes:\n", " c_comm, c_inc = colors[shape]\n", " m_comm, m_inc = markers[shape]\n", " d = results[shape]\n", " ax.plot(d['Ns_comm'], d['Fs_comm'],\n", " marker=m_comm, ms=3, ls='none', color=c_comm,\n", " label=r'%s comm. ($\\theta=0$)' % shape)\n", " ax.plot(d['Ns_inc'], d['Fs_inc'],\n", " marker=m_inc, ms=3, ls='none', color=c_inc,\n", " label=r'%s incomm. ($\\theta=%.1f°$)' % (shape, theta_inc))\n", "\n", "# Reference lines anchored to the data (not to F1s for the sqrt line,\n", "# to allow for the finite-size offset at small N).\n", "Nc = results['circle']['Ns_comm']\n", "ax.plot(Nc, F1s * Nc,\n", " '--k', lw=0.8, label=r'$F_s \\propto N$ (exact for comm.)')\n", "ax.plot(Nc, 2. * (results['circle']['Fs_inc'][0] / Nc[0]**0.5) * Nc**0.5,\n", " ':k', lw=0.8, label=r'$F_s \\propto N^{1/2}$')\n", "ax.plot(Nc, 6. * (results['circle']['Fs_inc'][0] / Nc[0]**0.25) * Nc**0.25,\n", " '-.k', lw=0.8, label=r'$F_s \\propto N^{1/4}$')\n", "\n", "ax.set_xscale('log')\n", "ax.set_yscale('log')\n", "ax.set_xlabel(r'$N$')\n", "ax.set_ylabel(r'$F_\\mathrm{s}$ [$\\epsilon/a$]')\n", "ax.legend(fontsize=7, ncol=2)\n", "ax.set_title('Structural superlubricity: scaling laws (theta_inc=%.1f deg)' % theta_inc)\n", "plt.tight_layout()\n", "plt.show()" ] }, { "cell_type": "markdown", "id": "9da57742", "metadata": {}, "source": [ "As expected:\n", " \n", " - **Commensurate**: $F_s$ grows **linearly** with $N$, matching $F_s = N \\cdot F_{1s}$ exactly.\n", " - **Incommensurate**: $F_s$ grows **sublinearly**. At large $N$ the data approach\n", " $\\sqrt{N}$; at small $N$ (cluster smaller than one Moiré tile) a steeper $N^{1/4}$\n", " behaviour is sometimes visible. Non-monotonic oscillations at intermediate $N$ are\n", " related to the number of complete Moiré tiles fitting in the cluster — each\n", " additional tile contributes a nearly-zero net force increment.\n", " \n", "### References\n", "[1] Dietzel et al., PRL 111, 235502 (2013) -- scaling laws of structural lubricity\n", "\n", "[2] Koren & Duerig, PRB 94, 045401 (2016) -- moiré scaling in twisted graphene\n", "\n", "[3] Panizon, Silva et al., Nanoscale 15, 1299 (2023) -- frictionless nanohighways\n", "\n", "[4] Dienwiebel et al., PRL 92, 126101 (2004) -- superlubricity of graphite\n", "\n", "[5] Cihan et al., Nat. Commun. 7, 12055 (2016) -- structural lubricity under ambient conditions\n", "\n", "[6] Vanossi et al., Nat. Commun. 11, 4657 (2020) -- structural lubricity review" ] } ], "metadata": { "kernelspec": { "display_name": "DRIFT2", "language": "python", "name": "drift" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.11.13" } }, "nbformat": 4, "nbformat_minor": 5 }