{ "cells": [ { "cell_type": "markdown", "id": "8218498b", "metadata": {}, "source": [ "# Phase equilibria\n", "\n", "Two basic approaches are implemented in teqp:\n", "\n", "* Iterative calculations given guess values\n", "* Tracing along iso-curves (constant temperature, etc.) powered by the isochoric thermodynamics formalism" ] }, { "cell_type": "code", "execution_count": 1, "id": "c0b4e863", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:44.218072Z", "iopub.status.busy": "2024-12-12T18:10:44.217612Z", "iopub.status.idle": "2024-12-12T18:10:44.683785Z", "shell.execute_reply": "2024-12-12T18:10:44.683305Z" } }, "outputs": [ { "data": { "text/plain": [ "'0.22.0'" ] }, "execution_count": 1, "metadata": {}, "output_type": "execute_result" } ], "source": [ "import teqp\n", "import numpy as np\n", "import pandas\n", "import matplotlib.pyplot as plt\n", "teqp.__version__" ] }, { "cell_type": "markdown", "id": "e57e532b", "metadata": {}, "source": [ "## Iterative Phase Equilibria" ] }, { "cell_type": "markdown", "id": "1baf38e1", "metadata": {}, "source": [ "### Pure fluid\n", "\n", "For a pure fluid, phase equilibrium between two phases is defined by equating the pressures and Gibbs energies in the two phases. This represents a 2D non-linear rootfinding problem. Newton's method can be used for the rootfinding, and in teqp, automatic differentiation is used to obtain the necessary Jacobian matrix so the implementation is quite efficient.\n", "\n", "The method requires guess values, which are the densities of the liquid and vapor densities. In some cases, ancillary or superancillary equations have been developed which provide curves of guess densities as a function of temperature.\n", "\n", "For a pure fluid, you can use the ``pure_VLE_T`` method to carry out the iteration." ] }, { "cell_type": "raw", "id": "f0ca0b22", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The Python method is here: :py:meth:`~teqp.teqp.AbstractModel.pure_VLE_T`" ] }, { "cell_type": "code", "execution_count": 2, "id": "2674227c", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:44.685960Z", "iopub.status.busy": "2024-12-12T18:10:44.685349Z", "iopub.status.idle": "2024-12-12T18:10:44.692303Z", "shell.execute_reply": "2024-12-12T18:10:44.691774Z" } }, "outputs": [ { "data": { "text/plain": [ "'guess:'" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/plain": [ "[12735.311173407898, 752.4082303122791]" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/plain": [ "array([12735.31117341, 752.40823031])" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Instantiate the model\n", "model = teqp.canonical_PR([300], [4e6], [0.1])\n", "\n", "T = 250 # [K], Temperature to be used\n", "\n", "# Here we use the superancillary to get guess values (actually these are more \n", "# accurate than the results we will obtain from iteration!)\n", "rhoL0, rhoV0 = model.superanc_rhoLV(T)\n", "display('guess:', [rhoL0, rhoV0])\n", "\n", "# Carry out the iteration, return the liquid and vapor densities\n", "# The guess values are perturbed to make sure the iteration is actually\n", "# changing the values\n", "model.pure_VLE_T(T, rhoL0*0.98, rhoV0*1.02, 10)" ] }, { "cell_type": "markdown", "id": "f8805ae1", "metadata": {}, "source": [ "### Binary Mixture" ] }, { "cell_type": "markdown", "id": "76bccf19", "metadata": {}, "source": [ "For a binary mixture, the approach is roughly similar to that of a pure fluid. The pressure is equated between phases, and the chemical potentials of each component in each phase are forced to be the same. \n", "\n", "Again, the user is required to provide guess values, in this case molar concentrations in each phase, and a Newton method is implemented to solve for the phase equilibrium. The analytical Jacobian is obtained from automatic differentiation.\n", "\n", "The ``mix_VLE_Tx`` function is the binary mixture analog to ``pure_VLE_T`` for pure fluids." ] }, { "cell_type": "raw", "id": "eef189fd", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The Python method is here: :py:meth:`~teqp.teqp.AbstractModel.mix_VLE_Tx`" ] }, { "cell_type": "code", "execution_count": 3, "id": "b12bd318", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:44.694017Z", "iopub.status.busy": "2024-12-12T18:10:44.693736Z", "iopub.status.idle": "2024-12-12T18:10:44.699081Z", "shell.execute_reply": "2024-12-12T18:10:44.698583Z" } }, "outputs": [ { "data": { "text/plain": [ "(,\n", " array([ 128.66049209, 12737.38871682]),\n", " array([ 12.91868229, 1133.77242677]))" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "zA = np.array([0.01, 0.99])\n", "model = teqp.canonical_PR([300,310], [4e6,4.5e6], [0.1, 0.2])\n", "model1 = teqp.canonical_PR([300], [4e6], [0.1])\n", "T = 273.0 # [K]\n", "# start off at pure of the first component\n", "rhoL0, rhoV0 = model1.superanc_rhoLV(T)\n", "\n", "# then we shift to the given composition in the first phase\n", "# to get guess values\n", "rhovecA0 = rhoL0*zA\n", "rhovecB0 = rhoV0*zA\n", "\n", "# carry out the iteration\n", "code, rhovecA, rhovecB = model.mix_VLE_Tx(T, rhovecA0, rhovecB0, zA, \n", " 1e-10, 1e-10, 1e-10, 1e-10, # stopping conditions\n", " 10 # maximum number of iterations\n", " )\n", "code, rhovecA, rhovecB" ] }, { "cell_type": "markdown", "id": "d72ef08e", "metadata": {}, "source": [ "You can (and should) check the value of the return code to make sure the iteration succeeded. Do not rely on the numerical value of the enumerated return codes!" ] }, { "cell_type": "markdown", "id": "f4e3f914", "metadata": {}, "source": [ "# Tracing (isobars and isotherms)\n", "\n", "When it comes to mixture thermodynamics, as soon as you add another component to a pure component to form a binary mixture, the complexity of the thermodynamics entirely changes. For that reason, mixture iterative calculations for mixtures are orders of magnitude more difficult to carry out. Asymmetric mixtures can do all sorts of interesting things that are entirely unlike those of pure fluids, and the algorithms are therefore much, much more complicated. Formulating phase equilibrium problems is not much more complicated than for pure fluids, but the most challenging aspect is to obtain good guess values from which to start an iterative routine, and the difficulty of this problem increases with the complexity of the mixture thermodynamics.\n", "\n", "Ulrich Deiters and Ian Bell have developed a number of algorithms for tracing phase equilibrium solutions as the solution of ordinary differential equations rather than carrying out iterative routines for a given state point. The advantage of the tracing calculations is that they can often be initiated at a state point that is entirely known, for instance the pure fluid endpoint for a subcritical isotherm or isobar." ] }, { "cell_type": "raw", "id": "e0097771", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The Python method is here: :py:meth:`~teqp.teqp.AbstractModel.trace_VLE_isotherm_binary`" ] }, { "cell_type": "markdown", "id": "63902dba", "metadata": {}, "source": [ "The C++ implementation returns a string in JSON format, which can be conveniently operated upon, for instance after converting the returned data structure to a ``pandas.DataFrame``. A simple example of plotting a subcritical isotherm for a \"boring\" mixture is presented here:" ] }, { "cell_type": "code", "execution_count": 4, "id": "49dcba2b", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:44.700829Z", "iopub.status.busy": "2024-12-12T18:10:44.700538Z", "iopub.status.idle": "2024-12-12T18:10:44.719017Z", "shell.execute_reply": "2024-12-12T18:10:44.718611Z" } }, "outputs": [ { "data": { "text/plain": [ "\"[{'T / K': 273.0, 'c': -1.0, 'drho/dt': [-0.618312383229212, 0.7690760182230469, -0.1277526773161415...\"" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
T / Kcdrho/dtdtpL / PapV / ParhoL / mol/m^3rhoV / mol/m^3txL_0 / mole frac.xV_0 / mole frac.
0273.0-1.0[-0.618312383229212, 0.7690760182230469, -0.12...0.0000102.203397e+062.203397e+06[10697.985891540735, 0.0][1504.6120879290752, 0.0]0.0000001.01.0
1273.0-1.0[-0.6183123817120353, 0.7690760162922189, -0.1...0.0000452.203397e+062.203397e+06[10697.985885357639, 7.690760309421386e-06][1504.6120866515366, 9.945415375682985e-07]0.0000101.01.0
2273.0-1.0[-0.6183123827116788, 0.7690760173388914, -0.1...0.0002032.203397e+062.203397e+06[10697.98585753358, 4.229918121248511e-05][1504.6120809026731, 5.469978386095445e-06]0.0000551.01.0
\n", "
" ], "text/plain": [ " T / K c drho/dt dt \\\n", "0 273.0 -1.0 [-0.618312383229212, 0.7690760182230469, -0.12... 0.000010 \n", "1 273.0 -1.0 [-0.6183123817120353, 0.7690760162922189, -0.1... 0.000045 \n", "2 273.0 -1.0 [-0.6183123827116788, 0.7690760173388914, -0.1... 0.000203 \n", "\n", " pL / Pa pV / Pa rhoL / mol/m^3 \\\n", "0 2.203397e+06 2.203397e+06 [10697.985891540735, 0.0] \n", "1 2.203397e+06 2.203397e+06 [10697.985885357639, 7.690760309421386e-06] \n", "2 2.203397e+06 2.203397e+06 [10697.98585753358, 4.229918121248511e-05] \n", "\n", " rhoV / mol/m^3 t xL_0 / mole frac. \\\n", "0 [1504.6120879290752, 0.0] 0.000000 1.0 \n", "1 [1504.6120866515366, 9.945415375682985e-07] 0.000010 1.0 \n", "2 [1504.6120809026731, 5.469978386095445e-06] 0.000055 1.0 \n", "\n", " xV_0 / mole frac. \n", "0 1.0 \n", "1 1.0 \n", "2 1.0 " ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model = teqp.canonical_PR([300,310], [4e6,4.5e6], [0.1, 0.2])\n", "model1 = teqp.canonical_PR([300], [4e6], [0.1])\n", "T = 273.0 # [K]\n", "rhoL0, rhoV0 = model1.superanc_rhoLV(T) # start off at pure of the first component\n", "j = model.trace_VLE_isotherm_binary(T, np.array([rhoL0, 0]), np.array([rhoV0, 0]))\n", "display(str(j)[0:100]+'...') # The first few bits of the data\n", "df = pandas.DataFrame(j) # Now as a data frame\n", "df.head(3)" ] }, { "cell_type": "code", "execution_count": 5, "id": "9aecca78", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:44.720525Z", "iopub.status.busy": "2024-12-12T18:10:44.720258Z", "iopub.status.idle": "2024-12-12T18:10:44.883572Z", "shell.execute_reply": "2024-12-12T18:10:44.883024Z" } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkAAAAG0CAYAAADacZikAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjkuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/GU6VOAAAACXBIWXMAAA9hAAAPYQGoP6dpAABnd0lEQVR4nO3dZ3RU1R6G8WfSAyQhENIgEHrvIr2ICNIEQSkqiGK7hiZWUERBDaKgNEEUCCpF6YiIIlKkS++hE1pCT0JC6pz7YTQaSZQAmclk3t9aWcs9s8/Jf87VzHvP3mdvk2EYBiIiIiIOxMnWBYiIiIhYmwKQiIiIOBwFIBEREXE4CkAiIiLicBSARERExOEoAImIiIjDUQASERERh+Ni6wLyIrPZzLlz5/Dy8sJkMtm6HBEREbkFhmEQHx9PcHAwTk7/fo9HASgL586dIyQkxNZliIiIyG04ffo0JUqU+Nc+CkBZ8PLyAiwX0Nvb28bViIiIyK2Ii4sjJCQk43v83ygAZeHPYS9vb28FIBERETtzK9NXNAlaREREHI4CkIiIiDgcBSARERFxOApAIiIi4nAUgERERMThKACJiIiIw1EAEhEREYejACQiIiIORwFIREREHI4CkIiIiDgcBSARERFxOApAIiIi4nAUgERERMSq9q9bxI3EBJvWoAAkIiIiVpGUlMSGz16g6q992DvtfzatxaYBKDw8nHr16uHl5YW/vz+dO3cmMjLyX4/54osvaNq0Kb6+vvj6+tKqVSu2bt2aqY9hGLz99tsEBQXh6elJq1atOHLkSG5+FBEREfkXhw5HcuSj+2h8YQ4Aac6eGOZ0m9Vj0wC0du1awsLC2Lx5MytXriQ1NZXWrVuTkJD9bbE1a9bQs2dPVq9ezaZNmwgJCaF169acPXs2o8/o0aMZP348U6ZMYcuWLRQsWJA2bdqQlJRkjY8lIiIif0g3GyxdNJtis1pRPf0A1ynA/qaTaPTi55icnG1Wl8kwDMNmv/0fLl68iL+/P2vXrqVZs2a3dEx6ejq+vr5MnDiR3r17YxgGwcHBvPzyy7zyyisAxMbGEhAQQEREBD169LjpHMnJySQnJ2e04+LiCAkJITY2Fm9v77vz4URERBzMmSvX2TD9DR6N/wYnk8Fpt7J49Z5N4RKVcuX3xcXF4ePjc0vf33lqDlBsbCwARYoUueVjEhMTSU1NzTjmxIkTREdH06pVq4w+Pj4+1K9fn02bNmV5jvDwcHx8fDJ+QkJC7uBTiIiIODbDMPhh815OjmtH9+tf42QyOB7SlRKvrM+18JNTeSYAmc1mBg0aROPGjalWrdotH/f6668THBycEXiio6MBCAgIyNQvICAg471/GjJkCLGxsRk/p0+fvs1PISIi4thiE1MZM30WtX7sRBPTbpJx51KrTyjTdzomtwK2Li+Di60L+FNYWBj79u1j/fr1t3zMqFGjmDt3LmvWrMHDw+O2f7e7uzvu7u63fbyIiIjAxqMX2TLnAwamzcTVlM5Vz5J49ZqNX3B1W5d2kzxxB6hfv34sW7aM1atXU6JEiVs65uOPP2bUqFH8/PPP1KhRI+P1wMBAAGJiYjL1j4mJyXhPRERE7p7ktHQ+/n4b12Y+xkvp0y3hp3QHfAdtxCUPhh+wcQAyDIN+/fqxaNEifv31V0qXLn1Lx40ePZqRI0eyYsUK7rnnnkzvlS5dmsDAQFatWpXxWlxcHFu2bKFhw4Z3tX4RERFHFxkdz6Bx39D198do57yVNJMLyQ+E49v7G3D3snV52bLpEFhYWBizZ89myZIleHl5ZczR8fHxwdPTE4DevXtTvHhxwsPDAfjwww95++23mT17NqGhoRnHFCpUiEKFCmEymRg0aBDvvfce5cuXp3Tp0gwbNozg4GA6d+5sk88pIiKS35jNBjM3nSRyxRQ+cZqGh1MqNwoE4dnza1xC6tm6vP9k0wA0efJkAFq0aJHp9RkzZtCnTx8AoqKicHJyynRMSkoKjzzySKZjhg8fzjvvvAPAa6+9RkJCAs899xzXrl2jSZMmrFix4o7mCYmIiIjFhbgkhnz3O61PfswolzUAJIe2xLPbNChw609y21KeWgcor8jJOgIiIiKOZMW+aCYv+IlR6R9T2SkKM06Y7huKqenL4GTbqcU5+f7OM0+BiYiISN6VkJzGu9/vJ27HQr5x/Rwvpxukefrh8ug0KNPC1uXlmAKQiIiI/KsdUVd5de7vPBY3nb5uPwJgDmmIy6MzwDvIxtXdHgUgERERyVJaupmJq48y/9ctjHMZR12XPzYWbzQAp/uHg7P9xgj7rVxERERyzanLCQz6dhdeZ9ay1HUSRUzXMdy9MT38OVRqZ+vy7pgCkIiIiGQwDIN5288wculenjHPo7/bIpwwIKgmpkdnQpFbW7Mvr1MAEhEREQCuJqQwZOFeft8fyWeuk2jqss/yRt2n4MFR4Jp/lpNRABIRERHWH7nEy/N2USJ+Dz+4jyfQdBXDtQCmDp9Cze62Lu+uUwASERFxYEmp6Xz0UyTT1h/nGeflDHGfgzNm8KuAqdtX4F/Z1iXmCgUgERERBxUZHc/AuTs5Gx3DFNfPedD5d8sb1R6BjuPAvZBtC8xFCkAiIiIOxmw2iNh4klErDlE+/Tg/eoyjBDHg7AYPhsM9fcFksnWZuUoBSERExIFciEvi5Xm7+e3IRXo6/8oI969wJRUKl4RHZ0LxOrYu0SoUgERERBzET/ujeWPBHpIS4/nUbQadnX6zvFGhLTw8GTx9bVugFSkAiYiI5HMJyWmMXHaAub+fpqzpLNMKTiA0PQpMznD/29BogM03MrU2BSAREZF8bPfpawycu5OTlxPp6LyRMe7TcEu/AYUC4JEZENrY1iXahAKQiIhIPpRuNpi85iif/nIEJ3MKHxeYyyPmH8EMhDaFrtPAK8DWZdqMApCIiEg+c/pKIoO/28XvJ69SwnSRb3wmE5p8yPJm05ehxVC73sj0bnDsTy8iIpLPLN55lmGL9xGfnEY79z186vYZbslxlgnOD0+FCq1tXWKeoAAkIiKSD8TeSGXY4n0s3X0OZ9IZU2QZXRO/hVSgeF14NMLyqLsACkAiIiJ2b8vxywz+bjdnr90gwCmWBcW+oETsDsub9z4Prd8DFzfbFpnHKACJiIjYqZQ0M5/+cpjJa49hGNCp8DE+No3HNfYiuBWChyZAtS62LjNPUgASERGxQ8cvXmfQt7vYcyYWE2Yml1rHgxe+xGSYwb8KdPsK/Mrbusw8SwFIRETEjhiGwdzfTzPi+wPcSE0nxCOJBYEz8Y9ea+lQ8zFoPwbcCti20DxOAUhERMROXElI4fUFe1h5IAaAXiUv8k7SxzhHnwYXD2j3EdTule83Mr0bFIBERETswLrDF3l53m4uxifj6gwRVXfT6OhYTOZU8C1tGfIKqmHrMu2GApCIiEgelpSazugVkUzfcAKA6sWc+KbYN/gcXmbpULkjdJoEHj42rNL+KACJiIjkUZHR8Qycu5ND0fEAvFIrjRcvjMTp+BFwcoEHRkKD/2nI6zYoAImIiOQxhmEQsfEk4T8eIiXNTNGCbnx9z1GqbH8X0m6Ad3HLRqYl69u6VLulACQiIpKHXIhP4tV5e1h7+CIArct7M85nDp5bZlk6lL0funwBBYvasEr7pwAkIiKSR/xyIIbXFuzhSkIK7i5OjGpRgM5H3sC0bx9ggvuGQtNXwMnJ1qXaPQUgERERG7uRks77yw/wzeYoACoFejG9/nmC1zwHyXFQwA+6fgll77NxpfmHApCIiIgN7Tsby6Bvd3H0wnUAnmtUgtdcZuPy0xRLh5IN4ZHp4B1swyrzHwUgERERGzCbDb5cf5yPfookNd3A38udCe2LUX/bK3Bmq6VTowFw/9vg7GrbYvMhBSAREREri45N4uV5u9hw9DIAD1QJYGydi3j90BluXLGs6dN5MlRqb9tC8zEFIBEREStasS+aNxbu4VpiKp6uzgzvUJHuCbMwzf8YMCCopmVVZ99QW5earykAiYiIWEFiShojvj/A3N9PA1C9uA8THipO6JqBcOKPjUzveRrahIOrhw0rdQwKQCIiIrlsz5lrDJq7i+OXEjCZ4IXmZRlc4RKu8x6E69HgWhA6joMaj9q6VIehACQiIpJL0s0GU9cdZ8zPkaSZDYJ8PBj7aA0aRs+Cr0eAkQ5+FaH711Csoq3LdSgKQCIiIrng3LUbDP5uF5uPXwGgXfVAwtuG4LOiPxz+0dKpejfo8Am4F7JhpY7JpktJhoeHU69ePby8vPD396dz585ERkb+6zH79++na9euhIaGYjKZ+PTTT2/q884772AymTL9VKpUKZc+hYiISGbL956n7bjf2Hz8CgXcnBn9SA0mtTDh89X9lvDj7GYJPl2mKvzYiE0D0Nq1awkLC2Pz5s2sXLmS1NRUWrduTUJCQrbHJCYmUqZMGUaNGkVgYGC2/apWrcr58+czftavX58bH0FERCRDQnIar87bzYuzdhB7I5WaJXxY3r8J3cwrME1vA9eioHAp6LvSMuFZu7jbjE2HwFasWJGpHRERgb+/P9u3b6dZs2ZZHlOvXj3q1asHwBtvvJHtuV1cXP41IImIiNxNu05fY9DcnZy8nIjJBGEtyjGwaSCuywfAvgWWTpU6QKdJ4FnYprVKHpsDFBsbC0CRIkXu+FxHjhwhODgYDw8PGjZsSHh4OCVLlsyyb3JyMsnJyRntuLi4O/79IiLiGNLNBlPWHuOTlYdJMxsE+3jwSfda1C8YA9Puh8tHwOQMD7wLDfvprk8ekWe2kzWbzQwaNIjGjRtTrVq1OzpX/fr1iYiIYMWKFUyePJkTJ07QtGlT4uPjs+wfHh6Oj49Pxk9ISMgd/X4REXEMZ6/doOcXm/noJ8tTXu1rBPHjwGbUj/sZvmhpCT9ewfDUcmjUX+EnD8kzd4DCwsLYt2/fXZmr07Zt24x/rlGjBvXr16dUqVJ899139O3b96b+Q4YMYfDgwRntuLg4hSAREflXy/acY+jCvcQlpVHQzZl3O1Wja/UimH58CXZ+belU5j7LLu4F/WxbrNwkTwSgfv36sWzZMtatW0eJEiXu+vkLFy5MhQoVOHr0aJbvu7u74+7uftd/r4iI5D/Xk9N4Z+l+5m8/A0DNkMKM616LUFM0TGsNMXsBEzR/HZq/Bk7Oti1YsmTTAGQYBv3792fRokWsWbOG0qVL58rvuX79OseOHaNXr165cn4REXEMu05fY+DcnZz6+0TnVuVxjfweFodBSjwU8IOuX0DZlrYuV/6FTQNQWFgYs2fPZsmSJXh5eREdHQ2Aj48Pnp6eAPTu3ZvixYsTHh4OQEpKCgcOHMj457Nnz7Jr1y4KFSpEuXLlAHjllVfo2LEjpUqV4ty5cwwfPhxnZ2d69uxpg08pIiL2Lt1sMHnNUT755QjpZoPihT0Z260m9Ut6wco3YfNnlo4hDeDRGeAdbNuC5T+ZDMMwbPbLs5kMNmPGDPr06QNAixYtCA0NJSIiAoCTJ09meaeoefPmrFmzBoAePXqwbt06Ll++TLFixWjSpAnvv/8+ZcuWvaW64uLi8PHxITY2Fm9v7xx/LhERyT/OXrvBS9/uYusJy4rOHWoE8f7D1fFJiYF5feDM75aOjfrD/cPB2dV2xTq4nHx/2zQA5VUKQCIiAtlMdK5THNPRVbDwWbhxBTx8oPNkqNTe1uU6vJx8f+eJSdAiIiJ5yT8nOtcKKcy4HrUo5esBq9+HdR8DBgTVgm4zwTfUluXKbVAAEhER+ZtsJzonXoSvusPJ3ywd6z0Drd8HVw/bFiy3RQFIRESEf1nRuUxROLke5j8N12PAtSA8NB6qP2LrkuUOKACJiIjDO/fHROctf0x0bl8jiA86V8fHwxl+Gwu/jgTDDMUqQ7evoFgFG1csd0oBSEREHNqPe8/zxsK9xN5IpYCbM+8+VJVH6pbAdOMqzHkBjvxk6VijB3QYC24FbVuw3BUKQCIi4pASU9J4d+kBvt12GoAaJXwY16M2pf0KwpntMO9JiD0Nzu7Q7iOo01t7eeUjCkAiIuJw9p6JZeDcnRy/lIDJBP9rXpaXHqiAq5MJtnwOP70J5lTwLW0Z8gqqYeuS5S5TABIREYdhNhtM/e04Y36OJDXdINDbMtG5YdmikBQHS/vDgcWWzpU7QqdJlnV+JN9RABIREYcQE5fE4O92seHoZQDaVgskvEt1Chdwg+h98F1vuHIMnFzggZHQ4H8a8srHFIBERCTf+3l/NK8v2MPVxFQ8XZ0Z3rEK3euFWLZk2vE1LH8F0pLAuzg8GgEh99q6ZMllCkAiIpJv3UhJ570fDjBrSxQAVYO9Gd+zNmWLFYKUREvw2TXL0rncA/Dw51CwqA0rFmtRABIRkXzpwLk4BszdydEL1wF4rlkZXm5dAXcXZ7h0xDLkdeEAmJzgvjehyWBwcrJx1WItCkAiIpKvGIbBjA0nGfXjIVLSzfh7uTOmW02ali9m6bBvASwdACnXoaA/PDINSjezbdFidQpAIiKSb1yMT+bV+btZE3kRgFaV/fmwaw2KFnKHtGT4aSj8/qWlc6kmlvDjFWjDisVWFIBERCRfWB15gVfn7ebS9RTcXZx4q31lnmhQyjLR+epJmNcHzu20dG76MrQYCs76GnRU+l9eRETsWlJqOh+uOMSMDScBqBToxfietakQ4GXpcGg5LH4BkmLB0xe6fAHlH7BdwZInKACJiIjdOhITT/85OzkUHQ9An0ahvNG2Eh6uzpCeCqtGwMbxls7F77E84l44xHYFS56hACQiInbHMAxmb41i5LIDJKWaKVrQjY8erUHLSgGWDnHnYP7TELXJ0m7wIrR6F1zcbFe05CkKQCIiYleuJqTw+oI9/HwgBoCm5f0Y060m/l4elg7HfoUFz0DiZXD3hk4ToUonG1YseZECkIiI2I2Nxy7x0re7iIlLxtXZxOsPVuLpxqVxcjKBOR3Wjoa1HwIGBFaHR2dC0bK2LlvyIAUgERHJ81LTzXyy8jCT1x7DMKBMsYKM71GbasX/2Kj0+kVY+AwcX2Np1+0DD44CV09blSx5nAKQiIjkaacuJzBgzk52n4kFoEe9EN7uWIUCbn98hZ3aaJnvE38eXAtAh0+hZnfbFSx2QQFIRETyJMMwWLTzLMMW7yMhJR1vDxdGda1Bu+pBlg5ms+UJr1UjwEgHv4rQ7Svwr2TbwsUuKACJiEieE5eUyrDF+1iy6xwA95YuwqfdaxFc+I8hrcQrsPhFOPyjpV29G3T4BNwL2ahisTcKQCIikqdsP3WVQd/u5PSVGzg7mRh0f3levK8czk4mS4ez2+G7PhAbBc7u0PZDy5wfk8mWZYudUQASEZE8Id1s8Nnqo3y66gjpZoMSvp6M61GbuqV8LR0Mw7KP14ohYE4F39LQbSYE1bRt4WKXFIBERMTmzl27wUvf7mLLiSsAPFQzmPceroa3h6ulQ1IcfD8A9i+ytCt1gM6fgYePjSoWe6cAJCIiNrViXzSvL9hD7I1UCro5M6JTNbrUKW7ZxBQgei989yRcOQZOLtD6Paj/goa85I4oAImIiE3cSEln5A8HmL0lCoAaJXwY36M2oX4FLR0MA3Z+DctfhbQk8C5h2csrpJ7tipZ8QwFIRESs7sC5OAbM3cnRC9cBeL55GV5+oCJuLk6WDikJ8MPLsHuOpV2+NTz8ORQoYqOKJb9RABIREasxDIOIjScJ//EQKWlminm580m3WjQp7/dXp4uRliGviwfB5AQt34LGL4GTk+0Kl3xHAUhERKzi8vVkXp2/h18PXQDg/kr+jH6kBkULuf/Vac88+H4gpCZAoQB4ZDqENrFRxZKfKQCJiEiu++3IRQZ/t5uL8cm4uTjxZrvK9G5Y6q+JzqlJsOIN2D7D0i7dDLpOg0L+tita8jUFIBERyTUpaWbG/BzJ5+uOA1DevxDje9amcpD3X52uHIfvelue9sIEzV6FFm+Ak7NtihaHoAAkIiK54sSlBAbO3cmePzYxfbx+Sd5qXwVPt78FmwNLYUkYJMdBgaLQZSqUa2WjisWRKACJiMhdZRgGC3ecZdiSfSSmpOPj6cqHXWvwYLXAvzqlpcDKt2HLZEs7pIFlvo9PcdsULQ5HAUhERO6a+KRU3vrbJqb1Sxfh0x61CPLx/KvTtdMwrw+c3WZpN+oP9w8HZ1frFywOSwFIRETuip1RVxk4dxdRVxKz3sQU4PBPsOh5uHHVso1F5ylQqZ3tihaHpQAkIiJ3xGw2mLLuGGN/Pkya2aB4YU/G9/zbJqYA6Wmw+j1Y/4mlHVwbHp0JvqVsU7Q4PJuuKhUeHk69evXw8vLC39+fzp07ExkZ+a/H7N+/n65duxIaGorJZOLTTz/Nst+kSZMIDQ3Fw8OD+vXrs3Xr1lz4BCIiji0mLoknpm1h9IpI0swGHWoEsXxg08zhJ+48zOz4V/i593l4+ieFH7EpmwagtWvXEhYWxubNm1m5ciWpqam0bt2ahISEbI9JTEykTJkyjBo1isDAwCz7fPvttwwePJjhw4ezY8cOatasSZs2bbhw4UJufRQREYez6mAMD366jo3HLuPp6szorjWY0LM2Pp5/m8tzbDVMaQJRG8HNy7KXV7vR4OKe7XlFrMFkGIZh6yL+dPHiRfz9/Vm7di3NmjX7z/6hoaEMGjSIQYMGZXq9fv361KtXj4kTJwJgNpsJCQmhf//+vPHGGzedJzk5meTk5Ix2XFwcISEhxMbG4u3tfVN/ERFHlpSazqgfDxGx8SQAVYK8mfBYbcoWK/RXJ3M6rPsI1owCDAioDt1mQtGyNqlZHENcXBw+Pj639P2dpzZWiY21rBVRpMjtb3aXkpLC9u3badXqr3UknJycaNWqFZs2bcrymPDwcHx8fDJ+QkJCbvv3i4jkZ0di4uk8aUNG+OnbpDSLwhplDj/XL8I3XWBNOGBAnd7wzEqFH8lT8kwAMpvNDBo0iMaNG1OtWrXbPs+lS5dIT08nICAg0+sBAQFER0dnecyQIUOIjY3N+Dl9+vRt/34RkfzIMAzmbI2i48T1HIqOp2hBN2Y8VY9hHarg7vK3hQ1PbYTPm8LxNeBawPKU10MTwNUz23OL2EKeeQosLCyMffv2sX79eqv/bnd3d9zdNR4tIpKV2MRU3li4hx/3Wf5PZNPyfozpVhN/L4+/OpnNsHEcrBoJRjr4VbQMeflXtlHVIv8uTwSgfv36sWzZMtatW0eJEiXu6Fx+fn44OzsTExOT6fWYmJhsJ02LiEjWfj95hYFzdnIuNglXZxOvtqnIM03K4PT3tX0Sr8CiF+DIT5Z2je7Qfiy4F8r6pCJ5gE2HwAzDoF+/fixatIhff/2V0qVL3/E53dzcqFu3LqtWrcp4zWw2s2rVKho2bHjH5xcRcQTpZoNxvxyh++ebOBebRGjRAiz4XyOea1Y2c/g5sw0+b2YJP87u0HEcPPy5wo/keTa9AxQWFsbs2bNZsmQJXl5eGXN0fHx88PS0jBf37t2b4sWLEx4eDlgmOR84cCDjn8+ePcuuXbsoVKgQ5cqVA2Dw4ME8+eST3HPPPdx77718+umnJCQk8NRTT9ngU4qI2Jdz124w6NtdbD1xBYAudYozolM1Crn/7SvDMGDzZMt+XuZUKFLGsrBhUA0bVS2SMzZ9DN5kMmX5+owZM+jTpw8ALVq0IDQ0lIiICABOnjyZ5Z2i5s2bs2bNmoz2xIkT+eijj4iOjqZWrVqMHz+e+vXr31JdOXmMTkQkP1mxL5rXF+wh9kYqBd2cee/hajxc+x9TE25cs+zgfmiZpV2ls2Wis4f+Xopt5eT7O0+tA5RXKACJiKNJSk3nvR8O8M3mKABqlPBhfI/ahPoVzNzx3C6Y9yRcPQlOrtDmA7j3Wcjm/9CKWFNOvr/zxCRoERGxncMx8fSfvZPImHgAnm9WhpdbV8TN5W/TRA0Dtk2DFUMgPQUKl7Ss6ly8rm2KFrlDCkAiIg7KMAxmb41ixPcHSE4z41fInbHdatKsQrHMHZPj4fuBsG+BpV2xHXT+DDx9bz6piJ1QABIRcUDXElN4Y8FeVuy3PHzSrEIxxjxak2Je/1gTLXqfZcjr8lFwcoFW70DDfhryErunACQi4mC2nrjCoLl/re3zWptK9G1SOvPj7YYBO7+B5a9AWhJ4F4dHZkDJW3uYRCSvUwASEXEQ6WaDib8eZdyqw5gNCC1agPE9a1OjROHMHVMS4IdXYPdsS7tcK3h4KhQsavWaRXKLApCIiAO4pbV9AC5Gwne94eIhMDlBy7eg8UvglGe2jhS5KxSARETyuZ/3R/Pq/P9Y2wdg97ewbBCkJkKhAHhkOoQ2sXq9ItagACQikk8lpabzwfKDfLXpFPAva/uk3oAfX4cdMy3t0s2h65dQyN/KFYtYjwKQiEg+dPRCPP1m7+RQtGVtn+ealeGVf67tA3DpKMzrAzF7ARO0eAOavQpOzlavWcSaFIBERPIRwzD4bttp3ll6gBup6RQt6MaYbjVpUTGLuzn7FsLSAZASDwWLQZcvoOx91i9axAYUgERE8om4pFSGLtzLsj3nAWhSzo+x3Wvi7+WRuWNaMvw0FH7/0tIu1Ri6TgPvICtXLGI7CkAiIvnAzqirDJi7k9NXbuDiZOLl1hV5vlmZzGv7AFw5YRnyOr/L0m4yGO57E5z1dSCORf/Gi4jYMbPZYMq6Y4z9+TBpZoMSvp6M71mbOiWz2Kbi4PewOAySY8GzCHSZCuUfsH7RInmAApCIiJ26EJ/E4G93s/7oJQA61Ajigy7V8fZwzdwxLQV+GQ6bP7O0Q+pbHnH3yeJReBEHoQAkImKH1kRe4OXvdnM5IQUPVyfefagq3e4JwfTPPbquRcG8p+DsNku7UX+4fzg4u958UhEHogAkImJHUtLMfPTTIb747QQAlQK9mPhYbcr5e93cOXIFLHoekq6BR2F4eApUbGvVekXyKgUgERE7cfJSAgPm7mTPmVgAejcsxdB2lfFw/ceaPemp8OtI2DDO0i5e17KRqW8pK1cskncpAImI2IElu84ydOFeElLS8fF0ZfQjNWhTNfDmjrFnYUFfiNpkadf/HzwwAlzcrFuwSB6nACQikoclJKcxfOl+5m8/A8C9oUX4tEctggt73tz5yC+w6DlIvAzu3tBpIlTpZOWKReyDApCISB514Fwc/ebs4PjFBJxM0K9leQa0LIeL8z+2s0hPgzUfwG9jLO3AGtBtJhQpY/2iReyEApCISB5jGAZfbz7Fez8cJCXNTIC3O592r03DskVv7hx3HhY8A6fWW9r1noHW74Orx819RSSDApCISB5yLTGF1+bv4ecDMQDcX8mfjx6tSZGCWczhObYaFj4LCRfBzQseGgfVulq5YhH7pAAkIpJH/H7yCgPn7ORcbBKuziaGtK3MU41Db17bx5wOa0fD2g8BAwKqW4a8ipa1Sd0i9kgBSETExtLNBp+tPsonvxzGbEBpv4JM6FmbasV9bu4cHwMLn4ET6yztun3gwVHgmsWkaBHJlgKQiIgNxcQlMWjuLjYdvwxAl9rFGdG5GoXcs/jzfGKdZb7P9RhwLQgdP4Ua3axbsEg+oQAkImIjqw9d4OV5u7mSkEIBN2dGdqpG17pZ7M9lNlue8FrzARhm8K8Cj86EYhWsX7RIPqEAJCJiZSlpZkavOMSX6y3bWVQJ8mbiY7UpU6zQzZ2vX7RMdD6+2tKu/QS0/QjcClixYpH8RwFIRMSK/rmdRZ9GoQxpVwl3F+csOm+wrOocfx5cPKHDWKj1mJUrFsmfFIBERKxkya6zvLloH9eT0yhcwJWPHqnJA1UCbu5oNsOGT+DX98FIB7+Klqe8/Ctbv2iRfEoBSEQklyWmpPHO0v18t+2v7SzG9axFkE8WT24lXLbs4H50paVdo4flzo9bQStWLJL/KQCJiOSig+fj6Dd7B8f+2M6if8vy9M9qOwuAqM0w/2mIOwsuHtDuI6jdC/65DpCI3DEFIBGRXGAYBrO2RDFi2YH/3s7CbIaN42HVCMuQV9Fylqe8AqtZv3ARB6EAJCJyl8XeSOWNBXv4cV80AC0r+fNxdttZJF6BRS/AkZ8s7eqPQodPwN3LihWLOB4FIBGRu2hH1FX6z97J2Ws3cHU28fqDlejbpPTN21kAnN4K856CuDPg7A5tP7Ss7KwhL5FcpwAkInIXmM0GU387zsc/RZJmNihZpAATetamZkjhmzsbBmyaCL+8A+Y0KFLW8pRXYHVrly3isG47AB04cICoqChSUlIyvf7QQw/dcVEiIvbk0vVkBn+3m3WHLwLQoUYQH3SpjreH682dE6/A4hfh8I+WdrWu0HGchrxErCzHAej48eM8/PDD7N27F5PJhGEYABm3d9PT0+9uhSIiediGo5cY9O0uLsYn4+HqxDsdq9K9XkjWQ15ntlmGvGKjLENeD4bDPU9ryEvEBrJ4DvPfDRw4kNKlS3PhwgUKFCjA/v37WbduHffccw9r1qzJhRJFRPKetHQzH/8UyRPTtnAxPpkKAYVY2q8JPe4teXP4MQzYNAmmt7GEnyJl4JmVUK+vwo+IjeQ4AG3atIkRI0bg5+eHk5MTTk5ONGnShPDwcAYMGJCjc4WHh1OvXj28vLzw9/enc+fOREZG/udx8+bNo1KlSnh4eFC9enWWL1+e6f0+ffpgMpky/Tz44IM5qk1EJDvnrt2gx9TNTFx9FMOAnveWZElYEyoEZDGMdeMqzH0cfhpqme9TpTM8txaCalq9bhH5S44DUHp6Ol5elv/I/fz8OHfuHAClSpW6pfDyd2vXriUsLIzNmzezcuVKUlNTad26NQkJCdkes3HjRnr27Enfvn3ZuXMnnTt3pnPnzuzbty9TvwcffJDz589n/MyZMyeHn1RE5GYrD8TQdtxvbDt1FS93Fyb0rE14l+p4umWxl9eZ7fB5M4j8AZzdoN3H8GgEeHhbvW4RySzHc4CqVavG7t27KV26NPXr12f06NG4ubkxdepUypQpk6NzrVixIlM7IiICf39/tm/fTrNmzbI8Zty4cTz44IO8+uqrAIwcOZKVK1cyceJEpkyZktHP3d2dwMDAHH46EZGsJaelE778EBEbTwJQs4QPE3rWoWTRLHZlNwzYMgV+HgbmVPANtSxsGFzLmiWLyL/IcQB66623Mu7QjBgxgg4dOtC0aVOKFi3Kt99+e0fFxMZadkcuUqRItn02bdrE4MGDM73Wpk0bFi9enOm1NWvW4O/vj6+vLy1btuS9996jaNEsVmAFkpOTSU5OzmjHxcXd5icQkfzoxKUE+s3ewf5zlr8NzzYtzattKuHmksVN9BvXYGk/OPi9pV35Ieg0ETx8rFewiPynHAegFi1akJaWBkC5cuU4dOgQV65cwdfXN+unHm6R2Wxm0KBBNG7cmGrVsl/+PTo6moCAzLsnBwQEEB0dndF+8MEH6dKlC6VLl+bYsWMMHTqUtm3bsmnTJpydb75NHR4ezrvvvnvbtYtI/rV451neXLSXhJR0ihR0Y8yjNbmvkn/Wnc/ugHl94NopcHKFNh/Avc9qorNIHnTLAejixYv07t2bX375BbPZTL169fjmm28oV67cv96xuVVhYWHs27eP9evX3/G5evTokfHP1atXp0aNGpQtW5Y1a9Zw//3339R/yJAhme4qxcXFERIScsd1iIj9SkxJY/iS/czbbtnBvUGZInzavTaBPh43dzYM2DoVfnrTMuRVuJRlrk/xOtYtWkRu2S0HoNdff51du3YxYsQIPDw8+Pzzz3n22WdZvXr1HRfRr18/li1bxrp16yhRosS/9g0MDCQmJibTazExMf8636dMmTL4+flx9OjRLAOQu7s77u7ut1e8iOQ7kdHxhM3ewdEL13EywYD7y9O/ZXmcnbK4k5MUC0v6wcGllnbljvDQRPAsbNWaRSRnbjkArVy5koiICNq0aQNAhw4dqFy5MsnJybcdHgzDoH///ixatIg1a9ZQunTp/zymYcOGrFq1ikGDBmWqrWHDhtkec+bMGS5fvkxQUNBt1SkijsEwDOb+fpp3lu4nOc2Mv5c743pks4M7WIa85j8FV0/+MeT1Ptz7nIa8ROzALQegc+fOUbPmX+tWlC9fHnd3d86fP09oaOht/fKwsDBmz57NkiVL8PLyypjH4+Pjg6enJwC9e/emePHihIeHA5aFGJs3b86YMWNo3749c+fOZdu2bUydOhWA69ev8+6779K1a1cCAwM5duwYr732GuXKlcsIbyIi/xSflMrQRfv4frdlaY/mFYoxpltN/Apl8X/wshzymgHF61q5ahG5XTmaBP3PCcTOzs4ZW2HcjsmTJwOWidV/N2PGDPr06QNAVFQUTk5/PWnRqFEjZs+ezVtvvcXQoUMpX748ixcvzpg47ezszJ49e5g5cybXrl0jODiY1q1bM3LkSA1ziUiW9p6Jpd+cHZy6nIiLk4lX21Tk2aZlcLqVIa9KHaDTJA15idgZk3GLCcbJyQkfH59MT3pdu3YNb2/vTAHlypUrd79KK4uLi8PHx4fY2Fi8vbVgmUh+ZRgGERtP8sHyg6SmGxQv7Mn4nrWpW8o36wPO7bQ85fXnkFfr96D+8xryEskjcvL9fct3gGbMmHHHhYmI5BXXElN4df4eVh6wPFTRpmoAo7vWxKdAFju4GwZs/QJ+fhPSU6BwyT+e8tKQl4i9uuUA9OSTT+ZmHSIiVrP91BX6z97Judgk3JydeLN9ZXo3LJX1WmYa8hLJl3K8EKKIiL0ymw0mrz3G2JWHSTcbhBYtwMTH6lCteDarNN805DUS6r+gIS+RfOCWA9Ct7vN1/Pjx2y5GRCS3XIxPZvB3u/jtyCUAOtUK5v2Hq1PIPYs/g1kNeT0SASU05CWSX9xyADp58iSlSpXisccew98/m2XgRUTyoI1HLzHw211cjE/Gw9WJEQ9V49F7SuRgyGsieGYzMVpE7NItB6Bvv/2W6dOnM3bsWNq2bcvTTz9Nu3btMj0BJiKSl6Slmxm/6ggTVh/FMKC8fyEmPV6HCgFeWR+gIS8Rh3HLj8H/6ezZs0RERBAREUFiYiK9evWib9++lC9fPrdqtDo9Bi9i/6JjkxgwdydbT1iW5uhRL4ThHavi6XbzhsgYBvz+Jfw01DLk5fPHU14a8hKxKzn5/s5xAPq7tWvX8s4777Bu3TouXbqEr2/+uEWsACRi31YfusDg73ZxNTGVgm7OfNClOp1qFc+68z+HvCq2h86TNOQlYodyZR2gv0tKSmL+/PlMnz6dLVu28Oijj1KgQIHbKlZE5G5JTTfz8U+RfL7O8jBGlSBvJj1eh9J+BbM+4J97eT3wLjR4UUNeIg4gRwFoy5YtTJs2je+++44yZcrw9NNPs2DBgnxz50dE7NfZazfoP3sHO6KuAdC7YSmGtquMh2s2Q16Z9vLSU14ijuaWA1DVqlW5cOECjz32GGvXrs20MaqIiC2tPBDDK/N2E3sjFS8PF0Z3rUHb6kFZd75xDZb2g4PfW9p6ykvEIeVoL7CCBQvi4uKS9aOjf9BeYCJiLSlpZkb9eIjpG04AULOEDxN61qFk0WyG5M9uh3lPwbVT2stLJB/SXmAiku9FXU6k35wd7DkTC0DfJqV5/cFKuLlksTSHYcCWKfDzsD+GvErBozO0l5eIA9NeYCJid37ce57X5u8hPjkNH09XPn60Jg9UCci6842rlqe8Di2ztCt3hIcmai8vEQenvcBExG4kpabzwfKDfLXpFAB1ShZmfM/alPDNZsjrzHaY3weuRYGzG7R+H+59VkNeIqIAJCL24cSlBPrN3sH+c3EAPN+8DK+0roirczZDXpsnw8q3LUNevqGWhQ2Da1u1ZhHJuxSARCTPW7LrLEMX7iUhJZ0iBd0Y060m91XMZk/CG1dhcRhE/mBpV+kED00Aj2x2fBcRh6QAJCJ5VlJqOu9+v585W08DcG9oEcb3rE2gj0fWB5zZZnnKK/aPIa82H0C9ZzTkJSI3ueUA1LRpUzp16sRDDz1EhQoVcrMmERGOXognbNZOImPiMZmg333lGHh/eVyyG/LaNAl+GQ7mNPAt/ceQVy1rly0iduKWt3J/9tln2bRpE3Xr1qVy5cq8/vrrbNiwgTvYSkxEJEsLtp+h44QNRMbE41fIna+frs/LrStmHX4Sr8Dcx+DnNy3hp+rD8Pw6hR8R+Vc53gw1OTmZVatWsWTJEr7//nvS09Np3749Dz30EG3atMHT0zO3arUaLYQoYhuJKWm8vWQ/87efAaBR2aJ82qMW/l7ZDHmd3grzn4bY0+DsDg+Gwz1Pa8hLxEFZbTd4sOwPtnTpUpYuXcqxY8do2bIlQ4YMoXHjxndyWptSABKxviMx8bw4awdHLlzHyQQD769Av5blcHbKIsyYzbBpAqwaYbnrU6SMZcgrSFv0iDgyqwagvzt27BhLly4lJCSERx555G6d1uoUgESsa96207y9ZD83UtMp5uXO+B61aVi2aNadEy7D4hfgyM+WdrWu0OFT8NB/qyKOzmYBKL9QABKxjsSUNIYt3s+CHZYhr6bl/RjbrRbFvNyzPuDUJsuQV/w5y5BX2w+hbh8NeYkIkEt7gYmI3E2H/xjyOvrHkNfgByrwYotyOGU35LXhE/j1fTDSoWg5eHQmBFazfuEiki8oAImIVRmGwbztZ3h7yT6SUs34e7kzvmdtGpTJZsjr+kVY9DwcW2Vp1+gO7ceCeyHrFS0i+Y4CkIhYTUJyGsMW72PhzrOAZcjrk+618CuUzZDXyfUwvy9cjwYXT2j3EdR+QkNeInLH7igA/Tl9yKQ/RiLyHyKj43lx1naOXUzAyQQvt67I/5qXzWbIKx1+GwNrwsEwg19Fy1NeAVWsXreI5E+3vBDi302bNo1q1arh4eGBh4cH1apV48svv7zbtYlIPmAYBt/+HkWnSes5djGBAG935jzbgLD7spnvc/0CfP0wrH7fEn5qPQ7PrVb4EZG7Ksd3gN5++23Gjh1L//79adiwIQCbNm3ipZdeIioqihEjRtz1IkXEPiUkp/HW4n0s+mPIq1mFYnzSrSZFsxvyOr4WFjwDCRfAtYBlrk+tnlasWEQcRY4fgy9WrBjjx4+nZ8/Mf5TmzJlD//79uXTp0l0t0Bb0GLzInTsUHceLs3Zw/GICzk4mXm5dgRea/cuQ19oPYe1owAD/KpYhr2IVrV22iNixXH0MPjU1lXvuueem1+vWrUtaWlpOTyci+YxlyOs0w5fuJznNTKC3BxMeq0290CJZHxAfbbnrc/I3S7tOb3jwQ3ArYL2iRcTh5HgOUK9evZg8efJNr0+dOpXHH3/8rhQlIvbpenIag77dxRsL95KcZqZFxWIsH9g0+/Bz7FeY0sQSftwKQZcv4aEJCj8ikutu6ymwadOm8fPPP9OgQQPAsh9YVFQUvXv3ZvDgwRn9xo4de3eqFJE87+D5OMJm7eD4JcuQ1yutK/J8szJZD3mlp1me8PptDGBAQHXLkJdfOWuXLSIOKscBaN++fdSpUwew7P0F4Ofnh5+fH/v27cvop0fjRRzDP4e8gnw8mNCzNvdkd9cn9qxlyCtqo6V9z9PQ5gNw9bRe0SLi8HIcgFavXp0bdYiIHfrnU14tKhZjbLdaFCnolvUBh3+2rOp84wq4ecFD4yybmYqIWJlWghaR2/L3hQ3/e8grFVaNgI3jLe2gmvDIDCha1rpFi4j8QQFIRHLsu22nM/by+s+nvK5FWXZwP/O7pX3v89B6JLhksxaQiIgVKACJyC1LTElj2OL9LNhxBriFhQ0PLoMlL0JSLHj4QKdJULmjFSsWEcnabW2FcbeEh4dTr149vLy88Pf3p3PnzkRGRv7ncfPmzaNSpUp4eHhQvXp1li9fnul9wzB4++23CQoKwtPTk1atWnHkyJHc+hgiDuFITDydJm5gwY4zOJng1TYViehTL+vwk5YMP74O3z5uCT/F68Lzvyn8iEieYdMAtHbtWsLCwti8eTMrV64kNTWV1q1bk5CQkO0xGzdupGfPnvTt25edO3fSuXNnOnfunOkJtNGjRzN+/HimTJnCli1bKFiwIG3atCEpKckaH0sk31mw/QwPTdzAkQvX8fdyZ/a/7eV15ThMaw1bpljajfrDUyvAt5R1ixYR+Rc53gojN128eBF/f3/Wrl1Ls2bNsuzTvXt3EhISWLZsWcZrDRo0oFatWkyZMgXDMAgODubll1/mlVdeASA2NpaAgAAiIiLo0aPHTedMTk4mOTk5ox0XF0dISIi2whCHdyMlneFL9/HdNsuQV9PyfnzSvRZ+2Q157VsISwdASjx4FoGHp0CFNlasWEQcWU62wrDpHaB/io2NBaBIkWwmU2LZeLVVq1aZXmvTpg2bNm0C4MSJE0RHR2fq4+PjQ/369TP6/FN4eDg+Pj4ZPyEhIXf6UUTs3tEL1+k8aQPfbbMMeQ1+oAIRT92bdfhJvQHLXoL5T1nCT8mG8MJ6hR8RybPyzCRos9nMoEGDaNy4MdWqVcu2X3R0NAEBAZleCwgIIDo6OuP9P1/Lrs8/DRkyJNMK1n/eARJxVIt3nmXoor0kpqTjV8id8T1r0aisX9adLx2BeX0gZh9ggqaDocVQcM4zf15ERG6SZ/5ChYWFsW/fPtavX2/13+3u7o67ux7JFUlKTefd7/czZ+tpABqVLcqnPWrh7+WR9QG7v7Xc+UlNgAJ+0GUqlLvfihWLiNyePBGA+vXrx7Jly1i3bh0lSpT4176BgYHExMRkei0mJobAwMCM9/98LSgoKFOfWrVq3d3CRfKR4xev8+KsHRyKjsdkggEtyzPg/vI4ZzXROSUBlr8Gu76xtEObQtcvwSvQukWLiNwmm84BMgyDfv36sWjRIn799VdKly79n8c0bNiQVatWZXpt5cqVNGzYEIDSpUsTGBiYqU9cXBxbtmzJ6CMimS3dfY6OE9ZzKDoev0JufNO3Pi89UCHr8HPhIHzR0hJ+TE6W4a7eSxR+RMSu2PQOUFhYGLNnz2bJkiV4eXllzNHx8fHB09OyMWLv3r0pXrw44eHhAAwcOJDmzZszZswY2rdvz9y5c9m2bRtTp04FLJuwDho0iPfee4/y5ctTunRphg0bRnBwMJ07d7bJ5xTJq5JS0xm57ACztkQB0KBMEcb3qI2/dxZDXoYBO7+23PlJuwGFAi13fUo3tXLVIiJ3zqYBaPLkyQC0aNEi0+szZsygT58+AERFReHk9NeNqkaNGjF79mzeeusthg4dSvny5Vm8eHGmidOvvfYaCQkJPPfcc1y7do0mTZqwYsUKPDyymccg4oBOXU7gxVk72H8uDpMJ+t1XjkGtsrnrkxxvmeuzd56lXbYlPDwVChWzbtEiIndJnloHKK/IyToCIvZoxb7zvDpvD/HJaRQp6MYn3WvRvEI2Yeb8bstTXleOg8kZ7h8GjQaCU55aRUNEJEff33liErSIWEdKmpnwHw8yY8NJAO4p5cuEx2oT5ON5c2fDgK1fwM9vQnoKeJeAR6ZDyfrWLVpEJBcoAIk4iDNXEwmbvZPdp68B8HzzMrzSuiKuzlncyblxFZb0g0N/rLhesZ1lI9MC2S9SKiJiTxSARBzAqoMxDP5uN7E3UvHxdGXMozVpVSUg686nf4f5T0NsFDi5QuuRUP8FMGUxN0hExE4pAInkY6npZj7+OZLP1x4HoGYJHyY+VoeQIgVu7mw2w6YJsGoEmNPANxQemQHF61i3aBERK1AAEsmnomOT6D9nB7+fvApAn0ahDG1XGTeXLIa8Ei7Bohfg6EpLu+rD0HEcePhYsWIREetRABLJh9Ydvsigb3dxJSGFQu4ujH6kBu2qB2Xd+eR6WPAMxJ8HFw94cBTU7aMhLxHJ1xSARPKRdLPBuFVHmPDrEQwDqgR589njdQj1K3hzZ3M6/DYG1oSDYQa/CvBoBARUtXrdIiLWpgAkkk9ciE9i0NxdbDx2GYDH6pfk7Q5V8HB1vrlzfDQsfBZOrLO0az0O7T4CtyyCkohIPqQAJJIPbDp2mQFzd3IxPpkCbs588HB1OtcunnXno6tg4XOQeAlcC0KHsVCzh3ULFhGxMQUgETtmNhtMXnuMMT9HYjagQkAhPnu8LuX8C93cOT0NVr8P68da2gHVLENefuWtWrOISF6gACRip64kpPDSt7tYe/giAF3rlGBk56oUcMviP+trp2FBXzi9xdK+py+0eR9cs1gBWkTEASgAidih7aeu0n/2Ds7FJuHu4sTITtXoVi8k686HfoDFL0LSNXD3hofGWx5zFxFxYApAInbEMAymbzhJ+PKDpJkNSvsV5LPH61A5KItN/9KS4edhsPVzSzu4jmUvryKlrVu0iEgepAAkYifiklJ5ff4eftwXDUD76kGM6lodLw/XmztfPgbzn7Ls5A7QsB/cPxxc3KxYsYhI3qUAJGIH9p+LJWzWDk5eTsTV2cRb7avQu2EpTFktVrhnHiwbBCnXwbMIPDwFKrSxes0iInmZApBIHmYYBnN/P83wpftJSTNTvLAnkx6vQ62Qwjd3TkmAH1+Dnd9Y2qUaQ5cvwCebx+FFRByYApBIHpWYksZbi/axcOdZAFpW8mdst5oULpDFMFbMAcuQ18VDgAmavwbNXgNn/ScuIpIV/XUUyYOOXojnxVk7OBxzHScTvNKmIi80K4uT0z+GvAwDdsyEH1+HtCQoFAhdv4DSzWxTuIiInVAAEsljluw6y5CFe0lMScffy53xPWvToEzRmzsmxcH3A2H/Qku77P3w8OdQqJh1CxYRsUMKQCJ5RHJaOiOXHeCbzVEANCpblHE9alPMy/3mzmd3WIa8rp4EJxe4/21o2B+cnKxbtIiInVIAEskDoi4nEjZ7B3vPxgLQv2U5BrWqgHNWQ16bJ8PKt8GcCj4lLWv7hNSzQdUiIvZLAUjExlYeiGHwd7uIT0rDt4Arn3SvRYuK/jd3TLxiWdH58I+WduWO8NAE8PS1bsEiIvmAApCIjaSlm/no50g+X3scgDolCzPxsToEF85if66TG2DBMxB/DpzdoM0HUO8ZyGodIBER+U8KQCI2cCEuiX6zd7L15BUAnm5cmjfaVsLN5R9zeMzpsO5jWDsKDDMULQePzICgGjaoWkQk/1AAErGyjccuMWDOLi5dT6aQuwsfPVKDttWDbu4Ydx4WPgsnf7O0a/aEdh+DeyHrFiwikg8pAIlYidlsMHntMcb8HInZgEqBXkx+oi6l/Qre3PnISlj0PCReBteC0H4M1Opp/aJFRPIpBSARK7iWmMJL3+5ideRFAB6tW4IRnarh6eacuWNaCvw6AjZOsLQDqsOjM8CvvJUrFhHJ3xSARHLZ7tPXeHHWDs5eu4G7ixMjO1WjW72QmzteOQEL+sLZ7Zb2vc/BAyPB1cO6BYuIOAAFIJFcYhgG32w+xchlB0lJN1OqaAE+e7wOVYN9bu68b6FlVefkOPAoDJ0mQeUOVq9ZRMRRKACJ5IKE5DTeWLiX73efA+DBqoGMfrQG3h6umTumJMJPQ2B7hKUdUh+6ToPCWdwhEhGRu0YBSOQuOxITzwvfbOfYxQRcnEy80bYSfZuUxvTPNXsuHIR5T8HFg4AJmg6GFkO1g7uIiBXoL63IXbR4p2Uj0xup6QR4uzPpsTrcE1okcyfDgB1f/bGD+w0o6A9dpkLZ+2xTtIiIA1IAErkLktPSGfH9AWZtsWxk2qScH5/2qIVfoX9sZJoUB8sGwb4FlnbZln/s4J7F1hciIpJrFIBE7tDpK4m8OMuykanJBP1blmfg/eVv3sj07zu4m5zh/mHQaKB2cBcRsQEFIJE78OuhGF76djexN1IpXMCVT7PayNRshs2fwS/v/G0H92kQcq9NahYREQUgkduSbjb4ZOVhJq4+CkDNkMJ89ngdiv9zI9OES7DoBTi60tKu/NAfO7gXtm7BIiKSiQKQSA5dup7MgDk72XjsMgC9G5bizfaVcXf5x6rOJ9bBgmfhejQ4u8OD4XDP09rBXUQkD1AAEsmBbSevEDZ7BzFxyRRwcya8S3U61SqeuVN6Gqz9ENZ9BBjgV9GynUVAVZvULCIiN7Pp7Mt169bRsWNHgoODMZlMLF68+D+PmTRpEpUrV8bT05OKFSvy1VdfZXo/IiICk8mU6cfDQ1sJyJ0xDINp60/QY+pmYuKSKVusIEvCGt8cfmLPwMwOsG40YEDtXvDcaoUfEZE8xqZ3gBISEqhZsyZPP/00Xbp0+c/+kydPZsiQIXzxxRfUq1ePrVu38uyzz+Lr60vHjh0z+nl7exMZGZnRvmkBOpEciE9K5fUFe1i+NxqADjWCGNW1BoXc//Gfz6EfYPGLkHQN3Lyg46dQ/RGr1ysiIv/NpgGobdu2tG3b9pb7f/311zz//PN0794dgDJlyvD777/z4YcfZgpAJpOJwMDAu16vOJ5D0XG8+M0Ojl9KwNXZxFvtq9C7YanMoTo1CVYOg61TLe3gOpanvIqUsU3RIiLyn+xqDlBycvJNw1menp5s3bqV1NRUXF0t+yxdv36dUqVKYTabqVOnDh988AFVq2Y/BJGcnExycnJGOy4uLnc+gNiVhTvOMHTRXpJSzQT7eDDx8TrUKembudOlI5a1faL3WtqN+kPLt8HFzfoFi4jILbOrFdjatGnDl19+yfbt2zEMg23btvHll1+SmprKpUuXAKhYsSLTp09nyZIlfPPNN5jNZho1asSZM2eyPW94eDg+Pj4ZPyEh2ojSkSWlpjN00V4Gf7ebpFQzTcv7sWxA05vDz6458HlzS/gpUBQenw+t31P4ERGxAybDMAxbFwGWYatFixbRuXPnbPvcuHGDsLAwvv76awzDICAggCeeeILRo0cTHR1NQEDATcekpqZSuXJlevbsyciRI7M8b1Z3gEJCQoiNjcXb2/uOP5vYj3+u6jygZXkG/HNV5+R4+OEV2DPX0i7dDB6eCt5BtilaREQAy/e3j4/PLX1/29UdIE9PT6ZPn05iYiInT54kKiqK0NBQvLy8KFasWJbHuLq6Urt2bY4ePZrted3d3fH29s70I45n9aELdJiwnr1nY/Et4ErEU/fy0gMVMoefc7ssd332zAWTE7R8C3otVvgREbEzdjUH6E+urq6UKFECgLlz59KhQwecstlPKT09nb1799KuXTtrlih2JN1s8Okvh5nw67+s6mwYsHkyrHzbsp2Fdwno+iWUamijqkVE5E7YNABdv349052ZEydOsGvXLooUKULJkiUZMmQIZ8+ezVjr5/Dhw2zdupX69etz9epVxo4dy759+5g5c2bGOUaMGEGDBg0oV64c165d46OPPuLUqVM888wzVv98kvddvp7MwLm7WH/UMofsyYaleLN9Fdxc/haoEy7Dkhfh8ApLu1IHy3YWBYrYoGIREbkbbBqAtm3bxn333ZfRHjx4MABPPvkkERERnD9/nqioqIz309PTGTNmDJGRkbi6unLfffexceNGQkNDM/pcvXqVZ599lujoaHx9falbty4bN26kSpUqVvtcYh92RF0lbNYOzscm4enqzKiuWazqfOI3WPgsxJ+3bGfR5n2o94y2sxARsXN5ZhJ0XpKTSVRifwzD4OvNpxi57ACp6QZl/AoypVddKgR4/dXppu0sKsAjMyCwms3qFhGRf5eT72+7nAMkcrsSktMYsnAvS3efA6Bd9UA+7FoDLw/XvzpdO2256xO1ydKu3QvafghuBW1QsYiI5AYFIHEYRy9c53/fbOfIheu4OJkY0q4yTzcOzbyq88HvYUk/bWchIpLPKQCJQ/hhz3lem7+bhJR0/L3cmfR4HeqF/m0Sc2oS/Pwm/P6lpV28LnSdBkVK26ZgERHJVQpAkq+lppsJX36I6RtOANCgTBHG96yNv9fftlS5GAnzn4aYfZZ2owHQcphWdBYRyccUgCTfiolLImzWDradugrAC83L8krrCrg4//GIu2HAzq/hx9chNREKFoOHp0C5VjasWkRErEEBSPKlTccu03/ODi5dT8HL3YWPu9WkTdXAvzokxcL3g2D/Qku7TAvLdhZeN2+nIiIi+Y8CkOQrhmEwdd1xRv8USbrZoFKgF1OeqEuo39+e4DqzzTLkde0UOLlYtrNoNBCyWU1cRETyHwUgyTfiklJ5dd5uftofA0CX2sV5/+HqeLo5WzqYzbDhU1j9PpjToHApeGQ6lLjHdkWLiIhNKABJvhAZHc8L32znxKUE3JydeLtjFR6vX/KvR9zjo2HR83B8jaVdrSt0+AQ8fGxWs4iI2I4CkNi9xTvPMmThXm6kphPs48FnT9SlVkjhvzocWQmLXoDES+BaANqOhtpPaDsLEREHpgAkdislzcx7Pxzgq02nAGha3o9xPWpTpOAfj6+npcCqd2HTREs7oJplyKtYRRtVLCIieYUCkNil87E3eHHWDnZGXQNgQMtyDGxVAWenP+7qXD5mmeh8fpelfe/z8MAIcPXI8nwiIuJYFIDE7mw4eon+c3ZyJSEFbw8XPu1Ri5aV/vb4+u658MPLkHIdPH2h02dQqZ3tChYRkTxHAUjshtlsMHntMcb8HInZgKrB3kx5oi4hRQpYOiTHww+vwJ65lnapJtBlKvgUt13RIiKSJykAiV2IvZHKy9/t5peDlkfcu91TghGdquHh+scj7ud2Woa8rhwHkxO0GAJNXwYnZxtWLSIieZUCkOR5B8/H8b9vtnPyciJuLk6MeKgqPe4taXnTbIbNk+CXd8GcCt4loOsXUKqRbYsWEZE8TQFI8rRFO88wZOFeklLNFC/syZQn6lK9xB9r91y/AIv/B0d/sbQrd4SO46FAkexPKCIiggKQ5FH/fMS9WYVijOteC98/H3E/9issfB4SLoCLBzwYDnWf0to+IiJySxSAJM/510fc01Jg9XuwYZylc7HK8OgM8K9su4JFRMTuKABJnrLx2CX6z97J5T8ecf+key3ur/zHI+5XjsOCZ+Dsdkv7nr7Q5n1w9bRdwSIiYpcUgCRPMAyDz9cdZ/SKQ5gNqBzkzZQn6lCq6B+7uO+ZB8tegpR48CgMnSZa5vyIiIjcBgUgsbn4pFRenbeHFfujAehSpzjvd/5jF/fk67D8Vdg929K5ZEPo8gUUDrFhxSIiYu8UgMSmjsTE8/zX2zl+KQFXZxPDO1b9axf3czthfl+4csyytk+z16DZq+Csf21FROTO6JtEbGbZnnO8Nn8PiSnpBPl4MOnxOtQp6WtZ22fTJPjlnT/W9iluuesT2tjWJYuISD6hACRWl5ZuZtSPh/hy/QkAGpYpyoTHauNXyP3mtX0qdYCHJmhtHxERuasUgMSqLsYn02/2DracuALA883L8Grrirg4O8HRVbDoBa3tIyIiuU4BSKxm+6krvDhrBzFxyRRyd+GjR2rQtnqQZW2fn0fAxgmWjv5V4JHpWttHRERyjQKQ5DrDMPhq0yne++EAqekG5fwLMeWJupTzLwSXj8GCvpYJzwD1noXWI7W2j4iI5CoFIMlVN1LSGbpoL4t2ngWgfY0gRnetQUF3F9g9F354GVKug6cvdJoEldrbuGIREXEECkCSa05dTuD5r7dzKDoeZycTQ9pWom+T0piS42HhK7DnW0vHUk2gy1TwKW7bgkVExGEoAEmu+PVQDIPm7iIuKQ2/Qm5MfKwODcoUhTPbLENeV0+CyRlaDIGmg8HJ2dYli4iIA1EAkrvKbDYYt+oI41YdAaBOycJ89nhdAr3c4LexsPp9MKeBT0no+iWUrG/jikVExBEpAMldcy0xhUHf7mJN5EUAejcsxVvtq+CWGANfPw8n1lo6Vu0CHT4Bz8K2K1ZERByaApDcFfvPxfLCN9s5feUG7i5OhHepTpc6JSDyR1j8Ity4Aq4FoO1oqP2E1vYRERGbUgCSO7ZwxxmGLNxLcpqZkCKeTHmiLlWLuVs2Md061dIpsIZlbR+/8rYtVkREBAUguQMpaWbe/+EAMzedAqB5hWKM61GLwtePwxdPw4X9lo4NwqDVcHBxt2G1IiIif1EAktsSE5fEi7N2sP3UVQAGtCzHwPvL47xjBvw0FNKSoGAx6DwFyreycbUiIiKZKQBJjm09cYWw2Tu4GJ+Ml4cLn3SrRatQV5jXCw4ts3Qq29ISfrwCbFusiIhIFhSA5JYZhsHMjSd574eDpJkNKgZ48XmvuoTG74DJz0H8OXByhVbvQIMXwcnJ1iWLiIhkyabfUOvWraNjx44EBwdjMplYvHjxfx4zadIkKleujKenJxUrVuSrr766qc+8efOoVKkSHh4eVK9eneXLl+dC9Y7lRko6L327i3e+P0Ca2aBjzWAWvVCP0N1jYWZHS/gpWg6e+QUa9VP4ERGRPM2m31IJCQnUrFmTSZMm3VL/yZMnM2TIEN555x3279/Pu+++S1hYGN9//31Gn40bN9KzZ0/69u3Lzp076dy5M507d2bfvn259THyvajLiXSZvJHFu87h7GRiWIcqjG9TmAKzOsJvHwMG1HoCnlsLwbVsXa6IiMh/MhmGYdi6CACTycSiRYvo3Llztn0aNWpE48aN+eijjzJee/nll9myZQvr168HoHv37iQkJLBs2bKMPg0aNKBWrVpMmTIly/MmJyeTnJyc0Y6LiyMkJITY2Fi8vb3v8JPZt9WRFxg4Z2fmLS0SVsOylyA5Dtx9oOMnUK2rrUsVEREHFxcXh4+Pzy19f9vVOEVycjIeHh6ZXvP09GTr1q2kpqYCsGnTJlq1yvzUUZs2bdi0aVO25w0PD8fHxyfjJyQk5O4Xb2fMZoPxq47wdMTvxCWlUbtkYZY9X4sGu9+y7OWVHAch9eGF3xR+RETE7thVAGrTpg1ffvkl27dvxzAMtm3bxpdffklqaiqXLl0CIDo6moCAzE8eBQQEEB0dne15hwwZQmxsbMbP6dOnc/Vz5HVxSak89/V2xq48jGHA4/VL8m0HNwLntIbds8HkBM1fhz7LwbeUrcsVERHJMbt6CmzYsGFER0fToEEDDMMgICCAJ598ktGjR+N0B5Nu3d3dcXfXIn0Ah2Pief7r7Zy4lICbixPvPVSFbimLIGKkZRNT7xLQZSqENrZ1qSIiIrfNru4AeXp6Mn36dBITEzl58iRRUVGEhobi5eVFsWLFAAgMDCQmJibTcTExMQQGBtqiZLvyw57zdJ60gROXEgj28WBxrzJ0O9gffhluCT9VOsH/1iv8iIiI3bOrAPQnV1dXSpQogbOzM3PnzqVDhw4Zd4AaNmzIqlWrMvVfuXIlDRs2tEWpdiEt3Uz4jwcJm72DxJR0GpUtyop2CVRZ/KBlB3fXAvDQBHh0Jnj62rpcERGRO2bTIbDr169z9OjRjPaJEyfYtWsXRYoUoWTJkgwZMoSzZ89mrPVz+PBhtm7dSv369bl69Spjx45l3759zJw5M+McAwcOpHnz5owZM4b27dszd+5ctm3bxtSpU63++ezBlYQU+s/ZwYajlwEIaxLMy3yD06IvLR0Ca0DXaVCsgg2rFBERubtsGoC2bdvGfffdl9EePHgwAE8++SQRERGcP3+eqKiojPfT09MZM2YMkZGRuLq6ct9997Fx40ZCQ0Mz+jRq1IjZs2fz1ltvMXToUMqXL8/ixYupVq2a1T6Xvdh7JpYXvtnO2Ws3KODmzJQHPGi2539w8aClQ8N+cP/b2sRURETynTyzDlBekpN1BOzV/O1nGLpoLylpZkKLePJtnX0EbHoP0pOhoD88PBnKaRNTERGxHzn5/rarp8DkzqWkmXnvhwN8tekUAJ3Ku/Gx22Rc1/9s6VC+NXT6DAoVs2GVIiIiuUsByIFciEvixVk72HbqKgBj617h4VMjMV2PAWd3aD0S7n0OTCYbVyoiIpK7FIAcxPZTV/jfNzu4EJ+MrwcsrriKUvunWd4sVsky0TlQ86RERMQxKADlc4Zh8M2WKEZ8v5/UdIOWfteY7DkZ98i9lg73PA2t3we3ArYtVERExIoUgPKxpNR03l6yj++2nQEM3i+5k8eufobpeiJ4FoFOE6FSe1uXKSIiYnUKQPnUuWs3+N8329l9JpbCpussKPEdZS/8YnmzdHN4+HPwDrJtkSIiIjaiAJQPbTp2mX6zd3A5IYX7PQ/zmecU3C9Gg5OLZV2fhv3hDvZOExERsXcKQPmIYRhM33CSD5YfxGROZVThH+ie9B2mRAOKlIWuX0LxOrYuU0RExOYUgPKJGynpvLFwD0t2naOkKYavCn9BaNIBy5u1n4AHPwT3QrYtUkREJI9QAMoHTl9J5Pmvt3PgfCxdXTYwyj0C16RE8PCBjuOg6sO2LlFERCRPUQCyc78duUj/OTtJT4xlimcEDxrrIR0o2Qi6TIXCIbYuUUREJM9RALJThmEwdd1xPlxxiNpE8lmByQSYL4DJGVoMgaaDwcnZ1mWKiIjkSQpAdigxJY3X5u/hxz1nGOCyiAEui3Eym8E3FLp8CSH1bF2iiIhInqYAZGdOXU7g+a+3cz3mGN+5fUZdp8OWN2r0gHYfgUf+3L1eRETkblIAsiNrD19kwJydtEhezfvuERQiEdy9of1YqPGorcsTERGxGwpAdsAwDCavPcaUn3byjksEXdzWW94IaWCZ6OxbyrYFioiI2BkFoDwuITmNV+fvJnrfOn5wnUSI00UMkxOm5q9D01fAWf8TioiI5JS+PfOwk5cS+N9XW3ng8jeMd1uIi8kMhUti6vIllKxv6/JERETslgJQHrUm8gKj5vzESPN46rn+MdG5ejdo/7FlgUMRERG5bQpAeYxhGHy25hiRv8zgO5dpeDvdwOxWCKcOn0CNbrYuT0REJF9QAMpDEpLTGPbtRpoc+ZAwV8tEZ3OJejh1/dKyxo+IiIjcFQpAecTJSwl8MuMbXo7/iJLOFzHjhFOL13HSRGcREZG7Tt+secDag+fY9+1wxhjzcXEyk1yoBO7dpkHJBrYuTUREJF9SALIhwzD4ZsU6Km96hTCnw2CCG5W64tn5E010FhERyUUKQDaSkJzGvOlj6BL9Cd5ON0hyKojzQ5/gWau7rUsTERHJ9xSAbOD0uWiOzHiOPqlrwQQXCtfC/8mvtKKziIiIlSgAWdnuDSvwW9mPllwkDSdiag2geMdhmugsIiJiRfrWtaLfvxtNnf0f4GwyiHYKxPXRLyleuamtyxIREXE4TrYuwJH4lG9IOk5s92mD78ubKarwIyIiYhO6A2RFFWo35VSBX6lToSYmk8nW5YiIiDgsBSArK1Wxlq1LEBERcXgaAhMRERGHowAkIiIiDkcBSERERByOApCIiIg4HAUgERERcTgKQCIiIuJwFIBERETE4dg0AK1bt46OHTsSHByMyWRi8eLF/3nMrFmzqFmzJgUKFCAoKIinn36ay5cvZ7wfERGByWTK9OPh4ZGLn0JERETsjU0DUEJCAjVr1mTSpEm31H/Dhg307t2bvn37sn//fubNm8fWrVt59tlnM/Xz9vbm/PnzGT+nTp3KjfJFRETETtl0Jei2bdvStm3bW+6/adMmQkNDGTBgAAClS5fm+eef58MPP8zUz2QyERgYeFdrFRERkfzDruYANWzYkNOnT7N8+XIMwyAmJob58+fTrl27TP2uX79OqVKlCAkJoVOnTuzfv/9fz5ucnExcXFymHxEREcm/7CoANW7cmFmzZtG9e3fc3NwIDAzEx8cn0xBaxYoVmT59OkuWLOGbb77BbDbTqFEjzpw5k+15w8PD8fHxyfgJCQmxxscRERERGzEZhmHYugiwDFstWrSIzp07Z9vnwIEDtGrVipdeeok2bdpw/vx5Xn31VerVq8e0adOyPCY1NZXKlSvTs2dPRo4cmWWf5ORkkpOTM9pxcXGEhIQQGxuLt7f3HX0uERERsY64uDh8fHxu6fvbrnaDDw8Pp3Hjxrz66qsA1KhRg4IFC9K0aVPee+89goKCbjrG1dWV2rVrc/To0WzP6+7ujru7e0b7z0yooTARERH78ef39q3c27GrAJSYmIiLS+aSnZ2dgew/bHp6Onv37r1pntC/iY+PB9BQmIiIiB2Kj4/Hx8fnX/vYNABdv349052ZEydOsGvXLooUKULJkiUZMmQIZ8+e5auvvgKgY8eOPPvss0yePDljCGzQoEHce++9BAcHAzBixAgaNGhAuXLluHbtGh999BGnTp3imWeeueW6goODOX36NF5eXphMprv6mf8cXjt9+rSG13KRrrN16Dpbh66zdeg6W09uXWvDMIiPj8/IBP/GpgFo27Zt3HfffRntwYMHA/Dkk08SERHB+fPniYqKyni/T58+xMfHM3HiRF5++WUKFy5My5YtMz0Gf/XqVZ599lmio6Px9fWlbt26bNy4kSpVqtxyXU5OTpQoUeIufMLseXt76z8wK9B1tg5dZ+vQdbYOXWfryY1r/V93fv6UZyZBO4qcTNCS26frbB26ztah62wdus7WkxeutV09Bi8iIiJyNygAWZm7uzvDhw/P9NSZ3H26ztah62wdus7WoetsPXnhWmsITERERByO7gCJiIiIw1EAEhEREYejACQiIiIORwFIREREHI4CUC6YNGkSoaGheHh4UL9+fbZu3fqv/efNm0elSpXw8PCgevXqLF++3EqV2recXOcvvviCpk2b4uvri6+vL61atfrP/13EIqf/Pv9p7ty5mEymf93gWP6S0+t87do1wsLCCAoKwt3dnQoVKuhvxy3I6XX+9NNPqVixIp6enoSEhPDSSy+RlJRkpWrt07p16+jYsSPBwcGYTCYWL178n8esWbOGOnXq4O7uTrly5YiIiMj1OjHkrpo7d67h5uZmTJ8+3di/f7/x7LPPGoULFzZiYmKy7L9hwwbD2dnZGD16tHHgwAHjrbfeMlxdXY29e/dauXL7ktPr/NhjjxmTJk0ydu7caRw8eNDo06eP4ePjY5w5c8bKlduXnF7nP504ccIoXry40bRpU6NTp07WKdaO5fQ6JycnG/fcc4/Rrl07Y/369caJEyeMNWvWGLt27bJy5fYlp9d51qxZhru7uzFr1izjxIkTxk8//WQEBQUZL730kpUrty/Lly833nzzTWPhwoUGYCxatOhf+x8/ftwoUKCAMXjwYOPAgQPGhAkTDGdnZ2PFihW5WqcC0F127733GmFhYRnt9PR0Izg42AgPD8+yf7du3Yz27dtneq1+/frG888/n6t12rucXud/SktLM7y8vIyZM2fmVon5wu1c57S0NKNRo0bGl19+aTz55JMKQLcgp9d58uTJRpkyZYyUlBRrlZgv5PQ6h4WFGS1btsz02uDBg43GjRvnap35ya0EoNdee82oWrVqpte6d+9utGnTJhcrMwwNgd1FKSkpbN++nVatWmW85uTkRKtWrdi0aVOWx2zatClTf4A2bdpk219u7zr/U2JiIqmpqRQpUiS3yrR7t3udR4wYgb+/P3379rVGmXbvdq7z0qVLadiwIWFhYQQEBFCtWjU++OAD0tPTrVW23bmd69yoUSO2b9+eMUx2/Phxli9fTrt27axSs6Ow1fegTTdDzW8uXbpEeno6AQEBmV4PCAjg0KFDWR4THR2dZf/o6Ohcq9Pe3c51/qfXX3+d4ODgm/6jk7/cznVev34906ZNY9euXVaoMH+4net8/Phxfv31Vx5//HGWL1/O0aNHefHFF0lNTWX48OHWKNvu3M51fuyxx7h06RJNmjTBMAzS0tJ44YUXGDp0qDVKdhjZfQ/GxcVx48YNPD09c+X36g6QOJxRo0Yxd+5cFi1ahIeHh63LyTfi4+Pp1asXX3zxBX5+frYuJ18zm834+/szdepU6tatS/fu3XnzzTeZMmWKrUvLV9asWcMHH3zAZ599xo4dO1i4cCE//PADI0eOtHVpchfoDtBd5Ofnh7OzMzExMZlej4mJITAwMMtjAgMDc9Rfbu86/+njjz9m1KhR/PLLL9SoUSM3y7R7Ob3Ox44d4+TJk3Ts2DHjNbPZDICLiwuRkZGULVs2d4u2Q7fz73NQUBCurq44OztnvFa5cmWio6NJSUnBzc0tV2u2R7dznYcNG0avXr145plnAKhevToJCQk899xzvPnmmzg56R7C3ZDd96C3t3eu3f0B3QG6q9zc3Khbty6rVq3KeM1sNrNq1SoaNmyY5TENGzbM1B9g5cqV2faX27vOAKNHj2bkyJGsWLGCe+65xxql2rWcXudKlSqxd+9edu3alfHz0EMPcd9997Fr1y5CQkKsWb7duJ1/nxs3bszRo0czAibA4cOHCQoKUvjJxu1c58TExJtCzp+h09A2mneNzb4Hc3WKtQOaO3eu4e7ubkRERBgHDhwwnnvuOaNw4cJGdHS0YRiG0atXL+ONN97I6L9hwwbDxcXF+Pjjj42DBw8aw4cP12PwtyCn13nUqFGGm5ubMX/+fOP8+fMZP/Hx8bb6CHYhp9f5n/QU2K3J6XWOiooyvLy8jH79+hmRkZHGsmXLDH9/f+O9996z1UewCzm9zsOHDze8vLyMOXPmGMePHzd+/vlno2zZska3bt1s9RHsQnx8vLFz505j586dBmCMHTvW2Llzp3Hq1CnDMAzjjTfeMHr16pXR/8/H4F999VXj4MGDxqRJk/QYvL2aMGGCUbJkScPNzc249957jc2bN2e817x5c+PJJ5/M1P+7774zKlSoYLi5uRlVq1Y1fvjhBytXbJ9ycp1LlSplADf9DB8+3PqF25mc/vv8dwpAty6n13njxo1G/fr1DXd3d6NMmTLG+++/b6SlpVm5avuTk+ucmppqvPPOO0bZsmUNDw8PIyQkxHjxxReNq1evWr9wO7J69eos/97+eW2ffPJJo3nz5jcdU6tWLcPNzc0oU6aMMWPGjFyv02QYuo8nIiIijkVzgERERMThKACJiIiIw1EAEhEREYejACQiIiIORwFIREREHI4CkIiIiDgcBSARERFxOApAIiIi4nAUgERERMThKACJiIiIw1EAEhG5RS1atGDQoEF3dA7DMHjuuecoUqQIJpOJXbt23ZXaRCRnFIBExC499dRTvPXWW7YuI8dWrFhBREQEy5Yt4/z581SrVs3WJYk4JBdbFyAiklPp6eksW7aMH374wdal5NixY8cICgqiUaNG2fZJSUnBzc3NilWJOB7dARKRDHPmzMHT05Pz589nvPbUU09Ro0YNYmNjb/u8JUqU4LPPPsv02saNGylQoACnTp3K8fk2btyIq6sr9erVy/L9Fi1a0L9/fwYNGoSvry8BAQF88cUXJCQk8NRTT+Hl5UW5cuX48ccfMx2XnJzMgAED8Pf3x8PDgyZNmvD7779nW4fZbCY8PJzSpUvj6elJzZo1mT9/frb9+/TpQ//+/YmKisJkMhEaGppRb79+/Rg0aBB+fn60adOGFStW0KRJEwoXLkzRokXp0KEDx44du+n3jx49mnLlyuHu7k7JkiV5//33b/Eqijg2BSARydCjRw8qVKjABx98AMDw4cP55Zdf+PHHH/Hx8bnt89avXz9TkDAMg0GDBvHSSy9RqlSpHJ9v6dKldOzYEZPJlG2fmTNn4ufnx9atW+nfvz//+9//ePTRR2nUqBE7duygdevW9OrVi8TExIxjXnvtNRYsWMDMmTPZsWMH5cqVo02bNly5ciXL3xEeHs5XX33FlClT2L9/Py+99BJPPPEEa9euzbL/uHHjGDFiBCVKlOD8+fOZrsnMmTNxc3Njw4YNTJkyhYSEBAYPHsy2bdtYtWoVTk5OPPzww5jN5oxjhgwZwqhRoxg2bBgHDhxg9uzZBAQE5PRyijgmQ0Tkb77//nvD3d3deO+99wxfX19j3759Ge917tzZKFy4sNG1a9ccnXP06NFG1apVM9ozZ840AgMDjfj4+Ns6b/ny5Y1ly5Zl+37z5s2NJk2aZLTT0tKMggULGr169cp47fz58wZgbNq0yTAMw7h+/brh6upqzJo1K6NPSkqKERwcbIwePTrjvAMHDjQMwzCSkpKMAgUKGBs3bsz0u/v27Wv07Nkz29o++eQTo1SpUjfVW7t27X/9zBcvXjQAY+/evYZhGEZcXJzh7u5ufPHFF/96nIhkTXeARCSTDh06UKVKFUaMGMGiRYuoWrVqxnsDBw7kq6++yvE5GzRowMGDB7l+/ToJCQkMHTqU9957j0KFCuX4vAcPHuTcuXPcf//9/9qvRo0aGf/s7OxM0aJFqV69esZrf94puXDhAmCZm5Oamkrjxo0z+ri6unLvvfdy8ODBm85/9OhREhMTeeCBByhUqFDGz1dffXXTUNWtqFu3bqb2kSNH6NmzJ2XKlMHb2ztjuCwqKgqwXIfk5OT/vA4ikjVNghaRTFasWMGhQ4dIT0+/aTilRYsWrFmzJsfnrFu3Lk5OTuzYsYNffvmFYsWK8dRTT93WeZcuXcoDDzyAh4fHv/ZzdXXN1DaZTJle+3P47O9DSjlx/fp1AH744QeKFy+e6T13d/ccn69gwYKZ2h07dqRUqVJ88cUXBAcHYzabqVatGikpKQB4enreVt0iYqE7QCKSYceOHXTr1o1p06Zx//33M2zYsLty3gIFClC9enUWLFjAxx9/zCeffIKT0+39+VmyZAmdOnW6K3X9XdmyZTPm4PwpNTWV33//nSpVqtzUv0qVKri7uxMVFUW5cuUy/YSEhNxRLZcvXyYyMpK33nqL+++/n8qVK3P16tVMfcqXL4+npyerVq26o98l4qh0B0hEADh58iTt27dn6NChGUMvDRs2ZMeOHdSpU+eOz9+gQQMmTJhAp06daNGixW2d48KFC2zbto2lS5fecT3/VLBgQf73v//x6quvUqRIEUqWLMno0aNJTEykb9++N/X38vLilVde4aWXXsJsNtOkSRNiY2PZsGED3t7ePPnkk7ddi6+vL0WLFmXq1KkEBQURFRXFG2+8kamPh4cHr7/+Oq+99hpubm40btyYixcvsn///ox6J06cyKJFixSSRLKgACQiXLlyhQcffJBOnTplfNHWr1+ftm3bMnToUFasWPGf54iIiOCpp57CMIws369Zsyaurq589NFHt13n999/z7333oufn99tn+PfjBo1CrPZTK9evYiPj+eee+7hp59+wtfXN8v+I0eOpFixYoSHh3P8+HEKFy5MnTp1GDp06B3V4eTkxNy5cxkwYADVqlWjYsWKjB8//qbgOGzYMFxcXHj77bc5d+4cQUFBvPDCCxnvX7p06bbmI4k4ApOR3V8rEZEsrFmzhokTJ9603s3w4cNZu3ZttnN57rvvPurUqcOYMWNydN6/e+ihh2jSpAmvvfbabdcvIgK6AyQiOdCqVSt2795NQkICJUqUYN68eTRs2BCAH3/8kYkTJ2bqbzabuXjxItOmTePIkSMsWbIkx+f9uyZNmtCzZ8+7/8FExOHoDpCI5Jo1a9bQsmVLKlWqxIwZM6hfv76tSxIRARSARERExAHpMXgRERFxOApAIiIi4nAUgERERMThKACJiIiIw1EAEhEREYejACQiIiIORwFIREREHI4CkIiIiDgcBSARERFxOApAIiIi4nD+D2RjvyA8QdRSAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "plt.plot(df['xL_0 / mole frac.'], df['pL / Pa']/1e6)\n", "plt.plot(df['xV_0 / mole frac.'], df['pL / Pa']/1e6)\n", "plt.gca().set(xlabel='$x_1,y_1$ / mole frac.', ylabel='p / MPa')\n", "plt.show()" ] }, { "cell_type": "markdown", "id": "a9c9fe55", "metadata": {}, "source": [ "Isn't that exciting!\n", "\n", "You can also provide an optional set of flags to the function to control other behaviors of the function, and switch between simple Euler and adaptive RK45 integration (the default)" ] }, { "cell_type": "raw", "id": "264c5123", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The options class is here: :py:meth:`~teqp.teqp.TVLEOptions`" ] }, { "cell_type": "markdown", "id": "d4193110", "metadata": {}, "source": [ "Supercritical isotherms work approximately in the same manner" ] }, { "cell_type": "code", "execution_count": 6, "id": "c5b925ad", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:44.885664Z", "iopub.status.busy": "2024-12-12T18:10:44.885258Z", "iopub.status.idle": "2024-12-12T18:10:45.007973Z", "shell.execute_reply": "2024-12-12T18:10:45.007414Z" } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjcAAAG0CAYAAADO5AZFAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjkuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/GU6VOAAAACXBIWXMAAA9hAAAPYQGoP6dpAABnVklEQVR4nO3dd3gU5QLF4d9uOqTQe+i9hQ6BoCgoAlKlSBdBLICg14ZdLNgLgiCKgCKg9A4CCkjvEHon9FDTIHXn/jEYRSCSkGSym/M+zz5mJ7O7J3OBnDvzzffZDMMwEBEREXERdqsDiIiIiKQnlRsRERFxKSo3IiIi4lJUbkRERMSlqNyIiIiIS1G5EREREZeiciMiIiIuxd3qAJnN4XBw+vRp/Pz8sNlsVscRERGRO2AYBlFRURQpUgS7PeVzM9mu3Jw+fZrAwECrY4iIiEganDhxgmLFiqW4T7YrN35+foB5cPz9/S1OIyIiInciMjKSwMDA5N/jKcl25eavS1H+/v4qNyIiIk7mToaUaECxiIiIuBSVGxEREXEpKjciIiLiUlRuRERExKWo3IiIiIhLUbkRERERl6JyIyIiIi5F5UZERERcisqNiIiIuBSVGxEREXEpKjciIiLiUlRuRERExKWo3GSQazFRVkcQERHJlrLdquAZ7tIRIn/qif/lXWx1lGVfq1k0LpePwDw5rE4mIiKSLajcpKeDy2BqN/yT4gCoZT9Eh1mhAFQrGkC3+sVpX7Mo3h5uVqYUERFxaboslV6SEmDeYLhebP5Sr1Qe7DYIPRXB0JmhhHz0ByN/P0jEtQSLgoqIiLg2lZv0Er4HIk+Cl/8Nm399MpiNrzXjtZaVKBLgzYXoOD797QBNPvmD8WuOEp/osCiwiIiIa1K5SS+Xj5v/zVf+pm/l8/XiiXtKs/Kl+/iySw3KFvDl8tUE3pm3hwe/WMnSPecyOayIiIjrUrlJL16+5n8Trt52Fw83O+1qFmXx4MZ80L4a+Xy9OHbxKk/8uJn+P27m9JVrmRRWRETEdancpJec+c3/Rp7+z13d3ex0q1+cFS824al7y+But/HbnnM88PlKvv/zCIlJulQlIiKSVio36SVfeXD3gdgrd/wSXy93XmlRkQXPNqZOidzExCfx3oK9dBi9loPnNE+OiIhIWqjcpBd3LyheP00vrVDIj1+fDObDDtXw93Zn58kIWn29mm9XHibJYaRzUBEREdemcpOeSjZO80vtdhuP1ivO0ufv5b4K+YlPdDB80T46f7uOoxdi0jGkiIiIa1O5SU+V2iR/ubPCYA53WJTqtyjo780Pj9Xlo0eq4evlzpbjl2nx1Sp+XHcMw9BZHBERkf+icpOe8peHYnUBqF6iAGWqN0zT29hsNrrULc7iIY0JLp2X2AQHb87ZTb+Jm7kYHfffbyAiIpKNqdyktxrdzP9u/xnu8kxLsdw5+Llffd5qXRlPdzvL94XT4qs/WXPoQjoEFRERcU0qN+mtSgdw9zZnLD615a7fzm630adRKWY/04gy+XMSHhVHj3Eb+HDRPhJ0y7iIiMhNVG7Sm08uqPqI+fWaL9PtbSsX8Wf+oMZ0rVccw4AxKw/Tccw6jl/UYGMREZF/UrnJCA2fNf+7dz6cP5Bub+vj6cbwDtUY3b0W/t7u7DhxhVYjVjNn+6l0+wwRERFnp3KTEQpUhAqtAAPWfpXub9+iWmEWDbmHeiXzEB2XyOCp23l9diixCUnp/lkiIiLORuUmo4Q8Z/53xy8Qkf5nVorm8mHyE/UZdH9ZbDaYtD6MjmPWEnbx9mtbiYiIZAcqNxklsC6UCAFHAqwblSEf4e5m538PVmD8Y3XJncODXaciafX1nyzZfTZDPk9ERMQZqNxkpL/O3mwelyFnb/7SpEIBFjzbmFrFcxEVm8iTP23h/QV7dDeViIhkSyo3GalsUygeDImxsGJ4hn5UkVw+/PJkMP1CSgHw3Z9HeXTses5EXMvQzxUREclqVG4yks0GDwwzv97+M4Tvy9CP83Cz8/rDlRnTozZ+3ubSDa1GrGbVgfMZ+rkiIiJZicpNRgusB5Vag+GAZW9nykc+VLUQ8weFUKWIP5di4uk9fiOj/jiktalERCRbULnJDE3fApsbHFgEx9dmykeWyJuTGU83pGu9QAwDPlmyn6cmbSEqNiFTPl9ERMQqKjeZIV85qNXL/Hrpm3e95tSd8vZwY3iH6gzvUA1PNztLdp+j3ag1HAqPzpTPFxERsYLKTWZp8gp45ICTm2Dv3Ez96K71ivPLkw0o5O/N4fMxtBu1RreLi4iIy1K5ySx+hSB4oPn1b69DQubexVSzeG7mDQqhfilzVuMnf9rCJ0v2keTQOBwREXEtKjeZqdFg8C8KV8Lgz88z/ePz+3kxqV99Hm9k3i4+6o/DPD5hE1euxmd6FhERkYyicpOZvHzhoevz3az5Ei4ezvQIHm523mxdma8erYG3h52VB87TeuRq9pyOzPQsIiIiGUHlJrNVagNl7oekeFj0UqYNLv63tjWKMvPpRgTm8eHEpWt0GL1Gq4uLiIhLsLTcvP3229hsthseFStWTPE106ZNo2LFinh7e1OtWjUWLlyYSWnTic0GLT8FN084tAz2zrMsSuUi/swbGMI95fMTm+Bg8NTtDJunZRtERMS5WX7mpkqVKpw5cyb5sXr16tvuu3btWrp27Urfvn3Ztm0b7dq1o127duzatSsTE6eDvGWg4bPm14uHQnyMZVFy5fBk/GN1GXhfWQB+WHOUHt9v4HxUnGWZRERE7obl5cbd3Z1ChQolP/Lly3fbfb/66iseeughXnzxRSpVqsS7775LrVq1GDlyZCYmTieN/wcBxSHyJKz6xNIobnYbLzSvwLc9a+Pr5c6Go5do/fVqtp+4YmkuERGRtLC83Bw8eJAiRYpQunRpunfvTlhY2G33XbduHc2aNbthW/PmzVm3bt1tXxMXF0dkZOQNjyzBMwe0+ND8eu1IOH/A2jxA8yqFmD2gEWXy5+RsZCydx6xj6sbb/+8hIiKSFVlaburXr8+ECRNYvHgxo0eP5ujRozRu3JioqKhb7n/27FkKFix4w7aCBQty9uztJ6QbPnw4AQEByY/AwMB0/RnuSoWWUK45OBJg4f8sG1z8T2UL+DJnYAgPVSlEfJKDV2aGMnRmKPGJGocjIiLOwdJy06JFCzp16kT16tVp3rw5Cxcu5MqVK/z666/p9hlDhw4lIiIi+XHixIl0e++7ZrNBi4/A3RuOroLtk61OBICvlzuje9TixeYVsNlgysYwun+/ngvRGocjIiJZn+WXpf4pV65clC9fnkOHDt3y+4UKFeLcuXM3bDt37hyFChW67Xt6eXnh7+9/wyNLyVPKXJoBYMmrEB1ubZ7rbDYbA+4ryw+96+Ln5c6mY5dpO3INu09HWB1NREQkRVmq3ERHR3P48GEKFy58y+8HBwezfPnyG7YtXbqU4ODgzIiXcYIHQaHqEHsFFr5odZob3FexALMGNKJ0vpycunKNR0avZcHOM1bHEhERuS1Ly80LL7zAypUrOXbsGGvXrqV9+/a4ubnRtWtXAHr16sXQoUOT9x88eDCLFy/ms88+Y9++fbz99tts3ryZgQMHWvUjpA83d2g7EmxusGc27J1vdaIblC3gy6wBjZLnwxkweSuf/bYfh9alEhGRLMjScnPy5Em6du1KhQoV6Ny5M3nz5mX9+vXkz58fgLCwMM6c+fssQcOGDZk8eTJjx44lKCiI6dOnM3v2bKpWrWrVj5B+CgdBo+tz3yz4H1y7Ymmcfwvw8WD8Y3Xpf09pAL7+/RBPTtpCdFyixclERERuZDOMLHCLTiaKjIwkICCAiIiIrDf+JuEajG4Elw5DjR7QbpTViW5p5taTvHL9DqryBX35vlddiufNYXUsERFxYan5/Z2lxtxkex4+0O4bwAbbJ8GBJVYnuqUOtYrx65PBFPDz4sC5aNqMWs3aQxesjiUiIgKo3GQ9xRtA8ADz67nPwrXL1ua5jRqBuZg3KISgwFxcuZpAzx82MnHtMbLZiUAREcmCVG6yovtfh7xlIfosLHrF6jS3VdDfm1/6N6BDzaIkOQzemrtbE/6JiIjlVG6yIg8faDcGbHbYORX2Zd2Vz7093PiscxCvtayE3QZTN52g23frtfCmiIhYRuUmqwqsCw0HmV/PexZisu6YFpvNxhP3lOaHx+ri5+3O5uOXaTtyNbtOacI/ERHJfCo3WVmTVyF/JYg5D/MGZ4m1p1LSpEIBZg9oROn8OTkdEUvHMWuZt+O01bFERCSbUbnJyjy84ZHvwO4B++bDtklWJ/pPZfL7MuuZRtx7fcK/QVO28ekSTfgnIiKZR+UmqytUzRxgDLD4Fbh01No8dyDAx4MfHqvLk9cn/Bv5xyH6/7SFqNgEi5OJiEh2oHLjDBoOguINIT4aZj0JjiSrE/0nN7uNoS0r8UWXIDzd7Szbe44O36zl+MUYq6OJiIiLU7lxBnY3aD8GPP3gxAZY/YXVie5Y+5rFmPZkMAX9vTgYHk2bkWtYfTDrDo4WERHnp3LjLHKXgJafmF+vGA6ntlibJxWCAnMxb2AINQJzEXEtgd7jNzJ+zVFN+CciIhlC5caZBD0KlduCIxGm94W4KKsT3bEC/t5M7d+AR2oVI8lh8M68Pbw8YydxiVn/EpuIiDgXlRtnYrNB668gIBAuH4WFL1qdKFW8Pdz4tFN1Xm9lTvj36+aTdPtuA+FRsVZHExERF6Jy42x8ckOH78zZi3dMgZ2/Wp0oVWw2G/0al2Z8n3r4e7uz5fhl2o5cQ+hJTfgnIiLpQ+XGGZUIhntfNr+e/7xT3B7+b/eWz5884d+Z6xP+zdWEfyIikg5UbpxV4xegeDDER8GMvpDkfHPIlM7vy+wBjbivQn7iEh08O2UbHy/epwn/RETkrqjcOCs3d/PylHeAeefU7+9ZnShN/L09+L53XZ66twwA36w4zBM/btaEfyIikmYqN84sVyC0HmF+veZLOLjM0jhp5Wa38UqLinz1aA283O0s3xdO+2/WcvSCJvwTEZHUU7lxdlXaQd1+5tez+kPEKUvj3I22NYoy7algCvl7cyg8mrYjV/PnwfNWxxIRESejcuMKHnwfClWHqxevj79JtDpRmlUvlou5AxtRs3guImMT6f3DRsat1oR/IiJy51RuXIGHN3SaYC7PELYO/njf6kR35a8J/zrWLobDgHfn7+HF6ZrwT0RE7ozKjavIWwbaXB9/s/pzpx1/8xcvdzc+6VidNx6ujN0G07ecpOvY9ZrwT0RE/pPKjSup2sFlxt+AOeFf35BSTLg+4d/WsCu0+XoNO09esTqaiIhkYSo3ruam8TfOf0v1PeXzM2dgCGUL+HI2MpZOY9YxZ7tzFzcREck4Kjeu5t/jb5a9bXWidFEqX05mPdOQ+ysWIC7RweCp2/lw0T6SNOGfiIj8i8qNK8pbBtp9Y369biTsnm1pnPTi5+3Bd73q8EwTc8K/MSvNCf8iNeGfiIj8g8qNq6rcBho+a349ZwCcP2BtnnTiZrfx0kN/T/j3+75w2o9aw5Hz0VZHExGRLELlxpU1fQtKNob4aPilB8S5TgFoW6Mo059qSOEAbw6fj6HdqDWsOqAJ/0REROXGtbm5Q8cfwK8wXNgPcweBC02GV61YAHMGNqJ2idxExiby2PiNfP/nEU34JyKSzancuDrfAtBpItjdYfdM2DDG6kTpqoCfN5OfqE/nOuaEf+8t2MsL03YSm6AJ/0REsiuVm+ygeH3zFnGA316H4+uszZPOvNzd+OiR6rzVujJudhsztp7k0bHrCY/UhH8iItmRyk12Uf9JqNoRHIkwrTdEnrY6Ubqy2Wz0aVSKiX3qEeDjwfYTV2g9cjU7TlyxOpqIiGQylZvswmaD1l9BgcoQfc4cYJzgemc2QsrlY86ARpQr4Mu5yDg6fbuO2ds04Z+ISHaicpOdePnCoz+Ddy44tQUWPO9SA4z/UjJfTmY+05BmlQoQn+hgyC/bGb5wryb8ExHJJlRusps8pc0ZjG122P4zbBxrdaIM4eftwdiedRhwnznh37erjtB34iYirmnCPxERV6dykx2VuQ8eeNf8evFQOPqntXkyiN1u48XmFfm6a028Peys2H+e9t9owj8REVencpNdBQ+A6l3ASDIHGF8JszpRhmkdVITpTzWkSIA3R87H0HbUGlbsD7c6loiIZBCVm+zqrwHGhYPMFcSndoP4q1anyjBViwYwZ2AIdUrkJio2kccnbOK7VZrwT0TEFancZGcePtDlZ8iRD86GwtyBLjnA+C/5/bz4+Yn6dKkTiMOA9xfu5X+/7tCEfyIiLkblJrvLFQidfzRnMN41A9Z8aXWiDOXl7saHj1Tj7esT/s3cdoouY9dzThP+iYi4DJUbgZKNoMVH5tfL3oF9C6zNk8FsNhuPNSrFj4/XI1cOD3acuELrr1ezNeyy1dFERCQdqNyIqW4/84EBM56AMzutTpThGpU1J/wrX9CX8Kg4Hv12PdM2n7A6loiI3CWVG/nbQx9C6SaQEANTukLUOasTZbgSeXMy85lGNK9SkPgkBy9O38k783aTmOSwOpqIiKSRyo38zc3DnOAvb1mIPGneQeWCSzT8m6+XO6O712Zw03IAjF9zjN7jN3I5Jt7iZCIikhYqN3Ijn9zQ7dfrSzRshjkDXPoOqr/Y7Taee6A8Y3rUIoenG2sOXaTtqDXsPxtldTQREUkllRu5Wd4y0OWn63dQTYdVn1qdKNM8VLUwM59pSGAeH8IuXaXDN2tYsvus1bFERCQVsky5+fDDD7HZbAwZMuS2+0yYMAGbzXbDw9vbO/NCZiel7oGWn5hf//Ee7J5taZzMVLGQP3MHhNCwTF5i4pN48qctfLXsIA4tvCki4hSyRLnZtGkT3377LdWrV//Pff39/Tlz5kzy4/jx45mQMJuq8zjUf9r8etZTcGqrtXkyUe6cnkx8vB6PNSwJwBfLDvDMz1uJiUu0NpiIiPwny8tNdHQ03bt357vvviN37tz/ub/NZqNQoULJj4IFC2ZCymys+ftQ9gFIvAZTHnXpNaj+zcPNztttqvDxI9XxcLOxePdZHhm9lhOXXHeZChERV2B5uRkwYACtWrWiWbNmd7R/dHQ0JUqUIDAwkLZt27J79+4U94+LiyMyMvKGh6SC3Q06/gAFq0L0Ofi5E1y7YnWqTNW5biBT+zcgn68X+85G0WbkatYeumB1LBERuQ1Ly83UqVPZunUrw4cPv6P9K1SowA8//MCcOXOYNGkSDoeDhg0bcvLkydu+Zvjw4QQEBCQ/AgMD0yt+9uHtb95B5VcEzu+DX3pAYva6Tbp2iTzMG9SI6sUCuHw1gZ4/bGTCmqNaeFNEJAuyGRb963zixAnq1KnD0qVLk8faNGnShBo1avDll1/e0XskJCRQqVIlunbtyrvvvnvLfeLi4oiLi0t+HhkZSWBgIBEREfj7+9/1z5GtnA2FH1pAfBRUfxTajzFXF89GYhOSGDozlFnbTgHQpU4gw9pVwcvdzeJkIiKuLTIykoCAgDv6/W1ZuZk9ezbt27fHze3vXwpJSUnYbDbsdjtxcXE3fO92OnXqhLu7O1OmTLmjz03NwZFbOLQMfu4MRhLc+zLc96rViTKdYRh8/+dRhi/ai8OAWsVzMaZnbQr46c49EZGMkprf35ZdlmratCmhoaFs3749+VGnTh26d+/O9u3b76jYJCUlERoaSuHChTMhsQBQthk8/IX59cqPYNska/NYwGaz8cQ9pRnfpx7+3u5sDbtCm6/XsPPkFaujiYgIFpYbPz8/qlatesMjZ86c5M2bl6pVqwLQq1cvhg4dmvyaYcOG8dtvv3HkyBG2bt1Kjx49OH78OP369bPqx8ieaveGxv8zv543GA7/YW0ei9xbPj9zBoZQtoAvZyNj6TRmHbO23X78l4iIZA7L75ZKSVhYGGfOnEl+fvnyZZ544gkqVapEy5YtiYyMZO3atVSuXNnClNnU/W9AtU7gSIRfe8G5lO9ac1Wl8uVk1jMNaVqxAHGJDp77ZQcfLNxLkib8ExGxjGVjbqyiMTfpKDEOfuoAx1eDf1Hotwz8i1idyhIOh8HnSw8w8o9DANxTPj9fP1qTgBweFicTEXENTjHmRlyAuxc8OgnylYfIU+ZA49jsOY+Q3W7jheYVGNmtJt4edlYdOE+7b9ZwKFwLb4qIZDaVG7k7Prmh+zTImR/OhcIv3c0zOtnUw9WLMOPphhTN5cPRCzG0G7WW5XvPWR1LRCRbUbmRu5e7pFlwPH3h6CqY2R8cSVanskyVIgHMHdiIeqXyEB2XSL8fNzPqj0Oa8E9EJJOo3Ej6KFITukwCuwfsmQ2LXoZs/Ms8r68Xk/rWp0eD4hgGfLJkP4OmbONqvBbeFBHJaCo3kn7K3AcdxgI22PQdrPrU6kSW8nS38167arzfvirudhvzd56h4+h1nLyshTdFRDKSyo2kr6odoMVH5td/vAebx1ubJwvoXr8Ek59oQN6cnuw5E0nbkWvYcOSi1bFERFyWyo2kv/pPQuMXzK8XPA9751mbJwuoVyoPcweFUKWIPxdj4un+/QYmrT9udSwREZekciMZ4/7XoVYvMBwwvS8cW211IssVzeXD9Kca0jqoCIkOg9dn7+LVWaHEJzqsjiYi4lJUbiRj2GzQ6guo0AqS4mBKV3NV8WzOx9ONEY/W4KWHKmCzweQNYfT4fgMXorPv7fMiIulN5UYyjps7dBwHxYMhLhImPQKXj1mdynI2m41nmpRlXO86+Hm5s/HYJdp8vZpdpyKsjiYi4hJUbiRjefhA1ylQoDJEnzOXa4gOtzpVlnB/xYLMGtCI0vlycjoilo5j1jJn+ymrY4mIOD2VG8l4PrmhxwwIKA6XDsNP7eHaZatTZQllC/gya0AjmlTIT2yCg8FTt/PBwr0kJmkcjohIWqncSObwLwK9ZkPOAnBuF0zqCHFadwkgwMeDcb3r8nSTMgCMXXWEPhM2ceVqvMXJRESck8qNZJ68ZcyC450LTm02BxknxFqdKktws9t4+aGKjOxWEx8PN/48eIHWI1ez90z2XIhURORuqNxI5ipYBXrMNNehOvYnTOsNSQlWp8oyHq5ehJnPNCQwjw8nLl2jwzdrmb/ztNWxREScisqNZL5itaHrVHD3hgOLYdaT2XqhzX+rVNifuQNCaFwuH9cSkhg4eRsfLtpHkiP7rtUlIpIaKjdijVKNofNPYHeHXTNg/pBsvdDmv+XO6cn4x+rS/57SAIxZeZg+EzYRcVVnuURE/ovKjVin/IPwyPdgs8PWH+G311Vw/sHdzc6rLSvx1aM18Paws+rAedqMWs2BcxqILSKSEpUbsVaV9tB6hPn1upGw8mNr82RBbWsUZcbTDSmay4fjF6/SbtQaFu86Y3UsEZEsS+VGrFerJzz0ofn1ig9g3Shr82RBVYoEMG9QCA3L5OVqfBJPTdrKp0v249A4HBGRm6jcSNbQ4Gm473Xz6yWvwuYfrM2TBeXJ6cmPj9ejb0gpAEb+cYh+P24m4prG4YiI/JPKjWQd97wAjQabX89/DrZMtDZPFuTuZueNhyvzRZcgvNzt/L4vnHaj1nAoXONwRET+onIjWYfNBs3egQbPmM/nDYZtk6zNlEW1r1mMGU83pEiAN0cvxNBu1Fp+233W6lgiIlmCyo1kLTYbNP8A6j0JGDBnIGyfYnWqLKlq0QDmDgqhfqk8RMcl0v+nLXyx9IDG4YhItqdyI1mPzQYtPoK6/QADZj8NO3+1OlWWlM/Xi0n96vNYw5IAfLX8IP1/2kJUrMbhiEj2pXIjWZPNBi0+gdp9AMOcxTh0utWpsiQPNztvt6nCp52C8HS3s2zvOdqNWsPh89FWRxMRsYTKjWRddju0+hxq9QLDATP7w+5ZVqfKsjrWLsa0J4MpHODN4fMxtBu5hmV7zlkdS0Qk06ncSNZmt8PDX0GNHmAkwfS+sGeO1amyrKDAXMwdGEK9knmIikuk34+bGbH8oMbhiEi2onIjWZ/dDm1GQFDX6wXncdg73+pUWVZ+P3McTq/gEgB8vvQAT03aQnRcosXJREQyh8qNOAe7G7QdBdU6gyMRpj0G+xdZnSrL8nS3M6xtVT56pBqebnZ+23OO9qPWcPRCjNXRREQynMqNOA+7G7QbDVUfAUcC/NITDiyxOlWW1qVucaY+2YCC/l4cDI+mzcjV/LEv3OpYIiIZSuVGnIubO7QfC5XbXS84PWD/YqtTZWm1iudm3sAQapfITVRsIo9P3MTXGocjIi5M5Uacj5s7PPI9VGoDSfFmwdkz1+pUWVoBf2+mPNGAbvWLYxjw2dIDPDlJ8+GIiGtSuRHn5OYBHX/4+xLVtMc0D85/8HS380H7anzYwRyHs3TPOdpqXSoRcUEqN+K83Dygw3cQ1M28i2pGP9j2s9WpsrxH6xXn16fM+XCOnI+h7cg1LN51xupYIiLpRuVGnNtfd1HVfgxzLapnYPN4q1NleTUCczFvUAgNSuchJj6JpyZt5aPF+0jSOBwRcQEqN+L87HZ4+Euo/5T5fP4QWD/GykROIZ+vF5P61qdfSCkARq84zGPjN3I5Jt7iZCIid0flRlyDzQYPfQiNBpvPF78Mq7+0NJIzcHez8/rDlfnq0Rp4e9j58+AFWo9cza5TEVZHExFJM5UbcR02GzR7B+592Xy+7C1Y8REYutTyX9rWKMqsZxpRPE8OTl6+xiOj1zJr20mrY4mIpInKjbgWmw3uexWavmk+X/EBLB+mgnMHKhX2Z97AEJpUyE9cooPnftnB23N3k5DksDqaiEiqqNyIa2r8P2j+gfn16s9hyWsqOHcgIIcH43rX5dn7ywIwYe0xun+3gfCoWIuTiYjcOZUbcV3BA6DVZ+bX60fBgv+BQ2ch/oub3cbzD1bgu1518PNyZ+OxS7T+ejVbjl+2OpqIyB1RuRHXVrcftBkJ2GDzOJgzAJK0OvadeKByQWYPbETZAr6ci4zj0bHrmLT+OIbOgIlIFqdyI66vVk/oMBZsbrBjMvzaCxJ0meVOlMnvy+wBjWhRtRAJSQavz97FyzN2EpuQZHU0EZHbUrmR7KF6Z+gyCdy8YP8C+LkjxEZancop+Hq58033Wrz8UEXsNvh180k6f7uOU1euWR1NROSWVG4k+6jYEnrOBC9/OPYnTGwNMResTuUUbDYbTzcpw8TH65Erhwc7T0bQ+uvVrD2s4yciWU+WKTcffvghNpuNIUOGpLjftGnTqFixIt7e3lSrVo2FCxdmTkBxDSVDoPc8yJEPzmyHH5rDlRNWp3IajcvlZ97AEKoU8edSTDw9vt/Ad6uOaByOiGQpWaLcbNq0iW+//Zbq1aunuN/atWvp2rUrffv2Zdu2bbRr14527dqxa9euTEoqLqFIDXh8CQQEwsVDZsE5f8DqVE4jME8OZjzdkA61iuIw4P2Fexk0ZRtX4zVQW0SyBpuRxv/LtWfPHsLCwoiPv3EdmjZt2qTqfaKjo6lVqxbffPMN7733HjVq1ODLL7+85b5dunQhJiaG+fPnJ29r0KABNWrUYMyYO1tLKDIykoCAACIiIvD3909VVnExEafgp/ZwYT/45IEeM6BoLatTOQ3DMPhp/XGGzdtDosOgQkE/vu1Zm5L5clodTURcUGp+f6f6zM2RI0cICgqiatWqtGrVKvnsSfv27Wnfvn2qww4YMIBWrVrRrFmz/9x33bp1N+3XvHlz1q1bd9vXxMXFERkZecNDBICAotBnERSpBdcumWNwjqy0OpXTsNls9AouyZT+Dcjv58X+c1G0Hrma3/edszqaiGRzqS43gwcPplSpUoSHh5MjRw52797NqlWrqFOnDitWrEjVe02dOpWtW7cyfPjwO9r/7NmzFCxY8IZtBQsW5OzZs7d9zfDhwwkICEh+BAYGpiqjuLiceaH3XCh1L8RHm3dR7Z1ndSqnUrdkHuYPCqF2idxExSbSd+Jmvlp2EIdD43BExBqpLjfr1q1j2LBh5MuXD7vdjt1uJyQkhOHDh/Pss8/e8fucOHGCwYMH8/PPP+Pt7Z3aGHds6NChREREJD9OnNDgUfkXLz/oPg0qtYakeHMenK0/WZ3KqRT092bKEw3o2aAEhgFfLDtA/582ExmbYHU0EcmGUl1ukpKS8PPzAyBfvnycPn0agBIlSrB///47fp8tW7YQHh5OrVq1cHd3x93dnZUrVzJixAjc3d1JSrp5krBChQpx7tyNp7zPnTtHoUKFbvs5Xl5e+Pv73/AQuYm7F3ScADV7guGAuQNh7ddWp3Iqnu523m1XlU86VsfT3c6yveG0HbmGA+eirI4mItlMqstN1apV2bFjBwD169fn448/Zs2aNQwbNozSpUvf8fs0bdqU0NBQtm/fnvyoU6cO3bt3Z/v27bi5ud30muDgYJYvX37DtqVLlxIcHJzaH0PkZm7u0OZraHj9DORvr5sLbmo9qlTpVCeQGU81pGguH45eiKHdqDXM23Ha6lgiko24p/YFr7/+OjExMQAMGzaMhx9+mMaNG5M3b15++eWXO34fPz8/qlatesO2nDlzkjdv3uTtvXr1omjRosljcgYPHsy9997LZ599RqtWrZg6dSqbN29m7Nixqf0xRG7NZoMH34UceWHZW7BuJESehvZjzLM7ckeqFQtg7sBGPDt1G2sOXWTQlG1sDbvMqy0r4eGWJWagEBEXlup/ZZo0aULz5s0BKFu2LPv27ePChQuEh4dz//33p2u4sLAwzpw5k/y8YcOGTJ48mbFjxxIUFMT06dOZPXv2TSVJ5K6FDIH2Y8HuAbtnwk8d4JpWxU6NvL5eTOxTj6eblAFg/JpjdB27nnORWtdLRDLWHc9zc/78eXr16sWyZctwOBzUrVuXSZMmUbZs2YzOmK40z42kypEV8EtPiIuE/BWh+3TIpTvuUuu33Wf53687iIpLJJ+vF6O61aR+6bxWxxIRJ5Ih89y8/PLLbN++nWHDhvHpp59y5coVnnjiibsOK5KllW5izoXjVwTO74Pvm8HZUKtTOZ0HqxRi7qAQKhby40J0HN20bIOIZKA7PnMTGBjI999/n3xJ6uDBg1SqVImYmBi8vJxnLILO3EiaRJyESR3h/F7w9IMuP0GZ+6xO5XSuxify2qxdzNp2CoCW1QrxcccgfL1SPfxPRLKZDDlzc/r0aYKCgpKflytXDi8vrxvGxIi4rIBi8PhiKNkY4qPMyf52TLU6ldPJ4enO552DGNa2Ch5uNhaGnqXtyNUcCtft4iKSflI1oPjft2e7ubnptLJkHz65zPWnqj4CjkSY9ST8+Rno70Cq/LVsw9T+wRTy9+bw+RjajFzD/J26XVxE0scdX5ay2+0EBARgs9mSt125cgV/f3/s9r870qVLl9I/ZTrSZSm5aw4HLH8b1nxlPq/zOLT4xJwnR1LlQnQcgyZvY92RiwD0DSnFKy0q6nZxEblJan5/3/G/xuPHj7/rYCIuwW6HB4aBfzFY9BJs/gEiz0DHceCpFbFTI5+vFz/1rcenvx1gzMrDjFt9lNCTEYzsVpMC/hm3LIuIuLY7PnPjKnTmRtLV3nkwox8kxkLR2tD1F/DNb3Uqp7R411lemLaD6LhE8vt58U33WtQtmcfqWCKSRWTIgGIRuYVKraHXXPDJDae2wPf3w7k9VqdySg9VLcTcgY0oX9CX81FxdB27nnGrj2pcn4ik2h2fubnTdaOOHDlyV4Eyms7cSIa4cBAmd4ZLR8xbxTtNgHLNrE7llK7GJ/LKjFDmXl+P6uHqhfnokerk1O3iItlaan5/p2pAcYkSJejWrRsFChS47X6DBw9OXdpMpnIjGebqJXM24+OrwWaHhz6E+k9ancopGYbBxLXHeG/BXhIdBmUL+DKmR23KFvC1OpqIWCRDys20adP44YcfWLFiBS1atODxxx+nZcuWN9wp5QxUbiRDJcbD/Odg+yTzed0nzJKjO6nSZMvxSzzz81bORcaR09ONTzsF0aJaYatjiYgFMqTc/OXUqVNMmDCBCRMmcPXqVXr27Enfvn0pV67cXYXOLCo3kuEMw7xNfNnbgAFl7jcvU3kHWBzMOZ2PimPg5K1sOGpOM9H/ntK81LwC7rpdXCRbydBy808rV67k7bffZtWqVVy4cIHcuXOn9a0yjcqNZJq982HmE5BwFfJVgG6/QJ5SVqdySolJDj5Zsp9vV5lj+uqXysPX3WpSwE+3i4tkFxl+t1RsbCyTJk3inXfeYcOGDXTq1IkcOXKkKayIy6r08PVFNwvDhf3wfVM4vs7qVE7J3c3O0JaVGN29Fr5e7mw4eomHR6xm87GsPWmoiFgjVeVmw4YN9O/fn0KFCvH555/ToUMHTp06xdSpU51q8UyRTFOkBjzxOxQOgqsX4cc2sH2K1amcVotqhZkzsBHlCvgSHhXHo2PXM36NbhcXkRvd8WWpKlWqEB4eTrdu3Xj88cdvWETTmeiylFgiPsZci2rvPPN54xfgvtfM2Y4l1WLiEnl5xk7m7zQX7m0TVIThHarpdnERF5Zht4LnzJkTd3f3G9aX+jetLSVyGw4H/P4urP7cfF65LbQbA566pJsWhmEwfs0xPlho3i5evqB5u3jp/LpdXMQVZUi5mThx4h19eO/eve9oP6uo3Ijltk+Guc+CIwEK14BHf4aAYlanclqbjl1iwM9bCY+Kw9fLnU86Vtft4iIuKNPulnJGKjeSJRxbA7/0gGuXIGd+6PwTlAi2OpXTCo+KZeDkbWy8fru4VhcXcT1aW0okqyvZCPr/AQWrQsx5mPgwbBpndSqnVcDPm8n96vPkPeYyMeNWH+XRses5GxFrcTIRsYLKjYhVcpeEvr9BlfbgSIQFz8O8wZAYZ3Uyp/TX7eLf9qyNn7c7W45fptWIP1l98ILV0UQkk6nciFjJMyd0HA/N3gZssGUCTGwNUWctDua8mlcpxPxBIVQu7M/FmHh6/rCBEcsP4nBkqyvwItmayo2I1Ww2CHkOuk8DrwA4sQHGNoGTW6xO5rRK5M3JzGca0rVeIIYBny89QJ8Jm7gcE291NBHJBHdcbho3bsynn37KgQMHMjKPSPZV7gFzHE6+ChB1BsY/BNt+tjqV0/L2cGN4h+p82ikIbw87Kw+cp9WIP9kWdtnqaCKSwe643DzxxBOsW7eO2rVrU6lSJV5++WXWrFmjmUFF0lPeMtBvGVRoBUnxMOcZWPQyJCVYncxpdaxdjFnPNKJUvpycjoil87frmLj2mP7tEnFhqb4VPC4ujuXLlzNnzhzmzZtHUlISrVq1ok2bNjRv3hwfH5+MypoudCu4OAWHA1Z9DCuGm89LNjZXFs+Zz9JYziwqNoGXZ+xkYag5nqn19VmNfTWrsYhTyNR5bjZs2MDcuXOZO3cuhw8f5v7772fo0KE0atTobt42w6jciFPZO99ctiE+GgKKw6OTzHWqJE0Mw+CHNccYfn1W4zL5czK6R23KF/SzOpqI/AfLJvE7fPgwc+fOJTAwkI4dO6bX26YrlRtxOuH7YGpXuHQE3H2g7UioljX/fjmLLccvMeDnbZyNjMXHw43hHarRrmZRq2OJSAo0Q3EKVG7EKV27AjP6wqFl5vN6T8KD74G7p6WxnNnF6DiG/LKdP6/Pg9O9fnHeeLgy3h5uFicTkVvRDMUirsYnF3T7FRr/z3y+8VuY0BIiTloay5nl9fViQp96DG5aDpsNft4QRqcx6zhx6arV0UTkLqnciDgLuxs0fRO6TgXvADi5Cb69Bw7/bnUyp+Vmt/HcA+UZ/1hdcufwIPRUBA9/vZrle89ZHU1E7oLKjYizqdAC+q80BxZfvQg/dYAVH5l3WEmaNKlQgPnPNqZGYC4iriXQd+JmPl68j8QkHVMRZ3RX5cYwDM0VIWKFPKXg8d+gVm/AgBUfwOROcPWS1cmcVtFcPvz6ZDCPNSwJwDcrDtNj3AbCo7T4poizSVO5GTduHFWrVsXb2xtvb2+qVq3K999/n97ZRCQlHt7QZgS0G23eRXVomXmZSss2pJmnu52321RhZLea5PR0Y/2RSzw8YjUbjly0OpqIpEKqy82bb77J4MGDad26NdOmTWPatGm0bt2a5557jjfffDMjMopISmp0M2c1zlMaIk7AD81h43egs6pp9nD1IswZGEL5gr6ER8XR7fsNfLvysM5UiziJVN8Knj9/fkaMGEHXrl1v2D5lyhQGDRrEhQsX0jVgetOt4OKyYiNgzgDYO898Xq0TPPwlePlaGsuZXY1P5LVZu5i17RQAD1QuyKedggjw8bA4mUj2k6G3gickJFCnTp2btteuXZvExMTUvp2IpBfvAOj8Ezz4PtjcIHQafN8Uzmux27TK4enO552D+KB9NTzd7Czdc47WX69m16kIq6OJSApSXW569uzJ6NGjb9o+duxYunfvni6hRCSNbDZoOBAemw++heD8PvjuPtg1w+pkTstms9GtfnFmPN2QYrl9CLt0lQ6j1zJlY5guU4lkUam+LDVo0CB+/PFHAgMDadCgAWCuLxUWFkavXr3w8Pj7dO3nn3+evmnTgS5LSbYRdc6c1fjYn+bz2o9B8+HgmcPSWM4s4moCz/+6neX7wgFoX7Mo77WrSk4tvimS4TJ0+YX77rvvjvaz2Wz8/nvWm1xM5UaylaRE8zbxPz8HDMhfCTqNhwKVrE7mtBwOg29XHeHT3/aTpMU3RTKN1pZKgcqNZEuH/4CZ/SEm3LxtvMVHUKuXeRlL0mTj0UsMmrKVc5FxeHvYea9dNTrWLmZ1LBGXpbWlRORGZe6Dp9dAmfsh8RrMexamP27eYSVpUq9UHhY825jG5fIRm+DghWk7eHHaDq7FJ1kdTSTbU7kRyS58C0D3GdDsHbC7w+6Z5qR/pzTpX1rl8/ViYp96/O+B8thtMG3LSdqNWsOh8Giro4lkayo3ItmJ3Q4hQ6DPYshVHC4fg3EPwtqvtTZVGtntNgY1LcekfvXJ5+vF/nNRtBm5mtnX58YRkcynciOSHQXWhSf/hEptwJEIv70OU7pATNaehDMra1gmHwsHhxBcOi9X45MY8st2hs4MJTZBl6lEMpvKjUh25ZMLOv8ID38B7t5w8DcY3QiOrrI6mdMq4OfNpH71ebZpOWw2mLIxjPbfrOXohRiro4lkKyo3ItmZzQZ1Hocnfod8FSD6LExsA398YN5GLqnmZrfx/APl+fHxeuTN6cneM5G0/no183eetjqaSLZhabkZPXo01atXx9/fH39/f4KDg1m0aNFt958wYQI2m+2Gh7e3dyYmFnFRBatA/z+gZk/AgJUfwcTWcOWE1cmcVuNy+Vk4uDH1SuYhOi6RgZO38cbsXcQl6jKVSEaztNwUK1aMDz/8kC1btrB582buv/9+2rZty+7du2/7Gn9/f86cOZP8OH78eCYmFnFhnjmh7Uh4ZBx4+kHYWvMyVeh0q5M5rYL+3kx+oj7PNCkDwE/rj/PI6LWEXbxqcTIR15blJvHLkycPn3zyCX379r3pexMmTGDIkCFcuXLljt8vLi6OuLi45OeRkZEEBgZqEj+RlFw6AjOegFObzedVO0Krz8xxOpImf+wP5/lftnP5agJ+3u580jGIh6oWsjqWiNNwykn8kpKSmDp1KjExMQQHB992v+joaEqUKEFgYOB/nuUBGD58OAEBAcmPwMDA9I4u4nrylIbHl0CToeYK47uma7DxXbqvQgEWPNuY2iVyExWbyFOTtvDOvN3EJ+oWfJH0ZvmZm9DQUIKDg4mNjcXX15fJkyfTsmXLW+67bt06Dh48SPXq1YmIiODTTz9l1apV7N69m2LFbj3tuc7ciNylk5th5hPm2Rxs0HAQ3P86uHtZncwpJSQ5+GTJfsauOgJAUGAuRnatSWAeLWgqkhKnWlsqPj6esLAwIiIimD59Ot9//z0rV66kcuXK//nahIQEKlWqRNeuXXn33Xfv6PO0tpRIGsRFw5JXYetE83nBatBhLBT877+ncmvL9pzjf9N2EHEtAX9vdz7rXIMHKhe0OpZIluVU5ebfmjVrRpkyZfj222/vaP9OnTrh7u7OlClT7mh/lRuRu7BvAcwdBFcvgpsXPPAO1HvSnPlYUu3EpasMnLKNHSeuAND/ntK82LwCHm46niL/5pRjbv7icDhuuIyUkqSkJEJDQylcuHAGpxIRACq2gqfXQbkHISkOFr8CkzpApOZwSYvAPDmY9mQwjzcqBcDYVUd4dOx6Tl+5ZnEyEedmabkZOnQoq1at4tixY4SGhjJ06FBWrFhB9+7dAejVqxdDhw5N3n/YsGH89ttvHDlyhK1bt9KjRw+OHz9Ov379rPoRRLIfv4LQ7Vfz7il3HzjyB3wTDLtnW53MKXm623mzdWXG9KiFn7c7W45fptWIP/ljf7jV0USclqXlJjw8nF69elGhQgWaNm3Kpk2bWLJkCQ888AAAYWFhnDlzJnn/y5cv88QTT1CpUiVatmxJZGQka9euvaPxOSKSjmw2qNsPnlwFhWtA7BWY1htmPQ2xkVanc0oPVS3MgkGNqVrUn8tXE+gzfhMfLd5HYpLuphJJrSw35iajacyNSDpLjDdnNF79ORgOc7XxdqOhZIjVyZxSbEIS7y/Yy0/rzQlK65XMw4iuNSkUoNnYJXtz6jE3IuJk3D2h6Rvw2EKz2FwJgwmtYNHLEK8FI1PL28ONd9tV5euuNfH1cmfjsUu0HPEnf+zTZSqRO6VyIyLpo0QwPLUGavU2n28YA2NC4Pg6a3M5qdZBRZg3KIQqRfy5FBNPnwmbGL5wLwm6TCXyn1RuRCT9ePtDmxHQYwb4FzUn/hvfAha/Cgm6Ayi1SuXLyYynG9I7uAQA3646Qudv13HystamEkmJyo2IpL+yzeCZdVCzB2DA+lHmWZwTG61O5nS8Pdx4p21VRnc376baFnaFViNW89vus1ZHE8myVG5EJGN4B0DbUdBtGvgVhouH4Ifm8NsbkBBrdTqn06JaYRY+25igwFxEXEug/0/m2lRxiUlWRxPJclRuRCRjlX/QPIsT1NW8m2rtCPi2sblmlaTKX5P+PdHYnPRv/JpjdBy9juMXNXBb5J9UbkQk4/nkhvZjoOtU8C0IFw7AuAdg2duQeGczkovJ093Oa60qM653HXLl8CD0VAQPj1jN/J2aJVrkLyo3IpJ5KrSAZ9ZDtc7mWZzVX8C398KprVYnczpNKxVk4bONqVMiN1FxiQycvI3XZoUSm6DLVCIqNyKSuXLkgUe+gy6TIGd+OL8Xvm8Gy4dpLE4qFcnlw9T+DRhwXxlsNvh5QxjtRq3h8Ploq6OJWErlRkSsUak1PLMBqj4CRhL8+Zk5Fkfz4qSKu5udF5tXZGKfeuTN6cm+s1G0/no1s7adtDqaiGVUbkTEOjnzQscfoPOPf4/FGf8QzH9ea1Sl0j3l87NocGOCS+flanwSz/2ygxen7eBqfKLV0UQyncqNiFivclsYsAFq9jSfbx4Ho+rDvoXW5nIyBfy9mdSvPs81K4/dBtO2nKTtyDUcOBdldTSRTKVyIyJZg09uaDsSes2F3KUg6jRM7QrTHoNorat0p9zsNgY3K8fP/RpQwM+Lg+HRtBm5ml83nSCbrZMs2ZjKjYhkLaXvNefFaTQEbG6wexaMrAvbJoF+Od+x4DJ5WTi4MfeUz09sgoOXZuzkuV+2Ex2ny1Ti+mxGNqvyqVkyXUQsdno7zB0EZ3eaz0vdC62/hDylrUzlVBwOgzGrDvPZbwdIchiUypeTkd1qUqVIgNXRRFIlNb+/deZGRLKuIjXgiT/ggWHg7g1HV8I3DWHNCEjSGYg7YbfbeKZJWX7p34AiAd4cvRBD+2/W8tO6Y7pMJS5LZ25ExDlcPAzzh8DRVebzwkHQZiQUrm5pLGdyOSaeF6fvYNlecwxTy2qF+PCR6vh7e1icTOS/6cyNiLievGXMwcZtR5mLcp7ZAWObwNK3IP6q1emcQu6cnnzXqw6vt6qEh5uNhaFnaTXiT3acuGJ1NJF0pXIjIs7DZoOaPWDAJqjczpz8b82X8E192L/Y6nROwWaz0a9xaaY91ZDAPD6cuHSNR0av5btVR3A4stWJfHFhuiwlIs5r30JY+CJEXp+Nt+LD8NCHkCvQ2lxOIuJaAkNn7mRh6FkAmlTIz6edgsjn62VxMpGb6bKUiGQPFVuak/81fBbs7rBvvjn535qvICnB6nRZXoCPB6O61eKD9tXwcrezYv95Wn71J2sPXbA6mshd0ZkbEXEN5/bAguch7PraVAUqQ6vPoUSwtbmcxP6zUQycvJWD4dHYbPBMkzI816w87m76/8CSNejMjYhkPwUrw2MLzQHHPnkgfI+5TtXsARBz0ep0WV6FQn7MHRhC13rFMQwY9cdhuoxdz8nLGqwtzkdnbkTE9Vy9BMvegq0/ms99ckOzd8y1q+z6/3T/ZcHOM7wyYydRcYn4e7vzccfqPFS1sNWxJJtLze9vlRsRcV1hG8xLVed2mc+L1YOHP4dC1azN5QROXLrKoCnb2H79NvEeDYrzeqvKeHu4WRtMsi2VmxSo3IhkM0mJsGEMrBgO8dHmelX1n4L7hoKXn9XpsrSEJAef/XaAMSsPA1CxkB8ju9WkbAEdN8l8KjcpULkRyaYiTsHiV2DvXPO5XxFo/j5UaW/OnyO3terAeZ7/dTsXouPx9rDzTpsqdK4TiE3HTTKRyk0KVG5EsrmDS2HhC3D5mPm8RAi0+AgKVbU0VlYXHhXL/37dwZ8HzdvEWwcV4f32VbV0g2QalZsUqNyICAnXYPWX5uzGibFgs0PdftBkKOTIY3W6LMvhMPh21RE++20/iQ6D4nly8HXXmgQF5rI6mmQDKjcpULkRkWRXwmDJa39fqvLJA03fgFq9wa6Bs7ezNewyz07ZxsnL13C323jpoQr0CymN3a7LVJJxVG5SoHIjIjc5sgIWvQLn95rPC1WHlp9A8QaWxsrKIq4l8OrMUBaEngHg3vL5+ayzlm6QjKNykwKVGxG5paQE2DQO/vgA4iLMbdU6wwPvgH8Ra7NlUYZhMGXjCd6Zt5u4RAf5/bz4onMNQsrlszqauCCVmxSo3IhIimIuwPJh1ycANMAjJ9zzAgQPAHedlbiV/WejGDRlKwfO/b10w5Bm5fHQ0g2SjlRuUqByIyJ35PQ2WPgSnNxoPs9T2lxxvHxza3NlUdfikxg2fw9TNoYBUKt4LkZ0rUmx3DksTiauQuUmBSo3InLHHA4I/RWWvgnR58xt5R6E5sMhX1lrs2VRC3ae4ZWZO4mKNZdu+OiR6rSopqUb5O6p3KRA5UZEUi0uClZ9Auu+AUcC2N3NW8fvfVm3jt/Cv5du6FqvOG8+XBkfT92BJmmncpMClRsRSbMLh2DJq3BwifncO5dZcOr2A3dPS6NlNQlJDj5fai7dYBhQtoAvIx6tSeUi+ndX0kblJgUqNyJy1w7/Dkteh/Dd5vM8peGBd6FiKy3l8C9rDl3guV+2Ex4Vh6e7nVdbVKR3w5JaukFSTeUmBSo3IpIuHEmwbRL8/h7EhJvbSjQy16sqUtPabFnMxeg4Xpq+k+X7zOPUtGIBPu5YnbyaE0dSQeUmBSo3IpKu4qJgzVew9mtzKQeAoK5w/xsQUNTabFmIYRj8uO447y/cS3yigwJ+XnzRpQaNympOHLkzKjcpULkRkQwRcdKcH2fnL+Zzdx9o9Cw0fBa8fK3NloXsPRPJoCnbOBRuzonz5D1l+N+DmhNH/pvKTQpUbkQkQ53aYq5XFbbOfO5byFyvKqir1qu67t9z4gQVC2BE15qUyJvT4mSSlancpEDlRkQynGGYi3EufRMuHzO3FaoGD74HpZtYmSxLWbzrDC9N30lkbCI5Pd14r31V2tcsZnUsyaJUblKgciMimSYxDjaOhZWf/L1eVZmm0OxtKFzd0mhZxakr13hu6nY2HrsEQPuaRRnWtgp+3h4WJ5OsRuUmBSo3IpLpYi7Cyo9g8w/mJIBgLsp5/2uQu6Sl0bKCJIfByN8P8dXyAzgMKJ4nByO61qRGYC6ro0kWonKTApUbEbHMpSPw+/uwa7r53O5hTgB4zwuQU3cNbT52icFTt3PqyjXc7Taef7A8T91TBrtdc+JI6n5/Wzo8ffTo0VSvXh1/f3/8/f0JDg5m0aJFKb5m2rRpVKxYEW9vb6pVq8bChQszKa2IyF3KUxo6joP+K6H0feZZnA2j4asa5qWr+BirE1qqTsk8LBzcmFbVCpPoMPh48X56/rCBc5GxVkcTJ2NpuSlWrBgffvghW7ZsYfPmzdx///20bduW3bt333L/tWvX0rVrV/r27cu2bdto164d7dq1Y9euXZmcXETkLhSpAb1mQ89ZUDgI4qPgj/dgRE3YNA6SEqxOaJkAHw9GdqvJR49Uw8fDjTWHLtLiqz9Zvvec1dHEiWS5y1J58uThk08+oW/fvjd9r0uXLsTExDB//vzkbQ0aNKBGjRqMGTPmlu8XFxdHXFxc8vPIyEgCAwN1WUpEsgaHA3bPhN/f/fvOqjxloOmbULlttl7O4VB4NM9O2caeM5EAPNawJK+0qIi3h26pz46c5rLUPyUlJTF16lRiYmIIDg6+5T7r1q2jWbNmN2xr3rw569atu+37Dh8+nICAgORHYGBguuYWEbkrdjtU6wgDNkGLTyBHPrh0GKb1hu+bwtE/rU5ombIFfJk1oCGPNyoFwIS1x2g3ag0Hz0VZnEyyOsvLTWhoKL6+vnh5efHUU08xa9YsKleufMt9z549S8GCBW/YVrBgQc6ePXvb9x86dCgRERHJjxMnTqRrfhGRdOHuCfX7w+Dt5krjHjnNCQEnPgw/dYBTW61OaAkvdzfebF2Z8Y/VJW9OT/adjaL1yNVMWn+cLHbhQbIQy8tNhQoV2L59Oxs2bODpp5+md+/e7NmzJ93e38vLK3nA8l8PEZEsy8sP7nsVnt1m3klld4fDy+G7+2BqdziXfv8+OpP7KhZg0ZDGNC6Xj9gEB6/P3kX/n7ZwKSbe6miSBVlebjw9PSlbtiy1a9dm+PDhBAUF8dVXX91y30KFCnHu3I2Dys6dO0ehQoUyI6qISObxKwitPoOBm6D6o4AN9s2H0Q1hel+4eNjqhJmugJ83E/vU4/VWlfB0s7N0zzmaf7mKPw+etzqaZDGWl5t/czgcNwwA/qfg4GCWL19+w7alS5fedoyOiIjTy1MaOnwLz6w3BxhjmPPkjKwLcwbClTCrE2Yqu91Gv8almTWgIWUL+HI+Ko6e4zby3vw9xCUmWR1PsghLy83QoUNZtWoVx44dIzQ0lKFDh7JixQq6d+8OQK9evRg6dGjy/oMHD2bx4sV89tln7Nu3j7fffpvNmzczcOBAq34EEZHMUaAidP7RnCOn3INgJMG2n+Dr2rDwRYi6/dhDV1SlSADzBobQo0FxAL5ffZT2o9ZyKFyDjcXichMeHk6vXr2oUKECTZs2ZdOmTSxZsoQHHngAgLCwMM6cOZO8f8OGDZk8eTJjx44lKCiI6dOnM3v2bKpWrWrVjyAikrmK1IDu0+Dx36BkY0iKN9ev+qoG/PYGXL1kdcJM4+PpxnvtqvFdrzrkyenJnjORtBqxmp802Djby3Lz3GQ0Lb8gIi7lyEpzjpyTm8znnn4QPACCnwHvAGuzZaLwyFj+N20Hfx68AECzSgX56JFq5PX1sjiZpBetLZUClRsRcTmGAQd/M0vO2VBzm3cuaDgQ6j0J3tnj3zqHw+CHNUf5ePF+4pMc5Pfz4vPOQTQul9/qaJIOVG5SoHIjIi7L4YC9c+CPD+DCAXNbNiw5u09HMHjqdg6FRwPQL6QULz5UAS93zWzszFRuUqByIyIuz5EEu2bCyo/g4kFzm3cuCB4I9bNHybkWn8QHC/fy0/rjAFQq7M+IR2tQrqCfxckkrVRuUqByIyLZhiMJds8yS84/z+QEDzRnQ84GY3KW7TnHSzN2cikmHi93O68/XJke9Ytjy8ZrdjkrlZsUqNyISLZzy5IT8I8zOa5dcm4ebFyAjx6prsHGTkblJgUqNyKSbSWXnI/hwn5zm3cANBgADZ5y6ZLjcBiMX3uMjxbtSx5s/FmnIO4pr8HGzkLlJgUqNyKS7d225DwD9Z8Cn1yWxstIe05HMnjqNg5eH2zcN6QUL2mwsVNQuUmByo2IyHWOJNgz2yw55/eZ2zz9oF4/82yOr2ue1YhNMAcb/7hOg42dicpNClRuRET+xeGAPbNg1acQfn3VcXcfqN0bGj4LAUWtzZdBlu89x4vT/zHYuFUlejQoocHGWZTKTQpUbkREbsPhgAOLzJJzequ5ze4BNbpCyHPmIp4uJjwqlhem7WTVAXNlcQ02zrpUblKgciMi8h8MA478Aas+g+OrzW02O1R9BBr/DwpUsjZfOnM4DCasPcaH1wcb5/P14rPOQdyrwcZZispNClRuRERS4fg6+PMzOLT0720VH4Z7XoAiNa3LlQH2nonk2Sl/DzZ+vJE52NjbQ4ONswKVmxSo3IiIpMHpbWbJ2Tvv721lmppncko0BBcZpxKbkMTwhXuZeH2wccVCfozoWpPyGmxsOZWbFKjciIjchfB9sPpzCJ0ORpK5rVg9CBkC5VuA3W5pvPTy+75zvDhtJxevDzYe2qIivRuW1GBjC6ncpEDlRkQkHVw6Cmu+hO2TISne3JavAjR6Fqp1BndPS+Olh/CoWF6avpMV+83BxveUz88nHatT0N/b4mTZk8pNClRuRETSUdRZWD8aNv8AcZHmNr/C5oSAtR9z+kU6DcNg0vrjvLdgL3GJDnLn8GB4h2o8VLWw1dGyHZWbFKjciIhkgNgI2DzeLDrRZ81tXgFQty80eBp8C1ib7y4dCo9iyC/b2XXKLHCd6xTjzdZV8PVytzhZ9qFykwKVGxGRDJQYBzt/gTUj4OJBc5ubF9ToBg0HQd4y1ua7C/GJDr5YdoAxKw9jGFA8Tw6+6BJE7RJ5rI6WLajcpEDlRkQkEzgcsH8BrP4STm02t9nsUKmNOetxsdqWxrsbG45c5Plfd3DqyjXsNhh4X1kGNS2Hh5trDKbOqlRuUqByIyKSiQwDjq8xS84/58opHgzBA6FCC7A73zwykbEJvD1nNzO3nQIgqFgAX3SpQen8vhYnc10qNylQuRERscjZXbBuJIROA0eiuS1PGQh+BoK6gWcOa/Olwbwdp3ltViiRsYn4eLjxxsOV6VovULeMZwCVmxSo3IiIWCzyNGz4FraMNwciA/jkMQcf130C/Apamy+VzkRc43+/7mDt4YuAuT7Vh49UJ5/Wp0pXKjcpULkREcki4qJh2yRY/w1cMWcExs0Tqnc2L1k50RpWDofBD2uO8vHi/dfXp/Lko0eq07SScxW1rEzlJgUqNyIiWYwjyVzWYd1IOLnp7+1lm5klp3QTp1neYe+ZSIZM3c7+c1EAdK9fnNdaVSKHp24Zv1sqNylQuRERycLCNsC6r2HvfOD6r6cClaH+U+YZHQ8fS+PdidiEJD5dsp/vVx8FoHS+nHzRpQZBgbmsDebkVG5SoHIjIuIELh0xJwTc9jMkxJjbfPJAnT5Qtx/4F7E23x1Yc+gC//t1B2cjY3G32xjctBxPNymDu24ZTxOVmxSo3IiIOJFrV2DbT7BhLESEmdvs7lC5nbnEQxafL+fK1Xhem7WLBaFnAKhdIjdfdK5B8bzOd2eY1VRuUqByIyLihJISYf9C82xO2Nq/txeray7vUKkNuHlYly8FhmEwa9sp3pyzm+i4RHJ6uvF2myp0rF1Mt4yngspNClRuRESc3OntsGEMhE4HR4K5zb+oeSt5rd6QM5+l8W7nxKWr/O/XHWw8dgmAFlUL8UH7auTO6fwrqGcGlZsUqNyIiLiIqHPmauSbx0HMeXObmydU6QD1+mfJS1ZJDoNvVx3m898OkOgwKODnxaedgrinfH6ro2V5KjcpULkREXExCbGwe6Y5MeCZ7X9vL1LTnBSwaocsd5dV6MkIhvyyjcPnzcHSjzUsySstKuLt4XxLUWQWlZsUqNyIiLgow4BTW2Djd2bZSYo3t/vkhpo9zctWuUtaGvGfrsUnMXzRXn5cZ05gWK6AL18+WoMqRQIsTpY1qdykQOVGRCQbiLkAW380L1tFnLi+0QblHoR6T0CZpmDPGrdk/7EvnBen7+RCdBwebjZeeLAC/RqXxs2uwcb/pHKTApUbEZFsxJEEBxabZ3OO/PH39jyloU5fqNndPLNjsYvRcbwyM5Sle84BUL9UHj7rHESx3Lpl/C8qNylQuRERyaYuHIJN38P2yRB3fcFOdx+o1tE8m1M4yNJ4hmHw6+YTvDNvD1fjk/D1cuet1pV1y/h1KjcpULkREcnm4mNg569m0Tm36+/txeqZd1lVbgPu1q3offxiDM//uoMtxy8D0LxKQT5oX4282XyVcZWbFKjciIgIYA5ADltnXrLaOxccieb2nPnN+XLq9IGAYpZES3IYjFl5mC+XHSAhySCfrxcfPVItW68yrnKTApUbERG5SdRZ2DIRtoyHKHOpBGx2KNccaj8G5R4Ae+bfpr3rVATP/bKdg+HRAHStF8jrrSqT0yv7rTKucpMClRsREbmtpATYt8C8ZHXsz7+3+xc1byev1TPTz+b8tcr4uDVHMQwonicHn3cOok7JPJmaw2oqNylQuRERkTty/gBsnWgOQL5mLpmAzQ5lHzAvWZV9ANwy7wzKusMXeWHaDk5duYbdBk/dW4Yhzcrj6Z41bmnPaCo3KVC5ERGRVEmIhX3zYcuEG8/m+BWBmj3Mszm5imdKlMjYBN6eu5uZW08BULmwP190qUGFQn6Z8vlWUrlJgcqNiIik2YVDsHWCeTbn6sXrG21Qtpl5Nqdc80w5m7Mo9Ayvzgrl8tUEPN3tvNS8Ao83KoXdhSf+U7lJgcqNiIjctcS4v8/mHF3193bfQuaZnJo9IXeJDI0QHhXLy9N38sd+c9HQBqXz8Gkn1534T+UmBSo3IiKSri4eNsfmbPsZrl64vtEGZZuad1qVfwjcPDLkow3DYMrGE7y3wJz4z8/LnbfbVKFDraIuN/Gfyk0KVG5ERCRDJMbD/gXm2ZwjK/7enjM/BD0KNXpAgYoZ8tHHLsTw/K/b2Rp2BYCHqhTigw7VyJPTM0M+zwoqNylQuRERkQx36Yi5cOe2nyEm/O/txeqag5CrdADv9P0dlJjk4NtVR/hi6QESHebEf590rM59FQuk6+dYJTW/vy29f2z48OHUrVsXPz8/ChQoQLt27di/f3+Kr5kwYQI2m+2Gh7e3dyYlFhERuQN5SkOzt+H5PfDoFKjQCmxucHITzBsMn5aHWU/BsdXmTMnpwN3NzoD7yjJ7QCPKFfDlQnQcfSZs4tVZocTEJabLZzgLS8vNypUrGTBgAOvXr2fp0qUkJCTw4IMPEhMTk+Lr/P39OXPmTPLj+PHjmZRYREQkFdw8oGJL6DoZnt8LD7wL+cpD4jXYMQUmtIIRNWHVJxBxKl0+smrRAOYNCqFvSCkAJm8Io+WIP5PXqsoOstRlqfPnz1OgQAFWrlzJPffcc8t9JkyYwJAhQ7hy5UqaPkOXpURExFKGASc3w7afYNdMiI8yt9vsUOZ+87JVhZbpsnjn2sMXeOHXHZyOiMVug6eblGFwU+ec+M9pLkv9W0SEuQR9njwpTykdHR1NiRIlCAwMpG3btuzevfu2+8bFxREZGXnDQ0RExDI2GwTWhTYj4IX90G4MlAgBwwGHlsG0x+CzCrDwJTi97a4uWzUsk49FQ+6hQ82iOAwY9cdh2n+zhgPnotLv58mCssyZG4fDQZs2bbhy5QqrV6++7X7r1q3j4MGDVK9enYiICD799FNWrVrF7t27KVbs5vU+3n77bd55552btuvMjYiIZCkXD5uTA26fDFGn/96ev5J5t1X1zuBfJM1vv/D6xH9Xrk/898KD5ekbUho3J5n4zynvlnr66adZtGgRq1evvmVJuZ2EhAQqVapE165deffdd2/6flxcHHFxccnPIyMjCQwMVLkREZGsyZEEh383S86+BZD01+8wG5RuAkFdodLD4Jkz1W8dHhnLKzND+X2feQdXvZLmxH/F82b9if+crtwMHDiQOXPmsGrVKkqVKpXq13fq1Al3d3emTJnyn/tqzI2IiDiNa1dgzxzYMRXC1v693dMXKrc1z+iUCAH7nY8yMQyDXzefYNi8PcTEJ5HD043XWlWiW73iWXriP6cpN4ZhMGjQIGbNmsWKFSsoV65cqt8jKSmJKlWq0LJlSz7//PP/3F/lRkREnNKlo7DzF/Muq8vH/t4eEGhesgrqCvnu/PfoiUtXeWHaDjYcNVc8v7d8fj56pDqFArLm9CpOU26eeeYZJk+ezJw5c6hQoULy9oCAAHx8fADo1asXRYsWZfjw4QAMGzaMBg0aULZsWa5cucInn3zC7Nmz2bJlC5UrV/7Pz1S5ERERp2YYcGKDedlq92yIi/j7e0XrmGdzqj4COVK+OQfA4TD4Yc1RPl6yn/hEB/7e7rzbriptgopkubM4TlNubnfgxo8fz2OPPQZAkyZNKFmyJBMmTADgueeeY+bMmZw9e5bcuXNTu3Zt3nvvPWrWrHlHn6lyIyIiLiPhGuxfZF62OrQMjCRzu90DKjwE1R+Fcg/8523lh8KjeP7XHew8aRalltUK8V67rLV8g9OUGyuo3IiIiEuKDofQaeZlq7Ohf2/3zmWOz6neGYo3vO34nIQkB6NXHGbE8oPJyzd82KEazSoXzJz8/0HlJgUqNyIi4vLO7jJLzq4ZEHXm7+3+Rc1LVtU6QaFq5pw7/7LrVATP/bKdg+HRAHSqXYw3W1fGzztjVja/Uyo3KVC5ERGRbMORZK5fFfor7Jl34/ic/BXNklOtI+QuecPLYhOS+GLpAcb+eQTDgKK5fPikY3Uals2Xufn/QeUmBSo3IiKSLSXEwsHfzKJzYAkkxf/9vcD6ZtGp0h5y/l1gNh69xAvTdhB26SoAjzUsycsPVcTH0y2z06vcpETlRkREsr1rV2DvPLPoHP0TuF4F7O7m+lbVOpsLfnrmJCYukfcX7mXyhjAASufLyWedg6hZPHemRla5SYHKjYiIyD9EnjHH5oROgzPb/97ukQPKNzfH6JR9gBVHInl5xk7ORcZht8EzTcrybNNymbYIp8pNClRuREREbuP8AbPkhE6Dy0f/3u7pBxVbEVOuDW/sys/MHecBqFzYn8+7BFGxUMb/PlW5SYHKjYiIyH8wDHNF8l0zzIkCI0/+/T3vXJwo1Iz3j1di6bXy2N3cee6B8jx5T5kMXYRT5SYFKjciIiKp4HDAyY2waybsngUx4cnfirTnYnZ8XeYnNSCpWH0+7VKLUvlSv6DnnVC5SYHKjYiISBo5kuD4GrPo7JkD1y4lf+uskZslRjD5g7vy0IMPY3dL37E4KjcpULkRERFJB0kJcHQl7JqJY8887PGRyd/a7V2Lyi//nq7rU6Xm93fmDHEWERER1+LmAWWbQbtvsL90CEeXyRwu1IIYw4urBWpauvCmu2WfLCIiIq7B3Qt7pVaUqdSKsLPnqelnbb1QuREREZF0U7xQfqsj6LKUiIiIuBaVGxEREXEpKjciIiLiUlRuRERExKWo3IiIiIhLUbkRERERl6JyIyIiIi5F5UZERERcisqNiIiIuBSVGxEREXEpKjciIiLiUlRuRERExKWo3IiIiIhLyXarghuGAUBkZKTFSURERORO/fV7+6/f4ynJduUmKioKgMDAQIuTiIiISGpFRUUREBCQ4j42404qkAtxOBycPn0aPz8/bDab1XHSVWRkJIGBgZw4cQJ/f3+r42QZOi63puNyazouN9MxuTUdl1vLqONiGAZRUVEUKVIEuz3lUTXZ7syN3W6nWLFiVsfIUP7+/vqLdgs6Lrem43JrOi430zG5NR2XW8uI4/JfZ2z+ogHFIiIi4lJUbkRERMSlqNy4EC8vL9566y28vLysjpKl6Ljcmo7Lrem43EzH5NZ0XG4tKxyXbDegWERERFybztyIiIiIS1G5EREREZeiciMiIiIuReVGREREXIrKjZMZNWoUJUuWxNvbm/r167Nx48Y7et3UqVOx2Wy0a9cuYwNaJDXHZcKECdhsthse3t7emZg286T2z8uVK1cYMGAAhQsXxsvLi/Lly7Nw4cJMSps5UnNMmjRpctOfFZvNRqtWrTIxceZI7Z+VL7/8kgoVKuDj40NgYCDPPfccsbGxmZQ286TmuCQkJDBs2DDKlCmDt7c3QUFBLF68OBPTZrxVq1bRunVrihQpgs1mY/bs2f/5mhUrVlCrVi28vLwoW7YsEyZMyPCcGOI0pk6danh6eho//PCDsXv3buOJJ54wcuXKZZw7dy7F1x09etQoWrSo0bhxY6Nt27aZEzYTpfa4jB8/3vD39zfOnDmT/Dh79mwmp854qT0ucXFxRp06dYyWLVsaq1evNo4ePWqsWLHC2L59eyYnzzipPSYXL1684c/Jrl27DDc3N2P8+PGZGzyDpfa4/Pzzz4aXl5fx888/G0ePHjWWLFliFC5c2HjuuecyOXnGSu1xeemll4wiRYoYCxYsMA4fPmx88803hre3t7F169ZMTp5xFi5caLz22mvGzJkzDcCYNWtWivsfOXLEyJEjh/H8888be/bsMb7++mvDzc3NWLx4cYbmVLlxIvXq1TMGDBiQ/DwpKckoUqSIMXz48Nu+JjEx0WjYsKHx/fffG71793bJcpPa4zJ+/HgjICAgk9JZJ7XHZfTo0Ubp0qWN+Pj4zIqY6dLyd+ifvvjiC8PPz8+Ijo7OqIiWSO1xGTBggHH//fffsO355583GjVqlKE5M1tqj0vhwoWNkSNH3rCtQ4cORvfu3TM0p1XupNy89NJLRpUqVW7Y1qVLF6N58+YZmMwwdFnKScTHx7NlyxaaNWuWvM1ut9OsWTPWrVt329cNGzaMAgUK0Ldv38yImenSelyio6MpUaIEgYGBtG3blt27d2dG3EyTluMyd+5cgoODGTBgAAULFqRq1ap88MEHJCUlZVbsDJXWPyv/NG7cOB599FFy5syZUTEzXVqOS8OGDdmyZUvyJZojR46wcOFCWrZsmSmZM0NajktcXNxNl7h9fHxYvXp1hmbNytatW3fDMQRo3rz5Hf+dSyuVGydx4cIFkpKSKFiw4A3bCxYsyNmzZ2/5mtWrVzNu3Di+++67zIhoibQclwoVKvDDDz8wZ84cJk2ahMPhoGHDhpw8eTIzImeKtByXI0eOMH36dJKSkli4cCFvvPEGn332Ge+9915mRM5waTkm/7Rx40Z27dpFv379MiqiJdJyXLp168awYcMICQnBw8ODMmXK0KRJE1599dXMiJwp0nJcmjdvzueff87BgwdxOBwsXbqUmTNncubMmcyInCWdPXv2lscwMjKSa9euZdjnqty4qKioKHr27Ml3331Hvnz5rI6TpQQHB9OrVy9q1KjBvffey8yZM8mfPz/ffvut1dEs5XA4KFCgAGPHjqV27dp06dKF1157jTFjxlgdLUsYN24c1apVo169elZHsdyKFSv44IMP+Oabb9i6dSszZ85kwYIFvPvuu1ZHs9RXX31FuXLlqFixIp6engwcOJA+ffpgt+tXbWZztzqA3Jl8+fLh5ubGuXPnbth+7tw5ChUqdNP+hw8f5tixY7Ru3Tp5m8PhAMDd3Z39+/dTpkyZjA2dCVJ7XG7Fw8ODmjVrcujQoYyIaIm0HJfChQvj4eGBm5tb8rZKlSpx9uxZ4uPj8fT0zNDMGe1u/qzExMQwdepUhg0blpERLZGW4/LGG2/Qs2fP5LNY1apVIyYmhv79+/Paa6+5xC/ztByX/PnzM3v2bGJjY7l48SJFihThlVdeoXTp0pkROUsqVKjQLY+hv78/Pj4+Gfa5zv8nMJvw9PSkdu3aLF++PHmbw+Fg+fLlBAcH37R/xYoVCQ0NZfv27cmPNm3acN9997F9+3YCAwMzM36GSe1xuZWkpCRCQ0MpXLhwRsXMdGk5Lo0aNeLQoUPJJRjgwIEDFC5c2OmLDdzdn5Vp06YRFxdHjx49MjpmpkvLcbl69epNBeavUmy4yHKFd/Pnxdvbm6JFi5KYmMiMGTNo27ZtRsfNsoKDg284hgBLly6943+f0yxDhytLupo6darh5eVlTJgwwdizZ4/Rv39/I1euXMm3Mffs2dN45ZVXbvt6V71bKrXH5Z133jGWLFliHD582NiyZYvx6KOPGt7e3sbu3but+hEyRGqPS1hYmOHn52cMHDjQ2L9/vzF//nyjQIECxnvvvWfVj5Du0vp3KCQkxOjSpUtmx800qT0ub731luHn52dMmTLFOHLkiPHbb78ZZcqUMTp37mzVj5AhUntc1q9fb8yYMcM4fPiwsWrVKuP+++83SpUqZVy+fNminyD9RUVFGdu2bTO2bdtmAMbnn39ubNu2zTh+/LhhGIbxyiuvGD179kze/69bwV988UVj7969xqhRo3QruNzs66+/NooXL254enoa9erVM9avX5/8vXvvvdfo3bv3bV/rquXGMFJ3XIYMGZK8b8GCBY2WLVu61DwU/5TaPy9r16416tevb3h5eRmlS5c23n//fSMxMTGTU2es1B6Tffv2GYDx22+/ZXLSzJWa45KQkGC8/fbbRpkyZQxvb28jMDDQeOaZZ1zql/hfUnNcVqxYYVSqVMnw8vIy8ubNa/Ts2dM4deqUBakzzh9//GEANz3+Og69e/c27r333pteU6NGDcPT09MoXbp0pswTZTMMFzmHKCIiIoLG3IiIiIiLUbkRERERl6JyIyIiIi5F5UZERERcisqNiIiIuBSVGxEREXEpKjciIiLiUlRuRERExKWo3IiIiIhLUbkRERERl6JyIyICNGnShCFDhtzVexiGQf/+/cmTJw82m43t27enSzYRSR2VGxHJcvr06cPrr79udYxUW7x4MRMmTGD+/PmcOXOGqlWrWh1JJFtytzqAiMg/JSUlMX/+fBYsWGB1lFQ7fPgwhQsXpmHDhrfdJz4+Hk9Pz0xMJZL96MyNSDYxZcoUfHx8OHPmTPK2Pn36UL16dSIiItL8vsWKFeObb765YdvatWvJkSMHx48fT/X7rV27Fg8PD+rWrXvL7zdp0oRBgwYxZMgQcufOTcGCBfnuu++IiYmhT58++Pn5UbZsWRYtWnTD6+Li4nj22WcpUKAA3t7ehISEsGnTptvmcDgcDB8+nFKlSuHj40NQUBDTp0+/7f6PPfYYgwYNIiwsDJvNRsmSJZPzDhw4kCFDhpAvXz6aN2/O4sWLCQkJIVeuXOTNm5eHH36Yw4cP3/T5H3/8MWXLlsXLy4vixYvz/vvv3+FRFMneVG5EsolHH32U8uXL88EHHwDw1ltvsWzZMhYtWkRAQECa37d+/fo3lATDMBgyZAjPPfccJUqUSPX7zZ07l9atW2Oz2W67z8SJE8mXLx8bN25k0KBBPP3003Tq1ImGDRuydetWHnzwQXr27MnVq1eTX/PSSy8xY8YMJk6cyNatWylbtizNmzfn0qVLt/yM4cOH8+OPPzJmzBh2797Nc889R48ePVi5cuUt9//qq68YNmwYxYoV48yZMzcck4kTJ+Lp6cmaNWsYM2YMMTExPP/882zevJnly5djt9tp3749Docj+TVDhw7lww8/5I033mDPnj1MnjyZggULpvZwimRPhohkG/PmzTO8vLyM9957z8idO7exa9eu5O+1a9fOyJUrl/HII4+k6j0//vhjo0qVKsnPJ06caBQqVMiIiopK0/uWK1fOmD9//m2/f++99xohISHJzxMTE42cOXMaPXv2TN525swZAzDWrVtnGIZhREdHGx4eHsbPP/+cvE98fLxRpEgR4+OPP05+38GDBxuGYRixsbFGjhw5jLVr197w2X379jW6du1622xffPGFUaJEiZvy1qxZM8Wf+fz58wZghIaGGoZhGJGRkYaXl5fx3Xffpfg6Ebk1nbkRyUYefvhhKleuzLBhw5g1axZVqlRJ/t7gwYP58ccfU/2eDRo0YO/evURHRxMTE8Orr77Ke++9h6+vb6rfd+/evZw+fZqmTZumuF/16tWTv3ZzcyNv3rxUq1YtedtfZzjCw8MBcyxMQkICjRo1St7Hw8ODevXqsXfv3pve/9ChQ1y9epUHHngAX1/f5MePP/540+WjO1G7du0bnh88eJCuXbtSunRp/P39ky9hhYWFAeZxiIuL+8/jICK3pgHFItnI4sWL2bdvH0lJSTdd4mjSpAkrVqxI9XvWrl0bu93O1q1bWbZsGfnz56dPnz5pet+5c+fywAMP4O3tneJ+Hh4eNzy32Ww3bPvrktY/L/OkRnR0NAALFiygaNGiN3zPy8sr1e+XM2fOG563bt2aEiVK8N1331GkSBEcDgdVq1YlPj4eAB8fnzTlFhGTztyIZBNbt26lc+fOjBs3jqZNm/LGG2+ky/vmyJGDatWqMWPGDD799FO++OIL7Pa0/dMyZ84c2rZtmy65/qlMmTLJY17+kpCQwKZNm6hcufJN+1euXBkvLy/CwsIoW7bsDY/AwMC7ynLx4kX279/P66+/TtOmTalUqRKXL1++YZ9y5crh4+PD8uXL7+qzRLIrnbkRyQaOHTtGq1atePXVV5MvhwQHB7N161Zq1ap11+/foEEDvv76a9q2bUuTJk3S9B7h4eFs3ryZuXPn3nWef8uZMydPP/00L774Inny5KF48eJ8/PHHXL16lb59+960v5+fHy+88ALPPfccDoeDkJAQIiIiWLNmDf7+/vTu3TvNWXLnzk3evHkZO3YshQsXJiwsjFdeeeWGfby9vXn55Zd56aWX8PT0pFGjRpw/f57du3cn5x05ciSzZs1SARK5BZUbERd36dIlHnroIdq2bZv8S7R+/fq0aNGCV199lcWLF//ne0yYMIE+ffpgGMYtvx8UFISHhweffPJJmnPOmzePevXqkS9fvjS/R0o+/PBDHA4HPXv2JCoqijp16rBkyRJy5859y/3fffdd8ufPz/Dhwzly5Ai5cuWiVq1avPrqq3eVw263M3XqVJ599lmqVq1KhQoVGDFixE2l8I033sDd3Z0333yT06dPU7hwYZ566qnk71+4cCFN439EsgObcbt/rUQk21mxYgUjR468aT6Xt956i5UrV9527Mx9991HrVq1+Oyzz1L1vv/Upk0bQkJCeOmll9KcX0QEdOZGRK5r1qwZO3bsICYmhmLFijFt2jSCg4MBWLRoESNHjrxhf4fDwfnz5xk3bhwHDx5kzpw5qX7ffwoJCaFr167p/4OJSLajMzcikiYrVqzg/vvvp2LFiowfP5769etbHUlEBFC5ERERERejW8FFRETEpajciIiIiEtRuRERERGXonIjIiIiLkXlRkRERFyKyo2IiIi4FJUbERERcSkqNyIiIuJSVG5ERETEpajciIiIiEv5P7W36PO67moxAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "Tc_K = [190.564, 154.581]\n", "pc_Pa = [4599200, 5042800]\n", "acentric = [0.011, 0.022]\n", "model = teqp.canonical_PR(Tc_K, pc_Pa, acentric)\n", "model1 = teqp.canonical_PR([Tc_K[0]], [pc_Pa[0]], [acentric[0]])\n", "T = 170.0 # [K] # Note: above Tc of the second component\n", "rhoL0, rhoV0 = model1.superanc_rhoLV(T) # start off at pure of the first component\n", "j = model.trace_VLE_isotherm_binary(T, np.array([rhoL0, 0]), np.array([rhoV0, 0]))\n", "df = pandas.DataFrame(j) # Now as a data frame\n", "plt.plot(df['xL_0 / mole frac.'], df['pL / Pa']/1e6)\n", "plt.plot(df['xV_0 / mole frac.'], df['pL / Pa']/1e6)\n", "plt.gca().set(xlabel='$x_1,y_1$ / mole frac.', ylabel='p / MPa')\n", "plt.show()" ] }, { "cell_type": "markdown", "id": "936c5f55", "metadata": {}, "source": [ "As of version 0.10.0, isobar tracing has been added to ``teqp``. It operates in fundamentally the same fashion as the isotherm tracing and the same recommendations about starting at a pure fluid apply" ] }, { "cell_type": "raw", "id": "81eacaf0", "metadata": { "raw_mimetype": "text/restructuredtext" }, "source": [ "The tracer function class is here: :py:meth:`~teqp.teqp.AbstractModel.trace_VLE_isobar_binary`" ] }, { "cell_type": "code", "execution_count": 7, "id": "109bc65b", "metadata": { "execution": { "iopub.execute_input": "2024-12-12T18:10:45.009816Z", "iopub.status.busy": "2024-12-12T18:10:45.009484Z", "iopub.status.idle": "2024-12-12T18:10:45.100496Z", "shell.execute_reply": "2024-12-12T18:10:45.099962Z" } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjsAAAG0CAYAAADU2ObLAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjkuMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/GU6VOAAAACXBIWXMAAA9hAAAPYQGoP6dpAABoXklEQVR4nO3dd3QU5dvG8e+mE0ghCSEEQu9VagBRRDpKt1BFQKyggq8FrGABe8HCTwQBBbGCgApSBFR6Cb2FXgIBQjqpO+8fo9EoYAJJZndzfc7Zc/LMzs7eO4bs5cxTbIZhGIiIiIi4KDerCxAREREpTAo7IiIi4tIUdkRERMSlKeyIiIiIS1PYEREREZemsCMiIiIuTWFHREREXJqH1QU4ArvdzqlTp/Dz88Nms1ldjoiIiOSBYRgkJSURHh6Om9vlr98o7ACnTp0iIiLC6jJERETkKhw/fpwKFSpc9nmFHcDPzw8wT5a/v7/F1YiIiEheJCYmEhERkfM9fjkKO5Bz68rf319hR0RExMn8VxcUdVAWERERl6awIyIiIi5NYUdERERcmsKOiIiIuDSFHREREXFpCjsiIiLi0hR2RERExKUp7IiIiIhLU9gRERERl6awIyIiIi7N0rCzevVqunfvTnh4ODabjfnz5+d63mazXfLx+uuv5+wTFxfHwIED8ff3JzAwkOHDh5OcnFzEn0REREQclaVhJyUlhUaNGvHBBx9c8vmYmJhcj+nTp2Oz2ejbt2/OPgMHDmTXrl0sXbqURYsWsXr1au69996i+ggiIiLi4GyGYRhWFwHmVZx58+bRq1evy+7Tq1cvkpKSWL58OQB79uyhbt26bNy4kWbNmgGwePFiunXrxokTJwgPD8/TeycmJhIQEEBCQoIWAhURESlACcmpnIr6mTptehX4sfP6/e00fXbOnDnDDz/8wPDhw3O2rV27lsDAwJygA9ChQwfc3NxYv379ZY+Vnp5OYmJiroeIiIgUILudjKivSXm7GXWWDeH3X1dYVorThJ2ZM2fi5+dHnz59cradPn2a0NDQXPt5eHgQFBTE6dOnL3usiRMnEhAQkPOIiIgotLpFRESKFcOA6GXYP74Jr/n3EJ59kvP4E+EeZ1lJThN2pk+fzsCBA/Hx8bnmY40dO5aEhIScx/HjxwugQhERkWLuxCaY2R0+74vb6W0kGSV4z347Rwf+TsXWt1lWlodl75wPv/76K/v27ePLL7/MtT0sLIzY2Nhc27KysoiLiyMsLOyyx/P29sbb27tQahURESl2zu6D5RNg7yIAsmyezMjswMdGL14f0p4mNcpYWp5ThJ1p06bRtGlTGjVqlGt7q1atiI+PZ/PmzTRt2hSAFStWYLfbiYyMtKJUERGR4iPhBKycCFFzwLBj2NzYHtyNB0504rQthMn9m9C2prVBBywOO8nJyURHR+e0Dx8+TFRUFEFBQVSsWBEwe1p//fXXvPnmm/96fZ06dejSpQsjRoxgypQpZGZmMnLkSPr165fnkVgiIiKST6lx8OubsGEqZKeb22rfyhelhjDut0wAJvVuwC0Ny1lY5F8sDTubNm2iXbt2Oe0xY8YAMGTIEGbMmAHA3LlzMQyD/v37X/IYs2fPZuTIkbRv3x43Nzf69u3Le++9V+i1i4iIFDsZKbDuQ/j9PUj/YyRzpTbQ4QU+O1GGZ7/fBcDT3erQr0VFCwvNzWHm2bGS5tkRERG5guxM2PoZrJwEyWfMbWXrQ4cXoHoH5ked4tEvowAYdXN1HutUq0jKyuv3t1P02RERERELGAbsng8rXoLzf3Q7CawENz8D9W8DNzeW7j7DY19vA+Du1pUZ07GmdfVehsKOiIiI/Nvh1bD0eTi1xWz7hkDbJ6DpUPDwAmBN9Dkemr2FbLtBnyblee7WuthsNguLvjSFHREREflLzHZY9gIcNJdmwqsUtBoJrUeCt1/ObluOXeCeWZvIyLbTpV4Yr/VtiJub4wUdUNgRERERgAtHYMXLsOMrs+3mCc2Gwo1PQKncw8d3n0rk7ukbSM3I5oYaIbzb/zo83B13nmKFHRERkeIs5RysfgM2fgJ2c9g49W+Dm5+GoKr/2v3Q2WTumr6exLQsmlYqzf8GN8Xbw72Ii84fhR0REZHiKCMF1n4Iv78LGUnmtqo3mSOswhtf8iUnLqQy6JP1nEvOoG45f6bf3RxfL8ePEo5foYiIiBScSw0jD2sIHcdDtZsv+7LYxDQGfbKeUwlpVC1TklnDWxBQwrOIir42CjsiIiLFgWGYa1ctGw/nD5jbAitB++egXh9wu3yfmwspGQyetoEj51OpULoEs++JJKSU86wxqbAjIiLi6o6tg6XPwfH1Zts32Ox43GwoeFw5tCSlZTLk0w3sO5NEqJ83c+5pSbmAEkVQdMFR2BEREXFVZ/ebw8j3/WC2PUqYQ8hbPww+/71iwMWMbIbP3MT2EwmU9vVk9j2RVAz2LdyaC4HCjoiIiKtJOm2uRr7lMzCyweYOTQZD26fAP2+Lc6ZnZXPf55vZcDgOP28PZg2LpEZZv/9+oQNS2BEREXEV6UnmIp1r34fMVHNbrVugw/NQJu/rVWVl23n4i62s3n+WEp7uTB/anAYVAgqp6MKnsCMiIuLssjNh8wxY9SqknDW3VWgOHV+ESq3ydSi73eDxb7azZNcZvNzdmHpXM5pXDir4mouQwo6IiIizMgzYs9DslxN30NwWVM28klOnB+RznSrDMHh6/k7mbT2Jh5uNDwY2oU2NkIKvu4gp7IiIiDijY+vg52fhxAaz7RsCNz0FTe8G9/zPf2MYBi/9sIcvNhzDZoO377yOjnXLFmzNFlHYERERcSbnomHZ8+acOQCevn8s1DkqTyOsLuftpfuZ9tthAF7t05DujcILolqHoLAjIiLiDJLPwqpJsOnTP0ZYuUHjQXDTuDyPsLqcj1Ye5L0V0QC80L0udzSPKIiKHYbCjoiIiCPLSIV1H8Bvf1vDqkZnc3mH0DrXfPgZvx/m1cV7AXiyS23uvr7KNR/T0SjsiIiIOCJ7Nmz7Ala8DEmnzG3lGkGnl6DKjQXyFl9uPMYLC3cD8PDN1XngpmoFclxHo7AjIiLiaKKXm8s7nNlptgMqmmtY1e97xTWs8uP7qJM89d0OAO5pU4XRHWsWyHEdkcKOiIiIozi90ww5B5ebbZ8AuOH/oMW94OlTYG+zeGcMY77ahmHAoJYVefqWOtjyOUzdmSjsiIiIWC0xBn55CbbOBgxw84QWI+DGx8G3YCf0+2VvLKO+2Eq23eC2phWY0KO+SwcdUNgRERGxTnoyrHkP1kz+a3mHer3NW1ZBVQv87X6PPsd9n28mM9vg1obleLVvQ9zcXDvogMKOiIhI0bNnw9bPzM7HKbHmtohIs/NxRItCecv1h84zfOZGMrLsdKxblrfvvA73YhB0QGFHRESkaB1YBkufhVhzFBSlq5jDyK9ieYe82nLsAsNmbCQt007bmmV4f0BjPN0LpqOzM1DYERERKQqnd5oh5+AKs+0TCG2fhOb3gIdXob3tjhMJDJm+gZSMbFpXC+Z/g5vi7eFeaO/niBR2REREClPSaVjxEmz9nJzOx5H3wY3/ByVKF+pb74lJZPD09SSlZdG8cmk+GdIMH8/iFXRAYUdERKRwZKTAmvfh93chM8XcVrcXdHgBggp/luIDZ5IY9Ml64lMzaRQRyPS7m+PrVTy/9ovnpxYRESksdjtsnwvLX/xr5uMKzaHTy1AxskhKOHQ2mQGfrOd8Sgb1y/sza2gL/HzyvxK6q1DYERERKSiHV8OSp+H0drMdUBE6vgD1+hRa5+N/OnY+lQFT13M2KZ3aYX58NiySAN/iG3RAYUdEROTanYs2Ox/v+9Fse/vDDY9B5P0FOvPxfzkZf5H+U9dxOjGNGqGl+PyeSEqXLLzOz85CYUdERORqpcbBqtdg41SwZ4HNHZoNhZvGQsmQIi0lJuEi/T9ex8n4i1QJKcnseyIJKeVdpDU4KoUdERGR/MrKgE3TYOUkSIs3t9XoZE4KWKZWkZcTm5jGgKnrORaXSsUgX+aMiCTUv+iuKDk6hR0REZG8MgzY9xP8/AzEHTS3hdaDzi9BtZstKelsUjr9p67j8LkUygeWYM6ISMoFlLCkFkelsCMiIpIXMdvh56fNTsgAJctAu6ehyV3gZs3cNeeT0xn4yToOnk2hXIAPc+9tSYXSvpbU4sgUdkRERK4k6QysePGvSQHdvaHVg9BmDPj4W1ZWfGoGg6ZtYP+ZZMr6e/PFiJZEBCnoXIrCjoiIyKVkpsG6D+DXtyAj2dxWr485KWDpSpaWlpCayaBp69kTk0hIKW/mjGhJ5ZCSltbkyBR2RERE/s4wYPd8+Pk5SDhmbivfFDpPLLJJAa8k4WImg6evZ+fJRIJLevHFiEiqlSlldVkOTWFHRETkTye3wJJxcGyt2fYLN1ckr38buFm/SnhiWiZ3Td/A9hMJBJX0Ys6IltQo62d1WQ7P0v9yq1evpnv37oSHh2Oz2Zg/f/6/9tmzZw89evQgICCAkiVL0rx5c44dO5bz/E033YTNZsv1uP/++4vwU4iIiNNLjIF5D8DUdmbQ8ShhzpUzahM0vMMhgk5SWiZDpm9g2/F4An09mX1PJLXCFHTywtIrOykpKTRq1Ihhw4bRp0+ffz1/8OBB2rRpw/Dhwxk/fjz+/v7s2rULH5/ccweMGDGCCRMm5LR9fdVBS0RE8iDzorlY529v/7VYZ8M7of3zEFDe2tr+Jjk9i6GfbmTrsXgCSphBp0456zpHOxtLw07Xrl3p2rXrZZ9/+umn6datG6+99lrOtmrVqv1rP19fX8LCwvL8vunp6aSnp+e0ExMT8/xaERFxAYYBO7+FZS9AwnFzW4Xm0GUSVGhmaWn/lJqRxbBPN7Lp6AX8fTyYfU8k9cIDrC7LqVh/Xe4y7HY7P/zwAzVr1qRz586EhoYSGRl5yVtds2fPJiQkhPr16zN27FhSU1OveOyJEycSEBCQ84iIiCikTyEiIg7n1FaY3gW+HW4GHf8K0HcaDF/qkEFn6Kcb2XAkDj8fDz4bHkn98go6+WUzDMOwuggAm83GvHnz6NWrFwCnT5+mXLly+Pr68tJLL9GuXTsWL17MuHHj+OWXX2jbti0AH3/8MZUqVSI8PJzt27fz5JNP0qJFC7777rvLvtelruxERESQkJCAv78uC4qIuKSkM7B8AkTNBgzw9IXrH4XWo8DL8bo/pGZkMWzGRtYdisPP24NZw1vQuGJpq8tyKImJiQQEBPzn97fDjsay2+0A9OzZk9GjRwNw3XXXsWbNGqZMmZITdu69996c1zRo0IBy5crRvn17Dh48eMlbXgDe3t54e2txNBGRYiErHdZ9BKvfgIwkc1uDO8z5chyoX87fXczIZviMTaw7FEcpbw9mKuhcE4cNOyEhIXh4eFC3bt1c2+vUqcNvv/122ddFRppzIERHR1827IiISDFgGLDvR1jyNFw4bG4LbwJdX4WIFtbWdgUXM7IZNmMjaw+dp9QfV3SaKOhcE4cNO15eXjRv3px9+/bl2r5//34qVbr8zJVRUVEAlCtXrjDLExERR3ZmNywZC4dWmu1SZc0rOQ37OcQw8su5mJHN8Jl/BZ2ZwxR0CoKlYSc5OZno6Oic9uHDh4mKiiIoKIiKFSvy+OOPc+edd3LjjTfm9NlZuHAhK1euBMyh6XPmzKFbt24EBwezfft2Ro8ezY033kjDhg0t+lQiImKZ1Dj45RXYNB2MbHD3glYj4YYx4O3Yc9JczMjmnlkbWXPwPCW93Jk5rDlNKynoFARLOyivXLmSdu3a/Wv7kCFDmDFjBgDTp09n4sSJnDhxglq1ajF+/Hh69uwJwPHjxxk0aBA7d+4kJSWFiIgIevfuzTPPPJOvjsZ57eAkIiIOKjvLDDi/vAxp8ea2Ot2h44sQVMXS0vLiYkY2I2Zt4rfoc5T0cmfW8BY0rRRkdVkOL6/f3w4zGstKCjsiIk7s0CpY/BTE7jbbZetDl4lQ5UZr68qjP6/o/B795xWdFjSrrKCTF04/GktEROSKLhyFn5+BPQvMdonScPMz0ORucHeOr7fUjCyGz9jE2kNm0JmhoFMonOO3QURE5E8ZqfD7O/D7u5CVBjY3aH6PuZaVr/MEhb/Po2N2Rm6uW1eFRGFHREScg2HArnnw87OQeMLcVvkGcyh52XrW1pZPKelZDJ2xkQ2H43JGXakzcuFR2BEREcd3egf89BQc/WOetYAI6PQS1O0JNpu1teVTSvrfloD4Y8JADS8vXAo7IiLiuFLjYMVLsPlTMOzg4QNtRkPrhx1yiYf/Yq5evoGNRy5oCYgipLAjIiKOJzvLDDgrXvprKHndXtDpRQisaGVlVy0pLZO7P93I5qMXchb1vC4i0OqyigWFHRERcSxHfoefnoAzO812aD2zX06VG6yt6xokpmVy9/QNbDkWj7+PB5/fE0nDCoFWl1VsKOyIiIhjSIyBpc/Cjq/Ntk+gOZS86VCnGUp+KQmpmdz16Qa2HTeDzux7WtKgQoDVZRUrzvvbIyIiriErA9Z/BKteg4xkwAZN74abn4WSwVZXd00upGQwePp6dp5MpLSvJ58Nj6R+eQWdoqawIyIi1oleDj89CecPmO0KzaHb6xDe2Nq6CsD55HQGTdvAnphEgkp6MfueSOqU0yz9VlDYERGRonfhKCwZB3sXme2SZaDDeGjU36FXJc+rc8npDJy6nn1nkggp5c2cEZHULOvYC5G6MoUdEREpOpkX4ff34Le3/pj92B0i74ObngIf17i9E5uYxoBP1hMdm0yonzdzRrSkemgpq8sq1hR2RESk8BkG7PvJXLAz/qi5rfIN0PU1KFvX2toK0OmENAZMXcehcymUC/BhzoiWVAkpaXVZxZ7CjoiIFK7zB81+OdFLzbZfOHR+Ger1drrZj6/kZPxFBkxdx9HzqZQPLMEXI1pSMdj5Jj50RQo7IiJSODJSYPUbsPZ9yM4AN09oPQpueAy8Xeu2zvG4VPpPXceJCxeJCCrBnHtaEhGkoOMoFHZERKRgGQbsWQCLx0LiSXNb9Y7QZRKEVLe2tkJw5FwKAz9Zz8n4i1QK9uWLES0JDyxhdVnyNwo7IiJScM4fhB//Dw6uMNuBlczZj2t2calbVn+Kjk1iwNT1xCalUzWkJHNGtCQswMfqsuQfFHZEROTaZaSaI6x+f9e8ZeXuDW0eNRft9HTNqxy7TyUyeNp6zqdkUKusH5/fE0kZP2+ry5JLUNgREZFrs/dHWPwkxB8z29U7mKOsgqtZW1ch2n4insHTNpBwMZP65f2ZNSySoJJeVpcll6GwIyIiV+fCEXOU1f7FZtu/AnSZCHW6u+Qtqz9tOhLH0E83kpSeReOKgcwY2oKAEp5WlyVXoLAjIiL5k5kGa96DX980JwZ084TWI+HGx8HLteeUWXPwHPfM3ERqRjYtqgQx/e7mlPLWV6mj038hERHJuwPL4KfHIe6Q2a5yI3R7E8rUtLauIrByXyz3fbaZ9Cw7N9QI4ePBzSjh5W51WZIHCjsiIvLf4o/DkrGwZ6HZLhVmTgxYv69L37L608+7TvPQnC1kZht0qBPK+wOa4OOpoOMsFHZEROTysjLMSQFXvw6ZqeZaVi0fgLZPgk/xWMF74bZTjP4yiiy7QbcGYbxzZ2O8PJx/sdLiRGFHREQu7dAqc86cc/vNdsXWcMsbULaetXUVoW82n+CJb7ZhN6B34/K8fltDPNwVdJyNwo6IiOSWdBqWjIOd35rtkmWg44vQqF+xuGX1p9nrj/L0vJ0A9G8Rwcu9GuDmVnw+vytR2BEREZM9GzZOgxUvQnoi2Nyg+T3Q7mkoEWh1dUXq49UHeeXHvQDc3boyz3evi60YBT1Xo7AjIiJwKgoWjYZTW8x2+aZwy1sQfp2VVRU5wzB4a+l+Jq+IBuD+ttV4skstBR0np7AjIlKcpSfBL6/A+ilg2MHbH9o/B82GgVvxGm1ktxtMWLSbGWuOAPBEl1o8eJPrLVxaHCnsiIgUR4ZhDiP/6UlIOmVuq98XOr8CfmHW1maBrGw7T323g282nwDgxZ71GNyqsrVFSYFR2BERKW4uHIWfnvhrmYfSleGWN801rYqh9KxsHp0bxU87T+PuZuP12xrSp0kFq8uSAqSwIyJSXGRnwtoPYNWr5pw5bp7myuQ3POayK5P/l4sZ2dz3+WZW7z+Ll7sb7/VvTJf6xe/KlqtT2BERKQ6OrTc7IMfuMtuV2sCtb0GZWtbWZaHEtEzumbGJDUfiKOHpzsd3NeWGGmWsLksKgcKOiIgrS42D5eNh8wyzXSLIXOahUf9iNWfOP8WlZHDX9PXsPJmIn48Hn97dnGaVg6wuSwqJwo6IiCsyDNj+lTk5YOo5c1vjQebkgL7F+0v9dEIag6atJzo2maCSXswa1oL65QOsLksKkcKOiIirORcNP4yGw6vNdkgtuPVtqHy9tXU5gGPnUxk4bR3H4y4S5u/D5/dEUj20lNVlSSFT2BERcRWZafDb2/DbW5CdAR4+0PYJaDUKPLysrs5y+88kMeiT9cQmpVMp2JfPh0cSEeRrdVlSBBR2RERcwaFVZgfkuINmu3oH6PYGBFWxti4HseNEAndNX8+F1ExqlfXjs+EtCPX3sbosKSKWLt26evVqunfvTnh4ODabjfnz5/9rnz179tCjRw8CAgIoWbIkzZs359ixYznPp6Wl8dBDDxEcHEypUqXo27cvZ86cKcJPISJioeSz8N29MKuHGXRKhcHtM2DgNwo6f1h/6Dz9p67jQmomjSoEMPfelgo6xYylYSclJYVGjRrxwQcfXPL5gwcP0qZNG2rXrs3KlSvZvn07zz77LD4+f/2Sjh49moULF/L111+zatUqTp06RZ8+fYrqI4iIWMNuh02fwvtNYfuXgA1a3AsjN0C93sV6pNXf/bIvlrumbyA5PYvIKkHMHtGS0iV1S6+4sRmGYVhdBIDNZmPevHn06tUrZ1u/fv3w9PTks88+u+RrEhISKFOmDHPmzOG2224DYO/evdSpU4e1a9fSsmXLPL13YmIiAQEBJCQk4O/vf82fRUSkUJ3dBwsfgWNrzXa5RnDrO1C+iaVlOZqF204x5qsoMrMNbq4dyocDm+DjWbzW+3J1ef3+tvTKzpXY7XZ++OEHatasSefOnQkNDSUyMjLXra7NmzeTmZlJhw5/TXFeu3ZtKlasyNq1ay977PT0dBITE3M9REQcXlYGrHwVprQxg45nSegyCe5ZoaDzD5+tPcLDc7eSmW1wa8NyTBnUVEGnGHPYsBMbG0tycjKTJk2iS5cu/Pzzz/Tu3Zs+ffqwatUqAE6fPo2XlxeBgYG5Xlu2bFlOnz592WNPnDiRgICAnEdERERhfhQRkWt3fAP870ZY+Yo50qpGJ3hoPbR8ANw11uRPhmHw9tL9PPv9LgwDBresxLv9GuPl4bBfd1IEHPZfiN1uB6Bnz56MHj0agOuuu441a9YwZcoU2rZte9XHHjt2LGPGjMlpJyYmKvCIiGNKS4TlE2DjJ4ABviHQ9VVzhXL1y8nFbjd4YeEuZq09CsCjHWrwSPsa2HSeij2HDTshISF4eHhQt27dXNvr1KnDb7/9BkBYWBgZGRnEx8fnurpz5swZwsIuv5Cbt7c33t7ehVK3iEiB2bcYfhgDiSfN9nUDodNLxX4G5EvJyLIz5qsoFm2PwWaD8T3qcVerylaXJQ7CYa/reXl50bx5c/bt25dr+/79+6lUqRIATZs2xdPTk+XLl+c8v2/fPo4dO0arVq2KtF4RkQKTHAtf3w1f3GkGndKVYfB86PWhgs4lpKRnMXzmRhZtj8HT3ca7/Ror6Egull7ZSU5OJjo6Oqd9+PBhoqKiCAoKomLFijz++OPceeed3HjjjbRr147FixezcOFCVq5cCUBAQADDhw9nzJgxBAUF4e/vz6hRo2jVqlWeR2KJiDgMw4Ctn8PPz0BaPNjcodVDcNNY8NJMv5dyISWDoTM2EnU8nhKe7vxvcFNurKmVyyU3S4eer1y5knbt2v1r+5AhQ5gxYwYA06dPZ+LEiZw4cYJatWoxfvx4evbsmbNvWloajz32GF988QXp6el07tyZDz/88Iq3sf5JQ89FxHLnD8KiR/9azyqsIfSYDOHXWVmVQzsVf5G7pm8gOjaZQF9PPr27OY0rlra6LClCef3+dph5dqyksCMilsnOhDWTYdWrkJUGHiWg3Tho+aBGWV1BdGwyd01bz6mENMoF+DBrWAtqlPWzuiwpYnn9/ta/JBERq5zcAgsehjM7zHbVm8zVyYOqWlqWo9t+Ip67P91IXEoGVcuU5LPhkZQPLGF1WeLAFHZERIpaRgr88gqs+xAMO/gEQpeJ0Ki/hpP/h98OnOO+zzaRkpFNwwoBfHp3c4JLaXStXJnCjohIUYpebvbNif9jQeP6t5mzIJdSp9r/8uOOGB6dG0VGtp3rqwfzv8HNKOWtrzH5b/otEREpCinnYck42D7XbPtXMG9Z1exkbV1OYvb6ozwzfyeGAd0ahPH2ndfh7aHlHyRvFHZERAqTYcDOb+GnJyD1PGCDyPvg5mfAWx1q/4thGLy/Ipo3l+4HYGBkRSb0rI+7m273Sd4p7IiIFJbEGHMG5H0/mu3QutD9PYhobm1dTiLbbvDc9zuZvd685ffwzdUZ3bGmln+QfFPYEREpaIYBUXNgyVhISwA3T7jx/6DNGPDwsro6p5CWmc3DX2zl591nsNnghe71GNK6stVliZNS2BERKUjxx80OyNHLzHZ4Y+j5AZStZ2lZziQ+NYPhMzex+egFvDzcePfO6+jaoJzVZYkTU9gRESkIhgGbP4Wfn4OMJHD3hnZjodUoTQ6YDycupDJk+gYOnk3B38eDT4Y0p0UVrQcm10b/AkVErlXcYVj48F9LPVRoYV7NKVPT2rqczO5Tidz96QZik9IpF+DDzGEtqKlZkaUAKOyIiFwtux02ToVlL0BmqrnUQ/vnzNFWbhoWnR9rDp7jvlmbSUrPolZZP2YMa065AM2KLAVDYUdE5Gqci4YFI+HYWrNdqQ30eA+Cq1lblxNauO0Uj321jYxsOy2qBDH1rmYElPC0uixxIQo7IiL5Yc+GtR/ALy+bC3d6lYKO46HpMHBzs7o6p/PJr4d46Yc9ANzSoBxv3tEIH09dFZOCpbAjIpJXsXvh+wfh5GazXbWdeTUnsKK1dTkhu91g4k97mPrrYQDubl2Z526ti5smC5RCoLAjIvJfsjPh93dh1auQnQHeAdD5ZWg8SAt3XoWMLDv/9/U2Fmw7BcDYrrW598aqmixQCo3CjojIlZzeAfMfhNPbzXbNLuaaVv7h1tblpJLSMrn/8838Hn0eDzcbr9/ekN6NK1hdlrg4hR0RkUvJyoBf34Bf3wR7FpQoDV1fgwa362rOVYpNTGPIpxvZE5NISS93pgxuyg01tNq7FD6FHRGRf4rZBvMegNhdZrtOd+j2JviVtbYuJxYdm8SQ6Rs5GX+RkFLezBjanPrlA6wuS4oJhR0RkT9lZ8Kvb8Hq18yrOb4hcMsbUK+31ZU5tTUHz3H/Z5tJTMuiSkhJZg5tQcVgX6vLkmJEYUdEBCB2D8y7H2KizHadHmbfnJIhlpbl7L7bcoInv91OZrZBs0ql+fiuZgSV1GKoUrQUdkSkeLNnw5rJ5rw52RngEwi3vAn1+6pvzjUwDIPJK6J5a+l+AG5pWI43b9ccOmINhR0RKb7OHzSv5pzYYLZrdDbnzfELs7YuJ5eZbWfcdzv4evMJAO5rW5UnO9fWHDpiGYUdESl+/lzTaunzkHURvPyg6yS4bqCu5lyjxLRMHvx8C79Fn8PNBhN61mdQy0pWlyXFnMKOiBQvF47C9w/BkV/NdpW25grlgRHW1uUCTsZfZOinG9h/JhlfL3c+GNCEdrVDrS5LRGFHRIoJw4Ats2DJOMhIBk9f6DgBmg3XmlYFYOfJBIbN2EhsUjqhft5Mv1tDy8VxKOyIiOtLPAULHobopWa7Yivzao5WKC8QK/aeYeScraRmZFOrrB/ThzanfGAJq8sSyaGwIyKuyzBg+1fw0+OQlgDu3tD+WWj5ILhpVFBB+GzdUZ7/fid2A9pUD+HDQU3w9/G0uiyRXBR2RMQ1JZ+FRY/C3kVmO7wx9P4flKllaVmuwm43mLR4Lx+vPgTAbU0rMLFPAzzddUtQHI/Cjoi4nl3z4YcxkHoe3Dyh7ZPQZjS4609eQUjLzOaxr7bxw44YAMZ0rMmom6tr1XJxWPqXLyKuIzUOfnwcdn5jtsvWh14fQbmG1tblQuJSMhgxaxObj17A093Ga7dp1XJxfAo7IuIaDiyF70dC8mmwuZlXcto+CR7eVlfmMg6dTWbYjI0cOZ+Kv48H/xvcjFbVgq0uS+Q/KeyIiHPLSIWlz8LGT8x2cA3oPQUqNLO2LhezJvoc939uLuZZPrAEM4Y2p0ZZP6vLEskThR0RcV6ntsK3I+D8AbPd4j7oOB48Ney5IH2x4RjPzt9Jlt2gccVAPh7cjDJ+umImzkNhR0Scjz0bfnsLVk4CexaUCoNeH0D1DlZX5lKy7QYTf9zDJ78dBqBHo3Beu62hFvMUp6OwIyLOJe4wzLsPjq8323V6QPd3wTfI2rpcTHJ6Fo98sZXle2MBGN2hJg+314grcU4KOyLiHAwDombDT0+ayz14+UG316BRfy3eWcBOxl9k+IyN7D2dhJeHG2/c3ogejcKtLkvkqinsiIjjSzkPCx/+a4LAiq3MTsilK1talivaeuwCI2Zt5lxyOiGlvPn4rqY0qVja6rJEronCjog4tgPL4PsHIfmMOUFgu3Fw/SNa7qEQLNp+ise+2kZ6lp3aYX5Mu1trXIlrUNgREceUkQpLn4ONU812SC3oOxXKNbK2LhdkGAaTV0Tz1tL9ALSvHcq7/RtTyltfEeIaLF3EZPXq1XTv3p3w8HBsNhvz58/P9fzdd9+NzWbL9ejSpUuufSpXrvyvfSZNmlSEn0JECtyprfBx27+CTov74L5VCjqFIC0zm0e/jMoJOve0qcLHdzVT0BGXYulvc0pKCo0aNWLYsGH06dPnkvt06dKFTz/9NKft7f3vuR0mTJjAiBEjctp+fproSsQp2bPht7dh5UQNKS8C55LTuXfWJrYci8fDzcaEnvUZEFnR6rJECpylYadr16507dr1ivt4e3sTFhZ2xX38/Pz+c5+/S09PJz09PaedmJiY59eKSCG5cAS+uw+OrzPbdXvCre9oSHkh2Xc6iWEzNnIy/iL+Ph58NKgp11cPsboskUJh6W2svFi5ciWhoaHUqlWLBx54gPPnz/9rn0mTJhEcHEzjxo15/fXXycrKuuIxJ06cSEBAQM4jIiKisMoXkf9iGLD1c/joejPoePlBrylw+0wFnULyy75Y+n60hpPxF6kc7Mu8h65X0BGXZjMMw7C6CACbzca8efPo1atXzra5c+fi6+tLlSpVOHjwIOPGjaNUqVKsXbsWd3dzJMZbb71FkyZNCAoKYs2aNYwdO5ahQ4fy1ltvXfa9LnVlJyIigoSEBPz9/QvtM4rIP6TGmUPK9yw02xVbQe//QelK1tblogzDYMaaI7y4aDd2AyKrBDFlUFNKl/SyujSRq5KYmEhAQMB/fn87dNj5p0OHDlGtWjWWLVtG+/btL7nP9OnTue+++0hOTr5k/55LyevJEpECdPhX+O5eSDqlIeVFID0rm+fm7+LLTccBuKNZBV7q1QAvD4e/wC9yWXn9/naq3/KqVasSEhJCdHT0ZfeJjIwkKyuLI0eOFF1hIpJ32Vmw4iWY2d0MOsHV4Z5lcMMYBZ1CcjYpnQFT1/PlpuO42eDpbnV4tW9DBR0pNpxqbOGJEyc4f/485cqVu+w+UVFRuLm5ERoaWoSViUieXDgK394DJzaY7caDoMur4F3K2rpc2I4TCdz72SZiEtLw8/Hg/QFNaFuzjNVliRSpAg87WVlZeHjk7bDJycm5rtIcPnyYqKgogoKCCAoKYvz48fTt25ewsDAOHjzIE088QfXq1encuTMAa9euZf369bRr1w4/Pz/Wrl3L6NGjGTRoEKVLa3pzEYey81tY+CikJ4K3P3R/B+r3tboql7Zg2yke/9qcEblqmZJ8clczqpZRsJTiJ1/XML/66qsrPp+VlcUdd9yR5+Nt2rSJxo0b07hxYwDGjBlD48aNee6553B3d2f79u306NGDmjVrMnz4cJo2bcqvv/6a0xfH29ubuXPn0rZtW+rVq8fLL7/M6NGj+fjjj/PzsUSkMKUnw/yH4JthZtCp0ALu/01BpxBl2w1eXbyXh7/YSnqWnXa1yjD/oesVdKTYylcHZR8fHxYuXEjHjh3/9Vx2dja33347a9euJSYmpkCLLGzqoCxSSE5FwbfD4Xw0YIMb/w/aPgXuTnUH3akkpWXyyNwoVuyNBeD+ttV4vHMt3N20Mry4nrx+f+frL86rr75Knz59WLZsGZGRkTnb7XY7d9xxB7///jsrVqy4+qpFxDXY7bDuQ1j2AtgzwS8c+nwMVW6wujKXdvhcCiNmbSI6NhlvDzde7duQXo3LW12WiOXyFXYeeeQR4uLi6NatG6tXr6ZevXpkZ2dz55138uuvv7JixQrq1atXWLWKiDNIjoX5D0D0MrNd+1boMVkTBBay1fvPMnLOFhLTsgjz9+Hju5rSsEKg1WWJOIR8X0seP348cXFxdOrUiV9++YVnnnmGVatWsXz5curXr18YNYqIs4heBvMegJRY8PCBzq9As2Fg0y2UwmIYBtN+O8wrP+7BbkCTioFMGdSUUH8fq0sTcRhXdeN88uTJXLhwgUaNGlGqVCmWL19Ow4YNC7o2EXEWWRmwfDysfd9sh9aFvtOgbF1r63Jx6VnZPD1vJ99sPgHA7U0r8FLv+nh7aL4ikb/LV9gZM2ZMzs+lS5fGMAyuu+46ZsyYkWu/Ky3VICIu5lw0fDsMYraZ7eYjoNOL4FnC2rpcXGxiGvd9vpmtx+Jxs8Ezt9Rl6PWVsekqmsi/5CvsbN26NVe7VatWZGVl5dquf2gixYRhQNQc+PFxyEyBEqWh5wdQ+xarK3N5247Hc+9nmziTmE5ACU8+GNCENjW0kKfI5eQr7Pzyyy+FVYeIOJO0BFg02pwoEKDyDeZoK/9wa+sqBuZvPckT324nI8tOjdBSTL2rGZVDSlpdlohD02QXIpI/J7fA13dD/FGwuZsLeLYZrXWtCllWtp1XF+9l6q+HAehQJ5S377wOPx9PiysTcXwKOyKSN4YBGz6GJU+bc+cEVoS+0yGiudWVubxzyemMmrOVtYfOA/BQu2o81rEWbpooUCRPFHZE5L9djIcFI2HPQrNd+1azf06JQCurKhaijsfzwOebiUlIo6SXO2/c3oiuDS6/GLKI/JvCjohc2d9vW7l5QqeXIPI+zZ1TBOZuOMZz3+8iI9tcyPPjwU2pHupndVkiTidfC4FOnz6dc+fOFVYtIuJIDAPW/w+mdTKDTmBFGL4EWt6voFPI0rOyGfvddp76bgcZ2XY61S3L9w9dr6AjcpXyFXY+//xzKlSoQOvWrXn11VfZs2dPYdUlIlZKS4Cv7oKfnjD759S+Fe77Fco3tboylxeTcJE7/reOLzYcx2aDxzvXYsqgpuqILHIN8hV2VqxYQUxMDA8++CCbN28mMjKSGjVq8Nhjj7F69Wrsdnth1SkiReXUVvjfjbBngXnbqsskuPNz9c8pAmsPnufW935j2/F4Akp4MmNoCx5qV10dkUWukc0wDONqX5yRkcGKFStYsGABCxcu5OLFi3Tr1o0ePXrQtWtXSpZ0jrkf8rpEvIhLMwzY+AksGQfZGeZtq9tmQAVdzSlsf65vNfGnvWTbDeqW8+d/g5sSEeRrdWkiDi2v39/XFHb+adOmTSxYsIDvv/+e2267jWeffbagDl2oFHak2EtLgAWjYPf3Zrv2rdDzfXNWZClUqRlZPPntDhZuOwVAn8blebl3A0p4ad4ikf9iSdj5u8zMTDw9neMes8KOFGunouDrIXDhyB+jrV6ESHVCLgpHzqVw32eb2XcmCQ83G8/eWpe7WlXSsjsieZTX7+9CG3ruLEFHpNj6522rgIpw+wzdtioiK/ae4ZG5USSlZVHGz5sPBzaheeUgq8sScUmaZ0ekOEpLgAUPw+75ZrvWLdDrA922KgJ2u8F7Kw7wzrIDADStVJoPBzahrL+PxZWJuC6FHZHi5lSUOUnghcPg5gEdJ0DLB3XbqggkXMxkzJdRLN8bC8DglpV49ta6eHnka2CsiORTvv6FTZgwgdTU1MKqRUQKk2HAhqkwraMZdAIiYNgSaPWQgk4R2Hs6kR7v/8byvbF4ebjx+m0NebFXfQUdkSKQrw7K7u7uxMTEEBoaWpg1FTl1UBaXl54MCx+Gnd+a7VrdzLWtfNVHpCjM23qCcd/t5GJmNuUDS/C/wU2pXz7A6rJEnF6hdFAupIFbIlKYzh2ALwfB2b3mbasOL0CrkbqaUwTSMrMZv3A3X2w4BkCb6iG8178xQSW9LK5MpHjJd58dDYkUcSJ7FsK8ByAjCUqFwR0zoWJLq6sqFo6eT+HB2VvYdSoRmw1G3VyDR9rXwF2zIYsUuXyHnZo1a/5n4ImLi7vqgkSkAGRnwYoJ8Pu7ZrvS9XDbp+BX1tq6ionFO0/z+DfbSErLIqikF+/ceR031ixjdVkixVa+w8748eMJCNC9ZhGHlXwWvhkKR341261Gmreu3DX3VWHLzLYz6ae9TPvtMGAOK39/QGPKBZSwuDKR4i3fYadfv34u10FZxGUc32iuVp50CjxLmks+1O9jdVXFwqn4i4ycs4Utx+IBGHFDFZ7oUhtPd422ErFavsKO+uuIOKg/Z0NePBbsmRBS01ypvEwtqysrFlbtP8ujc7dyITUTPx8P3ri9EZ3rhVldloj8QaOxRJxdRiosGg3b55rtuj3NYeXeftbWVQxk2w3eXbafyb9EYxhQL9yfjwY2pWKwVisXcST5Cjt2u72w6hCRqxF3CL4cDGd2gs0dOo7XsPIicjYpnUe/3Mrv0ecBGBhZkWdvrYuPp1YrF3E0Wi5CxFnt+wm+uw/SE6BkGXMRz8ptrK6qWFh/6DyjvthKbFI6vl7uTOzTgJ7Xlbe6LBG5DIUdEWdjz4ZfXoFf3zDbEZFm0PEPt7Ss4sBuN/jf6kO88fM+su0GNUJL8dGgJlQP1S1DEUemsCPiTFLOw7fD4dAvZrvFfdDpJfDQjLyFLT41g8e+2paziGfvxuV5uXd9fL30Z1TE0elfqYizOLnFHFaecBw8faH7e9DwdqurKha2HY/nwdlbOBl/ES8PN8b3qEe/5hEaoSriJBR2RJzB5pnw4/9BdgYEVTWHlZetZ3VVLs8wDGatPcpLP+wmM9ugUrAvHwxookU8RZyMwo6II8tKN0POlllmu9Yt0Psj8NGXbWFLuJjJuO928MOOGAA61yvL67c3wt9HM1GLOBuFHRFHlXTaHFZ+YgPY3ODmZ+H6R8FNM/IWti3HLvDwF1s5ceEiHm42nupam+Ftqui2lYiTUtgRcUQnN8PcQeayDz4BcNt0qN7B6qpc3p+jrd78eR9ZdoOKQb68178x10UEWl2aiFwDS/8XcfXq1XTv3p3w8HBsNhvz58/P9fzdd9+NzWbL9ejSpUuufeLi4hg4cCD+/v4EBgYyfPhwkpOTi/BTiBSwqC9gelcz6ITUghG/KOgUgdikNIZ8uoFXF+8ly25wa8NyLHq4jYKOiAuw9MpOSkoKjRo1YtiwYfTpc+nFCrt06cKnn36a0/b29s71/MCBA4mJiWHp0qVkZmYydOhQ7r33XubMmVOotYsUuOwsWPY8rH3fbNfqBr3/Bz7+1tZVDKzef5YxX0VxLjkDH09ztNUdzTTaSsRVWBp2unbtSteuXa+4j7e3N2Fhl15Qb8+ePSxevJiNGzfSrFkzACZPnky3bt144403CA/XJGviJFLj4Jthf82fc+MTcNNY9c8pZJnZdt74eR//W3UIgNphfrw/oLEmCRRxMQ7/l3TlypWEhoZSq1YtHnjgAc6fP5/z3Nq1awkMDMwJOgAdOnTAzc2N9evXX/aY6enpJCYm5nqIWCZ2D0xtZwYdT1+4fSbc/LSCTiE7dj6V26aszQk6g1tWYv5D1yvoiLggh+6g3KVLF/r06UOVKlU4ePAg48aNo2vXrqxduxZ3d3dOnz5NaGhortd4eHgQFBTE6dOnL3vciRMnMn78+MIuX+S/7VkE8+6DjGQIrAj9voCw+lZX5fIWbjvFuO92kJSehb+PB6/d1ogu9S99BVlEnJ9Dh51+/frl/NygQQMaNmxItWrVWLlyJe3bt7/q444dO5YxY8bktBMTE4mIiLimWkXyxW6H1a/DylfMduUbzCs6JYOtrcvFXczIZvzCXczdeByAZpVK806/66hQ2tfiykSkMDl02PmnqlWrEhISQnR0NO3btycsLIzY2Nhc+2RlZREXF3fZfj5g9gP6Z0dnkSKTngzz74c9C8125P3m+lbumqyuMO09ncjIOVuJjk3GZoOR7arzSPsaeLjrdqGIq3OqsHPixAnOnz9PuXLlAGjVqhXx8fFs3ryZpk2bArBixQrsdjuRkZFWlipyaXGHYe4AiN0N7l5w69vQeJDVVbk0wzD4fP0xXly0m4wsO6F+3rzT7zpaVwuxujQRKSKWhp3k5GSio6Nz2ocPHyYqKoqgoCCCgoIYP348ffv2JSwsjIMHD/LEE09QvXp1OnfuDECdOnXo0qULI0aMYMqUKWRmZjJy5Ej69eunkVjieA6thK/vhosXoFRZuHM2RDS3uiqXlpCayZPfbmfxLrMPX7taZXjj9kYEl9KVXZHixGYYhmHVm69cuZJ27dr9a/uQIUP46KOP6NWrF1u3biU+Pp7w8HA6derEiy++SNmyZXP2jYuLY+TIkSxcuBA3Nzf69u3Le++9R6lSpfJcR2JiIgEBASQkJODvrzlNpIAZBqz7CH5+BoxsKN/UXMjTX4G8MG06Escjc6M4GX8RT3cbT3bRkg8iriav39+Whh1HobAjhSYzDRaNhm1/THLZaIB568rTx9q6XFi23eCjldG8vewA2XaDysG+TO7fhAYVtHiqiKvJ6/e3U/XZEXEqSWfM/jknN5kLeXZ6GVo+ALqyUGhOxV/ksa+2sfaQOR9X78blebFXfUp560+dSHGmvwAiheH0TphzJySeAJ9AuH0GVPv3LVspOAu3neLpeTtITMvC18udF3vWp2/TClaXJSIOQGFHpKDtX2Iu/ZCRDMHVYcBXEFzN6qpcVmJaJi98v4vvtp4EoFFEIO/ceR1VQkpaXJmIOAqFHZGCktMR+Wkw7FDlRrhjFpQobXVlLmvD4ThGf2l2QnazwcibazDq5up4au4cEfkbhR2RgpCdCT8+Dps/NdtNhsAtb2qiwEKSkWXnnWX7+WjVQQwDKgb58vad19G0koKliPybwo7ItboYD18PMefRwWbOhtzqIXVELiTRsck8+uVWdp40F/C9vWkFnu9RT52QReSy9NdB5FrEHTI7Ip/bD54l4bZpUKur1VW5pD9nQn75h92kZdoJ9PVkYu8GdG1QzurSRMTBKeyIXK2ja2DuQLgYB/7lof9cKNfQ6qpc0tmkdJ74Zhu/7DsLwA01Qnjj9kaU9dd8RSLy3xR2RK5G1BewYBTYMyG8sRl0/C6/+KxcvaW7z/DUt9s5n5KBl4cbY7vWZkiryri56TahiOSNwo5Iftjt8MtL8OubZrtuT+g1Bbx8ra3LBaVmZPHioj18seEYALXD/Hi3X2NqhflZXJmIOBuFHZG8ykiF+ffD7u/N9g3/B+2eBjcNcy5oUcfjGf1lFIfPpWCzwYgbqvJYp5p4e7hbXZqIOCGFHZG8SDoNX/SDU1vBzRN6TIbr+ltdlcvJyrbz0cqDvLPcXNeqXIAPb97eiNbVQ6wuTUScmMKOyH85veOPpR9OQokg6DcbKrW2uiqXc+x8KqO/imLz0QsA3NKwHK/0akCAr+YqEpFro7AjciX7foJvhkNmCoTUhAFfQlBVq6tyKYZh8M3mE7ywYBcpGdn4eXswoVc9el1XHpvmKhKRAqCwI3I5az+EJeMAA6reBLfPhBKBFhflWs4np/PM/J38tPM0AC0qB/HmHY2ICFKHbxEpOAo7Iv9kt5vrW6370Gw3HQrdXtfSDwVs8c4Ynp63k/MpGXi42RjTqSb33VgNdw0pF5ECprAj8neZF+G7e2HPArPd8UVoPUpLPxSg+NQMnl+wi++jTgFQq6wfb97RiPrlAyyuTERclcKOyJ9S4+CL/nB8Hbh7Qa+PoMFtVlflUlbsPcNT3+4gNikdNxvc37Yaj3SooSHlIlKoFHZEAC4cgc9vg/MHwCcA+s2Bym2srsplJKZl8uLC3Xy9+QQAVcuU5M3bG9G4olYpF5HCp7AjcmorzL4DUmLBvwIM+gZC61hdlctYvf8sT367nZiENGw2uKdNFR7rVAsfT13NEZGiobAjxduBpfDVEHNoedkGMPBr8Ncq2gUhOT2LV37cw5z15nIPlYJ9eeP2RjSvHGRxZSJS3CjsSPG1ZRYsfBSMbKjaDu6YBT7+VlflEtYePM/j32zjxIWLAAxpVYknu9bG10t/ckSk6OkvjxQ/hgErJ8KqV812owHQ4z0NLS8AFzOyeXXxXmasOQJA+cASvH5bQy33ICKWUtiR4iU707yaE/W52b7xCWg3TkPLC8Dmo3H839fbOXwuBYD+LSoyrltt/HwUIkXEWgo7UnykJ8FXd8HBFWBzh1vfhqZDrK7K6aVlZvP20v1M/fUQdgPC/H149baGtK1ZxurSREQAhR0pLpJOw+zbzEU9PX3NpR9qdrK6Kqe3/UQ8j321jQOxyQD0bVKB57rXJaCEruaIiONQ2BHXd3YffN4XEo5DyTIw4Cso38TqqpxaRpadySsO8OHKg2TbDUJKeTOxTwM61i1rdWkiIv+isCOu7ega+KIfpCVAcHUY+A0EVbG6Kqe2+1Qij329jT0xiQB0bxTOhB71KF3Sy+LKREQuTWFHXNeueeY6V9kZEBEJ/eeCr+Z4uVqZ2XY+WnmQySsOkJltEFTSi5d61adbA81LJCKOTWFHXNOGqfDj44ABdbpDn6ngWcLqqpzWjhMJPP7NNvaeTgKgc72yvNy7ASGlvC2uTETkvynsiGsxDHP+nJUTzXbze6Dra+CmpQmuRlpmNm8v28/U1eZIq6CSXjzfvS49GoVj03B9EXESCjviOux2WPwkbPjYbN80Fto+qTl0rtL6Q+d56rsdOfPm9GgUzvPd6xKsqzki4mQUdsQ1ZGXA/Adg5zeADbq9Di1GWF2VU0pKy+TVxXv5fJ25plWYvw8v9apPB420EhEnpbAjzi8jBb4cDAeXg5sH9P4fNLjN6qqc0i97Y3l63g5OJaQBMCCyIk91rY2/ZkEWESemsCPOLTUO5twBJzaakwXe+RlU72B1VU4nLiWDFxftZt7Wk4C5QvnEPg1oXU1rWomI81PYEeeVcBI+7wNn90KJ0jDga4hobnVVTsUwDBZtj+GFBbs4n5KBmw3uuaEqozvUpISXOnWLiGtQ2BHndO4AfNbbnBXZLxwGz4PQ2lZX5VTOJKbx9LydLNtzBoBaZf149baGXBcRaG1hIiIFTGFHnM/JLeY6V6nnzVmRB8+DwIpWV+U0DMPgy43HefnHPSSlZeHpbmNkuxo8cFM1vDzcrC5PRKTAKeyIczm0EuYOhIxkKHcdDPoWSqpfSV4dO5/KU99tZ83B8wA0igjktb4NqRXmZ3FlIiKFx9L/jVu9ejXdu3cnPNycoGz+/PmX3ff+++/HZrPxzjvv5NpeuXJlbDZbrsekSZMKt3Cxxq75MPt2M+hUuRHuXqSgk0fZdoNPfj1Ep3dWsebgeXw83Xjmljp890BrBR0RcXmWXtlJSUmhUaNGDBs2jD59+lx2v3nz5rFu3TrCw8Mv+fyECRMYMeKvOVX8/PTH2+VsngkLH8Fc/qEH9P0EPDS5XV7siUlk7Hc7iDoeD0CrqsFM6tuASsElrS1MRKSIWBp2unbtSteuXa+4z8mTJxk1ahRLlizhlltuueQ+fn5+hIWF5fl909PTSU9Pz2knJibm+bVigbUfwpKx5s9N74Zb3tLyD3lwMSObd5bv55NfD5NtN/Dz9uDpW+pwZ/MILfUgIsWKQ/dGtNvtDB48mMcff5x69epddr9JkyYRHBxM48aNef3118nKyrricSdOnEhAQEDOIyIioqBLl4JgGLDq9b+CzvWPwK3vKOjkwcp9sXR8exX/W3WIbLtB1/phLB3Tln4tKiroiEix49AdlF999VU8PDx4+OGHL7vPww8/TJMmTQgKCmLNmjWMHTuWmJgY3nrrrcu+ZuzYsYwZMyannZiYqMDjaAwDlj0Pv79rtts9Azf+n9a5+g+xSWlMWLibRdtjACgfWILxPeppqQcRKdYcNuxs3ryZd999ly1btlzx/0T/HloaNmyIl5cX9913HxMnTsTb+9J9Ory9vS/7nDgAux1+ehw2fmK2O78CrR6ytiYHZ7cbzN14nEk/7SExLQs3Gwy7vgqjO9akpLfD/jMXESkSDvtX8NdffyU2NpaKFf+aPyU7O5vHHnuMd955hyNHjlzydZGRkWRlZXHkyBFq1apVRNVKgcnOggWjYNscwAbd3zH76chl7TudxLh5O9h89AIADcoHMLFPA+qXD7C4MhERx+CwYWfw4MF06JB7jaPOnTszePBghg4detnXRUVF4ebmRmhoaGGXKAUtKwO+uwd2fw82d3NBz4a3W12Vw0rLzGbyigP8b9UhsuwGJb3ceaxTLYa0roy7m273iYj8ydKwk5ycTHR0dE778OHDREVFERQURMWKFQkODs61v6enJ2FhYTlXbNauXcv69etp164dfn5+rF27ltGjRzNo0CBKly5dpJ9FrlHmRfjqLjjwM7h7wW2fQp1bra7KYf124BxPz9/B0fOpAHSsW5bxPeoRHljC4spERByPpWFn06ZNtGvXLqf9Z/+bIUOGMGPGjP98vbe3N3PnzuWFF14gPT2dKlWqMHr06Fz9eMQJpCfDF/3gyK/gUQL6zYbq7a2uyiGdS07n5R/25KxOHubvw/ie9ehcL+9TL4iIFDc2wzAMq4uwWmJiIgEBASQkJODv7291OcXLxXhzVuQTG8DLDwZ+BZVaW12VwzEMg682HeeVH/eScDETmw2GtKrMY51q4ufjaXV5IiKWyOv3t8P22ZFiIOUcfNYLTu+AEqXNda7KN7W6KocTHZvMuHk72HA4DoC65fyZ2KcBjbQ6uYhInijsiDUSY2BWTzi3D0qGwl3zoezlJ44sjtIys/lw5UE+WhlNZrZBCU93xnSsydDrK+Ph7tDzgYqIOBSFHSl6CSdhZneIOwj+FeCu7yGkutVVOZRf9sby/IJdHIszOyDfXDuUCT3rUaG0r8WViYg4H4UdKVrxx2HmrXDhCARWhCGLoHQlq6tyGCfjLzJh4S6W7DoDmB2Qn+tel671w7TMg4jIVVLYkaJz4agZdOKPQenKZtAJ1DIdABlZdqb9dpj3lh/gYmY27m42hrepwsPta1BKMyCLiFwT/RWVohF3CGb2gITjEFQNhiyEgPJWV+UQ1hw8x3Pf7yI6NhmAFpWDeLFXfWqF+VlcmYiIa1DYkcJ3/qDZRyfxJATXMIOOfzmrq7JcbFIar/ywh/lRpwAILunFuG516NOkvG5ZiYgUIIUdKVznDsCMWyH5NITUMoOOX/FegTsr287n647y5s/7SUrPwmaDQZGV+L9OtQjw1Zw5IiIFTWFHCk/sXvOKTkoshNaFuxZAqTJWV2WpLccu8My8neyOSQSgYYUAXupVn4YVAq0tTETEhSnsSOE4sxtm9YCUs1C2gTm8vGTwf7/ORV1IyeDVxXuZu/E4AP4+HjzRpTb9W1TUop0iIoVMYUcK3pld5hWd1PMQ1tAMOr5BVldliWy7wZcbj/P6kr1cSM0E4LamFXiqa21CSnlbXJ2ISPGgsCMF68zuv4JOeGMYPM9cCqIY2nz0Ai8s2MWOkwkA1Crrx0u969O8cvEMfiIiVlHYkYLzZx+d1PNQ7joYPB9KBFpcVNGLTUrj1Z/28e2WEwD4eXvwaMea3NWqEp5a5kFEpMgp7EjBOLv/j6Bzzrx1NXhesQs6mdl2Zq45wrvLDpCUngXA7U0r8ESX2pTx0y0rERGrKOzItTt3wJwZOSX2r87IxayPzu/R53hhwS4O/DExYMMKAbzQox5NKhbPW3giIo5EYUeuzbnoP+bROQNl6xe7oHMy/iIv/7CbH3ecBiCopBdPdK7FHc0icNMoKxERh6CwI1fv/EHzik7y6T/m0Sk+w8vTMrOZuvoQH6yMJi3TjpsN7mpVmdEdampiQBERB6OwI1cn7pDZRycpBsrUNicMLBlidVWFzjAMlu2J5cVFuzkWlwpAiypBjO9Rjzrl/C2uTkRELkVhR/LvwhGY8cdaV38uAVEMZkY+dDaZCYt2s3LfWQDK+nszrlsdejQK11pWIiIOTGFH8if++B9B58Rfi3qWCrW6qkKVlJbJ+yuimf77YTKzDTzdbdxzQ1VGtqtOSW/9ExIRcXT6Sy15lxhj3rpKOAZB1Vx+Uc9su8E3m4/z+pJ9nEvOAOCmWmV47ta6VC1TyuLqREQkrxR2JG+Sz8KsnnDhMARWhCELwL+c1VUVmg2H4xi/cBe7TpkLdlYNKcmzt9alXW3XvoolIuKKFHbkv6XGwWe94dw+8As3r+gEVLC6qkJxMv4iE3/cw6LtMQD4+XjwSPsa3NWqMl4emv1YRMQZKezIlaUlwud94cwOKBlqBp3Sla2uqsClZmQxZdUh/rfqIOlZdmw26N+iIo91rEmwFuwUEXFqCjtyeRkpMOcOOLUFSgSZ8+iEVLe6qgJlGAYLtp1i4o97OZ2YBkDLqkE8d2s96oZrKLmIiCtQ2JFLy7wIX/SDY2vBO8Bc66psXaurKlDbjsczYdFuNh+9AECF0iV4ulsdutQP01ByEREXorAj/5adCV8NgcOrwasUDPoWwq+zuqoCE5NwkdeX7OO7LScB8PVy56F21Rnepgo+nu4WVyciIgVNYUdys2fDd/fCgSXgUQIGfAkRza2uqkCkpGfxv9WH+Hj1QdIy7QD0aVyeJ7vWpqy/j8XViYhIYVHYkb8YBix6FHZ9B26ecOfnULmN1VVdM7vd4JstJ3hjyT5ik9IBaF65NM/eWpeGFQKtLU5ERAqdwo6YDAN+fga2zAKbG/T9BGp0sLqqa7b24Hle+mF3znw5FYN8Gdu1tvrliIgUIwo7Ylr9Oqx93/y5x2So18vScq7V4XMpTPxxDz/vPgOAn7cHo9pXZ0jrynh7qF+OiEhxorAjsO4j+OVl8+cuk6DxIGvruQYJqZm8u/wAs9YeIctu4O5mY0CLijzaoYbmyxERKaYUdoq7rZ/D4qfMn28aBy0fsLaeq5SZbefzdUd5d/kB4lMzAWhXqwzjutWhRlk/i6sTERErKewUZ7sXwIJR5s+tRkLbJ6yt5yoYhsGSXWd4dfFeDp9LAaBm2VI8c0tdbqxZxuLqRETEESjsFFeHVsK3w8GwQ+PB0OklcLIOu1uPXeDlH/aw6Y9JAYNLejGmU03ubBaBh7vWsRIREZPCTnF0cjPMHQjZGVCnO3R/16mCzrHzqby6ZC8//LFYp4+nGyNuqMp9batRylu/0iIikpu+GYqbs/vh89sgIxmqtIW+08DNOUYnxadmMHlFNLPWHiEz28Bmg9uaVGBMp5qUCyhhdXkiIuKgFHaKk/jj8FkvuBgH4Y2h32zwcPwRSulZ2cxac5TJKw6QmJYFwA01QhjbtY4W6xQRkf9kaceG1atX0717d8LDw7HZbMyfP/+y+95///3YbDbeeeedXNvj4uIYOHAg/v7+BAYGMnz4cJKTkwu3cGeUch4+6w2JJyG4Bgz8Frwde5SS3W7wfdRJOry1ipd/3ENiWha1w/yYOawFnw2PVNAREZE8sfTKTkpKCo0aNWLYsGH06dPnsvvNmzePdevWER4e/q/nBg4cSExMDEuXLiUzM5OhQ4dy7733MmfOnMIs3blkpMCc2+H8AfCvAHfNh5LBVld1Rb9Hn2PST3vZcTIBgLL+3jzWsRZ9m1bA3c15+heJiIj1LA07Xbt2pWvXrlfc5+TJk4waNYolS5Zwyy235Hpuz549LF68mI0bN9KsWTMAJk+eTLdu3XjjjTcuGY6KnT9XMD+5GUqUhsHfQUAFq6u6rN2nEpm0eC+r958FoJS3B/fdWJXhN1TB10t3XUVEJP8c+tvDbrczePBgHn/8cerVq/ev59euXUtgYGBO0AHo0KEDbm5urF+/nt69e1/yuOnp6aSnp+e0ExMTC754R2C3w/cjIXrpHyuYfw1lalld1SWduJDKWz/vZ17USQwDPN1tDIysxKibq2vmYxERuSYOHXZeffVVPDw8ePjhhy/5/OnTpwkNDc21zcPDg6CgIE6fPn3Z406cOJHx48cXaK0OadnzsH0u2NzhjpkQ0dzqiv4lPjWDD36JZubao2Rk2QG4tWE5Hu9ci0rBJS2uTkREXIHDhp3Nmzfz7rvvsmXLlgJfnXrs2LGMGTMmp52YmEhERESBvofl1n4Aa94zf+4xGWp2traef0jLzGbGmiN8+Et0zgirVlWDeaprbRpFBFpbnIiIuBSHDTu//vorsbGxVKxYMWdbdnY2jz32GO+88w5HjhwhLCyM2NjYXK/LysoiLi6OsLCwyx7b29sbb28XvjWy4xtYMs78ucML0HigpeX8XVa2ne+2nOTtZfuJSUgDoHaYH092rc1NNcsUeLAVERFx2LAzePBgOnTokGtb586dGTx4MEOHDgWgVatWxMfHs3nzZpo2bQrAihUrsNvtREZGFnnNDuHQKph3v/lz5ANw/aOWlvOnP9eweuPnfUTHmlMDhAf4MKZTLXo3Lq8RViIiUmgsDTvJyclER0fntA8fPkxUVBRBQUFUrFiR4ODcw6M9PT0JCwujVi2zk22dOnXo0qULI0aMYMqUKWRmZjJy5Ej69etXPEdind4BXw4CeybU6w2dX3GIZSDWHTrPq4v3svVYPACBvp48dFN1BreqhI+nc8zeLCIizsvSsLNp0ybatWuX0/6zH82QIUOYMWNGno4xe/ZsRo4cSfv27XFzc6Nv37689957hVGuY4s/DrNvh/REqNQGek0BN2sXw9x9KpHXluxl5T5zGHkJT3eGt6nCvW2r4u/jaWltIiJSfNgMwzCsLsJqiYmJBAQEkJCQgL+/E87Ke/ECTO8CZ/dCmTowbDGUCLSsnKPnU3h76X6+33YKwwAPNxv9WkTw8M01CPX3sawuERFxLXn9/nbYPjuSR1np5grmZ/eCXzgM+sayoHMmMY33lh/gy43HybKbGfrWhuX4v061qByiYeQiImINhR1nZrfD/Afg6O/g7W8GHQtmR45PzeCjVQeZueYIaZnmXDlta5bh8c61qF8+oMjrERER+TuFHWe2fDzs/BbcPODOz6Dsv2eZLkwp6Vl8+vth/rf6EEl/zJXTtFJpnuhci8iqjr32loiIFB8KO85q4yfw+zvmzz3eh6o3Fdlbp2dl88X6Y7z/SzTnkjMAc66cxzvX4ubaoZorR0REHIrCjjPavwR+fNz8ud3TcF3/InnbPycEfHf5AU7GXwSgUrAvYzrWpHvDcNw0V46IiDgghR1nE7MNvh4Khh0aD4IbHy/0t7TbDX7YEcPbS/dz6FwKAGX9vXm4fQ3uaBaBp7u1Q9xFRESuRGHHmSSchDl3QmaKedvq1ncKddJAwzBYsTeWN37ez54Yc2X40r6ePNSuOoNaakJAERFxDgo7ziI9CebcAUkx5lw6d8wC98KZmM8wDNYcPM8bP+/LmfXYz9uDETdWZVibKpTy1q+NiIg4D31rOYPsLPPW1ZmdUKosDPwafApnSPfGI3G8+fM+1h2KA8DH0427W1fh/rZVCfT1KpT3FBERKUwKO47OMGDxUxC9FDxKQP+5EBhR4G+z7Xg8by7dz+r95tIOXu5uDIisyIPtqhHqp1mPRUTEeSnsOLr1U2DjVMAGfadC+SYFevg9MYm8tXQ/S3efAcylHW5vFsGom6sTHliiQN9LRETECgo7jmzfYlg81vy54wSo073ADr3/TBLvLNvPjztOA+Bmg96NK/BI+xpUDPYtsPcRERGxmsKOozq9A74ZBhjQZAi0HlUghz14Npn3lh9gwR+LdIK5ftWjHWpSPbRUgbyHiIiII1HYcURJZ2BOP3OIeZUb4ZY3r3mI+dHzKby3PJp5W0/wxxqddKkXxqMda1A7zAlXehcREckjhR1Hk3kR5g6AxBMQXP2ah5gfj0tl8ooDfLvlJNl/pJwOdUJ5tENNLdIpIiLFgsKOIzEM+H4knNwEPoEw4CsoUfqqDnXiQiof/BLN15tOkPVHyGlbswyjO9bkuojAgqtZRETEwSnsOJLVb8DOb/5YxfxzCK6W70OcjL/4R8g5Tma2GXJuqBHCox1q0rTS1QUnERERZ6aw4yh2fw+/vGT+fMubUOWGfL381B8h56u/hZzW1YIZ3bEmzSsHFXS1IiIiTkNhxxHEbIN595s/t3wQmt6d55eeir/Ihyuj+WrjCTKy7QC0qhrMox1qEFk1uBCKFRERcS4KO1ZLOgNf9IfMVKjeATq+mKeX/Rlyvtz415WcVlWDeaRDDVoq5IiIiORQ2LFSVjp8OQgST0JwDeg7Ddyv/J/kZPxFPvzH7aqWVYN4pH1NWlVTyBEREfknhR2rGAYsGg0nNpiLevafCyUCL7v78bhUPlx5kG8260qOiIhIfijsWGXdRxA1G2xucNunEFL9krsdO28OIf92y19DyFtXC+aR9uqTIyIikhcKO1Y4uAJ+ftr8udNLUL39v3Y5fC6FD36JZt7WvyYDbFM9hEc61NDoKhERkXxQ2ClqcYfg66Fg2KHRAHP01d9ExybzwS/RfB91MmdZhxtrluGR9tVpWkkhR0REJL8UdopSehJ8MQDS4qF8M7j17Zw1r/adTmLyigP8sCMmZ4HOm2uHMurm6jSuqMkARURErpbCTlGx22H+A3B2D5QKM2dI9vRh58kEJq84wJJdZ3J27Vi3LA/fXIMGFbR2lYiIyLVS2Ckqv74JexaCuxfc+Tlb4n14/7uNrNgbC5gXeLrWD2PUzTWoU06rkIuIiBQUhZ2isH8J/PIyAAdbjOf5JXZ+i14DgJsNujcKZ2S76tQo62dllSIiIi5JYaewnYvG+PYebBgsKXEr9/1SCTiHh5uN3o3L82C76lQJKWl1lSIiIi5LYacQ2S8mkjrrTkqlJ7LRXpORF+7Ay92N25tV4P621YgI8rW6RBEREZensFNIDLud9e/0o1V6NKeN0ow2xjD4+hrce2NVwgJ8rC5PRESk2FDYKSQ2IxufgDJknnFnWf3X+b5rD4JLeVtdloiISLGjsFNY3D2pevcnpMUdYFCFulZXIyIiUmwp7BSiAF9P8FXQERERsZKb1QWIiIiIFCaFHREREXFpCjsiIiLi0iwNO6tXr6Z79+6Eh4djs9mYP39+rudfeOEFateuTcmSJSldujQdOnRg/fr1ufapXLkyNpst12PSpElF+ClERETEkVkadlJSUmjUqBEffPDBJZ+vWbMm77//Pjt27OC3336jcuXKdOrUibNnz+bab8KECcTExOQ8Ro0aVRTli4iIiBOwdDRW165d6dq162WfHzBgQK72W2+9xbRp09i+fTvt27fP2e7n50dYWFih1SkiIiLOy2n67GRkZPDxxx8TEBBAo0aNcj03adIkgoODady4Ma+//jpZWVlXPFZ6ejqJiYm5HiIiIuKaHH6enUWLFtGvXz9SU1MpV64cS5cuJSQkJOf5hx9+mCZNmhAUFMSaNWsYO3YsMTExvPXWW5c95sSJExk/fnxRlC8iIiIWsxmGYVhdBIDNZmPevHn06tUr1/aUlBRiYmI4d+4cU6dOZcWKFaxfv57Q0NBLHmf69Oncd999JCcn4+196eUZ0tPTSU9Pz2knJiYSERFBQkIC/v7+BfaZREREpPAkJiYSEBDwn9/fDn8bq2TJklSvXp2WLVsybdo0PDw8mDZt2mX3j4yMJCsriyNHjlx2H29vb/z9/XM9RERExDU5fNj5J7vdnuuqzD9FRUXh5uZ22Ss/IiIiUrxY2mcnOTmZ6OjonPbhw4eJiooiKCiI4OBgXn75ZXr06EG5cuU4d+4cH3zwASdPnuT2228HYO3ataxfv5527drh5+fH2rVrGT16NIMGDaJ06dJWfSwRERFxIJaGnU2bNtGuXbuc9pgxYwAYMmQIU6ZMYe/evcycOZNz584RHBxM8+bN+fXXX6lXrx5g3o6aO3cuL7zwAunp6VSpUoXRo0fnHEdERETEYTooWykhIYHAwECOHz+u/jsiIiJO4s8BRvHx8QQEBFx2P4cfel4UkpKSAIiIiLC4EhEREcmvpKSkK4YdXdnB7PR86tQp/Pz8sNlsBXbcPxOnrhgVLp3noqNzXTR0nouGznPRKMzzbBgGSUlJhIeH4+Z2+TFXurIDuLm5UaFChUI7voa3Fw2d56Kjc100dJ6Lhs5z0Sis83ylKzp/crqh5yIiIiL5obAjIiIiLk1hpxB5e3vz/PPPX3bZCikYOs9FR+e6aOg8Fw2d56LhCOdZHZRFRETEpenKjoiIiLg0hR0RERFxaQo7IiIi4tIUdkRERMSlKexcow8++IDKlSvj4+NDZGQkGzZsuOL+X3/9NbVr18bHx4cGDRrw448/FlGlzi0/53nq1KnccMMNlC5dmtKlS9OhQ4f//O8if8nv7/Sf5s6di81mo1evXoVboIvI73mOj4/noYceoly5cnh7e1OzZk39/ciD/J7nd955h1q1alGiRAkiIiIYPXo0aWlpRVStc1q9ejXdu3cnPDwcm83G/Pnz//M1K1eupEmTJnh7e1O9enVmzJhRuEUactXmzp1reHl5GdOnTzd27dpljBgxwggMDDTOnDlzyf1///13w93d3XjttdeM3bt3G88884zh6elp7Nixo4grdy75Pc8DBgwwPvjgA2Pr1q3Gnj17jLvvvtsICAgwTpw4UcSVO5/8nus/HT582Chfvrxxww03GD179iyaYp1Yfs9zenq60axZM6Nbt27Gb7/9Zhw+fNhYuXKlERUVVcSVO5f8nufZs2cb3t7exuzZs43Dhw8bS5YsMcqVK2eMHj26iCt3Lj/++KPx9NNPG999950BGPPmzbvi/ocOHTJ8fX2NMWPGGLt37zYmT55suLu7G4sXLy60GhV2rkGLFi2Mhx56KKednZ1thIeHGxMnTrzk/nfccYdxyy235NoWGRlp3HfffYVap7PL73n+p6ysLMPPz8+YOXNmYZXoMq7mXGdlZRmtW7c2PvnkE2PIkCEKO3mQ3/P80UcfGVWrVjUyMjKKqkSXkN/z/NBDDxk333xzrm1jxowxrr/++kKt05XkJew88cQTRr169XJtu/POO43OnTsXWl26jXWVMjIy2Lx5Mx06dMjZ5ubmRocOHVi7du0lX7N27dpc+wN07tz5svvL1Z3nf0pNTSUzM5OgoKDCKtMlXO25njBhAqGhoQwfPrwoynR6V3OeFyxYQKtWrXjooYcoW7Ys9evX55VXXiE7O7uoynY6V3OeW7duzebNm3NudR06dIgff/yRbt26FUnNxYUV34VaCPQqnTt3juzsbMqWLZtre9myZdm7d+8lX3P69OlL7n/69OlCq9PZXc15/qcnn3yS8PDwf/3jktyu5lz/9ttvTJs2jaioqCKo0DVczXk+dOgQK1asYODAgfz4449ER0fz4IMPkpmZyfPPP18UZTudqznPAwYM4Ny5c7Rp0wbDMMjKyuL+++9n3LhxRVFysXG578LExEQuXrxIiRIlCvw9dWVHXNqkSZOYO3cu8+bNw8fHx+pyXEpSUhKDBw9m6tSphISEWF2OS7Pb7YSGhvLxxx/TtGlT7rzzTp5++mmmTJlidWkuZeXKlbzyyit8+OGHbNmyhe+++44ffviBF1980erS5Brpys5VCgkJwd3dnTNnzuTafubMGcLCwi75mrCwsHztL1d3nv/0xhtvMGnSJJYtW0bDhg0Ls0yXkN9zffDgQY4cOUL37t1zttntdgA8PDzYt28f1apVK9yindDV/E6XK1cOT09P3N3dc7bVqVOH06dPk5GRgZeXV6HW7Iyu5jw/++yzDB48mHvuuQeABg0akJKSwr333svTTz+Nm5uuDxSEy30X+vv7F8pVHdCVnavm5eVF06ZNWb58ec42u93O8uXLadWq1SVf06pVq1z7AyxduvSy+8vVnWeA1157jRdffJHFixfTrFmzoijV6eX3XNeuXZsdO3YQFRWV8+jRowft2rUjKiqKiIiIoizfaVzN7/T1119PdHR0TpgE2L9/P+XKlVPQuYyrOc+pqan/CjR/BkxDy0gWGEu+Cwut63MxMHfuXMPb29uYMWOGsXv3buPee+81AgMDjdOnTxuGYRiDBw82nnrqqZz9f//9d8PDw8N44403jD179hjPP/+8hp7nQX7P86RJkwwvLy/jm2++MWJiYnIeSUlJVn0Ep5Hfc/1PGo2VN/k9z8eOHTP8/PyMkSNHGvv27TMWLVpkhIaGGi+99JJVH8Ep5Pc8P//884afn5/xxRdfGIcOHTJ+/vlno1q1asYdd9xh1UdwCklJScbWrVuNrVu3GoDx1ltvGVu3bjWOHj1qGIZhPPXUU8bgwYNz9v9z6Pnjjz9u7Nmzx/jggw809NzRTZ482ahYsaLh5eVltGjRwli3bl3Oc23btjWGDBmSa/+vvvrKqFmzpuHl5WXUq1fP+OGHH4q4YueUn/NcqVIlA/jX4/nnny/6wp1Qfn+n/05hJ+/ye57XrFljREZGGt7e3kbVqlWNl19+2cjKyiriqp1Pfs5zZmam8cILLxjVqlUzfHx8jIiICOPBBx80Lly4UPSFO5Fffvnlkn9z/zy3Q4YMMdq2bfuv11x33XWGl5eXUbVqVePTTz8t1BpthqFrcyIiIuK61GdHREREXJrCjoiIiLg0hR0RERFxaQo7IiIi4tIUdkRERMSlKeyIiIiIS1PYEREREZemsCMiIiIuTWFHREREXJrCjoiIiLg0hR0RkUu46aabePTRR6/pGIZhcO+99xIUFITNZiMqKqpAahOR/FHYERGHN3ToUJ555hmry8i3xYsXM2PGDBYtWkRMTAz169e3uiSRYsnD6gJERK4kOzubRYsW8cMPP1hdSr4dPHiQcuXK0bp168vuk5GRgZeXVxFWJVL86MqOSDH1xRdfUKJECWJiYnK2DR06lIYNG5KQkHDVx61QoQIffvhhrm1r1qzB19eXo0eP5vt4a9aswdPTk+bNm1/y+ZtuuolRo0bx6KOPUrp0acqWLcvUqVNJSUlh6NCh+Pn5Ub16dX766adcr0tPT+fhhx8mNDQUHx8f2rRpw8aNGy9bh91uZ+LEiVSpUoUSJUrQqFEjvvnmm8vuf/fddzNq1CiOHTuGzWajcuXKOfWOHDmSRx99lJCQEDp37szixYtp06YNgYGBBAcHc+utt3Lw4MF/vf9rr71G9erV8fb2pmLFirz88st5PIsixZvCjkgx1a9fP2rWrMkrr7wCwPPPP8+yZcv46aefCAgIuOrjRkZG5goNhmHw6KOPMnr0aCpVqpTv4y1YsIDu3btjs9kuu8/MmTMJCQlhw4YNjBo1igceeIDbb7+d1q1bs2XLFjp16sTgwYNJTU3Nec0TTzzBt99+y8yZM9myZQvVq1enc+fOxMXFXfI9Jk6cyKxZs5gyZQq7du1i9OjRDBo0iFWrVl1y/3fffZcJEyZQoUIFYmJicp2TmTNn4uXlxe+//86UKVNISUlhzJgxbNq0ieXLl+Pm5kbv3r2x2+05rxk7diyTJk3i2WefZffu3cyZM4eyZcvm93SKFE+GiBRbCxcuNLy9vY2XXnrJKF26tLFz586c53r16mUEBgYaffv2zdcxX3vtNaNevXo57ZkzZxphYWFGUlLSVR23Ro0axqJFiy77fNu2bY02bdrktLOysoySJUsagwcPztkWExNjAMbatWsNwzCM5ORkw9PT05g9e3bOPhkZGUZ4eLjx2muv5Rz3kUceMQzDMNLS0gxfX19jzZo1ud57+PDhRv/+/S9b29tvv21UqlTpX/U2btz4ip/57NmzBmDs2LHDMAzDSExMNLy9vY2pU6de8XUicmm6siNSjN16663UrVuXCRMmMG/ePOrVq5fz3COPPMKsWbPyfcyWLVuyZ88ekpOTSUlJYdy4cbz00kuUKlUq38fds2cPp06don379lfcr2HDhjk/u7u7ExwcTIMGDXK2/XkFJDY2FjD70mRmZnL99dfn7OPp6UmLFi3Ys2fPv44fHR1NamoqHTt2pFSpUjmPWbNm/et2U140bdo0V/vAgQP079+fqlWr4u/vn3PL69ixY4B5HtLT0//zPIjIpamDskgxtnjxYvbu3Ut2dva/boncdNNNrFy5Mt/HbNq0KW5ubmzZsoVly5ZRpkwZhg4delXHXbBgAR07dsTHx+eK+3l6euZq22y2XNv+vAX299tC+ZGcnAzADz/8QPny5XM95+3tne/jlSxZMle7e/fuVKpUialTpxIeHo7dbqd+/fpkZGQAUKJEiauqW0RMurIjUkxt2bKFO+64g2nTptG+fXueffbZAjmur68vDRo04Ntvv+WNN97g7bffxs3t6v7UfP/99/Ts2bNA6vq7atWq5fSZ+VNmZiYbN26kbt26/9q/bt26eHt7c+zYMapXr57rERERcU21nD9/nn379vHMM8/Qvn176tSpw4ULF3LtU6NGDUqUKMHy5cuv6b1Eiitd2REpho4cOcItt9zCuHHjcm6ftGrVii1bttCkSZNrPn7Lli2ZPHkyPXv25KabbrqqY8TGxrJp0yYWLFhwzfX8U8mSJXnggQd4/PHHCQoKomLFirz22mukpqYyfPjwf+3v5+fH//3f/zF69Gjsdjtt2rQhISGB33//HX9/f4YMGXLVtZQuXZrg4GA+/vhjypUrx7Fjx3jqqady7ePj48OTTz7JE088gZeXF9dffz1nz55l165dOfW+//77zJs3T4FI5BIUdkSKmbi4OLp06ULPnj1zvlQjIyPp2rUr48aNY/Hixf95jBkzZjB06FAMw7jk840aNcLT05PXX3/9qutcuHAhLVq0ICQk5KqPcSWTJk3CbrczePBgkpKSaNasGUuWLKF06dKX3P/FF1+kTJkyTJw4kUOHDhEYGEiTJk0YN27cNdXh5ubG3Llzefjhh6lfvz61atXivffe+1dIfPbZZ/Hw8OC5557j1KlTlCtXjvvvvz/n+XPnzl1V/yGR4sBmXO6vlYgUeytXruT999//13wyzz//PKtWrbps35t27drRpEkT3nzzzXwd9+969OhBmzZteOKJJ666fhER0JUdEbmMDh06sG3bNlJSUqhQoQJff/01rVq1AuCnn37i/fffz7W/3W7n7NmzTJs2jQMHDvD999/n+7h/16ZNG/r371/wH0xEih1d2RGRArFy5UpuvvlmateuzaeffkpkZKTVJYmIAAo7IiIi4uI09FxERERcmsKOiIiIuDSFHREREXFpCjsiIiLi0hR2RERExKUp7IiIiIhLU9gRERERl6awIyIiIi5NYUdERERcmsKOiIiIuLT/B90FZYXenSDAAAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "Tc_K = [190.564, 154.581]\n", "pc_Pa = [4599200, 5042800]\n", "acentric = [0.011, 0.022]\n", "model = teqp.canonical_PR(Tc_K, pc_Pa, acentric)\n", "model1 = teqp.canonical_PR([Tc_K[0]], [pc_Pa[0]], [acentric[0]])\n", "T = 170.0 # [K] # Note: above Tc of the second component\n", "rhoL0, rhoV0 = model1.superanc_rhoLV(T) # start off at pure of the first component\n", "p0 = rhoL0*model1.get_R(np.array([1.0]))*T*(1+model1.get_Ar01(T, rhoL0, np.array([1.0])))\n", "j = model.trace_VLE_isobar_binary(p0, T, np.array([rhoL0, 0]), np.array([rhoV0, 0]))\n", "df = pandas.DataFrame(j) # Now as a data frame\n", "plt.plot(df['xL_0 / mole frac.'], df['T / K'])\n", "plt.plot(df['xV_0 / mole frac.'], df['T / K'])\n", "plt.gca().set(xlabel='$x_1,y_1$ / mole frac.', ylabel='T / K')\n", "plt.show()" ] } ], "metadata": { "celltoolbar": "Raw Cell Format", "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.11.0" } }, "nbformat": 4, "nbformat_minor": 5 }