{ "cells": [ { "cell_type": "markdown", "id": "7b0fb08a-ad81-4af6-a4e3-def2d9b56dea", "metadata": {}, "source": [ "# Introduction to the EGM2008 model\n", "\n", "```{versionadded} 7.3.0\n", "\n", "```\n", "\n", "heyoka.py includes an implementation of the [EGM2008](https://en.wikipedia.org/wiki/Earth_Gravitational_Model#EGM2008) geopotential model, enabling precise calculations of Earth's gravitational potential and acceleration. This is essential for accurately modeling the trajectories of Earth-orbiting spacecraft, particularly those in [LEO](https://en.wikipedia.org/wiki/Low_Earth_orbit).\n", "\n", "Two functions are available: {py:func}`heyoka.model.egm2008_pot()`, for the computation of the potential, and {py:func}`heyoka.model.egm2008_acc()`, for the computation of the acceleration. Both functions require as input a Cartesian position in the Earth-centred Earth-fixed [WGS84](https://en.wikipedia.org/wiki/World_Geodetic_System#WGS84) frame, which, for astrodynamical applications, can be considered coincident with the [ITRS](https://en.wikipedia.org/wiki/International_Terrestrial_Reference_System_and_Frame).\n", "\n", "## Gravity on the reference ellipsoid\n", "\n", "For this tutorial, we will be computing the gravitational acceleration on Earth's [reference ellipsoid](https://en.wikipedia.org/wiki/Earth_ellipsoid). We begin with the introduction of symbolic variables representing the geodetic latitude and longitude:" ] }, { "cell_type": "code", "execution_count": 2, "id": "c93159ea-7a0b-4897-bd44-99a8e9b34530", "metadata": {}, "outputs": [], "source": [ "import heyoka as hy\n", "\n", "lat, lon = hy.make_vars(\"lat\", \"lon\")" ] }, { "cell_type": "markdown", "id": "2f2ba212-d3b5-4db1-87cc-c98e7f1f5445", "metadata": {}, "source": [ "Next, we formulate the analytical expression for the gravitational acceleration as a function of latitude and longitude via the {py:func}`~heyoka.model.egm2008_acc()` function:" ] }, { "cell_type": "code", "execution_count": 4, "id": "d51b4ac3-ae3c-48cf-a5be-4de83afe028c", "metadata": {}, "outputs": [], "source": [ "acc_model = hy.model.egm2008_acc(hy.model.geo2cart([0.0, lat, lon]), n=80, m=80)" ] }, { "cell_type": "markdown", "id": "9b504a2e-fb13-45b4-a146-4c88894c5a4b", "metadata": {}, "source": [ "There are several things to note here.\n", "\n", "To begin with, recall how {py:func}`~heyoka.model.egm2008_acc()` accepts in input a position in Cartesian coordinates. Here, however, we want to express the acceleration as a function of latitude and longitude, thus we use the {py:func}`~heyoka.model.geo2cart()` function to convert from geodetic coordinates to Cartesian coordinates. In the conversion, we set the geodetic height $h$ to zero, so that we will be computing the acceleration at \"sea level\".\n", "\n", "Additionally, in the invocation of {py:func}`~heyoka.model.egm2008_acc()`, we supply as ``n`` and ``m`` parameters the maximum degree and order of the expansion of the geopotential in [spherical harmonics](https://en.wikipedia.org/wiki/Spherical_harmonics). Higher values of ``n`` and ``m`` result in a more accurate calculation, at the price of increased computational complexity.\n", "\n", "We can now proceed to the construction of a [compiled function](<./compiled_functions.ipynb>) for the numerical evaluation of the acceleration:" ] }, { "cell_type": "code", "execution_count": null, "id": "9cbe5309-0cb7-4a0a-a54a-0805d476d32a", "metadata": {}, "outputs": [], "source": [ "cf = hy.cfunc(\n", " acc_model,\n", " [lat, lon],\n", " compact_mode=True,\n", ")" ] }, { "cell_type": "markdown", "id": "0eb6c848-6229-4626-abb0-674d06241d97", "metadata": {}, "source": [ "We are now ready to compute the acceleration on a grid of latitudes and longitudes. We will be using [sphere point picking](https://mathworld.wolfram.com/SpherePointPicking.html) to select uniformly-distributed locations on the reference ellipsoid:" ] }, { "cell_type": "code", "execution_count": 11, "id": "e4300e84-a11b-46bd-b4c7-6ad91e70b264", "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "\n", "# Generate a grid of latitudes and longitudes (500x500).\n", "nsamples = 500\n", "u = np.linspace(0, 1.0, nsamples)\n", "v = np.linspace(0, 1.0, nsamples)\n", "lons = 2 * np.pi * u\n", "lats = np.pi / 2 - np.arccos(2 * v - 1)\n", "lons, lats = np.meshgrid(lons, lats)\n", "cf_grid = np.ascontiguousarray([lats.flatten(), lons.flatten()])\n", "\n", "# Compute the acceleration on the grid.\n", "acc = np.linalg.norm(cf(cf_grid), axis=0)" ] }, { "cell_type": "markdown", "id": "71b4a91f-4733-4a30-a214-66592a06e994", "metadata": {}, "source": [ "We can now visualise the acceleration on the reference ellipsoid with the help of [cartopy](https://github.com/SciTools/cartopy):" ] }, { "cell_type": "code", "execution_count": 12, "id": "aa800276-f29e-48ee-9778-cccfab77a23f", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAw0AAAGiCAYAAAClNCryAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjkuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/TGe4hAAAACXBIWXMAAA9hAAAPYQGoP6dpAAEAAElEQVR4nOzddZhV1frA8e8+ndM9TBDD0N1Id4soYICgmBioP691veq10WtcA/WqmCiChIEiIWnQ3TXdfebM6f3748CBYYIBBmYG1ud5eHTO2bH2qb3eFe+SZFmWEQRBEARBEARBqIKirgsgCIIgCIIgCEL9JoIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAEQRCqJYIGQRAEQRAE4YpWUlLCrFmziIuLQ6/X06tXLzZv3lztPl9//TXt27fHYDAQGRnJ9OnTycvL8z2/aNEiunTpQkBAAEajkQ4dOvDll1+WO8bLL79M165dMZvNhIWFce2113Lw4MFy20ybNg1Jksr969GjR+1dfC0RQYMgCIIgCIJwRZsxYwYrVqzgyy+/ZPfu3QwdOpTBgweTlpZW6fYbNmxg6tSp3H777ezdu5cFCxawefNmZsyY4dsmKCiIp556ij///JNdu3Yxffp0pk+fzvLly33brF27lpkzZ/LXX3+xYsUKXC4XQ4cOpbS0tNz5hg8fTkZGhu/fsmXLLs0LcREkWZblui6EIAiCIAiCIFwKZWVlmM1mli5dyqhRo3yPd+jQgdGjR/PCCy9U2Of1119nzpw5HD161PfYO++8w+zZs0lJSanyXJ06dWLUqFE8//zzlT6fk5NDWFgYa9eupW/fvoC3p6GwsJAlS5Zc4BVeHqq6LoAgCIIgCIJw5bPZbDgcjlo5lizLSJJU7jGtVotWq62wrcvlwu12o9Ppyj2u1+vZsGFDpcfv1asXTz31FMuWLWPEiBFkZ2ezcOHCckHH2eVZvXo1Bw8e5NVXX62y3EVFRYC3l+JMa9asISwsjICAAPr168eLL75IWFhYlcepC6KnQRAEQRAEQbikbDYbsbFGcnI8tXI8k8mExWIp99gzzzzDs88+W+n2vXr1QqPRMG/ePMLDw/nmm2+YOnUqCQkJFeYYnLJw4UKmT5+OzWbD5XIxduxYFi5ciFqt9m1TVFREdHQ0drsdpVLJ+++/z2233Vbp8WRZZty4cRQUFLB+/Xrf4/Pnz8dkMhEXF8fx48d5+umncblcbN26tdIgqK6IoEEQBEEQBEG4pIqLi/H392fN32GYTNK5d6iGxSLTv3s2KSkp+Pn5+R6vqqcB4OjRo9x2222sW7cOpVJJp06daN68Odu2bWPfvn0Vtt+3bx+DBw/moYceYtiwYWRkZPDoo4/StWtXPvnkE992Ho+HY8eOYbFYWLVqFc8//zxLliyhf//+FY45c+ZMfv75ZzZs2ECjRo2qvL6MjAzi4uL49ttvue66687jlbm0RNAgCIIgCIIgXFKngoYte8MxmS8uD4+lxEOX1lkUFRWVCxpqorS0lOLiYiIjI5k0aRIWi4Wff/65wnZTpkzBZrOxYMEC32MbNmygT58+pKenExkZWenxZ8yYQUpKSrnJ0AD3338/S5YsYd26dTRu3Pic5UxISGDGjBk89thj53V9l5LIniQIgiAIgiBcFYxGI5GRkRQUFLB8+XLGjRtX6XZWqxWFonw1WalUAt5hRlWRZRm73V7u7/vuu49FixaxevXqGgUMeXl5pKSkVBmY1BUxEVoQBEEQBEG4oi1fvhxZlklMTOTIkSM8+uijJCYmMn36dACeeOIJ0tLS+OKLLwAYM2YMd9xxB3PmzPENT5o1axbdunUjKioK8K7B0KVLF5o2bYrD4WDZsmV88cUXzJkzx3femTNnMm/ePJYuXYrZbCYzMxMAf39/9Ho9FouFZ599lgkTJhAZGcmJEyd48sknCQkJYfz48Zf5VaqeCBoEQRAEQRCEK1pRURFPPPEEqampBAUFMWHCBF588UXfpOaMjAySk5N920+bNo2SkhLeffddHnnkEQICAhg4cGC5zEilpaXce++9pKamotfradGiBV999RWTJk3ybXMqgDh7jsPcuXOZNm0aSqWS3bt388UXX1BYWEhkZCQDBgxg/vz5mM3mS/iKnD8xp0EQBEEQBEG4pOrLnAbhwok5DYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEsEDYIgCIIgCIIgVEtV1wUQ6q8jR46wbNkyPB5PXRdFEARBOA8KhYKRI0fSrFmzui6KIAhXCBE0COXIsszKlSt5661XWblyDd17mTAYpLoullBDsgxSXb1dMiCdLofD7sHt9n6mZPmMp09uI0kSvqJK3ocVSlCpJJRK8ZkThIthtco8+ujDDB7cn1mzHmPw4MFIdfbjIAjClUAEDQIAVquVr776irfeeoWs7Awm36JmxbPBhEco67poV6yyMpn5X5cC4HSA3S5jK5MpKPCQneXB7ZYJCFRQVOgh6bgbq1XGbpdxOmVkD5zqANLqJJwOGafz9LEVClAqKffY2Vq2VnH9ZANOp4zDDk6njNt9+nmNBkLDlRTkecjJ9mCxePAPUBAUpGDbFgerV9jLne9SdkjFxSt5471AbGUy1jKZMqsHhwPUatDpJAxGifYdNWg0olIkCKdkZWqZ//UmJt84lvCwSGbNepxbbrkFg8FQ10UTBKEBkmRZluu6EELdyczM5N13/8ucOe8SFi4z5XYlo8fq0epE5as2uFwyixeUsX2rg0XflfkeV6shJFRBRnrNatp33GPEaFag1Z5qiQeF5G29t9lk1GoJtwsyM91s3eRg145qooUzKBSg00uo1aDRSKhU+FojrVYPhQUyZj+J0DAFJpM3gCnI91BcXP5n45kX/dBoJG+PwsnejvAIBSazAoWCk/+8x3W7ZSaOzatR+c5H/0FaPpgbVOvHPcXlkinI95Cf5yEkVEFQsIKD+10UFnhYuqiMUosHjVZCq5UIDVMiSZCd5cZolDCZvdt7PHDdRD0m08VPJ/u/Bwr4aYmt0udatVExb1EIunr+Pf7681Jmv1CM3e79PlhLZRRKiG6kpFGMksSWalq1UdO6rZrIKNGAcaHsNpmffijjy0/cZGdJ3Hvv/cyceT8RERF1XTThKlJcXIy/vz9b9oZjMl/cb6ClxEOX1lkUFRXh5+dXSyUUzkUEDVepffv28frrLzNv3rf07G3g1hkqevTWiO7rk/Jy3az93e5t8XfJTL7FQHDI6UrL+jU27phaQPMWKlq3VRMRqSQsXMHQETqcTlj1mw2VSsJm8/DycyXljv3QP8zodBIb19vYvsVJSUn1X8HgEAUbt4X7/r5rWj5rV3tb+RNbKsnJ8pCff/oYkgSBQRL5eeWPO2ioltZt1fgHKNAbJIaP0qHXS9W+5w6HXGXrvd0mc/SIi4hIBUHBtVOhczplLBYZlRK++qyUwwddDB+t44G7CsttZzBIJLRQERCgwM/fW1G/cYqR1m3VtVKOsy34xsrTjxXVyrE++yaIHr21F32cUQNzOHrEVeXzW/eHYzTW71wXreIzatxDdeMUAxnpbuw2mcAgBVHR3sCiYxcNiS0vzft+pZFlmT83OPjiExd//VHGTTdN5pFHHqdVq1Z1XTThKlDXQUNJSQlPP/00ixcvJjs7m44dO/L222/TtWvXKvf5+uuvmT17NocPH8bf35/hw4fz+uuvExwcDMCiRYt46aWXOHLkCE6nk4SEBB555BGmTJlS7jjvv/8+r732GhkZGbRu3Zq33nqLPn36+J6XZZnnnnuOjz76iIKCArp37857771H69atL+DVuXRE0HAVkWWZ9evX8+qr/2bVqjWMGW9k+h0amiZcPTdch13m2FEXhw+5sNu8w3F27XBw/WQDHTtr2PSXnakT8yvsN3iYlnf/d7oVu0VsRqXHb9VGRVi4kjWr7JU+f4peL2G3y5VWmJ5/1Z8RY3TYyrzDkw4dcHHimIt1a+wcPexi987TvQhxjZUkHXdXPMhJUdFKkODW24zcOsNYbZnqsx7tMykskHnoH2Ym3WIgIODyVobTUl0M6pVT4fF/v+JPl24aNqyz89KzxVXuL0lwzwMmwiOUjBqru+gb5pmsVg9ZmR7MZomQ0IbXGm+3yRw84GTFrzaOH3WRlur29WZZS899ezIaJcaM1zNhkp627TWXocRXhiOHnHz2Pwc/Lill0KD+PP74M1xzzTWi4Ui4ZOo6aJg0aRJ79uxhzpw5REVF8dVXX/Hmm2+yb98+oqOjK2y/YcMG+vXrx5tvvsmYMWNIS0vj7rvvJiEhgcWLFwOwZs0aCgoKaNGiBRqNhp9++olHHnmEn3/+mWHDhgEwf/58pkyZwvvvv0/v3r358MMP+fjjj9m3bx+xsbEAvPrqq7z44ot89tlnNG/enBdeeIF169Zx8OBBzGbzRb1WtUkEDVcBj8fDjz/+yMsvP8P+/fuZPEXLLdP0hIU3vArG+bJaPXz6YSnHj7lYu9qOpYpW/e49Ndw41cChA07ef7u03HMRkQo+/iqIiEilb1jJhnU2PvufldAwiawMD6mpbsqsMiaTRMs2apb9UHHYSGi4Aj+zRKs2an6sYlhJi1YqxlyrJzhEgdlPwavPF5OcdDoouGGynshoJQmJavr01/LzD2U89X9Vt4B/+FkgWq132JBaLaHVQeMmqlqttF4qySdcDO1bvqL+8n/8GX9D3Y3HTk1xsX6NneeeOh0gtGytYv9eF2HhCibfYmDQMB3jhuZW2PdAcuTlLGqNlRR7OHTARXq6G2upjNXqwVIic+SQi/17nb7P3+RbDDzzot9lqVRmZrgpLPBgtXqQJAgKVvKvx4rISHfjdMlkpHmjbY0GHI7T+xmMEj17ayi1eL/nkgRqjURUtJLYeCXxjVV07KwhMOjSf/7z89xkpHuw22VUKu+wwuAQBcEhClSq+lMxz85y89VnZXz7pZ1WrVrz+OPPMGbMGBSK+v8bITQsdRk0lJWVYTabWbp0KaNGjfI93qFDB0aPHs0LL7xQYZ/XX3+dOXPmcPToUd9j77zzDrNnzyYlJaXKc3Xq1IlRo0bx/PPPA9C9e3c6derEnDlzfNu0bNmSa6+9lpdffhlZlomKimLWrFk89thjANjtdsLDw3n11Ve56667zv2CXCYiaLiCuVwuvv32W1566Rny8zOYOkPFpJsMDaLCWFvS09wM7Jl9Xvs0b6Hi0IHKh33sOhyBRisxYVQue3efbvE3GCVfq2iLlioO7C+/f2JLFQdPPqbVgr2Sjogz62KVfSvnzgui5zUVh7VkZbq5eUIeqSlV9zicfR4/f4mRY/Q886J/jfapC0cOORk9+HTl+6apBu6caSIisu6DXYddZtcOJ3t3Ozl21IVWC7fdVb5ssuydYF7f5gcdPuhk1Qo7u3c42L/XRXpazT43YeEKVv8ZVqsV3r27naQkuzh80MWmvxwU5nsoLPCQk1P5mKXAIAWFBR6CQxR8uziY6Bgl+Xke5vzXQmSUkrIymbWr7UQ3UqJQArI3wUBqipuUJG8yAYDHnjYz/Q5TrV3H2Y4edjJqUMXA8RSDQSIsQkHnrhoSElV88amV2FglbdqpiYhSojdItGqtpmXry9cLbCnxMP9rK1984iIoKJInn3yOyZMno1KJfClC7bgUQUNKSkq5oEGr1aLVVrxPlpSU4Ofnx8qVKxk0aJDv8Z49e6LValmzZk2Fff744w8GDBjA4sWLGTFiBNnZ2UycOJGWLVvywQcfVNhelmVWr17N2LFjWbJkCUOGDMHhcGAwGFiwYAHjx4/3bfvggw+yY8cO1q5dy7Fjx2jatCnbtm2jY8eOvm3GjRtHQEAAn3/++YW+TLVOBA1XILvdzueff84rrzyH21PE7XerGH+9Ho22flVealtGupt5n5diLfNmISorkzEaJUqKPajUEkaTAn9/+Og96wUdf9x1el59KwCA//6nuEKPhMksVejJ6NZTw6Y/Hej1EB6p5MSxmlXQAMxm+HNnOCpVzX5c3W5v67DTCWoNqFUSCiUoFd5AweEEa6nM4YNOnnjE2zvRp7+WD+YG1ssUp3+st3PbzaeHitXXlvqG5PBBJ2OGeCuzPXppiIhS8tsym68y7e8v0biZirh4FTGx3jkD0TFKoqK9/y6ml8HjkSkqlJEUoFLC8l9svl4yPz+J7r21REZ659vodBI9r9Gi10s4HDJZmR5atlYR3ej8K7CyLPPUo0XlEhH06a/lrTkBl2zOh8MuM6RPNlmZ556w0ae/hvVrHJU+tz8p4rIPF3LYZRYvLOPjOS5USn8ef/wZbr311korYoJwPk4FDU/8ORyd6eICYpvFycs9f63w+DPPPMOzzz5b6T69evVCo9Ewb948wsPD+eabb5g6dSoJCQkcPHiw0n0WLlzI9OnTsdlsuFwuxo4dy8KFC1GrT5e/qKiI6Oho7HY7SqWS999/n9tuuw2A9PR0oqOj2bhxI7169fLt89JLL/H5559z8OBB/vjjD3r37k1aWhpRUVG+be68806SkpJYvnz5hbxEl4RoQriC2Gw2PvnkE15++Tn0hjLueVDJyLHmetUVfil4PN7JfQ/dW0BxsUxEpILwSCU6ncSxI97x3gX55W/eY8brmHSTgXlflLLsR2+zf/MWKvoO0FJY4KHMKrPpLwc52d79nn/Vj/E3GMjKdPPME0Vs31r+Jq/RQJduatb97ig3T2HTn97t7HZvBiGzn0TJGZmHAoMUyLKM3eat3BiMpycwl5TAuKF5TJluYNIthnNWHpRKqUYTQtu0UzNyjJ5ffy7jsYeK+PozK1Nvrx/zHYoKPTz/ryI0GqlcJa9DZzWyLIvx1hchL9ftCxge+oeZUeN03H5LPjq9xE1TDQwerqNdB7Uvy1V18vPcHDviJjBIIihEib+/hEIh4XLJuF3le1dkWWbSuLxqM3ot+z202vkYbdufx4WexVIil/ssAaxfY6dzyyz8/SVcbu+Ef4NeQm/wpu/t3UdLl+4atFoJjUbyBa/RMd5rBfD3V9CmnZq2HdT07qstF4BotBJrN4WTk+1m/14n+/a68LhlYmJVHD3i5IN3Ss8oi4Pho3VMnW4gJdlNUZFMYYGHYaN0dfJ512glJt1sYMIkmWU/lPLa64/w73//kyeffJbbb79dBA9CvVJZT0NVvvzyS2677Taio6NRKpV06tSJm266iW3btlW6/b59+3jggQf417/+xbBhw8jIyODRRx/l7rvv5pNPPvFtZzab2bFjBxaLhVWrVvHwww/TpEkT+vfv79vm7O9yZfezmmxT10RPwxXAZrPx8ccf8/LLz2H2s3HPAyqGjtTV6ObfkBXke1i0wMq3X1pJSXYT11jJ+Ov13H63id07nfzwvZXVK+3kZHsqHe5zSqs2Kk4cOz10Abwt9c4z4oLb7zbyyQenb/Sjxuno0UtbIaNOk2ZKsjI9vjHVADod2GwQG6ekU1cNSxaersC0bqsmvom3RffDd8v3XJzps2+D6NGr9m/W44blcHC/i8W/hFzWoRCV+e9/Snj/bYvvb5NZIjZOyctvBJDY4uqZrH8p5OW66d3JO0zvlmkGbp5mZOoNeej0Ep9+HURMXM3bjwoKPPTvllVuiJ3BKPHibH8+/aiU3Tud3HaXkXHX6dm62UFJscybs70ZxCTJO3ncZJKQJDCaJIKClZcs69UpsixTUiyTlen2DePb9JeDg/uc9O6nRaeTKLPKWK3eCvtPS8oqpBUGb8NCpy7eyda5OW5273SSlenB31/izvtM3H7XuYc8ZWe56du14pDJpglK/P29c5nMfhJmc/n/tmytpm37mgV1tcnjkfltmY3333ZhKdH5ggedTndZyyE0fJeip+FCUq6WlpZSXFxMZGQkkyZNwmKx8PPPP1fYbsqUKdhsNhYsWOB7bMOGDfTp04f09HQiIyvv/Z4xYwYpKSksX778ihueJHoaGjCHw8Enn3zCv555Ere7lHse0HPLNFOt31RKij3s2uHkxHEXsXFKuvfS1vkiWlarh54dsnx/v/leAIFBCk4cczFyQDYpyeV7FkaM1jFhsp7WbTR88K6FwgIPA4dosZbKrF5pZ9hIva9iM+sfZibdZKBnhywGDNbyn3cDOHrYVS5o6Npdw5FDLgKDJAryZSTJOw/h2BE3Pa/RMHa8Hj9/BTNnFGCzwZxPAxkw2HuTves+E99+aeW6iXpfz4DNJvuCBoMBGsV6J22mpbo5dMDFtMn5JLZU0aKVmhGjdfQfdPE37GefLPLNs/hpSRlxjZUYDJd/vovbLfP151Y+/6R80DRvUTDNE6+sYKG4yPu5VCggJ8dDQZ6Hxs1UBAZe2tf91edPp/3teY2WY0dc5OR4GH2tjqIimUY1bNFyu2U+/7gUux1efM2fxk1U5OV5+M8rxTw0s9C33acflvLph973MzBIQUiogtwcb/D+9GNFrPoj9IKGGl0oSZLw85fw81eQcPIzder7WJknnvHDWupdTNHh8M5NcbllmiWoKvy+piS5uP+uAl57sYR1v9sZNtLboPDcP4v4+w8HI8fqCA1TcvigE7XaOyn7TK3aqOg/SEdxkYeSEpmSYm/v6NHDLoqLZSzF3sc9Hu/aJz16a2nVWk2rtipatlL7xoafmkNjLfP2mtTWXBqFQmL4aD1DR8r89ouNd//7OC+88C+eeeZFbrvtNjQakbFKaFiMRiNGo5GCggKWL1/O7NmzK93OarVWmNOjVHq/v9W1t8uyjP1kq4pGo6Fz586sWLGiXNCwYsUKxo0bB0Djxo2JiIhgxYoVvqDB4XCwdu1aXn311Qu/0EtABA0NkMvl4quvvuKppx4lLz8fu81bETl8oOINraYyM9zs3e0kK9PNsSMujhxycfyYt0J5alxuZav+3vugCbXaWzFo3VZNaJgSk0nCaJLweLwrEl/oAlNOp0xpqYyfn8TqFXY8Hpkhw71d9mo1PPuSH6+9VEKpRfZVWBQK6NRVXSFo+OUn77hto6mMx/7pV26l6/E3GJBl70rLOTGJ+A+K5PE5h4AsSgIj+OJEM+5uv4+9xyP4+YcyXn2+mGefLCYgUKKwQCYsXMG4CXqUSlj4bRkRkUrG32DgHw8W+s6xdrWdbj01GI0KGjdR8cQz5VtGdDqJtZvCfAuBnelUK/HB/S4cDpmflpTRopWaFq1UDB+lo0//CwsgIs5YLOuTD0v55MNSbpxiuKyTo91umQ/ftfDf/3h7GDp0VrNjq3cYS/IJd4MOGmRZJiPdQ0mx97P40fsWfl5aedasqGildxjPzQbuus9U613St84wsmWTg8ICDzNnFBAT633vf1pi8y0QN39pMO07Vl8BXDDPygfvWBgyQsd1E/W+ci77oazS+Tr+/hJ/7vCuMWKxePh9pZ09u5y1tq7HpaJWS/gH1Ow9iIlT8fX3wfy0xMZvv5Tx0rPFuM7Ig5CV4WbPTicJiSpkGbZtcRAeocDthtwcD9dOMJxzeKDLJbNti4NVy+1s3+rg15/KfD09oaEK7Hbvb+WZK7rrdN55VN16aBg1Tn/ePZVnD41QKCSGj9IzZLiHrq2zmPXQTF544V+8+OJr3HzzzWLCtFDvLV++HFmWSUxM5MiRIzz66KMkJiYyffp0AJ544gnS0tL44osvABgzZgx33HEHc+bM8Q1PmjVrFt26dfPNPXj55Zfp0qULTZs2xeFwsGzZMr744otymZIefvhhpkyZQpcuXejZsycfffQRycnJ3H333YC3UWPWrFm89NJLJCQkkJCQwEsvvYTBYOCmm266zK9S9cTwpAZElmUWLVrEk08+gt2RS+++MiGhEuERKpomqGjTTn1e8xeyMt18+5WVZT+W+XL9q1QQF+89XpNm3ptAXLyS9h01xDVW8s4bFpwOmZ+WlpGV6UGhgIBA70rBZ96wzH4SpRZvC3y7DmpatlHTtJm3jKcqJh6PzInjbvbtcVJY4CEg0NtTsHunk8MHT2d1OTW8B+Da6/W+4T3demrQaKCoUKZ7Lw0jx+pplqBCo5GqXEcB4NbbjTzxjB/zi7pUuU1+SimfTN2IJdd7Z+53dwJqrZLdy9LISyqlcfcQDq/PJqqRklUbQ/F44L47Cvh9pXf753aPIXl7Pp9M3VjuuEqNArVOidagxC9CT3C8EZfNg7XQQfL2fBKaKeh5jZYHHzFX2lLodMos/NbKrh1O9uzyvk4fzA28oJ6Hx2YVsnRRGUNG6Fjxy+nK7PaDEej1l6cnadkPZTx8X2Glz/2wIqTBBQ2yLPPXRgdffVbKjm1O8nJPB6+BQQpmzjIREuqtMIaEKggIUHDooItHHyj0bdetp4Z27dUMG6WjZevz+07XpHy7djhZ8I0Va6l3Ib11v3s/s0uXh/h6vub8t4Qtmxw4o6Mwh+owh+oISzBj8Nfw1T1/kZdUysgn29JuVDS3NtoJwJa/Hdx2Sx6Ok5XZBT8GExevws//6snWBt7epH17nOze6cTsJzH5ltqfL+Ryedeb2b/HmxJXr/c21BiNEjq9d6hVUaFMUpKLjevsnDjmZuRYHT17awkNU3D0iHedGo8H/AMUhIR4J6CXlnrTw6anuX29U42bKgkKUuDn7/3c2u2n5ly50agVZKQrMZvDeOXlN7juuuvq3Rhsof6o6+FJ3333HU888QSpqakEBQUxYcIEXnzxRfz9vQ1l06ZN48SJE+UyKb3zzjt88MEHHD9+nICAAAYOHMirr77qW9fhn//8J/Pnzyc1NRW9Xk+LFi148MEHmTRpUrlzv//++8yePZuMjAzatGnDm2++Sd++fX3Pn1rc7cMPPyy3uFubNm0u6nWqbSJoaCB+//13Hnz4ftLTD/HAQwaum6hHrZaqrfie4nK4ObQ2i6Rt+XhcMrYSJ6m7CshPsaLRK2k7MpqmvUKJ7RiEMVhbobdikv+Wc57DZvNm5Tm1KFNaipuAQAUej8ymPx0cPugi6YQ3s0+X7hrCIxRs+sPhS62oVEm4XTJ6fzXRbQKIaOGPX7iOrBUHSWylZvRYPV99XsrK5TZfpQQg4ZowrIUO0vYUolBKhDYxEd7cj4NrsrCXnm7uG/Z/rXDa3ah1SgwBmtOvQbKV4iwbTrsbc6iWqFYBBETp2fljKoXpZWdfZjl6PzVlxU7uXtCXkMYmXur+Cx736fzw3W9uzNBHWlGUWcZPz+/m2N+5BMUYcJS6sFlcOMvcGIM1hDU1o/NT4x+hZ/viFOylLsKb+zF9bi/0fhV/WN0uD398foyM/UXsX5nBxBv1PPvS+fcOTJ2U55uofaZndoxGoZRq9L5fKJdLZt8eJ2tX23nvLUu55zp1UfPu/wLrfWv02bb87eDt/5Sw+S+Hb8hJm3ZqQk62BDdP9K7GXZmcbDdP/l8RWzc5aN1OzeGDTgoLvEF3QKCC4GAFMXFKevXRck0/LfGNvZmMzvz+W3LtSEowBFS/srvb5SF1VwFH/8xFa1TR/ebGSBJ43DLZh0s4uDaLtR8c8m1vCFBjLfT2/ijVCiJb+pG6qxAAtV5Jj1sa03ViPP4RevatSGf+w1sB0BpVNO0VypwX3eV69mrymwU1+92pSzW9jrPV5LqqOvb5vibzi7ogyzLbFqewce4R8pNKkWVvylejUQIJCgs95eZvVUZrUhHXORilSkKpVmAtdFCYZqU4w+rrVdFoJRISW/LOW+8yYMCA8yqncHWo66BBuHgiaKjndu7cySOPPsQff/5BmcVOUKyB+5YOQFnDNJyWXDuvDfgNgIBoPTqTGrVeSVTrAGI7BNK0Zyh6/5qNSb3Ym/i83E4cXJvF1oVJ2EpcxHUOokmPEKJaBaD3V2MvdaE1qqqs8LhdHjwuGWuhA1mWMYfqfK9DzrESkrblk763kOOb8shPrnpScVVMId7ue5fDg7PMRcI1YfSa3hRJkvju4S2U5JyOVnR+ahxWF+YQLbd/eQ3+EfpKexZGPdWWPz47SnFWGRNe6UTrYVFk7C9iwT+2kneilNFPt6XrxHgAUnbk8/NLe8jY751crdYp6HRdHMP+rxVK9en3O+dYCe+OW+P7u9OEWLreEMdDvY6c1/X+71Br3hu/BluJ964f3TaAif/pTEDk6cXTLvY9z8p08/MPZRw97B1alZHuIT3NTVaGG7fbO9m5Z28Nvftq6T9IVy/WYDhfu3Y4eOu1Ev5Y76BlaxX3P2xmwGDtBbe4zi/qgtvlIXl7PvkpVkrz7ZTm2ck6VEzytnzcLpnAaANNe4fSclAETXuG4nHL/Lvj6Yl8wfFGYtoFMuzR1hgCNFjy7BzZmM3h9dkc2ZiDreR0JiOdWeX7DACoNApcDm8wn9g/nJve6YbT5ibrUDEpuwo4/ncuB9ecnk8E0LRnCFM/6oksy/z2xn7++Oz0Ykidr49l+D9ao9GL4St1zWlzU5pvxy9C72sckmUZu8WFw+pCa/L+Dh/7K4fDG7LZ80s6ADe/143mfcMrHM/jlinOLiP3uIWUnQUkbc4n64CFXj178frsN2jf/iLSXglXHBE0NHwiaKin0tPTeeKpx/n222/pOimOa2Y0ZfviZH57Yz/Neocy5KFW6P3VGIM0qDRVV7RcTg8LH93K/lWZTH6rCy0H1U6u+5pWJi+0Ne5irP3wEKvfrTzn8pkatQvgpne6oTGqOPZXLr+8sgdJIVUacJhDtfS8tSm/vb4P8FbKmnQPZejDLdEYvJUhl9ND8tY87KUuvp11+vVp1juUEY+3ISTehMvpYclTO9j9SxqT3uxCq8He9+Pwhmy+uudv3z5hzcwExxvZvzKT6Z/2JL5rSLnyFKZb2fVzGnuXp5N9tASPS+b62Z1oMTCCm0O3VVthnV/UhaSteXx5z9/oTCoG3NeC+M5BBMUay+13MQGDwy7z4XsWPnrPgqSAxBZqtFqJ8EiFL+d/80Q17TqqUasvzXCGy/HZ2/VTKoue2k5oEzP972lOy8GRlzS7jd3q4sTmPI5s9Fb+85NLie8SzJQPu7N9SQo/Pb+73PYJ14RhDNKw44dUAKLbBNDsmjCa9wnj8IZstixIIrSpmdZDI9EaVJjDdMR2DCJpax5J2/PpfWtT3+f7TG6nh7S9hRz9M4fjf+fSa2pTWgyM8D1fkGrlrRGrfAHIgJmJ9L+7+SV7XYRLw+PxVg/O5zNtLXSw4eOjbJ6fxI033shLL7xcLve8cPUSQUPDJ4KGeqa0tJTXXn+N2a+9SsI1YQx4sBlBMUZkWWb+w1vYvzKz3PaB0QYe/GWgr7InyzKF6WVkHykh+3Axu39JR++n5sSWPGI6BDLjy2vq4rIuK1mWKUi1kp9SisvuQaNXojGo0Pmpvb0RewppMTCCxt1DUKkr9tjkp5SSuqsAY6AWQ6AG2SMT2MhA7olSPp22kclvdSWxX8VWt1OcdjefTNlAxv5iOoyLoe3IKApSrRz7K5ejf+RgL3XRakgkN7zWGcXJRdU2zD3Cijf2A97WXaVGwYHVmRiDtNz9XV9Mwd5eELfTQ+ruAiw5dmI7B2EO0eF2elj05Hb2/OptFWzWO5QpH/QAylf8z6xEvz1yFfkpVu5d1I/whKp/cM8VOJxdMc9LsqDZuIXvv7OSkuRmxj0mbr/LiNmv4ut85r4up4e03QW0LtlHSKiCbj28Q2xSk10YzQqMBomN6+1s2+xApZYYOkLnSxFbF4EpwI6lKSx5egcdxsUw9tn2vveytmQeLGLO9esA6DC2EWHN/QhPMBPWzA9zqBaPS2b1ewfZ8MkRukyMwxyq4/f3ygfLCqWExy0TEm9k+tzevt60S23ZK3v4++vjvqDhzm/7EN064LKcW6gf8lNK+f3tIxzekM1jjz7G//3foxiN9WM9GKFuiKCh4RP9xfWELMt8/fXX3P/gfZSWWjAEaRj3UhtfL4IkSRRleMfYD3u0FWWFTvb8mkZ+ipWjf+ZwYnMeabsLSdtbiN3iqvQcxsCrIzWeJEkExRgJiql4gwqJN9H5uthq969q35j2Gp7ZPvqc58/YX4TzZEarnT+msGNpCpLkHf7Ta1pTEvuFE9HCr1yrfu9pTel4bQypO/P5/okdKNUKBj/Ykq6T4nzDOjwemW8e3Mzh9d4c74HRBu6Ydw2GQA3txzTyBQ2thpzuTaqsQm23uihIsxLdNoDQJuZzXk9lx3I7PWQdKib7aApZh4rJPFhMYXoZ+cmlqHUKmvUO487ZiYQ392OZDBRVfVxLnp03Bq/A7TrdfhHXKQizq4Q9u7zDaPR6ibIy78J9mRke5vzXwrO7RtfZpMstC5P46d+76HRdLKP/1e6cLbHZR0vYujCJzAPFuF0eWg6OpMOYRhiDqq7En9mcc6qn4BS9nxpHmRu30/s52/JdEhrD6R7HAfcmojWpyDlWQkhjE62HRl22gAHg6EbvZ9Tl8ND/3uZENBc39atNUIyRCa+3J3l7PnNfn8Ocj+bwn9lvcNNNN4nJ0oLQQImehnpg27Zt3HbHdPbt3YfTXr7C//S2UajUCpK25fHjc7vIOXZ60mhAlJ6EPmHIsrfS0GJAONFtA4ls4UdYgh9up4cl/9xB4x4h+IXradYrFP8I/eW+vKuGJc/Ozy/tZt9v5TM3jX2mHe3GNEKtPV2ps5U4WfHWfvatyMDjkvG4PTis3mxRap2CR1YOqTDXxFbi5PVBK/CP0DP4wRb8+O9dqLRK/MJ1pOwooFG7AAY90JIm3csPZTqbLMt898hWDm/I5v4fBpz3Z+LEljx+fG4nuSe8w7gCovVEJPrjH6mncddgmvYKPa/x61sWJPHjv3f5/u45pQlpewrQ+2toP7oRTrubkmwbif0jCGtm9g0/63FLY0Y8dvkzS5zqFeo2OZ4RT7Q5Z8CQcaCIz277E7VeSUz7QFx2N4fWZRPdJoA7v+lT7b4ej8yRDdks+L+tdJ0cT4exMXz3yBZyjllI7B/um1sw9P9a0XNKE0qybRzflEvroVGodXU3P8RW4iT3uAW1Tkm4CBiueh6PzJ5f0lj95hESmiTywXsflVvESrg6iJ6Ghk/0NNSh3Nxcnnjqcb784ktM4WoGP5RIfJdgijJtzLtvEwCLn9qOy+7hwOpMotsGcOPbXTGFajEFa/GP9OZJX/L0DgBKCxzsXpbGrp9TsRU7fRN3i7NszPp1UF1d5hXLZnFSlF6GtcjB7p/T2PlzKpzMsNn/3uY4y9zsXpbGD8/t4ofndjF9bi/iuwTj8cj876b1vkr3mQKi9Vz/SqdKJ6frzGpufq8bix7fzvyHtzD0kVak7SnEkmvn5ve6kdAnrEYteJIk0W5UNPtWZOCwVt4rVRVroYMv7/6LyBb+TJ/bnsiW/miNF/czEt8luNzfxmAtt39RcRhdYbqVbYuTOfZ3LoCvlf1ysZe6WPHmfjbPP0HfOxIYeH9ijV7vQ+uysJU4uWt+H4ozy/j5pT2Ad17MuSgUEs37htPxuli2L07m6J85voaDMycjF2WUoVBI+Efo6TA25gKvsPbozGoatQus62II9YRCIbF5/gkKskrYnL2F7j26ceu0W3n5xVcICam+kUMQhPpD9DTUAY/Hw8cff8yjj/0fjTr4M/jR5gTHlq9AOKwu5j+8heRt+Sg1CoJijfSe1pTmfcMrtCAWZljZ+UMqSdvySd1ZUC7VqKSA+38cWOH4woWzFjn47fV97P4lDZfdW3HVGlW+1737TY0xBmn4e95xSvPL5zKc9mlPGncNYfV7B1n7wSHMoVq0JjU6swq7xUVxtg27xcV1L3Wg/Rhv5W/Xz6ls/T6Za//dgcBGBlwON2+PXE3LQZGMfOLCWto3fnaU1e8e4Mm/RtQ4E1dpgZ3vH9tG8vZ8Zv062DfPojYUZlh5c+gq39+thkQS2MiAOVRH9pESjm/KpSDViiRBVOsABsxMJOGasFo7f3XcLm/Q/uvsvZQVORjycCu639j4nPvZrS7Wf3yYjXOPEhRjJKq1P7t+SiOmfSCj/tmWyBY1T5NbmGFl8VM7OLE5DwCFSuKeBf049ncOHrdMh3ExGGqYBU0Q6sL2xcks+ddO398BkUbcVonXZ/+HGTNmoFBcXet5XI1ET0PDJ4KGy2znzp3ccfftHE85wtDHm5PYP6LKbYsyy3hjyMpyjzXvF87N73ZDlmWStuaz4ZMjpO4qALxfIvlk42tkSz86TYij68Q4MX60FsiyzLG/cknZUcCB3zMoTC+j17SmNO4agtaoYv4jm8k9drrnQKmWCIo1Ysm1Yy914Tk5Xv+Or6+ptgX22N+5fD7jTwDu/b4f4c39WP/xYVa+fQCAtiOiie0cxG+v76X3bc0YcE/iBV1P5sEiPpy0nt7TmzL4wZbVbluUWcbGuUfZtjgZhULixne60rhr7bcOHvg9k5QdBZQVOchPtVKYaqU4q4yQxibiu4bQuFsw8V2Ca5wi+GJlHy1hx9IUdv2USkmOnYRrwhj1z7YERhvOue+OpSmsfHs/ZUVOrrm9GWXFTjZ9c5wRj7eh66T4C86w5LS5OfZXDkhStZPxBeFyOHvV6Jqwl7rY9XMqf315jNwTpWgNalq1bMVnn35Bu3btLlFJhfpABA0NnwgaLhOLxcK/nnma995/j+43N6bvnZWnMjxb1qFijv2dy/FNucgyXHNbU5K25bPrx1RyjlkIb+5Hm+FRlObbCW1ipuXgCLRGVbVpWIWa83hk9v2WwfpPDpN5oBhjkIbARkZGPtGG6DYB5bYtyixjw6dHOLIxm/xkq+/xkU+0Yc+vacR3CcEUomXj3CM4ytxEtwkgtKkZW4mTnKMWek5pQuthkexYmsrKt/ZjybMT1dqfxt1C2Dj3aLlzxXUO4ub3u6OtwWeoKqeCkVPDpqry8ZQN5By10OOWxnSbHF/t5N3adiGVkto455oPDrHm/UPo/dW0G9WIDmMbEdnKv0ZlSdtbyEeT19N6aCRDHm5FYLSBl3r8QpsRUYx9RuStFxq+H/+9iy0LkghrZmbKhz3wCzu/FendTg9/zTvuS2EtKUChUHL//ffx/L9fwGQyXYpiC3VMBA0Nn5jTcBn88ssvzLjrNnRhMnd805uwZjXPWBPe3I/w5n50nRTHyrcP8OVdfwHQYmAkwx9rTdOeoRdUqXKUuVDrlKIXohqZB4v46fndpOwsoEmPEG79Xw8adw+p8jXzTlBuyfYlKQB0mRhH5wmxfHXP35TmO7AWOMlPKaXl4EhCm5pI31PIoXVZKFUKjMFaFjy6lYwDzRgyqyXtRkVzaF0Wu35O82VFAu9icS0HRWAM0l50is/etzXj729OcOSP7CqDBkeZi5JsG62GRDLg3gvr1bgYdfH5TNqaz5r3D9HvrgT63plw3gH4obVZ6P3UTHi1k2/oV+vhUez5JZ1BD7TAGHj5gi5BuBiH12dxfHMeGfuKcDk99J7WlMT+4RhOZuLLPlLCF3f+yfS5vWr8uf7r62P89eVxCjOsdJoQiyFAw/G/c9GaVSxd+y3ftZrPJx/NZfjw4Zfy0gRBuAAiaLiE8vLyeGDW/SxZuphBDzWn04TYCx6WsG1xCn9+cYz2Yxox4rHWFzREoyTXxoo39rPzR2/6xuH/aE3LwRGkbC+gKLOM9mMbYQ45vxajK9HK1KYc/3Q9qd9vw9AokPavTySgQwwngCbSoWr31RpV/HPTyHKP9ZjShKKMMhxWF7YSJ9c+377S7EK/vLqHDZ8cYe+v6fhF6Jj4n9OL8VkLHeSesBDZ0r9cFqaLoVBINGobwLbvk1GplUgK0Ptr6DQhFpVagSzLLPnnDkoLHHS/Kb5WztkQnEpd2rxv+AX12OWesKDSKSnNs+MX7s1M1fWGOLZ9n0xhepkIGoQGweX08PXMTRiDtTRqF4in2Mk3D2wmuk0AiQO8Q+OGPNSSjZ8dZdET231rw1TmVI9haYGd5a/vo1nvMCa+0ZmoVgHltvN4ZLYuTGLCxPGMv3Y8b7/5DsHBVfeCCoJweYmg4RKQZZmFCxdyz8y7CG9t5O5F15x3WsukbXns/CGVzEPFmIK1HFyTRUC0noH3JdZoWNOZso+WsPKt/eWyrQAgwTtjfsdl96BUSfz99XF6Tm1CbKcgbMVO/CP155XHvyFZnVX56rTOIit7n11A8f4MGk/vTaPru6BQKyvdb2B41QHEqe0Ghh+i74wEnDY3L/f6lf73Nq8yHWn/exIxBmvZNO842YdL+KMwgaHByQAYAjTYS5ys+u8BhjzUssaTl89l7LPtWfiPrWz+7gTI3uBk38oM+tzejCMbs9n7WwaT3uxCRGLNJ+02dCHx3qERucctF5QBaND9LZg7/Q8+mLiOkU+0pc3wKLRmNZIEX971F4MeaEGbYVHkHLOQtDXPO8k7rQxHqQu1Xkn3mxrTdVKcGGIo1CmlSkJSSPS/pzldJ8YD3jlXq985wOp3DhLXKYiuk+JZ//GRCuva5CVZKEizEhJvYsnTO0jfV0TCNWGk7SnE45K59vn2lQbPCoVE14nxNO8bzo/P/05iywTmvPch119/vegVF4R6QAQNtSwnJ4fb77qDVWuWM/KJVrQZHnVBP3Y7lqaybVFyuccK08p4c9gqek1ryrBHWlW7f35KKbuXpbF3eTpZh0vKPRfa1ER+spVfX91LwjVhXPdyRxxlLla+dYCVb+0vt8hWeIKZtiOj6XVrU5SVrJ58uVVV2a8NHoeLPU8voSy9kPavT8S/TfRFl2V1VnM8DhfHP9uI2+mheZ+qJ6/+WdYaR+94St45SLN7+qMy6Vid1Rx3mZOi979n9y9pALQbGU1ULa2uawjQMPWjnr6/T2zJ49sHN/PFnd5hcAPvS6TV4Miqdr8iaQwq/CP1ZB8tqfCcxy2TtDWP2E5BVQZuQTFG7vymD1/P3MSiJ7fTcnAEwbFGZi7pz8a5R/np+d389PzucvsolBJhzcwUpFn5dfZeDAFqX/asK9WZgXVDURu/P3V9vbYSJ7YSJ1qTGq1JhUIhkbq7gO2LUzAEaghsZKDtyGjUWiVak4qDa7Jo3C2EkHgTTbqH0Lhbb5K25hMcZyTvhAVbiZMO47yf1ffGryH7SPnvjV+4jm6T4zmyMZuCVO9cL/msbMkVXlcJbn5Xx55f05l211S+/OZrPvnwf4SGhl6y10UQhHMTQUMtWrp0KdNn3IauVRj3Lel7URNGB89qgSlYizFYS7vR0Wz7PpkVb+4H4I/PjmIrcdLlhjjcTg97fkmn161NCIjyZnWxlTh5e+RqwJv3/0yhTUzEdQqm68R4QuJNxHcNRqlWYAjQcP2rnSh9rDXF2TZ0ZjVZB4vZ82s6q/57AK1JTbfJ8Rd8PZW5lAHAhTjy3u8U78ug6cwB5wwYaqpoXzoHX/sVW0YRcbf2Yq9/Z/ZlSb6Kw9mvQfaag+CR8WsVBYA1OY99z/9EWUYRzR8ZytH/ruTE1rxaCxrOFt8lmFm/DqKsyInWrLoi0nieGhrh8chIUs3mSTRqG8DxTbnlHivOtvH949s4sTmPtiOjue7ljpUON1yZ2pSCrSfISbUTPrwta/NaeJ8wgnlmJ1q0PMCBl5f5tg+I0lOabyfzYDEoJGImd6P1sIaXu/5Cv8/17XfgUqtpb2VtsRY6+OmF3eQlWShKL6Os2Ol7TpJAa1Ljsrt98xRKcmxsX5xMj1uaYA7VcXh9NofXZ3PN7d75VpIk+eZA5Z3wrhmiNXmrEhGJfmQfKSEo1kjfOxIoK3bQcVwMen8Ngx9sSX5KKVmHvAklzvW+/56dCJ0Saf9Rd7b893cSW7Xgs0/mMnbs2EvxMgmCUAMiaKgFhYWF3Hv/TBb/sIT4mf0JG9iCv50SA7nwG4IxUMugB1r4/j68Ibvc89u+T2bb96d7IuI6B+F2yWz7PpktC04AMOWD7oQ392PDp0eIaR9EdNuAcukiV2c1Jym/4rkHtvCWOzDaQIuBEShUEr+/fxBZljEGalFpFPiF68674lrfKweuEhuSSsHR934nfGAL1P7nTq15Lkff/x17dgmd3r8FU5PTrWRVvRZKvRqFTs3upxZhiA2maFcq+ugAOr1zE8bGIeT9cYR9KzLoNbVpuf3OPt7FVEZ0ZjU688Vltqgvdv2Uyg//3oXb4cHjlulyQxxj/lU+rWNlLd4tB0Wy8LFtLHt5Dy0HR5B1qIR1Hx1CoVLQ/97mrHn/EM16hZLfo+KiiUlfbiB53t8YYoOIvaX8OG9Jkgjt15yDs39Fdntoend/Gl3fGdntwVVqR6lTo9CoWFdQ9TWdz3t7OVrz6/v3uj67HO/Pgd8z2fdbOp0mxNJmWBQB0QYOuhvjttpxWezESakgexM36ExqUnYW8M0Dm/juka0oVBKBjQyU5tspiGvN6qxm5cprODnEyFroXY/mupc7krankLguwXS8tnxP2eqs5qAB2kBG+dtZtTSBBhL+NZLsVfuZNOVGrht3He+/8y7+/lfPkElBqC9EytWLtG7dOibeNBmijTR+aCDa0KrnANRkDHxVrKkFWJPy6N+uiNXvHmD/qkzfc/Fdg8lPLqU4y4ak8LYeXfdKJ9oOj66VG7o9z8KhN1aQv+kYnPFpuevbPhUCh4Zegdj/8jIsR7Pp+vG0Wjle/pYT7H78e9q8cC3BPZqeewegLL2Qbfd+hSkhnIihrQnpk4BS563EZ/62l4Ozf6Xnd3ejCaqdBfvqerjEhajsc3b2dWQeLOKjGzfQvG+Y7/sS3TaADmNjsAwYWu3xZbeHpK/+JP2HnTiLypCUCoJ6NKHZvQMo2HqCox+sRWXS0u2z21Boyre9HHztV0qT8uj4zk0VejWcRWWk/7wLZ6GVtEXbAOi38pHzvn7h6lbZd7aq394Dr/xCaVIunedMqfGxnTY3TpubP8taI1WRvCNn3SH2/ftHAG56t5tv3ZBTC1cm9g/HFteU0L7NMcbVzmRme04Jx95YhZRu5bt539K3b99aOa5weYiUqw2f6Gm4QE6nk2eee5b/vPEfYm/vTdS4Ducc9nAxFWpNgIHifel8+dRBivdllnsuaVsBsts7SFT2eOv1Bxxx5GQ1ueDznUkbbKLti+PxOFx4XB623/c17jIH+wO7cCCr7uc51KbCHSkgwZY7v8AYF0zY4JYEdYlHUl7YdQZ2igMJ7LmWGu+jjwqg95L7Kn0uqGs8AAXbkggfXP28lpq6kM/lxQQa5zrf+VSIqtpmYPghNn+XhClEy7UvdGD/ql8ByDxqZdkre7mmzyAUqqonGktKBfG39iZmcjesyfmo/fWkLtzKzke+w5ZZhEKjwp5rIXn+ZuKn9Cy3r8vqQGXSVfg9cJc52XLXFzjyvYsA+rWOwnL0PJpcBeGkmn5nZVmmcGcyof1rni75zGNL1fzsSSdTPpuahdGs9+le1P73NMcUrGXjsiIsi7aRMn8TLZ8YSUjvhBqXoSraUDMtXhpH+tIdDBk+lEcefoTnnnkWtfrK6BkVhPpOBA0X4NixY1x/40RO5KbR7r+TMTa+NOOPPQ4X2WsOkrPmIPmbjld4PmJEW/zaRnNotrdCpIvwI2xgS8KHtMIQE1Tr5VFoVCg0ILtlgro3ueCKdH1mTgzHbXOij/CnaG862U8dwNwigtbPjkUbcv6ZpLJ/PwAytdbSZjmcDQoJZ2FZrRzvQlU2Lru2eplq4zirs5qz988NGBNi+KO0Ndf83BzZ5WH7g9+gCTBUGzCcSalVY04I58Arv5C95gCBXeKxZRYhqRREDW9N8ld/EdytMebECDwOFygkSo/l4N+u4iTm3A2HceRa6PrZdNIWbydj2W7MiVWvCC8IF8uanI89x0JAh9haP3ZI7wQCOsTgcbhYk90cSalgYPghFAoJy4ChtB8AbpuTA7N/Ze+zP9BkRl8aTexy0VmQJEki+tqO+LdrxAevfMqvK5bz/bcLaNy4cS1dmSAIVRFBw3n69ttvuf3OGQQPTKT1MxNRamuvhcPjclOw5QT5m09QejyXol2pFbYJH9qawM5xBHaKRRNoRPbIKNVKjPEhGOKDL3laOneZE9njwVFgPffGDVCbf1/r+39Zlinak8b+539i77M/EDOxKyHXJFTZXX8mj9NN6sItHJ+7kdC+zfFv26hWypf6/VbwyGhD68+KqfV1SFrYgBYkffEn+VtOENQlnn2v/oj1RB6B18Wd97FKk/MIH9SShIeGsH74W7itDqxphSh0ag7+5zc6/vdGNt36KY48b49SwqwhFY4hncy2pDbrSbh/EAn3V5wPIQi1KX/TcRQaFQHtL00mrtgbu7P7yUUcenMFzR8ZWuG3QKlT0+qpUawbdohj/1uHLsKP0H61s0ikqUkord+ZSNJHG2jTvi2ffPQxkydPrpVjC4JQORE01JDNZuOBWQ/y1Tdf0+TRwYT0alYrx3WW2Cjcnkzen0fJWrGv0m10Ef6E9k8kdnJXVKbyi69JComwAS0q3e9SSP52E7aMIlo+OeqynbOuSJJEQNtGtHhiJMc+Wsu+f/+If9tozC0jUeo1aIOMBHaJRxd+ejylu8xJzvpDJH3xJ7bsYmImdaXxtN61VqYWj4/k0Ju/se/Fn2lusRPUvQnakPoTQNQHJQczKcsoQhfunSi5+/Hv6fLJNEyJ4RTtSSNt0TZC+iQQcB6BnCbAQFlGERk/7wJAadBQuDUJgNJjOex8eD6OglKa3TcQhUZZaSXNlBAGQN7fx4gY2vpiL1MQzil/03ECOsT45kTVtsDOcST+3zAOvPoLSqOWpnf2rdADLSkVmBMjKDmYiVTDHr6aUmrVNLl/AH6dYrjtrhmsXvM777z9X7RasYCiIFwKYiJ0DRw7doyxE64l015As6eGo4u4uKwNskcm+dtNnPhsI3gqvvzhQ1oR0DEWpU6NIS4YQ2xQvVnYpmhPGjtmfUuTO/oSM6lrXRfnsirYnsyJuRtwFpVhzyvFY3MSNrAFLZ8cRcHWJFK+20zhrlRkp5vg3s1oPL03xvjaH7omuz3se/5HcjccAbw37tbPjkOpv7rH9bptTo59tJb0H3YCoDRqcZfaAYga256I4W0xNQtj+wPz8NhdtHpmDIZGNRvGl7P2IPue/4moMe1J/3EnkWPakfHzbvDI6KMDMTYNJahrPJEj2uJxuMj76xjZaw7iKrHh1yoKlUnLic834rG5aDyjD7GTu12y10EQwDu35o/r3sOvVRQhvZsR1DUeQ+ylWV05bfE2jsxZQ0C7RrR8clSFJA0ui52d/1iA5VAW0RM60+ye/rVeBltmEUde+JVIfRBLFy6mSZPamdMn1B4xEbrhEz0N57BkyRJuuXUKwYMSaXXHhAqZUi5EwbYkTny6wfe3JthI1NgOBPdsirFxSL0JEMBb1vQfdqDQqAjtl4g1JR+k062mV5PAjrEEdrwJgMNvryTj1z3EnZwEm/nbXgq2JtH07v4E92yCPvr8VxKuKUmpoPWz47DnlpC/+QSH/vMbBTuSCelZs+xMV6LiAxkceOUX7NklNLt/ECG9m7Jt5jxUJi2mJqHkbjhC5m/76Pj2ZBIeHMzOR75jy4zP6fLR1BpVpEL6NAeFhEKrQmnUkvHjyR4HowZbZhGmZqHoIv058t5qslbux1Viw5QQhtpPT/LX3oXyTE1DafnUKPSXYL6RIJxNkiQC2sdQllbA8U82cOyjdTS6oQt+rSJJmb8ZhVaFf+sooq/tROHOFPb9+8cLDmijx3fC2DiU/S/9zJa7vqDVU6MJ6HC6t01l0qJQK5FUCkxNL80CbboIf1q9MYETH62nfacOfPX5l4wbN+6SnEsQrlYiaKiC2+3mn/96mjfffosmDw0i7DyyT5xLQLtGtJt9PdowP/SR/vV2QnHpiVx2P7HIl5lJoVVhahYGkoShlib2NkRuu5OslfuIvam7b8K5IS4YhU5F9PiOl+391IaYfWlckz7biKlJaLmhUleTXf9YiC7Mj84f3IIhNhhHgRVHnoUmd/dDbdZhbhlJ8ld/kfnLHprdN5Ce8+/m76kfk7ZkOwkPDK5wPFmWKd6bTsHWJAp3pmA5mgMeGXt2CT3n38WBV34hf/NxgjrHY24RQfoPO8lZuwB1gJ6IEW2IGNbGN/m9+EAmziIr5sQINAEXv/aHINSEUq+m3avXA96kGsnzN5M872/kb90ABHVrTMqCLaQu3EpQD2+r/PGP12M5mk3i/w077/l6AR1i6PzhFPa/tIyd/1hA/LTexE7u5psDFtQ1nuK96bV4hRUpNCqa3DeA7DZRTLr5Rh56cBYvPv8CCkX9vMcKQkMjgoZKFBYWcsONk9i8dztt/zup1jLfnKLQqLypOOu5vL+P+QKGqHEdaHbvAI68u7rGmWeuVEV70nCXOQnte3rSX8GWEwR2jLvsAaA6QE+Lx0Zw4NVfyFqxj7izFhO7WpgTI0DC12vgLPZml1KZtBx8bTmSSoHs8mBuGQl4K1RRY9qTMn8z5haRBLRr5Bt26CiwsusfCyg9novKT4farCOwozf7TPOHh6DUqWn9bPlVaSNHtydv4xGCezdDZSi/irZfC5EhSahbCo2K+Ck9aTS+I8X7M/E4XYT0aoaj0Erqwq2kL93uXeBHlsn5/SBlyfm0emYs+qiA8zqPJtBIu1cmcOLLPzkxdwPFe9Jo8dgI1P56jI1DUWhVHPrPbwR3b1wri2dWJax/IsbGIcx59mO2bt/Gd/O+JSAg4JKdTxCuFiJoOMu+ffsYOXYUtlA1bd6ZWGHi8dXEY3MCEDaoJc3uG4jH5iT9510YYoNQaq/ej46ryFsh1YV5W/WtKfkU702n2X0DL3tZJEny9v4ApmaXptu/vpNlGUkhYcsuBryT0Y/OWYM60OCbJJ5w/yDCh7UuF/A2uq4z+X8f5+DsX1HoVFzz4wNIkkTa0u2UZRTRbvb1BHSIrVG2LJVBQ/iQ2lk3QxAuFZVJ51vrBbwT/JvM6EPMxC4censluWsPodAocRSVsf3Bb+jw5qQaz/s5RVIqaDytN/6to9n/8s9svedLWv1rLAVbT+Cxu4i+rhMqP30tX1lFxrhg2rwzkV0vL6dDl4788uMyWrZsecnPKwhXMtFnd4Yff/yRLt27IvWIJPG50Vd1wACgPjmUwtjEO89CoVMT0CEG64k8SpPz67h0dedU6oCju5Qc2a5k2yOLUIWHYGtee1mSalQOj0z2mgPs/udiVCZtg+i9uhSyVuyjYGsSfi0iOfLeav6c/CEF25JIfHQ4xz5aD3h7h87uIVOZtHR4cxJ+raNQGbRIkoTb7iTjx51EjmhDYKe4GgUMgtDQqf30tPrnaEL7NcfjcKMy60CGXY8upCy98IKOGdQ1ns5zpqAJMrLjoW/RRfhjahZG7obDlSYAuRRUJh2J/x6D1D2Czt268NNPP12W8wrClerqbS4+gyzLvPHmmzz+zydo/n9Day2PdEMXNbo9Ae1jfOP2i3amUrgtmeYPDcG/VVQdl+7SOJF67tZ6T7w/6qgtZL8xF9lmR2E0EPaPGSiMek6knm5Bi2+UU+vlc5c5yN90nMJdqd7FwvK8qwsb4oNrZZJ+QyS7PWiCTWSt2Ic6QE/UmPZEjW6HLsKfE3O9CQeyVuxDdnlo+VT5VMHOEhvFe9PxaxNFxs+7yPh1Dy6LnejxneriUq5Ip75Tl+L7UNl5Lse5alKGuizHhZAkieYPDSF3w2Fs6YV47C4Ats38mpZPjSKoS/x5H1MX7keHNyZx9MO1HPtoHdpQkzfznMuN8uRQzkv9+ZCUCuJu642+SQjjJ07g1Rdf4aFZs+pVwhFBaCiuzlrGGZxOJ3fPvJdvFi4g5NE7KG0SQ2nFNdVqpCHdIGpCUiowxofgKLSS9OWfWI7lICkVaEJMyG5PvZ3ADTW/EdUkSDibQqsh5K5J5H44H33XdpiHX4MqoOIE5HMd+3zLZt26h6JvluDItaD20/vG7QO1uhZEQxM5oi2RI9ritjuRlApfj4Ise1cutxzOBsCaVkDWqv04i8ooOZSJy2LHVWIDoHhPOsX7MvBvE0272def91juhq6mld3zrRSfuf2l/E5Wd4yaXsvZqitnTct4PuWoD/cPlUlHaL9E8v48irllBCX7M3GV2Nj9xPfE39qb2Bu7nffvvkKjIuH+QQR2juPQGyswNg7xTbK+kM/H+fKdo1koIf8Xyz+ff459B/Yz5933UKuv7jTVgnC+rup1GgoLCxk34Tq2HT2M/wO3oAoOuKjj1Ycf/XO5kBty3ueLsfz+NygkX7dycK+m5VZPPvv45/NaVFWmmh6jNioZdcVtKaXwu19wpGbhLi0DtxtJq0EdFoQqNAiFUY+7oBhnWhb2w0no27dA2zyewgW/+o6R+OgwIoa1qcOrqJ+Of7aR5K+86U71jQIpSy0AQFIrMTUNQxNowONy4ywoxePy0OE/k1D7X/qx1rXtclR+z6WyczTk72V9URf3FJfVwY6HvsVVYsO/XQzZK/dhah6O5XAWAR1iafH4CLTB3rlC5xtEuix2PA6Xbx2Hi/ntv9DPlyu3gKJ3vqZT0wSWfr9ITJC+jMQ6DQ3fVRs0pKamMmDoEHJ0CvzunIRCf/ErSNbWD3x9u9lmvzGXsl0H0bVpjiYhluLFKwEInTUNQ4dzr0Z99utS366vLsgeD2W7DpH/xWJkhxN9x1ZIOg2WFX8AoDDqkfRaZIcLVaA/ymB/jD07Ium05Lz1GZzMagWg9Dfjf91QWo6LRGWq+nN8scFZQ3Pojd/IWLab4N7NUPvpcDfrgMfuoPD73wi8YTj2w0lYNm5FdrgIuvVazP1O56e/HK9Jdd+DS9Uaf+q44jvYsF2y1viTIrXH2Dbza+96PAqJwu3JxE/tRcqCLSBDi8dGUBxZ9eKetVHpv5TBrqfMRvFH3xFq87BmxUqio6Mv6nhCzYigoeG7KoOGffv2MWDIEBzNY/CbMhZJWXspRC9Hq97l5rE7yHzuXZzp2YQ/eTdZr3+CQq1GdruJeOpuNDGRdV3EBkN2uyleto6SVX/gLiwBQOFnQlKr0LVpTumGLeUCgqiXH0Ed6f3clO09TM47X6FrHk/gLWMp27aPgm9/9m2rMOqJevVRlKbaTWVYF0HFxQ7tOuXYYTOSVoMrMxd3UQnOzFzyP18MgMJsxNS/G+YBPVAFVb/Ke01b0i/H0BtBuBys2/eR8/YXNH9kKOlLdyC7PNjzLL4hfdoWTYh4/M4q96/v90LZ7ab4yx/QHEphzcqVIrPSZSCChobvqgsaNm7cyPBRI1EP6IZ53CAxGaqGPFYbKfc+i6TVoG0ej233IQACrhuK/9jLn2q0IXKXlJL9xlwcx09PmlHHRqJtHIMjOR3H8VSCpl6LwqCj9O+dIEPIXZNR6LUULvqNoh9Wo2vdjND7p6DQeXsUPDY7jhNpZL3yEQDhT9yJLrFJnVwfnP/48UvdC+VMzybzpQ/wWKzlHg+4fjh+Q3sjacSYZkGoSs6cb7DtPUzIzJvIfv3Tcg0aKBXEfvg80jnW7anPPc2yLFOyZCWutVv49edl9OrVq66LdEUTQUPDd1VNhP7hhx+YOHkyxhuGYR54dS6CdaEUBh3GPl0oXb/FGzAoleDxYOzZoa6L1iDYT6SS+ey7vr8NPTrgN7gnmqaxFP+6HsvaTehaNsXQqTXKADPGHh182zozcyj6YTUAYY/chnTG6qYKnRZVWDCaJjE4jqXgKSlfOb7czrdCcKkrEO4SCx6LFfOgnpj6dcWZnoMzKxdT3y4iYBCEcwi8cRRps17CU2gh5I6J5H7wLQDmYdcQdOPoGh2jPgUJZ5MkCb/xQyjxNzNw8GC++/Zbxo4de+4dBeEqddUEDfPmzWPa7bcTMON6DF3EpNHzJcsy9sMn0LVsiiM5HU9pGcY+nVGFnt/CP1cbWZYpWb7BN4xIFRpE+ON3lpt0X7pxG4Zu7Qi5e3K5gOCUkpV/+v4//Yk3kJ1OZLsDVXgISpOBsl0HCbxxNI5jKVxlHYfnpI4MQxUeTMmqPyn9eyeBN44mYNygui7WZSG73TjTslCFBft6pgThfCj9vBOeZZcLU58ueGx28j9bTMnyDbhzCwi8ZSyqwOqH9jUE5oE9UPqZuH7SJD7/9FNuvPHGui6SINRLV0XQ8OFHH/HAQ7MIvO9m9G0S6ro4DZJt9yFcmbkET59A1ssfAhA4YXgdl6p+89jspNz9jO9vv9H9CbhuaLnAQPZ4cGbkoAzyx51fhCok0Pu4y03eJwswdGmDeWAPFEY9stsDbjey00XJyj9wHEvxHceyYSuSRo0m7spcP+NCKf1MRL3yfziS0ilZvp68/30HbjemvlVP4qzPZFk+ORdGrlBZk2UZd34RjhOpFC5egTM1Czi/VuEK5/N4sKzbjGXdFrRNYlCFBKIw6jH16XKxlyI0AJJCAZKE7PSu2WDu3x1nWhYlK/7AunUv1q17CZo+AXO/hvl9OpOhSxsknYZpM26nxGLhzjvuqOsiCUK9c8UHDf954w2e/NfTBD44FV1i47ouToMku90ULPgFbbNYtM3jCX1wKtomMSj9zXVdtHrNmZnr+//wx+9E16LiXAPZ5gC3G9uug6T936tIei2Gjq1xFRVj33sEhdGAoXMbAsYP8e3jKiz2znlwu9G1bIp1616cKRkE3jwGdVjwZbm2hkSSJLTx0WjunISk0ZD3xVI0cdENMsDKnfMN1k27QKkg5K7JaJvGogoOwFVQRN7HC7DtPVJhH0OnVhd0ruJf11O8fD3ugmJ0bZpj3boHT2kZssOJKiQQXcumF3s5QgOgiYvCun2fb0hv4MSR2I+m+Bot8ud+j6fUiv/IfnVZzFqhb9Mc6cGp3P/QLCwWCw8/9FBdF0kQ6pUrOmh4/sUXeHH2qwT9321oGzeq6+I0WCWr/sSZmkXEMzORJAlDxwurhFxttPHRxH32SvUbKRVIWjWy3QmAXGandNNOcLkBcBdbKuyiCvAj5p2nvdvLMvYjyeDxoG0eX6vlv9JIkkTQzWNwnEgl95MFRP37wVo7trvYgv3wCTTxjS56vZeqONOzsW7ahbFnB8r2HCb3/XkAGPt0QRMVhm3vEYLvmIgmJgIkBRlPv4V5UM/znhgvO5xYdx6gYP4yDF3a4DeqH9p47++nLMtkPP02lnWbRdBwldC3b0HJqtNDJCW1irCHp5F63/O+xwq/++WKCBoAdImNCXpkOk8++wyl1lKefuqfdV0kQag36u+Svhfp+Rdf4KXXXiPo0RkiYLgIzqxcChetwDSgu6/iIFw4V2ExpZt3U/TLOgoX/Ub2m58RdMs4tK2aeQMIjdoXMADnnDMiSRK6hDh0iY1FJrAakDRq9O1bVMimdDHK9h0h7bHXyXnnK9IeeYWSNX/jcThr7fin5H2+GFV4MEHTJxDzztPoO7cGwHEiDV0rbwW+aOkqMv71XzKefgsUCsyDep73eXI/WUjue1+jjgoj6NZrfd97yx/byX33KzyWUpzp2bV2XXVBm6Qp90+omrvIgvKslMRKk5FG//0n6lhvuu3Qh6bVQckuHW2TGIIevZ2XXnuNF156sa6LI9SSkpISZs2aRVxcHHq9nl69erF58+Zq9/n6669p3749BoOByMhIpk+fTl5enu/5//3vf/Tp04fAwEACAwMZPHgwmzZtKneM+Ph4JEmq8G/mzJm+baZNm1bh+R496l/Cniuyp+HFl1/ipdmzCXz0Nm+rm3BBPHYHOe98hdLfTOD1Yv7CxfJYbWQ++45vfYZTNHFRGDq2BNmD/eAJAHRtmxM4cYRYA+MScJeUoqjFtSwcx1KRHU6MPTtS+ud28j9bjG3vEUJn3lwrx5fdbnLe/Qr7weMETByB4mTWJ3VoEGWAMyUDj8OFMsgfV/bpm5mkVCDpzq9CfKo3I+iWcZgHnw44nNl55H+2CGWAH8qgAN/cmyvFqcDBHueo45LUP66c/Erfb6WfCW2TGJzJGeS8+RmmgT0ImDAUpbF214mpK5qYCAIfmcaLr76KQqHgycefqOsiCRdpxowZ7Nmzhy+//JKoqCi++uorBg8ezL59+ypd4G/Dhg1MnTqVN998kzFjxpCWlsbdd9/NjBkzWLzYu97PmjVruPHGG+nVqxc6nY7Zs2czdOhQ9u7d6zvm5s2bcbtPNwbu2bOHIUOGcMMNN5Q73/Dhw5k7d67vb42m/jVoXHFBwyuvvsrzr7xM0P/dJipcF0GWZfI/W4wrO4+If81EYdDVdZEaLGd6Npb1WyhZswnZeboFWtMkBmQZy+q/kGUZXcumBE0Zh75DS1SBIu/0peIpKUVpNp33frLDSeHSVbgLi/Ebdg2a2JNzIhQSuN0EXD8MWfbgKbFi3bwb644DNVox/ZzlLS2jbPt+lMEB5VatDpw8Cm1CPDnvfEnWqx+h0GoIuG4o2sTGOJLTKfjmZ9IefgVVRAih909BEx1+7nM5nCDLqBud3tZtKSX9H68BoG0aS8hdky76mi7W2b0DtVXZ1yZpROBwFmdGNsaeHSt9LvDG0Vi37MFjsVL6x3ZK/9yOsWdHDF3aYFmz6eTwuIY7l1ATG0XgI9P598svo1Qoeewf/6jrIgkXqKysjO+//56lS5fSt29fAJ599lmWLFnCnDlzeOGFFyrs89dffxEfH88DDzwAQOPGjbnrrruYPXu2b5uvv/663D7/+9//WLhwIatWrWLq1KkAhIaWTzv8yiuv0LRpU/r1Kz+kT6vVEhFRvxu6r6ig4Y033+TZl14k6JHbTt/QhQtiWf0XpX9uJ+TuyWga1e8PcX3ksTsoWb6B0i27cSZnoDAaMPXtgrFnBzKffRddq2Y4UjLA48Fv9ABvyj+zsa6L3aCV7TqI/XgqmpgIlH5mQEYZ6I8yyN83dEt2uXCkZKJtFntex3YkpZH/zc/YjyR7F9/buA1DlzaE3ncLxp4dKFq6iuw35xI8YyLq8GBS7nmWnLc+I+aD5y463al8soUq+NbxKIz6cs/pO7Yk+I6JODNy8Bva25ciU5fYGG1CPIXf/YJt/xEy/vkWARNH4D+ib7XnUkeEoG0eT/bbnxP+6Ay0TWJQGA2+dUAk/ZWRutUe56hyWFJVgcOZ218tgYXH7sCdX+Rblf5sCq2GwBtHkfe/BQRNvRZnejalG7ZgWf0XANZNu4h87n40ceVbcS9V0HcpaOKiCHpkGs+8+AIajYaHZs2q6yIJZyguLi73t1arRaut+Dvlcrlwu93odOUbQPV6PRs2bKj02L169eKpp55i2bJljBgxguzsbBYuXMioUaOqLI/VasXpdBIUVPnQYofDwVdffcXDDz9cYUjxmjVrCAsLIyAggH79+vHiiy8SFhZW5bnqwhWzIvSnn37KPQ/cLyY914LSv3aS++G3mAf3JOjmK2ehm/MdgnChNzaP1Ub2W5/hOJ6KvnNrjF3aom/fAkmtwn4kicwX5gBg6NKWoCljRRaqi+QutpD/1Q9YN+1C0mmRbfZyz0t6nbeHR/YO2ZFdbiKfmVmhIgNgP5aCdcseyvYcQqFRo24UgdJsouiXdSj9jARPuw5ty6YULvyVkuUbiH79MVQhgTiS08n93wLceQVEzf4Hruw8smZ/jEKnxTywB34j+yKpLqyNpmDhrxQvW0fUC7NQR9XsBnLmZ7eYE2TP/h+e0jL8rx2MeVDPagNUj81O1jMfoA4IIuSJKYC357Fw/i8Ur9xI+GN3oG0WV+UcmvOpXFdVca9uv9qscNbH+Qx1WYGWZZmiH1ejDg/B0K0dtt2HyH5jLuFP3FVlj0HeF0uwrP4LTdNYIp++F9njwX40BVWgH5nPv4/s8RD5/IOoAk73nlb2utfnwAHAfjyV/Nc/5YN33mX69Ol1XZwG6dSK0L2X3ofKeHENEK5SOxvHvVvh8WeeeYZnn3220n169eqFRqNh3rx5hIeH88033zB16lQSEhI4ePBgpfssXLiQ6dOnY7PZcLlcjB07loULF6JWV7446MyZM1m+fDl79uypEKAAfPfdd9x0000kJycTFXW6cXv+/PmYTCbi4uI4fvw4Tz/9NC6Xi61bt1YaBNWVKyJo+P7775l0882EPHRrpWktG5JzdY9f7E3u7GOffbzSIwdI+/pjjD3aEzzjhkoXG2soavJanev1OJ99T8l67RMcx1MIe/i2Slu0XYXFyFZbpRXACynz1S7n/XlYN+8m+I6JGHt2wJ1fhMdqA8k7HtuZlo2k0yAhIbtcqGMj0bdqVuE4ntIyUu77NwqTAX3b5iCDIzUDZ0Yuxu7tCL51PB6bnfQn30ATH41t/1H8hvQmcLK31cldVELaP17D2LsTQVPG4crIofi3jVjWb0HTKJyACcPQtW1+XhPWZVkm9YEXMPbqWOO1Fir7DJVFWyn87hdKftuIpNUQNGUcpms6V3mM4o9+pWjHZuLuegh1gLfFzON0kPTZWzhTM1FHhxM0bTy6hPhznhsu7ntWnYv9LtTHoOFsl/P7XrL6L/K/WAKAOiYSV1YumqYxhD9yO5JKWek+aY+9BpIC/7EDQZbRxEWhaRRB5gvvYz+SjCoylLCHp6M+mdShutf8fK/1fIOPiw04bfuPkvvWF3w3bx7XXXfdee0rXJqgISUlBT+/MwLSKnoaAI4ePcptt93GunXrUCqVdOrUiebNm7Nt2zb27dtXYft9+/YxePBgHnroIYYNG0ZGRgaPPvooXbt25ZNPPqmw/ezZs3nllVdYs2YN7dq1q7QMw4YNQ6PR8OOPP1Z7fRkZGcTFxfHtt9/Wq89agx+e9Ntvv3HTLbcQPuEWTPoWkFTXJbp4l/JGdq5jF2xcgy46luDbr68XAcOlvqlfzPErvWHlZGHbexhT/25VDoFRBfhBgN8Fn/tqHCZRHUOXNlg37UIbH40kSd6UpyeXq9A0ioAapgh25ReCLBP24K3l3jtZln0V/cJvfsJjsWLbc5iAG4ZTuOBXJK0G/2sHo/Q3E3D9MAq+/hFkmYDrhxE0dRymvl3I//pHst+Yi7ZFE8IemFrjOUKyzYHHYvUNOzrlzM+ALMs483LwOBzIHjeugCBUpvK9V/o0A4qbxuA3oh+FC34h75OFaJvFoY4IqfS8hht6Ydm/m9TP5xAz/T5Ufv4o1Bripz9M6bFD5G74hfzPFxM87ToUfqZzrg9Sl5XzhhAYVOdyfd8daVnkf7EEXZvm+A2/Bsu6LehaNSVgwrAqAwYAXYumWLftRdsslvQn3wCXG22zOJxZeZj6diFq0E1gpdbuzed6P8/n/T7fz4bW0BLlhFu48eab+emHHxgyZMi5dxIuKT8/v3JBQ3WaNm3K2rVrKS0tpbi4mMjISCZNmkTjxpX3or388sv07t2bRx99FIB27dphNBrp06cPL7zwApGRp+fNvv7667z00kusXLmyyoAhKSmJlStXsmjRonOWNTIykri4OA4fPlyja7tcGnTQsHnzZsaNH0/QqAmYElvXdXEaPEdeDtYTRwgbMR5dqv7cO1TjQm9uDfkGL3s8OLIzAbCs2URo5+Go/QIu6Tlr0nJ2IUNAGhJJ7f0ZUxgu7jPrLvCOjZXdrvLHP6NnwJVTcOpBzIN7gSxTuHA5rvwigqeMw29IbySNhvzPF2P5/W80cdFEPDOTiH/eg233IXLmzCPz5Q+JfPpeb3rdc1DotZj6dKb41/WY+nVDeUbWJ9njoWDj7xRu+RNXYf7pnVQqTH27ENpuMOrA0+NqtUkaiPMjeNp1lO0+hGXNJgInj6z0vKpAPxrdeg8pn75L6hcf0GjavahMZiSVClPzVliPHqLwr3VkvjAHVWgQkS/MQqG9/N/dq23i8oX+PtbkNbKs8aaJDGjRDT9zG/T3Nq+2HKeOGTBhKGV7DpH7/jzUESG4C0qQtBpvwoGy2htWUV/uDabE1nhGTWDstdeybs0aunZt+KthX22MRiNGo5GCggKWL19ebmLzmaxWK6qzhpUqld4A+sxBOq+99hovvPACy5cvp0uXLlWed+7cuYSFhVU7J+KUvLw8UlJSygUm9UGDDRqOHj3KkOHDMfcdgl+7qrvZhZqRZZmcX5eiMvvj1+HifwTryw/8pSR7PJTs2U7+2hU4crNRGoy4raW+55X6yz+x+WJb2eqyAnYhnxl7nAPp5NhSd0kpyoDzmx8iu93YDh7HfvgExT+vRR/bGLOyCYoqXhtDlzbY9hxCGRyAQqvBf/QAlIF+5M1djONoMvp2LTAP603ks/eR/cZn3gnUny0mrOtIAvzbIvcfQc4vS3AVFlfZOi+73aBQ+IIV86CeWNZtIfW+f6MM9COgfU/KTGYsB/ZgPXoQv47dMLdqh7uZDkkhUbbrIMXLN3B87SZMzVtjbtMBbWQjHHk5cFDGak/FY7WdMx2rOjCIRrfeTcpn75P62fvE3H4/Sr03aPFr15nCv9Z5X/eiErJf/4To8beDsX4FDlfD71B1ZI8HSaGo0XwuSalAodViiPcO3atpi749zkTYA1PIfPkj73wipRKFUU/c/Y9X6PE6l5pMQq8P/Np1xm2xMHTECLZu2kSTJg17WPTVYvny5ciyTGJiIkeOHOHRRx8lMTHRN0fliSeeIC0tjS+++AKAMWPGcMcddzBnzhzf8KRZs2bRrVs333yE2bNn8/TTTzNv3jzi4+PJzPQ2HJpMJkym0z3EHo+HuXPncuutt1YIRCwWC88++ywTJkwgMjKSEydO8OSTTxISEsL48eMvx0tTYw0yaMjNzWXgkCGoWrQlsOeVsQplXSv8ewOlh/YROWkainqYG7i+ceTncuLtl8o9pjT7ET7+RlQmP7QRUfVieNf5qo05HRdzzPOlTdKg1jYjG3Akp1e5LktVZchfv4rclT+j0Okwt2xP2OjrUairKO+mXPI/83Yru3MLcObkow4NwtS7M5rYKAoXLseyfgvWHfsIuG4oUa88QsnKPyj+aR3H/9hOyJDRqE72POnz/FCWna4gyU4XzvQsCpesxLb3CJJajdo/mNCho8lZ/I2vCO6CYgq3/om71IJSpyd6yp0YmyYCp98PTWwU5iG9sazdROmfO8hY8GX565Ak/Ib3wa8GK/hqQsKImXYvyR+/Tfr8z4kYNxF1YDC66Bg0oeE4crJQB4ZgP5xE9o8LiJpcNxNE61ulsj6wHNxL5uJvMCa0IHz09Si0umqDh5C2gyn+dT25q38hYlzN0+p6e7CiafT2U6Q+/ArKADNlO/aTk2cj6sbbz7vcDeW9DOzVj7ySQgYOGcLWTZsIDq5+iJ5Q94qKinjiiSdITU0lKCiICRMm8OKLL/omNWdkZJCcnOzbftq0aZSUlPDuu+/yyCOPEBAQwMCBA3n11Vd927z//vs4HA6uv/76cuc6e0L2ypUrSU5O5rbbbqtQLqVSye7du/niiy8oLCwkMjKSAQMGMH/+fMzm+pUopcFNhC4rK+Oafv04bnMRct3NDbJiVt/Ys9JJ+uANArr3IWz4uLouToNw6JmHff+vDgrBv2M3Anr0vWoDrrpsHSxLTSLlf2/TaPpMDPFNa7yfIzeb5P+9jT62MVE33nbO3xK3tZSjrz6Nf+ceWA7tw9SiLeGjJ5TbxqJNI2/uIuwHjxM0ZRxBMd05+tozvvU5tJGNsGekEj31Ll9l32O3c/ztF3GXWgDQRcdiaNKc/PUrK5TB3LYTkdffguzxgCwjKU+PNa8qiJO2l+DMy0YTGuFdU0IGtX9AtftA+fev9PB+0r76H4Cv7KVHD5L2xYdIag26RrGEDR+HNqJiRirh8rNnZ5L8wX/QxcRjS09FH9eYRrfcWe0+hX9vIHvZIiJvmIK5TeXrMlR7zjgHWa9/giu3kMBJI8id8w362KZETZ6OoopMMw2d7PGQu+hrmujVrF+zBr3+4oZIXukuxUTooqKiGs9pEC5eg6pxezwebrzlFo7mFRA8brIIGGqJq6QYPB4CuvWu66LUa7LHg9tmw3Wycmdq2Zbmz71B4wefJKjv4Ks2YABvBfPsf5eL9chBJI0WfUz8ee2X9/uvKPUGIibUrPFBaTCiDgxCoTegMvuDx11hG5M9mtgb78PcpgP5Xy4lbd4nyE4nxubeydiOnEy0UTGkff0xtoxU4OTn6uRnCgCFgpDBIwkdfq33vEYTERNuJnLSNEJHeB+TFIpyAQNUHaRpgoIxJrREHRCI2i/AFzCc2qcm75U+/nS2KUnyvlYqk/dGLbuchI++XgQM9UjOr0tQBwYTfcudqMx+qIxmHHk5ZP3wHUdfe4akD9+kaOtfuCze1elL9mwne9kiDE0SMLVqf0Hn1CZpCLp5LO7CYko3bid01q2UnThCzq9La/PS6hVJoSB43CSO5OZz4y230MDaYAXhvDWoWvezzz3HynXrCbl+6hXbclEXbGkpSCoVSp1oJalK6ZED5K9bwdGXn+TY7H+hNJkJHVm/xhpejVzFReSvX4m5dfsKlejqeJwOLAf24N+5x3l97hU6A5b9u7Gnp/jmUpxNkiQib5iKsUUbyk4cRdJoceRkEjnxVqJvnkGjqXej0OrIW/MbAEq9Hr+O3pWeQwaNpNGt93jL6PCuNxF+7WT82nXG3KodKmP1K1lfaPB2rv0klQoUCsJGTcDQJMG7T3gkUTfeBrKM9fgRUWGqB+w5WWT9/D1lScfw69gNhVqNOjCY4h2bSf3iQyyH9uHXthMKjYasH74j5dN3kWUZj937WfPr0PWiGuNMjmhCZtyAdctuPBYrgZNHUbTlD2zpqbV1ifWOQq0h+PqprFy3nmefe66uiyMIl1SDCRoWLFjAq6+9TsjEaSjPceMUzk/xjs2Y23ZCaRArElfGcmAPaV9+5KvkAcTcdt8lz4wknJstMw3Z5SK4/7Dz2k9SKFEaTNizMs5rP110LM68HCSl8pw9G0F9BgEgO+y4raWYW7fH0KQ5Sr2esOHjKD2wh6Kt3pVzTa28KfqcRQW+BhG/jl1RBwWTt2rZeZXxUrAeOQAeD6ozeikA9HFN0cc3JfunhRT+tb5uCif4FGxYTfH2zehjm/gyChqbeYfBua0WGk25i9Dh44iZPpOoydNx5uVgS03C1Lo9mrAIMhfNI2PhV+USOpyvwNBOqGMiKNt9ENPA7gDYUk5c9LXVZyqjiZCJt/LK7NdYuHBhXRdHEC6ZBjERetu2bUy59VZCxt+ENrx+pZ+6Eig02vNqpb2a2NJSSP/mU8BbYbSlJhH/4FNogsSkt/rgVDaf3BU/EnH9lBovnGbPSsdttUDN11kDIGzkeAxNEjA2b1Wut1N2u7Ec2INCp0NSqnAW5GHZv/v0fmMmljuOX/suWE8cJXvZIvSNm2Fo3Ax94wSKt29CFxOPf4eu3iFQsnceRF1zFRcBoGsUV+5xpV5PzPSZpHz6LrbUE0Dfy184wUdl9kN2OjA0buYL8Ewt26Ey+2NMbFVukr++sbfHqPTQPvQx8cTd+yglO7eQ/etSkj/+L9G33HnBv3PGqARK9+xFl6JDoTf4es2uZNrwKEKuncyUW2+ladOmdOx4/vNCBKG+q/c9DTk5OYwYNRq/3gPFWgyXiDoomLKkY3VdjHrJWZDn+39bahIhQ0aLgKEe0TWKw9S6AyV7dpSfF1AN2e2m8O8NKI1mwsfWPEsMgKRUYm7dvsLwyNIjB8j47nPSvviQ1LnvkbV0PqUH9wJgbNEGv7YVKxBhI8YjKZWU7N7urcy5XchuN1lL558sqAySROmRA+fdI1LbjC3aIGm05K/5zZsS9izq4FDsJ9coEepO8KCRBPUbQu6qZRx5+Sny1vyGOiCQstQkMhfNw3riKLLLuw6JQqtFFRBI/rqVHHrmYTIXzcPUsi2xdzwIskzKx29jS0+5oHK4rae/i0qjCeuxw1fF8DVTizaYew1gxKjR5OTk1HVxBKHW1eugwe12c90NN+AMjSDgmoF1XZwr18nKiVCRMfH0asL6+KYEic9hvSJJEn7tOgHeRdls6Sm4y8qq3N5tK+PE+69RvGMzfu0612huVFlqEpb9u33j9q0njpD2zadkLPwKa9IxXMVFyC5vdiSFVktQn0HE3vkQAH6duhNx7eRKj6vQaDAmtKR411achfmUJR8HwNDYO+lYUiiIue0+lHoDGd9/5c2YVEdURhNBfQZRuGkDR199mrw1v5UrjymxNY7sTMpSa2nZX+GCSJJEyMARxN33GP6de5K39jdSv/qIwj/XYktLJnXuexx+/h/krVmOJEnE3jEL88l1jkr2bCfl03eRVGpiZjyAOjCIlLnvYUtLPsdZvexZ6WT99D2lRw5SsmcHgb0GABA2fBzWY4co2bnlkl13fRJwzUAcIeFMmDgRdyUBtiA0ZPV6eNJT//wnOw8eJvy2+2o87EA4f25rKdqIqLouRr1kSzvd0nZqgqpQv5xahTvpvdneCZ2SRNjI6wjo1htnYT6lRw56J3fKMpaDe3GXFBMz4wHfUBuP3Yb1xFGc+bk48nNxW0rw2MpOLrKmpOz4Yd+5FAYjHmspmvBIHFkZlOzehso/kMYP/ZNG02eSv2EV+etXkb/xd8A7ZOHUEKrKBHTvQ8on75C3dgUAxsTWRN90Ore9yuxHxLWTSf7oLXJ+XYomLAJjQgvU/oGUHt6Pq6QYlZ8/hqaJvt/Ikr07saWnYGiSgKRSY4irnYWngvoMwtAkAcueHeStWY5l/260UY1QGoy+iqU9LQX9WUOYhMtPGxpO2IhrUQcEUrR9E+Y2HYiYcAt5v/9KwZ/ryPt9OaZW7dGGRRA54WYiJ9yMLT2F5A/f5MTbL9LsqVdodOu9pH4+h7R5nxB39yOozNWntbSlp1K0eSNFmzei8g/Ev5N3PoMxoSWmlm3J3/h7rSwcWt9JkkTwmBvY/uk7/PPpp3n5pZfOvZMgNBD1Nmj44YcfePPt/xJ5+/0otLq6Ls4VzV1mFUFDFXSR3jSSYaOvFyl+6yljYmtkjxuPw4GuUSzFO7ZQ8MfvqPwDyPz+azx2m29bSaUiYvxN6GPikd1uclYto2jLH3jsdu+CaoHBqMz+KPRG8LhxFhcSOfFW9PFNsaUmYUtPQRcdizGhJVlL51O8fROa0HAkScIQ3xRDfFMc+XkUbFhF0da/MMRXX2HXxcTj37UXRZv/8D5QSW+CLjqWgO59KPzbO9FYodUS0LMf+WtXeHsJAVPr9vi164yzII/cFT+BJFGwYTUACc/+p1YaXSRJQt8oDn2jOAwJLSjevhl7ZhpuSwnayEaEDBnta7W+2piTve9DSWz9adySlEqCrhno6x31uFzkb1iNOjAYc+sOqAPLD7NUnly9+dSK5AqlkqgbbyPpvdnkrlpWZY/ZKaeO59e+C+b2nb0Zt04ytWxL5qJ5OIsKUPsH1uZl1ksKrY7gCVN446236dWzJ2PGjKnrIglCraiXQcPx48e56ZYpBI++Hm1oeF0X58rn8VwV400vhCPXOy5VBFX1lzY8slyCBE1QKElzXid93icoTWaaPvw0Co13ISFJocCenUnuyp+9cwWyMwnq1R+/Tt1RBwRVGxiaEluXm1cV3G8oxds3YT6Z+ej0+YMJHzuRsFETqk0wcKqi6dduApkt21J2/Aia8EjKUpMqtNaHjRxP6PBxeBx2cpb/QP7aFWgjGxF9yx2UJR0lc/G3WPbuBIUCfUw8oSOuJWvpfBw5WbitpedM1Xq+jE0TfYvTCd5g4dT7WV8V/rkWPB4irrup0t4gtV8AERNuJvP7rylYv4qgvoNRmcwE9RtKzq9LCOzZF214+d9B2ePxDUs7NVfC1KZDhc+GMaElKj9/MuZ/TswdD14VIwe0YREEj57AjTffwu6dO2jcuHFdF0kQLlq9CxqcTifXTZyIrlU7zG061HVxrgraqJgaj1u9msiyTOGWP5BUqnILV51dOahPrYsC5TIiuS0lZH7/NbLswRDfDGSZvDXLUWh16GLiCR0x/oKH76gCAtFGxVC49U/8O/eoWIzzyEgWoW5OyeBE0r+d6826JEloI6Jxl5USc/v9qP0CkBQKArJ1BHSciNx+ApJCSYlRwtyqPdrwKCwH9lC4+Q/Kko6R9tX/iBh/I2lffoTbUlLrQYNQ0YX8DlyuHgpnUQG5K38m8JqB1Q4f82vXGXtWBrm//+pdw8RoIqBLTwo2rqZo29+EjTi9Nk3hlj/J/e0H35DAU71elSUkUBqMhAweReaieXjsdgKyT68GfCX/fprbdMSRfJwJEyfy9x9/oBbrSwkNXL0LGp56+mmOZmQRPn1mXRelzlzuFqsSZRDZOfsu6znPZE6W692Nw5GXQ8HG3yne9jdx19yAf4YKqPx9udj3q75de0N16n3IW7cWtcGP2AceI+e3H3BbSnzBgux04t+5B6Ejri2XfvJCSJKEOjAYd2lJbRQfc7KMURWCBUCWsZ9cMfr4f/5N59v/U/7cCqVvn5JYiZRP3/Ve50mRN0z1JTeQK1m5Wqh7Z/5unO9vyPn+ZrhKigHvECJZlqtt6Q/s2Y+Cjb9TsmcHfh27kvXDd7hKLb6gALwNKgV/rEHXKI6gvoPRRcXgsdspSznhW/xPlmVchfl47N51Sgo3bURSqyv05tXW739tBGBnvg+19bscOGQ0Rz59l6f/9S9eefnlWjmmINSVehU0rFq1ijfefItGdzxY5Q29vncBN0Rl+Rm4HWUUfreEqE7DsDYrv0JuZa95bf7IV3aOc3X3X+z5qzq222Hj+Np5FCXvRaHSENPzOkISK7Yi16YL/Uxf7cFGZa+brSiHvCNbiek+loAcHcpxp1OqelwuZIe91hYxtBzci2XfToL6DqmV4wE06jaaqE7DyN63HkdpESqtHnNEs2r3MSfLBMW0I2f/RpAk/Lv2Qh8T71ugy5Gdia4Gaz1cigrTpdTQevxq895Vk2Od+XpogkO9i/D9uICibX8Red3NaELCKt1PZTJjbtuJ3FXLKPhzLa5SCyGDRuJ/ctVygLw1y3Hm5RA28joUGh3p8z/HlnoCj82GX+ceRIydSP76VeUWJlQHBRM95S4UGg1nN8CcK3CoyfOV/f+Z+5zv638h71dlZVSoNQSNv4nX33iToUOGMHCgyMAnNFz1JmjIzc1l0o03Ed11LCG2cBDBwWUT0X4QSAqy964je+86mo+8B3PkuSsqp1RVwa/sB7SmP8Tn2q62g0eP20X2nnXkHd6Mw1pMfL+bCIxvi0J1ca3Rl9KlDKpqW20GnlVdt60wm2O/f4naYCYksUeF4ytUKlCd309eda9x9pGDIMuYW7U9r2OeUtX3RqFSE9Hu/CoWMT3HE5zQBaVGj7OttzJo2bcLFAoMTZqfc/+zy1HfA4jKXrf6OBm5Lhu5znw9lHoDjabdS9nxw2T99D3p33xK7J2zqkwyEjbiWrI9bhRaHX6dupcb0uS2lpK/5jcACjdtoPTQvnK9EMVb/8LUog1Fmzdiat2BwF79UGi0qINCvN/Bi7yeU9dU2ePV7XM5VHU+M+HIXccwcfKNHNy/j+BgsdaP0DBJcj2ZATt+wvVs3HGEmP41X9VVqF3F6Yc4/MuHKLUG2t/03FWTLcjjdnFo2ftYc1Pxi04kustI9EFX7srjl7NSda6b9vmUpapjeVxOMnauJGv3GjSmQCJumlpuDsqFOmfZ049w6NcPUShVtL7+cTRG/3Me82IC6fN1dOVc7AobMdPurXa7mp6/Ju/VxVxLbXwWLuR45wqQrrTebVthNvt/eAv/Ri1oPOAWJKni73xVr53H4SD7l8UUb/vb95hCq0N2u0ChRHaUn98Qe/fDVfZynU9DwpXwHsiyTMqaL7mmQwKLvl9Q18WpE8XFxfj7+9N76X2ojNpz71ANV6mdjePepaioCD+/6tMBC7WnXvQ0zJ8/n+W/rSBh7MMiYKhDflHNSRh2J4eXf0RJ5hH8os7dQtnQ2YvzSNq4gNKcFBJHzcQUHl/XRbrkLnZYR20Ps6jpsIPKyLLMifXfUnhiN+HtBmAaNfjk8IeLK1ONtotqRrsb/8Xeha+y9/tXCWzcHn1gBGqDP37RzVFpK67PcHZruOzxUJJxDIelAH1QJLqAcBTKi/9ZlmUZS9YJghO6VH8N5/FeXuqW/No+/oWMlb8SKqfnogsII77PZI6t/gKVfgkxPcZXuO9W9RthS02ieNvfGJomEjZyPFlL5/sWJQRv9iS/6ESCm3XGXpJfo2FxZ5+3PvUU1SZJkojsMZ5fl/6H+fPnM2nS+a1GLwj1QZ0HDRkZGdxx511EdrsWtd5c18W56mlMAXVdhMsmZ/8fpPz9A2q9iYShM66KgKEydV1RupjzFxzbTsGxHTQeMAV1/46X7bynqPVmWk/4BzkH/iTv0GbyDm8GWUZSqghq2onwNv3QB0ZUeu6icAfp8+diPXrI97hSo6f5yHswBF9cT4nDUoDLZsEc0bTKitiFXn9Vlfva+hydK6g930DnSq2EXozAxu2I7T2B5I0LUSjVRHcdXW2D3anX3KRsirNNP7L2rMWy/Hea9byZ3amvoFRrcNutNBk0Db+oBJSak8OeLiAQPPs9q+vfp9qk1puJ7D6eO+68i379+hERUfG3QRDqszoNGmRZZtr02zBEJBDYuH1dFkU4yW4pAEDnFwqAo7QItcHviuoBkmWZ1L9/IHvvOkJb9ia66yiU6ovrKhXqRsGJ3RjD4s8rYKjtSoja4EdUp2FEdhyC7PHgslnIP7KV7H0byT+6jTY3PFnp0KXChUuwJSfRbMjtmCKaUFaQQdL679i/5A38Y1tjDI0lpHk31Ibz73p32bxpL1UGb0PMpah4Xa7K3MWe51yBQ0NYY+FSCG3RE9ntIuWvJbjsVuJ6X+/LylUVSZKI7jYGtdGfzJ2rcJYVYwyOQuMXQlSHIWj9gis9xvkGb/VxfkptCWzcHkvqXqZNv41flv18Rd1bhStfnQYN8+fPZ8Mff5Iw7v/qshjCGVxlJysbOhOFSXs4unIuhuBGRHUejn9MyzouXe3I3LmK7L3riOk5nrBW19R1cYQL5HbaKU49UONJw5e6YihJCiSlAo0xgIj2gwht2Ytd3z5P1q7fadR9TIXKVP7xnYS1ugb/2FYAmMIb03zkPeQd3kJh8h4yd64i98Bf3oqYUkVAbGsCG7dHpfNmfyrLzyD30N84rSXoAsMJa9nb95wl0ztkROdfeYacq825KqFXa+AQ1roPSq2BE+u+RanSEtPz2nPuI0kS4W36IXvcZGxfQWDjdpTmpqILqP6zdiGBwJX6nkR2u5b1S1/nu+++E8OUhAalzoKG/Px87pl5H+FdxlQ69leoG7LsAbx54DN3rQbA5bBydOVcWl77cKVDLRoSh6WA9G3LiWg/SAQMDVxRyj48LgdBzTqjraYls64qHkqNnrDWfcjcsZLcw5sJbNyOkObdMYbFIUkS+sAIbIVZ5fZRG/yIaD+QiPYDseank7nz5HfQZiH5z8Wk/LWEgLg2GEJiyNy5CoVKjdYvhILjO8jYthxjaBx+0c3JPfQ3+qAo0YN2luoqrldr4OCdf5BH1q7fiew0tMb344DYNqRv/ZW8I9vQB4bX+HxX42t8NpXWQHjnMdxz70yGDh1KYGBgXRdJEGqkzoKGWQ89jNo/isDGHeqqCEIlTrWGelx2yvIziO46CpXORNKG71Bq9OfY+/JzO2w4rEWodSYUKg32kjzUerOvxfVsOQf/QqnWetPMCg3byQwtpyo59bEyEtVpOAGxbShK3kvekS3kHdqEMTSOqC4jcNlLUVazwJwhKIomA27x/e0sKyH/6DZyD/5NwfGdBDXtREzP61Bp9eQe/IukDQuQlCpyD29CpTPRuP/Nl+y6zrVAWH13JQ9/uRD24jw8LgeO0sIaBw26gDBajHmAtK2/ENVx6CUu4ZUnsEkHLEk7mPXQw3z+2dy6Lo4g1EidBA2///4787/7juYiW1K9cypoKE4/jMflwC+6BcWp+1Fq9DVKKXk5nVj3DXmHt1R4XKHW0WTALWj9QihM2k1Qk45oTN6WHFthNoaQGNECewXQGAMAcJQWotfUzx4wSZIwhsZgDI0hstNQitMOkbbpJw7/8gFKjZ5mg6fX+FhqvZnwNv0Ia90Xt6OsXOUuJLEHIYk9LnllXva4Kck8StKGhej8Q9GYAglt0QtDcNQlO+elVB8Dzbqg1psAsBVkog+MrPFnyBDSiIRhd1zKol2xJEkivOs4vv32P0yfdiv9+/ev6yIJwjld9qDB4XBw+x13EtZ+CFpz0OU+vXAOpxYzcztsACfXapCQPZ7L2rroLCsh/8g2ilL28f/svXd4XNWZ+P+504tmRr13WV1yb7jRTAsEEiAhkITQsuymkLJJfpsvmw2bkEaypG4IYSEJhBBCD6EFjG0w7rZsS7J6712aGU2fub8/xh5bqM1II49k38/zzGPr3nPPee/MLec9b3NZRzBmFIMooo1JJqF4E6LoY7j5GLr4DNLXfxS33YzP40IdFUfPsbdp/Of/BfKFDzcfo/iGryEIAj6340xmD4kljdMyDIDP44qwJMEhCDJM6UUYU/MZaTlBVFJ2QJkNrR9hytXg4aajdOx/hbJP/AdylRZR9GEf7g3rhN7a10LDG48C4B4fRfR5sQ12UPTRr1wwdV3OR1LXXIPTMkTLrqexj/SStvYjkRbpgkBtiCVx5ZXcdc/nqas5iVKpjLRIEhIzcs6f8r/61a8Ytjgkf/JFisaUAIKMkZbjyNU6BusOoDbG43M7AhlZFhprbzMnX/wZXUdeR6ZQYUgrYKT5GGOdtbTvfYH2vS9S+8ovEb0e4gs3YkjJIzZ3FfEFGzCkLiNz880Y04vI2HADedvvxD7UhbXPHxgqyOSIPu85OQ+JhcM21E373heIyV2JLj4j0uKEhCCTE5u3ak4Kw0y47dZA7ANA34md1Lz8PzjNQ2EbIyopB4D4gg2svvMhCq/7ErbBToYaDoZtDIlzj0yuIO/yO4gv3MhI87FIi3NBkViyhWGznV/96leRFkVCYlbOqdLQ09PDf333ARLXXDdrajeJyKAxJZC97VOYO2vRxiQz1HAIn8d5au/CWhm8bic9x96h/o3foYlOovyW/2TZlXeTtflmVnz6vyn7xLdJLN3GaHs1XreDqORcjKn5U55D/lWfJ7F0K6bMEpQ6I+aOGsCfFco9Prag5yGxsIy2V1P76i9RG+PJ2vJJycXxFAnFmwAYbjyCz+Oi51QQdTgtAIJMjimzFEtfM6LPR1RSDrF5a+g6/AYex3jYxpGIDMa0ApyWIdw2c6RFuWAQZHIS117Hf333AXp7eyMtjoTEjJxTpeEb3/wWhtSCKSd6EouH6KwyEAS8Tjtelx2P0w74s9UsFEMNh6h89kF6Kv5JYskWCq65d1KxP0EQyNh4A8s/9R3KPvFtCq/9Impj3MwdiyIIMrxuv+KjjUnBPtKDKEq+zEsNp2WIoYbDNL/7JKb0Yoo++mUpNuUsZHIFscvWAND2wQv43H4XQ1+YLWspKy7HOTbASOsJANLXX4co+mj74Pklc1+ZmpyYmpyzN7zA8HncgD/2S+LcYUwtIColn29881uRFkXiPMHtdtPR0UFdXR3Dw8Nh6/ecKQ0HDx7kueefJ3G15Cu5mBmsP0jda/+LIMhIXXM14F+5N2WU0HngFQZq94V9zJ7jO2h976+YMoop+8S3Sd9wfdgsUaNtVbjHR4kvWA/4U66qDLHS6vQSo7dyJ1V/+yGt7z1DVFIuOZd8OhB/I3GGzItuBGC48UyCgPrXHwnEKIUDfWIWhtQCeo69jc/jRqkzkrnpJkZbT9C+9wVsw91hG2uhOa08SEqEn86DryJX66bNPiexcCSu/gh/e+45Dh6UXP0k5obVauXRRx/lkksuwWQykZ2dTUlJCQkJCWRlZfH5z3+eQ4cOzWuMc6I0iKLIfV/5GgklW6Xg50VO95E3sQ91kXfFXegTsxFkchxj/WRuvhmvy4GluyGs4433t9F9+A2SV2wne9utYffzdllHkClUaOPSgPC6akicG6z9bXQd/AdJ5Zey/Nbvkn/1vyBTSAGDUyFXacjcdFPg75xLPo3XZafz4KthHSd1zdU4Rvs5+fL/4HU7ic1dSerajzDaWknNS/9DT8Xbk6wOXrcTt90SVjnCzUIpEEtBIfG6nXgcVjI23IA2NiXS4lxwqA2xxBdv4Stf/fqSsdhJLB5+/vOfk52dzWOPPcZll13Giy++yLFjx6irq2Pfvn1897vfxePxcMUVV3D11VfT0DC3udw5yZ702muvcaKqisKP/3/nYjiJeRBfuJ6eYzsYqNlL7qWfRReXzlj7SRJLt2FIWYbbZglrFqXB+oOoDbGkrr5qQVb/1aYEfB4X7vFRVFExyFVavC4HHqeNnmPv4LQM4XXYWHblXUHVoTj75T+WJ7nGnAt6jr6JJjqZtHUfQRAkpW824grW03PsbQwp+cTmrcbrdtL+wfNEZ5djSi8KyxhRiVkU3/A1av/+C5re+SMrsm8mZcXlJJdfQs+xd+g++iYO8yApq65AqYmi8+CrDDUeRvR5SV19NSkrt4dFjrkylqeedSL/4f1zvd9P97PYnx0u6wiAtLAXQZLKL+H4iz/h9ddf59prr420OBJLiL1797Jz507Ky8un3L9+/Xruuusufve73/H444+ze/du8vNDDxVYcKXB6/Xyta9/g/iyy6RUl0uA5BXbcZqHGG46StueZ0kqv5jmd5/EPtJDYtnFNL39OF2H/kH6+o+GZbyxzlpispcvmAVAfyqzzvhAO6qoGFy2MeQqDe0fvMBIy7FAu6Z3/si6wrtC6tvU5Azp5R/qauNinFicS0Sfl+6jb2Huqif3stslhSFIZHIFiWXb6D78BmnrPkJ84UZGW0/QtudvlN387bBZaXRxqSy78m5adzxNnesVylX+YnLRpkuILjJSX/93hhsPI5erAZG4vDUM1h+g59jbJBRdtORcYOYy6Z/ung/12XEucI2PAoTd2isRPHKVlriyy/ja17/B1VdfjVwuJYyRCI7nnnsuqHZqtZovfOELcx5nwZWGJ598kr6hUZZtvmihh7qgWChz91iemuxtn0JtjKen4p8IMv8lIiCQ6cqhifAFyXldDtzjo+ji0wPbTp9XuF6oSp0RVVQM1r4WorOXY+6sxZRZikpnQvR5EX1exjpOYhTi59T/QrodLPaVyYWmv3oPvSfeJXXNNURnL1+QMc7X7zih8CJ6j+2g9/i7ZG66kYxNN1H9/E8YqN1HUtm2sI1jTCskJXkNvX0VHDvxR+QyJaUlnyIleTWJCWUMDdUxZu4gLXUDxyv/CIDo9TBQu4/kFZfP27oYqcl3MFaIYC0Zi+W6c4+PAgJKnTHSolzQJBRvovHv+3nyySe5887giz9KSJwLFlRpcLlcfOs/7id++XZk8ogUn15SLAa/19MyCKuvYqTlBPZOf30DY5cPhzAKQKw8PMWibENd/r6H5JiEiec+l+9iupevMa2Qgdp9aGNTcVlHiMkqx5C6DPBbILzDw2RlXhzyeOeS83VyOx2iKNJXtZu4ZWsXzJVlMdxvC4VcpSF55eV0HXqN+KKN6GJTic9fR9fh14hKzglY4OaLqcmJRWPC5bIwPOyPVxgcqiExoQy5XEViYjmJieUMDFRjt/vrRRhSC+g+8ga9J3aii0v1Z45JzkUXlzbJGi2KPkSfD6/ThsM8iNM8SKItHoMhLSzyh4v5XEuh3ttTKRtTjR/qc8JpGUapMwSSUHhddmRKjZQ0Igjm+vtP9RvJ5Ariyi/nW/9xP5/5zGekgm8Ss2K32xkeHiYtbeJzsbq6mtLS0rCOtaAz+T/96U943ZAtFiMsQnNsMITbr3WpoKsdwTk2QEb6Jto73sdmHyQutoCE+FJaWt9Bvz8R9cY1c+5f9Hlp3/sCRmMmcXGFYZF5updvYulWBuv20/b+s5iMWUSl5GEf7qFn50vExxexft19YRn/XBGOCUIw/UfqfjU1OfF4nLhtYxjTChak/wuBpNJtDNbup+vQ6yy78m4yLvo41v5Wuo+8Rf5V98y7/9PfY3LSauobzgRaj4/3Q8KZdi6Xhdr6l9Fq47Dbh8iOWYeteBOOsQHGB9rordyF7+ibgIAmOhGl1oDTMoRrfAxE36RxO2QKVpTfQUxMbkCOhXQTPJeEIlu4LRnW3hb0idmA31Wp+vmfoE/IJG/7nUG5Fi9Gl6uZmO9EPxzX0XR9GCmhwrWTP/3pT9xzz/zvVYnzl+eff56vfe1rxMbGIooijz32GBs2bADgs5/9LEePHg3reAumNLhcLr773f8mJXlbwBd5uhtkMT5oZvJFvRAYGq5HFL3ExubT3vE+TucYgiCjrPRWjlT8nt7eCpY3lc35txtuPoZjtI/yNV9cEF/1s38no2jClfcRxm39ZGdegrbZRWXVG4yMNiGTyclI3xz28c810ym3M12vs61Uznath+qSEcoqav9AJQC6uNBWlBfD/TmfhYZg5J/ut/3wOIJMTtq6a2ne8Sccew+g3byRmJzl9FXuxut2hq3GhUKhpqT4Fk7WPOv/Wz6x39r6VwBYs+pfOFrxGJ1d+1kZfxdCVAFEbcaX5cVmG8Bs6cRi6cTtthETU44mJRqZTIEgyFAotOi0sahURqpOPsOJqidZu/oL6PWJge9iLi5CFwLB3Mdet5PxgbZArFp/9fsgyLD0NDLaVklc/rqg+l9oxSFc/YfLMrRQCIKM1JStfPe//pvbb78dlUpKLy0xNQ8++CBHjx4lISGBw4cP87nPfY7777+f2267bUGycC2Y0vCnP/0Jh91HUsHUkdxnM5+V03CaBaUXzBlionPR65M4dvwJABoaXyc1ZR2CICM5aSX1Da8yNtaGqSkrcEwoD/Oho+8RG5N/TtwMBEFGZsaWCdvEUyuYLrdtwcePBMFcy/O93kM9fiY3jLP3jY61UVf3Msb0YjTRSWGXY6Z+5jMhmWvQ63y+x9m2myhkPGkldQ1/Z70pm2zVSnq979L26hOkbf8EamN4YnmSk1agUuqorvkbDU2vMTrW6t8hCAwOnqSs5FZUKgP5+ddx/MQfaWl9h5zsyxEEGTKZnKioZKKikiFl7axjLy/7DAcO/ZLm1ncoL71txvOXmB1jo4MTPa+AIGDKKAZgrKOG2LxVjDRXIO8cwSQL33c7V1fLULNQLfXrISlxOd097/Pkk09K1gaJaXG73SQk+E27a9eu5b333uPGG2+ksbFxQVwLFyQdicfj4YEHvk9KytY5ryLPVnRnvrm0P9z/Un/AhBu12sjK5XegVOgA8PncjI62ApCWuh6jMYOTtc9PqDY71Xf64e/V1OREXd2HxdJJSvLqc3Y+HyY5aSUASqUuYjJcyEx3jXi9bk7W/A1dYiZ52+8Iqp9wyzWXY0JN3zmf8UKlYNlH8XqddHTuQauNoazkVuz2Iepf+V/U1X1hU1piY/MpzL8eALt9GK/XhdNpJiN9MwkJZQDExRaQk72d1radfLD3x5yseQ6zpSuk8eVyFdlZlzIwUMXho7/DbO4M6XiJiXR1H2Co/iCZm25CbYxHFEWclmE0pgSUMi0eT2iFAWe6nqZ6H4TjPX4+ErA2fPd7eDyeSIsjsUhJTEzkxIkTgb/j4uJ4++23qampmbA9XCyIpeGFF17AanVQsGx2K0OwnK8PhsWMWm1iWd411NS9AMDgUC0xMbkIgozC/I9y6Mj/MjRcT0J88Yz9fPi3GxqqQ0BGbGz4/dXBb0Ww2Ydw2IcZtw2gUkWhkGuIjc1HdirILzo6h6TE5aSnSVm9FhP9/SdwOEYovfbzMyZPWOisVaG4Ui12FAo1MkGBWm0CICGhBKMxg0OHf01L67sUF90YtnNJTCzn0oSyGVe4crIvIzZmGYNDNfQPVDFQUc26tV9Cpwve6pGSvAalQkdr27scPfZ7igpvIjlpRThO4YKif6CK+oZXSU+7iDz5CmhyMpItB9GHesCFXKbE63WH3O+H76GFziS1VO7FUElMKKez611efPFFPvnJT0ZaHIlFyFNPPYVCMfFdqVKpeOaZZ/jSl74U9vHCrjSIosj3v/8DkhI3BCZoEksXk+mM+9HZrkQGQxqGqFR6eg7PqjR8mIGhGkymLJTK2YuphYIo+ujuOUxr2y6cztFJ+7XaOEqKbsZozESliqK05FNhHV9i/vQPVhNtyp7SLelcTgymm8QsxcmJIMiIikqhq/sAiQllaLWxqNUGMjO20tTyFjnZl6HRRIdxvNlN4iZTJiZTJlmZl3D4yP9SVf0XVq28J2jLnyAIJCSUEBubT139y5yseRanc4yszPClk50rw8MNdPUcQqM2kZW5DZXKEGmRAvh8XkTRR3//CazjfXR27SMxsZz8ZWcKicW0ekmML6Or+yAKuQqfL3SlAeZutYPFGecYCWQyOUmJG/n+93/AJz7xCSmTlcQk0tPTp923eXP44zXD7p60c+dOmptbSQ3CN3UhGRltweEYjagM5wMaTQwpyWsoLfnUpJW8pKQVDI804PMFbzp1uawMDzeQlBTevPuiKNLQ+Bp19S9PUhgu2vhN1q35EgqFlmMn/sB7e/6bE5VPhSS3xMLj8/kDYtWa6EXjPhhOGaZyzThXlJZ+CkGQcaTi94yP9wGQmroemUxJT294s2uEgkKhpqz0NsZtAxw++kjIgXtyuZLioptJT9tEU/NbWK29CyRp8IyZ2xkYqKKr+yBV1X+d4MIZSXw+D4cO/5oP9v2YmroXGRioIjl5FSVFn5jkRpyVuQ2ncwzreG9Eiioulvt/MZCSvIampmZ27doVaVGWPBaLha9+9atkZWWh1WrZtGkThw4dmvGYp59+mhUrVqDT6UhJSeHOO+9kaGgosP+xxx5j69atxMTEEBMTw/bt2zl48OCEPh544AEEQZjwSU5OntBGFEUeeOABUlNT0Wq1XHLJJVRXV8/pPD/3uc/xhz/8IfB3W1sbb7zxBmNjY3Pq72zC/jT4wQ9+THLSOuTyyEX7W6w9VBx7jKrqZ+ZkWpU4g0wmp7joJpISJ0/yddoEfD4PLtd40P25XFZAJEqfEkYpoar6L3R27QNAJigoKf5kYIx9+39KXf3L5GZfTkryalJT1jE0XE9d/SsLkl1AYm40NL6G3T5MRtqmSIty3qHVxLB65b+gVGo5UfkUouhDoVCjUZtwu4O/fxeCqKhkUlPW4XHb53S8IAi4TyU0UCjCa70MBlH0MT7eT8Wx/+PdXf+PltYdACwv+yyjY61098w8KQmfHCJWay9er2vSvvHxfk5U/ZlxWz8ej52S4k+w6aJvUVx445QeAVFRKQHrg1SJPbIoFGqSk9bxgx/8ONKiLHnuuece3n77bZ566ikqKyu58sor2b59O11dU8dV7dmzh9tvv527776b6upqnnvuOQ4dOjQhMH3Xrl3ceuut7Ny5k3379pGZmcmVV145qc/S0lJ6enoCn8rKygn7H3roIR5++GF+85vfcOjQIZKTk7niiiuwWCwhn+dbb71FUVERACMjI6xevZobb7yRkpIS6urqQu7vbMLqntTQ0MB77+1iw7pvhLPbkNHrEkhJXoMhKpXd738XrTaetav/VQp6DTMKhT93t9NlRqMxBXWMShUFgN0xjMmUGfRYFks3g0O1jI21IQgy3B4bhqhUCgtuAAi8+FJS1rIs92qUSh1JictPpXLsorn5nxyv/BOCICczYytFhR+npvZ51GojWZkXR1TJlTjtWnaI7KxLMRqnN7cudSK5YqpWGygq+DhHKn7H6GgrMTG5aLVxjNsGIiYT+H/7kdEm4uKK5ux+kZ5+EYNDJzlZ+xwrym8P3M+jY20o5Gp/ZqYFoKn5Lbq6D0wKFk5L3UBs7DKSEpfT3rGHtNQNC+pa4vW6qar+C0PDdcgEBdExuWSmb2FktJmxsTbGxtpQa0yUl36G+PjioGRJS/Xneo+NzV8wuSWCIzVlI7t3/4zGxkaWLVsWaXGWJHa7nRdeeIFXXnmFbdv8bowPPPAAL7/8Mo888ggPPvjgpGP2799PdnY2993nr+WUk5PDvffey0MPPRRo8/TTT0845rHHHuP5559nx44d3H777YHtCoViknXhNKIo8otf/IL777+fG2+8EfBnIE1KSuIvf/kL9957b0jnOjY2FnBb+tvf/kZqaipHjhzh/vvv59vf/jYvvvhiSP2dTViXEP73f39LUlIZanVkfThlMgXFRTeRmOjP2GG3D/L+Bw/S1v6+tLIcRk6vQIWyQqhSRaFQaHE4gjOTOZ1mjlY8xqEjv6G9430EmRxBkGG19NDVfTCQOrW05FNcdskPKS68MaAcCoIMvT6JlOTVXLTxG5SV3IooehkZaSQ5aRVZmRfT2raT9z94kP6BqhDPXiKcuFxWRNG7YJM7CT9GYwZyuZqBQb/ZOzo6m7HRVuz24YjJ1NH5ATbbICkpc8+mZjJmsKL8c5jNHVRWPc3oaAtV1c9wtOJRTtY+H0ZpJzJmbkep1LNyxd2sW/MlLt76AJdd8sPAYkZy0iocjuEzKWgXAFH0cbLmb4yMNlFS/Elyc69keLie5pa3aWvfhSCTk59/HRvXf42EhJKglRdBEEhP24hOG7dgsksEh1ptICmpjP/9399GWpRFh9lsnvBxOqdemPF4PHi9XjSaiYUKtVote/bsmfKYTZs20dnZyeuvv44oivT19fH8889z7bXXTtkewGaz4Xa7iY2NnbC9oaGB1NRUcnJy+NSnPkVzc3NgX0tLC729vVx55ZWBbWq1mosvvpi9e/fO+h18mIyMDFpaWgB48cUX+dznPodKpeLzn/88H3zwQcj9nU3YlAabzcbjjz9OUkJkYxnORqUycPHWB1iWdx0ATc1vSH7sYcLrdVNZ9Wc0mliio3NCOlb0eYN+cTU1/5NxWz/lpZ9m6+b/ZEX57RQVfhy1xkS0KTto07lMpiAxsZzSkk9htnTy3p7/Rq9PZs2qe/H5PFRV/4Wq6r/g8Vy4vrORRKHQolabaG3duWh8wM9HBEGAsxZOUlPWo1TpqW94NSILKiMjTTQ2vUFiQikx0bnz6is6Oofl5bczZm7j6LHHGB1tISV5LVZr94IpRVH6FLxeF7ExeRgMqZMslrGxyzAaMqirfznk1KXB0tq2k8HBGspKbiU5aSWZGVuIjyvGbOkAICfrUtLTNiKTLVhZJolzQGLCWh5//HFstqVfW6i9O57WzoR5fdq7/dnWMjIyMJlMgc+PfvSjKcc0GAxcdNFFfP/736e7uxuv18uf//xnDhw4QE9Pz5THbNq0iaeffppbbrkFlUpFcnIy0dHR/PrXv5723P7jP/6DtLQ0tm/fHti2YcMGnnzySd566y0ee+wxent72bRpUyA2orfXH4uVlDQxAUhSUlJgXyjccccdfOlLX+Lb3/427777Lh/72McA8Hq9WK3WkPs7m7ApDc888wxqlQmTKTtcXYYFuVxFZsamU4XJ5AG/d4n5MTLajMttJTVlLXK5MqRjDcZ0BgaqZ52kiKKPwaEa0lLXk5BQGnBBOlnzHB6Pg+Kim0KWOzGhnKKCj+P1ujhZ8ywej5OiwhuJiy2kf6CKhqbXsNmGsNuHEUUfVmsPdvtIyONIhIZcrqS89Das4z0MDddHWpzzGq02NrB4olCoKcj/KEPDdVQc+78Fm9hOhc/n5UTVnwFInCJmai7ExuSxZdP9rF39BTZu+Popv3yB4ZGmsPR/NuPj/fT3nyAutnDaNoIgo7joJlxOC1XVz4RdIfZ63XR27SMtbQPxZ2WxyzirmGV/f5VkYT8PiDZlo1Qa+etf/xppURYVHR0djI2NBT7f/va3p2371FNPIYoiaWlpqNVqfvWrX3Hbbbchl0+d6fPkyZPcd999/Nd//RdHjhzhzTffpKWlhX/913+dsv1DDz3EM888w4svvjjBonHNNddw0003UV5ezvbt23nttdcAvwvS2Xx4MVUUxTm5NX7729/mlltuYe/evfz4xz8OuLQdOnSIzMzg3cKnImxLD7/+9W+Ji1u9aFOCmc2diKKXpuY3iYstkNwg5on2VIrGqKjQA5oz0jdRWfVnbLYB9PrEaduNjLbg8diJjTnjU2s2dzI80kBZya1otbHTHjsdgiDQN3AmAOl45R8D/5fL1aiUevYf/J9Jx8XELCM3+/IJKWglwotenwwIuFyhB35JBE9UVMqEgmoJ8SWsKL+DqpPPUH3yWZaXf/acBL+6XBa8Xicryj9HXNz0E+9QkcuVE+JitNpYbGGM2xBFke7ugzQ0vY5GbSIn+7IZ2+v1iZSV3cbxE3/kg70/oqTkk8SFqUbNwOBJ3G7bhHozDscolaeUMZ02ns7ufURHZ5OYGL66SRLnHkEQiI9bza9+9b/cddddkRZn0WA0GjEajUG1zcvLY/fu3YyPj2M2m0lJSeGWW24hJ2dqb4kf/ehHbN68mW9+85sALF++HL1ez9atW3nwwQdJSTkz//nZz37GD3/4Q9555x2WL595EUSv11NeXk5DQwNAINaht7d3Qp/9/f2TrA/BIAgC999/P/fff/+E7X19fdx2220h93c2YXkz1NTUUF1dSVLi4iyu4/W6sY6fMT85pUnJvNHpEtHrEmlr3x2yS8/prEYzuQy4XBaamt4gKiplwkTdcmqyY5hHsGxRwcdYteIeNqz/Ghdt+CarV93LxvVfZ9uW7xCl99+8AjLyl13H2tVfoKjwRkZGGjlS8SiNTW/MeVyJmfF6nYAoJSxYYGJilmG1duN0mgPb4uIKKCu5laHhOvr7K2c4Onyc/p3tjoWNp9DpErDZ+sPW3+DgSeoaXiE5aRXr1n4pqBoXsTHLWLv633B7bAwNzS97ydk4nWMoFNoJhfE6Oj9AQGDTxm9RWnorQCD2S2Jpk5S4gurqSmprayMtypJGr9eTkpLCyMgIb731FjfccMOU7Ww2GzLZxGnyaavE2da7n/70p3z/+9/nzTffZO3a2V30nU4nNTU1AQUhJyeH5ORk3n777UAbl8vF7t272bQp+GyC/+///b9J6V7P5pvf/Cbf+c53gu5vKsKiNDzxxB9ISioLe7GucOFwnHEvKS76BLExUvaB+SIIAvn512G19nC04tFAysNg0GhMaDSxtLS9y8DAycALzefzMDLSTH3Dqxw49EucTjOFBR8LWK+8XjcDgyeJ0qegUUfPSW5/qkktMTG56HUJaLUxRJuyEAQ5ff0niIpKITYmH5XaSFrqeozGdFJT1rJ6lT97QXvH+wEXN5fLitttk0z/YUKp1GMyZuKeY+pNieCIjytEQMbg0MSJR1xcAUZDOn39JxZcBpfLSsXxxwEWvJ6O0ZDG6FhbWOLZfD4vTc1vEROzjKLCj4WUdc1gSCMr82K6ug9MsPTMB4VCi8fjwGzuCGxTKvX4RC8qlSGwyJKQUDqvcXw+Lx6PA1EUpeddBFEqtSQllvLEE3+YvbHEJN56662Ai9Hbb7/NpZdeSmFhIXfeeSfgd+s5O+PRRz/6UV588UUeeeQRmpub+eCDD7jvvvtYv349qampgN8l6T//8z954oknyM7Opre3l97e3gmxA9/4xjfYvXs3LS0tHDhwgJtvvhmz2cznPvc5wD+f+upXv8oPf/hDXnrpJaqqqrjjjjvQ6XQhWQZ6enq47rrrSElJ4V/+5V947bXXpg0Mnyvzdk/yer388Y9/Ii3lI+GQZ0HQ6xNZteIe9PrEQMpPifkTG7OMNavu5fDRR2hseoP8ZdehUMxeyVMQZBQV3EBj81tUVv+Z+PgSnI4xrNYeRHyoVEaSEleSlXlxIBOXz+dh9/vfBaCo8MY5ucE5HGMcOvIb3O5xVCq/UpCddQmCIKOh8R8MDtVMaG+zDdDTW0F3zyG02liUSj1u9zjtHR/g8bhobnkLQZCjUkWxeuXn5+QuJXEGQRCwO0Yk96QFRqnUERO7jKbmt9Bp44iJyQvsS05eTUPDPxgf75/RdXC+tLa9i802SHHRzcTHhVZRPlQS4ktpad3B8EgT8fN0g+rpPYLNPkhpyS1zOj4n+3KGhhs4WfM31q354rxTPackr6Kn5zBV1c+wbu2XUSq1xMbk0dzyTw4e+iU2+yB6fdK8g6APH/kt1vEeBEGGUhnFujVfjHiWxAuV+LgV/OEPf+RHP/rhtL74ElNzOuahs7OT2NhYbrrpJn7wgx+gVPrjMnt6emhvbw+0v+OOO7BYLPzmN7/h3//934mOjuayyy7jJz/5SaDNb3/7W1wuFzfffPOEsb773e/ywAMPANDZ2cmtt97K4OAgCQkJbNy4kf3795OVdcaL4lvf+hZ2u50vfOELjIyMsGHDBv75z39iMAR/n/3hD39AFEX27NnDq6++yr//+7/T1dXFFVdcwfXXX891111HfHz87B3NgCDOc9ngrbfe4uabb2Xdmn+XisBcoLR37KGp+U1kMiXRpiyKCm8K+oXS1PwW7e3vExdfRFxsPkZDOlFRKZOuJa/XHVAaLt76wJxetsMjTRw7tbp5mtKST5GUuJzGpjdo73h/2mP1+iSMhgyGhmopKf4E7Z0fMHxWwO7ystuJjy8KWSaJiezd9xBJSSvIy70q0qKc17jddqpOPsPoaDOrV34+4ALo83nYf/AX6LRxLC+/fcrCX/PB43HS3PJPenqPYDJls3L5HWHtfypEUWT/wYeJNmXPKXnCaez2YQ4e/jUJ8aWUFN88+wHTMD7ex4FDv6Ko4GOkpq6bcz8BuRwjHDr8G6JN2ZSXfQZBEOjtO0ZH5wfodYnk5V497wn+e+9/j5jYZVgsXTgcI6xacfcEZVPi3OHzeTl89GGef/4ZrrpqaT0nzWYzJpOJjEceQKbVzH7ADPjsDjr+7QHGxsaCjmm4EKmpqeHVV1/llVde4fDhw2zYsIHrr7+eW2+9lbS0tJD7m7el4c9//gtxsaWSwnABk5mxhaTEcrp7jtDatpPqk8+wvPz2QPG3mcjNuYKszItnbTsy6s9+snYeq3OxMXkU5t9AU/NbCDIZen1yIAd5Xu7VREfnYojyB+NaLF2cqHqK2NgCVpTfPun6FhGJOrWCp9PGz6owqGo6cBVnzEnuCwVR9OF0WVCrpBfAQqNUallR/jkOHfkNre27WVHuN8nLZAoK8q+jsvLP1NQ+P+cV9elobXuXnp4jmKKzycm+PKx9T4cgCCTEl9LbexRR9M3pXeVyWTlR+SRKpY6C/OvmJY9en0RMdA59/cfDojRoNTGUFN3Miaqn6OzaS0b6ZpKTVpKctHLefZ/GJ3oxGTMpLf4ke/b+kL7+E5LSECFkMjlxsaU8/fTSUxokzj3FxcUUFxfzrW99i/7+fl599VX+/ve/A363qVCZl9Lgcrl4+eWXKMz/9Hy6kTgPUJ/KIhITncuJqic5ePg3FOZfT1zczFlCBEEWlHJhsXQhkymJ0oeeSeBs0tI2kJKyBhAmrKL6M1OccV1Qq42sWf1v6HUJU04y4mILQs6AIikOM2OzDSKKXjSSm9c5QSaTk5G+hdq6Fxi3DaDXJQAQH1dEXt7VNDa9QV7uVUEF+s7GyEgT7R17/OPok86JheFs4uMKae94D4ulO+SK46Lo49iJP+J221i16vNBPa9mlSe+mIbG1yekVBwarkejNqGfwzMuPr6YjPTNNDa9icmUjdEQ+griTGg00TgcI8hkCnKzr6C+8R9kZGwJXDMS55b4uDJefukvuN2PBVxrJCROMzo6yuOPP05vby85OTmsWrWK5cuXk5iYyN13383dd989577nZR7YsWMHMpkKo1GaCEn4iY7OZu3qL6DRRHO88o/sP/hzOjo+mHf2jtHRVozGjLAUKJLJFEG5XZiMGWGZIJyNqqZj9kYXKD29R/3+9vMs8iURPMlJK1Aq9fT0HJ6wPSV5NWq1kYbG1+bVf0/vUfYf/DkVxx/H6TSj1cRQVPixefU5F4zGTBQK7aTg72Bwu21Yrd1kZ18atkmygAxBkAUUhuHhBo6f+BOVVU/POWA7L/cqovTJVFb9OWyB1qeRy1VYx/1FpmJjCwAR11nZtyTOLUZjOghKduzYEWlRJBYhN954Iz/+8Y+prq7mt7/9LRdffDEmk4mCggJuuWV+1uN5zcCefvoZYqKLF21tBonIoNPFs2rF3QwN1QYKpg0O15GZvoXY2GV4vW46OvcwOFSLWmUkM3Mb0TPUP2hqfouR0Sbyl12L71Q16cGhWuQyFV6fm9iYvHkHFEpEHp/PjUoZFXKxQIm5I5MpiI8rYnColmV51wS2K5U60lI30Na+a84uPaOjLdTUPk9CfAk52ZeTmFAWMTdWv0tHPkNDdeTmbJ/9APxpnweH6mhpeQe5XB0oXDqV4h+qBdHrcyMT5Hg8TiyWLiqrnsZoTMds7qSr+yAZ6cGnWTyNTKagvOwzVFb/maNHH6Wo6CaSk+afBt1uHw5kYQIYM7edGk+6TyOFIMiIjSnm6aef4eqrr460OBKLjAMHDrB79+5A+len00l1dTXHjx/n+PHj8+p7zkqDx+Ph73//OwXLbp2XABLnJ4IgEB9fTHx8MUmJy2loep3jlX9Er0/C7RrH43GQkFh2KmXr7ykp/sSUPrg+n5e29t0ANDS+RkPja6fSDJ5Jy6lQaFm5/I4lY/GS3JSmxmLpRq02RVqMC474+BJ6eo/Q03PklOueH5UqCq/XhdfrmpPFbWi4HpmgoKz0tkUR85YQX0rVyWcYGWma0R//tFKwz/w37PYhEuJLWZb3EUytVmBqS2Go97ROG4/P5+GDfT9C9PkwRWezvOwznKx9nt6+Y3NSGsCfznr1yn+htu4lamqfR69LwGBInVNfpzn7nhRFkYbGfxAXVxSym5dEeImPK+WVV/6K1+uVsihJTKCsrGxCfQm1Ws3q1atZvXr1vPues9Kwf/9+vF5RenBIzEpcXCGxsQWMjrXQ2bkXdXQeGemb0WpjEEUfJ2ufp67uZUymLLSamAnHymRyyks/TXPrO8TG5KNWGfD6XMTELEN+aqWrruEVKo4/QXnZp5d8DY4LVaEYGW1mzNxGWen8qlVKhE58XBGpKeuorXsJn+glLXU94I8xAUJ2l7Hbh+no/IDOrn2YTNnnXGGYzgUwTTTRYcyirv4V1q+7b0pXx9PHenxu7PYhCqI2kitfBa3WSW3nQ0JCCRs3/jt9fSfw+dxkZV6MXK4kLraQ2roXcbttcy5yKJcrKS66EZttgKqTz7Bh3Vfm5dYpk8n9LoMBdyoBkzFrUSiCFzJGYwZer8j+/fvZvHlzpMWRWET85Cc/4Tvf+Q4vvPACGk14Xazn/CR59dV/EBtTID04JIJCEARionMn+asLgozC/OsZHWmmvv7vLC+/fZK7W0JC6YzFiVaU38EH+35MW/t7i1JpmGoSM5VycLrd6X8vFOXBau2hsuppDIZ0EuLnV4RKInQEQaCw4AYEQUZd/cs4HKPk5lwRCMgdt/UHVd9mzNxBS8s7DI80oFBoyM25isyMczuZmSlmSBAEypQb2Gt+nq4j/yAvas20bb24AdDKg09VGqrCr9XEkJ118YRt/ueXyPBII0mJy4Pu68PIZAqKim7k0OFfMzTcQEL8/GphJCevoab2OY6f+BMgYrWGN2ZCInT8Lkr5vPrqPySlQWICOTk5WCwWiouLufXWW9mwYQOrVq0iMzNz3n3Pecb/0kuvYDItvgmaxNJDodBQUHA9Q8N1dPdMXwJ9OpRKLUqlDkPU/MzwoaKq6ZhXYPPp46frZ779LxXMli48Hjsx0blSfFSEEAQZBfnXsyz3Gtrad9HQ+BrJSStQq00MDtbMerzF2sOx40/gclkpLryJzRf9B9lZF4clccFUfPjeCfZeMSjiyNKV02w7ild0T9tOLdOhkxsZcfWEU+wA08mq0ZiIikqlf6Bq3mMYolKI0ifT13ds3n0lJ60gKiqVoeE6PB4HA4M1tLTuoLnlHVyu8FphJILHZMrnpZdeibQYEouMm266iY6ODi699FIOHjzI3XffTU5ODnFxcVx22WXz6ntOT/T29naamhrYfNEn5jW4ROh8eBU61EnlYl29TogvIS11I3X1f2dg4CRRUclkZmxBpZp9pc9i6cLpNKPTza/S4dl8+Hud6fteaJeis8dcrL/ffEhJXsPISBP9/SfIy71Ssl5GCEEQyMzcikyuoL7hVYyGNGJjlgWCpKf6XQYHa2lpfQe7YxidNo5VKz8fVFX4uRIOJTpDW0Kr7TiDzk6SNDnTtotRpjLs7g6p72CeBWdbFKdqm5y4guaWt/F4HPPO3paUtJKW1nfm3ZcgyFi98h7e2/M9DIZ0xsf7aGndgYCMjs49ZGZsJTvrMknpP8fExuTzwb7n6OjoICPj/Hs3SMyNkydPsn//fpYvP2OtbG9vp6KigmPHjs2r7zkpDTt27CAhIQelUjuvwSVmZ7qX5Fxfnot5AlqQfx06bRzDI410dx+ip/coyUkryc66jNGxFjSaGAxRKZOO6+07hlplICV57kE+s32fwew/F9/n+ei6JAgCqSnr6Os/jsXStWQC2s9X0lI3YjZ3UVv/MrnZ2+npPYLNNohenwj40x9brF1YrX309lWg08aRmrKezIzNYVEYFtq6pldEo5ObGHJ1zKw0qJLpctTiFd3IheAzBU13j0634PDhtklJy2lsfpPOrgNkZ12MzTZIT+9RBEEgPW0TKpU+aFmSEpfT1PwWA4Mn5/V8BL/Lk16fhMXSGdgm4sPrddHSuoOY6Dyio7PnNYZEaCiVWhLis9mxYwd33HFHpMWRWCSsW7cOq3WiBTAzM5PMzExuuOGGefU9J6Xhn/98B61GCoAOF5FyQVlsCoQgyMjI2ExGxmbs9mFaWnfQ3XOYjs4PAm0u3vrfk1Jymi1dmKLnFnC50N/9QvUfivIwFxnO9fVwOthWoZhb8KdE+PDHOFyPxdpFY/MbgL8isl6fSHf3IWrrX0IQZGg0MWRnXUpW5rawuSGdq2dhrCqNIdfMfvlRCn+RQatnFJMy9PoMoZzL2YsOarWJzIwttLS+Q1JiOQPV79E27q+jodMlhpRGVaOJJkqfhNncPm+lQRRFxsf7Tv0lkJy0it6+o4H9ff3HFqXScL4nl9Bo0vnnP9+RlAaJAF/96ld54IEHePbZZ4mJiZn9gBAI+UkviiLvvvsuqckfCasgFyKLyV99NlnO9UNXq42lpPgTOJ0WBodqqKt/+dQecVJbr9eJcg6m9+nO2eYZA0CnWPzpP2d7Ic7XInWufneZ3P8oEkXvORlPYmbkchVlJbdytOL3uD02qk4+g1ymxOEcJS62kOXln13SbmRxqjQ67ScZ94yhn+Y+F8XTz5rJz5yF4Ox7NdebQ7v4Ps66alSyM8+21vo3GR6uJy/3KtRqY3D9qo1hiTuQy5WUlnyK6pN/BURMpixSklfT3vE+Q8N1dHUfRCZTEhWVMm8FJVzM5gp2PhATncuOHW9OqC4ucWFz0003AZCfn8/111/Pxo0bA1Wh1er5WYNDVhqampoYHByguHD+UdiR4lxOiBaTYjAfIuEWo6rpQAWoi9djiErlSMWjNDS+RlHhxwNtHI5RrNYeEhNCyzQy3e8y7Orm0MjfERG5KPYmBp0dxKnTiVYmzUn+c8FCjhPu3326/nTaOMCfgee0G0ykOR9dwUJBr09k7dovsm//T3G7x4lNWkl6VCppqetnVBh8Pi+DQzXotHFETeFOOB3n8lmZqM5CLdPRPH6EctPUgYGj7j5kyDEo4s6ZXKc5rSi4fA6Ugv8ln6ktxyO6GByoZbC/mrLlnyE2Ztms16lKGYXdPhgWuRITyulPqGJgoAq5XEVMTC7R0dm0tu2irX13wCqs0URHtLL7+fLeDQaTKYsTVf00NzeTlzd9/RGJC4eWlhaOHTvG8ePHOXbsGD/5yU9obW1FLpdTVFTEiRMn5tx3yErDe++9R5IhfclW4D37YRLO1YcL5SF1riZSH/6djMUZLMu7hobG18hI3xyYWI6MNgOQlrou5H4/TFvMAPW1ryGeWllstB5iwNVOv6uNDTE3IBNmLqBzPl8D053bXK0ckwL6VUaMhgyaW/5JYkLpvANA58tCPSdmGysUd7NgA27PJtTz6Onxu8UUFd5IasraGfsVRZGR0SY6OvcyNFSLXp/EhnVfWZT3hVxQkqNfRa1lLzn6lQFXpLNxeC1o5YZZ7/uFQCbI0cgMDDo7MCkTkAsKig2bEQQBt8/JkdHXaax+mc2xn4RTq8vTxlKoohgbaw2LXIIgUFL0Cdr1SbS0vsPJmmcpLvIX5jSZsjh2/HEAf4zDyoVXGhbjtXWu0db3kWhI57333pOUBgkAsrKyyMrKmhC/YLFYOHbs2LwUBpiD0rB3z140npg5Tx4XatIZzAsy2EC0YF+2F/IDa7rsQuHCVZwxaeKWLSbTpYun6uRfWbv635DLlQwPN6BU6idMMoP9Xbyihz5HMz2ORkSDlsHaWlQqA16XBQDPqZSMY+4+dvQ/weroa4hTS7E8ZzPfe+D0hFwQBPKXXcvRit+zd/9PWbXynimD3s8FkbQQzaRAfLj9TM/SmRIoBHuvqmo6kFlGAXB0NKEanWxts3ut9DtbaDvwHA7HKD7Rg1ptQquJxefzLOpnZIa2hHZbFdXm3ayP+dgk1w6VTIvb54iQdJCrX8VJy3t4RRd6eXRAPqVMTWHURg6MvEyfs5lkzcSJ4tm/sctlpa//BFrtZKUoVBoaX0enjSMtbQM52ZfT2bUPgJra54iLK0L0nXEtHB1twek0B+1CNRfCfW0FMzdZrNZHrSeGvR/s5c4774y0KBKLFIPBwNatW9m6deu8+glZadiz5wNMijOuSXOdYM/l5pvqhRrMamaofc+3rwuRhYiJ+LDiIBeUZKRvoq7+FfoHKklJXo3H40CjiQnZv9rlc3Bo5FUsHr/Z3uBKo7DgBqJNORw49AsA7F4zBVEb0cujabNXcmzsbbbF34aAgN1nQS+PRibI8YoezO4B3KITj8+Nw2cBUcawu5MUTT5p2sKQz/1C4vQkx2TKZMP6r1F98hkqjv0fq1feE5J7S6hjQnALC/MdY6GPD3WcUNKC5kdtQBDkNI0fJl1bhF4ezbC7B4t7kDHPAL2OJkAkNq6QtLSNGKJSMJmyaWh8jdHe2es7RBK5oKDIsImjo29g8Qxi/FCws1quwyU68IneiFgb5IL/9ez02dDJoyfsi1GlEKfKoNayj1hVKirZxEyGqpoOHF4rxxw78PnclBTdPGFfqM9jn89DR+ceANLSNgCQmbGNpuY3USp0DA3VTjpmdKx1XgXqwkkoFtGpvp/p5i+z9R0OgllIMCkS2fPeB0hcuLS3t4dUvK2rq4u0tLSQxwlJaRgfH6e+sZ5tcTO7gsxlsh6qm4M0iV9ahGuFRt/vAkAh91sWEhJKqa17CbfbhlKpC/q6qLXsxeG1si7meoyKeMRSf6FCl2s80GZb/GeQnVJGNPIo9g0/z46BJwL7TYpE1sRcS41lDz2OhinHUQgqSWkIAZ0ujpUr7ubQkd/Q0bWP4sIbw9r/TC//cPd9PiAIArn6lfQ6Gjk08iog4PSNI0OBXmGiIGoDGdoSfKUTV7uF4TF8+CIjdAjEqzJRybR0OeomKA2iKAayK7l8djTy2StihxuH14pCUOHwWklQZ0/aX2rcxt6h53h34I/EKFNYYboCjdyfjtXiHuLw6D8QkLEu+iMYWyyAJXBsqIqD2XzWxNVlRaWKQqvxZ2Vxe2yBfXK5Cq/XRWxMfqC6eySDkOfiwjfT9pnazqZoTBtzMkO7YJ9XJmUiJ5rewWazodNJGeguRNatW8f111/P5z//edavXz9lm7GxMf72t7/xy1/+knvvvZcvf/nLIY8TktJw9OhRojQGNLLg80QHy/n4wpWYTKgvkA9bG3ynsutUVf+F7OxLSUtdj1ympLnlnxQWfCyoPv0Tgk7StUXEqfyatuvUvpbWdwCQIefsrCma8pWUDijweJzI5UpU3aNUmt/l3YE/AJCuLWaZfh3Hx95mxH2mimwkJhtLHaVSS2xMPv0DlaSlrAu5bsO5fpYsxWfXdPfhVOciF5SsibmW5vGjyAUFKeplmJRJE9x5XB86ZtjVRZxq8bvyyQQZqZoCuh31FEVtClgsR9299DgaiFLE4vY5I3IfRyli8IiuwP8/jE5uZEvcLQy6OmiwHuLAyEtsjLkRpUzDCfO7KAU162KuRy2fehIZyrO4pXUH4M/Uo1T6+4uNLSAzYyu9fcdwnXLpjI7OY2ioBrXGhEy2sNaZD78bpto/Ewvl3jTb/tm8JOYil0YWhU4dxZEjR+btfiKxNKmpqeGHP/whV199NUqlkrVr15KamopGo2FkZISTJ09SXV3N2rVr+elPf8o111wzp3FCUhoqKiowqRKltF4S82I+AZoxqhRUKiMul5mW1h3ExRaSl3cVDQ2vEx9fQgp+M71X9FBt3o3ZPYBJmUR+1PrAKtyIuwenb5wYZfIkmYYHTqIQVHhEF06fDa3cX5FaEGQTTO2q4Q6ilUkMuTpRy/TEq/3yFxu2cGJsB+PeERSCimS1FJgWDB+ewOTmXIHF2s2JyqdYv+4+VKrgJm2SwhA8ociukxspM14SVFvFyTZsXjNZivI5SnZuSVLn0mo7zoi7l1hVKgA9jgbUMj0GeTz7R15kTfS1gX3nirNdjgSmfudq5FGka4sxKOLZN/w8Y55+xj2jWDyDXBR707QKQ2CMINxsRkZbAgknSktuCShWCoWauNhC2jvex2TMxGzporjwRkbHWqmqfprx8T5WLr8TFQtnbZhOcVhsMQdnsxDPDEEQiFYncezYMUlpuECJjY3lZz/7GQ8++CCvv/4677//Pq2trdjtduLj4/n0pz/NVVddRVlZ2bzGCVFpOIbSHQWRTWwicR4S7MqnWqZj7ep/o7HpdazjvRypeJQ1q/6FhIRSjp/4I2L57fjaOmiwHsQlOhBFL1bvCL3uFi6K/hhe0cPhkX+gFNQYpijYFKNKodtRj1ZuRCObeaKqkUeRpi2asM2ojGdL/C3YPGb2j7xIl6OWaFXo6VovdFQqPSvKP8v+gz+npfVdCguuj7RIE1jKysJCcPb3YfNaEPGhlS9cEGw4iVYmoZbp6XM0E6tKxem10W6vBsDqHcIrejky+hpb4j4VWEQ4F9i9/tV7uaA4tRAxRp5+6voHbp8TALN7kMbxw2TrVmBShpa6eCoFwu22UXHsMcAfw6BSGSa0NRakE6VPZszcjsmUhUqlD7h0ms0dnKx5juXln13Q+2UxKwjnEpU7imNHj0VaDIkIo9FouPHGG7nxxvC69p4mpOjRo4crMEyRmk5CIhyoajqCerloNCbKSm9l/dovgyjSP1AVqJJ6ovJJqsw7McRkwqnVuSh9Mmq1iaP2nRwe/gcGZTwXJ3wW7RQuBxnaEuJVGYH0hh+WLVgZdQojJkUiHfaTDLk6g/wGLmwm+e+qDGSkb6G751CgYnQox8+HqfyTQ/n9L2SG3d0AIU9aI4UgCCRrculy1DHuGaXKvMufklW3Cod3nGzdchSCijrL/nMql8M3jkJQsTb6OhLV2TRbjwSUgw9jUiagEFQ0jh8iSZ1DYdTGeY19+hpvbPJXBFepDGRnXTKpnba+j7KyT2OISkUh12A2d3Ki6ilAQK9PYnikEd+prErzuW+ke292ouSxHDlSEWkxJM5zglYafD4f9Q11ESl0I3FhEewLQiZToNaYaO94nxNVT5GedhEF+dezbNl1jI61oFLpWV72WazjvcRE55CcvIro+DzWRF+DQlBO2WeMKoW1Mf6X9HyJU6WhFNQcGnmVpvGjZ1WYlZiOD//uJmMGoujF6TSfMxmklcv5MejswKRInFDJeLGzTL8OlUzLvuEXGHC1kR+1nryoNbhFBwZFHMv06+h1NjLuGT1nMvlELx7RRZejjmLDFnz46HLUTdlWKVNTZNhMsnoZy02Xh6Vat6qmg57eIwCsX/vlQFrrDz+Xnc4xXG4bQ8N1HD76CDJBwaaN36So8EZ8PjcWS1dQY0nMjyhFLPUNtfh8iz8BgcTSJegnS3t7Oy6XE70iegHFkZCYnbNfMKtXfh6l0h+r0Nm1j/qGv9PWthOPx8GqFfdgMmUDoNXGkZ15MeVln5mUnnChyNav4LKEO8nTr6HBeoBm29FzMu75wGnlIbrbb2Gw1Bw552NLhI5P9DLo6gjE+CwVlDI1a6I/gkd0IcMf7M0pJV9AIFVbgFLQBNyWzgWnF+g67TUMuNqIU6Ux6Jz+ukzXFrEy+opAqtZwsCHm41yWcAdRTSOBbWcr1W6fk+PH/4hGbcJoyABE1q+7D40mGkNUKjKZkrGxtkD72VKkB5MpMRz35vl4f0cpYnA6nXR0nH/nJrF4CPrp0tjYiEkfG5F81RISH+b0Q18FbNn0bUBgYLCaquq/4Hb706aOjrWSkrwalTIKq7U3InIKgkB+1Hp8oo9G62GS1XmS4h0CGnkUieps6q37ia9KD2SxkawBi5NhVzce0UWSOifSooSMXhHNtrjbkAvKQACxWqZn3DuKXFCQoM5kyNkB5yisIUGdSaa2FLvXwknz+ySqsxlwteP2OVHK1OdEhhhV8pTbT99/A4M1+AY82OyDgedud89hluVdzfBwAzJBjrWrAaKzgx4zkilalzIyQY5RH0NjYyNZWVmRFkfiPCVoS0NzczM6xdIIbJO4sFDXdqGu7SQxoYzCgo8DkJi4nMSEMgRBICf7cqzjvYyOtgD+F97Zn3NBnn41Ij4axw+dk/HOJ0oNFyNDxt7h508VE5u8Knk+rhwuRfqczWhlBgyK+EiLMid0CtOEjEN6uSngkmRQxOHwjU9zZPiRCXJKjNtYFX01arkOt+hEhoyK0TcZcfXM3sE8GfeMcmz0n3TYTk77nHQ4/BYIAYGszIvRamKx2Qaw2QY5UfUkHq8Dh88661gzxRBNRTDpTWf6nK/oFSaam5sjLYbEeUxIlgaFRwfnZoFDQiJkVDUdpBWvIzVlzQSf3tTUdXR276ezaz/R0ZNXQIOpLj5f5IISvTyaHkcjCaosUrUFCzbW+YZaruOi2JupNO/k2Ng/WS5eHvj+wvmbiaKI0zeOT/QiIiIXlIE0vRKz4xU99DqaSdMWnjdpuQ3KeHocDXhFDy6fHfk0sVALiUyQk64tocF6kFzdSjrtNRwYeZky4yWka4sXbNwW2zF6nU30OpsQ+lNJTJycQjcpcQVmcwfZWZeg1yfRP1CFUqljeKQx0Ga6+LH5EkqtkQsFhUdHY2Pj7A0lznt27NjBjh076O/vnxTn8sQTT0xz1OwErTTUnqxDc67sshISc2S6F0a0NxqLc/TcCnMWgiCwJe5TVJrfpcq8Gx8+tHIDI65u3KKLOFU68ap0yf1vGtRyHWuiP0KleSdV5t2o5fpAYb75IIoig652uux1DLo6AsW0/AisNl1FombpudpEgj5HM27RsaAT2XNNpraMdlsl+4dfwuoZIlUTmeruubqVDDrbGXR1cmnCHVRb3qPavBtRFMnQlYR9PId3nB5HA8v06xhx99DW/h7x8cXIZBOnDCqVntKSWwJ/JyWuoKNjD4aoMzUtTMqFSzl9ISsIU6Ehitqa+kiLIRFh/vu//5vvfe97rF27lpSUlLAu4gStNLS0tKKVT+3fKCGx2HH5HIg2W0RlEASBUuPFOHzjVJl3BrZrZFG02U6gFDRk6ErI0pbPWpTpQsT//W3F6Rvn8Mg/KDdeRqo2f059iaLIgKuNRushzJ5BDIo4snUrMCrjkQtKBARabcc5bt7BZsUn0ClMYT6b848O+0lilalTVi9equgVJrJ1Kxh2d1Nk2EyGtnTBxhpzD+Dy2dDKjejkxgkLCIIgI0mTS71lP4IgUGLYAoictLyHXhEd1sJzPtHH8bG3UQhqMnWlqBPyqD75LPsP/pyszG2kJK+ZpDycJiN9MxZrN/WNrwKQo1tFnn7NhDZSzMLCoZUbaWluibQYEhHmd7/7HX/84x/57Gc/G/a+g1Ya+vr7yJXlhl0ACYlzgdUzTII68sFhckHB+pjrcfnsjLr7MSriUct0WD3DdNpraLNV0mmvYWPMx6WJ6hTIBSVroj9ClXkXleZ3UcrUJKgzJ7Rxem0Mu7sxuwdw+MbRyPQkqLPQyY1o5FG4fA4qx95lwNVGjDKFdTHXE6tMnbQaY1TEs3vwaTrsJyk0XHQuT3PJYXYPMOLuYYVpe6RFCTvn6rc/MvoaLp8dAAEZBkUsy02XE3WqNpJP9CIIAj7R6493MGxj3DNCtXk3W+JuCUuaVa/o4cTYO4y6e1kfcwMqmZakxAL0ukTq6l+mrv4VPB4HWZkXB445WwlQKrUsL/ssAwNV9DXuJUe/IlDsTWLhUct0tPT1RVoMiQjjcrnYtGnTgvQd1N3s8/kYHhlELZP8eyWWJlqFEacvspaGs1HJtCSqs9DI9QiCgEEZR7FxC9vib0MhqDg0+g+c3sUj72JCJsgpM15KvCqDY2P/pMtehyiKuH1OGqwH2T34NMfH3qbH0YTDa6XLXsvBkVfYNfhnDg6/wu7Bpxhx97A6+hrWx9xAnCptSvOtQqYiWZNLt6MBj+gGTrkzOTupMu+Sfp+zqLMeQC+PJkktLSzNlShFLEpBzZroaykybMYrejg29jZe0Z92OEGViVf00GmvBUAmyCgybGbcO0qvc/7Br6IocnzsbQac7ayKvpoYVQrgVwr0+kQsFn/RPoXijBU0kMXuLDchQRBITCynfNO95yy9tYQftUzH8MigVBPoAueee+7hL3/5y4L0HZSlYWhoCK/Xi1ouPQAkli4+0TvtvsXiG6uW6Vgbcx37h1/khHkHa6OvO2+CSsOJTJCxIvoKqk9ZHJrGj+D02RBFkSxdOdn65ahl/smNT/Qy7h1l1NVHj7ORbN1KMrWlQbmAZetW0mWv4/jo28gEGf3ONkT8QWUy5JQYty7oeS4FBp2dDLk6WGW6SorJmQeFURexf/hFBpxtlBi3EqtKYd/QC9RaPqDUeDEGZRxxqnQGXe1k6vxuUiZlIhpZFBb3ICmaZfMav9fZSL+zldXR10wqbmmxdIEgkBBfRmrKWmD+FZ4lF6Xwo5Lr8Hg9DA0NER+/NDOYScwfh8PB73//e9555x2WL1+OUjkxGcHDDz88576DUhp6e3vRqLQRyRwhIREOLO4hohcwIC+c6ORG0rXFNI8fpcVWQa5+daRFWpQoBCUrTFeQqimkx9GIRq6fMh5EJsgxKOIwKOJCDhrVK0yUGLfSNH4UAcjSlROvyqTf2Uy3o55UbcGSua4WAp/opcayh2hlMolLsDbDYsKkTKDYsJmTlveJViaRqi2g2LCFastu9PJosvUrMCji6HU04RN9AbcfpUyNe0IAf+h4RQ/1lgMkqXMmKQxun5OKiqfQ6xPJzd6OurZzXmNJLBwKQYlapaGnp0dSGi5gTpw4wcqVKwGoqqqasG++i5BBKQ3Dw8NoVFJgpsTSRSbIsMqnzhe+WKwMZ5OvX4/FPUS3o0FSGmYhQZ05Ka4hnKRriydlBIpSRDPs7uHwyD+4JP6zKGSqBRt/MdMyfgybd4xNppsli1gYyNCWMuruo9L8rv9vXQk27xi11r0oZGqSNXm02k7QYjtG3qnngsfnmncV6Obxozh846yJunbSvi5HHT7RyxrN5WjanSH16yrOmHOtBYm5oVXqGBkZmb2hxHnLzp07Z280R4J60pjNZlRyqUCDxNJk0NmJ02ejLPdjkRYlaARBIEmTS5V55zmtACsRHBp5FGujr2X34NN0OerJ0pVFWqRzzrhnjKbxI2TrVmBQxkVanPMCQRAoN14KQJV5F7GqNAqiNuIRXVSZd7I6+hqydctpsh4mR7fCb0VTxmHxDIc8lkd002NvYNjdRa+jmVz96ikzXw062zEo4qSaJUsEpUKD2WyOtBgSEWZ0dJTHH3+cmpoaf8a1khLuuusuTKb5JVgJKhB6bGwMhXBhrqRJLG3abFUcH3sbkzKJhPjJ6RIX82rXabeXMXd/hCWRmAqNPAqFoMQjhrb6ej7g9jmpGHsTtVzHsqg1sx8gETSCIKPYsAWZIKfFduxUitVtJKiyqDLvIl6VgQ8vrbYTAOjlMYy5+wKVq4NBFEX2D73ASct72DxmcvWrApYLV3HGhE+atgizZ4B+Z+sCnK1EuFEKKsbGxiIthkQEOXz4MHl5efz85z9neHiYwcFBHn74YfLy8jh69Oi8+g5aaZAjxTNILC08Phd1lr3EqFJYYdo+yX1iMSsMAHp5NEpBzai7N9KiSEyBy+fALTrRyY2RFuWc4hXdHB19A6d3nDXRH5Fi3RYApUxNjm4lbbYTtNuqAjVeRNFHo/UQJkUi9db91Fj2kKNfgUJQ0zQ++2TA4hnG6bPh8I1j9Y6wwnQFF8XdRH7UemSCfMrg5NiVlxCnSqPBenBe52T1jARSykosHHKUktJwgfO1r32N66+/ntbWVl588UVeeuklWlpauO666/jqV786r76Dck+yWCwIvvn5TEpInGt6nc348JKmKUQnN3J2qOBiVxjA76oQo0ql1XYCvSKaFM3cCplJLAynLUBGRUKEJTl3uHwOjoy+jtUzxNrojwZqCEiEn1z9alw+Oyct76MQ1KRq8yk3XU69ZR9u0QkItNkq8fjcpGkLaLWdwCtumza+Ydwzyr6h5wOJAcBfi2Q2BEFAKWhxiOOIohhU7IrLZWV4pImMoRgQBERRZM/QXwFYs+peTKYzNXOWwrN4KSGICiwWS6TFkIgghw8f5rHHHkOhOPMsUCgUfOtb32Lt2rXz6jsoS4PL5UIQpSA3iaWDT/TRfGrlLVqZFFhBU9V0LKmXVIlhK2qZjj6HVOVzsTHq7kUpqNHJL4wifDavmQPDL2H3jLE+5gZiVMmRFum8RhAEigybiVdl0G6vBCBRncWW+E9xcfxnSFRnISDQ5ailw1aNV3TP6KJUY9mDWq4nRbMMmSCj3Hhp0AUkUzTLGPeOUmP5IFADoM/RQrV596QxnU4Le/f+hJM1zzLm8SvWDt+ZJBRHKh7FZh8K/C2lXg0vgijgcs0vm9b5isVi4atf/SpZWVlotVo2bdrEoUOHZjzm6aefZsWKFeh0OlJSUrjzzjsZGjpz/T722GNs3bqVmJgYYmJi2L59OwcPTrTK/ehHP2LdunUYDAYSExP52Mc+Rl1d3YQ2d9xxB4IgTPhs3LhxTudpNBppb2+ftL2jowODwTCnPk8TlNLgdrtBUhoklghm9yAHhl/C5jVTZNiMUFYYaZHmjL84knDqI7GYcPrsgaJv5zsOr5WDw6/gw8eG2BsxKRMjLdIFgSAIpGkLGXX3YfOeCW6VCXJWmK5ALdNjVCSglkeRpM6ZUYEddfeRri2m1Hgx62KuJ01bNKnNdAsqCkGJTm6i3V5Jl6OOUVcvx8fepstex/tDz7Bz4E8c2vMwFR/8ikMHfoEPf02c0657No/fXWblirv9f9sG5/aFSMyOT/DP2SQmcc899/D222/z1FNPUVlZyZVXXsn27dvp6uqasv2ePXu4/fbbufvuu6murua5557j0KFD3HPPPYE2u3bt4tZbb2Xnzp3s27ePzMxMrrzyygl97t69my9+8Yvs37+ft99+G4/Hw5VXXsn4+PiE8a6++mp6enoCn9dff31O53nLLbdw99138+yzz9LR0UFnZyd//etfueeee7j11lvn1OdpgvI5kpQGiaVEp/0k495RNsTcgH7F+sD2pWRhOE3LeAXj3lHKTZdGWhSJD6EUVKhl2vM+1ajL5+DwyD8A2BBzAxp5VIQlurCIV2UCAkPOTnRn1RmRCwqWRa2jyryTEsO2QMG3qfCJXjyiK1DwMBTkJ1uoGHsLpUxDkjoXrczAgZFXEPGRpimk19GE02fD6ZtYIb3EsDVQEfp07ZSWlncAJhUBnCk1q0SIiJLSMBV2u50XXniBV155hW3btgHwwAMP8PLLL/PII4/w4IMPTjpm//79ZGdnc9999wGQk5PDvffey0MPPRRo8/TTT0845rHHHuP5559nx44d3H777QC8+eabE9r84Q9/IDExkSNHjgRkAVCr1SQnz9+C+7Of/QxBELj99tvxePwV5ZVKJf/2b//Gj3/843n1HZSlwePxSEqDxJLB5jUTpYglRpUSaVHmhSiKtNqOE6WIYdwzypCrSwokXESMuvsxnucr7qIoUm3ejdNnY23MdZLCEAGUMjUmZQKdjtpJVe3TNAVkaks5aXmPTnvtrH2Jom/WNqddOE9/hlydeEQXa6I/wqroq/CKrkBV9BF376lK1WmYlEkICMiQk6opJEN7RsHRyaMBGDO3AaBWT04eILkphYkLTGkwm80TPk7n1NnsPB4PXq8XjUYzYbtWq2XPnj1THrNp0yY6Ozt5/fXXEUWRvr4+nn/+ea69dnI9k9PYbDbcbjexsdPHe50OVP9wm127dpGYmEhBQQGf//zn6e+fW+ZElUrFL3/5S0ZGRjh27BgVFRUMDw/z85//HLV6funbpehmifMKm9fMoKuDcuNlS/4lJAgChYZN1Fv3B4o9qWRaVpiuIEoeM6nyscS5wyt6GHP3URC1IdKiLCi9zkb6nM2sNF05ZQ5/iXNDUdQmDo78nTrrfooNmwPbBUFGuraEdns1Kplm2uNlghy1TD8htiBYPD6/f7z2lMI46PK7XeToVpKz/pMBC4HLZ8fmGSNKEYdCNjGj1unq1QCrVt6DXj9R2fZ63TQ2vYHTM0pCQjlpQ8YLtmDihYCqXYVcM7/f1+vwK64ZGRPf89/97nd54IEHJrU3GAxcdNFFfP/736e4uJikpCSeeeYZDhw4QH7+1ElGNm3axNNPP80tt9yCw+HA4/Fw/fXX8+tf/3pauf7jP/6DtLQ0tm/fPuV+URT5+te/zpYtWygrO1Pf55prruETn/gEWVlZtLS08J3vfIfLLruMI0eOzHmir9PpKC8vn9Ox0xGU0qBQKEAQwzqwhMRC0ONoACBBnTVLy6VBuraIdG0RXtFNq+0EDdaDHBr5OwICCeossrTlxKrSznsXmcVGj6MBHz4S1NmRFmXBcHptnDS/T7I6j2RNXqTFuaCJUaVQaLiIWssHpGryJ8SUDDj9q/eyWV7nJmUCA842lunXhfS8cJ+qQ+ITfcgFGHZ1AjDk6iJhvI8u+VFG+moCMRfbE+6Zsp+CqI3UW/fT3v4+MdG5APT0HqW37xjmsXZEfERFpVBT+xx1yFkdfQ3x6qW98BMRBBGl8sJJg9zR0YHReMZyNdME+6mnnuKuu+4iLS0NuVzO6tWrue2226atXXDy5Enuu+8+/uu//ourrrqKnp4evvnNb/Kv//qvPP7445PaP/TQQzzzzDPs2rVrkkXjNF/60pc4ceLEJOvGLbfcEvh/WVkZa9euJSsri9dee40bb7xxxu8A4Otf/zrf//730ev1fP3rX5+x7cMPPzxrf9MRlNKgVColpUFi0eMTvTRYD2JUJKCSaZZcitWZkAtKNDJ/RdZMbRl6RQyd9pMcGn2VaGUyefo1xKsyJOXhHGH3Wk5lTppfJorFTJejDq/oocS4NdKiSABZ2jJaxo/Raa+ZoDSkaQvpczZzePRV1kZfN+1EO0NbwpHR1xl0tQe9qNLraKLWspckdQ4KQYXTa8PqHSFbt4I2exWHDv8apTIK11lB2tOlfM3UlmF2D5CYdsY61919iDFzO6kp68jI2Ixel4DDMcqJyifptNdISsNcuMCUBqPROEFpmIm8vDx2797N+Pg4ZrOZlJQUbrnlFnJycqZs/6Mf/YjNmzfzzW9+E4Dly5ej1+vZunUrDz74ICkpZ1ygf/azn/HDH/6Qd955h+XLl0/Z35e//GX+/ve/895775Genj6jrCkpKWRlZdHQ0BDUuVVUVATc0ioqKqZtN985gqQ0SJw3yAQ5SkFD4nm6+huvzqTYsIVMbSmCICNTW8qgq53G8cMcGX2NIsNmsnVTP6wkwkuiOpum8SP0O1tJ0uRGWpywY3YP0mY7QaI6JxDMKhFZ/K5IRbTaTlBo2ITiVFE9jTyKjbE3snPgSXqdTdNOtONVGcSq0jg6+gZb425Dp5h5oiWKIpVj75Kozg4UxzR7BgBIKb8S9UgOHZ0foNclMjB4kixd2YxWDIVMycroK6Ef6Pcv4izLu4aK44/T3XMQQYDCgo+h0URjMKRhG5y40GP1jFBvPYDbZ0cjN5CnXyO5zE2F7MJSGuaCXq9Hr9czMjLCW2+9NSGw+WxsNtuEWgcAcrk/iP906mGAn/70pzz44IO89dZbU9ZBEEWRL3/5y7z00kvs2rVrWiXlbIaGhujo6JigmMzEzp07p/x/uAkqEFqlUiFKSoPEIkcUReSCkmF394SgwKVuZTiNWqYjS1eOcMo/WBD8LkobY24kU1tKg+UADq8Vhzd0v2WJ0DApE4lWJtNqOz7h5XE+MObuZ//wi6hkOooMmyItjsRZpJ1yVeyxT1x9lAlyUjTLGHC2Me4ZnfKaFE7VZhARGfeOzjqWDy9ePCSpcwLZjqyeYeQo0GhMpKWuJyN9CwOD1WRlXkxh1EUoZaH5Xid0C1wWdzvZWZfR1X2I2rqXqat/heGRRpy+M+ko221V7B16DqtnGK3cyIirlw+GnqXesh9vkGmPXcUZEz7nK6IgolJJ8SBT8dZbb/Hmm2/S0tLC22+/zaWXXkphYSF33nknAN/+9rcDGY8APvrRj/Liiy/yyCOP0NzczAcffMB9993H+vXrSU1NBfwuSf/5n//JE088QXZ2Nr29vfT29mK1nnkPf/GLX+TPf/4zf/nLXzAYDIE2drs/sYnVauUb3/gG+/bto7W1lV27dvHRj36U+Ph4Pv7xj4d8nu3t7dO+l6aq3xAKQVkaDAYDoswzr4EkJBYaER8iXoZdXQy7uth0Ih6jcvaKp0sdQRBOFYCqZtfgUwCYFImkagtI0SxDKWjw4UFAPiEgUWJ+5OlXc2T0dfqczeeNz7/ZPcjhkdfQyg1sjP34tK4mEpFBJzeSqM6m1XacdG0xgiAEJsCZlWV02Wt5f+gZcnSrKDRMLgylEPyTSY84e/Gv0wHQ8lMWDZ/oo8NeQ0xcPoIgY3CwFq/HgUKuwedz4ynJnnWBxum14RLtCMjQyg3IBQUKmYplthxaeZfunjNFsQRkiKLIgKudk5b3ydCWUmS4CLmgxCt6aBk/RtP4EVpsx4lXZZCiWUaKJj9o94uzFQePx0nF8ccxRKVQ4l0RsOIsRUTBM+8CXucrY2NjfPvb36azs5PY2FhuuukmfvCDHwQsMz09PRMm1XfccQcWi4Xf/OY3/Pu//zvR0dFcdtll/OQnPwm0+e1vf4vL5eLmm2+eMNbZAdmPPPIIAJdccsmENn/4wx+44447kMvlVFZW8uSTTzI6OkpKSgqXXnopzz777Jx+y5ycHHp6ekhMnJhwYGhoiJycHLxe7zRHzk5QbwSTyYSXCyeFl8TSRCbIWR9zA+8PPQME92I8X1DI1KhletK1Rejl0fQ4Gqm17KXGcibYSi3TUW68nHj1zL6UEsGRoM4iQZVFrWUvcar0kFdZFyP9zlbcooP1puslhWGRkqtbxf6Rl+i015BxVt0GVflyLq7W0DxeQavtOJm6UrQfirnxKw0CjdZDdNvrUMo0lBi20e9sJUoRjVGZEGhr9vgLsJ12AbJ4hrB5x8hP+zjtHXtobPIXnsrJ3k5b2y4y0rfw4fVts3uAUXc/PrzEKFPYN/wCIAZkSdYso9iwiV5n86mxYsnUlmH3WmixVVBt2c2gs4PYmHyWLb8NryBwerqTSQ6plfkMuNrpcTRwwryDFttxCqM2EqdKD8l3e3y8D4ulE4ulkwHlSTabbkIj1wd9/GLCixuT6cKoUh8qn/zkJ/nkJz857f4//vGPk7Z9+ctf5stf/vK0x7S2ts467mzWaK1Wy1tvvTVrP8EiiuKU17/Vap02QDtYglYaLqQJmMTSRSOPwqiIRxDkxCiXdp2GUIhVpXJpwhmzaqq2AJfPzoCzDRERGXK6HHUcHn2VPP0a8vRrJatDGCg2bmHv0HNUmXey0nTVkg5EF0WRIVcnMcpUDMq4SIsjMQ3RqmTSNIXUWfeRoM6c6GNcWkDqiXFabcexey2TlAZBEFhuvIxBVwcun51uRz29jiZ8eBGQkajOIkNbRpwqjS57LSpBE6gyHaWIQSZTUN/wKlqNP7+8ISqNxMRyWlrfYczchpGJMQZV5l0B5SNWmcZphQH8izqd9pNEK5NIVufRq2pkyNWFRqZDrzChV5ioNr+HgIxVBTdMeW8pystIoYysmg5GXL3UWfdxePQfGBUJZOnKSdEsm1RIbipMpkzWrLqXIxWP4naPs3fob1yScHtQxy423KJLUhouUE5nTRIEge985zvodGfSsnu9Xg4cOMDKlSvnNUZQSoPRaMTlnbpghoTEYkEUfRwfe5txzyjrY6d+yVxIqGRa0rRFgb9TNPk0247SYD2E3Wuh3HjZBf8dzRed3Ei58TIqxt6kzV65pAPRR919jLh7WB19TaRFkZiFIsNmBlztNFgPUUjxhH0jiSIMg1ExtWtmqraAVG0BAA1WvztQqqaALnst/c5WKs07UMv0mD0DCAiBZ4RcUJCZsZXWtp04Hf7iVHb7EOPjfQDIOvtBO1FpSFRnB5SGGGUyBmUcouhDKVMz6u5jyNVJ6/gx0rVFZOtWMOjq4OiYv3puUdQmLoq9CY/oQqudvlAW+F2N9GSw4WQSg6522myVVJrfpc6ylzLjJSRqclDVdMwYy2AyZbF+/Vc5ePAXuEQHtZZ9lBi3zDjuYsTtcQSdTUji/OJ01iRRFKmsrJwQ26JSqVixYgXf+MY35jVGUEpDbGwsDpcNlqa1TuICYdDVQb+zldXR10xISSjhRxAE8vRr0MmNHB97B63cQH7U+kiLteRJ0uSQ6SqnwXKAZHXukq2aPOhqRy4oSVBlRloUiVlQytTk6FZRb91PyrG9GFeeCVi3WnvQauOCKo529v1fYNhIijafvUPPoZH5r+ESw7YJ7QvthbQJ7+ETPeh1iYiIdHTuRSs3UGl+d9LKfp5+LdGF61CrozE0j07oSxRFep2NCKdsJdHKJBSCCp/oJU1bSK11L+naElI0eXi9LuTy2c/HXZKJiUzW1mSdyra0n4qxf7JWdi1xqvQZFQefz8vRo48G/m63V1Js2LzkFlbsbhsxMVJWqQuR01mT7rzzTn75y18uiPIYlNKQnJyMw2XHK7oDQVESEouNAWcHGpmBBNX5UdhtoUjR5GP3Wqi3HkAuKMnWLV+SZvjFRH7UOnodDdRbD7DcdHmkxQkZs3uQ5vEKMrQlgexcEoubTF0pw65Ojoy+TvlRO/Gr/dedwzGKVhMbmByHkj3OoIjj0oTPoRQ0U06WvXjglH+2ShWFwuZlwNYa2N/jaJhg3XSXZHJ62uIqNkyQRRAEUjRnKvEqZCq2xt+KQlAhFxQYFPHUWj6g034S9cH3KFh2HfHxJUFN4l3FGUTVwErTlRwZfZ1jo2+zOe6TU8YpeL3+eE2fz4PHY5+wT0REQJjWR3yx4RHdOF2OoNN0Spyf/OEPfwD8xena29txuSaGF1x//fVz7jsopSEuLg65XI7Ta0enkJQGicWJSqbBLdrxiC6UwtIPSl1IcnSrcPkc1Fv30zJegVqmI1e/OuC2IBEaSpma/KgNVFt2k6tfveTyx/c4GlAKKinF6hJCLihYFX011eb3OGHewaoKJTGrtuF221Crz/i0h6o8zFSXY8jZiYgPvTwGhc1Lpq6UAVdbYL/bN7Mbs6s4Y0Y51DK/D7YoiqRri0jXFmH1jNBgPUBl9dPExRZSUHA9Ws3s99fpsVaYtrNr4Cl6HY1k61egqunAUZiCz+dBodCw78D/oFLp8fnOZIiMV6WTo18NwLsDfyRLW05e1JpZx4w0Lq8NhVxBXJwUk3Qh09LSwsc+9jEqKysRBCEQiH1a8Z1P9qSglpRkMhmxMfET8iZLSCw20rXF+EQfXY46gPM+H/d8EASBIsMmtsTdQpZuOVq5kUrzu9i9lkiLtmRJ0eYjIGPI1RlpUUKm39lGgjpbsjgtMWSC3O+zr872Bx2bO3G6zCgUakTRN6FtOJ6Hdp//+eDy2fyubOosig1n/P6jFOGZrFaa32XXwJMAGJXxrI7+CKtMV2Ed7eTk4T8F3Y+rOAOVTIteEU27vYpeRxOiKFJ38Cne2/M9Ptj3EC6XGau1B5vNX7hOLihYbtpOnCoNER8un52G8YP+dLGL/J3i9NmIjYlfElYRiYXjvvvuIycnh76+PnQ6HdXV1bz33nusXbuWXbt2zavvoO3QSYlJOH22eQ0mIbGQyAUFMkGO2+eYsH0xP+QjTZQilmVRayk3XYaIyIirJ9IiLVkUgpJoZRKDzvkVzznXOLzjjHtHSFBLsQxLEUEQKDVejCAIHD76W+z2Ibq6D7Bz93+y78D/MDzSOKF9MBPf6fZnaEvI1JaSqM4mV78qsC1Lt5zV0ddMSuf8YatCMNYOm9dMt6Mel+jAfqpQpSAIJGlyKTJsYszTj7uyKmjLias4g+XGy9HKTRwb+yf7hl8IFMB0OkdJTlqFUul3W0pS57Ih5uMBa4tcUARiLlptx4MaL5I4fTaSkpIiLYZEhNm3bx/f+973SEhIQCaTIZPJ2LJlCz/60Y+477775tV30EpDTk62tAopsajxim4EBLrt9VjcQxP2SYrDzKhk/tSKo+6+SIuypNHKDYx5BpZklejZXEskFi9qmY7NcZ+kzHApAKkp6ykqvBG12sSx409Q3/CPgOVhcLCGyqqnaY8dxl7gn2BOVS15KuVCLigoMW6j3HRZoHCmTJBTbNhMojp7StlUNR2BTzCMunoD/5eVTHSXTFBnISAErHnB9mlQxrE2+loKojZg9gyQoM4iW7cCGXJ6+yrwuh3k6lZNOK/TZGj92alabMdoa9uN0zkW1JiRwO41k5ObE2kxJCKM1+slKsqfzCA+Pp7u7m4AsrKyqKurm1ffQVfvKSop5MTu5nkNJiGxkGjkUWyK+wQVo29xaOTvrBpPRa8PLYvSbD635zPRymRG3b2zN5SYlnRtMd2OenocDUsmPkQt05GgyqLashu9IppYVWqkRZKYA2qZjr5TRdLc7nGSk1aSkryauvq/09m1l86uvSiVOrxeNzKZgoHBauQyFamp68ln+kWVUJ6Jcwm+/jA+/MpNTvZ25HLVhPEVghKlTDPBmjxbKtXTcg1V7KTeegDwV7kuMqwiR7cSi2eIKEXMtFnPMnXltNurAWhqeYumlre4NOFzgfiLxYQDK0XFS+O5I7FwlJWVceLECXJzc9mwYQMPPfQQKpWK3//+9+Tm5s6r76CVhmXLluFRSO5JEosbndzIupjrODj8CpXVT7Nh3VcD/p2zvfzOXmGD+b34liIxymS6HfV4RDcKKUvanIhVpRKtTKbP2bJklAZBECgzXsLOwT9h95oBSWlYqmTqSlHKNHQPVtPU/BYpyaspLLie2Jg83B4bbrcdj8dBWqo/1WpH5146OveQnLQSg2H6332mSXm4n5OJq69gi+silMrJk3KrZwRRFHGJE7McTaU4iKKIKPqQyU7F6WSlwQnYHPvJQPFCtVyHWj7z5D9KEUNZ6W1Yrb0MjzRhNrdRMfpmoKJ1+lnZoiKNR2Fj2bJlkRZDIsL853/+J+Pj/hjkBx98kOuuu46tW7cSFxfHs88+O6++g1YacnNzsXnMICWlkVjkqGRaio1bOTTyd8bGWomOPmOuDXXV7EJSHPxBsEvPrWax4Ve+GiItRkj0O1sBgQS1lK54KZOgziJelYFcUNDddZCOzg9QKvW43eOo1dEkJ60kM2MzKpW/UnRe7lX09R/n0JHfoFIZMBozyM/7yKRiasGmHA3X81Klmrjq7yrOgOp69g+/iFJQk6zOm/ZYUfTR3XOYjs4PcDhGSU1ZS2bGVoyGDFQqA5XmnayPuQGFLLiFEVdxBolkkJhQRmbGVurqX0Gh0NDVvZ9BVwdWzzCFURsXRaricc/YvFeSJZY+V111VeD/ubm5nDx5kuHhYWJiYuYdJB+SpWFsfBifzitl2JBY9MQqU9HpEmhoep01q+5FJgv6Up/AhaQ4jHtH0cj0kpVhnkQrk2ixHcPhtS6ZQm8e0Y0c+YzpNiWWBoIgo9S4jWLDZoZcXQwYLGg00YyP99PZtZf2jveJNmVhNPon0akp6zBbOonSJzEweJKDh39FYcHHSU5agcMxSkvrDnr7KkhL3UBa2kU4HaMMDFaj1ycRF1tAk2UPbbZK1sfcEHBtE0URm9eM1TPMiLuHZE0e0crZA3Rnsmi0jB/HI7rI1JZyaOw1oqJSSE1ZS2rK2kAbn8/LseOPMzrWRkJCKQnxpXR1H6C3t4KLNn6TFeWf4+ix31OvqWNZXuiVzxUKNaUlnwTAaFFSY3mfVttx1DIdXtFNmrYIrdww5bks9HvEJ3oxj49IloYLHLfbzZVXXsmjjz5KQcEZa3ds7MxV1YMl6JlUZmYmKpWacc9owLQnIbFYEQSBFZpt7Bt+kZ7eowFzPISuCFwoisO4ZxSdIjrSYix5ohT+h7PZM7hklAa36EAuKJZMESuJ2ZEJchLUmSS4ABe4ijeSm3MFvX0VDA830Nt3HJfLAqKIiIjF0oVel4jDMUJn1150uniOHH0kEEDd2bUPp9PMwGD1lONpZP4MRCOuXirN72LzngkYNrsHWB97w7SyBpOoYkQ5AkCzrQK9PgmlQkNt3YtoNTHExPgtDw7nKKNjrRQX3UxKsr/OQlrqevbuf4ih4TqSk1aSnLSK/v5K8nKvnte1nrLmGsSeBGrrXqTDfhKbdwyf6CV7/SeCOsdwv1OsnhHUajUZGVLSjwsZpVJJVVXVgj3Hg1YaZDIZBfmFWNqHJKVBYklgUiYSH1dEZ+deUlPWTbiJJMVhIqIoYvEMEqdKn72xxIx0O+qRC0pMitCC8CNJnCqd5vGjjLp7iVFJ1WTPR1Q1HaiAPDLIU2ZAtP++dxVn4HKZaW3bhXW8j6LCG0lKXEF1zbMTaj3I5Wo0UxRVUwpqsnTLUcm01Fv202o7gY8zxaOilcmUmy6bUqZQstqtXHFXIIjbbzkWqTj+ONU1fyM2Jh+tNhaPxx8grdOemaNoNNEYDGn09B4lKXFFwPpwtOL3lJR8MqhCcdORlLgCr8eJ2+OgtW0H/bJesoM8Ntyxc1bPMAX5RchkkXeTkogst99+O48//jg//vGPw953SD4bq9euYnfzsbALISGxUOT48jlkq2FktInYmPmZbc9+wQUTUL2UlIweRwM2r5lyzdII3l3MnK5IrpJpz7mLwlyJVaailunod7ZKSsMFhCAIqGs7UQPlrAYNMAqM9pLlycGgU6FXxKAvXEFr2y66ew4hl2sozL8efZ8DvTwajVyPKIocHX2DIVdnQGEwKRNZpl9LvCpz0qrnXFJgn1EWAtJTUvQJGhpfw24fZGCwClEUycq8BKNxYv85WZdxouop+vqPk5hQTlRUKmPmNtradlGQf/2ZYOkQ6e45REPTawG5rNYe3G7blEHc0xGuBSmrd5hL1qyadz8SSx+Xy8X//d//8fbbb7N27Vr0ev2E/Q8//PCc+w5JaVi1aiXvvLRnzoNJSJxrYpWpROmT6ezcO0lp+PDDOpjUfWcfe/qYqbYvJZxeGzWWD0hS50oTxjAQr8qgzVbJSJaaDzsnLVaLlSAIxKnSGXR1UCBulFyUJPyuTacK/jm0cVjHe/F6/bU86hpe4dLYz6CQqQB/IP2Aqw29PIZx7wgrTVeRrJk6IDecz0iNJprysk8D4PE4EUUfSuXkuJz4+GKSEldSV/8ySqWO8tLb2HfgZ3T3HCJKn0x6+kVzE+BUPRafz4Nen8z4eC9yuSrkbsKx0ORSWlm5euWcj5c4f6iqqmL1ar97Xn19/YR95ywQGmDVqlWMufoR1aH5vYbjRblYX7YSixtBEEhP30Rt3UvYbIPodBML98xHcTh9/Ez7Fvs1K4oiJ8zvIggCJcatkRbnvCCqfAPyve/Q33+CqJwrJu1frNdFsiaP7tF6WmzHAtV+JRYv0z17FuLa0tR1s2bVvdTWv0R//wm8XucEF6QBl78K+rh3hHz9+kkKw7lYTJHJ5Mhk06d3LCy4AWtFL8dP/BE4M3+pb3wVrTaWuLjCkMdMT78IrTYO63gfzS3vAH7lRaWae+INCP03FEWRUWcfK1eunNO4EucXO3fuXLC+Q3NPWr0aq8OCZYURjSZ60v6ZVl3no0lf6PnzJWZmtusiYySOZlUUbe27KS66adbjQ1UcljJtthMMuTpYE33toixWtBSYVDUXSEleQ0fnB2RlXoJcPjkb1WJUHBLV2SSr8+h3tkpKwxJmoa4tXUM/q9iELa4Uj+iekGkrS1uGTm4kUZ2FqnwFrhD69cdNCPNaAW1pfZeW1nfYtPFbU85NwJ/5aN2aLzJu68di6aK3t4LRsRYAtNq5xWkKgoz4+CLi44vo6j6A0znK/oMPs23Ld+Z6KsDUStbZv+mH9zsco4wfsLJmzZp5jStx/vD+++/z6KOP0tzczHPPPUdaWhpPPfUUOTk5bNmyZc79hqQ06PV6CvKLMJs7prwxz558TTfpCnXiP1U/kvIQfpaK7/VMTPeylAsKsjK20dj0BtlZl07KQT7V8eFSHBbzter2Oam3HiRTWx5wQ1hMBBtDcq4I5XqIiy2ks2sfbvc4cnl0SP1F8lxlwtxWSCXOLZFa1BAEAf0UGdYMyjgMyriQ5erqPkRd/UtoNbEUF900oaZOsLhcVlpa/av8I6PNgaxJUyGTyTFEpWCISiE5aRXjtn60mhgUCk3I436YspJP0dm9H70uYd59TcVM363Z3EFBQTE6nbTwIwEvvPACn/3sZ/n0pz/N0aNHcTr9boUWi4Uf/vCHvP7663PuO+Q3xOYtF7HjnXoSE8unbRPMgyOYyVQwpeEXw2TibM71JHG2FYnp2syl36n6XmxMd02kpq6jtX0Xre27KC68ccbjF0qus1kM32OvowkfXnL1KyMtygSmWygI9TsL13c+l2vC6/Ovtc6lPkgw9/RUbcNxTbl9DlSy+U+gJCLLVO+hcFwr091Tc31uyk/dH3bHMMcrn6S05Bai9Mmo1SbcbhsWaxcmYxYKxfRuR719xwL/93qcQY99WoEIFyZTJibT5MWX+X5HwWC2dHLlVZsWrH+JpcWDDz7I7373O26//Xb++te/BrZv2rSJ733ve/PqO3SlYfMmXn/tvXkNejbzDSg9ezIRyss23Hx47JmUh4VeJVrI/md78SzE6nCofU41wdTW95GVsY2mlrfIzrwUrXbuafbCQaQtEKIo0uWoI16VHrFaAqFep9MpDvMJXp9JOZ7PfWS19iKTKVEq9bM3DoJgF2Lmez358KIg9EBOR2EKx078gfi4IjIz/LExsyk6i0FxPt+ZzeJ/NrP9HjNZ/UPFbh+ho/MDFAoNarUJvS4Rr8/FiconAQJBxQCpKespKvzYlP14vW46OvaQmLicwYGTWE8ds5g4F4ubTlcPmzbduaBjSCwd6urq2LZt26TtRqOR0dHRefUdstKwbds2Bgbvxet1zSlLwHTM5wU9WzDqfJnLwzSc4y9WgrEEzZXpFMFgH8BTtUtL20B75x6amt+krPTWOcsWTmabCC+E4gUwOFjDaH8vq6Onr4q6kNamuV4b4XQZW6hxnE4LnV17SUleteSyECkEJV7RHXT7099TT/chRkdbGB1tISV5zbQpJ0OJcwvGUrTYXNjmSyQVqnP1rhoba+fYiT8gCDJE0YfJlEV21qWYjJmYzR0MjzQyOHiS3Jwr6Ozcx/BIw7R9DQ3V4nSZycm6FId9GI83eEvD+YLX62JwsHXKSaLEhUlKSgqNjY1kZ2dP2L5nzx5yc6fOahYsISsNeXl5xMXFMzbWTmzshVGu/Hye+C9WglHEZos/+PCEXC5XkZd7FTW1zzMysj5QRTTShHKuwRw/20q62dwZKIhkWr4FajuDlieY/WfLcSHeO3X1LyMIcnKyJ2dOigShWAeFfTJExJDHGB/vRanU4fE46es/QXraxjmvsIayOj7V/qWqPHxYoQpWmVpKeDxOqk4+Q5Q+iRXL75gUS3DaxSc761KGhutxe2zkZkx/Hw0N1yOTKVGro4mOzqW3r2KhT2FOLORvNjbWRlxcwrwngxLnD/feey9f+cpXeOKJJxAEge7ubvbt28c3vvEN/uu//mtefYesNAiCwGWXXcb+fc0XjNIgsTgJZlLy4eDm5KKVdHUfpL7xH6xfe9+SWQmebuI3F5eBru6DeL1OCgtuQBCEwHcU7hfbUp7czJWxsTYGh2ooK7kVlSo8rkmhMJu75ocn1h9u49YrkI97Jx033Vg+n5fR0WZ6+ypIiC9lfLyfMXM76Wkbpx072O1zZSkqDwtt/VostLXvxukcoyD/+hmDjzu79tHQ+A9iYpaRkb55yjY+n4f+gSoyM7agUKjx+dxhCWheaoyMNnP55ZctmXeZxMLzrW99i7GxMS699FIcDgfbtm1DrVbzjW98gy996Uvz6ntOqTKuvHI7O3f+cF4DS0iEg1B9vQVBRnxcMc0tbyGKXoQlmC1mvpMHi6WT1JR1E7JInQ8TksVAb98xNJoYEhJKIyrHXK1FXp8bhckY9DgjI40cr/wTIJCXezWNTa9jtw9NGCdSFqel4ro0H8teMPE55xq3205H5wcMDp4kLq6I3JwrEASBkdFm2tp3kZy0irjY/Bn7GBtrA2B52WenTSYwZu7A63WSEO+/1yyWLqL0yWF3nV7sOBydXHnlHZEWQ2KR8YMf/ID777+fkydP4vP5KCkpISpq/vGLc5oxXX755QwMfB73MvuU1RclJBYbZysO0dHZAAwO1ZKYUBZBqSKDIMgjLcJ5iSiKDA3XExdbiCDIIibHfCaPos+DIMiDjx1yWQFYUf45VCo9MpkCnzjRUrEYJrPhkmGuboLBHhPsvqlYDO6Avb0VnKx9DoC4uCLa2nehUKjJSN9Mbd1LmIxZFBfdNOP9MThUS/9AJUUFH5+yxslpnE4zABqNP6mFXK6mf6CSgcGTbNr4TdTq4JXfhcbn8yKThf+563bbGRhs5fLLLw973xJLH51Ox9q1a8Pa55yUhszMTPLy8hkeaSApcXlYBZKQWChOv1CjTVnEROfS2PQGRkP6tMWAzldUqqjAZE8ifLjdNhyOEWLmkGt+saBWm7CdshTMtFLvKs7AYummsflN4uOKiIsrAPxpLH3e4AOplxpTKVMzuV0tdFzCbO5o55KOzg9oaHwNgLLS20hMKKO+4VVa23ZhtfZitw+Rn/eRGRUGr9dNTc3zxMUWkpIy82RnbKwFnTY+4JJkMKQyMtqMKHpp79hD/rKPhO/k5oHDMcre/Q+Rv+zaaV2t5srwSAPLlhWQkRF5xVwisnz9618Puu3DDz8853Hm7Jvx8Y/fwF+e3iEpDRJLkqLCG6k4/jhHKh5lRfntRIUxX/diR6WKYny8L9JinHdYLF0A6PWJEZZkboiijzFzO3GxhZP2fXhCarF0U3Hs/9Bq4yg6VfdEFH34fF5stkHcbtu0GZSWOqFMzj+sOCzExD7SygLA0FA9jU1vkpa6gfxl1wZcirKzLsFi7WF4pBEAxSzXhM/nweNxoFabZvXRt9mHUamNgXa5OVeQm3MFdfWv0Nt7dNEoDafPvb1jT9iVhrGxBj79mRvC2qfE0qSiIrgkAPONfZmz0vDRj17Hb37zCKLoi6gpXkJiLmi1saxZ9S8cP/Ekh488Qv6ya0lNXX9BBJOpVFGMjDRFWozzjsGhGjSaGHS6pak0DA7W4HSOkZqyZto2Xq+bnt4jtLa+i1YXx6oVdwdWeltad9Ddc+hcibtkWAyT+oVkfLyfE1VPERebP0FhAFCpDKxZ9S9AcC46SqWWZcs+QkPjPxAEyMu9ZsrCbmNj7YyMNFJS/MnAttPzEKu1d9FkxgP/9wOErWbLaUTRx/BIPR/96C/D2q/E0mTnzp3nZJw5Kw0bN25ELhcwmzunrIIoIbHYUatNrFn9rzQ2vU5dwyuMjDZRWPDxCyBOR5hDUk2J2RgbayM2ZtmSVDx9Pg8tre8SHZ2L0Th5kmu19tLbV0Ff33FcLiv6qGSWl31mQrYau30E4FRRu/PTyiBxBlH0YTZ3cKTiUUCgrPS2GSugB+vTn552EYIgo6npTcbHB1i54s5J/ba170avS5zS08HtHid6EbkIZmdd5i9kJ5++qnWoeDxOGhr/gVwuY+PGyZnKJCQWijmbCBQKBddffz2DQ1XhlEdC4pwilyspLLiBspJbGRpu4MChX9I/UIUonr/Tao/HjlJxvitG5xZRFLHbh9Hp4iMtSsj4fF5q615m3NbPstyrp2zT3PI27R3vExOzjA3rv8r6tV9CrTZNaHM68HTj+q8tuMwSkcVuH+H9Dx48pTBAWur6GRWGUBAEgfS0jaxYfgdj5jYam97A5zsTXG8d72NwqJb09Ium9HLweB2LKvWqKHoAkdGxFrxe17z68vm81NS+yHt7/pue3qN89KMfRS6XEltITOb999/nM5/5DBdddBFdXX7X2aeeeoo9e/bMq995+RV9+tO3MjJagyj65iWEhESkSUwsZ8O6r2A0pFJV/RdOVD4ZWDk93/B47CiUi+elej4gil68PhdKxdJZYRdFkf7+Sg4c+gW9fRUUF96I0Zg+qZ3H42R0rIXsrMsoKb55SsXI63UzPNxATMyyCy6xwIWGx+PgROWf8HgcABTkX0/+suvCPk50dDZZmRfT2bWPIxWP4vH4qz03NL6GVhtLSvLUbnSCIMdsXjwpdlWqKDLSt5Cbc8W8FauxsVZ6eg8D/oXbz3720+EQUeI844UXXuCqq65Cq9VSUVGB0+m/dywWCz/84fzKJcxLabj88svx+VyYzZ2zN5aQWORoNNGUl32W8tJPY7X2cODQL7BaeyMtVtjxeJwo5JLSEE5kMgVyuRq32xZpUYLCYunmSMWjVJ18Bq02jnVrvkhy8qop23Z27cXn9ZA6QzabkdFmrOM95GRdtlAiSywSWlrfxeEcwxCVBvgnxQuRThQgJ3s7pSWfwmLppK7+FZpb3mZkpJFluddMOwE3GTNxe+wLIs9cyV/2EbKzLp1z/Gdv3zH27P0xxyufJD6umNzsK9BotFKqVYkpefDBB/nd737HY489hlJ5Jm3xpk2bOHr06Lz6npfaq1Kp+NjHPs7771VJcQ0S5wWCIJCQUEpMzDL24zHf2gAAYe5JREFUHfgpDU2vU1r8SVSq+RdFWSx4PHY0mtjZG0qEhEKuxuN1RlqMWRkf76fi+P+hVkezcvldxMYum7W90ThzamKbbQCZTCm9By4AhobrSEwsZ1nu1fh87kluauFEEASSEpdjsw3Q1rYbEZGszEuIjy+e9hivz31eFXczmzs5WfM3dNp4XC4zUVEpeH12Pv7xj0+YEEpInKauro5t27ZN2m40GhkdHZ1X3/NOe/SZz9zG0HD1BJ9DCYmljkKhpqToE1it3Rw49Eus1p5IixQW7PZhxswdRJuyIi3KeYdKFcWYuX1Rx8O4XFaOV/4JtcrImlX/MqvCAH7XK2EWt4qR0WYMUalSJr3znN6+Y9hsA8THFaNU6hZUYTibnOzL2XTRN9m08Zvk5V45bbIBl8vKyHAjcbEF50Suc4F13G/tdjjHkMmUpCSvZmi4mk9/+tYISyaxWElJSaGxsXHS9j179pCbmzuvvuf9hN++fTsajVJK4Shx3hEXV8j6tV9BJlNMSiU5lqdmLC982TDOFQ7HKCAuqpSE5wvZWZcxMtLI4GBNpEWZEut4H0ePPYbP52H58s8FHSwql6vxzODu4XbbGB6ul2r2nOf09Z/gZM3fSE5aRXzc5FoeC41KZZi1yrPH40DER5Q++RxJtfAkJa6gIP+jZGddyro1X8RmG0CjUbJ9+/ZIiyaxSLn33nv5yle+woEDBxAEge7ubp5++mm+8Y1v8IUvfGFefc9baZDL5dxxx+cYHD4+364kJBYdarWB2Nh8Rkaap9y/1BSH00kLBEHKuBFu4uOLiY3Jp7n17UiLMgmXa5wTJ/6EgMCqFXeh1cQEfaxOl4DNNjCtNbm/vxJEfzIBifMTu2OE+oZXSUwoo7jo5kVrUVKpDMhkSppb3sblskRanLAglytJT7uI7KxL0OsTGRw+zl133SFlTZKYlm9961t87GMf49JLL8VqtbJt2zbuuece7r33Xr70pS/Nq++w3Pl33XUnfX1VuN2LK/hIQiIcxETnMW7rn/YltJQUB6/Pn/Jvsb70lzKCIJCcvIrx8T5cLmukxQngdI5x7MQf8PrcrFh+B3p9UkjHR0Wl4PW6sDuGp9zf03eU2Nj88yruR8KfMMHhGKWj8wOOHHkEuVxF/rLrFnUdEoVCzaoVd2N3DHO04rFAhqfzBbfbTl9fNXfeeWekRZFY5PzgBz9gcHCQgwcPsn//fgYGBvj+978/737DMnMoLi6mtLScvn7J2iBx/hETnQ3A2NjiSeM3V3p6DqPTxqNWGyItynlJTLTfX3R4ZLI/aaTo7DqAwzHMyuV3zikd6vh4H4Ign9I6MT7ej9ncQUry6jBIKrEYsNuHqal9gff2/Dd79z9EY9ObmExZrF39b7O6By0GTKZMVq/8PC6XlZM1f1vwlPBWay/HT/yJzq79CzoOQF//cUpLyykqKlrwsc5HLBYLX/3qV8nKykKr1bJp0yYOHZq5iv3TTz/NihUr0Ol0pKSkcOeddzI0NBTY/9hjj7F161ZiYmKIiYlh+/btHDx4cFI/v/3tb8nJyUGj0bBmzRref//9CftFUeSBBx4gNTUVrVbLJZdcQnV19ZzO80c/+hFPPPEEOp2OtWvXsn79eqKionjiiSf4yU9+Mqc+TxO25cYvf/kLDA0dnTIIcKn6f0tIAKhURpRKPRZrd6RFmRcWSzeDQ7VkzSP1n8TMqNVGDFGpDA3VRlqUADZbP4aoVAyG1JCP9XictHfsISG+ZMoUl719FSgU2hmz2UgsbkRRZHColsamNzh2/A/sO/AzensryM66jLLS29iy6f9RXvbpJWVJ0uniycu9isGhWoaG6hZ0rNGxFoaG66hv+DsDAycXbBz/73SU++774oKNcb5zzz338Pbbb/PUU09RWVnJlVdeyfbt2wPFzz7Mnj17uP3227n77ruprq7mueee49ChQ9xzzz2BNrt27eLWW29l586d7Nu3j8zMTK688soJfT777LN89atf5f7776eiooKtW7dyzTXX0N7eHmjz0EMP8fDDD/Ob3/yGQ4cOkZyczBVXXIHFErqb3aOPPjqlYllaWsrvfve7kPs7m7DNHG699VZcLjNjY60TtkvKgsRSRxAEDFGpS15p6Oj8AK0mVgpYXWDi4ooYGqqfd/XXcKHTJWCxduP1ukM+trNrH273OHlTVIoWRR+9fcdISiwPWzVgiXPL+Hg/Ryp+x4nKJ+no3Mu4rZ+iwhvZsvl+cnO2k5hQhlK5NKvHa7R+y5hKtbBW1dP1KuRyNfWNry5YJsnRsVbcbjOf+tSnFqT/pYrZbJ7wOV3I7MPY7XZeeOEFHnroIbZt28ayZct44IEHyMnJ4ZFHHpnymP3795Odnc19991HTk4OW7Zs4d577+Xw4cOBNk8//TRf+MIXWLlyJUVFRTz22GP4fD527NgRaPPwww9z9913c88991BcXMwvfvELMjIyAuOKosgvfvEL7r//fm688UbKysr405/+hM1m4y9/+UvI30lvby8pKSmTtickJNDTM79MkGFTGnQ6HXfdfRd9A4enbSMpEBJLFYMhdca0q0vh2nY6xzAY0xesEJOEn5Tk1Xi8Tnr7Foe7ZkryGjweB8P/f3v3HR5FtT5w/Ls9m94LAdIIoUmTjhSxoKJYEFFUFAT14pWLBRX9qSgoyvXaULEhCiJXFBELlyYCIlKkCQSQkoQkpPe2fX5/rERCekiyCbyf5/GR7J45c2Z3dua8c1resXptl5efwMmEdfj6RGI0Vu6alF+QiNlcQEhI1YvCiZZJURw4HHbSM/b9tdKyiZ49JnH5sNkMGvAEbcL6tNpA4Wy+PlGoVBoKCk/Vnvg8ONcw8cPNzReLuYjUJuqmlJn1O/fddx/u7q1n1fnqeKUoeJ06z/9SnL1a2rVrh4+PT/l/c+fOrXKfNpsNu92Om1vFWeOMRiNbt26tcptBgwaRkpLC6tWrURSFjIwMvv76a0aNGlXtsZWWlmK1WvH3d66FZLFY2L17N1dffXWFdFdffTXbtm0DICEhgfT09AppDAYDw4YNK09TH+3atePXX3+t9Pqvv/5Kmzb1b3E+W6M+Hnrooam8++67RLa/BoPBq1VUpISoCzc3X8zmQhTFQWGH1ndDVRQHVmspBoOvq4tywTMa/QkM6ETSqU2EhvRw+UJTRmMAADk5RwkK7FLn7TIy9mN086f7JROqfd/NzQ8fb1nQrbVwOGxs2vJc+d+BAZ3p3OnW8iChJQ9yri+VSo2i2LFaS5p8P5ERl3Pk6DcYDD4kJP5EWFgftNrGq/+YzUVkZBzkoYeWN1qeF4rk5GS8vf8ea2MwVP25e3l5MXDgQGbPnk3nzp0JCQlh2bJl7Nixg9jY2Cq3GTRoEEuXLmXcuHGYTCZsNhujR49m/vz51ZbnqaeeIjw8vHxK3OzsbOx2OyEhFSegCAkJIT3duQbHmf9XlSYpKamWT6CyyZMnM336dKxWKyNGjADgp59+4oknnuCxxx6rd35na9SOzbGxsQwdOpzTaTsaM1shXE6ndQeUFtPlpD4UReGPA0soLskgIODCWfSoJesQcx1mUwEZLaC14UxFsL4TVRQVpeLrG4VGU3nVWYfDRmbWAUKCe1xQFc0L3YmT68r/3af3VLpfcvcF0apQlTPBQmLSz/x57Psm3VdYaG8CAzpjNhdgs5vIz696iu6GOp22nWHDLqdDh9oXY7zYeHt7V/ivuqABYMmSJSiKQnh4OAaDgbfffpvx48dXO31tfHw806ZN47nnnmP37t2sWbOGhIQEHnzwwSrTz5s3j2XLlvHNN99UatE49zqpKEql1+qSpi6eeOIJ7rvvPqZOnUp0dDTR0dE8/PDDTJs2jZkzZ9Y7v7M1+mjIZ555ivSMXa2yciVEdbR/3Vit1tJq07S0ljWbzURm1kFSUn8jJ/conTuNkfEMzcTdPQCDmw+lZTm1J24G3l7t8PWJrNc2iuJAra4cMADk5h7HZjMREtKjEUonmoOiKCSnOLssxMaMwtu7rYtL1LT0ek/atxsCOMfmNCWVSk276+8hMG4A4X2uw78RV6S22cykZ+zimWeearQ8L1YxMTFs3ryZ4uJikpOT2blzJ1arlaioqCrTz507l8GDBzNjxgy6d+/OyJEjee+99/jkk08qjQ147bXXePnll1m3bh3du/99nw0MDESj0ZS3JpyRmZlZ3rIQGupcjLCmNPWhUql49dVXycrKYvv27ezfv5/c3Fyee+652jeuRaMHDZdffjnR0ZGcTqt6bENLq1gJURc67V9BQw0r47YUBYXJ/HFgCb/8+hIHD33BseM/oNEYCPBv/lVcL2bOGaqc/W5dOYOcoihYrMW4uwfVazuDwbvS2gxnjiEz6yDu7kF41nPNB+E6xSXpgEL3SybQrt1gVxenyalUajrEXIubmx9twvo26b4KYgyotToiLhtLaI8rGnXcWFr6bmJiohk+fHij5Xmx8/DwICwsjLy8PNauXcuNN95YZbrS0lLU6orV5DOtEmfPFPrvf/+b2bNns2bNGvr06VMhvV6v59JLL2X9+oqLfq5fv55BgwYBEBUVRWhoaIU0FouFzZs3l6dpCE9PT/r27Uu3bt1qbIGpj0af8kKlUvHss89w75SpdHEMRdWAH09BjAGfE1WPgBfCFbRa5+Azq7W0xh+NK89dm81MRuY+jh3/EaMxgJjokQQFdcWg9wJUMgC6mSkOOxb/isGCK84Pq7UYkykPna5+U2Z6eYWTetrZ1fTsY8iL1JC9NZ624QMbtZyi6Tgcdo4d+x43gy/+fhdWF5eagnHF4cC0uQDCAhr02zs376q2r2r/9d3X2XmcvZ3DYScjczsfffSudANsBGvXrkVRFOLi4jh+/DgzZswgLi6ufLG8mTNnkpqayuLFiwG44YYbmDJlCgsWLGDkyJGkpaUxffp0+vXrVz6geN68eTz77LN88cUXREZGlrcWeHp64unpvOY++uij3H333fTp04eBAwfy4YcfcurUqfJuTiqViunTp/Pyyy8TGxtLbGwsL7/8Mu7u7owfP77exzl37lxCQkKYNGlShdc/+eQTsrKyePLJJxv2AdIEQQPAmDFjePhfj5KXsB//mMqL/lT3gzr35gpV/0iFaArVXbgBdDpn0FDob8G/Dvk093l7Om03fx77DofDRkjwJXSKG1NlX3RRN9VVROr6vRbEGHDsVqNUMf1ic58flr+61Ol19Zt1RaVSo0JV6bPITzqIzWYiWLq6tRqZmX+QX5BIr56TL6jpcWtrvVOp1fjH9Cbr6G+EdL+8Xr+98w0Garqf1HUfiarDGIx6brnlljrtU9SsoKCAmTNnkpKSgr+/P2PGjOGll15Cp3PeK9PS0iqsnXDvvfdSVFTEO++8w2OPPYavry8jRoyosEDae++9h8Vi4dZbb62wr+eff55Zs2YBMG7cOHJycnjxxRdJS0ujW7durF69moiIiPL0TzzxBGVlZUydOpW8vDz69+/PunXr8PKq/5TBH3zwQZVTtXbt2pXbb7+95QUNWq2Wl2bPYsYzs/CL6olKLQtJiZbt3Iv3uTcHjUaPWqXFZq5+TEN98j/X+VQi09L3cOToN4SF9iYq8ooGrfornGr7nupSaTiTh97DD0tJXr3yqcuTzerKeW5au91KSupvJJ3ahIdHaL0r+Xa7BVUVFcys+F/xDI1xadekunT1kgdOf8vOOYKHR0j5iuUXgrp29wvtMYLc47vJO7mPgNg+tf6G63oNqE93w/p2TSyIMaA4HOR9/zOvzX0RrfbCCfRc6bbbbuO2226r9v1PP/200msPP/wwDz/8cLXbJCYm1mnfU6dOZerUqdW+r1KpmDVrVnmgcT6acp2GJjsT77nnHp6b9QJ5Cfvq1NpQ3Y+qqh94Q5oM5Qbyt9oqHI3x+dZ0kaztglvbeVHXfdX1mOpy7qlUKnQ6d2ymuk3f15CbxLnlqkuepbmnObx5BW3C+hDX8cZWvdJzfT+zupwHZ9I05niC6p4gnrsPvacf5sKsOuVzPmnOTet93ERG5h8cO7UWW1kRgR370+bSkZS4edbrOujl1YZTyVvIObYL/w59KE4/ScbBzRRnnCR6xARw1DmrSmU8V1N/T1Xt62Ki03m4ughVquq619AuQdUx+obgEdye/KSDBMT2qff2VWmO8Ul5CftwN2iYMKHqKY+FqM6ZdRrOHeDd4tZpOJter+fFWc/X2NpQ16a+s9PVpcmwujTNpaluSk15DDXl3VSVm6Z6vy7pGvJZqj096xw0NFRdguSzpe1Zh97Tn9BrxlKo1rSIClFz/dYa+7xs7DIYvPwoTD3apPs/l6I42J+6ktzjv+MbcQnh/a7HzTuw/P26dJk4Q9+vF/7moyRu+S9p+3/CXJCF0b8N7QeNwTeyOwUqVb3O1Zo05/W5rvuqSytPS/i91aasLJfs4uM4lNpnNKzPw7bzaT2tqmW3prTnG1R6BkeSl/BHg7Z1BcXhIPfgz7w29wX0eteu8yJan6Zcp6FJ27zuuecenn/hxWpbG6BxKoD1yac5tKSyiPNz9g1L6+aJtaywWfZZF+bCHPKTDhAxZFz5hAONMdhPNA6jfzi2siIsJQXoPXyaZZ+l2ankHv+dtv1HE9JtWI1pa/ve1UDU8Lvwj+5N5uFfCetxJf4dLq0wIPNCPnfqG5S2tADiTNkOLP8ES5FzFqz8KG2FyUlqq4yf2xJQ34dD59tKcD7bnM1cnIvBO+C88mhOuSf34m5Qc88997i6KKIVeuKJJ8jNzWXq1KlYLM6HBW5ubjz55JPnvU5DkwYNer2el2a/yL8en4lvZHfUGumXJ1qfMzcsQ3oAJZn1X52xqVjLigDwCKq4Im9tXZ0u5IpeS+IeEA5AWe7pZgsatG7Owc6W4vxGy9OnfRd82td9JemLVV27XNZnu5rU5fdtt5rLAwag0myGjfXQrqbtmqLrWX3YrWaKTh8jqPNlLtl/fTnsNnIObOCt114pH6ArRH2cWafh2Wef5fDhwxiNRmJjYzEYDOzbt4+ePXs2OO8mr8VPmDCBl+e+SvaR3wjuOqSpdydEk3HzDiT3+O4Gr9LY2DQ6503YWlKA0S+00vsSHLiW3sMHVGosxVUPhm4KGr1zPZGqzgfhOk3xWzy3laPKfZw1l3y3255p9DLUhauvQ9lHtmO3WgiM6+/SctRV1uFthAb6y1gGcd7OrNNQUFDAwoUL+fjjj9m/fz92e+VZ/eqqyUdNajQa3nj9NbIPbsRuMTX17oRoMgbvQBw2C7a/nvC7mptfKHqvAPJPHXJ1UUQVVGoNBu8ASrNTmm2feQn7QaXCp13nZtuncL3qKuYavRt6rwCCuw7F4FXbZNEXHsVhJ23/BgLj+reK47dbysg5uJE3Xn+tfBExIRpq48aN3HXXXYSFhTF//nyuu+46fv+96oWX66pZploZNWoU3bt1I/PgpubYnRBNQu/hC9CsT45rolKp0Bm9cNhaVl9q8Te/yO7kJf6B3do831FB8mE8Q6LRuXs3y/5Ey6YoCg6rCbXu4hxMa7eUYTeX4h0e5+qi1EnGgU306N6d6667ztVFEa1USkoKc+bMITo6mjvuuAM/Pz+sVisrVqxgzpw59OrV67zyb5agQaVS8fZbb5AV/wvms/pXCtGaWE3FAGiN9V9spcmoQHEotacTLhEY1x/FYSdl5/dNvi9FcVCanYxncETticVFwVpaiM1UgtHXdWtquJKpMAcAraF+Cxu6grkol+zDW3nrzddbRPdX0fpcd911dOnShfj4eObPn8/p06eZP39+o+6j2SZ179evH2NvvZXMPauba5dCNCpLUS6o1M02qLU2dosJU34mei8/VxdFVMPgFUB43+vJPvIbJdnJTbqv7KM7sJYWyqBlUU7n7o3G4E5J5qnaE1+A0vasweAThGdIpKuLUqvMPau5bexY+vXr5+qiiFZq3bp1TJ48mRdeeIFRo0Y1SRe3Zl0J6j+v/Zui039SePpYc+5WiEahOOyoVCocdpuriwJA0tavUBwOAmL7urooogZBnQaiNXqRe2x3k+0j4+AWTv36NYFxA/AMiap9A3FROPPEWuvWMhd3a0o5x3dTmPon4X2uqzRrVEtTePpPitOO8dq/57m6KKIV++WXXygqKqJPnz7079+fd955h6ys6hcYbYhmDRpCQ0N58YVZZO7+AcXR8NHbQriCd7vOKA57sy/YVRVrWRF5Cftp23dUhcW7RMujUqvxj+lN7sm9TXLdy0v4g5Qdqwi5ZDjtB9/a6PmL1uvM+dZSHnQ0l7K8dJK2foV/h0vxjbjE1cWpkeKwk/n7D7z4wixCQ2XWM9FwAwcO5KOPPiItLY0HHniA//73v4SHh+NwOFi/fj1FRec/iUuzBg0A06ZNw9/Ljcz4rc29ayHOi5t3IG5+oRScind1UShIPgyAb2TLviEKp4CY3thMxRSkNG7AabeYOLXta3wjuxPe93rpCy0qKEg+jN1cil9Ud1cXpdk4bBZOblyMwSuAiMG3tvjfRGb8Vvy9jUybNs3VRREXCHd3dyZNmsTWrVs5cOAAjz32GK+88grBwcGMHj36vPJu9qBBr9fzyccfkbl/vQyKFq2OwSsAm6nE1cWg4NQhPILbo2tJg7JFtdwD22L0b0POnzsaNd+0feuxW0y06z+6xVeORPMryUxC5+FTvtDgxSD9j42Yi3KIHnE3am3LnjXKXJRL5r51LFr4sSzkJppEXFwc8+bNIyUlhWXLlp13fs0eNAAMHz6ccbfdRsau71AUmflFtCKKA8Xh2qZ+m6mYwtQ/8Wnf1aXlEPUTGDeA/FPxmItyGiW/wtPHyDiwmTa9r0HvKYPhRWU2Uwk648Uz/a65MIeMg1sI7jy4xS9wqCgK6bu+5fbbb2fYsGGuLo64wGk0Gm666Sa+++6788rHJUEDwJtvvI614DR5CftcVQQh6s0rLJaitBPYzKUuK0Py9lWo1BoCO8osG61JQGwfNHo3Mg5uPu+8SrNTOLHhU7zadCDkkuHnXzhxQVIUR6MPAjYX5ZJ56BcSNn3OsbUfcfi7tzi5cbHLF2912Kyc/HkJWjdPQnte5dKy1KaovYrThfuxF2bw5huvu7o4QtSZ1lU79vf3Z8G77zBh8v14h8e1inmUhfDv0JuUXT+Qd3IfbiMHA+B1qvlay7IObyP3xB4ih94hXZNaGY3OgF9kd4rTTp5XPqaCLI6t/Qg3nyBirrgXldplz35EC2cuzEbn3vDrhM1ciqUo1zkBRNpx8hP+oDQnBZVag3tQO3RuXhj9QslP/IP4lf/BMySCdgNuQuvm2YhHUTvFYefkz0soy0sn7vqH0BqMzbr/uipq7+xCaC8tIed/37Jk4Uf4+UkroWg9XBY0AIwbN44PFn7C7v0r6djvTlcWRYg60Rm98GnbiZzjuwn/K2goaq9qlsDBbjFx6reVBHUejH+HS5t8f6JxKQ4HpTmpqLUN77tsKsjiz9XvoXXzIHbkZDR6t0YsobiQlOWlU5yZSMTgsfXe1lJSQOahLWQd3obDZgFArdXj064zIZcMx6dd5wrnXmnnwWQd+Y38pIPEf/s6nUdPb7ZVye2WMpJ+XUFB8mE6XDUJj8B2zbLf+jgTLJyRt3YVlw+9jNtuu81FJRKiYVwaNKhUKr5Y/BkxHeM43fkP2nhVP8PD2T+65nyyK8S53IPakXV4W4XXagoczpy753vemvIzQHEQ2LGfDHpthWzmUkqzkwmMG9Cg7QtSjpC4eRkag5GO1z7YZE9zz63ggFxzW6PU33/E4BVQrwcMpvxMErf8l5KsJNQ6A8Fdh+AX1QNFcWD0Dal2YLF7YFsiLhuL3sOX03vWUJJ1Ct+IboBzutfj6xbiGRpFm15Xn/dxnX09zUv8g1O/rsBhsxA1fDw+7Tqfd/6N7dzfU9GhfVgTjvHZ/76X67hodVwaNACEhYWx6OOPuPu+yRinRuOXW/uNsLEqYUI0hEbnht1Wuf9ubS0O59sikZF7AJVWi5tvcIPzEK6jM3ri5uscnHluRaK688JuMVGYepSsI79RdPoY3m07/dU1rfrrZF3zrsu2Z7/eGPnINbvp2S0mMuO3UnAqnojLxqLW1O02bzOXcmLjZyh2G5HDxuPTrku9u/n4RHTl9J41KIoDAGtpIcfWfkhZbhooDuhV78Mpd+45dbrkAGkbF+PbvivtBt6M3sO34Zk3UpnAeY5Xd/7biovI+99KPvv4Y1mTQbRKLg8awNlN6b/Lv2LDd1+huWMi3sl12665uoWIi0ttrVoavRuK1Ypit6M6Z5n2c8/JSk+ZqjlnC9sqZHz/FT6XDsDYNqLC/ovaq7DkZJG3bRMBI66lJNpQbdlEy2YqyMBqKcKfv7slKA4Hp0sOkr/jVygowVKSj0dQeyzFec7WJcAYEE7U5XfhF9UDlarqMQznU+Gvbtua8qnLNnUtS13zqm3/tbX21ZTmQpG8/VtyT+whMG5ApdXirWXFpP+xEUtRDmX5GWgNHrj5BGEpyceUn4nDbiVu1EMNnnnI3b8NGr0bprwMMoo2kbJ7NdidC8xFDr29QXlW9T0XHzlI+ldL8IvsTtTwuxo0rqeuwXV9z/Pq0iuKQt7qb7h25EjpliRarRYRNAB89MH7dOrSlYLd2/EOrth8X9OP9uybUE3pakpT35tIXW5A5/O0TzSt6lqqqjo3qqrklDcpVzNd8JltaqrEndn/mX8rFiuFe3ZQuGcH3r364dY2Ap9e/Shqr0FRFPJ3bUOl1eI3YGiNZWtujVGRrK+GHHNduo/VJ+/ajrO6c0ttcHMOgty0Dp9LB5C94UdKThzFXlSIW7tI9O3D8aUrJdnJGP3DcA9og390b7zbda62K0NDKvz13b4haZsij9q2retnUR1X/57Ol6Io5J7ch+Kw027ATeUzJ1nLismM/4XsoztQHHY8AtvhE94Jc3EuZblp6L38neMVul9+XivMF2ckYreYSDu4EcVqLb9GevXog7mLP/rGaPU6tJ+0rxbj2fkSIvuNrxQw1PZ7ru26XFe2okLyd27Fq2sPDKG1r4NRsHs76sw0Ptr8U732I0RL0mKChsDAQL5c9gUjrxuFcUoUgaaQOm97vjeKpryJNcY+mkptwVZdgrEz6Vx1fHWp+FenPt/d2fsp9baBWg2a6qcyrO85qdbrcY/uSOnJPyncu5PCvTvJ+XkNADpff0wpSQQMH4lar6+UR3NVdBrrZttU5WiM7RrrWKrLJ/qx50lb8Tk5P68h5+c1aDy98O7ZF89O3TC2iyxPd+58N8VNWCbxt4Z+Ri0l2FCpVMRccQ8nfvqUhE1L8WnXCZuphLT9G3HYLHiHx/019sCnSfZf6FEEgGKxVHjdkuVsMTvfc7DsVALp3yzFq1tPYvtUDBiqe+gDTXOPyl7/A4X7f8dhNhF83S01pjVnppO//nv+9+OPBAQENGo5hGhOLSZoALjiiit47JHpvL1oMbr7/4VPmqyQ2JQa46ldfdI1heba99n7UVKsqHW6Rh/EFn73/eT+8hM5P69BHxiMe1QsKq0W0+lkwm6bgFfXnrWWTbRsar2e8DsmUXRwH+asdHz7DELrdfEsvnWhqu9vsCmDDJ92nYkcMo6krV+Rn3QAUBHYaQDhl16L1s2jyfZb1F6FF92xXXMjWWtWVXgv8Iprzzv/0pN/cvrLz3ALb0/ITXegOq0u329dytbYSpOcUyfn7/wVz07dcI/uWGU6h9VC7soveGT6dEaMGNHo5RCiObWooAFgzuzZrP/pJ06u+R71DWNcXRwhKlHsNmjkBZMAVGo1AcOuQrFaKdiznaDrbpbZNS5QXt16VmpNEBePph4o7h/TG/+Y3iiKA8VuP69pfuvi7OPxGzgMQ2g4JcePkLd1I+AcAHw+iv+M5/SyT3CPiiVs7ATUWi1F7c8ry/NmLylGHxKGJSMNjUf1v+b89T/SITyM2S++2IylE6JptLigQafTsWL5ci7p0ZOiqBi8uvV0dZGEqEitAYejybJ3ax9J7i8bSFu+mLCxd8viXUJcJBraYlF9EKJGpW349aO2MYDV7dc9qgPuUR3w6taTzB+/oXDPDnx69q0ybW1MqcmkLV+MR8cutLntnkqTT7iKYrdhyUhDFxiMISSsyjSFB/ZSFr+fb/7Yj04nPSdE69figgaAqKgoli39nFvH3Y4+OBRDsExNJloOlVrlbG1oIp4duxB6y3jSv/kCU8pQjO2jmmxfQojWq7G6mDbmPs9mCA3HXlZaYbxOfdiKizj9308whIQSdutdLSZgsObnlj84smZnUpp0EveI6AppzJnp5P74NSuWLycyMtIFpRSi8bXYR5g33HADj/xrGrkrPsdhrjwnvhCukrVmFYqt6YIGAH2gcyIApQlbNIQQoimVJZ7Amp2Jd/fe9d7WYbORsepLFLudsHH3otZVvbCcKxTH/1Hhb7Vej8NqLb8v2E0mclcs4dHp07n++utdUUQhmkSLbGk446U5c/j1t9849P1XBI65S/p3ixZBHxiMraQx5rOpXvGRA6iN7g1+QieEEK5mzjgNajXGiJgKryuKQlnCMWxFhXh1vxSVSoXicFCacIyiP/ZgLy3BnH4aW1EBYbdNQOft65oDqIKtqJDcX37CvUMcpcePovX2pTj+D3K3bACNhujHnifvfyvp2SmOObNnu7q4QjSqFh00aDQavvnqK7r16En+1o34DbnC1UUSAt8BQ8n8cQXWgnx0Pr6NmreiKJQcO0zB3p14xnVtMc3xQghRX27h7cHhIHnRu3h26oZP7/44zCbSVy2nLOEYAF7delGadILM1SuxZGWgDwpBFxCER8cu+A0Ygj6o7tOvNzWH2cTpZZ+AWk3oLXdiyc4kY+UyZ8CAc3rsot3b0WdnsOKntWjk+i0uMC06aAAICgrifz98z6DLLkMXHIpnXFdXF0lc5Lwu6U32xv+Rs2ktoTeOa3A+ptPJFO7dhc4/EIfZhCU3G0vGaczppzFGRBMw/OpGLLUQQjQvY7tIwu++n7xtm8j5eQ25v/yEYrdXGBOW/s0XFB3cizEyhnY3jMWtfVSL61WgKAqm5ESy1n2PJSeTtvdMRevhidbDk6CRo0lbsRS3Nm3x7TuYvNUr2LZ1K0FBQa4uthCNrsUHDQC9e/dmyWefcdc996Kb+FC1MxUI0Rw0bm4EDLuarDXf4j9oeIOfhJUlJZC/c6szTw9PdP6BGELaEDDiWjw6dmlxN04hhKgvjw6d8OjQCVtxEblbN6LWavEfciW2kiJOffgmpYnHCbnpdrx79m2R1zzFbuf08s8oOXIQnX8g4Xfdj87Xj7SvllBy/AgOUxmeXbrjN/hyMj7/iKWLP6NXr16uLrYQTaJVBA0AY8eO5cDBg7z+3gJCJj2M1sPT1UUSFzH3mI6gKNiKChscNPhc2p/sjf8jYOiV+EvXOyHEBUzr6UXwNTeW/603GIj855Oo9XrUeoMLS1azsuQESo4cxBDWlvb3TweHg9QvFmI6fQq/AUPR+fljjIkj89N3eeqJGdx6662uLrIQTabFzp5UlVnPP8+VQ4eQ8/ViHFZL7RsI0UTsfw2E1ng0fIVVtd6AztvHOX2fEEJcZLSeXi0yYLCbTNiKCgEwto/GI64r1twsrLnZnFr4NqUJxwgbcxcBl4/Es2sPcld8zpVDhzDr+eddXHIhmlarChrUajXLPv+cmAA/clZ9KdNRCpcxhLYBlQpTWmq9tnNYrWRv+JHiIwdR7Hbc2kVSuG8X9rLSJiqpEEK4Xv6ubSQvepe87Vuw5uVQsHs7OT+vpejgPlcXrZzDYiFn01pOzH2ahDfnUJaSRN62TZQcO4xis5H0wes4zGbaT/kXHrGdURwOcr79Lx2C/Fn2+ectsnuVEI2p1XRPOsNoNLJ29Wou7deP3PXfEzDyxto3EqKRadyM6AODMSUnla90ao6o3PplSKo4t3jRgT3k/vKT8w+1GhwOjFGxLfJpmxBCNAZLdiaZq7/BENKGrDWryPrft+Xv6QKC8OrW02VlO8NeVkbqFx9jOpUAgGKzkbxwfvkibgrg3aMPwdfchNrgvF7nrv8e79JC1m7aidFodFXRhWg2rS5oAAgMDOTnDRu4tF8/stp74n3NEKByBU2IpuTWNgJTSiJQdcBw9utnzk23Nu0A8Bs0HLWbEUNYOB7RHWVqVVGuLsHn+eZXm7P3V9v2ct0VtSk9eQxUKtrd90+K1cmo/izBfDqF3F9+ImT0WJeWTVEUCvfsIGvDj6A4CBt3L2lffup886+AQefnT+DI0Xh17l6+Xdbhn1COHmTjzp0EBAS4oORCNL9WGTQAREdHs+5//2Po8OFofLzwGNizTjdHucGJxuLWNoLCvTv58/lHce/XnYDJY1HrdVWmLQ8eaIMxKpaCPdudebSLxKNDp0rpzrgQz9eGVGKr4orPpj5lr0/5asv33OCzofk0VnmaYp+uVlWgdCH+/pqboigUJu5HH9EGawfQK20pOPATBb9uxG/gMNwjO7ikXOYIC44yMzmfraR0+z48LrsUz6F9KVi+CQCfPoMwxsSi1mjxiO2MJcqGGed5UbxtL2Xf/cwvmzcTHR3tkvIL4QqtNmgA6Nu3L999+y3Xjx6NymjAvWfnWrepS6WsoTfB+uTVGm5GZ5e9NZQXmq8CY0jSo+kbBt87/y7d+QfG7nF4XnZpjduZIywEPDKO/K/XUrxlF6XHjlBoOYZbXFS16c/sryVzRcWxpVdWm6J8Lf2YW7OqPtumChIvJgUl8ZgOHSdo2t1YM3PIW/o9ZfuP4DN6BIE9r2328pgjLCh2OyVb9pC/Yh0Os5mAKWNxFJeS+dpC1J7uBD9yL8Yefz/MsfD3uhKle+MpWryKH777jj59+jR7+UXDFRUV8eyzz7Jy5UoyMzPp1asXb731Fn379q12m6VLlzJv3jyOHTuGj48P11xzDa+99lp569KhQ4d47rnn2L17N0lJSbzxxhtMnz69Qh6RkZEkJSVVynvq1Km8++67ANx777189tlnFd7v378/27dvP8+jblytOmgAuOqqq/ji888Zf/ddqKfdjVvnmNo3Oktj3oQv5Kdzra28Tc0cYUFnDwGtFmzOG4ol6TTUEjQAaLw9QatBpdfhN+46DLERddpfYwa40PBKjpwLQlR2vt3AWoP6PAQ7k9a0/DgAluQ0st5eAoDv2GvwGTUcVZKmym3qk391qmo5cpgtlPy0m8I1v2DLysW9fw/ce3Um/5v12LJy8RoxAN9br0FtrHqMmenwCQo+XM6ypUu56qqraty/aHkmT57MwYMHWbJkCW3atOHzzz/nyiuvJD4+nvDw8Erpt27dyoQJE3jjjTe44YYbSE1N5cEHH2Ty5MmsXLkSgNLSUqKjoxk7diyPPPJIlfvdtWsXdru9/O+DBw9y1VVXMXZsxa5511xzDYsWLSr/W69vedeHVh80AIwZM4b3Cwt58OF/4v/4JAxRbV1dJHERUGk0BN5/G4U/bsaanoXl1Ok6b1uybQ9eVwzE64qBdd6msSvrUvkXwrUulN9gTcdhz3NOXVqwcgMeQy7FrWMUHn89XKnr8TfkczqzjeJwYDmRQunvByne8juO0jLc+3TD66pB2LLzyPt6LWqjgbAXp6FvV/3CseaTyeS9s5QP3lvALbfcUu/yiKZRWFhY4W+DwYDBUDnoKysrY8WKFaxatYqhQ4cCMGvWLL799lsWLFjAnDlzKm2zfft2IiMjmTZtGgBRUVE88MADzJs3rzxN3759y1sqnnrqqSrLeO7q4K+88goxMTEMGzasUtlDQ0NrO2SXuiCCBoCJEyeSX1DAzOefw/+xiegj2ri6SOIioI8Mx5KUikqvw3LqNIrDgUpd+0zG+rZh2HLym76AQgjhQv733ozniAFoPN3RhQXVvkEjURSFsn2HyV+xDmtKOmoPIx6XXYpbXBT5X6+ldNcBNP4+aAP98L35qhoDBktSKnlvLuaVOS9x7733NtsxXKi8E8xotec3Pa3NZgagXbt2FV5//vnnmTVrVhXpbdjtdtzc3Cq8bjQa2bp1a5X7GDRoEM888wyrV6/m2muvJTMzk6+//ppRo0Y1uNwWi4XPP/+cRx99tNIUvZs2bSI4OBhfX1+GDRvGSy+9RHBwcIP31RQumKAB4JHp07FYLDz/0mz8H5+Evr0EDqJpabw9QaNG7e2JxserTtuYE1NQGfRYU9KbuHRCCOFaaoMetzp0wWxMpqMnyf9qDebjpzB0iibkySkYOkZSuvMA2e//F21ooPO1uKhaH/JYkk6T+59PefH/nuWRc/qqC9dLTk7G29u7/O+qWhkAvLy8GDhwILNnz6Zz586EhISwbNkyduzYQWxsbJXbDBo0iKVLlzJu3DhMJhM2m43Ro0czf/78Bpf322+/JT8/v1Lwee211zJ27FgiIiJISEjg2WefZcSIEezevbvaY3KFCypoAHjyiSewO+zMnjsXvxmTanx6IMT5UrsZMES1RePnQ9BDd9aaviz+OJmvfQIOB4bOMuuGEEI0FkepidxlP1Dyy+/oI8MJfnwSbl1jUUwWcj5ZQcmve/AY2Av/e29Gbai9v7jl1Gny/rOIWU8/wxMzZjTDEYj68vb2rhA01GTJkiVMmjSJ8PBwNBoNvXv3Zvz48ezZs6fK9PHx8UybNo3nnnuOkSNHkpaWxowZM3jwwQdZuHBhg8q7cOFCrr32Wtq0qfhQe9y4ceX/7tatG3369CEiIoIff/yxRXWHu+CCBoCnn5qJw+HgpXnz8Ht8Evp2LbuPmGjd9DHtKd11oE5pi3/6DV14CKFPP4CqDjctIYQQtbNl55H2/HwUux3/e2/Bc1hfsNsp3b6f/JXrsRcWETBlLB6Detdp5WZLchp5//mU/3vqKZ568slmOALR1GJiYti8eTMlJSUUFhYSFhbGuHHjiIqqevbCuXPnMnjwYGb8FTB2794dDw8PhgwZwpw5cwgLq99D6aSkJDZs2MA333xTa9qwsDAiIiI4duxYvfbR1C7IoAHg/55+BkVRnIHDI/fI4GjRJBRFwXzkJLo2IXVKr3I34igpxZ5f1Kz9e4UQ4kJWtHknjpJSABxlJgpW/UTxzzuwFxTh1jWW4McmogsJrFNe5pPJ5L25mGeeeIJnZj7dlMUWLuDh4YGHhwd5eXmsXbu2wsDms5WWlqLVVqwma/5aiFVRlHrvd9GiRQQHB9dpTEROTg7Jycn1Dkya2gUbNAA8+8z/4W505+nnnsVv2t3VzoUvREMpVhvW05ko1tOkPPYKPqOvwGtY9XM+e18zBPOxRNKefxvva4Zg7N4JfUy7Oj35EkIIUTWfa4aiDw8hf9VP5H+5GpVeh8eg3nhdNQh9eN0e6oBzPETe258zd/YcHq1mCk3ROq1duxZFUYiLi+P48ePMmDGDuLg4Jk6cCMDMmTNJTU1l8eLFANxwww1MmTKFBQsWlHdPmj59Ov369SvvXmSxWIiPjy//d2pqKvv27cPT05MOHf5euNDhcLBo0SLuueeeSoFIcXExs2bNYsyYMYSFhZGYmMjTTz9NYGAgN998c3N8NHWmUhoSLrUyH3z4IdMemY7v1Dswduvo6uKIC4w9v4jirb+T//VafG6+Ct8br6gxvcNsIX/5/yj+dTeKyULwjPswdq16IJYQQoi6UxwOHKVlqHS6Oo1bOFvZgT/JX7CM+W++xf1TpjRRCS9ehYWF+Pj4MPSy59Bq3WrfoAY2m4ktW1+koKCgzmMali9fzsyZM0lJScHf358xY8bw0ksv4ePjAzgXWEtMTGTTpk3l28yfP5/333+fhIQEfH19GTFiBK+++mr5ug6JiYlVdm8aNmxYhXzWrVvHyJEjOXr0KB07VqyHlpWVcdNNN7F3717y8/MJCwvj8ssvZ/bs2ZVmh3K1iyJoAFi2bBn3TJqE7+Rbce/TzdXFEReYvBVrKd64nfA3nkat19VpG8Vq49T9z+J/z014De/fxCUUQghRndLfD5L/8dd89skn3HHHHa4uzgXJ1UGDOH8XdPeks91xxx14eHhw2+23Yy8sxmvEAFcXSVxANB7uKFYbKl3dflKKw0H2+8tAUdAG+jVx6YRo2ez5RdgLCtG1byNd9USzK9q4nZKv1rJi+XJuuOEGVxdHiBbrogkaAEaPHs1P69dzzajrKCwowuumK+UGJRqFPqINisWK6Y+jGHt0qjW9Jek0pbsPETB5rHSZExelyLZZAFhyS/jt3vcBULvpiLp3MLYB17iyaOIioSgKRSs3YNvyOxs3bGDQoEGuLpIQLdpFFTQADB48mJ2/bWf4lVdSmF+M94TRqP4aDS9EQxk6RWOIiyLv6zV0utYflbpiMJqYUnGmJMVmB0Avs3qJi8iZQOFsyV//DoDOx4i1oAxzbjExZ6U797cjRGNQ7HYKF6/CcDyFX7f9RufOnV1dJCFavIsuaADo3Lkze3buZPhVV5L1zlK87x+H2thyVtwTrY9KpcJv7DWkz1lA5k+HCbmqS4X3I9tmlVd+HGVm8r/8EbQaNF4elfKqqmJVF+dTuYps4ZW0hn4mZ2uu42qMsl5MNG46UKvQB3oSMKgDURMvq/B+c36eLfHcF43PUWai4IPlhFgUft6xs3xQqxCiZhfNQOiqFBQUMPqWm9lz4hg+0+5CG+Dr6iK1Ko11M2/sG/WZcp1vvjUdX1V5R7bNIv7F78nZcZIuz48moF/VU/we+iaN7Pe+IOTJKbh1jql1X+eym61oDHUbbF1dWeuyz4ZudyFKTAm6KI9bNJwEIHXTFL+rmj57W3YeBfM/59IOcaxa8U35zDmi6clA6Nbvog4aAKxWK//450N88dVyfB6+C0N03aa3kgqEqIrdbCX+he8pOJDC4G//iUqjrpSm+GQWe/7xOe1u70fUxMF1yrfoaDonP/4Fva87mZuPEnx5J6LvH4ohwLOxD0EI0UpIMFu9qgIH84lTFLyzlDtvG8eCd96tNF++aFoSNLR+F/0vRqfT8dH7H9ClU2dmPvM0MY9dSfDwOFcXS7RSGoOO8DG9yd2ZQPbWYwQNq3wueUYHEXH3QBIXb8MQ5EXYqEtqHJDvsNjY89DSCq9l/nQYvb87MQ8Mb+xDEEK0EhIwVO/sLqEAJTv/oHDRN7zy8lym/+tfMgmKEA1w0QcN4OyP/ugjjxDboQO33jGOspPZtL9nYJVPiYWojV+vCIKGx5GyYjfZv52g3W198Yz+++ZlKzZRcDAFHArH3lyPNa+EiLsHVp/hWYOq1W5aHCYbHjFBhFzVtSkPQwghWrXItlkodgcJi36l6PsDrFj+Fddff72riyVEq3XRd086V3x8PNeNHoUpUEuHmSPRep5fE5q4OCl2B6dX/8Hxt34CwLNjCHFPXItnZADb7/oIc3ohKp0GxWrHKy6U7v8ei9a9+tVLrQWlFB5Ox+eScGzFZtxCpDlWCCFqYi0yceKVdRhzbKz+7keZIcnFpHtS6yeP0s/RpUsX9v2+l+7+MRx8eDklSTmuLpJohVQaNW1G9UClc07nW/xnBrsnf8quyZ9iTi8EoOvzo+nx+jhKk3M5MHMF1sKyavPT+bgTMCAarYdBAgYhhKhFSWI2h6Ytp7t/NPt+3ysBgxCNQIKGKvj6+rLmh9VMnTCFA9O+JHPTUVcXSbRCKrWKof+bTreXbyl/rTTx7yDUEOSJb/e2dJ93KyWJOWwf/xHH3/uZ4hOZKHaHK4oshBCtXubPRzjwr+VMnTCFNT+slhmShGgk0j2pFqtWreKue+7G7/JYIu8fglovw0BE/dmKzaSs+J1T/92FYnUu7Nbjtdvw7BCM1tOAJbeE1FX7OLV0e/k24bf0psPUy11VZCGEaFUcFhuJH/xC3qZjfP7ZEm688UZXF0mcRbontX5SA67FjTfeyP49+7jx1puJf2wFHZ65BrdQeWoh6kfraSDynsG0HdOH9LUHSf12L/sfXw6AV1woQZfHofd3r7BNWWq+C0oqhBCtT1laASdeWkOY0Z/Ne/cTFVX1OjlCiIaToKEOoqOj2bVtB/96ZDpLpn5O9GNXEji4g6uLJVohraeBtmMuJfzm3hSfyKQkIZusTUdJWLgVxWrHu1s4gYNiUGnUBAyIdnVxhRCixcv+9TgnXlvPhDvv5u0338JgMLi6SEJckCRoqCM3Nzc+WPA+lw8bzn33T6ZwTzIR919Wr5V5hThDpVbhFRuCV2wIoVd3RbE7sJdZ0HgYZP5wIYSoA7vZStKHW8nZeJRFHy7k9ttvd3WRhLigyUDoerr99ts5uP8AgWkqDj28nOKTsriOOH8qjRqtp9sFETCMCPmzwn9CCNHYik9mcfCfXxKUrubg/gMSMAjRDCRoaICoqCh2bP2Nf4y/jwP/+pLUb/ci48mFK7WECnp1ZZAAQgjRWBRFIfXbvRz415dMvXMyO7b+JuMXhGgm0j2pgXQ6HS/NnsPIq67mtvG3c2RHEtGPXoEhyMvVRRMXgIZWsM/dbmNGx8YoTr33W1O6upSpqvya61has5b2uZ0pT01laMi5LufCxcmcVcTJ139CdbqU9WvWMXToUFcXSYiLiky52ggKCgqY+vBDfLNqJZFThxF8RecLopuJaD4X21P46ip9jfE5uKpCebF9hy1JVd95bd+HBB6th6IoZG44TOJ7mxlz8y28+/Y7svZCKyRTrrZ+EjQ0olWrVjFx8iTcugQTOe1y9L7utW/kQnV96luf/BpLY5WrtiedTVnR25jRUSqSNTj7O5HPSVyoWmpwcr6/ueY6LkteKYlv/4z5SBaLPv6E0aNHN8t+ReOToKH1k6ChkWVnZ3PfA1PY8PNPRPxjCEGXd2q0VgepWAkhhGjJGiuYUBSFrJ+PkLTgF668/AoWfvARgYGBjZK3cA0JGlo/GdPQyAIDA/n2629YsWIF9//jAWxb93DDs93wCTW6umhCCCFEk6rp4VZdA4pLlf18P/sg2fEWPvvwE8aMGSNdfoVoASRoaAIqlYpbb72Vyy+/nH898jDv37KSEdNjufTWCNRqufAJIYS4+NTWWu5wKOz+Ool33zjK2Ftu5a2V8/H392+m0gkhaiNBQxMKCAjg88VfsGbNGu67fyLxqzO59v86E9xBZlgSQgghzsg4Vsial45gylSx6uvvGDlypKuLJIQ4h6zT0AyuueYajsYf46bh4/nojl/Z8OYRLKU2VxdLXMAyTxRxYlsWNovd1UURQohqWUptbHjjCAvv/I2bho/naPwxCRiEaKEkaGgmnp6e/Oe119m143dMh7z54JZtHN2U7upiiQuI3eYgeX8eP8w5wLs3bWLxA9v5YNwvri6WEEJU6eimdN6/+VdM8T7s3L6L/7z2Op6enq4ulhCiGtI9qZl1796d7b/uZOHChTz+xGPs7ZHGVU90JKC9h6uLJlqp1EP5nNyezYY3D1d6L7ybb/MXSAghapCTVMz6fx8jdX8B//n360yaNAm1Wp5hCtHSSdDgAmq1milTpnDzzTfz9P/N5IMxi+l3ZyRD7o/B4C5fiai7nV8m8uOcA+V/G310DLonhrICC12vbkPb7n4uLJ0QLY/d6sBmdbTYa+3pQ/mYS2xE9bvwphc1l9r45cMT7FyayIR7JjB3xSsEBAS4ulhCiDpqmVfNi0RgYCAfvv8R/3hgKv946AEWjN7KiOkduGRUuEwvJ2qkKAp7Viaz+uUDdLkqDKO3jt0rTnHbf/oQ3f/Cq2yIxlWSZyZxVw6n4wvIPVVCcY4ZnZsG72A34oaH0GlEaIu6BplLbJiLrXiHNHzq6uzEYvZ8c4p9q5IpK7Qy4M4orn6sS5XHqSgKmceLSIsvwFRsJS2+AL92HgyZ3AGNtnGeiCuKgsOm4HAolOSaKcmx8MePKWz/PAGA5/aMQqO7MJ6+OxwKB1ensvHN48TFdOa3bdvp1auXq4slhKgnCRpagF69evHbrzv44osveOyJR/l9WSpXPh5L+14y1Zyo2pYPj7HxnaP0vLEd/W6P5L//2gXAd7P286/VI1pUhU+4nsOuUJpvoSTHzMZ3jnB0UwaKAt4hbgRGeeLbxojN7CB+fRp7v01mzCu96D6qrcvKa7c5MBVZUalVHNuSwTdP70Olgmd3j0KtVWEptaMzalCrVVhKbWgNGtSav895u9WBApgKrOz/MYV93yaTebwIo7eOHqPbojdq2fLRMVRqFdH9A1EcCsU5FgrSy8hNLiHp9xzyT5cBoNGqCIj0ZN93KZzYlklgtBf9bo8kNM4blUqF3eYg9UA+RzdlkLAzm5I8M35tPfAJNVKUZcJUaEXvrsG/nQc2i4PshGLyUkooK7BS3dKqQdGexK9Pw25zYDM7MJfaUBwKil1BAXRuGtx9dMQODcHdR9/0X8h5OLU3lw2vHcOUDfP/s4Dx48fL9UmIVkpWhG5hSktL+fdr/+bVea8QOziYy6d3wL+djHcQfyvMKOM/V26g+/XhXPdUN94bsxm1VkX/8VF0uTIM3zburi6icLG8lFKSduegUqtIO1zAgf+lUpxtBsDLW8WjT3ox4io3QkI1Fbb74HAXPr5rK0VZZtr19COidwCDJ8Zg9NGRebyIP35IxVJqo10PP/wjPPALd8fNS4dKo2rQGjRlhVb++DGF41szSTtSiFanxm51UJRtQnH8nU6tVeGwKahUoNGrsZkdqDUqNDoVVpMDjVZFcKw3pXlm7DaF0jwL4Hyar1ar6HxVGF2vCiN2aAg6g/OYf15wlE3vVVw3wMNfj28bd9r19KPDoCAi+waic3Om3770JIfWniYnqYSSXAt6dw0GTx1lBRZsZgfufnpiLwvGM9BAbnIJRZkmPAPdcPfTYy6ykptcikanJijaE7+27ngEGNDq1KjU4O5nwDPQgMOusHjKb5iK/p5dT6UCvYcWlUqF+q+vqzTfCkDskGDueq9/vT/35pCbXMLPbx7n2K+ZPPXEUzz++Azc3eXadDGTFaFbPwkaWqjTp0/z9P/NZNmyZfQdF8Flk2Nw923ZT5RE0xvn8ztmk8LN12VRkK/w9CxvHvtnPrfebuSpZ73x9FLzZUEfVxdTNIKyAgvpRwsJivbCM9BQbbqC9DIyjxWSnVRCxp+FJO7KIS+ltPz9wCA1117vRv9BBjy9VPTopcdorL6CbzYpfP5ZCYcPWdm4zozJpGAwqCgrU/D1U+HlpSb5VMWpfN399Dyy9gr0xpobr+02B8e3ZnLk5wxS/sgj60QRKo2KyD4BtL3ED7vNgVavxjvEyNVtk7FYFLr30NEuQsvBP6wc2G9hR2E7PP0NmIqtzsq6r56ibBOZx4vwCTGi1qnwCnRzBgwaFZ1GhOLhV/XnV5pvwWpyHouHvx6tXlNlurPZLHYSduaQdaIIU7ENo4+Otpf4Ed7Nt0JrR0PZrA7KCizo3bVodGrUVQRk+79P5pun9wFw/bOXEBrnQ2CkB8a/Wh1MRVZMxVZ8w5q/kl6SZ2brxyf5fXkS48eP5+U5cwkLC2v2coiWR4KG1k+ChhZu//79PP7Eo2z7bRsDJ0Yx4M5I9C10AJ9ofON8fi//95dLS/nmq1LiOum48RY37rw1l3+/5cvheCuLF5YA0KOXjpHXGRk73r1CxVACidZhnM/v2O0K331TxkuzCikuUtDp4b4HPDmha0f7nv607+VP+tFCti89+VfFuxgArUFNbAc1vfvqGTjYQL8BetQaMBpVaBpYmc3NsbNxvZmiIgexcTr69dejN6goKnRwKsnO6VQ7Wzaa+Oq/ZWz8LZhfPKt+6n2TYRe//Wrm3y8VceK4jQ6xWi7tp+eSHjqGDDdUavFoqIvlPHfYFfZ/l8y+71JI2pNT3irjHeKG3l1L7qkSFIfCiIc70eWqMLyD3Zr8vmEptbF9aSK/LUpg8KDB/PvV/9CjR48m3adoXSRoaP0kaGglNm3axKMzHiHh1AkueyCS3je3v2AGyV1Mzg4C6qt353RKSxT8A9SYzQ5KimHemz6MvsWdjHQ7G9eb2PKzma2bzYSEanjrfT+6XqKrMq/WXLk6n8+wPkpKHKz+zkRmhp2AADUBQRquHOl8Yr13t5Wd2y3YrAphbTRYLAqJCXaKCh3YrKDTQ0CgmiRjBO16+tOux9+zWJ0pv6Io5GQ7SDtt5499Vvb8biHhhI3UZDsFBQqjbzYy8X4PPv24hA1rTTjsUFam0K69htQUO+FtNQwYrGfwUAPde+oJCVU3ODioL4tZ4ad1Jn7eYGbLJhOBgRpWrglEp1Nx8A8rSxaVkHLKRmmpQnaWg6xMZ622Tz89Tz7rxSU9mr7V9HzO8fqcY67+LVlNdrITislOKCbjWCFWk52ASE/SDxewe8Wp8nQRl/ozcdGgRh9PYLc62LPyFL9+kEhk+xhe//cbDB8+vFH3IS4MEjS0fhI0tCKKorBy5UqemPk4haY8hjwYySXXtW2UJvGLXXNVRBvKbFKY9mAemzea+WJFAMuWlPD9tyZUKrjvAQ+mPeaF3uA8D5ISbTwyNY8/j9i4eayRGU974+3TsgLMmipaDf0uFEXBagW93vk5/HnUyoH9VuI66ejW3Rk8WcwK2dkOHA6FwCANBgNkZzlISrCTdtpOaoqNQwes/HnURsopO3Y7uHuoMJsU7HZnpVelhl3bLXh6qdBooCBfQa2G9pEafP3UaLUqrBaFrCwH2Vl2bFbo3ktHZoaD0hIHDjs4FLBYFCzOYQZotdCth47YjlrC2mj+CgR0FSp4DofC6u9N7N1toXMXHaNuNNbYzagpzXwsn5VfOQcK3zPZg/se8CAp0caKL8v49usyIqI09Oipx91DRUCgmjZtNXSM09Ktu65ZB8HWtULfWL//hu6vqQKP/LRS8lJKWfbwLswlNkI6ehMU44lKpUJn1KA4FPTuWryD3Wjfy5+2PfzqPDbFYVf448cUtn6QiLebH/PmvsbNN98sg5xFtSRoaP0kaGiFbDYbS5cu5dlZz2DVlDHkwSi6XB3WoIGIrZErKvi/bTWz5kcTsXFavL1VeHqqCQxWY7eDn7+ayKima/r/5MNi5s0pKv/7yWe9mDjFk59/MvHo1DzKyiAqRsObC/yI6/RX5diisPSzEt57q5joDloW/zcAg9uFe37s3G7m2ScLSEqw4+GpokdPHdu2Wsrfbx+hoaxMKX/iDaDTOb+7zIy/X/PxVdGpi45OXXREx2gZMFhPm3ANKhVs3mhmwdvF+AeouXGMkWtGuaFWw7ZfLHTqoiUgsHIXG7td4YvFpezdbSG8rQZPTzVaLaByBjehYRrahGuIiHK+1xocibdy0zXZAAwbYWDcne7Mf72Iw4dstAnXcO8UD8ZPcEerdf35VpfKeEt6YNDYwcPZx/bNV6WsW23CVKbgUKC0REGjcf4/NdVOaYmzKrBmc1CF69m5ZXI4FOLXpbH1/QS0diNzXniZ8ePHo9VKt1lRMwkaWj8JGloxi8XCokWLeP7F5/D0KuDh6UauvtYNtVrl8ibz+mhJN+2q7NphZvo/8snJdlSbxsdHxa97Q6qsKOXm2PnwvRL6D9TT81I93t7OQaXVVRIVRanwtO6hybn8tM5cIY1GA5266Dh0wFr+WqcuWr5dE1Qh3YH9Fu66NYeYWC0ffOpPUHDj9B13tbIyheN/OlsSftlkZvNGM7376Ln5NiN/HrGxa7sZvUHFG+/6sWGtidOpdoxGFeHtNAQHOz/3//1gwsNDRd8BeqJjtISFa/DwaB0Vd1ea/VwBSz8trfBa7z46pj3uRf+B+hb5pLm662FLv/acrbZrekOO5VSijf/9YGLt6jLiDzpnbPrq+4Aqu485HArrVpuY/1YZxUU+vPj8bCZOnIhOV3UXSCHOJUFD6ydBwwXAZDKxcOFCXn55Fl7eJh58WMvIUW4XRctD54i08rnOYzpoufYGNw4fsvDTOgsdO2np2ElLcZGzf/jlV7lRUuwgIFCN2QybNpjw8VVzaT99eZeWqiz8oJh/v1RU7fsAwSFqrhnlRna2g+QkO6kpdkqKHUREaTl53IbNVnmbDrFaOnTU4mZU4eunJjRMze6dFn7dYqHvAL1zqkW9iphYLSqVs9tRSYnCzm1mHArYbWD9O2bgnvs8mPl85YvnoQNW7rsrB51OxZrNQa26YqwoCk89WsCqFX/Noa9xDv4efbORW8a51/g9isZhsyls+dlM2mk7np4qOnbS0bmrVBxburIyhb2/WziVZCPttJ0tP5s5fMiG0ahi+BUGhl9h4Iqr3fD0qnh9cDgU1v5oYsHbNoqL3HjmmReYNGkSbm7nV+kTFx9XBw1FRUU8++yzrFy5kszMTHr16sVbb71F3759q91m6dKlzJs3j2PHjuHj48M111zDa6+9Vr6S+aFDh3juuefYvXs3SUlJvPHGG0yfPr1CHrNmzeKFF16o8FpISAjp6enlfyuKwgsvvMCHH35IXl4e/fv3591336Vr1651/ESah7QnXgDc3Nx46KGHmDx5cnnw8O6bxUyZqmHUjcYW0U2gNkWFDgoLHVitkJfr4M8jNvbvtRB/0Fkr9vJWU1ToQFEgKFiDTgemMgWDG5ic9UdOHLfxzhvFBIc4b3p/HrHx5xFnbX3/HguP/jMfAL0BbFZw/NVw4B+g5oqrDeTnK/j7q+ncRYt/oJrgEA1R0VrG3elOUJCaJ6YXVChzRJSG9hFa9HpISbbz6xYL/oFqYuO0jLjKgNFdReJJO70u1dO9p44evXScOG6jIF9Br4ed2y2kn7aTka6wb4+D9NN2LBZni4HD4eznXljgYOmnJRQWOiMjgwF6XqqnU2cdQcFq/ALUhIRqiO2orXIGmox0Oyu+LMVgUJGb46CsVMGjFS/7cfSwrTxgmP6EF/dO9sDtAu521RJptSpGXCUVxpbMblc4Em9j+69mftlk5lSSnYx05xgdtRoCAtT0G6TngYc8GTrCgLt75QcJNpvCj6vK+Og9O2WlRp555mUmTZqEwVD99L9CtGSTJ0/m4MGDLFmyhDZt2vD5559z5ZVXEh8fT3h4eKX0W7duZcKECbzxxhvccMMNpKam8uCDDzJ58mRWrlwJONfWio6OZuzYsTzyyCPV7rtr165s2LCh/G+NpuL9et68ebz++ut8+umndOzYkTlz5nDVVVdx9OhRvLy8GukTOH/S0nABMpvNLF68mLlzZ2GzF3DfgxpuGetePlDW1bb8bOL+e/IA58DSvFwHJ0/YKqyOqlZDXGfnwEmtVkVhoQOfvwbzZmXasdnAYFDhH6jGYlYoLlIoKHBgtSq4u6uxWRWKihQsFmcFvX2klsFDDXh7q8u7qlzaT4/ZpLDqmzJ++dmMl4+KkiKFE8dt5QFFXSz5yp++/c//RqooCmYzlSrBVqvCyRM23I0qQsI01T5Nz0i3k3DSRuIJG6mpdhJO2NjxmwW7HQYM0nP7Xe4MvbzlVvYcDufxWy0KhYUOsrMcJJ+yc/K4jaOHrRw+ZCXttPOLaR+h4b/fBuAfcGF0txLifGVl2lnzo4mN60zs32eltETBYICBlxmI66ylTbiWS/vqiIrR1jjLlsWs8M1XpSx8345W48PMmbOYMGGCBAvivDVFS0NycnKFlgaDwVDluVpWVoaXlxerVq1i1KhR5a/37NmT66+/njlz5lTa5rXXXmPBggWcOHGi/LX58+czb948kpOTK6WPjIxk+vTpVbY0fPvtt+zbt6/KY1EUhTZt2jB9+nSefPJJwFmPCwkJ4dVXX+WBBx6o8bNoTtLScAEyGAxMmTKFiRMn8uWXX/LSS8+x4K00JtznfGp+bvNzYziVaOO1uUUMv9JAXCcdpjJnxc9UphDXRUdomJqiQoXEBBtLFv3dH/r3nRbuutedSQ96EBLqbEHw9VXTPlLbbDPDdO6q46ln//7bZlMoKlRIO20nKcGGQwFvb2eFHaCkWCEv1+F8cl+m0LN340wfqVKpqKrFX6dT0SFWy8P357FxvZmgIDW+fmqCQtTodM7PKDfHwYH9zlYZjQZC22ho207D3RM9GHene6PNg98YbDZnYLZnl4XNG82kp9nJyXaQl+uoshtXULCz9WbUaCNdL9E5g78WNhuUEK6Sk21n3pwifvyuDJUKBgw2MHWaJz166+jeQ1/nCRCKixz89/NSFi+0ERAQxtyXX2TcuHEywFk0Ov3RVLTq87tvqh3OiS7atWtX4fXnn3+eWbNmVUpvs9mw2+2VutUZjUa2bt1a5T4GDRrEM888w+rVq7n22mvJzMzk66+/rhB01NWxY8do06YNBoOB/v378/LLLxMdHQ1AQkIC6enpXH311eXpDQYDw4YNY9u2bRI0iOah1Wq58847ueOOO/jhhx+YO3cWH7xziHF3Gbh7opHgkPOrSJ44ZmXtahNhbTScPGFj3f9MrPufqdbtPDxVXHG1gaIihSf/z7vatQRcRatV4eevws9fTZduLaNsdhtsXO8cDJ2V5SAry8GxP53vXXG1gfC2Gm6/251evfW0bV99a0RDnThm5dctFg78YSEn24Gvr5rwthq699IzZLih2i5CdrvCV1+UcuiAldOn7WSkOUg+ZcNsdrYm9Rugp2dvPYFBavwD1Hh5q9HpwMtLTUCgmrbtNE0S5ApxoXjsn/ls32ah16U63lvoj59//X4vmRl2liwq48vPzXTp2pWPPprF9ddfj1otvzvR8lXV0lAVLy8vBg4cyOzZs+ncuTMhISEsW7aMHTt2EBsbW+U2gwYNYunSpYwbNw6TyYTNZmP06NHMnz+/XmXs378/ixcvpmPHjmRkZDBnzhwGDRrEoUOHCAgIKB/bEBISUmG7kJAQkpKS6rWvpiZBw0VArVYzevRobrjhBrZu3corr7zAVZdt4oabPLh3ip4OHRtWMU5KsPP2f4rrtU1omLMi7mZ0Tlv61bJSNm909ssPCdPQtZuu3je9i4HeoOIf0zxZ8HbFz3vYCAPvfuxf47YOh7Pf19kD47Oz7GzdbCYmtuq5861Whdwc56JcP60z8eG7xdjtcGlfHf6BGrKzHezeZeGjBSV4e6tYuSaQ8LaVLydH4m3Meqaw0utBQc7Wg6gOWqY96oWPr3znQlQnL8/BhjUmft9pIeWUDZvd2Y0oL9dBepqzy97e3VZKSx11vn4e/9PKog8t/LCqlCuuGMaPP87isssua8rDEKLReXt713kg9JIlS5g0aRLh4eFoNBp69+7N+PHj2bNnT5Xp4+PjmTZtGs899xwjR44kLS2NGTNm8OCDD7Jw4cI6l/Haa68t//cll1zCwIEDiYmJ4bPPPuPRRx8tf+/c+/C5Mym2BBI0XERUKhVDhgxhyJANxMfH85//vMKYUcsYMNide+7TMvCy+k2XOHiogTfe9WXXTgu/bjGTlGAvf2/ygx5cf5ORhBPOxbI+fr8EgPQ0B+lp5uqyBGDh5/4MHir9Z8814T4PTqfYiT9o5difzn48W34207tTOh4eKvwC1AQEqAkJ0xAcoiYwUIOvn4qvlpWxa4el2nzfX+TH8Cv+brJVFIVh/TLJzak4sOO6GwyYTJCfZ6ewQIG/xqAUFip8ubSEcXd6cPK4jYICBcXhzOfMAO5zOVtLLGzbamHpp6UsXu5PvwHynQtxtrkvFvLZx85rp1rtnCQhJlaLTqfCYHC2hnr7qPDycq4VU1XgfjZFUdj2i4XFn9jYsa2MO++8g717n6Jz587NcThCuFRMTAybN2+mpKSEwsJCwsLCGDduHFFRUVWmnzt3LoMHD2bGjBkAdO/eHQ8PD4YMGcKcOXMICwtrUDk8PDy45JJLOHbsGAChoaEApKenV8gzMzOzUuuDq0nQcJHq0qULCxcu5uWX/80777zN4w/PJzjEwt2TNFx/o7FCP9izo12HQ8FicT7lsloVLr/KjWtvMFa7n05ddIwc5VYeNNREo4H+A/V07d4yugS1NH5+al590xdwjiE5/qeNrCwHZWUOiosU8vIc5GQ5SEywsfM3BznZdsw1x2doNPDGvKLyhetysh1kZzowugM5zjQqFSgKrP6++swWvl/Kh++WVnpdo3EumObj4+x65O6hwmhU4e6uwt1DhalMQadX0amzfOdCnGvv738H+zGxWrp00xHWRlPele/Svno8PFS4GVU1dkk0mxR+WFXGkk/sZGWpmfqPh/nyi4dbXIVEiObg4eGBh4cHeXl5rF27lnnz5lWZrrS0tNKYnjOzHp3PHEJms5nDhw8zZMgQAKKioggNDWX9+vX06tULcK7DtXnzZl599dUG76cpyOxJAnDOLLBkyRLefPMVMjLTUBQzFotzAbKcbAf6v55qZWdWrogajSp8fFV4easxGFTo9c6ZjdyMKozuKtzcVGi1KjQaKCpykHDCTmqKjfy8+p96nbtq+fqHwBpn/xBOZ2ZjMpsUTCaF0lKF4iIHhQXOACM3xzn4ODPDzuFDztaL0FBnK0VAoIbAvwZc+/io8fVzfv96g6o8iLBancHjmXx9fdV076XD10+NWu1MZzBUbnIVQtTdyeM2Vn9fxvE/bSSfspOeZic/z4H974ZdVCoIDFRjsSqYTQqKAp5eatyMoFapKCnRERIcxvTpT3H33XdjNFb/oEeIpnJm9qQrg+4774HQNoeFDVkL67VOw9q1a1EUhbi4OI4fP86MGTMwGAxs3boVnU7HzJkzSU1NZfHixQB8+umnTJkyhbfffru8e9L06dNRq9Xs2LEDcFbu4+PjAbjuuuu48847ufPOO/H09KRDhw4APP7449xwww20b9+ezMxM5syZw+bNmzlw4AAREREAvPrqq8ydO5dFixYRGxvLyy+/zKZNm2TKVdGyKYrChg0beO0/L7Pxpy20j9Tj7q5GUc5MX6pGr3dWCNVqFTabs8XBZlOwWRUcioLD7myRcDicA2EdDlAcCg6F8mlV1Wrn9mVldqyWup+CgUE6ul7i3kRHL4QQrUNRkR2LWcHhcF5/zWYFtZryByplZQpJCWZGjBjC448/w5VXXikBvHApVwcNy5cvZ+bMmaSkpODv78+YMWN46aWX8PHxAeDee+8lMTGRTZs2lW8zf/583n//fRISEvD19WXEiBG8+uqr5es6JCYmVtm9adiwYeX53H777WzZsoXs7GyCgoIYMGAAs2fPpkuXLuXpzyzu9sEHH1RY3K1bt24N/ISahgQNolrHjx9n9erVOOqzaIEQQgiXU6vVXHfddeVPO4VwNVcHDeL8yZgGUa0OHTowbdo0VxdDCCGEEEK4mMxzKIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghhKiRBA1CCCGEEEKIGknQIIQQQgghLmhFRUVMnz6diIgIjEYjgwYNYteuXTVus3TpUnr06IG7uzthYWFMnDiRnJyc8vcPHTrEmDFjiIyMRKVS8eabb1bKY+7cufTt2xcvLy+Cg4O56aabOHr0aIU09957LyqVqsJ/AwYMaJTjbkwSNAghhBBCiAva5MmTWb9+PUuWLOHAgQNcffXVXHnllaSmplaZfuvWrUyYMIH77ruPQ4cO8dVXX7Fr1y4mT55cnqa0tJTo6GheeeUVQkNDq8xn8+bNPPTQQ2zfvp3169djs9m4+uqrKSkpqZDummuuIS0trfy/1atXN97BNxKtqwsghBBCCCFEUykrK2PFihWsWrWKoUOHAjBr1iy+/fZbFixYwJw5cypts337diIjI5k2bRoAUVFRPPDAA8ybN688Td++fenbty8ATz31VJX7XrNmTYW/Fy1aRHBwMLt37y4vC4DBYKg28GgppKVBCCGEEEI0C5tiweY4z/8UCwCFhYUV/jObzVXv02bDbrfj5uZW4XWj0cjWrVur3GbQoEGkpKSwevVqFEUhIyODr7/+mlGjRp3X8RcUFADg7+9f4fVNmzYRHBxMx44dmTJlCpmZmee1n6agUhRFcXUhhBBCCCHEhctkMhEVFUV6enqj5Ofp6UlxcXGF155//nlmzZpVZfpBgwah1+v54osvCAkJYdmyZUyYMIHY2NhKYwzO+Prrr5k4cSImkwmbzcbo0aP5+uuv0el0ldJGRkYyffp0pk+fXm2ZFUXhxhtvJC8vj19++aX89S+//BJPT08iIiJISEjg2WefxWazsXv3bgwGQ+0fRjOR7klCCCGEEKJJubm5kZCQgMViaZT8FEVBpVJVeK2mCvaSJUuYNGkS4eHhaDQaevfuzfjx49mzZ0+V6ePj45k2bRrPPfccI0eOJC0tjRkzZvDggw+ycOHCBpX5n//8J3/88Uel1o1x48aV/7tbt2706dOHiIgIfvzxR2655ZYG7aspSNAghBBCCCGanJubW6UuQs0lJiaGzZs3U1JSQmFhIWFhYYwbN46oqKgq08+dO5fBgwczY8YMALp3746HhwdDhgxhzpw5hIWF1Wv/Dz/8MN999x1btmyhbdu2NaYNCwsjIiKCY8eO1WsfTU3GNAghhBBCiIuCh4cHYWFh5OXlsXbtWm688cYq05WWlqJWV6wmazQawNnKUVeKovDPf/6Tb775ho0bN1YbpJwtJyeH5OTkegcmTU1aGoQQQgghxAVt7dq1KIpCXFwcx48fZ8aMGcTFxTFx4kQAZs6cSWpqKosXLwbghhtuYMqUKSxYsKC8e9L06dPp168fbdq0AcBisRAfH1/+79TUVPbt24enpycdOnQA4KGHHuKLL75g1apVeHl5lY/p8PHxwWg0UlxczKxZsxgzZgxhYWEkJiby9NNPExgYyM0339zcH1ONZCC0EEIIIYS4oC1fvpyZM2eSkpKCv78/Y8aM4aWXXsLHxwdwLrCWmJjIpk2byreZP38+77//PgkJCfj6+jJixAheffVVwsPDAUhMTKyy5WDYsGHl+Zw77uKMRYsWce+991JWVsZNN93E3r17yc/PJywsjMsvv5zZs2fTrl27xv0QzpMEDUIIIYQQQogayZgGIYQQQgghRI0kaBBCCCGEEELUSIIGIYQQQgghRI0kaBBCCCGEEELUSIIGIYQQQgghRI0kaBBCCCGEEELUSIIGIYQQQgghRI0kaBBCCCGEEELUSIIGIYQQQgghRI0kaBBCCCGEEELUSIIGIYQQQgghRI3+H2AoUQAD7/srAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "%matplotlib inline\n", "import matplotlib.pyplot as plt\n", "import cartopy.crs as ccrs\n", "\n", "fig = plt.figure(figsize=(10, 5))\n", "ax = fig.add_subplot(1, 1, 1, projection=ccrs.EqualEarth())\n", "\n", "cf = ax.contourf(\n", " np.rad2deg(lons),\n", " np.rad2deg(lats),\n", " acc.reshape((500, 500)),\n", " transform=ccrs.PlateCarree(),\n", " cmap=\"viridis\",\n", ")\n", "ax.coastlines()\n", "ax.set_global()\n", "plt.colorbar(\n", " cf, ax=ax, orientation=\"vertical\", label=r\"Acceleration $\\left(m/s^2\\right)$\"\n", ");" ] }, { "cell_type": "markdown", "id": "eb0e78ff-5c00-4395-8b38-f81aa77f4267", "metadata": {}, "source": [ "We can see how, as expected, the acceleration is highest at the poles, due to Earth's flattening. Also, despite the fact that we are not using a high-resolution model, the [Indian Ocean Geoid Low](https://en.wikipedia.org/wiki/Indian_Ocean_Geoid_Low) is also clearly visible." ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.12.4" } }, "nbformat": 4, "nbformat_minor": 5 }