{ "cells": [ { "cell_type": "markdown", "id": "ed5d30fe", "metadata": {}, "source": [ "# Steady-state torsional vibration in the crankshaft of an internal combustion engine\n", "\n", "Calculation of the steady-state torsional vibration in the crankshaft of an internal combustion engine. Torque produced by each cylinder is calculated from the force produced by the pressure from ignition which is scaled according to the rotation speed of the crankshaft. Torsional vibration analysis based on the article \"Analysis of torsional vibration in internal combustion engines: Modelling and experimental validation\" [1]. The model is based on Fig.4 of the original article." ] }, { "cell_type": "code", "execution_count": 1, "id": "6c0e2f6f", "metadata": {}, "outputs": [], "source": [ "import matplotlib.pyplot as plt\n", "import numpy as np\n", "import scipy\n", "\n", "import opentorsion\n", "from opentorsion import Shaft, Disk, Assembly" ] }, { "cell_type": "code", "execution_count": 2, "id": "b938a552", "metadata": {}, "outputs": [], "source": [ "def pressure_curve():\n", " '''Load digitized pressure curve from csv and pass it to interpolator'''\n", " curve = np.genfromtxt('ICE_example_data/pressure_data.csv', delimiter=';')\n", " \n", " return scipy.interpolate.interp1d(curve[:, 0], curve[:, 1])" ] }, { "cell_type": "code", "execution_count": 3, "id": "11a64081", "metadata": {}, "outputs": [], "source": [ "def peak_pressures():\n", " '''Load digitized peak pressure from csv and pass it to interpolator'''\n", " curve = np.genfromtxt('ICE_example_data/peak_data.csv', delimiter=';')\n", " \n", " return scipy.interpolate.interp1d(curve[:, 0], curve[:, 1])" ] }, { "cell_type": "code", "execution_count": 4, "id": "73d3a5a0", "metadata": {}, "outputs": [], "source": [ "def scaled_cylinder_pressure(rpm, num_points):\n", " '''\n", " Scales the cylinder pressure curve based on the variation \n", " of peak pressure with engine speed.\n", "\n", " Arguments\n", " ---------\n", " rpm: float\n", " Rotating speed of the engine in rpm\n", " num_points: int\n", " Number of points in x-axis\n", "\n", " Returns\n", " -------\n", " angles: np.ndarray\n", " Angles for which the cylinder pressures are given\n", " scaled_curve: np.ndarray\n", " The scaled pressure curve scaled to the given rpm corresponding\n", " to the angles. (Pa)\n", " '''\n", " # Load the base pressure curve\n", " angles = np.linspace(0, 720, num_points)\n", " base_curve = pressure_curve()\n", " base_curve_sampled = base_curve(angles)/15.2\n", " # Scale with the peak pressure\n", " pressures = peak_pressures()\n", " scaling = pressures(rpm)\n", " scaled_curve = base_curve_sampled*scaling*1e6\n", " \n", " return angles, scaled_curve" ] }, { "cell_type": "code", "execution_count": 5, "id": "cdc2e94e", "metadata": {}, "outputs": [], "source": [ "def calculate_cylinder_torque(speed_rpm, num_points=500):\n", " '''\n", " Calculates the torque produced by a single cylinder at given rpm.\n", " Torque is based on the tangential force caused by the pressure in the cylinder and the tangential component\n", " of the oscillating inertial force from the crankshaft spinning.\n", " \n", " Parameters\n", " ----------\n", " speed_rpm: float\n", " Current rotation speed in rpm\n", " num_points: float, optional\n", " Number of points for the x-axis\n", "\n", " Returns\n", " -------\n", " M_t: ndarray\n", " The torque produced by the cylinder\n", " alpha: ndarray\n", " Angles for which the cylinder torques are given\n", " '''\n", " def alpha_to_beta(alpha, r, l_rod):\n", " return np.arcsin(r/l_rod * np.sin(alpha))\n", " \n", " l_rod = 0.207 # m\n", " d_piston = 0.105 # m\n", " r = 0.137 / 2 # m (piston stroke / 2) crankshaft radius\n", " m_a = 2.521 # kg\n", "\n", " alpha_deg, p_curve = scaled_cylinder_pressure(speed_rpm, num_points)\n", " alpha = alpha_deg/180*np.pi\n", " w = speed_rpm/60*2*np.pi\n", "\n", " F_g = p_curve * 0.25*np.pi*d_piston**2\n", " beta = alpha_to_beta(alpha, r, l_rod)\n", " F_tg = F_g * np.sin(alpha + beta) * 1/np.cos(beta) # Tangential gas load\n", " lambda_rl = r/l_rod\n", " F_ia = -m_a*r*(np.cos(alpha)\n", " + lambda_rl*np.cos(2*alpha)\n", " - lambda_rl**3*1/4 * np.cos(4*alpha)\n", " + 9*lambda_rl**5*np.cos(6*alpha)/128)*w**2 # Oscillating inertial force\n", " F_ta = F_ia * np.sin(alpha + beta) * 1/np.cos(beta) # Tangential inertial force\n", " F_t = F_tg + F_ta # Total tangential force\n", " M_t = F_t * r\n", " \n", " return M_t, alpha" ] }, { "cell_type": "code", "execution_count": 6, "id": "f7ab717a", "metadata": {}, "outputs": [], "source": [ "def calculate_dft_components(signal, t, num_harmonics):\n", " '''\n", " Calculates dft components and harmonics (0,0.5,1,...) for the given signal, to be used at stedy-state\n", " vibration calculations.\n", "\n", " Parameters\n", " ----------\n", " signal: ndarray\n", " Cylinder torque\n", " t: ndarray\n", " Crankshaft rotation angle\n", " num_harmonics: float\n", " Number of harmonics considered\n", "\n", " Returns\n", " -------\n", " complex ndarray\n", " The first num_harmonics complex components of the Fourier transform\n", " complex ndarray\n", " The first num_harmonics components of the harmonics\n", " '''\n", " dft = np.fft.rfft(signal)/len(signal)\n", " dft[1:] *= 2\n", " omegas = np.fft.rfftfreq(len(signal))*1/(t[1]-t[0])*2*np.pi\n", " \n", " return [dft[:num_harmonics], omegas[:num_harmonics]]" ] }, { "cell_type": "code", "execution_count": 7, "id": "083a1a75", "metadata": {}, "outputs": [], "source": [ "def crankshaft_assembly():\n", " '''\n", " A shaft-line Finite Element Model of a crankshaft based on model presented\n", " in https://doi.org/10.1243/14644193JMBD126 Fig.4.\n", "\n", " Returns\n", " -------\n", " assembly: openTorsion Assembly class instance\n", " The created opentorsion assembly\n", " '''\n", " J2 = 0.0170\n", " J3 = 0.0090\n", " J4 = 0.0467\n", " J5 = 0.0327\n", " J6 = 0.0467\n", " J7 = 0.0467\n", " J8 = 0.0327\n", " J9 = 0.0487\n", " J10 = 2.0750\n", "\n", " k2 = 1.106e6\n", " k3 = 1.631e6\n", " k4 = 1.253e6\n", " k5 = 1.253e6\n", " k6 = 1.678e6\n", " k7 = 1.253e6\n", " k8 = 1.253e6\n", " k9 = 1.976e6\n", "\n", " c_a = 2 # absolute damping\n", " shafts, disks = [], []\n", " disks.append(Disk(0, J2))\n", " shafts.append(Shaft(0, 1, None, None, k=k2, I=0))\n", " disks.append(Disk(1, J3))\n", " shafts.append(Shaft(1, 2, None, None, k=k3, I=0))\n", " disks.append(Disk(2, J4, c=c_a)) # cylinder 1\n", " shafts.append(Shaft(2, 3, None, None, k=k4, I=0))\n", " disks.append(Disk(3, J5, c=c_a))\n", " shafts.append(Shaft(3, 4, None, None, k=k5, I=0))\n", " disks.append(Disk(4, J6, c=c_a))\n", " shafts.append(Shaft(4, 5, None, None, k=k6, I=0))\n", " disks.append(Disk(5, J7, c=c_a))\n", " shafts.append(Shaft(5, 6, None, None, k=k7, I=0))\n", " disks.append(Disk(6, J8, c=c_a))\n", " shafts.append(Shaft(6, 7, None, None, k=k8, I=0))\n", " disks.append(Disk(7, J9, c=c_a)) #cylinder 6\n", " shafts.append(Shaft(7, 8, None, None, k=k9, I=0))\n", " disks.append(Disk(8, J10)) # flywheel\n", " assembly = Assembly(shafts, disk_elements=disks)\n", " \n", " return assembly" ] }, { "cell_type": "code", "execution_count": 8, "id": "05931ab7", "metadata": {}, "outputs": [], "source": [ "def relative_damping_C(assembly, d, w):\n", " '''\n", " Updates the damping matrix C of assembly when using frequency dependent relative damping.\n", "\n", " Parameters\n", " ----------\n", " assembly: openTorsion Assembly class instance\n", " The assembly of whose damping matrix is to be updated\n", " d: float\n", " Loss factor property, used to calculate relative damping\n", " w: float\n", " Angular frequency of the system, used to calculate relative damping\n", "\n", " Returns\n", " -------\n", " C: ndarray\n", " The damping matrix assembled with new component specific damping coefficients\n", " '''\n", " if w!=0:\n", " c_r = d/w\n", " else:\n", " c_r = 0\n", " C = np.zeros((assembly.check_dof(), assembly.check_dof()))\n", "\n", " if assembly.shaft_elements is not None:\n", " for element in assembly.shaft_elements:\n", " dof = np.array([element.nl, element.nr])\n", " C[np.ix_(dof, dof)] += c_r*element.K()\n", "\n", " if assembly.disk_elements is not None:\n", " for element in assembly.disk_elements:\n", " C[element.node, element.node] += element.C()\n", "\n", " if assembly.gear_elements is not None:\n", " for element in assembly.gear_elements:\n", " C[element.node, element.node] += element.C()\n", "\n", " # Build transformation matrix\n", " E = assembly.E()\n", " transform = assembly.T(E)\n", " # Calculate transformed mass matrix\n", " C = np.dot(np.dot(transform.T, C), transform)\n", " \n", " return C" ] }, { "cell_type": "code", "execution_count": 9, "id": "dd60e398", "metadata": {}, "outputs": [], "source": [ "def calculate_response(crankshaft, rpm):\n", " '''\n", " Calculates the crankshaft's response to excitation at given rpm.\n", " \n", " Parameters\n", " ----------\n", " crankshaft: openTorsion Assembly class instance\n", " Opentorsion Assembly of the crankshaft\n", " rpm: float\n", " Current rotation speed in rpm\n", "\n", " Returns\n", " -------\n", " sum_response: ndarray\n", " Array containing maximum vibratory torque at current rotation speed for each shaft\n", " '''\n", " dof = 9\n", " cylinder_torque, alpha = calculate_cylinder_torque(rpm)\n", " dft_parameters, harmonics = calculate_dft_components(cylinder_torque, alpha, 25)\n", " q = np.zeros([dof, len(harmonics)], dtype='complex128')\n", " M = crankshaft.M\n", " K = crankshaft.K\n", " d = 0.035\n", " for i in range(len(harmonics)):\n", " # build T vector\n", " offset = 2 # offset to the first cylinder\n", " phase_shift = 2/3*np.pi\n", " T = np.zeros(dof, dtype='complex128')\n", " T[offset] = dft_parameters[i]\n", " T[offset+1] = dft_parameters[i]*np.exp( 2.0j*phase_shift*harmonics[i])\n", " T[offset+2] = dft_parameters[i]*np.exp(-2.0j*phase_shift*harmonics[i])\n", " T[offset+3] = dft_parameters[i]*np.exp( 1.0j*phase_shift*harmonics[i])\n", " T[offset+4] = dft_parameters[i]*np.exp(-1.0j*phase_shift*harmonics[i])\n", " T[offset+5] = dft_parameters[i]*np.exp( 3.0j*phase_shift*harmonics[i])\n", " w = harmonics[i]*rpm*2*np.pi/60\n", " C = relative_damping_C(crankshaft, d, w)\n", " receptance = np.linalg.inv(-w**2*M+1.0j*w*C+K)\n", " q.T[i] = receptance @ T\n", " # Calculate the angle difference between two consecutive disks\n", " q_difference = (q.T[:, 1:] - q.T[:, :-1]).T\n", " shaft_list = crankshaft.shaft_elements\n", " shaft_ks = np.array([shaft.k for shaft in shaft_list])\n", " # Multiply the angle difference between two disks with the connecting shaft stiffness to get\n", " # the torque in the shaft\n", " q_response = (shaft_ks*q_difference.T).T\n", " for i in range(dof-1):\n", " shaft_i = q_response[i]\n", " sum_wave = np.zeros_like(alpha)\n", " for i in range(len(shaft_i)):\n", " this_wave = np.real(shaft_i[i]*np.exp(1.0j*harmonics[i]*alpha))\n", " sum_wave += this_wave\n", " sum_response = np.sum(np.abs(q_response), axis=1)\n", " \n", " return sum_response" ] }, { "cell_type": "code", "execution_count": 10, "id": "6c41d2d1", "metadata": {}, "outputs": [], "source": [ "def plot_results(rpms, vibratory_torque):\n", " '''\n", " Plots the vibratory torque in wanted shafts for each considered engine speed.\n", " \n", " Parameters\n", " ----------\n", " rpms: ndarray\n", " All considered engine speeds in rpm\n", " vibratory_torque: ndarray\n", " Matrix containing maximum vibratory torque at each rotation speed for each shaft.\n", " Each row correspond to an engine speed and each column to a shaft.\n", "\n", " Returns\n", " -------\n", " shaft_8: ndarray\n", " Array containing maximum vibratory torque at all considered engine speeds for shaft \n", " between flywheel and 6th cylinder\n", " shaft_1: ndarray\n", " Array containing maximum vibratory torque at all considered engine speeds for shaft\n", " between crankshaft pulley and gear train\n", " '''\n", " shaft_8 = [shaft[7] for shaft in vibratory_torque]\n", " plt.plot(rpms, shaft_8, c='red', label='calculated')\n", " plt.xlabel('Engine speed (rpm)')\n", " plt.ylabel('Vibratory torque (Nm)')\n", " plt.title('Torque between the flywheel and the 6th cylinder (Fig.31)')\n", " plt.figure()\n", " shaft_1 = [shaft[0] for shaft in vibratory_torque]\n", " plt.plot(rpms, shaft_1, c='red', label='calculated')\n", " plt.xlabel('Engine speed (rpm)')\n", " plt.ylabel('Vibratory torque (Nm)')\n", " plt.title('Torque between the crankshaft pulley and the gear train (Fig.32)')\n", " plt.figure()\n", " \n", " return shaft_8, shaft_1" ] }, { "cell_type": "code", "execution_count": 11, "id": "6d0b86a7", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAZIAAAEWCAYAAABMoxE0AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuNSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/xnp5ZAAAACXBIWXMAAAsTAAALEwEAmpwYAABASUlEQVR4nO3dd5hU5dnH8e+PImAFlBil20GjKIg1SuwolhgLxoKGiBo1JnajRqMmscWoeROViIpRMcRKsWFDRRABAUVUQFBARJSiolLv94/nGRnWnd3Z3TlzZnfvz3XNNTPnnDnnnnrPearMDOecc666GqQdgHPOudrNE4lzzrka8UTinHOuRjyROOecqxFPJM4552rEE4lzzrka8URSRJJM0lZpx1Eskk6V9FqB9iVJ90paJGmspB6S5hRi31nHmCXpgELuM+73akkPFHq/OY5VpeeQ1HPO89j3Sbou4WN8/9pLaifpa0kNq7mvgr5WkppIelfSZpVs97SkPoU6br4kPSqpZz7b1tpEEj8QmctqSd9m3T8x7fgKrdSTkKQOMcZGCR1ib+BAoI2ZdU/oGHVakX64z5M0U9JSSVMlbROXF+xPRXWZ2cdmtr6ZrUozjiz9gFfMbB58//4sL/PbdryZ9TSzgdU5gKSXJC2Q9KWkSZKOzFq3maQhkj6J390OZR5+A5DX56XWJpL4gVjfzNYHPgYOz1r2YD77SPBHzxVee2CWmS1NOxBXPkm/BvoChwHrA72Az1MNqgRU8DtzJvCfMstuzP5tM7P/1vDw5wGbmdmGhMT1QNYZ0GrgGeAX5T3QzMYCG0rqVtlBam0iySWeLt4as+wn8XaTuK6HpDmSLpH0KXCvpGbxn8CieJp5UXaRSdkzgbL/6iT1kjRR0mJJr0vasZIQD5X0oaTPJd0k6fv3QNKv4r+4RZKeldQ+Ln8lbjIp8y9F0khJv4jr94pxHhbv7y9pYmX7jeu2kzRC0kJJ70s6rsxz/aek4ZK+kvSGpC1zPK9MjItjjHtk7efmeOyZ2afKkjaSNEDSPElzJV1XXrGDpL7A3cAecd9/KrP+IkmPlll2u6TbJP1M0ttZy0dIejPr/quSjsp6aBdJkyUtkfRfSU2zts35XkvaXKEoYEF8nr/N8TqVfW4tJA2Lj1sUb7fJWv+ypGsljYrvwXOSNslaf7KkjyR9IenyCo7TDzgRuDi+hkNr+pzL7L8BcBXwezN714IZZrZQUifgTta8f4uzHtoiz88XkvaOMSyWNFvhLGdXSfOzPzeSjpY0qZzHr3XWXJPXVlIDSZdKmhHXD5bUssxx+kr6GHixnFjaAVsAb+R6vlnbvqyQpJHUUNLfFH4/Zko6RxWUBJjZZDNbmbkLNAbaxnXzzexfwJvlPTZ6mfDHoGJmVusvwCzggHj7GmAM8COgFfA6cG1c1wNYSThlawI0A64HXgVaxhf4HWBO1r4N2Crr/n3AdfH2zsBnwG5AQ6BPjKVJjjgNeCkeqx3wAfDruO5IYDrQCWgEXAG8XkEc1wD/iLf/AMwAbshad1tl+wXWA2YDp8V1OxP+QXbOeq5fAN3j+geBh3M8tw4xxkZZy04FVgCnx9fnLOATQHH948BdMY4fAWOBM3Ls/1Tgtaz7PTLvE7AZsBRoHu83iu9L1/gefwdsQvgSzQfmAhvEdd8CG2d9jsYCm8f3aCpwZmXvNeEP2Xjgj8A6hB+ID4GD42OvBh7I8bw2JvwjXDfG9D/giaz1L8f3dpsY78vA9XFdZ+BrYJ8Yxy2Ez/cBOY51H/GzW+a7U+XnXM6+28X3/zzCZ2om8CegQXnvXzU+X+2Br4AT4vu4MdAlrnsX6Jm17ePABWVfe8p8Rmvy2sbnOQZoE9ffBQwqc5z7CZ/tZuU8n8OAKZW9P1lxZn4nzozPtw3QAnieMt+7ch4/jPAdMMIZSIMy6xvFdR3Keez5wGOV/gZX98e7lC6snUhmAIdmrTuYUCQC4cdnOdA0a/2HwCFZ9/uRfyK5g5iksta/D+ybI04rc6zfAC/E208DfbPWNQC+AdrniGN/YHK8/Qzwa2BMvD8SOLqy/QLHA6+WifEu4Kqs53p31rpDgfdyPLfMl6dsIpmedX/duM2PgU2BZdlfMsKPxEs59n8qORJJ1vM8Pd7uBbybte5V4Ghgd+A5YDBwCPCzzGuY9Tk6Kev+jcCdlb3XhB/aj8usuwy4N96+mhyJpJzn2QVYlHX/ZeCKMp+ZZ+LtP5L1w0v40VpO1RNJlZ9zOfveM763w4Hm8fPwQdZ7stb7V43P12XA4znWXQI8GG+3JHy+Nyv72lN+IqnWa0tIuPtnrd+M8KepUdZxtqjgfT6R+H0t83p8ByyOl8+z4swkkhfJ+rMFHEAliSRu1xjoCZxfzrqKEsnpwIuVfW7rXNEW4Z/VR1n3P4rLMhaY2Xdltp9dZvt8tQcuiKfai+Mpe9syxyur7LEy27YHbsvaz0JAQOsc+xkNbCNpU8KPz/1A23hq3p01RU0V7bc9sFuZ+E8k/NBnfJp1+xtC2XdVfP94M/sm3lw/HrsxMC/r2HcRzkyqYyBwUrx9EmuXPY8kJJ594u2XCQlg33i/3HhZ+/lW9F63BzYvs+4PhGRZIUnrSrorFqF8SXjfmmvtIr5cMa312bVQf/RFZccsR3Wec1nfxusbzWyxmc0ivJ+HVvPYZbUl/EkszwPA4ZLWA44j/DmaV8lxKzt+Za9te+DxrNdlKrCKtd/z7O96WYsIZ6Bl3WxmzeNlk3LWl/29qugY3zOzFWb2NHCQpCPyeUy0ASGpVaguJpJPCG9yRru4LMPKbD+PWGaYtX22bwj/pDOyf2RnA3/OeuObm9m6ZjaogvjKHisT22zCP43sfTUzs9fL20n8UR5POMV+x8yWE4rxzgdmmNnneex3NjCyzLr1zeysCuLPpezrWpnZhDOSTbKOvaGZbV+NYwM8AewoaQfCGUl2g4uyiWQkuRNJRfHmeq9nAzPLrNvAzCr7EQW4ANgW2M1Cheg+cbnyeOxan11J6xKKfHKpznuU7+f7fcI/9uxj5LpdHbOBcutPzGwu4Y/V0cDJ/LACuzoqe21nE4rTsl+bpjGW70OrYP+TgY656jYqiatN1v22uTbMoRE5XsccOgE/qG8qqy4mkkHAFZJaxX/nfyT8Y8llMHBZrPRsA5xbZv1E4JexkusQwo9Pxr+BMyXtpmA9SYdJKu+fRsZF8VhtCUkg0yrjzhjH9vB9RfSxWY+bTyh7zzYSOIc1P4Yvl7lf2X6HEc5qTpbUOF52jZWjVbWA0AqkbIzliv8YnwP+JmnDWHm5paR9K3tsjv19BzwCPASMNbOPs1a/Tvix7h7XTSGejbHmzK0yFb3XY4GvFBpxNIuflR0k7ZrHfjcg/JtfHCtrr8ozHgjPt5dCJfQ6hLqxir7T5X2GKpL35zv+sfkvoTJ/g/hd6kf4jGWO3SbGWR0PAgdIOk5SI0kbS+qStf5+4GLgJ8Bj1TxGtspe2zuBP2tNg5hWympaWxkzm0Oou6xqU/bBwHmSWktqTijWK5dCQ5qe8TPZWNJJrPkzldmmKaGOB6CJshpaRPsSio0rVBcTyXXAOELGfxuYQMVtof9EKGKaSfhhK/tv5jzgcMLp3YmEf74AmNk4Qhni/xFOVacTyoIr8iThTGIioTx5QNzX44RGAA/HIo53CGWaGVcDA+OpdKZl1UjCD9ErOe5XuF8z+wo4COhNODP6lDUNEaok/pD8GRgVY9w9j4edQqicfpfw+j1CKGuuroGEH5K13sNYLDGBULm5PC4eDXxkZp/ls+OK3msL/RJ6EYoYZxIaLNwNbJTHrm8lVPR+Tqi8fSafeOJxpwBnE5LnvBhXRZ00BwCd4/vzRB77r+rn+xxCBfUnhNf3IeCeuO5FYArwqaQqNwmOfwwOJZzBLSR8f3bK2uRxYnFTVhFqteXx2t4GDAGek/QV4b3brYqHuYtwBlUV/yb8Tk0G3gKeIjQCWAUg6U5Jd8ZtRfjd+IzwR+884Hgzm5C1v28J7xnAe6wpoiT+EfraQjPgCmVaz7hIUg9C5VybSjZ1JUahSeV7wI/N7Mu043HFJWkGoRj3+bRjyYdCt4S3CJX2+dbplN1HT0LjiPaVblz1fT8KDDCzpyrb1jvkuTpBoR/D+YSWNp5E6hmFPlVGOX02SpWZLSM0M86bpGaE1obPESr2ryKcjRWcmZXbUbE8nkhcrRdb68wnFFEeknI4rsgkvUz4QT7ZzFanHE7SRCiO/y+hGGo4oR44VV605ZxzrkbqYmW7c865IqqTRVubbLKJdejQIe0wnHOuVhk/fvznZtaqqo+rk4mkQ4cOjBs3Lu0wnHOuVpFUlZE9vudFW84552rEE4lzzrka8UTinHOuRjyROOecqxFPJM4552rEE4lzzrka8UTinHOuRjyROOdK3+zZ8OSTaUfhcqiTHRKdc3XIkiVwwAHwwQcwfTpsWZUJ/lwx+BmJc650rV4NJ50EM+J07U88kWo4rnyeSJxzpetPf4Jhw+Dvf4cuXTyRlChPJM650vTEE3DNNdCnD5xzDhx1FIwaBZ/lNTuyKyJPJM650jN1Kpx8MnTrBnfeCVJIJGYwdGja0bkyPJE450rLkiUhaay7Ljz2GDRtGpbvuCN06ACPJzKzrKsBTyTOudJy772hhdbgwdC27ZrlEvz85/D88/DVV+nF537AE4lzrrS8+mo489h33x+uO+ooWLYMnn222FG5Cngicc6VDrNQob733uWv32sv2GQTL94qMZ5InHOlY8YMmD8/dyJp2BCOOAKGD4fly4sbm8vJE4lzrnSMGhWu99or9zZHHRUq5EeOLEpIrnKeSJxzpeO116B5c+jcOfc2BxwQWnR558SSkWgikXSepHckTZH0u7ispaQRkqbF6xZxuSTdLmm6pMmSdsnaT5+4/TRJfZKM2TmXolGjYM89oUEFP03NmsEhh4REsnp10UJzuSWWSCTtAJwOdAd2AnpJ2gq4FHjBzLYGXoj3AXoCW8dLP+COuJ+WwFXAbnFfV2WSj3OuDvnii9ARMVf9SLaf/xw++QTGjUs+LlepJM9IOgFvmNk3ZrYSGAkcDRwJDIzbDASOirePBO63YAzQXNJmwMHACDNbaGaLgBHAIQnG7ZxLw+uvh+uK6kcyDjssVLx78VZJSDKRvAP8VNLGktYFDgXaApua2by4zafApvF2a2B21uPnxGW5lq9FUj9J4ySNW7BgQWGfiXMuea+9Bo0bw667Vr5tixbQo4cnkhKRWCIxs6nADcBzwDPARGBVmW0MsAIdr7+ZdTOzbq1atSrELp1zxTRqFHTtGupA8rHvvqEobOnSZONylUq0st3MBphZVzPbB1gEfADMj0VWxOvMUJ5zCWcsGW3islzLnXN1xXffwZtv5lc/ktGpU7h+//1kYnJ5S7rV1o/idTtC/chDwBAg0/KqD5CZP3MIcEpsvbU7sCQWgT0LHCSpRaxkPyguc87VFePHhw6G+dSPZGQSydSpycTk8pb0VLuPStoYWAGcbWaLJV0PDJbUF/gIOC5u+xShHmU68A1wGoCZLZR0LfBm3O4aM1uYcNzOuWLKpyNiWVtvHSrcPZGkLtFEYmY/LWfZF8D+5Sw34Owc+7kHuKfgATrnSsNrr8E220BV6jfXWSfM3+6JJHXes905l67Vq0PT36rUj2R06uSJpAR4InHOpev990NnxKoUa2V06gTTp8OKFYWPy+XNE4lzLl2Z+pHqnJFst11IIh9+WNiYXJV4InHOpeu110LdyNZbV/2x3nKrJHgicc6la9SoUKwlVf2x220Xrj2RpMoTiXMuPfPnhzqO6tSPAGy4IbRu7YkkZZ5InHPpqU7/kbK85VbqPJE459IzenToD7LLLpVvm0unTvDee2G+d5cKTyTOufSMHh0GamzSpPr76NQJvv4a5swpXFyuSjyROOfSsXx5mJhqjz1qtp9Mhft779U8Jlctnkicc+mYNAmWLat5IvEmwKnzROKcS8fo0eF6991rtp9NN4XmzT2RpMgTiXMuHaNHQ5s24VITkrfcSpknEudcOkaPrnmxVoYnklR5InHOFd+8efDRR4VNJJ99Bgt9qqI0eCJxzhVfpn6kkIkE/KwkJZ5InHPFl+mIuPPOhdmfj7mVqqTnbP+9pCmS3pE0SFJTSR0lvSFpuqT/Slonbtsk3p8e13fI2s9lcfn7kg5OMmbnXBGMGRN6s9ekI2K2Dh3CvrwvSSoSSySSWgO/BbqZ2Q5AQ6A3cAPwdzPbClgE9I0P6Qssisv/HrdDUuf4uO2BQ4B/SWqYVNzOuYQVqiNitoYNYdtt/YwkJUkXbTUCmklqBKwLzAP2Ax6J6wcCR8XbR8b7xPX7S1Jc/rCZLTOzmcB0oHvCcTvnkjJpEnz3XWETCXjLrRQllkjMbC5wM/AxIYEsAcYDi81sZdxsDtA63m4NzI6PXRm33zh7eTmP+Z6kfpLGSRq3YMGCwj8h51xhFLqiPaNTJ5g1C779trD7dZVKsmirBeFsoiOwObAeoWgqEWbW38y6mVm3Vq1aJXUY51xNFaojYlmdOoURgN9/v7D7dZVKsmjrAGCmmS0wsxXAY8BeQPNY1AXQBpgbb88F2gLE9RsBX2QvL+cxzrnaZvTomg+LUh5vApyaJBPJx8DuktaNdR37A+8CLwHHxG36AE/G20PifeL6F83M4vLesVVXR2BrYGyCcTvnklLojojZtt4aGjTwRJKCRpVvUj1m9oakR4AJwErgLaA/MBx4WNJ1cdmA+JABwH8kTQcWElpqYWZTJA0mJKGVwNlmtiqpuJ1zCRozJlwnkUiaNoWOHT2RpCCxRAJgZlcBV5VZ/CHltLoys++AY3Ps58/AnwseoHOuuAoxI2JFOnWCd99NZt8upwqLtiTtIemfkiZLWiDpY0lPSTpb0kbFCtI5V0eMHl3Yjohl7bhjqGz/7rtk9u/KlTORSHoa+DXwLKG11WZAZ+AKoCnwpKQjihGkc64OSKIjYlldusCqVX5WUmQVFW2dbGafl1n2NaHOYwLwN0mbJBaZc65umTAhnCkk0WIrY6edwvXEickVn7kfyJlIyiYRSRtmb29mC8tJNM45V75hw8JQJgcckNwxttwS1lsv9J53RVNpZbukM4A/Ad8BFhcbsEWCcTnn6pqhQ2GvvaBly+SO0bAh/OQn4YzEFU0+/UguBHYwsw5m1jFePIk45/L30UcweTIcfnjyx+rSJZyRmFW6qSuMfBLJDOCbpANxztVhw4aF62Ikkp12giVLQvJyRZFPP5LLgNclvQEsyyw0s98mFpVzrm4ZOjT0PN922+SP1aVLuJ40KcxT4hKXzxnJXcCLwBjC6L2Zi3POVe6rr+Cll4pzNgKhjkTyepIiyueMpLGZnZ94JM65uum550IfkiOK1O1svfXC2Y+33CqafM5Ino5zfWwmqWXmknhkzrm6YehQaNEitNgqlp128jOSIsrnjOSEeH1Z1jJv/uucq9yqVTB8OPTsCY0SHdpvbV26wP/+B19+CRtuWLzj1lOVvrNm1rEYgTjn6qA33oDPPy9e/UhGpof75Mmw997FPXY9lDORSNqngseZmb2aQDzOubpk6NBwJnJIYpOjli/TcmviRE8kRVDRGclF5SwzYEfCjIUNE4nIOVd3DB0KP/0pNG9e3ONuvjlsvLFXuBdJRWNtrXUuKmkvwsi/nwLnJhyXc662+/BDmDIFfv3r4h9bCmclXuFeFJW22pK0v6SXgeuAW8xsdzMbmnhkzrnabWj8mSh2/UjGTjvBO+/AypXpHL8eqaiO5DDgcmAJcIWZvVa0qJxztd/QoWHGwi23TOf4XbqEYes/+AA6d04nhnqiojOSoUAbwjzpF0sakn2pbMeStpU0MevypaTfxX4oIyRNi9ct4vaSdLuk6XFGxl2y9tUnbj9NUp+aPmnnXMIWLYKRI9M7G4G1h0pxiaqosv1nNdmxmb0PdAGQ1BCYCzwOXAq8YGbXS7o03r8E6AlsHS+7AXcAu8XOj1cB3QiV/eMlDTGzRTWJzzmXoP79Q5HSL3+ZXgzbbRfmh584EU44odLNXfVVVNk+soDH2R+YYWYfSToS6BGXDwReJiSSI4H7zcyAMZKaS9osbjvCzBYCSBpBmPp3UAHjc84VyvLlcNttcOCBa/pzpKFxY9h+ez8jKYKK5mwfKulwSY3LWbeFpGsk/SrP4/RmzQ//pmY2L97+FNg03m4NzM56zJy4LNfysjH1kzRO0rgFCxbkGZZzruAGDYJ58+CCC9KOxIdKKZKK6khOB34KvCfpTUlPSXpR0oeEEYHHm9k9lR1A0jrAEcD/yq6LZx8FmX3GzPqbWTcz69aqVatC7NI5V1VmcPPNsMMOcNBBaUcT6knmz4dPP007kjqtoqKtT4GLCRXtHYDNgG+BD8ysKhNd9QQmmNn8eH++pM3MbF4suvosLp9L6OiY0SYum8uaorDM8percHznXLE891xocnvffaEvR9oyRWuTJsGPf5xuLHVYPqP/YmazzGy0mU2sYhKBMOhjdn3GECDT8qoP8GTW8lNi663dgSWxCOxZ4CBJLWILr4PiMudcqbn55tCrvFQqt7MTiUtMosNxSloPOBA4I2vx9cBgSX2Bj4Dj4vKngEOB6YSpfU8DMLOFkq4F3ozbXZOpeHfOlZCJE+H55+Gvfw2tpUpBixbQrp3XkyQs0URiZkuBjcss+4LQiqvstgacnWM/9wCV1sc451J0yy1hUqkzzqh822Lq1g3GjEk7ijotr6ItSc0kFWGyZedcrTRnTmit9etfh7OAUrLvvjBzJnz8cdqR1Fn5jLV1ODAReCbe75JPz3bnXD1y++2wejX87ndpR/JDPXqE65GF7BrnsuVzRnI10B1YDGBmEwGf7Mo5F8ybB3fdBcceCx06pB3ND+2wA7RsCS+/nHYkdVY+iWSFmS0ps6wgfT+cc7WcGfTrF3qzX3112tGUr0GDMCeKn5EkJp9EMkXSL4GGkraW9A/g9YTjcs7VBgMHwrBh8Je/hLGtSlWPHjBjRqjLcQWXTyI5F9geWEboD/Il8LsEY3LO1QazZ8N558E++4TrUrbvvuHaz0oSUWkiMbNvzOxyM9s1DkFyuZl9V4zgnHMlygz69oVVq+Dee0PxUSnbcccw3a8nkkRU2o9E0kuUUydiZvslEpFzrvTddReMGAF33AFbbJF2NJVr2DDUk3iFeyLy6ZB4YdbtpsAvCJNdOefqow8/hAsvDMPEl1rnw4r06BFmbfzkkzCMiyuYShOJmY0vs2iUpLEJxeOcK2WrVsFpp4V/+AMGlMbAjPnKricplbHA6oh8OiS2zLpsIulgYKMixOacKzW33AKvvBI6ILZtW/n2paRLF9hwQ68nSUA+RVvjCXUkIhRpzQT6JhmUc64ETZoEl18ORx8Np5ySdjRV5/UkicmnaMt7sTtX3333HZx8Mmy8cahor01FWtl69IDhw8NEVz4/ScHk02rr6IrWm9ljhQvHOVeSrrwS3n4bnnoKNtkk7WiqL7ue5Pjj042lDsmnaKsvsCfwYrz/M0LP9gWEIi9PJM7VZS+/DH/7G5x1FvTsmXY0NbPzzrDBBp5ICiyfRNIY6BxnKyROj3ufmZ2WaGTOufQtWQJ9+sBWW8FNN6UdTc01agR77+31JAWWT3fUtpkkEs0H2iUUj3OulPzudzB3LjzwQJi0qi7o0QOmToXPPks7kjojn0TygqRnJZ0q6VRgOPB8smE551I3YQLcdx9cdBF07552NIXj424VXD5jbZ0D3AnsFC/9zezcfHYuqbmkRyS9J2mqpD1if5QRkqbF6xZxW0m6XdJ0SZMl7ZK1nz5x+2mS+lTvqTrnquTyy8M8HpdemnYkhbXLLrD++vDSS2lHUmdUmEgkNZT0npk9bma/j5fHq7D/24BnzGw7QhKaClwKvGBmWwMvxPsAPYGt46UfcEeMoSVwFbAbYYKtqzLJxzmXkFdegWeeCUlkozrW/7hx4zC8y+OPh576rsYqTCRmtgp4X1KV60QkbQTsAwyI+1puZouBI4GBcbOBwFHx9pHA/RaMAZrHiv2DgRFmttDMFgEjgEOqGo9zLk9mcNllYTyqc85JO5pknHBC6EvixVsFkU+rrRaEya3GAkszC83siEoe15HQRPheSTsResifB2yaVXn/KbBpvN0amJ31+DlxWa7la5HUj3AmQ7t23hbAuWobPhxefx3uvBOaNUs7mmT06hWKtwYNgv18IPOayieRXFmDfe8CnGtmb0i6jTXFWACYmUkqyLS9ZtYf6A/QrVs3nwrYuepYvTrUjWy1FfzqV2lHk5xmzeCoo+DRR+Gf/4R11kk7olotn8r2kcB7wAbxMjUuq8wcYI6ZvRHvP0JILPNjkVWmT0qmDd5cIHsUuDZxWa7lzrlCe/hhmDwZrrkm1CXUZSecAIsWwbPPph1JrZfP6L/HAWOBY4HjgDckHVPZ48zsU2C2pG3jov2Bd4EhQKblVR/gyXh7CHBKbL21O7AkFoE9CxwkqUWsZD8oLnPOFdLy5WEolJ12qh+9vg88MIwdNmhQ2pHUevkUbV0O7GpmnwFIakXoR/JIHo89F3hQ0jrAh8BphOQ1WFJf4CNCcgJ4CjgUmA58E7fFzBZKuhZ4M253jZktzOPYzrmqGDAgTFo1bFjpT51bCI0bwzHHwH/+A0uX1p0OlymQWcXVCZLeNrOfZN1vAEzKXlZqunXrZuPGjUs7DOdqj6VLQ73IllvCq6/W3tF9q2rkyNDTfdAg6N077WhSJ2m8mXWr6uPy+dvxTDk925+u6oGccyXs1ltDc9gbb6w/SQTC/CStW3vxVg3lU9l+EXAXsGO89Dezi5MOzDlXJAsWwA03hFZMe+6ZdjTF1aBBOBN5+ulQ8e6qJZ/K9hvM7DEzOz9eHpd0QzGCc84VwZ//HIq2/vKXtCNJxwknwIoV8JjPiFFd+RRtHVjOslo+KYFzDgiV6//6V+gz0qlT2tGkY5ddYOutvXirBnImEklnSXob2DYOopi5zAQmFy9E51xirrwyzNFx9dVpR5IeKZyVvPRSqCdyVVbRGclDwOGE/h2HZ126mtlJRYjNOZekCRPgoYfCnCOtfzDqUP1ywgmhV//DD6cdSa1UafPf2sib/zqXh4MOgvHjQ/FWXRvhtzr22COckXzwQd3v1Z9Dks1/nXN1zYgR4XLFFZ5EMi6/HGbNggcfTDuSWsfPSJyrb5YuDRXMy5bB++9DkyZpR1QazKBrV/j66zAVb8OGaUdUdImdkUg61yeScq4OueiiUHxzzz2eRLJJ4Qxt2jT473/TjqZWyadoa1PgTUmDJR0i1adur87VMcOHwx13wAUX+Dwc5TnqKNh++9C3ZvXqtKOpNfLp2X4FYfrbAcCpwDRJf5G0ZcKxOecK6bPPQn+Rn/wk/FC6H2rQINSVvPtumIrX5SWvynYLFSmfxstKwqyJj0i6McHYnHOFYgannw6LF4fKZC/Syu2442CbbeC668Lr5iqVTx3JeZLGAzcCo4CfmNlZQFfgFwnH55wrhAEDYMgQuP76cEbicmvYEP7wB5g4MQyp7yqVzzDyVwP3mtlH5azrZGZTE4qt2rzVlnNZpk+HLl1g993huefqx1wjNbViBWy7LbRqBWPG1JsRkRNptSWpIdC7vCQCUIpJxDmXZdYsOPjgMCf5ffd5EslX48Zw2WUwdmxIvq5CFX6qzGwV8L6kdkWKxzlXKB98EObbWLgQnnkG2rRJO6LapU8f6NgRzj8/nKG4nPL5e9ICmCLpBUlDMpd8di5plqS3JU2UNC4uaylphKRp8bpFXC5Jt0uaHgeH3CVrP33i9tMk9cl1POdc9M47sM8+odPhyy9D9+5pR1T7rLMO3HZbaMF1661pR1PS8qkj2be85WY2stKdS7OAbmb2edayG4GFZna9pEuBFmZ2iaRDCXO8HwrsBtxmZrtJagmMA7oBBownDByZcxYaryNx9dqECWEcrSZN4Pnn6+/w8IVyxBHw4ovw3nt1/qwusZ7tMWG8B2wQL1PzSSIVOBIYGG8PBI7KWn6/BWOA5pI2Aw4GRpjZwpg8RgCH1OD4ziVv2jT429/C0OTFYgZPPBE6Gq6/PrzyiieRQrjtNli1Cn7/+7QjKVn5NP89DhgLHAscB7wh6Zg892/Ac5LGS+oXl21qZvPi7U8JPecBWgOzsx47Jy7Ltdy50jJ3LtxyC+y6a+iHcOGF4Uf95JNDZ8AkjRkTirJ+/vPwr/mVV2BL7zNcEB07hk6KjzziFe855FNHcjmwq5n1MbNTgO7AlXnuf28z24Uwo+LZkvbJXhk7Ohakx4+kfpLGSRq3YMGCQuzSufysXg2nnAJt24ahR8zg5pvDWcmVV4Zxmzp1CmNbFbqD2wcfwDHHhCHQp02DO+8M/R/aefuYgrroojCL4jnnhHont5Z8EkkDM8v+O/VFno/DzObG68+AxwlJaH4ssiJeZ/Y9F2ib9fA2cVmu5WWP1d/MuplZt1atWuUTnnOFcfPN8J//wLnnhtF0x40LCWWrreCaa8IPe+fO0Lcv9OgRKsJrwgxGjgyTMW2/fWiR9ac/hf4iZ5wRZjx0hdWkCfzf/4VkffPNaUdTesyswgtwE/AsYZytU4GngRvzeNx6wAZZt18n1G3cBFwal1+a2RdwWNy3gN2BsXF5S2AmofVYi3i7ZUXH7tq1qzlXFKNHmzVqZHbssWarV+febtUqs7vvNmvRwqxBA7MzzjD79NOqHWvxYrPbbzfr3NkMzJo3N/v976u+H1d9xxxj1rSp2cyZaUeSCGCcVfLbXt4lv43gaOCWePl5no/ZApgUL1OAy+PyjYEXgGnA85mkEBPIP4EZwNuE1l6Zff0KmB4vp1V2bE8krigWLjRr396sQwezRYvye8znn5v99rch+Wywgdlf/2r27be5t58zx6x/f7MjjzRr1ix8ZXfd1eyee8yWLi3Ak3BVMnu22brrmh13XNqRJKK6iSSf5r83mNkllS0rJd781yXOLNRNDBkCo0ZVvZ/G++/DxReHx7drFyZUWn99WG+9cG0GL7wQisUA2reHXr3gtNPCti49V18dihJHjYI990w7moKqbvPffBLJBAsV5tnLJpvZjlU9WLF4InGJu+MO+M1v4KabQuus6nrxxTCQ4qefhpn5MpcVK0IF+mGHhQTSuXO9Ge+p5C1dGire27WD0aPr1PtS8EQi6SzgN4QiqhlZqzYARpnZSdUJtBg8kbhETZoEu+0WmvYOG5bM+FVmdeoHqs65775wdjhoEPTunXY0BZNEItmIULn9V0KleMZXZrawWlEWiSeSErVoUWg+uc028Mc/1s4fym++CUVLS5aEhOItBOun1auhWzf44ovQ471Zs7QjKojqJpKc7QTNbAmwBDghHuBHQFNgfUnrm9nH1Q3W1UNTpoRpTKdPD/e/+SYU6dS2ZHLBBeGH4/nnPYnUZw0ahJEL9tsv9Hy/9NLKH1OH5dOz/XBJ0wjNbkcCswjNdJ3LzxNPhLkwvvoKXn0VzjoLbrwx9BauTTPQDRkSOvxdeCHsv3/a0bi0/exnYRyuv/wl+ZELSlw+hbvXEfp1fGBmHYH9gTGJRuXqhtWr4aqrwrAdnTrB+PGw996hY1e/fvDXv4b1tcG8eaFD4c47hylYnYPwh+jbb2vP5zgh+SSSFWb2BdBAUgMze4kwEq9zFevbN/Ts7tMnjP3UOg6R1qBBaPXUty9ce21oSlnKVq+GU08NrXUeesjnO3drbLttaL3Xvz9Mrb/z/OWTSBZLWh94BXhQ0m3A0mTDcrXeE0+Eli2XXgr33gtNm669vkGD8OU79dTQLv+GG4ofY75uvz0M1vf3v8N226UdjSs1V14ZPt833ph2JKnJpx/JesC3hKRzIrAR8GA8SylJ3morZYsXh34PrVqFcacaN8697apVcOKJYWDDBx+EX/6yaGHmZfLkMJrvIYeE5FjbGge44jj3XLjrLpg5c82Zdy2U5Jztw8xstZmtNLOBZnZ7KScRVwIuugjmzw+j3VaURAAaNoSBA8MQ6KedFgYjLBUzZoQOgS1bwt13exJxuZ1/fvhTdPvtaUeSinzmbF8d+5Q4V7kXXgg/uhdemP9QHk2awOOPwxZbhCbCpVDWPH067LtvqEh95hlv6usq1rEjHHtsaNX35ZdpR1N0+dSRfA28LWlAnFP9dkn1M+26ii1dGlpjbbVVqPeoipYt4amnQlLp2TMMGZKWadNCElm2LAxhstNO6cXiao+LLgpJpH//tCMpunwSyWOEiaxeIcyXnrk4t7Y//hE+/DCckVSnp2/HjmHIkQULwvhSX39d+Bgr8/77IYmsWBGmyd2xZIeUc6Wma9fQt+TWW2H58rSjKap85mwfCAwC3gImAIPiMufWeOON8AU688zwQ1xd3bqFive33oIddgjTmxar0+Lbb4eJp1atCklkhx2Kc1xXd1x0UZhy+eGH046kqPLp2X4oYdDG24H/A6ZL6pl0YK4W+eab0Ix3s80K04y3V6/wQ77RRqHceb/9wo98UpYtC31ZunYNSeull8LMg85V1SGHhD8gN99cu0ZtqKF8irZuAX5mZj3MbF/gZ8Dfkw3L1SqXXBLGn7rvPthww8Lsc599Qk/4f/0rNMHt0gXOPhvmzCnM/jNefz30Vr/66jC/yOTJoemyc9UhhbOSt9+GZ59NO5qiySeRfGVm07Pufwh8lVA8rrZ59tkw5Ml558EBBxR2340ahXG5pk0LvYfvuitM8HTUUaEl1erV1d/3F1+ExLT33qGRwPDhodf6j35UsPBdPdW7d+hLUo86KFY0jPzR8eaBQHtgMGDAscDHZvabokRYDd4hsUg+/xx+8hPYeGN4883kh9KeORP+/W8YMCAMktexI5x+ekgs221XeT8PszBUS//+oe5lxYqQAK+9NsxK6Fyh3HzzmjOTWlTXlsR8JPdW8Dgzs1/lGVhDYBww18x6SeoIPEyYu308cLKZLZfUBLgf6Ap8ARxvZrPiPi4D+gKrgN+aWYXnjJ5IiiAz1ezQoSGJFLOJ7PLloZf5nXeG+gyANm3gwAPhoIPCWcaqVaEp5pdfhrlDpkwJrck++CDUvZxyCpxxhteFuGTMnw8//nEY4PPyy9OOJm+JTbVbU5LOJwzyuGFMJIOBx8zsYUl3ApPM7A5JvwF2NLMzJfUGfm5mx0vqTGg11h3YHHge2CZ2liyXJ5IiyMwQd8MNYe7xtHz0URgH67nnQmfIRYtyb7vnniF5HHMMrLtu8WJ09VP37mHkhtGj044kb0mckVxsZjdK+gehSGstZvbbPIJqAwwE/gycDxwOLAB+bGYrJe0BXG1mB0t6Nt4eLakR8CnQijg7o5n9Ne7z++1yHdcTScJmzgxnILvsEn68GzZMO6Jg1apQQZ8pZttww3DZaKPw77B9+7QjdPXJNdeERhyfflpr6t4KPkMi8G68rskv8q3AxYR53iEUZy02s5Xx/hwgM8JZa2A2QEwyS+L2rVl7/pPsx3xPUj+gH0C7du1qELKr0PLlcMIJoT5i4MDSSSIQYunePVycS1uvXmGekqefDlMp1GEVJZKekhZVt/OhpF7AZ2Y2XlKP6uyjKsysP9AfwhlJ0serty6+OHQ+fOQR/4fvXEV23hk23zyM1lDHE0lFzX8/AG6WNEvSjZJ2ruK+9wKOkDSLULm+H3Ab0DwWXQG0AebG23OBtgBx/UaESvfvl5fzGFdMjz4a5qc+7zz4xS/Sjsa50iaF0aOffbbOD5mSM5GY2W1mtgewL+EH/R5J70m6StI2le3YzC4zszZm1gHoDbxoZicCLwHHxM36AE/G20PifeL6Fy1U4AwBektqElt8bQ2MreoTdTU0fTr86lew2271qn28czXSqxd89RW8+mrakSQqn7G2PjKzG8xsZ+AE4CigJuN8XwKcL2k6oQ5kQFw+ANg4Lj+fNZXsUwh9WN4FngHOrqjFlkvAt9+GoUoaNYLBg2GdddKOyLnaYf/9w4jWw4enHUmi8pkhsRHQk3BWsT/wMmHgxicrelyavNVWgfXrFzoCDh8Ohx6adjTO1S6HHhrO6D/4IO1IKlXwGRIlHSjpHkIrqdOB4cCWZta7lJOIK7C77w5J5LLLPIk4Vx29eoVhfmpBIqmuioq2LgNeBzqZ2RFm9pCZLS1SXK4UPP546MB30EGhTbxzruoOOyxcDxuWbhwJqqiyfT8zu9vMKugq7OqsF18Mg8917w6PPRbqR5xzVde+fRiTrj4mElePvfkmHHkkbLNNqBdZb720I3KuduvVK7TcWrw47UgS4YnErW3q1DBneqtWof17y5ZpR+Rc7derF6xcGcaEq4M8kbg1pk8P9SGNGoUP/Oabpx2Rc3XDbruF6RbqaPGWJxIXPPZYmGr2m2/CmchWW6UdkXN1R8OGodXjU0+FwUXrGE8k9d2KFXDBBWHIk223hQkTiju3iHP1xRFHhJk562Avd08k9dmcOdCjB9xyC5xzTviA+0CMziXj0ENDw5WHH047koLzRFLfLFoUZhc891zo0gUmTw4f7H/8Iwzl4JxLxrrrhrOSzDTPdYh3Dqhrli0L85IvWBDqO5YuDZcvvgjLJ0yA1avDhzpzNrLttmlH7Vz90Ls3DBoUJoQ75JC0oykYTyR1wcqVoQPhoEGhN/qSJT/cpnHj0HLkyivDQHK77eaDLzpXbAcfHGbs/O9/PZG4EjF3bhjSfdCgcAay4Ybw85/DccfBlluG8th11w3X66wT5kdwzqWnSRM4+ugwt8+dd9aZ4mRPJLXRggVw/fXwz3+GYqqjjgrT3/bsCU2bph2dc64ivXvDvffCM8+EESTqAE8ktcmSJfC3v8Hf/x7qP045JcwJ3aFD2pE55/K1336wySahkYsnElcps3BpUIPGcatWwUsvwQMPhNPhr78Ok0xdcw1st13hYnXOFUejRnDMMXD//aEhTB0Yy84TSaEsXw4jR8K77659+fZb2HVX2H132GOPcP3jH+fez+rVMG8ezJgBQ4fCQw/BJ5+E+o/jj4ezz4addy7e83LOFV7v3qGOZNiw8L2u5SqdIbE2KvoMibNmhX8Y48eH+y1bwvbbQ+fOoZJ77NjQ7DbTdrxlS2jRApo3Dy04mjcPRVUffggffRSa8EL459KzJ5x8Mhx+uNd/OFdXrFoF7dqFaRoefzztaL5X3RkSEzsjkdQUeAVoEo/ziJldJakj8DBhvvbxwMlmtlxSE+B+oCvwBXC8mc2K+7oM6AusAn5rZs8mFXeVPf00nHhiOJN44AE48MAwcm7ZFlLffQdvvQWjR4fBEZcsCUNKL1kC770HzZrBjjuGMtOOHcOlW7dQluqcq1saNgytK//1r/AbsNFGaUdUI0kWbS0D9jOzryU1Bl6T9DRwPvB3M3tY0p2EBHFHvF5kZltJ6g3cABwvqTNhvvjtgc2B5yVtY2bpjny2alWop7j22pAAHn00NLnNpWnTULS1xx7Fi9E5V7p694Zbbw0jTfTpk3Y0NZLYECkWfB3vNo4XA/YDHonLBwJHxdtHxvvE9ftLUlz+sJktM7OZwHSge1Jx52X+/DB95jXXhJZTr79ecRJxzrmyuncPLS7rwNhbiY61JamhpInAZ8AIYAaw2MxWxk3mAK3j7dbAbIC4fgmh+Ov75eU8JvtY/SSNkzRuwYIFCTwbQg/yf/wjzBz40ktw112hPfi66yZzPOdc3SWFs5IRI8IAqrVYoonEzFaZWRegDeEsIrH2qmbW38y6mVm3Vq1aFf4AY8aE1le//W0YXuTtt6FfP+8t7pyrvn79QheB229PO5IaKcrov2a2GHgJ2ANoLilTN9MGmBtvzwXaAsT1GxEq3b9fXs5jkvfhh/DrX4e6jQULYPDgMPHTNtsULQTnXB3VsWPoF3bXXfDll2lHU22JJRJJrSQ1j7ebAQcCUwkJ5Zi4WR/gyXh7SLxPXP+ihbbJQ4DekprEFl9bA2OTihsIRVhPPBEGVdtySxg4EC68MMxnfuyxfhbinCucCy8MSeTuu9OOpNqSbLW1GTBQUkNCwhpsZsMkvQs8LOk64C1gQNx+APAfSdOBhYSWWpjZFEmDgXeBlcDZibXY+uwzuOMO+Pe/w4CIrVvD1VdD377Qpk0ih3TO1XPdusG++4YWXOeeG0bqrmW8Q2K2t96CXXYJQz2fdVZomdXIO/875xI2bFjodPzgg/DLX6YWRnU7JHoiKWv2bGjbtvLtnHOuUFavDqNhNG0aRsFIqfi8uonEp9oty5OIc67YGjSACy6AiRPDJHW1jCcS55wrBSedBJtuCjffnHYkVeaJxDnnSkHTpqGy/Zln4J130o6mSjyROOdcqTjzzDBSRi07K/FE4pxzpWLjjeH00+E//6lVdSWeSJxzrpRcdx1suy2ccEKY1K4W8ETinHOlZP314ZFHwrTaJ5wQRtoocZ5InHOu1HTuHMbfeuUVuOKKtKOplCcS55wrRSedFEYHvuGG0PO9hHkicc65UnXbbbDzzmECvVmz0o4mJ08kzjlXqpo2DfUlq1dDz55hPMAS5InEOedK2RZbhGktliwJ0/Neey2sWJF2VGvxROKcc6WuR4/Q2/244+CPf4Q994R33007qu95InHOudqgZcswzPz//hfqS3bZBf7wB5g+Pe3IPJE451ytcswx4ezkiCPg+uth661hr72gf39YtCiVkDyROOdcbbPppjB4MHz8cUgmixfDGWfAZpuF4eiLLMk529tKeknSu5KmSDovLm8paYSkafG6RVwuSbdLmi5psqRdsvbVJ24/TVKfXMd0zrl6pU0buOSScIYyblxIJu3bFz2MJOeRXQlcYGYTJG0AjJc0AjgVeMHMrpd0KXApcAnQE9g6XnYD7gB2k9QSuAroBljczxAzS+cczjnnSo0EXbuGSwoSOyMxs3lmNiHe/gqYCrQGjgQGxs0GAkfF20cC91swBmguaTPgYGCEmS2MyWMEcEhScTvnnKuaotSRSOoA7Ay8AWxqZvPiqk+BTePt1sDsrIfNictyLXfOOVcCEk8kktYHHgV+Z2ZfZq8zMyMUVxXiOP0kjZM0bsGCBYXYpXPOuTwkmkgkNSYkkQfN7LG4eH4ssiJefxaXzwXaZj28TVyWa/lazKy/mXUzs26tWrUq7BNxzjmXU5KttgQMAKaa2S1Zq4YAmZZXfYAns5afEltv7Q4siUVgzwIHSWoRW3gdFJc555wrAUm22toLOBl4W9LEuOwPwPXAYEl9gY+A4+K6p4BDgenAN8BpAGa2UNK1wJtxu2vMbGGCcTvnnKsChWqKuqVbt242bty4tMNwzrlaRdJ4M+tW1cd5z3bnnHM1UifPSCQtIBSbVdcmwOcFCqfQPLbqKeXYoLTj89iqp5Rjg/Lja29mVW6tVCcTSU1JGled07ti8Niqp5Rjg9KOz2OrnlKODQobnxdtOeecqxFPJM4552rEE0n5+qcdQAU8tuop5digtOPz2KqnlGODAsbndSTOOedqxM9InHPO1YgnEuecczVSLxKJpHskfSbpnaxlJTFTY47YbpL0Xjz+45KaZ627LMb2vqSDs5YfEpdNjxOGFUR58WWtu0CSSdok3k/9tYvLz42v3xRJN2YtL9prl+N97SJpjKSJcaTq7nF5sV+3kp29tILYSuI7kSu+rPWpfScqii3x74SZ1fkLsA+wC/BO1rIbgUvj7UuBG+LtQ4GnAQG7A2/E5S2BD+N1i3i7RUKxHQQ0irdvyIqtMzAJaAJ0BGYADeNlBrAFsE7cpnNSr11c3pYweOZHwCYl9Nr9DHgeaBLv/yiN1y5HbM8BPbNeq5dTet02A3aJtzcAPoivT+rfiQpiK4nvRK74SuE7UcFrl/h3ol6ckZjZK0DZgR5LYqbG8mIzs+fMbGW8O4YwdH4mtofNbJmZzSQMcNk9Xqab2Ydmthx4OG5bYzleO4C/Axez9nwyqb92wFnA9Wa2LG6TmaagqK9djtgM2DDe3gj4JCu2Yr5uJTt7aa7YSuU7UcFrByl/JyqILfHvRL1IJDnUlpkaf0X4R1MysUk6EphrZpPKrCqF+LYBfirpDUkjJe1aQrH9DrhJ0mzgZuCytGNTCc9eWia2bCXxnciOr9S+E2Veu8S/E0kOI19rmJlJKrl20JIuB1YCD6YdS4akdQnTARyUdiw5NCIUF+wO7EqYsmCLdEP63lnA783sUUnHEebrOSCtYFRm9lJJ369L+ztRNras5SXxnciOL8ZTMt+Jct7XxL8T9fmMJJGZGgtF0qlAL+BEiwWaJRLbloTy1EmSZsVjTZD04xKJbw7wWCxKGAusJgxOVwqx9QEyM4X+j1CEQBqxqYizlxYotpL5TpQTX8l8J3K8dsl/J6paoVNbL0AH1q74vIm1KxZvjLcPY+3KsbG2pnJsJqFirEW83TKh2A4B3gValdlue9auHPuQUDHWKN7uyJrKse2Teu3KrJvFmorFUnjtziRMfgbhlH52jKfor105sU0FesTb+wPj03jd4nHuB24tszz170QFsZXEdyJXfKXwnajgtUv8O1HjL0ttuACDgHnACkJ27gtsDLwATCO0aGiZ9Wb8k9Bq4W2gW9Z+fkWokJoOnJZgbNPjmz0xXu7M2v7yGNv7xBZAcfmhhFYaM4DLk3ztyqzP/tKUwmu3DvAA8A4wAdgvjdcuR2x7A+PjF/MNoGtKr9vehArhyVmfsUNL4TtRQWwl8Z3IFV8pfCcqeO0S/074ECnOOedqpD7XkTjnnCsATyTOOedqxBOJc865GvFE4pxzrkY8kTjnnKsRTySuVpG0SmH03Myl2qO6Snq9kLEVg6QekoblWLezpAEJHbeXpGuS2Ler/bz5r6tVJH1tZuunHUdaJPUALjSzXuWs+x9wnZUZ70lSI1sz4GF1jytCH4S9zOybmuzL1T1+RuLqBEmzJP1J0gRJb0vaLi5vpTC3xhRJd0v6KGuuiK/jdQ9JL0t6JM7Z8GD84URS1zjQ3XhJz2aGEClz7GMlvSNpkqRX4rJTJT0Z9ztN0lVZ258kaWw8o7pLUsO4/CBJo+Nz+F8cMykzN8R7kiYAR+d4/hsAO2aSiKSrJf1H0ijgP7nikdQh7vs+SR/E536ApFFxu+4Qxt4CXiYMUeLcWjyRuNqmWZmireOz1n1uZrsAdwAXxmVXAS+a2fbAI0C7HPvdmTAAX2fCPAx7xXGL/gEcY2ZdgXuAP5fz2D8CB5vZTsARWcu7A78AdgSOldRNUifgeMI/+y7AKuDEmNyuAA6Iz2EccL6kpsC/gcOBrsCPc8TfjdBzOVvnuL8TcsUTl28F/A3YLl5+SeglfSFhMMKMccBPcxzf1WM++q+rbb6NP8DlyQxSN541/9z3Bn4OYGbPSFqU47FjzWwOgKSJhHGyFgM7ACPiCUpDwrAnZY0C7pM0OCsGCPNNfBH3+ViMZSUhIbwZ99mMMDji7oQf/lFx+TrAaMIP+0wzmxb38wDQr5wYNgMWlFk2xMy+rSSeJ+L+347LpwAvmJlJeju+DhmfAZuXc2xXz3kicXXJsni9iqp/tpdl3c48XsAUM9ujogea2ZmSdiMM0DdeUtfMqrKbxn0ONLPLsldIOpzwQ39CmeVd8oz/W6BpmWVLyzl+efezn/vqrPurWft1bBqP49xavGjL1XWjgOMg1EEQRlrN1/tAK0l7xMc3lrR92Y0kbWlmb5jZHwlnBZkhuA9UmAe9GWG2wVGEQRGPkfSj+NiWktoTZv3bS9JWcfl6krYB3gM6SNoy7nOtRJNlKqGIqiLlxVMV2/DD4jPnPJG4WqdsHcn1lWz/J+AgSe8AxxJm/vsqnwNZmGb0GOAGSZMIo6nuWc6mN8UK/neA1wmj+wKMJcwNMRl41MzGmdm7hLqQ5yRNJkyxupmZLQBOBQbF5aOB7czsO0JR1vBY2f4Z5TCz94CNYqV7Lj+IJ5/XIcvPgOFVfIyrB7z5r6vTJDUBVpnZynhmcUcFdSyFPO6phCHDz0n6WFnH/D3wlZndXeh4JG0KPGRm+9csSlcXeR2Jq+vaEaYWbQAsB05POZ4k3UE460pCO+CChPbtajk/I3HOOVcjXkfinHOuRjyROOecqxFPJM4552rEE4lzzrka8UTinHOuRv4fxAfAO40oFLwAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAaAAAAEWCAYAAAAgpUMxAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuNSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/xnp5ZAAAACXBIWXMAAAsTAAALEwEAmpwYAABAxklEQVR4nO3dd5xU1fnH8c+XbkcUUSliwd5B1Ch2o2I3SjQWUBNj7xpN4s+eYElsMSqWCIlREQvYNSh2QLAhWEBFBWmiFAEL8Pz+eM6GYbNldtmZO7v7vF+vee3MvXfufebuvfPMPefcc2RmhBBCCMXWJOsAQgghNE6RgEIIIWQiElAIIYRMRAIKIYSQiUhAIYQQMhEJKIQQQiYafQKSZJI2yDqOYpHUR9KrWcdRU5ImStqrhu+5TNK/lmGbV0n6WtLU2q6jBtv6b6ySOqfjslmht7ssanruNLZzrSqSfi/prmV4/z6SHqtmmU6SvpPUtLbbqQ1J7SR9IKlldcvWaQJKH7bssVjSgpzXR9fltkpBqZ9Q9eWLrBRJ6gScB2xqZmvW18SdFUnDJP066zgKoS4+m5n9ycyWZR1XA31zYjJJ83K+b2eZ2RdmtqKZLarpyiXtLmmMpFmSZkp6VFL7nPnXSxovaa6kDyUdl/PZpgEvAidVt506TUDpw65oZisCXwAH5ky7L591xJdlw1NP/6edgJlmNj3rQEJ2anP1UOjjXdJ2wCpmNrzcrK1yvm9bL+NmxgH7pPWsDYwHbsuZPw84EFgF6A3cJOlnOfPvA35b3UaKUgQnqaWkGyV9lR43ll2eSdpN0iRJv0tFHf+QtJykeyV9K2mcpAskTcpZ31JXHmnZq3JeHyDpnZS9X5e0ZTUh9pT0aSpuuU7Sf/eLpBPS5eS3kp6VtE6a/nJa5N30i+OXkl6S9Is0f6cU5/7p9Z6S3qluvWnexpKel/SNpI8k9Sr3WW+V9GT69TFC0vqVfK6yGGelGHfMWc/1adufSdovZ/oqku6WNEXSZHkxVIUnoaSm8qKET1IsoyV1TPNM0mmSxuMHL5JukvSlpDlp2R4567pM0kBJA9K6xkrqVsl2N0lxH5Ve/y7FOjftrz1zFm9R2TolXZQT+zhJh6bpewHPA2un/fYgcDuwY3o9q5K4hkn6s6SR6TMOltQmzdst9xhO0/IqVqzsfyKpRTpGtshZdg1J8yW1rWA960t6Qf6L9mtJ90lqXS6e8yW9J2m2pAcltcqZf0GK4StJJ1QR79VAD+BvaX/9LWf2XvJfzrPScayc91V6TlSwjeMkfZ4+yyW5+1JSk5z/7cx0XLXJee9Dkqamz/iypM1y5t0r6TZJT0maB+yez2er5fFevsi1t6Qv0v/mD5V9dmA/4KUq5pdtY6kSEEnrps87V9J/0v6vsIjazKaZ2Vc5kxYBG+TMv9TMPjSzxWY2AngF2DFn+RHAelX9D8tWVJAHMBHYKz2/AhgOrAG0BV4HrkzzdgMWAtcALYHl8EvLV4A2QEfgfWBSzroN2CDn9b3AVen5NsB0YHugKZ6dJwItK4nT8MvFNviv3o+BX6d5BwMTgE2AZsAfgderiOMK4Jb0/PfAJ8A1OfNuqm69wArAl8Dxad42wNd4UVDZZ50JdE/z7wMeqOSzdU4xNsuZ1gf4CfhN2j+nAF8BSvMfBe5IcawBjAR+W8n6LwDGABsBArYCVsvZN8+n/bpcmnYMsFqK+zxgKtAqzbsM+B7omeL6MzC8/PEEbItfXR+Qpm+U9tfaOZ95/TzXeQT+664J8Ev8V91aOcflpHL77dVqjvlhwGRg87T/Hgb+VdH6KjhHLstZdqn/W1X/E+DvpGMsvT4LeLyS+DYA9sbPs7b4D5Qby8UzMu2TNsAHwMlp3r7AtJzP9m/KHf8V7ItfV3CuPQG0xs+1GcC++Zxr5dazKfAdsDPQArgeP6b3ytkHw4EO6bPeAdyf8/4TgJXSvBuBd8p9l8wGdkrHRasafLaaHu/l/9934t9/WwE/AJtU8vkfAi6oYPsblJtW/jh6I+2rFmnfzSmLoZLtdAJmAYvT/u1TyXLLAVPK/pc5098DDqrynKlq5rI8WPrk+gTomTNvH2Bizon5Y+4/Gvg098PgZYn5JqDbSMktZ/5HwK6VxGnltnUqMDQ9fxo4MWdeE2A+sE4lcewJvJeePwP8mvSFh/9iOay69eJfhK+Ui/EO4NKcz3pXzryewIeVfLalDsA0rQ8wIef18mmZNYF26cBfLmf+UcCLlaz/I+DgKvbrHtUcI9/ixQbgJ+R/cuZtCiwodzxdDkwCdsuZvgH+g2MvoHm59Ve5zgrieafs81D7BNS33PZ+xJPfUuur4By5jAoSUHX/E/yH1hcs+QExCuiV5zl6CPB2uXiOyXl9LXB7en5Puc+2IbVLQDvnvB4IXJTPuVZuPf/H0gll+bSfy/blB8CeOfPXwr9Am1WwrtYprlVyzq8BefyfK/psNT3ey/+/O+QsOxI4spL1PE/6YVBu+3PwhDELuLnccdQJ/6G/fM57/kUVCShnuTbA74AdKpnfH/++U7nprwHHVbXuYrWCWxv4POf152lamRlm9n255b8st3y+1gHOS5f4s1JxScdy2yuv/LbKll0HL9ssW883+C/99lTsDWBDSe2ArYEBQEdJq+NXLGVFYlWtdx1g+3LxH40niDK5rbLmAytW8dkq8t/3m9n89HTFtO3mwJScbd+B/+quSEf8x0VlcvcrqXjng1T0MQsvP169orjwz9VKS5enn4z/Kh6WE/8E4Gz8hJ4u6QFJuf/rSteZinHeyfmsm5eLpzbKH0vNl3GdVf5PzIs/5gO7SdoYT8hDKlqRvHXSA6kYbw7+BVQ+tsqOrWU5J/NZf03OtaViScfwzJz56wCP5qzrA7wIqZ286LJvKp6bgyddWHo/LHXc1kBNj/fy8j2vv8Wv4Mrb1sxap8eZ5eatDXyTc77/T7yVMbNv8CQzuNz5iKTr8POml6Wsk2MlPBlWqlgJ6Cv8oCjTKU0rUz7wKfiXW+7yuebjv3rK5H45fwlcnfOPaG1my5vZ/VXEV35bZbF9iRd15K5rOTN7vaKVpH/uaLwI4H0z+xEvbjwX+MTMvs5jvV8CL5Wbt6KZnVJF/JUpv1+r8yX+a3v1nG2vbGabVbF8ZfVPS20/lX9fCPQCVjWv3JyNf8nk62Sgk6QbltqI2b/NbGf8GDO8OLdKqWz6TuB0vNiwNV7UW1k8+e7L8sfST3gR6jxyjll5vdr/1NNUIJ//SX+8uOdYYFC5H3O5/pQ+xxZmtnJ6T777v7pzsrzaHHv5nmtT8OI1ACQthxd15a5rv3LramVmk4Ff4cV9e+EJoXPZamoQe2Xz6/p4r8x7+BVoTUwB2kjK/d7sWNnCFWiG/+hZuWyCpMvx+qifm9mc3IVTotoAeLeqlRYrAd0P/FFS23Q18H/4r6/KDAQulrSqpA7AGeXmvwP8Kv2a2RfYNWfencDJkraXW0HS/pIq+sVQ5oK0rY548ngwTb89xbEZ/Lcy+Iic900D1iu3rpfwL7WySsJh5V5Xt94n8KuoYyU1T4/tJG1SRfyVmYGX35aPsUJmNgV4DviLpJXllbnrS9q1krfcBVwpqUva11tKWq2SZVfCiwBmAM0k/R85B3Oe5uJ1EbtI6gsgaSNJe8gbtXwPLMA/c3VWwL8wZqT1HI//kqvMNKCDpBbVrPcYSZumE/0KPCEswusWW6VjsTlex1HtfRJ5/k/+BRyKJ5QBVaxuJbzuZLa8Se0F1W0/x0CgT85nu7Sa5Ss6N6pS3bmWaxBwoKSfpf/HZSz9xX47cLWWNBhqK+ngNG8lPKHPxH8Q/KkGMZbJ57PVxfFemadY+juvWmb2OV48e5m88cqOeCu2Ckk6LJ1bTeQNWv6KF9d+k+ZfjCfzvcxsZgWr6I5Xs1R5pVysBHQV/uHfwyut30rTKnM5fon/GX7y/bPc/LPwnTcLL556rGyGmY3CK9j/hl+qTsDL76syGL9yeQd4Erg7retR/Nf0A+ly/X0845e5DOifLvXLWqq9hB98L1fyusr1mtlc4OfAkfiV2FSWNNCokXRFdjXwWopxhzzedhxeSTkO33+D8DL0ivwV/2J6Di9/vhuvkKzIs3g58cf4//Z7alHUYWaz8Ir0/SRdie+XvvhVxlT8V9rFeaxnHPAXvNh0GrAFXmZdmReAscBUSV9Xsdw/8XqEqUAr4My0vdl4/eJdeEOFeXh9Vj6q/J+Y2Zf4OWV4453KXI434piNH+eP5Ll9zOxpvML+BfyceqGat9wEHC5v0XZzHuuv7lzLXXYs/qP0AfyX/Xd4PeAPOdseAjwnaS7eIGH7NG8AfvxNxvdn+abM+cjns9XJ8V4RM3sL/xGxfbULL+1ovKXaTPz790GW7LOy+zjLWuq1x+Ofi39nL8Z/5JT5E34VPEFL7j36fblt3V5dQGUVlyVN0m54ZVmHahYNITOShuHHaa3vcF+Gbd8DfGVmfyz2trMmaUX8x2gXM/ss43CKQtLPgVPN7JBlWMeDeAOm6q5ma7reNfAf3ttUURwMRFc8IdR7kjoDh5Gu3BsDSQdKWl7SCnjT4jEsaVDQ4JnZczVNPqkof/1UrLYvXhf2WAFim25mm1SXfCASUAj1WiqGfB+4rrH8+k8OxouovwK64E2WS784J1tr4nXS3+HNtE8xs7ezDKheFMGFEEJoeOIKKIQQQiYK3WneWXiLNAF3mtmN8j6ZHsTb30/Eb2D6VpLw1iU98ft8+qTWHpVaffXVrXPnzoX7ACGE0ACNHj36azPL5z60gipYApK0OZ58uuPdZDwj6Qm8W52hZtZX0kXARXg3D/vhZbld8CaTt7Gk6WSFOnfuzKhRowr1EUIIoUGSVNueLOpUIYvgNgFGmNl8M1tI6gsNrzzsn5bpj/dHRZo+wNxwoLWkyu4/CSGEUM8VMgG9D/SQtFq6c7on3vVDu3R3N/jNeu3S8/YsfaPWJCroB0rSSZJGSRo1Y8aMwkUfQgihoAqWgMzsA/zO5ufwO2rfwTsEzF3GqGGfUWbWz8y6mVm3tm0zL8IMIYRQSwVtBWdmd5tZVzPbBe9C5GNgWlnRWvpbNuLkZJbuHK9DmhZCCKEBKmgCSl0yIKkTXv/zb7yPpt5pkd54P2yk6celTi13AGbnFNWFEEJoYAraDBt4OPWO/BNwmpnNSr0YD5R0It5JX1knnk/h9UQT8GbYxxc4thBCCBkqaAIysx4VTJuJjxxafroBpxUynhBCCKUjekIIIZSO11+HYcOyjiIUSaGL4EIIIT8//ghHHAHffQfjx8MalY0EHxqKuAIKIZSGQYPgq69gzhy47LKsowlFEAkohJA9M7jhBthwQzjtNLjjDhg7NuuoQoFFAgohZO+NN2DUKDjrLLj8clh5ZbjggqyjCgUWCSiEkL2bboLWreG442C11eCSS+Dpp+HZZ7OOLBRQJKAQQra++AIefhh+8xtYcUWfdtppsN56cN55sHBhtvGFgokEFELI1q23+t/TT18yrWVLuPZarwe6555s4goFFwkohJCdefOgXz847DDo1GnpeYcdBjvv7MVxc+ZkE18oqEhAIYTsDBgAs2bB2Wf/7zwJ/vpXmD4drr++2JGFIogEFELIxuLF3vhgu+1gxx0rXma77WCPPeDxx4sbWyiKSEAhhGw89xx89JFf/UiVL7fjjjBmDMyfX7TQQnFEAgohZOPee6FdOzj88KqX694dFi2Ct94qSliheCIBhRCyMWIE7LortGhR9XLdu/vfkSMLH1MoqkhAIYTimz4dJk5cklyqsuaa3kJuxIiChxWKq9Ajop4jaayk9yXdL6mVpHUljZA0QdKDklqkZVum1xPS/M6FjC2EkKE33/S/+SQggO23jyugBqhgCUhSe+BMoJuZbQ40BY4ErgFuMLMNgG+BE9NbTgS+TdNvSMuFEBqikSOhSRPYdtv8lu/e3a+Ypk8vaFihuApdBNcMWE5SM2B5YAqwBzAoze8PHJKeH5xek+bvKVXVNCaEUG+9+SZsthmssEJ+y0c9UINUsARkZpOB64Ev8MQzGxgNzDKzss6dJgHt0/P2wJfpvQvT8quVX6+kkySNkjRqxowZhQo/hFAoZp5I8i1+A+ja1a+YIgE1KIUsglsVv6pZF1gbWAHYd1nXa2b9zKybmXVr27btsq4uhFBsn30GM2fWLAGtsAJsvnk0RGhgClkEtxfwmZnNMLOfgEeAnYDWqUgOoAMwOT2fDHQESPNXAWYWML4QQhbKrmK2265m7ytriGBW9zGFTBQyAX0B7CBp+VSXsycwDngRKLvzrDcwOD0fkl6T5r9gFkdaCA3OyJHQqpVf0dRE9+7eb9yECQUJKxRfIeuARuCNCd4CxqRt9QN+B5wraQJex3N3esvdwGpp+rnARYWKLYSQoTff9NZvzZvX7H1lRXZRDNdgNKt+kdozs0uBS8tN/hT4n8JfM/seOKKQ8YQQMrZwIYweDb/9bc3fW9ZqbuRIOOaYuo8tFF30hBBCKJ6xY2HBgprX/wA0beqt4eIKqMGIBBRCKJ6yBgg1aQGXa/vt4Z134Icf6iykkJ1IQCGE4nnzTVh1VVh//dq9v3t3+PFHePfduo0rZCISUAiheMpuQK1tJyfbb79kPaHeiwQUQiiOefPg/fdrV/9TpkMH7x07ElCDEAkohFAcb7/tA8vVtv4H/Mqpe/doiNBARAIKIRRHbXtAKG/77eHjj+Hbb5c9ppCpSEAhhOJ4800fWG7NNZdtPWVXUGVjCoV6KxJQCKE4Ro5c9qsfWLKOKIar9yIBhRAK7+uv4dNPl63+p8wqq8Cmm8Ibbyz7ukKmIgGFEAqvpkNwV2fnneG117xRQ6i3IgGFEApv5Ehvwda1a92sr0cPmDMHxoypm/WFTEQCCiEU3ogRXmy20kp1s76dd/a/r75aN+sLmYgEFEIoLDNPQDvsUHfrXGcdvyn1lVfqbp2h6CIBhRAKa8IE+OabJd3o1AXJi+FeeSVGSK3HCpaAJG0k6Z2cxxxJZ0tqI+l5SePT31XT8pJ0s6QJkt6TtG2hYgshFFFZc+m6vAICL4abMgU++6xu1xuKppAjon5kZlub2dZAV2A+8Cg+0ulQM+sCDGXJyKf7AV3S4yTgtkLFFkIoouHDYcUVvQ6oLvXo4X+jGK7eKlYR3J7AJ2b2OXAw0D9N7w8ckp4fDAwwNxxoLWmtIsUXQiiUESP85tGmTet2vZttBq1bRwKqx4qVgI4E7k/P25nZlPR8KtAuPW8PfJnznklpWgihvlqwwAeQq8v6nzJNmsBOO0VLuHqs4AlIUgvgIOCh8vPMzIAa1SBKOknSKEmjZsyYUUdRhhAK4u23YeHCuq//KdOjB3z0EUyfXpj1h4IqxhXQfsBbZjYtvZ5WVrSW/pYdOZOBjjnv65CmLcXM+plZNzPr1rZt2wKGHUJYZsOH+99CXAHBknqg114rzPpDQRUjAR3FkuI3gCFA7/S8NzA4Z/pxqTXcDsDsnKK6EEJ9NGKE37OzrD1gV6ZrV2jZMuqB6qlmhVy5pBWAvYHf5kzuCwyUdCLwOdArTX8K6AlMwFvMHV/I2EIIRTB8eOGK38CTz/bbRwKqpwqagMxsHrBauWkz8VZx5Zc14LRCxhNCKKIpU+CLL+Dsswu7nR49oG9f+O47b+4d6o0qi+Ak7Sjp1nRj6AxJX0h6StJpklYpVpAhhHqo7AbUQtX/lNl5Z+8Vu6y+KdQblSYgSU8DvwaeBfYF1gI2Bf4ItAIGSzqoGEGGEOqhESOgeXPYZpvCbudnP/Mm2VEMV+9UVQR3rJl9XW7ad8Bb6fEXSasXLLIQQv02fDhstRUst1xht7Pyyr6duB+o3qn0Cqh88pG0curHrY2kNhUtE0IIgBeJvflmYRsg5Np5Z094P/1UnO2FOlFtM2xJv5U0FXgPGJ0eowodWAihHhs7FubNK14C6tED5s+Ht94qzvZCncinFdz5wOZxtRNCyFuxGiCUKbsh9aWXirfNsMzyuRH1E/y+nBBCyM/w4bDaarD++sXZ3pprwhZbwBNPFGd7oU7kcwV0MfC6pBHAD2UTzezMgkUVQqjfRozwKxGpeNs89FC48krvF26NNYq33VBr+VwB3QG8AAxnSR3Q6EIGFUKox+bMgXHjilf/U+bQQ3101CFDirvdUGv5XAE1N7NzCx5JCKFhePNNTwTFrovZaivo3BkefRR+/evibjvUSj5XQE+nIRDWKt8MO4QQ/scrr3jRW/fuxd2u5FdB//kPzJ1b3G2HWsknAR1FqgcimmGHEKrz5JNe/Na6dfG3fcgh8OOP8PTTxd92qLFqE5CZrVvBY71iBBdCqGemTIFRo+DAA7PZ/k47Qdu2XgwXSl6ldUCSdqnifWZm0fFSCGFpTz3lfw84IJvtN20KBx0EAwfCDz/4cA2hZFXVCOGCCqYZsCU+cmnTgkQUQqi/nngCOnWCzTfPLoZDD4W774YXXoD99ssujlCtqvqCOzD3gQ8k1xyYChxSpPhCCPXF99/D88/71U8x7/8pb889fVygxx7LLoaQl3z6gttT0jDgKuCvZraDmT2ez8oltZY0SNKHkj5I4wu1kfS8pPHp76ppWUm6WdKENP7Qtsv0yUIIxfXSS97/W1bFb2VatfIrn8GDvVPUULKqGg9of0mv433B/dHMdjez52u4/puAZ8xsY2Ar4APgImComXUBhqbXAPsBXdLjJOC2Gm4rhJClJ56A5ZeH3XfPOhIvhps2LQapK3FVXQE9DnQAFgIXShqS+6huxWnE1F2AuwHM7EczmwUcDPRPi/VnSXHewcAAc8OB1pLWqsVnCiEUm5knoL328iuQrPXs6YPhRWu4klZVI4Rl/RmzLjAD+IekrfD7h84C2pnZlLTMVKBdet4e+DLn/ZPStCk505B0En6FRKdOnZYxxBBCnRg3DiZOhN//PutI3CqreF3QY4/BdddlWycVKlVpAjKzl+pg3dsCZ5jZCEk3saS4rWwbJslqslIz6wf0A+jWrVuN3htCKJCyXqh79sw2jlyHHAInnwzvv+89ZYeSU1Ud0OOSDpTUvIJ560m6QtIJVax7EjDJzNLAIAzCE9K0sqK19Hd6mj8Zb95dpkOaFkIodU88AdtuC+3bZx3JEgcf7PcF3Xtv1pGESlRVB/QboAfwoaQ3JT0l6QVJn+I9ZI82s3sqe7OZTQW+lLRRmrQnMA4YAvRO03oDg9PzIcBxqTXcDsDsnKK6EEKpmjkTXn89+9Zv5a25JvTqBXfeCbNnZx1NqEBVRXBTgQvxBgidgbWABcDHZpbvAHVnAPdJagF8ChyPJ72Bkk4EPgd6pWWfAnoCE/AB8I6v8acJIRTfM8/A4sWll4AAzjsP7r8f7rrLn4eSIrP6W43SrVs3GzUq+kUNIVNHHQUvvghffQVN8unfuMh23x0mTIBPP/WWcQFJo82sW9ZxlODREkKoN376ya+A9t+/NJMP+JXPpEnw0ENZRxLKKdEjJoRQL7z6KsyaVZrFb2V69oSNNoLrr/f7lULJyCsBSVoupzFBCCG4O+6AlVeGvffOOpLKNWniV0Fvvw3DhmUdTciRT19wBwLvAM+k11vn0xNCCKGB++ILGDQITjrJO/8sZcce6+ME/eUvWUcScuRzBXQZ0B2YBWBm7+C9HIQQGrNbbvG/Z5yRbRz5aNUKTj/dR2v94IOsowlJPgnoJzMr34g+ClJDaMzmzoV+/eDww338n/rglFM8Ef31r1lHEpJ8EtBYSb8CmkrqIukW4PUCxxVCKGX/+AfMmQPnnpt1JPlr2xb69IEBA2Dq1KyjCeSXgM4ANgN+AO4H5gBnFzCmEEIpW7QIbrwRfvYz6N4962hq5txz/abZCy/MOpJAHgnIzOab2R/MbDsz65aef1+M4EIIJWjwYPjss/p19VOmSxe4+GL45z/hqaeyjqbRq7YnBEkvUkGdj5ntUaig8hU9IYSQgR49YPJkGD/eO/usb374wTtOnTMHxo71ZuSNTKn0hFDVeEBlzs953gr4BT5IXQihsRk50m8+veGG+pl8AFq2hHvu8SLECy+E22/POqJGq9oEZGajy016TdLIAsUTQihlN9zgVwwnVDUSSz2w/fZw9tneIu7II2G33bKOqFHK50bUNjmP1SXtA6xShNhCCKXk88+9P7Vf/7phFFtdeSWsv75/nvn5dvAf6lI+reBGA6PS3zeA84ATCxlUCKHEmMGpp3rx1VlnZR1N3Vh+eR+m4ZNP4JJLso6mUcqnCC56PQihsbvvPm81duON9efG03zsthv89rf+ufbfH/bIvG1Vo5JPK7jDqppvZo9U8d6JwFxgEbDQzLpJagM8CHQGJgK9zOxbSQJuwgelmw/0MbO3qtp2tIILoQimTYNNN/UepV95pf42PqjM3LleJzRtmjeyWH/9rCMquFJpBZdPEdyJwN3A0elxF3ACcCCQTx/su5vZ1jkf9iJgqJl1AYam1wD7AV3S4yTgtnw/RAihgE4/HebN85ZjDS35AKy0EgwZ4sWMBx3kzbNDUeSTgJoDm5rZL8zsF3ivCM3N7Hgzq01TmIOB/ul5f+CQnOkDzA0HWktaqxbrDyHUlUGD/HHZZbDxxllHUzgbbOCf86OP4OijvbeHUHD5JKCOZjYl5/U0IN9CYAOekzRa0klpWruc9U0F2qXn7YEvc947KU1biqSTJI2SNGrGjBl5hhFCqLGZM+G00/ymzfPPr375+m6PPeCmm+CJJ+CPf8w6mkYhnxtRh0p6Fu8HDuCXwH/yXP/OZjZZ0hrA85I+zJ1pZiapRj1rm1k/oB94HVBN3htCqIFzzoFvvoHnnoNm+XxVNACnngpjxkDfvrD55n41FAomn1Zwp0s6FNglTepnZo/ms3Izm5z+Tpf0KD6u0DRJa5nZlFTENj0tPhnomPP2DmlaCKHYHn/c+0u75BLYaqusoykeycc5+vBDOPFE2Gwz2HrrrKNqsKosgpPUVNKHZvaomZ2THnklH0krSFqp7Dnwc+B9YAjQOy3WGxicng8BjpPbAZhdrugvhFAMM2b4zZlbbdU4i6KaN/cbbldbzXtJ+O67rCNqsKpMQGa2CPhIUm0a/rcDXpX0LjASeNLMngH6AntLGg/slV4DPAV8CkwA7gROrcU2QwjLwswHbvv2Wx83p0WLrCPKRtu2fu/Txx97K8BQEPkU7K6KD0o3EphXNtHMDqrqTWb2KfA/1+5mNhPYs4LpBpyWRzwhhEL597/h4Ye9DmTLLbOOJlu77eZFkFdcAXvuCccem3VEDU4+N6LuWtF0M3upIBHVQNyIGkIdmjQJttjCbzp9+eWGec9PTS1c6K3j3nrLHxtumHVEdaLe3IiaEs2HwErp8UEpJJ8QQh0y80r3H3+E/v0j+ZRp1syL4lq29PqgH37IOqIGJZ/esHvhdThHAL2AEZIOL3RgIYQiuv12b259/fV+U2ZYomNH+Mc/4O23YyjvOpZPEdy7wN5mNj29bgv8x8wyb5sZRXAh1IFJk7yftx494OmnvSly+F9nnulNtF97zQezq8fqTREc0KQs+SQz83xfCKE++L//87qO22+P5FOVP/0J1l7bh6NYvDjraBqEfBLJM5KeldRHUh/gSeDpwoYVQiiKMWPg3nv9133nzllHU9pWXBGuuQZGjfJ6srDMqi2Cg/8OybBzevlKvjejFloUwYWwjHr2hOHDfVC2VVfNOprSZwY77QSffur3CNXTkWHrTRGcpGvM7BEzOzc9HpV0TTGCCyEU0NChXufzhz9E8smX5B2WTpsGV12VdTT1Xj5FcHtXMG2/ug4khFBEixd7i6511vEer0P+ttsOjj/eR1EdPz7raOq1ShOQpFMkjQE2kvRezuMz4L3ihRhCqHMPPOA3Vl59NbRqlXU09c+f/uT77bzzso6kXqu0DkjSKng3PH9myailAHPN7JsixFatqAMKoRZ++MGbXbdp4xXqTaJRa61cd51fRT7zDOyzT9bR1EjJ1wGZ2Wwzm2hmR5nZ5zmPkkg+IYRauvVW+Pxz/wKN5FN7Z57pN+2ee643Tgg1FkdfCI3JN9945fm++3oHm6H2Wrb0e6jGjYNXX806mnopElAIjcmVV8Ls2XDttVlH0jAcdhistJJ31RNqLJ9m2GdIijaaIdR3H38Mf/ubDza3xRZZR9MwrLAC9OoFAwfGwHW1kM8VUDvgTUkDJe0r1ayvjjSq6tuSnkiv15U0QtIESQ9KapGmt0yvJ6T5nWv8aUIIlbvwQm+5dcUVWUfSsPTpA/PmwSOPZB1JvZPPcAx/BLoAdwN9gPGS/iRp/Ty3cRbwQc7ra4AbzGwD4FvgxDT9RODbNP2GtFwIoS68+CIMHgy//z20a5d1NA3LTjt5Y4QohquxvOqA0milU9NjId48e5CkKguSJXUA9gfuSq8F7AEMSov0Bw5Jzw9Or0nz96zp1VYIoQKLFnlLrXXWgXPOyTqahkfyq6Bhw+Czz7KOpl7Jpw7oLEmjgWuB14AtzOwUoCvwi2refiNwIVDWdexqwCwzW5heTwLap+ftgS8B0vzZafny8ZwkaZSkUTNmzKgu/BDCgAHwzjs+zHbcdFoYxx7riSg6Ka2RfK6AVgUOM7N9zOwhM/sJwMwWAwdU9iZJBwDTzWx03YTqzKyfmXUzs25t27aty1WH0PB895339bbDDvDLX2YdTcPVqZM3a+/fP4ZqqIEqE5CkpsCRZvZ5RfPN7IOKpic7AQdJmgg8gBe93QS0ltQsLdMBmJyeTwY6pu02A1bBxx4KIdTWtdfClClwww0x1k+hHX88TJwIL7+cdST1RpUJyMwWAR9J6lTTFZvZxWbWwcw6A0cCL5jZ0cCLQNmQ3r2Bwen5kPSaNP8Fy2esiBBCxd591xPQkUf6FVAorEMO8eEZojFC3vItghsraaikIWWPZdjm74BzJU3A63juTtPvBlZL089l6f7nQgg1MX++J542beDmm7OOpnFYfnkv5hw0CObOzTqaeqFZ9YtwybJuxMyGAcPS80+B7hUs8z1wxLJuK4SAt3b76CN47jmIutLiOf54uPNOeOghOOGErKMpefncB/QS8CGwUnp8kKaFEErRI49Av35wwQWw115ZR9O47LADbLhhFMPlKZ9m2L2AkfjVSS9ghKTDq35XCCETX3wBJ54I3bp5v2+huCTv6ujVV+G9GDatOvnUAf0B2M7MepvZcXjx2TIXy4UQ6tiiRXDMMbBwIdx/P7RokXVEjdOJJ8Jyy3m/e6FK+SSgJmY2Pef1zDzfF0IopiuugFde8fF+Ntgg62garzZt4Oij4V//8uEvQqXySSTPSHpWUh9JfYAngacLG1YIoUb69vUEdNxxfld+yNbpp8OCBXDPPVlHUtLyaYRwAXAHsGV69DOzCwsdWAghD2Y+KNrFF8NRR8Fdd8UNp6Vgq62gRw+/Gl20KOtoSlY+jRCuMbNHzOzc9HhUUvRUHULWzHyIhSuv9Ca///wnNG+edVShzBlneM8ITz6ZdSQlK58iuL0rmLZfXQcSQqiBxYu9mOf66+G00/zek6ZNs44q5DrkEGjfHm65JetISlalCUjSKZLGABtJei/n8RkQ7QtDyMq0aX7H/d//7vf63HILNIl2QSWneXM45RT4z3/gg6q6zWy8qjpq/w0ciPfRdmDOo6uZHVOE2EIIuRYvhjvugI03hiFD4Jpr/BF1PqXrN7/x5vC33pp1JCWp0gRkZrPNbKKZHZV6w14AGLBibTonDSEsg3ff9ZE3Tz4ZttnGb3K88MJIPqVujTW8T77+/WHOnKyjKTn5NEI4UNJ44DPgJWAi0Qw7BDdjBjzwgA95UAhffeWV2V27wief+OByQ4fCRhsVZnuh7p1xho/LdO+9WUdScvIpOL4K2AH42MzWBfYEhhc0qhBKmRm89JI3e+7Qwf+uuy6cemrdDck8dap3KLr++nDbbd69y4cfLhl5M9Qf3brBjjvCjTfCTz9lHU1JyScB/WRmM4EmkpqY2YtAtwLHFULpMfOK/003hd12g2ee8SKxF16A3r3h7ruhSxdPEmPH1m4bX30F558P663njQuOOgo+/hhuv93vsA/10x/+4D9OBgzIOpKSkk8CmiVpReBl4D5JNwHzChtWCCXob3/zJs+rrOLFKZMnw003we67e+OATz+Fs87y3qg339yHaH7wQfjxx6rX+9NP8NhjcOCB0LGjj156xBF+xXPPPZ6MQv3Wsyd07w5XXVX98dCIqLpBRyWtgDdAaAIcjQ+VfV+6Kqrqfa3wpNUSH3dokJldKmldfIju1YDRwLFm9qOklsAAoCve39wvzWxiVdvo1q2bjRo1qtoPGcIyGzUKfvYz2HdfGDy46mKwmTM9Id15p9+I2LYt9OnjldGSV0bPmQOzZ3tjggEDvGn1Wmv5eDInnOBFb6FhefppT0T9+nnruAxJGm1mmZdkVZmAJDUF/mNmu9d4xZKAFczsO0nNgVeBs/DRTh8xswck3Q68a2a3SToV2NLMTpZ0JHComf2yqm1EAgpFMWsWbLutd6ny9tv5F4UtXuwDwvXr582mK+qSpWlTv/I58URPbs3yGSMy1EtmXhc0ZQqMH59pb+WlkoCqPNrNbJGkxZJWMbPZNVmxeWb7Lr1snh4G7AH8Kk3vD1wG3AYcnJ4DDAL+JklW3SVaCIVk5lckX37pPU3XpB6mSRNPKvvu63U7w4b5sM0rr+zFeCuv7M10V1mlYOGHEiLBZZfBfvv5gHW//W3WEWUun59b3wFjJD1PTt2PmZ1Z3RvTFdRoYAPgVuATYJaZLUyLTALap+ftgS/TuhdKmo0X031dbp0nAScBdOoUtyOFArvlFnj0UfjLX3y0y9pae2341a+qXy40bPvs48fR1Vd7sWzLlllHlKl8GiE8gg9A9zKeTMoe1TKzRWa2NdABH8hu49qFudQ6+5lZNzPr1jbGug+F9Oab3iLtwAO9SXQIy0qCyy/3K+oYqqH6KyAz6y+pBZ48DPjIzGrUjMPMZkl6EdgRaC2pWboK6gBMTotNBjoCkyQ1wxs7VNnQIYSCmTULevXyhgH33hv33oS6s/fe3qDlT3/y4t1GfBWUT08IPfGis5uBvwETJFXbG7aktpJap+fL4b1qfwC8CByeFusNDE7Ph6TXpPkvRP1PyISZ3/g5aZI3o477b0JdKrsKmjTJx29qxPKpA/orsLuZTQCQtD75jYq6FtA/1QM1AQaa2ROSxgEPSLoKeBu4Oy1/N/BPSROAb4Aja/xpQqgLt90GDz8M1123bPU+IVRmzz1h113h0ku9ef5qq2UdUSbyuQ/oTTPbLue1gJG507ISzbBDnXvnHdh+e9hrL3j88RjmIBTOmDHesewJJ3hT/SIqlWbYVY0HdJikw4BRkp6S1EdSb+Bx4M2iRRhCscyd6/U+q6/uvRdH8gmFtMUW3rjlzjvhjTeyjiYTVZ1hZeP/tAKmAbsCuwEz0rQQGg4z79ftk0/g/vs9CYVQaJde6h3annIKLFxY/fINTKV1QGZ2fDEDCQ3MiBGwySZ+s2V98I9/wL//DVdeCbvsknU0obFYcUXvT/AXv/C+Bs8+O+uIiqqqIrgL099bJN1c/lG8EEO9c9ttXnm//vp+I2epd744eLBf/eyxB1x8cdbRhMbm0EO9d4RLLvEObhuRqorgxqW/o1j6BtS8b0QNjdCQIXD66X6vwxZbwJln+vAFDz7oxVyl5qGH4PDDva+3hx/2vtlCKCbJr34WLoRzz806mqKqKgHtJ2knM+tf0aNoEYb6Y8QIb1Latat3XzN0KDz1lPd/duSR3rps3Ljq11Ms//qXx7XDDt5paOvWWUcUGqv11vMxgwYOhGefzTqaoqkqAX0MXC9poqRrJW1TrKBCPTRhAhxwgPcc8MQTsMIK/stuv/28B+l774XPP/cv+yeeyDpa7wbluOP8Xoxnnqk/dVWh4brgAth4Y2+W/fXX1S/fAFSagMzsJjPbEW/9NhO4R9KHki6VtGHRIgylb8YM7/HZzL/M11hj6flNm/qIoaNG+YihBx0E11yTTZGcmQ+NfOKJ8POfw5NPerIMIWstW3pDmK+/9uOzFIus61i1NzqY2edmdo2ZbQMcBRyCd6kTgo/medBBPtzAE094gqlMx44+pEGvXnDRRT509YIFxYv188+9N+JzzvGYH3sMlluueNsPoTrbbAPXXut1qX//e9bRFFw+fcE1k3SgpPvw7nc+Ag4reGShfrj5Zhg+3Jsx59NtzfLL+302V10F993nTZ6HDi3sr73Fi+H2232Y7Dfe8FZ6jz4KreJ2tlCCzjwT9t8fzjvPR8xtwCrtikfS3vgVT09gJD6M9mAzm1fhGzIQXfFk7IsvvIXbHntUP0x1RQYP9qKGmTNhww19gK4+feq2888JE3y9L7zg3evcdRess07drT+EQpgxA7bcElZd1Yuul1++Tldf8l3xABcDrwObmNlBZvbvUko+oQScdZZfXdx8c+2GKzj4YO8R+J//9J4HzjvPB2479lgYNMiHRKgNM3jtNW9evdFGPq5Pv37e0i2ST6gP2rb1Vpofftiwb041s3r76Nq1q4WMDBliBmZ9+9bdOt991+zUU81at/Z1N21qtvPOZldfbTZypNm8eVW/f8ECs/vvN+ve3d+/6qpmF19s9tVXdRdjCMV00UV+LD/wQJ2uFhhlJfAdXm1v2KUsiuAyMm8ebLaZtx57+21o0aJu179wodcrPf20P95+26dL0LmzF/ttuqkXT3z6qfffNmGCX02ZeUOIc87xZtbRwi3UZz/9BLvvDm+95Vf129TN3TClUgQXCSjU3MUXQ9++8NJLxek3bepUePVV+OADv5F13DgvmvjxRy+q2GCDJY/ttvOWbtGTdWgopk2Dbt38B9ibb0K7dsu8ykhAdSASUAbGjoWtt4ZjjvGWb1lZuNCbcK+0UnYxhFAsb70FO+/sXUa98MIylzqUSgIq2M9ESR0lvShpnKSxks5K09tIel7S+PR31TRdqaPTCZLek7RtoWILtbRwobcoW3llv1chS82aRfIJjce22/oPvtdeg9NOazA3qRaynGIhcJ6ZbQrsAJwmaVPgImComXUBhqbXAPsBXdLjJOC2AsYWauPii/0EuOkmL/oKIRTPL3/p/cXddRfcemvW0dSJgiUgM5tiZm+l53Px3hPaAwcDZZ2Z9sd7ViBNH5AaaQwHWktaq1DxhRoaOBCuv95/fR1zTNbRhNA4XXGF375w9tl+A3c9V5SaWkmdgW2AEUA7M5uSZk0FymrU2gNf5rxtUppWfl0nSRoladSMGTMKF3RY4v33vYPEnXaCv/4162hCaLyaNPH75nbaqUE0tKl0RNS6ImlF4GHgbDObo5wbFs3MJNWoMNPM+gH9wBsh1GWsoQKzZvmAWSut5GPn1HWT6xBCzay0EgwbVrubv0tMQROQpOZ48rnPzB5Jk6dJWsvMpqQitulp+mSgY87bO6RpISuLF3uvBBMn+gG/VpSIhlASGkDygcK2ghNwN/CBmeWW2wwBeqfnvYHBOdOPS63hdgBm5xTVhSxccYX3cH3jjX7JH0IIdaiQV0A7AccCYyS9k6b9HugLDJR0IvA50CvNewrv+HQCMB84voCxhercdRdcfrl3DnrqqVlHE0JogAqWgMzsVaCy68Q9K1jegNMKFU+ogcce8/t99t3XO/FsIJf7IYTSUv+bUYS69dJLcOSR3qXNoEHQvHnWEYUQGqhIQGGJd97xkULXWy+Gqg4hFFwkoOA+/dSL3FZeGZ59FlZbLeuIQggNXCSgAO++C7vu6l2/P/ccdOxY/XtCCGEZRQJq7J55xnvZNfNedjfZJOuIQgiNRCSgxuyOO+CAA3wcnREjYKutso4ohNCIRAJqjBYvht/9Dk4+2Qdve/llaP8/3e6FEEJBRQJqbMaN86Rz7bVwyikweHCMqxNCyEQkoMZi5kw44wzYcksYNQpuu83HFGlW8P5oQwihQpGAGroff4Sbb4YuXeDvf/ceDsaP9+K36OEghJCh+PnbEH32mbdue+YZH7Rq3jzYay+44QbYfPOsowshBCASUP01f74Pjz15sj+++sofY8f6FQ5A585w3HE+ns9ee8UVTwihpEQCqk/MPOnce68PkT137pJ5bdrA2mvDxhvD6ad7rwZdukTSCSGUrEhA9cE333iDgf794ZNPvI+2Xr2809AuXXyguFatso4yhBBqJBJQKVuwwBsQ/PnPMHs27LEHXHopHHZYdBQaQqj3Cjki6j2Spkt6P2daG0nPSxqf/q6apkvSzZImSHpP0raFiqteWLgQ7r7br24uugh69IAxY7xBwbHHRvIJITQIhWyGfS+wb7lpFwFDzawLMDS9BtgP6JIeJwG3FTCu0jVtmg9/vcUW8OtfQ4cOMGwYPP54tF4LITQ4BUtAZvYy8E25yQcD/dPz/sAhOdMHmBsOtJa0VqFiKykLFsCDD8L++3t3OOec41c4gwbBG294L9UhhNAAFbsOqJ2ZTUnPpwLt0vP2wJc5y01K06ZQjqST8KskOnXqVLhIq7NggV+dPPWU328zf743e+7cGdZZx/+usgq0aOGjirZoAU2bepPpjz6Cjz/2v+PHw/ff+9XOBRd4Edumm2b3uUIIoUgya4RgZibJavG+fkA/gG7dutX4/ctsyBDvRfqFFzxxLLecNw5o2xY+/9yvWgYO9HqcyjRt6qOObrgh7L039OwJu+3m00MIoZEodgKaJmktM5uSitimp+mTgdxR0DqkaaVj4UK4+GK4/nq/ujnpJE8cu+76v02gFy3ym0LnzvWucH76acnftdaCddf1K6IQQmjEip2AhgC9gb7p7+Cc6adLegDYHpidU1SXva+/9ntuhg6FU0/1Lm2qSiBNm8aooiGEUI2CJSBJ9wO7AatLmgRciieegZJOBD4HeqXFnwJ6AhOA+cDxhYqrxkaP9vtupk2Df/wD+vTJOqIQQmgQCpaAzOyoSmbtWcGyBpxWqFhqZdEiuP12OO88WGMNePVV6NYt66hCCKHBiOEYKvLKK9C1q/eptuuufhUUySeEEOpUJKBckyfD0UfDLrt4/2sDB3oT67Zts44shBAanOgLDvx+nLvu8lFCFy6ESy6B3/0uurwJIYQCarwJaP58ePhhuPNOL3Jr1swbG/z5z36PTgghhIJqnAnorrvg/PO9h+kNNoC+faF3b1hzzawjCyGERqNxJqBOneDAA73Dz112iUHbQgghA40zAf385/4IIYSQmWgFF0IIIRORgEIIIWQiElAIIYRMRAIKIYSQiUhAIYQQMhEJKIQQQiYiAYUQQshEJKAQQgiZkA/FUz9JmoEPbFcbqwNf12E4da2U44vYaqeUY4PSji9iq53KYlvHzDLv5r9eJ6BlIWmUmZXsID+lHF/EVjulHBuUdnwRW+2UcmwQRXAhhBAyEgkohBBCJhpzAuqXdQDVKOX4IrbaKeXYoLTji9hqp5Rja7x1QCGEELLVmK+AQgghZCgSUAghhEw0qAQk6R5J0yW9nzOtjaTnJY1Pf1dN0yXpZkkTJL0naduc9/ROy4+X1LuAsV0n6cO0/Ucltc6Zd3GK7SNJ++RM3zdNmyDpokLFljPvPEkmafX0uqj7rar4JJ2R9t9YSdfmTM9030naWtJwSe9IGiWpe5pe7GOuo6QXJY1L++isND3zc6KK2DI/JyqLLWd+pudEVfGVwjlRI2bWYB7ALsC2wPs5064FLkrPLwKuSc97Ak8DAnYARqTpbYBP099V0/NVCxTbz4Fm6fk1ObFtCrwLtATWBT4BmqbHJ8B6QIu0zKaFiC1N7wg8i9/su3oW+62Kfbc78B+gZXq9RqnsO+A5YL+c/TUso2NuLWDb9Hwl4OO0fzI/J6qILfNzorLYSuWcqGLflcQ5UZNHg7oCMrOXgW/KTT4Y6J+e9wcOyZk+wNxwoLWktYB9gOfN7Bsz+xZ4Hti3ELGZ2XNmtjC9HA50yIntATP7wcw+AyYA3dNjgpl9amY/Ag+kZes8tuQG4EIgt6VKUfdbFfGdAvQ1sx/SMtNz4st63xmwcnq+CvBVTmzFPOammNlb6flc4AOgPSVwTlQWWymcE1XsNyiBc6KK+ErinKiJBpWAKtHOzKak51OBdul5e+DLnOUmpWmVTS+0E/BfUSURm6SDgclm9m65WZnHlmwI9JA0QtJLkrYrofjOBq6T9CVwPXBx1rFJ6gxsA4ygxM6JcrHlyvycyI2tFM+JcvuulM+JCjUr5sayZmYmqeTanUv6A7AQuC/rWAAkLQ/8Hi8OKVXN8KKNHYDtgIGS1ss2pP86BTjHzB6W1Au4G9grq2AkrQg8DJxtZnMk/Xde1udE+dhypmd+TuTGlmIpqXOigv9rKZ8TFWoMV0DT0uUw6W/ZZelkvDy3TIc0rbLpBSGpD3AAcLSlAtsSiG19vKz4XUkT03bekrRmCcRWZhLwSCr2GAksxjteLIX4egOPpOcP4UUdZBGbpOb4l9R9ZlYWU0mcE5XEVhLnRAWxldQ5Ucm+K+VzomKFqlzK6gF0ZukK4etYusL12vR8f5auOBxpSyoOP8MrDVdNz9sUKLZ9gXFA23LLbcbSlYaf4hWGzdLzdVlSabhZIWIrN28iSypci77fKtl3JwNXpOcb4kUJKoV9h5fJ75ae7wmMzmLfpe0MAG4sNz3zc6KK2DI/JyqLrVTOiSr2XcmcE3l/lmJurOAfBu4HpgA/4b8GTgRWA4YC4/EWIm1y/om34q1AxgDdctZzAl5RNwE4voCxTUgHyTvpcXvO8n9IsX1EalGVpvfEW718AvyhULGVm597shV1v1Wx71oA/wLeB94C9iiVfQfsDIxOJ/QIoGtGx9zOeGX5eznHWM9SOCeqiC3zc6Ky2ErlnKhi35XEOVGTR3TFE0IIIRONoQ4ohBBCCYoEFEIIIRORgEIIIWQiElAIIYRMRAIKIYSQiUhAoSRJWiTvTbrsUeueeiW9XpexFYOk3SQ9Ucm8bSTdXaDtHiDpikKsO4Tyohl2KEmSvjOzFbOOIyuSdgPON7MDKpj3EHCVleuTTFIzW9KRZ223K/wekp3MbP6yrCuE6sQVUKhXJE2UdLmktySNkbRxmt5WPrbNWEl3Sfo8Z7yW79Lf3SQNkzQojZlyX/rCRVLX1IHjaEnPlnVVU27bR0h6X9K7kl5O0/pIGpzWO17SpTnLHyNpZLqCu0NS0zT955LeSJ/hodSnV9nYLB9Kegs4rJLPvxKwZVnykXSZpH9Keg34Z2XxSOqc1n2vpI/TZ99L0mtpue7gfcMBw/CucEIoqEhAoVQtV64I7pc58742s22B24Dz07RLgRfMbDNgENCpkvVug3cuuSk+DspOqV+tW4DDzawrcA9wdQXv/T9gHzPbCjgoZ3p34BfAlsARkrpJ2gT4JX4lsTWwCDg6JcU/AnulzzAKOFdSK+BO4ECgK7BmJfF3w+90z7VpWt9RlcWTpm8A/AXYOD1+hd9Vfz7e0WaZUUCPSrYfQp1pVL1hh3plQfrirkhZ54ujWXKlsDNwKICZPSPp20reO9LMJgFIegfvx20WsDnwfLogaop3r1Pea8C9kgbmxAA+5svMtM5HUiwL8UTyZlrncninnzvgCeO1NL0F8AaeED4zs/FpPf8CTqoghrWAGeWmDTGzBdXE81ha/5g0fSww1MxM0pi0H8pMB9auYNsh1KlIQKE++iH9XUTNj+Efcp6XvV/AWDPbsao3mtnJkrbHO58cLalr2azyi6Z19jezi3NnSDoQTxBHlZu+dZ7xLwBalZs2r4LtV/Q697Mvznm9mKX3Y6u0nRAKKorgQkPxGtALvI4F7304Xx8BbSXtmN7fXNJm5ReStL6ZjTCz/8OvQsq6st9bUhtJy+Gji76Gd/Z5uKQ10nvbSFoHH+VzJ0kbpOkrSNoQ+BDoLGn9tM6lElSOD/CitKpUFE9NbMj/FvOFUOciAYVSVb4OqG81y18O/FzS+8AR+Eifc/PZkPlwxIcD10h6F+9d+GcVLHpdavjwPvA63ts1wEh8bJb3gIfNbJSZjcPrep6T9B4+HPNaZjYD6APcn6a/AWxsZt/jRW5PpkYI06mAmX0IrJIaI1Tmf+LJZz/k2B14sobvCaHGohl2aBAktQQWmdnCdCVzWxV1SHW53T549/unF3pbOds8B5hrZnfVdTyS2gH/NrM9ly3KEKoXdUChoeiED0HcBPgR+E3G8RTSbfhVXiF0As4r0LpDWEpcAYUQQshE1AGFEELIRCSgEEIImYgEFEIIIRORgEIIIWQiElAIIYRM/D/3Qn62rBpg9QAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAY8AAAEJCAYAAABsc6siAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuNSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/xnp5ZAAAACXBIWXMAAAsTAAALEwEAmpwYAABOmElEQVR4nO3dd3xUVfrH8c+TnpCeQGih9xKqCIJCLNiwYUFUFkXFsupadl11d3V19be2XfuuFcuuiiIWYG2ooNjovSb0JJBQQ0IIac/vjxlwiKlTMpPwvF+veeXObfOdG24O955zzxFVxRhjjKmPIH8HMMYY0/hY4WGMMaberPAwxhhTb1Z4GGOMqTcrPIwxxtRbiL8DeEtycrJ26NDB7e0PHjxIs2bNvBfIiyybewI5GwR2PsvmvkDOV1W2xYsX71bV5vXemao2idegQYPUE3PmzPFoe1+ybO4J5GyqgZ3PsrkvkPNVlQ1YpG78zbXbVsYYY+rNCg9jjDH1ZoWHMSbgHDxcxj+/XM/5z39P1r4if8cxVWgyFebGmMavokKZviSLJ75YT17BYUKDhd9PW8471w0lKEjqtI/S0lKysrIoLi72cdqqxcXFsXbtWr98dk0iIiIQqdsxrAsrPIwxAWHB5r38bdYaVmbn0y81nn9fNZCNeQe5e/oKpvywmetO7lSn/WRlZRETE0OHDh28+seyrgoKCoiJiWnwz62JqrJnzx6vtgKzwsMY43ffbdjFb6YsoFVcBE+P68/5/VoTFCQMbJfAl2tyefyL9ZzSrTndUmr/o1xcXOy3giNQiQhJSUls377da/u0Og9jjF+pKk98sZ7UxEi+vmskFw5oc/QWlYjw97F9iQkP4Y73llFSVlGnfVrB8WvePiYeFR4ikiAivUWkk4hYQWSMqbcv1+SyMjuf207tSlTYr2+GNI8J55GL+rI65wDPfZPhh4SmKvX+gy8icSJyn4isBH4GXgLeB7aKyDQRSfd2SGNM01RRoTw1ewMdk5tx0YA21a53Vp+WXDywLS/MySRzf3kDJnTPI488Qu/evUlLS6N///7Mnz8fgA4dOrB79+567Wv8+PGkpaXx1FNP8fTTT1NUVHXrs+eff54uXbogIvX+DHe4U+fxAfAWcLKq7nddICKDgAki0klVX/NCPmNME/bZqp2s21nA0+P6ExJc8/9lHzi/Fz9v2sMrK4qZeF4FobWs7y/z589n1qxZLFmyhPDwcHbv3k1JSYlb+9q5cycLFy4kMzMTcBQ+V111FVFRUb9ad/jw4YwZM4ZRo0Z5Er/O6n30VfUMVf1P5YLDuWyxqt5uBYcxpjblFcpTX22ga4tozuvXutb1YyNC+fO5PcktUn7I9P3/rN2Vm5tLcnIy4eHhACQnJ9O69S/f77nnnmPgwIH07duXdevWAbBgwQKGDRvGgAEDOOmkk1i/fj0Ao0ePJjs7m/79+/Pggw+Sk5NDeno66em/vsEzYMAAPOnfr748am0lImlAB9f9qOqHHmYyxhwHZi7PITOvkBeuGEhwHZ/hOLVnC6JCYMayHEZ1b1Hr+g/OXM2anAOeRj1Gr9axPHBe72qXn3rqqTzxxBN069aN008/nXHjxjFy5Mijy5OTk1myZAn/+te/ePLJJ3n11Vfp0aMH8+bNIyQkhK+++or77ruP6dOnM2PGDMaMGcOyZcsAeP3115kzZw7Jycle/U7ucLvwEJEpQBqwGjjSBEIBKzyMMTUqK6/gma8z6NEyhrP7tKzzduEhwQxuGcIXq3dyqKScyLBgH6Z0T3R0NIsXL2bevHnMmTOHcePG8eijj3L11VcDMHbsWAAGDRrEhx86/lzm5+czceJEMjIyEBFKS0v9Fb/OPLnyGKqqvbyWxBhz3PhwaTabdx/k5QmD6vzk+BHDWoXwXVYxX6/LZUxazbe7arpC8KXg4GBGjRrFqFGj6Nu3L2+++ebRwuPI7azg4GDKysoA+Mtf/kJ6ejofffQRW7ZsabB6C094UuP0k4hY4WGMqZeSsgqe/TqDvm3iOKNXSr23754YREpsOJ8sy/FBOs9lZGSQkfFLk+Jly5bRvn37GrfJz8+nTRtHa7M33nij2vViYmIoKCjwSk5PeVJ4vIWjAFkvIitEZKWIrPBWMGNM0/Thkiyy9h3ijjO6uvXgWpAI56W1Zu76PPKLAu/2TmFhIRMnTqRXr16kpaWxZs0a/vrXv9a4zd133829997LgAEDjl6NVGXy5MmcddZZVVaYP/vss7Rt25asrCzS0tK47rrrPP0qNfLkttVrwARgJb/UeRhjTLXKyiv419yNpLWNI70OFd7VuaB/G179fjOfrdrB5UPaeTGh5wYMGMCPP/5Y5bItW7YcnR48eDBz584FYNiwYWzYsOHosocffhhwNM1dtWrV0fm33nort956a5X7vu2227jttts8TF93nlx57FLVGaq6WVW3Hnl5LZkxpsn5ZFkO2/YWcUt6F4+6y+jTJpZOyc0C9tbV8cCTK4+lIvIOMBM4fGSmNdU1xlSlvEJ5YU4mPVrGuFXX4UpEOL9/a575OoOd+cW0jIvwUkpTV55ceUTiKDRGA+c5X2O8EcoY0/T8b+UONu0+yK2nulfXUdn5/VqjCrNW/PrqwzE0t3Hl7WPi9pWHql7jzSDGmKarokJ5/psMurSIrtdzHTXp1DyatLZxfLIs55ixPiIiItizZw9JSUnWu67TkfE8ysu91y9YvQsPEbm/hsWqqn/zII8xpgn6cs1ONuQW8vS4/vV+rqMm5/drzcP/W8vGXYV0bh4NcLTF0a5du7z2OfVRXFxMRETg3UaLiIjg4MGDXtufO1ceVX16M+BaIAmwwsMYc5Sq8tw3mXRIimJMWiuv7vu8fq155NO1zFiWwx1ndAMgNDSUjh07evVz6mPu3LkMGDDAb59fk61bvdemyZ2OEf9x5AW8jKPu4xpgKlC3cSKNMceNb9blsTrnADend6m159z6SomNYFinJD5Zlm31HA3Mrd+kiCSKyMPAChxXLwNV9Y+qmufVdMaYRk1VefabTNomRNY4Xocnxg5sy5Y9Rfy8aa9P9m+q5s5gUE8AC4ECoK+q/lVV99Vj++4isszldUBEbq+0zigRyXdZp6Z6FmNMgPpgcRbLt+/nlvQuPht/Y0xaK+IiQ3l7vj1m1pDcqfO4C0cT3T8Df3JpzSA4Ksxja9pYVdcD/QFEJBjIBj6qYtV5qmpNf41ppHbmF/PQrDUM6ZDIZYNTffY5EaHBXDKoLW/9tIVdBYdpHhPus88yv3CnziNIVSNVNUZVY11eMbUVHFU4DdhoT6Yb07SoKvd9tJLS8goevyTNqy2sqnLFie0oLVfeX7Tdp59jfiHuVjKJyLWVRwwUkUdV9Z567GMKsERVn680fxQwHcgCcoDfq+rqKrafDEwGSElJGTR16tT6fo2jCgsLiY6Odnt7X7Js7gnkbBDY+TzN9kN2Ka+sLGF8jzDO7BDqxWTVZ3tswSHyipQnRkYS5MfnOxrb7zU9PX2xqg6u985U1a0X8Clwpcv7F4Ap9dg+DNgNpFSxLBaIdk6fA2TUtr9BgwapJ+bMmePR9r5k2dwTyNlUAzufJ9ly8w9p3wc+14v/9YOWlVd4L5RTddlmLc/R9n+cpd+szfX6Z9ZHY/u9AovUjTLAkxqsi4GrRWS8iLwJlKnqpHpsfzaOq47cygtU9YCqFjqnPwVCRcT/4y4aY2qkqtz30SoOlzluV9V1eFlvGN07heYx4fz3Z7sL3hDcaW2VKCKJOJ7vuA64G0fLqwed8+tqPPBuNZ/RUpw18SIyxJlzT32zGmMa1ozlOXy1Npffj+5Op+YNe+smNDiIcYNT+WZ9Hln7ihr0s49H7lx5LAYWOX/OAeKBc13m10pEmgFn4DLeuYjcKCI3Ot9eAqwSkeXAs8DlzssrY0yAWp2Tz18+XsWAdvFMGuGfJ7zHn9gOAaYusIpzX6t3U11V9fhfhaoexNGVieu8F12mnweer7ydMSYwrck5wJWvzic6PIRnLx/QoLerXLWJjyS9ewumLtzO707v6rNnS4x7t61G1LI8VkT6uB/JGNOYOAqOn4kMDWbq5GGkJkb5Nc9VQ9uzu/AwX67+VXWq8SJ3HhK8WEQeBz7HcatqFxABdAHSgfY4HiQ0xjRxa3c4Co6I0GCmTh5KuyT/FhwAp3RrTpv4SN76aQvnerkjRvMLd25b3eGsGL8YuBRoBRwC1gIvqer33o1ojAlE63Y6blWFhwTz7vVDaZ/UzN+RAAgOEiaN6MjfZq3h+4zdjOhqDTV9wa3BoFR1L/CK82WMOY6oKtMWZ/G3mWtoFh7C1MlD6ZAcGAXHEVcNbcfrP2zm75+tZWbnET5/wv14ZLVJxpg6yz1QzKQ3FnL3Byvo2SqWaTcOC7iCAyA8JJg/nNmd1TkHmLH818PUGs9Z4WGMqZWq8uGSLM7457f8tGkPD5zXi6mTh/q9crwm56W1pk+bWJ74Yj3Fpd4bftU4WOFhjKnRzvxirn9rEXe+v5xuKTF89rtTuGZ4x4C/FRQUJNx7dk+y9x/iPz/ZU+fe5ladB4CIROFoVdVOVa8Xka5Ad1Wd5bV0xhi/UVXeW7idRz5dS2l5BX8+tyfXDO/ot2c43DG8SzIjuzXn+TmZXDY4lbgo73bSeDzz5MrjdRzjegxzvs8GHvY4kTHG73YVVTDhtQXc8+FKerWK5fPfncJ1J3dqVAXHEfec3YMDxaX8a26mv6M0KZ4UHp1V9XGgFEBVi3AMCGWMacTeXbCNP/1wiKXb9vHwhX149/rAa01VHz1bxTJ2QFte/3EL2fsP+TtOk+FJ4VEiIpGAAohIZxxXIsaYRqi8Qnl41hru/XAlXeOD+PLOkVw1tH3A123UxV2juwHw6Gfr/Jyk6fCk8HgAx1PmqSLyNvA1jh52jTGNzKGScm5+ezGvfr+ZicPac9fgCNrER/o7lte0jo/klvQuzFyew/sLrdNEb3C78FDV2cBY4GocXasPVtW53olljGkoeQXFXP7yT3y5Jpf7x/TiwQv6+HUkPl/5bXoXRnRJ5i+frGJ1Tr6/4zR67nSMOPDIC0c/VjtwDBXbzjnPGNNIbMgt4KIXfmRDbiEvXTXIb12pN4TgIOHpy/sTHxXKzW8v4UBxqb8jNWruXHn8w/l6AZgPvIyjm5L5znnGmEZg7vo8Lv7Xj5SUV/DeDUMZ3bulvyP5XHJ0OC9cMZCsfYf4w7Tl2DBB7qt34aGq6aqajuOKY6CqDlbVQcAAHM11jTEBTFV544fNTHpjIW0To/jkt8NJaxvv71gNZnCHRO49uwdfrM7lte83+ztOo+X2Q4I4HghceeSNqq4SkZ5eyGSM8ZHS8goenLma//68jdN7pvDM5f1pFu7Jn4HG6doRHVm4ZS9//2wd/VLjOaFDfUbQNuBZa6sVIvKqiIxyvl4BVngrmDHGu/KLSrnm9YX89+dt3HBKJ16aMOi4LDgARIQnLu1HakIkk15fyMIte/0dqdHxpPC4BlgN/M75WuOcZ4wJMF+u3skZT33L/M17ePySNO49p2ejfFrcm2IjQnn7+qE0jwlnwmvzmbs+z9+RGhVPmuoWq+pTqnqR8/WUqhZ7M5wxxjO7Cg7z23eWMPk/i0lsFsb0m07issGp/o4VMNrER/L+jcPo3Dya699axKwV1n17XXnSMeJmnE+Xu1LVTh4lMsZ4TFX5aGk2D81aQ9Hhcn4/uhs3jOxMaLB1pF1ZcnQ4704eynVvLOLWd5dSUFzG+CHt/B0r4Hlyw3Owy3QEjiFp61TrJCJbgAKgHChT1cGVlgvwDHAOUARcrapLPMhqzHGhuLScmctzmPLDFtbuOMCg9gk8dnFfurSI8Xe0gBYbEcqbk4Zw89uLuffDlezML+a207oe97f2auJ24aGqeyrNelpEFgP313EX6aq6u5plZwNdna8TgX87fxpjqpBXUMzbP2/j7flb2V1YQreUaB6/OI1LBrVtEn1TNYTIsGBemjCYez9cyTNfZ7Bk2z6eGtef5Ohwf0cLSJ7ctnJ9mjwIx5WIt5puXAC8pY4neH4WkXgRaaWqO7y0f2M8VlxazqrsfDbuKiQzr5CNuw6ycVch8ZGhjOzWnJHdW9A/Nd4n/3stLi1nydZ9/Lx5L/M37WHJtn2Uliun9WjBpBEdOalzEtIEuxjxtbCQIJ68NI0TOiTwwIzVnPPMPJ65fADDOif5O1rAEXefsBSROS5vy4DNwD9UdX0dtt0M7MNRZ/KSqr5cafks4FFV/d75/mvgj6q6qNJ6k4HJACkpKYOmTp3q1ncBKCwsJDo62u3tfcmyucdX2VSVRbnlvL22hP2HHedPSBC0ahZEy2bCvmJl4/4KFGgWCr2TgumTHEzPxGCaR/1S51DXfKUVSk5hBdsOVLCtoIKtByrYtL+CMnWMgdA+NoheScGc0jaEls28U6dxPP5eK9teUMELS4vJLVIu6hrKmE6hderzq7Edu/T09MWVqw7qwpMrhWtVdZPrDBGpa8c4I1Q1W0RaALNFZJ2qflffAM5C52WAwYMH66hRo+q7i6Pmzp2LJ9v7kmVzjy+ybd1zkPs/Wc23G3bRq1Usj57WlZ6tYmibEHXMFcb+ohLmZexm7vpdfLthFwt2OkYraJsQybBOSQzrnETunvUM7pNGaHAQIUFCcJCQe6CYbXuL2LaniG17i9i6p4iNuwopq3AUUpGhwfRsFcu1fRMZ2imJQR0SiI3w/uh4x9vvtToXjS7jvg9X8uHyHPYFxfPM+AG1Hu/j5dh5Unh8AFTuCPEDYFBtG6pqtvNnnoh8BAwBXAuPbMC1PWFbrOsT40eHy8p55btNPPdNJiFBwv1jevGbYe0Jqab1UnxUGOf1a815/VqjqmTkFfJj5m5+2rSHL9fkMm1xlmPFhT9VuX1EaBDtEqNolxjF6b1a0LNVLL1axdI+qZlV4jag6PAQnrm8Pyd0TOTBGau58PkfePk3g+nSIjCvLBpSvQsPEekB9AbiRGSsy6JYHK2uatu+GRCkqgXO6dHAQ5VWmwHcIiJTcVSU51t9R+OTV1DMjGU5fLwsmyARJg3vyLlprRpdc9FFW/Zyz4crycwr5Ny+rfjLmF60jKv1n/pRIkK3lBi6pcRw9fCOlFcoG3ILmPvjQvqkpVFWrpSUV1BeobSICaddUhTNo8OtziJAiAgThranW4tobn57CRe98ANPX96f03qm+DuaX7lz5dEdGAPEA+e5zC8Arq/D9inAR84TIwR4R1U/F5EbAVT1ReBTHM10M3E01bUn1xuJ4tJyvli9kw+XZDMvYxcVCv3axlFUUs7t7y3jiS/Wc+2Ijow7ITXgu8YoKC7lsc/X8d+ft9EmPpLXrz6B9B4tPN5vcJDQs1UsuUnBnNy1uReSmoZwYqckZtw6ghv+s4jr3lrE70d35+ZRnY/bQr7eZ6+qfgJ8IiLDVLXqa+6at98E9Kti/osu0wr8tr77Nv61bPt+7nhvGZt3H6RNfCQ3jerMRQPa0qVFNBUVypz1ebz07SYemrWGZ77O4JrhHbh2REdifHDP3lNfrt7J/Z+sJregmEnDO3LX6G4BX9gZ32sTH8m0G07ing9X8MQX68nad4hHLuxzXDaHdue21d2q+jhwhYiMr7xcVW/zSjLTaJSVV/DCnI08+00GKTHhvH71CYzs1vyYEyooSDitZwqn9Uxh8dZ9vPjtRp7+KoM3ftzCTSM785thHYgMC/bjt3C0opq7YRevztvED5l76NEyhhcnDKJ/arxfc5nAEhkWzNPj+tM2IZIX5mxEVfm/i/oedwWIO/+VWuv8uajGtcxxIa+ogktf+oml2/ZzQf/WPHRBH+Iia76SGNQ+gVd+M5iVWfk8+eV6/v7ZOl79fjO3ntqFywanEhHasIXI4bJyPlmawyvzNpGRV0hKbDh/OqcnVw/v0OjqZ0zDEBF+P7o7wSI8+00mFao8OjbtuCpA3LltNdP5803vxzGNRUWF8u7CbTz0wyHCQkt55vL+XNC/Tb320bdtHG9OGsLCLXt54ov13P/Jah7531pO6JDI8C7JDO+SRO/WcT57yG7+5r3MWZfHrBU72F14mJ6tYvnnZf0Yk9aasBArNEzNRIQ7R3dHRHjm6wwqFB67OM3fsRqMO7etZlJFh4hHqOr5HiUyAW/z7oPcM30F8zfvpWdiEK9NPoXW8ZFu7++EDom8N3koP2/ay+w1ufy4cTePfb4OgLjIUPqlxtOzVQy9nM1VOyY3q7aJbHWKS8vJzCtk6fb9zF2Xx48b93CotJyI0CBGdGnONcM72FPZxi13nNENEXj6qwxU4dzmx8fQtu7ctnrS6ylMo1BWXsFr32/mn7M3EBYSxGMX96VF4UaPCo4jRIRhnZOOdgORV1DMTxv38GPmHlZm5/P6xj2UlFcAji4kUhMiaZMQRduESNrER9I6PgJBKCmroKS8gpKyClZtLGH6jqWs3XGAzbsPUu580C41MZLLBrdlVI8WDOuU1OC3yUzTc/vp3RCEp77aAF1DOTXd34l8z53bVt8emRaRMKAHjiuR9apa4sVsJoCsys7nvo9WsiIrnzN6pfDwhX1IiY1g7txNtW/shhYxEVzQv83RW2Gl5RVs3FXI2h0HWLejgG17i8jef4hV2fnsPVj9P7s28fvo2SqWs/u0pEfLWHq3jqV9UpRdYRiv+93pXVm38wCfrN7JTXkFTb4nY086RjwXeBHYiKOLnY4icoOqfuatcMb/8otKefLL9bw9fysJUWE8f8UAzu3bqsH/+IYGB9GjZSw9WsbCgGOXFZWUsSO/GMFxVRIWHERYSBALfvqB0acdB/8FNAHjwQt68936ndz9wQqm3XhSk+4NwJOG6//A0a16JoCIdAb+B1jh0QRUVCgfLM7i0c/Xsb+ohAlD23PnGd2Jiwq8ZzKiwkLo3PzX3UWEBTfdE9cEphYxEVzRI4xXVu7nzR+3MGlEXbv7a3w8KTwKjhQcTptwPGVuGrmfN+3hsc/XsXTbfga3T+DBC4bQu3Wcv2MZ0yic1DqEjJI4nvhiPaf3TKFdUpS/I/mEJ+0RF4nIpyJytYhMBGYCC0VkbKU+r0wjoKr8kLmby176ictf/pmsfYf4x6X9mHbjMCs4jKkHEeH/LupLcJBw70crcHfYi0DnyZVHBJALjHS+3wVE4ujvSoEPPYtmGkJFhfJtxi6e+zqDJdv2kxIbzgPn9WL8kHbWCskYN7WOj+Ses3vw549X8d7C7VzeBMdE92QYWuussJFSVVZlH2DG8mxmrdjBjvxiWsdF8LcL+3DpoLZWaBjjBVcMacfM5Tk88r+1nNKtuVeatAcST1pbdQRuBTq47sceEgwspeUV7MwvJnv/IXL2HyIzr5DPV+1k0+6DhAQJI7s1556ze3B2n1b2VLUxXhQUJDx2cRrnPjuPW95Zwns3DGtS3d14ctvqY+A1HHUdFV5JY+qkrLyCrXuL2OgcNzszr5DNuwspPFxGWblSWlFBuXOMiD0HS3C95SoCQzsmcf0pnTi7T0vio8L890WMaeI6JDfj0YvTuPXdpTzxxXruO6envyN5jSeFR7GqPuu1JKZaxWXKvIxdLNyyj4Wb97J0+z6KS38pr1vEhNOpeTM6xji67QgNEsfP4CCax4TTJj6CNvFRtI6PoHV8pN2WMqYBndevNfM37+Hl7zYxpEMip/dqGoNIeVJ4PCMiDwBfAoePzFTVJR6nakQOHi5j78ES9hWVsK+olH0HSzhcVk7zmHBaxETQIjacpGbhtT4spKoUl1aQvf8QmXkFZOQWkpFXyIbcAjbkFlGhCwgS6NkqlstPaEffNnF0bhFNp+bNfDKGtTHGe/58bi+Wbd/PXdOWM+vWEaQmNv7mu54UHn2BCcCp/HLbSp3vm7TdhYf5dOUOPlmWw+Kt+2pdPzhIiI8MJTzE8eRzWEgQ4SHBVKhSUFxGQXEpBcVllFUc26SvbUIkXVtE0zWqmEtGDWBgu/iAHDjJGFOziNBgXrhiIGOe/Z5b3l3KtBuGNfo6Rk8Kj0uBTsdLf1Z7D5YwZ10eM5bn8H3mbsorlB4tY7j99K60jo8kISqMhKhQ4qPCCA8JYlfhYfIOHCavoJjcA8XsKyp1dNrnfB0uK0dE6NIihJiIEGIiQomJCKFFTATdUqLp3Dz66Mh1c+fOZWQ3G67UmMasfVIzHr8kjZveXsLfP1vLA+f19nckj3hSeKzCMY55nneiBJaSsgqWbNvHvIxdfLdhN6ty8lF1XA3cOLIT5/drQ/eW1Xd81hQuS40x3nV231ZcM7wDr/+whbYJUVzbiLsv8aTwiAfWichCjq3zaHRNdQ8eLiNzXznbf97K2h0Hjvbceqi0nJAgYWC7BO48vRundGtOWts465HVGOO2+87pyc78Yv42aw0Cjbb/K08Kjwe8lsKPsvYVcfLjc5zNWVcRGxHiqJQeksrQTkmc1DnJ6hmMMV4TGhzEs+MHcOs7S3lo1hqgcRYgnjxh/q3rexEZAYwHvq16i6PrpQJvASk4KthfVtVnKq0zCvgE2Oyc9aGqPuRu1pq0jovkrjO6Ubp7K5eecRJt4iPtysIY41OhwUE8d8UAbnlnCQ/NWoMIXDO8cRUgnlx5ICIDgCtwVJ5vBqbXYbMy4C5VXSIiMcBiEZmtqmsqrTdPVcd4kq8ugoKEW07tyty52bRNsHoKY0zDCA0O4vkrBnLLO0t4cOYaKhQmDe/QaP7zWu+2YiLSTUQeEJF1wHPANkBUNV1Vn69te1XdceRZEFUtANYCbeqbwxhjGrvQ4CCeGz+QM3un8LdZa/jNlAVs2lXo71h1IvXtLlhEKoB5wLUuA0FtUtVO9f5wkQ7Ad0AfVT3gMn8UjquYLCAH+L2qrq5i+8nAZICUlJRBU6dOrW+EowoLC4mO/vWAQoHAsrknkLNBYOezbO5zJ195hfLNtjI+zCyhpBzO6hDK+Z1DCQ/x7lVIVdnS09MXq+rgeu9MVev1Ai4EpgLbgVeA04DNbuwnGlgMjK1iWSwQ7Zw+B8iobX+DBg1ST8yZM8ej7X3JsrknkLOpBnY+y+Y+T/LlHjikd763TNv/cZYO/b+v9INF27WguNSn2YBFWs+/36pa/9tWqvqxql4O9ADmALcDLUTk3yIyui77EJFQHFcWb6vqr8b9UNUDqlronP4UCBWR5PpmNcaYxqRFTAT/uKwfH9w4jISoMO6atpyBD81m4pQF/PfnreQeKPZ3xKM8aW11EHgHeEdEEnBUmv8RR19X1RJHbdBrwFpV/Wc167QEclVVRWQIjrqZPe5mNcaYxmRwh0Rm3jqCRVv2MntNLrPX5vLnj1fx549X0T0lhu4tYxxdF6XE0C0lmnaJUYQ0cHfvHrW2OkJV9wEvO1+1GY6jT6yVIrLMOe8+oJ1zXy8ClwA3iUgZcAi43Hl5ZYwxx4XgIOHETkmc2CmJP53bk4y8QmavyWXRlr0s3rqPGctzjq7bu3Us/7vt5AbN55XCoz5U9XugxlogdbTaqrXlljHGHA9EhG4pMXRL+aVLpIOHy8h09rwdEtzwzXsbvPAwxhjjuWbhIfRLjadfarxfPr9x9wlsjDHGL6zwMMYYU2/1fkgwUInILmCrB7tIBnZ7KY63WTb3BHI2COx8ls19gZyvqmztVbXeAwY1mcLDUyKySN15yrIBWDb3BHI2COx8ls19gZzPm9nstpUxxph6s8LDGGNMvVnh8Yu6PODoL5bNPYGcDQI7n2VzXyDn81o2q/MwxhhTb3blYYwxpt6s8DDGGFNvTbbwEJEpIpInIqtc5iWKyGwRyXD+THDOFxF5VkQyRWSFiAx02Waic/0MEZnow2xPiMg65+d/JCLxLsvudWZbLyJnusw/yzkvU0Tu8Ua26vK5LLtLRPRIF/mBcOyc8291Hr/VIvK4y/wGO3bV/F77i8jPIrJMRBY5e4n2x3FLFZE5IrLGeYx+55zv93OihmwBcU5Ul89lud/OiZqy+fyccGcQkMbwAk4BBgKrXOY9DtzjnL4HeEx/GXDqMxwdNg4F5jvnJwKbnD8TnNMJPso2GghxTj/mkq0XsBwIBzoCG4Fg52sj0AkIc67Ty1fHzjk/FfgCx8OYyQF07NKBr4Bw5/sW/jh21WT7Ejjb5VjN9dNxawUMdE7HABucx8fv50QN2QLinKguXyCcEzUcO5+fE032ykNVvwP2Vpp9AfCmc/pNHKMiHpn/ljr8DMSLSCvgTGC2qu5VR7fzs4GzfJFNVb9U1TLn25+Bti7ZpqrqYVXdDGQCQ5yvTFXdpKolOEZ3vMDTbNXlc3oKuBtwbWXh92MH3AQ8qqqHnevkuWRrsGNXTTbFMTImQByOYZWPZGvI47ZDVZc4pwuAtUAbAuCcqC5boJwTNRw78PM5UUM2n58TTbbwqEaKqu5wTu8EUpzTbXAMq3tElnNedfN9bRKO/7kETDYRuQDIVtXllRYFQr5uwMkiMl9EvhWREwIo2+3AEyKyHXgSuNff2USkAzAAmE+AnROVsrkKiHPCNV+gnROVjp3Pz4njtkt2VVURCbh2yiLyJ6AMeNvfWY4QkSgcA3bVaZhhPwjBcStgKHAC8L6IdPJvpKNuAu5Q1ekichmOUTRP91cYEYnGMQT07ap6QOSXcSD8fU5UzuYyPyDOCdd8zjwBc05U8Xv1+TlxvF155DovH3H+PHIpl43j3uURbZ3zqpvvEyJyNTAGuFKdNygDJFtnHPdHl4vIFudnLRHHcMGBkC8L+NB5m2ABUIGjA7hAyDYR+NA5PQ3H7QH8kU1EQnH8gXlbVY9kCohzoppsAXNOVJEvYM6Jao6d78+J+lbQNKYX0IFjKy+f4NjKwced0+dybAXXAv2lgmszjsqtBOd0oo+ynQWsAZpXWq83x1ZwbcJRuRXinO7ILxVcvX117Cot28IvlYOBcOxuBB5yTnfDcfkt/jh2VWRbC4xyTp8GLPbHcXN+zlvA05Xm+/2cqCFbQJwT1eULhHOihmPn83PC45MlUF/Au8AOoBRHKXwtkAR8DWTgaImQ6PILeAFHa4OVwGCX/UzCUamUCVzjw2yZzl/wMufrRZf1/+TMth5nyx3n/HNwtK7YCPzJl8euhhMlEI5dGPBfYBWwBDjVH8eummwjgMXOk3E+MMhPx20EjkrdFS7/xs4JhHOihmwBcU5Uly8Qzokajp3PzwnrnsQYY0y9HW91HsYYY7zACg9jjDH1ZoWHMcaYemsyz3kkJydrhw4d3N7+4MGDNGvWzHuBvMiyuSeQs0Fg57Ns7gvkfFVlW7x48W51Ywxzj1siBMpr0KBB6ok5c+Z4tL0vWTb3BHI21cDOZ9ncF8j5qsoGLFI3/ubabStjjDH1ZoWHMSYglZZXsCo7/8gzCCbAWOFhjAk4ZeUV3PrOUsY89z13f7CC4tJyf0cylTSZCnNjTNNQXqHcNW05n6/eyak9WjBtcRYZeYW8NGEQKbERtW5fWlpKVlYWxcXFDZD21+Li4li7dq1fPrsmERERuHaE6SkrPIwxAaOiQrnvw5V8siyHP5zZnd+md+HzVTu48/3ljHnue168alCt+8jKyiImJoYOHTp49Y9lXRUUFBATE9Pgn1sTVWXPnj1ebQVmt62MMQFBVXlw5mreW7Sd207twm/TuwBwVp9WfHTzcCJDgxn/8s/MyyqtcT/FxcUkJSX5peAIVCJCUlISwcHBXtunFR7GGL9TVR79bB1v/rSV60/uyB1ndDtmefeWMcy4ZThDOiby2qoS1u08UM2eHKzg+DVvH5M6FR4ikiAivUWkk4hYgWOM8arZa3J56btNTBjanvvO6VnlH7r4qDCev2IA4cHw0reb/JDSuKq2IBCROBG5T0RW4hg/+CXgfWCriEwTkfSGCmmMabpUlX/N3Ui7xCgeOK9Xjf9Djo8KY1TbEGYszyFrX1EDpqyfRx55hN69e5OWlkb//v2ZP98xqm6HDh3YvXt3vfY1fvx40tLSeOqpp3j66acpKqr6e1955ZV0796dPn36MGnSJEpLa76956mariI+wNGX/smq2l1VR6jqYFVNBR4FLhCRa32azhjT5C3cso9l2/dz/ckdCQmu/cbGmR1DCRJ4dd7mBkhXf/Pnz2fWrFksWbKEFStW8NVXX5Gamlr7hlXYuXMnCxcuZMWKFdxxxx21Fh7r1q1j5cqVHDp0iFdffdWTr1GraltbqeoZNSxbjGOAG2OM8ciL324ksVkYlwyq2x/YxIggLuzfhqkLt3HrqV1Iig6vdt0HZ65mTU7N9SP11at1LA+c17va5bm5uSQnJxMe7siVnJx8zPLnnnuOmTNnUlpayrRp0+jRowcLFizgd7/7HcXFxURGRvL666/TvXt3Ro8eTXZ2Nv379+eiiy4iJyeH9PR0kpOTmTNnzjH7Peecc45ODxkyhKysLC9+61+ra51HmoicLyJjj7x8msoYc1xYv7OAb9blcfVJHYgMq3tLoBtGdqK4tII3f9ziu3BuOvXUU9m+fTvdunXj5ptv5ttvvz1meXJyMkuWLOGmm27iySefBKBHjx7MmzePpUuX8tBDD3HfffcBMGPGDDp37syyZct44IEHaN26NXPmzPlVweGqtLSU//znP5x11lm++5LU4TkPEZkCpAGrcQyiDo5hDz+sdiNjjKmDl7/bRGRoMBOGtq/Xdl1axDC6Vwpv/rSVG0Z2pll41X/KarpC8JXo6GgWL17MvHnzmDNnDuPGjePRRx/l6quvBmDsWMf/vQcNGsSHHzr+jObn5zNx4kQyMjIQEY/qK26++WZOOeUUTj75ZI+/S03qcuUx1FnXMVFVr3G+JtVl5yJyloisF5FMEbmniuVPicgy52uDiOx3WVbusmxG3b+SMaYx2JF/iE+WZTPuhFQSmoXVe/sbR3Um/1Ap7y7Y5oN0ngkODmbUqFE8+OCDPP/880yfPv3osiO3s4KDgykrKwPgL3/5C+np6axatYqZM2e6/XT8gw8+yK5du/jnP//p+ZeoRV0Kj59EpFd9dywiwTgGgT8b6AWMr7wfVb1DVfuran/gOY69mjl0ZJmqnl/fzzfGBLbX5m1GgWtHdHRr+4HtEhjaKZFX522mpKyi9g0aSEZGBhkZGUffL1u2jPbta76yys/Pp02bNgC88cYb1a4XExNDQUFBlcteffVVvvjiC959912Cgnz/REVdPuEtHAXIehFZISIrRWRFHbYbAmSq6iZVLQGmAhfUsP544N067NcY08jlFzmuGMaktSI1Mcrt/dw4sjM7DxTz8bJsL6bzTGFhIRMnTqRXr16kpaWxZs0a/vrXv9a4zd133829997LgAEDjl6NVGXy5MmcddZZpKf/+kmJG2+8kdzcXIYNG0b//v156KGHPP0qNZLaujsWkUzgTmAlv9R5oKpba9nuEuAsVb3O+X4CcKKq3lLFuu1xPEvSVlXLnfPKgGVAGfCoqn5cxXaTgckAKSkpg6ZOnVrjd6lJYWEh0dHRbm/vS5bNPYGcDQI7n6+zzdpYwgcZpTx0UgTtYuvXZYZrNlXl/h+LKa1Q/m9EJEEixMXF0aVLF1/ErpPy8nKvdgPiTRkZGRw4cGzrs/T09MWqOrjeO6tttCjgJ3dGmQIuAV51eT8BeL6adf8IPFdpXhvnz07AFqBzTZ9nIwn6h2VzXyDn82W2QyVlOuhvs3XCa/Pd2r5ytlnLc7T9H2fp9MXbVVV1zZo1nkb0yIEDB/z6+TVZsmTJr+bhw5EEl4rIOyIyvp5NdbMB14bbbZ3zqnI5lW5ZqWq28+cmYC4woA6faYwJcB8szmJ34WFuPKWTV/Z3dp+W9GoVy9NfZVBaHjh1H01dXQqPSOAwMBo4z/kaU4ftFgJdRaSjiIThKCB+1WpKRHoACcBPLvMSRCTcOZ0MDAfW1OEzjTEBrKy8ghe/3Uj/1HiGdU7yyj6DgoQ/nNmdbXuLeH/RdgAbfbAK3j4mtT7noarXuLNjVS0TkVuAL4BgYIqqrhaRh3BcJh0pSC4Hpuqx36wn8JKIVOAo4B5VVSs8jGnkHH1SHeKB83p7tZfXUd2bM6h9As9+ncE74zqxZ88e65bdhTrH8ygv996IjNUWHiJyf81Z9G+17VxVPwU+rTTv/krv/1rFdj8CfWvbvzGm8aiocHSA2KNlDKf1aOHVfYsIvx/dnfGv/My8nApGBhWwa9cur35GXRUXFxMRUfuIhw0tIiKCgwcPem1/NV15VPUpzYBrgSSg1sLDGGOO+HJNLpl5hTxzeX+Cgrx/RTCscxInd03m2bmbueTudKKreerc1+bOncuAAYFZRbt1a42NZOul2joPVf3HkRfwMo66j2twPK/hnZouY8xxQVX519xM2idFcW7fVj77nLtGd2fvwRKmfB+YPe42JTVWmItIoog8DKzAcZUyUFX/qKp5DZLOGNMkfJ+5mxVZ+dw4snOdul13V//UeM7olcIr321if1GJzz7H1DwY1BM4WkwVAH1V9a+quq/BkhljmowX5mTSMjaCsQPb+Pyz7hrdjcKSMl600QZ9qqb/AtwFtAb+DOSIyAHnq0BEvNtBvjGmyVq8dS8/b9rL9ad0IjzE909e92gZy4X92zDl+821jnVu3FdTnUeQqkaqaoyqxrq8YlQ1tiFDGmMarxfmbCQhKpTxQ9wbTc8dfz63J7GRIdw+dRmHy7zXPNX8otabj1UNNSsij/omjjGmKflx426+WZfHpOEdiQpruNZPSdHhPHZxGut2FvDU7IzaNzD1Vpeaq4tF5Mojb0TkBcC7jbSNMU1OUUkZ90xfSfukKK47ueEbaJ7WM4XxQ1J56buNLNi8t8E/v6mrU+EBXO3s2+pNoEzrOBiUMeb49Y8vN7BtbxGPXZxWryFmvenP5/YiNSGKO99fRkGx+6PzmV+rqbVVoogk4ni+4zrgbhwtrx50zjfGmCot2baPKT9s5soT2zG0k3f6sHJHs/AQnhrXj5z9h3hopvVw5E01XXksBhY5f84B4oFzXeYbY8yvHC4r5+4PVtAqNoJ7zu7h7zgMap/ITaM6M21xFl+s3unvOE1GTa2tOqpqp0o/j7zsCXNjTJWe+zqTzLxC/m9sX2IiQv0dB4DfndaNPm1iufO9ZSzeavUf3lDTbasRNW0oIrEi0sf7kYwxjdXqnHz+/e1Gxg5sw6jugdOuJiwkiNcmnkCL2AgmTlnI4q32vLOnarptdbGI/Cgi94vIuSIyREROEZFJIvIfYBaO+hBjjKHwcBm/n7aChKgw7h/Ty99xfiUlNoJ3rx9KcnQYE6csYMk2K0A8UdNtqztwDPq0A7gURy+6dwJdgZdU9RRVXdggKY0xAa2opIxJbyxkQ24BT1ySRnxUmL8jVallXATvTh5KUnQYE19bwLLt+/0dqdGqsamuqu5V1VdU9WpVPVNVL1TVe1X1+4YKaIwJbMWl5Vz/1iIWbdnLU+P6k+7lsTq8rVVcJO9eP5SEZmFMeG2+XYG4yXfdWxpjmrzi0nIm/2cxP27cw5OX9uP8fq39HalOWsdH8u7koSREhTHupZ/499yNlFfY0LX1YYWHMcYtJWUV/PbtJXy3YRePjU1j7MC2/o5UL23iI/nkt8M5o1cKj32+jstf/onte4v8HavR8GnhISJnich6EckUkXuqWH61iOwSkWXO13UuyyaKSIbzNdGXOY0x9bO/qIQb/7uYr9fl8chFfbjshIbr9NCbEpqF8cIVA3lqXD/W7Sjg7GfmMW3RdlTtKqQ2dekYMUpE/iIirzjfdxWRMXXYLhh4ATgb6AWMF5GqmmC8p6r9na9XndsmAg8AJwJDgAdEJKHO38oY4zPfbtjFmU9/x3cbdvG3C/tw5Ynt/R3JIyLCRQPa8vkdp9CnTSx/+GAFv5mywLpzr0VdrjxeBw4Dw5zvs4GH67DdECBTVTepagmO4WsvqGOuM4HZzgr7fcBs4Kw6bmuM8YFDJeXc/8kqJk5ZQGxEKB//djgThjbugsNVm/hI3rluKA+c14sVWfmc88w87v5gOTvzi/0dLSBJbZdnIrJIVQeLyFJVHeCct1xV+9Wy3SXAWap6nfP9BOBEVb3FZZ2rgb8Du4ANwB2qul1Efg9EqOrDzvX+AhxS1ScrfcZkYDJASkrKoKlTp9bjqx+rsLCQ6Ohot7f3JcvmnkDOBoGdr3K2TfnlvLz8MDuLlNHtQ7ikWxhhwRIQ2XzyGSXKzE0lfL21jCCBMzuGcm7HUCJCav/Ojen3CpCenr5YVQfXd1916WC/REQiAQUQkc44rkS8YSbwrqoeFpEbgDeBU+u6saq+DLwMMHjwYB01apTbQebOnYsn2/uSZXNPIGeDwM53JFtxaTnPfJ3By/M30SImnLev68fwLskBkc3XxgDb9xbx+Bfrmbk8h0W7g/nr+b0Z3SsFkeoLkcbwe/WGuty2egD4HEgVkbeBr3H0sFubbMC1Fq2tc95RqrpHVY8URK8Cg+q6rTHGt5Zu28eY577n33M3cvHANnx++yl+LzgaWmpiFM+NH8D0m04iLjKUG/6zmOvfWkTWPmuVVWvhoaqzgbHA1cC7wGBVnVuHfS8EuopIRxEJAy4HZriuICKtXN6eD6x1Tn8BjBaRBGdF+WjnPGOMjxWXlvP++hIu/vePHDxcxhvXnMDjl/QjLjIwOjn0h0HtE5h16wj+dE5Pfty4hzP++R0vfruR0vIKf0fzm2pvW4nIwEqzdjh/thORdqq6pKYdq2qZiNyC449+MDBFVVeLyEPAIlWdAdwmIucDZcBeHAUUqrpXRP6GowACeEhVrStMY3wsa18R17+1mLU7Shk/JJV7z+lJbID0jOtvIcFBXH9KJ85Ja8WDM1bz6Gfr+GZtHs9fOYAWMRH+jtfgaqrz+IfzZwQwGFgOCJCGYzyPYdVsd5Sqfgp8Wmne/S7T9wL3VrPtFGBKbZ9hjPGOxVv3csN/FnO4rILbB4Zz+9g0f0cKSG3iI3n5N4P5ZFk290xfyZhnv+dfVw5kcIfja4y8mjpGTFfVdBxXHANVdbCqDgIGYPUPxjQp0xdnMf7l+USHh/DRzcPp36IubWmObxf0b8NHvz2JqLBgLn/5Z974YfNx9XBhXSrMu6vqyiNvVHUV0NN3kYwxDaW8Qnn0s3XcNW05gzsk8PFvh9OlRWA2Mw1EPVrG8sktIxjVvTl/nbmGO99fTkn58VGA1OW/FytE5FXgv873VwIrfBfJGNMQDhSXcsfUZXy9Lo8rT2zHX8/vTWiwdXdXX3GRobw8YTD/mpvJP2ZvIDslmDNO1Rqb8zYFdSk8rgFuAn7nfP8d8G+fJTLG+NzGXYVc/9Yitu0p4qELevObYR38HalRCwoSbjm1K8FBQTz2+Tqe/yaTW0/r6u9YPlVr4aGqxcBTzpcxppH7em0ut09dRlhIEP+97kSGdkryd6Qm48aRnfhuRQb/mL2BrinRnNWnVe0bNVK1Fh4ishnn0+WuVLWTTxIZY3xCVXlhjuPWSu/Wsbw0YTBt4m0kaW8SEa7pHU5RcBh3vLec1MQoereO83csn6jLDc7BwAnO18nAs/xS/2GMaQS27jnIb6Ys4MkvN3BBv9ZMu+EkKzh8JCxYeGXCIOIiQ7n+zUXsKvBWb06BpS5PmO9xeWWr6tPAub6PZozxVGl5Bf+am8nop75j6bb9/O2C3jw1rj+RYcH+jtaktYiN4JXfDGavc9yT4tJyf0fyurrctnJ90jwIx5WINQI3JsAt2baP+z5cybqdBZzZO4UHz+9Dy7jj70lof+nbNo5/XNqf376zhNveXcq/rhxISBNqzVaXQuAfLtNlwGbgMt/EMcZ4QlVZtHUfr/+wmc9W7aRlbAQvTxjE6N4t/R3tuHRuWivyCnrx4Mw13D19BU9e0o+goKbRhLcuhce1qrrJdYaIdPRRHmOMG0rKKvjfyhymfL+Fldn5xEWGctPIztyc3oXocLtR4E/XDO9IQXEZ/5y9gdiIUB44r1eTeAakLv+qPgAqd5L4Ab90n26M8YODh8v4aeMevt2wi89X72RXwWE6N2/Gwxf2YezANkSFWaERKG49tQv5h0p57fvNxEaGcucZ3fwdyWM19arbA+gNxInIWJdFsTg6SzTGuDhcVk5mXiERocF0TGrm9dsTJWUVrMrJZ8HmvXy3YRcLt+yltFyJCgtmeJdkrhranpO7JDeZ2yJNiYjw53N7UlBcyrNfZxAbEcJ1Jzfupx1q+q9JdxyDacUD57nMLwCu92EmYxqFjNwCvlmXx9odB1i7o4CNuwopq3A8EhUbEUK/1HgGpMYzoF0CfdrE0TwmvM77VlVy8otZm3OAxdv2sXjLPpZn7edwmWP8iB4tY5g0vCMjuzVnUIcEwkOs9VSgExH+PjaNwsNlPPy/tZSWKzeO7NRob2FVW3io6ifAJyIyTFV/asBMxgS0guJSnpqdwZs/baG8QmkVF0HPVrGc3qsFPVrGcqi0nKXb9rNs+36en5OJszwhOTqMnq1i6dEyhh4tY9m0s4z8ZdmUlFVQUl5BcWkFm3YVsn5nAetzCygoLgMgJEjo3SaOq4a2Z3D7BAZ1SDgux49oCoKDhKfHDSAkaDmPfb6OvIJi/nJur0Z5tVjTbau7VfVx4AoRGV95uare5tNkxgQYVWXG8hwe+d9adhUe5ooh7bj99G5VXlFcNtgxivLBw2WszM5nTc4BxxXKzgO8+dNWSpxXECxbdsx2MREh9GgZwwX9W9OjpaOg6dMmjohQu7JoKsJCgnh6XH+So8OZ8sNmdheW8OSlaY3u6rGm21ZHhoRd1BBBjAlkmXmF3P/JKn7cuIe0tnG88pvB9EuNr3W7ZuEhDO2UdEz/UWXlFWzZU8RP8xcw7MQhhIcEERYSRFhwEPFRoY32Noapu6Ag4S9jetIiNpxHP1vHvoMlvDhhUKNqGVfTbauZzp9vNlwcYwLLoZJynvsmg1fmbSIyNJiHL+zD+CHtCPbgNkNIcBBdWkSTFRNkY2ccx0SEG0d2Jjk6nD9OX8G4l37ixasGkZoY5e9odVLTbauZVNEh4hGqer5PEhkTIL5em8sDM1aTte8QFw9sy73n9CA5uu6V3sbUxSWD2pIUHcZt7y7lnGfn8fjFaZzdN/B7463pGulJT3cuImcBzwDBwKuq+mil5XcC1+F4cn0XMElVtzqXlQNHRjDcZoVV47UqO5+8gmIGtU8kLjLU33FqtXn3Qf7v07XMXpNL1xbRvDd5KCdat+XGh9K7t+DT207mlneXctPbS/jNsPbcd07PgK7rqum21bdHpkUkDOiB40pkvaqW1LZjEQkGXgDOALKAhSIyQ1XXuKy2FBisqkUichPwODDOueyQqvav5/cxAWRVdj5Pzd7A1+vyABCBXq1iObFjEid2SmRY5yRiIwKjMFFVftq0hynfb+brdXlEhARzz9k9mDS8I2EhTac/IhO4UhOjmHbDMJ74Yh2vzNvMoi37eP6KAXRqHpi3NuvSMeK5wIvARkCAjiJyg6p+VsumQ4DMI12biMhU4ALgaOGhqnNc1v8ZuKp+8U0g2pBbwFOzN/DZqp3ERoTwhzO7MyA1ngVb9jJ/017enr+VKT9sJjwkiDN7t+TiQW0Z0SXZo3oEdx0qKed/K3cw5fvNrNlxgMRmYdya3oWrhrW35rCmwYWFBPGnc3txYsckfv/Bcs599nvuOKMr1wzvGHBDBItqzYO1i8g6YIyqZjrfdwb+p6o9atnuEuAsVb3O+X4CcKKq3lLN+s8DO1X1Yef7MmAZjltaj6rqx1VsMxmYDJCSkjJo6tSpNX6XmhQWFhIdHZglfGPJVlKuTF1fwpxtZYQHw5kdQhndIZRmoccWCqUVyqb9FSzYWcbPO8o4WAoJ4cJJrUNIax5MakwQUaGeFyTVHbd9xRUsyytn2a5y1uwpp7QC2kQLo9uHMqx1CGHBDVOINZbfa6AJ5GzgvXx7iyt4a3UJy3aVkxoTxNW9w+gc79ltrKqypaenL1bVwfXdV13ahRUcKTicNuF4ytxrROQqHF29j3SZ3V5Vs0WkE/CNiKxU1Y2u26nqy8DLAIMHD9ZRo0a5nWHu3Ll4sr0vNYZsmXkF3PLOUtbtLOOa4R247dSuJDQLq3HbG3B06fH12jw+WJzF5xt28b/NpQCkJkbSs2UsPVvF0ql5M1ITo0hNiCI5OqzOTVnnzp3LyaeMZOOuQlZm5bMyO5+FW/ayOucgAO0So7hqWAtG92rJ0E6JDd5EtjH8XgNRIGcD7+a76Ezli9W5/HXGah6eX8yEoe35/Znd3b7d681sdSk8FonIp8D7OOo8LsVRfzEWQFU/rGa7bCDV5X1b57xjiMjpwJ+Akap6dMgtVc12/twkInOBAThunZkAoqq8v2g7D3yymqiwYN645gRGdW9R5+3DQ4I5p28rzunbij2Fh1mRlc+aHQdYs+MAa3MOMHttLq4Xx5GhwaQmRpLULJyEZqHER4WREBVKXGQoRSXlHDhUxoHiUg4cKmVTziGyv/6CQ86BeKLCgunTJo67z+rOGT1T6NIi2p6pMAFNRDirT0uGd0niH19u4M2ftvDpyp3ceUY3Lhvc1q/jg9Sl8IgAcvnlqmAXEImjvysFqis8FgJdnd23ZwOXA1e4riAiA4CXcNzeynOZnwAUqephEUkGhuOoTDcBpPBwGS+vPMxPOSsY1imJpy/vT0qs+/UESdHhpPdoQXqPXwqfQyXlbN9XxLY9RWzfV8T2vYfI2lfE3oMlrN9ZwP6iUvYfKqXc2QdIdHgIsREhxEaGEhIE405IJa1tHGlt4+iYHO2XehVjPBUTEcpfz+/N2IFteGjmGu77aCWv/7CZ+87pyajuzf3yn6BaCw9VvcadHatqmYjcAnyBo6nuFFVdLSIPAYtUdQbwBBANTHN++SNNcnsCL4lIBY7RCx+t1ErL+FFFhTJ9SRaPf7Ge3QXl3HF6N245tYtP/jBHhgXTLSWGbikxNeY5WFJGZGjwMf8Tc1yi9/Z6JmP8Ja1tPNNuHMYXq3fy6GfruOaNhQzvksR95/Skd+u4Bs1Sl9ZWHYFbgQ6u69fluQtV/RT4tNK8+12mT69mux+BvrXt3zS8xVv38uDMNazIyqd/ajw39hauPb2rXzMFBQkxAdLk1xhfc9zKasWpPVJ4e/5Wnvk6gzvfW87nt5/coFcgdblt9THwGjATqPBpGhOwtu8t4okv1jNjeQ4pseE8Na4fF/Rrw3fffVv7xsYYrwsLCeKa4R0ZO7AtuQeKG/zWVV0Kj2JVfdbnSUxA2raniBfmZDJ9SRbBQcKtp3bhxpGdadaIOnAzpimLiwz1S88NdfkL8IyIPAB8Cbi2hlris1TG77btKeL5ORlMX5JNcJBw1dD23DiyMy3j7ME5Y0zdCo++wATgVH65baXO96YJKa9QvsvYxXsLtjN7bS7BQcKEoe25aVRnj1pRGWOanroUHpcCnerSn5VpnLL3H+L9hduZtmg7OfnFJDUL47oRHZk0oqMVGsaYKtWl8FiFYxzzvFrWM41E4eEyFmzeww+Ze/ghczfrdhYgAiO6JPPnMb04vWeKdQZojKlRXQqPeGCdiCzk2DoP6yI9AKkqxaUVFJWUUXi4jJz9xWzfV0TW3iKy9h1i0+6DrMrOp6xCCQsJYnD7BP5wZnfO79e60QxCY4zxv7oUHg/4PIWpl+LScjbkFrAht5CM3ALW5xaQmVfI/qJSDpaUUVVfl0ECreIiaZsQyQ0jOzG8czID2ycE9HgBxpjAVZcnzI9pyC8iI4DxgDXwbyCHypS56/OYv3kv8zftYUWW48oBICw4iE7NmzGwXQLJ0eE0Cw8mKizk6M+WsRGkJkbSKi7SbkUZY7ymTo31nX1QXYGj8nwzMN2XoRozVaWkvIKw4KB6P7RTUlZB9v5DbHBeSRy5uli/s4gKXUhIkNC3bRzXndyJtLZxdEuJoUNSlF87RzPGHJ9qGsO8G44rjPHAbuA9HON/pDdQtoC3bU8R327IY+ueIrbtLWL7vkNs31tE4eEyQoKEmIgQoiNCiAkPJTIsmGARgoIgSITgIKGkrIL8Q6VHX0Ul5cfsv018JF1ToukSGcpl6QMY1D6BqDB7OM8Y4381/SVaB8zj2IGg7miQVAFs38ES/rdyBx8vzWbR1n0ARIQGkZoQRbvEKE7smEhydBhFJeUUHi6joNjxOlRaRkUFlKtSXlFBcakSEhREamIUfSJDiXc+JZoSF0G3lBi6tIgm2vkU99y5czm5a3N/fm1jjDlGTYXHWBzdqM8Rkc+BqTiGoT3ubN9bxLyM3XyzLo9vN+RRWq50bRHN3Wd1Z0zf1qQmRtq4EMaY40q1hYdz2NePRaQZjrHHbwdaiMi/gY9U9csGSdjAKiqUbXuLWJWTz48b9/B9xm627S0CoGVsBBOHdeDCAW3o3TrWCgxjzHGrLq2tDgLvAO84B2m6FPgjjr6uGr0DxaVMX5zFnFWHeWr1D2zYWXB05Lno8BCGdkpk0vAOjOiaTOfmNvKcMcZAHVtbHaGq+3CMGf6yb+I0PK2AB2euIToU+qYGc/mQVHq0jKFHy1h6tY4l1FoyGWPMrxz3TXfiokJZ+KfTWbXoR9LTh/o7jjHGNAr232qgeUy43Y4yxph6sMLDGGNMvVnhYYwxpt5Eq+pFrxESkV3AVg92kYzjSfpAZNncE8jZILDzWTb3BXK+qrK1V9V6P4XcZAoPT4nIIlUd7O8cVbFs7gnkbBDY+Syb+wI5nzez2W0rY4wx9WaFhzHGmHqzwuMXgfzgo2VzTyBng8DOZ9ncF8j5vJbN6jyMMcbUm115GGOMqTcrPIwxxtRbky08RGSKiOSJyCqXeYkiMltEMpw/E5zzRUSeFZFMEVkhIgNdtpnoXD9DRCb6MNsTIrLO+fkfiUi8y7J7ndnWi8iZLvPPcs7LFJF7vJGtunwuy+4SERWRZOd7vx875/xbncdvtYg87jK/wY5dNb/X/iLys4gsE5FFIjLEOb+hj1uqiMwRkTXOY/Q753y/nxM1ZAuIc6K6fC7L/XZO1JTN5+eEqjbJF3AKMBBY5TLvceAe5/Q9wGPO6XOAz3AMdjUUmO+cnwhscv5McE4n+CjbaCDEOf2YS7ZewHIgHOgIbASCna+NQCcgzLlOL18dO+f8VOALHA9jJgfQsUsHvgLCne9b+OPYVZPtS+Bsl2M110/HrRUw0DkdA2xwHh+/nxM1ZAuIc6K6fIFwTtRw7Hx+TjTZKw9V/Q7YW2n2BcCbzuk3gQtd5r+lDj8D8SLSCjgTmK2qe9XRHf1s4CxfZFPVL1W1zPn2Z6CtS7apqnpYVTcDmcAQ5ytTVTepagmOkR4v8DRbdfmcngLuBlxbWfj92AE3AY+q6mHnOnku2Rrs2FWTTYFY53QckOOSrSGP2w5VXeKcLgDWAm0IgHOiumyBck7UcOzAz+dEDdl8fk402cKjGimqusM5vRNIcU63Aba7rJflnFfdfF+bhON/LgGTTUQuALJVdXmlRYGQrxtwsojMF5FvReSEAMp2O/CEiGwHngTu9Xc2EekADADmE2DnRKVsrgLinHDNF2jnRKVj5/Nz4rgdz0NVVUQCrp2yiPwJKAPe9neWI0QkCrgPx22EQBSC41bAUOAE4H0R6eTfSEfdBNyhqtNF5DLgNeB0f4URkWhgOnC7qh4Ql6EI/H1OVM7mMj8gzgnXfM48AXNOVPF79fk5cbxdeeQ6Lx9x/jxyKZeN497lEW2d86qb7xMicjUwBrhSnTcoAyRbZxz3R5eLyBbnZy0RkZYBki8L+NB5m2ABUIGjA7hAyDYR+NA5PQ3H7QH8kU1EQnH8gXlbVY9kCohzoppsAXNOVJEvYM6Jao6d78+J+lbQNKYX0IFjKy+f4NjKwced0+dybAXXAv2lgmszjsqtBOd0oo+ynQWsAZpXWq83x1ZwbcJRuRXinO7ILxVcvX117Cot28IvlYOBcOxuBB5yTnfDcfkt/jh2VWRbC4xyTp8GLPbHcXN+zlvA05Xm+/2cqCFbQJwT1eULhHOihmPn83PC45MlUF/Au8AOoBRHKXwtkAR8DWTgaImQ6PILeAFHa4OVwGCX/UzCUamUCVzjw2yZzl/wMufrRZf1/+TMth5nyx3n/HNwtK7YCPzJl8euhhMlEI5dGPBfYBWwBDjVH8eummwjgMXOk3E+MMhPx20EjkrdFS7/xs4JhHOihmwBcU5Uly8Qzokajp3PzwnrnsQYY0y9HW91HsYYY7zACg9jjDH1ZoWHMcaYerPCwxhjTL1Z4WGMMaberPAwxhhTb1Z4GGOMqTcrPIzxAhGJdHZAF+x8Xy6OMTxWicjMSmNRuC6b5uw7DOeYEP91WS9ERHaJyCzn+zAR+c7Zb5ExfmWFhzHeMQlHX0LlzveHVLW/qvbB0U37b13WdV1WgqMrCYCDQB8RiXS+PwOX/oXU0VX218A4H34PY+rECg9jqiEi54vI9ErzbhKR56pY/Urgk2p29RPVd289D+ji8v5THH0jAYzH0eWJq4+dn2WMX1nhYUz1HgEeqDRvI9DTdYaIhAGdVHVL5R04b2OdBsyoYlkIcDaO/o+OmApcLiIRQBq/HtdiFY4uto3xKys8jKmCiPQDglR1lYi0F5GbnItCOXbUOHB0db2/0rxIEVnGLwMsza5i2SJgG44xPgBQ1RU4euYdj+Mq5BjO22IlIhLj1hczxkus4s2YqvXH0RsuOOoeujqnj4wB7eoQEFF5nqr2d1aGf4GjzuNZ12U1fPYMHKMOjsLR621l4UBxrd/AGB+yKw9jqhYERDtvO40FYpwV2VcD77iuqI7xqIOdt5qotKwIuA24qx6tpKYAD6rqysoLRCQJ2K2qpfX5MsZ4mxUexlTtU6ATznEkcAyiswh4WVWXVLH+lzjGVvgVVV2KY7yF8XX5YFXNUtVnq1mcDvyvLvsxxpdsPA9jvEBEBuIYq3yCjz/nQxwj/23w5ecYUxu78jDGC5xXI3OOPCToC85WXR9bwWECgV15GGOMqTe78jDGGFNvVngYY4ypNys8jDHG1JsVHsYYY+rNCg9jjDH1ZoWHMcaYevt/jxdLe4YwGUkAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "assembly = crankshaft_assembly()\n", "rpms = np.arange(1000, 2575, 25)\n", "# Vibratory torque for every shaft and every considered engine speed\n", "vibratory_torque = []\n", "for rpm in rpms:\n", " vibratory_torque.append(calculate_response(assembly, rpm))\n", "\n", "shaft_8, shaft_1 = plot_results(rpms, vibratory_torque)\n", "# Same plots but using openTorsions Plots class\n", "plots = opentorsion.Plots(assembly)\n", "plots.torque_response_plot(rpms, [np.array(shaft_8), np.array(shaft_1)], True)" ] }, { "cell_type": "markdown", "id": "6e0f0d44", "metadata": {}, "source": [ "[1] Mendes AS, Meirelles PS, Zampieri DE. Analysis of torsional vibration in internal combustion engines: Modelling and experimental validation. Proceedings of the Institution of Mechanical Engineers, Part K: Journal of Multi-body Dynamics. 2008;222(2):155-178. doi:10.1243/14644193JMBD126" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.8.10" } }, "nbformat": 4, "nbformat_minor": 5 }