{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Z peak example" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We will plot a invariant mass distribution from a LHE input file.\n", "\n", "First, some basic imports." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "import hist\n", "\n", "import pylhe" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "# Use an example LHE file from package scikit-hep-testdata\n", "from skhep_testdata import data_path\n", "\n", "lhe_file = data_path(\"pylhe-drell-yan-ll-lhe.gz\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Prepare a histogram to calculate the invariant mass of two particles." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [], "source": [ "mass_hist = hist.Hist.new.Reg(30, 50, 150).Int64()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Use the generator provided by pylhe to read the events." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [ { "data": { "text/html": [ "\n", "
\n", "
\n", "\n", "\n", "\n", "50\n", "\n", "\n", "150\n", "\n", "\n", "Axis 0\n", "\n", "\n", "\n", "
\n", "
\n", "Regular(30, 50, 150, label='Axis 0')
\n", "
\n", "Int64() Σ=13189080.0 (16780000.0 with flow)\n", "\n", "
\n", "
\n", "" ], "text/plain": [ "Hist(Regular(30, 50, 150, label='Axis 0'), storage=Int64()) # Sum: 13189080.0 (16780000.0 with flow)" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "events = pylhe.to_awkward(pylhe.read_lhe_with_attributes(lhe_file))\n", "mass_hist.fill(\n", " (events.particles.vector[:, -1] + events.particles.vector[:, -2]).mass,\n", " weight=events.eventinfo.weight,\n", ")" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjoAAAG1CAYAAADwRl5QAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjcuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/bCgiHAAAACXBIWXMAAA9hAAAPYQGoP6dpAAAsyklEQVR4nO3df3RU9Z3/8dfkJ4khifwwIRpAKiJRCGgCatUKpkQUXLEtlKMQgWXX3WGpxvqDFQEpP+rSRfpjhBUP0FYpyBERfyzURF20RkyCgWIEyTZCVBJQJENCJGFyv3/4ZbYxQDKTSe7MZ56Pc+YcZ+5n7rznI8l95X4+n3sdlmVZAgAAMFCE3QUAAAB0FoIOAAAwFkEHAAAYi6ADAACMRdABAADGIugAAABjEXQAAICxCDoAAMBYBB0AAGAsgg4AADAWQQcAABgr5IPO/v37NWzYMO8jLi5OW7ZssbssAAAQBBwm3dSzrq5O/fv318GDB3XBBRe06z3Nzc364osv1L17dzkcjk6uEAAABIJlWTpx4oTS0tIUEXHu8zZRXVhTp9u6datuueWWdoccSfriiy+Unp7eiVUBAIDOUlVVpUsuueSc220POjt27NCyZctUWlqqw4cP66WXXtKdd97Zoo3L5dKyZctUXV2tzMxM/fa3v9WIESNa7euFF17Q1KlTffr87t27S/q2oxITE/3+HgAAoOu43W6lp6d7j+PnYnvQqa+vV2ZmpqZPn6677rqr1faNGzcqPz9fq1at0siRI7VixQrl5uZq//79uuiii7zt3G633nvvPW3YsOG8n3fq1CmdOnXK+/zEiROSpMTERIIOAAAhpq1pJ7ZPRh47dqwWLVqkCRMmnHX78uXLNXPmTE2bNk0ZGRlatWqV4uPjtWbNmhbtXn75ZY0ZM0bdunU77+ctXbpUSUlJ3gfDVgAAmMv2oHM+jY2NKi0tVU5Ojve1iIgI5eTkqKioqEXbF154QZMmTWpzn3PmzFFtba33UVVVFfC6AQBAcLB96Op8vvzyS3k8HqWkpLR4PSUlRfv27fM+r62t1QcffKAXX3yxzX3GxsYqNjY24LUCAIDgE9RBp72SkpJUU1NjdxkAACDIBPXQVa9evRQZGdkqxNTU1Cg1NdWmqgAAQKgI6qATExOja665RoWFhd7XmpubVVhYqOuuu65D+3a5XMrIyFB2dnZHywQAAEHK9qGruro6VVRUeJ9XVlaqrKxMPXr0UN++fZWfn6+8vDxlZWVpxIgRWrFiherr6zVt2rQOfa7T6ZTT6ZTb7VZSUlJHvwYAAAhCtgedkpISjRo1yvs8Pz9fkpSXl6d169Zp0qRJOnr0qObNm6fq6moNGzZM27ZtazVBGQAA4LuMuteVP86c0amtreWCgQAAhIj2Hr+Deo4OAABAR4Rt0GEyMgAA5mPoiqErAABCTnuP37ZPRgYQOizLUkOTp0P7iIuObPMmfAAQKAQdAO3W0ORRxrztHdpH+cJcxcfwqwdA1wjbOToAAMB8/FkFwC8lc3MUHxPZrrYnGz3KWlTQyRUBQGthG3RcLpdcLpc8no7NNwDCVXxMJENQAIJe2A5dOZ1OlZeXq7i42O5SAABAJwnboAMAAMxH0AEAAMYi6AAAAGMRdAAAgLHCNuhwrysAAMwXtkGHVVcAAJgvbIMOAAAwH0EHAAAYi6ADAACMRdABAADGIugAAABjEXQAAICxwjbocB0dAADMF7ZBh+voAABgvrANOgAAwHwEHQAAYCyCDgAAMBZBBwAAGIugAwAAjEXQAQAAxiLoAAAAYxF0AACAscI26HBlZAAAzBe2QYcrIwMAYL6wDToAAMB8BB0AAGAsgg4AADAWQQcAABiLoAMAAIxF0AEAAMYi6AAAAGMRdAAAgLEIOgAAwFgEHQAAYCyCDgAAMFbYBh1u6gkAgPnCNuhwU08AAMwXtkEHAACYj6ADAACMRdABAADGIugAAABjEXQAAICxCDoAAMBYBB0AAGAsgg4AADAWQQcAABgryu4CAHQty7LU0OTx670nG/17HwDYhaADhJmGJo8y5m23uwwA6BIMXQEAAGNxRgcIYyVzcxQfE+nXe+Oi/XsfAHQlgg4QxuJjIhUfw68BAOZi6AoAABiLoAMAAIwVtkHH5XIpIyND2dnZdpcCAAA6SdgGHafTqfLychUXF9tdCgAA6CRhG3QAAID5CDoAAMBYBB0AAGAsgg4AADAWQQcAABiLoAMAAIxF0AEAAMYi6AAAAGMRdAAAgLEIOgAAwFgEHQAAYCyCDgAAMBZBBwAAGIugAwAAjEXQAQAAxiLoAAAAY0XZXQCA8HKy0ePX++KiI+VwOAJcDQDTEXQAdKmsRQV+va98Ya7iY/iVBcA3DF0BAABj8ecRgE4XFx2p8oW5Pr/vZKPH7zNAACARdAB0AYfDwbATAFswdAUAAIxlRNCprKzUqFGjlJGRoSFDhqi+vt7ukgAAQBAw4lzyvffeq0WLFunGG2/UsWPHFBsba3dJAAAgCIR80Pnoo48UHR2tG2+8UZLUo0cPmysCAADBwvahqx07dmj8+PFKS0uTw+HQli1bWrVxuVzq37+/unXrppEjR+qDDz7wbjtw4IASEhI0fvx4XX311VqyZEkXVg8AAIKZ7UGnvr5emZmZcrlcZ92+ceNG5efna/78+dq1a5cyMzOVm5urI0eOSJJOnz6td955R08//bSKior0xhtv6I033ujKrwAAAIKU7UFn7NixWrRokSZMmHDW7cuXL9fMmTM1bdo0ZWRkaNWqVYqPj9eaNWskSRdffLGysrKUnp6u2NhY3XbbbSorKzvn5506dUput7vFAwAAmMn2oHM+jY2NKi0tVU5Ojve1iIgI5eTkqKioSJKUnZ2tI0eO6Ouvv1Zzc7N27NihwYMHn3OfS5cuVVJSkveRnp7e6d8DAADYI6iDzpdffimPx6OUlJQWr6ekpKi6ulqSFBUVpSVLluimm27S0KFDNXDgQI0bN+6c+5wzZ45qa2u9j6qqqk79DgAAwD4hv+pK+nb4a+zYse1qGxsby/JzAADCRFCf0enVq5ciIyNVU1PT4vWamhqlpqbaVBUAAAgVQR10YmJidM0116iwsND7WnNzswoLC3Xdddd1aN8ul0sZGRnKzs7uaJkAACBI2T50VVdXp4qKCu/zyspKlZWVqUePHurbt6/y8/OVl5enrKwsjRgxQitWrFB9fb2mTZvWoc91Op1yOp1yu91KSkrq6NcAAABByPagU1JSolGjRnmf5+fnS5Ly8vK0bt06TZo0SUePHtW8efNUXV2tYcOGadu2ba0mKAMAAHyX7UHn5ptvlmVZ520za9YszZo1q4sqAgAApgjqOToAAAAdQdABAADGCtugw6orAADMF7ZBx+l0qry8XMXFxXaXAgAAOknYBh0AAGA+gg4AADAWQQcAABgrbIMOk5EBADBf2AYdJiMDAGC+sA06AADAfAQdAABgLIIOAAAwFkEHAAAYi6ADAACMFbZBh+XlAACYL2yDDsvLAQAwX9gGHQAAYD6CDgAAMBZBBwAAGIugAwAAjEXQAQAAxiLoAAAAY0XZXYBdXC6XXC6XPB6P3aUAfrEsSw1Nvv/7PdnIv3kA4SNsg47T6ZTT6ZTb7VZSUpLd5QA+a2jyKGPedrvLAICgxtAVAAAwVtie0QFMUjI3R/ExkT6/Ly7a9/cAQCgh6AAGiI+JVHwMP84A8F0MXQEAAGMRdAAAgLEIOgAAwFgEHQAAYCyCDgAAMFbYBh2Xy6WMjAxlZ2fbXQoAAOgkYRt0nE6nysvLVVxcbHcpAACgk4Rt0AEAAOYj6AAAAGMRdAAAgLEIOgAAwFgEHQAAYCyCDgAAMBZBBwAAGIugAwAAjBVldwEA0B4nGz1+vzcuOlIOhyOA1QAIFQQdACEha1GB3+8tX5ir+Bh+3QHhKGyHrrjXFQAA5gvbP3GcTqecTqfcbreSkpLsLgfAWcRFR6p8Ya5f7z3Z6OnQWSAAZgjboAMg+DkcDoacAHRI2A5dAQAA8xF0AACAsQg6AADAWAQdAABgLIIOAAAwFkEHAAAYi3WbgM0sy1JDk++3N+jILREAIFwQdACbNTR5lDFvu91lAICR/Bq6GjBggL766qtWrx8/flwDBgzocFEAAACB4NcZnU8//VQeT+vT5qdOndLnn3/e4aKAcFUyN0fxMZE+vy8u2vf3AEA48CnobN261fvf27dvb3GPKI/Ho8LCQvXv3z9gxQHhJj4mklseAEAA+fQb9c4775T07f1n8vLyWmyLjo5W//799Z//+Z8BKw4AAKAjfAo6zc3NkqRLL71UxcXF6tWrV6cUBQAAEAh+nSOvrKwMdB0AAAAB5/dkgMLCQhUWFurIkSPeMz1nrFmzpsOFAQAAdJRfQeeJJ57QwoULlZWVpT59+sjhcAS6rk7ncrnkcrnOunoMAACYwa+gs2rVKq1bt05TpkwJdD1dxul0yul0yu12t1g9BgAAzOFX0GlsbNT1118f6FoAoFP4e7uMuOjIkDxjDeD/+BV0/vEf/1Hr16/X448/Huh6ACDgshYV+PW+8oW5XNcICHF+/QR/8803euaZZ1RQUKChQ4cqOjq6xfbly5cHpDgAAICO8Cvo7NmzR8OGDZMk7d27t8U2TvMCCAZx0ZEqX5jr8/tONnr8PgMEIPj4FXTeeuutQNcBAAHlcDgYdgLg393LAQAAQoFff+6MGjXqvENUb775pt8FAQAABIpfQefM/JwzmpqaVFZWpr1797a62ScAAIBd/Ao6Tz311FlfX7Bggerq6jpUEAAAQKAEdI7OPffcw32uAABA0Aho0CkqKlK3bt0CuUsAAAC/+TV0ddddd7V4blmWDh8+rJKSEq6WjLBlWZYamny/1YC/tycAALTNr6Dz3ZtgRkREaNCgQVq4cKHGjBkTkMKAUNPQ5FHGvO12lwEA+Dt+BZ21a9cGug4AAICA69BlQ0tLS/Xxxx9Lkq688koNHz48IEUBoa5kbo7iYyJ9fl9ctO/vAQCcm19B58iRI/rpT3+qt99+W8nJyZKk48ePa9SoUdqwYYN69+4dyBqBkBMfE8ntBwAgCPi16urf/u3fdOLECX300Uc6duyYjh07pr1798rtdmv27NmBrhEAAMAvfv3JuW3bNhUUFGjw4MHe1zIyMuRyuZiMDAAAgoZfZ3Sam5sVHR3d6vXo6Gg1Nzd3uCgAAIBA8OuMzujRo/Wzn/1Mf/rTn5SWliZJ+vzzz/XAAw/olltuCWiBpvP32ivfFRcded4brQIAEI78Cjq/+93vdMcdd6h///5KT0+XJFVVVemqq67Sc889F9ACTReoa6+UL8xl8isAAN/h15ExPT1du3btUkFBgfbt2ydJGjx4sHJycgJaHAAAQEf4FHTefPNNzZo1S++//74SExP1wx/+UD/84Q8lSbW1tbryyiu1atUq3XjjjZ1SrOl8vfbKyUaPshYVdGJFAACENp+CzooVKzRz5kwlJia22paUlKR//ud/1vLly7s86PTv31+JiYmKiIjQhRdeqLfeeqtLPz9QuPYKAACB5dOqq927d+vWW2895/YxY8aotLS0w0X547333lNZWVnIhhwAABB4Pp0+qKmpOeuycu/OoqJ09OjRDhcF2KUjq+C4CzkABB+fgs7FF1+svXv36rLLLjvr9j179qhPnz4+FbBjxw4tW7ZMpaWlOnz4sF566SXdeeedLdq4XC4tW7ZM1dXVyszM1G9/+1uNGDHCu93hcOgHP/iBIiIidP/99+vuu+/2qQbgDO5ADgBm8Wno6rbbbtPjjz+ub775ptW2hoYGzZ8/X+PGjfOpgPr6emVmZsrlcp11+8aNG5Wfn6/58+dr165dyszMVG5uro4cOeJt8+6776q0tFRbt27VkiVLtGfPnnN+3qlTp+R2u1s8AACAmXw6ozN37lxt3rxZl19+uWbNmqVBgwZJkvbt2yeXyyWPx6PHHnvMpwLGjh2rsWPHnnP78uXLNXPmTE2bNk2StGrVKr322mtas2aNHn30UUnfnmmSpD59+ui2227Trl27NHTo0LPub+nSpXriiSd8qhHhyd87kEvchRwAgoVPQSclJUXvvfee/uVf/kVz5syRZVmSvh06ys3NlcvlUkpKSsCKa2xsVGlpqebMmeN9LSIiQjk5OSoqKpL07Rmh5uZmde/eXXV1dXrzzTc1ceLEc+5zzpw5ys/P9z53u93eix4Cf49VcAAQ+nz+Ld6vXz+9/vrr+vrrr1VRUSHLsjRw4EBdeOGFAS/uyy+/lMfjaRWeUlJSvBcqrKmp0YQJEyRJHo9HM2fOVHZ29jn3GRsbq9jY2IDXGkoCcdsJbjkBAAgFfv+5euGFF543UHSVAQMGaPfu3XaXEVICMeGWW04AAEJBUB+pevXqpcjISNXU1LR4vaamRqmpqTZVFZx8WdrMMmgAQLgI6qATExOja665RoWFhd4l583NzSosLNSsWbM6tG+Xy+WdQG0Cf28F4cuEW245AQAINbYHnbq6OlVUVHifV1ZWqqysTD169FDfvn2Vn5+vvLw8ZWVlacSIEVqxYoXq6+u9q7D85XQ65XQ65Xa7lZSU1NGvEbKYcAsAMJntR7iSkhKNGjXK+/zMiqi8vDytW7dOkyZN0tGjRzVv3jxVV1dr2LBh2rZtW0BXd4WquOhIlS/M7fA+/OHv8BeTmAEAXcn2oHPzzTd7l6mfy6xZszo8VGUih8Nh29kYf4ewmMSMUEKgB0IfRxwAOAcCPRD6wvYn0bTJyF3B36EyOyYx+3utIFakAYBZwjboMBnZd3YMlXUksLBCDP4IpUAPoG1hG3QQGribOLqanXPfAAQeP80wnr835+TGnPAXk5iB4EHQQcjoSGDh4IGuxCRmIHjwE4WQwcUNAQC+CtujBquuAAQSk5iB4BS2QYdVVwACiUnMQHCKsLsAAACAzsKfH+hSvq5G4QJ+AICOIOigSzEXAQDQlRi6AgAAxuKMDjqdv6tRzrYfAAB8EbZBh+XlXYfVKAAAu4Tt0JXT6VR5ebmKi4vtLgUAAHSSsA06AADAfAQdAABgLIIOAAAwFkEHAAAYi6ADAACMRdABAADGCtug43K5lJGRoezsbLtLAQAAnSRsgw7X0QEAwHxhG3QAAID5CDoAAMBYBB0AAGAsgg4AADAWQQcAABiLoAMAAIxF0AEAAMaKsrsAAIB/LMtSQ5OnQ/uIi46Uw+EIUEVA8AnboONyueRyueTxdOyXBADYpaHJo4x52zu0j/KFuYqPCdtDAcJA2A5dcWVkAADMR4wHAAOUzM1RfExku9qebPQoa1FBJ1cEBAeCDgAYID4mkiEo4CzCdugKAACYj/gPAEHiZKNviyN8bQ+EI4IOAAQJ5s0AgcfQFQAAMBZndADARnHRkSpfmBuQ/QBojaADADZyOByslgI6EUNXAADAWAQdAABgrLANOi6XSxkZGcrOzra7FAAA0EnCNuhwrysAAMwXtkEHAACYj6ADAACMRdABAADGIugAAABjcZUqAAhj/t4YNC46Ug6HI8DVAIFH0AGAMObvjUTLF+ZyRWeEBIauAACAsYjjABBm/L2R6MlGj99ngAC7EHQAIMxwI1GEE4auAACAsYj0AACf+btaS2LFFroWQQcA4LOOzNVhxRa6EkNXAADAWERqAEC7+LtaS2LFFuxD0AEAtIsdq7Usy1JDk//zgSTmBIU7gg4AIGg1NHmUMW97h/bBnKDwFrZzdFwulzIyMpSdnW13KQAAoJOEbcR1Op1yOp1yu91KSkqyuxwAQBtK5uYoPiayXW2ZE4QzwjboAABCS3xMJENQ8FnYDl0BAADzEXQAAICxCDoAAMBYBB0AAGAsgg4AADAWQQcAABiLoAMAAIxF0AEAAMYi6AAAAGMRdAAAgLEIOgAAwFgEHQAAYCzujgYA6FInGz2d0tZulmWpoalj9cZFR8rhcASoIkgEHQBAF8taVGB3CZ2iocmjjHnbO7SP8oW53KE9wBi6AgAAxiI2AgA6XVx0pMoX5nZ4H6GiZG6O4mPaV+/JRo+xZ7mCAUEHANDpHA5HWA3JxMdEhtX3DWYMXQEAAGMRdAAAgLEIOgAAwFgEHQAAYCyCDgAAMJYxQefkyZPq16+ffv7zn9tdCgAACBLGBJ3Fixfr2muvtbsMAAAQRIwIOgcOHNC+ffs0duxYu0sBAABBxPags2PHDo0fP15paWlyOBzasmVLqzYul0v9+/dXt27dNHLkSH3wwQcttv/85z/X0qVLu6hiAAAQKmwPOvX19crMzJTL5Trr9o0bNyo/P1/z58/Xrl27lJmZqdzcXB05ckSS9PLLL+vyyy/X5Zdf3q7PO3XqlNxud4sHAAAwk+3Xpx47dux5h5yWL1+umTNnatq0aZKkVatW6bXXXtOaNWv06KOP6v3339eGDRu0adMm1dXVqampSYmJiZo3b95Z97d06VI98cQTnfJdAABAcLH9jM75NDY2qrS0VDk5Od7XIiIilJOTo6KiIknfBpeqqip9+umn+tWvfqWZM2eeM+RI0pw5c1RbW+t9VFVVdfr3AAAA9rD9jM75fPnll/J4PEpJSWnxekpKivbt2+fXPmNjYxUbGxuI8gAAQJAL6qDjq3vvvdfuEgAAQBAJ6qGrXr16KTIyUjU1NS1er6mpUWpqqk1VAQCAUBHUQScmJkbXXHONCgsLva81NzersLBQ1113XYf27XK5lJGRoezs7I6WCQAAgpTtQ1d1dXWqqKjwPq+srFRZWZl69Oihvn37Kj8/X3l5ecrKytKIESO0YsUK1dfXe1dh+cvpdMrpdMrtdispKamjXwMAEKRONnqM+hz4xvagU1JSolGjRnmf5+fnS5Ly8vK0bt06TZo0SUePHtW8efNUXV2tYcOGadu2ba0mKAMAcDZZiwrsLgE2sj3o3HzzzbIs67xtZs2apVmzZnVRRQAAwBS2Bx0AAAItLjpS5Qtzbf18BIewDToul0sul0seD2OqAGAah8Oh+JiwPcTh7wT1qqvO5HQ6VV5eruLiYrtLAQAAnSRsgw4AADAfQQcAABiLoAMAAIxF0AEAAMYK2ynprLoCAAQbf6+uHBcdKYfDEeBqzBC2QYdbQAAAgo2/V3EuX5jLcvpzYOgKAAAYi/gHAICN/L2K88lGD/fxageCDgAANuIqzp2LoSsAAGAsgg4AADBW2AYdl8uljIwMZWdn210KAADoJGEbdLipJwAA5gvboAMAAMxH0AEAAMYi6AAAAGOxcB8AgBDn7z2yJPPvk0XQAQAgxHXkCsmm3ycrbIeuWF4OAID5zI1wbeDu5QCAUObvPbKk8LpPVtgGHQAAQhn3yGqfsB26AgAA5iPoAAAAYxF0AACAsQg6AADAWAQdAABgLIIOAAAwFkEHAAAYK2yDDldGBgDAfGEbdJxOp8rLy1VcXGx3KQAAoJOEbdABAADmI+gAAABjEXQAAICxCDoAAMBY3PYUAIAwdrLR49f74qIj5XA4AlxN4BF0AAAIY1mLCvx6X/nCXMXHBH+MYOgKAAAYK/ijGAAACKi46EiVL8z1+X0nGz1+nwGyC0EHAIAw43A4QmLYKRAYugIAAMYi6NigqqpKK1eu1OnTpwPaFgAAOwXj8S1sg46dN/WcO3euFixYoJMnTwa0LQAAdgrG41t4DNCdhdPplNPplNvtVlJSUpd97q5du/SHP/xBK1euVGJiok42njvJfrctAADBypdjVlce3xyWZVmd+glB7kzQqa2t7fTOtixLo0ePVk1Njfbs2aOoqCidbDytjHnbJbW8JsHZ2gIAYKdAHLMCdXxr7/Gbo2cXeu211/T222/r1VdfbfN/rC9tAQCwUzAf3zij00VndE6fPq0hQ4YoLS1NBQUF3stmny0dn6stAAB26ugxK5DHN87oBJnVq1dr//79Wr9+fZv/Y31pCwCAnYL9+Ba2q666ktvt1vz58zVlyhQNHz48YG0BALBTKBzfCDpd4Mknn9SJEye0ePHigLYFAMBOoXB8I+h0sqqqKi1fvlwPPvigLrnkkvO2/cyHtgAA2MmXY5Yvx8JAI+h0srlz5yoxMVGPPPJIm22fWDCv3W0BALCTL8csX46FgcZk5E709xdE6t69+3nbnqqu0PrnnmtXWwAA7OTLMcuXY2FnCPvl5bW1tUpOTlZVVVVAl5dblqVx48bpyJEjKioqOue1Ak42nlb2ogIdefEJ9Ys7rfffP3dbAADs5Osxq73HQn+43W6lp6fr+PHj573DQdgHnc8++0zp6el2lwEAAPxQVVV13nk/YR90mpub9cUXX6h79+4BW9N/+vRpXXvtterTp4+2bt3a5sWT2ts2GJxJ0IE+Axbu6NfOQb8GHn3aOUKlX4Pp+GZZlk6cOKG0tDRFRJx7ynHYj5FEREQEfAb4ypUrVVFRoY0bN7Z5w1Bf2gaTxMTEoP5hDFX0a+egXwOPPu0cwd6vwXZ8a89+WXUVYGcuiDR16tR2XzypPW0BALBTqB7fCDoBduaCSIsWLQpoWwAA7BSqxzeCTgD5ckEkOy+e1BGxsbGaP3++YmNj7S7FKPRr56BfA48+7RzB3q+hfHwL+8nIgTR16lRt375dFRUVbV4rwJe2AADYKZSPb2E/GTlQdu3apT/+8Y/tvnhSe9sCAGCnUD++cUYnACzL0ujRo1VTU6M9e/a0efGk9rYFAMBOJhzfgqOKEPfaa6/p7bff1quvvtrm/1hf2gIAYCcTjm+c0emgpqYmDR06VGlpaSooKDjvBZF8aQsAgJ1MOb6x6qqDnn32We3fv1+/+tWv2vwf60tbu33++ee655571LNnT8XFxWnIkCEqKSnxbrcsS/PmzVOfPn0UFxennJwcHThwwMaKg5/H49Hjjz+uSy+9VHFxcfre976nX/ziF/r7vzXo17bt2LFD48ePV1pamhwOh7Zs2dJie3v68NixY7r77ruVmJio5ORkzZgxQ3V1dV34LYLP+fq1qalJjzzyiIYMGaILLrhAaWlpmjp1qr744osW+6BfW2rr3+rfu+++++RwOLRixYoWr9vZp8Yc3yz4rba21urdu7eVl5cX0LZ2O3bsmNWvXz/r3nvvtXbu3Gn97W9/s7Zv325VVFR42/zyl7+0kpKSrC1btli7d++27rjjDuvSSy+1GhoabKw8uC1evNjq2bOn9eqrr1qVlZXWpk2brISEBOvXv/61tw392rbXX3/deuyxx6zNmzdbkqyXXnqpxfb29OGtt95qZWZmWu+//771zjvvWJdddpk1efLkLv4mweV8/Xr8+HErJyfH2rhxo7Vv3z6rqKjIGjFihHXNNde02Af92lJb/1bP2Lx5s5WZmWmlpaVZTz31VIttdvWpScc3gk4HNDU1WStXrrSqqqoC2tZujzzyiHXDDTecc3tzc7OVmppqLVu2zPva8ePHrdjYWOtPf/pTV5QYkm6//XZr+vTpLV676667rLvvvtuyLPrVH989eLSnD8vLyy1JVnFxsbfNf//3f1sOh8P6/PPPu6z2YHa+g/IZH3zwgSXJOnjwoGVZ9GtbztWnn332mXXxxRdbe/futfr169ci6NjZpyYd3xi66oCoqCjdd9997bogki9t7bZ161ZlZWXpJz/5iS666CINHz5cq1ev9m6vrKxUdXW1cnJyvK8lJSVp5MiRKioqsqPkkHD99dersLBQn3zyiSRp9+7devfddzV27FhJ9GsgtKcPi4qKlJycrKysLG+bnJwcRUREaOfOnV1ec6iqra2Vw+FQcnKyJPrVH83NzZoyZYoeeughXXnlla2229mnJh3fgmdaNILG3/72N61cuVL5+fn693//dxUXF2v27NmKiYlRXl6eqqurJUkpKSkt3peSkuLdhtYeffRRud1uXXHFFYqMjJTH49HixYt19913SxL9GgDt6cPq6mpddNFFLbZHRUWpR48e9HM7ffPNN3rkkUc0efJk7w0o6VffPfnkk4qKitLs2bPPup0+DQyCDlppbm5WVlaWlixZIkkaPny49u7dq1WrVikvL8/m6kLXCy+8oOeff17r16/XlVdeqbKyMt1///1KS0ujXxEympqaNHHiRFmWpZUrV9pdTsgqLS3Vr3/9a+3atSv4Ju8ahqErtNKnTx9lZGS0eG3w4ME6dOiQJCk1NVWSVFNT06JNTU2Ndxtae+ihh/Too4/qpz/9qYYMGaIpU6bogQce0NKlSyXRr4HQnj5MTU3VkSNHWmw/ffq0jh07Rj+34UzIOXjwoN544w3v2RyJfvXVO++8oyNHjqhv376KiopSVFSUDh48qAcffFD9+/eXRJ8GCkEHrXz/+9/X/v37W7z2ySefqF+/fpKkSy+9VKmpqSosLPRud7vd2rlzp6677rourTWUnDx5UhERLX/kIiMj1dzcLIl+DYT29OF1112n48ePq7S01NvmzTffVHNzs0aOHNnlNYeKMyHnwIEDKigoUM+ePVtsp199M2XKFO3Zs0dlZWXeR1pamh566CFt375dEn0aMHbPhkbw+eCDD6yoqChr8eLF1oEDB6znn3/eio+Pt5577jlvm1/+8pdWcnKy9fLLL1t79uyx/uEf/oFl0G3Iy8uzLr74Yu/y8s2bN1u9evWyHn74YW8b+rVtJ06csD788EPrww8/tCRZy5cvtz788EPv6p/29OGtt95qDR8+3Nq5c6f17rvvWgMHDgzrZdCWdf5+bWxstO644w7rkksuscrKyqzDhw97H6dOnfLug35tqa1/q9/13VVXlkWfBgJBB2f1yiuvWFdddZUVGxtrXXHFFdYzzzzTYntzc7P1+OOPWykpKVZsbKx1yy23WPv377ep2tDgdrutn/3sZ1bfvn2tbt26WQMGDLAee+yxFgcK+rVtb731liWp1ePMNTza04dfffWVNXnyZCshIcFKTEy0pk2bZp04ccKGbxM8ztevlZWVZ90myXrrrbe8+6BfW2rr3+p3nS3o0Kcdxy0gAACAsZijAwAAjEXQAQAAxiLoAAAAYxF0AACAsQg6AADAWAQdAABgLIIOAAAwFkEHAAAYi6ADwGgLFiyQw+GQw+HQihUrbKvj008/9dYxbNgw2+oAwg1BB0CH3XvvvXI4HLrvvvtabXM6nXI4HLr33nu7vrD/78orr9Thw4f1T//0Ty1e//DDDzVp0iT16dNHsbGx6tevn8aNG6dXXnlF7b1o/Pjx43Xrrbeedds777wjh8OhPXv2KD09XYcPH9aDDz7Y4e8DoP0IOgACIj09XRs2bFBDQ4P3tW+++Ubr169X3759baxMioqKUmpqquLj472vvfzyy7r22mtVV1en3//+9/r444+1bds2TZgwQXPnzlVtbW279j1jxgy98cYb+uyzz1ptW7t2rbKysjR06FBFRkYqNTVVCQkJAfteANpG0AEQEFdffbXS09O1efNm72ubN29W3759NXz48BZtt23bphtuuEHJycnq2bOnxo0bp//93//1bm9sbNSsWbPUp08fdevWTf369dPSpUslSZZlacGCBerbt69iY2OVlpam2bNn+1RrfX29ZsyYodtvv12vvfaaxowZowEDBmjw4MGaMWOGdu/eraSkJG/7vXv3auzYsUpISFBKSoqmTJmiL7/8UpI0btw49e7dW+vWrWvxGXV1ddq0aZNmzJjhU20AAougAyBgpk+frrVr13qfr1mzRtOmTWvVrr6+Xvn5+SopKVFhYaEiIiI0YcIENTc3S5J+85vfaOvWrXrhhRe0f/9+Pf/88+rfv78k6cUXX9RTTz2l//qv/9KBAwe0ZcsWDRkyxKc6//znP+urr77Sww8/fM42DodDknT8+HGNHj1aw4cPV0lJibZt26aamhpNnDhR0rdni6ZOnap169a1GO7atGmTPB6PJk+e7FNtAAIryu4CAJjjnnvu0Zw5c3Tw4EFJ0l/+8hdt2LBBb7/9dot2P/rRj1o8X7NmjXr37q3y8nJdddVVOnTokAYOHKgbbrhBDodD/fr187Y9dOiQUlNTlZOTo+joaPXt21cjRozwqc5PPvlEkjRo0CDva8XFxRo1apT3+YYNGzRu3Dj97ne/0/Dhw7VkyZIW9aanp+uTTz7R5ZdfrunTp2vZsmX6n//5H918882Svh22+tGPftTizBCArscZHQAB07t3b91+++1at26d1q5dq9tvv129evVq1e7AgQOaPHmyBgwYoMTERO/ZmkOHDkn6dnJzWVmZBg0apNmzZ+vPf/6z970/+clP1NDQoAEDBmjmzJl66aWXdPr06Q7XPnToUJWVlamsrEz19fXefe7evVtvvfWWEhISvI8rrrhCkrzDbVdccYWuv/56rVmzRpJUUVGhd955h2ErIAgQdAAE1PTp07Vu3Tr9/ve/1/Tp08/aZvz48Tp27JhWr16tnTt3aufOnZK+nZsjfTvfp7KyUr/4xS/U0NCgiRMn6sc//rGkbyc979+/X08//bTi4uL0r//6r7rpppvU1NTU7hoHDhwoSdq/f7/3tdjYWF122WW67LLLWrStq6vT+PHjvSHozOPAgQO66aabvO1mzJihF198USdOnNDatWv1ve99Tz/4wQ/aXROAzkHQARBQt956qxobG9XU1KTc3NxW27/66ivt379fc+fO1S233KLBgwfr66+/btUuMTFRkyZN0urVq7Vx40a9+OKLOnbsmCQpLi5O48eP129+8xu9/fbbKioq0l//+td21zhmzBj16NFDTz75ZJttr776an300Ufq37+/NwideVxwwQXedhMnTlRERITWr1+vP/zhD5o+fbp3ng8A+zBHB0BARUZG6uOPP/b+93ddeOGF6tmzp5555hn16dNHhw4d0qOPPtqizfLly9WnTx8NHz5cERER2rRpk1JTU5WcnKx169bJ4/Fo5MiRio+P13PPPae4uLgW83jakpCQoGeffVaTJk3S7bffrtmzZ2vgwIGqq6vTtm3bWtTudDq1evVqTZ48WQ8//LB69OihiooKbdiwQc8++6y3XUJCgiZNmqQ5c+bI7Xbbet0gAP+HMzoAAi4xMVGJiYln3RYREaENGzaotLRUV111lR544AEtW7asRZvu3bvrP/7jP5SVlaXs7Gx9+umnev311xUREaHk5GStXr1a3//+9zV06FAVFBTolVdeUc+ePX2qccKECXrvvfcUHx+vqVOnatCgQRo9erTefPNN70RkSUpLS9Nf/vIXeTwejRkzRkOGDNH999+v5ORkRUS0/BU6Y8YMff3118rNzVVaWppP9QDoHA6rvZf/BIAQtGDBAm3ZskVlZWV2lyIp+OoBTMcZHQDG++tf/6qEhAQ9/fTTttVw6NAhJSQktFimDqDzcUYHgNGOHTvmncTcu3dv265rc/r0aX366aeSvl3hlZ6ebksdQLgh6AAAAGMxdAUAAIxF0AEAAMYi6AAAAGMRdAAAgLEIOgAAwFgEHQAAYCyCDgAAMBZBBwAAGOv/AYS7zCP5sZ+7AAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "artists = mass_hist.plot1d()\n", "ax = artists[0].stairs.axes\n", "ax.set_yscale(\"log\")\n", "ax.set_xlabel(\"Mass [GeV]\")\n", "ax.set_ylabel(\"Count\");" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "`pylhe` also has helpful graph representation of events so you can view what the LHE events you are studying look like." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [ { "data": { "image/svg+xml": [ "\n", "\n", "\n", "\n", "\n", "\n", "%3\n", "\n", "\n", "\n", "0\n", "\n", "u\n", "\n", "\n", "\n", "2\n", "\n", "e\n", "+\n", "\n", "\n", "\n", "0->2\n", "\n", "\n", "\n", "\n", "\n", "3\n", "\n", "e\n", "-\n", "\n", "\n", "\n", "0->3\n", "\n", "\n", "\n", "\n", "\n", "1\n", "\n", "\n", "\n", "\n", "\n", "1->2\n", "\n", "\n", "\n", "\n", "\n", "1->3\n", "\n", "\n", "\n", "\n", "\n" ], "text/plain": [ "" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "events = pylhe.read_lhe(lhe_file)\n", "this_event = next(events)\n", "this_event" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false } }, "outputs": [ { "data": { "image/svg+xml": [ "\n", "\n", "\n", "\n", "\n", "\n", "%3\n", "\n", "\n", "\n", "0\n", "\n", "u\n", "\n", "\n", "\n", "2\n", "\n", "Z\n", "0\n", "\n", "\n", "\n", "0->2\n", "\n", "\n", "\n", "\n", "\n", "1\n", "\n", "\n", "\n", "\n", "\n", "1->2\n", "\n", "\n", "\n", "\n", "\n", "3\n", "\n", "μ\n", "+\n", "\n", "\n", "\n", "2->3\n", "\n", "\n", "\n", "\n", "\n", "4\n", "\n", "μ\n", "-\n", "\n", "\n", "\n", "2->4\n", "\n", "\n", "\n", "\n", "\n" ], "text/plain": [ "" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "this_event = next(events)\n", "this_event" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "You can also render these graphs into a PDF file." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "tags": [] }, "outputs": [ { "data": { "text/plain": [ "'z0-event.pdf'" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "this_event.graph.render(filename=\"z0-event\", format=\"pdf\", cleanup=True)" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3.7.2 ('_env': venv)", "name": "pythonjvsc74a57bd073182083f4e2a7a153b1c22bd64054bec714030b17650bd3a885c40287b848ca" }, "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.2" }, "metadata": { "interpreter": { "hash": "73182083f4e2a7a153b1c22bd64054bec714030b17650bd3a885c40287b848ca" } } }, "nbformat": 4, "nbformat_minor": 4 }