{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "cab3f50e",
   "metadata": {},
   "source": [
    "# rf308_normintegration2d\n",
    "Multidimensional models: normalization and  integration of pdfs, construction of\n",
    "cumulative distribution functions from pdfs in two dimensions\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": "0992c586",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:54.949582Z",
     "iopub.status.busy": "2026-05-19T20:30:54.949467Z",
     "iopub.status.idle": "2026-05-19T20:30:55.902581Z",
     "shell.execute_reply": "2026-05-19T20:30:55.901849Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "73d599db",
   "metadata": {},
   "source": [
    "Set up model\n",
    "---------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a2f8fb3b",
   "metadata": {},
   "source": [
    "Create observables x,y"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "401c3979",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:55.904363Z",
     "iopub.status.busy": "2026-05-19T20:30:55.904231Z",
     "iopub.status.idle": "2026-05-19T20:30:56.061671Z",
     "shell.execute_reply": "2026-05-19T20:30:56.060992Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", -10, 10)\n",
    "y = ROOT.RooRealVar(\"y\", \"y\", -10, 10)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "75e3b63e",
   "metadata": {},
   "source": [
    "Create pdf gaussx(x,-2,3), gaussy(y,2,2)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "0d0de452",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.074490Z",
     "iopub.status.busy": "2026-05-19T20:30:56.074351Z",
     "iopub.status.idle": "2026-05-19T20:30:56.201127Z",
     "shell.execute_reply": "2026-05-19T20:30:56.200448Z"
    }
   },
   "outputs": [],
   "source": [
    "gx = ROOT.RooGaussian(\"gx\", \"gx\", x, -2.0, 3.0)\n",
    "gy = ROOT.RooGaussian(\"gy\", \"gy\", y, +2.0, 2.0)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "53482bb7",
   "metadata": {},
   "source": [
    "gxy = gx(x)*gy(y)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "b7f83fc4",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.203219Z",
     "iopub.status.busy": "2026-05-19T20:30:56.203095Z",
     "iopub.status.idle": "2026-05-19T20:30:56.381330Z",
     "shell.execute_reply": "2026-05-19T20:30:56.380643Z"
    }
   },
   "outputs": [],
   "source": [
    "gxy = ROOT.RooProdPdf(\"gxy\", \"gxy\", [gx, gy])"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "ed4ed652",
   "metadata": {},
   "source": [
    "Retrieve raw & normalized values of RooFit pdfs\n",
    "--------------------------------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a99c332c",
   "metadata": {},
   "source": [
    "Return 'raw' unnormalized value of gx"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "a722a223",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.383522Z",
     "iopub.status.busy": "2026-05-19T20:30:56.383395Z",
     "iopub.status.idle": "2026-05-19T20:30:56.506925Z",
     "shell.execute_reply": "2026-05-19T20:30:56.506261Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "gxy =  0.4856717852477124\n"
     ]
    }
   ],
   "source": [
    "print(\"gxy = \", gxy.getVal())"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "94414777",
   "metadata": {},
   "source": [
    "Return value of gxy normalized over x _and_ y in range [-10,10]"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "895eb5ff",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.508655Z",
     "iopub.status.busy": "2026-05-19T20:30:56.508514Z",
     "iopub.status.idle": "2026-05-19T20:30:56.630125Z",
     "shell.execute_reply": "2026-05-19T20:30:56.629424Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "gx_Norm[x,y] =  0.012933200957206766\n"
     ]
    }
   ],
   "source": [
    "nset_xy = ROOT.RooArgSet(x, y)\n",
    "print(\"gx_Norm[x,y] = \", gxy.getVal(nset_xy))"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "ec726200",
   "metadata": {},
   "source": [
    "Create object representing integral over gx\n",
    "which is used to calculate  gx_Norm[x,y] == gx / gx_Int[x,y]"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "20cf3b2b",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.632035Z",
     "iopub.status.busy": "2026-05-19T20:30:56.631913Z",
     "iopub.status.idle": "2026-05-19T20:30:56.752558Z",
     "shell.execute_reply": "2026-05-19T20:30:56.751909Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "gx_Int[x,y] =  37.552326516436096\n"
     ]
    }
   ],
   "source": [
    "x_and_y = {x, y}\n",
    "igxy = gxy.createIntegral(x_and_y)\n",
    "print(\"gx_Int[x,y] = \", igxy.getVal())"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "4a81abcd",
   "metadata": {},
   "source": [
    "NB: it is also possible to do the following"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f96baaec",
   "metadata": {},
   "source": [
    "Return value of gxy normalized over x in range [-10,10] (i.e. treating y\n",
    "as parameter)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "ce16a86e",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.765480Z",
     "iopub.status.busy": "2026-05-19T20:30:56.765331Z",
     "iopub.status.idle": "2026-05-19T20:30:56.869043Z",
     "shell.execute_reply": "2026-05-19T20:30:56.868358Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "gx_Norm[x] =  0.1068955044839622\n"
     ]
    }
   ],
   "source": [
    "nset_x = ROOT.RooArgSet(x)\n",
    "print(\"gx_Norm[x] = \", gxy.getVal(nset_x))"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f66e75d3",
   "metadata": {},
   "source": [
    "Return value of gxy normalized over y in range [-10,10] (i.e. treating x\n",
    "as parameter)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "ee889ba3",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.879388Z",
     "iopub.status.busy": "2026-05-19T20:30:56.879255Z",
     "iopub.status.idle": "2026-05-19T20:30:56.985171Z",
     "shell.execute_reply": "2026-05-19T20:30:56.984591Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "gx_Norm[y] =  0.12098919425696865\n"
     ]
    }
   ],
   "source": [
    "nset_y = ROOT.RooArgSet(y)\n",
    "print(\"gx_Norm[y] = \", gxy.getVal(nset_y))"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "92572583",
   "metadata": {},
   "source": [
    "Integrate normalized pdf over subrange\n",
    "----------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "502ba548",
   "metadata": {},
   "source": [
    "Define a range named \"signal\" in x from -5,5"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "57d36bb3",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:56.987156Z",
     "iopub.status.busy": "2026-05-19T20:30:56.987035Z",
     "iopub.status.idle": "2026-05-19T20:30:57.096408Z",
     "shell.execute_reply": "2026-05-19T20:30:57.095695Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'signal' created with bounds [-5,5]\n",
      "[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'signal' created with bounds [-3,3]\n"
     ]
    }
   ],
   "source": [
    "x.setRange(\"signal\", -5, 5)\n",
    "y.setRange(\"signal\", -3, 3)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "82e62aeb",
   "metadata": {},
   "source": [
    "Create an integral of gxy_Norm[x,y] over x and y in range \"signal\"\n",
    "ROOT.This is the fraction of of pdf gxy_Norm[x,y] which is in the\n",
    "range named \"signal\""
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "825747ed",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:57.098140Z",
     "iopub.status.busy": "2026-05-19T20:30:57.098018Z",
     "iopub.status.idle": "2026-05-19T20:30:57.236131Z",
     "shell.execute_reply": "2026-05-19T20:30:57.235442Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "gx_Int[x,y|signal]_Norm[x,y] =  0.5720351351990985\n"
     ]
    }
   ],
   "source": [
    "igxy_sig = gxy.createIntegral(x_and_y, NormSet=x_and_y, Range=\"signal\")\n",
    "print(\"gx_Int[x,y|signal]_Norm[x,y] = \", igxy_sig.getVal())"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2a63188d",
   "metadata": {},
   "source": [
    "Construct cumulative distribution function from pdf\n",
    "-----------------------------------------------------------------------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "73d26e64",
   "metadata": {},
   "source": [
    "Create the cumulative distribution function of gx\n",
    "i.e. calculate Int[-10,x] gx(x') dx'"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "d21a6470",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:57.237974Z",
     "iopub.status.busy": "2026-05-19T20:30:57.237854Z",
     "iopub.status.idle": "2026-05-19T20:30:57.348159Z",
     "shell.execute_reply": "2026-05-19T20:30:57.347474Z"
    }
   },
   "outputs": [],
   "source": [
    "gxy_cdf = gxy.createCdf({x, y})"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "16fdcf42",
   "metadata": {},
   "source": [
    "Plot cdf of gx versus x"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 13,
   "id": "3d512374",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:57.349955Z",
     "iopub.status.busy": "2026-05-19T20:30:57.349832Z",
     "iopub.status.idle": "2026-05-19T20:30:57.607146Z",
     "shell.execute_reply": "2026-05-19T20:30:57.606524Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf308_normintegration2d.png has been created\n"
     ]
    }
   ],
   "source": [
    "hh_cdf = gxy_cdf.createHistogram(\"hh_cdf\", x, Binning=40, YVar=dict(var=y, Binning=40))\n",
    "hh_cdf.SetLineColor(\"kBlue\")\n",
    "\n",
    "c = ROOT.TCanvas(\"rf308_normintegration2d\", \"rf308_normintegration2d\", 600, 600)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "hh_cdf.GetZaxis().SetTitleOffset(1.8)\n",
    "hh_cdf.Draw(\"surf\")\n",
    "\n",
    "c.SaveAs(\"rf308_normintegration2d.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "afa1e1bc",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 14,
   "id": "052c994f",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:30:57.609295Z",
     "iopub.status.busy": "2026-05-19T20:30:57.609172Z",
     "iopub.status.idle": "2026-05-19T20:30:57.787677Z",
     "shell.execute_reply": "2026-05-19T20:30:57.787124Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222657777\" style=\"width: 600px; height: 600px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222657777() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(40549,'WkwI8j4AZZ4AeAHtveuOHUeSoPkqRKJ/7ABRsW7mVzuB+UFJpVHtSiWipGpR09sopMiTZK7ITE5mUqK6Ue8++Mw8TiYpqfqy1bvV2AGYhLmHR/jF7ubmfv757E93P705Xp2/Pp4dzr7++Pzqh/Pbb47ffXV1/ub25fXd2XJ28cery//x9vi7T84OaTm7+Ojy7jagL7/7v4/P7qg/o9mXb+4ur69m4f+8vHp+dsjL2cXpS4d//qW+fq2DnIvmUpazi88vr44fX7+6vjk7yCx+dffTq+N98ZvL53cvo/jp5atXszGDpbg3TslfP17cfXF+8+Ly6uyQVqnL2cUfLl+8fFjHJK/v7q5fP2i3nF18ff3m/YqnF5eMQpezi2/vwccB8sbT27vzO7oxo817pcdR4uVPb85fHz8cOHUfzPzU7v0ZnZru1fsnef3h0vDJj65vnh9vvrr8p7l8Dyq/uH5+DMQ+lbPDb2SVVKRUTcnayD6DWS/dVKzXPlKi/qnGJIeNUVOzWnypv9Wzg6w/b313/fi72yeX746vvj87ZNHVUkpJVMXR8fTu+l9+eHbQkdc8WiuarVlj1t++92UdLer2r71fcXb4jda6jiwj9d4kD4jtj+99AhL64/1wHhbPDtX4/t+/90LtIPTv7195v+Ls8Juo2Bfg7vrp97+81j7of/kxC59Srjp66SNrGqX25exBB99+/4tYmN/nKaj+EKX747PDb7wDS2PISNYlZ6b49NX1j7//5OOgl28fFp7+8Y0/cIQ8gL851X52gh5/d/vehx5/d/vetx5/d3v/2uPvbu/f/OO713Dwb9Lax2i9Si7VyhAkzh9/2p8NKzn32rpaLwl0/fHd6/N3LNovvhaPfv7W1y+Pd+dnhwwBPHl5OaHHt2+Oz+7+cH53eR3r8Pu3r7873gT89eWz7+kIwXH57PufAvz8+sWs/Pz6xX3dP8XTJ+fPn5xfXiEylrOLj2+ub29fnl/OD56KT66n/H3Iy+Akyvds/MX188uLy+Pzs8PF+avb43J28d9uLp+/e7/4033x8Xe3H19f3zxo/9vnl3fn3yHk7m7e8oFPL98dn7837/3TT24uX1/eXf5wvP2ZoP/88hY9suuYWTy/uTk7/MM/LmfXb+4A/rycXfz23fHZ7dnh6u2rV8vZxe9DKd1c5DT+dHV98/ry6u744ob1vtLnKI6vL+8Y3F9q8vu3r5+cvzre3e06hKX9/fHd3c9rP/ndV08+f/zt2eHs73ZwObv45Prtd6+OH729uNhR+4fj3fnlFes6V+Xp7eU/Hf94uz//9v2iP/3D8RyZRef++EH5m8ur59c/fn39BmZezu7L3z4sT0F+3+CzI2prksqPuzj6+OXZlDkfn9/d/QwVj+/uQsUzs6cfHe9+PB6vpgp7r+Qr++nN9euvr98gP5znn5/foSG88O1egOMeR0H+vJx9/8X1D8cv35z/j7cnqvn+D0dW5P3Ki88uX7z8nClMde0UfH737OW+rN9/9fL6x9/+cLy6++ru/O7t7YlQv3/89u4aUjm1/OJ49faj85soQ0iPn0GKpzcu/nA8f/7l1auf9jcuvrm8e3n99u4h1e6U/Nn57aTDveZhq3/4wIr5q9lKaJ9ftZW+OX7nkuHy6sWvGUxQxsevzm9vJ9vQLiy0hxVvINKzdNBal/m3ySEtaUmbei1QPiR/mrbiz7TWrZ7eSVubLanvsy3wOAxbVGQZZTPgURaVtgkdlmX+bSIHMV1k6CJdN9FDb4v/2yQfJOky/zYpB9G+zL9N6kFqXubfJu0gfSzzb5N+0FSW+bfJOKjaMv82sYMWXebfpskbi/F+2lTeL+pBLC9iskgbm+aDjLbIoJg3LdFvq4vUtGk9yMg8WaTWTRnVHGRJm/aD1LpI0UVy2nQcJFef4Gib2kF6XTQ1RrLldJCsi3fd8pblILkzgkVa2rIeJKdF6lj4cM7RzxjxLkul0U/NW64Hkb74+6VuuR0kFUbgQ879IMkWfyHrlsdBtPi6arIt2/2q17qVdJCeFp+uyVbkAIp9zLluRQ+qsvj6UGSpRgwyp62Ugy8xWNOylXrwNWUYYltpBz5jZbG2lR7dsFClb2UcvEfxYW7FDlL64pMQ3WryQTiJbVUOUtriNLtV9XWYhXywscPl0MsO10OVHW4H7TvcD3kHB4Q9P2OHtBQIuG4tHdLSRsCwjknAemBSzlZbg4MkWzyBiaAOf71SGD0KjT4kOmnOSXn24r07k8rW6N4L2rdO/15IZesMgIKMtHVG4IXatr7zMPjtjMCfpLF1RkBh1K37AGpdmmx95+Tct753L3nrdogRLmkb6RADBJZD9WkA62H4MgDDwiFd0jagSl8GCsG+9Je2AZu0vQDrxqjSNmBcXwYKsK2v9pI2m4IrDwoSEkkKhRBeYogvy/5EWqNQolCUQkgxiCptFmLM/MvdW3WhfjhcOrA5rHlJm6R0EJGlVBcKksQptOhSdJOkIRQM9MomKR9GX8SphsfF5yG1LZliPaiOJfellk1SO0hNS4dBbZOEyGjeVhqfgj11yXlRpciSLI1VlU1ctlZ4SDcRxJgsWRalpL6qCKNmmyBbgfnTugnCFWyy8HxXaow46FmQrkkWs0VS2kT6AQHHfGRsgnSF7potgkS3gwyjY9lE00EKMqss0pDvyLFYt0FRD/5a6suwTTQfNHs/aRNE66iohQU2F2Sr5MVFmdK4OaurjBgz0vWBCtAkBxNY0mRTUIJIjL9NE7IpL/NvU3AizMf/Nk015C1ifVOQMrKLgLRp6gd4s8vS+fI4aLYYlvRNkx14C2FdNwUlUAEvbypy6LLEy5uKulyeb28KUuiPlzeVwpvx8qZSeTNe3hR8INGl8/am0mO0e6fjNDnv1Xy8e7eoO8RRDHpTRwmjjX5BScxVRt7UUeKT1Vw21XKarfcMUuZsGbQixe9XVfs+33h5nCbMy3nnYWRbKlve2RhiQm/ByVOLRzmY2Z8n3TL8vD/38gObwMvB1bSvsmUYezb34jhZEF5EtoW5UWUr6d4wSltJYRPwOiU9vUnp3k6hVFAnaJStJFYjpBfzK87MLvt9PgVmpow2rnkrkE7S0/yKpCjzPOlWMJYYIhrUy0pfc3pFspdS8fGLj4OmzEZcyfGyl3ZRn7YiLupjXoKo5xXqUTS8kLai6JmYkk6UhM24FcyhOSVUWoFp55R8ytDHnFI8R5KoTzmeY33MKbEEOvYpRXPb5+Stc6JIB944ixcDa1vJu9LzV7NrvcD3VrIr3tOLrnlP77EcjrStQpKTRvhKNTmhOsr3qGdE1e6RH89LzCfmuFVIMvrxIVdrDDmwtlXre8m/5XOPCW7VfOqUkm4tgQLXiltLYaXzTkvM2ZG2tcSMHWlbS2GtYzFAhcFhPsSG+AoOm+XhU55o2xo6heeBtq25CHMKjC6hwqBAn1JDtTiHxeck++eYIgNEszgRzsaxHvu05KH70KbVjtTjzeCOOTm0CtQYy6FOCGBta5DkPj/eQ2LN+UElDZEVHDbLU/FP9DeXWTE/X2q30p3DZjnUWxBIU9unR99Y6TE7LwV/zpb53olqORYFi6RhncNjDsdyBOykAdh3QRjD3voUhHPZtw7VzSFW2TpEN4u8XE/z9Yf36Gd6fZo3PoStTwOHnrxxzC4ejl0GunQYyMC5rlW2kUISOFK2MYWgC4mRTjjnyRxOyJ4xqXGOdxs7NbqkGdDipD3v47TY2JghkaOLSYQ0piGqdIoVWroknJ+8F4Q8meIopjvQo7tATroN9GiQvTcep5L3cZJFGLlOgHRPYWqHWJ1toEJ3uc5nEZJQWwjjgYycRV6elvAc0qTGuQrbmNS49xS8Sokh4TPOjjC2p+yKQUGMc26wxcBlnEsEW4wcixR0s417AcmIs+sLxhRtp9Seinrk3T+Ixvcy2ztyF2V/tTz0UUZxic13/LvFJXYw3zaKe0qTTsousJlnuR8Ooyv3w4nPnIbDxMppONF4Hw4fqvtovDBXKFT0NnAXJ158IvVERv7dyoimdB7hNc4F3UaNEc0F30ZlRP6UUEj1Ed0XGdGUj8N9yKkPRguFFmjbRviRgbRtuCM5BecIRzJGu40W6j1IYRvtgZKlf3cn778a67PP2j3K06vhUu7P+gyF7APqU/vtQ+pB2Pui9CDskNijB2HvE+1hDJ1WoePZ7B726EHY+0LgaMZDKGoEYYf6GyPQtjcdgbVQGmM4zkJLjOE20BwOHqfTViAXn3PKTF+jEay/9xGsf+rDWX92YVP3+GtTUu8rMAX17H8aB66mxhTSsxBCccfmlNEnfE4ZvX/mpICcoKeM9g/ZFNGzcC+hWTd7IKKhYENIT12cdLMppqMXm3J6zsR2OT1RbVNOx2LatBn2UlgME9M2BfV8OKNrTlMMCneDhY+/zRCW0HH8baaOwdChpuVQM261ZF51Mi+4xn0zbYea3AMvOP39UAxHXQofHYcyFhzP0jZTO5S21LRIGZvldCiIz0WqbJblUMrigymbZT2UvLh46pvlfCi64H43IgnlUNJCxK3pZrkesi3ViIVtltshj6VFiNByP+S+4CJ3ehiHXBcc4E4PRtiktUV626ykQ85L64t026zIISvePL6xFT1kd/lk1M1KPqgt+MKjb1bKQYcHCYywRiVs1SvRws2IpLWl90WMHjoRgm6LGD3grC4DSpbNih1UlxF2nXkcbXHXr2+GcLRlVGJSmyEaxzLaopI3q0QPljEWlbpZLQdpixG9GptVdwgJhKlshlysi+kswfIenlC+iX1HyCVK7jZbi1LDrlqszxIksxh+PqEaYgbEIlRss5YpCYELL/pgnGq9WONpnk9bFCNwYa1HkUmOzdpgWpKYJUUCsguhCS92Qp8LoVKVvlkXLxIO8yIhn4WQUBR9hUR0FokSL0QvCIEbknEs4rEGikQXIl7oxR5F+uUpzvwiRDlYbSQjRYZRNxvJi0QKvChRZBg8JXbtQY8o5igyjLIZATmeMn2KBGQXUaZPsUWRYVDsUWQYeTMEJI0ZBkXzYmYYeTOLURF7Et0MJ2oQfZzFWCviTSKbISWJTjMMioHBzDAoBgYzwyAG5wQlmWFQdIry6KoXnaSkQN+2GS59WYQoViKWRtCOMgQ/KBOMQmDsZaLqixTCEf48H4QIdFs0ddoTwSQivZdrlMdeJiyzSDEi4rTvXq4Mx8toFcTOoqny3MOpUhkPZUkHm6rYi+JF38bgqR4IpkesRJJkLzIWf1q8yFAyT+uBkA4j8WKj2BgIEUcCeMQT9+LwIsPwp3YgmseqEGLU5O96dJCi0BHizFhTVQbJ9oexZG7qerhVjCXD1tVZZk5Yu8QgkU4M24N4s8xA3dyNbQ5jLISMGDnyywcT5O8yku5xwIxtm0UG/Wcsllmmf5wwtpZ8V4IyEcaIKQ7Gk4nFzjLjyTiHERclPpmy7wh4nDQiuR3BK4RnB6uRB2KZjRskuCTEO/1TZjwIePoPeU8aDQI/yvSPkKd/NAD9l4yyiDgm61EKqiTKjKdUFI0gzzv9l4YaksGa03/pKKkoe/8DFRa7VaxHMRRclFmPmlB/vv3WGE8VlKOHRRv9V0WRRpn+p03se1zgo5ZD5X0UOv1jF/dFBo4D/dd2QEd62fvrh8Z4eO79jQO042XmX+3AFhdxXw+Vt0QEMcr01+RAzJrnhf6abxhGmf5aPrA27NAV+mvlwNail+mv1QN9eZn+WvNAt5eZX+sHj9pjR3h/48DmH2UPxzc7GLt/lOmvJ3aVosz6ohYSA0be0SGKAQfHK+jRbWaWyCUiFeWAGokK+kQ9wDTDfJtQEgqCKCsV6r0iWOiWCu+WbQG6pcK7Zd8OHk2+2SYJRaFsLVIBZlEVkK5X0C3KwnccaQFtoS4gRnY+ifMnFAbUSAuhW1QG5HdfQUyeXUtasP6oDQjMKxgYiqPMcbBbkFAdkJi3YGAoD6wzr2CB9hicVzBSD8LN9RBG6vGQuR4u6V2JsKbsWDB0VyNzTV2Wo0ic8mjBSFEliCWvYKQok5AzbDVSgUEwcUuFbwK5JBiLIGZ9GwhW9wp2XfadIK/wHRr2HfZvsJWCTvGNiLEgT9kNcmqmzN6Kbz0wtbEgUH1DiIUYnY1S8Q0h30DvCwKWHSEx6K0vCFj2hMSYVl8QsGwLCfs4lH2XiDAOA9rLEanyMgP2bQimhIlHGc+F8TTfsWFriF1L2iOA2RtiS9jLjI/NISfmwSajsDvEJgvzQeCyQaRO22MZjI9YChKf7Vz6J7o3SR2B61tESHhLCwKXXSLfc6ZM/0T3kPgmS2e9iKc4I8iCwGWXSJHwlL1/dhihaV0QuEJ8DwlPmfUhqOJskhcErqBg2NqyvMRGmh7UuaYsCFxxBQPBlwWBKygYZyKMfsr1oBn6xyWg3A7sYXqZ/lEwzlJtQeAKCqbQX2NjViRH7gI4ReBKSQd1BusLAlfYkkeCU6Z/9uSd38aCABb8CGe3sSCABU+iuHewVPrHl3DuswUBLHgTbqC4iymCP+FxhrRU1gOPwsu6IJAFn8LLeSn0j1fh5bogkKXio2JatAWBLJUoqRtgCwJZKtFbLBFbCuOp7NliOIpvrUplfJTzgoAWtuy9XBYEtFTGx/O2IKClMj7Knc1skcr4KNuSGU9jfHgraUFeS2N8lHVBXEtjfJTLkhlPY3yU64Kwlsb4KPclsx4t1g8aYpNSGuNLi2IaM57G+CjLguSWxvgo50V9PIyPMtvmbLUyPsptQWwLIRkv4w1RZnw8xxKnDH5tgWYQ4tIDv+S6COPpdT6HBym3WcZCoYy/yPtYXZShP8ps7pJI9NHbC3LJ/s6T9s4+fXV9fpf1bDl75dljtS5nP5wd/sG0LKa4fm0xRbaMxdQWy2kxtnizLpbzYrkslutiuS2WkVFjsWyLlbRYkcWKLlbyYqUsVnAeyRFBto3Fii1W02JVFqu6WM2L1bIYe+C1LVaRiWOxaou1tJAZbE0Xa3mxRqJJXay1xRqydCzWbLGeFuuyWNfFel6sl8U6rmhbrPfFOnLYFhtpsSGLDV1s5MVGWWzgpLbFkKUuv20xS4uxB226mOXFSHCxupjhwvbFDLnv++rEYZDEySMhiR0Ij4wm9mIQwQmzISF8ybXGXeA/3kDAJuzqhGhNyNPkXjGSNCE+E1ZyQnAm7ISEiExYwgnhmJCICfMgIQsTAjAh9RLGbELeJYRcwnZNiLfkiUsIsoRFkBBhCbmVMD8TEitBcglrM0FsCamUILMUyTu8AVUlAjMJmZMqb3hUDqMyoeMT4ZeETZnQ8AlrMqHaE9GWhE5PRFoSHklCoyeCKwlVnjAlEzo84SAklHcihJKwIROqO2GpJ3R2wgZPJBAkLMiESZ1Q1QkDMqGjE6ZjwvZMGIwJKzC5anejAbWcMJYSCjlhNyZUcTLewGxMaOKEwZjAOVELjAb+I8kBnEdWEzgnEIG+5z88fnDuhiJRBnELMfIhwLnbhwLOIzEAnHu+FKEBtCr/eSIFfYBz/H7UI//RBzh3gzAyssC5m4MCznHW0V38xxvgHCccfUQggz7AuRt/ONQoFv7jDXDuCVtu5+Eeoxb4jzfAuQcDBZx7fpxnvAg49wQ33FgENf/xBjj3DDv3QiMzDJy7+eYZcm63uVPoBlukkIBzN9fcO3M7zf0u6eMf//znPy//UdmcnCv51WzOOKhBNvqvnKyZGZukct9c0XBPdI6aR6cPfHBChOKDUyAfHgC5PxsjHNv4+eGYL85vvj/ePDhsExUPPjkrTudHvj6+u3t89YIsbBJQKcbDtCbWwJ+/unxxdUZUJMoPvs/jT6/JeG+eYXz+7vLn+eOP7+4eU0/i9vPLHy5vL6+vbs8OVeiRJw8++Pn5d8f9ABD9eTl6KPTg5S8vLm6PfjAHOTsrT8POPu7LZ99/frx6wbmitCaSnB0H+6s+F7LPP3zt7tWeaH5qsndPYvK3/2lm6Oj8d8zwv/+nmeEJQf9GHH50fvPgWNVH5zc7UXgKNifVYNNXT74Knvjk5vzHOJoR5S/f3N0fA4nCPAkShXkY5Ms3d59Evn2crSMJHiZyNvryzd0UCUziyzd3n/pRtNn008t5LuBnGfQ08Mrnl3ccqNvLX19fv/IEeiriGMvH11d3129vbuephcd3czgfSMzHd3cwsQupvyAL9N8oDOAVph8HnJgkJU5ErF767dXz397cXM+TbHC2F705XX369urZFAs8pPhAilGcKOQpx2JmY+ZPcTaG7yk+wPfnxxfHq+cPz94wuqh9IGH50H3l3vd+BpFP7GLh1JD1m6S4nF18xlmI4+0HMnzWfvXm/BlnAbzv07nAB3M4HQqcdYzx1O790Zya7tV70w+69nYfzvtB5f3Bo88ubyHIh+Ohiu/N4bTEtPd2e8exOnvTvXY2/GA0tPri8ury9dvX//14c31/1IMH7x3VdBEfp16e3Bwvjjf/7fP71lH/YOGi4uE0GenD2vt5Ru0nx4vPzg41gbRTzTdnh/5+zdMz9gEeNPl2Vjw5f0h8T87foy06P1Xd9+xVPz+i+uT8+XtTZ+2enD//4MxrVP7Cqdcn588h9qf3yzNrvn2vBt04zyrR4eWz7+dJpSfnb+Js6dMpNE4V356xPXl28dWzm+Px6tPzZy59GAli7cHyU4QPHpAtVQ/xsb/1gH9oQvGee6jZCSjIKtrcvEaKnrW1+EEaKl1m2JhjcRGDcxpD+2Y/70TLzyiIz+Pu5vLNJ8dnl6/PX92eDhW5SJ7Wjp5Mhgez8wYfTM/rHs4PIvHKBxM8lU+qKowUDJA4C3166yQC/Z2Yjp1aM5/9+0wnhM/1766ujjd/YHq0hNX8s7dnh3/gINGjR/wn+khw+h+VR2Mvt0flkTRKmmhQZpP7klc9rNf0fvPTs5Hms5EelfSoPCre0vv+D/rvHxEux/Pnxxu0tZ+Z8mU7lT69vPt0J5o6icZPYoHF0xOnJV+1Z+ev/GWw/39cX15RuRsCH5+/eVj8+vL1yZzsY5iM4hLjd6/PXxz50EnAf3x+9fzV8ZuXl7ffH2/+cH71Yh7ejvqPrt/NusBe1PpIHhzl/PvL61eXV3vtPLgYTT++vHn26kNpPx9xBpVBP1CATzG5f/vuzdOHZs9e+e3Dym9/qeVe+V5LGn5x/u6Tyxd+uwBE+OXN3cvrj89fH2/Op/T5q7trt29vLjBc4rIClyv7tQcfmjmf6ae0/CWHDckxHbaXL//07PnFn/707k8/uUEEf54dzlBJ1y9uzl8/ur549EEbFndKiObop2KXXJMXHyz+B2r1r+HC/f7Z8RVCTHrjDOKv+GDTAfulBUCIzAV458b/LmDODmfvKPyn9toY/3eXOJyFiT6Ns+Wuyp/6gfEAownnmS8+vby5nQb45+c7xLUZXJPga/P6+Mnl7ZtX5w8OXyMNTvKERXNP9f4w9BfXzz8//26W/4In+a/D008f4snp9W8GT7FK/zbv+m8VT7/iD//r8PRPH+LJz0PfPvq//vdH/9ujtNZH7/z///K3xWUrItGF38lcmUGVvxAduccfYi/YLBiOWxm87v91HvvAt3/g+CeX1b+9uru5xEqSkN1f3759/eNkcsCd4YH3yyAcng++OH+HA8MVHOHXT4/mvuL31zevd3OZ1Zhe+YebMY1dmtiMSbEXk4ltmVkWzZKHHn8jecF51mJWW87WtGmjihoTbd1U2TYk3tWLjWZJuPWDY5RrotHQmusoUoZXNS3WhzQbyQPXaU29mPU8Uu4t9fj86MWaldy7plwqrczMWk1dWgdY0kp42lrqqbbM1k1a7wdazVSoYVRVrOdax2A2wqAKXVUZpVDDmErWIqWPxCiFIeUxei25t0ENI8qlJsvs9qRVGE9OUnKuzXRJq9K71t57llb5sLJKKoXAeRM+rAxHeipVcyWcjMsQSHc18N59SEQcuLfi59HNv3Q7xgexlifnPxzxQ3710iSiug8Mit1++H9gTnAVUFprbgqJ1FR7Ltwvw1VAaS252rA6+ujD/H6VeRUQy1hHV0vG8p35XUBpta7VWFut4WI9lbhEZuUqIS4Basky1wTQw3zUrPqz3stoyU0UnY9svmWpZ2iSfvZHOd7SPAp5C/v1KCGRWJnfXREmA/rq5fnz6x8fuJwfX9/gD50/v+TaCfjuhMzvbugAYRW3k3CNk+MDh+CXI+EIQx7OMDiS7kMf8BQGvw9ST3fk8+urF0eMCEGo3l9AFYHry6t/430rH1DU5+d3RzfR/gW7jgnv1uwHBuyvz3tNubkmuJ/8vhYPfeLT3Hk44wozoPCvpmXW9HMuoAlZGUjGB/vy5vLF5dW9GipFtNWivY1h42f7P/+e9ZhK4NF/fYQWeOSu6l9YlfepYS/+f7UgWksZlliXoiZ/lQX54nh+9ejdo//66FFeWyVi8L+W4yeWo6692P/fl+Oru+ePPjn+EARSMEj+14rEijiN6NrMU25+Xa7u/BT6ZC/+rUiQ041efovXr/3n1pLvW01P+cs37+9w/VLg7RdDuE/Ob45sIP/z2d/dHC/ODpL//P4YvKu4v2vqgUho8vu83FL76PKKT//25uZL7iLjOeUvfzjeXLy6/pEtMKqI7V08CF27KT+37gJ+z96fTx7f3ODn/1ryVERd3pyRoe1JVLvh7ga7rJK1pHL8DadKVu6yS9ULZZWhYwwvtLXnkUf2gqxJa048kUXWYrV4M1l0JRPGfQEKY0idhbyOXHoXf6es1pskXAZZ2ppTzZXhyNLXQaqJw7aWnnLduxEtwwtpkTWTa0QzCsWGDGZAoY1c614YTXJnOmnRNWlXiye6Su/j9CSTl8PQaFYU230WqjLVvWC9+RLQrNUx2v6kpzrG/oFeyHea7/TW6umdPkqr+3C69TT2ZiPV4uvOpwcZQqxHFIqW/Z0hLZ8+MKS3B82Y9+md0TurGx+wdJrcEIxsf+IH/1dN1lJztI41K9cwThz30tSbypJxu5oFWspqNUsN5HUs9zZxLPhdOtdB1lqG7bMQGVL3BRoyRoqZo797p/+0lLXlMmYrNJnMObTVGrdCequxahslx4dtbVBCFLhPVBSc2CI4T1U14NyG5llfLAlUQJuWuHRywjagQar7gNQCHlUrQ6LekpXWJ1xLyWPCYxgzNYhMysj+rq4pZ0eQ15dmcJnDtaXiQ9M1NZXaZn1rcnq3jQYZe/ueFOIImKypE9yg4VlvMRVdU1c5VasWHz7ottVKqZOls2N4cmRfR+UGTmc8WeH1tBO3qqlGIa/ZxPwd8GWlFyYiCRx1KROrttasNgKtAsJanuvWSlP428c8qrW5bjj96uuZV8k57XBXyXzVlrLmXGqsFTTYcvM1rGvldlafZFul9uGfbGtPrTvYV03FfEX62krxz401qTqxGpSfg5eAq+XeZ5tuTbov91itOkHaYitfjg/aKlVrzM7w6jOsTxutDHfCXGzqBGir2hgwDG0gXvFZA3exHVbJMQSDKYfMd7NWsf1drQYX+ne0ZScisJxXy1YhpBDPozhP84C0uwaS01LXgjAINhzcUpsmU8maS+6wvhM9wRHvXtfepIuvRl57b8VHW9EGDdKxpa1VpQROYdQ8VwaNkRs8MqCHQgpfwKUXMQm4Sxd4mTbWcgHXY0FGtQzLApdsZYKNXmeTUYtm/0xekw5nzcGE+0jibfKakxqIor6kUlhUh0epNuEKeWrUV1MF97RpWnub32nFCiLB61vtMZ68tt5aTCujhU6vjnEaWrNktn/GtIE//4zltg+zWREoO+pLafsQrDSQwm1zrvVKby1Eal1LzS6mXAVKKWnnqySF8dmSVy1VQuyUFRJH+ZgzhwZ9iUsv2AMM5NF78WWRdRRnJRCQrZV9lim3hhZnqNAGXxxLWWtLA3ocjEy6+mKh580166BTqQIpjWVgMwxkKbC1ovDyWGztLSszZDSpjEG3HTEvCt4dbH3ohOGZHNVa3bqghY5UUMDAOWWBPR3WlNDsDmdluA4Wcn53uCcQ7fWVdOEdzu3UptbMPKJNQ8tM0GCngLtCXxMukPKE3dba4fgkyC2rWcrOwWB6dDO0H0i0kgdLDq9ZKYIBAJuKjBb8ImsHoc5HupahrrhBUbGRg3DL2kfu8OwAL61J8CY8m1yqxJq3Dh0xQLXWoVbg2tUlFfBo3Q2X7jZVEBhwGV0CG7p2K7B4hzctCYQOnFutEDpwk1rLrB+5KdPqS3HREuuFOZqNNaE+ixaEO3CRkoI0yloMYqO21lTKbFFNRwymrE11YHDQppVc0ZQO157VMVPW1lqGCL0eIySAWln5gHtlxQIeY++0datBDQV4oCO8jZuX3AHpJqoMDSne1oSgm/aF9uEC2WVeFbdPXQpod/EOOxGdDqE31hRUDqokh6Bw9FTJyDDQkFopzQehayOE6jBqvNTAQ1k1VcGM6IhyGwVzsDO0WjvKqS99TZC80/BYU6sVE6S7IsxoAbpNSRMoacAmhbUClq4FnQOspbi1BpxTLSDF4crgJmzjVF+ylDTrC7dsz2+WXhuyiXfLsIJ967CZoAaAa6oDxeUwRw72esnK1KO+MtsAmyCaJzxQvwF33YdZCfXf1zfQPts0ZB83hIZdJBrGHgwihWGASEvFQi20dfTqDOXYS70hmVnGNmrGkHMmqkMxAuGQ0rKxFJBTzy34EpmaK5/vyNFujRVybORWQBjD05ENe9aHWluHuYBHsUCSrlJzRzI2WLdhsAXcW8/IrgbvYq1NOGctmGjU17CFALsNV43ANlpCejdYt6irSWDtKQUey5pxnLzbspYSHhltaqoN5nW4dHeaHO49x6zKWq13eJD6xkmA2VfThJb3am0pZgjf232TLBWy9TZZO4I94FyCrPBFfE8rqksKSpqal98RcNlZVywc0OJoVMRSrL/1Xpk9qJOBbAEuK7TNjJ2pelLIDEZis6xPTspNc0xB1iGKUQhectUak9TVBqd1qM9rM3WxwVBL7xlV2hiajASfNsxidIXDfU0jTwT01aoVDP2GCZhSir4Mc7zuZJLUnNmqM7YWaANYCsJlwubuHdUY1ZiVDteekG4Od4PCHDSc54Bz6oZ4oJ6tMCSXw6rosABbhmICDuvZ4ZxcIgfMVt1skzXBplGvbq5M2Nlswsi8CQ7Wsjp68zqqOJeDUvbwQkPkNZXcEdvIyVyNecF1DZN6x5zktgsRk+bUAOoa+39OZI6uFNxYVpPsdj4oyr32kF1jTfx4gKMRuZrGPnuV7PKQYbPRsC9/g3XmdEysQknV2TpiDMB5sOEY9dX6XGi0c3KXiDZskDpe8soJJ6ZSsRw1ogbAWUeGqoBLam7cOdyaYcEAV5HTu7WUHojMq0t5H0Je68gDZvH2VgouEXBLUtA6AdeB7xLwmEPGWkhYrF7NWaZ7uJ2ai7hwizaKlA/82jpyLJZbram5ooIde6+7hZmb+UKAO03VKQtWY184lEHhRxw6Zg24G9anTjF0lqtHl7zo4knoTUeDvVjn1LJbVsC1thaMgXlkLnp91IOtX9qXVccomLnAlnKDJKtbA8nZvboloRkLELjlWqJNX0WHYsNUgjUtZwQcsDFJbz9WLcMjfRX3NynCCrDj2XnzsY7aCrqIelx9pl6xAJii19vKQS2EG/VSiyHcHB6qyGGHjV3ugJWjZj4EXOfhQRHaqKrs76qWGrKC+moYsNGmdSRXwF0IEgSGFfS5owF7msJLsCqGtBviYFLRHq4TEZQK8yAD08hMEtwJjDF5isAMEUbwpZ1MAGBczZx34jbNgovoOOqcOwMmIjFyTLKvpTVXMCzi4DdvfAKyJsIl9fibQr8hPgGzkh0Q1YXfIJGAqw0XiLTpQSWAQzGRo4mJuFVKPYkWKI4C5ZWRkM4OW1OWE1iyB9gcbMLCOshYZmux6jRIvabq4VeHOermY9dVpblTFfXBZg5rcqYOWFwVTdjjKAFrRstOuGKqTtidTG5WJ2yRR/WYiZu38BToKqvCmlPDFcsN1nEWTN2tW0TmIBDokh5r2DyICRp7x7IDXZX8Dg/hgLrRzEKsjtUEC4w2sooSFWF4smYzN0mAmzX36Bh2kir7oufS4V6qQTqirkA90qBawJyLm/jALctAgANbEte9BSMqm8elgHMpLvyBq1S3BYAbdlzh3YL7nKNbNzsqJgv1xonY2YYwWix0WQ2d4wioa+J4qbfx8M5AsBTgXD267nARtwcDDmM+YKy02b5YYgWj3tx1dbgmxR6ZcEGCBILFJUcowrLmmlznNiSKdt0djJa6e/guQMfwmBCYtOyeMojspfqr8Bpe1q4fObALEYCw0nPGq3UuwYWe3FBUK3LesZS4biLglkfFBGAZs4rbnz6FXBKyFBi7z2fc1sJmhGO4r1Kloh8LmoYNcf/kWHPrFdFSEAl9KEKpeDiwjcCMQTQeEKXeRqrIT1KGUiJ+WCs5R8CkXJFcBFwzp3ADxvb1tKiEVh9jb99L89wo6nsfjRwoh41D8gEP0dp3WDn7POu1cpWitx868t7XUA6hz/rM2eYd5twvLKxrzuTThMQdnLx10wiHZNQwF/OqNfYCQKXyI0uwXl+tDw0TRVZpuaGGnPWydTwGlxhFd3azwW93UZ2R/in4CiXaCvZtQXH24XExMJMbKVjUj7VXsouAYwZiNU9YGmHmqM+SOMXubQrB1rlytfea5moRuxuzvtvgWLW3H4MoYMDW/CA39UoAp8S7LDlXEHg9h9djCCpZxxyCStXOL0/wqjRijxPurUp0q8QcPKWNNqNz80e0HzZaDF+x4rmV19tg6e/fseyZbrxqmd8xm01ybfvILDey6tz/Lb11DDHkbxut7DitxJxchjZsSpe4GLjZfQ3QWLUZnlkIftL3gNkSGB7zg/HUZOKurp0T8M6ofS1uAzjvWMrBjoyTU9065ooXv1TAqwlopUnCVklPjGlpRncGXBPRi4BJL5zkn5MYx8T5TtaSiwSyyHLJPQgjN462x/dzrzmn+H4ebV9wEhH7pB0CI30it6SBB+OfJ5TA8XW6IudwfqVIq3Uiokjv3CwRTUbjKoWAzc/UT7joHFmB1iazFyFlbu/KyLPb248pS4oQFArpTLBWnaMwolJlgcLMKSl7jDj0Y6+7IcKWJWYvGO0tjxClZSUUu2uJNlgV2nQiDjWYmmFj3/ado3KRvMs6nJM5ZeeKOWytVYdOBCQCPbFyuWYZk6NKyqNOxBeCLZFrmmoStgy931oaWizg4TtJDhM95LINxtaIiXFvNHC3XD2xMyXiGzlyYT0BlUtbaIPxyWUgDpMvOlkTRVUntntNBCeiTcVUCOLrtZEfOut7m5jvdTRumPBPVvOt/glXOVV37kCY1WMngl7Ztpwj8z0eZ17i/iPDdi54cwpfjJCeu3quKbO4gnbGFPVdYjc42Hl3edxWS75TBDf2rL6ZGGNgc892tuvdb29gbKo2dBdoxldj/JkJTjlaci/z1Yo0ngvbUt0ThlPD/d0R0Ymx+ddHtZ4m/rHXQkJKYt/eW3C1x1wabrjgUhVymaUVLjpx0LilxROcmY/jhSTeXp1kyPrjmhZv0PtUD+xUcLeF1xp23QRr5xdiSE7a+UzUsI+8MqfUA62S0/5+5t4VXgEpc3wIIosPZTaJIwKpaM9pxQ5i4jkCAmxtSyWwEDK1eUQmbLDsQVQ4sA0j3jDRJco1GmBIKvEkBxUfYbKKjl7HZMvcM7/z6G1K19FmfR21x4jdbtNJ1aOQ1xxUTeo3V+wwQdFcT6uqNXhGkJYaDRDcsy1Gbkhr6Vm5cIcvDBk1tJoMS1wSQy2/38iP+QBa5tAiyelJzS83B6y1RQNNvfTAlqYx+E1Kb2AEIx1EkcSclLtdYpGUq9b4LZ8VEMk0wcwtNFGb232tBUqVy058HYBc+/hG7m6fIlSxpdgucUFaylSTZWVsYXBiWO6ukuu9UbonkqdEEHu3y9BXtiujjCsTMqyUnLiUiHcrge6J966jtmleGBI2ZBW30GhUQ7q9BMdUydzAxBr37IcOYuVrDjGuUqfE0ly0hCDTKshbX6DW4WIHByhz5kEicEUOKb9QYgwtJ9IbHKFYu/zyKQ2k4aIHaK1L1Kr0GpI7a0ZIeAMtUvlNJlzw0kfQF9vsffamtWggNGutLUyorLXVmGUmxhkSG1DDzAMsvOYuK+InbCBZTfgZUdBIaqJ5fhNcl3qX8JQCdyVxEZDDZMRP2U/elU7VnRnD1FHFTx/4UqVWxtiNk1HEuKQMHCQSdQIzUExQvzSUVTCNpQmpSuOXTqHXmss0BHSwUx44kGI5PptL5kKpyMX230cFNH7I1WtJhhoh64q2wiVGnPXIVrmxCJCDII6BUjnr4XXEx0IBlDbSCFOroC34la01FTbrerTtJFVFZa092Kj0xtZO1HIQZAfJJ5m1bPtPcDfmSm82aQuwM64lgR1LHs7H32x9dwjr2khQc8vF1tFaI3rt00ZaWgmeT1pklBDTbNegzOmY3aDUpgXdPB/Mh5kIPJWgSq5e62HESBk9+0LJKFonO2nOI5oqJBmoRMVOsZlrznG6JWXCUsEjJVtqEyt9iIUQrpJ7fKAWLF9fVTQVd6Vx346xSeNDaKK9qE+1KVu9AeahgasGUkLstqp1zKaVVY+m7N+EQdqa9BzL01pO3LWUVsDBHXgOljytqEYCUZgDrZXOBW7RYHAXXIDc5DahHYOEZbEwMGvyWjEdpi9eyBxxD8J5jVhgELaL4jwxgvTgsjzHGkdPpr0JxqcVYhjKgRTRRAe0FhTKrnP6mO4YqRhcuQaDWcc/AcxVScYDLNzEtVM9VrC3xbZt4bhV9mljmrkRj+KtNnLuIc271jT9w97YdXdS75YGV3ylNZGPx11lgPBw4IotwRpYGWTHxRBGG36VGm177sFMo3fSA7xyYAUEJDawQlaUZmoxxzHybgWPQc6CU9NA/YfBM0ZBNcdrBfTsoDFw3zXXlNxlRHBWBJQjzdYyEBsntms97eaFEuAOKZ1KssK9cCCwIptiEgjIkid/ITFCkgvRk7BvhOjpVPlVeg+phIs49WnGpppSp5hMw6Vqlh59Vyy8iKewg13iu50zZHN5kxEr9yk3a9PhJPto2p7WhnZnBa7sy06H3JDHOabETg2r4dBAiDpkDAaI2/GcPFDG5sjnSrfhC8BFblzq6efjNMe5OFE2g6OuhBjhXregFiBO3fgbJAJMKI768bRhQ0SKUquFXTYszKqpkWDpFJfYQOXeNtChKRuXoQLjqHN5nvOZ2aizDQe9emCB5NqcvQmJH9P2A0mhO+Lsl3eTU989rjwSv8LOh0uTNtV/LWwGhYxhV8gH1DOWXRC2wj1ea/hpTpY2BhvCzDs1mXjB/PDxedzdeRTJ0Pxoo2hrIQREwS+vZiII/pHM2Tdf7owah3kk117j1UyeYzxtpB/7056yUxHymjsg/ehi5g7GgJr7N+RDhZHLwUbu5/WnJTunUZeZu6Oqq2caEA/juE7EisGC5IZr4bCymjO4lHFspgtQR695WiWDbMwwz0RTryGgfWM3RIyK9akDlK3WYER64YZiNwh648pFFE7RPZpF1tm0DoirTWdu1EESDG2tYS372iZP4GG20iQsAs8XmqusJVSM5MFlnpzyLGQhOjRG9ZlIVbJ3qaukI/jT2tmz97pRs5umGFKzriXkM085L+pCDj9pcl+DRJxv2Exy2S1NSljfgrnilgkQV8nyFakRawPyJXNUGSkcM+aVu2jE8h1XlV0EJ3KS3XX3uVBjGlI5NRK2pt9vPZdoDiH1UD1s28yDtZyYjXXBZmYP2PHD9ZuhF2saZVrdTVqp8d2eVblwEiVQdfQgARxl7seFbXz/H0gJzgV1i7XsTwthS1+oyq6RY7T2Ef4BDnkNDiKZMpa755AZJMSHXyh9FH6ZkrO4qDH/3CA73BkIvvbJyNBw+AiPhikmI8twR0FwrQJRcLwPbpCWukMxtpHDtiMDH7fKQ9x4FSH6Buzv+au+Ion9uj1ESMknnUhSrM0ZOKWOWxySj93tENBsKCGiWVZN5EtBkgnFsNt3WL9hPJVK9NY75EBHCXOlkc4cJtVw1ykYp/gBbUdMJzvKEYOB4KuWc+2xCAVxGehoOpeo1RGxVekljAwZuDuxzK3Mj1gaXKrKSehCWopDqGqnAhs1uyoUM7ILOPtMvMwFqibkDNPWJMR7A7LqMyVUzlWqPCUrMZ7qNOI0kY0CfQFxs65DHoqNnQgZc+9prGTOuGPmrMQOa1jbyFrOH3g9+XdcoUqbbuRpOBa4a3Saw9K0Ti9QOe4e8k2tcvsJVt7ARnOwjJa5HBXTzW1Cr+0ZS9LxPFqPaXItb3UjkLtMI3pJ6KjE8fIC7fjy1jaVOP6EW6bSRyyfDM83CiTkkA3qwXempqmS8eGQWyFAhCV89chTgRpUOMTjraTCdl7XcnfHUknHdRmg0rlslxfIjIhXe/bfSKbOFQ9A4RrggLicO6A8o826lkIAGKvBVoKO7Or6IiUxwhK+dilnVpcRJ2g95Vk/yOsNa0J0TsVpqIc9rZlYhSM2cwd76PfijhIfqyJxDUAiSTPsfM717hYaKsuFG5cFRzhaFPwzLMmIWqfvyg27DrVCEI2nvU0dhrydFxKQ/OEGqKZSQ1tB3cOVlEodEVpSsdZcQhF7jUApaRJx14BqS81VLQZkd+mq7A/56iqHb5yZdJTmfI2ZEpsWpMYWD06pjsFl0VxLgLnjONKBOJ91Dch5iD2b8KCM9NS5ie9MZNa5ABU4V2l75KKSBjE9q9HEuM4ZX4n8tsABkeESegODOQcIS7YIYRTBvXAkEX+IMEir+I9OB0g6rmzGMCi1ulzgUuTiRqIoARkXb0VqoJPt2Z1jUuGHxVAeOXX3bMQKbM3kUx7FvWhlx9hJdgTFIZ79l79ZcLLjeOjc7ACmsbNVLtOzon2E+zVXdku9XdOIp2omKTneJfsj3m2jumDnQEV4UkBhZwFNjOfG5eETQ6OQ3hDB+J7ut1eUmMwMwObKjeXOQt04QOAwoSR+ScHRM6qFk0ZeyNQ0Wgize+Os1sM1SEXJSHOs1oIM8AatYwc6UoaUPDde2Mjh4n73KMi0AdIxwicQzjmEYdvSjK0KNixflNG1hTrGTXN8EyeMuCviKGiAKMnOE5hCTs5ZWwgm7uRwx504+XB5oRwzCVYspFk7xxQhLQAUkfkQ6GBzzPWkFrpzFBXVWBEltODmOVB8uPBpx3Ph0zuKmnJDe+BodNC4CzoSWLgSw5monSK6bnNPm24QzA66l8yuoksvGbnnQK5W0RrGBQlgzrsE4TIZB3zZk7XCLcfX7zMaQBQzlKAN1J2jhgOJsfo5TXYTzod4Q2mqsU8ihAdC/g0b3LgP70gZzuyoltjMc3/AcajaCjffwzToX+etXC3cDMUpC4EFIkKcFR3cXw9K8hi+SlqKG+ta6jQM2OQL1uL82GxESBsS0gKPOUEUtuydZ0kdCeMDeY1mcUnXwKL7sLaSDUU6ZWgj9ggixMU2KHEHx1c1kolCABqbQZOTspDU5kyFTA9DTdkjjuAMm55zc449rxkRq7hDERppozWNUOogEhqerQ2uSQkkVekRJ8h4u85dpXOsBLbBp3RjA+OAi+4xzIit+or4KRmneF52ulJ2gWKVFBPYn2bUPu/CN9WtGyXO4tJAYdGwCvagIYsZmlcLiVK++j3tuOk5cE9y/+SV3npI2MJZpGAuooaBOK50CcT14XZhIAnEBCOZu+r3ODI2OicjkcfA2AkHcYA46i3VHKIfqzPNPQn21KcdwXUWOUK2ORMic/yXzIEq5zoS7Gfstg3iu85fnGnhsn+/miiTDseSsxe3+6RT5rORWd3IJ/GrOrGSiRixK1Y04oaKmbKbZH6/PUpHa8hFVSMp2VFTeuQMwD+xk6GFU83BNUWHywXHjbslyklfJxhWnd8RgEkwGyZUgrpJBI/4H67GNCpwn330QEFsWoaF9QPkuiTw1FCkO3KQfjpNBeXnI6bA64lkhYmnNgMbRMh7mTEjyZr3KMNo3A5Fa+WwTsS6MfLCA0yYk/z2BhKvjsBdM06YOvpHJpTnCCPxTx2j/FBAZFoIRy6dVzjZHGYBQYdwXbqVqcVMewkVnUrlhxOwagsRc0dSZrfGITN+WQP+qdnVIrqh70qEW494SFJRWHuYFOHhlF7LVCz4X/6RMvrUcWDIJYLfMRW2IDmRTrd+fVV0W0hLcgZmOzEGVax4YHPiaCbmEQjCzp5pYbj3Y24P5k7U1BedyCpJaI4uIzYQmy0sHL/WgREx2D6ZGMr83gIuEZmN4fyUkkbocTDU5h5EMyUPn7YjGyMBNGtsMjgntSkvJWsv4eQXkxZBzsYWOJ36CRl/wdQ8XwlG6iVcTikRusZaiC098h+cgdmQHoE/CGgKJWwJx2SpE7mldZ0rSi5f4G/UaWgUSxHyBRsh/qyGgV6s7aglOh2sRpb+DsX6+H1oSNWJnuKHwXxlHUEzYA+C2GZklUAQLrnD+Ah77oulYTMuDYLmFo4MNtycB7TWElQGhqyE/OLQoatmEGS7PjLWwIlgcPg4rEezLt3nQiAuOIu8NX4pxO9W0wipSquJHztxBPEKEDlnQbDQGowgxVqwhGZynpxz2F0MQVc7P/ICN9mIyWhRizAJHe2M0NjOdsbq3CzmEMbjXGl+wgOuI70/8GfthL/eQ6oWG2VH0i4YT9zMKbQHOCKI90DMQeCOC+Wk6Yz4uL6PrTdwxEart7FkEgQFjtpMNARHU04p0YAIlpKDOj1TzoZEEDFVP9fvlECuIT9UE1w0wysgKUdkxKOlvuQgKfYppJBAzyKBpDDnuvXY/3KP3OUtOHK+Z3s/cstAUguNklOOiD332NVp4Zmpe22Kuz+ZoGZ+0IbVx74OhIzE7724DopUG4TatAU45j5x1EMdgJk8DTszfvTIcWl5SlUjUjDrGmQ1Gal6EuqJkQo/pwfzgKRwSZyRbA+jcoJmRu4IQ85NNXC0e7Ck+cVig6OQBY6jGlb7vCiMXmo1yRE7BUc9FBOZ4hHSB0doUbAgjb0OIPZK3JYGR2FqgaOIVHfr3U1ucDQXOhXjB3PgpJpyGACaS4Q4QNLkkFzJJ6QdXD9t48z2jS9cjWQYcBQhcC1cxzix1fs0sE3DetFiSPJYc3LdAhp111C7lRj7aXBPhEji1sj8H/ArOt+9mlfs/SvuZeZKRG5K4g6pX7rD7q94IyLXAnEd4m/S6ndp9CK1kXEUtxVyH6L11MksaJwPN+6fi/sQUX/EtLFVMKPjPkRuRGmdu+hICm7efL/0UEbF0LPBBSx+9d/pPkTLBX1aMmTlvwZyug9xCDxKSLYoaVcP70OctygWdqk6jz68c/7feR/ijqx5H+Kd33m9X1H1i78MtD+M+6uUsYDCB/cv3t8J+Ov3ISrL/rdwH+Jfut37P8slXhyGyTr68APL48MLtPxGLX684fLN3e28mPKzyxcvX12+eHn38fXV1fHZ3f1PRHx6+e74PK7aujh/dXv88/8EUvDADQ==').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222657777', 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_1779222657777();\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
}
