{ "cells": [ { "cell_type": "markdown", "id": "f583c3b9", "metadata": {}, "source": [ "# High Frequency Moments of the Green's Functions and Self-energy\n", "\n", "\n", "We present here the derivation of the high-frequency moments of the Green's function and self-energy, as well as show case how these high frequency moments can be computed and utilized within [TRIQS/cthyb](https://triqs.github.io/cthyb/latest/).\n", "\n", "## Derivation\n", "The following derivation follows closely [Rev. Mod. Phys 83, 349 (2011)](https://journals.aps.org/rmp/abstract/10.1103/RevModPhys.83.349) and [Phys. Rev. B 84, 073104 (2011)](https://journals.aps.org/prb/abstract/10.1103/PhysRevB.84.073104).\n", "\n", "The high-frequency expansion of the Green's function takes the following form \n", "\n", "$$ \\mathbf{G}(i\\omega_{n}) = \\frac{1}{i\\omega_{n}}\\mathbf{G}_{0} + \\frac{1}{(i\\omega_{n})^{2}}\\mathbf{G}_{1} + \\frac{1}{(i\\omega_{n})^{3}}\\mathbf{G}_{2} + \\cdots \\ , $$\n", "where the coefficients are related to time derivatives of the Green's function:\n", "\n", "$$ \\mathbf{G}_{0} = -(G(0) + G(\\beta)) = \\delta_{ab}$$\n", "$$ \\mathbf{G}_{1} = \\partial_{\\tau}G(0) + \\partial_{\\tau}G(\\beta) = -\\langle \\{[H, d_{a}], d_{b}^{\\dagger} \\} \\rangle $$\n", "$$ \\mathbf{G}_{2} = -\\big (\\partial^{2}_{\\tau}G(0) + \\partial^{2}_{\\tau}G(\\beta) \\big ) = \\langle \\{[H, d_{a}], [d_{b}^{\\dagger},H] \\} \\rangle \\ , $$\n", "where $H$ is the Hamiltonian of the system. For the purpose of this cthyb tutorial we assume $H$ to be the general Anderson impurity hamiltonian, \n", "\n", "$$H = \\sum_{k,\\alpha} \\varepsilon_{k\\alpha}c^{\\dagger}_{k\\alpha}c_{k\\alpha} + \\sum_{k\\alpha b} \\big ( V_{k}^{\\alpha b} c^{\\dagger}_{k\\alpha}d_{b} + \\mathrm{h.c.}\\big ) + \\sum_{ab}E_{ab}d^{\\dagger}_{a}d_{b} + \\frac{1}{2}\\sum_{abcd}U_{abcd}d^{\\dagger}_{a}d^{\\dagger}_{b}d_{d}d_{c} \\ .$$\n", "\n", "In DMFT the coefficients for the Weiss field $\\mathbf{G}_{0}$ can be obtained from the same expressions as the Green's function with the hamiltonian replaced with only the non-interacting part. \n", "\n", "To derive the analogous coefficients for the self-energy, we start with the Dyson equation,\n", "\n", "$$\\mathbf{\\Sigma} = \\Big ((\\mathbf{G}_{0} )^{-1} - \\mathbf{G}^{-1} \\Big ) \\ .$$\n", "\n", "Using the high-frequency expression for the Green's function and Weiss field we obtain,\n", "\n", "$$ \\mathbf{\\Sigma} = \\frac{1}{\\frac{1}{i\\omega_{n}} + \\frac{1}{(i\\omega_{n})^{2}}(\\mathbf{G}_{0})_{1} + \\frac{1}{(i\\omega_{n})^{3}} (\\mathbf{G}_{0})_{2}}-\\frac{1}{\\frac{1}{i\\omega_{n}} + \\frac{1}{(i\\omega_{n})^{2}}\\mathbf{G}_{1} + \\frac{1}{(i\\omega_{n})^{3}} \\mathbf{G}_{2}} \\ .$$\n", "\n", "Expanding this expression for large $i\\omega_{n}$ and simplifying yields,\n", "\n", "$$\\mathbf{\\Sigma} = \\Big (\\mathbf{G}_{1} - (\\mathbf{G}_{0})_{1} \\Big ) + \\frac{1}{i\\omega_{n}} \\Big (\\mathbf{G}_{2} - (\\mathbf{G}_{0})_{2} +((\\mathbf{G}_{0})_{1})^{2} - \\mathbf{G}_{1}^{2} \\Big) + \\mathbf{O}\\Big ( (i\\omega_{n})^{-2}\\Big )$$\n", "\n", "We now identify the leading terms,\n", "\n", "$$ \\mathbf{\\Sigma}_{\\infty} \\equiv \\mathbf{G}_{1} - (\\mathbf{G}_{0})_{1} $$\n", "$$ \\mathbf{\\Sigma}_{1} \\equiv \\mathbf{G}_{2} - (\\mathbf{G}_{0})_{2} +((\\mathbf{G}_{0})_{1})^{2} - \\mathbf{G}_{1}^{2} $$\n", "\n", "Note, that $(\\mathbf{G}_{0})_{1}$ are just the crystal field levels inside of the impurity problem $\\mathbf{E}$. We can simplify these expressions by plugging in their definitions from above,\n", "\n", "$$\\Sigma_{\\infty}^{ab} = -\\langle \\{[H, d_{a}], d_{b}^{\\dagger} \\} \\rangle + \\langle \\{[H-H_{\\mathrm{int}}, d_{a}], d_{b}^{\\dagger} \\} \\rangle = - \\langle \\{[H_{\\mathrm{int}}, d_{a}], d_{b}^{\\dagger} \\} \\rangle $$\n", "\n", "$$\\Sigma_{1}^{ab} = \\langle \\{ [H_{\\mathrm{int}}, [H_{\\mathrm{int}}, d_{a} ]], d_{b}^\\dagger \\} \\rangle - (\\Sigma_{\\infty}^{ab} )^{2} $$\n" ] }, { "cell_type": "markdown", "id": "74d0c030", "metadata": {}, "source": [ "## Implentation in TRIQS/cthyb\n", "\n", "The leading two moments of $\\mathbf{G}$ and of $\\mathbf{\\Sigma}$ are automatically calculated and stored as members of the TRIQS/cthyb solver when ``measure_density_matrix = True``. The moments for the Green's function and self-energy are stored as ``Solver.G_moments`` and ``Solver.Sigma_moments``, with the matching block structure of the ``Gf`` objects. Additionally, the Hartree shift $\\mathbf{\\Sigma}^\\infty$ is stored as ``Solver.Sigma_Hartree``.\n", "\n", "The moments of the Green's function are used in the Fourier transform of $G(\\tau)$ to $G(i\\omega_{n})$ in the solver post-processing automatically if available. The moments of the self-energy are used in the tail fitting routine, where the first two ``known_moments`` now correspond to the calculated first and second moments of the self-energy." ] }, { "cell_type": "code", "execution_count": 1, "id": "0d5e94f5", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Warning: could not identify MPI environment!\n" ] }, { "name": "stderr", "output_type": "stream", "text": [ "Starting serial run at: 2023-02-02 09:38:30.738159\n" ] } ], "source": [ "from triqs.gf import *\n", "from triqs.operators import *\n", "from triqs_cthyb import Solver\n", "import triqs.utility.mpi as mpi\n", "\n", "import numpy as np\n", "from triqs.plot.mpl_interface import plt,oplot" ] }, { "cell_type": "code", "execution_count": 2, "id": "faa8b046", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "╔╦╗╦═╗╦╔═╗ ╔═╗ ┌─┐┌┬┐┬ ┬┬ ┬┌┐ \n", " ║ ╠╦╝║║═╬╗╚═╗ │ │ ├─┤└┬┘├┴┐\n", " ╩ ╩╚═╩╚═╝╚╚═╝ └─┘ ┴ ┴ ┴ ┴ └─┘\n", "\n", "The local Hamiltonian of the problem:\n", "-2*c_dag('down',0)*c('down',0) + -2*c_dag('up',0)*c('up',0) + 4*c_dag('down',0)*c_dag('up',0)*c('up',0)*c('down',0)\n", "Using autopartition algorithm to partition the local Hilbert space\n", "Found 4 subspaces.\n", "\n", "Warming up ...\n", "09:38:33 0% ETA 00:00:15 cycle 661 of 100000\n", "09:38:35 12% ETA 00:00:14 cycle 12541 of 100000\n", "09:38:38 30% ETA 00:00:10 cycle 30560 of 100000\n", "09:38:41 52% ETA 00:00:07 cycle 52654 of 100000\n", "09:38:45 75% ETA 00:00:03 cycle 75715 of 100000\n", "\n", "\n", "\n", "Accumulating ...\n", "09:38:48 0% ETA 00:01:09 cycle 723 of 500000\n", "09:38:50 3% ETA 00:01:08 cycle 15022 of 500000\n", "09:38:53 6% ETA 00:01:05 cycle 32974 of 500000\n", "09:38:56 11% ETA 00:01:02 cycle 55408 of 500000\n", "09:39:00 16% ETA 00:00:58 cycle 83197 of 500000\n", "09:39:05 23% ETA 00:00:54 cycle 117863 of 500000\n", "09:39:11 32% ETA 00:00:48 cycle 160877 of 500000\n", "09:39:19 43% ETA 00:00:40 cycle 215435 of 500000\n", "09:39:29 56% ETA 00:00:31 cycle 281909 of 500000\n", "09:39:41 71% ETA 00:00:21 cycle 355564 of 500000\n", "09:39:56 91% ETA 00:00:06 cycle 455301 of 500000\n", "\n", "\n", "[Rank 0] Collect results: Waiting for all mpi-threads to finish accumulating...\n", "[Rank 0] Timings for all measures:\n", "Measure | seconds \n", "Auto-correlation time | 0.144108 \n", "Average order | 0.0375342 \n", "Average sign | 0.0416869 \n", "Density Matrix for local static observable | 0.732204 \n", "G_tau measure | 0.0985052 \n", "Total measure time | 1.05404 \n", "[Rank 0] Acceptance rate for all moves:\n", "Move set Insert two operators: 0.0448675\n", " Move Insert Delta_up: 0.044915\n", " Move Insert Delta_down: 0.0448199\n", "Move set Remove two operators: 0.0448576\n", " Move Remove Delta_up: 0.04491\n", " Move Remove Delta_down: 0.0448052\n", "Move set Insert four operators: 0.00347959\n", " Move Insert Delta_up_up: 0.00269973\n", " Move Insert Delta_up_down: 0.00431061\n", " Move Insert Delta_down_up: 0.00423298\n", " Move Insert Delta_down_down: 0.002676\n", "Move set Remove four operators: 0.00348254\n", " Move Remove Delta_up_up: 0.00270066\n", " Move Remove Delta_up_down: 0.00431291\n", " Move Remove Delta_down_up: 0.00424282\n", " Move Remove Delta_down_down: 0.0026742\n", "Move Shift one operator: 0.0350781\n", "[Rank 0] Warmup lasted: 15.3982 seconds [00:00:15]\n", "[Rank 0] Simulation lasted: 74.6982 seconds [00:01:14]\n", "[Rank 0] Number of measures: 500000\n", "Total number of measures: 500000\n", "Average sign: 1\n", "Average order: 1.2027\n", "Auto-correlation time: 0.358059\n" ] } ], "source": [ "D, V, U = 1.0, 0.2, 4.0\n", "\n", "ef, beta = -U/2.0, 50\n", "H = U*n('up',0)*n('down',0)\n", "S = Solver(beta=beta, gf_struct=[('up',1), ('down', 1)])\n", "\n", "for name, g0 in S.G0_iw: g0 << inverse(iOmega_n - ef - V**2 * Wilson(D))\n", "\n", "S.solve(h_int = H,\n", " n_cycles =500000,\n", " length_cycle = 200,\n", " n_warmup_cycles = 100000,\n", " measure_density_matrix=True, # moments will be calculated automatically\n", " use_norm_as_weight = True,\n", " )" ] }, { "cell_type": "markdown", "id": "d6a0f954", "metadata": {}, "source": [ "Now, we can compare the high-frequency expansion of the self-energy to the real and imaginary parts of the self-energy obtained from the Dyson equation using the sampled Green's function and input Weiss field. Let us first have a look at the calculated moments and then plot the real part, corresponding to the overall Hartree shift, and the imaginary part of the measured self-energy together with the high-frequency expansion obtained from the moments:" ] }, { "cell_type": "code", "execution_count": 3, "id": "b8570aba-c879-423e-8c69-4f4a2d871ade", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "-------\n", "up:\n", "𝚺_0= 2.0033\n", "𝚺_1= 4.0000\n", "-------\n", "down:\n", "𝚺_0= 1.9967\n", "𝚺_1= 4.0000\n" ] } ], "source": [ "for block, moment in S.Sigma_moments.items():\n", " print('-------')\n", " print(f'{block}:')\n", " for i, value in enumerate(moment):\n", " print(f'𝚺_{i}= {value[0,0].real:1.4f}')" ] }, { "cell_type": "code", "execution_count": 4, "id": "233a2954-8d33-4928-8f11-bbfa64cf298d", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "-------\n", "up:\n", "Gimp_0= 0.0000\n", "Gimp_1= 1.0000\n", "Gimp_2= 0.0033\n", "-------\n", "down:\n", "Gimp_0= 0.0000\n", "Gimp_1= 1.0000\n", "Gimp_2= -0.0033\n" ] } ], "source": [ "for block, moment in S.G_moments.items():\n", " print('-------')\n", " print(f'{block}:')\n", " for i, value in enumerate(moment):\n", " print(f'Gimp_{i}= {value[0,0].real:1.4f}')" ] }, { "cell_type": "code", "execution_count": 5, "id": "69705321-21c3-4c51-8d57-63451e54bce8", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAABLkAAANvCAYAAAAx+whxAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjYuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/av/WaAAAACXBIWXMAAB7CAAAewgFu0HU+AAEAAElEQVR4nOzdeXhTZd7G8fukO1AEgbIVEBDEYlkFUYsgsggKDIIojuKGy7iPuLyOI4uvMioy+Cqjg8IMLgiCA8MiiKKCFnVYC5WKyCa0BWQRC13okrx/MImlbdJsTXKS7+e6uIzJyTlP0pyTc+78nucxbDabTQAAAAAAAICJWYLdAAAAAAAAAMBXhFwAAAAAAAAwPUIuAAAAAAAAmB4hFwAAAAAAAEyPkAsAAAAAAACmR8gFAAAAAAAA0yPkAgAAAAAAgOkRcgEAAAAAAMD0CLkAAAAAAABgeoRcAAAAAAAAMD1CLgAAAAAAAJgeIRcAAAAAAABMj5ALAAAAAAAApkfIBQAAAAAAANMj5AIAAAAAAIDpEXIBAAAAAADA9Ai5AAAAAAAAYHqEXAAAAAAAADA9Qi4AAAAAAACYHiEXAAAAAAAATI+QywTy8vI0f/58jR8/Xn369NH555+vc845R7GxsUpKSlLfvn310ksv6dixY37b5sqVKzVixAglJycrLi5OycnJGjFihFauXOm3bQAAAAAAAPiLYbPZbMFuBFxbvXq1BgwYUO1yDRs21HvvvadBgwZ5vS2r1aq7775bs2fPdrrMuHHjNHPmTFksZKQAAAAAACA0RAe7AXBPixYtdOWVV6p79+5q0aKFmjZtKqvVquzsbH344YdatGiRjh49qmHDhmn9+vXq3LmzV9t5+umnHQFX165d9cQTT6ht27bavXu3XnrpJW3ZskWzZs1So0aNNGXKFH++RAAAAAAAAK9RyWUCZWVlioqKcrnMv//9b40YMUKSNGLECC1atMjj7ezcuVMdO3ZUaWmpLr74Yn355ZdKSEhwPF5QUKA+ffpo48aNio6O1vfff6/zzz/f4+0AAAAAAAD4G/3NTKC6gEuSfve73+mCCy6QJH311VdebeeVV15RaWmpJOm11147K+CSpFq1aum1116TJJWWlmr69OlebQcAAAAAAMDfCLnCSGJioiSpqKjI4+fabDYtWbJEktShQwf16tWryuV69erlCNOWLFkiCgEBAAAAAEAoIOQKEz/88IMyMjIknQmpPLV3717l5uZKkvr06eNyWfvjOTk52rdvn8fbAgAAAAAA8DcGnjexgoIC5eTkaNmyZXrppZccXQ0feeQRj9eVlZXluF1dSFb+8e+//16tW7f2aFvZ2dkuHy8qKtKOHTvUuHFjNWrUSNHRfEwBAAAAAPC30tJSHTlyRJKUmpqq+Pj4ILfIN6QHJjNnzhzdfvvtTh//n//5H910000er7d88JScnOxy2RYtWjhuHzhwwONtlX8+AAAAAAAIvvXr16tHjx7BboZPCLnCRJcuXfTmm296/YE8efKk43adOnVcLlu7dm3H7VOnTnm1PQAAAAAAAH8i5DKZ3/3ud7r44oslSYWFhdq9e7cWLFigxYsXa8yYMXrllVd07bXXerze8oPVx8bGulw2Li7OcbuwsNDjbVVX/XXgwAFddtllks4kyU2bNvV4GwAAAAAAwLWDBw+qZ8+ekqRGjRoFuTW+I+QymXr16qlevXqO/+/Ro4duvPFGvfvuu7r11ls1fPhwzZ49W7fddptH6y3f77a4uNjlsqdPn3bcTkhI8Gg7UvXdIctr2rSpR8sDAAAAAADPhcN42MyuGCZuueUWXX/99bJarXrggQd0/Phxj56fmJjouF1dF8T8/HzH7eq6NgIAAAAAAAQCIVcYGT58uKQzIdTHH3/s0XPLV0tVN/th+e6GDCIPAAAAAABCASFXGCnff/ann37y6LkpKSmO2zt27HC5bPnHL7zwQo+2AwAAAAAAUBMIucJITk6O47an3Qhbt26tZs2aSZLWrl3rctkvv/xSktS8eXOdd955njUSAAAAAACgBhByhZGFCxc6bqempnr0XMMwHN0dd+zYoW+//bbK5b799ltHJdfw4cNlGIaXrQUAAAAAAPAfQi4TmDNnjoqKilwuM336dK1YsULSmaqs3r17n/X4mjVrZBiGDMNwOvPiI488oqioKEnSgw8+qMLCwrMeLyws1IMPPijpzKwLjzzyiBevBgAAAAAAwP/MPz9kBJg0aZLGjx+vkSNHKi0tTW3btlWdOnV08uRJZWZmau7cuVq3bp0kKTY2Vm+++aYjrPJE+/bt9fjjj+uFF17Qxo0bdfnll+vJJ59U27ZttXv3br344ovasmWLJOnxxx9Xu3bt/Po6AQAAAAAAvEXIZRLHjx/XW2+9pbfeesvpMsnJyfrHP/6h/v37e72d559/Xj///LP+8Y9/aMuWLbrxxhsrLXPnnXfqueee83obAAAAAAAA/kbIZQKrVq3SRx99pHXr1mnXrl06fPiwjh07poSEBCUlJalLly669tprNXr0aNWqVcunbVksFs2ePVsjR47Um2++qQ0bNujo0aNq2LChevTooXvuuUeDBw/20ysDAAAAAADwD8Nms9mC3QigvOzsbLVo0UKSdODAASUnJwe5RQAAAAAAhJ9wu/5m4HkAAAAAAACYHiEXAAAAAAAATI+QCwAAAAAAAKZHyAUAAAAAAADTI+QCAAAAAACA6RFyAQAAAAAAwPQIuQAAAAAAAGB6hFwAAAAAAAAwPUIuAAAAAAAAmB4hFwAAAAAAAEyPkAsAAAAAAACmR8gFAAAAAAAA0yPkAgAAAAAAgOkRcgEAAAAAAMD0CLkAAAAAAABgeoRcAAAAAAAAMD1CLgAAAAAAAJgeIRcAAAAAAABMj5ALAAAAAAAApkfIBQAAAAAAANMj5AIAAAAAAIDpEXIBAAAAAADA9Ai5AAAAAAAAYHqEXAAAAAAAADA9Qi4AAAAAAACYHiEXAAAAAAAATI+QCwAAAAAAAKZHyAUAAAAAAADTI+QCAAAAAACA6RFyAQAAAAAAwPQIuQAAAAAAAGB6hFwAAAAAAAAwPUIuAAAAAAAAmB4hFwAAAAAAAEyPkAsAAAAAAACmR8gFAAAAAAAA0yPkAgAAAAAAgOkRcgEAAAAAAMD0CLkAAAAAAABgeoRcAAAAAAAAMD1CLgAAAAAAAJgeIRcAAAAAAABMj5ALAAAAAAAApkfIBQAAAAAAANMj5AIAAAAARDyr1aaC4lJZrbZgNwWAl6KD3QAAAAAAAIIlKzdPs9L3aGXmIRWWlCkhJkqDU5toXFobpTSrG+zmAfAAIRcAAAAAICItycjR+AVbVVquequwpEyLNudoaUaupo3urOFdmgetfVarTUWlZYqPjpLFYgStHYBZEHIBAAAAACJOVm5epYCrvFKrTeMXbFW7pMSAV3RRXQZ4hzG5AAAAAAARZ1b6HqcBl12p1abZ6XsD1KIzlmTkaNiMdC3anKPCkjJJv1WXDZuRriUZOQFtD2AmhFwAAAAAgIhitdq0MvOQW8uuyDwYsMHo3a0uy8rNC0h7ALMh5AIAAAAARJSi0jJHlVR1CkvKVFTq3rK+CtXqMsAsCLkAAAAAABElPjpKCTFRbi2bEBOl+Gj3lvVFqFaXAWZCyAUAAAAAiCgWi6HBqU3cWnZIatOAzGwYqtVlgJkQcgEAAAAAIs64tDaKria8irYYujOtdUDaE4rVZYDZEHIBAAAAACJOSrO6mja6s9OgK9piaNrozkppVjcg7QnF6jLAbAi5AAAAAAARaXiX5lr6QJpGdkt2VFElxERpZLdkLX0gTcO7NA9oe0Ktugwwm+hgNwAAAAAAgGCxV3RNHdVJRaVlio+OClqVlL0t4xdsrXKWxUBXlwFmQ8gFAAAAADAlq9Xmt2DKYjFUKzb4l8jDuzRXu6REzU7fqxWZB1VYUqaEmCgNSW2qO9NaE3ABLgR/DwYAAAAAwANZuXmalb5HKzMPOUKgwalNNC6tTViEQKFUXQaYCSEXAAAAAMA0lmTkVOrOV1hSpkWbc7Q0I1fTRncO+FhaNSVUqssAs2DgeQAAAACAKWTl5jkdr0qSSq02jV+wVVm5eQFuGYBQQMgFAAAAADCFWel7nAZcdqVWm2an7w1QiwCEEkIuk9i4caOeffZZDRw4UMnJyYqLi1OdOnXUvn173X777UpPT/fLdiZNmiTDMNz6t2bNGr9sEwAAAACqY7XatDLzkFvLrsg8KGs1YRiA8EPnXhO44oor9NVXX1W6v7i4WD/++KN+/PFHzZkzR2PHjtVbb72l2NjYILQSAAAAAGpOUWmZCkvK3Fq2sKRMRaVljGcFRBj2eBPIzc2VJDVr1kzXX3+9evfurZYtW6qsrEzffPONpk2bppycHL3zzjsqKSnR+++/75ftZmZmuny8devWftkOAAAAAFQnPjpKCTFRbgVdCTFRio+OCkCrAIQSQi4Xjhw5oj179ujQoUPKz89XTEyM6tWrp5YtW+r8889XVFRgDpodOnTQlClTNHLkyErb7NWrl2655RZdfvnl2rlzp+bNm6d7771XV1xxhc/bveiii3xeBwAAAAD4g8ViaHBqEy3anFPtskNSm8piMQLQKvdZrTYVlZYpPjrK57b5c11m2jZQHUKucvLz87VkyRKtXLlSa9euVU6O84NnXFycunbtqoEDB2rEiBHq1KlTjbVr+fLlLh9v2LChpk2bpqFDh0qSPvzwQ7+EXAAAAAAQSsaltdHSjFyXg89HWwzdmRY6vU6ycvM0K32PVmYeUmFJmRJiojQ4tYnGpbVRSrO6QVuXp4K5bcBdhs1mi/jR+LZs2aLXXntNCxcuVEFBgSTJ3bfFMM4k1x07dtT999+vW265RbVq1aqxtjqTn5+vOnXqSJKGDBmijz76yKv1TJo0SZMnT5bk/nvgb9nZ2WrRooUk6cCBA0pOTg5KOwAAAACEniUZORq/YGuVQVe0xdC00Z01vEvzILSsMn+2NZiv20zvOTwTbtffEV3JtWXLFj3zzDNauXKlpN9CnSZNmqhnz57q3r27kpKSdO6556p+/foqLCzU8ePH9csvv2jnzp3asGGDtm3bppKSEn333Xe677779Mwzz+iJJ57Qgw8+qLi4uIC9ltOnTztuB6obJQAAAAAE2vAuzdUuKVGz0/dqReZBR1XRkNSmujOttVKa1Q2JLnVZuXlOgyFJKrXaNH7BVrVLSqy2Esqf6/JUMLcNeCpiQ67bb79d7777rqxWqySpW7du+v3vf6+RI0eqZcuWbq+nuLhYX375pebOnavFixfr6NGjevLJJ/X666/rnXfeUVpaWk29hLOsXbvWcfvCCy/0yzoHDhyojIwMnThxQvXq1VNKSoquvvpq3XPPPapfv77X683Oznb5+MGDB71eNwAAAIDwl9KsrqaN7qypozqdFWZl5ebp0QUZIdGlblb6HpfdKqUzAdHs9L2aNrpzwNblqWBuG/BUxHZXtFgsio2N1a233qrx48erffv2Pq/z9OnTWrhwoaZMmaIdO3Zo0qRJmjBhgh9a65rVatWll16q9evXS5I2btyo7t27e7Wu8t0VnalXr57mzJmj4cOHe7UNexdPd4RDuSQAAACAmhdKXeqsVps6Tlzl9kyQmRMHqthqrbLyzNN1bZ88yG/Va8HcNgKD7oph4r777tOTTz7p+GP6Q1xcnG6++Wb9/ve/18KFC1VWVv2BwB+mT5/uCLiuu+46rwMuu9TUVP3ud79Tz5491axZM5WUlOiHH37Q3Llz9cknn+jEiRMaOXKkli1bpsGDB/vjJQAAAACA10KtS11RaZlbwZAkFZaUKXXyJ04rzzxdV1FpmWrF+udSP5jbBrwRsZVc4WLt2rXq37+/SktLlZSUpMzMTCUlJXm9PnvXRGdmzpype++9V5LUrFkz7d69W/Hx8R5tw53uij179pQUHkkyAAAAgJr16IIMLdqcU+1yI7slB6RLnScVUFUpX3lGJRdqUrhVclmC3QB4b/v27RoxYoRKS0sVHx+vhQsX+hRwSXIZcEnSPffcozvvvFOSlJubq3/9618ebyM5Odnlv6ZNm3rTdAAAAAARyGq1aWXmIbeWXZF5UNZqxpfyB4vF0ODUJl4/3155lpWb59G6hqQ29WvIFMxtA96I6JDrj3/8ozIyMoLdDK/s3btXAwcO1C+//KKoqCjNnz9fV1xxRUC2fc899zhulx/wHgAAAAACzZsudYEwLq2Non0IfeyDubu7rmiLoTvTWnu9PWeCuW3AUxEdcv3f//2funfvrk6dOunll182zax+ubm56t+/v3Jzc2UYhv7xj394PQi8N1JSUhy3c3KqLwmuKVarTQXFpQH5JQYAAABAaIqPjlJCTJRbyybERCk+2r1lfWWfAdKXoMteeVbduuzdG2tivLFgbhvwVMSPCGez2bR9+3Y9+eSTeuqpp9SvXz/deuutGjFihBISEoLdvEqOHj2qAQMGaM+ePZKk1157TWPHjg1oGzyZHbEmZOXmaVb6npCYFhgAAABAcNm71LkzJlegu9QN79Jc7ZISNTt9r1ZkHlRhSZnioy0qKrW69fzyg7lXta6EmCgNSW2qO9Na1+i1UDC3DXgiogee//TTT/Xuu+9q8eLFys/Pl/RbgFOnTh2NHDlSt9xyi6688spgNtPh119/Vb9+/bR582ZJ0gsvvKAnn3wy4O3YuHGjevToIUkaN26c3nrrLb+u39XAd6E0LTAAAACA0JCVm6dhM9Kdzq4onbleWPpAWtACGavVpqLSMsVaLI7ZFKvjbDB3+7rio6MCPg5WMLcN/2Pg+TAyYMAAvfPOOzp8+LDeeecdDRgwQIZhyGaz6eTJk3r77bfVv39/tWrVSk8//bR27NgRtLYWFBTommuucQRcTz/9dFACLunMDIt2ffr0Cdh23Z0WOCs3L2BtAgAAABB8ZuhSZ7EYqhUbrehoi8+DudvXFYyQKZjbBqoT0SGXXa1atXTzzTdr1apVOnDggF566SV16tRJNptNNptNBw4c0AsvvKCOHTuqZ8+e+tvf/qZjx44FrH3FxcUaMWKE1q1bJ0l6+OGH9dxzz3m8njlz5sgwDBmGoUmTJlV6PDMzU7t27XK5jjfffFOzZs2SJDVp0kQjRozwuB3empW+x+UvM9LZgzMCAAAAiBzDuzTX0gfSNLJbsmOMroSYKI3slqylD6SFVI8PBnMHakbEj8lVUdOmTfXYY4/pscceU2Zmpt555x3NmzdPubm5kqRNmzZp06ZNGj9+vK6++mqNHTtWQ4cOVUxMTI21acyYMfrkk08kSf369dOdd96p7777zunysbGxat++vcfb2bRpk8aNG6crr7xSgwcPVmpqqho0aKDS0lLt2LFDc+fOdbQjKipKb775pmrXru3di/KQp9MCTx3ViV8WAAAAgAhjr+iaOqpTSHeps7ezuqFYGOsK8AwhlwupqamaOnWqXnrpJX322Wd69913tWjRIuXn56u4uFjLli3TsmXLVL9+fd1www265ZZb1KtXL7+3Y9GiRY7bn3/+uTp16uRy+VatWmnfvn1ebausrEyrV6/W6tWrnS7ToEEDzZ49W0OHDvVqG56yWm06XnDa42mBa8Xy8QYAAAAikb1LXShjMHfA/0J7rw8RhmGof//+6t+/v9544w0tWrRI7777rj7//HOVlZXp+PHjeuONNzRz5kyVlpYGu7leGzJkiGbPnq1vvvlGW7Zs0eHDh3Xs2DHZbDade+656ty5s66++mrddtttqls3MAfcycu2Kz33O7cDLimw0wIDAAAAgLfMUnkGmEVEz67oq+3bt+vGG29UVlaWbDabDMNQWZn7YQyqVn52h+Z/mKPoug09ev7IbsmaNrpzTTQNAAAAAICwEW6zK1LJ5aGSkhItW7ZM7733nlasWKGSkpJgNwnlMDgjAAAAAACRiZDLTevWrdO7776rhQsX6sSJE5IkexFcYmKiRo0apVtvvTWILQSDMwIAAAAAELkIuVz48ccf9e6772ru3LmOgdztwVZUVJT69++vsWPHasSIEYqPjw9iSyMbgzMCAAAAAABCrgqOHj2q+fPn691339XGjRsl/RZsSWdmXBw7dqx+//vfq0mTJsFqJv5r45+v0rm14hicEQAAAACACEfIJen06dNasmSJ3nvvPa1atcoxQ6I93GrcuLFuuukmjR07Vp07M6B5qEiIiSLgAgAAAAAAkiI85FqzZo3ee+89/etf/1JeXp6k34Kt+Ph4DRs2TGPHjtWgQYMUFRUVzKaiCkNSmxJwAQAAAAAASREecvXr10+GYTiCLcMwlJaWprFjx2r06NGqW5fxnUIVsygCAAAAAIDyIjrkks5UbrVt21a33HKLbrnlFrVuTXAS6phFEQAAhAqr1aai0jLFR0dRYQ4AQJBFdMh19913a+zYsbrsssuC3RQ4MSS1qb7KtaqwpIxZFAEAQMjIys3TrPQ9Wpl5yHGeMji1icalteE8BQCAIDFs5acOBEJAdna2WrRoIUk6cOCAmjVrzi+kAAAgZCzJyNH4BVtVaq18Gm2vOB/epXkQWgYAgGcqXn8nJycHuUW+iehKLnfs3r1b33zzjQ4dOqSCggLdd999atiwYbCbFVEsFkO1Yt3/qNJtAAAA1JSs3DynAZcklVptGr9gq9olJVLRBQBAgBFyObF582Y98sgjWrdu3Vn3jxo16qyQ629/+5smT56sc845R1lZWYqJiQl0U/FfdBsAAAA1bVb6HqcBl12p1abZ6Xs1bXTnALUKAABIkiXYDQhFy5cv1+WXX65169bJZrM5/lVl7NixKiws1J49e7R8+fIAtxR2SzJyNGxGuhZtzlFhSZkkqbCkTIs2n7l/SUZOkFsIAADMzmq1aWXmIbeWXZF5UNZqwjAAACqyWm0qKC7lO8RLhFwVHDx4UGPGjNHp06eVkpKilStX6uTJk06XT0xM1LBhwyRJK1euDFQzUY673QaycvMC3DIAABBOikrLHD+mVaewpEwFxaVcqAAA3JKVm6dHF2So48RVSpmwSh0nrtKjCzK4jvUQ3RUrmD59uvLz89WqVSt99dVXqlevXrXP6du3r+bNm6dNmzbVfANRCd0GAABAIMRHRykhJsqtoCvKMNTj+c8YQgEAUK2qJjSx90xampHLhCYeoJKrgo8//liGYWj8+PFuBVyS1KFDB0nS3r17a7BlqArdBgAAQKBYLIYGpzZxa9kym40hFAAA1aJnkn8RclXw008/SZJ69uzp9nPq1j3zi9ypU6dqpE1wztNuA0Wl7i0LAABQlXFpbRTt5ezNXKgAACrypGcSqkfIVUFpaakkyWq1uv2cX3/9VZJUp06dGmkTnLN3G3BHQkyU4qPdWxYAAKAqKc3qatrozj4FXVyoAAAkeibVBEKuCpo0OVOCvmfPHrefs379eklSy5Yta6RNcM6TbgNDUpvK4uUJKQAAgN3wLs219IE0jeyW7PixLT7aoig3zzO4UAFgBszyV/PomeR/DDxfQe/evbV3714tXLhQN910U7XLFxcXa+bMmTIMQ3379q35BqKScWlttDQj12WJZ7TF0J1prQPYKgAAEM7sFV1TR3VSUWmZrFabLpr0iVvPLSwp0/GC0zq3Vhw/wAEIOVm5eZqVvkcrMw8xeUYN82RCE3omuYdKrgpuu+02SdLSpUv16aefuly2uLhYY8eO1e7du2UYhu66664AtBAVVddtINpiaNrozhyQAQAwCTNVD1gshmrFRqtWbLTbQyhI0sXPfcb08ABCzpKMM5NkLNqcw+QZAUDPJP+jkquCvn376oYbbtAHH3ygoUOH6uGHH9bIkSMdj+/bt08nTpzQunXr9Oabb2rPnj0yDEP33nuvOnbsGMSWR7bhXZqrXVKiZqfv1YrMg45fHIakNtWdaa0JuAAAMAEzVw/YL1QWbXb/ApDp4QGEEndn+WuXlBjyx2QzoWeSfxk2my30fyILsNOnT2vkyJFasWKFDMN5Ump/66677jp98MEHioqidNAfsrOz1aJFC0nSgQMHlJyc7NHzrVabikrLFB8dRdINAIBJLMnIcXpxZa/KDvUQKCs3T8NmpFc7S1ZVoi2Glj6QxoUjgKB5dEGGW0H9yG7Jmja6cwBaFDmC+R3o6/V3qKG7YhXi4uK0fPlyzZw5U23atJHNZqvyX3Jysl5//XV9+OGHBFwhxN5tgIALAABzcLd6INS79fky8yKzLgIIJmb5C66qJjRJiInSyG7JWvpAWsj/yBNKqORyQ1ZWljZu3Kiff/5ZZWVlatCggbp27apu3bq5rPSCd8ItSQYAAK6FW/VAVm7eWUMouCshJkrbJw/ihzrA5MzYs6SguFQpE1a5vXzWs4NUK7by6EdmfO3B4Op9CvR7GG7X34zJ5YaUlBSlpKQEuxkAAABhx9PqgamjOoX8hVP5mRePF5zWxc995tbz7NPDV3XhCJhRpAUeZh5X0NdZ/sz82gPJnffJ3jMJ3uGdAwAAQNAUlZa5Xe1kthDIYjF0bq04podHWPAksIrEwKOqMZXMNLmEJ5NnVJzlz+yvPVB4nwKDMbkQMcw0HTkAAJHCXj3gDjOGQEwPD7PLys3Towsy1HHiKqVMWKWOE1fp0QUZTsfIW5KRo2Ez0rVoc44j3LVfyA+bka4lGe7PQFpeKJ/Lh8u4guPS2lQ7pmDFWf7C5bXXNN6nwInYkGvRokU1uv7c3Fx9++23NboNuMfTL2YAABA4kRACeXPhCIQCTwOrmriQN8O5/Kz0PdXOqmqGySWqmzzDPstf+Wq8cHntNY33KXAiNuQaNWqUunTpog8//NCv6z1w4IDuu+8+tW3bVp988olf1w3P1dQvSQAAwH/CPQTy5sIRCDZvAit/X8ib4Vw+3GYl9GSWv3B77TWF9ymwIjbkatu2rbZt26YbbrhBrVu31tNPP63t27d7ta78/Hy99957GjJkiNq2bau///3vKisrU9u2bf3caniCklAAAMwhEkIgpoeH2XgaWPn7Qt4s5/LejCsY6uzH5O2TBynr2UHaPnlQlcfgcHztNYH3KbDMMWpnDcjKytIrr7yil156ST/99JNeeOEFvfDCC2rXrp169eqlHj16qGvXrkpKSlL9+vVVv359FRYW6vjx4/rll1+0c+dObdiwQevXr9f69etVVFQkm+3MAfi6667TlClT1L59+yC/ysjmyRezGaYjBwAgnA3v0lztkhI1O32vVmQedAxWPSS1qe5Ma23qgMuu/KyLkTTjHMzHm1lP/T2JhFnO5X2dlTCUVTfLXzi/dn/ifQqsiA25YmJi9Pjjj+sPf/iDXn/9df3tb3/TgQMHtHPnTv34449699133VqPPdiKi4vTddddp4cfflg9e/asyabDDeE4HTkAAOEuUkIgpodHqPMmsPLnhXxNnst7MkukO3yZldDsIvm1e4L3KbAitruiXZ06dfTEE09o7969WrlypW6//Xa1atVKNput2n9xcXHq06eP/vrXvyonJ0dz584l4AoRlIQCAGBe9hCIE31EqmDPJOjNrKf+nESiJs7la3IA+3AfV9CVSH7tnuB9Chx+Qvovi8WiQYMGadCgQZKknJwcff3118rOztaRI0d0/PhxxcfHq1GjRmrUqJFSU1N18cUXKyYmJsgtR1UoCQUAIPz4uwIDCDVZuXmalb5HKzMPObrsDk5tonFpbQLaZdfbypNxaW20NCPXZTdDdy7k/X0uvyQjp9L4XvYB7Jdm5Gra6M4+jYtnr0J1NoZYOIwr6Ewkv3ZP8D4FDiGXE82bN9f1118f7GbAS5SEAgAQGvwRTIXKhb+ZRGogaObXXdNBjKe8Caz8dSHvz3N5dwewb5eU6NPxJBLGFXQmkl+7J3ifAsOw2QeVAkJEdna2WrRoIUk6cOCAkpOTvVpPVm6ehs1Ir/aLeekDaRxQAADwM38FU1Vd+NvZL5iZmfA31b3vZg6BXDF7EBqq563e7n9ZuXk+X8j76z15dEGGW2HZyG7JfhvA3mz7mT/ba7bXHiyh9D756/o7VBByIeT4cyfjxBgAAO/4cgLur+/fUL3wD1Wu3vcoQ+rWqr6+y8kzZQjkSjic7wUjiHGXL4GVrxfyvv5trVabOk5c5Xa3x+2TBwU9cAgks4fD8A9CLqCG+Xsn88cvSQAARApfL3r8GUxRgeE+d973qpglBHImHIJQswQxwdo3fDmXLyguVcqEVe5v69lBETPzaTiEw/CPcAu5ImMPRkSLlOnIAQDwlT/GBJqVvqfaoKXUatPs9L0ugymr1aaVmYfcaveKzIOaOqqTT9/vZq9ocOd9r4q/xiMKFn993oLJm5kEgxHE2Gc9DTRfzuWZjKpqgRqnDAgGS7AbAAQK05EDAOCcuxc9Wbl5TtfhaTBldRFOeHPh7wmr1aaC4lJZrTYtycjRsBnpWrQ5x7FNe7g3bEa6lmRUX00WLFarTaeKStx+36tiD4HMxp+ft2CyBzHuiKQgpiJvzuXtA9i7I5Imo/IkHIZvyn/XIDCo5AIAAIBfKmL8WZFSUxUYFSu24qItKi61ytkrD9WKhoqvw1f+qIYLtGBVQPm72x6zgtcsb2aJDGeBrpKNVGavDjYzKrkAAAAinL8qYvxZkVITFRhVVWyddhFw2YVaRUNVr8NX3lTDOROoyoVAV0Bl5ebp0QUZ6jhxlVImrFLHiav06IIMl9WN7hqX1kbR1XyGIymI8Sd7d0dn7699/KlICR5qukrWmUiqaDJzdXA4oJILcMGsg88CAOAJf1XE+LsixZ8VGNV1x6xOqFQ0+Po6nPFXCBTIyoVAVkD5Y7w6V+xBTHUDgUdKEONvw7s0V7ukRCajUuDHKYu0iibGOws+Qi6gCpF2MAYARDZPLnrioy2yWm2yWm1Vhgb+DKb8ceFv/8Fq1lfeDcxuF8wBv8vzdoD56oR6CORMILqiBeqi1Z0ghh9gvcdkVGeEUzgcisJhMgyzM2w2W/jXC8JUgj2FKdPpAgAi0aMLMty66ImyGCqz2lz+AOTv79Ks3DyPKzD8PWZVQkyUtk8eFNSLYqvVpo4TV/mti6JdtMXQ0gfSvA5osnLzNGxGerVBky/bcKWmz93c3TdGdkv220VrxTCrJn6AJTCLXIHYZ4N9XAgGT47RofCdYhfs629/o5IL0G9f8nuP5FNeCgCIGOUvct2piJGksv8+7urXeH93DfK0AsNV6OGtUBjw25NupXYJMVFKbX6ONu3/xfG3K88f3eCCXblQk13RgjVIt30mQcn/1TD0WEAguscG+7gQDMGaDANno5ILISeQSXLFL/kow1CZG7uEP3+pAwAg0Jxd5J7fqI7++ulOj8MhV7/GB7paxJ3qAU+FSrWBJ1UC8dEWbfxzf9WKjXZUAtVUCBRKlQv+/rwVFJcqZcIqt5fPenaQXy9a/V0NQ48FlBcpx4VAMevrppIrghw5ckR79uzRoUOHlJ+fr5iYGNWrV08tW7bU+eefr6go3wbhQ3BV9SXvTsAlhc7gswAAeMpVVUi0xdCjA9pr95F8x0WPOz8AlVptmvXVHj034qJK4UL5ipRA8PeYVaE04LcnY+lc06mZ6sTHOP6/psYjCrXKBX9/3gI9SHdF/qyGYUBs74Vr185IOS4ESiDHO4Nz5v8k+VF+fr6WLFmilStXau3atcrJcf7hjIuLU9euXTVw4ECNGDFCnTp1CmBL4StfZyYKp4MxACByuHOR+9dPd2rpA2maOqqTCopL1eP5z9y6WFm0JUeLtuQEteuTJ13LqmJIio226HSpNWRnXvN1oPVwC4FqWjAvWv3dVTISu4/5KlK6dnJc8J9ATIYB1yzBbkAo2LJli+644w41adJEt9xyi95//31lZ2fLZrM5/VdUVKRvvvlGzz77rLp27apOnTpp5syZKigoCPbLgRt8/ZW3/MxSAACYhScXuRaLIYvF8HgMKHtV2LAZ6VqSUX0w4E/ejFllF20x9MqNXfT9s1cr69lB2j55UMhUcJVnr7yIdhJmBLryzB4CucOslQvj0to4fb/tauKi1ZtqGGc8Ccw+2parU0UlEX+euyTjzHFs0eYcx98hmMc3M4mE44IzoXaMjkQRXYayZcsWPfPMM1q5cqUkyT48WZMmTdSzZ091795dSUlJOvfcc1W/fn0VFhbq+PHj+uWXX7Rz505t2LBB27ZtU0lJib777jvdd999euaZZ/TEE0/owQcfVFxcXDBfHpzw9VdeSSqx2nTRpE8q/ZoTrqXMAADz86YqxJNf4ysKRtcnb9pbVcVWqFdq1+RA694I98qFQAzSXRV/VsN4EpgVlVqrPM+NJHTt9F24HxdcCbVjdKSJ2IHnb7/9dr377ruyWq2SpG7duun3v/+9Ro4cqZYtW7q9nuLiYn355ZeaO3euFi9erLy8PBmGoVatWumdd95RWlpaTb2EsFXTA995OoCoO6IMqVur+vouJy+sS5kBAObl7QDajy7IcKurljOBnqzF3fZe17V5lWOImU2o/MAWCQOa19Qg3a64+3mubj/zZEDsisLl7+cJf73vkS4SjgvVCZVjtCvhNvB8xHZXfPvttxUdHa277rpLO3bs0MaNG/XHP/7Ro4BLkmJjY9W/f3/985//1OHDh/XOO+/oggsu0L59+/T555/XUOvhC/uvYv5UZpM27PuFUmYAQMjy5PuvfFWIO121XFmReTCg3Z7c7Vo2rncbx8yDZmYfSyfYr2N4l+Za+kCaRnZLdnzOEmKiNLJbspY+kBYWF7L2iq7tkwcFrEurv7pKetJ9rCJ71VJWbp5XzzcbT6teI71bpyuRcFyoTqgcoyNJxFZyPfDAA3ryyScdiaU/2Ww2LVy4UGVlZRozZozf1x/uApEku/vrTJTFUJnV5tbMUs54O+24GVJ/AIC5eFud4OrXeHfYq8ICheqB4OIcxr/89XnOys3TsBnpXu/HkVK15G3VK1zjuBB47r7n4VbJFbF744wZM2ps3YZhaPTo0TW2fvjO3T7i/77/cp3XoJbbM0tVxdNZaiJlFhcAQOB5O0ZKVeOLuCsYM2cxHkpw+Xumtkjnr89zdWOLVcedGRzDQSTPDFiTOC4ETqRfT0ZsJRdCV6CSZHd/FfPHGF4JMVHaPnlQtScF7rRpaKdm/AoCACYUKr9i+1oVYn8df178nRZtCf0xa/z1vofK3w+RzR+fw4pji3n03AipWjLzmFwcqyKbN9/xVHIBYcLdX8V8mVnKzj6ts6uTAndmcXlkfoae+HCbTpdaIy6RDxRODAD4W6j9ouprVYj91/hxvdto6dbQnznL1+qBUPv7IbL5oxrGXtE1dVQnFRSXut1jIZKqlsw4MyDHKjAr6BlUciHkBCNJri7Y8HVmKXcqubzdBmOL+AcnBv5DUAj8JtTHhvJ1fw311+ercH99gGTuqqWaZKb930xtRc3xdl8Ot0quiJ1d0ReHDx/WiBEjlJiYqNjYWF1wwQW64YYb9MILL+jjjz/W4cOHg91EeKi6WS98nVlqSGpTSWcGsqxqBhZPZnGpqPyMN1arrcptOLvfGU+XN7slGWdmwly0OYcZMn2QlZunRxdkqOPEVUqZsEodJ67SowsyImY2JphDII9v7v6iGsx9xNdZn8J55iwz/P2CLdLOF8KVv2Zw9IQZPjtmOb5xrILErKDlUcnlheuvv17/+te/1LJlS7Vu3Vp79uzRgQMHJJ0ZdF6SGjdurK5du+qjjz4KZlNNKVSTZG9nlooypL4XJOnr3ceqrBCyWm06XnBaFz/3mU/ta3luLR05efqsbfS7IEmf//Cz0+qkir/gR2I1kzszDXkzQ2akVTOFyi+Ikfa+w33BOL75Wh1hts+z2dpbHapbnIvE8wVfhfr+EajvcbN+dkL578exCpJvs4KG6vW3twi5vNCgQQO1bdtW69atU0xMjCTp+PHj2rx5s7Zs2aLNmzdr06ZN2rNnj0pLS4PcWvMJ5Z2s4kCdCTFRSm1+jjbt/0VlVZwU2L8Dq8pPogypW6v6+i4nz6fxvrxRcdsJMVG6qHldbd5/osrXYdYyZ3dOSNw9Mbiua3M9N+Kiak9uqjt5c9Umsw6OXFNBoadt8OdJcyifzMJzwQhhrVabOk5c5fY4N+W7tJv1IjCc+PL3q6n2hMoxKVR+1DCLQO7Pvn5OqjrP9eeMpHx2/M/TY1XmxIEqtlpD4lgC//LleyuUr7+9QcjlhXPPPVd33323XnjhBZfLnTp1SnXq1PHLNjdu3KgVK1YoPT1dWVlZOnLkiGJiYtSsWTNdfvnluvPOO5WWluaXbdnNmzdP//znP7Vt2zadOHFCjRs3Vu/evXX//ffr0ksv9eu2yjPDTlZVBVTFk4LL2jbQmp1HqgyNzCjaYujf91+uNo1qB+2L0dnJm6cVafblYy0WpU7+xKOQ0dXJqauTt6qCRft6JPnlBDhYF8bB/gXRnyfNhAvhJ1ghrLe/qHIRGBp8+UXcn0LtmBQKP2qYiVmro2oiVOWzUzM8PVbZJ9MK9rEENYMxuc4g5PLCtddeq8TERM2bNy8g27viiiv01VdfVbvc2LFj9dZbbyk2Ntan7RUWFmrUqFFasWJFlY9bLBZNmDBBEydO9Gk7zph5Jyt/UvDYh1t9Gqw+FEVZDJVZbVV+MbobQHnD2clbVd0xXVWkVQya4qMtKiq1etWmiien7py8VcVVtZ99G0M7Nav2PXT3RNrTv0d1ywe72sGfJ82EC+GpJkPY6iozPd03dhw6yUVgiAjWsa38Z2rZttyQOyYF+0cNMwlUqGOW7y4+OzXDk2NVVULpMxJKQql61hPeHnfMfP1dFf//5BQBnnjiCV1zzTX66aef1KpVqxrfXm5uriSpWbNmuv7669W7d2+1bNlSZWVl+uabbzRt2jTl5OTonXfeUUlJid5//32ftnfHHXc4Aq4rr7xSDz/8sJo1a6bMzExNmTJFu3fv1qRJk9S0aVPdfffdPr++cGIfwNeXgeRDmT04sg/KvjQjV48OaK9dR065FUC5022v4v1VnbzZt1/xZKmwpEwb9v3ivP02nfW4twGX9Nsgnm0b1VGbRrU166s9HgdcUtXhVvltPDI/Q098uE2nS61Of3VzZ8DRRz/I0NKMXKdjw1Xk7q/CRaVlbp9YFZaUqai0zO1qB3dOMGalV/++l1ptmp2+1+VJczCmXPbmBMqsJ10VBep1eDoI69RRndxqjzv7h8ViaHBqE7cu6oakNpXFYvjt82wmofqZ9ubv54uKn6m4aIuKS61y9mkIxjTwNbU/hatA7M/B+O7yBp+dmuPJsaoqofIZCRWhVj3rqZRmdTVtdOdqg28zvBZfUMnlhTfffFOfffaZtm3bpgULFig1NbVGt3fttddq7NixGjlypKKioio9fvToUV1++eXauXOnJGnt2rW64oorvNrW559/rquuukqSNHToUC1evPisbR49elTdu3fX/v37Va9ePe3Zs0f169f3alvOhEOS7GnpsCvVneiajbNue1WFYpe2baC1Id7l017dFkgVf3Vz99dRd9YlefarcE1UO7h7guHPbQfyF2ZvTqBC+aTLk5DClzHrvFETXc482T88+UW1Q5PEkBoDqqaF8mfaLhQqcaoTyKoXf+xPoRpq+lugKgH9Pa5oTQmV7r/hytveBOVRQWeeqkh3eDq+Xjhcf5dHyOUFi8UiwzBks9lkGIYuueQS9e/fX926dVP37t0dH5BAWr58uYYOHSpJevDBB/Xqq696tZ4hQ4Zo5cqVio6O1t69e6v8gM+fP19jxoyRJL300kt6/PHHvW94FcJhJ/O1dNhu45+v0rm14lx2WUBk8ubCuLp1pTSr6/ZFXfnx2dztmuvOCZQn3S49nZXU2UlzIAdt9eYEKlRPujwNKbwds86X8MAfF5rlL8q96U7ozt9vaKdmfvs8m0GofqarUtNt9fXiNJCBJ5MpuC8QoY4355rBet+DPbRBJPAlLJd438N1zDh3f1gIh+vv8sx5dhRk77//vrZu3aqMjAxlZGTo22+/1bfffivDOPPBadCggbp27aru3btrypQpAWnTlVde6bi9e/dur9Zx8uRJffbZmRPs/v37O/1wX3fddapbt67y8vK0ePFiv4dc4cDX0mHpzJfNubXiZLEYGt6ludolJVZK5Bslxmn/8QI/thxmYe/i8L+/6+hzmFq+u4S73SuG/22dY3y2S9s2qLaiLdpi6M601i7X6023S3clxEQpPrpyJazkeZdL+2QFnl4seNOtJFS7orjqRrw0I7fSWHI7Dp10+ToqdiWuuC5nIUJ1J2++dDmr6qK8UWKsx92PnB2/h6Q21ZUXNNLnP/ys//lXpt8+z6EuVD/Tzrj6+/ljxjl3jrmueNoN3BV3Zv8dfFETLdri2f7kzvEiVEJNf4mPjnIM8F0db/dnT7677IL1vge6+28kqupY5cnYs/48lphRuA4XYB9KJ9JE3iv2gxtvvFE33nij4/8PHz7sCLzs/z777DOtXr06YCHX6dOnHber6tLojg0bNqi4uFiS1KdPH6fLxcbGqlevXvrkk0+0YcMGlZSUKCYmxqttVqvguJQfXzPrrmH3dD9H6Rk7vD55vfbCZrIUHnP8f8o50rRrmmvq4N8uHH/8+ZRumf0fKrz8KD4mSkU+hkaB8vW2HSrr10jNYvJ9bvM3mT+odFBTfZu5U+fKjXXZz5lKpIwdv6i+IclwPoD+5OEdlXJOsZR/1Okq5635TnWtv7rerk3asiNPCZISqm+lQ8X9qbx4q82z97Dkv9sukdZs/kXpGTs0eXhHDb6oaZWL2y8O532xo/rXZ5Xmr9miSUNTPH7Os8M7utd+H/1w6KSeX/Af1XV23LFKz87/Un/58ExX6/iYKDWoE6O61iLPN2aVnl/wlTokXqILmiSe1Yb3/vOTPvv+ZxWVlCk+JkpXXZikmy9pddZyknvH4miLobu713V8Pld+d1ATl2xXqdX222etRDp1XDrXjWZ/k/mDrIObOS7Uqjp+r8o6pIkffHn2Ntzk6vMc6tzazwP8ma5OVX+/M39b18e06litNvePuU7Ex0Qp/vQvUon3oYCr/UnSWY/FRlvUQK6HUCi/P7lzvKhqHzc7i6RRF8bro20Hq13W2/3Z4++u8oLwvntzLIZnKh6rYi0WXfHyGrc+I/44lpiVJ8fiit/vYaPgeLBb4Fd0V6whhYWF2rZtmy655JKAbG/x4sW67rrrJJ0ZGP/FF1/0eB0zZszQgw8+6Fjf7373O6fLPvzww44ukdu3b1dKSorb28nOznb5+MGDB9WzZ09J0oE/1lFyXYvb6wYAAAAAAO7JzrOqxfRTkuiuCBcSEhICFnBZrVa98MILjv8fPXq0V+spHz5V98EuP+7YgQMHPAq5gjFmGQAAAAAACG+EXFXIzc3Vpk2btHv3bh08eFAlJSWqVauWWrZsqe7du6tr166yWEKnumj69Olav369pDPjZXXv3t2r9Zw8edJxu06dOi6XrV27tuP2qVOnvNoeAAAAAACAvxByVdC6dWvt37/f5TJ16tTRsGHDdNttt+mqq64KUMuqtnbtWv3P//yPJCkpKUlvvPGG1+sqKvptvJTY2FiXy8bFxTluFxYWerSdAwcOuHy8fHdFAAAAAAAAdxByVfDTTz9Vu8zJkyf1/vvv6/3339ell16q2bNn64ILLghA6862fft2jRgxQqWlpYqPj9fChQuVlJTk9fri438b4N0+AL0z5Qe6T0jwZLjc6rtCnuXedVLzZh6tH/73w6GTmvuf/Vr9/WHHoLT9L2ys8xrU0t/X7q7xge+jDOnydg21fu8vZ20/7fwGSt91rFK7fn9JS13QJLHSTFHOXod9eVfs69p/rEBj/7G+2oFT372z6sFc3Zm9yjEo9X8HvvZUXLRFnz56hWrFRMtiMc4aRLumxcdEKf2JK896bc4GNe59fkN9tevoWQMal5S6HtDY3e1PXLbdrQF/qxJtMfT2HT3VqkEtjwZtjYu26JNHemvQ/6WbZvICd5X/TD+z5Duv39tAsk8i4WpQ+oqsVpvSXvqiRv5+9vewXVIdt7dRcX9ytS/bJ3hwNgFCqPjh0MlqJ0xxdQwNR+78XQelNHFrGnh3FRSXKu3FL3xejySlP3ml09m73D1eXNup2VkTDZjts+7JZBg1yX4e8ZcVO/RRpufvuz+27elnlGNCYHn7d4oEZjvu+FVOrjTdPLNGVoeQq4KOHTvq4osvVseOHZWcnKyGDRtKkvLz87V//35lZGRo7dq12rNnjyTpm2++Uffu3bV8+XL17ds3YO3cu3evBg4cqF9++UVRUVGaP3++rrjiCp/WmZj42xdHdV0Q8/PzHber69rok1rnSrUb1tz64ZYL2jbUs21ba1IVX4y9Ui+ocmr1Ky9opC9+OHLW/anNz9Gm/b+ozIOwJdpiOKa5ruqLeWBPVdku6czsRrXcfB3Vsa+rQ33p6dFxlaZDr9jeC9pWPS13xTY5e2zwJY3UqkWrs95bt5VKljqNZPnvRcfgSxp6v67/ijIMlbkxT8nI1GRZEhs5/v/M1PHf//e9+m835xLp3W35eneb/ThS29FuX41MTZZqN9SH3xepUHW9W4lVmrX5pKaNbiNJ6pXa3q1pz6OshlJf3iLH6wknVunNTXmaNrq1xvTtpnnfpYf+rK4lv/333W35mvfd945jSVXsx5deF7XXoi3V/71bnVtLP588rcKSMvf2j3LvobufqfL7U1Zunh5cckClVicXeVbpwSUHtLRFK6U08/KzHwAXtG2op0fHen0MDUdVHaPt36V3prV2/D2dfXd4Iz7BpsKY+l59H5SXEBOl+HMaS06+S905XkRbDN3Yt6tU+8zrNNtn3dX3XHXHHX+zn0eMubKW5m337H33VlZunmal79HKzEOOz+7g1CYal9bGrb/PzE3Z+tnZ39qu3PETvnF1Hhrp3D0Wh6VaXsyAHcKYXdFL27Zt06uvvqq3335bZWVlqlu3rnbs2KEmTZrU+LZzc3PVu3dv7dmzR4ZhaM6cORo7dqzP6w3U7IrVyc7OdgxOHw6zO0QKZ78MVbw/KzfP7VAslL9UnL2Ommiv1WpTQXGpejz/mVsXJAkxUdo+eVCVAZ7VatPxgtO6+LnP3N7+xj9fpcO/ntbwv62r9oR56QNpjteflZunYTO8C0MMSbHRFp0utSohJkqXtW2gNTuPuAxI7ds/r2EtpUxY5fE2yyv/HvryOtwRH21RidXmUfgbDOXfkzMXdVWHFNUpHw55G357q+JnVKp8gRYXbVFxNRWF9vV0aJLo1b6549DJaj9TFdv66IIMt4Kx67o213MjLgr5X+kDeQw1k0BWWbj7mXJlZLdkTRvt+td/V8eL8j9kedoud7Zd09z5fqjquBMInr7vwdiG1WpTx4mrfD63CSYqo8JTpP1dw+36m0ouL3Xq1EmzZs3SrbfeqmuuuUYnT57U1KlTNW3atBrd7tGjRzVgwABHJdlrr73ml4BL0llB1Y4dO1wua388Ojpa7dq188v2YW4Wi1Fld4WK96c0q6tpoztr6qhOlb48ru3crMr7Q5Gr1+FvFouhOvExGpzaxK0T/yGpTZ22xWIxdG6tOCXERLl9UnlurTg1rBOvaaM7V3syW/4kflb6Hq+DIZuka1KbnnWx7s7JdEqzurJabW6/PmcKS8pUVFqmWrHRjr+1t6GOM+XDiMc+3Op2gPGnazp4FFL6S/n3ZHiX5mqXlHhWSOFuOPTGzd3VoUmiy/C7ppRabZqdvtdxYVzVZ+p0qdXlOip+1i0Ww+0229/D6j5TFbdhtdq0MvOQW9tYtCVHi7bkeFxNEWiBPIaaibPv0powLq2Nlmbken1ci7YYujOt+sqaqo4XzkJNTz7rKzIPauqoTkH93LjzPVfxuBMonrzv3sjKzXP5vVhqtWn8gq1ql5TodFtFpWUeHz8DtX9Ux9cKNoS2QB6LPRFp4Zu3qOTyg+nTp2v8+PFq166dfvjhhxrbzq+//qp+/fpp8+bNkqQXXnhBTz75pN/Wf/LkSTVs2FDFxcW6+uqrtXLlyiqXKy4uVqNGjZSXl6dLL71UX3/9td/aIIVfkgz4iz9/Mfb2l3J3qy88+XXWmap+tXV3+75WKLizbXe7cFbFm6q38tVDvr633nD2K3r5E65l23J9/lW/qLRMf178nVtdBr3hSTVVxYpCXz/rFd9Ddz/PBcWlXlcn+qtiA+HJ1Y8H9l29ql3E289VdRdonn7Ws54dFLQLUTNVIdXEhbE/Ku7M9B6WF4gqOaC8mg5Vw+36O/TiSRPq3bu3pOpnDfRFQUGBrrnmGkfA9fTTT/s14JLOjMl11VVXaeXKlVq9erWys7Or/IAvWrRIeXl5kqQRI0b4tQ0AnPO0+sMVd37Br+pXenerLzz5ddaZqn61dXf7vlYoVFUNV37bnnRRq6iqv5Onf1t3q/r8yVmFYPlfO32tHLCva1zvNlq61fXfr3wAFR9tUVE1FVh29s+VOxUYVVUUVtVmb6ss3f08x0dHeV2d6E41BSJXdfusJL9WAlVXHeHJZz0hJkrx0VEet8FfzFSF5O+qFH9V3Ply/AwWf1SwAZ6oKlQtLCnTos05WpqRS6haBUuwGxAOVqxYIcnzWQbdVVxcrBEjRmjdunWSzoyH9dxzz3m8njlz5sgwDBmGoUmTJlW5zGOPPSZJKi0t1f3336+ysrO/vI8ePeoI1+rVq6dx48Z53A4A3hvepbmWPpCmkd2SlRBz5uQ+ISZKI7sla+kDaW5/ydkvrqOdnDBWF5jZT5idnXDaL1R84eoCprrtV/f6XKmuC47FYnjURc2uur+TJ3/bcWltvHpt3nK3W5L023u/ffIgZT07SNsnD3I7fK24Dlefz1du7KLvn71aWc8O0neTBrn9eUuIiVKsxeL2BdrK7w5VW/3gzt/D1XtY3efZfiHorVKrTbO+2qOC4lJZQ3zsNwSeq33WH/uzJzz5rAc79PDke87V95l93E0z7ZveBHzO+Hr8DDRPuqgCvnI3VM3KzQtwy0IblVxuuuWWW7Rx40a1b99ezZs3V+3atXXq1Clt2LBBW7ZskWEYGjRoUI1se8yYMfrkk08kSf369dOdd96p7777zunysbGxat++vVfb6tevn2688UbNnz9fS5cu1YABA/TII4+oWbNmyszM1PPPP6/9+/dLkl588UXVr1/fq+0A8J6/xrKpyfE6PPl11hlfL2CcvT5XA527Ww3nSbVBfLRFG//c32WIYefu39aXccKiDOnKDklat+uYR++JfZB1dz9vvlYOuPv5tG/Dk2qAYqvVrxUY/qyydMbX6kSzjNWF4HG1z4baWGGhEHr4WoVk5jGd/FlxF4jjp7+Ybcw4mJ+7oeqsr/aYYtKZQGFMLjeNGTNGH3zwgQzj7A+NzWaTYRjq37+/5s6dq4YNG/p92xW3WZ1WrVpp3759le6fM2eObr/9dknSxIkTnVZzFRYWatSoUY4KtYosFoueeeYZp8/3Vbj1CQbMoCbG6/BlVkJ/z0bl7iyfnoR7oTADWFWvw53AaniX5h7NfPr5Dz8H9ULMnc9nTY1r5sk4MDU9Y6Avs1pWxLgxCGVmGfPI27EyzfL6XPH3d6AZZlw105hxMD9vxrf19hwt3K6/CbnctHHjRn388cdat26dvvrqKxUUFMgwDPXu3VtTp05Vjx49amzbgQ657N5//33NmTNHW7du1YkTJ9S4cWP17t1bDzzwgC699FKP2uSJcNvJgEjmzUV5IE/wfQn3QmnqeH+GeP4cSD7QPLlwrMmQsiZnP/LnTJSB+nwGAjNOhR8zhB6S54FVKH13+KKmXkco78tmHSgf5hTISWfC7fqbkMsLp0+f1rJly/S3v/1Na9euVXx8vGbMmKE77rgj2E0LC+G2kwGRzlWF0Bc/HAn5CxhXQv3X+HAJ8Tzh7oWxWV+fnb9moqzJSsNAMHOXL/zG1bEqlEMPO08CuVCoAvaXUP8OrAnh9PdDaPN1pnJPzmHC7fqbkKsCe/dDdy1YsEC33367ioqK9Omnn6pfv3412LrIEG47GYAznF2omOECxhWzVBt4yuwn8u58rsLhAs2XbsGSuasNwuHvF+kCGVIG4rumum2EYyVQuH4HOmP2H0hgLu6eiznj7jlauF1/E3JV0KVLF7322mvq3bu328957bXX9PDDD6tfv35avXp1DbYuMoTbTgYgMpg9rCsvHC/EnAmHCzRfx+oy47gxXGiaX6BCylCq9gvnMZ3C6TuwOgTsCJRA/ZAVbtfflmA3INRs27ZNffv21fDhw7Vp0ya3nmMfn2r9+vU12TQAQAizzz4WDif3/pwePtTZZ/baPnmQsp4dpO2TB4XMTF7uGt6luZY+kKaR3ZKVEON8FrOqVDfzWahyd8ap2el7A9QieCIrN89lMFtqtWn8gq3Kys3zaTtLMnI0bEa6Fm3OcRzTCkvKtGjzmfuXZHhfIeEN+6yE7jDbvhlO34HVqeqYmxATpZHdkrX0gTQCLviN/Rwl2sv9yuznaN4i5KrgwgsvlM1m0/Lly9WzZ0/16dNH77zzjvLynH/J2mchpCgOABAOwvlCzBmzX6BVDOuu6+reRdaQ1Kame81Wq00rMw+5teyKzIOy+mE2SvhXIELKQAVpnrBYDA1ObeLWsmbcNyNJOPxAAnOIxB+yfEXIVcHWrVv1/PPPKyEhQTabTenp6br99tvVsGFDXXLJJbr33ns1depUvfXWW5o+fbpGjBihZ599VoZhqEOHDsFuPgAAPuNCzLzsYd243m2q/eU32mLozrTWAWqZ/0RSpWE4ClRIGarVfuPSwnffjERm/4EE5hBJP2T5AyFXBdHR0Xrqqae0e/duPfDAA4qNjZXNZlNpaak2btyot956S//zP/+je++9V4899piWLl0qq9UqSfrDH/4Q5NYDAOAfXIiZW3VdHOzjxpix6iASKw3DSSBCSn8FaVarTQXFpX6tBgznfRNAzYqEH7L8gZDLicaNG+vVV19VTk6Opk6dqp49e8pischms1X6FxcXp4kTJ+qOO+4IdrMBAPALLsTML1zHjaHS0NwCEVL6GqRl5ebp0QUZ6jhxlVImrFLHiav06IIMv3VtDNd9E0BgcI7mGrMreuDkyZPavn27du3apV9//VWS1KJFC/Xp00fnnHNOkFsXPsJtdgcAMLNwmH0Q4TfzGbMrmtujCzK0aHP1g76P7JasaaM7e7x+X2aIDfTMeeG2bwL+wr5RPX+do4Xb9TchF0JOuO1kABAOONlEqAl0GAH/CURI6U2QRngKBF9Wbp5mpe/RysxDjuBmcGoTjUtrExL7XSieD/napnC7/o7Y7oovvPCCVq9erRMnTgS7KQAAhDwG10WoocuXefmjq01142V5M65gqA5WD0SKJRk5GjYjXYs25zgqMQtLyrRo85n7l2RUH1zXlJruxuwLztHOFrGVXBaLRYZhyDAMlZaWOu5fsGCBunTponbt2skw+JAEQ7glyQAAoGaF4i/rqJ43XW08qfLwpNrPly6OAHwXypWU4V45HG7X3xEbckVHR8tqtcowDJWV/fZlZg+/EhISlJqaqs6dO6tLly7q3LmzOnXqpNq1awex1ZEh3HYyAAAAOOduSOnNhaa7QVpBcalSJqxyu81Zzw5Srdhot5cH4FpNj9XnrVAO3/wl3K6/I/bIfOrUKW3ZskUbN26s9JjNZlNBQYH+85//aP369Y77DcNQmzZtHKGX/b9m/xAAAAAAwWLvauNKVm6e04BLOtONcPyCrWqXlHjWhaa9a+TUUZ1cBmn2WR/dreTyZtZHAFWzWm1amXnIrWVXZB7U1FGdPKqk9KXa15NuzIEM3+BcxIZc8fHxuvTSS3XppZeedf/333+vrVu3Ov5lZGQoNzdX0pnwa9euXdq9e7f+9a9/OZ5Tv359de7cWV27dtXLL78c0NcBAAAAhDtfLzSrC9IsFkODU5u4VUkyJLUpXRUBPyoqLXMrYJbOjNFVVFrmViVldd2bqwu/ajp8Q82I2JDLmQsuuEAXXHCBRo8e7bjv2LFjZ4VeW7du1ffff6+SkhJJ0vHjx/XFF19ozZo1hFwAAACAHwXqQnNcWhstzcittltS+cHqAfiuJiopq+rebB/EfsmWHHVrVV/f5eS5HNuvpsI31Cz+Am5o0KCB+vXrp379+jnuKykp0ffff+8Ivez/AAAAAPhPoC407V0bqxv3y6zj7gChyt+VlNV1by6zSRv2/eL4f3v4tTQj96yx/SKlG/Oq7QeD3QS/IuTyUkxMjDp16qROnToFuykAAABA2ArkhebwLs3VLinR41kfAfjGn5WU7nRvrkrFsf0ioRtzVm6e/nfZ98Fuhl9Zgt0AAAAAAHDGfqHpDn9caNorurZPHqSsZwdp++RBVHABNcy+30U72X/draT0pHtzVexj+9mNS2vjtE3l22bWbszeBoKhjJALAAAAQEgLxoWmfbB6M1ZnAGY0vEtzLX0gTSO7JSsh5kxFZkJMlEZ2S9bSB9Ic3QidsVptOl5w2u3uzc6syDwo63+DH3+Fb6HI10AwVEVsd8VFixbpuuuuq7H15+bmav/+/erVq1eNbQMAAACIBIyXBUQG+74+dVQnlzMflldxFkVfVRzbL1y7MXsy3qGZRGzINWrUKHXq1El//vOfNWrUKL+t98CBA/rLX/6if/7zn3rqqacIuQAAAAA/CNcLTQCV2Sspq1PVLIq+qmpsP2/Ct1BnH+/wZLAb4mcRG3K1bdtW27Zt0w033KCWLVvqpptu0k033aSOHTt6vK78/HwtXrxY77//vlavXq3S0lJFR0erbdu2NdByAAAAIDKF44UmAO9UN4uit1yN7edu+GYG9vEOF6w5HOym+FV4/HW8kJWVpVdeeUUvvfSSfvrpJ73wwgt64YUX1K5dO/Xq1Us9evRQ165dlZSUpPr166t+/foqLCzU8ePH9csvv2jnzp3asGGD1q9fr/Xr16uoqEg225md67rrrtOUKVPUvn37IL9KAAAAIPyE04UmAO/UxKDpZh5E3hvj0tpo0Zfbgt0MvzJs9mQmQp06dUqvv/66/va3v+nAgQOSJMNw/9cg+9sXFxen6667Tg8//LB69uxZI22NFNnZ2WrRooWkM90/k5OTg9wiAAAAAECosFpt6jhxlUdjSiXERCm1+TnatP8XlbkY26+6Ae7DzexVGzTu6jMZRjhcf0f8zx916tTRE088occee0yffvqpFixYoC+++EL79u2r9rnx8fG65JJLNHz4cI0dO1bnnntuzTcYAAAAAIAI5umg6Rv/fJXOrRUni8VQVm4eY/uVM6hj02A3wa8iPuSys1gsGjRokAYNGiRJysnJ0ddff63s7GwdOXJEx48fV3x8vBo1aqRGjRopNTVVF198sWJiYoLccgAAAAAAIod90HR3gq6EmChHwCUxtl+4I+Ryonnz5rr++uuD3QwAAAAAAFCOfdD0RZtzql3W2UDyjO0XnizBbgAAAAAAAIAnxqW1UXQ1FViRNpA8CLkAAAAAAIDJ2LsdOgu67APJR+I4W5Esomvznn32WUnSfffdp4YNG3q9nn379umOO+6QYRj67LPP/NU8AAAAAADgxPAuzdUuKZGB5OFg2Gy2ynNnRgiLxSLDMJSZmamUlBSv17N9+3alpqbKMAyVlbk/wwOqlp2drRYtWkgKjylMAQAAAAA1y2q1MZC8F8Lt+juiK7kAAAAAAID5MZA8JMbkqtaePXuozgIAAAAAAAhxhFwufPTRR+rSpYsuu+wy7dq1K9jNAQAAAAAAgBOEXE68+uqr+t3vfqf8/Hxt3LhRXbt21ezZs4PdLAAAAAAAAFSBkKsCm82mBx98UH/84x9VVlamxMRESVJ+fr7uvvtujRo1Sr/88kuQWwkAAAAAAIDyCLnKyc/P19ChQ/X666/LZrOpY8eOysjI0IoVK9S4cWPZbDYtXrxYnTp10hdffBHs5gIAAAAAEDGsVpsKiktltdqC3RSEKKYe+K/s7Gxde+21yszMlM1m08CBA7VgwQLVrVtX5513njIzM3Xrrbdq5cqVysnJ0YABAzR+/Hg9//zzwW46AAAAAABhKys3T7PS92hl5iEVlpQpISZKg1ObaFxaG6U0qxvs5iGEUMklacOGDbrkkkscAdd9992njz76SHXr/razNGzYUB999JFeeeUVxcfHy2q16uWXX1avXr20Y8eOILYeAAAAAIDwtCQjR8NmpGvR5hwVlpRJkgpLyrRo85n7l2TkBLmFCCWGzWaL2Do/i8UiwzBkGIasVqssFov++te/6qGHHnL5vMzMTI0ZM0ZZWVkyDEMWi0VlZWUyDENlZWUBan34ys7OVosWLSRJBw4cUHJycpBbBAAAAAAItKzcPA2bka5SF90Toy2Glj6QRkWXl8Lt+ptKLklWq1WJiYlaunRptQGXJKWmpmrTpk36wx/+IJvNJqvVGoBWAgAAAAAQOWal73EZcElSqdWm2el7A9QihLqID7lsNptatmypdevWaciQIW4/Ly4uTn/729+0ZMkSNWzYsAZbCAAAAABAZLFabVqZecitZVdkHmQwekiK8IHn7TMkXnTRRWrQoIFX6xg6dKi+++47ZWVl+bNpAAAAAABErKLSMscYXNUpLClTUWmZasVGdMQBRXjI1adPH7+sp1GjRn5bFwAAAAAAkS4+OkoJMVFuBV0JMVGKj44KQKsQ6iK+uyIAAAAAAAgtFouhwalN3Fp2SGpTWSxGDbcIZkDIBQAAAAAAQs64tDaKria8irYYujOtdYBahFAX0d0V3bF792598803OnTokAoKCnTfffcx0LwJFBcX69SpU8rPz1dxcTEzYAJhwmKxKDY2VrVr11adOnUUGxsb7CYBAEzOarWpqLRM8dFRVIIAISalWV1NG91Z4xdsrXKWxWiLoWmjOyulWd0gtA6hiJDLic2bN+uRRx7RunXrzrp/1KhRZ4Vcf/vb3zR58mSdc845ysrKUkxMTKCbinJsNpuOHj2qo0ePBrspAGqIPcQ+fPiwGjVqpAYNGsgwuCgBAHgmKzdPs9L3aGXmIRWWlCkhJkqDU5toXFobLpiBEDK8S3O1S0rU7PS9WpF50LG/DkltqjvTWrO/eqF8uB9uDJvNxjybFSxfvlzXX3+9iouLVf7tMQxDmZmZSklJcdx38uRJNWvWTAUFBfrwww81YsSIYDQ5rGRnZ6tFixaSpAMHDig5Odnt5+bm5urXX3896z7DMBQVFX47LxCJysrKVPFr65xzzlGzZs2C1CIAgBktyciptjJkeJfmQWgZAFeovPRNVeF+WjNDs+67WpLn19+hiEquCg4ePKgxY8bo9OnT6tixo15++WWlpaUpMTGxyuUTExM1bNgwzZ8/XytXriTkCqKioqKzAq4GDRqobt26iouLo8oDCBM2m02nT59WXl6ejh07Jkn69ddf1aBBA8XFxQW5dQAAM8jKzXMacElSqdWm8Qu2ql1SIhUiQIixWAzVig29GMMM4VtV4X5hSZlWZoZXL6jQ+3QE2fTp05Wfn69WrVrpq6++Ur169ap9Tt++fTVv3jxt2rSp5hsIp06cOOG4nZSUpAYNGgSvMQBqhGEYio+PV3x8vKKiovTzzz9Lkn755Rc1aeLe7DsAgMg2K32P04DLrtRq0+z0vZo2unOAWgXAjMzS7bm6cD+cMLtiBR9//LEMw9D48ePdCrgkqUOHDpKkvXv31mDLUJ2CggLHbXf/dgDMq/x+Xn7/BwDAGavVppWZh9xadkXmQVkj4IIQgHeWZORo2Ix0Ldqco8KSMklnKqMWbT5z/5KMnCC38DfuhPvhgpCrgp9++kmS1LNnT7efU7fumYT21KlTNdImuKes7MyBJTo6mjG4gAgQFRXl2Nft+z8AAK4UlZY5LkarU1hSpqJSvl8AVOZut+es3LwAt6wyT8L9cEDIVUFpaakkyWq1uv0c+zhQderUqZE2AQCqxnh7AABPxEdHKSHGvR9DE2KiwnLmMQC+86Tbc7B5Eu6HA0KuCuxjuuzZs8ft56xfv16S1LJlyxppEwAAAADfWSyGBqe6N4bjkNSmITuANIDgMVu3Z0/C/XBAyFVB7969ZbPZtHDhQreWLy4u1syZM2UYhvr27VuzjQMAAADgk3FpbRRdTXgVbTF0Z1rrALUIgJmYrduzJ+F+OCDkquC2226TJC1dulSffvqpy2WLi4s1duxY7d69W4Zh6K677gpACwEAAAB4K6VZXU0b3dlp0BVtMTRtdOeQmhkNQOgwY7dnd8L9cEHIVUHfvn11ww03yGazaejQoXryyScd3RElad++ffr66681depUdezYUQsXLpRhGLr33nvVsWPHILYcAAAAgDuGd2mupQ+kaWS3ZMfFakJMlEZ2S9bSB9I0vEvzILcQQKgyY7fn6sL9cGLYbLbImEfSA6dPn9bIkSO1YsUKl4Ma29+66667Th988AEz+vlJdna2WrRoIUk6cOCAkpOT3Xrejz/+qNLSUkVHR6tdu3Y12UQAIYL9HgDgK6vVpqLSMsVHR4XExSiA0JeVm6dhM9JdDj4fbTG09IG0kKoKzcrN0+z0vVqReVCFJWVKiIlS72YWvXXfIEmeXX+HKiq5qhAXF6fly5dr5syZatOmjWw2W5X/kpOT9frrr+vDDz8k4AJ0pgvvvHnzNHbsWHXo0EENGjRQTEyMGjZsqO7du+sPf/iDVq9e7Zi9tG/fvjIMw6d/c+bMcWzffp+74+PNmTOnyvXcf//9jvvnzp3r1rp27dqlWrVqyTAMde3a1TFT6759+1y2v3bt2mrdurWuv/56/fvf/3ZrWwAAwH8sFkO1YqMJuAC4zazdnu3t3j55kLKeHaTtkwdpwtCUYDfLr6KD3YBQdtddd+muu+5SVlaWNm7cqJ9//lllZWVq0KCBunbtqm7dujF9PfBfixYt0vjx47Vv375Kjx07dkzHjh3T5s2b9fe//13t27fXX//618A30k0vvPCCli9frv379+uRRx7RwIED1ahRI6fL22w23XXXXSosLFR0dLT+8Y9/KDravcNrQUGB9u3bp3379unDDz/UwIEDtWjRItWuXdtfLwcAAACAnw3v0lztkhIrVUYNSW2qO9Nah1zAVZ493A9H4fmq/CwlJUUpKeGVbgL+9L//+7+aMGGC4/8HDBigYcOGKSUlRfXq1dPx48f1ww8/aNmyZfr000+1c+dOPf3001q8eLHy8/OrXOef//xnLVmyRJK0atUqNWvWrMrlaqKcNjExUTNnztTgwYN19OhRPfTQQ5o3b57T5WfNmqU1a9ZIkh577DF17dq1yuWGDx+u55577qz7fv31V23atEmvvPKK9u7dq08++UT33HOP3nvvPb+9HgAAAAD+Z6+MmjqqE92eQwQhlx9t2rRJ3bt3r5F1//zzz1q/fr3Wr1+vDRs2aMOGDTp27Jgk6dZbbz2rq5UvJk2apMmTJ7u17BdffOF2tzCEr3/+85+OgCspKUkLFixQnz59Ki3Xv39/3X///fruu+/0xz/+UUeOHFHr1s6n5q5Xr57jdvv27XXeeef5u+kuXX311Ro7dqzeeecdzZ8/X7///e917bXXVlouNzdXjz/+uCTpggsu0MSJE52us169errooosq3X/55ZfrhhtuUGpqqo4cOaL3339fL730ktNgDwAAAEDoCOfKKLNhTC4/+PrrrzV48GBdcsklNbaNxo0ba+jQofrf//1fffzxx46ACwimnJwcPfDAA5Kk2rVra+3atVUGXOVddNFFWrVqlR577LFANNEn06dPV+PGjSVJ9957r/Ly8iotc9999+nXX3+VYRiaNWuW4uPjvdpW48aNNXbsWElnuj9u3LjR+4YDAAAAQAQiavTBZ599pueee05ffvllQLfbsmVLdejQQZ988kmNbiczM9Pl466qcHC2cJ21Z/r06SooKJAkPfvss+rQoYNbz7NYLLr55ptrsml+ce6552rGjBm6/vrrlZOTo8cff1wzZ850PL5w4UJHl8r7779faWlpPm2v/D51+vRpn9YFAAAAAJGGkEtnqiYWL16s1atX68CBA4qJidF5552nUaNG6bLLLqu0/Jo1a/SnP/1J//nPfxzPl6SBAwfWWBsnTJigHj16qEePHmrcuLH27dtX4yFTVd2q4Jms3DzNSt+jlZmHHAMRDk5tonFpbUJ6IEJ32Gw2vf3225LOVHHdddddQW5RzRg1apRGjBihxYsX66233tKYMWPUt29f/fLLL3rwwQclSa1atdJf/vIXn7f1008/OW63bNnS5/UBAAAAQCSJ+JDrp59+0vDhw6usWvq///s/XX/99Zo7d66ioqJ07NgxjRs3TkuXLpV05iLfMAwNHz5cTz/9tC6++OIaa6e742QhdCzJyNH4BVtVarU57issKdOizTlampGraaM7a3iX5kFsoW+2b9+uo0ePSpJ69+6txMTEILeo5rz++utas2aNfvnlF40bN06ZmZn64x//qMOHD0uSZs6cqTp16vi0jZ9//lnvvPOOJKlFixbq1q2bz+0GAAAAgEgS0SFXcXGxrr32Wm3fvt3pMgsXLlTLli314IMPqk+fPvrpp59ks9kUFRWl0aNH609/+pM6duwYwFbDDLJy8yoFXOWVWm0av2Cr2iUlmraia+vWrY7bNTXhgrfy8/P13XffVbtcTk6OW+tr0qSJpk2bpjvuuEO7d+/W0KFD9dlnn0k6M/HDoEGD3FrPiRMnKrUrLy9Pmzdv1v/93//p8OHDio2N1WuvvaaYmBi31gkAAAAAOCOiQ665c+dq+/btMgxDrVq10p///GelpqYqNjZW33//vaZOnaotW7bojTfe0DfffKN9+/ZJkkaOHKkpU6aoXbt2wX0BCFmz0vc4DbjsSq02zU7fq2mjOweoVf5VfvKDpKSkILakso0bNyo1NdWv67z99ts1b948ffrpp46Aq0mTJpo+fbrb61iyZIljDK+q3HDDDXr88cdDLjQEAAAAADOI6JBr0aJFkqTk5GRt27btrO5GnTt31ujRo3XFFVfo66+/1rp16xQVFaXZs2c7ZkALdwMHDlRGRoZOnDihevXqKSUlRVdffbXuuece1a9f3+v1Zmdnu3z84MGDXq87FFitNq3MPOTWsisyD2rqqE6mHIz+5MmTjtu1a9cOYksC56233tJFF12kU6dOSZJmzJjh075Q0bJly5SQkKDp06erXr16flsvAAAAAEQCS7AbEExbt26VYRh6/PHHqxxPx2Kx6Nlnn5UkGYahW265JWICLkn69NNPdeTIEZWUlOjIkSNau3atnnrqKbVp08ZlNUp1WrRo4fJfz549/fgqAq+otEyFJWVuLVtYUqaiUveWDTXlx+DKz88PYksq69Onj2w2W7X//vnPf3q03latWp1VZTVy5EiPnn/rrbdWakNRUZF++OEHTZkyRYZhaM6cObr88ssd430BAAAAANwT0SGXvbuVq1kEO3Xq5Lg9atSoGm9TKEhNTdUzzzyjZcuWadOmTfr222/19ttvO2aPPHHihEaOHKmVK1cGuaWhKT46SgkxUW4tmxATpfho95YNNQ0aNHDcJpDxXlxcnNq3b6+nnnrKUV2alZWl8ePHB7llAAAAAGAuEd1dsbCwUIZhuBxPqGHDho7bycnJgWhWUD3yyCOaNGlSpfsvueQSjR07VjNnztS9996rsrIyjRs3Trt371Z8fLxH2zhw4IDLxw8ePGjqai6LxdDg1CZatLn6Qc2HpDY1ZVdF6UyXXrvNmzcHsSXhY+DAgercubO2bt2qBQsWaObMmRHTFRQAAAAAfBXRlVyeio4O/0ywunGA7rnnHt15552SpNzcXP3rX//yeBvJycku/zVt2tSbpoeUcWltFF1NeBVtMXRnWusAtcj/Onbs6AiBv/rqK+Xl5QW5ReGhQ4cOkqSSkhLt2LEjyK0BAAAAao7ValNBcams1UzaBbiLkAseu+eeexy3165dG8SWhK6UZnU1bXRnp0FXtMXQtNGdldKsboBb5j+GYejWW2+VdGZMrlmzZgW5ReGhtLS0ytsAAABAuMjKzdOjCzLUceIqpUxYpY4TV+nRBRnKyuWHc/gm/EuT3PD666+77LLoyXITJkzwV7NCVkpKiuN2Tk71XfIi1fAuzdUuKVGz0/dqReZBFZaUKSEmSkNSm+rOtNamDrjs/vjHP+qNN95QQUGBJkyYoCFDhjgqkVyxWq2aN2+efv/73wegleZhs9m0adMmx/+3aNEiiK0BAAAA/G9JRo7GL9iq0nLVW4UlZVq0OUdLM3I1bXRnDe/SPIgthJkRckl64403XD5uGIZby0mREXLZ3w9Uz17RNXVUJxWVlik+Osq0Y3BVpXnz5poxY4buuOMO5efnq0+fPlqwYIH69Onj9DlZWVl65JFH9PPPPxNyVfD6669r3759kqQuXbqoWbNmwW0QAAAA4EdZuXmVAq7ySq02jV+wVe2SEsOiKACBF/Ehl83mv76/kRL+ZGVlOW5zEe4ei8VQrdjw3N1uv/12ZWdna8KECfr555/Vt29fDRw4UMOHD9eFF16oevXq6fjx49q5c6c++ugjffzxxyorKztr4PpIceLECX333Xdn3VdcXKx9+/bpww8/1Lx58yRJFotFL774YjCaCAAAANSYWel7nAZcdqVWm2an79W00ZF3vQDfhedVt5u++OKLYDfBlGbOnOm47apiB5HjmWeeUceOHTV+/Hjt27dPn3zyiT755BOny3fs2FEvvfRSAFsYGpYsWaIlS5a4XKZOnTp64403NHDgwAC1CgAAAKh5VqtNKzMPubXsisyDmjqqU1j1gkFgRHTIFWkBzZw5c3T77bdLkiZOnKhJkyad9XhmZqYSEhJ0/vnnO13Hm2++6RhgvEmTJhoxYkSNtRfmct111+naa6/Vhx9+qJUrV2rDhg36+eefdfLkSdWtW1fnnXeeevXqpVGjRqlv374RU/lYnZiYGJ1zzjm68MILNWDAAN15551USAIAACDsFJWWqbCkzK1lC0vKVFRaFra9YVBz+MSYRHp6unbt2uX4/6NHjzpu79q1S3PmzDlr+dtuu83jbWzatEnjxo3TlVdeqcGDBys1NVUNGjRQaWmpduzYoblz5zqqc6KiovTmm2+qdu3aXr0ehKfY2FjddNNNuummm3xe15w5cyp9rqvjaffj2267zeN9Zc2aNR4tf9555/m1WzQAAABgRvHRUUqIiXIr6EqIiVJ8dFQAWoVwQ8hlErNmzdLbb79d5WPr1q3TunXrzrrPm5BLksrKyrR69WqtXr3a6TINGjTQ7NmzNXToUK+2AQAAAACILBaLocGpTbRoc061yw5JbUpXRXiFkAsOQ4YM0ezZs/XNN99oy5YtOnz4sI4dOyabzaZzzz1XnTt31tVXX63bbrtNdesy0wUAAAAAwH3j0tpoaUauy8Hnoy2G7kxrHcBWIZwYNvrRIMRkZ2erRYsWkqQDBw4oOTnZref9+OOPKi0tVXR0tNq1a1eTTQQQItjvAQAAzGVJRo7GL9haZdAVbTE0bXRnDe/SPAgti0zeXn+HKiq5AAAAAABAQAzv0lztkhI1O32vVmQeVGFJmRJiojQktanuTGutlGb0GoL3CLkAAAAAAEDApDSrq2mjO2vqqE4qKi1TfHQUY3DBLwi5AAAAAABAwFkshmrFEkvAfyzBbgAAAAAAAADgK0IuAAAAAAAAmB4hFwAAAAAAAEyPkAsAAAAAAACmR8gFAAAAAAAA0yPkAgAAAAAAgOkRcgEAAAAAAMD0CLkAAAAAAABgeoRcAAAAAAAAMD1CLgAAAAAAAJgeIRcAAAAAAABMj5ALAAAAAAAApkfIBQAAAAAAANMj5AIAAAAAAIDpEXIBAAAAAADA9Ai5ACACFRYWKjo6WoZhaMqUKcFuDgAAAAD4jJALgNfWrFkjwzA8+jdp0qRgN7vGhNL7kZOT49jGli1bKj2+ZcsWlZWVSZIuvvjiGmkDAAAAAAQSIRcAhKHly5dLkpKTk9W1a9dKj2/cuNFxm5ALAAAAQDiIDnYDAJjXRRddpJYtW2r//v0yDEOzZ89Wjx49XD4nKSkpQK0LvFB6P5YtWyZJuvbaa6t8/KGHHtJDDz1UI9sGAAAAgGAg5ALgtYYNG2r58uW6/PLLdfLkSf3pT3/Sf/7zH7Vs2TLYTXPLpEmTNHnyZPXp00dr1qzxeX2h8n4UFBTos88+kyQNHTo0oNsGAAAAgGChuyIAn6SmpmrevHmyWCw6dOiQrr32Wp08eTLYzQqaUHg/Vq9eraKiItWqVUv9+vUL6LYBAAAAIFgIuQD47JprrtFf//pXSVJmZqZuuOEGx6DmkSjY74d9PK7+/fsrPj6+0uM7duxwDEo/f/78So+XlJSoY8eOMgxDnTp1crvtqampMgxDrVq18u0FAAAAAIAXCLkA+MXDDz+sP/zhD5KklStXRvx4T8F6P2w2myPkctZVcevWrY7bXbp0qfT4yy+/rKysLEnSq6++qqioKLe23b17d0nS/v37lZ2d7UmzAQAAAMBnhFwA/ObVV1/VwIEDJUmvv/66XnnlleA2KMiC8X5s2rRJBw8elGEYuuaaa6pcxh5yJSQkqF27dmc9dujQIT333HOSzlSk9e3b1+1td+zY8ax2AAAAAEAgMfA8Io/VKhUeD3YrAifhXMkSmDw7OjpaCxYs0GWXXaasrCyNHz9ebdu2jdjBz4PxfthnVezevbuaNm1a5TL2kCs1NbVSldaUKVNUUFAgSZowYYJH205OTnbc3r17t0fPBQAAAABfEXIh8hQel6a2DXYrAufx3VLthgHb3DnnnKPly5erS5cuysvL05gxY5SRkaHzzz8/YG0IJYF+P6rrqihJGRkZkqTOnTufdf+xY8f05ptvSpLS0tLUs2dPj7bdsOFvn7ODBw969FwAAAAA8BXdFQH43c6dOx0zCnbo0CHiByIP1PuRk5OjzZs3S3Iech09elS5ubmSKo/HNWfOHJ0+fVqSdMcdd/jUlpKSEp+eDwAAAACeIuQC4FdHjhzRrbfeKpvNplq1amnu3LmKiYkJWnvsswhW9W/y5MmSpLVr17pcbs6cOV5v35v347333tM999yjiy++WHFxcW63wV7FlZycrK5du1a5TPlB5ytWci1YsECSFBUVpREjRlR67oYNG3TjjTfq5ptvls1mq/R4fn6+43atWrWqbS8AAAAA+BPdFQH41W233abDhw9LkqZPn64LLrggyC0KLm/ejz//+c/66aef1LBhQzVt2lQ//fSTW9uyj8flbMB56beQyzAMderUyXF/Xl6eY7D4rl27ql69epWe+/HHH+uDDz5Q+/btZRhGpcf379/vuN2yZUu32gwAAAAA/kLIhciTcO6ZcaoiRcK5AdvUq6++qhUrVkiSfve73+nuu+8O2LadyczMdPrY66+/rjfeeEMXX3yx/vnPfzpdrvyA6p7w9v2YNWuW2rVrp1atWumFF17QU089Ve1zCgsL9fnnn0tybzyuNm3aKDEx0XH/9u3bVVZWJklOq8DWrVsnSWrfvr3LdUtnBr4HAAAAgEAi5ELksVgCOhB7pMjMzNSTTz4pSWrWrJlmzZoV5BadcdFFFzl9LCkpSZJUu3Ztl8t5w5f3o3///h5vb/Xq1SosLFStWrV01VVXOV3OXslVcTyuXbt2OW63bVt5YoaioiJ99dVXkqTmzZtXue61a9dKkurUqVNlUDZ//nyNGTNGb7zxhi688EK99NJL+vrrr2UYhvr166cZM2aoSZMmrl8oAAAAADjBmFwAfFZUVKQxY8aoqKhIhmHo7bffVoMGDTxaR2FhobKzs6tcd15enr+aGhD+eD88Ze+q2L9/f8XHx1e5THFxsb7//ntJlcfj+vXXXx23q+qquHz5chUUFEg6EwpWtHXrVu3Zs0fSmUqy6OjKv6HYK70+//xzDRkyRHXq1NFdd92lNm3a6F//+pfuueeeal4lAAAAADhHyAXAZ+PHj9f27dslSY8++qhHlUhWq1UPP/ywEhMT1aJFC11yySXKzs6W1WrVgw8+qDp16qhevXrq16+fDh48WFMvwa98eT+8YbPZHIPOu+qq+P333ztmPaxYyVXeqVOnKt03Y8YMxzhcVT3+2muvOW7ffvvtVa7XHnJt3bpVmZmZ+uCDDxzVXC1bttTq1aurHNAeAAAAANxByAXAJ8uWLdPrr78u6UxwMmXKFI+e//bbb+vrr7/Wzp07tWfPHlksFl1//fWaMmWKFi5cqBUrVig3N1ddunTRuHHjauIl+JWv74c3Nm3apIMHD8owDJeDzpcfM6tiJVeLFi0ct7/88suzHluwYIHWrl2rtLQ0SWdmWSzvq6++csz+2KtXLw0YMMDl9j/44AO1adPGcX9sbKxat26t06dPO207AAAAAFSHMbkAeO3nn3/WnXfeKUmyWCz605/+pJ07d7p8Tu3atdW6dWvH/8+dO1d/+ctfHKHHkiVLdOGFF+rbb7/VmjVr1KdPH0nSSy+9pBYtWuj48eM699zADabvCX+8H96wV3F1795dTZs2dbqcfTyu+vXrV5r9sE+fPoqLi9Pp06e1bNkyTZ48WUOGDNG3336rJ598UrGxsZo6daouvfRSbdmyxfF4enq6JkyYoLKyMiUkJGj27NlVbvvQoUM6fPiwevbsWWUV2Z49e9S6desqZ20EAAAAAHcQcgHwWlZWlo4cOSLpTLfD0aNHV/ucPn36aM2aNY7/z83N1fnnn+/4/6SkJN16662aN2+eo3JIkqKjo9WyZUvl5OSEbMjlj/fDG/bxuFx1VZR+C7kqVnFJZ8bheuqppzRp0iTZbDZNmjRJkyZNcjz+5ptv6pJLLlH//v316aefVno8MTFRixcvVkpKSpXbtldxVVXldezYMR04cECjRo1y2X4AAAAAcIXuigCCqnXr1vruu+8c///rr79q4cKFOnr0qGbMmOG4v6CgQLt27dJ5550XhFaGrpycHG3evFmS+yGXs/G4Jk6cqLfeekvdunVTrVq1VLt2bfXt21effPKJ7rrrLknSe++9pxtuuEH16tVTXFyczj//fD388MPasWOHy1kd7SFXt27dKj1mb39VMzICAAAAgLuo5ALgtb59+/o8UPi9996rhx9+WFFRUUpISNDTTz+tiy66SM8++6zuvfdexcTE6LLLLtPUqVM1ePBgJSYm+qn1qlSN5Ct/vB+esndVbN68ebUh0dGjR6td37hx41yOfZaUlKT58+d71kj9FnJ179690mNbtmyR5HowfAAAAACoDiEXgKAaOnSocnNzdf/99+vo0aMaPHiwZs6cqXr16unkyZN65plndPr0aV177bX6+9//Huzmhhx7yHXttdcGuSWuZWRkqEGDBmrVqlWlx6jkAgAAAOAPhFwAgu6ee+7RPffcU+n+hx56SA899JBsNltEDUg+a9YspaenS5IyMzMd99nH7kpLS3NUW/Xu3Vvdu3fXiBEjgtJWdxQUFOjHH39Uv379qnx88+bNaty4sctB8wEAAACgOoRcAEJeJAVckpSenq633377rPvWrVundevWOf7fHnI98cQTAW2bN7Zt2yar1VrleFwnT57Url27NHDgwCC0DAAAAEA4IeQCgBAzZ84czZkzJ9jN8JtevXo5HassMTFRVqs1wC0CAAAAEI6YXREAAAAAAACmR8gFAAAAAAAA0yPkAgAAAAAAgOkRcgEAAAAAAMD0CLkAAAAAAABgeoRcAAAAAAAAMD1CLgAAAAAAAJgeIRcAAAAAAABMj5ALYSMqKkqSVFZWJpvNFuTWAKhpNptNZWVlkn7b/wEAAABELkIuhI3Y2FhJZy58CwoKgtwaADWtoKDAEWjb938AAAAAkYuQC2Gjbt26jtvHjx+nmgsIYzabTcePH3f8f/n9HwAAAEBkIuRC2KhTp44Mw5AknTp1StnZ2crPzyfsAsKIzWZTfn6+srOzderUKUmSYRiqU6dOkFsGAAAAINiig90AwF8sFouaN2+unJwc2Ww2nTp1SqdOnZJhGIzXA4SJimPuGYah5s2by2LhNxsAAAAg0hFyIawkJiaeFXRJZyo/SktLg9wyAP5mD7gSExOD3RQAAAAAIYCQC2EnMTFR7du316lTp5SXl6fi4mLHDGwAzC0qKkqxsbGqW7eu6tSpQwUXAAAAAAdCLoQli8WiunXrMhg1AAAAAAARgp/AAQAAAAAAYHqEXCbx888/a/ny5ZowYYIGDx6shg0byjAMGYah2267rUa2OW/ePA0cOFBNmjRRfHy8WrVqpZtvvlnffPNNjWwPAAAAAADAW3RXNInGjRsHbFuFhYUaNWqUVqxYcdb9+/fv19y5czVv3jxNmDBBEydODFibAAAAAAAAXKGSy4RatmypgQMH1tj677jjDkfAdeWVV+rf//631q9fr9mzZ6tt27ayWq2aNGmS3nzzzRprAwAAAAAAgCeo5DKJCRMmqEePHurRo4caN26sffv2qXXr1n7fzueff6758+dLkoYOHarFixcrKipKktSjRw8NGzZM3bt31/79+/Xkk0/q+uuvV/369f3eDgAAAAAAAE9QyWUSkydP1rXXXlvj3RZffvllSVJ0dLRef/11R8Bl17BhQ7344ouSpBMnTmjWrFk12h4AAAAAAAB3EHLB4eTJk/rss88kSf3791dycnKVy1133XWqW7euJGnx4sUBax8AAAAAAIAzhFxw2LBhg4qLiyVJffr0cbpcbGysevXq5XhOSUlJQNoHAAAAAADgDGNywSErK8txu0OHDi6X7dChgz755BOVlpbqxx9/VEpKitvbyc7Odvn4gQMHHLcPHjzo9noBAAAAAID7yl9zl5aWBrEl/kHIBYfy4ZOzrop2LVq0cNw+cOCARyFX+edWp2fPnm4vCwAAAAAAvHPkyBGdd955wW6GT+iuCIeTJ086btepU8flsrVr13bcPnXqVI21CQAAAAAA1LzDhw8Huwk+o5ILDkVFRY7bsbGxLpeNi4tz3C4sLPRoO+W7I1Zl7969uuKKKyRJX3/9tUeVXzCfgwcPOir21q9fr6ZNmwa5RahJ/L0jC3/vyMLfO7Lw944s/L0jC3/vyHLgwAFddtllkqoftsgMCLngEB8f77htH4DemdOnTztuJyQkeLSd6rpClteiRQuPloe5NW3alL93BOHvHVn4e0cW/t6Rhb93ZOHvHVn4e0eW8pmAWdFdEQ6JiYmO29V1QczPz3fcrq5rIwAAAAAAQE0j5IJD+YTekxkQ6U4IAAAAAACCjZALDuVnSNyxY4fLZe2PR0dHq127djXaLgAAAAAAgOoQcsGhR48ejgHn165d63S54uJiffvtt47nxMTEBKR9AAAAAAAAzhBywSExMVFXXXWVJGn16tVOuywuWrRIeXl5kqQRI0YErH0AAAAAAADOEHJFkDlz5sgwDBmGoUmTJlW5zGOPPSZJKi0t1f3336+ysrKzHj969KiefPJJSVK9evU0bty4Gm0zAAAAAACAO6KD3QC4Jz09Xbt27XL8/9GjRx23d+3apTlz5py1/G233ebVdvr166cbb7xR8+fP19KlSzVgwAA98sgjatasmTIzM/X8889r//79kqQXX3xR9evX92o7AAAAAAAA/mTYbDZbsBuB6t122216++233V6+qj/rnDlzdPvtt0uSJk6c6LSaq7CwUKNGjdKKFSuqfNxiseiZZ55x+nwAAAAAAIBAo7siKklISNBHH32kuXPnasCAAUpKSlJsbKxatGihm266Senp6QRcAAAAAAAgpFDJBQAAAAAAANOjkgsAAAAAAACmR8gFAAAAAAAA0yPkAgAAAAAAgOkRcgEAAAAAAMD0CLkAAAAAAABgeoRcAAAAAAAAMD1CLgAAAAAAAJgeIRcAAAAAAABMj5ALAAAAAAAApkfIBQAAAAAAANMj5AIAAAAAAIDpEXIBAABAklRYWKjo6GgZhqEpU6YEuzkAAAAeIeQCAACAJGnLli0qKyuTJF188cVBbg0AAIBnCLkAAAAgSdq4caPjNiEXAAAwG8Nms9mC3QgAAAAAAADAF1RyAQAAAAAAwPQIuQAAAAAAAGB6hFwAAADQjh07ZBiGDMPQ/PnzKz1eUlKijh07yjAMderUyTFAfXVSU1NlGIZatWrl7yYDAACchZALAAAA2rp1q+N2ly5dKj3+8ssvKysrS5L06quvKioqyq31du/eXZK0f/9+ZWdn+95QAAAAJwi5AAAA4Ai5EhIS1K5du7MeO3TokJ577jlJ0jXXXKO+ffu6vd6OHTs6bm/atMn3hgIAADhByAUAAABHyJWamlqpSmvKlCkqKCiQJE2YMMGj9SYnJztu796928dWAgAAOEfIBQAAAGVkZEiSOnfufNb9x44d05tvvilJSktLU8+ePT1ab8OGDR23Dx486FsjAQAAXCDkAgAAiHBHjx5Vbm6upMrjcc2ZM0enT5+WJN1xxx0+baekpMSn5wMAALhCyAUAABDhyg86X7GSa8GCBZKkqKgojRgxotJzN2zYoBtvvFE333yzbDZbpcfz8/Mdt2vVquWvJgMAAFQSHewGAAAAILjsIZdhGOrUqZPj/ry8PMdg8V27dlW9evUqPffjjz/WBx98oPbt28swjEqP79+/33G7ZcuWfm45AADAb6jkAgAAiHD28bjatGmjxMREx/3bt29XWVmZpDMhV1XWrVsnSWrfvr3LdUtS9+7d/dBaAACAqhFyAQAARDh7JVfF8bh27drluN22bdtKzysqKtJXX30lSWrevHmV6167dq0kqU6dOpWCsvnz58swDP3973/X2rVrdc0116h+/fo699xzNWrUKB06dMjr1wQAACIPIRcAAEAEKy4u1vfffy+p8nhcv/76q+N2VV0Vly9froKCAklS7dq1Kz2+detW7dmzR5I0dOhQRUefPVKGvcrr888/15AhQ1SnTh3dddddatOmjf71r3/pnnvu8fp1AQCAyMOYXAAAABHs+++/d8x6WLGSq7xTp05Vum/GjBkyDEM2m63Kx1977TXH7dtvv73S4/aQa+vWrcrMzFSbNm0knQne2rVrp9WrV8tms1U51hcAAEBFVHIBAABEsPJjZlWs5GrRooXj9pdffnnWYwsWLNDatWuVlpYm6cwsi+V99dVXmjNnjiSpV69eGjBggNNtf/DBB46AS5JiY2PVunVrnT592uPXAwAAIhchFwAAQASzj8dVv379SrMf9unTR3FxcZKkZcuWafLkydqwYYNee+013XbbbYqNjdXUqVNlGIa2bNnieHz69OkaMmSIysrKlJCQoNmzZ1fa7qFDh3T48GH17NmzygqyPXv2qHXr1lRxAQAAtxFyAQAARDB7yFWxiks6Mw7XU089JUmy2WyaNGmSevbsqYceekiFhYWaMWOGLrnkEvXv31+SHI8/+uijOnXqlBITE7Vs2TKlpKRUWre9iquqCq9jx47pwIEDLrtPAgAAVETIBQAAEMGczaxoN3HiRL311lvq1q2batWqpdq1a6tv37765JNPdNddd0mS3nvvPd1www2qV6+e4uLidP755+vhhx/Wjh07dNVVV1W5XnvI1a1bt0qPbd68WZIqzcYIAADgCgPPAwAARLCjR49Wu8y4ceM0btw4p48nJSVp/vz5Hm3XHnJ179690mNbtmyR5HogfAAAgIqo5AIAAEDAZWRkqEGDBmrVqlWlx6jkAgAA3iDkAgAAQEAVFBToxx9/dBpibd68WY0bN1bTpk0D3DIAAGBmhFwAAAAIqG3btslqtVY5HtfJkye1a9cuuioCAACPGTabzRbsRgAAAAAAAAC+oJILAAAAAAAApkfIBQAAAAAAANMj5AIAAAAAAIDpEXIBAAAAAADA9Ai5AAAAAAAAYHqEXAAAAAAAADA9Qi4AAAAAAACYHiGXSRiG4da/vn37+mV78+bN08CBA9WkSRPFx8erVatWuvnmm/XNN9/4Zf0AAAAAAAD+ZNhsNluwG4HqGYbh1nJ9+vTRmjVrvN5OYWGhRo0apRUrVlT5uMVi0YQJEzRx4kSvtwEAAAAAAOBvhFwmYQ+5/vCHP+i+++5zulzt2rXVunVrr7czZswYzZ8/X5J05ZVX6uGHH1azZs2UmZmpKVOmaPfu3ZKkmTNn6u677/Z6OwAAAAAAAP5EyGUS9pBr4sSJmjRpUo1s4/PPP9dVV10lSRo6dKgWL16sqKgox+NHjx5V9+7dtX//ftWrV0979uxR/fr1a6QtAAAAAAAAnmBMLji8/PLLkqTo6Gi9/vrrZwVcktSwYUO9+OKLkqQTJ05o1qxZAW8jAAAAAABAVQi5IEk6efKkPvvsM0lS//79lZycXOVy1113nerWrStJWrx4ccDaBwAAAAAA4AohFyRJGzZsUHFxsaQzg9c7Exsbq169ejmeU1JSEpD2AQAAAAAAuELIZTILFy5USkqKatWqpcTERLVr10633nqrvvjiC5/Wm5WV5bjdoUMHl8vaHy8tLdWPP/7o03YBAAAAAAD8ITrYDYBnyodRkrRr1y7t2rVL77zzjn73u99pzpw5Ouecczxeb3Z2tuO2s66Kdi1atHDcPnDggFJSUrzeVlWKioq0Y8cONW7cWI0aNVJ0NB9TAAAAAAD8rbS0VEeOHJEkpaamKj4+Psgt8g3pgUnUqlVLw4YN01VXXaUOHTqoTp06OnLkiNauXau///3vOnbsmP79739r+PDh+vTTTxUTE+PR+k+ePOm4XadOHZfL1q5d23H71KlTnr0QnR2SAQAAAACA4Fu/fr169OgR7Gb4hJDLJHJyclSvXr1K9w8YMEAPPvigBg8erC1btmjt2rV644039NBDD3m0/qKiIsft2NhYl8vGxcU5bhcWFnq0HQAAAAAAgJpAyGUSVQVcdo0bN9aHH36oDh06qKSkRK+99prHIVf5kkT7APTOnD592nE7ISHBo+1IZ7o4Vvf4ZZddJkla/9ajanr1Hz3eBgAAAAAAcO3gwYPq2bOnJKlRo0ZBbo3vCLnCRJs2bTRgwACtWLFCu3btUm5urpo1a+b28xMTEx23q+uCmJ+f77hdXdfGqlQ35ld5TRue49HyAAAAAADAc+EwHjazK4aR8gPA5+TkePTc8kFSdQPDl6/EqvHxtWw1u3oAAAAAABAeCLnCiGEYXj+3fEC2Y8cOl8vaH4+Ojla7du283iYAAAAAAIC/EHKFkaysLMdtT7oqSlKPHj0cA86vXbvW6XLFxcX69ttvHc/xdBZHz1HKBQAAAAAAqkfIFSb27t2rTz/9VJLUtm1bNW/e3KPnJyYm6qqrrpIkrV692mmXxUWLFikvL0+SNGLECB9aDAAAAAAA4D+EXCawbNkylZaWOn388OHDGjlypGNWxPvuu6/SMnPmzJFhGDIMQ5MmTapyPY899pgkqbS0VPfff7/KysrOevzo0aN68sknJZ2Z7XHcuHHevBzP2KjkAgAAAAAA1TP/0PkR4MEHH1RJSYlGjhypSy+9VOedd54SEhJ09OhRrVmzRjNnztTRo0clSWlpabr//vu92k6/fv104403av78+Vq6dKkGDBigRx55RM2aNVNmZqaef/557d+/X5L04osvqn79+n57jc4RcgEAAAAAgOoRcplEbm6uXnvtNb322mtOlxk5cqRmzZqluLg4r7fzj3/8Q3l5eVqxYoW++OILffHFF2c9brFY9Mwzz+juu+/2ehsAAAAAAAD+RshlAm+//bbWrl2rb775Rnv27NHRo0eVl5enOnXqqEWLFrrssst066236tJLL/V5WwkJCfroo4/0/vvva86cOdq6datOnDihxo0bq3fv3nrggQf8sh230V0RAAAAAAC4wbDZSBEQWrKzs9WiRQtJ0oEPn1HyyGeD3CIAAAAAAMLPWdffBw4oOTk5yC3yDQPPI8SRwQIAAAAAgOrRXRFhqbi4WKdOnVJ+fr6Ki4tltVqD3SQAfmCxWBQbG6vatWurTp06io2NDXaTAAAAAIQIQi6ENg9709psNh09etQx2ySA8GMPsQ8fPqxGjRqpQYMGMgwj2M0CAAAAEGSEXAgrBw8e1K+//nrWfYZhKCoqKkgtAuBPZWVlKj+U5JEjR1RcXKxmzZoFsVUAAAAAQgEhF0Kc+5VcRUVFZwVcDRo0UN26dRUXF0eVBxAmbDabTp8+rby8PB07dkyS9Ouvv6pBgwaKi4sLcusAAAAABBMDzyO0edBb8cSJE47bSUlJSkpKUnx8PAEXEEYMw1B8fLxjH7f75ZdfgtgqAAAAAKGAkAtho6CgwHG7Xr16wWsIgIAov5+X3/8BAAAARCZCLoQ490u5ysrKJEnR0dGMwQVEgKioKMe+bt//AQAAAEQuQi4AgGnRHRkAAACAHSEXAAAAAAAATI+QCwAAAAAAAKZHyIUQ58H0igAAAAAAIGIRcgEAAAAAAMD0CLkQ2mxUcgEAAAAAgOoRciHEEXKZSXFxsebNm6exY8eqQ4cOatCggWJiYtSwYUN1795df/jDH7R69WpZrVZJUt++fWUYhk//5syZ49i+/b6+ffu61d45c+ZUuZ7777/fcf/cuXPdWteuXbtUq1YtGYahrl27qrS0VJK0b98+l+2vXbu2Wrdureuvv17//ve/3doWAAAAAKAyQi4AfrFo0SJdcMEFuummm/Tuu+/qhx9+0PHjx1VaWqpjx45p8+bN+vvf/64BAwbowgsv1EcffRTsJjv1wgsvqGXLlpKkRx55REeOHHG5vM1m01133aXCwkJFR0frH//4h6Kjo93aVkFBgfbt26cPP/xQI0aM0KBBg5Sfn+/zawAAAACASOPeVRgQLHRXNIX//d//1YQJExz/P2DAAA0bNkwpKSmqV6+ejh8/rh9++EHLli3Tp59+qp07d+rpp5/W4sWLnQY6f/7zn7VkyRJJ0qpVq9SsWbMql0tOTvb760lMTNTMmTM1ePBgHT16VA899JDmzZvndPlZs2ZpzZo1kqTHHntMXbt2rXK54cOH67nnnjvrvl9//VWbNm3SK6+8or179+qTTz7RPffco/fee89vrwcAAAAAIgEhFwCf/POf/3QEXElJSVqwYIH69OlTabn+/fvr/vvv13fffac//vGPOnLkiFq3bu10vfXq1XPcbt++vc477zx/N92lq6++WmPHjtU777yj+fPn6/e//72uvfbaSsvl5ubq8ccflyRdcMEFmjhxotN11qtXTxdddFGl+y+//HLdcMMNSk1N1ZEjR/T+++/rpZdechrsAQAAAAAqo7siQhyVXKEsJydHDzzwgCSpdu3aWrt2bZUBV3kXXXSRVq1apcceeywQTfTJ9OnT1bhxY0nSvffeq7y8vErL3Hffffr1119lGIZmzZql+Ph4r7bVuHFjjR07VtKZ7o8bN270vuEAAAAAEIEIuYAAsFptKiguldUaXqHd9OnTVVBQIEl69tln1aFDB7eeZ7FYdPPNN9dk0/zi3HPP1YwZMySdCfTsFVt2CxcudHSpvP/++5WWlubT9spXtp0+fdqndQEAAABApKG7IkKbycfkysrN06z0PVqZeUiFJWVKiInS4NQmGpfWRinN6ga7eT6x2Wx6++23JZ2p4rrrrruC3KKaMWrUKI0YMUKLFy/WW2+9pTFjxqhv37765Zdf9OCDD0qSWrVqpb/85S8+b+unn35y3LYPfA8AAAAAcA+VXAhx5g25lmTkaNiMdC3anKPCkjJJUmFJmRZtPnP/koycILfQN9u3b9fRo0clSb1791ZiYmKQW1RzXn/9ddWvX182m+3/2bvv+Kaq/3/gr3uTpru0FEopZUOZBSpDlI1QBAUEkSFbUYZbvh/5oB+WP0UEcSCiCAiCCIKibEWmFNGyKZYhlNEFlNndZtzfH21C0iZt9mhfz8ej9Obem3NOk6ZtXrzPuZgwYQLy8vLwxhtv4MaNGwCApUuXIiAgwKY+bt68idWrVwMAateujYceesjmcRMREREREVUmrOQicoDEtExM3XAKKhPTE1UaCVM3nELjsECPreg6deqUbrtt27YuHElpOTk5OHPmTLnnpaaaFzSGh4dj4cKFeO6553Dp0iX0798fe/bsAQCMHTsWffr0Maude/fulRpXZmYmjh8/js8++ww3btyAQqHA559/Di8vL7PaJCIiIiIioiIMuci9eeh0xeVxSSYDLi2VRsKKuMtYOLS1k0ZlX7dv39Zth4WFuXAkpR09ehTR0dF2bXP8+PFYt24dfv/9d13AFR4ejk8++cTsNjZv3qxbw8uYYcOG4T//+Y/bhYZERERERESegNMViexMo5GwM+G6WefuSEj32MXos7KydNv+/v4uHInzLFu2zGBa4uLFixESEmK39rdu3YrFixfj3r17dmuTiIiIiIiosmDIRWRn+Sq1bg2u8uQp1chXmXeuu9FfgysnJ8eFIymtW7dukCSp3I+VK1da1G7dunUNqqyefvppi+4/duzYUmPIz8/H+fPnMXfuXAiCgFWrVqFTp0669b6IiIiIiIjIPAy5iOzMRy6Dr5fMrHN9vWTwkZt3rrsJDQ3VbTOQsZ63tzeioqIwffp0bNq0CQCQmJiIqVOnunhkREREREREnoUhF7k3D1yTSxQF9I0ON+vcftE1IYqCg0fkGK1bP1hL7Pjx4y4cScURGxure1w3bNjgdhVyRERERERE7owhF5EDTOjcAPJywiu5KOD5zvWdNCL7a9GiBapVqwYAOHjwIDIzM108ooqhadOmAAClUolz5865eDRERERERESegyEXuTnPq+QCgOYRQVg4tLXJoEsuClg4tDWaRwQ5eWT2IwgCxo4dC6BoTa7ly5e7eEQVg0qlMrpNREREREREZWPIRe7NA6crag1sUwtbXu6Mpx+K1K3R5eslw9MPRWLLy50xsE0tF4/Qdm+88Qb8/PwAADNnzjS78kij0WDt2rWOHJpHkiQJx44d092uXbu2C0dDRERERETkWeSuHgBRRaat6FowpBXyVWr4yGUeuwaXMbVq1cLixYvx3HPPIScnB926dcOGDRvQrVs3k/dJTEzE66+/jps3b2LkyJFOHK37W7JkCa5cuQIAaNOmDSIiIlw7ICIiIiIiIg/CkIvcnOdWcukTRQF+ior5chs/fjxSUlIwc+ZM3Lx5E927d0dsbCwGDhyIZs2aITg4GHfu3MGFCxewfft2/Prrr1Cr1QYL11cW9+7dw5kzZwz2FRYW4sqVK/jxxx+xbt06AIAoivjwww9dMUQiIiIiIiKPVTHfdRORU82YMQMtWrTA1KlTceXKFezatQu7du0yeX6LFi0wf/58J47QPWzevBmbN28u85yAgAB8+eWXiI2NddKoiIiIiIiIKgauyeUhjh49infffRexsbGIjIyEt7c3AgICEBUVhfHjxyMuLs4u/cyePRuCIJj1sX//frv0WSYPXpOrshk8eDDOnz+PtWvXYtSoUWjSpAlCQkIgl8tRtWpVPPTQQ5gyZQr27t2LhIQEhjjFvLy8UK1aNXTp0gXvvvsuzp8/j1GjRrl6WERERERERB6HlVweoGvXrjh48GCp/YWFhfj333/x77//YtWqVRgzZgyWLVsGhULhglESAQqFAs8++yyeffZZm9tatWoVVq1aZdF9JAtD0XHjxmHcuHEW3cfScLdevXoWj4uIiIiIiIgsx5DLA6SlpQEAIiIi8Mwzz6BLly6oU6cO1Go1Dh8+jIULFyI1NRWrV6+GUqnE999/b5d+ExISyjxev359u/RDRERERERERGQrhlweoGnTppg7dy6efvppyGQyg2MdO3bE6NGj0alTJ1y4cAHr1q3DpEmT0LVrV5v7bdmypc1tEBERERERERE5A9fk8gDbtm3D0KFDSwVcWtWqVcPChQt1t3/88UdnDc3xOM2LiIiIiIiIiMzAkKuC6NGjh2770qVLLhyJvTHkIiIiIiIiIqLyMeSqIAoKCnTbpiq+iIiIiIiIiIgqKoZcFcSBAwd0282aNbNLm7GxsQgLC4NCoUBYWBi6d++OefPm4e7du3Zp3zys5CIiIiIiIiKi8nHh+QpAo9Fg3rx5uttDhw61S7u///67bjsjIwMHDhzAgQMH8OGHH2LVqlUYOHCgVe2mpKSUeTw9Pd2qdomIiIiIiIio8mLIVQF88skniI+PBwAMHjwYbdu2tam96OhoPPXUU+jQoQMiIiKgVCpx/vx5rF27Frt27cK9e/fw9NNPY+vWrejbt6/F7deuXdum8RERERERERERlSRIEi9f58kOHDiAXr16QaVSISwsDAkJCQgLC7O6vXv37iE4ONjk8aVLl2LSpEkAgIiICFy6dAk+Pj4W9SEIgtnnJq+ejMjRS8w6999//4VKpYJcLkfjxo0tGhMReSa+7omIiIiIrJeSkqIrRElOTkZkZKSLR2QbVnJ5sH/++QeDBg2CSqWCj48PNm7caFPABaDMgAsAJk6ciCNHjmDFihVIS0vDTz/9hJEjR1rUR3JycpnH09PT0aFDh6IbzGCJiIiIiIiIyAwMuTzU5cuXERsbi7t370Imk2H9+vXo2rWrU/qeOHEiVqxYAaCokszSkMuyZJghFxERERERERGVj1dX9EBpaWno1asX0tLSIAgCvvnmG6sXgbdG8+bNddupqalO65eIiIiIiIiIyBSGXB7m1q1b6N27N5KSkgAAn3/+OcaMGePUMViyppbNOF2RiIiIiIiIiMzAkMuD3L9/H3369EFiYiIAYN68eXjppZecPg5t/0DR4vNERERERERERK7GkMtD5Obm4oknnsDx48cBAO+88w6mTZvmkrEsXbpUt92tWzcH98ZKLiIiIiIiIiIqH0MuD1BYWIhBgwbh0KFDAIDXXnsN7733nsXtrFq1CoIgQBAEzJ49u9TxhIQEXLx4scw2vv76ayxfvhwAEB4ejkGDBlk8DiIiIiIiIiIie+PVFT3AiBEjsGvXLgBAz5498fzzz+PMmTMmz1coFIiKirK4n2PHjmHChAno0aMH+vbti+joaISGhkKlUuHcuXNYu3atbhwymQxff/01/P39rfuizMU1uYiIiIiIiIjIDAy5PMCmTZt023v37kWrVq3KPL9u3bq4cuWKVX2p1Wrs3r0bu3fvNnlOaGgoVqxYgf79+1vVBxERERERERGRvTHkIp1+/fphxYoVOHz4ME6cOIEbN27g9u3bkCQJVatWRevWrfH4449j3LhxCAoKctKoWMlFREREREREROVjyOUBJDtN2Rs3bhzGjRtn8nhYWBiee+45PPfcc3bpzy44XZGIiIiIiIiIzMCF54mIiIiIiIiIyOMx5CI3x0ouIkfIy8uDXC6HIAiYO3euq4dDRERERERkM4ZcRGS1/fv3QxAEiz5mz57t6mE7jDs9Hqmpqbo+Tpw4Uer4iRMnoFarAQDt2rVzyBiIiIiIiIiciSEXuTcWchFZZdu2bQCAyMhIxMTElDp+9OhR3TZDLiIiIiIiqgi48DwRWa1ly5aoU6cOrl27BkEQsGLFCrRv377M+4SFhTlpdM7nTo/H1q1bAQBPPvmk0eOvvvoqXn31VYf0TURERERE5AoMucjNsZTLnVWrVg3btm1Dp06dkJWVhbfffht///036tSp4+qhmWX27NmYM2cOunXrhv3799vcnrs8Hrm5udizZw8AoH///k7tm4iIiIiIyFU4XZGIbBIdHY1169ZBFEVcv34dTz75JLKyslw9LJdxh8dj9+7dyM/Ph5+fH3r27OnUvomIiIiIiFyFIRe5N4mVXJ7giSeewMcffwwASEhIwLBhw3SLmldGrn48tOtx9erVCz4+PqWOnzt3Trco/fr160sdVyqVaNGiBQRBQKtWrcwee3R0NARBQN26dW37AoiIiIiIiKzAkIvcHEMuT/Haa69h8uTJAICdO3dW+vWeXPV4SJKkC7lMTVU8deqUbrtNmzaljn/00UdITEwEACxatAgymcysvtu2bQsAuHbtGlJSUiwZNhERERERkc0YchGR3SxatAixsbEAgCVLluDTTz917YBczBWPx7Fjx5Ceng5BEPDEE08YPUcbcvn6+qJx48YGx65fv4733nsPQFFFWvfu3c3uu0WLFgbjICIiIiIiciaXLTyfkZGBpKQkXL9+HTk5OfDy8kJwcDDq1KmDRo0amV05QBWcIwq5NBog744DGnZTvlUB0Tl5tlwux4YNG/Doo48iMTERU6dORcOGDSvt4ueueDy0V1Vs27YtatasafQcbcgVHR1d6mft3LlzkZubCwCYOXOmRX1HRkbqti9dumTRfYmIiIiIiGzltJArJycHmzdvxs6dO3HgwAGkpqaaPNfb2xsxMTGIjY3FoEGD0KpVK2cNkyqDvDvAgoauHoXz/OcS4F/Nad1VqVIF27ZtQ5s2bZCZmYkRI0bg5MmTaNSokdPG4E6c/XiUN1URAE6ePAkAaN26tcH+27dv4+uvvwYAdO7cGR06dLCo72rVHnyfpaenW3RfIiIiIiIiWzk85Dpx4gQ+//xzbNy4UVcdIJWzmHh+fj4OHz6Mv/76C++++y5atGiBl156CaNHj4afn5+jh0xuhWtyeaILFy7orijYtGnTSr8QubMej9TUVBw/fhyA6ZDr1q1bSEtLA1B6Pa5Vq1ahoKAAAPDcc8/ZNBalUmnT/YmIiIiIiCzlsJDrxIkTmDFjBnbu3AngQbAVHh6ODh06oG3btggLC0PVqlUREhKCvLw83LlzB3fv3sWFCxdw5MgRnD59GkqlEmfOnMGUKVMwY8YMvPXWW3jllVfg7e3tqKETkQ0yMjIwduxYSJIEPz8/rF27Fl5eXi4bjyAI5Z5z4MCBMs9buXIlxo0bZ1X/1jwe3333HQ4ePIhjx44hISEBhYWFZo1BW8UVGRmJmJgYo+foLzpfspJrw4YNAACZTIZBgwaVuu+RI0ewcOFCyOVyrFmzptRjlpOTo9vmf0gQEREREZGzOSTkGj9+PNasWQONRgMAeOihhzBy5Eg8/fTTqFOnjtntFBYW4o8//sDatWvx888/49atW5g2bRqWLFmC1atXo3Pnzo4YPrmTcqr+yP2MGzcON27cAAB88sknaNKkiYtH5FrWPB7/+9//cPXqVVSrVg01a9bE1atXzepLux6XqQXngQchlyAIBlPBMzMzdYvFx8TEIDg4uNR9f/31V/zwww+IiooyGgpeu3ZNt23Jz3oiIiIiIiJ7cEjI9e2330KhUGDs2LGYOnUqoqKirGpHoVCgV69e6NWrF7766its3LgRc+fOxblz57B3716GXJWCA0Iu36pF61RVFr5VndbVokWLsGPHDgDAU089hRdffNFpfZuSkJBg8tiSJUvw5Zdfol27dli5cqXJ8/QXVLeEtY/H8uXL0bhxY9StWxfz5s3D9OnTy71PXl4e9u7dC8C89bgaNGiAwMBA3f5//vkHarUaAExWgR06dAgATP5M17YNFC18T0RERERE5EwOCbmmTJmCadOmoXbt2nZr09vbG6NGjcLIkSOxceNG3ZsxIouJolMXYq8sEhISMG3aNABAREQEli9f7uIRFWnZsqXJY2FhYQAAf3//Ms+zhi2PR69evSzub/fu3cjLy4Ofnx8ee+wxk+dpK7lKrsd18eJF3XbDhqUvzJCfn4+DBw8CAGrVqmW07QMHDgAAAgICjAZl69evx4gRI/Dll1+iWbNmmD9/Pv78808IgoCePXti8eLFCA8PL/sLJSIiIiIiMkF0RKOLFy+2a8ClTxAEDB06FCNGjHBI++RuOF3RE+Tn52PEiBHIz8+HIAj49ttvERoaalEbeXl5SElJMdp2ZmamvYbqFPZ4PCylnarYq1cv+Pj4GD2nsLAQZ8+eBVB6Pa779+/rto1NVdy2bZvu4iH+/v6ljp86dQpJSUkAiirJ5PLS/4eirfTau3cv+vXrh4CAALzwwgto0KABfvrpJ0ycOLGcr5KIiIiIiMg0h4RcRFS5TJ06Ff/88w8A4M0337SoEkmj0eC1115DYGAgateujYcffhgpKSnQaDR45ZVXEBAQgODgYPTs2RPp6emO+hLsypbHwxqSJOkWnS9rquLZs2d1Vz0sWcmlLzs7u9S+xYsX69bhMnb8888/122PHz/eaLvakOvUqVNISEjADz/8oKvmqlOnDnbv3l3u1XeJiIiIiIhMcVjI9cYbbxisz0JkHb7hdXdbt27FkiVLABQFJ3PnzrXo/t9++y3+/PNPXLhwAUlJSRBFEc888wzmzp2LjRs3YseOHUhLS0ObNm0wYcIER3wJdmXr42GNY8eOIT09HYIglLnovP7P5JKVXPrVt3/88YfBsQ0bNuDAgQO6dRCPHDlicPzgwYNYtWoVAKBjx47o3bt3mf3/8MMPaNCggW6/QqFA/fr1UVBQYHLsRERERERE5XHImlwA8Nlnn2HRokVo0aIFxowZg5EjR6JmzZqO6o6IXODmzZt4/vnnAQCiKOLtt9/GhQsXyryPv78/6tevr7u9du1afPDBB7rQY/PmzWjWrBn++usv7N+/H926dQMAzJ8/H7Vr18adO3dQtarzFtO3hD0eD2toq7jatm1b5s9Z7XpcISEhpa5+2K1bN3h7e6OgoABbt27FnDlz0K9fP/z111+YNm0aFAoFFixYgEceeQQnTpzQHY+Li8PMmTOhVqvh6+uLFStWGO37+vXruHHjBjp06GC0iiwpKQn169c3etVGIiIiIiIiczgs5AKKptD8888/mDZtGqZPn46ePXti7NixGDRoEHx9fR3ZNVUULORya4mJicjIyABQNO1w6NCh5d6nW7du2L9/v+52WloaGjVqpLsdFhaGsWPHYt26dQZXUJXL5ahTpw5SU1PdNuSyx+NhDe16XGVNVQQehFwlq7iAonW4pk+fjtmzZ0OSJMyePRuzZ8/WHf/666/x8MMPo1evXvj9999LHQ8MDMTPP/+M5s2bG+1bW8VlrMrr9u3bSE5OxpAhQ8ocPxERERERUVkcNl3xt99+w6hRo+Dn5wdJkqBWq7F7926MHj0a4eHheO6557Bv3z5HdU9EHqJ+/fo4c+aM7vb9+/exceNG3Lp1C4sXL9btz83NxcWLF1GvXj0XjNJ9paam4vjx4wDMD7lMrcc1a9YsLFu2DA899BD8/Pzg7++P7t27Y9euXXjhhRcAAN999x2GDRuG4OBgeHt7o1GjRnjttddw7ty5Mq/qqA25HnrooVLHtOM3dkVGIiIiIiIiczmskqt3797o3bs3cnNzsWnTJqxZswZ79uyBRqNBVlYWvv32W3z77beIjIzEqFGjMHr0aDRt2tRRwyGPxVIud9a9e3ebFwqfNGkSXnvtNchkMvj6+uKdd95By5Yt8e6772LSpEnw8vLCo48+igULFqBv374IDAy00+hRqhrJVvZ4PCylnapYq1atckOiW7duldvehAkTylz7LCwsDOvXr7dskHgQcrVt27bUsRMnTgAoezF8IiIiIvJcGo2EfJUaPnIZRJHLU5DjOHS6IgD4+flh1KhRGDVqFNLT07F27Vp89913OH36NAAgOTkZ8+bNw7x589C2bVuMHTsWw4cPR2hoqKOHRp6AV1qr8Pr374+0tDS89NJLuHXrFvr27YulS5ciODgYWVlZmDFjBgoKCvDkk0/iq6++cvVw3Y425HryySddPJKynTx5EqGhoahbt26pY6zkIiIiIqqYEtMysTwuCTsTriNPqYavlwx9o8MxoXMDNI8IcvXwqAISJBddrz0hIQGrV6/GunXrkJaWVjSY4gWHvby88Pjjj2PMmDHo378/vLy8XDFEcpGUlBTdld6SvxqGyInmVY38+++/UKlUkMvlaNy4sSOHSE4mSVKlWpB8+fLliIuLA1D0s/L48ePo1KmTbu2yzp0766qt5s+fj9zcXAwaNMjoWlvuIDc3F4GBgejZsyd+//33UsejoqKQmZmJ69evW9w2X/dERERE7mnzyVRM3XAKKk3pyEEuClg4tDUGtqnlgpGRPoP338nJiIyMdPGIbOPwSi5ToqOjsWDBAsyfPx979uzBmjVrsGnTJuTk5KCwsBBbt27F1q1bERISgmHDhmH06NHo2LGjq4ZLLsNKLkKlCrgAIC4uDt9++63BvkOHDuHQoUO629qQ66233nLq2Kxx+vRpaDQao+txZWVl4eLFi4iNjXXByIiIiIjIERLTMk0GXACg0kiYuuEUGocFsqKL7MpllVzG6K/ftXfvXqjVat0xURShUqlcODpyFsNKrqGInPiDWfdjRQdR5cPXPREREZH7eXPDSWw6nlrueU8/FImFQ91zNkJlUdEquRx2dUVraNfv+u2333Dy5Em0aNFCV8HhRlkcOROfdyIiIiIiIo+h0UjYmWDeMhQ7EtKhMVHtRWQNl01XNEapVGLr1q347rvvsGPHDiiVSlcPiYiIiIiIiIjMlK9SI0+pLv9EAHlKNfJVavgp3CqaIA/mFt9Jhw4dwpo1a7Bx40bcu3cPwIPKrcDAQAwZMgRjx4514QjJdZjqExEREREReQofuQy+XjKzgi5fLxl85DInjIoqC5eFXP/++y/WrFmDtWvX4sqVKwAeBFsymQy9evXCmDFjMGjQIPj4+LhqmEREREREREQOpdFIyFep4SOXQRQ9+6JLoiigb3S4WWty9Yuu6fFfL7kXp4Zct27dwvr167FmzRocPXoUgOFaW9HR0RgzZgxGjhyJ8PBwZw6N3BXX5CIiIiIiogrAWJCVmJaJ5XFJ2JlwHXlKNXy9ZOgbHY4JnRvY/aqDzgzSJnRugC0n00xeXREA5KKA5zvXd+g4qPJxeMhVUFCAzZs347vvvsNvv/2mu0KiNtyqUaMGnn32WYwZMwatW/OqClQSQy4iIiIiIvIcJcMkU0FWo+oB+Pj3CwZBUJ5SjU3HU7HlZBoWDm2N/q0ibA6mnBmkaTWPCMLCoa0xdcMpo0GXXBSwcGhrh/Wvzx2r5NxxTBWFw0Ku/fv347vvvsNPP/2EzMxMAA+CLR8fHwwYMABjxoxBnz59IJNxDq4lrl69ikWLFmH79u1ITk6Gt7c3GjZsiKFDh+Kll16Cn5+fXfrZuXMnvv76axw5cgQZGRmoXr062rdvjxdffBF9+/a1Sx9EREREREQVgbEwqWWtIBy/dg9qI0FWWVQaCa+vP4m3fjyNApXG6mBq88nUUkGTI4I0Lf3wZmCbWmgcFogVcZexIyFd95j0i66J5zvXd3jA5YpwzxVjYmBmSJAkx8wHE0URgiDogi1BENCpUyeMGTMGQ4cORVCQa76pPN3WrVsxatQoXXBYUlRUFLZv345GjRpZ3YdGo8GLL76IFStWmDxnwoQJWLp0KURRtLofU1JSUlC7dm0AQPIXgxA5ZZNZ9/v333+hUqkgl8vRuHFju4+LiNzPhQsXoFar+bonIiKqpFz5Bl+/762n00xWLdmTtgJqYJta5Y7p3PUsDFgcV+aYBAAKuWhTkAaUH97Y63kytx1j4Z5WeY+ho9h7TPYKzAzefycnIzIy0uz7uiOHhlwA0LBhQ4wePRqjR49G/fqcb2uLEydOoFOnTsjLy0NAQACmT5+OHj16IC8vD+vXr8eyZcsAFAVdR48eRWBgoFX9TJ8+HfPmzQMAxMTE4K233kLDhg1x6dIlzJ8/HydOnNCdN3fuXPt8cXqsDbmSkpJQUFAAoOgxYIUgUcWmVqtx4cIFAIC3tzcaNGjg4hERERF5Pk+pCnFllU7Jvr3lIgpVGqcttCIXBWx5ubPB12ns8ageqMC1O3lWtW9J4OKIQMncKZ/Gnu/EtMxywz1jj6Ej2XtM9nzMGXKZadKkSRgzZgweffRRRzRfKXXt2hUHDx6EXC7HH3/8gUceecTg+IIFC/DWW28BAGbNmoXZs2db3MeFCxfQokULqFQqtGvXDn/88Qd8fX11x3Nzc9GtWzccPXoUcrkcZ8+etalqzBhrQ67r16/j7t27AICwsDCEhobadVxE5F5u376NmzdvAgBCQkJ4wRIiIiIbuOPULlOcXaXjioqt8jz9UCQWDi1a07qsx8Na5gYu9g5vzJ3yqd92yef7zQ0nzbqyo/5jqM8RQa+tY9Jn78e8ooVc9p9rVuyrr75iwGVH8fHxOHjwIADg+eefLxVwAcDUqVPRrFkzAMBnn30GpVJpcT+ffvqp7uIAn3/+uUHABQB+fn74/PPPAQAqlQqffPKJxX04SnBwsG775s2buHnzJvLz8+GgHJeIXECSJOTn5+te41ohISEuHBUREZFn23wyFQMWx2HT8VTkKdUAHqzbNGBxHDafLP/NuSU0Ggm5hSporAhlEtMyywx0VBoJUzecwpnU+0b7sKTvxLRMvLnhJFrM+g3NZ/6GZjN/xevrT7o84AKAHQnp0Gikch8Pa6k0ElbEXTZ6TP8xXB6XVG7fZbWlz9T34ZErd40GXNq2p244hcS0TN3YdiZcL7cv4MFjqFXy+W4x6ze8ueGkrm1r2TImY+z5mFdEDr+6YlkuXbqEw4cP4/r168jNzcWUKVNQrVo1Vw7Jbf3yyy+67fHjxxs9RxRFjBkzBtOnT8e9e/ewb98+xMbGmt2HJEnYvHkzAKBp06bo2LGj0fM6duyIJk2a4Pz589i8eTMWL14MQXBQKbMFAZWPjw+qVKmC+/fvAyiq8rh9+zYEQeDURaIKQq1Wlwquq1SpAm9vbxeNiIiIyLOZGxo1Dgu0uaLLlmoxbXXN8oPmvcEf+MUhqDWSro+eTcKw9/xNs/s2Vh1VoNJY94U7QJ5SXfR4mBF4WGtHQjoWDGmlq2Yq+fz5yEUozey7ZFsl2RLWaQOdhUNbI1+l1gVk5clTqnEntwBV/byNVuiVXKC/rOrAsqq/LB1TvkoNP4XxqMbSwKysx7yicknIdfz4cbz++us4dOiQwf4hQ4YYhFxffPEF5syZgypVqiAxMRFeXl7OHqrbiIuLAwD4+/ujbdu2Js/r1q2bbvvQoUMWhVyXL19GWlpaqXZM9XP+/HmkpqbiypUrbrPeWs2aNaFQKJCRkaHbJ0mSrjqNiCqW6tWrc2oyERGRDcytCll+MAnvDWpp9RQua6/yVzJYMZe28kfbR8mpYmUFGI6qjrInXy8ZFKJoduBhjfJCoHwLQr/ywhtbwzptoOMjl8HXS2b290q79/aUu6ZaWUGvOcGtJWPy9ZLBR266QMOegVlF5fSvdtu2bXjmmWdQWFho8L/xxiqBxowZg//+97+4ffs2tm3bhkGDBjlzqG7l7NmzAIBGjRpBLjf9tDVt2rTUfcyVmJhotB1z+nGXkEsQBFSrVg1BQUHIzs5GTk4OCgsLodG4z/+6EJH1RFGEQqGAv78/AgICoFAoXD0kIiK3ZKqqwFMWFa/s7Pk8lfW9kFuoMjsk2XQiFZtOpFq1Vpc51WKvrz+Jt348bXCVv0bVA/Dx7xccGjbpBxhNwwPNrhZztb4tw3Evv9Ci4M8a5oRA5jAV3lj6fWiKfqDTNzrcrPWvtMyp0NOvFtMyJ7gd2KYWRFEwe0z9omuW+Zq3Z2BWUTk15EpPT8eIESNQUFCAFi1a4KOPPkLnzp1NXgUwMDAQAwYMwPr167Fz585KG3Ll5+fj1q1bAFDuInAhISHw9/dHTk4OkpOTLeonJSVFt11eP9qF6QDY1I8x6enpD25YuZ6WQqFA1apVUbVqVavuT0REROSJTFUVWDpVi1zDnou/m/u9YCljb+LLC+XMqdKR8CBs0PbhLCqNhEnfHUNGVoHDQyN9clHAm72jcCkjBzsS0s2+UqMAYHtCOjadcM5jZI9pmiXDG2sr9EzRD3QmdG6ALSfT7B5U7khIx4eDo1Go0eByRo5F03zNGZNcFPB857KLR+wZmFVUTg25PvnkE+Tk5KBu3bo4ePCgwULhpnTv3h3r1q3DsWPHHD9AN5WVlaXbDggIKPd8bciVnZ3tsH78/f1125b2ox+QEREREZF9lFVVYMlUrZJY/eUc1k7ns7Qte4VHKo2EN384iS0n0/DnpdsmQzlL1hBypWt3ch3WtkwA2tYNQUJqpu5x6hddE893rq97nBYMaWX2FRz1A0Fz1a3qh5vFIZ49KrMsUTK8ccTVIPtF1wQA5Baq0DQ8EAuHtrZ7H3lKNaLn7EKeUg2ZIEBdTkGG/jTf8sakvUpk84igcn/m2iswq6icGnL9+uuvEAQBU6dONSvgAh5Mi7t8uXJeGQAoquTSMmdqjnYB5ry8PIf1o7/Is6X9WMa9y4SJiIiIXEn7Zqi8qgJTbF1rhmxj7vNnajqfsefCmetJqSVgz7kHVzs2FsrdyXVudZSrCQAUclH3POmHWWWFF6Io6NZOGtimFhqHBWJF3GWLKrxMkYsCvhzVVjcd05wgzV70wxvAMd+fMgG4l1uIFrN+M/hZ9emwNth3PkP3GNqDtp3yAi6tktN8S45J/3sEAN7ccLLcn7nNI4LMDswqI6eGXFevXgUAdOjQwez7BAUVPTGWVgtVJD4+PrrtwsLCcs8vKCgAAPj6+jqsH20f1vRT3vTG9PR0i75HiIiIqHJgRdEDJQMoc6oKTLF2rRlzq4rK42nPqz3Ga83zZ2w6n7FKPEdebc9cJUO5ykIbLph6begHWeXRBhnaCq///XzGqumJJQOPsoI0W8kEAQq5aLJaDbD/96coFL02ygpbFwxphTu5BWj33h679WspY2PS/x4xd30vLWPPn6nHvLJxasilvcKdJYuA379/H4B50/QqKv01y8wJ+3JycgBY/phZ0o+2D2v6KW+9LwNW/rFGRERkD5Xxzbc7sqWiqCI+JsbeDFkbcGnpX2re2kXCtc+HuY+5u1SKOXu89nz+SlbiudPUQGum1NnClqDXGmVVbEGS4CcXAEkFqCRA0pj4kACUPF58u3i/KEnw0UhIPHMC9QU1RBQ9piIkCMUf+tva2z5yEZ0ahmJA63DUr5oCXL1W3OaD9ptLGixsK2FBjBey8iW8uu4EBGggABBLfIZePw/6g8H57euFYFi7WlCqNfASARHXgbTjQGpRfxqNBiFnEjFWptG1K+g+G7aPEv3IBA0iqvjgVlY+VBoNFCIQGeKD5Du5gKQpNS7t45T5k4R7F8IR6ivHPMU1aDSaUo+bKBT1Z/7XbdiPfn+i/tcjGD9f3AQUHPKHn5cISBrkFarROCMLW+Wlv27dWH+WULDXB94yQfe90VwCFkLCRyESJKm432sa4HtJ73tL//tJuw+G++6rHPEScRmnhlzh4eG4cuUKkpKS0LFjR7PuEx8fDwCoU6eOI4fm1nx8fBAaGorbt2+Xu2j73bt3dQGUpWtf6YdP5fWjX43FNbaIiMgVHBleOOLNt6XjteR8dwgLHPV8WFtR5KzHxNkhmqOmoulfmczaRcI3n0jFQ3VDcEZv3SFTj7mlVQuOUN73iP5za2xql1mVbZIEaNTF4YUaZ9PuY9aGw/ArfqMtgwYipKI3vtBABg0EoeR+6cExY/eRJOz/7Qaad6uPwkIlOqiPQxQ1BvcRdR+GYYAICaJg/Bz9fYKuzwe3dW/YdW/mDc+BXnsP+tMYaaN4v1D6HJgYu36bNYIUaFTND/mFKvyTetegX22IUnLsD7aNjEsoOc7S5xSFSEXT5LSBg3BOA5wtDqrsTATwqwyApRfKu1L8YUb7VQB8a+tFq9MAbAG8TRwWAcwQizeskVN8X+3978O8xySx6NNwW/q2t4wHm74Ampszrkzju4XiD6tVsGpLp4ZcXbp0weXLl7Fx40Y8++yz5Z5fWFiIpUuXQhAEdO/e3fEDdGPNmzfHwYMHcfHiRahUKsjlxp+6c+fO6babNWtmcR/G2rF3P5ZhJRcRkSdydQBlS//mvvk21UfJ/Za8kTbnfGvGa2rqTFmPkzOqW8rrw9qKokbVA/Dx7xfsGqBY8zw5Iti0ZaqPUBx6yPQCFe1nfzkg3U+DSpRwOuE06gpKyKGGWByq6J8rg7ros2DYjhwayK5p0B0ayEQNRLUG0kkNvj2twfB2tdC6ViBUKiVuZ+Yh4eBFjBfUkMn02ygKXERISP9pHW5droVQXxlUGjXkggRBW+miFxqVDJEkjQZqTdH4BEkqPqf0fW5l5SE3IwsjocEYQQNBURwenZEgO6NBskKAUqkCJA1kgoSHJA32e0m6r1UbeOjG/bMG+T8XBSEqaCATABlKT/9qBuCkrUGCMVcBrAZ8YIegwtPkFH34A+jgzACj8iwzRuRRBElyXk3n/v370bNnTwiCgF9//RW9e/cGAIiiCEEQkJCQoAtaCgsLMWbMGGzYsAGiKOLUqVNo0aKFs4bqdt5++2188MEHAIC//voLDz/8sNHz5s2bh+nTpwMAfvvtN8TGxprdhyRJiIyMRFpaGpo2bYqzZ8+aPLdZs2Y4d+4catWqheTkZAiC/d7ApKSk6KrDkj97ApGvbrNb20REZF+2BjRltWVsf1kL5WovxX4xI9vsQKmkxLRMDFgcV2aIIBOA7k3CSl1NrGeTMOw9f9Og75a1gnD82j2ojbQnE1Cq6qWs87XrqugHbJczcjDwi0NljrfklBpTY9U+TgDMfg7LukKW/niNPc7mBH/WrEFTFLaoIYO6KHQx+KyGTNDAW9BgyYhWiAxWwFuQIEINaFQlPor2Jd/Kwr6zaUhIvg2NSgVfuYTawQqk3cmGKKkg0/VX9NlL1ODRelVwPzcf1zIyIWlU8BYl1Kvqg6Y1/BHiK4OkUUGtUhUHMWpk5hXg2q0s3M7Kg6BRwUuUUN1fjppBCvh7CbqxSJIaF9LvQ5DUBsGSKJQOr/Q/dCGSwP88JCIiQymZGtT+pGi5ouTkZMuWF3JDTg25AGDEiBH44YcfoFAo8Nprr+Hpp59Gx44dIQgCtm7diuDgYBw6dAhff/01kpKSAACTJ0/G4sWLnTlMtxMfH68LtiZOnIivvvqq1DkajQYtW7bE2bNnERwcjJs3b8LLy8uifqZMmYIvv/wSAHD48GGj00r/+usvPPLII7rzv/jiC0u/nDIx5CIisow1VUu2VloZCyksDWjKC8ZKBjG2XFnKWKBkLLh5c8NJbDpu+cK+zlIyYLN8DRoJXlBDDlXx56IPL6ggF9RQCGp4QQ1R0jsuFB33FjR4oVNtdKgTCGhUSL2diSV7zkKQiu8Ptd59iu6vENR4qlUYgrwFyKGBoFEj+XYmTly5BbFk+AQNvAQ1qniLKCgshCCpTIZU+qGSrESoxSCHiCoTjSRAAoprCgG5TAZRFAEIgCACQvFnCMXz2kwdK9rOVapxJ1cFqUS7GhSdHxbkgwAfBSQJEEQRgrYd/TbLaD+nUI3TqZnQSIJuYixK9CMIAtrXD0Wgj6L4fqXblwQR206nQ1k8Tgli0fJTeNCu9kMURQzvULeoMEMQcTtXiTNpWbh8KxdKjQRRlKFB9QBER4agWqC3wZh/S7yJs9ezoJFEg8cDev1oIEAQBKgl7WTZEo8bDMelgQAvmQz/76loQBDw9s9nkK+C3gRd6Mb+4GsBFDIZPhsRU/T86j+uxZ+v3snDrC2JRcvA6bUhQdA9n6Io4oOnW6N+tQBoJODZ5X/j7p272LVkNgCGXFYpKCjA008/jR07dpRZ/aMd1uDBg/HDDz9AJrN0AnLF07VrVxw8eBByuRx//PGHLmjSWrBgAd566y0AwKxZszB79myD4/v370ePHj0AAGPHjsWqVatK9XHhwgU0b94carUa7dq1wx9//GFw9cS8vDx07doVR48ehVwuR2JiIho3bmzXr9Mw5OqHyFe327V9IqLyuOtC1bZOgwMsW4PGmmqqshirgCorGHMG/Sth5Raq0P79PWZfYUqApijIgRJe0IY6KnhBpdunKN7vJZTcp4KXNjzS3V9VfJ5+m0XnyY0GUw/CJIOwShdMFd1He25RkFWx1t0gqiy0gYC6eJUotW7SZNEKU15yGQJ8it+cizLkKjW4lavW3U//XO0baQ0E3Yd29am61QIRGlDcTnGgkF0oIfVePm5kK6HSFAUbIX7eyMhRQi3pt/GgH90bdMlwn3YMMlHEsx3rQhBkutDhdGoWLmbkolAjQSbK0DAsEG3qVEX1IF9AKGpDJQFymRyiNjwp9SGY2DZyXnFoEn/1Hr45dBUqqeTjAQiCiBe7N0KnRmGl2ygZ5FgzFpQ8Ryjdtt72ltPp+M+PCSjUaJdJf/B7uqzqWUslpmU69Kp91lYBl2Tuf0w9/VCkwdVjtcz5e8+cCm+5KOCXlzqhQXV/syuQ9cdk69ehZWk7uYUqNJ/5G1SZt5D65TgADLlssmzZMsyfPx+XLl0yejwyMhJvv/02Jk2a5OSRua8TJ06gU6dOyMvLQ0BAAN5++2306NEDeXl5WL9+Pb7++msAQFRUFI4ePWpwtUTAvJALAKZPn4558+YBAGJiYjBt2jQ0bNgQly5dwocffogTJ07ozps7d67dv06GXETkKs5cY8fWRcUtnQZnap2i8u5jz2oqW8ighjeUUEBZ9Fko+uwNFbxRCIWgMjwOJbwFJRTa41AV3y7ap9ALlbyhhLeghkzSD5qUUBQHQ4b7VLrAyEvggixUOWikB6GKCqJuW/+zCjKD89S6D1mpfRpJhBqC7phhew9CEbV2+XCpZLgjoH71IDzSqDriLt3F+Zs5eiFO8f0l0eA+LWsFo39MbRSqgfd3ni/VnsagH/2xlB6P8fsIZXwN+qGUWCKwKnosHqz0VbRfEEWoNALkcjnyVZKu2sUYuShgy8ud0TQ8sNR/hBgLKno0qY595zMsDjBK/t4qK6gojy2hgyM4OtSxJ2eO1dHra9r6dZgbQG15ubNNj40loZw1Y7I0SDO1jmaLWb+Z9R91vl4y/DOnDwCgxazfkHX7BkMue0pMTMTRo0dx8+ZNqNVqhIaGIiYmBg899JBd13mqKLZu3YpRo0YhM9P4pRWioqKwfft2NGrUqNQxc0MujUaDF154Ad98843JcTz//PP4+uuvi0th7csw5OqLyFd32L0PIqpczPkjrbw/YCxd68ncKXjlBWnWVk3ZmwANvKGEDwqLPoRC+GhvC0X7vHXHlLrb3mWETd7aD93x0iGW9n4yTkEjB1BJRSGMCg8+y+UKBPh5Iz1TBaX0IMApeZ4aMr37F4U1St15Mqglw5CndDAkFAdDhsf0PzQltlWSrMQ+wWiYBEGG959ujVoh/lB4eUOUyQBRBghF989TSxi45G/kKGEwngdBlmE/Nl63yy34yEUc/V8v+MhliJ6zy+yKTVco+WbWUWsR2iPAKBlUmPMfIfYIHRzFXSu5jfGksZbF1q/DXlVh5bEklLNmTGXdx5wlF7RVWWZ/Pe/2gZ9Cjjc3nMSG/acYcpFrXb16FZ999hm2b9+OlJQUKBQKNGrUCM888wxefvll+Pn5Gb2fuSGX1o4dO/D111/jyJEjuHXrFqpVq4b27dtj4sSJ6Nu3r72/LB2GXERUHntd1U7LnP9BM8XSxcNNMfZGpfw3CxK8oYQvCuCHAvgKxZ9RAD+h6LMvCuCrH0Lph1N6AZUPCuEt6O03OFYURFHlUSjJoIIcKsigLA5tlJBDJckM9qkgL9qW9M4pPq49Xz/4KRUWlQqZZHqhUlGoYzRc0n5IJu6n7aec9itCeFOSOW/qLH0z5Oslc2klpz35eslQPVCBa3fyXD0Uo0w9f6beYDes7m+yQteeb/DLY8kFQpw1Jqo83LG6zZoxGbtPdK0qOHbtrlnrnVpayaX9W7nf3F9w9YuxABhyETmEYcj1OCJf3eniERGRPVj6P8mWrCdl6VXt9NdhsvYKbuYSoYEf8nXBk58ueDIMpB5s58MXhSX25+tCK7/i0KpoO58VTm6sQJKjEF66wEcJOQql4s/Ft/X3afcXQg5l8T4VZCgs/qyCHEqDkEkvhJL0AqeSx6USgZTuPkZCrAoc/lREdav64WZWgcVv6ix9M5QwKxaFGo1bVZdWNNZMGTx3PcspU7Ws4UnT/6jicMfqNlsuDmTOFZT1X+OWrMm1YEgr3bhW/n4UEx7vAIAhF5FDMOQicn9l/cK2ddqeqf1lrSdlDQGAQi6iQFVyIW4JvihAAPLhL+TBH/m6bYN9xZ/9kYcAId/kPj+hwC7jJeMKJRkK4YUCeBV9lrTbchRAgUJJ/uBY8edC6UHIVGgkaNKFTeXs0wVVxWFWoV6YpYIMDIocpyJUFNlDyQDKkjd1tix0bGm1QXlkogC1RiqusvLGtTu5FrfhLPq/O+z1fXj0f4+hqp+3xW/K7bVYtSO5Y+hA5CksfY2bMzPB2EWAOkcIWD7lcQAMuUzatGkTBg8ebO9mddLS0nDt2jV07NjRYX2Q6xiEXJ8+jsjXGHIROZo9pv8BsGhxdEcQoUEAchEk5CKo+HMgirYDBWOfcxCgC6+Kwil/VkeZTSMJyIcC+fAq+iwpUACF3nbR/kLIUSAZBk0F8EJByQDKjICqKMQqPg45JNh/bUhjfL1keLRhKPZfyHDZlSDNoQ0LKmoINDimFt4b1NLsdYouZeToghiZIEBdAf9vV7uuisX3s8OCzcZ+R1izPpP++lPmVCfZU8lqOHOmBmmrgO1R2aY/bcgS1k5NIiLPYO1rvKy1vbQ/AkoeqmhXV7T8N6IZhgwZglatWuF///sfhgwZYrd2k5OT8cEHH2DlypWYPn06Q65KoeL9MUrkaI6qsjIWWOUp1dh0PBW/FE/10/+lmadU48iVuxaNXYAGgchFiJCNYGSjipBjNKAKEnIQiLxSnwMF91xjxZkKJDny4I1ceCNP8kYevItDp+IgCsVBlKR4EFDptouPS4bBlXZ/XokgqxByuKpaSb+awpy1aaylXazaTyG3+mpiMgFoWzcECXrrtpX1RtrY+eYEbJYsVm2KqT+ALeGoMEkuCpjQpYEu0BnYphYahwWWOyVKOyXDnGkfljD2PDnq+9AUXy8ZfOQyq+7bPCIIC4e2LnftpLKmlomiYBCwadvUnwZjzvpMLWtVMXtc9nYzq6BUNZw5U+3K+z40tyKtX3RNq8KnfJXa7MXz85Rq5KvUVoWhROQa1r7GTf1M8oT/qLMXh1RyNW7cGJcuXYIgCKhTpw6effZZPPvss2jRooXFbeXk5ODnn3/G999/j927d0OlUkEul2PlypUYOXKkvYdObsCwkqsPIl/71cUjIvIM7lRlpQ2rgoUchCALwUIOqiBbF14FC8UfyEawkINg3Tk5ECt4FZVGEooCKHgjV/I22M6DD/KgKN4uOpYr+SAP3iX2+xick6e3rYZ1b3jtxdqAxhLGqin0122z5zQqc6dq9YuuiR5NqmPf+QyTb4zLq3op73xrrtZkzVgBGNzH3NBKW2VlTpgkE4AeTcNw6OJtiypoTC1Wbe6UKGtCypLKe54Ay78PS75uzH3M7TENzRlrJ9ljAWZHMlUNZ+lUO2eul8VKLqKKzR6vcf2fSf/34ymTUx8rWiWXQ0IupVKJTz/9FPPnz8ft27chCEUPduPGjdGxY0e0b98eMTExCAsLQ0hICEJCQpCXl4c7d+7g7t27uHDhAo4cOYL4+HjEx8cjPz8f2mEOHjwYc+fORVRUlL2HTW7CIOT6JBaRr5t/9R+iisicqYRl/U+5PSozRGgQgiyECpkIFTJRFVkIFe4jVMhCKO6jqlB8DJmoKmSiCnIqzJQ/pSRDDnyQDV/kSD5F25IvcuCDHPgiu3hfjuSLbP3PxednwRe5ko8uzCqAFyrKWk3GqqmsCWhMVcSUFxqVxV7TqCyZqlXefnPHWhZrAwlrxmrtwreA+ZdPtzT4swdrw5OSVX3msDbYtOYxt5Uz1k6yZQFmcy8QMjimFt4d2ALt39/j8gDImmDaEp6wJhcRWc9er/HyAjOGXBbIzs7GkiVL8MUXXyA5ObmoQ8H8XyDaoXl7e2Pw4MF47bXX0KFDB4eMldwHQy6qrKydSmjt2jveKER14R5q4C7ChHsIFTJRTbhfIsAqCq1CkO1xFVZqSUAm/JEl+RZ/9kMm/JAFP2RK+p99kVMcWhUFVg/Cqhz4uHUoZek0OGPrFFlb6VRWNZUpllYuadnzzXdFucy9sxdztmcVmTUhpSNYGp7YMygw5+tzdEDiaSxdQ8xdAiBHBrf2WFeNiNyXvV7juYUqNJ9p+j01Qy4raDQa/P7779iwYQP27duHK1eulHsfHx8fPPzwwxg4cCDGjBmDqlWrOnqY5CYMQ67eiHx9l4tHRGQ/pqa12HMqoS/yESbcQxjuoYZQFGAVfdxFWHGgVUO4iyqC+169CgAKJC/cgz8yJf+icEobUpX4nGlw+0GolQtv2BpOuWKNHVMsqZoy502VufexpZrKFHe62hYvc28+e1eRuQt3Dgr4/WnIkuDP3Z5XR70OGIYSVWz2eI2zkssJUlNT8eeffyIlJQUZGRm4c+cOfHx8UL16dVSvXh3R0dFo164dvLy8nD00cgMMuciT2HpVwkbVA8wOTwRoUA2ZqCncRkTxR7hwpzi8Kgquqgv3EORmi68XSHLcQwDuSQEPPksBuIsA3JcCcA/+uCsF4n7x53uSP+4hAPnwdtmY9a/gZq81dsoLjcoK0qypmgJsmxrkyGoqd1TRvz57qoiPlbsHBRXxMbeWJcGfuz+v9sIwlKhis8drvKzqVoZcRA5mEHJ93AuRb/zu4hER2f+qhMYEIFcXXkUIt4vDrFuIwB1ECLcQLtyBt6By0lds2j3JH7elINxGEO5IQcXbgXrbQUUhlhSIe/BHnh0qqqxVMmhyxDpM9pyCxzcqRK7D159nMTf4q0zPK8NQoorNltd4WdWtDLmIHIwhFzmaOQstlxVmWTuVUA4VIoTbqC3cRJ3ij9pChu52iJBtzy/TbHmSAjekENxClaKQSgrEbVTBHSkQt6UquI3iz1Ig7iIQKrjmEuSm1pOy9Kp2jlqHyVHrRvGNCpFz8fVXMfF5JaLKzlR1K0MuIgczDLkeQ+Qbu108IqooTFVfTejcAADsEmb5oAD1hetoIKSjrnBdF2TVEW4iQrjl1CsO5kjeuCGFIAPBuCGF4Kb04PNNFH+WQpAFX9haaVXWFRwtnbZnzhpQtk6pq0z/s09EREREBBj/G7hLhIhlU/oAYMhF5BAGIdfCnoh8c4+LR0SeytzqnbICGmMEaFBLuI0GQlpxoJWGBkI6GojpqCXctuNXYJxSkuEGQpAqVcN1qSpuSCG4IQUjQwrBTb1AKwe+du/bVGD1fOf6AGDXK+e56+XsiYiIiIg8mf7fwGlpqQ/ef1eAkMs1c06IiByoZMVWeeswmQq35FChvnAdTYVraCymoGFxmFVfuA4fQemw8WdIQUiXQpEmVUO6VBWpUrXi20UfGQiGBqLZ7dmzyqq8wGrh0NZYMKSV0WOiKMBPUfrXjqX77ckZfRARERERuZOK/DdwxfyqiKhSKa9iq0ClKacFCeG4g6ZiMpoK19BETEZTIRkNhDSHLPSeKfkiWQrDteKPZN1HdaRK1VAAhdVtCwAUchEFKo3NVVZPto6wOLAq7xgREREREZGj8F0IEXmM8haFN+fKed4oRHPhKpqLV9FESC4OtK6hipBr17GmSqG4rAnHNamGLsjShlr34Q9HXG1Qu2h6/1YRDq+yIiIiIiIicjd850LujUvGEYwvGG9sUfiSFVsyqNFYSEUr8RJaC0loJV5CUyEZXoLaLuPKlHyRJNVEkhSBy5pwJEkRSJJq4opUA3nwsbl9a6cSAmCVFRERERERVTp8p0NEbqVktZaxS93mKdU4cuVuiXtKqCvc0IVZrcQktBSuwE8osHlM6VJVnNfUxgUpsijU0tTEZakmMlAFtlZkyUUBb/aOwqWMHLtNJSQiIiIiIqqMGHKRm2MlV0VV3tRDXy8ZHmkYigMXMgyqtbR8UIAY8SLaC+fRTjyPVmISgoUcm8aUJfnivFQb5zW1ca7483mpNu4jwKZ2gbKvSqitvuJUQiIiIiIiIuvx3RG5N2ZcFY65Uw/zlGrsPXdTd7sKstFBPIf24nm0F8+jpXDZ6mmHGknAJSkCZ6U6OFccZJ3T1EEqqsHcyix7X5UQYGhFRERERERkC76bIiKHKu/Kh8anHgJ+yEd78TweEf9BJ/EMWghXIQrWpZ7Jmuo4JTXAaU0DnJYa4oymHrLhpzsuCijKtow0X14Flj2vSkhERERERETW4zstcnMs5fJUll75UIAG0cJl9BBPopPsDGKEi1ZVamVIQTitaYjTmgY4JRV9voMgvX4AhVwEVBqDwAoAVsRdtrgCi1MJiYiIiIiI3IPL3oFlZGQgKSkJ169fR05ODry8vBAcHIw6deqgUaNGkMlkrhoaEVmhvIqtklc+BIBgZKGrmIDuspPoKp5GNSHToj5VkogzUj0c1TTBMU0UTmkaIg2hMDXlUC4KWDi0Nfq3ijAaWC0c2poVWERERERERB7Kae/YcnJysHnzZuzcuRMHDhxAamqqyXO9vb0RExOD2NhYDBo0CK1atXLWMMndSKzkcneWVmzVFa6jj3gEsbJjiBH+hcyCKYh5kgLHNY1xRGqCeE1TnNQ0Qi58jJ6rrdgqKFGxpV3k3VRgxTCLiIiIiIjIMzn8ndyJEyfw+eefY+PGjcjNzQUASOUEF/n5+Th8+DD++usvvPvuu2jRogVeeukljB49Gn5+fmXel4icZ/PJVDMqtiQ0F66ij+wo+ohH0FRMNrt9pSTDSakh/tS0wJ/qljghNUIhvMq9X3kVW0RERERERFTxOCzkOnHiBGbMmIGdO3cCeBBshYeHo0OHDmjbti3CwsJQtWpVhISEIC8vD3fu3MHdu3dx4cIFHDlyBKdPn4ZSqcSZM2cwZcoUzJgxA2+99RZeeeUVeHt7O2ro5FZYyeVO9KcknrueVSrg0hclJOMp2SE8KR5GHTHD7D4uaWrigKY1/tBE44imKXLga/Lc8haFB0xXbBEREREREVHF4pB3f+PHj8eaNWug0RRVdDz00EMYOXIknn76adSpU8fsdgoLC/HHH39g7dq1+Pnnn3Hr1i1MmzYNS5YswerVq9G5c2dHDJ+ISig5JdHXS4bqgYpSAVctZKC/7DAGyg6hmZkVW3mSAn9qWmC/pjX2a1ojWapR6hyZAPRoGoZDF29btCg8ERERERERVR4OCbm+/fZbKBQKjB07FlOnTkVUVJRV7SgUCvTq1Qu9evXCV199hY0bN2Lu3Lk4d+4c9u7dy5CrMuCaXC5R3iLyeUo1rt3JAwD4IR9Pyg5jiOwPdBDPm9V+hlQFv6vbYpemHQ5rmqMACpPnaqceDmxTy+IrHBIREREREVHl4ZB3hVOmTMG0adNQu3Ztu7Xp7e2NUaNGYeTIkdi4cSPUarXd2iZ3xpDLmSxZRL6lkIRnZXsxQPYnAoT8ctu+pqmO3zTt8Zu6HY5LUdBANDhuztRDhllERERERERkikPeLS5evNgRzQIABEHA0KFDHdY+UWVlziLyAcjFQNmfGCHbi5bilXLbTJOqYov6UWxVP4p/pLoouuZh0b/eJq58yKmHREREREREZA2WRJCbYyWXI2kDpcsZOWUuIh8pZGCc7FcMk+1HoJBXZpv3JT9sVz+MzerOiJeaQCpRsVXelQ9ZrUVERERERETW4DtJokqo5LREmSBAbWT9szbCRUyQb0dfMR4ywXTgqJYE7NO0wQZ1d+zXtEEhvACUXbEF8MqHREREREREZD98h0lujpVc9mZsWqJ+wCVAg1jxGCbIt6O9eKHMtlKkatig6o4N6m64jlCDY+VVbBERERERERHZk9uFXDdu3MCkSZOwe/duFBQUoH79+mjTpg1iYmJ0n2vUqOHqYRJ5FPOmJUqIFY/iDflPaCZeM9mWWhKwW9MW69Q98YemFTQQUbeqH3yzCkwuGM+KLSIiIiIiInI0t3vn+fLLL2Pz5s2oU6cO6tevj6SkJGzcuBEbN26EIBRVgdSoUQMxMTHYvn27i0dLDmdkCh2Zz7xpiRJ6iCfxpnwjostYTD5b8sEGdXd8o34cKVKYbr9cFPDlqLZoGh7Iii0iIiIiIiJyGbcLufbu3Yt27drh0KFD8PIqWtfnzp07OH78OE6cOIHjx4/j2LFj+O2331w8Uue6cuUKtm7div379+P06dNITU2FRqNBtWrV0K5dOwwfPhxDhgyBXG7bU7p//3706NHDrHNnzZqF2bNn29QfOU550xIBCV3F03hT/iPaiJdMtpMmVcUqVR+sV/dEJvwNjmmnJLJii4iIiIiIiFzN7d6RSpKEnj176gIuAKhatSp69eqFXr166fZlZ2e7YnguMWPGDLz//vuQjFQ1paamIjU1FZs3b8bHH3+MH3/8EXXq1HHBKMmdJKZllnm1xMZCCmbJv0Vn2T8m27igqYUvVAOxXdMRKsjLXUSeiIiIiIiIyJXcLuR69NFHcfXq1XLPCwgIcMJo3EN6ejokSYK/vz8GDRqExx57DI0bN4aPjw/Onj2LRYsW4ciRIzhy5Ah69eqF48eP2+Xx+eabb9C+fXuTx8PCwkwesxtOV7SIdu2t5QeTjAZcQcjB6/KfMEa2C3JBY7SNS5qa+Ez1NLZpOkIDEQAXkSciIiIiIiL353Yh11tvvYUnnngCV69eRd26dV09HLcQGhqKDz/8EJMnT0ZgYKDBsbZt22LEiBF49tlnsWHDBvz777/4+OOPMXPmTJv7rV+/Plq2bGlzO+R4JdfeKkmEBsNk+/B/8g0IFbKMtnFFUwOfqQZji+ZRQJRDA4mLyBMREREREZHHcLt3q+fOnUO/fv3w+OOPY8OGDYiOjnb1kFzuww8/LPO4TCbDkiVL8Msvv6CwsBA//vijXUIu98BKrvIYW3tLX4zwL/6f10q0NLGofJpUFZ+ohuBndWeoIIdcFPDLS53QoLo/K7aIiIiIiIjIY7hdyDVp0iQIggBJktCmTRs8/PDD6NWrFx566CG0bdsWtWvXdvUQ3VJoaChatWqFo0eP4tIl04uIU8VS1tpb3ijEG/If8aJsO0Sh9PECyQtL1U/gS9UA5MEHwINpiS1rVXH42ImIiIiIiIjsye1Cru+//x6nTp3CyZMncfLkSfz111/466+/IAhF1SShoaGIiYlB27ZtMXfuXBeP1r0UFBQAKKrsqjC4Jlcp2nW3fOQyLI8zvvZWSyEJH3t9iSgx1Wgbv6rb4z3VSKQLNaA2MS2RiIiIiIiIyJO4Xcg1fPhwDB8+XHf7xo0busBL+7Fnzx7s3r2bIZeemzdv4uzZswCAZs2a2aXNd955BykpKbh+/Tr8/PxQr149dO/eHZMnT0ZUVJRd+iDzlVx3y0cuQlki4PKCCi/Lf8FLsl+MLix/QVMLc1RjcEgTDbkoYDOnJRIREREREVEF4XYhV0k1atRAnz590KdPH92+vLw8nD592oWjcj8LFiyASqUCAAwdOtQubf7555+67cLCQl3IuGjRIsyYMQOzZs3SVdhZIiUlpczj6enperdYyQUYX3crX2UYYjURrmGh11dG194qkOT4VDUEy9T9dOtucVoiERERERERVSRuH3IZ4+vri4cfftjVw3Abf//9Nz799FMAQGRkJCZPnmxTezVr1sTgwYPRuXNnNGjQAHK5HNeuXcO2bduwevVqKJVKzJkzB4WFhVZV01m0rhozrjLX3dJ6RrYf78lXwltQljp2RlMPbyon44JUG75eMgzktEQiIiIiIiKqgARJcu6iR2lpaTh27BguXbqE9PR0KJVK+Pn5oU6dOmjbti1iYmIgiqIzh+TRbty4gXbt2iElJQWCIGD37t3o2bOn1e3l5ORAoVDAy8vL6PH4+HjExsbi/v37EAQBJ06cQOvWrS3qw5Lqr+T32yHy7SMWtV/RvLnhJDYdN762lhdUmCFfgzHy30sdU0oyfKEeiMWqpzAgpi7eG9SS0xKJiIiIiIhIJyUlRVeIkpycjMjISBePyDZOreSqX78+rl27VuY5AQEBGDBgAMaNG4fHHnvMSSOzD2um7pW0cuVKjBs3zqxzs7Ky8MQTT+im/82bN8+mgAsA/P39yzzeoUMHLF68GKNHj4YkSVi8eDGWLVtmUR/JycllHk9PT0eHDh2Kb1XeUi6NRkJuoQo7E64bPV4d97BE8SnaixdKHTuvicRU5SSckRpALgqY0KUB/BQeWbhJREREREREZBanvuu9evVquedkZWXh+++/x/fff49HHnkEK1asQJMmTZwwOs+Sn5+PgQMH4tixYwCA//u//8Nbb73llL6HDx+Ol156CZmZmThw4IDF9/f0ZNjRSi4wb0yM8C++VHyKcOFuqWMrVX0wTzUCBVDo1t7i1EQiIiIiIiKq6JwacrVo0QLt2rVDixYtEBkZiWrVqgEomiJ37do1nDx5EgcOHEBSUhIA4PDhw2jbti22bduG7t27O3OoVtFe3dAWNWvWLPcclUqFoUOHYt++fQCACRMmYMGCBTb3bS65XI6oqCgcPXoUqanGp9HZjXNn07qcsQXmSxoh24M58lVQCIYBWL7khZnSRGxQPQpfLxme5tpbREREREREVIk4NeRKSEgw67zTp09j0aJF+Pbbb5Gbm4unnnoK586dQ3h4uINHaJumTZs6vA+NRoPRo0dj69atAIBhw4Zh6dKlDu+3JHtMzSRD5S8wL+G/8nWYJN9W6kiKVA0bGnyAeaOfwWyVmmtvERERERERUaXjliu8t2rVCsuXL8fevXsREBCArKwsp1YqubOJEydi/fr1AID+/fvju+++c/pC/SqVChcuFK0DFRER4eDeKk8l1/K4JJMBlwAN3pN/YzTgilO3wCDl+3i89+MQRQF+CjkDLiIiIiIiIqp03DLk0urSpQvmzJkDSZKwbVvpN/eVzZtvvonly5cDAB577DFs3LgRcrnzFxP/4YcfcP/+fQBAt27dnN5/RaTRSCYXmJdBjYVeX2GUfE+pY0tVT+B59XT8b2gXTkskIiIiIiKiSs2tQy6gKOgCyr8iX0U3e/ZsfPLJJwCARx99FJs3b4a3t7dFbVy5cgWCIEAQBKNrnN29exf79+8vs434+Hi8/PLLAIqmLE6ePNmiMViskqzJla9SG11kXgElvvBahMGyOIP9aknAdM1kXGg1DT+/3A0D29Ry1lCJiIiIiIiI3JLzy4AstGPHDgCAr6+vi0fiOp9//jnmzJkDAKhVqxbmz5+Py5cvl3mfJk2awMvLy6J+7t+/jx49eqBVq1Z46qmn0LZtW9SsWRMymQzXrl3Dtm3bsGbNGhQWFgIouqJj27ZtrfuiSEejkaDRSPD1khkEXb7Ix1KvT9BVZriWXaEkg3rQcrzfahCnJRIREREREREVc4uQa/To0Th69CiioqJQq1Yt+Pv7Izs7G0eOHMGJEycgCAL69Onj6mG6zE8//aTbTk1NRefOncu9z+XLl1GvXj2r+jt9+jROnz5t8rhMJsOMGTMwc+ZMq9qnIolpmVgel4SdCdeRp1RDpreYfyBysUKxAB3E8wb3yZe8sLruXLzYZrCzh0tERERERETk1twi5FKpVDh//rxuMXMtSZIgCAJ69eqFRYsWuWh0lUdERAQ2btyIw4cPIz4+Hqmpqbh16xby8/NRpUoVNGnSBN27d8eECROsDtAsVkGnK24+mVrqSorq4q/VG4VYqZiPdqLh6yFL8sWLqv9gxuPDnTpWIiIiIiIiIk/gFiHX1KlT0aJFCxw6dAgHDx5Ebm4uBEFA165dsWDBArRv397VQ3Sp8tbJMle9evUglREaKRQKDBkyBEOGDLFLf2RcYlpmqYBLS4AGn3gtKRVw3ZP8MV41HeOeGcwF5omIiIiIiIiMcIuQq127dmjXrh0AoKCgAFu3bsUXX3yBAwcOoFu3bli8eDGee+45F4+SXKPiVXItj0syGnABwHT5OvSTxRvsy5CqYFXDT/F+r14MuIiIiIiIiIhMcOrVFcuqItLy9vbGkCFDsG/fPqxfvx6CIOCFF17A3r17nTBCIsfSaCTsTLhu9Nho2S68KN9usE/yDkTopB34zxhWcBERERERERGVxakhV0xMDA4ePGj2+UOHDsW8efMgSRLmzp3rwJGR26pga3Llq9QGV1DUekw8htnybw32KSUZCgZ/C7FmS2cNj4iIiIiIiMhjOTXkOn36NLp3746BAwfi2LFjZt3nkUceAQDEx8eXcyaRe9NoJGg0Eny9ZAb7o4UkfO61GDLBMNCbIU2EonFPZw6RiIiIiIiIyGM5NeRq1qwZJEnCtm3b0KFDB3Tr1g2rV69GZmamyfvs2LEDgHlTHaki8vznPTEtE29uOIkWs35Dy9m7UKjS6I5FChn4RrEAfkKBwX0+VQ2GsuVwiKLg7OESEREREREReSSnLjx/6tQpLFiwAO+//z5yc3MRFxeHuLg4TJgwATExMYiJiUHDhg0RHByM7Oxs/PHHH9i6dSsEQUDTpk2dOVRyG54dcm0+mVrqSorq4sDWF/n4xms+qgv3De7zo7orFmuGYEvn+k4dKxEREREREZEnc2rIJZfLMX36dDz33HN4//33sWzZMhQUFEClUuHo0aM4evRoqftIkgRBEDB58mRnDpXIZolpmaUCLn3/k69FlJhqsC9O3QIz1C9g4dA2XGieiIiIiIiIyAJOna6oVaNGDSxatAipqalYsGABOnToAFEUIUlSqQ9vb2/MmjULzz33nCuGSi7nuZVcy+OSTAZcvcRjGCnfY7DvglQbO5svwE8vd8fANrWcMUQiIiIiIiKiCsOplVwlVa1aFVOnTsXUqVORlZWFf/75BxcvXsT9+0XTt2rXro1u3bqhSpUqrhwmkcU0Ggk7E64bPVYd9/Ch19cG+yQvfzR68Re8X72RM4ZHREREREREVOG4NOTSFxgYiI4dO6Jjx46uHgq5Ew+94EC+So08pdrIEQnzvZYiVMgy2FvY+314M+AiIiIiIiIisppDpivOmzcPu3fvxr179xzRPJHb85HL4OslK7V/tOx39JCdMti3W2oPr7ZjnTU0IiIiIiIiogrJISHX22+/jT59+qBatWoG+zds2IALFy5A8tDqHHIFz/xeEUUBfaPDDfY1ElLwjnytwb6bUjD+aDIDoswly+MRERERERERVRgOma4oiiI0Gk2p/cOHD4cgCPD19UV0dDRat26NNm3aoHXr1mjVqhX8/f0dMRwip9FoJOSr1PCRyzChcwNsOZkGlUaCF1T4zOsL+AhKg/OnqSbhP91jXDRaIiIiIiIioorDISFXdnY2Tpw4gaNHj5Y6JkkScnNz8ffffyM+Pl63XxAENGjQQBd6aT9HRkY6YojkKTykkCsxLRPL45KwM+E68pRq+HrJ0Dc6HG/2jsLHv1/AVHEjWohXDe6zSv04nnpmDJpHBLlo1EREREREREQVhyA5ce7g+fPncerUKd3HyZMnkZaWZjggQTC4HRISgtatWyMmJgYfffSRs4ZKLpSSkoLatWsDAJJnNUHk7HMuHlHZNp9MxdQNp6DSlH4pyUUBH3bIx6CTz0HUS+zSvevj3sjf0KxODWcOlYiIiIiIiEjH4P13crLHFxo5NeQy5vbt2wah16lTp3D27FkolYbTugRBgFpt7Gp1VNF4UsiVmJaJAYvjjAZcRST8opiJNuKlB3tkCggv7AXCo50zSCIiIiIiIiIjKlrI5ZDpipYIDQ1Fz5490bNnT90+pVKJs2fP6kIv7QdVQm4+XXF5XFIZARfQXzxsEHABgNDjHQZcRERERERERHbm8pDLGC8vL7Rq1QqtWrVy9VCITNJoJOxMuG7yuDcKMc1rveHOqg2BjlMcPDIiIiIiIiKiykd09QCIyua+pVz5KjXylKan0I6V/YZI4Zbhzt5zALnCwSMjIiIiIiIiqnwYchFZyUcug6+XzOixEGTiZflmg31SnUeApk86Y2hERERERERElY5DQq5NmzY5olmdtLQ0/PXXXw7tg9yEa6+LUCZRFNA3OtzosdfkmxAk5BrsE/rMBUpcPZSIiIiIiIiI7MMhIdeQIUPQpk0b/Pjjj3ZtNzk5GVOmTEHDhg2xa9cuu7ZN7sp9Qy4AmNC5AeSiYXDVQEjDSNkeg333Gg0Caj3kzKERERERERERVSoOCbkaNmyI06dPY9iwYahfvz7eeecd/PPPP1a1lZOTg++++w79+vVDw4YN8dVXX0GtVqNhw4Z2HjWR+TQaCbmFKjQND8TCoa0Ngq7/ytfBS3iwVpdaVCD4yXddMUwiIiIiIiKiSsMhV1dMTEzEp59+ivnz5+Pq1auYN28e5s2bh8aNG6Njx45o3749YmJiEBYWhpCQEISEhCAvLw937tzB3bt3ceHCBRw5cgTx8fGIj49Hfn4+pOJpa4MHD8bcuXMRFRXliKGTu3Gz6YqJaZlYHpeEnQnXkadUw9dLhr7R4fh0WBvsO5+BjIQ9iJUdM7iP7NGXgOA6LhoxERERERERUeUgSJLjUoTs7GwsWbIEX3zxBZKTk4s6tGBNIu3QvL29MXjwYLz22mvo0KGDQ8ZK7iMlJQW1a9cGACT/ryEi/99FF4+oyOaTqZi64RRUmtIvGbkoYOEz0RgQPxJC+qkHB/yqAa+eAHyCnDhSIiIiIiIiovIZvP9OTkZkZKSLR2Qbh15dMSAgAG+99RYuX76MnTt3Yvz48ahbty4kSSr3w9vbG926dcPHH3+M1NRUrF27lgFXpeQelVyJaZkmAy4AUGkkHPjpS8OACwC6/5cBFxEREREREZETOGS6YkmiKKJPnz7o06cPACA1NRV//vknUlJSkJGRgTt37sDHxwfVq1dH9erVER0djXbt2sHLy8sZwyMq1/K4JJMBVxEJL4pbDHdViwLajnfouIiIiIiIiIioiFNCrpJq1aqFZ555xhVdk6dxgzW5NBoJOxOul3lOR/EsmorJhjt7zQZkLnmJEREREREREVU6Dp2uSFQR5KvUyFOqyzxnrOw3wx1VGwJRfR04KiIiIiIiIiLSx5CL3JzrK7l85DL4eslMHo/ALcSKRw13dngREPnyIiIiIiIiInIWh70Lf/fdd/Huu+/i1q1bNrVz5coV9OzZE4899pidRkYexfUZF0RRQN/ocJPHR8l3QyboDVQRALR51gkjIyIiIiIiIiIthy0YNHv2bAiCgCFDhqBatWpWt5OTk4P9+/dDEAQ7jo7IMhM6N8CWk2mlFp/3RiFGyPYantx6BK+oSERERERERORknE9Fbs4NSrkANI8IwsKhrSEXDcPWAbI/ESJkG57c4UUnjoyIiIiIiIiIABeHXElJSVCry17Qm8hdDGxTC1te7oynH4osXqNLwnj5LsOTGvQAqke5ZHxERERERERElZnLQq7t27ejTZs2ePTRR3Hx4kVXDcNjdO/eHYIgmPVhL3/++SdGjRqFunXrwsfHB+Hh4ejTpw/WrVtntz7KJblHJZeWtqLrnzl9cO6FYDQXrhie8PAkl4yLiIiIiIiIqLJzSci1aNEiPPXUU8jJycHRo0cRExODFStWuGIoZMLs2bPRpUsXrF27FteuXUNBQQFu3LiBXbt24dlnn8WTTz6J/Px8Vw/T4TQaCbmFKmhKrMUligJ8ji83PDmkHtC4t/MGR0REREREREQ6Dlt43hhJkvDqq69iyZIlkCQJQUFByMrKQk5ODl588UXs3LkTy5YtQ0hIiDOH5VHatWuHlStXOrSPpUuXYs6cOQCAhg0b4u2330Z0dDTS0tLw2WefYd++fdi+fTuee+45fP/99w4di6vW5EpMy8TyuCTsTLiOPKUavl4y9I0Ox4TODdA8IgjITAMStxjeqf0LgChzyXiJiIiIiIiIKjunhVw5OTkYNmwYdu7cCUmS0KJFC2zduhXnz5/HuHHjcOPGDfz888/4+++/sXr1avTo0cNZQ/Mo/v7+aNmypcPav3PnDqZNmwYAqFOnDv766y+Dq2M++eSTGDRoELZu3Yp169bhxRdfRPfu3R02HlfYfDIVUzecMriSYp5SjU3HU7HlZBoWDm2NgXdWApLeenJefkDMKBeMloiIiIiIiIgAJ01XTElJQadOnXQBV2xsLA4dOoR69eqhT58+SEhIQN++fSFJElJTU9G7d29MmzYNKpXKGcMjPcuXL8f9+/cBAB9++KFBwAUAMpkMS5YsgUxWVLG0YMECB4/IuZVciWmZpQIufSqNhOkbjkJ1pEQ1XathgG+w4wdIREREREREREY5POQ6cuQIHn74YSQkJECSJEyZMgXbt29HUFCQ7pxq1aph+/bt+PTTT+Hj4wONRoOPPvoIHTt2xLlz5xw9RNLzyy+/AACCgoIwePBgo+dERkaiV69eAIA9e/YgKyvLgSNybsi1PC7JZMCl1QeHIc+7Zbizw4sOHBURERERERERlcfhIdfzzz+P9PR0CIKATz/9FIsXL9ZVAZX06quv4u+//0bz5s0hSRJOnDiB4cOHO3qIVKywsBDx8fEAgEceeQQKhcLkud26dQMAFBQU4OjRo04Zn6NpNBJ2Jlwv97wx8t8Nd9TrAtRo7qBREREREREREZE5HB5yaTQaBAYGYsuWLXj11VfLPT86OhrHjh3D5MmTIUkSNBqNo4foUc6dO4eHH34YwcHB8PHxQWRkJAYOHIjVq1dDqVTa1PaFCxegVhetM9W0adMyz9U/fvbsWZv6LZPkvEqufJUaeUp1medE4BZixIuGOx+e6MBREREREREREZE5HLrwvCRJqFOnDrZt22bRYune3t744osv8Pjjj2PChAnIyMhw4Cg9y40bN3Djxg3d7dTUVKSmpmLLli348MMP8eOPP6JZs2ZWtZ2SkqLbjoyMLPPc2rVr67aTk5Ot7seY9PR0i9qzFx+5DL5esjKDrt6yY4Y7fEOAqL4OHhkRERERERERlcdhIde+ffsAAC1btkRoaKhVbfTv3x9nzpxBYmKiPYfmkURRxGOPPYZ+/fqhdevWCA0NRVZWFo4fP46lS5fi7NmzSExMRI8ePRAfH486depY3If+2loBAQFlnuvv76/bzs7Otqgf/YCsfM6r5BJFAX2jw7HpeKrJc2LFElMzo/oCMqddpJSIiIiIiIiITHDYu3Ptmk22ql69ut3a8mSbNm1CcHBwqf1dunTBlClT8MILL+Dbb7/FjRs38Prrr2PTpk0W95Gfn6/bLms9LqCo2k4rLy/P4r7c1YTODbDlZJrRxeerIBsPiyWmZjZ9wkkjIyIiIiIiIqKyOHxNrspEEASbP1atWmW0bWMBl5aXlxeWL1+OJk2aAAB+/vlnpKaarkYyxcfHR7ddWFhY5rkFBQW6bV9fX4v6SU5OLvNDu/g9AGdfXBHNI4KwcGhryEWh1LGe4gnIBb014uS+QMOeThwdEREREREREZnCeVYVhFwux/PPP4+33noLAHDgwAE8++yzFrURGBio2y5vCmJOTo5uu7ypjSWVt96XISenXAAGtqmFxmGBWBF3GTsS0pGnVMPXS4bxwf8AWXonNuwJKPycPj4iIiIiIiIiKo0hlx3Z4yqDNWvWtPq+zZs3121bU8mlHz6Vtzi8/mLzlq2x5Rm0FV0LhrRCvkoNH6kQ4kfjDU/iVEUiIiIiIiIit+GSkEulUmH79u04ePAgkpKSkJWVBbXa9BXtgKKpgHv27HHSCK3TtGlTl/YvCKWn2FkiKioKMpkMarUa586dK/Nc/ePWXs3RLM4v5DIgigL8FHLg/O+AMvfBAUEEoh533cCIiIiIiIiIyIDTQ64DBw5g3LhxuHbtmm6fJJlOMgRBgCRJNgc4lYH+VSgjIiIsvr9CoUCHDh1w+PBhHD58GIWFhSYXoD9w4ACAogXo27VrZ92APcm5bYa36zwK+Ft31VAiIiIiIiIisj+nhlwnT57E448/jsLCQkiSBB8fHzRu3BjBwcEQRa6BbwuVSoVvvvlGd7tr165WtfPUU0/h8OHDyMzMxKZNmzB8+PBS56SkpGD37t0AgMcee8xgLS/7c3EpFwBo1MD5nYb7OFWRiIiIiIiIyK04NeSaPXs2CgoK4O3tjY8//hjjx483uKIfGbdv3z7ExMSYvMKiUqnECy+8oFsTrH///kbXybpy5Qrq168PAOjWrRv2799f6pwJEyZg7ty5uH//Pv773/+id+/eCA19ULGkVqsxZcoU3fTS//znPzZ+dR4g+W8g97bhvqb9XDMWIiIiIiIiIjLKqSFXXFwcBEHAO++8g8mTJzuza4/27bffYsCAARgwYAC6d++OJk2aICgoCNnZ2Th27Bi+/vpr3VTFsLAwfPbZZ1b3VbVqVXz44YeYNGkSrl69iocffhjvvPMOoqOjkZaWhk8//RT79u0DAIwYMQLdu3e3x5doWhlTWZ3m3HbD2zWigZB6LhkKERERERERERnn1JArPz8fAPD441yw21LZ2dn4/vvv8f3335s8Jzo6GuvXr9dVa1lr4sSJSEtLw//7f/8Ply5dwnPPPVfqnH79+hlMj/R0Go1UdBVFuQyiqLf+mySVXo+LUxWJiIiIiIiI3I5TQ6569erh7NmzUCqVzuzW402bNg1t2rTB4cOHkZiYiIyMDNy5cwfe3t6oUaMG2rVrhyFDhmDQoEGQyWR26XPOnDno06cPvvjiCxw8eBA3btxAcHAwWrdujfHjx2PEiBF26ad8jq3kSkzLxPK4JOxMuI48pRq+XjL0jQ7HhM4N0DwiCLiZCNy9YngnhlxEREREREREbsepIddTTz2Fs2fP4o8//sAjjzzizK49WrNmzdCsWTO8/vrrNrVTr169Mq9kWdKjjz6KRx991KY+bebA6YqbT6Zi6oZTUGke9JGnVGPT8VRsOZmGhUNbY+C9ElVcVeoA4dEOGxMRERERERERWceplzR87bXXULNmTXz00Ue4cuWKM7smMpCYllkq4NKn0kiYuuEU8hK2GB5o+gQgCEbvQ0RERERERESu49SQq3r16tixYwd8fX3x8MMPY9myZbh//74zh0AexzGVXMvjkkwGXFo1NDfhe/uM4U5OVSQiIiIiIiJyS06drggArVq1wh9//IGHH34YkyZNwuTJk1GtWjX4+fmVeT9BEHDp0iUnjZIqMo1Gws6E6+We11t2zHCHbwhQh9NsiYiIiIiIiNyR00Oun376Cc8//zyysrIgSRIkScLNmzfLvZ/AKWKVkwPW5MpXqZGnVJd7Xqx41HBHVF9A5vSXDBERERERERGZwanv2A8fPozhw4dDrS4KGOrWrYtWrVohODgYoujUmZNUifnIZfD1kpUZdAUhGx3Ec4Y7OVWRiIiIiIiIyG05NeR67733oFarUaVKFaxduxb9+vVzZvfkkexfySWKAvpGh2PT8VST57QV/4Vc0DzYIfcBGva0+1iIiIiIiIiIyD6cWj519OhRCIKAOXPmMOAil5rQuQHkoukpsDHiv4Y7arUFFGWvG0dEREREREREruPUkCs3NxcA0LlzZ2d2S57MAWtyAUDziCAsHNraZND1kFjiIge12jpkHERERERERERkH04NuerXrw/gQdhF5EoD29TClpc74+mHIuHrJQMA+HrJMCQmAo/4XDU8ObK9C0ZIREREREREROZy6ppcgwcPxpkzZ/Dbb7+xmovcgraia8GQVshXqeEjl0G8fQE4m2l4YmQ71wyQiIiIiIiIiMzi1EquqVOnonHjxvj0009x9OhRZ3ZNnspB0xVLEkUBfgo5RFEAUkp8bwbVAoIinDIOIiIiIiIiIrKOU0OuwMBA7NmzBy1btkTXrl3xzjvv4PTp08jPz3fmMIjKlloi5OJ6XERERERERERuz6nTFWUymW5bkiTMmzcP8+bNM+u+giBApVI5amjktpxTyWUg5YjhbU5VJCIiIiIiInJ7Tg25pBJTz0reJnK5whzgRqLhPi46T0REREREROT2nBpyzZo1y5ndUYXg5CA07SQgqR/cFmRAzTbOHQMRERERERERWYwhF7k5J4dcJdfjqtECUPg5dwxEREREREREZDGnLjxP5PZKXlmR63EREREREREReQSGXET6SoZctRhyEREREREREXkCh01XfPfdd+3e5syZM+3eJpHO/VQgK81wHxedJyIiIiIiIvIIDgu5Zs+eDUEQ7NomQ65KSpIAO38vGVVyPS7vKkBoI8f3S0REREREREQ2c+jC85Jkv0XD7R2YEZVSaqriQ4DIGb1EREREREREnsBhIde+ffsc1TRVNk6r5DpmeJtTFYmIiIiIiIg8hsNCrm7dujmqaSKraDQS8lVq+MhlEMUSoZlaBaSdMNzHKysSEREREREReQyHTlcksg/bpr0mpmVieVwSdiZcR55SDV8vGfpGh2NC5wZoHhFUdNLNRECZa3hHXlmRiIiIiIiIyGMw5CL3Z8PabptPpmLqhlNQaR60kadUY9PxVGw5mYaFQ1tjYJtaQMoRwzuG1Af8Q63ul4iIiIiIiIici6tqU4WVmJZZKuDSp9JImLrhFBLTMo2sx8UqLiIiIiIiIiJPwpCLPIB1lVzL45JMBlxaKo2EFXGXS19ZkYvOExEREREREXkUTlekCkmjkbAz4bpZ5x5MuAjIzhvu5HpcRERERERERB6FlVzk/qxYkytfpUaeUm3WuVHqfw13yBRAeEuL+yQiIiIiIiIi12HIRRWSj1wGXy+ZWee2l18y3FGzNSD3dsCoiIiIiIiIiMhRGHKRB7C8kksUBfSNDjfr3N5ByYY7OFWRiIiIiIiIyOMw5KIKa0LnBpCLQpnnyEUgSlViPS5eWZGIiIiIiIjI4zDkIvdnxZpcANA8IggLh7Y2GXTJRQFfPlEV8vw7hgcYchERERERERF5HF5dkTyAdSEXAAxsUwuNwwKxIu4ydiSkI0+phq+XDP2ia+L5zvXR/PYuwzv4VQOC69o4XiIiIiIiIiJyNlZyeYDZs2dDEASLPmbPnm1VX/v373d4H86mrej6Z04fJL7bB//M6YOFQ1ujeUQQkFFiqmJEDCCUPcWRiIiIiIiIiNwPK7kqqCZNmrh6CPZj5XTFkkRRgJ+ixLf83cuGt0Mb2aUvIiIiIiIiInIuhlweYMqUKRgyZEiZ56jVanTt2hWZmZkICgrCU089ZXO/33zzDdq3b2/yeFhYmM19uNydJMPbVeu7ZhxEREREREREZBOGXB4gLCys3EBp586dyMzMBAA888wz8PX1tbnf+vXro2XLlja3Yzv7VHIZdadEJVfVBo7ri4iIiIiIiIgchmtyVRCrV6/WbY8ZM8aFI/Eg+feBvBJXVgxhJRcRERERERGRJ2LIVQFkZmZi8+bNAIqqr7p06eLiEdmZndbkKqVkFZcgAsF1HNMXERERERERETkUQ64KYOPGjcjLywMAjB49GgKvDmiekutxBUUCcoVrxkJERERERERENmHIVQE4aqriO++8g7p168Lb2xshISGIiYnBG2+8gQsXLtitD/M4qJKr5JUVueg8ERERERERkcfiwvMe7sqVKzh48CAAoFOnTmjYsKHd2v7zzz9124WFhTh58iROnjyJRYsWYcaMGZg1a5ZVVWMpKSllHk9PTzfc4azpigy5iIiIiIiIiDwWQy4Pt2bNGkjFIdDYsWPt0mbNmjUxePBgdO7cGQ0aNIBcLse1a9ewbds2rF69GkqlEnPmzEFhYSHmzp1rcfu1a9e2yzhtdveK4W0uOk9ERERERETksQRJclSZDDlDkyZNcOHCBfj4+OD69euoUqWKTe3l5ORAoVDAy8vL6PH4+HjExsbi/v37EAQBJ06cQOvWrS3qw5Lqr+Q3AhA5NwXwse3rMurj5kBm6oPbQ1cDzQfavx8iIiIiIiIiN5SSkqIrRElOTkZkZKSLR2QbrsllR4Ig2PyxatUqs/v766+/dOtjDRw40OaACwD8/f1NBlwA0KFDByxevBgAIEmSbtsSycnJZX7Ex8dbPX6zKfOBzDTDfVUbOL5fIiIiIiIiInIITlf0YI5acL48w4cPx0svvYTMzEwcOHDA4vtbnAw7otjw3lWUWtA+pJ79+yEiIiIiIiIip2DIZUdnz561uY2aNWuadV5hYSF++OEHAECNGjXQp08fm/s2l1wuR1RUFI4ePYrU1NTy7+COSi46718d8A50zViIiIiIiIiIyGYMueyoadOmTutr27ZtuHPnDgBg5MiRkMlkTusbsGxdLds5oJLrTpLhbS46T0REREREROTRuCaXh3LVVEUAUKlUurXAIiIiHN+hI6Yr3i1RycX1uIiIiIiIiIg8GkMuD3T79m3s2LEDANC6dWuLr25oqx9++AH3798HAHTr1s2pfdtNyemKVVnJRUREREREROTJGHJ5oHXr1kGpVAKwrIrrypUruqs4du/evdTxu3fvYv/+/WW2ER8fj5dffhlA0ZTFyZMnm92/WylZycXpikREREREREQejWtyeSDtVEW5XI6RI0fard379++jR48eaNWqFZ566im0bdsWNWvWhEwmw7Vr17Bt2zasWbMGhYWFAID/+7//Q9u2be3Wv9No1MDdq4b7WMlFRERERERE5NEYcnmYc+fO4ciRIwCA2NhY1KhRw+59nD59GqdPnzZ5XCaTYcaMGZg5c6bd+zbK3mty3U8BNErDfVyTi4iIiIiIiMijMeTyMGvWrNFt23vB+YiICGzcuBGHDx9GfHw8UlNTcevWLeTn56NKlSpo0qQJunfvjgkTJqBevXp27dupSk5VVAQCfqGuGQsRERERERER2YUgSY64dB2R9VJSUlC7dm0AQPIbAYj8f5cB/2r26+DoSmDb6w9uh0cDk+Ls1z4RERERERGRBzB4/52cjMjISBePyDZceJ4qHI1GQm6hChqNifz2TpLhbS46T0REREREROTxOF2R3J+ZxYaJaZlYHpeEnQnXkadUw9dLhr7R4ZjQuQGaRwQ9OLHkdEWux0VERERERETk8RhykQcoP+TafDIVUzecgkqveitPqcam46nYcjINC4e2xsA2tYoO3LlieGdeWZGIiIiIiIjI43G6Inm8xLTMUgGXPpVGwtQNp5CYlllUFVaykovTFYmIiIiIiIg8HkMucn/lTFdcHpdkMuDSUmkkrIi7DORkAIXZhgdZyUVERERERETk8RhykUfTaCTsTLhu1rk7EtKhuV1i0XmZAgiq5YCREREREREREZEzMeQiD2C6SitfpUaeUm1WK3lKNZS3LhnuDK4LiDJbBkdEREREREREboAhF3k0H7kMvl7mhVS+XjIo7l813MmpikREREREREQVAkMucn9lrMkligL6Roeb1Uy/6JoQuOg8ERERERERUYXEkIs83oTODSAXhTLPkYsCnu9cv/SVFas2cODIiIiIiIiIiMhZGHKRByj7yonNI4KwcGhrk0GXXBSwcGhrNI8IAu6UDLlYyUVERERERERUEchdPQCicpUxXVFrYJtaaBwWiBVxl7EjIR15SjV8vWToF10Tz3euXxRw5WcCubcM78jpikREREREREQVAkMuqjC0FV0LhrRCvkoNH7kMon51V8mpihCAkLpOHSMREREREREROQZDLvIA5Vdy6RNFAX4KI9/aJacqVokE5N42jIuIiIiIiIiI3AXX5KLKo9SVFeu5ZBhEREREREREZH8Mucj9mbEml1m46DwRERERERFRhcWQiyqPO0mGt7noPBEREREREVGFwZCLPICdKrnuXjG8XbWBfdolIiIiIiIiIpdjyEXuzx7TFVUFwP0Uw32crkhERERERERUYTDkosrh7lWUqgjjdEUiIiIiIiKiCoMhF3kAO1Rylbyyol8o4BNke7tEDhbkqwAAJ/1JREFURERERERE5BYYclHlUOrKilyPi4iIiIiIiKgiYchF7s8ea3KVrOTiVEUiIiIiIiKiCoUhF1UOJRedD6nnkmEQERERERERkWMw5CIPYIdKrtw7hrf9q9veJhERERERERG5DYZcVDnk3ja87VfVNeMgIiIiIiIiIodgyEXuzx5rcuWVqORiyEVERERERERUoTDkoopPoyk9XdEv1DVjISIiIiIiIiKHYMhFFV/BfUBSG+5jyEVERERERERUoTDkIvdn63TFklVcAODL6YpEREREREREFQlDLqr4SoZccl9A4eeasRARERERERGRQzDkIg9gYyUXF50nIiIiIiIiqvAYclHFl3vb8DZDLiIiIiIiIqIKhyEXuT97r8nFReeJiIiIiIiIKhyGXA6WnZ2NP/74Ax999BGGDh2K+vXrQxAECIKAevXqWdzemTNnMHHiRDRs2BC+vr6oXr06unTpgq+++goqlcquY1+3bh1iY2MRHh4OHx8f1K1bF6NGjcLhw4ft2o/Dlazk4qLzRERERERERBWO3NUDqOj69++P/fv326WtZcuW4eWXX0ZhYaFuX35+PuLi4hAXF4eVK1di+/btqFatmk395OXlYciQIdixY4fB/mvXrmHt2rVYt24dZs6ciVmzZtnUj/nsvSYXK7mIiIiIiIiIKhpWcjmYpDfVrmrVqoiNjUVAQIDF7ezYsQOTJk1CYWEhatSogUWLFuHvv//Gzp07MXjwYABAfHw8Bg0aBLVabdOYn3vuOV3A1aNHD/zyyy+Ij4/HihUr0LBhQ2g0GsyePRtff/21Tf2YzebpilyTi4iIiIiIiKiiYyWXgz377LOYOHEi2rdvj0aNGgEA6tWrh+zsbLPbUCqVeOWVV6DRaBAUFIRDhw6hYcOGuuOPP/44XnrpJSxZsgRxcXFYs2YNxo0bZ9V49+7di/Xr1wMoqkL7+eefIZPJAADt27fHgAED0LZtW1y7dg3Tpk3DM888g5CQEKv6chquyUVERERERERU4bGSy8FefPFFjBgxQhdwWePnn39GUlISAGD69OkGAZfWggULdGHTggULrO7ro48+AgDI5XIsWbJEF3BpVatWDR9++CEA4N69e1i+fLnVfZmPC88TERERERERUdkYcnmAX375RbdtqkLLz88PQ4cOBQAkJibiwoULFveTlZWFPXv2AAB69eqFyMhIo+cNHjwYQUFBAIoCOLdXauF5N688IyIiIiIiIiKLMeTyAHFxcQCAJk2aIDw83OR53bp1020fOnTI4n6OHDmiW9Rev62SFAoFOnbsqLuPUqm0uC+L2LImlyRx4XkiIiIiIiKiSoAhl5vLzs5GcnIyAKBp06Zlnqt//OzZsxb3lZiYaLStsvpSqVT4999/Le7LaQoyAY3KcB8XniciIiIiIiKqcLjwvJtLSUnRbZuaPqhVu3Zt3bY2GHNWX82bN7eqH2P0x56epQHS0gFlkNntG7h7DcjUGO67kw9klz0GIiIiIiIiooouPT1dt61Sqco40zMw5HJzWVlZuu2AgIAyz/X399dtW3L1Rmf3pR+QlafD8lxgeS+L2i/XJ03s2x4RERERERGRh8vIyEC9evVcPQybcLqim8vPz9dtKxSKMs/19vbWbefl5bl1X0RERERERETkPm7cuOHqIdiMlVwABEGwuY2VK1eavPKhLXx8fHTb2kXhTSkoKNBt+/r6um1f5U2lvHz5Mrp27QoA+PPPPy2q/CLPk56ejg4dOgAA4uPjUbNmTRePiByJz3flwue7cuHzXbnw+a5c+HxXLny+K5fk5GQ8+uijAMpfm9sTMORyc4GBgbrt8qYF5uTk6LbLm27oyr7KW+9LX+3atS06nzxbzZo1+XxXIny+Kxc+35ULn+/Khc935cLnu3Lh81256Be+eCqGXLDuSoQlOSrdrlWrlm7bkgXbral+0v/hlZKSgnbt2jmsLyIiIiIiIiIie2LIBfcuyQsMDETt2rWRnJyMc+fOlXmu/vFmzZpZ3Jf+FRLN7Usul6Nx48YW90VEREREREREZE9ceN4DdO7cGQBw/vx5XL9+3eR5Bw4c0G136tTJ4n7at2+vW3Bev62SCgsL8ddff+nu4+XlZXFfRERERERERET2xJDLAzz11FO67VWrVhk9Jzc3Fxs2bABQVJEVFRVlcT+BgYF47LHHAAC7d+82OT1y06ZNyMzMBAAMGjTI4n6IiIiIiIiIiOyNIZcHGDRoEBo0aAAA+OCDD3Dp0qVS5/znP//B3bt3ddvGrFq1CoIgQBAEzJ492+g5//d//wcAUKlUeOmll6BWqw2O37p1C9OmTQMABAcHY8KECVZ9TURERERERERE9sQ1uRzs4sWLiIuLM9invXJhdnZ2qcqsxx9/HOHh4Qb7vLy88Pnnn6N///7IzMxEp06d8L///Q8dOnTA3bt3sWzZMvz0008AiqY2jh492urx9uzZE8OHD8f69euxZcsW9O7dG6+//joiIiKQkJCA999/H9euXQMAfPjhhwgJCbG6LyIiIiIiIiIie2HI5WBxcXEYP3680WO3b98udWzfvn2lQi4A6NevH7766iu8/PLLuHHjBl555ZVS53To0AE///wzZDKZTWP+5ptvkJmZiR07dmDfvn3Yt2+fwXFRFDFjxgy8+OKLNvVDRERERERERGQvDLk8yAsvvIBHHnkEixYtwp49e5CWlgZ/f380a9YMI0eOxIQJEyCX2/6U+vr6Yvv27fj++++xatUqnDp1Cvfu3UONGjXQpUsXvPzyy3jkkUfs8BUZFxkZCUmSHNY+uRc+35ULn+/Khc935cLnu3Lh81258PmuXPh8Vy4V7fkWpIr01RARERERERERUaXEheeJiIiIiIiIiMjjMeQiIiIiIiIiIiKPx5CLiIiIiIiIiIg8HkMuIiIiIiIiIiLyeAy5iIiIiIiIiIjI4zHkIiIiIiIiIiIij8eQi4iIiIiIiIiIPB5DLiIiIiIiIiIi8ngMuYiIiIiIiIiIyOMx5CKHuHr1KqZOnYqmTZvC398fVatWRfv27bFgwQLk5ubarZ+dO3di0KBBiIyMhLe3NyIjIzFo0CDs3LnTbn2QaUePHsW7776L2NhY3XMQEBCAqKgojB8/HnFxcXbpZ/bs2RAEwayP/fv326VPKs3c56B79+526W/dunWIjY1FeHg4fHx8ULduXYwaNQqHDx+2S/tkWvfu3c1+vm157fG17Rw3b97Etm3bMHPmTPTt2xfVqlXTPa7jxo2zuD1n/e7Nzc3F/Pnz0b59e1StWhX+/v5o2rQppk6diqtXr9q1r4rEHs93bm4uNm3ahMmTJ6N9+/YICQmBl5cXQkND8cgjj2D27Nm4fv26XcZbr149s34G1KtXzy79VTT2eL5XrVpl9s/iVatW2WXct27dwsyZM9GqVSsEBQUhKCgIrVq1wsyZM3H79m279FER2fp8X7lyxeLf77a89vj6to2932tVmt/fEpGdbdmyRQoKCpIAGP2IioqS/v33X5v6UKvV0vPPP2+yDwDShAkTJLVabaevikrq0qVLmY+/9mPMmDFSQUGBTX3NmjXLrL4ASPv27bPPF0ilmPscdOvWzaZ+cnNzpX79+plsXxRFafbs2fb5osiobt26mf18a5+TlJQUi/vha9s5ynpcx44da3Y7zvzd+++//0qNGzc22U9QUJC0detWm/upiGx9vk+dOiUFBASU+5oMCgqS1q9fb/N469ata9bPgLp169rcV0Vkj9f3ypUrzf5ZvHLlSpvH/Ndff0nh4eEm+6hZs6b0999/29xPRWTr83358mWLfr8DkGJjY60eL1/f1rPne63K9vtbDiI7OnHiBIYNG4a8vDwEBARg+vTp6NGjB/Ly8rB+/XosW7YMFy5cwBNPPIGjR48iMDDQqn7eeecdrFixAgAQExODt956Cw0bNsSlS5cwf/58nDhxAsuXL0f16tUxd+5ce36JVCwtLQ0AEBERgWeeeQZdunRBnTp1oFarcfjwYSxcuBCpqalYvXo1lEolvv/+e7v0m5CQUObx+vXr26UfMm3y5MmYMmWKyeP+/v42tf/cc89hx44dAIAePXrgtddeQ0REBBISEjB37lxcunQJs2fPRs2aNfHiiy/a1BcZt3LlSuTk5JR5TmJiIoYNGwYAeOyxx1CrVi2b+uRr2znq1KmDpk2bYteuXRbf11m/e7OysvDEE0/g33//BQC88MILGD58OHx9fbFv3z588MEHyMzMxLBhw3Do0CG0adPG6r4qOmue78zMTGRnZwMAOnXqhCeffBLt2rVDaGgoMjIysGnTJixbtgyZmZkYOXIkgoKC0LdvX5vHOnDgQLz33nsmjysUCpv7qOhseX1r/fbbb4iIiDB5PDIy0uq2ASA5ORn9+/dHRkYG5HI53nzzTTz55JMAgG3btuHjjz9Geno6+vfvj2PHjtncX0VmzfNdq1atcn/fAsAHH3yg+9t97NixVo9Ri69vy9nzvVal+/3t0AiNKh1t4iyXy6U///yz1PH58+frUtxZs2ZZ1cf58+cluVwuAZDatWsn5ebmGhzPycmR2rVrpxuHrVVjZNwTTzwh/fDDD5JKpTJ6PCMjQ4qKitI93wcOHLC6L/1qD3IdW1+75tizZ4+un/79+5f6/srIyJDq1KkjAZCCg4OlO3fuOGwsVLa33npL91ytWbPGqjb42naOmTNnSlu3bpWuX78uSZLh/+SbW+nhzN+9M2bM0I1v/vz5pY4fOnRINxZbK0crIluf70OHDklDhw6V/vnnH5Pn/PLLL5IgCBIAqWHDhpJGo7F6vNpKD0uqCukBe7y+9Su5Ll++7LjBSpI0evRoXV8bNmwodfyHH36wePyViT2e7/KoVCopIiJCAiAFBgaW+nlvCb6+rWev91qV8fc3/6oku/n7779139QTJ040eo5arZaaNWume4NaWFhocT+TJ0/W9XP48GGj5xw+fFh3zpQpUyzug+xj69atuufhlVdesbodvhF2D84Iufr27av7JZucnGz0nHXr1pX5C5QcT61WS7Vq1ZIASAEBAVJOTo5V7fC17RrWvCly1u/ewsJCqUqVKhIAqVmzZianTkycOFHXV3x8vFV9VRaOeBMsSZL09NNP69o9duyY1e3wTbB9uXPIlZ6eLomiKAGQ+vTpY/K8Pn36SEDRVPj09HSHjacicMTr+9dff9W1OX78eJva4uvbscx5r1UZf39z4Xmym19++UW3PX78eKPniKKIMWPGAADu3buHffv2WdSHJEnYvHkzAKBp06bo2LGj0fM6duyIJk2aAAA2b94MSZIs6ofso0ePHrrtS5cuuXAk5AmysrKwZ88eAECvXr1MTlEYPHgwgoKCAAA///yz08ZHD+zZswepqakAgCFDhsDPz8/FIyJHcubv3n379uH+/fsAiqbIiKLxP1X1F1jmzwHX4O94stSWLVug0WgAmH6vADx4fWs0GmzZssUZQyM9q1ev1m3bY6oiOU55P4cr6+9vhlxkN9qrO/j7+6Nt27Ymz+vWrZtu+9ChQxb1cfnyZd38ZP12yuonNTUVV65csagfso+CggLdtkwmc+FIyBMcOXIEhYWFAMp+fSsUCt0v6SNHjkCpVDplfPSA/h/A2v+4oIrLmb979a8UVVZf7dq104Wrlv4tQfbB3/FkKXNf37a8VyDbZGVl6QoX6tWrh65du7p2QFSm8n4OV9bf3wy5yG7Onj0LAGjUqBHkctPXNGjatGmp+5grMTHRaDv27ofs48CBA7rtZs2a2aXN2NhYhIWFQaFQICwsDN27d8e8efNw9+5du7RP5du4cSOaN28OPz8/BAYGonHjxhg7dqzFlZklWfP6VqlUusUtyTmys7N1//NWt25ddO/e3S7t8rXtvpz5u9fcvuRyORo1amR1P2Q7e/+O/+OPP9CmTRsEBgbCz88P9evXx7Bhw/DLL7+wIt+Jxo8fj4iICCgUClSrVg0dO3bE//73P131ri20r+8qVaogPDzc5Hk1a9bUVWzz9e1cP/74I3JzcwEAo0ePhiAIdmmXr2/HKO/ncGX9/c2Qi+wiPz8ft27dAlD+VVdCQkJ0V19LTk62qJ+UlBTddnn91K5dW7dtaT9kO41Gg3nz5uluDx061C7t/v7778jIyIBSqURGRgYOHDiA6dOno0GDBrpyXHKsxMREnD17Fnl5ecjOzsbFixexevVq9OzZE4MGDdKVKluKr2/P8NNPP+muvDhq1Ci7/QHM17b7cuZrU9uXv78/goODzeorIyPD4H+zyfFOnTqF7du3AwCio6PtEnJdvnwZp06dQnZ2NvLy8nDlyhVs2LABgwYNQpcuXewSslD59u/fj/T0dCiVSty+fRt///033n//fTRq1AhLly61qW3t69ucKyZqX9/8He9cjqrU5uvb/sx5r1VZf3+bLrchskBWVpZuOyAgoNzz/f39kZOTo7tEtSP60QZpACzuh2z3ySefID4+HkDRGkplTWE1R3R0NJ566il06NABERERUCqVOH/+PNauXYtdu3bh3r17ePrpp7F161a7XMqcSvPz88OAAQPw2GOPoWnTpggICNCFEV999RVu376NX375BQMHDsTvv/8OLy8vi9rn69sz2PsPYL623Z8zX5vavsz9W0K/L29vb4v7I8sVFBRgwoQJUKvVAID333/fpvYUCgUGDBiA2NhYtGzZElWqVMG9e/dw+PBhfPnll0hOTsahQ4fQu3dvHD58GFWqVLHHl0ElNGjQAIMHD8YjjzyiewOalJSEn376CT/++CPy8/MxadIkCIKAF1980ao+rHl983e881y7dk1XGfToo4/qqm1swde345jzXquy/v5myEV2kZ+fr9tWKBTlnq/9Rs7Ly3NYP/ovFkv7IdscOHAA//3vfwEAYWFh+PLLL21q7/XXX8fs2bNL7X/44YcxZswYLF26FJMmTYJarcaECRNw6dIl+Pj42NQnlZaammr0f2Z69+6NV155BX379sWJEydw4MABfPnll3j11Vctap+vb/eXkpKC/fv3AyhapDQqKsqm9vja9gzOfG1q+7Lkbwlr+yLrvPzyyzh69CiAosWF+/fvb1N78fHxRn+3dO/eHS+//DKGDBmCXbt24ezZs5gzZw4+/vhjm/qj0gYNGoSxY8eWqsxt3749hg0bhm3btmHw4MFQKpV44403MGDAgDKnG5pizeubr23n+e6773RTB+1VxcXXt2OY+16rsv7+5nRFsgv9Nx3ahaPLoi1L9PX1dVg/+qWPlvZD1vvnn38waNAgqFQq+Pj4YOPGjQgLC7OpzfJKXidOnIjnn38eAJCWloaffvrJpv7IuLKehxo1auDHH3/UVW99/vnnFrfP17f7++6773RXxrLHFZf42vYMznxtavuy5G8Ja/siy33wwQdYvnw5gKIA5IsvvrC5zbJ+DgQGBmLDhg2oWrUqAODrr78263uDLFOlSpUyp54/+eSTmDlzJgAgNzcXK1assKofa17ffG07z5o1awAUBRDDhg2zS5t8fdufJe+1Kuvvb4ZcZBeBgYG6bXPKG7XruZhTzmhtP9o+rOmHrHP58mXExsbi7t27kMlkWL9+vdOuyjJx4kTdtv4ijOQ8DRo0QO/evQEAFy9e1F3NxVx8fbs/R/wBXB6+tl3Pma9NbV+W/C1hbV9kmaVLl+Ltt98GULSo8I4dOwymnDhKlSpVMHz4cABFz7m2ioyc68UXX9QFYdb+LLbm9c3XtnPEx8fj3LlzAIABAwaU+59Q9sLXt2Usfa9VWX9/M+Qiu/Dx8UFoaCgAwwXu/n97dx4bVdXGcfw3tClp2QplRyhQqAkIyFqRkrKDZamIC4gKsoREQYgGiRoFEyDAi1EEI2AQiKhl0ZACAVnEAhYMLWVR2UohSGgJU6EsRWjLef9oei2dTpnpTKcO8/0kk1x67tznDGfOPXOfOXNuWa5du2a9sUsucOeKkgvmPSxOyQXz3I0D912+fFkDBgzQ5cuXZbPZ9PXXXyshIcFn8du1a2dts3hl1fGkHejf/22pqanWnXOGDRumunXr+iQufbvq+bJvFse6ffu2rl+/7lKsBg0asB5XJfv+++/1xhtvSCq6q+quXbtUv359n8XnPFD1GjZsaH3Wr2gbFPfvh51HpH/7N2O8b1TWgvOuoH+7piLXWoE6fpPkgtcUn6AyMjJUUFDgdL/ibwkk9285XfIkWPI43o4D99jtdg0cOFCZmZmSin6q5usB0lt3eINnPGmHivTv4OBgtW3btsIx4bqSH4C98VNFV9G3q54vx15XYxUUFOjcuXMVjgPXJSUl6bXXXtP9+/fVpEkT7dmzx6W743kT54H/Bk/bobh/5+bmKjs72+l+WVlZunHjhiT6ty/k5+crMTFRUlEyc8iQIT6NT/9+uIpeawXq+E2SC14TGxsrqSh7m5aW5nS/klOce/Xq5VaMVq1aqWnTpg7HKcu+ffskSc2aNVPLli3digPX5ebmavDgwdYMjwULFujNN9/0eT2K40uy3iPwPU/aoXv37tZileX173v37unQoUPWc9y9iyPcV/IDcIMGDXx6l0P6dtXz5dhb/FniYbFSU1OtWeHufpaA6/bs2aMXX3xRBQUFioiI0K5duxQVFeXzenAeqHpXr16V3W6XVPE2cLV/e3KtAPdt27ZNOTk5kqSXX35ZwcG+vTcd/bt8nlxrBer4TZILXvPss89a26tXry5zn/v371uzAcLDw9W3b1+3YthsNmta5qlTp6wL3dIOHTpkZZATEhL4hqCS5OXlaejQoTpy5Igk6YMPPtCsWbOqpC4rVqywtuPi4qqkDoHu/Pnz2rVrlyQpKipKzZo1c+v5tWrVUv/+/SVJu3fvdjqt+scff7S+4R05cqQHNYartm/frqtXr0ry/Qdg+nbV8+XY26dPH+sW8mvXrrXu9FXamjVrrG3OA5UjJSVFCQkJunv3rurUqaOffvpJ7du393k9cnNzrSR7WFiYunXr5vM6oGhR8OL+WNFz8YgRI1StWtHlp7NrBenf/l2tWjWNGDGiQrHguqqaqS3Rvx/G02utgB2/DeBFvXv3NpJMcHCwSUlJcShftGiRkWQkmdmzZzuU79271yofN25cmTFOnz5tgoKCjCTTrVs3k5eX90B5Xl6e6datm1WPM2fOeOOloZS7d++aQYMGWe01ffr0Ch1n9erV5b4njh8/bs6ePVvuMVasWGEdo3HjxubWrVsVqgucS0pKMvn5+U7Ls7OzTefOna12+OSTTxz2eVhbG2PMnj17rH1GjBhhCgoKHii/evWqadGihZFkwsPDzd9//+3R64JrRo0aZbVLWlqaS8+hb/93nT9//qFjbWneGnvHjRtnxd67d2+Z+3z44YfWPosWLXIoT0lJMcHBwUaSiYuLc6n+gawi7Z2enm7Cw8ONJFOjRg1z4MCBCsWOi4uzYp8/f96hfPv27Q7vpZJu3rz5wGeNadOmVagegcTd9j5//rw5cuRIufts2bLFhISEGEkmNDTUXLp0qcz9Htbexhjz6quvWvts3LjRoXzDhg1uv18DWUX6d0k5OTlW23bo0MGt59K/K5e3rrUCcfz27VxEPPKWLFmiXr166c6dOxo0aJDef/999e3bV3fu3FFiYqJWrlwpSYqOjtY777xToRjR0dGaOXOmFixYoNTUVPXq1UuzZs1SVFSUzp07p4ULFyo9PV2SNHPmTNbrqSRjxozRzp07JUn9+vXTxIkT9fvvvzvdPyQkRNHR0W7HSUtL06RJk9S3b18988wz6tChgyIiIlRQUKBTp07p22+/teoRFBSklStX+uRuT4Fm2rRpys/P16hRo9SzZ0+1bNlSoaGhstvt+uWXX7RixQrrZwyxsbEV/slqv379NHr0aCUmJiopKUkDBw7UjBkz1LRpU504cULz5s3TxYsXJUkLFy702eLngezatWvaunWrJOmJJ55Qly5dvHJc+rbvHDhwQBkZGda/i/uqVLSOZslvVSVp/PjxDsfw5dg7c+ZMrV+/XmfOnNG7776rjIwMjR49WqGhodq7d6/mz5+vgoIChYaG6rPPPqtwnEeVp+197tw5DR482Fo4eO7cuapTp065Y3zDhg2d3sK+PAsWLNDYsWP13HPPKTY2VlFRUapZs6Zyc3OVkpKi5cuXW+f8xx9/XHPmzHE7xqPO0/a+cOGC+vbtq549e2r48OHq1KmT1ZaZmZnatGmTNm3aZM3KWLx4sdsztUuaN2+eduzYoatXr2rMmDFKTU3VsGHDJElbt27VJ598Iqnop/Fz586tcJxHlTfO5yUlJibq3r17krw/i4v+7RlvXWsF5PhdaekzBKykpCRTu3ZtK4tb+hEdHe3023tXZnIZY0xhYaGZMGGC0xiSzMSJE01hYWElvUqU939f1iMyMrLM4zxstkfJ8vIeERERZvPmzZX7ogNYZGSkS+0watQoc+3atTKP4cpMLmOKvlGKj493GqNatWrlPh/e9eWXX5b7rZwz9O3/jpLfvrrycMYbY68r3wQbY8zZs2dN27ZtncapXbu22bJliyf/LY8sT9vb1b5Z8uHsnPywmR4ly8t7xMXFOZ09FOg8be+Sn73Le4SFhZkVK1aUWxdXZnIZY8yhQ4dM48aNncZq3LixOXTokKf/NY8kb53Pi8XExBhJJigoyGRlZblVF/p35XL3POzsWsuYwBu/mckFrxs+fLiOHz+uJUuWaNu2bbp06ZJCQkLUpk0bvfDCC5o6darCwsI8ilGtWjWtWrVKo0aN0sqVK3X48GHZ7XbVr19f3bt315QpU3y6MDIqT3x8vFatWqWDBw8qPT1dV65cUU5Ojowxqlevnjp16qQhQ4Zo/Pjxql27dlVX95G1du1aJScn6+DBg8rMzJTdbteNGzdUs2ZNNW/eXE8//bTGjRunnj17ehwrNDRU27Zt03fffac1a9bo2LFjun79uho1aqTevXtr6tSpXokD13zzzTeSimZTjR071mvHpW/7H1+OvW3atFF6erq++OILbdy4URkZGbp3756aN2+u+Ph4TZ8+XZGRkV6JhaqzePFi7dmzRwcPHtTp06dlt9t1/fp1hYWFqWnTpoqJidGYMWM0aNAg1letJF27dtW6det08OBBpaamKisrS3a7XQUFBapbt67at2+v/v37a9KkSRWarVeWmJgYnThxQkuWLNHmzZt14cIFSUWLZCckJGjGjBmKiIjwSiw4d/bsWf3222+SpIEDB6px48ZePT79+78j0MZvmzFOVgQDAAAAAAAA/AR3VwQAAAAAAIDfI8kFAAAAAAAAv0eSCwAAAAAAAH6PJBcAAAAAAAD8HkkuAAAAAAAA+D2SXAAAAAAAAPB7JLkAAAAAAADg90hyAQAAAAAAwO+R5AIAAAAAAIDfI8kFAAAAAAAAv0eSCwAAAAAAAH6PJBcAAAAAAAD8HkkuAAAAAAAA+D2SXAAAAAAAAPB7JLkAAAAAAADg90hyAQAAAAAAwO+R5AIAAAAAAIDfI8kFAAAAAAAAv0eSCwAAAJKkO3fuKDg4WDabTfPnz6/q6gAAALiFJBcAAAAkSenp6SosLJQkdevWrYprAwAA4B6SXAAAAJAkpaamWtskuQAAgL+xGWNMVVcCAAAAAAAA8AQzuQAAAAAAAOD3SHIBAAAAAADA75HkAgAAgE6dOiWbzSabzabExESH8vz8fLVv3142m00dO3a0Fqh/mA4dOshmsykyMtLbVQYAAHgASS4AAADo2LFj1vaTTz7pUL548WL9+eefkqTPP/9cQUFBLh23a9eukqSLFy/q0qVLnlcUAADACZJcAAAAsJJcoaGhatu27QNl2dnZmjt3riRp6NCh6tOnj8vHbd++vbWdlpbmeUUBAACcIMkFAAAAK8nVoUMHh1la8+fPV15eniTpo48+cuu4jz32mLV97tw5D2sJAADgHEkuAAAA6OjRo5KkTp06PfD3nJwcrVy5UpIUGxurHj16uHXc+vXrW9tZWVmeVRIAAKAcJLkAAAACnN1u1+XLlyU5rse1Zs0a3b17V5I0YcIEj+Lk5+d79HwAAIDykOQCAAAIcCUXnS89k2vDhg2SpKCgII0cOdLhuYcPH9bo0aP1yiuvyBjjUH779m1rOywszFtVBgAAcBBc1RUAAABA1SpOctlsNnXs2NH6+40bN6zF4jt37qzw8HCH5+7YsUPr169XdHS0bDabQ/nFixet7RYtWni55gAAAP9iJhcAAECAK16Pq3Xr1qpVq5b19z/++EOFhYWSipJcZfn1118lSdHR0eUeW5K6du3qhdoCAACUjSQXAABAgCueyVV6Pa6MjAxrOyoqyuF5//zzj/bv3y9JatasWZnHTk5OliTVrFnTIVGWmJgom82m5cuXKzk5WUOHDlXdunVVr149Pf/888rOzq7wawIAAIGHJBcAAEAAu3fvnk6ePCnJcT2u3Nxca7usnypu3bpVeXl5kqQaNWo4lB87dkyZmZmSpOHDhys4+MGVMopnef3888+Kj49XzZo1NXnyZLVu3Vo//PCDpkyZUuHXBQAAAg9rcgEAAASwkydPWnc9LD2Tq6Rbt245/G3ZsmWy2WwyxpRZvnTpUmv79ddfdygvTnIdO3ZMJ06cUOvWrSUVJd7atm2r3bt3yxhT5lpfAAAApTGTCwAAIICVXDOr9Eyu5s2bW9v79u17oGzDhg1KTk5WbGyspKK7LJa0f/9+rVmzRpL01FNPaeDAgU5jr1+/3kpwSVJISIhatWqlu3fvuv16AABA4CLJBQAAEMCK1+OqW7euw90P4+LiVL16dUnSli1b9PHHH+vw4cNaunSpxo8fr5CQEP3vf/+TzWZTenq6Vf7pp58qPj5ehYWFCg0N1apVqxziZmdn68qVK+rRo0eZM8gyMzPVqlUrZnEBAACXkeQCAAAIYMVJrtKzuKSidbjee+89SZIxRnPmzFGPHj301ltv6c6dO1q2bJliYmI0YMAASbLK3377bd26dUu1atXSli1b1K5dO4djF8/iKmuGV05Ojv76669yfz4JAABQGkkuAACAAObszorFZs+era+++kpdunRRWFiYatSooT59+mjnzp2aPHmyJGndunV66aWXFB4erurVq6tNmzaaPn26Tp06pf79+5d53OIkV5cuXRzKjhw5IkkOd2MEAAAoDwvPAwAABDC73f7QfSZNmqRJkyY5LW/YsKESExPdiluc5OratatDWXp6uqTyF8IHAAAojZlcAAAA8LmjR48qIiJCkZGRDmXM5AIAABVBkgsAAAA+lZeXp7NnzzpNYh05ckSNGjVSkyZNfFwzAADgz0hyAQAAwKeOHz+u+/fvl7ke182bN5WRkcFPFQEAgNtsxhhT1ZUAAAAAAAAAPMFMLgAAAAAAAPg9klwAAAAAAADweyS5AAAAAAAA4PdIcgEAAAAAAMDvkeQCAAAAAACA3yPJBQAAAAAAAL9HkgsAAAAAAAB+jyQXAAAAAAAA/B5JLgAAAAAAAPg9klwAAAAAAADweyS5AAAAAAAA4PdIcgEAAAAAAMDvkeQCAAAAAACA3yPJBQAAAAAAAL9HkgsAAAAAAAB+jyQXAAAAAAAA/B5JLgAAAAAAAPg9klwAAAAAAADweyS5AAAAAAAA4PdIcgEAAAAAAMDvkeQCAAAAAACA3/s/m8kBZ3PR9EgAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots(2,1, sharex=True, dpi=200)\n", "\n", "# prepare highfreq expansion for plotting on grid\n", "iwn = np.array([complex(x) for x in S.Sigma_iw.mesh])\n", "high_freq = S.Sigma_moments['up'][0][0,0] + S.Sigma_moments['up'][1][0,0]/iwn\n", "\n", "# plot Re component\n", "ax[0].oplot(S.Sigma_iw['up'][0,0].real, 'o', ms=3, label='CTHYB')\n", "ax[0].plot(iwn.imag, high_freq.real, '-', label=r'$\\Sigma_{\\infty} + \\Sigma_{1}/i\\omega_{n}$')\n", "\n", "# plot Im component\n", "ax[1].oplot(S.Sigma_iw['up'][0,0].imag, 'o', ms=3, label='CTHYB')\n", "ax[1].plot(iwn.imag, high_freq.imag, '-', label=r'$\\Sigma_{\\infty} + \\Sigma_{1}/i\\omega_{n}$')\n", "\n", "ax[0].legend(); ax[1].legend()\n", "ax[0].set_xlabel(r'$i\\omega_{n}$'); \n", "ax[1].set_xlabel(r'$i\\omega_{n}$'); \n", "ax[0].set_ylabel(r'Re$\\Sigma(i\\omega_{n})$ (eV)')\n", "ax[1].set_ylabel(r'Im$\\Sigma(i\\omega_{n})$ (eV)')\n", "ax[0].set_xlim(0, 20);\n", "ax[0].set_ylim(1,3); ax[1].set_ylim(-10, 5)\n", "plt.subplots_adjust(wspace=0.5)\n", "plt.show()" ] }, { "cell_type": "markdown", "id": "2d925b8f", "metadata": {}, "source": [ "Additionally, we can verify the relationship between the first moment of the Green's function and first moment of the self-energy,\n", "$$ \\mathbf{\\Sigma}_{\\infty} = \\mathbf{G}_{1} - (\\mathbf{G}_{0})_{1}$$\n", "$$ \\mathbf{\\Sigma}_{\\infty} = \\mathbf{G}_{1} - \\mathbf{E}, $$\n", "where $\\mathbf{E}$ are the local impurity levels." ] }, { "cell_type": "code", "execution_count": 6, "id": "de27ad45", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "-2.3092638912203256e-14\n" ] } ], "source": [ "E = {key[-1][1][0] : coef for key, coef in S.h_loc0}\n", "print((S.Sigma_moments['up'][0]- S.G_moments['up'][-1] + E['up'])[0,0].real)" ] }, { "cell_type": "markdown", "id": "b6338b04-0800-4409-929c-dc3abdff9291", "metadata": {}, "source": [ "Hence, the above equation holds exactly up to machine precision!" ] } ], "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.10.7" } }, "nbformat": 4, "nbformat_minor": 5 }