{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Coffea-Casa Benchmark Example 8" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "%matplotlib inline\n", "from coffea import hist\n", "import coffea.processor as processor\n", "import awkward as ak" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "# This program plots the transverse mass of MET and a third lepton, where the third lepton is associated with a lepton pair\n", "# that has the same flavor, opposite charge, and closest mass to 91.2.\n", "\n", "import math\n", "\n", "class Processor(processor.ProcessorABC):\n", " def __init__(self):\n", " dataset_axis = hist.Cat(\"dataset\", \"MET and Third Lepton\")\n", " muon_axis = hist.Bin(\"massT\", \"Transverse Mass\", 50, 15, 250)\n", " \n", " self._accumulator = processor.dict_accumulator({\n", " 'massT': hist.Hist(\"Counts\", dataset_axis, muon_axis),\n", " 'cutflow': processor.defaultdict_accumulator(int)\n", " })\n", " \n", " @property\n", " def accumulator(self):\n", " return self._accumulator\n", " \n", " def process(self, events):\n", " output = self.accumulator.identity()\n", " \n", " dataset = events.metadata[\"dataset\"]\n", " \n", " # Keep track of muons and electrons by tagging them 0/1.\n", " muons = ak.with_field(events.Muon, 0, 'flavor')\n", " electrons = ak.with_field(events.Electron, 1, 'flavor')\n", " MET = events.MET\n", " \n", " output['cutflow']['all events'] += ak.size(events.MET, axis=0)\n", " \n", " # A few reasonable muon and electron selection cuts\n", " muons = muons[(muons.pt > 10) & (np.abs(muons.eta) < 2.4)]\n", " electrons = electrons[(electrons.pt > 10) & (np.abs(electrons.eta) < 2.5)]\n", " \n", " output['cutflow']['all muons'] += ak.sum(ak.count(muons, axis=1))\n", " output['cutflow']['all electrons'] += ak.sum(ak.count(electrons, axis=1))\n", "\n", " # Stack muons and electrons into a single array.\n", " leptons = ak.with_name(ak.concatenate([muons, electrons], axis=1), 'PtEtaPhiMCandidate')\n", " \n", " # Filter out events with less than 3 leptons.\n", " trileptons = leptons[ak.num(leptons, axis=1) >= 3]\n", " output['cutflow']['trileptons'] += ak.sum(ak.num(trileptons, axis=1))\n", " \n", " # Generate the indices of every pair; indices because we'll be removing these elements later.\n", " lepton_pairs = ak.argcombinations(trileptons, 2, fields=['i0', 'i1'])\n", " \n", " # Select pairs that are SFOS.\n", " SFOS_pairs = lepton_pairs[(trileptons[lepton_pairs['i0']].flavor == trileptons[lepton_pairs['i1']].flavor) & (trileptons[lepton_pairs['i0']].charge != trileptons[lepton_pairs['i1']].charge)]\n", " \n", " # Find the pair with mass closest to Z.\n", " closest_pairs = SFOS_pairs[ak.local_index(SFOS_pairs) == ak.argmin(np.abs((trileptons[SFOS_pairs['i0']] + trileptons[SFOS_pairs['i1']]).mass - 91.2), axis=1)]\n", " \n", " # Make trileptons and closest_pairs have same shape. First, fill nones with empty arrays. Then filter out events that don't meet the pair requirement.\n", " closest_pairs = ak.fill_none(closest_pairs, [])\n", " closest_pairs = closest_pairs[ak.num(closest_pairs) > 0]\n", " trileptons = trileptons[ak.num(closest_pairs) > 0]\n", " MET = MET[ak.num(closest_pairs) > 0]\n", " \n", " # Remove elements of the closest pairs from leptons, because we want the pt of the third lepton.\n", " trileptons_no_pair = trileptons[(ak.local_index(trileptons) != ak.flatten(closest_pairs.i0)) & (ak.local_index(trileptons) != ak.flatten(closest_pairs.i1))]\n", " \n", " # Find the highest-pt lepton out of the ones that remain.\n", " leading_lepton = trileptons_no_pair[ak.argmax(trileptons_no_pair.pt, axis=1)]\n", " output['cutflow']['number of final leading leptons'] += ak.sum(ak.num(trileptons_no_pair, axis=1))\n", " \n", " # Cross MET with the leading lepton.\n", " met_plus_lep = ak.cartesian({'i0': MET, 'i1': leading_lepton})\n", " \n", " # Do some math to get what we want.\n", " dphi_met_lep = (met_plus_lep.i0.phi - met_plus_lep.i1.phi + math.pi) % (2*math.pi) - math.pi\n", " mt_lep = np.sqrt(2.0*met_plus_lep.i0.pt*met_plus_lep.i1.pt*(1.0-np.cos(dphi_met_lep)))\n", " \n", " output['massT'].fill(dataset=dataset, massT=ak.flatten(mt_lep))\n", " \n", " return output\n", "\n", " def postprocess(self, accumulator):\n", " return accumulator" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "[########################################] | 100% Completed | 6min 59.8s\r" ] } ], "source": [ "from dask.distributed import Client\n", "import time\n", "\n", "client = Client(\"tls://localhost:8786\")\n", "\n", "fileset = {'SingleMu' : [\"root://eospublic.cern.ch//eos/root-eos/benchmark/Run2012B_SingleMu.root\"]}\n", "\n", "output = processor.run_uproot_job(fileset,\n", " treename = 'Events',\n", " processor_instance = Processor(),\n", " executor = processor.dask_executor,\n", " executor_args = {'schema': processor.NanoAODSchema, 'client': client}\n", " )" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAaIAAAEGCAYAAAAnhpGXAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjQuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/Z1A+gAAAACXBIWXMAAAsTAAALEwEAmpwYAAAiFklEQVR4nO3de3hV1bnv8e/PBAmK3KlF4lOgYkU0pXIxVms524q0csT2iOA+Kj2l0u2lm9bqUU8rtbbuI7Xeu/UULxt0Wy5aL2irgljUKsrFAoqogLBLKvWCStWKmvieP9ZMXAlrJSHJyiTJ7/M861lzvXOOMcecLPM6xxxrTEUEZmZmadkj7QaYmVnH5kRkZmapciIyM7NUORGZmVmqnIjMzCxVxWk3YHfRp0+fGDBgQNrNMDNrU1auXPlmRPRtTh1ORIkBAwawYsWKtJthZtamSPqv5tbhrjkzM0uVE5GZmaXKicjMzFLle0RmVjAff/wxFRUV7NixI+2mWDOVlJRQWlpKp06dWrzugiUiSbcC44DXI+KQJNYLmAcMADYDJ0fE28m6i4ApQBXwrxHxcBIfDswCugB/AKZFREjqDNwGDAe2ARMjYnNSZjLwk6Qpv4iI2YU6TjPLr6Kign322YcBAwYgKe3mWBNFBNu2baOiooKBAwe2eP2F7JqbBYytE7sQWBwRg4HFyWckHQxMAoYmZW6QVJSUuRGYCgxOXtV1TgHejogDgKuBGUldvYCfAocDo4CfSupZgOMzswbs2LGD3r17Owm1cZLo3bt3wa5sC5aIIuJx4K064fFA9dXJbODErPjciPgwIjYBG4BRkvoB3SJiaWSmCb+tTpnquu4CjlHm234csCgi3kquthaxc0I0s1biJNQ+FPLfsbUHK+wbEVsBkvfPJPH+wJas7SqSWP9kuW68VpmIqAS2A73rqWsnkqZKWiFpxRtvvNGMwzIzs6baXUbN5Uq1UU+8qWVqByNmRsSIiBjRt2+zfhhsZo0kidNOO63mc2VlJX379mXcuHEAzJo1i759+zJs2LCa1+rVq2uWe/XqxcCBAxk2bBhf+9rXWrx9Xbt2rfV527ZtNfv+7Gc/S//+/Ws+v/zyyxxyyCE565k+fTqPPPJIg/tbsmRJzbE3Jr6r7r33Xl544YVm11NIrT1q7jVJ/SJia9Lt9noSrwD2z9quFHg1iZfmiGeXqZBUDHQn0xVYAYyuU2ZJyx7Gp0YdeTRvvvVOznV9evVg2ZOPF2rXZm3S3nvvzfPPP88HH3xAly5dWLRoEf371+60mDhxIr/+9a9rxVatWgXAt7/9bcaNG8dJJ53UKu3t3bt3zb4vueQSunbtynnnnQfA5s2b85a79NJLc8arqqooKirKua4Q7r33XsaNG8fBBx/cavvcVa19RbQAmJwsTwbuy4pPktRZ0kAygxKWJd1370oqT+7/nF6nTHVdJwGPJveRHgbGSOqZDFIYk8QK4s233mHoWTfkfOVLUGYd3de//nV+//vfAzBnzhxOOeWUFqv7xBNPZPjw4QwdOpSZM2fWxLt27cqPf/xjvvjFL1JeXs5rr70GwKZNmzjiiCMYOXIkF1988S7vr6qqijPOOIOhQ4cyZswYPvjgAyCTMO+66y4gM4XYpZdeylFHHcWdd97JQw89xEEHHcRRRx3F3XffvUv7W7hwIUcccQSHHXYYEyZM4L333qvZxwUXXMCoUaMYNWoUGzZs4KmnnmLBggWcf/75DBs2jI0bN7Jq1SrKy8spKyvjm9/8Jm+//TYAo0ePril/4IEH8sQTT+zyuWiqgiUiSXOApcAXJFVImgJcDhwraT1wbPKZiFgLzAdeAB4Czo6IqqSqM4GbyQxg2Ag8mMRvAXpL2gCcSzICLyLeAn4OLE9elyaxJht15NEMGlKW81VVXNKcqs06pEmTJjF37lx27NjBmjVrOPzww2utnzdvXq2uueo/7o1x6623snLlSlasWMF1113Htm3bAHj//fcpLy9n9erVHH300dx0000ATJs2jTPPPJPly5fz2c9+dpePZf369Zx99tmsXbuWHj168Lvf/S7ndiUlJfzpT3/ixBNP5IwzzuD+++/niSee4G9/+1uj9/Xmm2/yi1/8gkceeYRnn32WESNGcNVVV9Ws79atG8uWLeOcc87hBz/4AV/+8pc54YQTuOKKK1i1ahWf//znOf3005kxYwZr1qzh0EMP5Wc/+1lN+crKSpYtW8Y111xTK15oBeuai4h8/4tzTJ7tLwMuyxFfAezUCRsRO4AJeeq6Fbi10Y1tQPVVj5m1jLKyMjZv3sycOXP4xje+sdP6XF1zjXXddddxzz33ALBlyxbWr19P79692XPPPWvuuQwfPpxFixYB8OSTT9Ykj9NOO40LLrhgl/ZXfb+qut583XUTJ04E4MUXX2TgwIEMHjwYgFNPPbXWlVt9nn76aV544QWOPPJIAD766COOOOKImvXVV5annHIKP/zhD3cqv337dt555x2++tWvAjB58mQmTPj0z+i3vvWtBo+jEDyzgpml4oQTTuC8885jyZIlNVctzbVkyRIeeeQRli5dyl577cXo0aNrfvvSqVOnmiHIRUVFVFZW1pRrztDkzp071ywXFRXlvXrbe++9m72/iODYY49lzpw5Oddn19uUfVQfS93zU2i7y6g5M+tgvvOd7zB9+nQOPfTQFqtz+/bt9OzZk7322osXX3yRp59+usEyRx55JHPnzgXgjjvuaLG25HPQQQexadMmNm7cCJA3qeRSXl7Ok08+yYYNGwD4xz/+wcsvv1yzft68eTXv1VdK++yzD++++y4A3bt3p2fPnjX3f26//faaq6M0ORGZWSpKS0uZNm1aznV17xE99dRTjapz7NixVFZWUlZWxsUXX0x5eXmDZa699lr+/d//nZEjR7J9+/ZdOoamKCkpYebMmRx//PEcddRRfO5zn8u77eLFiyktLa15bdiwgVmzZnHKKadQVlZGeXk5L774Ys32H374IYcffjjXXnstV199NZC5H3fFFVfwpS99iY0bNzJ79mzOP/98ysrKWLVqFdOnTy/4MTdEmYFmNmLEiMj3YLxBQ8qadI9o7Q1n8cq6Nc1tmlmbtW7dOoYMGZJ2MzqE6od79unTp2D7yPXvKWllRIxoTr2+R1RAVcUlDBpStlPcvy8yM/uUE1Fi3Ysv5UwaQJOHaJdNvSpnfO0NZzWpPjOzfFpzlFtLcyJKVFZVeYi2mVkKPFjBzMxS5URkZmapciIyM7NU+R6RmbUpXz12LFsqXm14w0bav3Q/Hlv0UIPbXXbZZfz2t7+lqKiIPfbYg9/85jfcdNNNnHvuuU2a2Xrz5s2MGzeO559/vt5tBg4cyE9+8hN+/vOfA5n55vr168f3vve9Jk+DtLtxIjKzNmVLxastOrCoMaNYly5dygMPPMCzzz5L586defPNN/noo4+4+eabW6wd+QwaNIgHHnigJhHdeeedDB06tOD7bU3umjMza8DWrVvp06dPzVxsffr0Yb/99mP06NFU/xA+32MmNm7cSHl5OSNHjmT69Ok7PXgPMo+SOP/88xk5ciRlZWX85je/qVnXpUsXhgwZUrOfefPmcfLJJ9esz37cRHU72honIjOzBowZM4YtW7Zw4IEHctZZZ/HYY4/ttE19j5mYNm0ay5cvZ7/99stZ/y233EL37t1Zvnw5y5cv56abbmLTpk0166sfm1FRUUFRUVHeetoqJyIzswZ07dqVlStXMnPmTPr27cvEiROZNWtWrW3qPmai+gemS5curXnUwj//8z/nrH/hwoXcdtttDBs2jMMPP5xt27axfv36mvVjx45l0aJFzJkzp+ZxEu2J7xGZmTVCUVERo0ePZvTo0Rx66KHMnj271vr6HjPRkIjg+uuv57jjjqsVr05me+65J8OHD+fKK69k7dq13H///TXbFBcX88knn9TU89FHHzXl8FLlKyIzswa89NJLta5QVq1aVe+s2dnKy8trHrxX/biJuo477jhuvPFGPv74YwBefvll3n///Vrb/OhHP2LGjBn07t27VnzAgAGsXLkSgPvuu6+mjrbEV0Rm1qbsX7pfi87XuH9pw/db3nvvPb7//e/zzjvvUFxczAEHHMDMmTM56aSTGix7zTXXcOqpp3LllVdy/PHH07179522+e53v8vmzZs57LDDiAj69u3LvffeW2uboUOH5hwtd8YZZzB+/HhGjRrFMcccU+sBfG2FHwOR6Nxlrxjzy4Wtsq81M8+lqHJHznWemdvaEz8GIvPwui5duiCJuXPnMmfOHO677760m9UkfgxEO5JvVm7wzNxm7c3KlSs555xziAh69OjBrbfemnaTdjtORGZmBfSVr3yF1atXp92M3ZoHK5hZQbn7v30o5L+jE5GZFUxJSQnbtm1zMmrjIoJt27ZRUtK0h4Q2xF1zZlYwpaWlVFRU8MYbb6TdFGumkpISSktLC1K3E5GZFUynTp0YOHBg2s2w3Zy75szMLFVORGZmlionIjMzS5UTkZmZpcqJyMzMUuVEZGZmqXIiMjOzVDkRmZlZqlJJRJJ+KGmtpOclzZFUIqmXpEWS1ifvPbO2v0jSBkkvSTouKz5c0nPJuuuUPB5RUmdJ85L4M5IGpHCYZmbWCK2eiCT1B/4VGBERhwBFwCTgQmBxRAwGFiefkXRwsn4oMBa4QVJRUt2NwFRgcPIam8SnAG9HxAHA1cCMVjg0MzNrgrS65oqBLpKKgb2AV4HxQPVD4GcDJybL44G5EfFhRGwCNgCjJPUDukXE0sjMqHhbnTLVdd0FHFN9tWRmZruXVk9EEfFX4FfAX4CtwPaIWAjsGxFbk222Ap9JivQHtmRVUZHE+ifLdeO1ykREJbAdqP2gd0DSVEkrJK2oqqxsmQM0M7NdkkbXXE8yVywDgf2AvSWdWl+RHLGoJ15fmdqBiJkRMSIiRhQVe/5XM7M0pPHX92vApoh4A0DS3cCXgdck9YuIrUm32+vJ9hXA/lnlS8l05VUky3Xj2WUqku6/7sBbBTqeFlVVXMKgIWU51/Xp1YNlTz7eyi0yMyusNBLRX4BySXsBHwDHACuA94HJwOXJ+33J9guA30q6iswV1GBgWURUSXpXUjnwDHA6cH1WmcnAUuAk4NFoI0/mKpt6Vd51a284qxVbYmbWOlo9EUXEM5LuAp4FKoE/AzOBrsB8SVPIJKsJyfZrJc0HXki2PzsiqpLqzgRmAV2AB5MXwC3A7ZI2kLkSmtQKh2ZmZk2gNnKhUHCdu+wVY365MO1m1GvtDWfxyro1aTfDzKyGpJURMaI5dXhmBTMzS5UTkZmZpcqJyMzMUuVEZGZmqXIiMjOzVDkRmZlZqpyIzMwsVU5EZmaWKiciMzNLlRORmZmlyonIzMxS5URkZmapciIyM7NUORGZmVmqnIjMzCxVTkRmZpYqJyIzM0uVE5GZmaXKicjMzFLlRGRmZqlyIjIzs1Q5EZmZWaqciMzMLFXFaTfAGq+quIRBQ8p2ivfp1YNlTz6eQovMzJrPiagNKZt6Vc742hvOauWWmJm1HHfNmZlZqpyIzMwsVU5EZmaWKiciMzNLlRORmZmlyonIzMxS5URkZmapSiURSeoh6S5JL0paJ+kISb0kLZK0PnnvmbX9RZI2SHpJ0nFZ8eGSnkvWXSdJSbyzpHlJ/BlJA1I4TDMza4S0roiuBR6KiIOALwLrgAuBxRExGFicfEbSwcAkYCgwFrhBUlFSz43AVGBw8hqbxKcAb0fEAcDVwIzWOCgzM9t1rZ6IJHUDjgZuAYiIjyLiHWA8MDvZbDZwYrI8HpgbER9GxCZgAzBKUj+gW0QsjYgAbqtTprquu4Bjqq+WzMxs95LGFdEg4A3gPyT9WdLNkvYG9o2IrQDJ+2eS7fsDW7LKVySx/sly3XitMhFRCWwHehfmcMzMrDnSSETFwGHAjRHxJeB9km64PHJdyUQ98frK1K5YmipphaQVVZWV9bfazMwKYpcTkaSeknaeArrxKoCKiHgm+XwXmcT0WtLdRvL+etb2+2eVLwVeTeKlOeK1ykgqBroDb9VtSETMjIgRETGiqNjzv5qZpaFRiUjSEkndJPUCVpPpVss9FXQDIuJvwBZJX0hCxwAvAAuAyUlsMnBfsrwAmJSMhBtIZlDCsqT77l1J5cn9n9PrlKmu6yTg0eQ+kpmZ7WYaexnQPSL+Lum7wH9ExE8lrWnGfr8P3CFpT+AV4H+RSYrzJU0B/gJMAIiItZLmk0lWlcDZEVGV1HMmMAvoAjyYvCAzEOJ2SRvIXAlNakZbzcysgBqbiIqT7rKTgR83d6cRsQoYkWPVMXm2vwy4LEd8BXBIjvgOkkRmZma7t8beI/oZ8DCwISKWSxoErC9cs8zMrKNo7BXR1oioGaAQEa809R6RmZlZtsZeEV3fyJiZmdkuqfeKSNIRwJeBvpLOzVrVDSjKXcrMzKzxGuqa2xPommy3T1b872SGRZuZmTVLvYkoIh4DHpM0KyL+q5XaZGZmHUhjByt0ljQTGJBdJiL+qRCNMjOzjqOxiehO4P8BNwNVDWxrZmbWaI1NRJURcWNBW2JNVlVcwqAhuaf/69OrB8uefLyVW2Rm1niNTUT3SzoLuAf4sDoYETtNJGqtr2xq/p90rb3hrFZsiZnZrmtsIqqeQPT8rFiQebaQmZlZkzUqEUXEwEI3xMzMOqZGJSJJp+eKR8RtLdscMzPraBrbNTcya7mEzCzZzwJORGZm1iyN7Zr7fvZnSd2B2wvSIjMz61B2+VHhiX+QeVKqmZlZszT2HtH9ZEbJQWay0yHA/EI1yszMOo7G3iP6VdZyJfBfEVFRgPaYmVkH06iuuWTy0xfJzMDdE/iokI0yM7OOo1GJSNLJwDJgAnAy8IwkPwbCzMyarbFdcz8GRkbE6wCS+gKPAHcVqmFmZtYxNHbU3B7VSSixbRfKmpmZ5dXYK6KHJD0MzEk+TwT+UJgmmZlZR1JvIpJ0ALBvRJwv6VvAUYCApcAdrdA+MzNr5xrqXrsGeBcgIu6OiHMj4odkroauKWzTzMysI2goEQ2IiDV1gxGxgsxjw83MzJqloURUUs+6Li3ZEDMz65gaSkTLJZ1RNyhpCrCyME0yM7OOpKFRcz8A7pH0P/k08YwA9gS+WcB2mZlZB1FvIoqI14AvS/pvwCFJ+PcR8WjBW2ZmZh1CY59H9EfgjwVui5mZdUCN/UGrtVFVxSUMGlKWc12fXj1Y9uTjrdwiM7PanIjaubKpV+Vdt/aGs1qxJWZmuaU2X5ykIkl/lvRA8rmXpEWS1ifvPbO2vUjSBkkvSTouKz5c0nPJuuskKYl3ljQviT8jaUCrH6CZmTVKmhOXTgPWZX2+EFgcEYOBxclnJB0MTAKGAmOBGyQVJWVuBKaSeWz54GQ9wBTg7Yg4ALgamFHYQzEzs6ZKJRFJKgWOB27OCo8HZifLs4ETs+JzI+LDiNgEbABGSeoHdIuIpRERwG11ylTXdRdwTPXVkpmZ7V7SuiK6BvjfwCdZsX0jYitA8v6ZJN4f2JK1XUUS658s143XKhMRlcB2oHfdRkiaKmmFpBVVlZXNPCQzM2uKVk9EksYBr0dEY2dmyHUlE/XE6ytTOxAxMyJGRMSIomKP2zAzS0Maf32PBE6Q9A0yc9l1k/SfwGuS+kXE1qTbrfpBfBXA/lnlS4FXk3hpjnh2mQpJxUB34K1CHZCZmTVdq18RRcRFEVEaEQPIDEJ4NCJOBRYAk5PNJgP3JcsLgEnJSLiBZAYlLEu6796VVJ7c/zm9Tpnquk5K9rHTFZGZmaVvd+qPuhyYn0yo+hdgAkBErJU0H3gBqATOjoiqpMyZwCwyM4E/mLwAbgFul7SBzJXQpNY6CDMz2zWpJqKIWAIsSZa3Acfk2e4y4LIc8RV8OgdednwHSSIzM7PdW5q/IzIzM3MiMjOzdDkRmZlZqpyIzMwsVU5EZmaWKiciMzNLlRORmZmlyonIzMxStTvNrGCtLN9jxP0IcTNrTU5EHVi+x4j7EeJm1prcNWdmZqlyIjIzs1Q5EZmZWaqciMzMLFVORGZmlionIjMzS5UTkZmZpcqJyMzMUuVEZGZmqXIiMjOzVDkRmZlZqpyIzMwsVU5EZmaWKs++bTvJ93gI8CMizKzlORHZTvI9HgL8iAgza3numjMzs1Q5EZmZWaqciMzMLFVORGZmlionIjMzS5UTkZmZpcqJyMzMUuVEZGZmqWr1RCRpf0l/lLRO0lpJ05J4L0mLJK1P3ntmlblI0gZJL0k6Lis+XNJzybrrJCmJd5Y0L4k/I2lAax+nmZk1ThozK1QCP4qIZyXtA6yUtAj4NrA4Ii6XdCFwIXCBpIOBScBQYD/gEUkHRkQVcCMwFXga+AMwFngQmAK8HREHSJoEzAAmtupRtlOe/sfMWlqrJ6KI2ApsTZbflbQO6A+MB0Ynm80GlgAXJPG5EfEhsEnSBmCUpM1At4hYCiDpNuBEMoloPHBJUtddwK8lKSKiwIfX7nn6HzNraaneI0q6zL4EPAPsmySp6mT1mWSz/sCWrGIVSax/slw3XqtMRFQC24HeOfY/VdIKSSuqKitb6KjMzGxXpJaIJHUFfgf8ICL+Xt+mOWJRT7y+MrUDETMjYkREjCgq9vyvZmZpSCURSepEJgndERF3J+HXJPVL1vcDXk/iFcD+WcVLgVeTeGmOeK0ykoqB7sBbLX8kZmbWXGmMmhNwC7AuIrJvOCwAJifLk4H7suKTkpFwA4HBwLKk++5dSeVJnafXKVNd10nAo74/ZGa2e0qjP+pI4DTgOUmrktj/AS4H5kuaAvwFmAAQEWslzQdeIDPi7uxkxBzAmcAsoAuZQQoPJvFbgNuTgQ1vkRl1Z2Zmu6E0Rs39idz3cACOyVPmMuCyHPEVwCE54jtIEpmZme3efIfeWky+3xj590VmVh8nImsx+X5j5N8XmVl9PNecmZmlyonIzMxS5URkZmapciIyM7NUORGZmVmqnIjMzCxVHr5tBednGJlZfZyIrOD8DCMzq4+75szMLFVORGZmlionIjMzS5XvEVmqPJDBzJyILFUeyGBm7pozM7NUORGZmVmqnIjMzCxVvkdkuy0PZDDrGJyIbLflgQxmHYO75szMLFVORGZmlip3zVmblO/+ke8dmbU9TkTWJuW7f+R7R2Ztj7vmzMwsVb4isnbFQ77N2h4nImtXPOTbrO1xIrIOw1dLZrsnJyLrMHy1ZLZ7ciIyw8PBzdLkRGRG/qulNTPPdXeeWYE5EZnVo77uPCcps5bhRGTWRE1JUk5QZjtr14lI0ljgWqAIuDkiLk+5SdZBNKWrrz5OYNaeKSLSbkNBSCoCXgaOBSqA5cApEfFCru07d9krxvxyYSu20Kzx1sw8l6LKHS1WnxObtRRJKyNiRHPqaM9XRKOADRHxCoCkucB4IGciMtud1dcN2BRNvTLb3TnBtk3tORH1B7Zkfa4ADs/eQNJUYGry8ZMHpn21spXaVjgReyB9knYzdls+P/m1g3OzKeITFeYYioG2//ehML7U3AracyJSjlitfsiImAnMbJ3mtA5JK+KTT5p1mdye+fzk53OTn6QVze1+aq8krWhuHe159u0KYP+sz6XAqym1xczM8mjPiWg5MFjSQEl7ApOABSm3yczM6mi3XXMRUSnpHOBhMsO3b42ItSk3qzW0q67GAvD5yc/nJj+fm/yafW7a7fBtMzNrG9pz15yZmbUBTkRmZpYqJ6I2TtJmSc9JWlU9jFJSL0mLJK1P3num3c7WIOlWSa9Lej4rlvdcSLpI0gZJL0k6Lp1Wt4485+YSSX9NvjurJH0ja11HOjf7S/qjpHWS1kqalsQ7/HennnPTst+diPCrDb+AzUCfOrFfAhcmyxcCM9JuZyudi6OBw4DnGzoXwMHAaqAzMBDYCBSlfQytfG4uAc7LsW1HOzf9gMOS5X3ITA12sL879Z6bFv3u+IqofRoPzE6WZwMnpteU1hMRjwNv1QnnOxfjgbkR8WFEbAI2kJkWql3Kc27y6WjnZmtEPJssvwusIzMzS4f/7tRzbvJp0rlxImr7AlgoaWUyZRHAvhGxFTJfJOAzqbUuffnORa4poOr7D6y9OkfSmqTrrrrrqcOeG0kDyExZ8wz+7tRS59xAC353nIjaviMj4jDg68DZko5Ou0FtRINTQHUANwKfB4YBW4Erk3iHPDeSugK/A34QEX+vb9McsXZ9fnKcmxb97jgRtXER8Wry/jpwD5nL4Nck9QNI3l9Pr4Wpy3cuOvwUUBHxWkRURcQnwE182oXS4c6NpE5k/tDeERF3J2F/d8h9blr6u+NE1IZJ2lvSPtXLwBjgeTJTGU1ONpsM3JdOC3cL+c7FAmCSpM6SBgKDgWUptC811X9kE98k892BDnZuJAm4BVgXEdnP2+jw351856alvzvtdoqfDmJf4J7Md4Vi4LcR8ZCk5cB8SVOAvwATUmxjq5E0BxgN9JFUAfwUuJwc5yIi1kqaT+b5VJXA2RFRlUrDW0GeczNa0jAyXSebge9Bxzs3wJHAacBzklYlsf+DvzuQ/9yc0pLfHU/xY2ZmqXLXnJmZpcqJyMzMUuVEZGZmqXIiMjOzVDkRmZlZqpyIrMOT1DtrFuG/1ZlVeM+029dSJIWk27M+F0t6Q9IDabbLzL8jsg4vIraRmaoESZcA70XEr6rXSyqOiMp0WvepFmjH+8AhkrpExAfAscBfW6Z1Zk3nKyKzHCTNknSVpD8CMySNkvSUpD8n719Itvu2pLslPZQ8t+aXSbwoqeN5ZZ4X9UNJQyQty9rHAElrkuXhkh5LJq99OGtqmSWS/k3SY8A0SROSOldLejxrX1dIWp5MQvm9eg7tQeD4ZPkUYE5We/Id41BJy5IrxDWSBiezevw+acfzkia21Lm3jsdXRGb5HQh8LSKqJHUDjo6ISklfA/4N+B/JdsPIzEr8IfCSpOvJzNTcPyIOAZDUIyLekbSnpEER8Qowkcwv9zsB1wPjI+KN5I/6ZcB3kvp7RMRXk3qeA46LiL9K6pGsnwJsj4iRkjoDT0pamEzDX9dcYHrSHVcG3Ap8JVn3Yp5j/Bfg2oi4I+mqLAK+AbwaEccn7erexHNs5kRkVo87s6Yn6Q7MljSYzLQmnbK2WxwR2wEkvQB8DlgLDEqS0u+Bhcm284GTyUwfMzF5fQE4BFiUTNdURGZG42rzspafBGYl06hUT845BiiTdFJWWwcDOyWiiFijzHT+pwB/qLM63zEuBX4sqRS4OyLWJwnxV5JmAA9ExBN192XWWO6aM8vv/azlnwN/TK5w/jtQkrXuw6zlKqA4It4GvggsAc4Gbk7WzwNOlnQgEBGxnszU+WsjYljyOjQixuRqR0T8C/ATMjMcr5LUOyn//azyAyNiIfktAH5FVrdcfccYEb8FTgA+AB6W9E8R8TIwHHgO+L+SptezP7N6ORGZNU53Pr2x/+2GNpbUB9gjIn4HXEzmMd1ExEYyyepiPr3SeQnoK+mIpGwnSUPz1Pv5iHgmIqYDb5JJSA8DZyZdfEg6UJnZ2PO5Fbg0Ip5rzDFKGgS8EhHXkUliZZL2A/4REf9JJqkd1tA5McvHXXNmjfNLMt1W5wKPNmL7/sB/SKr+n72LstbNA64ABgJExEdJt9p1yb2WYuAaMt17dV2RdJ0JWAysBtYAA4Bnlenbe4N6Hg8fERXAtbtwjBOBUyV9DPwNuBQYmbTlE+Bj4Mx8+zNriGffNjOzVLlrzszMUuVEZGZmqXIiMjOzVDkRmZlZqpyIzMwsVU5EZmaWKiciMzNL1f8Hh9+YaP1u5IgAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "hist.plot1d(output['massT'], overlay='dataset', fill_opts={'edgecolor': (0,0,0,0.3), 'alpha': 0.8})" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "all events 53446198\n", "all muons 991421725\n", "all electrons 87391040\n", "trileptons 5043516\n", "number of final leading leptons 1227754\n" ] } ], "source": [ "for key, value in output['cutflow'].items():\n", " print(key, value)" ] } ], "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.8.8" } }, "nbformat": 4, "nbformat_minor": 4 }