{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "c7d1d455",
   "metadata": {},
   "source": [
    "# rf103_interprfuncs\n",
    "Basic functionality: interpreted functions and pdfs\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:28 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "d1e10ae9",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:34.560571Z",
     "iopub.status.busy": "2026-05-19T20:28:34.560454Z",
     "iopub.status.idle": "2026-05-19T20:28:35.496305Z",
     "shell.execute_reply": "2026-05-19T20:28:35.495646Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "98d9bb17",
   "metadata": {},
   "source": [
    "Generic interpreted pdf\n",
    "------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "d8583891",
   "metadata": {},
   "source": [
    "Declare observable x"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "2d9af52c",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:35.498342Z",
     "iopub.status.busy": "2026-05-19T20:28:35.498208Z",
     "iopub.status.idle": "2026-05-19T20:28:35.660852Z",
     "shell.execute_reply": "2026-05-19T20:28:35.660131Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", -20, 20)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "920c2c9e",
   "metadata": {},
   "source": [
    "Construct generic pdf from interpreted expression\n",
    "------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0fe5dbeb",
   "metadata": {},
   "source": [
    "ROOT.To construct a proper pdf, the formula expression is explicitly normalized internally by dividing\n",
    "it by a numeric integral of the expression over x in the range [-20,20]"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "630c73c9",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:35.662895Z",
     "iopub.status.busy": "2026-05-19T20:28:35.662768Z",
     "iopub.status.idle": "2026-05-19T20:28:35.853877Z",
     "shell.execute_reply": "2026-05-19T20:28:35.853332Z"
    }
   },
   "outputs": [],
   "source": [
    "alpha = ROOT.RooRealVar(\"alpha\", \"alpha\", 5, 0.1, 10)\n",
    "genpdf = ROOT.RooGenericPdf(\"genpdf\", \"genpdf\", \"(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))\", [x, alpha])"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "565be49e",
   "metadata": {},
   "source": [
    "Sample, fit and plot generic pdf\n",
    "---------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "556bc2ae",
   "metadata": {},
   "source": [
    "Generate a toy dataset from the interpreted pdf"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "2e2bf05a",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:35.856021Z",
     "iopub.status.busy": "2026-05-19T20:28:35.855900Z",
     "iopub.status.idle": "2026-05-19T20:28:36.010648Z",
     "shell.execute_reply": "2026-05-19T20:28:36.010014Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n"
     ]
    }
   ],
   "source": [
    "data = genpdf.generate({x}, 10000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "7df8ada2",
   "metadata": {},
   "source": [
    "Fit the interpreted pdf to the generated data"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "7d449c1d",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.012699Z",
     "iopub.status.busy": "2026-05-19T20:28:36.012556Z",
     "iopub.status.idle": "2026-05-19T20:28:36.213697Z",
     "shell.execute_reply": "2026-05-19T20:28:36.213078Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(genpdf_over_genpdf_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- using generic CPU library compiled with no vectorizations\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 877.084 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_genpdf_over_genpdf_Int[x]_genpdfData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n",
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooFitResult object at 0x(nil)>"
      ]
     },
     "execution_count": 5,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "genpdf.fitTo(data, PrintLevel=-1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "449e9e75",
   "metadata": {},
   "source": [
    "Make a plot of the data and the pdf overlaid"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "aeaecad3",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.215673Z",
     "iopub.status.busy": "2026-05-19T20:28:36.215535Z",
     "iopub.status.idle": "2026-05-19T20:28:36.391287Z",
     "shell.execute_reply": "2026-05-19T20:28:36.390613Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55ce4accfc00>"
      ]
     },
     "execution_count": 6,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "xframe = x.frame(Title=\"Interpreted expression pdf\")\n",
    "data.plotOn(xframe)\n",
    "genpdf.plotOn(xframe)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "92fe289f",
   "metadata": {},
   "source": [
    "Standard pdf adjust with interpreted helper function\n",
    "------------------------------------------------------------------------------------------------------------\n",
    "Make a gauss(x,sqrt(mean2),sigma) from a standard ROOT.RooGaussian                                               #\n",
    "\n",
    "Construct standard pdf with formula replacing parameter\n",
    "------------------------------------------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "37f9e274",
   "metadata": {},
   "source": [
    "Construct parameter mean2 and sigma"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "8f0c317a",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.393483Z",
     "iopub.status.busy": "2026-05-19T20:28:36.393364Z",
     "iopub.status.idle": "2026-05-19T20:28:36.496893Z",
     "shell.execute_reply": "2026-05-19T20:28:36.496138Z"
    }
   },
   "outputs": [],
   "source": [
    "mean2 = ROOT.RooRealVar(\"mean2\", \"mean^2\", 10, 0, 200)\n",
    "sigma = ROOT.RooRealVar(\"sigma\", \"sigma\", 3, 0.1, 10)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a9ab0bcf",
   "metadata": {},
   "source": [
    "Construct interpreted function mean = sqrt(mean^2)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "20abc5c0",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.498784Z",
     "iopub.status.busy": "2026-05-19T20:28:36.498662Z",
     "iopub.status.idle": "2026-05-19T20:28:36.616215Z",
     "shell.execute_reply": "2026-05-19T20:28:36.615582Z"
    }
   },
   "outputs": [],
   "source": [
    "mean = ROOT.RooFormulaVar(\"mean\", \"mean\", \"sqrt(mean2)\", [mean2])"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8e231978",
   "metadata": {},
   "source": [
    "Construct a gaussian g2(x,sqrt(mean2),sigma)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "2dff2954",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.618337Z",
     "iopub.status.busy": "2026-05-19T20:28:36.618217Z",
     "iopub.status.idle": "2026-05-19T20:28:36.729773Z",
     "shell.execute_reply": "2026-05-19T20:28:36.729103Z"
    }
   },
   "outputs": [],
   "source": [
    "g2 = ROOT.RooGaussian(\"g2\", \"h2\", x, mean, sigma)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "7db531d7",
   "metadata": {},
   "source": [
    "Generate toy data\n",
    "---------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "6a2b4e28",
   "metadata": {},
   "source": [
    "Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian\n",
    "dataset with mean 10 and width 3"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "b2e97907",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.731940Z",
     "iopub.status.busy": "2026-05-19T20:28:36.731817Z",
     "iopub.status.idle": "2026-05-19T20:28:36.846301Z",
     "shell.execute_reply": "2026-05-19T20:28:36.845811Z"
    }
   },
   "outputs": [],
   "source": [
    "g1 = ROOT.RooGaussian(\"g1\", \"g1\", x, 10, 3)\n",
    "data2 = g1.generate({x}, 1000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "fa60cade",
   "metadata": {},
   "source": [
    "Fit and plot tailored standard pdf\n",
    "-------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "b8007c29",
   "metadata": {},
   "source": [
    "Fit g2 to data from g1"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "67e731b0",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.848036Z",
     "iopub.status.busy": "2026-05-19T20:28:36.847919Z",
     "iopub.status.idle": "2026-05-19T20:28:36.961560Z",
     "shell.execute_reply": "2026-05-19T20:28:36.961132Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(g2_over_g2_Int[x]) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 197.932 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_g2_over_g2_Int[x]_g1Data) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n"
     ]
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "\n",
      "  RooFitResult: minimized FCN value: 2551.39, estimated distance to minimum: 4.39288e-06\n",
      "                covariance matrix quality: Full, accurate covariance matrix\n",
      "                Status : MINIMIZE=0 HESSE=0 \n",
      "\n",
      "    Floating Parameter    FinalValue +/-  Error   \n",
      "  --------------------  --------------------------\n",
      "                 mean2    1.0010e+02 +/-  1.98e+00\n",
      "                 sigma    3.1172e+00 +/-  7.12e-02\n",
      "\n"
     ]
    }
   ],
   "source": [
    "r = g2.fitTo(data2, Save=True, PrintLevel=-1)  # ROOT.RooFitResult\n",
    "r.Print()"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "094f634e",
   "metadata": {},
   "source": [
    "Plot data on frame and overlay projection of g2"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "8f56bcf8",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:36.962954Z",
     "iopub.status.busy": "2026-05-19T20:28:36.962839Z",
     "iopub.status.idle": "2026-05-19T20:28:37.067319Z",
     "shell.execute_reply": "2026-05-19T20:28:37.066813Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55ce4ae96040>"
      ]
     },
     "execution_count": 12,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "xframe2 = x.frame(Title=\"Tailored Gaussian pdf\")\n",
    "data2.plotOn(xframe2)\n",
    "g2.plotOn(xframe2)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "d0886a7e",
   "metadata": {},
   "source": [
    "Draw all frames on a canvas"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 13,
   "id": "a3bd5b49",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:37.069009Z",
     "iopub.status.busy": "2026-05-19T20:28:37.068895Z",
     "iopub.status.idle": "2026-05-19T20:28:37.305798Z",
     "shell.execute_reply": "2026-05-19T20:28:37.299355Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf103_interprfuncs.png has been created\n"
     ]
    }
   ],
   "source": [
    "c = ROOT.TCanvas(\"rf103_interprfuncs\", \"rf103_interprfuncs\", 800, 400)\n",
    "c.Divide(2)\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",
    "xframe2.GetYaxis().SetTitleOffset(1.4)\n",
    "xframe2.Draw()\n",
    "\n",
    "c.SaveAs(\"rf103_interprfuncs.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "049231a6",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 14,
   "id": "b63300ed",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:28:37.307503Z",
     "iopub.status.busy": "2026-05-19T20:28:37.307331Z",
     "iopub.status.idle": "2026-05-19T20:28:37.489371Z",
     "shell.execute_reply": "2026-05-19T20:28:37.488817Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222517480\" style=\"width: 800px; height: 400px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222517480() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(44304,'WkwIqTgAEK0AeAHtvWtzHEeypvlXaLDzYccskBvu4XGr/ETdjnpXF5qkPqJG0yYrkQWyViTAAUCJOmP939cej6xCAQXqMq3u1ukRTSilR0ZmeNzcIzz89fxfJ99c//hqc75+uTlZnXzx7vr8+/XVl5tvPz9fv7p6fnF9Ek7O/ny+/Z+vN39672QVw8nZO9vrq3H16bf/3+bJNeknZPv01fX24nwh/t/t+dOTVQonZ/s3rf7XfWW9rYCUTJNZODn7aHu+effixcXlyUoW8vPrH19sbsgvt0+vnw/yg+2LF0tmmIXcZY7RH9+cXX+8vny2PT9ZxYmUz7bPnt9Jeufi+vri5e1sX1y8up3w+GwLExpOzr66uXw4Lnnx46vr9TWl9E6eW9TDQfHwB5frl5u7fJN2p+L7fLcrtM+6S969kscPW4ZXvnNx+XRz+fn2P5fWO0j8+OLpZvTrYxn//2r5/2MdbfvV8v/H1xcPv716tH2zefHd8sT1xX3kyar2QtVvPZCqN9nNI7cTTlanI+HPt56iN/9888whuSvmP249MF7yHzeP3E7YF7OryvXF46U2XpVjkp6MUTQXbWKxxpSlhZODF3z13Wio5QWQ8LkjT1anvEFLa1I1WooSa2agvLj44ZP33l1a/ZB4/OdXfoNO/erg+st96of7q4ffXt160cNvr2696+G3VzePPfz26ubJP7956dOBJv7x5vLNy/WbUaE//7i//OL55np9skpesefb5erh1avNk+vP1tfbi1GLT16//HZzOa6/2D757s3N5Y/j8qOLZ0viRxfPbtL+c9x9tH76aL09Z/qEk7N3Ly+urp6vt8sL9+Sji0UUHY5rRtegb4b0xxdPt2fbzdOT1dn6xdUmnJz9++X26Zvb5I835MNvr969uLg8yP/+0+31+lsm/PXla17wwfbN5umteu9e/ehy+3J7vf1+c3Uk8z7aXiFSd+J2IdeXlyerr/8STi5eXXPx13By9v6bzZOrk9X56xcvwsnZJ0M+X55JTN9sz683l68uz16fP3Hx+cX2Gr7ecveT1y8frV9srq93QpQG/WTz5vo49b0/ff7oo4dfnaxO/m13GU7O3rt4/e2LzTuvz852HfrZ5nq9Pac1l7Z4fLX9z82fr3b3v7pN+t3PNusXJyulcL89aIH+cnv+9OKHLy5eMcUO6a8O6UWU3WT4cIPgXgbIDzsh8O7zk2Wmv7u+vj7qgIfX10PHUbPH72yuf9hszhchfovyRv3g8uLlFxevTlYyMaoeP11fIyOd+GpHoOQeDkL+Gk6++/ji+82nr9b/8/V+rHz32YYWuZ149uH22fOPqMKir3zcrq+fPN8163efP7/44f3vN+fXn1+vr19f7Yfndw9fX18wQPY5P96cv35nfTlohs/DJwzA/RNnn23WTz89f/Hj7omzL7fXzy9eXx+O1d34/XB9tYy+Xcphrq/vqPHfbLGArn/rYuHLzbcuD7bnz962YmBkvPtifXW1TBbyjSXKYcIrBulJXGnOYfmbZRVDDHFWT+UqraLfjbP5Pc15zvtn4lyWnKTXJS/XbdV6UJHQbO5cNwsqZRYKtLD8zSIr6RqkaZCqs+iqluD/zZJWEjUsf7PYSrSG5W+WvJKcwvI3S1lJbWH5m6WuNFpY/mZpK9Uelr9Z+kpNw/I3a/TM0nk+ziq3SV1JT0G6BClt1rSSVoI0yDSrjXJLDpLjrHklLXEnSM6zwtXCpMVZ60pyDmIaJMVZ20pS9gq2MmtfSc1BY4GTOcWVJA1edElzkpWkCgdBSpyTriTFILkFXpzSKKe18SxNpaOcnOaUVyI1+POW51RWEg0OnOVUVxJ78AeSzqmtRM3bVWOfU79p9ZxniyupMXh1u8wmK7rYeU55Nl2pSvD2gaSp2mAyxdls5U1Mr6nNllfeprAhfbay4jXdQi+z1VEMDWV1trbyEsXZnK2vxGrwSojOOToTPsTmLCuxEnzMzlm9HRYirXrbXduq2u46r7LsrstK6+66rtLusjGwl9f0VQzGAM5ziasYShvXTJ0u41pXVMqn1VyYQZL6uMMkYnT44xmi1UEUypBRSPGZlJZSvHSfpDIXindC61wp34loc4UBCGlxrnDgRC5z3c1h+rfCgd+Jba5wANHyXJ2BnEORue5mcqpz3RUvaa4dPUWmOLe4GgxyLavs1eBaV82bgWum8JAucW6MSm8GiDF9l1cxTcqOYOruC2HiejPwDNPWWzvEuS+CKzUIGRJJDGIIL+mIr578jpQCYYMwhRhSjEEV5z7EWPc3V89VhfTm11a57n6tKcRZYlyJSLDsQkGi+Ag1DaazRB1CodO9MktMq1aD+Kjhtnk9JJeQIPNKtYVUQ7ZZYllJjqEyQfssEZFRPK8UXsX01JBSUIWkSUKhVWUWl62ZOaSzCGJMQpKgUOqtijAqfRZkK9f8aZ4F4Upv0vC8V/LgeIxnQbpGCb0HiXEWqSsEHPWRNgvSlXFXehAkel9J6xQss2hciSGzLEhBviPHRrs1SF35Y7GG1mfRtNLk5cRZEK0toxYC01yQrZKCizIlc/GprtIGz0jXAxWgUVZdmJJdZqVLEInjb9aIbEph+ZuVPhHq43+zxjzkLWJ9VjqlJRcBcdZYV8zNKqHy5rbS1AdbUmeNfcVTCOs8K13CKODhWUVWVcJ4eFZRl8vL07PSKZTHw7OK8eR4eFbJPDkenpX+QKJL5elZpQ5ud4W2feW81O787opF3SGOBtOzepfA7SiXLhl1lZZm9S7xymqyWdX2tfWS6ZSltjCtSPGbVtW6q+94uO0rzMNpN4eRbdHmtJvGDCb0FjN50eKDHpPZ70edE/N5d9/pgzWB02NWkz/LnJjYS3Yn234F4SSybSw3sswWbxZGcbY41gQ8DqX7J6Fu1ilQhjpBo8wWaY0hvaif+WR22e/1MSYzNNo4p9kYOlH39TOJg+Z+1NlYLMEiGtRppayleibJqWjOvzgfZKU24kqOh53aifo4m7ioH/USRD2PkI6i4YE4m6JnRpV06ZKxZpyN5dBSJVSaMWmXKnmVGR9LlcZ9JIl6lcd9Vh9LlWgCbbsqjex9VyfPnSIkBXjmJE6OXpst7ZSeP5pc643+ni254t0/6Jp3/xzN4Z02Z4bkMkZ4S+6y7+pB33Q9HOV+0/njvo36jDrOmSE5ynGWcy+wPHptzr3uKH+X131UcM7dqw4VdS6RLnCtOJc4Vuk8UyJ19k6bS6TG3mlziWO1zoqBUThmmLNYEF9jhi108yov3TYXdAr3R7fNxUWYj8BRJKNwjECvUkG1+Awbr5Pkr6OKMIhm8UG4ZB7tsauWHG4fyrJqR+rx5JgdS+XQKozG0RzqA4FemwtDclc/nkNiLfVjlBRE1phhC70o/qX7i8usUT9val+l+wxb6KHexgAp2nfVo2xW6aN2To35ueRMN5uokkajsCIprM6ZY349mmNc+9Dgsu4E4WB7rosgXJp9roy6hcUsc2XQLSQP5319/eZN91O9uixvnIW5LgscSvLMo3bjZtvJQJcODRm4tGuWucUhCbxT5rYIQRcSLe77nDsLO0P2tGU0LvzObTcaXdI0xuIy9ryMfWOzxhwSeRSxDEIykxFVuogVcrokXF55Iwi5s4ijUd2GHt0J5KhzQ4+OYe+Z257yMvayiEWuD0CKh1i0w2iduaFCd3Kd1yIkGW1DGDdk5ELy8LISXlhaRuPSCnNbRuOupDFXoWCJPeNSEIvtRXYNphiMS92YFo0t49JETIuWRiONcTO3GwEJx8n1BTyNvIvUXhR1S7v9wch8I7O9IN+i7B61wz1KM5fYvMffay6xx+Sbm/lOaRknthPY1NNu2IE7u2FnvGbPDhWzPTsj844dXpR33DixtNBQ0XNju7j0i1ck74eRvzfD0SKd29g1Lg06tzw4Whp8bhmO/C6mkOwc3ZBwtMjH5nvIRR+0MhTa6La5jX3k6LS5+UZyEZxtbCQHt3MrQ72PoTC3cqBkKd+3kzdvHe2zq7XvKPePji3l7l5dTCE7huqi/XYs1TGwd41Sx8AeErvVMbB3Fa1jMbRvhcrOZrfDbnUM7F1DsNEcNxlRbQzsof5aG922y9pGrw2l0Zr32dASrfkaaGGHHaePrdG57DkXmelt1MbU35Uxpv6+DJ/6SxF90T3+2CKpdy2wCOql/GVx4GqqLUJ6IYZQ3PXmIqP3/bnI6N1r9grIB/Qio/1FfRHRC3EjoWm3fiCiGcEdIb3o4qhzX8T0KKUvcnqpSd/J6aWr+yKnR2P2Zc2wo8aKYenpvgjq5eZiXfMxBVNsN2j48Td3hCXjePzNXb0Hhw7taquc2FZL4lEf5sbWuM5dyypH34Ebm/66ss5GXYyXtpW1wMbTyty1r6yEHINYm3uKK0N8Bsky9yQrs+DM2NyTriwFF0917imtTAPb74IlwVYWAxa3onNPeZV6yB1b2NxTWaUWyjAR9lRXqQa2yJUS2irlwAa4UkLHbFJKkFrmbnGVUig1SO1zN1klZTfP3ribrpJv+aTluVtaaQ/shVudu9lKmxsJOmaNjNmqZqyFc8eSVkKtQTolVCwEtQfplMBmNTRGsszd+ko1tLGu625HC771q3NHOPbQMjapuSMaW2glqKS5Z6wHobWgkueebSUldKxXbe7ZN4QYwlTmjlzMoetCMeXdPKG8k/UdJpdB+ba5l0EV1lWh14ViyITOPh9TDTYDbBEqfe4lQQmGCyedGR+1TuZxNy13yyCH4aKXOkgq2eZeGtWSSC0hMcgGTBNOVkyfAVOpSp17FScxhzmJySdgEhqkt5CILiRW4oD1AhN4RzK2IG5rgMS6MOyFTtZBUi532cwHwcpBayMZIWEjz71FJ7EUOCmDhA3uYrt2o8cg0yBhw+aOQY67VB8Sg2wQpfqQZZCwAVkHCRtp7ghIMsMGZHcywUaaex9cYXsSnTubqIb1cSFHW2FvEpk7UhLrNGxAjh5MsAE5ejDBBjY4H1CSYAPSR5RbV530ISXG+O5zZ0tvQbBiRWxpGO2gGfANGmMUAmNHY1UPYpgj/H5aCRboEjRW8mPBxCK9o/Og247GLBPEOhZx8lenM+w4jVZB7ASNmftuTpUMP9ASV31RxU6Kk36MwV1dYUwfthKJkpyEF79rTsJK4m5eYdKBEycLZIERLI4Y8LAn7sjmJGz43b7CmkerYGLU6M+6dRBSKAhx1mlTVZjk+KPTZL7UdXOrdJqMta4uNHVitYsNEukE227EW2gY9eXuOObo8ILJCM6RX87MGP4uIymeDVjn2CZIo/zEimWhKZ9NGEdLfioBjYVx2BQb/CRssQsNP4nN4bCLYp+MyU8E3E46LLkVwSuYZxutkRpimYMbJLhExDvlQ8MPAp7yh7yXiIinfGjKR8hTPhqA8i2hLIYdk/YwQ5UMGn4so2gEeV4p3wpqSBptTvlWUVKD9vIbKmycVtEe1lFwg6Y9ckT9+fFbgZ8sKEc3ixbKz4oiHTTlL2tiP+OiP7KtMs+j0CmfdXEN0tg4UH4uK3Sk015eXRX44b6X11aMHaepf+4rjriw+7qpvEQsiIOmvCIrbNbcN8orfmA4aMoraUXbcEJnlFdsxdGi05RX8oqynKa8UtzQ7TT1K3XlVnvWEV5eW3H4B+3m+NJXndM/aMqrkVOlQdO+qIUIw8g7CkQxsMHxBEr0NTNN5BKRBFuhRkYCZaIemDSt+zGhRBQEVlYS1EtFsFAsCV4sxwIUS4IXy7kdczT6YZtEFIVytEgCPYuqYOh6AsWiLPzEkRyMLdQFg5GTT+z8EYXBaCSHUCwqg+F3k4BNnlNLctD+qA0GmCfAGIrDFj44LYioDoaY54AxlAerM0+ggXY2OE+AUzfCLe0hcOr2kKU9XNK7EqFNObGAdVcjS5u6LEeR+MgjB5yiShBLngCnKJMhZzhqJIEFwdK3JPghkEuCFgQx68dATG1P4NRldxLkCX5Cw7nD7h0cpaBT/CCiBeQpp0E+mqE5W/GjB6rWAgLVD4RoiFY5KBU/EPID9BoQsJwISWe81YCA5UxIOtWqAQHLsZBwjgPtp0SYcWBoRw9LldMw7McQVIklHjQ7F/gpfmLD0RCnluRHAHM2xJGw0/DH4ZAP5sYho3A6xCEL9UHgckCkPrZbaPCHLQWJz3Eu5WPdW4Y6AtePiJDwPQYELqdEfuYMTflY95D4XUKlvbCn+ESQgMDllEiR8NBePieMjGkNCFzBvoeEh6Z9MKr4NEkBgSsoGI62egrjIE1X6rPGAgJXXMEw4C0gcAUF45OIRT90Xmli/LMlgC4rzjCdpnwUjE+pEhC4goIxyisczIqk4btAnyJwxeJKfYLVgMAVjuSR4NCUz5m8z7cWEMDCPsKnWwsIYGEnYb47CJny2Uv47OsBASzsJnyB4ltMEfYTbmeIIdMe7Cic1oBAFvYUTqdglM+uwukcEMiS2aOytCgBgSwZK6kvwAICWTLWW1YiPRj8ZM5sWTiKH61Khj/oFBDQwpG90xYQ0JLhj/slIKAlwx905TBbJMMfdA8Jfgr8sVuJAXktBf6gNSCupcAftIUEPwX+oHNAWEuBP+gaEu1RRvsxhjiklAJ/MShLY/gp8ActAcktBf6gU1DnB/6gOTbnqBX+oEtAbAsmGafZDUHDH/dZiUPTvz0wZhDiUkf/4usi8FPzcp85CF0WmhUJNPtFnmfVBc34g+ZwF0eid16f4Uv2b+6qd/LBi4v1ddKTcPLCvcdyDiffn6y+7mqhK1u/EroiW1ro2kNPMXSOeJOGnlLoyUJPOfRUQk/IqBZ66qFbDN0kdNPQLYVuFrqxecRHBNnWQrceeo6hZwk9a+g5hZ4tdM7Acwk9IxNb6LmHXmLoRUIvGnpJoRccTXLopYRekKUt9NJDrzH0KqFXDb2m0KuFXtmKltBrDb0ih3voLYbeJPSmobcUerPQG5vUEjqy1OV3D73H0DmD7hp6T6Hj4NJz6J0tbA29I/f9XB07DJI4uiUkcgLhltHIWQwiOLJsiAjfiMNOROxGZG1EwEbW1RHRGpGn0XfFSNKI+IyskiOCM7JOiIjIyEo4IhwjEjGyPIjIwogAjEi9yGI2Iu8iQi6ydo2It+iOSwiyyIogIsIiciuy/IxIrMiQi6w2I4MtIpUiwywO5x2eYFRFDDMRmRMzT7hVjkVlRMdHzC+RNWVEw0dWkxHVHrG2RHR6xNIS2ZFENHrEuBJR5ZGlZESHRzYIEeUdMaFE1pAR1Y3bMQssfngCB4LICjKypI6o6sgCMqKjI0vHyNozsmCMrAKjq3ZfNKCWI4uliEKOrBsjqjh2nmDZGNHEkQVjpM+xWrBo4AcnB/p8eDXR5xgi0Pf8sOOnz32hiJVBfIU4/CHoc18fCn0+HAPoc/eXwjSAVuXHHSkogz5n34965Icy6HNfEA6PLPrcl4NCn7NZR3fxwxP0OZtw9BGGDMqgz33xx4YaxcIPT9Dn7rDl6zy2x6gFfniCPndjoNDn7h/nHi9Cn7uDG9tYBDU/PEGfu4ed70KHZxh97ss395DzdZtvCn3BNlxI6HNfrvnuzNdpvu+S2v7y17/+Nfy9vDnxa3+rN+eAKvwEtGTx2MSB+/KcjDsf55HyYP+COxgJyAMcxF0IxA04RIBkHKNDPl5ffre5PECbjISDVy4JewTFF5s31w/Pn+GAjQMq5LgZp0gb+P0X22fnJ1hFBn3wfm5/cIGfe3EP4/Wb7bHX+MPr64ek47j9dPv99mp7cX51ssruPc2dgxd+tP52s0PAUJ7TowSjBKc/PTu72jg0BTm7JO7ZTs739sl3H23OnwGsiVPEydn7YPeo1wXH87uPXb/YOZrvs+yKxzH5q/8yNfTu/N+o4X//L1PDfQf9yj58Z315ACx6Z325GxTugg1Ui2n64tHnY068d7n+YQAyBv3pq+sb8McgFvzHIBYIyKevrt8b/vYDXIYTPJPIp9Gnr64XkUAlPn11/YGDsZasH2wXXMCRBz0ZPPHp9hpE2Y7+4uLihTvQkzDAK+9enF9fvL68WlALD68Xdu5IzIfX10xiF1I/IQv0VwoD5grVH7AjKgkFImJy6v3zp+9fXl4sWC5mtpOenaI+eH3+ZBEL3IQ8kGKQSxdyFzDMkpn6Qy6ZmfeQB/390ebZ5vzpIeIG7kbqgYTlRTeJu7J3IDxesRML+4y03zIUw8nZh2AhNld3ZPiS+vmr9ROwAF72Hhl3UIc9LG5Jg8d9vtvc7LPukndZ7xTt+e7W+yDxBm704faKAXnID0m8b2GnRKq9y7creLTOLusudcl4hxtyfbw93758/fK/by4vbqAe3LgFVnQRP1Avjy43Z5vLf//oJvdIP2i4kXBYTTg9TL2p50h9b3P24ckqRzptn/LlyareTnl8MuA++yxfLQmP1oeD79H61tii8H3STcmedAzSfLR+eqvqtN2j9dNj0Oej9dN7cJ+P1k8Z7I9vmmdJ+epWCrpxwSpR4PbJdwtS6dH61UBXPl6Exj7hqxOOJ0/OPn9yudmcf7B+4tIH9hBrB80PyTw4GLYkHfbH7qmD+UMWyJvZQ8puAI1hNfJcvkSKnpTJHEhDosuMDpQRykUMm9NBfrnDO3HvQwjxelxfbl+9t3myfbl+cbUHFblIXlY7ul8yHNTOM9ypnqcd1o9x5IkHFdzTe1U1FiksQAYYeP/UXgT6M6M6fZ+b+uzeT3WG8Ln40/n55vIzqkdOppq/9upk9TVAogcP+BF9IGz6H9iDtqPLA3sgBUojGWzJckN50mG6xtvZ9/daXO61+MDiA3tgntPL/jv9/AXhslk/3VyirR0z5c22pz7YXn+wGzR5GTSOxKIX93d8LHmrPVm/8Ifp/f/nYntO4m4h8O761SH5xfblfjlZW+vSzCXGn16un2140V7Av7s+f/pi8+Xz7dV3m8vP1ufPFvjySH/n4s2SNnpvpDonBwDO/9hevNie71IX4OLI+u728smLu9J+uQXyFKYPFOBjltzvv3n1+HDZs0v86jDxq/ty7hJv5STjx+s3722fObyeQfjp5fXzi3fXLzeX60X6HG/XHq2f/mb4O3rsrTu2R+unb9uvsZOihQ4kFuRO+Cx3D9rvjmY83oXdEdc0xj3C+g+Q/qm2Kfq/JLV1AUEJaP9UbWqp51a0WKq5gl3Vk5XmKfk/E20YkU4c06+appyLFsk1NZf8tyH+GtuUs8YYa5dqDJQDxL/GOOXO3V56KYd3T1Z1kpJ67z0n7SU6h7dh+kmmpqn32nrKR6EBtE1Y2Hrttdy6ebI6lclKLc5zVSuopdsxA+rUM2U31dgwEv5ECIHUdIqNxkw5N7LeCSdQ2qR+P7aESccz7KIepGJTLjwtpWinwXeRBwg0cHz3dtSBU+0TGxvMSWLYvPZRAwhKcF8/74IUuEb2vqYVJXXfHt28niAEqmWKPXZJUS1XeN/FJPC3H4+W3W0PWVCapFx70tgzluXDgAVTjLRw773WYk399k3YgXtv74IaxMn6GL5RNWWh1XZBDgiacffeiF4QJxvdEKMZQRjQYCN6wtSP7twKinAvN7fCJLwtx9uL3gdSuK/0JbDCqcPej0Mr7JI9zEJvU5VSak6xW2YE/MpQC0jaP0It7EPZ/OpQC3f29b/8+bcGZ/hGDk2Xx8Ebxv1fGqDht4P1L4uffwiq3zcru2BCd9vYN89vW1cgaH7bdQWBbcZcPAxtM2JSENvmeAoebkuoyaCXbejxcuxvXIt98OdP3qU1RvCkn2y5D+W9t7Vb9gotg/KMqD/fvPkm5ycbWz95cvYksqcYi/yT1cmfRiyRzfXm6YPNm1eXmyuMyw9ePT0j10HzD4PF37iw+/Xm9U+ebF6wweSo6OTsLebxxTb+tjhVS1O8cbvsQd3fQPyXNqjD/7dbzgI4VDs5e+zBe8YQf+xqxcfC45GHCDNnH2wvrxbj6Efr3RUxvXS3MX65eW979erF+iAcDju1/V7PBwanDDfhaT6+ePrR+tuF/gkr/y/rqB/vdpTHQLl68D/+7wf/14M42YP/9vvquInVou+a98aJ5QjlJ85CDrpu33GjC3dBl/7xvfaWk4tf1mv/ebfXfl+d9C/SRXcOWw5OYqJLgPfPry+3mK28wlevX/6wTGsud1Oc611ALr9ebny8foM5+V5NuJia2XaNE5hPLi5f7uyYDNzluGQEsjr7fBTnIudWYEBOHogTeHzK+csXXCMW1Y2tcxebykt7Z3uOqfL9y8tPiaYFa9Cffr+5PHtx8QPHOcRHurxEwN326Cm4+rhHj6sb4j65X8/dhcFvf27+isny82r/s4sLrPtkvk/ZubvKPljY82+ebc5fPeXcaM0TLqFOVie84eLZ5frlg4uzBzdZvnn14uL6m29cKR7ofcbR36D1af6f1vos8+4eqn/y6mJ7TpjJRa39VEfF4Xp1Kn1qgV/jN5xKm4r/ajiV6veq36vhVIrfK34v+73s93I4FfN75veS30t+L4VTUb+nfk/8HtL/VCScSvR7BAc7dVackx5OnY9GsnPhTNRw6iw4B86Al5/DqZfuhXvZXnQKp17wKNeLpVT2yLf+keQGj5t/JMkNyb8QJ4w1N/8kxMltF/t/EiRQBGVRKKWnAC8wBXewmQNMwz3VoD41UDuqSX2peA80A+3hLeON5I48vJ3Xe2N6u+LJRRHe6N7+3hW4dVGO95N3mfcePl4U5l0rFOcdjsMXJTIaPI7eV79g7OCkhKOXY3tA3fFXHOfj3j84+kDj3IOXj3sJgm/AH99d3geWA2cdnIxwtsKxCFcwHH1wGXJXIjysqDfefqDaiqPgcg6pBOA7NQAtIryJ/zZ+s4H9wzGvEE0MLzz8rBSUWq2g6EoC84crJB6LARgUrp5+3YhIk3IwAYtXjShjODG7q9aIhuU1x+vIBtdwC+d4MeFK5rWg9kQ8oTMA+/BHgB73qx9xX7jnLUen4iTlPlL4ivFu74zfQvjvxdc+OOKRatrrsZ2uuqWB3idK5s+PCo/yh9dO9J33+4+fb589//WPYQv7+ae+7lNtsZVSrY4hMlm1nHttxZr7z02tp1SaWCmVOSFTTFFjrard3Nns7iMyFas5W29WirvBTclSbBpLIZl3lGjZai/W3RNQpiYtSotN1H3mjl/RW4m54onn7pQydTMtrbeeGZg64c2fS7KiFc88napYrb1oKu6glyaVgltdtOhO6zpJllq1pNYr7o5HbOnUY+/RspZqeB3qVFJurcSYuuE6KZNotN5K1ZYZeWnC51NTw//THQsn7a1LoX6ZsSpTj5qSSdaWcRXUKUnW3KxoxONWZKopxV5waMTrUaZcJLXYS7JMhqNOiZM0KRpL7Gm4g07NaqqSrPeGZ+xUtGhLJaXScaduk5WcFdt3T3g21gkmtZtm/D5LaJMUa1WzWnZ3xjolz56zlm5JQ556ld5z7TH1pAV57KZQqbUUSxIyZvhaaY+ccyshTymnniQXa9nh0Xdz2FRLxCKpeI+m+16huVegLLUnK2iDO4XUKWkuFlORVmqv6IhYO+cBtcbaUBF3an9nHvRJrLSqBa9HA347tdKTacvVmuGyPOVqGkV71wa+QuIUm3bNOUuj7n3SKKkavrjNEu+otEGT2hh2JBy18J0+KBNNW7um0qoJ2k6lASnEG1VLP25Qm1qxIqoxNctag03Zcoul954Yd+jQbDlq6tnMerkvx5135OM2z6nmkpPmFnOFj9YqYqTG3DJif0q1tF5LbJEWrFOvPZpqis1wHm9Tbtos5WrWcjL0dm+55krFFNfuidbMVgrsttCnmEpriXLNcS9xMjWtsbXUG+pUptRqrKXkWDNeyUfSKU4t4oVaas7D8XkqRVLvrUZL7t17PGPvzsejGRy1Vk4GSnevXp26tVxL4RDK3XKPH7FaCgcQGnGTF737jnukkVbLHcdeRBJ1ayUyXc2q4tOsUyxdelcmJRpfJ02pNoT54tw8xZrNssZe26hs6dEimVJzP+Zj2XxHJRyLpzsS7h7xdKeJI4NDOicqTHRfnx12wghF/NUvU31fU2DSGKvRj7im65QRiKlakYSbukw9de2ZA7mKG7bSELnkln1a3POI0jJZNPZYKs7sytiSHmu36OstRUuVYrFZdqdvpUNSLlpj7MAEjt6R6OUWNZWccKRX9JqUKMaMHfojS6vaWjRgFpKmmk2SWZKScVK3SXPqnHpl6SilNEnXbl2qdkczHfGVGI6lWM9VDR/yNJWarZWU1PFcOknJ1nrMxbwmNonFJBmtWz1gKTpceo29xIIzO0PcUtXK0SGqME2cYPKvsBzTqdZcLaYuxXxVy8AquSLIDLDBUZ/IpCnWrFYjJ38+wHuyri33ZA6smEptaqmYdAeuTGinLsigBCajIWW0mhUtFcREn1QL+Ax8/auKi37c2FkX6BBMESGTtRn6E+3BoiDW3Gov2bVaiWqpuc6oFsrECaD0YqrVWkYq386Rp1ayJis9lwbS5ShHmRLSt0TRGDmgPsrRULfNetcCr0jQ0lIGmCLaGQpHtT+eCOy5JFtryAEQPlNXBkGrCajFmHlazJK2XgrDSSZJmlNNNSd1xAdiuGpNxlmv76pEWSESUp5QEcdNfKcTkP0loSpjb9rZqmlvnJ5ndLmm4ybNUy+Fs99SQUGgYAqaVXoVaSCf8lSZ3MWigZhgOXKU48477mn0opU+SK02837qVS32kluKSI46GcVmQZGDv+hTNCldrOXe2Jz0qRi4EmYX61vWErlEVm91rE3jlDjCjrnV2hw9NMWuqinlZAMaMVllnEqXSnR/0ckS0CJpgp4m4UiiNavWq6ig++i3GltSKzG27KFGj6fs0ZQ8msOx9ITwKuoiLk29VVUxyb0jfe6Z9lkKKhdfCXab6e477hFHiaUwggQt7S/1CVZYqoMMSqzmC4MuFvwvWLxXK62nnrWOBhKNRVqyxtqTd9SUNKbCioQd81GxR014j4i6I+TuEVF3GllYgfBVhtxbd/jQ7Rxosm9+xJbtZxcnq3tPFb7Zf/7gm+WTB43Fz17M93DyzfnF5cv/WL8gvHt08uX2fP3ine354uIUOQ/45vzz7bOXa7cgfrPZG2WxsISTby7XP+wNtacSTr65uNw+4y1fekx7t40eQg0PDZOLvesPm8U/wmbxmx/xvvjFht53X19+/9az8aw5HR7zDjPuN9jkv37zl0Nr76PLC77Sw6Hu3tzL/QMTb3Er7duNvHh0/aTP3s8bebEiv93Im91n7CesJGQAX3uqcTLwoX7Br1t8sfN2t8k2N442t8liDa5+r/q94veK3+OJ7Pey3zO/Z34Pe27ye8nvqd9Tv4ctWfweX4TA5kspzPhTidh9Cz+K0bfxY+G0Yfct/ChG38aPhdMSdnZfvzJIDcMC3PZmYONKw6mFnTm4cOWkhmEYbnvrMC1CMpyRcdimSbnXVNzu2IVP41SPk9ivHvx7i43ZDrJgYz6Nk91OI8k9AW/+edJxLr3JwT9y6XGuu1ZtDN3HuXbeaeNf85fFW9lQ2Wg5NsRqRbOUzankAH5vbxDHsQ8zuQvw/T9PuZ3pXuO6HuVJ+3fwj6fSUR47ymNHeW4l3GvIjxP74Zt/lLXzkhv/SNn52u1SMP9rYBBh/OfP/BCAucBBgC6HAWUcCCxHAm7SwPDhRwMsBDkewNjUAst3FteNNSdL6MAhAcYG7AnK6i04tBvjPdPKZ5hPtrGB85UEhnQYY36ya2bG+uT1eexTGpsjRfl896nvUoC1CuW5iHBpIRQ5go0jVvwnxrHl/alTA4TR/mjwa6HkpLmnmqS4lfAooU7WktVWVIbZtU+9mvSSSuSrVKI6YdXFRuEgY81T1pgt1Sjmn/9ok7VWpBk2qCBJpqxSLUdpbkNN+J+WmmK1VMCwpzplXpljKxXgssXJSiqF3WrHJmHsV43dWXe7veUJ/90arRcDlWx16jHH1DHVsQExVtRWldU04QkkC/ZOkRa7AFfOymqwYwvKxD8gIUtKlkruBeR51onPXjHfxDDIZGGFbSkaeykPcDS13JT5YA7Jtzpl05IRHg6iN2MvELNhaaWuFieW/MVyzY61T5kRXdlktwymO8XJEqy3NiqvNmk0y7Ekbb4JrggZLVZKz7436FNJiA/rWL0IJgCnPZWUayEwg04119gKi3IQ9S1NNeVMbAAtcFqNXkqp9tYVqH7JEw1XBMflWkL2zWYs1quV2nKwOhVrrZdazOMOpI75PJtpjIl3pjLFVvDlNUsEL0gsuCObxpQKw0v7VJpvFKWWRCzD4puHig92IVSE2lRqssRoY1gG1ak2Sb1iHCYghXJUwEFAaTkSUMAVb9fYWk5NOT7qnD4may1xFuA5JMfCbsqIlaNYnxQrqfVYiLqmtHCqPRtpBTbUVEuphb4xGG3NMvbLKjn3kNye3lPp0rr1HJJNbHmLdC2tthYSXrOxsfavuXXCp1CsKL1QsTMYJveI5bJIbk1CzhxyaLbSmipHaH2ShI2ql544W6l5anS6FNVCHIOmE1NVeoo5ltJCa/DRqpWuROcJPU4YYEvC0I0XdNepSM2WenNbYeg2gctvEt1gRICPqUtj02vZunpCw2pl0mNJhNIwNqzF1VUj7EeaLNcSzayZRR+CUrv2mFUrVelxSuzncVLHwNTq1JlYWVIVJmCzqZTEZPRWJirmpNLp5qytlBJqmnqOPRVRYAUcUk4tp1465w9E+IxTwaRdOXugp6tOqtJbNAauhNpw3paKRZHhxryQnnrli3cEM2h1ko7FC5t9iim0PtXamrVO64iFLhNfycsZG10hkkqaIgc+GEcLsS6674ixbFlUj6OCfb2UjAyHVVqvdo4n2N4T8oEcMZdM4AVOPrw9sRVjCkmJWCaCHE6dVyA7YcxytIaJrhAJpZWpYTBNuSVMKbRfF+pZEbO5hVIZCVpMskn275ug3JOoYEcRDdYmTd1aqqnQmcEYK7lpjZiMiblTp5yF0dEzLwspTZw6qZaWMJkzybW2IppirkQ90TJxEqFm3TRx+psmbGFNU3UjEXPaq4ZtRBqfMYlTxkRf6CsPFNE5PECtJM0cGCKKxWKPGhN6SLEXMwU4a9NKOMrJsBRg587YAxS7RLQuLRVrfKLFpq6spCQxPznYnsBPRI43OShkkiPRVJvkiApIKKZSsyQpiegoiXOj3hljRmQUi5O2HkGspJQ006JYiTqncoVgM7lMvSAmc9LK+WkpU0yKkUKlZDHmuBbkCuZdzrKbTVa6xdS6puRH9BPWS0Ku1DYCxEyxIySsYuokOAkd2XoqqEmCdNSpdquaaxePF+aSFou3tugBvhImcmtZe2xo0FSYO4aMUw8EY3FqUnr38ymUsDE2UqxJm78SFYFlTHJl1rtWxoDMOsASlo8sUyqpMaEjkYhcLRfrJUfObD3I4JQEqEhujEB/pBVWG9VaxR6e49SqcoSUciGWjHUWpJIScXI8XmidMiFRoreqBxCdtEarFSMbkUlMp2qK9byKRyNL6AzppUtqxU1UlbM6FUPiY9pPNkVUIcfFhSVHkikqk54zAiqnFb1jOVcTL0XzFI2HShfnFD2DFbClWo1INIJ9thLzJ1ZjuSUVYdL9XBNZLSwhkxFhR3txg/fPJPxuvBlum7Testk/tmrdXZfdtmr95laPh4//9PnnDz9+/+3Gj3+73Jzhd/ibF/3tC7BBv8i77otH6+8dnEz2+/zrfkOkIyYThyPgtSQ1tlxyz8T/GxA+jRje0RhimMTJDkQho5xNGoeHKTliCmCfMqFNOZboBGji5QsmSiqrdU6pUqrRn9jd6tqlVOXETkp0YJcuTzUmC1+4xXWkY5/5aner97HrtIoywpnyECVBxf50TvwNbD6fP18/vfjhABn67sUlQOv10y3fsyTz3jdn31OLv/417tN0BUjj+0Ps7G6O+DoKL3fB5fv4OjfRbxac80cX5882eMA7MP/m284jIs72/Fd+vvUufGh9vXFfyvsGEowvtfQKek1/DozxS5th1yqHsPt9K3BzCV2wxCw4MAT+tK8nrfsRX7YdXr/D/ZaR8akbr/e+7zn3XFth7aK5y1/337XlW7bEnPrLUdSpP2DMN8OPTjiKYfEv/a3x++CtDmPunO9qY01BiD4G4E/AmFueonJGiMRii+Ao5R0S9ruTVYnl9whjZnvFaX1WERb/d1HMhiq4wTH/fkHMvrcsbKC0i/mXzXco5d8AxNzyhI+hFgwIbMpvg5jvGSu3QMyIIpOOK1TGR+cWiDnvzMQOYnZtt/8K+72Q4BsI8zFM+e8LYT5m9XcCYK55yrXhcNlyIyjmrwYwu/r+41vxdz4Of4f8iW/F312B/OJvzS9LkXsAykS72QNH3nb/DwDz/duFfziA+Z4peLg0Z2XxXxzAvCGGyC0A8xfr7YuLy83TB/++fn11tV3/gV0+2G25q8td+iCq1j8zGChboD+wy//sKK5/YJeBRP6BXf5HBAX4x2CX71GCRwCx/4Owy3eb47d3afs7YZflp3HLfvsPzDJ+bH9glvn3fxBm+dXJKuvwiFQctIKEzG8MLXAg6Jhpvu4DpppPrvSQkqOAY8iOF841pMb3P62AHeYbooXPiCYNqfojfGzI35MCuO4eSpBgIQfx8iTI7+YM7PeH6D3oIEBs+MDkEvH9SDjW/UxCnJoS0BAgEkiGAECWU1Mxs1j4WMfdd4B2bQUQAY4bfLVw4nQGqG4BiIsHX25iNcJIwaMi4e2VY0zSmgJpTFPTnlqMDv7jy0UczkpLHNsD0LEpaVTpoLkyDjyEjcTamlIkzGTIx2DEI5Rp5TBbHHyZCsDUI7jiHehq5fC2YFRqWvk07F1s5t0MZXJnl8QnR5oWntBcI6f+VSofFSw4buUWe7fa+VjJPfBOAJNaUiRsJV+9On4Ex6eKq0X00JQZCIvgalLx6uGRmLOYliKtgzg+ah6bLJlZJUAlwLP7WvhOn9xtcp1S7T1Jwg0o13rczzhkJCvulYa7x/EjR+NNpi6WEv4nqsoXb+8OwHseuTOqj3L8bMKvwhrKJAaCzQrQ2aQBqL4LRTlZfQ08vLYM5txy5lOvvyABsKaWDLAztxzSlJpa1ZIVDHI/foe5B4nhjhYjKDjgfRlvQOmtAEJJU61RQTpGS8bHvSb8nZpDHlsTnGV7bjU1S1qVkKx50phwl2odJCi4a+NUt0gRa3wC06bSAee6lyGllGPw1BEurkXFP7QAkIz44h7Bq+6A7dqUu8VYFVihCUD1O2CyoxwVTLS1BO4HBw6Qc5jdcy0C0BVsu8dOJVQtSHYHu98GpLWJc+XEiC34Wd37SDR1h0p8ZAH9VXxhEv7SFclWEQA559gEB0Zm4h10WcZVLEtX7dk/0nNPI9/pl6NWT1PGx0myZVz82n2djZtrxKWXpmdI3XlEjwedFK3IEPc8JN7K0Si8M07vGdpHOX4m4X8TG1Xl1pBw1NPfho3626FRB5r3n7o0+s39WP5ewCH9edDQcz8FOPAT+KcDhvwTBD8FGCLDvyJgyLFC4IEGRgh4EHggBwUB/BlgIHBAewiQo3/uwfXcxeK8BQRzjEK5ix45RIqw8EQSsposATwIfw2t5kiQ+1EgHgfEw0U5GmTBgwDgBVtNSA8gIWxwAHAQZ2tBiOxBIh5MysNOgNkFc/q7gYswFG/gIninCyAjyzVq1g2xyu5NlCkB9Y64y3cxchJ6y3JKsbWWgBiAVGI3YAU0duzEzyl9c+qLdcK9aC4tAj/wNEAj+JkmfDylyubUVyCslDoup2DNPQ0YdGvgXpUIHQMQZZOxXm5K8IUaeTrDJRgL81gM1ikHpdWzKq7VMVuLsE4sMdaUeNGrEGtjc8rmFHBtyYQ1shpz9sQ+8VlCV53EpPc0xVUZ2AAu+Brb5tSDsNSuWbLknKo4bgvQUosdCHKJeNnCuxDMxvAS72Y1xWieKJPiNp+tWW/SeWf0nVRkgYXDdxTYjKFOgO7xhWtasiRPlKn12qqv8EvTujllbZZbbJo1SculdtJkiq0WllI0fW6bUxbU+KmXWCzx3QslLU8t5ZSbsCeTaqTJlFjCOAYbL+LNKSt81vsxeaAUEU8rEwOl19oyzsVkkwlPcULqx1oaI4ClfwejRNilarkl0soUK67UOYP000aaAMyoXVtjf5c2p6ygMy44SdleRjqaDSJu79ZyLPhVeL4+eTSVDhLATMgmeJEXUPx87LHnzSlDxDr+tq12SakbaYUAJ4QLIo5hMk/D/ULAGPXaQSmwJCQNXxlfVbeci6inEVgB4JC5f0jyNKA1iQASsKCexkQZm4PujvJ8UIswO6zUPKoUgLiYYpMaI77YmU9FOhyQL+kQw0WJnEQSIWw6C3xrJdcCQlBLTMry14TP2RJ+L0cgV8QgqhUWtJdm0WrrmoEFgZU04Dy1FPqbFEIDVW2F+c7nYCOhQ1IswCjEi5cJ90w8somhQWwQUCvE4yH8DzGvgBXgFN5SrMpXBdKUjdhDeJBaacjmynJTiSchuPwXYFGFGBFVCyFc2hQBIuVSi4ftQvQK+CapDGdC++hUCR8QPS6To76mIh4RAQdWoFNtatWUPiTGBSABAzbQUqWvHPijVqyZNj6g4Lgfj5ESe2bLA+ynSsu1pqbJP3KM0zfYDeAy1YIJLUYcjw7aIYEDSgqCB5SUWQpGQKZKSAkpkQ++Wpm0dmwRBFlT7F9TTnjzp5gSOBVrU0maLZcSM91rnY0Dob3McgXAFSclLkbNxDpJ0LWXRDwsgGY9gJCTok0Mu0psJIDiqOBGxW9raVEQ9KlWXtCkskdMgoM6CSmCouoJU0R1HhpO9g3jTHUmMXNEtodupqlT7rUUhA7jiYomH0ydPZIQT3AC+pUsAugCgJIcVg3sqkba12Sqhn6KROPQBM4KHj3SlLVmdAne/kbAukTooMS4zj5trJgCSEE1tNSLxdozvY4Lsm/dkHQMi05gG+cERsFvApVA/GNE8aHlW10fLP6F20mTEAGm5crXhtsUcwe1wbgu2QO8WC4ZPE0qrWNTIY4ZAaxSzqmwMsmFWGkGoqcIS5aStBLQxbEiqKpOdJcOZKD4t4EBVAiQK21FG7O6Z3yzQS0B3WDCEtqsjW9TJEMUmDUXx1i7ij8EshmEEzH4GOc02c+k/G6srD+FNPAdwhHK4PZQ/51ADEx/863Z3nH95wP4/nMgBnliLQcoS3DrXyAGtfE1oAr2hiB1eAc5xCBNwIKARhYhjM74SFCzKeckOfMN4eZYgR2MQLUpGodYeOZO+L8EYYCs0awsMtr4Js2/MMKAxr1x8f7nIQze6i31LwQu4CsHnz+53L4ijo7jKz7cPnv+Yvvs+fW7F+fnmyfXNx9z/GD7ZvN0IBfO1i+uNn/9/wGZuO7A').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222517480', 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_1779222517480();\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
}
