{ "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": "2025-10-15T23:08:21.569765Z", "iopub.status.busy": "2025-10-15T23:08:21.569628Z", "iopub.status.idle": "2025-10-15T23:08:22.048665Z", "shell.execute_reply": "2025-10-15T23:08:22.048202Z" } }, "outputs": [ { "data": { "text/plain": [ "'0.23.1'" ] }, "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": "2025-10-15T23:08:22.050522Z", "iopub.status.busy": "2025-10-15T23:08:22.050299Z", "iopub.status.idle": "2025-10-15T23:08:22.056191Z", "shell.execute_reply": "2025-10-15T23:08:22.055721Z" } }, "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": "2025-10-15T23:08:22.057547Z", "iopub.status.busy": "2025-10-15T23:08:22.057425Z", "iopub.status.idle": "2025-10-15T23:08:22.062085Z", "shell.execute_reply": "2025-10-15T23:08:22.061550Z" } }, "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": "2025-10-15T23:08:22.063225Z", "iopub.status.busy": "2025-10-15T23:08:22.063108Z", "iopub.status.idle": "2025-10-15T23:08:22.078134Z", "shell.execute_reply": "2025-10-15T23:08:22.077722Z" } }, "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": "2025-10-15T23:08:22.079418Z", "iopub.status.busy": "2025-10-15T23:08:22.079299Z", "iopub.status.idle": "2025-10-15T23:08:22.243590Z", "shell.execute_reply": "2025-10-15T23:08:22.243103Z" } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkAAAAG0CAYAAADacZikAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjcsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvTLEjVAAAAAlwSFlzAAAPYQAAD2EBqD+naQAAZ3dJREFUeJzt3Wd0VNUehvFn0gMkIRDSIBB67yK9iAjSBEEpKohiu4YmVlBEQQ2ioDRBFAgqRemIiCJSpEvvoRNaQk9CQuqc+2E0GkmUAJnJZN7fWlnLPbPPyX/O1cx7z95nb5NhGAYiIiIiDsTJ1gWIiIiIWJsCkIiIiDgcBSARERFxOApAIiIi4nAUgERERMThKACJiIiIw1EAEhEREYfjYusC8iKz2cy5c+fw8vLCZDLZuhwRERG5BYZhEB8fT3BwME5O/36PRwEoC+fOnSMkJMTWZYiIiMhtOH36NCVKlPjXPgpAWfDy8gIsF9Db29vG1YiIiMitiIuLIyQkJON7/N8oAGXhz2Evb29vBSARERE7cyvTVzQJWkRERByOApCIiIg4HAUgERERcTgKQCIiIuJwFIBERETE4SgAiYiIiMNRABIRERGHowAkIiIiDkcBSERERByOApCIiIg4HAUgERERcTgKQCIiIuJwFIBERETEqvavW8SNxASb1qAAJCIiIlaRlJTEhs9eoOqvfdg77X82rcWmASg8PJx69erh5eWFv78/nTt3JjIy8l+P+eKLL2jatCm+vr74+vrSqlUrtm7dmqmPYRi8/fbbBAUF4enpSatWrThy5EhufhQRERH5F4cOR3Lko/tofGEOAGnOnhjmdJvVY9MAtHbtWsLCwti8eTMrV64kNTWV1q1bk5CQ/W2xNWvW0LNnT1avXs2mTZsICQmhdevWnD17NqPP6NGjGT9+PFOmTGHLli0ULFiQNm3akJSUZI2PJSIiIn9INxssXTSbYrNaUT39ANcpwP6mk2j04ueYnJxtVpfJMAzDZr/9Hy5evIi/vz9r166lWbNmt3RMeno6vr6+TJw4kd69e2MYBsHBwbz88su88sorAMTGxhIQEEBERAQ9evS46RzJyckkJydntOPi4ggJCSE2NhZvb++78+FEREQczJkr19kw/Q0ejf8GJ5PBabeyePWeTeESlXLl98XFxeHj43NL3995ag5QbGwsAEWKFLnlYxITE0lNTc045sSJE0RHR9OqVauMPj4+PtSvX59NmzZleY7w8HB8fHwyfkJCQu7gU4iIiDg2wzD4YfNeTo5rR/frX+NkMjge0pUSr6zPtfCTU3kmAJnNZgYNGkTjxo2pVq3aLR/3+uuvExwcnBF4oqOjAQgICMjULyAgIOO9fxoyZAixsbEZP6dPn77NTyEiIuLYYhNTGTN9FrV+7EQT026ScedSq08o03c6JrcCti4vg4utC/hTWFgY+/btY/369bd8zKhRo5g7dy5r1qzBw8Pjtn+3u7s77u7ut328iIiIwMajF9ky5wMGps3E1ZTOVc+SePWajV9wdVuXdpM8cQeoX79+LFu2jNWrV1OiRIlbOubjjz9m1KhR/Pzzz9SoUSPj9cDAQABiYmIy9Y+Jicl4T0RERO6e5LR0Pv5+G9dmPsZL6dMt4ad0B3wHbcQlD4YfsHEAMgyDfv36sWjRIn799VdKly59S8eNHj2akSNHsmLFCu65555M75UuXZrAwEBWrVqV8VpcXBxbtmyhYcOGd7V+ERERRxcZHc+gcd/Q9ffHaOe8lTSTC8kPhOPb+xtw97J1edmy6RBYWFgYs2fPZsmSJXh5eWXM0fHx8cHT0xOA3r17U7x4ccLDwwH48MMPefvtt5k9ezahoaEZxxQqVIhChQphMpkYNGgQ7733HuXLl6d06dIMGzaM4OBgOnfubJPPKSIikt+YzQYzN50kcsUUPnGahodTKjcKBOHZ82tcQurZurz/ZNMANHnyZABatGiR6fUZM2bQp08fAKKionBycsp0TEpKCo888kimY4YPH84777wDwGuvvUZCQgLPPfcc165do0mTJqxYseKO5gmJiIiIxYW4JIZ89zutT37MKJc1ACSHtsSz2zQocOtPcttSnloHKK/IyToCIiIijmTFvmgmL/iJUekfU9kpCjNOmO4biqnpy+Bk26nFOfn+zjNPgYmIiEjelZCcxrvf7ydux0K+cf0cL6cbpHn64fLoNCjTwtbl5ZgCkIiIiPyrHVFXeXXu7zwWN52+bj8CYA5piMujM8A7yMbV3R4FIBEREclSWrqZiauPMv/XLYxzGUddlz82Fm80AKf7h4Oz/cYI+61cREREcs2pywkM+nYXXmfWstR1EkVM1zHcvTE9/DlUamfr8u6YApCIiIhkMAyDedvPMHLpXp4xz6O/2yKcMCCoJqZHZ0KRW1uzL69TABIREREAriakMGThXn7fH8lnrpNo6rLP8kbdp+DBUeCaf5aTUQASERER1h+5xMvzdlEifg8/uI8n0HQVw7UApg6fQs3uti7vrlMAEhERcWBJqel89FMk09Yf5xnn5Qxxn4MzZvCrgKnbV+Bf2dYl5goFIBEREQcVGR3PwLk7ORsdwxTXz3nQ+XfLG9UegY7jwL2QbQvMRQpAIiIiDsZsNojYeJJRKw5RPv04P3qMowQx4OwGD4bDPX3BZLJ1mblKAUhERMSBXIhL4uV5u/ntyEV6Ov/KCPevcCUVCpeER2dC8Tq2LtEqFIBEREQcxE/7o3ljwR6SEuP51G0GnZ1+s7xRoS08PBk8fW1boBUpAImIiORzCclpjFx2gLm/n6as6SzTCk4gND0KTM5w/9vQaIDNNzK1NgUgERGRfGz36WsMnLuTk5cT6ei8kTHu03BLvwGFAuCRGRDa2NYl2oQCkIiISD6UbjaYvOYon/5yBCdzCh8XmMsj5h/BDIQ2ha7TwCvA1mXajAKQiIhIPnP6SiKDv9vF7yevUsJ0kW98JhOafMjyZtOXocVQu97I9G5w7E8vIiKSzyzeeZZhi/cRn5xGO/c9fOr2GW7JcZYJzg9PhQqtbV1inqAAJCIikg/E3khl2OJ9LN19DmfSGVNkGV0Tv4VUoHhdeDTC8qi7AApAIiIidm/L8csM/m43Z6/dIMAplgXFvqBE7A7Lm/c+D63fAxc32xaZxygAiYiI2KmUNDOf/nKYyWuPYRjQqfAxPjaNxzX2IrgVgocmQLUuti4zT1IAEhERsUPHL15n0Le72HMmFhNmJpdax4MXvsRkmMG/CnT7CvzK27rMPEsBSERExI4YhsHc308z4vsD3EhNJ8QjiQWBM/GPXmvpUPMxaD8G3ArYttA8TgFIRETETlxJSOH1BXtYeSAGgF4lL/JO0sc4R58GFw9o9xHU7pXvNzK9GxSARERE7MC6wxd5ed5uLsYn4+oMEVV30+joWEzmVPAtbRnyCqph6zLthgKQiIhIHpaUms7oFZFM33ACgOrFnPim2Df4HF5m6VC5I3SaBB4+NqzS/igAiYiI5FGR0fEMnLuTQ9HxALxSK40XL4zE6fgRcHKBB0ZCg/9pyOs2KACJiIjkMYZhELHxJOE/HiIlzUzRgm58fc9Rqmx/F9JugHdxy0amJevbulS7pQAkIiKSh1yIT+LVeXtYe/giAK3LezPOZw6eW2ZZOpS9H7p8AQWL2rBK+6cAJCIikkf8ciCG1xbs4UpCCu4uToxqUYDOR97AtG8fYIL7hkLTV8DJydal2j0FIBERERu7kZLO+8sP8M3mKAAqBXoxvf55gtc8B8lxUMAPun4JZe+zcaX5hwKQiIiIDe07G8ugb3dx9MJ1AJ5rVILXXGbj8tMUS4eSDeGR6eAdbMMq8x8FIBERERswmw2+XH+cj36KJDXdwN/LnQnti1F/2ytwZqulU6MBcP/b4Oxq22LzIQUgERERK4uOTeLlebvYcPQyAA9UCWBsnYt4/dAZblyxrOnTeTJUam/bQvMxBSARERErWrEvmjcW7uFaYiqers4M71CR7gmzMM3/GDAgqKZlVWffUFuXmq8pAImIiFhBYkoaI74/wNzfTwNQvbgPEx4qTuiagXDij41M73ka2oSDq4cNK3UMCkAiIiK5bM+Zawyau4vjlxIwmeCF5mUZXOESrvMehOvR4FoQOo6DGo/aulSHoQAkIiKSS9LNBlPXHWfMz5GkmQ2CfDwY+2gNGkbPgq9HgJEOfhWh+9dQrKKty3UoCkAiIiK54Ny1Gwz+bhebj18BoF31QMLbhuCzoj8c/tHSqXo36PAJuBeyYaWOyaZLSYaHh1OvXj28vLzw9/enc+fOREZG/usx+/fvp2vXroSGhmIymfj0009v6vPOO+9gMpky/VSqVCmXPoWIiEhmy/eep+2439h8/AoF3JwZ/UgNJrUw4fPV/Zbw4+xmCT5dpir82IhNA9DatWsJCwtj8+bNrFy5ktTUVFq3bk1CQkK2xyQmJlKmTBlGjRpFYGBgtv2qVq3K+fPnM37Wr1+fGx9BREQkQ0JyGq/O282Ls3YQeyOVmiV8WN6/Cd3MKzBNbwPXoqBwKei70jLhWbu424xNh8BWrFiRqR0REYG/vz/bt2+nWbNmWR5Tr1496tWrB8Abb7yR7bldXFz+NSCJiIjcTbtOX2PQ3J2cvJyIyQRhLcoxsGkgrssHwL4Flk6VOkCnSeBZ2Ka1Sh6bAxQbGwtAkSJF7vhcR44cITg4GA8PDxo2bEh4eDglS5bMsm9ycjLJyckZ7bi4uDv+/SIi4hjSzQZT1h7jk5WHSTMbBPt48En3WtQvGAPT7ofLR8DkDA+8Cw376a5PHpFntpM1m80MGjSIxo0bU61atTs6V/369YmIiGDFihVMnjyZEydO0LRpU+Lj47PsHx4ejo+PT8ZPSEjIHf1+ERFxDGev3aDnF5v56CfLU17tawTx48Bm1I/7Gb5oaQk/XsHw1HJo1F/hJw/JM3eAwsLC2Ldv312Zq9O2bduMf65Rowb169enVKlSfPfdd/Tt2/em/kOGDGHw4MEZ7bi4OIUgERH5V8v2nGPowr3EJaVR0M2ZdztVo2v1Iph+fAl2fm3pVOY+yy7uBf1sW6zcJE8EoH79+rFs2TLWrVtHiRIl7vr5CxcuTIUKFTh69GiW77u7u+Pu7n7Xf6+IiOQ/15PTeGfpfuZvPwNAzZDCjOtei1BTNExrDTF7ARM0fx2avwZOzrYtWLJk0wBkGAb9+/dn0aJFrFmzhtKlS+fK77l+/TrHjh2jV69euXJ+ERFxDLtOX2Pg3J2c+vtE51blcY38HhaHQUo8FPCDrl9A2Za2Llf+hU0DUFhYGLNnz2bJkiV4eXkRHR0NgI+PD56engD07t2b4sWLEx4eDkBKSgoHDhzI+OezZ8+ya9cuChUqRLly5QB45ZVX6NixI6VKleLcuXMMHz4cZ2dnevbsaYNPKSIi9i7dbDB5zVE++eUI6WaD4oU9GdutJvVLesHKN2HzZ5aOIQ3g0RngHWzbguU/mQzDMGz2y7OZDDZjxgz69OkDQIsWLQgNDSUiIgKAkydPZnmnqHnz5qxZswaAHj16sG7dOi5fvkyxYsVo0qQJ77//PmXLlr2luuLi4vDx8SE2NhZvb+8cfy4REck/zl67wUvf7mLrCcuKzh1qBPH+w9XxSYmBeX3gzO+Wjo36w/3DwdnVdsU6uJx8f9s0AOVVCkAiIgLZTHSuUxzT0VWw8Fm4cQU8fKDzZKjU3tblOrycfH/niUnQIiIieck/JzrXCinMuB61KOXrAavfh3UfAwYE1YJuM8E31Jblym1QABIREfmbbCc6J16Er7rDyd8sHes9A63fB1cP2xYst0UBSEREhH9Z0blMUTi5HuY/DddjwLUgPDQeqj9i65LlDigAiYiIwzv3x0TnLX9MdG5fI4gPOlfHx8MZfhsLv44EwwzFKkO3r6BYBRtXLHdKAUhERBzaj3vP88bCvcTeSKWAmzPvPlSVR+qWwHTjKsx5AY78ZOlYowd0GAtuBW1bsNwVCkAiIuKQElPSeHfpAb7ddhqAGiV8GNejNqX9CsKZ7TDvSYg9Dc7u0O4jqNNbe3nlIwpAIiLicPaeiWXg3J0cv5SAyQT/a16Wlx6ogKuTCbZ8Dj+9CeZU8C1tGfIKqmHrkuUuUwASERGHYTYbTP3tOGN+jiQ13SDQ2zLRuWHZopAUB0v7w4HFls6VO0KnSZZ1fiTfUQASERGHEBOXxODvdrHh6GUA2lYLJLxLdQoXcIPoffBdb7hyDJxc4IGR0OB/GvLKxxSAREQk3/t5fzSvL9jD1cRUPF2dGd6xCt3rhVi2ZNrxNSx/BdKSwLs4PBoBIffaumTJZQpAIiKSb91ISee9Hw4wa0sUAFWDvRnfszZlixWClERL8Nk1y9K53APw8OdQsKgNKxZrUQASEZF86cC5OAbM3cnRC9cBeK5ZGV5uXQF3F2e4dMQy5HXhAJic4L43oclgcHKycdViLQpAIiKSrxiGwYwNJxn14yFS0s34e7kzpltNmpYvZumwbwEsHQAp16GgPzwyDUo3s23RYnUKQCIikm9cjE/m1fm7WRN5EYBWlf35sGsNihZyh7Rk+Gko/P6lpXOpJpbw4xVow4rFVhSAREQkX1gdeYFX5+3m0vUU3F2ceKt9ZZ5oUMoy0fnqSZjXB87ttHRu+jK0GArO+hp0VPpfXkRE7FpSajofrjjEjA0nAagU6MX4nrWpEOBl6XBoOSx+AZJiwdMXunwB5R+wXcGSJygAiYiI3ToSE0//OTs5FB0PQJ9GobzRthIers6QngqrRsDG8ZbOxe+xPOJeOMR2BUueoQAkIiJ2xzAMZm+NYuSyAySlmila0I2PHq1By0oBlg5x52D+0xC1ydJu8CK0ehdc3GxXtOQpCkAiImJXriak8PqCPfx8IAaApuX9GNOtJv5eHpYOx36FBc9A4mVw94ZOE6FKJxtWLHmRApCIiNiNjccu8dK3u4iJS8bV2cTrD1bi6calcXIygTkd1o6GtR8CBgRWh0dnQtGyti5b8iAFIBERyfNS0818svIwk9cewzCgTLGCjO9Rm2rF/9io9PpFWPgMHF9jadftAw+OAldPW5UseZwCkIiI5GmnLicwYM5Odp+JBaBHvRDe7liFAm5/fIWd2miZ7xN/HlwLQIdPoWZ32xUsdkEBSERE8iTDMFi08yzDFu8jISUdbw8XRnWtQbvqQZYOZrPlCa9VI8BIB7+K0O0r8K9k28LFLigAiYhInhOXlMqwxftYsuscAPeWLsKn3WsRXPiPIa3EK7D4RTj8o6VdvRt0+ATcC9moYrE3CkAiIpKnbD91lUHf7uT0lRs4O5kYdH95XryvHM5OJkuHs9vhuz4QGwXO7tD2Q8ucH5PJlmWLnVEAEhGRPCHdbPDZ6qN8uuoI6WaDEr6ejOtRm7qlfC0dDMOyj9eKIWBOBd/S0G0mBNW0beFilxSARETE5s5du8FL3+5iy4krADxUM5j3Hq6Gt4erpUNSHHw/APYvsrQrdYDOn4GHj40qFnunACQiIja1Yl80ry/YQ+yNVAq6OTOiUzW61Clu2cQUIHovfPckXDkGTi7Q+j2o/4KGvOSOKACJiIhN3EhJZ+QPB5i9JQqAGiV8GN+jNqF+BS0dDAN2fg3LX4W0JPAuYdnLK6Se7YqWfEMBSERErO7AuTgGzN3J0QvXAXi+eRlefqAibi5Olg4pCfDDy7B7jqVdvjU8/DkUKGKjiiW/UQASERGrMQyDiI0nCf/xEClpZop5ufNJt1o0Ke/3V6eLkZYhr4sHweQELd+Cxi+Bk5PtCpd8RwFIRESs4vL1ZF6dv4dfD10A4P5K/ox+pAZFC7n/1WnPPPh+IKQmQKEAeGQ6hDaxUcWSnykAiYhIrvvtyEUGf7ebi/HJuLk48Wa7yvRuWOqvic6pSbDiDdg+w9Iu3Qy6ToNC/rYrWvI1BSAREck1KWlmxvwcyefrjgNQ3r8Q43vWpnKQ91+drhyH73pbnvbCBM1ehRZvgJOzbYoWh6AAJCIiueLEpQQGzt3Jnj82MX28fkneal8FT7e/BZsDS2FJGCTHQYGi0GUqlGtlo4rFkSgAiYjIXWUYBgt3nGXYkn0kpqTj4+nKh11r8GC1wL86paXAyrdhy2RLO6SBZb6PT3HbFC0ORwFIRETumvikVN762yam9UsX4dMetQjy8fyr07XTMK8PnN1maTfqD/cPB2dX6xcsDksBSERE7oqdUVcZOHcXUVcSs97EFODwT7Doebhx1bKNRecpUKmd7YoWh6UAJCIid8RsNpiy7hhjfz5MmtmgeGFPxvf82yamAOlpsPo9WP+JpR1cGx6dCb6lbFO0ODybrioVHh5OvXr18PLywt/fn86dOxMZGfmvx+zfv5+uXbsSGhqKyWTi008/zbLfpEmTCA0NxcPDg/r167N169Zc+AQiIo4tJi6JJ6ZtYfSKSNLMBh1qBLF8YNPM4SfuPMzs+Ff4ufd5ePonhR+xKZsGoLVr1xIWFsbmzZtZuXIlqamptG7dmoSEhGyPSUxMpEyZMowaNYrAwMAs+3z77bcMHjyY4cOHs2PHDmrWrEmbNm24cOFCbn0UERGHs+pgDA9+uo6Nxy7j6erM6K41mNCzNj6ef5vLc2w1TGkCURvBzcuyl1e70eDinu15RazBZBiGYesi/nTx4kX8/f1Zu3YtzZo1+8/+oaGhDBo0iEGDBmV6vX79+tSrV4+JEycCYDabCQkJoX///rzxxhs3nSc5OZnk5OSMdlxcHCEhIcTGxuLt7X1TfxERR5aUms6oHw8RsfEkAFWCvJnwWG3KFiv0VydzOqz7CNaMAgwIqA7dZkLRsjapWRxDXFwcPj4+t/T9nac2VomNtawVUaTI7W92l5KSwvbt22nV6q91JJycnGjVqhWbNm3K8pjw8HB8fHwyfkJCQm7794uI5GdHYuLpPGlDRvjp26Q0i8IaZQ4/1y/CN11gTThgQJ3e8MxKhR/JU/JMADKbzQwaNIjGjRtTrVq12z7PpUuXSE9PJyAgINPrAQEBREdHZ3nMkCFDiI2Nzfg5ffr0bf9+EZH8yDAM5myNouPE9RyKjqdoQTdmPFWPYR2q4O7yt4UNT22Ez5vC8TXgWsDylNdDE8DVM9tzi9hCnnkKLCwsjH379rF+/Xqr/253d3fc3TUeLSKSldjEVN5YuIcf91n+T2TT8n6M6VYTfy+PvzqZzbBxHKwaCUY6+FW0DHn5V7ZR1SL/Lk8EoH79+rFs2TLWrVtHiRIl7uhcfn5+ODs7ExMTk+n1mJiYbCdNi4hI1n4/eYWBc3ZyLjYJV2cTr7apyDNNyuD097V9Eq/AohfgyE+Wdo3u0H4suBfK+qQieYBNh8AMw6Bfv34sWrSIX3/9ldKlS9/xOd3c3Khbty6rVq3KeM1sNrNq1SoaNmx4x+cXEXEE6WaDcb8cofvnmzgXm0Ro0QIs+F8jnmtWNnP4ObMNPm9mCT/O7tBxHDz8ucKP5Hk2vQMUFhbG7NmzWbJkCV5eXhlzdHx8fPD0tIwX9+7dm+LFixMeHg5YJjkfOHAg45/Pnj3Lrl27KFSoEOXKlQNg8ODBPPnkk9xzzz3ce++9fPrppyQkJPDUU0/Z4FOKiNiXc9duMOjbXWw9cQWALnWKM6JTNQq5/+0rwzBg82TLfl7mVChSxrKwYVANG1UtkjM2fQzeZDJl+fqMGTPo06cPAC1atCA0NJSIiAgATp48meWdoubNm7NmzZqM9sSJE/noo4+Ijo6mVq1ajB8/nvr1699SXTl5jE5EJD9ZsS+a1xfsIfZGKgXdnHnv4Wo8XPsfUxNuXLPs4H5omaVdpbNlorOH/l6KbeXk+ztPrQOUVygAiYijSUpN570fDvDN5igAapTwYXyP2oT6Fczc8dwumPckXD0JTq7Q5gO491nI5v/QilhTTr6/88QkaBERsZ3DMfH0n72TyJh4AJ5vVoaXW1fEzeVv00QNA7ZNgxVDID0FCpe0rOpcvK5tiha5QwpAIiIOyjAMZm+NYsT3B0hOM+NXyJ2x3WrSrEKxzB2T4+H7gbBvgaVdsR10/gw8fW8+qYidUAASEXFA1xJTeGPBXlbstzx80qxCMcY8WpNiXv9YEy16n2XI6/JRcHKBVu9Aw34a8hK7pwAkIuJgtp64wqC5f63t81qbSvRtUjrz4+2GATu/geWvQFoSeBeHR2ZAyVt7mEQkr1MAEhFxEOlmg4m/HmXcqsOYDQgtWoDxPWtTo0ThzB1TEuCHV2D3bEu7XCt4eCoULGr1mkVyiwKQiIgDuKW1fQAuRsJ3veHiITA5Qcu3oPFL4JRnto4UuSsUgERE8rmf90fz6vz/WNsHYPe3sGwQpCZCoQB4ZDqENrF6vSLWoAAkIpJPJaWm88Hyg3y16RTwL2v7pN6AH1+HHTMt7dLNoeuXUMjfyhWLWI8CkIhIPnT0Qjz9Zu/kULRlbZ/nmpXhlX+u7QNw6SjM6wMxewETtHgDmr0KTs5Wr1nEmhSARETyEcMw+G7bad5ZeoAbqekULejGmG41aVExi7s5+xbC0gGQEg8Fi0GXL6DsfdYvWsQGFIBERPKJuKRUhi7cy7I95wFoUs6Psd1r4u/lkbljWjL8NBR+/9LSLtUYuk4D7yArVyxiOwpAIiL5wM6oqwyYu5PTV27g4mTi5dYVeb5Zmcxr+wBcOWEZ8jq/y9JuMhjuexOc9XUgjkX/xouI2DGz2WDKumOM/fkwaWaDEr6ejO9Zmzols9im4uD3sDgMkmPBswh0mQrlH7B+0SJ5gAKQiIiduhCfxOBvd7P+6CUAOtQI4oMu1fH2cM3cMS0FfhkOmz+ztEPqWx5x98niUXgRB6EAJCJih9ZEXuDl73ZzOSEFD1cn3n2oKt3uCcH0zz26rkXBvKfg7DZLu1F/uH84OLvefFIRB6IAJCJiR1LSzHz00yG++O0EAJUCvZj4WG3K+Xvd3DlyBSx6HpKugUdheHgKVGxr1XpF8ioFIBERO3HyUgID5u5kz5lYAHo3LMXQdpXxcP3Hmj3pqfDrSNgwztIuXteykalvKStXLJJ3KQCJiNiBJbvOMnThXhJS0vHxdGX0IzVoUzXw5o6xZ2FBX4jaZGnX/x88MAJc3KxbsEgepwAkIpKHJSSnMXzpfuZvPwPAvaFF+LRHLYILe97c+cgvsOg5SLwM7t7QaSJU6WTlikXsgwKQiEgedeBcHP3m7OD4xQScTNCvZXkGtCyHi/M/trNIT4M1H8BvYyztwBrQbSYUKWP9okXshAKQiEgeYxgGX28+xXs/HCQlzUyAtzufdq9Nw7JFb+4cdx4WPAOn1lva9Z6B1u+Dq8fNfUUkgwKQiEgeci0xhdfm7+HnAzEA3F/Jn48erUmRglnM4Tm2GhY+CwkXwc0LHhoH1bpauWIR+6QAJCKSR/x+8goD5+zkXGwSrs4mhrStzFONQ29e28ecDmtHw9oPAQMCqluGvIqWtUndIvZIAUhExMbSzQafrT7KJ78cxmxAab+CTOhZm2rFfW7uHB8DC5+BE+ss7bp94MFR4JrFpGgRyZYCkIiIDcXEJTFo7i42Hb8MQJfaxRnRuRqF3LP483xinWW+z/UYcC0IHT+FGt2sW7BIPqEAJCJiI6sPXeDlebu5kpBCATdnRnaqRte6WezPZTZbnvBa8wEYZvCvAo/OhGIVrF+0SD6hACQiYmUpaWZGrzjEl+st21lUCfJm4mO1KVOs0M2dr1+0THQ+vtrSrv0EtP0I3ApYsWKR/EcBSETEiv65nUWfRqEMaVcJdxfnLDpvsKzqHH8eXDyhw1io9ZiVKxbJnxSARESsZMmus7y5aB/Xk9MoXMCVjx6pyQNVAm7uaDbDhk/g1/fBSAe/ipanvPwrW79okXxKAUhEJJclpqTxztL9fLftr+0sxvWsRZBPFk9uJVy27OB+dKWlXaOH5c6PW0ErViyS/ykAiYjkooPn4+g3ewfH/tjOon/L8vTPajsLgKjNMP9piDsLLh7Q7iOo3Qv+uQ6QiNwxBSARkVxgGAaztkQxYtmB/97OwmyGjeNh1QjLkFfRcpanvAKrWb9wEQehACQicpfF3kjljQV7+HFfNAAtK/nzcXbbWSRegUUvwJGfLO3qj0KHT8Ddy4oVizgeBSARkbtoR9RV+s/eydlrN3B1NvH6g5Xo26T0zdtZAJzeCvOegrgz4OwObT+0rOysIS+RXKcAJCJyF5jNBlN/O87HP0WSZjYoWaQAE3rWpmZI4Zs7GwZsmgi/vAPmNChS1vKUV2B1a5ct4rBuOwAdOHCAqKgoUlJSMr3+0EMP3XFRIiL25NL1ZAZ/t5t1hy8C0KFGEB90qY63h+vNnROvwOIX4fCPlna1rtBxnIa8RKwsxwHo+PHjPPzww+zduxeTyYRhGAAZt3fT09PvboUiInnYhqOXGPTtLi7GJ+Ph6sQ7HavSvV5I1kNeZ7ZZhrxioyxDXg+Gwz1Pa8hLxAayeA7z3w0cOJDSpUtz4cIFChQowP79+1m3bh333HMPa9asyYUSRUTynrR0Mx//FMkT07ZwMT6ZCgGFWNqvCT3uLXlz+DEM2DQJprexhJ8iZeCZlVCvr8KPiI3kOABt2rSJESNG4Ofnh5OTE05OTjRp0oTw8HAGDBiQo3OFh4dTr149vLy88Pf3p3PnzkRGRv7ncfPmzaNSpUp4eHhQvXp1li9fnun9Pn36YDKZMv08+OCDOapNRCQ7567doMfUzUxcfRTDgJ73lmRJWBMqBGQxjHXjKsx9HH4aapnvU6UzPLcWgmpavW4R+UuOA1B6ejpeXpb/yP38/Dh37hwApUqVuqXw8ndr164lLCyMzZs3s3LlSlJTU2ndujUJCQnZHrNx40Z69uxJ37592blzJ507d6Zz587s27cvU78HH3yQ8+fPZ/zMmTMnh59URORmKw/E0Hbcb2w7dRUvdxcm9KxNeJfqeLplsZfXme3weTOI/AGc3aDdx/BoBHh4W71uEcksx3OAqlWrxu7duyldujT169dn9OjRuLm5MXXqVMqUKZOjc61YsSJTOyIiAn9/f7Zv306zZs2yPGbcuHE8+OCDvPrqqwCMHDmSlStXMnHiRKZMmZLRz93dncDAwBx+OhGRrCWnpRO+/BARG08CULOEDxN61qFk0Sx2ZTcM2DIFfh4G5lTwDbUsbBhcy5oli8i/yHEAeuuttzLu0IwYMYIOHTrQtGlTihYtyrfffntHxcTGWnZHLlKkSLZ9Nm3axODBgzO91qZNGxYvXpzptTVr1uDv74+vry8tW7bkvffeo2jRLFZgBZKTk0lOTs5ox8XF3eYnEJH86MSlBPrN3sH+c5a/Dc82Lc2rbSrh5pLFTfQb12BpPzj4vaVd+SHoNBE8fKxXsIj8pxwHoBYtWpCWlgZAuXLlOHToEFeuXMHX1zfrpx5ukdlsZtCgQTRu3Jhq1bJf/j06OpqAgMy7JwcEBBAdHZ3RfvDBB+nSpQulS5fm2LFjDB06lLZt27Jp0yacnW++TR0eHs67775727WLSP61eOdZ3ly0l4SUdIoUdGPMozW5r5J/1p3P7oB5feDaKXByhTYfwL3PaqKzSB50ywHo4sWL9O7dm19++QWz2Uy9evX45ptvKFeu3L/esblVYWFh7Nu3j/Xr19/xuXr06JHxz9WrV6dGjRqULVuWNWvWcP/999/Uf8iQIZnuKsXFxRESEnLHdYiI/UpMSWP4kv3M227Zwb1BmSJ82r02gT4eN3c2DNg6FX560zLkVbiUZa5P8TrWLVpEbtktB6DXX3+dXbt2MWLECDw8PPj888959tlnWb169R0X0a9fP5YtW8a6desoUaLEv/YNDAwkJiYm02sxMTH/Ot+nTJky+Pn5cfTo0SwDkLu7O+7u7rdXvIjkO5HR8YTN3sHRC9dxMsGA+8vTv2V5nJ2yuJOTFAtL+sHBpZZ25Y7w0ETwLGzVmkUkZ245AK1cuZKIiAjatGkDQIcOHahcuTLJycm3HR4Mw6B///4sWrSINWvWULp06f88pmHDhqxatYpBgwZlqq1hw4bZHnPmzBkuX75MUFDQbdUpIo7BMAzm/n6ad5buJznNjL+XO+N6ZLODO1iGvOY/BVdP/jHk9T7c+5yGvETswC0HoHPnzlGz5l/rVpQvXx53d3fOnz9PaGjobf3ysLAwZs+ezZIlS/Dy8sqYx+Pj44OnpycAvXv3pnjx4oSHhwOWhRibN2/OmDFjaN++PXPnzmXbtm1MnToVgOvXr/Puu+/StWtXAgMDOXbsGK+99hrlypXLCG8iIv8Un5TK0EX7+H63ZWmP5hWKMaZbTfwKZfF/8LIc8poBxetauWoRuV05mgT9zwnEzs7OGVth3I7JkycDlonVfzdjxgz69OkDQFRUFE5Ofz1p0ahRI2bPns1bb73F0KFDKV++PIsXL86YOO3s7MyePXuYOXMm165dIzg4mNatWzNy5EgNc4lIlvaeiaXfnB2cupyIi5OJV9tU5NmmZXC6lSGvSh2g0yQNeYnYGZNxiwnGyckJHx+fTE96Xbt2DW9v70wB5cqVK3e/SiuLi4vDx8eH2NhYvL21YJlIfmUYBhEbT/LB8oOkphsUL+zJ+J61qVvKN+sDzu20POX155BX6/eg/vMa8hLJI3Ly/X3Ld4BmzJhxx4WJiOQV1xJTeHX+HlYesDxU0aZqAKO71sSnQBY7uBsGbP0Cfn4T0lOgcMk/nvLSkJeIvbrlAPTkk0/mZh0iIlaz/dQV+s/eybnYJNycnXizfWV6NyyV9VpmGvISyZdyvBCiiIi9MpsNJq89xtiVh0k3G4QWLcDEx+pQrXg2qzTfNOQ1Euq/oCEvkXzglgPQre7zdfz48dsuRkQkt1yMT2bwd7v47cglADrVCub9h6tTyD2LP4NZDXk9EgElNOQlkl/ccgA6efIkpUqV4rHHHsPfP5tl4EVE8qCNRy8x8NtdXIxPxsPViREPVePRe0rkYMhrInhmMzFaROzSLQegb7/9lunTpzN27Fjatm3L008/Tbt27TI9ASYikpekpZsZv+oIE1YfxTCgvH8hJj1ehwoBXlkfoCEvEYdxy4/B/+ns2bNEREQQERFBYmIivXr1om/fvpQvXz63arQ6PQYvYv+iY5MYMHcnW09YluboUS+E4R2r4ul284bIGAb8/iX8NNQy5OXzx1NeGvISsSs5+f7OcQD6u7Vr1/LOO++wbt06Ll26hK9v/rhFrAAkYt9WH7rA4O92cTUxlYJuznzQpTqdahXPuvM/h7wqtofOkzTkJWKHcmUdoL9LSkpi/vz5TJ8+nS1btvDoo49SoECB2ypWRORuSU038/FPkXy+zvIwRpUgbyY9XofSfgWzPuCfe3k98C40eFFDXiIOIEcBaMuWLUybNo3vvvuOMmXK8PTTT7NgwYJ8c+dHROzX2Ws36D97BzuirgHQu2EphrarjIdrNkNemfby0lNeIo7mlgNQ1apVuXDhAo899hhr167NtDGqiIgtrTwQwyvzdhN7IxUvDxdGd61B2+pBWXe+cQ2W9oOD31vaespLxCHlaC+wggUL4uLikvWjo3/QXmAiYi0paWZG/XiI6RtOAFCzhA8TetahZNFshuTPbod5T8G1U9rLSyQf0l5gIpLvRV1OpN+cHew5EwtA3yalef3BSri5ZLE0h2HAlinw87A/hrxKwaMztJeXiAPTXmAiYnd+3Hue1+bvIT45DR9PVz5+tCYPVAnIuvONq5anvA4ts7Qrd4SHJmovLxEHp73ARMRuJKWm88Hyg3y16RQAdUoWZnzP2pTwzWbI68x2mN8HrkWBsxu0fh/ufVZDXiKiACQi9uHEpQT6zd7B/nNxADzfvAyvtK6Iq3M2Q16bJ8PKty1DXr6hloUNg2tbtWYRybsUgEQkz1uy6yxDF+4lISWdIgXdGNOtJvdVzGZPwhtXYXEYRP5gaVfpBA9NAI9sdnwXEYekACQieVZSajrvfr+fOVtPA3BvaBHG96xNoI9H1gec2WZ5yiv2jyGvNh9AvWc05CUiN7nlANS0aVM6derEQw89RIUKFXKzJhERjl6IJ2zWTiJj4jGZoN995Rh4f3lcshvy2jQJfhkO5jTwLf3HkFcta5ctInbilrdyf/bZZ9m0aRN169alcuXKvP7662zYsIE72EpMRCRLC7afoeOEDUTGxONXyJ2vn67Py60rZh1+Eq/A3Mfg5zct4afqw/D8OoUfEflXOd4MNTk5mVWrVrFkyRK+//570tPTad++PQ899BBt2rTB09Mzt2q1Gi2EKGIbiSlpvL1kP/O3nwGgUdmifNqjFv5e2Qx5nd4K85+G2NPg7A4PhsM9T2vIS8RBWW03eLDsD7Z06VKWLl3KsWPHaNmyJUOGDKFx48Z3clqbUgASsb4jMfG8OGsHRy5cx8kEA++vQL+W5XB2yiLMmM2waQKsGmG561OkjGXIK0hb9Ig4MqsGoL87duwYS5cuJSQkhEceeeRundbqFIBErGvettO8vWQ/N1LTKeblzvgetWlYtmjWnRMuw+IX4MjPlna1rtDhU/DQf6sijs5mASi/UAASsY7ElDSGLd7Pgh2WIa+m5f0Y260Wxbzcsz7g1CbLkFf8OcuQV9sPoW4fDXmJCJBLe4GJiNxNh/8Y8jr6x5DX4Acq8GKLcjhlN+S14RP49X0w0qFoOXh0JgRWs37hIpIvKACJiFUZhsG87Wd4e8k+klLN+Hu5M75nbRqUyWbI6/pFWPQ8HFtladfoDu3Hgnsh6xUtIvmOApCIWE1CchrDFu9j4c6zgGXI65PutfArlM2Q18n1ML8vXI8GF09o9xHUfkJDXiJyx+4oAP05fcikP0Yi8h8io+N5cdZ2jl1MwMkEL7euyP+al81myCsdfhsDa8LBMINfRctTXgFVrF63iORPt7wQ4t9NmzaNatWq4eHhgYeHB9WqVePLL7+827WJSD5gGAbf/h5Fp0nrOXYxgQBvd+Y824Cw+7KZ73P9Anz9MKx+3xJ+aj0Oz61W+BGRuyrHd4Defvttxo4dS//+/WnYsCEAmzZt4qWXXiIqKooRI0bc9SJFxD4lJKfx1uJ9LPpjyKtZhWJ80q0mRbMb8jq+FhY8AwkXwLWAZa5PrZ5WrFhEHEWOH4MvVqwY48ePp2fPzH+U5syZQ//+/bl06dJdLdAW9Bi8yJ07FB3Hi7N2cPxiAs5OJl5uXYEXmv3LkNfaD2HtaMAA/yqWIa9iFa1dtojYsVx9DD41NZV77rnnptfr1q1LWlpaTk8nIvmMZcjrNMOX7ic5zUygtwcTHqtNvdAiWR8QH22563PyN0u7Tm948ENwK2C9okXE4eR4DlCvXr2YPHnyTa9PnTqVxx9//K4UJSL26XpyGoO+3cUbC/eSnGamRcViLB/YNPvwc+xXmNLEEn7cCkGXL+GhCQo/IpLrbuspsGnTpvHzzz/ToEEDwLIfWFRUFL1792bw4MEZ/caOHXt3qhSRPO/g+TjCZu3g+CXLkNcrrSvyfLMyWQ95padZnvD6bQxgQEB1y5CXXzlrly0iDirHAWjfvn3UqVMHsOz9BeDn54efnx/79u3L6KdH40Ucwz+HvIJ8PJjQszb3ZHfXJ/asZcgraqOlfc/T0OYDcPW0XtEi4vByHIBWr16dG3WIiB3651NeLSoWY2y3WhQp6Jb1AYd/tqzqfOMKuHnBQ+Msm5mKiFiZVoIWkdvy94UN/3vIKxVWjYCN4y3toJrwyAwoWta6RYuI/EEBSERy7LttpzP28vrPp7yuRVl2cD/zu6V97/PQeiS4ZLMWkIiIFSgAicgtS0xJY9ji/SzYcQa4hYUNDy6DJS9CUix4+ECnSVC5oxUrFhHJ2m1thXG3hIeHU69ePby8vPD396dz585ERkb+53Hz5s2jUqVKeHh4UL16dZYvX57pfcMwePvttwkKCsLT05NWrVpx5MiR3PoYIg7hSEw8nSZuYMGOMziZ4NU2FYnoUy/r8JOWDD++Dt8+bgk/xevC878p/IhInmHTALR27VrCwsLYvHkzK1euJDU1ldatW5OQkJDtMRs3bqRnz5707duXnTt30rlzZzp37pzpCbTRo0czfvx4pkyZwpYtWyhYsCBt2rQhKSnJGh9LJN9ZsP0MD03cwJEL1/H3cmf2v+3ldeU4TGsNW6ZY2o36w1MrwLeUdYsWEfkXOd4KIzddvHgRf39/1q5dS7NmzbLs0717dxISEli2bFnGaw0aNKBWrVpMmTIFwzAIDg7m5Zdf5pVXXgEgNjaWgIAAIiIi6NGjx03nTE5OJjk5OaMdFxdHSEiItsIQh3cjJZ3hS/fx3TbLkFfT8n580r0WftkNee1bCEsHQEo8eBaBh6dAhTZWrFhEHFlOtsKw6R2gf4qNjQWgSJFsJlNi2Xi1VatWmV5r06YNmzZtAuDEiRNER0dn6uPj40P9+vUz+vxTeHg4Pj4+GT8hISF3+lFE7N7RC9fpPGkD322zDHkNfqACEU/dm3X4Sb0By16C+U9Zwk/JhvDCeoUfEcmz8swkaLPZzKBBg2jcuDHVqlXLtl90dDQBAQGZXgsICCA6Ojrj/T9fy67PPw0ZMiTTCtZ/3gEScVSLd55l6KK9JKak41fInfE9a9GorF/WnS8dgXl9IGYfYIKmg6HFUHDOM39eRERukmf+QoWFhbFv3z7Wr19v9d/t7u6Ou7seyRVJSk3n3e/3M2fraQAalS3Kpz1q4e/lkfUBu7+13PlJTYACftBlKpS734oVi4jcnjwRgPr168eyZctYt24dJUqU+Ne+gYGBxMTEZHotJiaGwMDAjPf/fC0oKChTn1q1at3dwkXykeMXr/PirB0cio7HZIIBLcsz4P7yOGc10TklAZa/Bru+sbRDm0LXL8Er0LpFi4jcJpvOATIMg379+rFo0SJ+/fVXSpcu/Z/HNGzYkFWrVmV6beXKlTRs2BCA0qVLExgYmKlPXFwcW7ZsyegjIpkt3X2OjhPWcyg6Hr9CbnzTtz4vPVAh6/Bz4SB80dISfkxOluGu3ksUfkTErtj0DlBYWBizZ89myZIleHl5ZczR8fHxwdPTsjFi7969KV68OOHh4QAMHDiQ5s2bM2bMGNq3b8/cuXPZtm0bU6dOBSybsA4aNIj33nuP8uXLU7p0aYYNG0ZwcDCdO3e2yecUyauSUtMZuewAs7ZEAdCgTBHG96iNv3cWQ16GATu/ttz5SbsBhQItd31KN7Vy1SIid86mAWjy5MkAtGjRItPrM2bMoE+fPgBERUXh5PTXjapGjRoxe/Zs3nrrLYYOHUr58uVZvHhxponTr732GgkJCTz33HNcu3aNJk2asGLFCjw8spnHIOKATl1O4MVZO9h/Lg6TCfrdV45BrbK565Mcb5nrs3eepV22JTw8FQoVs27RIiJ3SZ5aByivyMk6AiL2aMW+87w6bw/xyWkUKejGJ91r0bxCNmHm/G7LU15XjoPJGe4fBo0GglOeWkVDRCRH3995YhK0iFhHSpqZ8B8PMmPDSQDuKeXLhMdqE+TjeXNnw4CtX8DPb0J6CniXgEemQ8n61i1aRCQXKACJOIgzVxMJm72T3aevAfB88zK80roirs5Z3Mm5cRWW9INDf6y4XrGdZSPTAtkvUioiYk8UgEQcwKqDMQz+bjexN1Lx8XRlzKM1aVUlIOvOp3+H+U9DbBQ4uULrkVD/BTBlMTdIRMROKQCJ5GOp6WY+/jmSz9ceB6BmCR8mPlaHkCIFbu5sNsOmCbBqBJjTwDcUHpkBxetYt2gREStQABLJp6Jjk+g/Zwe/n7wKQJ9GoQxtVxk3lyyGvBIuwaIX4OhKS7vqw9BxHHj4WLFiERHrUQASyYfWHb7IoG93cSUhhULuLox+pAbtqgdl3fnkeljwDMSfBxcPeHAU1O2jIS8RydcUgETykXSzwbhVR5jw6xEMA6oEefPZ43UI9St4c2dzOvw2BtaEg2EGvwrwaAQEVLV63SIi1qYAJJJPXIhPYtDcXWw8dhmAx+qX5O0OVfBwdb65c3w0LHwWTqyztGs9Du0+ArcsgpKISD6kACSSD2w6dpkBc3dyMT6ZAm7OfPBwdTrXLp5156OrYOFzkHgJXAtCh7FQs4d1CxYRsTEFIBE7ZjYbTF57jDE/R2I2oEJAIT57vC7l/Avd3Dk9DVa/D+vHWtoB1SxDXn7lrVqziEheoAAkYqeuJKTw0re7WHv4IgBd65RgZOeqFHDL4j/ra6dhQV84vcXSvqcvtHkfXLNYAVpExAEoAInYoe2nrtJ/9g7OxSbh7uLEyE7V6FYvJOvOh36AxS9C0jVw94aHxlsecxcRcWAKQCJ2xDAMpm84Sfjyg6SZDUr7FeSzx+tQOSiLTf/SkuHnYbD1c0s7uI5lL68ipa1btIhIHqQAJGIn4pJSeX3+Hn7cFw1A++pBjOpaHS8P15s7Xz4G85+y7OQO0LAf3D8cXNysWLGISN6lACRiB/afiyVs1g5OXk7E1dnEW+2r0LthKUxZLVa4Zx4sGwQp18GzCDw8BSq0sXrNIiJ5mQKQSB5mGAZzfz/N8KX7SUkzU7ywJ5Mer0OtkMI3d05JgB9fg53fWNqlGkOXL8Anm8fhRUQcmAKQSB6VmJLGW4v2sXDnWQBaVvJnbLeaFC6QxTBWzAHLkNfFQ4AJmr8GzV4DZ/0nLiKSFf11FMmDjl6I58VZOzgccx0nE7zSpiIvNCuLk9M/hrwMA3bMhB9fh7QkKBQIXb+A0s1sU7iIiJ1QABLJY5bsOsuQhXtJTEnH38ud8T1r06BM0Zs7JsXB9wNh/0JLu+z98PDnUKiYdQsWEbFDCkAieURyWjojlx3gm81RADQqW5RxPWpTzMv95s5nd1iGvK6eBCcXuP9taNgfnJysW7SIiJ1SABLJA6IuJxI2ewd7z8YC0L9lOQa1qoBzVkNemyfDyrfBnAo+JS1r+4TUs0HVIiL2SwFIxMZWHohh8He7iE9Kw7eAK590r0WLiv43d0y8YlnR+fCPlnbljvDQBPD0tW7BIiL5gAKQiI2kpZv56OdIPl97HIA6JQsz8bE6BBfOYn+ukxtgwTMQfw6c3aDNB1DvGchqHSAREflPCkAiNnAhLol+s3ey9eQVAJ5uXJo32lbCzeUfc3jM6bDuY1g7CgwzFC0Hj8yAoBo2qFpEJP9QABKxso3HLjFgzi4uXU+mkLsLHz1Sg7bVg27uGHceFj4LJ3+ztGv2hHYfg3sh6xYsIpIPKQCJWInZbDB57THG/ByJ2YBKgV5MfqIupf0K3tz5yEpY9DwkXgbXgtB+DNTqaf2iRUTyKQUgESu4lpjCS9/uYnXkRQAerVuCEZ2q4enmnLljWgr8OgI2TrC0A6rDozPAr7yVKxYRyd8UgERy2e7T13hx1g7OXruBu4sTIztVo1u9kJs7XjkBC/rC2e2W9r3PwQMjwdXDugWLiDgABSCRXGIYBt9sPsXIZQdJSTdTqmgBPnu8DlWDfW7uvG+hZVXn5DjwKAydJkHlDlavWUTEUSgAieSChOQ03li4l+93nwPgwaqBjH60Bt4erpk7piTCT0Nge4SlHVIfuk6DwlncIRIRkbtGAUjkLjsSE88L32zn2MUEXJxMvNG2En2blMb0zzV7LhyEeU/BxYOACZoOhhZDtYO7iIgV6C+tyF20eKdlI9MbqekEeLsz6bE63BNaJHMnw4AdX/2xg/sNKOgPXaZC2ftsU7SIiANSABK5C5LT0hnx/QFmbbFsZNqknB+f9qiFX6F/bGSaFAfLBsG+BZZ22ZZ/7OCexdYXIiKSaxSARO7Q6SuJvDjLspGpyQT9W5Zn4P3lb97I9O87uJuc4f5h0GigdnAXEbEBBSCRO/DroRhe+nY3sTdSKVzAlU+z2sjUbIbNn8Ev7/xtB/dpEHKvTWoWEREFIJHbkm42+GTlYSauPgpAzZDCfPZ4HYr/cyPThEuw6AU4utLSrvzQHzu4F7ZuwSIikokCkEgOXbqezIA5O9l47DIAvRuW4s32lXF3+ceqzifWwYJn4Xo0OLvDg+Fwz9PawV1EJA9QABLJgW0nrxA2ewcxcckUcHMmvEt1OtUqnrlTehqs/RDWfQQY4FfRsp1FQFWb1CwiIjez6ezLdevW0bFjR4KDgzGZTCxevPg/j5k0aRKVK1fG09OTihUr8tVXX2V6PyIiApPJlOnHw0NbCcidMQyDaetP0GPqZmLikilbrCBLwhrfHH5iz8DMDrBuNGBA7V7w3GqFHxGRPMamd4ASEhKoWbMmTz/9NF26dPnP/pMnT2bIkCF88cUX1KtXj61bt/Lss8/i6+tLx44dM/p5e3sTGRmZ0b5pATqRHIhPSuX1BXtYvjcagA41ghjVtQaF3P/xn8+hH2Dxi5B0Ddy8oOOnUP0Rq9crIiL/zaYBqG3btrRt2/aW+3/99dc8//zzdO/eHYAyZcrw+++/8+GHH2YKQCaTicDAwLterzieQ9FxvPjNDo5fSsDV2cRb7avQu2GpzKE6NQlWDoOtUy3t4DqWp7yKlLFN0SIi8p/sag5QcnLyTcNZnp6ebN26ldTUVFxdLfssXb9+nVKlSmE2m6lTpw4ffPABVatmPwSRnJxMcnJyRjsuLi53PoDYlYU7zjB00V6SUs0E+3gw8fE61Cnpm7nTpSOWtX2i91rajfpDy7fBxc36BYuIyC2zqxXY2rRpw5dffsn27dsxDINt27bx5ZdfkpqayqVLlwCoWLEi06dPZ8mSJXzzzTeYzWYaNWrEmTNnsj1veHg4Pj4+GT8hIdqI0pElpaYzdNFeBn+3m6RUM03L+7FsQNObw8+uOfB5c0v4KVAUHp8Prd9T+BERsQMmwzAMWxcBlmGrRYsW0blz52z73Lhxg7CwML7++msMwyAgIIAnnniC0aNHEx0dTUBAwE3HpKamUrlyZXr27MnIkSOzPG9Wd4BCQkKIjY3F29v7jj+b2I9/ruo8oGV5BvxzVefkePjhFdgz19Iu3QwengreQbYpWkREAMv3t4+Pzy19f9vVHSBPT0+mT59OYmIiJ0+eJCoqitDQULy8vChWrFiWx7i6ulK7dm2OHj2a7Xnd3d3x9vbO9COOZ/WhC3SYsJ69Z2PxLeBKxFP38tIDFTKHn3O7LHd99swFkxO0fAt6LVb4ERGxM3Y1B+hPrq6ulChRAoC5c+fSoUMHnLLZTyk9PZ29e/fSrl07a5YodiTdbPDpL4eZ8Ou/rOpsGLB5Mqx827KdhXcJ6PollGpoo6pFRORO2DQAXb9+PdOdmRMnTrBr1y6KFClCyZIlGTJkCGfPns1Y6+fw4cNs3bqV+vXrc/XqVcaOHcu+ffuYOXNmxjlGjBhBgwYNKFeuHNeuXeOjjz7i1KlTPPPMM1b/fJL3Xb6ezMC5u1h/1DKH7MmGpXizfRXcXP4WqBMuw5IX4fAKS7tSB8t2FgWK2KBiERG5G2wagLZt28Z9992X0R48eDAATz75JBEREZw/f56oqKiM99PT0xkzZgyRkZG4urpy3333sXHjRkJDQzP6XL16lWeffZbo6Gh8fX2pW7cuGzdupEqVKlb7XGIfdkRdJWzWDs7HJuHp6syorlms6nziN1j4LMSft2xn0eZ9qPeMtrMQEbFzeWYSdF6Sk0lUYn8Mw+DrzacYuewAqekGZfwKMqVXXSoEeP3V6abtLCrAIzMgsJrN6hYRkX+Xk+9vu5wDJHK7EpLTGLJwL0t3nwOgXfVAPuxaAy8P1786XTttuesTtcnSrt0L2n4IbgVtULGIiOQGBSBxGEcvXOd/32znyIXruDiZGNKuMk83Ds28qvPB72FJP21nISKSzykAiUP4Yc95Xpu/m4SUdPy93Jn0eB3qhf5tEnNqEvz8Jvz+paVdvC50nQZFStumYBERyVUKQJKvpaabCV9+iOkbTgDQoEwRxvesjb/X37ZUuRgJ85+GmH2WdqMB0HKYVnQWEcnHFIAk34qJSyJs1g62nboKwAvNy/JK6wq4OP/xiLthwM6v4cfXITURChaDh6dAuVY2rFpERKxBAUjypU3HLtN/zg4uXU/By92Fj7vVpE3VwL86JMXC94Ng/0JLu0wLy3YWXjdvpyIiIvmPApDkK4ZhMHXdcUb/FEm62aBSoBdTnqhLqN/fnuA6s80y5HXtFDi5WLazaDQQsllNXERE8h8FIMk34pJSeXXebn7aHwNAl9rFef/h6ni6OVs6mM2w4VNY/T6Y06BwKXhkOpS4x3ZFi4iITSgASb4QGR3PC99s58SlBNycnXi7YxUer1/yr0fc46Nh0fNwfI2lXa0rdPgEPHxsVrOIiNiOApDYvcU7zzJk4V5upKYT7OPBZ0/UpVZI4b86HFkJi16AxEvgWgDajobaT2g7CxERB6YAJHYrJc3Mez8c4KtNpwBoWt6PcT1qU6TgH4+vp6XAqndh00RLO6CaZcirWEUbVSwiInmFApDYpfOxN3hx1g52Rl0DYEDLcgxsVQFnpz/u6lw+ZpnofH6XpX3v8/DACHD1yPJ8IiLiWBSAxO5sOHqJ/nN2ciUhBW8PFz7tUYuWlf72+PruufDDy5ByHTx9odNnUKmd7QoWEZE8RwFI7IbZbDB57THG/ByJ2YCqwd5MeaIuIUUKWDokx8MPr8CeuZZ2qSbQZSr4FLdd0SIikicpAIldiL2Rysvf7eaXg5ZH3LvdU4IRnarh4frHI+7ndlqGvK4cB5MTtBgCTV8GJ2cbVi0iInmVApDkeQfPx/G/b7Zz8nIibi5OjHioKj3uLWl502yGzZPgl3fBnAreJaDrF1CqkW2LFhGRPE0BSPK0RTvPMGThXpJSzRQv7MmUJ+pSvcQfa/dcvwCL/wdHf7G0K3eEjuOhQJHsTygiIoICkORR/3zEvVmFYozrXgvfPx9xP/YrLHweEi6Aiwc8GA51n9LaPiIicksUgCTP+ddH3NNSYPV7sGGcpXOxyvDoDPCvbLuCRUTE7igASZ6y8dgl+s/eyeU/HnH/pHst7q/8xyPuV47Dgmfg7HZL+56+0OZ9cPW0XcEiImKXFIAkTzAMg8/XHWf0ikOYDagc5M2UJ+pQqugfu7jvmQfLXoKUePAoDJ0mWub8iIiI3AYFILG5+KRUXp23hxX7owHoUqc473f+Yxf35Ouw/FXYPdvSuWRD6PIFFA6xYcUiImLvFIDEpo7ExPP819s5fikBV2cTwztW/WsX93M7YX5fuHLMsrZPs9eg2avgrH9tRUTkzuibRGxm2Z5zvDZ/D4kp6QT5eDDp8TrUKelrWdtn0yT45Z0/1vYpbrnrE9rY1iWLiEg+oQAkVpeWbmbUj4f4cv0JABqWKcqEx2rjV8j95rV9KnWAhyZobR8REbmrFIDEqi7GJ9Nv9g62nLgCwPPNy/Bq64q4ODvB0VWw6AWt7SMiIrlOAUisZvupK7w4awcxcckUcnfho0dq0LZ6kGVtn59HwMYJlo7+VeCR6VrbR0REco0CkOQ6wzD4atMp3vvhAKnpBuX8CzHlibqU8y8El4/Bgr6WCc8A9Z6F1iO1to+IiOQqBSDJVTdS0hm6aC+Ldp4FoH2NIEZ3rUFBdxfYPRd+eBlSroOnL3SaBJXa27hiERFxBApAkmtOXU7g+a+3cyg6HmcnE0PaVqJvk9KYkuNh4Suw51tLx1JNoMtU8Clu24JFRMRhKABJrvj1UAyD5u4iLikNv0JuTHysDg3KFIUz2yxDXldPgskZWgyBpoPBydnWJYuIiANRAJK7ymw2GLfqCONWHQGgTsnCfPZ4XQK93OC3sbD6fTCngU9J6PollKxv44pFRMQRKQDJXXMtMYVB3+5iTeRFAHo3LMVb7avglhgDXz8PJ9ZaOlbtAh0+Ac/CtitWREQcmgKQ3BX7z8XywjfbOX3lBu4uToR3qU6XOiUg8kdY/CLcuAKuBaDtaKj9hNb2ERERm1IAkju2cMcZhizcS3KamZAinkx5oi5Vi7lbNjHdOtXSKbCGZW0fv/K2LVZERAQFILkDKWlm3v/hADM3nQKgeYVijOtRi8LXj8MXT8OF/ZaODcKg1XBwcbdhtSIiIn9RAJLbEhOXxIuzdrD91FUABrQsx8D7y+O8Ywb8NBTSkqBgMeg8Bcq3snG1IiIimSkASY5tPXGFsNk7uBifjJeHC590q0WrUFeY1wsOLbN0KtvSEn68AmxbrIiISBYUgOSWGYbBzI0nee+Hg6SZDSoGePF5r7qExu+Ayc9B/DlwcoVW70CDF8HJydYli4iIZMmm31Dr1q2jY8eOBAcHYzKZWLx48X8eM2nSJCpXroynpycVK1bkq6++uqnPvHnzqFSpEh4eHlSvXp3ly5fnQvWO5UZKOi99u4t3vj9AmtmgY81gFr1Qj9DdY2FmR0v4KVoOnvkFGvVT+BERkTzNpt9SCQkJ1KxZk0mTJt1S/8mTJzNkyBDeeecd9u/fz7vvvktYWBjff/99Rp+NGzfSs2dP+vbty86dO+ncuTOdO3dm3759ufUx8r2oy4l0mbyRxbvO4exkYliHKoxvU5gCszrCbx8DBtR6Ap5bC8G1bF2uiIjIfzIZhmHYuggAk8nEokWL6Ny5c7Z9GjVqROPGjfnoo48yXnv55ZfZsmUL69evB6B79+4kJCSwbNmyjD4NGjSgVq1aTJkyJcvzJicnk5ycnNGOi4sjJCSE2NhYvL297/CT2bfVkRcYOGdn5i0tElbDspcgOQ7cfaDjJ1Ctq61LFRERBxcXF4ePj88tfX/b1ThFcnIyHh4emV7z9PRk69atpKamArBp0yZatcr81FGbNm3YtGlTtucNDw/Hx8cn4yckJOTuF29nzGaD8auO8HTE78QlpVG7ZGGWPV+LBrvfsuzllRwHIfXhhd8UfkRExO7YVQBq06YNX375Jdu3b8cwDLZt28aXX35Jamoqly5dAiA6OpqAgMxPHgUEBBAdHZ3teYcMGUJsbGzGz+nTp3P1c+R1cUmpPPf1dsauPIxhwOP1S/JtBzcC57SG3bPB5ATNX4c+y8G3lK3LFRERyTG7egps2LBhREdH06BBAwzDICAggCeffJLRo0fjdAeTbt3d3XF31yJ9AIdj4nn+6+2cuJSAm4sT7z1UhW4piyBipGUTU+8S0GUqhDa2dakiIiK3za7uAHl6ejJ9+nQSExM5efIkUVFRhIaG4uXlRbFixQAIDAwkJiYm03ExMTEEBgbaomS78sOe83SetIETlxII9vFgca8ydDvYH34Zbgk/VTrB/9Yr/IiIiN2zqwD0J1dXV0qUKIGzszNz586lQ4cOGXeAGjZsyKpVqzL1X7lyJQ0bNrRFqXYhLd1M+I8HCZu9g8SUdBqVLcqKdglUWfygZQd31wLw0AR4dCZ4+tq6XBERkTtm0yGw69evc/To0Yz2iRMn2LVrF0WKFKFkyZIMGTKEs2fPZqz1c/jwYbZu3Ur9+vW5evUqY8eOZd++fcycOTPjHAMHDqR58+aMGTOG9u3bM3fuXLZt28bUqVOt/vnswZWEFPrP2cGGo5cBCGsSzMt8g9OiLy0dAmtA12lQrIINqxQREbm7bBqAtm3bxn333ZfRHjx4MABPPvkkERERnD9/nqioqIz309PTGTNmDJGRkbi6unLfffexceNGQkNDM/o0atSI2bNn89ZbbzF06FDKly/P4sWLqVatmtU+l73YeyaWF77ZztlrNyjg5syUBzxotud/cPGgpUPDfnD/29rEVERE8p08sw5QXpKTdQTs1fztZxi6aC8paWZCi3jybZ19BGx6D9KToaA/PDwZymkTUxERsR85+f62q6fA5M6lpJl574cDfLXpFACdyrvxsdtkXNf/bOlQvjV0+gwKFbNhlSIiIrlLAciBXIhL4sVZO9h26ioAY+te4eFTIzFdjwFnd2g9Eu59DkwmG1cqIiKSuxSAHMT2U1f43zc7uBCfjK8HLK64ilL7p1neLFbJMtE5UPOkRETEMSgA5XOGYfDNlihGfL+f1HSDln7XmOw5GffIvZYO9zwNrd8HtwK2LVRERMSKFIDysaTUdN5eso/vtp0BDN4vuZPHrn6G6XoieBaBThOhUntblykiImJ1CkD51LlrN/jfN9vZfSaWwqbrLCjxHWUv/GJ5s3RzePhz8A6ybZEiIiI2ogCUD206dpl+s3dwOSGF+z0P85nnFNwvRoOTi2Vdn4b94Q72ThMREbF3CkD5iGEYTN9wkg+WH8RkTmVU4R/onvQdpkQDipSFrl9C8Tq2LlNERMTmFIDyiRsp6byxcA9Ldp2jpCmGrwp/QWjSAcubtZ+ABz8E90K2LVJERCSPUADKB05fSeT5r7dz4HwsXV02MMo9AtekRPDwgY7joOrDti5RREQkT1EAsnO/HblI/zk7SU+MZYpnBA8a6yEdKNkIukyFwiG2LlFERCTPUQCyU4ZhMHXdcT5ccYjaRPJZgckEmC+AyRlaDIGmg8HJ2dZlioiI5EkKQHYoMSWN1+bv4cc9ZxjgsogBLotxMpvBNxS6fAkh9WxdooiISJ6mAGRnTl1O4Pmvt3M95hjfuX1GXafDljdq9IB2H4FH/ty9XkRE5G5SALIjaw9fZMCcnbRIXs377hEUIhHcvaH9WKjxqK3LExERsRsKQHbAMAwmrz3GlJ928o5LBF3c1lveCGlgmejsW8q2BYqIiNgZBaA8LiE5jVfn7yZ63zp+cJ1EiNNFDJMTpuavQ9NXwFn/E4qIiOSUvj3zsJOXEvjfV1t54PI3jHdbiIvJDIVLYuryJZSsb+vyRERE7JYCUB61JvICo+b8xEjzeOq5/jHRuXo3aP+xZYFDERERuW0KQHmMYRh8tuYYkb/M4DuXaXg73cDsVginDp9AjW62Lk9ERCRfUADKQxKS0xj27UaaHPmQMFfLRGdziXo4df3SssaPiIiI3BUKQHnEyUsJfDLjG16O/4iSzhcx44RTi9dx0kRnERGRu07frHnA2oPn2PftcMYY83FxMpNcqATu3aZByQa2Lk1ERCRfUgCyIcMw+GbFOipveoUwp8NgghuVuuLZ+RNNdBYREclFCkA2kpCcxrzpY+gS/QneTjdIciqI80Of4Fmru61LExERyfcUgGzg9Llojsx4jj6pa8EEFwrXwv/Jr7Sis4iIiJUoAFnZ7g0r8FvZj5ZcJA0nYmoNoHjHYZroLCIiYkX61rWi378bTZ39H+BsMoh2CsT10S8pXrmprcsSERFxOE62LsCR+JRvSDpObPdpg+/Lmymq8CMiImITugNkRRVqN+VUgV+pU6EmJpPJ1uWIiIg4LAUgKytVsZatSxAREXF4GgITERERh6MAJCIiIg5HAUhEREQcjgKQiIiIOBwFIBEREXE4CkAiIiLicBSARERExOHYNACtW7eOjh07EhwcjMlkYvHixf95zKxZs6hZsyYFChQgKCiIp59+msuXL2e8HxERgclkyvTj4eGRi59CRERE7I1NA1BCQgI1a9Zk0qRJt9R/w4YN9O7dm759+7J//37mzZvH1q1befbZZzP18/b25vz58xk/p06dyo3yRURExE7ZdCXotm3b0rZt21vuv2nTJkJDQxkwYAAApUuX5vnnn+fDDz/M1M9kMhEYGHhXaxUREZH8w67mADVs2JDTp0+zfPlyDMMgJiaG+fPn065du0z9rl+/TqlSpQgJCaFTp07s37//X8+bnJxMXFxcph8RERHJv+wqADVu3JhZs2bRvXt33NzcCAwMxMfHJ9MQWsWKFZk+fTpLlizhm2++wWw206hRI86cOZPtecPDw/Hx8cn4CQkJscbHERERERsxGYZh2LoIsAxbLVq0iM6dO2fb58CBA7Rq1YqXXnqJNm3acP78eV599VXq1avHtGnTsjwmNTWVypUr07NnT0aOHJlln+TkZJKTkzPacXFxhISEEBsbi7e39x19LhEREbGOuLg4fHx8bun72652gw8PD6dx48a8+uqrANSoUYOCBQvStGlT3nvvPYKCgm46xtXVldq1a3P06NFsz+vu7o67u3tG+89MqKEwERER+/Hn9/at3NuxqwCUmJiIi0vmkp2dnYHsP2x6ejp79+69aZ7Qv4mPjwfQUJiIiIgdio+Px8fH51/72DQAXb9+PdOdmRMnTrBr1y6KFClCyZIlGTJkCGfPnuWrr74CoGPHjjz77LNMnjw5Ywhs0KBB3HvvvQQHBwMwYsQIGjRoQLly5bh27RofffQRp06d4plnnrnluoKDgzl9+jReXl6YTKa7+pn/HF47ffq0htdyka6zdeg6W4eus3XoOltPbl1rwzCIj4/PyAT/xqYBaNu2bdx3330Z7cGDBwPw5JNPEhERwfnz54mKisp4v0+fPsTHxzNx4kRefvllChcuTMuWLTM9Bn/16lWeffZZoqOj8fX1pW7dumzcuJEqVarccl1OTk6UKFHiLnzC7Hl7e+s/MCvQdbYOXWfr0HW2Dl1n68mNa/1fd37+lGcmQTuKnEzQktun62wdus7WoetsHbrO1pMXrrVdPQYvIiIicjcoAFmZu7s7w4cPz/TUmdx9us7WoetsHbrO1qHrbD154VprCExEREQcju4AiYiIiMNRABIRERGHowAkIiIiDkcBSERERByOAlAumDRpEqGhoXh4eFC/fn22bt36r/3nzZtHpUqV8PDwoHr16ixfvtxKldq3nFznL774gqZNm+Lr64uvry+tWrX6z/9dxCKn/z7/ae7cuZhMpn/d4Fj+ktPrfO3aNcLCwggKCsLd3Z0KFSrob8ctyOl1/vTTT6lYsSKenp6EhITw0ksvkZSUZKVq7dO6devo2LEjwcHBmEwmFi9e/J/HrFmzhjp16uDu7k65cuWIiIjI9Tox5K6aO3eu4ebmZkyfPt3Yv3+/8eyzzxqFCxc2YmJisuy/YcMGw9nZ2Rg9erRx4MAB46233jJcXV2NvXv3Wrly+5LT6/zYY48ZkyZNMnbu3GkcPHjQ6NOnj+Hj42OcOXPGypXbl5xe5z+dOHHCKF68uNG0aVOjU6dO1inWjuX0OicnJxv33HOP0a5dO2P9+vXGiRMnjDVr1hi7du2ycuX2JafXedasWYa7u7sxa9Ys48SJE8ZPP/1kBAUFGS+99JKVK7cvy5cvN958801j4cKFBmAsWrToX/sfP37cKFCggDF48GDjwIEDxoQJEwxnZ2djxYoVuVqnAtBddu+99xphYWEZ7fT0dCM4ONgIDw/Psn+3bt2M9u3bZ3qtfv36xvPPP5+rddq7nF7nf0pLSzO8vLyMmTNn5laJ+cLtXOe0tDSjUaNGxpdffmk8+eSTCkC3IKfXefLkyUaZMmWMlJQUa5WYL+T0OoeFhRktW7bM9NrgwYONxo0b52qd+cmtBKDXXnvNqFq1aqbXunfvbrRp0yYXKzMMDYHdRSkpKWzfvp1WrVplvObk5ESrVq3YtGlTlsds2rQpU3+ANm3aZNtfbu86/1NiYiKpqakUKVIkt8q0e7d7nUeMGIG/vz99+/a1Rpl273au89KlS2nYsCFhYWEEBARQrVo1PvjgA9LT061Vtt25nevcqFEjtm/fnjFMdvz4cZYvX067du2sUrOjsNX3oE03Q81vLl26RHp6OgEBAZleDwgI4NChQ1keEx0dnWX/6OjoXKvT3t3Odf6n119/neDg4Jv+o5O/3M51Xr9+PdOmTWPXrl1WqDB/uJ3rfPz4cX799Vcef/xxli9fztGjR3nxxRdJTU1l+PDh1ijb7tzOdX7ssce4dOkSTZo0wTAM0tLSeOGFFxg6dKg1SnYY2X0PxsXFcePGDTw9PXPl9+oOkDicUaNGMXfuXBYtWoSHh4ety8k34uPj6dWrF1988QV+fn62LidfM5vN+Pv7M3XqVOrWrUv37t158803mTJliq1Ly1fWrFnDBx98wGeffcaOHTtYuHAhP/zwAyNHjrR1aXIX6A7QXeTn54ezszMxMTGZXo+JiSEwMDDLYwIDA3PUX27vOv/p448/ZtSoUfzyyy/UqFEjN8u0ezm9zseOHePkyZN07Ngx4zWz2QyAi4sLkZGRlC1bNneLtkO38+9zUFAQrq6uODs7Z7xWuXJloqOjSUlJwc3NLVdrtke3c52HDRtGr169eOaZZwCoXr06CQkJPPfcc7z55ps4Oekewt2Q3fegt7d3rt39Ad0Buqvc3NyoW7cuq1atynjNbDazatUqGjZsmOUxDRs2zNQfYOXKldn2l9u7zgCjR49m5MiRrFixgnvuuccapdq1nF7nSpUqsXfvXnbt2pXx89BDD3Hfffexa9cuQkJCrFm+3bidf58bN27M0aNHMwImwOHDhwkKClL4ycbtXOfExMSbQs6fodPQNpp3jc2+B3N1irUDmjt3ruHu7m5EREQYBw4cMJ577jmjcOHCRnR0tGEYhtGrVy/jjTfeyOi/YcMGw8XFxfj444+NgwcPGsOHD9dj8Lcgp9d51KhRhpubmzF//nzj/PnzGT/x8fG2+gh2IafX+Z/0FNityel1joqKMry8vIx+/foZkZGRxrJlywx/f3/jvffes9VHsAs5vc7Dhw83vLy8jDlz5hjHjx83fv75Z6Ns2bJGt27dbPUR7EJ8fLyxc+dOY+fOnQZgjB071ti5c6dx6tQpwzAM44033jB69eqV0f/Px+BfffVV4+DBg8akSZP0GLy9mjBhglGyZEnDzc3NuPfee43NmzdnvNe8eXPjySefzNT/u+++MypUqGC4ubkZVatWNX744QcrV2yfcnKdS5UqZQA3/QwfPtz6hduZnP77/HcKQLcup9d548aNRv369Q13d3ejTJkyxvvvv2+kpaVZuWr7k5PrnJqaarzzzjtG2bJlDQ8PDyMkJMR48cUXjatXr1q/cDuyevXqLP/e/nltn3zySaN58+Y3HVOrVi3Dzc3NKFOmjDFjxoxcr9NkGLqPJyIiIo5Fc4BERETE4SgAiYiIiMNRABIRERGHowAkIiIiDkcBSERERByOApCIiIg4HAUgERERcTgKQCIiIuJwFIBERETE4SgAiYiIiMNRABIRuUUtWrRg0KBBd3QOwzB47rnnKFKkCCaTiV27dt2V2kQkZxSARMQuPfXUU7z11lu2LiPHVqxYQUREBMuWLeP8+fNUq1bN1iWJOCQXWxcgIpJT6enpLFu2jB9++MHWpeTYsWPHCAoKolGjRtn2SUlJwc3NzYpViTge3QESkQxz5szB09OT8+fPZ7z21FNPUaNGDWJjY2/7vCVKlOCzzz7L9NrGjRspUKAAp06dyvH5Nm7ciKurK/Xq1cvy/RYtWtC/f38GDRqEr68vAQEBfPHFFyQkJPDUU0/h5eVFuXLl+PHHHzMdl5yczIABA/D398fDw4MmTZrw+++/Z1uH2WwmPDyc0qVL4+npSc2aNZk/f362/fv06UP//v2JiorCZDIRGhqaUW+/fv0YNGgQfn5+tGnThhUrVtCkSRMKFy5M0aJF6dChA8eOHbvp948ePZpy5crh7u5OyZIlef/992/xKoo4NgUgEcnQo0cPKlSowAcffADA8OHD+eWXX/jxxx/x8fG57fPWr18/U5AwDINBgwbx0ksvUapUqRyfb+nSpXTs2BGTyZRtn5kzZ+Ln58fWrVvp378///vf/3j00Udp1KgRO3bsoHXr1vTq1YvExMSMY1577TUWLFjAzJkz2bFjB+XKlaNNmzZcuXIly98RHh7OV199xZQpU9i/fz8vvfQSTzzxBGvXrs2y/7hx4xgxYgQlSpTg/Pnzma7JzJkzcXNzY8OGDUyZMoWEhAQGDx7Mtm3bWLVqFU5OTjz88MOYzeaMY4YMGcKoUaMYNmwYBw4cYPbs2QQEBOT0coo4JkNE5G++//57w93d3XjvvfcMX19fY9++fRnvde7c2ShcuLDRtWvXHJ1z9OjRRtWqVTPaM2fONAIDA434+PjbOm/58uWNZcuWZft+8+bNjSZNmmS009LSjIIFCxq9evXKeO38+fMGYGzatMkwDMO4fv264erqasyaNSujT0pKihEcHGyMHj0647wDBw40DMMwkpKSjAIFChgbN27M9Lv79u1r9OzZM9vaPvnkE6NUqVI31Vu7du1//cwXL140AGPv3r2GYRhGXFyc4e7ubnzxxRf/epyIZE13gEQkkw4dOlClShVGjBjBokWLqFq1asZ7AwcO5KuvvsrxORs0aMDBgwe5fv06CQkJDB06lPfee49ChQrl+LwHDx7k3Llz3H///f/ar0aNGhn/7OzsTNGiRalevXrGa3/eKblw4QJgmZuTmppK48aNM/q4urpy7733cvDgwZvOf/ToURITE3nggQcoVKhQxs9XX31101DVrahbt26m9pEjR+jZsydlypTB29s7Y7gsKioKsFyH5OTk/7wOIpI1TYIWkUxWrFjBoUOHSE9Pv2k4pUWLFqxZsybH56xbty5OTk7s2LGDX375hWLFivHUU0/d1nmXLl3KAw88gIeHx7/2c3V1zdQ2mUyZXvtz+OzvQ0o5cf36dQB++OEHihcvnuk9d3f3HJ+vYMGCmdodO3akVKlSfPHFFwQHB2M2m6lWrRopKSkAeHp63lbdImKhO0AikmHHjh1069aNadOmcf/99zNs2LC7ct4CBQpQvXp1FixYwMcff8wnn3yCk9Pt/flZsmQJnTp1uit1/V3ZsmUz5uD8KTU1ld9//50qVarc1L9KlSq4u7sTFRVFuXLlMv2EhITcUS2XL18mMjKSt956i/vvv5/KlStz9erVTH3Kly+Pp6cnq1atuqPfJeKodAdIRAA4efIk7du3Z+jQoRlDLw0bNmTHjh3UqVPnjs/foEEDJkyYQKdOnWjRosVtnePChQts27aNpUuX3nE9/1SwYEH+97//8eqrr1KkSBFKlizJ6NGjSUxMpG/fvjf19/Ly4pVXXuGll17CbDbTpEkTYmNj2bBhA97e3jz55JO3XYuvry9FixZl6tSpBAUFERUVxRtvvJGpj4eHB6+//jqvvfYabm5uNG7cmIsXL7J///6MeidOnMiiRYsUkkSyoAAkIly5coUHH3yQTp06ZXzR1q9fn7Zt2zJ06FBWrFjxn+eIiIjgqaeewjCMLN+vWbMmrq6ufPTRR7dd5/fff8+9996Ln5/fbZ/j34waNQqz2UyvXr2Ij4/nnnvu4aeffsLX1zfL/iNHjqRYsWKEh4dz/PhxChcuTJ06dRg6dOgd1eHk5MTcuXMZMGAA1apVo2LFiowfP/6m4Dhs2DBcXFx4++23OXfuHEFBQbzwwgsZ71+6dOm25iOJOAKTkd1fKxGRLKxZs4aJEyfetN7N8OHDWbt2bbZzee677z7q1KnDmDFjcnTev3vooYdo0qQJr7322m3XLyICugMkIjnQqlUrdu/eTUJCAiVKlGDevHk0bNgQgB9//JGJEydm6m82m7l48SLTpk3jyJEjLFmyJMfn/bsmTZrQs2fPu//BRMTh6A6QiOSaNWvW0LJlSypVqsSMGTOoX7++rUsSEQEUgERERMQB6TF4ERERcTgKQCIiIuJwFIBERETE4SgAiYiIiMNRABIRERGHowAkIiIiDkcBSERERByOApCIiIg4HAUgERERcTgKQCIiIuJw/g9kY78gPEHUUgAAAABJRU5ErkJggg==", "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": "2025-10-15T23:08:22.244908Z", "iopub.status.busy": "2025-10-15T23:08:22.244769Z", "iopub.status.idle": "2025-10-15T23:08:22.340148Z", "shell.execute_reply": "2025-10-15T23:08:22.339588Z" } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjcAAAG0CAYAAADO5AZFAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjcsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvTLEjVAAAAAlwSFlzAAAPYQAAD2EBqD+naQAAZ1ZJREFUeJzt3Xd4FOUCxeHfbjqk0HvovYUOgaAoKAJSpUgXQSyAoNeGXSzYC4IgioAioPQOAgpI7xB6J/RQ0yB15/4xGEUgkpBkspvzPs8+ZiezuydzgZw788332QzDMBARERFxEXarA4iIiIikJ5UbERERcSkqNyIiIuJSVG5ERETEpajciIiIiEtRuRERERGXonIjIiIiLsXd6gCZzeFwcPr0afz8/LDZbFbHERERkTtgGAZRUVEUKVIEuz3lczPZrtycPn2awMBAq2OIiIhIGpw4cYJixYqluE+2Kzd+fn6AeXD8/f0tTiMiIiJ3IjIyksDAwOTf4ynJduXmr0tR/v7+KjciIiJO5k6GlGhAsYiIiLgUlRsRERFxKSo3IiIi4lJUbkRERMSlqNyIiIiIS1G5EREREZeiciMiIiIuReVGREREXIrKjYiIiLgUlRsRERFxKSo3IiIi4lJUbkRERMSlqNxkkGsxUVZHEBERyZay3argGe7SESJ/6on/5V1sdZRlX6tZNC6Xj8A8OaxOJiIiki2o3KSng8tgajf8k+IAqGU/RIdZoQBUKxpAt/rFaV+zKN4eblamFBERcWm6LJVekhJg3mC4Xmz+Uq9UHuw2CD0VwdCZoYR89Acjfz9IxLUEi4KKiIi4NpWb9BK+ByJPgpf/DZt/fTKYja8147WWlSgS4M2F6Dg+/e0ATT75g/FrjhKf6LAosIiIiGtSuUkvl4+b/81X/qZv5fP14ol7SrPypfv4sksNyhbw5fLVBN6Zt4cHv1jJ0j3nMjmsiIiI61K5SS9evuZ/E67edhcPNzvtahZl8eDGfNC+Gvl8vTh28SpP/LiZ/j9u5vSVa5kUVkRExHWp3KSXnPnN/0ae/s9d3d3sdKtfnBUvNuGpe8vgbrfx255zPPD5Sr7/8wiJSbpUJSIiklYqN+klX3lw94HYK3f8El8vd15pUZEFzzamToncxMQn8d6CvXQYvZaD5zRPjoiISFqo3KQXdy8oXj9NL61QyI9fnwzmww7V8Pd2Z+fJCFp9vZpvVx4myWGkc1ARERHXpnKTnko2TvNL7XYbj9YrztLn7+W+CvmJT3QwfNE+On+7jqMXYtIxpIiIiGtTuUlPldokf7mzwmAOd1iU6rco6O/ND4/V5aNHquHr5c6W45dp8dUqflx3DMPQWRwREZH/onKTnvKXh2J1AaheogBlqjdM09vYbDa61C3O4iGNCS6dl9gEB2/O2U2/iZu5GB33328gIiKSjancpLca3cz/bv8Z7vJMS7HcOfi5X33eal0ZT3c7y/eF0+KrP1lz6EI6BBUREXFNKjfprUoHcPc2Zyw+teWu385ut9GnUSlmP9OIMvlzEh4VR49xG/hw0T4SdMu4iIjITVRu0ptPLqj6iPn1mi/T7W0rF/Fn/qDGdK1XHMOAMSsP03HMOo5f1GBjERGRf1K5yQgNnzX/u3c+nD+Qbm/r4+nG8A7VGN29Fv7e7uw4cYVWI1YzZ/updPsMERERZ6dykxEKVIQKrQAD1n6V7m/folphFg25h3ol8xAdl8jgqdt5fXYosQlJ6f5ZIiIizkblJqOEPGf+d8cvEJH+Z1aK5vJh8hP1GXR/WWw2mLQ+jI5j1hJ28fZrW4mIiGQHKjcZJbAulAgBRwKsG5UhH+HuZud/D1Zg/GN1yZ3Dg12nImn19Z8s2X02Qz5PRETEGajcZKS/zt5sHpchZ2/+0qRCARY825haxXMRFZvIkz9t4f0Fe3Q3lYiIZEsqNxmpbFMoHgyJsbBieIZ+VJFcPvzyZDD9QkoB8N2fR3l07HrORFzL0M8VERHJalRuMpLNBg8MM7/e/jOE78vQj/Nws/P6w5UZ06M2ft7m0g2tRqxm1YHzGfq5IiIiWYnKTUYLrAeVWoPhgGVvZ8pHPlS1EPMHhVCliD+XYuLpPX4jo/44pLWpREQkW1C5yQxN3wKbGxxYBMfXZspHlsibkxlPN6RrvUAMAz5Zsp+nJm0hKjYhUz5fRETEKio3mSFfOajVy/x66Zt3vebUnfL2cGN4h+oM71ANTzc7S3afo92oNRwKj86UzxcREbGCyk1mafIKeOSAk5tg79xM/eiu9Yrzy5MNKOTvzeHzMbQbtUa3i4uIiMtSucksfoUgeKD59W+vQ0Lm3sVUs3hu5g0KoX4pc1bjJ3/awidL9pHk0DgcERFxLSo3manRYPAvClfC4M/PM/3j8/t5MalffR5vZN4uPuqPwzw+YRNXrsZnehYREZGMonKTmbx84aHr892s+RIuHs70CB5udt5sXZmvHq2Bt4edlQfO03rkavacjsz0LCIiIhlB5SazVWoDZe6HpHhY9FKmDS7+t7Y1ijLz6UYE5vHhxKVrdBi9RquLi4iIS7C03Lz99tvYbLYbHhUrVkzxNdOmTaNixYp4e3tTrVo1Fi5cmElp04nNBi0/BTdPOLQM9s6zLErlIv7MGxjCPeXzE5vgYPDU7Qybp2UbRETEuVl+5qZKlSqcOXMm+bF69erb7rt27Vq6du1K37592bZtG+3ataNdu3bs2rUrExOng7xloOGz5teLh0J8jGVRcuXwZPxjdRl4X1kAflhzlB7fb+B8VJxlmURERO6G5eXG3d2dQoUKJT/y5ct3232/+uorHnroIV588UUqVarEu+++S61atRg5cmQmJk4njf8HAcUh8iSs+sTSKG52Gy80r8C3PWvj6+XOhqOXaP31arafuGJpLhERkbSwvNwcPHiQIkWKULp0abp3705YWNht9123bh3NmjW7YVvz5s1Zt27dbV8TFxdHZGTkDY8swTMHtPjQ/HrtSDh/wNo8QPMqhZg9oBFl8ufkbGQsncesY+rG2//vISIikhVZWm7q16/PhAkTWLx4MaNHj+bo0aM0btyYqKioW+5/9uxZChYseMO2ggULcvbs7SekGz58OAEBAcmPwMDAdP0Z7kqFllCuOTgSYOH/LBtc/E9lC/gyZ2AID1UpRHySg1dmhjJ0ZijxiRqHIyIizsHSctOiRQs6depE9erVad68OQsXLuTKlSv8+uuv6fYZQ4cOJSIiIvlx4sSJdHvvu2azQYuPwN0bjq6C7ZOtTgSAr5c7o3vU4sXmFbDZYMrGMLp/v54L0RqHIyIiWZ/ll6X+KVeuXJQvX55Dhw7d8vuFChXi3LlzN2w7d+4chQoVuu17enl54e/vf8MjS8lTylyaAWDJqxAdbm2e62w2GwPuK8sPvevi5+XOpmOXaTtyDbtPR1gdTUREJEVZqtxER0dz+PBhChcufMvvBwcHs3z58hu2LV26lODg4MyIl3GCB0Gh6hB7BRa+aHWaG9xXsQCzBjSidL6cnLpyjUdGr2XBzjNWxxIREbktS8vNCy+8wMqVKzl27Bhr166lffv2uLm50bVrVwB69erF0KFDk/cfPHgwixcv5rPPPmPfvn28/fbbbN68mYEDB1r1I6QPN3doOxJsbrBnNuydb3WiG5Qt4MusAY2S58MZMHkrn/22H4fWpRIRkSzI0nJz8uRJunbtSoUKFejcuTN58+Zl/fr15M+fH4CwsDDOnPn7LEHDhg2ZPHkyY8eOJSgoiOnTpzN79myqVq1q1Y+QfgoHQaPrc98s+B9cu2JpnH8L8PFg/GN16X9PaQC+/v0QT07aQnRcosXJREREbmQzjCxwi04mioyMJCAggIiIiKw3/ibhGoxuBJcOQ40e0G6U1YluaebWk7xy/Q6q8gV9+b5XXYrnzWF1LBERcWGp+f2dpcbcZHsePtDuG8AG2yfBgSVWJ7qlDrWK8euTwRTw8+LAuWjajFrN2kMXrI4lIiICqNxkPcUbQPAA8+u5z8K1y9bmuY0agbmYNyiEoMBcXLmaQM8fNjJx7TGy2YlAERHJglRusqL7X4e8ZSH6LCx6xeo0t1XQ35tf+jegQ82iJDkM3pq7WxP+iYiI5VRusiIPH2g3Bmx22DkV9mXdlc+9Pdz4rHMQr7WshN0GUzedoNt367XwpoiIWEblJqsKrAsNB5lfz3sWYrLumBabzcYT95Tmh8fq4uftzubjl2k7cjW7TmnCPxERyXwqN1lZk1chfyWIOQ/zBmeJtadS0qRCAWYPaETp/Dk5HRFLxzFrmbfjtNWxREQkm1G5yco8vOGR78DuAfvmw7ZJVif6T2Xy+zLrmUbce33Cv0FTtvHpEk34JyIimUflJqsrVM0cYAyw+BW4dNTaPHcgwMeDHx6ry5PXJ/wb+cch+v+0hajYBIuTiYhIdqBy4wwaDoLiDSE+GmY9CY4kqxP9Jze7jaEtK/FFlyA83e0s23uODt+s5fjFGKujiYiIi1O5cQZ2N2g/Bjz94MQGWP2F1YnuWPuaxZj2ZDAF/b04GB5Nm5FrWH0w6w6OFhER56dy4yxyl4CWn5hfrxgOp7ZYmycVggJzMW9gCDUCcxFxLYHe4zcyfs1RTfgnIiIZQuXGmQQ9CpXbgiMRpveFuCirE92xAv7eTO3fgEdqFSPJYfDOvD28PGMncYlZ/xKbiIg4F5UbZ2KzQeuvICAQLh+FhS9anShVvD3c+LRTdV5vZU749+vmk3T7bgPhUbFWRxMREReicuNsfHJDh+/M2Yt3TIGdv1qdKFVsNhv9GpdmfJ96+Hu7s+X4ZdqOXEPoSU34JyIi6UPlxhmVCIZ7Xza/nv+8U9we/m/3ls+fPOHfmesT/s3VhH8iIpIOVG6cVeMXoHgwxEfBjL6Q5HxzyJTO78vsAY24r0J+4hIdPDtlGx8v3qcJ/0RE5K6o3DgrN3fz8pR3gHnn1O/vWZ0oTfy9Pfi+d12eurcMAN+sOMwTP27WhH8iIpJmKjfOLFcgtB5hfr3mSzi4zNI4aeVmt/FKi4p89WgNvNztLN8XTvtv1nL0gib8ExGR1FO5cXZV2kHdfubXs/pDxClL49yNtjWKMu2pYAr5e3MoPJq2I1fz58HzVscSEREno3LjCh58HwpVh6sXr4+/SbQ6UZpVL5aLuQMbUbN4LiJjE+n9w0bGrdaEfyIicudUblyBhzd0mmAuzxC2Dv543+pEd+WvCf861i6Gw4B35+/hxema8E9ERO6Myo2ryFsG2lwff7P6c6cdf/MXL3c3PulYnTcerozdBtO3nKTr2PWa8E9ERP6Tyo0rqdrBZcbfgDnhX9+QUky4PuHf1rArtPl6DTtPXrE6moiIZGEqN67mpvE3zn9L9T3l8zNnYAhlC/hyNjKWTmPWMWe7cxc3ERHJOCo3rubf42+WvW11onRRKl9OZj3TkPsrFiAu0cHgqdv5cNE+kjThn4iI/IvKjSvKWwbafWN+vW4k7J5taZz04uftwXe96vBME3PCvzErzQn/IjXhn4iI/IPKjauq3AYaPmt+PWcAnD9gbZ504ma38dJDf0/49/u+cNqPWsOR89FWRxMRkSxC5caVNX0LSjaG+Gj4pQfEuU4BaFujKNOfakjhAG8On4+h3ag1rDqgCf9ERETlxrW5uUPHH8CvMFzYD3MHgQtNhletWABzBjaidoncRMYm8tj4jXz/5xFN+Cciks2p3Lg63wLQaSLY3WH3TNgwxupE6aqAnzeTn6hP5zrmhH/vLdjLC9N2EpugCf9ERLIrlZvsoHh98xZxgN9eh+PrrM2Tzrzc3fjokeq81boybnYbM7ae5NGx6wmP1IR/IiLZkcpNdlH/SajaERyJMK03RJ62OlG6stls9GlUiol96hHg48H2E1doPXI1O05csTqaiIhkMpWb7MJmg9ZfQYHKEH3OHGCc4HpnNkLK5WPOgEaUK+DLucg4On27jtnbNOGfiEh2onKTnXj5wqM/g3cuOLUFFjzvUgOM/1IyX05mPtOQZpUKEJ/oYMgv2xm+cK8m/BMRySZUbrKbPKXNGYxtdtj+M2wca3WiDOHn7cHYnnUYcJ854d+3q47Qd+ImIq5pwj8REVencpMdlbkPHnjX/HrxUDj6p7V5MojdbuPF5hX5umtNvD3srNh/nvbfaMI/ERFXp3KTXQUPgOpdwEgyBxhfCbM6UYZpHVSE6U81pEiAN0fOx9B21BpW7A+3OpaIiGQQlZvs6q8BxoWDzBXEp3aD+KtWp8owVYsGMGdgCHVK5CYqNpHHJ2ziu1Wa8E9ExBWp3GRnHj7Q5WfIkQ/OhsLcgS45wPgv+f28+PmJ+nSpE4jDgPcX7uV/v+7QhH8iIi5G5Sa7yxUInX80ZzDeNQPWfGl1ogzl5e7Gh49U4+3rE/7N3HaKLmPXc04T/omIuAyVG4GSjaDFR+bXy96BfQuszZPBbDYbjzUqxY+P1yNXDg92nLhC669XszXsstXRREQkHajciKluP/OBATOegDM7rU6U4RqVNSf8K1/Ql/CoOB79dj3TNp+wOpaIiNwllRv520MfQukmkBADU7pC1DmrE2W4EnlzMvOZRjSvUpD4JAcvTt/JO/N2k5jksDqaiIikkcqN/M3Nw5zgL29ZiDxp3kHlgks0/Juvlzuju9dmcNNyAIxfc4ze4zdyOSbe4mQiIpIWKjdyI5/c0O3X60s0bIY5A1z6Dqq/2O02nnugPGN61CKHpxtrDl2k7ag17D8bZXU0ERFJJZUbuVneMtDlp+t3UE2HVZ9anSjTPFS1MDOfaUhgHh/CLl2lwzdrWLL7rNWxREQkFbJMufnwww+x2WwMGTLktvtMmDABm812w8Pb2zvzQmYnpe6Blp+YX//xHuyebWmczFSxkD9zB4TQsExeYuKTePKnLXy17CAOLbwpIuIUskS52bRpE99++y3Vq1f/z339/f05c+ZM8uP48eOZkDCbqvM41H/a/HrWU3Bqq7V5MlHunJ5MfLwejzUsCcAXyw7wzM9biYlLtDaYiIj8J8vLTXR0NN27d+e7774jd+7c/7m/zWajUKFCyY+CBQtmQspsrPn7UPYBSLwGUx516TWo/s3Dzc7bbarw8SPV8XCzsXj3WR4ZvZYTl1x3mQoREVdgebkZMGAArVq1olmzZne0f3R0NCVKlCAwMJC2bduye/fuFPePi4sjMjLyhoekgt0NOv4ABatC9Dn4uRNcu2J1qkzVuW4gU/s3IJ+vF/vORtFm5GrWHrpgdSwREbkNS8vN1KlT2bp1K8OHD7+j/StUqMAPP/zAnDlzmDRpEg6Hg4YNG3Ly5Mnbvmb48OEEBAQkPwIDA9Mrfvbh7W/eQeVXBM7vg196QGL2uk26dok8zBvUiOrFArh8NYGeP2xkwpqjWnhTRCQLshkW/et84sQJ6tSpw9KlS5PH2jRp0oQaNWrw5Zdf3tF7JCQkUKlSJbp27cq77757y33i4uKIi4tLfh4ZGUlgYCARERH4+/vf9c+RrZwNhR9aQHwUVH8U2o8xVxfPRmITkhg6M5RZ204B0KVOIMPaVcHL3c3iZCIiri0yMpKAgIA7+v1tWbmZPXs27du3x83t718KSUlJ2Gw27HY7cXFxN3zvdjp16oS7uztTpky5o89NzcGRWzi0DH7uDEYS3Psy3Peq1YkynWEYfP/nUYYv2ovDgFrFczGmZ20K+OnOPRGRjJKa39+WXZZq2rQpoaGhbN++PflRp04dunfvzvbt2++o2CQlJREaGkrhwoUzIbEAULYZPPyF+fXKj2DbJGvzWMBms/HEPaUZ36ce/t7ubA27Qpuv17Dz5BWro4mICBaWGz8/P6pWrXrDI2fOnOTNm5eqVasC0KtXL4YOHZr8mmHDhvHbb79x5MgRtm7dSo8ePTh+/Dj9+vWz6sfInmr3hsb/M7+eNxgO/2FtHovcWz4/cwaGULaAL2cjY+k0Zh2ztt1+/JeIiGQOy++WSklYWBhnzpxJfn758mWeeOIJKlWqRMuWLYmMjGTt2rVUrlzZwpTZ1P1vQLVO4EiEX3vBuZTvWnNVpfLlZNYzDWlasQBxiQ6e+2UHHyzcS5Im/BMRsYxlY26sojE36SgxDn7qAMdXg39R6LcM/ItYncoSDofB50sPMPKPQwDcUz4/Xz9ak4AcHhYnExFxDU4x5kZcgLsXPDoJ8pWHyFPmQOPY7DmPkN1u44XmFRjZrSbeHnZWHThPu2/WcChcC2+KiGQ2lRu5Oz65ofs0yJkfzoXCL93NMzrZ1MPVizDj6YYUzeXD0QsxtBu1luV7z1kdS0QkW1G5kbuXu6RZcDx94egqmNkfHElWp7JMlSIBzB3YiHql8hAdl0i/Hzcz6o9DmvBPRCSTqNxI+ihSE7pMArsH7JkNi16GbPzLPK+vF5P61qdHg+IYBnyyZD+DpmzjarwW3hQRyWgqN5J+ytwHHcYCNtj0Haz61OpElvJ0t/Neu2q8374q7nYb83eeoePodZy8rIU3RUQyksqNpK+qHaDFR+bXf7wHm8dbmycL6F6/BJOfaEDenJ7sORNJ25Fr2HDkotWxRERclsqNpL/6T0LjF8yvFzwPe+dZmycLqFcqD3MHhVCliD8XY+Lp/v0GJq0/bnUsERGXpHIjGeP+16FWLzAcML0vHFttdSLLFc3lw/SnGtI6qAiJDoPXZ+/i1VmhxCc6rI4mIuJSVG4kY9hs0OoLqNAKkuJgSldzVfFszsfTjRGP1uClhypgs8HkDWH0+H4DF6Kz7+3zIiLpTeVGMo6bO3QcB8WDIS4SJj0Cl49ZncpyNpuNZ5qUZVzvOvh5ubPx2CXafL2aXacirI4mIuISVG4kY3n4QNcpUKAyRJ8zl2uIDrc6VZZwf8WCzBrQiNL5cnI6IpaOY9YyZ/spq2OJiDg9lRvJeD65occMCCgOlw7DT+3h2mWrU2UJZQv4MmtAI5pUyE9sgoPBU7fzwcK9JCZpHI6ISFqp3Ejm8C8CvWZDzgJwbhdM6ghxWncJIMDHg3G96/J0kzIAjF11hD4TNnHlarzFyUREnJPKjWSevGXMguOdC05tNgcZJ8RanSpLcLPbePmhiozsVhMfDzf+PHiB1iNXs/dM9lyIVETkbqjcSOYqWAV6zDTXoTr2J0zrDUkJVqfKMh6uXoSZzzQkMI8PJy5do8M3a5m/87TVsUREnIrKjWS+YrWh61Rw94YDi2HWk9l6oc1/q1TYn7kDQmhcLh/XEpIYOHkbHy7aR5Ij+67VJSKSGio3Yo1SjaHzT2B3h10zYP6QbL3Q5r/lzunJ+Mfq0v+e0gCMWXmYPhM2EXFVZ7lERP6Lyo1Yp/yD8Mj3YLPD1h/ht9dVcP7B3c3Oqy0r8dWjNfD2sLPqwHnajFrNgXMaiC0ikhKVG7FWlfbQeoT59bqRsPJja/NkQW1rFGXG0w0pmsuH4xev0m7UGhbvOmN1LBGRLEvlRqxXqyc89KH59YoPYN0oa/NkQVWKBDBvUAgNy+TlanwST03ayqdL9uPQOBwRkZuo3EjW0OBpuO918+slr8LmH6zNkwXlyenJj4/Xo29IKQBG/nGIfj9uJuKaxuGIiPyTyo1kHfe8AI0Gm1/Pfw62TLQ2Txbk7mbnjYcr80WXILzc7fy+L5x2o9ZwKFzjcERE/qJyI1mHzQbN3oEGz5jP5w2GbZOszZRFta9ZjBlPN6RIgDdHL8TQbtRaftt91upYIiJZgsqNZC02GzT/AOo9CRgwZyBsn2J1qiypatEA5g4KoX6pPETHJdL/py18sfSAxuGISLanciNZj80GLT6Cuv0AA2Y/DTt/tTpVlpTP14tJ/erzWMOSAHy1/CD9f9pCVKzG4YhI9qVyI1mTzQYtPoHafQDDnMU4dLrVqbIkDzc7b7epwqedgvB0t7Ns7znajVrD4fPRVkcTEbGEyo1kXXY7tPocavUCwwEz+8PuWVanyrI61i7GtCeDKRzgzeHzMbQbuYZle85ZHUtEJNOp3EjWZrfDw19BjR5gJMH0vrBnjtWpsqygwFzMHRhCvZJ5iIpLpN+Pmxmx/KDG4YhItqJyI1mf3Q5tRkBQ1+sF53HYO9/qVFlWfj9zHE6v4BIAfL70AE9N2kJ0XKLFyUREMofKjTgHuxu0HQXVOoMjEaY9BvsXWZ0qy/J0tzOsbVU+eqQanm52fttzjvaj1nD0QozV0UREMpzKjTgPuxu0Gw1VHwFHAvzSEw4ssTpVltalbnGmPtmAgv5eHAyPps3I1fyxL9zqWCIiGUrlRpyLmzu0HwuV210vOD1g/2KrU2VptYrnZt7AEGqXyE1UbCKPT9zE1xqHIyIuTOVGnI+bOzzyPVRqA0nxZsHZM9fqVFlaAX9vpjzRgG71i2MY8NnSAzw5SfPhiIhrUrkR5+TmAR1/+PsS1bTHNA/Of/B0t/NB+2p82MEch7N0zznaal0qEXFBKjfivNw8oMN3ENTNvItqRj/Y9rPVqbK8R+sV59enzPlwjpyPoe3INSzedcbqWCIi6UblRpzbX3dR1X4Mcy2qZ2DzeKtTZXk1AnMxb1AIDUrnISY+iacmbeWjxftI0jgcEXEBKjfi/Ox2ePhLqP+U+Xz+EFg/xspETiGfrxeT+tanX0gpAEavOMxj4zdyOSbe4mQiIndH5UZcg80GD30IjQabzxe/DKu/tDSSM3B3s/P6w5X56tEaeHvY+fPgBVqPXM2uUxFWRxMRSTOVG3EdNhs0ewfufdl8vuwtWPERGLrU8l/a1ijKrGcaUTxPDk5evsYjo9cya9tJq2OJiKSJyo24FpsN7nsVmr5pPl/xASwfpoJzByoV9mfewBCaVMhPXKKD537Zwdtzd5OQ5LA6mohIqqjciGtq/D9o/oH59erPYclrKjh3ICCHB+N61+XZ+8sCMGHtMbp/t4HwqFiLk4mI3DmVG3FdwQOg1Wfm1+tHwYL/gUNnIf6Lm93G8w9W4LtedfDzcmfjsUu0/no1W45ftjqaiMgdUbkR11a3H7QZCdhg8ziYMwCStDr2nXigckFmD2xE2QK+nIuM49Gx65i0/jiGzoCJSBanciOur1ZP6DAWbG6wYzL82gsSdJnlTpTJ78vsAY1oUbUQCUkGr8/excszdhKbkGR1NBGR21K5keyhemfoMgncvGD/Avi5I8RGWp3KKfh6ufNN91q8/FBF7Db4dfNJOn+7jlNXrlkdTUTkllRuJPuo2BJ6zgQvfzj2J0xsDTEXrE7lFGw2G083KcPEx+uRK4cHO09G0Prr1aw9rOMnIllPlik3H374ITabjSFDhqS437Rp06hYsSLe3t5Uq1aNhQsXZk5AcQ0lQ6D3PMiRD85shx+aw5UTVqdyGo3L5WfewBCqFPHnUkw8Pb7fwHerjmgcjohkKVmi3GzatIlvv/2W6tWrp7jf2rVr6dq1K3379mXbtm20a9eOdu3asWvXrkxKKi6hSA14fAkEBMLFQ2bBOX/A6lROIzBPDmY83ZAOtYriMOD9hXsZNGUbV+M1UFtEsgabkcb/y7Vnzx7CwsKIj79xHZo2bdqk6n2io6OpVasW33zzDe+99x41atTgyy+/vOW+Xbp0ISYmhvnz5ydva9CgATVq1GDMmDtbSygyMpKAgAAiIiLw9/dPVVZxMRGn4Kf2cGE/+OSBHjOgaC2rUzkNwzD4af1xhs3bQ6LDoEJBP77tWZuS+XJaHU1EXFBqfn+n+szNkSNHCAoKomrVqrRq1Sr57En79u1p3759qsMOGDCAVq1a0axZs//cd926dTft17x5c9atW3fb18TFxREZGXnDQwSAgKLQZxEUqQXXLpljcI6stDqV07DZbPQKLsmU/g3I7+fF/nNRtB65mt/3nbM6mohkc6kuN4MHD6ZUqVKEh4eTI0cOdu/ezapVq6hTpw4rVqxI1XtNnTqVrVu3Mnz48Dva/+zZsxQsWPCGbQULFuTs2bO3fc3w4cMJCAhIfgQGBqYqo7i4nHmh91wodS/ER5t3Ue2dZ3Uqp1K3ZB7mDwqhdoncRMUm0nfiZr5adhCHQ+NwRMQaqS4369atY9iwYeTLlw+73Y7dbickJIThw4fz7LPP3vH7nDhxgsGDB/Pzzz/j7e2d2hh3bOjQoURERCQ/TpzQ4FH5Fy8/6D4NKrWGpHhzHpytP1mdyqkU9PdmyhMN6NmgBIYBXyw7QP+fNhMZm2B1NBHJhlJdbpKSkvDz8wMgX758nD59GoASJUqwf//+O36fLVu2EB4eTq1atXB3d8fd3Z2VK1cyYsQI3N3dSUq6eZKwQoUKce7cjae8z507R6FChW77OV5eXvj7+9/wELmJuxd0nAA1e4LhgLkDYe3XVqdyKp7udt5tV5VPOlbH093Osr3htB25hgPnoqyOJiLZTKrLTdWqVdmxYwcA9evX5+OPP2bNmjUMGzaM0qVL3/H7NG3alNDQULZv3578qFOnDt27d2f79u24ubnd9Jrg4GCWL19+w7alS5cSHByc2h9D5GZu7tDma2h4/Qzkb6+bC25qPapU6VQnkBlPNaRoLh+OXoih3ag1zNtx2upYIpKNuKf2Ba+//joxMTEADBs2jIcffpjGjRuTN29efvnllzt+Hz8/P6pWrXrDtpw5c5I3b97k7b169aJo0aLJY3IGDx7Mvffey2effUarVq2YOnUqmzdvZuzYsan9MURuzWaDB9+FHHlh2VuwbiREnob2Y8yzO3JHqhULYO7ARjw7dRtrDl1k0JRtbA27zKstK+HhliVmoBARF5bqf2WaNGlC8+bNAShbtiz79u3jwoULhIeHc//996druLCwMM6cOZP8vGHDhkyePJmxY8cSFBTE9OnTmT179k0lSeSuhQyB9mPB7gG7Z8JPHeCaVsVOjby+XkzsU4+nm5QBYPyaY3Qdu55zkVrXS0Qy1h3Pc3P+/Hl69erFsmXLcDgc1K1bl0mTJlG2bNmMzpiuNM+NpMqRFfBLT4iLhPwVoft0yKU77lLrt91n+d+vO4iKSySfrxejutWkfum8VscSESeSIfPcvPzyy2zfvp1hw4bx6aefcuXKFZ544om7DiuSpZVuYs6F41cEzu+D75vB2VCrUzmdB6sUYu6gECoW8uNCdBzdtGyDiGSgOz5zExgYyPfff598SergwYNUqlSJmJgYvLycZyyCztxImkSchEkd4fxe8PSDLj9BmfusTuV0rsYn8tqsXczadgqAltUK8XHHIHy9Uj38T0SymQw5c3P69GmCgoKSn5crVw4vL68bxsSIuKyAYvD4YijZGOKjzMn+dky1OpXTyeHpzuedgxjWtgoebjYWhp6l7cjVHArX7eIikn5SNaD437dnu7m56bSyZB8+ucz1p6o+Ao5EmPUk/PkZ6O9Aqvy1bMPU/sEU8vfm8PkY2oxcw/ydul1cRNLHHV+WstvtBAQEYLPZkrdduXIFf39/7Pa/O9KlS5fSP2U60mUpuWsOByx/G9Z8ZT6v8zi0+MScJ0dS5UJ0HIMmb2PdkYsA9A0pxSstKup2cRG5SWp+f9/xv8bjx4+/62AiLsFuhweGgX8xWPQSbP4BIs9Ax3HgqRWxUyOfrxc/9a3Hp78dYMzKw4xbfZTQkxGM7FaTAv4ZtyyLiLi2Oz5z4yp05kbS1d55MKMfJMZC0drQ9RfwzW91Kqe0eNdZXpi2g+i4RPL7efFN91rULZnH6lgikkVkyIBiEbmFSq2h11zwyQ2ntsD398O5PVanckoPVS3E3IGNKF/Ql/NRcXQdu55xq49qXJ+IpNodn7m503Wjjhw5cleBMprO3EiGuHAQJneGS0fMW8U7TYByzaxO5ZSuxifyyoxQ5l5fj+rh6oX56JHq5NTt4iLZWmp+f6dqQHGJEiXo1q0bBQoUuO1+gwcPTl3aTKZyIxnm6iVzNuPjq8Fmh4c+hPpPWp3KKRmGwcS1x3hvwV4SHQZlC/gypkdtyhbwtTqaiFgkQ8rNtGnT+OGHH1ixYgUtWrTg8ccfp2XLljfcKeUMVG4kQyXGw/znYPsk83ndJ8ySozup0mTL8Us88/NWzkXGkdPTjU87BdGiWmGrY4mIBTKk3Pzl1KlTTJgwgQkTJnD16lV69uxJ3759KVeu3F2FziwqN5LhDMO8TXzZ24ABZe43L1N5B1gczDmdj4pj4OStbDhqTjPR/57SvNS8Au66XVwkW8nQcvNPK1eu5O2332bVqlVcuHCB3Llzp/WtMo3KjWSavfNh5hOQcBXyVYBuv0CeUlanckqJSQ4+WbKfb1eZY/rql8rD191qUsBPt4uLZBcZfrdUbGwskyZN4p133mHDhg106tSJHDlypCmsiMuq9PD1RTcLw4X98H1TOL7O6lROyd3NztCWlRjdvRa+Xu5sOHqJh0esZvOxrD1pqIhYI1XlZsOGDfTv359ChQrx+eef06FDB06dOsXUqVOdavFMkUxTpAY88TsUDoKrF+HHNrB9itWpnFaLaoWZM7AR5Qr4Eh4Vx6Nj1zN+jW4XF5Eb3fFlqSpVqhAeHk63bt14/PHHb1hE05nospRYIj7GXItq7zzzeeMX4L7XzNmOJdVi4hJ5ecZO5u80F+5tE1SE4R2q6XZxEReWYbeC58yZE3d39xvWl/o3rS0lchsOB/z+Lqz+3HxeuS20GwOeuqSbFoZhMH7NMT5YaN4uXr6gebt46fy6XVzEFWVIuZk4ceIdfXjv3r3vaD+rqNyI5bZPhrnPgiMBCteAR3+GgGJWp3Jam45dYsDPWwmPisPXy51POlbX7eIiLijT7pZyRio3kiUcWwO/9IBrlyBnfuj8E5QItjqV0wqPimXg5G1svH67uFYXF3E9WltKJKsr2Qj6/wEFq0LMeZj4MGwaZ3Uqp1XAz5vJ/erz5D3mMjHjVh/l0bHrORsRa3EyEbGCyo2IVXKXhL6/QZX24EiEBc/DvMGQGGd1Mqf01+3i3/asjZ+3O1uOX6bViD9ZffCC1dFEJJOp3IhYyTMndBwPzd4GbLBlAkxsDVFnLQ7mvJpXKcT8QSFULuzPxZh4ev6wgRHLD+JwZKsr8CLZmsqNiNVsNgh5DrpPA68AOLEBxjaBk1usTua0SuTNycxnGtK1XiCGAZ8vPUCfCZu4HBNvdTQRyQR3XG4aN27Mp59+yoEDBzIyj0j2Ve4BcxxOvgoQdQbGPwTbfrY6ldPy9nBjeIfqfNopCG8POysPnKfViD/ZFnbZ6mgiksHuuNw88cQTrFu3jtq1a1OpUiVefvll1qxZo5lBRdJT3jLQbxlUaAVJ8TDnGVj0MiQlWJ3MaXWsXYxZzzSiVL6cnI6IpfO365i49pj+7RJxYam+FTwuLo7ly5czZ84c5s2bR1JSEq1ataJNmzY0b94cHx+fjMqaLnQruDgFhwNWfQwrhpvPSzY2VxbPmc/SWM4sKjaBl2fsZGGoOZ6p9fVZjX01q7GIU8jUeW42bNjA3LlzmTt3LocPH+b+++9n6NChNGrU6G7eNsOo3IhT2TvfXLYhPhoCisOjk8x1qiRNDMPghzXHGH59VuMy+XMyukdtyhf0szqaiPwHyybxO3z4MHPnziUwMJCOHTum19umK5UbcTrh+2BqV7h0BNx9oO1IqJY1/345iy3HLzHg522cjYzFx8ON4R2q0a5mUatjiUgKNENxClRuxClduwIz+sKhZebzek/Cg++Bu6elsZzZxeg4hvyynT+vz4PTvX5x3ni4Mt4ebhYnE5Fb0QzFIq7GJxd0+xUa/898vvFbmNASIk5aGsuZ5fX1YkKfegxuWg6bDX7eEEanMes4cemq1dFE5C6p3Ig4C7sbNH0Tuk4F7wA4uQm+vQcO/251MqflZrfx3APlGf9YXXLn8CD0VAQPf72a5XvPWR1NRO6Cyo2Is6nQAvqvNAcWX70IP3WAFR+Zd1hJmjSpUID5zzamRmAuIq4l0HfiZj5evI/EJB1TEWd0V+XGMAzNFSFihTyl4PHfoFZvwIAVH8DkTnD1ktXJnFbRXD78+mQwjzUsCcA3Kw7TY9wGwqO0+KaIs0lTuRk3bhxVq1bF29sbb29vqlatyvfff5/e2UQkJR7e0GYEtBtt3kV1aJl5mUrLNqSZp7udt9tUYWS3muT0dGP9kUs8PGI1G45ctDqaiKRCqsvNm2++yeDBg2ndujXTpk1j2rRptG7dmueee44333wzIzKKSEpqdDNnNc5TGiJOwA/NYeN3oLOqafZw9SLMGRhC+YK+hEfF0e37DXy78rDOVIs4iVTfCp4/f35GjBhB165db9g+ZcoUBg0axIULF9I1YHrTreDismIjYM4A2DvPfF6tEzz8JXj5WhrLmV2NT+S1WbuYte0UAA9ULsinnYII8PGwOJlI9pOht4InJCRQp06dm7bXrl2bxMTE1L6diKQX7wDo/BM8+D7Y3CB0GnzfFM5rsdu0yuHpzuedg/igfTU83ews3XOO1l+vZtepCKujiUgKUl1uevbsyejRo2/aPnbsWLp3754uoUQkjWw2aDgQHpsPvoXg/D747j7YNcPqZE7LZrPRrX5xZjzdkGK5fQi7dJUOo9cyZWOYLlOJZFGpviw1aNAgfvzxRwIDA2nQoAFgri8VFhZGr1698PD4+3Tt559/nr5p04EuS0m2EXXOnNX42J/m89qPQfPh4JnD0ljOLOJqAs//up3l+8IBaF+zKO+1q0pOLb4pkuEydPmF++677472s9ls/P571ptcTOVGspWkRPM28T8/BwzIXwk6jYcClaxO5rQcDoNvVx3h09/2k6TFN0UyjdaWSoHKjWRLh/+Amf0hJty8bbzFR1Crl3kZS9Jk49FLDJqylXORcXh72HmvXTU61i5mdSwRl6W1pUTkRmXug6fXQJn7IfEazHsWpj9u3mElaVKvVB4WPNuYxuXyEZvg4IVpO3hx2g6uxSdZHU0k21O5EckufAtA9xnQ7B2wu8Pumeakf6c06V9a5fP1YmKfevzvgfLYbTBty0najVrDofBoq6OJZGsqNyLZid0OIUOgz2LIVRwuH4NxD8Lar7U2VRrZ7TYGNS3HpH71yefrxf5zUbQZuZrZ1+fGEZHMp3Ijkh0F1oUn/4RKbcCRCL+9DlO6QEzWnoQzK2tYJh8LB4cQXDovV+OTGPLLdobODCU2QZepRDKbyo1IduWTCzr/CA9/Ae7ecPA3GN0Ijq6yOpnTKuDnzaR+9Xm2aTlsNpiyMYz236zl6IUYq6OJZCsqNyLZmc0GdR6HJ36HfBUg+ixMbAN/fGDeRi6p5ma38fwD5fnx8XrkzenJ3jORtP56NfN3nrY6mki2YWm5GT16NNWrV8ff3x9/f3+Cg4NZtGjRbfefMGECNpvthoe3t3cmJhZxUQWrQP8/oGZPwICVH8HE1nDlhNXJnFbjcvlZOLgx9UrmIToukYGTt/HG7F3EJeoylUhGs7TcFCtWjA8//JAtW7awefNm7r//ftq2bcvu3btv+xp/f3/OnDmT/Dh+/HgmJhZxYZ45oe1IeGQcePpB2FrzMlXodKuTOa2C/t5MfqI+zzQpA8BP64/zyOi1hF28anEyEdeW5Sbxy5MnD5988gl9+/a96XsTJkxgyJAhXLly5Y7fLy4ujri4uOTnkZGRBAYGahI/kZRcOgIznoBTm83nVTtCq8/McTqSJn/sD+f5X7Zz+WoCft7ufNIxiIeqFrI6lojTcMpJ/JKSkpg6dSoxMTEEBwffdr/o6GhKlChBYGDgf57lARg+fDgBAQHJj8DAwPSOLuJ68pSGx5dAk6HmCuO7pmuw8V26r0IBFjzbmNolchMVm8hTk7bwzrzdxCfqFnyR9Gb5mZvQ0FCCg4OJjY3F19eXyZMn07Jly1vuu27dOg4ePEj16tWJiIjg008/ZdWqVezevZtixW497bnO3IjcpZObYeYT5tkcbNBwENz/Orh7WZ3MKSUkOfhkyX7GrjoCQFBgLkZ2rUlgHi1oKpISp1pbKj4+nrCwMCIiIpg+fTrff/89K1eupHLlyv/52oSEBCpVqkTXrl1599137+jztLaUSBrERcOSV2HrRPN5wWrQYSwU/O+/p3Jry/ac43/TdhBxLQF/b3c+61yDByoXtDqWSJblVOXm35o1a0aZMmX49ttv72j/Tp064e7uzpQpU+5of5UbkbuwbwHMHQRXL4KbFzzwDtR70pz5WFLtxKWrDJyyjR0nrgDQ/57SvNi8Ah5uOp4i/+aUY27+4nA4briMlJKkpCRCQ0MpXLhwBqcSEQAqtoKn10G5ByEpDha/ApM6QKTmcEmLwDw5mPZkMI83KgXA2FVHeHTsek5fuWZxMhHnZmm5GTp0KKtWreLYsWOEhoYydOhQVqxYQffu3QHo1asXQ4cOTd5/2LBh/Pbbbxw5coStW7fSo0cPjh8/Tr9+/az6EUSyH7+C0O1X8+4pdx848gd8Ewy7Z1udzCl5utt5s3VlxvSohZ+3O1uOX6bViD/5Y3+41dFEnJal5SY8PJxevXpRoUIFmjZtyqZNm1iyZAkPPPAAAGFhYZw5cyZ5/8uXL/PEE09QqVIlWrZsSWRkJGvXrr2j8Tkiko5sNqjbD55cBYVrQOwVmNYbZj0NsZFWp3NKD1UtzIJBjala1J/LVxPoM34THy3eR2KS7qYSSa0sN+Ymo2nMjUg6S4w3ZzRe/TkYDnO18XajoWSI1cmcUmxCEu8v2MtP680JSuuVzMOIrjUpFKDZ2CV7c+oxNyLiZNw9oekb8NhCs9hcCYMJrWDRyxCvBSNTy9vDjXfbVeXrrjXx9XJn47FLtBzxJ3/s02UqkTulciMi6aNEMDy1Bmr1Np9vGANjQuD4OmtzOanWQUWYNyiEKkX8uRQTT58Jmxi+cC8Jukwl8p9UbkQk/Xj7Q5sR0GMG+Bc1J/4b3wIWvwoJugMotUrly8mMpxvSO7gEAN+uOkLnb9dx8rLWphJJicqNiKS/ss3gmXVQswdgwPpR5lmcExutTuZ0vD3ceKdtVUZ3N++m2hZ2hVYjVvPb7rNWRxPJslRuRCRjeAdA21HQbRr4FYaLh+CH5vDbG5AQa3U6p9OiWmEWPtuYoMBcRFxLoP9P5tpUcYlJVkcTyXJUbkQkY5V/0DyLE9TVvJtq7Qj4trG5ZpWkyl+T/j3R2Jz0b/yaY3QcvY7jFzVwW+SfVG5EJOP55Ib2Y6DrVPAtCBcOwLgHYNnbkHhnM5KLydPdzmutKjOudx1y5fAg9FQED49YzfydmiVa5C8qNyKSeSq0gGfWQ7XO5lmc1V/At/fCqa1WJ3M6TSsVZOGzjalTIjdRcYkMnLyN12aFEpugy1QiKjcikrly5IFHvoMukyBnfji/F75vBsuHaSxOKhXJ5cPU/g0YcF8ZbDb4eUMY7Uat4fD5aKujiVhK5UZErFGpNTyzAao+AkYS/PmZORZH8+KkirubnRebV2Rin3rkzenJvrNRtP56NbO2nbQ6mohlVG5ExDo580LHH6Dzj3+PxRn/EMx/XmtUpdI95fOzaHBjgkvn5Wp8Es/9soMXp+3ganyi1dFEMp3KjYhYr3JbGLABavY0n28eB6Pqw76F1uZyMgX8vZnUrz7PNSuP3QbTtpyk7cg1HDgXZXU0kUylciMiWYNPbmg7EnrNhdylIOo0TO0K0x6DaK2rdKfc7DYGNyvHz/0aUMDPi4Ph0bQZuZpfN50gm62TLNmYyo2IZC2l7zXnxWk0BGxusHsWjKwL2yaBfjnfseAyeVk4uDH3lM9PbIKDl2bs5LlfthMdp8tU4vpsRjar8qlZMl1ELHZ6O8wdBGd3ms9L3Qutv4Q8pa1M5VQcDoMxqw7z2W8HSHIYlMqXk5HdalKlSIDV0URSJTW/v3XmRkSyriI14Ik/4IFh4O4NR1fCNw1hzQhI0hmIO2G323imSVl+6d+AIgHeHL0QQ/tv1vLTumO6TCUuS2duRMQ5XDwM84fA0VXm88JB0GYkFK5uaSxncjkmnhen72DZXnMMU8tqhfjwker4e3tYnEzkv+nMjYi4nrxlzMHGbUeZi3Ke2QFjm8DStyD+qtXpnELunJ5816sOr7eqhIebjYWhZ2k14k92nLhidTSRdKVyIyLOw2aDmj1gwCao3M6c/G/Nl/BNfdi/2Op0TsFms9GvcWmmPdWQwDw+nLh0jUdGr+W7VUdwOLLViXxxYbosJSLOa99CWPgiRF6fjbfiw/DQh5Ar0NpcTiLiWgJDZ+5kYehZAJpUyM+nnYLI5+tlcTKRm+mylIhkDxVbmpP/NXwW7O6wb745+d+aryApwep0WV6AjwejutXig/bV8HK3s2L/eVp+9SdrD12wOprIXdGZGxFxDef2wILnIez62lQFKkOrz6FEsLW5nMT+s1EMnLyVg+HR2GzwTJMyPNesPO5u+v/AkjXozI2IZD8FK8NjC80Bxz55IHyPuU7V7AEQc9HqdFlehUJ+zB0YQtd6xTEMGPXHYbqMXc/JyxqsLc5HZ25ExPVcvQTL3oKtP5rPfXJDs3fMtavs+v90/2XBzjO8MmMnUXGJ+Hu783HH6jxUtbDVsSSbS83vb5UbEXFdYRvMS1XndpnPi9WDhz+HQtWszeUETly6yqAp29h+/TbxHg2K83qrynh7uFkbTLItlZsUqNyIZDNJibBhDKwYDvHR5npV9Z+C+4aCl5/V6bK0hCQHn/12gDErDwNQsZAfI7vVpGwBHTfJfCo3KVC5EcmmIk7B4ldg71zzuV8RaP4+VGlvzp8jt7XqwHme/3U7F6Lj8faw806bKnSuE4hNx00ykcpNClRuRLK5g0th4Qtw+Zj5vEQItPgIClW1NFZWFx4Vy/9+3cGfB83bxFsHFeH99lW1dINkGpWbFKjciAgJ12D1l+bsxomxYLND3X7QZCjkyGN1uizL4TD4dtURPvttP4kOg+J5cvB115oEBeayOppkAyo3KVC5EZFkV8JgyWt/X6ryyQNN34BavcGugbO3szXsMs9O2cbJy9dwt9t46aEK9Aspjd2uy1SScVRuUqByIyI3ObICFr0C5/eazwtVh5afQPEGlsbKyiKuJfDqzFAWhJ4B4N7y+fmss5ZukIyjcpMClRsRuaWkBNg0Dv74AOIizG3VOsMD74B/EWuzZVGGYTBl4wnembebuEQH+f28+KJzDULK5bM6mrgglZsUqNyISIpiLsDyYdcnADTAIyfc8wIEDwB3nZW4lf1noxg0ZSsHzv29dMOQZuXx0NINko5UblKgciMid+T0Nlj4EpzcaD7PU9pccbx8c2tzZVHX4pMYNn8PUzaGAVCreC5GdK1Jsdw5LE4mrkLlJgUqNyJyxxwOCP0Vlr4J0efMbeUehObDIV9Za7NlUQt2nuGVmTuJijWXbvjokeq0qKalG+TuqdykQOVGRFItLgpWfQLrvgFHAtjdzVvH731Zt47fwr+XbuharzhvPlwZH0/dgSZpp3KTApUbEUmzC4dgyatwcIn53DuXWXDq9gN3T0ujZTUJSQ4+X2ou3WAYULaALyMerUnlIvp3V9JG5SYFKjcictcO/w5LXofw3ebzPKXhgXehYist5fAvaw5d4LlfthMeFYenu51XW1Skd8OSWrpBUk3lJgUqNyKSLhxJsG0S/P4exISb20o0MterKlLT2mxZzMXoOF6avpPl+8zj1LRiAT7uWJ28mhNHUkHlJgUqNyKSruKiYM1XsPZrcykHgKCucP8bEFDU2mxZiGEY/LjuOO8v3Et8ooMCfl580aUGjcpqThy5Myo3KVC5EZEMEXHSnB9n5y/mc3cfaPQsNHwWvHytzZaF7D0TyaAp2zgUbs6J8+Q9Zfjfg5oTR/6byk0KVG5EJEOd2mKuVxW2znzuW8hcryqoq9aruu7fc+IEFQtgRNealMib0+JkkpWp3KRA5UZEMpxhmItxLn0TLh8ztxWqBg++B6WbWJksS1m86wwvTd9JZGwiOT3deK99VdrXLGZ1LMmiVG5SoHIjIpkmMQ42joWVn/y9XlWZptDsbShc3dJoWcWpK9d4bup2Nh67BED7mkUZ1rYKft4eFieTrEblJgUqNyKS6WIuwsqPYPMP5iSAYC7Kef9rkLukpdGygiSHwcjfD/HV8gM4DCieJwcjutakRmAuq6NJFqJykwKVGxGxzKUj8Pv7sGu6+dzuYU4AeM8LkFN3DW0+donBU7dz6so13O02nn+wPE/dUwa7XXPiSOp+f1s6PH306NFUr14df39//P39CQ4OZtGiRSm+Ztq0aVSsWBFvb2+qVavGwoULMymtiMhdylMaOo6D/iuh9H3mWZwNo+GrGualq/gYqxNaqk7JPCwc3JhW1QqT6DD4ePF+ev6wgXORsVZHEydjabkpVqwYH374IVu2bGHz5s3cf//9tG3blt27d99y/7Vr19K1a1f69u3Ltm3baNeuHe3atWPXrl2ZnFxE5C4UqQG9ZkPPWVA4COKj4I/3YERN2DQOkhKsTmiZAB8PRnaryUePVMPHw401hy7S4qs/Wb73nNXRxIlkuctSefLk4ZNPPqFv3743fa9Lly7ExMQwf/785G0NGjSgRo0ajBkz5pbvFxcXR1xcXPLzyMhIAgMDdVlKRLIGhwN2z4Tf3/37zqo8ZaDpm1C5bbZezuFQeDTPTtnGnjORADzWsCSvtKiIt4duqc+OnOay1D8lJSUxdepUYmJiCA4OvuU+69ato1mzZjdsa968OevWrbvt+w4fPpyAgIDkR2BgYLrmFhG5K3Y7VOsIAzZBi08gRz64dBim9Ybvm8LRP61OaJmyBXyZNaAhjzcqBcCEtcdoN2oNB89FWZxMsjrLy01oaCi+vr54eXnx1FNPMWvWLCpXrnzLfc+ePUvBggVv2FawYEHOnj172/cfOnQoERERyY8TJ06ka34RkXTh7gn1+8Pg7eZK4x45zQkBJz4MP3WAU1utTmgJL3c33mxdmfGP1SVvTk/2nY2i9cjVTFp/nCx24UGyEMvLTYUKFdi+fTsbNmzg6aefpnfv3uzZsyfd3t/Lyyt5wPJfDxGRLMvLD+57FZ7dZt5JZXeHw8vhu/tganc4l37/PjqT+yoWYNGQxjQul4/YBAevz95F/5+2cCkm3upokgVZXm48PT0pW7YstWvXZvjw4QQFBfHVV1/dct9ChQpx7tyNg8rOnTtHoUKFMiOqiEjm8SsIrT6DgZug+qOADfbNh9ENYXpfuHjY6oSZroCfNxP71OP1VpXwdLOzdM85mn+5ij8Pnrc6mmQxlpebf3M4HDcMAP6n4OBgli9ffsO2pUuX3naMjoiI08tTGjp8C8+sNwcYY5jz5IysC3MGwpUwqxNmKrvdRr/GpZk1oCFlC/hyPiqOnuM28t78PcQlJlkdT7IIS8vN0KFDWbVqFceOHSM0NJShQ4eyYsUKunfvDkCvXr0YOnRo8v6DBw9m8eLFfPbZZ+zbt4+3336bzZs3M3DgQKt+BBGRzFGgInT+0Zwjp9yDYCTBtp/g69qw8EWIuv3YQ1dUpUgA8waG0KNBcQC+X32U9qPWcihcg43F4nITHh5Or169qFChAk2bNmXTpk0sWbKEBx54AICwsDDOnDmTvH/Dhg2ZPHkyY8eOJSgoiOnTpzN79myqVq1q1Y8gIpK5itSA7tPg8d+gZGNIijfXr/qqBvz2Bly9ZHXCTOPj6cZ77arxXa865MnpyZ4zkbQasZqfNNg428ty89xkNC2/ICIu5chKc46ck5vM555+EDwAgp8B7wBrs2Wi8MhY/jdtB38evABAs0oF+eiRauT19bI4maQXrS2VApUbEXE5hgEHfzNLztlQc5t3Lmg4EOo9Cd7Z4986h8PghzVH+XjxfuKTHOT38+LzzkE0Lpff6miSDlRuUqByIyIuy+GAvXPgjw/gwgFzWzYsObtPRzB46nYOhUcD0C+kFC8+VAEvd81s7MxUblKgciMiLs+RBLtmwsqP4OJBc5t3LggeCPWzR8m5Fp/EBwv38tP64wBUKuzPiEdrUK6gn8XJJK1UblKgciMi2YYjCXbPMkvOP8/kBA80Z0POBmNylu05x0szdnIpJh4vdzuvP1yZHvWLY8vGa3Y5K5WbFKjciEi2c8uSE/CPMzmuXXJuHmxcgI8eqa7Bxk5G5SYFKjcikm0ll5yP4cJ+c5t3ADQYAA2ecumS43AYjF97jI8W7UsebPxZpyDuKa/Bxs5C5SYFKjciku3dtuQ8A/WfAp9clsbLSHtORzJ46jYOXh9s3DekFC9psLFTULlJgcqNiMh1jiTYM9ssOef3mds8/aBeP/Nsjq9rntWITTAHG/+4ToONnYnKTQpUbkRE/sXhgD2zYNWnEH591XF3H6jdGxo+CwFFrc2XQZbvPceL0/8x2LhVJXo0KKHBxlmUyk0KVG5ERG7D4YADi8ySc3qruc3uATW6Qshz5iKeLiY8KpYXpu1k1QFzZXENNs66VG5SoHIjIvIfDAOO/AGrPoPjq81tNjtUfQQa/w8KVLI2XzpzOAwmrD3Gh9cHG+fz9eKzzkHcq8HGWYrKTQpUbkREUuH4OvjzMzi09O9tFR+Ge16AIjWty5UB9p6J5Nkpfw82fryROdjY20ODjbMClZsUqNyIiKTB6W1mydk77+9tZZqaZ3JKNAQXGacSm5DE8IV7mXh9sHHFQn6M6FqT8hpsbDmVmxSo3IiI3IXwfbD6cwidDkaSua1YPQgZAuVbgN1uabz08vu+c7w4bScXrw82HtqiIr0bltRgYwup3KRA5UZEJB1cOgprvoTtkyEp3tyWrwI0ehaqdQZ3T0vjpYfwqFhemr6TFfvNwcb3lM/PJx2rU9Df2+Jk2ZPKTQpUbkRE0lHUWVg/Gjb/AHGR5ja/wuaEgLUfc/pFOg3DYNL647y3YC9xiQ5y5/BgeIdqPFS1sNXRsh2VmxSo3IiIZIDYCNg83iw60WfNbV4BULcvNHgafAtYm+8uHQqPYsgv29l1yixwnesU483WVfD1crc4WfahcpMClRsRkQyUGAc7f4E1I+DiQXObmxfU6AYNB0HeMtbmuwvxiQ6+WHaAMSsPYxhQPE8OvugSRO0SeayOli2o3KRA5UZEJBM4HLB/Aaz+Ek5tNrfZ7FCpjTnrcbHalsa7GxuOXOT5X3dw6so17DYYeF9ZBjUth4ebawymzqpUblKgciMikokMA46vMUvOP+fKKR4MwQOhQguwO988MpGxCbw9Zzczt50CIKhYAF90qUHp/L4WJ3NdKjcpULkREbHI2V2wbiSETgNHorktTxkIfgaCuoFnDmvzpcG8Had5bVYokbGJ+Hi48cbDlelaL1C3jGcAlZsUqNyIiFgs8jRs+Ba2jDcHIgP45DEHH9d9AvwKWpsvlc5EXON/v+5g7eGLgLk+1YePVCef1qdKVyo3KVC5ERHJIuKiYdskWP8NXDFnBMbNE6p3Ni9ZOdEaVg6HwQ9rjvLx4v3X16fy5KNHqtO0knMVtaxM5SYFKjciIlmMI8lc1mHdSDi56e/tZZuZJad0E6dZ3mHvmUiGTN3O/nNRAHSvX5zXWlUih6duGb9bKjcpULkREcnCwjbAuq9h73zg+q+nApWh/lPmGR0PH0vj3YnYhCQ+XbKf71cfBaB0vpx80aUGQYG5rA3m5FRuUqByIyLiBC4dMScE3PYzJMSY23zyQJ0+ULcf+BexNt8dWHPoAv/7dQdnI2Nxt9sY3LQcTzcpg7tuGU8TlZsUqNyIiDiRa1dg20+wYSxEhJnb7O5QuZ25xEMWny/nytV4Xpu1iwWhZwCoXSI3X3SuQfG8zndnmNVUblKgciMi4oSSEmH/QvNsTtjav7cXq2su71CpDbh5WJcvBYZhMGvbKd6cs5vouERyerrxdpsqdKxdTLeMp4LKTQpUbkREnNzp7bBhDIROB0eCuc2/qHkrea3ekDOfpfFu58Slq/zv1x1sPHYJgBZVC/FB+2rkzun8K6hnBpWbFKjciIi4iKhz5mrkm8dBzHlzm5snVOkA9fpnyUtWSQ6Db1cd5vPfDpDoMCjg58WnnYK4p3x+q6NleSo3KVC5ERFxMQmxsHumOTHgme1/by9S05wUsGqHLHeXVejJCIb8so3D583B0o81LMkrLSri7eF8S1FkFpWbFKjciIi4KMOAU1tg43dm2UmKN7f75IaaPc3LVrlLWhrxn67FJzF80V5+XGdOYFiugC9fPlqDKkUCLE6WNancpEDlRkQkG4i5AFt/NC9bRZy4vtEG5R6Eek9AmaZgzxq3ZP+xL5wXp+/kQnQcHm42XniwAv0al8bNrsHG/6RykwKVGxGRbMSRBAcWm2dzjvzx9/Y8paFOX6jZ3TyzY7GL0XG8MjOUpXvOAVC/VB4+6xxEsdy6ZfwvKjcpULkREcmmLhyCTd/D9skQd33BTncfqNbRPJtTOMjSeIZh8OvmE7wzbw9X45Pw9XLnrdaVdcv4dSo3KVC5ERHJ5uJjYOevZtE5t+vv7cXqmXdZVW4D7tat6H38YgzP/7qDLccvA9C8SkE+aF+NvNl8lXGVmxSo3IiICGAOQA5bZ16y2jsXHInm9pz5zfly6vSBgGKWREtyGIxZeZgvlx0gIckgn68XHz1SLVuvMq5ykwKVGxERuUnUWdgyEbaMhyhzqQRsdijXHGo/BuUeAHvm36a961QEz/2ynYPh0QB0rRfI660qk9Mr+60yrnKTApUbERG5raQE2LfAvGR17M+/t/sXNW8nr9Uz08/m/LXK+Lg1RzEMKJ4nB593DqJOyTyZmsNqKjcpULkREZE7cv4AbJ1oDkC+Zi6ZgM0OZR8wL1mVfQDcMu8MyrrDF3lh2g5OXbmG3QZP3VuGIc3K4+meNW5pz2gqNylQuRERkVRJiIV982HLhBvP5vgVgZo9zLM5uYpnSpTI2ATenrubmVtPAVC5sD9fdKlBhUJ+mfL5VlK5SYHKjYiIpNmFQ7B1gnk25+rF6xttULaZeTanXPNMOZuzKPQMr84K5fLVBDzd7bzUvAKPNyqF3YUn/lO5SYHKjYiI3LXEuL/P5hxd9fd230LmmZyaPSF3iQyNEB4Vy8vTd/LHfnPR0Aal8/BpJ9ed+E/lJgUqNyIikq4uHjbH5mz7Ga5euL7RBmWbmndalX8I3Dwy5KMNw2DKxhO8t8Cc+M/Py52321ShQ62iLjfxn8pNClRuREQkQyTGw/4F5tmcIyv+3p4zPwQ9CjV6QIGKGfLRxy7E8Pyv29kadgWAh6oU4oMO1ciT0zNDPs8KKjcpULkREZEMd+mIuXDntp8hJvzv7cXqmoOQq3QA7/T9HZSY5ODbVUf4YukBEh3mxH+fdKzOfRULpOvnWCU1v78tvX9s+PDh1K1bFz8/PwoUKEC7du3Yv39/iq+ZMGECNpvthoe3t3cmJRYREbkDeUpDs7fh+T3w6BSo0ApsbnByE8wbDJ+Wh1lPwbHV5kzJ6cDdzc6A+8oye0AjyhXw5UJ0HH0mbOLVWaHExCWmy2c4C0vLzcqVKxkwYADr169n6dKlJCQk8OCDDxITE5Pi6/z9/Tlz5kzy4/jx45mUWEREJBXcPKBiS+g6GZ7fCw+8C/nKQ+I12DEFJrSCETVh1ScQcSpdPrJq0QDmDQqhb0gpACZvCKPliD+T16rKDrLUZanz589ToEABVq5cyT333HPLfSZMmMCQIUO4cuVKmj5Dl6VERMRShgEnN8O2n2DXTIiPMrfb7FDmfvOyVYWW6bJ459rDF3jh1x2cjojFboOnm5RhcFPnnPjPaS5L/VtEhLkEfZ48KU8pHR0dTYkSJQgMDKRt27bs3r37tvvGxcURGRl5w0NERMQyNhsE1oU2I+CF/dBuDJQIAcMBh5bBtMfgswqw8CU4ve2uLls1LJOPRUPuoUPNojgMGPXHYdp/s4YD56LS7+fJgrLMmRuHw0GbNm24cuUKq1evvu1+69at4+DBg1SvXp2IiAg+/fRTVq1axe7duylW7Ob1Pt5++23eeeedm7brzI2IiGQpFw+bkwNunwxRp//enr+SebdV9c7gXyTNb7/w+sR/V65P/PfCg+XpG1IaNyeZ+M8p75Z6+umnWbRoEatXr75lSbmdhIQEKlWqRNeuXXn33Xdv+n5cXBxxcXHJzyMjIwkMDFS5ERGRrMmRBId/N0vOvgWQ9NfvMBuUbgJBXaHSw+CZM9VvHR4ZyyszQ/l9n3kHV72S5sR/xfNm/Yn/nK7cDBw4kDlz5rBq1SpKlSqV6td36tQJd3d3pkyZ8p/7asyNiIg4jWtXYM8c2DEVwtb+vd3TFyq3Nc/olAgB+52PMjEMg183n2DYvD3ExCeRw9ON11pVolu94ll64j+nKTeGYTBo0CBmzZrFihUrKFeuXKrfIykpiSpVqtCyZUs+//zz/9xf5UZERJzSpaOw8xfzLqvLx/7eHhBoXrIK6gr57vz36IlLV3lh2g42HDVXPL+3fH4+eqQ6hQKy5vQqTlNunnnmGSZPnsycOXOoUKFC8vaAgAB8fHwA6NWrF0WLFmX48OEADBs2jAYNGlC2bFmuXLnCJ598wuzZs9myZQuVK1f+z89UuREREadmGHBig3nZavdsiIv4+3tF65hnc6o+AjlSvjkHwOEw+GHNUT5esp/4RAf+3u68264qbYKKZLmzOE5Tbm534MaPH89jjz0GQJMmTShZsiQTJkwA4LnnnmPmzJmcPXuW3LlzU7t2bd577z1q1qx5R5+pciMiIi4j4RrsX2Retjq0DIwkc7vdAyo8BNUfhXIP/Odt5YfCo3j+1x3sPGkWpZbVCvFeu6y1fIPTlBsrqNyIiIhLig6H0GnmZauzoX9v985ljs+p3hmKN7zt+JyEJAejVxxmxPKDycs3fNihGs0qF8yc/P9B5SYFKjciIuLyzu4yS86uGRB15u/t/kXNS1bVOkGhauacO/+y61QEz/2ynYPh0QB0ql2MN1tXxs87Y1Y2v1MqNylQuRERkWzDkWSuXxX6K+yZd+P4nPwVzZJTrSPkLnnDy2ITkvhi6QHG/nkEw4CiuXz4pGN1GpbNl7n5/0HlJgUqNyIiki0lxMLB38yic2AJJMX//b3A+mbRqdIecv5dYDYevcQL03YQdukqAI81LMnLD1XEx9Mts9Or3KRE5UZERLK9a1dg7zyz6Bz9E7heBezu5vpW1TqbC3565iQmLpH3F+5l8oYwAErny8lnnYOoWTx3pkZWuUmByo2IiMg/RJ4xx+aEToMz2//e7pEDyjc3x+iUfYAVRyJ5ecZOzkXGYbfBM03K8mzTcpm2CKfKTQpUbkRERG7j/AGz5IROg8tH/97u6QcVWxFTrg1v7MrPzB3nAahc2J/PuwRRsVDG/z5VuUmByo2IiMh/MAxzRfJdM8yJAiNP/v0971ycKNSM949XYum18tjd3HnugfI8eU+ZDF2EU+UmBSo3IiIiqeBwwMmNsGsm7J4FMeHJ34q052J2fF3mJzUgqVh9Pu1Si1L5Ur+g551QuUmByo2IiEgaOZLg+Bqz6OyZA9cuJX/rrJGbJUYw+YO78tCDD2N3S9+xOCo3KVC5ERERSQdJCXB0JeyaiWPPPOzxkcnf2u1di8ov/56u61Ol5vd35gxxFhEREdfi5gFlm0G7b7C/dAhHl8kcLtSCGMOLqwVqWrrwprtlnywiIiKuwd0Le6VWlKnUirCz56npZ229ULkRERGRdFO8UH6rI+iylIiIiLgWlRsRERFxKSo3IiIi4lJUbkRERMSlqNyIiIiIS1G5EREREZeiciMiIiIuReVGREREXIrKjYiIiLgUlRsRERFxKSo3IiIi4lJUbkRERMSlqNyIiIiIS8l2q4IbhgFAZGSkxUlERETkTv31e/uv3+MpyXblJioqCoDAwECLk4iIiEhqRUVFERAQkOI+NuNOKpALcTgcnD59Gj8/P2w2m9Vx0lVkZCSBgYGcOHECf39/q+NkGTout6bjcms6LjfTMbk1HZdby6jjYhgGUVFRFClSBLs95VE12e7Mjd1up1ixYlbHyFD+/v76i3YLOi63puNyazouN9MxuTUdl1vLiOPyX2ds/qIBxSIiIuJSVG5ERETEpajcuBAvLy/eeustvLy8rI6Spei43JqOy63puNxMx+TWdFxuLSscl2w3oFhERERcm87ciIiIiEtRuRERERGXonIjIiIiLkXlRkRERFyKyo2TGTVqFCVLlsTb25v69euzcePGO3rd1KlTsdlstGvXLmMDWiQ1x2XChAnYbLYbHt7e3pmYNvOk9s/LlStXGDBgAIULF8bLy4vy5cuzcOHCTEqbOVJzTJo0aXLTnxWbzUarVq0yMXHmSO2flS+//JIKFSrg4+NDYGAgzz33HLGxsZmUNvOk5rgkJCQwbNgwypQpg7e3N0FBQSxevDgT02a8VatW0bp1a4oUKYLNZmP27Nn/+ZoVK1ZQq1YtvLy8KFu2LBMmTMjwnBjiNKZOnWp4enoaP/zwg7F7927jiSeeMHLlymWcO3cuxdcdPXrUKFq0qNG4cWOjbdu2mRM2E6X2uIwfP97w9/c3zpw5k/w4e/ZsJqfOeKk9LnFxcUadOnWMli1bGqtXrzaOHj1qrFixwti+fXsmJ884qT0mFy9evOHPya5duww3Nzdj/PjxmRs8g6X2uPz888+Gl5eX8fPPPxtHjx41lixZYhQuXNh47rnnMjl5xkrtcXnppZeMIkWKGAsWLDAOHz5sfPPNN4a3t7exdevWTE6ecRYuXGi89tprxsyZMw3AmDVrVor7HzlyxMiRI4fx/PPPG3v27DG+/vprw83NzVi8eHGG5lS5cSL16tUzBgwYkPw8KSnJKFKkiDF8+PDbviYxMdFo2LCh8f333xu9e/d2yXKT2uMyfvx4IyAgIJPSWSe1x2X06NFG6dKljfj4+MyKmOnS8nfon7744gvDz8/PiI6OzqiIlkjtcRkwYIBx//3337Dt+eefNxo1apShOTNbao9L4cKFjZEjR96wrUOHDkb37t0zNKdV7qTcvPTSS0aVKlVu2NalSxejefPmGZjMMHRZyknEx8ezZcsWmjVrlrzNbrfTrFkz1q1bd9vXDRs2jAIFCtC3b9/MiJnp0npcoqOjKVGiBIGBgbRt25bdu3dnRtxMk5bjMnfuXIKDgxkwYAAFCxakatWqfPDBByQlJWVW7AyV1j8r/zRu3DgeffRRcubMmVExM11ajkvDhg3ZsmVL8iWaI0eOsHDhQlq2bJkpmTNDWo5LXFzcTZe4fXx8WL16dYZmzcrWrVt3wzEEaN68+R3/nUsrlRsnceHCBZKSkihYsOAN2wsWLMjZs2dv+ZrVq1czbtw4vvvuu8yIaIm0HJcKFSrwww8/MGfOHCZNmoTD4aBhw4acPHkyMyJnirQclyNHjjB9+nSSkpJYuHAhb7zxBp999hnvvfdeZkTOcGk5Jv+0ceNGdu3aRb9+/TIqoiXScly6devGsGHDCAkJwcPDgzJlytCkSRNeffXVzIicKdJyXJo3b87nn3/OwYMHcTgcLF26lJkzZ3LmzJnMiJwlnT179pbHMDIykmvXrmXY56rcuKioqCh69uzJd999R758+ayOk6UEBwfTq1cvatSowb333svMmTPJnz8/3377rdXRLOVwOChQoABjx46ldu3adOnShddee40xY8ZYHS1LGDduHNWqVaNevXpWR7HcihUr+OCDD/jmm2/YunUrM2fOZMGCBbz77rtWR7PUV199Rbly5ahYsSKenp4MHDiQPn36YLfrV21mc7c6gNyZfPny4ebmxrlz527Yfu7cOQoVKnTT/ocPH+bYsWO0bt06eZvD4QDA3d2d/fv3U6ZMmYwNnQlSe1xuxcPDg5o1a3Lo0KGMiGiJtByXwoUL4+HhgZubW/K2SpUqcfbsWeLj4/H09MzQzBntbv6sxMTEMHXqVIYNG5aRES2RluPyxhtv0LNnz+SzWNWqVSMmJob+/fvz2muvucQv87Qcl/z58zN79mxiY2O5ePEiRYoU4ZVXXqF06dKZETlLKlSo0C2Pob+/Pz4+Phn2uc7/JzCb8PT0pHbt2ixfvjx5m8PhYPny5QQHB9+0f8WKFQkNDWX79u3JjzZt2nDfffexfft2AgMDMzN+hkntcbmVpKQkQkNDKVy4cEbFzHRpOS6NGjXi0KFDySUY4MCBAxQuXNjpiw3c3Z+VadOmERcXR48ePTI6ZqZLy3G5evXqTQXmr1JsuMhyhXfz58Xb25uiRYuSmJjIjBkzaNu2bUbHzbKCg4NvOIYAS5cuveN/n9MsQ4crS7qaOnWq4eXlZUyYMMHYs2eP0b9/fyNXrlzJtzH37NnTeOWVV277ele9Wyq1x+Wdd94xlixZYhw+fNjYsmWL8eijjxre3t7G7t27rfoRMkRqj0tYWJjh5+dnDBw40Ni/f78xf/58o0CBAsZ7771n1Y+Q7tL6dygkJMTo0qVLZsfNNKk9Lm+99Zbh5+dnTJkyxThy5Ijx22+/GWXKlDE6d+5s1Y+QIVJ7XNavX2/MmDHDOHz4sLFq1Srj/vvvN0qVKmVcvnzZop8g/UVFRRnbtm0ztm3bZgDG559/bmzbts04fvy4YRiG8corrxg9e/ZM3v+vW8FffPFFY+/evcaoUaN0K7jc7OuvvzaKFy9ueHp6GvXq1TPWr1+f/L17773X6N27921f66rlxjBSd1yGDBmSvG/BggWNli1butQ8FP+U2j8va9euNerXr294eXkZpUuXNt5//30jMTExk1NnrNQek3379hmA8dtvv2Vy0syVmuOSkJBgvP3220aZMmUMb29vIzAw0HjmmWdc6pf4X1JzXFasWGFUqlTJ8PLyMvLmzWv07NnTOHXqlAWpM84ff/xhADc9/joOvXv3Nu69996bXlOjRg3D09PTKF26dKbME2UzDBc5hygiIiKCxtyIiIiIi1G5EREREZeiciMiIiIuReVGREREXIrKjYiIiLgUlRsRERFxKSo3IiIi4lJUbkRERMSlqNyIiIiIS1G5EREREZeiciMiAjRp0oQhQ4bc1XsYhkH//v3JkycPNpuN7du3p0s2EUkdlRsRyXL69OnD66+/bnWMVFu8eDETJkxg/vz5nDlzhqpVq1odSSRbcrc6gIjIPyUlJTF//nwWLFhgdZRUO3z4MIULF6Zhw4a33Sc+Ph5PT89MTCWS/ejMjUg2MWXKFHx8fDhz5kzytj59+lC9enUiIiLS/L7FihXjm2++uWHb2rVryZEjB8ePH0/1+61duxYPDw/q1q17y+83adKEQYMGMWTIEHLnzk3BggX57rvviImJoU+fPvj5+VG2bFkWLVp0w+vi4uJ49tlnKVCgAN7e3oSEhLBp06bb5nA4HAwfPpxSpUrh4+NDUFAQ06dPv+3+jz32GIMGDSIsLAybzUbJkiWT8w4cOJAhQ4aQL18+mjdvzuLFiwkJCSFXrlzkzZuXhx9+mMOHD9/0+R9//DFly5bFy8uL4sWL8/7779/hURTJ3lRuRLKJRx99lPLly/PBBx8A8NZbb7Fs2TIWLVpEQEBAmt+3fv36N5QEwzAYMmQIzz33HCVKlEj1+82dO5fWrVtjs9luu8/EiRPJly8fGzduZNCgQTz99NN06tSJhg0bsnXrVh588EF69uzJ1atXk1/z0ksvMWPGDCZOnMjWrVspW7YszZs359KlS7f8jOHDh/Pjjz8yZswYdu/ezXPPPUePHj1YuXLlLff/6quvGDZsGMWKFePMmTM3HJOJEyfi6enJmjVrGDNmDDExMTz//PNs3ryZ5cuXY7fbad++PQ6HI/k1Q4cO5cMPP+SNN95gz549TJ48mYIFC6b2cIpkT4aIZBvz5s0zvLy8jPfee8/InTu3sWvXruTvtWvXzsiVK5fxyCOPpOo9P/74Y6NKlSrJzydOnGgUKlTIiIqKStP7litXzpg/f/5tv3/vvfcaISEhyc8TExONnDlzGj179kzedubMGQMw1q1bZxiGYURHRxseHh7Gzz//nLxPfHy8UaRIEePjjz9Oft/BgwcbhmEYsbGxRo4cOYy1a9fe8Nl9+/Y1unbtettsX3zxhVGiRImb8tasWTPFn/n8+fMGYISGhhqGYRiRkZGGl5eX8d1336X4OhG5NZ25EclGHn74YSpXrsywYcOYNWsWVapUSf7e4MGD+fHHH1P9ng0aNGDv3r1ER0cTExPDq6++ynvvvYevr2+q33fv3r2cPn2apk2bprhf9erVk792c3Mjb968VKtWLXnbX2c4wsPDAXMsTEJCAo0aNUrex8PDg3r16rF3796b3v/QoUNcvXqVBx54AF9f3+THjz/+eNPloztRu3btG54fPHiQrl27Urp0afz9/ZMvYYWFhQHmcYiLi/vP4yAit6YBxSLZyOLFi9m3bx9JSUk3XeJo0qQJK1asSPV71q5dG7vdztatW1m2bBn58+enT58+aXrfuXPn8sADD+Dt7Z3ifh4eHjc8t9lsN2z765LWPy/zpEZ0dDQACxYsoGjRojd8z8vLK9XvlzNnzhuet27dmhIlSvDdd99RpEgRHA4HVatWJT4+HgAfH5805RYRk87ciGQTW7dupXPnzowbN46mTZvyxhtvpMv75siRg2rVqjFjxgw+/fRTvvjiC+z2tP3TMmfOHNq2bZsuuf6pTJkyyWNe/pKQkMCmTZuoXLnyTftXrlwZLy8vwsLCKFu27A2PwMDAu8py8eJF9u/fz+uvv07Tpk2pVKkSly9fvmGfcuXK4ePjw/Lly+/qs0SyK525EckGjh07RqtWrXj11VeTL4cEBwezdetWatWqddfv36BBA77++mvatm1LkyZN0vQe4eHhbN68mblz5951nn/LmTMnTz/9NC+++CJ58uShePHifPzxx1y9epW+ffvetL+fnx8vvPACzz33HA6Hg5CQECIiIlizZg3+/v707t07zVly585N3rx5GTt2LIULFyYsLIxXXnnlhn28vb15+eWXeemll/D09KRRo0acP3+e3bt3J+cdOXIks2bNUgESuQWVGxEXd+nSJR566CHatm2b/Eu0fv36tGjRgldffZXFixf/53tMmDCBPn36YBjGLb8fFBSEh4cHn3zySZpzzps3j3r16pEvX740v0dKPvzwQxwOBz179iQqKoo6deqwZMkScufOfcv93333XfLnz8/w4cM5cuQIuXLlolatWrz66qt3lcNutzN16lSeffZZqlatSoUKFRgxYsRNpfCNN97A3d2dN998k9OnT1O4cGGeeuqp5O9fuHAhTeN/RLIDm3G7f61EJNtZsWIFI0eOvGk+l7feeouVK1feduzMfffdR61atfjss89S9b7/1KZNG0JCQnjppZfSnF9EBHTmRkSua9asGTt27CAmJoZixYoxbdo0goODAVi0aBEjR468YX+Hw8H58+cZN24cBw8eZM6cOal+338KCQmha9eu6f+DiUi2ozM3IpImK1as4P7776dixYqMHz+e+vXrWx1JRARQuREREREXo1vBRURExKWo3IiIiIhLUbkRERERl6JyIyIiIi5F5UZERERcisqNiIiIuBSVGxEREXEpKjciIiLiUlRuRERExKWo3IiIiIhL+T+1t+jzuu5qMQAAAABJRU5ErkJggg==", "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": "2025-10-15T23:08:22.341433Z", "iopub.status.busy": "2025-10-15T23:08:22.341288Z", "iopub.status.idle": "2025-10-15T23:08:22.411825Z", "shell.execute_reply": "2025-10-15T23:08:22.411313Z" } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjsAAAG0CAYAAADU2ObLAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjcsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvTLEjVAAAAAlwSFlzAAAPYQAAD2EBqD+naQAAaF5JREFUeJzt3Xd0FOXbxvHvphNIIQkhBELvVWoAUUQ6SrdQRUCsoIKvBaxgAXvBwk8EAQWxgoAKUgRUegm9hV4CAUI6qTvvH6PRKGACSWZ3c33O2XPyzM7O3juG7OXMU2yGYRiIiIiIuCg3qwsQERERKUwKOyIiIuLSFHZERETEpSnsiIiIiEtT2BERERGXprAjIiIiLk1hR0RERFyah9UFOAK73c6pU6fw8/PDZrNZXY6IiIjkgWEYJCUlER4ejpvb5a/fKOwAp06dIiIiwuoyRERE5CocP36cChUqXPZ5hR3Az88PME+Wv7+/xdWIiIhIXiQmJhIREZHzPX45CjuQc+vK399fYUdERMTJ/FcXFHVQFhEREZemsCMiIiIuTWFHREREXJrCjoiIiLg0hR0RERFxaQo7IiIi4tIUdkRERMSlKeyIiIiIS1PYEREREZemsCMiIiIuzdKws3r1arp37054eDg2m4358+fnet5ms13y8frrr+fsExcXx8CBA/H39ycwMJDhw4eTnJxcxJ9EREREHJWlYSclJYVGjRrxwQcfXPL5mJiYXI/p06djs9no27dvzj4DBw5k165dLF26lEWLFrF69WruvffeovoIIiIi4uBshmEYVhcB5lWcefPm0atXr8vu06tXL5KSkli+fDkAe/bsoW7dumzcuJFmzZoBsHjxYrp168aJEycIDw/P03snJiYSEBBAQkKCFgIVEREpQAnJqZyK+pk6bXoV+LHz+v3tNH12zpw5ww8//MDw4cNztq1du5bAwMCcoAPQoUMH3NzcWL9+/WWPlZ6eTmJiYq6HiIiIFCC7nYyor0l5uxl1lg3h919XWFaK04SdmTNn4ufnR58+fXK2nT59mtDQ0Fz7eXh4EBQUxOnTpy97rIkTJxIQEJDziIiIKLS6RUREihXDgOhl2D++Ca/59xCefZLz+BPhHmdZSU4TdqZPn87AgQPx8fG55mONHTuWhISEnMfx48cLoEIREZFi7sQmmNkdPu+L2+ltJBkleM9+O0cH/k7F1rdZVpaHZe+cD7/++iv79u3jyy+/zLU9LCyM2NjYXNuysrKIi4sjLCzsssfz9vbG29u7UGoVEREpds7ug+UTYO8iALJsnszI7MDHRi9eH9KeJjXKWFqeU4SdadOm0bRpUxo1apRre6tWrYiPj2fz5s00bdoUgBUrVmC324mMjLSiVBERkeIj4QSsnAhRc8CwY9jc2B7cjQdOdOK0LYTJ/ZvQtqa1QQcsDjvJyclER0fntA8fPkxUVBRBQUFUrFgRMHtaf/3117z55pv/en2dOnXo0qULI0aMYMqUKWRmZjJy5Ej69euX55FYIiIikk+pcfDrm7BhKmSnm9tq38oXpYYw7rdMACb1bsAtDctZWORfLA07mzZtol27djntMWPGADBkyBBmzJgBwNy5czEMg/79+1/yGLNnz2bkyJG0b98eNzc3+vbty3vvvVfotYuIiBQ7GSmw7kP4/T1I/2Mkc6U20OEFPjtRhme/3wXA093q0K9FRQsLzc1h5tmxkubZERERuYLsTNj6GaycBMlnzG1l60OHF6B6B+ZHneLRL6MAGHVzdR7rVKtIysrr97dT9NkRERERCxgG7J4PK16C8390OwmsBDc/A/VvAzc3lu4+w2NfbwPg7taVGdOxpnX1XobCjoiIiPzb4dWw9Hk4tcVs+4ZA2yeg6VDw8AJgTfQ5Hpq9hWy7QZ8m5Xnu1rrYbDYLi740hR0RERH5S8x2WPYCHDSXZsKrFLQaCa1Hgrdfzm5bjl3gnlmbyMi206VeGK/1bYibm+MFHVDYEREREYALR2DFy7DjK7Pt5gnNhsKNT0Cp3MPHd59K5O7pG0jNyOaGGiG82/86PNwdd55ihR0REZHiLOUcrH4DNn4CdnPYOPVvg5ufhqCq/9r90Nlk7pq+nsS0LJpWKs3/BjfF28O9iIvOH4UdERGR4igjBdZ+CL+/CxlJ5raqN5kjrMIbX/IlJy6kMuiT9ZxLzqBuOX+m390cXy/HjxKOX6GIiIgUnEsNIw9rCB3HQ7WbL/uy2MQ0Bn2ynlMJaVQtU5JZw1sQUMKziIq+Ngo7IiIixYFhmGtXLRsP5w+Y2wIrQfvnoF4fcLt8n5sLKRkMnraBI+dTqVC6BLPviSSklPOsMamwIyIi4uqOrYOlz8Hx9WbbN9jseNxsKHhcObQkpWUy5NMN7DuTRKifN3PuaUm5gBJFUHTBUdgRERFxVWf3m8PI9/1gtj1KmEPIWz8MPv+9YsDFjGyGz9zE9hMJlPb1ZPY9kVQM9i3cmguBwo6IiIirSTptrka+5TMwssHmDk0GQ9unwD9vi3OmZ2Vz3+eb2XA4Dj9vD2YNi6RGWb//fqEDUtgRERFxFelJ5iKda9+HzFRzW61boMPzUCbv61VlZdt5+IutrN5/lhKe7kwf2pwGFQIKqejCp7AjIiLi7LIzYfMMWPUqpJw1t1VoDh1fhEqt8nUou93g8W+2s2TXGbzc3Zh6VzOaVw4q+JqLkMKOiIiIszIM2LPQ7JcTd9DcFlTNvJJTpwfkc50qwzB4ev5O5m09iYebjQ8GNqFNjZCCr7uIKeyIiIg4o2Pr4Odn4cQGs+0bAjc9BU3vBvf8z39jGAYv/bCHLzYcw2aDt++8jo51yxZszRZR2BEREXEm56Jh2fPmnDkAnr5/LNQ5Kk8jrC7n7aX7mfbbYQBe7dOQ7o3CC6Jah6CwIyIi4gySz8KqSbDp0z9GWLlB40Fw07g8j7C6nI9WHuS9FdEAvNC9Lnc0jyiIih2Gwo6IiIgjy0iFdR/Ab39bw6pGZ3N5h9A613z4Gb8f5tXFewF4sktt7r6+yjUf09Eo7IiIiDgiezZs+wJWvAxJp8xt5RpBp5egyo0F8hZfbjzGCwt3A/DwzdV54KZqBXJcR6OwIyIi4miil5vLO5zZabYDKpprWNXve8U1rPLj+6iTPPXdDgDuaVOF0R1rFshxHZHCjoiIiKM4vdMMOQeXm22fALjh/6DFveDpU2Bvs3hnDGO+2oZhwKCWFXn6ljrY8jlM3Zko7IiIiFgtMQZ+eQm2zgYMcPOEFiPgxsfBt2An9PtlbyyjvthKtt3gtqYVmNCjvksHHVDYERERsU56Mqx5D9ZM/mt5h3q9zVtWQVUL/O1+jz7HfZ9vJjPb4NaG5Xi1b0Pc3Fw76IDCjoiISNGzZ8PWz8zOxymx5raISLPzcUSLQnnL9YfOM3zmRjKy7HSsW5a377wO92IQdEBhR0REpGgdWAZLn4VYcxQUpauYw8ivYnmHvNpy7ALDZmwkLdNO25pleH9AYzzdC6ajszNQ2BERESkKp3eaIefgCrPtEwhtn4Tm94CHV6G97Y4TCQyZvoGUjGxaVwvmf4Ob4u3hXmjv54gUdkRERApT0mlY8RJs/ZyczseR98GN/wclShfqW++JSWTw9PUkpWXRvHJpPhnSDB/P4hV0QGFHRESkcGSkwJr34fd3ITPF3Fa3F3R4AYIKf5biA2eSGPTJeuJTM2kUEcj0u5vj61U8v/aL56cWEREpLHY7bJ8Ly1/8a+bjCs2h08tQMbJISjh0NpkBn6znfEoG9cv7M2toC/x88r8SuqtQ2BERESkoh1fDkqfh9HazHVAROr4A9foUWufjfzp2PpUBU9dzNimd2mF+fDYskgDf4ht0QGFHRETk2p2LNjsf7/vRbHv7ww2PQeT9BTrz8X85GX+R/lPXcToxjRqhpfj8nkhKlyy8zs/OQmFHRETkaqXGwarXYONUsGeBzR2aDYWbxkLJkCItJSbhIv0/XsfJ+ItUCSnJ7HsiCSnlXaQ1OCqFHRERkfzKyoBN02DlJEiLN7fV6GROClimVpGXE5uYxoCp6zkWl0rFIF/mjIgk1L/orig5OoUdERGRvDIM2PcT/PwMxB00t4XWg84vQbWbLSnpbFI6/aeu4/C5FMoHlmDOiEjKBZSwpBZHpbAjIiKSFzHb4eenzU7IACXLQLunocld4GbN3DXnk9MZ+Mk6Dp5NoVyAD3PvbUmF0r6W1OLIFHZERESuJOkMrHjxr0kB3b2h1YPQZgz4+FtWVnxqBoOmbWD/mWTK+nvzxYiWRAQp6FyKwo6IiMilZKbBug/g17cgI9ncVq+POSlg6UqWlpaQmsmgaevZE5NISClv5oxoSeWQkpbW5MgUdkRERP7OMGD3fPj5OUg4Zm4r3xQ6TyyySQGvJOFiJoOnr2fnyUSCS3rxxYhIqpUpZXVZDk1hR0RE5E8nt8CScXBsrdn2CzdXJK9/G7hZv0p4Ylomd03fwPYTCQSV9GLOiJbUKOtndVkOz9L/cqtXr6Z79+6Eh4djs9mYP3/+v/bZs2cPPXr0ICAggJIlS9K8eXOOHTuW8/xNN92EzWbL9bj//vuL8FOIiIjTS4yBeQ/A1HZm0PEoYc6VM2oTNLzDIYJOUlomQ6ZvYNvxeAJ9PZl9TyS1whR08sLSKzspKSk0atSIYcOG0adPn389f/DgQdq0acPw4cMZP348/v7+7Nq1Cx+f3HMHjBgxggkTJuS0fX3VQUtERPIg86K5WOdvb/+1WGfDO6H98xBQ3tra/iY5PYuhn25k67F4AkqYQadOOes6RzsbS8NO165d6dq162Wff/rpp+nWrRuvvfZazrZq1ar9az9fX1/CwsLy/L7p6emkp6fntBMTE/P8WhERcQGGATu/hWUvQMJxc1uF5tBlElRoZmlp/5SakcWwTzey6egF/H08mH1PJPXCA6wuy6lYf13uMux2Oz/88AM1a9akc+fOhIaGEhkZeclbXbNnzyYkJIT69eszduxYUlNTr3jsiRMnEhAQkPOIiIgopE8hIiIO59RWmN4Fvh1uBh3/CtB3Ggxf6pBBZ+inG9lwJA4/Hw8+Gx5J/fIKOvllMwzDsLoIAJvNxrx58+jVqxcAp0+fply5cvj6+vLSSy/Rrl07Fi9ezLhx4/jll19o27YtAB9//DGVKlUiPDyc7du38+STT9KiRQu+++67y77Xpa7sREREkJCQgL+/LguKiLikpDOwfAJEzQYM8PSF6x+F1qPAy/G6P6RmZDFsxkbWHYrDz9uDWcNb0LhiaavLciiJiYkEBAT85/e3w47GstvtAPTs2ZPRo0cDcN1117FmzRqmTJmSE3buvffenNc0aNCAcuXK0b59ew4ePHjJW14A3t7eeHtrcTQRkWIhKx3WfQSr34CMJHNbgzvM+XIcqF/O313MyGb4jE2sOxRHKW8PZiroXBOHDTshISF4eHhQt27dXNvr1KnDb7/9dtnXRUaacyBER0dfNuyIiEgxYBiw70dY8jRcOGxuC28CXV+FiBbW1nYFFzOyGTZjI2sPnafUH1d0mijoXBOHDTteXl40b96cffv25dq+f/9+KlW6/MyVUVFRAJQrV64wyxMREUd2ZjcsGQuHVprtUmXNKzkN+znEMPLLuZiRzfCZfwWdmcMUdAqCpWEnOTmZ6OjonPbhw4eJiooiKCiIihUr8vjjj3PnnXdy44035vTZWbhwIStXrgTMoelz5syhW7duBAcHs337dkaPHs2NN95Iw4YNLfpUIiJimdQ4+OUV2DQdjGxw94JWI+GGMeDt2HPSXMzI5p5ZG1lz8DwlvdyZOaw5TSsp6BQESzsor1y5knbt2v1r+5AhQ5gxYwYA06dPZ+LEiZw4cYJatWoxfvx4evbsCcDx48cZNGgQO3fuJCUlhYiICHr37s0zzzyTr47Gee3gJCIiDio7yww4v7wMafHmtjrdoeOLEFTF0tLy4mJGNiNmbeK36HOU9HJn1vAWNK0UZHVZDi+v398OMxrLSgo7IiJO7NAqWPwUxO4222XrQ5eJUOVGa+vKoz+v6Pwe/ecVnRY0q6ygkxdOPxpLRETkii4chZ+fgT0LzHaJ0nDzM9DkbnB3jq+31Iwshs/YxNpDZtCZoaBTKJzjt0FERORPGanw+zvw+7uQlQY2N2h+j7mWla/zBIW/z6NjdkZurltXhURhR0REnINhwK558POzkHjC3Fb5BnMoedl61taWTynpWQydsZENh+NyRl2pM3LhUdgRERHHd3oH/PQUHP1jnrWACOj0EtTtCTabtbXlU0r635aA+GPCQA0vL1wKOyIi4rhS42DFS7D5UzDs4OEDbUZD64cdcomH/2KuXr6BjUcuaAmIIqSwIyIijic7yww4K176ayh53V7Q6UUIrGhlZVctKS2Tuz/dyOajF3IW9bwuItDqsooFhR0REXEsR36Hn56AMzvNdmg9s19OlRusresaJKZlcvf0DWw5Fo+/jwef3xNJwwqBVpdVbCjsiIiIY0iMgaXPwo6vzbZPoDmUvOlQpxlKfikJqZnc9ekGth03g87se1rSoEKA1WUVK8772yMiIq4hKwPWfwSrXoOMZMAGTe+Gm5+FksFWV3dNLqRkMHj6enaeTKS0ryefDY+kfnkFnaKmsCMiItaJXg4/PQnnD5jtCs2h2+sQ3tjaugrA+eR0Bk3bwJ6YRIJKejH7nkjqlNMs/VZQ2BERkaJ34SgsGQd7F5ntkmWgw3ho1N+hVyXPq3PJ6Qycup59Z5IIKeXNnBGR1Czr2AuRujKFHRERKTqZF+H39+C3t/6Y/dgdIu+Dm54CH9e4vRObmMaAT9YTHZtMqJ83c0a0pHpoKavLKtYUdkREpPAZBuz7yVywM/6oua3yDdD1NShb19raCtDphDQGTF3HoXMplAvwYc6IllQJKWl1WcWewo6IiBSu8wfNfjnRS822Xzh0fhnq9Xa62Y+v5GT8RQZMXcfR86mUDyzBFyNaUjHY+SY+dEUKOyIiUjgyUmD1G7D2fcjOADdPaD0KbngMvF3rts7xuFT6T13HiQsXiQgqwZx7WhIRpKDjKBR2RESkYBkG7FkAi8dC4klzW/WO0GUShFS3trZCcORcCgM/Wc/J+ItUCvblixEtCQ8sYXVZ8jcKOyIiUnDOH4Qf/w8OrjDbgZXM2Y9rdnGpW1Z/io5NYsDU9cQmpVM1pCRzRrQkLMDH6rLkHxR2RETk2mWkmiOsfn/XvGXl7g1tHjUX7fR0zascu08lMnjaes6nZFCrrB+f3xNJGT9vq8uSS1DYERGRa7P3R1j8JMQfM9vVO5ijrIKrWVtXIdp+Ip7B0zaQcDGT+uX9mTUskqCSXlaXJZehsCMiIlfnwhFzlNX+xWbbvwJ0mQh1urvkLas/bToSx9BPN5KUnkXjioHMGNqCgBKeVpclV6CwIyIi+ZOZBmveg1/fNCcGdPOE1iPhxsfBy7XnlFlz8Bz3zNxEakY2LaoEMf3u5pTy1lepo9N/IRERybsDy+CnxyHukNmuciN0exPK1LS2riKwcl8s9322mfQsOzfUCOHjwc0o4eVudVmSBwo7IiLy3+KPw5KxsGeh2S4VZk4MWL+vS9+y+tPPu07z0JwtZGYbdKgTyvsDmuDjqaDjLBR2RETk8rIyzEkBV78OmanmWlYtH4C2T4JP8VjBe+G2U4z+Moosu0G3BmG8c2djvDycf7HS4kRhR0RELu3QKnPOnHP7zXbF1nDLG1C2nrV1FaFvNp/giW+2YTegd+PyvH5bQzzcFXScjcKOiIjklnQaloyDnd+a7ZJloOOL0Khfsbhl9afZ64/y9LydAPRvEcHLvRrg5lZ8Pr8rUdgRERGTPRs2ToMVL0J6ItjcoPk90O5pKBFodXVF6uPVB3nlx70A3N26Ms93r4utGAU9V6OwIyIicCoKFo2GU1vMdvmmcMtbEH6dlVUVOcMweGvpfiaviAbg/rbVeLJLLQUdJ6ewIyJSnKUnwS+vwPopYNjB2x/aPwfNhoFb8RptZLcbTFi0mxlrjgDwRJdaPHiT6y1cWhwp7IiIFEeGYQ4j/+lJSDplbqvfFzq/An5h1tZmgaxsO099t4NvNp8A4MWe9RjcqrK1RUmBUdgRESluLhyFn574a5mH0pXhljfNNa2KofSsbB6dG8VPO0/j7mbj9dsa0qdJBavLkgKksCMiUlxkZ8LaD2DVq+acOW6e5srkNzzmsiuT/5eLGdnc9/lmVu8/i5e7G+/1b0yX+sXvyparU9gRESkOjq03OyDH7jLbldrArW9BmVrW1mWhxLRM7pmxiQ1H4ijh6c7HdzXlhhplrC5LCoHCjoiIK0uNg+XjYfMMs10iyFzmoVH/YjVnzj/FpWRw1/T17DyZiJ+PB5/e3ZxmlYOsLksKicKOiIgrMgzY/pU5OWDqOXNb40Hm5IC+xftL/XRCGoOmrSc6Npmgkl7MGtaC+uUDrC5LCpHCjoiIqzkXDT+MhsOrzXZILbj1bah8vbV1OYBj51MZOG0dx+MuEubvw+f3RFI9tJTVZUkhU9gREXEVmWnw29vw21uQnQEePtD2CWg1Cjy8rK7OcvvPJDHok/XEJqVTKdiXz4dHEhHka3VZUgQUdkREXMGhVWYH5LiDZrt6B+j2BgRVsbYuB7HjRAJ3TV/PhdRMapX147PhLQj197G6LCkili7dunr1arp37054eDg2m4358+f/a589e/bQo0cPAgICKFmyJM2bN+fYsWM5z6elpfHQQw8RHBxMqVKl6Nu3L2fOnCnCTyEiYqHks/DdvTCrhxl0SoXB7TNg4DcKOn9Yf+g8/aeu40JqJo0qBDD33pYKOsWMpWEnJSWFRo0a8cEHH1zy+YMHD9KmTRtq167NypUr2b59O88++yw+Pn/9ko4ePZqFCxfy9ddfs2rVKk6dOkWfPn2K6iOIiFjDbodNn8L7TWH7l4ANWtwLIzdAvd7FeqTV3/2yL5a7pm8gOT2LyCpBzB7RktIldUuvuLEZhmFYXQSAzWZj3rx59OrVK2dbv3798PT05LPPPrvkaxISEihTpgxz5szhtttuA2Dv3r3UqVOHtWvX0rJlyzy9d2JiIgEBASQkJODv73/Nn0VEpFCd3QcLH4Fja812uUZw6ztQvomlZTmahdtOMearKDKzDW6uHcqHA5vg41m81vtydXn9/rb0ys6V2O12fvjhB2rWrEnnzp0JDQ0lMjIy162uzZs3k5mZSYcOf01xXrt2bSpWrMjatWsve+z09HQSExNzPUREHF5WBqx8Faa0MYOOZ0noMgnuWaGg8w+frT3Cw3O3kpltcGvDckwZ1FRBpxhz2LATGxtLcnIykyZNokuXLvz888/07t2bPn36sGrVKgBOnz6Nl5cXgYGBuV5btmxZTp8+fdljT5w4kYCAgJxHREREYX4UEZFrd3wD/O9GWPmKOdKqRid4aD20fADcNdbkT4Zh8PbS/Tz7/S4MAwa3rMS7/Rrj5eGwX3dSBBz2X4jdbgegZ8+ejB49GoDrrruONWvWMGXKFNq2bXvVxx47dixjxozJaScmJirwiIhjSkuE5RNg4yeAAb4h0PVVc4Vy9cvJxW43eGHhLmatPQrAox1q8Ej7Gth0noo9hw07ISEheHh4ULdu3Vzb69Spw2+//QZAWFgYGRkZxMfH57q6c+bMGcLCLr+Qm7e3N97e3oVSt4hIgdm3GH4YA4knzfZ1A6HTS8V+BuRLyciyM+arKBZtj8Fmg/E96nFXq8pWlyUOwmGv63l5edG8eXP27duXa/v+/fupVKkSAE2bNsXT05Ply5fnPL9v3z6OHTtGq1atirReEZECkxwLX98NX9xpBp3SlWHwfOj1oYLOJaSkZzF85kYWbY/B093Gu/0aK+hILpZe2UlOTiY6OjqnffjwYaKioggKCqJixYo8/vjj3Hnnndx44420a9eOxYsXs3DhQlauXAlAQEAAw4cPZ8yYMQQFBeHv78+oUaNo1apVnkdiiYg4DMOArZ/Dz89AWjzY3KHVQ3DTWPDSTL+XciElg6EzNhJ1PJ4Snu78b3BTbqyplcslN0uHnq9cuZJ27dr9a/uQIUOYMWMGANOnT2fixImcOHGCWrVqMX78eHr27Jmzb1paGo899hhffPEF6enpdO7cmQ8//PCKt7H+SUPPRcRy5w/Cokf/Ws8qrCH0mAzh11lZlUM7FX+Ru6ZvIDo2mUBfTz69uzmNK5a2uiwpQnn9/naYeXaspLAjIpbJzoQ1k2HVq5CVBh4loN04aPmgRlldQXRsMndNW8+phDTKBfgwa1gLapT1s7osKWJ5/f7WvyQREauc3AILHoYzO8x21ZvM1cmDqlpalqPbfiKeuz/dSFxKBlXLlOSz4ZGUDyxhdVniwBR2RESKWkYK/PIKrPsQDDv4BEKXidCov4aT/4ffDpzjvs82kZKRTcMKAXx6d3OCS2l0rVyZwo6ISFGKXm72zYn/Y0Hj+reZsyCXUqfa//LjjhgenRtFRrad66sH87/BzSjlra8x+W/6LRERKQop52HJONg+12z7VzBvWdXsZG1dTmL2+qM8M38nhgHdGoTx9p3X4e2h5R8kbxR2REQKk2HAzm/hpycg9Txgg8j74OZnwFsdav+LYRi8vyKaN5fuB2BgZEUm9KyPu5tu90neKeyIiBSWxBhzBuR9P5rt0LrQ/T2IaG5tXU4i227w3Pc7mb3evOX38M3VGd2xppZ/kHxT2BERKWiGAVFzYMlYSEsAN0+48f+gzRjw8LK6OqeQlpnNw19s5efdZ7DZ4IXu9RjSurLVZYmTUtgRESlI8cfNDsjRy8x2eGPo+QGUrWdpWc4kPjWD4TM3sfnoBbw83Hj3zuvo2qCc1WWJE1PYEREpCIYBmz+Fn5+DjCRw94Z2Y6HVKE0OmA8nLqQyZPoGDp5Nwd/Hg0+GNKdFFa0HJtdG/wJFRK5V3GFY+PBfSz1UaGFezSlT09q6nMzuU4nc/ekGYpPSKRfgw8xhLaipWZGlACjsiIhcLbsdNk6FZS9AZqq51EP758zRVm4aFp0faw6e475Zm0lKz6JWWT9mDGtOuQDNiiwFQ2FHRORqnIuGBSPh2FqzXakN9HgPgqtZW5cTWrjtFI99tY2MbDstqgQx9a5mBJTwtLoscSEKOyIi+WHPhrUfwC8vmwt3epWCjuOh6TBwc7O6Oqfzya+HeOmHPQDc0qAcb97RCB9PXRWTgqWwIyKSV7F74fsH4eRms121nXk1J7CitXU5IbvdYOJPe5j662EA7m5dmedurYubJguUQqCwIyLyX7Iz4fd3YdWrkJ0B3gHQ+WVoPEgLd16FjCw7//f1NhZsOwXA2K61uffGqposUAqNwo6IyJWc3gHzH4TT2812zS7mmlb+4dbW5aSS0jK5//PN/B59Hg83G6/f3pDejStYXZa4OIUdEZFLycqAX9+AX98EexaUKA1dX4MGt+tqzlWKTUxjyKcb2ROTSEkvd6YMbsoNNbTauxQ+hR0RkX+K2QbzHoDYXWa7Tnfo9ib4lbW2LicWHZvEkOkbORl/kZBS3swY2pz65QOsLkuKCYUdEZE/ZWfCr2/B6tfMqzm+IXDLG1Cvt9WVObU1B89x/2ebSUzLokpISWYObUHFYF+ry5JiRGFHRAQgdg/Mux9iosx2nR5m35ySIZaW5ey+23KCJ7/dTma2QbNKpfn4rmYEldRiqFK0FHZEpHizZ8Oayea8OdkZ4BMIt7wJ9fuqb841MAyDySuieWvpfgBuaViON2/XHDpiDYUdESm+zh80r+ac2GC2a3Q2583xC7O2LieXmW1n3Hc7+HrzCQDua1uVJzvX1hw6YhmFHREpfv5c02rp85B1Ebz8oOskuG6gruZco8S0TB78fAu/RZ/DzQYTetZnUMtKVpclxZzCjogULxeOwvcPwZFfzXaVtuYK5YER1tblAk7GX2TopxvYfyYZXy93PhjQhHa1Q60uS0RhR0SKCcOALbNgyTjISAZPX+g4AZoN15pWBWDnyQSGzdhIbFI6oX7eTL9bQ8vFcSjsiIjrSzwFCx6G6KVmu2Ir82qOVigvECv2nmHknK2kZmRTq6wf04c2p3xgCavLEsmhsCMirsswYPtX8NPjkJYA7t7Q/llo+SC4aVRQQfhs3VGe/34ndgPaVA/hw0FN8PfxtLoskVwUdkTENSWfhUWPwt5FZju8MfT+H5SpZWlZrsJuN5i0eC8frz4EwG1NKzCxTwM83XVLUByPwo6IuJ5d8+GHMZB6Htw8oe2T0GY0uOtPXkFIy8zmsa+28cOOGADGdKzJqJura9VycVj6ly8iriM1Dn58HHZ+Y7bL1odeH0G5htbW5ULiUjIYMWsTm49ewNPdxmu3adVycXwKOyLiGg4she9HQvJpsLmZV3LaPgke3lZX5jIOnU1m2IyNHDmfir+PB/8b3IxW1YKtLkvkPynsiIhzy0iFpc/Cxk/MdnAN6D0FKjSzti4Xsyb6HPd/bi7mWT6wBDOGNqdGWT+ryxLJE4UdEXFep7bCtyPg/AGz3eI+6DgePDXsuSB9seEYz87fSZbdoHHFQD4e3IwyfrpiJs5DYUdEnI89G357C1ZOAnsWlAqDXh9A9Q5WV+ZSsu0GE3/cwye/HQagR6NwXrutoRbzFKejsCMiziXuMMy7D46vN9t1ekD3d8E3yNq6XExyehaPfLGV5XtjARjdoSYPt9eIK3FOCjsi4hwMA6Jmw09Pmss9ePlBt9egUX8t3lnATsZfZPiMjew9nYSXhxtv3N6IHo3CrS5L5Kop7IiI40s5Dwsf/muCwIqtzE7IpStbWpYr2nrsAiNmbeZccjohpbz5+K6mNKlY2uqyRK6Jwo6IOLYDy+D7ByH5jDlBYLtxcP0jWu6hECzaforHvtpGepad2mF+TLtba1yJa1DYERHHlJEKS5+DjVPNdkgt6DsVyjWyti4XZBgGk1dE89bS/QC0rx3Ku/0bU8pbXxHiGixdxGT16tV0796d8PBwbDYb8+fPz/X83Xffjc1my/Xo0qVLrn0qV678r30mTZpUhJ9CRArcqa3wcdu/gk6L++C+VQo6hSAtM5tHv4zKCTr3tKnCx3c1U9ARl2Lpb3NKSgqNGjVi2LBh9OnT55L7dOnShU8//TSn7e3977kdJkyYwIgRI3Lafn6a6ErEKdmz4be3YeVEDSkvAueS07l31ia2HIvHw83GhJ71GRBZ0eqyRAqcpWGna9eudO3a9Yr7eHt7ExYWdsV9/Pz8/nOfv0tPTyc9PT2nnZiYmOfXikghuXAEvrsPjq8z23V7wq3vaEh5Idl3OolhMzZyMv4i/j4efDSoKddXD7G6LJFCYeltrLxYuXIloaGh1KpViwceeIDz58//a59JkyYRHBxM48aNef3118nKyrriMSdOnEhAQEDOIyIiorDKF5H/Yhiw9XP46Hoz6Hj5Qa8pcPtMBZ1C8su+WPp+tIaT8RepHOzLvIeuV9ARl2YzDMOwuggAm83GvHnz6NWrV862uXPn4uvrS5UqVTh48CDjxo2jVKlSrF27Fnd3cyTGW2+9RZMmTQgKCmLNmjWMHTuWoUOH8tZbb132vS51ZSciIoKEhAT8/f0L7TOKyD+kxplDyvcsNNsVW0Hv/0HpStbW5aIMw2DGmiO8uGg3dgMiqwQxZVBTSpf0sro0kauSmJhIQEDAf35/O3TY+adDhw5RrVo1li1bRvv27S+5z/Tp07nvvvtITk6+ZP+eS8nryRKRAnT4V/juXkg6pSHlRSA9K5vn5u/iy03HAbijWQVe6tUALw+Hv8Avcll5/f52qt/yqlWrEhISQnR09GX3iYyMJCsriyNHjhRdYSKSd9lZsOIlmNndDDrB1eGeZXDDGAWdQnI2KZ0BU9fz5abjuNng6W51eLVvQwUdKTacamzhiRMnOH/+POXKlbvsPlFRUbi5uREaGlqElYlInlw4Ct/eAyc2mO3Gg6DLq+Bdytq6XNiOEwnc+9kmYhLS8PPx4P0BTWhbs4zVZYkUqQIPO1lZWXh45O2wycnJua7SHD58mKioKIKCgggKCmL8+PH07duXsLAwDh48yBNPPEH16tXp3LkzAGvXrmX9+vW0a9cOPz8/1q5dy+jRoxk0aBClS2t6cxGHsvNbWPgopCeCtz90fwfq97W6Kpe2YNspHv/anBG5apmSfHJXM6qWUbCU4idf1zC/+uqrKz6flZXFHXfckefjbdq0icaNG9O4cWMAxowZQ+PGjXnuuedwd3dn+/bt9OjRg5o1azJ8+HCaNm3Kr7/+mtMXx9vbm7lz59K2bVvq1avHyy+/zOjRo/n444/z87FEpDClJ8P8h+CbYWbQqdAC7v9NQacQZdsNXl28l4e/2Ep6lp12tcow/6HrFXSk2MpXB2UfHx8WLlxIx44d//VcdnY2t99+O2vXriUmJqZAiyxs6qAsUkhORcG3w+F8NGCDG/8P2j4F7k51B92pJKVl8sjcKFbsjQXg/rbVeLxzLdzdtDK8uJ68fn/n6y/Oq6++Sp8+fVi2bBmRkZE52+12O3fccQe///47K1asuPqqRcQ12O2w7kNY9gLYM8EvHPp8DFVusLoyl3b4XAojZm0iOjYZbw83Xu3bkF6Ny1tdlojl8hV2HnnkEeLi4ujWrRurV6+mXr16ZGdnc+edd/Lrr7+yYsUK6tWrV1i1iogzSI6F+Q9A9DKzXftW6DFZEwQWstX7zzJyzhYS07II8/fh47ua0rBCoNVliTiEfF9LHj9+PHFxcXTq1IlffvmFZ555hlWrVrF8+XLq169fGDWKiLOIXgbzHoCUWPDwgc6vQLNhYNMtlMJiGAbTfjvMKz/uwW5Ak4qBTBnUlFB/H6tLE3EYV3XjfPLkyVy4cIFGjRpRqlQpli9fTsOGDQu6NhFxFlkZsHw8rH3fbIfWhb7ToGxda+tycelZ2Tw9byffbD4BwO1NK/BS7/p4e2i+IpG/y1fYGTNmTM7PpUuXxjAMrrvuOmbMmJFrvyst1SAiLuZcNHw7DGK2me3mI6DTi+BZwtq6XFxsYhr3fb6ZrcficbPBM7fUZej1lbHpKprIv+Qr7GzdujVXu1WrVmRlZeXarn9oIsWEYUDUHPjxcchMgRKloecHUPsWqytzeduOx3PvZ5s4k5hOQAlPPhjQhDY1tJCnyOXkK+z88ssvhVWHiDiTtARYNNqcKBCg8g3maCv/cGvrKgbmbz3JE99uJyPLTo3QUky9qxmVQ0paXZaIQ9NkFyKSPye3wNd3Q/xRsLmbC3i2Ga11rQpZVradVxfvZeqvhwHoUCeUt++8Dj8fT4srE3F8CjsikjeGARs+hiVPm3PnBFaEvtMhornVlbm8c8npjJqzlbWHzgPwULtqPNaxFm6aKFAkTxR2ROS/XYyHBSNhz0KzXftWs39OiUArqyoWoo7H88Dnm4lJSKOklztv3N6Irg0uvxiyiPybwo6IXNnfb1u5eUKnlyDyPs2dUwTmbjjGc9/vIiPbXMjz48FNqR7qZ3VZIk4nXwuBTp8+nXPnzhVWLSLiSAwD1v8PpnUyg05gRRi+BFrer6BTyNKzshn73Xae+m4HGdl2OtUty/cPXa+gI3KV8hV2Pv/8cypUqEDr1q159dVX2bNnT2HVJSJWSkuAr+6Cn54w++fUvhXu+xXKN7W6MpcXk3CRO/63ji82HMdmg8c712LKoKbqiCxyDfIVdlasWEFMTAwPPvggmzdvJjIykho1avDYY4+xevVq7HZ7YdUpIkXl1Fb4342wZ4F526rLJLjzc/XPKQJrD57n1vd+Y9vxeAJKeDJjaAsealddHZFFrpHNMAzjal+ckZHBihUrWLBgAQsXLuTixYt069aNHj160LVrV0qWdI65H/K6RLyISzMM2PgJLBkH2RnmbavbZkAFXc0pbH+ubzXxp71k2w3qlvPnf4ObEhHka3VpIg4tr9/f1xR2/mnTpk0sWLCA77//nttuu41nn322oA5dqBR2pNhLS4AFo2D392a79q3Q831zVmQpVKkZWTz57Q4WbjsFQJ/G5Xm5dwNKeGneIpH/YknY+bvMzEw8PZ3jHrPCjhRrp6Lg6yFw4cgfo61ehEh1Qi4KR86lcN9nm9l3JgkPNxvP3lqXu1pV0rI7InmU1+/vQht67ixBR6TY+udtq4CKcPsM3bYqIiv2nuGRuVEkpWVRxs+bDwc2oXnlIKvLEnFJmmdHpDhKS4AFD8Pu+Wa71i3Q6wPdtioCdrvBeysO8M6yAwA0rVSaDwc2oay/j8WVibguhR2R4uZUlDlJ4IXD4OYBHSdAywd126oIJFzMZMyXUSzfGwvA4JaVePbWunh55GtgrIjkU77+hU2YMIHU1NTCqkVECpNhwIapMK2jGXQCImDYEmj1kIJOEdh7OpEe7//G8r2xeHm48fptDXmxV30FHZEikK8Oyu7u7sTExBAaGlqYNRU5dVAWl5eeDAsfhp3fmu1a3cy1rXzVR6QozNt6gnHf7eRiZjblA0vwv8FNqV8+wOqyRJxeoXRQLqSBWyJSmM4dgC8Hwdm95m2rDi9Aq5G6mlME0jKzGb9wN19sOAZAm+ohvNe/MUElvSyuTKR4yXefHQ2JFHEiexbCvAcgIwlKhcEdM6FiS6urKhaOnk/hwdlb2HUqEZsNRt1cg0fa18BdsyGLFLl8h52aNWv+Z+CJi4u76oJEpABkZ8GKCfD7u2a70vVw26fgV9bauoqJxTtP8/g320hKyyKopBfv3HkdN9YsY3VZIsVWvsPO+PHjCQjQvWYRh5V8Fr4ZCkd+NdutRpq3rtw191Vhy8y2M+mnvUz77TBgDit/f0BjygWUsLgykeIt32GnX79+LtdBWcRlHN9orlaedAo8S5pLPtTvY3VVxcKp+IuMnLOFLcfiARhxQxWe6FIbT3eNthKxWr7CjvrriDioP2dDXjwW7JkQUtNcqbxMLasrKxZW7T/Lo3O3ciE1Ez8fD964vRGd64VZXZaI/EGjsUScXUYqLBoN2+ea7bo9zWHl3n7W1lUMZNsN3l22n8m/RGMYUC/cn48GNqVisFYrF3Ek+Qo7dru9sOoQkasRdwi+HAxndoLNHTqO17DyInI2KZ1Hv9zK79HnARgYWZFnb62Lj6dWKxdxNFouQsRZ7fsJvrsP0hOgZBlzEc/KbayuqlhYf+g8o77YSmxSOr5e7kzs04Ce15W3uiwRuQyFHRFnY8+GX16BX98w2xGRZtDxD7e0rOLAbjf43+pDvPHzPrLtBjVCS/HRoCZUD9UtQxFHprAj4kxSzsO3w+HQL2a7xX3Q6SXw0Iy8hS0+NYPHvtqWs4hn78blebl3fXy99GdUxNHpX6mIszi5xRxWnnAcPH2h+3vQ8HarqyoWth2P58HZWzgZfxEvDzfG96hHv+YRGqEq4iQUdkScweaZ8OP/QXYGBFU1h5WXrWd1VS7PMAxmrT3KSz/sJjPboFKwLx8MaKJFPEWcjMKOiCPLSjdDzpZZZrvWLdD7I/DRl21hS7iYybjvdvDDjhgAOtcry+u3N8LfRzNRizgbhR0RR5V02hxWfmID2Nzg5mfh+kfBTTPyFrYtxy7w8BdbOXHhIh5uNp7qWpvhbarotpWIk1LYEXFEJzfD3EHmsg8+AXDbdKjeweqqXN6fo63e/HkfWXaDikG+vNe/MddFBFpdmohcA0v/F3H16tV0796d8PBwbDYb8+fPz/X83Xffjc1my/Xo0qVLrn3i4uIYOHAg/v7+BAYGMnz4cJKTk4vwU4gUsKgvYHpXM+iE1IIRvyjoFIHYpDSGfLqBVxfvJctucGvDcix6uI2CjogLsPTKTkpKCo0aNWLYsGH06XPpxQq7dOnCp59+mtP29vbO9fzAgQOJiYlh6dKlZGZmMnToUO69917mzJlTqLWLFLjsLFj2PKx932zX6ga9/wc+/tbWVQys3n+WMV9FcS45Ax9Pc7TVHc002krEVVgadrp27UrXrl2vuI+3tzdhYZdeUG/Pnj0sXryYjRs30qxZMwAmT55Mt27deOONNwgP1yRr4iRS4+CbYX/Nn3PjE3DTWPXPKWSZ2Xbe+Hkf/1t1CIDaYX68P6CxJgkUcTEO/5d05cqVhIaGUqtWLR544AHOnz+f89zatWsJDAzMCToAHTp0wM3NjfXr11/2mOnp6SQmJuZ6iFgmdg9MbWcGHU9fuH0m3Py0gk4hO3Y+ldumrM0JOoNbVmL+Q9cr6Ii4IIfuoNylSxf69OlDlSpVOHjwIOPGjaNr166sXbsWd3d3Tp8+TWhoaK7XeHh4EBQUxOnTpy973IkTJzJ+/PjCLl/kv+1ZBPPug4xkCKwI/b6AsPpWV+XyFm47xbjvdpCUnoW/jwev3daILvUvfQVZRJyfQ4edfv365fzcoEEDGjZsSLVq1Vi5ciXt27e/6uOOHTuWMWPG5LQTExOJiIi4plpF8sVuh9Wvw8pXzHblG8wrOiWDra3LxV3MyGb8wl3M3XgcgGaVSvNOv+uoUNrX4spEpDA5dNj5p6pVqxISEkJ0dDTt27cnLCyM2NjYXPtkZWURFxd32X4+YPYD+mdHZ5Eik54M8++HPQvNduT95vpW7pqsrjDtPZ3IyDlbiY5NxmaDke2q80j7Gni463ahiKtzqrBz4sQJzp8/T7ly5QBo1aoV8fHxbN68maZNmwKwYsUK7HY7kZGRVpYqcmlxh2HuAIjdDe5ecOvb0HiQ1VW5NMMw+Hz9MV5ctJuMLDuhft680+86WlcLsbo0ESkiload5ORkoqOjc9qHDx8mKiqKoKAggoKCGD9+PH379iUsLIyDBw/yxBNPUL16dTp37gxAnTp16NKlCyNGjGDKlClkZmYycuRI+vXrp5FY4ngOrYSv74aLF6BUWbhzNkQ0t7oql5aQmsmT325n8S6zD1+7WmV44/ZGBJfSlV2R4sRmGIZh1ZuvXLmSdu3a/Wv7kCFD+Oijj+jVqxdbt24lPj6e8PBwOnXqxIsvvkjZsmVz9o2Li2PkyJEsXLgQNzc3+vbty3vvvUepUqXyXEdiYiIBAQEkJCTg7685TaSAGQas+wh+fgaMbCjf1FzI01+BvDBtOhLHI3OjOBl/EU93G0920ZIPIq4mr9/floYdR6GwI4UmMw0WjYZtf0xy2WiAeevK08faulxYtt3go5XRvL3sANl2g8rBvkzu34QGFbR4qoiryev3t1P12RFxKklnzP45JzeZC3l2ehlaPgC6slBoTsVf5LGvtrH2kDkfV+/G5XmxV31KeetPnUhxpr8AIoXh9E6YcyckngCfQLh9BlT79y1bKTgLt53i6Xk7SEzLwtfLnRd71qdv0wpWlyUiDkBhR6Sg7V9iLv2QkQzB1WHAVxBczeqqXFZiWiYvfL+L77aeBKBRRCDv3HkdVUJKWlyZiDgKhR2RgpLTEflpMOxQ5Ua4YxaUKG11ZS5rw+E4Rn9pdkJ2s8HIm2sw6ubqeGruHBH5G4UdkYKQnQk/Pg6bPzXbTYbALW9qosBCkpFl551l+/lo1UEMAyoG+fL2ndfRtJKCpYj8m8KOyLW6GA9fDzHn0cFmzobc6iF1RC4k0bHJPPrlVnaeNBfwvb1pBZ7vUU+dkEXksvTXQeRaxB0yOyKf2w+eJeG2aVCrq9VVuaQ/Z0J++YfdpGXaCfT1ZGLvBnRtUM7q0kTEwSnsiFyto2tg7kC4GAf+5aH/XCjX0OqqXNLZpHSe+GYbv+w7C8ANNUJ44/ZGlPXXfEUi8t8UdkSuRtQXsGAU2DMhvLEZdPwuv/isXL2lu8/w1LfbOZ+SgZeHG2O71mZIq8q4uek2oYjkjcKOSH7Y7fDLS/Drm2a7bk/oNQW8fK2tywWlZmTx4qI9fLHhGAC1w/x4t19jaoX5WVyZiDgbhR2RvMpIhfn3w+7vzfYN/wftngY3DXMuaFHH4xn9ZRSHz6Vgs8GIG6ryWKeaeHu4W12aiDghhR2RvEg6DV/0g1Nbwc0TekyG6/pbXZXLycq289HKg7yz3FzXqlyAD2/e3ojW1UOsLk1EnJjCjsh/Ob3jj6UfTkKJIOg3Gyq1troql3PsfCqjv4pi89ELANzSsByv9GpAgK/mKhKRa6OwI3Il+36Cb4ZDZgqE1IQBX0JQVaurcimGYfDN5hO8sGAXKRnZ+Hl7MKFXPXpdVx6b5ioSkQKgsCNyOWs/hCXjAAOq3gS3z4QSgRYX5VrOJ6fzzPyd/LTzNAAtKgfx5h2NiAhSh28RKTgKOyL/ZLeb61ut+9BsNx0K3V7X0g8FbPHOGJ6et5PzKRl4uNkY06km991YDXcNKReRAqawI/J3mRfhu3thzwKz3fFFaD1KSz8UoPjUDJ5fsIvvo04BUKusH2/e0Yj65QMsrkxEXJXCjsifUuPgi/5wfB24e0Gvj6DBbVZX5VJW7D3DU9/uIDYpHTcb3N+2Go90qKEh5SJSqBR2RAAuHIHPb4PzB8AnAPrNgcptrK7KZSSmZfLiwt18vfkEAFXLlOTN2xvRuKJWKReRwqewI3JqK8y+A1Jiwb8CDPoGQutYXZXLWL3/LE9+u52YhDRsNrinTRUe61QLH09dzRGRoqGwI8XbgaXw1RBzaHnZBjDwa/DXKtoFITk9i1d+3MOc9eZyD5WCfXnj9kY0rxxkcWUiUtwo7EjxtWUWLHwUjGyo2g7umAU+/lZX5RLWHjzP499s48SFiwAMaVWJJ7vWxtdLf3JEpOjpL48UP4YBKyfCqlfNdqMB0OM9DS0vABczsnl18V5mrDkCQPnAErx+W0Mt9yAillLYkeIlO9O8mhP1udm+8QloN05DywvA5qNx/N/X2zl8LgWA/i0qMq5bbfx8FCJFxFoKO1J8pCfBV3fBwRVgc4db34amQ6yuyumlZWbz9tL9TP31EHYDwvx9ePW2hrStWcbq0kREAIUdKS6STsPs28xFPT19zaUfanayuiqnt/1EPI99tY0DsckA9G1Sgee61yWghK7miIjjUNgR13d2H3zeFxKOQ8kyMOArKN/E6qqcWkaWnckrDvDhyoNk2w1CSnkzsU8DOtYta3VpIiL/orAjru3oGviiH6QlQHB1GPgNBFWxuiqntvtUIo99vY09MYkAdG8UzoQe9Shd0sviykRELk1hR1zXrnnmOlfZGRARCf3ngq/meLlamdl2Plp5kMkrDpCZbRBU0ouXetWnWwPNSyQijk1hR1zThqnw4+OAAXW6Q5+p4FnC6qqc1o4TCTz+zTb2nk4CoHO9srzcuwEhpbwtrkxE5L8p7IhrMQxz/pyVE81283ug62vgpqUJrkZaZjZvL9vP1NXmSKugkl48370uPRqFY9NwfRFxEgo74jrsdlj8JGz42GzfNBbaPqk5dK7S+kPneeq7HTnz5vRoFM7z3esSrKs5IuJkFHbENWRlwPwHYOc3gA26vQ4tRlhdlVNKSsvk1cV7+XyduaZVmL8PL/WqTweNtBIRJ6WwI84vIwW+HAwHl4ObB/T+HzS4zeqqnNIve2N5et4OTiWkATAgsiJPda2Nv2ZBFhEnprAjzi01DubcASc2mpMF3vkZVO9gdVVOJy4lgxcX7Wbe1pOAuUL5xD4NaF1Na1qJiPNT2BHnlXASPu8DZ/dCidIw4GuIaG51VU7FMAwWbY/hhQW7OJ+SgZsN7rmhKqM71KSElzp1i4hrUNgR53TuAHzW25wV2S8cBs+D0NpWV+VUziSm8fS8nSzbcwaAWmX9ePW2hlwXEWhtYSIiBUxhR5zPyS3mOlep581ZkQfPg8CKVlflNAzD4MuNx3n5xz0kpWXh6W5jZLsaPHBTNbw83KwuT0SkwCnsiHM5tBLmDoSMZCh3HQz6FkqqX0leHTufylPfbWfNwfMANIoI5LW+DakV5mdxZSIihcfS/41bvXo13bt3JzzcnKBs/vz5l933/vvvx2az8c477+TaXrlyZWw2W67HpEmTCrdwscau+TD7djPoVLkR7l6koJNH2XaDT349RKd3VrHm4Hl8PN145pY6fPdAawUdEXF5ll7ZSUlJoVGjRgwbNow+ffpcdr958+axbt06wsPDL/n8hAkTGDHirzlV/Pz0x9vlbJ4JCx/BXP6hB/T9BDw0uV1e7IlJZOx3O4g6Hg9Aq6rBTOrbgErBJa0tTESkiFgadrp27UrXrl2vuM/JkycZNWoUS5Ys4ZZbbrnkPn5+foSFheX5fdPT00lPT89pJyYm5vm1YoG1H8KSsebPTe+GW97S8g95cDEjm3eW7+eTXw+TbTfw8/bg6VvqcGfzCC31ICLFikP3RrTb7QwePJjHH3+cevXqXXa/SZMmERwcTOPGjXn99dfJysq64nEnTpxIQEBAziMiIqKgS5eCYBiw6vW/gs71j8Ct7yjo5MHKfbF0fHsV/1t1iGy7Qdf6YSwd05Z+LSoq6IhIsePQHZRfffVVPDw8ePjhhy+7z8MPP0yTJk0ICgpizZo1jB07lpiYGN56663Lvmbs2LGMGTMmp52YmKjA42gMA5Y9D7+/a7bbPQM3/p/WufoPsUlpTFi4m0XbYwAoH1iC8T3qaakHESnWHDbsbN68mXfffZctW7Zc8f9E/x5aGjZsiJeXF/fddx8TJ07E2/vSfTq8vb0v+5w4ALsdfnocNn5itju/Aq0esrYmB2e3G8zdeJxJP+0hMS0LNxsMu74KozvWpKS3w/4zFxEpEg77V/DXX38lNjaWihX/mj8lOzubxx57jHfeeYcjR45c8nWRkZFkZWVx5MgRatWqVUTVSoHJzoIFo2DbHMAG3d8x++nIZe07ncS4eTvYfPQCAA3KBzCxTwPqlw+wuDIREcfgsGFn8ODBdOiQe42jzp07M3jwYIYOHXrZ10VFReHm5kZoaGhhlygFLSsDvrsHdn8PNndzQc+Gt1tdlcNKy8xm8ooD/G/VIbLsBiW93HmsUy2GtK6Mu5tu94mI/MnSsJOcnEx0dHRO+/Dhw0RFRREUFETFihUJDg7Otb+npydhYWE5V2zWrl3L+vXradeuHX5+fqxdu5bRo0czaNAgSpcuXaSfRa5R5kX46i448DO4e8Ftn0KdW62uymH9duAcT8/fwdHzqQB0rFuW8T3qER5YwuLKREQcj6VhZ9OmTbRr1y6n/Wf/myFDhjBjxoz/fL23tzdz587lhRdeID09nSpVqjB69Ohc/XjECaQnwxf94Miv4FEC+s2G6u2trsohnUtO5+Uf9uSsTh7m78P4nvXoXC/vUy+IiBQ3NsMwDKuLsFpiYiIBAQEkJCTg7+9vdTnFy8V4c1bkExvAyw8GfgWVWltdlcMxDIOvNh3nlR/3knAxE5sNhrSqzGOdauLn42l1eSIilsjr97fD9tmRYiDlHHzWC07vgBKlzXWuyje1uiqHEx2bzLh5O9hwOA6AuuX8mdinAY20OrmISJ4o7Ig1EmNgVk84tw9KhsJd86Hs5SeOLI7SMrP5cOVBPloZTWa2QQlPd8Z0rMnQ6yvj4e7Q84GKiDgUhR0pegknYWZ3iDsI/hXgru8hpLrVVTmUX/bG8vyCXRyLMzsg31w7lAk961GhtK/FlYmIOB+FHSla8cdh5q1w4QgEVoQhi6B0Jaurchgn4y8yYeEuluw6A5gdkJ/rXpeu9cO0zIOIyFVS2JGic+GoGXTij0HpymbQCdQyHQAZWXam/XaY95Yf4GJmNu5uNoa3qcLD7WtQSjMgi4hcE/0VlaIRdwhm9oCE4xBUDYYshIDyVlflENYcPMdz3+8iOjYZgBaVg3ixV31qhflZXJmIiGtQ2JHCd/6g2Ucn8SQE1zCDjn85q6uyXGxSGq/8sIf5UacACC7pxbhudejTpLxuWYmIFCCFHSlc5w7AjFsh+TSE1DKDjl/xXoE7K9vO5+uO8ubP+0lKz8Jmg0GRlfi/TrUI8NWcOSIiBU1hRwpP7F7zik5KLITWhbsWQKkyVldlqS3HLvDMvJ3sjkkEoGGFAF7qVZ+GFQKtLUxExIUp7EjhOLMbZvWAlLNQtoE5vLxk8H+/zkVdSMng1cV7mbvxOAD+Ph480aU2/VtU1KKdIiKFTGFHCt6ZXeYVndTzENbQDDq+QVZXZYlsu8GXG4/z+pK9XEjNBOC2phV4qmttQkp5W1ydiEjxoLAjBevM7r+CTnhjGDzPXAqiGNp89AIvLNjFjpMJANQq68dLvevTvHLxDH4iIlZR2JGC82cfndTzUO46GDwfSgRaXFTRi01K49Wf9vHtlhMA+Hl78GjHmtzVqhKeWuZBRKTIKexIwTi7/4+gc868dTV4XrELOpnZdmauOcK7yw6QlJ4FwO1NK/BEl9qU8dMtKxERqyjsyLU7d8CcGTkl9q/OyMWsj87v0ed4YcEuDvwxMWDDCgG80KMeTSoWz1t4IiKORGFHrs256D/m0TkDZesXu6BzMv4iL/+wmx93nAYgqKQXT3SuxR3NInDTKCsREYegsCNX7/xB84pO8uk/5tEpPsPL0zKzmbr6EB+sjCYt046bDe5qVZnRHWpqYkAREQejsCNXJ+6Q2UcnKQbK1DYnDCwZYnVVhc4wDJbtieXFRbs5FpcKQIsqQYzvUY865fwtrk5ERC5FYUfy78IRmPHHWld/LgFRDGZGPnQ2mQmLdrNy31kAyvp7M65bHXo0CtdaViIiDkxhR/In/vgfQefEX4t6lgq1uqpClZSWyfsropn++2Eysw083W3cc0NVRrarTklv/RMSEXF0+ksteZcYY966SjgGQdVcflHPbLvBN5uP8/qSfZxLzgDgplpleO7WulQtU8ri6kREJK8UdiRvks/CrJ5w4TAEVoQhC8C/nNVVFZoNh+MYv3AXu06ZC3ZWDSnJs7fWpV1t176KJSLiihR25L+lxsFnveHcPvALN6/oBFSwuqpCcTL+IhN/3MOi7TEA+Pl48Ej7GtzVqjJeHpr9WETEGSnsyJWlJcLnfeHMDigZagad0pWtrqrApWZkMWXVIf636iDpWXZsNujfoiKPdaxJsBbsFBFxago7cnkZKTDnDji1BUoEmfPohFS3uqoCZRgGC7adYuKPezmdmAZAy6pBPHdrPeqGayi5iIgrUNiRS8u8CF/0g2NrwTvAXOuqbF2rqypQ247HM2HRbjYfvQBAhdIleLpbHbrUD9NQchERF6KwI/+WnQlfDYHDq8GrFAz6FsKvs7qqAhOTcJHXl+zjuy0nAfD1cuehdtUZ3qYKPp7uFlcnIiIFTWFHcrNnw3f3woEl4FECBnwJEc2trqpApKRn8b/Vh/h49UHSMu0A9Glcnie71qasv4/F1YmISGFR2JG/GAYsehR2fQdunnDn51C5jdVVXTO73eCbLSd4Y8k+YpPSAWheuTTP3lqXhhUCrS1OREQKncKOmAwDfn4GtswCmxv0/QRqdLC6qmu29uB5Xvphd858ORWDfBnbtbb65YiIFCMKO2Ja/Tqsfd/8ucdkqNfL0nKu1eFzKUz8cQ8/7z4DgJ+3B6PaV2dI68p4e6hfjohIcaKwI7DuI/jlZfPnLpOg8SBr67kGCamZvLv8ALPWHiHLbuDuZmNAi4o82qGG5ssRESmmFHaKu62fw+KnzJ9vGgctH7C2nquUmW3n83VHeXf5AeJTMwFoV6sM47rVoUZZP4urExERKynsFGe7F8CCUebPrUZC2yesrecqGIbBkl1neHXxXg6fSwGgZtlSPHNLXW6sWcbi6kRExBEo7BRXh1bCt8PBsEPjwdDpJXCyDrtbj13g5R/2sOmPSQGDS3oxplNN7mwWgYe71rESERGTwk5xdHIzzB0I2RlQpzt0f9epgs6x86m8umQvP/yxWKePpxsjbqjKfW2rUcpbv9IiIpKbvhmKm7P74fPbICMZqrSFvtPAzTlGJ8WnZjB5RTSz1h4hM9vAZoPbmlRgTKealAsoYXV5IiLioBR2ipP44/BZL7gYB+GNod9s8HD8EUrpWdnMWnOUySsOkJiWBcANNUIY27WOFusUEZH/ZGnHhtWrV9O9e3fCw8Ox2WzMnz//svvef//92Gw23nnnnVzb4+LiGDhwIP7+/gQGBjJ8+HCSk5MLt3BnlHIePusNiSchuAYM/Ba8HXuUkt1u8H3USTq8tYqXf9xDYloWtcP8mDmsBZ8Nj1TQERGRPLH0yk5KSgqNGjVi2LBh9OnT57L7zZs3j3Xr1hEeHv6v5wYOHEhMTAxLly4lMzOToUOHcu+99zJnzpzCLN25ZKTAnNvh/AHwrwB3zYeSwVZXdUW/R59j0k972XEyAYCy/t481rEWfZtWwN3NefoXiYiI9SwNO127dqVr165X3OfkyZOMGjWKJUuWcMstt+R6bs+ePSxevJiNGzfSrFkzACZPnky3bt144403LhmOip0/VzA/uRlKlIbB30FABauruqzdpxKZtHgvq/efBaCUtwf33ViV4TdUwddLd11FRCT/HPrbw263M3jwYB5//HHq1av3r+fXrl1LYGBgTtAB6NChA25ubqxfv57evXtf8rjp6emkp6fntBMTEwu+eEdgt8P3IyF66R8rmH8NZWpZXdUlnbiQyls/72de1EkMAzzdbQyMrMSom6tr5mMREbkmDh12Xn31VTw8PHj44Ycv+fzp06cJDQ3Ntc3Dw4OgoCBOnz592eNOnDiR8ePHF2itDmnZ87B9Ltjc4Y6ZENHc6or+JT41gw9+iWbm2qNkZNkBuLVhOR7vXItKwSUtrk5ERFyBw4adzZs38+6777Jly5YCX5167NixjBkzJqedmJhIREREgb6H5dZ+AGveM3/uMRlqdra2nn9Iy8xmxpojfPhLdM4Iq1ZVg3mqa20aRQRaW5yIiLgUhw07v/76K7GxsVSsWDFnW3Z2No899hjvvPMOR44cISwsjNjY2Fyvy8rKIi4ujrCwsMse29vbG29vF741suMbWDLO/LnDC9B4oKXl/F1Wtp3vtpzk7WX7iUlIA6B2mB9Pdq3NTTXLFHiwFRERcdiwM3jwYDp06JBrW+fOnRk8eDBDhw4FoFWrVsTHx7N582aaNm0KwIoVK7Db7URGRhZ5zQ7h0CqYd7/5c+QDcP2jlpbzpz/XsHrj531Ex5pTA4QH+DCmUy16Ny6vEVYiIlJoLA07ycnJREdH57QPHz5MVFQUQUFBVKxYkeDg3MOjPT09CQsLo1Yts5NtnTp16NKlCyNGjGDKlClkZmYycuRI+vXrVzxHYp3eAV8OAnsm1OsNnV9xiGUg1h06z6uL97L1WDwAgb6ePHRTdQa3qoSPp3PM3iwiIs7L0rCzadMm2rVrl9P+sx/NkCFDmDFjRp6OMXv2bEaOHEn79u1xc3Ojb9++vPfee4VRrmOLPw6zb4f0RKjUBnpNATdrF8PcfSqR15bsZeU+cxh5CU93hrepwr1tq+Lv42lpbSIiUnzYDMMwrC7CaomJiQQEBJCQkIC/vxPOynvxAkzvAmf3Qpk6MGwxlAi0rJyj51N4e+l+vt92CsMADzcb/VpE8PDNNQj197GsLhERcS15/f522D47kkdZ6eYK5mf3gl84DPrGsqBzJjGN95Yf4MuNx8mymxn61obl+L9OtagcomHkIiJiDYUdZ2a3w/wH4Ojv4O1vBh0LZkeOT83go1UHmbnmCGmZ5lw5bWuW4fHOtahfPqDI6xEREfk7hR1ntnw87PwW3Dzgzs+g7L9nmS5MKelZfPr7Yf63+hBJf8yV07RSaZ7oXIvIqo699paIiBQfCjvOauMn8Ps75s893oeqNxXZW6dnZfPF+mO8/0s055IzAHOunMc71+Lm2qGaK0dERByKwo4z2r8Efnzc/Lnd03Bd/yJ52z8nBHx3+QFOxl8EoFKwL2M61qR7w3DcNFeOiIg4IIUdZxOzDb4eCoYdGg+CGx8v9Le02w1+2BHD20v3c+hcCgBl/b15uH0N7mgWgae7tUPcRURErkRhx5kknIQ5d0Jminnb6tZ3CnXSQMMwWLE3ljd+3s+eGHNl+NK+njzUrjqDWmpCQBERcQ4KO84iPQnm3AFJMeZcOnfMAvfCmZjPMAzWHDzPGz/vy5n12M/bgxE3VmVYmyqU8tavjYiIOA99azmD7Czz1tWZnVCqLAz8GnwKZ0j3xiNxvPnzPtYdigPAx9ONu1tX4f62VQn09SqU9xQRESlMCjuOzjBg8VMQvRQ8SkD/uRAYUeBvs+14PG8u3c/q/ebSDl7ubgyIrMiD7aoR6qdZj0VExHkp7Di69VNg41TABn2nQvkmBXr4PTGJvLV0P0t3nwHMpR1ubxbBqJurEx5YokDfS0RExAoKO45s32JYPNb8ueMEqNO9wA69/0wS7yzbz487TgPgZoPejSvwSPsaVAz2LbD3ERERsZrCjqM6vQO+GQYY0GQItB5VIIc9eDaZ95YfYMEfi3SCuX7Vox1qUj20VIG8h4iIiCNR2HFESWdgTj9ziHmVG+GWN695iPnR8ym8tzyaeVtP8McanXSpF8ajHWtQO8wJV3oXERHJI4UdR5N5EeYOgMQTEFz9moeYH49LZfKKA3y75STZf6ScDnVCebRDTS3SKSIixYLCjiMxDPh+JJzcBD6BMOArKFH6qg514kIqH/wSzdebTpD1R8hpW7MMozvW5LqIwIKrWURExMEp7DiS1W/Azm/+WMX8cwiulu9DnIy/+EfIOU5mthlybqgRwqMdatK00tUFJxEREWemsOModn8Pv7xk/nzLm1Dlhny9/NQfIeerv4Wc1tWCGd2xJs0rBxV0tSIiIk5DYccRxGyDefebP7d8EJreneeXnoq/yIcro/lq4wkysu0AtKoazKMdahBZNbgQihUREXEuCjtWSzoDX/SHzFSo3gE6vpinl/0Zcr7c+NeVnFZVg3mkQw1aKuSIiIjkUNixUlY6fDkIEk9CcA3oOw3cr/yf5GT8RT78x+2qllWDeKR9TVpVU8gRERH5J4UdqxgGLBoNJzaYi3r2nwslAi+7+/G4VD5ceZBvNutKjoiISH4o7Fhl3UcQNRtsbnDbpxBS/ZK7HTtvDiH/dstfQ8hbVwvmkfbqkyMiIpIXCjtWOLgCfn7a/LnTS1C9/b92OXwuhQ9+iWbe1r8mA2xTPYRHOtTQ6CoREZF8UNgpanGH4OuhYNih0QBz9NXfRMcm88Ev0XwfdTJnWYcba5bhkfbVaVpJIUdERCS/FHaKUnoSfDEA0uKhfDO49e2cNa/2nU5i8ooD/LAjJmeBzptrhzLq5uo0rqjJAEVERK6Wwk5Rsdth/gNwdg+UCjNnSPb0YefJBCavOMCSXWdydu1YtywP31yDBhW0dpWIiMi1UtgpKr++CXsWgrsX3Pk5W+J9eP+7jazYGwuYF3i61g9j1M01qFNOq5CLiIgUFIWdorB/CfzyMgAHW4zn+SV2foteA4CbDbo3Cmdku+rUKOtnZZUiIiIuSWGnsJ2Lxvj2HmwYLClxK/f9Ugk4h4ebjd6Ny/Ngu+pUCSlpdZUiIiIuS2GnENkvJpI6605KpSey0V6TkRfuwMvdjdubVeD+ttWICPK1ukQRERGXp7BTSAy7nfXv9KNVejSnjdKMNsYw+Poa3HtjVcICfKwuT0REpNhQ2CkkNiMbn4AyZJ5xZ1n91/m+aw+CS3lbXZaIiEixo7BTWNw9qXr3J6TFHWBQhbpWVyMiIlJsKewUogBfT/BV0BEREbGSm9UFiIiIiBQmhR0RERFxaQo7IiIi4tIsDTurV6+me/fuhIeHY7PZmD9/fq7nX3jhBWrXrk3JkiUpXbo0HTp0YP369bn2qVy5MjabLddj0qRJRfgpRERExJFZGnZSUlJo1KgRH3zwwSWfr1mzJu+//z47duzgt99+o3LlynTq1ImzZ8/m2m/ChAnExMTkPEaNGlUU5YuIiIgTsHQ0VteuXenatetlnx8wYECu9ltvvcW0adPYvn077du3z9nu5+dHWFhYodUpIiIizstp+uxkZGTw8ccfExAQQKNGjXI9N2nSJIKDg2ncuDGvv/46WVlZVzxWeno6iYmJuR4iIiLimhx+np1FixbRr18/UlNTKVeuHEuXLiUkJCTn+YcffpgmTZoQFBTEmjVrGDt2LDExMbz11luXPebEiRMZP358UZQvIiIiFrMZhmFYXQSAzWZj3rx59OrVK9f2lJQUYmJiOHfuHFOnTmXFihWsX7+e0NDQSx5n+vTp3HfffSQnJ+PtfenlGdLT00lPT89pJyYmEhERQUJCAv7+/gX2mURERKTwJCYmEhAQ8J/f3w5/G6tkyZJUr16dli1bMm3aNDw8PJg2bdpl94+MjCQrK4sjR45cdh9vb2/8/f1zPURERMQ1OXzY+Se73Z7rqsw/RUVF4ebmdtkrPyIiIlK8WNpnJzk5mejo6Jz24cOHiYqKIigoiODgYF5++WV69OhBuXLlOHfuHB988AEnT57k9ttvB2Dt2rWsX7+edu3a4efnx9q1axk9ejSDBg2idOnSVn0sERERcSCWhp1NmzbRrl27nPaYMWMAGDJkCFOmTGHv3r3MnDmTc+fOERwcTPPmzfn111+pV68eYN6Omjt3Li+88ALp6elUqVKF0aNH5xxHRERExGE6KFspISGBwMBAjh8/rv47IiIiTuLPAUbx8fEEBARcdj+HH3peFJKSkgCIiIiwuBIRERHJr6SkpCuGHV3Zwez0fOrUKfz8/LDZbAV23D8Tp64YFS6d56Kjc100dJ6Lhs5z0SjM82wYBklJSYSHh+PmdvkxV7qyA7i5uVGhQoVCO76GtxcNneeio3NdNHSei4bOc9EorPN8pSs6f3K6oeciIiIi+aGwIyIiIi5NYacQeXt78/zzz1922QopGDrPRUfnumjoPBcNneei4QjnWR2URURExKXpyo6IiIi4NIUdERERcWkKOyIiIuLSFHZERETEpSnsXKMPPviAypUr4+PjQ2RkJBs2bLji/l9//TW1a9fGx8eHBg0a8OOPPxZRpc4tP+d56tSp3HDDDZQuXZrSpUvToUOH//zvIn/J7+/0n+bOnYvNZqNXr16FW6CLyO95jo+P56GHHqJcuXJ4e3tTs2ZN/f3Ig/ye53feeYdatWpRokQJIiIiGD16NGlpaUVUrXNavXo13bt3Jzw8HJvNxvz58//zNStXrqRJkyZ4e3tTvXp1ZsyYUbhFGnLV5s6da3h5eRnTp083du3aZYwYMcIIDAw0zpw5c8n9f//9d8Pd3d147bXXjN27dxvPPPOM4enpaezYsaOIK3cu+T3PAwYMMD744ANj69atxp49e4y7777bCAgIME6cOFHElTuf/J7rPx0+fNgoX768ccMNNxg9e/YsmmKdWH7Pc3p6utGsWTOjW7duxm+//WYcPnzYWLlypREVFVXElTuX/J7n2bNnG97e3sbs2bONw4cPG0uWLDHKlStnjB49uogrdy4//vij8fTTTxvfffedARjz5s274v6HDh0yfH19jTFjxhi7d+82Jk+ebLi7uxuLFy8utBoVdq5BixYtjIceeiinnZ2dbYSHhxsTJ0685P533HGHccstt+TaFhkZadx3332FWqezy+95/qesrCzDz8/PmDlzZmGV6DKu5lxnZWUZrVu3Nj755BNjyJAhCjt5kN/z/NFHHxlVq1Y1MjIyiqpEl5Df8/zQQw8ZN998c65tY8aMMa6//vpCrdOV5CXsPPHEE0a9evVybbvzzjuNzp07F1pduo11lTIyMti8eTMdOnTI2ebm5kaHDh1Yu3btJV+zdu3aXPsDdO7c+bL7y9Wd539KTU0lMzOToKCgwirTJVztuZ4wYQKhoaEMHz68KMp0eldznhcsWECrVq146KGHKFu2LPXr1+eVV14hOzu7qMp2Oldznlu3bs3mzZtzbnUdOnSIH3/8kW7duhVJzcWFFd+FWgj0Kp07d47s7GzKli2ba3vZsmXZu3fvJV9z+vTpS+5/+vTpQqvT2V3Nef6nJ598kvDw8H/945LcruZc//bbb0ybNo2oqKgiqNA1XM15PnToECtWrGDgwIH8+OOPREdH8+CDD5KZmcnzzz9fFGU7nas5zwMGDODcuXO0adMGwzDIysri/vvvZ9y4cUVRcrFxue/CxMRELl68SIkSJQr8PXVlR1zapEmTmDt3LvPmzcPHx8fqclxKUlISgwcPZurUqYSEhFhdjkuz2+2Ehoby8ccf07RpU+68806efvpppkyZYnVpLmXlypW88sorfPjhh2zZsoXvvvuOH374gRdffNHq0uQa6crOVQoJCcHd3Z0zZ87k2n7mzBnCwsIu+ZqwsLB87S9Xd57/9MYbbzBp0iSWLVtGw4YNC7NMl5Dfc33w4EGOHDlC9+7dc7bZ7XYAPDw82LdvH9WqVSvcop3Q1fxOlytXDk9PT9zd3XO21alTh9OnT5ORkYGXl1eh1uyMruY8P/vsswwePJh77rkHgAYNGpCSksK9997L008/jZubrg8UhMt9F/r7+xfKVR3QlZ2r5uXlRdOmTVm+fHnONrvdzvLly2nVqtUlX9OqVatc+wMsXbr0svvL1Z1ngNdee40XX3yRxYsX06xZs6Io1enl91zXrl2bHTt2EBUVlfPo0aMH7dq1IyoqioiIiKIs32lcze/09ddfT3R0dE6YBNi/fz/lypVT0LmMqznPqamp/wo0fwZMQ8tIFhhLvgsLretzMTB37lzD29vbmDFjhrF7927j3nvvNQIDA43Tp08bhmEYgwcPNp566qmc/X///XfDw8PDeOONN4w9e/YYzz//vIae50F+z/OkSZMMLy8v45tvvjFiYmJyHklJSVZ9BKeR33P9TxqNlTf5Pc/Hjh0z/Pz8jJEjRxr79u0zFi1aZISGhhovvfSSVR/BKeT3PD///POGn5+f8cUXXxiHDh0yfv75Z6NatWrGHXfcYdVHcApJSUnG1q1bja1btxqA8dZbbxlbt241jh49ahiGYTz11FPG4MGDc/b/c+j5448/buzZs8f44IMPNPTc0U2ePNmoWLGi4eXlZbRo0cJYt25dznNt27Y1hgwZkmv/r776yqhZs6bh5eVl1KtXz/jhhx+KuGLnlJ/zXKlSJQP41+P5558v+sKdUH5/p/9OYSfv8nue16xZY0RGRhre3t5G1apVjZdfftnIysoq4qqdT37Oc2ZmpvHCCy8Y1apVM3x8fIyIiAjjwQcfNC5cuFD0hTuRX3755ZJ/c/88t0OGDDHatm37r9dcd911hpeXl1G1alXj008/LdQabYaha3MiIiLiutRnR0RERFyawo6IiIi4NIUdERERcWkKOyIiIuLSFHZERETEpSnsiIiIiEtT2BERERGXprAjIiIiLk1hR0RERFyawo6IiIi4NIUdEZFLuOmmm3j00Uev6RiGYXDvvfcSFBSEzWYjKiqqQGoTkfxR2BERhzd06FCeeeYZq8vIt8WLFzNjxgwWLVpETEwM9evXt7okkWLJw+oCRESuJDs7m0WLFvHDDz9YXUq+HTx4kHLlytG6devL7pORkYGXl1cRViVS/OjKjkgx9cUXX1CiRAliYmJytg0dOpSGDRuSkJBw1cetUKECH374Ya5ta9aswdfXl6NHj+b7eGvWrMHT05PmzZtf8vmbbrqJUaNG8eijj1K6dGnKli3L1KlTSUlJYejQofj5+VG9enV++umnXK9LT0/n4YcfJjQ0FB8fH9q0acPGjRsvW4fdbmfixIlUqVKFEiVK0KhRI7755pvL7n/33XczatQojh07hs1mo3Llyjn1jhw5kkcffZSQkBA6d+7M4sWLadOmDYGBgQQHB3Prrbdy8ODBf73/a6+9RvXq1fH29qZixYq8/PLLeTyLIsWbwo5IMdWvXz9q1qzJK6+8AsDzzz/PsmXL+OmnnwgICLjq40ZGRuYKDYZh8OijjzJ69GgqVaqU7+MtWLCA7t27Y7PZLrvPzJkzCQkJYcOGDYwaNYoHHniA22+/ndatW7NlyxY6derE4MGDSU1NzXnNE088wbfffsvMmTPZsmUL1atXp3PnzsTFxV3yPSZOnMisWbOYMmUKu3btYvTo0QwaNIhVq1Zdcv93332XCRMmUKFCBWJiYnKdk5kzZ+Ll5cXvv//OlClTSElJYcyYMWzatInly5fj5uZG7969sdvtOa8ZO3YskyZN4tlnn2X37t3MmTOHsmXL5vd0ihRPhogUWwsXLjS8vb2Nl156yShdurSxc+fOnOd69eplBAYGGn379s3XMV977TWjXr16Oe2ZM2caYWFhRlJS0lUdt0aNGsaiRYsu+3zbtm2NNm3a5LSzsrKMkiVLGoMHD87ZFhMTYwDG2rVrDcMwjOTkZMPT09OYPXt2zj4ZGRlGeHi48dprr+Uc95FHHjEMwzDS0tIMX19fY82aNbnee/jw4Ub//v0vW9vbb79tVKpU6V/1Nm7c+Iqf+ezZswZg7NixwzAMw0hMTDS8vb2NqVOnXvF1InJpurIjUozdeuut1K1blwkTJjBv3jzq1auX89wjjzzCrFmz8n3Mli1bsmfPHpKTk0lJSWHcuHG89NJLlCpVKt/H3bNnD6dOnaJ9+/ZX3K9hw4Y5P7u7uxMcHEyDBg1ytv15BSQ2NhYw+9JkZmZy/fXX5+zj6elJixYt2LNnz7+OHx0dTWpqKh07dqRUqVI5j1mzZv3rdlNeNG3aNFf7wIED9O/fn6pVq+Lv759zy+vYsWOAeR7S09P/8zyIyKWpg7JIMbZ48WL27t1Ldnb2v26J3HTTTaxcuTLfx2zatClubm5s2bKFZcuWUaZMGYYOHXpVx12wYAEdO3bEx8fnivt5enrmattstlzb/rwF9vfbQvmRnJwMwA8//ED58uVzPeft7Z3v45UsWTJXu3v37lSqVImpU6cSHh6O3W6nfv36ZGRkAFCiRImrqltETLqyI1JMbdmyhTvuuINp06bRvn17nn322QI5rq+vLw0aNODbb7/ljTfe4O2338bN7er+1Hz//ff07NmzQOr6u2rVquX0mflTZmYmGzdupG7duv/av27dunh7e3Ps2DGqV6+e6xEREXFNtZw/f559+/bxzDPP0L59e+rUqcOFCxdy7VOjRg1KlCjB8uXLr+m9RIorXdkRKYaOHDnCLbfcwrhx43Jun7Rq1YotW7bQpEmTaz5+y5YtmTx5Mj179uSmm266qmPExsayadMmFixYcM31/FPJkiV54IEHePzxxwkKCqJixYq89tprpKamMnz48H/t7+fnx//93/8xevRo7HY7bdq0ISEhgd9//x1/f3+GDBly1bWULl2a4OBgPv74Y8qVK8exY8d46qmncu3j4+PDk08+yRNPPIGXlxfXX389Z8+eZdeuXTn1vv/++8ybN0+BSOQSFHZEipm4uDi6dOlCz549c75UIyMj6dq1K+PGjWPx4sX/eYwZM2YwdOhQDMO45PONGjXC09OT119//arrXLhwIS1atCAkJOSqj3ElkyZNwm63M3jwYJKSkmjWrBlLliyhdOnSl9z/xRdfpEyZMkycOJFDhw4RGBhIkyZNGDdu3DXV4ebmxty5c3n44YepX78+tWrV4r333vtXSHz22Wfx8PDgueee49SpU5QrV477778/5/lz585dVf8hkeLAZlzur5WIFHsrV67k/fff/9d8Ms8//zyrVq26bN+bdu3a0aRJE9588818HffvevToQZs2bXjiiSeuun4REdCVHRG5jA4dOrBt2zZSUlKoUKECX3/9Na1atQLgp59+4v3338+1v91u5+zZs0ybNo0DBw7w/fff5/u4f9emTRv69+9f8B9MRIodXdkRkQKxcuVKbr75ZmrXrs2nn35KZGSk1SWJiAAKOyIiIuLiNPRcREREXJrCjoiIiLg0hR0RERFxaQo7IiIi4tIUdkRERMSlKeyIiIiIS1PYEREREZemsCMiIiIuTWFHREREXJrCjoiIiLi0/wfdBWWF3p0gwAAAAABJRU5ErkJggg==", "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.14" } }, "nbformat": 4, "nbformat_minor": 5 }