{ "cells": [ { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "# Turbulent Pipe Flow (Reynolds 180)\n", "\n", "The simulation of a periodic turbulent pipe flow is used for the validation of the immersed boundary method effectiveness to delineate a curved boundary in a turbulent flow. \n", "As for the turbulent channel case, the external force density term represents the pressure gradient $F_{x}=-\\mathrm{d}p/\\mathrm{d}x$, the regularized no-slip BC is employed at $y=0$ and $y=H$, and periodicity is considered in the remaining boundaries. \n", "\n", "> **_Note:_**\n", " A precursor turbulent flow is used to start the simulation and assure a turbulent flow. \n", " This field can be generated through a precursor poiseuille flow simulation with an Lagrangian body placed to generate perturbations and removed after turbulence is generated." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from nassu.cfg.model import ConfigScheme\n", "\n", "filename = \"tests/validation/cases/05_turbulent_pipe_flow.nassu.yaml\"\n", "\n", "sim_cfgs = ConfigScheme.sim_cfgs_from_file_dct(filename)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "An extra spacing of 2 lattices at each side of the cylinder is kept to assure a complete interpolation-spread procedure." ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "sim_cfg = next(iter(sim_cfgs.values()))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Functions used for processing of turbulent pipe" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "from nassu.cfg.schemes.simul import SimulationConfigs\n", "import numpy as np\n", "\n", "from vtk.util.numpy_support import vtk_to_numpy\n", "from tests.validation.notebooks import common\n", "\n", "stats_export = sim_cfg.output.stats[\"default\"]\n", "last_step = stats_export.interval.get_all_process_steps(sim_cfg.n_steps)[-1]\n", "reader = stats_export.read_vtm_export(last_step)\n", "reader_output = reader.GetOutput()\n", "\n", "\n", "def get_macr_compressed(\n", " sim_cfg: SimulationConfigs, macr_name: str, is_2nd_order: bool\n", ") -> tuple[np.ndarray, np.ndarray]:\n", " global reader, reader_output\n", "\n", " macr_name_read = macr_name if not is_2nd_order else f\"{macr_name}_2nd\"\n", " ds = sim_cfg.domain.domain_size\n", "\n", " p0 = np.array((ds.x // 2, 0, ds.z // 2))\n", " p1 = np.array((ds.x // 2, ds.y - 1, ds.z // 2))\n", " n_points = ds.y\n", " pos = np.linspace(p0, p1, num=n_points, endpoint=True)\n", " norm_pos = (pos[:, 1] + 0.5) / (len(pos) + 1)\n", "\n", " # Sum 0.5 because data is cell centered in vtm\n", " line = common.create_line(p0 + 0.5, p1 + 0.5, n_points - 1)\n", "\n", " probe_filter = common.probe_over_line(line, reader_output)\n", " probed_data = vtk_to_numpy(probe_filter.GetOutput().GetPointData().GetArray(macr_name_read))\n", "\n", " return norm_pos, probed_data" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "## Results" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Load values for comparison" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "import pandas as pd\n", "import numpy as np\n", "import os\n", "\n", "comparison_folder = \"tests/validation/comparison/Turbulent_pipe/Re_tau=180\"\n", "files = [\"ux_avg\", \"uang_rms\", \"ur_rms\", \"ux_rms\"]\n", "get_filename_csv = lambda f: os.path.join(comparison_folder, f + \".csv\")\n", "\n", "\n", "df_cp = {f: pd.read_csv(get_filename_csv(f), delimiter=\",\") for f in files}" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Friction velocity and y+ to use" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "yp = 2.5\n", "u_fric = 0.003401361 * 0.925" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Load simulation velocity fields" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "((152,), (152,))" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "pos, ux_avg = get_macr_compressed(sim_cfg, \"ux\", is_2nd_order=False)\n", "pos, ux_2nd = get_macr_compressed(sim_cfg, \"ux\", is_2nd_order=True)\n", "ux_rms = (ux_2nd - ux_avg**2) ** 0.5\n", "\n", "pos, uy_avg = get_macr_compressed(sim_cfg, \"uy\", is_2nd_order=False)\n", "pos, uy_2nd = get_macr_compressed(sim_cfg, \"uy\", is_2nd_order=True)\n", "uy_rms = (uy_2nd - uy_avg**2) ** 0.5\n", "\n", "pos, uz_avg = get_macr_compressed(sim_cfg, \"uz\", is_2nd_order=False)\n", "pos, uz_2nd = get_macr_compressed(sim_cfg, \"uz\", is_2nd_order=True)\n", "uz_rms = (uz_2nd - uz_avg**2) ** 0.5\n", "\n", "ux_avg.shape, ux_rms.shape" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [], "source": [ "wall_pos = 5 # From 4 to 148\n", "mid_pos = ux_avg.shape[0] // 2\n", "\n", "ux_vals = ux_avg[wall_pos:mid_pos] / u_fric\n", "x_vals = np.arange(len(ux_vals)) * yp" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAiIAAAHLCAYAAAD1F/P9AAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAX5pJREFUeJzt3XlYlPX6P/D3ALIIgqkJosiQuDtCmYgUCskJTdERMVNza7FySRO/KZxM2yAT1EpPWplmJaE4UpFLSqCkqKVi4vEoeFhcAFc2F4SZ5/eHP57jCKMMzMwzwPt1Xc91nXm2uacDzs3nuT/3RyYIggAiIiIiCVhIHQARERE1X0xEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJWUgdgrjQaDS5evIhWrVpBJpNJHQ4REVGjIQgCysrK4OrqCguLB495MBHR4eLFi3Bzc5M6DCIiokbr3Llz6NSp0wPPYSKiQ6tWrQDc/Y/o6OgocTRERESNR2lpKdzc3MTv0gdhIqJD9eMYR0dHJiJERET1UJfSBharEhERkWSYiBAREZFkmIgQERGRZFgj0kBqtRqVlZVSh0GNQIsWLWBpaSl1GEREZoWJSD0JgoDCwkIUFxdLHQo1Iq1bt4aLiwt70xAR/X9MROqpOglp3749WrZsyS8WeiBBEHDz5k1cunQJANChQweJIyIiMg9MROpBrVaLSUjbtm2lDocaCTs7OwDApUuX0L59ez6mISICi1XrpbompGXLlhJHQo1N9c8M64qIiO5iItIAfBxD+uLPDBGRNiYiZDamTp0KpVIpdRgNIpfLsXLlSqnDICJqNJiISEitViM1NRVxcXFITU2FWq02+ntOnToVMpmsxjZ06FCjv/fDfPrpp9iwYYPUYQC4O3KRmJgodRhERE0ei1UlolKpEB4ejtzcXHGfXC5HbGwsQkNDjfreQ4cOxfr167X22djYGPU9H0StVkMmk8HJyUmyGIiISBocEZGASqVCWFgYFAoF0tPTUVZWhvT0dCgUCoSFhUGlUhn1/W1sbODi4qK1PfLII0hNTYW1tTXS0tLEcz/55BO0b98eRUVFAICAgADMmjULs2bNgpOTE9q1a4dFixZBEATxmoqKCsyfPx8dO3aEvb09BgwYgNTUVPH4hg0b0Lp1a/z888/o1asXbGxskJ+fX+PRTEBAAGbPno25c+fikUcegbOzM7766ivcuHED06ZNQ6tWreDp6YkdO3Zofb7MzEwMGzYMDg4OcHZ2xqRJk3DlyhWt+7755pt4++230aZNG7i4uGDJkiXicblcDgAYPXo0ZDKZ+Prs2bMYNWoUnJ2d4eDggP79+2PPnj0N/H+DiKh5YyJiYmq1GuHh4RgxYgQSExPh6+sLBwcH+Pr6IjExESNGjMD8+fNN8pjmfgEBAZg7dy4mTZqEkpISHDt2DIsWLcLXX38NZ2dn8bxvv/0WVlZWOHz4MD799FMsX74cX3/9tXh81qxZSE9Px48//oi///4bY8eOxdChQ5GVlSWec/PmTSxduhRff/01Tp48ifbt29ca07fffot27drh8OHDmD17Nt544w2MHTsWfn5+OHr0KJ599llMmjQJN2/eBAAUFxfjmWeeweOPP46//voLO3fuRFFREZ5//vka97W3t8ehQ4fwySef4P3338fu3bsBAH/++ScAYP369SgoKBBfl5eX47nnnkNycjKOHTuGoUOHIiQkBPn5+Qb4r09E1HhVVVUhMzNT6w/ZOhOoViUlJQIAoaSkpMaxW7duCf/+97+FW7du6X3flJQUAYCQnp5e6/EDBw4IAISUlBS9710XU6ZMESwtLQV7e3ut7aOPPhIEQRAqKioEb29v4fnnnxd69eolvPrqq1rXDx48WOjZs6eg0WjEfQsWLBB69uwpCIIg5OXlCZaWlsKFCxe0rhsyZIgQEREhCIIgrF+/XgAgZGRk1Iht1KhRWu/19NNPi6+rqqoEe3t7YdKkSeK+goICrf+eH3zwgfDss89q3ffcuXMCAOH06dO13lcQBKF///7CggULxNcAhG3btun4r/g/vXv3Fj7//HPxtbu7u7BixQqd5zfkZ4eIyJz85z//EaKjo4WAgADBzs5OACB069ZNEIQHf4fejzUiJlZQUAAA6NOnT63Hq/dXn2cMgYGB+OKLL7T2tWnTBgBgbW2NH374AX379oW7uztWrFhR43pfX1+taagDBw5EbGws1Go1Tpw4AbVajW7dumldU1FRodX8zdraGn379n1orPeeY2lpibZt20KhUIj7qkdqqjuWHj9+HCkpKXBwcKhxr7Nnz4px3f/eHTp0EO+hS3l5OZYsWYJff/0VBQUFqKqqwq1btzgiQkTNyoYNG7Bq1SocOXJEa3+rVq3QoUMHrUf1dcFExMSqW3tnZmbC19e3xvHMzEyt84zB3t4enp6eOo8fOHAAAHDt2jVcu3YN9vb2db53eXk5LC0tceTIkRqdQ+9NDuzs7OrUU6NFixZar2Uymda+6ntoNBrx/UNCQrB06dIa97r3v2lt962+hy7z58/H7t27ERMTA09PT9jZ2SEsLAx37tx56OcgImoqUlJScOTIEVhZWSEoKAghISEIDAxE9+7dYWGhf8UHExET8/f3h1wuR1RUFBITE7X+T9NoNIiOjoaHhwf8/f0lie/s2bN466238NVXXyE+Ph5TpkzBnj17tOI8dOiQ1jUHDx5E165dYWlpiccffxxqtRqXLl2S5DM88cQT2Lp1K+RyOays6v/j3aJFixp1Ovv378fUqVMxevRoAHeTnntnPRERNTVlZWWIiYnBpEmTxD9g33rrLXh5eWHSpEl49NFHG/wekherRkdHo3///mjVqhXat28PpVKJ06dPa51z+/ZtzJw5E23btoWDgwPGjBkjzuLQRRAEvPvuu+jQoQPs7OwQFBSkVSwpFUtLS8TGxiIpKQlKpVJr1oxSqURSUhJiYmKMug5JRUUFCgsLtbYrV65ArVbjxRdfRHBwMKZNm4b169fj77//RmxsrNb1+fn5mDdvHk6fPo24uDh8/vnnmDNnDgCgW7dumDhxIiZPngyVSoWcnBwcPnwY0dHR+PXXX432marNnDkT165dw/jx4/Hnn3/i7Nmz2LVrF6ZNm6ZXAbBcLkdycjIKCwtx/fp1AEDXrl2hUqmQkZGB48ePY8KECQ8dRSEiaozUajW++uoreHp64v3338cHH3wgHvP29sa8efMMkoQAZpCI7N27FzNnzsTBgwexe/duVFZW4tlnn8WNGzfEc9566y388ssv2LJlC/bu3YuLFy8+tNfGJ598gs8++wxr1qzBoUOHYG9vj+DgYNy+fdvYH+mhQkNDkZCQgBMnTsDPzw+Ojo7w8/NDZmYmEhISjN5HZOfOnejQoYPW9vTTT+Ojjz5CXl4e1q5dC+Duo4wvv/wS77zzDo4fPy5eP3nyZNy6dQs+Pj6YOXMm5syZg+nTp4vH169fj8mTJyM8PBzdu3eHUqnEn3/+ic6dOxv1cwGAq6sr9u/fD7VajWeffRYKhQJz585F69at9RoyjI2Nxe7du+Hm5obHH38cALB8+XI88sgj8PPzQ0hICIKDg/HEE08Y66MQEUni0KFD8PHxwfTp03Hp0iV07drVuF2vjVtTq79Lly4JAIS9e/cKgiAIxcXFQosWLYQtW7aI55w6deqBM080Go3g4uIiLFu2TNxXXFws2NjYCHFxcXWKw1izZu5VVVUlpKSkCJs2bRJSUlKEqqqqBt3PFAYPHizMmTNH6jAaLc6aISJzVVJSIsycOVOQyWQCAMHJyUlYsWKFcOfOnXrdS9d36P3MrkakpKQEwP9mcRw5cgSVlZUICgoSz+nRowc6d+6M9PT0Wgs+c3JyUFhYqHWNk5MTBgwYgPT0dLzwwgs1rqmoqEBFRYX4urS01GCfSRdLS0sEBAQY/X2IiIgeJjY2FqtXrwYATJo0CcuWLdPqIWUskj+auZdGo8HcuXPx1FNPidNYCwsLYW1tjdatW2ud6+zsjMLCwlrvU73//v+AD7omOjoaTk5O4ubm5tbAT0NERNR4vP3223j22WeRnJyMjRs3miQJAcxs1szMmTORmZmJP/74w+TvHRERgXnz5omvS0tLmYzU4t5W7URE1HipVCrExcUhPj4eFhYWsLe3x65du0weh9mMiMyaNQtJSUlISUlBp06dxP0uLi64c+cOiouLtc4vKiqCi4tLrfeq3n//zJoHXWNjYwNHR0etjYiIqKm5fv06Jk6ciDFjxiAhIQHffvutpPFInogIgoBZs2Zh27Zt+P333+Hh4aF1vF+/fmjRogWSk5PFfadPn0Z+fj4GDhxY6z09PDzg4uKidU1paSkOHTqk8xoiIqKmbvfu3VAoFNi0aRMsLCwQGRmJCRMmSBqT5I9mZs6ciU2bNuGnn35Cq1atxBoOJycn2NnZwcnJCS+//DLmzZuHNm3awNHREbNnz8bAgQO1ClV79OiB6OhoccXUuXPn4sMPP0TXrl3h4eGBRYsWwdXV1bhTkIiIiMzQ7du3ERERgZUrVwK42/Np48aNGDBggLSBwQwSkeo1T+6fPbJ+/XpMnToVALBixQpYWFhgzJgxqKioQHBwMP71r39pnX/69Glxxg1wt+jmxo0bmD59OoqLi/H0009j586dsLW1NernISIiMjeTJ0/Gli1bAAAzZszAsmXL0LJlS4mjuksmCHquTtNMlJaWwsnJCSUlJTXqRW7fvo2cnBx4eHgwsSG98GeHiKSQkZGBkSNHYs2aNXjuueeM/n4P+g69n+QjIkRERGRYV69exR9//IFRo0YBuNuWPTs7G9bW1hJHVpPkxapERERkOKmpqfDy8kJYWBj+/PNPcb85JiEAE5FmZ+rUqZDJZPj444+19icmJkImk0kUFRERNVRlZSUWLVqEZ555BhcuXMBjjz1mtsnHvZiINEO2trZYunSpuKosERE1bmfPnsWgQYPw4YcfQhAEvPTSSzh69Ci8vLykDu2hmIg0Q0FBQXBxcUF0dHStx69evYrx48ejY8eOaNmyJRQKBeLi4rTOSUhIgEKhgJ2dHdq2bYugoCBxxeTU1FT4+PjA3t4erVu3xlNPPYW8vDwAd0dk7p9CPXfuXK65Q0RUT99//z28vb1x8OBBODk5IT4+HuvWrYO9vb3UodUJi1UNrPrLuDaWlpZaMyUedK6FhQXs7OweeG59f8gsLS0RFRWFCRMm4M0339TqZAvcndnRr18/LFiwAI6Ojvj1118xadIkdOnSBT4+PigoKMD48ePxySefYPTo0SgrK0NaWhoEQUBVVRWUSiVeffVVxMXF4c6dOzh8+DAf+xARGcmVK1dQXl6OQYMGYePGjXB3d5ckDrVajbS0NBQUFOjVnZyJiIE5ODjoPPbcc8/h119/FV+3b98eN2/erPXcwYMHa63rIpfLceXKFa1zGjLzevTo0fD29sbixYuxbt06rWMdO3bE/PnzxdezZ8/Grl27sHnzZjERqaqqQmhoqPgDr1AoAADXrl1DSUkJRowYgS5dugAAevbsWe84iYiopuLiYnEx2DfffBPt2rXD+PHjYWlpKUk8KpUK4eHhyM3N1ftaPpppxpYuXYpvv/0Wp06d0tqvVqvxwQcfQKFQoE2bNnBwcMCuXbuQn58PAPDy8sKQIUOgUCgwduxYfPXVV2K9SZs2bTB16lQEBwcjJCQEn376KQoKCkz+2YiImqLr169j8uTJePLJJ8WRcgsLC7z44ouSJiFhYWFQKBRIT09HWVkZ9uzZU+frmYgYWHl5uc5t69atWudeunRJ57k7duzQOjc3N7fGOQ01aNAgBAcHIyIiQmv/smXL8Omnn2LBggVISUlBRkYGgoODcefOHQB3H+3s3r0bO3bsQK9evfD555+je/fuyMnJAXC3K256ejr8/PwQHx+Pbt264eDBgwDu/sLcP5JTWVnZ4M9CRNTU/fzzz+jduze+++475OTkaK2nJhW1Wo3w8HCMGDECiYmJ8PX1hYODA/r371/ne/DRjIHpU7dhrHP18fHHH8Pb2xvdu3cX9+3fvx+jRo3Ciy++CADQaDQ4c+YMevXqJZ4jk8nw1FNP4amnnsK7774Ld3d3bNu2DfPmzQMAPP7443j88ccRERGBgQMHYtOmTfD19cWjjz6KzMxMrRgyMjLQokULo3w+IqLG7tKlS3jzzTcRHx8PAOjevTvWr19vkkVc76376NChA/z9/bVGXtLS0pCbm4u4uDhYWNRvbIMjIs2cQqHAxIkT8dlnn4n7unbtit27d+PAgQM4deoUXnvtNRQVFYnHDx06hKioKPz111/Iz8+HSqXC5cuX0bNnT+Tk5CAiIgLp6enIy8vDb7/9hqysLLFO5JlnnsFff/2FjRs3IisrC4sXL66RmBAR0d06wG+//RY9e/ZEfHw8LC0t8fbbb+PYsWMmSUJUKhU8PT0RGBiICRMmIDAwEJ6enlCpVOI51Y/e+/TpU+/3YSJCeP/996HRaMTX77zzDp544gkEBwcjICAALi4uWlNuHR0dsW/fPjz33HPo1q0b3nnnHcTGxmLYsGFo2bIl/vOf/2DMmDHo1q0bpk+fjpkzZ+K1114DAAQHB2PRokV4++230b9/f5SVlWHy5Mmm/shERGZPJpMhMTER165dg5eXFw4dOoSlS5dqzag0ltrqPtLT06FQKBAWFiYmIx06dACABv1ByUXvdOCid2QM/NkhogepqKjArVu3xBkx58+fx48//og5c+aY7BG2Wq2Gp6cnFAoFEhMTtR65aDQaKJVKZGZmIisrCwBqPVefRe84IkJERGQGUlNT0bdvX8yePVvc16lTJ8yfP9+kdXTVdR+RkZE16j4sLCwQERGBnJwcpKWlwdLSErGxsUhKSoJSqRRHTw4fPlzn92OxKhERkYSKi4sxf/58sadTaWkprl69irZt20oSz8PqPqr3V58XGhqKhIQEhIeHw8/PT+/344gIERGRRH7++Wf06tVLTEJee+01nDp1yqhJiFqtRmpqKuLi4pCamgq1Wq11/GF1H9X7q88D7iYj2dnZSElJwaZNm5CUlFTneDgiQkREZGLFxcWYPXs2vv/+ewBAt27d8PXXX8Pf39+o71tbB1S5XI7Y2FiEhoYCAPz9/SGXyxEVFVVrjUh0dDQ8PDxqxGppaSmuG1ZaWlrnmDgiQkREZGIajQbJycmwsLDA22+/jYyMDJMkIXWZCaOr7iM9PR1KpRJJSUmIiYkxWCdXzprRoS6zZuRyuUmmUVHTcevWLeTm5nLWDFEzVFFRAWtra3ER0N9//x0tW7aEr6+v0d9bn5kw1QlGbaMnHh4eiImJEUdPdOGsGSOrrl7WtWAdkS7VPzPsJEvUvBw7dgze3t744YcfxH3PPPOMSZIQQL+ZMNXur/tISUlBVlbWQ5MQfbFGpB4sLS3RunVrXLp0CQDQsmVLLnNPDyQIAm7evIlLly6hdevWki1ORUSmpdFosHz5ckRGRqKyshJRUVGSrJKr70yYavfWfRgLE5F6cnFxAQAxGSGqi9atW4s/O0TUtBUWFmLy5MnYvXs3AGD06NH46quvJPlD5N6ZMLWNwtQ2E8ZUWCOiQ12fb6nVaq4eS3XSokULjoQQNRM7d+7E5MmTcfnyZdjZ2eHTTz/FK6+8ItnoeX1qRBpCnxoRjog0kKWlJb9ciIhIdObMGQwfPhwajQZ9+/bFjz/+KC78aWgPWx23WvVMmLCwMCiVSkRERKBPnz7IzMxEdHQ0kpKSkJCQIMn3GRMRIiIiA+rWrRvCw8Nx48YNxMbGGm2GXF16gtxLVwdUDw8PJCQkGLwIta74aEYHfYaViIioedu1axd69uyJzp07A7hboG7MxzDVPUFGjBiByMhIcXQjKipKHN3QlVjUdRSlIfT5DmUiogMTESIiehi1Wo0lS5bgo48+woABA7B3715YW1sb/T1NWe9RH+wjQkREZGSXL19GcHAwPvzwQwiCAG9vb5jib/v69AQxZ6wRISIi0tPBgwcxduxYnD9/Hi1btsSXX36JiRMnmuS969sTxFxxRISIiKiOBEHA6tWrMWjQIJw/fx7du3fH4cOHTZaEAPVbHdecsUZEB9aIEBHR/SoqKuDj44O///4bYWFh+Oabb9CqVSuD3LuuRaSsESEiImqmbGxsoFKpsHLlSmzevNlgSYhKpYKnpycCAwMxYcIEBAYGwtPTU1wR916mXh3X2CRPRPbt24eQkBC4urpCJpMhMTFR67hMJqt1W7Zsmc57LlmypMb5PXr0MPInISKipmjHjh1YuXKl+LpLly6YM2eOwabnVk/FVSgUWkmFQqFAWFhYrclIdU+QEydOwM/PD46OjvDz80NmZqakPUHqQ/JHMzt27MD+/fvRr18/hIaGYtu2bVAqleLxwsLCGue//PLLyM7OxmOPPVbrPZcsWYKEhATs2bNH3GdlZYV27drVOS4+miEiat40Gg2ioqLw7rvvAgBSU1MxaNAgg75HQx+zmKInSH00qhbvw4YNw7Bhw3Qev3+BsJ9++gmBgYE6k5BqVlZWXFyMiIjqpaysDJMnTxZH6V977TUMGDDA4O9TPRU3Li5O51RcPz8/pKWl1boKrilWxzU2yR/N6KOoqAi//vorXn755Yeem5WVBVdXVzz22GOYOHEi8vPzH3h+RUUFSktLtTYiImp+zpw5gwEDBiAxMRHW1tb46quvsGbNGtjY2Bj8vZraVNz6aFSJyLfffotWrVo99NnXgAEDsGHDBuzcuRNffPEFcnJy4O/vj7KyMp3XREdHw8nJSdzc3NwMHT4REZm5X3/9Ff3798epU6fg6uqKffv24ZVXXjHa+zW1qbj10agSkW+++QYTJ0586AJCw4YNw9ixY9G3b18EBwdj+/btKC4uxubNm3VeExERgZKSEnE7d+6cocMnIiIzd+HCBZSWluLpp5/GkSNHjPI45l7+/v6Qy+WIioqCRqPROqbRaBAdHQ0PDw/4+/sbNQ4pSV4jUldpaWk4ffo04uPj9b62devW6NatG7Kzs3WeY2NjY5RhNyIiajymT5+OVq1aYcyYMQ1aM6auRaTVU3HDwsKgVCoREREhLmAXHR0tLmBnDgWoxtJoRkTWrVuHfv36wcvLS+9ry8vLcfbs2SY9tEVERPrLzs7GyJEjce3aNXHf+PHjG5SE6NMTBGhaU3HrQ/JEpLy8HBkZGcjIyAAA5OTkICMjQ6u4tLS0FFu2bNH5nG7IkCFYtWqV+Hr+/PnYu3cvcnNzceDAAYwePRqWlpYYP368UT8LERE1Hjt37kT//v3xyy+/YM6cOQa5Z316ggB3k5Hs7GykpKRg06ZNSElJQVZWVpNPQgAAgsRSUlIEADW2KVOmiOesXbtWsLOzE4qLi2u9h7u7u7B48WLx9bhx44QOHToI1tbWQseOHYVx48YJ2dnZesVVUlIiABBKSkrq87GIiMhMaTQa4eOPPxZkMpkAQBg4cKBw8eLFBt+3qqpKkMvlQkhIiKBWq7WOqdVqISQkRPDw8BCqqqoa/F7mTp/vUMkbmpkrNjQjImp6bt26hVdeeQWbNm0CALzyyitYtWqVQWoEU1NTERgYiPT0dPj6+tY4np6eDj8/P6SkpDT63h8P06gamhEREZlCQUEBRo0ahT///BOWlpb47LPP8MYbbxisVTt7gtSP5DUiREREptCiRQtcvnwZbdq0we7duzFjxgyDJSEAe4LUFx/N6MBHM0RETc/Jkydha2uLLl26GPzeDV03pinR5zuUIyJERNQkCYKApUuXYv369eK+3r17652EqNVqpKamIi4uDqmpqVCr1bWeV90TJCkpCUqlUmvWjFKpRFJSEmJiYpp8EqIv1ogQEVGTU1lZiTfeeAPr1q2DlZUVnnrqKXTr1k3v+6hUKoSHhyM3N1fcJ5fLERsbW+vU2uqeIOHh4fDz8xP3e3h4NIueIPXBEREiImpSSktLMXz4cKxbtw4WFhZYvnx5vZMQ9gQxPtaI6MAaESKixufixYt47rnncPz4cdjb2yM+Ph7Dhw/X+z6s92gY1ogQEVGzc/LkSfj6+uL48eNwdnbG3r1765WEAHfXN8vNzUVkZKRWEgIAFhYWiIiIQE5ODtLS0gwRerPGGhEiImoSfvrpJ5w7dw7du3fHjh074OHhUe97sSeI6TARISKiJiEiIgJWVlZ4+eWX0bZt2wbd696eILV1SWVPEMNhjYgOrBEhIjJ/mzdvRkhICOzs7Ax6X9aINAxrRIiIqEkTBAEREREYN24cXnjhBZ29Pe5V134gAHuCmBIfzRARUaNSWVmJV155BRs3bgQA+Pj41CgovZ++/UAA9gQxFY6IEBFRo1FWVoYRI0Zg48aNsLS0xLp16/DPf/7zgWvG1LcfCMCeIKbAGhEdWCNCRGReioqKMHz4cBw5cgQtW7bEli1b8Nxzzz3wGtZ6SIM1IkRE1KQIgoBRo0bhyJEjaNeuHVJSUh6ahADsB9IYMBEhIiKzJ5PJ8Omnn4qPV3x8fOp0HfuBmD8mIkREZLauXLki/u8BAwYgIyMDnp6edb7+3n4gtWE/EOkxESEiIrO0atUqdOnSBUeOHBH3PWx2zP38/f0hl8sRFRUFjUajdUyj0SA6OhoeHh7w9/c3SMykPyYiRERkVjQaDf7v//4Ps2fPRmlpKRISEup9L/YDMX/sI0JERGbj9u3bmDJlCjZv3gwA+OijjxAREaF1jlqtRlpaGgoKCtChQwf4+/s/MJFgPxDzxum7OnD6LhGRaV25cgVKpRL79+9HixYtsG7dOkyaNEnrnPo0JqumbwJD9cfpu0RE1KgUFBTAz88P+/fvh5OTE3bu3FlrElLfxmTA3cc0AQEBGD9+PAICApiEmAmOiOjAEREiItOpqqrCqFGjcPLkSWzfvh29evXSOs7GZI2LPt+hrBEhIiLJCIIAmUwGKysrxMfHo7y8HC4uLjXOq25MFhcXp7MxmZ+fH9LS0hAQEGCi6MkQ+GiGiIgksWLFCrzxxhuoHph3cHCoNQkB2JisKWMiQkREJqVWqzFnzhzMmzcPa9euxa5dux56DRuTNV1MRIiIyGRu3ryJsWPH4rPPPgMAxMTEIDg4+KHXsTFZ08VEhIiITOLy5csYMmQItm3bBmtra8THxyM8PBwymeyh17IxWdPFYlUiIjK6//73vwgODkZ2djYeeeQR/PTTT/D399ertwcbkzVNTESIiMjoTp8+jZycHLi7u2Pnzp3o0aNHvZqThYaGYtSoUWxM1oTw0QwRERndsGHDsGXLFqSnp4tJSH2bk7ExWdMieSKyb98+hISEwNXVFTKZDImJiVrHp06dCplMprUNHTr0ofddvXo15HI5bG1tMWDAABw+fNhIn4CIiGqTmJiIs2fPiq9Hjx6NDh06QK1WIzw8HCNGjEBiYiJ8fX3h4OAAX19fJCYmYsSIEZg/fz7UarWE0ZOpSJ6I3LhxA15eXli9erXOc4YOHYqCggJxi4uLe+A94+PjMW/ePCxevBhHjx6Fl5cXgoODcenSJUOHT0REtfjiiy8QGhqK4OBgXL16VetYdXOyyMhInc3JcnJykJaWZsqQSSKS14gMGzYMw4YNe+A5NjY2Opvc1Gb58uV49dVXMW3aNADAmjVr8Ouvv+Kbb77BwoULGxQvERHpJggCPvroIyxatAgAEBQUhNatW2udw+ZkdC/JR0TqIjU1Fe3bt0f37t3xxhtv1Miu73Xnzh0cOXIEQUFB4j4LCwsEBQUhPT1d53UVFRUoLS3V2oiIqO4EQcDChQvFJGTx4sX44osvatRwsDkZ3cvsE5GhQ4di48aNSE5OxtKlS7F3714MGzZM57PDK1euQK1Ww9nZWWu/s7MzCgsLdb5PdHQ0nJycxM3Nzc2gn4OIqCnTaDSYNWsWPvnkEwB3R6aXLFlSa48QNieje5l9IvLCCy9g5MiRUCgUYtOaP//8E6mpqQZ9n4iICJSUlIjbuXPnDHp/IqKm7MMPP8S//vUvyGQyrF27Fm+99ZbOc9mcjO5l9onI/R577DG0a9cO2dnZtR5v164dLC0tUVRUpLW/qKjogXUmNjY2cHR01NqIiKhuXn/9dfTp0wffffcdpk+f/tDzq5uTnThxAn5+fnB0dISfnx8yMzPZnKyZaXSJyPnz53H16lWdzw6tra3Rr18/JCcni/s0Gg2Sk5MxcOBAU4VJRNTkVa+aCwBt27bFihUrYGFhgdTU1DpNvQ0NDUV2djZSUlKwadMmpKSkICsri0lIMyP5rJny8nKt0Y2cnBxkZGSgTZs2aNOmDd577z2MGTMGLi4uOHv2LN5++214enpqLZI0ZMgQjB49GrNmzQIAzJs3D1OmTMGTTz4JHx8frFy5Ejdu3BBn0RARUcPcvn0boaGhCAsLQ+vWrfXukFqtujkZNWOCxFJSUgQANbYpU6YIN2/eFJ599lnh0UcfFVq0aCG4u7sLr776qlBYWKh1D3d3d2Hx4sVa+z7//HOhc+fOgrW1teDj4yMcPHhQr7hKSkoEAEJJSUlDPyIRUZNy69YtITg4WAAg2NraCgCEkJAQIT09XSgrKxPS09OFkJAQQSaTCVu3bpU6XJKAPt+hMkG4Z2yNRKWlpXByckJJSQnrRYiI/r/bt29DqVRi165daNmyJRwdHdG/f38kJiZqNSfTaDRQKpXIzMxEVlYWC0+bGX2+QxtdjQgREUnj/iQkKioKhYWF7JBKDcJEhIiIHur+JGTHjh1o3749AHZIpYZhIkJERA8VFxcnJiHbt2/HoEGD2CGVDII1IjqwRoSI6H8EQcB7772HwMBADB48GACgVqvh6ekJhULBGhHSwhoRIiJqsMrKSlRUVAAAZDIZlixZIiYhADukkmEwESEiohrUajUmT56MkSNH4ubNmzrPY4dUaig+mtGBj2aIqLnSaDR49dVX8c0338DKygrJyckYNGjQA69Rq9VIS0tDQUEBOnToAH9/f46ENGP6fIdK3lmViIjMR1VVFZ5//nls27YNFhYW+P777x+ahADskEr1x0czREQEANi6dSvatm2Lbdu2Abg7MrJw4UKoVCqJI6OmjIkIERFBpVIhLCwMpaWlAIBPP/0U6enpUCgUCAsLYzJCRsMaER1YI0JEzYVarYZcLkdhYSGqqqqwfPlyvPXWWwA4DZfqhzUiRERUZ2lpaTh//jzWr1+Pa9euiUkI8L9W7X5+fkhLS2MdCBkcExEiombs9u3bYgv2sLAwODg41DiHrdrJmFgjQkTUTO3btw9dunTBtWvXALBVO0mDiQgRUTN0/PhxhISE4OLFi0hLS4NcLkdUVBQ0Go3WeRqNBtHR0fDw8IC/v79E0VJTxkczRERN3P3NxlxdXREcHIzS0lIMGjQI69evx44dOxAWFgalUomIiAj06dMHmZmZiI6ORlJSEhISElioSkbBRISIqAlTqVQIDw9Hbm6uuM/KygpVVVXw8vLCzz//DDs7O7FVe3h4OPz8/MRzPTw82KqdjIqPZoiImqjq3iAKhQLp6ek4d+4cunbtiqqqKgDA7Nmz4eTkJJ4fGhqK7OxspKSkYNOmTUhJSUFWVhaTEDIq9hHRgX1EiKgxU6vV8PT0hEKhQGJiIiwsLDB//nzExsbC2dkZvXv3Rk5ODnuDkFGwjwgRUTOXlpaG3NxcxMXFwcLi7uD3hx9+iEuXLmHevHm4desWe4OQWWAiQkTUBFX3/KjuAQIAtra22LhxIwCgrKxM6zwiqbBGhIioCaru+bF48WJERkbWmJbL3iBkLlgjogNrRIioMVOr1XB1dcWlS5cAANu2bYNSqQTA9WPI+FgjQkTUzB07dgwlJSUAgM6dO6N9+/YoKytjbxAyO0xEiIiamLy8PIwYMQIVFRXw9vbG9evX8dRTT4nH2RuEzAkTESKiJqSkpATDhw9HUVERFAoF9u7dC3t7e63Oqv7+/hwJIbPBRISIqInQaDQYO3YsTp48CVdXV/z666/i83lO0SVzxUSEiKiRuX/tmOoRDgsLC0yZMgVHjhxBUlIS3NzcpA6V6KE4a0YHzpohInNU29oxcrkcsbGxYs1HaWkp/90iSenzHco+IkREjcT9a8eUlZUhPT0dLi4uGDNmDFQqFQAwCaFGhSMiOnBEhIjMSW1rxwDA8ePHMXDgQAiCgHbt2iE3N5eFqCQ5jogQETUx1WvHREZGiknIlStXMGrUKNy6dQteXl44f/480tLSJI6USD+SJyL79u1DSEgIXF1dIZPJkJiYKB6rrKzEggULoFAoYG9vD1dXV0yePBkXL1584D2XLFkCmUymtfXo0cPIn4SIyHjuXzumsrISzz//PPLy8tClSxds3rxZ6zyixkLyROTGjRvw8vLC6tWraxy7efMmjh49ikWLFuHo0aNQqVQ4ffo0Ro4c+dD79u7dGwUFBeL2xx9/GCN8IiKTqF4TpnqNmPnz5yMlJQUODg746aefcOHCBa3ziBoLyafvDhs2DMOGDav1mJOTE3bv3q21b9WqVfDx8UF+fj46d+6s875WVlZwcXExaKxERFLx9/eHXC5HVFQURo0ahc8++wwA8N1336Fnz55QKpXw8PCAv7+/xJES6UfyRERfJSUlkMlkaN269QPPy8rKgqurK2xtbTFw4EBER0c/MHGpqKhARUWF+Lq0tNRQIRMRNZilpSViY2MxZswY7N27FwCwcOFCODs7Q6lUcu0YarQkfzSjj9u3b2PBggUYP378A6twBwwYgA0bNmDnzp344osvkJOTA39/f5SVlem8Jjo6Gk5OTuLGRkBEZG5CQ0OxdetW8Q+xjz/+GH5+fsjMzOTaMdRomdX0XZlMprVU9b0qKysxZswYnD9/HqmpqXpNqS0uLoa7uzuWL1+Ol19+udZzahsRcXNz4/RdIjI7ujqrEpkLfabvNopHM/dWh//+++96JwatW7dGt27dkJ2drfMcGxsb2NjYNDRUIiK9PSyxEAQBM2bMQI8ePfDmm2/C0tKSa8dQk2H2j2aqk5CsrCzs2bMHbdu21fse5eXlOHv2LKvJicjsqFQqeHp6IjAwEBMmTEBgYCA8PT3FLqkA8MUXX2DNmjWYN28eTpw4IWG0RIYneSJSXl6OjIwMZGRkAABycnKQkZGB/Px8VFZWIiwsDH/99Rd++OEHqNVqFBYWorCwEHfu3BHvMWTIEKxatUp8PX/+fOzduxe5ubk4cOAARo8eDUtLS4wfP97UH4+ISCddLdsVCgXCwsKgUqmwb98+zJkzB8DdWra+fftKHDWRgQkSS0lJEQDU2KZMmSLk5OTUegyAkJKSIt7D3d1dWLx4sfh63LhxQocOHQRra2uhY8eOwrhx44Ts7Gy94iopKREACCUlJQb6pERE/1NVVSXI5XIhJCREUKvVWsfUarUQEhIiuLm5Ce3btxcACOPHjxc0Go1E0RLpR5/vULMqVjUnXGuGiIwpNTUVgYGBSE9Ph6+vb43jaWlpGDRoEACgb9++SE9PR8uWLU0dJlG9cK0ZIiIzd3/L9vv98ssvAABbW1ts2bKFSQg1WUxEiIgkcH/L9vvdunULwN2at27dupksLiJTYyJCRCSBe1u2azQarWMajQZ5eXno1KkTlixZIk2ARCbCRISISALVLduTkpKgVCqRnp6O69ev4/fffxdbtn/66adsVEZNXqNoaEZE1BSFhoYiISEB4eHh8PPzE/d37NiRLdup2eCICBGRhEJDQ5GdnY2lS5eK+2JiYpiEULPBRISISGIFBQX45JNPAAAzZszACy+8IHFERKbDRzNEREZSl8XpKisr8cILL+Dq1at44oknEBsbK1G0RNLgiAgRkRHUZQ0ZAFi0aBH2798PR0dHbN68Gba2thJFTCQNJiJERAZWlzVkAGDHjh1ibcg333yDLl26SBk2kSTY4l0HtngnovpQq9Xw9PSEQqFAYmIiLCz+9/eeRqOBUqlEZmYmsrKyUFRUhAkTJqBv37747LPPJIyayLD0+Q5ljQgRkQGlpaUhNzcXcXFxWkkIAFhYWCAiIgJ+fn5IS0tDQEAA9uzZU6OhGVFzwkSEiMiAHraGTPX+6vOsrPjPMDVvrBEhIjKgh60hs2HDBgDAjz/+yJEQIjARISIyqAetIXP58mW8/fbbAO4mLPc/uiFqjvhbQERkQLWtIVNWVoYDBw7Ay8sLt2/fhqurK5YvXy51qERmgQ8niYgMTNcaMsDdROXnn39Gy5YtJYqOyLzoPSJy8OBBXLhwAcDdYqv09HSDB0VE1NhVryGTkpKC5cuXi43KPvroI/Tr10/i6IjMh96JyI0bNxAeHg4AmDdvHm7dumXwoIiImgJLS0sMGjQI8fHxuH37NgYPHoz58+dLHRaRWdE7ERkyZAjatm2Ld955B23atMEzzzxjjLiIiMyWWq1Gamoq4uLikJqaCrVarfNcCwsLvP322+jSpQu+++67GmvNEDV3enVWDQwMhEwmQ2lpKY4ePYp+/frBwcEBMpkMv//+uzHjNDl2ViWi2qhUKoSHhyM3N1fcJ5fLERsbi9DQUJ3XVVVVsWcINRv6fIfqNSKSkpKC33//HQMGDMDChQvh4+Mj7iMiaurquoYMcPcf4uqmZQAblxHpovejmeTkZFy5cgVRUVG4du0akxAiahbUajXCw8MxYsQIJCYmwtfXFw4ODvD19UViYiJGjBiB+fPni49pZs2aBYVCgV27dkkcOZF50ztFt7OzQ2xsLAAgNjZWa3iSiKip0mcNmaKiInz33XewsLCAg4ODRBETNQ56JyL3zol3dXWFq6urQQMiIjJHdV1D5sSJE3j33XcBAP/85z/x1FNPmSZAokaKnVWJiOrgYWvIVO9ft24diouL4ePjg0WLFpksPqLGql6JyLVr1wwdBxGRWXvQGjIajQbR0dF45JFHcPz4cdjb2+OHH35AixYtJIqWqPGoVyLSrl07uLm5YcSIEfjnP/+JzZs34/Tp09BjJjARUaOiaw2Z9PR0KJVK/PLLLygtLQUAfPbZZ/D09JQ4YqLGQa8+ItVOnjyJjIwMHD9+HBkZGTh27BiuXbsGW1tb9OnTB4cOHTJGrCbFPiJEVJva+oh4eHjgww8/RGpqKq5du4YtW7ZAJpNJFySRxPT5Dq1XInI/QRCwc+dOzJ49G88//zyioqIaekvJMREhIl3UajXS0tJQUFCADh06wN/fX+yYysZlRPp9hxrkt0Umk2HYsGH4/vvv8eWXXxrilkREZsvS0hIBAQEAgJycHK3RDyYhRPox6KwZX19fpKSk6HXNvn37EBISAldXV8hkMiQmJmodFwQB7777Ljp06AA7OzsEBQUhKyvrofddvXo15HI5bG1tMWDAABw+fFivuIiIHubKlSsYOHAggoODcfnyZanDIWqU6pWIODg4YODAgXj99dfxr3/9C/v378fVq1exc+dOlJWV6XWvGzduwMvLC6tXr671+CeffILPPvsMa9aswaFDh2Bvb4/g4GDcvn1b5z3j4+Mxb948LF68GEePHoWXlxeCg4Nx6dIlvWIjouZDn4XsgLt/JM2YMQNFRUW4ePEiWrVqZaJIiZoYoR527NghREdHC+PGjRO6desmWFpaChYWFoKlpaUQFRVVn1sK/79WRdi2bZv4WqPRCC4uLsKyZcvEfcXFxYKNjY0QFxen8z4+Pj7CzJkzxddqtVpwdXUVoqOj6xxLSUmJAEAoKSnR70MQUaOzdetWQS6XCwDETS6XC1u3btV5zaZNmwQAgpWVlfDXX3+ZMFoi86fPd2i9RkSGDh2KhQsX4scff8Tp06dRWlqKv//+G+fPn0dERITBkqScnBwUFhYiKChI3Ofk5IQBAwYgPT291mvu3LmDI0eOaF1jYWGBoKAgndcQUfOlz0J21S5cuIAZM2YAABYtWoR+/fqZOmyiJkOvROTdd9/FkSNHauxv2bIlevfuDRcXF4MFBgCFhYUAAGdnZ639zs7O4rH7XblyBWq1Wq9rAKCiogKlpaVaGxE1bfouZAfcfSTz8ssvo7i4GE8++aRB//giao70SkTOnz+PYcOGoVOnTnjjjTewY8cO3Llzx1ixmVR0dDScnJzEzc3NTeqQiMjIqheyi4yM1LmQXU5ODtLS0sT9a9euxa5du2Bra4uNGzeyeypRA+mViHzzzTcoLCxEXFwcWrVqhblz56Jdu3YYM2YMNm7caPDW79UjLEVFRVr7i4qKdI6+tGvXDpaWlnpdAwAREREoKSkRt3PnzjUweiIyd3VdyK76PADo168fevTogY8//hg9e/Y0fpBETZzeNSIWFhbw9/fHJ598gtOnT+PQoUMYMGAA1q5dC1dXVwwaNAgxMTG4cOFCg4Pz8PCAi4sLkpOTxX2lpaU4dOgQBg4cWOs11tbW6Nevn9Y1Go0GycnJOq8BABsbGzg6OmptRNS01XUhu+rzAKB///44evQoZs+ebfwAiZqBBi9617NnT7z99tvYv38/8vPzMWXKFKSlpSEuLq5O9yovL0dGRgYyMjIA3C1QzcjIQH5+PmQyGebOnYsPP/wQP//8M06cOIHJkyfD1dUVSqVSvMeQIUOwatUq8fW8efPw1Vdf4dtvv8WpU6fwxhtv4MaNG5g2bVp9Pi4RNVF1WcjOw8MD/v7+uH79unjMzs6uxqMcIqqn+kzLkclkQqdOnYThw4cLkZGRQnx8vHD69GlBo9Hofa+UlBStKXPV25QpUwRBuDuFd9GiRYKzs7NgY2MjDBkyRDh9+rTWPdzd3YXFixdr7fv888+Fzp07C9bW1oKPj49w8OBBveLi9F2i5mHr1q2CTCYTQkJChAMHDgilpaXCgQMHhJCQEEEmkwlbt24Vjh8/LrRq1UpYunSpoFarpQ6ZyOzp8x1q8EXvFAoFDh48aNBkSQpca4ao+dC1kF1MTAyGDx8OHx8f/P333xg5ciQSExO5oB3RQ3DROwNgIkLUvOhayC4yMhLR0dFo164dMjMza7QGIKKauOgdEZGe7l3Irlp6ejqWLl0K4O60XSYhRIYn+aJ3RETm6MaNG5g8eTI0Gg0mTZqE0NBQqUMiapLqNSLi4OAAhUIBLy8v9O3bF15eXujRowf+/PNPvRe9IyIyFl2PW+piwYIFyM7ORqdOnfDZZ58ZOVKi5qteiUhCQoI45fbTTz/F2bNnIQgCZDIZPvjgA0PHSESkt9oKUOVyOWJjY+s0uuHu7g4bGxt88803aN26tfECJWrmDFKsevPmTeTk5KBt27YGX29GKixWJWq8qheyGzFiBCIjI9GnTx9kZmYiKioKSUlJSEhIqFMyUlhY2GT+TSMyJZPPmmmKmIgQNU5qtRqenp5QKBRITEzUajym0WigVCqRmZmJrKysWh/TVFZWcv0YogbS5zuUrQGJqEmpz0J21VQqFfr27Yu//vrLVOESNXtMRIioSanPQnbA3YUxX3vtNfznP//Btm3bjBskEYmYiBBRk1KfhewEQcD06dNx5coVeHl5YfHixcYPlIgAMBEhoiZGn4Xsqn3zzTf4+eefYW1tjY0bN8La2trUYRM1W0xEiKhJsbS0RGxsLJKSkqBUKpGeno6ysjKkp6dDqVQiKSkJMTExYqHq2bNnMWfOHADABx98gL59+0oZPlGzY5AW70RE5iQ0NBQJCQkIDw+Hn5+fuN/Dw0Nr6m5VVRUmTZqEGzduYNCgQQgPD5cqZKJmi4kIETVJoaGhGDVq1AM7q5aVlcHBwQGOjo7YuHFjnbuuEpHhsI+IDuwjQtQ8aDQanDlzBj169JA6FKImg31EiIge4N4iVgsLCyYhRBJiIkJEzc6cOXMwffp0lJeXSx0KUbPHGhEialZ27tyJVatWAQDGjRuHIUOGSBwRUfPGEREiajauXr2Kl156CQAwe/ZsJiFEZoCJCBE1C4Ig4LXXXkNBQQF69OiBjz/+WOqQiAhMRIiomfjuu++wdetWWFlZ4YcffkDLli2lDomIwBoRImoE1Gr1A/uBPExubi5mzZoFAHjvvffwxBNPGCtUItITExEiMmsqlQrh4eHIzc0V98nlcsTGxoodUh8mLy8PNjY26Nu3LxYsWGCkSImoPvhohojMlkqlQlhYGBQKhdaaMQqFAmFhYVCpVHW6z+DBg3HixAls2rSJ3VOJzAw7q+rAzqpE0lKr1fD09IRCoUBiYiIsLP73d5NGo4FSqURmZiaysrJ0JheCIEAmk5kqZCL6/9hZlYgavbS0NOTm5iIyMlIrCQHudkONiIhATk4O0tLSar3+9u3bCAwMxNatW00RLhHVExMRIjJLBQUFAIA+ffrUerx6f/V594uIiMDevXsxc+ZMlJWVGSdIImowJiJEZJY6dOgAAMjMzKz1ePX+6vPulZycjJUrVwIA1q1bh1atWhknSCJqMCYiRGSW/P39IZfLERUVpbVIHXC3RiQ6OhoeHh7w9/fXOnb9+nVMmTIFAPD6669j+PDhJouZiPTHRISIzJKlpSViY2ORlJQEpVKpNWtGqVQiKSkJMTExNQpVZ86ciQsXLqBr166IiYmRKHoiqiv2ESEisxUaGoqEhASEh4fDz89P3O/h4YGEhIQafUTi4uIQFxcHS0tLfP/997C3tzd1yESkJyYiRGTWQkNDMWrUqDp1Vj1+/DgAYNGiRfDx8TF1qERUD42ij4hcLkdeXl6N/TNmzMDq1atr7N+wYQOmTZumtc/Gxga3b9+u83uyjwhR45Samoqnn34aVlb8O4tIKvp8hzaK39Q///wTarVafJ2ZmYl//OMfGDt2rM5rHB0dcfr0afE1mxoRNQ8BAQFSh0BEemgUxaqPPvooXFxcxC0pKQldunTB4MGDdV4jk8m0rnF2djZhxERkKidOnMBzzz2Hc+fOSR0KEdVDo0hE7nXnzh18//33eOmllx44ylFeXg53d3e4ublh1KhROHnypAmjJKK6UKvVSE1NRVxcHFJTU7VGPuuioqICL774Inbs2MHF7IgaqUaXiCQmJqK4uBhTp07VeU737t3xzTff4KeffsL3338PjUYDPz8/nD9/Xuc1FRUVKC0t1dqIyHhUKhU8PT0RGBiICRMmIDAwEJ6ennVeyA4A5s6di7///hvt2rXDihUrjBgtERlLo0tE1q1bh2HDhsHV1VXnOQMHDsTkyZPh7e2NwYMHQ6VS4dFHH8XatWt1XhMdHQ0nJydxc3NzM0b4RATDrKq7ceNGrFmzBjKZDN999x0fvxI1Uo1i1ky1vLw8PPbYY1CpVBg1apRe144dOxZWVlaIi4ur9XhFRQUqKirE16WlpXBzc+OsGSIDM8Squn///Td8fX1x69YtLF68GEuWLDFR9ERUF0129d3169ejffv2erdsVqvVOHHiRK1rUlSzsbGBo6Oj1kZEhtfQVXVLSkowZswY3Lp1C8HBwVi0aJEpwiYiI2k0iYhGo8H69esxZcqUGv0BJk+ejIiICPH1+++/j99++w3//e9/cfToUbz44ovIy8vDK6+8Yuqwieg+DV1V9/r167C1tUXnzp3x/fff6xw1IaLGoVH0EQGAPXv2ID8/Hy+99FKNY/n5+Vp/WV2/fh2vvvoqCgsL8cgjj6Bfv344cOAAevXqZcqQiagW966q6+vrW+P4g1bVBe42ODx48CDOnz+Pdu3aGS9QIjKJRlUjYkrsrEpkHPWtEblx4wbXjiFqJJpsjQgRNX71WVW3oKAA3bt3x8cffwyNRiNh9ERkaExEiMjkqlfVPXHiBPz8/ODo6Ag/Pz9kZmbWWFW3srIS48aNw4ULF/DDDz9ozW4josav0dSIEFHTUtdVdSMjI5GWloZWrVph69atsLOzkyhiIjIGJiJEJBlLS8sHLlKnUqkQExMD4O6q2t26dTNRZERkKkxEiMig1Gr1Q0c56uLMmTPiUg7z58/XelxDRE0Ha0SIyGAMsX4MANy+fRtjxoxBWVkZ/P39ER0dbaSIiUhqTESIyCAMsX5MNRsbG8yYMQOdO3dGfHx8jSaGRNR0sI+IDuwjQlR3hlg/pja3b9+Gra2tMUImIiNiHxEiMqmGrh9T7d///jeKi4vF10xCiJo+JiJE1GANXT8GAK5evYphw4ahX79+OHPmjOGDJCKzxESEiBrs3vVjavOw9WPUajUmTpworhvl7OxsnECJyOwwESGiBvP394dcLkdUVFSNFuwajQbR0dHw8PCAv79/rdd/+OGH2LVrF+zs7LB161Y4OTmZImwiMgNMRIioweqzfky1nTt34r333gMArF27Fn379jV1+EQkIc6JIyKDqF4/Jjw8HH5+fuJ+Dw+PGuvHVMvLy8PEiRMhCAJee+01TJo0yZQhE5EZ4PRdHTh9l6h+9OmsOmbMGKhUKjz55JNIS0vjLBmiJkKf71COiBCRQT1s/Zh7rV27FtbW1oiOjmYSQtRMMREhIsm0a9cOcXFxUodBRBJisSoRmdSJEyewYcMGqcMgIjPBEREiMpnS0lKMGTMGWVlZKC8vx6xZs6QOiYgkxkSEiGqlT9FpXQiCgGnTpiErKwtubm544YUXDBgtETVWfDRDRDWoVCp4enoiMDAQEyZMQGBgIDw9PfVaQfd+y5cvh0qlQosWLZCQkIB27doZMGIiaqyYiBCRFpVKhbCwMCgUCq3GZAqFAmFhYfVKRvbt24cFCxYAAFauXAkfHx9Dh01EjRT7iOjAPiLUHKnVanh6ekKhUCAxMVFrJV2NRgOlUonMzExkZWXV+TFNYWEhHn/8cRQWFmLixIn47rvvIJPJjPURiMgM6PMdyhERIhKlpaUhNzcXkZGRWkkIAFhYWCAiIgI5OTlIS0ur8z3T09NRVFSE3r17Y+3atUxCiEgLi1WJSFRQUAAA6NOnT63Hq/dXn1cXo0ePRlxcHB5//HHY29s3PEgialI4IkJEog4dOgAAMjMzaz1evb/6PF3Ky8tx6dIl8fW4cePQrVs3A0VJRE0JExEiEvn7+0MulyMqKgoajUbrmEajQXR0NDw8PODv76/zHjdv3kRISAgGDx6MixcvGjtkImrkmIgQkcjS0hKxsbFISkqCUqnUmjWjVCqRlJSEmJgYnYWqt2/fxujRo5GamooLFy7gwoULJv4ERNTYsEaEiLSEhoYiISEB4eHh8PPzE/d7eHggISEBoaGhtV53584dhIWF4bfffoO9vT22b9+O/v37mypsImqkOH1XB07fpeZOn86qlZWVGDduHLZt2wZbW1ts374dgYGBJo6YiMyFPt+hHBEholpZWloiICDgoeep1WpMnjwZ27Ztg7W1NRITE5mEEFGdsUaEiBrkypUrOHz4MFq0aIGtW7ciODhY6pCIqBEx+0RkyZIlkMlkWluPHj0eeM2WLVvQo0cP2NraQqFQYPv27SaKlqj5cXZ2xr59+5CYmIgRI0ZIHQ4RNTJmn4gAQO/evVFQUCBuf/zxh85zDxw4gPHjx+Pll1/GsWPHoFQqxbbURM2VWq1Gamoq4uLikJqaCrVa3aD7CYKAjIwM8XXHjh3x3HPPNTBKImqOGkUiYmVlBRcXF3F70Kqdn376KYYOHYr/+7//Q8+ePfHBBx/giSeewKpVq0wYMZH5MPRKuoIg4P/+7//w5JNP4scffzRwtETU3DSKRCQrKwuurq547LHHMHHiROTn5+s8Nz09HUFBQVr7goODkZ6e/sD3qKioQGlpqdZG1NgZeiVdQRDwzjvvIDY2Fmq1GmVlZUaKnIiaC7NPRAYMGIANGzZg586d+OKLL5CTkwN/f3+d/wAWFhbC2dlZa5+zszMKCwsf+D7R0dFwcnISNzc3N4N9BiIpqNVqhIeHY8SIEUhMTISvry8cHBzg6+sr1nPMnz9fr8c0H3zwAaKiogAAq1atwquvvmqs8ImomTD7RGTYsGEYO3Ys+vbti+DgYGzfvh3FxcXYvHmzQd8nIiICJSUl4nbu3DmD3p/I1Ay9ku7SpUuxePFiAEBsbCxmzpxp8JiJqPlpdH1EWrdujW7duiE7O7vW4y4uLigqKtLaV1RUBBcXlwfe18bGBjY2NgaLk0hqhlxJd+XKlVi4cCEA4KOPPsK8efMMFCURNXdmPyJyv/Lycpw9e1bn6p8DBw5EcnKy1r7du3dj4MCBpgiPyGwYaiVdQRCQlZUFAHj33XcRGRlpwCiJqLkz+xbv8+fPR0hICNzd3XHx4kUsXrwYGRkZ+Pe//41HH30UkydPRseOHREdHQ3g7vTdwYMH4+OPP8bw4cPx448/IioqCkePHtX5l2Ft2OKdGju1Wg1PT08oFAokJiZqPZ7RaDTitPasrCydrdurCYKAX3/9FcOHD4dMJjN26ETUyOnzHWr2IyLnz5/H+PHj0b17dzz//PNo27YtDh48iEcffRQAkJ+frzW07Ofnh02bNuHLL7+El5cXEhISkJiYqFcSQtQUNHQl3b1796KyshIAIJPJMGLECCYhRGRwZj8iIhWOiFBToVKpEB4ejtzcXHGfh4cHYmJidK6ku2XLFrzwwgsICQnB5s2bYW1tbaJoiagp4KJ3RCQKDQ3FqFGj6ryS7k8//YQJEyZAo9Ggbdu2sLLiPxNEZDz8F4aoGajrSro7duzA2LFjUVVVhRdffBFffvlljam/RESGxH9hiAgAsGfPHowePRqVlZUYO3Ys1q9f/9AiViKihmIiQkTYu3cvRo4ciYqKCowaNQo//PADH8kQkUkwESEiCIIACwsLDBs2DPHx8WjRooXUIRFRM8E/eYgIAQEB+OOPP9C9e3d2GCYik2IiQtRMZWRkwMrKSuyx4+3tLW1ARNQsMREhMkNqtbrO023rIzMzE0FBQZDJZEhNTUXv3r0Ndm8iIn2wRoTIzKhUKnh6eiIwMBATJkxAYGAgPD09oVKpDHL/06dPIygoCFevXoVcLkenTp0Mcl8iovpgIkJkRlQqFcLCwqBQKLRasisUCoSFhTU4GTl79iyeeeYZFBUVwdvbG7t27YKTk5OBoici0h9bvOvAFu9kaoZcpK42eXl5GDRoEPLz89G7d2+kpKSIazYRERlSk1r0jqi5SEtLQ25uLiIjI2t0M7WwsEBERARycnKQlpam970vXLiAZ555Bvn5+ejWrRv27NnDJISIzAKLVYnMRPUq0rpWiq7ef+9q03XVqlUruLq6AgB+//13uLi41DNKIiLD4ogIkZno0KEDgLszWmpTvb/6PH04Ojpix44d2Lt3Lzp27Fj/IImIDIyJCJGZ8Pf3h1wuR1RUFDQajdYxjUaD6OhoeHh4wN/f/6H3EgQBGzZswBtvvCHuc3Bw4AwZIjI7TESIzISlpSViY2ORlJQEpVKpNWtGqVQiKSkJMTExDy1UzcrKwpAhQzBt2jSsWbMGycnJJvoERET6YyJCZEZCQ0ORkJCAEydOwM/PD46OjvDz80NmZiYSEhIQGhqq89o7d+4gKioKCoUCKSkpsLOzw7JlyzB48GATfgIiIv1w+q4OnL5LUtK3s2p6ejqmT58u1pEEBwfjiy++gIeHh6lCJiIS6fMdykREByYi1BDGbtF+r8rKSnTt2hV5eXl49NFHsXLlSowfPx4ymcwo70dE9DD6fIdy+i6RgalUKoSHhyM3N1fcJ5fLERsb+8BHK/oSBAEymQwtWrTAqlWrsHXrVsTExKBt27YGew8iImNjjQiRARm7RTsAnD9/HqNHj8a6devEfSNGjMD69euZhBBRo8NHMzrw0Qzpy9gt2tVqNdasWYOIiAiUlZWhffv2yMvLg62trSE/BhFRg7HFO5EEjNmi/cSJE3j66acxa9YslJWVwdfXF8nJyUxCiKjRYyJCZCDGaNF+69Yt/POf/8QTTzyBgwcPolWrVli9ejX++OMPne9DRNSYMBEhMhBjtGg/efIkoqOjUVVVBaVSiVOnTmHGjBlGm4FDRGRqrBHRgTUipC9D1YhUVlaiRYsW4usPPvgAffr0wejRo40aPxGRobBGhEgCDW3RLggCfvjhBzz22GM4ffq0uH/RokVMQoioyWIiQmRA9W3R/t///hdDhw7Fiy++iPPnz2PZsmUmjpyISBp8NKMDH81QQ9S1s2plZSVWrFiBJUuW4NatW7CxscG7776L+fPnw9raWoLIiYgajp1ViQygIW3aLS0tERAQ8MBz/vzzT7z66qs4fvw4ACAwMBBr165F165dGxo6EVGjwUczRLVQqVTw9PREYGAgJkyYgMDAQHh6ehqkM2q15ORkHD9+HG3atMH69euRnJzMJISImh0mIkT3MWab9pKSEvF/h4eHY+HChTh16hSmTp3KReqIqFky+xqR6OhoqFQq/Oc//4GdnR38/PywdOlSdO/eXec1GzZswLRp07T22djY4Pbt23V+X9aINE/GatNeUFCAOXPm4OTJkzh27BjrP4ioSWtS03f37t2LmTNn4uDBg9i9ezcqKyvx7LPP4saNGw+8ztHREQUFBeKWl5dnooipMTN0m3aNRoMvv/wSPXv2xJYtW3D69Ol6tXgnImqqzL5YdefOnVqvN2zYgPbt2+PIkSMYNGiQzutkMhlcXFyMHR41MYZs037q1ClMnz4df/zxBwDgySefxFdffQVvb2/DBEtE1ASY/YjI/aqfsbdp0+aB55WXl8Pd3R1ubm4YNWoUTp48+cDzKyoqUFpaqrVR82OINu1VVVVYsmQJvLy88Mcff8De3h4rV67EwYMHmYQQEd3H7GtE7qXRaDBy5EgUFxeLf2XWJj09HVlZWejbty9KSkoQExODffv24eTJk+jUqVOt1yxZsgTvvfdejf2sEWm86jP91hA1IoIgICgoCL///juGDx+Of/3rX+jcubNBPxsRkTnTq85SaERef/11wd3dXTh37pxe1925c0fo0qWL8M477+g85/bt20JJSYm4nTt3TgAglJSUNDRsksDWrVsFuVwuABA3uVwubN26tU7XymQyISQkRDhw4IBQWloqHDhwQAgJCRFkMlmt97h27ZrWz8qZM2eEzZs3CxqNxqCfi4ioMSgpKanzd2ijeTQza9YsJCUlISUlReeohi4tWrTA448/juzsbJ3n2NjYwNHRUWujxqmh02/1adMuCALi4+PRs2dPREREiPu7du2KsWPHckouEdFDmP2jGUEQMHv2bGzbtg2pqan1avikVqvRu3dvPPfcc1i+fHmdruH03cbJkNNvH/ZoJy8vDzNmzMD27dsBAL169cJff/0FOzs743w4IqJGokm1eJ85cyY2bdqEn376Ca1atUJhYSEAwMnJSfwHf/LkyejYsSOio6MBAO+//z58fX3h6emJ4uJiLFu2DHl5eXjllVck+xxkGtXTb+Pi4nROv/Xz80NaWtpDW7DratNeVVWFzz//HO+88w5u3rwJa2trREZGYuHChbCxsTHgpyEiavrMPhH54osvAKDGF8L69esxdepUAEB+fr7Wl87169fx6quvorCwEI888gj69euHAwcOoFevXqYKmyRiyOm3tTlz5gwmTJiAI0eOAAD8/f3x5ZdfokePHvW6HxFRc2f2iUhdnhylpqZqvV6xYgVWrFhhpIjInN07/dbX17fG8bpMv32Q1q1b47///S9at26NZcuW4aWXXqox8kJERHVn9jUiUmGNSONkjBbtR48exRNPPCG+TklJQc+ePdkwj4hIhybV4p2aBrVajdTUVMTFxSE1NRVqtdoo72NpaYnY2FgkJSVBqVRqzZpRKpVISkpCTExMnZKQS5cuYeLEiejXrx8SExPF/YGBgUxCiIgMhIkIGZ1KpYKnpycCAwMxYcIEBAYGwtPTs0Gr2D6IPtNvayMIAtavX4+ePXti06ZNsLCwwL///W+jxEpE1NwxESGjamhPj/oKDQ1FdnY2UlJSsGnTJqSkpCArK+uhSciZM2cwZMgQvPTSS7h27Rq8vb1x6NAhREZGGiVOIqLmjjUiOrBGpOGMUa9hTGvXrsWcOXNQUVEBOzs7vP/++5g7dy6srMy+ppuIyKywRoTMQnVPj8jISJ09PXJycpCWliZRhNrc3d1RUVGB4OBgnDx5EvPnz2cSQkRkZPxXlozG2D099CEIAoqKipCbmytueXl50Gg0WLt2LQBg6NChSEtLw1NPPcXW7EREJsJEpAmpz2qzxmTsnh730mg0KCwsRG5uLq5du4YRI0aIx0aNGoVdu3ahoqKixnX29vZYs2aNmHg8/fTTDY6FiIjqjolIE6FSqRAeHo7c3Fxxn1wuR2xs7EMLNI3F398fcrkcUVFRtdaIREdHw8PDA/7+/g+9lyAIWqMU69atw8GDB8XRjfz8fNy5cwfA3eSirKxM6/yKigpYWFigU6dOkMvlcHd3h1wuh1wuh0ajMYsaFSKi5oiJSBNQPTNlxIgRiIuLQ58+fZCZmYmoqCiEhYXVacqqMVT39AgLC4NSqURERIQYW3R0NJKSkpCQkCAmAQUFBcjOzhYfm9z7GKW4uBiXL18Wk4uffvoJv/zyi9b7WVhYwM3NDe7u7rh16xZatmwJAIiNjcXKlSvRqVMntGjRwrT/EYiI6IE4a0aHxjJrpjHMTFGpVJg3bx7y8vLEfa1bt4aPjw927twpJhcjR46skVzc6+rVq2jTpg0A4IcffkB2drbW6EbHjh2ZaBARmYEmtfouPZghV5ttiMrKSpw/fx55eXkYPHiwmFxERkYiLi4O58+f1zq/uLgYv/32G65fvy4mF926dcNjjz0mPjK5/xFK69atxesnTpxotM9CRESmw0TEQKQqFJViZsru3bvFBKj6Mcr58+eh0WgAANeuXcMjjzwCACgrKxPrVqytrdG5c2etROPeOo6YmBjExMQYLE4iIjJ/TEQMQMpCUUPNTKmoqEB+fn6t9Rl5eXk4fvy4OHLx888/Y9WqVTXuYWNjA3d3d1y9elVMRGbOnInx48dDLpfDxcWFK9USEZEWJiINJHWhaF1npvTv3x9nzpzRSjAWLlwoPrsLDw/H6tWrdb5PXl6emIgEBgaiqqpK67GJXC5H+/btayQaPXr0MMKnJiKipoLFqjrUpdDGXApFVSoVxowZg4CAACxZsgRPPPEEMjMz8cYbb+D48eN45JFHcP369RrXHTt2DN7e3gCAZcuWYcmSJTXqMqo3hUIBOzs7o30GIiJqOvQpVmUiokNd/iOmpqYiMDAQ6enptT4WSU9Ph5+fH1JSUgxWKJqZmVmjPiM3NxdFRUU1zr0/AbG3t9dKLubOnQtPT08Ad4tNrays2FGUiIgajLNmTMSQhaLl5eU66zPWr1+P3r17AwC2b9+OBQsW1HqPVq1aYfHixXB1dUWHDh3Qvn17nDp1Skw82rRpozPR4LRXIiKSAhORBtCnULS0tFQryRg9ejQ6deoEAFi5ciXeeustne+TnZ0tJiKPP/44Ro0aVesjlNatW9dINHr16mWQz0pERGQMTEQa4N5C0Q0bNqBly5awtbUFAOzZsweTJk2CtbU1QkNDa9RouLu7i4mIq6srgLuPUmrrnzFw4EDxun/84x/4xz/+YaJPSEREZFxMRPRQUFCgtb5JdafQX375BW3btsUnn3yC119/HZmZmViwYAEKCwsBQFwDpW3btmJyUT29FQBCQkLMvoMrERGRMbBYVYfqQptff/0Vzz33HAAgPj4eL7zwQp2u79SpE4KDg6FUKsURjlatWhkzZCIiIrPAYlUDyszMFBORrl27wsfHp9bprR07dsSxY8dM3lmViIioMeOIiA7V2dzhw4fRv39/qcMhIiJqNPQZEWG/7Yfo3r271CEQERE1WUxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDKNJhFZvXo15HI5bG1tMWDAABw+fPiB52/ZsgU9evSAra0tFAoFtm/fbqJIiYiIqK4aRSISHx+PefPmYfHixTh69Ci8vLwQHByMS5cu1Xr+gQMHMH78eLz88ss4duwYlEollEolMjMzTRw5ERERPUijaPE+YMAA9O/fH6tWrQIAaDQauLm5Yfbs2Vi4cGGN88eNG4cbN24gKSlJ3Ofr6wtvb2+sWbOmTu+pT3taIiIi+p8m1eL9zp07OHLkCIKCgsR9FhYWCAoKQnp6eq3XpKena50PAMHBwTrPB4CKigqUlpZqbURERGRcZp+IXLlyBWq1Gs7Ozlr7nZ2dUVhYWOs1hYWFep0PANHR0XBychI3Nze3hgdPRERED2T2iYipREREoKSkRNzOnTsndUhERERNnpXUATxMu3btYGlpiaKiIq39RUVFcHFxqfUaFxcXvc4HABsbG9jY2DQ8YCIiIqozsx8Rsba2Rr9+/ZCcnCzu02g0SE5OxsCBA2u9ZuDAgVrnA8Du3bt1nk9ERETSMPsREQCYN28epkyZgieffBI+Pj5YuXIlbty4gWnTpgEAJk+ejI4dOyI6OhoAMGfOHAwePBixsbEYPnw4fvzxR/z111/48ssvpfwYREREdJ9GkYiMGzcOly9fxrvvvovCwkJ4e3tj586dYkFqfn4+LCz+N7jj5+eHTZs24Z133kFkZCS6du2KxMRE9OnTR6qPQERERLVoFH1EpMA+IkRERPXTpPqIEBERUdPFRISIiIgkw0SEiIiIJMNEhIiIiCTDRISIiIgkw0SEiIiIJMNEhIiIiCTDRISIiIgkw0SEiIiIJNMoWrxLobrhbGlpqcSREBERNS7V3511ad7ORESHq1evAgDc3NwkjoSIiKhxKisrg5OT0wPPYSKiQ5s2bQDcXVDvYf8RiUg6paWlcHNzw7lz57guFDVpjelnXRAElJWVwdXV9aHnMhHRoXo1XycnJ7P/P5yIAEdHR/6uUrPQWH7W6/pHPItViYiISDJMRIiIiEgyTER0sLGxweLFi2FjYyN1KET0APxdpeaiqf6sy4S6zK0hIiIiMgKOiBAREZFkmIgQERGRZJiIEBERkWSYiBAREZFkmIjUQ1JSErp3746uXbvi66+/ljocItJh9OjReOSRRxAWFiZ1KERGc+7cOQQEBKBXr17o27cvtmzZInVIeuGsGT1VVVWhV69eSElJgZOTE/r164cDBw6gbdu2UodGRPdJTU1FWVkZvv32WyQkJEgdDpFRFBQUoKioCN7e3igsLES/fv1w5swZ2NvbSx1anXBERE+HDx9G79690bFjRzg4OGDYsGH47bffpA6LiGoREBCAVq1aSR0GkVF16NAB3t7eAAAXFxe0a9cO165dkzYoPTS7RGTfvn0ICQmBq6srZDIZEhMTa5yzevVqyOVy2NraYsCAATh8+LB47OLFi+jYsaP4umPHjrhw4YIpQidqVhr6u0rUWBjyZ/3IkSNQq9WNauX4ZpeI3LhxA15eXli9enWtx+Pj4zFv3jwsXrwYR48ehZeXF4KDg3Hp0iUTR0rUvPF3lZoLQ/2sX7t2DZMnT8aXX35pirANR2jGAAjbtm3T2ufj4yPMnDlTfK1WqwVXV1chOjpaEARB2L9/v6BUKsXjc+bMEX744QeTxEvUXNXnd7VaSkqKMGbMGFOESdRg9f1Zv337tuDv7y9s3LjRVKEaTLMbEXmQO3fu4MiRIwgKChL3WVhYICgoCOnp6QAAHx8fZGZm4sKFCygvL8eOHTsQHBwsVchEzVJdfleJmoK6/KwLgoCpU6fimWeewaRJk6QKtd6YiNzjypUrUKvVcHZ21trv7OyMwsJCAICVlRViY2MRGBgIb29vhIeHc8YMkYnV5XcVAIKCgjB27Fhs374dnTp1YpJCjU5dftb379+P+Ph4JCYmwtvbG97e3jhx4oQU4daLldQBNEYjR47EyJEjpQ6DiB5iz549UodAZHRPP/00NBqN1GHUG0dE7tGuXTtYWlqiqKhIa39RURFcXFwkioqI7sffVWoumsPPOhORe1hbW6Nfv35ITk4W92k0GiQnJ2PgwIESRkZE9+LvKjUXzeFnvdk9mikvL0d2drb4OicnBxkZGWjTpg06d+6MefPmYcqUKXjyySfh4+ODlStX4saNG5g2bZqEURM1P/xdpeai2f+sSz1tx9RSUlIEADW2KVOmiOd8/vnnQufOnQVra2vBx8dHOHjwoHQBEzVT/F2l5qK5/6xzrRkiIiKSDGtEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIioydiwYQM2bNggdRhEpAcmIkRERCQZJiJEREQkGa6+S0SN2p07d+Dj4wMAuHbtGgCgTZs2AIDDhw/D2tpastiI6OGYiBCRWerUqRMiIyMxY8YMcd+BAwcQFBSEU6dOwd3dvcY11fUhU6dONVGURNRQfDRDRGZpwIAB+PPPP8XXgiBg7ty5eOutt2pNQoiocWIiQkRmydfXVysR+e6773Du3DlERERIGBURGRofzRCRWUpLS0NAQABKSkogk8nQvXt3vPfee3j55ZelDo2IDMhK6gCIiGrTr18/WFhY4OjRo9izZw8effRRTJs2TeqwiMjAmIgQkVlq2bIlFAoFtm7diq+++grbt2+HhQWfJhM1NfytJiKz5evri88//xzBwcEICAiQOhwiMgImIkRktry8vNCiRQssW7ZM6lCIyEhYrEpEZiswMBBPPPEEYmNjpQ6FiIyENSJEZFY0Gg0uX76MdevWISsrCz/99JPUIRGRETERISKzsm/fPjzzzDPo0aMHtm7dCkdHR6lDIiIj4qMZIiIikgyLVYmIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgyTESIiIhIMkxEiIiISDJMRIiIiEgy/w931VhUlNGIBAAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "import matplotlib.pyplot as plt\n", "\n", "fig, ax = plt.subplots()\n", "fig.set_size_inches(6, 5)\n", "\n", "ax.plot(df_cp[\"ux_avg\"][\"y+\"], df_cp[\"ux_avg\"][\"u/u*\"], \"ko\", fillstyle=\"none\")\n", "ax.plot(x_vals, ux_vals, \"k--\")\n", "ax.legend([\"Experimental\", \"Nassu\"])\n", "\n", "ax.set_xscale(\"symlog\")\n", "ax.set_xlim((1, 170))\n", "\n", "ax.set_ylabel(\"$u/u*$\")\n", "ax.set_xlabel(\"$y^+$\")\n", "\n", "plt.show(fig)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "The average flow velocity shown below presents an excellent agreement with experimental results. Below, the root mean squared velocity $u_{\\alpha,\\mathrm{rms}}$ is also presented against experimental data." ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAgQAAAHGCAYAAAALypO0AAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAmEBJREFUeJzs3XdYFFf3B/DvAgIKYkERsIFiB8XYUSKoiQ0Uib7WqDHxTdEUSwpGkxhfQY1GY/nFaKImKiYRECMaOxgi2EHBFlTsKFZAEJDd8/vjZheWvn0Xzud59gFmZmfu7C47Z245V0JEBMYYY4xVa2aGLgBjjDHGDI8DAsYYY4xxQMAYY4wxDggYY4wxBg4IGGOMMQYOCBhjjDEGDggYY4wxBsDC0AWoDJlMhrt376J27dqQSCSGLg5jjDFmMogIWVlZcHZ2hplZ2fUAJhEQ3L17F02bNjV0MRhjjDGTdevWLTRp0qTM9SYRENSuXRuAOBk7OzsDl4YxxhgzHZmZmWjatKniWloWkwgI5M0EdnZ2HBAwxhhjaqioyZ07FTLGGGOMAwLGGGOMcUDAGGOMMZhIHwLGGGO6J5PJkJ+fb+hiMBXVqFED5ubmGu+HAwLGGGPIz89HamoqZDKZoYvC1FC3bl04OjpqlKuHAwLGGKvmiAhpaWkwNzdH06ZNy01ew4wLESEnJwfp6ekAACcnJ7X3xQEBY4xVcwUFBcjJyYGzszNq1apl6OIwFdWsWRMAkJ6eDgcHB7WbDzgMZIyxak4qlQIALC0tDVwSpi55IPfixQu198EBAWOMMQAVJ65hxksb7x0HBIwxxhjjPgSVJZVKERsbi7S0NDg5OcHb21srwzwYY4wxY8A1BJUQEREBNzc3+Pr6Yty4cfD19YWbmxsiIiIMXTTGGGNMKzggqEBERARGjhwJDw8PxMfHIysrC/Hx8fDw8MDIkSM5KGCMsX9JpVLExMRg27ZtiImJUXRW1KXJkydDIpGUeAwaNEjnx65qJEREhi5ERTIzM1GnTh1kZGTodbZDqVQKNzc3eHh4IDIyUmlsrkwmQ0BAAJKTk5GSkgIA3KTAGDNJubm5SE1NhaurK6ytrdXaR0REBGbNmoXr168rlrm4uGDZsmUIDAzUUklLmjx5Mu7fv4+NGzcqLbeyskK9evV0dlxjU957WNlrKNcQlCM2NhbXr1/HnDlzSiTqMDMzQ1BQEFJTU7Fw4UKVmhSeP3+OzMxMfZwCY4zpnKFrUq2srODo6Kj0qFevHh48eABHR0cEBwcrto2Li4OlpSUOHToEAPDx8cH06dMxffp01KlTBw0aNMC8efOg6r3y3r17YWNjo5TpMTk5GRKJBA8fPtTOieoYBwTlSEtLAwC4u7uXul6+/KuvvqrUP0Jubi6Cg4PRsGFD2NvbIyAgAJGRkRqNG2WMMUOSSqWYNWsW/Pz8EBkZiZ49e8LW1hY9e/ZEZGQk/Pz8MHv2bL00HxTXsGFDbNiwAV999RVOnTqFrKwsvP7665g+fTr69++v2O7nn3+GhYUFTpw4ge+++w7ffvstfvzxR5WOlZCQAHd3d6Wbx8TERDg7O6NBgwZaOyedIhOQkZFBACgjI0Ovx42OjiYAFB8fX+r62NhYAkBeXl4klUqV1kmlUvL39ydXV1cqKCggIqKbN29SzZo1CYDSo2HDhvTRRx/RpUuXdH5OjDFW3PPnz+nChQv0/PlzlZ9b0fdkXFwcAaDo6GgNS1m6SZMmkbm5OdnY2Cg9Fi5cqNjmvffeo9atW9O4cePIw8ODcnNzFev69u1L7dq1I5lMplj26aefUrt27VQqx3/+8x+aOnWq0rLZs2fT4MGD1Twz1ZT3Hlb2Gso1BOXw9vaGi4sLgoODS0z4IZPJ8MknnwAAlixZUm6TQmxsLACgadOm+Oabb7B161YkJyfj448/hqOjIx48eIAVK1bgxIkT+jkxxhjTksrWpMq30wVfX18kJiYqPd555x3F+qVLl6KgoADbt2/H1q1bYWVlpfT8nj17KiX26dWrF1JSUlSq1UhISEDHjh2VliUmJqJTp05qnpX+cUBQDnNzcyxbtgxRUVEICAhQahKQ/w2g1Df88ePHik4uhw8fViyfNm0axo0bhw4dOmDJkiW4desWoqKiMG7cOKWONzk5OTo+O8YY05x8Mp3k5ORS18uXazLpTkVsbGzg5uam9Khfv75i/dWrV3H37l3IZDKlTo/akp2djatXrypdC2QyGRISEpSW9e7dG8ePHwcAvPnmm1i+fLnWy6IJDggqEBgYiLCwMCQlJcHLywt2dnbw8vJCcnIy5s+fD6DkP8KLFy/w8ssvY/369QDEh7EsFhYWGDp0KLZu3QobGxsAwL1799C5c2csWbJER2fFGGPaUVFNakhICFxdXeHt7W2Q8uXn52PChAkYPXo0FixYgLfeeksxM6Cc/CItd+zYMbRq1arSI8Xk00a3bdtWsWzfvn149OiRUkAwb948LFq0CN9++y3MzMwwY8YMDc5MB3TVnqFNhupDUFRBQQFFR0dTaGgoRUdHU0FBARUUFJCLiwv5+/sr9SFYu3YtASBLS0tydHRU9CGorNWrVyv6F8yfP1+pbYsxxrRNkz4ERETh4eEkkUjI39+f4uLiKDMzk+Li4sjf358kEgmFh4drucSFJk2aRIMGDaK0tDSlx4MHD4hItOO7uLhQRkYGSaVS6tOnDw0dOlTx/L59+5KtrS3NmDGDLl26RKGhoWRjY0Nr165VbLNq1Srq169fmWW4e/cuSSQS2r17NxERxcfHU4sWLcja2rrE9/9LL71E/fr1o/z8fG2+DFrpQ8ABAZV+sa+s4v8I9+7dowYNGigu6Or+IyxcuFCxjzlz5nBQwBjTGU0DAiLxXeji4qLUYdrV1VWnwQCRCAhQrKM2AGrTpg1FR0eThYUFxcbGKrZPTU0lOzs7+r//+z8iEgHBe++9R++88w7Z2dlRvXr1Snznfvnll9S8efNyy/G///2P6tWrR82aNaNJkybRp59+Sl27dlXa5sSJE+Ti4kIjRozQ3gvwLw4ItKC0D7GLi4tKH+LS9mFhYUG//vqrRmVbtmyZYn8zZszgoIAxphPaCAiINLu5MpS+ffvShx9+qPPj3L59mzp27EjXrl2jzp07U1JSklb3z6MMNKStZBqBgYG4cuUK/vjjD8Wc1D/99BNGjx6tUflmzpyJNWvWAACWL1+OadOmlWijY4wxY2Fubg4fHx+MHTsWPj4+nK31X8+fP8eoUaOwatUquLq6IigoCAsWLDB0sUqotqmLVUlLXNkPtUwmw++//47ffvsNYWFhWvtn+OmnnzB16lS0bNkSx44dg729vVb2yxhjgHZSF5sqHx8feHp6YsWKFYYuika0kbq42k5/LE9LvG3btjJzCHh5eSE2NhY+Pj6V2qeZmRnGjBmDMWPGaLWsb775JurUqYPu3btzMMAYY1oUExNj6CIYjWobEGg7mUZBQQEsLHT3co4cOVLp79zc3GoXyTPGGNOdatuHQJvJNP755x+4urri+++/r3BCDG1MD7phwwa0adMGt27dUvm5jDHGWGmqbUCgzWQa8+bNw+3bt7Fnzx6l9JfFRUREqDQrYmny8/OxcuVK3Lx5EyNGjMDz588r/VzGGGOsLNU2IKgoLXFUVBSWLl1aYcfA06dP4/fff4dEIsHChQvL3E5bIxosLS0RGRkJe3t7nD59Gm+//bbK03QyxhhjJWh1IKSO6DsPgSrJNF599VUCQBMmTChzm7IyGhKVPitiZRw6dIjMzc0JAC1fvrzSz2OMseK0lYeAGY428hBU22GHRUmlUsTGxiItLQ1OTk7w9vau1JDBw4cPo3///qhRowYuXbqEFi1alLpdTEwMfH19ER8fj549e5ZYHx8fDy8vL0RHR1d6RAMAfPfdd/joo49gbm6Offv2Kc3vzRhjlVWdhx1WFTzsUEvkyTRUJW8iePvtt8sMBgDdTQ/6wQcfICEhAT///DP+85//ICkpCc7OzirtgzHGGANU7EMQEhKCbt26oXbt2nBwcEBAQAAuX75c7nM2bdoEiUSi9KgKEejTp09x5MgRACKjYHl0NT2oRCLB2rVr0b17d3z44YdwdHRU6fmMMcaYnEo1BEeOHMG0adPQrVs3FBQUYM6cOXj11Vdx4cIFxdS9pbGzs1MKHMrriW8qateujSNHjuDYsWNwdXUtd9uiIxpKy4qoyfSg1tbWiI2NhaWlpcrPZYwxxuRUCgj27t2r9PemTZvg4OCA06dP4+WXXy7zeRKJpMrdvZqbm6N3797o3bt3pbZdtmwZRo4ciYCAAAQFBcHd3R3JyckICQlBVFSURqmOiwYDOTk5iI6OxtChQ9XaF2OMsepJo2GHGRkZAID69euXu92zZ8/QvHlzNG3aFMOHD8f58+fL3T4vLw+ZmZlKD1MXGBiIsLAwJCUlwcvLC3Z2dvDy8kJycjLCwsIQGBio8TGysrLQu3dvDBs2DPv379dCqRljTAVSKRATA2zbJn6qkXhNVZMnTy7RLC2RSDBo0CCdH7uqUTsgkMlk+Oijj9C7d+8yO8sBQJs2bbBhwwbs3LkTW7ZsgUwmg5eXF27fvl3mc0JCQlCnTh3Fo2nTpuoWUycuXbqEd999F3/++adKz5PPihgdHY3Q0FBER0cjJSVFK8EAANja2qJLly6QyWQYM2YMrl69qpX9MsZYhSIiADc3wNcXGDdO/HRzE8t1bNCgQUhLS1N6bNu2TefHrXLUHfP4zjvvUPPmzenWrVsqPS8/P59atmxJc+fOLXOb3NxcysjIUDxu3bqlszwE6liyZAkBoMGDBxu6KCXk5uZSjx49CAB5eHhQVlaWoYvEGDNyGuchCA8nkkiI/P2J4uOJsrLET39/sbySeV3UMWnSJBo+fHip69LT06lRo0a0cOFCxbKjR49SjRo16ODBg0RE1LdvX5o2bRpNmzaN7OzsyN7enubOnUsymUylcvz5559Uq1YtpVwzSUlJBIAePHig+ompSBt5CNSqIZg+fTqioqIQHR2NJk2aqPTcGjVqoHPnzrhy5UqZ21hZWcHOzk7pYUzkfSkGDx5s4JKUZGVlhfDwcDg6OiIpKQlvvPEGZzJkjOmOVArMmgX4+QGRkUDPnoCtrfgZGSmWz56tl+aD4ho2bIgNGzbgq6++wqlTp5CVlYXXX38d06dPV8rb8vPPP8PCwgInTpzAd999h2+//RY//vijSsdKSEiAu7u7UqfxxMREODs7o0GDBlo7J11SKSAgIkyfPh07duzA4cOHK+xdXxqpVIqkpCSVh9gZi2fPniE2NhYAdNZGpekESI0bN0ZYWBhq1KiBsLAwLF68WCflZIwxxMYC168Dc+YAxaaSh5kZEBQEpKaK7XQkKioKtra2So/g4GAAwJAhQzB16lSMHz8e77zzDmxsbBASEqL0/KZNm2L58uVo06YNxo8fj/fffx/Lly9XqQyJiYno1KmT0rKzZ8+WWGbMVAoIpk2bhi1btiA0NBS1a9fGvXv3cO/ePaUJdiZOnIigoCDF319//TX279+Pa9eu4cyZM5gwYQJu3LiBt956S3tnoUfR0dF48eIFWrRoATc3N63vXxsTIAFA7969sWrVKgDAihUrqkTHTMaYEZInVCurL5l8uYqJ11Th6+uLxMREpcc777yjWL906VIUFBRg+/bt2Lp1K6ysrJSe37NnT6Xh8L169UJKSopKN2MJCQno2LGj0rLSggQ5dWa61TWVAoLvv/8eGRkZ8PHxgZOTk+Lx22+/Kba5efOmUsa9J0+eYOrUqWjXrh2GDBmCzMxMxMXFoX379to7Cz2SdyQcNGiQ1vMpaGsCJLm3334bixYtwvHjx42u2YUxVkXIa3vLSLymWK7DWmEbGxu4ubkpPYqOfrt69Sru3r0LmUyG69eva/342dnZuHr1qtLFXyaTISEhQWnZsGHD8N5776Fbt27YuHGj1suhMd10b9AubU5uVFBQQNHR0RQaGkrR0dEqTSgkk8kUEyHt2rVL47IUL5e2J0AqTU5OjkbPZ4xVPRp1KiwoIHJxER0Ii313kVQqlru6iu10oLxOhUREeXl51KlTJ5o0aRIFBweTg4MD3b9/X7G+b9++1L59e6XnfPbZZ9SuXbtKl0HeeTA9PV2xbM+ePQSALly4oFjm4uJCK1asqPR+VaGNToXVKiAobWZDFxeXSs9s+PDhQ3JzcyNLS0ut996Pjo4mABQfH1/q+ri4OAJA0dHRah8jMjKSnJ2dKTk5We19MMaqHq2OMoiLI8rMFD/1NMpg0KBBlJaWpvSQ9+yfPXs2ubi4UEZGBkmlUurTpw8NHTpU8fy+ffuSra0tzZgxgy5dukShoaFkY2NDa9euVWyzatUq6tevX5lluHv3LkkkEtq9ezcREcXHx1OLFi3I2tpacROXmZlJTZs21cVLQEQGHGVgirRRHW9vb4+UlBSkpqbC1tZWq+XT1QRIcjKZDEuWLMHdu3fxyiuv4Nq1a+oVlDHGigsMBMLCgKQkwMsLsLMTP5OTxXIt5Vopy969e5WasZ2cnNCnTx/ExMRgxYoV2Lx5M+zs7GBmZobNmzcjNjYW33//veL5EydOxPPnz9G9e3dMmzYNH374If773/8q1j98+LDcvC5OTk5YsGABJkyYgObNm2Pt2rUYNWoU3N3dFRloz58/Dy8vL929CFpQLaY/lkqlcHNzg4eHR6lzCQQEBCA5ORkpKSlqpw/WlK6mSC7q8ePH6Nu3L5KTk+Hi4oK///4bjRs31rDkjDFTp7Xpj6VSMZogLU30GfD2Bgz0nVpZPj4+8PT0xIoVK3R6nPXr1yM9PR2ff/65TvavjemPq0UNQWxsLK5fv445c+YoBQMAYGZmhqCgIKSmpiqGE5amoKAA+fn5Oitj0QmQZDKZ0jpNJ0CSq1+/Pvbv34+WLVvi+vXreOWVV/DgwQNNi84YY4K5OeDjA4wdK34aeTCgT0lJSfDw8DB0McpVLQICbVTHHzp0CA0aNMB7772n/QKicAKkqKgoBAQEKDVrBAQEICoqCkuXLtW4BsPJyQkHDx5EkyZNcPHiRQwaNAhPnz7Vzkkwxhgr1cqVKzFs2DBDF6Nc1SIgkCdBSi5jWIx8eXnJkv78809kZWXhxYsX2i/gv/QxARIAuLi44MCBA2jYsCHOnDmD9evXa2W/jDFmauT9DBj3Iah0H4K2bdvi8uXLCA8P19qFubzyxsbGIi0tDU5OTvD29tZJ34Zz585h3bp1WLlyZYmmFMZY9aG1PgTMYLTRh8BC14U0BvLq+JEjRyIgIABBQUFwd3dHcnIyQkJCEBUVhbCwsDIvuqmpqbh8+TLMzc2V8l/rsrzqdhxURceOHbF69WrF33l5eThx4oRG/RQYY4yZpmpzW6hJdfy+ffsAAF5eXqhTp46+iqxXMpkMb7zxBnx8fPDDDz8YujiMMcb0rFrUEMgFBgZi+PDhKlfHF01XXFXJZDJYW1tDJpPhnXfeQWpqKoKDg7kpgTHGqolqFRAAqlfH5+fn49ChQwCMb7pjbfY1sLCwwE8//QRXV1d88cUXWLx4MW7cuIGNGzdymyJjjFUDfPtXgby8PAQFBcHf39+oprHU1qyIRUkkEsybN08xN/ivv/6KQYMG4cmTJ1osOWOMMWPEAUEFateujc8//xx//PGH0VSfa3tWxOImTpyIvXv3ws7ODkeOHIG/vz9MYDAKY4wxDVSLYYdViT7TMJ87dw4jRozATz/9pJdRD4wxw+Bhh6aPUxfrGBHhzz//xK1bt4zmDlkbaZgrq2PHjrh48aJSMPDs2TON98sYY8z4cEBQjlu3bmHIkCFo0aKFTjMUqkLXsyIWZ2lpqfj9/PnzaNGiBbZu3aqVfTPGGDMeHBCU49y5cwBElsKiF0ZD0kYaZnWtX78eDx48wIQJE7B48WKt758xZtqkUiAmBti2TfyUSnV/zMmTJ0MikZR4VOVh4rrCAUE5kpKSAIiqc2Ohj1kRy/Ltt99i5syZAIDPPvsMCxcu1PoxGGOmKSICcHMDfH2BcePETzc3sVzXBg0ahLS0NKXHtm3bdH/gKoYDgnLIawiMKSDQ16yIpTEzM8OyZcsUtQNz587F8uXLtX4cxphpiYgARo4EPDyA+HggK0v89PAQy3UdFFhZWcHR0VHpUa9ePTx48ACOjo4IDg5WbBsXFwdLS0tFfhkfHx9Mnz4d06dPR506ddCgQQPMmzdP5X5je/fuhY2NjdKNWnJyMiQSCR4+fKidE9U1MgEZGRkEgDIyMvR63Pbt2xMA2r17t16PWxnh4eHk4uJCABQPV1dXCg8P18vx58+frzju2rVr9XJMxphuPH/+nC5cuEDPnz9X+bkFBUQuLkT+/kRSqfI6qVQsd3UV2+nCpEmTaPjw4WWu3717N9WoUYNOnjxJmZmZ1KJFC5oxY4Zifd++fcnW1pY+/PBDunTpEm3ZsoVq1apF69atU6kcwcHB1L17d6VlmzdvJmdnZ5X2o67y3sPKXkO5hqAMeXl5uHz5MgDjqiGQCwwMxJUrVxAdHY3Q0FBER0cjJSVF5zMxys2bNw+ffPIJAGDr1q2Q6qOxkDFmdGJjgevXgTlzgOKpWszMgKAgIDVVbKcrUVFRsLW1VXrIawWGDBmCqVOnYvz48XjnnXdgY2ODkJAQpec3bdoUy5cvR5s2bTB+/Hi8//77Ktd+JiYmlkhed/bsWaVl+fn5mDp1Ktq3b49evXrh8ePHap6xblS71MWVdfHiRUilUtSrVw+NGzc2dHFKpa9ZEUsjkUiwaNEiNG/eHJMmTdJJEwVjzPjJBzSVMfBJsVxLA59K5evri++//15pWf369RW/L126FO7u7ti+fTtOnz4NKysrpW179uwJiUSi+LtXr15YtmwZpFJppb/bEhIS8MEHHygtS0xMRNeuXRV/z58/H3369MH69esxd+5c/Prrr3jvvfcqfZ66xgFBGZo2bYpffvkFWVlZSh8UVkgikZT4MF+7dg0tWrQwUIkYY/omH9CUnAz07FlyvXxAlA4GPinY2NjAzc2tzPVXr17F3bt3IZPJcP36dXh4eGj1+NnZ2bh69apSbYBMJkNCQgLefPNNAEBGRgb++usvRWdsV1dXXLx4Uavl0BQ3GZTB3t4er7/+ulFFb8aMiDB//ny0b98eBw8eNHRxGGN64u0NuLgAwcFAsYFPkMmAkBDA1VVsZwj5+fmYMGECRo8ejQULFuCtt95Cenq60jbHjx9X+vvYsWNo1apVpWsHUlNTIZPJ0LZtW8Wyffv24dGjR4og4eDBg7hy5Qo8PT3h6emJuXPnKtViGAMOCJhWyGQyJCYmIi8vD8OHD0dCQoKhi8QY0wNzc2DZMiAqCggIUB5lEBAgli9dKrbTlby8PNy7d0/pIe/Z//nnnyMjIwMrV67Ep59+itatW2PKlClKz7958yZmzpyJy5cvY9u2bVi1ahU+/PBDxfrVq1ejf//+ZR7f3t4eEokEJ0+eBCACiunTp8Pa2hqtW7cGIPoTLFmyBImJiUhMTISnp6dRTZgHcEBQph9++AExMTHIz883dFFMgrm5OX799VcMGDAAOTk5GDZsmNayJTLGjFtgIBAWBiQlAV5egJ2d+JmcLJbruq/z3r174eTkpPTo06cPYmJisGLFCmzevBl2dnYwMzPD5s2bERsbq9TnYOLEiXj+/Dm6d++OadOm4cMPP8R///tfxfqHDx/i6tWrZR7fyckJCxYswIQJE9C8eXOsXbsWo0aNgru7u6KW4enTp4q+C5mZmUhISICvr6+OXhH18ORGpXjw4AEcHBwAAFlZWbC1tdX5MauKp0+folevXrh06RK6deuGI0eOoGbNmoYuFmOsHNqa3EgqFaMJ0tJEnwFvb93WDGiDj48PPD09sWLFCp0eZ+XKlbh58yaWLl2K2bNno379+pgzZ47W9s+TG1WCVCpFTEwMtm3bhpiYmEoNj5NnKGzZsmWVCQbUeR3UUbduXezatQv169fHyZMnMXny5BIZFRljVZO5OeDjA4wdK34aezCgT+PGjUNMTAxatWqFvLw8fPrpp4YuUglVOiCIiIiAm5sbfH19MW7cOPj6+sLNzQ0RFaTNkmco1HZPVENR93VQl3zfNWrUwO+//46YmBidHIcxxkxFgwYNcOrUKaSkpGDVqlVGOVS7ygYEERERGDlyJDw8PJTS+3p4eGDkyJHlXgyNMWWxujR5HTTRt29frF+/Hlu2bEG/fv10cgzGGNOUvJ8Bq6J9CKRSKdzc3ODh4YHIyEiYFUmfJZPJEBAQgOTkZKSkpJQapXXr1g2nTp3C9u3bMXLkSK2eiz5p+jowxqoHbfUhYIbDfQjKEBsbi+vXr2POnDlKF0FATNATFBSE1NRUxJaSS1MqlSqmEDb1GgJNXgdtS0tLg7+/P27cuKHzYzHGGFNdlQwI5MPd3MvIpSlfXtqwuCtXriA3Nxc1a9ZEy5YtdVdIPdDkddC2t99+G1FRURg2bBhycnJ0fjzGGGOqqZIBgdO/OTLld/rFyZc7lZJLs3nz5oiLi8Mvv/xi8tXomrwO2rZmzRo0atQI586dw/Tp03V+PMYYY6rhPgQmftEvj7G9DtHR0RgwYABkMhk2bNiAN954Q+fHZIxVjPsQmD7uQ1AGc3NzLFu2DFFRUQgICFDqXR8QEICoqCgsXbq0SgcDgPG9Dr6+vpg/fz4AYNq0aYp8D4wxxgyvSgYEABAYGIiwsDAkJSXBy8sLdnZ28PLyQnJyMsLCwhBYRi7Nr776Cj/++CMyMzP1XGLdUPd10JU5c+Zg0KBBeP78OUaOHFllXmfGGDN1VbLJoCipVIrY2FikpaXByckJ3t7eZd4RP3v2DLVr1wYApKeno2HDhhqX3Vio8jro2sOHD9G5c2fY2dkhKioKrq6uBikHY0zgJgPTp40mAwtdF9LQzM3N4ePjU6lti3ayq0rBAKDa66BrDRo0wP79+9GsWTPY2NgYujiMMcZQDQICVVSlDIXGrl27dkp/5+bm8p0JY4wZUJXtQ6AOeSe3qjKHgSmQyWRYsmQJ3N3d8fjxY0MXhzFmYiZPngyJRIJFixYpLY+MjIREIjFQqUwTBwRFcA2B/mVnZ+OHH37A1atXMXnyZJhAlxbGmJGxtrbG4sWL8eTJE0MXxaRxQPAvIuKA4F/6mioZAGrXro2wsDBYWlpi165dWLNmjc6OxRirmgYMGABHR0eEhISUuc3evXvRp08f1K1bF/b29vDz88PVq1cV68PCwuDh4YGaNWvC3t4eAwYMQHZ2dqXWu7i4lJggydPTE1999VW55bGxsVGaHj45ORkSiQQPHz5U8RXQDg4I/nXnzh08ffoU5ubmaNu2raGLYzD6nioZADp37oxvvvkGADB79mxFYMYYM6zs7OwyH7m5uZXe9vnz55XaVl3m5uYIDg7GqlWrcPv27TLPZebMmTh16hQOHToEMzMzjBgxAjKZDGlpaRg7diymTJmCixcvIiYmBoGBgYoay4rWqyMhIQHu7u5KCeMSExPh7OyMBg0aqL1fTXCnwn81adIE9+/fx5UrV2BlZWXo4hiEfKpkPz8/bNu2De7u7khOTkZwcDBGjhyp07wF77//Pvbv34/du3djzJgxOHXqFGrVqqWTYzHGKsfW1rbMdUOGDMHu3bsVfzs4OJQ5T0nfvn0RExOj+NvFxaXUu2BNLrAjRoyAp6cnvvzyS/z0008l1r/22mtKf2/YsAENGzbEhQsXkJ+fj4KCAgQGBqJ58+YAlPuSpaWllbteHYmJiejUqZPSsrNnz5ZYpk9cQ1CEg4MDvLy8DF0Mg5BKpZg1axb8/PwQGRmJnj17wtbWFj179kRkZCT8/Pwwe/ZsnTUfSCQSbNy4EU5OTrh48SJmzJihk+MwxqquxYsX4+eff8bFixdLrEtJScHYsWPRokUL2NnZwcXFBQBw8+ZNdOrUCf3794eHhwdGjRqF9evXK/VHqGi9OhISEko0T5cWJOgTBwQMgHFMldywYUNs3rwZ1tbWaN++PXcwZMzAnj17VuYjPDxcadv09PQyt/3zzz+Vtr1+/Xqp22nq5ZdfxsCBAxEUFFRinb+/Px4/foz169fj+PHjOH78OAAgPz8f5ubmOHDgAP7880+0b98eq1atQps2bZCamgoAFa43MzMr8X314sWLMsuZnZ2Nq1evKl38ZTIZEhISlJb17t1bUc4333wTy5cvV/OVqRxuMvjXu+++izp16uCDDz6As7OzoYujd8YyVXL//v1x/fp1NGrUSKfHYYxVTJXEYbraVlWLFi2Cp6cn2rRpo1j26NEjXL58GevXr4e3tzcA4O+//1Z6nkQiQe/evdG7d2988cUXaN68OXbs2IGZM2dWuL5hw4ZK342ZmZmKYKE0qampkMlkSv3V9u3bh0ePHikFBPPmzcOiRYvg7e0NMzMzndeccg0BgLy8PPz4449YvHixTnvUGzNjmiq5aDDw7NkzFBQU6PyYjLGqwcPDA+PHj8fKlSsVy+rVqwd7e3usW7cOV65cweHDhxUXegA4fvw4goODcerUKdy8eRMRERF48OCBIoFaRev79euHzZs3IzY2FklJSZg0aVK5qeHt7e0hkUhw8uRJAMCxY8cwffp0WFtbo3Xr1ortBg0ahJs3b2L37t34v//7P62+TqXhgADApUuXUFBQgLp166JJkyaGLo5BeHt7w8XFBcHBwUrDYABRlRUSEgJXV1dFdK0Pp0+fRufOnbFw4UK9HZMxZvq+/vprpe8xMzMz/Prrrzh9+jTc3d0xY8YMxcgmALCzs8Nff/2FIUOGoHXr1pg7dy6WLVuGwYMHV2p9UFAQ+vbtCz8/PwwdOhQBAQFo2bJlmeVzcnLCggULMGHCBDRv3hxr167FqFGj4O7urhRInDx5Eo8fP0adOnVQo0YNbb9MJZEJyMjIIACUkZGhk/1v3ryZAJC3t7dO9m8qwsPDSSKRkL+/P8XFxVFmZibFxcWRv78/SSQSCg8P12t5tmzZQgDIzMyM/vrrL70em7Hq5Pnz53ThwgV6/vy5oYvC/nX79m3q2LEjXbt2jTp37kxJSUnlbl/ee1jZa2iVqyFQJ6mOvDq8uqcsNrapksePH4+JEydCJpNhzJgxuH//vl6PzxhjhvD8+XOMGjUKq1atgqurK4KCgrBgwQKdH7dKdSqMiIjArFmzcP36dcUyFxcXLFu2rNyLmTxbVdG2m+oqMDAQw4cPN5qpktesWYOTJ0/i4sWLGDNmDA4cOAALiyr1sWWMMSU1a9ZEXFyc4u9Ro0Zh1KhROj9ulakhkCfV8fDwQHx8PLKyshAfHw8PDw+MHDmy3Ex78gBCPi61upNPlTx27Fj4+PgYLBgARGKU8PBw2NraIiYmBnPnzjVYWRhjrCqrEgGBpkl17t27BwBwdXXVZ7FZJbVr1w4bNmwAIBKP7Nq1y8AlYoyxqkelgCAkJATdunVD7dq14eDggICAAFy+fLnC523fvh1t27aFtbU1PDw8sGfPHrULXBpNk+rcvHkT9+7dUwwhYcZn1KhRmDFjBl599VX06tXL0MVhjLEqR6WA4MiRI5g2bRqOHTuGAwcO4MWLF3j11VfLnZQiLi4OY8eOxZtvvomEhAQEBAQgICCgzPHu6tA0qY5EIkGjRo30M6yDqW3JkiXYs2ePwSb+YIyxqkylgGDv3r2YPHkyOnTogE6dOmHTpk24efMmTp8+XeZzvvvuOwwaNAgff/wx2rVrhwULFuCll17C6tWrNS68nDEl1WG6Y2FhodSfITo6mtMbM8aYlmjUhyAjIwMAUL9+/TK3iY+Px4ABA5SWDRw4EPHx8WU+Jy8vD5mZmUqP8miSVGf79u147bXXsGnTpnKPwSqmzpBPdb377rvo16+fXrJ3McZYdaB2QCCTyfDRRx+hd+/eZVbVA6LDXvG89I0aNVJ05CtNSEgI6tSpo3g0bdq03LKYm5tj2bJliIqKQkBAgNIog4CAAERFRWHp0qWl9pY/fvw4IiIicO7cuQrOmJUnIiICbm5u8PX1xbhx4+Dr6ws3N7dyR3dows3NDQAwY8aMEjnJGWOMqU7tgGDatGlITk7Gr7/+qs3yABBpIDMyMhSPW7duVfgcdZPq8JBDzWky5FNdM2fOxMiRI/HixQsMGTKk3BonxhhjFVMrw8v06dMRFRWFv/76q8Lc/46OjiUyzN2/fx+Ojo5lPsfKygpWVlYql0udpDrygICHHKqn+JBP+SgP+ZDPgIAAzJ49G8OHD9dqPgOJRIJNmzbh4cOHiImJwcCBA7F37154eXlp7RiMMVadqFRDQESYPn06duzYgcOHD1fqItqrVy8cOnRIadmBAwd0NnRM1aQ6XEOgGU2HfGrCxsYGUVFR8PX1RVZWFgYOHIijR49q/TiMMVYdqBQQTJs2DVu2bEFoaChq166Ne/fu4d69e3j+/Llim4kTJyIoKEjx94cffoi9e/di2bJluHTpEr766iucOnUK06dP195ZqCkrKwuPHj0CADRv3tzApTFNmg751JQ8KOjXrx9ycnIq1bzEGKs6Jk+eDIlEgkWLFiktj4yMhEQiMVCpTJNKAcH333+PjIwM+Pj4wMnJSfH47bffFNvcvHlT6cvfy8sLoaGhWLduHTp16oSwsDBERkaW2xFRX+S1A/Xr14ednZ1hC2OijGHIZ61atbBr1y7s2bMHY8aM0dlxGGPGydraGosXL8aTJ08MXRSTpnKTQWmPyZMnK7aJiYkpMYRv1KhRuHz5MvLy8pCcnIwhQ4Zoo+waS09Ph5WVFTcXaECTIZ/aVKtWLQwcOFDx9507d7j5gLFqYsCAAXB0dERISEiZ2+zduxd9+vRB3bp1YW9vDz8/P8XEdgAQFhYGDw8P1KxZE/b29hgwYIBS0r3y1ru4uGDFihVKx/P09MRXX31VbnlsbGyUvjeTk5MhkUjw8OFDFV8B7agScxmoq3///sjJycHhw4cNXRSTpcmQT125f/8++vXrh1dffRW///673o7LWJWTnV32Ize38tsWaVYud1s1mZubIzg4GKtWrcLt27fLOJVszJw5E6dOncKhQ4dgZmaGESNGQCaTIS0tDWPHjsWUKVNw8eJFxMTEIDAwUJH4rKL16khISIC7u7tS36vExEQ4OzsbLBtrtZ9H1szMDHXq1DF0MUyafMjnrFmzlHr5u7q6ljvkU1fq1KkDV1dX/PPPPxg9ejTCw8OxZs0aTnnMmKpsbcteN2QIsHt34d8ODkBOTunb9u0LxMQU/u3iApR2F6zBBXbEiBHw9PTEl19+iZ9++qnE+tdee03p7w0bNqBhw4a4cOEC8vPzUVBQgMDAQEV/Mg8PD8W2aWlp5a5XR2JiIjp16qS07OzZsyWW6VO1riFg2hMYGIgrV64gOjoaoaGhiI6ORkpKit6DAUC0J/7xxx/44osvYG5ujt9//x0dOnRAZGSk3svCGNOfxYsX4+eff8bFixdLrEtJScHYsWPRokUL2NnZKZqKb968iU6dOqF///7w8PDAqFGjsH79eqX+CBWtV0dCQgI6duyotKy0IEFOl5lf5ap1QDBhwgSMHDmy1A8PU52qQz51ydLSEvPnz8fx48fRoUMHpKenY8SIEXj99deVRsUwxsrx7FnZj/Bw5W3T08ve9s8/lbe9fr307TT08ssvY+DAgUoj3eT8/f3x+PFjrF+/HsePH8fx48cBAPn5+TA3N8eBAwfw559/on379li1ahXatGmD1NRUAKhwvZmZWYnmgxcvXpRZzuzsbFy9elXp4i+TyZCQkKC0bNiwYXjvvffQrVs3bNy4Uf0XppKqdUCwZ88ehIeHl+gMx6qOLl264PTp0/j0009hZmaGe/fuwdra2tDFYsw02NiU/Sj+f1TetjVrVm5bLVi0aBF27dqllL300aNHuHz5MubOnYv+/fujXbt2Je7wJRIJevfujfnz5yMhIQGWlpbYsWNHpdY3bNhQaXRdZmamIlgoTWpqKmQyGdq2batYtm/fPjx69EgpIEhKSkKbNm1w8uRJvPXWW+q/KJVUbfsQZGRkKD4QnIOgarOyssKiRYsQEBAAZ2dnxdjktLQ03L9/H56enoYtIGNMazw8PDB+/HisXLlSsaxevXqwt7fHunXr4OTkhJs3b+Kzzz5TrD9+/DgOHTqEV199FQ4ODjh+/DgePHiAdu3aVWp9v379sGnTJvj7+6Nu3bqK5sqy2NvbQyKR4OTJkxgyZAiOHTuG6dOnw9raGq1btwYg8uRIpVJ8+OGHuniZSlVtawjkOQgaNGgA2/I6zrAqo2fPnmjWrJni76+++gqdO3fG2LFjceXKFQOWjDGmTV9//bVSza+ZmRl+/fVXnD59Gu7u7pgxYwa++eYbxXo7Ozv89ddfGDJkCFq3bo25c+di2bJlGDx4cKXWBwUFoW/fvvDz88PQoUMREBCAli1bllk+JycnLFiwABMmTEDz5s2xdu1ajBo1Cu7u7opA4vz583pPxS4hE5hQPjMzE3Xq1EFGRobWEgjt3LkTAQEB6Nq1K06ePKmVfTLTQUSYNGkSNm/eDACwsLDAm2++iS+++ALOzs4GLh1j+pWbm4vU1FS4urpyk5qRWL9+PdLT0/H5559Xavvy3sPKXkOrfQ0BJyWqniQSCX755RckJCRgyJAhKCgowA8//ICWLVvik08+MVhiEMYYA0T/AU2HNqqq2gcEPMth9ebp6Yndu3fjr7/+Qu/evZGbm4tvvvmmRNYxxhjTp5UrV2LYsGF6PWa1DQjy8vI4bTFT8Pb2RmxsLKKiotC3b1/MmDFDsS4lJUUxCRZjjFVV1bYPASDGfRYUFMDS0lJr+2RVCxHBx8cHCQkJ+OCDDzBz5kzUr1/f0MViTKu4D4Hp4z4ERUilUsTExGDbtm2IiYmpVFYnMzMzDgaMnDrvqzY9ffoUGRkZyMrKwsKFC+Hi4oJPPvlEZ9M5M8aYoVSJgCAiIgJubm7w9fXFuHHj4OvrCzc3N0RERBi6aEwDxvC+1qtXD2fOnEFERAQ6duyIrKwsfPPNN3B1dcU777xTbvIRxhgzJSYfEERERGDkyJHw8PBQmmnPw8MDI0eOLPXikZycjG7duuG///2vAUrMKkOd91VX5LOiJSYmYteuXfDy8kJeXh5++OEHxMXF6a0cjDGmSybdh0AqlcLNzQ0eHh6IjIxUmkZSJpMhICAAycnJSElJUcoatWPHDgQGBqJ79+6KfNbMeKj7vuoLESE2NhYbN27E+vXrYWEhEn5GRETA1tYWr7zyiiIbImOmgPsQmL5q34cgNjYW169fx5w5c5QuGoC4qwsKCkJqaipiY2OV1nEOAuOm7vuqLxKJBC+//DI2btyoCAYKCgrw0UcfYeDAgWjfvj3WrFmDrKwsg5SPMcbUYdIBgbxjl7u7e6nr5cuLdwDjHATGTd331ZBycnIwYsQI1K5dG5cuXcL06dPRpEkTfPTRR5wWmTFmEkw6IHBycgIg+gSURr5cvp0c1xAYN3XfV0Oys7PDd999hzt37mDVqlVo3bo1MjMz8d1336F169ZYunSpoYvIGGPlMumAwNvbGy4uLggODi4xhbFMJkNISAhcXV3h7e2ttI4DAuOm7vtqDGrXro3p06fj4sWL+PPPPzFkyBBFLgO5J0+eID8/33CFZIyxUph0QGBubo5ly5YhKioKAQEBSr3RAwICEBUVhaVLlyp1PCMixVAxDgiMkzrvq7ExMzPDoEGDsHv3bty4cQNdu3ZVrAsKCoKrqysWLVpUYk52xphqJk+eDIlEgkWLFiktj4yM5M69KjLpgAAAAgMDERYWhqSkJHh5ecHOzg5eXl5ITk5GWFgYAgMDlbbPycmBk5MTrKys0Lx5cwOVmlVE1ffVmBWdcrmgoAAHDhzA3bt3ERQUhCZNmuCTTz7hwIAxDVhbW2Px4sX8f6Qhkw8IAHHxuHLlCqKjoxEaGoro6GikpKSUetGwsbHB5cuXkZOTg5o1axqgtKyyVHlfTYWFhQUuXLiAn3/+GR07dkROTg6++eYbtGzZEsuXL0deXp6hi8iYyRkwYAAcHR0REhJS5jZ79+5Fnz59ULduXdjb28PPzw9Xr15VrA8LC4OHhwdq1qwJe3t7DBgwANnZ2ZVa7+LiUmJCNE9PT3z11VfllsfGxkapWTQ5ORkSicRgs61WiYAAENXMPj4+GDt2LHx8fCqsTi4+nI0ZJ1XfV1NgZWWFiRMnIjExEVFRUejQoQOePHmCmTNnlvuFxhgrnbm5OYKDg7Fq1Srcvn271G2ys7Mxc+ZMnDp1CocOHVIkHJPJZEhLS8PYsWMxZcoUXLx4ETExMQgMDIQ8TU9F69WRkJAAd3d3pWtRYmIinJ2d0aBBA7X3qwkLgxyVMQaJRIKhQ4di4MCB2LRpE7799lu8//77ivXPnz/nWixmUEVukEswNweK5r8pb1szM6DoR7msbW1sVCtfUSNGjICnpye+/PJL/PTTTyXWv/baa0p/b9iwAQ0bNsSFCxeQn5+PgoICBAYGKpqSPTw8FNumpaWVu14diYmJ6NSpk9Kys2fPKi3Lz8/HtGnTcPToUdSpUwe7d+/W6eRq1e42ee7cuejWrRu2bt1q6KIwPTL0JEnlsbCwwFtvvYXz58/D3t4egOj86u/vj/Hjx/PUy8xgbG3LfhS7vsLBoextBw9W3tbFpfTtNLV48WL8/PPPuHjxYol1KSkpGDt2LFq0aAE7OztFp/KbN2+iU6dO6N+/Pzw8PDBq1CisX79eqT9CRevVkZCQgI4dOyotKx4kzJ8/H3369MGFCxfQv39//PrrrxodsyLVLiA4e/YsTp06hWfPnhm6KExPjGGSpMoo2iP64sWLir4THTp0wM6dOw1YMsZMw8svv4yBAwciKCioxDp/f388fvwY69evx/HjxxVp6/Pz82Fubo4DBw7gzz//RPv27bFq1Sq0adNGMSKtovVmZmYlmg9evHhRZjmzs7Nx9epVpYu/TCZDQkKCYllGRgb++usvTJo0CYBIpHft2jUNXp2KVbuAgHMQVC/GNEmSKtq3b4+4uDi0bdsW9+/fR0BAAF5//XU8fvzY0EVj1cizZ2U/wsOVt01PL3vbP/9U3vb69dK304ZFixZh165diI+PVyx79OgRLl++jLlz56J///5o165diTt8iUSC3r17Y/78+UhISIClpSV27NhRqfUNGzZUypyamZlZ7kyoqampkMlkaNu2rWLZvn378OjRI0VAcPDgQVy5cgWenp7w9PTE3LlzddpcAFSzgIBzEFQvUqkUs2bNgp+fHyIjI9GzZ0/Y2tqiZ8+eiIyMhJ+fH2bPnm1UzQdF9ejRAwkJCfjkk09gZmaGLVu2oEOHDti1a5ehi8aqCRubsh/F50Aqb9viXWHK2k4bPDw8MH78eKxcuVKxrF69erC3t8e6detw5coVHD58GDNnzlSsP378OIKDg3Hq1CncvHkTERERePDgAdq1a1ep9f369cPmzZsRGxuLpKQkTJo0qdwO0Pb29pBIJDh58iQA4NixY5g+fTqsra3RunVrAKI2e8mSJUhMTERiYiI8PT1L9DnQOjIBGRkZBIAyMjI02s+DBw8IAAGg58+fa6l0zFhFR0cTAIqPjy91fVxcHAGg6Oho/RZMDfHx8dS2bVsCQK1ataK8vDxDF4lVIc+fP6cLFy6Y5PfipEmTaPjw4UrLUlNTydLSkope4g4cOEDt2rUjKysr6tixI8XExBAA2rFjB124cIEGDhxIDRs2JCsrK2rdujWtWrVK8dyK1mdkZNDo0aPJzs6OmjZtSps2baJOnTrRl19+WWa5//e//1G9evWoWbNmNGnSJPr000+pa9euivXvv/8+/fbbb4r9N2rUiLKzs8vcX3nvYWWvoSY9/bGqTp06hW7dusHJyQl3797VYgmZMdq2bRvGjRuHrKws2JbSYykrKwt2dnYIDQ3F2LFjDVBC1Tx//hxffvklhg8fjt69exu6OKwK4emPjc/KlStx8+ZNLF26FLNnz0b9+vUxZ86cMrev9tMfq4r7D1QvpjhJUnlq1qyJJUuWKAUD69atw7fffqvReGjGmPEZN24cYmJi0KpVK+Tl5eHTTz/V+TGrVR4CmUwGNzc3RRsNq9qKTpIUGRmplADE2CdJqoxr167hgw8+QF5eHqKjo7Fp0ybFsEXGmGlr0KABTp06pddjVqsagv/85z9ISUnBpk2bDF0UpgdVYZKk8ri6umL58uWwsrJCVFQUPD098ffffxu6WIwxE1WtAgJW/VSlSZKKk0gkePfdd3Hs2DG0bt0at2/fho+PT6nTRjPGWEWqVadCVn1JpVLExsYiLS0NTk5O8Pb2NtmagdJkZWXh3XffVWTg9Pf3L9FMwlhZuFOh6dNGp8Jq04eAiNCiRQs0bNgQf/zxBxwdHQ1dJKZH8kmSqqratWtj8+bN6NevH6ZNm4bu3btzMMAYU0m1CQgePnyI69ev48aNG6hXr56hi8OY1kkkEkyZMgV9+vSBm5ubYnl2djZstJX1hTFWZVWbWwj5kENnZ2dYWVkZtjCM6VDr1q0VtQM5OTnw8vLC9OnTkZ+fb+CSMWNnAi3IrAza6DdUbWoIOGUx0zVj7Kewb98+nDt3DufOncOZM2ewfft2NG7c2KBlYsanRo0akEgkePDgARo2bKg00RYzbkSE/Px8PHjwAGZmZrC0tFR7X9UmILh58yYAKOayZkybIiIiMGvWLEVNFCCCz2XLlhl0JMOIESOwa9cuvP7664iPj8dLL72EX375BQMHDjRYmZjxMTc3R5MmTXD79m2lzzAzHbVq1UKzZs006jtUJQKCytyZyWeicnZ2NkQRWRUmn1HRz88P27Ztg7u7O5KTkxEcHIyRI0cafHijn58fTp06hddeew1nz57FoEGDMHPmTAQHB3PzGVOwtbVFq1atyp22lxknc3NzWFhYaF6zU+5MB0aivIkZwsPDycXFRTFpEQBycXGh8PBwpe3Gjh1LAGjZsmX6KjarBgoKCsjFxYX8/f1JKpUqrZNKpeTv70+urq5UUFBgoBIWysnJoWnTpin+T9577z1DF4kxpgeVndzIpDsVqjLXvYODA1q1asVNBkyrYmNjcf36dcyZM6dEVZ2ZmRmCgoKQmpqK2NhYA5WwUM2aNbF69Wrs3LkTrVu3RlBQkKGLxBgzIiabmEgqlcLNzQ0eHh6l5qkPCAhAcnIyUlJSDN6xi1VdpjqjolQqVfq/2LRpE4YPH85Dchmrgqr8bIemdGfGqi5TnVGxaDDwxx9/4I033kCnTp1w+PBhA5aKMWZIJhsQyDsJuru7l7pevly+HWO6UHRGxeLjgE1lRkUnJye4ubnh1q1b6N+/P2bMmIHnz58buliMMT0z2YBAlTuza9euwcnJCX369NFb+Vj1UBVmVOzWrRsSEhLw9ttvAwBWrFiBrl274syZMwYuGWNMn0w2IFDlzuzu3bu4d+8e1xYwnagKMyra2tpi7dq12L17NxwdHXHhwgX06NEDK1euNHTRGGN6YrIBgSp3ZvJAwNjacVnVERgYiCtXriA6OhqhoaGIjo5GSkqKSQQDRQ0ZMgRJSUl47bXXUFBQgKZNmxq6SIwxPTHpxETyO7NZs2bBy8tLsdzV1VXpzowDAqYPup5RUV+pkRs0aIDt27fj6NGjSs1sx44dg4eHB0+UxFgVZdIBASCCguHDh5f7RckBATN1+k6NLJFIlIKBBw8eYPDgwahduzaWLl2KUaNGcb57xqoYk20yKEp+ZzZ27Fj4+PiUmbaYAwJmilRJwKUrN27cQJ06dXDr1i2MHj0avr6+OHfunM6PyxjTnyoREFSEAwJmqqRSKWbNmgU/Pz9ERkaiZ8+esLW1Rc+ePREZGQk/Pz/Mnj0bUqlUp+Xo2rUrLl68iPnz58Pa2hpHjhxB586d8f777+Px48c6PTZjTD9UDgj++usv+Pv7w9nZGRKJBJGRkeVuHxMTA4lEUuJx7949dcussqZNm6JNmzZo1qyZ3o7JmDYYUwKumjVr4osvvsClS5cwcuRIyGQyrF69Gq1bt+a8BYxVASr3IcjOzkanTp0wZcoUldouL1++rJQy0cHBQdVDq+3HH3/U27EY0yZjTMDVvHlzbN++HYcPH8bMmTPRpUsX1KxZU7H+3Llz8PDw4D4GjJkYlQOCwYMHY/DgwSofyMHBAXXr1lX5eYxVZ0UTcPXs2bPEekOmRu7Xrx8SEhKQk5OjWJaYmIjOnTujW7dumDFjBl577TVYWlrqvWyMMdXprQ+Bp6cnnJyc8Morr+Do0aPlbpuXl4fMzEylB2PVkbGnRpZIJErDEBMSEmBlZYWTJ09i3LhxaNiwIUaPHo0tW7bg0aNHBikjY6xydB4QODk5Ye3atQgPD0d4eDiaNm0KHx+fctOihoSEoE6dOoqHJslRTp06BUdHR7VqNRgzNFNLjfzGG2/g5s2bmD9/PhwdHZGZmYnff/8dr7/+OhwcHJCUlGToIjLGyqDR9McSiQQ7duxAQECASs/r27cvmjVrhs2bN5e6Pi8vD3l5eYq/MzMz0bRp0wqnbizNzp07ERAQgG7duuHEiRMqPZcxY1FaHgJXV1csXbpUK3kIdJH0SCaT4cSJE9i1axd27dqlSB8u329QUBDu37+P/v37w9fXF87OzhqfB2OspMpOf2yQxETdu3fH33//XeZ6KysrWFlZaeVY8tEMjo6OWtkfY4ZQmQRc6tJV0iMzMzP07NkTPXv2xMKFC/HkyRNFeYkIW7duxa1bt7Bx40YAQLt27dCvXz/069cP/fv3R506dTQ6L8aYagyShyAxMVFvnaA4BwGrKipKwKUOfSY9qlevnuJ3IsL69evx8ccfo0uXLpBIJLh48SLWrFmD1157jWcmZcwAVK4hePbsGa5cuaL4OzU1FYmJiahfvz6aNWuGoKAg3LlzB7/88gsAMZWqq6srOnTogNzcXPz44484fPgw9u/fr72zKAcHBIyVrnjSI3meA3nSo4CAAMyePRvDhw/Xeh8FMzMzDBw4EAMHDgQAPH78GEeOHMGhQ4dw4MABvPrqq4pt8/Ly0KVLF7z88ssYOnQoXnnlFR65wJgukIqio6MJQInHpEmTiIho0qRJ1LdvX8X2ixcvppYtW5K1tTXVr1+ffHx86PDhwyodMyMjgwBQRkaGqsUlf39/AkBr165V+bmMVWXy/+X4+PhS18fFxREAio6O1m/BiCg/P1/x+/79+5W+a+rXr09vv/02xcbGklQq1XvZGDM1lb2GatSpUF8q2yGiNN26dcOpU6ewc+dODBs2TEclZMz0bNu2DePGjUNWVhZsbW1LrM/KyoKdnR1CQ0MxduxYA5RQyMnJQXR0NHbv3o3IyEilJEwuLi744YcflGoUGGPKKnsNrfJzGbRs2RJt27bled0ZK6Zo0qPSGDLpUVG1atXC0KFD8X//93+4desWDhw4gEmTJsHW1hbXr19HkyZNFNump6cjNzfXgKVlzHRV+RoCxljppFIp3Nzc4OHhodSHABBDBgMCApCcnIyUlBSjyXNQlLzmYOjQoYpl48ePx/79+/H222/j3XffRePGjQ1YQsaMA9cQMMbKZWpJj4qT1xzI5efn4/jx43j48CEWLlwIFxcXjBkzBnFxcTCB+x7GDI4DAsaqscDAQISFhSEpKQleXl6ws7ODl5cXkpOTERYWppWkR/piaWmJS5cuISwsDC+//DIKCgrw22+/oXfv3ujevTuioqIMXUTGjFqVDgj279+PRo0aYfTo0YYuCmNGKzAwEFeuXEF0dDRCQ0MRHR2NlJQUkwoG5CwsLPDaa6/hyJEjSEhIwJQpU2BlZYVTp07hwoULhi4eY0bNIJkK9eXOnTtIT09HRkaGoYvCmFGTJz3SF12kSi7O09MTP/30ExYtWoS1a9fiv//9r2Ldrl27EB8fj/fff9/gnSYZMxZVuoaAkxIxZnwiIiLg5uYGX19fjBs3Dr6+vnBzc9NqVsSiGjZsiHnz5immXyciLFiwACEhIXBxccHUqVNx9epVnRybMVPCAQFjTG/0mSq5PHPmzIGXlxfy8/Px448/ok2bNpg0aRIuX76sl+MzZoxMOiCQSqWIiYnBtm3bEBMTA6lUqrSeJzZizHgUT5Xcs2dP2NraKlIl+/n5Yfbs2SX+j7VNIpEgICAAR48exdGjRzF48GBIpVL88ssvaNeuHb766iudHp8xY2WyAUFlqh25hoAx4xEbG4vr169jzpw5SjkPADG3QVBQEFJTUxEbG6u3Mnl5eWHPnj04efIkhg0bBiJCp06d9HZ8xoyJSQYEla125ICAMeMh/390d3cvdb18edHUxPrStWtX7Ny5E0lJSRg+fLhi+YoVKzBhwgRcu3ZN72ViTN9MLiBQpdqxXbt2aNu2LWcrY8wImEKqZHd3d0XtRW5uLhYuXIitW7eibdu2eP/993H//n2DlY0xXTO51MVnzpyBr68v4uPj0bNnzxLbxsfHw8vLC9HR0XodRsUYK58ppko+ffo05syZo5iu3cbGBjNnzsTs2bM5jTozGVU2dbExVzsyxsqmz1TJFXU4rqwuXbpg3759OHToELp3747s7GwsWLAALVq0wB9//KFxORkzJiYXEJhCtSNjrHT6SJWsizwH/fr1w7FjxxAeHo42bdrgyZMnaNGihcZlZcyYmFxA4O3tDRcXFwQHB0Mmkymtk8lkCAkJgaurK+7evQsHBwe89dZbBiopY6w0ukyVrMs8BxKJBIGBgUhOTsbhw4eVaimXLVuGmJgYjcvPmEGRCcjIyCAAlJGRQURE4eHhJJFIyN/fn+Li4igzM5Pi4uLI39+fJBIJhYeH09KlSwkAjR071sClZ4zpQ0FBAbm4uJC/vz9JpVKldVKplPz9/cnV1ZUKCgq0etyLFy+Subk5AaChQ4dSUlKSVvfPmKaKX0PLYnI1BEDlqh15yCFj1Yuh8hzY29vj3XffhYWFBXbv3o1OnTphypQpuHnzplaPw5iumWRAAFRc7cgBAWPVi6E6HDds2BCrVq3ChQsXMHLkSMhkMmzcuBGtWrXCjBkz8OTJE60ejzFdMdmAACicoW3s2LHw8fFR6p3MAQFj1YuhOxy3atUK27dvx7Fjx+Dr64v8/Hxs3LgRZPwjuxkDYOIBQXk4IGCseqlsh2Nvb2+19l/ZoYw9evTAoUOHsG/fPqxatQr169cHIGZZ/O2335CXl6fW8RnTtSobEMgnNuKAgLHqQZd5DlQdyiiRSPDqq6/i9ddfVyzbu3cvxowZg3bt2iE0NLRE0MKYoVXJgEAqleKll15Cu3btOCBgrBrRRZ4DbQ1lzM3NhZOTE1JTUzF+/Hh07doVBw4cULk8jOmKyaUu5nShjLGKSKVSxMbGIi0tDU5OTvD29larZkDb6Zazs7OxYsUKLF68GFlZWQCAV155BYsXL0bnzp1VLh9jlVHZaygHBIwxVoaYmBidzJ3y8OFD/O9//8P//d//4cWLF2jXrh2Sk5NLDJdkTBuq7FwGjDGmL7oaytigQQOsWLECly9fxrhx4xAcHKwIBvLy8vDw4UMNSs2YeqpkQPB///d/cHBwwMyZMw1dFMaYCdP1UEZXV1ds3boVAQEBimXff/89WrZsiZCQEOTk5Ki1X8bUUSUDgjt37uDBgwfIz883dFEYYyZM10MZgZLDGaOiopCZmYk5c+agdevW+Omnn9SerZExVVTJgIBzEDDGtEHXUzaXNpzxypUr+OCDD9C8eXPcuXMHb731Fjp16oT9+/dr+ewYU8YBAWOMlUNXUzaXNZyxY8eOWLVqFRYtWoRly5ahXr16OH/+PAYOHIj58+dr+ewYK1QlRxl4enri7Nmz2LNnDwYPHqyHEjLGqjptDWWU76uywxkzMzMxf/58rF27FidPnoSHh4e2TolVE9V62GGjRo2Qnp6OhIQEeHp66r6AjDGmAnWGMz548AANGzZUbPP111/DxcUFEyZM4OGKrFzVdthhQUEBHjx4AICbDBhjxkmd4YxFg4ELFy7g66+/xqRJk9C7d2+cOnVKh6Vl1UWVCwiys7PRr18/eHh4oEGDBoYuDmOMlaDpcMaWLVti4cKFsLW1xbFjx9C9e3dMnTpVcTPEmDqqZJMBY4wZM22lRL579y4+/fRTbNmyBQBQt25dfP3113j33XdhYWFR4pja6gPBTEu1bTJgjDFjp63hjM7Ozti8eTP+/vtveHp64unTp5g7dy4eP36stJ2qszWy6qnKBQQmUOHBGGNaHc4o70fw/fffY9myZXBwcFCs27Bhg1Zma2RVn8k2GZRV/bVgwQKsXLkS06dPx5dffmngkjPGWPl0WZW/a9cuDB8+HK1atcKpU6dQu3ZtxTp1ZmtkpqmyTQYWZa4xYhEREZg1axauX7+uWObi4oJly5bh7t27ePjwYYk0o8x4EAF37wLnzwPJyeLnhQtAjRrAX38VbhccLJbXqgW0aAF06iQeTk6ARGK48jOmTebm5irNlKiKtWvXgojwzz//oH379liyZAnGjBkDiUQCMzMzBAUFwcvLC7GxsTorAzMdJhcQyLN7+fn5Ydu2bXB3d0dycjKCg4MxcuRIdO3aFQDg6Oho4JKy0kyeDOzcCTx9WnKdvb3y3wcPAtHRJbdr0ADo3BnYswewMLlPMGP6M378eOzZswfNmjXDzZs3MW7cOKxevRrfffcdunbtqvZsjaxqMqmvU6lUilmzZsHPz0+pZ27Pnj0RGRmJgIAARb5vzkFgeBcuALt3A7NnF97R5+aKYMDcHGjVCnB3Bzp0EI/69ZWfP20a4OcHZGUBly4BiYnAP/8ADx8CV68qBwPTpgEyGTBgAODrW3JfjFVHzs7OAICff/4ZR48eRXBwMOLi4tC9e3d88cUXGDhwIAD+vmSCSfUhiIqKgp+fX4XZvQDg2LFj6NGjh76LWu3dvg2EhorH2bNi2cmTwL8VNzh/HpBKgTZtACsr1fefkyP28fQp8MorYtmLFyIAePZM/C2RAF26iOBg0CDAy0s0RzBW3RQf3piWlobPPvsMW7ZsQUREBDZu3Mh9CKqBKpm6+Mcff8Rbb72FrKws2NralrkdANy4cQPNmjXTd1GrJSIgPh5YsQKIiBAXfEBchAcPBr78EnjpJd0d/8ULYO9e0cRw8KComShq4ECxnrHqqGgza1BQENzd3bFz50789ttv2L17N8LCwvDkyRPk5+fjrbfeQg2OnqucKtmpUN4vIDk5udQagqNHj5bYtrqSSoHYWCAtTXTC8/YW1fS6eN6hQ4V364B4zuuvA6+9pp+q+xo1AH9/8QBEh8VDh4D9+0UgULSv1MOHwKuvikffvkCvXkDdurovI2OGIh/eOGvWLEUNKgC4uroiLCwMvr6+aNmyJZ48eYLvvvsOixcvxrBhwyBRo+cuJz8ycWQCMjIyCAA9fvyYXFxcyN/fn6RSqdI2UqmUBgwYQNbW1uTl5WWgkhqH8HAiFxcice8uHi4uYrk2npeeTvTXX4V/FxQQtW9PNGUK0dmz4u/oaKLQUPGzoEDbZ1h5UilRTk7h31u3Kp+fRELUsSPRe++J8t6/b7iyMqZLBQUFFB0dTaGhoRQdHU0F//5j5ufn0+rVq6lBgwYEgACQl5cX7d27l2QyWaX3Hx4eTi4uLop9ACAXFxcKr+iLh+mc/BqakZFR7nYmFRBkZGRQeHg4SSQS8vf3p7i4OMrMzKS4uDjy9/cniURS7T984eHiIufvTxQfT5SVJX76+4vlZb08FT0vLEwEAWPHEllaEjk6EuXmFj7/xYvC/agTjOjLw4dEW7YQTZ5M5OamXE5ArJM7fZpoyRKiiAiipCTlwIKxqiYjI4PmzJlD1tbWigt6t27d6NSpUxU+t+j3cnx8PGVlZVF8fDx/LxuJKhsQEJUeibq6ulapD506d9kFBeLi6+8v7oyLkkrFclfXkvsq73lPnhC5uxPVqKF84ezalejaNeVt1Q1GtHX+6khLE8HORx+Jc7p+vXDdokUlAwYHB6K2bYm8vIgSEwu3PXWKaPlyonXrRJl37iQ6eJDo2DERTDx7ppvyM6Ztd+/epRkzZlDNmjVJIpHQxYsXy92+oKCg3Jpbf39/cnV1VdRIMP2r0gEBUenVX6pUbxkzde+yo6PFtvHxpa+PixPro6Mr97ytW4lsbArLYGVF9OabRCdPlty3usFIaYylliEyUtSIdO1KVKdOyeAgIaFw2+DgkuuLPg4dKtz2t9+IevUiGjaM6K23iBYuFEFJcrJyrQtjhnT//n3aUrTKjIi+/PJL2rJlC+Xl5SmWRUdHEwCKL+OLJy4ujgBQdPEvHqY3VT4gKM2nn35K9vb2FBISoqeSaZ8md9mhoeLik5VV+vrMTLE+NLT05x05Ii5KcmfPiuVt2oif69eXfWx1g5HitFnLoE0ymeg7cfasOIfwcPF6yoWFEY0ZI8rZrx9Rjx5EHTqIQKZBA+XX9X//KztwMDNTDh7u3iW6ckUcnzFD+ueff8jMzIwAkKOjI82fP5/u3btHoaGhBICyyvjiyczMJAAUWvyLh+lNtQwIJk2aRABMNiDQ9C5b1YuyTCaqvceNK7wgTZig/JyTJ4mOHq34Yq5uMFKUNmsZjKljY3FXroiAYu1aovnziV5/nahbNyI7O/EaXb1auO2CBWJZgwZEfn7i7wMHiJ4+NVz5WfX05MkTWrBgATk5OSmaai0tLenVV1/lGgIjVy0DgsGDBxMA+umnn/RUsoqpcmHS9C67shfUQ4eIZswgatGi5N3plCllP0+XZdfWPoiMp8lBVTKZ6NNQ9L2bNUt04ixekyCREHXpQnTjhuHKy6qnvLw8Cg0NpR49eij143rppZe4D4GRqpYBQZcuXQgA7dq1S08lK5+qFyZt3GUXrXKPixPPiYlRrnLv1Em5X8CIEUQzZ4q/iz4vLq7yVfXauLvX9vkbU5ODJnJzRefEFStEnwZXV/E62NoS5ecXbvftt0RffCHe7yJNvIzpzLFjx2jcuHFka2tLABSjv/bt20ebN2/WySiDsoZPsrJV2YCgvA9D06ZNCQAdP37cgKUV1LkwafMOuWlT5UCkefPCYy5bRjRxovi76MW3tADG1bXyF9HSghFVggp91ZBo+v1hDM0Rd+4o9zUgImrduvB9s7ERTQwrVxJdusR9EJhuZWdnlzr6y8rKit577z16pqVhNpzrQD1VMiDYvHlzmR8GmUxGVlZWBIBSU1MNWl5dDP+rzAXt3j3RLv3KK0QWFsoX9p07K192TS52mgQV+u5DoQ5jbY6QyYh+/FH0B3FwKNnE8PLLhi0fqx4KCgrozz//pG7dupG5ubnie9rW1pbGjRtHERERlKNmQg/OdaC+KhkQyKukSvsw/PLLL4ptsrOzDVpeTS5M6t5lh4aKPgBFLwIdO4pOa8nJ+r1D1CSo0KSWQRtNDpUtmzE3R0ilYkjk4sVE/fuLPgjvvFO4/sULsfyLL0SHUXlSKca06d69e7Ro0SJyc3NTuombULznciVwrgPN6CwgOHLkCPn5+Sl6mu7YsaPC50RHR1Pnzp3J0tKSWrZsSRs3blTpmPKTGTRoUJkfhmbNmlH//v2pR48eKu1bFzS9MFV0l33pkhj3vnt34XOuXStMGBQSQvTPP9o9J31St5ZBlzUE+mqO0IXsbOWUzH//rfza1q1L9NprYljpzZuGKyermqRSKR09epRmzpxJzZo1o7CwMMW65ORkGjZsGC1dupTi4uIot4xEHJzrQDM6Cwj27NlDn3/+OUVERFQqILh27RrVqlWLZs6cSRcuXKBVq1aRubk57d27t9LHlJ/MwYMHS11vbB8GbVyYit5l799PdPgw0Zw5Ys4A+Rf5sGHKz7lzR0snYAT0mamxMvTRHKEvDx6I5oVRo4jq1SvZvPDdd4YuIauqZDKZ0l38l19+qVR7YGlpSb169aJZs2ZRWFgYPf13fC3nOtCMXpoMKhMQfPLJJ9ShQwelZaNHj6aBAwdW+jjyk7lTxhVP3x+Gii5W2rowyWREAQHK2QIBkUZ40CCiTZu0eVZVg6YdG8ui6+aIsui6A2NBgRjBMH++SMdsZqYc9PzxB9HgwSJI4M6JTNvOnz9PISEhNGzYMKXJleSPY8eOEVFhDcGnn36qmMOmKGO7KTQ2RhMQeHt704cffqi0bMOGDWRnZ1fmc3JzcykjI0PxuHXrltHUEFS2U5mqF6b790Wq4K++Ul7ep484hoOD6DC2ZQsnpamIpqMlSmOIGgJDdGB8/Fg56Jg6teTx335bTPjEn0OmTTKZjFJSUujnn3+md955h7p27apIkVxQUKAY2ih/NG/enIYMGUKzZ88mT09Pat68uUZ9CKrycEajCQhatWpFwcHBSst2795NAMrsbVq8Gkn+KK8PQd26dal+/fo0d+5cTU6pXKp2KivvwvT8ucg49/HHRJ6eyglnHj0q3MfRo6KDWPGaBlY+bd9Z67sPgbY7MKr7ely4IGZ8lHdOLPpZNjcXiZQY04epU6eSfChjadeHzZs3K7b9/vvvacaMGfTDDz9QTEwMpaWllTvXTVUfzmjSAUFZNQRFE18Un/ZYnqXw888/1+SUyqTJUMLiX8SffUZkbV2y7dbTk+iTT8TwQWZ8dNUcUZy2gw9t1TQ8eyY6sr7/vsh50Lq18vrRo0XzwjffiKmjq9ANFjMSpV24a9euTb169VLarn///iUChjp16lD37t3pjTfeUOq8WB2GMxpNQKBOk0Fx5eUhkE97/J///IcA0Hc66hGlSZVxRoby9Lc//yy2dXYmmjxZNBUU7QXOjJcumiOK02bzhC6HShZtMsjPF5kTi74u9eqJ5Ej/+x9RbKz6x2GsqMpU7YeGhtKMGTNo8ODB1KJFC5JIJIprhoODg6K2QN4U0ahRI5o7dy5FRETQ9evXSSaTVanhjEYTEHzyySfk7u6utGzs2LFqdSosL1Nh3759CQBt27ZN5fOoDHU6lV28SDR9uviiXL68cHlurpg1jztomSZdd/TTVgdGbdQ0VPZcpVIxUda334ogoHZt5eDg1VeVt9+6lej4ceUZIxnTlefPn1NSUhJt376dNhXpjS3vrFj8Ua9ePerfvz9NmTKlSnRWrGxAYAEVPXv2DFeuXFH8nZqaisTERNSvXx/NmjVDUFAQ7ty5g19++QUA8M4772D16tX45JNPMGXKFBw+fBi///47du/ereqhAQDm5ubw8fEpsTw9PR0A4ODgoNZ+K+LkJH4mJwM9e5Zcn5xcuF1GBjB1KrB9e+H6gweBjz4Sv1tZAR076qSYTA/MzYFSPoJao8pnrTyxscD168C2bYCZmfI6MzMgKAjw8hLblXY+ERHArFliH3IuLsCyZUBgYMn9deokHjNmAAUFwOnTQHw8cOyYOI7c/fvA+PGFfzdpArRrB7RvL3726QN06FD+uTGmCmtra7i7u8Pd3V1p+d27dwEA33zzDS5evIiEhAQkJyfjyZMnOHToEHr06AEASEtLAwC89dZbaNCgAbp3745u3bqhSZMmkEgkpR5TKpUiNjYWaWlpcHJygre3N8zNzXV4llqgaqRRVkQ1adIkIhJTEPft27fEczw9PcnS0pJatGihdmKi8qIbe3t7AkBJSUkqnlHlVPZuKympMKe8REI0fLjoPMi1AayytNWHQJOaBl02NVy4QNSvH5GjY8l+NIDoaCt3966Y1XH4cKJp00RCrp9+EsMhjx3jpjammdISHuXm5tLp06dp/fr19PXXXytqCLKyssjMzEzputeoUSPy8/OjoKAgpdw6xtZJsUqmLi7rZF68eKF40e/r8Buiok5ln39e2I7apAnRiRM6Kwqr4rTRgVHdvgj6bGp4/FiMpPnxRzHj5uDBRL/8Urj+2LHSgwb547PPCre9fZuoZ0+Ru+Pdd4m+/ppo3TqiXbuIzpzhYZKsJFVSImdlZdG6deto6tSp5OnpqTRXAwCaPHkyEYlgAAA5OTnR5MmTaePGjbR582YaMmSIwTopVquAICMjg1555RXy9PTUeeeP8jqVnT1LVKsWUd++fOfCNKdpB0Z1L+yadmrUZv6Ex49FbcD334uAe/JkoiFDiLp1EzN4rl5duG18fPnBw+zZhdump4v5HUJCiH79lejkSeXhvqz6KDrKoLQRbGVdwLOzs+no0aO0cuVKevvtt+nXX39VBBg+Pj4latHNzMzIxsaGatWqpZRETx5s6FJlAwIJEZE+miY0kZmZiTp16iAjIwN2dnaGLg6kUtHumpYGODiI9ld509CpU6IdtUYNgxaRVRFFP2tOToC3d+FnrTIiIoCRIwE/P9FnwN1d9EEICQGiooCwsJL9AbZtA8aNA7KyAFvbkvvMygLs7IDQUGDs2LKPN2dO4fGCg8s+nrbO99Ej4K+/xHPv3RM/5b/fvg3MmwdMmya2jYsDevcuuY+6dYGWLYEPPgAmThTL8vPFvh0dgTKai5mJi4iIwKxZs3C9SIcZV1dXLF26FIHlfWCLiYmJga+vL6KiopCSkoLz58/j4sWLOH/+PJ4+farYbsqUKfjpp58AAMnJyfDw8ICdnR0aN26MJk2awNHREfXr14e9vT369+8Pr3874eTk5OC3337Dw4cP4eDggK5du8LCwgISiQQSiQT29vaoX78+AOD58+f4559/IJPJkJGRAV9f34qvoToNS7SkstGNJtTpOX7+vJhbIC5OZ8ViTGOq1jQYqqlB15kZi/bjuXJF1DhMmEDUu3fJ/gxFax7kzRa2tiJXyH/+I567ebNoilBzNl9mZLSRqbCsORdkMhndvXuX/vjjDwKglJtn//79pfbLkz9CQkKISNRkyCcVLOvx9ddfK/Z77ty5Euu1PsrAGBFRmT09K0OV3tRyjx8DAweKO49Zs4CjR/nugRmnwEBg+PDK33l7e4vPf3AwEBmpPEJBJhO1C66uYruiNBnVULRmYds25ZqFkSPLr1mobK1C0f/Pli2B//1PeX12NpCaCly9Cnh4FC6/fVuU/9kzIDFRPIr69lsxsgIA7twBDh8GOncG2rYFLKrEN2z1UNYINlU4/Tv0Jzk5GT2LDBGSSCRwcnJCgwYNAAC9evVSrHvllVeQmZmJO3fu4Pbt27hz5w7S09Px+PFjPH78GF26dEFERARGjhyJXr16IScnB2ZmZigoKEBubi5evHiBWrVqwdzcHDWKVE3XqFEDjo6OMDMzg0QiwZ07dyo+AZVDIAOoqIZg+fLlVK9ePZoxY4bK+1anN7VMRhQYKO4aWrUSs8cxVpWo06lR3VENmtQs6Gu+h9xckVfkjz+Ili0T8zm8/DKRvT3Rvn2F28lfA4CoZk0xYdQHH4hkZOfPc/bGqk6VTor63Ge16lT42WefEYASGRErou4X0bp1hbMOnjql0iEZMxn6ampQ93m6HBpZWTKZclPEH38QeXuXzNoof0REFG779ClRdrbuy8j0S91OimUpbWhkUZWZ3K9aBQTybFILFy5Uab/qfBFduCAif0DkbGdGSJNUgpqmIdR1GkM9U+V01A2w1alZ0FauBl29XVKpqFHYvJnoo4/ErKW1ahHdvFm4zcKFRBYWRD16iBEQO3cSPXmineMzwyotD4E8zb6qyuqXIJeZmUkAlEYuFFetAoKhQ4cSAFq/fr1K+1X1iyg3l6hTJ7HslVd4BkKdqcy3tExW+hunSR2ypvXPhpiv2Mio09SgTmCujfketPl2VeYjW1CgXJswblzJGgQzM5GI6eOPOW+CqdPWdMpcQ1BMRSfTrVs3AkA7d+5Uab+qfqlkZYl/4gYNRAY1pgOlfUvXrEnk5kbUti1R48YiUb5EItYV/ScbMEAss7ISM0e1by++XRs3Fst//rlw23v3xEMe1Wla/6zN+mtNb1uzsoiuXRNd4A8fFvXUP/0kGr+/+koM7pfbupVoxAiiQYNE+kBfXyIfH5FMo29fon/+Kdw2LIxo2DCikSOJxo8nmjJFDOafPVtkAUpNVbwU3ZrcpW44Ti64RjWRXW5Tgzp3+5rO96DNt0vdwEImEy/Z5s1E//0vUZs2hc+3sSHKyyvcdudOkfGUmxiqH+5DUExFkxs1b9683AiqLOp8EclkIiMa06LLl4mWLBHjvwAxO478W/qVV8rPNiO/jSooEHWy5W3bvHnhm/n222KZubkIGCwtiRo1Esu/+qpw1h2pVMzM06yZqCIqjTbnKy7r6rJyJVFkJNGGDaKtKiiIaOpUcTH39lbOqvP+++W/DhcvFm47d275254+Xbjt//5X/rZFpjSULv9OaZ3M1lYEdb17E732GlFCQuF+MzPpj1+eqFSzoEkNgbbfLm32Y7hzh2jLFjFJVFHyYKFGDdH8MHcu0aFDPOSxutC0X0KVTEy0efNmzJs3Tyl5RPPmzZGWlob8/HxcvXoVLVq0UGnflUncMnAgUKsWDyvUqgsXxIsbFgYkJSmvu3wZaN1a/H7gAHDtGvDTT8CtW8DOnUC9ekDt2uIhf2NiYgBfXzGjVJs2QE4OkJkJpKeL2XQSE4HNm4HoaDHmbeJEYMsWcbkqTWam2D8ADBsG7NoljmNnV/ioXVv8nD5dbBMfD9y4AZw4IcabmZuLn2lpwI8/imw/69YBNjZivz//DOzbJ8azZWYCN2+KcW+WlmJ9crLIiBMcLI5fnosXxTg3APjyS2DJEvE61a0rxszVqAHY2wNubiJDT+PGYtsTJ8QsRDVriuNKJMqPV14R+wHEa3jypMjUk5dX+MjOFuWfMwdo3lxs+8MPwMKFwIMHQG5uyfLGxopZjADg+++B995Dvm09XMpvgUv5LXAR7XAeHZDZuD3e/rYNRvxHOdOXVCpOxcOj9KGRAQHi5UtJKTkEUf5RiY8vffKo+HgxPFL+USmLJmUoa3+lDZ988QJ4801Rntu3lZ9jaQmMGSM+Sqxq0yR5UqWT+2kvhtEdeXQDgPz9/Sk+Pp6ysrIoPj5e0X/A1dVV7fSP5fWmlkpFbeqwYTy8sEKVrer++mvlF9vCQuSiBcSt0bNnJZ9TUaOwOnXIL16I6p4FC8S6FSuI5s0T1eBFG3v/85/y74x/+qnw2FOmlL9tWlrhfqdPL3/bS5fEdlIpkbu7aAoZNEhU17//vngd16wh+u035d5oRZtB9NlIXhqZjCgjQ9QC/fUX0fbtRKtWidzBcl9+Wf7rUPQ9P3dO1JTcvEnhYTK15nvQ1vTS2ujHIFeZt0omEwmV1q8XHwFnZ7Hdm28WbpOfL/o5TZ4sRkOdOyc+5qxqULdfQpVsMhg0aJDWxnYWV9b33erV4p/O2pooOVn9c6jyyvtGi40Vbdpyf/0lquj9/Ig2bRLV3Zp+S2vy7VzZ50ZEiIv0yZOivjYyUszEc+hQ4fPDwog+/ZRo1iyiDz8UU/QNHy7WjxypfOE+fJho+XLxDT9vnthm7VoxKUZqqvh2r0z5S2MMjeSqyMoSU4Xu3Cn6Orz5ppipqHZt5YlBPv64sAwNG9I9z4G0uk4QvYbt1BypBMgqnO9BWxdybQUW6r5VMpno4nH5cuGykydLxlM1axL16iViyMp+fFjVUiUDgoMHD5a6vjK9LNWRlSU6EALipoaVobRvtLg4ou7dC7+Vit7GFBSU7EKt6be0Jg3DmjYqa6NRWltXF22VR07bjeSq1jQUnzd86VJxC2xhUWptQvyWK4W7TE9XDqqKFEEbL482AgttvlVE4u354w+Rkrl9+8Ih0vJHkYy5dP266My4Zg3R338XdpthVU+VDAju3LlTYp1MJqvUOEx1fPON+Cdq2ZKr3cpU/ButoIBo27bC8Znyx3vvlfxyL28/RVX2m1GTOYM1nW9Y0+drs/5ZW/vS9tVKmzUNz56JKRCnTCEaOlSMJmnSRPkzNnKkuCL27St64R0+TPT8uaIomk4vrY2XR5tvu/y8ir/Ezs6ismrGDOUp2bdvLxlTtWghpo+eN49rRKuSKhkQlFZDsG7dOqpdu7bWawiys4kcHMQ/yYYNWtut8VP1Dq7oN9qhQ2JooPzbxcaGaOxY1RpSNf2W1mTOYE3nG9bk+dq8+BprI7mumzB+/115O3f3klc8a2sxPHXJEgoPk2n0dhc/LXU+stqsGFL1JT57luizz8R00vKRuUUfYWGF28bEiMEhX3wh4v3ERB7hYEqqZEBQWh+CBQsWEACytbXVqA9BccuXF37PlFLrWDWpcwdX9Btt9mzxu7090fz5om+AKt9oZZVB1W9pU81UqI2AiMj4GskN1YQhlYrUovJeeEWnNOzenaiggAoORlP03IN0ZNQqitt6Va1MhZp8ZI2pMufhw8JuLW+9pdztJzi4ZMAgkYh9DhmiPIqUE7YZnyoZEMhHGRQdh+nq6koAKDAwsNx9qPI9LZOJGkhA9NStFtS5g3vxQiyXf6NlZYnx8UU7zqla50lU5dL/qkRbAZGxNJJrcz+anpdMJmYYWrlS1KEXf50Bcav80UdEe/eqdAus7kfW2N6q0oSHF45okD/MzJT/Tkws3H75chF79esnBtL83/8RHTminA+L6VeVDAg2b95cIj+0jY0NAaDly5eX+Xx1bnxzcoh++EE5W1iVpc630tGjRB07ip7gzZtrr52ZaScgMpZGciLja8IoHvyeOSM+y/Lsl0WbF159VYx8qIgG75k23iptNj2UVbai9wl+fmJ/X38tvif/7ZpBRKJ2oXicJX80a6acF+vZM65R0IcqGRCUlqnQx8eHANDWrVtLfa4xzIhm9FT5os3OLszyBxDVqyeGYGijqptplzZqG7RxtTKmJozygpxHj0TVoK2t6KAof83WrCnc5u5d0Veh6O2uFjpLavpW6aKGQN14MDNTdF7ctInok09E8NC8eWEzw+7dhXHTtGliZGn//kSff060a5dyigqmHVU2ICiuQ4cOBIAOHDhQYp06H+irV6vhzWxlv2iXLCHq0KHwG+vNNwuzNWnj4sO0T1u1DZq8t8ZUL17ZfRw+LLrZL1tGdONG4frvvy+sM+/WTXTJB4gGD9b4jqOggCj6YAGFzj1P0XMPUsHB6Eq/X9oeEEKk/SDjl1+Uu3HIK2FKq0lo0ULMG8MdF7Wj2gQEDRs2JAB09uzZEutU/UDn5YkqrbZtC5PEVQsVvVBHj4r1Vlbip6OjGFFQXHVu+6/qNH1vjaUJQ9Naho0bxQD/4lcwc3PRfFY0u6Seh2Vqq0+qnD5GQAwdKvbx9ttEb7xB1K5d4ak7OyuPIl29WtQ63Lql2nmwahIQyGQyGjx4MHXu3JnuF81m9i9VP9Dr1xde70w6MlX1y7uiL9qhQ0VmQUC0qZbyWjNWIWNowtDWbe+tW0Rz5hReueS1BkUTbsknmZo5U/S6Ky+ZiZbaNkt9iZ1yKHy76sG5ofqCPnlCtH+/yMhddDt5kjhATPY0bZrIFfPLL3z/UZFqERBURJUPdH6++FACJWcaMynq3mVU9EX73XeiyYB7ADFNGLoJQ1f5Hq5fL9n5sGdP5TLa2BC9/LIIFNatK7z91XJ9f8Hv4RTdaDSFYgxFoy8VwEytBFAFeQXk4phD/i/dIumhaKXj67ulJztbxF/du5cc4aBmt41qhQMCUu3/bONG8aFycDDhOcc1vcvgfgDMVJhCvoelS8X6zp1Fz7mi/1fOziX3M3GiyIJ25IiYdEv+pWWIBFD/fheEYwRJICV/7KQ4xxGUuWWnwUdA/PJL4YiF+vXF75MnF5YpNFT0c+Zp6gtVi4BAVl4q3H9V5n//xQsxVbu835xJ0tZdRn4+0ahRhVWgRQcYM1ZV6DvfQ0GB6KS4caNoQggKKtxWfsUs/rC2Fv0VJk5UvmKePClmNCp+hdXWd0CxoCJ8Sw65OD1X+6VSiptKCeRUiXeKn6JUSnTsmJgLTH6K8gyzgJjUadkyUYFT0X6rcvenahEQhIaGUp06dWjChAnlPr+i//0tW8Qye/uyo1ijp416uZwckZ9U/iIFB5c//wBjpsxY8j3IZ8scMYLolVdEF3tz88L/Q09P5f/dolmCatcWdzPduhF17VryO2DLFtET77ffxN0OIBIznTpFlJKiXI47d8SIiiZNRF+he/fEKKJHj6jgcQZFH5JS6EvfULTjGCp4miW+L3JzxU3EixdlNicqLuJd75K0uavSF7G0uSv5d71b6RaRygYXpWWt7tqVaNGiksMa9TGZp6FVi4BgxYoVBIBGjx5d4T7K+98fP158CBYs0G659UrTern0dBFOA6IDYRl5HRhjxehiWGZ+vrhg79kj/i/lV0yZTNQa2NqWXqtQ/Dug6O1y8UeXLsrlKC1zo/zRrp3YRn7FlScWKP6QSMRscEUNHkzZtR3oH7hRqkVLymregV5070VP2vWiI1avkARSCv/43yBm0yZxI7JypWg++f13on37RDXAxYuFX3NbIku9imdu2an4mrt9WzQd9O2r3O/gn38Ki7Zpk1imaGF5WkDxa06T/0u3SCKRqdUZs+jbaiy1DtUiIJgzZw4BoPfff1+j/ctk4nNu0gkxNKkh+OefwjaTevVEGyZjrPIMMSwzM1MMc/zrL9Gh8dNPS34HjBtHNHCg6MwoH9Pn7CwegwYp78/NTblmouijffvCY8r3UVbw0KKF8n47dy5z23TzRhTeNbgw4OnTp+z9WlsXfs2hp6jJaNZM3Pr36kXUuDHdhSPNxzz6Z9oKpffg3sXHtPa7XHpjcmGNZ0EBUa1a4v5n8mSiHZ/EUUYzUbUghYT8sZNcLW5Qwe+qVxWEbxcdMtWuddByNFEtAoI333yTANACk7611xJN2g/XrCn8xBbNK8oY0x9jSABV9MZCJiuc0lw+ZFJ+Y7F3r6iJyMgQ4wQfPRJ3VGlp4lGUvBfg2rVUsO8gnQ2Oor8+CqPzc7aQdOs25ZuVb74RCc9GjxbDnb29RVrp5s2J2rQRIx8sbpJ/o2Mk61128CCzslJu7hw2rLBfRtOmRF26UHr3IVQfD5WeaoYC6tLyCc188wktm35VFAs+qmWd/Dhe0REzHj0oCzYU7xhA/l3vVq4VqZw2jArjhDI2qBYBgb+/PwGgH374Qa395ueb8IiC0mjSnrlihWgzZIwZjqETQOki5aE2hxlER4uRDxIZ/df3H0pac4RytkVSypwN9GP7ZfQ/zKGrcBUBQFF9+5YaOOTDgg7VeJU+qv0judncUVrdF9GUhkb0uKYTUc2atLXP/9HJSaso+2mR6W8fPFC6iBT8Hk4uuEb+jY6R9KjyKA8pzCruL1HOKJFwBJKLw7Oyax3KCSSqRUDQo0cPAkA7duxQa7/btxPZ2RHNm6eFQhqLyt5lREebcA9KxliZNK1p0HbKQ23mQP43uAjfklP6KW7JKT24kMnE9921a2KihagouvjJBvoEi+jm0LcV5bt1i2hrx0X0X/xAa/CuYucPUb9INwkZtWol+oB+0Xob/Y6R9I9tZ6J27eiRZSPahIl0e8Knoh+E/MovlRL5+VGcU2DZp1pOMBa+XSpqHWodoPi/C0qOJv04vtzhpps/OFj1AwL51MdHjx5Va7+DBok3eM4cbZRSyzS5U6jouZGRRBYWor3u2TPtlZkxZhy0UdOgrZwk2qx1KBJclHqKKgQXimK9dIukkCjdIEnzXtDrr6TREOcEKtgeTlfQgnwdkqlBLeU7dPnjI3xb2C8CDcgHh2mi+RaaEySjZctE/8jI9kEUgz60DwPoTqdBRFOmiIvPd9+JUSDy0SbFAidFOfs8IikkVHAwmh4/FkMtz5wh8uolIwezdDr98oeK1/fWLTGpVPfuMnKwfELA06ofEAQGBtJLL71EqampKu/zxo3C2U6vXNFSQbVFl+Ngdu4kqlFD7HPMmPLTqTLGqi9tdmzTVq2Dlps0RLFkIvHSuqTSi1UkyJDJRMvqgQNEy5cTTXlDRj27vaD1X98lmjOHLqAtvY5NZfaLBIg+wSLFH5fRiszxgmrhGdWzySUn3CGX5lJqYZNGzcxv0cJG39FjTx/ajtfoeL9PyBo5Ze53dPerIkq4fZvunk0vtr4KNhlERUUppj0u0LDX5VdfiRfK11dLhdQWXc7X/McfHAwwxgxDW7UOWm7SCN8uOiqWWiwVczQXwIyaNnpOnTuL0ZPTp4th7YMHE/Vo84Tq4jF1r3mWpOt/IlqwgM6OCS43eJiJpYo/bqCp0jpra6JGjYhaW6ZSNxynbzBLsfIFzGkd3qIIsxEEEP0M/6oXEBR9uLi4UHixN76yAW1BgRitAhjZcHtddOiR27WrMBgYPZqDAcaY/mmr1kHLadYLfg+naPiIxEurk6ngifqzcYZ3DRa1DkXjlb+l5N/omMi5UCS3QX6+yJdw9SrRhaQCSnAaTMf6zKK4Hffo5OaLdGtLDJ2fG0rTsIquNfWm67Yd6P6Uzyg3t/Cw6V0H0Q00pfyadkQ1a4rmYHlgUMOaAKI7qFX1AoKDBw9SVlYWxcfHk5+fH0kkEkVQoEot+759Yn3dukY2q6G2JyCX27OncLbCUaM4GGCMmT5tZ/7R4myc4V2DS6Z7xtXCBEwVPL9oNFEQG0cuNe+RP3aSdLtyWaRSIn8/mciX4De88EZSJiN68YKO7BF9Hg428qt6AUHRk4mMjCQLCwuqVasWbd9eoFItuzxV//Tpej6Rimh7FhC5xESRl3nkSBGSMsYYK0mLs3EWwIyi0VfMOuk4pvIJjkoJTMId3ilZ6xBXxiiDIhsU+A0nF1yjQZ0vV+2A4Mcff1Q0Hzg65qhUy37lCtHnnxOdO6enE6gsXdUQEImT5mCAMcZ0T9PAopTnV1iBUcYG4R/HV7pToYSICEYuMzMTderUQUZGBuzs7AAAISEhmDNnDgBHAGmIjwd69iz53Ph4wMsLiI4GfHz0WWo1SKWAmxvg4QFERgJmZoXrZDIgIABITgZSUgBz8/L39c8/wNOnQPfuOiwwY4wxfZFKgdhYIC0NcHICvL2LXQrK2GDLlky8/rryNbQ0Fro/Bd24f//+v7/VAwC4u5e+nXx5Wpruy6Qxc3Ng2TJg5Ehx8Q8KEieQnAyEhABRUUBYWMXBQFoaMHAgkJ4O7NkD9O2rl+IzxhjTHXPzCm5sy9hg2LDK7d+s4k2MkzwgqF3bCoC4ZpZGvtzJCThxQlxn9+7VQwHVFRgoLvpJSaJqw85O/ExOFssDA8t/fkYGMHgwcP064OwMtGunl2IzxhgzbSYVEJw4cQJZWVmIj4/HwYMHAQCTJr0MFxcgOFjUqhclk4kba1dXUXOyfj2wcyewbZv+y66SwEDgyhXRzhEaKn6mpFQcDOTmiojn7FmgUSNg3z7AwUEvRWaMMWbaTKoPQVE1atTAixcvsG/fPjx79ipGjgT8/MquZX/lFVFLkJ0N/PWXCBCqFKkUGDNGnGzt2sCRI0DnzoYuFWOMMQMrrR9eaUyqD0FUVBQyMzPh5OSE9evX4/Lly2jSpAnatxfXwVmzRO26nKtrYS37pk0iGGjdGujTx2CnoBtEwIcfipOtUUN0SORggDHGmApMKiDw9vZWRDc+xTpOBAYCw4eX3QNzxw7xc/x4QCLRY6H14cUL4MYNcWKbNwP9+hm6RIwxxkyMSQUEFSmrB2ZODnDggPh9+HC9FqkS40S0wNJSRDxHjgD9+2t334wxxqoFk+pUKKdqt4eDB4Hnz4HmzYGOHXVUqNJERIi8Ar6+wLhx4qebm1iuDVeviuYCALCw4GCAMcaY2kwyIIiJiYGdnR0GDRpUqe3NzUV+noAAPTYXRESIfAIeHiI7UlaW+OnhIZZrGhScPSv6CUydKpoMGGOMMQ2YZJNBeno6srKykJubW6nthw4Vj+LDEnVGKhU9HP38lDMO9uwp/g4IAGbPFu0X6jQf3LwJDBkigowrV/R4Yowxxqoqk6whSE9PBwA4qDjG3kxfZxsbKxIDzZlT8qBmZmJsZGqq2E5VT56IxEN37wIdOoi+A1ZWWik2Y4yx6sskAwJ5lsLKBATnzonkfXolz5Os7XzKeXmiduHCBZGFcM8eoF49tYvJGGOMyZlkQKBKDcGYMUDDhsDhw7ouVRFOTuJnZfIpV5ZMBkyaJLIq1a4tgoFmzTQrJ2OMMfavKh0QpKQAFy+KjvgvvaSPkv3L2xuVzqdcWWfOAOHhYjTBjh1Ap05aLTJjjLHqrUoHBH/8IX76+AB16+q2TErksxZGRYkq/qKjDAICxPKlS1XrUNi1q+iQ+MsvPLyQMcaY1qkVEKxZswYuLi6wtrZGjx49cOLEiTK33bRpEyQSidLD2tpa7QIDQIcOHdC1a1c8eNAM27YBMTGiY39x8oCgslM/apWmsxbKPXtW+PvQocDYsbopL2OMsWpN5WGHv/32G2bOnIm1a9eiR48eWLFiBQYOHIjLly+XecduZ2eHy5cvK/6WqJkMYPv27WjZsiUGDlyLgwfN8c47hetcXMRNufw6+/Ah8Pff4neDBARAxfmUK7JmjahJOHgQaNlSt2VljDFWralcQ/Dtt99i6tSpeOONN9C+fXusXbsWtWrVwoYNG8p8jkQigaOjo+LRqFEjtQr71ltvwdd3JUaNkqBBg7Ry8/3s2SOa6zt1EhkKDUaeT3nsWPGzssHAzz8D06eL4YvaymzIGGOMlUGlgCA/Px+nT5/GgAEDCndgZoYBAwYgPj6+zOc9e/YMzZs3R9OmTTF8+HCcP3++3OPk5eUhMzNT6QEAN2/egaPjVjRqdBKnTjXB3bsRsLUtzPfj5yfy/UilwM6dYl86rR2QSkV7RXntFuoIDwemTBG/f/SROCnGGGNMh1QKCB4+fAipVFriDr9Ro0a4d+9eqc9p06YNNmzYgJ07d2LLli2QyWTw8vLC7du3yzxOSEgI6tSpo3g0bdoUAHDunC3u3UtCZuYA1K9fF7Nnz4b034tw8Xw/K1YAq1eLKQR0QlfzFOzbJ2oTZDLgzTeBb7+tgtMzMsYYMzY6H2XQq1cvTJw4EZ6enujbty8iIiLQsGFD/PDDD2U+JygoCBkZGYrHrVu3AAAi5riN58+fwcnJCampqYgtku2vaL6fpk2BadOAtm11cFK6mqfgwAFgxAgxN8F//gP88AMHA4wxxvRCpYCgQYMGMDc3V2QKlLt//z4cHR0rtY8aNWqgc+fOuHLlSpnbWFlZwc7OTukBAOIQIjho3bo1ACCtSLY/dfL9qKz4PAU9e6LMdgtVEAFffCGmZfTzAzZv1v40yYwxxlgZVAoILC0t0aVLFxw6dEixTCaT4dChQ+jVq1el9iGVSpGUlAQnNa7aYvTebUVZACj2I8/34+ICrF8vbq5zclQ+RMV0NU+BRCLGSc6eLfoQ/Ht+jDHGmD6o3GQwc+ZMrF+/Hj///DMuXryId999F9nZ2XjjjTcAABMnTkRQUJBi+6+//hr79+/HtWvXcObMGUyYMAE3btzAW2+9pXJhzc0BDw8REOze3QhOToHw9PRWyvfz7rtAaKi4rurkBlvb8xRculT4e8OGwDffcDDAGGNM71TOQzB69Gg8ePAAX3zxBe7duwdPT0/s3btX0dHw5s2bMCty5/zkyRNMnToV9+7dQ7169dClSxfExcWhffv2Khc2KysLz55dBAA8e9Yaz559p5jbx9VV5PtJSBB/Dxqko0kAi85T0LNnyfWqtFv88APw3nvA998D//2v9srIGGOMqUhCRGToQlQkMzMTderUUVr2v/8Fo3fvoBL5fjp1EjMc/vIL8PrrOiiMVCpGE3h4iD4DRZsNZDJRVZGcLCZSKK+KYsUKYMYM8fsHHwDffaeDwjLGGKvu5NfQjIwMRZ+80qhcQ2BIP/74IzZv3oxnz55h8uSJaNxYeX1qqggGzM1Fll+dkM9TMHKkuPgHBYlmguRk0YkhKkpUVZQXDISEiD4IAPDJJ8CiRToqLGOMMVY5JhUQjBo1Cm+++WaZ6w8eFD+9vID69XVYEPk8BbNmiYPJydstypqngEgEEIsXi7+/+kqMLOChhYwxxgzMpAKCisjnLujbVw8HU3WeAiLR41Gef+GbbzgDIWOMMaNhUgGBVCoFEZU5OVJuLmBhAfTpo6cCyecpqAyJBHB2Fj/XrQPUGGXBGGOM6YrOMxVq0/r162Fra4v333+/1PW//QZkZIgswkZp3jzgzBkOBhhjjBkdkwoI7ty5g5ycHFhYlF2xUauWlobxa2PiosxMMTnRs2fib4kE8PTUQuEYY4wx7TKpJoO7d+8CAJo0aVJinVSqxUREERGiw+D164XLXFzE6IKyOgwW9/ChSIZw+jRw+7bobMgYY4wZKZOrIQBKDwh69QK6dBE18hrRxsRF6emib8Hp00CDBoVDDBljjDEjZVKJiZo1a4abN29i1aqjsLf3UnTsz8kB6tYVeYFu3QJKiRcqRxtJhx4+BPr1A5KSRCfCQ4d0NOUiY4wxVrEqmZjozh3RZPD++4VXfBcXYOJEcb1u3lyDYAAonLho27ayJy7y8hLblTa64NEjYMAAEQw4OQHR0cC/szIyxhhjxsykmgyk0gIAEsTGOinV5H/9tViv8XBDTScuGjMGOHsWaNQIOHyYgwHGGGMmw6QCgoYNh+DVV19Fnz41YGsr5haKjBTN9IBy0kC1FJ24qDQVTVy0ZAnQrp0IBriZgDHGmAkxqYBg27Zt2Ldvr9IyqRTIzha/16yp4QG8vUUbRHCwaIMoSiYTcxC4uortStO5s2guUGMmR8YYY8yQTCogaNeu5LKzZ4Hnz8XvGucfkE9cFBUlOhAWHWUQECCWL11a2KEwK0sMLYyLU94HY4wxZmJMKiAIDU2AtFiCoBo1gFdeEb8Xn/1QLfKJi5KSRBuEnZ34mZysPHFRfr6Yy2DfPmDsWPE3Y4wxZqJMatghUAP169fB+vU/IPDfC3NlRwOqTCote+IiIuCNN4CffwZq1xZ9Brp21dKBGWOMMe2pksMOgReQSj/Ea68txZYtFmjRYhhCQkRNfliYlmvry5u4aOFCEQyYmwPbt3MwwBhjzOSZWEAA1Kzpg4yMuZgwQfzduDHw+++VzyissdBQMUkRAKxZAwwcqKcDM8YYY7pjUn0IAODo0aZYvfo8gLEYO/YG7twBfvxRTwc/eVI0FQBiroO339bTgRljjDHdMrmAoGlTZ0yc2AzAr7h37wUAoFMnPR3c3R0YNgwYMULkHGCMMcaqCJNqMmjUqBFq1KiBU6dOAQCuXhUJgjTOUFhZNWsCv/0mRhQUT23MGGOMmTCTuqo1btwYMpkMISEhaNbsJdy8aQNACxkKy5OXJ9ok5IMxzMwAa2sdHpAxxhjTP5MKCNzc3BAQEICoqCiMG7cagEhWZG+vw4O+9x4wdSrw7rs6PAhjjDFmWCbVZPD777/D1dUVYWFhOHasFwAVmwvKyy1Qmo0bgQ0bRK2A3oYxMMYYY/pnUgFBVFQUBg0aBHNzcyxdKpZVOiCIiBAjA65fL1zm4iJSFZd2sT93TtQOAGI6xVdf1aDkjDHGmHEzqSYDb29vmP97R//VV8CcOWXnDlISEQGMHCnmSi46P4GHh1geEaG8fWamWJ6bCwweDAQFaftUGGOMMaNiUqmL58w5hFde6VdhTb8SqRRwcxMX/8hI5dEBpeU9JgJGjxYZCJs2BRISdNxJgTHGGNOdyqYuNqkaguDgBvD1Fdf34jf1ZYqNFc0Ec+aUHCpoZibu/lNTxXaAaCrYsUPMmvT77xwMMMYYqxZMKiBITXVGfDxQty7w2mvA1q2VeFJamvjp7l76evly+XadOong4IcfgJ49NS0yY4wxZhJMKiCoX98S3bqJG3oA+OQT0SJQLieRvAjJyaWvly+XbweIQECeopgxxhirBkwqIADE9TsjA6hVC7h7t7Cmv0ze3mI0QXCw6DNQlEwGhISI9du3i+YCxhhjrBoyuYAgOlr8lGcnlNf0l8ncXAwtjIoSHQiLjjIICBDL+/QB/u//AF9fsY4xxhirZkw2IGjdWvwsWtNfpsBAICwMSEoSkYSdnfiZnCxqDn79VWy3eDFQu7ZOys0YY4wZM5Madvj4cQZcXe2QkSFaAm7fLhwtWCnFMxV26gR06SI6JYweDWzbBkgkOj0XxhhjTJ8qO+zQpDIVxseL/gMWFuK6Hh6uQjAAiI2LZjKaMEEEA82bA2vXcjDAGGOs2jKpJoOhQ8VPS0sRDGg0vcDmzWLcork5EBoqxjIyxhhj1ZRJ1RCsXy/yBXl4AC+9pMGOiESfAkDkQNbp/MmMMcaY8TOpPgQVtX+oRCoFtmwRzQYqtTswxhhjpqNK9iGQVpiFSAXm5sCkSdrbH2OMMWbCTKoPQePG+/HllxVlIirHwYPA9OnA8+faKxRjjDFWBZhUDcHz54Px9defoVOnBwgs3qOw+JDC4lMipqcDr78O3LsHODgAX3yh38IzxhhjRsykaggAwMeHMHv2bOXmg4gIMQWiry8wbhxKTIkolQLjx4tgoH17YPZswxSeMcYYM1ImFRA0agQsXDgRqampiJVPYhARAYwcKYYeFE1L7OEhlkdEAPPni+aCWrWA334TPxljjDGmYFJNBi+/DHh4iOmK09LSxJ3/rFmAnx8QGQmY/Rvf9Owp/g4IAN57D7h/Xyxft67saZAZY4yxasykaghefhlI/ne6YicnJ9Fn4Pp1YM6cwmBAzswMmDy5MBiYNk00GzDGGGOsBJMKCHr3liEkJASurq7w9vYunOqwrLv+mjXFz5YtxYyHjDHGGCuVSTUZBAWNxf79+xAWFgZzc/PCqQ6Tk0UzQXHydMTBwYCVld7KyRhjjJkak6ohuHTpIsLCwgqHHHp7Ay4u4oIvkxVumJcn/g4JAVxdgddeM0h5GWOMMVNhUjUECQkJqFevXuECc3PRFDBypOhAGBQklg8fLgKBkyfFnAWcmpgxxhgrl0kFBOalXdgDA8VFf9Ys5UmKnj0Dtm/XcEpExhhjrHowqSaDMgUGAkeOAM2bi78bNRKjD7ipgDHGGKuUqhEQpKSIMYk3bgCNGwOHD4v0xIwxxhirFNMKCLZvB2JiREIiuYQEoE8fEQy0agUcPSrSEzPGGGOs0kwrIHjrrZLzFGzcKCYu6twZ+PvvwmYDxhhjjFWaWgHBmjVr4OLiAmtra/To0QMnTpwod/vt27ejbdu2sLa2hoeHB/bs2aNWYXHnTsl5Cr79FvjqKyA6mpsJGGOMMTWpHBD89ttvmDlzJr788kucOXMGnTp1wsCBA5Genl7q9nFxcRg7dizefPNNJCQkICAgAAEBAYoUxCqxtRVZBydMAIYMEbMWSiTAl18Cdeqovj/GGGOMAQAkRESqPKFHjx7o1q0bVq9eDQCQyWRo2rQp3n//fXz22Wclth89ejSys7MRFRWlWNazZ094enpi7dq1lTpmZmYm6tSpg4wGDWD38KFY6OcHREWJmgEfH1VOgTHGGKs2FNfQjAzY2dmVuZ1KNQT5+fk4ffo0BgwYULgDMzMMGDAA8fHxpT4nPj5eaXsAGDhwYJnbA0BeXh4yMzOVHgCAhw9FjUCLFoXzF8jnM2CMMcaY2lQKCB4+fAipVIpGjRopLW/UqBHu3btX6nPu3bun0vYAEBISgjp16igeTZs2FSuio4GsLODqVWDYMLFMPp8BY4wxxtRmlKMMgoKCkJGRoXjcunVLrHjpJcDGRnmeAm9vwxaWMcYYqwJUSl3coEEDmJub4/79+0rL79+/D0dHx1Kf4+joqNL2AGBlZQWr0mYnzMoCzp8XwUBUFM9TwBhjjGmJSjUElpaW6NKlCw4dOqRYJpPJcOjQIfTq1avU5/Tq1UtpewA4cOBAmduXq0kTMV9BcrIIBnieAsYYY0wrVJ7caObMmZg0aRK6du2K7t27Y8WKFcjOzsYbb7wBAJg4cSIaN26MkJAQAMCHH36Ivn37YtmyZRg6dCh+/fVXnDp1CuvWrVO9tD/+KIYdentzzQBjjDGmRSoHBKNHj8aDBw/wxRdf4N69e/D09MTevXsVHQdv3rwJM7PCigcvLy+EhoZi7ty5mDNnDlq1aoXIyEi4y0cJqGLUKKCcIROMMcYYU4/KeQgMobJjKBljjDGmTCd5CBhjjDFWNXFAwBhjjDEOCBhjjDHGAQFjjDHGwAEBY4wxxsABAWOMMcbAAQFjjDHGwAEBY4wxxsABAWOMMcagRupiQ5AnU8zMzDRwSRhjjDHTIr92VpSY2CQCgqysLABA06ZNDVwSxhhjzDRlZWWhTp06Za43ibkMZDIZ7t69i9q1a0MikRi6OBrLzMxE06ZNcevWrSo3NwOfm2niczNNfG6mSd/nRkTIysqCs7Oz0uSDxZlEDYGZmRmaNGli6GJonZ2dXZX7oMvxuZkmPjfTxOdmmvR5buXVDMhxp0LGGGOMcUDAGGOMMQ4IDMLKygpffvklrKysDF0UreNzM018bqaJz800Geu5mUSnQsYYY4zpFtcQMMYYY4wDAsYYY4xxQMAYY4wxcEDAGGOMMXBAwBhjjDFwQKAzISEh6NatG2rXrg0HBwcEBATg8uXLStvk5uZi2rRpsLe3h62tLV577TXcv3/fQCVW36JFiyCRSPDRRx8plpnyud25cwcTJkyAvb09atasCQ8PD5w6dUqxnojwxRdfwMnJCTVr1sSAAQOQkpJiwBJXjlQqxbx58+Dq6oqaNWuiZcuWWLBggdKEJ6Zybn/99Rf8/f3h7OwMiUSCyMhIpfWVOY/Hjx9j/PjxsLOzQ926dfHmm2/i2bNnejyLspV3fi9evMCnn34KDw8P2NjYwNnZGRMnTsTdu3eV9mGs51fRe1fUO++8A4lEghUrVigtN+Vzu3jxIoYNG4Y6derAxsYG3bp1w82bNxXrDfndyQGBjhw5cgTTpk3DsWPHcODAAbx48QKvvvoqsrOzFdvMmDEDu3btwvbt23HkyBHcvXsXgYGBBiy16k6ePIkffvgBHTt2VFpuquf25MkT9O7dGzVq1MCff/6JCxcuYNmyZahXr55imyVLlmDlypVYu3Ytjh8/DhsbGwwcOBC5ubkGLHnFFi9ejO+//x6rV6/GxYsXsXjxYixZsgSrVq1SbGMq55adnY1OnTphzZo1pa6vzHmMHz8e58+fx4EDBxAVFYW//voL//3vf/V1CuUq7/xycnJw5swZzJs3D2fOnEFERAQuX76MYcOGKW1nrOdX0Xsnt2PHDhw7dgzOzs4l1pnquV29ehV9+vRB27ZtERMTg3PnzmHevHmwtrZWbGPQ705iepGenk4A6MiRI0RE9PTpU6pRowZt375dsc3FixcJAMXHxxuqmCrJysqiVq1a0YEDB6hv37704YcfEpFpn9unn35Kffr0KXO9TCYjR0dH+uabbxTLnj59SlZWVrRt2zZ9FFFtQ4cOpSlTpigtCwwMpPHjxxOR6Z4bANqxY4fi78qcx4ULFwgAnTx5UrHNn3/+SRKJhO7cuaO3sldG8fMrzYkTJwgA3bhxg4hM5/zKOrfbt29T48aNKTk5mZo3b07Lly9XrDPlcxs9ejRNmDChzOcY+ruTawj0JCMjAwBQv359AMDp06fx4sULDBgwQLFN27Zt0axZM8THxxukjKqaNm0ahg4dqnQOgGmf2x9//IGuXbti1KhRcHBwQOfOnbF+/XrF+tTUVNy7d0/p3OrUqYMePXoY/bl5eXnh0KFD+OeffwAAZ8+exd9//43BgwcDMO1zK6oy5xEfH4+6deuia9euim0GDBgAMzMzHD9+XO9l1lRGRgYkEgnq1q0LwLTPTyaT4fXXX8fHH3+MDh06lFhvqucmk8mwe/dutG7dGgMHDoSDgwN69Oih1Kxg6O9ODgj0QCaT4aOPPkLv3r3h7u4OALh37x4sLS0V/8ByjRo1wr179wxQStX8+uuvOHPmDEJCQkqsM+Vzu3btGr7//nu0atUK+/btw7vvvosPPvgAP//8MwAoyt+oUSOl55nCuX322WcYM2YM2rZtixo1aqBz58746KOPMH78eACmfW5FVeY87t27BwcHB6X1FhYWqF+/vkmdKyDanD/99FOMHTtWMXOeKZ/f4sWLYWFhgQ8++KDU9aZ6bunp6Xj27BkWLVqEQYMGYf/+/RgxYgQCAwNx5MgRAIb/7jSJ6Y9N3bRp05CcnIy///7b0EXRilu3buHDDz/EgQMHlNq+qgKZTIauXbsiODgYANC5c2ckJydj7dq1mDRpkoFLp5nff/8dW7duRWhoKDp06IDExER89NFHcHZ2Nvlzq65evHiB//znPyAifP/994YujsZOnz6N7777DmfOnIFEIjF0cbRKJpMBAIYPH44ZM2YAADw9PREXF4e1a9eib9++hiweAK4h0Lnp06cjKioK0dHRaNKkiWK5o6Mj8vPz8fTpU6Xt79+/D0dHRz2XUjWnT59Geno6XnrpJVhYWMDCwgJHjhzBypUrYWFhgUaNGpnsuTk5OaF9+/ZKy9q1a6foBSwvf/Fev6Zwbh9//LGilsDDwwOvv/46ZsyYoajlMeVzK6oy5+Ho6Ij09HSl9QUFBXj8+LHJnKs8GLhx4wYOHDigqB0ATPf8YmNjkZ6ejmbNmim+W27cuIFZs2bBxcUFgOmeW4MGDWBhYVHh94shvzs5INARIsL06dOxY8cOHD58GK6urkrru3Tpgho1auDQoUOKZZcvX8bNmzfRq1cvfRdXJf3790dSUhISExMVj65du2L8+PGK30313Hr37l1ieOg///yD5s2bAwBcXV3h6OiodG6ZmZk4fvy40Z9bTk4OzMyU/+XNzc0Vdy6mfG5FVeY8evXqhadPn+L06dOKbQ4fPgyZTIYePXrovcyqkgcDKSkpOHjwIOzt7ZXWm+r5vf766zh37pzSd4uzszM+/vhj7Nu3D4DpnpulpSW6detW7veLwa8LOu+2WE29++67VKdOHYqJiaG0tDTFIycnR7HNO++8Q82aNaPDhw/TqVOnqFevXtSrVy8Dllp9RUcZEJnuuZ04cYIsLCxo4cKFlJKSQlu3bqVatWrRli1bFNssWrSI6tatSzt37qRz587R8OHDydXVlZ4/f27Aklds0qRJ1LhxY4qKiqLU1FSKiIigBg0a0CeffKLYxlTOLSsrixISEighIYEA0LfffksJCQmKXvaVOY9BgwZR586d6fjx4/T3339Tq1ataOzYsYY6JSXlnV9+fj4NGzaMmjRpQomJiUrfL3l5eYp9GOv5VfTeFVd8lAGR6Z5bREQE1ahRg9atW0cpKSm0atUqMjc3p9jYWMU+DPndyQGBjgAo9bFx40bFNs+fP6f33nuP6tWrR7Vq1aIRI0ZQWlqa4QqtgeIBgSmf265du8jd3Z2srKyobdu2tG7dOqX1MpmM5s2bR40aNSIrKyvq378/Xb582UClrbzMzEz68MMPqVmzZmRtbU0tWrSgzz//XOkiYirnFh0dXer/16RJk4iocufx6NEjGjt2LNna2pKdnR298cYblJWVZYCzKam880tNTS3z+yU6OlqxD2M9v4reu+JKCwhM+dx++ukncnNzI2tra+rUqRNFRkYq7cOQ350SoiJpyhhjjDFWLXEfAsYYY4xxQMAYY4wxDggYY4wxBg4IGGOMMQYOCBhjjDEGDggYY4wxBg4IGGOMMQYOCBhjjDEGDggYY4wxBg4IGGM6smnTJmzatMnQxWCMVRIHBIwxxhjjgIAxxhhjAE9uxBjTmvz8fHTv3h0A8PjxYwBA/fr1AQAnTpyApaWlwcrGGCsfBwSMsQo1adIEc+bMwf+3a4euicYBGMefya0s2MxGV2SgQaMTYf/GxLa0xcXlsWQcS8uGlSWTgkWwLiwuLvkP6KUbHHdD4XTv3fH5tPeFF57ywhd+v4uLi493s9ksvV4vLy8vqVarv3zz4/7A+fn5F60E/oQjA2CjVquV+Xz+8bxer3N5eZmrq6vfxgDw7xEEwEbtdvunIHh8fMzb21uur68LXAXskiMDYKPpdJpOp5PlcpmDg4PUarXc3NxkMBgUPQ3YkW9FDwD+fs1mM6VSKYvFIuPxOJVKJf1+v+hZwA4JAmCjo6Oj1Ov1jEaj3N/f5/n5OaWSE0f4n/ijga202+0Mh8OcnZ2l0+kUPQfYMUEAbOXk5CSHh4e5vb0tegqwBy4VAls5PT1No9HI3d1d0VOAPXCHAPjUarXK+/t7Hh4e8vr6mqenp6InAXsiCIBPTSaTdLvdHB8fZzQapVwuFz0J2BNHBgCAS4UAgCAAACIIAIAIAgAgggAAiCAAACIIAIAIAgAgggAAiCAAACIIAIAIAgAgyXevD0soa8s0OwAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "import matplotlib.pyplot as plt\n", "\n", "fig, ax = plt.subplots()\n", "fig.set_size_inches(6, 5)\n", "\n", "get_prof_plot = lambda arr: arr[wall_pos:mid_pos] / u_fric\n", "\n", "ax.plot(df_cp[\"ux_rms\"][\"y+\"], df_cp[\"ux_rms\"][\"u/u*\"], \"ko\", fillstyle=\"none\")\n", "ax.plot(df_cp[\"ur_rms\"][\"y+\"], df_cp[\"ur_rms\"][\"u/u*\"], \"ro\", fillstyle=\"none\")\n", "ax.plot(df_cp[\"uang_rms\"][\"y+\"], df_cp[\"uang_rms\"][\"u/u*\"], \"bo\", fillstyle=\"none\")\n", "ax.plot(x_vals, get_prof_plot(ux_rms), \"k--\")\n", "ax.plot(x_vals, get_prof_plot(uy_rms), \"r--\")\n", "ax.plot(x_vals, get_prof_plot(uz_rms), \"b--\")\n", "ax.legend(\n", " [\n", " r\"Exp. $u_x$\",\n", " r\"Exp. $u_r$\",\n", " r\"Exp. $u_{\\theta}$\",\n", " r\"Nassu $u_x$\",\n", " r\"Nassu $u_r$\",\n", " r\"Nassu $u_{\\theta}$\",\n", " ]\n", ")\n", "\n", "ax.set_xlim((1, 170))\n", "\n", "ax.set_xlabel(\"$y^+$\")\n", "\n", "plt.show(fig)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Also, an excellement agreeement is obtained for all directions in cylindrical coordinates. In general the results confirm the solver capability of solving a pipe flow turbulence, confirming the IBM as adequate to represent a curved boundary and the D3Q27 velocity set as capable of axissymetric turbulence." ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "## Version" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Version: 1.6.33\n", "Commit hash: fbc0edb5260d2734f0a290e1806c26ac6d865ff4\n" ] } ], "source": [ "sim_info = sim_cfg.output.read_info()\n", "\n", "nassu_commit = sim_info[\"commit\"]\n", "nassu_version = sim_info[\"version\"]\n", "print(\"Version:\", nassu_version)\n", "print(\"Commit hash:\", nassu_commit)" ] }, { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "## Configuration" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
simulations:\n",
       "  - name: periodicTurbulentPipe\n",
       "    save_path: ./tests/validation/results/05_turbulent_pipe_flow/periodic\n",
       "\n",
       "    n_steps: 400000\n",
       "    # u* = 0.0034013 / R = 72 / ETT = R/u* = 18,820\n",
       "    # Re_tau = 180\n",
       "    # y+ = 2.5\n",
       "\n",
       "    report:\n",
       "      frequency: 1000\n",
       "\n",
       "    domain:\n",
       "      domain_size:\n",
       "        x: 456\n",
       "        y: 152\n",
       "        z: 152\n",
       "      block_size: 8\n",
       "      bodies:\n",
       "        cylinder:\n",
       "          lnas_path: fixture/lnas/basic/cylinder.lnas\n",
       "          transformation:\n",
       "            scale: [72, 72, 72]\n",
       "            translation: [-4, 4, 4]\n",
       "\n",
       "    data:\n",
       "      divergence: { frequency: 50 }\n",
       "      instantaneous:\n",
       "        default: { interval: { frequency: 200000 }, macrs: [rho, u] }\n",
       "      statistics:\n",
       "        interval: { frequency: 100, start_step: 200000 }\n",
       "        macrs_1st_order: [rho, u]\n",
       "        macrs_2nd_order: [u]\n",
       "        exports:\n",
       "          default: { interval: { frequency: 50000 } }\n",
       "\n",
       "    models:\n",
       "      precision:\n",
       "        default: single\n",
       "\n",
       "      LBM:\n",
       "        tau: 0.504081632653061\n",
       "        F:\n",
       "          x: 3.21368E-07\n",
       "          y: 0\n",
       "          z: 0\n",
       "        vel_set: D3Q27\n",
       "        coll_oper: RRBGK\n",
       "      initialization:\n",
       "        vtm_filename: "../nassuArtifacts/macrs/turbulent_pipe.vtm"\n",
       "\n",
       "      engine:\n",
       "        name: CUDA\n",
       "\n",
       "      BC:\n",
       "        periodic_dims: [true, false, false]\n",
       "        BC_map:\n",
       "          - pos: N\n",
       "            BC: RegularizedHWBB\n",
       "            wall_normal: N\n",
       "            order: 1\n",
       "\n",
       "          - pos: S\n",
       "            BC: RegularizedHWBB\n",
       "            wall_normal: S\n",
       "            order: 1\n",
       "\n",
       "          - pos: F\n",
       "            BC: RegularizedHWBB\n",
       "            wall_normal: F\n",
       "            order: 2\n",
       "\n",
       "          - pos: B\n",
       "            BC: RegularizedHWBB\n",
       "            wall_normal: B\n",
       "            order: 2\n",
       "\n",
       "      IBM:\n",
       "        forces_accomodate_time: 1000\n",
       "        body_cfgs:\n",
       "          default: {}\n",
       "
\n" ], "text/latex": [ "\\begin{Verbatim}[commandchars=\\\\\\{\\}]\n", "\\PY{n+nt}{simulations}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZhy{}}\\PY{+w}{ }\\PY{n+nt}{name}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{periodicTurbulentPipe}\n", "\\PY{+w}{ }\\PY{n+nt}{save\\PYZus{}path}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{./tests/validation/results/05\\PYZus{}turbulent\\PYZus{}pipe\\PYZus{}flow/periodic}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{n\\PYZus{}steps}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{400000}\n", "\\PY{+w}{ }\\PY{c+c1}{\\PYZsh{} u* = 0.0034013 / R = 72 / ETT = R/u* = 18,820}\n", "\\PY{+w}{ }\\PY{c+c1}{\\PYZsh{} Re\\PYZus{}tau = 180}\n", "\\PY{+w}{ }\\PY{c+c1}{\\PYZsh{} y+ = 2.5}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{report}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{frequency}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{1000}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{domain}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{domain\\PYZus{}size}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{x}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{456}\n", "\\PY{+w}{ }\\PY{n+nt}{y}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{152}\n", "\\PY{+w}{ }\\PY{n+nt}{z}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{152}\n", "\\PY{+w}{ }\\PY{n+nt}{block\\PYZus{}size}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{8}\n", "\\PY{+w}{ }\\PY{n+nt}{bodies}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{cylinder}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{lnas\\PYZus{}path}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{fixture/lnas/basic/cylinder.lnas}\n", "\\PY{+w}{ }\\PY{n+nt}{transformation}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{scale}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{[}\\PY{n+nv}{72}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{72}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{72}\\PY{p+pIndicator}{]}\n", "\\PY{+w}{ }\\PY{n+nt}{translation}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{[}\\PY{n+nv}{\\PYZhy{}4}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{4}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{4}\\PY{p+pIndicator}{]}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{data}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{divergence}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{n+nt}{ frequency}\\PY{p}{:}\\PY{+w}{ }\\PY{n+nv}{50}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZcb{}}\n", "\\PY{+w}{ }\\PY{n+nt}{instantaneous}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{default}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{n+nt}{ interval}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{n+nt}{ frequency}\\PY{p}{:}\\PY{+w}{ }\\PY{n+nv}{200000}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZcb{}}\\PY{p+pIndicator}{,}\\PY{n+nt}{ macrs}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{[}\\PY{n+nv}{rho}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{u}\\PY{p+pIndicator}{]}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZcb{}}\n", "\\PY{+w}{ }\\PY{n+nt}{statistics}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{interval}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{n+nt}{ frequency}\\PY{p}{:}\\PY{+w}{ }\\PY{n+nv}{100}\\PY{p+pIndicator}{,}\\PY{n+nt}{ start\\PYZus{}step}\\PY{p}{:}\\PY{+w}{ }\\PY{n+nv}{200000}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZcb{}}\n", "\\PY{+w}{ }\\PY{n+nt}{macrs\\PYZus{}1st\\PYZus{}order}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{[}\\PY{n+nv}{rho}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{u}\\PY{p+pIndicator}{]}\n", "\\PY{+w}{ }\\PY{n+nt}{macrs\\PYZus{}2nd\\PYZus{}order}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{[}\\PY{n+nv}{u}\\PY{p+pIndicator}{]}\n", "\\PY{+w}{ }\\PY{n+nt}{exports}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{default}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{n+nt}{ interval}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{n+nt}{ frequency}\\PY{p}{:}\\PY{+w}{ }\\PY{n+nv}{50000}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZcb{}}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZcb{}}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{models}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{precision}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{default}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{single}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{LBM}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{tau}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{0.504081632653061}\n", "\\PY{+w}{ }\\PY{n+nt}{F}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{x}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{3.21368E\\PYZhy{}07}\n", "\\PY{+w}{ }\\PY{n+nt}{y}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{0}\n", "\\PY{+w}{ }\\PY{n+nt}{z}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{0}\n", "\\PY{+w}{ }\\PY{n+nt}{vel\\PYZus{}set}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{D3Q27}\n", "\\PY{+w}{ }\\PY{n+nt}{coll\\PYZus{}oper}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{RRBGK}\n", "\\PY{+w}{ }\\PY{n+nt}{initialization}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{vtm\\PYZus{}filename}\\PY{p}{:}\\PY{+w}{ }\\PY{l+s}{\\PYZdq{}}\\PY{l+s}{../nassuArtifacts/macrs/turbulent\\PYZus{}pipe.vtm}\\PY{l+s}{\\PYZdq{}}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{engine}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{name}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{CUDA}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{BC}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{periodic\\PYZus{}dims}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{[}\\PY{n+nv}{true}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{false}\\PY{p+pIndicator}{,}\\PY{+w}{ }\\PY{n+nv}{false}\\PY{p+pIndicator}{]}\n", "\\PY{+w}{ }\\PY{n+nt}{BC\\PYZus{}map}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZhy{}}\\PY{+w}{ }\\PY{n+nt}{pos}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{N}\n", "\\PY{+w}{ }\\PY{n+nt}{BC}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{RegularizedHWBB}\n", "\\PY{+w}{ }\\PY{n+nt}{wall\\PYZus{}normal}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{N}\n", "\\PY{+w}{ }\\PY{n+nt}{order}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{1}\n", "\n", "\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZhy{}}\\PY{+w}{ }\\PY{n+nt}{pos}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{S}\n", "\\PY{+w}{ }\\PY{n+nt}{BC}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{RegularizedHWBB}\n", "\\PY{+w}{ }\\PY{n+nt}{wall\\PYZus{}normal}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{S}\n", "\\PY{+w}{ }\\PY{n+nt}{order}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{1}\n", "\n", "\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZhy{}}\\PY{+w}{ }\\PY{n+nt}{pos}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{F}\n", "\\PY{+w}{ }\\PY{n+nt}{BC}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{RegularizedHWBB}\n", "\\PY{+w}{ }\\PY{n+nt}{wall\\PYZus{}normal}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{F}\n", "\\PY{+w}{ }\\PY{n+nt}{order}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{2}\n", "\n", "\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZhy{}}\\PY{+w}{ }\\PY{n+nt}{pos}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{B}\n", "\\PY{+w}{ }\\PY{n+nt}{BC}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{RegularizedHWBB}\n", "\\PY{+w}{ }\\PY{n+nt}{wall\\PYZus{}normal}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{B}\n", "\\PY{+w}{ }\\PY{n+nt}{order}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{2}\n", "\n", "\\PY{+w}{ }\\PY{n+nt}{IBM}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{forces\\PYZus{}accomodate\\PYZus{}time}\\PY{p}{:}\\PY{+w}{ }\\PY{l+lScalar+lScalarPlain}{1000}\n", "\\PY{+w}{ }\\PY{n+nt}{body\\PYZus{}cfgs}\\PY{p}{:}\n", "\\PY{+w}{ }\\PY{n+nt}{default}\\PY{p}{:}\\PY{+w}{ }\\PY{p+pIndicator}{\\PYZob{}}\\PY{p+pIndicator}{\\PYZcb{}}\n", "\\end{Verbatim}\n" ], "text/plain": [ "simulations:\n", " - name: periodicTurbulentPipe\n", " save_path: ./tests/validation/results/05_turbulent_pipe_flow/periodic\n", "\n", " n_steps: 400000\n", " # u* = 0.0034013 / R = 72 / ETT = R/u* = 18,820\n", " # Re_tau = 180\n", " # y+ = 2.5\n", "\n", " report:\n", " frequency: 1000\n", "\n", " domain:\n", " domain_size:\n", " x: 456\n", " y: 152\n", " z: 152\n", " block_size: 8\n", " bodies:\n", " cylinder:\n", " lnas_path: fixture/lnas/basic/cylinder.lnas\n", " transformation:\n", " scale: [72, 72, 72]\n", " translation: [-4, 4, 4]\n", "\n", " data:\n", " divergence: { frequency: 50 }\n", " instantaneous:\n", " default: { interval: { frequency: 200000 }, macrs: [rho, u] }\n", " statistics:\n", " interval: { frequency: 100, start_step: 200000 }\n", " macrs_1st_order: [rho, u]\n", " macrs_2nd_order: [u]\n", " exports:\n", " default: { interval: { frequency: 50000 } }\n", "\n", " models:\n", " precision:\n", " default: single\n", "\n", " LBM:\n", " tau: 0.504081632653061\n", " F:\n", " x: 3.21368E-07\n", " y: 0\n", " z: 0\n", " vel_set: D3Q27\n", " coll_oper: RRBGK\n", " initialization:\n", " vtm_filename: \"../nassuArtifacts/macrs/turbulent_pipe.vtm\"\n", "\n", " engine:\n", " name: CUDA\n", "\n", " BC:\n", " periodic_dims: [true, false, false]\n", " BC_map:\n", " - pos: N\n", " BC: RegularizedHWBB\n", " wall_normal: N\n", " order: 1\n", "\n", " - pos: S\n", " BC: RegularizedHWBB\n", " wall_normal: S\n", " order: 1\n", "\n", " - pos: F\n", " BC: RegularizedHWBB\n", " wall_normal: F\n", " order: 2\n", "\n", " - pos: B\n", " BC: RegularizedHWBB\n", " wall_normal: B\n", " order: 2\n", "\n", " IBM:\n", " forces_accomodate_time: 1000\n", " body_cfgs:\n", " default: {}" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from IPython.display import Code\n", "\n", "Code(filename=filename)" ] } ], "metadata": { "kernelspec": { "display_name": "nassu", "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.9" }, "orig_nbformat": 4 }, "nbformat": 4, "nbformat_minor": 2 }