{ "cells": [ { "cell_type": "markdown", "id": "e8258b96", "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/dielectrics.ipynb)" ] }, { "cell_type": "markdown", "id": "5515278d", "metadata": {}, "source": [ "# Dielectrics" ] }, { "cell_type": "markdown", "id": "1fb8776a", "metadata": {}, "source": [ "We calculate MDMD." ] }, { "cell_type": "code", "execution_count": 1, "id": "ff3dba1f", "metadata": {}, "outputs": [], "source": [ "import time\n", "import sys\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "from pathlib import Path\n", "import glob" ] }, { "cell_type": "markdown", "id": "0921ba1a", "metadata": {}, "source": [ "Some configurations if we are in a google colab environment:" ] }, { "cell_type": "code", "execution_count": 14, "id": "926253a1", "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/" ] }, { "cell_type": "markdown", "id": "95da3a41", "metadata": {}, "source": [ "Next we import MDAnalysis and functions needed from newanalysis." ] }, { "cell_type": "code", "execution_count": 8, "id": "a47b041a", "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.correl import correlateParallel\n", " from newanalysis.helpers import dipByResidue\n", " from newanalysis.functions import atomsPerResidue, residueFirstAtom\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": 9, "id": "d664f2b4", "metadata": {}, "outputs": [], "source": [ "base='../../test_cases/data/emim_dca_equilibrium/'\n", "psf=base+'emim_dca.psf'\n", "#Check PSF:\n", "if np.array_equal(MDAnalysis.Universe(psf).atoms.masses, MDAnalysis.Universe(psf).atoms.masses.astype(bool)):\n", " print(\"Used wrong PSF format (masses unreadable!)\")\n", " sys.exit()\n", "\n", "u=MDAnalysis.Universe(psf,base+\"emim_dca.dcd\")\n", "\n", "skip=20\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": 10, "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", "mass_cat1 = sel_cat1.masses\n", "charge_cat1 = sel_cat1.charges\n", "com_cat1 = np.zeros((n,ncat1,3),dtype=np.float64)\n", "mdcat1 = np.zeros((n,3),dtype=np.float64)\n", "apr_cat1 = atomsPerResidue(sel_cat1)\n", "rfa_cat1 = residueFirstAtom(sel_cat1)\n", "\n", "\n", "print(\"Number EMIM = \",ncat1)\n", "\n", "sel_an1 = u.select_atoms(\"resname DCA\")\n", "nan1 = sel_an1.n_residues\n", "mass_an1 = sel_an1.masses\n", "charge_an1 = sel_an1.charges\n", "com_an1 = np.zeros((n,nan1, 3),dtype=np.float64)\n", "mdan1 = np.zeros((n,3),dtype=np.float64)\n", "apr_an1 = atomsPerResidue(sel_an1)\n", "rfa_an1 = residueFirstAtom(sel_an1)\n", "\n", "print(\"Number DCA = \",nan1)" ] }, { "cell_type": "markdown", "id": "791aeff3", "metadata": {}, "source": [ "## Analysis" ] }, { "cell_type": "code", "execution_count": 11, "id": "5ef3bf49", "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": [ "md = np.zeros((n,3),dtype=np.float64)\n", "\n", "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,u.trajectory.n_frames), \"\\tElapsed time: %.2f hours\" % ((time.time()-start)/3600))\n", " \n", " # efficiently calculate center-of-mass coordinates\n", " coor_cat1 = np.ascontiguousarray(sel_cat1.positions,dtype='double')\n", " coor_an1 = np.ascontiguousarray(sel_an1.positions,dtype='double')\n", " com_an1 = sel_an1.center_of_mass(compound='residues')\n", " com_cat1 = sel_cat1.center_of_mass(compound='residues')\n", " mdcat1[ctr] += np.sum(dipByResidue(coor_cat1,charge_cat1,mass_cat1,ncat1,apr_cat1,rfa_cat1,com_cat1),axis=0)\n", " mdan1[ctr] += np.sum(dipByResidue(coor_an1,charge_an1,mass_an1,nan1,apr_an1,rfa_an1,com_an1),axis=0)\n", "\n", " #Alternatively, use\n", " #com_an1 = centerOfMassByResidue(sel_an1,coor=coor_an1,masses=mass_an1,apr=apr_an1,rfa=rfa_an1)\n", " #com_cat1 = centerOfMassByResidue(sel_cat1,coor=coor_cat1,masses=mass_cat1,apr=apr_cat1,rfa=rfa_cat1)\n", " #mdcat1[ctr] += np.sum(dipoleMomentByResidue(sel_cat1,coor=coor_cat1,charges=charge_cat1,masses=mass_cat1,com=com_cat1,apr=apr_cat1,rfa=rfa_cat1),axis=0)\n", " #mdan1[ctr] += np.sum(dipoleMomentByResidue(sel_an1,coor=coor_an1,charges=charge_an1,masses=mass_an1,com=com_an1,apr=apr_an1,rfa=rfa_an1),axis=0)\n", "\n", " md[ctr] += mdcat1[ctr]+mdan1[ctr]\n", " ctr+=1" ] }, { "cell_type": "markdown", "id": "d9d5ef3c", "metadata": {}, "source": [ "Calculate " ] }, { "cell_type": "code", "execution_count": 12, "id": "c908acd9", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "calculating \n" ] } ], "source": [ "print(\"calculating \")\n", "\n", "# total\n", "md = np.ascontiguousarray(md.T)\n", "md0mdt = np.zeros(n)\n", "correlateParallel(md,md,md0mdt,ltc=1)\n", "\n", "#f1 = open(\"md0mdt.dat\",'w')\n", "#for i in range(len(md0mdt)):\n", "# f1.write(\"%5.5f\\t%5.5f\\n\" % (i*skip*dt, md0mdt[i]))\n", "#f1.close()\n", "\n", "# cations\n", "mdcat1 = np.ascontiguousarray(mdcat1.T)\n", "md0mdt_cat = np.zeros(n)\n", "correlateParallel(mdcat1,md,md0mdt_cat,ltc=1)\n", "#f1 = open(\"md0mdt_emim.dat\",'w')\n", "#for i in range(len(md0mdt)):\n", "# f1.write(\"%5.5f\\t%5.5f\\n\" % (i*skip*dt, md0mdt[i]))\n", "#f1.close()\n", "\n", "\n", "# anions\n", "mdan1 = np.ascontiguousarray(mdan1.T)\n", "md0mdt_an = np.zeros(n)\n", "correlateParallel(mdan1,md,md0mdt_an,ltc=1)\n", "#f1 = open(\"md0mdt_dca.dat\",'w')\n", "#for i in range(len(md0mdt)):\n", "# f1.write(\"%5.5f\\t%5.5f\\n\" % (i*skip*dt, md0mdt[i]))\n", "#f1.close()" ] }, { "cell_type": "markdown", "id": "51b054b7", "metadata": {}, "source": [ "## Plot the results" ] }, { "cell_type": "code", "execution_count": 13, "id": "0b89b9e7", "metadata": { "tags": [ "nbsphinx-thumbnail" ] }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAisAAAGxCAYAAACju/aQAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/bCgiHAAAACXBIWXMAAA9hAAAPYQGoP6dpAACb4UlEQVR4nOzdd1yV5fvA8c85h40MAdm4UBHEhZo7NbdmzrRyNmxoQ61s2NC+v7Rsl2lplprlyK2ZM7VyL5Qh4kIRRQSUPc95fn8cQQkwxjkcxvV+vc7rS89zP8+5jl+Fi/u57utWKYqiIIQQQghRSalNHYAQQgghxP1IsiKEEEKISk2SFSGEEEJUapKsCCGEEKJSk2RFCCGEEJWaJCtCCCGEqNQkWRFCCCFEpSbJihBCCCEqNUlWhBBCCFGpSbIiRDW3ZMkSVCoVKpWKvXv3FjqvKAqNGjVCpVLRvXv3AufyrlOpVGg0GmrXrk3Lli157rnnOHToULHv+ffff2Npacnly5cLHD9x4gS9evWiVq1aODo6MmzYMC5evFhgTGRkJBYWFpw4caLQfSdMmIBKpcLOzo7U1NRC5y9fvoxarUalUjFz5sz843v37i3wWSwsLKhTpw6dO3dmxowZheK81wcffEBAQAA6nQ6A9PR0Zs6cWeSf5eLFi/Hy8iItLa3Y+wkhSk+SFSFqCDs7OxYvXlzo+L59+7hw4QJ2dnZFXjdixAgOHjzIP//8w8qVKxk3bhyHDh2iY8eOvPLKK4XGK4rClClTmDhxIvXq1cs/HhERQffu3cnOzmb16tX8+OOPREZG0rVrV27evJk/rkmTJowePZqpU6cWGY+5uTm5ubmsWrWq0Lmffvqp2M8BMHv2bA4ePMiePXtYvHgx3bt358cff8Tf359ffvml0Phr164xd+5cPvjgA9Rq/bfL9PR0Zs2aVWSyMn78eGxtbZk7d26xMQghykARQlQ7165dU5KSkhRFUZSffvpJAZRnnnlGsba2zj+eZ8yYMUrHjh2VZs2aKd26dStwDlAmT55c6P65ubnKU089pQDK/PnzC5zbunWrAigREREFjj/66KOKi4tLgfePiopSzM3NlenTpxcYe+zYMQVQ9u/fX+D4+PHjFVtbW+Wxxx5TOnXqVOCcTqdT6tWrp0ycOFEBlPfffz//3J49exRA+e233wp9loSEBKV169aKmZmZcvr06QLnpk+frnh5eSlarTb/2M2bNwvd/16ffvqp4uDgoKSlpRU4HhkZqeh0uiKvEULcn8ysCFFN3Lp1i8WLF9OrVy+8vb0LPV55/PHHAVixYkX+saSkJNauXctTTz1VqvfSaDTMmzcPFxcXPvnkkwLnFixYQLt27fDz88s/lpuby5YtWxg+fDj29vb5x+vVq0ePHj1Yv359gXu0adMGf39/vvvuuyLf/6mnnuLAgQOcPXs2/9iuXbu4fPkyTz75ZKk+i5OTE99//z25ubl88cUX+cezs7NZvHgxTzzxRP6sSlRUFHXq1AFg1qxZ+Y+VJkyYkH/d6NGjSU5OZuXKlQXeJ2+mafr06Zw8ebJUMQpR00myIkQVlp6ezqpVqxg8eDDu7u689NJLODo6smrVKgICAgqMtbe3Z8SIEfz444/5x1asWIFarWbUqFGlfm9ra2t69erFpUuXuHr1KqD/Ab9r1y569OhRYOyFCxfIyMigRYsWhe7TokULzp8/T2ZmZoHj3bt3548//kApYmP4Xr16Ua9evQKfZfHixTz44IM0bty41J+lXbt2eHh48Ndff+UfO3z4MAkJCQU+i4eHB9u2bQPg6aef5uDBgxw8eJB33303f4y7uztNmzbl999/L/Ae8+fPZ8yYMaxdu5agoCD8/f354IMPOHfuXKnjFaKmkWRFiComJyeH33//ndGjR+Pq6sq4cePQ6XT8+OOPxMXFsWbNGkaMGIGFhUWha5966imOHDlCWFgYAD/++COPPvrofes87ievJuXatWsABAcHk5GRQVBQUIFxCQkJgH4W49+cnJxQFIVbt24VOB4UFER8fHyB2ZM8ebMZy5YtIzc3l8TERDZs2FDqGaJ71a1bN/9zABw8eDA/jjyWlpa0adMGAG9vbzp06ECHDh3w9fUtFPv+/fsLHAsICGD27NlcuHCBI0eOMGDAABYtWkSTJk1o27Ytn3/+OTExMWWOX4jqTJIVIaqQ06dP4+7uzuDBg4mPj+frr7/mxo0bbN68mdGjR1OrVq37Xt+tWzd8fX358ccfCQkJ4ejRo+X6Af/vWY+8H/aurq5FjlepVMXe69/n8u5R3A/wJ598khs3bvDHH3/wyy+/YGFhwaOPPlri2P+tqM+iUqlwcXEp9b1cXV2Ji4sjNze3yPPt2rXjs88+48qVK/z111906NCBjz/+mLp16/Lqq6+WKX4hqjMzUwcghCg5c3NzHBwcSExMJCkpiaSkJFJTU3F0dCzR9SqViieffJKvv/6azMxMmjRpQteuXcscT96SX09PTwAyMjIAsLKyKjDO2dkZuDvDcq/ExERUKlWhz5B3j7x7/lu9evXo2bMnP/74I1FRUTz22GPY2NiQnp5eps9y5cqV/M+R977m5uZoNJpS38vKygpFUcjMzLxvApmTk0NSUhK3b98mIyMDCwuLMs9yCVGdycyKEFWIv78/Fy9e5MCBA7Rr146PPvqIunXr0rVrV+bNm0dsbOx/3mPChAnEx8fz3XfflboY9V4ZGRns2rULX19fvL29AfJnIRITEwuM9fX1xdrampCQkEL3CQkJoVGjRoUSnLx73G9m46mnnmLTpk0EBweXa4boyJEjxMbGFugz4+LiQnZ2dpl6piQmJmJpaVlkopKbm8v27dt56qmncHNzY8iQIcTFxfHVV19x48aNAv1hhBB6kqwIUQV17NiRb775hmvXrrFt2zZ8fX2ZMWMGXl5ePPTQQ3z//ffFzkh4eXnx+uuvM2jQIMaPH1+m99dqtbz44oskJCTwxhtv5B/39/cH9AW19zIzM2PQoEGsW7eOlJSU/ONXrlxhz549DBs2rNB7XLx4EbVaXWBV0b8NHTqUoUOH8tRTT9GhQ4cyfZbExESef/55zM3NC/R2adq0aZGfxdLSEih+xicv9n8XOB88eJDnn38eDw8P+vXrx5kzZ5g1axYxMTHs2LGDJ598ssBKKSHEXfIYSIgqTKPR0KdPH/r06cN3333H77//zq+//sqUKVNo3749rVq1KvK6jz76qMTvcePGDQ4dOoSiKKSkpBAaGsqyZcs4deoUU6dOZeLEifljvb29adiwIYcOHeLll18ucJ9Zs2bRrl07Hn74Yd58800yMzN57733cHFxKbJO49ChQ7Rq1YratWsXG5uVlRVr1qwp8Wc5d+4chw4dQqfTkZCQwOHDh1m8eDHJycksW7aMZs2a5Y/Nm2U5dOhQgVVMdnZ21KtXj40bN9KzZ0+cnJxwcXGhfv36AOh0Oo4cOcLTTz9d4L3feust4uPjmTJlCk888QQNGjQocdxC1HimbPIihDCOpKQkJTU1VVGUu03hjh49et9rimsKl/dSq9WKvb290rx5c+XZZ59VDh48WOR93n33XaV27dpKZmZmoXPHjh1TevbsqdjY2Cj29vbKkCFDlPPnzxcal5KSotjY2CifffZZgeN5TeHup6imbXlN4fJeZmZmirOzs9KxY0fl7bffVqKiooq8V9euXZUBAwYUOr5r1y6ldevWiqWlpQIo48ePzz+3e/duBVCOHz9e4Jpr167dN24hRPFUilJEEwMhhCija9eu0aBBA5YtW1am/i2g75nyyiuvEB0dfd+ZFWNbu3Yto0aN4vLly3h5eZXomrFjx3Lx4sVCS5eFEGUnyYoQwuDeeOMN/vjjD4KDg/O7v5ZUbm4uAQEBjB8/nhkzZhgpwpJRFIVOnTrRpk0b5s2b95/jL1y4gL+/P3/++SddunSpgAiFqBmkwFYIYXDvvPMOw4cPL1OTs+joaMaMGVMp+o2oVCoWLVqEp6dn/q7L93PlyhXmzZsniYoQBiYzK0IIIYSo1GRmRQghhBCVmiQrQgghhKjUJFkRQgghRKVW5ZvC6XQ6rl27hp2d3X03SRNCCCFE5aHcaTTp6en5n6sGq3yycu3aNXx8fEwdhhBCCCHKIDo6On9/seJU+WQlb4fS6Oho2VdDCCGEqCKSk5Px8fEp0U7jVT5ZyXv0Y29vL8mKEEIIUcWUpIRDCmyFEEIIUalJsiKEEEKISk2SFSGEEEJUalW+ZkUIIYQoDUVRyM3NRavVmjqUak2j0WBmZmaQtiKSrAghhKgxsrOzuX79Ounp6aYOpUawsbHBw8MDCwuLct1HkhUhhBA1gk6n49KlS2g0Gjw9PbGwsJBmokaiKArZ2dncvHmTS5cu0bhx4/9s/HY/kqwIIYSoEbKzs9HpdPj4+GBjY2PqcKo9a2trzM3NuXz5MtnZ2VhZWZX5XlJgK4QQokYpz2/4onQM9Wct/48JIYQQolKTZEUIIYQQlZokK0IIIYTIN2HCBIYMGWLqMAqQZEUIIYSo5Lp3786UKVOMfk1lJclKMcKvJfPuhlA2BseYOhQhhBCiRpNkpRj7Im/y86HLLD902dShCCGEMBJFUUjPzjXJS1GUEsU4YcIE9u3bx1dffYVKpUKlUhEVFcW+fft44IEHsLS0xMPDgzfffJPc3Nz7XqPVann66adp0KAB1tbW+Pn58dVXXxnzj9ggpM9KMYa29uKT7REcjbpFVHwa9V1sTR2SEEIIA8vI0RLw3naTvHf4B32xsfjvH8NfffUVkZGRBAYG8sEHHwCg1WoZMGAAEyZMYNmyZURERDBx4kSsrKyYOXNmkdfUqVMHnU6Ht7c3q1evxsXFhQMHDvDss8/i4eHByJEjjfp5y0OSlWK4O1jRtXEd9kXeZO2Jq7zax8/UIQkhhKiBHBwcsLCwwMbGBnd3dwBmzJiBj48P8+bNQ6VS0bRpU65du8Ybb7zBe++9V+Q1oN+vZ9asWfn/3aBBAw4cOMDq1aslWamqhrfxZl/kTdadiGFqryao1dKWWQghqhNrcw3hH/Q12XuX1ZkzZ+jYsWOB7QI6d+5MamoqV69epW7dusVe+9133/HDDz9w+fJlMjIyyM7OplWrVmWOpSJIsnIffQLcsLMyI+Z2BocuJtCpkYupQxJCCGFAKpWqRI9iKhtFUQrta5RXA3O//Y5Wr17N1KlT+eyzz+jYsSN2dnZ88sknHD582KjxlpcU2N6HlbmGQS09AVhz/KqJoxFCCFFTWVhYoNVq8/87ICCAAwcOFCjSPXDgAHZ2dnh5eRV5DcDff/9Np06dmDRpEq1bt6ZRo0ZcuHChYj5EOUiy8h9GtPEGYGvodVIyc0wcjRBCiJqofv36HD58mKioKOLj45k0aRLR0dG89NJLREREsHHjRt5//32mTZuWvx/Pv6/R6XQ0atSIY8eOsX37diIjI3n33Xc5evSoiT/df5Nk5T+09nGkYR1bMnN0/BESa+pwhBBC1ECvvfYaGo2GgIAA6tSpQ05ODlu3buXIkSO0bNmS559/nqeffpp33nmn2GuuXLnC888/z7Bhwxg1ahTt27cnISGBSZMmmfCTlYxKKelC70oqOTkZBwcHkpKSsLe3N8p7zN97nrnbzvJAfSdWP9/RKO8hhBDCuDIzM7l06RINGjTAysrK1OHUCPf7My/Nz2+ZWSmBoa29UKngSFQilxPSTB2OEEIIUaNIslICHg7WdLmzEmjtCWm/L4QQQlQkSVZKKK/Qdu3xq+h0VfrJmRBCCFGlVFiyMmfOHFQqVYEdIBVFYebMmXh6emJtbU337t0JCwurqJBKpW8zd+ws7/RcuZRg6nCEEEKIGqNCkpWjR4+ycOFCWrRoUeD43Llz+fzzz5k3bx5Hjx7F3d2d3r17k5KSUhFhlYqVuYaHpeeKEEIIUeGMnqykpqYyevRoFi1aRO3atfOPK4rCl19+yYwZMxg2bBiBgYEsXbqU9PR0fv31V2OHVSYj2ugb7fwREktqVq6JoxFCCCFqBqMnK5MnT2bgwIH06tWrwPFLly4RGxtLnz598o9ZWlrSrVs3Dhw4UOz9srKySE5OLvCqKEF1a9PAxZaMHC1/hFyvsPcVQgghajKjJisrV67kxIkTzJkzp9C52Fh9gzU3N7cCx93c3PLPFWXOnDk4ODjkv3x8fAwb9H2oVKr8Qlt5FCSEEEJUDKMlK9HR0bzyyissX778vs13itqI6X6bML311lskJSXlv6Kjow0Wc0nk9Vw5fCmRKwnpFfreQgghRE1ktGTl+PHjxMXF0aZNG8zMzDAzM2Pfvn18/fXXmJmZ5c+o/HsWJS4urtBsy70sLS2xt7cv8KpIno739lyR2RUhhBBVk0qlYsOGDaYOo0SMlqz07NmTkJAQgoOD819t27Zl9OjRBAcH07BhQ9zd3dm5c2f+NdnZ2ezbt49OnToZKyyDyO+5ckJ6rgghhKjcZs6cSatWrQodv379Ov3796/4gMrAzFg3trOzIzAwsMAxW1tbnJ2d849PmTKF2bNn07hxYxo3bszs2bOxsbHhiSeeMFZYBtEnwJ1almZcvZXB4UuJdPR1NnVIQgghRKm4u7ubOoQSM2kH2+nTpzNlyhQmTZpE27ZtiYmJYceOHdjZ2ZkyrP9kbaHh4RYegDwKEkKIKk1RIDvNNK9S7COs0+n4+OOPadSoEZaWltStW5cPP/wQgDfeeIMmTZpgY2NDw4YNeffdd8nJyQFgyZIlzJo1i1OnTqFSqVCpVCxZsgQo/BgoJCSEhx56CGtra5ydnXn22WdJTU3NPz9hwgSGDBnCp59+ioeHB87OzkyePDn/vYzJaDMrRdm7d2+B/1apVMycOZOZM2dWZBgGMaKNNyuPRrM15DqzHmmGrWWF/lEKIYQwhJx0mO1pmvd++xpY2JZo6FtvvcWiRYv44osv6NKlC9evXyciIgLQP8lYsmQJnp6ehISEMHHiROzs7Jg+fTqjRo0iNDSUbdu2sWvXLgAcHBwK3T89PZ1+/frRoUMHjh49SlxcHM888wwvvvhifnIDsGfPHjw8PNizZw/nz59n1KhRtGrViokTJ5b/z+M+5CdsGbWpV5v6zjZEJaTzR2hsfh2LEEIIYUgpKSl89dVXzJs3j/HjxwPg6+tLly5dAHjnnXfyx9avX59XX32VVatWMX36dKytralVqxZmZmb3fezzyy+/kJGRwbJly7C11SdQ8+bNY9CgQXz88cf5C19q167NvHnz0Gg0NG3alIEDB7J7925JViqrvJ4rn+6IZM3xaElWhBCiKjK30c9wmOq9S+DMmTNkZWXRs2fPIs+vWbOGL7/8kvPnz5Oamkpubm6pV8qeOXOGli1b5icqAJ07d0an03H27Nn8ZKVZs2ZoNJr8MR4eHoSEhJTqvcpCdl0uh6FB3qhUcOhiItGJ0nNFCCGqHJVK/yjGFK/79BS7l7W1dbHnDh06xGOPPUb//v3ZsmULJ0+eZMaMGWRnZ5fqj+F+Pc7uPW5ubl7onE6nK9V7lYUkK+Xg5WhNpzsrgdadiDFxNEIIIaqjxo0bY21tze7duwud279/P/Xq1WPGjBm0bduWxo0bc/ny5QJjLCws0Gq1932PgIAAgoODSUtLK3BvtVpNkyZNDPNBykGSlXLKb79/Ilp6rgghhDA4Kysr3njjDaZPn86yZcu4cOEChw4dYvHixTRq1IgrV66wcuVKLly4wNdff8369esLXF+/fn0uXbpEcHAw8fHxZGVlFXqP0aNHY2Vlxfjx4wkNDWXPnj289NJLjB079r6NWiuKJCvl1LeZvudKdGIGR6MSTR2OEEKIaujdd9/l1Vdf5b333sPf359Ro0YRFxfH4MGDmTp1Ki+++CKtWrXiwIEDvPvuuwWuHT58OP369aNHjx7UqVOHFStWFLq/jY0N27dvJzExkXbt2jFixAh69uzJvHnzKuoj3pdKUUqx0LsSSk5OxsHBgaSkpApvvZ/njTWnWXUsmkfbePPJoy1NEoMQQoj7y8zM5NKlSzRo0OC+e9YJw7nfn3lpfn7LzIoBjGirfxS0NeQ66dm5Jo5GCCGEqF4kWTGAtvVqU8/ZhrRsLdtCY//7AiGEEEKUmCQrBqBSqRgedKfQ9ri03xdCCCEMSZIVAxkW5AXAgQsJXL0lPVeEEEIIQ5FkxUC8a9tIzxUhhBDCCCRZMaD8nivHr1LFF1kJIYQQlYYkKwbUL9AdWwsNVxLTORp1y9ThCCGEENWCJCsGZGNhxoDmHgCslUJbIYQQwiAkWTGwvEdBv0vPFSGEEMIgJFkxsHb1najrZENqVi7bw6TnihBCiMphyZIlODo6mjqMMpFkxcDUaum5IoQQovIZNWoUkZGRpg6jTCRZMYJ7e67E3M4wcTRCCCEEWFtb4+rqauowykSSFSPwcbKhQ0MnFAXWn5DZFSGEqKwURSE9J90kr9K2uNi2bRtdunTB0dERZ2dnHn74YS5cuABAVFQUKpWKdevW0aNHD2xsbGjZsiUHDx7Mv76ox0ALFizA19cXCwsL/Pz8+PnnnwucV6lU/PDDDwwdOhQbGxsaN27Mpk2byvaHXQ5mFf6ONcSINj4cupjImuNXmdyjESqVytQhCSGE+JeM3Aza/9reJO99+InD2JjblHh8Wloa06ZNo3nz5qSlpfHee+8xdOhQgoOD88fMmDGDTz/9lMaNGzNjxgwef/xxzp8/j5lZ4R/369ev55VXXuHLL7+kV69ebNmyhSeffBJvb2969OiRP27WrFnMnTuXTz75hG+++YbRo0dz+fJlnJycyvX5S0NmVoykf6A7NhYaohLSOX5Zeq4IIYQon+HDhzNs2DAaN25Mq1atWLx4MSEhIYSHh+ePee211xg4cCBNmjRh1qxZXL58mfPnzxd5v08//ZQJEyYwadIkmjRpwrRp0xg2bBiffvppgXETJkzg8ccfp1GjRsyePZu0tDSOHDli1M/6bzKzYiS2lvqeK2uOX2XN8au0rV9xGagQQoiSsTaz5vATh0323qVx4cIF3n33XQ4dOkR8fDw6nQ6AK1euEBAQAECLFi3yx3t46Pt+xcXF0bRp00L3O3PmDM8++2yBY507d+arr74qcOzee9ra2mJnZ0dcXFypYi8vSVaMaEQbb9Ycv8qW09d5f1AzrC00pg5JCCHEPVQqVakexZjSoEGD8PHxYdGiRXh6eqLT6QgMDCQ7Ozt/jLm5ef7XeeUHeUlNUf5doqAoSqFj994z75r73dMY5DGQET1Q3wnv2tbSc0UIIUS5JCQkcObMGd555x169uyJv78/t26Vr8TA39+ff/75p8CxAwcO4O/vX677GoPMrBhRXs+Vr3afY+2Jqwxp7WXqkIQQQlRBtWvXxtnZmYULF+Lh4cGVK1d48803y3XP119/nZEjRxIUFETPnj3ZvHkz69atY9euXQaK2nBkZsXI8hrE/XM+nmvSc0UIIUQZqNVqVq5cyfHjxwkMDGTq1Kl88skn5brnkCFD+Oqrr/jkk09o1qwZ33//PT/99BPdu3c3TNAGpFJKu9C7kklOTsbBwYGkpCTs7e1NHU6RRn1/kMOXEnm9rx+TezQydThCCFEjZWZmcunSJRo0aICVlZWpw6kR7vdnXpqf3zKzUgHyNjdcc/xqqZsACSGEEDWdJCsVoH9zD6zNNVyKT+PEFem5IoQQQpSGJCsVoJalGf2buwOw5niMiaMRQgghqhZJVipI3qOgLaeukZmjNXE0QgghRNUhyUoF6dDAGS9Ha1Kk54oQQpiU1A5WHEP9WUuyUozguGBe2v0S68+tN8j91GoVw+8ptBVCCFGx8jqxpqenmziSmiPvz/rfXXBLS5rCFePUzVPsvbqX6JRohjQaYpBdk4cHefH17nP8cz6e60kZeDiUbl8IIYQQZafRaHB0dMzf18bGxsYg39tFYYqikJ6eTlxcHI6Ojmg05dtuRpKVYgxrPIwFpxZwIekC/8T8Q1fvruW+Zz1nWx6o78SRqETWn4xhUnfpuSKEEBXJ3V2/2KGiN+KrqRwdHfP/zMtDkpVi2FnYMazxMH4O/5ml4UsNkqyAvtD2SFQia45f5YVuvpLVCyFEBVKpVHh4eODq6kpOTo6pw6nWzM3Nyz2jkseoycqCBQtYsGABUVFRADRr1oz33nuP/v37A/ppolmzZrFw4UJu3bpF+/bt+fbbb2nWrJkxwyqxMf5j+PXMrxy+fpiziWfxc/Ir9z0HtPDg/U1hXLyZxsno2wTVrW2ASIUQQpSGRqMx2A9SYXxGLbD19vbmo48+4tixYxw7doyHHnqIwYMHExYWBsDcuXP5/PPPmTdvHkePHsXd3Z3evXuTkpJizLBKzLOWJ73r9QZgWfgyg9yzlqUZ/QPzeq5Ioa0QQgjxX4yarAwaNIgBAwbQpEkTmjRpwocffkitWrU4dOgQiqLw5ZdfMmPGDIYNG0ZgYCBLly4lPT2dX3/9tdh7ZmVlkZycXOBlTOMCxgGw9dJW4tIN84wzr+fKZum5IoQQQvynClu6rNVqWblyJWlpaXTs2JFLly4RGxtLnz598sdYWlrSrVs3Dhw4UOx95syZg4ODQ/7Lx8fHqHE3r9OcINcgcnW5/Hqm+CSqNDo0vNNzJTOXHeE3DHJPIYQQoroyerISEhJCrVq1sLS05Pnnn2f9+vUEBAQQG6tvjObm5lZgvJubW/65orz11lskJSXlv6Kjo40aP8C4ZvrZldWRq0nPKf/6fLVaxbAgLwDWyqMgIYQQ4r6Mnqz4+fkRHBzMoUOHeOGFFxg/fjzh4eH55/+9GkZRlPuukLG0tMTe3r7Ay9i6e3enrl1dUrJT2HB+g0HuOTxI/yjo73M3iU3KNMg9hRBCiOrI6MmKhYUFjRo1om3btsyZM4eWLVvy1Vdf5a+7/vcsSlxcXKHZFlPTqDWMCRgDwM/hP6PVlb/OpL6LLe3q10anwPqTsrmhEEIIUZwKb7evKApZWVk0aNAAd3d3du7cmX8uOzubffv20alTp4oO6z8N9h2MvYU9V1Ovsid6j0HuOSK//X607FUhhBBCFMOoycrbb7/N33//TVRUFCEhIcyYMYO9e/cyevRoVCoVU6ZMYfbs2axfv57Q0FAmTJiAjY0NTzzxhDHDKhMbcxtG+Y0CDLeMeUBzD6zM1Vy4mUZw9G2D3FMIIYSoboyarNy4cYOxY8fi5+dHz549OXz4MNu2baN3b33vkunTpzNlyhQmTZpE27ZtiYmJYceOHdjZ2RkzrDJ7vOnjmKnNOBl3ktM3T5f7fnZW5vRrpn8ctvaEFNoKIYQQRVEpVfz5Q3JyMg4ODiQlJVVIse2Mf2aw6cIm+tTrw2fdPyv3/f45F8+YxYextzLjyIxeWJlLR0UhhBDVX2l+fld4zUpVl9ckbteVXcSklr8wtqOvM54OViRn5rLrjPRcEUIIIf5NkpVS8nPyo6NHR3SKjuXhy8t9P41axbCgvEJbeRQkhBBC/JskK2Uwvtl4ANadW0dydvnb/Q+/syror8ib3EiWnitCCCHEvSRZKYNOnp1o5NiI9Nx01kauLff9GrjY0qae9FwRQgghiiLJShmoVKr82pVfzvxCji6n3PfM67my9vhV6bkihBBC3EOSlTIa2HAgzlbO3Ei/wY6oHeW/XwsPLM3UnItL5fTVJANEKIQQQlQPkqyUkYXGgsebPg7A0rCl5Z4Nsbcyp1+gvueKFNoKIYQQd0myUg4j/UZipbHiTOIZjt04Vu775T0K2nTqGpk55d9/SAghhKgOJFkph9pWtXnE9xEAloWVvwV/J18X3O2tSMrIYfeZuHLfTwghhKgOJFkpp7EBYwHYe3Uvl5Iulete+p4rXoC03xdCCCHySLJSTvUd6tPduzsAP4f/XO775fVc2Rd5kzjpuSKEEEJIsmII45rplzFvurCJxMzEct3Lt04tguo6otUpbAiWnitCCCGEJCsG0NatLQHOAWRps1h9dnW57zeijQ+gXxUkPVeEEELUdJKsGIBKpWJ8gL4F/4qIFWRps8p1v7yeK5E3UgmJkZ4rQgghajZJVgykd/3euNm4kZiZyO8Xfy/XvRyszenTTN9zZa30XBFCCFHDSbJiIOZqc8b4jwH0y5jL+/gmr+fKxlPXyMqVnitCCCFqLklWDGh4k+HYmttyIekC+6/tL9e9ujRywc3ektvpOfwpPVeEEELUYJKsGJCdhR3DGg8D9C34y0Pfc0U/uyLt94UQQtRkkqwY2Gj/0ahVag5dP8TZxLPlutfwO8nK3sibxKVIzxUhhBA1kyQrBuZVy4ve9XoDsCy8fC34G7nWopWPvufKxpPXDBGeEEIIUeVIsmIEecuYt17aSlx6+epN8gptpeeKEEKImkqSFSNoXqc5Qa5B5OpyWRGxolz3GtTCEwszNWdvpBB2LdlAEQohhBBVhyQrRjIuQN+Cf/XZ1aTnpJf5Pg425vQJcAOk0FYIIUTNJMmKkXT36Y6PnQ/J2clsOL+hXPfKexS0IThGeq4IIYSocSRZMRKNWsPYgLEALD+zHK2u7ElG18Z18nuu7ImQnitCCCFqFklWjGiw72DsLeyJTolmb/TeMt9Ho1YxpLUXAGuOy07MQgghahZJVozIxtyGkX4jAVgaXr4mcSPu9FzZczaOmynl2yhRCCGEqEokWTGyx5s+jpnajJNxJzl983SZ79PYzY6WeT1XgmV2RQghRM0hyYqRudq4MqDBAKD8TeKk54oQQoiaSJKVCpC3jHnn5Z3EpJZ9VuSRFp5YaNRExErPFSGEEDWHJCsVwM/Jjw4eHdApOpaHLy/zfRxszOktPVeEEELUMJKsVJDxzfQt+NedW0dydtlnRfIeBW06dY3sXJ1BYhNCCCEqM0lWKkhnz840cmxEem466yLXlfk+XRu7UMfOksS0bPaclZ4rQgghqj9JViqISqXKr11ZfmY5ObqcMt3HTKNmWH7PFXkUJIQQovqTZKUCDWg4ACcrJ26k32BH1I4y32f4nUdBeyLiiE+VnitCCCGqN0lWKpClxpLHmz4OwNKwpWVeftzEzY4W3g7k6hQ2Bl8zZIhCCCFEpWPUZGXOnDm0a9cOOzs7XF1dGTJkCGfPni0wRlEUZs6ciaenJ9bW1nTv3p2wsDBjhmVSo/xGYaWx4kziGY7dOFbm++QV2q6VR0FCCCGqOaMmK/v27WPy5MkcOnSInTt3kpubS58+fUhLS8sfM3fuXD7//HPmzZvH0aNHcXd3p3fv3qSkpBgzNJOpbVWbR3wfAWBZWNmbxA2603Ml/HoyYdeSDBWeEEIIUekYNVnZtm0bEyZMoFmzZrRs2ZKffvqJK1eucPz4cUA/q/Lll18yY8YMhg0bRmBgIEuXLiU9PZ1ff/3VmKGZ1JiAMQDsvbqXS0mXynSP2rYW9ApwBWCtbG4ohBCiGqvQmpWkJP0MgJOTEwCXLl0iNjaWPn365I+xtLSkW7duHDhwoMh7ZGVlkZycXOBV1TRwaEB37+4A5WoSl/coaENwjPRcEUIIUW1VWLKiKArTpk2jS5cuBAYGAhAbGwuAm5tbgbFubm755/5tzpw5ODg45L98fHyMG7iRjGumX8a88cJGbmXeKtM9HmxcB5da+p4re6XnihBCiGqqwpKVF198kdOnT7NixYpC51QqVYH/VhSl0LE8b731FklJSfmv6Ohoo8RrbG3d2hLgHECWNotVZ1eV6R5mGjVDW3sCsPaEFNoKIYSoniokWXnppZfYtGkTe/bswdvbO/+4u7s7QKFZlLi4uEKzLXksLS2xt7cv8KqK7m0StyJiBVnasvVLyeu5svtMHAnSc0UIIUQ1ZNRkRVEUXnzxRdatW8eff/5JgwYNCpxv0KAB7u7u7Ny5M/9YdnY2+/bto1OnTsYMrVLoU78PbjZuJGYmsvXi1jLdo6m7Pc299D1XNp2SnitCCCGqH6MmK5MnT2b58uX8+uuv2NnZERsbS2xsLBkZGYB+dmHKlCnMnj2b9evXExoayoQJE7CxseGJJ54wZmiVgrnanDH++pVBy8KXlblJXF6hrbTfF0IIUR0ZNVlZsGABSUlJdO/eHQ8Pj/zXqlV3azSmT5/OlClTmDRpEm3btiUmJoYdO3ZgZ2dnzNAqjeFNhmNjZsP52+fZf21/me7xSEtPzDUqwq4lE36t6q2OEkIIIe7H6I+BinpNmDAhf4xKpWLmzJlcv36dzMxM9u3bl79aqCaws7BjWONhgL4Ff1nUtrWgZ1N9jY8U2gohhKhuZG+gSmBMwBjUKjWHrh/ibOLZ/76gCHmPgjYGx5CjlZ4rQgghqg9JVioBr1pe9K7XG9DXrpRFN786uNSyID41m31nbxoyPCGEEMKkJFmpJPKWMW+9tJW49NI3eDPXqBnSyguQQlshhBDViyQrlUSLOi1o7dqaXF0uKyIKN84rifyeKxE3SEzLNmR4QgghhMlIslKJjA8YD8Dqs6tJz0kv9fX+HvYEetmTo1XYFCybGwohhKgeJFmpRLr7dMfHzofk7GQ2XthYpnsMD9LPrqw9IcmKEEKI6kGSlUpEo9bkN4n7OfxntDptqe8xuJUX5hoVITFJRMRKzxUhhBBVnyQrlcyQRkOwt7AnOiWavdF7S329k60FDzV1BWCtFNoKIYSoBiRZqWRszG0Y6TcSKPsy5hFtfABYf/Ka9FwRQghR5UmyUgk93vRxzNRmnIg7QcjNkFJf392vDs62FsSnZvFXpPRcEUIIUbVJslIJudq4MqDBAACWhpe+Bb+5Rs1g6bkihBCimpBkpZLKaxK38/JOYlJLv7Inr/3+7jNx3JKeK0IIIaowSVYqKT8nPzp4dECn6PjlzC+lvj7A054AD3uytTo2n75mhAiFEEKIiiHJSiU2vpm+Sdy6c+tIyU4p9fV5syvyKEgIIURVJslKJdbZszO+Dr6k5aSxNnJtqa8f3MoTM7WK01eTOBtb+mRHCCGEqAwkWanEVCoV45rpa1eWn1lOji6nVNc717K823PlhMyuCCGEqJokWankBjYciJOVEzfSb7Azamepr8/b3HD9yRhypeeKEEKIKkiSlUrOUmPJ400fB/TLmBVFKdX1PfxccbK14GZKFn+fizdGiEIIIYRRSbJSBYz0G4mlxpLwhHCO3ThWqmstzNQMbuUJSKGtEEKIqkmSlSrAycqJR3wfAWBZWOlb8OetCtoZfoPb6dJzRQghRNUiyUoVMTZgLAB7r+4lKimqVNc283TAP6/nyinpuSKEEKJqkWSlimjg0IDu3t0B+Dn851JfPzzoTvv9E6XvhiuEEEKYkiQrVUjeMuaNFzZyK/NWqa4d0toLM7WKU9G3OXdDeq4IIYSoOiRZqULaurXF38mfLG0Wq86uKtW1LrUs6e6n77myRnquCCGEqEIkWalCVCpVfgv+FREryNJmler6vELb9Sek54oQQoiqQ5KVKqZP/T642biRmJnI1otbS3XtQ01dqW1jTlxKFn+fl54rQgghqgZJVqoYc7U5o/1HA7AsfFmpmsTpe67cKbSVnitCCCGqCElWqqDhTYZjY2bD+dvnOXDtQKmuvbfnSlJ66fYaEkIIIUxBkpXiKArcumzqKIpkb2HPsMbDAFgatrRU1zbztKepux3ZuTo2n5aeK0IIISo/SVaKc3k/fNUSVjwOl/7WJy+VyJiAMahVag5eP8jZxLMlvk6lUuXPrsijICGEEFWBJCvFuXIQUODsVlj6MHzfFYJXQG7pVuAYi1ctL3rV7QXoa1dKY3ArLzRqFcHRtzkfJz1XhBBCVG6SrBTnwddh8lFo+xSYWUNsCGx4Hr5sDvvmQprpV9PkLWPeemkrN9Nvlvi6OnaW9PCrA8Ca49LRVgghROUmycr91GkCD38B08Kh53tg5wGpN2DPh/BFM9j0EsSdMVl4Leq0oLVra3J1uayIWFGqa4cH3em5cvIqWl3lesQlhBBC3EuSlZKwcYKur8KUEBj2A3i2htxMOLEM5neAn4fCuV0mqWsZH6CfXVl1dhXpOeklvu4hf1ccbcy5kZzFP9JzRQghRCUmyUppaMyhxaMwcQ88uQ38B4FKDRf+hF+Gw7ft4diPkF3ypKG8uvt0x8fOh+TsZDZe2Fji6yzNNAxu6QlIoa0QQojKTZKVslCpoF5HGLUcXj4JHSaBhR3En4UtU/WPiHZ/AMnXjR6KRq1hjP8YAJaHL0er05b42hFtfADYHhZLUob0XBFCCFE5GTVZ+euvvxg0aBCenp6oVCo2bNhQ4LyiKMycORNPT0+sra3p3r07YWFhxgzJ8GrXh35z9HUtfeeAYz3ISIS/P9MX4657Fq4FGzWEIY2GYG9hz5WUK+y9urfE1wV62ePnpu+5skV6rgghhKikjJqspKWl0bJlS+bNm1fk+blz5/L5558zb948jh49iru7O7179yYlpQoup7Wyh46T9DMtI3+Guh1BlwOnV8HCbvDTADizGUox81FSNuY2jPQbCcCysJIvY1apVAxvI+33hRBCVG4qpTSby5TnjVQq1q9fz5AhQwD9rIqnpydTpkzhjTfeACArKws3Nzc+/vhjnnvuuRLdNzk5GQcHB5KSkrC3tzdW+GUTcwIOLYCwdaDL1R+rXR/aPw+tx4ClncHeKi49jr5r+5Kry+XXAb/SvE7zkl2XnEnHj/5Eq1PY/Wo3fOvUMlhMQgghRHFK8/PbZDUrly5dIjY2lj59+uQfs7S0pFu3bhw4UPx+N1lZWSQnJxd4VVpeQTB8kX4VUZdpYF0bbkXBtjfh8wDYPsNgLf1dbVwZ0GAAULomca72VnRrou+5slZmV4QQQlRCJktWYmNjAXBzcytw3M3NLf9cUebMmYODg0P+y8fHx6hxGoS9J/R6H6aGw8DPwbkxZCXDwXnwdStYPQ6uHCr30udxAeMA2Hl5J9dSS16Dktd+f92JGOm5IoQQotIx+WoglUpV4L8VRSl07F5vvfUWSUlJ+a/o6Ghjh2g4FjbQ7mmYfARGr4GGPUDRQfhG+LEvLHoIQtaAtmwrc/yc/Ojg0QGtomX5meUlvq6nvysO1ubEJmeyX3quCCGEqGRMlqy4u7sDFJpFiYuLKzTbci9LS0vs7e0LvKoctRoa94ZxG+CFg9B6LGgs4doJWPu0fgPFf76A9MRS3zpvdmXduXWkZJesUNnSTMPgVtJzRQghROVksmSlQYMGuLu7s3Pnzvxj2dnZ7Nu3j06dOpkqrIrnFgCD58HUMOgxA2xdITkGds3U92vZMg3iz5f4dl28uuDr4EtaThrrzq0r8XV57fe3h8WSnCk9V4QQQuhl5Rp+FWtpGTVZSU1NJTg4mODgYEBfVBscHMyVK1dQqVRMmTKF2bNns379ekJDQ5kwYQI2NjY88cQTxgyrcqpVB7pNh6mhMGQBuDWHnHQ4thjmtYFfRsLFvf9Z16JSqRjXTD+7svzMcnJ0JUs8Wng70Ni1Flm5On4/bfxmdkIIISq/kKtJdP14j8l7cRk1WTl27BitW7emdevWAEybNo3WrVvz3nvvATB9+nSmTJnCpEmTaNu2LTExMezYsQM7O8Mt6a1yzCyh1RPw/N8wfjM06Q+o4Nx2WDYYvusCJ5dDTmaxtxjYcCBOVk7EpsWyM2pnsePupVKp8gtt5VGQEEKIHK2O19ecIi4li22hxS98qQgV1mfFWCp1nxVDSbig79cS/It+tgXAtg60fVpfsFvLtdAlC04tYH7wfAKcA1g5cOV9i5bzxCVn0mHObnQK/PlqNxpKzxUhhKixvtl9js92RlLbxpyd07rhUsvSoPevEn1WRCk4+8LAT/Ut/Xt/APbekHYT9n2kr2vZMBluFNymYJTfKCw1loQnhHP8xvESvU2BnisnZHZFCCFqqnM3UvjmT3295PuDmhk8USktSVaqEuva0PkVeCUYRvwIXm1Bmw3By2FBJ1j6CERuB50OJysnHvF9BICl4UtL/BbDpeeKEELUaFqdwvS1p8nW6hjY2IbBUf+DlBsmjUmSlapIYw6Bw2Hibnh6JwQMAZUaLu2DX0fCt+3gyCLGNB4OwL7ofUQlRZXo1r383bC3MuN6UiYHLyQY7zMIIYSolJYciOLkldvUstTwmeVCVKdWwOqx5W5cWh6SrFR1Pg/AyKXwyino9BJYOkDCedj6Gg1/GEA3S3cUFH4O/7lEt7My1/BIfs+VKtRwTwghRLldSUjn0+1nAVjifwyr81tBYwH9P4YS1D4aiyQr1YVjXejzfzAtDPrPhdoNIPM246NOAbApcg23Lu4t0a1GtNFvYbBNeq4IIUSNoSgKb60/TUaOltHeN2kT+aX+RN/Z4NnapLFJslLdWNpB++fgpePw2AraurfDPyubTHSs3jQWFveBsA2gzS32Fi29HWjkWovMHB1bpeeKEELUCKuORrP/fAJu5unMyvoElS5HX2bQ7hlThybJSrWl1kDTAagm/M64tq8AsMLejqyrh+G38fB1azgwDzKTCl2qUqnyO9pKzxVhSvGpWVyKTzN1GEJUe7FJmXz4+xlAYZXbz5ilXNXP0D/ytUkf/+SRZKUG6NtqIm42biRoNGxtPRxsnCHpCuyYAZ8HwB9vQOLFAtcMbe2FWgXHLt8iSn5YCBOIvJFCj0/30vvzfRyLKv0+WUKIklEUhXc2hJKSlct7znuoH79Pv1/dyKVg5WDq8ABJVmoEc7U5o/1HA7BMuYUyJRQGfQ11mkJ2Khz+Dr4OgpWjIWo/KAruDlZ0bSw9V4QJZCaRGLqTXQvf5GPtp/xPvZA3Vx0lLav4R5dCiLLbcvo6u87coJ3ZOZ7MWKI/2G8OeLQ0aVz3kg62NURydjK9f+tNem463/X6js5enfXL0C78CYfmw/lddwd7tIQOk9mia8+Lq8LwdLDinzceQq02/VSgqGZyMuFGKMScgJjjcO0ESvw5VBT8tvRrbg9C2/yP2UObmyhQIaqnxLRsen++D21aAn85vId91g19a4zhi43++Kc0P7/NjBqJqDTsLewZ1ngYy88sZ2nYUn2yolJBo576182z+qTl1Eq4fgrWP8vAWu5ctOrOT0ndOXgxgc6NXEz9MURVptPq/57dSUqIOaHvvPyvzTZVwFXFhQhVIx4IaoPdifk8YbaHN48tYU/AG/TwK7y9hBCibGZtDiMxLZOVtRbpExUnX3j4y0pRp3IvmVmpQa6mXGXg+oHoFB1rBq3Bz8mv8KC0BDj+ExxZBKn6jasyFXNO1O5Hp9HvQp0irhHi3xQFbkXdTUpiTuiT4Jwi6p9snMGrDYpna5ZeduKbCHvSzGuz8tmOtPJxhL8/g90fkKWY8ZzZ//HFtKepbWtR0Z9IiGpn95kbPL30GM+bbeZNsxVgZgXP7AL3ipnBLM3Pb0lWaphX977Kjss7GOw7mP/r8n/FD8zNhrD1pP/1NTYJoXePN+oFHSaB70OVLvMWJpQap09Irt15nBNzAjKKKIq1qAUercCrNXi1Ac8gfY8glYpv95znk+1nUavg+7Ft6R3gpr9GUdCuHIPm7BauKU587fsDH43rWaEfT4jqJjkzhz6f/4V3SjCrLT9EjRYGfQVtJlRcDJKsiOKcvnma0VtHY6Y2Y8fwHdSxqXPf8YpOx9RPFtAvdT19Ncfu1hLU8YcOL0CLkWBuXQGRi0ojMxmuB99NSq6dhKQiuh2rzcE98G5S4hUELk30y+r/ZcPJGKasCgbgg8HNGNexfqH3zFzQHaukCxzQBhA/dCWPBNUz+EcToqZ4a10I24+Est16BnWUBGg+EoYtrNBfQqVmRRSrRZ0WtKrTiuCbwayIWMHLQS/fd7xKraZJ+748v60+D7tlMc/3CJz8GW6egc0vw+5Z0PYpfdMgO/cK+hSiwtxbAJs3axJ/Dvj37zgqfSLi1UaflHgG6RMVs//eqfXAhXheX6PvtDyxa4PCiQqAlT1WY1aQ/V13OhHOkk0ziPX9AXcHq3J/RCFqmgMX4ll5JIol5vP1iYpzY3j4i0o9Wy4zKzXQrsu7mLp3KvYW9uwcsRMbc5v7jo9NyqTTR7vRKbD3te7Ur5ULJ36Gw9/r+7WA/rfo5iP0j4g8WlTApxAGl1cAm19ncrzIAlgAHOrqH+V4BukTFI+WYFX6f3+RN1IYvuAAKZm5DGjuzrzHg+676iw3dANma8YD8I3zDF588XVUlfgbrBCVTUa2ln5f/cXA2yuYbr4KzKz1m+K6NavwWOQxkLgvrU7LoA2DiE6JZkb7GTzW9LH/vGbcj0f4K/ImLz/UiGl97hTZanMhYgscWgDRh+4OrtcFOk6CJv2KnPIXlYCiwO3L9ywZPgnXgosvgM1LSvJmTWrd//FhScQlZzJ0/gFibmfQtl5tlj/THivz//77cmvj29Q++S1piiV/dvmVQb17lTsWIWqKD38P59Q/W1lh+SEadPDIPAgaa5JYJFkR/+nXM78y58gc6trVZdOQTWj+I6nYGBzDKyuD8XK05u/pPQr/9nv1uH7pc/gG0N1p3lW7gb6updVosKxlnA8iSib1ZsHi12snID2h8DhzW/2GZffOmtwpgDWktKxcRi08SGhMMg1cbFn3QqeSr/DR5hLz7QC8Eg8TpXjAs3uo7+Vh0PiEqI6Co28zcf5Wtli8hZvqNrR8HIYsMNnjH0lWxH9Kz0mn15pepGSn8GWPL+lZ9/6rKzJztLT7v12kZOXy68T2dPItpudKUgwcWQjHl0Dmbf0xSwdoMw4eeA4cfQz6OUQR8gtgT9x9pHO/Ath7Z02KKYA1pFytjonLjrHn7E2cbS1YN6kT9ZxtS3UPXWo8CV90pI42jiMWHQia/jtmZlKCJ0RxsnK1DP56H2/feo8HNSH6DuYT/wSL0v3bMyRJVkSJfHn8SxaHLibINYil/Zf+5/i31oWw4sgVhgV58fnIVvcfnJ0Gp1boHxElnNcfU2nAfxB0nAw+D5T/AwjIzYLY0IKzJvGRFF8AG3R3dU4JC2ANSVEU3l4fyoojV7AyV7NiYgda161dpnvFnT2Ew68PY6nK4VD9F+gw4SMDRytE9fHFzkh0ez/mVfM1KOY2qCbuAdemJo1JkhVRInHpcfRd25dcXS4rBq4g0CXwvuOPX77F8AUHsDbXcPSdXtSyLMFvsjodnN8JB7+FS/vuHvdqq69r8R8MGvmNuER0Wn0icu+jnNjQYgpgfe7Wl3gF6XublKEA1tDm7z3P3G1nUanguzFt6NusfCvIjq3/hran3kGnqLjSfwn1OwwxTKBCVCMRscn83zffs9TsQzQqRf/op9UTpg5Lli6LknG1cWVAgwFsurCJpWFL+aTbJ/cdH1TXkYYutlyMT2NryHVGti3BIx21Gpr01b9iQ/UzLSGrIeYYrHkK7L3hgYnQZjxYl+037GpJUeD2lXta05/UP9rJTi08Nr8A9p5ZEwMUwBraxuAY5m47C8D7DweUO1EBaDPkRfZdOky35M04b5tEVsNALF0blfu+QlQXuVods1fv43OzeWhUCkqr0agqQaJSWjKzUsOdTTzLiM0j0Kg0bB22Fc9anvcdn9dl9IEGTqx+rmPZ3jQ1Do79CEd/gLSb+mPmtvpMv8ML4OxbtvtWZfkFsHc39Cu+ALZVwVkTx3qVuj8CwKGLCYxbfIRsrY5nujTgnYcDDHbvhNvJxHzZkxZEcsO6EW5T/waL+y/HF6KmWLQvkoBdE+isCSPHuSnmz+2pNP8+5DGQKJVndjzD4euHGRcwjtfbvX7fsdeTMuj00Z8oCvz1eg/qOpfjL31OJoSugYPzIS7szkGVfslzx0lQv2ul/yFcJlkp+mXC986a5PWruVeBAtg7syYVUABraOfjUhg2/wDJmbn0D3Tn2yfu30ulLP46dgr/zQ9TR5VMfIPBuIxbWj3/7ghRCpfi09jy9Uu8pF5LjsYa8+f/gjpNTB1WPklWRKn8ffVvJu2ehK25LTtH7MTOwu6+48cuPszf5+J5pWdjpvY2wF98RdHXsxycD+e23z3u1lw/09J8RIUXghpMgQLYvA6w/1EAm7c6x60ZmFftDq1xKZkM/VbfSyWoriO/TuxQol4qZbFg6VImXpyCmUpHVu85WHaeZJT3EaIq0OkUZs+bz9sJM1CrFJShC1G1HGXqsAqQZEWUiqIoDN04lAtJF3it7WuMbzb+vuPzeq5417bmr9eL6LlSHvHn4fACCP4VctL1x2xd9e382z0NtsUsma4M8gtg72lNf78CWM/Wd5cMV5ICWENKy8rlsYWHCIlJor6zDesmdcbJiLslp2Tm8OOn03kl90e0aNA8uQXqdTLa+wlRma3de4wH9wyljiqZ1GajqfXofFOHVIgkK6LU1kauZebBmbjburN12FbM1ebFjs3I1vLAh/qeKysmdqCjr7PhA0pPhBNL4fBCSLmmP6ax1G+c2GESuBmu5qFM8gpg85cM36cA1tqpYPdXryCo5VrhIVekXK2O534+zu6IOJxsLVj3Qifquxi/n8OhC/HcWDKWwZoDZFm5YDnpH7CXhnGiZrmWmELMV31opwonsVYTnF75q1JuOCvJiii1LG0Wfdb0ITEzkbkPzqV/g/73Hf/WutOsOBLN8CBvPhvZ0niBaXMgfKN+6fO1E3ePN+yhT1oa9dKvODK2ewtg8/43Pb7wuLwC2HtnTapAAawhKYrCuxtDWX7oCpZmalY824GgMvZSKYu5m44z6Nh4/NXR5Hi2xfypP8DMeDM6QlQmiqKw+cvJPJL0Cxkqaywm/Y2mTmNTh1UkSVZEmSw4tYD5wfNp5tyMFQNX3HeDuGNRiYz47iA2FhqOzuiFbUl6rpSHokD0YX1L/zObQdHpj7s0gfbP69tGG6rCPa8A9t5Zk+IKYN2aFZw1qeNX5QpgDe27fRf46I8IVCpYMLoN/QIrdjfuzBwtz371G9+kTMVBlY7S9hlUD39WoTEIYSr7t6+i44HnUKsUYnt/i3vnMaYOqViSrIgyScxMpM+aPmRps/ip70+0dW9b7FhFUejx6V6iEtL59NGWjGjjXXGB3rqsb+l/YhlkJeuPWdeGNk/qe7bY33/5dQG5WXAj9E7x651Zk5tnKboAtvHdPibVpADW0DadusbLK04C8N7DATzVpYFJ4giNSeKL+fNYZPYJ6krUBEsIY0q4FoVqYRecSCHEYzjNn/vR1CHdlyQrosxmHZzFmsg1dPfpzjcPfXPfsfP+PMenOyLp0NCJlc+WsedKeWSlwMlf9AW5t6L0x9Rm0GyYfhWRV1DB8TotxJ8r2Jr+RihoswvfO78A9k5iUg0LYA3t8MUExt7ppfJU5wa8N8i0dUXf7D5H7p9zmGq+FkVjherp7fpHdEJUR9pczn/ag0YZp7mgaUjd6fsxt6wc/VSKI8mKKLOLSRcZvGEwKlRsGrKJ+g71ix0bczuDLh/re678Pb0HPk4m+oeh08LZP/SPiC7vv3u8bkdo/qg+kbl2Uv8qtgD2nu6vNaAA1tDOx6UyfMEBkjJy6NfMnW9HB6ExcC+V0srV6nh0wX4m33iXXpqTKA4+qJ77C2ycTBqXEMZwceV0GkZ8T4pizfVR22kSYMRaQgORZEWUy4u7X2Tf1X2M8hvFOx3eue/YMT8c5p/z8Uzp1ZgpvSpBs6FrwfqkJXQt6HILny9QAHsnQalhBbCGdjMli6Hz93P1Vgat6zqywoi9VErr4s1URn29jd9Ub1NffUNfmD1mbY2vKxLVS1roH9iueQyADY1nM2T0ZBNHVDKSrIhyORp7lKe2P4WVxoqdI3biaOVY7NgNJ2OYsioYHydr9r1m4J4r5ZF8HY4ugugjBWtNpADWoNKz9b1UTl9Nop6zDete6IRzrcrVwG/ZwSh+2fQH6y3ex0aVBV2mQa/3TR2WEIaRdJW0rzthq01ig/kA+k1fXml+Wfgvpfn5XQFrPkVV09atLf5O/mRqM1l1dtV9x/Zt5k4tSzOiEzM4GpVYQRGWgL0H9HwPJmyBh7+A1mP0vVkkUTGYXK2Ol349yemrSdS2MWfJkw9UukQFYEz7erg2CuKNnIn6A/98rl9RJkRVp80h6eex2GqTOK1rgM9jnxs0UdHqtETeimT9ufUExwUb7L5lIcmKKESlUjGu2TgAVkSsILuoAtQ7rC00DGyub7q15vjVColPmJ6iKMzaHM7uiDgszdT8ML4dDSqg6VtZqNUq5o5owT6LB1mUO0B/cP0LcDPStIEJUU7ZO2bhEH+CZMWaPwM/po1v2RsgKopCdEo02y5t49OjnzL+j/F0XNGR4ZuG896B9/j94u8GjLz0KkWyMn/+fBo0aICVlRVt2rTh77//NnVINV7f+n1xtXElITPhP/+SjmirX7b8e8h10rKKqBMR1c6ivy/y86HLqFTw5ahWtKlXcU3fysLDwZr/DQnko9zHOazzh+wUWDVav6JMiKro7DYsDutXbH5s8RITB/cs1eXxGfHsi97Ht8Hf8vyu53lw1YMMWDeA1/96naXhSzkRd4KM3AxszGxo596Oho4NjfEpSszInbz+26pVq5gyZQrz58+nc+fOfP/99/Tv35/w8HDq1q1r6vBqLHO1OaP9R/PF8S9YFr6MIY2GFNskrm292tRztuFyQjrbQmMZXpE9V0SF23L6GrO3RgDwzsAA+jc3Xjv7jec3EpEYQePajfFz8qORYyMsNWV71PRIS092hN1gcsjL/GH9DnXiI2HDCzDyZymwFlXL7Why1z6LGfBTbl/6PvHsfRtzpmanEp4QTmhCKKHx+tf1tOuFxpmpzWhauynNXJoR6BJIc5fm1Levj6YSPD43eYFt+/btCQoKYsGCBfnH/P39GTJkCHPmzPnP66XA1niSs5Pp/Vtv0nPT+a7Xd3T26lzs2G92n+OznZF0bOjMimc7VGCUoiIduZTImB8Ok63VMaFTfWY+0swo76MoCt+c/IZFIYsKHNeoNNS3r08Tpyb41fajqVNT/Jz8cLEu2QaXiWnZ9P3yL7xSw1hj9T/MlBzoNRO6TDXCpxDCCHKz0f3UH3XMMYJ1DVnRbBEfj7rbwDNbm83ZxLMFEpNLSZdQ/tXoUoWKhg4NCyQmTWo3wUJTcVtTlObnt0lnVrKzszl+/DhvvvlmgeN9+vThwIEDRV6TlZVFVlZW/n8nJycbNcaazN7CnmGNh7H8zHKWhS+7b7IyNMiLz3ZGcvBiAldvpeNdu3I3IxKld+FmKhOXHSNbq6NPgBvvPmycpm86Rcecw3NYeXYlAP3q9+NW5i3O3jrL7azbXEi6wIWkC/xx6Y/8a5ysnPCr7Yefkx9NajehqVNT6jvUL7Qhp5OtBR8Pb85TS7J4N3scc8wXw+4P9E3/fHsY5fMIYVC7Z6GOOUaSYsO7Fq/ywYO12HB+Q35icvbWWXKLaNvgaetZIDHxd/KnlkUtE3yAsjFpshIfH49Wq8XNza3AcTc3N2JjY4u8Zs6cOcyaNasiwhPAaP/R/BrxKweuHSDyViRNahfdS8W7tg2dfJ05cCGBdSdieLln5dw4S5TNzZQsJvx0hKSMHFr5OPLVY62N0vQtR5fDe/vfY8vFLahQ8U6HdxjpNxLQz7bEpcdx9tZZIm9FEpEYwdnEs1xOvkxiZiIHrx/k4PWD+fcyV5vTyLERTWo3wc/JLz+ZeaipG48/4MOKIw/Ryeoyg7S7YM1T8Nw+cJRHz6JyUhSFa6eXE3pqMaG1HVlr0ZDsWl8zbkdGobGOlo4EugTmJybNnJvhbO1sgqgNx+Q1K0ChWghFUYqtj3jrrbeYNm1a/n8nJyfj4+Nj1PhqMm87b3rW7cnOyztZFraM/+vyf8WOHdHGmwMXElhz/CovPdTovhshiqojPTuXZ5YeJToxg3rONiwe3xZrC8M/w87SZvHavtfYG70XM5UZH3b5kAENB+SfV6lUuNm64WbrxoPeD+Yfz8jN4Pyt85y9dZaziWfzk5m0nDTOJJ7hTOIZuHD3fdxt3fF1aEIdHytev90KF6srtE+LRL1qDDy1HcytDf7ZhCitxMzE/NmS0PhQwm6eJjE7Cdzq3BkRDzqwNrMmwDmAQOdAAusEEugciFctr2r3/dekyYqLiwsajabQLEpcXFyh2ZY8lpaWWFpWvl4O1dn4ZuPZeXknv1/6nVeCXqGOTZ0ix/ULdOfdDaFcSUznaNQtHmggbc2rOq1O4eUVwZy6moSjjTk/TWhnlF4qaTlpvPznyxyJPYKlxpLPu39eICG5H2sza5rXaU7zOs3zj+kUHTGpMUQmRhZIYmJSY4hNiyU2LRZqgXkteBaw1vnQODsWvzWP4Nf6afycmtKkdhNszOVxpjC+tJw0fQFsXmKSEEZMakyhcWaKgneWmivZD/Bq97509G5FQ4eGlaIA1thMmqxYWFjQpk0bdu7cydChQ/OP79y5k8GDB5swMnGvlnVa0qpOK4JvBrMiYgUvB71c5DgbCzMGtvBg9bGrrDkeLclKFafvpRLGrjM3sDBT88O4tjSsY/hn3Lczb/PCrhcITQjF1tyWbx76hnbu7cp1T7VKjY+dDz52PvSsd3dJZ0p2CpG3IjmbqJ992XMpmITsy2SoczltZcnp7Fg4/GH++Lp2dfPrYPIeI3nYelS731pFxcnR5hB5K5KQ+JD8xOTC7QuFCmABGjg00M+YxF0kMHI3blnmDMuczXuP9mJ485q16tLkq4FWrVrF2LFj+e677+jYsSMLFy5k0aJFhIWFUa9evf+8XlYDVYxdl3cxde9UHCwd2DF8R7G/cR6+mMCohYewtdBw9J1e2FhUiieNogwW/XWRD7eeQaWCb58IYoARlijHpcfx3M7nOH/7PI6WjnzX6zuauRhnhVFRsnK1DPrmL87fimKsx25qZ//FWUtLIh3cicu+XeQ1dhZ2+YmLX20/mjg1KdeSalF96RQdUUlRhCaEEnIzhLCEMCISI8jR5RQa627rrk9M7tSaBDgHYGdhB+GbYPVYAJ7Kfg1to74sebJdtUiYq8xqIIBRo0aRkJDABx98wPXr1wkMDGTr1q0lSlRExenh0wPvWt5cTb3KpgubeKzpY0WOa1ffibpONlxJTGd7WCxDW9es7L+6+P30dT7cegaAGQP8jZKoRKdEM3HHRGJSY3C1dmVhn4X4Ovoa/H3ux9JMw5ej2jD42wyWXn6Mvxtk4XN9B6RpSJywibNZCfkzMWdvneXi7YukZKdw7MYxjt04ln8fjUpDA4cGhYp5S7qkWlR9iqIQmxarnzFJCCUsPoywhDDSctIKjXWwdCiQmAS6BBb9dyXxEmx8EYDvch/msFk7dgxrXi0SldIy+cxKecnMSsX55cwvfHTkI+ra1WXTkE3FPif9atc5vtgVSedGzvzyjPRcqWqORSXyxA+Hyc7V91J5f1CAwb85nrt1jud2PsfNjJv42PmwsPdCvO1Ml9jO33ueudvO4mqZw37n/8M88RzU6wzjNoLm7vLnbG02F5Mu3i3kvVMTczvrdpH3dbJy0veCuTMD41fbr8gl1aLquZV5S19jcicxCYkPITGz8P5o1mbW+Dv5301MnAPxtvP+739TuVmwuA9cD+aE0oSRWe/w/uCWjO1Y3zgfyARk12VhFOk56fRa04uU7BS+6vEVD9V9qMhx0YnpdJ27B5UK/p7eQ3quVCEXb6YybMEBbqfn0DvAje/GtDH4EuWQmyG8sPsFkrKSaFy7Md/3+r7You2KotUpjPz+IMcv32Jo3TQ+vz0NVXYKdJgE/e7fnFJRFG6k3ygwA5O3pLqoOoS8JdX3zsA0qd0EB0sHY308UU7pOemEJ4QTlhBGaHwoIfEhRRfAqsxoXLsxzVya5S8Z9nX0xUxdhocYv78GRxeRoranT/qH+NRvzMpnO1Sene0NQJIVYTRfHv+SxaGLCXINYmn/pcWOe3zhIQ5eTODV3k14SXquVAnxqVkMm3+AK4nptPRxZOXEDgZfonz4+mFe/vNl0nPTaVGnBfN7zq80P6Sj4tPo/9XfZORo+eGBWHqdvtMiYdgP0OLRUt8vb0l1xK2I/ILevCXVRXG3dS9QC+Pn5IePnQ9qVaXYwq3GyNHlcO7WufyVOSHxIVxMuohO0RUaW9++foHEpKlTU6zMrMofRNh6+G0CABOyp3NQHcQfr3Q1SoG7KUmyIozmRtoN+q3tR66Sy4qBKwh0CSxy3NrjV3n1t1PUd7Zhz2vda+Qz1qokI1vLY4sOcSr6Nj5O1qyf1BkXAy9R/vPKn7y+73Wyddm092jP1z2+rnRLg5cfusw7G0KxMFNz6IEDOJ34Gsys4Zld4F703/XSyFtSfe8MTOStyCJ/Swf9I4TGtRvnby3QpHYTWVJtQDpFx+XkywX6mUQkRpCtK7zTvKuNK4HOgTSvo09Mmrk0w97CCD9zEi7A990gO4UfGML/ZY7kzf5Neb5bxdZzVQRJVoRRvf3322y+uJl+9fvxSbdPihyTlpVLuw93kZ6tZc3zHWlbX5YxV1ZancILy4+zI/wGjjbmrH2hE74G/g1u84XNvLv/XbSKlp51ezL3wbkVugdJSSmKwoSfjrIv8iYtPWuxvvZXqC/shtoN4Nk9YG2c3aXzllRHJEbkP046f/s8WdqsQmNVqPCx8ymwtYBfbT/cbd3ll4L7yHtcV6DRWkIYqTmphcbaWdgVKoB1tXE1fpA5mbC4F8SGcM4qkH633yDAy4n1kzphpql+M2ySrAijOpt4lhGbR6BRadg6bCuetTyLHPfab6dYc/wqj7Xz4aPhLSo4SlES+l4q4Sw5EIWFmZpfnmlPOwMnlisiVjD78GwAHvF9hFmdZpXtGX4FuZGcSZ8v/iIpI4fpD7oyKfJpuH0ZGveBx1eBumJ+aOTqcrmSfIWzt87qtxa4U9B7M+NmkePtLewLrUbydfStsUuqk7KSCiQmoQmhxGfEFxpnqbEsWADrEkhdu7qmSfy2TIVjP5JlUZsHk/9HgtqZTS92IcCzev5sk2RFGN0zO57h8PXDjAsYx+vtXi9yzKGLCTy28BC1LM04OqOXUVq0i/L54e+L/N/v+iXK3z4RxMAWhluirCgKi0IW8c3JbwD9PlPT202vEjUYm09d46UVJ9GoVWx91B6/34dDbiZ0exN6vGXS2BIzE/MfH+UlMZduXyJXKbx5XU1ZUp2Rm8GZhDMFEpPolOhC4zQqDY0cGxVITHwdfSvH6qyQNbD2aRRUvKSewZb0AF5+qBHT+viZOjKjkWRFGN1fV/9i8u7J2JrbsnPETn3zon/R6RS6fbqH6MQMvhzViiGtvUwQqSjO1pDrTP71BIoCbw9oyrMPGu6ZuKIofHbsM5aG64uwX2j5Ai+0fKFKPaZ4acVJNp+6RsM6tmzrfg2LzS/oTzy+Cvz6mTa4f/n3kuq8/03KSipyvLOVc4Gmdk1r63eprswzXnlydDmcv3We0IS7syYXbl9Aq2gLja1rV7dAYtLUqSnWZpVw76f487CwG2SnstN5LBNj+tPYtRZbXu6CpVn1/SVPkhVhdDpFx9CNQ7mYdJHX2r7G+Gbjixz35a5Ivtx1ji6NXFj+TPsKjlIU5/jlRB5fpO+lMq5jPWY90sxgiYRWp+WDQx+w7tw6AKa3m87YgLEGuXdFup2eTd8v/+JGchYTOtVnptkSOLIQLB309SvOlbvgsbRLqi3UFvg6+laqJdU6RceV5CsFEpOIxIgia3nqWNfJX5kT6BxIM5dmlWal2X3lZMAPveBGKLdcH6DtlZfQqTSsfaETQXWNUyNVWUiyIirE2si1zDw4Ew9bD7YO21rkb2X39lzZ/8ZDeDpWwt9qaphL8WkMm7+fW+k59PJ35fuxbQ3WSyVHm8Obf7/Jjss7UKvUzOw4k6GNh/73hZXUvsibjP/xCAC/TGhN5/1PQvQhcA3QrxCysDVxhKWXnpPO+dvnC6xGOpt4lvTc9CLHe9h6FGhqZ8wl1TfSbhRoshaWEEZKdkqhcXbmdgS4BOQnJoEugbjZFr35baW36WU4sRSdTR0G584hJNmGp7s04N2HA0wdmdFJsiIqRJY2iz5r+pCYmcjcB+fSv0H/IseN+v4ghy8l8npfPyb3aFTBUYp7JaRmMWzBAS4npNPS24EVz3Yw2P5N6TnpTNs7jf3X9mOuNmfug3PpVa+XQe5tSu9sCGH5oSt4OFix/Zkm2C/tCak3IHA4DF8MVejRVnF0io6YlBh9AlPCJdX3bu5YliXVSVlJhCWE3U1M4sOIy4grNM5CbUFT56b5vUyauzSnrn3dKlH79J9Or4Z1EwEVP/l+wawwV+o62bBtStcasa+aJCuiwiwIXsD8U/Np5tyMFQNXFPko4bdj0by+5jQNXGz589VuVapuoTrJyNby+KJDBN/ppbLuhc7UsTPMSpHk7GRe3P0iJ+NOYm1mzZc9vqSTZyeD3NvU0rNz6f/V31xOSGdYay8+75gJSwaCLhf6zoaOk00dotEkZyfnbymQNwNz7ta5IvuQqFBR175ugSQmb0l1ljaLiMSI/CZrYQlhXE6+XOgeapUaX0ffAolJo9qNKkcBrKHdjISF3SEnjastXqLLkY4A/PpMezo1ql4F0MWRZEVUmMTMRPqs6UOWNosl/ZbQxq1NoTH39lxZ+0JH2tSTnisVTatTmPTLcbaH3cDBWt9LpZGrYXqpJGQk8Pyu54lIjMDOwo75PefTyrWVQe5dWRy/nMij3x1Ep8CC0UH0T98Mf7wOKo1+/6AGXU0dYoXJ1eVyOfny3TqY/1hSbWduR3puepEFsN61vPWJyZ1ak6ZOTWtGw7vsdPihJ8SFo63/IH1uTuFCQiaPP+DDnGE1p81Dldp1WVRtTlZODPIdxJrINSwNW1pksmJraUb/QA/WnrjKmuNXJVkxgf/7PZztYTew0KhZNK6twRKV66nXeXbns0QlR+Fk5cTC3gvxc6p+Sy3b1HPi+W6+zN97gbfXh9BmylhcY47D6ZX6tujP/QUONWO1m5naDF9HX3wdfRnAgPzjCRkJBTZ3zFtSnZKjrzlxtnIukJg0c26Go5WjiT6Fif3xOsSFg60r39Z+gwsRSbjZW/LWAH9TR1ZpycyKKLeLSRcZvGEwKlRsHrqZevb1Co05eCGBxxcdws7SjKPv9MLKvPoux6tsFv9zif9tCQfgm8dbM6hl0U38SutS0iWe3fkssWmxeNh6sKjPoiL/v68usnN1DP52P2euJ9OzqSs/PBGA6se+EBsCXm3gyT/ArGY2YCtOtjabS0mXsLewlw67eYJXwIbnQaXmQv9f6L1eQafAD+Pa0iugihYJl1Fpfn5XgwolYWoNHRryoPeDKCj8HP5zkWPaN3DCy9GalKxctofFVnCENde20Ov83+/6ROWt/k0NlqicSTjDhG0TiE2LpYFDA5b1X1atExUACzM1X4xqiYVGze6IOFafSoBRy8HKEWKOwx9vmDrESsdCY4Gfkx8etTwkUQGIi4Df9Rtk5j74JpMP1EKnwKCWnjUuUSktSVaEQYwP0PdZ2Xh+I7czbxc6r1arGN7GG4A1x69WZGg11vHLt3hlZTCKAmM61OXZBxsa5L4nbpzgqe1PkZiZiL+TP0v6LcHd1t0g967smrrb82qfJgB8sDmcaMUVRiwGVHD8JzixzLQBisorOw1+Gw856dCwOwu0g4mITcHJ1oKZg6r/MuXykmRFGEQ793b4O/mTqc1kdeTqIscMD9I/0//nfDzXkzIqMrwaJyo+jYnLjpGVq6NnU1dmDjJM07d/Yv7huZ3PkZqTSpBrEIv7LsbJqmbVID3TtSHt6tcmLVvLq6tPoW3YEx6aoT/5+2v6WRYh7qUo8PurcDMCarlzoeuXfL3nIgDvDwrA2cA7nFdHkqwIg1CpVIxrNg7Qb1yXrS28tLGesy0PNHBCUWDdiaL7N4jyS0zLZsJPR0hMy6a5lwPfPNHaIDu2bo/azkt/vkSmNpOuXl35rvd3RW6zUN1p1Co+e7QVthYajkQlsvifi9DlVfAbCNosWDUO0gpvmCdqsOBf4NQKUKnRDvuBV7deI0er0LOpK48Y6NFsdSfJijCYvvX74mrjSnxGPL9f/L3IMSPuPApae/wqVby2u1LKzNHyzNKjRCWk413bmsUT2hqkudS6c+uY/td0cnW59Kvfj696fFU591ipIHWdbfI7jH66PZKzcWkwdAE4N4Lkq7DmSdAW3lhQ1EA3wvUzbgA9ZvBTjBfB0bexszTjw6HNpZanhCRZEQZjrjZntP9oAJaFLysyGRnQ3ANrcw0X49M4GX27giOs3rQ6hSkrgzlx5Tb2VmYsebIdrnZW5b7v0rClvH/gfXSKjhFNRvBR148w11TDJl2lNKqdDw81dSVbq2PqqmCyzexg1C9gbguX/oLds0wdojC1rFR9nUpuBvj25HLAc3y64ywAbw/0x92h/P8+awpJVoRBjWgyAhszG87fPs/BawcLna9laUb/QH0xphTaGtbsrWfYFhZ7Ty+V8j2iURSFr098zafHPgXgqcCneK/De2jUsuwc9I8+PxrenNo25oRfT+br3efAtSkMma8fcOBrCFtv2iCF6SgKbJkK8ZFg54ky9HveXBdGZo6Ojg2deaydj6kjrFIkWREGZW9hz7DGwwBYGr60yDF5j4I2n7pGZk7hrpai9H7af4nF/1wC4JNHW9C+oXO57qdTdMw+PJtFIYsAeCXoFaa2mSpT1v/iamfFh0ObAzB/73lOXLkFzYZAp5f1AzZM1i9XFTXPiaUQslrf5XjEj6wMz+DgxQSszNV8NFwe/5SWJCvC4Eb7j0atUnPg2gEib0UWOt+hobO+50pmLjvCb5ggwuple1gsH9xp+vZGv6YMblW+Tqo5uhxm/DODlWdXokLFO+3f4Znmzxgi1GppQHMPhrb2QqfAtFXBpGfnQs/3ocGDkJMGq0ZDZpKpwxQVKTYEtk7Xf93zXa47tmL272cAeK2PH/Wcq95u3aYmyYowOG87b3rW7QlQZJM4tVqVv4xZHgWVz4krt3h5xUkUBUa3r8vz3crXSyVLm8W0vdPYcnELGpWGOV3nMKrpKANFW33NfKQZ7vZWRCWkM2drBGjMYMRPYO8NCedh/Qug05k6TFERslL0WzBos6BxH5ROL/PO+lBSsnJp5ePIk50bmDrCKkmSFWEU45vpm8T9fvF34jMKL+McFqR/FPTPuZvEJmVWaGzVxeWENJ5Zqu+l8lBTV2Y9Ur5eKmk5aUzaNYm90Xux1FjyVY+vGNhwoOECrsYcrM355FH9BnQ/H7rMvsibYOsCo34GjSWc/R3++czEUQqjUxTY/Io+QbX3gqHfs+l0LLsj4jDXqJg7ogUatTz+KQtJVoRRtKzTklZ1WpGjy+HXM78WOl/fxZZ29WujU2D9Sem5Ulr6XipHSUzLJtDLnm8eL18vlduZt5m4YyJHYo9ga27Lgl4L6ObTzYARV39dG9dhfEf9lgPT15wiKT0HvIJg4J0k5c8P4fwuE0YojO7YjxC6FtT6mbUEnS2zNusf0b7YozFN3GpeXyJDkWRFGE1ek7jVkatJz0kvdH5Efvv9aOm5UgqZOVomLjvGpfg0vByt+XFCO2wty95LJS49jie3P0lIfAiOlo4s7rOYdu7tDBhxzfFmf38authyIzmL9zaF6g8GjYU2EwAF1jwNt6JMGKEwmuunYNtb+q97vg912zNrcziJadk0dbfjhe6+po2vipNkRRjNQz4P4V3Lm6SsJDZf2Fzo/IDmHliZq7lwM41g6blSIjqdwrTVwRy/fMsgvVSiU6IZ98c4zt8+j6u1K0v6LaGZSzMDRlyzWFto+HxUKzRqFRuDr7Hl9DX9if5z9TszZ96GVWMgu3DyLqqwzCRYPV5fp9KkP3R6iV3hN9h06hpqFcwd0QILM/lxWx7ypyeMRqPWMCZgDAA/n/kZnVKwwNDOypz+gR6AFNqW1Jw/zrA1RN9LZeG4tjQux7Ty+VvnGf/HeGJSY/Cx82Fp/6X4Ospvf+XVyseRyXd+i35nQyg3kjPBzBJG/gw2LvqVIlum6usbRNWnKLDpZbh1CRzqwpD5JGflMmNDCAATuzakhbejaWOsBiRZEUY1tNFQ7CzsuJx8mb3Rewudl54rJbdk/yUW/X23l0qHcvRSCbkZwoTtE7iZcZNGjo1Y2m8p3nbehgq1xnvxocYEetlzOz2HN9ae1j/mdPCCR5fo+26cXglHFpk6TGEIR3+A8A36OpVHfwIbJ+ZsPcON5CzqO9swtXcTU0dYLUiyIozKxtyGR5s8Cujbtv9bx4bOeDpYkZyZy64z0nOlODvCYpl1p5fK6339ytVL5cj1Izyz4xmSspJo4dKCJf2WUMemjqFCFYCFmZrPR7bCwkzN3rM3WXEkWn+iQVfo8z/919vfgiuHTBekKL9rJ2H72/qve38A3m05cD4+///vj4e3wMpcOj4bgiQrwuieaPoEZiozTsSdIDQ+tMA5tVqVv4xZHgUVLTj6Ni+v1PdSefyBukwqR6Henit7eGHXC6TnptPeoz2L+izCwdLBgNGKPE3c7Jje1w+A//s9nMsJafoTHSZBs2Ggy4XV4yAl1oRRijLLuH2nTiUbmj4MHSaRnp3Lm+v0j3/GdKhb7k7S4i5JVoTRudm60b9BfwCWhS0rdH74nUdBf0Xe1D/fF/muJKTz9JKjZObo6OFXh/8NLnsvlc0XNjN171Syddk85PMQ3/b8FhtzGwNHLO71VOcGtG/gRHq2lmmrT6HVKaBSweB54BoAqTf0P/Bys00dqigNRYGNk+H2ZXCsq///U6Xisx2RXElMx9PBijf6NTV1lNWKJCuiQuQtY95xeQfXU68XONfAxZa29aTnyr/dSstmwk9HSLjTS2XeE0Fl7qWyImIFb//zNlpFyyO+j/BZ98+w1FgaOGLxb2q1ik8fbUktSzOOX77Fwr8u6k9Y2MKo5WDpANGHYMcM0wYqSufw9xCxBdTm+jok69qcuHKLH/fra8o+HNYcOyvZmdyQJFkRFaKpU1Pau7dHq2hZfmZ5ofN3e65clZ4r3O2lcjGvl8r4svVSURSFhacXMvvwbEC/b9P/Ov8PM3XZ+7KI0vFxsuG9QQEAfL7zLGeuJ+tPOPvCsIX6r48shOAVJopQlMrV47DjHf3XfT8ErzZk5Wp5Y81pFAWGtfaih5+raWOshiRZERUmb3Zl7bm1pGanFjg3oIUHlmZqzselcvpqzd70TadTeHX1KY5dvoWdlRk/PdkOV/vS91JRFIXPj3/ONye/AeCFli/wRrs3UKvkn31Fe7SNN7383cjRKkxdFUxW7p2Vb379oNub+q+3TNE3FhOVV8Yt/b4/uhzwfwQeeBaAb/88z7m4VFxqWfDuwwGmjbGaMup3rQ8//JBOnTphY2ODo6NjkWOuXLnCoEGDsLW1xcXFhZdffpnsbHl+Wx118epCQ4eGpOWksfbc2gLn7K3M6RfoDkih7UfbIvg95DrmGhXfj21TphbdWp2WWQdnsSRsCQDT201nUqtJsi29iahUKuYMa46zrQURsSl8sfPc3ZPd3oDGfSA3U98wLj3RdIGK4ikKbJgMSVegdv38OpUz15OZv/cCALMeCaS2rYVp46ymjJqsZGdn8+ijj/LCCy8UeV6r1TJw4EDS0tL4559/WLlyJWvXruXVV181ZljCRNQqNWMDxgLwy5lfyNXlFjif9yhoUw3uubLsYFR+XcMnI1rSydel1PfI0eYw/a/prD23FrVKzQedPsj/cxemU8fOkg+HNgfg+78ucDTqTlKiVusfB9WuD7evwNpnQFcz//5Xage/1W9IqbGAR5eClQO5Wh3T15wmV6fQt5kbA5q7mzrKasuoycqsWbOYOnUqzZs3L/L8jh07CA8PZ/ny5bRu3ZpevXrx2WefsWjRIpKTk40ZmjCRQb6DcLJy4nradXZe3lngXCdfFzwcrEjKyGH3mTgTRWg6O8NvMHNTGACv9WnCkNal76WSkZvBS3teYsflHZipzfi026cMbTzU0KGKMuoX6M7wIG8UBV5dfYq0rDsJu3VtGPULmFnDhd2wZ7ZpAxUFRR+FXe/rv+47GzxbAfDDP5cIiUnC3sqM/w0OlJlLIzLpw+uDBw8SGBiIp6dn/rG+ffuSlZXF8ePHi7wmKyuL5OTkAi9RdVhqLHnM7zFA3yTu3mJajVrFsCD9D+g1x6NNEp+pnIq+zUsrTqBT4LF2Pkzu0ajU90jOTua5nc+xP2Y/1mbWfPvQt/Su19sI0YryeP+RADwdrLiSmM6HW8/cPeEeCI/o64v4+1M4s8U0AYqC0hPv1Knk6vvjtHsGgIs3U/liZyQA7zwcUKa6MlFyJk1WYmNjcXNzK3Csdu3aWFhYEBtbdKOkOXPm4ODgkP/y8fGpiFCFAY1qOgoLtQVhCWGciDtR4NzwOw3i9kXeJK6G9FyJTkzn6aX6XirdmtThf0NK/xtaQkYCT29/mpNxJ7GzsGNh74V08upkpIhFedhbmfPpoy0B+PXwFfacvWcWscWj+qZxAOufh/hzRdxBVBidDja8AMlXwakhDPoKVCp0OoU314aQlauja2MXHm0jW1UYW6mTlZkzZ6JSqe77OnbsWInvV9Q3ZUVRiv1m/dZbb5GUlJT/io6uWb+BVwdOVk4M8h0EFG7B37BOLYLqOqJTYENw9e+5cjs9m/E/HSE+NZsAD3u+HR2EeSl7qVxPvc6EbROISIzAycqJn/r+RCvXVsYJWBhEp0YuPNm5PgDT15zmVto9iwp6fwD1OkN2ir7gNivFNEEKOPgNRG4DjeWdOhV7AH45fJkjUYnYWGiYPbS5PP6pAKVOVl588UXOnDlz31dgYGCJ7uXu7l5oBuXWrVvk5OQUmnHJY2lpib29fYGXqHryljHvjd7L5eTLBc6NaKOfLavuPVfye6ncTMPTwYqfnmxHrVL2UolKimLctnFEJUfhYevBsv7L8HPyM1LEwpDe6NcU3zq23EzJ4p0NoXf/rmvuNBqz84CbEfpOqdX430GldeUQ7Jql/7r/R+DRAoCY2xl89EcEANP7+uHjJF2gK0KpkxUXFxeaNm1635eVVcme3XXs2JHQ0FCuX7/b0XTHjh1YWlrSpk2b0oYmqpCGDg150PtBFBR+Dv+5wLmBd3quRN5IJSSmevZc0ekUXvvtFEejbmFnacZPTz6AWymfeZ9JOMP4beOJTYulvn19lvVfRj37ekaKWBialbmGL0a1QqNW8XvIdTadunb3ZC1XGLlM3yE1fCMc+Np0gdZEaQnw25OgaCFwBLR5EtDP+r+9LoS0bC1t69VmXMf6po2zBjFqzcqVK1cIDg7mypUraLVagoODCQ4OJjVV3xCsT58+BAQEMHbsWE6ePMnu3bt57bXXmDhxosyY1ADjA8YDsPH8Rm5n3s4/7mBtTt9m1bvnysfbI9hy+m4vFT/30vVSORl3kqe3P01iZiL+Tv4s7b8Ud1tZNlnVtPB25KWH9MXU724IJTbpnjotnweg/8f6r3fNhIt7Kzy+Gkmng/XPQco1cG4Eg77U7+eEfjuQfZE3sTBT89HwFqjV8vinohg1WXnvvfdo3bo177//PqmpqbRu3ZrWrVvn17RoNBp+//13rKys6Ny5MyNHjmTIkCF8+umnxgxLVBLt3NvR1KkpmdpMVkeuLnAur+fKxuBrd7t9VhM/H7rM9/v0vVQ+Ht6CTo1K10vln5h/eHbHs6TkpBDkGsTivotxsnIyRqiiAkzu0YgW3g4kZ+by+ppTBR99tn0KWo0BRQdrnoLbUqNndPu/hPM7wcxKX6diqf9F4mZKFh9sCQfglZ6NaeRay4RB1jxGTVaWLFmCoiiFXt27d88fU7duXbZs2UJ6ejoJCQl88803WFrKBms1gUqlYlyAvnZlRcQKsrV3iww7N3LB3V7fc+XPatRzZfeZG7y/MRSAV3s3YVhQ6VYRbI/azkt/vkSmNpMuXl34rvd32FmUvsOtqDzMNWo+H9kKSzM1f5+LZ/mhe2q4VCoY+Cl4tIL0BFg9FnJqxio5k7h8AP78P/3X/efql5Pf8f6mUG6n5xDgYc+zDzY0UYA1l2wSIkyqX4N+uNq4Ep8Rz9ZLW/OPa9Qqhub3XKkej4JOX73Ni7+eRKfAqLY+vPhQ6XqprDu3jul/TSdXl0u/+v34usfXWJtZGylaUZEaudbizf5NAfhw6xkuxafdPWluDaN+BmsnuHYStr4qBbfGkHpTP3ulaKHFKAgal39qW+h1tobEolGrmDuiRalX7Inykz9xYVLmanNG+48GYFn4sgJT4Hk9V/ZG3iQupWr/NhmdmM5TS46SkaPlwSZ1+L+hpeulsjRsKe8feB+domNEkxF81PUjzDWyBX11Mr5jfTr5OpOZo2Pa6mBytbq7Jx3rwojFoFLDyeVwfInJ4qyWdDpY/yykXAeXJjDw8/w6laT0HN7dqO8s/Xy3hgR6OZgy0hpLkhVhcsMbD8fazJpzt85x8NrB/OONXGvRuq4jWp3CxpPX7nOHyu12ejYT7vRS8few59snWpf4NzNFUfjm5Dd8ekxfx/Vk4JO81+E9NGqNMUMWJqBWq/jk0ZbYWZpx8sptvr+zR1Q+34eg53v6r7e+rm8BLwzj78/gwp/67Q4eXQqWd+tR/vd7ODdTsvCtY8tLDzU2YZA1myQrwuQcLB0Y1ngYoJ9duVdeoW1V7bmSlavl2Z+Pc+FmGh4OVvw0oR12ViWbEdEpOuYcmcPC0wsBeCXoFaa1mSYNqKoxL0drZj7SDIAvdkYS+u+l+52ngP8joMuB1eMgtfrUc5nMpb9h7529mAZ+Bm4B+af+irzJmuNXUalg7ogWWJnLLwmmIsmKqBTG+I9BrVKz/9p+zt2622L84RaeWJipOXsjhdCYqrUPlL6XymmOXEq800ulHe4OJeulkqPLYcY/M1gRsQIVKt5p/w7PNH/GyBGLymBYkBf9mrmTq1OYtjq44A7kKhUMma9/VJFyTd8LRJtb/M3E/aXGwdqn9autWj4BrUfnn0rLyuWtdSGA/hFdm3qy4s6UJFkRlYK3nTc96/YECs6uOFib0ydA38147YmqVWj7yY6zbD51DTO1iu/GtqGpe8l6B2Vps5i2dxpbLm5Bo9Iwp+scRjUdZeRoRWWhUqn4cGggLrUsiLyRyud3NsvLZ2mn36HZwg4u/3N3N2BROjotrH0GUm9Anab6VVf3+GT7WWJuZ+Bd25rX+0pXaFOTZEVUGnnLmH+/+DvxGfH5x/MeBW0IjqkyPVeWH7rMgr0XAPhoeAs6l7CXSlpOGpN3TWZv9F4sNZZ81eMrBjYcaMRIRWXkXMuSOcP07d0X/X2RwxcTCg6o0wSGLtB/fXAehKyp4Airgb8+gUv7wNxGX6diYZt/6lhUIksPRgEwZ1hzbEu5DYYwPElWRKXRyrUVLeu0JEeXw4qIFfnHuzaug5u9JbfTc9gTUfmf0e8+c4P37vRSmdqrSX6y9V9uZ95m4o6JHI49jI2ZDQt6LaCbTzdjhioqsd4Bboxs642iwKu/nSI161+Pe/wHQZdp+q83vQQ3wio+yKrq4l7Y+5H+64e/ANem+acyc7RMX3saRYFH23jTtXEd08QoCpBkRVQq45vpW/CvPruajNwM4E7PldZ3C20rs5CrSfm9VB5t483LPUvWSyUuPY4ntz9JSHwIjpaO/Nj3R9q5tzNytKKye/fhALxrW3P1Vgb/2xxeeMBD70DDHpCTDitHQ8btCo+xykm5AWsnAgq0HgstHytw+uvd57h4M406dpa8MzCg6HuICifJiqhUHvJ5CK9aXtzOus2m85vyj49oo28Qt+fsTW6mZJkqvPuKTkznqaX6XipdG7swe1jJto6PTolm/B/jOX/7PK7Wrizpt4RmLs0qIGJR2dlZmfPpoy1RqWDVsWh2hd8oOECtgRE/gkNduHVJv6eNTlf0zcSdOpWnIS0OXAP0XWrvERqTlL9k/P+GBOJgI72MKgtJVkSlolFrGBswFoCfz/yMTtF/423kakdLnzs9V4JjTBlikZLSc3hyyVFupmTR1N2O+aODStRL5fyt84z/YzxXU6/iY+fD0v5L8XX0rYCIRVXRoaEzz3RpAMCb606TkPqvZN3GSd/h1swKIrfpazFE0fZ+BFF/g7ntnToVm/xTOVod09ecRqtTGNjcI38zVVE5SLIiKp2hjYZiZ2HH5eTL7Ivel3+8svZc0fdSOcb5uFTc7a346cmS9VIJuRnChO0TuJlxk0aOjVjabynedqXbK0jUDK/28aOxay3iU7OZsT608N9/z1b62guAvXMgcnuFx1jpnd99N5Eb9JW+SPkeC/+6SPj1ZBxtzPN73YjKQ5IVUenYmNvwaJNHAVgavjT/+CMtPLHQqImITSHsWuXouaLTKbz+22kOX0qk1p1eKh4O/71fz5HrR3hmxzMkZSXRwqUFS/otoY6NFPKJolmZa/hiVCvM1Cq2hcWy/mQRs4utnoB2zwAKrJsICRcqPM5KK/k6rHsWUKDNBGjxaIHT5+NS+GqXvr/T+4MCqGMnm+lWNpKsiErpiaZPYKYy4/iN44TF61c5ONiY07uZvudKZSm0/XTHWTbd6aWyYEwQ/h7/3Utlz5U9vLDrBdJz02nv0Z5FfRbhYCn7jYj7C/Ry4JWe+nbv728M49rtjMKD+s4B7wcgMwlWjYXstMJjahptrn6DwvR4cGsO/T4qeFqnMH3NabK1Orr71WFIKy8TBSruR5IVUSm52brRr0E/QL+JX568R0Ebg2PIzjVtIeGvh68w/04vlTnDmpdoiePmC5uZuncq2bpsHvJ5iG97fouNuc1/XicEwAvdfWnl40hKVi6vrzmFTvevx0FmFjByGdi6QlwYbHpZdmjeOxuuHACLWjByqX4X63ssOxjFiSu3sbXQMHtoyYriRcWTZEVUWnnLmHdc3sH11OsAdG3kQh07S26l5/CnCXuu7ImI4907vVRe6dmYR9v6/Oc1KyJW8PY/b6NVtDzi+wifdf8MS41MN4uSM9Oo+XxkS6zM1ew/n8CyO43LCrD30P9QVptB6Bo4/F2Fx1lpnNul36QQ4JGvwblg8Xp0Yjpzt50F4M0B/ng6/vcjXGEakqyISqupU1Pau7dHq2j55cwvgP6b9bDW+mlaU7XfD41JYvKvJ9DqFIYHeTOl1/13YlUUhUWnFzH7sH6ztCeaPsH/Ov8PM7V0xRSl17BOLd4e4A/AnD8iOB+XWnhQvU7Q50P919tnQNQ/FRhhJZEUo6/dAWj7NAQOL3BaURTeWhdCRo6WBxo4MfqBuiYIUpSUJCuiUhvXTN+Cf+25taRm678pD7/zKGhPRBzx/17GaWRXb6Xz5JKjpGdr6dLIhTn/0UtFURQ+P/45X5/8GoDnWz7Pmw+8iVol//RE2Y1pX4+ujV3IytXx6upgcrRFPBJt/xw0HwmKFn6bAMnXKjxOk9Hm6OtUMhLBvQX0nV1oyG/HrvLP+XgszdR8PLwFarU8/qnM5DumqNS6eHWhgUMDUnNSWXtuLQBN3Oxo6e1Ark5hY3DFfQNOysjhyZ/u6aUyJggLs+L/CWl1WmYdnMWSsCUAvN72dSa3mizPxEW5qdUq5o5ogb2VGaeuJjF/TxErf1Qq/RJdt+aQdhNWj4PcytlQ0eD+/B9EHwJL+zt1KgV3O7+RnMn/ftd3BJ7WuwkNXGyLuouoRCRZEZWaWqXO3+DwlzO/kKvT749yb8+VipCVq+W5n49x7p5eKvb36aWSo81h+l/TWXtuLWqVmg86fZA/SySEIXg4WPO/IYEAfPPnOU5fvV14kIWNvmGclQNcPQrb3qzYIE0hcjvs/0r/9SPfgFPDAqcVReGdDaGkZObSwtuBp+803BOVmyQrotIb5DsIJysnrqddZ9flXfpjLfU9V85cTybsWpJR319RFN5Yc5pDF/W9VH6ccP9eKhm5Gby05yV2XN6BmdqMT7t9ytDGQ40ao6iZHmnpycDmHuTqFKatPkVmThG7kjs1gOGLARUc+xFOLq/wOCvM7Wj9lgMADzwLzYYUGvJ7yHV2ht/ATK3i4+EtMCtBp2lhevL/kqj0LDWWPOan32xsadhSFEXB0caCXgGuAKw9btz2+5/tiGRD8DU0ahXzRwcR4Fl8L5Xk7GSe2/kc+2P2Y21mzbcPfUvver2NGp+ouVQqFf8bEkgdO0vOx6XyyfazRQ9s3Bt6zNB/vWUaXDtZcUFWlPw6lVvg2Rr6/F+hIbfSsnl/o75v06QejUrUF0lUDpKsiCphpN9ILNQWhCaEciLuBHD3UdAGI/ZcWXHkCvP2nAdgztDmPNik+F4qCRkJPL39aU7GncTO3I6FvRfSyauTUeISIo+TrQUfD28OwOJ/LnHgQnzRA7u+Ck36gzZL3zAuLaECo6wAu2bC1SNg6QCPLgGzwm0BPtgSTkJaNk3cajG5h+zBVZVIsiKqBGdrZwb5DgJgWdgyAB5sXAeXWpYkpmWz96zhe67sPRvHOxv0vVRe7tmYke2K76VyPfU6E7ZNICIxAicrJ37q9xOtXFsZPCYhivJQUzcev7P09vXfTpOcmVN4kFoNw74HJ19IioY1T+q7u1YHEVvh4Dz910O+hdr1Cw3ZExHH+pMxqFXw8fAWWJppKjZGUS6SrIgqI6/Qdk/0Hi4nX9b3XAnS91wxdKFtaEwSk3/R91IZFuTF1Pv0UolKimLctnFEJUfhYevBsv7L8HPyM2g8QvyXdwb6U9fJhpjbGXywObzoQVYOMGq5ftfhS/v0q2aquluXYcPz+q87TAL/QYWGpGTm8Pb6EACe6tyA1nVrV2SEwgAkWRFVRkPHhnT16oqCws/hPwMwPEj/KOjPiDgSDNRzJeZ2Bk8tOUpatpZOvs58NKxFscuNIxIjGL9tPLFpsdS3r8+y/suoZ1/PIHEIURq2lmZ8NrIlKpU+ed8eFlv0QLcAGHxnFmL/lxC+scJiNLjcbP0MUWYSeLWBXrOKHPbRHxFcT8qkrpMNr/aRXySqIklWRJWS14J/4/mN3M68jZ+7Hc299D1XNp0qf88VfS+VI8SlZOHnZsd3Y9sU20vlZNxJntr2FImZifg7+bO0/1Lcbd3LHYMQZdWuvhPPPqhfqvv2upDimyYGDoNOL+m/3jAJbhZTmFvZ7XofYo7rZ4xG/KTfG+lfDl5I4JfDVwD4aHhzrC3k8U9VJMmKqFIecH+Apk5NydRm8lvkb4Dheq5k5+p4/ufjRN5Ixc3e8r69VPbH7OfZHc+SkpNCkGsQi/suxsnKqVzvL4QhTOvdhKbudiSkZfPWuhCU4jYy7DkT6neF7FRYORoykys0znI7sxkOzdd/PeQ7qF14RjMjW8tb604D8PgDdenk61KREQoDkmRFVCkqlSq/duXXiF/J1mbzSEtPzDUqwq4lE36tbN9wFUXhzbWnOXgxAVsLDT9OaFfspmY7onbw4p8vkqnNpItXF77r/R12FnZl/kxCGJKlmYbPR7bCXKNiZ/iN4pN4jZl+NsLeCxLOwYYXQGfancxLLPESbJis/7rji9B0QJHDvtgVSVRCOu72Vrw1oGkFBigMTZIVUeX0q98PV2tX4jPi2XppK7VtLejl7waUfXPDL3ZGsu5kjL6Xypg2NPN0KHLcunPreP2v18nV5dKvfj++7vE11mayU6uoXAI87ZnauwkAszaHc/VWetEDa9WBkT+DxgIitsD+LyowyjLKzdLvdZSVBN7toNfMIoedir7ND39fBODDoYH37TgtKj9JVkSVY64x5wn/JwBYFr4MRVHu9lw5GVP0pm73seroFb7+U99LZfbQQLoV00tladhS3j/wPjpFx/DGw/mo60eYa+QboKicnnvQlzb1apOalctrv51CpyvmcZB3Gxjwqf7r3f+D87srLsiy2PEuXA8G69r6maEi/g1m5+p4Y+1pdAoMbuVJzzu/zIiqS5IVUSWNaDICazNrzt06x8HrB3mwib7nSkJaNnvP3izxffZF3uTt9fpeKi891IhR7QpvE68oCt+c/IZPj+m/oT8Z+CTvd3wfjVoK9UTlpVGr+OzRllibazh0MZEf918qfnCb8RA0HlBg7dP65cCVUdgGOPK9/uuh34Nj0b2P5u89T0RsCk62Frz3cEDFxSeMRpIVUSU5WDowrPEwQN8kzlyjZkgrTwDWlrDQNuxaEpOWH9f3UmntxbQ70+b30ik65hyZw8LTCwF4JegVprWZJjsniyqhvostMwb6AzB3+1nO3UgpfvCAT8AzSN+uftUYyMmooChLKPEibLqzgqnzK9Ckb5HDzsam8O2drtMzH2mGc63CnWxF1SPJiqiyRvuPRq1Ss//afs7dOsfwO4+CdkfcIDEt+77XXrunl0rHhs58NLxwL5VcXS4z/pnBiogVqFDxTvt3eKb5M0b7PEIYw+j2denWpA7ZuTqmrg4u/jGpmaV+h2YbF4g9DVumQnEriSpaTiasHg9ZyeDTAR56t8hhWp3C9LWnydEq9PJ3Y1ALjwoOVBiLJCuiyvKx86Fn3Z4A/Bz+M/4e9gR62ZOjVdgUXPzmhsmZOTz501FuJGfRxK1Wkb1UsrRZTNs7jS0Xt6BRaZjTdQ6jmo4y6ucRwhhUKhVzR7TAwdqc0JhkvrlTn1UkB2949CdQqeHUCjj6Q8UFej/b39YnUNZOMOLHIutUAH785xKnom9jZ2nG/w0JlBnQakSSFVGl5S1j3nJxC/EZ8Yy409F2TTGrgrJzdbyw/Dhnb6TgamfJT08+gIN1wW98aTlpTN41mT3Re7BQW/Bljy8Z2HCgcT+IEEbkZm/F/w0JBODbPecJjr5d/OAGD0LvD/Rfb3sTrhw2foD3E7oWji3Wfz1sITh4FTksKj6Nz3bqm9vNGOiPu4NVRUUoKoAkK6JKa+XaipZ1WpKjy2FFxAoeaeWFuUZFaEwyEbEFe64oisKb606z/3wCNnd6qXj9q5dKUlYSE3dM5HDsYWzMbPiu93d09+legZ9ICOMY1NKTQS090eoUpq0KJiNbW/zgji9Cs6Ggy4XV4yDlRsUFeq+EC7DpFf3XXaZB495FDtPp9P+2M3N0dPJ1ZtR9Nh0VVZPRkpWoqCiefvppGjRogLW1Nb6+vrz//vtkZxesJbhy5QqDBg3C1tYWFxcXXn755UJjhLifvNmV1WdXY22p5aGmrkDhQtsvdp1j3Ql9L5VvRwcR6FWwl0pcehwTtk0gJD4ER0tHfuz7I+3c21XMhxCiAvxvcDPc7C25GJ/Gx9siih+oUsEj86COP6TGwm/jQVvETs7GlJOhr1PJToF6naHHjGKHrjwazaGLiViba+67l5eouoyWrERERKDT6fj+++8JCwvjiy++4LvvvuPtt9/OH6PVahk4cCBpaWn8888/rFy5krVr1/Lqq68aKyxRDfWs2xOvWl7czrrN5gubGdFG/1vV+pPX8osJVx+L5uvd5wD4vyGB9PBzLXCP6JRoxv8xnvO3z+Nq7cqSfkto5tKsYj+IEEbmaGPB3BEtAVhyIIp/zsUXP9iyFjz2C1jaw5WDsOOdCoryjm1vwo0QfcHv8MX6jrtFuJ6UweytZwB4ra8fdZ1tKjJKUUGMlqz069ePn376iT59+tCwYUMeeeQRXnvtNdatW5c/ZseOHYSHh7N8+XJat25Nr169+Oyzz1i0aBHJyUW3Tc/KyiI5ObnAS9RsGrWGsQFjAX2h7YNNnHG2tSA+NYu/Im/yV+RN3l6n3x5+cg9fHn+gYC+V87fOM/6P8VxNvYp3LW+W9l+Kr6NvhX8OISpCtyZ1GNNB/2/g9TWnSMq4z4yJs6++nwnA4e/g1KoKiBA4/RscXwKo9HUq9kWv6lEUhRnrQ0nNyqV1XUcmdKpfMfGJClehNStJSUk4Od3d7O3gwYMEBgbi6emZf6xv375kZWVx/PjxIu8xZ84cHBwc8l8+PvJsUsDQRkOxM7cjKjmKA9f+ZkhrfRHe13+eZ9IvJ8jVKQxu5clr/9oePjQ+lAnbJ3Az4yaNHBuxrP8yvO28TfERhKgwbw/wp56zDdf/v717j4qy3PcA/p0LM8M9FOUieEMCFUGFzPGK0ka8hlvNVh2y086T5Y3yEtpOXbU9aLF0185taR73zXPQ3Yha5tlQCuRJMxBFvN9FgdBEQLABhuf8MUIRw8AgA/MO389a71o4/J53nt/6sZpf7/vO85T+hDX7TpsPDp4EjFlu/PnzxUBhrnUnd/uC8X0AYMxSoF9Uk6H7Thbg4LliqBRyvDcjFAo5b//Yq3ZrVi5fvow//elPmDdvXv1rRUVF8PJquAyyh4cHVCoVioqKTJ5nxYoVKC0trT/y8/OtOm+SBicHJ8wMmgkA+OuZv9Yvv38y/x7u62swvG8XvDez4b3sY4XH8Lt//Q6l+lKEeobiLzF/QTcn00vtE9kTJ5USG54ZDLkMSMm5hS9PFZofEJkA9PsNUPPAuGBc5V3rTKyq0vh8THWFcUfoyBVNht65r69vtBaO74dAL24mas8sblbWrFkDmUxm9sjKymowpqCgADExMZg1axZefrnholqmHoQSQjT5gJRarYabm1uDgwgAngt+DkqZEtk/ZKPWIR8DfY1/G4HdXfDJv0VArfx5efz0/HS8+tWrqKypxJM+T2Jr9Fa4q01vXkhkj8J7eeDVSOPtzrdSTqG4/Kemg+UK4+0Yj97AvevA7rlArZlvE7XWgeVA8RnAuTsw41Pj+zZhzb7TKKmsRrC3K+ZF8ratvbO4WVmwYAHOnj1r9ggJCamPLygowLhx46DVarFly5YG5/L29m50BaWkpATV1dWNrrgQNcfb2RsxfWIAGK+uvPP0QMyO8MdfXhoGd6ef11L5/PLniD8Uj6raKozzH4dNUZvg5MCH8qjzWRz1OPr7uKGkshordKcgzK1Y69QFmP0PQOkIXPoKSF/XtpM5mQzk/B2ADJixFXD1bjI09XQRvsgthFwGvD8zDA4KrsJh7yyusKenJ4KDg80eGo1xMZ5bt24hMjISQ4cOxfbt2yGXN3w7rVaLvLw8FBb+fAkyNTUVarUa4eHhj5gadUZ1X2NOvZYK3656rJ8Z2mAtleRzyVh5eCUMwoBpAdOwIXID1AruHUKdk0opx8bZYVAp5Pj6XDF2ZTVzW917EDD1A+PPme8B575sm4kUnzMu7w8Ybzn1jWwytPRBNX6/x7j56NwxfTHIj1dEOwOrtaMFBQWIjIyEv78/kpKScPv2bRQVFTW4khIdHY0BAwYgLi4OOTk5+Prrr7F06VLMnTuXt3eoVfp37Y9h3sNgEAbsOLuj/nUhBLbmbsXa79YCMN4yenfku1DKTX8dkqizCPZ2w5Jo4yae73x+Bjd+rDQ/IGw28OTDZw9TXgHumFm+vyWqKh4+p1IJ9BkLjFlmNjzxy7MoLtejj6czXn+q8eajZJ+s1qykpqbi0qVLOHjwIPz8/ODj41N/1FEoFNi/fz80Gg1GjhyJZ555BrGxsUhKSrLWtKgTmDNwDgBAd1GH+1X3IYTAhuwN+DDnQwDAvLB5SBiWALmMl46JAODl0X0xrHcXVFQZsPSfJ2GobWYDw+g/AD1HGDcW3Pk8oL/f+jf/chlw+xzg4tXscyr/d+kOkr83Xv1ZPyMUGoemY8m+yITZm5S2r6ysDO7u7igtLeXVGAIA1IpaxO6NxdXSq1gSvgTXyq5Bd1EHAFgWsQwvDHyhg2dIZHtu/FiJiR9koqLKgJWTgvEfY5p5aLX8B+CTMcYVbgdOB2ZuN658a4mcHcDe14wbJ76wD+gzusnQyqoaTPhjJvLvPkDc8F54NzakyViSBks+v/m/lmR35DJ5/bMrG49vhO6iDnKZHO+MeIeNClETenZ1wttTBgAAkv51AeeLys0PcPUCZv8dkDsAp1OAIx9Z9oY/nAH2P1ytPHKl2Ualbk75dx+gx2OOeHNisGXvRZLHZoXs0pS+U9BF0wW1ohZKuRJJY5MwPXB6R0+LyKbNfsIf44O7o8pQi9d3nkBVTa35Af7DgJhE489pq4ArGS17I/1943MqNQ+AgPHAaPNbrGRfL8H2b68CANZOD4GLms+adTZsVsguaZQaxA+NR4B7ADaN34Tf9DK9WysR/Uwmk2HdjEHwcHLAmcIyfPD1heYHPfEyEPYcIGqBz/4dKL1pPl4IYP8bwJ0LgKsPMH0LIG/6o0hfY8CbulwIAfx2aA9E/mpfL+oc2KyQ3ZoeOB17YvdgRI8RHT0VIsno7qrB2umDAACb0y8j+3qJ+QEyGTBlA+AdClT+COyMA6rNLDB3/G9A7k7jcyoztgEu5leN/ujgJVwqvg9PFxVWPbxNRZ0PmxUiImpg0iAfTB/SA7UCWLLrBCqraswPcHA0Lhjn6AEUHAcONPH146I84yq1ADD+90DvkWZPe6agDJvTLwMA3nk6BI85qSxNhewEmxUiImpkzbSB8HbT4NqPlUj88lzzAzx6ATP/y3jF5PjfHu6a/Av68ofPqfxk3Gdo5OtmT1djqMVy3UnU1ArEDPTGpEGmd16mzoHNChERNeLu6ID3Z4UCAP5+9DoyLtxuflDAeOMVE8C4fsrNbOPPQgCfxwM/XgLcegDTPzH7nAoAbP3mKvJulcFNo8Q7sQMfIROyB2xWiIjIpNGB3TBH2wsAsPyzk7hXWdX8oFFvAMFTAEMVsCsOuH8byN4O5H0GyBTGqy/OXc2e4vLt+9j4lfHh3renDEB3V80j50LSxmaFiIialDCxP/p6OuOHMj1W7T3d/ACZDIjdDHQNBMpuATtmAgcSjL+LWgX0HG52eG2tQIIuF1U1tRjzeDfMDPdrgyxI6tisEBFRkxxVCmyYPRgKuQz7Thbg85MFzQ/SuAHP7gBULkDhCcCgBwInACMWNTv0H99dx/fXSuCkUuA/p4dAZumquGSX2KwQEZFZg/0fw/xI4/L7b+/Nww9lZr6aXKdbkPEKC2SAe09g+sfNPqdys6QS6w8YH+Z9MyYYfh5Ojzp1shNsVoiIqFkLxgcipIcb7lVWP1ykrQXbyg2YBizMBl77FnDqYjZUCIGVKXmoqDIgopcH4ob3aqOZkz1gs0JERM1SKeXY+MxgqJRypJ+/jf8+dqNlA7sGAGrXZsN0x28h88JtqJRyrJ8ZCrmct3/oZ2xWiIioRQK9XLF8QhAA4A9fnMW1OxVtct7i8p/w7hdnAADxTwUioJtLm5yX7AebFSIiarGXRvbBk3264EG1AUv+eRKG2hbcDmrG6r2nUfqgGiE93PAfo/u2wSzJ3rBZISKiFpPLZUiaFQYXtRLZ10uwJfPKI53vwKlCHMgrglIuw/oZoVAq+LFEjfGvgoiILOLfxQmrpho3FdyQdh5nCspadZ57lVV4++HaLfPGBmCgr3ubzZHsC5sVIiKy2KxwPzzV3wvVBoE3dp2AvsZg8Tne/eIs7tzXI6CbMxZG9bPCLMlesFkhIiKLyWQyrJsxCF2dVThXVI6NaRctGp9+vhi64zchkwHvzQyDWqmw0kzJHrBZISKiVvF0UWPt9EEAgE8yL+P7a3dbNO6+vgZvpeQBAF4c0RvhvTysNkeyD2xWiIio1WJCvDFjqB+EAJbsOon7+ppmx7z3v+dw694D+Hk4Yml0UDvMkqSOzQoRET2S1dMGwNddgxt3K7F2/1mzsceu3sXfjlwHAKz7bSic1cr2mCJJHJsVIiJ6JG4aByTNCgMA/M+xGzh0rthk3E/VBiTocgEAsyP8MSrQs93mSNLGZoWIiB7ZiH6eeGlkHwDAcl0uSiqqGsX88auLuHKnAt1d1Vg5uX97T5EkjM0KERG1ieUxQQjo5ozb5Xr8fk9eg80OT90sxdZvjAvI/SE2BO6ODh01TZIgNitERNQmNA4KbJw9GAq5DPtPFWLfyQIAQLWhFst1uTDUCkwJ9UH0QO8OnilJDZsVIiJqM6F+j2HheOMCb2/vyUNh6QN8nH4ZZwvL4OHkgDXTBnbwDEmK+Bg2ERG1qfnj+uHguWLk3izFvH8cx9mHy/GvnjoQni7qDp4dSRGvrBARUZtyUMix4ZnBUCvlOJl/D1WGWowL6oanB/t29NRIotisEBFRm+vX3QUJE4MBAC5qJdZOHwSZTNbBsyKp4m0gIiKyijna3nBRK9Gvuwt8H3Ps6OmQhLFZISIiq5DLZZgV4d/R0yA7wNtAREREZNPYrBAREZFNY7NCRERENo3NChEREdk0qzYr06ZNQ8+ePaHRaODj44O4uDgUFBQ0iLlx4wamTp0KZ2dneHp6YtGiRaiqarwBFhEREXVOVm1Wxo0bh127duH8+fPQ6XS4fPkyZs6cWf97g8GAyZMno6KiAocPH0ZycjJ0Oh2WLFlizWkRERGRhMjEL7fFtLJ9+/YhNjYWer0eDg4OOHDgAKZMmYL8/Hz4+hpXNkxOTsaLL76I4uJiuLm5NXvOsrIyuLu7o7S0tEXxRERE1PEs+fxut2dW7t69ix07dmDEiBFwcDBuDX7kyBGEhITUNyoAMGHCBOj1emRnZ5s8j16vR1lZWYODiIiI7JfVm5U333wTzs7O6Nq1K27cuIG9e/fW/66oqAheXl4N4j08PKBSqVBUVGTyfImJiXB3d68//P254BAREZE9s7hZWbNmDWQymdkjKyurPn7ZsmXIyclBamoqFAoFXnjhBfzyzpOpvSKEEE3uIbFixQqUlpbWH/n5+ZamQERERBJi8XL7CxYswLPPPms2pnfv3vU/e3p6wtPTE48//jj69+8Pf39/HD16FFqtFt7e3vjuu+8ajC0pKUF1dXWjKy511Go11GpuMU5ERNRZWNys1DUfrVF3RUWv1wMAtFot1q5di8LCQvj4+AAAUlNToVarER4e3qr3ICIiIvtitY0Mjx07hmPHjmHUqFHw8PDAlStXsGrVKgQEBECr1QIAoqOjMWDAAMTFxeH999/H3bt3sXTpUsydO5ff7CEiIiIAVmxWHB0dsXv3bqxevRoVFRXw8fFBTEwMkpOT62/jKBQK7N+/H6+99hpGjhwJR0dHPPfcc0hKSmrx+9RdreG3goiIiKSj7nO7JSuotOs6K9Zw8+ZNfiOIiIhIovLz8+Hn52c2RvLNSm1tLQoKCuDq6trkN4haq6ysDP7+/sjPz7fL21LMT/rsPUfmJ332nqO95wdYL0chBMrLy+Hr6wu53PyXk612G6i9yOXyZjuyR+Xm5ma3f4QA87MH9p4j85M+e8/R3vMDrJOju7t7i+K46zIRERHZNDYrREREZNPYrJihVquxevVqu12EjvlJn73nyPykz95ztPf8ANvIUfIP2BIREZF945UVIiIismlsVoiIiMimsVkhIiIim8ZmhYiIiGwamxUiIiKyaZ26Wfnzn/+MPn36QKPRIDw8HN98843Z+IyMDISHh0Oj0aBv3774+OOP22mmrWdJjunp6ZDJZI2Oc+fOteOMWy4zMxNTp06Fr68vZDIZ9uzZ0+wYKdXQ0vykVr/ExEQ88cQTcHV1Rffu3REbG4vz5883O04qNWxNflKr4ebNmxEaGlq/sqlWq8WBAwfMjpFK/QDL85Na/X4tMTERMpkM8fHxZuM6ooadtlnZuXMn4uPj8dZbbyEnJwejR4/GxIkTcePGDZPxV69exaRJkzB69Gjk5ORg5cqVWLRoEXQ6XTvPvOUszbHO+fPnUVhYWH8EBga204wtU1FRgbCwMHz00UctipdaDS3Nr45U6peRkYH58+fj6NGjSEtLQ01NDaKjo1FRUdHkGCnVsDX51ZFKDf38/LBu3TpkZWUhKysL48ePx9NPP43Tp0+bjJdS/QDL86sjlfr90vfff48tW7YgNDTUbFyH1VB0UsOGDRPz5s1r8FpwcLBISEgwGb98+XIRHBzc4LVXXnlFDB8+3GpzfFSW5njo0CEBQJSUlLTD7NoWAJGSkmI2Roo1rNOS/KRcPyGEKC4uFgBERkZGkzFSrmFL8pN6DYUQwsPDQ3z66acmfyfl+tUxl59U61deXi4CAwNFWlqaGDt2rFi8eHGTsR1Vw055ZaWqqgrZ2dmIjo5u8Hp0dDS+/fZbk2OOHDnSKH7ChAnIyspCdXW11ebaWq3Jsc6QIUPg4+ODqKgoHDp0yJrTbFdSq2FrSbV+paWlAIAuXbo0GSPlGrYkvzpSrKHBYEBycjIqKiqg1WpNxki5fi3Jr47U6jd//nxMnjwZTz31VLOxHVXDTtms3LlzBwaDAV5eXg1e9/LyQlFRkckxRUVFJuNrampw584dq821tVqTo4+PD7Zs2QKdTofdu3cjKCgIUVFRyMzMbI8pW53UamgpKddPCIE33ngDo0aNQkhISJNxUq1hS/OTYg1PnToFFxcXqNVqzJs3DykpKRgwYIDJWCnWz5L8pFi/5ORkHD9+HImJiS2K76gaKq12ZgmQyWQN/i2EaPRac/GmXrclluQYFBSEoKCg+n9rtVrk5+cjKSkJY8aMseo824sUa9hSUq7fggULkJubi8OHDzcbK8UatjQ/KdYwKCgIJ06cwL1796DT6TBnzhxkZGQ0+YEutfpZkp/U6pefn4/FixcjNTUVGo2mxeM6ooad8sqKp6cnFApFoysMxcXFjTrGOt7e3ibjlUolunbtarW5tlZrcjRl+PDhuHjxYltPr0NIrYZtQQr1W7hwIfbt24dDhw7Bz8/PbKwUa2hJfqbYeg1VKhX69euHiIgIJCYmIiwsDB988IHJWCnWz5L8TLHl+mVnZ6O4uBjh4eFQKpVQKpXIyMjAhx9+CKVSCYPB0GhMR9WwUzYrKpUK4eHhSEtLa/B6WloaRowYYXKMVqttFJ+amoqIiAg4ODhYba6t1ZocTcnJyYGPj09bT69DSK2GbcGW6yeEwIIFC7B7924cPHgQffr0aXaMlGrYmvxMseUamiKEgF6vN/k7KdWvKebyM8WW6xcVFYVTp07hxIkT9UdERASef/55nDhxAgqFotGYDquhVR/ftWHJycnCwcFBbNu2TZw5c0bEx8cLZ2dnce3aNSGEEAkJCSIuLq4+/sqVK8LJyUm8/vrr4syZM2Lbtm3CwcFBfPbZZx2VQrMszXHjxo0iJSVFXLhwQeTl5YmEhAQBQOh0uo5Kwazy8nKRk5MjcnJyBACxYcMGkZOTI65fvy6EkH4NLc1PavV79dVXhbu7u0hPTxeFhYX1R2VlZX2MlGvYmvykVsMVK1aIzMxMcfXqVZGbmytWrlwp5HK5SE1NFUJIu35CWJ6f1Opnyq+/DWQrNey0zYoQQmzatEn06tVLqFQqMXTo0AZfKZwzZ44YO3Zsg/j09HQxZMgQoVKpRO/evcXmzZvbecaWsyTH9evXi4CAAKHRaISHh4cYNWqU2L9/fwfMumXqvib462POnDlCCOnX0NL8pFY/U7kBENu3b6+PkXINW5Of1Gr40ksv1f/3pVu3biIqKqr+g1wIaddPCMvzk1r9TPl1s2IrNZQJ8fDJGCIiIiIb1CmfWSEiIiLpYLNCRERENo3NChEREdk0NitERERk09isEBERkU1js0JEREQ2jc0KERER2TQ2K0RERGTT2KwQERGRTWOzQkRERDaNzQoRERHZtP8HGM2A3vFW5UgAAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "time = np.arange(ctr*dt)\n", "\n", "plt.title(\"\")\n", "plt.plot(time,md0mdt, label=\"total\")\n", "plt.plot(time,md0mdt_cat, label=\"cation\")\n", "plt.plot(time,md0mdt_an, label=\"anion\")\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 }