{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "920dc66e",
   "metadata": {},
   "source": [
    "# rf208_convolution\n",
    "'ADDITION AND CONVOLUTION' RooFit tutorial macro #208\n",
    "One-dimensional numeric convolution\n",
    "(require ROOT to be compiled with --enable-fftw3)\n",
    "\n",
    "pdf = landau(t) (x) gauss(t)\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": "d47cc015",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:07.417409Z",
     "iopub.status.busy": "2026-05-19T20:30:07.417291Z",
     "iopub.status.idle": "2026-05-19T20:30:08.385254Z",
     "shell.execute_reply": "2026-05-19T20:30:08.384586Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "e131f095",
   "metadata": {},
   "source": [
    "Set up component pdfs\n",
    "---------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "ee50b928",
   "metadata": {},
   "source": [
    "Construct observable"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "e6b37953",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:08.387954Z",
     "iopub.status.busy": "2026-05-19T20:30:08.387820Z",
     "iopub.status.idle": "2026-05-19T20:30:08.548303Z",
     "shell.execute_reply": "2026-05-19T20:30:08.547572Z"
    }
   },
   "outputs": [],
   "source": [
    "t = ROOT.RooRealVar(\"t\", \"t\", -10, 30)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "05666d7b",
   "metadata": {},
   "source": [
    "Construct landau(t,ml,sl)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "d8f274e5",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:08.550429Z",
     "iopub.status.busy": "2026-05-19T20:30:08.550303Z",
     "iopub.status.idle": "2026-05-19T20:30:08.677024Z",
     "shell.execute_reply": "2026-05-19T20:30:08.676261Z"
    }
   },
   "outputs": [],
   "source": [
    "ml = ROOT.RooRealVar(\"ml\", \"mean landau\", 5.0, -20, 20)\n",
    "sl = ROOT.RooRealVar(\"sl\", \"sigma landau\", 1, 0.1, 10)\n",
    "landau = ROOT.RooLandau(\"lx\", \"lx\", t, ml, sl)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "de4517a6",
   "metadata": {},
   "source": [
    "Construct gauss(t,mg,sg)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "d2779cd6",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:08.678711Z",
     "iopub.status.busy": "2026-05-19T20:30:08.678559Z",
     "iopub.status.idle": "2026-05-19T20:30:08.793881Z",
     "shell.execute_reply": "2026-05-19T20:30:08.793113Z"
    }
   },
   "outputs": [],
   "source": [
    "mg = ROOT.RooRealVar(\"mg\", \"mg\", 0)\n",
    "sg = ROOT.RooRealVar(\"sg\", \"sg\", 2, 0.1, 10)\n",
    "gauss = ROOT.RooGaussian(\"gauss\", \"gauss\", t, mg, sg)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "b75b0b42",
   "metadata": {},
   "source": [
    "Construct convolution pdf\n",
    "---------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "6e9fff45",
   "metadata": {},
   "source": [
    "Set #bins to be used for FFT sampling to 10000"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "9e94ee90",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:08.795939Z",
     "iopub.status.busy": "2026-05-19T20:30:08.795811Z",
     "iopub.status.idle": "2026-05-19T20:30:08.904315Z",
     "shell.execute_reply": "2026-05-19T20:30:08.903638Z"
    }
   },
   "outputs": [],
   "source": [
    "t.setBins(10000, \"cache\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "94e46ad1",
   "metadata": {},
   "source": [
    "Construct landau (x) gauss"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "29ec0c78",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:08.906439Z",
     "iopub.status.busy": "2026-05-19T20:30:08.906216Z",
     "iopub.status.idle": "2026-05-19T20:30:09.022400Z",
     "shell.execute_reply": "2026-05-19T20:30:09.021719Z"
    }
   },
   "outputs": [],
   "source": [
    "lxg = ROOT.RooFFTConvPdf(\"lxg\", \"landau (X) gauss\", t, landau, gauss)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "c38ddbc1",
   "metadata": {},
   "source": [
    "Sample, fit and plot convoluted pdf\n",
    "----------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8141027b",
   "metadata": {},
   "source": [
    "Sample 1000 events in x from gxlx"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "dd7bc672",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:09.024526Z",
     "iopub.status.busy": "2026-05-19T20:30:09.024402Z",
     "iopub.status.idle": "2026-05-19T20:30:09.163345Z",
     "shell.execute_reply": "2026-05-19T20:30:09.162649Z"
    }
   },
   "outputs": [],
   "source": [
    "data = lxg.generate({t}, 10000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f951d02b",
   "metadata": {},
   "source": [
    "Fit gxlx to data"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "66e79f4d",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:09.165088Z",
     "iopub.status.busy": "2026-05-19T20:30:09.164961Z",
     "iopub.status.idle": "2026-05-19T20:30:09.512396Z",
     "shell.execute_reply": "2026-05-19T20:30:09.511740Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Eval -- RooRealVar::setRange(t) new range named 'refrange_fft_lxg' created with bounds [-10,30]\n",
      "[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x5558a2748920 with pdf lx_CONV_gauss_CACHE_Obs[t]_NORM_t for nset (t) with code 0\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(lxg_over_lxg_Int[t]) 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 34.3745 ms\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_lxg_over_lxg_Int[t]_lxgData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooFitResult object at 0x(nil)>"
      ]
     },
     "execution_count": 8,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "lxg.fitTo(data, PrintLevel=-1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "e54c1620",
   "metadata": {},
   "source": [
    "Plot data, pdf, landau (X) gauss pdf"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "d0327244",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:09.514391Z",
     "iopub.status.busy": "2026-05-19T20:30:09.514266Z",
     "iopub.status.idle": "2026-05-19T20:30:09.700812Z",
     "shell.execute_reply": "2026-05-19T20:30:09.700036Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x5558a29b6fe0 with pdf lx_CONV_gauss_CACHE_Obs[t]_NORM_t for nset (t) with code 0\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x5558a27f3bb0>"
      ]
     },
     "execution_count": 9,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "frame = t.frame(Title=\"landau (x) gauss convolution\")\n",
    "data.plotOn(frame)\n",
    "lxg.plotOn(frame)\n",
    "landau.plotOn(frame, LineStyle=\"--\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "91d52b08",
   "metadata": {},
   "source": [
    "Draw frame on canvas"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "b08ee9bd",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:09.702338Z",
     "iopub.status.busy": "2026-05-19T20:30:09.702214Z",
     "iopub.status.idle": "2026-05-19T20:30:09.914893Z",
     "shell.execute_reply": "2026-05-19T20:30:09.914113Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf208_convolution.png has been created\n"
     ]
    }
   ],
   "source": [
    "c = ROOT.TCanvas(\"rf208_convolution\", \"rf208_convolution\", 600, 600)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "frame.GetYaxis().SetTitleOffset(1.4)\n",
    "frame.Draw()\n",
    "\n",
    "c.SaveAs(\"rf208_convolution.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "cf06f2ff",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "fff862f9",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:09.916683Z",
     "iopub.status.busy": "2026-05-19T20:30:09.916537Z",
     "iopub.status.idle": "2026-05-19T20:30:10.097299Z",
     "shell.execute_reply": "2026-05-19T20:30:10.096588Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222610087\" style=\"width: 600px; height: 600px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222610087() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(31459,'WkwIty8A43oAeAHtnVmPHEeSoP8KkZiHGcAZ62Z+WsQTda16VwchqVvk9DaEJJlF5qpYxakqSlQP+r8vPvPIrGSROnqnezG9GEKVSvNwDzO/7HbPf998d/PTq93F9uVuM2+++XB78cP2+tvdk68vtq+uX1zebMLm7PcX+397vfvdR5s5hs3ZB/ub6/Htyyf/e/f0hvIN1b58dbO/vFiB/7m/eLaZU9icHd80//v7cP0cgpSyppzD5uyz/cXuw8vzy6vNLCv49c1P57tb8Nv9s5sXA/xkf36+VoZYwEPlGL357uzm8+3V8/3FZo6TlLA5+2r//MVpGZ28vLm5fHlSL2zOvrl89XbBo7M9VGjYnD2+/fpgfKXFo+ub7Q1ozKjzFvRgQDT+5Gr7cneXcMru9PxY7+0eHaseig+vpPnp0PDKDy6vnu2uvt7/eR2+k8LPL5/txsQ+ks18X/oU/V+S1k2YiceUW51yaSWapZ5a937qZk5lSv4vi/ZYfVB0M/faJ5FWklhRzZ1Rubl88OT64f7N7vz7zaxRJikxRo0ltcKKeXRz+RuebmaRSVo2M5PakvHux2+9u0ieGs9rsVhHF27f/QtPN/P9OBVLMTYtTXMtQod+/9bbWV+/v33dKbiZi9WwOfvDWw1K4yV/uG3ydsFmvj8KDqNzc/no+/dPhA/Prz9m6cVuuZdcTJtVYxZvX//4+/dO0Pp2nr5vtg+PIW2qvSaNMaYUi7LKH51f/vjFRx+OhfT4FHj0+1f+gBXz+OT7t8fST4/fHjy5futFD55cv/WuB0+ub5s9eHJ92/L3b16yte+Lz8dPDvjXNy+3bzZzWov53ppMtcaoGiX1xhZ/sbvZrpUevtiv3x5cv9o9vflqe7O/HN364vXLJ7ur8f2b/dPv39x+/Wl8/ezy+Vr42eXz27I/j6cPt88ebvcXsIawOfvw6vL6+sV2v77wCD68XPns6Z5lBQ34drt+fvlsf7bfPdvMZ9vz613YnP33q/2zN2+DP92CD55cf3h5eXVS/+Nn+5vtE5jZzdVrXvDJ/s3u2Vv9Prz64dX+5f5m/8Pu+h2G/tn+GnlxkCUruL262sx//FPYXL664ctfwubs4ze7p9eb+eL1+XnYnH0xhM/Vmcb+3dPLix8uz1+7KGFG9jeQ9f6HX7x++XB7vru5OcgHhvOL3Zubd0s/+t3XDz978Hgzb/7p8DVszj66fP3kfPfB67Ozw3R+tbvZ7i8Yy3UkHl3v/7z7/fXh+eO3QX/61W57DiNzHr//8yn87f7i2eWP31y+Yi+GzS38+BRemfRthU93iKR1efx44CYfvtisLOPD7c3NO8P/4OZmiG969uiD3c2Pu93FKp7egnxMP7m6fPnN5avNLBNr6tGz7Q3c34HHBwBu/GAA8pew+f7zyx92X77a/tvr40r5/qsdXX678OzT/fMXn9GFVRT7qt3ePH1xGNbvv35x+ePHP+wubr6+2d68vj4uzu8fvL65ZHkca36+u3j9wfZqwCyeB09ZfscWZ1/tts++vDj/6dDi7Nv9zYvL1zenK/Wwej/dXq9r71ByWuuPdzSUv5keBNv9WT3o290T5wb7i+c/pwyxMj48315fr1uFekP7Oi14xSLdxFlLCevfInMMMcRFvZRvaY7+NC7Zn2kpSzm2iUtda1Le1rp873O3oCKh58X43nNQqYuAMIf1bxGZxTRI1yBNF9G51eD/LZJmiRrWv0XyLNrC+rdImaWksP4tUmdpPax/i7RZYw7r3yJ9VrWw/i1is2YN69+i0SuL0T4uKm+DOoulICZBal80zdJrkA6YFs0Dby1BSly0zNITT4KUsihUrUTmuGibpZQgWYOkuGifJRXvYK+L2iytBI0VSpYUZ0kaHHVNS5JZUoOCIDUuSWdJMUjpgRenNPD0PtoyVDrwlLSkMou04O1zWVKdJWYocJJTmyVa8AZJl9Rn0ezjqtGWZLejXsqS4ywtBu+uyZJlZoqd5lSWrLOqBB8fQIaqDyJTXHKefYiZNc1LLrOPKWSILbnOvMZysLrkNtAwULktuc+OUZzMJdssuQXvhOhSohPhS2wpMkuuwdfsUtTHYQXSbP3wPc8tH76Xucjhe521Hb63OR2+dhb2+hqbY8gs4LLUOMdQ+/jO1jEZ33WmU76tlsoOkmTjCZuI1eHNC0BvA6jgkIGk+k5KKxbH7ptUlgp6B7QtDfwOxLw0CACQHpcGBQ6UurTDHmZ+GxT4k9iXBgUAvSzNCSglVFnaYSentrQDeklLs3lQGOLS4zwI5LvMxbvBd527DwPf2cKDu8Slsyp9GADG9gVfXDrbpB4Atu6gKi6djevDAMC29dEOcbGVcaUOIIMjSQYYzEsM9mXJn0itAHkAWQEGF2NRxcUGGzN/c/NaTSjv/j03vpt/1xTiIjHOIhJycaYgUXyFZg1ZF4k6mIIxvbJITHNvQXzV8Dh7P6TUkADLrNpDaqHkRWKdpcTQ2KC2SIRlVK8rlVexPTWkFFQBGZJQGVVZxHlrYQ/pIgIbk5AkKJD6qMKMqi0Cb+U7f1oWgbkymww875UyKB7rWeCuUYJZkBgXkTbD4OiP9EXgrqy7akHg6DZLNxDLIhpnyfCsHKTC3+FjY9w6oM7eLLbQbRFNsybHExeBtfaCWAhsc4G3SgrOypTK1be6Sh80w11PRIBGmU3YkiaLMiWwxPG3aIQ3pbD+LcqcCP3xv0VjGfwWtr4ok9KTs4C4aGwze7NJaLy5z5pskCVt0WgzrWDWZVGmhFVA40VF5iZhNF5U1Pny2npRJgV8NF5UMi1H40Wl0HI0XpT5gKNLo/Wi0ga1B6T92DnHak7vAS3iDnY0iF7UpwRqB16mZPRVelrUp8Q7qykvqvnYW8fMpKy9hWiFi9+OqrZDf0fjfuwwjdNhD8PbYl7SYRuzmJBb7ORVig94bGZ/HnVJ7OfDc4dPdAKHx66mfpElsbHX6g72owbhILxtqBtFlhxvFaO45Dh0ApoD6bEl0K2eApQRJ0iUJUdGY3Av+pd9Mzvv9/5kNjMw0rikJbN0oh77lyUOmOdRl4yyBIlIUIcVXGv3siSHYnb6xemgKr0RF3I0dujA6uOSxVn96JfA6mlCOYKGBnHJipwZXdJ1SobOuGTUobVLiLTMpl275F1mfaxdGs/hJOpdHs/RPtYuMQTaD10a1e3QJ6+dIiAIvHISB8esLTkdhJ43TS71xnwvObngPTZ0yXtsx3D4pC2FJbmuEd5STI5TPeDbqYeiYreTP57n0Z/Rx6WwJAceJ7lYheQxa0uxdoD8Xd730cGlmHcdKOpSI1PgUnGpcWjptKmRPvukLTXSY5+0pcahraMxsArHDnMSK+xr7LAV7t7lddqWikzh+Zi2pToL8xU4ULIKxwr0LlVEi++w8TpJ/jq6CIFIFl+Ea+UxHoduyan5UFetHa5Hy7E71s4hVViNYzjUFwKztlSW5KF/tINjrf1jlVRY1thhK7wK/nX6q/Os0T8fatfSfYet8BBvY4FUtUP3wI2WPnrn0Nifa810a0TVNAYFjaSinbPH/PsYjvHdlwZf24ERDrKXtjLCddiXxqpbSSyyNBbdCtK4HPvrD2+nn+61Vb1xEpa2Kjhg8sqjd+NhP/BA5w4dHriOa5Glx8EJfFKWvjJBZxI9HuecJys5g/f0dTWu9C79sBqd03TW4rr2HMdxsNExB0ceKNZFSGUqIkpXtkJN54TrK28ZIU9WdjS625GjB4YcdenI0bHsvXI/Qo7jyItQcn0Bgh5glQ5jdJaOCD3wdV4Lk2S1DWbc4ZErSONVE15JWlfjOgpLX1fjAdPYq0CQhM24IkLZXnnXIIrFuPaNbdExGdchYlv0NAZprJul3zJIKE4uL6Bp1F259iqoezrYB6PyLc92RG6iHJrmUxulZ+fYvMffm51jj8239OyW0rpO8oFh0898Sw7U5VtyxmuO5NCxfCRnVD6Qw4vKgRoH1hEaInrpmIvrvHhHynEZ+XsLFK3cuQ+rcR3QpZdB0TrgSy9Q5E9xhRSn6BaEopU/drchV3nQ6xBoY9qWPuzIMWlLd0NyZZx9GJKD2qXXId7HUlh6PRGy4Hdz8vatY3wOvXaL8th0mJSHZ211hRwIaqv0O5DUxsI+DEobC3tw7N7Gwj50tA1l6DgKDcvmYGH3Nhb2YSAwNMdDVlQfC3uIv97HtB2q9jFrQ2j07nM2pETvrgOt5GBx+toak4vNufJMH6M+tv4Bx9j6Rxy+9VcUtsoeb7Zy6sMIrIx6xb8qBy6m+sqkV2AwxcNsrjz6OJ8rjz685iiAfEGvPNpfZCuLXoFbDs242QmLZgUbTHqVxVEXW9n0wGIrn157Ygc+vU61rXx6DKatOsMBGhrDOtO2Mur14epd8zUFUZgbDPz4WwxmyToef4upz+CQoaZ5LgmzWhJNfZlnTOO2mNa5RLfAM0Z/m7NhqEvmpX3OPWB45rqY2pxrKDFI7oulOGfYZ5AiiyWZcw5OTF4s6ZxTcPbUFktpzhowvyuehDznGPC4VV0slTlZKIYvbLFU59RDHS5CS21OLWAiNzD0OZWAAdzAYLhNag3S6mI5zimF2oI0WyzLnBRrHtvYss7JTT7pZbGcZrWALdzbYjnP2t1JYLg1Cm6rVvAWLoYnrYbWghgYGh6CZkEMDBirobOSZbFss2roQ68z96MFN/3aYjBHC73gk1oM1thDr0ElLVbwHoTeg0pZrORZajC8V32x4gYhjjCVxeCLJZiuEFve3RPKO9HvcLkMyM1mqwOq6FXB2gqxZIJh5+OqwWeAL0LFFqsJSHBcOOjE+Kp1sIynaX1aBzgcF1bbAOlkX6x2uiWRXgLikA24JhxsuD4DrlKVtlgTB3GHOYjLJ+ASGqCPkIiuIF7igPcCF7jBGXsQ9zUA4l0Y/kIH2wDBy1OM+SB4ORhtOCMgZJTFenQQT4GDMkDI4Cm+a3d6DDANEDLyYjjkeEr3AXHIBlG6D1gHCBmAbYCQkRaDQVIZMgDNwQQZaTEbVOF7El0MI6rjfVzBMVb4m0QWg0vinYYMwDGDCTIAxwwmyMAH5wtKEmQA+opy76qDvqQks75tMUz6HAQvVsSXhtMOmAXfgXFGwTAOMF71IBl3hD9Ps+CBrkFjoz4eTDzSB7gMuB9g3DJBsuERp35zuECOw0gV2E7QWHju7lQp0AMscbZVFDsoDnoYg6c640wfvhKJkhyEFn+aHYSUxNMy49KBEgcrYIUQPI448PAnHsDuIGT4U5vx5jEquBg1elv3DgIKiGBnxpiqQiThD2PIXNV1d6sYQ4auqytMn9B28UHCnSDbnXgrDKGu7o4wh0ELLiMoh385MWP5O48EPQaYEbYJ0sGf0FhWGPwYYYSWPCoBjIdx+BQ79CR8sSsMPQnjcPhF8U/G5BEB95MOT26D8Qru2c5opA5bJnADB5cIewc/MPTA4ME/+L1EWDz4gcEPkwc/EgD8OSEshh+T8cgZUTJg6MkFQSPw8wb+XBFD0hlz8OeGkBqw4++IsBGtYjyyIeAGzHiUiPjz8FuFniIIR3eLVvAXRZAOGPyrTuwxLuaj5LnQHoEOfvTiFqRjOIC/1BkZ6bDja3OFHp47vj6zdhym/8VmQlz4fd1VXiMexAGDr8qMz5rnGXzVA4YDBl9NM2NDhC6Dr+aZ0KLD4KtlBpfD4KvVHd0O07/aZvfao0c4vj4T/AN2d3y12Yj+AYOvRaJKA2Z8EQsRguF3IEQwYOB4ARhdZ2aInCNSkGfEyCgAJ+KBTdPNw4QSERB4WSlQxwpjAS0FjpawAGgpcLTE7dij0YNtEhEUSmiRAmYWUcHS9QLQIiw84kgN1hbigsVI5BM/f0RgsBqpIaBFZLD8bgvwyRO1pAbjj9hggXkBhCE48koH0YKI6GCJeQ0IQ3ignXkBA3TwwXkBlLoTbh0PgVL3h6zj4ZzehQhjSsQC0l2MrGPqvBxB4iuPGlCKKIEteQGUIkwGnyHUSAEKwTq3FHgQyDlBDwKb9TAQW9sLiLocIkFe4BEa4g6HdxBKQaZ4IKIH+CnRIF/NwMRWPPRA13qAoXpAiIHojUCpeEDIA+gtwGCJCImx3lqAwRITEqNbLcBgCQsJcRxgjxLhxoGgAzw8VQ5DsIch6BIqHjCWC/RUj9gQGiJqSX0YMLEhQsIOQx/BIV/MnSCjEB0iyEJ/YLgEiNTXdg8d+vClwPEJ54If79661GG4HiKCw1sMMFyiRB5zBgY/3j04vklojBf+FN8IEmC4RIkUDg/s+IkwsqY1wHAF/x4cHpjxwani2yQFGK4gYAhtWQojkKaz+q7JAYYrLmBY8DnAcAUB45sIpR+4zJpY/5gEwHUmhukw+BEwvqVqgOEKAiaDrxKYFUkjd4E5heFKjrP6BmsBhiuE5OHgwOAnJu/7rQcYsGBH+HbrAQYsWBLZrYNQwI8t4bvPAgxYsCZcQXETUwR7wv0MMRTGA4vCYQ0wZMGmcDiFDH6sCodLgCFLwUZFtagBhiwFL6krYAGGLAXvLZqIhQw9hZgtiqN4aFUK9AGnAIMWQvYO5wCDlgJ9PK8BBi0F+oAbwWyRAn3AFhL0VOjDWokBfi0V+oA1wK6lQh9wDgl6KvQBlwCzlgp9wC0kxqOO8WMNEaSUCn0xKKox9FToA5YA55YKfcApqNMDfcCEzQm1Qh9wDbBtwSXjMNYQMPTxHE0cmPm1wJqBiUsb80uui0BPK+tz9iBwXWE0EmDsRdqjdQGz/oAJ7pJI9MHrM3LJ/skT9TafnF9ub5Juwubcs8dKCZsfNvMfTXMwxfSrwRTe0oOpBUsxGCHepMFSCpZysFSCpRoswaN6sGTBcgyWJVjWYDkFyzlYxngkRwTe1oNlC1ZisCLBigYrKVjJwYiBlxqswBN7sGLBagxWJVjVYDUFqySalGC1Bqvw0h6sWrAWgzUJ1jRYS8FaDtYwRWuw1oI1+LAF6zFYl2Bdg/UUrOdgHSO1BoOXOv+2YBaDEYM2DWYpGAkuVoIZJmwLZvB9j6vjh4ETR/eERCIQ7hmNxGJgwRG1IcJ8Iwk7EbYb4bURBhvRqyOsNcJPo1vFcNII+4xoyRHGGdETIiwyoglHmGOEI0bUgwgvjDDACNeLKLMRfhdhchHdNcLeoicuwcgiGkGEhUX4VkT9jHCsyJKLaJuRxRbhSpFlFkfyDi1YVRHHTITnxEIL98qhVEZkfMT9EtEpIxI+ok1GRHvE20LmOEofH7RAokecKxFRHlElIzI8YiBEhHfEhRLRISOiO6KpR2R2RAePJBBENMiISh0R1REFMiKjI6pjRPeMKIwRLTC6aHelAbEcUZYiAjmiN0ZEcTRaoDZGJHFEYYzMOV4LlAY+SHJgzkdWE3OOIwJ5zwcWP3PuiiJeBnENceRDMOeuHwpzPhIDmHPPl8I1gFTlwxMpwMGcY/cjHvkAB3PuCuHIyGLOXR0U5hxjHdnFBy2Yc4xw5BGODHAw5678YVAjWPigBXPuCVuu52EeIxb4oAVz7s5AYc49P84zXoQ59wQ3zFgYNR+0YM49w86t0JEZxpy7+uYZcq63uVHoCttIIWHOXV1z68z1NLe7pPU//eUvfwl/r2xOzoz8bDbnOITxC6dm1oxN0revLqh4SHEeJfeOL7hz+gPw5ITH3cMdt+dehDT8dw++fL69+n53dXKQZhScvHItOJ4N+Wb35ubBxXPyr0lABRwP4xQZA39+vn9+wVGMFT55P48/uSTL3U+DPNq+2b+bM/7g5uYB5SRuP9v/sL/eX15cb+bi6fs8OXnhZ9snu8PhHvA5PDBkMowd/vLs7Hrnh27gs2vhkezkdO+ffv/Z7uI5Z4biFEly9jk4NPW+kHd+t9nN+SHR/FjlgJ7E5Mf/MD306fy/6OG//sP08DhBf+UcfrC9Ojky9cH26rAoPAWbU2hs0/OHX4898dHV9sdxHGPAX766uT36MYD19McA1gMgX766+Wjk249zcyTBs4l8G3356mZlCXTiy1c3n/gxs7XqJ/v1XMA7GfRU8MJn+xsOyx3gby4vzz2BnoJxdOXDy4uby9dX1+uphQc3Kzl3OOaDmxs2sTOpX+AF+lcyA/YK3R/nk+gkECciJoc+vnj28dXV5XpKjZ3toFcH1SevL56ubIGHgCdcDHCdQp5yFGatTP8B18rse8CT+f5s93x38ez0vA3UjdITDsuLbgsPuA/nC3nFgS0cKzJ+61IMm7NPOQuxu77Dw9fSr19tn3IWwHEfz/yd9OF44G8tg8ZjvbepOVY9FB+q3kHt9e72+6Tw9rDRp/trFuQpPRTxvpWcGun2od4B8RidQ9VD6VrxDjXU+nx/sX/5+uW/7q4ub4968OCtY5jO4sepl4dXu7Pd1X//7Lb2KD8ZuFFw2k0oPS297eco/Wh39ulmLpFJO5Z8u5nb2yWPNsQBTqo8Xgsebk8X38PtW2sL5MeiW8xe9O7x04fbZ291nbF7uH125zzrKHzPidaH22cs9ke3w7OWPH6rBNm4nlUC4f7p9+tJpYfbV+Pc6KOVaRwLHm8IT27Ovn56tdtdfLJ96twHSmBrJ8MPyD44WbYUnc7HodXJ/qEK4O3uoeSwgMayGnWuXsJFN3XKfpCGQucZfjIUyFkMxukg7dvDeSeefQog3o+bq/2rj3ZP9y+359fHQ0XOkldtR48qw0nvvMKd7nnZaf9YJF540sEjfBRVQ0lBARnnnI+tjizQ24zu2LE2/Tm8n+4M5nP5u4uL3dVXdI+abDV/7fVm/iMHie7d40P0nmD038v3+gGu9/I9qUAaqZDXKreQF52Wa3y7+vFZj+uzHu/leC/fy17Tcf+dPv4Ec9ltn+2ukNZ+ZsqH7Qh9sr/55LBoyrpo/CQWs3h84mvJR+3p9twbM/v/43J/QeFBEfhw++oU/Gb/8qhOtt5NenaO8buX2+c7XnRk8B9uL56d7759sb/+fnf11fbi+Xowe5R/cPlmLRuzN0qdkpPjm3/YX57vLw6l68HFUfXD/dXT87vcfn3EuVOIPhGAj1C5P37z6tGp2nMofHxa+Ph9NQ+Fb9Wk4ufbNx/tn/vNASzCL69uXlx+uH25u9qu3OfvaK45TzlcZ3BXxXEZ93PmGifKGaAThgV44D3r05PhuyMY3zXCxsl6diZn6fm/n533An3fieRT5gHCAa/C4m8+aJ/8/osPGYxxe8MvDtyn8tHPDVvxo7CrnXvGrQPf3XxXSulbbWfpyRN2/tiKm3lzvr14tn1975/f/Mu959vX19f37hwAPhn/oVj8B2fgrzeDv3i6O0cQ4NLZ/JwZu9qwP3dVxjoYb9x+Oum969b/0IYvNsWTPTY7zq/N2aOT4/ePbo/cPxp1OAd+9sn+6no1Yj7bHr5xrYgeBNjL3Uf761fn25ND63DUI09mCN3avz1E/vnls8+2T1b4F6zx3zZRP92dKD+rfH3vf/23e/98L0753r9Awn+eiZs4YezS7ahErK6OX/BZnEzdceLGFHI7gpf9P5+1n/Ew/LZZ+/PdWfvPNUn/n0zRHafIicckOgf4+OLmao96yQL75vr1yx/Xbc3Xwxbn++HaDP++Pvh8+waz732icLUIudtjOEq+uLx6eTA3WLerV2PcNnH29cDmHOetq4lwEHC9xLvOyN9+gcW4MuLWJDlcIeHYPthfYFF8fHX1JVdeQBrwlz/srs7OL39cR+bB1RX87e3AWyUi54E3lzZcz+DhtztXlfztvduv2Cq/Lva/urzEBqfy+0SdB5WOF3q8+O78zfOPtjdbqjtz2swbml8+v9q+vHd5dm99/t2r88ub777z157IezbMf0DaM+6/LO3R7+46vb94dbm/4IarVZz90gzFERq9b1MP923K4b6F+32qfGi43yhuFLdwv1JcKS4UF4pLuJ8pzhQnihPFKdxXipVioRj+fl+CX0r01j+K5K0SC3HiQp7bfxLiVG/B6DHJwCt5N0jAlgK4IQJqIKsEiIRayIb+FugN3aJ/dNQC3ab/wtUhEifCNkH87bxeeL+AgGgmKAQcAhIBC6FN8AiIBEwCKuKcIBOwCegEfAQ9wSigVFAqKD09diL9lGLvEijJiwWlglJBqaAkPxaUCkoFpYKSxFhQKigVlApKMmJBqTZ1v7zm8a8sCTatjZWRggaCXlBCqrl6JrpHBUfOIVkTZHMkDnZXDclayJpDLj0UGhULuVvIRnZECbloSK2ElFtIpKmSpU329yHfjMgaYUuCsR7qJdLWSN/uLbRIFnnL5LxXJbc++5fECznfzX8MV+CIPifSGVlikGTnHv7TwIz4DDGRAdQ5kE3gXxtfCZzSc7IIg0Qft78F3z1yj+PlQe9IhaMEOYiJt5j/x1wr9SsT+MPqFyKuFSO65MePXuyfv/jrm3Fn1a+0OlktMtWkLebaKyHmFmRSky651CidZaNT6ZJbpKRWMlQmraU0yzFa8XP4U6uxSa1KYDyxda2z5VurNXP8YJLeeiot5phaZUsnS6VZKVote1LAJEmTJS3RskfXp55bapKyWScyniZpXbVmWnkiyZQycf+soPctXVtqajWaZbIZ2pSSNbEquXiw3aYorUitPXdPjrDJpOaWe5Pql0zEqfQkuWuNiWC+CqRKigSLK9d06NSrdIvdmpEQpDJlqblbqq0VblHQiVwDE4t+U5kXWI1ditVMCk6cUk0pW2w5NU/vn7RoTiXnmKsnlk2pVO3WcqqZiHmfQDoIJ92Eglx6lN60dlIF2hTFSlPT3hpB8To1TSrarUthfPKUm/XULXlQP03VjOB9k6qEw3UqVrWkbrlVJiVNrYlJ7ZEp4JaHSYq0pjV1a2RSyJR6i63WElvxbIkpa9YWe0/WPVth6pHoeG2lkJBhU1HpKcdeavEkmilFy2LdSozZLPQp5VpqTU0sk1nSp1ZKr11aB3kLfSpde06l5dwLFwtwcR73A0bS0jsrrhdtotzBl0u2UKeeCOEn1mhtyBaLiSlMtTdWyztNytR7qy23FksvWSiIjGeOlkruGsrE9YG1JC2dGX9fwZ0md/dKnnrNVVRj6rloe0+FWnKJmqzknK2GPBUmvppZ0gQPnay0Hsm1yTF1QZhG7b2psOwi77zT5LcU3EGrU6s99RpLzVUb4vodLLEztd0sZQHrnQosuFpiTNK7MuRpIgfGck6WKnkY79TQyXLKlfwgLWTYvFPjnXfkO2S85xXvYL1T8B6WxwKsmlLMLgqnrpZ6jNa5lrK8Q8W4ce7xb+PgMkluveVcmeGkwaUA4lw28x91gh/30pP2jqScamqdnZpyKWQQTq1F7bFpzCmT8ghDKhqLJqk1JhZuLZpytVI7CWptSuzi2EpvVktFtyot5dKqWuyJTZgrWTFVU1Vmk4svuzT4VhlJQFO3lE17YbqHeiVVrGuVlOEFZUqtWM5dNBaSlurUYqks1pIKLKpPqfbSY0oaSeDVOMWcOgyhlsK5sThZlgoHgRbntjW2XGMpauQQq47eSi2dTeiKV08N/muNNC5qJOvFGRtpypomdrEaWUUwB02TSq9aq0ni2k24foc39FI4y6aTtmKdddt74xDLlLpKToibkenr7DjnzBIn4cum2rTm2GMtieynPsXcrdWSpZHYJ23iitLU2Myw+DIVgVBrzbPF8tRSqVlytpZI40pTTcVab9IKaY+Spx67Fi3NOgOeJjG1bAJV5FTRczL5pAsMxHXkplpFDDnoDLznlq2JCqjRpgtCB7lckmfvIaycsBYLCZRTjjXXVnpKUslknLqotV6j1s4Ru6lmEtKylVZLgX/XYvTEulTM6I7UTZVFEJsa/Nti1m4afd0WChBBZhqtU+OdJnWypjky2immWuDwvUuMJcXUiwjWA5qBaOqt54oMeKfgTpN3NkqZrFZLGmsjaQ0Of2crlalVSbXmmElwQyjUwoBbE+mkspYp9hRrKTVHqcp2RArFllRTHy99p8mvF9xBmybLRbIUsV4ttvdh0YhMTd16T/09dOSJseyxlVJ7F2QNeyJbq2SLJljLOzVIO4xZUqu11fK+GnfeUe7Skd99x50m79DxHp5XrZXGhqyR7ufJSm+p56RNyfC8+w5Y9Hc/4eVzr+5mfq+/9bvjba7f/eQuypJZaEf+0sLmu4vLq5d/2J5zQWV08OX+Ynv+wf5ijdJEPKXfXXy9f/5y686V73YHdxXuK9pcbX88urDuS9h8d3m1f85bvvVbOd1tdJosfeqzieFUof8v8+9g/v3NI1Xnv9lf9eHrqx9+NsJXtKTTaNX5m+ff4Vb8482fTp1WD68uuep8f3nhXqs1ZPVoDVlR88RtVd3z9POOK8LHvxgo/HXHFZ6xn3dcqbvuf8HupAI5/fclTpmcdP/CJ26syofiw+p85HDf3VgVX5biw+p85HC/4sbyD8WH1fnI4T4eLIoTxUqxUuy+rIovS/Fa+b3ct/8oSrcg/97n2ypTKaRGiybR0vruvtT3OrzyHYfXuy4wjKHTfxJwiunqGKvDORYgnT7QmeEgQ+/tIcPD0O7cVUbnES5IHFxmFdkU6uo4y+48Y9wYwDbVgAMNOYrAVoT08Jjg2Ir4u5wOodhNwYz2gEMH8hhbXC3gRutAWeEqGor5j2JwCcgEbOgW4BMQChjRiMbld3QPlApKBSVOH+81KNHG+MBLhfPMP7ixkWJQcrgclApKBaWCksMyoFRQKigVlCmO//A0jqX3qx411urRDf5HmaJpTZatR9T98E5BnFrFHss1WTROdcWpSOcG9xwVO5OSlBr+kqxqfuIiTioaxWKUUrLlECfRZq1rq61q477IKbasvav2pC33xHtiltxSllqRdrFTpJpSRcxZreRGktUoXbXhuai5cT/lFKNbs9J6ixo51EJZKd16yr3lWA9lvXQUtRJr6eQu8bqeNTWsjtIwdZ0ORK6k2JPiPIxTtNRqLk2KaW6ZLqHY5xxbiTUl4VVZLGNTpCSJ8HqcesI5kEoySwXn3VSjYTZUzZFDaDpZa6mqgQ3drkwofYZdiFItwSZUP1xSvWAosThjbpwgqCVnFliarKfYMICaX0JkU6/0qjTcHD2kijGekklrramFXOhRxvzCtNVQ8lQhIUrWVEVDLVOLgkadW+QUR+sTWnDOWlPiQn5DB9cae8Je9NMnfVLNrRYuyee8CIs69dakoA/7HZ9oqtY1p+anb2qZOJ/D6OBY4eDkhGqMWWMNGyelqUuOveWipbAZEkNkWrP6mZCUdOriS7Ja7S2kUqfaFY9Ux4AMqdmUsRNRllNqIUf4lWV8JMJCy273MHNRI26WnNrU+SUErSVW05BLmmLDbxc7PsWQGaFUUSF75ExVbmVqYqULhr/WkLuSH5eLNQwl/MNtSpivVcV/1SH3Ppllw92VaNDwChhegc7BktzzlLT1imuG8xu59YlfkCgRNwM3hbY4iRiOsdozV7wUQ7PPJQpLrYacO14AruWPBWMqJ6yr0qpZ6T4YeTLJrbiRk0NqOqWYO346LU1CSn1qucfUsKL9JF+ZKsszacN3F7TlSRKGZ/NfsgiaXZ9Wy6V3DtCo9In9lGOvkji4ZXmShvzJJXcOxzSdDFMehpTd65knFZPeo41jT8wJy1M148jkpPZUcJT1rI2jcCIyYfr3ItjdhFCI6DRcqfj5crA4VWO3NYuKB7frxFKkr9bZGq1MsUrPRvdwg9Y+Vc34RyR3HH2Vac0dl2ROWVIobeod6zSLdK7SSUjfnpk3fksjZJuSiUrHLBIWE+KSzYUjFRVep9JSxCOnCcsy2VQrTineW1sKqU0xpYKfs+dCkCNPNTOYJUWID0mn3KP2EjH+LSEUCl5pDHMcpR6BweLBoZswz7ThlVC4aFNOBWmZiuaMvxdLnIONsG4pJVW8xO4DxhGWBYeXX4YCGVZrylnxeuJozq3WJC2acogMx2gtPVuJ3ZiUPpXWtYjilMV4b/jPoka8p4mgS506v0EC2eonourEUBUc0yWtLphupgXPXebkWp5Yizm22KRxpipP0WLEyZzUT5inKZu1blGYBHci4MPBSVlYEO5VUI29mLLw0BxM8a7GWmOmgky5xpi7lFyTnxqboiTN1lqN2c+/4cElgtA4n+iRPKlMilXctvDzhjmL0CgsL5uSxqS9dKjAhdB5marUlLTgE85wjKw9m+EHgjOl3iVVohC14vRKtZYmmW3ff0PBf5oQ0tt26q9aBu8YrXc46Nsm6z+OnfQbzKQ3v2YYkSCCnbTa5n93wyiRJfVLhhEV/gEMo5+xgt61ed6xcND9sQgwBCYL2DeypgAQXcTWKZ4OoBOaHukB5qkBB+sHFzcWEKFHrCJzS8jjHCFPabWIcD3hscL1Q0KBuXWEqnhrIWEdYRl1TyvALsL3CCfo/2Uf/Xb7iCX7ahzJYOV6iLTjjq85t1h39xFoqMrJcsz4TnvT3X2uJNFJe1NMgIazPO/uc46YII5rNKXUZLq7j9N9Smj01XLPWVvb3eck99RUzLpZr5V3cpAY8wfG362gr+3uk8QwZRX84i6RW93dJ7wwWTfF+EpiVfvuvrH4XEdGn7Tcd/eHtRMt1ZylV0PNj1NNnDi3jM5ZK+sU2UyYONeS0GYmReiZtl5SUSxsQiMxNzzhCCQc/KpdMDqyX5iWp4oF1EuyjsLY6HJuhAEQuJgbnLDueDJxZKJ1FandJJmO885REJTac7OWObiedGqtJyWQ4xK/dv+RNBzU2SPasU1EfVIRV11KnDA3XR8nbSNlnVKOnI9X6dyYkcqUo/YaUUlrDrm78YcuVGJtNRRlsloyDv1jfpYqU1FDdZGKyVyMHACt5jHainbY3XLI6Hsejy11ysSsWieyUkNtKHc6RoRoUTVUTLFS0Kn9Onk0yJ7EYm8cAW/oYbXkpKUkzo63xA+95Uq0KXHau6EPps49AtJL5roWV+tTrLVZ5zdjkk5qhPtaVE7PN9Q/DxRHaUoNPBmJ0AJhLNXQyEOKseXYsZWg1CZIEo5hcLFLbZ41kGOK2YgG1awoiIJVG9VqKK1Nmghqo9pJD0XiZFGTSI+55OIWSqkkSWRphFchQ8lFaDknQqOmUyrFpBqHP0pItU6N4dGo0QvAisodUa64nIRMJmzXXlP2fKCISwDD2xImd9Captgzil5Lys+MJJl8eRK/YK415gmDsGDC+kUQXTDGcWD06gfmK7OA2pmNECOX5U3WsqIv1sydAWqTe0Z6j92we8jKqo2tVvACcGfDFFu0asbW5G6OicyG6OZ+ai30OjVRdONCVCU0m6SmEnM0aTnm0HTKllWJcJbkV+mQ2sIabTW7ySJoma1WKUZCSMBEVtcn4WEUEG207tsvFpOQycPQRPytiYfny6R4gVSKNreAdZJKfCQRoUyG0ZJw8lguhDlJwiII10h3IbqYMVoKfgiNZp2bGTBakhEGLEqCixstUjz1oPr1GEqSBFdXwBz5XUJtk0flBM6Q2nCdpRYJwTGoxZPTauTuiYIfwsOc1htRP/rDBYaC/Z9ySqkSwcJqIc+npeGB8sS4hoPJNLbsN0xNrRD3JA5VuLcCZtM8sqqlYXy2KZJJ1Fvyhc9Ed6w7rb0kAs1lKi16CK1jLbnD0Rr2qJLtA6y1ZjQPwlnuo6wtlRhJpUjcQ5EmLn2w5r+zR6CZ1JZonlDDZFMQSQuqsem4jGJiyRO0M7//YpLes8ZIfM5v3iAVpljJDROESzkIoQvB3Ob3X7BvyArByQVzJwpqhGs15ZS5hGZKwg7i7gofK3dqeH6TeUYRcVKLPSZsHA8zY9MojqBE8k8mTC9KiB47CR/Trxf8Y5gwaO3vGC3vsOa/s9ny4NHvvv76wecf/3yU55+udmebWcrf3GJ6cs6vUP6mbOhvHm5/8COfVH9fPjQ3cZxEhv4jCc0Ehfz0mE2tkqiAb7d3LBqOkvWYp0rKX1Ft6gdUOVk25OXguo0bfzZn41dZ06SWK9kh+D/93f7jm/wkLm6RglvRY+i8f32E/7XGnhP73+8qeKTro54LKYXVSH8o/lOsh0dGMllEmaod5nc4xHabPf+7Cy41YHS+frF9dvnjyXm7Dy+vOL26fbbnRwIZhGM653Gi1sNVN5x1YSY4vvn+e0sOD8elJeoHFN4+kDyOfnJpye2VIuvh0c8uL57vOK6ka974+hvC45qR/cVf+YuYd+5a+Gx7s3OD+X3rCMLXXnoHvae/fnbutw7EYVxOTzMfL2/h4XoifD0K/ptX9LDvX+5vximNMeFY+V96RP14VKlIc3ZMPkPR/pfjj4XyA6Fc5cOpsa+fXu1fEX33IeAXJs/HL0xeXOye3tweYvefLh2v9t94/Mv/AbHJkv4=').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222610087', 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_1779222610087();\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
}
