{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Contact analysis: number of contacts within a cutoff\n", "\n", "We calculate the number of salt bridges in an enzyme as it transitions from a closed to an open conformation.\n", "\n", "**Last updated:** June 29, 2020 with MDAnalysis 1.0.0\n", "\n", "**Minimum version of MDAnalysis:** 0.17.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", "* [Write your own contacts analysis method](contacts_custom.ipynb)\n", "* [Q1 vs Q2 contact analysis](contacts_q1q2.ipynb)\n", "* [Fraction of native contacts over a trajectory](contacts_native_fraction.ipynb)" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "ExecuteTime": { "end_time": "2020-07-29T05:08:24.953983Z", "start_time": "2020-07-29T05:08:21.329652Z" } }, "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": [ "## Background\n", "\n", "Quantifying the number of contacts over a trajectory can give insight into the formation and rearrangements of secondary and tertiary structure. This is closely related to native contacts analysis; where the fraction of native contacts refers to the fraction of contacts retained by a protein from the contacts in a reference frame, the number of contacts simply counts how many residues are within a certain cutoff for each frame. No reference is necessary. Please see the [Fraction of native contacts](contacts_native_fraction.ipynb#Background) for an introduction to native contacts analysis." ] }, { "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": [ "## Defining the groups for contact analysis" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We define salt bridges as contacts between NH/NZ in ARG/LYS and OE\\*/OD\\* in ASP/GLU. It is not recommend to use this overly simplistic definition for real work that you want to publish." ] }, { "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": [ "## Calculating number of contacts within a cutoff" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Below, we define a function that calculates the number of contacts between `group_a` and `group_b` within the `radius` cutoff, for each frame in a trajectory." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "def contacts_within_cutoff(u, group_a, group_b, radius=4.5):\n", " timeseries = []\n", " for ts in u.trajectory:\n", " # calculate distances between group_a and group_b\n", " dist = contacts.distance_array(group_a.positions, group_b.positions)\n", " # determine which distances <= radius\n", " n_contacts = contacts.contact_matrix(dist, radius).sum()\n", " timeseries.append([ts.frame, n_contacts])\n", " return np.array(timeseries)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The results are returned as a numpy array. The first column is the frame, and the second is the number of contacts present in that frame." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "(98, 2)" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ca = contacts_within_cutoff(u, acidic, basic, radius=4.5)\n", "ca.shape" ] }, { "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", "
Frame# Contacts
0069
1173
2277
3377
4485
\n", "
" ], "text/plain": [ " Frame # Contacts\n", "0 0 69\n", "1 1 73\n", "2 2 77\n", "3 3 77\n", "4 4 85" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ca_df = pd.DataFrame(ca, columns=['Frame', \n", " '# Contacts'])\n", "ca_df.head()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Plotting" ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "scrolled": true }, "outputs": [ { "data": { "text/plain": [ "Text(0, 0.5, '# salt bridges')" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYIAAAEGCAYAAABo25JHAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4yLjEsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy+j8jraAAAgAElEQVR4nOyde5xcdX33P9+57szOzM7sLbvZTQgBkmxCQgIBEZSLAgIqCFTFp1aoVaq1tU+fKoq1WqxPtZUHUWutaBXEVi2SIFoUBBEo5Zp7yG5IQhKym93sbPYyMzv3md/zxzm/s2dmztw2c9md+b5fr33tzuyZmd+Z2f19z/f2+ZIQAgzDMEzzYqr3AhiGYZj6woaAYRimyWFDwDAM0+SwIWAYhmly2BAwDMM0OZZ6L6AUOjs7xYoVK+q9DIZhmEXFtm3bJoQQXcWOWxSGYMWKFXjllVfqvQyGYZhFBREdLeU4Dg0xDMM0OWwIGIZhmhw2BAzDME3OosgRMAyz+EgkEhgeHkY0Gq33UhqelpYW9Pf3w2q1zuvxbAgYhqkKw8PDcLvdWLFiBYio3stpWIQQOHnyJIaHh3H66afP6zk4NMQwTFWIRqPo6OhgI1BliAgdHR2n5HmxIWAYpmqwEagNp/o+syFgFhw73pjC7uHpei+DYZoGNgTMguOOLXvwD48O1nsZTANxxx134Pe//z0efvhhfPWrX8173I9+9COcffbZWLduHdauXYu77rprXq+3c+dOPProo/NdLu655x6Ew+F5P75c2BAwC4pYMoWD4yFMhxP1XgrTQLz44ot405vehKeffhpvfetbDY/59a9/jXvuuQePP/44Xn31VWzfvh1tbW3zej02BAxzChwcDyGZFghGk/VeCtMAfPrTn8aGDRvw8ssv481vfjO+//3v4+Mf/zi+9KUv5Rz7la98BXfddReWLl0KQCnJ/OhHPwpA2dgvvPBCbNiwATfccAOmpqYAAJdddhk+85nP4IILLsCqVavw7LPPIh6P4wtf+AJ+9rOfYePGjfjZz36Gl156CRdddBE2bdqEiy66CPv37wcApFIpfOpTn8L69euxYcMGfOtb38I3v/lNHD9+HJdffjkuv/xypFIp3HrrrTj77LOxfv16fP3rX6/4+8Tlo8yCYmg0CAAIRNkjaCTu/OWr2Hc8UNHnXLvUgy++e13BY772ta/hve99Lx544AHcfffduOyyy/Dcc88ZHrt3716cd955hr/70Ic+hG9961u49NJL8YUvfAF33nkn7rnnHgBAMpnESy+9hEcffRR33nknnnjiCXzpS1/CK6+8gn/+538GAAQCATzzzDOwWCx44okn8LnPfQ4PPfQQ7r33Xhw+fBg7duyAxWLB5OQk2tvbcffdd+Opp55CZ2cntm3bhpGREezduxcAMD1d+fwZGwJmQTE4qmwWoVgS6bSAycRVJ8ypsWPHDmzcuBFDQ0NYu3Zt2Y+fmZnB9PQ0Lr30UgDALbfcgve+973a72+88UYAwHnnnYcjR47kfY5bbrkFBw4cABEhkVAudJ544gl87GMfg8WibMXt7e05j125ciVef/11/MVf/AXe+c534qqrrir7HIrBhoBZUAyOKYZACCAUT8LTMr9OSWZhUezKvRrs3LkTt956K4aHh9HZ2YlwOAwhBDZu3Ijnn38eDocj4/h169Zh27ZteNvb3lbW69jtdgCA2WxGMmkc0vzbv/1bXH755di6dSuOHDmCyy67DIDSDFas9NPn82HXrl147LHH8O1vfxv/+Z//iR/84AdlrbEYnCNgFgxCCAyOBmG3KH+WgQiHh5j5s3HjRuzcuROrVq3Cvn378La3vQ2PPfYYdu7cmWMEAKWy6Pbbb8fY2BgAIBaL4Zvf/Cba2trg8/nw7LPPAgAeeOABzTvIh9vtRjAY1G7PzMygr68PAHDfffdp91911VX413/9V82ATE5O5jx+YmIC6XQaN910E/7+7/8e27dvn+c7kh82BMyCwR+MYXI2jo3LvADACWPmlPH7/fD5fDCZTEVDQ9deey0+8YlP4IorrsC6detw3nnnaRv0/fffryWed+7ciS984QsFX/fyyy/Hvn37tGTx7bffjjvuuAMXX3wxUqmUdtxHPvIRLF++HBs2bMA555yD//iP/wAA3Hbbbbjmmmtw+eWXY2RkBJdddhk2btyIW2+9FV/5ylcq8M5kQkKIij9ppdm8ebNYbINpJkIx/Ml9L+MbN2/Cis7Wei9n3qTTArfe9zI+dOFpuGLtkpIe89zBCXzjyQP40YcvQIvVXPJr/X7/OG794cv45NvOxDd/dxA/u+1CvGllx3yXztSZwcFBDAwM1HsZTYPR+01E24QQm4s9lj2CKrF7eBq7hmfw0pHJei/llDg2FcYzr/nx3KGJkh/z/KGTeOnwJH6770RZrzWoVgxdcLqy+QfYI2CYmsCGoEoMT0Uyvi9WZBXPRChe8mMmQjEAwNYdI2W/Vp/XgX6fEr8Ncgkpw9QENgRVYkQ1ACOL3hAoV+n+YOnKhv6gYgiefs2v/VwKQ2MBDPS64XEolUKcLF78LIbQcyNwqu8zG4IqMTytGoLp2rWJV4P5egT9PgdSaYFf7jpe0mOiiRQO+WexpscDd4tS1czJ4sVNS0sLTp48ycagysh5BC0tLfN+Du4jqBKaRzC9yD2CMWkISr+ynwjFceHKDnidAWzZMYwPv6X4sIyD4yGk0gIDvR5YzSY4rGbuLl7k9Pf3Y3h4GH6/v95LaXjkhLL5UlVDQER/CeCjAAjA94QQ9xBRO4CfAVgB4AiA9wkhpqq5jnogDcDodBSptIB5EXbIBqMJHJuMoNVmxnQ4gXgyDZulsBMphIA/GEOn24Z1S/vxpV/tw2sngli1xF3wcftUz2OgVznO3WJhj2CRY7Va5z0xi6ktVQsNEdHZUIzABQDOAfAuIjoLwGcBPCmEOAvAk+rthiKaSMEfjKHH04JkWmC8jPj6QmL/mJIfePMZnQCAk7PFvYJANIl4Ko0ulx3XbVwKs4mwZXvxpPHgaAAOqxmndSilth6HlT0ChqkR1cwRDAB4QQgRFkIkATwN4AYA1wO4Xz3mfgDvqeIa6sLojLLxX3C6ohtS78ohfzA2rxCVzA9csqpTex4948EoxmYyjZwMIXW57eh02XHpqi78YucIUuncOPHOY9P4zd4x/GbvGF46PIlVPW7Nc/K0WBCINLdHcGwyjKnZ0nMzDDNfqmkI9gK4hIg6iMgJ4FoAywAsEUKMAoD6vdvowUR0GxG9QkSvLLYY4/CUkiCWhqDelUN3bNmDj/94W9mP2zcahKfFgvV9iiZ7dp7gMz/fjU/+ZEfGfdJYdLoU/ZUbNvVhdCaK5w+dzDhuZDqCG//lOXzsx9vwsR9vw6vHA9ikdhQDgLvF2vTlox++72V87fH99V4G0wRULUcghBgkon8E8FsAIQC7AJR8iSeEuBfAvYDSWVyVRVYJufG/SRqCOieMdw9PIxxPlSRwpUcp5/Sgy61s6hPBzKvTQ/5ZhOOZH6k0FtIQXLl2Cdx2C7bsGMZbzurUjnt4xwjSAvj3j7wJPqcNRMAZXS7t9x6HFW9MLu6Kq1Pl+HQE44HSk/QMM1+qWj4qhPg3IcS5QohLAEwCOADgBBH1AoD6fbyaa6gHI9MRmAhY0dmKjlZbXUNDJ0MxjAdjCMWSZYVa0mmB/WNBDPR6tE3dr/MI0mmB0ZkIJkJxRBNz2ikTwbnQEAC0WM24dn0vfrN3TDMaQghs2T6MC1a04+IzO7F2qQcDvZ6MRLSSLG5ejyCeTGM2nsJsrLnDY0xtqKohIKJu9ftyADcC+AmARwDcoh5yC4BfVHMN9WBkKoIeTwusZhP6fI66egRDY3MKiMNl9DQcnQwjHE9hba8HLVYz3HZLRo5gPBhDIqU4avrz84diMJsIXsecfPSN5/YhHE/hsVcVVcfdwzM45J/Fjef25X19T4u1qXME0xHF+5qNN+97wNSOajeUPURE+wD8EsAn1DLRrwK4kogOALhSvd1QDE9H0KfKJPR5HRiZql+IQyZ8gfJyFUPq49ao5ZxdbnuGR6BvlNN7PBPBODpabRkDZc5f0Y4+r0OrHtq6YwQ2iwnXrO/N+/ruFgviqXSGt9FMyJnNIS6hZWpAVfsIhBA5U6KFECcBvL2ar1tvRqYiOH+FD4BiCJ7aP152fL5S7BsNoNVmxmw8VVaIanA0ABNBq//vdNm1sA+QufnrDcxEKKaFhSQmE+HGc/vw7acOYngqjEd2HceVa5egzZF/6IwmMxFNlKVg2ijIaqEQh4aYGsASExUmmUpjLBCd8wh8DkQTaZysUxng4GgQ561oh8NqLitEtW80iJVdLm0T7nTbMqqGpCEwUaZ34A/FtJyCnhs29SEtgE89uAuTs3HcuCl/WAhQykcBNG14aFrVWeIcAVML2BBUmLGA0knc73MCgPa9HiWkiVQaB8eDGOh1K7mKMj2CgV6PdrvLZc/IEYxMR+BzWrHUm/m8E8FcjwAAVna5cM4yL154fRIdrTZcsqqr4OvLEZXNmjCeDsscQQppgx4MhqkkrDVUYeSm2Od1ZHwfmY7gHF2d/HyZCSfw2S27tZCB1WzC564dwJndrpxjD/lDSKQE1vZ6MDQazPEI9o8F8Y+/GUIilc49j+kI/teblmu3O112BKJJxJIp2C1mjEwpeZBWm0V7XiEEJkJxQ48AAG46tw+7jk3j3ecshdVc+BrE41A9ggrGyH/0/JGMGQkXn9mJj116RsWev5JMhecM4Gw8CTfPbmaqCHsEFUZuivrQEFA5j2DPyAx+vXcMJwJRBKJJ/G5oHL/fb1yBO6jp93gMq5d+ues4fr9/HKFYMufrTae34x3r5iaSdcpeAlWFdGQ6gj6vI8PTCEQUeYlOl81wPddv7MO163vwxxevKHqe7ip4BD9+4Sj2jMwgFEvitRNBfPupgwtWGXMqPBdKnI01Z8KcqR3sEVSYbI+gzWGF226pWAmprMW/+30bsW6pB6s//5uMah49Q6NB2MwmrOxsRZ/XgcnZOMLxJJw25WMfHA3gzG4Xtv7ZxUVft8slm8piWNrWgpGpCC45qwsuuxljgSgSqbS2DqPQEKC8F//yh+eVdJ4yNFTJHEEwmsRVa5fgn/7gHPz4haP4/MN7MTId0cJ3C4kZnUcQiiUAzF9imGGKwR5BhRmZjqDTZcuodOnzOSrWVBaOK1eHTpsZRIROly2n41eybzSAs5a4YDGbtKlfes8kOw9QiDmPIIapcAKRRAr9PsUjSAtgbCaq5RC68oSGykGGhirpEQQiCc3TkOctB+8sNPQeQYg9AqbKsCGoMDJkoqfP69D0h04V2WDUalc2yuz6fj2Do0Ftw5NrkgNzZsIJHJ+JYk1PiYZADff4g7E5r8fnQJ9XuZoenorMyUvk8QjKwWE1w2yiiimQJlNKp670NFb3KGWx+j6LhcRUOAGL2ovBvQRMtWFDUGGGp+aaySSV7C6O6DwCILe+X+IPxjARimmGILt6SQ6ckfr/xZAJ4IlQTDNq+vnCI9NzhqASHgERVVSBVCbX5fQzl92C0zqcGBpbmIZgJpxAr1cJB3EvAVNt2BBUkHRaGMac+30OBKPJilzdysShjPN3uuyG08OGsjb6brcdVjNpBkleCa8tMTTUYjXD3WLBRCiuPUe/z6FtViNTEfiDMVhMVLBRrBwqqUAqDYpHt7aBHs+CDg31q94W9xIw1YYNQQWZmI0hnkwbhIYq10sQjidht5g03f4utx0nZ+M5teZaxZAa+jGZCL1tcxU+Q6NBtLfa8iZ2jZC9BMNTytSyNocVdosZ3W47hqfCmAjF0OHKlJc4FTwOS8XKR6URlh4BoOQJjpzMVVCtN0IITIcTmrfFHgFTbdgQVJDsiiFJJUtIw/GUlh8AlNh9Ki0ykouAkh/o8bTA1zpXyqnPVQyOBTDQ6y5L9qJTzUeMqFpK8rEy9DURipdlWIrhts/PI4gmUjllodIQeHT1+Gt63RAiU5hvIRBJpBBPpTXPslaGIJUWhj0lTOPDhqCCyMlkSw2SxQBwfObUDcFsPKnlBwCgy62EZmR9v2RoLKgJxmnrUDfsZCqtSEyXmCjWXkvNR4xMZYa/+rzK8/qDxvIS88XjKD9HEIgmcP7/fQL/tWc04345/1hWIwFzYbGFljCWzWRLPHZYTFSz0NDXHtuP93/3+Zq8FrOwYENQQaRiZHtrZkOVz6lchZ4MnbreUDiWyjAE+moePcOTYaxQ5/9K+rwOjAdjeO1ECLFkGmtKzA/oX0vzCHTGrs/nwOh0FOPBaGUNwTxyBPuOBxCMJnFofDbj/kAk1yPo9zngtlswtMDyBFJwzuu0odVuqZlHcHA8iAPjoZq8FrOw4IayCmIUhwYAi9kET4sFM5EKJIt1DWFAZn2/ZCaSQDCWNAxRCQE8pXYil1oxJOly27Ura31lVL/XgXgqjRMBY52h+eJusZadI5BX99mhMvk8ekNARFjT615wHoG8oPA5rXDV0BBMhRMIRpNIptKwFJEAYRoL/rQrSCCSgNlEGVfsEl+rLWdzmg+ReAqtdn1oKNcQ6Ov89fSrhuGJwROwmMhQn6gQ+qv9bI/A6JhTxeNQNkGjwff5kJv6dNZ7LT0LV5aRXtPjwdBYcEEJu8mhNL5WG1x2S81CQ/Lvc7oCFyzM4oINQQUJRpPwtFgME7BehzVDSGy+zMZTcFjnNjO33QKbxZSjDArkJq1lXH/nsWmc2e2C3VKezn+GIdB7BLp8QT6dofkgu4DLaaiS5aDZm1kgkoTLbtGqrSQDvR6EYsm6jhPNRv6deJ1WtNrNNfMIpKzFdAX+TpnFBRuCChKIJvKqRHqdNsxUwCMIx5MZHgERKWWdGbMC1IavLI+gp60FRIAQwJqe8sJCQKaGUL/eI9D9XMnQkDaToMQ8QTKVxmsnFEOQbXSD0UROyA6YC48NLqDGsmmZI3DY4Gqx1kRiQgihGc9sb4ppfNgQVJBgNJlRlaLH56yMRxCOpzJyBICSJ9BXDY1MRdBiNaEjK2lts5iwRK0yKlVjKPt15PPovYNWuwVeNSFeia5iiTSqpRqCIydnEUumYTVTzmYWiCYy8gOS1T1uEC2syqHpSAKtNjNsFhNcdnNNQkOB6FwIrhJ/p8zigg1BBQlEEnDb83sElcgRhGPJnBxEl8uWExrq8zoMQ1TSS5iPIZCGpc/ryGkak15BpXMEQOkKpPvUsNCmZT6t8kaSz0g7bRas6GhdUIZgKhyH16m81602S020hvSGkz2C5oMNQQUp5BF4nVatIkPPM6/58ZVfD2pf+sEp2QghEE6k0JplCLJlJpSGL2NpZblhZ/cYlEKL1QxPiyUn9yCf12qunLwEkDulbHgqjJ9vG857/OBoABYT4bwVvowrXKBw2G6g110VqYmjJ2fx2KtjZT9uOpzQPCxXS+Fk8Xgwii3b878nkp9vG8aY2ueS7zWNfq4mU7Nx/PuLRxfsTIhmgg1BBckXfgAAn3qFl11C+g+PDuJ7z7yO+547gu8/exh3bNmd9/mjiTSEAJz2TGPT5bZjcjaubXwjU7kKqJK3nNWJi8/sQLd7fvr2l6/pxqUGYyYvWdWFy1Z3V0xeAtDNJFCviO96bD8+9eAuLQ+QzZA6X6FbDWHp3+tAJKnlHLI5s8uFY1PhHCN9qvzbfx/GJ3+yo+yNbioc1/5eXHYLQvFk3ue4/3+O4P/85y5DvSnJyVAMn3pwF+7+7f6Cr2n0czX5ztOH8Ddb9+KNycoo8zLzhw1BBdHr3Wcjr/Cy46/+YAw3X7Ac+798Df7m2gFMhOIYDxpfuUkJ6uzQUKfLrslMhONJnJyNazo12bxv8zL8+0cuLOu89Hzj5k346CUrc+7/4IWn4Xsf2jzv5zViLjSUQCiWxGOvKt7Slu0jhsdL2W25ieo3tGABj6DL0wIhgMnZym6A44EYYsl02b0QMzqPoNVugRBzcyiykZ5Moaon+btf7xlDNGH8PBkeQQ3KR1NpgYd3KJ9jdjMkU3vYEFQITe8+b2hI2Zz08ddkKo3J8NyM32LDUsJZyqMSvUT08Tylo4sRl10Op0niN3vHtGE4D+8YyektmJqNYywQxUCvW9tE5eYmhECgQNiuSy15Ha/whiQrucrd6JQcgRoaUt+DfOEhmdsopGMly4mDsWTe0KP8u+xotdUkR/DcwQnt/WZDUH/YEFSIOb37fKGhzM0JUK5AhZjbiLRSxjyJy3BCHUqTnSxWQyFSGRTILR1djFjMJrTazAhEE9i6YxjL2534zNVrMBaI4oXXT2YcK9+zNT1zHoHc0CKJFFJpkd8jMGjKqwTy+cp53nRaYCaSyAgNAcomns10OK7pW41M5w+vSCPR3mrLm0+QnuryDiemZqvvEWzdMQKbRdl+Kv2+M+XDhqBCaKJmeeLQRuGK7Bm/XqcNvW0tGMpjCOQsAkdOaEh57glVBwhoDI8AUAzrayeC+J9DJ3HDpj5cuXYJ3HYLHsra0AZVBdGBXk9OGE6bRZDHEMx5VJW9EpZXuuVsdMFoEmkx50EW8gj0nmMxj8Blt+D95y/DMwcmDK/Ap8NxuFss6Gi1Vz1HMBtTPLwbNvbBROwRLATYEFQImZjMd9XZZuARyH8AfcnlQG/+YSnhrDGVEu2KNhjHyFQEFhNhiacxhp17HBY8e2ACQgA3ntuHFqsZ79zQi9/sHcuYIzA4GkCny44utz0nDCerjvKFhuT7X8kNKRxPanH9cp5XbsI+XY4AMJaill5Qp8tecALesFo8cOOmPqTSAo/sOp5zzLTqhfic1opoYhVChvn+YHM/2ltt8FfYADPlw4agQgSKbDZuuwUWE2VcbckrUH037kCvG4f8IcSSuUm9uelkmR6By26B3WKCP6SEhnq9LTlSCosVeRV/3mk+nKaqqd6wqQ/heCqjNHNwNKCF1jwtipSENLpzYoDGRrrVboHTZq5oiGIiqP+cyzcE2TkCo14CxfjZcE5/W5FkcRh9PgfOWuLG+r42bN2RGx6aCifgc1orpolViC07hrGs3YHNp/nyTthjagsbggoRNFC31ENE8DqtGRUZ2rD3LI8gmRY4aCAHHNFyBJnGhoi02cXZEtGLHSkLceO5fdp9569oR7/PoVUPJVJpHDgR0uYLECn9DHJDmwsN5Rfb7VSnr1UKf2iu8ktvFIoxrekMqaEhdc2zBlPUhsaCWNPjQX+RmdjK+FTlb+KGTX3YOxLIKcGdDsfR5rShzWFFNJHOW110qozORNQwX78ij+Ku7PvOzA82BCXyd4+8in/8zVDe3xvp3WfT5rBmVGT4gzE4rOaMUM+anvyVQ/k8AkDxKvyhmNpDYNxMthhpc1hhM5vwrvVLtftMJsINm/rw7IEJrP/iY9h45+OIp9IZTXJep7VkjwBQ3r/5Xpl+5ue78a9PH8q4z69u/lYzZehAFUNTHpWdxaquVLbeUDKVxv4TQQz0utGnzsQ2CukEooq0tLw4uG7jUphNhEd2ZoaHlN4Fqy7RXn546MP3vYyfvPRGwWP+a/cohFAMEpB/5jZTW3geQQmcDMXw4xeOYqDXg89cbXyM0QSsbHxOW0ZFxkQoV7//9M5WtFhNhpVDMiae3VAGKP9QhydCOBGMNkTFkOSjl6zE1Wf3aDkWya0XrUA8lUYiqZSRtlhNuGJgifZ7n9OmbaqBEj6bTpcNhydm8/4+H6m0wMM7R3BOvxcfu/QM7X65+Z/Z7S4vNKT+fXgdhUNDhydmEU+mMdDr0VRkR6YiOZ3d2ZLknS47zuhqxf4cj2AuRwAohqGnrfQ8UzSRwu+GxnHIH8LN5y/LOwJ1/1gQ3W47Tu9UwnzSAAshyhqbylQWNgQl8Mtdx5FMC21jMUJedboMNmmJ12nLcOEnQrEc2WazibB6ifGwFJl8dFiNPAIbnhyahRCZyqCLnXVL27BuaVvO/R0uO+64ZiDv47wOq1ZaGTSYV5xNp8uOlw5Plr2+wxOK0F12aGYiGAMRsHqJCy+W8bzT4TiIAI+6oTusZpgot2pon65cNq52RI9MR7B2aaaGlNEc7T6vI6PKKJlKIxhNwuu0aga33DyBPP+jJ8PY/sYUzjutPe9xmfMrbIgm0gjFkgU9Nqa6cGioBLaqHZDTBeqrg9EkWm3mgpOdlHBFZmjISLZ5TY8Hg6OBHFmBcDyFFqvJMBHc5bJDHt5IHsF88TrnGqMCkSRsZhPslvyfTZfbjqlwouzh7dJgj85EMh47EYrB57Shp82Bk6F4yTIT05EE2hxW7TMmIsNxlUNjQVjNynAhucmPTOX2EmjlxL7MQUJ6wyXzVl7HXGhopszQkD5Zna/zWx6nN0rVKt1lyoMNQREOjoewa3gGnS4bgrFk3o2ikLyERJGizqwaMlLrHOh1YyqcyOl0nY0lcxLFkk79rAA2BBmy34FoAh6H8cAgifwcypWZkIYgLZAh6uYPKt5ep8uGeCpdsoLqVDihhYUkRuMqB0cDOKPLpUqC22C3mAwTxsNTYeWYVv3fh1MZZxrNHETja7Xp+l3KMwTSw7hgRTt+tXvUsOotnRYYnYlkDDLSN0My9YMNQRG27hiGiYAPXLAcQP4kmtxsCuF12rSKjEQqjalwPkOguPf7ssJDkXgKTrvxVDH5PERAbxsbAl+rDZFECtFECsFo8bDDfHsJhsZ0TV1ZYb8ut31uowvlV/7UM62ToJYYjatUymXnqqSyr/L1a+rPkg3XPAj1+GmtZNWma8YrNzQUhtlE+PhlZ2AmksBTQ+M5x4wHY0ikhOFoU04Y15eqGgIi+isiepWI9hLRT4iohYj+johGiGin+nVtNddwKqTTAg/vOI63ntWlzfedyZMnUMZUFt5s9Bo4mryEUWhI/Qcfyqocmo0n4bTm0ctRn6fbbdda95sZmTSdiSQQiCQKlo4CuivTMjekwdEAzlnmBZDZ3Su9vS7NwJS2scrqHT3ZoaHJ2ThOBGJa3wSQG/eXjExFckKF8rY8XhuN6bCixWpGi9VUdlPZyFQEPZ4WvPWsTnS57XjIIDwkZTD6DUNDbAjqSdV2DCLqA/BJAJuFEGcDMAO4Wf3114UQG9WvR6u1hlPlxcOTGBPCiVMAACAASURBVJmO4MZz+4q6zIE8oxD16GUmjLqKJW0OK/q8jpyEcbgEj6CReghOBf17XUh5VCI37IkyPAKp9fP2Nd0AMj0CJTRk10J2pW50yiyCTI/A3ZJpCKQEiX64UL5eAqO+kv4sj2Cum9mmfc8e7FMMmQS2mE14z8aleGpoPCfMZqSD1d5qY5mJBUC1Lx0tABxEZAHgBJDb277A2HZ0Eo/sOo5Hdh3Hvc8cQqvNjKvW9sxtLHn+QZShNKV5BFPhuLYxdLmNh70rw1IyDUHBHIFafZRvIE2zoRf5K6Q8KulUP4dyPALZ63HOMq8i86BudLOxJCKJlBIaKiHk9Mxrfu1v7mRoTnlU0mrLDA3tMzAEfV4HJkJxRHRy1dFEChOheI4h6HTZYTObtPXKxLC3VXldpRmvfI9AGpgbNvUjmRb41e7Mf3cjHSyzidDeWriXYGgsoOUzmOpQtfJRIcQIEd0F4A0AEQCPCyEeJ6KLAPw5EX0IwCsA/loIMZX9eCK6DcBtALB8+fJqLTODaCKF93/3BSR1EscfuGA5HDbzXFgnj8usJItL8whmwgntCq/LZVyrvWqJG0/t9yOVFloFSTieQkeeUZAuuwVL21qwbmn5IygbkTltJ9UjyDNCVOK0WdBqM5fVBTyobcjujBi93ttrc1hhMVHejW7b0Sl86AcvZdy3vD3TmLfaM8dVDo4GFW9D97eghXumI1oY06hiCFAa8pZ6WzCs8wjMJoJbLX32OW15Q6BGJFJpjAXm+lfWLvXgrG4XfrvvBD705hXaccNTEXid1hytrM6sUat60mmBm/7lf/DBC0/DHdfmLxdmTo2qGQIi8gG4HsDpAKYBPEhEHwTwHQB/D0Co3/8fgA9nP14IcS+AewFg8+bNNZllF4olkUwL/MXbzsT1G5XOx9M6lH9Kr25jMVhrWTmCqXBCi8F25vEIlnhatGEz8h8+HM8dUykhIjzx15dqzUXNjj6UF4gU9wgApfKqnFi11Prpdreg3+fAvuOKYZjQqcqaTFSwe/ah7cNosZqw5eMXw2ZRSoNXdGQagpzQ0FggIz8AQOsm1xuCYYMeAkm/z5mRI/A6rFpVla/VitdO5Eqc5GNsJoq0yKxW27jMi6f2ZyaMR6YihhVtSle8seGZjiQwG09hz8hMyethyqeaoaErABwWQviFEAkAWwBcJIQ4IYRICSHSAL4H4IIqrqEs5OCX5e1OnNntwpndLljVvgCXJhqX6xFEEikk06JoaEgft54IxdBqM+cMmZEYaeSH40nDrmKJ02ZpGLG5U0W+1/5gDJFEqqiRBpQ8QTmx6sGxgCYJ0q8ma9NpodOQUtbQ6Ta+4o0mUvjVruO4el0P1i714MxuF07vbM0pc221mzEbT0EIkaOrJOnPSgDrf+5vzw0X9nnnPJiZSGY4qs1R3nCauZDP3OsM9Hpypu3l08HqUnWyjJDvpVFfDVM5qmkI3gBwIRE5SfnLfjuAQSLq1R1zA4C9VVxDWWiDXww2W0U0zmZYPiprxIuFhmRFxrSaLO40qBiSGJUzhuMpOA26iplcWqwm2CwmHFPn4Rb7bIDydG+SqTReOxHSrsz7fA7EU2lMhGLaZybzA8rz5m6sTw2NIxBN4oZz+wu+VqvdglRaIJpI43X/LOKpdEZ+AFA8SIuJMgbUyJLOJQZ/Z30+B/zBGKKJFKZm5wbhAEp+ZTqcKHnjNUoCZ0/bE0Lk1cHqVHWyjF5PvpdT4QROBDihXC2qZgiEEC8C+DmA7QD2qK91L4B/IqI9RLQbwOUA/qpaayiXfINfJNmdwZJSJAy053AoxkSRlyhkCOaGzQBKrFSpGmJVkFIgIvicVm0wejFvDVCv3Es0BHqtH2Au/DI8HYE/pMhEtLcqn2FXHgOzZccIut12XHxGR8HXcutmEgyNqdISWaEhs4nQ09aS4xH0eFoMu93lekdnohmjMQHFm0qmheEMBCPka/bqtImkgZQVTlPhBCKJlGHXe5fLjngybTiFTf++DY4ZD2xiTp2qVg0JIb4ohFgjhDhbCPFHQoiY+n29EGKDEOI6IcRoNddQDtrglzzhmuzOYMmcumXxTdqrdrxOhGLaFaMR+mEzgBJ+UtbGHkGp+Jw2zRCUomPT5WrBdIkyE3qtHyCzNt8fjKHdadM24E6dsJpkcjaOp4bGcf3GpQVlSYA5D3U2lsS+0QBsZhPO6HLlHKcP9wC5uj4Zx6r3D0+FMRPJLFk1GqJUiJHpMLrcdrTovFU5bS97prJRaEjmyYzCQ3qPON8IV+bUKWoIiOhiImpVf/4gEd1NRKdVf2m1R4q6Gck8A8gfGtLULYtvNj5VA0cJDRknioHMYTOlrI3Jxeu0YiygxKiLNZQBcxvSyRJ0b/RaP0Bmt262t9fpsiOREhlNWr/arQgZ3lgkLARkTikbHA1m5K709PkcOR5BPgHCOX2iSE4Tm9FY1ULki/3rp+0NqzpIhslitXLOKI/iD8VgM5vUvhrjyX3MqVOKR/AdAGEiOgfA7QCOAvhRVVdVJzSZ53yGwGHNkyOQoaHSPIKJUAxT4UTe0lEgc9hM5to4NFQqXodNE+IrzSMoXWZCr/Ujn7/NYcXIVCRHXtwo8f/Q9hGs6XHnxPqNcGcYgkDex/T7nBgLRJFIpbWSzny6U71tyhS7Q/4Qool0hkfgK9cjyFMNpJ+2Jz0Vo+M0j8DAAE8E4+hy2w37apjKUYohSArFp70ewDeEEN8A4C7ymEWJzBEYJYsB5B3jFygynUyP12nDMfWqrZBHAMwNm8lcG3sEpeJrnfs8Si0fBYp3AQshMDgayKnckaEZKTinPa/6sxQRPOQPYdex6Yypa4WQf49vnAzDH4zllI5K+r0OpAXw/u8+j/d/93mkRX4lWovZhB5PC15VS171OYJCekNTs3H85U93aAJ76bTA8WnjGRhreuam7Q1PRdBqM+fMSwD0hRG5ekx+Vap9oNeD1/2hqk1Oa3ZKMQRBIroDwB8B+C8iMgNoSOHwSNHQkBWxZDqjexPQD0cvJTRkRUptWCuULJa/l1enckwlewSlo7/KLeWz0TyCIoZgz8gMTgRiOP/0TM39Pp8Dw1PhHI+gWzMwysa6dfsITAStV6UY0hC8clSZa5DPI7j4rE5cuqoLVrMJVrMJl67qwsVnduZ93j6vA3vV+nx91ZB834w8gp9vG8Yvdh7HAy8cUc8phngqbRiC0lcOyXyFkQKsz2mD2UR5PAIlzLamx4O0AA6U0d/AlE4pu8r7AfwvAB8WQowR0XIAX6vusurDbJHwi9eh/oNE4nDY5v7wA5EkrGYqqHcv0f/DGQnO6ely27Hz2LSytgJjKhljpJwzEeAqwYCWqkC6ZfsIbBYTrl3fm3F/n9eB3+8fRyIlcnIE8nnTaYGtO0Zw8ZmdWOIpbQKYLEJ4+YjSgJ/PEPR5Hbj/w6W35fT5HHjpiGJc9NLX8mcjQ7BFnc3x8I7j+OsrV2vdyUYegX7a3shU/lnaisyEzdATmwjFsKG/TfOCBkcDWN+fO6iIOTWK7lxCiDEADwGQf9kTALZWc1H1IhxPwW4xHvwCzMVOp7IG1ASjCXharCWN2tOPXCxUNaT83obJ2RhSacE5gnkgja7LbsmQYc6Hw2aGy24pGBpKpNJ4ZNdxXDmwJCfM0e9zIJHK9fbaHFZYzYrMxEtHFCHDm0pIEkukR3B4YhZLPHatLPVU0W/Meu/JYjbBbbfkhIYGRwMYHA3g/BU+jExH8OLhSV33cm5/gJy2NzQWKFjBBGR6v5J0WuDkrNJZf1pHKxxWM5eQVolSqoY+CqUf4LvqXX0AHq7moupFOJ7Mmx8A9C5z5j9IIJosqXQUyPQIioaG3HakhVJqyDmC8pGx7lJyNxJlhm7+apmn9/sxORvXhq/r0W+sem9Pn/jfun1EETJctyTn8flwWs2Q1xiyXLUS6DdmfT4FUATosv/Ot+4YgcVE+Pr7N6LVZsbWHcM5M5GzGej1YNexGcxEEobGQtJlIO8xFY4jlRbodCmho9U9nDCuFqXkCD4B4GIAAQAQQhwA0F3NRdWLcCxVMPSST3guGE2UFIMG5rwKl92St3FNokkjh2IIJwo3uzG5+NQr51KNNCAF0PIPkdm6YwTtrTZcuror53f6yVvZRr7TZcexqTAe3TOKq8/uLcuzM5lI620ppcqoVPQVPL4s6Wuf05bxd55KCzy8YwSXre5Gv8+Jq8/uxaN7xnDIH4LXac07q3tNj1trTCs0Oa/TZcsxwDJX0+VWQmhK5VCQpSaqQCmGICaE0D4hVVK6IT+J2XiyoCHIV19divKoRBqT7KH1RsgqFn8whnCscLMbk4uMdZfjEeSTgwCUITe/HTyB685ZmreOX3uerIqwTpcNLx6eRDCWxE0lVgvpkZ5gvoqh+SA9GLvFlNEMBuRKUT93cALjwZi29pvO7UMolsSvdh8vOAMjQyq7gCGQOk/6TV42U8r/lYFeD2YiCYzOlDbtjSmdUgzB00T0OShzBa4E8CCAX1Z3WfUhHE8VvFLz5qmvDpSgPDr3HKrsQJFEMZA5vWlWrVRysNZQycj3upTSUYlRiELy6J5RxJNpw7AQoHh7DqsZJgI6WjM/3y63HUIo9fsXriwsKWGEvOLOLlk9FZaqG3i2NyDv04eGtu4YgafFgrcNKMGAC1d2oLetBdFEuqAhWKMfnlPguC63XZntrJPb1sT71P8VaVSGOE9QcUoxBJ8F4IeiF/SnAB4F8PlqLqpehOOpgjF4vWicHpksLgV5lVosPwBkNiJF4kllk2F10ZKZT46g02XHdDiB/WNBHBwPZXw9+MoxnNHVig15qlbk7OD2VltOwYH8vN+zqW9en6HLboHNYsLpna1lPzYfLVYzOl32nEE4gGLUJmfjODgewuBoAL/ZO4Z3nbNUkzk3mQjvUQ1ioSt9OW3PZjYV0dbKrdjSxPvU/4PVPbJyKLfDOJUWOOSf+6zqOfoyGE0gnV5cQZOil0o6uejvVX859WU2loTPWXjCl89py5GiDkRKTxZb1H+IUkZKttoUw+MPKh4BJ4rLw2o2odNlK8n7ksir5Hfc84zh7z/9jtUFq8NWdrbiRCD3c5Kb5Y15vIlidLrsWNtLRXWJymVlV6thOLTb04JgNIkr7n5auy87pHXTuX34zu8PFTVOG/rbcMgfKmgA5Wc0Hohqsh0ToRhsFpPWWe1psSpzHwwSxl//7Wv456cOarftFhO2/e2VeXMX1WI2lsRFX/0d/ubaAdx8QW0GalWCou8SEe1Bbk5gBsp0sS8LIU5WY2H1IJIovtl6s1zmRCqt6N2XmCwGgJ/e9qaC8hISrdokFIcQgktH58FPb7uwpPda8q4NvXDZzYincq/orCbC5WsK10l8+T1nI5bMFa276dx+rO9rw1lL5hfj//INZ6MaF5lff/9GGO3PH3rzaVjR0YqUGrP3Oqw477TMBrozu9145M8vxqoi53Tn9etymjCzkZv//hNBXKQ2wfmDijCj3vAO9Ho0RVM9245O4YyuVvzlFauw7cgk7n/+KPzBWM0NwdBYEMFoEi8dmWwsQwDg1wBSAP5DvS0H0AcA3Afg3ZVfVn2YjRXOEQBzWu2SYLS0WQR6zuwufTOQMWuH1czNZPOgnPcaUMIlV5/dW/zAPHTnaRJrsZqxod877+ftbSvuQc6HfJ6pu8WKd24o/j6Uck7d7uKGuNut9Efoy0P9odyZHQO9Hjw5eALRREpLcAshMDQWwDvW9eC6c5bCbbfg/uePqhdslQullYJc/2ITyCtl97pYCHGx7vYeInpOCHGxOnqyYQgXqRoClLjz/rG5D7mcWQTzodNlx7HJMDpddjYETMNCRBjodWNI9781EYqjz5tpRNb2upEWwP6xIM5ZphihEwFFxFEmk8uV0a4k0hAcHA8inkxrooQLnVJW6SKiN8kbRHQBACmGXtrkikWAHPxSTO8/W4paTicrJzRUDnJq1myRZjeGWeys6fFg/1gQSXUehD+YO7xJNtTpPQf5szQE5cpoV5LB0QCIgERKSV4vFkoxBB8B8H0iOkxERwB8H8BH1RkFX6nm4mpJNKlq+RTZbH1OK6Yjc2P8yhlKMx+63HacnI0jGE1y6SjT0Az0ehBLpnHk5CxSaYHJ2VhOon95uxOtNnOG5yCTx7KqSJOCqbFHkE4LDI0FcZE6cW4xdUGXUjX0MoD1RNQGgIQQ07pf/2fVVlZjShV18zpsSKUFgjGld6DaoaEul6KpPzIVwfo+FttiGpc5YbkgvE4b0iK3zNqkSk3oK4eGxoLo8zo07SdF9wuYqbFHcGwqjHA8hWvX9+LlI1MZxmqhk9cQENH/yXM/AEAIcXeV1lQXShV105rKZpXegVIH188X+Y8QSRSWv2CYxc6Z3S5YTITB0YBWRWTUe7Cm14Nf7joOIQSIKGdYj8lEOZ3RtUB6AGcvbcPqJYtLF6nQ7iXLLVYDOB/AI+rtdwMwLrJexMhRkMVyBPr44/IOpxYaqlaOQO8asyFgGhm7xYwzu10YHA3gzWp4xagHZKDXg/948Q0cn4mio9WG1/0hXHt2T8YxSr9PbT2CfaNBmAhYtcSNNT1uPLV/POP3vxs6gWcPTGi3z13uw7vPWVrTNeYjryEQQtwJAET0OIBzhRBB9fbfQZGZaCikR1BM1C1beG7PyAzaHFat6aXS6K+IuI+AaXTW9LjxwuuTWlexkSbXWhlCOh5At0dR6M0W4/M6rRkzomvB4GgAKzpb4bCZMdDrwYPbhjEejKLb3YJEKo1PPbgboVgSdosJsUQaP982jHdt6C1Jvr7alJIsXg5Ab1rjAFZUZTV1pNiYSoleijoUS+KxV8fwrg29VZN+0F8RcWcx0+gM9HowFoji4LhScWPkEazWVQ7J8MuabEPgsNbcIxgamwtRabpIaj+BlC//zh+eiz1/9w787bvXIhhN4vgCEdArxRA8AOAlIvo7IvoigBfRgMPrw0XGVErmhtPE8es9o4gm0rixjCEj5dJqt2jVQuwRMI2O3ECfPTABu8Vk2BnssluwvN2JwbEABkeDcNrMOK09UxrG57TlDJCqJsFoAscmI5oooH6iGqCI9nW02nDJKkW+XO/VLARKmVD2fwH8MYApANMA/lgI8Q/VXlitkaGhYjLPsjJhOpLAlu0jWNHhxLnL598xWgpS0phzBEyjIw3B3uMz6MySl8g8zo2h0SAGRwNY3ePO8ci9TltNQ0OyQkgaAK/Tht62FgyOBjT58nfr5MtXG/RD1JO8hoCIPOr3dgBHoHgGDwA4qt7XUMyW6BFYzCa4WyzYdzyAFw6fxA2b+qse45MDatgjYBqdLrcdnWrJdCGxwIFeDw6fnMWrxwOGw3p8TitCsSTiBrpP1UALUekmyA30ejA4GtTky2/Uifa57Bac1uFcMKM3C3kEUltoGxSBOfklbzcUEVk+WkLS1+e04YnBExACebXpK4lMGHOOgGkG5MZeSLZ6TY8HQgChWBIDPbl6UnNFHbXJEwyOBtHmsKK3bU4SY6DXjUP+EH768jGc2e3K6QNa0+PWcgj1Jq8hEEK8i5RL3UuFECt1X6cLIVbWcI01QSaLS+ne9TqtSAvg/BU+LO8oLFtdCaTwFoeGmGZAGoIud/4pfvoBPUYegSzqmKlRL4HSy+DOiA6s6fEgmRbYdWwaN2zqy4kcSK9GhqXrScEcgVB0FLbWaC11JRxPosVqyhkoYoT8I7thU/WSxHo4NMQ0EzLO3lXAI+j3ObREcnbFEKDv96m+IUilBfaPBXMMkv72ewwiBwO9ilezfwF0IJeys7xAROerUhMNiyI4V9pG29Fqg81iwjvXz1+uuBy6Pco/RK211RmmHsx5BPkNgclEWNPjxniemQNeTW+oOqGh3cPT+MPvv6jMnhBAPJXOMQSnd7aixWrCpmU+Q7lv6dUMjQWxabmvKusslVJ2lssB/CkRHQUwC4CgOAsbqrqyGhOOp+AsMQb/sUvPwDvX92pyt9Xm3ecshcVE6C8wEpBhGoXVS9z4x5vW4+p1hS+0Pv+utXnDKnPzxatjCP774ASC0SRuu2QlTESwW0y4Jqu72WwifOPmTXknuPV5Fa9mIVQOlWIIrqn6KhYAs7EknNbSrrhX97g1pcNa4Gmx4v3nL55pRwxzKhBRSX/vG5flL9v2aY2f1QkNDY0qQnefu3ag4HHvWNeT93fSq1kIhqCUPoKjADoAXA/gOgAd6n0NRSRRukfAMMzCxmkzw2Y2VS1HIJPDp4oyejOoydrXi6KGgIi+AOB+KMagE8APiejz1V5YrZmNJUvOETAMs7AhIrQ5rVUJDUUTKbw+MWtYrVQua3rdCMaSGJ6KVGBl86eUne8DADYJIaIAQERfBbAdwJerubBaE46n0FGgSoFhmMVF9nzxSnHgRAiptMhoHpsv0pgMjgawrL36pej5KEVr6AgA/eBQO4BDpTw5Ef0VEb1KRHuJ6CdE1EJE7UT0WyI6oH6vb7pcZTaeLCpBzTDM4sFbJSlq2Q1cidDQmh43iOo/7L6QxMS3iOibAGIAXiWi+4johwD2Aig6jJOI+gB8EsBmIcTZAMwAbgbwWQBPCiHOAvCkervuROKpkrqKGYZZHHgd1fEIBkcDcFjNOK3DuBqoHJw2C1Z0tGKozlIThXY+KSOxDZlNZb8v8/kdRJQA4ARwHMAdAC5Tf3+/+nyfKeM5q8JsLAUnzwRmmIbB57Rh1/B08QMB7B2ZyRgtuXGZV5uSlo0Uuiul+bQU8lUOTc3G8dt9J/DWVZ3obatu6XihwTT3n8oTCyFGiOguAG8AiAB4XAjxOBEtEUKMqseMElG30eOJ6DYAtwHA8uXVLZ1Mp4VaNcQeAcM0Ct5WZVylHGlZiD99YBtGpucStucu92LLn12cc5wQyoD67J6BU2HVEjd+vXcMsWQKdsvcxeiB8RBuf2g3fvThC6puCErJEcwLNfZ/PYDTASwF0EpEHyz18UKIe4UQm4UQm7u6uqq1TABK6ShQfEwlwzCLB5/Thngyrf1/5yOWTGFkOoI/ecvpePb2y3HdOUvxxmTY8NixQBTT4URFKoYkSzxKCvZkKDOfMRGSU9qqX8RSNUMA4AoAh4UQfiFEAsAWABcBOEFEvQCgfh8v8Bw1YVYbXM+GgGEaBa+cHVIkTzA6rUwJW9PjxrJ2J87qdmEiFEfUwIDIEE4lDYEcxyk3fom8XUhqo1KU0kfw3lLuM+ANABcSkVNVMX07gEEAjwC4RT3mFgC/KH251SEck7MIODTEMI2CVxOeK1w5JENCfaqEi/yuDxVJZHVPJZUF5EafbQj8wRhMBLS35ldhrRSleAR3lHhfBkKIFwH8HErPwR71te4F8FUAVxLRAQBXqrfrihxTyXr/DNM4+JyleQQjajNXv1ep45cCcSMGTV6DowH0+xzwtFROZ0yGfvzBXI+gvdVWsaR0IfJeAhPRNQCuBdCnlpFKPABKEtAWQnwRwBez7o5B8Q4WDGEtNMQeAcM0Ct4S9YaGpyMgAnrUoTKFPQLjiWinwpxHkOm5+IPxmuQHgMIewXEopaNR9bv8egTAO6q/tNpR6phKhmEWD74SpahHpiJY4m6BzaJshz2eFphNlOMRRBMpHJ6YNZyIdiq0WM1w2y05HoE/FKtJfgAoXD66C8AuIvqxEKL+I3SqSIQ9AoZpONpKlKIemQ5nSLxbzCb0eFpyPILXTgSRFpVNFEs63Xb4s5PFwRjOyCNhXWkKhYb2ABDqzzm/b6R5BHJMJecIGKZxsFvMcNrMmgJpKJbEZx7ajU9ftRordBvs8FQE552WqXTT53NgeCqzhLQaFUOSLpcdEzqPQAiBiVBMG1NbbQpdAr+rJitYAMgcgYNDQwzTUPicNi1H8Ktdx/Ffu0exvq8NH7v0DADKmMmxmWjOBLF+rwMvvH4y477B0SCcNjOWV0EcrtNty+hsDsaSiCXTWmlptSkUGmq4mQP50KqGODTEMA2FVydFvWX7CABkyDmcCESRTAstQSzp8zkwFogikUrDalZyB/tUaQlTFap4Ol12TAQntNvSO6hVjqCUPoILiehlIgoRUZyIUkRU/5E6FUQmix2sNcQwDYVPVSA9NhnGS0cmQaRMF5NoPQRZHkGf14G0AMZmlGYzIQSGqlAxJOly2RGIJhFLKnuRTBwvhKohyT9DmUlwAIADwEcAfKuai6o14VgSDqu5KpaeYZj60ea0YjqSwNYdijfwno19OOQPaRuu1kPgywz3yNvSUByfiSIQTVbNEHRmlZDK7wvJEEAIcRCAWQiREkL8EMpA+4YhnEhxophhGhA5nGbrjhFcuLIdbx/oRjItcOCEoqSf1yPwZTaVDR5XgiBrKzCDwIgudcOXIaFayksApRmCMBHZAOwkon8ior8CUJuaphoRjiW5dJRhGhCf04bJ2TgOT8zixnP7taliMk8wPBVBR6stp1CkV20uk4ZCHr+6AlPJjOjMkpmQ8hI+Z22SxaUYgj9Sj/tzALMAlgG4qZqLqjWz8RQ3kzFMA9KmCs/ZLSZcc3YPTu9sRYvVpGkGjUxHchLFgNLk1eW2z3kEYwEsb3fCVSWpelkd5Nd5BB0ue03kJYASZhbrqoeiqtTEMjVU1DCE40k2BAzTgMgr6qvW9cCt6gOtXuLWJoINT4WxeolxuKfP68DwtNJLMDQarMhoynzIXID0CCZCsZrlB4DSqoZ+T0QeImoHsAvAD4no7uovrXaE4ym08lAahmk4ZIjnD87r1+4b6PVgcDQAIQSOT0dy8gOSPp8DI1MRhONJHD45W7VEMaDKTLRYtCSxPxirWQ8BUFpoqE0IEQBwI4AfCiHOgzJroGEIx1JcOsowDcibz+jAI39+MS5dNTfcak2PG1PhBPaNBhBNpDPkJfT0+xw4Ph3F0FgQQkDLL1SLLrddFxqK1yxRDJRmCCzqAJn3AfhVlddTF2bjpf80WgAADZdJREFUSfYIGKYBISJs6Pdm3Cev7J/Yp8zE6vMZdwr3ex2Ip9J49jWl0WttFT0CQAkP+UMxCCEUwbmFFBoC8CUAjwE4KIR4mYhWQukpaBginCxmmKZhjTQEgycA5JaOSmQS+YnBE3DZLXk9h0oh9YYC0STiyfTCyhEIIR4UQmwQQvyZevt1IUSDVQ1xsphhmoU2hxV9Xgf2jMwAgGHVEAD0qYNq9ozMYE2VpCX0dKkKpLXuIQCqO7N4UZBKC0QTae4jYJgmQoaH3HaLVmKajd5AVDNRLOl02RCMJrWS1QXlETQ6kQRLUDNMsyFLQfN5AwDg0hmJNVUsHZXIjV+Wtna6F1bVUEMTjvFQGoZpNuQVfrG4v/x9LTwCGQqSoni1TBYX3f2I6PNCiC+rP9uFELFij1nozIQT+ObvDiCaSCGkGQL2CBimWZAbe75EsaTP68C+0QDWVHg8pRHSI9g3GoDZRDWTlwAKTyi7HcCzAP4AwJfVu58HcG4N1lVVnj7gx7/992H4nFaYTYQ+r6MmFp9hmIXBae1OXLKqC5eu7ip43BVrl8BpM9ckYiA9goPjIXS02mqqhlzo7PYDeC+AlUT0LIBBAB1EtFoIsb8mq6sSsmnjd399GXyttbO6DMMsDEwmwo8+fEHR4963eRnet3lZDVYEdKidxMm0qGmiGCicI5gC8DkABwFcBuCb6v2fJaL/qfK6qspEKAaLifJWCzAMw9Qau8UMT4tybV6rWcWSQobgagD/BeAMAHcDuADArBDij4UQF9VicdViIqgIOvEgGoZhFhIyPFTLRDFQwBAIIT4nhHg7gCMAfgwljNRFRP9NRL+s0fqqgj8Uq2lpFsMwTCnIkFCt96dSMiCPCSFeBvAyEX1cCPEWIuqs9sKqSa0lXhmGYUqhc6F5BBIhxO26m7eq901Ua0G1YCIYr/kbzTAMUwy5L9VSXgIos6FMCLGrWgupFem0UDyCGr/RDMMwxZAGYCFVDTUkM5FEXcqzGIZhitGtGoLuGl+oNp2ugr8Oyn4MwzCl8M4NvTAR4cxuV01ft+kMwYTaTFbLMXAMwzCl4LRZcJNurGataLrQkPQIau16MQzDLFSazxBoHgEbAoZhGKAJDcFEKA6rmeUlGIZhJFXLERDRagA/0921EsAXAHgBfBSAX73/c0KIR6u1jmxkMxkRy0swDMMAVTQEqkLpRgAgIjOAEQBbAfwxgK8LIe6q1msXwh/krmKGYRg9tQoNvR3AISHE0Rq9Xl4Uj4ArhhiGYSS1MgQ3A/iJ7vafE9FuIvoBEfmMHkBEtxHRK0T0it/vNzpkXkyEYtxDwDAMo6PqhoCIbACuA/Cgetd3oEhbbwQwCuD/GT1OCHGvEGKzEGJzV1fhKUKloshLxDk0xDAMo6MWHsE1ALYLIU4AgBDihBAiJYRIA/gelDkHNWE6kkCK5SUYhmEyqIUh+AB0YSEi6tX97gYAe2uwBgBzPQQcGmIYhpmjqhITROQEcCWAP9Xd/U9EtBGAgDL05k8NHloVJkLcTMYwDJNNVQ2BECIMoCPrvj+q5msWYkITnOOqIYZhGElTdRZroSFXS51XwjAMs3BoLkMQisFmNsHjaDrRVYZhmLw0lSGYCMbR6bKxvATDMIyOpjIEfh5RyTAMk0NTGYIJ1hliGIbJobkMQSiGLjYEDMMwGTSNIUinBU7OxtHJpaMMwzAZNI0hmArHWV6CYRjGgKYxBBOhOACWl2AYhsmmaQwBzypmGIYxpmkMAesMMQzDGNM0hmBOXoINAcMwjJ6mMQQTLC/BMAxjSNMYAr86q5jlJRiGYTJpHkMQ5FnFDMMwRjSNIeBZxQzDMMY0kSFgj4BhGMaIpjAEqbTAyRALzjEMwxjRFIZgKhxHWgCdLtYZYhiGyaYpDIHWQ+DmEZUMwzDZNIUhmOsqZo+AYRgmm+YyBJwsZhiGyaEpDMFcaIgNAcMwTDZNYQgmQnHYLCa47SwvwTAMk01zGIKgMqKS5SUYhmFyaQpD4A/FOD/AMAyTh+YwBMEYurhiiGEYxpCmMAQsL8EwDJOfhjcEqbTA5CwLzjEMw+Sj4Q3B5KyUl2BDwDAMY0TDGwLuIWAYhilMwxsCHlrPMAxTmKYxBOwRMAzDGNPwhkCGhlhwjmEYxpiqGQIiWk1EO3VfASL630TUTkS/JaID6ndftdYAKB6B3WKCi+UlGIZhDKmaIRBC7BdCbBRCbARwHoAwgK0APgvgSSHEWQCeVG9XDTm0nuUlGIZhjKlVaOjtAA4JIY4CuB7A/er99wN4TzVfmIfWMwzDFKZWhuBmAD9Rf14ihBgFAPV7t9EDiOg2InqFiF7x+/3zfuEJnlXMMAxTkKobAiKyAbgOwIPlPE4Ica8QYrMQYnNXV9e8X1+GhhiGYRhjauERXANguxDihHr7BBH1AoD6fbxaL5xMpTEZjrPgHMMwTAFqYQg+gLmwEAA8AuAW9edbAPyiWi88GY5DCO4hYBiGKURVDQEROQFcCWCL7u6vAriSiA6ov/tqtV5/roeADQHDMEw+qlpcL4QIA+jIuu8klCqiqjMRigPgofUMwzCFaOjOYk1wjj0ChmGYvDS0IdAE59gjYBiGyUtjG4JgDC1WE1pt5novhWEYZsHS0IbgzG4XrjtnKctLMAzDFKChldhuvmA5br5geb2XwTAMs6BpaI+AYRiGKQ4bAoZhmCaHDQHDMEyTw4aAYRimyWFDwDAM0+SwIWAYhmly2BAwDMM0OWwIGIZhmhwSQtR7DUUhIj+Ao/N8eCeAiQouZzHB5958NOt5A3zuRud+mhCi6IjHRWEITgUiekUIsbne66gHfO7Nd+7Net4An/upnDuHhhiGYZocNgQMwzBNTjMYgnvrvYA6wufefDTreQN87vOm4XMEDMMwTGGawSNgGIZhCsCGgGEYpslpaENARFcT0X4iOkhEn633eqoFES0joqeIaJCIXiWiv1Tvbyei3xLRAfW7r95rrRZEZCaiHUT0K/V2U5w7EXmJ6OdENKR+/m9uhnMnor9S/9b3EtFPiKilkc+biH5ARONEtFd3X97zJaI71H1vPxG9o9jzN6whICIzgG8DuAbAWgAfIKK19V1V1UgC+GshxACACwF8Qj3XzwJ4UghxFoAn1duNyl8CGNTdbpZz/waA3wgh1gA4B8p70NDnTkR9AD4JYLMQ4mwAZgA3o7HP+z4AV2fdZ3i+6v/+zQDWqY/5F3U/zEvDGgIAFwA4KIR4XQgRB/BTANfXeU1VQQgxKoTYrv4chLIZ9EE53/vVw+4H8J76rLC6EFE/gHcC+L7u7oY/dyLyALgEwL8BgBAiLoSYRhOcO5Qxuw4isgBwAjiOBj5vIcQzACaz7s53vtcD+KkQIiaEOAzgIJT9MC+NbAj6ABzT3R5W72toiGgFgE0AXgSwRAgxCijGAkB3/VZWVe4BcDuAtO6+Zjj3lQD8AH6ohsW+T0StaPBzF0KMALgLwBsARgHMCCEeR4OftwH5zrfsva+RDQEZ3NfQtbJE5ALwEID/LYQI1Hs9tYCI3gVgXAixrd5rqQMWAOcC+I4QYhOAWTRWOMQQNRZ+PYDTASwF0EpEH6zvqhYUZe99jWwIhgEs093uh+I+NiREZIViBP5dCLFFvfsEEfWqv+8FMF6v9VWRiwFcR0RHoIT/3kZEP0ZznPswgGEhxIvq7Z9DMQyNfu5XADgshPALIRIAtgC4CI1/3tnkO9+y975GNgQvAziLiE4nIhuU5MkjdV5TVSAighInHhRC3K371SMAblF/vgXAL2q9tmojhLhDCNEvhFgB5TP+nRDig2iOcx8DcIyIVqt3vR3APjT+ub8B4EIicqp/+2+Hkhdr9PPOJt/5PgLgZiKyE9HpAM4C8FLBZxJCNOwXgGsBvAbgEIC/qfd6qnieb4Hi+u0GsFP9uhZAB5RqggPq9/Z6r7XK78NlAH6l/twU5w5gI4BX1M/+YQC+Zjh3AHcCGAKwF8ADAOyNfN4AfgIlH5KAcsX/J4XOF8DfqPvefgDXFHt+lphgGIZpcho5NMQwDMOUABsChmGYJocNAcMwTJPDhoBhGKbJYUPAMAzT5FjqvQCGWQgQUQrAHt1d7xFCHKnTchimpnD5KMMAIKKQEMKV53cE5X8lbfR7hlnscGiIYQwgohWqvv+/ANgOYBkRfYeIXlF18O/UHXuEiP6BiJ5Xf38uET1GRIeI6GO64z5NRC8T0W794xmm3rAhYBgFBxHtVL+2qvetBvAjIcQmIcRRKN3pmwFsAHApEW3QPf6YEOLNAJ6Foh3/B1BmQ3wJAIjoKiit/hdA6QY+j4guqcWJMUwxOEfAMAoRIcRGeUOV8z4qhHhBd8z7iOg2KP83vVAGHu1Wfyd1rPYAcAllLkSQiKJE5AVwlfq1Qz3OBcUwPFOd02GY0mFDwDD5mZU/qOJdnwJwvhBiiojuA9CiOzamfk/rfpa3LVCkgb8ihPhuVVfMMPOAQ0MMUxoeKIZhhoiWQBmBWg6PAfiwOjMCRNRHRI0+OIVZJLBHwDAlIITYRUQ7ALwK4HUAz5X5+MeJaADA80oREkLA/2/nDmgAgEEgiCEU/zrQsV1r4hLyYXb+/5nPA8xHAeKchgDihAAgTggA4oQAIE4IAOKEACBOCADiDtsdk5dws2zCAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "ca_df.plot(x='Frame')\n", "plt.ylabel('# salt bridges')" ] }, { "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 }