{ "cells": [ { "cell_type": "markdown", "id": "71990519", "metadata": { "papermill": { "duration": 0.002836, "end_time": "2024-11-19T22:31:12.550591", "exception": false, "start_time": "2024-11-19T22:31:12.547755", "status": "completed" }, "tags": [] }, "source": [ "# GE Discovery MI (Listmode Reconstruction; With Time of Flight)" ] }, { "cell_type": "markdown", "id": "7c513fdf", "metadata": { "papermill": { "duration": 0.00225, "end_time": "2024-11-19T22:31:12.555794", "exception": false, "start_time": "2024-11-19T22:31:12.553544", "status": "completed" }, "tags": [] }, "source": [ "This tutorial demonstrates how to reconstruct data listmode data in HDF5 data obtain from GE scanners. In this tutorial, we'll reconstruct using time-of-flight listmode. The data can be obtained from https://zenodo.org/records/8404015" ] }, { "cell_type": "code", "execution_count": 1, "id": "ea48e0e5", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:31:12.563036Z", "iopub.status.busy": "2024-11-19T22:31:12.562415Z", "iopub.status.idle": "2024-11-19T22:31:16.653254Z", "shell.execute_reply": "2024-11-19T22:31:16.652547Z" }, "papermill": { "duration": 4.096959, "end_time": "2024-11-19T22:31:16.655415", "exception": false, "start_time": "2024-11-19T22:31:12.558456", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "from pytomography.metadata import ObjectMeta\n", "from pytomography.metadata.PET import PETLMProjMeta\n", "from pytomography.projectors.PET import PETLMSystemMatrix\n", "from pytomography.algorithms import OSEM, BSREM\n", "from pytomography.priors import RelativeDifferencePrior\n", "from pytomography.likelihoods import PoissonLogLikelihood\n", "from pytomography.transforms.shared import GaussianFilter\n", "import matplotlib.pyplot as plt\n", "from pytomography.io.PET import clinical" ] }, { "cell_type": "markdown", "id": "fc81208e", "metadata": { "papermill": { "duration": 0.002419, "end_time": "2024-11-19T22:31:16.660738", "exception": false, "start_time": "2024-11-19T22:31:16.658319", "status": "completed" }, "tags": [] }, "source": [ "We first need to get the required PET geometry information dictionary from the scanner name. Currently, we only support `discovery_MI`, but feel free to make a commit to support other scanners (add to the `src/data/pet_scanner_info.txt` file)" ] }, { "cell_type": "code", "execution_count": 2, "id": "b9b22ab3", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:31:16.667007Z", "iopub.status.busy": "2024-11-19T22:31:16.666609Z", "iopub.status.idle": "2024-11-19T22:31:16.684871Z", "shell.execute_reply": "2024-11-19T22:31:16.684247Z" }, "papermill": { "duration": 0.023246, "end_time": "2024-11-19T22:31:16.686263", "exception": false, "start_time": "2024-11-19T22:31:16.663017", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "scanner_name = 'discovery_MI'\n", "info = clinical.get_detector_info(scanner_name)" ] }, { "cell_type": "markdown", "id": "fc6ba139", "metadata": { "papermill": { "duration": 0.003787, "end_time": "2024-11-19T22:31:16.692708", "exception": false, "start_time": "2024-11-19T22:31:16.688921", "status": "completed" }, "tags": [] }, "source": [ "Now we obtain the required data from the listmode files." ] }, { "cell_type": "code", "execution_count": 3, "id": "2d2bd0c3", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:31:16.699534Z", "iopub.status.busy": "2024-11-19T22:31:16.699042Z", "iopub.status.idle": "2024-11-19T22:31:55.254439Z", "shell.execute_reply": "2024-11-19T22:31:55.253351Z" }, "papermill": { "duration": 38.561221, "end_time": "2024-11-19T22:31:55.256256", "exception": false, "start_time": "2024-11-19T22:31:16.695035", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "# Get listmode events\n", "detector_ids = clinical.get_detector_ids_hdf5('/disk1/ge_discovery_tof_mi_pet/dmi_nema_lm/LIST0000.BLF', scanner_name)\n", "# Get scatter/random correction term at each event\n", "additive_term = clinical.get_additive_term_hdf5('/disk1/ge_discovery_tof_mi_pet/dmi_nema_lm/corrections.h5')\n", "# Get multiplicative weights (attenuation/normalization) at each event\n", "weights = clinical.get_weights_hdf5('/disk1/ge_discovery_tof_mi_pet/dmi_nema_lm/corrections.h5')\n", "# Get ALL valid detector-pairs and corresponding multiplicative weights\n", "detector_ids_sensitivity, weights_sensitivity = clinical.get_sensitivity_ids_and_weights_hdf5('/disk1/ge_discovery_tof_mi_pet/dmi_nema_lm/corrections.h5', scanner_name)" ] }, { "cell_type": "markdown", "id": "f025ba63", "metadata": { "papermill": { "duration": 0.005228, "end_time": "2024-11-19T22:31:55.264302", "exception": false, "start_time": "2024-11-19T22:31:55.259074", "status": "completed" }, "tags": [] }, "source": [ "First we'll get the required object space and projection space metadata.\n", "\n", "* The `object_meta` specifies the size of each voxel, and the number of voxels in each direction\n", "* The `proj_meta` specifies the list of all measured listmode events (`detector_ids`), the time-of-flight metadata, and the detector IDs / weights used for computing the sensitivity image." ] }, { "cell_type": "code", "execution_count": 4, "id": "0bc888ff", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:31:55.271104Z", "iopub.status.busy": "2024-11-19T22:31:55.270513Z", "iopub.status.idle": "2024-11-19T22:31:55.409643Z", "shell.execute_reply": "2024-11-19T22:31:55.409044Z" }, "papermill": { "duration": 0.144441, "end_time": "2024-11-19T22:31:55.411319", "exception": false, "start_time": "2024-11-19T22:31:55.266878", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "# Define object space reconstruction matrix\n", "object_meta = ObjectMeta(\n", " dr=(2.78,2.78,2.78), #mm\n", " shape=(192,192,71) #voxels\n", ")\n", "# Get time-of-flight metadata and define projection space metadata\n", "tof_meta = clinical.get_tof_meta(scanner_name)\n", "proj_meta = PETLMProjMeta(\n", " detector_ids,\n", " info,\n", " tof_meta=tof_meta,\n", " detector_ids_sensitivity=detector_ids_sensitivity,\n", " weights_sensitivity=weights_sensitivity\n", " )" ] }, { "cell_type": "markdown", "id": "344fe1e5", "metadata": { "papermill": { "duration": 0.002602, "end_time": "2024-11-19T22:31:55.416741", "exception": false, "start_time": "2024-11-19T22:31:55.414139", "status": "completed" }, "tags": [] }, "source": [ "Now we'll define the system matrix for the PET system. In this case we'll use 4.5mm PSF modeling." ] }, { "cell_type": "code", "execution_count": 5, "id": "9daf97ba", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:31:55.422945Z", "iopub.status.busy": "2024-11-19T22:31:55.422422Z", "iopub.status.idle": "2024-11-19T22:32:00.599612Z", "shell.execute_reply": "2024-11-19T22:32:00.598593Z" }, "papermill": { "duration": 5.18235, "end_time": "2024-11-19T22:32:00.601404", "exception": false, "start_time": "2024-11-19T22:31:55.419054", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "psf_transform = GaussianFilter(4.5) # 4.5mm gaussian psf\n", "system_matrix = PETLMSystemMatrix(\n", " object_meta,\n", " proj_meta,\n", " obj2obj_transforms = [psf_transform],\n", " N_splits=8,\n", ")" ] }, { "cell_type": "markdown", "id": "561328db", "metadata": { "papermill": { "duration": 0.002514, "end_time": "2024-11-19T22:32:00.606717", "exception": false, "start_time": "2024-11-19T22:32:00.604203", "status": "completed" }, "tags": [] }, "source": [ "Now we'll define the likelihood function (`PoissonLogLikelihood` for PET systems). The additive term contains scatters/randoms and needs to be normalized by the `weights`." ] }, { "cell_type": "code", "execution_count": 6, "id": "e7ff8d16", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:32:00.613616Z", "iopub.status.busy": "2024-11-19T22:32:00.613230Z", "iopub.status.idle": "2024-11-19T22:32:00.691825Z", "shell.execute_reply": "2024-11-19T22:32:00.691129Z" }, "papermill": { "duration": 0.083695, "end_time": "2024-11-19T22:32:00.693545", "exception": false, "start_time": "2024-11-19T22:32:00.609850", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "likelihood = PoissonLogLikelihood(\n", " system_matrix,\n", " additive_term=additive_term/weights,\n", ")" ] }, { "cell_type": "markdown", "id": "bb8c61cc", "metadata": { "papermill": { "duration": 0.002528, "end_time": "2024-11-19T22:32:00.698885", "exception": false, "start_time": "2024-11-19T22:32:00.696357", "status": "completed" }, "tags": [] }, "source": [ "We can now reconstruct using a variety of reconstruction algorithms:" ] }, { "cell_type": "code", "execution_count": 7, "id": "b83a54d6", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:32:00.705301Z", "iopub.status.busy": "2024-11-19T22:32:00.705056Z", "iopub.status.idle": "2024-11-19T22:32:50.069910Z", "shell.execute_reply": "2024-11-19T22:32:50.068854Z" }, "papermill": { "duration": 49.370579, "end_time": "2024-11-19T22:32:50.071796", "exception": false, "start_time": "2024-11-19T22:32:00.701217", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "recon_algorithm = OSEM(likelihood)\n", "recon_OSEM = recon_algorithm(n_iters=2, n_subsets=34)" ] }, { "cell_type": "code", "execution_count": 8, "id": "85c06fa9", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:32:50.080010Z", "iopub.status.busy": "2024-11-19T22:32:50.079536Z", "iopub.status.idle": "2024-11-19T22:42:32.549816Z", "shell.execute_reply": "2024-11-19T22:42:32.548836Z" }, "papermill": { "duration": 582.476715, "end_time": "2024-11-19T22:42:32.551609", "exception": false, "start_time": "2024-11-19T22:32:50.074894", "status": "completed" }, "tags": [] }, "outputs": [], "source": [ "# This is probably very similar to \"Q.Clear\" is ;)\n", "prior = RelativeDifferencePrior(beta=25, gamma=2)\n", "recon_algorithm = BSREM(likelihood, prior=prior)\n", "recon_BSREM = recon_algorithm(n_iters=20, n_subsets=34)" ] }, { "cell_type": "markdown", "id": "94c2c796", "metadata": { "papermill": { "duration": 0.002541, "end_time": "2024-11-19T22:42:32.557084", "exception": false, "start_time": "2024-11-19T22:42:32.554543", "status": "completed" }, "tags": [] }, "source": [ "Plot axial slice of reconstruction:" ] }, { "cell_type": "code", "execution_count": 9, "id": "da80f06d", "metadata": { "execution": { "iopub.execute_input": "2024-11-19T22:42:32.563539Z", "iopub.status.busy": "2024-11-19T22:42:32.562994Z", "iopub.status.idle": "2024-11-19T22:42:32.694752Z", "shell.execute_reply": "2024-11-19T22:42:32.694141Z" }, "papermill": { "duration": 0.136674, "end_time": "2024-11-19T22:42:32.696090", "exception": false, "start_time": "2024-11-19T22:42:32.559416", "status": "completed" }, "tags": [] }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAhwAAAD4CAYAAACntD/BAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjkuMiwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8hTgPZAAAACXBIWXMAAA9hAAAPYQGoP6dpAAEAAElEQVR4nOz9ebxtSVkfjH+r1lp7OOOd7+2mZ2hAoQURIhC1mfwZNYo4oCgRNFFRf3nzKppEowzJa4g/jdHwagRRcQBRwZhINApRcACVyQihoWnoufv2Hc+4p7VW1e+PqqfqqVq19t7nDn2H3k9/bp9z9l6rVq3pqe/zfSahtdZYyEIWspCFLGQhC7mIIi/1BBaykIUsZCELWcjVLwvAsZCFLGQhC1nIQi66LADHQhaykIUsZCELueiyABwLWchCFrKQhSzkossCcCxkIQtZyEIWspCLLgvAsZCFLGQhC1nIQi66LADHQhaykIUsZCELueiyABwLWchCFrKQhSzkossCcCxkIQtZyEIWspCLLvm8Gwox96YLWchCFrKQhSzkMSJaV3Ntt2A4FrKQhSxkIQtZyEWXBeBYyEIWspCFLGQhF10WgGMhC1nIQhaykIVcdFkAjoUsZCELWchCFnLRZQE4FrKQhSxkIQtZyEWXBeBYyEIWspCFLGQhF10uGuDo9Xp41au+B+9+93/HAw/ci+FwBxsbp3HHHZ/AL/3Sm/BlX/ZlM8f46q/+avz+7/8eHnjgXoxGu9jcPIO7774LH/jAX+IXf/EX8KpXfQ96vV6wz2tf+xpoXc31b3193e13++23N74/ffoEut1ucm5vfOPPNbb/1V/95fO7aAtZyFUkN954Y+u7Nxhs47777sYf/dH/wHd91z9DURSt49x22214y1vejDvvvAO7u1sYDLZx//334KMf/RB+8zd/HT/0Q6/GDTfcEOyTep/p387OJj7zmU/h7W//TTz/+c9PHnPa/vG/F7/4a4N9U9s861nPSh7nG77h6xvb3n33XXu80kZ+4Af+bzfGr/zKW4LvfvVXfzk5r6oa48yZk/jQh/4ab3jDv8exY8eSY6f2H48HOHPmJO688w784R++Gz/6oz+Ca665pnV+r3jFt7dew9FoF/fddzfe9a7fxVd91Vc19r3llltQliNoXeHBB+9Dv98/p2u0kEssek4Bsrn/Pfe5X6rvv//+mWP+1//6+3p1dV9yjF/8xTfNNa8bb7wl2O+1r339vKek19cPuP1uv/0FyW2+7dv+SWNu/f6KPnv2bGPbX/3Vt+7pOi3+Lf5dzf9uvPGWud/FD3zgg7rbXWqM8e3f/ko9mUxm7v+KV3xHsF/b+5ySH/ux1zSOu5f9X/zilwT7puSXfuktyWv0J3/ynsa2d999956v9fr6AX3q1CmttdZ1XesnPenzg+9/9VffOte5nD59Wt9229Mb48+7/3g81j/xE2/QUhaNMV7xiu+Yawyt07r0t37rHe77f/NvfvySP9+Lf9Of+ZRc8Gpez372s/Gnf/qegBk4fvw4PvKRj2J1dRXPfvYXo9PpAAC+7utejPe+90/wpV96OyaTidv+pS/9JnzP93y3+3tnZwcf+tCHsbm5iYMHD+KpT30K9u/fP9d8PvShD+Hee+9LfseP2Sbf/d3fhbe97e3BZ9/yLd+Mffv2zXX8hSxkIUZOnjyJ97//zwEAR44cwZd8yT+ElIZkfc5zno3v/u7vwhvf+P+67a+//nq86U3/xbEfVVXhIx/5KI4fP47l5WU8+clPwnXXXTfXse+55x58+MMfQa/Xw9Oe9gW4/vrr3Xevfe2P421vezvuvvvumfun5MEHH5p5/G/5lm/GD/zAq7Gzs+M+u+WWW/DCF75grvnPkh/6oVfj4MGDAID/8T/+EJ/+9Kenbv/+9/85Tp48iYMHD+I5z3m2Y4oPHDiAn/u5/4QXvOBFM/c/deoU9u3bh6c//Wnu2J1OBz/6o/8aT3rSE/GN3/jSqWN88pOfxCc/eQeEEHjiE2/Fbbfd5r575StfgY997O/wn//zG91n//E//id8y7d8MwDgX/7LH8LP//wvYGNjY+oxFnKZybzIZB6U0+n09b333hvs9zM/87M6z7tumxtvvEV//OMfD7Z5wxt+MhjnD/7g3e67z372s3rfvoPB90Lk+jnP+RL9pje9WV9zzXXBdzHDEVs+bf+mWTSxtfDXf/03ye0WDMfi3+Kf/xczHH/2Z+8Lvn/Vq74v+P5d7/q94PtXv/qH3Xd1XeunP/2LGse49dYn69e85nX6q77qa4LP4/eZv5tF0dPvf/+fB99/13d9z9z7z/rXJt/zPd8bbPcf/sP/L7ndXhmOoujpRx55xO3/zd/8ssY2MUNx++0vcN/ddtvT9Xg8dt9VVaV7veW595ey0K985Xfq7e3tYJsf/MEfCsaIGY7Xvvb1wfc//uOvDb6/8847G+dx5513uu//+T//F9H3Oft36Z//x9K/eeWCxnC8/OXfFvhS//zP/wI/+IOvRlX5sqf33nsvvuEbXoqyLN1n3//934u1tTX39xOe8Hj3+9///ccbKFZrjQ9+8IP4nu/5Xjz88MMX8hScPPSQt1q++7u/y/3+tKc9DV/8xf+gsc1CFrKQvQmxHSRxPBbXA5ubm/i7v/u7xhif+cxn8G//7b/DH/7hH8593LIs8Xu/91+Dz8hCv9DC9dN3f/c/c78XRYHv+I5XAACGwyHOnDlzzsd4yUu+DkeOHAEA7O7u4g/+4N172v/jH/84PvnJT7q/syzbE4OrlMJb3/pr+NZvfXnw+Y/8yL/aU6zFG97wHwLW+dZbb8Xq6mqwze/8zjvd79/1Xf80GkFCIIMJTRRzH3chj55cUMDxtV/7NcHfP//zv5Dc7s4778R73/u/3N+rq6t43vNud39zMPLVX/1VeN3rXounPOUpF3KqM+W3f/t3sbu7CwD49m9/uXMDvepV3tXzK7/y1kd1TgtZyNUkt98eBo7/3d/97+Bvrgf279+Pt7/9N/ElX/Il7l08HxEiXJAuluHy6U/fib/8y78CADzjGc/AF33RFwEAvv7rX+JAwrve9XvY2to652PwoNUPfejDGAwGex6DX4+qqnDy5Mk9j/EHf/Bu/M3f/K37+9ChQ61BuSmpqqphXK6srAR//9mfvc/9ftttt+Hmm2/e8zwXcunkggKOZzzjC4O/P/jBv27dNv7uGc94hvv9Ax/4oPu9KAq89rU/jk984n9jY+M03vveP8GP/uiPzP2gff/3fy9+93d/u/Hvh3/4h6but7m5iXe847cBmBfn67/+JVheXsa3fuvLAABbW1v4rd96x1xzWMhCFgI85Smf796/973vT/HzP+/98/feey9+7uf+c7A91wMA8LKXfQv+4i/eh62ts/jQh/4aP/3TP4VnP/vZe55HURT4+q9/ift7PB7jT/7kPVP3ed7zbk/qkbe+9VdmHu9Nb3qz+52sch6j9qY3/dJeTyGQL/3SL3G/t8WZTJMv+IIvwOd//ue7v//n//xj1HV9TnN5z3veG/z9rGc9c+59r7nmGgfCAAM4T58+HWzz4Q9/GEop93cIWmto1ABqAHov017IoyQXNGj00KFDwd/Hjx9v3Tb+7vBhv+9P/uRP4aUv/aZGYOj6+jpe+MIX4IUvfAFe//rX4ud//hfw6lf/8NSX41nPelYyJS3PZ5/6m970S/in//Q7ARi3yurqqnP9vO1tb3cMyEIWspDZcvjwYXzjN35D4/Pd3V288pX/FCdOnAg+f+c734Uf/MH/u/H+drtdPPOZz8Qzn/lMvPrVP4D3vOe9ePnLv72xPxcCDN1uF09/+tOCoNEf/uF/NZPhuOmmm3DTTTc1Pp8naPF3f/ed+Nmf/RkcPHgQL3vZt+DNb34Lnv/85wEA7rjjDvzlX/7lzDHa5NChQ8G53HHHHXPt9/rXvxYnT57EgQMH8NznPscF5j788MP4wR+cboxNkwceeCD4+8iRw3Ptd+utt+K//JefDz77i7/4y0Zg/+bmJh566CEXLPyMZ3wh3vrWX2NbLIDG5SyXrPBXTGlq7R+Uu+++G895zpcEbpdY8jzHv/gX/xde97rXXrQ5fuhDH8LHPvYxAMDzn/88/OiP/mv33Zvf/Jb0TgtZyEL2JMvLy3jve/8Y3/RN3xh8XlUVXvSir8Cb3/xLGI/Hrft/+Ze/CO961+9OPcZNN92Eb/zGb8DXfM0/dgv0mTNn8OVf/o+CzJiLIePxGL/+678JAFhbW8M73/nb7rvzZTc4IwBg7liQ22//MnzjN34DXvCC57vYmfe+93/h6U//InzmM5855/lM0+uxvO51vmbSnXfeEWTslGWJ17zmdcn9+DkePXr0nOe6kEdfLijgOHXqVPB3WxEZoPmgnDoVUmef/vSn8eVf/hV4/OOfiO/7vv8v3vGO38YjjzzSGOf7vu9VjYecyytf+Z0QIm/8e8lLmpZWSrhCIAvnb//2Q8kAtoUsZCHt8r73vd+9f2tr+/HN3/wyjEYjACZQ8Rd+4f9tBBlubW3he77ne3HttdfjZS/7NrzpTW9Opnx+yZf8Qzz96U/f03wOHDiAN77xZ6fqKZK3vvXXknpk//5DM/cFgDe/2esRcgcPh0P8+q//xp7mHEsc3Lm9vZPecA550YteGBhV5yJxAbaTJ0+1bNkup0+fxktf+i34q7/6q+T3W1vb7vdFeYIrSy4o4PjoRz8W/P2c57T7V+PvPvrRjya3+9znPof/8l9+ES972bfh2LHH4Wu+5uuCoKgDBw7g8OH5aLtzkbe97e3Y3t4OPuM+2YUsZCF7l+3tbfzO7/xuUOPm0KFDrTEZZ86cwTve8dt41au+D09+8lNw221Pb7gPnvjEW1uPR4Dhppsej3e96/fc509+8pPxG7/xa637XSj51Kc+1cjKeec734WzZ8+e17ixS2d1dSW9YSTPe94L0en08bznvTBwJ/2Lf/F/4du+7VvPeT5f/uUvDP7+27/9UOu2n/zkJ/HOd74L73znu/C2t70dP/uzP4dv/daX44Ybbsbv//5/a91vfd1nNC7qcFxZckEBR5yO9b3f+6rkdrfeeite9CL/YO7s7OB973u/+3saTfbud78bf/EXoc+Tp91eaNnZ2QmCQ3kw6UIWspDzE26tAsDRo95FELsLuHziE59wbgqSefTAvffei2/91pfjc5/7nPvsRS96Ib7yK79y3imfs8SGyoVwy8ZxKwcOHJh737Is8f73vx/f/u3fEXz+hjf8xDllAr3kJV+HZz7TB4meOnUK73vf+1q3/53feSe+6Zu+Gd/0Td+Ml7/82/EDP/Bq/NZvvWNmlg0/x3PJplnIpZMLCjh+4zd+E/fff7/7+3nPux0//dM/hSzL3Gc33HAD3vnO3w56J/zCL/xikBb29rf/Jn7/938PX/mVXxnsS/vzbJgTJ06cVw77PPKmN/0STp06hVOnTuGXf/lXzyntbCELWUgoBw4cwNd9XdiH5Phx7zb93u99FT7ykb/Fd3zHKxv1GLrdLv7RP/qK4LNPfWp6dU2SyWSC17/+3wWf/diP/ehepn5O8q53/R7uuusunDp1Cn/zN397XsGiJKdOnQoCNXm2ybzy3ve+N0g3vf766/HKV75i7v2llPjO7/yOBlP0hjf8JIbD4Z7nM03W19eDfi0f+9jfXdDxF3Jx5YJmqUwmE7zsZS/Hn/7pexxCfvWrfwDf+q3fgo985KNYWVnBc5/7nAA9f+QjH8FrXhMGfkop8eIXfy1e/OKvxe7uLv73//57nDx5Evv27cOzn/3FQdn03/iNt02d0/d///fiH//jr05+95M/+VP48Ic/PPO8PvrRj+Lw4dl+3oUsZCHtQmmxgKm988Vf/A8CH/zDDz+MD3zgA8E+z3jGM/Arv/IWvPnNv4hPfOITeOCBB1EUBb7oi54RZMV95CMfCYpXzZK3v/238G//7etw4403AgCe+9zn4Pbbb8f73//+5PaU5ZKSd7/7f+DXfu3XZx5zMpng1lufPPcc55U///O/cOn6z3zmF53TGG94w0+6zBnAlA5/y1t+OUhB5fL6178Wp06dwvr6Or7wC5/eKJz2u7/7TvzMz/ync5rLNHnWs57lyuEDzeJxC7m85YL3Uvmrv/orvOhFX4F3vONtuPbaawGY/OrUov8Hf/Bu/JN/8opGBDqPbF5eXsZzn/uclmN9AK997eumzqctLRYAfvM334Y58MZCFrKQCyBtabEAMBgM8IpXfGeQBsn1QJ7nePrTn54MDH344YfxT/7JK/c0l6qq8NM//TN44xt/zn32Iz/yr1oBR1taLGBSQX/t4oeBtMp/+2//3QGOZz3rmVheXp4zZV+A0kjf85734CMf+YgrTPb4xz8eL33pN7W6j+OibSSTyQQ/9VP/ET/+46/Z83nMI89//vPc73fccQfuuuvcOusu5NLIRUmL/Yu/+As8/vFPxPd//z/HH/3R/8RDDz2E8XiM7e1t3HnnnfiVX/lVPP/5L8LXfu3XYXNzs7H/i1/8Enzt134d/vN/fiM+8IEP4nOf+xx2dnYwmUxw/Phx/PEf/wn+2T/7btx++/MXtTAWspArUOq6xubmJj760Y/ip3/6Z/D5n38b3vOesPjWv//3b8CXfunz8PrX/zv8r//1p/jUpz6FM2fOoCxLnDlzBh/84F/jNa95HZ7ylC+Yu/4El7e85ZeDGIiv+Ir/T1CA8EqR//pff9+dx9LSEr7ma/7xOY3zhjf8ZPD3v/7X/3Lq9mVZYmNjA3fddRf++I//BK95zetw881PwI/92I9PTYc9H3npS33q9PmmFC/k0Reh53wyhLjgZMhCFrKQhSzkAsj/8//8O/ybf/MjAIyL52u+5sVTtuZlBK6cQlnPetaz8Ld/a6rPbm1t4aabHn/eWT4LuTCi9XyJG5es8NdCFrKQhSzkwshP/dRPuzLgX/VVX4knP3larIhm/64cefWrf8D9/lM/9R8XYOMKlAXgWMhCFrKQK1w2NzfxEz/xBgAm6P5f/asfvsQzurByyy234Bu+4esBmC7dFyMgdSEXXxYulYUsZCELWchCFnLOsnCpLGQhC1nIQhaykMtGFrTFQiJp70szW64sn/BCFrKQhSzk0ZMF4HhMSQpMNEkuQduJBAGm04WA3NcB6EhtuwAlC1nI5SNcJ8jom0hfxPpghi5wm03VCQt98FiSBeC46qSNoZD22xhM2M8TfzfidqyC0VZpaKZwhJDQWtnRFftOBYrJKx+ueBZKZyELubiSNjZCfSCtHmA6gEAI/90K6YFYL8RidIEKf3f7cH2w0ANXuywAxxUvaQslsE4iJeIVC4wS4b/TSCK3oMN/prUHEjoAFX4qXAm5bUSwCTwgUXY7jYXls5CFnK/EoGJv+oB0gYj0A9cDAqa3lUZtfupQD6RAh7a6AFpBC64TjIHCAYnXB8BCJ1x9sgAcV6SQ8kixFi1Wiv3JgYRA5hWJ/RmzFiKiUfn3SlVO8fDvtFUs4Yy94tJaQekKSlfQ7l+b5QN4xbNQOAtZSFO8PmjTBVwPmD1kwGDEeoCABX//Q11glo4AbATAI9QLsU6Ij0s6wRgpkU4AWoDIQh9cabIAHFeURIolZitE3mqlSJk7gCEEBx5hN16tayib4hQrhWAK8G4UznoAgBZe2ZhdMnN8puCUrqBUhVqNPfBgzIe3fIB2JmShcBbyWJZp+sC841wXcNDPjQ0gBhPtn4VHl4E+CFnPPGRDmU7geogDDgM6SihVNXWCAMRUZnShC64EWQCOy1qa7hKjWHILGjoGZBB7EYGKYCSuXGKXCpOAKrUvs0AEYCxIkQAUKrc97U/WER3XKznpQI5EAWSA1ksB+ODH5fPQuoJSE4TKhmShdBbyWJEmyDCxVhJSdlrBBRAylrPAxOxZNPWB1rVxq0bvMG3vfo/0AZ2W1goSObRUjj1NuWzMZx6UwBpITZ2w0AeXmywAx2UpKcvF/i1yy07kyLMepCjcS9vGWpBoXQduDyEkNJRTBv5l9tZGwIogbwIZyKTfdpZykyKHZBZSJSfIVMcqEeX241ZPDWnYF8EDURUMcuEWz0LRLORqE4FWfSA7kCJHJjtJgyM2Lui9IuHvb1tgKH8ng1mJrMGEKpTmFUzEqbYZOkJkAXCRIk/GiXHdVKuJAx5hcLrRCVrXWOiDy0sWgOOyEgY0ROYtF0aNGsXSRZZ1kMkupMjtIl0HYCEGHu579hILSAi7rdlGuQVf6dJtb7bLYPCBhESoIARkQzFw1w0A559VylgjUuSQsgitHDR9xgR0lKggZW4tn8iNY+cN+28BPBZy9UjGWE3mlmQMQSY7kKKwDGcTMKRcozHoSElbAOgs94tjOLU9fgvIoPH8GJn7KYQBDPEcKA5E6jyIB+NjeZ0wgVYTqw9CN+9CLo0sAMdlISHQIHpUio5jLYjJkDJHLvvIZReZ6AIAKjVAqYdJMCFlERxJqZBBkJL5WrVq/ONT1NozIWZ/oxgItJi02DBORIrCKoDKHVujhkIFoSbGKpsSmOrGyXJo3YHOVDBfAJ4BUSMoPbHAQyx8vAu5gsUyGiKDEB3HYhCA53EY9A5JYd51owdq+z42mQTAGCAp0EEAJXZr0s9prhgJCUUAh1Ra9NqlgtP539xICg0mr9fc+WeyoQvod8OA5KghAT2BaeCxiP+61LIAHJdUmoyGcZl0kMkO8qwfKhhLXxLYkEJCaQVlgUZdT1wMBlk8RrGkXSx+Ft6PCnjFkrRaov2EkNCi9r9HFosQMni3XYwIam8BRVRvKhjVgCPPiEgKPqVzV2NUdReVGkKpCZSugiAzQ6+aPRbKZiGXr3CgYXSBd5l0mwwgZJMp1Pw9Cxdqv2/G3olQYjcrfy+niWqpwxHM1x53HoZl2jiAD1qNz8UYNR3IuoAUOWo9SWbDLZjQR18WgOOSSOyPtQFfIkcue8hkF3nWQyH7YUAWvM9Uo0ZtrRilK9RqYix9XZkXUQJCe0vD05VhVDj/TEBCobJAZUaEusiQidy87IKDgzp4fyUklGCWCTKnDEPlGQai8vkpVJAwbEkuu+a6MACkUaNSY5RyCFnnqNXY+nd9wKn/Z0akPReykEsvTcMjBhrkLkm6MxKuivYjhSzltG3od+4Wib+nQFHo+QNRp4KdOQwkGiMwPuwlJNCjtTKMkMqRqa7TjTzg1LtgaS4LfXCxZQE4HjUR9v9ZGFnuWA2vWDJJ8RmGIlW69NaDJh+mt/CVKmcc2TMYUuRuQedR5p5lyMx4omkdeVcKzdvslwkbWa4rE+sBAjR1EDiWCl4TkIHy0PSdlD49dw4ly904WigI6YGNUpWzcijTZcF4LOTSi2jVB0FcRpzREYkvthf9fb6zc27SKRlvc4CDeJ94/LbvaXxiRPhnpHf43xoKmQ0WVQImE05a3aAkpC6MO1dVELoyMR5aWuCxYDweDVkAjosuTaDB3SakVFJBlAQ0gqDQRER5AARACsoHYAFAJnLjfmDuDz6GEBKZKMzLmuWNtFjyCwuRGXAkDF1pGAsFof1in7wKdm4EZOi4DTeJkG6uPhjM08O1LiG1GYfYEhXFnPhUXJMmqDPl0m4rMbKgw0S2LxiPhTz60ozPyETHpbXPAhmp+AcD1uske8CNieAzFmcRGxdtGWbzAoxZkgIS07Z1vyfScd1n2m0EqWGZ1cwEuUvSYQpSFDbuK3S34IIxoDw9Z6FTuCwAx0UVsmBy5zKhwC8K8iKqVMrCxWvQYq8AKFUG6ar0Pa9vQX9LFMyFYq0l9oIa94YJ8iRp5MuLzL2sRJNSMBjs/DjYMHEkFiSIsqFAKMaCKw1ibpJXzG6Xo4C2AKi2tT7cnLVCjRICLGhUlwHwMldfuqBZB860MsrHpdP5lLpFkOlCLq54feB0QRSfQboBSMcmAC0MRsvjGsZH1Raoz3KXNgFKaty2+JB4rLZYjbnZkQgESXY8/l637iv87xq2iJi2gfhK2jpAE2iRM/CxYDwuhiwAxwWX2B/bQSZ7KLJlZFkn9H9SyqljNMzCqlTpimFR8Rs3uovBqA3AAFcQYaliCtqUYMFa8ECFskaUKiFkGO9BoIP2zRBZGfZvReAHmclIEZ4ZCVPWfDArV0Ia1GfBHk/bgFdRIBOFZTQMoKjp+lB9EA44ONhg6XicSaHvTSXDslFefZHdspCLI57RkLLXYDdT7y0Xzja2BlpaIDFPIOb5uFzi+I54nsl95oztaDteDDSC7JUowNUfMzScGrFsOgxQV7ITVDg9vzT7hc5okwXguGASBoJy10me9VDkSy5FNCjARQwAf/gt2GjUwpgjawRgYCNheYTulhafMIGDOL6CKRSlyxA40D6yCM9P+8ql/HjOxcNjPRwN6i2njKX6iUhRpqoZmqwZFfxtNzagQwEik5A6Z+lzkwh4qEjR2AEWspC5Ja0PMtlzMVox0CAh12XsGrkoCzot3vr8gAiXafU7+Dlxl8i0fQGvi+Y9Ht9PJcAY6QkHPCzjQSn2IgYeLu5rwXicjywAxwURsmAKUOU/8slScS4hTKwDsRepQlnQPnuEwIYZnQIiC8dOcImrhrqx2SIuW/LhXewIY06ILZApZWjjKOJqoP5KMN8wmx+xGFoQeMpsgJehOE2RH18WXUBCMXAkRd7I1KH5kFDWTCy19kG1sYvFMR3EesQBpotCYgvZsxjXibC1dETKjcrABoDmYsiAelslX3+09qDOGOhPywIJy4+3L+5t1UiDz1KuoDlSa9sqF08bP87ko89c5puL7ZBh6i5l1lCAKawL28Z8Ubq9Zz0o3mMBPM5VFoDjvMX7ZTPZM2AjCgYlocZocVGtkClQjbRRYiXIJePqVLj4CNlQAnQ8DjrCWTeBwSwhdka5lw72pW1aKn4OLQrTnkcmauMqIdChayiYyAqhZcBqOAXCrimvOUABrT62xAA3CYk6qvfhmBOdo9aV7d9QQskKUuWo1cSAHl0BemJdLTUWSmYh04WBDdlDLnsNQ6HtfSWZh8lIMSOz9vHjh6DDL9i86nA7MJkFCNpSXuNjAM3zSJ271ippNE07Hp9/UDKdDA2W6ss/kxQwb/WBlCbNXqkK0DIyRBb6YK+yABznLJ4yJeVSZEsosuXg5Ykr4bVFf8cvGgVWCiFdYCkv/EPoGwKOQTGf183xojz5NmXH2Y3wMyPEbFAga2MsFpDWFuBGbpJcdg1Qk4DQGWpdGptB12Z8212SB6hmDnBIPx+Tz8LGz5CJAhkK1KJ0rIkvAmbpWQJuYFSqVTZ51je1POqJrelhK5iqyULJLGSKZEHcVp71kGf9mUGaexEfgzU74DLFMgCeQYnHaOqkNHCY1qsp/Hs+F81U1y7LpKnRZEhk4toqKEg0zzEFPPjnFKdGBkwtTSC+FDlqMYHU3BBZ6INzkQXgOCeJqgFaS6bIltHJlgEAlZ4A8IGgTlFEzEIqICxu30wlwl2hLcozl4UJ+OTxGs49EebmN9PbwvgLBQAifBx47Q9eRjxoqkRBWYBzG/E01maqqnRBpxkKaGHiPSTM4q9RO/dLrcauf4tgYMO7XYgL8cyQASrNLJgg/gTKqRtv5RSQgO/VYDMHqjqHqCUqJVFrBeEYlYWSWQhJmObKwUYmTfuBlCEwT8qrP0I6Pit2abbFT8wrs7bl7p64V9PexokKd1mp4XVSMG6iVLqwLhK+P+msNtAxjX0JrikUcmEMISFtkL02ho8pl77QB+ciC8CxJ/GKJa4OWuTLTsGYICQTE0CSalSmtYIWyhfaCopxMYACOEue04pxWWMTi1A2FvrGWbCcdR60Sou7mx8DD8GcE8yNaaxWuhiJuNMjgQaJ8LxIFFRjHxpfofTxKZDIqAR7oiGT1goVRqhhMlsqPQ5dQDCArgKQwV87l15s55WJwsw1s/esAqghFBZ+3IUAiA0PHiRObQliQ6DNmk9lpvijNMFG4JZgrGccA9HmomhjORrbxgCm1YDJ0AaqmuzulHTbGEDF849KXHDQEbOzHHTwz7hwUMJZW1dITJhiYUKYBpZOXyi50AfnIAvAMZfM7nmSyS5y0UUmClR67PbkMRZxUKMpT145ABArkRRAScVicL8r78jKX3CeHltrw4qAxZNQ3xVyi1CA67Q4EWIVBDIoVbkxYoURnlcGKhamYNJqfZ2RqGCRVQK17Ysi4dkNfv4Uia6gAJsyC1SodYlajZMghk4lE0WS8jYdenNXORW53U+xvBVXJGihZB57ku554qsE5y52AwAoNTMV7BiMuodeIzz24FxlarwFmu9xkIY6dxnyhNGTYDemzzPMcnEAgYGOlDRdxBEosfpURiwtgQ47WfO9yMGnq7UCtLSFBCcLfTCHLADHTOHpbb5gjy/iVQSLsZpFJbKgT/4iKpSQ7Ha0+zWbfkyy4Cnbgo9BBcX4wmqsAT+e+0efCgTZNHTMZD8H3piNwANCFqRxDqhR2/0qNUCtq2RWDmd3zNjSxmK0KzoDOnymD3dphdfRkLfkykkrxQxSAJmooWUXmTRuMqFNgaBF1PpjUdL6wHdzpQDRvRW1ig/RBN/zjwmgAfz5ODP3bVmkOehIxoTNEYQ+k1GZ0/0TsBJsLimmY9qYdE4uM49VK3XMknWDc9ChoZDpLoSWqCGhdR7pg3Yg91iWBeBoFcZqWOXCq4Uad0beXBhRN15Y+o6LJKrOBjRKFEH9CB/XwXyZLJ88TnVzhcKYi4NQeliG3BxTq9otymYMW78iBhCJl1WK3JyrUAHASZ1n4zrYhV5YYFJrU3I8LtlOrEzsnzZptHlgkQChtdeWVpw6LyVMPEe77zxDJrtQSpksISHd2LUySkap0cK6uepluj5IGR97Gj2RpZEyTLjMw2yk0lhbt21xv8wrsfuo8X1cijyS1PmkmJ55XUKzYltSVUqp2CHtE6cKSw3j1rWxOUL5IopCc32w0AUpWQCOpKSL9gjuk7VCVr2phBlmTyRHZimZSgMu6DGRzdG0NMLSxG20axwHQkIBn5XNLecAhbtOSIyLxPtfSanSJfJxFc19U/PgVKyGyYiJg1BD91CYPsfdQsE1YeIL/WSMnvaKzlQg9Z/HlK7rGWPZDVcHBOb8IeGYKOpCWda5UTJqgkXNjqtNmkBDMJ0Qtw8g4SnpQDMeKsyemA4IUgvsNCDS2H8voCMR/9F4rxPn22A950y5nwc0xcaC/7xZ3JB3izVTmR6E2zYnqZnbBl7PUDqt1srpAyoaRkZIDZPFol0Wi5nJQhaAIyGhb5Y6OMIuuCTe2jbdSCFhWnII5ZRL0mpH5hcxAdSu2FcY+T1v1T8eAOr+Bl+4lV3czUtDjEJdTwIWwGSLSFaYy2fU8O6McfqtyvIGQxH3NOHWn4oUWjK+gmqHoAyOnVEqLYxbhtwl4fXIbHn0Dnj/F/d9xBDFZZIp9Zfqf0hWEtkAFHsMkQOZcWcJZCghbfT6wpd79UiL4QH6PYyzcu8cGFNIAc4pi5/tS8C3Tiz49LuLYYgyNvjx98JIcEm+hwgX+7hmRswuOJYxAlp7mcO5VDudF3jNA7q4i4VY1FhnubGEDekQ/rss66Cqc1T1iOmDhcuVZAE4Amk2V+KKpbE1S9E02RQVUwyKbcOsdmHiF6SQdnGrGYr2MRT+GKkqemGwlaaaHbIZoa6FD7qM4xporNarwRZ76m8CwLqNfDtoSJ/+q3XtfKrkd6YYEvq+1lVg8SWrkxJrQgGrQiK3NTYA2BiQKlAGFKcipWUgwOJbOPhLnDOBDdre39ewzHTcLVdrBRRmzKrOXc2OBdtxpUs6VoMMD57dlGozEICOeORp8UcXWFJZK/Ms0ClWgz6PC/21ZcKlsmZSDLH7PXhH5y/tzs8pyfhGx2yyUdMDTqexJY2SBGRs2dIGSsWxHY9tfbAAHE7Czq5SdILOjUCM9sN+JBws8PRKEkprFXbxV5oWOd96PpwNs5aQVlK+YI09hlCNsZQq3TxjC4J3beXuF/4CZyJHLpdcMS0AmOhBQNNSCqlJbZVuHlQ7hMCKCRb1MSkB08D8qalgVanbO8x6hVBDyhy56JnCXyhRYWTOWIdumwCwwV83iocRWkGL2nbXDa08Ahu56FkGSCIXHUxkB1XdRVnvLtiOq0DaYrc40IjZDYCxEVFA4zRxixt/FqOFei9uFJpTUhel2IuWbePx4hog8TnsJevE/d7qfm6CjnCMMK4tHo9nmkw3rNJxIny8VAadO0akT4T0zFctchbrNYHWJR7L+mABOACQNQMhA7BhAi8z54fjBbmSo0SLp4tvcC4WiRq+NkesVFJUrNku9P1qbdNjWfR4qikclzgYLcs60RVo0qVBwBR4SfXQXUQMQntWig9O5fU2uELmvml/ToYt4XEyfsy6oZgNDapc8a82VoqfHz+n4PjsXpg8/JQVG94XV+QJlJ4MLNLlrkQhfeBjt6RjPCXiwnztozQBfCzNha4JNoAwDoNnaMSL7bSxZ8k0yz/OkjnXQNY4niSOkfDj7yEjp+WYqdiVGBAFXaZbxphVOj0WYyRFBpJ1t5smcPaDx6A+WACOIGbDWDLUCyXcKnNZKa4ZW0skdtB4iblYFKhGBMJtEqCjDbG78VQY7DmP9cPnLkWOOIWsHTDUUCjNvGAYiCqO20iUUa71GFqaToycYYnrewTXowW4GOaBFwALj60JaPAgrkgRcX+5j1PJnPKGrYnS5jvn7i3FztMFmlpFk+mui+tRdJ8ew0rmyhNeRbjTCjYopmsqnZ9y30WFpvjnXB/ECz59xp/HZEZcS0psm8Qukr3GUaQYk2luiGQwfGPM9tLr4fs5pXIoA2VAugz6NJnnOtD4yTgPy86Sa1vLjjVCqO39Y69C6QJwkJ8WXonwgj2GZgyVDH2v2cIXv6yx4nALbpQVQpJSLrGYDAs/TjhebSrhTbG8BKRzc1AAJg++pMwRPn9yAwHETtQuC4cryBjw0N+1GgfBcR502GMmUnH9fLNAsVMsCt82ZktqNTaAQvpzi7vbumuZ8u0muuomm8bBdpqN7rGwQcZS5gYUio6b8QJ0XAkizH88/Z11euUulNjqnybx+51eZNPVdt3+ifLes441S1KNF1MZcmbc6VkyDbamxSCZNudQh7T1ckmBuCnnzLtmU2DrHGAizmqbNvdW1pvWFdIHElBW5T0Wy6I/xgGHCILCeBtp59PkCsZZPr66ZUAP8tiJRGCmENLFJzTrV9T8jyAQjbt0YuH7kQvCzzWiQSGRyy4yYaqiSmSoMGrMl7toal2bzBVRupc1frl5ei1/8eIiZLGFQtchBipu8SaAZ8u3O3bEzq9ZGMyzQzU/H7sfWaga0gSUAg7ATLMKeYYOV8aKHZcHmgK2zopQLnvJVPtYuFcub/GdXkNXiu/2Gr9bYZ2GZv2clCQBRcRqzFxA2TbnU/eD75+KZUjNxbkYpzzCKbAxDxBK6RE/52Y6/DQmhc/PsczxtWPHCXqysFg8Ah3zzJvmwvWJWS8Kpws46DAK4bGjDx7DgCPKSJGdAGzEDzxXNpQNQVUvAR7R7FvQ0378p2MmMJ32JDaD9rFTtj/Crob+jLzbJ1aI0+hBCpak7BE6D+4akbIwAXR28aUAUZO6KgNFGfxNc2QpfG1ui/h6xdQ1XTu6vm0WlkIFfrqOgYGE0Mo0vQPl2qeZIIq9aVPmnN2hwFwz78xZNAB8wKliVVCd5nlsKJkrQ0KwEcRxTQEbDVZvXqu44VaYw7UZfc7ZOnMG8wOPwJ0ZxYLEC2brHJw+agbIm22bbta91ARxh5lSAp4fr+28gm3idN0pQb3O6GxhfID0ucfXzrtym6DDx3Q8NiqTPkYBR7tyyWQXmTQBlWGalqFUiSGQQkLZuAyqZ6F0GYAN2o+YkTi2o02BJGdMCzGLjI8XcR6jQeIribKX0sY6SCFRqbFrcEYLZ2xtEQVsrpOx9jL7vdI526bdyoiVtslkSbuRuGKXPAOG9YrRqCFRTL1mBDRoe4EMWdYx9yG6TmYBUZ6l0Gx/2zxOIg8o2fh+S1gWhOp1SGniV2yzOCigFsZFtYAal5MICFGA90ThbhSq2gt4F1uqhoZiLo95QAfJLLDRZvW3fT6vyyEFHsznnDVNxKUwcNJoOjcFbPBt5gUdMUM7bV4xO+QYp6hPSjA/627xho51n1INjojpmCbTQIkQtpKzBRtCySjG67FhgDzGAIfPrReUjdLogRA+eCSuPTxjOCBKY62yxYn/i632Pc00aXlP77IYVAIFggUxsF5kF1DmHKn2BPVOSQmP55B2HlKYjBvzu6myBxGDFP87d4/4ccPiXTE97c6D+bWJHaLfU9eh4aIhNonmxywzDt6kzUyg1N4guFUr20Oh6f7hDFRATwvlGI4MHsQq9dj0316ekoECxnmAaOBGicAGF07Dp/p4zJK9xFucj7jYr3h+LY9fclu3SzPThC/IrezIHq5NzGDOA1A4I5wGXS0gSDTHD2vvnHtBtWAeGshkx+mEIMbL6bWrWx88hgBHmI0iZQeZ6JgaD6JoLlLBYkk1Nzy7oKFcN9JSDRulwmlfGo//JJ+eSd9sp0ob46CeATr8/AhoxBVFpc6hVAUlTSv5ODKeK1gOohQqk50CIIOyFVJ95dIY2DjrL8rsoewYwPQxEToLXBJ0XLI6XFwGVMAmGKXSjGLnx9I2/kSimbHCj0P7ZqKL3AbS1mqMSvomdrAN7cjq9fcxNy6c1P0jyynSfQIStRotQMclFdIHhUuFp+y0NrdqSnjtCR7TA7QvvFNnFbtREeqDtt+b48Tt7tOLLQ9Eb1vUp7llebZMW+bd+UicrTLP4j8NqMSuHdIjFN9BJc1JOIh0rtaoY21cBoHPIWZRta4h4bvOCkhA22aarhT61asPHiOAowk2ctnzdKkIH5iUTzQAJNbCVrpCpSeo1cTR9oSu47TaOIjUvDiywQz4GTeDRt3vOvycz5O/IDqi/TUUtPIVUaWuGgqOsyQmC8MDCaVKY+VHvl4Httj5x75RoqF5LxRqhJSyvjgY4HRvSoK4D/jqpsqeTwwIGgoXxu0RtKqXQF03A2klclb+fXrFRafI7XxEZo5T0wKlJ48JJXP5iY/fymTPtZVPsRpuD/eehF1KAQSBiOcDNmbOegZjGjMPwTlEwZKpGCoeN5YaI7V4k7uVf5ZK643PYR63SoPliICECV71DHPbPGddr7hHVbrcQehGMvc3kboLOVVf+WaTtg8LTGHBCjawXI2uan3wGAAcYbxGJntuUeXdXgPQQYuzs7R5dLRCDcMYVGpsrHtVuVgJM15YICxOdXPHlNIt1sqOyfenHG4T0Vy7hV+I6EUOXAS+oyEdi7+QSle226vpneLmIyil08SpAECtJYTOIuVUW4AUpsRq1BZ0UN0NwyyQEDtA37vPdbsrx5+edMpByhy6NiAtFfchZeFYFJ1QaMTW1GoMSFsOXUjkosvKr4euIOhQKXMFRQxYnDLo77tXXlLkkHmOXPdQ1V3DjmEHvt/CYyNw7NIK6QOjBwzgYPoA8zEbQGgl09/AbLCRWpBj4B+7JIOfUTxBWyxH8JltqCaYSzEJOrgrOAGk25iX1HjJc2+JG2nMewoY4UA+aI8QbTOPcOCoIR143EvNjhTQm4fhISNEarMO1SI3ptFVrA+ucsDhKwZSvEbKkqG0S4Ah5Qg1u9gBQS4XU4iqYdkjXMC5W8J9x47JGQE3aws2eBv4Rs0LNJWQs8wEIKBcS/ogy4PHmqBpQWQiRya6bswKTYXXABs6BFT0EptCZxa0iHRMS0qBpZQNXTOta2ip0mCJgQ0+pgOR9h66eiU6s4yLSRHOUEChdtVgY0UeLyyUicLdP/Giw9kyITIUsg8AyLM+qnoIDYWqBqBH0BTBtpCLJN74iFlO3sZgmjuC/k41WUy5FOL3dB53QApAxKwL1y3++Cqw1I3LzwJiqrEjzTyJ9TsXl09bjATQdLvMC1BS38+aB4Go9m2m30f3OY9BiVir9rmGAcHT2KVUIK1zOVsjJ5Md1GoCAFe1PrjKAYevGEgVRDPZDViNuBNqrUsoAHxRhgKkLILCU1qrZNxCHNVOMs0KSb14mewiFx3f6KxWQCJWIA52BcgtYhfeHKjVxDExbVHnLvYDCrUeh24a26Qs1UiuTULLQQXvDbE48blw4NCYowNw0rl1wsXcK/RYWcTBoW47VUJluVXGBmjUunRZO3x/DzSVK1bmGKkUO6Z9gTMhPM1Kjd8K9FGJrjuHCmrhXrnoYowPIXITvxWBjVS7c8o2Afxznyxzv4daE7EFPGuhjVnYtuPELgcJrh+64IXwKHYqdjdeTOEsSmruKReK68yMqE4GB3QzXK7J47SwRG6sxCsY685ptTnazjNkT1m2kz3HXPYti3516oOrGHCEFQOdG4VlSrjUS/tiSuvDLyPWgdPo3gdKFixDreRK4XEJMbvBXiilfRYE/SMFIW01UKkLt41yrdN9yqhbqEXmzoPcA1ooVNkIpRpiUu+gVhNjVcfunWBOvsZFzBRMvdozXrBWtxIDGvQ5Z2T4sU0WiQUWwte/CBQxY4CSgbAtc6wwhrKAw1Un5QFgLCbFxd1YIEVgVevapdHRdqnjSRSmAZwsUGXkkqtQawU8xgoBPXpi9IHrlSRzBzY4w9lwbUyhx2cB75RLoo3B8GM23Rx+39mAn5g00g+57CJHz71TJvDb6IQK0vT4aDFC2uaYmhvfrzU9lDGc04BZrN943BeJ1M005Hju5iuuc8j1Ob2pndnPs1jzukhcPFjM5Or0NQa87jAnBejMZtVBMX1w9bhWrlLAwTq/OldK1zEcGStR7YIZRQEpjAVN32sRolGSoIocsgBhcwBgvk+/XJRySfEfwWJsrRMamy/MptZDGCdCL04mu+jKFSxjP5bVGiQkalFhkO1gS57AWG1hiDMmfqFFyFUkhIRqseQ4u+A/Dylkuhb8OyB0WcX0r7sPEQ1M+9csGDYOVIutO065miqgxhKipnW0Dc+8UVqCso94XZK28yEQ6r9vxs/wf/wa0DWU8O4zJXI81ps7eeFm6/leh7jujs9MI7ARL578XrsYoj00FJsVg8BPadoCnWI3UvOkY9L2mQUbPbGOPtawpFZAAeUjMcB2dhoTuYNRvZV81ttknjgNmtc010sKdKRco86oYsUOCfBLbZiO+HoGY0bgIwYdbednvme1OuaQeVJ33Xkm2DTJgFaoD64eXXAVAg4Wt+F8tX3kWQ+Z7KKQfc8+gLMFlKmAoBocECofvlgYhZVD65YIb5FBtLx0BDYopXT2WfmGZG2WRSGWsIKDOKaO4li3j44UqDWwUx7Cw+oATmUPQ+sa43obdT0Jzsm5VIjhgASieYkEW8FZnhjZO7eJMIqCXDrxi04vWCa7rrAaAGhRQ+gKNSjTpk6Cl9h36tgZIU2J8yioNrgPQfl2dp91HYzPgZ3L3KFAXKEa6XSkBF3FUfhYkhqlPafSlWCX9nlVSsHk5C/Ev4gXYBxrfGTCGB3OlcKfFyu8Pw5JI8YL3mLmkgLZfHvzXd0ANPHYexUP3DNkFlT1xDpW9QEc0gdxIO+hk5kXYaescUKtYVOewWaWYaJ2UCFM+UxJHMPEjxvPhafbmnnN5zYisFHIvgFN6AXXit4dbc0PDjpmxoTETMgUF1USeMDfw3nBRTiflBuZFWlkjGpTH1wdrpWrDHDw9FfWC4GVJ24TpRUQKZPYnuF0ud+uGdjF0StfqMwMZfLloO/csRCWP/fgpmlhCGF6pCyJ/TiqjuDWlWU8aR3YV5jvT44L3LGxH/cNuyhz0ztlgl0oXTWuySyfaNLvGV2PVMS/1gpKGoDF553JLvKsj262ip5cRwdLKNBFhRITDFDqIUo1QAWgjpV7wCRkjcXDxb9oH9Qbz4muMe/k235uYTEgvq0SMIohYn+Skf50TB12zZUihxY27XZRowPnDzbCWhuC64IE+DRHbIKNcMQQHLR3jW4aH/HvMfCg/eaKA0lY9Zwl7IpVrOnDOIaDuGG5hxtXBFZyDSGAs5MM92zvx8OjJdwngC15AgI7qPQYegbwaDCbcbBkNH9vfDTdKbHuJBdQIfvoiXX0sIKO7kNqU2RwLIYoMcYYO6CsvgbTMYfMWwEWaHcd7VVilqzpzmkC16tRH1xlgMODDSrs5YJE7Y1W2iwsvg4/NTxL1+snaTwQ/IFhKZqe0jT9RgQya6dVDsEDcMwHuUn4g23qZ1TO9RFb2nFNDykKdOUaDqtjuHlpGf/goMLtjzuBA4d3kRUaJx5cxfseOoKlM8sYbd+CB3IF4ARKNWQvv21MJsICPvz4cTpfSrlIUTjXVS67QYdVXiTNxKvk6GTL6GX7cQg34qDahxXZQS4FdusKZ/UOtuUGtuQJ42sWFhgkct+5r5cKkvEaJBxsUPorVUqNY21mSezTVboEBZHGYMM/F8YF5jJ4LLvh4lBsNpEWHd9n4TEdz3EhztmDDV7Yi9feIYmBxrS4hsbncy5gjf1deqpP7Wxkjs1wzXC9QeyGcaWs4ag6gCes9fAF+zS+6MAG9i8PkUmNE1vL+FhvDXds9VFvXAtIYMuSSRXGzhVJ0tbZ1s0hil3g50K/B27mxPubWTdXT65jSezDQXUE+8QSelkGCYGRqrGpRtgS24AAxnoHUhjQTqCDZK/VXuk8+HnF5zo/mGlnTWJg5pthetcxd9VfjfrgKgIcAq6wF0WhM2vG3FTjXxeRcjF/NLu7pnKxUwxFw9K2P8mfbyqEhi+jlAVr4lO5qqM0L6VK1DIMGHP+R8hgEc1lF32xjqNyFU9Yk/gHh0/hxq9SkE+8BcgyLP3dvXjm/9jAbrUfx4er2KmPmmBFKA9qongEN5eEC4RfJyGagWdmUTfsUi56PsYkK0CZQACQiy6WsoM4qm7EU5cO4glrEgc6CoXQeGjUx52bXTw46juF6PrCcFdYxCiZ72sHNnjHWjoXSnt15400TTzLCuH3KxUU6+4bn5uuoeEb5fH9yKoB4PosLOI5zlVIHyS6vkbsBr1LYYBzbHFm4QIagc6pM0n462kc+4uXxG1uY0s4Q+Kq+IrCuFb1PhzpdXHzCvAF+7bx5KecROfaAuhIHLrnBPAJQGMNp0d9DEYHMZYDKFlCK/t+MPdDIwV0BjB3QCjBQKYks+7UQixhRRzEIXUEN/ZWcLSfYTkHpAA2JjkeGeQ4OcmNS0VQH6vwmgDnBz6mAQbzWR08C24/tAOW6cdo71vjWA4oSAs6tC5xJeuDqwRw2KAw0fGFvaKS5V6Z1MG9irMoOAXIA8qkDcSiCGJuzXL/KaU9kmuCGoJxgODYD2ruRTEl7B91PI0XLhdrYdPZpCzQkSvYpw7j+vUCz9g3xJNeuA393d+GOrPxA5/3BDzp4T/Exgc7uGd3BVsbhzCQ26jkyCzADnSE8RFtlj+5d6iEuIBMxqGQi0IyV0UhgUybLJqeWMMxdT2evLQP33D9BF9y2/3o31pA5BInP6DxJ5+9Dh86s4xy6whKOUYpB6jrygMYG/tBvV0A4xqjINx47vE1jc85PsdZDFcyNofd51RAIsVs8EXNfa8BKWEKozkwOrnilcyjL74/SiZ7AbPB620ATdZuVnrrXgJHaftpGR9uYd5DamdqLhRoLkWBLlawpldw7VKGz1sd4cmPP4n+C64Frj8G3SnQv+chfD7uh/6EwIPDdWyVK9hR+1BiiFqUhiGcc6FuzVZpcxtH2/rMuh6WxH4cUUdxQ28FTzsg8HmrI6x3jHHwwKCHz3QLfHarj+FwP0phXCvkWuKuHWJdTBlxavw4DxAI3/c2VqMNjDYMlCmxLnHMUNJVrwEhFDJ0oCB9W/srWB9cBYDDR6Bnsoc86wW9LkLrMk2XxosTt9x5ypw/YvPBos+CVFJmKdMCRtShQJhSOy0fPgYb5KahQlI9sYY1vYwjPeCGtW1ktz3OgQ0A0EevQfeLDuOmT23iyOkVrGYddHTfuRZS1yS4NlEsCQdYPKCyEXiFGqaFPXcVKUhhKnv2sIIDchk3rkg8+9YHsfKqZ6B+2tOgpcTRW/4Iz/qFMzg1OYR7d3o4hW50TaRjUDJb0VShBkQ6s4S7q7j15bZDSA03LDKRVjT8OrnrI/w14RlRvOR8mwi7WFo/XMR0XD3pcRdPWDq8bLpRYiDJZWovkCltzFtnEgcRJ7cJ60JMA7nBfvHzSzEQoosCXSyJDtYK4FBvjP41Grj+GNSN1wP9PpBl6N59Ekfv28H+U2tYzjN0xqbiqtCZGy8Ostyr7MUNkYsuurqPVdnF4Z7EzcsTPPnQGayujKE10D+5D4N6DWfHGfrDAhkSxRoRsSuJd7Z9DnFA5x5qDrFslrYA4JTbqX0uiSqzkBDIYapMX7n64CoAHNKCjSUU2TKKfAlAaKHME4jFsy20rpEJXyDMb2N9bSwqOkb4CgqZS+EyaY9KwLS0h7RpuQWUVsaNYt0Rwio5ntVA8R6NwjciB0SOXHTRkStY1QexL+9iX0djfWUI7L++eX5PfQIOXf8BHLxPYbmQKMZd8Fb3qevFM1GS1yvhhvEMjWnDDiDIMlFQrgZFjhz7ujluXKqx/qI11F/4hW6c+mu/Erf8zS/j8Rvr+GuZmevEFoIMBmwUou8YlEqPTdEuRIt/is5mzBHd2xSQ8icG8GwcXpBIAUHJaArc4zn2QaQ7o9Jpv5jWl6KAyCRqRZ9N7M8rz6p5dMXrA17or624VyrFmkvKgk9Zs20yzyITuyBmFZhKgg33zBUodBeFkMglkEsF5ML4JTod6P4SxPISxHoP/f4uOhIoJPEjLMZiir6cN7CVz7/pnmzGqhToYjnPsd4BjvbGOHh0F91DZipHxjs4sLuElSJDkQD688wvzjBq7Je6p1GcTXzM1LNxIYSznhCm94qAtO0wpGU+rzx9cIUDDmEVfwd51nOpr5R10L5Xk1IFokCo6IGdFr2u3aJRu5gNt/iIDJkNUiWFYNJvM9RtUeZRISKJcHHMYF7gQi6hI5bQVV3kwqTADoZd4KFTzXMej6FKs42a8ozy4EpiXXjMSOO80f4CurgZFsQJAFoodKSy10GY37uJ6qK5uYgKGhVKtzj462AsHWn610KhDkABt0rblLZnPsIgN5K45Xh8HahnS0DtJoL+4tRh951kfWaiZ0va4m2Z7KFWwNWUHndxxOsDAhs8W6o1lqIl6HPaohrqjVSWyl5dL20uiOnj8DgmSfrDnuO4BrYmBcqTCsWJMxCrK8BkAmxtQ2+OMBoVmCigVBrGREicR8y8JJ7xvUjgttUKCN4lc3ekAAR/xpWG0gJam08VtI9/ilyTbS6eWUAgmU3GmKfU/Pnf09wvqXghLmSUOn0Vl2xnf0uwGlG6uuL0wRUMOHzcRp4ZdiOT3fato1gI77s1fQVSbgGtfRlr9xkDHg0WQJrxa10hQ0ypUpaMgiTUHL3gBDao6E0zWKt2889lF4Xoo4cVAMB2VeGz2zn+6uHD6L/9QRyb/DfgOU+DuvEmyM99FqO3fQQf/uTjcNe2xKnxABM5NDEDxKQof52cS8p2s6Wf/FqSn5SuDfdDxy+UcopBOeWlUaMUY4xqjYdGGc789zM4sPpe6JuvA+oK6o8+ig+//yg+ttHBI9UGRmILFStYJkRmAYaBAzV8iml8z/n15axD8BxQrwkguKe82BhlvDjFIEOgwp8NArxtCwZ3Z5H7hOp7UOE1aM90AEClK+AqiFS/OCIgRAHqBh0GiMroXfSsHl9Q44JQ6ZgL3qTv3CnteOGJ3b1JEBMxIE5fRSxeiTF21AQPDwvcsbWEfZ8+hCf170bnnpOQB1egjm/i1AcUPnniMB4cCGxMSkzEKBiXFsHA/cgYvGkxLknwFrFCqb9LTDCoa2xMMtw/6GPtvv1YPTlGVUvcu7GGuwcFjg80tvUQYzFwwed8DPuLAyKz4jfizJHwXOZncqaCUxZoSsdTlv2Oj+dS81tAhxAm7s8UKCR9cOW4Vq5gwCEB56ftOuoUaEezLpIbLJBIGH1vcp3DqqEAWtwG04rj1JYdCItB0Tim9gIFOFYmWyGi24M26QLBiyVtMSmBDAV6yHSOUpTY1APcuyNR6wKlvg7Pfetp3PCX70XvyX2MPzfEX/3tdfjg6SV8dnuCk2IDI+yAyvcKHbpvXNpoBDhS1zWO3G8LtuSFtAzjUaIUY2yVJR4Y9PBnn7keX/DGM+gWJzEYd3DnxrX48NkO7tgocUoex6jebDAQGj7jRGsVlSOnFtJxU73MBZqaMdKVC+PzVcoHotLxOF3OP+eKpdahUuFBpYAHHlrUgdvGzcmmGZvnZWLjOS4MbXu1CcVt8No7wfeReyN2ZRBL2Ra/4N1e8/n1/bufLkzVdh+ngY1gu4RrR0OhFGPs6CEeGXRxV1agK9ehPiZw7J4tLK2fwM6ZLv7++DF8YquHh3ZrbOgBhnK31T3hzoEZFfwcYhZxGuiIxehJ8+6OxRDb9QSnRjnuGRQA1rG0pVBqgQeGBT67DTwyLLElNlHqoX13vaFA2Sn0Dk3LUJkGNM5V2tgV8x3LbrEGCwcdtCYF19hVSg6rqZpYQAUtKxtEeuUYIFco4OB9EfLA2tZQrmz4tOhwTmOZhVcGNxWAH3PKg0Tje2vKK5kahvbyA1bBfMiabQQc8gVV+/gHANCigNQ++6YWFWpUGIldDNUI21tr2CmX8NDwMK47cQBHP1ri1LjAR89muGd7gnv1CWyJk5joAZRVIpnIXfYEZX6YqyydpU5l3t2ptJQ499/7Jmru/JhrSOkKI72Fk3oTxZbEpC7wic2j0AB2S+DkSOGR0QDHxUls46Rzx3AwYwBGbeJhWGMqii+JXT46WnD4PnzeABrl5nnWDonSlQNa8f7EeEEhyJbiQouQeRapD0+YckfxB5nsQoqOZdwWrpVYDNvp01/NZwn/ehT8mVpQOeiIZfqiMo0ZaX9X/DZhhgyBmwDATxlDaQPAJ2KIbbmBh2uJcmsFg7qDR8b7ceD0OlZzhc1S4r6BxIO7Cg9MdrAhz2KCgavgmaxBkghunZbJMY+QjjYsZYkJhtgU28iHEkAPJ4YFMmlcQ2fHCifGY5zCBnbFBio9amTGxQAjZo3cuVxAsJEKLp4mMQhtgA5qEMtABxV8bIsXUXUF/Sg24TsfuUIBhzQ0ou0AC5AFHbZNh2iPP4jFuDPCSp5ASNPF1mxgOUexF6lcfr4vH19r5VJMqZ4+LZrU5yDwI6IbWPclzMu3g9PYkF2cHh3G/cN9WM06WMoK7FY1HtQncFacwBBnUdXjkL5DpIRZ0CtVXG28zBaMtFptrLOyZx38MWpdYqy3cUYex0SPcGJ7BcV2jiFG2JVbGGOAsTDVDydqJwAtpIhrjFFrYkyq4Fip+6eFYkrOXNtKjdPMTbRApOJ7+LmmQIdSFYSUzi1HgLKRySSM+4qeQbMuhuWUpciRZz0LQifAedD5V58Ipw98cGjMLIULd2OEwJUyu99GWzxTY9sZsR7TxphHbwGWlYMF3Sgx0QPUokQpxxjoVZzdXscDu30sZRk6mcCkVtgoR9jALjblaYywg1IPA7CeYgdSoKNNOMvR5oYy25h519pUFt4Up1Ghxu5wBfcPc2hoTFBjIAbYERsYYwelHpj3NrqesYGQnNecYGOvLKJjsaaAG7NdPXVdagMdMbMkZY4MHWvYjQB9ZbSyvwIBh097M9YlxUZUgU8+ZcUAjHlACAocs8Cs8nCRV437KUVu/WphuVpXByJyAQRj0TGiRSu04CsHOtx2EqiUaR+vZQ2JAqYAjnlxJ3qAiRhgS55Errsoqi4mGGKAs6ZjrLXmTVpZc1Gto2Bbuj6pXhOpc+Pnxd0r9JNXHlW6wlCfRYkhNqUpcFPqASbVThBk6o7FrAMlAG2LejWa302hNY1S9ZQrMWKcFo4tXCHI8qXtwyJxjqq31ohSVUNhxcG/BPjc/SZWXxOoDRWSlAUy3TXN/kRuz+HKoVIvnrDYDVZjgxscscwDGPYCOvbiGkmxIKmFct4sGAegta+coSxjYNiOAQZiG5t6GUXVRVEWqEWFoRxggiHGegc1SlR6NBVkk8Sgg0sjRiMCHW1CGW0VlS0XNUZiBxIZKnsepR6i0qOA8Y2Dufn1S/U/aVy7PQDHmNFpeybMtu0MWOq7OJ6D+kAZ/ZSu/kzuVikqCKcPLn8j5AoDHLY3ArWYjtqYxyiwze/J3QE87oCPp4UK+n7Q+NwPHBcWc4FKFiC4Byyh+PiiFDMlbp5Qrlqm8+3XKgA2PPLeHH+MCXYwwGkzR+sbpHNx58roY2WPQ+dB4I3PjY/lLkcUFEvnzi62+46uW24bs5nS3iV0XaMUQwCG9ajVGFU9Cq6NzzTIArBGzE9wre19d/OPmAelKxdDw91vnMGKz4vfc36d/DZUYTaLwIufQ1wRlYqVBSyTlKi02c7dL3YemexAyQqVyiFEE9Q8NsUyGxS7kXBfzeMGmeaq2Evw4LRjTPuszY0RZ9akFnoC4FKbfkPGZSkBGCZzjB0MbFxYJgv7HtSBoRLoL6hgfvMyxPHfnAVJg6/aMaGOnbGunZG7LrVjIclYmsUmtYGbNvA5r6QqvJIbf6/SllrdqJjKgQcHcYwxNdV0SR9c3gbIFQY4CBjYSPSWpmyplxLwD2LsdoHNBiCUSdHZnIFoRjCHqU68b8e0/doUi5TUlp0WId8dlRfegrCLngK0VJAIM1p4D5Hw3G1glb3l1OG0qVx8jxmiPF0AZMz4AKFCpABcur7RSyXZwk2xMTUUtJ64c63VxFU+NfsU7j7xFtocaNA1EjoLngdS2HGmjz8338MiqdyZCyRgntAEom3CQUwQZArTpTajdDchARQwDallkC7LY1/MM59fUVbNxZPMXpOeaWXQEivDxTMCzWd5Hqs4GGvGYjUNzLSl5p+LEGCuQcGFzGK260/dkhJM+3PWL7kNuQxazim1gLbp4ZQoKMAGffMWBtPcpeH82utkxOXI265DG9vRyA5KuocS3Wrj/ejvREAxsRyecffsZ4YMlEXn1iw7nmE5OowFubz1wRUEOMiVQhVF+25BaLMW2lAk0HRnpCjAYBukQIdnGmihbD7sfjHmx+b7uxgQyliw0depOQOGHYENSIQEIHJkyCxT4XuIxJZJm2vJTI+nCkcvdqL4jR8zc4sypwA5OHDnKDJn0ZviQlXAMihVOfdIDCSJneB/x2nFVLY+YHJskTQeeJs6F85wkaT6tLhjs2ukdBUALX7O/Gc4d3+thTCZSfSv1iWEzpxzi4M0s20HWlRXjFVzcUQEYCPP+o4J28sC3mYVpxawedwd8aKccq3Rd21zjbNp9nYeZtxGafIoYD3VJ8pt2nK+e3FbcNeLHaSVeaCxqeMy/2zeWLi9iI70XpvEa0vKZWT/mPuZm+bGcYARvtaSmweki9OpdWVYczuW7yqbXfb64AoBHMyVIjso8mXkWR8AuQQiFwaaN7aN2kta9wmgEQeG8uPQAm9cH5FSmtEjoeFGoYdPV2HsBptHOH9TbIyqebq5OPrRn69E4QIYedfVWf5iZ/1ErAlfoIXIHAVY68p3xUUYKMmrufr5p0vLx5VQQ6s0lTrYtEp5dVbeejv2ycbXgIMNmguPKXHHZ11i4wWGPjeARDL6WCFMdTaMRibMOZhKtKVbkBoLgKVRDdok5Hn5KpmLIxImbsOAjUx2ZtLcsaXLPzc/p/vgY5l3Oz5uQy9FoCNk2Hzw9l5kOhPgAck0d+ms8ZIAjZf2TgCCWTU6YtapEbcwY06N4wXu1NlBw+crs1xv3ACOr5efJ0u9bolB4schl7PSXB9cvizHFQI4yALtoMiW0MlXUcg+tFao9BhQQK0nXjEn/OsuPgM21TOR7hr3ugip+WZRLqr7QC6VYMGJIpE5qIgfGDoORYnHLEUwF2RBjQGtlTufeJ94EVcobX+Y3K1PSQUcMUNu4YzmIyRPA/aLOPUPiQtscWXEi4HF94yocSrzzt0mdA2EpRnd3wmFIRhokICJy2F1UIzlMMV3b8EGlaI3dHVzO3M+Ifvi5i18518BT3UT48IpVS0UctltAKfYyqPnycxHPOaghmc3Ora6sC9hzq/5PMLTp/1n7TELcWzTtEU0BsmpcRq/i3CRDPbZY/XSNuGZW3TMZhZam2slBAENFpnR/bMkle47DWTwbae7u0Pd5XVlCDrmDcydJvNUYA3YY+v+msZ0GbdJbus4MV0V6UqnB0CuZOOQvZxZjisAcBh2g4LDMtk1DYrEEiqMIXXp8FwjWC9x86cpAl58hYMEMwtfNjgubJWiTZ1lL3mqHmNZZMSu2AcxBht8UXUxKy2Ua/A5s/Lih5x/Fovi/tKoq25QMwSm4p2peicbdKjLyGDX0FjtviBXSlG7CpFEFWrlX2Q2H/rbsSGt9zqssUIdf92c7DVqy8IJxks8N+TGiZ+5QDnYgm8SBatNEqZwC51BW2AjA/db3VCsXMFoceVVGzw/sWwnq7khZeFArquO2QIY4hTNae6DNkqdgw3+k8acOvsp7+8s4f79eE78HOYeb4rLmR8v/Cx2q/iYAm7EUZxUm6S+mxdspObUChBnAMdp29B3bS57PudzkZhx4/fXXFfZqqdpO+5Wka6e0uXLclwBgMNeWFtzI5Md5NQAjaVwxtR4/FDQX0qx4E67WLiGaSJzATpmPxagI6Tr9BqP7+YAv/gJYRu12TLLcSOvigVKunmyl5iDDep2aYo/5d5Kjuh2qhrqYgRqOPYjBcJSLgD6SdcEQOAyIkDC3SlaNFkbLQpb3t28NK4seRQAlgIO3HWlhXJlvxv3mgETKYokDUs1CrxrpzBVPe02gi003u8c3ltfWIxVMmXPCY1FtL6GSrv5mPIIwAazwLSooag087RnDLRw5TCFwK6MPPwLIcIGymU2cJxq2PAgYESgjUscPB58Fy3CMYCcBRTa/Pr82ebGR2pe4WfhGMGxYuAS3f5ziSvgcRNtC21saTcBsf05J/Dbq0wDJHwOgeuqBaTNA/wCdiEBNuaJJ0mxGG1sFrlTqFo1j6FTTF+RIeb/5VYfVJetPrjMAYdnN3zKG6vkqcvgBpjFpwheGFLyXMHzGIc4joGUlYisaooadrRtZHmaRlshxU+fu5b0wtBf1CW2wiR91swiB+DARiH7riy3ycqoXN0MHzNiLT1a2GteNj3dZbFNccTCF0HKlBEwrdTJfeReJNnc188jtBABkz1Dyp4UcyZy62qozfci/WLTfYErBezjQ/gefCHiTA2ds0Bm9/fKhQCIibuoXInz4BoyX6skIByxQzzILAaWwXXSdKyy5Rkz30tt/LYCirEcly+VeuHEvhv2frnnbopFl4rfSj3/ZtvQ/cc/T405jWFIspCiJRAZYdB6PK/UcdtcCfzZ1mj+npLwmW8yP8E2kUFkxs8CVsjcIwL17YxzfG5BymlifvQOp1iNGPDwcWedd3Mus5mfeYBGcv/Avdy83/w+UrpzshCbvRYmaFRBigpa5KhFftmynpc54ABckS9XI6NGqYamkp4auhRKKQvkNjWOFHatxs4q9y9QlH0i+LFYASdQuenM9RfJqIYEKheDAJgbLrM8YAzI8nYKRkjXKRYoXKxfhUnDMpOAbUNcOaYkFx105Aoy2D4nKFFhFFgKQVSzVqgIfESFqJzCYNfEX4Nm9dEgEIwpImUuhtuWZ5gYi181gAcvm+4P1MwicdvDsCWU9tVmaeiEQvUsB6N6VdNiNI6ywikyXosgePkZwxGwYNF8UnMgIBSwHQlRFkTXauzqoVDsEPVdIXZDawWtFGM5rnYxmWq+r6gHwXHq5DQgwJmNeReN1KLVlgbKWY5Yv/jxGog8ubintufjzOMKjOcXHja+NnXr9YmBBtchQjDDj4EPmrtnD6eDj7aYmIBBjRgsfn587NaCXy3XKY5bi90pKQZ1nrHaJNYv4ee+70qccRQ/C1rUhunTOTQ61ti8PGM5LmPAYdmNyEevtaHnhSgDpez9qaGv3y0UsVU61wx80KBXcLZhEFsgKSPDjB/65rU21TuFjQ/JrHWT2UVUQqHmAYFE7VtXgqdh/fmFL25o8RBgyGQBrRRqi35TNGi8aNJlb6TBRdaMf8lskTAtAxYJ2t46G+PBLbo4BoZLDBoce2XvaTozpRlTE4zJQZVuWrNCZIAwTAEpfQIbYVBrCFTjnifJ84mUVRwTFAeZ0XF5TRfP9vgeN7ndJZMmUFqL6rHDcjAgT0IUdOADT7AUZtv2Rb31kHug/pOgA+mFibtGaV8OTBuuCtE0BuaV2c9q080UGCUJoMFdgqSTzN8+04KDD2I9uHUfFy6Mz5PPrQE6ImkDGik2qTl21MMI6VotdP7xZ22SigGZR9oyqmhugI85JDc/MaCXcyzHZQw40LBmAE/l0+9KVY4GhQSU8pR/irJO0XCBtRQtSNyCqLWp3887vAoHIHK7VEUvpDIvopIVctEB5BJy0YVEgUJmvpZGLFEkc61LlMr0SDCfedATRnubl1y1PORtDz//XKE0FU3ty+IXQB/DYcQ81JUKI60FJJSqgiZaHLjFL6G7B8JbEbxVfMNFZsekOA8XMDjDxjPnwmJBICG0cdnEcShxUz0+VxqLA2H6GZdZD47FXD/m2QpfvxhsaK0alpxxIyhk1FvFnpMBHFczyyEsmCZ3ShxzEFuIYeBiQF2jeV/bqfU02JhmwcagIzU20eX+OHtbmNxzmuje6n5vWVjjz/h+cVA9AY0YZLjxtYLL+mI/heuSXHmjSZvAd9jwAq7fRKTn+fyD+JzEdeJAo+2+pPZrq8vTttC3SUpPJOcwxa01LdYl3C5kjyQoxTljTd5M5V3fWfry0QqXMeDw7AYpGPPgS1AAIcViCEgPQlCyhbLZQc/FUdiFuS1wx8yAA50ysECblKGpQ+Faz7Ny5AIZpK6ADMhg3DKuBbwokTGqzpxD+KJRLEIFQLL6DUFbZvYSU3wKr/UAhAwBTy11cSnETgDud7egObARUawABPLwkRZ5qLiEZTbgi5txP6xAhsyyRMoCAKixq5oaKzofH5MHlj/UGDpS9u5aRQ3waP503Xxw7XzWbBjf05I1xP3Z0eeuJgp9FtG58X7Gis/9MyUyEx+UdQzT91jIWBE2MI7YvBb3lLe0JVKLanLoiKaetnjMQ5fHi0vMGrY9G+a7sI8LP5e4YmarAeEAe3tBMz4f2mcW0CCQwY/hs9PofnjwEDDOAqBCiPMwNLPiToJto0q+6W0knejUcc8n8+RCSDqgOGShoT3YdEYMbFakyu35SVxuLMdlCjhYgzbWARKIX5CwmQ6xHRw9+7+9S8K8vKGP1weReuRIUgeLXugykSBXQGhBadSuyZegTAtRQkuzyJnmaRK1sC3W0QQIXAh0KPa3Gdu4BAh0GAXh+6Pw+JVkd1jmqkiDDoXQJwggsXD7F93UmzBshS3qheieJQp00fykABTdL5dNpBAvwkJkLq5G0jMgvDskvnZm/FAh+XgMCg4Ni30R4xIOBucKcc9KI+05ZG/4PDR7rnhNF1KulTLnAaZMyM0DVvKdvqPW9SHLcXlGqJ+PhO3nWRB2YoEhi7ntPiRdLfPQ43MAjTaJQYIrvJfQVXF2Sww6/JjtdSVSWRGz5tYGNpStCIwGmxS9ZwHTEYIPAcPouWMJr2diF9i5uiEAlo7fEguiyI0ZZX/MI42g1pbaG81YkinAMAKdDZampXAk6enY9UagQysFdRlmrFymgAPw7aZ9UyaS1AJKnzdvYEh7O0tDxGjeKgGEJa5p4W4LplKoIETWIPMDHyexHjp3C7CARAZTyrrSpncKHzelAGM/q2AvNQR1jGzGAQSXtY26i4BBE2jE1gwflOYtIZxbhP+roYRycRJ+++lC84mDXt39hAwaoHFrqk2JEDDQunYKIwAy1kXDFVZNCjc1FnxQpzluWOit7Zx4Vgt/tvluQTyLrqDQVGaUMi1VfhWzHKb9PPVQymQ3qElzYY4wuwjcNEmxGVO3Y1Y/Bx1ujMSCmXqX9xq3kppL7I4KgAZjNTxj3AT19LlnSy0LpeE6p2tY9pJiqdg6GLp65gMbjbioyLWSYsAo4yNmXsxxZ2c6hYUOEwAxej/53GgtSDGtqe3t5IM5xucf4whikRUqE8sh5EXQB6zwyh7lMgQcglkztu6G7Zsy60HkvkcgBBipgKIUiCAfmS/QZBfMqEto+JJ6y7XhT44W2VqbVtAlClRijEqPQYGoKRbFX5V2ujBYRNkc+L/ktbJWDDEh0DGgYmBj2rWn74R9tqEM22BjbNxmNk5CRoGk1MgM8ExSpcbNmilmY3fNTcdZy6g07kmbrzQq1sQsHf6MUA8WY41a95JzZzQLQvHGfwKmIA93GzXmxF05gKswKoVEJrrgMTpUs4WuldmdBeJalqPWE1x9dTlIHxiwkcs+sqzT+j7Mu1jF+3LDwG8z3b3WBlBaP2eWqDsmy0aYtk88LmdwuLtlz3PlRlekE8xnFdNp7eyGH8/Oj9hWcq/YTDQoQEog6DHFFtVUdlEqDZWzk4hcoQ2dD369jb4wLkpf6NExSNEx43HjuJIAuNB2jO301yVMDGhet7R7jcenpYzRNkOIyp0L0QH06AJqAipVAZwLiLkMAYd01owpW9y1BZXSFBmXeIEj65N+ByJFMyUdk6Th14wAh1t8WMCj2zexANZqDAmJsQSEziwTUSbjFfj+VI48fomEyFx8BFi1z5T1EroTap8yDKtYOJPh9g3Bhk48usKZCOQKIEVl40ns7kqwNE/YNE87HI+R0VAObFAsjDsPp8hsKqTyFg3vdZLykTctH1/mPr5GtPjX7HvuMmnzvRvXVu66vcaWCX9+SMHlAGotkYsMOXqAABRqKFGajCztr4FTeFoFfX2kzJGpDstYuVoYDul6KGWygyzruEJfe/HxA+lnIv4uuV9iEeQSBxzyz/n77D9XbnF1+87oudSwitvYjJbPw/k29dL8YGMK08nOIwQeAFCBYruo+STv7uwCqBMF29pAlztk9I4DxBx69hOAq5OknJ6Cf6+Ed0XH4CN2VdE8mqCCXK2JHlVR/Z82oYwrMmJoLN4x1jOzUxgZG58oXV2OGhfKACF9fy5GzWUGOHxHWKJOqVKneQA8aODNhpyLgSkT51ePXr44LdNE+VZuTDN+5h4+Ol78YprZ2ohh6UFAEJjJ4km4VWMKSIVzT6XTxorB0MjZTPAVp3C6udifvGooDwYNX4aEYp3z4XIKChWU9pkuwioAqStk0rRUJgRPgItnaHAlGI9PBbAUFDIwJQNAMxcFFUrjysewB2bs2gan0j0gJsvcy0Rpa1snxd0r+PifafZwvBA4V45WKK0ihARg43ukOyeFTER9eiIhRsY/7/RP4cpmOWzPFNFBJjq+0i5vQ89T4NkiMM2l5n8Pt019Nw1spCuDppmGFKAJ/fDN/QKAm6LiI+u8bZ7h/C4s2CCdEBgd5gM7tgcdAspR/SbezsZNaV88kOvwlEUfC78/ZHjlsoscPcP82XcpXvyN7jAxZjVKW0SyDPtZJYy2afOZKxB2jngOYsyp6FfAygY6vb1mCYENJWwWp5a4UG4Vvw5cDS4Vm9JDrhRe8IiCpFyaJcxNMg8KS79C2NXVARIW1MStXBozaA8PCaV9kB5/4T27QQjdZ7t4ZBmWQOcPo9K+eRzfhsbmiy2NyedB47kYByjAdYdNuxUaFLJlIQh0JN0mQjLGZK8PFzEVE0D4NDkhJISS0JmClBWUtVbbGI2UReqvUW0XI9Nfh9JuabHOUKBAF4XuuIJpyiqYUo4xwRATMUCpB6jV2BVhq3UdgCGaA4Fcxa6tef7Csus0tzg+hH/vzoPR6dTADXLJKUpJz5f2RelaaVRRmHQ4XcEEi135Qil+1BGWKveSBGmhLNh5FigPa1qE2yYt1CkSu0rC4/h71rbY8IUmrkuRCkRuY2liN0BrTZr4p2U7DXg2rtVp8VspXcA/ExCI3axCwHZNzqFQOT1tnmvlaw4ljMTUNUsFagqRIZdddMQKelhBV/fR0T3kyJEhoyO6fRQ0SlSoUGEiRhiLISZiiEqPMdEDxz5Tfx5zntPB3DT31rRKqrR/wK5oGB3MRIL6VqWf89jIlSKHgoSAOCdGIjFLnA9wuYwAh2/KlNmW07m1akj0HAqAwEYuOgDC6ntkEVPwIndP8C6dtTbBoLRwuLEtdU1D0g2VIkcuuh6QiNpkoyjvo/fWtV9gkoFmU3pA8Hl4pVFbO7YZQ5KikAm9e3o+Ru/MYkkqGNpeBvvEbhWK5fBuHgbwIFGrsZlnZtwDbRb8tLiVQvbRzdaxLo5hVe3Dsu6hQIbM7tOREr0sQ0cKdDMzv1JpTGqN7arCthpjQ2xhW57GLk5DKYVKeTZBSzYf+/Ib0OGvJc2Zgoc93RkuQj6VM+y14V1tpsy9lNJkK7W0puZC9Cu/Jg5ku2CxK5nlsPrAGiCpnkSc6WzE1gT1G+ZNd26WDk+DhOlKl9iwOLhwmhUcu/xcnAZvjjZjjHhu06xpM6QHG1obo4P0wrxgozF2iz4g0EEsHM9gofOVKNoHjubPY64ykaOQffTEOvbpo9inV7Equ1guMnQygUIKZAKwasAaDVYfKI1BpTCsa2zqAbbEJnblBsZ6B6UaADHoSDJbkauUuUZSkorx8MaKB2vUT2Xa89vKdIAM6cun3PllBDjMIidlB0W+jCJbRidbgYTtjWL9V5z6A2jBZb5t4VMFzQtUomb+Lg3jlzf0m3/BBTJkKKBciqih2zkzYgIJWXlua1XmootCLkGjhtS+uZoQsnF/uQKJFyFiUZQqzQKmmxQrd42kXC/cJ2radIfBrlSdVcc1SrQKlElascQulxB00D30Soaq3bGqmhZ81HZ3oXxNBfqnNWOTGNvh3GW2t0wv249DuBE3iSO4ZrXA/o7AagHXbbUjNZYzjaVMoZeZMUa1xG4tcXLcxclRFw/tLuGRahXHZYFNPBiUw6fS8nFlSzOvZjyPUiWkLFzND+fTZYybvUj+fnPQoitUeuLueyYK31OFgxj4hTBeUKTMITQ1cVJXLNQgdwp1iKaOsCnFy0FHY5QWRT3d/92+oLdR2Px4fAHlgDA2NDhbSW4xXtUYALSwfXXg441S+6fnOt0qj8EGd6M4sNEKNFJjztIH3r3CDREylgz7UQaxHfG5cP1HP0nf98V+HFDHcH12AEeXchzuCawXGkuZRi9TKIR2gEPDpN9PlMBYCWxVBTYmBR4edHBytIRHdA9nxQlAAhO144AYuaunSQw6uLSyZyJ8trwhA7iyB4h3aQax+uNwnUqtQei+XjqtcBkBDmkurLVocmH8cEJI+8L5ACLuX6PPpglX1pT6atC1coGLWtRQ9gZ6cKMcQGj4jq1kojDda1EAsM3lpGVWFKBlaO1yK5RbNBRbUuvSud+DwlDsuCmrLpUiJiAB6RcmipFQgX82JU0Woy0qvXmtraKJlAzArE8bWEpl0cmqSYEpsCwSwMRpZLKDXrYPR3Azbs4O45mHCnz+2hjX9odY749QK4lxlUMIjUIq5FmNPCNQIDGpMpwd9XBy1MWndzq4e3sNxU5m7lWmMNRnGuxFOK/p8TPmZwg2qOYGBX1BwPmrm1R37ZgxrVVQ2bZ5vSOwqTPbzMk2u7tiWQ4DVHlsSuxuIGnztbfJuRZ2mslqBG6aKMaEBQ023lORuT5QhsXxVn5tCxtWGINixqYdPxVoGZ5DCORTYGM6qzHtGk9hPxuggywOmNRZZjTGzH8MPtzn8AZfT67jgDqGx2X7cctagRuWNa7plTjYmWApr9DJFKTQyAS90wK1kqiVwERl2C5znJkU2NfJsTboodg5CFFLnJGmOCPFVbQ1CUxlELXNmc4pGiB4Rugnz5pMAetZgdNSFFCiisqdXzp9cJkADl/oK7PBornsBjdFaYVaj11TNqfQE5HCLs0TlE7ELEQN8EWffPW1Ng+/0jxwSgXHMd1aWT0QW8QqRw+5bezmFKJQBnjo0L/PFyCqksmDGis9MvEE0pdn57EpJNICNA1pmQyeoREGSSmw9vJqEiiW6Up62ncpq8orHA1tQQe9fPweAdSBNcX42D981pGLeJfIZAfdbA0H5PV4Yn4Uzzwk8XU3HseNXzZE9pRrgMPHgM0d6ONnobfH0IMKelRDlxqQgOxmQC4AtYPqdIXP/8Qq/v7UASyfWoHcuBn3ZCbGZoxtd83oOaF7we+Lfz4Y6xCDALB+PJaejZ9Df218zRJ63uOidC6oOSrfTy6/DB3zuVC4UmM5TKnmjnNZZrbCKOCVeww05gnq9N+1u+9SmR7x2G3MCTcgABZjEi0ojq0jg8XGIXXEkos3AoBSjFDqIaCAkvZrLMbNeBRuiKTOPfzXHiDqwcZ8YK65rUyCDjNXxnZYAyTIVknUJ3HnTK5zacDGfn0tbikO4pa1HLetV7h1dQdH13ewsjZG1lGQdljrNYLWwv4EdCVRlhK7u11ct7mKe3t99LMOsq0DkJWElgoDAErtNEDfrCypdDAwMTr2GWGMVQxoHEuiEZRqCMZjoIN/J6UJzs1kxz5z+SXXB5cJ4ABckzZk7kZobUL8Kj22qaPMQtfKPZxBKpoFGFS1k7blxaOSwUmMUSCA0ghMEtJVtTRj0661a9AGsIVTh7QqsSS5bTVPQCVH153zRAwxEptAvYNa+KBVyoaQKBza5QFEUhZOwfJgWVexk+o6QDXdKQ1XySyZBjbodxn8Hrp9/JYc1JnvfLojj1gHzAuUZ30sZQdxRB3FLasZblsb4rqnbSP7sqdCPeEW6PV1iFOnIPc9BHFmE9jchd4eQW9PjEG13AGWuxCrfWTDCW7sPILO39fYro5ht1rCeOd6jLMdmFDk3eC5oWA3LjFrReeUkuSiZd02NBZJ/Pz6Y2XumRJsQYldP7SQ+MI/VxLL4eO5DKuYKvXcbDY2r6T28W6q6fUs9ipJq5Tdx0L2zaIp1tHDCnp6GYWmtgcVhmKAIbGmvFvxnGxOSgJ2gwN/DjYCF+u5H8u9/zyug+kEnsEC4YMl24AGwECdjdtYFYdxVB3Adcs5bl5WuHV1BzcdO4vlYyXyAxlEbtLUdaWBWpufSsPaQ9Cqhp5UWB5M0HukRP/MGkq1iloXqLfWMVZjKGnq/hDoO5fnbpqrjjOdcSwSAHddZh6D9D65qWzMYSZJH+SXVB9cJoCDlEvHXCAAlRpbKkuZ4ljMJcAtGhftLKTrusobvBHjYR6QsHgUjUH0tYJyCz8F28T+9hphmXMAUKJ0ed8AHAvDXRcENjrZMnpyHcvYj2W1hmXdRxc5ciGhAWzqAU7LE9jKjqNUQ1Rq5M9XGqtVygIapiImBx2xv9Ccf5giG/9sWjKz71VTZr14KarVWt/C1M8QWjrAl1n6jxZyrU1KcC776Gf7cUg/Djf0VnDbeoXbrj2B4sn7oPetA3UFceoUxKnTwJlNYHsIPZxAD0uoHdsrZVRDVgpY6gCrfRS3Hcbjls/iuX9+CrU+jFqvYbh7M5AB2/o4Sr0bgCLOaPCYCt4xmNxJdJ3JDx9ed7h9+X3jWQkunoddC7ONdEpZwj+jmeyyKqnE9E1wJcZyxACKZ2KF3Xybza2mSWCdMgVPf6fiBvjYs9wqtDDETdXi+VGgeU+uY0nsw351GGvoYzkr0MnMHEqlsFVNcFr0sSkkallCqVlN6NPz5IwGzVNbcz9wo1xEcawnYJ9ZBjqsPoB9non5DOJc2O+Z7KIjV7AijuCYOooblvt44prCU9Z3cdOxs1h7fIXscM8YGLZfgp7UQK2ASkNXzfPNRjX290fo9StIodGRq8hEB3rjCKSWqKSJ6SuVajwb851/mMacCguIGSrz3ezOsVyPxLVgyK3uMtiSzDan4i6eXCaAA+D+WlIuSvvYg9DfqJoKGN5iBEJXBuBRfRjgFS+AcOwBT72jbU2juCzorWJ2ngQBXxre707HzGQX3WwNq/II9qnDOCLWsa9bYK0jsZQJ5JYpPTHsoDfsmGAluYNaTfz5q9otyq6gmSs2JZ0iBtiCxd1DbUolAh/T7tGehcwINj7PYHF1L3TYPZXYLmHBBlkzy+Igjuj9uGZJ4sblHaw/bgyxehC6riBOn4UYjoATZ6BPbQKDErqsjVtlrKCV0T16XEFMaqALYP8K5A0K1958Cp8/6OPkeA0nR/sxUgOMs21/7SMmhubo58ro+MhCoUyilBB44DR+w7rhANvl5Nvrw4LQTKaVsakzaZq6USG9K6sQWBi/AVjDgYM4BjZCy7Bm1yhceBuLbSRtwZfzSFxinLsBOTClzwhsrOMo9qsDuKZYxoFehrVCoGeSnTCuM5waZSiGGaQWKMUQShgrO+jdtAdWptGQrbHoXUh2g4QznunPvWtV+fch8X4ZIEpxG2tY0/txtNPHtUsCNy6NcN2+Law8rkR2zRLEeh+iWwBKQVcKolMDkxqolKenrQgpoJcUil4F2S1xfbUBpQV26zVslQVG2/uwLddNdWiEPY3ORVLPX5MBbV6zOH07LvfQEMuESplDKMliOS6NXAaAw8RvSJGbi2IRoBJEJzer9XGfubMYYNIVHR3PL74ABDL3k3zdvGIh96v56m7MmtLKZruYWg1xnxIhJJQtTASESk2K3MUd3KSvx3XLXdy8KnBtr8aBzgSrRYnl3DTZuWt7GX+3sYLu2Zuwk53GpN5FpYeM1jdsgNASSijkKNjczfe8iFZdTwIqllvpe0Ppbcq4fYzAd8vFVSM1VKpCBanJ39i1R7Npc3ZBLqRxpRxWx3DDShc3LSscWRkgWxbQu2OIh04YJTIcQ5/cRn1yaPKb6ZAdCVlIiJUcol9ATyqI7aH5sq7RPZbh1rOncXrcxdlJF5Oz12BXbqCWY5eyGrvi4sBfe9I2tiftV21co2kLnw63c7E6EkFMD7EarpGdApT0GTY+WOzKAR1hXIz3YwPebTrtus4Kpms/bhbc31b3GFmdDSbBl6CnNOn4vHLRRT/bj336KK7DUTxuuYsnrAlc26twqFtiKTMAdVDluHfQwV3bPfS2D2GoB4bl0CVS7E58HnYiwXPIWQ7ObpyPi2Yv0tAJUWA5tIIWTH/ygHGb6ZVbdmNdH8ZR7Mf1KxluWa5xw+oO9h8dIL+mC7FvCWKlC2QZoDREWQGVhJYSUAqCAIf0BxBKQxQZRC6xWo1xkziLYZVjVC9jVPVwenQQY7GDSgxMI84UWznz/Nma0fKMTkvnDgODo+Z9bUGqEtDKgI5a+0zAUB88OhzoZQA4rAUmWOtyNMEGFx40JMAqhMY/meLmNDQpYp7XH0eOBy4YhEouZA18m3u6uYBXPM4ylwdxjboGn7evhyeuKTxpdRfXrOxidXmE3kqJ/iEF2RM49ok+MnEUu2UX9+wewa48CSFCWr9xBUWYu66gTLdAVbnUXiro4675eVsu57G/vadEpQrNlWENIclqV6A27rlcwqo+iENyBdctC1zXn2B5aWzIk50xdLZtAccE6uwI9WYNXWsIAciegOhnEL3MUKydHKhr6N0KqBT0pIZcLbB64xiPP7OFB4YHcXLYwYOjAxjIs5BiF7WeBKdAlhYPEATgeqlwPyow2+fenq/PitdZBcxz/On5Nc+zDJ5FKQqEhYMuROGfiy2i8T5qtjiZv8P3Mnx3oyJ5LN0w5ScPrOgIbMT7xNJWWMvv59MiCZRmokAhl7CM/TikD+Km1S5uWQWeujbEDWvb2L8+QHepgtbAeLfAgZP7UMhVlLrAic2DGMkdlGIAoavg3JJXMgE66JrF7IYR1TQOLpgwNiOI5/DfxW4UIGSv6HnPRBd9sY51tQ8Hex1c29e4bmmMA/t20TkEyPWuARtd606p/HmZbFxpdIUUEAxwaKWBPIOQArnSWFFjXLezhTOTDo6POrh3uIJNuYKhOGueFbY2zQM2YpkGiPcaRxQHopKoAODQ833uTN75yiUGHHSzwwBOxQLmSPiDJ7mSZYyEgGr0TgHgymBzxiGTJio8TkMzjcTSbeJjt47r3cEqv7lqgaym/1J2CI9TN+Hz1pbxgqMjPOXwaew/MkCxbpR/1hfIP+8QxLUH8bjrHsTtf/YgtqsbcMf2MWzkD6Ksdxkt66ub8ofLtWiHp+hdzAt8oOj0hW9eIHLhaFbHtsSgEkY5E+jryhUcVAdx7UqBW1dK3LS6jW6/gpoAareEzCWQS2+xKA09AZQCdK1N0FReQ1QKkIZa1eMKemD7sSxlyA4WOHZ0C5836OPEeAV3D/djU65hJDcaLjqqiRI8O1HQLpc2sBgvkHGGA3UvVqybMGer+Pg1SgBlYq4S+nLHGVaocaOw2VskqdiLNoXNg6aBpmvFbceCEPl+5yIpN1g8lrRuop5YwyF1BDf0l/DUfRpPWRvgideewvr1E+RHO5CrSwDMc7386RPofFZhVO/Dw7srGNSHMZJbULpEFYFgfv6Ad8H568XiuS4YsEg9WC0BLMFePHOFZ7JJZ5BohAwHgfyOWMKy3odD2RKuWZJ4XH+Ca5Z3sby/RLaeA/0OUBTGsCBR9llSGiZDzopkYBXK6M9eAak08kph36kBjm0PcahbYC3roqeWkYmuc6vEhurM856z+q2fU6RHWp55rVUjDIC72Pcyx4spl5zhoHTYkEJtKpdgnwhsuIhcNPutAKZPh4vst8xIZnPeqSQ2SaXHAEobTNpONfK6IDSnwAK1n5kI9DUczpbw+FXgi244joP/QEAuL0PXCurMGKIjIa47BPXUz4M4vB+H1Mfx1Ee2caSzgnur/diVJ13Jb5cRwQEVagvYzQJYI8xuSTVmCy+otGhdYzromKakooIDrX5btgdTyjGjZEbyxWv6WMOBrI9r+gJPWN3B467ZRGe1NkbZqIbulC5ATHQyyF4FXWuoMaArQI00RK6gx4bVwLiGHtWot83inHUkxFKO/rESN5zdwjXbS9iX9dDTK8hFB7UYN1gm6t1AQDOwJloslJiNMtuGjZrcthrOrSSR+2J3zKqnSoQatbsFQSE8a9GY0saXuwjzLEaF1niwHf8s2DNiKOiz2Lo3Rwkbc9G+JAFjNQe4TgGgtqyXTBToYQX7ZR/XLEk8cXWAW4+dxv7Pq5DftA/i0BqwugTkGbLhGEv9h3FzeRonhj3c0e/h9PYqzsh+g+Xgx+HB9I1UYcZunL87pe2Jaolg36NO4EIGlhQFCtHHslrC/m6Gwz2Na/pj7F8boNgHiKUcopMb44OLJFbDzkElFn5pajNpJYBeDrHaQW//Lg6cGGJ/sYLlPENvvGRi9aLA8TZpBHDOsR0whfGMWDsuPOmBy7nWnbkYcskBB6XD8gCxILVHhKCDXiQONrylkoFalcc5zpQ2GhYPCis40mdKV0FgaNDcyWUHFKDsFi6+dK+XAl0s5xkOdRT2Pb5E9nk3AsMxcHIbutTQ4wpZWQH9PvT110M+bRs3vf8zOHLvOta3DmNT3o9KjAKXEA+u9edtMmUkq+9hlI65ErNpfcFAh9mHndnUfc9daHFUrKFcbeN5fMXNnl7BWjfDtX2Fm687g9XbMkBJ6IlV9qMaoqiNolg2ForoVBA7CrrW0LVGvaMATMwzVWqokYYam6sjl2uIXED2BFZXxlgvaqNgJivGdQLpmQOe1gvlnhUuKWVPJZhjIQBN50s+evfsCuNWEppnyfhsljgktfHMXlFCGWvN6q4pSSnmvewXsxxACDLi4N+9ShDkKzJIFCh0F/0iw3oHONIfYO3aMfLr1iCO7QcOrEMvG8CB8QRyd4il40McuW+Eg70+Vne76GAJI7EFYUFwPL9URsNsUNH2/bm6X9tAxzx7WiCE0AgxBdK6KNDFErpY60gcKGrs742wvH8CuZpB9AuXleLrFjBgRN+RKGXBiNUj9F2WGcOlL7Dcn6CfafQyiQ68G3/ep+FCAI1Z+5HEz24qQPpS6oRLDDis9RVUE2zGX8z73NKDoJiy5hHiIaAxarrCKGhbXOsxKjV27pLYL+usTVQBwDd0pQUCdmHiWRdSAL1MITvYAQ6uAae3ALkDPdFQpUYxmgB1Bb22DnXzjVi//v/gyN8D+zZX8ZBcQiZ3fXpktNhoPYHJ6Qek6DrqUcrcxQXNVDjCUJkedND9IdkLuzGvMIvHBkVqTR1g4ZinXHTRV32sdSSOdsdYeZJAdtvjgM1dqNO70NsT6ImCLmuIXg4sdyGXOsBgAnF2BLVbo95VUGNtslYUoEpAVwJaCcgc0MMaOheABDrLFVbyGr0sQ1d3G6wELzwltY83cpcy0SfBXB1fwpqLQBZkOJmS0x7E0LuhLQCnOfhgwGb590Y8wR7vzKWRuBJvu4tjHiq7LSA3jN1o95WHWS3nbiVyl41EhhwZCmlK7/c7FfJ1CbG+BKyvQK+uAMtLZhHMc2BtBdn+Dtb7I6wVGst5hm7ZD+oB8fm62JErUNrYDcBeO5FDIkOhu+jJHKsFsFbUWFsaoVjTkGsF0Cs8u1HZ66C0AxQAXOyGAxdKmd9jYJJLyL5Et1ehnyl0sgwF8oARn31O6Tif+DxnAec2/R0D55hta2RyuWfxQjVz25tccobDuwgSfUOI4Ui88KTs+TVLofo4qIduTG07agpROheEggmudI3Fonm64FM374yN36x4SvUXJhhiu6xxYpxj/Lkxlm46CRQ5sNpDdniCrFJAr+NfkG4X2VqGfYU2fkO9hkm2g6oeuutF50t9O0IXiwoVtnum6JwSL7V7UM/lAWwHGs0AMS7TkTYpaNP1NUcvA5byGnJfB1hfMf7Y4QR6t/QKpcghVvtG8QzHkJ0M4uwIUGPUuxr1GFATgXoioZWAkKYIjotLLCSyboVcaGQiPe/YGm60tU8ozqApmK0sm75etQNddF+du8VSqT5zKszgCu812PGkBZQal3+mii8ASNIWW9EezMlcdQmgPA1spNw2qYyYaVbiNCtVoUaFGpNaY1QLDCc56l2FvKrNW5JnPq6gqoCyhJ7UqFUHAoAQ5MrLEDcSjOc/7RpcjpJahHm/K8DG3+kCmTDN2AqpkeUKoiOATNp/WRDL5X7yeygly9YPwQgAiDyDVhqiI5F1anSkQiaMc9LMNQviY/w5pOOFZmVUpc4/3pe7ysJjNsFGW2E8t34IiUuRKn8JAYewC1HTkiGwEXeGrCkAEgrQQA1TAjrtt/KFkwCfUuXSSyMa3FmtyhcJC5SesJ1iE+uxYz5sWi+5gbRWqNUYQ72Jk9UAd++s4+47D+DW9eMobjsEsX8FWSc3QKPfhdjdNSlbGxvQE4WO1OhlGXrVCoay786LfvKW7nQOjah9eIvR1b+IFdB5gY3zFfYCMDEK1QZmokAOiVwAhVAQpFQ6uQEYuSk/jEwa/+2qsRYxnkAUBZBvQ1YKWpWoxxp1JVBNJOpaoujUyCUgcgHRy4BKQ3Yq5JIcUZ7WDSxj+Jecpxq6+gGsDoM7J6c0JXL07Pimb4px79WulLnL2IkDG4WE0mN37Dg9G3GxJEHXVsKAjcs7U4UD6ibr6avQtkkq3iN9nMhAaBgrs3uSTLNSUzEl9H5OxBA7VYWzkxwnhku49pEtdM4OII6OgPHYgI5JCbGxCX1yC9XJEjvjfZgoAVcWPKXzohiSBvjQ6f0ujpybO8Xv3VyESR9bCO6OwAOihRSWHcos2FBGv2oVulkkn5+N6aCsFSmBXEIgh+hlyLsVMqEhrQGSymicByzwc3Ou0eg848q00wKkUwzItJTxAGyA64RHTy4bhsP9HvUFcApH5BA6C4AHQMCiebNVVJyLp6il+lDEkdxmAvTDuzDiBYTSD6muhxSFOTZrvFbqATbFBh4ZruBTm+vYf8cQxw5vG0t8uefqR4g7PwsxGEM/cAYP//0Sjo8kdqsxapCLRro2ya4Rm/IZCfS7qwXizkMCmnyxaTfJxQEb5niO5QheShksiESZeqvMlIvPdIFalFDQqDQwrDOo7Qmy3aFRJDJSbPHfygARudqxH5SQtpFbVmnkvRrFmkZ2tA95aAV6OEF2doJCqMYV8cCSpcMKCSXKQMHTc8djkZoFu6SNC/HuPIof8o35mj0Y/DGazzBtw3838ScTu+BcCYGjDCgFn0VdWM0vgTXJKzLOAiaxTG05HrGk5yO1LjEWA2zrMU6POrhv0MWRR9aw9Lkz6CydhJhUQL8L7I6gz2yj+tQZnLxrCfftLOPUCNiuKpRiEsVuNFNcG6ynVlaP8fbw58t8pMDrHoFGYqFOpSZzUVCotdEHEyVRVZmJxaptgS+lmvEbBDZsOmxDT5BYsIEiN9lsuQSkhhTtb07MmDUKzs24Nxxs0M8AdMwAG3y7aWADgGvuqJw+eHQNkEsPOBIMBwlZeOQzE8hcoKZjO3SNOLshRHvkb698WV2Elim3YGMfW8NVwmgtM7euLVpWuI6PSocLvtIVduQGjo8O4eOby1jODmPtjvuxvLoBcXgVEBL6gdMo79nB1udy3PPIAfzt6TV8cqPGw/oMhtgMUh3jEu80X1qQeAqlEKZImDkf3qIYSAeEtbldUtueq9B8fMaN68Zrj6t0ZcrIC1NEa4gxtkuNjTJH9cgE2ZktowwCylQbEFKWwHAE7A6B3ZH5bqUL2cshlibItkrI0yVUqZCvCmSHe5A3HwGOHYQ4cQby4V1kluEwDdDCa0wpuzTXTCpTbAtNEGsAFQFe2/unBSDzLsixC49YPV75lI7BS/b7oFK/vxQ56ktIo84ntoFjzM5FTEFMHzeriYZloFu77M4qUZ4IAm6deWNOiYJwIHdYiVIPsSk2cHzYxWe2e+jKfejfUeHYcAOdY9sQSznUZonxcY2H71/D/zmzH5/c6uCB3Qqn1A5Gcqc1dT84h5aFyoDjEMCe+/t9HnFbqdHIAGl1OZr071IrjGpgUEkMxwXqXQE9qCCWS6CsDMMBmHe/rkN2A2i6XEgIbHRyoJKGTQWgtECtDEdlzprYgrbzSFe5decYSRtLP4+ktktlTEFw1pC6ST9mGA6rAFmsAQXnEfVIJVhr7X26zkJMRN/GLpKk8k9sy//mqDqVwhgGnFlESt1fpUmxrRQgRQXfYrzGWO/grNjCPdt9rOQ9HLjzCJ6wexr9A1sAgM2Hu3jg9BHctb2Mu3Zy3LWlcPfkDM7I45jUO0Evjrihl3PlROcdBh3mACqE7aEvA+HuKKa4tXM1SAiMMRID7JQKG2WG0WmJztldE2gHQFfKdIOlwNFRaSyVkSltDqVMieNuYVwwxRhZpSBrbcDGtfugb74O+sB+yN0hdKkwriVKpmCkLEy9DiFd51J3/7WEts9KCnRImQeKyd0j95yXtjlhuuSwyUZSzvKZZZnSs0ssm9IlhKLqgo8+jTq3sPfK/2vGRs2SvTZgm+l+aVkgU4znNFHWBVxhhB2xgTNqBffvFihkjl52AJvDLvbdO0KR19gd9XFq0Me9u318djfHPdsKj5QDbMkNTPTAgdOgE/aUeXBwQbS6cGnXMRv66OsGMkDizwKdbdnnWlQYqwqDSmOnlhiWBcqRQHesTOnySWVcrnlmGE5hWAonCXbDFQCj/fLMFAWsFeqxxERJ1BpOH7g5s7Tr2H02jzhDQrSAjhZ2pLGeMXZjVn0aKXMInUOICs2KoxdXLoMsFYCKG1HKILEXLj1QGWXLAUDQnp1R0DGlxF0lQnvQwoXKnvO/yUVi9o1TZ5WzygGwPiqZPZcCmfDlh2tdodQDnJHH8blRH+NTKzg7WceR02voSY1aC5yeCJwYapwY1ni42sAZeRIb4kFMqp2gwA9vxkbWKzXt4kGvvG4HZcuYa1XZbJdpivYi+Xm1SlKogL/mdG/MuYztYq4wyLZxZlzi5LiHzdN9LD8yQNYrDHqqNNRAQRQlsv4YKOxiMy5Dt0suITpdE30+tLE91x+EvuVxUI+/BSg6wF2fQ3myxplJgVGtbRdgm1ki6aXN3JyBsLYG3RceB8QBoNIVKoxskS7Ap2FXLjOKB5i62hu6bqgFUh46SJc183PvkpT2fpOCMSNffnEcZIUxl1MUA0EyG3ClwcZeQPYs8DAtSHPa9jVKA1xlhtPyEWSjDLvlMjYmHRzudbCca2QCGNYCmxPg1EjjkeEEp9QOTssTGGILJQGOGT0xYrcSd19Tpp1nOygW4mKDjdQ1ilgtyHDxZu9OLUqUGGOECbYmS9iYSGyNOxjtFFgaTCDLGihr0zeF4jGgAJk3A0cb05A2O8U+P1UNvVVisFlgu8owqDRKVM5dKjVrnIkm07VXSaW0mt/DMufunrYkS/D7nkr/vpQVRy8jl0oYoBkoaZSGRWgJEOVphhrSMAu8hwC5TBLdYv08MrdAuHgMW7VUqRIQCOIlKDhvnu6UCqZ77Ehs4oR8COPyILbPrKOfZZAQqLXGZj3GptjGpjiNXXEKk2oHpRo2Hl7OYNACnVl3BHSZVBcEOpRjOKYF1Z3rQ3iOlKxWDVqSAzyKz5nA+r3HPZzYXsHh47vIDk+A1Z6fsoZLbwuGnNT+815hfpc2qn11CVhZMVHr4xHw0GlsPNTD6UmOUW3iKaigWptw1o0HjJop+b4qvkAYIFquVeAqYAqHK6AgloHijJJWbWa89haYqkugYPYmMjy3SGa5SYAm2Jinw+u8cq5jNeMtapRqgIGUOCUlxvV+DLbW8NBugU4mIAUwqTUGdY1tNcaG2MJAbmGgz6LSYxfnE3eCDubKrqMLZIbx4UPAuFmnxHKE6fGXRoJKs854K1GKEYZihJ2yxlYpcbYssDPoYn13jGxUQ1TMhSIRulc4RosDSUmUMmzp9hDV6RIb2/uxWUqMa4UKTZCXfPcgG67+tm2BJlvmP29v0ncuMTiO5aA4jkdZLhHgsP5aTqFTfwjb/ZQWGlp8VF25oEIere786LQgCOOCqdXY097ag41mTIb/3FtW3lcPGDqdsl1ipsP8lMwFEAb+0T6VngD1DnYyYCx3sIElZCpHjQolRijFEKUaoFRD1/8kjmQmV0pQI2TeGiU2lqONMj0/BSPZz1RMCAsctZ85q4rur6X9IQGw4ml0Lyo9xpbYxonhCu7aXsbR+5Zx7Q1jyNWeqywqCgmR20yVJdMEDrsjU110WwH5GOgXMGXPFUTfzrssIT9zF8TDj2D4l8fxyePX4f6BxLCuoKSC1AWk8O4xelZ5JhOlrypbZ4PHfpj9TFZRLaXZQ/iX38SEeGqcgw7aN2DZon15TQ6aoxIKme0mK5kL6HKtOMr1QeM7EQbLzsMqcBdjquU3fce3bZ3btFodUQBgvH3KYKgJMCgFJUsM5Ra2sYZu3UdWGz1Wo0QpJhjLASqMUeoBKqsXasZuNIvLJdImGejwsT0UWMtjOizLYRqOsBEuZPwWzakt+KEZ0+fABioIPUaphxiILWzUqzgzznF8VODQ7jIOnh0g319BrlZAnwI+c890UPqrUmFch+0mC6khxiYGRG8PoR7eweb9Hdy3vYLjI4mdaoLS6gHq1Cx1yEwA/nnlabOB6y2K7dmLxIA3BsE8UD2ZtRitb4+2PrjkDAeQcHEIXzHUKGr7cgnvRnD7kPJlF5UUv2D0No8BoXFjIcagEUHMKnwG/nnhbyxV/YSm2gzevwoYlqQUpobCBDsYWKVUqbFvsKZ8sGdwXVx2jV+Q+DWjY9es4Z2bP2XYTPPTNUDI9MCu6TJDOSXcKnyhpWvgI7dzC9hGGMgtnJ0cwv3DHh7YWMPR08chj1YQ/RzZQUD0c4j1PrB/BVhZNoNv7Bhf7KgGcsNqCClMGqwUQF1DnN0Ejp9G/X+O475P78OdO32cGGrsYoIaJSSBUVBsifKgWHhAoKCcMpka9KYrBzokbEC0yFxwL1Wr5QtFeE8T8Q2aMUNuTr6vz7kquEdHWE+lPc6zTaHz36f5tGfObAbYmDo3puwDQMied61rVGKMsdhhwfGsXD7FM9l3m7dbaAS2tsRu2IOGTAfkHlmOCwE6ZgG76S4pUxjQsBwTMcQ2Bjgz7uHkuMAjoy6u2eyhd3YHcr2EXCoh8gwodCJzzQKNSkFPKoOybA0kXVbApIY6OcTo/hIPnTqM+wZdnB5rDHWJUvoaTcQwEuiIs0vceSVcHOZ8m/FJqbVwGqt2btVvebr8oyuXDnDYIKF44eSxPXzxNwFDJg6DAw/znXSLv9sXMgALJBRg6V966dwTcQVIUhiZBThaK0A2KVJCsooDg4jhoCqhbjGx41Fqa1DCN9G0yimoREaOYVTg6nLQPtTNFLCxU5HV3S5tD+K8D+iU7VIPuVUmSpee2bDxNpBmUVa6wkhv4bText3bXRzuLOHxdxc4eHgAeXgJ8uiaSTFeX4E+sB96eRny+CPAyU3rlxWmZ81Kx5QtXqqN1bM7BE5vYfK/T+Khjy/jIycP4M4tgUeGEwzEwIAIYWqBOEZDE7hQTpGHpfDbgwmJHSHQYbJwsuB5Nfv41Fli+yj7qBHfIJrPJAWM8owqIXJctpkqiWrD06SxQOrpCxZnOdwuLW5FHxSeaKbFx5njXUot5LyFvQkanHiqe8q74wwZKMRWrfs9NpbY7944qiB05vUBfediOGSC5QDa3+tZQCRmXMgPyNcAH7vRBto5yzHWO9iWZ3G66uOB3Rz9rMDhzVX0Hi6xvz9B0ZG2uJcA0LEJepTFVhmgYcGFntQmpXZiXCn1doXd+zI89MgBfHJjFXftSDw8KLGNAUqMoVqMt7i6NWfVmqm+nhkNzzFixGYwGqnv2lzm7pha+YzAOZ/j6TJ/ptIlAhz+QeOLLy3art28ZSlInAuBgmUYyAiKJIkM0KGlTL9Tp1liH0iJ85iN1E3OkJniUgquBLeIzoHmmFYytKA2mRaKvqZYgHjesfDteHEzTru7aGS6PoI3dIond/50aUyRxu6Z8HuvZMyXCdbJgTcJpUpUkJiIATblaTw8WMfneh3cd2I/lj/7CPodaWqadDvQ/b4pD722Dr29bdwrhSkCJLo5RN8EjepRCYxLqAfOon54iIc+voy/P3EId+7keGC3xGk1wFgOXXBncG4sANCV0o/YJf5MpBQnBAFs1WpFx11jhbu3WfDexPNyCxKL7+Apx5ejECBKLzRNtwT/PAXo5mE5UseZtiBMYzSmgZSUtAWd1oljT2Nl2+5pqPea/Tq4MRbPIWTT5nG1ptiPNlcXBxvRfBPulFjo2a70CCOxgy2xgxPDPvp5gaPdPpbOrKHbO4uVzgiZdZeISgEdWwhsUpnMtVFpXK0TBT2poYY11ECj2gKGGznuP7EPn9tZwWd3czywU+N0NcRA7qDCeK7YiThY1/+eBhrzyjSwQT+nxjgxMHcpWM9LADioNKxRMFL6KShdGoqZ03524VTK7eovWODnbxZJIiXt6HlbbVQIYzUqlIayF8qlthp3jF+4qaW0QAYhDdtBWTS0vZt/y4NIoMK4CyoPUmygKp8vPz/6nX4q1YxtCK+sBxl0TaTIvQJRgBAT0wTsIvhlzSSIoUiMzSwaIagqa2jVkITZHgasVWqAgdjAcb2J5a0D+GB3DYOPZ7j1kTPYf9ND6Nx0FvLGQ4AU0JmNSs+laegGAD3bRbJS0Ce2Ud4/wtl7TDryHZsr+PR2hnu2azxS7WJbbKLEGDVK1Lp0fnPeZ8fUhfFxBaRkCNSmzgmYTcUDzJprsVbcJSULyi4gPluFAdcZ1v+lF2bhzgAP/rtzCJhjyjjVxTP1Pk0bJ0V/8/3mcbmY44duM+cKti6xtvU+9WzEi1vjHAQYcA3dKrAxHHEsB4GE6cBjioEUW8DM2OC6IGA30O4GVFAQukKph9gWZ/BI3UO2s4blPEep11AqieuHm1g7M0LnyATZegeib9nwsjZ9l0Y1dKmgRxr1UGOyJTHa6eD05jJODPr43G4f9+xK3LtT4/hkgLPiLIbYwkQPQFWBfRB4xGy0ufBasq5mSVtA6SyJ3Td+7fDs6YWR+SmSS8JwuJLmEYIH4OMh4n1aqB9qWzwL6fMFwHT6U4GCMK4TC0isr57qIlDKayYK1LqE0GNoXYfpkFDGNRAgW2+VBg8jnYe1XONrwAFV8LmcXh2QW7J0nibboYaWyneyDawbBIqFlEpbUNdUpcNfCmGKjPkmJQn6NKLRmwo/C0BHpScYqU2cyR5BMSnQPbOGjXIFDw97ePyZXdzwwFnsf+QB5NtDyCcMgN0h9PYIelACY5Muh0pB704wvGOMh+9bwydOH8Bndjq4Z0fj4cEEp+oBzsrTKDF21oxCaQHPOExHtNczXjgAAzpy20gPsAG/Nv5jLipem+P6KrfNBc5fJ1qc0t/L4PrSv8vFreIDyGPmjr8nc482g2Foc5M0ZzVdGfOWC/zvmFkBWtwqDfo9Yvf24CJKxrEkiqXZA/trbd0qLnjU6oK2ANK96oRg++BcuEGVh7ogAhuphZEMuEqPMMQmzsouRCVRbK1iVGco1Sq2ywLX7e7g4Jld9PeNka+MQJ1iVQnzbyJQjTKMRxm2d3s4PejjwWEfDw1z3LMr8NBuhUfKAU6L0xgIAzYUvPERx9RMu1fmerQxnk02pOlamdLYLpGZxceI15ZwbAn3wD0KclkEjXJXSRyZzyV4CNlFdP50N17t/WauPLm/kTWawCAo5mRfNBpTaYXM9XcxpcshvNLgD2BsYTbmy194lhUTnB8BBsFKaCOqrhpdH1d8LEGdBywKWTUiEZkOxkzEDzdTOq2gwwWE0r4yGCewZqyCke5nEWwHoBFTo1SJEkMM5FmckBIYApuTZRwf5rhzZx2P21jFjQ+McPPHtnD4mo9DdrRp1FYJ1OMMWgNK19jZ6eGzZw/jnkEXd25J3LdT4XQ1xIbYxkjuBn7a2gKNSpvg3rjomnf9sHshCuSii1x2IVEYsIQR6pb0VufLD+6ZX3Ap8LPte7eNfY9oXLIEKZhVulik9O27ZMJYrmnSVrY53mZWn5RHW1LAY1owa+w6pnvq2CxQMHE709K6uEVUPwWPmrTpqh10JCfarhOaYCOeT8h0cj3AKf/QNRmyhqamzRi74qwZfgIM6iUMqwLHR0u4dtDD0Y11HOyOsa83RievIYRGrQRqJVHVEqMqx26V49S4g5PjHA8OJU4Mfe2TDemZjb2CjWANaLkf4fZNhm1asP80cB23CInHnTdW6kLLJQAc0311VKeAfm9DaTx2g0pGc3rLHY25FFJ5/KmX01nVNqbELAihQqQxFFRjIZLI3QJqAINx22gov7BHLxUXevGovgaNJW1FygphLwVauGmRA9K9OYLj2sqj3v9KGQ70N48RsJ/vgcZu3tvpCobHm/DPzbVWrqS30iXG9SZqWWIsd3BKr+G+nXWsbC/jQNbHkX4f155ewvUPH8b+osZybtiIiZIYK4lBJXGmzHDfrsDDgxr3jzdxSp7AQG6i0mNI22tAIoOyQICDjTiLxgPMJljM0TOK3Bb7AkKF74I79ezFNFa0fKzAzRY8o6YMdDhHiUvVmjot7UovBlUpaVOaFwtkzJsVkEyXnQIy5hVuuHDQMXVbmhNjUzx48cGjUuR2fsaAowKBQGwVR0YJc7uEEwjBBf+dmFaest2mC8zPposSMG74Ug+wA4VSjrGtVrGxvR8PDnq4t5vjYK+Lg90uDnYUljKFjqS1QKBUAoNaYrcSODkWODvWOD4srfGxhV25gQkGttCaL0Q5C2w4QIhwzTJnn27O1gSNs1yps4vipQAPPx79ezSrjV5yhiMu7sIVTBvLEaNes20c6MlusOAvW/uF5ZX7qNy0QOY+p86hVHyHFgm+EPFiXEYZ5M41Ewd1pgLSGn43IV0rcym6psZIHRafCfz1gHs5lC5R6Ql4oy86hrRVWgHYRY8CaA3YgL16QQVCuuZapVkOWtAssIjPhxSM65uSsGSEyJCLjisTr7RCrceoFFBrk5JW6QmUMp+PxCZ2xGl05BJOq304sbMfx4dLuHcnx1qng9XCtPUuFVBpmHLIpcaJ0Rin9TZOyYewq06jth1YM9FFDnNsrVn2iQ6zT/i8KTaIPs9FF1JIw5JoGHZDjd294PeWegF59iodGMiZLkQg2j07LYuP2bdAJo2L53KAGaF4JiwGb1xS1HN6u/bS5qlFv40Vid/PWemJvFhbPEbo+sgaDOxehZjYWaDDzy+d3aK1yVgw7jtfu0GhMoYS6LnnY/GRW4yRBtgIgYYHGeRWzYKyBLF+55/Re0AxPpUaQwkDridigLEYYkOt4MxgFSdGPezr5FgrJHpZjk7myaNaG0/rsNLYnFTYriqc1tvYlhsYYBOlHrp1AEASbLTVcmljNFJgo+1+8fMPxo6M5DDIt92NF49hdHAOgcmjZoBcMsARBAiJ+MVu3pSmsm+P9m1LLUsxJpx+BjsW+eqFNq9eJmr3YtdqjFpN/EIUFRWTwip3ZCbtUZvFhUCKOX4WzIkvaI4R0YZZUVDIRIEMhaFXE1adtFZIDR/cGDMvNEcpPW1LD6zSlQm5QLiwGvBRMSXTpFG5deOt6DBGw9833103ee9sqfZCLKEjlkwpY01Ayrdtr3WNGmMIkaEUQ4zFFgbiLDbFIzih17Ay3If+oIceOpAQNnrClCYeCtPLYoQtDOuztqKrXXQkIInhQu2tGkTPH2cyZDdQiBTzo1Ci0mNUtqAbAVcD+Ar3nFJ/FKXKQIHFbjV6ZqTtvWAeoHAxTCkoCQkhC1MM7zJOjZ0GIvYqMYuZGjOu4si3BaYr7XmCVgmABItPVA8jpeNjFsxO3o1DVjIfb9YcuSuWiymNryAp7Vt4tw0ZIAJh4TkRuBEtI9rCaLYBjfj55uCDrp05bcZ6pix6l2Zu4qsqMUIpBthBgS25hFNYQW+8jKXREjrIUCCDEAJaa1BVnQlqDMQAI7GLgTBAo9KjQF+TzAU2ImBE+7XdmzZJZ1aFrnr+c6/AwxmGIn/U9MElZTjaXAoAu5kJasksiGEQKEn8UNIDErerp1gKUgiBxchBAJSp0AkArsdLhRhs8HMK55gFdWd0fK4B+2J/WrChhAJ0aedpCrsL2CBZq2wa10x7diZuXx+wCTJEyLUyFr4JLA0VlWc6rIKJotfNQEa5wLlLcgceYlAZWzHuvmhzfKUrp2Cb8SgJFkzUUML0WVDCWDkDcRa57CJDAXKP0P6m26QBAaUaurRiApSK3EgIlUoMeE1MT25YEdGFZEGuLrvFMhv0zFD6tWfA6FmrPSCJlHdg0SFUXrHlnBKaq2di6N/lCTr24l+eN7gTiBdy3yGaH/t8yqDPk03g3lvnWk0wtKxHRuwiaxs/HifF4iTnwv+Gr5DpMvgcg8aBB+1PY7ecZwAumgHiPIaNV5Hm80nNM5bAla7hDBEqqLYrutgURhfkoPoTsExO7ao9V3ps2euwc3ObzAM0SLgrpQ1spAqH7eV5TK2Vbe+Hu7Z4dN0qjzLg8BUF3SdCJutfcJ/TtLx7P06Y5saRqPf/szgLmbMXO02p+v2r4GUnsGHm6XuwNJga+txWqcxEAVep0j5cMdjw52DnC+lqf+TcO8HyrgUkagYUaNFuPHzBy+1pucB6V9LVGTFgCxYURZYNp1Gdcskd2MikdYvITuByiBuUccWvLKAzC/kYuehOtzLZcyG0AqRtlaDHqDAIjmPulV/IajU2JaYd2OD3LnzxqHEb/5sYBylydMQSCtG34EZiJHachUTMhlKVKS4nC+cyEshQY2yALHum7GQDZcGLPtG5UIM2qng6j3i/7XzbP5qSYjBJ4iJKbT7yeYQ/czGbMO8iFyv32HCZtn/KWHC6ZcaiBDQBJ2c5YvdKq5stunZUKI6ew0x0vaFGRpYABDPUYuARii3MaEGGBxQh0JiXtZ52DWLGQQEQWqIWpk5OSSwKMkj70+3LgsMds83etXmeMf45ryHl59dsUx8zXxykpEIJUtKm39u+D/bFbHB8MeQSVxplf0ZWzawXj0dnxxVKab9URkewwPMiYpH1EFvgsYsrsCx5/Qz4ssR8IfCKoXYPl4RNR4tevNQ5kigWh+FefKLU2bm0xRtQZ9lGSi+NL6oQoDHrxljgLMbDKRp2DkK6GA0pPcPhrzFZ9aHSpzkrVaHGGLU01gY/5/i+cAVP14sYJNdkid0zKXLjloJ3Y7jzt9vVqGy9Eg/ITC2WDhuHZdSgySpoKJfdEgSawl+b1H7+PBi4hXIKkFgoAA6km6J0YZVD//z5Xi8mJuXyYzRIZi30ewmOnCap546P0wYcpgGgacdqZ51S7oH5rFOAFYWb4Vppi21LnU/sZjXHKRpsBzdEzL7Nc0sFhabAhr/v89eoaKt7ESzs1KmZWAOdBgP8/GNQT2O3MRkk06rDpuRCBDQ3YhsTLH9b5Ws3htMzjx7wuOQuFcAvIvGNq2058JSfk4QeugpV48HgYCO2IAOkyRZEsv4Bk23i24OrEDwgTWGRUH+ToGw1msjXnHdhWAVrZVA6qCvQIkL2xLycvtopxWrEcyPXiTsfkaOQfeSihwwG9Vd6jFqXNnYhRy3GUDo3FrkunGtGw/QaNQ8uXRN+PbwlYyz53Fk10jIwPE6hTbkQ0KnUGFIMHIPElXcAloQHlXQ/w/GavvlmMKB36WihoKWChE8xDq9p5oJ4a22a9dXUndjev0qPXRM+ntUihK31Iqanu/FrYY7jg5NJSZhGfBbEyNyUSdemdge5cMwglXc5qUn0vF7qTJWwiWPz29DqB/w711wE2ixK/9xMZTZoXJFiT5oFvVJ1MPh+87hXgDRwAFsoUmPELmRagKdlwcxyVVB7CJ7aT0ynYViZLmVuFrBr5sbmrhNQ1klLCQDSb3MwS7OARtr1bq6pcRPZ+enQzTwNuNiTSuqraUAjBpttWS3caJ5HWo0G0gWIYr2C5zxkWhqs6qMglzxLBWAWJnv4tLWqY1dFmgVJp8+mrBj+XVJYzjshew460kzEbBrVfR9ZIpRFI6W3QBsBU5DONZOJzFZkjbNv2JhxyqkdJ5Nd9OQ6cph4A0r7rFE6kCOlAQc1xi6QFDB9PHxxoLzJCAmmaOCBktm3GRSZujYxRU3xJN7Kb+uKmH4J2wIpm6CEPSMSNgZOOmvP3IOQjqVy9BrKdPIUErWN9C/VgLFDzTmFnZCb59TcJ8yE4uDVpE17pcVjjpxLECb42fTseXQVzHTR0NDgGQ6OymYBl9MX0qb7pXGUaEFKsaczQUf0+zQLcy9C7tXGWDw8asb5mW3ajZ+YPaDnJwhqJOzJwQ7YO8DrcbiCYABarlUb2IhZ5GlzcnNLvUctYCMAp2hek9h9HjOk/PNp+lwGc09nSqZirtoCTVPzobGb4zZdKbG4Zyhi+vk8Hm2wATzqgKNpTRkavXQXLYO3zp3LIXoBw/3r1uBIGp9/nh4jBAzkGwesxavCWAvuOuHH4Q8pIU16+SjFM35RDCWegUqqx3MydH6BzBaTAgAta1QKgC5BVUtSbhMOXKQo0BFLWMZ+FLqDSlSYYGizWkwcA6XKuhdI+XOl4kDcxQKEVh23UuIYnFqNG/fC3Q+rwKjMvb+mPmWUV/gUIrNWS7OXSBsV7ufDrFwWHEzXAAoQmXTXw2cehWCj1mWQHm22MQBP6So4X36urjW5biqh4NpAuvNzQcosVsgADc/C+PvC7oUAOz/bIPASKJl5JKaxU2mm0yRmQ5LADc33C5htKMTb8cV4HteK2yfabh7rNpUWOnueaTYjtUiSSA3ULfOI31M6D+deQfO8POiwz6tspskTWOQxUakg8Qw2mDtyeaQW0Xhf9sdckgKhsUwDG/zdTq1BfG4+SLd9O34Mfhw+Bp/XPK4pz+o++i7WS8RwMH+z9jnfCqV9MYyQ6yDFDKToM5eeqkP67lwloG8Zm8ARfMrnRxkmkpUnpjgE+t6VTafzFFmTAiPLHjWAMG4gdjV5i8IXDOMWWyYKFKKPQpvksNK2Xq/0CKUamm1Apd3D1LWws6S38FJBVcE/tsjHYCC2TDmo5C+NaYxWBvu3KXluPXgFzxYulk7K3WxBp14CbRZsEHtArjHDvCgHKribQsoKHGiZz0JF648bpUfHwJFfHx2CYq7M3ZgsC4e2o/uVYsFMWvPlIJaaRwgYeRDkuUqb33qaxAxqSnmnnr3k80h/C/9M7m0u6fblSYuf6ab4XDjN7r5DeJ6zJAAdZgAPPBILc9J9kgAbcQyHK3LI9CYVXoQuHas7q1BePKdYN1wo4QAgTpnlrr9pjFibHpt63Hm24a7mKeyemc+jow8uAeBQADJoa7VRwzZlYzCEzgzJSA8eJLQtTNNQRvQisViPeVFbkoaay+KQrow4YG+csi4HUpgS0LVpCJeDFiszJi+NS+fi0m7BHwj6XrpKlTyimiqIOpZFIwgIdXEOAo4doYDJEmMMsIlBfRqjagNVPbIBnp0AqND5GhBYOD8uzZPHPwTsBnPnpBRzrICkzJNKEQC0oufEsFhOyYkpyiZ+qXj9C1kk96EFPJMd5FnfFR6Lx6XaHDwg1D1zKgIM0fG0rhtupfgpNNcud5Vma1RQgjX8EzKYm5kPzatZ5yUGHSa76PyA+AUXrQz1yxYvHsi7F9CRcrHNu8CkdIKZz/QS0v6aTwcdND9gvgUjJfOcC7EZnGGdNYZy73OzbUI8NsAWbtE8b+6GjN/1BvNojSFi66RNY3dWP2yKOXVt1mVyZiljNDWn1vNqAZbxNWyL22iCjZoGnnmcc3kWUuMoYoPseCpaS2i/FKP6aBkfl4bhsDdGacNsAPwhMeyAgGcWeJBPTFcBEd0YNT1qc61wC5yLqwIYLbo0Li0IgcIXpVOOxkVkAZACaimNHSq4pc9qLiTOh8/V1O+3FSmFf6ApjYvmTAsULUYUbKm0efEd2BATlBhjqM5iVG1gUm2bWAkloWSFLOs0ephQcSCl7PWxBYIoJqAZa+KDgF1rdWbtEFBMKRs6P3JbCGHT8mhhZ4Wu0tHmCeaLXWMXFxA8P5k7T8MQtWeRhCWOm9aDU0yMZaLAPsXmEyysNI/o2jiGD9IXaGLFkcw1ooDSsG6A979nzmXZ9sxfWtH2/8rdd61ns5MNKp8F9tE7fK6SYgmmbj8n6NiLi8jsxwJe0bSWm4GzHmiQO5D3VnKAmRYiJCzgOcBG8JljE1nn66jrdYPpYEHsEsQomlYAua30Swt7LSrDxGKMCQbOQIsZvdnXcroROSu2Yl6J47KmFeSj7ecRN69E7yb+fSPegzGeqWsQM+UXWy5Z0CihY/rJP6fF2NWqiChX/uA6K1LBR1UjQuJTXqA2qjF2n8AWinGLow2cpCAkfkPppwnW87EpqWpztEjF9Bt/8Gs1hhIZoCfNc7D0fya7LvuExDBFgIU8qDDGGDso9QCjagNlvWtdAhWEyM1cazNZxz64Wg+WudEqYJL4dc5kF1nWCdkBe14U/0HsRy67JmvGVhTN0UWOAsoqxAkGKDGEEhVKeOXYBjJS91gpn+JLc9VaJaucxguyshko5tyaeft0TiKTpnYJfy4hHQB0PXWc1VEG81ZRC/oQePqMAQK8vGZNnYhYD90mxCAVJh6HBe9dHu6Udgl92806BkDzPYmZNgKnbcxFSs51wZmWtuuP2cyY2YuEVnwd6E2KeSN2M7fF6KjwHQAom1GlYBZxDRPoTLVoTMxU00XScJdABnor+DzQz02QG8e1mR5QPRSijx5W0NMrKHSBHDkkhJmpKG0V0BwTvQMAqJR9l1iRtLmvY9u9bzFIuKSCfOeVmOUJnuUZwCQ1V/976DJuBCGzjCpzMHtvdFz35OJnrV3SLJXYuqPPTIAQEPSMCG5OlC2iWZaHaD5QSRqJW+UzFBJ9nyF31kNK+IMUWAxCNcBG2zGCubEXmvd5cT58Ybrc5nLJFZ8qYAJLa1RQEbVMFTZH9RbKaoCqHkE5EEPNmqQrUhX4V0UGrdg5aIDiOhzzw1Jh42suhSnLHgONPtbQ1X10dRcFcuI2MBRD7IoNaFljInZcv4fmNWxZjMj3ywAsNFwWiuvbQIyEPV9ijwBAagJdvnAQv1+Z7JpgUwnHahEoiF1bPrXXWxTxfO3NnUnv2rvl3oM2BRgomClK9FKLhg7eGWA6s9AIDkwwCIFLqWW/NomVeRsIiZ/xaamynjU9v0C9NuuZng8pCnTlCjpYQgd9dHUfEj7XqxRj1KLEBEMT9AwL4qNCgQ3QEbM+iViO4Ht4xs79HrlRMtF1emtF70NPL2EFffREjkIa7qPWGmNdY0ctoRAd7IjCGj6+wzIZa/PFNYTriTP2ZjFqc2zjrpOI3Bgi1BvTj9NkbvcKglXiurRl5oRjX8W9VEhooUoF0ejoIefsg6H4S1tpm6ra5Z45Sfgj26yjVDl0vk1YBS7sZlsn6u3z83CppbIF+DAroC2+IAAvbpHJXc+WrlzBMvZjVe1D3wKOClS0t7KAR2MshjABujvu/FISWijtFpmfO097Y4GMjGaWkMjlEnpyDcvYjzW1D0u6i57I0csy9DKJwtaAr7XGoFrDRrUPp8QqUABDYRiZVv8si1WI77NB/woQFP9RA/AFiUhp0nkq+hkU/yrs4sMCfIVp9gYFCNlMD2xeryzJhnHl50ALV9Tw1jOPbQjGtctKzWOdEkoyvj6XXkwtjr3tMd/i4oAyPbv8uiUMk3hfOpb5298TqpnBWUkzvEKK6UhlJMwz91jiFNCAUYN0C3gfa1hV+7GKJfRFgdxuV2mFiaoxRoWBGGAgdrArJMYKqEXZAA7+UjXBKj//mEXi8+J/A3AMcSH76IgVpwv2i2WsdQqsFRLLhUAhgUyYpoujWmOn7OLseAmn1CpOygKDrMBYbZlYKuqInJgLv15tf8fn2nYPZt2/tuu3V5nGdvDg+vS+oTE+C4w5/fgoyaUt/IXQ508vq/NzC+ZXE7BdSkPmgFNFEhJUwbS2/Tjcttxvnnhwgps75aGhJm/cgkmhe670ybXi/I6xNYcwlZWnWKloLC65bXK2pg/joD6Ig3kf+zoZpBColEalNSa1+Tmqa+yqCWpZYiIHkHITopYAexiDlw0tQZx8gWSxD47VsHUmasZGkDWzLA/iiLoOB8Qyji11sNaR6EigkEAvAzoSyIRGrSV2qwynxzkeHnSR1QVOFw9gWxxHVQ/dSxLPscFwJO4jp4yp5Do9E7zCqtm4YsXLeOS8uZe0AGlRuGdx2rULr22TViUJAmgTtKtz4zGWTopiZn2WK0XoWUoB3tRCT0LWYeM7EW83/ZqkUlHbim2Fh5mu2IOYk+B4U3TRDOHsRi666Ik17FeHcUAs42C3g7WOcEC+VBrjGhhUCjtlH5tqGadFji1pXIhKVODxcvE1aDM+AldwApi4eQbAaAWr+iAO60M4mPdxpJ/hYE9gf0djOVPoZxoCGqUWGNUCG6XEiVEHy4N1iInEGVlg07aad+4VMwk3Fz6HWTEc0/6ODYnYrZJ8Fqa46s+XaZwXcKckrsb6aOuHy4LhIL+++cArBFefQzQfcsNwVC7MPxNd9yCE5X6ntZT27o4MXlHxnickbZXi+FjTvlNR4Ktb9GwzL3I1AIbZIaucL4hcpCzQlWs4iOtwk7gGN651cOOyxnX9CoXUKJVAqQQGtcRuLXB2UuDkqIvu9g3IZIEqH5jGYvUI0BVbvDzt6efra500rCt7Dnx7xSqsZra66Yo4giPqWtzaX8PNqxJPWi1xuDtELoyGyKVCLjVyqSCFxqDMcXbSwWd3u9h35hjuHqzinrzAFh7EpN4NQJJ5RnKT4aRKc08VZT1ZIBQUUpKOFaPaJhXGgAtIDQuBUeaU0va50KWvgktjBop5ihVECymaSo4vsHTd2yTwycIqDw7G2fUJY6Aq8w8ajwaFOq9w91b6++kxEoHLInJL8YC7wK3FU92j553caEAzrVTrugHsgmO1BPa1MWBti1DMngTfcdbALuJ9sY796jCuy9dx3UqOa/rA4W6Nnu37UGmBnUpgu8pwdpLh1KhAb1Cg0F2ckhZwixoKIejg5z1NHMMQsdX8vSBGdp8+hsN6P67r93HNksSNywqHuyUOdSZYLkrkmYKERqkyjKoMG5MOTvYKPNDN0N1aw8OjLu4TwJY8AWALlULwnDeZp3ScFL+eMZPTFrTbJntdvP0a1Z4JNc2dR/Oath3fJpWBZPQd6YO9yt5jPi454AASaDCinPlPLgQ6hJS2KqVsFK8RIoNoo5ahXNleSCAXXW+pW4o6VXAmeDCjzqWpOVJcBN8mXrApeIpKZYesSeayWWifQvaxLo7hOn0UT9nfwdP2lfiCg2dx7XWbyPsaqgbqscRkN8P2Tg8Pbq7is7t9dOQSsq3HYTc/axZurVCr8CWcB2xwV1CDnYKC0Mo2b+tiWRzCteo63NxfwTMPCTx93zaeeutx9G+0L9xEQdlQEtkBkAtAAdWmwpM/u4r1/AjWNlYx2bgZdV6i1g+hUkM3F4rEV5CopbSKU7n6IbxarGPVBA+yk9BKocY4yYZR5hRlp1AhMoHQmk5Zw6luwqbjLxq0Pl1P3huIjh8/YyHY8LU9aG5+UKtoQPUMFHjtiytJ4uC3FFXuvmth7Zr+9bSiTrlTuT+cMy7xHOaVIPMk0lHTXXMhayBEho5Ywqo+gGNyHbes5XjyWo3HLw9xbHUHvU4FKTQmVYbBuIPNcRcPD7t4aFSgI3sodg+i1GNU2cjoujpalKK5zeuSiD83tYBM8cHDej8e1+/jiesSNy1VeMLqDo6s7WJ1fYRiSUFk5gWoxwL1RGK4W2Bzp4+j2yvoZT0sbfcw2T4GSGBDlu7ZhwjBaexSSAGPlAHgz7W9H848kjIsUtukGItp+3qWeQrQiNjROLzAGyDnow/2BjouIeBI+45iHxW3NqfRlRRpDcAtgMHfEeqmMWs9MS4BnUGLAjwIkoAHLzhD4/EGOc5iYuMH8zsH2oo/IDxACjB+0K5cw3XqGjx5Xw//8NAIz7rxYex/GpA98RjQtWyRMmWjD57ewePufBg3faqPnjyGbtbD1pknYFRsmXErM0dyjTTn4gtkxY3qQtrb1+aAMBVb+3I/rlHX4ymrq3jqPo3nHjmDxz/xNJa+7Ahw0zFgUgGbO9BntqE3R+aAmYRYKpDnGa4/tIEXfvQhrBZHsVuuQ+8+AffkY+xWJxvzNFVbzfPgOrwSoCTCQvo22CbILu2bjs/f1EKhy8prwqRpfJdFkKiIagBmJ9jeFDirXIwGB87TFjh3LMCBjQAYQjrq1IHGy9DNMi22JKTs56x7AXvdpsUptTQ8C8aKnw8GOnjDSJVY0B4NMYDbuFKO6P24fqXAk9dqfMG+LdxwdAMr11eQSwJCCuiJQrW9g/HZDEdPL+PQ1goklpHLLkabRzCWQ9TS1viJ6sWc7znRPLtiBfvVQVzT6+OGFYlblkvcsrKLG4+exfLREsXB7P/f3rs1SXJj6YEf4O4ReU/WlWQ3mz3T092jXUkmk2nX1mR6077uv92Hfdmb2b5oHiTtmqSZ7hmym2yyilWsS97i5hcA+wAc4AAO94jMqqzMIv2YkZUZ6QGHu8OBD9/5zjkQh+7d6LTdjKw19tYbHF422PuxhcADlGIPy+4Qev0pWllDSw3olU+oyNks/j6kjIs9NPl9hE3b1SWyq9Az7QOdY5fvjbnnYpd/kpGVbarsfJBntG7D7gRw+PoJzm+udRuVMPduBPdi0+CP9AMJ2ovoNAYYpLRiPwl43z/5671pe87OuQAkKr/z9W4VY8PIAESx7QLaF16j5F/2GnMlofPx6f7hC+WzWEYZNtm1S1l6duM3R/v4Hx52+Hf//fc4+V8+g/n9X0E/fQKU7rHu78PM94C2QfX1n/DL//QH/M//63fY//Mv8WJ1grfNr9AVKx8+6sWffteV9/Gl9CMXNHJ2oJBznJon+M3+Mf7tY4X/8ekb/NX/dIXqX/8S5l//c5inT4G6hnj9GuKfvgE2r2GWLaAVxP4h8MUTlL96gs9/8xr//v/8C942v0GjP8F5/Qu0coXO1B4IcZDK3RY2f4WEkOHec5GyzwCb7HpTBoEEyvS8+L/0vKUriEdjRjLanguUbXTRHDwteatXaOl9MJQ3g4Pmvk+YgK9mYzpickTISBveEX1D+nRXu2loHQOrgy3HivvBjLPZMRvfy17E0ADw2DUMMhKUeiZli1akB1LzbqFB0ORYg5k4wIl5gF8d7OGfnQD/6sElfve3r7D/uz3ILx7ZDYgUgFIoFzXm5yscfHeFT75bofr+CWbyEJfNPlbrz1DLBVq9tiD9PQ0TO/9Kx8I8wpPiCF8eSfzuSOF3Jwt88eQcJ7/pUH6+D3G6BzGvYDoNtB3QKMi6AzqDaqMwO1lgVirMiwdYqSNoc4Tl5gkaubLzgCadHfouNHvTsvc/ZXW55VgOPi5STc+uYINbNMZEeKeB3dw56TXQ96Jso8k40u/FvfrRuFTYQhrt5gJ7QBVGt/nfhtwuHD0KIZ1vv40mYU+bgwLHbNIxCxjGHzK1wwFO6FMqrOoXMCIRmXF0vU3wpVk20v69KeQMe8UDPNKP8eWhwN+eXuHoX8xg/tlvYD77FJjPgboG6g1QFsB8D6hmMJ99CvnbS3zyL8/xt+fn+OLtU3yzeYqL4gVqdQWjtkyOCYuTu9/0Gd2LShzg2Jzi6b7Eb4+v8Iu/vkD5u8fAZ4/swRcXEKslxNkFcLGEvqihr9yCub+GeNoBp0cQv9vD0dsl/sXXV3i+PsU36ye4kC+gla0lo1yRt9yLT889Yqd8llgFBVvyXVElVg6ker7OmGXg484DZEaRKzeWZAKAiFoufQhz66IEYh0JH325jIdcWJqCvdxziVmE2wAd14s2IeMbkJ4rETEwj+nx4YgQOt4eF1Lc87bTDIzb+7l95ymZW3cXFmXX8/XmE584z7Ibp/oYn+4L/Pqgxi8eX2Dvr2eQXz4Anj4A5jNASkB1EKsa4miOqpA4kkv8anmBRVvim4N9vK4P8QYnWIszaJGZf26wkPL+VmIfJ/oYD2YlPtsz+MV+jacnCxw8alE8mkMczyH258CsgOg0jBSAtOwMAIhDg6oUeKDX+LIr8MN6D5dtiZebY1zhBJ2obaABqzHV68cAyNjW922W0/TkXH/8/UxrsnCASVFwEVhg1bmH1roh64WS+w3xu7hTrj+H3AHgMHA8t/cf0c2jKI0gPtLuhVUR69FrMXK7KA8CgMxOUCjvx++7buKQV57ZMz5foMvJOItio1KorHpgXEKeiqB3AKw4tpOAdnkueIrq0L4FG/PiBKfiM3xeHuI3hy1+/fkZ5F89hT49AboO4vwM4tVb4PU5sFfBPPwEODwAygJmfx/l7x/ji2cv8JuXj/FPF4f4QT/AQryAQgPK8kj3hav1e37cgecQ+mrp01Ps4dEc+Px0gb1fVxBPTuwx//gniKsVzNUa+mID9XKD5q2Bqu059roFqtlLiL8WMJ89gfybT/G7L77Hd6sD/P35IZ53D9CIBTrT2HtlYJOj0bPruRXY74JVsXXZUolFowyhdI08KyzdEz8x8MgU/nnCuEGUkECUkGkPRyhQWXYDNTpRQ4na60hs/7powSW9Cj2Toefg+5poSYzpcLv6DV5S9DrfsbSuNg0kZiCXHLDd583deaHFeJEcAmG57+aMNDW7LDxcsO5vRWZeHgOFuxixG6WY48h8gsfVHn65r/Hl8QInXzSQnz8CHp7AnB5bwAEAnYKoLNshOo1SaTy8WOHLzQy/XO7h5brCi9UDnIvn1wIXQ5oa+pvtqx3zp3IPT/YlfrHf4PPDJY4fbVA9kZCHFVCVlonRBugUQJugUgJVCVFKyL0SMynwpF7ibxZLXHbH+Mtijrf1A2zkAp3Y+Pk/t/hve4bXzZPCgSuBjtw9GAId1EaU28lX401BQj7SKQUwQxFU1EZwkX941+rdZRp1uxptuBgxn2SH/O1An+bnoVh2B9RnEny6XwMbwuhcE5THqkehZ3y33LyoT/T7kzId9vfgqvAVRaktt5gRZZ+jwGjRKeQcM3mEI32Ck/0CD2dr7D3ogPkMommA5RLix9fAi7cwr67sxHJ6ARwfAKeHdpdzvI/5ZwWezFs8mM+wtzrsvYR8l0iCVSpClnUnZCaaQpSoMMdcFtiTBrNZBzFz51muYf7yCvrlCuqig1oatCuJdlVCaQEpDIrzFsXLJYrjC+D0GJACe486fDqvcVId4qj9BAv5o6V/fbhxRqDpkw2lWfWQ9VV7N4S9EeCZRfmOW0JGkSqprsePaVC6aQnpAANlf5SGUr9LSAQBYHSfQeJFeN1HOjlJSO/uS91d1NaH0xXcjDmh+SBymw5Q1vS364Q+jtk2BnVX4wmXhtrftR/pnBIBaBHmk1LMsa/3cFBKfFJpfHKwQXlaQOzPgFkJFIV976UEtIGREqIsgHkFsV+hPBU4OdzgQaVwMiuxv7LpxZudd//xXDcEOqSQqMwce0WBwxI4rToc7deYfWIgj2fArLBgo3Pvcqf9z5ACopTAXgVRFpCdxuxxg8fPlni03MfJrMJxfYALHKARCxfe22fB3kX8yS0HGOznIVSe27ZqwME9H8pGWIbTsmXRvJbU5CG7afbTD213HKXi3BimzfrNaPDSBGsRpPYLBQDvL+29lHxhgQt7daJCODEhlSGndrKJVjwtGkcKaIFodxwdn3yPaPxCzlHKOesXIAh0QcNoHU04HDgVco5K7mMujnBo9rBfCuwVCmIGYLmBePkjcLWCeXkOc76CPm8AbSDebCD2LiFO5xDH+8D+DPK4wsNZg0d7cxwtT1CKGVpHO3P0Sy9DwfpM16sZOxVVW0Wo+ligCsm8Ogm9aCAvVvb3H5aon3VoFwWauoBWAsYICGEgS8t0dK9aiINzyFLCLGrICjiatTiZCRwtDzGTR2jlGq3WvfETgY1kAab7S+wFv9c0zlIBJrUTpSt3k1maidTWNZHez+sThrHqrpbVaP3vQkiUco7WRd8A8JFNQkgbiQWNEsPF53y2XcS7923aiLs3x3KAQvQoqsaCraEJO6fuj9iPnfzeSdpn9Dc+ucmcotfiz/Isi20zD2qG8owMtcv7JEWFCns4wBwHpcBBqTCbdTZjlmcJFNBJAArQGqJT9m/aAkMxE5jNFA5Lhf3CFhiwdY1k5BK6iSsi6isqVKbCrBDYL4DDssXBYYPitLDsxsxtaLSxzEbT2T5KYV3DZWGZmspAdArytMLxUY2HZy1OqhkO5Qxzs48VKgjnJk01P/RZznIC2ZwmIm2Du+qGAE3qHhkynqdp8Bia50R+XkttKE/RXdkdAY68W4XXB9BswgmTtg1VJVYgHQB8sefx0zDWT86zN9LknIIU2qnGaalD+l9/BUykFw3QNAtkwrh415Dpexn7fnbfCAD43XEBGdjaBjDrGmK5toDjYg2zaGFaDXQGptYQrdMDSGFZBglUUtukWwhFxvp9D/eDT4ytE5TRPYi+418+6VODA4DSEkYZmLoFtIZZK+hGQCvnnxWALGwODlkYGA3oGtDnDcTpCug0xEzgcNZgr4CdGOESdyG4gXjf+ZhIE/VE1KIJOVFogTNI2TfG+DhwYcWqNpyZWAoNhQ52zBSIWRUqSGegrWbDfSe+f/ae8yq0wyF7GR87dxMmO+67oFCvZYkWBQjAe4g1CALNYT3GdRbINLotzBnbqPjt93isH9vyjETHEvMjJApXc6SQ7KqVgekUhFIQTWtfa2IPVGf/bd2CzrpbSqDwcxRlac2DjaEkYLnrsu5mKxwthEAhDGaFRjHTEHslcFBBVJaJMZ0FRt6ksJNDKS3o0BooC4i9AuVejYOyw14BVEKiMKHSrDSZFAkZ5o/bkEg+NVpT0rolufd0l3eO6zMGj8m4h1IbY1Lyjf6MXCpkxmgo3UCbzov/6HMy8nH5SR7BF2UPgNdG2F9DsixS/dt2gkvDFhOTLFupHailsCwElUmmBUSZOqY82cRIO1EeVpoOClJqD6Hm9JpJi0Kf22uxC1WNDuvO4Kqt0Cwk9ho22UoBU0rIOWAqAyEFxF4BcTq3avVNC33RYqMKtIlLKVchtQe+hELn6q/wTHWcIehfO1CWCqJwk0dRQB6WmH3Sojp24NG68aFqC0KETN6HWYniQYXDvQYzCQa52HkYQOL3nRtV2SWWBrDPrZAzz+T4vzOGgT+nUE228rQ26TFscawW2kg/7mjHrgwgYMeRMIUHkIATjlJU0sBEEKfYDyHTY0npesl+7sEuJ28GVFMlsGd9rRQQbyzGFog+q3WzOiapoK/X8wRs3OQ8N9EOeEYNBq022GiJ9aaCWjQo1g2wqq27FYi1EW0L03Qw6xZ6pbFe72GjJGoFKKMhOc2fmcdyvw+WPgBz8bnnKAVQSI2iMhAzq8/ArAwCUSlgCg0YAwgBQa6hxGRh4A6HGNBO3NDDt33suHGpBfzGgn93VzaB1p6cTvA6ACJ156VjcRfW7EPYHRdvM/BlhXQHLVvwEuWpSSGJCYx9vYIxIGxHrtBaNwVCyXhRFCBRUU/c5doqsRf85KQMdkwLiYJyvlbfRsa/HDQA/UyS+TYyGT9dUrAaDdZK421T4vLtPo6vGghtgP05xIMDyL0GIBAyKyAO94DTAxvX/uIMzUuF82aGZQcol1tDitL6PjNGPliBwi2asTgyFGwL+gdjNJTs0GoDRfdwJiEOZkBRQJysUGgLiGyiLwPTaKgrjW7hXGl7sL5bmmx8rRWgpcqXJmgw+LPP6hgQh49yZsNW27XF2Dp0PTcEgQ7rhtEoyN3kWacqSudOfdK5hd5Y0KJNC+NADoVDk3tlbLIMY6+flC6XkXeQObuXxoC86exIGyhoBwxP7u+ixei19R7946n//1r9YPoNMg2NFgqbzuCiLXC+2sfj8yWqixriuAEOGqvfKArnXlE2F3jbwaw6NOcCi/UcV12BZWfg4r7s+RL26Hop1/PXqIyBMtZ9CiCAjMJpOIR0n2nLZkhp2Q3p2iNtR2eglUCnhU03lKkvBMaO7XrPcwBraJ4Oc10cAbXre0YuV17KYqxfsW5O9a5tbB4A4k0yuSw/tN0h4LAjxN6A4C/3kSHMV2V3OpT2m0WsDPg26fscFNDPlPPDK5kRGA5/DJ0DEhoKlM46Oo8oYYRCIWeRa0fKKtJ0kPmFzTEFtLRzSj9lCAioUPVR61st0IoW502HPy3nePLqIR599R0OvlwCXzyB+cVTiKYFlitgVsEcH8IcHcHs70M+ew71n7/H8z+f4OtlhZdrhaVYBRZDxnoGHybsmIAiTQ8uJKiIHIlK6Zq0aVFjhaVqsewq1E0JSA08/gTY34OsW4hiYSeTWeEnRPlmDfnG+WDnEuKkgiglzLqBetXgzeIhzhvgCiu0Zh0XzxOhmq8FUQSSKEOoBShKNxErRfqYUuyhMxtobafx9HnYMNoOwtQwUDYzrGOdDDQ6U6MxK3QuAVEsvg2TArnUDBRa7XQKpvUTDoXPcobFjkvNGEALNiiiaWiS0/7d4jH399XCZKhNC2lK785MdQRD9UyAcdp/aNEYs/TcN2MwhhmoseNy2hJJbkRYZmyNGhfNPp6tK3xydYjT5xt8fnyF2ewSsrRiUcxndgFvO6BuYBYN1FmH81cHeLY4xPO1xNuNwhLrPujNuCv7GpW8FoZfVytaNMpg1Qms2hLtqoBpNAxtPGalc/MUQKEsSCcdBwlI6wZmWUNftVhfVbhqKyxaYKU7tKKGZvczKnMx8jxyC+8u4Cpyrbj7lDJh2xhIPoaHjuH3Nu03j1rZxmrQZ5EU4ANvQu5cNMp3NHzxB2hyLvyil6P8gb6/fsjic/QV4ULYFOmUt8G4/QMl4kpNpkW7+M56YF7n0R/Zv9NkSWyNYKXOHftTizUu9AbfLiocV3v4m29O8Otv3qJ4+gDmV49gpLC+2/kcOD62yb9UB3zzLZZfafzD2wf4dinwqlljIxeB8cm5fIx2LgHSNrCoDaqc6sJ9gTDolelQmwUW2OCyPcDVZg69XKGYz2AenEKcX9owPZpo5jPHtRYoSqfrmJfAgc0hYK5qLL8v8GxxiLe1wlJeolUrD1T98xSWHeLJ2aSQUG7HTGCDsqaSELQUVp1vWaQ6CwC5pqMQyiYTE/CTXGtWaPUaSofQVp7QLmqLMR89LYYbPzZpWTwh0f33zAYTvqbnIHcfsQXBYX+fQQeQ5uPwH7OFPxcRMrTYpWHoUZsju187nsb95TdhjYbOOQQ2UlcxZzkUWmzEChdqH8+WJQ7LCg/PT3DwlwYP5mtUlc1pgcM9+37VjQ1FP6uxfinww8URvlvN8XJt8LZtsBGrMHaJ0U3z7wwCvXxyOgDO2VhjpTosuxKXbYX1ssLxsgHq1r7n2rEckoR8Albsauy/SsEsa5irGt2ZwsXiGG+bCletwco0aEXdYyaH7v+YRekACPAOAUPmWpEMnHGXSE5rlTITvl226eDH08a7x2oMgOFcmzHbeTcbkDvXcNiy4UxJL8INBhxTAZueWaONHkg0wSDQ6R5VMqq9R6+nTArgwUhnNq6EufX5h0RcLt5ahEiBMaUwBzRK17Z/kmdJdAInVm019LfwjIhk7IaGQoM1LsQ5/rIqIcUBTssn+Hf/+2v89eqPKLWG+fIXMI8fA5Xd1chn30P84Ss0//ef8Xd//AJ/92aOry43+FG8wgYLv0BrGa6R+63tTrpGgTj5VVGEME2qVup1J27HvxBX+HH9Cb6+PMYv/niJR//8hfcpi4OZdf10GtgX1iX0RNhoGgIhbQvz6hLdsw2+fvYU/3A1x4vNCgvxxobE9nyTdpxIBty00U4fwTQncEBJWn2K9zV7V1bRGzN+jGgN7QCWdpogY+z1dmrtAYAFDKo3PlLqk18D9YeS1EEwoAztQ3HTSYwzU1HMPRSUaWAMhQred7cKuVLCfJCjrMf85O+D+me9QYF8ivT0GfDzj7Eq/LPcosM/T8GGTMapgcZaLHGBGX5YV5Bihj15AOAxftue4fH6CtVnG8hH+1YbsemgX2+w+JPBd88e4O8vjvDHqwLPVg3emEus5CVUTwi+XZOQYz38tbqNWy3WWOoGZ/Uc36/neHh+gpPvf8R+tYQEII611XP4L2rvQjHaWDfQxQbdDzUun8/wl8tjfLcucVZ3WIglGqx831MNxCDDlGFq6Hq2XRedh4OOuLozmzMy2UrHWJYUSPTa9ONnOH2/vyZIv75atr5z8wFtQj4c6LhzwJGNv3c7PL+bM3kXxRizQbs/Mk7527/n/W4KnQug4cW7dJSvwS4m1SDYiHQaUNDa9rH1BcccawMSr2oPMnx/XSVTy3CwMurQaLBCKzZoRYN6+RTACRbqKf79/1Pg9/gnlP+2g9nfB4oa4u1biL//Gov/7Qf8v//wOf6vH/fxX842+EZ8h0vzY9jNiBKlmEGL/g7BRwXxhVFUvr9SlN6VpUzt768yNRbyHK82n+KfFnv49Plj/Jv//AJV00EcOYFm0wHGWJ3JrAL293ySMswqm1fkTz/i6tsS/+3iCH+8AH4Ur7HRF71xERZk60YR0E730xdjkouLV7qlBTo8x3iyIRcHABhJ9y0AAKVrKN3Ek8JAHiwCGnFSMTaeM8mLKEuuX2xYGDN9j7sbA8NhJ5b7ViG2byZsQBjzE+3gBmjknI0JPbeyoey8HISnupk4Jf6ONDzr39hCHuVo4Pok1ndlWjRihUsBCC2hlqc2vNocoNESv21LPDxbYP/pFUQJ6Npg/bLEtz88xB8vjvGHqwJ/vuzwg7rAhXyD1ti8NsQuU72YsbBOvjgPmTYanaixEEu8qa3757Q6wKPnR3haXmFPLiEbZZmOQlrxO4XIdgpoFEzdQb2qsXou8eLNCb5ZzfFsBbzp1ljJS1vqgFVNBsbZjG1go/888gm8CHQYo7NjLrAOQSCajml+HO/LNssBGd5f6icZzQd3VTH6jgEHXSzXcdhdc8Eoeksbh0iTfjhsX0jVCxekn6Nda15YQ77ymHGQnp4WwlYX5W4PPiEAQCHtdVCEl9YdIIHONJDuttvvVJbVESxxFEKefrtDtW6BgFDdQooOnWwhlwK1OoIxj7D4Pyr8/us/4fhffgdRCnSvarz5wwz/4btf4z+dzTzYODfP0Oq1r+lRiQOUxZ7XOZCmIAJuSdKZNNdJTmndmjXemiW+vpxjvzjC/t89xl99/xaHX2rIBzOIUljXyfGe/VJZwFAtmE0N/PgWm39Y4OsfPsU/XhX4y2qNC/EKnalBIXfxjl6j0zUKERYJlQhLuVA0RJwwNis3dohxIxBMHhEOTEd2eH0mhrMc8W63933DdD6QUeXg6Lqc3odb7Ku97+yGu79+QlRZ0HZTsVtKR+c1Hv2NCIU/5nIF0TH9c40XmaPrGArhHQQb0UJl2++MTel9IW1fisVDNGqGjdrHWVPhl5dHePzDGoXUaLoCr9b7+Gqxj2+WEn++7PC8u8SFfIONuUSHDTTdI5bcMJe6e+j+5sxAoTM11mKF867G98sSM1nhsHiAui3x2eoC+4+vUJwWEAeljV6RwgrdNwqm1lCXCqsXJZ6/PMFXl0f4ZinxYtXhUiywwcIy4ORmvOkY2SIQzoGO9Drtcfkw57GolN6xQyyHd8fyfiWp+pMNLAfw+g7ZzjtnOCgfB/nXKTTRvlzaUlUObGi0EeWdsyFBVzoRp5NFjz5nqZVtTLebfJjbgGh5I2wpdiErn8uBxonfaZrW1ysxQnltBvmIfS4MAJJR6LZ/nY8O4buqTtRoxAqdrHFZf4rLF4/xzfIUvz87xu+/2kDA4G1T4dtVhf9yZvDN+hJ/kV/jSr1AoxYAgHlhwzrnwqba7kSN2izQmQ2ga1CNVE7hpfeL/52bhkZjVjiTr/D1co5GH2DZPcBvr47w3724wKePrnDwuEX1pEBxUNk6CpsmrC/rGt1/fYlv//AA/9/5Mb666PAML7BSb6BN5/UX9LJ1pvGCXS14MrbCjwHhcoMQ2LAsh4R2ESY8JDO6Rs84sDBTxoTxiYgzFTlmZcyCCzGOLvEp+92Yi/UKAWyQlkaDxt7d7WZuZgakS9G6s+4lo5PJMxeC2Z8PxpiNVOPB2013ogQ6wvn7YJHnUqBxwtvdBXRk2VK2EYm/p1zEnm3Txqq06ESN1jS4Wj7AZXuIF+sST/eO8Hh+gEIAtRY4awT+sjB4uW7wXJ/hXL7C2lx4hiC6Vwx05CyXAynn5jSQUKLFBguciX08W5XQZgYp9nDRllg0Mzx9u8DJJxvMjmsU+4AoBUxnoNZAu5JYXe7h5fkxvr46wj9eVfj2SuFZu/DMDN9Y9J5N0l/6eehv/PnknglvP/15m+2Su2RbH+jzlIXj3+ebJcHcKXfJdt4DwGGcW6Xrid8kpA9RdIcC4ItejOAo0Qstyj2/aLI7SMVbvD3uC6SSx0LY/Akhv4KG0QqQgDASxpSgaBq7IyiiF5DOaWABTZSJklHl3ieY7Lb5tdrJwL5gK31mJxpd4/L8MV5vDvDN8gBSAFetweuNxp+bt3glv8eVeoG6u4wGq0Rh05CbfbSighY2goIKm3Gwwd0nY8I96bIVGiiscYkfxSu064fYqCP8uKnwun6IL6+O8MtXKzz+fonT5+eYPT637XWA6Qy6hcDzb0/wH18/xB8uJZ5357gQL512o6+NoGeS9kki6G4IHHBRqXYvIwENEorxcUNRK7ZGTnie6Ts7yFKwMWx1O5XbifVDkfv+XDaWBWBMMskj1m7Y61fRBPMxsBvBbH/5fDBOHe9eHyN959N2eixG4sLJts/mnVQQeBPLhXVH52I6qfB5SFgHAbSiQds8wlV7gFfrCicziUIArQYWrcKrZoMzcYkz+aPbYNQRixz1hzR1A/PprteqTIsWayzEOV7pAlifQoo5ll2JpTrG55s5nl5tcLq/wd6sQ1FodF2BTVNi2cxwVs/xw3qGPy1LfL80eNmscSHOvHaDdHaRdmNk3O/KlA2CrTFGh20+hvQ54ffh8TsW5ur/brY/g7CBCS7Wu7B7ADgAntY4lGQnX2nwU/Oy4qmvzIOAHsXVp6U4DW3/TTNChjwAaeIkmti1S2suTMg/oUwHafKpp1NLQ6oMJArHwtgEgF3kK840ALidVGMWaLDARlzgTD7Hi+Ypvto8hIRAgxYrucCFeIFNd4ZGLaOEV3RNBSrMzB4ECnSiRYegxaBjSLzqdwbJ7p3uJT9OoEBrVjjDcyzlOV53R/ju/BGeL0/w6cEevjjYw+dnp3j4osVh2UEZiY2SWCuJZVfgRV3gz1cG3yzXeCG/w1qd2SgQaA8iQsKsgUgENwmTm45P6N4X78aaYpoKAlfWLRMq+/bGUwK4eos/c6lwQbB1qttnHWkWOItkAksWThKDjdS0A+9WT9IvBPgxGE2ONs/KzC3kw6xkvo28cDNl67iFiBZGjWdB//ACwHed1xGw5kwmINKeo7+btfoBiihTaLDCRi7wBvvYa4+w3+xDQqCDQi1qrOQlaizQmhU6Xfvxb/NX5O9tGrUypJ3r/e7c4AYKHTbYGMBIjcZssFo+wJvNPt5sSjzb38PT+QyfVEc4KhUKYdBoiZWSuGglzluJ1xvg2UrhVbPGC/kSS5yh1gsoU0dgY2jeHIyw2cJuZBm0DDhN28rORwPjd5tAdZc28scGvZjWrtjlHbGd9wRwwC1elkZVuraFzBitzK0/4QdQkHWRgJWlTyyXnZTa5PVCqE3a5QphAQO9SME/ZoGPYmAl9DsMWt8fV+W0BKBMoOria0rrhFCGQVeS3dWW6bBGI5eo5SUu5JEfwJ3eoFFLKF17sMF9fLSQKsds5EwIlxlTzJ0Q015XRyJR9rJ5Jb0gIWlrc1RggZU4w1KeYd19gbcXpziv9/FsXuC4KjCTQGeAWgEbBaw6g8umw+t2g9fiDdbmzE+O/rl6cBgyh9L1RYsMQoRHKLzHnjUCwATgE4HxiCQlagvS9PhOmd+P+B4GMBIKvRWO7hyYCDPumhTQROc1dqwr1XxcobCRBcaTXCt2kx02BGPh5dsiRPhxY6Ajbav//WGqm7vgxo5NF6d0THE9F/0en4z93YEB+07aCsQbcYEV9mxEloty06ZDZzZZVsCHdrLh0hdXEjM60q/0GtymSsFGXmmjoESLTnRYqxPUi2OcNxXe7BU4qST2ihKFABoNrDpg2RksWo3zpsNrtcKFOMMSZ2jMwtfi2gY2xmxIL2Gvty8m5RtQ/24mm+AxPcjQZncIdAxpQrLXkuhEglv/7tnOewI4gltFmQad2gR6En3AASABB7oHDtKoFL7A8Hh2/jtv107cgXEhgacXHaKCBgsDdRN9J2NESWwNr6nigZRm/n+houqjQ2DD/s2qogXRqMznbdQSrVhiI88jEJUuxhQGaqDRYYMVLqCFfWkVCz+mEE3LWthFuBA2jffGWNDE27btx/VGbNr6AAgavUBbrHAhH+BV/QhHmyPMUaISBVqjsEGDBi0asUEt1mjkGrVZoFGLKDLFTjKIisnRs6drp0WLi3ztIh8v5LxQG7lcKnngU493QlqwY3Qm4klHY7HvqgsgIXUHStgEcpy251S13bkWPtGidLqT1I9sjA09BAClG5cL5GMJhe2bnQ8aaF1CyxmgE51AskDkJuqhhSfnIsml7M4Bl9xiMrRI5ADENlEid50RaBkqEkf3QEF7IM03DLQZkFh4Vo/6qdh4jzsoo/kltVzZAC8qHRhnXgNjrAbP3he3CRErrOUlFjjBQX2E4/oAe6JCJSSEAFqt0RiFGh1WYoUaG6zkJVqzRqNXPYE73RfqF+8n/6x32VvCS3suri3aFt6PtJ1dPtvFRjc8rF+W7awZ23l3m497AjgALx7VDTphQ7MKOd9KHaU7hFBxb6RQD3uZ+Pf5S8OBTvSCi9CuROUWi4J9ZxxsUHu+OJcTAhqpmQBxOF47+jx5CWwIbgeNDkL366NY7UDpBYa2bYVO1z6JFV2HPd4Kd604s4IUoYy6p/UzGoZI+AibA4OzKxot1kajlWts5CXO5Bwl5pAooESLxqyg0aLTNbs2Fz3EGAHuhhjfyYbU9lKWEEb7bLPRLoAtJsIpW/z9TdPhs37lWLgc6KD7Pbbz5eMk+j1DY/cn2MD2UO6Nj0csmpqGMTRZNoC07xtZzPrEiZoAZMeFfwYDu9le+G0KOtgCsssulvqxi/H2cqCDt8UXVFuNunBMaXxe7262X0r6lRcZxhFyQ6zbODhLv+u1Up7pCNonLVovfl/KOS5wgApzlKaEMAVaUUOJ1lZXRh2yBSfsDL8325iA3DuXAye7Wk5EOnjuzNiL3HB8jCUMnL3PiRxghOmk9uzGudt5LN6m3SPAYdxAbGxlUb/b6ydgomiVwFSwZEdGJ5NRKtCJdyb0GfnqifqTooSR2uXQoIcXcmqEhSFO2mWPC4sQ9Z8Wei44DLthKz41QvsFkYAHBy/8GoLg0+bDkKKEFSCyOiAussWe2/XDRXXQ95XLomq0TdFN6cANlI+4sXqJ0jM8EkWvEiO/n0AQvNKkoHQTPx9tF8ZWr3tCVMrUSeeliJJK7kObDi0QpY/nz9w/e0qY5RZg+ywkjKJ7LAPAckJRY2LAp6EgTSisZu99P9SMnic3Ol+UHdctIMrUPopmCDCFtMxFlOacxhcHS7wv2rTo1Mb5au9OHPbuZuB1XW7sWJwb38/UzeIFjugvANvcYNy4wv+6u9id2h/YDKUCee7eiIBGunhk8kBsc+dQP/JJqfqfpULq9PwwYZGL2nLjFUBv3rCsbg0lajSQ2IhLPwdRP7xGK1mMgZBxN3fe9Fj7c+bYJLood61Dz2ss/DkHbPj94Pc+d230cxrtRO/+mJbM6xB1h06vofQGxjS4+XzAY9NvvoG5R4DDmjHKTjAS0LqElPED9wuD0T6ElkLnpCxt1AgCik53jMAArSXYzh2OxhalDXHQjHJnk5wUQXTZby/epfJicMa1n+5M6VyUwyPcEx25EoSxixCNAR8SyUS3uf5EgkU6b8TkKLe42+8UovKMT+ESkBVIw4tjIEQLJfU5dXXx72jYbHcUDeOvlS2uKIACc8+yCMeCGBFXqk3Nsx/peV2UR/p5b3GCS2uPNnIH0RjhICW3sPHzcbcbZ0xSy7Fb/G9p/1LzdWJMA2PuVhz2foy0HI11qUD6+cAfwViOIeZiG0ua2tCOdZuobxdqfIyp7TEoHHTQ54mI3R3MG4psyK3EvzukGxlafG9ifXcXz1fB36H1MDNNfci974kb5UZ9Gwhrpn6kgCL0Z3vNHv4ZH18c5FzHIpdeZm6geVTpGlo3dk01GjebD2hQuQXqHeyeAQ67q7GJfySUm2hkkVnkqNggrHK6kHOfF4G7VsiGVMg+EgE2gRcVhzOighKlpXRFDaFlvOizBbwvDg0uB77I+ERVThwKbd0LvQWKjmUaDP5SWn++ZSooD4X9vvKLNR/cnN0oxQyUDIu+47oCiNLXrKEok5K0BqJEhT3nM4535737arR3oaQLcgoAe24r9uyIxfLF5VBBCG01JEL7OjbRri4BAcQGkb6BzkF5PPw5k8lAmRZGBbrWgynGnETfS3be6TO1/RwGEylLxnduAeT0aWp/3W6CocnFPtePld0IZqAA3UBLacemofoW8W4RCKCDjLMd3MbCOXtusLR4WcIspW0OWf5ccVGuUdBBn7F3zi9YyY55zNL3Ivd3rhsh4HFdV8OYiyLdwQOI52szzA4MRdHQvdsmrhxiIbguaIgJ2wreRs6TtuPH1Mj6z92wQ20PAViag5WxGxC4RHrXNwFsGdvXsXsGOKyRgNSYDnpsAFjS3i8EUlZBwGjiiXzMz8YzlQoUfldKlV2NUBG4C22SEKs/IIYGru+LdINFxaCjdxzbHacUOhkxFoWcQxcdoOAZgOzuxYT2eDsaThTmwVzhr8WWYLfnUYD3pWaFWaAQrDDIQ8I22fvXL95RWvGge7GsjHXl6OTeRAsP0HsRqX2JKjrHEOvi74XpoNiLLlFByAoFil6q994OI/VnU9uJ+p+75nJjNB+hxNxpTodDLihNbMAd1Ei4PQuuVmNKaCOjyCF7RAw6+L9DO9MYWMSLYhR9Qc+MFc+67m0dAhvpuccWqKGFMgc6eucaYGZSRogvuNTvNEom27/MOB08dgc3T3ZRT57vmF0noiP0K7yjYSPC3BkDQCN3bwP4GXexpOfL9WmQERth62huC2GwNA/f/XxwDwEH+W4BrRsIaZkFiirgkStatxYIOJaCtBeUTImO5xU7/a5YxKWejdHIJWGSju0gVwQNgjTNMFn0NxFeRvqsYG4YYexn2oSqoL6Qm99RM8EpAx/et+zdKoVjAkqgAJRq3I6PUnGHRdIv8gmAkYAvv+4aRWGqaHDbu1SjM7UL/dW9vwfaNwY9uZ18Tr8Q6V5IT2E0tEtXLkVlC8mxyqmWxYBN+81cTULY50Bsi9Zw/YqBEB8v/J7Q59Re9K+ghb+MxqVE5dgW1ZtwhOBAKgCNsICEZxp2f31gQ2PavwsfTQn6m5itsWLBVJgPgARwZkCHYAuJ/SefUGvMskzHUE9HFtL4uJxGYhjM7N5uX6zM2ZOhnTfXEnDmDwi5grblJImZzL6bYaf+5+aDiH3oJ/zj10htbAUmA0wVfT+/mdle62T3+z1c8XWMHRkzIWSY5/h88N5cq+/OmN5DwAEQ6OCx+DyhFl8wKSICfFFDf4GL/WZxJj9qJ6KxRciVQCAhV9ish24TYMFDXQmMABTvXkJI6bKU6t4xQJggbXhqvIBbkKQh3SJrQZFjZWQQjtLCrJwwNe4vBwuh8BnrgJ9srVi0dcmCVt5tIkHpxccRefpsIjqTuSai5+YnkbgoE9dS8Oii+BwxKLTXl0SIiPicuxofVyTEJV2H9ung84mqwnfiZ53uWv3PkNGklJa8t75aok41fjrsRjCqsaJ1Ay1KyAQIp5ZO+ClI3NVFEI0tRttn+2j6Kcx3HVPcvWLbGp/c/XXt0H5uEYz7HbQEObeUMTqUX8/Yrpk9+fmyn+/w3Z0YjgFGh7fxLsBw6/np/pJe7xoAIsf8jo7zTMhvYDt5Rdibmkm+/27zyj0FHNZ4LL6CdOXEyyx12Cu45agqXqPE/y09j7E5MLjgUsrKMyf0UEtUniHh4MS3zRYh+k7B+pYrACWEtMm0XL/4sUIU9vty5qvK5qj/DvBizoItZCE/RdCfpItZj6FwrJEx2uYGccyODVGz4cqNWti6JTwpE0V9JP7X8Zcl7Bz57o4zAPR8fNIgwMfe86ygoc1+OW+d+MIpCoeAHM+rkdsd2nthdw3EUFB0ED+XBp+8c647mfzHABHXfZBuBPG4pe9w0a/WFHL8sYfBbjPj3h3KolpHk/OYXXfHOJbEapD+ZvqB6zIo2fMMtDEEHHK6lDS6YZd7kLpIOdAduy/b2r7Jjp3bGDMVsYjpOpA7nu7dAGt1HWFsdL4Mq7MN/PDz8XWDP7M+a9UX8XKXchDrvy9XyvubT+4x4CCWQ3gqtdDzKOlRL0QRIS+Dd6WwB5VbCAk48DwR7g8QkjEnCKJSDYDXFUhfptTnSeAhF8pFok8f8cKyXRJwonwkOWqPkl8Zo31bhXOvKCm9kBa6Dzb8pZpQFM+IkJXVL7DJ4k0TPmk0LKDpIBEzEX5hT/Qb3OjFT3dsXiPDn59zUdi6Dza7Kgl/c0CDs0mU0dOfl7lt+AKeMxvL3kJJCUE+bREL6+xxsjexR8wV08T0WDZ37VEItIiBBv1cCFuHpTW1TfLlwMbHmuRrd1OgXD3kavXJ+JIFqc887raI7AKWs5+PUPS7Wi40NXfuFIyO9ittdyA6KrVclMYuli6AQ+fKAvLM4swZy7h/40xz7ryBRZUBaCTPbQxsDGlD0mtOxcB0ft7vnFsoai8Jx+aufDomasOESErSblh2435tPu4x4AAsx+FU6pDoxBoAUBSzsEBHPveYdvfuFgS3CRlNVN7flbhT0l3oMLpmfnbjBKDQoBL0tBilg4j3RzDKlsz7+qMdsPI/8xwX5PfX0Cgd02FzTRR+EbN9jRfCMdVyj953izZdc0646RN76X7Ss0LOBndfqXGXCkWoUPs8iZg/j3TRRk5EycEGkLBfjDlJ3SHh2jPZXaH9AkflyndN7ZwyG/4cxMLB6o+4oDW9V1Sfhp9P6QadWlthmG7e027mflvI1ePGnYSNZBtgPv33ehN/JvdExpVA33kfrAW3myWYCtqUm4RS2u9uj3DYtX8mM2/57291a4wAjwHmJqeby+WvyGpjkGe6OLO6i/Xa3oENGboXKZjKXV+2z2xTBbgNrPuPwmCNaXAftVz3HHAAFnQYpAnBUABFb0fAXA5UcXUghJPCq3gYZUCRzp3BKoqSv5+MKtnaHsZUOEwHjdZWBIUMbg5Iz46QsFKZOlooff8gPTMg6bxsJy4B+FgEDzpKv+OnnBX2nNK5SUIMP/lsuQ3tEHv0PjEHAym+CcRRe9YVxu7ljr5aP8m4l47ARqc2cb4RDYgiFpQFQWVcEyNmUPp5L8Z8puRaQZLVMe03cuHAiAEkaYl8tWDY3CxIwCExMIVgpedNcCtQQp/7OLncjhkYqlek3bOUef1QTswY64aGRYBji+ltWx6EJzU3xPCCno75flu7Af9s2+8BRKSfx+/JAMuRC69NNpH+84wANn/uPhNxHRuOLuHjKg/whn72vyep09PNb8qiG6OgVOPKGmzc924GSm/TPgLAAfCoFUhA6bCAc7U+kAz2JB8Cf/BDxdzCV/sRBL2aIY61kNQ+YsCjdYtCzt1AKeIMesKCDY9KHUDxCyHR9UKC4qA1bztTvyT0OwlRZVRf9OLqwHSkYIPfN7pnErE7i989Dshykz+/l6OIn1gHxCXkfahtZjGncxroniiWX0POvJuL30/DXnDBMwnmKeZ0cQrPsO9Xj6MA+izH2ARJ19jqNTq1dsDrpysUHTar59C6gYCEovtc9Herg26yHRaYm7IQ+UVwgEq/DghPx9NIgbVd+8p/3prDI/fOZtjfyO2TuAZyIscc6OBt5T7P9i9yGQ8V5hufg3axXUFo6sq68flG1inKzKx04zOKvrtQ9PbsIwEcFCrbQmsNxQaLEGFQp7vynAuBL/j0wu6krIb20RuRoBPa0uFeTMoWRLYLpiydWihIt7vSovVuAcAuerTgSklgpgJVeRTCXkOn28BWRG4iRvE7liEtPe0znrpDte7CfUjvV7qos9LsvYmGhQDbdpF12fDJJZc1UYrShr2Kki3cBQBXN8EV1OPt8Wdk76vuTQnpefhnqfuIwAYlhdMJiOFVZeOLY3QnVZVFnPuDjqN7GK7DsRyirwMJ523RmQZNe4VWLaH06mfjSumbgjGA4vsLuu9DoaVDugLsHvGRbZdpA3gUy7acE6koODc+R339iQ8/1RAMnTvVAtiv6+R7u+Wg8Cwm4uJt8c5+QPsQgfl3c+9sa2Os7euwHNlNxIC+jp/D/hHRmpBzp2zTcpDxbM6d2rgNyApGb+71fPCRAA4ggA74VMfKuSokc1n0adV+UikOCPgA8jkw+K7DhIXMAwMWIsYFhNZdwqqycv0G6PjCuT0qSNFBiwBQjNFRynIp4u/R60CUm49mIA0FNARCGnTtk3NppysBKPIGlLYdiECPbx/9nQfdO05XEr1HL5Vi0TAcyPA8CeE8Re+8kXaGMULRcxJMp8LvEQlw2X330Siyip5hzoZ89ZGv1MRjIQJLvg3L/ggpodwuPO+HLfx1wtjz0ChN+0Lp1Tu1DrURfrZgg8xtQkwJmAZCBxfeNqaIjFwTqWbDH+vGbW5e2WbvorUYaz+l0zV0MifkNUgAkuscX8jHQkLjOS0POvz5GAjrtfMO9+h9fH+43d3Ay1awkV5/bw68HkhK1yTScXV6Betavd/zwUcEOAAPOnQDLYHOTTBlsY9C2Poju0wMQ7QmX+wARPkooh2yAZD4FXnRN2I5cuei80lRoZROsJpxA+T6K40lyqweoPDjSpvORmC4U/FoEouCY62JEAVKVNCQMMKxMrTUDdSf4cbFq3Qt0lWTLcwcjeuTUvYFMMa6sLiOY8xXS8nY/L0VlNStslVDPdAMwt+0jaDVKHzOE0M7Bj5BZHYd6bV7oMGAZ9btJC3TJUUJuLT49HxC34JQlVxHvK5OdB2g6Kawk2nV0u5kPoLJ5faN5oONpZb5mOGhziNjjgB7uiBzrQ8dt4tFInbGwl5XI5DTWZBR9s8oCivDQuRA8ZBLjy+MqYB6SCiZe0/Sa9jFfZFzWw6Bphwr5HMP7ehGSPs1FM48pBvJbSCGcu5QP/15RS57cL9mS85CwIHT/zmw0arFvWc2yD4ywAEAceSKYgtE6jdPJ46xh8lFOf5MGT89P9a3u2WgkzuGStlztkW6uiA8OyqZNhoyqaEhXZ0TqvnCz0EhsnQ0dzOR4NGGr0q/yPbEjQS6/Ms4lFMin2k1PqYvTE372/vcuURgqCKwhk7oYh61YUWxrZ+EQztxn4mJIlOmC0AjSTcOhEmCgw0Kw+Xt0/3j3/OVh2FDk4UJkUIcbEg3dnnG1Nz9MLA7GetGcSLR95Y98GM3EpWHIm9kQ8D5JgBgm6UUOxBT4jkAEYGcHRfLXV0/o7qlhLXZBgai37e88/3vb7/XW7UjDIzkQEfu5xhYbb9nN9HrpMzGoOsscX1xgJvtOwOKnE2le2nnIst2BtH4xzEffISAw1qIXAm7zULOvao/9VECebYhuAb6+SY4q8H/xpOJEcWuDCsLn57DKHRmFX1G/StECVFIXw+D/w0IUTaABRuAS/IlACVsmXNa9JRpshMf/1frFp0EjA4o2V5f4e8dRc2kGgt/D9j1SzdZKlNDGQcYmFuJEmulUSr8uFjs5URQyb2ivtD5iQUhtwvPgkiLNF/AJWJAKYTNqEr6G/peLmqFV73NGXcFFU6DYrVFNjRZGQlh4vHId6j8HGkeDluQLUSkaK9Av5+isLsx5cY2IraQxl4EPGicDEQZbYtsSENngQxjysay/b3Ijq1cm/z4bZYeky5QWzdCA66VbW4UvssXYrd6Kztfw4DgNtq8JBqoQbA3Ypzl8N/dFfRlXCjbQJW/P6KfiHDoHOn8S/M8zQedXrFCjfcfbAAfLeAIeg7ARq/wmA0hY3dHT5dxozMOlL62zovejp2L/3yFUTDxJkIW0wIFIGMxovZajAoQdrevo/PmJyRe2p0zP/zfkIE06D2kC7uUsvI1YyD6uwMbLlwipgHD4s6FrORLJzcC//vQAk7onY7t9T9xgWTFlRmwGOpBBCai8Inc+iJf31fmmoqYjJRB80AsaE+IhfJFALNg1MXP82RqHCBSyXndJGDj45hgPpxRJtImCEklbNkA0ddh7MoqbAMg3NJswfw8fSA7Ii5EHqTwvvBFPsuM7cC6Rq6V9PwD1zzWbk5bsiuTlNXPcIY51UANhL0OMSpjLm7/3R3XiD5Iut6awqtOZ9sfcd36DZmhqtAt7mP465B9pIADiEFHA61zpdeCpa6VXmvOnwkDKMQRDQEEAJRIDCZMCtQuRLwbD9/XkXjT6y1EiRJpcTTaoYeS6vacrV/EqL1Wrz1Dwa+z5yKJ3E6ufzr0CbAgTcrKl683UGj1OrSD8PdSHnhBqzIt1EC0Ser6MNBozTpK2pXS3jxpGJ8keGZQHqVi71n8PCkLqb8fIFdW2I15tkMAhasRk3Ohad3G/YjAVxhTxG7Q+OKjMdoJuXaV6UBZSwlc8nECQ+Omc7QpgdEJbAxbmA9U8ornoleyeouI2h7P58A1Xbmdfm4XHX1vxPr5Y/p+/gjkZ3bEY5a9HhOKTFL/qf0x8yLtAf3FdTZ5Q27spEEYoz0Tex0LbMntMYQ3EoLy343y7z6BjDDvu3pChsDGx8V0fsSAA8gxHQqZBUHEL0QuFpx+jkqkJ+DB/yvdJGI6P3hjjUcRvbj0XY9WfaXXvrDNHytc3IvXMki/iFHEgt35dtmFsK/LiOt2eLYAYSGVzB2gaLI1OgIbhZz7EF+FFoJ2bzw6h8RsokQp5h4wKKpCi/wk5ilDBs4AuErBtgqr1p1Na40AKNL2iMGhayQ6lgOPkvUxBof5UN7sc/LuuJhRS0Oo+ecRLewSiHGqNJpooEPWwAls7GDJJsTIrKYjt9jm6Pkxv35f5BiDAH7O3Xf5I26XBPykC2duR5xez3WM3MdASHKXgqehe3XdyJFtLqDe/RO0cQybxG3XmWo6bmJDwtromJzb3rHcQ/2Kj42fo58L/bzQJbk2Pq754CMHHEA6ySjNdivk33ehptmXhVP9iVuEI02i+YWQfjHrxU0PZIJTDJiklKF2zAWBCBImat3ahVVQBIzybXl/v6tlQi9cbvftf0fsx5YuLDboIoKIkUfBkPaiEKUHGxb8BHBGfwdiGjWlWNP7n1KFQ+ZDwYS9VqpQy3Op0HFpPpb0fvOdEanbIQBpZPTqc3A66FdOrk+ZDjzlvF0gYt0KjTPpzgsKcYPdxYQdjJ1QJmbjusbC512iQHqwBKa5i28bwBhiAuy/tFAq/76kbfLfI33BQH6FIbp+aCGz5+ciw/zY38UtlM5fNPaJhePvLNdw3URw2buW67pfqIo10AMdg26ZHd1jqbstK0z1IEz15rLo/IyFSgMaQh+TlA3kEia9hkvu5wuyeWbj45sPfgKAA+i7V2zab7tjtwmYpCh94ir+QgH2gQ/51EhT4Gl8WtwdkKGQMm8CKFBGoVrCLZgeIbNFSLkSwlyYaBfF0h+b7ppSmm1n3yt7KaWsrH/bvSCFnIcFOGmLhJCl2PPXYBfHWANioCP3BOUFsd/hRd3yE4LtJ0sElCRmC+xQCC9LJ+Jowk2ecXptlPnVQEKLrvf+ppMLbyc3aRD4464zYj9yAjBq2yc005wutcXYjBvbH+PkcndGhR87aA0POozQvsYQJBjo6AvMt54hWWB4iOkueoHrnCP83nfxpELRnGUZnR4oCgUWY7dQH3SQ3VQoOmZ9rcb7PUc6Bw2BkHSu7vVjJHNoNCbiRqPxkXt2ntkgobjXath1ws4HH49mI7WfCOAA+qBDw4gSWpSQwrITBWwkC49KAeB2uABV+wwt6qjCaO+MbNHwXQCiTJfpgEp9xLxwHN/p25+V7xN/4aWs/ATKa8X4axEx4LF9CAsiENw+PokW65el8ZIc/tF9UZ6R4dfFdwZUgdYYjUIoP3EFF8Qw1ZwKM/lay9O2D36XMT5pDDwQ8ngIFDbiBxWUaD0wytnQLsTdkPyuSqAH4vhixNkzCnHThgqxEdBwJ/hZG9/O7mo2fF4Y+OgVI0poE6oDQwJUZHE4MdUwOOafp6Bj9PidNBy77faJLXENOwDb/+6oWHJHwNVjLEfcP0OfDbuzx6/3OsAjmstGxJlA/5o5G5X+jYt0U1YDiMFMDjh5F5CzoSgUYzQ6vQ6bj9588PHaTwhwABHoMBpGdBCihCG3CS0mbhwV6Y6e0YQUupTz+5KA0Z+VAw9H+RtXQpwiWKj9dJDmhEH2OKcNYenEQ2iVc4dAZ8cgRZIAId8EzyGRhgGnReoo8iZnGm3QkDh2I72P3OcoIG01VxeVk7Pgmik8oPAvPrR1pUTXF0cdhfBYmyAsTaTGxZ3RfQTb5RoWXZKwGLlQwxQg5ibElK5PdzIAXMElG/JKvtkAen7uQIPspvfB5ewxNk+HzUwczwfEengXF67HdAD9Zz9G3edCKsfazIKHLGNRjI5FOiY99y7ZRsfaGjp2KEKlF3nDFtz+ORIN3g7PhQtYr6Mj4eehn4fuZy6jahRkMMDkEuhIGU++JmjTOnZjw6JQfjoM508McAAedABuotF+AaLdZCEbSFmikLMopwSQpOkVdhEs5MxqBsYWl4R6i1J8D/Y0XrTStoMmQfqS6L6PorT6jwxTwF0Jgs4x8FKnWT25jiVMejahWIeNv4/xItwv/kZASghpq7k6PUrOUpCRgrx+obr8xCMhIWQVsS+pK8n30VggQ4uDFPky8Ybdv1RzwoFSrl9+TBA9S6DIBDcdVXf8KU4u98OS+cCBaaM1hClRiBmElj6HDxcBc3Emt215K3bVS/DFZ8jSSIysDo2LvxFE7HxOSkWdY0bRNmP9Gup7biEey7XhXafZEPzw/oyBqD5LsYUF3YHJSd1H3FIwlbpJ+WaJn5ODjt5cYlgIfBTy+tOaD36CgAMIDyhMNFp1dpEWdoIXwgKOQs4toHDiyNzCLmRfRJR7kaLsoynNxtLnpgtYDsHvQiFGO30To+rQx7ifcYhsP7rC59Bg9CHpTAoxkhGTtBrkvuHuIRG/nDm6MeoTZJQPRGK4Dgrd01xCNlJ4WxZn7r8Tas44IEY7WxdtwzNB8pwl4Zxp2G4RXUPu3vjv9cDGKqFMfzqTy/0xuqeh4KEQHQAJLRpbjsB0kLoMwIPpb2KXXpxX47psCLexhXGs/SFRqkQV5qkoAq1fqG3I6H1J3/EcOzkGNtJNztCmavS6GPjj/c9dw7b5sicE3eIqy7EYPAzZ/tuPYByy9LmkbvT+fPChk3ll0i7fgv1EAQdZH3hYalVCiBmMy+yodO2Bh/X1211FbxFMXjqFLhrEvrroAKuQ7giihTZB8T0EbRAxBHZiKKzGRIQ+UhvatJE7h5+Dh/Jxn2KvD9CQLi+JvYva3QWWA8DtGCmSJkelRiCHfk5U5UKEGiMU9UL5QMgt1BfRZXafYFVXPfNQsN2Fc2U47UsDO7l0uo5q50TPKKFG0xwhIbttAv6YcIxcTFRKmiJRpgJsH9qs6M4+Iwk4QK1Ntx14JDYUdcU/i85s8nVFhlgO/t7kFsBeoi2Suohh1mUUbPhNQ/+92pZoLGonAufDgla+CUuHP2eZ6PeksyPnj9029Gy2gcMhsJUeY09/Mz3FoAuFWI07mw8+zPl+4oCDrA88jFnDiALaNBCihNKNZzyk7HrUapSyGm6wIBY00i4bpkWHkDU0FajyLtmFtvLtE8tAqdLDFeR9ornkNz76BeGFDy9umh9EMRdGkg2RoAbVAWGumgIFOpdXwxgdxerzfz2YcN9VLgyYitz5CZWXpU9zhgxcu38+7F/+LHzBNoRJJ2UbOtR2kmGhrGnUkUiAWbzjDexGLi36sOp8yq1xt8aBB2BM6+cDKWaWWTNlSI/OFuG+S2V4IcsJDcd26tR+LpquwHBhM38uETY8NzEOOrZdW/b7AwvyYDTIAMtDc22OWUlFm/l+7JYnI7eJ2WbXARu59vmmRbkIFJvc76fnQkntZwI4yAK1Cgg7GIwVl9owy87vcPyAl2UcPpdGJCBO/uQX7gxi9+wHeOKqwoMO6VNtM+3FSDIgPikBGf/pCI3po0mgwQulpbv0VGBG54z+zoIIohwojLUgIyAAaSlg/7nLSGorqJJLQ8WLuwi7r96kzCYm/+xMiOKJ6d646iLvt0QMkjS0d6lZLUtn+5ewTdHOBTFVal0onRtfdoKZwMZ9sdjdAqOhJSKNRwQ8GHBPF+feQuRYyfcdPpqKP9PIDw48wuf5RTLnFk6TZPG5Y2yOGWt7V9dTGp7bB/FDUWQxO5EyQmPAKQcgtglrhxKtUV/4v/QzzxZq54cwL3xsKcpvaj8zwMHNuP8rwCjAaCghPeMhXA4PolcLaQdP4RbEoFMIoXD2vwJAXAsEAKSYuaJjgBEanQ5ghRJuSSFdmnDWNgJFmAvh4os0X0ABG1HCNSQAo1bpHEJBiwLCuORdbP1LX9gQDZIuuHG4FwE1Hj0S7gPdg5iZCSnQA9jgVWLpGni6dH49EbshJAQqGKGgDIXWttAUrZMBBuHLHFi5OjcuV4ZCEMJy0KI1EKUyh/IgwwIMO7kA9EwnoHE/zUa1GF1DoLGLtphBmpkdcyaIzAlMp2HUETPAssm+U69GNh2p8eM4I5KGhw4v3CHXBxD0aKmrYRDAXAOQDPWbg7rx76S5cMLPCn23UApCYsYzTebG5p2Mi4mDjSEdRyQ+90CDQIZN5vVz0279jAEHWR948HBa7fy6lJhLMV99FIkhAR0t5qEaqkTpso323RkQzj/qXDTkqulnDSx6FKNfGFH4RRaCFmvttR+pLzT7kgp4N8zwnXIuJBG/cGl/yA1BzM2gmDIzkfpS7bQjYOcgAFOKuWeD6Bp8CC40gAqFKGwbiZaG9BwpC+GNtF1CgeM7Oi6XMwVMy+PPxyaYkLDn5zOxfLzmXC0AYIR9jqAcHp1PJkjAe0gsLCAhjM1quwsL0DsmBStmuCJrjulI3YJcizXclzh0E4hdLJz93MYK7Hqe1IbydAB90Wav3QzY4fNeLrqE3/9UXxH6FGvO+N+5i4T/Htofmg9+nrqtCXB448DDIFayS3TKpUMWIS0yWSGtz9dIBSksOOnUxi9k2i1I2rQgoRSBFUoxTmAlpC0PoZ30chZyPzp3GNgKcAxLIewE14kNFLk5mE6D6NZUy0DtcHCSU5lL2JdPQ0d9lAiuE4oy4eLPcJcTjYrLwFoIBS266Nx0HN030rlU8gAAoAygjBV70ottQUkd3FcZX7TXczCXR3AxdRCadrCb6Pp7lKibSMJExScr+tsEND5e6zMeygs54/kgDkENmYiN3kETkCyggqBFwkpo9IWqQHDZRPV7klIO8fmGq53mQAf/W7pIR3/nwH5Ap7I1KRq5TF117KGQfaDPgPSSLArp57Scizi+Ju1LRXAwwYXttv9FBC54G/5nN7fRPPDh54MPE3FyE5sAR88odazTeEDCGEBAQLmJRrt/gXjxMkZDys4LBHn+CCocFu8SbOIhA+00DSGSoUf1Zwq+0QtIGT1LCUhTQBPDwdwzGi0kKibuVHH/kt2AzwHCXibK+0HHW8ET0a2dTSGfiFgpgoV+9ou8y4lC7SlRepDCr82zF45FKgF0ZgMqBtfp2rcXbhVvI9SKCfddQ6kmTifvJgTQbgYSQoRKvNGuLgEZ/QnFXm0YT5N9vEaMh7AMKABjZDIfEOMhIwBC7/xuZ2HJ79xYH8pyzBc/AdljLqj+Evi7zlwHW9kKgwh0xBuPWLswBBwGNROJ9on/zCPZYLoeYOKLumUdugwg0fwL/rj0nNH1E3vJNhP+fCbvtu25XiLQ0WXmhA+16RDJz/dr/hHGmJ16JMTPEZvknK/0UlDyBuk1H6T74Cg3tGQnI/qZi02lCGF4HGmH79rdUlnsR6JNvigLIX1xNcBODJ1p0Kl1D6lzgMSBUBB5BpcR5YzITSz85c25T7jug/rLc1BE3/fCW6pvolgf+b2qUBQzf34+WaR9ywEO+pvyPtSYndhm8XmGQAZw3170yd630dwgo7kAoPnA/TvAMqTibPqMj/OhRZu/y2Pt++RltGFg2Yb5uXsMi38P+wt+Oj9tAxzbANcQo8Db6H0nAQ1c95W7lkE3CWvPizcj0JJemxz4PMdq8uNuey4YEgl9uDnIbHHHk02A40bGH7CbcNwEQ7H9uYEfm6NDGS2b7sLD2ShaJgjWUtAAwBVgC8+JsyWpcYSeTnIc+HhwMLIgCwSgwidSLzD1/XHMhldmpxOHBWVRroLMvRAMFKXXwr8TgEQCOKCj84N9b1ebmIzJYksAiB+juVDQAZDAxr9Mxnju+70yACl4YXMF/T1iFdimJgb/gYlNmQASRHO2YaiPOZF7er5I55DeC/a9uN86OiZ3/dxM7v3m99YBDuTAxo5zwt3VPgrjrg+EJsDxE7R0opHDgzQ7gdDEEO+GIrEko2uzvlE6X+7vyWI8GJ7GWQbwLKUp6k++xsECXeO1Xux+W/7wkZ1GfjJPNBTunOgdy+9LbvjT+XIvcWoT0JgsNYGI+Yj+NPD+RYtnv/MH2QAAA1hJREFUnmFI27lJltNs3o8BYDDEcHDt0pjt0r88ozD0vWFGZtv3huYv3+7gfLD7JuTDzwXp+Lq7uWgCHB/cxpDm8IvnQUrvuNxAlz20v7tpINkR5PoY7xbSFzy04fu/hfbM9uOdmATb15jGpoP5LmhsWN/k/tlWJ5tsN9slFjZ5l5K5IMcIRMcPgpfrWN6FkQMLeXF03nbp2/CGIndwDhDE89HN5oS7YwXe3e6POHQCHHdm1w26TzQhzoYWTAGxZTEdt+xLOWZj4GCsjUF687p9z+kiOLjb9p0hu/uXdLKfu+2gEWOWvjtZBmXkO4PHR+/xtnlhB1fDLnNLqrna6X0c0kjFLu7t3821Mdm72AQ47r0NTzZ51P2O2YOy57mObaM7r9vOdWxoiO5yT6ZJZbKPxa47nm8iFhz6znXf65hxfDfbdU7Y5V3edg+n+eA2bAIcH5V9SJXx+wQu72rbJs9dr//++DInm+zjtfswN7yPd3eaDz607Qo4JhRxL+xDvhDv81y77sh2mQDepV/ThDLZZO9uHyIZ1Yc47zQf3FebAMdk72C7vtjTBDDZZD9vm+aAyd7d+TbZZJNNNtlkk0221SbAMdlkk0022WST3bpNLpXJdrBd470nsdZkk0022WR5mxiOyXawCThMNtlkk032bjYBjsl2tOuyG5NNNtlkk00WbAIck70nM8nPEysy2WSTTTZZsEnDMdl7tAlkTDbZZJNNlreJ4Zhssskmm2yyyW7dJsAx2WSTTTbZZJPduk2AY7LJJptssskmu3WbAMdkk0022WSTTXbrNgGOySabbLLJJpvs1m0CHJNNNtlkk0022a3bBDgmm2yyySabbLJbtwlwTDbZZJNNNtlkt24T4Jhssskmm2yyyW7dJsAx2WSTTTbZZJPduk2AY7LJJptssskmu3WbAMdkk0022WSTTXbrNgGOySabbLLJJpvs1k0YY6YSn5NNNtlkk0022a3axHBMNtlkk0022WS3bhPgmGyyySabbLLJbt0mwDHZZJNNNtlkk926TYBjsskmm2yyySa7dZsAx2STTTbZZJNNdus2AY7JJptssskmm+zWbQIck0022WSTTTbZrdsEOCabbLLJJptsslu3CXBMNtlkk0022WS3bv8//rIlYQLi4VgAAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fig, ax = plt.subplots(1,2,figsize=(5.5,3), gridspec_kw={'wspace': 0.0}) \n", "plt.subplot(121)\n", "plt.imshow(recon_OSEM[:,:,20].T.cpu(), cmap='magma', interpolation='gaussian', vmax=0.5)\n", "plt.axis('off')\n", "plt.text(0.03, 0.97, 'OSEM', color='white', fontsize=15, fontweight='bold', transform=plt.gca().transAxes, ha='left', va='top')\n", "plt.xlim(35,157)\n", "plt.ylim(45,152)\n", "plt.subplot(122)\n", "plt.imshow(recon_BSREM[:,:,20].T.cpu(), cmap='magma', interpolation='gaussian', vmax=0.5)\n", "plt.axis('off')\n", "plt.text(0.03, 0.97, 'BSREM (RDP)', color='white', fontsize=15, fontweight='bold', transform=plt.gca().transAxes, ha='left', va='top')\n", "plt.xlim(35,157)\n", "plt.ylim(45,152)\n", "fig.tight_layout()\n", "plt.show()" ] } ], "metadata": { "kernelspec": { "display_name": "pytomo_install_test", "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.11.10" }, "papermill": { "default_parameters": {}, "duration": 682.085087, "end_time": "2024-11-19T22:42:33.818138", "environment_variables": {}, "exception": null, "input_path": "t_GE_HDF5.ipynb", "output_path": "t_GE_HDF5.ipynb", "parameters": {}, "start_time": "2024-11-19T22:31:11.733051", "version": "2.6.0" } }, "nbformat": 4, "nbformat_minor": 5 }