{ "cells": [ { "cell_type": "markdown", "id": "e6e99643-fa2f-40c5-b748-19ae6baae4f1", "metadata": {}, "source": [ "# Inverting Kepler's equation in ODEs\n", "\n", "In the [previous tutorial](<./Comparing coordinate systems.ipynb>) we considered the numerical solution of the Stark problem, whose Hamiltonian, in Cartesian coordinates, reads\n", "\n", "$$\n", "\\mathcal{H}_\\mathrm{cart}\\left(v_x, v_y, v_z, x, y, z \\right) = \\frac{1}{2}\\left( v_x^2+v_y^2+v_z^2 \\right) - \\frac{1}{\\sqrt{x^2+y^2+z^2}} - \\varepsilon z.\n", "$$\n", "\n", "We saw how a reformulation of the Stark problem using Delaunay elements, instead of Cartesian coordinates, can help reducing the total number of timesteps per unit of integration time. The Delaunay elements are the Hamiltonian version of the [Keplerian orbital elements](https://en.wikipedia.org/wiki/Orbital_elements), to which they are related via the following equations (valid in adimensional units):\n", "\n", "$$\n", "\\begin{aligned}\n", "L & = \\sqrt{a}, & l & = M, \\\\\n", "G & = \\sqrt{a\\left( 1 - e^2\\right)}, & g & = \\omega, \\\\\n", "H & = \\sqrt{a\\left( 1 - e^2\\right)} \\cos i, & h & = \\Omega.\n", "\\end{aligned}\n", "$$\n", "\n", "The Hamiltonian of the Stark problem in Delaunay elements reads:\n", "\n", "$$\n", "\\mathcal{H}_\\mathrm{Del} \\left( L, G, H, l, g, h \\right) = -\\frac{1}{2L^2}-\\varepsilon L\\sqrt{1-\\frac{H^2}{G^2}}\\left[ L\\left( \\cos E - \\sqrt{1-\\frac{G^2}{L^2}} \\right)\\sin g + G\\sin E \\cos g \\right].\n", "$$\n", "\n", "In this Hamiltonian, $E$ is the eccentric anomaly, which is a function of $L$, $G$ and $l$ implicitly defined by [Kepler's equation](https://en.wikipedia.org/wiki/Kepler%27s_equation):\n", "\n", "$$\n", "l = E - \\sqrt{1-\\frac{G^2}{L^2}} \\sin E.\n", "$$\n", "\n", "Because Kepler's equation is trascendental, it cannot be inverted to yield an explicit expression for $E$ in terms of elementary functions. In the [previous tutorial](<./Comparing coordinate systems.ipynb>) we worked around this issue by leaving $E$ as an unspecified function in the Hamiltonian; we then augmented Hamilton's equations via the derivatives of $E$ (which can be formulated explicitly) and numerically solved the augmented system of differential equations.\n", "\n", "Here we will employ a different approach, and we will use instead the function ``kepE()`` provided by heyoka.py's expression system. As the name suggests, ``kepE()`` represents symbolically the inversion of Kepler's elliptic equation. That is, ``kepE(e, M)`` is the bivariate function $E\\left(e,M\\right)$ implicitly defined by Kepler's elliptic equation:\n", "\n", "$$\n", "M = E - e \\sin E.\n", "$$\n", "\n", "For the evaluation of ``kepE(e, M)`` heyoka.py uses an iterative numerical scheme. The high-order derivatives necessary to implement Taylor's integration method are built using automatic-differentiation techniques starting from the evaluation of ``kepE(e, M)``.\n", "\n", "Let us begin by writing down the Hamiltonian:" ] }, { "cell_type": "code", "execution_count": 1, "id": "c82ace35-08ea-4c28-b04c-81b468c5193b", "metadata": {}, "outputs": [], "source": [ "import heyoka as hy\n", "\n", "# Numerical value for the eps constant.\n", "eps = 1e-3\n", "\n", "# Symbolic variables for the Delaunay arguments.\n", "L, G, H, l, g, h = hy.make_vars(\"L\", \"G\", \"H\", \"l\", \"g\", \"h\")\n", "\n", "# Define E as the solution to Kepler's equation.\n", "E = hy.kepE(hy.sqrt(1 - G**2 * L**-2), l)\n", "\n", "# The Hamiltonian.\n", "Ham_del = -0.5 * L**-2 - eps * L * hy.sqrt(1.0 - H**2 * G**-2) * (\n", " L * (hy.cos(E) - hy.sqrt(1.0 - G**2 * L**-2)) * hy.sin(g)\n", " + G * hy.sin(E) * hy.cos(g)\n", ")" ] }, { "cell_type": "markdown", "id": "a6a91824-ca5b-4140-8485-2434733ad532", "metadata": {}, "source": [ "Note how, in terms of Delaunay arguments, $e=\\sqrt{1-\\frac{G^2}{L^2}}$ and $M=l$.\n", "\n", "Unlike in the [previous tutorial](<./Comparing coordinate systems.ipynb>), and thanks to the fact that ``kepE()`` supports symbolic differentiation, we can now use the standard formulation of Hamilton's equations:" ] }, { "cell_type": "code", "execution_count": 2, "id": "a1b17f85-4efc-4ec0-8b20-fb1aa8a3bddb", "metadata": {}, "outputs": [], "source": [ "# Define the system of ODEs.\n", "sys = [\n", " (L, -hy.diff(Ham_del, l)),\n", " (G, -hy.diff(Ham_del, g)),\n", " (H, -hy.diff(Ham_del, h)),\n", " (l, hy.diff(Ham_del, L)),\n", " (g, hy.diff(Ham_del, G)),\n", " (h, hy.diff(Ham_del, H)),\n", "]" ] }, { "cell_type": "markdown", "id": "72fe3a74-e035-440c-b420-3646aabb4516", "metadata": {}, "source": [ "We can now proceed to the creation of the integrator object. The initial conditions are taken from the [previous tutorial](<./Comparing coordinate systems.ipynb>):" ] }, { "cell_type": "code", "execution_count": 3, "id": "56aedf2a-767a-4333-a202-778a99c8430f", "metadata": {}, "outputs": [], "source": [ "ta = hy.taylor_adaptive(\n", " sys,\n", " [\n", " 1.0045488165591647,\n", " 0.9731906288081488,\n", " -0.9683287292736491,\n", " 2.776991035843252,\n", " 4.314274521695855,\n", " 3.3415926535897924,\n", " ],\n", ")" ] }, { "cell_type": "markdown", "id": "5cfa5b10-aec7-41d6-bb2f-4a31ac51c6a3", "metadata": {}, "source": [ "Let us now integrate up to $t=250$ on a time grid. As usual, we will provide a callback to periodically reduce the $l$ angle modulo $2\\pi$:" ] }, { "cell_type": "code", "execution_count": 4, "id": "b0c52e16-f32c-49d4-9bb3-0cd07e7ae1af", "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "\n", "\n", "# Callback to reduce l to the\n", "# [-pi, pi] range.\n", "def mod_cb_del(ta):\n", " l = ta.state[3]\n", " if l < -np.pi or l > np.pi:\n", " ta.state[3] = (l + np.pi) % (2 * np.pi) - np.pi\n", "\n", " return True\n", "\n", "\n", "# Define a time grid for the integration.\n", "t_grid = np.linspace(0, 250, 1000)\n", "\n", "# Integrate.\n", "_, _, _, _, _, out_del = ta.propagate_grid(t_grid, callback=mod_cb_del)" ] }, { "cell_type": "markdown", "id": "26208487-5c97-4412-9f6a-9efad835a9c2", "metadata": {}, "source": [ "Let us take a look at the time evolution of the Delaunay elements:" ] }, { "cell_type": "code", "execution_count": 5, "id": "77b9c512-ad7d-43ab-b11b-9a15b00d6d28", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAABJcAAAJTCAYAAAC1oPdwAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguMiwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8g+/7EAAAACXBIWXMAAA9hAAAPYQGoP6dpAACgBUlEQVR4nOz9eZxcZZ3/f79r6+ot3Um6yWZCEhDZwiIJIgKyqBkYcGCYLy434w2ofO/MgIL4YyTiSHRGgwPDOKMDCs6g3qgwiygzbuCwCYiThAARlR0SlhCydXe609VdVef3R3VVn+Wq65yqVHdVdb+ej0c9JH36nLpOVaUuzzuf63NijuM4AgAAAAAAAKoQr/cAAAAAAAAA0LwIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFC1pgqXstmsPve5z2np0qVqa2vTAQccoC9+8YvK5/P1HhoAoEEwVwAAwjBXAEBtJes9gEp85Stf0Te+8Q195zvf0eGHH67169froosuUnd3ty677LJ6Dw8A0ACYKwAAYZgrAKC2mipc+vWvf62zzz5bZ555piRpyZIl+sEPfqD169fXeWQAgEbBXAEACMNcAQC11VTh0oknnqhvfOMbeuaZZ/S2t71NTzzxhB566CF99atfLbtPJpNRJpMp/Tmfz2vnzp3q6elRLBabhFEDwORzHEcDAwNasGCB4vGmWgG9z5grACAa5grmCgAIE3mucJpIPp93rrrqKicWiznJZNKJxWLOl7/8Zes+11xzjSOJBw8ePKblY8uWLZP0Dd04mCt48ODBo7IHcwVzBQ8ePHiEPcLmipjjOI6axO23364rr7xS1113nQ4//HA9/vjjuvzyy3XDDTfoggsuMO7j/xeGvr4+7b///tqyZYu6uroma+gAMKn6+/u1aNEi7d69W93d3fUezqRirgCAaJgrmCsAIEzUuaKpwqVFixbpqquu0iWXXFL62d/+7d/qtttu0x/+8IdIx+jv71d3d7f6+vqYBABMWdP5u465AgCimc7fdcwVABBN1O+6plpcPTQ0FFjjl0gkuGUoAKCEuQIAEIa5AgBqq6kaer///e/Xl770Je2///46/PDDtXHjRt1www366Ec/Wu+hAQAaBHMFACAMcwUA1FZTLYsbGBjQX//1X+vOO+/Utm3btGDBAn34wx/W5z//ebW0tEQ6BuWrAKaD6fxdx1wBANFM5+865goAiCbqd11ThUu1wCQAYDrgu27f8PoBmA74rts3vH4ApoMp2XMJAAAAAAAAjYVwCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFVrunDp1Vdf1Z//+Z+rp6dH7e3tOvroo7Vhw4Z6DwsA0ECYKwAAYZgrAKB2kvUeQCV27dqlE044Qaeeeqp+9rOfac6cOXr++ec1c+bMeg8NANAgmCsAAGGYKwCgtpoqXPrKV76iRYsW6dZbby39bMmSJfUbEACg4TBXAADCMFcAQG011bK4u+66SytWrNB5552nOXPm6O1vf7tuueUW6z6ZTEb9/f2eBwBg6mKuAACEYa4AgNpqqnDphRde0E033aSDDjpIv/jFL7Rq1Sp98pOf1He/+92y+6xdu1bd3d2lx6JFiyZxxACAycZcAQAIw1wBALUVcxzHqfcgomppadGKFSv0yCOPlH72yU9+UuvWrdOvf/1r4z6ZTEaZTKb05/7+fi1atEh9fX3q6uqa8DEDQD309/eru7t7Wn7XMVcAQDTMFcwVABAm6lzRVJVL8+fP12GHHeb52aGHHqrNmzeX3SedTqurq8vzAABMXcwVAIAwzBUAUFtNFS6dcMIJevrppz0/e+aZZ7R48eI6jQgA0GiYKwAAYZgrAKC2mipc+tSnPqVHH31UX/7yl/Xcc8/p+9//vm6++WZdcskl9R4aAKBBMFcAAMIwVwBAbTVVuHTsscfqzjvv1A9+8AMtW7ZMf/M3f6OvfvWrOv/88+s9NABAg2CuAACEYa4AgNpqqobetTCdGxcCmD74rts3vH4ApgO+6/YNrx+A6WBKNvQGAAAAAABAYyFcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFC1pg6X1q5dq1gspssvv7zeQwEANCjmCgBAGOYKANg3TRsurVu3TjfffLOOPPLIeg8FANCgmCsAAGGYKwBg3zVluLRnzx6df/75uuWWWzRr1qx6DwcA0ICYKwAAYZgrAKA2mjJcuuSSS3TmmWfqve99b+jvZjIZ9ff3ex4AgKmPuQIAEIa5AgBqI1nvAVTq9ttv12OPPaZ169ZF+v21a9fqC1/4wgSPCgDQSJgrAABhmCsAoHaaqnJpy5Ytuuyyy3TbbbeptbU10j6rV69WX19f6bFly5YJHiUAoJ6YKwAAYZgrAKC2Yo7jOPUeRFQ/+tGP9Kd/+qdKJBKln+VyOcViMcXjcWUyGc82k/7+fnV3d6uvr09dXV0TPWQAqIvp/F3HXAEA0Uzn7zrmCgCIJup3XVMti3vPe96jTZs2eX520UUX6ZBDDtFnPvOZ0AkAADD1MVcAAMIwVwBAbTVVuDRjxgwtW7bM87OOjg719PQEfg4AmJ6YKwAAYZgrAKC2mqrnEgAAAAAAABpLU1Uumdx///31HgIAoMExVwAAwjBXAED1qFwCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAULWmCpfWrl2rY489VjNmzNCcOXN0zjnn6Omnn673sAAADYS5AgAQhrkCAGqrqcKlBx54QJdccokeffRR3XPPPcpms1q5cqUGBwfrPTQAQINgrgAAhGGuAIDaijmO49R7ENV68803NWfOHD3wwAN697vfHWmf/v5+dXd3q6+vT11dXRM8QgCoD77rxjFXAIAZ33XjmCsAwCzqd11yEsdUc319fZKk2bNnl/2dTCajTCZT+nN/f/+EjwsA0DiYKwAAYZgrAGDfNNWyODfHcXTFFVfoxBNP1LJly8r+3tq1a9Xd3V16LFq0aBJHCQCoJ+YKAEAY5goA2HdNuyzukksu0U9+8hM99NBDWrhwYdnfM/0Lw6JFiyhfBTClUapfwFwBAOUxVxQwVwBAeVN6WdwnPvEJ3XXXXXrwwQetE4AkpdNppdPpSRoZAKBRMFcAAMIwVwBAbTRVuOQ4jj7xiU/ozjvv1P3336+lS5fWe0gAgAbDXAEACMNcAQC11VTh0iWXXKLvf//7+vGPf6wZM2Zo69atkqTu7m61tbXVeXQAgEbAXAEACMNcAQC11VQ9l2KxmPHnt956qy688MJIx2BtOYDpYDp/1zFXAEA00/m7jrkCAKKZkj2XmigHAwDUCXMFACAMcwUA1Fa83gMAAAAAAABA8yJcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNWS9R5As1j9wyc1PJqv9zCAuovVewAToNFvRnz8gT36wIpF9R4GInjmjQHl8tE/UZXeCdup4tM60XfbbsRzqPQZqrkleeXPMdHPMBmvU4U7qPLXNuy3ww4X+vnat83W5w977vCxhz23/TcWzGzT2+bOCDkKAAATg3Apov9+4nUNZLL1HgaAaai9JUG41CT+7MZHmCsA1MX5x+2vL/3pEfUeBgBgmiJciujTK9+m0Vyj1zcAqFasgUuyDpnXVe8hIKLeGWm1tiQq3q+aj181n9lYFc9U3fNUsc8k/SWs6nwm6bWu/rmqeZ4qPgtVPE81O4XtEjb28P3Dnj/k+PvwUZ3IsS+Y2Vb5gAAAqBHCpYguPGFpvYcAAGhw9/0/p9R7CAAAAMCko6E3AAAAAAAAqka4BAAAAAAAgKoRLgEAAAAAAKBqhEsAAAAAAACoGuESAAAAAAAAqka4BAAAAAAAgKoRLgEAAAAAAKBqhEsAAAAAAACoWlOGSzfeeKOWLl2q1tZWLV++XL/61a/qPSQAQINhrgCAxuU4Tr2HIIm5AgBqpenCpTvuuEOXX365rr76am3cuFEnnXSSzjjjDG3evLneQwMASPrpT3+qo446qvTn/+f/+X9088036ze/+Y2GhoYmZQzMFQDQ2Do7O3XCCSfok5/8pL7zne/ot7/9rfL5/KSOgbkCAGon5tTwnw2eeeYZXXTRRXr44YdrdciA4447Tsccc4xuuumm0s8OPfRQnXPOOVq7dm3o/v39/eru7lZfX5+6urombJwAUE/1/K77kz/5E73nPe/RZZddJkmaMWOGcrmchoeHFY/HdfDBB+s3v/mNOjs7J2wMzBUAEK6e33U33nijHnvsMW3YsEFPPfWUcrmcWltbdeSRR2r58uVavny5jjnmGM8/VtQacwUAhIv6XVfTyqXR0VE9+uijtTykx8jIiDZs2KCVK1d6fr5y5Uo98sgjxn0ymYz6+/s9DwDAxHnyySf1zne+0/OzTZs26YUXXtCdd96p1tZW3XrrrRP2/MwVAND4/vIv/1Lf+ta3tHHjRg0MDMhxHH3qU5/S29/+dq1fv16XXHKJjjnmmAl7fuYKAKitploWt337duVyOc2dO9fz87lz52rr1q3GfdauXavu7u7SY9GiRZMxVACYtrZu3aoFCxaU/pxMJhWLxbRkyRK9//3v15VXXqnbb799wp6fuQIAmks6nZYkfeADH9CNN96oRx99VAMDA9q4ceOEPSdzBQDUVkXh0qpVq3TLLbdo/fr1GhkZmagxhYrFYp4/O44T+FnR6tWr1dfXV3ps2bJlMoYIANNWb2+vXn755dKft27dqsWLF5f+fPTRR+t3v/vdhI+DuQIAmlcikdCRRx454c/DXAEAtZGs5JeffPJJfe9739Pg4KBSqZQOO+wwHXPMMaU10fH4xBZC9fb2KpFIBP41Ydu2bYF/dShKp9Olfw0BAEy80047Tf/6r/+qE088UZIC38HxeFyjo6MT9vzMFQCAMMwVAFBbFaVBjzzyiPr7+/XUU0/pX//1X3XaaafphRde0NVXX613vetdgR4btdbS0qLly5frnnvu8fz8nnvu0bve9a4JfW4AQDRXXnmlvve97+mrX/2qcfvDDz+sAw44YMKen7kCABrfxRdfrG984xtav369MpmMpGAV0URirgCA2qqockkqfOkfeuihOvTQQ3X++eeXfv78889rw4YNevzxx2s5voArrrhCH/nIR7RixQodf/zxuvnmm7V582atWrVqQp8XABDNEUccodtuu03nn3++vve97+kzn/mMjjvuOCUSCT300ENavXq1Pv3pT0/oGJgrAKCxPf300/q3f/s3DQwMKJksXJJ84Qtf0CmnnKJjjjlGRx99tNrb2yd0DMwVAFA7FYdL5Rx44IE68MAD9YEPfKBWhzT64Ac/qB07duiLX/yiXn/9dS1btkw//elPPf08AAD1dd555+mtb32rPvWpT+kDH/hA6V+jHcfR2WefrU996lMT+vzMFQDQ2B588EFJ0rPPPqsNGzboscce04YNG/T5z39eu3fvViKR0Nve9jY99dRTEzYG5goAqJ2Y4zhOvQcxmfr7+9Xd3a2+vj51dXXVezgAMCEa6btu8+bN2rRpkwYGBrRs2TItW7asruOJopFePwCYKI36Xffiiy9q/fr12rhxo7785S/XezhlNerrBwC1FPW7rmaVSwAAmOy///7af//96z0MAECTWLp0qZYuXarzzjuv3kMBAEQ0sbd3AwAAAAAAwJRGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBoNvQEAAABgEjiOo72juXoPA8A01pZKKBaL1fy4hEsAAAAAMAn2juZ02Od/Ue9hAJjGfvfFP1J7S+2jIJbFAQAAAAAAoGpULgEAAADAJGhLJfS7L/5RvYcBYBprSyUm5LiESwAAAAAwCWKx2IQsRwGAemNZHAAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqFrThEsvvfSSPvaxj2np0qVqa2vTgQceqGuuuUYjIyP1HhoAoEEwVwAAwjBXAEDtJes9gKj+8Ic/KJ/P65vf/Kbe+ta36re//a0uvvhiDQ4O6vrrr6/38AAADYC5AgAQhrkCAGov5jiOU+9BVOu6667TTTfdpBdeeCHyPv39/eru7lZfX5+6uromcHQAUD98141jrgAAM77rxjFXAIBZ1O+6pqlcMunr69Ps2bOtv5PJZJTJZEp/7u/vn+hhAQAaCHMFACAMcwUA7Jum6bnk9/zzz+trX/uaVq1aZf29tWvXqru7u/RYtGjRJI0QAFBvzBUAgDDMFQCw7+oeLq1Zs0axWMz6WL9+vWef1157TaeffrrOO+88ffzjH7cef/Xq1err6ys9tmzZMpGnAwCYAMwVAIAwzBUAUD9177m0fft2bd++3fo7S5YsUWtrq6TCBHDqqafquOOO07e//W3F45XlY6yNBjAdTLXvOuYKAKi9qfZdx1wBALXXND2Xent71dvbG+l3X331VZ166qlavny5br311oonAABAc2KuAACEYa4AgPqpe7gU1WuvvaZTTjlF+++/v66//nq9+eabpW3z5s2r48gAAI2CuQIAEIa5AgBqr2nCpbvvvlvPPfecnnvuOS1cuNCzrc4r+wAADYK5AgAQhrkCAGqvaeo/L7zwQjmOY3wAACAxVwAAwjFXAEDtNU24BAAAAAAAgMZDuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgak0ZLmUyGR199NGKxWJ6/PHH6z0cAEADYq4AAIRhrgCA2mjKcOmv/uqvtGDBgnoPAwDQwJgrAABhmCsAoDaaLlz62c9+prvvvlvXX399vYcCAGhQzBUAgDDMFQBQO8l6D6ASb7zxhi6++GL96Ec/Unt7e6R9MpmMMplM6c/9/f0TNTwAQANgrgAAhGGuAIDaaprKJcdxdOGFF2rVqlVasWJF5P3Wrl2r7u7u0mPRokUTOEoAQD0xVwAAwjBXAEDt1T1cWrNmjWKxmPWxfv16fe1rX1N/f79Wr15d0fFXr16tvr6+0mPLli0TdCYAgInCXAEACMNcAQD1E3Mcx6nnALZv367t27dbf2fJkiX60Ic+pP/6r/9SLBYr/TyXyymRSOj888/Xd77znUjP19/fr+7ubvX19amrq2ufxg4AjWqqfdcxVwBA7U217zrmCgCovajfdXUPl6LavHmzZ13za6+9pj/6oz/Sf/zHf+i4447TwoULIx2HSQDAdDBdv+uYKwAguun6XcdcAQDRRf2ua5qG3vvvv7/nz52dnZKkAw88MPIEAACY2pgrAABhmCsAoPbq3nMJAAAAAAAAzatpKpf8lixZoiZZ0QcAqBPmCgBAGOYKANh3VC4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqFrThUs/+clPdNxxx6mtrU29vb0699xz6z0kAECDYa4AAIRhrgCA2knWewCV+M///E9dfPHF+vKXv6zTTjtNjuNo06ZN9R4WAKCBMFcAAMIwVwBAbTVNuJTNZnXZZZfpuuuu08c+9rHSzw8++OA6jgoA0EiYKwAAYZgrAKD2mmZZ3GOPPaZXX31V8Xhcb3/72zV//nydccYZeuqpp6z7ZTIZ9ff3ex4AgKmJuQIAEIa5AgBqr2nCpRdeeEGStGbNGn3uc5/Tf//3f2vWrFk6+eSTtXPnzrL7rV27Vt3d3aXHokWLJmvIAIBJxlwBAAjDXAEAtVf3cGnNmjWKxWLWx/r165XP5yVJV199tf7sz/5My5cv16233qpYLKZ///d/L3v81atXq6+vr/TYsmXLZJ0aAKBGmCsAAGGYKwCgfurec+nSSy/Vhz70IevvLFmyRAMDA5Kkww47rPTzdDqtAw44QJs3by67bzqdVjqdrs1gAQB1wVwBAAjDXAEA9VP3cKm3t1e9vb2hv7d8+XKl02k9/fTTOvHEEyVJo6Ojeumll7R48eKJHiYAoI6YKwAAYZgrAKB+6h4uRdXV1aVVq1bpmmuu0aJFi7R48WJdd911kqTzzjuvzqMDADQC5goAQBjmCgCovaYJlyTpuuuuUzKZ1Ec+8hHt3btXxx13nO69917NmjWr3kMDADQI5goAQBjmCgCorZjjOE69BzGZ+vv71d3drb6+PnV1ddV7OAAwIfiu2ze8fgCmA77r9g2vH4DpIOp3Xd3vFgcAAAAAAIDmRbgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqjVVuPTMM8/o7LPPVm9vr7q6unTCCSfovvvuq/ewAAANhLkCABCGuQIAaqupwqUzzzxT2WxW9957rzZs2KCjjz5aZ511lrZu3VrvoQEAGgRzBQAgDHMFANRW04RL27dv13PPPaerrrpKRx55pA466CBde+21Ghoa0lNPPVXv4QEAGgBzBQAgDHMFANRest4DiKqnp0eHHnqovvvd7+qYY45ROp3WN7/5Tc2dO1fLly8vu18mk1Emkyn9ua+vT5LU398/4WMGgHopfsc5jlPnkUwu5goAiI65grkCAMJEniucJvLKK684y5cvd2KxmJNIJJwFCxY4GzdutO5zzTXXOJJ48ODBY1o+tmzZMjlf0A2EuYIHDx48KnswVzBX8ODBg0fYI2yuiDlOff+pYs2aNfrCF75g/Z1169Zp+fLlOuecczQ6Oqqrr75abW1t+ta3vqW77rpL69at0/z58437+v+FIZ/Pa+fOnerp6VEsFos8zv7+fi1atEhbtmxRV1dX5P2ayXQ4R2l6nCfnOHVUe56O42hgYEALFixQPN40K6DLYq5oHNPhHKXpcZ6c49TBXFHAXNE4psM5StPjPDnHqWOi54q6h0vbt2/X9u3brb+zZMkSPfzww1q5cqV27drleSEOOuggfexjH9NVV101oePs7+9Xd3e3+vr6puwHbjqcozQ9zpNznDqmy3mGYa5oHNPhHKXpcZ6c49QxXc4zDHNF45gO5yhNj/PkHKeOiT7Puvdc6u3tVW9vb+jvDQ0NSVIgKYvH48rn8xMyNgBAY2CuAACEYa4AgPppmvrX448/XrNmzdIFF1ygJ554Qs8884yuvPJKvfjiizrzzDPrPTwAQANgrgAAhGGuAIDaa5pwqbe3Vz//+c+1Z88enXbaaVqxYoUeeugh/fjHP9ZRRx014c+fTqd1zTXXKJ1OT/hz1ct0OEdpepwn5zh1TJfzrBXmiok3Hc5Rmh7nyTlOHdPlPGuFuWLiTYdzlKbHeXKOU8dEn2fdey4BAAAAAACgeTVN5RIAAAAAAAAaD+ESAAAAAAAAqka4BAAAAAAAgKoRLgEAAAAAAKBqhEsR3XjjjVq6dKlaW1u1fPly/epXv6r3kKq2Zs0axWIxz2PevHml7Y7jaM2aNVqwYIHa2tp0yimn6KmnnqrjiMM9+OCDev/7368FCxYoFovpRz/6kWd7lHPKZDL6xCc+od7eXnV0dOhP/uRP9Morr0ziWdiFneOFF14YeF/f+c53en6n0c9x7dq1OvbYYzVjxgzNmTNH55xzjp5++mnP7zT7exnlHKfCezldMVcwV9Qbc0VBs7+XzBVTG3MFc0W9MVcUNPt72WhzBeFSBHfccYcuv/xyXX311dq4caNOOukknXHGGdq8eXO9h1a1ww8/XK+//nrpsWnTptK2v/u7v9MNN9ygr3/961q3bp3mzZun973vfRoYGKjjiO0GBwd11FFH6etf/7pxe5Rzuvzyy3XnnXfq9ttv10MPPaQ9e/borLPOUi6Xm6zTsAo7R0k6/fTTPe/rT3/6U8/2Rj/HBx54QJdccokeffRR3XPPPcpms1q5cqUGBwdLv9Ps72WUc5Sa/72cjpgrmCsaAXNFQbO/l8wVUxdzBXNFI2CuKGj297Lh5goHod7xjnc4q1at8vzskEMOca666qo6jWjfXHPNNc5RRx1l3JbP55158+Y51157belnw8PDTnd3t/ONb3xjkka4byQ5d955Z+nPUc5p9+7dTiqVcm6//fbS77z66qtOPB53fv7zn0/a2KPyn6PjOM4FF1zgnH322WX3abZzdBzH2bZtmyPJeeCBBxzHmZrvpf8cHWdqvpfTAXMFc0WjYa6YOu8lc8XUwVzBXNFomCumzntZ77mCyqUQIyMj2rBhg1auXOn5+cqVK/XII4/UaVT77tlnn9WCBQu0dOlSfehDH9ILL7wgSXrxxRe1detWz/mm02mdfPLJTXu+Uc5pw4YNGh0d9fzOggULtGzZsqY67/vvv19z5szR2972Nl188cXatm1baVsznmNfX58kafbs2ZKm5nvpP8eiqfZeTnXMFcwVzXTeU+37hbli6ryXUx1zBXNFM533VPt+Ya6Y+PeScCnE9u3blcvlNHfuXM/P586dq61bt9ZpVPvmuOOO03e/+1394he/0C233KKtW7fqXe96l3bs2FE6p6l0vlHOaevWrWppadGsWbPK/k6jO+OMM/S9731P9957r/7+7/9e69at02mnnaZMJiOp+c7RcRxdccUVOvHEE7Vs2TJJU++9NJ2jNPXey+mAuUKlPzfr+U6175dyptr3C3PF1HkvpwPmCpX+3KznO9W+X8qZat8vzBWT814m9/00podYLOb5s+M4gZ81izPOOKP030cccYSOP/54HXjggfrOd75Tau41lc63qJpzaqbz/uAHP1j672XLlmnFihVavHixfvKTn+jcc88tu1+jnuOll16qJ598Ug899FBg21R5L8ud41R7L6eTqfTdyVxR0KzfL+VMte8X5oqCqfBeTidT6buTuaKgWb9fyplq3y/MFQUT/V5SuRSit7dXiUQikNpt27YtkHI2q46ODh1xxBF69tlnS3d3mErnG+Wc5s2bp5GREe3atavs7zSb+fPna/HixXr22WclNdc5fuITn9Bdd92l++67TwsXLiz9fCq9l+XO0aSZ38vpgrmioJnPdyp9v1Simb9fmCu8mvm9nC6YKwqa+Xyn0vdLJZr5+4W5wmsi30vCpRAtLS1avny57rnnHs/P77nnHr3rXe+q06hqK5PJ6Pe//73mz5+vpUuXat68eZ7zHRkZ0QMPPNC05xvlnJYvX65UKuX5nddff12//e1vm/a8d+zYoS1btmj+/PmSmuMcHcfRpZdeqh/+8Ie69957tXTpUs/2qfBehp2jSTO+l9MNcwVzRbOedzN+vzBXmDXjezndMFcwVzTreTfj9wtzhdmEvpcVtf+epm6//XYnlUo5//Iv/+L87ne/cy6//HKno6PDeemll+o9tKp8+tOfdu6//37nhRdecB599FHnrLPOcmbMmFE6n2uvvdbp7u52fvjDHzqbNm1yPvzhDzvz5893+vv76zzy8gYGBpyNGzc6GzdudCQ5N9xwg7Nx40bn5Zdfdhwn2jmtWrXKWbhwofPLX/7Seeyxx5zTTjvNOeqoo5xsNluv0/KwnePAwIDz6U9/2nnkkUecF1980bnvvvuc448/3nnLW97SVOf4F3/xF053d7dz//33O6+//nrpMTQ0VPqdZn8vw85xqryX0xFzBXNFI2CuKGj295K5YupirmCuaATMFQXN/l422lxBuBTRP//zPzuLFy92WlpanGOOOcZze79m88EPftCZP3++k0qlnAULFjjnnnuu89RTT5W25/N555prrnHmzZvnpNNp593vfrezadOmOo443H333edICjwuuOACx3GindPevXudSy+91Jk9e7bT1tbmnHXWWc7mzZvrcDZmtnMcGhpyVq5c6ey3335OKpVy9t9/f+eCCy4IjL/Rz9F0fpKcW2+9tfQ7zf5ehp3jVHkvpyvmCuaKemOuKGj295K5YmpjrmCuqDfmioJmfy8bba6IjQ0KAAAAAAAAqBg9lwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNcAgAAAAAAQNUIl4AaWLNmjY4++uh6DwMA0MCYKwAAYZgr0KxijuM49R4E0MhisZh1+wUXXKCvf/3rymQy6unpmaRRAQAaCXMFACAMcwWmMsIlIMTWrVtL/33HHXfo85//vJ5++unSz9ra2tTd3V2PoQEAGgRzBQAgDHMFpjKWxQEh5s2bV3p0d3crFosFfuYvX73wwgt1zjnn6Mtf/rLmzp2rmTNn6gtf+IKy2ayuvPJKzZ49WwsXLtS//uu/ep7r1Vdf1Qc/+EHNmjVLPT09Ovvss/XSSy9N7gkDACrGXAEACMNcgamMcAmYIPfee69ee+01Pfjgg7rhhhu0Zs0anXXWWZo1a5Z+85vfaNWqVVq1apW2bNkiSRoaGtKpp56qzs5OPfjgg3rooYfU2dmp008/XSMjI3U+GwDARGCuAACEYa5AMyBcAibI7Nmz9U//9E86+OCD9dGPflQHH3ywhoaG9NnPflYHHXSQVq9erZaWFj388MOSpNtvv13xeFzf+ta3dMQRR+jQQw/Vrbfeqs2bN+v++++v78kAACYEcwUAIAxzBZpBst4DAKaqww8/XPH4eH47d+5cLVu2rPTnRCKhnp4ebdu2TZK0YcMGPffcc5oxY4bnOMPDw3r++ecnZ9AAgEnFXAEACMNcgWZAuARMkFQq5flzLBYz/iyfz0uS8vm8li9fru9973uBY+23334TN1AAQN0wVwAAwjBXoBkQLgEN4phjjtEdd9yhOXPmqKurq97DAQA0IOYKAEAY5grUAz2XgAZx/vnnq7e3V2effbZ+9atf6cUXX9QDDzygyy67TK+88kq9hwcAaADMFQCAMMwVqAfCJaBBtLe368EHH9T++++vc889V4ceeqg++tGPau/evfyLAwBAEnMFACAccwXqIeY4jlPvQQAAAAAAAKA5UbkEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACq1lThUjab1ec+9zktXbpUbW1tOuCAA/TFL35R+Xy+3kMDADQI5goAQBjmCgCorWS9B1CJr3zlK/rGN76h73znOzr88MO1fv16XXTRReru7tZll11W7+EBABoAcwUAIAxzBQDUVlOFS7/+9a919tln68wzz5QkLVmyRD/4wQ+0fv36Oo8MANAomCsAAGGYKwCgtppqWdyJJ56o//mf/9EzzzwjSXriiSf00EMP6Y//+I/rPDIAQKNgrgAAhGGuAIDaaqrKpc985jPq6+vTIYccokQioVwupy996Uv68Ic/XHafTCajTCZT+nM+n9fOnTvV09OjWCw2GcMGgEnnOI4GBga0YMECxeNN9e8I+4y5AgCiYa5grgCAMJHnCqeJ/OAHP3AWLlzo/OAHP3CefPJJ57vf/a4ze/Zs59vf/nbZfa655hpHEg8ePHhMy8eWLVsm8Vu6MTBX8ODBg0dlD+YK5goePHjwCHuEzRUxx3EcNYlFixbpqquu0iWXXFL62d/+7d/qtttu0x/+8AfjPv5/Yejr69P++++vLVu2qKura8LHDAD10N/fr0WLFmn37t3q7u6u93AmFXMFAETDXMFcAQBhos4VTbUsbmhoKFCGlUgkrLcMTafTSqfTgZ93dXUxCQCY8qZjmT5zBQBUhrmigLkCAMoLmyuaKlx6//vfry996Uvaf//9dfjhh2vjxo264YYb9NGPfrTeQwMANAjmCgBAGOYKAKitploWNzAwoL/+67/WnXfeqW3btmnBggX68Ic/rM9//vNqaWmJdIz+/n51d3err6+Pf2EAMGVN5+865goAiGY6f9cxVwBANFG/65oqXKoFJgEA0wHfdfuG1w/AdMB33b7h9QMwHUT9rpte9xwFAAAAAABATREuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKha04VLr776qv78z/9cPT09am9v19FHH60NGzbUe1gAgAbCXAEACMNcAQC1k6z3ACqxa9cunXDCCTr11FP1s5/9THPmzNHzzz+vmTNn1ntoAIAGwVwBAAjDXAEAtdVU4dJXvvIVLVq0SLfeemvpZ0uWLKnfgAAADYe5AgAQhrkCAGqrqZbF3XXXXVqxYoXOO+88zZkzR29/+9t1yy23WPfJZDLq7+/3PAAAUxdzBQAgDHMFANRWU4VLL7zwgm666SYddNBB+sUvfqFVq1bpk5/8pL773e+W3Wft2rXq7u4uPRYtWjSJIwYATDbmCgBAGOYKAKitmOM4Tr0HEVVLS4tWrFihRx55pPSzT37yk1q3bp1+/etfG/fJZDLKZDKlP/f392vRokXq6+tTV1fXhI8ZAOqhv79f3d3d0/K7jrkCAKJhrmCuAIAwUeeKpqpcmj9/vg477DDPzw499FBt3ry57D7pdFpdXV2eBwBg6mKuAACEYa4AgNpqqnDphBNO0NNPP+352TPPPKPFixfXaUQAgEbDXAEACMNcAQC11VTh0qc+9Sk9+uij+vKXv6znnntO3//+93XzzTfrkksuqffQAAANgrkCABCGuQIAaqupwqVjjz1Wd955p37wgx9o2bJl+pu/+Rt99atf1fnnn1/voQEAGgRzBQAgDHMFANRWUzX0roXp3LgQwPTBd92+4fUDMB3wXbdveP0ATAdTsqE3AAAAAAAAGgvhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgaoRLAAAAAAAAqBrhEgAAAAAAAKpGuAQAAAAAAICqES4BAAAAAACgasl6DwAAgKlizV1PKRGPaV5Xq+Z2t2q/zrQcORrJ5jWSzSsWi2m/GWnNmZHWfjPSSsZjymTz2juSkyR1taWUiMfqfBYAAABAZQiXAACoAcdx9G/rt2hoLCiKIhGPKZd3PD/rak1qZnuLZranFI/FNJjJKhGPaW5Xq+Z3t2p+d5vaWuIaGM5qYDirvOOot7MQWM3tblVvR1qDI1ntGc5qNJdXRzqpBTPbNL+7VR3pwrSfyea0ZzirdCqhjpaEYjECLQAAAFSvqcOltWvX6rOf/awuu+wyffWrX633cAAADWiy5oq8I3165cF6o39YW/uGtbV/WNv3ZJSIxdSSjKslGVc+7+jNgYze3JPRaM4JBEuS1D+cVf9wVpt3en/+h60D+zzG9paEsjlHI7l86WetqbjmdrVq7oxWtbYk1Ld3VMMjOeUcRz0dLVows029nS2a2d6iRDymvOPIcQphWioR15yutObMaNVoLq/+4axGsnkl4iocs6tVPR0t6kwn5UjK5Qv75h1HbamE4lRpAWgQXFcAwL5p2nBp3bp1uvnmm3XkkUfWeygAgAY1mXNFIh7Tx05cGul383lHu4ZGNJpz1J5OqD2VkCOpb++odg+Nqm/viHYNjirnOJqRTmokl9cb/cN6vW9Yr+8e1kgur850UjNak4rHYnpzIKNtA4XtOwdHCtvaUmpJxNS3d1Sv7tqrwZGcsapqeDSvl3cM6eUdQ4Ftz+3ri2KRjBeWCHa3pTQwnNXgSFaj2bw6W5Oa19Wq3s60ZnW0KJWIKZ+Xco6jvOMopph6O1tKywr79ma1dzQnx3HU09miuV2tmtXeou62VCnQyjuFIC+djGu/GWn1dKQ1PJpT395RZfOOYjGpM51URzqpdkIvYNrhugIA9l1Thkt79uzR+eefr1tuuUV/+7d/W+/hAAAaUCPPFfF4TD2d6cDPezvT6jX8fF85jqPBkZy2D2TUkoyrszWpzpakhrM5vTmQ0Rv9Gb3RP6zh0Zy621LqSCcVk/Tmnoxe2z2snYMZ7RoaleNI8ZgUj8UUjxeCqW0Dw9rWn1EqEVdXW1KtqYRGsnlt7R/Wm/0ZDWSyxjFl804hLOsb9vx8cCSnN/ozNX8NoorFpI6WpFKJmIZGcnIcaXZHi3pntGh2R1ptqbj69o5qaCSnXL5QvVUIphLqSCfVmU4qm3c0lMmqNZVQT2eLejrSmtWRUi5fCBBHsnk5ctSaTKgznVR7OqGOlqTSqbgyo3llsjl1taXU25nWzPaUulpT2pMpLHXM5gu9uzpakqVgMpko3J+lWAlH3y4gukaeKwCgmTRluHTJJZfozDPP1Hvf+97QSSCTySiTGf8/qf39/RM9PABAA2CuGBeLxdQ5Fny4tbcktbgnqcU9HRP23JlsTkOZXCmQisdiisUKIcu2/oz69o6qqy2lznRCqUQhuClWYO0cHFEu7ygRL+yTiMWUcxxtHxjRtoFhOU6hCXp7S0IxSdv3FIKyXUMjGhguhFqJeEyJeEzxmLR3JKfte0ZKywI7WgrBTD7vaHAkq7wjOY60xxeIbe0vLHNsVC3JuOKxQtgnSTPbU5o91rcrmYir3xWGtSTjam9JqL0lobaWQqXWaC6vwZGs2luSmtXeolntKXW3pTSSy6t/76gyY83oW1PxUoVXRzpZek2HRwthWDFEi8cKvcT2juQ0NFZV1pZKqL0lqbaWuFpTCUmF8WZzebUk40onE2P/W3jkHWk0l1cqEVdbKqHWlrhaEvFSf7DRXF4xqRSsAdVirgCA2mi6cOn222/XY489pnXr1kX6/bVr1+oLX/jCBI8KANBImCsaRzqZUDqZCPy8vSWp+d1txn2OXDhx43EcR3vGqopSrmDCcRwNj+YLFUKZ8WbokrRjT0Y79oxox+CI9o5Vd3WmE4rHYhrNORoc22dw7JGIFwKcvaM5bd+T0c7BEe0aGlEiHld3W0qtycLzDmfzpX2GRnLaO5pTWyqh1Nhyxp2DI9q9t1AxJhWW7iUTMeXzjoZGcsqOVSqNZPOec9w9VFheOdUk4jGlk3GNZPOlc4/H5AmnUvGYRvOOHMdROplQa6oQZhV7fGVGc8pk80omYmpJxMf6oSXUMvZZGM3lFY9JranE2KNw7LxTCMtyjlPaL5UoPFoSMWXzTqmSrbhPa6qwPZPNa3i0ELCmEjG1JBNj/1sIy7K5vPKOSsFaOpXwfEZGXOMtPGdMqURcecfRaC6v1lRCs9pbNLujpRTcIRrmCgConaYKl7Zs2aLLLrtMd999t1pbWyPts3r1al1xxRWlP/f392vRokUTNUQAQJ0xV8AmFotpRmvK+PO2loTaWhLab4Z3aeJbZppDsMmQyzsaGqsq8i93G8nmNTSS1eBITvm8o450UnnH0a6xqq9dQ6PK5vPqbkuV9h/NFQKtvSM5DY4FWi2JmNpakto7ktWOwRH1jYVT6VQhDEsn48rlpb2jOQ2NjAdpeafQJD6dLNy9cMeeEe0cGpHjOErG42ptSagtFVdMMe0dLVQ4FY5R6P3VlkooGY9pJJdXZixEyWRzY03hY2PjHW98nxsL1dzyTqECqli1NZ392TEL9fcfOKrew2gazBUAUFtNFS5t2LBB27Zt0/Lly0s/y+VyevDBB/X1r39dmUxGiYT3X2zS6bTS6dr3rwCaTf/wqLoMF1TAVMNcgakkETeHYZLGqm5aNLPd+/OJ6NtVT6O5fCGcGslpeLSwjK41VajsGQ+lCsFUsXooFpMyo/lSqDU8Wqj0aksVKpyyeUcjY/uOZPMacS2zyzuOMqO5sdAqp+FsTol44TkTsUJAN5Ib3380V6gsSo6FYZns+L6juXyhei8Vlxwpk8trdGyfYnCWSsYVk0rh2vBY3y1HGqtkiyuXd8aed2zfbGG5aHKsN9juoRHNameOrwRzBQDUVlOFS+95z3u0adMmz88uuugiHXLIIfrMZz4TmACAMNlcoY/EVG9++u2HX9QX/vt3+uafL9fKw+fVezj7JJ939OKOQR3Q21HqvTHRBjNZXX3nJp155AK977C5k/KcqB5zBTC1FJef8Q8k5TmOU1oqiGiYKwCgtpoqXJoxY4aWLVvm+VlHR4d6enoCPwfCjGTz+qOvPqjezhb9+6p31Xs4E2r9y7vkONLTWweaPly69ZGX9Df//TutPfcIffgd+0/Kc97/9Jv60eOv6bW+YWO4lM3lxxoOT+2QslkwVwCYbmJj/ZwQHXMFANQWt9iYJIW72lT+L0rPv7lHP930esX75vKOPvejTfrRxlcrfs7p4tXde/Xi9kFteHlXTY+byea05q6n9PBz22t63H2xa2hEkjTZ/6aZm4B/Rf3da4U7s2zeOVTzY5ezY7BwZ5hsLtjTY+9ITidfd78+/p31kzYeAAAAAGgkTVW5ZHL//ffXewih7n5qq/7v/3+D/ur0g/WXp7y1on2v+Lcn9MSW3frpJ0/SYQu6Iu/3+Jbduu3RzXro2e065+1vCWzP5x05UsXLwXJ5Rxte3qVlb+lSe0tlH58NL+/Uw8/t0F+ccqDnDj31snNwYgKXXz2zXd9+5CX97rV+nfDW3podN5939PQbA3rb3BkVv287Bwt3Daoi36zabY++rLU//b2+/dF36Ngls2t23GJQNplsn5Xn39yjV3fvVf+w+c5Mz74xoD2ZrN6+/6wJHCHCNMNcAQCoL+YKAKhe/a/wp4EnX+mTJD2/bbDifV8Zq87YvbeyC+rtewqVFqO54OVwPu/o7H9+WH9648MVV0T94qmt+sA3f62/+/nTFe0nSX/z37/XDfc8o3Uv7ax4X5s9maz6qrjl8q5iYFDjwGXbQOG1zxiqXPbF7eu26Ix//JW+9asXKt5351jljTOJtUsPPPOmBkdyemLL7sC2kWxe9z29TYOZbMXH3TlB71u1z1ncZnppHcfRh2/5jT74zUc1UCZ8AgAAAIBmR7g0CXaWliRVdjWczzvjVRoVXkgXgxOTHYMj2vRqn558pa/iW/c+t22PpMIyv0q90V/YZ3g0F9j289++rvf8/f367at9FR3TcRyd+U+/0ntuuF8j2crOZecEVcCMv2e1TT+e3lpYDvbq7r2BbZlsTv/4y2eNr5/jONplqVz6ryde0we++Wtt7av8PbWxhXf/+dgruujWdfr6fc9VfNzxKqLJS5d2Wv4+2ZYcDo7ktH1PRiO5vPYYgrRf/u4Nfe5Hm5TJBv9OAACAqenHj7+qXzy1VZte6dObAxnlx9oIOI5TugMhADSbpl8WN1kcx9FvX+3XW+d0qq2lsrtH7LJUNtj07R1VtS1rSoGW4crevayo0gt020W2jeM41uqP/3ridT3/5qAefm67lr2lO/Jx+4ezenlHobprTyar2ckWz/ZNr/TpyVd36//zjv0DzZbd5+I4jmf7aC6vT93xuI47oEcfeefiyONxH7fW0ccOy2v/P7/fpn/45TN6fMsu3XrROzzbBkdyGhn7PymmMd326Mv63xd36pHnt+vcYxYGtm8bGNacGa0Vj9cWAr24vVDFt32syqsStr9P/75+i/75vuf0rQtW6K1zZlR87LLPaQmQxj/Xhr9rIX9frr/7af1h64D++Ij5eteBwSWU/s8lAABofp//8VPq2zte0Vxsd+DuVTm7o0Xzulo1r7tVbamE9o7mNDyaUywmzWpv0eyOltL/JhMxDWVySiVimtfdpgUzW/WWmW3qSCc1mMlqMJOTI0f7zUhX3NYCAKLi2yWiHz3+qj51xxN65wGzdfv/Pb6ifXdUGTbs9IRAQZlsTpls3nhr3p17LJVLe9yhSmVjKl1kV7jf0EhhrGXHVFq2VeF4fAGR31U/fFJPvdavI98yU0cs9IZW3n0l9zX841t267+ffF1PvtJnDJdGsnk5cpROBoNGWwA3msvrgaff1LFLZqu7vbJbKtte+2LV0V5DVZjts+A+rinI/P5vNuuzd27Sdf/nSJ23YlFF491pGe+OkDGVM5LNa2CsAsj0WfmvJ1/XSzuG9JsXdwbCpVze0R+29uuQeV0V96wqjddwMrYwceeg/e/a9rHjmpav3rFus7780z/o1ouO1TH0awIAYErI5vI6/oAevd4/rNd379WbezLGG6DsHBzRzsER/e71/po+f1drUvO727TfjLSGRgrB02g+r5ltKc0fC6bmdbcpn3c0OJLV0EhOramE5sxIa15Xq+bPbFVHS1L9w6PKZPNyHGnOjLTmz2w1/v9iANMH4VJE//LQi5KkR18I9gtyHEf/fN9zOnxBt049ZE5g+y5LZcPDz23X3/7k91p77hE6etFM436FfYNj+j83/Vqbdw7p16tPC/wrxE5LpcWukNDKJqzfTS7vGC/cwy6ybcu2bNzVPKZdXxtbQmZakmQLgYphgmmyz+Udnf7VBxWLSfd86mTF4+aKKNO5/HTT67rs9sf1kXcu1t+cU9ltbktNuQ1naguePMv/bOGIYVuxX9Jzb+4JbMtkc/q39a/olLftp0Wz2z3bsrl86V/kbJ/BisPEkKWMtqV4t6/brKvv/K2uOuMQrTr5wKqe1165VH6baV/3slfTa/8/v9+mvr2jeuzlXYRLAABMEclEXN/4yPLSn0dzee3YM6JYTErGY0om4srm8nqjP6M3+oe1tX9YmdGc2loSak0llB9rd7BraKQUQOXyjtpbEspk83qtrxBabXNVh7e3JOQ4hX+E7B/Oqn94QE+/MWAY3e6qz6tYUTWzPaVU3Nt5JZ2Ka25Xq+Z1tSqbdzQwPKpszlEqGS9VWfV0pDWjtXA94/5/RcUwbGZ7SoOZQtiVyzvqbktpZnuKCm+ggRAuRVQMP0yeeq1f19/9jJb2dpjDJcuF6Y8ff1W/f71f9/7+jUC45A1OvHtnc3ltGuuvs2PPiNpne99K20W2fzmYWz7vaPUPN+ngeTP00ROXlt3XFHD8/LdbdcW/Pa4bPnC0Tl82zzuekECr2r5UtgAul3e0e2/5QMZ24W8LMXYMZvTC2LKuTDYfWCa5y3IuxeVgxUott9FcXmvuekonHdSr05fNN4x3rLrLVAlk6wkUGnCUD4F2Wnp+/eKpN/TXP/qt/uSoBfqnD7/ds61v76g1KAwLKTPZXGhVmCmQsVURPbO18H+iXivTs+qm+5/Xew+dG1iWGdazyvZZsY13YDhbCi/Dqp4AAMDUlErENa872H6gpzNd0Z2i/TLZnEayeXW0JBWPx+Q4jvZkstraN6zX+oa1fSCjjnRCnemUEvGYdg2NaGvfsF7v26vX+4aVjMfUnk6qoyWhwZGctvVntLV/r7b2DWtoJKfutpRaUwk5jqOt/cMaHs2Xwi6zyvqqRtGaimteV6tmdbSoJREv/P+psf9TFY9LvZ2FaqtUMq6B4VFlxnrN7jcjrfkz29Tb0aLutpQUK+znuI7b05FW74y00sm4cnlHo7m8cnlHralCyAcgiHApItsFZLG5talRtefi3Rj0WKo7LMHJLtfd0WwBUqWhytNvDOiO9VvU29liDJdsodVDz72poZGcNry8MxAu7bBcZBcu3ssfd8PLO/Xz327Vp973trIVWoVz8e4cGnAM+cc0/i8ftmqenZbQTxqvejJXaJX/HP3vizv1vd9s1mObdwfCJU/AYdh3lyVU2WH5HLkDDhPbcV/dtXfsGMHgNaxSzfb5/Lf1W/S5O3+rb35keSCsDavmG79zm+F9CelZ9dVfPqsntuwO9Kzak8m6elZZ3m9LRZlpSGEN5W3LCgEAAGzSyYTnH+pisZhmtKY0ozWlg+bWri+lVPj/qTsGR7R9T0a7BkcD//95aCSn1/uH9UbfsFKJuGa0JpVKxjU8ktOru/fq9b692jU4qoFMtvT/xmOxwv8H2jU0om0D40sH08m4kvGYBkdyGh7N66UdQ3pprP/qZJmRTqqnsxBMjeacUquStlRCM9tTmtneopltKbUk49o7ktNwthDGze1qVW9nWjPbU0rGY3KcQluKvOMoGY9pdkeLejpb1JFOKptzlMs7yubzSsTjmt3eoq62JJVaaGiESxENjZS/m9NOS1VD//Co/eI96nKmMvtJZQIk23EjLAfLGsZcnDgKz2kYk626wxJo9Q9njc9XdP0vntGvX9ihY5fM1srDvaGV51z8F+/u6qAKx7QvfXTs72n512j7nsJ4TXcIcTflNhkPVQzjsYRhO0JeoyiN4UNfI1PoYhnvI89t10gurydf6QuES7ZlkHtHcqWeU9aleBX2rLJVLYYd1xZo2T67heed/LviAQAAVCoWi6m3M63ezvSEHD+XL1RddbQklEwUltwNj+a0tW9Y2wYy2jU0ouxY/8pYrPDPxSO5vN4cKCwtzOYdzWhNqTUVl+NI2/oL1Vs7B0fUv3dUxX9iLuwb09BoVtsHRoz/v1CSBjLZUg/QyZSMxzRzrGfraM4pXTN0t6U0q71FszpSSicTGh5r/p53pBmtSXW0JNXZmlRnOqlEPKbh0ZwS8Zh6OtKa3VkIwjpbk5JTeK3zjqO8I3Wmk5rZntKsjhal4jGN5PIazTnK5fPqTBeWJfqruLgZzfRGuFQDtotAW5WQe7tx3z3h4YcU0sPINN4I1RSmY+4drb4pt228YXfUKlaGmZ7bHhCNlt1WGG91Y9plOe7waK4URJpfB8vStj3lg5ydIU3YrZ+jfQkpLRVltqbcts9YWFPuUgBX60qgGvSsqmX1oS2A8yzpJFsCAADTWCIeKyxfc2lNJbSkt0NLejsm7HkHM1llc44SiZiS8ZgS8ZiGRnLasSej7XsKwVRLMq50Ml6qxNq9t9APa/fQqEayebW3JJROxrVraFRv9A9rx+DI2OoKRzHFCoFWrBAW7Rqr/hocySkZjymViCuZiCmbK4Rr2bxTuhmM29BITq+P/UPpZEsn4+pMJ5XJ5pXJ5jSac9SaKvysI51Ue0tS6WRcmWxe+byjWR0p9XSk1dWWUldrUopJuZyjnOMoP7bssLs9pZltLWpNxZXNOWOhVl7pZEKdrUnNaE2qqzUpxymcezaf18z2Fu3XWThuMUQczeULFWCOo46WpFpTcYKvCUa4VCHT59Ha7yaw9MrL2gzYsq8ttMpkc8YG1qZ9/Tvbmo/vS1Nua7VUSB8da0+giJUh/sOO5vIaGB5/jYIVPbbKr/FqH/94wxpORwo4TPt5jlthOGIJpjwBnG9b5KbcFVbs7EtTblsgE75krrqeVd7327vN3ZTbWPkVtReWb9/dQyOln5EtAQAATL6OdPBSubstru62lA7Yb3LHksnmtHOwEFoVmr/Hx/pMOdo9VAi0dg2NaDTrKJ2KqzWVUEzS4EhWezI57RnOak+m0Ei9rSWhbK6wImXnYEa7h0Y1MJxVPC4lYjHF4zHFVGgNsWtoVLuHRpR3pFQiplQ8rng8pj2ZQmuNQqjk/f+7w6N5DY+OGIOwekrGY+psTaotldBINq+9YxVcvZ1pze5oKQWBxaqwbM5RLFao/upMF6q/kvF4qTKsszWp2R1pzWxLqa2l8HqP5h2NZvPK5sfCsLGQrSOdUC7vlCrKZrYXQrbO1qRaEnHlncJzjuYKlWMdLYV9ij3T8mOvdSoRK1XvNSLCpQrNbm8J/Mx2MbzDUn3kvng3ibJsq/C8/ovsaP2YCsc1hyrm8bjHalmKF1Jx4t83rCm37Y5atgDJGkSEBRyWc7FVRHneb1NQFjE4Kbefad/ITbkN2z3BiW9b5KbcoeO1BKMVNuWOHiZ6996XnlW299vblDusctH3OlgqosI+nwAAAJg+0smE5ne3aX53W2Db4p7JH0+xQfzuoVENjmSVTibUmoqXwpc9mawGM1kNjuSUGc0pnUooHisUDezYk1H/3kLYJUnxsUArESss2du9d1S7h0aVyebUkoiXKrgy2bwGhgtB2MBwVol4TK2phBLxwnXq9j0Z6yqbbL4QxO2W9/p7YDhbuulSI2pJxD3tUQrLIMeXJI5k88pk8xrJ5dWSiCudio/1XCtU1cViUj4vtbUk1NPRohPe2qs/W75wQsZKuBSB+0J1VkcwXNppqTjxXCQGLiDty7aiVBiZ9vVegFv63aj8cjDzMrLyFRzuptwm1hDIEn5U1JTbEiDZgifb9krvDmbb5jhOtOAkpLKm2qbcpn1tlUuRm3IbA6LqQlPP8+5TmOhVs55V/qAsYlPuSse7I2QZJAAAAFAv7gbxjcJxHI3mHO0dzSkeK9yJMZWIKx4rLJ/bk8lqYHhUg5mc0qm42lIJjeby2rGncJfDoZGcRnN5JRPxQpVWIq5s3imEZJlCoJXN59XeUqg2GshktWNPRn17RzU8mpfkKBkvLJFMjvW2KoVsmUKVVLHCaefQ2HNmCtcoxfG2JOLS2HiL13f+a5i+vaPWApUw7ekE4VI9xWIx/dOH365P/mCjejsN4ZK1GsUdIEUPPyR7cBK1UbApBLI1Rbbden6XJQRyN+UOr5Yqv822n0nUnjaVHtcWAtmCHttrHzXgCK0EsjTlrvTudrYwLHJTbuO+5YNI2+sXtSl36Hj978s+9KyqdsmcVEkI7NsvpBcWAAAAgHGxWEwtyZhaksFlYx1jy9PmdrUGtr11TuBHk8rUBN1xCsvg9mSyymTzSicLSx0zozntGhrRzsFCT6+84xQqlFIJpRKxwt0Lx/okZ7L50t3s47GYBkey2jU4okPnd03YuRAuRWRr/TV+AWm6sM8Efja+Lezi3VL9EbGxsf+oYU25d+4pH3BErUYxHjdi8+xKq2dsQZntwt7/2rqP7WnKHRIY+A+8Lz2BrMv/IjflDtppqZ6LGjzZmnKb7LRU5dnGG1YJZFt2GP245Zdthr7f/uNalsy5m3Ibx2tdOmhf2goAAACg+ZmajMdihSV//rvxdaaT6pmguzLWQuN2g2owxffcVmEUWrlkq6zxbfM35bZX5VRy0epfDma+QDdXwJSvkLE1Hy8ct3wFl7UvjWVboCm373ltjcL9oZ/72OFNuasLkHZY3jP3cY0BUQ2acpvYwrvojeoN26usKIsawEn2SiB7kOPdr6KeVdZg1LvR3ZTbdGxrBaEllAYAAACARkO4FFFsrHbJfxEY2pS7gr4/nv0slTX+fSvp5RI4rm+7vSdQ+SqNnZbG5e6m3CbWqqZ9aHpsDyLKv7479kQLgUxjslYYRQw4TKptyr3b17PKGshYxxs99AvbbguBbPu5m3IX/lz+Of1sS/H6h0ej96yyLE+1VWgZj2v9bihf8QQAAAAAjYZwKaJStZrvOs998R7an8dW5eLbL1Bh5NtebS+XHf6KHXfw4DiRlwdZl+L5tvmbcleyjM++nC56ABdspl7+9bVd9DuOE328FWzz3HWs4kog1+tgqdgpjt+z3bNkrvzSrErGU9i3uvCzFj2rQsdbQdgVNibbZ2XHHv9rbxuv5X0BAAAAgAZHuLSPbNUdwe3ebfZlW/ZQwHoXqgqqfdwX2gOZrEZz5askola52PYL3TfkdXDzB2W25626Ibpv29BITiO2nlW2pVmWpYOeptymSqCId8WzjafwvL7tloozW0VU1KbcJrY7KFqDJ39YU8ndAWvUs6qSIDf4vozvPZL1L+ksf1wKlwAAAAA0uqYKl9auXatjjz1WM2bM0Jw5c3TOOefo6aefnpTnHi9csi3jCbI2la7gYti9eWgkO3a7Q7OojYL9B7ZVd/jHVMlyprCgzNrbpwZNuf3b/Mf1jynqsi3T9l2WpYP70sPIVoW1wxIQ2Sq0KmrKbavY8b/flucMjrf6gMheEeUbk6Vn1Q7Ltkp6VgW3lf987g55jcK+V2BXz7kCANAcmCsAoLaaKlx64IEHdMkll+jRRx/VPffco2w2q5UrV2pwcHDCn7tcQ29bhUEmm7NfvNesEih64GVrFGyrRqlkTJX0wgk05baOKSzEGN9eSVNu//Puy3LFHZYQyH7c6ptyW5dmWSqXAgGHdbzRK8oCnxNrBZd9vLbjunf196yqZOmlrXIp0LOqgrDRFuyF3RWPyqV9U8+5AgDQHJgrAKC2kvUeQCV+/vOfe/586623as6cOdqwYYPe/e53T/Czmxt626p5dg/ZewJVGwrsS6+hQDWFezyW8CMf0pTbdoEedT8prELGK7BMyvULlfW78W63XfjbqmfcPauM47UGcOVDv0DPKttxLX29/IKfMVsvIstzhr22vufwf37LjTdsKaP7ed09q/zPGTiu7e9EWM8q/3ENn4diEG37DNo+u/7jmqoIYVffuQIA0AyYKwCgtpoqXPLr6+uTJM2ePbvs72QyGWUy4xfu/f39VT3XeOVS9Iv3YHVH+aU6fsHjju8btSl38BlNy8HKjNe3Y9/eUeU9z+Pdbqt6sgY5FVRhhTXltm3zBxyBAKnMmEKX07l27PcHHNalg/Jti770qrIQyB8muiprQoIT+5KvCqq7XL+wJ5P1NOWuqAn7BPXfqqRnVSCsDfSBinhcS9Vd2JJOVG4y5woAQHNirgCAfdNUy+LcHMfRFVdcoRNPPFHLli0r+3tr165Vd3d36bFo0aKqnq/MzeKstyIPbaZcSd8ad+VSjZpy+/e19qyxPKfkv/APWTpU5jmlCvvoWMIRWyBjaspd9nUIazBuCcr8bFVutmbfFTVE9z/nWOjSkoiPHds8Hv82f1Nua8VOBcsgA1V3tp5Lvn2DFWWWJWiWv2uV9KwqnkspXHZtM/WsKte7q7CvJdir4HOEykz2XAEAaD7MFQCw75o2XLr00kv15JNP6gc/+IH191avXq2+vr7SY8uWLfv0vNZKlgqWB+0dyXmacodXyLiPW35ZXNit522BgnU5WAUhkJ/tQjoYNozzV3CEj8l2YR/cL5WIycTWlLuiQCtwXEvlTSVNuV2/4A84gqFK4Vxmd7QEntNamRQS1kRtyu3fXlFTbv9nd+w500lDUGYJrfalZ1VxvDPbUoEx+XtWFZ7XPSbLss19+ByhMvWaKwAAzYO5AgD2XVMui/vEJz6hu+66Sw8++KAWLlxo/d10Oq10Or3PzxmLmYMIzx21fNuCoUD5pW2B40aoRomyX+iyuDLbKlnaNprLq992W/UqA63QapSI5xJ4zrFtszta9EZ/JrDd2pTbNl7Lc+byjnbXqCm3bTzl3u+ezhZt7R+OHMDZlmVWOl5r+On676hNuXs6WvRa37DvPS2/VLQWPat6OtOesUnj55lKxIzVgra/M1HvVhjYERWpx1wBAGguzBUAUBtNVbnkOI4uvfRS/fCHP9S9996rpUuXTtpzl1sW576o9VdaFKs7WpLBlznYlNu8FKoznRzb7t4W3p9nVnux0mJ8m6kpt/ei1hKUVRJwVBACVdSUu4LnjdKUe3bH+P85KO4a1rPK2qR57DxnFl971967h0Z8vbF8Y6qgKXfU5X/u7abKpUoCjkBfrwhNuYuVYe4R1aIpd09nemy7+bPrty89q3b4Xj9TAOd+baMuX905FC1o9D8noqnnXAEAaA7MFQBQW00VLl1yySW67bbb9P3vf18zZszQ1q1btXXrVu3du3fCn7tUuBRo6B1eudRjuDAN9jBy/bcr4Bi/qK2sX4tpP3dT7mQ8WIllDcrGjtuWShjGY++jUxxTV2uwUM52AW4LjxzHMVT0WKo/3NvGAqLeTncoUNjub8pdyXh3+t+ziEGDf3u5ptzjy8HKfxb8iq+D6TNoq4iyhUuDIzlvU+4y+85qN7wONWjKHen1jVihFTamUoBkOBdzSFn4BfeSzo6WRHDfSqoPSZcqVs+5AgDQHJgrAKC2mipcuummm9TX16dTTjlF8+fPLz3uuOOOCX9uU0NfKdodtUwXw7aLS3dTbvO+4YFMz9gFr+lieEZrUsliVYnnAj28csk0Htvd69z7lipODBf+iXjwBbadp7spt6kyzLokyRcCuZ+2kiqXwHEtQYSt4kmy3zmw+FnpLVXsuI5rGY+7Kfdsw+eh+DqUXnrLa28ba7nXyFhhNPac7WOBiyKGie59ewxVWP7X160WPat6OoNhrT+4c+/v7uvVaQhVbX3Twu6giHD1nCsAAM2BuQIAaqupei75L2InU0zBMCZwR60yS516DKFAlKbc7S0JtabiwX0tF8v+ZTzG8bh6DZm2G8+leGHf2aJXd+/1LQ8qv8RneDSnwZFiwNGiF7cPGpd1zWpv0fY9mcgNkYvb0sm42lsSGsnmIzcK91fWuI/tvjuY44SHiaagp6fTsAStkoDD/5yuptyv7t7r+YXicbtak56+V4X9ggGHaSnZrPYW7RgcsVeN2Zpy+8brr9jzjGnP+OdzaGRvVU25x4Oe4HMm4jHl8o61ebZprOXOpRjs9VhCXtvftVntLaXvDs/zupZQ7h4aLXu3QlSnnnMFAKA5MFcAQG01VeVSI7BfZJep4ChdmFqWthkuLt0XpvZ9Fdg223ABXjpuR8t4JdbYL4yG3FErUKVhqIgyjmdsWzIe04xiBYfhAt0URNiWK7kv7IvN1k1BT/E5wyqX/M9pWtJlasptel/cy6TGx1v+tfXfdaxc5ZdpqaOtSsi9n2lVZyk4MYRhxeO2GO7MZluW6W7KPdvyWTEuFbVUwO0ea8odi0kzTUvU/O+b6zj70rOqGFLaPiumnkvez6d3m+M4kZb4FZeu8v99AQAAADQ6wqWoDDeLs11kS/YQoxh+2JYk9XSOX5gWuZtym5eDjVV3GC7s3cu2xhuUF35ht6XBsGdfU3ASocn1rI4WxUshkGEZX2cwgLOFVu7+Rv7gxN2zqscQyOzYY7jwH9seDD/G93PfdazYrNotsK8TPM9S/x33eELuzOb+PASe07r00hVS+lYdupty2/attFdT//BoqWeVrQm2eYlkeJPw7rZUaQmlKWTr7TSdS+G4pSrAKnpWze4MLivc4XuNPOM1fT7HxuvuWWWqiHJ/Pgv7AQAAAEBjI1yKyFT5scN1wSt5LwJDqxMsd6EyhQJF7qbc47193BfZvmU8hhBjlrvaZ2xz8WLY/XymoMfUe2aXa/mVf5s7pDCFQIG7cdn6/hirhIIBkbsp96yIF/7FIUfpk+UJOMqEgq5DShoPDGZZPguG0/QeN8LSLNPnqHCe3g+SO+AI27fceGcU72Ro2G9GOumqegp+HozBU4Sm3O5g1BicGqu7in8nqu9ZZQuQikGZ+3l32v6ujW1rTcXV1lK+yb3ptQcAAACARkS4FJH/AlFyVWEYqoT2uJpyG3vP+Bpvu3nDD+/zuptyp5LB5WClJUmGSovAeF1K4YehD5Hk7ZUT2ObvYRRSPVPkbsptCzhK4zFsKzyn7zUaHK8SShuWdXnH5H0NTUudxvcbr/Qx9dGJstTJ2H9rqDieVOCYUZtym6qlbEuziu9nayqutlT5pYO2ptzFpZee0G/IFaoosDlSU25jVZhpSefYNnfPKuPn07JMz92zKvCcrtB0fHllsNLKtizO0+zbdy7eCkL3uRaO22t47QEAAACgEREuRWS+CPRVo3j63RQuENtSCbUZlkLZ775mCHrG9t1lqEYxVrIYes/sNF34F8djuYOavyl3udfBFJx4z9Mb5HibchsCjrHXsFQZZqnYcW8fX8rkeo2M4w0GILbqo/FQKmXsWTUw7As4bHcWi1p9FLEpt63JtadCK7D8Lx2o/DKO1z0my7ItT4WWLwRyN+W2BXC218EU7LmbeXe1BQO6XYH31HCelp5Vnt5n7jEZelb5P9vesNZb1TS70x36Fba5l72awjAAAAAAaESESxUy9dGxVSaZwg/Je4E5tnX8uIaeS8V9d3iqbuTZd9RwRy3PmNzBk++i1hzWFLYV+zEl4jF1taY829z7Gl+HUgA3ftFvbXpsWMZn62Hk3re0zV0Z4jtPd1NuU3g3PqZ0cUPgXGZ3pAN9dIrPGY9JM9si9EYyBjnBPkRRm3KbemEZQ0rf0qxZHanAcStuyl2mksofArmbcpsab/tfe+MyUsN76g79Eobm7lHu6mYLecs15TY1+y6Oyfb59Dbsl2e8A64lnabxAgAAAEAjIlyKyH+BKAVvPW++6A/2XPFUJxjvfGUIP3yhgKmHUTEE8ly8hwUVpee0LFcaO5dZ7S3jDchd24M9YgxLqAwVMqam3OPnG6zgKFch4w96PBVahp5VxeHNdFUg+cdr6lk1Ph5DUDY0HhjE495tUrS74kXtv1Xc7g44bA3co4eUBaaAw1SNNt7DKFhJ5Rmv77Pr7VkVDNl6DRVG7n5hRf6ljJ5zsYSfUXtW2ZZImppyu/c3N40fO64nePJuLH6OOt09qwLPDgAAAACNhXApIlOPHXsoMHaR3RGsjukfdjXltjSOnt2RGr/4LG4bCl5kl8Yztm1mW2o84Cizrz/wMlYulc5zPMDwj8fYlNs9pmJ1R3uwQsbYlHtsW39IBYenQsbfT8iydNDdlDuViAe225Yrmpo0G7cZllDZmjSXXj9jSFl+OVjUptze3lLeSitThVExTJyRTqolYWnKbaiOcwdwgdDPUH3k2TdCE3FPkOMPygyhaS16Vnn7PI2di7tn1diyV/PrEPwM+v++eI7rDqVL50m8BAAAAKCxES5F5A8wpIjLmQyhSvHickaruzrBXSFTXG6TDlYYeSoivCGGt5rH1IDc1ZDat6TOtMQvUJVjWELlbsptvK16seopdHmQr4+Ouyl3KhE4bpR9Z7cHw49AE2ZfVUlgKZ5hvKZGzOMhWjDsGh7NacjSs6oWTbnbi329LHdm87wOhuou/5IuU/Ns977G97v0+XTdQc2yDLLI2JTb85zjYW2gCbsnaDQ3wG9JxI09q4rjtfWssv2dmO37/I2fa7CCS75gz3PnO/n2c/esIlsCAAAA0OAIlyLyV2FIIbeed/fnKdcTqExPlfGL8FTgwt+0lKxUdWOo5inKZHPa47t4Nz2n6bbqpvH6l9Olk3G1G25NP171lA5U9Nj6PLmreWzN1G0NvY0VJ67+PJI81SGeptyd9qAn2LPKvQzS+yoV9yvcdczUs8q3tM0lSlPu2YYlc4F9ywRTsw1LL01Vd0XRm3KXX3JoCj89TblbyzfldldElY5rXDLn/eyaglHveO09q/yf3V3uz6d7QE6ZptwK+Xz6K+faU8YKOAAAAABoRIRLURmqCHa5+uwEtnkuhs2hiqkvjbspt7fhtO+4lqbcs3wVO+7ePIWL92TZu6SZLoZNjY3le85yy3g8VU++fY3hh+E8bU25Tc3AjbeBD5xnIUxwv/7F8ZRrym0KesbfU8MySH+YUCYEKgVThrvtRWvKPX4nviJ/U+7AZ3AoGI6Mn8vYXdAMlUC7XT2rbE25wyrKbKHfeD+m4HE9vbv8n0FjULbvPatM4Zw3eBrnyPH0rPKEWr5z6XG9b/7KOVPwBAAAAACNinCpQsXrPPfFu/HObKam0r4L3h7DRb+7KfdM94WrqerBNyZTw26pcHHq778TWFpkqaYyhwK+8zT00SlcvJcPkIxNuQNL8YJLhwJNuX3bbb2lTE25S/uNjXWmuym3a7tnvJbj+kMVUwBX5LnrWHv5gKPHEALZ7mbWPzxaCjhmtgcbkHuXDvrH61rSVdzPN56Z7fam3N4m4t5qNFMAt8sY1hiqxgxjsi23s/Ws2pPJRupZZVxG6tnmXRbnbsqdTiasfdMClUuWO1ACAAAAQKMiXIpofBlP4SrQffFuvcW552LYu5zJeJHtasqdiMfKXpia7kK3o8wFr3s8xYtWd5jjbsrtDsrMQY+5AXaPYfla/3BW2WIFh6nyxn0hXWa8s9uD51msrOlqTRaacvuDngh3xZtVeh3Gn9fWCDww3vGzkWReBunvCWR6zzx3HTMsr4xUCWRsDF3YNqMUcJi32wIO713xglVCtqbctmAqtCm3bzzuptyz3E3uA+M1VGhZArhiqNeWSpSacpftd+UP4AzL3orPu7PUZ8y/9NL7Ong/n4b3xbcNAAAAABoV4VJE5Soi3HfUksYvTo0X/qaLS99Fv7s6Rgr2etrlb0jt3mbszzN2wetanuY+H0eO9o7mlMl6Kzg85+pqIu4vvTHdkUy+8+xoSag1lQgEdDuNQUT512/8OYsVY94+OfaKHvMSKs+yOENT7iJ3U25Tc3JTT6DAeZrCrpCm3J7ePWX6b3leP1/g4v8cFbmbcvtfe3dT7uCSuWDQWBRoyh1orm1qyu397IY25R7r62Uarz0os/SsMrwv7n1Ny1fdPauCPaDGG/JL8rxv/iWdpREZKspYFgcAAACgWRAuRRSowjAsvZKCF8Sm5WDuO5b5l735q2PcF5ieptzGi1rzEjX/8jT3GbmXzLUkzRfv5l5DwZ5L/v0Cy6BsS5b8F++WCo5AU27XxbunKbehZ5V/vO5j+ytOSvs6jrcpd7p8zypTbx93xU6gZ5XhrmNll2b5z2VP+UAmEFK6dnY35Q6rXAosy7Q0CS/XlHv8fSs2Pd+HptyxmHUpWSAoi9SzKth8vNyy18By0A7L3zV/03hJu4dGxpd0tgWf1xOMBl4lAAAAAGhMhEsR+atyduwJhirFzdlcvtQ7ydRsuXSRbaiQCTblHg8cPE2525JlezkV7rDmrdPYYa0ECgYckjk48W+z3ZHMvwzKPdyyFRz+MRluyW66E1dxX3dT7q62YIjh39cdyrgrTjyVX44/PAr2rBpfdpj27CeVCTiK4zEsBysKNOUuU9Fj6ksVCCldnyN3U+6Zbe4QyHsunmo033iNyz1doV88HgsEcLvcr2+ZMNEdCBb5P3/e92a8Z5Wpv1m0nlVpucNWKdizyl+N5g0MvX/XTMtei8cuvi/dbSklE/HgXeiMf58oXQIAAADQ2AiXIrL20XH9nuOMhyaxmLc6YXw5k+ti2NIjxs0YcJS7PXrgbnHBbe6L+3INkZ2xMZeCik5T/6PyAYf/Itu9b7Aptzk4CeujU+5cZra3eHtWlanocZ+rp2LHH6r4Q6nifs5YZZNn2aF3+Z8tgDOFfsVt7ruOeZpyF8dkbWTtveOg+zPobsrtDjiMdwD0L9uM2JTb/ZzjYypfpRWliq1UQeT6HHmaclsq1ex3r0sFzrO4X2egZ5XvuJa/pz3+z72c0ufPX5nof+1nGQJXAAAAAGhUhEsVGu+j47qjVpn+O+PVCWP7+rZbewL5gxM5rm3+5WCF//U0GfaERE4w6HFtK1sJ5DgayGQ1mnPGz9V3Lp6KkzLB0+z2YOVScdt4U27PU1tvAx+1Kbf/PCVTRc/4uZa/vbyhesb1vgyN5DTi6llVLogwLXWyNeXeMVbhVmrK7RprueOObxtfgubmPhd/0GgaU7D6yPX6la0EClZLBZtyy7NzlKbc/solR+NVS8Wm3OUa2Vt7VpnOM/B3wnfcstVzjuHzOS7QC8v1Go5kx5d09hjGBAAAAACNinApomCo4u4fU2YJWruvcsEfuhiqE8pVyEj2C9qhkayGR4MBR2lf3zI+d2BjreYZG2t7sSm3rW+SazySOTgpPKfjqkZJ+84z2Cun7HED4YjjacrtP09/U273E7vP1b2Eyj/e2Z3+cxl/DQpNuZO+GKJMwOFrKm1qGl0uiCgd1/B5GP98+ppKG8YbDBorbMrtD8oCSw5dz1mjptzljhusKDNUjflfI18Vkfs5y/Wschxfzyr/Z9BSJSjHG0oXxjv++u52L+lsTQXfcAAAAABoUE0ZLt14441aunSpWltbtXz5cv3qV7+a8OcMVs+4Lt49QYQlVFGwKbd7v8Jxy1TIGI7r7Rc0dvGeLNx1zN8zKBhajf+G905d3vP29wvyV8h4ewKVORdDOFK+Kbc8Tbl7fK+hZ0yGJWrBptzjVTD+ptzufb1jCvasCi45DL72/jBRvs+K8a5uht49/oDI9Pr5Aw5bU27PmGQIZFyfo7Cm3J674pWpugtW7DgKNOX2vw4RmnKb/j4F3m/X61BRz6pyyyDbg83d3T2rZhmadtsquMpWlLmXVra3eHtWNeHCuE2bNimbzdZ7GHWZKwAAzYW5AgBqo+nCpTvuuEOXX365rr76am3cuFEnnXSSzjjjDG3evHmCn7lclYb3jk+2JWhynFKj70Q8phmtyfLLeDqDF6amxtCl8bgqdtwNp/3HNYVW3vDDFQk43moev1zecVXXGHoCWaqwbI21A025ywQnwZ42ESt2XI2qTQGdsd+VpYqo/JJDeXpWmQIOd1AWrHjKeM/TXeViCDgKz+kLegxLEsv25nIv6SrTlHv8NUp7junZZuj5Ve7zp+JrFKEpdyCscYLvt/vvS+SeVcY7L/qPO76vf9lr4XyCYa25J5jv8+nKIW3Vhc3mqKOOUmdnp4455hhddNFF+sd//Efdf//92r1796SNoX5zBQCgEuvWrdN73vMeHXnkkTr33HP1xS9+UXfdddekfF8zVwBA7TRduHTDDTfoYx/7mD7+8Y/r0EMP1Ve/+lUtWrRIN91004Q+b9nGxqYL6T3ei2z3hXRpuU2gOqGg/JIax9vvRv5qiWIFRzAEyjuWvkpyvEv8fM2GbKFUnyfgaJG/KXeggiNKkOMKOEpNuf39bizNlIMVO+PbTa+Be0mdtWdV2dfeCYQC4+Nx1L93POAwBnCmptyGiif3WB1XJVCpr1eZqrFiKOiuMvJXRI2Pt8Km3IEKLVuo4g8px0Ogck25/c9pGpO7Kbf/OSvqWRWowvJWGLnP1fTZ9VbPlVtm6gQ/n+7zDISU3kq2ZvLQQw9p9uzZWrp0qTKZjL797W/rtNNOU09Pjw4++GD99V//9YQHTfWaKwAAlfnIRz6iRCKhVatW6YADDtADDzygiy66SEuWLFFPT8+EPjdzBQDUTjL8VwrOPvtsHX300aXH0qVLJ3JcRiMjI9qwYYOuuuoqz89XrlypRx55xLhPJpNRJpMp/bm/v7+q5y67nKmz/F3SZvsv7FU+/PA3GTYuASrXDNh9RzdDldDA8HhTbmMVkfsW8Z7ndIKBgWFblKbc5fa13W0rGFL4QgF/VY6vKXfZ5zQEGHsy3qbc2bz3kr7cEj/jkkN3NcrYts5yAYexkbW/X5B/CZWpcbmvss5SuRSoiDKM1xRw+Jtyb+0fLm3zjte3PFC2kNLQlLtMQ+/gZ9spv0RSTkU9q4r9ysouSSwdt0y4NPYLtiWdhcDL/zq4KvYslV/N5tJLL9WNN96oc845p/SzBx54QB//+Md1wQUX6O6779Ztt92m//3f/9V+++1X8+ev51wBAKjMli1b9JOf/EQHHnig5+cvv/yyHn/88Ql7XuYKAKityJVLBx10kB5++GH93//7f3XggQdq5syZOvnkk3XZZZfp1ltv1caNGzU6OjqRY9X27duVy+U0d+5cz8/nzp2rrVu3GvdZu3aturu7S49Fixbt0xhMt10P9GsJhB/Bih3TRb+/Kbfkq/4ot2xG9v48xWqpYlPuwnbXmAxL2/zbShfShmqUYGBlXprl3rd8U25Fbspdaq5tqGQJVuwElxy6dy0GLulkXG0pQ88q35jcgu+p4f02VFK59/U05S7TVNq9vfxdBQ1NueWrvBnyByfBECgQPLk+C8Wm3P7Qr3xj+GDlnCesKbuMzBuyGYM933jdIvessjQR978Oco3XFFIWw6NSU275lsz5P5+u12GH7/PprqprNn/4wx902GGHeX528skn6x/+4R/02GOP6b777tOKFSv02c9+dkKevxHmCgBANCeccIK2bNkS+PnixYt19tlnT9jzMlcAQG1FDpeuv/56/fKXv9Sbb76pzZs367bbbtP73vc+vfrqq/rSl76kFStWqLOzU0cdddREjleSAmGO4ziBnxWtXr1afX19pYdp8qr0Of0X7/4gwl9NUdrmXl7l72HkCgyKTbklc+PtQLWKq5F1aXmQITAwXQz7KzHc5+npz2MJIkwX/aO5vPrLVXAoWtNjU1Pu4muQjAebcrvPJRD0yFwZUjzf7a47h5l6VpUdryXocb9G/p5AkqUp99j2cr2lpPIBh60pd+HYjqcpt3u8tjvJydKUu6hcIOP+O2HsZxV4v8dfB39T7sKYXJ8VS/AUpWdVLCbNbAv2rLI1sve/fu7tpXMZW/Za2DZ+PuX+LnqO2+5/7ZvPscceq9tuuy3w88MPP1x33323YrGYrrzySv3yl7+c0HHUY64AAIQ7++yz9fnPf17/+Z//qVWrVumLX/yiduzYUZexMFcAQG1EXhbntnDhQi1cuFBnnXVW6Wd79uzRxo0b9eSTT9ZscH69vb1KJBKBf03Ytm1b4F8ditLptNLptHFbJdzLmfwX78PZXOn33FUapVDFdPHuXxYjR/6m3N59g8u63D14Ag2I3ZVLxYtsd8XO2PZc3indAt0flEmW5UwyhAmGEKjYlNu9XbIsJTM1PXbx3rJ+7OLd0svJ3bPKFPoVN+8qG6qYe1Z5X/vy/aNsPYH6fE25i8cJ3H3NsMSvXEAkQ1Nu73ZT4GUab/Sm3KXG5YGlZIZgz7hkzv/ZHX8/+4dHjU25i6/TTn+Vm4LHjdKzyn2eMoy3tF3uoMz1+Rw79g5D/62Y4e+wMYj0V5T5qqmayY033qjjjz9ezz33nD7/+c/rkEMO0cjIiP7hH/5Bs2fPliTtt99+euONNybk+es5VwAAwh100EF65JFHdNNNN5VCpYMPPlhnn322jj/+eL397W/XEUccoZaW4D+M1gpzBQDUVs0aend2duqkk07SJZdcUqtDBrS0tGj58uW65557PD+/55579K53vWvCnlcyX5yP31HLVd/hWkIV6Ndi6c8jjYdApkbBecO+bqZqiuB4gxe8fXtHlC8FHC2e8Zguhk3bzFUYhdeg2JTb/Zy2Chlj02P3Mih/cOd+Xn9TbtnDMPf2cucpydOUe2Z7sBGzNSgreye58devFHCELJmzBXDuqhx/vyXP87oDuEDTeFMT9uB4Ak25HXtT7rDloIEAzrWtuF+xKbfnF2QKrYrPGfwcmcZjWpZpP67579r48kpz1aJUWNJZ7FkVWB5oqrTy3VmwmRx++OH69a9/rddff12HHXaY2tra1NHRoVtuuUXXXnutJGnjxo1asGDBhDx/PecKAEA4/4qIu+66S5dffrn6+vp07bXX6h3veIc6Ozt15JFHTtgYmCsAoLaqqlyqpyuuuEIf+chHtGLFCh1//PG6+eabtXnzZq1atWpCn9dUhVFumVkpJPJVJ0nl+xSZgif3vram3LY7VJme030+xUqLGWNNuf3nElySNL4tcFxXkGMLyjzVPpaAI2pT7uL2QV9T7nKvg/c18i6Ls4V+HWV6VpW9c6DK9wSyNeUuLYuLsHTQdFc84xK0sV/YO+ptyh0Yb5mll8Ylh55KIG9Tbv+YyvVyKldRVtgUbMpdGNP4eANhomvbeBWR918YzeHm+HmaelaNv0bOeIWRIaz1f67d24rnmUrE1Flc0mmo4DItp21Ghx9+uO677z69/PLLeuKJJ5RIJLR8+XLNmzdPUqFyqRg0TYR6zRUAgMrUa0WExFwBALXUdOHSBz/4Qe3YsUNf/OIX9frrr2vZsmX66U9/qsWLF0/K85v6x7gvoj1NuQ1Ne4OBjHtJV/DicrwiohAmuJtye/a1VkSVr8rxX/QHznVP+Vuyl7sbnGRugF3cdXg0H2zKXXxORWzKbaw+8jblLozJXdFja8o93ivHvV9h24hxrOXHO/4aBYITd9VNuTsDOuWaco8/qT/g8FQCmUJK37kUm3J7xivDZ9BYCWQINy3PabtzmymA84ZShr8TY79QWNLpH6/r/Q5UwLnO01JtZupZZezl5GkMPxbWGpvGF7YVG/K7bwIQdUliM1u8eLHx+/mkk06a0Oet91wBAKhecUUEcwUANI+mC5ck6S//8i/1l3/5l5P6nJ4qDcvyNPfFe4evgsPYN0njx/Uvi3EzVccYlxYZGm+PB0TlKy1MFU/Z/HhTbn8vJ3fvmfHAIBh2zbI0PfY25Q6GOdam3KYqLFdIYepZZQrZyi2LM/YoMoRSubyj3XvLBT2mJs3j5xK8Lf34cW1NuaVgoFXk+YwZznO8Z1UqGHA4wc+gpxKowqbc3sCmTFNu69LB4FIx9/biefqbcsu3PdL7Unrty/WsGn9zdg0GK5eKSn/XLJVL3ooyjb8OZe4W1+TZUl3VY64AADQX5goAqI2a9VyaLoxLagxBhPvi3S1KI2vbcptyTblLxzUFFabjjv3vDkNw4n/OWKzQF8i9n+dcDNUzO8ssSXIf19uUu/g6WJoeO04giPCMd49h29j/lutZVdzuXxbntsN03LEddw+NlCpLZrYHl6iVr5AxNOU2Pac74Ci+DIblge4gwhaGGftkefoJ+ZtKuyu//OcSrHhyv9/FAKlv73hT7uDSQdMSyfGTifJ3wt2U2/QamSrgyi/FC+lZpWAQaRqTufqw/GfM1LMqeD8+AAAAAGhMTVm5VA+mi3f/XdIk9wVk8CK73LKu4nEjNQo2XPB6mnIblqgZ+x+NnZCxqsl3njPbUoam3E75SosyzZT91TPmptzmJValbb7mz27mi/7Ck7rvOmaqrrFVo5lfP++5dLelSj2rbO+3N6Qw3yXNfS6mUMX8ObJ/xopj8jeUdm+zNuU2jMm7DDL4voy/RoXXb0Y6qZZkXP6dyzbPNjR39xzXEoyaQiJTBZzp73C54ElSqaeXe7zu09lh+DtcHJQ5ePJW3bWm4oGeVc2+LA4AAADA1EflUkS2ZVtu402EgxfDA8Ojpabctt4+pv4y9r5JhZDC3ZTbHVQYl4P5ntMUyNgqLaqpkAlUd7QHn3Mwkw005XaHGOOBVrAiKsp5drrvOqbx17/YA8q0RM3YpLm4xM8SSnkCGX+1j2spmbFyqbTMzH3c8YAocNcxY7AXfO2LjEHZ2Fg9Tbld+5Rryu1ecmgKZEzbXMVo5ZtyO+WbcksyNm8vvi+Re1ZVeFfBIndT7sJ2/9+Z8p9P490KS4GrO5Qujol0CQAAAEBjI1yKyNTw19Y827YEzXRHrbCmyPYQKFhZ47a71O/GXU3hG6+pCsO4VKywdXg0Nx7I2PoUGfY1V9Z4x+Npyu1ZmmVqFF4M4LxNud3GzzNY8eRWWkLl6VllWYrn6s/j35bLOxoo9qwyVQL5m3IbXj9jEGG665j7XC3LA4uMr72lN5fpzm22SiA3W2haaVNu9/MWmd6XsJ5V5YK9sj2r/M/pasrt3r/UqN70939P8PMZqEw0fT7JlgAAAAA0OMKlijmB0MVb3WFoyj32CxljpcX4cU3NquUKKsrtawqBTN1abP2PbM3JbdVHyXhMXa3eptyhd3WzLUFzLZnz96xy5BhDotIFuul18FWG+EMef0WKKegxVpz4xjvbUHFSFI+5elYZloP5m3JL9t49O0x3HYu5gp4IlTemAG5XKeSJ2JTbEIyaAqRi4FKLptymczEFOdtDelbtCixJLG4rU6lmqfwyibr0crxqzBA8ucYEAAAAAI2McCkie+VSSJWLL24wXVyamnKbmC4+bb1n3H/ubgsGMsZKC98Fr/eW6wWl5X2mptxl7qjnXu5U3Nf/nLuMVS6F/3U35fYsH/K/DhU0U/YHQaWm3K6fmcZUtMNXAWM66Kz2lmDAYbxzoGEpo+GzUOqLZAo4DE25/ceWzI3hxyuX7E25Z/r6XZmacrvPtcj02S2+n1GbcpuY3pdSlVVYz6rS9pBg1Hf8wGvv+wXT8rbi3xnTMj5bnyyHpksAAAAAGhwNvSPyLgHy9Wtx/Z6tEqjIdNG/e+9ooCm3e3uR7cLUVDVS5G7KbdoefE7HUyFTbj9T+OG5G5e14iQYdvn75IyPR+ovE3AUj2usDPMvmfOHS65BeZpyh/WsitCkucgUlBnvOuauljIuUbN9FgoGR7Klz1HkkLL4+paW8Nmbcrt7VhWZqu6C1T7ByqRS4BKxKXdhuy8oi3jnwOKr5O5Z5V+GVq5nlV/gc+T6b3dTbtN4Tcf1V4W5zwUAAAAAGh2VSxEVL/T2DBtuGW7olWMLVUx9aYoX0e6m3KZ9TUvmTP15/MIqdmzHtV9ku4OIwo7upty2ihN3U25rcOIbT0dLQq2pYFPu0pgsPW1sgUu58+zba7jwH/tf053v7EuoxgKOkfJNuT3jtVTPlKvukrx9vUw720KMyE25Dcvi7NU+wYoo43iKTblzwabcpuOaQtXifraeVS2J+HjPKk+wZwimQv4Ouz+/gaWX/vFGXOJXRN0SAAAAgEZHuBRR8RqwGCwFLt7H2PrzFHn63ViWdJn2NVV/bN9jrsrx7mfvNRS5ibjvuKblP6am3IXt5aueovQEKr1GlsDKP175Qhf/vu7nnWWopCo3Xn8AZxqvaT//crCyTbktYU5RjyGQKT1nBWFiYJvhXIpVYaYKrVzeKQVwtuo50x0Ui2xLJMs15TYd18/as6ojNb6kc2xb2Z5VlgDTfz5hn0/TkkTjeF2NzQEAAACgkREuRWS/lfv4Rv8dtfzbC/saloMNBys0TM9rqsopLQezXPAGLoZdF8sJV1Nu935RAg7T3a3coZTnjlqWfaOEH+WqjypptmyvXAr2cSpyN+V2j8lffeTeVmQK/Uo9qwxNuSX73eJKx7WFc773xbo0y7+ts/xn1zSe3XtH5Tgq25TbOF7ftqhNuQ3Dtb5vtp5VxqV25XpWhXzG3JttwZPk71lV/vUtIlsCAAAA0OjouVSlcre0L1YZ2HvPWJaDWUKgWEyaWcEdrGIavzC1Lm1rD96ZTXItBzP05zGNt3iIUpPwkOqZHsPd60rHtSxnsvW7kXz9mCz9j9zHLjxncIlfaT9XU+7w8VoCGct+bruHwptyW0M032tre40C+1bYlLvYB8vdlNu9vcizRNL/2TVU+5iachd2tXwGrc3zfe+LoVF9uZ5VUQJD03P6x9Tp61llrVwa+18aegMAAABodFQuRea/aPVdvEdY1jW+b/mKHdtFa6Apd2hw4uoDYwlVAv1jrBfoKr8tbGmWv8m4O+Cw9lwqbCvbaNm1b1clPat8Y7a99rbXNjBey76V7CdV1yOqsJ+/cmn8v/1Nue3j9T9ntFDKeNxahZQRqufGx1tZUFauZ5UtMPRvt91JzvaeBcY7XkwFAAAAAA2NcCki253OJO+F64x0Ui3JuHGbVFmzb8/Stgr6JvmfN2pT7ijH9e5bWd+komBT7vL7jldwjI3HEmL0dNordioK2dyhQFhwEjGQqWS/WjblruRzZO/lFK0pt+nAUcdbSf8od1Pu0OPuQzDqF6ymcu1r2RZWdRcWNgIAAABAIyJciiiswshduVBpCGQ/7vh/28IP03GtPZdsVU2u/25JxtVuaFxu2jfqki/J3lhbsi+hsgYnloodKeTCv4IqIfdxU4mYZrgDDlsgU1GQU0ElUEhljWe/kEq1yEFZ2Pvi+u+KmnJXFJSlPOceXHoZfRmp9znLh8emfW3hnWe/kM+nu2eVuw8UAAAAADQywqWIKrl4D160VrDMLBBwuEKrCppyB8bkD3Ms4/XfVt3WlDvqkjn/dltvKcke9FTy2ruP62/KXRhTtJDNFsiU61lVOq5lOdi+NOU29QyKMt6wkLIWTbn9xw005fbva7jrYLnjlhuP/zml6BVRtr5o5uP6X6Nor6/tfSnXs8ohXQIAAADQ4AiXIgq/eB9nC1ViMf9dx0KqVSL2azEFHO5j2y78bdU8trAmdEwV9Dey3VErtLdU1PFU0JS7MCZ38FS+4sTWY6cw3vLVPrVqyh1aoWWrYvMd0xNw+J7TE2j5n9PQlLvsc1oa2VvvbOcbVU8gNC1/3LAxeffzVy6NH9fflLswZtdxqwz2yoWm9PMGAAAA0OgIl6pU7cV72B21qg2tjBemUYOpSiqB/OP1nKu9ussWyPhfo4qackdckmg6T3uQ5t5W/s52YdUz9gCpuiVUpoAj+nNW0JS7gvfbVo1Wq6bc/n1tr31YzyrbkrnAZ8W3FM+mkiVz9uq4AsIlAAAAAI2OcCmiSu465q+msO1XyfZKesQUxlR+36iVN7awpt3flNsaRPjGU0l/I9/2ipopWyq0/DvbztUayFiWHLam4mpvsSxXrCTgKLNf4TnDKuDcQU/5JV/B59yHSrWIYWKgKXfYcSOOJ/Rz5HntKwlGg0Fj5L9rlVRSje3IsjgAAAAAja5pwqWXXnpJH/vYx7R06VK1tbXpwAMP1DXXXKORkZFJef7wO2ON/6e1eXbY8qAaNeV272tqym2tKnH9d+AiO+J+xu0Rg5ywptyVVETJMl73Vn9Tbr9q39NKKoH8nzHbax+2rNAeyFTQd8q1o78pt63Zt/85K2rKHfp+299T93HL7effNyzQ8uzXbghyXQfw96yyV3BZlsyVHQFs6j1XAAAaH3MFANRe+avpBvOHP/xB+Xxe3/zmN/XWt75Vv/3tb3XxxRdrcHBQ119//YQ//z5dvFv2czM15fZcmFoqZEz9Y4r7+ptyS9FDq8ruoBa9ysVWeRPWlLvL0nC60rt4FcdcvmeVU2ZM9uOW2xYIkCw9jKyVQIGeQN7fi9qU289+N73oTbn947Xeva6Cptz+562oQivked2sr4MxpCz8hn/Za+C4VVR+sSyuMvWeKwAAjY+5AgBqr2nCpdNPP12nn3566c8HHHCAnn76ad10002TPgmEXbxXcpHtD4hsIZCtCqayu9f5qpP8oZVlOVjU3lKSt+F0YV9LxYm1+sh9zBYl4qYQaGxMFYR+8myrrGdV1Pc0rNeQt3G5b0zWRuu+pVmujaaAwzNey/ttDcoqaMpdGJItpCy/jHRfmnJ7l0iW/+zOSCfVkqyur5fp71pxs22bVN3SS7KlyjTSXAEAaEzMFQBQe00TLpn09fVp9uzZ1t/JZDLKZDKlP/f39+/z84ZdvAcrl6ItZzJemLr+u9Km3MVjm4OTqNVJ5as/bIFWV2vS05Tbv28llUChd8yzBT0Re+WYK78s440cGJTfVnFT7oifMVMYYwt67Hczi74MclKacsv+dyZyUFZpzyrLvrbj+tmXmfr/rhW23fO7N3TfH7bp1EPmWI+N8uo1VwAAmgdzBQDsm6bpueT3/PPP62tf+5pWrVpl/b21a9equ7u79Fi0aFFVzxfWlNd68W6rPvIET6ZeLtEu7m3VSbYLXn9Tbv9xg31r3NssQVmn4e5qEZfiWSt2QkIgexPs8r1y/EsO3dvSybjaUpYQqKKKHfd4Kuy/ZQkpbfsVjm3fXu64UQM2f1Nu/86VfnZt26Lefa2mPatC//6P7WfaNrZzPCZPz6rAcS1/1+75/RuB4yKayZ4rAADNh7kCAPZd3cOlNWvWKBaLWR/r16/37PPaa6/p9NNP13nnnaePf/zj1uOvXr1afX19pceWLVuqGmfU5tmS1GO5bX1Fd21zMTblDh1TbGw85S+GbfuZt1vCLk94ZLp7nX1pUdnndFdZWQIiU1PusKDH3Zeq3LaejnL9mMzjtS+Zi7aEz7h9Hyp2ipv9Tbn9+1qXQVrPJWV4jdzHjV4R5Wa+C6Ltta+uSjCwpLOCvl7u5zUGmMXxtLd4elb5VXt3xemiWeYKAED9MFcAQP3UfVncpZdeqg996EPW31myZEnpv1977TWdeuqpOv7443XzzTeHHj+dTiudLh/2RBVWHeO+eJ/RGu3W8/7jWis4jA2nXWOquJqi/HNWW1Viu+h3P6epKXfUJXPmpW2F3wjrWWUK/YrbK13qVHUlkKcCxn5XvEqacodWdxXP09eUO7BvjZpyF/aNFqRV0pTbv3P1zcnLLzmc2W5Y9hry2baNN+rfNf/nM+pnbLpolrkCAFA/zBUAUD91D5d6e3vV29sb6XdfffVVnXrqqVq+fLluvfVWxeOTV3gVtiymyHjxHrU/T6UhUMSLz7BqCttxbUv8bBfoFTfljvoaVdiU27bEz83WiHlfjltJUOYOVcwBh+W47kDGsMSvNNaQ97uSJX62JuGB8daoKbd7u78pt39MlfSscrP9PZTC+ptZAsHQJX6WPlqES00zVwAA6oe5AgDqp+7hUlSvvfaaTjnlFO2///66/vrr9eabb5a2zZs3b8KfP6zfzfiyGHs/JutFtmVf88VlyMVn6YLXVE1RfslccccZNWzK7d4evqywunAkrHqm4obotv3Gjmtqyh21oXflAYdlOVjIvgtmtkmS3jZ3RvC4IX2/yj2nm+19qXVT7uLrYHvPTOON0jy73HHjsZgS8ZhyeUe9tgo4SwWh7T1NJWLBnlUuYY3CMa7ecwUAoPExVwBA7TVNuHT33Xfrueee03PPPaeFCxd6tjnO5N6s29bDKLQyJFDhEW0pWaVNuT1jslUuWY5rruaJuDzIVoUV1hDdEgoYK2SKx7VU7LSm4mpvCX7c53S16g9bB3RAb4fhuOPL7YLjLY61sp5ACgkxbNuihlamfY9cOFM/v/wkLZrVXva4YU25K72Ln73yK+L7bdh3SU+7EvGYjljYHdjm7c8VPdgLO5eWZFx/feahyuYddRuqqUrjNX4+xwJiS3XXbENfL0+lFT2XImukuQIA0JiYKwCg9pqm/vPCCy+U4zjGx2QIu3iX5UK6qCURV4elgsO2/Mq2bKtc35mWsYqa/Qx3bouy3K7iypCwBuPF3kiVNuUO6bFja8pdGk+Zbf/wgaP076uO17K3BIMKa+WS5bjxmP31Lap2KZ6pKXfYcSXpkHld6rBWx9ibcgcrjKJWhZV/v037hv1dW9zTod989j36xw8eHTyu6xi17FklSReesFQfP+kA47bi37ElPYaQMuS4pvEExmQJTuFV77kCAND4mCsAoPaapnKp3ryhSvlKINNF9tyuVrUk4jpk/oyKm3Kf9Lb99JNNr+u0Q+aUfU5TtYQkXXXGIXr+zT1665xO4/Zy47Uv4xn/3+4KmnJLKvWiMjZpHmNsyh02XlsYFlL51dOZVo8hfHM/caVVOclEXJ9+39s0PJpXb2f5Js3lAgVH1TTljlYRZVPLptySSue+tLf856/c85aOW2ab/3UtjWlsUN1t9p5V/r8ztvOM4sbzl+uVXUM6YL/guVo/g5bP0dBIrvTfVC4BAAAAaGSES1F5Ki3K9zAyXSTO7mjR/VeeYryLXFj1x8lv20+/Xv0e85BCKpf+z/KFxp9L4xfnxkoLywXv3K5Wze9u1ZKejoqackvSGcvmaf1LO63j2qfG5cbmz7aeVdHYArhyQcSlpx1k/Hl4UBaTHMd43EWz29WaiuuohTOtB640iLD3C4t2XNNn8IQDe/Xdj77DXBU2dmBjU273c1b5voX1rLJVLtl6fpUzr7tV87pbjdveeUCPfvPiTh23dHZwTKVqvuBz7hoaKf23v2cVAAAAADQSwqWIIjflLrOt2FA5eNxxlV5I70tw8qU/PULnH7dY7zyg/AWv6bitqYTuv/IUpQx31Ah7jQ6d36XvX/xO43i6xqqg3mJ4ncICuMU97frVs9JhC8qHGNW8RvY+WuFL8YzHDF06WGAKOObMaNVvPvteY+Nnd9Bna8ptUlyqaQpHqm3KLRUq1d79tv2Mz7mkp0Pzulq1Ysmsss8pVf93otKeVW61bp798ZMO0EdPWBqoNpOklkThZ6alq7sGRwI/AwAAAIBGRLgUke2W4e7t1VZalGvKHWVM1Tzn7I4WnXiQ+VatYUvJ/HdH8+8nVd4j5sS39urac4/QOw/oKfs76WRcbYbXaM37D9f/790HatHsYLPq4vW8belVOa2phPqHs5rbZahUG/vfSs8zaphY7rX3L0Us6kwndf5x+6slGdcMSz8mkz8+cr6278norCMXlP2dSptyh+lIJ/XwVacFqt+KR672uMXXZ+EsQ0g59r+mnlVR7+JXLVOwJEl/esxCbdm1V+cft39g286h0ZqPAwAAAAAmAuFSRN1tKb1lZptmtCaNlSNRmj+bhC1ti2JflnzV8rjF1yAZDzblDpNKxPWhdwQvsN16THfUUqHHkSlYkqT3HjpXDz23Q2cdOb+i8UjSF89eppd3DBr76BRVXrlU+N9yTbnHA8PKA44v/ekRFe8jSV2tqbLL+BaO3V3ukHkzAtv2pcJIUplgydvXy9aU2+S9h83R3/2fI3WSITiN2rPK1hOs1pb2dugfDI3JJemkt/bqv554TbMsd6cDAAAAgEZAuBRRKhHX/3z6ZCXiMWPAkRhbJlbpRXYxqDItBwtTDFRMF/77IsrdrWxmlQmBqh5PSCWVzXEH9Ohnl51U1fOevmxe2W3zZ7ZJL+/SQXPtzaoD+3W36YDeDh04p9NYzRIba+k9mQGHzdLeDj145anab0b5Ow5K+xaOBo479r+mptxh0smEPrBikXHb/rM71JqK6+hFwaV4+9KzaqL82fKFmt3RoqMWzaz3UAAAAADAinCpArZla6tOPkBPvdavQ+d3VXTMFUtm62/OPlwrlgR7H4X5qz86WB9YsVAHWiprqlE8T1OgYNM7o3h3sGCT8H1Ri6bctXbtuUfo/550gI5YGOzzZNOSjOuXV5xcdplUKhHTSK66ptITZf8ec1XYfp1pnXLwftqvMx1oyl0LtQ559puR1v9e/V51tNS2Z9VEScRjeu9hc+s9DAAAAAAIRbhUIx8/6YCq9kvEY/rI8Uuq2jeZiOutc2pbtSRJn/3jQ/XkK7t1ZIXBydvmztAPLn6nlvSaw4hqJccu/BspcOlIJysOlorKBUuStPqPD9Ub/cNll/k1klgspm9f9I6aH3dOV6GxuG05YrVMSxGlfetZBQAAAADTHeESAt572NyqKyaOP7B8Q+5q/dHh8/S/L+3UR45fXPNjN5o/f+fUP8cwRy3s1n+sOr7mFXlhqu1ZBQAAAADTHeESGt7+Pe265f+7ot7DwCSJxWJVLRMFAAAAANRH7RulAAAAAAAAYNogXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVK0pw6VMJqOjjz5asVhMjz/+eL2HAwBoQMwVAIAwzBUAUBtNGS791V/9lRYsWFDvYQAAGhhzBQAgDHMFANRG04VLP/vZz3T33Xfr+uuvr/dQAAANirkCABCGuQIAaidZ7wFU4o033tDFF1+sH/3oR2pvb4+0TyaTUSaTKf25v79/ooYHAGgAzBUAgDDMFQBQW01TueQ4ji688EKtWrVKK1asiLzf2rVr1d3dXXosWrRoAkcJAKgn5goAQBjmCgCovbqHS2vWrFEsFrM+1q9fr6997Wvq7+/X6tWrKzr+6tWr1dfXV3ps2bJlgs4EADBRmCsAAGGYKwCgfmKO4zj1HMD27du1fft26+8sWbJEH/rQh/Rf//VfisVipZ/ncjklEgmdf/75+s53vhPp+fr7+9Xd3a2+vj51dXXt09gBoFFNte865goAqL2p9l3HXAEAtRf1u67u4VJUmzdv9qxrfu211/RHf/RH+o//+A8dd9xxWrhwYaTjMAkAmA6m63cdcwUARDddv+uYKwAguqjfdU3T0Hv//ff3/Lmzs1OSdOCBB0aeAAAAUxtzBQAgDHMFANRe3XsuAQAAAAAAoHk1TeWS35IlS9QkK/oAAHXCXAEACMNcAQD7jsolAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFUjXAIAAAAAAEDVCJcAAAAAAABQNcIlAAAAAAAAVI1wCQAAAAAAAFVrunDpJz/5iY477ji1tbWpt7dX5557br2HBABoMMwVAIAwzBUAUDvJeg+gEv/5n/+piy++WF/+8pd12mmnyXEcbdq0qd7DAgA0EOYKAEAY5goAqK2mCZey2awuu+wyXXfddfrYxz5W+vnBBx9cx1EBABoJcwUAIAxzBQDUXtMsi3vsscf06quvKh6P6+1vf7vmz5+vM844Q0899ZR1v0wmo/7+fs8DADA1MVcAAMIwVwBA7TVNuPTCCy9IktasWaPPfe5z+u///m/NmjVLJ598snbu3Fl2v7Vr16q7u7v0WLRo0WQNGQAwyZgrAABhmCsAoPbqHi6tWbNGsVjM+li/fr3y+bwk6eqrr9af/dmfafny5br11lsVi8X07//+72WPv3r1avX19ZUeW7ZsmaxTAwDUCHMFACAMcwUA1E/dey5deuml+tCHPmT9nSVLlmhgYECSdNhhh5V+nk6ndcABB2jz5s1l902n00qn07UZLACgLpgrAABhmCsAoH7qHi719vaqt7c39PeWL1+udDqtp59+WieeeKIkaXR0VC+99JIWL1480cMEANQRcwUAIAxzBQDUT93Dpai6urq0atUqXXPNNVq0aJEWL16s6667TpJ03nnn1Xl0AIBGwFwBAAjDXAEAtdc04ZIkXXfddUomk/rIRz6ivXv36rjjjtO9996rWbNm1XtoAIAGwVwBAAjDXAEAtRVzHMep9yAmU39/v7q7u9XX16eurq56DwcAJgTfdfuG1w/AdMB33b7h9QMwHUT9rqv73eIAAAAAAADQvAiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFSNcAkAAAAAAABVI1wCAAAAAABA1QiXAAAAAAAAUDXCJQAAAAAAAFStqcKlZ555RmeffbZ6e3vV1dWlE044Qffdd1+9hwUAaCDMFQCAMMwVAFBbTRUunXnmmcpms7r33nu1YcMGHX300TrrrLO0devWeg8NANAgmCsAAGGYKwCgtpomXNq+fbuee+45XXXVVTryyCN10EEH6dprr9XQ0JCeeuqpeg8PANAAmCsAAGGYKwCg9pL1HkBUPT09OvTQQ/Xd735XxxxzjNLptL75zW9q7ty5Wr58edn9MpmMMplM6c99fX2SpP7+/gkfMwDUS/E7znGcOo9kcjFXAEB0zBXMFQAQJvJc4TSRV155xVm+fLkTi8WcRCLhLFiwwNm4caN1n2uuucaRxIMHDx7T8rFly5bJ+YJuIMwVPHjw4FHZg7mCuYIHDx48wh5hc0XMcer7TxVr1qzRF77wBevvrFu3TsuXL9c555yj0dFRXX311Wpra9O3vvUt3XXXXVq3bp3mz59v3Nf/Lwz5fF47d+5UT0+PYrFY5HH29/dr0aJF2rJli7q6uiLv10ymwzlK0+M8Ocepo9rzdBxHAwMDWrBggeLxplkBXRZzReOYDucoTY/z5BynDuaKAuaKxjEdzlGaHufJOU4dEz1X1D1c2r59u7Zv3279nSVLlujhhx/WypUrtWvXLs8LcdBBB+ljH/uYrrrqqgkdZ39/v7q7u9XX1zdlP3DT4Ryl6XGenOPUMV3OMwxzReOYDucoTY/z5BynjulynmGYKxrHdDhHaXqcJ+c4dUz0eda951Jvb696e3tDf29oaEiSAklZPB5XPp+fkLEBABoDcwUAIAxzBQDUT9PUvx5//PGaNWuWLrjgAj3xxBN65plndOWVV+rFF1/UmWeeWe/hAQAaAHMFACAMcwUA1F7ThEu9vb36+c9/rj179ui0007TihUr9NBDD+nHP/6xjjrqqAl//nQ6rWuuuUbpdHrCn6tepsM5StPjPDnHqWO6nGetMFdMvOlwjtL0OE/OceqYLudZK8wVE286nKM0Pc6Tc5w6Jvo8695zCQAAAAAAAM2raSqXAAAAAAAA0HgIlwAAAAAAAFA1wiUAAAAAAABUjXAJAAAAAAAAVSNciujGG2/U0qVL1draquXLl+tXv/pVvYdUtTVr1igWi3ke8+bNK213HEdr1qzRggUL1NbWplNOOUVPPfVUHUcc7sEHH9T73/9+LViwQLFYTD/60Y8826OcUyaT0Sc+8Qn19vaqo6NDf/Inf6JXXnllEs/CLuwcL7zwwsD7+s53vtPzO41+jmvXrtWxxx6rGTNmaM6cOTrnnHP09NNPe36n2d/LKOc4Fd7L6Yq5grmi3pgrCpr9vWSumNqYK5gr6o25oqDZ38tGmysIlyK44447dPnll+vqq6/Wxo0bddJJJ+mMM87Q5s2b6z20qh1++OF6/fXXS49NmzaVtv3d3/2dbrjhBn3961/XunXrNG/ePL3vfe/TwMBAHUdsNzg4qKOOOkpf//rXjdujnNPll1+uO++8U7fffrseeugh7dmzR2eddZZyudxknYZV2DlK0umnn+55X3/60596tjf6OT7wwAO65JJL9Oijj+qee+5RNpvVypUrNTg4WPqdZn8vo5yj1Pzv5XTEXMFc0QiYKwqa/b1krpi6mCuYKxoBc0VBs7+XDTdXOAj1jne8w1m1apXnZ4cccohz1VVX1WlE++aaa65xjjrqKOO2fD7vzJs3z7n22mtLPxseHna6u7udb3zjG5M0wn0jybnzzjtLf45yTrt373ZSqZRz++23l37n1VdfdeLxuPPzn/980sYelf8cHcdxLrjgAufss88uu0+znaPjOM62bdscSc4DDzzgOM7UfC/95+g4U/O9nA6YK5grGg1zxdR5L5krpg7mCuaKRsNcMXXey3rPFVQuhRgZGdGGDRu0cuVKz89XrlypRx55pE6j2nfPPvusFixYoKVLl+pDH/qQXnjhBUnSiy++qK1bt3rON51O6+STT27a841yThs2bNDo6KjndxYsWKBly5Y11Xnff//9mjNnjt72trfp4osv1rZt20rbmvEc+/r6JEmzZ8+WNDXfS/85Fk2193KqY65grmim855q3y/MFVPnvZzqmCuYK5rpvKfa9wtzxcS/l4RLIbZv365cLqe5c+d6fj537lxt3bq1TqPaN8cdd5y++93v6he/+IVuueUWbd26Ve9617u0Y8eO0jlNpfONck5bt25VS0uLZs2aVfZ3Gt0ZZ5yh733ve7r33nv193//91q3bp1OO+00ZTIZSc13jo7j6IorrtCJJ56oZcuWSZp676XpHKWp915OB8wVKv25Wc93qn2/lDPVvl+YK6bOezkdMFeo9OdmPd+p9v1SzlT7fmGumJz3MrnvpzE9xGIxz58dxwn8rFmcccYZpf8+4ogjdPzxx+vAAw/Ud77znVJzr6l0vkXVnFMznfcHP/jB0n8vW7ZMK1as0OLFi/WTn/xE5557btn9GvUcL730Uj355JN66KGHAtumyntZ7hyn2ns5nUyl707mioJm/X4pZ6p9vzBXFEyF93I6mUrfncwVBc36/VLOVPt+Ya4omOj3ksqlEL29vUokEoHUbtu2bYGUs1l1dHToiCOO0LPPPlu6u8NUOt8o5zRv3jyNjIxo165dZX+n2cyfP1+LFy/Ws88+K6m5zvETn/iE7rrrLt13331auHBh6edT6b0sd44mzfxeThfMFQXNfL5T6fulEs38/cJc4dXM7+V0wVxR0MznO5W+XyrRzN8vzBVeE/leEi6FaGlp0fLly3XPPfd4fn7PPffoXe96V51GVVuZTEa///3vNX/+fC1dulTz5s3znO/IyIgeeOCBpj3fKOe0fPlypVIpz++8/vrr+u1vf9u0571jxw5t2bJF8+fPl9Qc5+g4ji699FL98Ic/1L333qulS5d6tk+F9zLsHE2a8b2cbpgrmCua9byb8fuFucKsGd/L6Ya5grmiWc+7Gb9fmCvMJvS9rKj99zR1++23O6lUyvmXf/kX53e/+51z+eWXOx0dHc5LL71U76FV5dOf/rRz//33Oy+88ILz6KOPOmeddZYzY8aM0vlce+21Tnd3t/PDH/7Q2bRpk/PhD3/YmT9/vtPf31/nkZc3MDDgbNy40dm4caMjybnhhhucjRs3Oi+//LLjONHOadWqVc7ChQudX/7yl85jjz3mnHbaac5RRx3lZLPZep2Wh+0cBwYGnE9/+tPOI4884rz44ovOfffd5xx//PHOW97ylqY6x7/4i79wuru7nfvvv995/fXXS4+hoaHS7zT7exl2jlPlvZyOmCuYKxoBc0VBs7+XzBVTF3MFc0UjYK4oaPb3stHmCsKliP75n//ZWbx4sdPS0uIcc8wxntv7NZsPfvCDzvz5851UKuUsWLDAOffcc52nnnqqtD2fzzvXXHONM2/ePCedTjvvfve7nU2bNtVxxOHuu+8+R1LgccEFFziOE+2c9u7d61x66aXO7Nmznba2Nuess85yNm/eXIezMbOd49DQkLNy5Upnv/32c1KplLP//vs7F1xwQWD8jX6OpvOT5Nx6662l32n29zLsHKfKezldMVcwV9Qbc0VBs7+XzBVTG3MFc0W9MVcUNPt72WhzRWxsUAAAAAAAAEDF6LkEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RIAAAAAAACqRrgEAAAAAACAqhEuAQAAAAAAoGqESwAAAAAAAKga4RJQA2vWrNHRRx9d72EAABoYcwUAIAxzBZpVzHEcp96DABpZLBazbr/gggv09a9/XZlMRj09PZM0KgBAI2GuAAD8v+3cvUocexwG4HfNiRAxfuzGj0LsgizRxlyBXYqAwcYijVhtZ5UiFkGblBbWopWQvYmga1Kk8AICipIPkVyAERE9nQdPDgxnQXeV54GBmX8x8/tXL7wzu0VkBfeZcgkKHB8fX53X6/W8e/cuX79+vVp79OhRent7WzEaAG1CVgBQRFZwn/lZHBQYHh6+Onp7e1Mqlf5Y+/fnq3Nzc3n16lXev3+foaGh9PX1ZXl5Oefn53nz5k3K5XJGRkayvr5+7Vk/f/7M7Oxs+vv7U6lUMj09ncPDw9vdMAD/m6wAoIis4D5TLsEN+fjxY46OjtJoNLKyspKlpaW8fPky/f39+fLlS2q1Wmq1Wr5//54kOTk5ydTUVLq7u9NoNPLp06d0d3fnxYsXOTs7a/FuALgJsgKAIrKCu0C5BDekXC5ndXU1Y2NjmZ+fz9jYWE5OTrK4uJinT5/m7du36ezszOfPn5MkHz58SEdHR9bW1jIxMZFqtZqNjY18+/YtW1tbrd0MADdCVgBQRFZwF/zV6gHgvnr27Fk6Ov7pb4eGhjI+Pn51/eDBg1Qqlfz69StJsru7m729vTx+/PjafU5PT7O/v387QwNwq2QFAEVkBXeBcgluyMOHD69dl0ql/1y7uLhIklxcXOT58+fZ3Nz8414DAwM3NygALSMrACgiK7gLlEvQJiYnJ1Ov1zM4OJienp5WjwNAG5IVABSRFbSC/1yCNvH69es8efIk09PT2dnZycHBQba3t7OwsJAfP360ejwA2oCsAKCIrKAVlEvQJrq6utJoNDI6OpqZmZlUq9XMz8/n9+/f3jgAkERWAFBMVtAKpcvLy8tWDwEAAADA3eTLJQAAAACaplwCAAAAoGnKJQAAAACaplwCAAAAoGnKJQAAAACaplwCAAAAoGnKJQAAAACaplwCAAAAoGnKJQAAAACaplwCAAAAoGnKJQAAAACaplwCAAAAoGl/A4movCMkDfyBAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "%matplotlib inline\n", "from matplotlib.pylab import plt\n", "\n", "\n", "def plot_t_evol(out, labels):\n", " fig = plt.figure(figsize=(12, 6))\n", "\n", " max_abs = 8\n", "\n", " ncoord = out.shape[1]\n", " ncols = 3\n", " nrows = ncoord // ncols + (ncoord % ncols)\n", "\n", " for i in range(0, ncoord):\n", " ax = fig.add_subplot(nrows, ncols, i + 1)\n", " ax.plot(t_grid, out[:, i])\n", " ax.set_xlabel(\"Time\")\n", " ax.set_ylabel(labels[i])\n", " ax.set_ylim(-max_abs, max_abs)\n", "\n", " plt.tight_layout()\n", "\n", "\n", "plot_t_evol(out_del, [\"$L$\", \"$G$\", \"$H$\", \"$l$\", \"$g$\", \"$h$\"])" ] }, { "cell_type": "markdown", "id": "6f5585c7-2429-4f26-9cc8-8e8269e254d5", "metadata": {}, "source": [ "It can be confirmed by direct numerical comparison that the integration of the Stark problem using ``kepE()`` produced results matching those from the [previous tutorial](<./Comparing coordinate systems.ipynb>) 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.13" } }, "nbformat": 4, "nbformat_minor": 5 }