{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Calculating the RDF atom-to-atom\n", "\n", "We calculate the site-specific radial distribution functions of solvent around certain atoms.\n", "\n", "**Last executed:** Feb 06, 2020 with MDAnalysis 0.20.2-dev0\n", "\n", "**Last updated:** February 2020\n", "\n", "**Minimum version of MDAnalysis:** 0.19.0\n", "\n", "**Packages required:**\n", " \n", "* MDAnalysis (Michaud-Agrawal *et al.*, 2011, Gowers *et al.*, 2016)\n", "* MDAnalysisTests\n", "\n", "**Optional packages for visualisation:**\n", "\n", "* [matplotlib](https://matplotlib.org)\n" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import MDAnalysis as mda\n", "from MDAnalysis.tests.datafiles import TPR, XTC\n", "from MDAnalysis.analysis import rdf\n", "import matplotlib.pyplot as plt\n", "import numpy as np\n", "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Loading files" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The test files we will be working with here feature adenylate kinase (AdK), a phosophotransferase enzyme. (Beckstein *et al.*, 2009)" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "u = mda.Universe(TPR, XTC)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Calculating the site-specific radial distribution function" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A radial distribution function $g_{ab}(r)$ describes the time-averaged density of particles in $b$ from the reference group $a$ at distance $r$. It is normalised so that it becomes 1 for large separations in a homogenous system. See [the tutorial on averaged RDFs](average_rdf.ipynb) for more information. The `InterRDF_s` class ([API docs](https://docs.mdanalysis.org/stable/documentation_pages/analysis/rdf.html#MDAnalysis.analysis.rdf.InterRDF_s)) allows you to compute RDFs on an atom-to-atom basis, rather than simply giving the averaged RDF as in `InterRDF`.\n", "\n", "Below, I calculate the RDF between selected alpha-carbons and the water atoms within 15 angstroms of CA60, *in the first frame of the trajectory*. The water group does not update over the trajectory as the water moves towards and away from the alpha-carbon. \n", "\n", "The RDF is limited to a spherical shell around each atom by `range`. Note that the range is defined around *each atom*, rather than the center-of-mass of the entire group.\n", "\n", "If `density=True`, the final RDF is over the average density of the selected atoms in the trajectory box, making it comparable to the output of `rdf.InterRDF`. If `density=False`, the density is not taken into account. This can make it difficult to compare RDFs between AtomGroups that contain different numbers of atoms." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "ca60 = u.select_atoms('resid 60 and name CA')\n", "ca61 = u.select_atoms('resid 61 and name CA')\n", "ca62 = u.select_atoms('resid 62 and name CA')\n", "water = u.select_atoms('resname SOL and sphzone 15 group sel_a', sel_a=ca60)\n", "\n", "ags = [[ca60+ca61, water], [ca62, water]]\n", "\n", "ss_rdf = rdf.InterRDF_s(u, ags,\n", " nbins=75, # default\n", " range=(0.0, 15.0), # distance\n", " density=True,\n", " )\n", "ss_rdf.run();" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Like `rdf.InterRDF`, the distance bins are available at `ss_rdf.bins`." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([ 0.1, 0.3, 0.5, 0.7, 0.9, 1.1, 1.3, 1.5, 1.7, 1.9, 2.1,\n", " 2.3, 2.5, 2.7, 2.9, 3.1, 3.3, 3.5, 3.7, 3.9, 4.1, 4.3,\n", " 4.5, 4.7, 4.9, 5.1, 5.3, 5.5, 5.7, 5.9, 6.1, 6.3, 6.5,\n", " 6.7, 6.9, 7.1, 7.3, 7.5, 7.7, 7.9, 8.1, 8.3, 8.5, 8.7,\n", " 8.9, 9.1, 9.3, 9.5, 9.7, 9.9, 10.1, 10.3, 10.5, 10.7, 10.9,\n", " 11.1, 11.3, 11.5, 11.7, 11.9, 12.1, 12.3, 12.5, 12.7, 12.9, 13.1,\n", " 13.3, 13.5, 13.7, 13.9, 14.1, 14.3, 14.5, 14.7, 14.9])" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ss_rdf.bins" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "`ss_rdf.rdf` contains the atom-pairwise RDF for each of your pairs of AtomGroups. It is a list with the same length as your list of pairs `ags`. A result array has the shape `(len(ag1), len(ag2), nbins)` for the AtomGroup pair `(ag1, ag2)`. " ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "scrolled": true }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "There are 1041 water atoms\n", "The first result array has shape: (2, 1041, 75)\n", "The second result array has shape: (1, 1041, 75)\n" ] } ], "source": [ "print('There are {} water atoms'.format(len(water)))\n", "print('The first result array has shape: {}'.format(ss_rdf.rdf[0].shape))\n", "print('The second result array has shape: {}'.format(ss_rdf.rdf[1].shape))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Index the results array to get the RDF for a particular pair of atoms. `ss_rdf.rdf[i][j][k]` will return the RDF between atoms $j$ and $k$ in the $i$th pair of atom groups. For example, below we get the RDF between the alpha-carbon in residue 61 (i.e. the second atom of the first atom group) and the 571st atom of water." ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "scrolled": true }, "outputs": [ { "data": { "text/plain": [ "array([0. , 0. , 0. , 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0. ,\n", " 0.41218412, 0. , 0. , 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0.19906709,\n", " 0.18622866, 0. , 0.1640152 , 0. , 0. ,\n", " 0. , 0.13003857, 0. , 0. , 0. ,\n", " 0. , 0. , 0.09591515, 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0. ,\n", " 0. , 0.05427213, 0. , 0. , 0. ,\n", " 0. , 0. , 0.04435226, 0.04296637, 0. ,\n", " 0. , 0. , 0. , 0. , 0. ,\n", " 0. , 0. , 0. , 0. , 0. ])" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ca61_h2o_571 = ss_rdf.rdf[0][1][570]\n", "ca61_h2o_571" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYIAAAEWCAYAAABrDZDcAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8li6FKAAAgAElEQVR4nO29eZxkdXnv//50VVevswA9IMyAg4ALiwuOuKC4KxgF43bFDaO5/O7rSsSoUdRcoyQxaowxUcyN2wWNxihqJIpxR9xhAIUZZZlBkBnAmWGZmV6ma3t+f5zvqT5dXcuppquruup5v1796qpzTp16qpfzOc/yfR6ZGY7jOE7/MtBpAxzHcZzO4kLgOI7T57gQOI7j9DkuBI7jOH2OC4HjOE6f40LgOI7T57gQOC0j6QpJf9ppO5zW8d+dUwsXgh5A0m2SZiRNSrpb0sWSxhP7L5aUl7Q/fG2R9HeS1iSOea2kUjhH/PWxNti6Ii5Ekl4haXP4Odwl6VuSnlx1zGslmaT/UeP1o5I+LmmPpL2Srkzse7qkH4btty3Dx0mFpPeEz3N+1fbzw/b3hOc3JT+zpFOrfw5h235J2fB8naQvhM98n6TPJ479oKQ7JO2TdLukd1a9f0bS30i6M5zzOklra9j//WBHdsl+KH2CC0Hv8AIzGwceDTwGeEfV/g+a2SpgHfAnwBOAn0oaSxzzczMbT3ydtyyWdxmS3gx8BHgfcBhwFPBx4KyqQ88B7gVeU+M0nwAOBh4Rvv95Yt8U8BngL5bU8KXhZhZ+nnPC9pgrgdMSz08Dbqyx7edmVgzPvwrcTfSzPBT4UOLYTwMPN7PVwJOAV0p6UWL/e8P2JwKrgVcDB5IGSnolMJjuIzrVuBD0GGZ2N/BtIkGotf+AmV0NnAkcQiQKi+EYSVeFu7ivSzo43iHpCZJ+Jul+Sb+W9LSw/W+BpwAfiz0OSe+V9NGwf1DSlKS/D89HJB2Iz13vvGHfGkmfDnfvO8MdZCbse62kn0j6ULgb/Z2kM2p9qOAlXQi8wcy+amZTZlYws/8ys79IHPdg4KnAucBzJT0ose/h4ed7rpntNrOSmV2T+B1cZWafA25N84OW9OXg6e2VdKWkExL7LpZ0kaRvhrvlX0o6JrH/2ZJuDK/9GKAmb3c1MBq/R/g+HLbHVAvBU4AP1Nh2ZTjHc4Ajgb8ws73h53ld4udxk5lNJV5bBo4Nrz0IeBPwP83sdovYYmYVIQi/s78C3tbkszl1cCHoMSRtAM4AtjU6zsz2A98l+oddDK8BXgccDhSBfw7vvx74JvA3RHfCbwW+Immdmb0L+DFwXsLj+BHwtHDOxxHdNcYXlCcCN5nZvY3OG469ONhxLJFH9BwgGYJ6PHATMAF8EPi0pFoXxScSXfi+luLzbzazrwC/BV6Z2HcKcDvw3hAaukHSi5ucrxHfAo4jupO+Fvh81f6XE901H0T0e/9bAEkTRHfif0n0ubcDp6Z4v88x5xWcE54nuRI4QdLBkgaATcB/AGsT204Nx0Hkfd4EXCLpHklXS3pq8oSSLpA0CewAxoAvhF0nEf1eXxLE8GZJb6iy533AvxD97TiLwIWgd/hPSfuBO4BdRHdIzbiT6KIa84Rwtx1/PaHBaz8X7symgP8DvCzcgb8KuNzMLjezspl9F9gMPK/OeX4OHCfpECIB+DSwXlGO46lEQkGj80o6LJz/TeEOfhfwj0QXyJjbzeyTZlYCLiESsMNq2HMIsCcR0qjHa5i7WH2B+eGUDcCJwF7gCOA8oovgI5qcsyZm9hkz229ms8B7gEcpkd8Bvha8jCKRSMTe4POArWZ2qZkViMJdaS6W/wacLWmQ6Gf4b1X23A78nugm4lHALWY2A/w0sS0H/DK8ZAORMP8QeBDwD8DXg1DF53w/sAo4mUh49iZeuwZ4KHA08BLgPZKeDSBpE5HofDTF53Lq4ELQO7ww5ACeBjyc6A6wGeuJYtwxvzCztYmvXzR47R2Jx7cTxWcngAcDL00KCvBkogvvAsIFZDPRRf80ogv/z4j+uZNC0Oi8Dw7vf1di378S3UHHVC6AZjYdHo6zkHuAiUYJR0mnEl2Uvhg2fQE4SVJ8AZ4BCsDfmFnezH5EdBF8Tr1zNnivjKT3S9ouaR9wW9iV/P0mL+7TzH2uI0j8nizqMJn8vdXEzH5P5Fm8j+giX+s1cXjoNCIvD+AniW1XBeGC6Odxm5l9OoSFvhjsmOedhLDPdeH49yZeC3Chmc2Y2fVEP/fnBc/j48D5KYTbaYALQY8RLjoXMz8Zt4Bwx/0s5v6JW+XIxOOjiC58e4j+wT9XJShj4Y4PoFa72x8BzyAK6Vwdnj+XKMQShxcanfcOYBaYSOxbbWYnLHyrpvw8nOuFDY45hyjW/itJdzN353tO+H59jdcsts3vK4iS1M8iujPeGLY3i/UD3EXi9xRCYUfWP3wenwXeEr7XIhaCpzD3N/TjxLYrE8dez8LP3+jnkQXiPEf8s0weHz9eTQhLhd9DnMfYIWmxIc++xIWgN/kI8GxJj6reIWlI0mOB/wTuA/7fIt/jVZKOlzRKlFy9NIRd/g14gaTnhrvZYUlPC7kLgD8AD6k614+IQiu/MbM8cAVRfP93ZrY7HFP3vGZ2F/Ad4B8krZY0IOmY6jh0GsxsL/Bu4CJJL1RUBjoo6QxFZY7DwMuIksSPTnz9GfCK4ElcSRQ6eYekbPAgnk6UxCfYN0zkxSh8llwdk1YRCdM9wCjRXXpavkkUy39RsOuNRKGZNPwHkQfzpTr7ryQS7tOIQkIANxB5Sk9nvhB8DThI0jnhd/cSopDPT8PP4v+TdJAiTgHeAHwfwMy2EwnMu8Lf7iOIwlXfYC70Fv8O4vDjY5kTZycFLgQ9SLh4fpboghbztpBDuCfsuwZ4UlW1Rit8jsjzuJsoufrG8N53EN3BvhPYTXS3/hfM/a39E1Hi7z5J/xy2/QwYYe7i8Rui8sDKxSTFeV9DFJf+DZHAXUqdcFQzzOwfgDcTJVnj9zqPSDxfSBSu+KyZ3R1/EZWDZoHTQzz+LKIL017gk8BrzOzG8BanhXNcTuRNzRAJWS0+SxR62xk+W6NwXfXn2AO8FHg/0e/9OOYu2s1eO2Nm3wuhu1r7byb62dxtZveHbWXgKqI79Z8ljr2XqIrqrUQ/jwuAs4J9AH9MlMjeTyT4H2V+zP9sovDfPUTi9n/M7PshlJT8HcQ3DX8INxROSmQ+mMZxHKevcY/AcRynz3EhcBzH6XNcCBzHcfocFwLHcZw+Z8V16ZuYmLCNGzd22gzHcZwVxTXXXLPHzNbV2rfihGDjxo1s3ry502Y4juOsKCTdXm+fh4Ycx3H6HBcCx3GcPseFwHEcp89xIXAcx+lzXAgcx3H6HBcCx3GcPseFwHEcp89xIXBawsz40uY7mC2WOm2K4zhLhAuB0xJb79zH2y69np/csqf5wY7jrAhcCJyWiD2BA4Vyhy1xHGepcCFwWqJQigYZFcsuBI7TK7RVCCSdLukmSdskXdDguBdLMkmb2mmP88ApBiGIBcFxnJVP24RAUga4CDgDOB44W9LxNY5bBZyPD5teERRK5XnfHcdZ+bTTIzgF2GZmt4ZB0l8kGuhdzV8DHyAaVu50ObEAFF0IHKdnaKcQrAfuSDzfEbZVkHQycKSZfbPRiSSdK2mzpM27d+9eekud1BTLHhpynF6jY8liSQPAh4G3NDvWzD5hZpvMbNO6dTXnKjjLhIeGHKf3aKcQ7ASOTDzfELbFrAJOBK6QdBvwBOAyTxh3N3NVQ+4ROE6v0E4huBo4TtLRknLAy4HL4p1mttfMJsxso5ltBH4BnGlmPn6siym6R+A4PUfbhMDMisB5wLeB3wJfMrOtki6UdGa73tdpLx4acpzeo60zi83scuDyqm3vrnPs09ppi7M0VEJDnix2nJ7BVxY7LRGvKPaqIcfpHVwInJYoVFYWe2jIcXoFFwKnJSoLyrzXkOP0DC4ETkt4ryHH6T1cCJyW8Kohx+k9XAiclvCqIcfpPVwInJaYqxpyj8BxegUXAqclvMWE4/QeLgROS3iOwHF6DxcCpyW815Dj9B4uBE5LFMqeLHacXsOFwGmJQtE9AsfpNVwInJbwCWWO03u4EDgt4S0mHKf3cCFwWmKuasg9AsfpFVwInJYoevdRx+k5XAiclvCqIcfpPVwInJbwqiHH6T1cCJyW8F5DjtN7uBA4LVH0XkOO03O4EDgtkfcWE47Tc7gQOC2RnFBm5l6B4/QCLgROSyQXkpU8POQ4PYELgdMS+eKcEPiiMsfpDVwInJYolg0pelzwNhOO0xO4EDgtUSwZo4OZymPHcVY+LgROasyMfKnMSC4LeOWQ4/QKLgROauLk8Egu+rNxIXCc3sCFwElNvIhsdDDyCDw05Di9gQuBk5p4MdlILsoRuEfgOL2BC4GTmtgDGBmMhcA9AsfpBVwInNQUgwcwGjwCn1LmOL2BC4GTGg8NOU5v4kLgpMZDQ47Tm7gQOKmJQ0GV0JALgeP0BC4ETmryxXgdgS8oc5xewoXASU3sEcyFhlwIHKcXcCFwUhPnBOKVxT6lzHF6AxcCJzWFStWQh4Ycp5doqxBIOl3STZK2Sbqgxv7/JekGSb+S9BNJx7fTHueBESeHR71qyHF6irYJgaQMcBFwBnA8cHaNC/0XzOwkM3s08EHgw+2yx3ngxPMHRipVQ+4ROE4v0E6P4BRgm5ndamZ54IvAWckDzGxf4ukY4LeYXUyh6AvKHKcXybbx3OuBOxLPdwCPrz5I0huANwM54Bm1TiTpXOBcgKOOOmrJDXXSMdd91ENDjtNLdDxZbGYXmdkxwNuBv6xzzCfMbJOZbVq3bt3yGuhUKFS1mPBeQ47TGzT1CCStA/4nsDF5vJm9rslLdwJHJp5vCNvq8UXgX5rZ43SO2AMYzblH4Di9RJrQ0NeBHwPfA0otnPtq4DhJRxMJwMuBVyQPkHScmd0Snv4RcAtO1xInh4d9QZnj9BRphGDUzN7e6onNrCjpPODbQAb4jJltlXQhsNnMLgPOk/QsoADcB5zT6vs4y0ch5Ahy2QEG5L2GHKdXSCME35D0PDO7vNWTh9dcXrXt3YnH57d6TqdzxFVDucwAg5kB9wgcp0dIkyw+n0gMDkjaH772NX2V03PEyeFsRQjcI3CcXqCpR2Bmq5bDEKf7iS/82QGRzcirhhynR0i1jkDSmcBp4ekVZvaN9pnkdCtxKGjQQ0OO01M0DQ1Jej9ReOg34et8SX/XbsOc7qNYMgYEmQExOCAPDTlOj5DGI3ge8GgzKwNIugS4DnhHOw1zuo9CuUw2E907ZDMD3mvIcXqEtCuL1yYer2mHIU73UygauSAEgxn3CBynV0jjEfwdcJ2kHwIiyhUsaCnt9D7FcplsRgCeI3CcHiJN1dC/S7oCeFzY9HYzu7utVjldSaFkZAfi0JB8Qpnj9Ah1Q0OSHh6+nwwcTtQ9dAdwRNjm9BnFUpmcewSO03M08gjeTNT6+R9q7DPqtIx2epdCaS5ZPDjgQuA4vUJdITCzc8PDM8zsQHKfpOG2WuV0JYWyVXIE2YzIF10IHKcXSFM19LOU25weJwoNxVVDA5UmdI7jrGzqegSSHkQ0ZWxE0mOIKoYAVgOjy2Cb02UUSpaoGlKlCZ3jOCubRjmC5wKvJRookxwqvx94ZxttcrqUQqk8VzU0MOC9hhynR2iUI7gEuETSi83sK8tok9OlFEuJBWXZAZ9H4Dg9QpoFZSdKOqF6o5ld2AZ7nC6mUCqTy8ZVQyLvVUOO0xOkSRZPAlPhqwScQTS/2OkzoqqhxIIy9wgcpydIs7J43joCSR8iGj/p9BnVC8o8R+A4vUHapnNJRokSyE6fkUwWD2YGfB2B4/QITT0CSTcQrSSGaAj9OsDzA31IMVE+mh3wXkOO0yukSRY/P/G4CPzBzIptssfpYgrlslcNOU4P0jQ0ZGa3A4cAZwEvAk5qt1FOd1IoJhaUhaohMxcDx1nppBlV+W7gEiIxmAAulvSX7TbM6T6KVRPKAEoeHnKcFU+a0NArgUfFjefCDONfAX/TTsOc7qOQXFAWvhfLRjbTSascx3mgpKkauhNIdhsdAna2xxynm4mqhuZ6DQG+qMxxeoBGTec+SlQttBfYKum74fmzgauWxzynm4iqhuJeQ6pscxxnZdMoNLQ5fL8G+Fpi+xVts8bpWswsVA0FjyC0mii6R+A4K55mTeccB4iSwmbMm1AGHhpynF6gUWjoS2b2sqoFZRXM7JFttczpKuLFY8kJZeChIcfpBRqFhs4P35/f4BinT4jnEy+sGnKPwHFWOo1CQ3dJygAXm9nTl9EmpwsphDv/BVVDRfcIHGel07B81MxKQFnSmmWyx+lS4qTwXNWQewSO0yukWVA2CdwQyken4o1m9sa2WeV0HfGg+mSvIZjzFBzHWbmkEYKvhq8k/t/fZ8SD6pO9hmAud+A4zsoljRCsNbN/Sm6QdH69g53eJA4BVfca8qohx1n5pGkxcU6Nba9dYjucLicOAc1NKAsegecIHGfF02gdwdnAK4CjJV2W2LUauLfdhjndRRwCSk4og7mQkeM4K5dGoaGfAXcRtZ5Ozi3eD1zfTqOc7qNSPlq9oMzbUDvOiqfROoLbgdslPQuYMbOypIcCDwduWC4Dne6gWGdBmSeLHWflkyZHcCUwLGk98B3g1cDFaU4u6XRJN0naJumCGvvfLOk3kq6X9H1JD27FeGf5mPMI5vca8vJRx1n5pBECmdk00ZjKj5vZS4ETmr4oWpV8EXAGcDxwtqTjqw67DtgU+hZdCnywFeOd5aNQnl8+OtdryD0Cx1nppBICSU8kmlT2zbAtzUyqU4BtZnarmeWBLxLNPa5gZj8MIgPwC2BDOrOd5SYuE409gWylaqg7PIItO/dS7hJbHGelkUYI3gS8A/iamW2V9BDghyletx64I/F8R9hWj9cD36q1Q9K5kjZL2rx79+4Ub+0sNXEuYDAbCUCui6qGbt09yfM/+hM+8r2bO22K46xImgqBmf3IzM40sw+E57cudXsJSa8CNgF/X8eGT5jZJjPbtG7duqV8aycl1eWj2S7qPnrvVB6Ai67YzpadeztsjeOsPOoKgaSPhO//Jemy6q8U594JHJl4voEas45DVdK7gDPNbLY1853lohIainMElRYTnQ/HTOdLlcdv/fKvyXeBl+I4K4lG6wg+F75/aJHnvho4TtLRRALwcqIFahUkPQb4V+B0M9u1yPdxloFKaKgLy0djITj/mcfx4e/ezEd/cAtvec7DOmyV46wcGq0juCZ8/9FiTmxmRUnnAd8mSi5/JuQYLgQ2m9llRKGgceDLkgB+b2ZnLub9nPZSqJpQlhkQA+qOXkMzhSIAz3/k4dx2zxQfv2I7zzn+QZy0wbunO04aGrWYqDmiMibNqEozuxy4vGrbuxOPn5XOTKfTxGWicdUQRHmCbug1FHsEo7ksf/X8E/jptj289cu/5rI/O5WhbJoCN8fpbxoli58PvAD47/D1yvD1Laou7k7vM1c1NPcnk8sMUOiCCWUzQQhGchnWjA7yvj8+iZv+sJ+v/+rODlvmOCuDZi0mkPRsM3tMYtfbJV0LLFgp7PQu1aMqIQoTdUPV0JxHEN39P/m4CQB27/faA8dJQ9oFZacmnjwp5eucHmKuaigRGhoY6JqqocGMKrblMgNkB8R0vthhyxxnZZBmMM3rgc8k5hbfD7yufSY53UihVGZAUZI4JpdRV1QNzeSLjAzO5QIkMZrLMDVbavAqx3FimgpBqB56VCwEZuYrdvqQQrlcWUQWk80MdEWvoel8idHc/D/l0VzWPQLHSUkajwBwAeh3iiWrzCmOyWbUFb2GpgulSn4gZnQow1TePQLHSYPH+p1UFErleRVDEFcNdd4jmMmXGKkSgrFctlJN5DhOY1wInFQUSlbpMxQTVQ11gUeQLy70CHIZpmY9NOQ4aWi0oOxFjV5oZl9denOcbqVYKlf6DMVEVUPd4RGsGc3N2zY2lPXyUcdJSaMcwQsa7DPAhaCPKJZtXukoRKGhbmgxMZ0vcfia+R7BSC7DlCeLHScVjRaU/clyGuJ0N/lSudJnKCabUVd0+oyqhqpzBBmmG5SPbr7tXtaODnLsoavabZ7jdD2pqoYk/RHReMrheJuZXdguo5zuo1gqz+szBFH5aDdU5swUFiaLm5WPvu0r1/OIw1dz0StObrd5jtP1NE0WS/q/wP8A/gwQ8FLAh8z3GcWSVaaTxeQy6pJ1BEXGhubf04wNZZjOlzCrHbraN1Ng30xhOcxznK4nTdXQk8zsNcB9ZvZe4InAQ9trltNt5EvlhVVDXZAsLpeNA4XyvJXFEHkExbKRr2Pf5GyRSa8qchwgnRDMhO/Tko4ACsDh7TPJ6UaKJVtYNZRRx5PFM4X5Dedi4ue18gTFUpkDhTKTB1wIHAfSCcE3JK0lGiJzLXAb8O/tNMrpPorlcs2qoU7PI6juPBozFlpO1KocinsQ+ToDx4lI02vor8PDr0j6BjDs7Sb6j3zJGMktXFDW6XkEcUJ4pLrX0FAkDLVWF++fjXIDHhpynIhGC8qeYWY/qLWwTJIvKOszoqqh6tDQQMfnETT3CBYKQewRTM4WMTPCmFTH6VsaeQRPBX5A7YVlvqCsz4hyBDVCQx3OEUwnppMlGankCBbe9ceeQNmIEs05H2fp9DeNFpT9VfjuC8scCrUWlA10fh5BHPoZHUzvESRDQvtnCy4ETt/TKDT05kYvNLMPL705TrdSqJEsznZBi4k4R7BgHkHIEdRaVJZMEk/NlsAXFzt9TqPQUPzv8TDgccBl4fkLgKvaaZTTfdQqH81lRKFc7micPS4frdWGGuZCR0km5wmBJ4wdp1Fo6L0Akq4ETjaz/eH5e4BvLot1TtcQhYYWegRmUCrbgrDRclEvWRx7BLUu9Mn1A/t9LYHjpFpHcBiQTzzPh21OH1GoM6EM6OhMgrpCMJiZtz/JlHsEjjOPNE3nPgtcJelr4fkLgUvaZ5LTjUTzCBZWDUHkLQwPdibhOlNZRzD//bOZAXLZgZoLypKhIV9L4DjpFpT9raT/Bp4cNv2JmV3XXrOcbqNQsoWhoeAhdLKEdDpfIjOgiiglGctlai4om5wtIoGZC4HjQMo21GZ2jaQ7CG2oJR1lZr9vq2VOVxFVDS1cUAZ0tAPpdL7E6GCmZrJ6NJetLB5LMjVb5JCxIfZMzi46NGRm3DuV55DxoUW93nG6iTRtqM+UdAvwO+BH4fu32m2Y0z2UyoYZ9UNDHcwR1BpcHzOay9QsH52cLTIxnkNavEfw7a1/4Inv/wH3TeWbH+w4XU6aZPFfA08Abjazo4FnAb9oq1VOVxEvGqs1oQyg0MEpZdOFhdPJYkaHsnUXlK0eHmQ8l120ENx+zxT5Ypk/7D+wqNc7TjeRRggKZnYPMCBpwMx+CGxqs11OFxELQa0JZUBH+w3N5IsLGs7FROMqa3sEY0MZxoayiw4N7Q1DbfZO+3AbZ+WTJkdwv6Rx4Erg85J2AVPtNcvpJuLVw7UWlEHnk8V1PYJclvunZxZsn5otsfGQLGNDmUV7BBUh8ClnTg+QxiM4C5gG/hz4b2A7tRvROT3KXGho4YSy5P5O0EgIonGVtT2CVcNZxocHmWww4L4R+8JCtH2+IM3pAZoKgZlNmVnZzIpmdgnwMeD09pvmdAtxMrjWhDLorEcwky8tGFMZM5rL1M4RHCgylssyPpRh8sDi7ujdI3B6ibpCIGm1pHdI+pik5yjiPOBW4GXLZ6LTaeLy0HpVQx0tHy0UG4aGqnMEpbIxUygxPpxlfKh2eWkaXAicXqJRjuBzwH3Az4E/Bd4JCHihmf1qGWxzuoS6oaHKyuJOl482SBYXSvOa4sUrjceHsowNLb5qaF8QgH0uBE4P0EgIHmJmJwFI+hRwF3CUmXm9XJ8RX+jr9Rrq5NzihsnioSxWNXwmbjg3NhR5BJ4sdpzGOYLKX7iZlYAdLgL9yVzVUL3QUGc8ArMozFM3WRy2J/sNxeWi40NxaCgaV9nq+7oQOL1EI4/gUZL2hccCRsJzAWZmq9tundMV5JstKOtQjuBAoYzZwoZzMXHIaHq2BOPRtv2z80NDxbIxW2ytad5UvkQpJNBdCJxeoNE8Ap/f5wD1k8WdLh+tTCercxFv6BEMZ1k1HP35T84WWxKC5MXfhcDpBdKsI1g0kk6XdJOkbZIuqLH/NEnXSipKekk7bXEWT7HcnaGhuVkEte9nRocWTimLhWAsl52ba9xiniBOEK8ZGfRksdMTtE0IJGWAi4AzgOOBsyUdX3XY74HXAl9olx3OA6dbQ0P1xlTGxB5BclFZPJEsDg0lt6Ul9gKOOnjUPQKnJ2inR3AKsM3MbjWzPPBFolXKFczsNjO7Huhc2YnTlEqyeEGvobhqqNMeQb0cQTyucqFHkAwNteoRJIVgtljmQGFxaxEcp1topxCsB+5IPN8RtrWMpHMlbZa0effu3UtinJOeYh2PIBaGTi0oq+QI6q4jyM47DubaTsdN54CaU8waEQvBkQePAr6WwFn5tDVHsFSY2SfMbJOZbVq3bl2nzek78nWSxYPZziaLZ5p4BJUB9okcweRsicGMGMpmGA/7Ww0N7Ut4BNC+hPHnfnE7F/1wW1vO7ThJ2ikEO4EjE883hG3OCqNe99FOj6psFhqKPYKZqqqh8eAJjA8Nhm2thXb2zhSQYP1BI5Xn7eA/r9vJV67Z0ZZzO06SdgrB1cBxko6WlANeDlzWxvdz2kQ8b6C6xcRgh6uGYo+g7jqCwYU5gmgWQSQEY8EjmJxt7UK+d6bA6uFB1o4MVp63gz2Ts+yenG3LuR0nSduEwMyKwHnAt4HfAl8ys62SLpR0JoCkx0naAbwU+FdJW9tlj7N48nU8gsyAGFAXrCOokyMYGBAjg5kFOYLYI4g9hlZbUe+dKbBmZJDV7RaC/bPsP1D0ZLTTdlINr18sZnY5cHnVtncnHl9NFDJyuphinQllEHkJneo1NF1oHBqC6K5/qmodQSwEAwNiLJdZVNXQmpFB1rRRCGbypYrd907lOWLtyJK/h+PErKH4wegAABbASURBVIhksdNZ4tBPddUQRI3oOhkakmAoW//PeDSXrYSQYH5oCKLmc5OLSBavHsmyOpSf7ptZ+uE0exIhoT0eHnLajAuB05R6VUMQVQ51LjRUYnQwU2kxXYvRqjv+ydki48NzQjA+nGVyEeWja0YGyWYGGB/KtsUjcCFwlhMXAqcp9bqPQtRvqHNVQ/UH18eMDWXntZiYPFBkPPGa8UUMsN87U6yEhdaMDLZJCPJzj/fnGxzpOA8cFwKnKcVyGSlKDlczmFEHF5TVb0EdE42rrCofHZ4vBK2EhswshIYiIVjdNiGY8wK8cshpNy4ETlPypXJNbwAiL6GjoaEUQjAdqoLKZWMqX1qYI2jBIzhQKJMvlRMeQbYtK4v37I8u/kPZAQ8NOW3HhcBpSrFkC6aTxWQz6livoWhMZWMhGMtlmS5EF/q5MZVzr2l1StneROfR+Hu7PILVw1ketGZ4XpjIcdpBW8tHnd6gWCovWEwWMzgw0NFeQ009gqE5jyBeWBavKI4et5YjWD4hyDMxPsRBY7mKd+A47cI9Aqcp+ZLVDw1l1dEWEyODje9lRnPZiicQryAeS3gEY0PZllpMVAvB6uH2eQQT40McMpbjnikXAqe9uBA4TSmWygtWFcdEVUOdm0eQJkdwoFCmVLbKCuLxoWSyOEO+VGa2mE4M9tXwCGYKJfLFpf0Z7JmcZWJVjolVQx4actqOC4HTlGLZai4mg7hqqHMeQTMhSLaiTg6uj4kfp/UKFoSGRqPv+w4srVcQh4Ymxoe4bzrfsfCb0x+4EDhNKTSpGip2qMVEmmRx3Ip6Jl+qtJuurhqC9MNpYiFYPTznESS3LwX5Ypm9MwUmxodYN57DLGoz4TjtwoXAaUqhVK7ZZwiiXkP5DngEZpYqWVyZS5wvVS72qxLrCOLHaWcSVIQgsY4guX0piHMCsUcAvpbAaS9eNeQ0pVhqEBoa6MyCstlimbLV7zwaMzeuslhJGtf0CFK2mdg7U2DVULayuK4dHkG8knhiPMdBY7lom+cJnDbiQuA0pVBuUDWUGehIjqAyi2AwbY6gNG9wfWV/eJx2dXFyVTHMCcFSLiqLF5BNrBrioNEgBF5C6rQRFwKnKYVig6qhjDpSNZSmBTXM5QjiZHF2QPO6la6KhaCFHMGaGkKwlB5BHAaaGBvioLHo/L662GknniNwmlIsl8nWyREMdmgeQTx+Ms3KYog8grgFdbJb6WKSxUkhiJPGe6fb4RHkGB/K9nSbiel8kVd/+pf89q59nTalr3EhcJpSKFllUH01nSofnZtX3GxB2VyOIDmdLCZuQLdYjyCXHWBkMLO0yeLJPKO5DKO5SLQmxoe4p0dzBFvv3MePb9nDD27c1WlT+hoXAqcpUdVQvdBQZxaUNRtcHxPvnw5VQ9VCMDeuMmWO4MB8IYAoPLSU6wjiVcUxE+O5nq0a2rZrct53pzO4EDhNaVY11IkWE80G18ckq4Ki0ND84zNhrnHaZPHemUJlEVnMUvcbioQgV3k+Md67q4u3BwHYvtuFoJO4EDhNKZSbLCjrYo9gKDvAgCLhmJwtMT48uOCY8eFsqvLR2WKJA4VyZURlzJILwf58lUcw1LM5gm1BALbvmsSsfTcUrQ4f6je8ashpSqOVxVFoqBM5gugfe7RJ0zlJjOWixnKTBwqsXzu84JioFXXzFhPV7SViVo8MsvP+mbSmN2XP5CyP3XhQ5fnEqhz3TuUpl42BOiG6lcr23ZMMKFrwd/e+Axy+ZmTJ3+NLm+/gbZdez0Gjgxx76DjHHjrOSevXcvYpRzYcc9pPuEfgNKVYMrJ1LkCDGXWmaqiQLjQUHxOVj5YqOYEkY0MZJlPE+PdVrSqOWTMyuGTrCIqlMvdOL/QISmXjvuneCg8dKJTYcd8Mpxx9MADbd0215X1+sf0e1owMcvqJDwLgW1vu5p1fu4Ff79jblvdbibgQOE1pXDU0gBmUlnk4TdrQEIRW03GyeHihEIynbEVdzyNYytDQvdN5zFiQI4DeW1186+4pzOC5J0QX6G279rflfbbeuY/HPvgg/u5Fj+TL/+tJfPONTwHghh33t+X9ViIuBE5TGlcNqXLMcjKdcmUxhLnFs0Um8wurhiD9lLL6oaHo9UuRK5lrLzHfI4DeW1QW5wee8JBDWDWcZfvupfcIDhRKbNs9yQlHrK5sO2LNMAeP5bhhp3sEMS4ETlOaTSiD5ReCmXyRkcFMqpj5WC7LPZOzmNEWIai0mUhZedSIymKyhBCsW5Wbt69X2L4ryg8cPTHGMevG21JCeuPd+ymVbZ4QSOLE9Wu4YacvYotxIXCa0rjXUHQhXu5FZWlmEcSMDmXYFXr1jNUQgrGU4yr3zUTH1BWCJQgPzXUe7f3Q0Pbdkxx58CjDgxmOPXS8LSWkW8Jd/wlHrJm3/ZHr13DLH/ZzoJB+Ol0v40LgNKXQaEJZplMeQfNZBDGjuQy7gxDU8wj2t+AR1EoWJ/c/ECqhoVVzHsHq4UGyA+o5j2DbrkmOWTcOwDHrxtm1f3bJB/xsvXMfa0YG2XDQ/GqkE9evoVg2b20RcCFwGlIqG2Y06DUUcgQdSBan9ghyWYrBvnpCkC+Wm4rZ3pkCY7nMAu9oSYVgcpZcdqDSDA9gYEAcMt5bQ+xLZeN3e6Y4Zt0YAMceGgnC9iUOD229cy8nHLF6QZnoSRsiD2GL5wkAFwKnCfHFcTBbr3w0+hNa7kVl04USI036DMWM5eYPq1+wP2Xjueo+QzFLKQS7J2dZNz604MLVa4vKdt43w2yxXBGAWBCWMmFcKJW58e798/IDMZ4wno8LgdOQihA0mFCWPG65mMkXGU1RMQQwWmNGcZJ4W7MpZXurZhHELK1HkJ+XH4jptTYTcT4gDg0ddfAogxktacJ4265J8sUyJ65fs2CfJ4zn40LgNCROAjfqNQQs++rilkJDCcGouY5gON2UsnpCsJTjKvfsn+WQRMVQTK95BPEFPxaCbGaAjYeMLWnCeOud0UW+lkcAcNL61Z4wDrgQOA2JVw036jUEy1811FKyeN5EsoWvSRsa2lcnNDQ8mCGXHViSqqHqhnMxE6ty3DOZb2s/nuVk++5JDhmbG8UJUZ5gKXMEW3buZWQww9ET4zX3nxQSxjfe3Z6FbCsJFwKnIfGdfqMJZQD5DiwoS+sRJHMEq4ZqNJ1rITRUSwhgaVYXl8vGPVPz20vErBsfIl8qV0pYVzrbdk1yzKHzL9DHrBvn9nunyReX5m/pN3fu4xGHr6rMl64mDhl5nsCFwGlCnARuNKEsedxyMZ0vNh1KExN7BAOC4cGFn2O84hE0DhHU8whgaWYS7J0pUCpbTSGorCWY6o3w0Pbdc6WjMcceOk6pbNx+zwNPGJfLxtY799bMD8SsXzvCQaODbPGeQy4ETmMqHkGDXkNApTxzuZgppA8NxR5B9ZjKyv4QLpqcrX8hL5TKTOVLbfUIkkPrq6kIQQ+UkN4zOct904VKxVBMLAxLkSe4/d5ppvKluvkBiBLGJ21Y6x4BLgROE+aqhronNBTV/FvqqqFYMFbVqBiKtkcX90atqPfVaS8RsxRCUBlaXydHAL2xujguEY1LRmMesoQlpFvvrL2iuJqT1q/mZk8YuxA4jZmrGmrca2g5k8Vpp5PFxK2na60hiLbPzTWuR70+QzFL4xFEF/l1NUJDh4z1TuO56oqhmLGhLEesGV6SEtItO/cxmBEPPWxVw+PihPFNfZ4wdiFwGjJXNVSnfDQb9xpaPo9guhCG0qRdUBYu9LVKRyESueHBgYaN51IJwfQDFIIQ9qlVPnrwWI4B9YYQbN89yfDgAOvXLhxCc8wS9RzaeudeHnrYKnJ1QpoxcQ7h+j4PD7kQOA0pFBuXj8ZJ5OUMDbUyiyA6LhKAWovJYpp1IK3XZyhm9cgg+2eLlB9ArmTP5CyZAbG2xntkBsTBY7meEIJtuyZ5yMR4zc6xx6wbf8BjK82MrXfua5gfiPGEcURbhUDS6ZJukrRN0gU19g9J+o+w/5eSNrbTHqd14iRwowllsEJCQw08iGYdSOc8gtrnWD2cxax5CWoj9kzOcshYrm5r7YnxIXbv74UcweSCRHHMsYeOV8ZWLpa79x3g3ql8w4qhmLkVxi4EbUFSBrgIOAM4Hjhb0vFVh70euM/MjgX+EfhAu+xxFsdcr6FmVUPd6xHEglEvNATBI2hwEa83pjJmKdpMRO0lFoaFYnphdfFMvsTO+2cW5Adi4u0PJE+wZWfjFcXVnLR+Td8njNs5vP4UYJuZ3Qog6YvAWcBvEsecBbwnPL4U+JgkWRuWT37p6jv45I9vXerT9jzxXXL9XkPR3euHvnMzn/rx75bFplaFIJcdYDCjhqGhsaEsP9t+D8/+8I9q7o/nBTfKEQC86tO/ZKhJXLoed9w3zSlHH1J3/8R4jqt+d29dG1cChVIZMxp6BABv/fKvWT1c+2fdjPtnCkjwiMPTC0GxbJz+kSvrhkC7hTc+8zhe8Kgjlvy87RSC9cAdiec7gMfXO8bMipL2AocAe5IHSToXOBfgqKOOWpQxa0cHOe6w2n98TmNOG67/s1s3PsTrTj2au/fNLKtNTz52guMPb+76x7zreY9g08aD6+5/7ZM28o3r72x4jmPXjTOUrS0+pxx9MC8+eQMzhcWHho47bJwXn7yh7v6Xn3IUhZJhrOw2E4998ME8+diJmvsmxnP876cdw20PcFHZCUesSV1M8JSHruOlj93QtNdUN1DvRuSBonb1LpH0EuB0M/vT8PzVwOPN7LzEMVvCMTvC8+3hmD21zgmwadMm27x5c1tsdhzH6VUkXWNmm2rta6cftBM4MvF8Q9hW8xhJWWANcE8bbXIcx3GqaKcQXA0cJ+loSTng5cBlVcdcBpwTHr8E+EE78gOO4zhOfdqWIwgx//OAbwMZ4DNmtlXShcBmM7sM+DTwOUnbgHuJxMJxHMdZRtqZLMbMLgcur9r27sTjA8BL22mD4ziO05jurpVyHMdx2o4LgeM4Tp/jQuA4jtPnuBA4juP0OW1bUNYuJO0Gbm/xZRNUrVbuQtzGpWMl2Ok2Lg1uY3oebGbrau1YcUKwGCRtrreirltwG5eOlWCn27g0uI1Lg4eGHMdx+hwXAsdxnD6nX4TgE502IAVu49KxEux0G5cGt3EJ6IscgeM4jlOffvEIHMdxnDq4EDiO4/Q5PS8Ekk6XdJOkbZIu6LQ91Ug6UtIPJf1G0lZJ53fapnpIyki6TtI3Om1LLSStlXSppBsl/VbSEzttUzWS/jz8nrdI+ndJw522CUDSZyTtCsOi4m0HS/qupFvC94O60Ma/D7/v6yV9TdLabrMxse8tkkxS7fFsHaSnhUBSBrgIOAM4Hjhb0vGdtWoBReAtZnY88ATgDV1oY8z5wG87bUQD/gn4bzN7OPAousxWSeuBNwKbzOxEovbs3dJ6/WLg9KptFwDfN7PjgO+H553kYhba+F3gRDN7JHAz8I7lNqqKi1loI5KOBJ4D/H65DUpDTwsBcAqwzcxuNbM88EXgrA7bNA8zu8vMrg2P9xNdvNZ31qqFSNoA/BHwqU7bUgtJa4DTiGZcYGZ5M7u/s1bVJAuMhIl8o0DjQcnLhJldSTQTJMlZwCXh8SXAC5fVqCpq2Whm3zGzeNjwL4gmIXaMOj9HgH8E3gbdOXC614VgPXBH4vkOuvAiGyNpI/AY4JedtaQmHyH6Qy532pA6HA3sBv5fCF99StJYp41KYmY7gQ8R3RXeBew1s+901qqGHGZmd4XHdwOHddKYFLwO+FanjahG0lnATjP7dadtqUevC8GKQdI48BXgTWa2r9P2JJH0fGCXmV3TaVsakAVOBv7FzB4DTNH5UMY8Qoz9LCLROgIYk/SqzlqVjjBCtivvZgEkvYsozPr5TtuSRNIo8E7g3c2O7SS9LgQ7gSMTzzeEbV2FpEEiEfi8mX210/bU4FTgTEm3EYXXniHp3zpr0gJ2ADvMLPamLiUShm7iWcDvzGy3mRWArwJP6rBNjfiDpMMBwvddHbanJpJeCzwfeGUXzjw/hkj4fx3+fzYA10p6UEetqqLXheBq4DhJR0vKESXmLuuwTfOQJKK49m/N7MOdtqcWZvYOM9tgZhuJfoY/MLOuupM1s7uBOyQ9LGx6JvCbDppUi98DT5A0Gn7vz6TLEtpVXAacEx6fA3y9g7bURNLpRCHLM81sutP2VGNmN5jZoWa2Mfz/7ABODn+vXUNPC0FIIp0HfJvoH+5LZra1s1Yt4FTg1UR32b8KX8/rtFErlD8DPi/peuDRwPs6bM88grdyKXAtcAPR/19XtB+Q9O/Az4GHSdoh6fXA+4FnS7qFyJt5fxfa+DFgFfDd8L/zf7vQxq7HW0w4juP0OT3tETiO4zjNcSFwHMfpc1wIHMdx+hwXAsdxnD7HhcBxHKfPcSFwuhJJpVAOuEXSf7XaVVLSeyS9NTy+UNKzlsCmEUk/Cs0M20ropPq/23j+50u6sF3nd1YWLgROtzJjZo8OXTrvBd6w2BOZ2bvN7HtLYNPrgK+aWWkJztWMtUBNIQgN6x4o3wReEFogOH2OC4GzEvg5oVmgpHFJ35d0raQbQkMvwr53SbpZ0k+AhyW2XyzpJeHxbXE/eEmbJF0RHj81saDvOkmratjxSsLq2np2SNoYZiF8Mswd+I6kkbDvcaFv/q9CH/0tYfsJkq4K26+XdBzR4q1jEsc+TdKPJV1GWDEt6c3BY9oi6U2J978xfOabJX1e0rMk/VTRXIFToNI76Aqi1gxOv2Nm/uVfXfcFTIbvGeDLwOnheRZYHR5PANsAAY8lWq07CqwO298ajrsYeEl4fBswER5vAq4Ij/8LODU8HgeyVfbkgLsTz+vZsZGo+dmjw74vAa8Kj7cATwyP3w9sCY8/StQnJ36fkXCeLYn3expRI72jw/P4844Fe7cSda6N3/8kohu9a4DPBNvOAv4zcc5XAh/t9O/avzr/5R6B062MSPoVc+2Pvxu2C3hfaCPxPSJP4TDgKcDXzGzaou6trfaU+inwYUlvBNbaXI/7mAkgOd+gnh0QNZb7VXh8DbAx5DhWmdnPw/YvJM71c+Cdkt4OPNjMZurYeJWZ/S48fjLR550ys0miBnZPSbz/DWZWJhKI75uZEQnHxsT5dhF1QXX6HBcCp1uZMbNHAw8muujGOYJXAuuAx4b9fwBaGfdYZO7vvvI6M3s/8KdEd+M/lfTwanuq3qeRHbOJ40pE3kNdzOwLwJnhPS6X9Iw6h041Ok+C5PuXE8/LVbYMh/d0+hwXAqersaij5BuBt4Qk6Rqi2QgFSU8nEgqAK4EXhsqeVcAL6pzyNqKwCsCL442Sjgl30R8g6lo7TwjM7D4go7kZw/XsqPc57gf2S3p82FQZUSnpIcCtZvbPRDmIRwL7iZqp1ePH4fOOKhrA88dhWys8lChc5fQ5LgRO12Nm1wHXA2cTDR7ZJOkG4DXAjeGYa4H/AH5NNKXq6jqney/wT5I2E92tx7wpJF2vBwrUnnT1HaKQDPXsaMLrgU+GkNcYsDdsfxmwJWw/Efismd1D5JlskfT31ScKn/di4CqiiXafCj+nVng6UfWQ0+d491HHSYmkk4E/N7NXL/L14yGej6QLgMPN7PyltLEFWw4DvmBmz+zE+zvdxVLUIztOX2Bm10r6oaSMLW4twR9JegfR/93twGuX1MDWOAp4Swff3+ki3CNwHMfpczxH4DiO0+e4EDiO4/Q5LgSO4zh9jguB4zhOn+NC4DiO0+f8/82Py3iCipyzAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "plt.plot(ss_rdf.bins, ca61_h2o_571)\n", "w570 = water[570]\n", "plt.xlabel('Radius (angstrom)')\n", "plt.ylabel('Radial distribution')\n", "plt.title('RDF between CA61 and {}{}'.format(w570.name, w570.resid));" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "If you are having trouble finding pairs of atoms where the results are not simply 0, you can use Numpy functions to find the indices of the nonzero values. Below we count the nonzero entries in the first `rdf` array." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "4374 4374 4374\n" ] } ], "source": [ "j, k, nbin = np.nonzero(ss_rdf.rdf[0])\n", "print(len(j), len(k), len(nbin))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Each triplet of `[j, k, nbin]` indices is a nonzero value, corresponding to the `nbin`th bin between atoms $j$ and $k$. For example:" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "0.04893740480237007" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ss_rdf.rdf[0][j[0], k[0], nbin[0]]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Right now, we don't care which particular bin has a nonzero value. Let's find which water atom is the most present around the alpha-carbon of residue 60, i.e. the first atom." ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "The water atom with the highest distribution around CA60 has index 568\n" ] } ], "source": [ "# where j == 0, representing the first atom\n", "water_for_ca60 = k[j==0]\n", "# count how many of each atom index are in array\n", "k_values, k_counts = np.unique(water_for_ca60, \n", " return_counts=True)\n", "# get the first k value with the most counts\n", "k_max = k_values[np.argmax(k_counts)]\n", "print('The water atom with the highest distribution '\n", " 'around CA60 has index {}'.format(k_max))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "You can also calculate a cumulative distribution function for each of your results with `ss_rdf.get_cdf()`. This is the actual count of atoms within the given range, averaged over the trajectory; the volume of each radial shell is not taken into account. The result then gets saved into `ss_rdf.cdf`. The CDF has the same shape as the corresponding RDF array." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "(2, 1041, 75)\n" ] } ], "source": [ "cdf = ss_rdf.get_cdf()\n", "print(cdf[0].shape)" ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "scrolled": true }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYIAAAEWCAYAAABrDZDcAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8li6FKAAAgAElEQVR4nO3deZxcVZ338c83nYSEpBOWhARIQgJEICCLRhZBQUUFBwiOwIArjiPzPCNuKDMyzsMAjo77qMAooMiiiIgbIoqKsopCWIUKS1gTukMghFSThWy/5497Ktw03dU3TVdXVdf3/XrVq+9W9/6quvv+7j3nnnMUEZiZWesaVu8AzMysvpwIzMxanBOBmVmLcyIwM2txTgRmZi3OicDMrMU5EVghkq6X9E/1jsM2naTHJR1a7ziscTkRNKn0z71S0guSFkm6SNLY3PqLJK2W1JVe90n6b0njc9ucKGld2kfldU4NYm2KJCLp3ZLmpu+hU9JvJB3UbZsTJYWkf+jh/ZtL+l9Jz0paJunG3DpJ+pKkJen1JUkajM9VTfo7+a9uy6anzzhc0gHp76ctt/6CXpZ9J00fJ+nPklZIur6HY7ZJ+i9JHWk/d0naIq07XtKD6ftbLOliSeO6vf94SfMkLZf0iKQ39HCM09NncAIswImguR0ZEWOBvYF9gNO6rf9yRLQDE4EPAvsDt0gak9vm1ogYm3udPCiRNxhJpwDfAL4ATAKmAf8LzOm26QeA54D397Cb84GtgN3Sz0/m1p0EHA3sBewJHAn888B9gpqZS3aeeE1u2RuAhd2WvRGoJL7nyL7LL/ayzzOB1wMHAOOA9wGr0rpbgAMjYjywIzAc2JCoJL0V+BLZ33N7Ou6j+Z1L2gk4Fugs/jFbmxPBEBARi4BryRJCT+tXRcTtwFHA1mT/RP2xk6TbJJUl/VLSVpUVkvZPV4HPS7pH0iFp+efJThznVO44JJ0p6ey0fkS6svtKmh8taVVl373tN60bL+l76er9qXSV2ZbWnSjpZklflbRU0mOSDu/pQ6W7pLOAj0TEzyJieUSsiYhfRcSpue12AA4mO6m/XdLk3Lpd0/d7UkQ8ExHrIuKO3GE+AHwtIhZGxFPA14ATe4lnS0lXS3omxX61pCm59ddL+pykW9IV9e8kTcitf5+kJ9Kdx2d7+2UWERFrgL+QnXCRtA0wErii27JXkRJBRPwhIq4AOnr6bMAngA9HxBORuS8iVqX3LoiIZ3NvWQfsnJs/EzgrIv4SEesj4qn0feadC/wbsPqVfPZW4kQwBKSTxOHA/GrbRUQX8HuyE3N/vB/4R2BbYC3wrXT87YFfk125bQV8GvippIkR8VngJuDk3B3HDcAhaZ+vAxaRTipkV4kPRsRz1fabtr0oxbEz2R3R24B8EdR+wIPABODLwPd6KY45ABgF/LzA558bET8F5gHvya3bF3gCODMVDf1N0rty63cH7snN35OW9WQY8H1gB7I7k5VA9yK7d5Ml9MqJ+dMAkmYB3ya7yt6OLPFP4ZW5kZd+P28Ebk6v/LLHImJhgX29mux3doyyIs2HJH0kv4GkgyQtA7qAd5HdXZCS/GxgoqT5khamC4vRufceC7wYEdf098O2IieC5vYLSV3AAmAx8J8F3tNBdlKt2D9dbVde+1d576Xp6m058P+A49I/53uBayLimnSV9nuyIoV39LKfW4GZkrYmO4l8D9heWR3HwWSJgmr7lTQp7f8T6Qp+MfA/wPG54zwRERdExDrgYrIENqmHeLYGno2ItVU+O2SJ4LI0fRkbFw9NAfYAlpGdgE8GLpa0W1o/Nq2rWAaM7SkxRcSSiPhpRKxIyfvz6XvJ+35EPBQRK8muzit3g8cAV0fEjRHxItnvaX0fn+vT+b8B4N5u628ADkqxvoEssd9K9rdTWXYDxUwBxpPdQcxI8Z6Rinwqn//mVDQ0BfgK8HhaNQkYkd7zBl4qEv0PAEntZEV7Hy8YiyVOBM3t6FQHcAiwK9mVb1+2JyvDrfhLRGyRe/2lynsX5KafIPunnEB25Xpst5PJQWQn3pdJJ6+5ZCe3N5KdRP4MHMjGiaDafndIx+/MrTuP7Aq5YlHumCvS5FhebgkwQdLw3j64pAPJTlyXp0WXAa+WVDkBrwTWAP8VEasj4gbgT2R3KQAvkJWHV4wDXogeen1UVul8XireKZNdkW+hXOVs/rMBK3Kfaztyv6eUtJf09rmSr+b/BsjqMPL+kva/B9nv66aIeCEdp7LsRopZmX6eFRErI+Jesu/0ZRcNqcjnt7z0nVfee3ZEdKYipK/n3nsG2cXK4wVjscSJYAhIJ52LgK9W2y5dcR9KdkXXH1Nz09PITnzPkp0QLu2WUMZERKWysKcubm8A3kx2RXd7mn87WRFL5aRSbb8LgBeBCbl14yKit+KWam5N+zq6yjYfAATcLWkR8Nfccnj5VTRs/LnvJ6sortgrLevJp4BdgP0iYhwvFcEUecqok9zvSdLmZHc8/ZbK728nq+DeNiIeSKtuSsv2pHgiqHxP+e+mWhfIw4GdUhxLySqpe3vvW4CPpSKnRWTfwxWS/q1gbC3LiWDo+AbwVkl7dV8haTNJrwV+ASwlK3/uj/dKmpVOLmcBV6Zilx8AR0p6u7JHA0dJOiRXwfk02RMgeTeQFa2UImI1cD1Z+f5jEfFM2qbX/UZEJ/A74GuSxkkaJmknSd2LUPoUEcuA04FzJR2drshHSDpc0pcljQKOI6sk3jv3+ijw7nQncSPwJHCasscuDwTeRFaJD3AJcIqk7SVtR3ayv6iXkNrJrn6fV1ZpXqTIr+JK4IhUzj6S7Pc0EP/nN5IVufw5t+zmtKwzIh6pLKz8rshO4sPS720EQNruJuCz6e9yN7LivKvTe98jaVqa3oGsWOy63DG/D3xU0jap4vmTlfeSJYI9eOn300H2ZNa5A/D5hzQngiEinTwvITuhVfxrqkNYktbdAbw+FRf0x6VkJ69FZJWrH0vHXkD2mOW/A8+QXa2fykt/X98kqxxcKulbadmfgdG8dCVZInuEcMOVZYH9vp+sorREluCupJfiqL5ExNeAU8jKmyvHOpkseR5NdmK+JCIWVV7AhWQnu8PS0zVzyIoplgEXAO/PXT2fB/wK+BtwH1kl+Hm9hPMNsu/mWbJimd9uwue4H/gIWdFVJ9n3UqQSty83kBW73ZxbdnNa1v0O831k39e3ycryV5J9HxUnkBXtLSH7Hv5fRFRO9rOAP0taTvYo6YPAh3Pv/RzZ3clDZBX2d5Eli0rdSv73sw5YmoqxrAr1UERpZmYtxHcEZmYtzonAzKzFORGYmbU4JwIzsxbXawOaRjVhwoSYPn16vcMwM2sqd9xxx7MRMbGndU2XCKZPn87cuXPrHYaZWVOR9ERv61w0ZGbW4pwIzMxanBOBmVmLcyIwM2txTgRmZi3OicDMrMU5EZiZtTgnAjOzBrduffCFa+Zxz4Lna7J/JwIzswb32LPLOf/GR3l4cW2GVnAiMDNrcKXOMgCzth3Xx5b940RgZtbgSh1lRrSJnbcZW5P9OxGYmTW4UmeZmdu0M3J4bU7ZTgRmZg2u1FFm1na1KRYCJwIzs4a2uGsVz77wYs3qB8CJwMysoc3r7ALwHYGZWasqdWRPDO3mOwIzs9ZU6iwzZcvRjB89ombHcCIwM2tgpY5lNa0fACcCM7OGtWL1Wh59dnlN6wfAicDMrGE9uKiLiNq1KK5wIjAza1AbupbwHYGZWWsqdZQZN2o4228xuqbHcSIwM2tQpc6sRbGkmh7HicDMrAGtWx880NlV0/YDFU4EZmYN6PEly1m5Zl3NK4rBicDMrCFVWhTXuqIYnAjMzBpSqTMbg2DmNu01P1ZNE4GkwyQ9KGm+pM/0sH6apD9JukvSvZLeUct4zMyaRamjzM41HIMgr2ZHkNQGnAscDswCTpA0q9tm/wFcERH7AMcD/1ureMzMmkmpszwo9QMAw2u4732B+RHxKICky4E5QCm3TQCVTzoe6KhhPGZmDWHRslVcfOvjrFsfPa5fs249z3S9OCj1A1DbRLA9sCA3vxDYr9s2ZwC/k/RRYAxwaE87knQScBLAtGnTBjxQM7PB9KPbnuTb1z/C6BFtvW6z5eYjOHDnrQclnlomgiJOAC6KiK9JOgC4VNIeEbE+v1FEnA+cDzB79uyeU6iZWZModZbZeZux/OGUg+sdClDbyuKngKm5+SlpWd6HgCsAIuJWYBQwoYYxmZnVXamjPCgNxYqqZSK4HZgpaYakkWSVwVd12+ZJ4C0AknYjSwTP1DAmM7O6WrZiDU89v3LQKoKLKFQ0JOn1wPT89hFxSbX3RMRaSScD1wJtwIURcb+ks4C5EXEV8CngAkmfJKs4PjEiXPRjZkPWYPUouin6TASSLgV2Au4G1qXFAVRNBAARcQ1wTbdlp+emS8CBmxCvmVlT25AImuyOYDYwy1fqZmavXKmjzMT2zZjYvlm9Q9mgSB3BfcDkWgdiZtYKBrOhWFFF7ggmACVJtwEvVhZGxFE1i8rMbAhavXY98xd3ccguE+sdykaKJIIzah2EmVkreHhxF2vWRfPdEUTEDZImAa9Li26LiMW1DcvMbOgZzK6lN0WfdQSSjgNuA44FjgP+KumYWgdmZjbUlDrLjB7RxvStx9Q7lI0UKRr6LPC6yl2ApInAH4AraxmYmdlQU+oos+u27bQNq+0YxJuqyFNDw7oVBS0p+D4zM0signkN+MQQFLsj+K2ka4Efpfl/oFsjMTMzq+6p51dSXrW24eoHoFhl8amS3sVLLYDPj4if1zYsM7OhZUNFcZPeERARPwV+WuNYzMyGrFJnmWGCXSc3USKQdHNEHCSpi6xvoQ2rgIiIxvs0ZmYNqtRRZsaEMYwe2ftgNPXSayKIiIPSz/bBC8fMbGgqdZbZZ9qW9Q6jR0XaEVxaZJmZmfVs2co1LFzaWGMQ5BV5DHT3/Iyk4cBraxOOmdnQM68BxyDI6zURSDot1Q/sKamcXl3A08AvBy1CM7Mm18hPDEGVRBAR/53qB74SEePSqz0ito6I0wYxRjOzplbqLDNhbGONQZBX5PHR30h6Y/eFEXFjDeIxMxtySh1ldm/QYiEolghOzU2PAvYF7gDeXJOIzMyGkNVr1/Pw4i4ObrAxCPKKtCw+Mj8vaSrwjZpFZGY2hMxf/AJr1gW7NWj9APSv87iFwG4DHYiZ2VDUiIPVd9fnHYGks3mpZfEwYG/gzloGZWY2VJQ6yowaMYwZExprDIK8InUEc3PTa4EfRcQtNYrHzGxIKXUuY9fJ4xpuDIK8InUEF0saCexKdmfwYM2jMjMbAiKCUkeZI/bart6hVFWkaOgdwHnAI2Qdzs2Q9M8R8ZtaB2dm1sw2jEHQwPUDUKxo6OvAmyJiPoCknYBfA04EZmZVNOpg9d0VeWqoq5IEkkeBrhrFY2Y2ZJQ6y0iw6+TG7sS52ngEf58m50q6BriCrI7gWOD2QYjNzKypVcYg2HxkoTHA6qZadPmGZE8DB6fpZ8haGJuZWRWlzjJ7T92i3mH0qdrANB8czEDMzIaSyhgE795vWr1D6VO1oqF/jYgvd2tQtkFEfKymkZmZNbEHmqBFcUW1oqF56efcKtuYmVkPSg0+GE1etaKhX0lqA14dEZ8exJjMzJpeqSMbg2Cb9savUq36+GhErAMOHKRYzMyGjFJnuSnuBqBYg7K7JV0F/ARYXlkYET+rWVRmZk1s9dr1PPz0C7xhZuOOQZBXpEHZKGAJ2UA0R6bXEUV2LukwSQ9Kmi/pM71sc5ykkqT7JV1WNHAzs0b1yDMvsHrd+iF1R/Dd7r2NSuqzuCjVL5wLvJVsDIPbJV0VEaXcNjOB04ADI2KppG02KXozswb00mD1jd2iuKJIIjgbeE2BZd3tC8yPiEcBJF0OzAFKuW0+DJwbEUsBImJxkaDNrHX84C9PcM+C5+sdxia5f8MYBGPrHUoh1doRHAC8Hpgo6ZTcqnFAW4F9bw8syM0vBPbrts2r0rFuSfs8IyJ+20MsJwEnAUyb1viNM8xsYKxZt56zflVisxHDaN+ssbtp6O6d+0xp6DEI8qp9syOBsWmb/P1NGThmAI8/EzgEmALcKOnVEbFR+o+I84HzAWbPnv2yxm1mNjRVytq/fMyeHL3P9vUOZ8iq1o7gBuAGSRdFxBMAkoYBYyOiXGDfTwFTc/NT0rK8hcBfI2IN8Jikh8gSgzu1M7Om6ca52RV5aui/JY2TNAa4DyhJOrXA+24HZkqakUY4Ox64qts2vyC7G0DSBLKiokeLBm9mQ9u8zjIjhw9jxwYe73coKJIIZqU7gKPJBqOZAbyvrzdFxFrgZOBasu4qroiI+yWdJemotNm1wBJJJeBPwKkRsaQfn8PMhqBSZ5ldJ7czvK3Iqcr6q0jtywhJI8gSwTkRsUZSoXL6iLgGuKbbstNz0wGckl5mZhtUxvt9++6T6x3KkFckzZ4HPA6MIavM3YGswtjMrGYWlVexdMUa1w8Mgj7vCCLiW8C3couekPSm2oVkZpZvlOVEUGvV2hG8NyJ+0K0NQd7XaxSTmdmGRLCrE0HNVbsjqFTTN0cbaTMbUkqdZaZvvTljm6whWTOq1o7gvPTzzMELx8wsU+oss7vrBwZFtaKhb/W2DjxUpZnVTteqNTyxZAXHvnZKvUNpCdWeGrojvUaRdTD3cHrtTdb9hJlZTTywqAtwi+LBUq1o6GIASf8XOCg1EEPSd4CbBic8M2tFLz0xNL7OkbSGIu0ItiTrcbRibFpmZlYT8zrLbDVmJJPGbVbvUFpCker4LwJ3SfoTIOCNwBm1DMrMWlups8ysbcchNUc3zs2uzzuCiPg+2TgCPwd+BhxQKTYyMxtoa9et54FFXezWJKN7DQWFHtCNiEXAL2sci5kZjz67nNVrm2e836HAXfqZWUNxRfHgcyIws4ZSqoxBMNFjEAyWQolA0kGSPpimJ0qaUduwzKxVlTrK7DKpnREeg2DQ9FlHIOk/gdnALsD3gRHAD4ADaxuamb0Sy1au4Qu/nseKNevqHcomufPJpRy553b1DqOlFKksfiewD3AnQER0SHJ1vlmDu/7Bxfx47gKmbbU5w4c1z2OY220xmnfsuW29w2gpRRLB6oiIyqhkaexiM2twpc4yI9uGcd2nDnYxi1VV5K/jCknnAVtI+jDwB+CC2oZlZq9UqaPMqyaPdRKwPhUZoeyrkt5KNjzlLsDpEfH7mkdmZv1WGe/3LbttU+9QrAkUqSw+BfixT/5mzeOZrhdZsnw1u3l0LyugyD1jO/A7STdJOlnSpFoHZWavzP2dHu/XiivS19CZEbE78BFgW+AGSX+oeWRm1m+V1rm7uZsGK2BTapEWA4uAJYALHs0aWKmzzNStRjNu1Ih6h2JNoM9EIOlfJF0PXAdsDXw4IvasdWBm1n/zOsouFrLCirQjmAp8IiLurnUwZvbKrVi9lseWLGfO3tvXOxRrEtUGrx8XEWXgK2l+q/z6iHiuxrGZWT88sKiLCI/3a8VVuyO4DDiCbAD7IBudrCKAHWsYl5n104ZunJ0IrKBqg9cfkX66p1GzJlLqLDN+9Ai2Gz+q3qFYkyhSWXxdkWVm1hhKHR7v1zZNr4lA0qhULzBB0paStkqv6YBrocwa0Lr1wQOLyi4Wsk1SrY7gn4FPANuR1RNULi/KwDk1jsvM+uGxZ5ezas16dy1hm6RaHcE3gW9K+mhEnD2IMZlZP5XctYT1Q5HeR8+WtAcwCxiVW35JLQMzs01X6igzok3svM3YeodiTaToUJWHkCWCa4DDgZsBJwKzBlPqLDNzm3ZGDvcYBFZckb+WY4C3AIsi4oPAXsD4IjuXdJikByXNl/SZKtu9S1JIml0oajPrUanDFcW26YokgpURsR5YK2kcWedzU/t6k6Q24FyyO4hZwAmSZvWwXTvwceCvmxK4mW1scdcqnn3hRdcP2CYr0tfQXElbkA1PeQfwAnBrgfftC8yPiEcBJF0OzAFK3bb7HPAl4NSiQZsNBd/4w0Pcn1oBD4RlK9YAblFsm65IZfG/pMnvSPotMC4i7i2w7+2BBbn5hcB++Q0kvQaYGhG/ltRrIpB0EnASwLRp0woc2qyxrVi9lm9e9zDbtG/GVmM2G7D9HrTzBPaassWA7c9aQ7VO515TbV1E3PlKDixpGPB14MS+to2I84HzAWbPnh2v5LhmjaDSMdzn5uzB23afXO9wrMVVuyP4WpV1Aby5j30/xcZ1CVPSsop2YA/g+tQUfjJwlaSjImJuH/s2a2ruGM4aSbUGZW96hfu+HZgpaQZZAjgeeHdu/8uACZX5NPjNp50ErBWUOsuMGzWc7bcYXe9QzAq1I3h/T8v7alAWEWslnQxcC7QBF0bE/ZLOAuZGxFX9CdhsKKg85umO4awRFHlq6HW56VFkbQrupECDsoi4hqwRWn7Z6b1se0iBWMyaXqVjuHfvu0O9QzEDij019NH8fHqU9PKaRWQ2xFU6hnP9gDWK/rRDXw54sBqzfnLHcNZoitQR/IrsKSHIEscs4IpaBmU2lLljOGs0ReoIvpqbXgs8ERELaxSP2ZDnjuGs0RSpI7gBIPUzNDxNbxURz9U4NrMhqdRR5pBdJtY7DLMNihQNnQScBawC1pONVBbAjrUNzWzoccdw1oiKFA2dCuwREc/WOhizoc4tiq0RFSmkfARYUetAzFpB5YkhjylsjaTIHcFpwJ8l/RV4sbIwIj5Ws6jMhqhSR5kpW45m/OgR9Q7FbIMiieA84I/A38jqCMysn+Z1ll0/YA2nSCIYERGn1DwSsyFuxeq1PPrsco7ca7t6h2K2kSJ1BL+RdJKkbSVtVXnVPDKzIebBNAaB7wis0RS5Izgh/Twtt8yPj5ptog1dS/iJIWswRRqUuV8hswFQ6vAYBNaYajYegZltrNTpMQisMdV0PAKzVrVy9TpOvfIelq1cs2HZ/U+Vee/+HoPAGo/HIzCrgTueWMrV93ay6+R2Ro9sA2CvqeM5cq9t6xyZ2csVuSPozuMRmPWh1LkMgMs+vD9bjRlZ52jMqvN4BGY1UOoos+34UU4C1hQ8HoFZDZTcgtiaSJFE8CTQGRGrACSNljQ9Ih6vaWRmTWrVmnU88sxy3r775HqHYlZIkZbFP2HjPobWpWVm1oOHnu5i3frwHYE1jSKJYHhErK7MpGkXfJr1wmMOWLMpkgiekXRUZUbSHMCD1Jj1otRZZuxmw5m65eb1DsWskCJ1BP8H+KGkc9L8QuB9tQvJrLmVOsrstm07w4a5BbE1hyINyh4B9pc0Ns2/UPOozJrU+vXBvM4yx7x2Sr1DMSuscIMyJwCzvj353AqWr17n+gFrKkXqCMysoHmVrqa3HV/nSMyKcyIwG0ClzjJtw8TMSWPrHYpZYb0WDUn6+2pvjIifDXw4Zs2t1FFm54ljGTWird6hmBVWrY7gyCrrAnAiMOum1Flm/x23rncYZpuk10QQER8czEDMmt1zy1fTuWyVWxRb0yn01JCkvwN2JxuYBoCIOKtWQZk1o3kek9iaVJ+VxZK+A/wD8FFAwLGAh1ky66bStcRuviOwJlPkqaHXR8T7gaURcSZwAPCq2oZl1nxKnWUmj/MYBNZ8iiSClennCknbAWuAQuPtSTpM0oOS5kv6TA/rT5FUknSvpOsk+U7Dmlapo+xiIWtKRRLB1Wmc4q+QDVr/OPCjvt4kqQ04FzicbFSzEyTN6rbZXcDsiNgTuBL4cvHQzRrHqjXrmP/MC64otqZUpK+hz6XJn0q6GhgVEcsK7HtfYH5EPAog6XJgDlDK7ftPue3/Ary3aOBmm+L0X963oTK3Fl5cuz4bg8B3BNaEqjUoe3NE/LGnhmWSijQo2x5YkJtfCOxXZfsPAb/pJZaTgJMApk2b1sdhzTb2/IrVXHLrE+w4cQyTx43q+w39MKJtGG+dNYkDd5pQk/2b1VK1O4KDgT/Sc8OyAW1QJum9wOx0zJcfLOJ84HyA2bNnx0Ad11pDKd0JnHHk7rzxVRPrHI1Z46nWoOw/08/+Nix7Cpiam5+Slm1E0qHAZ4GDI+LFfh7LrFfzOrsAP9Zp1ptqRUOnVHtjRHy9j33fDsyUNIMsARwPvLvbMfYBzgMOi4jFhSI220SljjLbtG/GxPbN6h2KWUOqVjTUnn7uArwOuCrNHwnc1teOI2KtpJOBa4E24MKIuF/SWcDciLiK7EmkscBPJAE8GRFH9bpTs34odfqxTrNqqhUNnQkg6UbgNRHRlebPAH5dZOcRcQ1wTbdlp+emD930kM2KW712PfMXd/GmXVw3YNabIu0IJgGrc/Or0zKzhvfw4i7WrPNjnWbVFOl07hLgNkk/T/NHAxfXLiSzgeP+f8z6VqRB2ecl/RY4KC36YETcVduwzAZGqbPM6BFtTN96TL1DMWtYhbqhjog7JC0gdUMtaVpEPFnTyMwGQKmjzK7bttM2TPUOxaxhFemG+ihJDwOPATeknz22ADZrJBGRPTHkYiGzqopUFn8O2B94KCJmAIeS9Qtk1tAWLl1J16q1rig260ORRLAmIpYAwyQNSx3Fza5xXGavWKVrCd8RmFVXpI7geUljgRuBH0paDCyvbVhmr9y8zjLDBLtOdiIwq6bIHcEcYAXwSeC3wCP03BGdWUMpdZSZMWEMo0e21TsUs4bWZyKIiOURsT4i1kbExcA5wGG1D83slcm6lhhf7zDMGl6viUDSOEmnSTpH0tuUORl4FDhu8EI023TLVq5h4dKVrh8wK6BaHcGlwFLgVuCfgH8HBBwdEXcPQmxm/VYZjcxPDJn1rVoi2DEiXg0g6btAJzAtIlYNSmRmr8BLXUu097GlmVWrI1hTmYiIdcBCJwFrFqXOMhPGbsY27bUZmtJsKKl2R7CXpMpo3wJGp3kBERG+57aGVerwGARmRVUbj8DP3Nmg+PldC7nw5scHdJ8PPt3l8YnNCirU6ZxZLf349gU89fxK9p66xYDt883jt+HofbYbsP2ZDWVOBFZXEUGpo8yRe23H59/56nqHY9aSirQsNquZp55fSdkdw5nVlROB1VXlMU83/DKrHycCq6uSO4YzqzsnAqsrdwxnVsVpDacAAAoSSURBVH9OBFZXpc6yB5Y3qzMnAqubDR3DuaLYrK6cCKxu5nkEMbOG4ERgdbPhiSHfEZjVlROB1Y07hjNrDE4EVjfuGM6sMTgRWF2sXruehxd3uX7ArAE4EVhdzF/8AmvWhe8IzBqAE4HVRclPDJk1DCcCq4tSR5lRI4YxY8KYeodi1vKcCKwuSp3L2HXyONqGqd6hmLU8JwIbdJUxCFw/YNYYnAhs0G0Yg8D1A2YNoaaJQNJhkh6UNF/SZ3pYv5mkH6f1f5U0vZbxWGNwi2KzxlKzRCCpDTgXOByYBZwgaVa3zT4ELI2InYH/Ab5Uq3iscczr7EKCXSe31zsUM6O2YxbvC8yPiEcBJF0OzAFKuW3mAGek6SuBcyQpImKgg7ni9gVccNOjA71b64eny6uYMWEMm4/0kNlmjaCW/4nbAwty8wuB/XrbJiLWSloGbA08m99I0knASQDTpk3rVzBbbD6CmZPG9uu9NrBmThrL22ZNrncYZpY0xSVZRJwPnA8we/bsft0tvG33ybxtd598zMy6q2Vl8VPA1Nz8lLSsx20kDQfGA0tqGJOZmXVTy0RwOzBT0gxJI4Hjgau6bXMV8IE0fQzwx1rUD5iZWe9qVjSUyvxPBq4F2oALI+J+SWcBcyPiKuB7wKWS5gPPkSULMzMbRDWtI4iIa4Brui07PTe9Cji2ljGYmVl1bllsZtbinAjMzFqcE4GZWYtzIjAza3Fqtqc1JT0DPLGJb5tAt9bKDcgxDpxmiNMxDgzHWNwOETGxpxVNlwj6Q9LciJhd7ziqcYwDpxnidIwDwzEODBcNmZm1OCcCM7MW1yqJ4Px6B1CAYxw4zRCnYxwYjnEAtEQdgZmZ9a5V7gjMzKwXTgRmZi1uyCcCSYdJelDSfEmfqXc83UmaKulPkkqS7pf08XrH1BtJbZLuknR1vWPpiaQtJF0p6QFJ8yQdUO+YupP0yfR7vk/SjySNqndMAJIulLRY0n25ZVtJ+r2kh9PPLRswxq+k3/e9kn4uaYtGizG37lOSQtKEesRWzZBOBJLagHOBw4FZwAmSZtU3qpdZC3wqImYB+wMfacAYKz4OzKt3EFV8E/htROwK7EWDxSppe+BjwOyI2IOse/ZG6Xr9IuCwbss+A1wXETOB69J8PV3Ey2P8PbBHROwJPAScNthBdXMRL48RSVOBtwFPDnZARQzpRADsC8yPiEcjYjVwOTCnzjFtJCI6I+LONN1FdvLavr5RvZykKcDfAd+tdyw9kTQeeCPZGBdExOqIeL6+UfVoODA6jci3OdBR53gAiIgbycYEyZsDXJymLwaOHtSguukpxoj4XUSsTbN/IRsJsW56+R4B/gf4V6Ahn84Z6olge2BBbn4hDXiSrZA0HdgH+Gt9I+nRN8j+kNfXO5BezACeAb6fiq++K2lMvYPKi4ingK+SXRV2Assi4nf1jaqqSRHRmaYXAZPqGUwB/wj8pt5BdCdpDvBURNxT71h6M9QTQdOQNBb4KfCJiCjXO548SUcAiyPijnrHUsVw4DXAtyNiH2A59S/K2EgqY59DlrS2A8ZIem99oyomDSHbkFezAJI+S1bM+sN6x5InaXPg34HT+9q2noZ6IngKmJqbn5KWNRRJI8iSwA8j4mf1jqcHBwJHSXqcrHjtzZJ+UN+QXmYhsDAiKndTV5IlhkZyKPBYRDwTEWuAnwGvr3NM1TwtaVuA9HNxnePpkaQTgSOA9zTgmOc7kSX+e9L/zxTgTkmT6xpVN0M9EdwOzJQ0Q9JIsoq5q+oc00Ykiaxce15EfL3e8fQkIk6LiCkRMZ3sO/xjRDTUlWxELAIWSNolLXoLUKpjSD15Ethf0ubp9/4WGqxCu5urgA+k6Q8Av6xjLD2SdBhZkeVREbGi3vF0FxF/i4htImJ6+v9ZCLwm/b02jCGdCFIl0snAtWT/cFdExP31jeplDgTeR3aVfXd6vaPeQTWpjwI/lHQvsDfwhTrHs5F0t3IlcCfwN7L/v4bofkDSj4BbgV0kLZT0IeCLwFslPUx2N/PFBozxHKAd+H363/lOA8bY8NzFhJlZixvSdwRmZtY3JwIzsxbnRGBm1uKcCMzMWpwTgZlZi3MisIYkaV16HPA+Sb/a1F4lJZ0h6dNp+ixJhw5ATKMl3ZA6M6yp1JPqv9Rw/0dIOqtW+7fm4kRgjWplROydeul8DvhIf3cUEadHxB8GIKZ/BH4WEesGYF992QLoMRGkDuteqV8DR6YuEKzFORFYM7iV1FmgpLGSrpN0p6S/pQ69SOs+K+khSTcDu+SWXyTpmDT9eKU/eEmzJV2fpg/ONei7S1J7D3G8h9S6trc4JE1PYyFckMYd+J2k0Wnd61K/+XenfvTvS8t3l3RbWn6vpJlkjbd2ym17iKSbJF1FajEt6ZR0x3SfpE/kjv9A+swPSfqhpEMl3aJsXIF9YUPfQdeTdc1grS4i/PKr4V7AC+lnG/AT4LA0PxwYl6YnAPMBAa8la627OTAuLf902u4i4Jg0/TgwIU3PBq5P078CDkzTY4Hh3eIZCSzKzfcWx3Syzs/2TuuuAN6bpu8DDkjTXwTuS9Nnk/WTUznO6LSf+3LHO4SsI70Zab7yecekeO8n67m2cvxXk13o3QFcmGKbA/wit8/3AGfX+3ftV/1fviOwRjVa0t281P3x79NyAV9I3Uj8gexOYRLwBuDnEbEist5bN7VPqVuAr0v6GLBFvNTHfcUEID++QW9xQNax3N1p+g5geqrjaI+IW9Pyy3L7uhX4d0n/BuwQESt7ifG2iHgsTR9E9nmXR8QLZB3YvSF3/L9FxHqyBHFdRARZ4pie299isl5QrcU5EVijWhkRewM7kJ10K3UE7wEmAq9N658GNmW4x7W89He/4X0R8UXgn8iuxm+RtGv3eLodp1ocL+a2W0d299CriLgMOCod4xpJb+5l0+XV9pOTP/763Pz6brGMSse0FudEYA0tsh4lPwZ8KlWSjicbG2GNpDeRJQqAG4Gj05M97cCRvezycbJiFYB3VRZK2ildRX+JrNfajRJBRCwF2vTSGMO9xdHb53ge6JK0X1q0YYhKSTsCj0bEt8jqIPYEusg6U+vNTenzbq5sAJ53pmWb4lVkxVXW4pwIrOFFxF3AvcAJZAOPzJb0N+D9wANpmzuBHwP3kI1SdXsvuzsT+KakuWRX6xWfSJWu9wJr6Hmkq9+RFcnQWxx9+BBwQSryGgMsS8uPA+5Ly/cALomIJWR3JvdJ+kr3HaXPexFwG9mIdt9N39OmeBPZ00PW4tz7qFlBkl4DfDIi3tfP949N5flI+gywbUR8fCBj3IRYJgGXRcRb6nF8aywD8TyyWUuIiDsl/UlSW/SvLcHfSTqN7P/uCeDEAQ1w00wDPlXH41sD8R2BmVmLcx2BmVmLcyIwM2txTgRmZi3OicDMrMU5EZiZtbj/DxPvvUXIj0wKAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "plt.plot(ss_rdf.bins, ss_rdf.cdf[0][0][568])\n", "w568 = water[568]\n", "plt.xlabel('Radius (angstrom)')\n", "plt.ylabel('Radial cumulative distribution')\n", "plt.title('RDF between CA60 and {}{}'.format(w568.name, w568.resid));" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## The site-specific RDF without densities" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "When the `density` of the selected atom groups over the box volume is not accounted for, your distribution values will be proportionally lower." ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [], "source": [ "ss_rdf_nodensity = rdf.InterRDF_s(u, ags,\n", " nbins=75, # default\n", " range=(0.0, 15.0), # distance\n", " density=False,\n", " )\n", "ss_rdf_nodensity.run();\n", "ss_rdf_nodensity.get_cdf();" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAZQAAAEWCAYAAABBvWFzAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8li6FKAAAgAElEQVR4nO29eZhcZZn3//l2VXenlyRAEgKEJQhRBBXUiAsuKCjgqGEcVHDDkXl53+sHI46OIzjv6ygjDjoKM+My76vCgAgigzhGxVFEARUEwiIkIBg2SViyAEk66XTXcv/+OM+pPl1dW3eq6nRX35/rqqtOPWep+/RS37qX535kZjiO4zjOrtKVtgGO4zhOZ+CC4jiO4zQFFxTHcRynKbigOI7jOE3BBcVxHMdpCi4ojuM4TlNwQXFSRdINkv4qbTucyeO/O6ccFxRnHJIelTQsaUjSU5IukTSY2H+JpFFJ28JjtaR/kjQ/ccyHJBXCNeLHV1tg64z4QJP0Xkmrws/hSUk/lfTasmM+JMkkvafC+f2Svi5pk6Qtkm5K7HujpF+F8UfbcDsNIekz4X7OKhs/K4x/Jrx+IHnPko4q/zmEsW2SsuH1IklXhHt+VtLliWO/KOlxSVslPSbpU2Xvn5H0OUlPhGveJWm3CvZfH+zINu2HMgtwQXEq8XYzGwSOAF4KnFO2/4tmNhdYBPwl8Crgt5IGEsfcYmaDiceZbbF8miHpY8C/AJ8HFgP7A18HVpQdeirwDPDBCpf5BrAH8MLw/DeJfduBi4FPNNXw5vAgE+/n1DAecxPw+sTr1wN/qDB2i5nlw+trgKeIfpZ7Al9KHHsRcIiZzQNeA7xP0jsT+z8bxl8NzAM+AOxMGijpfUB3Y7foJHFBcapiZk8BPyMSlkr7d5rZ7cA7gAVE4jIVDpJ0W/hW+UNJe8Q7JL1K0s2SnpP0e0lHh/HzgNcBX409IEmflfSVsL9b0nZJ/xxe90naGV+72nXDvvmSLgrexPrwjTYT9n1I0m8kfSl8O35E0gmVbip4becCZ5jZNWa23cxyZvYjM/tE4rgDgDcApwPHSdorse+Q8PM93cw2mlnBzO5I/A5uM7PLgIcb+UFL+s/geW6RdJOkwxL7LpH0NUk/Cd/eb5V0UGL/myX9IZz7VUB13u52oD9+j/A8J4zHlAvK64AvVBi7KVzjLcB+wCfMbEv4ed6V+Hk8YGbbE+cWgYPDubsDHwX+h5k9ZhGrzawkKOF39g/A39W5N6cCLihOVSTtC5wArK11nJltA64j+sefCh8EPgzsDeSBfwvvvwT4CfA5om/mfwt8X9IiM/t74NfAmQkP6Ebg6HDNVxB9i40/mF4NPGBmz9S6bjj2kmDHwUQe2luAZGjtlcADwELgi8BFkip9uL6a6AP0Bw3c/yoz+z5wP/C+xL4jgceAz4aQ172S/qLO9WrxU2AZ0Tf7O4HLy/afTPQtfnei3/t5AJIWEnkG/5vovh8Cjmrg/S5jzEs5NbxOchNwmKQ9JHUBy4HvAbslxo4Kx0HkDT8AXCpps6TbJb0heUFJZ0saAtYBA8AVYdeLiX6vJwVRfVDSGWX2fB74d6K/HWeSuKA4lfgvSduAx4ENRN/Y6vEE0YdzzKvCt//48aoa514WviluB/4P8O7gEbwfuNbMrjWzopldB6wC3lrlOrcAyyQtIBKSi4AlinJAbyASHGpdV9LicP2PBo9iA3Ah0QdtzGNm9k0zKwCXEgnh4gr2LAA2JUI11fggYx96VzA+TLQv8CJgC7APcCbRh+kL61yzImZ2sZltM7MR4DPA4Urkv4AfBK8nTyQ2sXf6VmCNmV1tZjmiMF4jH7rfAU6R1E30M/xOmT2PAX8i+jJyOPBHMxsGfpsY6wFuDafsSyTwvwL2Ar4M/DAIXnzN84G5wMuIBGxL4tz5wPOBA4GTgM9IejOApOVE4vWVBu7LqYALilOJE0OO5GjgEKJvpPVYQpQDiPmdme2WePyuxrmPJ7YfI4pfLwQOAN6VFCbgtUQf4BMIH0SriMTj9UQCcjPRh0RSUGpd94Dw/k8m9v0/om/0MaUPUjPbETYHmchmYGGtxK6ko4g+3K4MQ1cAL5YUf5APAzngc2Y2amY3En2YvqXaNWu8V0bS+ZIekrQVeDTsSv5+kyKxg7H72ofE78mirrLJ31tFzOxPRJ7O54nEotI5cdjr9UReJ8BvEmO3BQGE6OfxqJldFMJdVwY7xnlLIZx1Vzj+s4lzAc41s2Ezu4fo5/7W4Al9HTirgS8AThVcUJyqhA+vSxif9JxA8ACOZezDYLLsl9jen+gDdBPRB8VlZcI0EL6BAlRqlX0j8CaiUNXt4fVxRKGjOGxS67qPAyPAwsS+eWZ22MS3qsst4Von1jjmVKJcxN2SnmLsm/ip4fmeCudMtUX4e4mKAY4l+qa+NIzXy4UAPEni9xRCfPtVP3wc3wY+Hp4rEQvK6xj7G/p1YuymxLH3MPH+a/08skCcB4p/lsnj4+15hHBb+D3EeZ51kqYayp11uKA49fgX4M2SDi/fIalX0suB/wKeBf5jiu/xfkmHSuonSmJfHcJJ3wHeLum48O16jqSjQ24H4GngeWXXupEoZHSfmY0CNxDlPx4xs43hmKrXNbMngZ8DX5Y0T1KXpIPK4/SNYGZbgE8DX5N0oqLy325JJygqb50DvJsoGX9E4vHXwHuDZ3MTUUjoHEnZ4NG8kahYgmDfHCKvSuFeeqqYNJdI4DYD/UReQ6P8hCjX8c5g10eIQk6N8D0ij+qqKvtvIvoC8HqiUBfAvUSe2xsZLyg/AHaXdGr43Z1EFMr6bfhZ/E9JuyviSOAM4HoAM3uISKj+PvztvpAoDPdjxkKK8e8gDqu+nDGRd+rgguLUJHwIf5vogzHm70KOZXPYdwfwmrLqmslwGZEn9BRREvsj4b0fJ/pG/SlgI5H38AnG/m7/lSjB+qykfwtjNwN9jH0I3UdUFlr6UGrguh8kitvfRySUV1MlzFYPM/sy8DGiZHb8XmcSifCJRGGYb5vZU/GDqAw4Cxwf8hUriD7gtgDfBD5oZn8Ib/H6cI1riby7YSJBrMS3iUKK68O91QpDlt/HJuBdwPlEv/dljH341zt32Mx+EUKSlfY/SPSzecrMngtjReA2Is/h5sSxzxBVvf0t0c/jbGBFsA/gz4kKBrYRfXH4CuNzIqcQhTU3E4nk/zGz60OILPk7iL98PB2+mDgNIPMFthzHcZwm4B6K4ziO0xRcUBzHcZym4ILiOI7jNAUXFMdxHKcpzOpOmgsXLrSlS5embYbjOM6M4o477thkZovKx2e1oCxdupRVq1albYbjOM6MQtJjlcY95OU4juM0BRcUx3Ecpym4oDiO4zhNwQXFcRzHaQouKI7jOE5TcEFxHMdxmoILiuM4jtMUXFCc1NiZK/Cfqx7HO147TmfgguKkxg0PbOQTV9/Dg08PpW2K4zhNwAXFSY2RfAGIPBXHcWY+LihOauQKUagrXyymbInjOM3ABcVJjXwhEpJYWBzHmdm4oDipkSsJinsojtMJuKA4qVEKebmH4jgdgQuKkxpx7sQ9FMfpDFxQnNSIPRTPoThOZ+CC4qRG7Jl4lZfjdAYuKE5q5N1DcZyOwgXFSQ2v8nKczsIFxUmNsSovFxTH6QRcUJzUGKvy8pCX43QCLihOanjIy3E6CxcUJzXGenm5h+I4nYALipMaefdQHKejcEFxUmNsYqMLiuN0Ai4oTmqUJjZ6Ut5xOgIXFCc14tyJV3k5TmfgguKkhrdecZzOwgXFSQ0vG3aczsIFxUkN7+XlOJ2FC4qTGrmit15xnE7CBcVJjVzeW684TifhguKkhq/Y6DidhQuKkxp5b73iOB2FC4qTGqNe5eU4HYULipMaeW+94jgdRUsFRdLxkh6QtFbS2RX290r6Xth/q6SliX3nhPEHJB0XxvaT9CtJ90laI+msxPF7SLpO0h/D8+6tvDdn14lzKN56xXE6g5YJiqQM8DXgBOBQ4BRJh5YddhrwrJkdDFwIfCGceyhwMnAYcDzw9XC9PPBxMzsUeBVwRuKaZwPXm9ky4Prw2pnGjOY95OU4nUQrPZQjgbVm9rCZjQJXAivKjlkBXBq2rwaOkaQwfqWZjZjZI8Ba4Egze9LM7gQws23A/cCSCte6FDixRfflNAnv5eU4nUUrBWUJ8Hji9TrGPvwnHGNmeWALsKCRc0N47KXArWFosZk9GbafAhZXMkrS6ZJWSVq1cePGyd2R01TGqrzcQ3GcTmBGJuUlDQLfBz5qZlvL95uZARW/9prZN8xsuZktX7RoUYstdaphZokqL/dQHKcTaKWgrAf2S7zeN4xVPEZSFpgPbK51rqRuIjG53MyuSRzztKS9wzF7AxuadidO0ykk5p54DsVxOoNWCsrtwDJJB0rqIUqyryw7ZiVwatg+Cfhl8C5WAieHKrADgWXAbSG/chFwv5ldUONapwI/bPodOU0jOZnRq7wcpzPIturCZpaXdCbwMyADXGxmaySdC6wys5VE4nCZpLXAM0SiQzjuKuA+osquM8ysIOm1wAeAeyXdHd7qU2Z2LXA+cJWk04DHgHe36t6cXWc04ZW4h+I4nUHLBAUgfNBfWzb26cT2TuBdVc49DzivbOw3gKocvxk4ZhdNdtpE7JVILiiO0ynMyKS8M/OJW9b3d2e8l5fjdAguKE4qxCGvvp6seyiO0yG4oDipEIe8+nq6yBWMqBbDcZyZjAuKkwrxZMb+7iiNV/Cwl+PMeFxQnFQYzcceSgbwyY2O0wm4oDipEHsofd1BULz9iuPMeFxQnFSIPZL+4KH45EbHmfm4oDipkCtVeWXGvXYcZ+biguKkQqnKq9sFxXE6BRcUJxXinImHvBync3BBcVIhlx+b2AjuoThOJ+CC4qRC3G5lLOTlHorjzHRcUJxUGEvKR3+Cvmqj48x8XFCcVMiVWq94yMtxOgUXFCcVkt2GwUNejtMJuKA4qZArjm+94lVejjPzcUFxUmGsysvnoThOp+CC4qTCWLdhFxTH6RRcUJxUGEvKh5CXt693nBlP3TXlJS0C/gewNHm8mX24dWY5nU6pbNg9FMfpGOoKCvBD4NfAL4BCa81xZgv5gtEl6M16lZfjdAqNCEq/mX2y5ZY4s4pcsUg200U2I2CsjNhxnJlLIzmUH0t6a8stcWYVubzRk+miOxP9CXrIy3FmPo0IyllEorJT0rbw2Npqw5zOJl8sks2I7uCheMjLcWY+dUNeZja3HYY4s4tcwch2dZHNeC8vx+kUGsmhIOkdwOvDyxvM7MetM8mZDeQKRXrcQ3GcjqJuyEvS+URhr/vC4yxJ/9Rqw5zOJl+IkvLdXZ5DcZxOoREP5a3AEWZWBJB0KXAXcE4rDXM6m1zRyGZEV5fokvfycpxOoNGZ8rsltue3whBndpHLF+kJ+ZPuTJd7KI7TATTiofwTcJekXwEiyqWc3VKrnI4nHzwUiAXFPRTHmek0UuX1XUk3AK8IQ580s6daapXT8eQKRbIhf5LNyKu8HKcDqBryknRIeH4ZsDewLjz2CWOOM2XyBSsLebmH4jgznVo5lI+F5y9XeHypkYtLOl7SA5LWSpoQJpPUK+l7Yf+tkpYm9p0Txh+QdFxi/GJJGyStLrvWZyStl3R3ePjs/mlMrlAcC3l1yXMojtMBVA15mdnpYfMEM9uZ3CdpTr0LS8oAXwPeTOTZ3C5ppZndlzjsNOBZMztY0snAF4D3SDoUOBk4DNgH+IWk55tZAbgE+Crw7Qpve6GZNSR2TrrkikZ/Jg55dXkvL8fpABqp8rq5wbFyjgTWmtnDZjYKXAmsKDtmBXBp2L4aOEaSwviVZjZiZo8Aa8P1MLObgGcaeH9nGpMPExsBujMqLQnsOM7MpaqHImkvYAnQJ+mlRBVeAPOA/gauvQR4PPF6HfDKaseYWV7SFmBBGP9d2blLGnjPMyV9EFgFfNzMnm3gHCcFkkn57kxXaUlgx3FmLrWqvI4DPgTsC1yQGN8GfKqFNk2Vfwf+EbDw/GVgwiJgkk4HTgfYf//922mfkyBfGCsbjqq83ENxnJlOrRzKpcClkv7CzL4/hWuvB/ZLvN43jFU6Zp2kLNGkyc0Nnltu79PxtqRvAhX7jZnZN4BvACxfvtw/xVIiV/SJjY7TaTQysfFFkg4rHzSzc+ucdzuwTNKBRGJwMvDesmNWAqcCtwAnAb80M5O0ErhC0gVESfllwG213kzS3mb2ZHj558DqWsc76ZLLJyY2drmgOE4n0IigDCW25wBvA+6vd1LIiZwJ/AzIABeb2RpJ5wKrzGwlcBFwmaS1RIn2k8O5ayRdRdSMMg+cESq8kPRd4GhgoaR1wD+Y2UXAFyUdQRTyehT4nw3cm5MS+bBiI0Qhr1HPoTjOjKeRmfJfTr6W9CUikaiLmV0LXFs29unE9k7gXVXOPQ84r8L4KVWO/0AjNjnTg1zZxMbto4WULXIcZ1dptDlkkn6inIbjTJmoyitRNuweiuPMeOp6KJLuJQojQRS6WgTUy584Tk2iKq8Q8urq8l5ejtMBNJJDeVtiOw88bWb5FtnjzALMLFR5BQ8l2+XroThOB9BIDuWx0AzytUSeym+IFthynClRKBpmlDyU7i4x6lVejjPjaWQJ4E8TtUdZACwELpH0v1ttmNO5xJMYx01sdA/FcWY8jYS83gccHjeIDGvM3w18rpWGOZ1LPOckWeXlORTHmfk0UuX1BNH8k5he6sxad5xaxGufjFV5dfk8FMfpAGo1h/wKUc5kC7BG0nXh9ZupM2vdcWoRt6ofq/LyXl6O0wnUCnmtCs93AD9IjN/QMmucWUHcqr4U8vIqL8fpCOo1h3ScphNPYkyu2DhaKGJmRMvhOI4zE6kV8rrKzN5dNrGxhJm9pKWWOR1LnIDPJlZshKicOBYZx3FmHrVCXmeF57fVOMZxJk2clB9bsTESlHzRyGZSM8txnF2kVsjrybAu/CVm9sY22uR0OHHZ8NiKjZGwjBaKzOl2RXGcmUrNsuHQMr4oaX6b7HFmAaWy4XhiYygf9sS848xsGl0P5d5QNrw9HjSzj7TMKqejyZdPbMx2jRt3HGdm0oigXBMeSfyrpDNlxjyUuJdX9Oz9vBxnZtOIoOxmZv+aHJB0VrWDHaceueL4suH42UNejjOzaaT1yqkVxj7UZDucWUS+UDaxsVTl5R6K48xkas1DOQV4L3CgpJWJXfOI1n93nClRqvLKjK3YCDCadw/FcWYytUJeNwNPErWsT64rvw24p5VGOZ1Nedlw/OweiuPMbGrNQ3kMeEzSscCwmRUlPR84BLi3XQY6nceEkFeo8sp5DsVxZjSN5FBuAuZIWgL8HPgAcEkrjXI6mwkhrzAPJedVXo4zo2lEUGRmO4B3Al83s3cBh7XWLKeTyU1YsTGeh+IeiuPMZBoSFEmvJlq58SdhzPtjOFMmnsAYzz+JhSU3DXIohaKxev0WzFzcHGeyNCIoHwXOAX5gZmskPQ/4VWvNcjqZOLQV507iXEpuGqzaePNDm3jbV37Df65al7YpjjPjqDux0cxuBG5MvH4Y8LYrzpQpXwK4NLFxGqzauHloFIB//PF9vHbZQvbZrS9lixxn5lDVQ5H0L+H5R5JWlj/aZ6LTacS5ku7M+LLh6ZCU3zFaiJ5zBc6+5l4PfTnOJKjloVwWnr/UDkOc2UOuUKRLkAkeSinkNQ2S8jtG8wB89JhlfPm6B7lq1eO85xX7p2yV48wMas1DuSM831jtGMeZCrlisVTZBcleXul7KMPBQzn9Dc/jtw9t4nM/vp/XLVvkoS/HaYBaIa97Jd1T7dFOI53OIl+w0twTSFZ5TQMPJVcg2yV6sxm++BeHUzDz0JfjNEitkFe89O8Z4TkOgb0fb1/v7AK5QrFU4QXTq8preLRAX09UFb//gn4+efwh/MPKNdzy8GZec9DClK1znOlNvdYrSHqzmb00seuTku4Ezm61cU5nkitYKREPiYmN02Aeyo7RPP09Y9OsXrssEpGN20bSMslxZgyNTmw8KvHiNQ2e5zgVyReKpQ7DMFY+PD2S8gX6e8a+Zw2E7bj6y3Gc6jQiDKcBX5f0qKRHga8DH27k4pKOl/SApLWSJng0knolfS/sv1XS0sS+c8L4A5KOS4xfLGmDpNVl19pD0nWS/hied2/ERqf95ArFUskwjJUPT4ey4eHRAn3dYx5Kf2+0vX0kn5ZJjjNjqCsoZnaHmR0OHA4cbmZHmNmd9c6TlAG+BpwAHAqcIunQssNOA541s4OBC4EvhHMPBU4m6hl2PJGgxf/ll4Sxcs4GrjezZcD1eEhu2pIrWikRD1H5cJemRy+vyENJCEoQF/dQHKc+DYeuzGyLmW2ZxLWPBNaa2cNmNgpcCawoO2YFcGnYvho4RpLC+JVmNmJmjwBrw/Uws5uovMBX8lqXAidOwlanjeQLxVIfr5hspmta9PLakRtLykNkV0+2i+2j7qE4Tj1amQtZAjyeeL0ujFU8xszywBZgQYPnlrPYzJ4M208BiysdJOl0Saskrdq4cWMj9+E0mVzB6M5q3FhPpovcNFixcbgsKQ8w0JMpzU9xHKc6HZlct2jSQMVPJzP7hpktN7PlixYtarNlDkS5kuwED0XTpMprfFIeoL8ny/YRFxTHqUetNeXfWetEM7umzrXXA/slXu8bxiods05SFpgPbG7w3HKelrS3mT0paW9gQ53jnZTIF2xclRdE/bymQ5VXch5KzEBvptSSxXGc6tSa2Pj2GvsMqCcotwPLJB1IJAYnA+8tO2YlcCpwC3AS8Eszs9B88gpJFwD7AMuA2+q8X3yt88PzD+sc76REvji+ygugJ6Np0Xplx2ihlIiP6evJst1DXo5Tl1oTG/9yVy5sZnlJZwI/I1qQ6+Kwnsq5wCozWwlcBFwmaS1Rov3kcO4aSVcB9wF54AwzKwBI+i5wNLBQ0jrgH8zsIiIhuUrSacBjwLt3xX6ndYwWjL6eCkn5lAWlWDSGc4WKOZQdNcqGb3hgAy/Yay57z/d+X87spu56KACS/oyohHdOPGZm59Y7z8yuBa4tG/t0Ynsn8K4q554HnFdh/JQqx28Gjqlnk5M+UZVXWcgro9R7ee3MR15IX4UcynM7hiueY2acftkd/OVRSznnhBe23EbHmc7UTcpL+r/Ae4C/BkQkAAe02C6ng4lyKOUhr67UQ17xXJOB3sZzKCP5IqP5IluHcy23z3GmO41Ueb3GzD5INAHxs8Crgee31iynk8kViuMmNkLwUFJOyselwX1lOZT+nkzVHMpQCIUNeRWY4zQkKLGvv0PSPkAO2Lt1JjmdTq5CUj6q8poeHkqlsuFqOZShnUFQdrqH4jiN5FB+LGk34J+BO4kqvL7VUqucjqZS2XAU8krXQ4nDWhWT8rkCZkbUyGGM2EPxeSqO04CgmNk/hs3vS/oxMGeSLVgcZxxRyGvixMbRlNdDiT2U8nko/b1ZzGBnrjhh31jIy+epOE6tiY1vMrNfVprgKKmRiY2OU5Fc2YqNEJUNpz3XYyzkNdFDAdg+mp8gKNtdUBynRC0P5Q3AL6k8wbGRiY2OU5F8YXpObKwW8orLiHeMFGBw/DljIS8XFMepNbHxH8LzLk1wdJxycgWbGPKaBkn5UpVXWVI+6aGUEwvKNhcUx6kZ8vpYrRPN7ILmm+PMBqIqr4llw+kn5UPIq7xsuLf6qo2xZzKaL05YOMxxZhu1Ql5zw/MLgFcQ9cqCKARWr6+W41SkUDTMqDixMe31UIZzlZPysYdSaXJjcv7J9pE8u/X3tNBCx5ne1Ap5fRZA0k3Ay8xsW3j9GeAnbbHO6TjisFbFiY0pr4eyYzRPl6A3O17s4nkplUqD43koANt2uqA4s5tG/PPFwGji9ShVFq9ynHrEglJpxca010OJ10Ipn2vSX8NDSSbjfVVHZ7bTyMTGbwO3SfpBeH0iY0vtOs6kiPMklSY2TofWK+XhLoD+3jgpX8FDSQhK0ltxnNlIIxMbz5P038Brw9BfmtldrTXL6VTGQl7lVV5Kvcor8lAmCspACHkNV6nyksDM56I4TkPt683sDkmPE9rXS9rfzP7UUsucjiRuUT+xyms6tF4pTGgMCWPNIivlULaP5Fkw0MumoZEpt1/JFYqe0Hc6gkba179D0h+BR4Abw/NPW22Y05nEkxcrTWzMFYuYpScqw7l8RQ+lq0v0dVduYT80kmfxvN6wPbUGkZfe/CjHXnBjqvfuOM2gkaT8PwKvAh40swOBY4HftdQqp2OpGvLKdGEWlRWnRZyUr8RAb+UW9kMjefaaNydsT81DeWzzDjYNjabeesZxdpVGBCUXVkPsktRlZr8ClrfYLqdDiRPvlVZsBMinKCjVkvJQvYX90EiePYOgTLX9ypawONcWX6TLmeE0kkN5TtIgcBNwuaQNwPbWmuV0KmNVXhMnNkLkwcypkMdoB9WS8hCVDlebKT+/r5vebNeUk/IlQdmRY8luvi69M3NpxENZAewA/gb4b+AhKjeMdJy6jFab2Bg8ljRLh2sJykBvdoKgjOQL5ArG3DlZ5s7JTllQtobFubb6Il3ODKeRsuHYGykCl0rqAk4BLm+lYU5nUi0pH+dU0uw4PDyap6+78r9Ef09mgmDE804GejIM9GanPA/FQ15Op1DVQ5E0T9I5kr4q6S2KOBN4GHh3+0x0Ool8sU7IK6UcipmxI1cn5FWWdI/LhAfndDPYm51yDmWrC4rTIdTyUC4DngVuAf4K+BQg4EQzu7sNtjkdSNWQV3idS2nVxpF8EbOJjSFjBnqy7MiVeShBQAZ7g4cyBUExs5KQbHVBcWY4tQTleWb2YgBJ3wKeBPY3s51tsczpSEpJ+Qq9vIDU+nlVW60xpr93oocSC8hAb5bB3ixPb538v8ZwrlDKG7mH4sx0aiXlS3/dZlYA1rmYOLtKKYeSLe/llW5SvtpqjTEDPdkJzR+3lzyU7JRDXkkRcUFxZjq1PJTDJW0N2wL6wmsBZmbzWm6d03GUQl7lHkrXWNlwGlRbrTGmryfDzlyRQtHIhIq0bQlBiUJek5+Y6ILidBK11kNJZzKA09FU6zacTd1DqbxXjIsAABqQSURBVLxaY0zcIHLHaJ65c7qBhIdSKhuevCBs2eGC4nQOvl6p01biHEm1Kq+0yoYbyaHAmCcDY4Iy0JtloCfLzlxx0vZvDaXG8/u6PSnvzHhcUJy2Mho8kIlVXnHIKx0PZThUcNWq8oLxa6JsK81DyTLQW70jcS1ir2T/PfrdQ3FmPC4oTlvJV12xMYS8Uq/yqp5DgbIVGkei7sSZLjF3TnTe0CRXbRwvKL6eijOzcUFx2kq+iocSC0xaa6LUC3mN5VDGPJChkTwDvdF4/DzZSq9YUPbdvY+twzlvYe/MaFxQnLYyWqX1SlxGnHaVV70cSrJ0eGgkz2CZoGybZPuVrcM55s7Jslt/D6OFIjtzrbn/C657kKtWPd6SaztOjAuK01aqdRtOu2y4XshrbBng8Un5WFDm7oKHMr+vm/l93aXXreB7t/+JH/3+iZZc23FiWiooko6X9ICktZLOrrC/V9L3wv5bJS1N7DsnjD8g6bh615R0iaRHJN0dHke08t6cqZEvFpEozeWIicuI0wp5DY9Ga8PP6a78L9FfIYcShbyi8dhDmWz7lXYISrFobB4aZdPQaNOv7ThJWiYokjLA14ATgEOBUyQdWnbYacCzZnYwcCHwhXDuocDJwGHA8cDXJWUauOYnzOyI8PB+Y9OQ0UJxgncCYx5Lmh5KX3cGSRX3x4IyPodSYLA3EoLBXRSUeX3Z0utms2U4R75obBoaafq1HSdJKz2UI4G1ZvawmY0CVxKtrZJkBXBp2L4aOEbRf/QK4EozGzGzR4C14XqNXNOZxuQLNmG1RkhWeaWUlK/RaRgSSffR8VVeg8FDGZzGIa/N2yMheWb7KMUUV8R0Op9WCsoSIJkFXBfGKh5jZnlgC7Cgxrn1rnmepHskXSipt5JRkk6XtErSqo0bN07+rpxdIl8oTlhPHpJVXukl5avNQQHozXbRpfE5lEpVXpNdE2XrcI55c8YEpRWTGzdui0JdhaLxnM91cVpIJyXlzwEOAV4B7AF8stJBZvYNM1tuZssXLVrUTvscoomNFUNe2bRDXnn6qyyuBSApahA5Ml5QBsP8k55sFz3ZrinNQ5nf31oPJRnq8rCX00paKSjrgf0Sr/cNYxWPkZQF5gOba5xb9Zpm9qRFjAD/QRQec6YZ+UJxQh8vSH8J4B11PBQILeyDYIzmi4zmiwwmqsIm23F4Z67ASL7I/L7uUn+wlgvKNhcUp3W0UlBuB5ZJOlBSD1GSfWXZMSuBU8P2ScAvLZrZtRI4OVSBHQgsA26rdU1Je4dnAScCq1t4b84UyRdtwqRGGEvKpzmxsVYOBaKS4rj1SrIxZMzgJJcBjsNb8/q6S7PtWy0oG91DcVpI3TXlp4qZ5cOSwT8DMsDFZrZG0rnAKjNbCVwEXCZpLfAMkUAQjrsKuA/IA2eENVmodM3wlpdLWkTUXv9u4H+16t6cqVOtyivTJbqUbpXX7v3dNY+JlgGOBCO5uFbMZFvYx+IRh7ta1SBy07ZR5nR3sTNX9NJhp6W0TFAAzOxa4NqysU8ntncC76py7nnAeY1cM4y/aVftdVpPvlCc0McrJpvpSq2X1/BovupaKDEDPdlS2fBQYi2UmMHezKRa2FcSlFZ5KEsXDLB2w5DnUJyW0klJeWcGkC9UDnkBdHcp3ZBXlbVQYpI5lO0VBSU7qW7D7RSURXN7WTDY4zkUp6W4oDhtpVrIC6JKrzR7edVNyvdkSjmUbVVCXpNJypcLyrw5rRKUURYN9rJgoJfN2z3k5bSOloa8HKecfMEqVnlB1M8rjSovM6s7sRGipHycQ6nmoWybhKBsrZRD2dlcQTGLZsgvnNvLwu2jHvJyWop7KE5byReLE9aTj+nOKJWJjaOFaK34eoIykPBQqlV5Tc5DiY6dF64xv7/5HsrQSJ6RfJGFgz0s9JCX02JcUJy2kitYaRJjOd2ZLvIptAaJZ7/XS8r392ZLx8Zt6pPzUAZ6o6R9ocF72DKcY6AnU+ocML+vm525IiP5ya36WIu4qmvhYC+LBnvZNDTqa644LcMFxWkruUKxYi8viPp5jabgodRbXCtmoCfDaCGa0Bgn3+Nuw0Bp1cbtDc6Wj/t4xcxrwWz5OMS1cLCXhYO9jBaKpXXsHafZuKA4baV2lVdXKiGvRgWlL7EmyvbRPHO6u8b1JZvsqo1bhnMlEQFa0s8rDnEtHOxl4dyeaMzzKE6LcEFx2kquWKvKK52y4VLIq07Z8EDP2KqN23bmxyXkYfINIreWeSit6OdV8lDm9rBwMOqX6nkUp1W4oDhtJVejbDjb1ZVSyCsSgGqrNcb0946tK59crTFm7iTXRCkPebVCUDYOjSLBHv0JQfHZ8k6LcEFx2kq+YKVGkOVEVV7t91B25OKkfGMeyo7R/LjW9aX9pZBXY0n1CTmUOc1fZGvT0Ai79/eQzXQlBKUzPZSnt+7kPf/vFjZs3Zm2KbMWFxSnrdSv8mq/hzLccA4lXga4ELWuLxOUsVUbGxOEqh7KjuYJyuahERYORrmT3fu7kaKxTuTmhzZx6yPPcMvDm9M2ZdbiguK0ldpVXl2MpuGhNFzlFYe88hVDXmOCUt9DGc0XGc4VKlZ5NbMKa9PQaMkzyWa62KO/h40dGvJau2Fo3LPTflxQnLZSbcVGiHt5peGhRB/gdUNevXFSvlAl5BXtH2pgtns8I35+osNxd6aLgZ5M00NesaBAVO3VqSGvhzZsj543uqCkhQuK01ZyxcorNkIIeaXqodRJypfKhoOHMqfMQynNQ6nvocSiMW/O+Jb5zW4QuWlbmaDM7elYQVkbhCQWllZgZpPqhjDb8F5eTlvJVVmxEaKJjWk0h9zRcNnwWNK9UtlwbzZDd0YNVXmVN4aMmddEQYnmyxRK808g8lDu/NOzTbn+dCJXKPLY5u10CR7ZtL2mJ7wrXPiLP/Jv1/+RRXN7OXjRIMsWD3LkgXvwtpfs0/T3mom4h+K0jULRMKNGL6901kMZzhXozXaRqZLbiYlDYlt35hjJF0sCk2SgwVUbtyRWa0zSTA8lOUs+ZuFgL5u2dV4O5U/P7CBXMI48cA9GC0XWPTvckvf53UOb2Xf3Pt7w/EUM5wpcc+d6zrziLjb63B7ABcVpI7H30Z2dZmXDo/m6CXmAnmwX3RmxIXx4lIe8oPEGkeWdhmOauWrjxpKgjPdQhnOFjgvbPBQS8ccdthfQmsR8sWjc9+RW3nTInnzpXYfzX2ccxUWnLgdg9fotTX+/mYgLitM2SoJSa8XGlEJe9fInMf09WTZsDYLSO1GEBnuz0ybklWy7EhOLS6flUeL8yVuCoLQiMf/YMzsYGslz2D7zSmOHLZmPBPe6oAAuKE4bib2PWis2prEeyvBo/bVQYgZ6MmzcFk2cG+yduAZ9w4Kyo7qH0ryQ11in4ZiFcztzcuNDG7azeF4vS3brY+Fgb0s8lDVPRKJx2D7zS2ODvVkOXDjgghJwQXHaRpwfqV3llZaH0pig9PdmSyGvgQoeSqOrNm7dmaOvO0NP2STP+X3d7BgtNMVTiycwLkiEvBZ1aPuVhzYOcfCegwAcvOdASzyU1eu30p0RyxYPjht/yZL5HvIKuKA4bSP2PqpXeaWzYmMjy//G9PdkSgnY8iqveKyRVRvLZ8nHNLPj8KahEebNydKbHbu3BR0Y8jIzHtowxEGLog/6gxYNsnbDUNPXfVnzxBaW7Tl33M8T4EVL5vPklp2emMcFxWkjsfdRa8XGNKq8duTyk8ihZEqLgO1KUr6eoDQj7LVpaLQU4opZMBB3HO4cD2XjthG2jeRLgnLwnoNs3ZlvqhdmZqx5YisvWjJvwr4XL4lCYO6luKA4baTkodTo5WVGwyseNosdk/BQkqXC1cqGG2kO2Q5B2Vg2Sx6iSrX5fd0d5aHE+ZI45BULSzPDXk9t3ckz20fH5U9iPDE/hguK0zbGqryqT2xMHtcuhkcL9NeZ1BjTnwhzVQ55ZRgayVOsI4pbhvPM65t4fjNXbdw0NFLKmSRZONhZs+Vj4Uh6KNDc0uHV67cCVPRQPDE/hguK0zbGqryq9fKKxtstKJNKyieEp7yXF4yFweKW+NXYWrZaY0xTQ17bRsYl5GM6rZ/X2g1DDPZmWTwvEs+958+hvyfTVA9lzRNbkOCQvSYKCkRhLw95uaA4bWSsyqv6xEag7ZMbo6R8gzmUUNnVk+2aUKEFjS8DXC3kFXstu5qUH8kX2LozPyHkBVHpcCdVeT20cTsHLRpAiv5+JJUS881i9fqtPG/hQMUvERAJypNbdnaUUE8FFxSnbeTytcuGY8+lnR5KvlBktFCcxDyU6ANlbpUPljgMtq1G+5V8ocjQSL6lOZTNFeagxCwa7O2oZYDXbhjioD3Hl/IetGiAhzc2r0nkfU9sqZg/iXlRSMzP9rCXC4rTNuLqqForNkLUkbhdxKGpxuehRMdV+6Y62ICHEotNJUHpzWaY0921y2uijAlKpZBXD9tG8uysE5abCQyN5Hlq685S/iTm4D0HWf/ccFNazDyzfZQntuysmD+JiWfPr17nguI4bWGsl1f1Ki+grZMb49UaJ1vlVU1QBhpYV75a25WY+X3du7xqY6kx5NxKSflobPP2mR/2eqiswismFphHNu26l1Jphnw5c+d087xFnph3QXHaRqlsuEYvr+i49glKo6s1xsTCUy/ktcuCsoshr7gxZOUqr3guyswPe5VXeMXEAtOMxPyaJ6IKr2QPr0p4Yt4FxWkjpYmNNXp5AW2dLR+HRPq6G0vKj3kolQWokZBXOwRlU4W2KzGdNFt+7YYhsl3igAX948YPWDBApktNScyvXr+FJbv1sVv/xJ9lkhcvmc8Tszwx74LitI04N1Krlxe0t8preIo5lME5lcUgLhtO20PZtG2U/p5MxQ4AJQ+lAz74Hto4xAEL+if8TfVkuzhgj/6meCj3VZkhX44n5l1QnDYyVuVVe2Lj6DQOecUeSqXW9dH4rgtKM1rYl68ln2TR3M5pELk20cOrnOc1oXR4aCTPw5u218yfxHhivsWCIul4SQ9IWivp7Ar7eyV9L+y/VdLSxL5zwvgDko6rd01JB4ZrrA3XrO2fOm0nX4xDXtMpKR9CXpNoDgmV264ApZUfGwl5VZrYCNE687s6DyUSlMr/AnO6Mwz2Zmd8M8No2d8dExLyMQfvOcijm3bs0t/T/U9WnyFfztw53Txvls+Yb5mgSMoAXwNOAA4FTpF0aNlhpwHPmtnBwIXAF8K5hwInA4cBxwNfl5Spc80vABeGaz0bru1MI+p1Gy4JSjvLhkseSuPNIaFyY0iIJtUN1lkGeOtwjp5sF3OqtHuZ39fNtpH8LvU0q+WhQGe0X3ls8w7yRavqoRy0aIDRQpHHd2E54DjJ3oiHAlHYazYn5hv7L5oaRwJrzexhAElXAiuA+xLHrAA+E7avBr6qaLrrCuBKMxsBHpG0NlyPSteUdD/wJuC94ZhLw3X/vRU39pXr/8jK3z/Rikt3NM+Fb93Vq7wiofm7q+9pOAS1q8TeQsMhr9445FX9X2ewN8sP7lrPzQ9trrh/w7aRquEuGAuFveXCG+lS7XXuq/HIpu0sX7pH1f0LB3v55R828OYLbpzS9acD8ZeBWh4KwPu/deuU/57iBpt7Vii/rsSLl8xn5e+f4NgLbmRqv7n28fl3vphX1PgbmQqtFJQlwOOJ1+uAV1Y7xszykrYAC8L478rOXRK2K11zAfCcmeUrHD8OSacDpwPsv//+k7ujwKK5vRMW2XEaY5/5fezWX/nD9IV7zeM9y/dj20hzVixslL3m9TX8gbHn3F4+euyy0trllfhfRx/ELQ9tqrp/2eJBjqzxj3zsCxfz+3XP7VL59PP3mstJL9+36v7TXnsgP7pn5n8pOvoFizi0Sjnvi5bM532v3J9nd0w9V7Rs8SCvW7ao1NalHn/2kr1Z88SWtuYBp0pfgw1RJ4OavQhN6cLSScDxZvZX4fUHgFea2ZmJY1aHY9aF1w8RCcRngN+Z2XfC+EXAT8NpE66ZOP7gML4f8FMze1EtG5cvX26rVq1qzg07juPMEiTdYWbLy8dbmZRfD+yXeL1vGKt4jKQsMB/YXOPcauObgd3CNaq9l+M4jtNCWikotwPLQvVVD1GSfWXZMSuBU8P2ScAvLXKZVgInhyqwA4FlwG3VrhnO+VW4BuGaP2zhvTmO4zhltCyHEnIiZwI/AzLAxWa2RtK5wCozWwlcBFwWku7PEAkE4biriBL4eeAMMysAVLpmeMtPAldK+hxwV7i24ziO0yZalkOZCXgOxXEcZ/KkkUNxHMdxZhEuKI7jOE5TcEFxHMdxmoILiuM4jtMUZnVSXtJG4LFJnrYQqD4NenrgNjaPmWCn29gc3MbGOcDMFpUPzmpBmQqSVlWqbphOuI3NYybY6TY2B7dx1/GQl+M4jtMUXFAcx3GcpuCCMnm+kbYBDeA2No+ZYKfb2Bzcxl3EcyiO4zhOU3APxXEcx2kKLiiO4zhOU3BBmQSSjpf0gKS1ks5O255yJO0n6VeS7pO0RtJZadtUDUkZSXdJ+nHatlRC0m6Srpb0B0n3S3p12jaVI+lvwu95taTvSpqTtk0Aki6WtCEsoBeP7SHpOkl/DM+7T0Mb/zn8vu+R9ANJu003GxP7Pi7JJC1Mw7ZquKA0iKQM8DXgBOBQ4BRJh6Zr1QTywMfN7FDgVcAZ09DGmLOA+9M2ogb/Cvy3mR0CHM40s1XSEuAjwPKwMmmGsPzDNOAS4PiysbOB681sGXB9eJ0mlzDRxuuAF5nZS4AHgXPabVQZlzDRxnhF2rcAf2q3QfVwQWmcI4G1ZvawmY0CVwIrUrZpHGb2pJndGba3EX0ILknXqolI2hf4M+BbadtSCUnzgdcT1tQxs1Ezey5dqyqSBfrCSqX9wLRYJN7MbiJa3yjJCuDSsH0pcGJbjSqjko1m9nMzy4eXvyNa+TU1qvwcAS4E/g6YdhVVLiiNswR4PPF6HdPwwzpG0lLgpcCt6VpSkX8h+ocopm1IFQ4ENgL/EcJy35I0kLZRScxsPfAlom+pTwJbzOzn6VpVk8Vm9mTYfgpYnKYxDfBh4KdpG1GOpBXAejP7fdq2VMIFpQORNAh8H/iomW1N254kkt4GbDCzO9K2pQZZ4GXAv5vZS4HtpB+iGUfIQawgEr99gAFJ70/XqsYIS3ZPu2/XMZL+nih8fHnatiSR1A98Cvh02rZUwwWlcdYD+yVe7xvGphWSuonE5HIzuyZteypwFPAOSY8ShQ3fJOk76Zo0gXXAOjOLvburiQRmOnEs8IiZbTSzHHAN8JqUbarF05L2BgjPG1K2pyKSPgS8DXifTb9JegcRfYH4ffj/2Re4U9JeqVqVwAWlcW4Hlkk6UFIPUQJ0Zco2jUOSiOL+95vZBWnbUwkzO8fM9jWzpUQ/w1+a2bT6Zm1mTwGPS3pBGDoGuC9FkyrxJ+BVkvrD7/0YplnhQBkrgVPD9qnAD1O0pSKSjicKxb7DzHakbU85Znavme1pZkvD/8864GXh73Va4ILSICFZdybwM6J/3KvMbE26Vk3gKOADRN/67w6Pt6Zt1Azlr4HLJd0DHAF8PmV7xhG8p6uBO4F7if6Xp0VbDknfBW4BXiBpnaTTgPOBN0v6I5F3df40tPGrwFzguvC/83+noY3TGm+94jiO4zQF91Acx3GcpuCC4jiO4zQFFxTHcRynKbigOI7jOE3BBcVxHMdpCi4oTscjqRDKQFdL+tFku8hK+oykvw3b50o6tgk29Um6MTQdbSmhc/L/18Lrv03Sua26vjNzcEFxZgPDZnZE6Mr7DHDGVC9kZp82s180waYPA9eYWaEJ16rHbkBFQQmNJXeVnwBvD61BnFmMC4oz27iF0NRT0qCk6yXdKene0HiPsO/vJT0o6TfACxLjl0g6KWw/Gq9HIWm5pBvC9hsSE0vvkjS3gh3vI8wWr2aHpKVhLZZvhnVPfi6pL+x7RVi34+6wjsfqMH6YpNvC+D2SlhFNIjwocezRkn4taSWhA4CkjwUPbrWkjybe/w/hnh+UdLmkYyX9VtG6JkdCqTfXDUQtS5zZjJn5wx8d/QCGwnMG+E/g+PA6C8wL2wuBtYCAlxPNPu8H5oXxvw3HXQKcFLYfBRaG7eXADWH7R8BRYXsQyJbZ0wM8lXhdzY6lRE0Kjwj7rgLeH7ZXA68O2+cDq8P2V4j6UMXv0xeuszrxfkcTNbw8MLyO73cg2LuGqFN1/P4vJvryeQdwcbBtBfBfiWu+D/hK2r9rf6T7cA/FmQ30Sbqbsbbp14VxAZ8P7VV+QeS5LAZeB/zAzHZY1K15sj3bfgtcIOkjwG42tsZGzEIgub5KNTsgagB5d9i+A1gackBzzeyWMH5F4lq3AJ+S9EngADMbrmLjbWb2SNh+LdH9bjezIaJGk69LvP+9ZlYkEprrzcyIBGhp4nobiLoeO7MYFxRnNjBsZkcABxB9eMc5lPcBi4CXh/1PA5NZRjfP2P9Q6TwzOx/4KyLv4LeSDim3p+x9atkxkjiuQOTNVMXMrgDeEd7jWklvqnLo9lrXSZB8/2LidbHMljnhPZ1ZjAuKM2uwqIPsR4CPh2T0fKK1WXKS3kgkOAA3ASeGSqy5wNurXPJRonARwF/Eg5IOCt/qv0DUpXqcoJjZs0BGY2vAV7Oj2n08B2yT9MowVFr6V9LzgIfN7N+IcjQvAbYRNT2sxq/D/fYrWkjsz8PYZHg+URjOmcW4oDizCjO7C7gHOIVoAaXlku4FPgj8IRxzJ/A94PdEq/bdXuVynwX+VdIqIu8h5qMhuX0PkKPyyn8/Jwo1Uc2OOpwGfDOE8gaALWH83cDqMP4i4NtmtpnIU1ot6Z/LLxTu9xLgNqIVPr8Vfk6T4Y1E1V7OLMa7DTtOCkh6GfA3ZvaBKZ4/GPIdSDob2NvMzmqmjZOwZTFwhZkdk8b7O9OHZtSgO44zSczsTkm/kpSxqc1F+TNJ5xD9Dz8GfKipBk6O/YGPp/j+zjTBPRTHcRynKXgOxXEcx2kKLiiO4zhOU3BBcRzHcZqCC4rjOI7TFFxQHMdxnKbw/wO37NTRu6sDhgAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "plt.plot(ss_rdf_nodensity.bins, \n", " ss_rdf_nodensity.rdf[0][1][570])\n", "plt.xlabel('Radius (angstrom)')\n", "plt.ylabel('Radial distribution')\n", "plt.title('RDF between CA61 and {}{}'.format(w570.name, w570.resid));" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## References\n", "\n", "[1] Oliver Beckstein, Elizabeth J. Denning, Juan R. Perilla, and Thomas B. Woolf.\n", "Zipping and Unzipping of Adenylate Kinase: Atomistic Insights into the Ensemble of OpenClosed Transitions.\n", "Journal of Molecular Biology, 394(1):160–176, November 2009.\n", "00107.\n", "URL: https://linkinghub.elsevier.com/retrieve/pii/S0022283609011164, doi:10.1016/j.jmb.2009.09.009.\n", "\n", "[2] Richard J. Gowers, Max Linke, Jonathan Barnoud, Tyler J. E. Reddy, Manuel N. Melo, Sean L. Seyler, Jan Domański, David L. Dotson, Sébastien Buchoux, Ian M. Kenney, and Oliver Beckstein.\n", "MDAnalysis: A Python Package for the Rapid Analysis of Molecular Dynamics Simulations.\n", "Proceedings of the 15th Python in Science Conference, pages 98–105, 2016.\n", "00152.\n", "URL: https://conference.scipy.org/proceedings/scipy2016/oliver_beckstein.html, doi:10.25080/Majora-629e541a-00e.\n", "\n", "[3] Naveen Michaud-Agrawal, Elizabeth J. Denning, Thomas B. Woolf, and Oliver Beckstein.\n", "MDAnalysis: A toolkit for the analysis of molecular dynamics simulations.\n", "Journal of Computational Chemistry, 32(10):2319–2327, July 2011.\n", "00778.\n", "URL: http://doi.wiley.com/10.1002/jcc.21787, doi:10.1002/jcc.21787." ] } ], "metadata": { "celltoolbar": "Raw Cell Format", "kernelspec": { "display_name": "Python (mda-user-guide)", "language": "python", "name": "mda-user-guide" }, "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.8.3" }, "toc": { "base_numbering": 1, "nav_menu": {}, "number_sections": false, "sideBar": true, "skip_h1_title": false, "title_cell": "Table of Contents", "title_sidebar": "Contents", "toc_cell": false, "toc_position": {}, "toc_section_display": true, "toc_window_display": false } }, "nbformat": 4, "nbformat_minor": 2 }