{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "64a4a0d2",
   "metadata": {},
   "source": [
    "# rf109_chi2residpull\n",
    "'BASIC FUNCTIONALITY' RooFit tutorial macro #109\n",
    "Calculating chi^2 from histograms and curves in ROOT.RooPlots,\n",
    "making histogram of residual and pull distributions\n",
    "\n",
    "\n",
    "\n",
    "\n",
    "**Author:**  Clemens Lange, Wouter Verkerke (C version)  \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:29 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "e2e2d044",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:18.850238Z",
     "iopub.status.busy": "2026-05-19T20:29:18.850120Z",
     "iopub.status.idle": "2026-05-19T20:29:19.832937Z",
     "shell.execute_reply": "2026-05-19T20:29:19.832173Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "38bae8bf",
   "metadata": {},
   "source": [
    "Set up model\n",
    "---------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "3f6875ce",
   "metadata": {},
   "source": [
    "Create observables"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "8984c1ac",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:19.841487Z",
     "iopub.status.busy": "2026-05-19T20:29:19.841324Z",
     "iopub.status.idle": "2026-05-19T20:29:20.001409Z",
     "shell.execute_reply": "2026-05-19T20:29:20.000718Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", -10, 10)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0399b256",
   "metadata": {},
   "source": [
    "Create Gaussian"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "356ec63f",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.003757Z",
     "iopub.status.busy": "2026-05-19T20:29:20.003582Z",
     "iopub.status.idle": "2026-05-19T20:29:20.135122Z",
     "shell.execute_reply": "2026-05-19T20:29:20.134377Z"
    }
   },
   "outputs": [],
   "source": [
    "sigma = ROOT.RooRealVar(\"sigma\", \"sigma\", 3, 0.1, 10)\n",
    "mean = ROOT.RooRealVar(\"mean\", \"mean\", 0, -10, 10)\n",
    "gauss = ROOT.RooGaussian(\"gauss\", \"gauss\", x, mean, sigma)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "aa838862",
   "metadata": {},
   "source": [
    "Generate a sample of 1000 events with sigma=3"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "cb6540fe",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.137004Z",
     "iopub.status.busy": "2026-05-19T20:29:20.136821Z",
     "iopub.status.idle": "2026-05-19T20:29:20.276522Z",
     "shell.execute_reply": "2026-05-19T20:29:20.275915Z"
    }
   },
   "outputs": [],
   "source": [
    "data = gauss.generate({x}, 10000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "6d48d6b3",
   "metadata": {},
   "source": [
    "Change sigma to 3.15"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "b7d46f7d",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.278641Z",
     "iopub.status.busy": "2026-05-19T20:29:20.278484Z",
     "iopub.status.idle": "2026-05-19T20:29:20.386518Z",
     "shell.execute_reply": "2026-05-19T20:29:20.385849Z"
    }
   },
   "outputs": [],
   "source": [
    "sigma.setVal(3.15)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "775663b1",
   "metadata": {},
   "source": [
    "Plot data and slightly distorted model\n",
    "---------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "3b6901f3",
   "metadata": {},
   "source": [
    "Overlay projection of gauss with sigma=3.15 on data with sigma=3.0"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "5d7569ca",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.388702Z",
     "iopub.status.busy": "2026-05-19T20:29:20.388537Z",
     "iopub.status.idle": "2026-05-19T20:29:20.587820Z",
     "shell.execute_reply": "2026-05-19T20:29:20.587380Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x56131acc8b10>"
      ]
     },
     "execution_count": 6,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "frame1 = x.frame(Title=\"Data with distorted Gaussian pdf\", Bins=40)\n",
    "data.plotOn(frame1, DataError=\"SumW2\")\n",
    "gauss.plotOn(frame1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "414616bf",
   "metadata": {},
   "source": [
    "Calculate chi^2\n",
    "------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "b7d1fefd",
   "metadata": {},
   "source": [
    "Show the chi^2 of the curve w.r.t. the histogram\n",
    "If multiple curves or datasets live in the frame you can specify\n",
    "the name of the relevant curve and/or dataset in chiSquare()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "5c35030f",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.600491Z",
     "iopub.status.busy": "2026-05-19T20:29:20.600342Z",
     "iopub.status.idle": "2026-05-19T20:29:20.717039Z",
     "shell.execute_reply": "2026-05-19T20:29:20.708525Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "chi^2 =  2.63197794362232\n"
     ]
    }
   ],
   "source": [
    "print(\"chi^2 = \", frame1.chiSquare())"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "db952825",
   "metadata": {},
   "source": [
    "Show residual and pull dists\n",
    "-------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "38e1c128",
   "metadata": {},
   "source": [
    "Construct a histogram with the residuals of the data w.r.t. the curve"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "af180581",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.725451Z",
     "iopub.status.busy": "2026-05-19T20:29:20.725304Z",
     "iopub.status.idle": "2026-05-19T20:29:20.842952Z",
     "shell.execute_reply": "2026-05-19T20:29:20.842309Z"
    }
   },
   "outputs": [],
   "source": [
    "hresid = frame1.residHist()"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "19c7fc47",
   "metadata": {},
   "source": [
    "Construct a histogram with the pulls of the data w.r.t the curve"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "69bda4d9",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.853496Z",
     "iopub.status.busy": "2026-05-19T20:29:20.853341Z",
     "iopub.status.idle": "2026-05-19T20:29:20.962508Z",
     "shell.execute_reply": "2026-05-19T20:29:20.961860Z"
    }
   },
   "outputs": [],
   "source": [
    "hpull = frame1.pullHist()"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "6f84d312",
   "metadata": {},
   "source": [
    "Create a frame to draw the residual distribution and add the\n",
    "distribution to the frame"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "00e90b18",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:20.964288Z",
     "iopub.status.busy": "2026-05-19T20:29:20.964157Z",
     "iopub.status.idle": "2026-05-19T20:29:21.072670Z",
     "shell.execute_reply": "2026-05-19T20:29:21.072022Z"
    }
   },
   "outputs": [],
   "source": [
    "frame2 = x.frame(Title=\"Residual Distribution\")\n",
    "frame2.addPlotable(hresid, \"P\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0e9468fe",
   "metadata": {},
   "source": [
    "Create a frame to draw the pull distribution and add the distribution to\n",
    "the frame"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "cf6190b2",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:21.074270Z",
     "iopub.status.busy": "2026-05-19T20:29:21.074142Z",
     "iopub.status.idle": "2026-05-19T20:29:21.293940Z",
     "shell.execute_reply": "2026-05-19T20:29:21.293288Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf109_chi2residpull.png has been created\n"
     ]
    }
   ],
   "source": [
    "frame3 = x.frame(Title=\"Pull Distribution\")\n",
    "frame3.addPlotable(hpull, \"P\")\n",
    "\n",
    "c = ROOT.TCanvas(\"rf109_chi2residpull\", \"rf109_chi2residpull\", 900, 300)\n",
    "c.Divide(3)\n",
    "c.cd(1)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "frame1.GetYaxis().SetTitleOffset(1.6)\n",
    "frame1.Draw()\n",
    "c.cd(2)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "frame2.GetYaxis().SetTitleOffset(1.6)\n",
    "frame2.Draw()\n",
    "c.cd(3)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "frame3.GetYaxis().SetTitleOffset(1.6)\n",
    "frame3.Draw()\n",
    "\n",
    "c.SaveAs(\"rf109_chi2residpull.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "018a7b55",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "6f955aee",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:21.295656Z",
     "iopub.status.busy": "2026-05-19T20:29:21.295497Z",
     "iopub.status.idle": "2026-05-19T20:29:21.479528Z",
     "shell.execute_reply": "2026-05-19T20:29:21.479007Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222561470\" style=\"width: 900px; height: 300px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222561470() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(48708,'WkwI/S0ARL4AeAHtfVtzHDey5l9h9M7DORFgbWbimlVPsmSNZ1e2FJY9lo53QtESm1KvqW6eZsuWZ8L/feNLoPrCJmV5TNuylg43VYlCFRJAIpHIW/1r8mz94/lsMX09m/STr+5OF99PL76ZPX+8mJ5fvFquJ25y+vVi/t9vZn+7N+nJTU4/ma8v6tXD5/939mKN8gmqPTxfz5eLBvzv+eJk0ns3Od28qf/XVW1d14D3QXwIbnL6YL6Y3V2eLVeTnhv4eP3j2WwLfjM/Wb+q4P352VmrDGQBjpWJ7PHZ6frz6erlfDHpqUPJl/OXry4VfbJcr5ev96t9tTzfL3hyOgcS4ianT7eXd+olXvzkYj1doxVV1NmD7lQID99fTV/PLuONsksd39Tb79Cm6lg8vhKP744MXvnJcnUyWz2e/7ON3k7h58uTWZ3XJ1z/fdr+fSJ1bJ+2f5+sl3eeXzyav52dfdeeWC+vAid90YSu7z0g2YZs+8h+waQ/rgVf7z2F2fx6+8wuODbz970H6kv+vn1kv2DTzNiV9fJJ64115RDETBIxc6LMQUrMHEDg2xc8/a4OVHsBQOA5gpP+GG/wKaeQKZYiPgpI48nZ8ocv7t1to74LPPn63G5gUp/uXH+zKf1sc3Xn+cXei+48v9h7153nF9vH7jy/2D759dvXthwwxD9uL9++nr6tHfr6x83lV69m6+mk99axV/N2defifPZi/eV0PV/WXnzx5vXz2apefzV/8d3b7eWP9fLB8mUrfLB8uS37Z737aHryaDpfYPm4yend1fLi4tV03l64AR8tGyvapWtQV4W3JP358mR+Op+dTPrT6dnFzE1O/7qan7zdB3/cgneeX9xdLlc79T89ma+nz7Hg16s3eMH9+dvZyV6/x1c/Ws1fz9fz72cXBzzvwfwCLHVktw2crlaT/tt/uMnyfI2Ln9zk9NO3sxcXk37x5uzMTU6/qPx5dcqkz168mstqdjE/OX9zdgb++dV8DcSuu/3Fm9ePpmez9XpkoxjSL2Zv14el9/72+NGDO08n/eQv46WbnN5bvnl+NvvkzenpOKVfztbT+QLj2UbjycX8n7OvL8b7T/dBu/vlbHo26QWN2+0Kpy7hv+wmp9/MFyfLH75anmOp7cJPd+HG0rYVPpuBgTdC+WFkBndfTdqKvztdrw8m4s56Xfc69O/JJ7P1D7PZojHzPcjG9v5q+fqr5fmk5w7U9eRkugavNODpCBgvqAD/5Cbffb78fvbwfPrfbzY0892XM4zLfuHpZ/OXrx6gC23fMvqdrl+8Ggf3u8evlj98+v1ssX68nq7fXGzI9Ls7b9ZLEMqm5uezxZtPpqsKg4zuvAAhbp44/XI2PXm4OPtxfOL0m/n61fLNepdmRzr+bHrRqHAs2a317aXt/MaEBuz51woN38yeG1+YL15eJzmAMu6eTS8u2qJBvSqq7Bacg1Qn1EuMrv0G7smRo0GsFFe+J7tLQ7B7EuMQN8/QkFpNlOdWF9elL+qE2ZUwKK5LcMJpYDQYXPsNzD2rOC7iOMvA0ufk7P+Bfc8krv0GDj1Ldu03cOw5etd+A6eec3HtN3DuhYJrv4FLL6Ku/QbWXoK49huErDIrnqdBeB+UntU7VnacyiC+55IcF4B+kFDbTdFxpEFiz8XjjuMYBwFWDclAg+SeY3QcxLGnQUrPPloHSxpEe87RCSVgMnjq2YuzppMfPPfsMzBwnGjw0rMnx7E4vNj72k4p9VkMldR2oh987Jmzs+dDHHzqmQIwMJR97pnU2QNeBl96lmDjKqSD1+2oxzgE6jmTs+4qD4F7TLHh7OMQpBdhZ+MDEENVKpKehhB6G2LMmoQhxN7GFGiwDiH1eI0Gp2kIuTaDgQp5CKW3FtnQHIL2HLKzTrAMkQwJI7Ehcs8hOaPZIYqNQwN8r2W8Dn0O43XsI4/XqZc8Xufej5cFhN1eoz25AAKOQ6KeXCr1GktHuV5Lj07ZshoSVhB7rXewiEAd9ngEUHIFEtrg2kiyleRbK9a6LVIeEpo3QPKQ0b4BFIYMBABwoSEDAwNiGvK4hjG/GRjYHSpDBgYAShyyIRCjSzzkcSX7POSxefZD1r5i6Ggo1FcEcc19tG7gWvpiw4BrLOHKXWgooEobBgB1+aI9GgqWSRoBLN2KFQ0FC9eGAQCWrY22o0Eb4/IFAFeOxAFAZV6sYF/q7Q6nBCBUIAiAysVAVDRoZWNqb85WKzPKi12HjGu1a/GOBibqmdmFaEyBiY1Cg7ggA5NUpqCYXh6YfF+yY6Ma3A7WD47JeYCxFynOZxfDwJR6juQyFqgOTGAZyepywquwPMV570QAYkhcwqjywMZbI9aQDMxgY+w8OwEkNqpgRkkHBm/FNX4SBwZzxWxi4PFejhXjSs8M7krsVB0TDcy5B4NDf7gMDO4KukvqGBxdey6KhnlgoZ4DeFZwnMDfwcfquBWA0ttjlF3RgcX34q0dGhistURsCw7LnMFb2TtjZYLKyZa6cKk4g7vubAFC3CtjSSoPgikBS6y/QQi8ybv2GwRzwuiP/QahWPkt2PogmJTijQXQIJR7rM3MLuPNpRevFS3Og5D2eArMOg6CKQEV4OFBmPvMrj48CIvx5fb0IJgUtIeHB+GAJ+vDg3DEk/XhQTAf4Oic8fQgnCu2Y6Nl0zlrVQ3fsVlsd2BHFelBbEqAbW0XU1L7ysUPYlNinRUfBpGw6a21jElpvQXSAi6+HVXJY3/rw2XTYTzsxzUM3kZh8OMyBjFh38JKbrt4hetitvskg8d6Hu8bvCMTGFxXNepHHjwWdqtuYNlIEAaCt1VxI/IQaCsY0RCoygR4HJBsngS0lVMABWwn2FGGQBiNyr3Qv2CL2Xi/9SdgMQPGbhz9EEA6JJv+BaYK4z7JECAsAUXsoAYL2mrdC+wNomD4s+GBqugN2yaHhw0aWT0NgY3V134xWD0eQTk2GjxAQxDsM7VL0qakyoxDgDjUuoQtLWDRti5Zl0EfrUv1PjiJWJfrfUgfrUsYAiljl2p1HftktT0BRANW2bOBddaG4MdNzx71tuvV+R6Ct41386DtvJvnMBw2aUMESTYawVui8maqK7ydemAUdTv59X6o/al9HCJIsrZjKEdNQLnO2hA1j5C9y/peOzhEta4DIhkSYQpsVxwSVSkdzyRCn23ShkTosU3akKhK65AYQIV1hRmKCeyrrrAGF+tym7YhYU/B/TptQzIWZhRYmwQVVgq0LiVsLbbC6uvY2+vQRSCIncWIsFWu4zF2i3ePD6lJ7eB6eLKujtY57CqgxjocYoSAWRsSSHLsH54Dx2r9A5UksKy6whrcNv42/cl4Vu2fDbVJ6bbCGly3t0ogSXTsHtqGlF57Z1Bdn62m3x6ikq+DAokkQTrHGrPrOhz12kgDl3lkhBXtITdG2IZ9yKC6hmLkIYPoGoiH46a/dnM7/ehebuKNoTDkJuCgJatce1dvlpEHGnco4IFtXCMPhSonsEkZSmOCxiQKbeYcdxo6lfeURo0N36GM1GicpoAWG+1ZG5vBhoxZOXJtohEhKqMittLGVlDTOGF75ZYR4k5jR7W7BfvoyJBJhoJ9tJK9VS4byNrY8CIIuUaAaB5A2x3q6AwFW+jI1/FaMElQW2XGBTyygXi4ScINpUaNbRSG0qhxbKmuVUBACWfG1hCE7ca7KlIgxtY3LIuCI2MbIiyL4usgVboZypZBAmNv+wVwqnUb124bdfHj+aBW3vJsa8iOKOOjYfeMUoJxbLzH3huMY9fFN5RgJ6VGJ2Fk2Ohn2KID7MIWnfqaDTroWNigUyuP6OBFccTGgDZCdYseCo6LbV6sI3FDRvbeCIwady711NgGdCixYtQGfCgRGNldqEKiYbQFgVHjj8XOkG0/KKluaHXahlLPkXXShmIHycY4Sz1IVmyHkur2XklhKGlnk0X7dpzcvrWOz9hrO1FuHq1HyvFebqqQEaHcdr8RpVwJexyUXAm7cuySK2GPHc1VGNqMQsbJZjxhl1wJexwIHDTrTVBUqYRdt79S6rSNVUudtbpplGJzVneJUkwGaujgxGm0VScXZ87GM22MSl36Yxt16W/asKXfmtC299hjjVOPI9AYdWu/CQe2TZXGpBtQmeI4m41Hb+az8ejxNZsNyAi68Wh7kTYW3YAth8a46Q6LBgUrmHTbi0kGbWy6tqKNT7ee6Min21Rr49N1MLXJDCNUJYY209oYdbvZtGtGU0AKxw0MfP0NCmYJOq6/QcVmsO6hKqGPHsdq9njUyDzgaJwHldRHshN4wKE/90FxUOeAl5Y+FIeDZ0iDivYhuUiOQxnUUx/APh1HHtRzH4IzZMKgXvrgnbGnPKj3fRCH43eCJiH0gRw0bkkG9bH36qJCFzaoT70vLlUVofrc++xwRM5oofQ+OhyAM1pQqE1ScpzToIF6713KjrMOGrj3gtM8zsYapPd25OMSBw2+F3U4C5c8aAi9FFMSKNQaEWqrHKEtHBSatORydqxoIUNDkNWxogUcVl0BJfOgQXsRV6pcp6ZHc3b0y4OCOaorETqpQcEaiyvJCftBI7QHrhQnHAeNoefkFNqrMmi0AyEUYcKDgi9Gp9IgLHlTTwjeCfkOKpcK2bFZU4US5CqnuUEgGac450NVA50BdBHCOmjygBiKCwMNGaNaA2O969vdVMGquNCUK4hOlkFTQbeY0EuAUMg6qCYMzFB9OqhKhfOgmQ2EOsxAqHwcVEIVtBFilgZCS+ygvYAKXMEZi2PTNQCEdqHqCw3MFUS7uIvDvGNoOTDa4IwAgUYctJCB0BQYyBUEGrgL3bUpPSroKwg0wqBQyOEuug8QClnHgu4DTBUEGgBzBYGGHxQMEpWBBkA10AMNP6hWrKB7YhkUh6gC7WMD61hB38Q8KLgktNNAA2CdQQ80ANYZ9EADOjgjKPZAA6BRlGlXDTSS4gD61kFxpA+OocUi6NKgtAMMgi+AoYwCwxhhaNUdB6gj7L7vGRro5IQy6kODCY30CMcKlxGGWsZxUGjEUT8bHIGOwdhVwHacUMR9U6dyBD6AmXptW7GBbKCZMXBXeijTq66Eib2BwMXuBgOBisfd2EOlA0wMTAATEIHGEQo86BNHsBgINOyu9tDmYVSgYhSyZ007CJDRENiZYkxFgCTMH4ohM1HX1K2sGDLIutJg9AnSLnSQ4E5A25R4DQaiJu5WM4cCF6iMgDn4lyFTyd94JJrHAUxhtnFc0L6HxNJgtI9DGExLZpUADA1j1SkW4OOhi20w8PE4HFa9KPST5M0iYHrSqsnNYLwM9WzBaPgCtgzDDTg4E9g72gcMfMDg0X7l90xg8WgfMNoHk0f72AHQfvDYLKoeE+MRAraSCgOfELHRMPh5RvshYRvigjFH+yFjk6qwtV+whVVrFcYjKDa4CmM8ImH7M/NbAj6RsTmaWjSh/SjYSCuM9ptMbDYuzEcMfcTz2NDRPuTi7Ljg4ID2Y+qxRxps7eU+AR/ct/ZKD9oxGP2P2sPEBb2vqcoTQYNYYbSXuIfOGvcD2ktmMKww2ku+x9jAQhfQXgo9TIsGo70Ue7RlMNpLyRTdBqN/KfemtYccYe2VHsY/wKaOT9orrH+A0V4mWJUqjPHFtkBAGPwODWJjwAHHCtCiycwYIuOIKAg9tpFagDaxPWDRFDUzIRM2CGhZUSDWKhgLmkWBNQuzAJpFgTULux3WKJmxjQkbhcC0iALMLLYKkK4VoFlsFmZxRA3QFrYLECMsn9DzEzYMUCNqMJrFlgHy2xZAJw+rJWpg/LFtgMCsAIhh4wgND1gLCFsHSMxqADFsHpDOrAADNOrgrACYmhKujQcDU9OHtPEwTm+bCMYUFgugbttIG1Pj5dhIjPJQA5hiKwFbsgJgis2k8hmYGlEAgaDNLQrMCGScoDgGmzUzEJa2FcDqMlqCrMAsNLA7jO+AKQV7ihkiigM/hTXIqBkwbCtmekDXigNDNYMQBqJkGErZDEJmQM8ODBYWIVbQW3ZgsLAJsaJb2YHBwizEsOMANisR1DhAaISrpspgIGxmCHQJIh5gnFyATzKLDUxDsFqiPhgwbEMwCRsM/GAcMmIuMDIyrEMwsqA/YLgwEInRdnEF+EGXAo4Pcy7ah3avkToYrpmIwOGVHBgurERmcwaM9qHdA8dXdhnjBX2KLQR2YLiwEgk4PGBrHxZG0LQ4MFyGfg8cHjDGB0oVWybegeEyNhiYttS7akiTXmzVBAeGy7bBgOCDA8NlbDC2iCD0A469eNA/jgSAUw8bpsFoHxuMLankwHAZG0xAewmGWWZffRcwp2C4HKgXW2DZgeEyTPLg4IDRPmzytt6KAwNmnCNsuRUHBsw4SQQ7HbiI9nGWsNWnDgyYcZowAcWOmMw4T5iegVzEeOBEYbA4MGTGmcJg7wLax6nC4OjAkDnijArRIjkwZI7QkpoA5sCQOUJ7C0lEXQA+ETZbCI5splWOwA+wd2DQDJO9wcGBQXMEfrifHBg0R+AHOMOYzRyBH2B1Hvgk4IfTCjnwa07AD7A4sGtOwA9wcB74JOAHODowa07AD3B2HuOR6viBhmCk5AT8yAlEY+CTgB9gduDcnIAfYO/E8AF+gGE2h6kV+AFODmyboZIxGKchwMAP9yGJA8b8qgPNgIlzrvMLXxcGPjm2+1iDgFODIZEAxnkRz0PqAgz6AwzjLhyJPnlzCl+yv5jL3uT+2XK69jJxkzPzHovRTb6f9N+qBKeCo19yKuAtxamoU09OYeL14tR7pz449dGpT049eFRx6tVpIKeBnQZxGrzTEJwGHB7hIwLeVpwGdRrJaWSnUZxG7zQGp7CBx+Q0gicWp1GdJnKa2GkSp8k7TXA0iU5TcprAS4vTpE4zOc3sNIvT7J3m4DTjKJqc5uw0gw+r00JOCzst4rR4pyU4LTikJqfgpca/1amSU9igVZyqdwoHF41OFUfY7FTB982uDj0MODGZJoRggTDNKMEWAxZMEBsIzJfgsENguwReS2CwBLmawFoJ/JTsVAxOSmCfBCmZwDgJcgKBRRIkYQJzJHBEgnhA4IUEBkjgegRhlsDvCEyOILsS2BuZ4xIYGUEiILAwAt8iiJ8EjkUgOYK0SSA2AlcikBlV5x08AaoiKGYIPIcinjCtHIRKwh5PUL8QZErCDk+QJglbO0HbQtjTCZoWwomEsKMTlCuErZwgShL2cMIBgbB5E1QoBBmSsHUTJHXCnk2QwQkOBAQJkiBSE7ZqggBJ2KMJoiNB9iQIjAQpkGxrN6EB2zJBWCJsyAS5kbAVk+IJiI2EnZggMBLmHFoLCA34AycHzHn1asKcQxGB/R5/cOLHnJugCC0Dm4RY/SEw5yYfMua8OgZgzs1fCqoB7Kr4Y44UaANzjnM/tkf8QRuYcxMIq0cW5tzEQcac47COvQt/8ATmHIdw7EdQZKANzLkJfzhQY2PBHzyBOTeHLZPzcDzGtoA/eAJzbspAxpybf5x5vDDm3BzccIwFo8YfPIE5Nw87O4VWzzDMuYlv5iFncpsdCk1gqy4kmHMT1+x0ZnKanbs4l3/89NNP7rfy5ozv8uasIQvvCDFpHptw5F4tUHF0da4lR5sXXIqVALgTD3E5FGIbJMLwvz+MEvl8uvputtqJOqkFO69sBZtIiq9mb9d3Fi/hhw0HVID1JnWEMbD7Z/OXiwm0IhXeeT9u31/C3z2Zh/H07fzQe/zOen0H5XDfPpl/P7+YLxcXkz4yWsSdnRc+mD6fjZEwaM/g2kJACwY/PD29mFmICvhsK9yg7Q3v+YvvHswWLxFgQx3BydnmYHzU+gL/88uPrc9Gd/NNlbF5OCY//dP00Kbz3+jhf/1periZoF84h59MVzsBRp9MVyNRmAs2QrawTM8ePa5r4t5q+kMNzKjww/P1NgikAi0OpAItFOTh+fpe9bevQWZwgscismX08HzdWAI68fB8fd+CslrV+/MWF3DgQY8KVngyXyOybIS/Wi7PzIEeBTWI5e5ysV6+WV202IU764bOJY55Z73GIjYm9Q5eIL+QGWCtoPs1/AidBISIiM6gTxcnn65WyxbThZVtoFVHU/ffLF40toCbAHe4GMA2hbiLoJhWGf0H2Cpj3QPcme8Hs5ezxclu5A2wq6U7HBYv2haObY/BeHjFyBY2FTF+jRTd5PQzxELMLi7x8Fb6+Hz6ArEA1vYmQm6nD5vwuFYGHDf19rHZVB2Lx6qXmrZ6l/u9U7gNO/psfgGC3MUHRXhfQycRuj3WGxuuozNWHUtbxUvYoNbn88X89ZvX/zVbLbehHrixF7RoLL5GvTxazU5nq78+2Nau5TsDVwt2uwlMd0u3/ayl92ann036SJi0Tck3kz7vlzyZwA6wU+VpK3g03SW+R9M92kLjm6Jty1Z0GKz5aHqy13WM3aPpyWHw56PpyRXxn4+mJyD2J9vhaSVP90qwN7ZYJTQ4f/Fdi1R6ND2vUZZPGtPYFDydwDw5OX38YjWbLe5PXxj3AXpgazvDDxDrYIdsUbQ7H+NTO+sHVQBuVw9KRgKqZFXrrF6Di05SFyyQBoXGM7Q0XIzF4HBaUftmjHdCzc8AsPVjvZqf35u9mL+enl1sgoqMJTdpRzYiw07vrMKl7lnZbv9AJFa408ENvNmqqpACAaQGBW+e2rBAe6Z2Rze10Z/x/ehOZT7Lvy0Ws9WX6B5qYqnZay8m/bcIJDo6wh+WI8ah/ygclRFOR+GIEyAhVAityhayot1yof3qm3uF2r1CR4GOwlGwmtb2b/TnH2Aus+nJbIXd2mKmbNg20P35+v5INLERjUViYRY3d4yWbNReTM/sYcz+/1rOFygcBYG70/Nd8Kv56404mUtRLsE4xt9eT1/O8KINg787XZyczb55Nb/4brb6crp42cKYa/kny7etrM5eLTVMdgI5/z5fns0XY2kLX6xV785XL84uc/t2CxGoQHpnA3wCkfvTt+dPdsWesfDpbuHTq2qOhXs1UfHz6dt785cWZg8ifLhav1renb6eraaN+xwe1x5NT24s/g4zdm383aPpyXXnNZykMEI7HAvgyHza3Z3xu7QzHp7CLrFrDMYVzPo2WP+YQ0f2H8eiISKYFsH7xypdoEiRxXuf7NxnwfxSg25TphSg35pYbH/x3KUQAyUSiRnsfT/Un2PqfBQiiiVCb28Vxsh/Drnziru5xIiXjmkBsL93UbKqqk+xkCG4H60fpEviVUuIWmxr2c0Q4LVTwt0olBgUOt61oHopJSQlTZR8MLF4P3lA6TSi8ayapXhQ6rW5BKRQlwNCkillQTcu5RVI2kmx0U6UFV3ZTTKQUhcT7rJwVGA63kXGgcO72+wBSD+AmSwROj3NofgETHfSEVw10eNtOxFo0BIKeSKfS0Yk7fb9SEdQAneUQ1CBGjHsvB53r6KXneQFvgspxRxL8CEle3ibcKAjwghjkHMKRSKoY5t/4MrbY3YD6jwIFP/FmIrRzpjtANkzKm2TiI9sgeI1WULnuT2VC/mEBmtOA+q0ThBRCLnqMfayI1yJzV6+hOtqXN/0JqPCVa23DAvHJvhekWMB4kdLspC9dhLIx5CLwng4+aU5F0CztzkXNjltfnHOhUsH+/d//vosDc94V3l5RRqHWuF9MzXcXGR/k39+l8B+0OW1goWdn68TLbC0b1a0QI6buhr3stxgHSLNzRWLcPdkgp5UuJ1EDyWyXymO3f/6i7sYjZpH6Z0j9xnfu27comXbaGR5igRAz94+i4k9T1+8KM8Zx4oq50/6yb3penr0w3z96uhkfrFertazk6O/Tt9cXMyni6Pzk1PU3ZmEqrn4lRLeL9ezf/FidoaTZhVmrlaxNh35dXmr2ni8tYd3BuAtgD+1Yh34P5/DJhBAyU8sl08l8yc1g48V1ypIOHN6f766aDrSB9PxCim+ZDwfv57dm1+cn013suPgwLY58hlZwNiwzVbz+fLkwfR5g9+h7H+/efrx8jxZKpSLo//zP4/+44i6ePSfH9a8dRCc7fC80VE0S8o7TCLbmcNarxNXp3DMwfT7z9o1Boz3m7V/Xp61D2uSPpIpumRz2THIkKlQPl2sV3Nor6zDF29e/9CWNS7HJY7rMT+XXbcbn0/fQqt85W7YNM7YRKsh5ovl6vWozgThNqtJzWt1+rg2ZyxnL08gDBBIG3ho7Hx/sasmptqqPMdEVdbaJ/MFNJafrlYPkVwLqAF++P1sdXq2/AFWHZySViswuH3HngSPH3PswWaD7E/m3XNZNrh56/k51srP7/xfLpfQ8aPyVVudOa1sUoe9evYSezk2eTxg/GnST/CC5cvV9PXR8vRoU+PZ+dly/eyZbYg7Wz6I6Fds+Bj7d2/4kPMuG9a/OF/OF0g5Wbe0d0wSVeerY+1ydMfaSXTHxa6LXWe7znad7DrZdbTraNfBroNde7v2di12LXbNds12TXZNuG5/4FID0GrZA/asvcbeaC+3dqxJa90QMZwMPcPUkDb8rSvolWWGe/rzQ+BddMkh4CY4Ly4KApvgeGPBH8jBoU6UnU/eBe8d5MSIoK0Av8vgUkwu5eJSYhcLfCfFRcoupOyglPDw34M/jDmemL+QK3DJdD5Y6hO4AbngKsI3scY3dLpJiXfAgTbsamRJe4zmU6gqfnbkbA4nCxAb5JZPn7yav3z1i5+C2uNnH/qWu+yFIhXKMRXEk3XiE6WsOQfoqpx0IWgomoPk4rk43yWKMbLkGJJXdaErWjRrKDGm5JOLXYIGJUjQoFDR5U7gNyUxkkqGv10HLyvvY1YO8C3siipq+5SaE1GnnHyOUnKxiKCOcw4hlVxyhItY7ELUJKF4nwPctnJHsWQR+H6ZO3BHUWIspXiJGUE6XaGSolBIqSBiSLpCPlKJxFrgLeo7zl4ShSixJkrqAqfoiWKQjHeEDvo7ZVUf4TApsUssQZMGzR7hZKkjX4L37CUWkGjsMqkKSQoeCisJnYSS4eQVNSHyyXeZUvTKRSQidZB0kVOCbqx4uFMKdwlOhJJL0AzfLu2E1GfJLAKnRy4dqw8UvZaUMByx86KZNHPKEe5b0pUSSEuWlOEvy9QplcwC9zBv7n1QnXoKJcQAnzPtvGSflaKnUoq43KUCj7IcckkF67IrnjRyUTzkXexigurNJykeXpS+85yShKzkYwzOd5xEck6JOBV40B5QXM08+fS9iP6Wfm/p90Ok32c/4oRiJ9JJf+VZ8dkmx+2zmtc2U9jjHG7ybLFcvf779Ay5O8nA1/PF9OyT+aLZr8zC+GzxeP7y9dQEw2ezjagNdbWbPFtNf9iI38fsJs+Wq/lLvOUbS1hqEu+uH/mOuNkkmQ9oG79xBZclrn0vGffum9X312oGo0S/q+QyEfYZDiPfvv3HrqD7aLVEtvL5crGRdHF7R7pNJqBeL9/ClvNOk+XPy7cQoK+Vb82A8Q75FvcRXHAMYxr8Ne0CfyHsQuKFuAuZ1x1n/JvdRsyNuIQ8nNxG3I24RFl0G7HXLlEW3LEH5N2x4N8q+YrbCL9o2Gqy2wjBEZeQhEdZGEIx5OFRIh6l4lEujs47NBNGoTg2wXgUjXGRIcs2ATk2Idmhg8Whu+rQeYsGwL6KSxuen5OVMZybM923oYO3upaSqGCPdgcFqcN2notq4oCgHaYuU4QLvCZWuLCnLvnMELSoBMRJx07Ywx6aOCnkN0j3JQfRoJwgWPsuFY2SQinqo7jkuxJ89CVGFaSJzbELgcXjmeCRlKB0sUTVHJKUbDESvvMpq/cSJJhDuXYlKrzxPUd4ovvSwXNfCyP9OHLxlI59DioCmQsSnHaQYYovQaMlYvUmfHpfOGZIW1K6qD4mdNgjKDzGTrwE5aApIPS/+I5L8hojxrI4H6QjjKiq18TJBZIONtqSsxKEq5CkKzFEjVExsC4od15TKuw1REb6ttIF8ZSy5wz3+Bh8F4lT1JAK4lliSp1wTJCQCW74saSOcyTVRAHBDYl8R5JzNCtjEpc4dSFpChQoIa4lCegrUMg5auTgkpeOVSmnkgOc9JMPkPIkFUkIw7uixsE7Dlq5jMcBpod9Oejt5fE4GLHDMb086gfzcjBzh3N7MPuX6eOAgg5p7IAKD+j0MiUf0PrBajhYLwcr6mDNHazKg3V7sLIP1v6WO4j6gui+Df+4psDY0YdwEt6Xj67eFg8kpAPy35eQblw6uPPkb48f3/n8U2zRVwsJf1nNTic96403/fwM3gXXtrv3uZRH0+/NixHVr1LB3aBLFISLZrQcfX1EOXj4eMCEmVPokqSsookQCjUxV5/SPH0Sjn0I8mquPpI7ksixYC8zh5on3LwqzM/q6Qip5JAFHBxKBshAT6RVLHDIejpCqtH8NEJOxdsnPnYtpkD/bwu440MGevxqerL8YcdR7O5yBb/L6ckc6e1ReaO82cxHM9utYUbBgMPx8OqIm/FmDbcRc6HYd6WtTosIt9kGwzS3xwfLxcsZLGEez20/+VIDZOaLX/hVh8vOBNP1zNSqV5ELEG+9tA5aT9/PMPu+gzGOza4v7ib0CDebP3NzZN4Rj9+t/MVYPcBnL6oNoCrjQS0P7dCzsYRFiB0SE8SITOGnzUcv8KELBKL94yAU7da3cUuEmIQDx/aP+kNEV7m8mW9jTp0vnlgiQdoDL3qHbyMzd9DqEVSICNC97NsYYuyISzIFndbXbb8d9Kt8GzmxffYnpCLRWwjd6L2IDwbF0knxuGuiMnhqaxjWtS5wjlADh+zNMW7fs9FTBttPKYUE0fvD9WwMvhPWEkvIrKlgBkbXxeb4WN38ti6s4+338mxk9l3MKlx8KojDHx83x8YriGXHsdEOciFEohBjQVTw7jeZOh+amyGCl6tr6uYbTVf6CW79GlO0r9ukVJSQ5mf7Faffxq/xENXf26vxOOYu5pKSek8lVJJtX5LC/tK8HBUHv0TwFC3szVn0F35Zynb12y9LXfqU1CXwHV+WuiyYvPeXqZqEcpUTI4JiNqblayvcejlefVr4rb0cr1qatmViXcLr8YpFuSvDg2X/ub0en3Oa+ZM9r0d8kuvkzfTs6N78Yr2aP39jH7G8rBmu4u2vCmZ5t+cDxvZAM9xcHZEeY3J66+tYOctOMGDzdUQikdFn7qNzdgQpfjg+qbe+jXCauvVt/D2chn8f38Yr9rytZ9EVsuxH7edoO83GKHZ5aG7eCPxbODral1Gf7bg7VsfHq2zCm71/eXq0cTDbmIfhKjk6Qh5NFydHV+vK319Hhj3+VxqQbx0kb8ZB8li6wJ59USmRQlGYuzPMll5DYnOlO85dykgKplQi55BQohIpq1eGJTCjRHwsOQaKwSNB3DH7LuTAHmoOpHQ7ZumEQogacgwxmNldNYpS4phigiVfpBMtMDYSh2xPJagzAsOjLJeYYYIvMSWBxt6HwsEdo06CbbhoQcK6Yzir+RxJskQpnM0xQAKsyykKHO+og4teLJ4iQRkFZ9DSaREJPmafkWvuOHaKlHSJAxzlkNqfukQcVTwLleRd6UKAWSCGIir4CGTocinJS4xwQoOr6OUCGLK3z2Qne++0LKQdwZvNl1gSUh1yvISZD/vom5PCXv+4dD5J0AyzBVLLHgtdGiSf9wfSHXu9PNaYs0vzcTBncX9a2R0zX5p62qcOKZiQyyTkD8gsHJAibchVYyCpjno3bE28pBD4jO9D9rzKThEtKnXUA+AA9Qfx2z84Bu3fPZgJBTDyNoAfdRza/tmMOxggWyhaAw6D0UwPYsFoFj//gUejWdjk8dGL0Uvtgzms4VB8E2FoxynAFTsxZ5VkAfh1BlW7pIEoS4GrjE3t4WwCjd8wtPA2SM3m+AqtyEYnUrUjf1wc4Xsf5Ixi/90gtSuocXuQOyThj/Act01EzSQ/IX5lOwAtSu9Xxsi8X0yM1bqNgbmNgbmNgbmN4bqN4foDYrj2fTwnbnToPI7IbUVZUklK2SPSpYbAlATvFxbx0H0gQ9avCYE5lo4zPKVzkEBR5UbCYW6VNrdKm1H79edU2ty4s/L7+kkHyIOXFDy/MkPRxi/3av/sP9xPOnexiCRNhBhZsCDzG2SizmckjIv45oKdkmTSJwQGq0+MWJwMx1U4RDAiAYRYOatnS427cZOWQCnhGdM1mtvexmc6VFc2r6EUgjpj4zOdo3qPBH+mo0XrH68DtSVP3Pqu/nEO1Bs7z2Ufjz+L13RIGc71nEi9xFuv6ZfzhWWghpvnbUbY5epkNiaiA6+x1HbXZIRNneU+QEYJ9S2cY9JfkxA2dAWOnDD9RETC4eV7WVuzRVtF75FAIWZjpjflNC3UaRSECXpf1Fxnd17NBcGHuKshqfHRrdO0IFYPCWEJcYWQci95TSdKnSi8etl7fEXpg80Hm7mLZrxT5Le1XWx0i74Jr+nQIaMGZ59DKvgS0r7b9CG17HlNE1JVsC8C0yY+h7XnNo2QI/sPcat1g3xft+nDlK+/bTrYQ1R/f7fpjqSUVIizcoQ4cYXXtO+QhQTpSWKMsEj/4tywWAu3XtMflte0/zmvaatw6zV9tTHarIU3mnZ+Pzfs8eHK3HOavmJNflxO09PTQNnvOU0/enN26zBtbvPINbMTx2ppZy7DO7a5P/Kra/A4qObZfaP8jnnQLn93E+41n4J7vzyjB9lhbx2mMYs37jRxa2f/eOzsOC7923b2K/a7rZn5cK/8CO3s20Rd+/7Sl0fm5tW9v4W/9Pmbs7P3c9+zXf/WVdqyXv3/k0v2GBmZvEgKyKgUJCHxV4oaS4gllORLghsue84SU8ZHpOEqjfxJhYRjSJnxSReEsueM1JNFoDMQ1PGF1TJkJXzwGC9WwReIOIUUfEY+MXwvKJUEJ2hKwjUVGbywY+GSM/ypgSG8ubMP0XtOloaMOCg+C5588NEqUWaJWTyjJ3Cf7iSVwFKyF/Xi8ViAPiQETuozPqaMAO8QfAmWorVYvs9OMmmMIWZWJJyypGeRYo6aUw6xeIeeBe+RTop8VORB88gTpMjpWXIJwVGnyFnlg2pCSldBxmBoeRI+K6QBabq6nCwhFaOrimfwBe+SzUEaycCQVzVKQEJYJBSzEiQe1YRMXbBaeMcdXLIpMYLXfbJB5eJDiQVJSCM+6o4EBSoKf29CFhkb08CSMuVAyVPAuKsvEgKe5Qjnd98lKeKVM9J1+YKncDfAWx7u0zY0jAroE4Lg65TmHJMoPqntkUntmDEWWXPJRJHVJpCC+ISsB8FbOlO0rxLxaRzxhOzDaC0hyZkQPuuEd0uHVKeZNXuyjLnH1HEJObEPiYUy598iEdMly9r7uk7/Qaz31msanyaumYjsQPaBH9A+Rq/p/9h1m/5PfMvjf1wgl+qzf51M19OfPqhDHM5TN+JIHTvL/JgKVOK52p7NLzd0jISPMUShXKqd+Xc/ht+e8D6eE56R7L/rSX0FNe6e8A5o+CM84v35XKl597MF75fBfeeRax0GkVRwT8QbHQb9ZQn3VzoM5i4lpP+XQOpDvhmHwdujy+3R5aM6uty4F9/7OhCmfONNf+gOhKmjhG+R5JgZ+oHmQBg6kZig/ZCQ2T6eDHtg6LKYmyCFCKeZ5kEYuuJxSE1ERaEqgX9Oy58qmhPCpaOPSExo72+33uVBSOJDEEQ9a/3A8MfrQWjBhh+CB+GVBs+PyHsQn4J8/GI1P8dnKcyE9tn85auz+ctX67vLxWL2Yt0+PW8pLN7OTmpC19Pp2cXsp/8H1+aFaw==').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222561470', 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_1779222561470();\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
}
