{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Write your own native contacts analysis method\n", "\n", "The `contacts.Contacts` class has been designed to be extensible for your own analysis. Here we demonstrate how to define a new method to use to determine the fraction of native contacts.\n", "\n", "**Last updated:** June 29, 2020 with MDAnalysis 1.0.0\n", "\n", "**Minimum version of MDAnalysis:** 1.0.0\n", "\n", "**Packages required:**\n", " \n", "* MDAnalysis (Michaud-Agrawal *et al.*, 2011, Gowers *et al.*, 2016)\n", "* MDAnalysisTests\n", "* [matplotlib](https://matplotlib.org)\n", "* [pandas](https://pandas.pydata.org)\n", "\n", "**See also**\n", "\n", "* [Fraction of native contacts over a trajectory](contacts_native_fraction.ipynb) (pre-defined metrics and a general introduction to native contacts analysis)\n", "* [Q1 vs Q2 contact analysis](contacts_q1q2.ipynb)\n", "* [Contact analysis: number of contacts within a cutoff](contacts_within_cutoff.ipynb)" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import MDAnalysis as mda\n", "from MDAnalysis.tests.datafiles import PSF, DCD\n", "from MDAnalysis.analysis import contacts\n", "\n", "import numpy as np\n", "import pandas as pd\n", "import matplotlib.pyplot as plt\n", "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Loading files\n", "\n", "The test files we will be working with here feature adenylate kinase (AdK), a phosophotransferase enzyme. (Beckstein *et al.*, 2009) The trajectory ``DCD`` samples a transition from a closed to an open conformation." ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "u = mda.Universe(PSF, DCD)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Background\n", "\n", "Please see the [Fraction of native contacts](contacts_native_fraction.ipynb#Background) for an introduction to general native contacts analysis." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Defining salt bridges" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We define salt bridges as contacts between NH/NZ in ARG/LYS and OE\\*/OD\\* in ASP/GLU. You may not want to use this definition for real work." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "sel_basic = \"(resname ARG LYS) and (name NH* NZ)\"\n", "sel_acidic = \"(resname ASP GLU) and (name OE* OD*)\"\n", "acidic = u.select_atoms(sel_acidic)\n", "basic = u.select_atoms(sel_basic)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Define your own function" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Any function you define *must* have `r` and `r0` as its first and second arguments respectively, even if you don't necessarily use them:\n", "\n", " - `r`: an array of distances between atoms at the current time\n", " - `r0`: an array of distances between atoms in the reference\n", "\n", "You can then define following arguments as keyword arguments.\n", "\n", "In the function below, we calculate the fraction of native contacts that are less than `radius`, but greater than `min_radius`." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "def fraction_contacts_between(r, r0, radius=3.4, min_radius=2.5):\n", " is_in_contact = (r < radius) & (r > min_radius) # array of bools\n", " fraction = is_in_contact.sum()/r.size\n", " return fraction" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we pass `fraction_contacts_between` to the `contacts.Contacts` class. Keyword arguments for our custom function must be in the `kwargs` dictionary. Even though we define a `radius` keyword in my custom analysis function, it is *not* automatically passed from `contacts.Contacts`. We have to make sure that it is in `kwargs`. " ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "ca = contacts.Contacts(u, \n", " select=(sel_acidic, sel_basic),\n", " refgroup=(acidic, basic),\n", " method=fraction_contacts_between,\n", " radius=5.0,\n", " kwargs={'radius': 5.0,\n", " 'min_radius': 2.4}).run()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "One easy way to post-process results is to turn them into a dataframe." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
FrameContacts from first frame
00.01.000000
11.00.988764
22.00.943820
33.00.943820
44.00.943820
\n", "
" ], "text/plain": [ " Frame Contacts from first frame\n", "0 0.0 1.000000\n", "1 1.0 0.988764\n", "2 2.0 0.943820\n", "3 3.0 0.943820\n", "4 4.0 0.943820" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ca_df = pd.DataFrame(ca.timeseries, \n", " columns=['Frame', \n", " 'Contacts from first frame'])\n", "ca_df.head()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Plotting\n", "\n", "We can plot directly from a dataframe (below), or you could use it with other plotting packags such as [seaborn](seaborn.pydata.org/)." ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "scrolled": false }, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYQAAAEGCAYAAABlxeIAAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4yLjEsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy+j8jraAAAgAElEQVR4nOy9eXhc1ZWv/a4qlaTSXCXLgySjEsbBE7ZsJHUCBHAMBDLg0LdzG8IFQkLT5EJu0t03NyRfp8PtdOfjy006QGfgyQQklzSQgQ6dOAMhjJ1BkgfANjbYlmxL8iCrNE8lVe3vjzqnVFUqSVWapVrv8/iRap+zz9m7bJ/fWXvttZYYY1AURVEUx3wPQFEURVkYqCAoiqIogAqCoiiKYqGCoCiKogAqCIqiKIpFxnwPIBWWLVtmfD7ffA9DURRlUbF79+5zxpiSyc5bVILg8/loaGiY72EoiqIsKkTkeDLn6ZKRoiiKAqggKIqiKBYqCIqiKAqwyHwIirKQGR4eprm5mcHBwfkeipKmZGdnU15ejsvlmlJ/FQRFmSGam5vJz8/H5/MhIvM9HCXNMMbQ3t5Oc3MzlZWVU7rGpEtGIvI9ETkrIvvHOS4i8pCIHBGR10RkW9Sxa0XksHXs3qh2r4g8KyJvWT89Uxq9oiwgBgcHKS4uVjFQ5gURobi4eFoWajI+hEeBayc4fh2w1vpzJ/BNa3BO4OvW8Q3ATSKywepzL/CcMWYt8Jz1WVEWPSoGynwy3X9/kwqCMeYlwD/BKTuB75swfwSKRGQVUAscMcYcM8YEgCesc+0+j1m/PwZ8YKoTSIbnD5/lGy8cmc1bKIqiLHpmYpdRGXAy6nOz1TZeO8AKY8wpAOvn8vEuLiJ3ikiDiDS0tbVNaYB/ONrOA8++xeBwcEr9FWWxcPr0aW688UbWrFnDhg0beM973sObb745pWs98MAD9Pf3T6nvvn372LVrV0p9HnroIdavX8/NN988pXsmy8svv8zGjRupqqqipaWFv/iLv0ip/xe/+MVxj/3oRz9i/fr1bN++fbrDnBdmQhAS2ShmgvaUMMZ8yxhTbYypLimZNPI6ITU+L4FgiNeau6bUX1EWA8YYbrjhBq688kqOHj3KwYMH+eIXv8iZM2emdL25FoRvfOMb7Nq1i8cffzymfWRkZEpjGI/HH3+c//k//yf79u2jrKyMH//4x2POmeieEwnCd7/7Xb7xjW/w/PPPJ329hcRMCEIzsDrqcznQOkE7wBlrWQnr59kZGMe4VFeEfdb1TROtfCnK4ub555/H5XJx1113Rdqqqqp45zvfiTGGT33qU2zatImLLrqIJ598EoAXXniBK6+8kr/4i79g3bp13HzzzRhjeOihh2htbWX79u2Rt92PfexjVFdXs3HjRj7/+c9H7lFfX88ll1zCli1bqK2tpauri3/4h3/gySefpKqqiieffJIXX3yRqqoqqqqq2Lp1Kz09PTFjv+uuuzh27BjXX389X/3qV7nvvvu48847ueaaa7j11ls5fvw4O3bsYPPmzezYsYMTJ04A8OEPf5iPfexjbN++nfPPP58XX3yRj3zkI6xfv54Pf/jDY76j73znOzz11FP84z/+IzfffDNNTU1s2rQJgEcffZQPfvCDvP/97+eaa67h1KlTXH755VRVVbFp0yZefvll7r33XgYGBqiqqhpjyfzjP/4jr7zyCnfddRef+tSnxlyvt7eXHTt2sG3bNi666CJ+9rOfAdDU1MS6deu444472LRpEzfffDO//e1vufTSS1m7di11dXUA9PX18ZGPfISamhq2bt0a6T+jGGMm/QP4gP3jHHsv8EvCFsHbgTqrPQM4BlQCmcCrwEbr2P8B7rV+vxf4UjLjuPjii81UueorL5hbv/unKfdXlMk4ePBg5Pf7ntlv/uvDv5/RP/c9s3/C+z/44IPmk5/8ZMJjP/7xj81VV11lRkZGzOnTp83q1atNa2uref75501BQYE5efKkCQaD5u1vf7t5+eWXjTHGVFRUmLa2tsg12tvbjTHGjIyMmCuuuMK8+uqrZmhoyFRWVpq6ujpjjDFdXV1meHjYPPLII+buu++O9H3f+95nXnnlFWOMMT09PWZ4eHjMGKPv9/nPf95s27bN9Pf3R/o/+uijxhhjvvvd75qdO3caY4y57bbbzF/+5V+aUChk/v3f/93k5+eb1157zQSDQbNt2zazd+/eMfe57bbbzI9+9CNjjDGNjY1m48aNxhhjHnnkEVNWVhaZ55e//GXzT//0T5E5d3d3G2OMyc3NHffv4IorrjD19fUJrzc8PGy6urqMMca0tbWZNWvWmFAoZBobG43T6YwZ9+233x6Zkz3Xz3zmM+YHP/iBMcaYjo4Os3btWtPb2ztmDNH/Dm2ABpPEMzaZbaf/BvwBuFBEmkXkoyJyl4jYryG7rAf/EeDbwH+3hGYEuAf4NfAG8JQx5oDV537gahF5C7ja+jyr1FR62XO8g2BIa0gr6ccrr7zCTTfdhNPpZMWKFVxxxRXU19cDUFtbS3l5OQ6Hg6qqKpqamhJe46mnnmLbtm1s3bqVAwcOcPDgQQ4fPsyqVauoqakBoKCggIyMseFNl156KX/7t3/LQw89RGdnZ8Jz4rn++utxu90A/OEPf+BDH/oQALfccguvvPJK5Lz3v//9iAgXXXQRK1as4KKLLsLhcLBx48Zx5zIeV199NV6vF4CamhoeeeQR7rvvPl5//XXy8/NTulb89YwxfPazn2Xz5s1cddVVtLS0RJbzKisrY8a9Y8eOyJzsOfzmN7/h/vvvp6qqiiuvvJLBwcGIpTRTTPq3Yoy5aZLjBrh7nGO7CAtGfHs7sCPJMc4ItT4vP/zTCd441c2mssK5vLWShnz+/Rvn/J4bN25MuB4O2NZ8QrKysiK/O53OhOvdjY2NfPnLX6a+vh6Px8OHP/xhBgcHMcYktdXx3nvv5b3vfS+7du3i7W9/O7/97W9Zt27dhH1yc3PHPRZ9T3v8DocjZi4OhyPltfvoe15++eW89NJL/OIXv+CWW27hU5/6FLfeeuuUr/f444/T1tbG7t27cblc+Hy+SMxA/Lij52TPwRjDT37yEy688MKUxpAKaZPLqKYyrNIN6kdQlijvete7GBoa4tvf/nakrb6+nhdffJHLL7+cJ598kmAwSFtbGy+99BK1tbUTXi8/Pz+y1t/d3U1ubi6FhYWcOXOGX/7ylwCsW7eO1tbWiLXR09PDyMhITF+Ao0ePctFFF/HpT3+a6upqDh06lNLcLrnkEp544gkg/GC97LLLUuo/FY4fP87y5cv5q7/6Kz760Y+yZ88eAFwuF8PDwylfr6uri+XLl+NyuXj++ec5fjypjNQR3v3ud/Ov//qvEXHfu3dvymOYjLQRhLIiN2VFbuqbOuZ7KIoyK4gITz/9NM8++yxr1qxh48aN3HfffZSWlnLDDTewefNmtmzZwrve9S6+9KUvsXLlygmvd+edd3Ldddexfft2tmzZwtatW9m4cSMf+chHuPTSSwHIzMzkySef5OMf/zhbtmzh6quvZnBwkO3bt3Pw4MGIU/mBBx5g06ZNbNmyBbfbzXXXXZfS3B566CEeeeQRNm/ezA9+8AMefPDBKX9PyfLCCy9EnOA/+clP+MQnPgGEv5fNmzenvD325ptvpqGhgerqah5//PFJLaR4Pve5zzE8PMzmzZvZtGkTn/vc51LqnwwykSm50KiurjbTKZDzySf28p9H26n77A6NKFVmnDfeeIP169fP9zCUNCfRv0MR2W2MqZ6sb9pYCBBeNmrrGeJ4+9T2ViuKoixl0koQan1hP0Kd+hEURVHGkFaCcMHyPDw5LuobVRCU2WExLcEqS4/p/vtLK0EQEap9Xo1YniMGh4McOt0938OYM7Kzs2lvb1dRUOYFY9VDyM7OnvI10q5Aztbzinj24Bl6BofJz55aVSElOZ5qOMkXfn6Qhr+/mkL30v+uy8vLaW5uZqpJGBVlutgV06ZK2gnCstxwwEfXgArCbHOme5DhoKGlYyAtBMHlck25UpWiLATSaskIID87rIE9g4sj++Bixv6OWzoH5nkkiqIkQxoKQvhNVQVh9rG/41YVBEVZFKSdIBS4wxZC90DqoedKatjfsVoIirI4SDtBiFgIQyoIs40uGSnK4iINBUF9CHNF96BlIXSoICjKYkAFQZk11IegKIuLtBOErAwnWRkO9SHMAbaFcLZniKGR4DyPRlGUyUg7QYCwH6FbLYRZJRQy9A6NsKowHDV5umtwnkekKMpkpKUgFGRn0DOoFsJs0hcYwRhYtzJcdlD9CIqy8ElLQcjPzlAfwixjf78XriwAdKeRoiwG0i51BdhLRlOzEAaHgzFO0uUF2eRljf0ax8uV1DM4TF5WRtIFeoIhQ2AkhDvTOaXxTnTdE/7+SCK2QreL4rysCfv0DY3gdjlxOGLHbvsHsjJGx2h/vxeuzAMmFoT+wAjZGYmvKwiZGWn53qIoc05S/9NE5FoROSwiR0Tk3gTHPSLytIi8JiJ1IrLJar9QRPZF/ekWkU9ax+4TkZaoY++Z2amNT4F76hbCPT/cw7u+8mLkz3/5xu/HnPPKW+e4+Au/5URcIZ7O/gBv/+Jz/KihOen7ff8PTVzxf54nFJrZDJpfffZNtn/5hcg83nH/72jvHRr3/JFgiMu/9Dzf+8/GMcf+5sl93P34npg2+/stzs1ieX7WuDuNBoeDXHL/7/jR7pNjjv2Pf9vL3zy1L5VpKYoyDSYVBBFxAl8HrgM2ADeJyIa40z4L7DPGbAZuBR4EMMYcNsZUGWOqgIuBfuDpqH5ftY8bY3ZNfzrJkZ/lmrIP4dDpHmp9Xh68sYprN67kSFsvI8FQzDkHWrsIBEP859FzMe31TR30BYK88ObZpO935GwvZ3uGODfBw3oqvPxWG+tW5vPgjVX8r2svJDASmjAt+JmeIdr7Ajx/OHbsI8EQLx5u41hbX0y7/f3mZ2dQWuQe10Jo6Rygs3+YA61j02Tvb+nmYIJ2RVFmh2QshFrgiDHmmDEmADwB7Iw7ZwPwHIAx5hDgE5EVcefsAI4aY45Pc8zTZqo+hGDIcLprkGqfh51VZbzzbcsIhgxne2If1vbbcHwhHvuBW9fYkXTO/I7+ADCza/B9QyPsb+1mx/rl7Kwq447Lzifb5aCusWPcPvac9hzvZDhKAA+e6qYvEMRvjdPG/n7zs12Uedy0dibeZWRfN96CCIYMp7sHaekc0PoCijJHJCMIZUC0Pd9stUXzKvDnACJSC1QA8Um5bwT+La7tHmuZ6Xsi4kl61NMkP9tFfyAY82BLhrM9g4yEDGUeNwBlReGf8Q9r+3N8qc46SyDO9Q7RlGRd5/bemReEvSc6CYYMNVZJ0cwMB1Wriya0EOxdQgPDwZi3eXtOnf3DMZaSHedRkJ1BmWUhJHqw29dtjtuFdKZ7MOI/OdcbGNNPUZSZJxlBSOT9jP+ffT/gEZF9wMeBvUDkFVxEMoHrgR9F9fkmsAaoAk4BX0l4c5E7RaRBRBpmqvCIneCuN0UrwX6LLS2KFYT4t9sW6224uWOAU13hY/2BEfa3dHH1hrDhlGwZT9tCmMlo37omPw6BiytGNbjW5+VAaxe9Q4m/k2hBih57tIh0RgX72XEeBW4XZUXucR/s41kI0Z810llR5oZkBKEZWB31uRxojT7BGNNtjLnd8hXcCpQA0d7H64A9xpgzUX3OGGOCxpgQ8G3CS1NjMMZ8yxhTbYypLikpSWpSkzHVFNj2W2y5JQS2MMS/3bZ09LP1vCJg9A1634lORkKGm2pX48lxjbEexsPfN/P5gOob/axfVRCzC6qm0kvIwJ7jiZeNWjoH8OS4qCjOiYzdGENDU0dkl1VH3+gDv2dwBJdTyMpwRL6nRA/2Zqute3Akxq8TLUC6ZVVR5oZkBKEeWCsildab/o3AM9EniEiRdQzgDuAlY0y0N/Am4paLRGRV1McbgP2pDn6q2PmMUt16aq+D2w+43KwMinJcMQ+6nsFhugdHuGr9CnIznaN+gyY/InBxhTfpus6hkInyIcxMpG9gJMTekx2R5SKbbed5cDpk3HG1dg5Q5nFT4/PS0OQnFDIcbeujvS/A9nXLAWiPEYTwtlsRGXdpzb7u6O+jc2xRC0FR5pxJBcEYMwLcA/waeAN4yhhzQETuEpG7rNPWAwdE5BBha+ATdn8RyQGuBn4ad+kvicjrIvIasB34m2nPJkmmKggtnf0U5bjIjYo7KCtyJ3yorfbmsK3CQ73lqK1v8rNuZQGFbhe1Pi/H2/s52z3xQ75ncISgtd10ph6K+1u7GBwOUVsZKwi5WRlsLC2IWDTxtHYOUFroptbnpaN/mKNtvRHxePfG8DJYtIXQPTgS+Z7HW1oLtw2y2jv2eEtH2CLJzXSOscAURZkdkgpMs7aE7oprezjq9z8Aa8fp2w8UJ2i/JaWRziAFU1wyau0cpLTQHdNWWuTmeHtf1Dnhh1dZUfjh+ZVn3+Rc7xB7jnfyX6vDfvYa62Fc1+TnfZtLx72fvXMnPztjxpZN7PX/eAvBbvu/fzzO0EgwJsjMmHBd5EsvWBYz9t1NHSzLy4z4IqJ3GvUMDke+5wJ3BnlZGWMe7KGQ4VTXAO/fUspJf0tk+QhGLZLASEgtBEWZI9IyBHSqgtDSMRDZYWRTVuSmpWN0B01zlCDYD8/v/76JgeFg5PPG0gLcLuekjmV/X3g766bSQroGhsd1+KZCfZOfymW5lOSPjUqu8XkYGgmxv6Urpr17YIS+QJCyIje+4hyW5WVR3+inrslPdYUXb254tdDfG+tDsC0EEaG0KHvMg72td4jhoKFqdREup8RaCJZFMlEMg6IoM0taCsJoTYRUfQgDkeUPm7IiN32BIN0Do7n/XU5heX5W5EH36O+bgPBOHgCX08G2iiLqmsbf9w+jDuVNZQWRa0+HUMhQ39RBjS/xDt9qa3zx8QjNneEtsmVFbkSE2koPvzt0luaOAWoqvWRlOMnLyhhjIdjfs903/sFuWwyrPTmsKnRHHOe2RVLmcY9ZklMUZfZIS0HIs30IA8m/cXcNDNMzNDJWEDyxDtOWjgFWFmbjcAjZLieby4voHhyhojiH5QXZkX41Pi+HTnfTNUFdBntNflNZYcw9pspbZ3vpGhhOuFwEsCwvi/NLcsc4luOd6TU+b2RbqS1ynlzXmF1G0buYShM82KO38UZbENEWSWmRm47+YfoDmoxQUWabtBQEl9OB2+VMyUKIj0GwKY3bQRNvRdgP3/iHcK3Pi5lgmyeM7tq5yBaEaTpX7e2i8Q7l+HHZu4hsWjosC8EzKggAuZlO1q8Kp7f25mbh74+KQxgY9SHYfeMf7C2R7zSbsqKcUVGN+q7LPeM7pBVFmVnSMtsppJ7gzn4YlxZlx7TH76Bp6RzgHWtGfei1lR4efpExyzRbz/OQ4RAajvsj2zbj6egPkO1yUFGcS4Yjdo298VwfN3zjP+kfSr4S2UgoxPL8LM7z5ox7To3PyxP1Jzl8pof1q6ylqq5BsjIcFFu+gvWrCsjPyqDqvCIynOF3Cm+OizYr31IwZOgLBMcsGQGc9A9woVUjobVzgILsjHB6i6JsznQPMhwMRQShrMhNwIp+bu4Y4ILl+UnPVVGU1ElbQcjPdtEzlIKFYEUcxzuVi3Mzycxw0NI5wHAwxJnuwUjgGsDla0v4ws6N7KyKzfbhznSyoiCbUxPEF/j7AnhzMnE6hJWF2TFLRr87dJbO/mHuuKwSVwrpoWsrvROm3rath/omf0QQWjoGIv4DAKdD+NcPbWVF1BKYJzeTN8/0AqMR4NGCsMG61qsnOyOCEPYThMWpzOMmZMKV1aKtMTu9yHi5kBRFmTnSWBAyUvIhtHQMkOl0sCw3dneOwyERh+nprkFCJnZZKcPp4JZ3+BJe05ubOSYpXDT+vgDevPBbebxztb7RT7nHzd+/Lz7x7PQo97hZWZBNXaOfW61xt3QOjFkqu/LCWKumODcTv7XEZcd3RC8ZXbA8LxKh/V9rVkeuay8JRUczt3QOkJnhYFleJsGQwRlnHSmKMjukpQ8BLAshBR9C+KGYPaaICxBxiEZiEOKsiPHw5GbGOGLj8fcF8OREC0L4LdkYQ32TP+LQnUlEhJrKcCS1vZU20e6qeDy5mQwMBxkIBCOCEG0hiMiYCO3o60aW3rrCgmBbJBlOBysLsnXrqaLMAWkrCAUppsBO9JZsY8citIzjeB6P4kkshI7+QGSPf5nHzenuQUaCIY6dC6eMqJnAOTwdan0eznQPcdI/wNBIkLM9Q5POyWsJl78/EPleC9yxFeOiI7TtFB/2dSPO+Y6ByBKVTaItq4qizDxpKwjhMprJC8JEb8mlRW7O9gzRdC4csTzZ27SNJyczJpgrHn9vlCAUuQmGDGd6hiaMNp4JoqORT3eFrZLJrB57nB19gahaCLErktHXta0d+7rZLifL8jJp6Qz7EKKd92Ue94wm91MUJTFp60MoyM5IOpdRYCQ04VuyLQC7T3RQnJtJtiu5+sfeXBd9gSCDw8ExfQIjIXqGRiJv3tFv0HVNfopzM1lTkpvUfVLlbcvzKXS7qG/0U1qYbd0/e8I+tiC09wWiqqXFWgjREdo5Vo3o6O+0tMhN07k+zvYMUVaUE9WezWmrPoIzwZKdoigzQxpbCBkERkKRAvETcbprEGPGf0u2BWHP8c6k/QcQ3rsP4eIy8XRaS0meqCUjCFsq9U1+qn2eCXcLTQeHQ6iu8FDf5I8s1ZQXjb9VNXqcHX2BSHGceAshOkLbzt5aHrc0tO9kJxArQGVFOWHraJJkgIqiTI+0FQR7fTsZP0L0vvhE2A/rgeHgmOR3E+HNDY+hvW9svWTbt2Dv/bevu/t4Byf9A7O2XGRTU+nl2Lk+XmvuQgRWFk5iIdg+hAmWjGA0QvvQqe7wrq280V1bpUVuBobDAl0WYzmE7607jRRldklbQRjNZzR9QVhZmI39sp6KhWDvIOroG2sh2L4F+83bnemkODeTX7x+Cpg42ngmsAXnF6+fYnl+FpmTxDoUul04JOwI7xkaISvDEZMx1caO0P7l/tOsitu1FeNIjvoey+PSgyiKMjukryBkhd/OuyfIJWRjv5mO95acleGkxHrTTXaHEYyuuyfaaWS32efY1/b3BcjNdEYCvWaLi8oKyXY58PcFkpqTwyF4cjIjPoR4/4GNHaHt7wskTCVuE/1dx6cHURRldkhfQUjFQugYoCQ/a0Jnsf1Gm+wOI4gShN6xS0Z2fIJtRURfe1uFJ5IyYrbIzHBQtboo5r6TYcdVdA+MUJBguQjClo6drC/emrItgeX5WTHWRU5mBp4cl+40UpRZJo0FwfYhJGEhdI0fg2BjH09FEArdLkSISQpn0x4RhNiMoTB7203jsQPfkp2T14pW7h4cJt+d2EKA0eWu8RIFJvquE2VLVRRlZklbQShwp2YhlE/yULSPp+JDyHA6KHS7EkYrd/QFKHS7YiyB+Gyjs40dN5DsnLw5mRGn8ngWAoyOP/479eS4cLucCe9XVuTmj8f8XP+1V7j+a69wzw/3xGRkhfDOrL/+QQOnumZPOAaHg9zzwz0cPt0za/dQlPkibQXBthAmi0XoHhymsb2PtSvyJjxvZ1UZd29fE/NGnwzeqBxA0fj7h2P8BwDXbFjBhy/xRUpWzja1lV5ue0cFV61fkdT5ntzMsFM5rjhOPJddsIxb3l7BletKYtpFhL+75m18qPa8MX1uqj2Pt5/vpTg3E2Pg56+d4khbb8w5vzt0ll8fOMMvXjuV1Hinwu7jHfz8tVP8ZE/zrN1DUeaLtA1My8uyiuRMYiHsPt6BMUyaN2hDaQEbSlN39Npv1fH4+4bGCMJqbw73Xb8x5XtMlawMJ/9756akzy/OzaSjfxiQiNM+Ee5MJ1/4QOLr3vHO8xO2b1+3PJImvPFcH9u//AJ1jX7etmI0JbadJ6m+yT/udaZLnRUlXjdJ+VNFWYykrYXgdAh5WRmT+hDqG/1kOISq84pmZRz2W3U8/r7hGIfyYsCTG85Oeq53aEILYbrYdZ0b4iq72Q/phqaOSGK+mcYWnf0tXVrFTVlyJCUIInKtiBwWkSMicm+C4x4ReVpEXhOROhHZFHWsSUReF5F9ItIQ1e4VkWdF5C3r59ysg0SRTIK7+iY/G8sKycmcnQdccW5mxIEcTUdfIBK4tliIHm98YruZxK7rXB9Vk7q9d4ijbX2cvyyX9r4AR9v6Zvy+w8EQe090cv6yXEZChn0nOmf8Hooyn0wqCCLiBL4OXAdsAG4Skfgk/J8F9hljNgO3Ag/GHd9ujKkyxlRHtd0LPGeMWQs8Z32eU/KzXRPGIQwOB3n1ZBe14xSlnwnsrZrRb7TGmHAthLjaCwud6PHOpoUAYcd0S+dohllbHO66co31eeaXdPa3dDEwHOTOy89HZLQkqaIsFZKxEGqBI8aYY8aYAPAEsDPunA2EH+oYYw4BPhGZzBO5E3jM+v0x4ANJj3qGyJ/EQnituYtAMDSru3q8OZmMhAw9Q6Pj6AsECQRDi89CiFriGi8wbaaw/07szK/1TX4yMxzsrCplWV5mpH0msUXmXeuXs25lwayIjqLMJ8kIQhlwMupzs9UWzavAnwOISC1QAZRbxwzwGxHZLSJ3RvVZYYw5BWD9TFhYWETuFJEGEWloa2tLYrjJk5+dMWEZTfs//KwKQiQ4bXTZKFFQ2mLAEyVgs20h2HWd66IcyVWri8jKcFJd4Z2Vt/f6pg58xTksz8+m1udhz/HOSIlPRVkKJCMIiVJqxnvs7gc8IrIP+DiwF7BfeS81xmwjvOR0t4hcnsoAjTHfMsZUG2OqS0pKJu+QAgVu14QWQl2jn7XL8yL5hGaDROkrbJ9Ccd7iEoTiqCWjglm2EJwOYVuFh/pGP31DIxxo7Y7sBKup9NLcMTCj8QihkKGhyR95Oaip9DIwHORAa/eM3UNR5ptkBKEZWB31uRxojT7BGNNtjLndGFNF2IdQAjRax1qtn2eBpwkvQQGcEZFVANbPs9OYx5QI11VObCEEQ4Y9xztmrSqZTXTaaJvFaiG4M51ku8L/pGbbQppfkZ0AACAASURBVIBwnMRbZ3t57tBZgiET+buyhWEmt4Yebeulo394zD1mY2lKUeaLZAShHlgrIpUikgncCDwTfYKIFFnHAO4AXjLGdItIrojkW+fkAtcA+63zngFus36/DfjZ9KaSOuG6yiMJtyi+caqbnqGRWalbHE1xVGEZGzsuIT4OYTFg+xFm20KA0aW8b75wFIfANmtr8PpV+eRmOmd0jd9egrL/PSwvyKaiOEcdy8qSYtLXOGPMiIjcA/wacALfM8YcEJG7rOMPA+uB74tIEDgIfNTqvgJ42irkkgH80BjzK+vY/cBTIvJR4ATwwZmbVnLkZ2cwEjIMDodwZ8Ymrov4D+bBQljUgpCXSWvX4JxYCJvLC8l0OnjjVDebygoijuwMp8NaTuqY5ArJU9/opyQ/i4ri0UJBNT4vz71xhlDIxKTxVpTFSlL/a40xu4BdcW0PR/3+B2Btgn7HgC3jXLMd2JHKYGeagqgEd4kEoazInVKyuqmQm+kk0+mI8SH4+wO4nBKJpl5M2MtccyEI2S4nW1YXUt/UMcbxX+vz8pVn36SzP0DRDCy91Td1UOvzxlSpq/V5+fHuZo629bI2KmJaURYri++JM4PYD63PP3NgzMP35bfOsWNdwo1PM4qI4LViEWw6+gJ4cjJnrUTmbOLNzSQn0znr6bltanzeyMM6pt2y7BqaOrhqQ+wO6F/tP0VuVgbvXBu7SeFgazeP/b6JUNwS4kjI0NI5wF+9szLhPT73s/2s9oQth51VZVy2dtn0J6Yo80BaC8KmskLOL8nl1ZNjI04L3S52VsXvrp0dPHEJ7o629aaUNXUh8a51y8lOUClttnj/llLqGv1csib2IVy1ugiXU6hv8scIgjGGv//3A5TkZ/HLT8QKwndePsYzr7ayPH9sQOCaklx2xCX58xXncMXbSnjrTA8n2vvx9wc40targqAsWtJaENaU5PG7v7tyvoeBN9cVEQQ7OvrDl/rmd1BTZGdV2ZwJKYTjEX78sUvGtGe7nGwuLxrj9G1q7+dc7xDtfUN0DQxTGJVio67Jz1XrV/DwLRcndW8R4bGP1EY+/3+/OsS3XzrGQCA4ZglSURYDaZvcbiHhybGzhM5NdHS6UOPz8npzFwOBYKTN3iZqDOw5Pup0PtU1QHPHwLQ2EdT6vIyEDHtPzpwzW1HmEhWEBUBxbibtVhlNe3dT9RzVPFjK1Pg8Yx7QdU3+cOEhh8RYD3bMwnS2GW+r8CDCjO5uUpS5RAVhAeDJzaR7cIThYMjK8T+70dHpQnWFF5GwY9mmvsnPn1V62VRWGBNU1tDUQW6mk/Wrpr5bqNDt4sIV+ZrjSFm0qCAsACLpK/oC4ehoXS6aEQpzYh/QZ7sHOd7eT22ll9pKL681dzE4HF5Oqm/ys63CM+3dUbWVXvac6GBEcxwpixAVhAWALQi/P3ouHB09y8Fw6USNz8ue4+EHdF1UssIan5dAMMSrJzvp6h/m8JmeGYlKr/F56Q8EOXhKcxwpiw8VhAWAne7hV/tPA7ObXTXdqKn00mc9oOsb/bhdTjaUFkR8NPVNfhqO+zFmZqLSbTHXEpvKYkQFYQFg+wtefLONsiI3pbMcHZ1ORCe6q2vqYFtFES6nA09uJm9bkUddUwd1TX5cTqFq9fTLpK4oyOY8b476EZRFiQrCAsBOcDc4HNLlohlmZWE2q71ufnfoLIdOd8dYX/Zy0h+P+dlcXkS2a2ZiB2p83lmt66wos4UKwgIgOteOLhfNPDU+L78/2o4xsdtKayu99A6N8OrJzhn93msrPbNW11lRZhMVhAVAZoaDfCuXUm2lxh/MNLYIZDiEreeNfr81MeIwc997pLynLhspi4y0Tl2xkPDmZeLKcLCmJG++h7LksJ3Fm8oKY1JKlFrZbFu7Bri4YuYshMpluSzLy+Qnu5vpt6KkL1iexxVvi82dNDQS5Me7mxkcHrtF1SHw3s2rWJ6fPWPjWow898YZ3rGmmJxMfVTNBfotLxA2lRZS4M5YlBlOFzrnL8tl7fI8rtm4Ysyxd29cyaHT3TE5jaaLiLBj3QqebDhJg5UeI9Pp4NXPXxMjSL8+cIb/5+n9412Gk/4B/uH9G2ZsXIuNw6d7+OhjDfzzDZu4+c8q5ns4aYEKwgLh6zdvUyfkLCEiPPu3VyQ8NlsP3Pv/y0V89r3rAXjlrXPc/cM97DvZyTvWFEfOqWtsJy8rg5f/1/YxBXb+6vsN1DW1z8rYFgt/agzP/4S/f55Hkj6oD2EBodbB0kFEKHS7KHS7uOyCZeEcR3E+hfrGDrZVePDkZkbOtf+8vdLLwdZuegYT1/xOB+xYjtbOwXkeSfqggqAos0x8Cg2Azv6AFR2d2JldU+klZGDPibG1OtIBY0zk+2rpUAthrlBBUJQ5IDqFBowm3Btvu+vW8zw4HRKTgC+dOOkf4Ez3EJlOh1oIc4gKgqLMAdU+D32BIG+c6gHCy0cup7BlnOjovKwMNqwqGFPgJ12w571j/XLO9AwSGNFkgXNBUoIgIteKyGEROSIi9yY47hGRp0XkNRGpE5FNVvtqEXleRN4QkQMi8omoPveJSIuI7LP+vGfmpqUoC4tIjiPrQVfXNHl0dI3Py76TnQyNBMc9Z6lS3xiuW3HF20owBs50q5UwF0wqCCLiBL4OXAdsAG4SkfitGZ8F9hljNgO3Ag9a7SPA3xlj1gNvB+6O6/tVY0yV9WfXNOeiKAuWVYVuyj1u6hv9DASCvN7cNWl0dG2lh8BIiP0tXXM0yoVDfZOfGp+Hck8OAC2dA/M8ovQgGQuhFjhijDlmjAkATwA7487ZADwHYIw5BPhEZIUx5pQxZo/V3gO8AcxdwV1FWUDU+rzUN/nZe6KDkZCZNDq6OpKYL70qsLX1DHHsXB81Pi9lnnCix5YOFYS5IBlBKANORn1uZuxD/VXgzwFEpBaoAMqjTxARH7AV+FNU8z3WMtP3RCTh/w4RuVNEGkSkoa2tLYnhKsrCpKbSS3tfgKcaTiLCpNHRy/KyOL8kN+1SYDTYdSsqvawqDEdqt6qFMCckIwiJNsfHR1DdD3hEZB/wcWAv4eWi8AVE8oCfAJ80xtiVQ74JrAGqgFPAVxLd3BjzLWNMtTGmuqSkJNEpirIosJeInnm1lQtX5CcVHV3r89LQ5CcUSp+gxbomP9kuB5tKC8l2OVmWl6VLRnNEMoLQDKyO+lwOtEafYIzpNsbcboypIuxDKAEaAUTERVgMHjfG/DSqzxljTNAYEwK+TXhpSlGWLGtKcinOzSRkSDrNeY3PS/fgCIfP9Mzy6BYO9U1+qlYXkZkRfjyVFWWrIMwRyaSuqAfWikgl0ALcCHwo+gQRKQL6LR/DHcBLxphuCYfefhd4wxjzL3F9VhljTlkfbwDGT+qiKEsAEaHa5+HXB84knW7bFo6f7mlm+4XLk75XfraLi8oLpzTO+aRncJiDrd3cs/2CSFuZx82h0+kjiPPJpIJgjBkRkXuAXwNO4HvGmAMicpd1/GFgPfB9EQkCB4GPWt0vBW4BXreWkwA+a+0o+pKIVBFefmoC/nrmpqUoC5PL1pbw/KE2/ixJC6Hc42a11823X27k2y83pnSvn3/8MjaVLS5RONjaTcjAtopRl2JpYbjAkTFG07vMMkklt7Me4Lvi2h6O+v0PwNoE/V4hsQ8CY8wtKY1UUZYAH6o9jyvfVsLyguTSWosIP/rrSzjennyxnd6hET76WAN/PNa+6AThXG8ACG/TtSnzuBkcDuHvC1CclzVfQ0sLNNuposwhToew2puTUp+VhdmsLEytLkJFcQ51jX7ueOf5KfWbb/x9QwB4c0erCNo1xls7B1UQZhlNXaEoS5Aan5eG44uvrrO/L5zdtShndAdWmSUILZ2a5G62UUFQlCVIjc+Dvy/A0bbe+R5KSnT0ByjIzsDlHH00jQqCpq+YbVQQFGUJUrNIo5zbE/gJinJcuF1OjVaeA1QQFGUJYtd1XmxRzh19ATw5sQF7IkKZx63RynOACoKiLEFEhBqfN1J1bLHg7wvEOJRtSovcGpw2B6ggKMoSpcbnpaVzYFG9WY8nCGVFaiHMBSoIirJEsaOcF8uykTEGf38AT0JByKa9L8BAIP1qQ8wlKgiKskRZv6qAvKyMRSMI/YEggZEQ3pwEgmClwW7tUithNlFBUJQlitMhbKvwUL9Idhr5+8JRygl9CIVaF2EuUEFQlCVMrc/D4TM9dPYHku5jjKG9dyjyZ3A4tWUaY8yUAuImEgTbQnjrbG9kXOmUEnyu0NQVirKEsauu7T3RyfZ1yWVL/edfvMF3XhlNpFeSn8UfP7MDp2M0LdlJfz9Xf/VFfvqxS9lQWhDT/90PvMQHtpbx36+8gFTwW6KVyIewoiAbl1P4ws8P8oWfHwTg5j87j3++4aKU7qFMjAqCoixhLlieB5BScrzfHTrLRWWFfLC6nL0nOnl6bwtnugcjOYUADrR2MTgc4s0zPTGCMBwM8eaZXn5z4EzqgmAltitOIAgup4PvfbiGxnPheTy9t4XnD51N6frK5OiSkaIsYYpzM8nKcNDalVzaB7ue8fs2r+LWd/jYWVUKjC1haaeRsJd5bDr7w7mI9rd00R8YIRU6JrAQAN65toRb3+ELj2tLKa1dgzR3aH6jmUQFQVGWMCJCWZE7aWdsdD1jCNdjAMYEhdnX64jzTdifR0KGfSc6Uxqrvy9AhkPIz5p84aJmkW2pXSyoICjKEieVKN/oesZ2XxgrCLbF0B5nIbT3jn6uS/Fh3WHFICRTBGfdygLyszIWXa6mhY4KgqIscUpTqElc3+Rn62pPpJ5xTmYGnhzXGAvDvl5HX2ILwe1ypvz23t4bSOg/SITTIVzs86iFMMOoICjKEqesKIe2niGGRibePmrXM66JK+9ZmiBthP053odgf77ywhL2HO9kOBhKepwd/QE8CYLSxqPG5+XI2d4xY1CmjgqCoixxSovC1dZOTVJPYM+JTkIGan2xglAWt+Q0EAhGloriH8a2xXDNxhUMDAc50Nqd9DjHy2M0HostNcdiQAVBUZY4kbQPkywb1Tf6cTqErecVxbSXWk5pO9jMTh+Rm+kc41Ru7wuQn5XBpWuWRa6ZLKkKwubyQjIzHCndQ5mYpARBRK4VkcMickRE7k1w3CMiT4vIayJSJyKbJusrIl4ReVZE3rJ+emZmSoqiRGNXHGueRBDqmvxsKi0gN26XT7nHTV8gSPdAeBup7U/YWFpIR/9wTMSw7RheXpAdruuc5Nt7MGToHBged8tpIrIynFSVF6mFMINMKggi4gS+DlwHbABuEpENcad9FthnjNkM3Ao8mETfe4HnjDFrgeesz4qizDArC7MRmdhCGBoJsu9kZ6TSWjTxO43s62wqKyQYMnQPDkfOjX7Lr/F5aWjyJ5ViomtgGGPAG1ccZzJqKj3sb+2mbyi1mAclMclYCLXAEWPMMWNMAHgC2Bl3zgbCD3WMMYcAn4ismKTvTuAx6/fHgA9MayaKoiQkK8NJSV7WhLEIrzd3ERgJRVJdRBMvCC2dAzgE1q/KB2L9CNGCUOvz0tE/nFRdZ3/fEADeuPKZk1Hj8xIMGfamGPOgJCYZQSgDTkZ9brbaonkV+HMAEakFKoDySfquMMacArB+JpdoRVGUlCnzuCdMHW0v7dT4xq7c2ktOrVGCsLIgm5L88MM72o8QLoEZFoRq61rJLBv5+8JWRqLU1xNxcYUHh6Qe86AkJplcRomiROJtwPuBB0VkH/A6sBcYSbLvxDcXuRO4E+C8885LpauiKBalRW4OtHTFtN31g93sPhEO7OoZHGZNSe6YAvcQTn+RmeEYtRA6BijzuCnODZ8bHYzm7w9QnBd+qNt1nf/p52/wwG/fGnPdtcvzePyOP0NEIlaGJze1JaP8bBfrVxUkdCz/8Vg7//KbN/n+R2vJdjlTum66kowgNAOroz6XA63RJxhjuoHbASQcZtho/cmZoO8ZEVlljDklIquAhJmqjDHfAr4FUF1drfluFWUKlBe5efbgGUIhg8MhnOsd4lcHTlPr87LGSoB39YbERrrDITFbT1u7Bth2nify8LYthIFAkMHhUMRCEBH+9/WbeOXIuTHXbO7o5+W3znHsXB9rSvIigmCLTCrU+Lw8UX+C4WAIl3N00eNn+1qoa/Kz50QHl1i7npSJSUYQ6oG1IlIJtAA3Ah+KPkFEioB+y09wB/CSMaZbRCbq+wxwG2Hr4jbgZzMwH0VRElBa5CYwEuJc3xDL87MjOYs+fd2FXFwx1m8wtn82LR0DBEOGU52DlG12R3wF9nJPu+0HiHrLf+/mVbx386ox1ztytper/uVF6hv9rCnJi4hKUYpOZQjHIzz6+yb2t3Sx9bzRJa86y2qob1RBSJZJfQjGmBHgHuDXwBvAU8aYAyJyl4jcZZ22HjggIocI7yj6xER9rT73A1eLyFvA1dZnRVFmgVE/QDg4ra6xg6wMBxeVFU3ULaZ/a+cAbT1DjIQMpUVucjIzyHY5Ig7hDksYkok2XlOSS3FuZmTt398XIDfTOaWlHXtnVPT20/beIY629Y1pVyYmqXoIxphdwK64toejfv8DsDbZvlZ7O7AjlcEqijI1IjuFOgaoWh3eu1+1uiiSsyiZ/md7hiL1COxgN29OZsRCsAvc2D6EiRARqqNyEfn7AniT6JeIkvwsKpflUtfYwZ2Xh9vqm8K+kY2lBew50cFIMESGU+NwJ0O/IUVJA6KjlXuHRjjQ2hVJ/ZBUf0tQdh/3x3z25GZGlntsSyHZfEQ1Pi8n/QOc7hoMC0KKO4xir+Wh4fhozEN9k5/MDAcfvayS/kBqKTTSGRUERUkDCrIzyMvKoKVzgD3HOwgZEgahjYctAHXWm7dtcXhzM6PyGg1H2pLBFqS6Jn8kwnmq1Pi8dPYPc8SKebAtoEsvWBb5rEyOCoKipAGRQjmdA9Q3+XEIbKtIPluMbWHsbvJT6HaRZ6W38OZmRhLadfQFcDqEguzkHMMbVhWQm+mkvtFPe29qeYziiYhLo5++oREOtHZT6/OyoiCb87w5EQezMjFaU1lR0gR7p1D3wDAbSwsjD/VkWFkYzpjaFwiyflVupN2TMyoI7X0BPDkuHI7JC9wAZDgdbKsI+xE6+qe3ZHSeN4fl+VnUN/mpKM4hGDKRNN41Pi/PHz6LMSap4jvpjFoIipImlHncnPT3j5uzaCKyMpwstyKT7eUjCAet9QyNMDQSjIlSTpYan5dDp3voDwSntWQkItRUeqlv9FPf1BG2gKysrbWVHvx9gaRSaKQ7KgiKkiaUFrmth3eI2srUkwvbfoMyq74CEHmId/YP4+9PfdknWpims2QE4dxJrV2D/MerrWwoLSDfWrqy76HlNidHBUFR0oToN/tESewm7W/5EeyfQFRwWoCOFOsZAFStLsLllJhrTRX7wd94ri9GaOwUGupYnhwVBEVJE2xBOH9ZLstSzCoa3b80SljsJSJ/XwB/X+o7hdyZTjaVFQLTF4QLV+aTnx32i0RXfRMRanxedSwngTqVFSVNsN/sU/UfRPpHloyifAhWMNm53iE6+gMUT+GhXuvzsvdEZ8r+h3icDuHiCg8vHG4bYwFV+7z8cv9pWjsHYgQtnm+9dJS1y/PZvi655MvNHf3c98wBhkbCtaOzMpz8750bY76jxYRaCIqSJqzIz+am2tXc9GdTyxr8rnXL+UBVKetXFUTa7If48fZ+Qib5oLRoPlhdzg1by6gozpnSuKK57RIfd1xWGUnNbVObIL1FPIPDQf7Prw/z+J9OJH2/X7x2it++cZaewRF6Bkf47Rtn2PXaqakNfgGggqAoaYLDIfy/f76ZqtXJ5S+KZ7U3hwdu3BqTb8hORnfkbHgHz1SWfS5Yns9X/7IqJlPpVNl+4XL+/n3xBR3DxXxyM50TCsK+k50MB82ktaejqW/yU7ksl3+/+1L+/e5L8aVQNnQhooKgKMqUcTkdFLpd0xKEuSAS8zDBTiO7pkJLkoIQChnqmzpiigqlUjZ0IaKCoCjKtPDmZnLs3MIWBAgvGx0+00NnVIW3aOw3+66BYXqTqNH81tleugaGY3wyNZXJlw1diKggKIoyLTw5LgaHw07V6QSXzTZ25HJD01grYSQYYs/xjohTPJllI1tAopME2r6KxbpspIKgKMq08EZVOZtO+onZxo55SORHeONUD32BIO+zivkks2xU3+hneX4W53lHneEVxTmU5GclLOm5GFBBUBRlWtgV0twuJ+7MhVu7ONvlZHN5UcK3d7tt59YyIFw3YiKMMdQ3+amp9MbkRxIRan3eSD2GxYYKgqIo08JeJlrI/gObGp+X15u7GAgEY9rrG/2Ue9xsKS8iwyGTLhk1dwxwqmswJgBu9B4eWjoHknZOLyRUEBRFmRbFi0gQais9jIQMe0+OvsHbb/u1Pi9Oh7CyMHvSh7m97JQoyM/2VSzGZSMVBEVRpoUdjLaQHco2F1d4ESFm++mxc3209wUiD3K7fvRE1Df5yc/O4MKV+WOOrVtZQH5WxqJ0LKsgKIoyLWzLwJuTXGGc+aTQ7eLCFfkxjmX7Td5+2y8rck/qQ6hr9FNd4cGZoPaD0yFc7POohaAoSvoREYTc1BPmzQe1lV72nOhgJBjeKlvX5Kc4N5M1JeHCP2UeN6e7ByPH42nvHeJoW1/EokhEjc/LW2d7I8WDFgtJCYKIXCsih0XkiIjcm+B4oYj8h4i8KiIHROR2q/1CEdkX9adbRD5pHbtPRFqijr1nZqemKMpcMCoIC99CgPDDuj8Q5LNPv84//+IgL73ZRrXPE9ktVFrkJmTgdPdgwv72DqKJkgTWJJE7aSEyabZTEXECXweuBpqBehF5xhhzMOq0u4GDxpj3i0gJcFhEHjfGHAaqoq7TAjwd1e+rxpgvz9BcFEWZB1YWZrOlvJCLK6aWRXWuufSCZZQWZvNzKwmdQ4T3bi6NHLczlbZ2DlLuGZtwr77JT2aGg83lhePeY3N5IZkZDuqb/FyzceUMz2D2SCb9dS1wxBhzDEBEngB2AtGCYIB8CUtsHuAH4mO/dwBHjTHHpz1qRVEWDFkZTn52z2XzPYyk8eZm8vvP7Bj3eGlEEBL7ERqa/FSVF5GVMX7MRbbLyZbywkUXj5DMklEZcDLqc7PVFs3XgPVAK/A68AljTPwC3I3Av8W13SMir4nI90QkYU0/EblTRBpEpKGtrS2J4SqKokwd20JItPW0b2iE/a3d1CRRgrTG52V/Sxf9gcnzIi0UkhGEsW70sEUQzbuBfUAp4SWir4lIJGm6iGQC1wM/iurzTWCNdf4p4CuJbm6M+ZYxptoYU11SUpLEcBVFUaaOO9OJNzczoSDsPdFJMGSSKjJUU+llJGTYd6JzNoY5KyQjCM3A6qjP5YQtgWhuB35qwhwBGoF1UcevA/YYY87YDcaYM8aYoGVJfJvw0pSiKMq8M97W07omPw6BiysmtxAurvAgsrgS3SUjCPXAWhGptN70bwSeiTvnBGEfASKyArgQOBZ1/CbilotEZFXUxxuA/akNXVEUZXYoLcpO6EOob/SzobSA/OzJd1QVZLtYv7JgUe00mlQQjDEjwD3Ar4E3gKeMMQdE5C4Rucs67QvAJSLyOvAc8GljzDkAEckhvEPpp3GX/pKIvC4irwHbgb+ZkRkpiqJMk7KiHFo6BzBmdHU8MBJi78mOlGpS11Z62XO8k+FxYhoWGsnsMsIYswvYFdf2cNTvrcA14/TtB4oTtN+S0kgVRVHmiNKibPoDQboGhimyUnPsb+1icDiUMKHdeNT4vDz6+yYOtHZPuXTpXKKRyoqiKHGUe8I7jZqj/Ah2KorqVATB2o20WNJYqCAoiqLEkSgWob7Jz/nLcinJTz5Fx/L8bHzFOYvGsayCoCiKEkd8LEIoZKhvSs1/YFPj89LQ5CcUit+tv/BQQVAURYnDm5tJtssRsRDeOttL18DwhAntxqOm0ktH/zBH23pnepgzTlJOZUVRlHRCRCgtctNwvIOf7G6m4Xg4BUUqDmUbu88jv2/i4vPGxi9cuDKfTWXj50WaS1QQFEVRErB+ZQG/eP0Ue61IY19xDqu97pSvU1GcQ0VxDj/80wl++KcTY46Xe9y88ul3TXu8M4EKgqIoSgK++pdVfPra0YQLxXmZkRTZqSAi7Pof76S9d2xthO/9ZyPf/0MTw8EQLuf8r+CrICiKoiQgM8PBecVj019PhdysDHKzxj5u16/KD9de6BpktXdm7jUd5l+SFEVR0pTJUm3PNSoIiqIo88REqbbnAxUERVGUeUItBEVRFAUIV1Zblpe49sJ8oIKgKIoyj5QWuWnpHJzvYQAqCIqiKPNKuBhP/3wPA1BBUBRFmVdKi9y0dg7G1F6YL1QQFEVR5pGyIjcDw0E6+ofneygqCIqiKPPJQtpppIKgKIoyjyQqxjNfqCAoiqLMI2ohKIqiKAB4cly4Xc4FEYuQlCCIyLUiclhEjojIvQmOF4rIf4jIqyJyQERujzrWJCKvi8g+EWmIaveKyLMi8pb1c2yicEVRlCVOuPZC9uKwEETECXwduA7YANwkIhviTrsbOGiM2QJcCXxFRDKjjm83xlQZY6qj2u4FnjPGrAWesz4riqKkHWWenEVjIdQCR4wxx4wxAeAJYGfcOQbIl3Cy8DzAD4xMct2dwGPW748BH0h61IqiKEuIssViIQBlwMmoz81WWzRfA9YDrcDrwCeMMSHrmAF+IyK7ReTOqD4rjDGnAKyfyxPdXETuFJEGEWloa2tLYriKoiiLi9JCN+d6AwwOB+d1HMkIQqISQfEhde8G9gGlQBXwNREpsI5daozZRnjJ6W4RuTyVARpjvmWMqTbGVJeUlKTSVVEUZVFQ5lkYO42SEYRmYHXU53LClkA0twM/NWGOAI3AOgBjTKv18yzwNOElKIAzIrIKwPp5dqqT3/TEegAACS5JREFUUBRFWcyULpC6CMkIQj2wVkQqLUfxjcAzceecAHYAiMgK4ELgmIjkiki+1Z4LXAPst/o8A9xm/X4b8LPpTERRFGWxUrZAYhEmralsjBkRkXuAXwNO4HvGmAMicpd1/GHgC8CjIvI64SWmTxtjzonI+cDTVmHqDOCHxphfWZe+H3hKRD5KWFA+OMNzUxRFWRSsLMzGIdAyz9HKkwoCgDFmF7Arru3hqN9bCb/9x/c7BmwZ55rtWFaFoihKOuNyOlhRkD3vdRE0UllRFGUBEC6UM791EVQQFEVRFgBlVl2E+UQFQVEUZQFQWuTmVNcAw8HQ5CfPEioIiqIoC4BNZQUMBw0HWrvnbQwqCIqiKAuAWp8XgPpG/7yNQQVBURRlAbC8IJuK4hzqm1QQFEVR0p4an5eG4x0YE58daG5QQVAURVkg1Pq8+PsCHG3rnZf7qyAoiqIsEGoqw36EusaOebm/CoKiKMoCwVecw7K8rHnzI6ggKIqiLBBEhBqfh7p52mmkgqAoirKAqPF5aekcmJfMpyoIiqIoC4hay48wH8tGKgiKoigLiPWrCsjLypiXZSMVBEVRlAWE0yFsq/DMi4WQVD0ERVEUZe6o9Xn48m/auPpfXoy0ffHPL6LGSm8xW6ggKIqiLDBu2FbOkbO9BKIyn7pdzlm/rwqCoijKAqOsyM0DN26d8/uqD0FRFEUBVBAURVEUi6QEQUSuFZHDInJERO5NcLxQRP5DRF4VkQMicrvVvlpEnheRN6z2T0T1uU9EWkRkn/XnPTM3LUVRFCVVJvUhiIgT+DpwNdAM1IvIM8aYg1Gn3Q0cNMa8X0RKgMMi8jgwAvydMWaPiOQDu0Xk2ai+XzXGfHlGZ6QoiqJMiWQshFrgiDHmmDEmADwB7Iw7xwD5IiJAHuAHRowxp4wxewCMMT3AG0DZjI1eURRFmTGSEYQy4GTU52bGPtS/BqwHWoHXgU8YY2IqRYuID9gK/Cmq+R4ReU1EvicinkQ3F5E7RaRBRBra2tqSGK6iKIoyFZIRBEnQFl/O593APqAUqAK+JiIFkQuI5AE/AT5pjLErSH8TWGOdfwr4SqKbG2O+ZYypNsZUl5SUJDFcRVEUZSokIwjNwOqoz+WELYFobgd+asIcARqBdQAi4iIsBo8bY35qdzDGnDHGBC1L4tuEl6YURVGUeSKZwLR6YK2IVAItwI3Ah+LOOQHsAF4WkRXAhcAxy6fwXeANY8y/RHcQkVXGmFPWxxuA/ZMNZPfu3edE5HgSY07EMuDcFPsudnTu6Ue6zht07onmXpFMZ0mmmLO1JfQBwAl8zxjzzyJyF4Ax5mERKQUeBVYRXmK63xjzf0XkMuBlwn4F26fwWWPMLhH5AeHlIgM0AX8dJRAzjog0GGOqZ+v6Cxmde/rNPV3nDTr36cw9qdQVxphdwK64toejfm8FrknQ7xUS+yAwxtyS0kgVRVGUWUUjlRVFURQgvQThW/M9gHlE555+pOu8Qec+ZZLyISiKoihLn3SyEBRFUZQJUEFQFEVRgDQRhMmytS4VxssuKyJeEXlWRN6yfiZME7IUEBGniOwVkZ9bn9Ni7iJSJCI/FpFD1t//O9Jh7iLyN9a/9f0i8m8ikr2U522l+TkrIvuj2sadr4h8xnruHRaRd092/SUvCFHZWq8DNgA3iciG+R3VrGFnl10PvB2425rrvcBzxpi1wHPW56XKJwgnUbRJl7k/CPzKGLMO2EL4O1jScxeRMuB/ANXGmE2E46RuZGnP+1Hg2ri2hPO1/u/fCGy0+nzDeh6Oy5IXBJLL1rokmCC77E7gMeu0x4APzM8IZxcRKQfeC3wnqnnJz93KG3Y54awAGGMCxphO0mDuhGOp3CKSAeQQTquzZOdtjHmJcDbpaMab707gCWPMkDGmETjCJCmC0kEQksnWuuSIyy67wo4Ct34un7+RzSoPAP+L0ah4SI+5nw+0AY9Yy2XfEZFclvjcjTEtwJcJp845BXQZY37DEp93Asabb8rPvnQQhGSytS4pxskuu6QRkfcBZ40xu+d7LPNABrAN+KYxZivQx9JaJkmItVa+E6gknGk5V0T+2/yOakGR8rMvHQQhmWytS4ZxssueEZFV1vFVwNn5Gt8scilwvYg0EV4WfJeI/F/SY+7NQLMxxq418mPCArHU534V0GiMaTPGDAM/BS5h6c87nvHmm/KzLx0EIZKtVUQyCTtZnpnnMc0KE2SXfQa4zfr9NuBncz222cYY8xljTLkxxkf47/h3xpj/RnrM/TRwUkQutJp2AAdZ+nM/Afz/7d07i1VnFIfx5w+CUQZJYiEiARsJaSQxKgqigcCAVilEUlgpBD9AUlmZQr9BxE60sFTERktFDDqoGRGiMmEkHyBIDCYEXSned2ArcS7COMyZ5webs8/Z7z7sdW6LfTlr7Uqytn/2v6adNxv1uN/0tngvA98mWd2rVW8Bbs/6TFU18hNwAHgMTAHHl3p7FjHOPbRdwklaw6L7Pfb1tKsPnvTbj5d6Wxf5dfgKuNLnV0TstMrBE/29vwR8tBJiB04Av9LK558HVo9y3MAF2vmSf2l7AEdnixc43n/3HgH753p+S1dIkoCVcchIkjQPJgRJEmBCkCR1JgRJEmBCkCR18+qpLK0USV4CDwYPfVNV00u0OdJ75WWn0kCS51U19pZloX1nXv3fcmm585CRNIskm3t/gZ+Au8AnSU4nmeh1+E8Mxk4nOZnkVl++LcnVJFNJjg3G/ZDkTpLJ4frSUjMhSK9bk+R+ny72xz4FzlXVF1X1lPZv9+3AVmBfkq2D9X+vqt3ADVrt+oO03hQ/AiQZp5UQ2En7d/GXSfa+j8CkuXgOQXrdi6r6fOZOLyP+tKp+How5lOQ72vdnI63x0mRfNlMn6wEwVq0vxZ9J/k7yITDep3t93BgtQVxfnHCk+TMhSHP7a2amFwn7HthRVX8kOQt8MBj7T799NZifub+KVpL4VFWdWdQtlt6Bh4ykhVlHSxDPkmygtWZdiKvAkd6zgiSbkox6AxctE+4hSAtQVb8kuQc8BH4Dbi5w/WtJPgNutYuWeA4cZvRr9msZ8LJTSRLgISNJUmdCkCQBJgRJUmdCkCQBJgRJUmdCkCQBJgRJUvcf2D8D7FzXTIEAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "ca_df.plot(x='Frame')" ] }, { "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": { "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": true } }, "nbformat": 4, "nbformat_minor": 2 }