{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "8606bf17",
   "metadata": {},
   "source": [
    "# rf301_composition\n",
    "Multidimensional models: multi-dimensional pdfs through composition, e.g. substituting\n",
    "a pdf parameter with a function that depends on other observables\n",
    "\n",
    "`pdf = gauss(x,f(y),s)` with `f(y) = a0 + a1*y`\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:30 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "8b9d3877",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:28.384748Z",
     "iopub.status.busy": "2026-05-19T20:30:28.384618Z",
     "iopub.status.idle": "2026-05-19T20:30:29.362149Z",
     "shell.execute_reply": "2026-05-19T20:30:29.361560Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "c639ceca",
   "metadata": {},
   "source": [
    "Setup composed model gauss(x, m(y), s)\n",
    "-----------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "cd4c362f",
   "metadata": {},
   "source": [
    "Create observables"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "f384920f",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:29.364233Z",
     "iopub.status.busy": "2026-05-19T20:30:29.364104Z",
     "iopub.status.idle": "2026-05-19T20:30:29.521721Z",
     "shell.execute_reply": "2026-05-19T20:30:29.520617Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", -5, 5)\n",
    "y = ROOT.RooRealVar(\"y\", \"y\", -5, 5)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "16fb2805",
   "metadata": {},
   "source": [
    "Create function f(y) = a0 + a1*y"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "fcf8cee8",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:29.523403Z",
     "iopub.status.busy": "2026-05-19T20:30:29.523277Z",
     "iopub.status.idle": "2026-05-19T20:30:29.704287Z",
     "shell.execute_reply": "2026-05-19T20:30:29.703630Z"
    }
   },
   "outputs": [],
   "source": [
    "a0 = ROOT.RooRealVar(\"a0\", \"a0\", -0.5, -5, 5)\n",
    "a1 = ROOT.RooRealVar(\"a1\", \"a1\", -0.5, -1, 1)\n",
    "fy = ROOT.RooPolyVar(\"fy\", \"fy\", y, [a0, a1])"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "7309fdcb",
   "metadata": {},
   "source": [
    "Creat gauss(x,f(y),s)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "c743b07d",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:29.706185Z",
     "iopub.status.busy": "2026-05-19T20:30:29.706057Z",
     "iopub.status.idle": "2026-05-19T20:30:29.832568Z",
     "shell.execute_reply": "2026-05-19T20:30:29.831700Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#0] WARNING:InputArguments -- The parameter 'sigma' with range [-inf, inf] of the RooGaussian 'model' exceeds the safe range of (0, inf). Advise to limit its range.\n"
     ]
    }
   ],
   "source": [
    "sigma = ROOT.RooRealVar(\"sigma\", \"width of gaussian\", 0.5)\n",
    "model = ROOT.RooGaussian(\"model\", \"Gaussian with shifting mean\", x, fy, sigma)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "af8db8b6",
   "metadata": {},
   "source": [
    "Sample data, plot data and pdf on x and y\n",
    "---------------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "4be02076",
   "metadata": {},
   "source": [
    "Generate 10000 events in x and y from model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "1edcbe10",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:29.834149Z",
     "iopub.status.busy": "2026-05-19T20:30:29.834016Z",
     "iopub.status.idle": "2026-05-19T20:30:29.987409Z",
     "shell.execute_reply": "2026-05-19T20:30:29.986926Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n"
     ]
    }
   ],
   "source": [
    "data = model.generate({x, y}, 10000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "47b6eee6",
   "metadata": {},
   "source": [
    "Plot x distribution of data and projection of model x = Int(dy)\n",
    "model(x,y)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "86e1df5c",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:29.989504Z",
     "iopub.status.busy": "2026-05-19T20:30:29.989371Z",
     "iopub.status.idle": "2026-05-19T20:30:30.155037Z",
     "shell.execute_reply": "2026-05-19T20:30:30.153918Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on x integrates over variables (y)\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[y]_Norm[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55c61076af00>"
      ]
     },
     "execution_count": 6,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "xframe = x.frame()\n",
    "data.plotOn(xframe)\n",
    "model.plotOn(xframe)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2a0e870a",
   "metadata": {},
   "source": [
    "Plot x distribution of data and projection of model y = Int(dx)\n",
    "model(x,y)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "67c9e77c",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:30.156690Z",
     "iopub.status.busy": "2026-05-19T20:30:30.156541Z",
     "iopub.status.idle": "2026-05-19T20:30:30.263259Z",
     "shell.execute_reply": "2026-05-19T20:30:30.262139Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on y integrates over variables (x)\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55c61089de80>"
      ]
     },
     "execution_count": 7,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "yframe = y.frame()\n",
    "data.plotOn(yframe)\n",
    "model.plotOn(yframe)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "07934802",
   "metadata": {},
   "source": [
    "Make two-dimensional plot in x vs y"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "f262b571",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:30.264985Z",
     "iopub.status.busy": "2026-05-19T20:30:30.264826Z",
     "iopub.status.idle": "2026-05-19T20:30:30.430920Z",
     "shell.execute_reply": "2026-05-19T20:30:30.429840Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)\n"
     ]
    }
   ],
   "source": [
    "hh_model = model.createHistogram(\"hh_model\", x, Binning=50, YVar=dict(var=y, Binning=50))\n",
    "hh_model.SetLineColor(\"kBlue\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "b0777a56",
   "metadata": {},
   "source": [
    "Make canvas and draw ROOT.RooPlots"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "7014e837",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:30.432557Z",
     "iopub.status.busy": "2026-05-19T20:30:30.432427Z",
     "iopub.status.idle": "2026-05-19T20:30:30.697114Z",
     "shell.execute_reply": "2026-05-19T20:30:30.696711Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf301_composition.png has been created\n"
     ]
    }
   ],
   "source": [
    "c = ROOT.TCanvas(\"rf301_composition\", \"rf301_composition\", 1200, 400)\n",
    "c.Divide(3)\n",
    "c.cd(1)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "xframe.GetYaxis().SetTitleOffset(1.4)\n",
    "xframe.Draw()\n",
    "c.cd(2)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "yframe.GetYaxis().SetTitleOffset(1.4)\n",
    "yframe.Draw()\n",
    "c.cd(3)\n",
    "ROOT.gPad.SetLeftMargin(0.20)\n",
    "hh_model.GetZaxis().SetTitleOffset(2.5)\n",
    "hh_model.Draw(\"surf\")\n",
    "\n",
    "c.SaveAs(\"rf301_composition.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "367d2c0c",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "12364b57",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:30.699144Z",
     "iopub.status.busy": "2026-05-19T20:30:30.699015Z",
     "iopub.status.idle": "2026-05-19T20:30:30.885529Z",
     "shell.execute_reply": "2026-05-19T20:30:30.884764Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222630876\" style=\"width: 1200px; height: 400px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222630876() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(83007,'WkwIKkEAP0QBeAHtnW2THLWSqP/KRMd+uDdCU1fvL1WxH4zBC3cNdmBY7GUJR9vTY/dlPO3taYN9NvjvN55UVXdV94wxYMBw2udMU6mSSikplUqlMlP/M3u8efNycTl/sZi1s69uzy9/mF99s3jy4HL+8ur5ajNTs/OvL5f//Wrx2cezVqvZ+UfLzVV9uvfk/y2ebkifke3ey81yddkD/768PJu1Ts3Ot19q/+e6um6qwDlvnfdqdn53ebm4vbpYrWet6cEHmzcXix34zfJs87yCd5YXF31mkAUcMmstxRfnm8/n62fLy1mrG1K+XD57vpf00WqzWb2YZvtq9XKa8PB8CRJWzc4f7R5v1Uc+/PBqM99QSynkmUC3KkThO+v5i8U+3qTtNXybb9qgbdYhefgkxcc9wyc/Wq3PFusHy3/0vTdK/Hx1tqjj+tDU/z7q//vQ1r591P/34WZ168nV/eXrxcX3fYnN6jpw1hpTIm2flHBJ+mxXZpowa09rwteTUgzn17syY3Bbz39MStSv/MeuzDRhW8/QmM3qYd8eacwhOGtz46Kx2mvjXErFLU6h0N0HHn1fu6r/ACCIDuCsPdWN1jbmbJLV3mmjU4BULlY/fvHx7b7fx8DDr1/KC4b10ej5m23qp9unW0+uJh+69eRq8q1bT652xW49udqV/Pr1C5kQ9PGb3ePrF/PXtUFfv9k+fvV8sZnPWicNe77sn25dvVw83Xw53yxXtRVfvHrxZLGuz18tn37/evf4pj7eXT3rE++unu3S/lHf3p+f3Z8vL5lAanZ+e726uno+X/Yf3IL3Vz0zGlM25FXhHVF/vjpbni8XZ7P2fH5xtVCz839bL89eT8E3O/DWk6vbq9V6lP+Ts+Vm/oQpv1m/4gN3lq8XZ5N2D5++v16+WG6WPyyuDrje3eUVTHVguD04X69n7bffqdnq5YaHn9Ts/JPXi6dXs/by1cWFmp1/UTn0+txp8/jp6sXL1dVS+K2anX+13IDW9S+/ePXi/vxisdkMTJTu/GLxenOY+vFnD+7fvfVo1s7+ZXhUs/OPV6+eXCw+enV+Pgznl4vNfHlJX/Y98fBq+Y/F11fD+0dTUN5+uZhfzFpL5fK6wrGJ/Etqdv7N8vJs9eNXq5cP64hv4UdjuGdofKZm+HQB++6J5MctJ7j9fNZP99vzzeZgFG5tNnWpo4EPP1psflwsLntePoGka++sVy++Wr2ctaaBtB6ezTewSgEeDQBr3a0KmJ/U7PvPVz8s7r2c//erLcF8/+WCjpkmnn+6fPb8Lm3oly0h3vnm6fOhd79/8Hz14yc/LC43DzbzzaurLY1+f+vVZgWVbHN+vrh89dF8XWFo6NZTqHBb4vzLxfzs3uXFm6HE+TfLzfPVq82YYAci/nR+1ZPgkDLO9e3eav7eZAYY6o0ywzeLJ8IUlpfPbhIchFtczK+u+hlDviqpjBNeQqsz3doQVP/XmVYrrXRnJZUn12p5qzsv72wIXdiW0V3sc5Ke+rw85zYXZY1R2XeF5+yVNbEzVOhV/9cZ05pilclWmWQ7Y9sUlfy/M6412qr+rzO+NTap/q8zoTXBqf6vM7E1Kav+rzOptdqr/q8zubW2qP6vM6W13qr+r7NaMptCed1ZMwVta4pTphhlYu6sa02OymRA11lf641BmaA7G1qTHW+UCaGzYNUj6XVnU2tCUMZbZZzubG6NC9LAHDtbWpOCsjqCSed0a5xVUnV0nTOtcQkMlIm6c7Y1TisTsuLDztV6cq5l6Spb6wmuc6E1Jikp70PnYmu0BwNB2aXW6KKkgLOdy62xXvrV6tK5suv1EDqvW5O0kuYW03nTMsSCswudt621Rkn/ANJVuSLpdOd9K13MqFnf+dBKn4KGKZ2PLZ8pXpXY+VSroaN86nxupUYjaHa+tMYnJY0wtgtakBAS64JpjY9KaLYLVvqhB1xb8vDs2+SH59AGMzzH1qbhObVueMwQdv+Z0mrlIeDQRd1qFXN9ZuoUU59tS6NkWnWRGWRcqW+YRFCHFA8AOVUgUoeplUSZSa6vRWqXSWq6SPUC2NQl6hdA+y6BAIDJuktgIECIXRrmMOObwEDe6NwlMADIoUuCQAgqmi4NM9mlLg3VG9el0lYMle6ybiuCPJs2SDN4tm2WbuCZKVy5i+4yVCndAFCnL/XpLjNN4gAwdStWustMXOkGAKat9LbSXekZl8sApnIk4wEq8zIF9lWcvDExAvgKeAtQuRhEpbtS2ViRLyfJlQzpWZ594rnIs3VKd0br1hijfBCmYLQRCvVWedsZbStTKAyv6Yx2bU7KCNXw2ks7TIjKAYbW2qxcUsF3RsfWBK0SE7R0RsMyouQ1kU8xPa1yTlkLSJeoSK+azghvDcwh2xkDGzPKGWWBrPQqzCiWzsBbeebPhs7AXBlNOp7vmlAxrvRs4K7aqFKU0bozJrUwONpjcmfgrtBdLMrA0UtrcqFi0xmrW+PhWV6ZCH+Hj9V+y4C2lWI6qVw6Y11rndSjOwNrzYFlQTHNDbzVOCWszJI5ylS3Jlec4a6jJcBq0xbDlCymswwJLLH+dVbDm5zq/zrLmBjaI3+d1aHyW9h6ZxmU7IQF6M7q1DI3k1GJL+fWulLRMqmzurSUglmHzjIkUAGFO2tMm4yqhTtrrPDlvnRnGRTqo3BnjadkLdxZEyhZC3eW8YCjm0TpzppUsR0qzdvGSa1F8B2qZbmDHVWkOytDAra1XoakttVk11kZEmmsdb6z1m9bKzUzKH1rQdrCxXe9atPQ3lo4bxtMYTfMYXib9p0bpjHExLrFTO5X8QrXySzvte0c83l4L/BIJhC4zmryB9M5JnafXcC8lSAEhLdVcSOYzuudYKQ7r6tMQHEguy0JtJNTgDzLCStK5zW9UbkX7fMymYX3S3s8kxmY1Ti4zkM62m7b542uMO+17TzCEiiyggpsqatvnjdOIO0FfyN4kJXWGFnkKCzQwOp1542w+touA6unCOksNBTQnbesM7VJth+SKjN2HnGobxJLmmfS9k2SJkMffZPqeziJlSbX90gffZPoApuHJtXsZWiT5HYakAokszMC1lHrvBsWPSnqZNWr4915JwvvtqCsvNtydIcMWhcgyZ5G+EooZjvUFd4NPRiFshv8+t7X9tQ2dgGSrPUIyqFEUK6j1oWSBki+JW2vDexCkaYDadtFzRDIqthFXaV0ykRNm2XQuqhpsQxaF3WV1pEYoMI6wwTFCPuqM6yHszS5H7Yusqbwvg5bF4WFCQXWKqHCSoHSpMjSIjOsfs44+RxNBEFWFiHCPnPtj6FZZrx9iL3UDtejZJ0dfeNYVaDG2h1WCIFR6yIkObSPcnCsvn1QSYRl1RnWw/3C3w9/FJ5V2yddLVK6zLAerstbJZBoy9A86kZKr60TqM7PPqfbbaKiq52CRBKRzplj8ly7oz4LafCYBkZY0e5Szwj7bu8SVNejGEyXILoepHDYtlde7oaf5qVevBEUutQLONQkmWvr6ss88EDhDhke2PdrMF3WlRPIoHS5Z4LCJLLejjlvenQq78k9Nfb4dnmgRuE0GVrsaU/q2HY2MmblyLWKngjJTEaW0p6tkFM4Yf/JHSPkTc+OanMz6+jAkLXtMutoJXvJnLeQ1LHlRQi5QoBUD9CvDrV3uswSOvB1PguThNoqM87wyB6kcC8J9yj11Nj3Qpd7ahxqqnMVCJTYM/YVIWz3vKsiBTH2bWNaZLaMfRcxLbKrnVTppss7BgnGTtYLcKp5e67dL9TZDfuDmnnHs6Ui2aIMRf14j5K9cGy+I9/1wrHr5Ouyl51STyd+YNi00+/QATu/Q6d+ZosODfNbdGrmAR0+FAZsBOh7qC7RXWa72I+LNCRsyUi+G8Co58657hr7Du1yqBj1Hd7lAEbyFlVIEIx2IBj1/DHLHrJfD3KsC1odti7XfWQdtC7LRrJnnLluJCu2XY51ea+k0OU4WmSpX7aTu6/W/hlaLTvKbdG6pRzepV4VMiCU+tVvQClVwh46JVXCrhw7p0rYQ0NTFYa2vZDY2Qw77JwqYQ8dwUazvoSiciXsuvzlXIdtyJrrqNVFI2cZs7pK5CwyUI8OO06hrTq47Dl7nil9lOvUH+qoU39bh0z9vorSrz1SrOfUQw/0jLqvvxcOZJnKPZPugcoUh9HsefR2PHsePXxmuwAJQfc8Wj5UehbdAzsOTb+VEYuGggtMul+Lte1Kz6ZrLaXn031LysCn+6EuPZ+unVl6mWGAqsTQj3TpGXX/steuCU2BFNsNOr7+dQVmCR3Xv65YGcG6hhbr2+DYVhtHUSFzz9Y4dcXGNmjZgXs2/an1hY268Xw0tz4rNp4+dsWW1kcVtDI+d8Xp1sM+lQmmK8603itBxnfF2dY7JewpdcW51lvF9juiSfCt1wqNW7RdcaF1RYWCLqwrLrYuq1hVhMWl1iXFFjlRQ25dUGyAEzUU1CYxKpNiV7xunVMxKZNKV7xpnWU3z964eNs62fKZHLriXWuLYi+cU1e8b20WJUFBrRFQW6WAtrAraNKiSkmZQg0JDUEqyhRqYLOqMpRsuuJLa63KVa4rokdTsvVLXYE5FpUDOqmuwBqzylFZ47oS0B6onJU1oSvBtyaqgvYqdyXIhhBFmDVdgS8GVWwPMeVFPWH5JvIdKpcKyba5xApF5CpVUg9BMqqwz0dVg84AXYQ1pSvRARkUFwIKMkK1Aob61vVvYwWr4qLEVEEambsSM80ymlYCopBVqCYETKg+FapSa1JXkhEQdZiAqHwUKqEKSg8ZY3sQLbFCe4EKvMAZszKiawBEu1D1hQKmClIvb9nMK4OWg96GMwKCRuhK1gKiKRDQVBA0eIvuWpQeFXQVBA3fFRRyvKX5gChklbE0HzBWEDQAUwVBw3UFBklm0AAsAjrQcF0pFSt0T8Z2hU1URvvYg7Wv0DcZ0xW4JNpp0ACsI+hAA7COoAMNdHBCUMaBBqBQlGhXBRSSMh76Ll1hS++VQYul0aWhtAOG4DMwyigYxgCjVVfGo46Q9641aKCjsjqRHw0mGukBDhXOA4xaRhlf0IiTPwkcQEdgVhXYjrI68F7UqSaAD7DRbemXYgGNgHKMwVvbokyvuhKjjRMQXOStFxBUHG9Di0oHTASMgBFE0DiiwEOfOIBZQNCQt6VFm0evoGK0WsqKdhDQUBHsrNCn1oIkxx+FLhNRV9StptBlyLq2h2kT0i46SLgTaIsSr4dBVMTdesxRwAWVEZjDvwSZSv7CI6meDVjh2EaZTP0OiaWHqZ9NGEdLcioBjIax6hQz+Dh0sT0MPo7NYdWLop/UTk4ERE9aNbkJxmtQz2Z6w2XYMgc3cHCjYe/UDww+MHjqr/zeaFg89QNTP0ye+lkBqN87Fouqx6Q/vGcpqTD4+MBCY+Dnifp9ZBkymT6nfp9YpCos9WeWsHpaRX/4wgJXYfojaJY/OX6L4BMMi6OoRSP1B8tCWmHq72ViOeNiPIJvA+VZ0KkfuTgpk9k4UH+ILWukwFJfaiP48F7qyy20IzDtD6XliAu9r6jKo0aDWGHqi6ZFZ817T31RDgwrTH3RtfQNJ3Se+qJvOVoUmPpiaKlLYOqLURTdAtO+mFrR2iNHSH255fAPWNTxsbSF0z9g6kuaU6UK078sCxqE4XdUyMLABkcSqFFkZrpIOCIJvmUZqQnUyfLApMlFjgmNZoFAy0qClVphLFRLglTLsQDVkiDVcm7HHNVy2GY0C4XlaJEERpalAtKVBKplsZATR3JAWywXECMnn+j5NQsG1EgOQ7UsGZDfLgGdPKeW5KD/WTYgMEkAMRYO3+PBaYFm6YDEJAeIsXggnUkCHTTo4CQBTEUJ1/eHAVPRh/T9IZxeFhH6lBMLUJdlpO9T4eUsJEJ55ABTlhLYkiSAKYtJ5TMcNZKAQNCPLQlyCCScICsDm5VjIKa2JHDqMpwESYKc0HDuMHyDoxTWFDmIyAp+ymmQUDMwZyty9EDTsoKhyoEQHZETB6VGDoTkAD0pGCwnQqZAb0nBYDkTMoVmJQWD5VjIcI4DLKdEqHFAaICrpkpgEJZjCJqEiAfMzgV8opzYcDTEqSX5YcCcDXEkLDD4cTgkxJw5ZDScDnHIQntguBwQWaHtrDL4oUuB43OcS/1o93pSh+HKEREcvmgFw+WUSM6cgakf7R4cvxiV6C/0KTIRjILhckpk4fDAUj8njNC0VTBcg34PDg9M/6BUkWniFAzXsMBwtFWcqgdptrUya7yC4RpZYCB4r2C4hgVGJhFCP3BorYP+2RIAx5YzTIGpnwVGplRUMFzDAuOpL3Iwa4yrtguMKQzXeN1amWBJwXANR/JwcGDq50xe5ltWMGDDPkKmW1YwYMNOwsvuQAXqZy8hs68oGLBhNyECimwxjWE/IXoGrQL9wY5CYKtgyIY9hcBOeepnVyFwUDBkE9ijIlpEBUM2AS2pCGAKhmwC2lskkaI8+ATObBEcjRytmgB+wE7BoA1H9gJ7BYM2Afx4HxUM2gTwA04cZhsTwA+4KAc+EfzYrWgFvzYR/ICtgl2bCH7AXjnwieAHHBTM2kTwA07K0R+x9h80xCGlieCnlUU0Bp8IfsBGwblNBD9gp6zgA37AHJtz1Ap+wFHBtg0qGYHZDQGDH++RxIEZ36KgGZi4SXV8sXUx4JNC/545CBx7GIkEmP0i5ZG6gKE/YA53MST66NU5tmT/IvZ6szsXq/nG2ZmaXYj1WAhq9sOs/bZYr4pl6xdVsfCWrIotqjitCke8zqrinCrOq+KCKi6q4uBRWRVXVPFaFW9U8VYV71TxXhXP5hEbEXhbVsUXVYJWJRhVglUlOFWCV4Uz8BBVCfDErEooqkStSjSqRKtKdKpEDE2CKjGqEuGlWZVYVElalWRUSVaV5FRJXpXEVjSqkpIqCT5cVMlalWxUyVaV7FTJXpXMJjWqAi8V/l1UKVoVzqCLVaU4VTBwKUGVwhY2qVLg+3Kujh4GTqxFE6I5gRDNqOYsBhasERs0zFdjsKNhuxpeq2GwGrlaw1o1/FTLrhhOqmGfGilZwzg1coKGRWokYQ1z1HBEjXig4YUaBqjhehphVsPvNExOI7tq2JsWwyUYmUYi0LAwDd/SiJ8ajqUhOY20qSE2DVfSkJmuxjuUgKo0ihkNz9GBEqKVQ6jUrPEa9YtGptSs8BppUrO0a7QtmjVdo2nR7Eg0K7pGuaJZyjWipGYN12wQNIu3RoWikSE1Sze2xwhY/FACAwKNBKkRqTVLtUaA1KzRGtFRI3tqBEaNFKhlaRehgWVZIyxpFmSN3KhZinWhBGKjZiXWCIyaMUdrgdDAD0YOjHm1amLMUUSw3vPDjp8xF0ERLYMRCbHaQzDmIh8axrwaBjDmYi+FaoBVlR8xpKAOxpx9P8sjP9TBmItAWC2yGHMRBw1jzmadtYsfSjDmbMJZj1BkUAdjLsIfG2oWFn4owZiLwZbIeWyPWRb4oQRjLspAw5iLfZxYvBjGXAzc2MbCqPmhBGMuFnayC62WYYy5iG9iISdym2wKRWCrJiSMuYhrsjsTOU32XSbl73766Sf1e1lzYtx+ozVn9Vh4i4dJb7GJFff6koyDpXNNOdl+YM9VAnDkDrHvCbHzETF4Zhw6iXw+X3+/WI+cTmrC6JN9wtaR4qvF682ty2eYYWOAClhf6kbTB/L+YvnsEgvlHh59n9d3Vhi7R7Ewnr9eHpqO39psbpGO/fbZ8ofl1XJ1eTVrg6FG3ow+eHf+ZDE4wlCfwLUGTw0C3zs/v1qIhwp8tk/cou0E7+XT7+8uLp/hX6MbjZGzjMFQVNqC+fl+sc3FYG++zTJUj2Hyo79MC2U4f0UL//Mv08LtAP3CMfxovh75F300Xw9EISbYeGwxTS/uP6hz4uP1/MfqlVHhey83Ow+QCvROIBXo/UDuvdx8XO3tq48ZRvBMIplG915uepZAI+693NwRn6w+651l7xhwYEFPBkk8W25wLBvgr1arCzGgJ6F6sNxeXW5Wr9ZXvfPCrU2Pzh7HvLXZMImFSb2FF9hfyAyYKzQfp4hGHNiAcImo0CeXZ5+s16vepYuZLaBkp6o7ry6f9myBl4AjLgbYDyFv8YjpM9N+wD4z8x5wNN53F88Wl2djtxuGoKaOOCwf2iUOdQ++eHxiYAvbjPRfT4pqdv4pvhCLqz0e3qc+eDl/ii+A1L11kBu1Yesd16eB4zbfFJtt1iF5yLpXteTbb/cocedz9OnyCoIc40MS3+vRiZpmD/mGimvvDFmH1D7jHjbk+nx5uXzx6sV/LtarnasHLyY+i8Liq9fL/fXifLH+t7u73DV91HE1YdxMMB2n7tpZUz9enH86a4Nm0LYp38zaNE15OOMcYJTlUZ9wfz4mvvvzCW1R+TZpV7MkHfpq3p+fTZpO392fnx36ft6fn13j/nl/fgaxP9x1T5/yaJLC2tg7K1Hh8un3vavS/fnL6mT5sGca24RHM44nZ+cPnq4Xi8s786fCfUAPtjbqfkDmwYhsSRqPx1BqNH/IAribPaQMBFTJquZZv4CLzmLjxZGGROEZJfe4CIthc1pR+2bwdyLnpwBG2rFZL19+vHi6fDG/uNo6FQlL7qUduxUZRq2TDHvNk7Rx+yASSRw1cAtvl6oqpCCAVJ/gbaktC5QytTllm5v2DN+nOZX5rD67vFysv6R55GSqyWevZu23OBKdnPBj7Ilh03/iT/IAxxN/YiKQ1WTwfZYdJEnjdKun2bfvsu7fZX3i9Yk/8ZJT6v6dfr6DuSzmZ4s1q7X4TEm3baE7y82dgWhCTzTiicUobt8ILUmvPZ1fSGFG//+ulpckDoLA7fnLMfjV8sVWnEw5F5O9cIzPXsyfLfjQlsHfnl+eXSy+eb68+n6x/nJ++az3Yq7pH61e92l19GqqYDLy4vyP5epieTmk9v6LNevt5frpxT6371/hfgrSowXwISL3J69fPhyLPUPio3Hio+tyDomTnGT8fP764+Uz8bKHCO+tN89Xt+cvFut5z30Ot2v352fvzf+OEbtxx3Z/fnbTfo2dFD004liAA/Pp3476b29lPNyF7bFrOuMaZn301T9Nja7/UvHJZvG4N7i32yYH74LL0fscYO348sfG1X/BWY3OYSau/baEJuZiYjAlxcCKO/X0tzo3SVutdcqGMwzJMDj+mxKb5K3WRvscZLSGsACz1uYG7VUp2elYILG9kADONNm6UlIuLhzEB5DivC0pTl7O2lPToF6yPuai0fJIlIDBIx93fNMUai7J2Zhcn2HAGsY/RBXAY9g3yeOQbDX6n9l+RIGYG5ulr7NDoSMZhm+56JsQeWtitAWn5iH4AINx+HbAskYaOE2NtyE4U3LyJnvGpQ8cQFyC60Z5eC3bARdCKqGYFHWUxXn3eeIQ2FIa7U0MWusSErj3xXl7Ha0MryVqAZ2rQzQ+QWKM4C7UQKM1HUwfp+izFUrbRR649vUQ10A3zg7kGGIOMIohzgGRM3rKttYFI27iNUxC40xfKmXtIq2p0Qx0U+oAae19qkqMSVyEa7GZREq4KcfNVW9jKVxXex9b4ZTxPAyu0KcSZ8FG22TnM1pHh1Jv9ktjLdB5x1gL22g2vzjWwt6e/t3L3xSd4bEZay0PgjfU1+8aoeH9OfT3Ys8f4s8v25QhmtB+D8u2+SaJgkn9fiUKItvIPByHttmuf4fTb7wdoR0V7refh2LYb5TB7nz9xW36osZOemu/fWo+vqnXgsTY6AnynKA/j18/DuFp5KRnfq7ZS1ThftbObp18uVrdv1htTlbnJ/81e/1fIvqPOr0qKH6jIPfL1elfPF1csKHkaGh2foM6vNeF3xSequ+C16KHHbX5NcBfWoEO/k+W6P45RJudP5SIPULYD2UlgaQf1hwElTm/s1xf9arQu/PhiUBedtgGv1h8vLx6eTEfRcBhX7bd2dFncmawi0jz+ers7vxJD79Fp/9uw/Rmf5gk4snVyX/9n5P/daIbc/K/P6xha4hTInvkrSqiPzB5y8nHaOC2w1YHcIiz9MeP2g3nFO82av/YH7UPa5D+JkO0d7QyOnfRMv8/udyslyippMFXr1782E9rHocpzvMQg0ue+xefz1+jPL5W/OwVy+yy6nnLF6v1i0FrCeH2hyM1dtX5g1qdsJxJNEDOGQhWdXim+e4iVg1AtdNsDgGppLaPlpcoJj9Zr+8RQAvUgO/9sFifX6x+5PCGM9L1GgY3td+JGPaI/Y4sNkR5EiuefXHg/Z+Sv2Sy/Pxi/+VqhS6fzNctdWKcso0P9vzxC87cPp5v5hQQBjVrZ3xg9Ww9f8E6v83x+OXFavP4sSyIozUfKvoNKz6d//YVH8Fu/wD9i5er5SWRJfsl7W3DpKuZ1alvSlCnvsnym+Q3ym+QXy+/Tn6t/Br51UGdOinrpKyTsk7KOinrpKyTsk7KOinrpKyVsihX1KmVslbKWilrpayVslbKWilrpSzaCHVqpKyRskbKGilrpKyRskbKGilrpKyWslrKaimrpayWslrKaimrpayWsrrRvt+e139JYZTQb6zrP6ckrxSTL8jH5LtShdQmFQsOgo5gJkgKvoK6tEIaJG2TZkqLpfHSD9Il0jvSUdJn0n3Sk9Kp0r/S1dLrMgAyFjIsMkIyWDJuMoQymp5O8rTZ03xPT3g6xdM/nq7y9JoQDLQjYfYe/QylvZy1sdKbUU55ZVTkF386hV2kwe8uaomIgqmOGDmJu4nY72C3hYVVdfoQ5wWJs4LNEyZUGE+JAwcmiphmYu6Kvw8m3tWqW0yMMZ7CsQyrKYnIJW5leGeJOSkmn5g2S9wk/LLw3KIO8uEBJs+Y5orPTI2pg/mXOB/g14PTgcHEHo81bKtSxFOP7vMK+yivxNKu7wKtjHTf+2DtW+60jXZ4sPBsV6lhJZqsL5+gi/qZcfxBlGRBzTDB0bKZ/uTh8+Wz57+iHLqtnym2IxvdZJusj9Zq9KqB+eFs0j7miHlcYoYY71wxOVprc2amTYvYxqVSnHGp2BIw1tsv4hrs8or3rriIbZhrYkHZ50zO1hsL2UcfjbXaZR9sUqEJLoUYnA1Zh1RUbLLDXsxZ9orJqtQk651JOmUbojYqN15HG7IzzqOkUrnJsThvc0g+e7ExbLLGrCymEKolY+Ott0nn7ErGENA20YWco9aueAwe4QHRlxiLt57JcpAAC3C5WMshJlaGtim6FO2Djckzh1xTTHTZ+mJ9xMrONS7GkFPMtiTM83wTTXK6JJNMwCrRN1brokNIGB7iWNvEaEMMRK3Fbta4JrkQ0OWm6Got2oaQoi8+JzB1jU/ZeQ02CWNB12Tvncu6YDaJJ04TjHeWsSkFW0bfaOtSCMbbWDByhTMZn0JOPkW6kKErmCcmEy0Gfwd4+Mankl0uTuwcDz7pmqyzSb44UY4a3xinyWu8zlhBuiZn7XXSztuAveQ1eOfgs6V/U8YyFN7pfCnRaqiCBOcxYPU2RTEl9U3MdFfw3kh3+UYXD/mV7BwmiKGxOQoqLoi96kFDHP0AyfmC3eY1Y3TQv9ZFnTDkdP6GInukcdBdvrEJs0wmSyVIn2JE6Ww5W6B3rInYgGqvxcOCRcg6h02qi/BV2zgTbMg+Wo2dtrGNtimhho5FTEttY5MPBUNQ6FYZ09jCeYvJ2gTMwJvAAYqxpdgszF830XhXOBVNGIDnJpYcUuAoo3DQw7GK89pGa2LyzqjYxOxMKs57OjiToFNxoptnGqvQ2FASrkfkiuGQUbjGm2BjKMU7jT31AW854D4/z7D0QZ4a5PbRu7Fg0xjI0Pvovcd+HjbOumxm7be2saWAMb0bcNtuiskmh+xszg5eaKJNwZlkbcQs+aCEa0J20ZvgQ8k+5MMivrEwqJIiFuiuwE+ddlknzi2yiSo0sDBndUxYCsNPI1PDWJdT9rGo1BTtbS5Wx2Cx8M5NykVbh8m2s5jlNz6G5EJikmDlXZpiQ7KeSaGzmDI32SdfkrEGhgodedrF7E79XI4p+Byds+Ie5JpQovPJFGsTgshBgm9sMdEn54LHKh926nOMvlA3Zs7IS9Z5Y73LEYNj33gmg4/Z46tgAlwseJ991uLlFZiYxnuTrRb7aC9srZgQTciVI6dkChXlIEV8o6MrxZYI8xTWyDdjigVjeuHZ8K1Mh1pcT0xogs9ZR5Y2h4wSGh2Kd5mTShwPQCyHCPNLKWJl7uHq0RvvS3LYmft9PILw7JJLSmIAf/BNDxfyJjgOK7FhD42JIVkWXBMql6dDQ2bqiZG6P8Q85xyw7M8xiA9pA4WzXIcYaWyAGxbvM0dvCIChSaZEl+lW6IlqocGoY0hiFR4b1vuYiovOi2x30BSPawC0VTJ2+NcM00EX22Sstclm7ZE9rymyRxwHHRYap6OFrTibEAJcE0xMOkUOdetQ2+AK54cBzBiVkAsHohoCkmXRhZD5F/HHcA1HyVoXE60Y9TtYsAtRs0LhOWEb550pSRfS8ENooo2e1bjEWKdOjCY4raPNCY+VJlEeX4US8FphWHXxzkYEM4NEFBJyh9be4GuTYLgB9xRjS9BRxcYFRofzUh21P+QTbEdKNkWH5AqeCAes5VfwqwMGuc/ivlOzHwjgf8hISY0/qdnjN2j9RMs7a6/Vvz7exoZ/XOPB28BJ7I7g1Ozx5Wr94j/mF4S91gK+WF7OLz5aXvamH3LpwuPLB8tnL+aia3m82Kqv2Iyq2eP1/MetSuvUqNnj1Xr5jK98I8G+RYs0dsEaq3C0Gonexx0bvjjjHdt7P7S6eGcl1u1X6x9uPOsLFmOUnR5LdFSPP7vcfPvmu8doHb99rd58N1Zp3V+vuH5kubpEp/Vv81dXV8v55cmPy83zk6vny3NihJ+8WMzFb2Sk3YqioLpZv4Vpw1tNk35ev4UC7Wb9ltiovE2/RQbcCE9DIzab/FediuaqiLqLn8RP5EdeeJ4cP5Yfo0696Lj4yfwkflCOOUoNWi6BUZM5vjDouxyPlh+jTh06rsJPFm0XP3zH8p1B4yUw37F8Z9B9SVELbNSpRetVRPXFT+In8iPqMM+T40fyGXVqFPoufjI/iZ/IT+DH8+P4sfwM9idbDZdtUgohJGSNaHNZnJqIjgsVFxouFFzot1Bvod1COkW3pQycFLUWWq2tZmur3dpquLZaLpRcCMTKot+C7bJb32q4tlquraZrq+1C2YWoqhwMHzUXWi7Eb8RR5HI0XCi40G8pVhLDWoFmC8UWCwlqLbRarDLotFTgf43hT+uqonmbhgta26q8v8X0JaaYYoy4J9qi6bBr09jsJBudc9EEYos32jkTEEyTsa4QWbpBqi3ZmlQ0ayID4JK30aQYrcN5Ujc2luKTFkHHZQYj4DnI/sTmWLUhCNVZZCAbE25wTXTFiliQnEOl1aSsvbPFW5dwpfVN8CXbgjGXxh82NSazTSw5xFwdPbELy7Ho4KpPX2xMMA4JtsQk4SIaFzVbPGu8w8/ZNdq44otz0SXi1afGIGQjwwdLvOkmORucccFGW6JyscmoSjDrSYT992hyrM3e4Q5qjfK5ccWmUGzU2DSpgHqSrZvWmg2MirkpJiSUKsUGQg34xiATsjOv/q2Mkda5aMRonAFDk4PVsmnAKdc40+jsfZSLZ9id+tAYj0bCghpSsmGj6IrFC1ZkzdREW7Q3PoYoGkzXmJx18dbYjI8mxnCx6Bijx5CQQDJNjrFoV3IWRWiKTcnRZFM8WxNCszTRa201Iq7EFkHvkJMoI3BJNhkrNY3AmSz+0abYJpvgsk3BRzbnJTQ+aOOQFpyoPVNjTfGJ5njwKFzTo2PM0Wf8kQ2mb8jExSJvESCgND4Xg1JOlGUkJKdNdsH60OfI0Sbriq8+saVAhin4EIuoJkhgs21McaKRLQ3EFU3UjJd8sxTH5t7qyD6DHCX5JNt8PFolISddcJWXtpGjJJdzsbLfqQmlxJx8IYRWP9bReefRyuxn6BP2P7GthB69Fo19RA+ast/W/c446K6DDj3o8r1BORi2g4E9GPp94jggn0MC2yfBAyI9IOMDQj+YCgeT5WA6TSccu5n9KTmdtOxW9qb1dN6jr51yhrjPOw64CxH69/jPAYfa42H7TI4QcFM2CEVNGKXa46Sy8Zry2j1uzBI24deslhOOLuvtlOfr6bIQWIMOl47p8sJ6f7gETZcpLlXYX8zkyOyaJW+yWmZuUOB6L9SyzpkUtM4JUUNrlDOxlGQT+nI5iLLoV4yP8H7uETNN1sm4KOq7pOVusdjEwOqTHbomkxan9Ax73GTZZEZtHWm50eh3PJom5+ziVKQMF7GnjuygfSYtN86mBOuNJnnrF6csncFpVvocNfrLX5D2wRwATXeuv2RncLij3WOw0w3te98+3Xr42YMHtz7/5OZd1L+sF+ez1pT3XvWTC0yp38kE4av78x/EX4vs1xkhvEfnD7ZXYqfpkOo5GMopcdqBYbiZtTaFRs7RggvGJmrGqcE2qHUy6ll095iKc1+hLb7RKI2LiZkzHT7em5DbYoPDKcCUwqEI3+9fFd+fi3N8oan6oe1fJY2HQHL4ARBCgXqGUqWeqvsUc71rcGxASsM+u8Qlmf3hg+fzs9WPI2eZ26s1vmfzsyVXfJF5e8S5HanepHGDjZnshW/0NMb2g5c15IBsOff97bYhB3YBAXrXr7ury2cLzATFjWB362UNErC8/IXX2u1bVc83CzE5uY6QQLxvpTRQWnqTneq7Nn/ojbEH4rb1vOy9OHv3zZGy4O2mMNgU3eWmv2oSVW2ToKN7oq/aGgbC48UrRBdng/9pe88fd/sRfuO7gwAcR4+uHdkxCAfuvH/r21ev8/URjy6TcJdh/5S18eJi9BaPLhNCU3QyyZqCpAClT25IjTo1LuYYk/cpiyfY4LOF59Sf5NFlG47js0OQ4rCiumgN3kzfz1pfFfkxcljkMMQdXLj6+1YH8M/26LJsqjnPso5jxZ3L1ftw6JJNR/LZJpM14bqmDl3XkMrEocungKtYNtn7A4cuse7gHwGbZBV4Z3euGORKT85XNNZFv7c71yGqH4Qzl3G54ZQzcODJLuAXO3OxuhyduT4kZy6c/rcmtYfOXPL66Mx1/Q7hD3bmumb6jWVxRIq/qDPXm96ZK5ezRX6bM9ebfwZnrjdMyaMz19GZi4BYdX9K5Lvt1uvdgtgdnblwEjk6c/0RXpJ/jDPXNevfgU39P5Ez1353vH87mA/FmUsWxHfWYCEH/UZjl6Mz19GZi3jtmP988M5cEp0dAwnLnV0SHlpXhydNSHRCU3PuyktilBOinPDNXNslRrz1HnCJek1o9v5ScEVsc0K6c6GXvCT6tUTPJkY24au5wEvCWcsrAkPzGicAKpCA3ETrxnSj5sPBrMYzl7wSJV0wkrsiFKa9EjS94icBuinNV3lFkHjBQf7jFGYjfdOCtE4wxuK1b1wtSOB2o7LhzZCf+PG1snrTGbYONEpaQjB5+kl6VHD6YA4FP0SvMKE/TrZF55rqWDUhJpNj5twKk2McRlKyEVNpiRvfFFdwWnLRRwn83rgUTfFZ64gN0XWOIxgCYXZSHPHoD+uc+mRhusLRRMbKGfepAyz2vlBwFZPgYiFhc80XnNMYxXuXxdtRY1vtCqYnHsv80kQajrOE1VFC7e/lMHLujlsDNjZclLD/UQwNbEnBuowJkFgW6MxBXggGo6rDWq75xkFb8bLKLjoXC3dbHPRnOcyRcFLA8icyoa5ByzjrOO7RxGQ7bJrZH8TrEnyyxuKSZhxXLDS4BuAwYrxPYtrLeX4uvliXs1xycFitxdSfEGkuSXj+RsK7YQplk9xesNeBRjcZA7dsOOYUdrTnhQSm+9VaPA6yLrl6J+om6Yx9vEQ0u7bEtG0H1HJAC2Yf84OBNXoPjWuIdI98SqOLrVRLpH2Hwb+POeFNE703iVqmM8Mc+DPuk/HBRw8QO2iL0XsT9qCPr5s92k0ny+H8mhJhaZwuHu8K3Bdg7I030dpQMlTGLR/7rS37hL7HnK6ZGXtObNfUsYfndU07INJ9fnYwkMEaHDFziDg7qdyUbFzGG1Z78Wo+bOvPN36Ppx0wwUPaSAE/tVCK9dze8rP9dd20n44atDHhJQdY7LFiqGu6ohx0+BTvX+SJ9y0Mwlmtk4e7c/OIgUEETrlwGJPpHqzB2ilgS4sA1JTsC1atMXjMKw0LEP6noeQiV7scuuZYp1Ow2Lg6uQLlsNo9vzfDjHA4+HHRBgLVPh4H39Dw1OTEEznLBSr0jdd4KeFYBKZJZ2d91DoHuTqFKYFbsU2YJF+TA8s+pxNXoegoF8Hsf1RDe5izOaxJEPEY4+CSSwFfI1Dfq+XgG9e0NqZsPW6MBQfFeko6dnfShzlSFvNm7YyVm2SuQSybFJ3zAWPxw7YdjOQ1CeLRlaKPJct9MswrU7xLvl4AgzyDHbUrWHxLb+xValMJOhsXcXYlQ4jJmcKCHXH1Oug/g9BkS7A6azHmNnuuXuaaWjHqcU685pHXiWHhi2EVK9zRdF2RvaYdUMsBLRygfs3A7uFxDZHukY/Gq0vjiBlCnU1EdsX12tHJYL43NQ6cRg8w3/+m2UfroCVmf8oedPEBVR90zzU5prSgG0wkvUOuDHK3EjzWeWtDdF6u/tlvqz6k8j3+dM3E2PMUPKjlANGDhJ/vn4PZplk9XbSRwAhyB1WjXdLMHGu59us6PH6++Xss7YALHvQYtqjGaXzcU5aLqX62x66Z9tOBM3vwARYHZH7NurLHFPcQZyGb2oJeF0XuwOYTW+zR+hT3jD7/HC/Go6oCbQqKmN9RVfHejWp/d3fI19/93dwhMRF6qzuk/93cIXFxHDk//qzbozg87pwcefqLuzfidiARiQjYdaMrI16MaHfxKvitbos1KNc2MJd4w/dejPjo92G6osTpwo8vSWIFa/AufB0lMZM9SyLhOLhy0zfF1QR5rl/kDkvfFHL+UvdIP3aPHDxyUgraEKX+bQmE2xC9qfjmpRQwMxcV7JCQg5MVfZtA6KLJRwshE8a1lCxBu0ZFSg0Us/tGSSIAjhK40nNXJBKhnz3akIMEuStylDDGImot0cUm7+XizHEB0edNEsZt5Rt7df5yeK8V1hliGhgJ59X3qdRj1OzSzLhKXiIejIatfz1tyuF39zP8Qvig5Yd9I5ryEc64MU56bwoejM90BLneYDLG4Roq2KeTA0o6oLV9ajyg1wOKPqD5/WnCvQYIsUPbbcyyfdsmiKfUqMNjsH1Ysb4IjrBWtuF9grcxyyXiQ4LJBOjZ9WBM0QWJOdLnIMBHCnI3a/2G59KJeglrTbCFb3I7rvhj4h4SE/doDAkYzxaX0fDUHNzJiXMul9ziwpmJ4OKD3LRbE3RxCU0jh0bVyTNlZyNbmJrAzQk++8BJlOQIxHAJnpB9NYEtnksucu2u5ED7bKz23DdcE3DqRXnRw5bYWEUTDbBmMNYU43Bl6xNwUfMxc8Wy5Mg548jmcUWuCck7bb2VC4upJBNULWZCxvU5CBfF9gsFouRgJ2+DJihYTdDGZC9BaWpCSglXZF3skBCDZtii4bpg3GQJD2NikPuLJcGl6DUBdfpacKPOzlqiNdUcOrnocdnpGydRg5ILRXMXcVPemvDBnFtNNzC/yZltoHgCH8Ev9zY27136fVdvNm/fe9VbH6mfD6j7IXmzGU+Eu6RxvXfRVE+zG73ZcOjQngusdSCCEZLz4LL2Nm82xyTJoTgdPdvbozcb5oJEyf8QvNl6Q92/ijdbcAb5whOZzWVz9GZbXkJOuGscrydbrc8WwwUpwp3khi5Own1E1ok7R15CWxPSxJQUUiYeRvXkNY02XhOg0gfsBqoj73WZJ75sRpsgIQm5+CkE0SpstndzWR0bLpWPyaWCILe7uww319g4HYgEg8DESE6vJzP1TrXd50y2+xeBDW9pcCR4tQ3W2xhSAZWvJ7hmXRoEKud0jHrwcBP3rO+r7++H4sEWubeNkC0pWu3kTqrei0xc2ND47w3t8LpaZEfUyCaEbIJO0tDBj28j944Rjg6xLeniAyLmnhPbIYlMnNi0jsZxOMOVZxqTkcm9ZBLJiH8eCxQI7J3d2A5vHvt9byU7RPWPcGPjZIRIzcb5UHw29H9/RxmBt4p3RIhG6e41ru5fv5abSm4oxiUm15X66vliM5+1Do5w//myf7p19XLxdPPlfLNcVZqv923PxIP66Pj2ITm+ubc7vsnrv6Pj29Wr9TlNf4c9xaf2Djmvi2nAYtPHNHje38/x+PHraqk9BDj4dHxBx0GukT33z4Uv/FmLbjK8/XqOXY7qRgM8XMhF3Hp4rNz7s3+b2zteGvS3vpEriNjzd7yQ68NytaObf/31Wx/SKL3fC7imB9725PWJbuwHdnkah0m/ZfRgSHWK1cn2512e9s7OXGhDafOvvZxruPbk3Z24brrliii7u4ifpgmZ25y5fETr6AlE5hqTQrQ2ZcJf2kyab5KEHERxrk3pg5MRgckHBERtQyJfwuxKB67P0NpZQ1ppAnr+7ILO2uVAGnp27hRjQ4gKX/saO83YpK0vhjicaOQln7dFOyJsRqslzBtB1ido13zJx0BMJ12IpFbT2NTJ/sJEnSSMm8YoD+MoXbAXrfHZrLW2OGTgQIg7qrUum2SSztziUAhMqrUNXhNX1Ir2ui9KPEh02zVweM2XQ9ZRBx2sx22Cso5jHI9ZMRE0reib3/ONaXuRlojaxTVoVzdJJuyIRrIF8wrwN9zZTuQuDLW5v8ZymhOcbLsk2RMFpTByGLFXPSTZS2SnkTKR6up2mdBdusGAtBDb3QZi5o2UnegLRhGV2NgNetDYB+FKyWfCr4+Unf3N0Kbo5Hzd0m9Dd7kausu67Otufz9cwK8M3bUmNthOALwahuNdtX1GwoAsXg/BqxijXfCqt4Tuott3MZT+PGXngTB7c8ub/iZvyVEjlwnT3LV+AO+sLvtr/n7n0F02GQLVZm+5GyMfxOn6NcHM+tj2J/96wqpwcnICfdzcK8PLD6RDtPU1tDYW3u+nQz5fzC9PXp/868nJCfG9j90xvzx5c/KvJwReN3lxquM/e5c82JydfLz4oRKJaYJNxx6pPQKdnNgm5780jUwiIirCIl73IyLTvZebOwTRhC/ee7lBvpm1/fW0d5ab7Y3VofHPIBIybBNjn3h/vl6wgvxPH141up+mOEhV7+HK2QdP5xeL8/nTzWo9up73TUVfdgNvJhf09m+uvanWEWBKbqqtWpiXsza4ahsZmqiDKYjvWPRZl51FHuc+AieuUotTV+TCTML1L04d1/8ZH7WLi1MXVGy0TinnxSkh6AmwmjVluEiASOrOL04tH4hRO2ITW+wAubQjhMUpIZ6x7i4lLU4tH8CdMubFqeU2BG1s9HFxarWKDXbt7B9wVyTisfZczoClH1KetXJTQ24CjgxpcYr9W+OjtxGAexUMlzcB8OkYvNN8gJZmQxz1xSn3X3K3oc5mcUo1JRhiJ3NPg005EeI4cfBQXObZNdxS4AMhknFZyDpEoiqz6bHaFp5Nw81vUTYv3F2UY6n7nSBunts9jUt6CANtLE69/X7AxGDrFofA0NHWvYPTwRjdh432Be8TyU90Zhw6ZfujuUZLLm7zxdZqLZeRmYpNZoRT8DU2dESVL3Giib3P/qq2iruUagsdDp+ZdQUa0CVktn2ZHkmp8J3CRXOGMxvpy9xw2QAFjFFaaW6B5K5KOoArMpILlg7zmm9bLIwWpw4zTZ8Z7cWp43K0mDARWpw6rziLc/UN7aJ5xi5OneZqi+ICVGDlbgxbcoKMJP431y5Cep5419E6IT0uzHDsxygjt3dEbkYUAnNN4LjaL05NxiMwaJcgMLk/0ycfyCY3bXDfGIBXoUk66J7auN3JmAC52mk/THpo3HXjLh13NTvtYQjGQ2NHQzYdyt0Qh2Y39GOSkOcRqYxJaExaY5KzI1Ick2hPuWlE0CM6H5P/dFrspst4Go2n127SjWbieILuzdzRnIbajuzsyM5OWaF+D3Y2meCV2o7s7MjOILjfgZ1NhBeo7SidjdeE3UIxXj7G43CUzt4unU3Y2USSgdqO0tlROpNNznh67SbdL5bOJuxsskur1Ca37x03m8fN5q/dbN7Mzlwz2rJBbcfN5nGziUbjd2FnE3UU1HbUnR11Z2NFz3jHNNEMTXRn78jOTDPSTVVqO+rOjrqznRb73XRnN0tnE3Y20btDbcejgONRwDsdBfwadjZRwldqOx4FDKc0x6OAtxwF/Bp2NjlgrCvp8WQT3eb0COd4ssnJ5m9mZ5PTxl5uO55s1sPo48kmyv7RyeZvZmcTS4pKbUdDjaOhxtZQ4/2yM9OMzCqgtqOhxs7m5mio8Z7Z2cRkDGo72p1hYDayxxmb6YgV0N/d7ux3ZGcT+zGozTZHuzMhuH9au7P3K51N2NnENrZS29GM9sjOtiank5PNyVHARHc22WxOpLMJO5sYykJthME/mtH2DqO9ZfXf3oz2j2JnEyeASm1HrwC8lY/SmdjWv192NmFlUNuE2U3Y4B6DPHoFQJV/Wa+AP4WdTbydoLaJIDcR8SaquckZxOSwdXIMOzGfmzRwT214dHIaH5P9AV4BEzPa31E6m7CzyRYUajuys8Ef7m/n5DSZ7RM+MOEQE94x4SoTfjPhRBMeNZHOJuxsol6D2iYfmXx+UvEEpQmyk2ZMGnhkZ6f4L/6RPpt/lO5sIn/dzM4mRwdQ21E6691+/xY+m5PZPuEDEw4x4R0TrjLhNxNO9GvY2cR/HWqbfH5S8QSlCbKTZkwaeGRnR3Y2clmfOLNXue242TxuNrchL94vO5sE6qgr6TGixl8tosaHLJ1N2NkkagfUNjmQmOj2JtvkyQZ60tqJd8TEb2JsgTx2tDg6Of0xTk4TUWgiJE3Ep9+RnU3cEKC2CUpH6exDDRA0meCTqT9hChN2MWEkv+Zk8x03mzezs0l4orqSjpywJphP2jRp7ZGdVbORqXX57xogaLIXm+zSJhxiwjv+fHY2cR+tK+kx3tmg3x1O5cfGR39evLPJBJ9M/QlT+JDZ2cQ1HmqbYD5p06S1R3Z2ZGcjl4LJZnMinU3YmRuH/agr6TF84wcTvnEywSdTf8IU/jLsbBLSCGqbtGnS2iM7O7Kz38jOJuHaoLaJ+Hk8CjgeBdx0FDDRnU02mxPpbMLOJpF1obYjO+ujLv9J0Wgns/2vudmcSGcTdjYJs1t528Rd9Wh39k9ndzYx1HhH3dk7srNJCHGo7SidHaWz9ymdTdjZ5HqEytuO7OzIzra3Dvxm6WzCznxjtcsFD8rQ87bRbQrHk80/Irj231o6O2BnbESP9538/hcETJRLfwOF2WSH+XYexu5zdz/M0SSDW3V+9zDaf5kzzHfcVl7DuLD+Pl7PNLobaGzhNJ5m4+k3MZWa7JwmKqJ/anbFUcDucqtxr457e7xojAXiqYPGMcDPP5Hu61omJddQHa+O+01Xxx15U90OooCQy1i21/AdOZLcyDfc0jjezNx8N9w/jzb+Bo4kIvkxEP9knT9KSb9sU1d50ZH/HPnPW04Dx/znu59++kn9z+zx5s3LxeX8xWLWzr76ZvHkweX85dXz1YYbl7++XP73q8VnH9dLjj9abq7625qf/L/F0w3pM7Lde7lZri5n7ezJxRcf3ybl35eXZ/WW5O3n2r2q7s9/WHy1eH1jPUXNzu8uLxe3VxfDjcuADzZvLhb104DfLM82zyt4Z3lx0WfmRmnAPjPgQzNrT3VTtHHecbWtyToFNTt/ZGatbkooIejinNVBa9IfcptzE1PUOtloTTa+GPLbWWsao2PyNkUuVI5Ovk/TW92Y+tEeKraYmKzVOpmos3y3f5Ul4wCVEjT/fIrZJatm5x+t1meL9YPlPxa10z+75LJsUHjwfH62+nHUMbdX68vF+sv52fLVMEIHQ/JFHeLNcnOxYIjo+1uXz+hM+gdwV5e8vFg+u5y1FlyAR/UB3llx7bbn7d35k8VFTwt3V5fPFlebWWujmp1/Pl8/W17SLdKnjNjVASXcXV5BBQMN9uB8vZ613+4Rzd35ZvH6bZTZt1IaKC1tZ58urzarZ+v5i5PV+cnz549frM4WF48fv3785pd0w9BHIwrb9YIQWO3HR/U/70y50n/LF8vNHblavA4BBHVvvXy2vOyHpNHBx+Rztt6bFIrZv+P8p5++I+n8wdP18iXTVLrg0+Wz5xfLZ883t1eXl4unm1l7Pr+4WsjseL04qx+XpJ/+Py7jT88=').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222630876', 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_1779222630876();\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
}
