{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "d2b93495",
   "metadata": {},
   "source": [
    "# df014_CSVDataSource\n",
    "Process a CSV file with RDataFrame and the CSV data source.\n",
    "\n",
    "This tutorial illustrates how use the RDataFrame in combination with a\n",
    "RDataSource. In this case we use a TCsvDS. This data source allows to read\n",
    "a CSV file from a RDataFrame.\n",
    "As a result of running this tutorial, we will produce plots of the dimuon\n",
    "spectrum starting from a subset of the CMS collision events of Run2010B.\n",
    "Dataset Reference:\n",
    "McCauley, T. (2014). Dimuon event information derived from the Run2010B\n",
    "public Mu dataset. CERN Open Data Portal.\n",
    "DOI: [10.7483/OPENDATA.CMS.CB8H.MFFA](http://opendata.cern.ch/record/700).\n",
    "\n",
    "\n",
    "\n",
    "\n",
    "**Author:** Enric Tejedor (CERN)  \n",
    "<i><small>This notebook tutorial was automatically generated with <a href= \"https://github.com/root-project/root/blob/master/documentation/doxygen/converttonotebook.py\">ROOTBOOK-izer</a> from the macro found in the ROOT repository  on Tuesday, May 19, 2026 at 08:09 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "8d5e21e9",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:41.514369Z",
     "iopub.status.busy": "2026-05-19T20:09:41.514227Z",
     "iopub.status.idle": "2026-05-19T20:09:42.462108Z",
     "shell.execute_reply": "2026-05-19T20:09:42.461099Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2e54bc9e",
   "metadata": {},
   "source": [
    "Let's first create a RDF that will read from the remote CSV file.\n",
    "The types of the columns will be automatically inferred."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "585fcd01",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:42.463844Z",
     "iopub.status.busy": "2026-05-19T20:09:42.463696Z",
     "iopub.status.idle": "2026-05-19T20:09:43.075757Z",
     "shell.execute_reply": "2026-05-19T20:09:43.075237Z"
    }
   },
   "outputs": [],
   "source": [
    "fileUrl = \"http://root.cern/files/tutorials/df014_CsvDataSource_MuRun2010B.csv\"\n",
    "df = ROOT.RDF.FromCSV(fileUrl)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "021f90e0",
   "metadata": {},
   "source": [
    "Now we will apply a first filter based on two columns of the CSV,\n",
    "and we will define a new column that will contain the invariant mass.\n",
    "Note how the new invariant mass column is defined from several other\n",
    "columns that already existed in the CSV file."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "9ebe405d",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:43.080310Z",
     "iopub.status.busy": "2026-05-19T20:09:43.080188Z",
     "iopub.status.idle": "2026-05-19T20:09:43.295830Z",
     "shell.execute_reply": "2026-05-19T20:09:43.295285Z"
    }
   },
   "outputs": [],
   "source": [
    "filteredEvents = df.Filter(\"Q1 * Q2 == -1\") \\\n",
    "                   .Define(\"m\", \"sqrt(pow(E1 + E2, 2) - (pow(px1 + px2, 2) + pow(py1 + py2, 2) + pow(pz1 + pz2, 2)))\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "586b0649",
   "metadata": {},
   "source": [
    "Next we create a histogram to hold the invariant mass values and we draw it."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "71cd4cc0",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:43.314401Z",
     "iopub.status.busy": "2026-05-19T20:09:43.314258Z",
     "iopub.status.idle": "2026-05-19T20:09:46.182904Z",
     "shell.execute_reply": "2026-05-19T20:09:46.181941Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file df014_invMass.png has been created\n"
     ]
    }
   ],
   "source": [
    "invMass = filteredEvents.Histo1D((\"invMass\", \"CMS Opendata: #mu#mu mass;#mu#mu mass [GeV];Events\", 512, 2, 110), \"m\")\n",
    "\n",
    "c = ROOT.TCanvas()\n",
    "c.SetLogx()\n",
    "c.SetLogy()\n",
    "invMass.Draw()\n",
    "c.SaveAs(\"df014_invMass.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "4348b38c",
   "metadata": {},
   "source": [
    "We will now produce a plot also for the J/Psi particle. We will plot\n",
    "on the same canvas the full spectrum and the zoom in on the J/psi particle.\n",
    "First we will create the full spectrum histogram from the invariant mass\n",
    "column, using a different histogram model than before."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "f3073135",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:46.184682Z",
     "iopub.status.busy": "2026-05-19T20:09:46.184523Z",
     "iopub.status.idle": "2026-05-19T20:09:46.288487Z",
     "shell.execute_reply": "2026-05-19T20:09:46.288016Z"
    }
   },
   "outputs": [],
   "source": [
    "fullSpectrum = filteredEvents.Histo1D((\"Spectrum\", \"Subset of CMS Run 2010B;#mu#mu mass [GeV];Events\", 1024, 2, 110), \"m\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a68ed3a9",
   "metadata": {},
   "source": [
    "Next we will create the histogram for the J/psi particle, applying first\n",
    "the corresponding cut."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "2d10b629",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:46.291267Z",
     "iopub.status.busy": "2026-05-19T20:09:46.291144Z",
     "iopub.status.idle": "2026-05-19T20:09:46.410938Z",
     "shell.execute_reply": "2026-05-19T20:09:46.410545Z"
    }
   },
   "outputs": [],
   "source": [
    "jpsiLow = 2.95\n",
    "jpsiHigh = 3.25\n",
    "jpsiCut = 'm < %s && m > %s' % (jpsiHigh, jpsiLow)\n",
    "jpsi = filteredEvents.Filter(jpsiCut) \\\n",
    "                     .Histo1D((\"jpsi\", \"Subset of CMS Run 2010B: J/#psi window;#mu#mu mass [GeV];Events\", 128, jpsiLow, jpsiHigh), \"m\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "b98fa322",
   "metadata": {},
   "source": [
    "Finally we draw the two histograms side by side."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "1a94ef22",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:46.414884Z",
     "iopub.status.busy": "2026-05-19T20:09:46.414753Z",
     "iopub.status.idle": "2026-05-19T20:09:47.843173Z",
     "shell.execute_reply": "2026-05-19T20:09:47.842337Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file df014_jpsi.png has been created\n"
     ]
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "Saved figures to df014_*.png\n"
     ]
    }
   ],
   "source": [
    "dualCanvas = ROOT.TCanvas(\"DualCanvas\", \"DualCanvas\", 800, 512)\n",
    "dualCanvas.Divide(2, 1)\n",
    "leftPad = dualCanvas.cd(1)\n",
    "leftPad.SetLogx()\n",
    "leftPad.SetLogy()\n",
    "fullSpectrum.Draw(\"Hist\")\n",
    "dualCanvas.cd(2)\n",
    "jpsi.SetMarkerStyle(20)\n",
    "jpsi.Draw(\"HistP\")\n",
    "dualCanvas.SaveAs(\"df014_jpsi.png\")\n",
    "\n",
    "print(\"Saved figures to df014_*.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "847a7cd0",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "e14f4334",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:09:47.844561Z",
     "iopub.status.busy": "2026-05-19T20:09:47.844435Z",
     "iopub.status.idle": "2026-05-19T20:09:48.036904Z",
     "shell.execute_reply": "2026-05-19T20:09:48.036405Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779221388024\" style=\"width: 700px; height: 500px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779221388024() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(23124,'WkwIBCEAVFoAeAHlnGuTXLeRpv9KR40/Qh1I3HEQ84EX0fIuKTJEyiLH4dgosqvJWnZ39VQXJcoO/feJJ4FTVd0kFbbHjrF3ZVfzJC4HCSDx5gXA+fPi/+x+vl5dLS9Xi2nx4sHy6sflzQ+r18+vltc37za7hVmcf3+1/s8Pq989XEzWLM7vr3c3/enp6/+7erMjfUGxp9e79eZqEP97fXW2mLxZnO/fNP35c219qQHvg/MhmMX54/XV6sHmYrNdTDLI57ufL1YH8of12e5dJx+tLy5GYZiFnAtbq9VX57sny+3b9dVisqekfLd+++5O0v3Nbre5vF3sxeb6dsLL8zVMOLM4f3V4vNcfefHLm91yRyu1UuYWda9TVH60XV6u7vJN2p2O78vd7tC+6Jw8v5LqxyPDK+9vtmer7fP1n8boHSU+2Zyt+ry+FJi2xYcSUvAulWKrdlMW01f2tCSXXYjFR+uLj/TULSZ36mKpvlifY/RM3Cu3mMJpjsE5W5O1tSSSX+42917fPFt/XF28X0xfuXSabY1SbExOJPcSf0n2YvJST6v33uYkqUZYeXXr7b7W02iDleBTrSpOr3b7l/9K7mL6qoRT66lobUo+J1bCrZcjYN8f3nZMLqZUqfD7WxVCZhR/f6hyO2ExfdUT5uHZbV6+/8JU6Pj8Bfla3XpxMUUJrsRUCxJyaOLV+8/O0miA3M9O+Zyv2VakWutt8TFahP3lxeanbx8+6PL06ph4+f21ZsDEq6PnH/ap3+yf7r2+ufWie69vbr3r3uubQ7V7r28ONb//eNkXuLdiXa01puR7z7//WbO+sp/L+3i5/Igo2yC+ulSixOKYtO9/1pxwKtVFl3x2LpZMV1+8W+2Wi8kz/c/ercfTvZvr1Zvdd8vdetMH4dsPl69X2/78Yv3m/cfD48/98fHm7ceBaZu3P++f/tRzny3Pni3XV+CJWZw/2G5ubt4t1+OFe/LZZmDz8UKH/04f1viTzdn6fL06W0zny4ublVmc/3a7Pvt4m/z5QN57ffNgs9kelf/6bL1bvgYBd9sPvODR+uPq7Fa/51c/264v17v1j6ubT5TA4/UNOmbWP4NcbreL6Q9/NIvN9Y6HX8zi/OuPqzc3i+nqw8WFWZx/2xXWG0F/vFjv4GNQ3364fLa8WO12s9ZgwL5dfdx9mvrwd8+fPb73ajEtfjM/msX5w82H1xer+x/Oz+cJ+261W66vGK3R15c36z+tvr+Z81/dJjX3u9XyYjE5GtfsTks4dSVmpP+H9dXZ5qcXm2sW6TH96pgeCH4o8M0KfTXE4KcZaB68Www0ebDc7T4Z5nu7XVft9O/l/dXup9XqauiuW5QO5aPt5vLF5noxySmy8/JsuUM1KPFqJtDt9zohv5jF+yebH1dPr5f/+WEvEe+/WzEutxPPv1m/ffeYLgw1rdK53L15Nw/u++fvNj99/ePqavd8t9x9uNkL4ft7H3YbxGBf8snq6sP95bbTCMm9N4jZvsb5d6vl2dOri5/nGuc/rHfvNh92xxI5S+k3y5shY3PKcak/3LFe/m42EirxizbSD6vXuurXV2+/ZCghGQ8uljc3Y0lQrltmxwnXiOrCTi5GM35NJmussc1pKk9+spprW9A8F2OL+zq2pVGS9DzK8lymUo0TMSW0ynMJxklqQoPBjF8TmaQ6I8UZya6Jm3Iy+v8mfhLrzPg1CZO4bMavSZwkejN+TdIkuZjxa5InZ4MZvyZlcq6a8WtSJxecGb/mrBaWSn3bnNwm3STVG6liJJXm/CQlGSmQvrnQ203RSLTNxUmKJ8dIjM3B1WAy2ObyJDEaCc6It82VSXzUDpbUXJ0kR+NsgpPm7STeGW06+eZlEp/hwEiyzbtJvDUSi+HF3vd2Sul1GSrX24m++ThhRGn9EJtPk9gAB8qyz5PYarSCd82XSVzQcXW2Nl8Pox5jC3aSbI12t0oLMjHFyrOPLbjJOTE6PpAMVelMettCmHSImTUXWoiTjilsSG0hTbymBlNTC7k3w0CF3EKZtEVRNluok4RstBPiWrTKhIpYizJJSEZltkWn4zAIP9UyP4cph/k5TlHm5zS5PD/nyc+PBcEer6mTNQEBji3ZyZpU+jNLp0p/dhOd0mXVEitIfO05LCKkQ6tHiJI7kWhDeiNJV5IfrWjrukilJZpXwuWWaV8JG1qGAQgptmU4UCKmluc1zPxmONAcW1qGA4gSW1YGYjRJWp5Xss8tz82Lb7lOnUNjW7FTZ5BnmaJ2g2c3FR0GnlnCHV1sK0ilDgNEX760Z1thmaSZYOl2rmwrLFwdBgiWrY62sa0O4PIFQjoiSYDo4CUV+KpecyQliNCJ4CA6iiFUttUOY1XfnLVUFtKLPofMc9Vn541tYu0kIiZEBQWxohIanAmuiXUdFCrTK02sn0o2olJDdtB+SEzGQ8bJuWJ8NjE0sWmSaE1mgdYmFshIWlYSr2J5OuO9cQ6SITGJUZUmiq2RNeSaCDAmxotxUE5HFTBKtQnYyjM/F5sArswmA897JXaOuzwL6GrF1GrE2iaSJwCO/khpAroid6kaAdHrJKXSsDRxdpIAZgUjCXwHx/q4FUg3aTWbTalNnJ+c13ZsE6C1RNSCYZkL2CreKJQ5Cidd6k5K5xl0PVIBzspUhSVZpTmmBEjsv+Ys2OTN+DXHnAj90V9zNna8BdabY1KKVwiwzdk8sTazmMyby+R87WxJbs7WiVqAdWyOKUEKqNycyJTF9MrNiVNcHrWbY1Joj8rNSaBmr9ycRGr2ys0xHyC6ZGo3J7lzOzda9p3TVqvyOzeLugOOOtPN6ZTAbW+XKel9leKb0ynRzjofmnNh31ttmUkZvYVpB4ofRtXlub+9ctl3mMp+XsNgmw3Nz8sYYUJvsZKHFu90X8yab13zrOc5X+kjm0DpvqopH6V5FvYormTZWxBKgm3d3IjSgj0YRrYF220CqkO5fU2og50CFVAnaJQWLKPR0Yv+BV3Miv3an8BihkYbR98ComPdvn9BbKfJt64FjCVYRIMq7WhrdC+IV8oG5V+UD4rSG1ElR2WlZqi3LYhCfe+XAPVUIR1FQwXbgkPP9C65MSXdZmwBc2h0CZUWWLSjS9pl5GN0qeeDJE673POxPkaXGAJX5i714nXuk5b2FpIGtLAXJfusteBnpadVvWq9Pt8teFW8+4qqeff1GA6dtBYRySEjvCVW2U91pw9TD0exHia/54fen97HFhHJ3o6yHGuC5T5rLdY8U/ou7XvvYItVuw5lXUuWKVCt2JLtVjp1kqXPOmktWXqsk9aS7dY6FgNS2FeYspiAr77CBl20y2PaWkKnkN+nrSWFMJXA3iRS2CVQu5RQLbrC+uvE6+voIgyiWVQIR+E+HnO35Nh9SMNqB/Wo2VfH6BxaBWnsw+FUEJi1lhDJuX/UA7FG/5CSBGT1FTboofjH9CfFrN4/HWq10nWFDbqrty4gydW5e7SNld57p1Rfn6OkPzhRyfdBwSJJWOesMX3uw9GfVTR4zDMQdrZbHkA4hr1lpG6wGKVlhG6QVI77/mrmYfrpXh7mjbLQ8jBwaEkL9971zDJjoKJDAQPHuEZpxXYk0ElpZYCggkSx+zknZ7DTsacMaRz8tjJLoyJNQRaH7Gkb+8HGxuyI3JsYQkhhCqJKB6xQUpFwvPIAhOQMOOrdLejRGZCtawU92sVeC5c9pW3ssQgjVwWQ5iGGduij0woqdMZ1XgtIIm0djAsYOUgqD0t4sDSkcYxCK0Ma55b6WoWCJXzG0RDG9sCuzhTCOPrGsii4jGOIWBbF90HqctPKASDh2Ku+gKdedqD2UNTFz/5BL3zAbG1IXZS5ajj2UUpQxOY9+t6giN0XXytBPaUhJ2EGbPoZDuzAXTiw01+zZ4eOhT07vfDMDi+KMzdKjBHqKroV3MUxL9qRuBcjfW+Eo4HOpXuNY0BbiZ2jMeCtRDjSXEIhUTk6kHA08LGoDzn0QUldofVpa6X7kX3SWlFHcgBn6Y5k57aV1NV7F4VW0pGSpX11Jw9v7eMz91o9yn3V7lLOeXmEQmaG8tB+M0u5C/Y8KLkLdkfskrtgzx3N3Rjaj0LGs5k97JK7YM8DgaPZM5Go0gW7q79S+rTNRUufta40StE561qiFLWBBjt4nCpbfXLxOQdm6hiVvvTnNvrS37ehS380UYfu0WoDqecRGEA92h/GgaqpMkB6EB0U59kcGL2fz4HR82v2CkgFemC0vqgOiB7EAaEZt3oE0UhwBaSHLrau1QHTvZU6cHr0pM44Paa6Dpzug1mHzTBT3WIYM10HUI/MEV1TmYIp3A0Gvv9aBSyR4/5r1ekMdh1aXZiix60WT1UV84BrnFt1aYpWPfCA05+nUHHUJfDSMoVicDxDatXVKSQTrZFQWvV2CsCnkSiteplCMMpMaNW7KXij8JRb9X4KzuB+JyIJYQrWEHFLrlUfJ19NrMTCWvVp8sWkHiKsPk8+G1zkTAtl8tHgAGdaqIRNUjKSU6vBTt6blI3k2mqQyTu8eXzjGtzk1eWTElsNfnLV4AuX3GoIkysaJKiENSJhqxyJFrZKJC2ZnI1UWshECHI1UmkBZ9UUJFlaDXVyzpRu11WNoxl1/XKrgGM1JRKTahVoLKYk48S3GokemFKMk9hqDJMkU4lelVajOoQEwpy0Ci5GU92gWPIannC8E/uOkEun1G2uqVMJu8rUPChExlT8fEI1xAyIRTiprSYPJQQulFRmVGqVjD3Xj9zUyR64qCl3kk6WVlOhW2LpJSQBWUNoQslM6NMQKnWSW82iJOEwJQn5GEJCndQREnGDJEpsiF4QAq8gYzGisQZIogs9Xqhk7iTtkoszb4QoB6MNMkLCRmy1WCWJFCgpnYQNcolda9Cjk76TsBFaJSBHLt2HJCBrxNF9yNRJ2IDMnYQN3yoASWHYgKxKetjwrdbOFbEnca3iRBWij4PsY0W8SaRVUJLoNGxA9hn0sAHZZ9DDBjE4FSjxsAGpEqXRVSVVpCQg37VVXPpghCiWJZZG0A4agS/QBKMAjJkmqm4kEI7QfD8JEehknM2UJ4JJRHqmY6fLTBOWMRIqEXHKZ6Uj7CiNVgF2jLORfA2nSoQfaLFTHapYSVFStzHIdRPB9B4rESteSXjR3KAkrHhy40RIB06UTJAJRog4EsAjnjiTRUnY0Nw6Ec1jVAgxOqt1NToIKTQEnFXG1DmYZPujMmRq6mq4VSpDhq3rBk2fsHaJQYJOsK1BvEHDqJq7fZujwgshIzgHv5SZLv6KkTSPA1bZtjFSaN9jsQya9nHC2FrSXQloIow9pljgxxOLHTT8eJzDHhclPmm97ghonLRHcjPAK4RnC6PhC7DMxg0ILhZ4p31o+AHgab/jvVggnvahaR+Qp300AO0Hj7LocUzGIwRUSafhJ0QUjYDnmfZDQg1JYcxpP2SUVKe1/YIK67tVjEeoKLhOMx7Rov50+y3BTxSUo4ZFE+1HhyLtNO0Pm1j3uJiPGKZIfRQ67WMXZyMFx4H2Y5rQkUpre3lK8EO+tlcmZEdp+h/rxBYXcV8NlSdLBLHTtJdkImZNfqC9pBuGnaa95CfGhh26QHspTGwtKk17KU60pTTtpaSBbqXpX8qTRu2xI7S9MrH5B63h+FSnyu4fNO1ly65Spxlf1IKFYfCOBlEMODiaQItqMzNEiogkhAk10hNoE/XAoilVtwnFoiCIspLgtFWAhWZJ0GbZFqBZErRZ9u1Yo1Y328SiKBxbiyQws6gKRFcTaBZloTuOlEC2UBcIIzufxPktCgNppITQLCoD8TskEJNn15ISjD9qAwHTBBhDcYTBB7sFFtWBiGkJGEN5YJ1pAgM0x+A0AU41CDfGQ+BU4yFjPBTpVYkwpuxYwLqqkTGmiuUoEpU8SsApqgRY0gQ4RZl0nGGrkQQMgjG3JOgmkCJBMQLM6jYQS1sT2HWZd4I0QXdo2HeY38FWCjpFNyKKAU/ZDVJphmZvRbce6FoxAKpuCDEQJbNRKrohpBvo2QCw7AhJRd6yAWDZE5JKt7IBYNkWEvZxoHWXiDAODM10j1QpDcO6DUGXMPGg8VzgJ+mODVtD7FpSHgBmb4gtYaXhj80hFebCJqOwO8QmC/0BcNkgcirbxRT4I5YC4rOdS/tE94aoA7i6RQTCV2sAXHaJdM8ZmvaJ7oH4VUxmvIin6EIQA+CyS+RAeGhtnx1GZNoZAFeI74Hw0IwPQRVdJt4AuIKCYWuretM30tzkdNUEA+CKKhgEPhgAV1Awuogw+qHj5Dzyj0sAnSb2MJWmfRSMLqlkAFxBwQTaS2zMivh+doE5BXAl2MnpAssGwBW25EFwaNpnT17XWzEAsOBH6HIrBgAWPImg3oGJtI8voauvGgBY8CbUQFEXUwR/QuMM1kTGA49CaWcAZMGnUNqbQPt4FUpHAyBLxEfFtEgGQJZIlFQNMAMgSyR6iyVSTYCfyJ4thqPo1qpE+IP2BoAWtuyVDgaAlgh/5CcDQEuEP+jMZrZIhD/oajz8JPjDW7EGvJYEf9DOANeS4A86GA8/Cf6gowGsJcEfdDae8Uh9/JAhNiklwZ81DtMYfhL8QYsBuSXBH7Q3TvmBP2i2zdlqhT/oZIBtISSjNN4QNPyRjyUOzfxWg8wA4pL7/HLWReAnx5HPGoROg8YigcZfpD5WFzTyB83mLgeJ7n845yzZb/RA3uLRxWa5825hFhd6eowDtz8upj9UF0x1uH7JVAe2FFNdNdVbU9ni9c5U7031wVQfTfXJVA9GFVN9NTVYU4OYGpypwRsO6taA88gZEbCtmBqqqdGaGsXU6EyN3tQYDKd+a0ymRjCxmBqrqcmamsTU5ExN3tTEQZNoakqmJrC0mJqqqdmamsXU7EzN3tQcTM24osnUnE3N4HA1tVjD4dFanKnFm1qCqQUnNZkKlip+V1OrNZU96OpMrd5UDrjUaGrFhc2mVnBf99WJw4DEViMhlh0IjYxa9mKAYIvZYAFfy4EdC+xasNYCsBa72gKtFjy16hWDpBb4tFjJFuC02AkWiLRYwhZwtCCixTywYKEFAC2oZzFmLXhnATmL7WqBN6sHlwAyi0VggTALblnMTwtiWUTOYm1ahM2CShYxs/3wDjWQKktgxoI5NlJDo3IYlRYdbwm/WGxKi4a3WJMW1W6Jtlh0uiXSYvFILBrdElyxqHKLKWnR4RYHwaK8LSEUiw1pUd0WS92isy02uOUAgcWCtJjUFlVtMSAtOtpiOlpsT4vBaLECrap2NRpQyxZjiaPxGCX8oUalBmajRRNbDEbLnBO1wGjgD4ccmPN+qok5JxCBvucPHj9zroYiUQZRC7Gfh2DO1T4U5rwfDGDO9bwUoQG0Kn/0IAVtMOf4/ahH/tAGc64GYT+RxZyrOSjMOc46uos/1GDOccLRRwQyaIM5V+MPhxrFwh9qMOd6YEvtPNxj1AJ/qMGcazBQmHM9H6cnXoQ51wNuuLEANX+owZzrCTv1QvvJMOZczTc9Iad2mzqFarD1IyTMuZpr6p2pnaZ+l+Tyx19++cX8o05zcgXhi6c5+w2NX7lRM05sckx7e0XB+WRzTznZv+DO1RDIo+sfd29+HO7ECCfWP70U82S5fb/aHl2y6QlHrxwJ+4sjL1Yfd/eu3nLsmgOokD3TnlrGQPMv1m+vFkRFOn30frIfbTjNnvSE8fLj+tOz4fd2u3ukc3z7bP3j+ma9ubpZTFFokZyjFz5evl7NF39oT+neQqAFpZ+en9+s9EYOODsS92zrTRaO5j9eXb3lPpE9tRxy1jmYq2pfOG5+t9ruYj5uvi8yN8/B5Ff/Mj3U6fwbevgf/zI93E/QXzmH95fbo/tU95fbWSj0CDY31FimF8+e9zXxcLv8qV+76PTT612/4kHxToybH50YFz2eXu8e9vP2vRiH4FlEuoyeXu8GJNCJp9e7R3oHbRR9tB73Aj45QU8BTTxb77hIN9MvNpsLPUBPQr+i8mBztdt82N6Muwv3doOdO4h5b7djEStI/QoWuL8SDFgrdL9fXKKTUNyIOFXq66uzr7fbzbjCxspWUovT1KMPV28GLJAJeYRikGMKyeXKyyhM/yFHYdY95NF8P169XV2dHd+rgbueeoSwvOiQOLc93z3kFTMs7AsyfkMUzeL8G+5CrG7uYPhIfX69fMNdAG17fyHwqA/724AjDR735W5zsy86J89F7zSt5e72+yjxcKnom/UNAnnMD0m8b7CTLN2ey80N99GZi86po+Adbij1ZH21vvxw+R+r7eZw1YOMW3c0FeL7rZdn29X5avvbx4fSPf1o4HrCcTfh9Dj10M+e+nB1/s1iipZJ26f8sJjy7ZSXC/YBjoq8GgnPlsfC92x5S7ZofJ90aFmTPr2b+mx5dqvrjN2z5dmnd12fLc8+c9312fIMYX95GJ6R8upWCrpx3FWiwfWb9+Om0rPldb9U+nKAxj7h1YLtycX58zfb1erq0fKNog/sAWtHww/JOjgSW5KO52OudbR+KAJ5WD2kzALUxaqX2V6Coot0GvQiDYmKGbUMXhRicE47az/M950o+Q2EaD922/X1w9Wb9eXy4mZ/qUgheVg7bm8yHPVOC9zpnqYd9w8h0cSjDu7pvapSPFX13O9A72vtIVBze3fqvjT9md9Pdzr4bH53dbXafkf3KMlS09feLKY/cJHo5IQ/4k4Ep/8knJSZTifhRBKUsxQIo8iB0qTjdGdvF9/nFTvyij0J9iScBC2pbf+D/vwRcFktz1ZbtLXemdJh21OP1rtHs9DEITR6E4tZ3OeoLOmovVleaGVm/39t1lckzobAg+X1Mflifbk3J3MpVUpQxPjd5fLtihftAf7B8ursYvXDu/XN+9X2u+XV23Fru6ff33wcaX32eqpycnRN8/frzcX6ak4d1xd70Qfr7ZuLu2g/srhfCtNHCvAlJvfXH69fHps9c+Kr48RXnys5J94qScEny48P12/1qwII4dPt7t3mwfJytV0O9PkHumuKKfOnDu6aOKrjvuSugRsM0BFgQc7YM3KPhu+OYvzUCevX7j+9sDxu3X+a0W/df3JVuV+7/+Se8jHUwF6nh2r5nxrib+Thrw3w8IfXVz8+Wd6oI6irdDEtHjx5fvL0enV1ttwtp5N/u/zwb5cfTi5HoaN56QbHf3Nm/nr3+Ns3qwsURBT94sHnXaPh237p8xqj8x+18ozdi2lx1NmTP/x29fs/kvkv7SDD/+t19+0hXuoN/f7EfXwCYYvzl70I18LPH623N8PXebycn/gyiZv13OXq4frm+mJ5dIcd4N1DN4OmQYHDnfInm7PHy9eD/hWn/S+bt5/vzpteaf7nimbsTYK9eTGCIL8SzThMFijSp0qnp0/V/8hEfSH28JdN1J/uTtT/Q+vpn2eK7oRLjmIpVq2fr6922zWGZ/Js+yzOX9x8uOTzBp4dokG6O/RH6OxcOrVJgoawR8mPbjGJL956l08zn6lhKJ4sP+I4LqavRjxleJKHhG8328vZTUGqRzSkf43i/PmHS3hQCLr1uSMCC3z96NMg5q994OJOSOXZ8scV7oZixOe0AuNwpNvo0X9Ts2FzyGku1vviXSyOLY2Fmhz+NFQbbHEhuZQTbY/P/EgMrlo2INh0mL/z40sVqVG8xOQVu6V/AYbX858kWz2fAuD9c1aOmhdrsqF/HciNrDpqVZt9KEjEq32W77WcL4GzCbMx01GMDvzuilAYT8/fLc82Px1ZaQ82W3ye5dmaT0swx/uZfL2lbQCuf12EbzR1q+OLES7kjcwR6qbfd/28faj7EIgeLsfjzdXbFdpLsCYPX5fqwen11V/5vZQ74vR4uVt9pDefEyUYH72kyGxWHVtaX+zzqfUJho86Po/Dsc+77zeZI24wAgZ/sRAzno/5eExfk32C8bGebtdv11cHtZVy5T/JOcRYP9nf+VvGYiDSyb+fKCSpK9r7/Nl9j3+yAYkllyC15MR3t/4uA/Jktbw6Ofn3k5N6WhO+///Xo/F8d3bycPXjiY6I2NNCiORfdkT231rS7yvd/qOfX9K9hmFyP72+vSvxuWDJZ8Nuz5bbFZt+f178Zrs6X0wSf7ndsjbVv7w01nY/hKJfYlK1e399xau/3m6f8m0o8qGf/rjanl9sfmLbgt3B7RbL//bJlcSRFj25om4ZJ1c4aM7hGj3rwLERySkbVwI3JfTYhO41SzZVxKTsTQrOJFtMspy5EhOrNTFxXquaGLKJyZtYuCySTLJikq0GfRhzz09iTfLVpGBNEj7AkfXfWKA5m8Kh8RxNjsEUKyb7bLK+IJgYudPNizjAJSZy0KpGE2IwIVkTcjAhJBNc1TRfgwkcfsre+Mi/YnzOxsdovOXQUjU+eOP1gBOHnILx3IHgAHiu/d8UDYfW9Hs6euafw0n9kFX/fo2eoNcT5nrYUI870wmO/Y7j1tAcFOAYAicJOOzBM2cJOGui34XhCIUeEtAPw4jhZE3iVkrm+CHXWmrhUkyufNoiRb7GEjMXeKLjAg4DUElhBDixxo2dwNExrhFxDSjw4Y0gfFHDJ+Oj8UG/mREN3+Rwhu5VQ2c51Gf0aCJ/naMA1zA4SGV0BPQ5Gg7dk0BZjn/oA/cGNFFra3I1nGswyRTDiRCexZtkOC0CwfkHctxI4VwuWZxO4bH0J86VUIKDJyYa2o8mUdubYLyp4281TilngsnaZoYPw4dwgikmaTm+ncOTGGpT32kNmKAmKeSQT55oPs+9NdqIXJUzfpShVC/f02CSVuY39HxKkE9r87vseKeYuH9bf1cY5XpO54L3if4vaj5056O/l6GlDeoE5bGXn/tAGjX6O6hDbdqjRO9np3nmM0F8NqiXr6MX5M+jx2X1qGXobzR5/Oap6jOps1ZNYb74PonJfd47Sa1sPBQvKzpdac8aXaRJ2IMVJofnPjw9vefBzKFreo9e//ZJOnSll+pvZRo7fZzf681v0C8n/QOOz7y+GHZ3/2qsmu1fCsriJGFafsmw/Tv6SGq36rdQs/jgSvBOiuuuEEHZcOqLDzb7XKLP+m1PPCQ5Ta7UGlwoGI7q7/QvobokPrlYJQVPJ/cekis1plJTjJI5zHnkIVXva/E2xsSnKLsXNjtPHLjjhF4oLsfuic3OU+0eUsipeOXsOPxLx/5GD2meqeE77IgU/TV2j35Mk/k78sgOnsKXPSQ9R/DP4CH9auCZnvwreAdRSgrBZj6Laq27a4OpUcae7fqabz2rh8hnGy/6ZxuvrlZvdoedYf3uZ7fW9MOJv/wXfs6pmQ==').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779221388024', obj, '');\n",
       "});\n",
       "\n",
       "      }\n",
       "      const servers = ['/static/', 'https://root.cern/js/7.11.0/', 'https://jsroot.gsi.de/7.11.0/'],\n",
       "            path = 'build/jsroot';\n",
       "      if (typeof JSROOT !== 'undefined')\n",
       "         execCode(JSROOT);\n",
       "      else if (typeof requirejs !== 'undefined') {\n",
       "         servers.forEach((s,i) => { servers[i] = s + path; });\n",
       "         requirejs.config({ paths: { 'jsroot' : servers } })(['jsroot'],  execCode);\n",
       "      } else {\n",
       "         const config = document.getElementById('jupyter-config-data');\n",
       "         if (config)\n",
       "            servers[0] = (JSON.parse(config.innerHTML || '{}')?.baseUrl || '/') + 'static/';\n",
       "         else\n",
       "            servers.shift();\n",
       "         function loadJsroot() {\n",
       "            return !servers.length ? 0 : import(servers.shift() + path + '.js').catch(loadJsroot).then(() => execCode(JSROOT));\n",
       "         }\n",
       "         loadJsroot();\n",
       "      }\n",
       "   }\n",
       "   process_root_plot_1779221388024();\n",
       "</script>\n"
      ],
      "text/plain": [
       "<IPython.core.display.HTML object>"
      ]
     },
     "metadata": {},
     "output_type": "display_data"
    },
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779221388031\" style=\"width: 800px; height: 512px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779221388031() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(33907,'WkwIPigAc4QAeAHtnetzWzeS6P8VFWc+Irx4P86p/eBHPJldO3bFycTe3NQWbVE2N5KopajEman877d+3TgkJdm+zmwek5nYInkar24AjQbQ6Mb52+y/tt9fLM8XZ8vZMPv83uL828Xll8sXT88XF5ev19uZmZ18cb76n6vln+/PBmtmJ3dX20t9evziv5cvt4TPSPb4Yrtan3fgP1bnx7MhmNnJrqThb2/D9S4EIUQfYjSzk4er8+W99el6MxtcB59uvz9d7sEvV8fb1wo+WJ2e9sQQCzgltlayL0+2jxabV6vz2WDnhHy2evX6RtDd9Xa7Prue7PP1xfWAZycriPBmdvJ8/3hHHyn42eV2sQVLa6S5Bt1RiMwPNouz5U26CbtR8V266xXaJZ2CpyLJftgyFHl3vTlebp6u/tpb7yDw0fp4qf36zOnv8/77zGvbPu+/z7brOy8un6zeLE+/6Tm267eBs6G0TNWvZYiVXn2+z3I9YDZ8pAFfXMtFb36xz3MITmj+ci2DFvKXfZbrATs0U1W262e9NlKV2yA9aa3zKfvqoi02JFfN7KCA599oQ/UCAKFzAmfDR5Tgbc7OpWJ9LC3QFs9O1999ev9eb/VD4NkXFxJBpz4/eP5yF/rJ7unOi8trBd15cXmtrDsvLvfZ7ry43Of84s2ZDAea+Pv945uzxRut0Bff7x4/f73cLmZDkIq9XvWnO5cXy5fbzxbb1Vpr8enV2YvlRp8/X7385s3+8Xt9fLh+1QMfrl/tw/6qsU8Wx08Wq3OGj5md3NusLy9fL1a9wB34ZN1F0SFfM6YU3rP0o/Xx6mS1PJ4NJ4vTy6WZnfxpszp+cx38fg/eeXF5b73eHKT/+Hi1XbxgwG83VxTwYPVmeXyt3lPRTzars9V29e3y8pbMe7i6RKRO4raDi81mNnz1tZmtL7Y8/GBmJx+/Wb68nA3nV6enZnbyqcrn+1eLUxXRiM3PV1vouRH66dXZk8XpcrudhCYN+OnyzfZ26P0/P33y8M7z2TD74/RoZif311cvTpd3r05Opg78bLldrM5pvV73Z5ervy6/uJzin18HJfaz5eJ0NniQS7TCzs8ZMl+uzo/X332+vmBQHcLPD+EuvPYJPlkiqjtLfDcN+3uvZ31s31tst7ea/M52u2+yZ3eX2++Wy/Mutq9B0pwPNuuzz9cXs8HN4aNnx4stUlGA5xPAtHZHAfeDmX3zaP3t8vHF4n+udtzxzWdL2uR64Mknq1evH1KFPkMJpy62L19PDfvN09fr7z7+dnm+fbpdbK8udwz5zZ2r7RqW2KV8tDy/urvYKAzD3HkJy+1ynHy2XBw/Pj/9fspx8uVq+3p9tT3kzoljP1lcdn6bQg5TfXVj4v7JlgfIvncuD75cvhAJsDp/9a41Apxx73RxedmHB+l0UXIYcAGbzuzgUzL9M7rBGmvs6CWUpzBYibVjlDif0ph2eeyYe0rCS0/Lcx1qM945U+PYeK7ReJdHB8Jo+md0bnDNG1e9ccWPzg8lG/kbXRic9aZ/RhcH54vpn9GlwaVg+md0eXClmv4ZXRm8jaZ/RlcH75vpn9G1wUdv+mf0VhK7Rn47encd9INrwbjmjMt19GFwNRtXAcPoo+LNybhkR58GVwMxxqU0eqjqREY7+jK4lIyL3rhgR18HF5JUsObRt8GVZLzNUDIGO7jgjaDOYQxucKFAgXHZjsEPLljjUjUUHILiqVXz0lRe8aQwhjQ4V4zkj2kMeXA2QoGQHMrgbDOSIfgx1MH5KO3qbRtD27d6SmO0gyvWSHWbG6Mb6GKhOaQx+sF7Z6R9AGmqqkQGO8Y4SBPTaz6OMQ3SppDh2hjzQDEtmpbHWBQNDRXLGOsgGJ2QOcY2uFiMVML5MVkhQlhsTG5wMRvh2TF5aYcOhKHV6TkOJU7PaUhues6DL9NzGcL0WGHsXkwbrIkwcBqzHazJVZ8ZOs3psx+olAyrMTOCXGgawyCCOyR7AqhFgQwOp0iyjKTQsQh2GaRuzKAXwJexgF8AG8cCAQCu2rFAgQApj2Uaw/RvgQKJsXUsUABQ01iEgJRMdmOZRnIoY5nQuzCWNiiFxo7VDkogz25IUg2e/VClGXhmCKt0sWOFK6UZAHT4gs+OlWGSJ4Chq1TZsTJwpRkAGLbS2saOrQuuUAGcSiQXAVR4uYb4akFiXM4AUYHoAVSKwVR2bCrGmpRcJFVxhFd5joXnJs8+GDs6awfnnIlJhIKzTjg0ehP96KxXodDoXjc6G4ZajBOuITpKPVzKJgCmwftqQjEpjs7mwSVrCgO0jc4iMrKkdZmiGJ7ehGC8B6RJTKZV3ehEtibGkB+dQ4w5E5zxQF5aFWGU2+iQrTzz8Wl0CFd6k4anXJeUYuVnh3S1zrRmnLWjc2VAwFEfV0eHdIXvcjMOid4GVxuI3ei8HVxEZkXjMvIdOabtVgH9INlsMbWNzofBB8FjR4dorYlpwTDMHbLVBSOizJM4y1D3rirNSNeDKcBbNzTHkGxu9HQJIlE/o7fIpmD6Z/T0iaM+8hm9TSpvEeujp1NqEBFgR2/LwNgszhRKroMPTclyZfS2DeRCWKfR0yVwAZlH79xQnNHMo3de5HLPPXo6BXxkHr2L5NTMo3eJnJp59PQHEt0Vco/eFaV2Qlp3lROsTeid0DLdIY6U6NFLl0Ct4qVLtK6uhtFLl0hlfYij93FXW8FMp/TaQrRHiu9b1Zepvpq57ipM5jCNYWSbjWOYhjHMxLzFSO6zuMI6mCXe+jEwnqd4gQ/WBALrqCZ9cmNgYPfkAtbdCkJAZJsuN5Ibo90vjOwYra4JyA7kdzmB9usUoMh0wowyRktrqPSiflEGs8h+qU9kMAMzG6cwRljH+l39orMKE2/9GFksQSIzqMAeXL160QWBbBT6ndBBUmrjZJIjs0CTqLdjdCLqtV4OUU8WwployGDH6JlntEq+d4muGcfIcqhXiSktMmh7laTK8EevksYjSbxUWeNZffQq0QS+TlXS5G2qk6QOFhAEkjg4AbXXxhimSU+yBpn1tL/HGGTi3WWUmXeXj+aQThsTLNl5hFJSc7uuVnjf9VCU2r7zNT5qfbSOY4IlFY+QnFqGZO21MbUyQVKW1F0rOKYmVQeyfsyWLpBZccxWV+nkyZY6S6eN2VJj6bQxW12ts2KAC3WECYkZ8aUjrMNVqty7bczMKcRrt41ZRJhwoKKEC5UDpUqZqUVGmBbnghRHFSGQmUWYsCfW9piq5Q63D7mv2pF65NTR0SvHrAI3anN4YQR6bcyw5FQ/8iGxev3gkozI0hHW4T7x9+7PIrO0ftLUskqXEdZhnd6UQbJvU/XAzSpdayeQjs+eMuw3UTloo7AiyazOGWPyrM2hz8IaPJZJECrZY+mCsDf7WOC6TmJyY4HpOkjmtKuvRO67n+qVvrwREsbSFzhgksRaO42skwwU6VCRgb1dkxurVUkgnTLWLgRFSFS763NiOjkqe2rnxk7vWCduFElT4cXOe4Jj19isMVUiK4rOhCQmIVNpFyukFEnYi9wLQmK6ONLqVubRSSBbP1bmUWV7SVx3kODYySIWucKAoAfos4O2zliZQie5TrEISbhNhXFFRnaQzH0l3Enq3NhbYaydGydMOlaBIIk9Y0fEYrvLLiUKZux1Y1hUtoy9iRgWNWgjKd+MdS8goTjIfAFNmrZL7T5R1zDtDzTxXmYLItmiTFnj4R6lRpHYlCPlRpHYOvjGGmWn1PkkTgKbesY9OVAX9+RoMTtyqFjckaOJJ3IoKE3UCNBbSKfosbJd7P0iFUk7NpJyExR16Vx119gbdKxJKeoNPtYERRKLKiQJRXsQirp8rLKH7PNBzTqhabeNVfeR2mljlY1kF5xVN5JK7VizTu/KCmPNB5Ms+GU7uS9V22eqtewod1l1SznFla4KmQgqffabSCrK2FOjFGVsldi1KGNPFS26GNq1QmFnM+2wa1HGnhqCjaZGwlFVGVunv1q126akVXtNJ41apc90lqhV1kCdHHacwlvauew5u8yUNqo69CccOvR3OGTodxStzz2SrUvqqQW6oO74++JApqnahXQHVChOvdll9K4/u4yeitlNQMLQXUZLQa2L6A7sJTTt1g5ENBzcENJ9LrZ+bF1MK5bW5XSvSZvkdO/q1uW0Nmbra4YJ0hVD7+nWBXWP7No14SmIYrtBw+tnbAhL+Fg/Y/PSgzqHNh+HFNhWu0BWYfPI1riMzechWdmBRzb9ZYiNjbqLFFqHWA0bz5jH5tsQs0nWuFjHFuwQEZ/GJTe24IYYjRATxxb8EIMR8VTGFsIQvWH7ndEkxCFag8Yt+7GFNIRmUkMXNraQh1BNVhVhC2UIxbBFLmCoQ0iGDXABQ0NtkrNxJY8t2iEEk4txpY0tuiF4dvPsjVv0Q5Atn6tpbDEMvhn2wrWMLcbBV1ESNNQaCbVVSWgLx4YmLZtSjGtgKGgISjOugYHNqqlwshtbbIP3puq6rokezcjWr4wN4dhMTeikxoZorKZm410YW0J7YGo13qWxpTi4bBraqzq2JBtCFGHejQ25mEzzHWLIi3rCUybrO1QuCsm2uWWFMusq00qHYBnT2OejqkFngC7Cuza2HIAcigsBhRjhWgGTxoYemxVUxUXLRUEqWceWK9VylloCopA1qCYELKg+DapS78rYihMQdZiAqHwMKiEFpYWc8x1ES2zQXqACb0jGapzoGgDRLqi+UMCiIHiJZTNvHFoOWhvJCAgZaWzVCoimQECnIGQQi+5alB4KBgUhI44NhRyxVB8QhaxxnuoDZgUhA7AoCBlhbAhIEkMGYBMwQEYYW1Oq0D05PzY2URXtYwe1rdA3OTc2pCTaacgA1B4MkAGoPRggAx2cMJQLkAEoHCXaVQGFpVyEv9vY2NJH49BiWXRpKO2AYfgKjDIKgTHBaNWNi6gjJD4MDg10Nt4W0qPBRCM9wUnhOsGoZYyLDY046YvACXIEZlZB7BhvE/GiTnUJeoCdHVqfigV0AsoxBrF+QJmuuhJnXRAQWiQ2Cggpgdg0oNKBEgEzYIYQNI4o8NAnTmAVEDIktg1o82gVVIzeSl7RDgI6ECHOGm3qPURy/NFoMlnqirrVNZqMta7vMHVitYsOEukE2aLE6zCEynJXjzkatKAygnLklxCj7C8yEvRswBrHNsZV8AdWLB0GP5swjpbkVAIYDaPqFCv0BHSxHYaewOZQ9aLoJ22QEwHRk6omtyB4HerZSmuEiljm4AYJ7iziHfzA0IOAB7/Ke2cR8eAHBj9CHvzMAOCPgclC9Zi0R4xMJQpDT0xMNA55XsAfM9OQq7Q5+GNhklJY8FemMD2toj1iY4JTmPZIlulPjt8y9CTH5Chq0Qz+5JlIFQZ/XxPLGRf9keKQyM+EDn7WxcW4ysYB/CkPzJECC74yZOghXvDVAd4RmPqnNnDEhd5XVOXZokFUGHzZDeisiY/gy3JgqDD4chhoG07oIvhyHDhaFBh8OQ3gEhh8OYuiW2Dql8sgWnvWEYKvDhz+AYs6PrehcfoHDL5iOVVSmPZlWrAQjLwDIRMDGxwJAKOsmWkikYgExIFpRAPAyfTAoKlNjgmdZYJAy0qAF6wIFtASIGg5FgAtAYKWczvGqJXDNmeZKDxHiwTQs0wVsK4EgJbJQk4cSQFvMV3AjJx8oue3TBhwIykcaJkyYL99ADp5Ti1JQfszbcBgEgBhTByx08FpgWXqgMUkBYQxebA6kwAaaNLBSQCUihKut4eDUtGH9PYQSS+TCG3KiQWkyzTS21RkOROJcB4poJSpBLEkAVDKZKJyhqNGAlgQ9L4lQA6BRBJU4xCzcgzE0JYATl2mkyAJkBMazh2mMjhKYU6Rg4hqkKecBgk3A3O2IkcPVK0aBKocCNEQtXBQ6uRASA7Qi0HAciLkGvxWDAKWMyHXqFYxCFiOhRznOMBySoQaB4ImWDVVAkOwHENQJZZ4wOxcoCfLiQ1HQ5xakh4BzNkQR8ICQx+HQ8LMlUNGx+kQhyzUB4HLAZEX3q6mQh+6FCQ+x7ngR7vXWR2BK0dESPhmDQKXUyI5cwYGP9o9JH5zptBe6FNkIDiDwOWUyCPhgQU/J4zwtDcIXId+DwkPTPugVJFhEgwC1zHBcLTVgtGDND94GTXRIHCdTDAwfDQIXMcEI4OIRT9wGnyA/9kSAOeBM0yBwc8EI0MqGwSuY4KJ4MsczDoX1HaBPkXgumgHLwOsGASu40geCQ4Mfs7kZbxVgwB27CNkuFWDAHbsJKLsDkwCP3sJGX3NIIAduwlZoMgW0zn2E6JnsCbRHuwoBPYGgezYUwgcTAQ/uwqBk0Egu8QelaVFNghkl9CSygLMIJBdQnvLSqSZCD2JM1sWjk6OVl2CPuBgENCOI3uBo0FAuwR9xGeDgHYJ+oALh9nOJegDbiZAT4Y+divWIK9dhj5gbxDXLkMfcDQBejL0ASeDsHYZ+oCLCbRH1vaDhzikdBn6rPEsjaEnQx+wM0hul6EPOBgv9EAfMMfmHLVCH3A2iG2HSkZgdkPA0Ec8K3Fg+rcZeAYh7or2L7YuDnpK6vGMQeDcYVYkwOwXyc+qCxj+A+ZwF0Oiu1cn2JL9UYzzZg9O14tt8DMzOxXrsZTM7NvZ8FXz0TTP1i+b5pEt1TTfTAvWNI54gzctBIOhZwvJtJBNC8ioalpopkVrWnSmRW9aDKbFaFpk84iNCLKtmhabacmalpxpyZuWgmkpmsYZeMqmJWRiNS0107I1LTvTsjctB9MyhibJtJxNy8jSalpuphVrWnGmFW9aCaaVaFphK5pNK8W0ghxuplVrWnWmVW9aDabVaFplk5pNQ5aK/G6mNWsaZ9DNm9aCaRi4tGRaYwtbTGvIfTlXRw+DJLaiCbGcQIhm1HIWgwi2LBsswtdisGMRuxZZaxGwlnW1RbRa5KmVXTGS1CI+Latki+C0rBMsItKyErYIR4tEtCwPLLLQIgAtUs+ymLXIO4uQs6xdLeLNiuESgsyyIrCIMIvcsiw/LRLLwnKW1aaF2SxSycJmVo13yAFXWRQzFpljEzlEK8ei0jLHW9QvljWlZYa3rCYtU7tF22KZ0y2aFsuOxDKjW5QrlqncspS0zOGWDYJl8raoUCxrSMvUbVmpW+ZsyxrcYkBgWUFaltSWqdqygLTM0Zalo2XtaVkwWlaBVqZ2WTQwLVsWS5YJ2bJutEzFtpGDZaNlJrYsGC19jtaCRQNfGDnQ52rVRJ+jiGC+54sdP30uC0W0DE5WiGoPQZ/L+tDR52oYQJ+LvRSqAWZVvsSQAhz0Oft+pke+wEGfy4JQLbLoc1kOOvqczTpzF1/koM/ZhDMfocgAB30uiz821EwsfJGDPheDLVnnsT1mWuCLHPS5KAMdfS72cWLx4uhzMXBjG4ug5osc9LlY2MkuVC3D6HNZvomFnKzbZFMoCzY1IaHPZbkmuzNZp8m+y5X69Q8//GB+LmvO9D5rTnVOeI8zSbfYxGR7c07CybpZQ452BdzwigA88Hy46fSwdwdxOGHc9gd5tNh8s9wc+JdowEGRPWDnM/H58s32zvkrTK8xQAXUSDu3tIHEn65enc/Qiih8UD7RD9ZYtmexMF68Wd22E7+z3d4hHNPt49W3q8vV+vxyNiQHRmIOCny4eLGcfF7AJ7BiiGAQ+PHJyeVSnFGQsz1wR3YQulcvv3m4PH+FK42dW4ycpQ+mrFIXTM5vZtueTqbmuyQTegyTn/9maijd+XfU8D9/MzXcddCP7MO7i82BK9HdxWZiCjHBxjmLYXr65KmOifubxXfqgqHw44vt3t1Dge7xoUB3+nh8sb2v9vbqToYRPINIhtHji20XCVTi8cX2gbhf9aQPVt0v4JYFPQkk8Hi1xYdsgj9fr0/FgJ4AdVe5tz7frq82l91v4c62k3NDYt7ZbhnEIqTeIwv8jxQGjBWqr45GVBIIj4i5QB+fH3+82ay79xYjW0BJDqoHV+cvu1ggEvBAigH2LiQW95eemPoD9sSMe8CD/n64fLU8Pz70sYE6DT2QsBS0D5xwT253FDGJhV1C2q+zopmdfIIvxPLyhgzvoU8vFi/xBRDcO1+4gzrsHOF6GDTu0l2nZpd0Cp6S3kAt6W7W+yBw72D0yeoShjykhyDK6+RkS7WndBNibZ0p6RTaE96ghlSPVuers6uz/1xu1ntXDyKuuSeKiFevlyeb5cly86eH+9QaftBwGnBYTSg9DN3XU0PvL08+mQ3J0mm7kC9nQ7ke8mzGOcBBkuc94MnikPmeLK7xFsh3QXvMEnTbLfPJ4vha1Wm7J4vj226eTxbHb/H0fLI4htmf7Zunhzy/FsLc2H2VQLh6+U33VHqyuFB/ymddaOwCns84npydPH25WS7PHyxeivSBPMTaQfMDMg4O2Jagw/6Ych2MH5IA7kcPIRMDKVtpms0ZUnSW51EcaQgUmdHwxAISEcPmVMEvJ38n4j4BcFKP7WZ1cX/5cnW2OL3cORWJSO6rHb9bMhzUThLcqJ6EHdYPJpHAgwru4N1UpYsUFiDq/rvLtROBkker03apqc9UPtVR4bP+8/n5cvMZ1SMlQ02KvZwNX+FIdHTEl/NHjk3/UTyqE5yP4pHLQN6SIPYke0iCDsO9vZ58F1dtj6v2KNqjeBQlpeD+mb6+RrgsF8fLDbO1+ExJs+2gB6vtg4lpUmca8cSiF3cxwkvSai8Xp5KZ3v/39eqcwGkhcG9xcQh+vjrbLSdLrc3VKBLjz2eLV0sK2gn4e4vz49Pll69Xl98sN58tzl91h2UNv7t+08O09zRUKDlw2fzLan26Op9Cu+uiJr232rw8vSntexS+phB9MAE+Y8n98ZuLZ4fLninw+WHg87elnAKvpSTho8Wb+6tX4lAPEz7ebF+v7y3OlptFlz63t2tPFsc/mf8dPfZO/7sni+N37dfYSdFCBxILcBI+Pfag/W7MjLd3YbfF9Vtk9e9e+XZua4g15hh8rhX9y0zc9D+y85pyiC223LILmU7Ab9/PfaotVBtKSuJsjhd/nOeWa0yuJuuyp5TrTv0f5XmxPtuYak4FXYOkmJz8P3JxnjPx0dasbvDTHQAz3PvmOfiQQs3epwKN113zbZ5bl2ryzpckdThw07du7gqRPuAkr7kF82z4qKZ5dKG24lsI0VH29XsCyryl1lqr3tuKmvA91waE6ue2WmttSEluJbhxhUBpc5eJR0NWIeXwPoES58ET610OmQE8xXK5wO3YGzcNzG2oLvvcaratVEtd+lUBchPB2/p6ipdNgU25hYKGp9So9O9RcPlAnJcUErrJlCoKwR0CYt/KMx1Bv6sA3VLJLefUgu7C9jcMzK2lnVtrpeRYfYJH9hcOvDV6us7AzmOj4Wg7H5J043S9Addl3IzTewvsPEpnWWtj5PoFZjK9N2HebsVcuw7hrdRcuyDhXSnejXp3hcLbsE9XKsyDRcPcWso5oEjfXbHw0Vvj9M4FP7fRheZzTS5VGZ/9CoY4d9Eln3INxVXh2R95IwMkvO9GBhHfciNDf/r9RoZrNzLcUAZ8+I0Ot+5w+C93qOfc3+2g4R96j8NP5/vfV0i/iOs/vPXOpYfssN+1+EDM/LSLD7nv5i0Dlftv3jpKdWK9NUR1Zr01Pg83OlRb4b6xvb3A+1+u7tAF0HR6AdN7m/kTd/99jdz59Sk3u2yuzg559enVi8vl9mh9cnTv0dOjz67Oj7x19i5JDrpGNR7/y5Xhj9fPf/pyecoOVQ76ZifvULB37fq77rbqlX8jmt1p9zgbZn84u/rD2dXR2eLy8uirPy3/8jWRv2kVPfS/WHG6IGeks5NncgMQ4c/01h/R4TzTNNxSc/Jgtbns6taHi+mJe8H8tNU+W95fXV6cLg6u1GHvt9s9CqdwbrG/4ubR+vjh4kWH33Nu8GEd9/3NjpNbVf6xDlR2WomdhqOfw7znQOWgt3ZdRTm9q3j6xTvqHccfH9ZRf73ZUf9MA+ofpotunNgcHOdYUcB8fL7drNB95YDlyezk88urs+8EDOz2BfQ34DfAxfs8t9lFOUXvKd9wU1yowQZf5iXUkGiKR4s36K5nw0f9SKcrs/cBn643Z5OmFK7uBzJ6OdbJ06szaBARdO2yQc42uHvw9jnq+1ZnNxZyTxbfLtF4iox427Tw02s9WHi4eak2hBp8qh6rCt3Qh3mMJZQcXcwFs4vddt6l6JvFBgK7h1lfdYSQai0phuZLQq/zzPWdUdGdkcu2BW4jovwpqiTZaaWWbcxyE52/vp9yzZYQq1zit4sKmsuHGjGPnFY0KsWowJ/POY3j6enrxfH6uwM90b31BrXr4njF7Vb08a4nX2zAjYDTy864IVE6A73j2w/c4Tci+2k70vGmqnl32r4/C+9az4fr81dLZi+5xW9/t6Oej6/Of+T1bTfY6eFiu3xDbd7GShDea0kS0b/OhutLrXdWem6D6GMOaj41xKHefVdxIvvZRT+0OFij0UnvXqHRoA+5zE4HpfYw3PV4s3q1Ot/PW7mIJsCVElNqt2xM/p7G6CLp6N+ORCaJOvyDWeHXbhCfXEje14hnRfA/SYM8Wi7Oj47+7eiozVvm/OFfujWebo+P7i+/PZIWcXZeOab5zbbI7u5Hue/x+pdcByn2Dn3N/fjiumXE2w5smMl2q+3p6O/JYrPE8Ohvsz9uliezwbUfrmMWVHrzYx/baggrN0HKvHt3dU7RH282j7mrknjgx98uNyen6+8wncBCabNh6X/dejZjVivWs7ozw3xWDAuxGhRjSLEAFZvMaE3z1lSs40s1FVNCTGNLcNxz5Eyq1jAv8cnBmRSyCSWZELBndiZgp8zNZC2agM205YPtcTIBz8fiTcDXkNv7SjG+YvecjMfbE3cyCcsGfL5UE3Bx5Gq/VrWcxpVTVeNwXSvOBFuN55q0wq1QXtO6ZAI21dhhe/J4yc+zb5SfjcfHE3/UBs1N8xXq2EzFfp66Wuy2o0HPHKlrrSbEZkJqJpSg5eAVhw8q/m8463GNUoG2ZDxOcMUaj7sbjoHYbWecDZPaVHP7GLbX5MM1UnwL8VLFH04cRzVefCDwaFTbdvGRwBoV+1LMRjEExcKUXsU8VDxL8FLC+w1vKBwH1b8P3wgxdyS/3CTYL8jB+pF0mMRibYktJ/mxy8RAE2Nd2ATrTLHUxfoS20nCMT8V81VsVnuYx21XzT8xxcbFFgtVgMQtjDXh8wuM6TM8h39wqXJDo3gOYwbtcS9uDjfmgpcCvs5FLsbjIkdxks54M5CRB/3jwhECI44G/OVgIoxLeJQ79nDFlrvCoviAR2tCxv07JhOqiQ5HbaJgEUACrdyNhnG//NGjJsDbRvxcuZYMX1Ks6eXBGbHUl2+609CB9D9lipupXPvFBZAUJa4jeheg0cLx4p1isb/Fbl6u7cJRU0rDZlWesW8lPX/4dkkgXuly+RZOGPJsDfbH/NEXpvAjoHhSyrcnLGk6OtRkiiMVRq6mmGqSiSSt8j/xSEHVYJpsiskmk8BKtJDUJJJkwUSTTDbBFKKkqGiCST2O3MV4E42X0CSpNUU2yQRTBfISAw2ClAeoInkxyXghYSp0SqxIIYLC+QQp1ApCkDTD3Y08UW/iHQh6OghWGEzcGOSl7IlU0pOXOx9JKTfxCEwF5epIyUU6jdNqTumJ1f9goXQ+lKd49YlvPtCltBELRg3RfJpH8/N9WMJU7hQr1z/1NPtytDzyQavm11IpjTzUTVtRLofqeScqpvxKKzgJmconNRAlaH00VmlVahSPptBnLVvzkkZL06cpp9JKy044lOoJs8YrxUkouBlPnGKaMEy0ar9puXxDG//5nXJRvuLQXy2L2mk6pWvfYlPqKVxxTm0CbYclaiuRhvKI41k5diqT3+nZCa8qF/KslJEnSW4NoSTFRLyGUYLWLwrng0N7m3FAmmSqcDwDk/9Rnhj8PPOf8cCAZUzqyCBtllFDOOsQiWA4UTV+VUpADhem8ksxGg9REK/f5NAGmWAaZPqvTUiKqQF1AJFfK7uvsqaYwrWEfb6bpU84aSDSKh37XHsKpifSkHbCDeMd/p/onn6J2z9Dp0Jad43TsrQUZ9zP4Gnw4rTrBz7gaANlDjvgd23Af0JdjmyveVPCPEabasvNhYJDiCpy4jxgHBCSzxhUcMgv71OYN+tq8K22mPDLmRQ5OTdbmg/V1iw2+DtFjsuc/FpckkLxqiiatDVRT8xDi7VaMdaYtDU1pBpLqy3HzMIUPFOupoqcWDhNZt1+eFRFvf5ORc7UUV3FsUWh/WN2Z3LgTfcdKI72+ox3K3LEAPEfQpHz7vMx6vEhyiya/6AJJlANJIXndPf1c6l0Yhbl4KQ2vLlP/OGHr28pNX43Tbv2XplbZsn/1C+M8fPmfBRzIG9Vcsn5eYjzjJjxSYydRDXtZ0OY+6raZ19aaGIjRrBz83qQHGl5zXTsI5+ymydfWkq1pCrydG889pFvucxjJTq3IkbAUzQHnK7Oa2mtppJzFBF6rfQY3Lz40lrxqYUiFk37wqPP80r2EnJJYmI3xXKEMQ+58i83nP3Ie90yLdp0aJv2D2uY5ufJRZt9bK66iKPzZBWGWdpb+nlvlNbm1UUM2nxLxduy/IhDhetGacGleahiVBhaqAfFi1HaW/jlmk1aCTa5FhterjVI9+9emDNP10zS4LUPNki7bXT28xqk3Sb1FzBH83MZEbdf8hPm4q3Ujcx8KfLCkx9pWSbHL7+/6+fGy31ugO9518/NE6QPfldQX2cdWJCpMlzWXYfvB9Lw3y3L3n4qh0j+4KMxJo//r1k7+4I+5g5fptZHGxPeNNQOV94U/fMbiT1hRf4BW6kPsxL774vLlSzxO8+9YwU8HP37//0/f7i4XB19J6+AIstBk/8CFmPyRqqbLt07kzF5rdnvFmN6In24+6JbD5yzJosxWeZ0gzGdW9QQqXP4L2+L9A5n8w+zRfrdaOyXs+773WjsXWOKPc/uJYy/pl3fBxuNcWMPmqVuNHYNxBLsMPqNbPJymkefXImperYIkhOrsZZDC/OmHivcIPe72ditV9iyqghzb1V7V1Lj4ifVNvoQ2KImG2JN3LWm2kZ226ohbDXEMu22fW3z4Elti+d09Xezsf6S4Q/RNor13V7l8+uZje3WXu9ULf7LmYyJyPmNmYxFn2yyzbsSqku3tKt/jw3dzmQszG37bRlI/QytcWgyxl1NLocfcyoxqeB/LSPLm03y69iMRfcr2YyxNeO9vHrvIiYq3EaJdUg3V+HFiob74LA+mgxYsELqNi+8xg8DK2xksPSKYkeDSRX2OEFeqoChDRdsYo+TislWLH3Uxsdi7NMylkFyoR1lc2sdt/dhGcNVb1wVh+WSvDUVQxoMYrgsHuOZILddcsMpN5lyW6ncaMmtmlzFzu2W3D4uFljy3kNNy+V53LbHmwd45iY6vSKwWz1xK6S+l89y6aRYIzmTk9SmYajEGyPEhChE/uSOT2opBmnyLcY9GAbpH8ZaJMNCiGu9ScNrErqpj1xtLg2txj9q8MNV6xKoBj//gmfAft5yLZyCVBujHGOi9/HNzS1+5CFHy32Ju1VZcJ5VnGsp6WmDnIHYOXf32YgHu3hO7M+AWaPt7Pffd+wrCf9pT3oDurp/hLXXj9BzvXt51o96uxvDBP5a80yo3uVakyvNV+evzzPYdXxNEPcarS62l92Zglebn+qrzc/Ply+3+9uTHqzeLI/VmlheLv7D/wNPdZJd').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779221388031', obj, '');\n",
       "});\n",
       "\n",
       "      }\n",
       "      const servers = ['/static/', 'https://root.cern/js/7.11.0/', 'https://jsroot.gsi.de/7.11.0/'],\n",
       "            path = 'build/jsroot';\n",
       "      if (typeof JSROOT !== 'undefined')\n",
       "         execCode(JSROOT);\n",
       "      else if (typeof requirejs !== 'undefined') {\n",
       "         servers.forEach((s,i) => { servers[i] = s + path; });\n",
       "         requirejs.config({ paths: { 'jsroot' : servers } })(['jsroot'],  execCode);\n",
       "      } else {\n",
       "         const config = document.getElementById('jupyter-config-data');\n",
       "         if (config)\n",
       "            servers[0] = (JSON.parse(config.innerHTML || '{}')?.baseUrl || '/') + 'static/';\n",
       "         else\n",
       "            servers.shift();\n",
       "         function loadJsroot() {\n",
       "            return !servers.length ? 0 : import(servers.shift() + path + '.js').catch(loadJsroot).then(() => execCode(JSROOT));\n",
       "         }\n",
       "         loadJsroot();\n",
       "      }\n",
       "   }\n",
       "   process_root_plot_1779221388031();\n",
       "</script>\n"
      ],
      "text/plain": [
       "<IPython.core.display.HTML object>"
      ]
     },
     "metadata": {},
     "output_type": "display_data"
    }
   ],
   "source": [
    "from ROOT import gROOT \n",
    "gROOT.GetListOfCanvases().Draw()"
   ]
  }
 ],
 "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.12.12"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 5
}
