{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# 01. Quickstart" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "# disable GPU. Remove this if you've compiled HOOMD for GPU\n", "import os\n", "os.environ['CUDA_VISIBLE_DEVICES'] = '-1'\n", "\n", "# import the hoomd, htf packages\n", "import hoomd\n", "import hoomd.htf as htf\n", "import tensorflow as tf" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Build the SimModel\n", "\n", "We prepare the computations that will be executed at each step during the simulation. We can have access to the neighbor list, positions, types, box dimensions of the simulation, but here we only use the neighbor list. Here we define a piece-wise repulsive potential, compute an RDF, and use the auto-differentiation tool to compute forces. " ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "class WCAPotential(htf.SimModel):\n", " def setup(self):\n", " self.avg_rdf = tf.keras.metrics.MeanTensor()\n", " def compute(self, nlist):\n", " # Use Weeks-Chandler-Anderson (WCA) repulisve potential\n", " r12 = htf.nlist_rinv(nlist)**12 # nlist_rinv is neighbor 1 / r^12\n", " # make it so anything above 2^1/6 is 0\n", " r = tf.norm(nlist[:,:,:3], axis=2)\n", " pair_energy = tf.cast(r < 2**(1/6), tf.float32) * r12\n", " particle_energy = tf.reduce_sum(pair_energy, axis=1) # sum over neighbors \n", " forces = htf.compute_nlist_forces(nlist, particle_energy)\n", " # compute rdf\n", " inst_rdf = htf.compute_rdf(nlist, [0, 3.5])\n", " self.avg_rdf.update_state(inst_rdf)\n", " return forces" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Running the simulation\n", "\n", "Now we run the simulation using the usual hoomd-blue syntax. We can specify things like how often the model is called, how often it is saved, etc. in the `attach` command. This simulation is 144 particles in 2D, whose forces are the ones we defined above. " ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "HOOMD-blue 2.5.2 DOUBLE HPMC_MIXED TBB SSE SSE2 SSE3 \n", "Compiled: 04/30/2019\n", "Copyright (c) 2009-2019 The Regents of the University of Michigan.\n", "-----\n", "You are using HOOMD-blue. Please cite the following:\n", "* J A Anderson, C D Lorenz, and A Travesset. \"General purpose molecular dynamics\n", " simulations fully implemented on graphics processing units\", Journal of\n", " Computational Physics 227 (2008) 5342--5359\n", "* J Glaser, T D Nguyen, J A Anderson, P Liu, F Spiga, J A Millan, D C Morse, and\n", " S C Glotzer. \"Strong scaling of general-purpose molecular dynamics simulations\n", " on GPUs\", Computer Physics Communications 192 (2015) 97--107\n", "-----\n", "HOOMD-blue is running on the CPU\n", "notice(2): Group \"all\" created containing 256 particles\n", "notice(2): Force mode is FORCE_MODE.tf2hoomd \n", "notice(2): Starting TensorflowCompute \n", "notice(2): completed reallocate\n", "notice(2): Setting flag indicating virial modification will occur\n", "notice(2): -- Neighborlist exclusion statistics -- :\n", "notice(2): Particles with 0 exclusions : 256\n", "notice(2): Neighbors included by diameter : no\n", "notice(2): Neighbors excluded when in the same body: no\n", "** starting run **\n", "Time 00:00:03 | Step 1000 / 1000 | TPS 290.852 | ETA 00:00:00\n", "Average TPS: 290.807\n", "---------\n", "-- Neighborlist stats:\n", "53 normal updates / 10 forced updates / 0 dangerous updates\n", "n_neigh_min: 58 / n_neigh_max: 68 / n_neigh_avg: 62.7969\n", "shortest rebuild period: 13\n", "-- Cell list stats:\n", "Dimension: 3, 3, 1\n", "n_min : 27 / n_max: 31 / n_avg: 28.4444\n", "** run complete **\n", "** starting run **\n", "Time 00:00:05 | Step 2000 / 2000 | TPS 488.064 | ETA 00:00:00\n", "Average TPS: 487.919\n", "---------\n", "-- Neighborlist stats:\n", "52 normal updates / 10 forced updates / 0 dangerous updates\n", "n_neigh_min: 58 / n_neigh_max: 66 / n_neigh_avg: 62.5703\n", "shortest rebuild period: 14\n", "-- Cell list stats:\n", "Dimension: 3, 3, 1\n", "n_min : 26 / n_max: 32 / n_avg: 28.4444\n", "** run complete **\n" ] } ], "source": [ "########### Hoomd-Sim Code ################\n", "hoomd.context.initialize('--mode=cpu')\n", "# this will start TensorFlow, so it goes\n", "# in a with statement for clean exit.\n", "#\n", "# if calling initialize() without params, \n", "# it will be throw error of using unexpected parameter 'f'\n", "# ref: https://github.com/glotzerlab/hoomd-blue/blob/master/hoomd/context.py#L204\n", "\n", "L = 16\n", "N = L * L\n", "model = WCAPotential(64)\n", "tfcompute = htf.tfcompute(model)\n", "\n", "# create a square lattice\n", "system = hoomd.init.create_lattice(unitcell=hoomd.lattice.sq(a=1.2),\n", " n=[L,L])\n", "nlist = hoomd.md.nlist.cell()\n", "# NVT ensemble with starting temperature of 1\n", "hoomd.md.integrate.mode_standard(dt=0.005)\n", "hoomd.md.integrate.nvt(group=hoomd.group.all(), kT=0.5, tau=0.5).randomize_velocities(seed=1)\n", "tfcompute.attach(nlist, r_cut=5)\n", "#equilibrate\n", "hoomd.run(1000)\n", "# reset rdf statistics\n", "model.avg_rdf.reset_states()\n", "hoomd.run(1000)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Analysis\n", "\n", "Now we'll plot RDF" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[]" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAWoAAAD4CAYAAADFAawfAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuMCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy86wFpkAAAACXBIWXMAAAsTAAALEwEAmpwYAAAhKklEQVR4nO3deXRcZ53m8e9bm/bNlmR5t+x4iZM4cSKchOzbkI1tQiDQhAayMcMauoemGXo43YfT0/R0gBCWJCR0WBMg0EwGAnQgOwl2lI0sdrwm3i3JtlxSSbW/80eVZEWWrJJct+5V3edzjs6RVOWqn8vWU69+912MtRYREfGugNsFiIjI0SmoRUQ8TkEtIuJxCmoREY9TUIuIeFzIiQdtbm62ixYtcuKhRUTK0rPPPttjrW0Z6zZHgnrRokV0dnY68dAiImXJGPPGeLep9SEi4nEKahERj1NQi4h4nIJaRMTjFNQiIh6noBYR8TgFtYiIxymoBYCN+/p4est+t8sQkTEoqAWAW/+4iU/d97zbZYjIGBTUAkDvQJLuvgT7+xNulyIioyioBYDoYBqA1/b2uVyJiIymoBYADg2mANigoBbxHAW1ABCN54JaI2oR71FQC9ZaokMj6n0KahGvUVAL/Yk0WQvhoGHj3j6yWZ1ML+IlCmohGs9dSDxpbgODqQzbDwy4XJGIjKSgluG2x5r2mYAuKIp4jYJahmd8dCxswhhdUBTxGgW1DI+o2xoqWTijmtf2RV2uSERGUlDLcI+6vjLM8rY6NuzRiFrESxTUMjyibqgKs7ytntf3x4inMi5XJSJDFNQy3KOurQyxoq2OrIVN+/pdrkpEhiiohWg8RV1FiGDAsLytDoANe9WnFvEKBbUQHUxTXxUGYNHMGipCAc38EPEQBbVwaDA1HNTBgGHZrDrNpRbxEAW1EI2nqK8MDX+9vK2OjdrzQ8QzFNRCdDBFQ35EDTCnsYru/gTpTNbFqkRkiIJaiI5ofQC01FVgLeyPJV2sSkSGKKiFaDxNfeXhoG6tqwCgu0/Hcol4gYLa59KZLP2J9JtaHy0KahFPUVD7XN/Q8vGqwxcTW2pzQd3VF3elJhF5s4KC2hhzszHmFWPMy8aYe40xlU4XJqUxdATXyNaHRtQi3jJhUBtj5gKfAjqstScCQeAapwuT0hg6fXzkxcTKcJD6ypCCWsQjCm19hIAqY0wIqAZ2O1eSlNKhERsyjdRSV0GXglrEEyYMamvtLuDfgO3AHuCQtfY/R9/PGHOjMabTGNPZ3d1d/ErFEcOtjxE9aoDWukqNqEU8opDWRxPwTqAdmAPUGGM+OPp+1to7rbUd1tqOlpaW4lcqjhja4nRkjxpyI+rufgW1iBcU0vq4GNhmre221qaAXwJvdbYsKZWjtj6iCazVieQibiskqLcDZxhjqo0xBrgIWO9sWVIq0XiKYMBQHQm+6futdRUMpjLEkjpAQMRthfSo1wL3A88BL+X/zJ0O1yUlEh1MU18ZIvcefJim6Il4R2jiu4C19kvAlxyuRVxwaNSGTEOGgrorGqe9uabUZYnICFqZ6HPR+Js3ZBrSWpdb06QLiiLuU1D7XHQwdcSMD1DrQ8RLFNQ+F42nx2x9NFaFCQWMFr2IeICC2udyx3AdeakiEDC5udQKahHXKah9brzWB6CgFvEIBbWPxVMZEunsmBcTIbfdqVofIu5TUPvY4X0+xg7q1nqNqEW8QEHtY8NbnFaOPZ2+pbaCA7EEmayWkYu4SUHtY0Mj6rFmfUCuR521sD+mUbWImxTUPja0IdO4Per8opeuqIJaxE0Kah8bb4vTIcOLXrQ6UcRVCmofi+YPth2v9dGq1YkinqCg9rGhEXXdeBcTFdQinqCg9rHoYIqKUIDKcHDM2yvDQep0yK2I6xTUPhaNj73F6UhanSjiPgW1j+X2+ZggqGsr6OqLl6giERmLgtrHooPpcfvTQ1rrdRq5iNsU1D7Wn0hTW3H0oNZ+HyLuU1D72EBy4qCe3VDJQDIzvDhGREpPQe1jsUSG6sjRg3peUxUAOw8OlKIkERmDgtrHcq2PsafmDZnXVA3AzoODpShJRMagoPYpay2xRJqaCVofh0fUCmoRtyiofSqZyZLO2gmDurE6TE0kyI4Dan2IuEVB7VOxRAaAmsjRWx/GGOY1VWtELeIiBbVPxRK5DZkmGlEDzJ9RpYuJIi5SUPtUfz6oJ5qeB7kLirsODmKtTnoRcYOC2qcGkrmgri4oqKvoS6SHj+4SkdJSUPtUf75HPdH0PDg882OH2h8irlBQ+9RketSaSy3iLgW1Tw0H9QQrE0GrE0XcpqD2qcmMqBuqwtRWhDSiFnGJgtqnYsn8POoCetS5udSaoifiFgW1T/Un0oSDhorQxEENaNGLiIsU1D41kEhPuHPeSLkRteZSi7ihoKA2xjQaY+43xmwwxqw3xpzpdGHirP5EpqDFLkPmNVXRn0hrX2oRFxQ6or4V+J21dgVwMrDeuZKkFHI75xXW9gBN0RNx04RBbYxpAM4F7gaw1iattb0O1yUOiyUn3uJ0JE3RE3FPISPqdqAb+HdjzPPGmLuMMTWj72SMudEY02mM6ezu7i56oVJcsUS6oDnUQ+ZrRC3imkKCOgScCnzHWrsaiAGfH30na+2d1toOa21HS0tLkcuUYoslMpNqfdRXhajTXGoRVxQS1DuBndbatfmv7ycX3DKN9RdwustIxhjmai61iCsmDGpr7V5ghzFmef5bFwGvOlqVOC6WnFzrA2D+DM2lFnFDoT+pnwR+bIyJAFuBjzhXkpTCQCIzqRE15C4oPrW5B2stxhiHKhOR0Qr6SbXWvgB0OFuKlEoynSWZyRa0xelI85qqiSUzHIglmVlb4VB1IjKaVib60GQ2ZBppcUtuss+W7ljRaxKR8Smofah/ElucjnRcSy0Am7v6i16TiIxPQe1DA8M7500uqOc2VlEVDiqoRUpMQe1DwyPqSfaoAwHD4pYaNncrqEVKSUHtQ1PtUQMc11rLFo2oRUpKQe1DQyeQT7ZHDbk+9a7eweGwFxHnKah96PAJ5FMbUQNs1cwPkZJRUPtQbIo9ajgc1Ju7+4pak4iMT0HtQ/3H0KNeOLOGYMBo5odICSmofWggmSYYMFSEJv/PHwkFWDizWkEtUkIKah+KJTLURIJT3q/juJZaBbVICSmofag/kZ7ShcQhx7XW8sb+AVKZbBGrEpHxKKh9KJZIU32MQZ3OWt7Yr5kfIqWgoPahWHLyW5yONDzzQ+0PkZJQUPtQLJGe9BanIy3R5kwiJaWg9qFYIk31FFYlDqmpCDGnoVJBLVIiCmofOtaLiQBLWmu1OZNIiSiofWggObkTyMeS25wpRjZri1SViIxHQe1Dkz2BfCxLW+sYTGXY1avDbkWcpqD2mVQmSzKdndLOeSMtm5W7oLipS3t+iDhNQe0zx7IX9UhLZ9UBsHGf+tQiTlNQ+0wsObTF6bH1qBuqwrTVV7Jxr0bUIk5TUPvM0Ij6WKbnDVnWVsdr+xTUIk5TUPvM0Banxzo9D2BZa25zpoxmfog4SkHtM8XqUUNuRJ1IZ9l+YOCYH0tExqeg9plY/hiuY51HDbA8f0HxNfWpRRyloPaZ4RF1EXrUS/NT9DaqTy3iKAW1z8SSxWt9VEdCzJ9RpaAWcZiC2mdix3AC+ViWz6pTUIs4TEHtM7FEmoCBynBx/umXzapja3eMZFqnvYg4RUHtM/2JNDWR0JTPSxxt2aw60lnL6zrtRcQxCmqfiRVhQ6aRlmnmh4jjFNQ+U4wtTkda3FJDMGDUpxZxkILaZ4qxxelIleEgi2ZWK6hFHFRwUBtjgsaY540xv3ayIHFWLN+jLqZls+q0i56IgyYzov40sN6pQqQ0ij2ihlxQv74/RjyVKerjikhOQUFtjJkHXAHc5Ww54rS+eJr6yuIG9fK2OqyFTRpVizii0BH114HPAeNOljXG3GiM6TTGdHZ3dxejNnFAXzxFXZGDekVbbubH+r3Roj6uiORMGNTGmCuBLmvts0e7n7X2Tmtth7W2o6WlpWgFSvFYa+lPpKmrDBf1cRfOrKE6EuTV3QpqEScUMqI+C3iHMeZ14D7gQmPMjxytShwRS2bIWoo+og4GDCva6nh1j4JaxAkTBrW19u+ttfOstYuAa4CHrbUfdLwyKbq+eAqg6CNqgJVz6lm/J4q1OkRApNg0j9pH+uK5nfOKPaIGWDm7gb54mp0HB4v+2CJ+N6mgttY+aq290qlixFmHR9QOBPWcegC1P0QcoBG1j0SHR9TFb30sn1VHwKALiiIOUFD7yFDro9jzqAGqIkHam2s0ohZxgILaR5y8mAiwck6DRtQiDlBQ+4iTFxMBVs6uZ1fvIIcGUo48vohfKah9pC+eIhgwVEeKt83pSLqgKOIMBbWP9MXT1FYU73SX0VbOVlCLOEFB7SN98bRjbQ+AlroKWuoq1KcWKTIFtY/kNmRy5kLikJWzcysURaR4FNQ+EnV4RA1w/Ox6NnX16VRykSJSUPtIdDDlyBzqkVbOqSeVsWzq0tFcIsWioPaRXI/a2dbHKfMaAeh8/aCjzyPiJwpqH3Hi0IDRFsyspr25hkde63L0eUT8REHtE4cPDXA2qAHOX97C01v2M5jUGYoixaCg9onDhwY42/oAuGB5K4l0lj9v3e/4c4n4gYLaJ5zc4nS0Ne0zqAoH1f4QKRIFtU/0ObjF6WiV4SBvXTKTR1/r1okvIkWgoPaJoRG109Pzhpy/opXtBwbY2hMryfOJlDMFtU84eWjAWM5fljuJ/pENan+IHCsFtU84eWjAWObPqGZpay2PvtZdkucTKWcKap9w+tCAsVywopV12w4QS6RL9pwi5UhB7RNOHxowlvOXt5DMZHlyc0/JnlOkHCmofcLpQwPG8pZFM2ioCvP7V/aW7DlFypGC2iecPjRgLOFggIuPn8UfXt2n3fREjoGC2iecPjRgPJed2EY0nuapLWp/iEyVgtonSnFowFjOXtpMbUWI372s9ofIVCmofaIUhwaMpTIc5MIVrfz+lb2kM2p/iEyFgton+uLpks2hHu3yk9o4OJBi3bYDrjy/yHSnoPYJt1ofAOcta6UqHOTBl/e48vwi052C2ifcupgIUBUJcv7yFn7/yj6yWW3SJDJZCmofKOWhAeO57KTZdPcl6HxDR3SJTJaC2gcGkhkyWeta6wPgwhWt1ESCfO/Jba7VIDJdKah9wI3l46PVVoS44dzF/O6VvTy3XaNqkclQUPuAGxsyjeX6cxbTXBvhK7/doAMFRCZBQe0DUQ+MqCE3qv7khUtZu+0Aj27U9qcihZowqI0x840xjxhjXjXGvGKM+XQpCpPiKfXpLkfz/jULWDCjmq/8doNmgIgUqJARdRr4G2vtSuAM4OPGmJXOliXFVMrzEicSCQX427ctZ8PePh54cbfb5YhMCxMGtbV2j7X2ufznfcB6YK7ThUnxeOFi4khXnjSb5bPquP2xLepVixRgUj1qY8wiYDWwdozbbjTGdBpjOru71X/0Eq9cTBwSCBiuP6edDXv7eGKTdtUTmUjBQyxjTC3wC+Az1tro6NuttXcCdwJ0dHRomOQhffE0AQM1JTw0YCLvOGUO//r71/juE1s5N38QrsAfXt3HLQ9tJJXJEg4GCAcNFaEAleEg1ZEg73vLfC5cMcvtMqXEChpRG2PC5EL6x9baXzpbkhRbXzxV8kMDJlIRCvLhty7iiU09rN9zxPu+76QzWb7yuw1c/4NOMtksy2bVMq+pihk1EUKBALFEmr/sPMRH7+nkunueYfv+AbdLlhKacERtcj/ddwPrrbVfdb4kKbbcPh/eaHuM9FenL+CbD2/mrie2cct7T3a7HNds7urni796iT9vPcD71yzgS29fSWX4yN9+kuks9zy1jVv/sImLv/YYn3vbcq47u91Tb8DijEJG1GcB1wIXGmNeyH9c7nBdUkRu7UU9kcbqCFd3zOOBF3exLxp3u5ySe/aNA9zwg04u/upjvLjjELdcfTL/+7+eNGZIQ27GzI3nLuHhvz2f85e18OXfrOemHz7LocFUiSuXUitk1seT1lpjrV1lrT0l//FgKYqT4uiLp6j34Iga4Lqz20lnra/mVWeyln/41ctc9Z2neeb1A3zqoqU8+XcXcNVp8wr687PqK7nj2tP4hytX8vCGLq687QmeeV17fZcz7w2zpOj64mlmN1S6XcaYFs6s4ZMXHMc3Ht4MwL++ZxWhYPkumB1MZvjkvc/zh/X7uOGcdm6+ZBnVkcn/GBpjuO7sdlYvaOSTP3meq29/mitXzebzl61gXlO1A5WLmxTUPhCNp1g2q9btMsZ18yXLCAcD3PLQRuLpDF9/32oiofIL6/39Ca77ficv7uzln955Ah86c9ExP+apC5p46LPncsdjW7nj8S089Oo+bjhnMR87fwm1FfrxLhfl99Mgb5LNWrqiCWZ5dEQNudHhJy9ayhevOJ4HX9rLZ376fNkthNnVO8jVtz/N+j1Rbv/gaUUJ6SHVkRA3X7KMh//mfC49sY1vPrKZ8//Po9y3bjsZn7STyp2Cusz1xBIkM1nmNla5XcqErj9nMX9/2QoefGkv33pks9vlFM2mfX1c9e2n6O5P8KPrT+dtJ7Q58jxzGqu49ZrV/OrjZ7FoZjWf/+VLXPGNJ3hSi4qmPQV1mdvdm5tNMafB+0ENcOO5i3nXKXO45aGN/HH9PrfLOWYv7zrE1Xc8TcZafnbTmbxl0QzHn/OU+Y38/GNn8q0PnEosmeaDd6/lunueYVtPzPHnFmcoqMvc7t5BIDfamg6MMfzLVas4YU49n7nvBbZ097td0pQdGkxx0w+fpSYS4v6Pncnxs+tL9tzGGK5YNZuHbj6Pz1+2gnXbDvDub/9Ji4umKQV1mdt1MBfUc5umR1ADVIaD3HFtB5FQgI/e8wxdfdNvjrW1li/8x0vsjcb55gdWs3BmjSt1VIaDfOy8JfzmU+dQFQ7yV3etZeO+PldqkalTUJe5Xb2D1FaEPLEX9WTMbaziu3/dQXdfgg/dvY5DA9NrUcf9z+7kN3/Zw2cvWcbqBU1ul8OCmdX85IYzCAUMH/ju2mn9m4ofKajL3O7eQeY0Vk7LZcanLmjizms72Nod48P3rCOWSLtdUkG29cT40gOvcMbiGXzsvCVulzOsvbmGn9xwOmB517f+xI/XvuGbRUbTnYK6zO0+NDht+tNjOXtpM994/2pe3NHLtXev9fwFsY37+rj27rWEgwG+9r5TCAa89QZ5XGsd93/srZw4p4H/+R8vc/UdT7Nhr/rWXqegLnO7e+PTYmre0Vx6Yhu3vf9UNnX1c+nXH+fOx7eQzmTdLusIT2zq5qpvP0UineWH161htkdn2izKj6z/7eqT2drdz+W3PsHNP32B1z3+Juhn06txKZMykExzIJac1iPqIVesmk3Hoia++KuX+ecHN3Dfuh38lxPauGB5C6ctbHJ12fmOAwP8/NmdfOuRzSxtreV7H36L519zYwzvOW0eF61o5fbHtvD9p1/ngRd3896OefzdpStorI64XaKMYJxYAdbR0WE7OzuL/rgyOZu7+rn4q4/x9fedwrtWl8fpadZafvvyXn705zdYt+0A6aylubaCG85p54NnLKSmRMumE+kMv3xuFz/v3MFz23sBeNsJs7jlvadMy6XbXX1xvv3IFn745zdoqo7w5XedyKUnOrMwR8ZmjHnWWtsx5m0K6vL1+MZuPvS9dfzspjNZ0+78QotS64uneHJTDz9Zt50nNvXQVB3mpvOWcMM5ix3rDcdTGX7WuYPvPLqFPYfiLJ9Vx7tWz+XtJ88ui82QXt51iM/d/xde3RPl8pPa+Pylx7Ng5vT/e00HRwvq6ffWLwU7vNjFu/t8HIu6yjCXnTSby06azXPbD3LbHzfxL7/dwIs7evna+04Zd1/nqdjVO8i9a7dz3zPb6elPctrCJr5y1SrOWdo8LWfUjOfEuQ3830+cxZ2Pb+W2hzfx0Kv7+MCaBXziwqW01FW4XV7JxVMZegdSHBxIsrt3kG09Mbb1xIjG04QDhmDAUFMRYkZNhJm1EVrrKrlkZfGPSlNQl7HdvYMEDLTVl2dQj3Tqgib+/SNruOuJrXz5N+vp7lvLdz/UQVPN1Hut2azlT1t6+OHTb/CH9fuwwEUrWvnoWe2cuWRmWQX0SOFggI9fcBzvOW0et/5xEz9au5171+3glPmNvKW9iY6FM1jcUsOcxirC02xL2oFkmoc3dPHoa92EAoYZNRGaqiMkM1migykODabYG42z6+Agu3sHiSUzRzxGY3WYpuoIqUyWTNbSn0jTF89NHW2tq1BQy+Ts7B2krb6yrPd3Hu36cxbT1lDJZ3/6Ipfd+gTtzTVUhgNUR0I010Zora9kRk2EdCZLPJVlMJUhk7Vkbe4jaAyBgCGdsfzmpT1s64kxsybCTect4QNrFjB/hn/aALPqK/nnd5/EDecs5t5121m77QC3P7aVTHYLAAGTW/F69nHNXLJyFm9d0lzU32Kmqqc/QU9/ggOxJAdiSfZFE+yLxtnWE+OJTd3EU1maqsOEggEOxpKk83PJK0IB6irDzKqvoL25hrOXNtNaV0lDVZjG6jCz6itZ3Fwz5pt/Ip3hYCxFf8KZhVkK6jKWW+zi7dkHTrhy1Rxa6yr5zqOb6U+k6elPE0sO0NOXIBofe9GMMRAwhqy1DF22OW1hE5+5eCmXnthGRcj9AHJLe3MNX7j8eABiiTQv7zrE9gMD7DgwwMZ9/Tzwwm7uXbeD6kiQVfMaWDWvkZPmNrBoZg1tDZXMrIkQKPI1A2stWcvwtYhM1vLQq/v43pPbWDfGaTeRUIA5DZVcfdp8Lj9pNmvaZxAMGKy19CXSRIKBY3qTqQgFaWsIAs789qqgLmO7e+OcMr/R7TJcsaZ9Bmva1xzx/Xgqw8GBJKFAgKpIkIpQgFDAvKmNkc2PsP30m0ihaipCnL54Jqcvnjn8vUQ6w1Nb9vPohi5e2HmIe/70OskR89wjwQCLW2o4YU4Dx8+uo74qTCqTJZXOEg4FqK8MU1cZIpnO0tWXoCsa5+BAilgiTSyZJpHOkrW5cB5MZnL36YuTSGeZUR2hpa6CvniaXb2DzG2s4n+8bTntzTU0VUeYURNhVn0FDVXhMVtVxhjPHlM3koK6TGWzlj2HBrli1Wy3S/GUynBwwoUogYAhQHn2n51QEQpywfJWLljeCuROS9/U1cfOg4PsPRRnV+8gG/f18fimbn7x3M4JHy9goKEqTE1FiJpIiIpwgIAxBEzuuVYvaKSltoLqSJCeWJKuaAJrLV+84nguWTmrLN9gFdRlqqc/QSpjfdn6EHdFQgFOmNPACXMajritpz/BYDJDJP+bTCpj6YuniMZThIMBZtXnWiXlGLbHQkFdpnbmp+bNLdOpeTI9NdceOcWvzcPHxHmF3rbK1HQ7MEBExqegLlO7h0fUCmqR6U5BXaZ298apqwxRNw2uaIvI0Smoy9TQVCURmf4U1GVq10F/LnYRKUcK6jKUzmTZcXBAI2qRMqGgLkNPbu6hL57mrONmTnxnEfE8BXUZ+sVzu2isDnPBila3SxGRIlBQl5loPMV/vrKXd5w8x9cbCYmUEwV1mfnNX/aQSGe56tR5bpciIkWioC4zv3h2J0taalg178h9FkRkelJQl5E39sfofOMgV502r2xPHxHxo4KC2hhzqTHmNWPMZmPM550uSqbmF8/twhh4d5mcOC4iORPunmeMCQLfAi4BdgLPGGMesNa+6nRxMjFrLVu6+3lqy37uXbeds49rnnC/ZRGZXgrZ5nQNsNlauxXAGHMf8E6g6EH99tueJJ468jBJeTMLpDJZkukssUR6+HipeU1VfPqipe4WJyJFV0hQzwV2jPh6J3D66DsZY24EbgRYsGDBlIpZ0lLzpiN8ZHzhYGD4nLcT5tRz1nHNvjp4VcRPinZwgLX2TuBOgI6ODjuVx/j6NauLVY6ISNko5GLiLmD+iK/n5b8nIiIlUEhQPwMsNca0G2MiwDXAA86WJSIiQyZsfVhr08aYTwC/B4LA96y1rzhemYiIAAX2qK21DwIPOlyLiIiMQSsTRUQ8TkEtIuJxCmoREY9TUIuIeJyxdkprU47+oMZ0A29M8o81Az1FL8Y5qtdZqtdZqtdZU6l3obW2ZawbHAnqqTDGdFprO9yuo1Cq11mq11mq11nFrletDxERj1NQi4h4nJeC+k63C5gk1ess1ess1eusotbrmR61iIiMzUsjahERGYOCWkTE40oe1BMdlGuMqTDG/DR/+1pjzKJS1ziqnonq/bAxptsY80L+43o36szX8j1jTJcx5uVxbjfGmG/k/y5/McacWuoaR9UzUb3nG2MOjXht/1epaxxVz3xjzCPGmFeNMa8YYz49xn088xoXWK9nXmNjTKUxZp0x5sV8vf84xn08kQ8F1lq8bLDWluyD3DapW4DFQAR4EVg56j7/Hbg9//k1wE9LWeMU6v0w8E23ahxVy7nAqcDL49x+OfBbwABnAGs9Xu/5wK/dfl1H1DMbODX/eR2wcYz/D555jQus1zOvcf41q81/HgbWAmeMuo8n8qHAWouWDaUeUQ8flGutTQJDB+WO9E7g+/nP7wcuMsaYEtY4UiH1eoa19nHgwFHu8k7gBzbnz0CjMWZ2aao7UgH1eoq1do+19rn8533AenJnio7kmde4wHo9I/+a9ee/DOc/Rs928EQ+FFhr0ZQ6qMc6KHf0f5zh+1hr08AhYGZJqjtSIfUCXJX/Nfd+Y8z8MW73ikL/Pl5yZv7Xy98aY05wu5gh+V+5V5MbSY3kydf4KPWCh15jY0zQGPMC0AU8ZK0d9/V1Ox8KqBWKlA26mHjs/h+wyFq7CniIw+/2cuyeI7f/wcnAbcCv3C0nxxhTC/wC+Iy1Nup2PROZoF5PvcbW2oy19hRyZ7OuMcac6GY9R1NArUXLhlIHdSEH5Q7fxxgTAhqA/SWp7kgT1mut3W+tTeS/vAs4rUS1TcW0OqjYWhsd+vXS5k4ZChtjmt2syRgTJhd6P7bW/nKMu3jqNZ6oXi++xvlaeoFHgEtH3eSlfADGr7WY2VDqoC7koNwHgL/Of/4e4GGb78y7YMJ6R/Uf30GuD+hVDwAfys9MOAM4ZK3d43ZR4zHGtA31H40xa8j9f3XthzJfy93AemvtV8e5m2de40Lq9dJrbIxpMcY05j+vAi4BNoy6myfyoZBai5kNBZ2ZWCx2nINyjTH/BHRaax8g9x/rh8aYzeQuNF1TyhqnUO+njDHvANL5ej/sVr3GmHvJXcVvNsbsBL5E7iIH1trbyZ17eTmwGRgAPuJOpTkF1Pse4L8ZY9LAIHCNi2/aAGcB1wIv5XuTAF8AFoAnX+NC6vXSazwb+L4xJkjuDeNn1tpfezQfCqm1aNmgJeQiIh6ni4kiIh6noBYR8TgFtYiIxymoRUQ8TkEtIuJxCmoREY9TUIuIeNz/Bxdskh9+wq1iAAAAAElFTkSuQmCC\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "import matplotlib.pyplot as plt\n", "import numpy as np \n", "\n", "rdf = model.avg_rdf.result().numpy()\n", "plt.plot(rdf[1, :], rdf[0, :] / rdf[0,-1])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can also examine the final positions" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "[]" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYYAAAD6CAYAAAClF+DrAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjMuMCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy86wFpkAAAACXBIWXMAAAsTAAALEwEAmpwYAAAlvklEQVR4nO2df7BdV3XfP0sSpjMOBWODfwvj4npCMoVgjTAtSXFthNHQuMlAcGCI+ZEqUNyWEiZAPMMw5h9I6rR0LDCKYQqpwaYJDi51sOTGGdIZRP2eBwPGNhbGGksYWwZhTEkDQqt/3PPM1eXed+97Z+29195nfWbevPvuPe+effZZe33XWvucfURVCYIgCIIVNpRuQBAEQeCLEIYgCILgGEIYgiAIgmMIYQiCIAiOIYQhCIIgOIYQhiAIguAYTIRBRD4mIo+IyNfG3nu6iOwRkfu63yfM+N/Lum3uE5HLLNoTBEEQrB+xuI9BRH4N+CHwCVX95e69PwK+p6rvF5F3ASeo6jsn/u/pwBKwBVBgGThPVQ+vtr+TTjpJzzrrrN7tDoIgGBLLy8uPquoz5m23yWJnqvoFETlr4u1LgJd0rz8O/A3wzoltXgbsUdXvAYjIHuBi4FOr7e+ss85iaWmpX6ODIAgGhojsX2S7lHMMJ6vqQ93r7wAnT9nmdODBsb8PdO8FQRAEhcgy+ayjelWvmpWI7BCRJRFZOnTokFHLgiAIgklSCsPDInIqQPf7kSnbHATOHPv7jO69n0NVd6nqFlXd8oxnzC2RBUEQBOskpTDcBKxcZXQZ8Nkp29wCbBORE7qrlrZ17wVBEASFsLpc9VPAF4FzReSAiLwJeD/wUhG5D7io+xsR2SIi1wJ0k87vA27vfq5cmYgOgiAIymByuWputmzZonFVUhAEwdoQkWVV3TJvu7jzOQgGyPL+w+y8bR/L+1e9Zaja/QX9MLmPIfDD8v7D7L3/u5x/9omc96ypN5sHA2d5/2Fee+1efnzkKMdt2sB1v3t+UlvJvb+gPyEMDREDMFiEvfd/lx8fOcpRhZ8cOcre+7+b1E5y7y/oT5SSGmLaAAzapE9p5vyzT+S4TRvYKPCkTRs4/+wTE7TQfn8W5ahaS1q52x0ZQ2ZSlnpWBuBPjhzNMuCDMvTNDM971glc97vnZys5WuzPIhsukVFbjPcS7Q5hmCCl4059gnMP+KAMFqWZ8551Qlb76Ls/i2POXdKyGu8lSnEhDGOkdtw5TnDuAW9BTJivjSFmhhbHnLvfrMZ7ifMdwjBGasfd2oCuNU0uTd9+G2JmaHHMufvNaryXON8hDGOkdtwtDWhPaXJNGYdVv5XKDEv2tcUx5+w3y/Ge+3yHMIyRw3HXWOqZhpc0ubaMo+ZLN2vraw/UOt5DGCao9UTmxkuaXJujrbmcWFtfB+snhCFYF17S5Nocbc3lxNr62pKaypUWxCJ6TklhiK0at9UkeIt9Y80Q+6mlEtqii+hFxuCQFIbYknFP0rf813LfWOO91JpCuIZYQoslMQywvl09xdIWsVzGbKJvfNB3HK0I/FW77+W11+41G4+5lxDxQGQMPUkRbaao5Q65PjyP6JvyWIyjVJF9zfNC6yWEoScpjDGFIXo27tJ1a899MxQsxlFKgU9RQitt96uRVBhE5FzghrG3zgbeo6r/eWyblzB6HvS3urc+o6pXpmhPihORyhhTGKLH+rCX+r7HvlnBswOxwmIc1STwXux+FkmFQVXvBZ4PICIbgYPAjVM2/VtVfUXKtqQ6ETUZo0eGOLG3Frw7ECusxpFngR/Hu93nLCVdCHxTVfdn3OcTpDwRtRijR6K+vzreHYglQxpH3u0+pzBcCnxqxmcvEpE7gW8D71DVu6x37v1EDBXLjKvFkouV3bbYNzXjvdKQ5QY3ETmOkdP/JVV9eOKzfwgcVdUfish24IOqes6U79gB7ADYvHnzefv3rz3xiMHRLi2XXPrabar7YloZSy0dyzy83eD2cuCOSVEAUNUfjL2+WUQ+JCInqeqjE9vtAnbB6M7n9TRiSKlqTjwMrJZLLn3t1rpvWhLhlo7Fklw3uP02M8pIInKKiEj3emvXprjDaEFKP8M21U1Fa8XyJqTSfWqN9Q1aLd0Q2NKxWJI8YxCR44GXAr839t6bAVT1GuCVwFtE5Ajwd8ClWuMCTgtgHVl7iHa8ROpWNdtUfVr6OQaW9eyW5utaOhZLkguDqv5f4MSJ964Ze301cHXqdpQmhcPx4JQ9DSyLUmGKPvUg4JZlVO8Tp2uhpWOxJO58zkQKh+PBKbc2sFL0qQcBt8ZKaDzMT8Xc488TwjCFWu6Qtiyf9H2WbisDK4XQeRBwj3jIpKxILXC5BTSEYYLa7pBuZclpD5HjCtZC50XAvdFKJpV6DJUYoyEMEwztDmkPg9OLOKWkFQG3pJVMKvUYKjFGQxgmaMVYF8XD8XoQJ++02EetzE+lHkMlxmgIwwStGOuieDheD+LknVJ9lLp85TGLXiupx1CJMRrPfA5cYOmAWqvFr5D7uGopX7V6vlPgbUmMIFgVy8sfa3Bm6yF3dF1D+arl812SeOZz0BSxxIEdNTzrOM53GiJjcEikxusn5ivs8DD/NA/r8x1jb0TMMRhiYVSWqfFQjXyoxz1UrM53zrJUKRuNOYbMWBmVVV13yAJjWYuv7diHiNX5zjWnUsO8SAiDEVZGZZUaexOYGh1sjgFcY7/0wfPx5ipD1jCpPzhhSGWYVkZlVdf1JDA1REjTSD2Ahya63u0g15xKDfNggxKGlIZpaVQWqbEngakhQppG6gE8NNGtwQ5yXBJcw6T+oIQhtWF6u4vTi8DUECFNI/UAHpro1moHKfDmKyYZlDCEYa6PvkZcQ4Q0i5QDeGiiW7MdDI3BXa5aSz02CBYlbDpYFDeXq4rIA8DjwE+BI5ONEhEBPghsB34EvF5V70jVHu8pXBCslbDpwJpcpaQLVPXRGZ+9HDin+3kh8OHudzCFiA6PJfojH6X6Os5xfjzMMVwCfEJHNa29IvI0ETlVVR8q3TBv1HQFSg6iP/JRqq/jHJchxyJ6CuwWkWUR2THl89OBB8f+PtC9dwwiskNElkRk6dChQ4ma6ptUC4Yt7z/Mztv2sbz/sMn35SIWUMtHqb623m+ttp6bHBnDi1X1oIg8E9gjIveo6hfW+iWqugvYBaPJZ+tG1kCKK1BqjshquiKndkr1teV+Pdq61zJZcmFQ1YPd70dE5EZgKzAuDAeBM8f+PqN7L5ggxeV+NV0HP0mrlz+mchZ9vrdUX1vu19LWvS2YaU1SYRCR44ENqvp493obcOXEZjcBl4vI9YwmnR/zMr/gUc2tr0CpPer2tmBe3+9I5SwsvrfU1U9W+7WydW8LZqYgdcZwMnDj6IpUNgGfVNXPi8ibAVT1GuBmRpeq7mN0ueobErdpITyruSWtRt1rxeJ8W3xHKmfh2QnlwsrWvS2YmYKkwqCq9wPPm/L+NWOvFXhrynash5QDyVsmEtfB25xvi+9I5Sw8O6GcWNi6twUzU+DhclWXpBpIQ8lEasPifFt8RypnUdoJeQuG+uBtwcwUDG5JjLWQwph33raPq3bfy1GFjQJv33Yub73gOSbfXQOeHYSHOYYWiWDID26WxKiZFGqeI6X36py8Owir1Wg9HZMHYn6jPkIYMpM6pffsfC2fKudR+ILpxPxGfYQwFCBlVOk5OrNwEJ6FL5hO6fmNYO2EMDRGyknzvgPbwkF4Fr5gNlFiq4vBCUPrZYgU0ZlllN7XQURZIrCkdX+wXgYlDJYOzrNBWUdnnqL0KEvMJoVNerbzvkRZcjaDEgbLyc8hGZS3KN1K+FpyeilssnU79xTweGNQwmDl4IZmUC1G6V6cnpU4pbDJ1u3ceuXWlsbHoITBysF5i6Bz0NrkoQenZylOKWyydTu38gdeggxLBiUMYHcTU+oIurUIxBsenJ6lOKWwSY+ZovW4sPAHHoIMawYnDFakjKBbjEC84cHpWYtTCpv0lCl6HRceggxrQhgckioCyZmF1JDxlHZ6HsSpJrxG5pYlKS+2EMLgkNof4ek1svNIaXHKhYXT8xyZ9z2P3sZMCINDUkSSOaMtr5FdUAYrp9dyhuVtzIQwOMU6kswZbQ3lMkDPbfOE9SR7i33tLRtKJgwicibwCUaP91Rgl6p+cGKblwCfBb7VvfUZVZ18JnRgQM5oawiXAXpumze8Ob0VPAm7t2woZcZwBPh9Vb1DRJ4CLIvIHlX9+sR2f6uqr0jYjqAjZ7TV+mWAntvmDW9OD3wKu6dsKJkwqOpDwEPd68dF5G7gdGBSGLLhKUII5uM10gTfbVuE3GPBk9ODEPZ5ZJljEJGzgF8BvjTl4xeJyJ3At4F3qOpdM75jB7ADYPPmzWtug8cIIVgdj5HmCp7bNo8YC/ULe2qSC4OI/ALwF8DbVPUHEx/fATxLVX8oItuBvwTOmfY9qroL2AWjZz6vtR0RIdSJt0hznBRtyxHJx1hItzx9jYHCNJIKg4g8iZEoXKeqn5n8fFwoVPVmEfmQiJykqo9at6WmCKElAwsWJ1ckX9NYSImlsLeWhaW8KkmAjwJ3q+qfzNjmFOBhVVUR2QpsAL6boj21pP6tGViwOLki+VrGQk20loWlzBj+GfA64Ksi8uXuvT8ENgOo6jXAK4G3iMgR4O+AS1V1zWWiRfFcllihNQMLFidnJF/DWKiJ1rKwlFcl/W9A5mxzNXB1qjZYk6PEk8vASperSu/fY3sikq+X1s6dJAzQk7FlyxZdWlrKus/caw2lXtK7ZLmq9P69tydYGx5E3YIcxyEiy6q6Zd52sSTGguQs8aRO80uXq0rv33t71korjnE9tCLq3o5jQ7E9V8ZKiWejUH0NsfSxWOx/ef9hdt62j+X9h120pxQrDuWq3ffy2mv3mvRHTUwT9RrxdhyRMSxISzXE0sfSd//W0VXp/uhD7mzHMjtpfSnuteDtOEIY1kBLV3KUPpY++0/hDEv3x3rJ6VAsBdnjUtwlS3LegpNBC8OQa7M14y26mkdKO8vpUCwF2dtS3B5q/J6Ck8EKgwdD6MtQhc1bdLUaOewsl0OxFGRv4l77BQjWDFYYajcEL8JWSpw8RVerUbudjWMpyN7E3Vqoag/aBisM3iKWteLB4XgRJ8/UbmeTWAqyJ3G3nquofVwMVhi8RSxrxYPD8SBO3qndzoaElVC1MC4GKwzQzxBKp4oeHI4HcaoBT5ExlLfd1mlhXMSSGOughVRxhb5OIpxMHqz6uSXb9YzXcRFLYiSkhVQRbJyEt2h4Fl4H6iJYOvNWbHcepc93LeNiFiEM66CFVBGG5SRqjpItz1MrtrsatZ9vKC9sIQzrwEN934IhOAmoXwAtz5Nn27VyhrWfbw/CFsKwTmpPFcG3k7CkdgG0Pk8ebdfSGdZ+vj0IWwjDwPHoJKwpJYCW5QBv58m61GG9REbNAY8HYUsuDCJyMfBBYCNwraq+f+LzJwOfAM5j9LznV6vqA6nbVTOl648W5D6G3I7VQzkgFSmOzdoZehPSteBB2JIKg4hsBHYCLwUOALeLyE2q+vWxzd4EHFbV54jIpcAHgFenbFdJLC4Prd3htHAM8/BQDkhFqtVtSztDT5QWttQZw1Zgn6reDyAi1wOXAOPCcAnw3u71nwNXi4hojTdYzMHCIbbgcFo4hnl4KAekItWxlXaGs2ghQ18rqYXhdODBsb8PAC+ctY2qHhGRx4ATgUcTt20qKY3AwiG24HBaOIZ5tBwBt3xskwwhu51GNZPPIrID2AGwefPmJPtIbQQWDrHkoLQSzaE4Fq8RsAUtH9s4Q8hup5FaGA4CZ479fUb33rRtDojIJuCpjCahj0FVdwG7YLQkRorGpjYCK4dYYlCmeJzmEAaYNUMsa5RkCNntNFILw+3AOSLybEYCcCnwmoltbgIuA74IvBL461LzCzmMoFaHOJTIybPjHWpZY71YnMuhZLeTJBWGbs7gcuAWRperfkxV7xKRK4ElVb0J+CjwZyKyD/geI/EoQstG0HeQDCFy8u54hyLOFliey1qDuT4kn2NQ1ZuBmyfee8/Y6/8HvCp1OxalRSOwWiyvtGimjua9O94hiLMV3s+ld6qZfA7Wj9UgKSmaOaJ5747XWpy9l80iwy13jkIYBkALgyRHBOghK5qHlTh7Lpu1kuH2peQ5CmEwxlLh4/LQn5FL3FosJU4jhdB6Wx219nNZshwWwmCIpcLH5aHH0oK4gZ/yjbXQxuqo9pTshxAGQywVPibPfp7axc1T+cZaaGN1VHtK9kMIgyGWCt9S1OQlSi6NN7G3FNohro6aw65L9UMIgyGWCt9K1OQpSi5NS2I/SSv2uiit23UIgzGWCu8hauobFXmLkkvSuvP0YK+5aN2uQxiCmVhERS1HyethSM6zZVq36xCGYCYWUVErUbK3eRJv7SlNiScCtmDXswhhCGZiFRWVvmO67+D1Vk/21p7SlOqPlrO/EIaG6esUa4+KrByGt3qydXtqzz68nZ8WGJQw1D4A1oKVU6w5KrJyGN7qyZbt8ZJ99Bmb3s7PPGrwQ4MRBi8DIBcRRdmWwjxlTpbt8WAnfcemt/OzGrX4ocEIg4cBkJPaoqgUWN9X4slerNrjwU6sLnLwdH5mUYsfGowweBgAixAL5/0MqydweTh2r+UDCzuJJbIXp5ZjlUJP0ezFli1bdGlpac3/ZzU4Uw3yWtLMHLTUFy0dyyRWx+ZVOFNQ8lhFZFlVt8zbLknGICJ/DPxL4MfAN4E3qOr3p2z3APA48FPgyCIN7oNF9JhykNeSZuagpb5o6VgmiSWy104Nx7oh0ffuAX5ZVf8J8A3g3atse4GqPj+1KFgxbSBYsZJmbhRcp5nzWN5/mJ237WN5/+F1f0crfQFtHcskLR/bkEmSMajq7rE/9wKvTLGfEqSsEbYyL2B1mWztfbFCS8cySUvHNqRy1jySzzGIyP8AblDV/zbls28BhwEFPqKquxb5zvXOMVjRugH1Ob6dt+3jqt33clRho8Dbt53LWy94TqKWBoENLc8DjZN8jkFEbgVOmfLRFar62W6bK4AjwHUzvubFqnpQRJ4J7BGRe1T1CzP2twPYAbB58+b1NtuEGmqE66XvAPFw1UXrwr0IQ+kDb48TbYV1C4OqXrTa5yLyeuAVwIU6Iy1R1YPd70dE5EZgKzBVGLpsYheMMob1tjtYnb4DpHRpIVXkV5OjHUr02/rjREvaXKqrki4G/gD456r6oxnbHA9sUNXHu9fbgCtTtGcoWBiSxQApmVGlesh9TY52KNFvy48TLW1zqW5wuxp4MqPyEMBeVX2ziJwGXKuq24GTgRu7zzcBn1TVzydqT/PEpO+IFJFfbY7WY/SbgpYfJ1ra5lJdlTR1tlFVvw1s717fDzwvxf49kjottI6evAyQtZJC2GpztLWL+6K0fJylbW5Qdz4vQgoHniMtXNnHiiF5L3fURk1zDBZ4XyVgCKTou6J3PtdKKgeeIy1sIXry7ERqzqLWiuUyFzXNzVhj8TyUUv0VwjBGKgeeKy2s2XkN3Yl4wmoclK6Tl6R2ew5hGCOVA28hmk/NkJ3ILEplUFbjoHSdvCS123MIwxgpHXjN0XwOhuxEplEy4rQaB0MOiGq355h8DnoTE5X2xNIi9ePRnmPyOciCZWSbK6vyOGAnqT3itKLkuap58rgvIQwV4NmR1VZLrWVScMhlmBVKnqvSdlJ6zIcwOCelgXpZQiMnNQlZzRGnBSXPVcl9lxYlCGFwTyoDHeoSGrUJ2ZApea5K7ttD8BLCYIx1CpjKQIe6hEZtQjZkSp6rkvv2ELzEVUmG1LTkcwtLaJSuw66HGtsc5CeVncRVSQVIlQKmiMhrj5w91GHXSi2BQ4hXeUpn4SEMhnhIAddCKeOzcDwe6rBrpYZnRdQouIE9IQzYRUi1R+E5sHI8tYkw1PGsiBoFN7Bn8MJgHSGVTgG9Y+V4ahThGp4VUaPgTiPKYf0YvDBEhJQXS8dTowhbt9labGoU3EmiHNafZMIgIu8F/jVwqHvrD1X15inbXQx8ENjI6LGf70/Vpmm0EiHVQguOxxspxKbm81Iq2GspS0mdMfwnVf2Psz4UkY3ATuClwAHgdhG5SVW/nrhdT+DNUbVkXLPw4HiG0M81Uuvd+K1lKaVLSVuBfd3znxGR64FLgGzCAD4cFaQxLouB1poTbW0Qt0LNd+O3VpJOLQyXi8jvAEvA76vq4YnPTwceHPv7APDCxG1yi7VxWQy0Fp1oa4O4FWq+G7+1kvSGPv8sIreKyNem/FwCfBj4R8DzgYeAq3rua4eILInI0qFDh+b/Q4WsGNdGwcS4pg20Et/hDet+tmJ5/2F23raP5f2T8dMw8HpeFmElS3n7tnObCJ56ZQyqetEi24nInwKfm/LRQeDMsb/P6N6btq9dwC4YLYmxtpbWgXUKbBHFtBYJgb95Jag3M7MsM3o8L2vBS0nagmRrJYnIqar6UPf6PwAvVNVLJ7bZBHwDuJCRINwOvEZV71rtuz2uleS1Dh9zDHWQ+4ltVnZRo5gNGQ9rJf2RiDwfUOAB4Pe6hp3G6LLU7ap6REQuB25hdLnqx+aJgkc8DxCLKKZEJDQ0McqZmVnZa8zVtEsyYVDV1814/9vA9rG/bwZ+7v6G3PRxRDFAbPEstNOwELGcZRQre22xzOiF0oFR6ctVXdDXEcUAsaUmoa3xmddW9lr7nIBXPARGIQz0d0QxQGypSWhrErEVLO21pQlXL3iwqRAGbBxRDBA7ahLamkRsnLBXv3iwqXiCW0fpmt40PLYp+HniPLVH6XNa+gluIQxO8VBnHAKlHUDgj5bH3qLC0OvO5yAdLd5x7I0VB3DV7nt57bV7q7njuOU7pD0cW4y9mGNwS8o6Y0TJIzxM8q2VlqNZq7W9alyd1RshDE5JNQHbsmNZK6WWZ+5zTmsUs0Xpe2yeVmetPfgKYXBMiitHWnIsfQdf7qufLBxXy9Fs32PzsjprC8FXCMPAaMWxWEaHuQathePydCmvdVTc99i82Haq4CtnFhLCMDA8OZY+1Jj5WN5xXPpYU0XFfY7Ni22nEKjcWUgIwyrUXiechQfH0hcv0eFa8OK4LPAqzB5sO8V5zt3fIQwzaKFO2DIlnKzVYnktlH9qFOacWJ/n3P0dwjCDFAqdKgNpNbOZR04n20KgYL3gXyvZTw3k7u8QhhlYK3Qqx1Kjw6pRyLyWTtaC9TGUzn6GRs7+DmGYgbVCp3IstTmsGoUM2iidtHAMQR5CGFbBUqFTDcraBnttQrZCC6WTFo4hyEMsopeR2uYYUnzvSsawImQ5M4YaS1hBYEnR1VVF5Abg3O7PpwHfV9XnT9nuAeBx4KfAkUUaDPUKQ02kLPmUcNA1lbCGJGAej9Vjm6xYVBiSlJJU9dVjDbkKeGyVzS9Q1UdTtCNYPylLPiUmLWspYdUkYH3xeKwe21SCpMtui4gAvwV8KuV+AntW5i42ClXMXcyjluPxuORzqqWwPR6rxzaVIPXk868CD6vqfTM+V2C3iCjwEVXdNeuLRGQHsANg8+bN5g0NjsXjRGWfFN/j8UzD28UEKSNob8fqtU0lWPccg4jcCpwy5aMrVPWz3TYfBvap6lUzvuN0VT0oIs8E9gD/VlW/MG/fFnMMLdcRW2RIKb4n29x52z6u2n0vRxU2Crx927m89YLnmH2/p2NdwWObrEg+x6CqF81pwCbgN4HzVvmOg93vR0TkRmArMFcY+mLtZFo2JC/UMkdggdUcTA0PrfF4k5zHNuUmZSnpIuAeVT0w7UMROR7YoKqPd6+3AVcmbM8TWDqZIUWyJSmZ4tco/J4eWhPUR0phuJSJSWcROQ24VlW3AycDN47mp9kEfFJVP5+wPU9g6WSsRKZG55OTUg6qVuH38tCaIdLCWE4mDKr6+invfRvY3r2+H3heqv2vhqWTsRCZWp1PbuIy18UZyiSqNyfcylge7JIYVk7GQmRyOp8cA8nbYO1DisUUc/TNEEpAHp1wrYHEJIMVBkv6ikyu6C7HQPI2WD09Fzp337ReAvLohFvJ1EIYHJArussxkDwNVm/PhfbUNzlJlSV5dMKtZGohDE7IEd3lGEieBqs3R+ypb3KRMkuycsLWwtVCphbCMCByRDMW+2j18ZOtRJNrIbU493XC3kqfXghhGBg5opk++2j98ZMtRJNrwZs4T+Itq/RCCEPginj8ZFt4FOdxvAtXKUIYAlfEQF2MlJe9Dqnm7l24ShFPcDOgpev2PRD9uTqpH6IUNfd2KfqgniERA8kezxGmB1LWxb3X3CNoyMPghMHasFIMpDD+YDVSlts8l/JqD8JqGteDEoYUhpViyYSajT9IT8q6uOeau/dsZjVqG9eDEoYUhmU9kGo2fqgrKqqZlOU2r6U8z9nMPGob14MShlSGZTmQan7uQKmoKJUYWX9viGY/PGcz86hN1AYlDDUYVqk2Wjj1ElFRKjFK8ZS/mkoJXimRzVgIeg2+Z5xBCQP4TZPHKdFGC6deIipKJUbW31tbKSEYYX0nfi3nfHDCEEzHwqmXiIpSiZH199ZWSghGDFXQe93gJiKvAt4L/CKwVVWXxj57N/Am4KfAv1PVW6b8/7OB64ETgWXgdar643n79XaDWyvUWgOPOYYgFSsZw4qg114CXPQGt77C8IvAUeAjwDtWhEFEnsvoec9bgdOAW4F/rKo/nfj/TwOfUdXrReQa4E5V/fC8/YYwBEGQi5YEPcudz6p6d7ezyY8uAa5X1b8HviUi+xiJxBfHGijAvwBe0731cUbZx1xhCIIgyEVNcwNWbEj0vacDD479faB7b5wTge+r6pFVtnkCEdkhIksisnTo0CHTxgZBEAQ/Y27GICK3AqdM+egKVf2sfZOmo6q7gF0wKiXl2m8QBMHQmCsMqnrROr73IHDm2N9ndO+N813gaSKyqcsapm0TBEEQZCZVKekm4FIReXJ35dE5wP8Z30BHs963Aa/s3roMyJaBBEEQBNPpJQwi8hsicgB4EfA/ReQWAFW9C/g08HXg88BbV65IEpGbReS07iveCby9m5w+Efhon/YEQRAE/YkH9QRBEAyELPcxlEJEDgH7x946CXi0UHMWxXsbo3398d5G7+0D/22svX3PUtVnzPuSKoVhEhFZWkQFS+K9jdG+/nhvo/f2gf82DqV9qSafgyAIgkoJYQiCIAiOoRVh2FW6AQvgvY3Rvv54b6P39oH/Ng6ifU3MMQRBEAR2tJIxBEEQBEZUIwwi8ioRuUtEjorIlonP3i0i+0TkXhF52Yz/f7aIfKnb7gYROS5hW28QkS93Pw+IyJdnbPeAiHy12y7rjRki8l4ROTjWzu0ztru469d9IvKujO37YxG5R0S+IiI3isjTZmyXvQ/n9Ul3x/8N3edfEpGzcrSr2/eZInKbiHy9Gy//fso2LxGRx8bO/Xtyta/b/6rnTEb8l67/viIiL8jcvnPH+ubLIvIDEXnbxDZZ+1BEPiYij4jI18bee7qI7BGR+7rfU5eAFZHLum3uE5HLFtqhqlbxw+hhQOcCfwNsGXv/ucCdwJOBZwPfBDZO+f9PA5d2r68B3pKp3VcB75nx2QPASYX6872MnqGx2jYbu/48Gziu6+fnZmrfNmBT9/oDwAc89OEifQL8G+Ca7vWlwA0Z23cq8ILu9VOAb0xp30uAz5Wwu0XOGbAd+CtAgPOBLxVs60bgO4yu/y/Wh8CvAS8Avjb23h8B7+pev2vaGAGeDtzf/T6he33CvP1VkzGo6t2qeu+Uj5549oOqfgtYefbDE4w9++HPu7c+DvyrhM0d3+9vMXpoUY1sBfap6v06erLe9Yz6Ozmqult/tiT7XkaLLHpgkT65hJGNwcjmLuxsITmq+pCq3tG9fhy4m1WWs3fKJcAndMReRottnlqoLRcC31TV/XO3TIiqfgH43sTb43Y2y6e9DNijqt9T1cPAHuDiefurRhhWwfzZD4b8KvCwqt4343MFdovIsojsyNCeSS7vUvWPzUhDF+nbHLyRUQQ5jdx9uEifPLFNZ3OPMbLBrHQlrF8BvjTl4xeJyJ0i8lci8kt5Wzb3nHmxOxhlfLMCu5J9CHCyqj7Uvf4OcPKUbdbVl72e4GaNOHn2wyIs2NbfZvVs4cWqelBEngnsEZF7usggeRsZPSnvfYwG6fsYlbzeaLXvRVikD0XkCuAIcN2Mr0nah7UiIr8A/AXwNlX9wcTHdzAqjfywm1v6S0YrIOeiinPWzUP+OvDuKR+X7sNjUFUVEbNLTF0Jg1b07Id5bRWRTcBvAuet8h0Hu9+PiMiNjMoUZgNk0f4UkT8FPjflo0X6dt0s0IevB14BXKhdwXTKdyTtwyks0icr2xzo7OCpjGwwCyLyJEaicJ2qfmby83GhUNWbReRDInKSqmZZA2iBc5bU7tbAy4E7VPXhyQ9K92HHwyJyqqo+1JXaHpmyzUFG8yErnMFonnZVWigleX32w0XAPap6YNqHInK8iDxl5TWjydavTds2BRM129+Yse/bgXNkdEXXcYzS6psyte9i4A+AX1fVH83YpkQfLtInNzGyMRjZ3F/PEjZrurmMjwJ3q+qfzNjmlJU5DxHZysgPZBGuBc/ZTcDvdFcnnQ88NlYyycnMjL9kH44xbmezfNotwDYROaErF2/r3ludXLPqfX8YOa8DwN8DDwO3jH12BaMrRe4FXj72/s3Aad3rsxkJxj7gvwNPTtze/wq8eeK904Cbx9pzZ/dzF6PySc7+/DPgq8BXOgM7dbKN3d/bGV3Z8s2cbezO04PAl7ufaybbV6oPp/UJcCUjEQP4B52N7ets7uyM/fZiRuXBr4z13XbgzSv2CFze9dedjCb2/2nG9k09ZxPtE2Bn179fZewqxIztPJ6Ro3/q2HvF+pCRQD0E/KTzg29iNG/1v4D7gFuBp3fbbgGuHfvfN3a2uA94wyL7izufgyAIgmNooZQUBEEQGBLCEARBEBxDCEMQBEFwDCEMQRAEwTGEMARBEATHEMIQBEEQHEMIQxAEQXAMIQxBEATBMfx/NToOGuYgC8wAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "pos = tfcompute.get_positions_array()\n", "plt.plot(pos[:,0], pos[:,1], '.')" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "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.7.8" } }, "nbformat": 4, "nbformat_minor": 2 }