{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Map (series) inversion\n", "(by Dario Izzo)\n", "\n", "In this notebook we show how to use the map inversion algorithm to locally invert systems of non linear equations. " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Importing stuff" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from pyaudi import gdual_double as gdual\n", "from pyaudi import sin, cos\n", "from pyaudi import invert_map" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## The Double Pendulum end-point" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "Let us consider, in 2D, the end point position of a double pendulum.\n", "\n", "$$\n", "x = l_1 sin(\\theta_1) + l_2 sin(\\theta_2) \\\\\n", "y = -l_1 cos(\\theta_1) - l_2 cos(\\theta_2) \n", "$$\n", "\n", "We will develop a model for the inversion of this system of equations, so that $\\theta_1$ and $\\theta_2$ can be expressed as a function of $x$ and $y$\n" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# We show an image for clarity\n", "from IPython.display import Image\n", "from IPython.core.display import HTML \n", "Image(url= \"http://www.physicsandbox.com/assets/images/pendulum.jpg\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Define our variables" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "# Lengths for the pendulum arms\n", "L1 = 1.\n", "L2 = 1.3\n", "# Nominal values for thetas\n", "th1n = 0.3\n", "th2n = -0.2\n", "# Gduals (high order)\n", "th1 = gdual(th1n, \"\\\\theta_1\", 11)\n", "th2 = gdual(th2n, \"\\\\theta_2\", 11)\n", "# Equations\n", "x = L1*sin(th1)+L2*sin(th2)\n", "y = -L1*cos(th1)-L2*cos(th2)" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "\\[ 1.27409{d\\theta_2}+2.63265e-06{d\\theta_1}^{9}-0.000410445{d\\theta_1}^{6}+3.51104e-06{d\\theta_2}^{9}-0.159223{d\\theta_1}^{3}+0.0372501+0.000358709{d\\theta_2}^{6}-0.000252795{d\\theta_2}^{7}-0.212348{d\\theta_2}^{3}-2.39332e-08{d\\theta_1}^{11}+7.32937e-06{d\\theta_1}^{8}+0.00796114{d\\theta_1}^{5}-0.0107613{d\\theta_2}^{4}-0.14776{d\\theta_1}^{2}+0.0106174{d\\theta_2}^{5}-8.14374e-08{d\\theta_1}^{10}-6.40551e-06{d\\theta_2}^{8}-0.000189551{d\\theta_1}^{7}+0.0123133{d\\theta_1}^{4}+0.955336{d\\theta_1}+\\ldots+\\mathcal{O}\\left(12\\right) \\]" ], "text/plain": [ "1.27409*d\\theta_2+2.63265e-06*d\\theta_1**9-0.000410445*d\\theta_1**6+3.51104e-06*d\\theta_2**9-0.159223*d\\theta_1**3+0.0372501+0.000358709*d\\theta_2**6-0.000252795*d\\theta_2**7-0.212348*d\\theta_2**3-2.39332e-08*d\\theta_1**11+7.32937e-06*d\\theta_1**8+0.00796114*d\\theta_1**5-0.0107613*d\\theta_2**4-0.14776*d\\theta_1**2+0.0106174*d\\theta_2**5-8.14374e-08*d\\theta_1**10-6.40551e-06*d\\theta_2**8-0.000189551*d\\theta_1**7+0.0123133*d\\theta_1**4+0.955336*d\\theta_1+..." ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Lets visualize the Taylor polinomial expansion for x\n", "x" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "\\[ -0.25827{d\\theta_2}+8.14374e-07{d\\theta_1}^{9}+0.00132686{d\\theta_1}^{6}-7.11723e-07{d\\theta_2}^{9}-0.0492534{d\\theta_1}^{3}-2.22942+0.00176956{d\\theta_2}^{6}+5.12441e-05{d\\theta_2}^{7}+0.043045{d\\theta_2}^{3}-7.4034e-09{d\\theta_1}^{11}-2.36939e-05{d\\theta_1}^{8}+0.00246267{d\\theta_1}^{5}-0.0530869{d\\theta_2}^{4}+0.477668{d\\theta_1}^{2}-0.00215225{d\\theta_2}^{5}+2.63265e-07{d\\theta_1}^{10}-3.15994e-05{d\\theta_2}^{8}-5.8635e-05{d\\theta_1}^{7}-0.0398057{d\\theta_1}^{4}+0.29552{d\\theta_1}+\\ldots+\\mathcal{O}\\left(12\\right) \\]" ], "text/plain": [ "-0.25827*d\\theta_2+8.14374e-07*d\\theta_1**9+0.00132686*d\\theta_1**6-7.11723e-07*d\\theta_2**9-0.0492534*d\\theta_1**3-2.22942+0.00176956*d\\theta_2**6+5.12441e-05*d\\theta_2**7+0.043045*d\\theta_2**3-7.4034e-09*d\\theta_1**11-2.36939e-05*d\\theta_1**8+0.00246267*d\\theta_1**5-0.0530869*d\\theta_2**4+0.477668*d\\theta_1**2-0.00215225*d\\theta_2**5+2.63265e-07*d\\theta_1**10-3.15994e-05*d\\theta_2**8-5.8635e-05*d\\theta_1**7-0.0398057*d\\theta_1**4+0.29552*d\\theta_1+..." ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Lets visualize the Taylor polinomial expansion for y\n", "y" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Compute the inverse map" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "# And let us invert the relationship\n", "res = invert_map([x,y])\n", "th1_map = res[0]\n", "th2_map = res[1]" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "\\[ -33.3889{dp0}^{2}{dp1}^{2}+2.97638{dp0}^{2}{dp1}-1.53985e+08{dp0}^{2}{dp1}^{8}+7.49742{dp0}^{4}{dp1}-1.18854{dp0}^{8}-450.043{dp1}^{4}-0.461981{dp0}^{2}-60579.9{dp0}{dp1}^{6}+13499.7{dp0}^{7}{dp1}^{3}-19.2519{dp0}^{3}{dp1}^{2}+4.30881e+06{dp0}^{6}{dp1}^{5}+1.03732{dp0}^{3}{dp1}-2.89428e+08{dp0}^{3}{dp1}^{8}+4.30674{dp0}^{5}{dp1}-0.913559{dp0}^{9}+7.34785e+07{dp1}^{9}-5141.9{dp0}^{3}{dp1}^{4}-0.0793933{dp0}^{3}-45079.9{dp0}^{5}{dp1}^{4}-47904.3{dp1}^{6}+\\ldots+\\mathcal{O}\\left(12\\right) \\]" ], "text/plain": [ "-33.3889*dp0**2*dp1**2+2.97638*dp0**2*dp1-1.53985e+08*dp0**2*dp1**8+7.49742*dp0**4*dp1-1.18854*dp0**8-450.043*dp1**4-0.461981*dp0**2-60579.9*dp0*dp1**6+13499.7*dp0**7*dp1**3-19.2519*dp0**3*dp1**2+4.30881e+06*dp0**6*dp1**5+1.03732*dp0**3*dp1-2.89428e+08*dp0**3*dp1**8+4.30674*dp0**5*dp1-0.913559*dp0**9+7.34785e+07*dp1**9-5141.9*dp0**3*dp1**4-0.0793933*dp0**3-45079.9*dp0**5*dp1**4-47904.3*dp1**6+..." ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Lets visualize the Taylor polinomial expansion for th1\n", "th1_map" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "\\[ 25.7695{dp0}^{2}{dp1}^{2}-2.34715{dp0}^{2}{dp1}+1.18575e+08{dp0}^{2}{dp1}^{8}-5.78329{dp0}^{4}{dp1}+0.917051{dp0}^{8}+347.27{dp1}^{4}+0.343531{dp0}^{2}+46667.1{dp0}{dp1}^{6}-10399{dp0}^{7}{dp1}^{3}+14.7909{dp0}^{3}{dp1}^{2}-3.3184e+06{dp0}^{6}{dp1}^{5}-0.872053{dp0}^{3}{dp1}+2.22845e+08{dp0}^{3}{dp1}^{8}-3.30875{dp0}^{5}{dp1}+0.704389{dp0}^{9}-5.65823e+07{dp1}^{9}+3962.12{dp0}^{3}{dp1}^{4}+0.00876767{dp0}^{3}+34726{dp0}^{5}{dp1}^{4}+36914.5{dp1}^{6}+\\ldots+\\mathcal{O}\\left(12\\right) \\]" ], "text/plain": [ "25.7695*dp0**2*dp1**2-2.34715*dp0**2*dp1+1.18575e+08*dp0**2*dp1**8-5.78329*dp0**4*dp1+0.917051*dp0**8+347.27*dp1**4+0.343531*dp0**2+46667.1*dp0*dp1**6-10399*dp0**7*dp1**3+14.7909*dp0**3*dp1**2-3.3184e+06*dp0**6*dp1**5-0.872053*dp0**3*dp1+2.22845e+08*dp0**3*dp1**8-3.30875*dp0**5*dp1+0.704389*dp0**9-5.65823e+07*dp1**9+3962.12*dp0**3*dp1**4+0.00876767*dp0**3+34726*dp0**5*dp1**4+36914.5*dp1**6+..." ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Lets visualize the Taylor polinomial expansion for th2\n", "th2_map" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Compute $\\theta$ from $x$" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [], "source": [ "# First we extract the x,y position around the nominal thetas\n", "xn = x.constant_cf\n", "yn = y.constant_cf" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "x nominal is: 0.037250076627759976\n", "y nominal is: -2.22942304031922\n" ] } ], "source": [ "print(\"x nominal is: \", xn)\n", "print(\"y nominal is: \", yn)" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [], "source": [ "# Lets assume some desired (close to nominal) values for the end point\n", "xd = 0.04\n", "yd = -2.21\n", "# And compute the change with respect to the nominal position\n", "dx = xd - xn\n", "dy = yd - yn\n", "# We now compute the thetas\n", "th1d = th1n + th1_map.evaluate({\"dp0\": dx, \"dp1\": dy})\n", "th2d = th2n + th2_map.evaluate({\"dp0\": dx, \"dp1\": dy})\n", "# Let us check that indeed they are producing the desired end point position\n", "xdi = L1*sin(th1d)+L2*sin(th2d)\n", "ydi = -L1*cos(th1d)-L2*cos(th2d)" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Error in x: -9.838373171699999e-12\n", "Error in y: 1.7152856912616699e-10\n" ] } ], "source": [ "print(\"Error in x: \", xdi-xd)\n", "print(\"Error in y: \", ydi-yd)" ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "-2.2099999998284714" ] }, "execution_count": 13, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ydi" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## The Kepler's Equation\n", "\n", "Let's consider the long standing problem of inverting Kepler's equation:\n", "\n", "$$\n", "M = E - e sin(E)\n", "$$\n", "\n", "and face it using the map inversion functionality of pyaudi" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Define our variables" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [], "source": [ "# Nominal values (expansion points) for eccentric anomaly and eccentricity\n", "En = 0\n", "en = 0\n", "# gduals\n", "E = gdual(En, \"E\", 11)\n", "e = gdual(en, \"e\", 11)\n", "# The equation\n", "M = E - e*sin(E)" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "\\[ {dE}-0.00833333{dE}^{5}{de}+0.000198413{dE}^{7}{de}-{dE}{de}-2.75573e-06{dE}^{9}{de}+0.166667{dE}^{3}{de}+\\mathcal{O}\\left(12\\right) \\]" ], "text/plain": [ "dE-0.00833333*dE**5*de+0.000198413*dE**7*de-dE*de-2.75573e-06*dE**9*de+0.166667*dE**3*de" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Lets visualize the Taylor polinomial expansion for M (mean anomaly)\n", "M" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Invert the map" ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [], "source": [ "# We need to have the same symbol set in the variables\n", "e.extend_symbol_set([\"de\", \"dE\"])\n", "# We may now create the inverse map\n", "res = invert_map([M, e])\n", "E_map = res[0]\n", "e_map = res[1]" ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/latex": [ "\\[ -20{dp0}^{3}{dp1}^{8}+0.00833333{dp0}^{5}{dp1}-3.33333{dp0}^{3}{dp1}^{4}+2.8{dp0}^{5}{dp1}^{4}+0.133333{dp0}^{5}{dp1}^{2}+{dp0}{dp1}^{10}-1.66667{dp0}^{3}{dp1}^{3}+0.758333{dp0}^{5}{dp1}^{3}+0.000705467{dp0}^{9}{dp1}^{2}+{dp0}{dp1}^{8}-0.000198413{dp0}^{7}{dp1}+2.75573e-06{dp0}^{9}{dp1}-1.07937{dp0}^{7}{dp1}^{4}-0.0126984{dp0}^{7}{dp1}^{2}+{dp0}{dp1}^{6}+{dp0}{dp1}^{2}-14{dp0}^{3}{dp1}^{7}+8.05{dp0}^{5}{dp1}^{5}-0.162698{dp0}^{7}{dp1}^{3}-5.83333{dp0}^{3}{dp1}^{5}+\\ldots+\\mathcal{O}\\left(12\\right) \\]" ], "text/plain": [ "-20*dp0**3*dp1**8+0.00833333*dp0**5*dp1-3.33333*dp0**3*dp1**4+2.8*dp0**5*dp1**4+0.133333*dp0**5*dp1**2+dp0*dp1**10-1.66667*dp0**3*dp1**3+0.758333*dp0**5*dp1**3+0.000705467*dp0**9*dp1**2+dp0*dp1**8-0.000198413*dp0**7*dp1+2.75573e-06*dp0**9*dp1-1.07937*dp0**7*dp1**4-0.0126984*dp0**7*dp1**2+dp0*dp1**6+dp0*dp1**2-14*dp0**3*dp1**7+8.05*dp0**5*dp1**5-0.162698*dp0**7*dp1**3-5.83333*dp0**3*dp1**5+..." ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Lets visualize the Taylor polinomial expansion for E as a function of M (p0) and e (p1)\n", "E_map" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Compute $E$ from $M$" ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAArMAAAH+CAYAAACPyYY0AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDIuMS4yLCBo\ndHRwOi8vbWF0cGxvdGxpYi5vcmcvNQv5yAAAIABJREFUeJzsvXmMJPd95fkiI++so6uqq7ur2Xd3\nVfXJPtlsWoLgFaU/1jQkEyvLHq9kEQMZomRrB1ofohcSQJsrm9oZWzBMeGAdBijbGIxoiToMj0RL\nXo0woxVFtURCFLsr677vKyPvuPaP6oiKyIzIjIiMyIjM/H4AQc2srMzIrMiMFy/e730ZWZZBEARB\nEARBEK1IyO8NIAiCIAiCIAinkJglCIIgCIIgWhYSswRBEARBEETLQmKWIAiCIAiCaFlIzBIEQRAE\nQRAtC4lZgiAIgiAIomUhMUsQBEEQBEG0LCRmCYIgCIIgiJaFxCxBEARBEATRsoRt3p/GhREEQRAE\nQRBew1i9IzmzBEEQBEEQRMtCYpYgCIIgCIJoWUjMEgRBEARBEC0LiVmCIAiCIAiiZSExSxAEQRAE\nQbQsJGYJgiAIgiCIloXELEEQBEEQBNGykJglCIIgCIIgWhYSswRBEARBEETLQmKWIAiCIAiCaFlI\nzBIEQRAEQRAtC4lZgiAIgiAIomUhMUsQBEEQBEG0LCRmCYIgCIIgiJaFxCxBEARBEATRspCYJQiC\nIAiCIFoWErMEQRAEQRBEy0JiliAIgiAIgmhZSMwSBEEQBEEQLQuJWYIgCIIgCKJlITFLEARBEARB\ntCwkZgmCIAiCIIiWhcQsQRAEQRAE0bKQmCUIgiAIgiBaFhKzBEEQBEEQRMtCYpYgCIIgCIJoWUjM\nEgRBEARBEC0LiVmCIAiCIAiiZSExSxAEQRAEQbQsJGYJgiAIgiCIliXs9wYQBNF5yLIMSZIgiiJY\nlkUoFALDMH5vFkEQBNGCkJglCKKpSJIEnuchCALK5bIqYlmWBcuyCIfDCIVCYFmWBC5BEARRF0aW\nZTv3t3VngiAIBUmSIAgCRFEEsOfO8jyPUCgEWZbV/ynIsoxQKIRwOKwKXUXgksglCIJoeyx/0ZOY\nJQjCU2RZhiAIEAQBAFQxKoqiKmbNfk/5f0mSdAJW6+JSTIEgCKItITFLEIS/mIlYhXpittbjVrq4\nAFAoFBCPx5FIJNSoArm4BEEQLYvlL2/KzBIE4Sr1RGyjGD2eLMtYXl5GX1+feh/l/0OhECKRiOrg\nkotLEATRXpCYJQjCFWRZhiiKEAQBsiw31RVVnkvJ2Gq3SZIkFItFdZsAVGVxycUlCIJoXUjMEgTR\nEIqIXV1dRSgUQl9fn6XogNvCkWGYquiBmYsryzLK5XLVfStzuOTiEgRBBB8SswRBOEJxPQVBgCRJ\nyGaziEQi6O/v92V7jMSs2f2MBC4A1VnWQpVhBEEQwYbELEEQtlG6YpWWAeUSv80FpYFBm7HVoh3u\noHVyGYapiimQi0sQBOEPJGYJgrBMZVes1uVkGAaSJPm2bVadWbuPaebiCoIAnud1PyMXlyAIovmQ\nmCUIoi5WGgr8FrMAmuIM23VxFVFLgx8IgiC8gcQsQRCm2KnZCoVCVXnTZuK3OKzl4vI8rxvdS5Vh\nBEEQ7kFiliCIKpx0xXpxmd8Ofj+/Ecp7xrKs7naqDCMIgnAPErMEQag00hVrR0yWy2VMTk5iZ2cH\n3d3d6OrqQldXF1KpVJXws0oQxawZdivDgD2hG4/HycUlCIKogMQsQRCGItbumFkrmVlBEDA7O4uV\nlRWcOnUKx48fRy6XQzabxeLiInK5HGRZRjKZRFdXlyp0o9FoIy+vJagVU1hfXwfHcTh16pT6s8rF\nZsr/CIIgOg0SswTRwVR2xToRsQq1nFFJkrCwsID5+Xk89NBDeOyxxwDsZUn7+vrUMbTKffP5PLLZ\nLLa2tjA7Owue5xGLxVQHt6urC8lkUif+WsmZtYrWlVWGOgBUGUYQBKGFxCxBdChGXbGNiB4jMSnL\nMlZXVzE1NYXBwUE8+uij6rhZMxc3FAqpglX7OOVyGRzHIZvNYn19Hfl8XnffQqGAZDLpePtbCaoM\nIwiC2IfELEF0GFoRC1hb3GWFyqEJm5ubGB8fR3d3N27evIlYLFb1O9rFT7VgGAaxWAyxWAwHDx5U\nbxdFEdlsFtlsFhzHYXNzE0tLS0gkErosbiwWa3nxZmUBntH9qDKMIIh2h8QsQXQIThoK7KBkZjOZ\nDMbHx8GyLK5cuYJUKmV6/0afn2VZ9Pb2ore3V73t6NGjKBQK4DgOu7u7WFxcRKlUQiQSUcVtd3c3\nkslky2RMG4lPOKkMC4fDuiwuCVyCIIIMiVmCaHO8FrEK5XIZ6+vryGazGBkZ0QnMZsIwDJLJJJLJ\nJA4fPqzbPsXFnZ2dRT6fBwCkUimdixuJRHzZ7mZSrzKsVCqhWCxSZRhBEC0BiVmCaFOaKWInJyex\nubmJVCqFGzdu+CJy6rUpRKNR9Pf3o7+/X71NFEXk83lwHIf19XVMT09DEATE43Fdm0I8Hu8I4Wa3\nMkzJ4dJiM4Ig/ITELEG0GY10xdpBW7N1+vRpDA0NYWFhwTcx46TNgGVZdHd3o7u7W71NlmUUi0XV\nxV1eXkaxWFTv60Ynbithtv9o9zMtVBlGEESzITFLEG2CG12xVtDWbB07dgyPPfYYQqEQOI5ri2os\nhmGQSCSQSCQwODio3s7zvCpwO70TFzB3cakyjCCIZkNiliBaHDe7Yus9j1nNFuB/z6vXzx+JRCx3\n4kajUZ2LW9mJ265QZRhBEH5AYpYgWhQly8jzPN544w08/PDDnl3OVWq2enp6TGu2rEwA85pmi2mn\nnbjd3d1IpVK6k4F2xUpl2M7ODtbX13Hy5EmqDCMIwjbt/01KEG1I5cCDXC7nycE+k8kgnU4jHA7X\nrNkCguHMBgErnbgrKyvIZrMQRdFyJ25QXp9baAWqcgKgnIzVqgyrbFQgCIIgMUsQLYQSJxBFEUC1\nIHDr4J7P5zE+Po5yuWy5ZisIYjbImV2jTlxZli114vrteHtN5ULFWpVh2v2cKsMIggBIzBJES1Cv\nZksRco0exJWarZ2dHQwPD2NgYMDyY/otJv1+fidY7cTd3d2FLMvI5/Nt2Ylbb9+tVxmm/X2qDCOI\nzoPELEEEGKtdsaFQCJIkOc7MCoKAmZkZrK6u4vTp0zh//rztg3/lOFvCOZWduGtra8hmsxgcHGzL\nTlwnJ2JUGUYQhAKJWYIIIHa7YhUxaxezmi0n+O2M+v38XhMKhdq2E9fNiAxVhhFE50FiliAChNOu\nWLuuqCzLWFlZwfT0tGHNlhM6sc3Ab5x24ipOblA6cd0Us0ZQZRhBtDckZgkiACjOEc/zjgYe2BGS\nVmq2nOC3M0piY596nbjb29uYm5sLTCeu12LWCCuVYVoXlyrDCCK4kJglCB/RdsU2MvDASszATs2W\nE4IgZjvNmbVDkDtxg/R3q+fiaivDAOgWm1FlGEH4A4lZgvCJyq7YRlyeWmLWSc2WE0hMth5edeI6\n3ZagUtl3q0CVYQQRDEjMEkSTqdUV6xQjIVlZs6UVK17g94G63cV0M9/fWp242Wy2ZiduMpm0fXXB\nj5iBG9ipDFNaSQ4cOECLzQjCZUjMEkSTsFqz5QStM6ut2Tpz5oyjmq1WpJ3FbBBel7YT99ChQ+rt\nlZ24+XweAJBKpSx34raqmDXC7HPNcRxWVlYQj8d1t1NlGEE0DolZgvAYL0WsQigUgiAImJubc6Vm\nq1m0i4DpZCo7cYG9mEI+n7fcidtOYtYMWZZV4aq9jSrDCKJxSMwShEfY7Ypt5Hny+TzeeustHD16\n1JWarVaknZ3ZVkPpubXaiRsKhRCJRJDJZALfiesUo6EmVBlGEO7QeUc8gvAYp12xTlBqtiRJwrlz\n53D06FFPnqcVIDEbbGp14s7OzqJQKLREJ65TlIWe9aDKMIKwD4lZgnCJRrti7aDUbEUiETz88MNY\nWlpqSTerEy4vE7WJRCKqyH3ooYcABLsT1ymyLDf0fUCVYQRhDolZgmgQt7piraDUbPE8j+HhYXW1\nud0JYO0IObOtS6XQq9WJm81m1SyuH524TjGKGTQKVYYRxB7B+8QTRAvhZldsLUqlEiYnJ5HJZHDu\n3Lmqmq0gjJL1m3YXs+0sOKz83bSduAMDA+rtfnTiOsELMWuGncowhmF0Li4tNiNaERKzBOEASZKQ\ny+WwsbGBI0eOeCZilZqttbU1nD59GhcuXDB8HisTwIKG8p5R1KA+7SzSFZzuA83uxHWK3/u52XeU\nNuOvxcjFDXo7CtG5kJglCBtoa7bK5TLW1tY8WXQlSRLm5+exsLCA48eP486dOzUPJK0oZt2m3Z3Z\ndsZtoedlJ65TJEkKZPyhnotbeV9ycYkgErxPFkEEEKOu2HA47LqAlGUZKysrmJqawuHDhy3XbFHM\nYA8Ss61Js1xLJ524iourdOI6xWqbQRCotdjMyMWlyjDCb0jMEkQNag08cNsN3djYwMTEBHp7e3Hr\n1i3EYjHLvxsKhao6Kf3Cr8updPBsXfy8BG+lE3dlZUXtxNU6uHY6cZuZmfUCq5VhmUwGPM9jcHCQ\nKsOIpkFiliAMsDLwIBQKQRTFhp9rd3cX6XQa0WgUDz/8MJLJpO3HCErMwO9MIDmzrYnfedJKzDpx\nBUFQ2xTsduI2Ws0VVCq/G4vFonpiTZVhRLMgMUsQGuwMPGhUQObzeaTTaQiCgNHRUfT09Dh+rKBU\ncylxh1bsvCX8I2hi1oxwOIwDBw7gwIED6m1WO3FFUWxLMVuJ8vlXRK7VyrDKLC65uIQdSMwSBJx1\nxToVkNqareHhYV3NkFOCkpn10x0lZ7Z1aRUxa4TVTtzt7W1sbW2ht7c38J24jSCKoum0tlpZXJ7n\ndS4uLTYj7NBenyKCcECzumIFQcD09DTW19dr1mw5ISgxA7sCP5vNIhwOIx6PN/zc7S5m2/kg3spi\n1gijTtyxsTE1Rxr0TtxGsOtAK6+z8moOVYYRdiAxS3QskiRBEAQ19+qViLVbs+WEIMUMrGxHNpvF\n2NgYJElSXSztZdnu7m4kEomWPaC7TRD+tl7SbmLWCEmSEIlE0N3dHfhO3EZwK2ZElWGEHUjMEh1H\nrYYCt5/HSc2WE4IUM6i1HaVSCRMTE8hmsxgZGUFXV5d6/3K5DI7jwHEc1tbWUCgUdKvHlcuyZgf0\ndndm25lO+LuZVXNZ6cSdm5tDLpcD0JxO3EbwsrWBKsMIM0jMEh1Ds0QssFezNT4+jgMHDuCRRx4x\nzZC5RVBiBmaCUjvJ7MyZM7h48SIYhlFXPWsvy2pH9fI8rx7Q5+fnkcvlwDBM1QE9HA6TmG1xWsF1\nbAS7bQZGnbjK5EGvO3EbQRTFpi4AtVoZpqCIWqoMay9IzBJtTzNFrLZm6+rVq45qtpwQ1JiBLMtY\nXFzE7Owsjh07ZjtiEYlE0NfXh76+PvU2URTVA/rq6iomJychiiLi8TiKxSI2NzdN65GIYBKEfddr\n3HAsQ6FQ3U7c1dXVqqsadjtxGyEorQ21XFyqDGs/SMwSbYuVrli3niefz2N8fNyVmi0nBClmoBww\n1tfXMTExgf7+fty+fdu1y6Esy6Knp0f3HsuyDI7j8NZbb+nqkWKxmC6m4KdjRZjTKZlZL0SeF524\njRDkaj6ti2tWGaaFKsNaBxKzRNthpyvWDd566y1wHOdazZYTghIzCIVCyGazuHfvHqLRKK5du4ZE\nIuH58yrRg0gkgnPnzgHY2w9KpRI4jtNNcQqHwzqB2yoLa9qZThCzzX6NjXTiJpNJx9vaipPOqDKs\n9SExS7QNytm1IAiWu2KdotRs5XI5nDhxQs2A+kUQYgbFYhE7OzvI5XK4ePGibsW2HzAMg3g8jng8\nrnOseJ5XBe7s7Czy+TwYhtG5VV1dXYF1l9qRThCzQRB5Vjtx8/m87r52OnGbnZn1CiuVYfPz8zhw\n4AC6u7upMsxnSMwSbUGzumIra7b6+vowODjo+4HYz5iBtj83kUhgeHjYFyFrdQFYJBKpWlgjiqKa\nOVxeXkY2m1UvyWpd3KCtHG8XSMz6h1EnLqD/TNjpxA3q63QL7bGlWCyqopUqw/yFxCzR0jSrK1aW\nZSwvL2N6ehpHjhxRa7Y2NzcDc3m/2dshSRIWFhYwNzeHEydO4M6dO7h3715Tt6ESp+40y7Lo7e3V\niXDlkizHcdjc3MTMzIxu5biyEKdZBfftfADsBDELtNbf0OgzYaUTt7Ieq51RXGiqDPMfErNES9LM\nrtjNzU3Tmq1QKKQKaT9pppiVZVld3DU4OIg7d+6olx/9HmfrJmaXZIvFotqHu7S0pB7MK3O4bm6P\n3xESr+kUMdvqWOnELZfLuHv3LoDgd+I2iiAIppEKO5VhsizrqsLu37+Pa9euebvxbQaJWaKl8KNm\nKxaLmdZsBWXhVbNEpPKexONx3Lhxo2oMrd3tcPNv1yx3VFk5XnkwV3K4GxsbyOfzYFm2KnPYDllC\nL2h3sd7uaDtxV1dX8cgjj7REJ26jiKJoexCOlcqwp556Cj/72c9a9n3xAxKzREugiNilpSUwDIND\nhw559kHP5XKWa7aCJGa9pFAoIJ1Og+f5mu9JEBai+UE0GsXAwIAuc6hUI2WzWV01UiqV0i00aze3\nyintnLPsRFqhE7dR3FrspnVxzRxdojYkZolAU1mzxfO8Z5ckS6USJicnkclkLNdsBUXMegXP85ia\nmsLW1hbOnTunawUwIih9t0HArBrJyK1KJBK6HG40Gu2og1knngC1K7X+lkHrxG0UL45F9FlwBolZ\nIpCYdcWyLFtVbN0o2tX4Z86cwYULFyx/QbEsG4jMrNtIkoS5uTksLi7i5MmTGBkZsfSe0EjZ2pi5\nVdpFNQsLCyiXy4hGo+pBvFwut4RT5ZROyMx2wudCaZOxi1+duEEkn88jlUr5vRktB4lZIlDU64p1\nUzxqBdvx48dtj1oF2u+yuizLWF1dxdTUFA4dOqS2NljFiZjtBCFTC7NFNdqBD+vr6yiXy1hbW9Mt\nNEulUm1xeb7T94F2wc3pX83oxG0UL/bZTCbT9AmS7QCJWSIwWOmKdUPMmtVsOSEobQZusLOzg7Gx\nMXR1deHmzZuIxWK2H4OcWfdQuj8PHjyIaDQKSZJw5MgRVeDOz88jl8sBQNXAh2YcyN2k3cVsu78+\nBa87Zt3uxG0Er77nOI4jMeuA1vrGI9oSrYgFajcUsCzrOJMpyzI2NjYwMTFhWLPlhHbIzObzeaTT\naYiiiEuXLumcELu0m1MdNMLhMPr6+tDX16feJooicrmcuqBmcnISoiiqeUNtDjeotLvYa/fXp+DX\n9C+nnbiNjLP26rVyHKeLIRHWIDFL+IaTmi2nTqiVmi0nNCKu/aZcLmNychK7u7uWF7zVg5xZ76h1\ngtfT06Nzc2RZVgc+aPOGsVhMJ3CDUovU7mLPaZa01QjS9C8rnbhzc3Pq1Q27nbheiVmKGTiDxCzR\ndBrpirUbM1BqtkRRxPnz510/4w2FQuB53tXH9BpRFDE3N4elpSWcPn0a58+fd+1AS20G3mD3BIFh\nGKRSKd1CElmWUSqV1LzhysoKisUiwuFw1cCHZguSdhezsiwHRuR5iV/OrB20nbgKSsvIz07+LwCA\nru9+vm4nbq2BCY1AYtYZJGaJpuHGwAOrYrZUKmFiYgIcx7nmOhoRpJiB4oqavafarPDQ0BDu3Lnj\n+pexn0MTiNowDIN4PI54PI6DBw+qt/M8r+ZwZ2dnkc/nVTGsOLhdXV2BFylBJkiOpZeIotiSrzMU\nCqlCFgBu3LhRtxNXORa5LeApM+sMErOE5xjVbDkVMfXEbGXN1sWLFz0VTEETs2aribe2tpBOp9HT\n0+NKVrjWNlDMoLWIRCJVTpWSw+U4DsvLy8hms7reT0XkujnwoZ1PbDpFzLrZZtAsftB/S/ff79j6\nCYD6nbhra2soFov42c9+5monbiaTqdvnTVRDYpbwDLOu2EYwE4/amq0TJ044qtlyc3v8QNkW7cEk\nm80inU6DYRhcuXLF8/5Cihm0B0Y5XG3v5+bmJmZmZnSXYhUHNyg53CDR7jEKhVYT7WZCthZKJy7P\n8wiHwzhz5oyrnbjZbBZnz55t6HV1IiRmCdep1xXbCJXOrJs1W25sj59omwSUaWYcx2FkZES3+t3r\nbVBiJER7oe3yPHLkCAD9eFKO47C0tKRbMa44uO1WbG+XVhN5TmmFzCwAvHblbSgsNjZ8RxRF9Vjj\nZicuZWadQWKWcBUrXbGNoDyWFzVbTgiSM8swDHiex/z8PFZWVmxPM3NrGyhm0DmYXYrVHsQ3NjZQ\nKBQMD+KtIHzcoFPaDFohM1spZMM9ezLol2Z+ZOtx6gl3O524r776Kr797W/j0qVLWFpaQqFQsOXm\nv/TSS3j22Wdx7949/PjHP8atW3uO8+bmJt73vvfhtddew1NPPYUXXnjB8Pe3trbwG7/xG5iZmcGp\nU6fwla98pWkGiFuQmCVcwU5XbKOIooif/OQniMViuHbtGhKJhCfPY4WgiFllpfpPf/pTHD9+HI89\n9pgvBxWKGRCA8YpxQRDUHO7i4iJyuZyaNSyXy9je3rZUidSKdEqbgSRJgR7YoRWyiogF7AtZYO84\n5MRAMerEvXbtGt7xjnfg7t27uHv3Lj73uc/h2WefRX9/P65du4Zr167hsccew/DwsOFjXr58GV/7\n2tfwkY98RHd7PB7Hc889hzfffBNvvvmm6TY9//zzePzxx/HMM8/g+eefx/PPP4/Pfvaztl+bnwR3\nryNaAiVOoFxq91LE5nI5pNNplEolXLt2LRDF0kEQs5ubm+rQgytXrujmmzcbcma9o9WdvXA4XHUQ\nVyqRdnd3sb6+junpaQiCgEQiUTXwoZVffyfFDII6nOO1K29T/60Vsk4RBMHVvvLLly/j8uXL+Na3\nvoUvfelLOHbsGDY3N/HGG2/g9ddfx6uvvmoqZi9cuGB4eyqVwtvf/nZMTEzUfP5vfOMb+P73vw8A\n+NCHPoRf/uVfJjFLdAZu1GxZpVgsYnJyEtlsFsPDwygUCoEQsoC/QxM4jkM6nQbLsrh69SomJiZ8\nv2xLYtYb2nUBUSgUUlsRRkZGANSe3KStCmulHG6niNkgthloRSzPiVVC1okrC3g7AUw54RsYGMA7\n3/lOvPOd73T9ebSsrq5iaGgIADA0NIS1tTVPn88LSMwStmimiK1VsxWUg7sfzmyxWMTExATy+TxG\nRkZUJzYILjGJWaJRzCY3aQc+rK2tqZ2fWgc3lUoFUjQG5fvKa4Im2iuFbCVOhSzg3dCEQqFg6Pi+\n613vwsrKStXtn/nMZ/De977X9e1oNUjMEpZwsyu2HvVqtpRV+0E4ODgdr+sEQRAwMzODtbU1nD17\nFocOHdK9B0EQs9pGBSsE4W9ItAZGi2mUzk+O4zA/P284mrS7u9v3HGfQRJ5XBKnN4BfvehwAEO7a\n2x6tmI10N76NXr5Wo33lu9/9rifPBQCHDx/G8vIyhoaGsLy8rDuJbBVIzBI18aIrttZzaWu2zCZU\nKXVYQTg4NCNmIMsyFhYWMDc3h2PHjpl26AbBFaUFYEQzUTo/tTlxJYerOLhTU1MQRRGJREIXU4jF\nYk3bzk4Ss0F4nb941+PIrxZVIVtYLgPQi9hHfv4/G3oObTWXW/j1/f2e97wHL774Ip555hm8+OKL\nLen0kpglDFG6YpeWltDb24toNOqpiFVqtvr6+urWbCluaBBWPXsp3pT3ZXx8HAcPHsTt27drvuYg\nOLNBENRBce0Jf1ByuNpcvSzLulL7+fl5lMtlxGIx3cCHRCLhyb7TKdVcQcjMGglZN5zYSryIGSjf\nnXb3lZdffhkf//jHsb6+jieeeALXrl3Dd77zHQDAqVOnkMlkUC6X8fWvfx2vvPIKLl68iA9/+MN4\n+umncevWLTzzzDN4//vfjy996Us4ceIEXnrpJVdfVzMgMUvokGUZsiyrNVvLy8tIJpOeuRg7OzsY\nHx+3VbMVpEEFXh2gMpkMxsbGEIvFcP36dUvvC4lZgjCGYRikUimkUikcPnwYwH6pPcdx4DgOq6ur\nKBQKCIfDOoHrRg63k6q5/BKzSqygnONVIQsYRwoadWUBb15rsVhEPB63/XtPPvkknnzyScOfzczM\nGN7+xS9+Uf33wMAAvve979l+3iBBYpZQMRp4EA6HPRGOSs2WJEk4f/68rXaCIIlZtykUChgfH0ep\nVMLo6KitSTAkZvefvxNcMKIxtKX2Bw8eVG/neV7N4c7NzSGfz6tiWLvQzM4l5qD3r7qFXzEDRciy\nkf3nDsfDKG5UT/lyQ8gC3lwBymQygWnqaTXa/9NF1KVWVyzLsq6OJ62s2dKWqlvFzzosrxAEAVNT\nU9jY2MDw8DAOHjxo+4vSbyEJ2F8A5jZBeA+I1iYSiaCvr083AUkURTWHq0xtEkURyWRSl8M1i0d1\nUsyg2WJ26tefALAnZAs7DwYiGAjZaG8Y5V33jmVe/D05jqNRtg4hMdvBWKnZcssF5Xke09PT2NjY\nwNmzZ3U1W3ZpZoOA10iSpC7uOnnypOniLiuQM9vetKsYaoX9hWVZ9PT06ISGJEkoFArgOA6bm5uY\nnZ0Fz/OIxWI6gRuPxzsqZtDM1zn160+gsF2sErJaor37/+2WK+sVHMeRM+sQErMdiJ2u2EZjBvVq\ntpzQDjEDWZaxtraGyclJDA4O4s6dOw1fhgyKmLWzDUq1UldXlysHwXYV0+34mhRaNRYSCoXUHO6R\nI0cA7L2WYrGoxhSWl5dRLBYhCAK6urogSZI68KEdxW0z/5ZTv/4EytmyGi3QiljFldUKWTfx6vOY\nyWTImXUIidkOwklXrNOYgSzLWFpawszMDIaGhkxrtpzQ6mJ2d3cXY2NjSCaTuHHjhqPAvxFBEbNW\nvuiVurHZ2Vkkk0kUCgUwDFNVgO9kn2ln4deOtKqYNYJhGCQSCSQSCQwODqq3j42NIR6Po1QqYXNz\nU83haheadXV1+d4E0CooQlaBL+qPB0Yi9ur/+O+uPb+X079IzDqDxGwH0EhXLMuy4Hne1nPZqdly\nQtBiBoobWe89zefzGB8fB8/V2kneAAAgAElEQVTzuHDhguuXk4LQ8WpFzG5ubiKdTmNgYAC3b99W\n90lRFFVHa3FxEdlsFgB0ArfeAb9dRFEn0U5i1oxQKISenp6qHK6yvy8tLSGXy0GWZTWHq+z3Qagg\ntEoz/o5aIcvGwiju7mdjo6kwijvVi77cRhAETxb0kTPrHBKzbYzSFcvzvOOBB+FwGMVi0dJ9d3Z2\nkE6nkUgkLNdsOSFoC8AUR9TsveV5HpOTk9je3lYXd3m5HX5SawFYLpfD2NgYQqEQrl69imQyqe6f\nwN7ftbe3V51LDuzFVLSXbLPZLGRZNp3w1K4xg3amE8Ss0feD2f6ez+fVHO7MzAwEQUA8Hq8a+BDE\n98zrz97Urz+BUHjvZJaNVWRjawhZN11ZgJzZIEJitg2p7IptZGqXlZhBNpvF+Pg4JEnyxHE02qZS\nyfuzb6uYiUhtXvjUqVMYHR319ADkd5MAYOwO8zyPiYkJ7O7uYnR0VOdO1UNxtCoX3uTzeWQyGayv\nr6sTnpLJJIrFInZ2dtDX19dSjlYn4/c+2wysLowKhUJq5EBByeFyHIfd3V0sLi6iVCohEonoTuiS\nyaSvAtfrkxJFyBZ3CqqQLe6WEE0Zy5hIcu/zf/EV98fAeilmT5w44frjdgIkZtsMo67YRr5gauVT\ni8UiJiYmkMvlHNdsub1NflDpFMuyjNXVVUxNTeHw4cOu5oVrEbSYgSRJmJ+fx8LCAk6dOoXz589X\n7YtO9k2zA34+n8cvfvELbG1tYWFhAYIg6EaYdnd3ux55IdyhHRdDaWlE6GlzuIcOHVJv1w582NjY\nQD6fB8uyVTncZr23XjYZKEKWL/CmQlbryipC1iu8mP4FkDPbCCRm24RaXbGNYNRmUFmzdenSpaY6\nAkETs9oM7/b2NtLpNLq6unDz5s2mzn8PQsxAEbPr6+sYHx/H4OAgHn30Uc8L45VS+3g8jlOnTql1\nSEp10vb2Nubm5sDzvO6SrSJwg3jJtlMgZ9YZ0WgUAwMDGBgYUG9T2kEqc+dmsRw38UrMzvy796hC\nViGSiEASqo8BlSLWC1cW2HNmvXgPScw6h8Rsi2OnZssJ2piBKIqYm5vD0tJSw52ojRAE0aYlFAqp\nE81kWcalS5d0rmEzt8Pv90Upll9eXna1qcEO2vnmyWQSyWRSN8LU6JJtLBbTLTSLx+MkcJtEp2Zm\nvSAcDuPAgQM4cOCA7rmVz6U2lqNctVD2+0ZPvL249D7z794DAAjHwqqYjSQiKHH6dRxiWfLcjdU9\nH2VmAweJ2RbFaxGroIjZxcVFT2q2nG5TUJzZcrmMTCYDjuNw4cKFpkUtjPBz8VO5XMbExAQ4jkMs\nFsPDDz/sy3bU+wwYXbKVZVl3yXZlZQXFYlHNJCr/SyQSvoqudhV8nSBm/RyaEAqF1H1Yuz3KVYvd\n3V0sLCygXC4jGo3qFprZ2efdFuxrv/ubAPaEbGGnAGBPyGqJJCK6NgMtXrmygHcxg0wmo1sQSFiH\nxGyLodRsKSvAvRKxynNtb29jZ2cHPT09ntRsOSEIYlYURczOzmJ5eRnxeBxnzpzxVcgC/jizkiRh\ndnYWS0tLOHPmDM6fP48f/ehHTd0GLU4EPcMwiMViiMViuqYJrcBdW1tDoVBAOBzWCdxmLbpp50vx\nnSBmgzbO1uyqhXafX11dRaFQAMuyOgc3lUoZilY33cq13/1NlHNlVchqRaziylYK22YiiqInV504\njiMx6xASsy2Ck4EHjaCt2UokEhgdHfXsueziZzWXLMtYXl7G9PQ0jh49ijt37mBiYiIQYqOZYlY7\nwayZi9yaiVEmked5w0U3Vg72hDGdImaDvk+YndTxPK/mcOfn55HL5dSMunahmVuvURGyCpVCVvvf\nla5srHvPbDn71X9peDtq4VXMIJfLIZVKuf64nQCJ2YDjRlesHZSaLVmW1ZqtH/7wh549nxP8cmY3\nNzcxPj6O3t5enUsdlCEOzYoZZDIZjI2NIZFINH2RWz28fg8ikQj6+/t1LrwgCOA4DtlsVnewrxz2\nEHQx4xedImZb9TVGIhH09fVVDXxQcrirq6uYnJxEuVwGwzCYnZ1V9327V/K0QpaNsihxerFq5sYq\nIrZZeBUzkGW57UyBZkFiNqC42RVrBW3N1sjIiK0u0GbTbPGYzWbVsv8rV65UnTkHYeFVM7ajVCph\nfHwchUIBo6OjgV2o0GyXPBwOGx7sG51m1il0gpgF2qt+jGXZqv7n9fV1bG5uIhaL6dpDYrGY7sqF\n2eLKtd/9zf2BCBVCNhwLo5zTC9vibslQxHrtygLetBkE4epeK0NiNoC43RVbC79rtpzQLGe2VCph\nYmIC2Wy2psAPykQyr8SsKIqYmZnB6uoqzp49i0OHDgV2HwnKdhlNd9K6WUbTzJQDvtc1ZkGjU8Rs\nuyNJEhKJBI4cOaLeJssySqWSeuVCWVypZM8VgZv/5If3BiLsFsBG9Sd44Vj15yGSiECWqsVfM4Qs\n4F3MAAjOd1ir0VnfmgHHq65YI+zWbAXpgOP1dmjF25kzZ3Dx4sWazxkkZ9bNs3tZlrGysoKpqSk1\nHxx0dynI42yN3CxtbZJyuVaZZqZdaNbOBOm7hXCOKIpV3w8MwyAejyMej2NwcFC9XcmeZ7NZQyFb\n4ko6Eat1ZSOJiC5T6wdexAxKpVIgFli3KiRmA4BSszU+Po7Tp0973lCwtLRkq2ZLcR7b/ZKoLMtY\nXFzE7OwsHnroIcviLUiZWbdE9e7uLu7fv4/u7u7AtFhYJahi1giz2qR8Pg+O47C5uYmZmRnk83nE\nYjEUi8W2nGZGYrb1sXOMULLn7H/6Q/BhFnyhrHNkjYRsvfaCZrmygDcxA47j2v7E1UtIzPpIZVfs\nysoKzp4969lzra+vY3JyEv39/bh9+7bl2fVK12w7i9mNjQ2Mj4/bfm+APUGiVKX5iRsOcbFYRDqd\nRrlc9m34QyO0gyhSVoqnUin1ku3s7Ky62rzdppm10skHYY4kSbYFHp8rgo2GwRf2nNZwLFx30Vel\nKxtNNf+kzot2CqrlagwSsz5gNvDAq4PQ9vY2xsfHkUgkcO3aNSQSCVu/H4ReV6/gOA5jY2OIRCK4\nevUqksmk7cdoh5iBIAiYnp7G+vo6hoeHdZcEW4kgxwwaQRGyhw8fNp1mphTfKwtulExi0KeZUcyg\nPbCbI934D/872GgYxd08AGMhK1d8r2qFrFbELv/e/4WNn/1Ml8P1ugPa7cfOZDKBXVTbCpCYbSJW\numLdPOMzqtlyQjgcDqSYbeQgWCwWdSvzGzkjDsoCMCdCThs7OXbsWEvkYjsRo79rvWlmmUwGy8vL\ngZxmpoXEbHtg59hlJGS1hON7biyfr57uZeTE3r59WzfwQdsBrfTgKh3QQb3CmMlkKGbQACRmm4CR\niDX60CuisVEx4XbNlhIzCBKKC2n3IKh1IM+dO4fBwcGGD6RBcWbtsr29jbGxMfT29tqOVgSVdnVm\nrdIK08wqaXcx2yn7o1VnVhGygF7EKq6sXSF7/Mtf3/uZwZATQRCQzWaRzWaxuLiIXC6nNogoArer\nqysQ330cx5Ez2wAkZj1EGXggCIKlrthwOAxBEBx/sHiex9TUFDY3N10TakAwYwbKoiurwl+SJCwu\nLmJubg7Hjx931YEMygIwq+TzeaTTaUiSZNib28p0upg1w+40M63A9dqp90PMhr/x1/vPL5qcqJvt\nRw9OXMVf/wNLz9XuYl3Byvfx9h9+SBWyQlGffVVEbCXK7ULR/rqEcDiMAwcO4MCBA+pt2gaR9fV1\nTE1NQRRFJBIJXQ+0Wf7cK+OCYgaNQWLWI5x0xSpi1i6VNVsjIyOufnkGMWagCOx6wl9Z+DYxMYGD\nBw/i0UcfdX0Vaqs4s4IgYHJyEltbWxgZGdEJG7fx8wBOYtYataaZcRyH2dlZ5PN5MAyjyyK6Pc3M\nq30l8s9/o38eg+8wp0IWANiX/pPhXSpFbqeI2XptBrk//aj6b62QDcejKHFF3X0VV7aekFVcWTuY\nNYgUCgVks1ld/jwajer2/UQi4UmTAUDObKOQmPUAQRDU1e12FnbZFbPaKimlB9SLPFAQYwZW3OLd\n3V2k02nE43HcuHED8Xjcs20JspiVZRkLCwuYm5vDiRMnXD/ZqURxR60+h5vb0gmiwUusTjNTGhfc\nmGbmltiL/Le/VTa44cdqBJ3IDYXAv+fjHZFDr5WZzf3pR1Hm8ghFHkiOB2LWTMiaubRewTAMkskk\nksmkmj8HoBv4oMRzGIYBz/NYWlpSc7hu/H05jsPQ0FDDj9OpkJj1CCftBFbFbCM1W04IYsygloAs\nFAoYHx9HqVRqytjVIDuzm5ubSKfT6O/v98SVNsLPS/3tHDPwS6h7Pc3M6d9LFa8WMHJlbePwMx7/\n5l/jFgBM/b8Q/7f/s/HtCChmmdlKIVvmHiz6ilfnX6OpaFWDQaUrqwjdoc+/5Mp218Iof76zs4Pp\n6WmIooj5+XnkcjkA0C00czLJj5zZxiAx6wFOa7asiFmlZiuZTOL69eueuY2V21UsFuvfsYkYCWwl\nM7y1tYVz587h4MGDTREAQczM5nI5jI2NIRQKOa4cc0o7C0q/CNr7aWeamSJwlQN95Ym3LMuWnK3I\nv/5d/Q2z8Tm0HTGwi8FrYr/6lwAAhmUh/Np/cOd5AoJRZtZIyGpFrNaVjaaqXVotzXZrzVBc3OPH\nj6u3KSd32Wy2apJfZQ7XDBKzjUFi1gOcCqhaYjabzSKdTgMALl682NQy+yDGDLQCUpIkzM/PY2Fh\nwZPMsJVtCYozK8sy7t27h93dXYyOjjbcZOEEt8fq2oGEtH8YZRElSaqaZiYIAhKJhHrfWgNHIv/9\nH/b+UXY+vtQVV9YDwl//KwBoG1FbmZlVMrJqtADmbqwZQpE3FLHNcGXNMHKgjU7utJP8tINOYrGY\nKnCVTuhQKERitkFIzAaIcDiMUklfR+J2zZYTghozEARBPQs+dOhQ0y6jVxIEMasI+lwuh5MnT+L8\n+fO+XZa2KyjdXiBDYjY4hEIh9fKrgrLYRjnIb2xsQBAE7OzsqAJ36Bff3n8QIyEreDRxz8LCr0Zh\nKoRQu4ha7ec496cfRSgaRjmzFylgo2HwOb3rKkuSTshWurLRVMzjLXaG1Qoy7SQ/BVmWdTncf/qn\nf8IXvvAFNZbwve99DyzL4vz585aigy+99BKeffZZ3Lt3Dz/+8Y9x69YtAHvxsve973147bXX8NRT\nT+GFF14w/P1nn30WX/jCF9QhOX/2Z3+GX/mVX7HyNgQOErMe0Igzq+RvvKrZcrpdQROz5XIZ6XQa\nfX19ni7usoLfYnZ9fR3j4+MYHBxEV1cXjh496utCKIZhfHs/aAFY8NEutjl8+DCSySREUcSZ+f8P\n2HbuwLoSMfCR8Nf/Ckw4Av5XP+b3pjSEkZCthI1FIBSqe2SBfRFbzhn/3E9XFthb4O3UNGEYBvF4\nHPF4HIODg/jYxz6Gj33sY1hdXcVv//ZvI5fL4bOf/Szu37+PSCSCK1eu4AMf+ADe8Y53GD7e5cuX\n8bWvfQ0f+chHdLfH43E899xzePPNN/Hmm2/W3KZPfOIT+IM/sFYzF2RIzAYIlmVRLpcxPT3tWc2W\n0+0KSsxA6UhVVn6eO3fO703y7bJ6NpvF/fv3EY1Gcf36dSQSCWxtbXkyN9wOTieRubGfU8yg9Tg5\n8z/2/sEbCNkG4gWAjxEDB58/JrznxCmVYq0mahmGUYWsgiJkFVeWjT2o2qoQsiWuGFgnthK7Y3ut\ncPjwYZRKJfzxH/+x6sjm83n8/Oc/rzkV7MKFC4a3p1IpvP3tb8fExISr2xlkSMwGBFmWsbW1hZWV\nFZw5c8azmi0nBCFmUC6XMTU1he3tbYyMjCCfz3esaCmXy5iYmADHcRgdHdUVgvvtEgP2BKWyWNLN\nv2Wn7hetRuSHX3H2i0YRAze+n3yIGNQi8s9/01KC9toP/yuEB+JVLJYtC1k2FkVU0r/3Rq5srDuO\n/r/4e9e32y5K5tttJEnSOb7JZBKPPvqo689TyQsvvIAvf/nLuHXrFv7iL/7ClyijG5CY9QA7DpO2\n1L+npwd9fX04ffq0h1tnHz9jBpIkYXZ2FktLSzh16hRGR0fBMAxKpVJVvrjd0b4XZ86cwYULF6r2\nNT8XXyn4Xc1FBJvwa18HY+TCegijXTArijDcSypEq5TLerpNVmgVQcv/1R+oIlWsmOyliNhK2Nhe\nXtYsbqAQ6/YvQmaEF0MTan1fvutd78LKykrV7Z/5zGfw3ve+t6Hn/ehHP4pPf/rTYBgGn/70p/H7\nv//7+Lu/s9AaEkBIzHqElQP69vY20uk0UqkUbty4AYZh8POf/7xJW2gdP2IGsixjZWUFU1NTGBoa\nqnKqg+AWNwtZlrG2tobJyUkcPny4pmvvZ15VgdoMCCPCr9WY1mQ1YmC28Ctl0O4imYhWi4SMHhOA\nxGXq/KLziIERkX/+GyAWB//uf2/7cZuBVsiykbBOzEqC/jtavV/MvMGgnCsZCtgguLKANzEDBaOT\n8e9+97uePBewF29Q+J3f+R386q/+qmfP5TUkZn1AqdliGAaXLl1SV/qKohiYbKqWZosTReR3d3fj\nkUceMezmC2K3qxdkMhmMjY1ZnmLWajEDL2hXMduqrnOliHXFlU3pnVbLmN3Xxmcm1P2gPumBaJV2\ndyz9np2IgRGRf/27wAna8ud+H2Jp7+/JRsLqQAT2QQWXlC3o7l8pYitd2XAiBiYU7P1cEATXxawX\nj2mF5eVlderYyy+/jMuXLzd9G9yCxGwTKRaLGB8fR6FQwPDwcFU2JQhCxE9yuZzapXv58mVdpUkl\nQR8h2yilUkndV+xMMQtKzIDaDNzF77+pE9jXvw2GtzBsxUjcRqJ6Z1YRr1YruSSPT3Q17muo94Du\nR3Wd2wYImqANRcJVJwKKkOU1QjaSjKGcrTEQIbG3+EsRxpUExZUFvHFmM5lMzYVeZrz88sv4+Mc/\njvX1dTzxxBO4du0avvOd7wAATp06hUwmg3K5jK9//et45ZVXcPHiRXz4wx/G008/jVu3buGP/uiP\n8Prrr4NhGJw6dQp/+7fWp+oFDRKzHqF1p7STqc6ePWtas9WuB+J6KAuaMpkMRkZG0N/fX/d3ghgz\ncGNFviiKmJmZwerqKs6ePYtDhw7ZeswgxAz8zsy2ovBrN9jXv214uy1XtvIyv1fdsm7ChBDq2RO3\nUsaaY2uXoAha4YVPQtKIT5EXVCGrJZI0bilQXFlFyLYKXmRmM5mMo4EJTz75JJ588knDn83MzBje\n/sUvflH999//fXBOEhqFxKyHiKKIubk5dfFSEGq2goQoipidncXy8rLpgiYzgiZmFVfd6Rm7NiN8\n9OhR3Llzx1G9VhDcfb8FJYlZ/9CKWCuurBRPIaQVuIkHV2OsCtcmRwzsYFXU1srLAgBi1dEivwWt\nkZDVwmcLOhFr5MpWithWcGUBbyIBHMc5cmaJfUjMesTy8jImJiZUYWJn53d7IpIbKI6fG/2lsixj\naWkJMzMzOHr0KB577DHbjxu0mIGyPU6+5HZ3d3H//n10dXWZZoStEpSYgdVt2N3dxfj4OCKRCLq7\nu9HT04NUKuX4YBG0z03HcO8HYEt5y3eX4nuiNSSU9wWsgpGQtePKNjFiYOnuPQeAcBjy7rarm+GX\noBVe+CRCkYgqZkORiE7MVgrZSiLJGPh8/SaaoDq2XvR4O3VmiX1IzHpEMpnE7du3LY2k06I4jn6M\nZa2Fsl2Nfog3NzfVyV2NCLegOrN2KBaLSKfTKJfLuHjxoitn5kGIGVgR1KVSCel0GqVSCadPn4Yk\nSchms1hcXEQul4Msy+jq6kJPTw+6u7vR1dVlad8LwuvvOCwKWYYvqyLWFQL0+QcAMLX3T6Z3b42E\nm6I28j//K/i3/YZrj1cPRcgK2RxCD45tfH4/G8vGolVOt+LK1hK4WldWK2J7/u8vGt3dV5RubDfh\nOI7EbIMESzG1EQcOHHAktsLhcEPj8rxCEY92xblCNpvF2NgYWJbF1atXkUwmG9qeoLUZ2NkeQRAw\nMzODtbU1DA8P4+DBg659OQYlZmC2Ddqu3HPnzuHQoUMQBAGSJKG3txcPPfQQgL0ISi6XQyaTweLi\nIrLZLBiGQVdXF7q7u20JXMJD7v3A8GZtxEARsAyr/+4ICTaqtxqlyRGDWmhFbd2IQYAwErIKSkuB\nWKiOE1SKWDNXNqhObDMgMds4wVJMbYRTcaKI2aDhdHCCsio/l8tVTatqhFZ0ZrXximPHjjnOxdbb\njqDGDNbX1zE+Po4jR47UjZawLIuenh7dF7woishms+A4zlDgKvf1+/V3BBoRa+bKal3Yhuq47Ajc\neBLIc/v/nXxwtUN+8NnMNtg04OTzamJMML19QL3BDAZ52Uqa4c4qQhaATsjy+YJpZywbi9YdiAAA\nTIjpaCELUMzADUjMBgw/p23Vwu7gBK37eObMGVy6dMnVSzNBy0bWy/Bub29jbGwMvb29juInVgnC\nZfZKMZvL5XD//n1EIhFLXblmsCyL3t5e9Pb2qrdpBe78/Dx2dnYgyzJyuZwaUUilUm3h4AZmnzdx\nYxUYvuhunMCIuMmVHVHcF7BGdFUIBrMTH27X3vbUiRiY0jcAbG86+93k/nvspaDVClntoi9JFHVC\nVuvKmglcrSsb7YpX3aYliBEDSZI8+Rxms1mcOnXK9cftJEjMBgw/pm1ZwaoTKssyFhcXMTs765n7\nGETMnNl8Po90Og1JknDlypWa3blebkczUcSsIAiYmJjAzs4ORkdHPZn5XSlwV1ZWUCwW0dfXh0wm\ng/n5eeRyOTAMo8YT2kngNhtx4jXUWponxlMIWVi8ZydiIBuIU8bOIi/Z5PNgJmQZBuipuILkUc0W\ngMYErQYvBC3z959BKBJBKBKGkN134EOxKKR8oer+WhFr5soqIhYwF7JBxavpXxzH6U7SCfuQmPWI\nTowZrK+vY2JiAv39/Z66j0GkMjMrCAKmpqawubmJkZERDAwMNG07giBmNzY2MDk5iZMnT2J0dLRp\nrqKyOMPIweU4TnVwKwVuT08PkslkYAWu39EJceI1AADLV0x0ehAxEGs4sXYjBlJiX7wyoo2IgZdX\ntHoO7LuvnA1hW2vtQ1Rzab3vwfdDg6L2Jz/5CVKplDvNIH//GYj5wt5gBA2hWBRCpZCVpJojaoG9\nOIFWyNYiiK4sAM/Ws1BmtnFIzAaMoIrZWo5xJpNBOp1GNBrFtWvXkEgkmrx1/qOISFmWsbCwgLm5\nOZw4cQJ37txp6uVhvztet7e3MT8/j1Qq5csJjdnrZ1kWBw4c0GW2BUFQIwqzs7PI5/MIhUK6FoUg\nC9xmoQhZw59ViNiQlYlfRo+T6AJj8XvPlivrBd0a15bbcR4xMELr0lrIy1byWGkSmyP/q5orV5pB\nFIGrLJysJ8gqhayQzSNkIlbDiRiEnF7cal3ZcGrveFC5OMzIlY10JXSTw4KGl84sidnGIDEbMIIs\nZiudWWU8b7FYxMjIiC+XSYLSyRsKhbC7u6s6048++qgvjRR+tTwoNWM8z2NoaAjJZDLwznw4HDYV\nuJlMRidwKyMKQdjnmoFWyFa6srDwHhi5snIkCghliIn9CV9GQtYVV9ZJxMAO3QcAlgXc7JC1EjtI\nmjvhAz//b+h522+ozSCSJCGXy4HjOKyurmJiYgKSJOkEbnd3t/p9Ff7GX0NiWVMhq7iyyqKtSiGr\noIhYwLjlQEuka/++QXVlAYoZBBkSsx7RSMygVApejigcDoPn9w4uyiX0jY0NnDt3znQ8r9coLpzf\nwiKXy2FlZcW12rFGCIVC6t+pGUiShJmZGaysrKhVW3Nzcy07ztZM4CoRhdnZWeRyObAsWxVR8Hs/\ndBszR1aM7gkplte3GFhxZaV4CozI64RsS1DPfX1Qt2VZ1EbrrN7vGwDyOWuPVQftyZiCJEnI5/Pg\nOA7r6+uYmpqCKIr4peVXIeVykIr7f0sjR7Zm+0AopBOyRiiurFbEtgJeTP8CqM3ADUjMBowgO7OF\nQgFzc3OYn5/H8ePHfV/c5dYgB6fwPI/JyUns7OxgYGBAFTV+0qzMrCzLakZ6aGhIty/4HXVw+7nD\n4TD6+vp0C9i0And6ehr5fF4VuNqIQqsK3OLcW9D66oorqwhZK4jxboT5PYdR7Zo1cFutxguAAEQM\n6qEVtQ1emZH7D4HZWnP0u6E3vwfp8uPmP38Qp+nq6sLQ0NDe88kypP/yb6qQZWJRiNmK2jVZ0glZ\nrSurCNjKhV+G3bMmIjbxzF/XeFX+49VAI4oZNA6J2YARRDEry7K6cOb48eO+XUKvpNFBDk6RJAkL\nCwuYn5/HqVOnMDo6ioWFhUBUqjVDSGazWdy/fx+xWAw3b95ELKZ3aey6w61Y2WYkcHmeVwXuxsYG\n8vk8wuGw7lJuKwjc4txbiJT1rmCliDVzZcV4t+42pxVdLRExAPYiBkb09gE5zvhnFpAbqDYTew86\n+j32v3wWTHjv9TBG+ViT91PrwtbqlWUT+xlg7cSvVsKrmIEgCA2NMSdIzHpGu7QZ7O7uqpO7BgcH\nMTw87PcmqfgxOEEp/h8cHNSJ+mZf3jfDS2eW53lMTEwgk8nUHIDhpzPr53NHIhH09/ejv79fva2e\nwO3p6UEikQiMwC3OvVV1Gx9NIVKZl61AK2K9wjdX1uECL6l3r6EgtOu8ocDQna2Rl3UK+5X/CIRZ\niNmcKmQVV5Z94MSKFQ0GQq5QN04gFoo6EQuYC9mgu7KANzEDv5tK2gUSsx7i5MAaFDFbKBTUBT0X\nLlyALMuYnZ31e7N01BtU4CbKON5wOIzr169XNTawLBuIrLMXYlbb0HDq1CmcP3++pvhqt5hBI9QS\nuJlMBuvr6ygUCqrAVSIKRgLXa8GrCFmtK2skZCtdWSGaAivqBYpRfrbRiIFZzpYx+3sbfA5C+QYn\ngFlECu+7bFLvgF7Q1mHNKFIAACAASURBVMvLVuA0blAvaqDAfuU/AoBOyKo/MxGyoXi8SjxUurLt\nONVLFMWqK1FuEZQT2laFxGzA8FvM8jyPqakpbG1tYXh4GAcP7l2yyuVygbiMrqUZzmy5XMbExAQ4\njqvpRgah31XZDjfF3Pb2Nu7fv2+rocFvZzbomAncTCYDjuOwtramClxF3JbLZU/fUzNHthaCshBM\ndPeSsRSNA6iupWIkg8YDG0IWAKRkdS4xlNuF1N2HUNbmYAQbDl2VoLVJI/nZWihClolGwcQ0Jxui\nqApZLaEH0/u0C8QqUVsODCIHrezKAt7EDPxc99FOkJj1ECcHdT8unQN7OdC5uTksLi7i5MmTGBkZ\n0QmDII7Z9bKGSvt+nDlzBhcuXKgplPyqxKrErXG2xWIRY2NjEEURDz/8sK3JZU72e7daKfx2hZ0S\niUQwMDCgG65RLpfViML29ja2trYwPz+vc3Dj8XjD71ulkNWKWCNXVrCxCEzBLANr5LYyooFoNRCy\nbiGl9iqRpK69E1WdqHWxQ9aqoG0kL1tJLXdWK2TFnCYjXfE9priyoRpjqIVCyZETG0nuXeG6e/cu\nurq6Aj+hz4uhCdlsFl1dLdbuEUBIzAaMZjtLsixjdXUVU1NTOHz4MO7cuWN45hnEMbteCH/tKv1a\n70clQXJmG9kOURQxMzOD1dVVDA8PY3Bw0NE2tKKgDBrRaFQVuKIo4sCBA+ju7lYjCsro3kgkoopb\nuwJ3Z2XewAM1ho8mIVc8rpErW6+iS4w9aDYIaDOBImoBIJRzN5Yg9Q4gVMg6+l25/xCYYu26LquL\nvyL/8p8hQS9kQw8WIEmF6kiBlkpX1krkANh3ZRUBq+Xq1avqABPthL5KgevF4is7eOHMUi2XO5CY\n7WC2t7eRTqfR1dVluCpdi1+OcS3czsxyHIf79+8jHo/jxo0biNdwIioJkph1IiSVk5rJyUk89NBD\nDdWuueUOO33udhbSWoGrUC6X1YiCHYG7szKPeHl/xX05kkRUk4fVurJ8tLHKOUXA1sPIlTW9r82I\ngSk19nOxa8+1ZbO7+h/UEDTavKzhY/YeBLu7YX37NPAHjyGysWDrdyrd2ci//GdI2SwYzep5IyHL\nRKOQirXXAYTi8ZqRAy1GIhYAwv/H/wMAVf3OoiiqAtdsmll3d3dTBa4XYpZqudyBxKyHBDW/l8/n\nkU6nIYoiLl26ZOkSRxBfi1uX9kulEsbHx5HP53H+/HlHXyzNXIxWCydCUhHxiUQCt27daniBg9+Z\n2XYWs0ZEo1EcPHhQzbcDe/u0ElFQBG40GlUjChFRL0AqhayCVsSG6zQaAHpXVojtNRyExGpBZMeV\n9TJiYBVTUev08RoQtI0Q/tYLkAoFVciKuZwqZLUwJjVRinCtFznQPWcyAclBFRfLsujt7dVNxjKa\nZibLMpLJpOE0M7fxos2AxKw7kJgNIIogcTszVC6XMTk5id3dXQwPD+vcnVakUbdYFEXMzs5ieXlZ\nnV7lVLQHJTNrxyHWLm47f/68a+MUO1FQBo1YLIZYLGYocMPCnihVXNlyxNh1refGGkUMFAHrBDuu\nrCkunlBWRioAjah1GBXQPZaBoLWSl3XizgJ7QpZhWeCBUGXCEZ2Q1YpcAIaubL3Ige75TJxY3X0e\nuLJWsTPNrFLgutFH7sXQBIoZuAOJ2QCiNBq4VaIsiiLm5uawtLSE06dP161WahWc1mFpL6kfPXoU\njz32WMMnDkHJiVrZDlmWMT8/j/n5eZw+fbru4ja72BWzfj53JxGLxVDY1q+IV4Ss1pUtR5KICnoX\n1oory0dTYIXmlOGbRgzs4vBzz3cPIMJVL+aqFzGoxE2HtlZeNvTm98CwLOTy3t+HCUcgFfROvJkb\nCwCheKxu5AAA2GSiatGcE1fWDmbTzBSBu7m5iZmZGQiCgEQiga6uLjWCY/cYK8uy6yYTObPuQGLW\nQxodnNComJVlGcvLy5ienlZHjjZ6icStVedu4MSZ3d3dxf3799HV1eXKJXWFoDiz9WIGW1tbGBsb\nw8DAgGeT3PwU9iRmzVle24bWKzNyZM1c2qr7RbuQKGzpmg+MhGw7RAwqkcJ73xl8996VLSNRa4Yc\nrb48b1XQCl2aaXM23Vl2adpUyDLRWJWwVYRrKG7+/ah1ZdkHLqxYYwJYJXZdWTswDINUKoVUKoUj\nR44A2Dt2FQoFZLNZ7OzsYG5uDjzPIxaL6RzcWCzW1GMcObPuQGI2gLjRNbu1tYV0Oo3e3l488sgj\nrri8yiVsv1eUKtgRkMViEel0GuVyGRcvXtRdpnKDoGRmzWIGhUIBY2NjkGUZV69eRTLZ2IKeWpCg\n9IZGD7AJaf/SeEiu/tyoLq1g7sKWo3v5+rBYqttDawcrEQMhvn/AZ7C3f4ULGc3PuxB2aSiCUcTA\nDDOX1g5i70GESvXdbyewS9Pqv5mw/lI7U2OAg1bImrmyrIUoQVBgGAbJZBLJZBKHDh0CsCdwtfny\npaUllEolNV/upCHELplMBsePH/fksTsJErMe4sdI22w2i3Q6DYZhcOXKFVv9oPVQnNCgiFkrAlIU\nRUxPT2NtbQ3nzp3D4OCgJ19KQWoz0G6HKIqYmprCxsYGRkZGmpKTpjYD92n0NWld2TIbR1zYr3iK\n8vmajmyYL6gi1gyrrqwZQsxsulf1fqQIWQAQEhpHS5YgVAxFcEvc1kNxadkCV+eeLj6nBXdWFbLl\nsipkFRdWEbKVrixQ25EFACbEVAlZI1e2MmIQ7vLuJNoJDMMgHo8jHo/ragiNFlCGw2EUi0Wsrq6a\nTulzQjabdd1c6URIzAYQJ2K2VCphcnISHMdhZGQEfX199X/JwXYF4VK6Qq2YgTZicezYsYaqpqzg\np4Cr3A5ZliHLMlZWVjA1NYVjx47h0UcfbVoJud+Csh3FbCMsr22rrmylkC2ztZ21UjgJidGfvIZt\niNRKBBPRHDIadWsgZG0/X7JHl+EM5zQDERx8HpSIgRliotuxoC0dGEJsZ9nR7xqhFbJazNxYs+le\nWldWEbBW67gAEwH775+1/Pt+YLSAslAo4I033kChUNBN6evu7lb7cJPJpO3vWY7jXFt828mQmA0g\ndsSsUnK/srKCs2fPur6YR0vQBieYidnt7W2MjY25GrGoR1ByxAzDQBRFvPbaa0ilUk17/ZXbQONs\ng8Hy2rb67zKrz2uW2QSiov7StjZiUAonEWlAuCrwmoYDo7yskZD1CiG132Mqg0GkUF21ZSdiUIkU\njkHqjiHCVWdgjfKyjaK4s5WLv1QhG47oxKxc8X2puLL1xtRqXVij+xi5sqFIGKFI+0gMlmURj8dx\n6tQp9Tae59Uu3NnZWeRyObAsa2uaGS0Ac4f22dMCiJcxA1mWsbS0hJmZGTz00EOurMivR9AGJ1TG\nDJRcqCiKuHz5cseNCCyXyxgfH0exWMTVq1d9+4L0W8ySM6snIWVVIau4srUc2VLY3qVgo4iBGI5C\ntLmyX4uZK6uNGOgwc3HrjKPlE3uOmJGobQS++6ChoK2HmTurXfxVD3UxWTgC5B+48MoJrWYoglTI\n1+yLBQAmFLKVi2U1Lqxs1mIQcFfWDKOIXSQSQV9fn+5KqCAItqaZkTPrDsEbfkzUFbMbGxv40Y9+\nBI7jcPv2bZw6daopl5CDFjNQFoAJgoB0Oo3XX38dx44dw82bNztKyEqShNnZWbz22mvo7+9HMpn0\n9Uzf75oyErN7KK6skSMLoMqVlcHUFbK1IgblaJf6PyOCNL5Wht5o4BO9qrB1C77b2mjZveffd4xL\nB4ZsPY/WlWW3VoActydkFR4IWVk73SsWqz+mNmHtpEYslMB2Ja0J2RbG6sCEcDiMAwcO4Pjx47h4\n8SIeeeQRXL9+HUNDQ5BlGYuLi/jLv/xLPPLII/it3/otrK2tYWxsDLlc7VHFWl566SVcunQJoVAI\nP/nJT9Tb//Vf/xU3b97ElStXcPPmTfzbv/2b4e9vbW3h3e9+N4aHh/Hud78b29vbhvdrJUjMBhAz\nMctxHO7evYvFxUVcvXoV58+fd6UI2ipBixmEQiHk83m8+uqriMfjePTRR3UZp05gc3MTr776Ksrl\nMu7cuYOhoSHfL7X77cwS+0KWZfY/r3EhZ+rIlthq4WIlYlCK1Baw9WhmxMAKfKIX5ZTxIsl6eVmj\nn9sRtE6RIhpRWizsC9l8zlTIykXz5oRQIolQIgm5VB0nqBa8CZ2IrUuLurJAYwMTlGlmx44dw4UL\nF/CHf/iH+MEPfoCPfexjyOfz+OY3v4l3v/vduHnzJj74wQ/ic5/7HCYnJ00f7/Lly/ja176Gd7zj\nHbrbDx48iG9961v4+c9/jhdffBEf/OAHDX//+eefx+OPP47x8XE8/vjjeP755x29riBBMQMPcStm\nUCwWMTExgXw+j5GREd386mYSpJjB5uYm0uk0BEHAY4891lRRHwTy+TzGxsbAMAyuXbuGRCI4FTl2\nF8PR0ARvKMtRJDRi1o6QNUIMRVRnthRRKrqqHThWqr7NVrdskyIGRkjs3veIImijucYqt4A9QRsu\n2ZsY5mQxGLu7sS9kBX4/WvAApkantlQsWnZigT0R22m43eSTSCTw9re/HbIs4/Of/zwYhgHP87h3\n7x5++tOfYmdnx/R3L1y4YHj79evX1X9funQJxWIRpVKpqk/9G9/4Br7//e8DAD70oQ/hl3/5l/HZ\nz3628RflIyRmPcbJwVURs4IgYHp6Guvr6zh79mxD41bdIAgxg1wup4q4q1ev4vXXXw+UkPV6qIQg\nCJiamsLm5mbTqrbsQoLSX7SLvhTKcgxx7LtxUbGgE7GVkQMzV1YRsa1MZcSgFuXUgCuCthGs5GXZ\n3Q1AcVIj0T0x+wC5UNAJ2UpXloknEKoQ/maubKWIlcrVJy7tGDEA3BezQPXxIhKJ4OGHH8bDDz/c\n8GN/9atfxfXr1w0HA62urqrT0oaGhrC2tlZ1n1aDxGwACYVC4DgOr776Ko4fP+55rZRVnI6PdQOe\n5zE5OYmdnR2MjIygv7/fl+2ohSLivBCz2qqx48ePN7Vqyy5+xwzaVUhb2a/Gl8voYh+4sqG9leqV\nQhaw7saq9w8nEZb0kQCrrqwZQYsYmKG4tHbdVS35niEkM/acVqvubGRjQS9kC5rsZSwOxuQqCRPf\nE6ayhe90JpFAyMH3WkjbgvCbn7T9+0FCEARPJiYCxp/td73rXVhZWam6/TOf+Qze+9731ny8X/zi\nF/jkJz+JV155xfa2BGnKpx1IzAYIWZaxsbGhXj7/pV/6pUC5jn7EDCRJwsLCAubn53Hy5EmMjo4G\n9oOmtCu4LTKVEbzd3d2+VG3ZhcbZuo/V19TF5lCW9/ePslztylQ6k5WurLaSS1kQVilk7WAWMSjH\nqhcpMjCJDGhef6ykGYTgQcTAjFKyH7H8lvHv1snTAuaCVrv4q5J6rmxkYwFgHxzGIxXfC7E4UDEQ\nQXFlFSFrhNaVZR44sbJBHVctVzbUQpPBrNJIZtaMWp/r7373u44ec2FhAU8++SS+/OUv4+zZs4b3\nOXz4MJaXlzE0NITl5WV1Ihqw9x26srKCe/fuIZfLIZFIYHBwEEeOHMHAwEBghiZVQmLWY6weXDOZ\nDMbGxhCLxXD9+nW88cYbgRKyQPNjBuvr6xgfH8fg4CAeffRRz86K3cLtKWClUgnj4+MoFAq2RvB6\n6RDbeX6/aEcxa4Xx5TKioT1BkwjlVSGrdWVLiCMKcydOEbJ267lqYbQ4LNRAs0HJQAQDQKxonjF0\ni1qC1gvEcBysYNz9qgrZfHZfyCqubMy4cqtSxJq5soyDTGwokYBsInRa3ZUF9pxZo0v2jZDL5Vxt\n3tnZ2cETTzyBP//zP8fb3vY20/u95z3vwYsvvohnnnkGL774our0Liws4B//8R9RKpWwsLCA2dlZ\nrK+vg2VZDA8P44knnsCv/dqvuTpZ1C2CrQ46gEKhgPHxcZRKJYyOjqqVSkE8IDerzSCbzWJsbAzh\ncBjXr1+vu7gpKJdFlKqwRpEkCXNzc1hcXMTZs2dx+PBhW69PWYDl5xk0tRk0lz0hu++eVjqyJdTu\nE60lXo1cWaOIAR+OgS2XUYzqT7pYqTkNKKW4xuFkGMSKBgMRbORlTZ/HpqAta4ZGOIkbGKEKWb6s\nF7JaEat1ZWPx/ShCDZhYDAjp36N6rqySo5UNnNp2wgtndnd319Eo25dffhkf//jHsb6+jieeeALX\nrl3Dd77zHbzwwguYmJjAc889h+eeew4A8Morr+DQoUP48Ic/jKeffhq3bt3CM888g/e///340pe+\nhBMnTuCll14CsLewOJvN4vbt2/jIRz6Cw4cPq8959+5dfPnLX8Y3v/lN/Mmf/AnOnz/vzpvgEiRm\nfUJZyLOxsYHh4WEcPHgw8Adir2MG5XIZExMT4DgOo6Ojllob/HYhtbjhzCpu9OHDh3Hnzh1HglTZ\nDr/ErJ2/Bc/zGB8fx87ODnp6etDd3Y2enh4kk0lHf1O/XeEgwDL7n9E4CjohW+nKykz9blkrFCMp\nsLLQkJC1EjGw9kB7+00pvtcbayRqG6XZDq0WNp8B2DDkcAQMrxGQJm6s2e1aV1ZZICbz1gQptRm4\ng9PpX08++SSefPLJqts/9alP4VOf+pTh73zxi19U/z0wMIDvfe97VfcZGRlRRXAlN2/exM2bNyFJ\nErJZ5/lxryAx6zGVB2RJkjA/P4+FhQWcOHEiMIu7rOBVzEDrRJ4+fdrWSF5FYAfhPWxEzCotDSzL\nWnKj621H0AWddoLdiRMncPLkSWSzWWQyGUxPTyOfzyMcDusEbjweD8RJS9DQurIlKYokq48VmGEU\nObA6vrYYcX6ZsZGIgRMUUQsA0WLG8D718rIiW51TVwStlbysFq07WysvCwBb3ccBAIPb4+ptoVJ+\nT8hqF3pVvqeFvF7Emriyteq6jAilUqYjb81oh4gBYH1ogh04jnPkzHqFYgzdvXsXd+/exdDQEJLJ\nJAYGBtDX14eenh709vYGcvwuidkmIcsy1tbWMDk5iUOHDlnKgAbFcVRwO2YgyzLW19cxMTHh2IlU\nxGwQ8sWV43WtIAgCJicnsbW1hfPnz+vGIjrFbs9rs8lkMrh37x56enpw+/ZtsCwLnuerxkKWy2Vw\nHIdMJoPV1VUUCgXEYjGdwK3MsHWiM1uWwoiG+GohK0URD+2LU0W4KgI3ytQXrpURg1I4CdYgdsDK\nwRmmYobEsCgm9vaveMGdiUelZD8iZeuTmxScxg0i2yt6IatEDEoVQxDMXFrNzyuPLEaurBIxCNXJ\nSLZ7xADwxpnNZDKBEobKFb0f/vCH+MQnPoFz584hGo1iaWkJ5XIZDz30EH7v934PTz/9dOAWgpGY\n9RiGYbCzs4N0Oo1kMokbN24gXmceNrDfNRsEkabgZsyA4zjcv38f8Xjc8nvi9TY1ip3MbKUzOTIy\n4tqJi9sL0dxCiZFks1ndgjazbY1GoxgYGNB16ZZKJWQyGWQyGSwuLqJUKiGRSKgCN5VKdZSYHV8u\noytcLS4rhax6ew2ntpYr6ySK4EvEwCLFRJ9rgjbbdRhd2VXDn2nzspXUc2UrUYTs/g0GQjYaB4r6\nBgMdisi1kJ8F6ovYerSLKwt4k5nlOA69ve6OUW4E5Rg0MDCAz3/+83jqqacQCoWwtLSEf/iHf8DI\nyAi++tWvore3Fx/4wAd83lo9JGY9ZmlpCfPz87hw4YKtywlBFLNuiCRlhX4+n8fo6GjDH+QgCTer\n26JUbSnOpNt/4yC9J8CecFdWxp45c8ZWjKSSWCyGwcFBDA4Oqo9dLBaRyWSwvb2N2dlZZLNZvPnm\nm6rA7e7uDnwThhWM3jOtK6ug/bdCFKW6i8CMqBSxRq6sHTyPGNjYr9x2ae3CRxKI/P/svXmMJXd9\n9vvUcvale3p69n3tnvFgj+3xzDgmQSFxIJiIa4QSroBwAwng3DeRCCigK65khMCygPyBCAKFvLoY\nsiACThSIAFvhjcQb7LENjuNx7zPT0z09S09vZ6u96v5x+lddVedXdarqVJ1TPXMeqaWZ7rPUqVPL\np556vt+v4j5W1ipOqpsgywj11jZcQBNkaZLEtk6t1ZVlsuvfudJ6cUOLGDhdWSsAJ+co1Lniyswm\nLWYAAD/84Q/xyCOPmPG93bt347//+79x7tw57N+/H4qSvB7Rm/8In3Dt2rXL1sPNr5wjbZOgTpxD\nTdMwOzuLGzdu4PDhw4Er9N2UNGfWCyIlScLk5CQkScI999wTaUsWq5J0q31tbQ1jY2MYHBzE+fPn\nI4dKhmGQy+WQy+WwY8cOGIaBl156CYcOHUKlUsHi4iIuXboETdNQLBZNwC0Wi4m7TeYl2vc5dV0G\nwJvwmucE899WV1bUMzAcFertIgYiW0AK/m4db5aIgZtqxWbFdl6gF3TR8rJE0nrbMS93tlPtvPEq\nGLkJkQafBpOygIQk2CGW5so6QdYtP5vtrBDQ6eKqv/c/Onq9pMkwjMhrMyqViq1jQK9FPt9jjz2G\n5557Dt/85jcxMjKCq1ev4saNG9i5cyfW1tYS5SYT9WE2ZoVt15REmA0jwzBw8+ZNzMzMYPfu3ZEX\nvCUJZt0ys7quY3Z2FgsLCzh69GjsY4mT4Mzquo7XX38dgiDg1KlTnuAex7ooFAooFArmyEZd11Gv\n11GpVLCwsIBarQaGYVAsFs38baFQSEQhoV/JOo8Uu7Ht0xxZUfdX4KMZPDQ2utNBkIiBmCogp1Rb\n/9DlC7JGbsgVaP0oKNCuZbdjWJl1/ftyaV8LyLKiI5/r5sYCQDbvHTkg4ngwnOO79+nKsvkCjATd\nPdxsqtVqOHr0aK8XwxQ5Fr/vfe9DPp/Hd7/7XXzrW9/Ctm3b8NRTT2FkZAQf+tCHcPDgwd4uKEV9\nmE2outXTNU6tra1hYmIChUIBZ86cibzhNJAsmKVB5K1btzA9PY2dO3fi4Ycf7gos9RJmDcPA3Nwc\nGo0Gjhw5gp07d3a1iNHtvViWNSMHe/bsAdC8W0A6KMzNzaFer9se10mLsLh18ZoOSUuZMGttx0VE\nQDbD2h1WP4VfNFc2bMRASDVvo7IGfT9lDN18jFM5md59gP5CnX9PUQKtV16W6HbpAIardKDdeeNV\n6FwKHEQYvP1CRU9nwTqLvgi4ejmsVlc2v+6khhhRzq4/128rr77oqlQqvlpQ9kIPPPAADh8+jEKh\ngG3btplu7AMPPNDjJaOrD7MJVbenbQVRuy4Loiiat9ODZoWDKkwHgbhkdeHJ4IdUKtVRgVsY9Spm\nsLKygvHxcQwNDdlc0aSK4zgMDAzYbpmpqmp2UHC2CCMRhV63CHv5Co9cSkYx1QQTSUshz1t6y7KS\nb0dW0rPgmfAXzdaIgchv3GZm3cbMBpSQbq30DgS4baRSIgRWoPWKGMQpRlebICtUTZAlrqzezo0l\ncnNl8+GKutgAz7sTIwZxKGyf2bj1wx/+EM888wxmZ2dRr9cxODiIz372s/it3/qtXi+aq/owG7PC\nnvSSGjMgTigt+6hpGi5fvoxbt27h6NGj2LZtW+wn/aimbkUhlmUhyzLGx8exurqKkZGRSFpthVmO\nbgK+NQv8pje9CcViEUtLS4FeIylt6Hie92wRduPGDYiiiHQ67dkiLE7lUrLpyjZB1n771wqyTldW\n1DNIcxIkvfOLK5EtgKOAcBCQZQJCr8GwaGQ2nKy81H6ErVde1k2N3BAycrjG8LXiDqQVH7f310Vz\nZzm54QmyLa4s4O3IAgDHtYIszZV1RgyyeUBufdzd5Mrquh7LXbWkwSw5Dn/2s5/F008/jbe+9a0A\nmtO/nnjiCfzkJz/pyTnNj/owm1DxPJ/MisF1x9gKs4Zh4Pr167h8+TL27NnT1UEQSYkZGIaB1dVV\nLC4u4tixYxgZGekZnHULZq0DQDrNAicBZN3kp0WYLMvIZrM2wI2jE8nLV3iwrGG6slaJWhpZjg4Y\nZuSAkduCLC1ioCANhbW7lDSQdZNbxKBTEbA1wKAgtwfbIJLSRU+gJcVfUYtfH+xgjRa4ubHG+u9J\nrtaU1ZXNrQOs8zFe6qAY7E5zZYF4OhkAyW3NRYwJogcffBCSJHX1DmNQ9WE2ZnXizAqCv7Yt3RTJ\n8hInamVlBRMTExgYGMBDDz2EdLq7t+U4joMUIvMVpcg6SKfT2LNnD/bu3dvT5enGBLDl5WVMTExg\neHjYddhFUtzWOOTVImx5eRmzs7NQVRX5fD6yFmEMw5iuLACbKytq7vudV+TAT8RAMPK+owhRxQvC\nqJ5ugm2UUFvJbUdZuBX4eXIq7+nOrmXpHW4IyLJWd9TRzoy4sq4gS5SzuLC0x9COmxwHcA6Qpbiy\nd5tUVY2lvV+lUkkUzBqGAcMw8Pjjj+PJJ5/EW9/6VgwODuKFF17AsWPHOppMGbf6MJtQJT1mIAgC\nJiYmoGla22r1bixPL0SywYqi4NSpUxAEASsrvelZaVWcE8Csn/nee+9FwaWpOsnt3qkw65SzRRjQ\nPDE0Gg1bizBd11EoFEK1CHv99nbkUiqKKdEEWmADZJ2urAHGHjlg/N0WFozOWjT5VdCIgV/5hVpa\nXpamMEBb5bZga8Cowc7F/wZgB1k9nW3pYGB4ZWaBpiubC5iLNYvB/Lm3d1PEAIjPmRVFMVFuJ8Mw\nZszgqaeewte//nXU63U8/PDD+M53vtPrxfNUH2a7oDAFOUmFWZZlceXKFdTrdRw7dgzDw8M9X55u\nF4BpmoYrV67g5s2b5u11oHlrJgnFaHGsE2t7sWPHjrXtnZykXre9EsMwgVqEWaeY0WI6uVTzeEAD\nWadELY0MFzympBqtpwSaK9uLiIHBBIsu1dODUJkUSnK47gSSxaEMA7RL2T3YKl7z9djh6uxGrGAd\nZmkgyzj2Kacrq2dyYP1cQBJX1pqh9QmyXlLe+X+3jMm9ExQXzAJITDtASZLwwx/+EAMDAygUCvjD\nP/xD/Pmf/znS20dsRAAAIABJREFU6TRUVTWjVElVH2YTqqTBrGEYuHbtGm7evImdO3fi3LlzidgJ\nu+nMGoaBxcVFTE1NUXvmJqG/K1mOKEFyaWkJExMT2L59u2ukwKk+zNIVtkXYIvcm5KDa+sqyDH1b\n8wuyBFJFS4a2k84GvYoYGG3wqZoeCg20booyLztcnQWnNkGSE6rUfKwZKaAVfqEJsQDA0uIEzt9l\n84Czr6ybfBR+MYWNbjUvvfQSeJ43t91SqZTY9nZBpKpq5DBLjo9JWTfVahU/+tGPUCqVIMsyFEWB\npmlgGAaqqmLfvn34/Oc/n9g7bn2Y7YI2uzO7tLSEyclJDA0NYe/evSiXy4kAWaB7MFur1TA+Po5M\nJuPaMzcpbcKiihmIoojx8XHouo7Tp08jn/d/+zkpYL8Z5KdFWK5chqTySKWb27qkNttzEZGIgZtT\nS4sYiBF0NAiruCIGbqqmhwCgI6iNy50lIAvYC72IK2u4dDBgZNGEWF/qcMKX6eKuNbclK8QSnT17\nFoqioFqtolqt4vbt22Z7O+vFWS6XSyQQucmtg08nStrF/tatW/GFL3wBS0tLaDQaphurKArq9brZ\ndSGp31sfZhOqJMBso9HAxMQEAOC+++5DPp/H7OxsIroHEMUNs4qiYHp6GpVKBSMjI54NrpMCcCzL\ndtQJQ9d1XLlyBTdu3MCxY8fMIqcg6juzncnaIux/z2SQw8axwAmygHcBGE1+4wSdRgxE5JBHuBZX\nVgWNGACAyrR2kyAurd+8rFOV3HZkNP9Z2HbasTJh/pvR7Ptsu2ysL5CVxVaIpeVdaREDmRJFAB1i\nAUB56wcBAKlUCkNDQxgaGtr42zrgkvy4IAgm4JJ4TZIBN46YgSAIgQyCuLWwsIA33ngDjz76qOtj\nFEWJpUtLFOrDbELVSzBSFAUzMzNYWVnByMiI7aCUtMlkca0nwzAwPz+Pq1ev4uDBgxgdHW17oE1K\nz9tOYga3b9/G5OQkduzY0VGLtT7MRqdcSoWk8iim/VWVOyMGusGC3IkX1t3YVAdxAprc2n0xjIEG\n03pLnmHs20Zep4yz9aF2EQOaqukh5DRvwJacFf2Ov3kBbZWz9+F0urOkk8FQY+N3GpcGb4FZxqWD\nAbABsc68rDNioGfy4e+g5YsAzYkMeex3A9xKpYJqtYpbt26ZgGvt/pEUwFVVNXKIq1QqsQ4UCipZ\nlvHcc8/hl7/8Jfbs2YMjR45gYGAAHMeB53ksLCzgwoULOHv2LH7913+914vboj7MdkFhdsZe7MC6\nrmN+fh5zc3M4cOAAtVcqx3GQ5eRUssbhzFonWZ07d8737aWkOLNhYgaCIGB8fBwAcP/993fcgqUP\ns9HoZxN5ZFMb27fTlRXVNLJ8+/1RaBMp8OvKNvQceMroXBadbfcN1n5SDwu3fiVwxbZA66a6UUQG\nnbmzQ41rSEnN99e4NHh5o9BL59Pg5NZsrNWJdW3HhSbEAgCrUB7TzpXNh8sCE1c2iFKpVEv/ZjKg\nxAq4qVTKFlHoxQQ+TdMib0uVtIEJhw4dwoc//GF85zvfwU9/+lMzSlepVDA7O4s9e/bgk5/8JB5+\n+OEeLyldfZjtC8CGIzc8POwJcEkbsxslzIqiaLYb82o75aakwGyQ5bB2Zjh+/Hhk3Sn6MBuNBnKS\n6cpKqn2fFNXWW+VOV1ZU00hzdigN4so2dPsJnAaybnK6r+bv0X67sMKtbrAoIHq4rXBDKGvhMrTL\n7HYM6cH7zwLNgi+d481RuQRk9fVuBk6QZSWhbaSAlUUTYgOLBrA0WO7CHTnagBLrBL6bN29CFEUT\ncLs1YjqOmEG1Wk2UMwsAIyMj+NznPgcAeOONN7C2toatW7fi+PHjPV6y9urDbMIVd+VgrVbDxMQE\neJ735cglLWYQBcxagS5sRpQsy2aCWdKZYdeuXZFPbUsK2G9mWV1ZN5D1cmVpIOtHxMXlfIJrp66s\nH9XRPOn7hVpaXtYqBc311wnQBhGJGqQUATrHg1Xt35vO0zO8Wrr1eNzajisP+MkUO1xZI9u8WGec\n42sTJjfAJREFMmI6lUrZIgpRAm4cQxOSNjCBiHzWkydP2n6f1C4GRH2Y7YLCbgAECOLobyfLMqan\np1GtVnH8+HHf85aTMj6WqJOdyzAM3Lx5EzMzM5GM4U1KZradK9poNDA+Pg6O4/DAAw/E0juw78x2\nLgKy1nZcuZRMdWSBJrxmOMX1725SDB6KoyDML8h2WwRqAURSWAa0Aq1XXtaqoO6sxqXB6YoNZHm5\nbgNZqyurpXPUuIFVeiZPbcdFjRisi0As4AKyQcbeApg6/BYcDPSMzpVOpzE8PGy7kyRJkhlRIICb\nTqdtEYVMJhPqnBGXM5ukmAERgXYnvCYZZIE+zCZapKNBlDuRruu4evUqrl27hkOHDuHEiROBNtKk\nxQzCqlqtYnx8HLlczrXVVlAlZWd3c0U1TcOlS5ewuLiIkZERm9MRtfow25mIKysqHFKZjf3NCqpW\nV9YvyDojBg0tY4NlL3UrYuBXNZRRRCWS14rKoXUWf1mVVurQ2BRYNL83nU9D11u7jtDcWCLiygaN\nFOjZAhiugwImt7txuUJswwSCKpPJIJPJUAG3Uqng+vXrEEURmUzGFlHwA7h3E8yOj49j69atoe9Q\n9kp9mO2CwkIOgdkoQIs0/J+eng7U/N6ppMUMgsrqSI+OjibyNk+ncsKsYRi4desWpqenI3Gg/Sgo\nzPbqQuA7/9ueWaMtxvt+Ld5iJJqIK1vKbAArbdm84NUrYtDQOj+mAPFHDHTDezutoQkDUUCtH6Ct\nG/Z8qV93dnvjCjQ2hZQlH2st+gKarqwVZOlFYHaIbefK6maUoDWO0lG8YH1c7pX9DyMZHcfpcgKu\nYRi2iMLCwgIkSTIBl0BuOp22HZPiGJpQqVRs3R16LV3XwbIs/vqv/xrvf//7+zDbV3SKqtcscSGz\n2WzHt5WTFjPwK8MwMDc3h7m5uVCO9GaSFWbr9TrGx8eRSqUic6D9KOnO7P/8X0WkU/6+/7/7zxK0\ndWb7wzfHD7Y/vljAQE61ubKCwiOfth8LAscJdA4K7Cdkmiub1IiBl6wurd+8LE1+IwZBRANZp9SU\n9/tq6RzYAPCpZ4MVr5qiRQxSmeYPbbk0Del0uJ69vRDDMMhkMti2bZsJa4Zh2CIKCwsLEEUR2WzW\nhFtFUWKB2YMHD0b6mlGoXC7jlVdeweHDh5HJZJBKpcCybNfOHWHVh9kEq1OYlSQJU1NTaDQaGBkZ\nicSFTCrMeoXTl5eXMTExga1btwZqtbVZxTAMNE3D5OQklpaWMDo66jsTHZWiHqkbpf7n/wo/ivSZ\nn5diB1oCskSC0n57dXYxcLqygpoBz4Y/liQtYqBT/EDi0mbhnTP10rI6hCE+WNzAy50tq8vQ2A24\nJiBrdWXVVB68Ym/1RVxZr8gBVRwHnXMMOaC14/IpI1/yfL6maYmZBhlWDMMgm80im81SAXdtbQ2i\nKOLll1+mRhTCKmkxA3L+LBQK+MY3voHXXnsNBw8eRDabhaZp+KM/+qNEOclO3dln9YSo05hBUJEJ\nTtevX8eRI0ewY8eOyFzIJEIKcQGdn1EQBExMTMAwDHOC2Z0uwzCwtLSEpaUls09uL042QXvdRl0p\nS3s9K8TSXFm3t9ccHyNOoP3xxQIyfHP/KmVkE2Strqyg8Mil/B0XBDU+N0XUMshz4cGxndpFDNwk\nGjlkme4CLU1ldRkptbkcKUehF+DPjSWiubLWiIG27sRyHoVfVtEiBgafAiOLMPL+2kUpD78H+sxM\nYjKzUcoJuMvLyzhz5gxEUTQBd35+HrIsmw6uNaLgR0mF2dHRUXzmM5+BJElYWVlBtVrF1atXE9+d\npg+zCVZQmLVW55N2S3figcYp4hYTaOtmoVOSVKvVMDY2hlQqhYGBARw4cKBny9LLmAENijtxY2mK\nA2i/91IJe7fKpivrdGTJ/50g63RlgVaIpbmyYSIGDdUCWIyOhkZpHeXizBa46MbAtlNFK6PMRVMc\n5kfL7HakmI3vYVBZBKc1HU2NS4PjN/7Gy3UbyDpdWY33HwPTwsYJXOQE2XauLslZ3g1iGAa5XA65\nXA7btzenuBmGYQPcubk5KIpiiyiUSiUq4CYNZone/e534+LFi7hx4wZGRkZw4MAB8/MmWX2YTbCC\nwOza2homJiaQz+e7mo1MggjM8jyPGzdu4NKlS9i7d29XCp1o6nY/PlVVzfHDo6OjyOfzeO2117r2\n/jT1Gmat34ETZIO4sl6KGmhLed0EWZ7dcEHyadVX1AAAGmoaGS6aGBABYo4x4PdGtRvIAkBda3Ui\nDYNBka9THu0uWsTAKllv3tanAa1XXtYqmjvrLP7ykhVk04r987k5sgRiObW9u2qwXAvI0lxZGowa\nHG86s9ZsLaMGjyPEUeWfNHkdx9oB7srKCq5evWoC7u3bt7G6uorz58+jWq0Giv5973vfw5NPPomx\nsTFcuHABZ86cAQA899xz+PSnPw1ZlpFOp/HFL34Rb33rW1ue/+STT+Jv/uZvzCjFF77wBbzjHe9o\nedwzzzyDZ599FktLS6jX69iyZQuefPJJvPnNb/a9rL1QH2a7oE5iBpLkHfoXRRFTU1MQRREnTpxI\n3ESRbojjOKytrWF2dhaFQgEPPfRQz4oSyOCEbhzgDcMw4X3//v04fvw4GIaBoig9vyUUBmbjuAj4\nyo/yKHZgXjkjBk5FBbQ/m8ib8QJRYVC0XIsGAVm/ormygprxPyiBiW77qqnNLygo1PpRtx1aAChq\nqwBgTvgiUrkseNUefyCurJcba40Y6OvRgzDgaZWzQCzo6ykPv6f5OneBMxs0F+wGuIIgYGVlBT/7\n2c/wla98BXNzc/jTP/1TPPzww3jwwQfx4IMPek5gPHXqFH7wgx/gox/9qO33w8PD+Nd//Vfs3r0b\nr7/+Ot72trfh2rVr1Nf4+Mc/jk9+8pOey//FL34R//iP/4h77rkHAPAf//Ef+Iu/+AtcuHDB9zro\nhfow2yWFObl7ObOapuHy5cu4desWjh49im3btnXNDSR5yCQcxEibFUEQcM899/T8tg0ZnBA3zFar\nVYyNjVHhPQnTt4JkqxmGidTJ7aYrLAg6vvFcAR99tDMQk1TGhNmipa+sc5d2Rgw0g2kLsbSIAS1L\n2+suBlaoDZuXpckv0NZVu2MaJjtb1FaRUeomyBJXVuXosOqEWDdXVm9TCObHlSVOrt+uCH4Kx+4G\nZ5bc9etEDMMgn8/jLW95C97ylrcAAN785jfjS1/6El599VU8//zzePrpp7G0tIQvfelLePTRR1te\n48SJE9TXvv/++81/33PPPRBF0Ww3FlSCIEDTNBw7dsz83YMPPtjWVEuC+jCbYNEGFBiGgevXr+Py\n5ctd6xnqtly9hFld1zE3N4f5+XnkcjkcPXq05yALxA+SiqJgenoalUrFtU9uEor0et2ayzCMdVfW\nvo3SIgaCZCCfbf19O1c2KlmLvqyurKhyyKWCAaafiIGgpJHi/H04ziM2EJdqasGE2ajc2opWRo4L\nNtkKCAa0Xo4sgBZXljF8fAcM2wKyflxULVMAvw6j1jgCDWTDurLA3ePMxgHsDMPg+PHjGBkZwR/8\nwR8AaK7PTjoYff/738f999/vCrJf/epX8cwzz+DMmTP48pe/3NLlhud5vO9978MnPvEJ/MZv/Aay\n2Sx+8YtfmJGGJOvO3go3uZzO7OrqKi5cuIDV1VU89NBDOHjwYE8OJL0enLC0tIQXX3wRsizj/Pnz\nKJVKPXciieKCWcMwcO3aNVy4cAGlUglnz551zVv1GiTJMvTqO2EYBl/7SeuFTUPYWJ6GaJg/1v+H\n1TeeC59lEOXmPtwE2SaMimrrydPqyjZkHprR/k6MojVfW1DS5o9fkHWTW8TAKy9Lk+Fj+YlbS+Q3\nL0vTshTfgJQt+iIAgLdM9GJ1jerIanwGGk9xxi2urJrKt+120PK6mYL5AzQhNmyBGM2V1XNF6Dl7\nbvhucGbjGJjgFqn6nd/5HTzwwAM4deqU7edf/uVf2r7mxYsX8alPfQrf+MY3qH9/4oknMDMzg1df\nfRW7du3CJz7xCerj3v72t6NYLOKrX/0qnnrqKaRSKdfXTJL6zmyX1EnMQBAETE5OQlVV3HPPPSgW\no63MDqpe9ZptNBqYmJgAwzA4ffo0crlcT5eHJpKZjVKVSgVjY2MmxKZSHYyl7JJ6CdQ/v94sVCCu\nLIHYdJptC6zWv2fS3rAlCJ1/zz++2IQNUdl4LwKybq5sQ6Yftt1cWUFpn6XtdcTAKmfEIM5MrR8t\nq0PIcO7uZY5p5l6tIKtwGWT0jS4FxJUlEMsr9NZhti4HlMgB1UVlORNgvR4XZOgCkRNeW/5+lziz\ncfUmdwLt888/H+p15ufn8fjjj+OZZ57BkSNHqI/ZsWOH+e8/+ZM/wTvf+c6Wx0xNTeHv/u7v8Fd/\n9VehlqOX6sNsgmUYBiqVCn71q1/h+PHjnuHwbooWf4hTqqri0qVLWFpawvHjx1tabcUBkGEVpTOr\nKAqmpqZQrVZx4sSJRMQo/CoIzOq6jhs3biCdTqNcLkd64mi0AU63mLmqGlBVA4W8/xN1mOysKLPI\npDaWkebICusdDoJcGzSU5gUPrdDLr3oRMfBSTS1A1nkMpsMX3C1LAxjKrLX83pmXDaKsVofKpk2Y\nVbgMMi3ttrzzi5wqBndi0wXzuWFFg14tVwLL08FXO/122//vFpiN2pkNm2mlaXV1FY899hieeuop\nPPLII66Pu379Onbt2gUAePbZZ3Hq1Cnq41ZWVnD9+nUMDg6C53lwHLcpvuM+zHZJQYqzyC3lK1eu\nAEDPWky5qVsxA2s+eN++fa4DAEjRVRIUxbKQ7392dhYHDx7clKN3/cLs8vIyxsfHMTg4CF3Xcfny\nZei6jmKxiHK5jHK5jEKh4Hv7/8qPNoCgHcj6Ub2hBwLaICKuLBFvOV8SV5aAbJa3b1MZimtLAJaI\nBrK0iEEQV7abEQM3rcqlWIA2jLboi1DZNLJq8yKGBrIMvLdDJVWAwbTfxgh4EogFOgNZq7TcRhec\nMA7unaw4YgaVSiVw56Fnn30Wf/Znf4bFxUU89thjOH36NH7yk5/gq1/9Kqanp/G5z30On/vc5wAA\nP/3pT7F9+3b88R//MT72sY/hzJkz+Mu//Eu8+uqrYBgGBw8etEUHSOyBFJY/8cQTePTRR5HP58Fx\nHEZHR3H27NlI10HU6sNswrS8vIzJyUkMDg7i3LlzeOmllxIFskB3buuvra1hfHwcpVKpbautJMUM\nOnVmyecul8ubevRuu/UgSRImJiagKApOnz4NnufNA6qu66jVaqhUKpibm0O9XgfLsmYT8nK5jFwu\n5wn4zl0mnW7dh/xeH9CA1i1iEMSdtbqyosyimLO/pqD4P4FqerKOEXFI1jf2BRrQeuVlO9VsZSsO\nlJdsvyMgCzQh1ill3Y1Nt7Tjav5fSbnnWWkRAyvEeslvxEDLlTpu8bXZLrLDKI6YQaVSCXyn7fHH\nH8fjjz/e8vvPfOYz+MxnPkN9zje/+U3z39/+9rddX5t8j6Io4tSpU9i5cycuXrwISZJw/fp1/PZv\n/zbOnj2b6Iz05jxT3oEieVAAeNOb3oRCIdrJLlEqzpiBJEmYmpqCIAg4efKkr6tXjuMS0zokLMzK\nsoypqSnU63XfnzvJcisAMwwDc3NzmJubw9GjR7F9+3YwDANZ3jipsixrQiuRqqqoVquoVCqYmZmB\nIAhmLIH8fP25QQCAIGoo5MMfcFW11WUM4tD6Adp//lUJkgJk1vnLC2S9XFmSn83w0cZsRJUOhl7s\nUkp3Pto2SEuuIA6tsxVZp+5slhFNkOWMjbtUxJV1A1mgFWKdXQ5aHm/GCToDT8DuwAYFWWfEAPAe\nKHCnSNO0yPuWh4HZOEWMhNu3b+Od73wn3v721u8aQGJBFujDbNfkdgWrKAouXbqE5eVlah40iYoj\nZqDrOq5evYpr167hyJEj2LFjh++r/iRlZoMui2EYmJ+fx9WrV3H48GGcPHnyjnA7aDGDtbU1jI2N\nYcuWLYFdZ57nsWXLFlsrGUmSUKlUzDGSwCNUkKW5smFEgDaKwi8AKOfdXVkvuRWAWeU3YiCqPGin\nAY4SJ2i3WVbl1l6oUQCulzqJHPgF2ht1+4VllhGRMpoXz1aQJVI88rGGjwZCvCqaAEtEA1laxMAN\nUK0Q66V+xKBVqqqaxcZRqVqtJtKwWFtbw/T0tCvMJll9mO2RdF3H/Pw85ubmcODAAXN6k1NJGlBA\nxHEcFKV1HnxYLS4uYmpqCjt27MD58+cDX/1t1szs6uoqxsfHQ8Fd0mWFWVLIVqvVcOrUKWo3jjDd\nDzKZDLZt24Zt27bhKz/KQxC913u9vvH3YjGcw1Bv6GA7vNawurJOkG1IDPKZgOshoCtbkzdcV56N\n11lzAm4pLXSUl6VpVS4hz4fLjgYp/CJRAwKyKpOywSxraDaQtbqyMt9cDynVDos0V9YJsmGlZpr7\nGRmta1Wn8QKgeQ67Ey682ykpMYM4RT5jo9HA97//fdy6dctse5lKpfDAAw+Y08ySqjvn7LmJdPv2\nbUxNTWHr1q1tIYa05+rVeFaaeJ6HKHZeeFCv1zExMQGO43D//feHvvrdbJlZWZYxMTFh5pPiarUW\nx3hYvyIXYQsLC7h8+XLPCtkEQYOi0GID9CgCLWIQRk9+m8WTH6BvB5LSdGVJf1kiGshaIwYNhQsd\nJ6grKfA+p0DQXNmoVJVz0AwGAxTH1itiYM3L0rQsFjGUrQVeHklLIcP5vzBPMQpgNEE2o298BoXN\nIKM1qM9xA1mnlPVuBrxmf5xfVxbYANioRYsYdGtsd68VR060Wq0mCmYJgwwODuKRRx6BIAj4+7//\ne4iiiCtXruDLX/4y3vWudyXOWLOqD7NdEsMwqNVqJrxZ+6R6KYkw22nMQFVVzMzMYHl5GSMjIxga\nGup4eTYDzBqGgatXr2J+fj5wlCKoiNPZK5iVJAnXrl3Dtm3butIbl7iyBFAFofl/nqcUfq3/qt7Y\n2Gb8ZmyJu1vycHav33C/tf7Pv2reWqSBrJcaLsVgNLiVNRYpVkNdsTqwrY+L25X10tq6Y0uD2rBy\nA1ra6N4wOjpwA1mjAZWxb8s0kE2rggmxbiKurGLtLauFu82vZIpgadBLcWXd5BYxUPNl0LbOJINN\nlIqjm0HSYJYMUnj3u9+Nd7/73a6PS/L3ndwlu8NUqVRw8eJFHD582DfIAq1TwJKgsPBIWk69+OKL\nyOfzOH/+fMcg28nyxCG3zOzKygpeeOEFSJKEc+fOYefOnbGCZtxjdd2kaRomJiZw9epVDA0N4eTJ\nk7GD7NPf34AVQdBMkA2iekPD2lr8+5mkAGmLhVDM6VSQbUgMdKMJsW4gS1Ndaj62rkS7zqPaVJ1T\ny9bknAm2vdJ8pT1UOEE2owtQ2AwUlg7KTpClubJKKm8D2TBSMkUoAd1YPxEDNV82f9yU5Mr2KBVH\nzKBarbpOcOyFvva1r+Hy5csAgCeffBL/9V//Zf7t05/+NH75y1/2atF8qw+zXVK5XMbZs2dbZiG3\nUxJhNkw3AzKKt1Kp4OzZs9i3b19kMNcrcKPJuSySJOG1117DpUuXcO+99+L48eNdycb2Yp3cunUL\nL7zwAnK5HEZHR7t+N8GZZaW5su1kzdV6qVoLDszfvdBa8GEFWRIxcHNpvSIGdYkzQbYTxRkx8NKa\nnMOqRIfadhEDq5bF6G+zHx24AQA2R9YKsVZXVuZykDlvOJf5LGS+dcwtzZWlRQx0jm+BWJorG0R6\nKuMKsMzxX6M/5y5xZu+GmMHU1JR5vvjpT39qixG+/PLL5t+S3L2iHzPokhiGCQVvSYTZIDEDSZIw\nOTkJSZJiG8WbJGeWFIBZuzNYW1B1czm6deBpNBoYHx8Hz/M4c+YMMpkMlpaWuvL+xJUVBRX5fPvD\nmY/e9KjXNRQKrScvJ+hWa1pL3MAaMaDlZkXJQJpvbgcs5TZ/u7iBVarGQNXs75/iW18zaREDmkif\n3FUph8FMsOhBsyPDhoLkZxfr3uBJRujy2MjVsi5DEAjEpjX68lsBNhUiTqCkmy6u3+hAu8dZ87Ws\nHryg925yZuOA2SQ5s+l0GqdPnwYADAwMYHR01Pwby7Lmsia54K8PswlXUmG2HTzquo7Z2VksLCzE\nDnNJg9l6vY4XX3wRw8PDobozRCG3Pq9RikzsunnzJkZGRmxt5cJ0JwgrvyDrJdVRJOYGtEElyzLS\n6TS+e6EEUTJQLjT3AUECCg6OsoJsNm1fHqcrK8gsUlx31m+3z1/EoQ0KtVYFLQibr5Sxt1yh/k01\nePBME/ZUpJCGHUTbObEpVaI6se2kszw4yCbERiGNSwNcNHdM7hZnFoge4pLmzP7iF7/Ahz70Idx7\n7734+c9/jm9+85s4evQo9uzZg5mZmdiKlKNUH2a7pLA7Q5JAjahdzODWrVuYnp4O3WorqJKyjkRR\nxKVLlyAIAh566CHk89GdhIIq7pjB0tISJiYmsHPnTuq45W7ArDUr61SYiIFTfoCW5s5adfHixfU2\ndr9l/q4Jsobl/wxyPttxCbL75/Lryrop7oiBMy/bTmFcWqs6Lfwq8nWoBo8s01wGGsgysK9zpysr\ncXkYjvIpmivrVvjlBFma29ouYtAujuDlyrpFDIC7x5mNQ0lzZr/0pS9hamoKly5dwu///u/jJz/5\nCb773e+iXq9D07REgbeb+jDbRYU5wfM8n5jpVkRuvVRJtwbSly6bDe5GhFGvb31YXeg9e/agXq/3\nFGSB+GIGkiRhfHwcmqZ5tlPrVmZXFNzvWjQa9r8VisEPdwRo/WZpnbr//vttrqzg2JVpIOt0ZYFW\niO3ElRUVDsVM63qrSnTHjmEMlDLR9ZUOolUph3wq3HsLago53v9zre6sFWRVtBbTKWiuKyfcEklc\nc/93ixzfYNU3AAAgAElEQVR4SUo14TOlhW9/qPNp6HznDuz8/DzK5TKKxWLLBevd5MxGraTB7Ec+\n8pFeL0LH6sNswsXzPOp1f7PeuyUnPCqKgpmZGayurmJkZCRwkdtm1u3btzE5OWm60I1GA5UK/XZl\nNxV1zMDaVuzYsWNtG2gHvXALCr+f/oaGXKEJZNaIQaOhgudZyHLrawkNDTlab1lKH1qr2oEscWdp\nLbme/DaLQwebr08DWT+K8pqgJvHgWQM1qfXQTxsGwTDNZa9K7t0ROgVdkpd103Iji6F8532t/Yph\nDKhGc/1YQZaAqxvIpjXBhFg3OV1Zic+D1xVobMrxuNbP6ycrq3gMSghaJFbb9SYwlQquXbuGWq0G\nhmFQLBbN8dFxtKxKmuK6u1Sr1TbFrfvNpD7MJlxJzMwSkVZbs7OzOHDgAEZGRnruknZLgiBgfHwc\nAGwOZVI6K0S5HGRS2dDQkO/YSFCYDXrSyBXsrpPThXWKXSc1Yb23LA1qO5Fbd4NDB5twk+I39otC\nzvAFsg2xCXnZdPvvkRYxEBUWfhvWdDLVzAq6hsGgnG2FpqARA6fcgNZZ/OXUQiWP3WX6MAOaOFZD\nimluSzxDGVULd7fTCbJerqzER3PnRuPSzRxsSHlFDEqlEkqlEvbs2dN8L01DrVZDpVLB7OwsVldX\nwbIsRFFEuVxGqVRCPp+/o84BcUUp7paBE91UH2a7qLAxgyTCrKqqePHFFzE4ONiVpvhJkbXo6fjx\n4xgeHrb9PUkw26mroCgKJicn0Wg0Ak8q62YBWDuQpYlArSxrKBS8t91arQlnxWIwaNi/rwksomgg\nVdw4wVtB1hoxECQGmRTQEL1hwCtiUBM3TpA85XHd6GJQEZvriQa1nShqh5bWyYCArGLwNphlYNhA\n1urKysggBX+fNSqIBZpxBF7vbFBCEHEch4GBAfP2+NzcHBiGQaFQQKVSweLiIgRBQCqVQqlUMgE3\nm81uWsCNA2aT3N5qM6sPswlX0mBWFEVMTk5ClmU88MADibpVEvfEq8XFRUxNTbkWPQHuQxO6rU5i\nBoZhYGFhAVeuXMGhQ4dw8uTJwOs1Tpj97HfWp3zV5RaHNorCrzAiwGsVOQeWihvL5BYxdHNq/biy\ndYkFH8MhgkQMwioOqA0CtCuN5kWKX3d2T7kKoAmyOXbjPWQ9jQxLz8bK8FdkJvF5euEXxRmlRQxU\nNg1Ok808bfO5AaZ7BYgYyJkyMgdOtX2cpmnI5XLYsmWLLVomyzKq1SoqlQpu3LgBURSRyWRMwC2X\ny4maaOklVVVj6wu+WQE/qerDbMKVFJjVNA1XrlzBzZs3cfToUQiC4HuKWTcU5/hW0keVZdm2hW1J\ncmbDLEe1WsXY2BhKpVJHjnsYmA3y/dFA1k1sm3vo9brS1p0FmsDq5s6urQgY2LKxP+za2fy31ZUV\nRAOFfKsrS0A242NVyypjOrN1qUnGPMU4ormybuokYmCV4RIlqIhpqDqDgRwlftAmLysGmH7WifaU\nq8hxIhTDfkqkgWwaUluIJRGDqJ3YdgriyqqpHNJSE6blTPBqdbdb5el0Glu3bjVb9RmGAUmSUK1W\nsba2hrm5OSiKglwuZ3Nwk3h3Lw5nVpblRH7Wza4+zHZRm3FogmEYZqut3bt3m47k1atXE1UAQNpz\nRVldq2kaLl++jFu3brX0UXWTW6eHbitozEBVVczMzGBlZQUnTpzouNI2Lqj/7Hc4CHX6CTuoKyvL\nG9+TG9DSHFc/ymYZiKJhurJOkAXaF4C5ubIEYoOqV4MSVL35OdeE5oUADWqDiLiz7fKyVnm5s06Q\nzbEiZN39QskJsrSIgRNiw7iyotWF7SA6QHNl5UwZrK6Eglgiv90MGIZBNptFNpvFtm3bADTPK4Ig\noFqtYmlpCVeuXIGqqigUCjbA7fX5Ja6BCaVS6zTAvjpTH2YTrm5OcnKqVqthfHwcmUzGnOxElJTe\nrkRkeaK64r116xampqawZ88e10gBTd3MirZbDj8waRgGbt68iZmZGezfvx/Hjx+PxN1Oynpo58pa\n5cehpbmz1+bWbP8fGspCFC1ZWLF1PTQEA/lcgGVbLwajFXrFoU4jBl5aE9KRAG0+3flFvjVaQGQF\nWasrKxtppBnv5ZaRQarNY7wkUtxXGsgGiRgQ+QVXPxEDoDPQYxgG+Xwe+XweO3bsANA8FtXrdVSr\nVdy8eRPT09MwDAOFQsGMJ9BahMWpOAybSqWyKfq2bjb1YbaL2iwZGUVRMD09jUqlgpGREQwODrY8\npt3ghG4rKrhuNBoYGxtDKpXCgw8+GLhXblK+Yz/OKPms6XQaDz30UKQ5tjhg1urKOiMGqqpDVVs/\nr6EbKFCiAVZX1ior0IZ1ZQfKrAmzfkHWLWJAIBZwGYqQ8IiBm4hLS+t3G5es7qy1+ItnN5aBcxtV\na9D3Dasr6zc/66Zmm67wIEyLGIjpMngu2KAEv4q6zyxp/VUsFrFr1y7zPer1OiqOFmGk00K5XEY+\nn48NcDVNizwzm7TpX3eK+jDblynDMDA/P4+rV6/i4MGDGB0ddYUzjuMSkeUl6vS2tqZpuHTpEm7f\nvo2RkREMDQ1FuHTdl9f60HUdly5dwuLiYmyfNQ6YpcULRKF5UuZT7u5J3QKlNLBteXwbh9YtO7u2\nIuDQkS0myPLWdlzrEYOG4G+daLodZDtVryIG7bRUS2FrkQ5W7fKy11cz2DVIL8wixV9e2lOuQtE5\nE2YVLQWO23i9DCvZINbNlW0Hse0iBiSSQAPZMPECMb0BSp3EE9qpGxPAWJZt2yKs0WiYjyMObi6X\ni8RYiOMzViqVfswgBvVhdpMo7kr9lZUVs5fouXPn2l6NJjVmEFTWTPDevXtx7ty5O2KqjRtMkiEP\nu3btivWzBoXZdtv2p7+x8d0yLGNCrJcMvfX9Cdim0t4nqHrdn3PljBhk0s3PIcmGDWatENsuXtAQ\ngczmKPZuK5KX9ZIX0EYt4s4SkM3x64MQtBSynB063dxYIgNMa37WZ8Sgk8IwGvSKqSI4rnNzwW/E\nAOjdBDBnizCgGQeoVCqoVquYmZmJrEWYqqq2eF0UStr0rztFfZjtosLCKAG1OFqECIKAiYkJ6LqO\ne++9F4VCwdfz7oSYgVcmeLPLWYgmiiLGx8dhGEZXRg3HVQAmSSqyObvr5uXK0qQqWluYrVUkFErB\niPLQkWZ7Ikk2UMhvnOS9dvt6w8DQwLpr69Fxym/EQJQZFHMGaoIdMjgKc7jlYp3Rg1Iu/v08KNCu\nNZof3suddVOaU1tA1ipZTyFDuTVPJOnN40SabQ+uTldWZAu+e9K2kzVfy+n+QTaKiAHQHWfWr3ie\nx9DQkO0uk7VF2PXr1yFJUuAWYf2YweZRH2Y3gUhHgyh3Kmul/rFjx8wqU79KWswgCMySyv3l5WWM\njo7ekeN3WZaFoijQdR1Xr17FwsJCqO85rKKMGRBXlgaybqK5slHKmacd2lowXVmrBFFHPkd3r+qN\n5jI6ITaMK0sGLbAsfIFsEFWFJrCQVTqQt+9nQfOybgrr0AYB2uGCBEVvfh4rxBJXVtbp2xeJGBCQ\npcnLlRVZb5PAb8RAdBmUQH3NGCMGQO+cWb+KokVYXN0M+jAbvfowuwkUZXsua/V60Ep9q5IWM/Dj\nBBqGgRs3buDSpUvYt28fzp8/H1t0I+5YSDuxLItGo4EXX3wRw8PDOHfuXFddlKg/uyTRt/8wrizg\nr09tvSp7urNrKxvjSrO55nJYXVlBbN0eScSAgGw2E249yQoDuTt35ltEXFGgFWyDSFJaP3vckYMM\nT0bVtn43biALtEKsH1cWaA+xfkVc2E4KxNopSMQASD7MOhWmRVij0Yj8WFatVrF9+/ZIX7OvPsx2\nVWF3iqhgtlqtYnx8HLlcruPb6jzPQxSjGy3ZqdrBda1Ww9jYGHK5XOSV+06Rdmq9gllZljE3NwdB\nEPDggw/6jo4kUZ/+hmaCbNSurBvQ1irBblsfOrIFmTQDSba242rCktOVJRDrJporSyIGtcbG9pSm\nrIogXOE3YkDktkrXGhx0ncFgoXXf85OXpWmplkIhEyyiYnVn3Yq/9gzUAaBZ9MVtvD7L6DaQdUYM\n/DjPTldW1HNIMZQ+spSIgRug0tp00dSLiAHQ+wv2KNSuRVij0cDk5KQ5tjeKFmH91lzxqA+zm0Cd\nwqwsy5ienka1WsXo6Ggk4fPNEjNQVRXT09NYXV2NZBiAHxGXuNuuhWEYuHbtGmZnZ7Ft2zaUy+VN\nDbKAuyMLBHdlaWrn0LZzZwFAkjQAzWUp5FmqIwsAYVMXVogNqk4jBu2krwPrar35+WlQmwQN5poQ\np+gcctwG0Mk6j6xLPlbSmt97mrUDoNOVFfUMUuuvIeqdTUUkTi6PzqDTT8RATNsr6oNaG5sdZN1k\nbRF269YtjI6OIpVKoVaroVqtdtwirB8ziEd9mO2iuu3M6rqO+fl5zM3N4dChQzhx4kRkB6CkxQw4\njoMkbThqhmHg+vXruHz5Mvbv34+RkZGuHXxJ8VVcM71pqlareOONN1Aul3Hu3Dmsra3h9u3bXXv/\nOPTxr9AdUtLJgHPpFUvkdHJJxMCpIKNxiRauLAEAjt2zE5nMBlRbQZa4sg2h+btc1n6i84oY1Ne7\nH/Bc62NormxStFrnIgHam6scdgzSX8cac7DKLTs7mFNQSMtmVpaIBrIZTjYhFmgFWafE9fhBJxDr\nN4rQacSgkS67vsbArgOBXy8JQ1HiFhmawLKs6cq2axFGsrduLcL6MBuP+jC7CRQGZpeWljA5OYmt\nW7f6arUVZpmSBrMkM1utVjE2NoZisRh7pKDdssQt4jyvra3hxIkT5kHS7wSwzSQ/7bhoj/cbTwDo\nEQM/7iwACKKGQt4OTARk/UjVDKhC+8fR1IuIgZtW6xxkFRgquX92Wl42CjkjBsMFCRleNUGWuLKy\nTj8eWkHWS6IlQ5tiWo/NfiMGKjq7MmkXMWhYe87GmLe9U6XrumutQZAWYbdv38bCwgJ+7dd+DZVK\nJdAdwu9973t48sknMTY2hgsXLuDMmTMAgAsXLuAjH/kIgOaFxZNPPonHH3+85fmXL1/Ge9/7Xiwv\nL+OBBx7At7/97a6fE7uhPsx2WWGqvHmeh6L4O5ELgoDx8XEAwH333Yd8Pnw/Qy8lLWbAsixkWcbY\n2BgqlYoN7LqtboCktZDvwIEDLc5zXK2x4pLTvbC6srKoIJ21n/S5APfPzcEKvPtzaAMZ2snqyjpB\nNp9jbSDbzpWtC0bkvWXjjhi003KV9QRaN9WE5rrxcmf9StVZsBZH1gqxVldW0lLIcO2PsToYG8iG\nkQjL9DG0HkNpEQO/ICqkiuB0FTIXb+u9OyEv61dBPqdbi7DXXnsNzz//PP7hH/4B09PT+OhHP4pH\nHnkEZ86cwUMPPWTmdWk6deoUfvCDH+CjH/1oy+9ffvll8DyP69ev47777sPv/d7vtRhXn/rUp/Dx\nj38c733ve/Gxj30Mf/u3f4snnnjC92faLNo8pYh3sfw4s5qmYWpqCr/61a+wb98+3H///bGBLJCs\nmIFhGFheXsbCwgLK5TLOnj3b09s4cTuz9Xodr7zyCm7fvo2HHnoIe/fubTngkiK0zS4ayIaRIqlt\ngbW65m6L1qsbzyURg2ZWtgmyTvl1ZOuCYUYKnEp6xED3UeC1XO3sFHNzNVguWpQ33m8wpyDNb3w3\nnIsjLWku7bgsEQNRy0DU/EEszZXVwEFEzgayUUlIFc2fMAoTMdhsnQx6qXQ6jTNnzuDpp5/Gj3/8\nY+zZswdf//rXcd999+E///M/8f73vx/33Xcffvazn1Gff+LECYyMjLT8Pp/Pm+AqiiIVug3DwL//\n+7/jPe95DwDggx/8IP75n/85wk+XHPWd2S4rrDPrBrPWdlN79+4N3WorqJISM6hUKhgbG0Mmk8G2\nbdvMPFMv5RxYEJWsI3fb9cfdzDED4srKIt0pC+LKOhUmH0vTsXt2mjALwObKiqKOrMWJtbqyBHI1\n3f4ZOnFlo9jd3SIGUSisQxtGl+dUnDjCtoCsrHLIpTaOocSVJSDr5spaATZNKRSjRQxaXkPPIkXJ\n3vp1ZWkS2CI4RiV1h20VdcSgD7Phpes6Dh8+jCNHjpixAMMwQt3pfPHFF/GhD30Is7Oz+Pa3v93i\nyi4tLWFwcND8/d69e3Ht2rXOP0QC1YfZTSA3mK1UKhgfH0ehUOh6NrTXMQNZljE1NYV6vY6TJ0+C\nYRjMzMz0bHmsiuMW/+LiIqamprB7925fY2g3W8xgeXkZKysr+Nz/V0ImZynAiciVjVJbdw5AkjSU\ny2mbKyuuF39ls63fjdWpzWSihYAapd0X50mnrX8bKNIvsMPMnpApq5s4tGGg1ho3cCv+smowp0DV\nWaShrS+PHWQBdzeWKM0qvp1YmkQ9utv8CtLgIUNg2zuvQdp0daIkTf+KS3EcP92MrEcffRQ3btxo\n+f3nP/95vOtd73J9vXPnzuHixYsYGxvDBz/4Qfzu7/6ubboj7f3u1HhIH2Y3gZzgaAW50dHRntxS\n79VtbMMwMD8/j6tXr+LQoUMmyDYajcTAW5QgKYoixsbGwDBMoDG0myVmIMsyxsfHoSgKtmzZYgNZ\nmjpxZYlo7iyJGDRqIvJF+jquV2XcXljBkXt2NV9nHWQLec4EWZraRQ5oriwtYiDJBiSKwRYFU6zV\nNt7PDWyt8hMxoGm5yqKQDb5d+s3PPnxKg6qzyKfXC73U1pXDwPv9RS0N3Wi/ndFcWc3goBn296S5\nsn4lGHnwjOoLZMMoTMQAuDuc2biAnWGYFqB8/vnnO3rNEydOoFAo4PXXXzcLxABgeHgYq6ur5gTR\n+fl57N69u6P3Sqr6MNtlhbkqIrf0dV3H3Nwc5ufncfjwYRPk7hatra1hbGwMg4ODLR0akpThjSIz\nq+s6Zmdncf36dRw/fhzDw8OBnp/0mIG1J+7Ro0exbds2/F//7y1kLJFC4sqKjQ2C84LZbL6VCt1c\n2bBxg5MP7IUkaUil7bECcxksrqwo6i0ubFBX1pqn5Snn1TjMMSvYarqBwZLHgyMWKf6iyY8rq1ri\nGwRkiSsrac3/Zzn7NkEiBqJHJwNaxMAqUhTmJ3YA0CMGipECzygQjPa1DpzP9wHi6WJwNzizcbRX\n1DQtsouAy5cvY9++feB5HrOzs5iYmMDBgwdtj2EYBr/5m7+Jf/qnf8J73/tefOtb3/J0ejez7uxL\nqztEPM9DEAS88MILkGUZ58+fx65du+4akJVlGa+//jomJydx6tQpjI6OthxkkgSznTqzKysrePHF\nF6HrOs6fPx8YZKNYhjhVq9Xw0ksvoVqt4ty5c2Ylr9WVlQQZYkP2DbIAWh7fTm4FYY2a+2Q7SdJs\nWVm32/leTm071RuGZ2FYt6StZwxWq82fqHRrObrXsmq4tN61gtVbHFk3kCUSbb1l/YOiqGc67m4A\nAMJ6n1onyPIBoLVbEQOg78yGVbVaRakU7Orw2Wefxd69e/GLX/wCjz32GN72trcBAH7+85/jvvvu\nw+nTp/H444/ja1/7mnmueMc73oGFhQUAwNNPP42/+qu/wtGjR7G0tIQPf/jDkX6mpKjvzCZcjUYD\n4+PjEEURjzzyCHK56KthkyrDMDA3N4e5uTkcOXIEO3bscAX4pMFsmGWRZRkTExOQZbnjtmpJgVlr\nCx9N0zAzM4Pl5eWWaWzv+9Q1FAeazeMlQUYqEz4rS4CW5tQ6JdRlqKq/72pgSxaSpKG43nNWFFQU\nChvLSVxZArJ+XFhVNZBJN9cPGXXL863bOM2VdZN3Xja8CNB6ObW0vCxNt5aB7UPtH2eVKANZl6+U\ngKyiseBZS0u0lGqCLE0GGE9HFqC7stQ4Aa3fLCVioOgp8OvALHQ4MawTNbgyws5DvFuc2ag/Y5hR\nto8//ji1f+wHPvABfOADH6A+59/+7d/Mfx8+fBgXLlwItqCbUH2Y7bL8uqmqquLSpUtYWlrC8ePH\nIYpi4kCW3MqO4wp9ZWUF4+PjGBoa8jX0IUkudVCQtOaA20F7XMsQh0jnDoZhcPv2bUxOTpoFbOTz\nGYaBD37mpg1kaQqTlRUbcuiMrVd2FmgFWaC9G+sGt3VKAZcf9Yolbq8AQyEpyPpZgwLtjUUVB/e0\nHgey6eZ6VzTWVuglaRxY1n3dCmoKWd4OoO1cWUFturBpF4fXr5wQG8iBDfBYmcmCh4wG1wpQu7a7\nd0Npp7vBmSU50yhVqVQCO7N9+VMfZhMm6xjWffv2+apc75WIGxrl8kmShMnJSUiShHvvvReFgr9R\nj0lSEGeWtBaj5YA7UZgWcFGLZVnz+9R1vaWAzTAM6LoORba7V524sk4JdRG5gjuUhnFlRaEVJkRR\nRTbr/7trCM33Tafs+w7NlU2iltea29bQQHTL65WXdVMxqyGb0qFo9vUoaa0dDEjEQFCb25cTZL1E\nIBagg6yfrGxDy/l+LNBZxEBgmsdNDhoVZDtV35kNp/4o2/jUh9kuy8txW1tbw/j4OEqlErXVVtKm\nrhCYTaU6hw9d13H16lVcu3YNR48exfbt2xP1WYOI47i2E9sURcH09DSq1SpOnjwZ+dV6r9edYRiQ\nJAmvvPIKjh8/ju3bt9v+pus6dF3HH/4/N5BqMwmg0w4GXkAr1ETkPBxYAOBTrJmT9QuyVheWgKvm\n6HPlBFnX9+9BxMC5rG5aXjM6Atqg7uyVa3Z31gqyzkIvpwjEuonmymoGYwPZIFIMHopm3y78gmxQ\nEXjtlu4WZzYJMYO+/KkPswmQJEmYmpqCIAiuYENuGyfpajiqwQnLy8uYmJjA8PAwzp8/n6jPGEZe\nzqzVeT948CBGR0d7Dp5Rq1qt4o033oCu6zhz5owt+2sYBjRNg2EY+D/+ZByFwZINZp2urCzaowJp\nt9AkRbK4EVlo59DSRKIGfKq5PVrH4RYKKYiiN5gQiAWi7y2bpF2EAK3fvKxTYfKzQDNeQHNkAbS4\nsoYRbB8T1vvQBikGs6qhZZBi/R0bw0YMBD1v+X3re3Fwf/9OIgbA3ePMRh0z6Duz8akPsz2U1Y1s\nl5UkgxOSdADpdHCCKIqYmJiAqqodFzwRJcG9dsur1mo1jI2NIZ/Pd33IRTfkLPCampoyt1erGws0\nnePCYAlZC2AqstISjXC6sgRQ9XXnMJsPX0leXa0D8OfOSqIKvrjxfVlBlubKWkGWJporuxkiBm67\nO4kdFAvun8ErGxwkYnDlmorRQywy/MbreRV6SQqHNG/fH50RA0Hhkc6oJsQCdJClRQwUnUOKU9Gw\nDFnwC7JuogFuXc9TnV0ayMYtXdcjB72kKaq7jlb1ndn4dGdvjQkUAS0y0WnHjh2+3EgCs5lM521g\nolLYDgKkh+rCwgKOHTtmuwXdiawFR72UE2atkDc6OorBwcEeLl08WlxcxOTkJPbu3WsWeJECQfJD\nvhuGYfD+Ty+0gCyfCn44EhvN0bdOqLW6skRh3NlsPgVJVFGwgCxHGWhgvkdDhabb9+VOXFlR1FEs\nsKjV7TBGu8PLsgzgMRCgXPS3HH4jBlGqIRjI59zX643FDYgbPcRC1RgTZiWVActurHMzbqA0f+cE\nWasEhV9/jGYD2cDL72NaWJiIQV3v/AI/Dt0Nzmwc5lGtVsP+/fsjfc2+murDbJelKAp++ctfguM4\n3H///b47FLiNtO2lwsQMlpaWMDExge3bt0ceKYijIC3schCYvXXrFqampmyQdydJkiSMj49D13U8\n+OCDtgIvhmGgaRo0TQPDMOb38nsfeh3FLQNIpVNm8RcNZN2ysjoFtsSG5MulJUBLXFkvDW1vjftI\nkop83tqOq7ncQiO6fbNet/Sx5RhfIOtHlZr9dfzCbRAtLmnYtjX4Pr28oiKfaw+TBGQLmeZnkVQG\nuXQrrBKQpSnLqybEdiKSpeVDRhEAugPb0HLgKM5uEBiOM2IA3B2Z2X7MYHOpD7NdFs/zOHLkiK3H\npt/nJQ1mg8QMRFHE+Pg4DMPA6dOnI4kU0JYnjltDQcWyLGRZNi9azpw5kyhHPQpZ24k53XUSKSiX\ny3j11VeRz+cxMDCAwcFBlEolFLcMgEtxLV0MOhVxadk2RVBCvXUogjNqkEo3D41WV9YJsuZzLSCb\nzrYHOVrEQJJ0SJL9d14OcBRarTSBZ7AcrfsUFmjnFxTs3e297zpB1ikWhg1kaa6sH5ClF4NxAFRb\nQRgNZMNEDEinAwBUkHVTmIjB4uIiyuVyR8eku8GZ7Xcz2Fzqw2yXxbJsYJAFkjUUgMjPMum6jitX\nruDGjRs4duwYtm3bFuvy9Lq3qq7rWFhYwNLSEk6fPo2tW7f2dHniECnwGhgYaGknRgq8dF3HwYMH\ncfDgQTQaDaytreH69et44nNVKJIMLrVx8u7UlXVKEmVkPArFqis15Ired0RKgzlIllysRBmL22go\nMNpsbq69Zev217MWl4VRO4B3yroeCdQC3mDb7rpVsPTZDQu0Xjq6f+P1CMhaXVlRYZFPuR+PRJVD\nhmv9e5pvf1xtqGlkOC10ZwOaq1pXsx1na8OoWq1iYWEBkiQhm82iXC6bP36dyLvFmY0DZsOc//tq\nrz7M9kBheoAm0Znled6zBRVplL9z506cP38+9oNf2MlbUWl5edkc9LBly5Y7DmQ1TcP09DRWVlZw\n8uRJm8PgLPBiWdaMVBQKBRQKBXziyxoUSUa2DUgqkgKVZTyB1PW5625vO6AVakIL0BJ3NpXmTZAt\nFNMmyFpd2UZDQS5rdxH9uLL1mgolHX4/iJsflleb+89AufM3sgJt2MEQRARkCxmd6sh6gayohgeS\nhtrchmgQHEZ11eL+BwDZSCMG6zEDwzAgiiIqlQqWlpZw+fJl6LqOQqFgwm2xWKQet+8GZzaOoQl9\nmI1PfZjdJEoizHIcB1Gk3LIVBIyPj4NhmEC54CiWpxcwK0mS2ZXh9OnTYFkWFy9e7PpyOBVlQRyt\nwArsDkwAACAASURBVAuAeVFG2m2RAi+aaCBr6AYUyX5BxKy7jJKjiKtdP9oolErzSK9Dqiwq4B2R\ngEbDXzRCaChIpThkMizqNUsMgQKyNFc27ohBO61V9MiB1o9oUYNibgOErSCbS+sQFfdlTPN6W5Cl\nubKqzkLV219I+Y0Y1JUMUqy/U22QiEGnYhgGuVwOuVwOO3bsANB0XOv1OiqVCq5du4ZarQaGYVAq\nlUzAzefzfWc2pPoxg/jUh9lNIp7nIQhCrxfDJic8apqGK1eu4ObNmxgZGem6M9ltmDUMA3Nzc5ib\nm7PlRhVF6XncAdgYN9zJAZlknQG0FHhZ3VhrgRdN7/kfM+a/rfAapIOBtUMBrd+sM4NLc2erKzXP\n91BkFelsCrKoIF/auKWcz6dsIOvmygoO2LWCbFzqJGLgpbVKcxseKLNtIwZe6tSVJeI5++t4gWxD\n4aDDvl78uKuCwiPFdb7v1pXo2nS5qV1etqYFH6TAsixKpRJKpRL27NkDoHlcr1arqFQquHz5MhqN\nBkRRxNzcHLZs2dJx/japigPY+zAbn/ow2wPdSTEDAo+kan/37t1diRTQ1E2YXVtbw9jYGIaGhlpy\no72OO1iXIyzMuoE6+ZuzZ6yX+0tAlkvxbUGWcQEzwwFgsij7GqDgFTdwRg3KQ0UTZJ1q58g6IZYm\nv66sm3plhK1VdBTy3m9uzcs61RB05HPuz19esR/XrO4scWUVlbHBrHNzs0YMGgqHbCoYkJKiMBrI\n0iDYrYOBFWK9FEfEoKauj7CNuO8sx3EYHBy0tRR8+eWXMTAwgGq1imvXrkGWZeRyOZuDeyf0oY26\n+4wsy3ck+CdBm39ru0uURJglMQNSte907rott2EFUUpRFExNTaFer+PUqVMoFos9WQ4/Ylk28EUT\nsFHgNTg46FngZc3FuskKslaF6SnrFHFq09m0Z2cEArRermx5qAjFMsKKuLLiOqRmLXlZqysrrINv\nOm2/YEilwrvhvYoYeG2y9UZ7oHXT2pqCfC74CbyYM6CQQi/SwUBp/j+bbt2uGx7tuGhAqhlMJC26\nCMB26sAGiRhUlTx4tvUL8wLZY7uiG9JiGIZt5LhhGBAEIVT+9m5RmGNxX/7Vh9keKMzVXtJgVtM0\nzM/PY3l5Gffffz+GhkLMo4xYcTqz1jG0hw4dwokTJ1y/x6T0kiUxA79SVRUzMzNYXV3FiRMnfBd4\ntZOmaS0wS11en66sU/W1OtK5cCdq4s4qsmq24yKigaz5PIt76wTZbqpWU1Eu+88S+40Y0LS8omBo\nS7jc8vUbEnbtDAa0iqPQyw1k8ynNE2RpEhTeVxcDN8XhwHqpqthbGdJAtptyZvEZhkE+n0c+n8fO\nnTsBuOdvy+Wy6eDm8/nEHC+7pbvt83ZLfZjdJEoKzBqGgVu3bmF6eho7duzAwMBAIkAWiA9ma7Ua\n3njjDRSLRZw9e7bnfWz9KohDTAq89u3bh+PHj4cq8HLq8ScmAQDprP3E71X0BQCpjP/1S15HFmRP\noF26vtSyHESkgwFRvpQxQZYmgRJDsIrmytIiBpKkmTGDen3jNWltybxWeaViX54gcBtUnQCtX+3d\nnYJ1E8hldBNknRJlBoB9fbeLGLiBrJ+IQU1OdZyp9Qu4dSUHjtEgInljr/0cA2j5W1VVUa1WUa1W\nzfxtKpWytQdLwm14UgeQ9Nfsa0N9mN0kSgLM1ut1jI+PI5VK4cyZM+A4Dq+88kpPl8kqjuMgOTvP\ndyDiVK6srODEiRObrqWKH5iNqsDLKSfIasrGtsu2yfA6QdctkuB8XDugpWloZ+s0JCvIWl1ZWVSR\nztqXxY8rKzQUpNMbJ2gCrjzP2iAWcO+vG0RWuPULtu2ueSRpA8CcQOuVl7XKjzu7d3cKqgobzFpB\n1urKijJDjRtYJchsc+wt529Ygptqsvd67NSBFbU0PDpq+ZJXxCDDykACoJjneWzZsgVbtmzsd7Is\no1KpmA4uyd8SuC2VSl3P38bRyaBWq6FUap0q2Fc06sNsD7TZYgaqquLSpUtYWlrC6OioeSAi+cmk\nKMqs6s2bNzE9Pd3iVG4meRUaWgu8jh8/bhtmEbTAyykCsrqq2SAWoIOsW7wAsLu4fhxbGtDWVqrN\nv4kS1Z0lrqwiq8hQYFgWg+931oIwPtUKrd1SGLD1o6gcWmvx18F9KYgSUFivyZMUIJehjC6W3bcX\n4soKsv2iICzItoNYL/kB3JrSvHgMEhvodcQgaqXTaQwPD2N4eBiAPX97+/ZtXLp0qev5W1VVI4fZ\nSqXSh9kY1YfZTaJewJRhGLh58yZmZmawb98+nD9/viUnlSRFETNoNBoYGxsz3eck3PIKKze4r1Qq\nGBsbw+DgIM6fP287aJMLlKCRAiIryKYydjBs58g65czKKpJiAq3TlQ2joZ1boCoaUmneVvxFlM2n\nbCDrdGWdEhsK+BQHVWkPG0G6GES1m62sNAvmBgaigdrlFQW5nPd3urZm/5683FknyDqVTRs2kHVz\nZa0gm+FbH9MuYtCQ3TsbdBoxkLQUJM3f+o+6K8FmUdD8LfnJ5XKRnZP6PWY3n/ow2xdVtVoNY2Nj\nyOVyeOihh5BO9/4WVTt1ArO6ruPy5cu4efMmRkdHE5MD7kROmFVVFdPT01hbW8PJkydtLkGnbizQ\nBFldba5/J8i6KWjRlx+ItbqzxJU1/7buzhKQzRUyJshaXVlZVNr2b02nOc9sbVBFETHwIwKYUUCt\nIGhtgdaPSLwA2ABZqyvr5cYSGUarI+tXBGDDqJ0Day0Wi8tVbR8xiE69yH565W8rlQpmZmYgCEJk\n+VtN02KZ/tV3ZuNTH2Z7oKQ5mlZZc6Kjo6O23oJJV1iYXVpawsTEBHbt2hVZj1zSSaCXrWisMQPS\nB3jfvn0YGRmxbYOdurEA8K4/GTP/7QRZVWnestN0+u16PkDBFwAocvPknPK4wGqXn1UVbf21WpdJ\nFhVk8/bnOl1ZsaFAV+1gwlMKv5wTxIDeuLI0+YVaa17WqWpVRrUKbN/uf8qf050lIFvI0R3ZdtlY\nQWquz2w6HCgqWvT7qKzxkDV/p9fNFjHo9XGNiJa/lSTJBNxO8rdxOLOVSqXvzMaoPsxuIsUJSNbW\nU/v37w+UE41qZGqn4jguUGZWFEVMTExA1/XIx+4SV7SXB32WZSGKIn71q1+BYZiW2EQnBV5WWUEW\naMKrVe1OCqqDYPy08QKaUOsFtJXFVbB863uXtzZPKFb4zOTS1GEJThEntpftuMJI09xh8PZtEUND\n3e0PTYCW5sgCG66smyNL4JaALE3tIgZB4wTtIgZ1eWNbpD02CIgmOWIQB+hFpUwmg0wmQ83fLi4u\n+s7fxpGZrVarm66IeDOpD7ObSGTiVhwj9sbGxlAoFAJHCjqZMhW1/E7e0nUdc3NzmJ+fb5luFeWy\n9HJwgmEYqFaruHHjBk6ePBlpgZdVb3/fL5HJ56Cvr3fO0bYs6HahGzp0D+eVuLLW/3sBrVPDe5rr\nwemiWkGW5spGGSegqVsRAzctL4sAEBpqb90SArmzQLP4K72+vTjNsnaxAifE+nVlO4kT0CSrHGT0\n/thHU1Vpfh+ZTPQxgyQ4s37kN39LYgwEcFVVjSVm0Hdm41MfZnugsOBAOhpE1edUURRMT0+jUqlg\ndHQ01FUjubWfBJj1EzNYXV3F2NgYhoeHW4qfolQvYbZSqeCNN94AwzDUTgWdRgqI3v6+XyKVSbuC\nbKdyxgmcIGt9nBNoxVoDQLMQjbiz+VJz3KciKSbMKpLiORZXbMgtgwZormwcEYN6TUaxFG9WXVHs\n+8vycnQurbP4y6md25vbi6wY4HlrYan9cc52XAFTKaZUvXVbD+XAShsL0GmRWBSxgYbqvo08sF/D\n1Zsdv4VNSTneh5Wf/G2lUgHP81AUxYwndFoMXKlUsHv37ig+Ql8U9WF2E4njuEjacxmGgYWFBVy5\ncgUHDhzA6OhoR4CdlPZcXjAryzKmpqbQaDTwpje9iTqGNupl6TbMWgu87rnnHiwuLpp/i9KNBYC3\n/Z+vIG0ZUUoD2TCuLE1+3Nd2j8mXCsiX81AkBbliE9ZoIEtcWbFBRuWGI6dGXUZ5sBUK67VWIHcr\ngmMZBrWq/fFxwy1gB1qvvKxTQdzZA/uary8rBvK5jc8vyQayGa8BCa2iubLWiEF93cXNpMJPP5NV\nDrIaHuCigFZRpW+LLOv9ufbviPZYt5mcWb9y5m+vXLmCVCqFTCaDSqWC+fn5jvvf9vvMxqs+zPZA\nnTqznYi0ZSqVSpFMs4oKsKMQLWZgGAauXbuG2dlZHD58GCdPnuxKvtdv5CEqkQKv/fv3mwVeS0tL\n0HU9UjeWQKwVZBmGNbsYmL9jGahtYJ63bHtuIEskVGvN5/jokkBcWaLilhKy+aytE4JbVwQCsUHE\npzg06jLld63vwXGdbXth4NYrL+um5WURhYL3saFabV1XfoB2aKi5zLKysVyS3Py3E2Sd7bj8qu6R\npQ2imtQ8Raa51nXYaZsuN1nzsjU5/k4IQbXZnVk/0jQNpVIJW7dubZu/LRaLJuAWCgVX0K9UKv3M\nbIzqw+wmUicwqygKpqamUK1W/3/23j1KjrrM/3939XX6Mj0zuU4yuV8nCUnITJgk8HVXkGUXD1nZ\nRTesiisbiQgs7KoIohBXPYDwU3ZlWdbVVcBlQTgiuqtEEDAokgQiYiYzk7klmUsymVvf71X1+6NT\nNVXVn+qu6q7qy+TzOscjmanurrl1v/qp9/M8aG1tNSy7Y9YK2WJQPokIWWCfz4eOjo6ybpEpV8wg\nkUigq6sLDMPkNHhZLBak02lkMpmSGrwEBJFlJeJqI7wZyrcEQUomPSN7pEYt8bhkSvbfakKbTqXA\nkqYTxFPii6/NbhVFVlqVTacyspW2akgjBvHzAmslRAxI6BFZRuMbDqnccjyP+nrj5iJHo+mCQquX\npiYHnA6LKLLuOososiSUIqslYpBhc793pKosSUaTGQbJTPF/J6VEDMJJB/Hj+chXld221Jzn5dlY\nmVVCEna1/G0kEkEoFMLQ0BCi0agsf+t2u8UGM5qZNRcqsxUi33YmNYqRWWllcsWKFWhtbTW0MllN\nMQMB6eV2I8VdD2bLLM/zOH36NIaHh7Fu3TqxeiB8juM4+Hw+9Pb2Ynh4GD6fD36/H36/X/f2nKuu\nz64sttptBUW2qK+F40UJtWoQynxCq3wRcnndqPNmK4WkXCtpNJcAKWIQV1RgtYqsGlrFvxDc+eeS\nUCi7zlmL1CrzsiQmx2OYM8+t61zOnYvD6cz9OTY1OZDJ8HA65LECAWlVNpFUjxtIkUYMoonS4gSR\nhBV20gQEQlVWD2pyKshrKfdRCS6EyqzWBjCGYcSqrPS2Qv721VdfxTe+8Q3MmTMHNpsN77zzDpYt\nW4YFCxYUvO/nnnsO+/fvR1dXFw4fPoz29nYAwOHDh3HTTTcByD7X79+/H9dee23O7f/u7/4Ov/71\nr8Vq8A9+8ANs3bpV09dfi1CZrSH0ymwwGER3dzf8fr9plclqihnwPI90Oo1Dhw7JLrdXAjMzs0KD\nV1NTU94NXvX19Whvbxe7dwOBAE6fPo1IJAKbzSbKrd/vV21ukIqsFvTKmXI5glJqpVXZQqTiWZET\nXmylIptJZ2Qya2EYmcjmq8om49lzsOWpHkshSXOp8QK96JFaNYyc4CCIrMc98yZKzYcSSbI8Squy\nsezwBbBc6RXCSEK/mBUTMQgn5W+ObISqajVJqxoXamVWK9L87bJly/C3f/u3OHXqFG6++WZ0d3fj\nYx/7GM6dO4e1a9fikksuwYc//GEsX7485342bdqEH//4x9i3b1/Ox99++23YbDacOXMGW7ZswTXX\nXEN8fX/ooYdw3XXXFfV11BpUZitEsZXZZDJZ8Dhps5Ny05PRVEvMQFhDm8lksHPnzopvLDOjMpvJ\nZNDb24tQKISNGzdq3uAl7d5dsmQJgOzvSDAYRDAYxNDQENLpNDwejyi3H/l0n3jfVrstRzoZW+7P\nnWEs4PP8KuhZZ0uKC0hRVmcFkRUQRFaYeetyz0hdOpnJu1AByH4/BYktN1ojBlqYmowDAPyEZjSt\nkKqzpLyslOGTU2hZPrNFTyqyQjVW2vjlclpkEqtWlRUkFgBIxXlSVZb0sWTagmRa/vtIqsrqgSS4\n4YQNNkJllySytcKFUJk1+mtctmwZUqkU7r//frjdbnAch97eXhw+fBjxeJx4m9bWVuLH3e6Zv8VE\nIlEVM96rASqzNYTNZkM0GlX9PM/zGBoawtDQUNmanSotsyzLYnBwEOPj41i/fj26uroMG11WCkY3\ngI2NjaGvr484faKYBi+Hw4F58+aJY7t4nkc0GsVff6oHwMwUBMZqJYqskkKrXwGII7zU7leJkKm1\nqlRE1eIGHr/v/O3lQpxOZv+tFFlpVTaVOD8OTBHOJFVlSREDUlVWDaMjBvkIBhKGC61W6utnvpck\nkQXUq7ECTrtcZIslEs8KtVZx1RMxCCfMeTlVq9ZWIi8LZCuz5ew/qARmLE1IJpNwubJ/gwzDYN26\ndVi3bl1R93Xo0CHceOONOHXqFJ566inVn8c999yDf/7nf8YVV1yBBx54oOTxYtXM7L5WMMvIFzMI\nBAI4dOgQ4vE4Ojo60NzcXJZ3bEZMWCiWiYkJHDp0CFarFR0dHWhsbCyq4m0GRlVm4/E4jh49irGx\nMbS3t6OlpUX8uQoSK3z/GYYp+mdusVjOi+wMpEpqsSKrJJPO5JV9aXMYm1E/LhmNIxmNgz//vRZE\nViqjVqtVFFkSqURK/J/ytkZQ7oiBQCYt//0LBhIIBrJGqCUvq2RyPFb4IAWCyHrcDLHRK5Eo/Lca\nT/BFi2w6Y0Ekzoj/MwJlBTacsCKsElUgVWXVqIWIAXBhVGZ5njc0SiG8Jimfnz/wgQ9g06ZNOf97\n8cUX895fR0cHOjs7ceTIEdx///1IJHL/QO6//350d3fjyJEjmJqawoMPPmjY11ONzO63V1VMMdJB\nEsdUKoWenh4kk0ls2rTJ9PmpSipRmU0kEuju7gYAbNu2TXy3Kz2fSme6SpVZnudx6tQpjI6OYu3a\ntcQGLyNnxsoem+PA2KzgJeOyBLHlCV8Te/5DahVUJdIlBMLvjvTFUSqy4nEZNuf+lZLr8fvAsSyc\nbvnvgxRpVTYVT8FOaFRSorUqm05liA1l3vrcimginkadJ7eqbGTEQI1gIAF3gSkFanlZPRVaUkUW\nmKnKqomsEDGI5xFdLRGDaMICUhSaVJXVEzEgiatdl7RW/s12KVwImVmg+BGaeu7zlVdeKen+Wltb\n4fF4cOzYMbFBTKC5uRlAdsXvJz/5STz88MMlPVa1Q2W2hpDKrHQl66pVq7BgwYKKZGfKOc2A4zic\nPn0ao6OjWLNmjWyzlYAgs5WOGpTSABYMBtHV1YU5c+ago6NDtcGrVIkFgD/7myOyfzM2KyyM/MVa\na9aVVEFVCqhym5Z42/O/QyRZlt6/cH/Kx/I2+mUxBiFiIP3eCSKbypOF1VKVTcRTsBEG6DOElbRW\nK5MzAQHIRgxIHxfweLVdDtQSMSARi6YLCq0ahfKy42dCWLl25s2XdLsXIJdYlyv391cpsQ67vt/x\naMKY50FpxCCSmPnZ6qm2mkWlIgbAhVGZNRrh+doIBgcHsWTJEthsNpw6dQo9PT3EBrIzZ86gubkZ\nPM/jJz/5CTZt2mTI41crVGZrCEFmp6en0d3djTlz5pi6klUL5ZpmIHzN8+bNyxE85flUao2sFIZh\nkCZUGPMhNHiFw+GcKrvR1VilxAIqEYISf7ek0mkpUM1hpXNnVR5XKbE8x8M3p0GyUtcmiqzTJRdC\npcRqrcom4rkLEcwmGpE3tGmVWz0EJmNomKM/Bzs2HMCCloa8x7hc59cFpzlZZbtQIqXYhFCGzZ0v\nq7Uqq4ZUYAV0xQYMkN5qjB5cKJVZI4nFYvB4PLpu88ILL+C2227D+Pg4PvjBD2Lr1q04cOAAfvOb\n3+CBBx6A3W4HwzB47LHHxCt3V199Nb773e9i0aJF+OhHP4rx8XHwPI+tW7fi8ccfN+NLqxqozFaI\nYkSEZVmEQiEMDAxg8+bNuv84zMDsmEEqlcKJEyeQSCQ0fc3l3ryV7zz0SLXRDV75IImshWFyGrIs\njEVW7ZRSjORKt4TlW5IAZJvFSI+hPEdPg08msgKCyKZTWUG2OwpXIYWqrHSKQUZjdIJUlVWjmPFl\nkdBMJo4UW1CizMsqSSWy35dihTYfGzYvQPr849fVzXz/kikOdU7y90mo1joVUwxIVVllxCAaJ8cO\ntKIUXCFb69AovnoiBmroldZQIt8XTO6ON4rZXpnlOM7wq5yhUEj3VKFrr72WOD/24x//OD7+8Y8T\nb/Pzn/9c/O9XX31V30nWOFRmawDh8vrIyAhsNhu2bdtWNeM4zIoZSJc96IlRVHq6goBWmY3H4+jq\n6oLNZsP27dtlI8Wk1VgjNnipSSwJNeGyWLLHkyYR5JM05fGC2ApSyxKq2Eqhld6Ho84J+3mDYVkW\nDkUVVpBYNZRV2WQilROBIGWA9VRlrToEVw+C2LIsB19D/tWxWlAKbSnzZd1ue47IJlPn/60QWZfL\nIoscKEW2EFGDnI3UGKZVZPWiNy8biOdKa6VjDrO9MmvGJAO6/ct8qMxWCK0yOjU1hZ6eHsybNw87\nduzAoUOHqkZkAXNiBqFQCF1dXUUte6gWmS0Ud5Dmf9etW4c5c+aInytHpCAfhURWDVJll/RxKVye\nSQUAQWh5Dg6XK7s17PzPWSqybDojSq6Asiprd9qQTKTyHqMHPVVZowkH4qYIbSFIUQO32w6HQ/69\nEESWhJZJBkoyLJBRSCypKqslYiBkazUW302NGAQTub9/Vov+78+frjW3KgvM/sosy7KGjx4rpjJL\n0QeV2SolkUigp6cHLMtiy5YtskHJ1YSR8ihdClDssodqkdl8ldlgMIjjx49j7ty5ZW/wEtBbkS0G\nnuPFiQhqIiyNMag9tvQYh8sFnudy5suy53Oy+URWqNZyfOGKudbJDKq3VxHcUjekkQgHsgIjSG2h\niIEaWoR2eiJC/LhSZOvqrDKRlVZl40n1uIEUZcQgnuBh19kMJoXUGFbij7kognH5y26lK616me2V\nWTNkPRQK0cqsyVCZrTI4jsOpU6dw5swZ1Y59IzsjS8WIeao8z+Ps2bMYGBggZkbLfT5GQMruShu8\nLrroorI3eAGlSSxPkMB81Vrp8SSpVeZxeY5XPQ9BYpUia7VaRZEloYwbkEZqaa3KkiIGlarKsmzu\nz0KQ2jpP/mYxIS9LQk/EQKjOut12ZDKcKLOpFAdGZbZuPEn+21RGDGJxXpTZfCO6tBKNA1oHnOiJ\nGGjJywZj2d+bSoizGVwIlVkzYgZ+v9/Q+6TIoTJbIUiiMjExgRMnTmDBggXYsWMH8d2vUHmslg0s\npUp1NBpFV1cXXC5XTma0GKqxMsvzPMbGxtDf31+xBi+g8DSBYpAKq1RUSeIr/bhaxVH4uFRq7U6H\neDub0yFGE3KqsOf/LcyptRU5no1UlU3FU8ikyaO3SKh9r92+4rdwaSU0FUV9U3HNoRNng5i7UPuL\nriCybnf2+SiV4uCqk3//6pyMTGLzVWVj8fONYA5LjsSSqrL5IgbSTG2pk/q0Vk+DMaupldZ89z0v\n/TaGh/2or6+H1+s1rXpaTcUUM8hkMoa/vobDYRozMJnqMKILnHg8Li4BuPjii1FXp56BE8ZzVYvM\nFgvLshgYGMDExARaW1vR0JB/1I9Wqk1mK93gBWiTWDW51HNZnOc5ooxqfSwSUmFlrFZVkbUwFuKy\nBSmFqrLCBjBSVYax6ZlWkEfWwtnmLS1Sq+f7JCBsOitFaLWyYu08ZDIzkppSyceqVWMFnE6LKLFG\noLUxzIhKqVB1LRY1OS0mL7t06VKEQiGMjIwgEomAYRjU19eL/3O5XIZJ6GyWWbNiBqSrrBTjqG0j\nqnE4jsPg4CDGxsZytjypUcn1sUYxPj6O3t5eLFq0CB0dHYZWEKxWK5LJZOEDTcZisSASieD3v/99\nxRu8hEUExVRm9UwtkB4r/W8tjWCk+7Y7HeKEAZtk7JZ0ioGAjVHkEAuU4tLJNGx2myiwAnpexEqZ\nViBIrRStVVtSxIBEaCoKADKpzRcxENBanRVE1u22iSIrrcomkxyciiqssiobi3NgucLfR1JVNpPh\n4XRkPx6NSSZdlJCrLRQxCMWkyxO03285IgY+nw8+nw+LFy8GkK0whkIhhEIhjI2NIZFIwOVyyQS3\n1osiZmDGNINIJIJVq1YZep8UOfQ3uULwPI9Dhw5h4cKFqpECEtUqs1ouPQkVaIvFkrOG1iiqITMb\nCARw/PhxsCyLXbt2yX625YoUkChFauX3o09U9Uqs9byIkkRWOsVAQPr57L9zRZbjOHAFKoRqGFWV\nLYRUcHmOg9tX+pQCoLgqLUlopc1fS1Zm35zZbAyxIpss8L2OxdU/LwiqFqQSC5BFttSIQSxpAVB7\nlUibzYampiY0NTUByD73JBIJhEIhTE5OYnBwEBzHwePxwO/PxhM8Hs+sbu7SglnTDGgDmLlQma0Q\nDMPgkksu0f0OsBpltlCOV9rUprUCXeq5VIJ0Oo3e3l5Eo1Fs3LgR3d3d4gtDpcdtCRidm+U5vqAk\nq62oFY4niayAzW6TybDVoRTX3N85ZdyAtHyBdLtyVWW1Egtnr5cbIbVahDYSjGm6L0FkUykWNono\nu+qsMolVVmUFpCKrdkw+orHs7fVIr5J8ldJQVH6/erzGiLysWsQg331fuiIElmXyxpUsFgvq6upQ\nV1eHBQsWAMg+N0ciEYRCIQwNDSESicBms+XEEy4kWJaF02nsxj06Z9Z8qMxWEJvNBl7n/sZqlFlh\ncQJJZoU5ufPnz9dVgS6WSsistMFr+fLlaG1tFeVV+HylqrFmIxVV6X8Loqomstlj5N8HQWSF7WU0\nAwAAIABJREFU29jsdk0im5FMMyBJaimUqyqbD6nU5osYCHlZNUJTUbjc2hssSdVZqci63fI3HoWq\nsTyXvyILkAU1lc6O5BIkVu24UuIFkRgPm9WCRFLbfRgRGzCyUUx6RUp4/hOea/IJrjRXK5BOp8V4\nwpkzZ5BMJlFXV4f6+nr4fD7dr1m1Bl2aUJtQma0xzFhSUCrCOUnfzaZSKXR3dyOdTpd1Tm65ZTYW\ni6GrqwsOh0PW4GWxWESJNbvBSw8kuSwuS1tIXNQ/r1xlK1ROpbdRRgUEkZUuWVD+FZBEtharsmrf\nu1g4jkyahcdf3N9SMpbUJbNKmuZnu7FTKfnfVzLJwlUn/55KK66JOHv+Y1bVYwohFVk9qEUMIjFt\nQmbE1eZy5GWF5x2O48Q30tI31MKbaavVWvC5yG63Y86cOWLOn+d5xONxhEIhjI+PIx6P48iRI/B6\nvbJ4wmxpCjMjZhAOhw1rcqaQoTJbQQTh0UM1VmalAsnzPIaGhjA0NITVq1dj/vz5ZX2SK7R5yyik\n0Qm1Bq9MJoOuri40NDSgoaGh6G5isyuxegW3kMjmQ01kpQgiK/xO2ey2nE1hyvshPpYGQU0nzy9S\nsJFm6FoAQr+U6oY0ycddbvMuzUaDsaKFNjgRhn+u9hFBQnXWwlhgd1hFkRWqssmk+htHQWKBXJHV\ngiCwpSxKUBKK8Chh2VtezBzJFY5b0eglP+9PR2wAso2MgqRK35wJ8Sap4Oqp3losFrjdbrjdbjQ1\nNSEej+Oiiy4S4wmnTp1CNBqF3W6XxROMvlRfLuic2dqEymyNYbPZqqJbX4oQMwgGg+ju7kZDQ4Pu\nNbRGQVpWYDSBQABdXV3iimFlg5fwwtHR0YFwOIxAIICenh7E43G43W74/X40NDTA5/MVfNKsVKSA\nFBlQflwPJPmcqcjysuNkUwpI1VbCfZGOIy1TYKxWUWAF9Gz70iKyAJCI5U4rKFVwM+mZ70v0fL5V\nkNpCEQMpJKHNl5cVRFaKVGKVVVme42UiW4hY7PxaYgeDTEbLbFntEYNQRGMFVmXJA/lYzYfmuQ/y\neYVjDBo8LAJR+YNYGUFac7lue/6sM8MwsucoZbWWFE8QjlcKrrD9y2q1wu/3ywQtlUqJ8YTh4WGk\n02m43W5RbrU831UDZshsPB7PO3KTUjpUZitIMVW6aqzMWiwWDAwMgGXZotfQGoWZMQNpg9fmzZvh\n8cw01JAavGw2GxobG9HY2CgeE4vFEAwGMTo6inA4LGbWGhoa4Pf7xWaLas3FFoOFYXImGjBWa+7H\nFJZQSGSFai1js8pys4DKrFgdL1BmXE2QCi7P8ajz5r646f1ea63SJmPyN8BaK7TzFjUgk+ZkMpsv\nZpGIZwhxgtzvu9PJiBJrBuEoB5ut8pe91cQ3ECF/D60MckTWaEiiqownCM+hynhCPtFzOByYO3eu\n2OArPN8Jo8H6+voAZEeICYLrdrurLp5gxhx3nucv+CkRZkNltsaoJpnleR5nzpzB6OgoFi5ciNbW\n1oo/MZkhs6QGr2I2eFksFng8Hng8HixatAhA9okzGAwiEAhgZGQEyWQS9z0yO/4sSVEFklCqVVqV\n626lSxPy3VZX/tWEqmwhBImPR2am+5PEVivRYAwOV2mb80gIIlvnyV6XT6dYuOrk1+iFqmwirv05\nKRbLgGULf9/1xAvCUW1vAvREDIyYYhCMkL+GaitQao0nhMNhWCwWZDIZTfEE4fmuubkZQFaOw+Ew\nQqEQBgcHEYvFYLfbxextfX19yVsgS8Xoyuxsb5irFmbHq+YFRLXIbCQSQVdXFzweD5YtW2bodplS\nMHrOrNDg5XQ6TdngZbPZxGaLbDV2dvxJkkTWwlhy1txmK7Tyj1ltVqLIlgLp9moiW4nfY0FsXW5y\nzlAaMSCRSqR0C22+6ixJZEkoJTZfNjYWy6ge43Bo+9tRRgwiUa7kTK0REQM1adUTSVAreOfrNywU\nMSgFaTyB53kMDg5iYmIC69evBwBd8QQBq9Uq9hAIJJNJhEIhBINBDA0NIZ1Ow+PxyOIJ5axqGl1F\nFWS2Gl4fZzOz45WzRqnFmAHLsujv78fU1BRaW1vh9/sxPDxcFYINZJ9EjXgnzHEcTp48ibNnz2L9\n+vXi4HHA+JmxQHXFCkpBrXGMVL3UI5gktFZljYgXGFWVzQencbMXiVQihVQiBW+DV9ftghPhnDW/\ngsgCcolVVmW14HRaRYk1gojGCqxZEYNwlPzcwuj4Pai2qmwhEokEOjs74ff70d7erhpPkDaWCf/P\nMPln3wKA0+nEvHnzxHWvPM8jGo2Ko8FOnDgBi8Uiay6rq6urGTlMJBI0L1sGqMzWGJWU2XPnzqG3\ntxdLlixBR0eH+GRitVqRSqUK3Lp2mJ6eRnd3N3E2rrIaSyV2BrNElnQsz/FgU/K/A8ZmRYbL/duw\nEGInVmtu9VftXKU4XKV1aBfahlbqooRIIJIjtMq8bD5cbqcosjY7+eeZTAoVVuU4LvnPKZHIgGX5\nvMcA5KqsstoajbKaq7dqFBsxkDaR6RljXI6RXFo2L5bC2NgYBgYGsH79ejH7L0UZT8jXXCZUPLXE\nE7xeL7xeryyOJcQT+vv7EY/H4XQ6ZYJrL3XVm0mEQqGK9pFcKFCZrTEqsRQgHo+jq6sLNpsN7e3t\nOSNXhGkGtU46ncaJEycQi8U0NXhVYoNXNaJHYgHAYsltCLMwFrAZgnQSIgcWi55FBtp/RlqOTSVm\nxFB5Ho664rN+yqpsLBwXhbZQxEAJSWjzEZoIon6uHy63E3We7N92JsPKZNZVZxclFsgVWSmJRPY4\nu710m4tG83/tpUYMYnEe9d6Z+whFzou8jsqunqqs0bz11ltwOp3iZAGjMqeZTAY9PT1gWRbt7e2a\nRVFLc1kxs2+VzbQAxNW809PTOHXqFDKZDLxeryi3Xq9Xd1zAjHwrXZhQHqjMVpBiZKicl1akl9qV\ns1SlVOMiBz3wPI+zZ89iYGAAK1aswIYNG4pq8NLKhS6yWo8lVWrVRJYUOVC9X5Ou86bi2asTQi7Y\nqRjFVagqq0So0hbKwqYSuVdF9AqtILKZ828olJGCZIHRX06nVZTYUolGM7CrVIW1oDdiIAhsIYxY\nLqf2q1d8XnYnEokEgsEgpqencfLkSVHqBMHVK3XBYBBdXV1YunQpmpubS36+U2sukxYIpBEF4TaF\nBNflcsHlcmH+/PnifQrxhJGREUQiEdmGM2E1b76vRxg9ZiThcJhWZssAlVkKkcnJSfT09GDhwoUF\n19BWolpsFLFYDMePH4fL5TKlwUvKbJJYAbWFCzzH58hkOUVWDTWR1Zt9VTsXaYNbUjKKSym2UvJl\nZTPpNDLpNNw+j+oxakQCEdgLXFuPBiOYv3S+qsimUywcLmWcQP7vRDydU9EiVWWVEYN4LAN/w8zf\nWzQqVHRzv7elRgwSCQ4OO4NQWP69JomvEXnbckQMBASpW7BgAYAZqQsEAjh9+jQikYhsLqx0BKAU\nnudx8uRJTExMYPPmzaZubdQTTwC0NZcxDAOfzwefz4fFixcDyFaYhdm3Y2NjSCQScLlcMsGVjuEy\nY8ZsKBSildkyQGW2glRjgD2ZTKKnpweZTAYXX3yxpuB6tcmsxWIp+A5bWnVubW2VXb6iDV7FoazS\nCpVIQRSVkwxIcQPheFLkwMLkHsswFrBpDlbFZVCj4wWlkpCM4contmrEwlEAKEpq8zF/abaqlVF8\nv9WmF0hJxLPLJ5wufS8jcUlDmCCwxaA1YhCJsLDZLJpEVi+VjBioIZW6JUuWAMhGqILBIILBoDgC\nUFjg4vf74XA40N3dDb/fj7a2trLPRC1l9m2heEJTU5PYwMvzvBhPmJycxODgIDiOE+MJTqfTlMos\nlVnzoTJbg2iRNb1I19CuWbNGvHSjhUpPWFAiyLXa94c2eBmPEXGDvMeTZtZKjmXTM1u9suO+yOep\nnFMrvW+7U1suUEtVNh9CxbZYqZUKLSliICU0FUR9E3mNpsfvBstycLlnqqOuOrtMZElVWUFi9cCx\nvExiAfUGs1JIJMxZZW1ExEAv+SIGerHb7cSFBsFgEAMDAwgEAnC5XHC73Th79iz8fn/FFxqYtZq3\nrq4OdXV1skq2sJp3ZGQEwWAQR48ezYknFEsoFKKrbMsAldkaRJBHo4ZLCxmppqYm7NixQ/dllmqr\nzArno2xaEBq84vE4tmzZIruMRhu8ikNNYvPexgCRVUPvPFrpfSvX3EoRRFdP85mAWlY2GUuAZVm4\n3OSrH5k0+XyUQlsIktCSRDadzOTd7pVKZDTlfqURAyFH69DYDFZsxCASyT6OzZZ7rJ4KrJlbw/Tm\nZfNR6nxZi8UCp9OJ6elp2O12vO997wPDMOL67b6+PsRiMVOay0qBtJpXrbkM0B5PkDaNnTt3DitW\nrBDjCWfOnEEymURdXZ1s9q3WLWHhcJjKbBmgMltBipUlo2RW2r2/adMmeL36ZlQKVJvMKhcnCJvK\nBgcHaYOXgWiRTNLUAlIFU11w1V6AylsxIomudEyX1qosiURMWJigfRxXLByFrUCHeTySKzyMjYHD\n5QAryemmzzd3ScUWmKnKpoTpBA5CDpYQMTCqESwfySSLpPaJY6ag9jtYzrxssag1eSkXGkibywYH\nB8GybEnNZUajpblMz+xbITNrt9vFZTZA9jUiHo8jFAphfHwc/f394HletprX4/EQXzvC4TCWLl1q\n+NdOkUNltsJYLBbd40BKvazP8zxGR0dx8uRJotzpxahFBUYhlWtpg9cll1wiq9bSBq/SEJq/tFZO\nC0UOiOKr0mDGnT9WKhSlVGWLQTqmSy2ikK+SqXwDKJVataqslFgoAne9tjegoakgGuY3yiYjuNwO\nUWRJpHRKaSyaJsotqSpLihgUmmAgzdeSG7e0/zz1VGDNihgEgtnfbUZlA5kZ79eETV6Tk5OamryK\naS5zOp01E09Qzr7NZDLEK5MWiwVutxtutxsLFy4EkP37FeIJp06dQjQahd1uF7O3HMdh6dKlNDNb\nJqjM1iClyGw4HEZXVxe8Xm+O3M0WrFYr0uk0+vv7ce7cuZyB37TBy1gKSW2+BqtSIwecRBY5wsIE\nQXaVoluqyCpRVm615m9JJGJx2DQalFahVYqs1crIRFZalU0m0nC6Cp+/IK6xqP4MrVbSaQ7pdPFV\nbzNjAySCIRZWFTklxQnKKbLxeBydnZ1obGwsusmrmOay+vp608bh6TnvfJvLhNeCaDQqFkMKFTik\nEi+QSqUQCoXQ1dWF++67D1NTU/D5fMhkMqivr8e2bdsKNlU/99xz2L9/P7q6unD48GG0t7fLPn/6\n9Gls2LAB+/fvx+c+97mc2w8ODmLPnj2YmprCtm3b8NRTT1U8HlIOLDoratVTfpslpFIp3VXN3t5e\n+P1+XU1amUwG/f39CAQCWL9+veEZnjfffBO7du0y9D6L5d1330UoFEJLSwuWL19OG7zKQFHZWcai\na1yX3sfQGkVQTkIwinxCWiiWk0lmG7tcXnLlLBWXX2NXCq00YuCf2wCnOxuJYNPnM6yKCi2QlVgB\npcwqIwaJeBp2R+7Xp6zMqmVl1SqzUYUYOwjRBjVB1ZOX1TOSS/ljDISyPzsr4fdLj8gCxcnsTR+I\nq39ShbNnz+LkyZNYt24dcZOXkUiby4TcKc/zqK+vF+Wv0s1lSjiOQ19fH8LhMFpbW2G1WnNel7XM\nvlXCsiw+8YlPYPv27RgZGcHRo0dhsViwbds27Nu3D1u2bMm5TVdXFxiGwb59+/Dwww/nyOxf//Vf\ng2EYdHR0EGX2Ix/5CP7qr/4Ke/bswac//Wls2bIFN998s+ZzrjI0/5LQymyFMTtmwPM8zp07h76+\nPixZsgRr166tqicRI0mlUjhx4gSCwSBWrVqFlpYW8XO0GmseeiVTKquFxnUJx5olsoB8EoKAEYIr\nFU5HnfY1uILIAkAiElMVWi34587kH0kiC8glFsgVWSnCJAMtIquGUmRjsTRsNgYpDePAVO/T4IhB\nODzz/EoSTpLIViOZTAbd3d3geR5tbW1luRJnsVjg8Xjg8XjEdbQsy4rNZf39/YhGo1XTXBaPx3Hs\n2DHMnTsX27ZtE18XtMy+LTQ9QVj1fsMNN4i52VgsJk5LINHa2qp6rj/5yU+wcuVK2XZKKTzP49VX\nX8XTTz8NAPjEJz6B/fv317LMaobKbA2iVWZjsRi6urrgcDiIa2iNxIxxYVqRNnitXLkSLpdL1mlK\nG7zModhqrOrnCFVaQWx5jiw6ehYm6EEpuHrlVjkCTBBbR51Td7OkFqFVxg2kEmuzW4kim0mzsOdZ\nSyvAshzYeOkjr6KR7PfA31iHWKxwNIFUldWD1ohB6Ly42gnHq1VOawGhyWvZsmVobm6u6LlYrVZi\nc5mwjtaIzWXFMDY2hsHBQaxfv152boD21bzC64ra7FvlNAO3243LLrtM97lGo1E8+OCDePnll/Hw\nww8Tj5mcnERDQ4P4+tfS0oKRkRHdj1WLUJmtQWw2G+Jx9UtNHMdhcHAQ586dw7p168SB0WafU77Z\nrmYRjUZx/PhxuN1uMQN88uRJ2btp4R01bfAyhlLypspFCuJ96hx5JS5iIDSIAQAYYyW3VLkVSMWT\nYFlWtVIrrcpKkQqtMmIgIAitIKwsy8JJWIWbSavLtLQqm4hlz0UpvFqrsg67VRRY8bZOmyaRVaOU\nDKwwvoskqCSRNQojIwZa0NvkVSnU1tEGAgEMDQ0hHA6b1lzGsix6e3uRTCZ1Vay1NpdlMhmxcBIM\nBnMqqR/4wAdw9uzZnPv/+te/jr/8y78kPvZ9992Hf/zHf8w7dYh0lXe2XolVQmW2whTzi5avMius\noW1ubkZHR0fZ5NJqtSKTyZStoUwq7MoGL4ZhkE6nxe8RrcYai1Igi5Fbte1gWqRWy8Yu0gYxgLwm\ntxjSEum0O+XCSFrMIJ7X+aqstFKrlUIVWkFkhcdQiixjtcpEVq0qK0hsvmPyEY9mb5/WeOmfFBHQ\nU5VVixgI8irFiEqrWsRALS9rJIXyskY0eVUKaXOZQDqdRigUQiAQyGkuEzK4epvLotEoOjs70dzc\njHXr1pX82kBqLgOyEY/vfe97mJiYyJHMV155RffjHDp0CM8//zzuvPNOBAIBMAwDl8uFW2+9VTxm\n7ty5CAQCyGQysNlsGB4eFqMesx0qszUISWYTiQR6enrAcZzmNbRGUs5Zs1NTU+ju7sbChQtzhJ3n\nebhcLnR3dyMajYqXtord4EIltjBq47MKQZJStXmtguSWunqWJLl6BZdTjNsSxFYptVpIxZOi0KpV\nZQUyqQwiUyGiAOcT2UyazcnJKkknM7DZrTKRJaFWlRUEVkBPhlUralXZSDi30mvEdrFaixgITV6k\nS+a1CmneaywWQygUEntBhHmvDQ0NBZvLRkdHcfr0aWzcuFEmzUbCMAwCgQBuvfVWzJkzBwMDA4YU\ned544w3xv/fv3w+v1ysTWSBbuHn/+9+P559/Hnv27METTzyhWumdbVCZrTClVmY5jsPp06cxOjqK\nNWvWYN68eUafouZzMltmhQavZDKJrVu3qm7wamxsREdHh/iOfnR0FMlkEh6PR5RbLXksKrLFw3Nc\nXqHVK6WC5PKKXzG9s2VJKAW32OptOpkSxZ600EDt70Oo0mptWpMKMJAVWas9981kvjiBUHEVxnPZ\nCFMHClVlk+fFl2O15WlJ91eM9IbDcnG2l1htNzNiUA4q0eRVKaTNZUIOWK25TDo9gWEYdHd3AwDa\n29s1b+8qhnfeeQe33XYb7rzzTlx//fW6X+NfeOEF3HbbbRgfH8cHP/hBbN26FQcOHMh7m6uvvhrf\n/e53sWjRIjz44IPYs2cPvvSlL+Hiiy/G3//935fy5dQMdDRXhWFZVvfM2Hg8ju7ubqxYsQJdXV2Y\nO3cuVq5cWdFZft3d3Zg3b574DtpIlA1eCxcu1L3Bi+d5MY8VCAQQiURgt9vh9/vFd/TCiwCV2NLR\nEz0otdqqpFTBlU5TsCkqkcqqrPx2uVInldp8b/aETK5adTeTkj9HOOqccNd7wXGcbKYtl+FyZtwW\nqspmz1ObzHKZ3K9Rq6RqPU4ZMQiHk+r3qSKypMqsWqVVTWbVjq/USC4gN2ZQTU1e1YTQXBYMBjE5\nOYloNAqfz4fm5mZxHa3REQyO4/D444/jueeew5NPPol169YZev8XKHQ012yG53lxf/bmzZtVx3SU\nE7NiBkKDl8fjybvBC8jf4GWxWOD1euH1esWRXalUCoFAAFNTUxgcHMQ95AZRismobckqVnI5wu+h\nVsFVnotUIouZnCBs8tIq9+lkqmBcwempg9Ptkq1sFiQzn8imU8JYLfkxhUQ2dX50l43w9ZcismoI\n8lroPtUwM2JQyZFcUpGVNnlt2bKl7LGyasflcsHpdCKZTIozWXmeRzAYxPDwsNhcVl9fLxYzSmku\nm56exmc+8xk0Nzfj9ddfpz+PCkBltsLo+ePheR4jIyM4deoUGIZBW1tb1XQqGh0z4DgOAwMDGB8f\nR2trqywDJh2JUsq4LYfDgfnz52P+/Pm0GmswxeZo5fdB6Mw1SHCLqd5KG7uUYqs6VUF4bJZVnYCg\nnJSgFFqpULt8bjhcMxEDq9VKrJaK95VSbCZzaLsEnUqYt9VLIJnIwOZ1IBpR5G0NkFGt1GLEoJab\nvMpFOp3G8ePH4XQ60dbWJl619Pl8YjFD2lw2OjqKRCJRVHPZ4cOHcfvtt+Puu+/G3/zN31TNa/KF\nBpXZGkFYQ+vz+dDR0YHDhw9X1R+NMM3ACIQGL9JEBuUGLzpuq3YwQ3CNlFu1CjHx9hKx1XoOgrRq\nGeulFFqpxArnblcMmReqspnzM2W1rNSVVmVTCaGRLfd2pVRl49GU7NioEBuwW0sSWT1Z2Vpr5lJj\nNjZ5GU0gEEBXVxdWrlyJBQsWqB6nt7msvr4eHo9HtlTh0UcfxYsvvojnn38ea9asKcvXRyFDZbbK\nyWQy6OvrQzAYRGtrq+rWkEpjtVqRJmxS0kMqlUJPTw9SqVTeBi+AjtuaDZQyr1bAqOqtVG71zLzN\nVl2z/20lrK8lRR7YdFoUWtL2MYF0MoW6ei9sTgccLqfsvkgiK0isGqSqrCCwRqKcbGC1W0WBLYZy\nRwzUiITT8PtzIyChUBqNjeSPkxB+P5sa9TVqHTt27IJo8ioWnudx8uRJTExMYOvWrbov9edrLgsG\ngxgYGMCdd96JRCKBzZs3o6urC62trXjttdeKnpZDMQ7aAFZheJ5HKpX7gsLzPMbGxtDf349ly5Zh\n8eLFMnl78803sWvXrnKeal7GxsYQiUSwatUq3bfleR6jo6M4efIkVq1ahQULFuhu8NIDldjKYoTE\n6nu84n9fCoktSVYFqSV9Tgt2pwM2SVXWZlcuLZj5HMuyOZ/P3keu7Ngddpm8krOy+quy8Wgy77FW\nwuOQHjv78dzvt5rMkiqzajKrp/krTBj1Bag3eak9H6keX8Tv4+VLf0ObvPKQTCbR2dkJn8+HVatW\nmRq9+OlPf4pvf/vbmDdvHgKBAILBIDZt2oQdO3bgsssuw5YtW0x77AsQ2gBWy0SjUXR1dcHlcmH7\n9u2qO6sFuasGio0Z6GnwoiI7Oyg0tsv4xysumsBzPHhIcrKK/JyarAqrY/VKixArKCSyhbLpSiFN\nJ9Ow2W2GVWET8RQyhHm9WkVWDbMiBtFICj7/TM44HDw/Ck2tyUttBEEFoU1e6kxOTuLEiRNYu3at\nKdN0BFiWxb/8y7/gF7/4BZ588kmxcJPJZNDZ2Ym33noLr7zyCpXZCkFltsJI5YxlWQwODmJ8fDxn\nq5USYdZstVxu0jvNoBwNXlKoxFYXRuRnpfej57ZSuVUTTlJ8QZBXxmrVVHVVW90rxWJhxPFfNsUU\nA0FkpQ1erJUlHqMknUznPUZPVTYRV+ZaSxt9puf2eiIGsVgK9X4XQsGE7ONWKyMKrNGUoyoLgIos\nAeE1JBgMYtu2bXA6tW/T08vExAT27duHdevW4dVXX5U9ls1mw5YtW6jEVhgqs1XCxMQETpw4gUWL\nFmlaQ1ttMqtnmkG+lbu0wevCRe+aXOnxxa7YJYltoUYwac5Vy+MoH0MqsAJSkeUyLGx2e86UAptT\nOVJr5j6EvKzNbpOJbDEkJfKa0VgBJVVl9aA360oa32W1WnJEVi9qVdlyrKqlaEeY6DBnzhxs27bN\n1CuUv/3tb/HZz34WX/nKV/ChD32oaq6GUuRQma0wPM/jD3/4A3iex7Zt2zQHyUkrbSuJlphBKpVC\nd3c30ul0zspdGimgKMlXvc03Ckv5eT1iW6jSq3xcLZVhZWOYcsMYY7XKpiOQtocpRZYlNHuRK7CF\nq7LpZAZWu1UmsaTz1IueiIEasUhuNIIxIAZQS9MNvviR6nmerwbOnTuH/v7+glcvS4VlWXzzm9/E\nK6+8gp/+9KdYvny5aY9FKR0qsxXGYrFgzZo1ss59LVSjzKpVZqXzcWmDF6UUCkmslttoEVWSDOd7\nbFIOuJDEArk5XDWRVcoraXKCVoQVtsVAigjoqcqqRQxipEUJOmRYdQNXFeZf9XL48GFx7mlDQwNc\nLtcFWR3kOA4nTpxAIpFAe3u7qVclz507h5tuugkXXXQRfvWrX6n2rVCqByqzVYDX65Vt89FCtcms\nWswgEong+PHj8Hq96OjokO3EptVYSiUoRlTVEG4vzUKSFjJoEVkBjp95PKvVqklk1aqyysgB6ThS\nBdWIqqx0ysHM4xMeS4dwGlGVVUNvxKBcedlt27aJq1l7enoQj8dRV1cHv98v/q+Sq8zLQTQaRWdn\nJxYuXIh169aZKvMHDx7EnXfeia997Wu45pprLsg3DrUIldkapdpklmEYmcyyLIuBgQFMTExgw4YN\n8Pv94ueExi6hyYtKLKVSFFPpBdRX22oVWYuFycnmWh02mciSyCeyyqavUrOzuY9D/pqr5lGrAAAg\nAElEQVSVDWIAwKQII7YMiB0YQS1FDIDsc31TUxOampoAZJ8/4/E4gsGgbLi/dDVrXV3drJGwM2fO\n4NSpU9iwYYOpc9ZZlsVDDz2EX//61/jf//1fLF261LTHohgPldkapdpkVvrEKTR4kZrZaIMXpRZQ\niyKoVdeEebRSQRWEl1NKK0l4HYSKaZ5qm7Raq3VzmRFV2URMXm1Vm6ZQKjRikIWUl7VYLHC73XC7\n3bLh/kL1tre3F7FYDC6XS5Tb+vp62VWxWoBlWXR3d4PjOLS3t5t6/mfPnsVNN92Ebdu24ZVXXqma\nxmqKdmrrt3uWUsw7aJvNhmTSnHEzxcJxHN577z1kMhna4EWpWYoVWSlqlVs9IquM7ZAiB6RKsNZG\nsHykFNVWxsYU3DA2c07aJbJaIgZ6qbaqp9VqRWNjo9gQxfM8EolETvVWWM3q9/vhdrur7usQCIfD\nOH78OFpaWrBo0SJTz/P111/HF77wBTzwwAO4+uqrq/Z7QskPldkaxWazIRaLVfo0AMw0eMViMaxZ\nsyZnHzZt8KLUAvkmEpBEVm07WCGR5SVRAuUkA+H2JJHNeZwSc5JWuzVnkUJWmLWN2NMjyLpmy5oY\nR9C7KEHvSC6j87LFYrFYUFdXh7q6OixcuBDAzGrWQCCAvr4+sXorzd5WunorvJaMjIxg48aN8Hq9\npj1WJpPBAw88gDfffBO/+MUv0NLSYtpjUcyHymyNUuzGLaMRGrx8Ph88Hg/mz58vfo5WYym1gF6J\nBdRF1sJYiDlcq80qk1iALKMkEdbT3EMSTJZlic2ZNlbb0z+jY3GBWVVZ9fuYvREDo7FarWhoaJAt\nqEkkEggEApiYmMDAwABYlhUnJ/j9fng8nrJVKtPpNLq6umC329He3m5qU9vZs2exd+9e7NixAy+/\n/DKNFcwCqMxWAcXGDCops6QGr0OHDoHjODAMQxu8KDWB0SKrRC1/qlVk1ZDeXogeMDYG6RRpLqu2\nKEKp0qx+bOlSEo8m4akvvAUrEsouTVCTWfWfqTFSLNx/05zSN3aVY76sy+XCwoULxeotx3Fi9nZg\nYADRaBROp1NWvTVD/ILBILq6urBixYqcK3tGwvM8Xn31VXzxi1/EQw89hKuuusoQWX/ppZdw++23\ng2VZ7N27F3fddZfs8wcPHsQdd9yB9957D8888wyuu+468XNPPPEEvva1rwEAvvSlL+ETn/hEyedz\nIUJltkappMyqbSuzWq1Ip9Ow2+20wYtS1RQjsQJChVUqtaWKrBRp7MDCMMhwitFcVitYTmt+tTSR\nNLoqG4/MbOjy+nNnaycTaVF+o+GZYxnGgmgoXtRjCpRLZAFgajL3XJXMmadvtng5YBiGWL0NBoOY\nmprC4OAgWJaFz+cT596WUr3leR6nTp3C+Pg4tmzZYura3kwmg69//et4++23ceDAASxatMiQ+2VZ\nFrfccgtefvlltLS0YPv27di9ezc2bNggHrN06VL84Ac/wMMPPyy77dTUFL7yla/g7bffhsViQVtb\nG3bv3m3qMojZCpXZKqBWKrPJZBI9PT2qDV4Mw2B8fBxz586F3W6n1VhKVVKKyErheY44xQDIVlmV\nUwyArJRxhEv+pOys1s1l2cfTfqxZVdlkjNyQmiZIodVGnkPLWBmZxOZDTWRrJWIwOU7ueag2yXW5\nXHC5XGLFlOM4hMNhWfXW4XCIcqu1eptKpdDZ2QmPx4O2traSCx/5GB0dxd69e/G+970PBw4cMDQb\nfPjwYaxevRorV64EAOzZswcvvviiTGaF7WHKr/HAgQO48sorxbFrV155JV566SVcf/31hp3fhQKV\n2RqlnDIr3eC1evVq1Qav5cuX4+zZsxgeHgYA1NfXo7GxEX6/X/OaXoBKLMU8Cgkiz/GahVZvAxij\ncr+k6qnaeRIbwVREVmu8QA3hfpWTDdSmGpCb1AxYPVvm5qlKMyO51bl1imEYMXIgzGJNJpMIBoOY\nnp7Oqd76/X54vV5ZcWNqago9PT1Ys2YN5s6da9q58jyPl19+GV/+8pfxzW9+E1deeaXhjzEyMoIl\nS5aI/25pacGhQ4eKvu3IyIjh53ghQGW2RmEYRvfWsGKQNngV2uAlNA4A2UsvwWAQgUAAIyMjSCaT\n4liYfJemqMhSzEBPlVNaZdXbhW6EyKqhp3qa735TCXlFlLFaiYKq5/F0HasavzCvMmd2xMBo/r9b\nqlNk1XA6nZg/f77YAMxxHCKRCAKBAE6ePIloNAq73Y76+nrE43Ekk0lcfPHFuoocekmn0/jqV7+K\nP/zhD/jlL38pzuQ1Gp7PvQKj9apkKbelyKEyWwUU88tr9i88y7Lo7+/H1NQUWltbdW/wslqtOVtr\nhCc34dKUMNS7oaEBf/PpPlO/HsqFDWnCQLFxA7WKLICcqACgHjkAAIbgdcJ5Ke+L9DUwNiu4NGF6\nAuHxLIwFKULEgYSanJo12UAv1RYxKPforWqHYRjU19fLNnaFQiEcO3ZMjKC9++678Hq9suqtUVGD\n4eFh7N27F1dccQVeeuklUycjtLS0YGhoSPbYWvO4LS0teP3112W3/dM//VODz/DCgMosJQehwWvx\n4sXo6OiQiWqxG7wsFgt8Ph98Pp94WSUejyMQCFCRpZQdM0SWRKFKrTI/m11zmyucpHPSs5hB9fwM\neJE3QhTU5PRCixjMVsbHx9HX14f169eLBQ6O4xCNRhEIBHD69GlEIhHYbDZRbhsaGuBw6KtQ8zyP\nl156CV/5ylfwyCOP4PLLLzfjy5Gxfft29Pb2YnBwEIsXL8YzzzyDp59+WtNtr7rqKnzxi1/E9PQ0\nAOCXv/wl7r//fjNPd9ZCZbbGESqjRpBMJsX1gdu2bZNdAjJjZuxf/t2xkm5PoRSDltws+XYWWdNX\nIfRGDvSM/FJDTSyNuA8jqrL5VuVSZh8cx4nrddva2mRyyjBMToEjlUohGAwiGAxiaGgIqVRKrN42\nNDTkrd6mUins378fXV1dePnll00d8SXFZrPh0UcfxVVXXQWWZXHjjTdi48aNuPfee9He3o7du3fj\nyJEjuPbaazE9PY2f/exnuO+++9DZ2YmmpiZ8+ctfxvbt2wEA9957ryj7FH1YSJmNPOg6mKKdYlbT\nHjp0CG1tbSV3ZvI8j+HhYZw+fRpr1qyRLT4QPs+yrDhDlk4poNQqerKzM7fRX6k1W2T1VmWJlV2V\nY/XIrOqxOmWWdLx6zlhfxKAceVkjYwa1lpdVIxaL4dixY1iwYAGWLl1a1OuGEE8TBDccDovLH3p7\ne7Fx40YsW7YMp0+fxt69e/EXf/EXuOuuu0yNFVDKiuZfGlqZrRIsFgsxDJ4PYaJBKTIr7MD2+/0F\nG7yoyFJqlWIkNns747KzFsYCVqVn08LwOS/AZoqsGkZUZVXv22ZFPJI7f1Xt/Lx+T8mPWQ5oXjaX\ns2fP4uTJkzn9FnqRxtOEdbPpdBrBYBBvvfUWHnzwQUxMTCAWi+EjH/kIrrzySrAsS2X2AoTKbA1T\nykpbaYPXhg0bZEF9QapZlqUbvCg1T7EiC6hPNtCbnc0rxefPT7lyluHJtyHFICyMRbX5DISeL8bC\ngMvkmjXLkJ9PyOPDLEgjnfNxI6QakC9Y0HI/wnMUaRkDpTywLCvOIm9razNlW5jdbsfcuXNxzz33\nIJlMoq+vD3fffTdOnDiB73znO3j33XdRV1eHP/uzP8OXv/xlwx+fUp1Qma1hbDYbced6IcbHx3Hi\nxAm0tLQY1uCVDyqylEqid5KB+v3w529rEbeAifenU24LnYfeOIIeGL2RBhOrXGqPqVeIpc9hkeDM\nMgK9DWQWhkF9Y21UhKuNSCSCzs5OLF68GIsXLzZ14s7Jkyexd+9e7N69G4888gisVive9773Ye/e\nvQAgTs2hXDhQma0SSokZaCWZTKKrqws8z6Otrc30Bi8qsZRqpRjBzVdZVMqtFDUZNEpk9cYRjMAo\nCa1GQtPRvJ+nsitHWKozMjKCjRs3wuv1mvpYP/3pT/HAAw/g3/7t33DZZZcRj2toaMC2bdtMOw9K\n9UFltobRKrM8z2NoaAhDQ0O0wYtCKYJiZExt1a2AdASXIJ5mi2w1VWX1oreRyyzUZNc/p3SJq7Xm\nr0wmg+PHj8Nms6G9vd3UrGoymcQ999yDoaEh/OpXvzJ1cxil9qAyWyUU84SsRWZpgxeFUphixnUZ\nOY+Wy7CqDWLZc5PHicrR4KJeUTYoGqDzfvRSTMSgFIKTEfG/eY5Dw7z6PEfXPsFgEF1dXVi+fDkW\nLlxo6mMNDAxg7969uO666/Doo48atlyBMnugMlvD5JNZlmXR19eHQCCA1tZW2uBFoRiMWsW1mGqm\n3sovy7JgGAt4QmTeYmHApnKfFyyMBYw9VwJqKQJQSwTGQ7J/zxa55Xkep0+fxtjYGDZv3gy327yG\nO57n8cILL+Chhx7C448/jp07d5r2WJTahspsDWOz2YjzaYUGryVLlmDt2rW0wYtCKYBRTWKk5jDx\nc0UsRTA6V5tJpxUfZ4jTDiyMJW9zqcPllP3b7KpstUQMSiEwHpKtNW5aUPzIqkqRSqXQ2dkJt9uN\n9vZ2Uyuk8Xgcd999N86dO4fXXnuNLhOg5IXKbJVgRMwgkUigu7sbAGiDF4VSJGbMpAXyNIlxDPG2\nlZh0oJVUQv4mWm8WV/i40+0ifr5Uyh0xkEJ6Y0Riaiwo/nfTAn/V52Wnp6fR3d2N1atXY968eaY+\nVm9vLz71qU/h+uuvx+OPP05jBZSCUJmtYQSZlTZ4rV27NueJRmjwMipSAFCRpcxOzBJZ1dvlaRLT\n0yBWzLmpfa16v5ZiRRYAkrGZWbJaZdxdX5tzZDmVWAogiK25glgsPM9jYGAA09PTuPjii2VFEjMe\n67nnnsMjjzyC73znO7jkkktMeyzK7ILKbA1js9kQj8dx+PBhNDQ0FGzwotVYCiU/QlVNr9SWo0EM\nAHiV+7SpDKfXK7K1QCyUnSNrYSxw+3LFthjhp5BJJBLo7OwUR12ZWSGNxWL4whe+gEAggNdeew2N\njY2mPRZl9kFltkrQK5mZTAanTp1CIBBAR0cHfD6f7PNGV2OpxFIuJIzI0BrZIJZ9fPW/YWUWduY2\nM+dssxd+ui9nVdYIYuFYzsdqZQ2ukv7+fvj9fvj9flM2Z+llYmICvb29WLdunel51e7ubuzbtw83\n3HADbrnlFhoroOiGymwNIjR4LV68GF6vVyaytMGLQjEHIwS3mAYx4XZ6UZ5bJp3N1zOMBSA0d9mc\nlctsGrkeOBJUX3pAEt1K5GWVPH5PPYJBHlNTUxgcHATLsqivr4ff70dDQwPcbnfZmt04jkNfXx+i\n0Sja2trgcJj3e8HzPJ555hk8+uij+O53v4u2tjbTHosyu6EyW0MoG7ycTidGR0cB0EgBhVIJ9Ahu\nsQ1ixY360tc8BgCZZIp8XxLRtDpmXjL0VmWrATXR9TX6iB83inx5WQBwOp2YP3++uNCG4ziEw2EE\nAgH09fUhFovB5XKJcltfXy+LlBlFPB7HsWPHMG/ePKxZs8ZUgY7FYvjc5z6HeDyO119/HX6/MdMd\nXnrpJdx+++1gWRZ79+7FXXfdJft8MpnEDTfcgHfeeQdz5szBs88+i+XLlyOdTmPv3r04evQoMpkM\nbrjhBtx9992GnBPFfKjMVgn5njSEuX7Dw8O0wYtCqVLyVfhKETxOZUSW3jW5xaCsmErn15LOyp6n\nuqs+Xqvyl5TD02HZv82W20IwDCNGDpYtWwae55FIJBAMBnHu3Dn09fWB53nU19ejoaEBfr8fdXV1\nJT33j42NYXBwEK2trYaJpRpdXV246aabsHfvXuzbt8+wWAHLsrjlllvw8ssvo6WlBdu3b8fu3bux\nYcMG8Zjvfe97aGxsRF9fH5555hl84QtfwLPPPovnnnsOyWQSf/zjHxGLxbBhwwZcf/31WL58uSHn\nRjEXKrNVTigUwvHjx9HY2IgdO3bI5jAKEhsIBODxeOgGLwqlQhi9Qazg41kY1UyumhuWozFKLbsr\n4Khz5v28FspR+S1GbouNGGjBYrGgrq4OdXV14rYtlmURCoUQCAQwNjaGeDwOt9stVm99Pp+mbWos\ny+LEiRNIpVJoa2szNa/L8zx++MMf4vHHH8d//dd/4eKLLzb0/g8fPozVq1dj5cqVAIA9e/bgxRdf\nlMnsiy++iP379wMArrvuOtx6661iISgajSKTySAej8PhcMiWDVGqGyqzVYJSQjOZDPr6+hAMBrFx\n40bVBq8VK1bg5MmTiEajcLlcaGxsFC9D0RA9hVIeeI4zfQKCeEyBSqYw+SDn4wBsDu1P+WZUTFPx\n3CUvAk53neGPZxRSueU5HvVzjJOc7/9zcSO5rFYrGhsbxa5/nucRj8cRCARw5swZ9PT0wGKxiHLr\n9/tzxmpFIhF0dnZi0aJFaGlpMTVWEIlE8NnPfhYsy+LXv/61KaI4MjKCJUuWiP9uaWnBoUOHVI+x\n2Wzw+/2YnJzEddddhxdffBHNzc2IxWL41re+RRc11BBUZquQc+fOobe3F0uWLMG6devybvBqbm5G\nc3OzeBlqenoao6Oj6O7uhs1mQ0NDAxobG+H3+zVlrH757Hbxv2mVlkLRjrIyV8zl/kKCW6pgZhRr\nboX7ttrkFbxiHqfUimkyFpfc18zjS7eNVTqPK/x8QpMzq2oLiW2hvKxRWCwWuN1uuN1uLFq0CEC2\nKBIMBhEIBDAyMoJkMgmPxwO/3490Oo3x8XFs2rQpp1hiNJ2dnfj0pz+Nffv2Ye/evaYVWoRV7VKU\ngq52zOHDh2G1WjE6Oorp6Wn8v//3//CBD3xArPJSqhsqs1VEMpnE8ePHwTAM2tvb4XTOPIkXavCS\nXoYSnshSqRQCgQAmJibQ398PAOK79IaGBtn9K4nH43jgC1a4XC6sWbNGvPREBZdC0YZRK3Klgsuf\nT6oWO95Lfi4zzx9sTjU3+2+1+bXlRLltTMCs7WF6EcSW5znUN1XXZWmbzYY5c+Zgzpw5ALKvI6FQ\nCN3d3Uin07DZbOjp6RHzuQ0NDYZOL+A4Dk899RT+8z//E9///vexZcsWw+6bREtLC4aGhsR/Dw8P\ni6+HymNaWlpE2W9qasLTTz+NP//zP4fdbsf8+fNx6aWX4u2336YyWyNQma0SeJ5Hd3c3WlpaDGvw\ncjgcsg5Z6bv04eFhpNNp+Hw+MZpQV1cnNpudPXsW69atyxlcLa3cAlRuKRQ9GCW4epvCch9T23OI\nWgZWKblGVky1fj+k28N4joPLW7nNYMIkitBUVmyrTWoFwuEwurq6sGzZMjQ3NwMA0uk0AoEAgsEg\nhoaGkE6n4fV6Rbn1er1FxQ/C4TDuuOMO2Gw2HDx4EF6v1+gvJ4ft27ejt7cXg4ODWLx4MZ555hk8\n/fTTsmN2796NJ554Ajt37sTzzz+Pyy+/HBaLBUuXLsWrr76Kj33sY4jFYnjrrbdwxx13mH7OFGOw\nkErueSjP9ZILlHQ6LVZeAXPGbUkRxr9MT08jEAggEokgk8mgvr4eq1atQn19ve7Ho3JLoZSGUdMI\nBMFURgbMulSfb5KBHoqrXpObr4wUXLWGO0B9rBoAeBvUL+EXm5fVi7Dy/OzZs9i4cSM8HvXFEhzH\nIRKJiIWPSCQCh8Mhy94WahL74x//iJtvvhm33HILbrzxxrLNyAWAn//857jjjjvAsixuvPFG3HPP\nPbj33nvR3t6O3bt3I5FI4OMf/zh+//vfo6mpCc888wxWrlyJSCSCT37ykzh+/Dh4nscnP/lJfP7z\nny/beVOIaP7FoTJbRUhl1oxxW2pkMhn09/cjFAph2bJlSKVSmJ6eRiQSgdPplDWVaemOlULllkIp\nnWKWM+TDiJiCVvRIrpEiC1ReZoVFFQDQMC93PWs5ZDadTqOzsxMulwtr164tKq+aTCZFuQ0Gg7Kl\nDgzDYP78+bBareA4Dt///vfxxBNP4Pvf/z4uuugiE74iygUEldlaJJ1Og2VZWYOX2e9ox8fH0dfX\nhyVLlmDx4sU5jyd0xwpPYlarVczcNjQ0FDXGhQouhVIaZsy0LZfg5p9Fa6zMSilVbNVkNl9VViqz\nAlKpNVtmp6en0d3djVWrVolxMyOQLnX493//d/zf//0f5s+fD57n4ff78eSTT2LBggWGPR7lgoXK\nbC1y9OhRLF++HA6Hw3SRTSaT6OnpAQCsW7cubzOYlFQqhWAwiOnpaQSDQXAcJ2sqU45+0QKVWwpF\nO2YtZyBRzgquo07/c0cxs10FKa2rV7/Unu92OR/XKbNSXvj3tbrOQSs8z2NwcBBTU1PYuHEj6urM\nHXv2+9//Hp/97GfR1tYGq9WKI0eOgOM4tLe3Y9euXfjQhz6UN9pAoahAZbbW4HkeN910Ew4fPow5\nc+Zg586duPTSS7F9+3ZDnwR4nsfIyAiGhoawevXqnGYzvbAsK15+mp6eRiqVkjWVFbNTnMothZJL\nMVVLM/Kx5RJcLXJbisxKKSS2xeZl88msWSKbTCbR2dmJ+vp6rFy50tR54xzH4Xvf+x7++7//Gz/4\nwQ9kywlisRjefvtt/O53v8PNN99MFxBQioHKbK3C8zxGR0dx8OBBvPHGGzh8+DAcDgc6Ojqwa9cu\n7Ny5E42NjUVVbaPRKLq6uuDz+bBq1SpTdntLLz8FAgHEYjG43W5x3q3X69X95ErllkLJYkRzWK1W\nb5VyW+zGrUIVVrc/t2Gr1LwsCTNkdnJyEidOnMDatWvFcVxmEQwGceutt6KxsRH/+q//Cre7ctMk\nKLMWKrOzBZ7nMT09jd/+9rc4ePAgfvvb3yKRSGD79u3YuXMnLrvsMjQ3N+eVW47jMDg4iMnJSaxb\nt870vdvK84/FYuLEhHA4DIfDIVZu/X4/bSqjUIrkQpbbYqYn6JVSQWyrPWLAcRz6+/sRDoexceNG\nzbGxYjl69ChuvfVWfO5zn8NHP/rRsk4roFxQUJmdzUSjURw6dAgHDx7Eb37zG5w7dw6bN2/Grl27\nsGvXLqxevVqsfh44cACjo6N4//vfj6VLl1bFittEIiFrKmMYBn6/XxRc2lRGoRSH0WO9jMJowVXK\npaNOm7wVI6Ucx8NDqNYWul25ZDYej+PYsWOYO3culi9fbqpYchyH//iP/8CPfvQjPPHEE1i/fr1p\nj0WhgMrshUU6ncbRo0dFue3v78fy5csRiUQQiUTw2GOPYdOmTZU+TVWEod3C/1iWFZvKGhsbaVMZ\nhVICegRX7dI9YzNWRkuV23wVVjWxLTb3Kl1Hq5TaSkcMzp07h/7+frS2tqKhocGQ+1Rjenoat9xy\nCxYuXIhvfetbpjeVUSigMnvhwvM8fvSjH+Hee+/FpZdeing8js7OTjQ3N4uV27a2tqIEsVywLItQ\nKCQ2lSWTSfh8PnFigsfjoU1lFEqRqMmtngyqEXIrlUurXXt+P5+UKpGKbTEyy+W5DZdh4VXZ9GW2\nzLIsi97eXiSTSWzYsKGoq1l6OHLkCP7hH/4Bd911F/bs2UNjBZRyQWX2QmR4eBif+cxnMG/ePHzj\nG9+Q7eM+efKkWLl955134PF4sGPHDlx22WXo6OiAz+er2iconudlTWXRaBR1dXVi5dbn89GmMgql\nSCwMU3QzlUAxcqsml4XEVo/MSlHL2GqtyuZ8LpNdKUwSWjNlNhqNigWKlpYW02MFjz32GH784x/j\nySefxNq15kxgoFBUoDJ7IdLf34/Tp0/j/e9/f97jeJ7H+Pg43njjDRw8eBCHDh0Cx3HYvn07Lr30\nUuzatQvz5s2rarmNx+OypjK73S5b5pCvqSyZTKK7uxs2mw1r164VqxpUcCkUYygkt3qEVCq3xYqs\nUlgdkitTpcqsgFRqzZLZ0dFRnD59Ghs3boTPp74m1wimpqZw8803Y8mSJfjmN79Z1VfzKLMWKrMU\n7QiVzzfffFOcmBAKhbB161ZxYsKyZcuqVm6BrKAKsYRgMAiLxSJrKnM4HOLYM2HG7ty5c/PeJ5Vb\nCqV0lGJbrJACxWVt88lqvokIekRWwNtUj2ggDABwesiZ0mJENpPJoLu7GwCwfv16U8YqSjl06BDu\nuOMO3HPPPfjwhz9c1c/9lFkNlVlKaSSTSRw+fBhvvPEG3njjDQwPD6O1tRW7du3CpZdeivXr1+se\nqVVOMpmMrKkslUohk8nA4/Fg7dq18Hq9NHdLoVQAoyYuaBXbfDIrRSm2xcgsIJ8EQRJavTIbDofR\n2dmJpUuXYtGiRbpuqxeO4/Dtb38bP/vZz/Dkk09i9erVpj4ehVIAKrMUY2FZFu+99564zKG7uxvL\nly8Xm8q2bt0Kh0P/3Eez4XkeQ0NDGBkZQUtLiyi5iUQCXq9XjCVQuaVQyo/ZYqtVZKUIUluMzKqN\nNJNKrVaZ5Xkew8PDOHPmDDZu3Gj6OtjJyUl8+tOfxqpVq/DQQw+ZPquWQtEAlVmKuXAch76+PlFu\n3333XVPX8BZDNBrF8ePH0dDQgJUrV8oqyTzPIxKJiNGEaDQKl8slxhLq6+tpUxmFUmZKlVul1BYj\nswJWlQkBWquySgSh1SKz6XQax48fh8PhwNq1a02/Cvbmm2/in/7pn7B//35ce+21hsQKXnrpJdx+\n++1gWRZ79+7FXXfdJft8MpnEDTfcgHfeeQdz5szBs88+i+XLlwMA3nvvPezbtw+hUAgMw+DIkSM0\ns3thQmWWUl5Ia3jtdjt27NhR8hpevXAch1OnTmF8fBzr16/XtBOc53kkEgmxqSwUCsFms4kTE/x+\nf1E5NSq4FEpxlCq2Ri1+UEptvqqsWNVlWbCE437+5JaCjxcIBNDd3Y0VK1ZgwYIFRZyxdliWxSOP\nPIIDBw7gqaeewooVKwy737Vr1+Lll19GS0sLtm/fjv/5n//Bhg0bxGMee+wxvPfee3j88cfxzDPP\n4IUXXsCzzz6LTCaDbdu24amnnsKWLVswOTlZsKmXMmuhMkupLMo1vG+++Sbi8Y6dymgAACAASURB\nVLiuNbzFEAqF0N3dLW7DKWXjWSqVkjWVARCXOTQ0NBR1GY7KLYWSpdBIMKMiCDP3V9pzjdVuB3P+\nPvLldTlWLrGC1BYSWWGE4sTEBDZt2mT6UoLx8XHs27cPra2tePDBBw2Nif3ud7/D/v37ceDAAQDA\n/fffDwC4++67xWOuuuoq7N+/Hzt37kQmk8HChQsxPj6OX/ziF3j66afxwx/+0LDzodQsmv9ozW2J\npFywWCwWNDU14ZprrsE111wDQL6G94c//CHOnTuHiy66SBwHJl3DqxeWZTEwMIBAIIANGzbA6/WW\n/DU4HA7Mnz8f8+fPB5BtKgsGgwgEAhgeHkY6nYbP5xOjCXV1dQXl/JfPbpf9m8otZTaST1SlkqpX\nWPMJqcWS/75IkQOplErPhVQFLCTDSokFgMe/mr2qUyhylUwm0dnZCZ/Ph7a2NtPXjr/xxhv4/Oc/\nj69+9avYvXu34UWFkZERLFmyRPx3S0sLDh06pHqMzWaD3+/H5OQkTpw4AYvFgquuugrj4+PYs2cP\n7rzzTkPPjzL7oDJLKRsejweXX345Lr/8cgDyNbz33Xcf+vv7sWbNGnFiwqZNmzRd2hcuyy1atAjt\n7e2mRRlsNhvmzJkjLqPgOA7hcBjT09M4ceIE4vE4PB6PGE3Q0lRG5ZYymyhFVPPfb3ESK72dBWRx\nlaJ2KVuryP7sv7Jrw4XnhmAwiIGBATGTL1zVqa+vFx9LELg1a9YUHBdYKizL4uGHH8Zrr72Gn/3s\nZ1i2bJkpj0O64qt8LlQ7JpPJ4De/+Q2OHDkCt9uNK664Am1tbbjiiitMOVfK7IDKLKVi2O12dHR0\noKOjA5///OfBcRyOHz+OgwcP4tFHH8WxY8fQ3NwsNpUp1/AGAgEcPHgQixcvxubNm+F2u8t6/gzD\nwO/3w+/3A8g+OUejUQQCAZw8eRKRSAROp1PWVFYo96WUW4AKLqU2MDoWMHO/+kU2720MFtkX/7M1\n52PS54alS5eKmfxAIICxsTGx+sjzPDiOw8aNGzVl+0thbGwMN910E7Zu3YpXXnnF1OkzLS0tGBoa\nEv89PDycM1ZMOEaYMhMMBtHU1ISWlhb8yZ/8iSj2V199NY4ePUpllpIXmpmlVC351vB6PB489dRT\nuO222/CpT32qaod6x+NxcdZtMBj8/9u797io6vx/4K/hfhsQvCGMSoaoM1qijHKxtof2CFeL3cqH\nyj7SWqWbWWbfAqnWKCsyazNXrbZ0Q13FUotqS/FhimLmpeihAgYqpoypiDPcB5gz5/eHvzkxyGUY\nZoYZfD3/Q87MOUeRec1nPu/3G+7u7maTyqyZqc5wS87EmUJsR4/r6Drd2nlMR3tj2wqxltLr9Thx\n4gR8fX3h6+uLqqoqNDY2drtdYHvy8vKQmpqKzMxMTJ8+3e6/Lw0GA6KiorBnzx6Eh4dDrVZj8+bN\nUKlU0jFr1qzBiRMnpAKwHTt24LPPPoNWq8WUKVOQn58PLy8vTJ06FYsXL8b06dPtes3klFgARr2P\nKIooLS3FE088gfLycgwcOBCNjY1Qq9XS1gRnHsMLXC8qq6qqkorKjEajVFQWHBzMojJyKTdrkAWs\nD7NXrlzBmTNnMHLkSAQHB0t/3rJdoE6nQ11dHby8vKTfD9Z0VDEYDFi+fDny8/OxadMms32s9vbt\nt9/i2WefhSAImDdvHl566SUsXboUMTExSEpKgl6vx5w5c1BQUICQkBBkZ2dj2LBhAIBNmzYhMzMT\nMpkM06ZNw9tvv+2w6yanwjDbG3SnT19mZibWrVsHd3d3rFq1ComJiT1wB7b15Zdf4pVXXsGLL76I\nmTNnArg+HefQoUPIy8vDwYMHUVVVhejoaGlrwtChQ+1eTNEdgiBIRWVarRZNTU1mRWV+fn4c5kBO\nz1ahtrOP9Dsr8urKc7p7tB9W7RFkjUYjSkpKoNfroVQqLfqYX6/XS78fqqqqIIqitH2hT58+8PHx\naff3w6VLl5CSkoKJEyfitddes+pTIKIexjDr6rrTp6+oqAjJyck4cuQILl68iLvvvhslJSUu3aev\nuroaqampeP311zssknD1MbymwhHT6kx9fT38/PzMiso6C+emzg5VVVVQKpX4698LHXT1RNYHW0cG\nWWtWaoHrIdeaIFtXV4fCwkKEhoZi8ODBVn96JAgCqqurpd8Per0e/v7+cHd3x9WrVxEbGwtvb2/s\n3bsX6enpWL58OaZOnerUn1YRdYBh1tV1p0/fW2+9ZXZsy+NuNi3H8Obn56O4uNglxvCaiKKI+vp6\naZhDbW0tvLy8pH11QUFBZuG8qqoKxcXFCAsL6/BFk6u3ZGvdWZ21tgeso7YcmJi6FXTF77//jt9+\n+w1KpdLmRV6motPCwkK89957KC4uhru7OxoaGrBs2TLce++9CAkJsek5iRyIfWZdXXf69Gk0GsTG\nxpo9VqPROObCnYy7uzuio6MRHR2NRYsWmY3h/c9//oNffvkFISEh0sqtM4zhbUkmk8Hf3x/+/v5Q\nKBQAIFVFX7lyBaWlpVLldENDA/R6PcaMGdPpPbAlGNmSvfbOdnjOzvrKGv9YezEFW0cGWUEQcOrU\nKRiNRsTExFg1QbAzMpkMAQEBmDhxIt577z2kpKRg3LhxmDRpEn788Ud88sknqKurw7hx4zBr1izc\nfffdNr8GImfAMOukutOnz5LH3qzc3NwQFRWFqKgopKSkmI3h/eabb/DKK6/Aw8MDEydOREJCgkPH\n8FrKx8cHoaGhCA0NBQBcu3YNRUVF8PHxgZubG06cOGFWVGbJTHOGW+qO1gMSuhpuWwbPP57D+gEJ\nbT2/zE3W4cQxuLW//airQbampgZFRUVQKBQICwuz6+8PURSxZ88evPzyy3jnnXdwzz33AADuv/9+\nANe3XhUUFNjt/ETOgGHWSXW3T19nj6XrZDIZwsPDkZycjOTk5BvG8K5cuRINDQ2IiYlBfHy83cbw\nWsNoNOLs2bPQarWIjo6WVmNb7qsrKipCY2Mj5HK5tDXB39+fwxzIrtoKjbYIuNefp+v/9yzpO2ts\n53z/+3SMxecRRREajQYajQYqlcomkwg70tzcjNdffx0///wzdu7c2ebveW9vb7NP6oh6I+6ZdVLd\n6dNXWFiIv/3tb1IB2JQpU1BaWurUhU/OrOUY3vz8fGkMr2lrQnfG8FqruroaxcXFGDhwIIYOHdph\nODW1/DHtu62rq4Ovr6+0ciuXy626fgZcsoY999a2tWJrzQAFk64E2ebmZhQXF8PDwwMjRoyw++9b\njUaDlJQU3HXXXfjHP/5hl20MRD2MBWC9QXf69L3xxhtYv349PDw8sHLlSvz5z3/u4bvpPVqO4c3P\nz5fG8JragY0ZM8ZuLyxGoxFlZWW4du0aRo0aZdXKjyiK0jAHrVaLmpoaeHp6mg1zsOaFmOGWOtOT\nRWKtH2/JtVgaZk2FlxEREdL2H3sRRRG5ublYunQpVq5cyclY1JsxzBI5itFoRHFxMfLy8pCfn282\nhjc+Ph4xMTEW7VvtTE1NDYqLi9G/f3+b989tbGw0m1QGAEFBQVK/W2s6PjDcUkvO1u2gNbdWfWct\nCbKiKOK3335DRUUFVCqV3UdqNzU14bXXXsPJkyexYcMGuwdnoh7GMEvUU0wvcPv378eBAwfMxvAm\nJCQgNjYWcrnc4n23RqMR586dw9WrV6FUKu2+Dw+4vs3FFG51Oh0MBgMCAwPNiso4zIGs5YghC9b0\npTX5btPYTo9pampCYWEh/P39HbLV6MKFC0hJScE999yDF198kdvG6GbAMEvkLERRREVFBfLz85GX\nl4fDhw9DEARMmDBB6nc7YMCANsOhaTW2X79+iIiI6LFpZkajEdXV1dK+W71e3+058gy31JIt23t1\nNsGrI5YE2WvXruHXX39FZGQk+vfvb/W5LCGKIr799lssW7YMq1atwl133WXX8xE5EYZZImcliqI0\nhte077b1GN5Bgwbhtddeg0ajwerVqyGXy3v6ss20nCOv1WpRV1cHHx8faVtCYGAgi8qoW2w9ScyS\nldrOgqwoijh79ix0Oh1UKpVNtg91pKmpCa+88gpKSkqQlZWFAQMG2PV8RE6GYZbIlbQcw/vdd9+h\nuLgYt912G6ZPn44777zT6cfwiqIIvV4vrdxWV1fDw8ND2pYQFBRkVVEcwy211FnAtapt1/8PuZ0F\nWb1ej8LCQvTp0wfDhg2ze3u+3377DSkpKbj33nuRmprq1P//ieyEYZbI1QiCgHfffRc7duzAmjVr\n4ObmZjaGd+jQoVI7MGcfwwtcX1Uyrdy2LCozbU3w9vbu8nMy3N687BFkTXb+d1yH36+oqMDp06cx\nYsQIu4+HFUURX3/9Nd58802sXr0ad955p13PR+TEGGbJue3cuROLFi2CIAhISUnBkiVLzL7f2NiI\nuXPn4qeffkLfvn2xdetWREREYPfu3ViyZAmamprg5eWFFStWYPLkyT10F7YjCALuueceTJo0CS+9\n9NINQbXlGN4DBw5IY3jj4uIwadIkpxvD2xbTYA9TUVlzczPkcrm0NcHX17dLq12VlZVIXnDWjldM\nzqKrWw66Emw7CrJGoxGlpaWor6+HSqWy+xvIxsZGvPzyyzh37hw+/fRTu+/H7UmbNm3CqlWr0NTU\nhIkTJ2Lt2rVcfabWGGbJeQmCgKioKOzevRsKhQJqtRpbtmyBUqmUjlm7di2OHz8uDYT44osvsHXr\nVhQUFGDgwIEICwvDyZMnkZiYCI1G04N3YzsVFRUWv3iJoojff/8deXl5OHDgAI4ePSqN4TUVlTnb\nGN7WjEYjampqpK0JDQ0N8Pf3l7YmtFdUJgiCFDCUSuUN+xa5etv72Kv7QUdBtr6+HidPnsTAgQMx\nZMgQu/9fKisrQ0pKCu6//348//zzNiv2tHbhwOT8+fNQKpXIyMjA888/b5NrKi4uRmpqKnbs2AFP\nT08sWLAAsbGxmDt3rk2en3oNhllyXocOHUJGRgZ27doFAMjMzAQApKenS8ckJiYiIyMDcXFxMBgM\nCA0NRUVFhdkLiiiK6NevHy5evGjVR9a9SesxvD/88IM0hte0emvvGfHdJYoi6urqpK0JtbW18Pb2\nNisqq62tRXFxMcLDw6FQKCy6H4bb3scW4XbXlvHtfu/SpUs4d+4cRo0ahaCgoG6fqyOiKOLLL7/E\n22+/jbVr1yIhIcFmz92dhQOTBx98EG5ubpg4caLNwuzq1avx5ptvSgVtDQ0NSE5ORkZGhk2en3oN\ni1+wOP+OHE6j0WDw4MHS1wqFAocPH273GA8PDwQFBaGyshL9+vWTjtm+fTuio6Nv+iALADKZDCEh\nIbjvvvtw3333ATAfw/vf//7XbAxvfHw8hg8f3mOtvtoik8kQEBCAgIAAKBQKAJCKyn7//XecOHFC\nemPj6+sLg8EAT0/PTp83d6va7GuGW9cnGo03/FlXAm57QVYQBPz6668wGAwYP368RT9f3aHX65Ge\nno5Lly7h+++/R9++fW36/EeOHEFkZKQ0GXL27NnIyckxC7M5OTlSiJwxYwYWLlwIURQhk8nw5Zdf\nYtiwYTbfwiSKIh5++GFpIYOouxhmyeHa+jSg9QpbZ8cUFhYiLS0Nubm5tr/AXsLf3x+TJ0+W9hQ3\nNzejoKAAeXl5yMjIcOgYXmv5+PggMDAQ5eXlUCgUCA8PR3V1Na5du4aysjIYjUapqCw4ONiiNzYM\nt72TpQG3vSBbW1uLwsJChIeHIzw83O6fYpw+fRqPPvooZs2ahQ8++MAubyy7s3Dg6+uL5cuXY/fu\n3XjnnXdsel1TpkzBX/7yFyxevBgDBgzAtWvXUFNTg6FDh9r0PHTzcK5XLropKBQKXLhwQfq6vLwc\nYWFhbR6jUCikwiFTFXF5eTnuv/9+bNiwAbfeeqtDr92VeXp6YsKECZgwYQJeeOEFaQzv/v37sWbN\nGpw4cQKDBg2SVm5tNYbXWqIoory8HBqNxuzj3v79+0t7iwVBkIrKNBoNmpqazIrK/Pz8Og0lDLe9\nV+uAm/3BrRAEwazQSBRFaDQaaDQaqFQqu0/YE0UR27dvx7vvvouPPvoIsbGxdj1Xa5YuHLzyyitY\nvHixXf4+lEolXn/9ddxzzz0wGo3w9PTEmjVrGGbJagyz5HBqtRqlpaUoKytDeHg4srOzsXnzZrNj\nkpKSkJWVhbi4OGzbtg2TJ0+GTCaDTqfD9OnTkZmZadO9ZTcjNzc3qFQqqFQqPPnkk2ZjeD///HMs\nWbIEfn5+0sptV8fwdoder0dRURH8/f2hVqvbrXJ2d3dHSEiI9EbHVFSm0+lw+vRp1NfXw8/Pz6yo\nrLMVsNbhFmDA7Q02vD8YV65cQWlpKWQyGfr06QO5XI7Lly/D09MTMTExdq+mb2hoQFpaGq5du4a9\ne/favc1XdxYODh8+jG3btiE1NRU6nQ5ubm7w8fHBwoULbXJts2bNwqxZs2zyXEQsAKMe8e233+LZ\nZ5+FIAiYN28eXnrpJSxduhQxMTFISkqCXq/HnDlzUFBQgJCQEGRnZ2PYsGF4/fXXkZmZieHDh0vP\nlZuby8k4dtByDO/+/fvx448/QhAEqNVqJCQkdDiGtzsuXbqEsrIyREVFdXsPoSiKqK+vlzom1NbW\nwsvLS+p1GxQUZFWAYbh1La3foBgMBly8eBFlZWXw9PSEm5sb5HK59HNhyYp+V5WUlODRRx/FQw89\nhKefftoh+9UNBgOioqKwZ88ehIeHQ61WY/PmzVCpVNIxpk9lTAVgO3bswGeffWb2PBkZGQgICLBZ\nARiRhdjNgIhsq70xvGPHjpU6JgwdOtTqF+nm5macOnUKMpkMI0aMsFvxjV6vl3rdVlVVwc3NTQox\nffr0seq8DLfOq3WQFUUR58+fx+XLlzF69Gj4+fmZrejrdDppRT8oKAjBwcGQy+VW/1yLoojPPvsM\n77//Pj7++GOo1Teu/NuTtQsHLTHMUg9hmCUi+2s5hjc/Px/nz5/HqFGjEB8fj0mTJlk8hreyshIl\nJSUYNmwYBg4c6IAr/0Nzc7MUYnQ6HQRBMCsqs2bfMMOtc2gdZJuamlBYWAg/P78Ou3mIooiGhgbp\nZ6LleGbTir4lb3rq6+uRmpqK6upqfPLJJ+jTp49N7ovoJsEwS0SOJwgCjh8/3uYY3vj4eERHR5tN\nUaqpqcHevXsxePBgKJVKp2izJggCqqurpX63jY2NZh9B+/v7W/URNAOuY7UOslqtFqdOnUJkZKRV\nk7VM45lNK/pGoxGBgYHSz4Wvr6/Z8adOncJjjz2GRx55BAsWLHCqNnhELoJhloh6XkdjePv164cP\nP/wQjz76KBYsWOC0Ax1EUURtba2077aurg6+vr7Syq21H0Ez3NpPyyAriiLOnj0LrVaL0aNH26xD\nR8s3PTqdDitXrkRVVRUmTpwINzc35OTkYN26dRg/vv3hDETUIYZZInI+pv2Kzz33HA4fPoyIiAgY\nDAbExsYiPj4ecXFxCAkJcdpgC5h/BK3ValFTUwNPT0+zfbcsKusZrVdj9Xo9CgsLERQUhGHDhtl1\ndVQQBPzwww9YtWoVTp8+DU9PTygUCiQkJGDSpEmYMGGCzYcPEPVyDLNEPckZ56E7g19//RXz58/H\ntGnTkJqaCnd3d5cfwwtc//ds+RE0ALNw23JrhaUYbrumdZC9evUqSktLMWLECLu3wAKAoqIiPP74\n43j00Ufx2GOPwc3NDeXl5Th48CAOHjyIw4cP45NPPsGYMWPsfi1EvQTDLFFPcdZ56D1NFEUkJycj\nLS0N0dHR7R5XV1eHI0eOIC8vD/n5+bhy5QpGjx4ttQNztjG8bTEYDGZFZQaDAYGBgdIwBx8fH4sC\nem1tLYqKihAaGoqUFy474MpdU8sga9raUltbi9GjR1v1RqIrRFHExo0b8e9//xvr16/H2LFj7Xo+\nopsIwyxRTzl06BAyMjKwa9cuAJDmj6enp0vHJCYmIiMjA3FxcTAYDAgNDUVFRYU0D/3gwYPw9/dn\nOxz8MYbXtO/29OnTGD58OOLj4512DG9rRqMR1dXV0r5bvV6PgIAAad9t66Iy0/SzixcvQqlUQi6X\nt/m8XL01D7INDQ04efIk+vfvj6FDh9p9Rb+2thbPPvss3Nzc8MEHH7T770REVrH4P7BzvwIQuSBn\nnYfuqlqO4X3++efbHcNrmlTW02N429Kyly3wR1GZTqfD2bNnUVdXBx8fHynYlpeXw9fXt9OpVDfz\nKN7W93758mWUlZVh5MiRDmmBdfLkSTzxxBNYsGAB5s+f7/RbYYh6M4ZZIhtz1nnovUVHY3i3bduG\ntLQ0+Pv7Iy4uDvHx8YiNjUVgYKBThQ2ZTAa5XA65XI7BgwdDFEXo9XpcuHBBmkolCALKysoQHByM\noKAgi1afb5Zw2/I+BUFASUkJmpqaMH78eLsN2zAxGo3IysrC+vXrkZWVxT2wRE6AYZbIxpx5Hnpv\nJJPJEBERgYiICMydO9dsDG9eXh6WL1/ukDG83WHq8lBfX4/4+Hh4e3tLfU2vXr2KM2fOAIA0zKFP\nnz4W9eTtbeG29f3U1dXh5MmTCAsLg0KhsPu/aU1NDZ555hl4e3sjLy+PbzqJnAT3zBLZGOehO5fO\nxvAmJCQgIiKix4rKWhZ5DR48uN1AZnrTYyoqa25uhlwul4rKfH19e/Uwh9a9Y3///XecP38eKpXK\nIXtVjx8/jieffBLPPPMMHnnkEad6M0TUS7EAjKgncR66c+toDG9CQgJGjRplVa/YrrC0yKs9RqMR\nNTU1UlFZQ0MD/P39paKygICAXhNuWwZZg8GA4uJiuLm5WTwuuTuMRiPWr1+PDRs2ICsry+xNKRHZ\nFcMsEZGlujqGt7uamppQVFQEHx8fDB8+3CaBTBRF1NXVScMcamtr4e3tLa3cBgUFudykstbbCqqr\nq1FUVIShQ4di0KBBdj9/VVUVnn76aQQGBmL16tXw8/Oz+zmJSMIwS0RkrbbG8AYHB0srt92Z5lRZ\nWYmSkhJERkaif//+Nr5yc3q9Xlq5ra6ulroqmIrKrCmWclS4bb2t4MKFC7h06RJUKpVDJmkVFBTg\nqaeewnPPPYc5c+ZwWwGR4zHMEhHZimmPZl5eHg4cOICjR4/C3d0dEydOREJCgkVjeI1GI0pLS1Ff\nXw+lUmlRAZetNTc3Syu3VVVVMBqNUlFZcHCwVddkj3DbMsg2NzejsLDQpqvYHTEajfj444+xZcsW\nZGVlYdSoUXY9HxG1i2GWiMheRFGETqfDwYMHkZeXJ43hHT9+POLj428Yw3vs2DH8/PPPmDp1aodF\nXo4mCIJUVKbVatHc3IyAgABpa4Kfn59D99223lag0+lQXFyMW2+9FQMGDLDqObtCp9Nh4cKF6Nu3\nL1atWgVfX1+7n5OI2sUwS0TkSKYxvKatCVeuXIFKpYKbmxuOHTuGDz/8EBMnTuzpy+yQqajM1DGh\nvr4efn5+ZkVl9tp323pbwblz53D16lWMHj3aIaHy2LFjePrpp5GWlobk5GSbvuHYuXMnFi1aBEEQ\nkJKSgiVLlph9v7GxEXPnzsVPP/2Evn37YuvWrYiIiMDu3buxZMkSNDU1wcvLCytWrMDkyZNtdl1E\nTo5hloioJ2k0GiQnJwMAgoODpTG8pnZgt912m9OP4RVFEfX19dK+29raWnh5eUm9boOCgqz62L91\nuG0ZZBsbG1FYWIjAwEAMGzbM7i3TjEYjPvjgA2zfvh1ZWVkYMWKETZ9fEARERUVh9+7dUCgUUKvV\n2LJlC5RKpXTM2rVrcfz4calV3xdffIGtW7eioKAAAwcORFhYGE6ePInExERoNBqbXh+RE2OYJSL7\nsna1Cbjes/Pxxx+XipKOHj3qdCNou2Pnzp1IS0tDZmYmpk2bBgBmY3jz8/Nx/PhxaQxvfHw8YmJi\nXOJjbb1eL63cVlVVmY3q7dOnT7cmcJmK46KiotC3b18bXnXbtFotFixYgLCwMLz33nt2+Rk8dOgQ\nMjIysGvXLgBAZmYmACA9PV06JjExERkZGYiLi4PBYEBoaCgqKirMVodFUUS/fv1w8eLFHtlvTdQD\nLA6zzr0sQEROSRAEPPXUU2arTUlJSWarTevWrZNWJLOzs5GWloatW7fCYDDgoYcewsaNG3H77bej\nsrLS7iNIHclgMOB///sfcnNzMXDgQOnPOxrDu337dqSnp8PPz8+px/ACgI+PD0JDQxEaGgrgj6Iy\nnU6Hc+fOQRAEs6IySwKi0WjEmTNnUFNTg3HjxjkkrB05cgSLFi3Ciy++iJkzZ9rt71mj0WDw4MHS\n1wqFAocPH273GA8PDwQFBaGyshL9+vWTjtm+fTuio6MZZInawDBLRF125MgRREZGSoMeZs+ejZyc\nHLMwm5OTg4yMDADAjBkzsHDhQoiiiNzcXNx22224/fbbAcAhK3CO5OHhgX/961+dHtfWGN6rV6/i\nwIEDLjOGFwA8PT3Rv39/qc2YIAiorq6GTqdDUVERGhsbIZfLpZVbf39/s3toaGjAyZMn0a9fP0RH\nR9v9/oxGI1avXo2cnBxs27YNw4cPt+v52vr0s/U9dnZMYWEh0tLSkJuba/sLJOoFGGaJqMu6s9pU\nUlICmUyGxMREVFRUYPbs2UhNTXXo9TsjmUyG/v3744EHHsADDzxwwxjejz/+GFqtFtHR0U4xhrc9\n7u7uCA4ORnBwMG655RaIooja2lpotVqcPXsWdXV18PX1RZ8+fQBAmoBm+tqeKisr8eSTTyIiIgJ7\n9+51yNYWhUKBCxcuSF+Xl5cjLCyszWMUCoU0tjgkJEQ6/v7778eGDRtw66232v16iVwRwyyRi3N3\nd8eYMWOkr2fPnn3D/lVb685qk8FgQH5+Po4ePQo/Pz9MmTIF48ePx5QpU+x2va5IJpMhMDAQiYmJ\nSExMBGA+hjc1NRXnz5/HyJEjkZCQ4LAxvF0lk8kgl8shl8sxZMgQKdyW93u/ngAACqJJREFUlJSg\nvr4eHh4eOHPmjNm+W3vcw6FDh7B48WIsXboUDz74oMNWuNVqNUpLS1FWVobw8HBkZ2dj8+bNZsck\nJSUhKysLcXFx2LZtGyZPngyZTAadTofp06cjMzMTCQkJDrleIlfEMEvk4nx9ffHLL7849JzdWW1S\nKBT405/+JO0HnDZtGn7++WeGWQt4e3vjjjvuwB133AHAfAzvihUrUFxcjCFDhkiTymw9htcW6uvr\nUVxcjNDQUKnnbmNjI3Q6Ha5evYozZ84AgFm47c49CIKA999/H9999x2++OILh69uenh4YPXq1UhM\nTIQgCJg3bx5UKhWWLl2KmJgYJCUlYf78+ZgzZw4iIyMREhKC7OxsAMDq1atx+vRpLFu2DMuWLQMA\n5ObmOqTnLpErYTcDIhcXEBCA2tpah57TYDAgKioKe/bsQXh4ONRqNTZv3gyVSiUds2bNGpw4cUJq\nN7Rjxw589tln0Gq1mDJlCvLz8+Hl5YWpU6di8eLFmD59ukPvoTcyFVKZJpWZxvDGxcVh0qRJ3RrD\nawsXL17E+fPnoVQqERgY2O5xBoNBKirT6XQwGAwIDAyUhjn4+PhYtLJaUVGBJ554AiNGjMDy5ctZ\nPEXkWtiai+hm0XqbQXp6OmbNmmX383777bd49tlnpdWml156yWy1Sa/XY86cOSgoKJBWm0wFY5s2\nbUJmZiZkMhmmTZuGt99+2+7XezMyjeHdv38/9u/fbzaGNz4+HvHx8Z2O4bUFg8GAU6dOAQBGjhzZ\n5f66RqMR1dXVUr9bvV6PgIAAqWNC66IyAMjPz8fzzz+PV199FX/961+drnCOiDrFMEt0s+iJlVly\nTS3H8O7fvx8HDx6UxvCaVm/Dw8NtGvxqampQWFiIIUOG3LAVxVqmfbemMbxarRaZmZmIiYnBnXfe\niZ9++gn79u3Dxo0bpd7GRORyGGaJbhYMs9Qd9fX1OHz4sDSG9/LlyxgzZoy0chsVFWVVxwRRFFFe\nXo6LFy9i9OjRdt3eYDQaUVpaim+++QZff/21VBhn2l8cFxcHuVxut/MTkV0wzBLdLBhmyZaam5tR\nUFAgTSorLS3t8hje5uZmFBUVwcvLC1FRUQ7psLB//3688MILeOONN3Dffffh6tWrOHjwIA4cOIBD\nhw7By8sLe/fu5XYDItfBMEt0s2i9Z3bq1Kl46623evCKqDdpawxvaGiotHLbegzvnj17UFVVhYSE\nBLMJaPYiCAJWrFiBvLw8bNy4EUOGDGnzuMbGRhaAEbkWhlkiIrK9lmN4W/YLjouLw+XLl/HLL79g\nw4YNGDFihN2v5dKlS3jssccwbtw4vPHGG71qLDIRMcwSEZEDiKKI4uJiPPTQQ/Dx8YHRaIQgCIiJ\nicGkSZMQFxeHgQMH2vzj/X379iEtLQ1vvfUWpk2bxu0DRL0PwywRkS3s3LkTixYtgiAISElJuWG6\nWmNjI+bOnYuffvoJffv2xdatWxEREYHm5makpKTg559/hsFgwNy5c5Gent5Dd2E/e/bswf/93/9h\n+fLlSExMvGEM78GDB6HVajF27FipY0J3xvAaDAa89dZb+OGHH7Bp0yYoFAob3xEROQmGWSKi7hIE\nAVFRUdi9ezcUCgXUajW2bNkCpVIpHbN27VocP35cGg7xxRdfYOvWrdi8eTO++uorZGdno76+Hkql\nEvv27etVraJqamowb948rFq1CoMGDWr3uJZjePPz86VuA/Hx8Zg0aZLFY3gvXbqElJQUxMbG4tVX\nX+W2AqLezeIwy3G2RETtOHLkCCIjI6VhD7Nnz0ZOTo5ZmM3JyUFGRgYAYMaMGVi4cCFEUYRMJkNd\nXR0MBgMaGhrg5eXV4dQrVySXy/H55593elxbY3hPnDiBvLy8G8bwxsfHY9y4cWYjbEVRxPfff48X\nX3wRK1asQGJiIrcVEJGEYZaIqB0ajQaDBw+WvlYoFDh8+HC7x3h4eCAoKAiVlZWYMWMGcnJyMGjQ\nINTX1+O9995DSEiIQ6/fWbm7u2Ps2LEYO3YsFi1aZDaGNysrC4sXL5bG8MbGxmLfvn04fvw4du3a\nZbPBC0TUezDMEhG1o61tWK1XBNs75siRI3B3d8fFixeh1Wpxxx134O6775ZWeekPbm5uGD58OIYP\nH46UlBSzMbwbN25Ec3Mzdu3a1eUxuER0c+BvBiKidigUCly4cEH6ury8/IaVQdMxCoUCBoMBVVVV\nCAkJwebNmzF16lR4enpiwIABSEhIwLFjxxhmLSCTyRAWFobZs2dj9uzZPX05ROTkrCsnJSK6CajV\napSWlqKsrAxNTU3Izs5GUlKS2TFJSUnIysoCAGzbtg2TJ0+GTCbDkCFD8P3330MURdTV1eHHH3/E\nyJEje+I2iIh6NYZZIqJ2eHh4YPXq1UhMTMSoUaMwc+ZMqFQqLF26FF999RUAYP78+aisrERkZCT+\n+c9/StPXnnrqKdTW1mL06NFQq9X4+9//jttuu60nb4eIqFdiay4iIurVrO0VDACZmZlYt24d3N3d\nsWrVKiQmJvbAHRDdlCxuWcKVWSIi6rUEQcBTTz2F7777DkVFRdiyZQuKiorMjlm3bh2Cg4Nx+vRp\nLF68GGlpaQCAoqIiZGdno7CwEDt37sSCBQsgCEJP3AYRdYBhloiIeq2WvYK9vLykXsEt5eTk4OGH\nHwZwvVfwnj17IIoicnJyMHv2bHh7e+OWW25BZGQkjhw50hO3QUQdYJglIqJeq61ewRqNpt1jWvYK\ntuSxRNTzGGaJiKjX6k6vYEseS0Q9j2GWiIh6ra70CgZg1ivYkscSUc9jmCUiIgDXq/5HjBiByMhI\nqcVYS/v378e4cePg4eGBbdu2mX0vKytLmuJl6rvrDLrTKzgpKQnZ2dlobGxEWVkZSktLMWHChJ64\nDSLqACeAERGRVPW/e/duKBQKqNVqJCUlQalUSscMGTIEn376Kd555x2zx167dg2vvvoqjh07BplM\nhvHjxyMpKQnBwcGOvo0btOwVLAgC5s2bJ/UKjomJQVJSEubPn485c+YgMjISISEhyM7OBgCoVCrM\nnDkTSqUSHh4eWLNmDdzd3Xv4joioNfaZJSIiHDp0CBkZGdi1axeA6/1VASA9Pf2GYx955BHce++9\nmDFjBgBgy5Yt2LdvHz766CMAwOOPP4677roLycnJDrp6IuqF2GeWiIgs153KfVb9E1FPYpglIqJu\nVe6z6p+IehLDLBERdatyn1X/RNSTGGaJiMiiqv/2JCYmIjc3F1qtFlqtFrm5uUhMTLTzFRMRXccw\nS0REZlX/o0aNwsyZM6Wq/6+++goAcPToUSgUCnz++ed4/PHHoVKpAAAhISH4xz/+AbVaDbVajaVL\nlyIkJKQnb4eIbiLsZkBEREREzobdDIiIiIio92OYJSIiIiKXxTBLRERERC6LYZaIiIiIXBbDLBER\nERG5LIZZIiIiInJZDLNERERE5LIYZomIiIjIZTHMEhEREZHLYpglIiIiIpfFMEtERERELsuji8db\nPCeXiIiIiMjeuDJLRERERC6LYZaIiIiIXBbDLBERERG5LIZZIiIiInJZDLNERERE5LIYZomIiIjI\nZTHMEhEREZHLYpglIiIiIpfFMEtERERELothloiIiIhc1v8D8OEBO44SPlMAAAAASUVORK5CYII=\n", "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "# We now investigate how accurate this map is around the nominal point\n", "from mpl_toolkits.mplot3d import Axes3D\n", "import matplotlib.pyplot as plt\n", "from matplotlib import cm\n", "from matplotlib.ticker import LinearLocator, FormatStrFormatter\n", "import numpy as np\n", "\n", "N = 50\n", "M_grid = np.linspace(0,0.1,N)\n", "e_grid = np.linspace(0,0.1,N)\n", "# Make data.\n", "X, Y = np.meshgrid(M_grid, e_grid)\n", "Z = np.zeros([N,N])\n", "for i in range(N):\n", " for j in range(N):\n", " res = E_map.evaluate({\"dp0\": X[i,j], \"dp1\": Y[i,j]}) + En\n", " Z[i,j] = np.log10(np.abs(X[i,j] - res + Y[i,j]*sin(res)) + 1e-14)\n", "\n", "# Prepare the plot\n", "fig = plt.figure(figsize=(12, 9))\n", "ax = fig.gca(projection='3d')\n", "\n", "# Plot the surface.\n", "surf = ax.plot_surface(X, Y, Z, cmap=cm.coolwarm,\n", " linewidth=0, antialiased=False)\n", "ax.set_xlabel('E')\n", "ax.set_ylabel('e')\n", "ax.set_zlabel('Error (log)')\n", "plt.show()\n" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "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.6.5" } }, "nbformat": 4, "nbformat_minor": 2 }