{ "cells": [ { "cell_type": "markdown", "id": "1f3a34f2", "metadata": {}, "source": [ "[![Open in Colab](https://colab.research.google.com/assets/colab-badge.svg)](https://colab.research.google.com/github/cbc-univie/mdy-newanalysis-package/blob/master/docs/notebooks/msd.ipynb)" ] }, { "cell_type": "markdown", "id": "5515278d", "metadata": {}, "source": [ "# Mean squared displacement" ] }, { "cell_type": "markdown", "id": "1fb8776a", "metadata": {}, "source": [ "We calculate the mean squared displacement (msd) as well as the msdmj to obtain diffusion and conductivity." ] }, { "cell_type": "code", "execution_count": 15, "id": "ff3dba1f", "metadata": {}, "outputs": [], "source": [ "import time\n", "import sys\n", "import numpy as np\n", "import matplotlib.pyplot as plt" ] }, { "cell_type": "markdown", "id": "96bb3bd0", "metadata": {}, "source": [ "Some configurations if we are in a google colab environment:" ] }, { "cell_type": "code", "execution_count": 24, "id": "7e80d48f", "metadata": {}, "outputs": [], "source": [ "IN_COLAB = 'google.colab' in sys.modules\n", "HAS_NEWANALYSIS = 'newanalysis' in sys.modules\n", "if IN_COLAB:\n", " if not'MDAnalysis' in sys.modules:\n", " !pip install MDAnalysis\n", " import os\n", " if not os.path.exists(\"/content/mdy-newanalysis-package/\"):\n", " !git clone https://github.com/cbc-univie/mdy-newanalysis-package.git\n", " !pwd\n", " if not HAS_NEWANALYSIS:\n", " %cd /content/mdy-newanalysis-package/newanalysis_source/\n", " !pip install .\n", " %cd /content/mdy-newanalysis-package/docs/notebooks/\n", "#if IN_COLAB and not HAS_NEWANALYSIS:\n", "# !pip install \"git+https://github.com/cbc-univie/mdy-newanalysis-package.git@master#egg=newanalysis&subdirectory=newanalysis_source\"" ] }, { "cell_type": "markdown", "id": "95da3a41", "metadata": {}, "source": [ "Next we import MDAnalysis and functions needed from newanalysis." ] }, { "cell_type": "code", "execution_count": 16, "id": "20845533", "metadata": {}, "outputs": [], "source": [ "try:\n", " import MDAnalysis\n", "except ModuleNotFoundError as e:\n", " if not IN_COLAB:\n", " !conda install --yes --prefix {sys.prefix} MDAnalysis\n", " import MDAnalysis\n", " else:\n", " print(\"Probably the cell above had a problem.\")\n", " raise e\n", "try:\n", " from newanalysis.diffusion import unfoldTraj, msdCOM, msdMJ\n", "except ModuleNotFoundError as e:\n", " if not IN_COLAB:\n", " print(\"Install newanalysis fist.\")\n", " print(\"Go to the newanalysis_source folder and run:\")\n", " print(\"conda activate newanalysis-dev\")\n", " print(\"pip install .\")\n", " else:\n", " print(\"Probably the cell above had a problem.\")\n", " raise e" ] }, { "cell_type": "markdown", "id": "6f8c1f63", "metadata": {}, "source": [ "## Preprocessing" ] }, { "cell_type": "markdown", "id": "d9a4e344", "metadata": {}, "source": [ "Next we create our MDAnalysis Universe:" ] }, { "cell_type": "code", "execution_count": 17, "id": "d664f2b4", "metadata": {}, "outputs": [], "source": [ "base='../../test_cases/data/emim_dca_equilibrium/'\n", "psf=base+'emim_dca.psf'\n", "\n", "skip=20\n", "u=MDAnalysis.Universe(psf,base+\"/emim_dca.dcd\")\n", "boxl = np.float64(round(u.coord.dimensions[0],4))\n", "dt=round(u.trajectory.dt,4)\n", "\n", "n = int(u.trajectory.n_frames/skip)\n", "if u.trajectory.n_frames%skip != 0:\n", " n+=1" ] }, { "cell_type": "markdown", "id": "fe28f64e", "metadata": {}, "source": [ "Now we do the selections" ] }, { "cell_type": "code", "execution_count": 18, "id": "998f92e6", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Number EMIM = 1000\n", "Number DCA = 1000\n" ] } ], "source": [ "sel_cat1 = u.select_atoms(\"resname EMIM\")\n", "ncat1 = sel_cat1.n_residues\n", "com_cat1 = np.zeros((n,ncat1,3),dtype=np.float64)\n", "print(\"Number EMIM = \",ncat1)\n", "\n", "sel_an1 = u.select_atoms(\"resname DCA\")\n", "nan1 = sel_an1.n_residues\n", "com_an1 = np.zeros((n,nan1, 3),dtype=np.float64)\n", "print(\"Number DCA = \",nan1)" ] }, { "cell_type": "markdown", "id": "f5d1b632", "metadata": {}, "source": [ "Loop throguh the trajectory and get the center of masses" ] }, { "cell_type": "code", "execution_count": 19, "id": "d78f9648", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "\u001b[1AFrame 0 of 100 \tElapsed time: 0.00 hours\n", "\u001b[1AFrame 20 of 100 \tElapsed time: 0.00 hours\n", "\u001b[1AFrame 40 of 100 \tElapsed time: 0.00 hours\n", "\u001b[1AFrame 60 of 100 \tElapsed time: 0.00 hours\n", "\u001b[1AFrame 80 of 100 \tElapsed time: 0.00 hours\n" ] } ], "source": [ "ctr=0\n", "start=time.time()\n", "print(\"\")\n", "\n", "for ts in u.trajectory[::skip]:\n", " print(\"\\033[1AFrame %d of %d\" % (ts.frame,len(u.trajectory)), \"\\tElapsed time: %.2f hours\" % ((time.time()-start)/3600))\n", " com_an1[ctr] = sel_an1.center_of_mass(compound='residues')\n", " com_cat1[ctr]= sel_cat1.center_of_mass(compound='residues')\n", " ctr+=1" ] }, { "cell_type": "markdown", "id": "a0542745", "metadata": {}, "source": [ "## Using newanalysis functions" ] }, { "cell_type": "markdown", "id": "d509c255", "metadata": {}, "source": [ "We now use the unfoldTraj function to undo all jumps due to PBC. This is done in-place." ] }, { "cell_type": "code", "execution_count": 20, "id": "76748098", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "unfolding coordinates ..\n" ] } ], "source": [ "print(\"unfolding coordinates ..\")\n", "unfoldTraj(com_cat1,boxl); #; to suppress output\n", "unfoldTraj(com_an1, boxl);" ] }, { "cell_type": "markdown", "id": "a79ea958", "metadata": {}, "source": [ "Next we can calculate msd and msdmj" ] }, { "cell_type": "code", "execution_count": 21, "id": "7fce42ea", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "calculating msd ..\n", "calculating msdMJ ..\n" ] } ], "source": [ "print(\"calculating msd ..\")\n", "msd_an1 = msdCOM(com_an1)\n", "msd_cat1 = msdCOM(com_cat1)\n", "\n", "#fan1=open(\"msd_an.dat\",'w')\n", "#fcat1=open(\"msd_cat.dat\",'w')\n", "\n", "#for i in range(ctr):\n", "# fan1.write(\"%5.5f\\t%5.5f\\n\" % (i*skip*dt, msd_an1[i]))\n", "# fcat1.write(\"%5.5f\\t%5.5f\\n\" % (i*skip*dt, msd_cat1[i]))\n", " \n", "print(\"calculating msdMJ ..\")\n", "msdmj = msdMJ(com_cat1,com_an1)\n", "#f1=open(\"msdMJ.dat\",'w')\n", "#for i in range(ctr):\n", "# f1.write(\"%5.5f\\t%5.5f\\n\" % (i*skip*dt, msdmj[i]))\n", "#f1.close()" ] }, { "cell_type": "markdown", "id": "51b054b7", "metadata": {}, "source": [ "## Plot the results" ] }, { "cell_type": "code", "execution_count": 22, "id": "0b89b9e7", "metadata": { "tags": [ "nbsphinx-thumbnail" ] }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAhYAAAGxCAYAAAA+tv8YAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/bCgiHAAAACXBIWXMAAA9hAAAPYQGoP6dpAABUy0lEQVR4nO3dd3hUddrG8e+kk5AEQgklAUIChN6b0kGUoqCiWHBRd10L9l1d29pd9F13LatiA+yCShUUQYEAolKk9xRIKCGEkl5nzvvHCQmRljKTKbk/15XLySkzz8kBcnvOeX4/i2EYBiIiIiJ24OXsAkRERMRzKFiIiIiI3ShYiIiIiN0oWIiIiIjdKFiIiIiI3ShYiIiIiN0oWIiIiIjdKFiIiIiI3ShYiIiIiN0oWIgIAB999BEWiwWLxcLKlSvPWm8YBjExMVgsFoYMGVK6/Pjx4zz++ON06NCBoKAgQkNDiY2N5ZZbbmHr1q3nfH+LxUJAQABNmjRh6NChTJ06lbS0tBo4ShFxNB9nFyAiriU4OJjp06eXCw8AcXFxJCQkEBwcXLosOzubfv36kZ2dzSOPPELXrl3Jy8tj7969zJ07l82bN9OlS5dy7zNz5kxiY2MpKioiLS2NNWvW8Morr/Dqq68ye/ZsRowYUROHKSIOomAhIuVMnDiRzz//nLfffpuQkJDS5dOnT6d///5kZmaWLvv666+Jj49n+fLlDB06tNz7PPzww9hstrPev1OnTvTq1av0+2uvvZaHHnqIAQMGcM0117Bv3z7Cw8MdcGQiUhN0K0REyrnxxhsB+PLLL0uXZWRkMGfOHG6//fZy2x4/fhyApk2bnvO9vLwq9k9MixYt+M9//kNWVhbvvfdeVcoWERehYCEi5YSEhDBhwgRmzJhRuuzLL7/Ey8uLiRMnltu2f//+APzpT39i/vz5pUGjKkaPHo23tzerVq2q8nuIiPMpWIjIWW6//XbWrVvHjh07AJgxYwbXXXdduecrAC699FKef/55tmzZwtVXX03Dhg1p3bo1d999d7kHNysiKCiIhg0bcvjwYbsdh4jUPAULETnL4MGDiY6OZsaMGWzbto3169efdRvktH/+858kJyczY8YM7rzzTurWrcu7775Lz549y91OqQjDMOxRvog4kYKFiJzFYrFw22238dlnn/Huu+/Stm1bBg4ceN7tw8PDue2223j33XfZunUrcXFx+Pn58cADD1T4M3Nycjh+/DjNmjWzxyGIiJMoWIjIOd16662kp6fz7rvvctttt1Vq30GDBjFy5EiOHTtW4fEpFi9ejNVqPavNVUTci9pNReScmjdvziOPPMLu3buZPHnyObc5evQojRo1Oqv7w2q1sm/fPgIDA6lXr95FPys5OZm///3vhIaGcuedd9qjfBFxEgULETmvl19++YLrP/30U9577z1uuukmevfuTWhoKAcPHuTDDz9kx44dPP300/j5+ZXbZ/v27RQXF1NcXExaWhqrV69m5syZeHt7M2/ePBo1auTIQxIRB1OwEJEqGzNmDKmpqXz33XdMmzaNkydPEhwcTJcuXfj000+ZNGnSWfucvq3i5+dHvXr1aN++Pf/4xz/4y1/+olAh4gEshh7DFhERETvRw5siIiJiNwoWIiIiYjcKFiIiImI3ChYiIiJiNwoWIiIiYjcKFiIiImI3NT6Ohc1m4/DhwwQHB2OxWGr640VERKQKDMMgKyuLZs2anTXa7plqPFgcPnyYyMjImv5YERERsYOUlBQiIiLOu77Gg0VwcDBgFhYSElLTHy8iIiJVkJmZSWRkZOnv8fOp8WBx+vZHSEiIgoWIiIibudhjDHp4U0REROxGwUJERETsRsFCRERE7KZSz1g8++yzPPfcc+WWhYeHk5qaateirFYrRUVFdn3P2szX1xdvb29nlyEiIrVApR/e7NixIz/++GPp9/b+hZWdnc3BgwfRbO72Y7FYiIiIoG7dus4uRUREPFylg4WPjw9NmjRxRC1YrVYOHjxIYGAgjRo10gBadmAYBseOHePgwYO0adNGVy5ERMShKh0s9u3bR7NmzfD396dv377861//onXr1ufdvqCggIKCgtLvMzMzz7ttUVERhmHQqFEj6tSpU9nS5DwaNWrE/v37KSoqUrAQERGHqtTDm3379uWTTz7hhx9+4IMPPiA1NZVLLrmE48ePn3efqVOnEhoaWvpVkVE3daXCvvTzFBGRmmIxqvEwQ05ODtHR0Tz66KM8/PDD59zmXFcsIiMjycjIOGuArPz8fJKSkoiKiiIgIKCqZckf6OcqIiLVlZmZSWho6Dl/f5+pWiNvBgUF0blzZ/bt23febfz9/fH396/Ox4iIiIibqNY4FgUFBezatYumTZvaqx4RERFxY5UKFn//+9+Ji4sjKSmJ3377jQkTJpCZmcnkyZMdVZ9buPXWW7FYLFgsFnx9fQkPD+eyyy5jxowZ2Gy2cttu2rSJ6667jvDwcAICAmjbti133HEHe/fuPet9R44cibe3N7/++mtNHYqIiEi1VCpYHDx4kBtvvJF27dpxzTXX4Ofnx6+//krLli0dVZ/buOKKKzhy5Aj79+/n+++/Z+jQoTzwwAOMHTuW4uJiABYtWkS/fv0oKCjg888/Z9euXXz66aeEhobyz3/+s9z7JScn88svv3Dvvfcyffp0ZxySiIi4kfwiK5/+sp+/fbXFqXVU6hmLWbNmOaqOczIMg7wia41+5ml1fL0r1U3h7+9fOr5H8+bN6dGjB/369WP48OF89NFH3HTTTdx2222MHj2aefPmle4XFRVF3759OXXqVLn3mzlzJmPHjuXuu++mT58+vP766wQFBdnl2ERExHPkFhbzxW/JvLcqkWNZZrPETX1b0LNlfafUU+PTpldGXpGVDk//4JTP3vn85QT6Ve/HM2zYMLp27crcuXNp0KAB6enpPProo+fctl69eqWvDcNg5syZvP3228TGxtK2bVu++uorbrvttmrVIyIiniMzv4hPfznAh6sTOZlrToPRLDSAu4ZE07HZ+bs2HM2lg4UniI2NZevWraWdM7GxsRfd58cffyQ3N5fLL78cgEmTJjF9+nQFCxER4WROITN/TmLm2v1k5Zu32ls2CGTKkBjGd2+On49z5xd16WBRx9ebnc9f7rTPtgfDMLBYLJWa+2T69OlMnDgRHx/z9Nx444088sgj7Nmzh3bt2tmlLhERcS9pWflMX53Ep78eILfQfEygTeO6TBkaw9guTfHxdo0Jy106WFgslmrfjnC2Xbt2ERUVRdu2bQHYvXs3/fv3P+/2J06cYP78+RQVFTFt2rTS5VarlRkzZvDKK684vGYREXEdh0/l8f6qRL5cl0xBsdlp2LFZCPcNi2FkhyZ4ebnW6Mru/VvbxS1fvpxt27bx0EMPMXLkSBo2bMj//d//lXt487RTp05Rr149Pv/8cyIiIpg/f3659T/99BNTp07lpZdeKr2SISIiniv5eC7T4uL5ZuNBiqzmVe9ukfW4f3gMQ9s1dtnpGvQbyk4KCgpITU3FarVy9OhRlixZwtSpUxk7dix/+tOf8Pb25sMPP+S6667jqquu4v777ycmJob09HS++uorkpOTmTVrFtOnT2fChAl06tSp3Pu3bNmSf/zjHyxevJhx48Y56ShFRMTR4tOyeGdFAgu2HMZqMwNFv9Zh3DesDZdEN3DZQHGagoWdLFmyhKZNm+Lj40P9+vXp2rUrb775JpMnT8bLy7zvNW7cONauXcvUqVO56aabSudNGTZsGC+++CIbN25ky5YtfPDBB2e9f3BwMCNHjmT69OkKFiIiHmjn4UzeXhHPd9uPcPqxvMFtG3HvsBh6twpzbnGVUK1JyKriQpOYaLIsx9DPVUTEdW1OOcVby/fx46600mUjO4Rz77AYukTUc15hf1Ajk5CJiIhI1fyWeJy3VsSzel86ABYLjO3SjClDo4lt4rxxKKpLwUJERKSGGIbB6n3pvLU8nnX7TwDg7WXh6u7NuXtINNGN6jq5wupTsBAREXEwwzD4cVcaby3fx5aDGQD4eXtxXa8I7hocTWRYoJMrtB8FCxEREQex2gy+336Et5bHszs1C4AAXy9u6tOSvw5qTZNQz3vuTcFCRETEzoqsNhZuPszbK+NJPJYDQJCfN3+6pBV/HhBFw7r+Tq7QcRQsRERE7KSg2MqcjYeYFhdPyok8AEICfLh9QBS3XtKKeoF+Tq7Q8RQsREREqimv0Mqs9cm8F5dIamY+AA2C/PjLwNZM6teC4ABfJ1dYcxQsREREqii7oJjPfjWnLk/PLgQgPMSfOwdFc2OfFtTxs8+Elu5EwUJERKSSMnKL+Gjtfmb8nERGXhEAEfXrcPeQaCb0jMDfp/YFitMULERERCroeHYB09ck8ckvB8guKAagdcMg7hkaw7huzfB1kanLnUk/ATu49dZbsVgsZ31dccUVALRq1QqLxcKsWbPO2rdjx45YLBY++uij0mWtWrXi9ddfL/d9ZfYXERH7OpqZzwuLdjLglRW8szKB7IJi2oUH878bu7Ps4cFM6BmhUFFCVyzs5IorrmDmzJnllvn7l7UTRUZGMnPmTG644YbSZb/++iupqakEBQVd9P2ru7+IiFTewZO5vBuXwFcbDlJYbAOgS0Qo9w6NYUT7cLy8XHumUWdw7WBhGFCU65zP9g00B26vIH9/f5o0aXLe9TfffDOvvfYaKSkpREZGAjBjxgxuvvlmPvnkk4u+f3X3FxGRiktKz+GdFfHM23SI4pKpy3u1rM99w9swqE1Dl5+63JlcO1gU5cK/mjnns584DH72uxIQHh7O5Zdfzscff8xTTz1Fbm4us2fPJi4urkLBoLr7i4jIxe1JzeLtFfEs2nqYkjzBgJiG3Dsshr5RYQoUFaAbQnayaNEi6tatW+7rhRdeKLfN7bffzkcffYRhGHzzzTdER0fTrVu3Cn9GdfcXEZFz23Ywgzs/3cDlr69i4RYzVAyPbczcey7hs7/0pV/rBgoVFeTaVyx8A80rB8767EoYOnQo06ZNK7csLCys3PdjxozhzjvvZNWqVcyYMYPbb7+9Up9R3f1FRKS8jQdO8tbyfazYcwww74CP6tSEe4bE0Kl5qJOrc0+uHSwsFrvejnCkoKAgYmJiLriNj48Pt9xyC8888wy//fYb8+bNq9RnVHd/ERExZxr9JfE4by2PZ23CcQC8LDCuW3PuGRJNm/BgJ1fo3lw7WHig22+/nVdffZWJEydSv379Gt9fRKS2MgyDlXuP8dbyeDYeOAmAr7eFa3uYU5e3auge/yPr6hQs7KSgoIDU1NRyy3x8fGjYsGG5Ze3btyc9PZ3AwMrdarHX/iIitY3NZrB051HeWrGP7YcyAfDz8eLG3pH8dXA0zevVcXKFnkXBwk6WLFlC06ZNyy1r164du3fvPmvbBg0aVOuzqru/iEhtYLUZLNp6mLdXxLP3aDYAgX7eTOrXkr8MiKJxSICTK/RMFsMwjJr8wMzMTEJDQ8nIyCAkJKTcuvz8fJKSkoiKiiIgQCfcXvRzFZHapMhqY96mQ0xbmUBSeg4Awf4+3HppK267NIqwIM+futwRLvT7+0y6YiEiIh4hv8jK1xsP8u7KBA6dygOgfqAvfx4QxS39WxFap/ZMXe5MChYiIuLWcguL+eK3ZN5flUhaVgEADev6c+eg1tzUtwVB/vpVV5P00xYREbeUlV/EJ78cYPqaJE7kFALQLDSAu4ZEc32vSAJ8a+/U5c6kYCEiIm7lVG4hM37ez0c/J5GZb05d3iIskClDo7m6ewR+PhpU2plcMljU8POkHk8/TxHxBMeyCvhwTSKf/XKAnEIrADGN6zJlaDRXdmmGj6YtdwkuFSy8vc3LVoWFhdSpo75ieyksNC8Rnv75ioi4kyMZebwXl8iX65IpKJm6vEPTEO4bFsPlHZto6nIX41LBwsfHh8DAQI4dO4avry9eXkqf1WWz2Th27BiBgYH4+LjU6RYRuaDk47lMi0vgm40pFFnNK6/dIutx37AYhsU21qRgLsqlftNYLBaaNm1KUlISBw4ccHY5HsPLy4sWLVroL6GIuIX4tGzeWRnPgs2HsZbMXd43Koz7hrXh0hjNMurqXCpYAPj5+dGmTZvSy/dSfX5+frr6IyIub9eRTN5aEc93245w+tGwQW0bce/QGPpEhV14Z3EZLhcswPw/bI0QKSJSO2xOOcVby+P5cdfR0mWXdQjn3qExdI2s57zCpEpcMliIiIjnW5d0gv8t38fqfekAWCwwpnNTpgyNoX3T8w8ZLa5NwUJERGqMYRisiU/nf8vjWZd0AgBvLwvjuzXnnqHRRDeq6+QKpboULERExOEMw+CnXWn8b0U8W1JOAeDn7cWEXhHcPTiayLBA5xYodqNgISIiDmO1GSzZnspbK+LZdSQTgABfL27s04K/DmpN01CNWeRpFCxERMTuiq02Fm45zNsr4kk4Zk5dHuTnzS39W/HnAVE0CvZ3coXiKAoWIiJiNwXFVub+fohpKxNIPpELQEiAD7ddGsVtl7aiXqCfkysUR1OwEBGRassvsjJrXTLvrUrkSEY+AGFBfvxlYBS39GtJcICvkyuUmqJgISIiVZZTUMznvx3g/VVJpGcXANA42J87B0dzY59IAv30a6a20RkXEZFKy8gr4pO1+5n+cxKncosAaF6vDncPiWZCzwgCfDXpYW2lYCEiIhV2IqeQGWuS+HjtfrIKigGIahjEPUOiGd+9Ob6aurzWU7AQEZGLSsvM54PViXz2azJ5RVYA2oUHM2VYDGM6N8VbU5dLCQULERE5r0On8ngvLoFZ61MoLLYB0Ll5KPcOi+Gy9uF4KVDIHyhYiIjIWfan5zBtZQJzfj9IccnU5T1b1ue+YTEMbttIU5fLeSlYiIhIqX1Hs3h7RTwLtxymJE9waUwD7h3ahn6twxQo5KIULEREhK0HT/HOigSW7EgtXTYstjFThsbQs2V9J1Ym7kbBQkSkljIMg18TT/DOyvjSqcsBRnVqwpShMXRqHurE6sRdKViIiNQyp2cafXtlPJuSTwHm1OXjujbj7iHRtAkPdm6B4tYULEREaoliq43F244wbWUCu1OzAPDz8WJir0j+Oqi1pi4Xu1CwEBHxcAXFVuZsPMS7cWUTg9X192FSv5bcPqAVjYMDnFyheJJqDZE2depULBYLDz74oJ3KERERe8kpKOaDVYkMfGUFT8zbRvKJXMKC/Pj7yLb8/NgwHhsVq1AhdlflKxbr16/n/fffp0uXLvasR0REqulkTiEfrd3Px7/sL53Ho2loAHcMbM0NmhhMHKxKf7qys7O5+eab+eCDD3jxxRcvuG1BQQEFBQWl32dmZlblI0VE5CKOZubz4epEPv8tmdxCc9jtqIZB3D3YnMfDz0fzeIjjVSlYTJkyhTFjxjBixIiLBoupU6fy3HPPVak4ERG5uAPHc3g3LpE5Gw9SaDWH3e7QNIR7hkYzqpPm8ZCaVelgMWvWLH7//XfWr19foe0ff/xxHn744dLvMzMziYyMrOzHiojIH+xOzeSdFQks2lo2SmafVmHcPTSaIRp2W5ykUsEiJSWFBx54gKVLlxIQULEHfvz9/fH3969ScSIicrbfk0/yzop4ftyVVrpsSLtG3DMkhj5RYU6sTAQshmEYFd14/vz5XH311Xh7e5cus1qtWCwWvLy8KCgoKLfuXDIzMwkNDSUjI4OQkJCqVy4iUosYhsGa+HTeXhHPr4knALBYYHTnptw9OFqjZIrDVfT3d6WuWAwfPpxt27aVW3bbbbcRGxvLP/7xj4uGChERqRybzWDpzlTeXpHAtkMZAPh6W7i6e3PuGhxN60Z1nVyhSHmVChbBwcF06tSp3LKgoCAaNGhw1nIREam6IquNhZsPMy0ugfi0bAACfL24sU8L7hjYmmb16ji5QpFzUzOziIgLyS+y8tWGFN6LS+TQqTwAggN8uPWSVtx6SSsa1NUza+Laqh0sVq5caYcyRERqt8z8Ij779QAz1iSRnl0IQMO6/vx5QBST+rUgOMDXyRWKVIyuWIiIONHx7AJm/myOkpmVXwxA83p1uGtwa67rFUmAr55dE/eiYCEi4gSHT+Xx/qpEZq1PJr/IHNQqpnFd7hkSzZVdm+HrrVEyxT0pWIiI1KCEY9m8uzKBeZsOUVwyqlWXiFCmDI3hsvbheGmUTHFzChYiIjVg+6EMpq1M4LvtRzg9elD/1g2YMjSGS2MaaJRM8RgKFiIiDrQu6QRvr4gnbu+x0mUj2odzz9BoerSo78TKRBxDwUJExM4Mw2DlnmO8szKe9ftPAuBlgSu7NuPuIdHENtGow+K5FCxEROzEajP4btsR3lmZwK4jmQD4eXsxoVcEdw5qTcsGQU6uUMTxFCxERKqpsNjGvE0HeTcukaT0HAAC/byZ1K8lfx4QRXhIxSZtFPEEChYiIlWUW1jMl+tS+HB1Ikcy8gGoF+hbOkpmvUA/J1coUvMULEREKikjt4iPf9nPzJ+TOJlbBEDjYH/+Oqg1N/ZpQZC//mmV2kt/+kVEKigtK5/pa5L47JcD5BRaAWjZIJC7BkdzTY/m+PtolEwRBQsRkYtIOZHLe6sS+GrDQQqLzVEyY5sEc/eQaMZ0boqPRskUKaVgISJyHvuOZjFtZQILthzGWjJKZo8W9ZgyNIZhsY01qJXIOShYiIj8weaUU7yzIp6lO4+WLhvYpiFThsbQNypMgULkAhQsREQwB7X6JeE4b6+M5+f44wBYLHBFxybcMySGzhGhTq5QxD0oWIhIrWazGfy46yjvrExgc8opAHy8LIzr1py7h7QmpnGwcwsUcTMKFiJSKxVbbSzaeoR3Vsaz92g2AP4+XtzQO5I7BrUmon6gkysUcU8KFiJSq+QXWflm40HeW5VAyok8AIL9fbilf0tuuzSKRsH+Tq5QxL0pWIhIrZBdUMznvx7gwzVJHMsqAKBBkB+3D4hiUr+WhNbxdXKFItVQlAeJK2H3IjiyFe5cZT4k5AQKFiLi0U7mFDJz7X4+XrufjDxzlMxmoQH8dVBrJvZuQR0/DWolbirvJOxdaoaJ+J+gKKds3eHfoXlPp5SlYCEiHik1I58PVifyxW/J5BWZo2S2bhTEXYOjGd+tOX4+GtRK3FDmYdi92AwT+9eArbhsXUgExI4xv5p0cVqJChYi4lH2p+fwblwCc34/SJHVHNSqU/MQ7hkSw+Udm+DtpTEoxM0c22MGiV2LzCsRZ2rUHtqPNcNE025Ou/1xJgULEfEIOw9nMi0ugcVbD1MySCZ9osKYMjSGQW0aalArcR82GxzaaIaJ3Yvh+L4zVlogsk/JlYmx0CDaaWWej4KFiLi1jQdO8PaKBJbvTitdNiy2MfcMiaZXqzAnViZSCcWFsH91SZj4DrJTy9Z5+ULrIWaYaDcagsOdVmZFKFiIiNsxDINV+9J5e0U865JOAOBlgdGdm3L3kGg6NtMomeIGCrIg/kfzqsTepVCQUbbOLxjaXGbe5oi5DAJCnFdnJSlYiIjbsNkMftiRytsr49l+KBMAX28L1/aI4M7B0UQ1DHJyhSIXkX0M9nxnhonElWAtKFsX1BhiR5u3OKIGgY97jqmiYCEiLq/IamP+pkNMi0sg8ZjZUlfH15ub+rbgLwOjaBpax8kVilzAiaSyTo7kXwGjbF1YazNIxI6FiN7g5f7dSgoWIuKy8gqtzF6fzPurEjmckQ9ASIAPt14axa2XtCIsyM/JFYqcg2FA6lYzTOxaBGk7yq9v2q2kk2MsNIp1iU4Oe1KwEBGXk5lfxKe/HGDGmiSO5xQC0LCuP3cMjOKmvi0IDtAomeJirMWQ/EvJlYnFkJFcts7iDa0uNYNEu9FQL9J5ddYABQsRcRnp2QXMWJPEp78cIKvAHPgnon4d7hoczYSeEQT4apRMcSFFeZCwwrzFsed7yDtRts6nDsQMN8NE28shsPZ0KClYiIjTHTyZywerEpm1PoWCYhsAbcPrcs+QGMZ2aYqPt/vfdxYPkXcS9v5wxjDauWXr6tSHtqPM2xyth4Jf7ZwhV8FCRJwmPi2bd+MSmL/pEMUlo1p1jazHlCHRjGgfjpdGyRRXkHGo/DDahrVsXWhk2TDaLS4Bb/1a1U9ARGrctoMZvLMyniU7UjFKHpC/NKYBU4bE0D+6gUbJFOcyjLJhtHcvgsObyq9v3KGkk2MMNO3qcQ9fVpeChYjUCMMw+C3pBG+viGf1vvTS5SM7hHPP0Bi6RdZzXnEiNhsc2nDGMNrxZ6y0QGTfsisTLjiMtitRsBARhzIMg+W703hnZQIbD5wEwNvLwlVdm3H3kGjahgc7uUKptYoLYf8qsyV0z3eQfbRsnbdf+WG06zZ2WpnuRsFCRByi2Gpj8bYjTFuZwO7ULAD8fLy4vlcEdw6KJjKsdj7YJk5WkAX7lplXJfYthYLMsnX+IdBmpBkmYka41TDarkTBQkTsqqDYytzfD/FuXAIHjptPzAf5eTOpf0v+fGkUjUMCnFyh1DrZaX8YRruwbF3dcPOKRPux0Gqg2w6j7UoULETELvKLrHz26wE+WJ3I0Uxz/oP6gb7cdmkUk/u3IjRQg1pJDTqRWDbyZcpvlB9GO7ps5MvmvTxiGG1XomAhItW2Zl86T83fxv6SKxRNQgK4Y1BrbuwTSaCf/pmRGmAYcGRLWVto2s7y65t1L3n48kpo1E6dHA6kv/EiUmXHsgp4cfFOFmw+DEDjYH8euqwt1/Rojr+PRskUBysdRrukkyMjpWydxRtaDShpCx0NoRHOq7OWUbAQkUqz2Qy+WJfMK0t2k5VfjMUCk/u34m8j22oeD3GswlxIXGHe4tj7vTkS5mm+gWXDaLcZWauG0XYlChYiUik7D2fyxLxtbE45BUCn5iH86+rOdImo59S6xIPlnigbRjth+R+G0Q6DdqPMMNF6SK0dRtuVKFiISIXkFBTz+o97mfHzfqw2g7r+PvxtZFv+1L8V3hp6W+wt4yDs/g52fwv7f/7DMNotzOcl2o+FyH4aRtvF6GyIyEUt3ZHKswt3cDgjH4DRnZvw9NiONAlV66jYiWHAsd3mVYldi+DI5vLrG3cs6eQYA0266OFLF6ZgISLndehUHs8u3MGyneaIhBH16/DCuE4MjdUohGIHNhscXF/28OWJhDNWWqBFv7JhtMNaO61MqRwFCxE5S5HVxkc/7+e1H/eSW2jFx8vCXwe15r5hbajjp24PqYbiAkhabd7i2P0d5KSVrfP2M6cbjx1jPjehYbTdkoKFiJTze/JJnpi7rXQY7t6t6vPS1Z01p4dUXX4mxJcMo713KRRmla07PYx2+7HmMNr++nPm7hQsRASAjNwi/u+H3XyxLhnDgHqBvjwxqj0TekbgpYczpbKy00oGq1oMSXF/GEa7iTm2ROzpYbT9nFen2J2ChUgtZxgGC7cc5oVFO0nPNv/xv7ZHBE+MjqVBXc2bIJVwPKFs5MuUdZQbRrtBTMlgVWOheU8No+3BFCxEarGk9Bz+OX87a+LTAYhuFMSL4zvTP7qBkysTt2AYZvfG6SsTZw2j3aNsTo5G7ZxSotQ8BQuRWqig2Mq0lQm8szKBwmIb/j5e3DcshjsGtdZQ3HJh1mJIXmu2hO5eDJkHy9Z5+ZQNo91uNIQ2d16d4jQKFiK1zNr4dJ6av53E9BwABrZpyIvjO9GyQZCTKxOXVZhrjni5exHsXXKOYbRHmGGi7UioU995dYpLULAQqSXSswt4afEu5m06BECjYH+eHtuBsV2aYtFgQ/JHuSfMELF7McT/BMV5ZevqhJlXJNqXDKPtW8dpZYrrUbAQ8XA2m8HsDSm8/P1uMvKKsFhgUt+W/P3ydoTW0YRhcoaMg2aQ2PUtHFh79jDap0e+1DDacgH6kyHiwXanZvLE3G38nnwKgA5NQ/jXNZ3pFlnPqXWJCzm5H3YuhF0LzVEwzxTeqaSTYww06axhtKVCFCxEPFBuYTFv/LiPD9ckYbUZBPl58/DIdkzu3xIfb7X51XrHE2DnfDNQlJuT4/Qw2iVhIizKSQWKO1OwEPEwP+06ytMLdnDolHlP/IqOTXjmqg40DdV98Fotbbd5VWLnAji6vWy5xcvs5Gh/FbS/EoKbOK9G8QiVChbTpk1j2rRp7N+/H4COHTvy9NNPM2rUKEfUJiKVcCTDnDDshx3mhGHN69Xh+XEdGd4+3MmViVMYBhzdYQaJnQsgfU/ZOi8fiBoEHcaZVyeCGjqvTvE4lQoWERERvPzyy8TExADw8ccfM27cODZt2kTHjh0dUqCIXFix1cZHa/fz2rK95JRMGPbngVE8MLwNgX66KFmrnB6w6nSYOJFYts7LF6KHmWGi3SgIDHNameLZLIZhGBff7PzCwsL497//zZ///OcKbZ+ZmUloaCgZGRmEhIRU56NFar3NKad4Yu42dh7JBKBny/q8dHUnYpvo71atYbPBoY3mMxO7FsKp5LJ13v7Q5jLzNke7KyAg1Gllivur6O/vKv/vjNVq5euvvyYnJ4f+/fufd7uCggIKCgrKFSYi1ZOZX8S/l+zhs98OYBgQWseXx0bFMrFXpCYMqw1sVkj5zbwqsetbyDxUts430JwttMNV5n81W6jUsEoHi23bttG/f3/y8/OpW7cu8+bNo0OHDufdfurUqTz33HPVKlJETIZh8O3WI7ywaCfHsszAfk335jwxpj0NNWGYZzs9lPbpMJF9tGydXzC0vdy8zREzAvwCnVen1HqVvhVSWFhIcnIyp06dYs6cOXz44YfExcWdN1yc64pFZGSkboWIVNL+9Bz+uWA7q/eZE4a1bhjEi+M7cUmMHrzzWNYiSFplhondiyD3eNk6/1Bz6vEO46D1UPANcF6dUitU9FZItZ+xGDFiBNHR0bz33nt2LUxETAXFVt6PS+R/K+IpLLbh5+PFlCEx3DVEE4Z5pOICSFxZEiYWQ/6psnV16ptdHB3Gm10dPn5OKlJqI4c/Y3GaYRjlrkiIiP38knCcp+ZvI+GYOWHYgJiGvDC+E1ENNWGYRynKM+fj2LnAnJ+j4Ixn0YIameNLtL/KHG/CW8Owi2urVLB44oknGDVqFJGRkWRlZTFr1ixWrlzJkiVLHFWfSK10PLuAf323mzm/m1NSN6zrxz/HduCqrs00YZinKMyBfUvN0S/3/gBFOWXrgpuaQaLDVdCiP3jpypS4j0oFi6NHj3LLLbdw5MgRQkND6dKlC0uWLOGyyy5zVH0itYrNZvD1xhSmfr+bU7nmhGE39WnBo5fHEhqo/1N1e/mZJWFiPuz7sfyMoaGRJWFiHET0Bi8NvS7uqVLBYvr06Y6qQ6TW23s0iyfnbWP9/pMAtG8awktXd6JHi/pOrkyqJe8k7Fli3uZI+AmshWXr6rcyg0SHcdCshyb5Eo+gYflEnCyv0Mqby/fxwapEim0GgX7ePDSiLbdd2koThrmrnOOwZ7EZJhLjwFZUtq5BjPnwZYdxmjFUPJKChYgTrdidxj8XbOfgSfOS+GUdwnn2qo40r6cJw9xOdpo5vsTOBbB/DRjWsnWNO5hBov1V0Li9woR4NAULESdIzcjn+UU7+G5bKgDNQgN49qqOjOyomSXdSubhkjCxEA78DJzRvd+ki/nwZftx0Kit00oUqWkKFiI1yGoz+OSX/fxn6V6yC4rx9rJw+6WteHBEW4L89dfRLZxKKZt+POW38uua9Sh5ZuIqCGvtnPpEnEz/konUkK0HT/HEvG1sP2SOUdC9RT1eGt+ZDs00UJzLO5FoXpXYuQAO/15+XWTfktscV0K9Fs6pT8SFKFiIOFhmfhH/+WEPn/xqThgWEuDDP0bFcmPvFpowzJWl7zPbQncugNRtZ6ywQMtLS8LEWAhp5qwKRVySgoWIgxiGwXfbUnnu2x2klUwYNq5bM54a04FGwZowzOUYBqTtKpnkayGk7SxbZ/GGqIHmw5exYyE43Hl1irg4BQsRB0g+nsvTC7ezcs8xAFo1COTF8Z0Z0EYThrkUw4DUrWW3OY7vK1vn5Quth5jPS7QbA0ENnFamiDtRsBCxo8JiGx+sTuTNn/ZRUGzDz9uLu4dEc/eQaAJ8NSyzSzAM8zmJnQvMr5P7y9Z5+0H0cPM2R7srzEm/RKRSFCxE7GRd0gmenLeNfWnZAFwS3YAXxnciulFdJ1cm2GxwcH3ZbY6MlLJ1PnWgzQhz0Ko2IyFAD9OKVIeChUg1ncgpZOp3u/h6ozlhWIMgP54a257x3ZprwjBnslkh+ZeSMPEtZB0pW+cbBG0vN69MtLkM/DRbrIi9KFiIVJFhGHyz8SD/+m4XJ3PNIZtv7BPJP66IpV6gn5Orq6WsxbB/tRkmdi+CnGNl6/xDoN0o8wHMmOHgq9FNRRxBwUKkCuLTsnhi3nbWJZ0AILZJMC9d3YmeLcOcXFktVFwISXElYWIx5J0oWxdQD2LHmFcmWg8BH3XjiDiagoVIJeQXWfnf8n28vyqRIqtBHV9vHhzRhtsHROGrCcNqTlE+JK4oCRPfQUFG2brABmZLaIdxEDUIvDXdvEhNUrAQqaCVe9J4esEOkk/kAjA8tjHPjetIRP1AJ1dWSxTmQvyPZpjYuwQKs8vW1Q03R77sMA5aXALe+qdNxFn0t0/kItIy83lu0U4WbzUf/msSYk4YdnnHcD2c6WgF2bDvBzNM7FsGRbll64KblczLMQ4i+4CX2nlFXIGChch5WG0Gn/92gH8v2UNWQTFeFrjt0igeuqwtdTVhmOPkZ8CeJWZbaPyPUJxftq5eC/Phyw7joXlP8NLtJxFXo38dRc5h+6EMnpy3jS0HzXv3XSNCeenqznRqHurkyjxU7gnY8715ZSJxBVgLy9aFtS67MtG0G+gqkYhLU7AQOUN2QTH/WbqHj9fux2ZAsL8Pj17Rjpv6tsRbE4bZV0662RK6cwEkrQJbcdm6hu3KwkR4R4UJETeiYCGCOSbFku2pPPftTlIzzUvvV3Ztxj/HtKdxSICTq/MgWanmYFU7F8CBn8Gwla0L71QyY+hV0DjWeTWKSLUoWEitl3Iil2cW7mD57jQAWoQF8sL4Tgxu28jJlXmIjINlYSL5V8AoW9e0W9mViQbRzqpQROxIwUJqrSKrjQ9XJ/HGT3vJL7Lh623hrsHRTBkaownDquvkfnPG0F0LzTk6zhTRu+QBzKugfitnVCciDqRgIbXShv0neGLeNvYeNcdC6BsVxktXdyamsSYMq7KsVNj8hXll4sjmM1ZYoEW/ktscV0JohLMqFJEaoGAhtcqp3EJe/n43s9abs1uGBfnxxOj2XNtDE4ZViWFAyjpY954ZKE4/gGnxgpaXloWJ4CbOrVNEaoyChdQKhmEw9/dDvPTdLk7kmK2ME3tF8tioWOoHacKwSivKg23fwLr3IXVr2fLIvtDtJmg3BurqGRWR2kjBQjxefFo2T83fxq+J5uRUbRrX5aWrO9MnShOGVdrJA7BhOvz+CeSdNJf5BEDnCdD7DmjWzanliYjzKViIx8ovsvLOinimxSVQZDUI8PXi/uFt+MuA1vj5aMTGCjMMc9CqdR+Yg1id7uqo1wJ6/wW63wKBCmkiYlKwEI+0et8xnpq/nQPHzbklhrZrxPPjOhEZpgnDKiw/E7Z8aQaK4/vKlrceCn3+Cm0v1/wcInIWBQvxKGlZ+by4aBcLtxwGIDzEn2eu7MioTk30cGZFHdtjhoktX5bNIOoXbD470fsv0Kitc+sTEZemYCEewWYz+HxdMv+3ZDdZ+eaEYX/q34q/jWxLcICvs8tzfTareZtj3fuQFFe2vGE76HMHdL0B/IOdV5+IuA0FC3F7Ow5n8OS87WxOOQVAl4hQXhrfmc4RmjDsonKOw6ZPYP10yDBbcLF4QbvRZqCIGqx5OkSkUhQsxG3lFBTz32V7mflzEjYD6vr78Mjl7ZjUTxOGXdThTebtjm3fgLXAXFYnDHpOhl63mw9miohUgYKFuKUfdqTy7MIdHMkwJwwb07kpT1/ZgXBNGHZ+xYXmIFbr3oeD68qWN+0Kfe6ETteAbx3n1SciHkHBQtzKwZO5PLtwJz/uOgpAZFgdnh/XiaHtGju5MheWeRg2zISNH0GOOdEaXr7Q8WqzuyOil253iIjdKFiIWyiy2pj5cxKvLdtHXpEVX28Lfx3UmnuHtqGOn1oez2IYkPyLeXVi17dlQ20HNzVvdfSYDMHhzq1RRDySgoW4vI0HTvLkvG3sTs0CoE+rMF66uhNtwtWlcJbCHNj2tfn8xNHtZctbXmo+jBk7FrzVJSMijqNgIS4rI7eIl5fs5st1yQDUD/Tl8dHtua5nhMak+KMTiWZnx6ZPIT/DXOZTB7pcbwaKJp2dW5+I1BoKFuJyDMNgwebDvLh4J+nZ5oRhE3pG8MTo9oRpwrAyNhskLDdvd+xbSulQ2/VbmfN2dL8Z6tR3ZoUiUgspWIhLSTyWzT8XbOfn+OMAxDSuy4vjO9GvdQMnV+ZC8jNg0+ew/gPzSsVpMSPMhzFjRmiobRFxGgULcQmGYfD2inje/CmeQqsNfx9zwrA7BmrCsFJHd5phYstsKMoxl/mHQPdJ5lDbDaKdW5+ICAoW4iI++y2ZV5fuBWBQ20a8MK4jLRsEObkqF2Athj2LzYcx968uW96oPfT9K3S+HvzrOq8+EZE/ULAQp4tPy+LFRTsBeOTydtwzJFoPZ2Yfg98/hg0zIPOQucziDbFjzNsdrQZo7AkRcUkKFuJUBcVW7vtyMwXFNga1bcTdg2t5qDi0EX57H3bMBav54CqBDaDnreb4E6ERTi1PRORiFCzEqf69ZA+7jmQSFuTHqxO64FUb5/goLoAd88zujkMby5Y372lenegwHnw1VLmIuAcFC3Ga1fuO8eGaJAD+PaELjWvbPB8ZB81bHRs/htx0c5m3H3S8pmSo7Z7OrU9EpAoULMQpTuQU8revtgBwS7+WDG9fS4aXNgzYv8a8OrF7MRhWc3lI87Khtus2cm6NIiLVoGAhNc4wDB79ZitpWQXENK7LE6PbO7skxyvIhq2zze6OY7vKlrcaaF6daDcavPXXUUTcn/4lkxr3xbpkftx1FD9vL968obtnTyJ2PAHWf2gOaFVQMtS2byB0vcEcHTO8g3PrExGxMwULqVHxaVm8UNJa+ugV7ejQLMTJFTmAzQbxy8zbHfE/li0Pa22GiW43QZ16TitPRMSRFCykxhQUW7n/y83kF9kY2KYht18a5eyS7CvvZNlQ2yf3lyy0QJuR5u2O6GHgpVFERcSzKVhIjXn1hz3sLGkt/c91XT2ntTR1u3l1YutXUJxnLgsIhe63QO8/m1cqRERqCQULqRFr9qXzwWqztfSVaz2gtdRaBLsXmYNZJa8tW964Y9lQ236BzqtPRMRJFCzE4U7kFPLwV5sBuLlvCy7r4MatpVlHy4bazjpiLrN4Q4erzNsdLfprqG0RqdUULMShDMPgH3PKWkufGuOGXRCGAQc3mLc7dswDW5G5PKgR9LwNet0GIc2cW6OIiItQsBCH+nJdCst2HsXX28IbN3Rzr9bSojzYPtcMFEc2ly2P6A197jSvUvj4O608ERFXpGAhDhOfls3zi3YA8OjlsXRsFurkiiroVDKsnw6/fwJ5J8xl3v7QeQL0/gs07+Hc+kREXJiChThEQbGVB2ZtKm0t/fMAF28tNQxIijNHxtzzHRg2c3lopNnZ0f1PENTAuTWKiLgBBQtxiP8u3cuOw5nUD/TlVVduLS3Igi2zzECRvqdsedTgkqG2R4GXG92+ERFxskoFi6lTpzJ37lx2795NnTp1uOSSS3jllVdo166do+oTN/RzfDrvrUoEzNbScFdsLU3fZ4aJzV9AYZa5zK8udL3RvN3RONa59YmIuKlKBYu4uDimTJlC7969KS4u5sknn2TkyJHs3LmToKAgR9UobuTkGa2lN/VtwciOTZxb0JlsVti3FH57DxJXlC1vEGNeneh6IwR44BDjIiI1qFLBYsmSJeW+nzlzJo0bN2bjxo0MGjTonPsUFBRQUFBQ+n1mZmYVyhR3YBgGj83dytHMAlo3CuKfrtJamnsCNn1qTgZ2KrlkocW8zdHnDogaoqG2RUTspFrPWGRkmLM1hoWFnXebqVOn8txzz1XnY8RNzFqfwg87zNZSl5i19MgWs1V02zdQnG8uC6gHPf5kPpBZv5UzqxMR8UgWwzCMquxoGAbjxo3j5MmTrF69+rzbneuKRWRkJBkZGYSE6LKzp0g4ls3YN9eQV2TlidGx/HVQtHMKKS6EXQvNQJHyW9nyJp3NsSc6XauhtkVEqiAzM5PQ0NCL/v6u8hWLe++9l61bt7JmzZoLbufv74+/vwYR8mSFxTYemLWJvCIrl8Y04C8DnDDpVuYR2PgRbJwJ2UfNZV4+0GG8+fxEZB8NtS0iUgOqFCzuu+8+Fi5cyKpVq4iIiLB3TeJm/rNsD9sPZVIv0Jf/XNet5lpLDcO8KrHufdi5AGzF5vK64dDrduh5KwS70MOjIiK1QKWChWEY3HfffcybN4+VK1cSFeXigx6Jw62NT+f9M1pLm4TWQGtpYS5s/8YMFKnbypZH9jNnFo29Enz8HF+HiIicpVLBYsqUKXzxxRcsWLCA4OBgUlNTAQgNDaVOnToOKVBcl9laugXDgBv7tOByR7eWntxvDrW96VPIO2ku8wmAzteZ3R1Nuzr280VE5KIq9fCm5Tz3qGfOnMmtt95aofeo6MMf4toMw+Duz35nyY5UWjcKYtF9Awj0c8BArjabOebEug9g7xKg5I9rvRbQ+w7oPgkCz9+VJCIi9uGQhzer2EAiHmj2+hSW7EgtbS21e6jIz4QtX5q3O47Hly2PHmY+jNlmpIbaFhFxQZorRCot8Vg2z327E4C/j2xHp+Z2nLU0bTes/8Ccv6Mw21zmFwzdbzaH2m7Yxn6fJSIidqdgIZVitpZuJq/IyiXRDbhjoJ1aS/MzYc5fYN8PZcsatjOfneh6A/gH2+dzRETEoRQspFL+u2wv2w5lUC/Ql/9eb6fW0uJCmD3JnLbc4gXtRpu3O6IGaewJERE3o2AhFbY2IZ33ViUA8PI1dmotNQxYeK8ZKvzqwuSF0Lxn9d9XREScQjMvSYWcyi3k4dlma+kNvSO5opOdWkt/eg62zgaLN1z/sUKFiIibU7CQizIMg8fnbiM1M5/WDYN4+ko7zVq67gNY85r5+qr/QcwI+7yviIg4jYKFXNTXGw7y/XaztfQNe7WW7loE3z1ivh76pNn1ISIibk/BQi4oKT2HZ7/dAcDfRrajc4QdWktT1sGcPwMG9JgMgx6p/nuKiIhLULCQ8zo9a2luoZX+rRvwV3u0lqbHwxcToTgf2l4BY/6rzg8REQ+iYCHn9dqPe9l6MIPQOr78d2LX6reWZh2Fz66BvBPQrAdMmAHeakwSEfEkChZyTr8kHOfduNOtpZ1pGlrNSeYKsuGL6+HUAagfBTd9BX5BdqhURERciYKFnOVUbiEPf7UZw4CJvSIZ1blp9d7QWgRfT4YjmyGwAUyaA3Ub2aVWERFxLQoWUo5hGDwxbxtHMvKJskdrqWHAogch/kfwqWNeqWgQbZdaRUTE9ShYSDlfbzzId9tS8fGy8MYN3Qjyr+YzECtfhk2fmUN1X/cRRPSyS50iIuKaFCykVFJ6Ds8uNFtLHx7Zli4R9ar3hhs/hriXzddj/gPtrqje+4mIiMtTsBAAiqw2HixpLe3XOow7B1XzdsXepbDoIfP1wL9Dr9urX6SIiLg8BQsB4PUf97LldGvp9d3wrk5r6aGN5sOahhW63gTDnrJfoSIi4tIULIRfE4/zzkqztXTqNZ1pVq8araUnEuHz66EoF6KHwVVvagAsEZFaRMGilsvILeKh2WZr6fW9IhhdndbSnHT4bALkpkOTLnD9J+Dta79iRUTE5SlY1GJntpa2ahDIM1d2rPqbFeaaQ3WfSIDQFnDz1+AfbL9iRUTELShY1GLfbDzI4m1HSlpLu1e9tdRabE4qdmgDBNSDSd9AcBO71ioiIu5BwaKW2p+ewzMlraUPXdaWrpH1qvZGhgHfPwJ7vgNvf7hpNjRqZ79CRUTErShY1EJFVhsPzN5MbqGVvlFh3DW4Gq2la/4LG2YAFrj2Q2jRz251ioiI+1GwqIXe+HEfW1JOERLgw2sTq9FaumUW/PS8+XrUK9DhKvsVKSIibknBopb5LfE4b6+MB2DqNV2q3lqasBwWTDFfX3I/9L3TThWKiIg7U7CoRTLyylpLr+sZwZguVWwtPbIVZv8JbMXQ6VoY8Zx9CxUREbelYFFLGIbBk/O2cTgjn5YNAnnmqiq2lp5Khs8nQGEWtBoI46eBl/4YiYiISb8Raok5vx9i0day1tK6VWktzT1hDoCVfRQad4CJn4GPv/2LFRERt6VgUQscOJ7DMwu2A2ZrabeqtJYW5cOsmyF9D4Q0h5u/gTpVeB8REfFoChYershq44FZm8kptNKnqq2lNhvM+yskrwX/UDNUhDa3f7EiIuL2FCw83Js/7WNzdVpLDQN+eAJ2LgBvP7jhMwjv4JhiRUTE7SlYeLB1SSd4e4XZWvqvazrTvCqtpb+8Bb9NM1+PnwZRg+xYoYiIeBoFCw91urXUZsC1PSIY26VZ5d9k2zew9Cnz9WUvQOcJ9i1SREQ8joKFBzIMg6fmb+fQqTxaNgjkuXFVaC1NWg3z7zZf970LLrnPvkWKiIhHUrDwQPM2HeLbLYfx9rLw+sRulW8tPbrT7ACxFkL7q+Dyf4GlisN+i4hIraJg4WGSj+fy9AJz1tIHh7ehe4v6lXuDjEPmAFgFGdCiP1zzAXh5O6BSERHxRAoWHsSctXQT2QXF9GkVxj1DYyr3BvkZ8Pl1kHkIGraFG74A3wDHFCsiIh5JwcKD/G95PJuSTxEc4MN/J3atXGtpcaF5+yNtB9QNN8eqCAxzXLEiIuKRFCw8xPr9J3hr+T4AXrq6MxH1Ayu+s80GC+6B/avBr64ZKuq3dFClIiLiyRQsPEBGXhEPzjJbS6/p0ZyrulaytfSnZ2Hb1+DlAxM/haZdHFKniIh4PgULN2cYBv8saS1tERbIc5WdtfS39+HnN8zXV70F0cPsX6SIiNQaChZubv7mQyw83Vp6QzeCA3wrvvPOhfD9o+brYU9BtxsdU6SIiNQaChZuLPl4Lv+cb7aWPjC8DT0q01qa/CvMvQMwoOdtMPDvjilSRERqFQULN1VstfFgSWtp71b1mVKZ1tJje+HLG6A4H9qOgtGvagAsERGxCwULN/W/5fH8nnyKYP9KzlqadRQ+vxbyTkLzXjBhBnhXcmROERGR81CwcEMb9p/gfyWtpS9e3aniraUFWeaomqeSIaw13DQb/CrRlioiInIRChZuJjO/iAdOt5Z2b864bs0rtqO1CL6aDKlbIbAhTJoDQQ0dW6yIiNQ6ChZu5umS1tLIsDoVn7XUMGDh/ZDwE/gGws1fmVcsRERE7EzBwo3M33SI+ZtPz1raveKtpSv+BVu+AIs3XPcRNO/p0DpFRKT2UrBwEykncnlq/nYA7h/Whp4tK9haumEmrPo/8/XY16Dt5Q6qUERERMHCLZitpZvJLiimV8v6TBkaXbEd9yyBxQ+brwf/A3pOdlyRIiIiKFi4hbdWxLPxwMnS1lIf7wqctoMb4ZvbwLBBt0kw5HHHFyoiIrWegoWL23jgBG/+VNZaGhlWgfbQ4wnwxXVQlAsxI+DK1zUAloiI1AgFCxeWdUZr6dUVbS3NSYfProXc49C0K1z3MXhXYv4QERGRalCwcGFPL9jBwZN5RNSvYGtpYQ58cT2cTIJ6LeCmr8G/ruMLFRERKaFg4aIWbD7EvE2H8Pay8MYN3Qi5WGuptRi+uR0ObYQ69WHSXAgOr5liRURESihYuKCUE7k8Nc9sLb1vWAw9W4ZdeAfDgO/+BnuXgE8A3DgbGrapgUpFRETKq3SwWLVqFVdeeSXNmjXDYrEwf/58B5RVexVbbTw0ezNZBcX0bFmfeysya+mqV2HjR4AFrv0QWvR1dJkiIiLnVOlgkZOTQ9euXXnrrbccUU+t9/aKBDaUtJa+XpHW0k2fw4oXzdej/w3tr3R8kSIiIudR6fmyR40axahRoxxRS6238cBJ3iyZtfSF8RVoLY3/Cb6933x96YPQ5w7HFigiInIRlQ4WlVVQUEBBQUHp95mZmY7+SLeUlV/Eg7M3YbUZjOvWjPHdL9JaengzfPUnsBVD5+th+DM1UqeIiMiFOPzhzalTpxIaGlr6FRkZ6eiPdEvPLNhBygmztfSF8Z0uvPHJA2ZbaWE2RA2GcW+Dl57DFRER53P4b6PHH3+cjIyM0q+UlBRHf6TbWbD5EHM3HcLLAq9PvEhrae4JcwCs7KMQ3gkmfgo+fjVXrIiIyAU4/FaIv78//v7+jv4Yt3XwZNmspfcOa0OvVhdoLS3Kgy9vgOP7ICQCbv4aAkJrqFIREZGL0/VzJyptLc0vpkeLetw/7AKtpTYrzL0DUn4zw8SkbyCkWc0VKyIiUgGVvmKRnZ1NfHx86fdJSUls3ryZsLAwWrRoYdfiPN20lQms33+Suv4+vD6x+/lbSw0DljwOu74Fbz+44Qto3L5mixUREamASgeLDRs2MHTo0NLvH374YQAmT57MRx99ZLfCPN3vySd5vWTW0ufHdaRFgwu0lq59E9a9Z76++j1oNaAGKhQREam8SgeLIUOGYBiGI2qpNbLyi3hw1masNoOrujbj6gu1lm79GpY9bb4e+RJ0uqZmihQREakCPWPhBM8s3EHyiVya1zNbSy0Wy7k3TIyD+Xebr/vdA5fcW3NFioiIVIGCRQ1buOUwc38vaS29oRuhdc7TWnp0B8yeBLYi6DDevFohIiLi4hQsatDBk7k8OW8bAPcOjaH3+VpLMw7CZxOgIBNaXmo+V6EBsERExA3ot1UNsdoMHp69haz8Yrq3qMf9w88zrXneKTNUZB2Ghu3ghs/BN6BGaxUREakqBYsaMm1lPOv2nyDIz/v8s5YWF8Csm+HYLqjbBCbNgTr1a75YERGRKlKwqAGbkk/y2o+nW0s70bJB0Nkb2Www7y44sAb8gs0BsOppXhUREXEvChYOll1QzAMlraVXdm3GNT3O01r649OwYy54+ZjzfzTpXLOFioiI2IGChYM9e0Zr6Yvnay399V1Y+z/z9bh3IHro2duIiIi4AQULB1q09TDfbDyIlwVem3ie1tKdC2DJY+br4c9A14k1W6SIiIgdKVg4yKFTeTw+12wtnTI0hj5R52gtPbAW5twBGNDrzzDgoZotUkRExM4ULBzAajNKZy3tFnme1tJje+DLG8FaAO3GwOh/w/lG4BQREXETChYO8G5cAuuSzNbSN27ohu8fW0szj8Bn10L+KYjoDdd+CF7eTqlVRETEnhQs7GxzyileW7YXgOfO1VqanwlfXAcZKRAWDTfOBr8LzGwqIiLiRhQs7CinoJgHZm2i2GYwtktTrv1ja2lxIXz1J0jdBkGNzAGwgho4p1gREREHULCwo2cX7uDA8VyahQbw0vjO5VtLDQO+vR8SV4BvENz0FYRFOa9YERERB1CwsJPFW4/w9ZmtpYF/aC1d/gJs+RIs3nD9x9C8h3MKFRERcSAFCzs4fCqPx+duBeCeITH0bf2H2xvrp8Pq/5ivr3wD2lxWwxWKiIjUDAWLajrdWpqZX0zXyHo8MOIPraW7v4Pv/m6+HvI49Lil5osUERGpIQoW1fRuXAK/nW4tnfiH1tKU9fDN7WDYoPstMPgfzitURESkBihYVMOWM1pLn72qI60antFaejwBvpwIxXkQcxmMfU0DYImIiMdTsKiinIJiHpy9mWKbwZjOTZnQM6JsZXYafHYN5B6Hpt3guo/A+xzzhIiIiHgYBYsqeu7bHSSl59AsNIB/XX1Ga2lBNnxxPZzcD/Vbwc1fg39dZ5YqIiJSYxQsquC7bUf4asNBLBb475mtpdZi+OY2OLwJ6oTBzXOgbmPnFisiIlKDFCwq6fAZs5bePTiafqdbSw0DFj8E+5aCTx1zAKyGMU6sVEREpOYpWFSC1Wbw8FebycgromtEKA9d1rZsZdz/we+fgMULJsyAyN7OK1RERMRJFCwq4f1VifyaeIJAP29ev6F7WWvp75/Cyn+Zr0f/G2JHO69IERERJ1KwqKCtB0/xn6V7AHj2yo5EnW4t3bcMvn3AfD3gYej9FydVKCIi4nwKFhVgzlpqtpaO7tyE63qVtJYe3gRfTQbDCl1ugOFPO7dQERERJ1OwqIAXFu0kKT2Hpme2lp7cD59fD0U50HoIXPU/DYAlIiK1noLFRXy/7Qiz1qeYraXXd6NeoB/kHIfProWcNAjvDNd/Cj5+zi5VRETE6RQsLuBIRh6PlbSW3jU4mv7RDaAwF768AY7HQ2ikOQBWQIiTKxUREXENChbnYbUZPDx7Cxl5RXSJCOWhEW3BZoW5d8DBdRAQCjd/AyFNnV2qiIiIy1CwOI8PVifyS+Jx6vh68/rEbvh5W+D7R2H3IvD2hxtnQeNYZ5cpIiLiUhQszmHbwQxe/aGktfSqDrRuVBd+fh3WfwhY4Jr3oeUlTq1RRETEFSlY/EFuYTEPzNpEsc1gVKcmXN8rErZ+BT8+a25wxVToON6ZJYqIiLgsBYs/eGHRThLTc2gSEsDUazpjSYqD+feYK/vfC/3udm6BIiIiLkzB4gxLtqfy5bqS1tKJXamXuQdmTQJbEXS8Bi57wdklioiIuDQFixKpGfk8NncrAHcOiuaSBnnw+XVQmAUtB8DV74KXflwiIiIX4uPsAlyBrWTW0lO5RXRuHsrDAxrDJ6Mg6wg0ag83fA4+/s4uU0RExOXpf8ExW0vXJpitpW9MiMXvm0lwbDcEN4VJ30Cdes4uUURExC3U+mCx/VAGr5bMWvrM2Fhar/k7HPgZ/EPMAbBCI5xcoYiIiPuo1cEit7CY+2dtoshqcHnHcCaefA92zAMvX5j4GTTp5OwSRURE3EqtDhYvLNpF4jGztfS1Fmux/PqOuWL8NGg92LnFiYiIuKFaGyzM1tJkLBb4uO8hAlc8ba4Y8Rx0uc65xYmIiLipWhksjmaWtZa+0DWDdmv/BhjQ569w6QPOLU5ERMSN1bpgcWZr6ejwk9y8/zGwFkLsWLjiZbBYnF2iiIiI26p1weLDNYn8HH+clr4ZvGF9CUt+BkT2hWs/BC9vZ5cnIiLi1mpVsNh+KIN//7CHYHKZX+81fLMPQ4M25hTovnWcXZ6IiIjbqzXBIq/QygOzNoG1iK/qv0P9rL0Q1NgcACswzNnliYiIeIRaM6T3C4t3knAsm2mBH9I+73fwqws3fw31Wzm7NBEREY9RK65YLN2Ryhe/JfOoz2xG2VaBxRuu/xiadXN2aSIiIh7F44PF0cx8/jFnK5O8l3GPz0Jz4VX/g5gRzi1MRETEA3l0sLDZDP721RZ656/led+PzIVDn4TuNzu1LhEREU/l0c9YzPg5idyEtXzo9xZeGNBjMgx6xNlliYiIeCyPDRY7Dmfw1ZIVzPJ7lQBLEbS9Asb8VwNgiYiIOJBH3grJK7TyzOfL+dB7KmGWbIxmPWDCDPD22BwlIiLiEjwyWPz72w08nfUcLbyOYa3XCstNX4FfkLPLEhER8XgeFyx+3H6QQZsfoYtXEkX+9fG+ZS7UbeTsskRERGqFKgWLd955h6ioKAICAujZsyerV6+2d11VkpaRR/ac+xjivYUiL398b/kGGkQ7uywREZFao9LBYvbs2Tz44IM8+eSTbNq0iYEDBzJq1CiSk5MdUV+F2WwGa6b/nfHGcqx4wXUfQUQvp9YkIiJS21gMwzAqs0Pfvn3p0aMH06ZNK13Wvn17xo8fz9SpU8/avqCggIKCgtLvMzMziYyMJCMjg5CQkGqUXt7qWa8ycPcLAKQNfpnGQ++223uLiIjUdpmZmYSGhl7093elrlgUFhayceNGRo4cWW75yJEjWbt27Tn3mTp1KqGhoaVfkZGRlfnICjl2KJG+u8xQs631HQoVIiIiTlKpYJGeno7VaiU8PLzc8vDwcFJTU8+5z+OPP05GRkbpV0pKStWrPY9GzVuza8Ab/Fr/KjpN+j+7v7+IiIhUTJUGdrD8YZApwzDOWnaav78//v7+VfmYSul62SS4bJLDP0dERETOr1JXLBo2bIi3t/dZVyfS0tLOuoohIiIitU+lgoWfnx89e/Zk2bJl5ZYvW7aMSy65xK6FiYiIiPup9K2Qhx9+mFtuuYVevXrRv39/3n//fZKTk7nrrrscUZ+IiIi4kUoHi4kTJ3L8+HGef/55jhw5QqdOnfjuu+9o2bKlI+oTERERN1LpcSyqq6J9sCIiIuI6HDKOhYiIiMiFKFiIiIiI3ShYiIiIiN0oWIiIiIjdKFiIiIiI3ShYiIiIiN0oWIiIiIjdKFiIiIiI3VRpdtPqOD0eV2ZmZk1/tIiIiFTR6d/bFxtXs8aDRVZWFgCRkZE1/dEiIiJSTVlZWYSGhp53fY0P6W2z2Th8+DDBwcFYLBa7vW9mZiaRkZGkpKR47FDhnn6MOj735+nHqONzf55+jI48PsMwyMrKolmzZnh5nf9Jihq/YuHl5UVERITD3j8kJMQj/7CcydOPUcfn/jz9GHV87s/Tj9FRx3ehKxWn6eFNERERsRsFCxEREbEbjwkW/v7+PPPMM/j7+zu7FIfx9GPU8bk/Tz9GHZ/78/RjdIXjq/GHN0VERMRzecwVCxEREXE+BQsRERGxGwULERERsRsFCxEREbEbBQsRERGxG7cKFu+88w5RUVEEBATQs2dPVq9efcHt4+Li6NmzJwEBAbRu3Zp33323hiqtmsoc38qVK7FYLGd97d69uwYrrrhVq1Zx5ZVX0qxZMywWC/Pnz7/oPu52/ip7jO52DqdOnUrv3r0JDg6mcePGjB8/nj179lx0P3c5j1U5Pnc6h9OmTaNLly6lIzL279+f77///oL7uMu5O62yx+hO5+9cpk6disVi4cEHH7zgdjV9Ht0mWMyePZsHH3yQJ598kk2bNjFw4EBGjRpFcnLyObdPSkpi9OjRDBw4kE2bNvHEE09w//33M2fOnBquvGIqe3yn7dmzhyNHjpR+tWnTpoYqrpycnBy6du3KW2+9VaHt3e38QeWP8TR3OYdxcXFMmTKFX3/9lWXLllFcXMzIkSPJyck57z7udB6rcnynucM5jIiI4OWXX2bDhg1s2LCBYcOGMW7cOHbs2HHO7d3p3J1W2WM8zR3O3x+tX7+e999/ny5dulxwO6ecR8NN9OnTx7jrrrvKLYuNjTUee+yxc27/6KOPGrGxseWW3XnnnUa/fv0cVmN1VPb4VqxYYQDGyZMna6A6+wKMefPmXXAbdzt/f1SRY3Tnc2gYhpGWlmYARlxc3Hm3cefzWJHjc/dzWL9+fePDDz885zp3PndnutAxuuv5y8rKMtq0aWMsW7bMGDx4sPHAAw+cd1tnnEe3uGJRWFjIxo0bGTlyZLnlI0eOZO3atefc55dffjlr+8svv5wNGzZQVFTksFqroirHd1r37t1p2rQpw4cPZ8WKFY4ss0a50/mrLnc9hxkZGQCEhYWddxt3Po8VOb7T3O0cWq1WZs2aRU5ODv379z/nNu587qBix3iau52/KVOmMGbMGEaMGHHRbZ1xHt0iWKSnp2O1WgkPDy+3PDw8nNTU1HPuk5qaes7ti4uLSU9Pd1itVVGV42vatCnvv/8+c+bMYe7cubRr147hw4ezatWqmijZ4dzp/FWVO59DwzB4+OGHGTBgAJ06dTrvdu56Hit6fO52Drdt20bdunXx9/fnrrvuYt68eXTo0OGc27rruavMMbrb+QOYNWsWv//+O1OnTq3Q9s44jzU+bXp1WCyWct8bhnHWsottf67lrqIyx9euXTvatWtX+n3//v1JSUnh1VdfZdCgQQ6ts6a42/mrLHc+h/feey9bt25lzZo1F93WHc9jRY/P3c5hu3bt2Lx5M6dOnWLOnDlMnjyZuLi48/7idcdzV5ljdLfzl5KSwgMPPMDSpUsJCAio8H41fR7d4opFw4YN8fb2Puv/3tPS0s5KYqc1adLknNv7+PjQoEEDh9VaFVU5vnPp168f+/bts3d5TuFO58+e3OEc3nfffSxcuJAVK1YQERFxwW3d8TxW5vjOxZXPoZ+fHzExMfTq1YupU6fStWtX3njjjXNu647nDip3jOfiyudv48aNpKWl0bNnT3x8fPDx8SEuLo4333wTHx8frFbrWfs44zy6RbDw8/OjZ8+eLFu2rNzyZcuWcckll5xzn/79+5+1/dKlS+nVqxe+vr4Oq7UqqnJ857Jp0yaaNm1q7/Kcwp3Onz258jk0DIN7772XuXPnsnz5cqKioi66jzudx6oc37m48jn8I8MwKCgoOOc6dzp3F3KhYzwXVz5/w4cPZ9u2bWzevLn0q1evXtx8881s3rwZb2/vs/Zxynl02GOhdjZr1izD19fXmD59urFz507jwQcfNIKCgoz9+/cbhmEYjz32mHHLLbeUbp+YmGgEBgYaDz30kLFz505j+vTphq+vr/HNN9846xAuqLLH99prrxnz5s0z9u7da2zfvt147LHHDMCYM2eOsw7hgrKysoxNmzYZmzZtMgDjv//9r7Fp0ybjwIEDhmG4//kzjMofo7udw7vvvtsIDQ01Vq5caRw5cqT0Kzc3t3Qbdz6PVTk+dzqHjz/+uLFq1SojKSnJ2Lp1q/HEE08YXl5extKlSw3DcO9zd1plj9Gdzt/5/LErxBXOo9sEC8MwjLffftto2bKl4efnZ/To0aNcG9jkyZONwYMHl9t+5cqVRvfu3Q0/Pz+jVatWxrRp02q44sqpzPG98sorRnR0tBEQEGDUr1/fGDBggLF48WInVF0xp9u6/vg1efJkwzA84/xV9hjd7Rye69gAY+bMmaXbuPN5rMrxudM5vP3220v/fWnUqJExfPjw0l+4huHe5+60yh6jO52/8/ljsHCF82gxjJKnOERERESqyS2esRARERH3oGAhIiIidqNgISIiInajYCEiIiJ2o2AhIiIidqNgISIiInajYCEiIiJ2o2AhIiIidqNgISIiInajYCEiIiJ2o2AhIiIidvP/IgAmMmylH+IAAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "time = np.arange(ctr*dt)\n", "\n", "plt.title(\"MSD\")\n", "plt.plot(time, msd_an1, label=\"DCA\")\n", "plt.plot(time, msd_cat1, label=\"EMIM\")\n", "plt.legend();" ] }, { "cell_type": "code", "execution_count": 23, "id": "1ce03cef", "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjoAAAGxCAYAAABr1xxGAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/bCgiHAAAACXBIWXMAAA9hAAAPYQGoP6dpAABWV0lEQVR4nO3deVxVdf7H8ddlFRCuLAKSaFrkEmqKG7Zoo6ElOS2TThq5ZbaoMWmLzdSU81PbppqyzMotdXKWsmkx08o0xwVFyQXXREEFcYHLIrLd8/vD8c4gqKDg4V7ez8fj/sE5n3Pv59sR7rvvPed7LYZhGIiIiIi4IDezGxARERGpKwo6IiIi4rIUdERERMRlKeiIiIiIy1LQEREREZeloCMiIiIuS0FHREREXJaCjoiIiLgsBR0RERFxWQo6InJFzZs3D4vFgsVi4ccff6y03zAMrr32WiwWC3369HFsP3HiBJMnT6Z9+/b4+flhtVpp27YtCQkJbN26tcrnt1gsNGrUiPDwcG699VamT59OdnZ2pdd88cUXsVgsuLm5sX///kr7CwsLCQgIwGKxMGLECMf2AwcOYLFYeP311y/rv4mI1B0FHRExhb+/P7Nnz660fdWqVfzyyy/4+/s7thUUFNCzZ0/mzZvHQw89xBdffMGiRYt4+OGHSUtLIyUlpdLzzJ07l3Xr1rFixQreffddbrjhBl555RXatWvHd999V2VPjRs3Zu7cuZW2/+Mf/6C0tBRPT89LH7CImMLD7AZEpGEaMmQIixYt4t133yUgIMCxffbs2cTGxpKXl+fY9o9//IN9+/bxww8/cOutt1Z4nieffBK73V7p+aOjo+natavj53vvvZff/e533HTTTdxzzz3s3buXsLCwSj3Nnz+fl156CTe3//5/4OzZs7n77rv54osvLnvcInJlaUZHRExx//33A/DJJ584ttlsNj799FNGjRpVofbEiRMANGvWrMrn+t9QciEtWrTgz3/+M/n5+cyaNavS/lGjRpGRkcGKFSsc2/bs2cOaNWsq9SQizkFBR0RMERAQwG9+8xvmzJnj2PbJJ5/g5ubGkCFDKtTGxsYC8OCDD/L55587gs+luOOOO3B3d2f16tWV9kVFRXHzzTdX6GnOnDlcffXV9O3b95JfU0TMo6AjIqYZNWoUSUlJ7NixAzgTKu67774K1+cA3HjjjUyZMoWff/6Zu+++m5CQEFq3bs2jjz5a4ULk6vDz8yMkJIQjR46ct6d//etfnDx5kvLycj7++GNGjBiBxWK5tEGKiKkUdETENL179+aaa65hzpw5bNu2jY0bN573I6Lnn3+e9PR05syZw9ixY2ncuDHvv/8+MTExFT7+qg7DMM6777777sPLy4tFixaxdOlSsrKyKtxpJSLORRcji4hpLBYLI0eO5O233+b06dNcd9113HzzzeetDwsLY+TIkYwcORKA1atXc/vtt/PEE084rvm5mMLCQk6cOEGHDh2q3O/n58eQIUOYM2cOLVu2pF+/frRs2bLmgxORekEzOiJiqhEjRnD8+HHef/99R4CprltuuYW4uDiOHTtW5fo4Vfn6668pLy+vsEbPuUaNGkVKSgpffvmlLkIWcXKa0RERU1111VU89dRT7Nq1i+HDh1dZc/ToUZo2bVrp7qry8nL27t2Lr68vTZo0uehrpaenM2nSJKxWK2PHjj1vXWxsLKNGjcJms3H33XfXaDwiUr8o6IiI6V5++eUL7l+wYAGzZs1i6NChdOvWDavVyqFDh/joo4/YsWMHL7zwAl5eXhWO2b59O2VlZZSVlZGdnc1PP/3E3LlzcXd3Z8mSJTRt2vSCr1nVYobnowuVReovBR0RqfcGDhxIVlYWS5cuZebMmeTk5ODv70/Hjh1ZsGABDzzwQKVjzn4M5uXlRZMmTWjXrh3PPPMMDz300EVDTnWdOnUKAG9v71p5PhGpfRbjQrcfiIjIeS1ZsoR77rmHr7/+mjvuuMPsdkSkCgo6IiI19Msvv5CSksJzzz1HXl4eaWlpNGrUyOy2RKQKuutKRKSG/vSnP5GQkEBERATffPONQo5IPaYZHREREXFZmtERERERl6WgIyIiIi5LQUdERERcVoNeR8dut3PkyBH8/f214JeIiIiTMAyD/Px8IiIiKq2Yfq4GHXSOHDlCZGSk2W2IiIjIJcjIyKB58+YXrGnQQcff3x848x8qICDA5G5ERESkOvLy8oiMjHS8j19Igw46Zz+uCggIUNARERFxMtW57EQXI4uIiIjLUtARERERl6WgIyIiIi6rQV+jUx2GYVBWVkZ5ebnZrbgkT09P3N3dzW5DRERclILOBZSUlJCZmcmpU6fMbsVlWSwWmjdvTuPGjc1uRUREXJCCznnY7XbS0tJwd3cnIiICLy8vLSpYywzD4NixYxw6dIioqCjN7IiISK1T0DmPkpIS7HY7kZGR+Pr6mt2Oy2ratCkHDhygtLRUQUdERGqdLka+iIstLS2XR7NkIiJSl/QuLiIiIi5LQUdERERcloKO1Ko+ffqQmJhodhsiIiKAgo7LysrKYvz48bRu3Rpvb28iIyO58847+f77781uTURE5IrRXVcu6MCBA9x44400adKEV199lY4dO1JaWsq3337L448/zq5du8xuUUREXJRhGBw8cYotGTlsSc+l1zUhDIgON60fBZ0aMAyDotIrv0Kyj6d7je5Oeuyxx7BYLCQlJeHn5+fYfv311zNq1CgA0tPTGT9+PN9//z1ubm4MGDCAd955h7CwMABefPFFPv/8cyZOnMjzzz9PTk4Ot99+Ox9++CH+/v4AFBYW8uijj/LZZ5/h7+/PpEmTanHUIiLiDPJPl7L1kI0t6WeCzZaMXE4Wljj2nyopV9BxFkWl5bR/4dsr/rqpU/rj61W9U3Xy5EmWLVvG1KlTK4Scs5o0aYJhGNx11134+fmxatUqysrKeOyxxxgyZAg//vijo/aXX37h888/56uvviInJ4fBgwfz8ssvM3XqVACeeuopVq5cyZIlSwgPD+e5554jOTmZG264oTaGLSIi9YzdbvDLsYL/BJozwWb30XwMo2Kdl7sb0VcF0LlFILdc19ScZv9DQcfF7Nu3D8MwaNu27XlrvvvuO7Zu3UpaWhqRkZEALFiwgOuvv56NGzfSrVs34Mzq0PPmzXPM4CQkJPD9998zdepUCgoKmD17Nh9//DG33XYbAPPnz6d58+Z1PEIREblSck+VsCUj90ywSc8hJSOX/NNlleqaB/rQuUUgnSOb0LlFE9pHBODtUT8WgVXQqQEfT3dSp/Q35XWry/hPrL7QR107d+4kMjLSEXIA2rdvT5MmTdi5c6cj6Fx99dWOkAPQrFkzsrOzgTOzPSUlJcTGxjr2BwUF0aZNm2r3KiIi9UdZuZ3dR/P/E2rOBJv9xwsr1fl4utOxufVMsGnRhM6RTQgNaGRCx9WjoFMDFoul2h8hmSUqKgqLxcLOnTu56667qqwxDKPKIHTudk9Pzwr7LRYLdrvdUSsiIs7rWH7xmetqMnLZfDCHrYdsVV6H2jrEjxtaNHHM2LQN98fD3Xlu2q7f79pSY0FBQfTv3593332XCRMmVLpOJzc3l/bt25Oenk5GRoZjVic1NRWbzUa7du2q9TrXXnstnp6erF+/nhYtWgCQk5PDnj176N27d+0OSkRELktJmZ3UzDy2pOew+T+zNYdyiirV+Xt7nAk1kWeCzQ2RTQj08zKh49qjoOOC3nvvPXr16kX37t2ZMmUKHTt2pKysjBUrVjBz5kxSU1Pp2LEjw4YN46233nJcjNy7d2+6du1arddo3Lgxo0eP5qmnniI4OJiwsDB+//vfV/pusMmTJ3P48GE+/vjjuhiqiIicwzAMMm2n2ZKey+b0HLak57D9SB4lZfYKdRYLXBfqf+bjpxZN6NIikGuaNsbNzbW+g1BBxwW1atWKzZs3M3XqVCZOnEhmZiZNmzYlJiaGmTNnYrFY+Pzzzxk/fjy33HJLhdvLa+K1116joKCAQYMG4e/vz8SJE7HZbBVqMjMzSU9Pr83hiYjI/zhdWs62wzY2H8xx3A11NK+4Ul2gr6fj46cuLQPp2NyKfyPPKp7RtViMBnyxRV5eHlarFZvNRkBAQIV9p0+fJi0tjVatWtGoUf29yMrZ6b+ziEj1GYZB+slT/zNbk8vOzDzK7BXfyt3dLLRr5k/nyEC6tGxC58hAWgb71mhNtvrsQu/f59KMjoiISD1VUFzG1oz/hppzF+M7q6m/N13+8/FT5xaBdLjKio9X/bi922wKOiIiIvWA3W6w/3gBmw9efDG+668K+E+oOXPRcIS1kcvM1tQ2BR0RERETVHcxvqua+NClZf1cjM8ZKOiIiIjUsUqL8WXksP+Y8y/G5wwUdC6iAV+rfUXov6+IuKL/XYxvS/qZxfhOlbjeYnzOQEHnPM6uCnzq1Cl8fHxM7sZ1lZScuajO3V3TsCLinP53Mb6zszUZJxvGYnzOQEHnPNzd3WnSpInju518fV3ntrz6wm63c+zYMXx9ffHw0D9FEXEOR3KLHNfVbMnIZdth20UX4+vcIpBrXXAxPmegd5cLCA8PB3CEHal9bm5utGjRQiFSROqls4vxOWZr0nPJyjtdqe5/F+Pr3CKQTpENYzE+Z6CgcwEWi4VmzZoRGhpKaWmp2e24JC8vr0pfGyEiYob/XYzv7GxN6pELL8Z39qsTXGkxPldT46CzevVqXnvtNZKTk8nMzGTJkiXn/ZbssWPH8sEHH/Dmm2+SmJjo2F5cXMykSZP45JNPKCoqom/fvrz33ns0b97cUZOTk8OECRP44osvABg0aBDvvPMOTZo0cdSkp6fz+OOP88MPP+Dj48PQoUN5/fXX8fKq3c883d3ddQ2JiIiLObsY39kLhrek53LiAovxnZ2x6di8iRbjcyI1DjqFhYV06tSJkSNHcu+995637vPPP2fDhg1ERERU2peYmMiXX37J4sWLCQ4OZuLEicTHx5OcnOwIFEOHDuXQoUMsW7YMgIcffpiEhAS+/PJLAMrLyxk4cCBNmzZlzZo1nDhxguHDh2MYRo2/s0lERFybYzG+9P+uW7PnaD728yzG55itaanF+JxdjYPO7bffzu23337BmsOHDzNu3Di+/fZbBg4cWGGfzWZj9uzZLFiwgH79+gGwcOFCIiMj+e677+jfvz87d+5k2bJlrF+/nh49egDw4YcfEhsby+7du2nTpg3Lly8nNTWVjIwMR5j685//zIgRI5g6depFv/tCRERcV+6pElL+sxjf5ossxtfZ8dUJWozPFdX6NTp2u52EhASeeuoprr/++kr7k5OTKS0tJS4uzrEtIiKC6Oho1q5dS//+/Vm3bh1Wq9URcgB69uyJ1Wpl7dq1tGnThnXr1hEdHV1hxqh///4UFxeTnJzMrbfeWum1i4uLKS7+7ze65uXl1dawRUTEJGXldvYcLWBLRo7j6xOqWoyvkacbHZv/N9RoMb6GodaDziuvvIKHhwcTJkyocn9WVhZeXl4EBgZW2B4WFkZWVpajJjQ0tNKxoaGhFWrCwsIq7A8MDMTLy8tRc67p06fz0ksv1XhMIiJSfxzLLybF8UWX51+Mr1WIn+PWbi3G13DVatBJTk7mL3/5C5s3b67x55mGYVQ4pqrjL6Xmf02ePJknn3zS8XNeXh6RkZE16lNERK4cLcYnl6tWg85PP/1EdnY2LVq0cGwrLy9n4sSJvPXWWxw4cIDw8HBKSkrIycmpMKuTnZ1Nr169gDPr1xw9erTS8x87dswxixMeHs6GDRsq7M/JyaG0tLTSTM9Z3t7eeHt7X/Y4RUSkbp0sLGHSP35mzb7j1VqM75qmjXHXYnxShVoNOgkJCY4LjM/q378/CQkJjBw5EoCYmBg8PT1ZsWIFgwcPBiAzM5Pt27fz6quvAhAbG4vNZiMpKYnu3bsDsGHDBmw2myMMxcbGMnXqVDIzM2nWrBkAy5cvx9vbm5iYmNocloiIXEFFJeWMmreRlIxcoPJifB0jrQRoMT6pphoHnYKCAvbt2+f4OS0tjZSUFIKCgmjRogXBwcEV6j09PQkPD6dNmzYAWK1WRo8ezcSJEwkODiYoKIhJkybRoUMHR0hq164dAwYMYMyYMcyaNQs4c3t5fHy843ni4uJo3749CQkJvPbaa5w8eZJJkyYxZswY3XElIuKkysrtjP9kMykZuVh9PJk/qjudmlt1e7dcshpflbVp0yY6d+5M586dAXjyySfp3LkzL7zwQrWf48033+Suu+5i8ODB3Hjjjfj6+vLll19WWJRv0aJFdOjQgbi4OOLi4ujYsSMLFixw7Hd3d+frr7+mUaNG3HjjjQwePJi77rqL119/vaZDEhGResAwDJ7/13a+25mNt4cbs4d35YbIJgo5clkshmEYFy9zTXl5eVitVmw2m2aBRERM9vb3e3ljxR7cLPDesBgGRIeb3ZLUUzV5/9Z9diIiYrq/b8zgjRV7AHjp19EKOVJrFHRERMRUK3dlM3nJNgAev/UaEnq2NLkjcSUKOiIiYpqfM3J5bNFmyu0G93S5iklxbcxuSVyMgo6IiJjiwPFCRs3bSFFpObdc15RX7u2oC4+l1inoiIjIFXe8oJjhc5M4UVhC9FUBvDesC576egapA/pXJSIiV1RhcRmj5m3k4IlTRAb5MGdENxp71/pXL4oACjoiInIFlZbbefyvm9l6yEagryfzR3Yn1F/fIC51R0FHRESuCMMweO6zbfy4+xiNPN2YPaIbrZs2NrstcXEKOiIickW8uWIP/0g+hJsFZtzfhS4tAi9+kMhlUtAREZE6t2jDQd7+4cz3JE69uwP92oeZ3JE0FAo6IiJSp1akHuX5z7cD8ETfKO7v3sLkjqQhUdAREZE6k3wwh/GfbMZuwJCukST2izK7JWlgFHRERKRO/HKsgIfmb+R0qZ1b2zRl6t3RWhBQrjgFHRERqXXZ+acZPieJnFOldGpu5d1hXfDQgoBiAv2rExGRWpV/upSRczdyKKeIq4N9mT2iG75eWhBQzKGgIyIitaakzM6jCzez40gewX5ezB/VnZDG3ma3JQ2Ygo6IiNQKwzB45tOtrNl3HF8vd+aO7EbLYD+z25IGTkFHRERqxavf7mbJlsO4u1l4d1gXOjZvYnZLIgo6IiJy+eavPcDMH38B4OV7OnBrm1CTOxI5Q0FHREQuy7Ltmbz45Q4AJsVdx31dI03uSOS/FHREROSSbTxwkgmLUzAMGNajBY/feq3ZLYlUoKAjIiKXZO/RfB6av4mSMju3tQ9jyq+1IKDUPwo6IiJSY1m2MwsC2opK6dKiCW//tjPubgo5Uv8o6IiISI3knS5lxNwkjthO07qpH7OHd8PHy93stkSqpKAjIiLVVlxWztiPk9mVlU9Tf2/mj+xOoJ+X2W2JnJeCjoiIVIvdbjDpH1tZt/8Efl7uzB3RjcggX7PbErkgBR0REamW6d/s5Mufj+DhZuH9hBiir7Ka3ZLIRSnoiIjIRX30034+/CkNgFd/05Gbo5qa3JFI9SjoiIjIBX358xH+7+udADwzoC33dGluckci1aegIyIi57XulxNM/PvPAIzodTWP9G5tckciNaOgIyIiVdqVlcfDCzZRUm7n9uhwno9vrwUBxeko6IiISCVHcosYMWcj+afL6H51EG8OuUELAopTUtAREZEKbKfOLAiYlXeaqNDGfPhgVxp5akFAcU4KOiIi4nC6tJwxCzax52gBYQHezBvVHauvp9ltiVyyGged1atXc+eddxIREYHFYuHzzz937CstLeWZZ56hQ4cO+Pn5ERERwYMPPsiRI0cqPEdxcTHjx48nJCQEPz8/Bg0axKFDhyrU5OTkkJCQgNVqxWq1kpCQQG5uboWa9PR07rzzTvz8/AgJCWHChAmUlJTUdEgiIgKU2w2e/HsKSWkn8ff2YN7I7lzVxMfstkQuS42DTmFhIZ06dWLGjBmV9p06dYrNmzfz/PPPs3nzZj777DP27NnDoEGDKtQlJiayZMkSFi9ezJo1aygoKCA+Pp7y8nJHzdChQ0lJSWHZsmUsW7aMlJQUEhISHPvLy8sZOHAghYWFrFmzhsWLF/Ppp58yceLEmg5JRKTBMwyDP32VytJtWXi6W5j1YAztmgWY3ZbI5TMuA2AsWbLkgjVJSUkGYBw8eNAwDMPIzc01PD09jcWLFztqDh8+bLi5uRnLli0zDMMwUlNTDcBYv369o2bdunUGYOzatcswDMNYunSp4ebmZhw+fNhR88knnxje3t6GzWarVv82m80Aql0vIuKq3v9xn9Hyma+Mls98Zfwr5fDFDxAxUU3ev+v8Gh2bzYbFYqFJkyYAJCcnU1paSlxcnKMmIiKC6Oho1q5dC8C6deuwWq306NHDUdOzZ0+sVmuFmujoaCIiIhw1/fv3p7i4mOTk5Cp7KS4uJi8vr8JDRKSh+3zLYaZ/swuAPwxsx6BOERc5QsR51GnQOX36NM8++yxDhw4lIODMFGhWVhZeXl4EBgZWqA0LCyMrK8tRExoaWun5QkNDK9SEhYVV2B8YGIiXl5ej5lzTp093XPNjtVqJjIy87DGKiDizNXuP89Q/zywIOPqmVjx0sxYEFNdSZ0GntLSU3/72t9jtdt57772L1huGUWEhqqoWpbqUmv81efJkbDab45GRkVGdoYiIuKQdR2w8sjCZ0nKD+I7N+P0d7cxuSaTW1UnQKS0tZfDgwaSlpbFixQrHbA5AeHg4JSUl5OTkVDgmOzvbMUMTHh7O0aNHKz3vsWPHKtScO3OTk5NDaWlppZmes7y9vQkICKjwEBFpiDJOnmLE3I0UFJfRs3UQfx7cCTctCCguqNaDztmQs3fvXr777juCg4Mr7I+JicHT05MVK1Y4tmVmZrJ9+3Z69eoFQGxsLDabjaSkJEfNhg0bsNlsFWq2b99OZmamo2b58uV4e3sTExNT28MSEXEZOYUlDJ+bxLH8YtqG+zMroSveHloQUFyTR00PKCgoYN++fY6f09LSSElJISgoiIiICH7zm9+wefNmvvrqK8rLyx2zLkFBQXh5eWG1Whk9ejQTJ04kODiYoKAgJk2aRIcOHejXrx8A7dq1Y8CAAYwZM4ZZs2YB8PDDDxMfH0+bNm0AiIuLo3379iQkJPDaa69x8uRJJk2axJgxYzRTIyJyHqdLy3no403sP1ZIhLUR80Z2x+qjBQHFhdX0lq6VK1caQKXH8OHDjbS0tCr3AcbKlSsdz1FUVGSMGzfOCAoKMnx8fIz4+HgjPT29wuucOHHCGDZsmOHv72/4+/sbw4YNM3JycirUHDx40Bg4cKDh4+NjBAUFGePGjTNOnz5d7bHo9nIRaUjKyu3GmPkbjZbPfGV0+OMyY3dWntktiVySmrx/WwzDMExJWPVAXl4eVqsVm82mWSARcWmGYfD8v7azcH06Xh5uLBzdg+6tgsxuS+SS1OT9W991JSLSALz34y8sXJ+OxQJvDblBIUcaDAUdEREX949NGbz27W4A/hjfnjs6NDO5I5ErR0FHRMSF/bg7m2c/2wbA2N6tGXFjK5M7ErmyFHRERFzUtkM2Hlu0mXK7wV03RPBM/7ZmtyRyxSnoiIi4oPQTpxg5L4lTJeXcdG0Ir/5GCwJKw6SgIyLiYk4UFDN8bhLHC0po3yyAmQ90wctDf+6lYdK/fBERF3KqpIxR8zeRdryQq5r4MG9kN/wbaUFAabgUdEREXERZuZ3xf93Czxm5NPH1ZP6o7oQGNDK7LRFTKeiIiLgAwzD4w+fb+X5XNt4ebswe3pVrQxub3ZaI6RR0RERcwF++38vijRm4WeDt+zsT01ILAoqAgo6IiNNbnJTOW9/tBWDKr6Ppf324yR2J1B8KOiIiTuz7nUf5/efbARh367U80LOlyR2J1C8KOiIiTmpLeg6P//XMgoD3dmnOxLjrzG5JpN5R0BERcUJpxwsZPX8Tp0vt9L6uKS/f2wGLRQsCipxLQUdExMkcyy9m+JwkThaW0OEqK+8N64Knu/6ci1RFvxkiIk6ksLiMUfM2kn7yFC2CfJkzoht+3h5mtyVSbynoiIg4idJyO48t2sy2wzaC/LyYP6o7Tf29zW5LpF5T0BERcQKGYTD5s22s2nMMH0935ozoRqsQP7PbEqn3FHRERJzAn5fv4Z/Jh3B3s/DusM7cENnE7JZEnIKCjohIPbdg/UFmrNwHwNS7ovlV2zCTOxJxHgo6IiL12Lc7svjjv84sCJjYL4rfdm9hckcizkVBR0Sknko+eJIJn2zBbsBvu0XyRN8os1sScToKOiIi9dC+7AJGz99EcZmdX7UN5f/uitaCgCKXQEFHRKSeyc47zfA5SeSeKqVTZBNmDO2MhxYEFLkk+s0REalH8k+XMmLuRg7nFtEqxI85w7vi66UFAUUulYKOiEg9UVJm59GFm0nNzCOksRfzR3YnuLEWBBS5HAo6IiL1gN1u8PQ/f2bNvuP4erkzd0R3WgT7mt2WiNNT0BERqQde+XYXn6ccwcPNwnvDutChudXslkRcgoKOiIjJ5v47jVmr9gPw8r0d6dMm1OSORFyHgo6IiImWbstkylepADzVvw2/iWluckcirkVBR0TEJBv2nyDxbykYBjzQswWP9bnG7JZEXI6CjoiICfYczWfMx5soKbMT1z6MlwZpQUCRuqCgIyJyhWXaihg+J4m802XEtAzk7fs74+6mkCNSFxR0RESuIFtRKSPmbCTTdpprmvrx0YNdaeTpbnZbIi6rxkFn9erV3HnnnURERGCxWPj8888r7DcMgxdffJGIiAh8fHzo06cPO3bsqFBTXFzM+PHjCQkJwc/Pj0GDBnHo0KEKNTk5OSQkJGC1WrFarSQkJJCbm1uhJj09nTvvvBM/Pz9CQkKYMGECJSUlNR2SiMgVUVxWztgFm9h9NJ9Qf2/mj+pOoJ+X2W2JuLQaB53CwkI6derEjBkzqtz/6quv8sYbbzBjxgw2btxIeHg4t912G/n5+Y6axMRElixZwuLFi1mzZg0FBQXEx8dTXl7uqBk6dCgpKSksW7aMZcuWkZKSQkJCgmN/eXk5AwcOpLCwkDVr1rB48WI+/fRTJk6cWNMhiYjUObvdYOLff2b9/pM09vZg7shuNA/UgoAidc64DICxZMkSx892u90IDw83Xn75Zce206dPG1ar1Xj//fcNwzCM3Nxcw9PT01i8eLGj5vDhw4abm5uxbNkywzAMIzU11QCM9evXO2rWrVtnAMauXbsMwzCMpUuXGm5ubsbhw4cdNZ988onh7e1t2Gy2avVvs9kMoNr1IiKX6k9f7jBaPvOVce1zXxtr9h4zux0Rp1aT9+9avUYnLS2NrKws4uLiHNu8vb3p3bs3a9euBSA5OZnS0tIKNREREURHRztq1q1bh9VqpUePHo6anj17YrVaK9RER0cTERHhqOnfvz/FxcUkJydX2V9xcTF5eXkVHiIide2jn/bz0Zo0AF77TSduvDbE5I5EGo5aDTpZWVkAhIWFVdgeFhbm2JeVlYWXlxeBgYEXrAkNrbwyaGhoaIWac18nMDAQLy8vR825pk+f7rjmx2q1EhkZeQmjFBGpvi9+PsL/fb0TgMm3t+WuzleZ3JFIw1Ind12duxaEYRgXXR/i3Jqq6i+l5n9NnjwZm83meGRkZFywJxGRy7F233Em/j0FgBG9rubhW1qb25BIA1SrQSc8PByg0oxKdna2Y/YlPDyckpIScnJyLlhz9OjRSs9/7NixCjXnvk5OTg6lpaWVZnrO8vb2JiAgoMJDRKQu7MzMY+yCZErLDe7oEM7z8e21IKCICWo16LRq1Yrw8HBWrFjh2FZSUsKqVavo1asXADExMXh6elaoyczMZPv27Y6a2NhYbDYbSUlJjpoNGzZgs9kq1Gzfvp3MzExHzfLly/H29iYmJqY2hyUiUiOHc4sYMTeJ/OIyul8dxBuDb9CCgCIm8ajpAQUFBezbt8/xc1paGikpKQQFBdGiRQsSExOZNm0aUVFRREVFMW3aNHx9fRk6dCgAVquV0aNHM3HiRIKDgwkKCmLSpEl06NCBfv36AdCuXTsGDBjAmDFjmDVrFgAPP/ww8fHxtGnTBoC4uDjat29PQkICr732GidPnmTSpEmMGTNGMzUiYprcUyUMn5PE0bxiokIb86EWBBQxV01v6Vq5cqUBVHoMHz7cMIwzt5j/8Y9/NMLDww1vb2/jlltuMbZt21bhOYqKioxx48YZQUFBho+PjxEfH2+kp6dXqDlx4oQxbNgww9/f3/D39zeGDRtm5OTkVKg5ePCgMXDgQMPHx8cICgoyxo0bZ5w+fbraY9Ht5SJSm4pKyozfzPy30fKZr4weU78zDuecMrslEZdUk/dvi2EYhok5y1R5eXlYrVZsNptmgUTkspTbDcb9dTPfbM/Cv5EH/3gklrbh+rsiUhdq8v6t77oSEblMhmEw5csdfLM9Cy93Nz5I6KqQI1JPKOiIiFymWav3M3/dQQDeGNKJ2GuCTe5IRM5S0BERuQyfbT7Ey9/sAuD5+PbEd4y4yBEiciUp6IiIXKKf9h7j6X9uBWDMza0YfVMrkzsSkXMp6IiIXILth208siCZMrvBnZ0imHx7O7NbEpEqKOiIiNRQxslTjJy3kcKScmJbB/P6fR1x04KAIvWSgo6ISA2cLDyzIOCx/GLahvsz68EYvD20IKBIfaWgIyJSTUUl5Tw0fyP7jxcSYW3EvJHdCWjkaXZbInIBCjoiItVQVm5n/Cdb2Jyei9XHk/mjuhNubWR2WyJyEQo6IiIXYRgGL3yxg+92HsXLw42PhnclKszf7LZEpBoUdERELuLdlfv464Z0LBZ4+7c30O3qILNbEpFqUtAREbmAv2/K4PXlewB48c7rGRDdzOSORKQmFHRERM5j5e5sJn+2DYBH+1zD8F5Xm9uQiNSYgo6ISBV+zsjlsYWbKbcb3NP5Kp7u38bslkTkEijoiIic4+CJQkbN20hRaTk3R4Xw8r0dsVi0IKCIM1LQERH5H8cLinlwThInCku4PiKAmQ/E4OWhP5Uizkq/vSIi/3GqpIzR8zZy8MQpmgf6MHdkNxp7e5jdlohcBgUdERHOLAj4+KLN/HzIRhPfMwsChvprQUARZ6egIyINnmEY/H7JdlbuPkYjTzdmD+/GNU0bm92WiNQCBR0RafDe+m4vf9uUgZsF3rm/CzEtA81uSURqiYKOiDRonySl85fv9wLwp7uiua19mMkdiUhtUtARkQbru9Sj/H7JmQUBJ/zqWob1aGlyRyJS2xR0RKRB2pyew7hPNmM34L6Y5vzutuvMbklE6oCCjog0OPuPFTB63kZOl9rp06Yp0+7poAUBRVyUgo6INCjZ+acZPjeJnFOldGxu5d2hXfB0159CEVel324RaTAKissYNW8jGSeLaBnsy5wR3fDTgoAiLk1BR0QahNJyO48uTGb74TyC/byYP7I7IY29zW5LROqYgo6IuDzDMHjm0638tPc4Pp7uzB7RjatD/MxuS0SuAAUdEXF5ry/fzWebD+PuZuG9YV24IbKJ2S2JyBWioCMiLm3BugO8u/IXAKbf3YFb24aa3JGIXEkKOiLispZtz+KFL3YA8ORt1zG4W6TJHYnIlaagIyIuadOBkzyxeAuGAfd3b8H4X11rdksiYgIFHRFxOTsz8xg9fxPFZXb6tQvlT7++XgsCijRQCjoi4lLSjheSMDsJW1EpXVo04Z37u+ChBQFFGqxa/+0vKyvjD3/4A61atcLHx4fWrVszZcoU7Ha7o8YwDF588UUiIiLw8fGhT58+7Nixo8LzFBcXM378eEJCQvDz82PQoEEcOnSoQk1OTg4JCQlYrVasVisJCQnk5ubW9pBExEkcyS3igY82cLygmHbNApg7sjs+Xu5mtyUiJqr1oPPKK6/w/vvvM2PGDHbu3Mmrr77Ka6+9xjvvvOOoefXVV3njjTeYMWMGGzduJDw8nNtuu438/HxHTWJiIkuWLGHx4sWsWbOGgoIC4uPjKS8vd9QMHTqUlJQUli1bxrJly0hJSSEhIaG2hyQiTuB4QTEPzN7A4dwiWof48fGo7lh9PM1uS0RMZjEMw6jNJ4yPjycsLIzZs2c7tt177734+vqyYMECDMMgIiKCxMREnnnmGeDM7E1YWBivvPIKY8eOxWaz0bRpUxYsWMCQIUMAOHLkCJGRkSxdupT+/fuzc+dO2rdvz/r16+nRowcA69evJzY2ll27dtGmTZuL9pqXl4fVasVmsxEQEFCb/xlE5AqyFZVy/wfrSc3MI8LaiH882ourmviY3ZaI1JGavH/X+ozOTTfdxPfff8+ePXsA+Pnnn1mzZg133HEHAGlpaWRlZREXF+c4xtvbm969e7N27VoAkpOTKS0trVATERFBdHS0o2bdunVYrVZHyAHo2bMnVqvVUXOu4uJi8vLyKjxExLmdKilj9LyNpGbmEdLYi4UP9VDIERGHWv82u2eeeQabzUbbtm1xd3envLycqVOncv/99wOQlZUFQFhYWIXjwsLCOHjwoKPGy8uLwMDASjVnj8/KyiI0tPLCX6GhoY6ac02fPp2XXnrp8gYoIvVGcVk5Yxcks+lgDv6NPPh4VA9aN21sdlsiUo/U+ozO3/72NxYuXMhf//pXNm/ezPz583n99deZP39+hbpzb/U0DOOit3+eW1NV/YWeZ/LkydhsNscjIyOjusMSkXqmrNxO4uIUx/dXzRvZjfYR+ghaRCqq9Rmdp556imeffZbf/va3AHTo0IGDBw8yffp0hg8fTnh4OHBmRqZZs2aO47Kzsx2zPOHh4ZSUlJCTk1NhVic7O5tevXo5ao4ePVrp9Y8dO1Zptugsb29vvL31bcUizs5uN3j2s218sz0LL3c3PnywKzEtg8xuS0TqoVqf0Tl16hRubhWf1t3d3XF7eatWrQgPD2fFihWO/SUlJaxatcoRYmJiYvD09KxQk5mZyfbt2x01sbGx2Gw2kpKSHDUbNmzAZrM5akTE9RiGwZSvUvln8iHc3Sy8fX9nbooKMbstEamnan1G584772Tq1Km0aNGC66+/ni1btvDGG28watQo4MzHTYmJiUybNo2oqCiioqKYNm0avr6+DB06FACr1cro0aOZOHEiwcHBBAUFMWnSJDp06EC/fv0AaNeuHQMGDGDMmDHMmjULgIcffpj4+Phq3XElIs7pze/2Mm/tAQBevbcjA6LDzW1IROq1Wg8677zzDs8//zyPPfYY2dnZREREMHbsWF544QVHzdNPP01RURGPPfYYOTk59OjRg+XLl+Pv7++oefPNN/Hw8GDw4MEUFRXRt29f5s2bh7v7fxf/WrRoERMmTHDcnTVo0CBmzJhR20MSkXrio5/28/b3ewF4adD13BvT3OSORKS+q/V1dJyJ1tERcR5/25jOM59uA2BS3HWM+1WUyR2JiFlMXUdHRKS2fbX1CM9+dibkjL2lNY/fqm8iF5HqUdARkXpt5e5sfve3FAwD7u/egmdvb6tvIheRalPQEZF6a8P+EzyyIJnScoM7O0Xwf3dFK+SISI0o6IhIvbTtkI3R8zdRXGbnV21DeWNwJ9zdFHJEpGYUdESk3tl7NJ8H52ygoLiMHq2CeG9YFzzd9edKRGpOfzlEpF7JOHmKB2ZvIOdUKZ2aW/loeFcaebpf/EARkSoo6IhIvZGdd5phH23gaF4xUaGNmTeyO/6NPM1uS0ScmIKOiNQLOYUlPDB7A+knT9EiyJeFD/Ug0M/L7LZExMkp6IiI6QqKyxgxN4k9RwsIC/Bm0UM9CAtoZHZbIuICFHRExFSnS8t5aP5Gfj5kI9DXk4WjexAZ5Gt2WyLiIhR0RMQ0peV2Hl+0mfX7T9LY24P5o7oTFeZ/8QNFRKpJQUdETFFuN5j495/5flc23h5uzB7elY7Nm5jdloi4GAUdEbniDMPg+X9t54ufj+DhZuH9B2Lo0TrY7LZExAUp6IjIFWUYBi8v28VfN6RjscCbQ27g1rahZrclIi5KQUdErqj3fvyFWav2AzD97g7c2SnC5I5ExJUp6IjIFfPxugO89u1uAH5/Rzt+272FyR2JiKtT0BGRK2LJlkO88K8dAEz41bWMuaW1yR2JSEOgoCMide7bHVlM+sdWAEb0uprf3XadyR2JSEOhoCMiderf+44z/q9bKLcb3NulOS/Et8disZjdlog0EAo6IlJnkg/mMObjTZSU2xlwfTiv3NsBNzeFHBG5chR0RKRO7MzMY+TcJE6VlHNzVAh/uf8GPNz1J0dEriz91RGRWpd2vJCE2UnknS4jpmUgsxJi8PZwN7stEWmAFHREpFYdyS3igY82cLygmPbNApgzohu+Xh5mtyUiDZSCjojUmuMFxTzw0QYO5xbROsSPj0d3x+rjaXZbItKAKeiISK2wFZXy4Owk9h8v5KomPix8qAchjb3NbktEGjgFHRG5bKdKyhg1byOpmXmENPZiwejuRDTxMbstEREFHRG5PMVl5YxdkEzywRwCGnmwYHQPWjdtbHZbIiKAgo6IXIaycjtPfJLCT3uP4+vlzrxR3WnXLMDstkREHBR0ROSS2O0Gz362jWU7svByd+ODhK50aRFodlsiIhUo6IhIjRmGwZSvUvln8iHc3Sy8M7QzN0WFmN2WiEglCjoiUmNvfreXeWsPAPDabzrS//pwcxsSETkPBR0RqZGPftrP29/vBWDKr6/nni7NTe5IROT8FHREpNoWJ6Xzf1/vBOCp/m14MPZqcxsSEbkIBR0RqZavth5h8pJtAIzt3ZrH+lxjckciIhdXJ0Hn8OHDPPDAAwQHB+Pr68sNN9xAcnKyY79hGLz44otERETg4+NDnz592LFjR4XnKC4uZvz48YSEhODn58egQYM4dOhQhZqcnBwSEhKwWq1YrVYSEhLIzc2tiyGJNGgrd2WTuDgFw4ChPVrw7IC2WCwWs9sSEbmoWg86OTk53HjjjXh6evLNN9+QmprKn//8Z5o0aeKoefXVV3njjTeYMWMGGzduJDw8nNtuu438/HxHTWJiIkuWLGHx4sWsWbOGgoIC4uPjKS8vd9QMHTqUlJQUli1bxrJly0hJSSEhIaG2hyTSoG3Yf4JHFiZTZjcY1CmCP/06WiFHRJyGxTAMozaf8Nlnn+Xf//43P/30U5X7DcMgIiKCxMREnnnmGeDM7E1YWBivvPIKY8eOxWaz0bRpUxYsWMCQIUMAOHLkCJGRkSxdupT+/fuzc+dO2rdvz/r16+nRowcA69evJzY2ll27dtGmTZuL9pqXl4fVasVmsxEQoEXORM617ZCN+z9cT0FxGX3bhvJ+Qgye7vrEW0TMVZP371r/i/XFF1/QtWtX7rvvPkJDQ+ncuTMffvihY39aWhpZWVnExcU5tnl7e9O7d2/Wrl0LQHJyMqWlpRVqIiIiiI6OdtSsW7cOq9XqCDkAPXv2xGq1OmrOVVxcTF5eXoWHiFRt79F8HpyzgYLiMnq2DuLdYV0UckTE6dT6X639+/czc+ZMoqKi+Pbbb3nkkUeYMGECH3/8MQBZWVkAhIWFVTguLCzMsS8rKwsvLy8CAwMvWBMaGlrp9UNDQx0155o+fbrjeh6r1UpkZOTlDVbERWWcPMUDszeQc6qUTs2tfDS8G4083c1uS0Skxmo96Njtdrp06cK0adPo3LkzY8eOZcyYMcycObNC3bmf8RuGcdHP/c+tqar+Qs8zefJkbDab45GRkVHdYYk0GEfzTjPsow0czSvmurDGzBvZncbeHma3JSJySWo96DRr1oz27dtX2NauXTvS09MBCA8/s4LqubMu2dnZjlme8PBwSkpKyMnJuWDN0aNHK73+sWPHKs0WneXt7U1AQECFh4j8V05hCQmzN5B+8hQtgnxZOLoHgX5eZrclInLJaj3o3HjjjezevbvCtj179tCyZUsAWrVqRXh4OCtWrHDsLykpYdWqVfTq1QuAmJgYPD09K9RkZmayfft2R01sbCw2m42kpCRHzYYNG7DZbI4aEam+/NOlDJ+bxJ6jBYQFeLPooR6EBjQyuy0RkctS6/PRv/vd7+jVqxfTpk1j8ODBJCUl8cEHH/DBBx8AZz5uSkxMZNq0aURFRREVFcW0adPw9fVl6NChAFitVkaPHs3EiRMJDg4mKCiISZMm0aFDB/r16wecmSUaMGAAY8aMYdasWQA8/PDDxMfHV+uOKxH5r9Ol5Tw0fxNbD9kI9PVk4egeRAb5mt2WiMhlq/Wg061bN5YsWcLkyZOZMmUKrVq14q233mLYsGGOmqeffpqioiIee+wxcnJy6NGjB8uXL8ff399R8+abb+Lh4cHgwYMpKiqib9++zJs3D3f3/14QuWjRIiZMmOC4O2vQoEHMmDGjtock4tJKy+08tmgzG9JO0tjbg49H9SAqzP/iB4qIOIFaX0fHmWgdHWnoyu0GiX9L4cufj+Dt4cbHo7rTo3Ww2W2JiFyQqevoiIhzMAyDP3y+nS9/PoKHm4X3E2IUckTE5SjoiDRAhmHw8je7+CQpHTcLvPXbG7i1TeV1qUREnJ2CjkgD9N6PvzBr9X4Apt/TgfiOESZ3JCJSNxR0RBqY+WsP8Nq3Z5aA+MPAdgzp1sLkjkRE6o6CjkgD8mnyIf74xQ4AJvSN4qGbW5vckYhI3VLQEWkglm3P4ulPtwIwotfV/K5flMkdiYjUPQUdkQZgzd7jTPhkC+V2g9/ENOeF+PYX/W45ERFXoKAj4uKSD+Yw5uNNlJTbuT06nJfv6YCbm0KOiDQMCjoiLiz1SB4j5yZRVFrOzVEhvPXbG/Bw16+9iDQc+osn4qL2HyvgwTkbyDtdRteWgcxKiMHbw/3iB4qIuBAFHREXdDi3iAc+2sDxghLaNwtg9ohu+HrV+lfbiYjUewo6Ii7mWH4xCR9t4IjtNK2b+vHx6O5YfTzNbktExBQKOiIuxFZUyoNzkth/vJCrmviwcHQPQhp7m92WiIhpFHREXMSpkjJGzdvIzsw8Qhp7s/ChHkQ08TG7LRERUynoiLiA4rJyxi5IJvlgDgGNPFgwujutQvzMbktExHQKOiJOrqzczoRPtvDT3uP4erkzb1R32jULMLstEZF6QUFHxInZ7QbPfLqNb3ccxcvdjQ8f7EqXFoFmtyUiUm8o6Ig4KcMwmPJVKp9uPoS7m4UZQztz47UhZrclIlKvKOiIOKk3V+xh3toDALx+X0firg83tyERkXpIQUfECX24ej9v/7APgD/9+nru7tzc5I5EROonBR0RJ/NJUjpTl+4E4Kn+bUiIvdrchkRE6jEFHREn8uXPR3huyTYAHul9DY/feq3JHYmI1G8KOiJOYuWubH73txQMA4b1aMEzA9qY3ZKISL2noCPiBNbvP8EjC5Mpsxv8+oYI/vTraCwWi9ltiYjUewo6IvXc1kO5PDR/E8Vldvq1C+X1+zrh5qaQIyJSHQo6IvXY3qP5DJ+TREFxGbGtg5kxtAue7vq1FRGpLv3FFKmn0k+cYthHG8g5VUqnyCZ8OLwrjTzdzW5LRMSpKOiI1ENH807zwOwNZOcX0ybMn/kju9HY28PstkREnI6Cjkg9k1NYwgMfbSD95ClaBvuyYHR3mvh6md2WiIhTUtARqUfyT5cyfG4Se7MLCA9oxMLRPQgNaGR2WyIiTktBR6SeOF1azkPzN7H1kI0gPy8WPtSdyCBfs9sSEXFqCjoi9UBJmZ3HFm1mQ9pJ/L09+HhUd64N9Te7LRERp6egI2KycrvBk39P4Ydd2TTydGP2iG5EX2U1uy0REZegoCNiIsMw+MPn2/hqayae7hbefyCG7q2CzG5LRMRl1HnQmT59OhaLhcTERMc2wzB48cUXiYiIwMfHhz59+rBjx44KxxUXFzN+/HhCQkLw8/Nj0KBBHDp0qEJNTk4OCQkJWK1WrFYrCQkJ5Obm1vWQRGqFYRhM/2YXnyRl4GaBt4Z0pk+bULPbEhFxKXUadDZu3MgHH3xAx44dK2x/9dVXeeONN5gxYwYbN24kPDyc2267jfz8fEdNYmIiS5YsYfHixaxZs4aCggLi4+MpLy931AwdOpSUlBSWLVvGsmXLSElJISEhoS6HJFJr3l25jw9W7wfg5Xs6MrBjM5M7EhFxQUYdyc/PN6KioowVK1YYvXv3Np544gnDMAzDbrcb4eHhxssvv+yoPX36tGG1Wo3333/fMAzDyM3NNTw9PY3Fixc7ag4fPmy4ubkZy5YtMwzDMFJTUw3AWL9+vaNm3bp1BmDs2rWrWj3abDYDMGw22+UOV6RG5v07zWj5zFdGy2e+Mj5c/YvZ7YiIOJWavH/X2YzO448/zsCBA+nXr1+F7WlpaWRlZREXF+fY5u3tTe/evVm7di0AycnJlJaWVqiJiIggOjraUbNu3TqsVis9evRw1PTs2ROr1eqoOVdxcTF5eXkVHiJX2qfJh/jjF2c+qn2ibxQP3dza5I5ERFxXnawpv3jxYjZv3szGjRsr7cvKygIgLCyswvawsDAOHjzoqPHy8iIwMLBSzdnjs7KyCA2tfD1DaGioo+Zc06dP56WXXqr5gERqybLtWTz1z58BGHnj1ST2izK5IxER11brMzoZGRk88cQTLFy4kEaNzr+iq8ViqfCzYRiVtp3r3Jqq6i/0PJMnT8ZmszkeGRkZF3w9kdr0095jTPhkC3YD7otpzvMD21/037yIiFyeWg86ycnJZGdnExMTg4eHBx4eHqxatYq3334bDw8Px0zOubMu2dnZjn3h4eGUlJSQk5NzwZqjR49Wev1jx45Vmi06y9vbm4CAgAoPkSsh+WAOD3+cTEm5ndujw5l+Twfc3BRyRETqWq0Hnb59+7Jt2zZSUlIcj65duzJs2DBSUlJo3bo14eHhrFixwnFMSUkJq1atolevXgDExMTg6elZoSYzM5Pt27c7amJjY7HZbCQlJTlqNmzYgM1mc9SI1AepR/IYOTeJotJybrmuKW/99gY83LWElYjIlVDr1+j4+/sTHR1dYZufnx/BwcGO7YmJiUybNo2oqCiioqKYNm0avr6+DB06FACr1cro0aOZOHEiwcHBBAUFMWnSJDp06OC4uLldu3YMGDCAMWPGMGvWLAAefvhh4uPjadOmTW0PS+SS7D9WwINzNpB3uoyuLQN5/4EueHu4m92WiEiDUScXI1/M008/TVFREY899hg5OTn06NGD5cuX4+//3+/2efPNN/Hw8GDw4MEUFRXRt29f5s2bh7v7f98kFi1axIQJExx3Zw0aNIgZM2Zc8fGIVOVwbhEPfLSB4wUlXB8RwJyR3fD1MuVXTkSkwbIYhmGY3YRZ8vLysFqt2Gw2Xa8jtepYfjFDZq1j//FCrmnqx9/HxhLc2NvstkREXEJN3r91oYBILbOdKuXBOUnsP17IVU18WPhQD4UcERGTKOiI1KJTJWWMnJfEzsw8Qhp7s+ihHjSz+pjdlohIg6WgI1JLisvKGbsgmc3puVh9PFn4UHeuDvEzuy0RkQZNQUekFpSV25nwyRZ+2nscXy935o3sRttwXfclImI2BR2Ry2S3Gzz96Va+3XEULw83PnqwK51bBF78QBERqXMKOiKXwTAMpnyVymebD+PuZuHdoV3odW2I2W2JiMh/KOiIXIY3Vuxh3toDWCzw5/s6cVv7qr9+REREzKGgI3KJPlj9C+/8sA+AKb+O5q7OV5nckYiInEtBR+QSfJKUzrSluwB4ekAbEnq2NLkjERGpioKOSA19+fMRnluyDYBH+1zDY32uNbkjERE5HwUdkRr4YddRfve3FAwDHujZgqf76wtkRUTqMwUdkWpav/8Ejy7cTJnd4Nc3RDBlUDQWi8XstkRE5AIUdESqYeuhXB6av4niMjv92oXx+n2dcHNTyBERqe8UdEQuYs/RfB6ck0RBcRmxrYOZMbQznu761RERcQb6ay1yAeknTvHARxvIPVXKDZFN+HB4Vxp5upvdloiIVJOCjsh5HM07zbDZ68nOL6ZNmD/zRnajsbeH2W2JiEgNKOiIVOFEQTEPfLSBjJNFtAz2ZcHo7jTx9TK7LRERqSH976nI/yi3GyzemM6fl+/hZGEJ4QGNWDi6B6EBjcxuTURELoGCjsh/rP3lOFO+TGVXVj4A14Y25v0HYogM8jW5MxERuVQKOtLgpZ84xdSlqXy74ygAVh9PftcvimE9W+ruKhERJ6egIw1WQXEZM37Yx5w1aZSU23F3s/BAjxYk9ruOQD9djyMi4goUdKTBsdsN/pl8iFe/3c3xgmIAbo4K4fn49lwX5m9ydyIiUpsUdKRB2XjgJC99uYPth/MAaBXixx8GtuNXbUP1dQ4iIi5IQUcahEM5p5j+zS6+3poJgL+3B0/0i+LB2Kvx8tB1OCIirkpBR1zaqZIyZv74Cx+s3k9xmR03Cwzp1oKJcdcR0tjb7PZERKSOKeiIS7LbDf7182Fe/mYXR/POXIfTs3UQL8RfT/uIAJO7ExGRK0VBR1zOlvQcXvoylZSMXAAig3z4/R3t6H99uK7DERFpYBR0xGVk2U7zyrJdLNlyGAA/L3ce/9W1jLqxlb6IU0SkgVLQEad3urScD1bvZ+aPv1BUWo7FAr/p0pyn+rfRVzeIiDRwCjritAzD4Kutmbz8zS4O5xYB0LVlIH+883o6NLea3J2IiNQHCjrilLYdsjHlqx1sPJADwFVNfHj29rbEd2ym63BERMRBQUecSnb+aV5btpt/bj6EYYCPpzuP9rmGh29pretwRESkEgUdcQqnS8uZ8+803v1hH4Ul5QDc3fkqnh7QhmZWH5O7ExGR+kpBR+o1wzD4dkcWU5fuJOPkmetwbohswgt3tqdLi0CTuxMRkfqu1te+nz59Ot26dcPf35/Q0FDuuusudu/eXaHGMAxefPFFIiIi8PHxoU+fPuzYsaNCTXFxMePHjyckJAQ/Pz8GDRrEoUOHKtTk5OSQkJCA1WrFarWSkJBAbm5ubQ9JTJJ6JI/7P1zPIws3k3GyiLAAb94c0onPHu2lkCMiItVS60Fn1apVPP7446xfv54VK1ZQVlZGXFwchYWFjppXX32VN954gxkzZrBx40bCw8O57bbbyM/Pd9QkJiayZMkSFi9ezJo1aygoKCA+Pp7y8nJHzdChQ0lJSWHZsmUsW7aMlJQUEhISantIcoWdKChm8mfbiH/nJ9bvP4m3hxsTfnUtKyf14e7OzXFz08XGIiJSPRbDMIy6fIFjx44RGhrKqlWruOWWWzAMg4iICBITE3nmmWeAM7M3YWFhvPLKK4wdOxabzUbTpk1ZsGABQ4YMAeDIkSNERkaydOlS+vfvz86dO2nfvj3r16+nR48eAKxfv57Y2Fh27dpFmzZtKvVSXFxMcXGx4+e8vDwiIyOx2WwEBOhrAcxWUmZn/toDvP39XvKLywCI79iMZ29vS/NAX5O7ExGR+iIvLw+r1Vqt9+86/9pmm80GQFBQEABpaWlkZWURFxfnqPH29qZ3796sXbsWgOTkZEpLSyvUREREEB0d7ahZt24dVqvVEXIAevbsidVqddSca/r06Y6PuaxWK5GRkbU7WLkkhmHwXepR+r+1mqlLd5JfXEb0VQH8fWwsM4Z2UcgREZFLVqcXIxuGwZNPPslNN91EdHQ0AFlZWQCEhYVVqA0LC+PgwYOOGi8vLwIDAyvVnD0+KyuL0NDQSq8ZGhrqqDnX5MmTefLJJx0/n53REfPsOZrPn75K5ae9xwEIaezN0/3b8JsYfUQlIiKXr06Dzrhx49i6dStr1qyptO/cRd0Mw7joQm/n1lRVf6Hn8fb2xtvbuzqtSx3LKSzhre/2sHBDOuV2Ay93N0bd1IrHb70G/0aeZrcnIiIuos6Czvjx4/niiy9YvXo1zZs3d2wPDw8HzszINGvWzLE9OzvbMcsTHh5OSUkJOTk5FWZ1srOz6dWrl6Pm6NGjlV732LFjlWaLpP4oLbezaP1B3vxuL7aiUgD6Xx/Gc3e0o2Wwn8ndiYiIq6n1a3QMw2DcuHF89tln/PDDD7Rq1arC/latWhEeHs6KFSsc20pKSli1apUjxMTExODp6VmhJjMzk+3btztqYmNjsdlsJCUlOWo2bNiAzWZz1Ej98uPubG7/y0+8+GUqtqJS2ob789cxPZiV0FUhR0RE6kStz+g8/vjj/PWvf+Vf//oX/v7+jutlrFYrPj4+WCwWEhMTmTZtGlFRUURFRTFt2jR8fX0ZOnSoo3b06NFMnDiR4OBggoKCmDRpEh06dKBfv34AtGvXjgEDBjBmzBhmzZoFwMMPP0x8fHyVd1yJeX45VsDUr3fyw65sAIL8vJgYdx2/7dYCd12HIyIidajWg87MmTMB6NOnT4Xtc+fOZcSIEQA8/fTTFBUV8dhjj5GTk0OPHj1Yvnw5/v7+jvo333wTDw8PBg8eTFFREX379mXevHm4u//3+4wWLVrEhAkTHHdnDRo0iBkzZtT2kOQS2YpKefv7vcxfe4Ayu4GHm4URva5mfN8orD66DkdEROpena+jU5/V5D58qb6ycjuLN2bwxoo9nCwsAaBv21B+P7AdrZs2Nrk7ERFxdjV5/9Z3XUmtWrvvOFO+SmVX1plVrqNCG/N8fHtuua6pyZ2JiEhDpKAjteLgiUKmfr2T5aln7oSz+njy5G3XMaxHCzzc63xdShERkSop6MhlyT9dyoyV+5i75gAl5Xbc3Swk9GxJYr8omvh6md2eiIg0cAo6cknK7Qb/TM7gtW/3cLzgzPeH3XJdU54f2I6oMP+LHC0iInJlKOhIjSWlneSlL3ew40geAK1D/PhDfDtubRN60dWtRUREriQFHam2jJOnePmbXXy9LRMA/0YePNE3igdjr8bLQ9fhiIhI/aOgIxdVWFzGzB9/4YOf9lNSZsfNAr/t3oKJt11HcGN9d5iIiNRfCjpyXna7wecph3ll2S6O5p25Die2dTAv3Nmeds207pCIiNR/CjpSpc3pObz0ZSo/Z+QC0CLIl+fuaEf/68N0HY6IiDgNBR2pINNWxMvf7OJfKUcA8PNyZ9yvohh109V4e7hf5GgREZH6RUFHACgqKeeD1fuZuWofp0vtWCxwX0xzJvVvQ6h/I7PbExERuSQKOg2cYRh8uTWTl5fu5IjtNADdrg7khfjr6dDcanJ3IiIil0dBpwHbeiiXKV+msulgDgBXNfFh8h1tGdihma7DERERl6Cg0wBl553m1W9388/kQwD4eLrzWJ9rGHNLaxp56jocERFxHQo6Dcjp0nJmr0njvZX7KCwpB+Cezlfx9IC2hFt1HY6IiLgeBZ0GwDAMlm3PYto3O8k4WQRA5xZNeCG+PZ1bBJrcnYiISN1R0HFxO47YmPJlKhvSTgIQHtCIZ29vy6BOEbi56TocERFxbQo6Lup4QTF/Xr6bxRszMAzw9nBjbO9reKR3a3y9dNpFRKRh0DueiykpszNvbRrvfL+P/OIyAOI7NuPZ29vSPNDX5O5ERESuLAUdF2EYBt/tzGbq16kcOHEKgA5XWXnhzvZ0uzrI5O5ERETMoaDjAnZn5fOnr1JZs+84AE39vXmqfxt+06W5rsMREZEGTUHHieUUlvDmd3tYtCGdcruBl7sbo29uxeO3Xktjb51aERERvRs6odJyOwvWHeSt7/aQd/rMdTgDrg/nuTva0SJY1+GIiIicpaDjZFbuzub/vkrll2OFALQN9+ePd15P7DXBJncmIiJS/yjoOIl92QX839ep/Lj7GADBfl5MjGvDkG6RuOs6HBERkSop6NRztlOl/OX7vXy87gBldgNPdwsjel3N+L5RBDTyNLs9ERGRek1Bp54qK7fzycYM3li+m5xTpQD0axfK7we2p1WIn8ndiYiIOAcFnXro3/uOM+XLVHYfzQcgKrQxz8e355brmprcmYiIiHNR0KlHDhwvZOrSnaxIPQpAE19PnrztOoZ2b4GHu5vJ3YmIiDgfBZ16IP90KTN+2Mecf6dRWm7g7mYhoWdLEvtF0cTXy+z2REREnJaCjonK7Qb/2JTB68t3c7ygBIBbrmvKC/HtuDbU3+TuREREnJ+Cjkk27D/BlK9S2XEkD4DWIX78Ib4dt7YJxWLR7eIiIiK1QUHnCss4eYrp3+xk6bYsAPwbeZDY7zoSerbEy0PX4YiIiNQmBZ0rpLC4jPd+3MeHP6VRUmbHzQL3d2/Bk7ddR3Bjb7PbExERcUlOP4Xw3nvv0apVKxo1akRMTAw//fST2S1VYLcb/DP5ELe+/iPvrvyFkjI7va4J5usJNzP17g4KOSIiInXIqWd0/va3v5GYmMh7773HjTfeyKxZs7j99ttJTU2lRYsWZrdH8sGTTPkylZ8P2QBoGezLc3e0I659mK7DERERuQIshmEYZjdxqXr06EGXLl2YOXOmY1u7du246667mD59eqX64uJiiouLHT/n5eURGRmJzWYjICCg1vo6klvEy9/s4oufjwDQ2NuDcb+6lpE3Xo23h3utvY6IiEhDlJeXh9Vqrdb7t9N+dFVSUkJycjJxcXEVtsfFxbF27doqj5k+fTpWq9XxiIyMrJPe5q89wBc/H8FigSFdI/lhUm8e6X2NQo6IiMgV5rQfXR0/fpzy8nLCwsIqbA8LCyMrK6vKYyZPnsyTTz7p+PnsjE5te+zWa0k7XsiEvlFEX2Wt9ecXERGR6nHaoHPWude6GIZx3utfvL298fau+4t/rT6efPBg1zp/HREREbkwp/3oKiQkBHd390qzN9nZ2ZVmeURERKRhctqg4+XlRUxMDCtWrKiwfcWKFfTq1cukrkRERKQ+ceqPrp588kkSEhLo2rUrsbGxfPDBB6Snp/PII4+Y3ZqIiIjUA04ddIYMGcKJEyeYMmUKmZmZREdHs3TpUlq2bGl2ayIiIlIPOPU6OperJvfhi4iISP3QINbREREREbkYBR0RERFxWQo6IiIi4rIUdERERMRlKeiIiIiIy1LQEREREZeloCMiIiIuS0FHREREXJZTr4x8uc6ulZiXl2dyJyIiIlJdZ9+3q7PmcYMOOvn5+QBERkaa3ImIiIjUVH5+Plar9YI1DforIOx2O0eOHMHf3x+LxVKrz52Xl0dkZCQZGRku+fUSGp/zc/UxanzOz9XH6Orjg7obo2EY5OfnExERgZvbha/CadAzOm5ubjRv3rxOXyMgIMBl/wGDxucKXH2MGp/zc/Uxuvr4oG7GeLGZnLN0MbKIiIi4LAUdERERcVkKOnXE29ubP/7xj3h7e5vdSp3Q+Jyfq49R43N+rj5GVx8f1I8xNuiLkUVERMS1aUZHREREXJaCjoiIiLgsBR0RERFxWQo6IiIi4rIUdERERMRlKehcovfee49WrVrRqFEjYmJi+Omnny5Yv2rVKmJiYmjUqBGtW7fm/fffv0KdXrqajPHHH3/EYrFUeuzatesKdlx9q1ev5s477yQiIgKLxcLnn39+0WOc6RzWdHzOdv6mT59Ot27d8Pf3JzQ0lLvuuovdu3df9DhnOYeXMj5nO4czZ86kY8eOjhVzY2Nj+eabby54jLOcP6j5+Jzt/J1r+vTpWCwWEhMTL1hnxjlU0LkEf/vb30hMTOT3v/89W7Zs4eabb+b2228nPT29yvq0tDTuuOMObr75ZrZs2cJzzz3HhAkT+PTTT69w59VX0zGetXv3bjIzMx2PqKioK9RxzRQWFtKpUydmzJhRrXpnO4c1Hd9ZznL+Vq1axeOPP8769etZsWIFZWVlxMXFUVhYeN5jnOkcXsr4znKWc9i8eXNefvllNm3axKZNm/jVr37Fr3/9a3bs2FFlvTOdP6j5+M5ylvP3vzZu3MgHH3xAx44dL1hn2jk0pMa6d+9uPPLIIxW2tW3b1nj22WerrH/66aeNtm3bVtg2duxYo2fPnnXW4+Wq6RhXrlxpAEZOTs4V6K52AcaSJUsuWOOM5/Cs6ozPmc+fYRhGdna2ARirVq06b40zn8PqjM/Zz6FhGEZgYKDx0UcfVbnPmc/fWRcan7Oev/z8fCMqKspYsWKF0bt3b+OJJ544b61Z51AzOjVUUlJCcnIycXFxFbbHxcWxdu3aKo9Zt25dpfr+/fuzadMmSktL66zXS3UpYzyrc+fONGvWjL59+7Jy5cq6bPOKcrZzeKmc9fzZbDYAgoKCzlvjzOewOuM7yxnPYXl5OYsXL6awsJDY2Ngqa5z5/FVnfGc52/l7/PHHGThwIP369btorVnnUEGnho4fP055eTlhYWEVtoeFhZGVlVXlMVlZWVXWl5WVcfz48Trr9VJdyhibNWvGBx98wKeffspnn31GmzZt6Nu3L6tXr74SLdc5ZzuHNeXM588wDJ588kluuukmoqOjz1vnrOewuuNzxnO4bds2GjdujLe3N4888ghLliyhffv2VdY64/mryfic8fwtXryYzZs3M3369GrVm3UOPersmV2cxWKp8LNhGJW2Xay+qu31SU3G2KZNG9q0aeP4OTY2loyMDF5//XVuueWWOu3zSnHGc1hdznz+xo0bx9atW1mzZs1Fa53xHFZ3fM54Dtu0aUNKSgq5ubl8+umnDB8+nFWrVp03DDjb+avJ+Jzt/GVkZPDEE0+wfPlyGjVqVO3jzDiHmtGpoZCQENzd3SvNbGRnZ1dKqmeFh4dXWe/h4UFwcHCd9XqpLmWMVenZsyd79+6t7fZM4WznsDY4w/kbP348X3zxBStXrqR58+YXrHXGc1iT8VWlvp9DLy8vrr32Wrp27cr06dPp1KkTf/nLX6qsdcbzV5PxVaU+n7/k5GSys7OJiYnBw8MDDw8PVq1axdtvv42Hhwfl5eWVjjHrHCro1JCXlxcxMTGsWLGiwvYVK1bQq1evKo+JjY2tVL98+XK6du2Kp6dnnfV6qS5ljFXZsmULzZo1q+32TOFs57A21OfzZxgG48aN47PPPuOHH36gVatWFz3Gmc7hpYyvKvX5HFbFMAyKi4ur3OdM5+98LjS+qtTn89e3b1+2bdtGSkqK49G1a1eGDRtGSkoK7u7ulY4x7RzW6aXOLmrx4sWGp6enMXv2bCM1NdVITEw0/Pz8jAMHDhiGYRjPPvuskZCQ4Kjfv3+/4evra/zud78zUlNTjdmzZxuenp7GP//5T7OGcFE1HeObb75pLFmyxNizZ4+xfft249lnnzUA49NPPzVrCBeUn59vbNmyxdiyZYsBGG+88YaxZcsW4+DBg4ZhOP85rOn4nO38Pfroo4bVajV+/PFHIzMz0/E4deqUo8aZz+GljM/ZzuHkyZON1atXG2lpacbWrVuN5557znBzczOWL19uGIZznz/DqPn4nO38VeXcu67qyzlU0LlE7777rtGyZUvDy8vL6NKlS4XbPocPH2707t27Qv2PP/5odO7c2fDy8jKuvvpqY+bMmVe445qryRhfeeUV45prrjEaNWpkBAYGGjfddJPx9ddfm9B19Zy9lfPcx/Dhww3DcP5zWNPxOdv5q2psgDF37lxHjTOfw0sZn7Odw1GjRjn+vjRt2tTo27evIwQYhnOfP8Oo+fic7fxV5dygU1/OocUw/nMlkIiIiIiL0TU6IiIi4rIUdERERMRlKeiIiIiIy1LQEREREZeloCMiIiIuS0FHREREXJaCjoiIiLgsBR0RERFxWQo6IiIi4rIUdERERMRlKeiIiIiIy/p/2vn8mFPGNb8AAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "plt.title(\"MSDMJ\")\n", "plt.plot(time, msdmj, label=\"Cond.\")\n", "plt.legend();" ] } ], "metadata": { "kernelspec": { "display_name": "newanalysis-dev", "language": "python", "name": "newanalysis-dev" }, "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.10.0" } }, "nbformat": 4, "nbformat_minor": 5 }