{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "ebd51f04",
   "metadata": {},
   "source": [
    "# rf203_ranges\n",
    "Addition and convolution: fitting and plotting in sub ranges\n",
    "\n",
    "\n",
    "\n",
    "\n",
    "**Author:**  Clemens Lange, Wouter Verkerke (C++ version)  \n",
    "<i><small>This notebook tutorial was automatically generated with <a href= \"https://github.com/root-project/root/blob/master/documentation/doxygen/converttonotebook.py\">ROOTBOOK-izer</a> from the macro found in the ROOT repository  on Tuesday, May 19, 2026 at 08:29 PM.</small></i>"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "fe657f13",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:44.376872Z",
     "iopub.status.busy": "2026-05-19T20:29:44.376750Z",
     "iopub.status.idle": "2026-05-19T20:29:45.337759Z",
     "shell.execute_reply": "2026-05-19T20:29:45.337226Z"
    }
   },
   "outputs": [],
   "source": [
    "import ROOT"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "b1f64d1a",
   "metadata": {},
   "source": [
    "Set up model\n",
    "---------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "9c75c764",
   "metadata": {},
   "source": [
    "Construct observables x"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "a09b5af7",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:45.339941Z",
     "iopub.status.busy": "2026-05-19T20:29:45.339793Z",
     "iopub.status.idle": "2026-05-19T20:29:45.497519Z",
     "shell.execute_reply": "2026-05-19T20:29:45.496721Z"
    }
   },
   "outputs": [],
   "source": [
    "x = ROOT.RooRealVar(\"x\", \"x\", -10, 10)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "c419d6f9",
   "metadata": {},
   "source": [
    "Construct gaussx(x,mx,1)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "71c96863",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:45.499244Z",
     "iopub.status.busy": "2026-05-19T20:29:45.499120Z",
     "iopub.status.idle": "2026-05-19T20:29:45.629703Z",
     "shell.execute_reply": "2026-05-19T20:29:45.629018Z"
    }
   },
   "outputs": [],
   "source": [
    "mx = ROOT.RooRealVar(\"mx\", \"mx\", 0, -10, 10)\n",
    "gx = ROOT.RooGaussian(\"gx\", \"gx\", x, mx, 1.0)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "c498be56",
   "metadata": {},
   "source": [
    "px = 1 (flat in x)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "d8c1900a",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:45.631723Z",
     "iopub.status.busy": "2026-05-19T20:29:45.631557Z",
     "iopub.status.idle": "2026-05-19T20:29:45.743395Z",
     "shell.execute_reply": "2026-05-19T20:29:45.742914Z"
    }
   },
   "outputs": [],
   "source": [
    "px = ROOT.RooPolynomial(\"px\", \"px\", x)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8942f7d6",
   "metadata": {},
   "source": [
    "model = f*gx + (1-f)px"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "346e1287",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:45.747423Z",
     "iopub.status.busy": "2026-05-19T20:29:45.747298Z",
     "iopub.status.idle": "2026-05-19T20:29:45.921647Z",
     "shell.execute_reply": "2026-05-19T20:29:45.920877Z"
    }
   },
   "outputs": [],
   "source": [
    "f = ROOT.RooRealVar(\"f\", \"f\", 0.0, 1.0)\n",
    "model = ROOT.RooAddPdf(\"model\", \"model\", [gx, px], [f])"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2c7a7b12",
   "metadata": {},
   "source": [
    "Generated 10000 events in (x,y) from pdf model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 6,
   "id": "3c72bb48",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:45.923258Z",
     "iopub.status.busy": "2026-05-19T20:29:45.923134Z",
     "iopub.status.idle": "2026-05-19T20:29:46.068932Z",
     "shell.execute_reply": "2026-05-19T20:29:46.068481Z"
    }
   },
   "outputs": [],
   "source": [
    "modelData = model.generate({x}, 10000)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8db97e5b",
   "metadata": {},
   "source": [
    "Fit full range\n",
    "---------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "66abfde8",
   "metadata": {},
   "source": [
    "Fit pdf to all data"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "id": "b6d9e1fa",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:46.081434Z",
     "iopub.status.busy": "2026-05-19T20:29:46.081299Z",
     "iopub.status.idle": "2026-05-19T20:29:46.274122Z",
     "shell.execute_reply": "2026-05-19T20:29:46.273415Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) 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 812.744 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n"
     ]
    }
   ],
   "source": [
    "r_full = model.fitTo(modelData, Save=True, PrintLevel=-1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f3dea233",
   "metadata": {},
   "source": [
    "Fit partial range\n",
    "----------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2fcf987b",
   "metadata": {},
   "source": [
    "Define \"signal\" range in x as [-3,3]"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "id": "9c53aa89",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:46.276070Z",
     "iopub.status.busy": "2026-05-19T20:29:46.275942Z",
     "iopub.status.idle": "2026-05-19T20:29:46.391032Z",
     "shell.execute_reply": "2026-05-19T20:29:46.384617Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'signal' created with bounds [-3,3]\n"
     ]
    }
   ],
   "source": [
    "x.setRange(\"signal\", -3, 3)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "3f433df4",
   "metadata": {},
   "source": [
    "Fit pdf only to data in \"signal\" range"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "id": "9fc9bdf7",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:46.392728Z",
     "iopub.status.busy": "2026-05-19T20:29:46.392583Z",
     "iopub.status.idle": "2026-05-19T20:29:46.506798Z",
     "shell.execute_reply": "2026-05-19T20:29:46.506011Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData' created with bounds [-3,3]\n",
      "[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data\n",
      "[#1] INFO:Fitting -- Creation of NLL object took 493.991 μs\n",
      "[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level\n",
      "[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only\n"
     ]
    }
   ],
   "source": [
    "r_sig = model.fitTo(modelData, Save=True, Range=\"signal\", PrintLevel=-1)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "9043bb7e",
   "metadata": {},
   "source": [
    "Plot/print results\n",
    "---------------------------------------"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "add6370e",
   "metadata": {},
   "source": [
    "Make plot frame in x and add data and fitted model"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "id": "5e1a5fef",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:46.508737Z",
     "iopub.status.busy": "2026-05-19T20:29:46.508509Z",
     "iopub.status.idle": "2026-05-19T20:29:46.710303Z",
     "shell.execute_reply": "2026-05-19T20:29:46.709490Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following\n",
      "\t- Clear the automatic fit range attribute: <pdf>.removeStringAttribute(\"fitrange\");\n",
      "\t- Explicitly specify the plotting range: Range(\"<rangeName>\").\n",
      "\t- Explicitly specify where to compute the normalisation: NormRange(\"<rangeName>\").\n",
      "\tThe default (full) range can be denoted with Range(\"\") / NormRange(\"\").\n",
      "[#0] ERROR:Plotting -- Range 'Full' not defined for variable 'x'. Ignoring ...\n",
      "[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'Full'\n",
      "[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'\n"
     ]
    },
    {
     "data": {
      "text/plain": [
       "<cppyy.gbl.RooPlot object at 0x55cdd0932b70>"
      ]
     },
     "execution_count": 10,
     "metadata": {},
     "output_type": "execute_result"
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following\n",
      "\t- Clear the automatic fit range attribute: <pdf>.removeStringAttribute(\"fitrange\");\n",
      "\t- Explicitly specify the plotting range: Range(\"<rangeName>\").\n",
      "\t- Explicitly specify where to compute the normalisation: NormRange(\"<rangeName>\").\n",
      "\tThe default (full) range can be denoted with Range(\"\") / NormRange(\"\").\n",
      "[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'\n",
      "[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'\n"
     ]
    }
   ],
   "source": [
    "frame = x.frame(Title=\"Fitting a sub range\")\n",
    "modelData.plotOn(frame)\n",
    "model.plotOn(frame, Range=\"Full\", LineColor=\"r\", LineStyle=\"--\")  # Add shape in full ranged dashed\n",
    "model.plotOn(frame)  # By default only fitted range is shown"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "a9012ef5",
   "metadata": {},
   "source": [
    "Print fit results"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "id": "b0fa3c36",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:46.711949Z",
     "iopub.status.busy": "2026-05-19T20:29:46.711823Z",
     "iopub.status.idle": "2026-05-19T20:29:46.821257Z",
     "shell.execute_reply": "2026-05-19T20:29:46.820383Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "result of fit on all data \n",
      "result of fit in in signal region (note increased error on signal fraction)\n"
     ]
    },
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "\n",
      "  RooFitResult: minimized FCN value: 25939.4, estimated distance to minimum: 3.77183e-06\n",
      "                covariance matrix quality: Full, accurate covariance matrix\n",
      "                Status : MINIMIZE=0 HESSE=0 \n",
      "\n",
      "    Floating Parameter    FinalValue +/-  Error   \n",
      "  --------------------  --------------------------\n",
      "                     f    5.0441e-01 +/-  6.32e-03\n",
      "                    mx   -2.1605e-02 +/-  1.77e-02\n",
      "\n",
      "\n",
      "  RooFitResult: minimized FCN value: 10339.5, estimated distance to minimum: 0.000279216\n",
      "                covariance matrix quality: Full, accurate covariance matrix\n",
      "                Status : MINIMIZE=0 HESSE=0 \n",
      "\n",
      "    Floating Parameter    FinalValue +/-  Error   \n",
      "  --------------------  --------------------------\n",
      "                     f    4.8979e-01 +/-  1.62e-02\n",
      "                    mx   -2.1518e-02 +/-  1.79e-02\n",
      "\n"
     ]
    }
   ],
   "source": [
    "print(\"result of fit on all data \")\n",
    "r_full.Print()\n",
    "print(\"result of fit in in signal region (note increased error on signal fraction)\")\n",
    "r_sig.Print()"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2cd533fb",
   "metadata": {},
   "source": [
    "Draw frame on canvas"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "id": "f1aeb073",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:46.823003Z",
     "iopub.status.busy": "2026-05-19T20:29:46.822877Z",
     "iopub.status.idle": "2026-05-19T20:29:47.034879Z",
     "shell.execute_reply": "2026-05-19T20:29:47.034063Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "Info in <TCanvas::Print>: png file rf203_ranges.png has been created\n"
     ]
    }
   ],
   "source": [
    "c = ROOT.TCanvas(\"rf203_ranges\", \"rf203_ranges\", 600, 600)\n",
    "ROOT.gPad.SetLeftMargin(0.15)\n",
    "frame.GetYaxis().SetTitleOffset(1.4)\n",
    "frame.Draw()\n",
    "\n",
    "c.SaveAs(\"rf203_ranges.png\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "2b849cb5",
   "metadata": {},
   "source": [
    "Draw all canvases "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 13,
   "id": "82eff388",
   "metadata": {
    "collapsed": false,
    "execution": {
     "iopub.execute_input": "2026-05-19T20:29:47.036622Z",
     "iopub.status.busy": "2026-05-19T20:29:47.036485Z",
     "iopub.status.idle": "2026-05-19T20:29:47.216430Z",
     "shell.execute_reply": "2026-05-19T20:29:47.215820Z"
    }
   },
   "outputs": [
    {
     "data": {
      "text/html": [
       "\n",
       "\n",
       "<div id=\"root_plot_1779222587207\" style=\"width: 600px; height: 600px; position: relative\">\n",
       "</div>\n",
       "\n",
       "</div>\n",
       "<script>\n",
       "   function process_root_plot_1779222587207() {\n",
       "      function execCode(Core) {\n",
       "         Core.settings.HandleKeys = false;\n",
       "         \n",
       "Core.unzipJSON(32285,'WkwISy8AHX4AeAHtneuPHDeS4P8VobAf7gAqj8F3ZH6SZes8d34IlmcsrW8glKRqqc6tLm13yZZn4f/98AtmVld3tWzPzcxh57CCO53BRzJIBuPFIOvfV8/3P7/bXKzfblbj6tuH64sf11ffbV48uVi/u3qz26/c6uyPF9t/e7/5w6er0bvV2Sfb/VV/+/rF/9683JO+otjX7/bb3cUM/M/txavVGN3q7PCl8d/vautjDcSYQkzJrc6+2F5sHu7Od5erUWbwyf7n8801+N321f5NBx9tz8/nwiALuBT23qpvzvZfri9fby9Wox8ku9XZN9vXb47T6ORuv9+9PSrnVmff7t7dTHh6tgWL4FZnz65fH/RXajy92q/3NKNKmRvQgw5R+dHl+u3mNuKk3er5odzNHh2KLsnLJ6l+PDR88pPd5avN5ZPtX+bhO0r8cvdq0yf2qazG+5IGb/8kN0250gPSSx6KLz6H2mrO0uhnWI0ShmL/qi/JR5B4FlZjVj+knFLM2mIVuvt0v3vw4urx9sPm/IfVGCUMSjsSgth0PN3vfjtzNYYwxKSq2mJMVZnIZze+nCUNlQIlqy8CKT27/vSv5K7G+35oVVr1LeQiNZXCMrjxdajrj9efOwbpNRX+dKNCrnT+T9dVbiasxvs9YRmb/e7pD3dPg43Ob2dDeD6lmkINWqK2xBxef/7ZD3dOz/x1cu+a6yUb1AZJRXNrWlKOgS4/Pd/99NWnDzsZPTsGnv7xnWUYYRy9f3dI/fzw9uDF1Y0PPXhxdeNbD15cXVd78OLquuYfP7xlYd8Xm4+fDbDXD2/XH1bjksx7DmloqcUWRYuA/LdvNvv1aoxUePxmO789uHq3ebn/Zr3f7nq3vnr/9sXmsr9/u335w4fr15/76xe713PiF7vX12l/6bmP168er7cXMAa3Ont4ubu6erPezh88gI93M5c9XrFQUIevF+uXu1fbs+3m1Wo8W59fbdzq7L9fbl99uAn+fA0+eHH1cLe7PCr/2avtfv0CVra/fM8HHm0/bF7d6Pfy6ceX27fb/fbHzdUJO/9ie4W0WCTJDK4vL1fj9392q927PS+/uNXZZx82L69W48X783O3Ovuqi57Ls+Dj88v1xeuNyYRvt3swOkn/6v3bx+vzzX6/yAQG8avNh/1p6qd/ePL4iwfPVuPqX5ZXtzr7dPf+xfnmk/dnZ8skfrPZr7cXjODc/6dX279s/ni15D+7CVruN5s1PIjGLfsI/m578Wr307e7d6xAt7qGnx3DM2O+LvD5BjE0E8VPCw95+GY1M4qH6/3+ZNAf7PddZNOzp59s9j9tNhezSLoB2XA+uty9/Xb3bjXKYIz41XoPxzfg2QIgsx90QH5xqx++3P24+frd+t/eH+jjh282jMjNxLPPt6/ffEEXZvFrtLrev3yzDOsPT97sfvrsx83F/sl+vX9/dSDJHx683+8gikPJLzcX7z9ZX3YYknnwEqI71Dj7ZrN+9fXF+c9LjbPvtvs3u/f7Y/pcaPbz9dVMcUvKcanvb2klfzfdB3nzUd3nu80L4wHbi9cfU4CgjIfn66ureYFQrmtcxwnvINKVH0PObv6bZPTOOz8FS+Utjt5y/ZQsL+Q85UMdP5W5JOl1Lst7G5u6IOJampT3llyQMgkNJjf/TSKjaHDSgpMaJgljLc7+mySO4oOb/yZJo4Tq5r9J8ig5uvlvkjJKbW7+m6SOwSc3/03SxhDUzX+T6BhScPPfFLwVFqW+n4LcBMMoGp2oOCltCnGUVpw0wDiF1Nst2Un2U8ijtEiOk5ynAFYzkslPoY6Ss5MUnEQ/hTZKzNbBVqago9Tsgi9gMkU/SgzOmi5xijJKrGDgpPgphlGid5Kb48Mx9nZa63UZqtDbyXGKeRSpzuqnPMUyik9gYCjHOopXZxVimGIbJSQb1+B1ino96jlPyY9SvbPuqkxJRqbYcI55SmEMQZyNDyBD1TqS0U8pjTbEzFpIU8qjjSloiE6pjHxGk9MypdqbYaBSnVIbrUUxNKeko6TqrBMSpuwNCSOxKcsoqTij2SkHG4cZiKO25T2NNS3vecyyvJcx1OW9jnF5bRD2/BkdvUsQcJ6KH70rrb+zdFT6exjplC2rqbCCJGrPYRFBHVY9A7TagUIb0hsptpLi3Iq1botUpkLzBoQ6Vdo3wKepggCAND9VMDAgl6kua5j5rWBgOb5NFQwAWp6qIZCzKzLVZSXHOtWleYlTVfQgCvmp+bEjyLuM2brBexibDQPvLOHOXfzUoEobBoC+fOdPsUzKArB0D42wcG0YqMOytdF2ftKZccUGIJ0jSQLozEsU9qXRcqQUgNSBFAA6F4Oo/KSdjal9uVqpKqQ3e0+Vd7X3EJ2fxPtRRFzKxhTEi1FoCi6FSXzoTEGZXpnEx7FVJ0Y1ZCfrh+TiImAeQ2guVpfTJL6Mkr2rLFCdxMMyipWVwqdYnsHF6EIAZEhcYVRlEuOtmTUUJhHYmLgoLgAFG1WYUdFJ4K288xfyJDBXZpOB57uSO8adngXu6sWpOvF+EqkjDI7+SJsE7grdFXUCR9dRmtKwTBL8KAmelZwU+Dt8rI9bAwyjVfPVNZ0kxDFEa8dPAmttGbHgWOYCb5XojJUFChdb6kFaxxnueiQCgpdRhSWpMgWmBJbY/6bg4U3RzX9TYE6E/tjfFHzu/Ba2PgUmpUVjAX4Kvo6szSqu8uU2hqgdLalT8DpSC2adp8CUQAVUnoLIWMX1ylOQYHx5rj0FJoX2qDwFSdTslacgmZq98hSYDzi6VGpPQWrHdmm0HTpnrarhuzSLuIMddaSnYFMCtr1dpqT3VVqcgk2JdTbENIWQDr21lpmUubcgHeDi16Ma6tLfXrkdOkzluKxheJtPU1yWMcSE3GIlz1K8w30xW74PU2Q9L/kGH+kEBvdVTfksU2Rhz8UNbAcNwkB4W1c3skzJXytGfkq+6wRUBwqHmkDXegpQQpwgUabkGY3OvehfssVsvN/6k1jMwEjjHKcE6fhw6F8S32HyfZgSyhIoIkENDrQ1dy9JNMgnw18MD4rSGzEhR2WDFlbvpyTG6nu/BFZPFdIRNFTwUwrImd6lME9J1xmnhDo0dwmRlli0c5esy9DH3KWeDycJ1uWej/Yxd4khCG3pUi+uS5+sdPSANGCFoxjYZ21KcRF6VjWa1OvzPaVogvdQ0STvoR7DYZM2ZUhyphG+klUOU93h66kHo6zXk9/zU+9P7+OUIcnejqGctYByn7Upa10g+5b1vXdwympdB/JhKp4pMKk4Fd+1dOoUT59t0qbi6bFN2lR819bRGKDCvsIMxQL76itshpt1eZ62qSBTyO/TNhVjYUaBvUmosFOgdakgWmyF9c9JtM/RRRBEshgRzoX7eCzdkmPzocxaO1yPmn11zJ1DqkCNfTiCEQKzNhVIcukf9eBYc/+gkgLL6itshmfBP09/MZ7V+2dDbVq6rbAZ7uKtE0gJunSPttHSe+8M6utzLhmvjagS+6CgkRS0c9aYvffh6O9GGrzWhRF2tKc6M8J52KcK1c0oZpkqRDeDVM6H/lrm9fTTvTqrN4bCVGcFh5ascO9dz2wLDzTu0OCB87hmmZrvnMAmZWozEzQm0fxhzsmZ0em8p83UOOM7tYUajdM0aHGmPWvjMNjomJ0j9yZmIqQwBRGlM1uhpHHC+ZPXjJCcmR317jbk6MKQfZgacrSTvRVuB8jaOPAilFwjQJoHmKVDH52pIUIXvs5nYZJQW2fGDR45g1SeNeEZpZka51GY2kyNS0t9rQKBEjbj3BDK9sy7OlIQ49w3lkXDZJyHiGXRYh+kTjdTu2aQYBxNXoBTLztz7VlQt7jYB73wNc+2hsxEWaqmYxulJePYfMe+m4xj98U3tWSW0kwnaWHY9DNdowN26Rqd/pkDOnQsHdDphRd0+FBesDFgHqEuoqeGuTjPi3UkH8jIvpvBaObOrVuN84BOLXeM5gGfWgYjy8UVkg2jaxCMZv7YzIac5UErXaD1aZtatyP7pE3NDMmZcbZuSHZsp1a6eO+kMLVyJGRp38zJ66/28Vl6bRbloWo3KZe8OrtCFoTqLP0WlGon7GVQaifszrFb7YS9dLR2ZegwChXLZrGwW+2EvQwEhmbPhKJaJ+wu/lrr07YUbX3WutBozeasS4nWTAea0cHiNNrqk4vNOfNMG6PWl/7SRl/6hzZs6c9N6Cx7rNrMqZcRmBn13P6sHJiYajOTnoHOFJfZnHn0YT5nHr185iCAjKBnHm0f0plFz8A1h2bc9IhFQ8EKk55lsQ+Tzmy6t6Izn557ogufnqdaZz7dB1NnnWGBusYwz7TOjHrOnL1rRlMghbnBwPe/SWGW0HH/mzTYDHYZqiGNOWJWS6SqkXnCNK6ThjJmbxZ4wuivY1IMdUl8tI2pOQzPVCYNOqbisneS2qTRjwn26STLpFHGlJwhkyaNYUzRGXuqk8Y4puAwvwuehDQm7/C4lTBpzGNUlxVf2KSxjLG50l2EGusYq8NErrTQxpgdBnClBcVtUoqTWiZNfozRleqk6qRJxhiw5rGNNYUxmsknLU+a4hjUYQu3OmlKY2jmJFDcGhm3Vc14CyfFk1ZcrU6UFioegqpOlBYwVl2DkmXSpGMIrnW9Ts2P5sz0q5PCHNW1jE9qUlhjc624IHHSjPfAteaC5ElzGqU4xXvVJs1mEOIICzIpfDE7DTPEkjf3ROCb6He4XDpkZrOWDhX0Kqd1hiAZp9j5uGrwGeCLCKKTlggkOC4MNGSMag3MPTfOuaWD3XGhpXaQTrZJS6Nb4uklIA5Zh2vCwIrr0+EqDVInrWIg7jADcfk4XEIdtBESCTOIl9jhvcAFrnDG5sR8DYB4F7q/0MDaQdolF2PeCV4ORhvOCAgaedLmDcRTYKB0EDTIxXdtTo8Oxg6CRpoUhxy5dB8Qh6yTQPcBSwdBA7B2EDTipDBICoMGoBoYQSNOqh0rfE8SJsWIangfZ7CPFf4mkUnhkninQQOwz2AEDcA+gxE08MEZQUkEDUCjKPOuGmgkJQn61kkx6ZMTvFgeXxpOO2AIvgHjjIJhLDBedScJd4Tlx1HwQBcXfKU8Hkw80gucO9wWGLeMk6R4xClfDc6gYzBSBbbjgs/kmztVMvgAix91FsUGioG2jUFuGHGmd1+JeIkGgovlJgNBJZKbR1w6YGJgASwggscRBx7+xAVsBoKG5eqIN49RwcUYvNU17yCg0BDsTBnTEECS7Q9lyEzVNXerKEOGrhtmmD6h7eKDhDuBtjnxZhhETd3t2xwKLriMwBz+Zch08jceSfMYYMq2jZNG+xGNZYZpHyOMrSXblQDGw9h9ig18Ir7YGQafiHHY/aL4J320HQHzk3ZPboXxCu7ZxmjEBltm4wYOLh72TvvA4AODp/3O78XD4mkfmPZh8rSPBKD9FBEW3Y/JeKSEKOkw+KSMoBH4eaX9VBBD0hhz2k8VIdVha78hwvpuFeORFAHXYcYje8Sfbb8V8MmCcDS3aKH9HBCkHab9WSe2PS7mI6cxUx+BTvvoxdVJw3Cg/VxGZKTB1l4dC/iQb+21EdoxmP5nHdniwu9rrvLi8SB2mPaKjPisyU+0V2zDsMO0V+LI2LBDl2ivpJGtRYNpr+SRtgymvVLM0W0w/St1NK89eoS110Y2/4DNHV90VHb/gGmvenaVOsz4IhY8CMPvaBDBgIFjCbRoOjNDZByRhDQiRnoCbSIeWDRNbZtQPAICLysJwVqFsdAsCdYs2wI0S4I1y74da9TbZpt4BEVga5EEZhZRAelaAs0iLGzHkRLQFuICYmTnEz+/R2BAjZQQmkVkQH7XCfjk2bWkBOOP2IDALAHEEBxpxoPdAo/ogMSsBIghPNDOLIEBWnxwlgCm5oSbx0PA1Pwh83gYpzchwpiyYwHqJkbmMTVejiAxyqMEmCJKYEuWAKYIk85n2GokAYVgnlsSbBPIOEFzApu1bSCWtiWw67LsBFmC7dCw77B8g60UZIptRDQHP2U3yKgZmL0V23qga83BUG1DiIFolY1SsQ0h20CvDgbLjpAo9FYdDJY9IVG6VR0Mlm0hYR8H2HaJcOOA0AJ3T5XBIGzbEHQJFQ8YywV8iu3YsDXEriXlYcDsDbElbDD4sTlkxNzYZBR2h9hkoT8wXDaIgtF2cw388KXA8dnOpX28ezOpw3BtiwgOr97BcNklsj1nYNrHuwfHV3GV8cKfYgtBHAyXXaIAhwe29tlhhKaDg+EK/j04PDDjg1PFlkl0MFxBwLC1pdH1jbQwBls1ycFwxQQMBJ8cDFcQMLaIUPqB8xgi9I9JAFxG9jANpn0EjC2p4mC4goBJtFfYmBWJPXaBOYXhSvJjsAVWHQxX2JKHgwPTPnvytt6agwELdoQtt+ZgwIIlkcw6cJn2sSVs9amDAQvWhCkoZmKKYE+Yn8G7zHhgURgcHAxZsCkMji7RPlaFwdnBkCVjo6JaFAdDloyX1BQwB0OWjPcWTURdAp/Mni2Ko9jWqmTwA44OBi1s2RucHAxaMviRXxwMWjL4AVc2s0Uy+AGri+BTwA9rxTv4tRTwAw4Odi0F/ICTi+BTwA84O5i1FPADri4yHqWPHzTEJqUU8PMuoBqDTwE/YHFwbingBxxdMHzAD5htc7ZawQ+4ONi24JIxGGsIGPzIRxMHZn7VQTMwcal9fol1EfCpec5nDQKXGUYjAcZepD5aFzD0B8zmLoFEn7w/I5bsXyw8b/XofLfex7Byq3OLHsvZrX5cjd9rSE4Dpl9xGuAtzWlQp9E7ZYs3BqcxOo3JacxOY3Ea4VHNaVSnyTtN4jQFpyk6TclpwngkRgTe1pwmdZq90yxOc3Cao9OcnLIHnovTDE9sTrM6Ld5pEaclOC3RaSHQJDstxWmBlzanRZ1W77SK0xqc1ui0JqcVU7Q4rdVphQ+r0+adNnHagtMWnbbktGGkFqfwUuPf6lS9U/agNTjV6JQAF81OFRO2OlX4vu2r44eBE3vzhHh2IMwz6tmLgQV71AYP8/UE7HjYrofXehisR6/2sFYPP/VmFcNJPezToyV7GKdHT/CwSI8m7GGOHo7oUQ88vNDDAD1cz6PMEpKNJs2DGrA3b4FLMDKPRuBhYR6+5VE/PRzLQ3IebdNDbB6u5CEz34N3qAFVeRwzHp7jMzXMK4dS6ZHxHveLR6f0SHiPNukR7R5vi0emezwtHovEI9E9zhWPKPeokh4Z7jEQPMLb40Lx6JAe0e3R1D0y26ODewIIPBqkR6X2iGqPAumR0R7V0aN7ehRGjxboTbSb0oBY9ihLHoHs0Rs9otgrNVAbPZLYozB65hyvBUoDD4IcmPMe1cSc44hA3vPA4mfOTVHEyyCmIfZ4CObc9ENhzntgAHNu8VK4BpCqPCyQgjaYc+x+xCMP2mDOTSHsEVnMuamDwpxjrCO7eFCDOccIRx7hyKAN5tyUPwxqBAsPajDnFrBleh7mMWKBBzWYc3MGCnNu8XEW8SLMuQW4YcbCqHlQgzm3CDuzQntkGHNu6ptFyJneZkahKWw9hIQ5N3XNrDPT08zuktr+/Msvv7h/VDQnxws+Gs3ZD178ykmZOWKToO3LCwou0c095d7hA7dOfAAeneq4faDj+qyLcMDk9LDLl+vLHzaXR4dnesLRJ+eEw3mQbzcf9g8uXhN6TQAqYM/0g2cMLP98+/pihVekw0ffJ/vRjtj2YhHG6w/b00jxB/v9A9IJ3H61/XF7td1dXK3GbNH55Bx98Iv1i81yoIf2DO4tJFow+Ouzs6uNHbSBz86JB7Sj4b19+cMXm4vXnBPygyfI2eZgqWp9IeT8drX9+RJofiiyNE9g8rN/mh7adP5f9PBf/2l6eJigv3IOP1lfHh2T+mR9uRCFhWBz8oxlev74SV8Tn16uf+qHMDr89bv99YGPDsxnPjowH/v4+t3+0x5v38/KEQTPIrJl9PW7/cwS6MTX7/aP7GjZXPTRdj4XcBJBTwFLfLXdc0Bugb/d7c4tgJ6EfmDl4e5iv3t/eTWfWniwn9G5xTEf7PcsYmNSv8ILwl/JDFgrdL+fSqKTQJyIGAz67OLVZ5eXu/lkGivbQCtOU4/eX7yc2QKZgEdcDHCeQnI5ADMXpv+Ac2HWPeDRfH+xeb25eHV8ygbseuoRh+VD14lL28uZQj6xsIVDQcZvJkW3OvucsxCbq1s8fE598m79krMA1vbhnN9RHw6H/OY0cDyUu4nNoeiSvBS91bSVu93vo8TrI0afb68gyGN8SOJ7MzrF0+2l3NJwH52l6JI6F7yFDaW+3F5s375/+6+by931UQ8ybhy9NBbfT708vtycbS7/+xfXpXv60cD1hONugulx6nU/e+qnm7PPV2P2TNoh5bvVWG+mPJ2Pkx2KPJsTHq+Pie/x+gZt0fgh6bplSzo9cvp4/epG1xm7x+tXt86w9sQ7TrE+Xr+C2J9eD8+c8uxGCrJxPqtEg9uXP8wnlR6v3/Wzok9npnFIeLZie3J19uTl5WZz8Wj90rgPmMDWjoYfkHVwRLYkHc/HUuto/VAE8Hr1kLIQUCerXubyLVx0VYZkB2lINJ6hnEsFMhaDcdrB75bzTuR9DmBH/57sL7fvPt283L5dn18dDhUZS561nXBQGY56ZwVudc/SjvsHHVniUQcP8EFUdSUFBaSfbT7UOrBAq9O7Y6ddDaY/y/fpTmc+uz9cXGwuv6F7fJelZp+9Wo3fc5Do3j0eEu4JRv+9dK8tcLmX7kkBCp4CaS5yDVnScXrwN4sf8pqf85q/l/y9dC9ZSWv7H/T4M8xls361uURa25kpG6YD9Gi7f7QQTZ6Jxk5iMYuHHKMlG7WX63OrzOz/j932gsRFEXi4fncMfrt9e1Ana2sqLRnH+MPb9esNHzow+Ifri1fnm+/ebK9+2Fx+w5HHzvh7+ie7D3Nan72eapgcHdr803Z3vr1YUueDi73ow+3ly/Pb3H7O4rQpSB8JwKeo3J99ePf0WO1ZEp8dJz67q+SSeKMkBb9cf/h0+9puC4AIv77cv9k9XL/dXK5n7vMPNNeMpyxXGNxWcUzGfcxcg28wQEcMC3DhPXPu0fDdEoynRlg/Tc/K5Pw8/7fz8pbACfmTc8jHzIMGOzwLi7/7oD3641cPGYx+Y8OvDtzn8unHhi3bUdjZzj3jpoHnH57n/PLVK68xvKis/L4UV+Pq0XbPycd763tX71/cszO/ZB8Ne9cn/saB/+ut369ebs7h/3hyVmcfsV5n0/Vjt2LMY/DBzKajTn8A+Ke2d8H/xRZTHZ/X6uzp0Vn7p9fn65/2Mhz6Pnu0vbyabZcv1ssbN4iERW693Xy6vXp3vj46oQ4jPbBiIwycANcnxr/cvfpi/WKGf8UI/30T9fPtibIjylf3/td/u/df7vkh3Puv/7EmbuBgsQm1g+4wezh+xVVxNHWHietTaNcikPb/fNY+4lj4fbP2l9uz9h9rkv4/maJbvpAjR4k3DvDZxf5yi1ZpHb56//aneVnzuixx3pc7Mux9zvhy/QFr7y4JOBuCXOTR/SNf7S7fLlYGdDs7M/rVEmdPemvGcW7cQoRfgEuJTn2Qv/+2in5TxLUlstwcYa19sr3AkPjs8vJr7rcANeCvf9xcnp3vfsLZwu0Fl5fwt5v7bYWNONtvM2nDrQy263ZLH/j7O7XfsVR+W9p/s9thelP4LlFne0mH2zvePH+Li+zT9X5NBWNPq3HFB3avL9dv7+3O7h1KPH93vts/f24C8UjmQ0N/g8Rn7H9d4qPa3fZ3f/Vut73gQqtZpP3aLPm+K3pfB3X3dag8Mo/IQ9z9RkYjo5HRyGhkVDIqGZWMSkYlo5BRyChkFDIKGZmMTEYmI5ORyUhkJDISGYmMREYkI5IRyYhkRDICGYGMQEYgI5AhZAgZQoaQIWT4oXFv0vU/kioPu3Dp+h9J4RrkH0n+ZlJ2foBhXP8T54d4DXL7k/NDPkmpJyl2VdThnzjwBXF6QFfoE52jl3SXfjMAjARDwtgwSIwWw8b4MZCMKEPLGDPYjDrDzzwwIcwMU8RcMWnMHtPIfDKxzDBTzZwz+VAB5ABdQCBQCiSjg9odOM9+B6ElIdw9N0Lne5h6VALnc7BntBcLbefIeiX+PWd7EZeTK56o96xEyRNVnwmlL0SgOELgA2UsQp9a1RPAbpuk7H1yzomIz9ZjECpxFOJiTi5FYjSIC0ku1epSFhcbR9TVBWK9a7I4C9vbY5cwqquNOPNaaK4IiKUEPkWsX9nlyjuHBJphZbH8KVtuMeT7CDTe6Y6F/ufIeXuy7JQAZwUsJc+FbaD/HiLgwMYOlxadCKiDMFsk1g059BnXWf32jHscfOyseY9a+9nTN9vXb/76atyV9du1vi9DrKVpLb75JJCwTzXG0kJuuUiDtjXW7NUTbNxYBCmnFn1lq7s2VkSQRiab2KGwLErWyFVutWa2209KlKH4qlFVa/W1BVZSqy3m6pOPtbC8TkqEXH3JkqpUomjvqHILsTqk4uEVmoNCxneg3iJ7yjEQOVDB46SVqDFXzTkUTZzpOK0iUjUGTSGH1hIcIWSuvyvSSlXGtIYUxXqai4dn3ELsZMTq7Y+WgTkKJfoYtBCyNpQWaTil7Eumc+qjVg2xtJokwIl8TD6UIKUmDqAMtbVaNUnOgcvyyhC0wXdrLaWXuDUNd/TlpLcnCbdGrA0xlVxKrKKJIJU2SEmthhxStlgGHXwsrcVcckwWc8jubAwSNRclyIiEFjRkLhsEUwlDlBxySyV4Qom4QTC1XEspUYkvkjxkicmXFtVb6EMZUguhtRJraBaJOKRcJdbYSkgW3Dw0X7yv6hOX2bnghxZSSa2l6rMFEw5NcvNeSyuELwQZaomtaM3CkiGhxcD6qTkRDxlkoFHNpYZU7DDPkAltCaVKrkRi1KFlCZJzycmCYeoQQlQlaKISO1KGnEJLrWQiTTgLNWjOKWuorVhARhhqTj7FFEoufDMMoqElLT73mB4Zas0lxpJUKzGaQ5NQW5MSaiS8Q4dYgtYcYsueCJE2FG255sqqDq2cUuUJRZ2Q+sk6PqHbk2+ccIt6m6LqkEMt0nJr4H76iZNlfLI47iDsk1V7ixWcYH4HorfW9R2r9haPuuMbJ/z21sq/gzfcYgUnJe4YD5+zpFCKNM0cajtFTGKVEglyagGBfZu/3JGQo7SYvM9BYkFnutXKHTz7RHbc7H6/K/DZ75OB39chwSCy1EQUnWlj0kJIxZdUY0ItK1lzSaJNCmZlHXIM0rywuHImIXBLaUg5ZK2csYN955JaSjVXj7p3q0RlAjJRfRI0++IatBNTriWobyg/pyViyCq5FsmBALM7qtxCrEH1rWpKUghkvgNz9Sk0Db7kQMhdO20kFcKsSogldHF/UkVakyCVESDYrLH4WlINhVootLWpD5HAvBgK7PwmXifj1U6/6VMtMcVUSiG+ug41w3Sz90kI+qwQfQyqwWsLRFoOUrJHAhe4KXi1HDSHGFuKzaRsrDUUX3OrWvIdk3BHV046e5Jwa8B0SL6kUnOLUQohpEMIhRhN4v2qRWQOXkMIMeaYemza4E1f8swMoZQySAw51li5MZUjsUPMufGvEJwZB201BEmSVRf+n0tJGrVH59UhSyjZI0ECwYs6pJZC9kVbtbBENTEkTUq1UEvkUkIfkQCixBQPLXHsB92Gs5BhqK20EipRgHbQc2g5Fx99iloCJz8HZHZpIUoLHPf0Q44FDcpLrnZaZ2hVQ9LQfFFiGCHr6FOuKBqCqMs1etXccq2EXuZBGySboQSEchzQpmoV1FBiMOPAJNs9wyEQDBqGJk1C8miZRO0OrWiOvjWfU7JQziFJiRoY9kSs8FBz8YWoTc1qqtFtKr9FTSdUfrKCT2j2Doq8xSfabWpqQ9bkPUtNUhJW7Aket/jEycq4g6pPVuzJN06W2wmmt9f0yYq9xZ1OmOIJp223v3HCF+4Y0pN5OumKxJpz9k2klWa2+Ak/F674zc2nStD4yQjekVB8C7WFmmMOhY/eauUObn0ygjcFDYLs+c84tm0jYzXeucXw/HBb8fOfzSufVG+sQrd6frG7fPun9TlXsXoD324v1uefbC/m/Ui7yvX5xZPt67dr8yc+3xw8tDhN3Or55fqng9f2vrjV893l9jVf+c7unzVP6fGxgGM35ewA+0+/BN6Yv8Uv8XffxD3/3T7dh+8vf/zo5ncOOR5v5JrL9jnO9+8//Pm5RSV8/+j9+fmfLe0IPnb6Pr7c8bsA293FwetL9pGnF0l3w9e7gDMZE2zxq9vqv+3rxZn8cV+vtf9rvl4KcAKGi8Tv69B4FB6JR3D3zdHbcO8WHolHcPcbfl17FN4Sj+Dum3+34dotPBKP4O4XfLqNh70l3oK7b27dhlu38Eg8gruf8Oc2HoWHgcHdN29uw6VbDh7edHDzhtnXG3Dtmru38WY+38KbOX4Tb+b9pZx5e2+6cC2JRz5Nv8MnjHvn6B8V+139h38k3eE5Djc/Rik51LF/lnRa6sTljBv6ZrGEUzdzyEJClJBr29yX4ghsPziU+UfSqbtabhTCXR1u1LKUkzKnTu54UiadfOdGEftyOSlz6hrHnXP9TxyO8TA7x9PsIC+zk7x1R/nsKuf/uMvT7DIvLgzN4SwPDkKD4iC95CBEKBLSTENzOMmDywNjWxzEXLBrsSkxmnGZObxgiF50fRwJzeGBQBPDgEBZbw5HOdoc2n5Bq3N2Luo3XeXzOs0yaEna8AvxSxPmRS1u9aOsxu9vZeJAvlH6VgJeuhsl8Ewf18DjdiOhtZsJibNuxyWU34M4Skga5bhEwZqpRyXUXD7XCejZzdz1M+Y1JvwfhxJZtGgz57iVUO9Tk9S4KnLQkj3Kfk0F0d0TcqjBh1bnVrKohNRSXVDPGVswL30tBZsWTX/GU2NJyQefWscCa9arx0vV+6ohxphiLJyMzWHwQWNs1WcMJhLMakqKD85nEkITj7GuyQ6+hgHXnI8hJY+ZnOMQWsWjhydV2OPAFMq1lYQJnlzGT+NTCkULLjdXPOutefN01RRcwe8ZY/CiIWZ0RuyFnDzWdSiY63lIrWhRxZcQsmsylJxES+GwVOSO28GXTJ3qc+H4Yh1ayaXmjAbfj5J17RVXQOYIWfD8VkbFo1LYNZEYB1wrueIM554DHUIN0orXWDjUVXBZ8tMjudAfbkdAuddqxg+33fg8mHum4UrlnvaAEu2rFJ9C5X6T7AfcnTnzIx0tuKBtyDXn1BT/a3Mx5UFTzCnWglfAxdIGyfg0sWtLdbG1oZRWik98V13yZdAkNrGNpZZCwLnKwcFQMjtUQQffWqsxe1wFbBmhyPpQJcdGlagDnletOE69uJTigP+/JcHvyP5WGRgOTakGpUAbalFjerX46lIyKzdIEY9f3BLUS1JJEjg9SBX2GHz1oXF6M6XKxJUWNOaYqZIGaVJVasZKdSkhp0rFO2oHNVPMQ27qPb4XTtjTN00xavBeVbJLIQ4p46nAu1zVJSmMIO7eohx+japDbIkFlThkHBsObXxEPmiuzcXc8M/G3JIXfK0x5kHpeUu1sqhDg55w4XqPl9aFxMxGzeZBsruNBmlFcI+EiDWuBeJoNTIJdllJHVjzMWcNnH+V3IZWpJbM8HCeMoUh+RqKMROOGUJPIQXvi8eXzvk/DPSUY0mx2qnUiESNPsYWOBCqkZ3j6FWapBqDa3mokrz57pWDkxUzt4BjC5hursZBaquhMYg5iittyBFOURo7Uc2VNNSKh67UGHEcZB0inkkvbOFIYd1zUtMHdgSij3AGDL8QqtZahFPxA+dFWwySfYTphYG1mgq4BY7J42qIAYdrbLWS4Bu//IPPkAPBJAg062tRjrbCSJs58aJGjqaTUJOkjEckVvYoBy0ts/Arw9wTSpTUcuJXkXpCTgXnYZjZeZbUvJaYuUy8M+tYW5BSDhW897mwrq9lREu1wQh7lcSv/UiTzOjwDYRj1pQ5Tn1ICDVE3GWHBGkR/nuUULxwevq6RE6M8nFCqDelGx7Ta7xMJBuDX9BARh8E0SK0WdM3S9wE820xnVX/7FY/EkB7/OEu/MNqrOyg/kfY/L1p+d9tHJ3a/klvcrObtv8/s/14tt0/vzg/7wFD12FDxyblx4r8lVammYkfjyj6h1uZPYj413biFyvTbDg1q62ZKdfMbqtm0VngTrHMbJnZMpNldkvPnsEyQ7f4LNOCf7wZd2pmn5oZ2Mzw63agZRZ7L5aZj6zCbhZaZrTMYJn2KpbpLdMHzC0t9uTREs/a7GmZxTLtkS0zWWayzGjpkdLBMsUyxTLZlT36Z3FLxB+Ezu+zlIMZd8tS9HzB8zH7rh8if5h3/HYEGHiQ8QNJ/LbEUEmvpDfSlT9MJuJh6Cu/eDZbVZy5Z1Rk4AT9wEF5C0HipgOGUhhVGbiewKwtHN6MP0LDBWalhykJxhYeduaPn7Sb7TDuq2OmA5MemP8eysRPbEAeAUph4xs77TftpIXSTjTkTJBH0Kal+lqxKHDeVlX20is+96xD06wxF/YvExE8Qwzqawhs4gR1JQzaJCUUq1CbIrFbZYPWs72g6gqxJS1rqihWJLARrkViwoGbqqthiCH44DM7E1pO9fCqg0/sOGgkaqG6FgdPHWXvLrPrU4dUS2mN+zC47iMMiSsWshJ0INz5caKox4GND6wUtdswxA8JQ6c0TXbVQh2KaI41sIVP8EEe0FNTTjma5n5DkyeSIDE1vtZSW+4XpIYhV8yFnM2KEwbM+1pTiilyJ0WNA/I7hyYJSXii7IvWIfumpXnRfpslP14ZEwklc6dniAMhAdl0SC58jKWbGIQjNBJOzIGSBvW1oCgl+9GMihXSTIkMkJfGoabGjh+6e3LRtyEIu7zqS4T2Qhhy0VxCloZNFmMZqg8phdQ0KtotUYOlNm2RO6diiUODcKSFGjI/mVKGrNgfhEZYNNctkwMNWmMkRiLXHIpL4gc2eAq6ebEEHYqkVNm35qKIFOrQ0OdqrRIJjzuxQZLndywD6lzl8peU8kCUTfQlZYw2DAhfmoiGUEoyC6Oo+pBijdDsiZCmBkoqcVWRqCPMGB8qFgbaJ3aNDKox5KxauD4rxTIUCRV9ny2vOwyMMGC8YxknLslKEgcfWqsSCfAJLvkwaK7CLq2y0qPKQLiN0JXYMOJumBxcmgQ5leJrDVJ8JLJviESNtFprZh4ZniTs/ebG9SgxME1BMb3ZgHT8cCfltRBFUi36b2gpqLQYM79m0+KAosd3pcLiSiN2hZ8BFcVpEnIawJs90hi4tPamXWM3MjEHavZps+tEI1uAPviQa4AZejE2gGNPPdy3YSpGSbVi6bDxeNv0KRl6qkmzauPumJyHFFQrpguapiRchaXGUlJtmOOxDDEQNVFbheLuMI7YaWylVAK1uMdGYD+hpKYtRbtuJ2NKK3FxFa+SKpOAeSQVL+cd5lMbWvDN4yUo3A7T0hCS16gwNGL0GkKOQCCF9lVdJZxDFJue8C+9bWAFXB6hVi+xEQeSYc4sbp8zu5HEyuUhlZgSY8M9RiUOSQXjuiX8PkiAWgNY1JK4sKmwHn2MPrK1KwEpokSXVHwG3CeV6wBpqSf6gKt6chk0wOpZF2ZBHBlxtd2y6kgwIfefGr0dMvjtH2F+8PQPT548+PKzj28M/cvl5mw1Sv67GxMvzvkl1991yODbx+sf7QA1xe86ZsC9Nkc7R3/LGQE2jewsJq6lmoWwQ8Lk7LeZZTXmRKyUmiYkxc57c1ATt1sitomAghDtd5BJbnUg1iiWUhoX+/Fx+wVbllfA+GblZoKhOPk5Z2nEj9ZS1KCEO3IYdM5CxWjSqoqi51BryVLtWwIJwWK/I3x8RpSO/eGCO0Ko9OTN+tXup6Pjqw93lxwGX7/a8pubFD7EJh9maj60uOcMmVlZH736g9MdZPY7gAI7d7cPwB/uALq+oWc+i/3Fjp93BVEG5frXtPutPduLv/JnZW9dXfLFer+xQyV3ERKIz720DlpPP3oU9ff2fxmO4zsBDt0nc75XYb5Q4XdTct8Qfbvd90NP/fQRBuvXtlt/OPmH5Ih48WG9of1y+KFdflyXC7E4hPnk5eX2HTv71nN+p/W8/07rxcXm5f76Kgj72d/+aful1F/+D47PR28=').then(json => {\n",
       "   const obj = Core.parse(json);\n",
       "   Core.draw('root_plot_1779222587207', 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_1779222587207();\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
}
